www.mooseframework.org
EqualValueConstraint.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 "EqualValueConstraint.h"
11 #include "SubProblem.h"
12 #include "FEProblem.h"
13 #include "Assembly.h"
14 
16 
19 {
21  params.addClassDescription(
22  "EqualValueConstraint enforces solution continuity between secondary and "
23  "primary sides of a mortar interface using lagrange multipliers");
24  params.addRangeCheckedParam<Real>(
25  "delta", 0, "0<=delta<=1", "The coefficient for stabilizing terms");
26  params.addParam<MaterialPropertyName>(
27  "diff_secondary", 1, "The diffusivity on the secondary side");
28  params.addParam<MaterialPropertyName>("diff_primary", 1, "The diffusivity on the primary side");
29  return params;
30 }
31 
33  : ADMortarConstraint(parameters),
34  _lower_secondary_volume(_assembly.lowerDElemVolume()),
35  _lower_primary_volume(_assembly.neighborLowerDElemVolume()),
36  _delta(getParam<Real>("delta")),
37  _diff_secondary(getADMaterialProperty<Real>("diff_secondary")),
38  _diff_primary(getADMaterialProperty<Real>("diff_primary")),
39  _stabilize(_delta > TOLERANCE * TOLERANCE)
40 {
41 }
42 
43 ADReal
45 {
46  switch (mortar_type)
47  {
49  {
50  // The sign choice here makes it so that for the true solution: lambda = normals_secondary *
51  // diff_secondary * grad_u_secondary
52  auto residual = -_lambda[_qp] * _test_secondary[_i][_qp];
53 
54  if (_stabilize)
55  residual += _delta * _lower_secondary_volume *
58 
59  return residual;
60  }
61 
63  {
64  // The sign choice here makes it so that for the true solution: lambda = -normals_primary *
65  // diff_primary * grad_u_primary
66  auto residual = _lambda[_qp] * _test_primary[_i][_qp];
67 
68  if (_stabilize)
69  residual -= _delta * _lower_primary_volume *
72 
73  return residual;
74  }
75 
77  {
78  auto residual = (_u_primary[_qp] - _u_secondary[_qp]) * _test[_i][_qp];
79 
80  if (_stabilize)
81  {
82  // secondary
83  residual -= _delta * _lower_secondary_volume * _test[_i][_qp] *
85 
86  // primary
87  residual -= _delta * _lower_primary_volume * _test[_i][_qp] *
89  }
90 
91  return residual;
92  }
93 
94  default:
95  return 0;
96  }
97 }
const Real & _lower_secondary_volume
The secondary face lower dimensional element (not the mortar element!) volume.
const VariableTestValue & _test_secondary
The shape functions corresponding to the secondary interior primal variable.
const VariableTestGradient & _grad_test_primary
The shape function gradients corresponding to the primary interior primal variable.
const Real _delta
The stabilization parameter.
MortarType
Definition: MooseTypes.h:683
EqualValueConstraint(const InputParameters &parameters)
const ADVariableGradient & _grad_u_secondary
The primal solution gradient on the secondary side.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const Real & _lower_primary_volume
The primary face lower dimensional element volume (not the mortar element!)
unsigned int _i
Definition: Constraint.h:35
const ADVariableGradient & _grad_u_primary
The primal solution gradient on the primary side.
const ADVariableValue & _lambda
The LM solution.
const VariableTestValue & _test
The shape functions corresponding to the lagrange multiplier variable.
Constrain the value of a variable to be the same on both sides of an interface.
registerMooseObject("MooseApp", EqualValueConstraint)
static InputParameters validParams()
std::vector< Point > _normals
the normals
const ADMaterialProperty< Real > & _diff_primary
The diffusion coefficient on the primary side.
DualReal ADReal
Definition: ADRealForward.h:14
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const bool _stabilize
whether to perform stabilization.
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...
const VariableTestGradient & _grad_test_secondary
The shape function gradients corresponding to the secondary interior primal variable.
static InputParameters validParams()
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object...
void addRangeCheckedParam(const std::string &name, const T &value, const std::string &parsed_function, const std::string &doc_string)
ADReal computeQpResidual(Moose::MortarType mortar_type) final
compute the residual at the quadrature points
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.
const ADMaterialProperty< Real > & _diff_secondary
The diffusion coefficient on the secondary side.
unsigned int _qp
Definition: Constraint.h:36
const ADVariableValue & _u_secondary
The primal solution on the secondary side.