https://mooseframework.inl.gov
ADPeriodicSegmentalConstraint.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  "ADPeriodicSegmentalConstraint enforces macro-micro periodic conditions between "
20  "secondary and primary sides of a mortar interface using Lagrange multipliers."
21  "Must be used alongside EqualValueConstraint.");
22  params.renameCoupledVar("scalar_variable", "epsilon", "Primary coupled scalar variable");
23  params.addRequiredCoupledVar("sigma", "Controlled scalar averaging variable");
24 
25  return params;
26 }
27 
30  _kappa_aux_ptr(getScalarVar("sigma", 0)),
31  _ka_order(_kappa_aux_ptr->order()),
32  _kappa_aux(coupledScalarValue("sigma"))
33 {
35  paramError("sigma",
36  "Must assign auxiliary scalar variable to sigma, rather than nonlinear variable");
37 }
38 
39 ADReal
41 {
43  ADRealVectorValue kappa_vec(_kappa[0], 0, 0);
44  Moose::derivInsert(kappa_vec(0).derivatives(), _kappa_var_ptr->dofIndices()[0], 1);
45  if (_k_order == 2)
46  {
47  kappa_vec(1) = _kappa[1];
48  Moose::derivInsert(kappa_vec(1).derivatives(), _kappa_var_ptr->dofIndices()[1], 1);
49  }
50  else if (_k_order == 3)
51  {
52  kappa_vec(1) = _kappa[1];
53  kappa_vec(2) = _kappa[2];
54  Moose::derivInsert(kappa_vec(1).derivatives(), _kappa_var_ptr->dofIndices()[1], 1);
55  Moose::derivInsert(kappa_vec(2).derivatives(), _kappa_var_ptr->dofIndices()[2], 1);
56  }
57  ADReal r = -(kappa_vec * dx);
58 
59  switch (mortar_type)
60  {
62  r *= _test[_i][_qp];
63  break;
64  default:
65  return 0;
66  }
67  return r;
68 }
69 
70 ADReal
72 {
73  // Stability/penalty term for residual of scalar variable
75  ADReal r = -dx(_h) * _lambda[_qp];
76 
77  RealVectorValue kappa_aux_vec(_kappa_aux[0], 0, 0);
78  if (_k_order == 2)
79  {
80  kappa_aux_vec(1) = _kappa_aux[1];
81  }
82  else if (_k_order == 3)
83  {
84  kappa_aux_vec(1) = _kappa_aux[1];
85  kappa_aux_vec(2) = _kappa_aux[2];
86  }
87 
88  r -= dx(_h) * (kappa_aux_vec * _normals[_qp]);
89 
90  return r;
91 }
static InputParameters validParams()
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 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
virtual ADReal computeScalarQpResidual() override
Method for computing the scalar part of residual at quadrature points.
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 ADVariableValue & _lambda
The LM solution.
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.
const VariableTestValue & _test
The shape functions corresponding to the lagrange multiplier variable.
This class enforces a periodic boundary condition between a microscale and macroscale field...
registerMooseObject("MooseApp", ADPeriodicSegmentalConstraint)
std::vector< Point > _normals
the normals
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
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.
ADPeriodicSegmentalConstraint(const InputParameters &parameters)
Interface class ("Veneer") to provide generator methods for derivative material property names...
Moose::VarKindType kind() const
Kind of the variable (Nonlinear, Auxiliary, ...)
virtual ADReal computeQpResidual(Moose::MortarType mortar_type) override
Method for computing the residual at quadrature points.
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
const VariableValue & _kappa_aux
The controlled scalar variable.
const MooseVariableScalar *const _kappa_aux_ptr
(Pointer to) the controlled scalar variable
unsigned int _qp
Definition: Constraint.h:36
unsigned int _h
Used internally to iterate over each scalar component.