https://mooseframework.inl.gov
ADPenaltyPeriodicSegmentalConstraint.C
Go to the documentation of this file.
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 
11 
13 
16 {
18  params.addClassDescription(
19  "ADPenaltyPeriodicSegmentalConstraint enforces macro-micro periodic conditions between "
20  "secondary and primary sides of a mortar interface using a penalty approach "
21  "(no Lagrange multipliers needed). Must be used alongside PenaltyEqualValueConstraint.");
22  params.renameCoupledVar("scalar_variable", "epsilon", "Primary coupled scalar variable");
23  params.addRequiredCoupledVar("sigma", "Controlled scalar averaging variable");
24  params.addParam<Real>(
25  "penalty_value",
26  1.0,
27  "Penalty value used to impose a generalized force capturing the mortar constraint equation");
28 
29  return params;
30 }
31 
33  const InputParameters & parameters)
35  _temp_jump_global(),
36  _tau_s(),
37  _kappa_aux_ptr(getScalarVar("sigma", 0)),
38  _ka_order(_kappa_aux_ptr->order()),
39  _kappa_aux(coupledScalarValue("sigma")),
40  _pen_scale(getParam<Real>("penalty_value"))
41 {
43  paramError("sigma",
44  "Must assign auxiliary scalar variable to sigma, rather than nonlinear variable");
45 }
46 
47 // Compute the stability parameters to use for all quadrature points
48 void
50 {
52 }
53 
54 // Compute the temperature jump for current quadrature point
55 void
57 {
59 }
60 
61 ADReal
63 {
64  // Compute penalty parameter times x-jump times average heat flux
66  ADRealVectorValue kappa_vec(_kappa[0], 0, 0);
67  Moose::derivInsert(kappa_vec(0).derivatives(), _kappa_var_ptr->dofIndices()[0], 1);
68  if (_k_order == 2)
69  {
70  kappa_vec(1) = _kappa[1];
71  Moose::derivInsert(kappa_vec(1).derivatives(), _kappa_var_ptr->dofIndices()[1], 1);
72  }
73  else if (_k_order == 3)
74  {
75  kappa_vec(1) = _kappa[1];
76  kappa_vec(2) = _kappa[2];
77  Moose::derivInsert(kappa_vec(1).derivatives(), _kappa_var_ptr->dofIndices()[1], 1);
78  Moose::derivInsert(kappa_vec(2).derivatives(), _kappa_var_ptr->dofIndices()[2], 1);
79  }
80  ADReal r = _tau_s * (kappa_vec * dx);
81 
82  switch (mortar_type)
83  {
85  r *= _test_secondary[_i][_qp];
86  break;
88  r *= -_test_primary[_i][_qp];
89  break;
91  return 0;
92  default:
93  return 0;
94  }
95  return r;
96 }
97 
98 ADReal
100 {
101  // Stability/penalty term for residual of scalar variable
104 
105  r *= -dx(_h);
106 
107  ADRealVectorValue kappa_vec(_kappa[0], 0, 0);
108  RealVectorValue kappa_aux_vec(_kappa_aux[0], 0, 0);
109  Moose::derivInsert(kappa_vec(0).derivatives(), _kappa_var_ptr->dofIndices()[0], 1);
110  if (_k_order == 2)
111  {
112  kappa_vec(1) = _kappa[1];
113  Moose::derivInsert(kappa_vec(1).derivatives(), _kappa_var_ptr->dofIndices()[1], 1);
114  kappa_aux_vec(1) = _kappa_aux[1];
115  }
116  else if (_k_order == 3)
117  {
118  kappa_vec(1) = _kappa[1];
119  kappa_vec(2) = _kappa[2];
120  Moose::derivInsert(kappa_vec(1).derivatives(), _kappa_var_ptr->dofIndices()[1], 1);
121  Moose::derivInsert(kappa_vec(2).derivatives(), _kappa_var_ptr->dofIndices()[2], 1);
122  kappa_aux_vec(1) = _kappa_aux[1];
123  kappa_aux_vec(2) = _kappa_aux[2];
124  }
125 
126  r += dx(_h) * _tau_s * (kappa_vec * dx);
127  r -= dx(_h) * (kappa_aux_vec * _normals[_qp]);
128 
129  return r;
130 }
131 
132 void
134 {
135  // Example showing how the penalty could be loaded from some function
136  _tau_s = _pen_scale;
137 }
138 
139 // Compute temperature jump and flux average/jump
140 void
142 {
144 }
static InputParameters validParams()
const VariableTestValue & _test_secondary
The shape functions corresponding to the secondary interior primal variable.
registerMooseObject("MooseApp", ADPenaltyPeriodicSegmentalConstraint)
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:435
This class enforces a periodic boundary condition between a microscale and macroscale field...
This Constraint adds standardized methods for assembling to a primary scalar variable associated with...
MortarType
Definition: MooseTypes.h:771
const MooseVariableScalar *const _kappa_var_ptr
(Pointer to) Scalar variable this kernel operates on
const Real _pen_scale
Input property from user as the value of the penalty parameter.
Real _tau_s
The stability parameter for the method.
const MooseVariableScalar *const _kappa_aux_ptr
(Pointer to) the controlled scalar variable
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
void renameCoupledVar(const std::string &old_name, const std::string &new_name, const std::string &new_docstring)
Rename a coupled variable and provide a new documentation string.
const MooseArray< Point > & _phys_points_primary
The locations of the quadrature points on the interior primary elements.
unsigned int _i
Definition: Constraint.h:35
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
const unsigned int _k_order
Order of the scalar variable, used in several places.
const ADVariableValue & _kappa
Reference to the current solution at the current quadrature point.
virtual ADReal computeScalarQpResidual() override
Method for computing the scalar part of residual at quadrature points.
std::vector< Point > _normals
the normals
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
virtual ADReal computeQpResidual(Moose::MortarType mortar_type) override
Method for computing the residual at quadrature points.
const MooseArray< Point > & _phys_points_secondary
The locations of the quadrature points on the interior secondary elements.
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Interface class ("Veneer") to provide generator methods for derivative material property names...
void precalculateMaterial()
Compute concentration jump before quadrature loop.
ADPenaltyPeriodicSegmentalConstraint(const InputParameters &parameters)
const VariableValue & _kappa_aux
The controlled scalar variable.
virtual void initScalarQpResidual() override
Put necessary evaluations depending on qp but independent of test functions here. ...
Moose::VarKindType kind() const
Kind of the variable (Nonlinear, Auxiliary, ...)
ADReal _temp_jump_global
the temperature jump in global and interface coordinates; TM-analogy: _displacement_jump_global, _interface_displacement_jump
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
Definition: ADReal.h:21
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
const ADVariableValue & _u_primary
The primal solution on the primary side.
const VariableTestValue & _test_primary
The shape functions corresponding to the primary interior primal variable.
unsigned int _qp
Definition: Constraint.h:36
const ADVariableValue & _u_secondary
The primal solution on the secondary side.
unsigned int _h
Used internally to iterate over each scalar component.