Line data Source code
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 :
10 : #include "ComputeFrictionalForceLMMechanicalContact.h"
11 : #include "DisplacedProblem.h"
12 : #include "Assembly.h"
13 : #include "MortarContactUtils.h"
14 : #include "WeightedVelocitiesUserObject.h"
15 :
16 : #include "metaphysicl/metaphysicl_version.h"
17 : #include "metaphysicl/dualsemidynamicsparsenumberarray.h"
18 : #include "metaphysicl/parallel_dualnumber.h"
19 : #if METAPHYSICL_MAJOR_VERSION < 2
20 : #include "metaphysicl/parallel_dynamic_std_array_wrapper.h"
21 : #else
22 : #include "metaphysicl/parallel_dynamic_array_wrapper.h"
23 : #endif
24 : #include "metaphysicl/parallel_semidynamicsparsenumberarray.h"
25 : #include "timpi/parallel_sync.h"
26 :
27 : #include <cmath>
28 :
29 : registerMooseObject("ContactApp", ComputeFrictionalForceLMMechanicalContact);
30 :
31 : InputParameters
32 566 : ComputeFrictionalForceLMMechanicalContact::validParams()
33 : {
34 566 : InputParameters params = ComputeWeightedGapLMMechanicalContact::validParams();
35 566 : params.addClassDescription("Computes the tangential frictional forces");
36 1132 : params.addRequiredCoupledVar("friction_lm", "The frictional Lagrange's multiplier");
37 1132 : params.addCoupledVar("friction_lm_dir",
38 : "The frictional Lagrange's multiplier for an addtional direction.");
39 1132 : params.addParam<FunctionName>(
40 : "function_friction",
41 : "Coupled function to evaluate friction with values from contact pressure and relative "
42 : "tangential velocities");
43 1132 : params.addParam<Real>("c_t", 1e0, "Numerical parameter for tangential constraints");
44 1132 : params.addParam<Real>(
45 : "epsilon",
46 1132 : 1.0e-7,
47 : "Minimum value of contact pressure that will trigger frictional enforcement");
48 1132 : params.addRangeCheckedParam<Real>(
49 : "mu", "mu > 0", "The friction coefficient for the Coulomb friction law");
50 1132 : params.addRequiredParam<UserObjectName>("weighted_velocities_uo",
51 : "The weighted tangential velocities user object.");
52 :
53 566 : return params;
54 0 : }
55 :
56 283 : ComputeFrictionalForceLMMechanicalContact::ComputeFrictionalForceLMMechanicalContact(
57 283 : const InputParameters & parameters)
58 : : ComputeWeightedGapLMMechanicalContact(parameters),
59 283 : _weighted_velocities_uo(getUserObject<WeightedVelocitiesUserObject>("weighted_velocities_uo")),
60 566 : _c_t(getParam<Real>("c_t")),
61 283 : _secondary_x_dot(_secondary_var.adUDot()),
62 283 : _primary_x_dot(_primary_var.adUDotNeighbor()),
63 283 : _secondary_y_dot(adCoupledDot("disp_y")),
64 283 : _primary_y_dot(adCoupledNeighborValueDot("disp_y")),
65 283 : _secondary_z_dot(_has_disp_z ? &adCoupledDot("disp_z") : nullptr),
66 283 : _primary_z_dot(_has_disp_z ? &adCoupledNeighborValueDot("disp_z") : nullptr),
67 566 : _epsilon(getParam<Real>("epsilon")),
68 840 : _mu(isParamValid("function_friction") ? std::numeric_limits<double>::quiet_NaN()
69 831 : : getParam<Real>("mu")),
70 575 : _function_friction(isParamValid("function_friction") ? &getFunction("function_friction")
71 : : nullptr),
72 566 : _has_friction_function(isParamValid("function_friction")),
73 283 : _3d(_has_disp_z)
74 :
75 : {
76 566 : if (parameters.isParamSetByUser("mu") && _has_friction_function)
77 0 : paramError(
78 : "mu",
79 : "Please only provide friction either as a function or as a constant value, but not both.");
80 566 : else if (!parameters.isParamSetByUser("mu") && !_has_friction_function)
81 0 : paramError("mu", "Please provide a value or a function for the coefficient of friction.");
82 :
83 566 : if (!getParam<bool>("use_displaced_mesh"))
84 0 : paramError("use_displaced_mesh",
85 : "'use_displaced_mesh' must be true for the "
86 : "ComputeFrictionalForceLMMechanicalContact object");
87 :
88 466 : if (_3d && !isParamValid("friction_lm_dir"))
89 0 : paramError("friction_lm_dir",
90 : "Three-dimensional mortar frictional contact simulations require an additional "
91 : "frictional Lagrange's multiplier to enforce a second tangential pressure");
92 :
93 283 : _friction_vars.push_back(getVar("friction_lm", 0));
94 :
95 283 : if (_3d)
96 122 : _friction_vars.push_back(getVar("friction_lm_dir", 0));
97 :
98 283 : if (!_friction_vars[0]->isNodal())
99 0 : if (_friction_vars[0]->feType().order != static_cast<Order>(0))
100 0 : paramError(
101 : "friction_lm",
102 : "Frictional contact constraints only support elemental variables of CONSTANT order");
103 283 : }
104 :
105 : void
106 0 : ComputeFrictionalForceLMMechanicalContact::computeQpProperties()
107 : {
108 0 : }
109 :
110 : void
111 0 : ComputeFrictionalForceLMMechanicalContact::computeQpIProperties()
112 : {
113 0 : }
114 :
115 : void
116 67481 : ComputeFrictionalForceLMMechanicalContact::residualSetup()
117 : {
118 67481 : }
119 :
120 : void
121 34759 : ComputeFrictionalForceLMMechanicalContact::post()
122 : {
123 : const auto & dof_to_weighted_tangential_velocity =
124 34759 : _weighted_velocities_uo.dofToWeightedVelocities();
125 :
126 : const std::unordered_map<const DofObject *, std::pair<ADReal, Real>> & dof_to_weighted_gap =
127 34759 : _weighted_gap_uo.dofToWeightedGap();
128 :
129 : // Enforce frictional constraints
130 :
131 146166 : for (const auto & [dof_object, weighted_velocities_pr] : dof_to_weighted_tangential_velocity)
132 : {
133 111407 : if (dof_object->processor_id() != this->processor_id())
134 4042 : continue;
135 :
136 : const auto & [weighted_gap_pr, normalization] =
137 107365 : libmesh_map_find(dof_to_weighted_gap, dof_object);
138 107365 : _weighted_gap_ptr = &weighted_gap_pr;
139 107365 : _normalization_ptr = &normalization;
140 107365 : _tangential_vel_ptr[0] = &(weighted_velocities_pr[0]);
141 :
142 107365 : if (_3d)
143 : {
144 0 : _tangential_vel_ptr[1] = &(weighted_velocities_pr[1]);
145 0 : enforceConstraintOnDof3d(dof_object);
146 : }
147 : else
148 107365 : enforceConstraintOnDof(dof_object);
149 : }
150 34759 : }
151 :
152 : void
153 48319 : ComputeFrictionalForceLMMechanicalContact::incorrectEdgeDroppingPost(
154 : const std::unordered_set<const Node *> & inactive_lm_nodes)
155 : {
156 : const auto & dof_to_weighted_tangential_velocity =
157 48319 : _weighted_velocities_uo.dofToWeightedVelocities();
158 48319 : const auto & dof_to_weighted_gap = _weighted_gap_uo.dofToWeightedGap();
159 : // Enforce frictional complementarity constraints
160 204262 : for (const auto & [dof_object, weighted_velocities_pr] : dof_to_weighted_tangential_velocity)
161 : {
162 : // If node inactive, skip
163 307761 : if ((inactive_lm_nodes.find(static_cast<const Node *>(dof_object)) !=
164 155943 : inactive_lm_nodes.end()) ||
165 155943 : (dof_object->processor_id() != this->processor_id()))
166 4125 : continue;
167 :
168 151818 : _weighted_gap_ptr = &dof_to_weighted_gap.at(dof_object).first;
169 151818 : _normalization_ptr = &dof_to_weighted_gap.at(dof_object).second;
170 151818 : _tangential_vel_ptr[0] = &weighted_velocities_pr[0];
171 :
172 151818 : if (_3d)
173 : {
174 46884 : _tangential_vel_ptr[1] = &weighted_velocities_pr[1];
175 46884 : enforceConstraintOnDof3d(dof_object);
176 : }
177 : else
178 104934 : enforceConstraintOnDof(dof_object);
179 : }
180 48319 : }
181 :
182 : void
183 46884 : ComputeFrictionalForceLMMechanicalContact::enforceConstraintOnDof3d(const DofObject * const dof)
184 : {
185 : using std::max, std::sqrt;
186 :
187 46884 : ComputeWeightedGapLMMechanicalContact::enforceConstraintOnDof(dof);
188 :
189 : // Get normal LM
190 46884 : const auto normal_dof_index = dof->dof_number(_sys.number(), _var->number(), 0);
191 46884 : const ADReal & weighted_gap = *_weighted_gap_ptr;
192 46884 : ADReal contact_pressure = (*_sys.currentSolution())(normal_dof_index);
193 : Moose::derivInsert(contact_pressure.derivatives(), normal_dof_index, 1.);
194 :
195 : // Get friction LMs
196 : std::array<const ADReal *, 2> & tangential_vel = _tangential_vel_ptr;
197 : std::array<dof_id_type, 2> friction_dof_indices;
198 : std::array<ADReal, 2> friction_lm_values;
199 :
200 : const unsigned int num_tangents = 2;
201 140652 : for (const auto i : make_range(num_tangents))
202 : {
203 93768 : friction_dof_indices[i] = dof->dof_number(_sys.number(), _friction_vars[i]->number(), 0);
204 93768 : friction_lm_values[i] = (*_sys.currentSolution())(friction_dof_indices[i]);
205 : Moose::derivInsert(friction_lm_values[i].derivatives(), friction_dof_indices[i], 1.);
206 : }
207 :
208 : // Get normalized c and c_t values (if normalization specified
209 46884 : const Real c = _normalize_c ? _c / *_normalization_ptr : _c;
210 46884 : const Real c_t = _normalize_c ? _c_t / *_normalization_ptr : _c_t;
211 :
212 : // Compute the friction coefficient (constant or function)
213 : ADReal mu_ad = computeFrictionValue(contact_pressure,
214 : _dof_to_real_tangential_velocity[dof][0],
215 46884 : _dof_to_real_tangential_velocity[dof][1]);
216 :
217 : ADReal dof_residual;
218 : ADReal dof_residual_dir;
219 :
220 : // Primal-dual active set strategy (PDASS)
221 46884 : if (contact_pressure < _epsilon)
222 : {
223 13160 : dof_residual = friction_lm_values[0];
224 13160 : dof_residual_dir = friction_lm_values[1];
225 : }
226 : else
227 : {
228 : // Espilon to avoid automatic differentiation singularity
229 : const Real epsilon_sqrt = 1.0e-48;
230 :
231 33724 : const auto lamdba_plus_cg = contact_pressure + c * weighted_gap;
232 : std::array<ADReal, 2> lambda_t_plus_ctu;
233 67448 : lambda_t_plus_ctu[0] = friction_lm_values[0] + c_t * *tangential_vel[0] * _dt;
234 67448 : lambda_t_plus_ctu[1] = friction_lm_values[1] + c_t * *tangential_vel[1] * _dt;
235 :
236 33724 : const auto term_1_x = max(mu_ad * lamdba_plus_cg,
237 33724 : sqrt(lambda_t_plus_ctu[0] * lambda_t_plus_ctu[0] +
238 33724 : lambda_t_plus_ctu[1] * lambda_t_plus_ctu[1] + epsilon_sqrt)) *
239 : friction_lm_values[0];
240 :
241 33724 : const auto term_1_y = max(mu_ad * lamdba_plus_cg,
242 33724 : sqrt(lambda_t_plus_ctu[0] * lambda_t_plus_ctu[0] +
243 33724 : lambda_t_plus_ctu[1] * lambda_t_plus_ctu[1] + epsilon_sqrt)) *
244 : friction_lm_values[1];
245 :
246 33724 : const auto term_2_x = mu_ad * max(0.0, lamdba_plus_cg) * lambda_t_plus_ctu[0];
247 :
248 67448 : const auto term_2_y = mu_ad * max(0.0, lamdba_plus_cg) * lambda_t_plus_ctu[1];
249 :
250 33724 : dof_residual = term_1_x - term_2_x;
251 33724 : dof_residual_dir = term_1_y - term_2_y;
252 : }
253 :
254 46884 : addResidualsAndJacobian(_assembly,
255 93768 : std::array<ADReal, 1>{{dof_residual}},
256 93768 : std::array<dof_id_type, 1>{{friction_dof_indices[0]}},
257 : _friction_vars[0]->scalingFactor());
258 46884 : addResidualsAndJacobian(_assembly,
259 93768 : std::array<ADReal, 1>{{dof_residual_dir}},
260 93768 : std::array<dof_id_type, 1>{{friction_dof_indices[1]}},
261 : _friction_vars[1]->scalingFactor());
262 46884 : }
263 :
264 : void
265 212299 : ComputeFrictionalForceLMMechanicalContact::enforceConstraintOnDof(const DofObject * const dof)
266 : {
267 : using std::abs, std::max;
268 :
269 212299 : ComputeWeightedGapLMMechanicalContact::enforceConstraintOnDof(dof);
270 :
271 : // Get friction LM
272 212299 : const auto friction_dof_index = dof->dof_number(_sys.number(), _friction_vars[0]->number(), 0);
273 212299 : const ADReal & tangential_vel = *_tangential_vel_ptr[0];
274 212299 : ADReal friction_lm_value = (*_sys.currentSolution())(friction_dof_index);
275 : Moose::derivInsert(friction_lm_value.derivatives(), friction_dof_index, 1.);
276 :
277 : // Get normal LM
278 212299 : const auto normal_dof_index = dof->dof_number(_sys.number(), _var->number(), 0);
279 212299 : const ADReal & weighted_gap = *_weighted_gap_ptr;
280 212299 : ADReal contact_pressure = (*_sys.currentSolution())(normal_dof_index);
281 : Moose::derivInsert(contact_pressure.derivatives(), normal_dof_index, 1.);
282 :
283 : // Get normalized c and c_t values (if normalization specified
284 212299 : const Real c = _normalize_c ? _c / *_normalization_ptr : _c;
285 212299 : const Real c_t = _normalize_c ? _c_t / *_normalization_ptr : _c_t;
286 :
287 : // Compute the friction coefficient (constant or function)
288 : ADReal mu_ad =
289 424598 : computeFrictionValue(contact_pressure, _dof_to_real_tangential_velocity[dof][0], 0.0);
290 :
291 : ADReal dof_residual;
292 : // Primal-dual active set strategy (PDASS)
293 212299 : if (contact_pressure < _epsilon)
294 60540 : dof_residual = friction_lm_value;
295 : else
296 : {
297 151759 : const auto term_1 = max(mu_ad * (contact_pressure + c * weighted_gap),
298 303518 : abs(friction_lm_value + c_t * tangential_vel * _dt)) *
299 : friction_lm_value;
300 151759 : const auto term_2 = mu_ad * max(0.0, contact_pressure + c * weighted_gap) *
301 303518 : (friction_lm_value + c_t * tangential_vel * _dt);
302 :
303 151759 : dof_residual = term_1 - term_2;
304 : }
305 :
306 212299 : addResidualsAndJacobian(_assembly,
307 424598 : std::array<ADReal, 1>{{dof_residual}},
308 424598 : std::array<dof_id_type, 1>{{friction_dof_index}},
309 : _friction_vars[0]->scalingFactor());
310 212299 : }
311 :
312 : ADReal
313 259183 : ComputeFrictionalForceLMMechanicalContact::computeFrictionValue(const ADReal & contact_pressure,
314 : const ADReal & tangential_vel,
315 : const ADReal & tangential_vel_dir)
316 : {
317 : using std::sqrt;
318 :
319 : // TODO: Introduce temperature dependence in the function. Do this when we have an example.
320 : ADReal mu_ad;
321 :
322 259183 : if (!_has_friction_function)
323 252607 : mu_ad = _mu;
324 : else
325 : {
326 : ADReal tangential_vel_magnitude =
327 6576 : sqrt(tangential_vel * tangential_vel + tangential_vel_dir * tangential_vel_dir + 1.0e-24);
328 6576 : mu_ad = _function_friction->value<ADReal>(0.0, contact_pressure, tangential_vel_magnitude, 0.0);
329 : }
330 :
331 259183 : return mu_ad;
332 : }
|