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  // compute new normal pressure for each node
136  for (const auto & [dof_object, wgap] : _dof_to_weighted_gap)
137  {
138  auto penalty = findValue(_dof_to_local_penalty, dof_object, _penalty);
139 
140  // If _use_physical_gap is true we "normalize" the penalty parameter with the surrounding area.
141  auto gap = _use_physical_gap ? adPhysicalGap(wgap) / wgap.second : wgap.first;
142 
143  const auto lagrange_multiplier =
145 
146  // keep the negative normal pressure (compressive per convention here)
147  auto normal_pressure = std::min(0.0, penalty * gap + lagrange_multiplier);
148 
149  // we switch conventins here and consider positive normal pressure as compressive
150  _dof_to_normal_pressure[dof_object] = -normal_pressure;
151  }
152 }
153 
154 void
156 {
158  selfFinalize();
159 }
160 
161 Real
163 {
164  const auto it = _dof_to_normal_pressure.find(_subproblem.mesh().nodePtr(node->id()));
165 
166  if (it != _dof_to_normal_pressure.end())
167  return MetaPhysicL::raw_value(it->second);
168  else
169  return 0.0;
170 }
171 
172 Real
174 {
175  const auto it = _dof_to_lagrange_multiplier.find(_subproblem.mesh().nodePtr(node->id()));
176 
177  if (it != _dof_to_lagrange_multiplier.end())
178  return -MetaPhysicL::raw_value(it->second);
179  else
180  return 0.0;
181 }
182 
183 void
185 {
187  for (const auto qp : make_range(_qrule_msm->n_points()))
188  _contact_pressure[qp] = 0.0;
189 
190  for (const auto i : make_range(_test->size()))
191  {
192  const Node * const node = _lower_secondary_elem->node_ptr(i);
193 
194  for (const auto qp : make_range(_qrule_msm->n_points()))
195  _contact_pressure[qp] += (*_test)[i][qp] * _dof_to_normal_pressure[node];
196  }
197 }
198 
199 void
201 {
202  // Let's not clear the LMs for improved performance in
203  // nonlinear problems
204 
205  // reset penalty
206  for (auto & dof_lp : _dof_to_local_penalty)
207  dof_lp.second = _penalty;
208 
209  // clear previous gap
210  for (auto & dof_pg : _dof_to_previous_gap)
211  dof_pg.second = 0.0;
212 
213  // save old timestep
214  _dt_old = _dt;
215 }
216 
217 void
219 {
221 }
222 
223 bool
225 {
226  Real max_positive_gap = 0.0;
227  Real min_negative_gap = 0.0;
228 
229  // Get maximum gap to ascertain whether we are converged.
230  for (auto & [dof_object, wgap] : _dof_to_weighted_gap)
231  {
232  const auto gap = physicalGap(wgap);
233  {
234  // Check condition for nodes that are active
235  if (gap < 0 || _dof_to_lagrange_multiplier[dof_object] < 0.0)
236  {
237  max_positive_gap = std::max(max_positive_gap, gap);
238  min_negative_gap = std::min(min_negative_gap, gap);
239  }
240  }
241  }
242 
243  // Communicate the extreme gap values in parallel.
244  std::vector<Real> recv;
245  if (this->_communicator.rank() == 0)
246  recv.resize(this->_communicator.size());
247  this->_communicator.gather(
248  0, -min_negative_gap > max_positive_gap ? min_negative_gap : max_positive_gap, recv);
249 
250  if (this->_communicator.rank() == 0)
251  {
252  min_negative_gap = *std::min_element(recv.begin(), recv.end());
253  max_positive_gap = *std::max_element(recv.begin(), recv.end());
254 
255  // report the gap value with the largest magnitude
256  const Real reported_gap =
257  -min_negative_gap > max_positive_gap ? min_negative_gap : max_positive_gap;
258  if (std::abs(reported_gap) > _penetration_tolerance)
259  {
260  mooseInfoRepeated("Penetration tolerance fail max_gap = ",
261  reported_gap,
262  " (gap_tol=",
264  "). Iteration number is: ",
266  ".");
267  return false;
268  }
269  else
270  mooseInfoRepeated("Penetration tolerance success max_gap = ",
271  reported_gap,
272  " (gap_tol=",
274  "). Iteration number is: ",
276  ".");
277  }
278 
279  return true;
280 }
281 
282 void
284 {
285  // Loop over all nodes for which a gap has been computed
286  for (auto & [dof_object, wgap] : _dof_to_weighted_gap)
287  {
288  const Real gap = physicalGap(wgap);
289  // Store previous augmented lagrange iteration gap
290  _dof_to_previous_gap[dof_object] = gap;
291  }
292 }
293 
294 void
296 {
297  for (const auto & [dof_object, wgap] : _dof_to_weighted_gap)
298  {
299  auto & penalty = _dof_to_local_penalty[dof_object];
300  if (penalty == 0.0)
301  penalty = _penalty;
302 
303  const auto gap = getNormalGap(static_cast<const Node *>(dof_object));
304  auto & lagrange_multiplier = _dof_to_lagrange_multiplier[dof_object];
305 
306  const auto possible_normalization =
307  (_use_physical_gap ? libmesh_map_find(_dof_to_weighted_gap, dof_object).second : 1.0);
308 
309  // Update penalty (the factor of 1/4 is suggested in the literature, the limit on AL iteration
310  // caps the penalty increase)
311  // Before we were updating the LM before adapting the penalty factor
312  if (lagrange_multiplier + gap * penalty / possible_normalization <= 0)
313  lagrange_multiplier += gap * penalty / possible_normalization;
314  else
315  lagrange_multiplier = 0.0;
316 
317  const auto previous_gap = _dof_to_previous_gap[dof_object];
318  Real eval_tn = 0;
319 
321  {
322  if (std::abs(gap) > 0.25 * std::abs(previous_gap) && std::abs(gap) > _penetration_tolerance)
323  penalty *= _penalty_multiplier;
324  }
326  bussettaAdaptivePenalty(previous_gap, gap, penalty, eval_tn);
327 
328  // Heuristics to bound the penalty factor
329  if (penalty < _penalty)
330  penalty = _penalty;
331  else if (penalty > _max_penalty_multiplier * _penalty)
332  penalty = _max_penalty_multiplier * _penalty;
333  }
334 }
335 
336 void
338  const Real gap,
339  Real & penalty,
340  Real & eval_tn)
341 {
342  // Positive gaps means no contact
343  if (previous_gap > 0.0)
344  {
345  penalty = _penalty;
346  eval_tn = 0.0;
347  }
348  else
349  {
350  if (previous_gap * gap < 0)
351  eval_tn = 0.0;
352  else
353  eval_tn = penalty * previous_gap;
354 
355  adaptiveNormalPenalty(previous_gap, gap, penalty);
356  }
357 }
358 
359 void
361  const Real gap,
362  Real & penalty)
363 {
364  const bool condition_one = std::abs(std::abs(previous_gap / gap) - 1.0) < 0.01;
365  const bool condition_two = std::abs(gap) > 1.01 * std::abs(previous_gap);
366 
367  if (previous_gap * gap < 0)
368  {
369  if (previous_gap > _penetration_tolerance)
370  penalty = std::abs(penalty * previous_gap / gap * (std::abs(gap) + _penetration_tolerance) /
371  (gap - previous_gap));
372  else
373  penalty = std::abs(penalty * previous_gap / 10.0 / gap);
374  }
375  else if (std::abs(gap) > _penetration_tolerance)
376  {
377  if (std::abs(gap - previous_gap) >
378  std::max(gap / 10.0, std::max(previous_gap / 10.0, 5.0 * _penetration_tolerance)))
379  penalty *= 10.0;
380  else if (condition_one && gap < 10.0 * _penetration_tolerance)
381  penalty *= MathUtils::pow(std::sqrt(std::abs(gap) / _penetration_tolerance - 1.0) + 1.0, 2);
382  else if (condition_two)
383  penalty *= 10.0 * previous_gap / gap;
384  else
385  penalty *= std::sqrt(std::abs(gap) / _penetration_tolerance - 1.0) + 1.0;
386  }
387 }
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
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 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)
SubProblem & _subproblem
bool isParamValid(const std::string &name) const
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
void paramError(const std::string &param, Args... args) 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.
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)
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)
registerMooseObject("ContactApp", PenaltyWeightedGapUserObject)
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
virtual void finalize() override