https://mooseframework.inl.gov
PenaltyWeightedGapUserObject.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 
12 #include "MooseVariableFE.h"
13 #include "SystemBase.h"
14 
16 
18 
21 {
23  params.addClassDescription("Computes the mortar normal contact force via a penalty approach.");
24  params.addRequiredParam<Real>("penalty", "The penalty factor");
25  params.addRangeCheckedParam<Real>(
26  "penalty_multiplier",
27  1.0,
28  "penalty_multiplier > 0",
29  "The penalty growth factor between augmented Lagrange "
30  "iterations if weighted gap does not get closed fast enough. For frictional simulations, "
31  "values smaller than 100 are recommended, e.g. 5.");
32  params.addParam<bool>(
33  "use_physical_gap",
34  false,
35  "Whether to use the physical normal gap (not scaled by mortar integration) and normalize the "
36  "penalty coefficient with a representative lower-dimensional volume assigned to the node in "
37  "the contacting boundary. This parameter is defaulted to 'true' in the contact action.");
38  params.addRangeCheckedParam<Real>(
39  "penetration_tolerance",
40  "penetration_tolerance > 0",
41  "Acceptable penetration distance at which augmented Lagrange iterations can be stopped");
42  params.addRangeCheckedParam<Real>(
43  "max_penalty_multiplier",
44  1.0e3,
45  "max_penalty_multiplier >= 1.0",
46  "Maximum multiplier applied to penalty factors when adaptivity is used in an augmented "
47  "Lagrange setting. The penalty factor supplied by the user is used as a reference.");
48  MooseEnum adaptivity_penalty_normal("SIMPLE BUSSETTA", "SIMPLE");
49  adaptivity_penalty_normal.addDocumentation(
50  "SIMPLE", "Keep multiplying by the penalty multiplier between AL iterations");
51  adaptivity_penalty_normal.addDocumentation(
52  "BUSSETTA",
53  "Modify the penalty using an algorithm from Bussetta et al, 2012, Comput Mech 49:259-275 "
54  "between AL iterations.");
55  params.addParam<MooseEnum>(
56  "adaptivity_penalty_normal",
57  adaptivity_penalty_normal,
58  "The augmented Lagrange update strategy used on the normal penalty coefficient.");
59  params.addCoupledVar(
60  "aux_lm",
61  "Auxiliary variable that is utilized together with the "
62  "penalty approach to interpolate the resulting contact tractions using dual bases.");
63  params.addParamNamesToGroup("penalty_multiplier penetration_tolerance max_penalty_multiplier",
64  "Augmented Lagrange");
65 
66  return params;
67 }
68 
70  : WeightedGapUserObject(parameters),
72  _penalty(getParam<Real>("penalty")),
73  _penalty_multiplier(getParam<Real>("penalty_multiplier")),
74  _penetration_tolerance(
75  isParamValid("penetration_tolerance") ? getParam<Real>("penetration_tolerance") : 0.0),
76  _augmented_lagrange_problem(
77  dynamic_cast<AugmentedLagrangianContactProblemInterface *>(&_fe_problem)),
78  _lagrangian_iteration_number(_augmented_lagrange_problem
79  ? _augmented_lagrange_problem->getLagrangianIterationNumber()
80  : _no_iterations),
81  _dt(_fe_problem.dt()),
82  _use_physical_gap(getParam<bool>("use_physical_gap")),
83  _aux_lm_var(isCoupled("aux_lm") ? getVar("aux_lm", 0) : nullptr),
84  _max_penalty_multiplier(getParam<Real>("max_penalty_multiplier")),
85  _adaptivity_normal(
86  getParam<MooseEnum>("adaptivity_penalty_normal").getEnum<AdaptivityNormalPenalty>())
87 {
88  auto check_type = [this](const auto & var, const auto & var_name)
89  {
90  if (!var.isNodal())
91  paramError(var_name,
92  "The displacement variables must have degrees of freedom exclusively on "
93  "nodes, e.g. they should probably be of finite element type 'Lagrange'.");
94  };
95  check_type(*_disp_x_var, "disp_x");
96  check_type(*_disp_y_var, "disp_y");
97  if (_has_disp_z)
98  check_type(*_disp_z_var, "disp_z");
99 
100  if (!_augmented_lagrange_problem == isParamValid("penetration_tolerance"))
101  paramError("penetration_tolerance",
102  "This parameter must be supplied if and only if an augmented Lagrange problem "
103  "object is used.");
104 }
105 
106 const VariableTestValue &
108 {
110 }
111 
112 const ADVariableValue &
114 {
115  return _contact_pressure;
116 }
117 
118 void
120 {
121  for (auto & dof_pn : _dof_to_normal_pressure)
122  dof_pn.second = 0.0;
123 }
124 
125 void
127 {
129  selfInitialize();
130 }
131 
132 void
134 {
135  using std::min;
136 
137  // compute new normal pressure for each node
138  for (const auto & [dof_object, wgap] : _dof_to_weighted_gap)
139  {
140  auto penalty = findValue(_dof_to_local_penalty, dof_object, _penalty);
141 
142  // If _use_physical_gap is true we "normalize" the penalty parameter with the surrounding area.
143  auto gap = _use_physical_gap ? adPhysicalGap(wgap) / wgap.second : wgap.first;
144 
145  const auto lagrange_multiplier =
147 
148  // keep the negative normal pressure (compressive per convention here)
149  auto normal_pressure = min(0.0, penalty * gap + lagrange_multiplier);
150 
151  // we switch conventins here and consider positive normal pressure as compressive
152  _dof_to_normal_pressure[dof_object] = -normal_pressure;
153  }
154 }
155 
156 void
158 {
160  selfFinalize();
161 }
162 
163 Real
165 {
166  const auto it = _dof_to_normal_pressure.find(_subproblem.mesh().nodePtr(node->id()));
167 
168  if (it != _dof_to_normal_pressure.end())
169  return MetaPhysicL::raw_value(it->second);
170  else
171  return 0.0;
172 }
173 
174 Real
176 {
177  const auto it = _dof_to_lagrange_multiplier.find(_subproblem.mesh().nodePtr(node->id()));
178 
179  if (it != _dof_to_lagrange_multiplier.end())
180  return -MetaPhysicL::raw_value(it->second);
181  else
182  return 0.0;
183 }
184 
185 void
187 {
189  for (const auto qp : make_range(_qrule_msm->n_points()))
190  _contact_pressure[qp] = 0.0;
191 
192  for (const auto i : make_range(_test->size()))
193  {
194  const Node * const node = _lower_secondary_elem->node_ptr(i);
195 
196  for (const auto qp : make_range(_qrule_msm->n_points()))
197  _contact_pressure[qp] += (*_test)[i][qp] * _dof_to_normal_pressure[node];
198  }
199 }
200 
201 void
203 {
204  // Let's not clear the LMs for improved performance in
205  // nonlinear problems
206 
207  // reset penalty
208  for (auto & dof_lp : _dof_to_local_penalty)
209  dof_lp.second = _penalty;
210 
211  // clear previous gap
212  for (auto & dof_pg : _dof_to_previous_gap)
213  dof_pg.second = 0.0;
214 
215  // save old timestep
216  _dt_old = _dt;
217 }
218 
219 void
221 {
223 }
224 
225 bool
227 {
228  using std::max, std::min, std::abs;
229 
230  Real max_positive_gap = 0.0;
231  Real min_negative_gap = 0.0;
232 
233  // Get maximum gap to ascertain whether we are converged.
234  for (auto & [dof_object, wgap] : _dof_to_weighted_gap)
235  {
236  const auto gap = physicalGap(wgap);
237  {
238  // Check condition for nodes that are active
239  if (gap < 0 || _dof_to_lagrange_multiplier[dof_object] < 0.0)
240  {
241  max_positive_gap = max(max_positive_gap, gap);
242  min_negative_gap = min(min_negative_gap, gap);
243  }
244  }
245  }
246 
247  // Communicate the extreme gap values in parallel.
248  std::vector<Real> recv;
249  if (this->_communicator.rank() == 0)
250  recv.resize(this->_communicator.size());
251  this->_communicator.gather(
252  0, -min_negative_gap > max_positive_gap ? min_negative_gap : max_positive_gap, recv);
253 
254  if (this->_communicator.rank() == 0)
255  {
256  min_negative_gap = *std::min_element(recv.begin(), recv.end());
257  max_positive_gap = *std::max_element(recv.begin(), recv.end());
258 
259  // report the gap value with the largest magnitude
260  const Real reported_gap =
261  -min_negative_gap > max_positive_gap ? min_negative_gap : max_positive_gap;
262  if (abs(reported_gap) > _penetration_tolerance)
263  {
264  mooseInfoRepeated("Penetration tolerance fail max_gap = ",
265  reported_gap,
266  " (gap_tol=",
268  "). Iteration number is: ",
270  ".");
271  return false;
272  }
273  else
274  mooseInfoRepeated("Penetration tolerance success max_gap = ",
275  reported_gap,
276  " (gap_tol=",
278  "). Iteration number is: ",
280  ".");
281  }
282 
283  return true;
284 }
285 
286 void
288 {
289  // Loop over all nodes for which a gap has been computed
290  for (auto & [dof_object, wgap] : _dof_to_weighted_gap)
291  {
292  const Real gap = physicalGap(wgap);
293  // Store previous augmented lagrange iteration gap
294  _dof_to_previous_gap[dof_object] = gap;
295  }
296 }
297 
298 void
300 {
301  using std::abs;
302 
303  for (const auto & [dof_object, wgap] : _dof_to_weighted_gap)
304  {
305  auto & penalty = _dof_to_local_penalty[dof_object];
306  if (penalty == 0.0)
307  penalty = _penalty;
308 
309  const auto gap = getNormalGap(static_cast<const Node *>(dof_object));
310  auto & lagrange_multiplier = _dof_to_lagrange_multiplier[dof_object];
311 
312  const auto possible_normalization =
313  (_use_physical_gap ? libmesh_map_find(_dof_to_weighted_gap, dof_object).second : 1.0);
314 
315  // Update penalty (the factor of 1/4 is suggested in the literature, the limit on AL iteration
316  // caps the penalty increase)
317  // Before we were updating the LM before adapting the penalty factor
318  if (lagrange_multiplier + gap * penalty / possible_normalization <= 0)
319  lagrange_multiplier += gap * penalty / possible_normalization;
320  else
321  lagrange_multiplier = 0.0;
322 
323  const auto previous_gap = _dof_to_previous_gap[dof_object];
324  Real eval_tn = 0;
325 
327  {
328  if (abs(gap) > 0.25 * abs(previous_gap) && abs(gap) > _penetration_tolerance)
329  penalty *= _penalty_multiplier;
330  }
332  bussettaAdaptivePenalty(previous_gap, gap, penalty, eval_tn);
333 
334  // Heuristics to bound the penalty factor
335  if (penalty < _penalty)
336  penalty = _penalty;
337  else if (penalty > _max_penalty_multiplier * _penalty)
338  penalty = _max_penalty_multiplier * _penalty;
339  }
340 }
341 
342 void
344  const Real gap,
345  Real & penalty,
346  Real & eval_tn)
347 {
348  // Positive gaps means no contact
349  if (previous_gap > 0.0)
350  {
351  penalty = _penalty;
352  eval_tn = 0.0;
353  }
354  else
355  {
356  if (previous_gap * gap < 0)
357  eval_tn = 0.0;
358  else
359  eval_tn = penalty * previous_gap;
360 
361  adaptiveNormalPenalty(previous_gap, gap, penalty);
362  }
363 }
364 
365 void
367  const Real gap,
368  Real & penalty)
369 {
370  using std::abs, std::max, std::sqrt;
371 
372  const bool condition_one = abs(abs(previous_gap / gap) - 1.0) < 0.01;
373  const bool condition_two = abs(gap) > 1.01 * abs(previous_gap);
374 
375  if (previous_gap * gap < 0)
376  {
377  if (previous_gap > _penetration_tolerance)
378  penalty = abs(penalty * previous_gap / gap * (abs(gap) + _penetration_tolerance) /
379  (gap - previous_gap));
380  else
381  penalty = abs(penalty * previous_gap / 10.0 / gap);
382  }
383  else if (abs(gap) > _penetration_tolerance)
384  {
385  if (abs(gap - previous_gap) >
386  max(gap / 10.0, max(previous_gap / 10.0, 5.0 * _penetration_tolerance)))
387  penalty *= 10.0;
388  else if (condition_one && gap < 10.0 * _penetration_tolerance)
389  penalty *= MathUtils::pow(sqrt(abs(gap) / _penetration_tolerance - 1.0) + 1.0, 2);
390  else if (condition_two)
391  penalty *= 10.0 * previous_gap / gap;
392  else
393  penalty *= sqrt(abs(gap) / _penetration_tolerance - 1.0) + 1.0;
394  }
395 }
const Real _max_penalty_multiplier
Maximum multiplier applied to the initial penalty factor in AL.
virtual MooseMesh & mesh()=0
ADReal adPhysicalGap(const std::pair< ADReal, Real > &gap) const
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
virtual void augmentedLagrangianSetup() override
static const unsigned int _no_iterations
virtual Real getNormalLagrangeMultiplier(const Node *const node) const
std::unordered_map< const DofObject *, ADReal > _dof_to_normal_pressure
Map from degree of freedom to normal pressure for reporting.
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)
V findValue(const std::unordered_map< K, V > &map, const K &key, const V &default_value=0) const
Find a value in a map or return a default if the key doesn&#39;t exist.
const Real _penalty
The penalty factor.
void bussettaAdaptivePenalty(const Real previous_gap, const Real gap, Real &penalty, Real &eval_tn)
Adaptive, local penalty for AL.
void gather(const unsigned int root_id, const T &send_data, std::vector< T, A > &recv) const
AugmentedLagrangianContactProblemInterface *const _augmented_lagrange_problem
augmented Lagrange problem and iteration number
const bool _use_physical_gap
Use scaled or physical gap.
const Real & _dt
Current delta t... or timestep size.
virtual void updateAugmentedLagrangianMultipliers() override
virtual void initialize() override
const MooseVariable *const _disp_z_var
The z displacement variable.
auto raw_value(const Eigen::Map< T > &in)
static InputParameters validParams()
void mooseInfoRepeated(Args &&... args)
processor_id_type rank() const
Interface class for user objects that support the augmented Lagrange formalism as implemented in Augm...
User object for computing weighted gaps and contact pressure for penalty based mortar constraints...
const bool _has_disp_z
For 2D mortar contact no displacement will be specified, so const pointers used.
const Parallel::Communicator & _communicator
virtual Real getNormalContactPressure(const Node *const node) const override
ADVariableValue _contact_pressure
The contact pressure on the mortar segument quadrature points.
void adaptiveNormalPenalty(const Real previous_gap, const Real gap, Real &penalty)
See Algorithm 3 of &#39;The adapted augmented Lagrangian method: a new method for the resolution of the m...
const libMesh::QBase *const & _qrule_msm
void addRequiredParam(const std::string &name, const std::string &doc_string)
auto max(const L &left, const R &right)
SubProblem & _subproblem
virtual Real getNormalGap(const Node *const) const
Elem const *const & _lower_secondary_elem
virtual const ADVariableValue & contactPressure() const override
const Real _penetration_tolerance
penetration tolerance for augmented Lagrange contact
PenaltyWeightedGapUserObject(const InputParameters &parameters)
OutputTools< Real >::VariableTestValue VariableTestValue
const VariableTestValue * _test
A pointer to the test function associated with the weighted gap.
virtual const Node * nodePtr(const dof_id_type i) const
std::unordered_map< const DofObject *, Real > _dof_to_local_penalty
Map from degree of freedom to local penalty value.
unsigned int n_points() const
virtual const VariableTestValue & test() const override
enum PenaltyWeightedGapUserObject::AdaptivityNormalPenalty _adaptivity_normal
virtual bool isAugmentedLagrangianConverged() override
void addCoupledVar(const std::string &name, const std::string &doc_string)
Creates dof object to weighted gap map.
const MooseVariable *const _aux_lm_var
The auxiliary Lagrange multiplier variable (used together whith the Petrov-Galerkin approach) ...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::unordered_map< const DofObject *, Real > _dof_to_previous_gap
Map from degree of freedom to previous AL iteration gap values.
const Real _penalty_multiplier
penalty growth factor for augmented Lagrange
const MooseVariable *const _disp_y_var
The y displacement variable.
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
Class to provide an interface for parameters and routines required to check convergence for the augme...
virtual const FieldVariablePhiValue & phiLower() const override
std::unordered_map< const DofObject *, Real > _dof_to_lagrange_multiplier
Map from degree of freedom to augmented lagrange multiplier.
IntRange< T > make_range(T beg, T end)
void resize(unsigned int size)
void addClassDescription(const std::string &doc_string)
Real physicalGap(const std::pair< ADReal, Real > &gap) const
Compute physical gap from integration gap quantity.
void addRangeCheckedParam(const std::string &name, const T &value, const std::string &parsed_function, const std::string &doc_string)
bool isParamValid(const std::string &name) const
const unsigned int & _lagrangian_iteration_number
const MooseVariable *const _disp_x_var
The x displacement variable.
T pow(T x, int e)
AdaptivityNormalPenalty
The adaptivity method for the penalty factor at augmentations.
std::unordered_map< const DofObject *, std::pair< ADReal, Real > > _dof_to_weighted_gap
A map from node to weighted gap and normalization (if requested)
auto min(const L &left, const R &right)
registerMooseObject("ContactApp", PenaltyWeightedGapUserObject)
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
virtual void finalize() override