https://mooseframework.inl.gov
ADMortarConstraint.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 
10 #include "ADMortarConstraint.h"
11 
12 // MOOSE includes
13 #include "MooseVariable.h"
14 #include "Assembly.h"
15 #include "SystemBase.h"
16 #include "ADUtils.h"
17 
20 {
22  return params;
23 }
24 
26  : MortarConstraintBase(parameters),
27  _lambda_dummy(),
28  _lambda(_var ? _var->adSlnLower() : _lambda_dummy),
29  _u_secondary(_secondary_var.adSln()),
30  _u_primary(_primary_var.adSlnNeighbor()),
31  _grad_u_secondary(_secondary_var.adGradSln()),
32  _grad_u_primary(_primary_var.adGradSlnNeighbor())
33 {
35 }
36 
37 void
39 {
40  unsigned int test_space_size = 0;
41  switch (mortar_type)
42  {
45  test_space_size = _test_secondary.size();
46  break;
47 
50  test_space_size = _test_primary.size();
51  break;
52 
54  mooseAssert(_var, "LM variable is null");
56  test_space_size = _test.size();
57  break;
58  }
59 
60  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
61  for (_i = 0; _i < test_space_size; _i++)
62  _local_re(_i) += raw_value(_JxW_msm[_qp] * _coord[_qp] * computeQpResidual(mortar_type));
63 
65 }
66 
67 void
69 {
70  std::vector<ADReal> residuals;
71  std::size_t test_space_size = 0;
72  typedef Moose::ConstraintJacobianType JType;
73  typedef Moose::MortarType MType;
74  std::vector<JType> jacobian_types;
75  std::vector<dof_id_type> dof_indices;
76  Real scaling_factor = 1;
77 
78  switch (mortar_type)
79  {
80  case MType::Secondary:
81  dof_indices = _secondary_var.dofIndices();
83  scaling_factor = _secondary_var.scalingFactor();
84  break;
85 
86  case MType::Primary:
87  dof_indices = _primary_var.dofIndicesNeighbor();
89  scaling_factor = _primary_var.scalingFactor();
90  break;
91 
92  case MType::Lower:
93  mooseAssert(_var, "The Lagrange Multiplier should be non-null if this is getting called");
94  dof_indices = _var->dofIndicesLower();
96  scaling_factor = _var->scalingFactor();
97  break;
98  }
99  test_space_size = dof_indices.size();
100 
101  residuals.resize(test_space_size, 0);
102  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
103  for (_i = 0; _i < test_space_size; _i++)
104  residuals[_i] += _JxW_msm[_qp] * _coord[_qp] * computeQpResidual(mortar_type);
105 
106  addResidualsAndJacobianWithoutConstraints(_assembly, residuals, dof_indices, scaling_factor);
107 }
108 
109 void
111 {
112  computeJacobian();
113 }
User for mortar methods.
MooseVariableField< Real > & _secondary_var
Reference to the secondary variable.
ADMortarConstraint(const InputParameters &parameters)
const VariableTestValue & _test_secondary
The shape functions corresponding to the secondary interior primal variable.
void computeResidualAndJacobian() override
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
virtual void haveADObjects(bool have_ad_objects)
Method for setting whether we have any ad objects.
Definition: SubProblem.h:767
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
unsigned int number() const
Get variable number coming from libMesh.
MortarType
Definition: MooseTypes.h:770
const std::vector< dof_id_type > & dofIndicesLower() const final
Get dof indices for the current lower dimensional element (this is meaningful when performing mortar ...
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const std::vector< Real > & _JxW_msm
The element Jacobian times weights.
unsigned int _i
Definition: Constraint.h:35
virtual ADReal computeQpResidual(Moose::MortarType mortar_type)=0
compute the residual at the quadrature points
virtual void computeResidual() override
Method for computing the residual.
const libMesh::QBase *const & _qrule_msm
The quadrature rule on the mortar segment element.
const VariableTestValue & _test
The shape functions corresponding to the lagrange multiplier variable.
static InputParameters validParams()
virtual const std::vector< dof_id_type > & dofIndicesNeighbor() const =0
Get neighbor DOF indices for currently selected element.
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
void addResidualsAndJacobianWithoutConstraints(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the provided d...
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
unsigned int n_points() const
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
MooseVariable *const _var
Pointer to the lagrange multipler variable. nullptr if none.
void prepareVectorTagNeighbor(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual the according to active tags for DG and interface kernels...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
ConstraintJacobianType
Definition: MooseTypes.h:796
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
virtual void computeJacobian() override
Method for computing the Jacobian.
static InputParameters validParams()
void prepareVectorTagLower(Assembly &assembly, unsigned int ivar)
Prepare data for computing the residual according to active tags for mortar constraints.
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual according to active tags.
const MooseArray< Real > & _coord
Member for handling change of coordinate systems (xyz, rz, spherical)
const VariableTestValue & _test_primary
The shape functions corresponding to the primary interior primal variable.
unsigned int _qp
Definition: Constraint.h:36
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
MooseVariableField< Real > & _primary_var
Reference to the primary variable.