https://mooseframework.inl.gov
ComputeWeightedGapLMMechanicalContact.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 #include "DisplacedProblem.h"
12 #include "Assembly.h"
13 #include "MortarContactUtils.h"
14 #include "WeightedGapUserObject.h"
15 #include "metaphysicl/metaphysicl_version.h"
16 #include "metaphysicl/dualsemidynamicsparsenumberarray.h"
17 #include "metaphysicl/parallel_dualnumber.h"
18 #if METAPHYSICL_MAJOR_VERSION < 2
19 #include "metaphysicl/parallel_dynamic_std_array_wrapper.h"
20 #else
21 #include "metaphysicl/parallel_dynamic_array_wrapper.h"
22 #endif
23 #include "metaphysicl/parallel_semidynamicsparsenumberarray.h"
24 #include "timpi/parallel_sync.h"
25 
27 
28 namespace
29 {
30 const InputParameters &
31 assignVarsInParams(const InputParameters & params_in)
32 {
33  InputParameters & ret = const_cast<InputParameters &>(params_in);
34  const auto & disp_x_name = ret.get<std::vector<VariableName>>("disp_x");
35  if (disp_x_name.size() != 1)
36  mooseError("We require that the disp_x parameter have exactly one coupled name");
37 
38  // We do this so we don't get any variable errors during MortarConstraint(Base) construction
39  ret.set<VariableName>("secondary_variable") = disp_x_name[0];
40  ret.set<VariableName>("primary_variable") = disp_x_name[0];
41 
42  return ret;
43 }
44 }
45 
48 {
50  params.addClassDescription("Computes the weighted gap that will later be used to enforce the "
51  "zero-penetration mechanical contact conditions");
52  params.suppressParameter<VariableName>("secondary_variable");
53  params.suppressParameter<VariableName>("primary_variable");
54  params.addRequiredCoupledVar("disp_x", "The x displacement variable");
55  params.addRequiredCoupledVar("disp_y", "The y displacement variable");
56  params.addCoupledVar("disp_z", "The z displacement variable");
57  params.addParam<Real>(
58  "c", 1e6, "Parameter for balancing the size of the gap and contact pressure");
59  params.addParam<bool>(
60  "normalize_c",
61  false,
62  "Whether to normalize c by weighting function norm. When unnormalized "
63  "the value of c effectively depends on element size since in the constraint we compare nodal "
64  "Lagrange Multiplier values to integrated gap values (LM nodal value is independent of "
65  "element size, where integrated values are dependent on element size).");
66  params.set<bool>("use_displaced_mesh") = true;
67  params.set<bool>("interpolate_normals") = false;
68  params.addRequiredParam<UserObjectName>("weighted_gap_uo", "The weighted gap user object");
69  return params;
70 }
71 
73  const InputParameters & parameters)
74  : ADMortarConstraint(assignVarsInParams(parameters)),
75  _secondary_disp_x(adCoupledValue("disp_x")),
76  _primary_disp_x(adCoupledNeighborValue("disp_x")),
77  _secondary_disp_y(adCoupledValue("disp_y")),
78  _primary_disp_y(adCoupledNeighborValue("disp_y")),
79  _has_disp_z(isCoupled("disp_z")),
80  _secondary_disp_z(_has_disp_z ? &adCoupledValue("disp_z") : nullptr),
81  _primary_disp_z(_has_disp_z ? &adCoupledNeighborValue("disp_z") : nullptr),
82  _c(getParam<Real>("c")),
83  _normalize_c(getParam<bool>("normalize_c")),
84  _nodal(getVar("disp_x", 0)->feType().family == LAGRANGE),
85  _disp_x_var(getVar("disp_x", 0)),
86  _disp_y_var(getVar("disp_y", 0)),
87  _disp_z_var(_has_disp_z ? getVar("disp_z", 0) : nullptr),
88  _weighted_gap_uo(getUserObject<WeightedGapUserObject>("weighted_gap_uo"))
89 {
90  if (!getParam<bool>("use_displaced_mesh"))
91  paramError(
92  "use_displaced_mesh",
93  "'use_displaced_mesh' must be true for the ComputeWeightedGapLMMechanicalContact object");
94 
95  if (!_var->isNodal())
96  if (_var->feType().order != static_cast<Order>(0))
97  mooseError("Normal contact constraints only support elemental variables of CONSTANT order");
98 }
99 
100 ADReal
102 {
103  mooseError("We should never call computeQpResidual for ComputeWeightedGapLMMechanicalContact");
104 }
105 
106 void
108 {
109 }
110 
111 void
113 {
114 }
115 
116 void
118 {
119 }
120 
121 void
123 {
124 }
125 
126 void
128 {
129 }
130 
131 void
133 {
134  // During "computeResidual" and "computeJacobian" we are actually just computing properties on the
135  // mortar segment element mesh. We are *not* actually assembling into the residual/Jacobian. For
136  // the zero-penetration constraint, the property of interest is the map from node to weighted gap.
137  // Computation of the properties proceeds identically for residual and Jacobian evaluation hence
138  // why we simply call computeResidual here. We will assemble into the residual/Jacobian later from
139  // the post() method
140  computeResidual(mortar_type);
141 }
142 
143 void
145 {
146  parallel_object_only();
147 
148  const auto & dof_to_weighted_gap = _weighted_gap_uo.dofToWeightedGap();
149 
150  for (const auto & [dof_object, weighted_gap_pr] : dof_to_weighted_gap)
151  {
152  if (dof_object->processor_id() != this->processor_id())
153  continue;
154 
155  const auto & [weighted_gap, normalization] = weighted_gap_pr;
156  _weighted_gap_ptr = &weighted_gap;
157  _normalization_ptr = &normalization;
158 
159  enforceConstraintOnDof(dof_object);
160  }
161 }
162 
163 void
165  const std::unordered_set<const Node *> & inactive_lm_nodes)
166 {
167  const auto & dof_to_weighted_gap = _weighted_gap_uo.dofToWeightedGap();
168 
169  for (const auto & [dof_object, weighted_gap_pr] : dof_to_weighted_gap)
170  {
171  if ((inactive_lm_nodes.find(static_cast<const Node *>(dof_object)) !=
172  inactive_lm_nodes.end()) ||
173  (dof_object->processor_id() != this->processor_id()))
174  continue;
175 
176  const auto & [weighted_gap, normalization] = weighted_gap_pr;
177  _weighted_gap_ptr = &weighted_gap;
178  _normalization_ptr = &normalization;
179 
180  enforceConstraintOnDof(dof_object);
181  }
182 }
183 
184 void
186 {
187  const auto & weighted_gap = *_weighted_gap_ptr;
188  const Real c = _normalize_c ? _c / *_normalization_ptr : _c;
189 
190  const auto dof_index = dof->dof_number(_sys.number(), _var->number(), 0);
191  ADReal lm_value = (*_sys.currentSolution())(dof_index);
192  Moose::derivInsert(lm_value.derivatives(), dof_index, 1.);
193 
194  const ADReal dof_residual = std::min(lm_value, weighted_gap * c);
195 
197  std::array<ADReal, 1>{{dof_residual}},
198  std::array<dof_id_type, 1>{{dof_index}},
199  _var->scalingFactor());
200 }
LAGRANGE
virtual const NumericVector< Number > *const & currentSolution() const=0
const libMesh::FEType & feType() const
Order
void addResidualsAndJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
void paramError(const std::string &param, Args... args) const
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
void mooseError(Args &&... args)
unsigned int number() const
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
T & set(const std::string &name, bool quiet_mode=false)
OrderWrapper order
ComputeWeightedGapLMMechanicalContact(const InputParameters &parameters)
const Real _c
This factor multiplies the weighted gap.
const WeightedGapUserObject & _weighted_gap_uo
The weighted gap user object.
DualNumber< Real, DNDerivativeType, true > ADReal
virtual void computeResidual() override
void addRequiredParam(const std::string &name, const std::string &doc_string)
virtual void computeQpProperties()
Computes properties that are functions only of the current quadrature point (_qp), e.g.
void suppressParameter(const std::string &name)
SystemBase & _sys
virtual void computeQpIProperties()
Computes properties that are functions both of _qp and _i, for example the weighted gap...
registerMooseObject("ContactApp", ComputeWeightedGapLMMechanicalContact)
ADReal computeQpResidual(Moose::MortarType mortar_type) final
const ADReal * _weighted_gap_ptr
A pointer members that can be used to help avoid copying ADReals.
unsigned int number() const
const std::unordered_map< const DofObject *, std::pair< ADReal, Real > > & dofToWeightedGap() const
Get the degree of freedom to weighted gap information.
Assembly & _assembly
void incorrectEdgeDroppingPost(const std::unordered_set< const Node *> &inactive_lm_nodes) override
Copy of the post routine but that skips assembling inactive nodes.
void addCoupledVar(const std::string &name, const std::string &doc_string)
bool isNodal() const override
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
MooseVariable *const _var
Creates dof object to weighted gap map.
const bool _normalize_c
Whether to normalize weighted gap by weighting function norm.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void computeJacobian() override
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
static InputParameters validParams()
virtual void enforceConstraintOnDof(const DofObject *const dof)
Method called from post().
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
processor_id_type processor_id() const
Computes the weighted gap that will later be used to enforce the zero-penetration mechanical contact ...
void scalingFactor(const std::vector< Real > &factor)