Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 3 : //* 4 : //* All rights reserved, see COPYRIGHT for full restrictions 5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT 6 : //* 7 : //* Licensed under LGPL 2.1, please see LICENSE for details 8 : //* https://www.gnu.org/licenses/lgpl-2.1.html 9 : 10 : #pragma once 11 : 12 : #include "MortarScalarBase.h" 13 : #include "DerivativeMaterialInterface.h" 14 : #include "MooseVariableScalar.h" 15 : #include "Assembly.h" 16 : 17 : /** 18 : * This class enforces a periodic boundary condition between a microscale and macroscale field. 19 : * Target example for the Diffusion equation with isotropic, unitary diffusivity. Coupling is 20 : * made between a scalar macro-gradient variable and the temperature/concentration field within 21 : * the periodic domain. Primary and secondary surfaces are on opposing sides of the domain. Only 22 : * the macro to micro coupling terms are handled here. The micro-micro coupling terms are 23 : * handled using the PenaltyEqualValueConstraint applied to the same primary/secondary pair. 24 : * 25 : * The applied macroscale conjugate gradient is applied as `kappa_aux` vector as an auxillary 26 : * scalar. The computed macroscale gradient `kappa` is equal to this value for isotropic-unitary 27 : * diffusivity. The volume integral of the gradient of the primary field will be equal to these 28 : * imposed values. 29 : */ 30 : class PenaltyPeriodicSegmentalConstraint : public DerivativeMaterialInterface<MortarScalarBase> 31 : { 32 : public: 33 : static InputParameters validParams(); 34 : PenaltyPeriodicSegmentalConstraint(const InputParameters & parameters); 35 : 36 : protected: 37 : virtual void precalculateResidual() override; 38 : virtual void precalculateJacobian() override; 39 : virtual void initScalarQpResidual() override; 40 : 41 : /** 42 : * Method for computing the residual at quadrature points 43 : */ 44 : virtual Real computeQpResidual(Moose::MortarType mortar_type) override; 45 890880 : Real computeQpJacobian(Moose::ConstraintJacobianType /*jacobian_type*/, 46 : unsigned int /*jvar*/) override 47 : { 48 890880 : return 0; 49 : }; 50 : 51 : /** 52 : * Method for computing the scalar part of residual at quadrature points 53 : */ 54 : virtual Real computeScalarQpResidual() override; 55 : 56 : /** 57 : * Method for computing the scalar variable part of Jacobian at 58 : * quadrature points 59 : */ 60 : virtual Real computeScalarQpJacobian() override; 61 : 62 : // using MortarScalarBase::computeOffDiagJacobianScalar; 63 : 64 : /** 65 : * Method for computing d-_var-residual / d-_kappa at quadrature points. 66 : */ 67 : virtual Real computeQpOffDiagJacobianScalar(Moose::MortarType mortar_type, 68 : const unsigned int jvar) override; 69 : 70 : /** 71 : * Method for computing d-_kappa-residual / d-_var at quadrature points. 72 : */ 73 : virtual Real computeScalarQpOffDiagJacobian(Moose::MortarType mortar_type, 74 : const unsigned int jvar) override; 75 : 76 : // Compute concentration jump 77 : void precalculateMaterial(); 78 : // Compute penalty parameter 79 : void precalculateStability(); 80 : 81 : protected: 82 : /// the temperature jump in global and interface coordinates; 83 : /// TM-analogy: _displacement_jump_global, _interface_displacement_jump 84 : ///@{ 85 : Real _temp_jump_global; 86 : ///@} 87 : 88 : /// The stability parameter for the method 89 : ///@{ 90 : Real _tau_s; 91 : ///@} 92 : 93 : /// (Pointer to) the controlled scalar variable 94 : const MooseVariableScalar * const _kappa_aux_ptr; 95 : 96 : /// Order of the homogenization variable, used in several places 97 : const unsigned int _ka_order; 98 : 99 : /// The controlled scalar variable 100 : const VariableValue & _kappa_aux; 101 : 102 : /// Input property from user as the value of the penalty parameter 103 : const Real _pen_scale; 104 : };