www.mooseframework.org
MortarConstraint.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "MortarConstraint.h"
11 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "MooseVariable.h"
15 #include "SystemBase.h"
16 
19 {
21  return params;
22 }
23 
25  : MortarConstraintBase(parameters),
26  _lambda_dummy(),
27  _lambda(_var ? _var->slnLower() : _lambda_dummy),
28  _u_secondary(_secondary_var.sln()),
29  _u_primary(_primary_var.slnNeighbor()),
30  _grad_u_secondary(_secondary_var.gradSln()),
31  _grad_u_primary(_primary_var.gradSlnNeighbor()),
32  _phi(nullptr),
33  _grad_phi(nullptr)
34 {
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) += _JxW_msm[_qp] * _coord[_qp] * computeQpResidual(mortar_type);
63 
65 }
66 
67 void
69 {
70  std::size_t test_space_size = 0;
71  typedef Moose::ConstraintJacobianType JType;
72  typedef Moose::MortarType MType;
73  std::array<JType, 3> jacobian_types;
74 
75  switch (mortar_type)
76  {
77  case MType::Secondary:
78  test_space_size = _secondary_var.dofIndices().size();
79  jacobian_types = {
81  break;
82 
83  case MType::Primary:
84  test_space_size = _primary_var.dofIndicesNeighbor().size();
86  break;
87 
88  case MType::Lower:
89  test_space_size = _var ? _var->dofIndicesLower().size() : 0;
91  break;
92  }
93 
94  auto & ce = _assembly.couplingEntries();
95  for (const auto & it : ce)
96  {
97  MooseVariableFEBase & ivariable = *(it.first);
98  MooseVariableFEBase & jvariable = *(it.second);
99 
100  unsigned int ivar = ivariable.number();
101  unsigned int jvar = jvariable.number();
102 
103  switch (mortar_type)
104  {
105  case MType::Secondary:
106  if (ivar != _secondary_var.number())
107  continue;
108  break;
109 
110  case MType::Primary:
111  if (ivar != _primary_var.number())
112  continue;
113  break;
114 
115  case MType::Lower:
116  if (!_var || _var->number() != ivar)
117  continue;
118  break;
119  }
120 
121  std::array<size_t, 3> shape_space_sizes{{jvariable.dofIndices().size(),
122  jvariable.dofIndicesNeighbor().size(),
123  jvariable.dofIndicesLower().size()}};
124  std::array<const VariablePhiValue *, 3> phis;
125  std::array<const VariablePhiGradient *, 3> grad_phis;
126  std::array<const VectorVariablePhiValue *, 3> vector_phis;
127  std::array<const VectorVariablePhiGradient *, 3> vector_grad_phis;
128  if (jvariable.isVector())
129  {
130  const auto & temp_var = static_cast<MooseVariableFE<RealVectorValue> &>(jvariable);
131  vector_phis = {{&temp_var.phiFace(), &temp_var.phiFaceNeighbor(), &temp_var.phiLower()}};
132  vector_grad_phis = {
133  {&temp_var.gradPhiFace(), &temp_var.gradPhiFaceNeighbor(), &temp_var.gradPhiLower()}};
134  }
135  else
136  {
137  const auto & temp_var = static_cast<MooseVariableFE<Real> &>(jvariable);
138  phis = {{&temp_var.phiFace(), &temp_var.phiFaceNeighbor(), &temp_var.phiLower()}};
139  grad_phis = {
140  {&temp_var.gradPhiFace(), &temp_var.gradPhiFaceNeighbor(), &temp_var.gradPhiLower()}};
141  }
142 
143  for (MooseIndex(3) type_index = 0; type_index < 3; ++type_index)
144  {
145  const auto jacobian_type = jacobian_types[type_index];
146 
147  prepareMatrixTagLower(_assembly, ivar, jvar, jacobian_type);
148 
150  if (jvariable.isVector())
151  {
152  _vector_phi = vector_phis[type_index];
153  _vector_grad_phi = vector_grad_phis[type_index];
154  }
155  else
156  {
157  _phi = phis[type_index];
158  _grad_phi = grad_phis[type_index];
159  }
160 
161  for (_i = 0; _i < test_space_size; _i++)
162  for (_j = 0; _j < shape_space_sizes[type_index]; _j++)
163  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
164  _local_ke(_i, _j) +=
165  _JxW_msm[_qp] * _coord[_qp] * computeQpJacobian(jacobian_type, jvar);
167  }
168  }
169 }
User for mortar methods.
MooseVariableField< Real > & _secondary_var
Reference to the secondary variable.
static InputParameters validParams()
const VectorVariablePhiValue * _vector_phi
The current shape functions for vector variables.
std::vector< std::pair< MooseVariableFieldBase *, MooseVariableFieldBase * > > & couplingEntries()
Definition: Assembly.h:1243
const VariableTestValue & _test_secondary
The shape functions corresponding to the secondary interior primal variable.
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
const VariablePhiValue * _phi
The current shape functions.
unsigned int number() const
Get variable number coming from libMesh.
virtual Real computeQpJacobian(Moose::ConstraintJacobianType jacobian_type, unsigned int jvar)=0
compute the jacobian at the quadrature points
MortarType
Definition: MooseTypes.h:683
const VariablePhiGradient * _grad_phi
The current shape function gradients.
const std::vector< dof_id_type > & dofIndicesLower() const final
Get dof indices for the current lower dimensional element (this is meaningful when performing mortar ...
const VectorVariablePhiGradient * _vector_grad_phi
The current shape function gradients for vector variables.
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.
This class provides an interface for common operations on field variables of both FE and FV types wit...
unsigned int _i
Definition: Constraint.h:35
MortarConstraint(const InputParameters &parameters)
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
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()
unsigned int _j
Definition: Constraint.h:35
virtual bool isVector() const =0
virtual const std::vector< dof_id_type > & dofIndicesNeighbor() const =0
Get neighbor DOF indices for currently selected element.
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
unsigned int n_points() const
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
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...
ConstraintJacobianType
Definition: MooseTypes.h:709
void prepareMatrixTagLower(Assembly &assembly, unsigned int ivar, unsigned int jvar, Moose::ConstraintJacobianType type)
Prepare data for computing the jacobian according to the active tags for mortar.
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
virtual void computeJacobian() override
Method for computing the Jacobian.
const FieldVariablePhiValue & phiFace() const override final
Return the variable&#39;s shape functions on an element face.
virtual const std::vector< dof_id_type > & dofIndicesLower() const =0
Get dof indices for the current lower dimensional element (this is meaningful when performing mortar ...
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
virtual Real computeQpResidual(Moose::MortarType mortar_type)=0
compute the residual at the quadrature points
MooseVariableField< Real > & _primary_var
Reference to the primary variable.