LCOV - code coverage report
Current view: top level - src/constraints - ComputeDynamicWeightedGapLMMechanicalContact.C (source / functions) Hit Total Coverage
Test: idaholab/moose contact: #31730 (e8b711) with base e0c998 Lines: 173 188 92.0 %
Date: 2025-10-29 16:50:33 Functions: 14 17 82.4 %
Legend: Lines: hit not hit

          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 "ComputeDynamicWeightedGapLMMechanicalContact.h"
      11             : #include "MortarContactUtils.h"
      12             : #include "DisplacedProblem.h"
      13             : #include "Assembly.h"
      14             : #include "AutomaticMortarGeneration.h"
      15             : 
      16             : #include "metaphysicl/dualsemidynamicsparsenumberarray.h"
      17             : #include "metaphysicl/parallel_dualnumber.h"
      18             : #include "metaphysicl/parallel_dynamic_std_array_wrapper.h"
      19             : #include "metaphysicl/parallel_semidynamicsparsenumberarray.h"
      20             : #include "timpi/parallel_sync.h"
      21             : 
      22             : registerMooseObject("ContactApp", ComputeDynamicWeightedGapLMMechanicalContact);
      23             : 
      24             : namespace
      25             : {
      26             : const InputParameters &
      27         146 : assignVarsInParamsDyn(const InputParameters & params_in)
      28             : {
      29             :   InputParameters & ret = const_cast<InputParameters &>(params_in);
      30         146 :   const auto & disp_x_name = ret.get<std::vector<VariableName>>("disp_x");
      31         146 :   if (disp_x_name.size() != 1)
      32           0 :     mooseError("We require that the disp_x parameter have exactly one coupled name");
      33             : 
      34             :   // We do this so we don't get any variable errors during MortarConstraint(Base) construction
      35         292 :   ret.set<VariableName>("secondary_variable") = disp_x_name[0];
      36         146 :   ret.set<VariableName>("primary_variable") = disp_x_name[0];
      37             : 
      38         146 :   return ret;
      39             : }
      40             : }
      41             : InputParameters
      42         292 : ComputeDynamicWeightedGapLMMechanicalContact::validParams()
      43             : {
      44         292 :   InputParameters params = ADMortarConstraint::validParams();
      45         292 :   params.addClassDescription(
      46             :       "Computes the normal contact mortar constraints for dynamic simulations");
      47         876 :   params.addRangeCheckedParam<Real>("capture_tolerance",
      48         584 :                                     1.0e-5,
      49             :                                     "capture_tolerance>=0",
      50             :                                     "Parameter describing a gap threshold for the application of "
      51             :                                     "the persistency constraint in dynamic simulations.");
      52         584 :   params.addCoupledVar("wear_depth",
      53             :                        "The name of the mortar auxiliary variable that is used to modify the "
      54             :                        "weighted gap definition");
      55         584 :   params.addRequiredRangeCheckedParam<Real>(
      56             :       "newmark_beta", "newmark_beta > 0", "Beta parameter for the Newmark time integrator");
      57         584 :   params.addRequiredRangeCheckedParam<Real>(
      58             :       "newmark_gamma", "newmark_gamma >= 0.0", "Gamma parameter for the Newmark time integrator");
      59         292 :   params.suppressParameter<VariableName>("secondary_variable");
      60         292 :   params.suppressParameter<VariableName>("primary_variable");
      61         584 :   params.addRequiredCoupledVar("disp_x", "The x displacement variable");
      62         584 :   params.addRequiredCoupledVar("disp_y", "The y displacement variable");
      63         584 :   params.addCoupledVar("disp_z", "The z displacement variable");
      64         584 :   params.addParam<Real>(
      65         584 :       "c", 1e6, "Parameter for balancing the size of the gap and contact pressure");
      66         584 :   params.addParam<bool>(
      67             :       "normalize_c",
      68         584 :       false,
      69             :       "Whether to normalize c by weighting function norm. When unnormalized "
      70             :       "the value of c effectively depends on element size since in the constraint we compare nodal "
      71             :       "Lagrange Multiplier values to integrated gap values (LM nodal value is independent of "
      72             :       "element size, where integrated values are dependent on element size).");
      73         292 :   params.set<bool>("use_displaced_mesh") = true;
      74         292 :   params.set<bool>("interpolate_normals") = false;
      75         292 :   return params;
      76           0 : }
      77             : 
      78         146 : ComputeDynamicWeightedGapLMMechanicalContact::ComputeDynamicWeightedGapLMMechanicalContact(
      79         146 :     const InputParameters & parameters)
      80             :   : ADMortarConstraint(assignVarsInParamsDyn(parameters)),
      81         146 :     _secondary_disp_x(adCoupledValue("disp_x")),
      82         146 :     _primary_disp_x(adCoupledNeighborValue("disp_x")),
      83         146 :     _secondary_disp_y(adCoupledValue("disp_y")),
      84         146 :     _primary_disp_y(adCoupledNeighborValue("disp_y")),
      85         146 :     _has_disp_z(isCoupled("disp_z")),
      86         146 :     _secondary_disp_z(_has_disp_z ? &adCoupledValue("disp_z") : nullptr),
      87         146 :     _primary_disp_z(_has_disp_z ? &adCoupledNeighborValue("disp_z") : nullptr),
      88         292 :     _c(getParam<Real>("c")),
      89         292 :     _normalize_c(getParam<bool>("normalize_c")),
      90         146 :     _nodal(getVar("disp_x", 0)->feType().family == LAGRANGE),
      91         146 :     _disp_x_var(getVar("disp_x", 0)),
      92         146 :     _disp_y_var(getVar("disp_y", 0)),
      93         180 :     _disp_z_var(_has_disp_z ? getVar("disp_z", 0) : nullptr),
      94         146 :     _capture_tolerance(getParam<Real>("capture_tolerance")),
      95         146 :     _secondary_x_dot(adCoupledDot("disp_x")),
      96         146 :     _primary_x_dot(adCoupledNeighborValueDot("disp_x")),
      97         146 :     _secondary_y_dot(adCoupledDot("disp_y")),
      98         146 :     _primary_y_dot(adCoupledNeighborValueDot("disp_y")),
      99         146 :     _secondary_z_dot(_has_disp_z ? &adCoupledDot("disp_z") : nullptr),
     100         146 :     _primary_z_dot(_has_disp_z ? &adCoupledNeighborValueDot("disp_z") : nullptr),
     101         292 :     _has_wear(isParamValid("wear_depth")),
     102         187 :     _wear_depth(_has_wear ? coupledValueLower("wear_depth") : _zero),
     103         292 :     _newmark_beta(getParam<Real>("newmark_beta")),
     104         584 :     _newmark_gamma(getParam<Real>("newmark_gamma"))
     105             : {
     106         146 :   if (!useDual())
     107           0 :     mooseError("Dynamic mortar contact constraints requires the use of Lagrange multipliers dual "
     108             :                "interpolation");
     109         146 : }
     110             : 
     111             : void
     112      801520 : ComputeDynamicWeightedGapLMMechanicalContact::computeQpProperties()
     113             : {
     114             :   // Trim interior node variable derivatives
     115             :   const auto & primary_ip_lowerd_map = amg().getPrimaryIpToLowerElementMap(
     116      801520 :       *_lower_primary_elem, *_lower_primary_elem->interior_parent(), *_lower_secondary_elem);
     117             :   const auto & secondary_ip_lowerd_map =
     118      801520 :       amg().getSecondaryIpToLowerElementMap(*_lower_secondary_elem);
     119             : 
     120      801520 :   std::array<const MooseVariable *, 3> var_array{{_disp_x_var, _disp_y_var, _disp_z_var}};
     121             :   std::array<ADReal, 3> primary_disp{
     122      801520 :       {_primary_disp_x[_qp], _primary_disp_y[_qp], _has_disp_z ? (*_primary_disp_z)[_qp] : 0}};
     123      801520 :   std::array<ADReal, 3> secondary_disp{{_secondary_disp_x[_qp],
     124      801520 :                                         _secondary_disp_y[_qp],
     125      801520 :                                         _has_disp_z ? (*_secondary_disp_z)[_qp] : 0}};
     126             : 
     127      801520 :   trimInteriorNodeDerivatives(primary_ip_lowerd_map, var_array, primary_disp, false);
     128      801520 :   trimInteriorNodeDerivatives(secondary_ip_lowerd_map, var_array, secondary_disp, true);
     129             : 
     130             :   const ADReal & prim_x = primary_disp[0];
     131             :   const ADReal & prim_y = primary_disp[1];
     132             :   const ADReal * prim_z = nullptr;
     133      801520 :   if (_has_disp_z)
     134             :     prim_z = &primary_disp[2];
     135             : 
     136             :   const ADReal & sec_x = secondary_disp[0];
     137             :   const ADReal & sec_y = secondary_disp[1];
     138             :   const ADReal * sec_z = nullptr;
     139      801520 :   if (_has_disp_z)
     140             :     sec_z = &secondary_disp[2];
     141             : 
     142             :   std::array<ADReal, 3> primary_disp_dot{
     143      801520 :       {_primary_x_dot[_qp], _primary_y_dot[_qp], _has_disp_z ? (*_primary_z_dot)[_qp] : 0}};
     144             :   std::array<ADReal, 3> secondary_disp_dot{
     145      801520 :       {_secondary_x_dot[_qp], _secondary_y_dot[_qp], _has_disp_z ? (*_secondary_z_dot)[_qp] : 0}};
     146             : 
     147      801520 :   trimInteriorNodeDerivatives(primary_ip_lowerd_map, var_array, primary_disp_dot, false);
     148      801520 :   trimInteriorNodeDerivatives(secondary_ip_lowerd_map, var_array, secondary_disp_dot, true);
     149             : 
     150             :   const ADReal & prim_x_dot = primary_disp_dot[0];
     151             :   const ADReal & prim_y_dot = primary_disp_dot[1];
     152             :   const ADReal * prim_z_dot = nullptr;
     153      801520 :   if (_has_disp_z)
     154             :     prim_z_dot = &primary_disp_dot[2];
     155             : 
     156             :   const ADReal & sec_x_dot = secondary_disp_dot[0];
     157             :   const ADReal & sec_y_dot = secondary_disp_dot[1];
     158             :   const ADReal * sec_z_dot = nullptr;
     159      801520 :   if (_has_disp_z)
     160             :     sec_z_dot = &secondary_disp_dot[2];
     161             : 
     162             :   // Compute dynamic constraint-related quantities
     163      801520 :   ADRealVectorValue gap_vec = _phys_points_primary[_qp] - _phys_points_secondary[_qp];
     164      801520 :   gap_vec(0).derivatives() = prim_x.derivatives() - sec_x.derivatives();
     165      801520 :   gap_vec(1).derivatives() = prim_y.derivatives() - sec_y.derivatives();
     166             : 
     167      801520 :   _relative_velocity = ADRealVectorValue(prim_x_dot - sec_x_dot, prim_y_dot - sec_y_dot, 0.0);
     168             : 
     169      801520 :   if (_has_disp_z)
     170             :   {
     171      550400 :     gap_vec(2).derivatives() = prim_z->derivatives() - sec_z->derivatives();
     172      550400 :     _relative_velocity(2) = *prim_z_dot - *sec_z_dot;
     173             :   }
     174             : 
     175      801520 :   _qp_gap_nodal = gap_vec * (_JxW_msm[_qp] * _coord[_qp]);
     176      801520 :   _qp_velocity = _relative_velocity * (_JxW_msm[_qp] * _coord[_qp]);
     177             : 
     178             :   // Current part of the gap velocity Newmark-beta time discretization
     179             :   _qp_gap_nodal_dynamics =
     180      801520 :       (_newmark_gamma / _newmark_beta * gap_vec / _dt) * (_JxW_msm[_qp] * _coord[_qp]);
     181             : 
     182             :   // To do normalization of constraint coefficient (c_n)
     183      801520 :   _qp_factor = _JxW_msm[_qp] * _coord[_qp];
     184      801520 : }
     185             : 
     186             : ADReal
     187           0 : ComputeDynamicWeightedGapLMMechanicalContact::computeQpResidual(Moose::MortarType)
     188             : {
     189           0 :   mooseError(
     190             :       "We should never call computeQpResidual for ComputeDynamicWeightedGapLMMechanicalContact");
     191             : }
     192             : 
     193             : void
     194     2703840 : ComputeDynamicWeightedGapLMMechanicalContact::computeQpIProperties()
     195             : {
     196             :   mooseAssert(_normals.size() == _lower_secondary_elem->n_nodes(),
     197             :               "Making sure that _normals is the expected size");
     198             : 
     199             :   // Get the _dof_to_weighted_gap map
     200     2703840 :   const DofObject * dof = _var->isNodal()
     201     2703840 :                               ? static_cast<const DofObject *>(_lower_secondary_elem->node_ptr(_i))
     202           0 :                               : static_cast<const DofObject *>(_lower_secondary_elem);
     203             : 
     204             :   // Regular normal contact constraint: Use before contact is established for contact detection
     205     5407680 :   _dof_to_weighted_gap[dof].first += _test[_i][_qp] * (_qp_gap_nodal * _normals[_i]);
     206             : 
     207             :   // Integrated part of the "persistency" constraint
     208     5407680 :   _dof_to_weighted_gap_dynamics[dof] += _test[_i][_qp] * _qp_gap_nodal_dynamics * _normals[_i];
     209     5407680 :   _dof_to_velocity[dof] += _test[_i][_qp] * _qp_velocity * _normals[_i];
     210             : 
     211     2703840 :   _dof_to_nodal_wear_depth[dof] += _test[_i][_qp] * _wear_depth[_qp] * _JxW_msm[_qp] * _coord[_qp];
     212             : 
     213     2703840 :   if (_normalize_c)
     214       89920 :     _dof_to_weighted_gap[dof].second += _test[_i][_qp] * _qp_factor;
     215     2703840 : }
     216             : 
     217             : void
     218         640 : ComputeDynamicWeightedGapLMMechanicalContact::timestepSetup()
     219             : {
     220             :   // These dof maps are not recoverable as they are maps of pointers, and recovering old pointers
     221             :   // would be wrong. We would need to create a custom dataStore() and dataLoad()
     222         640 :   if (_app.isRecovering())
     223           0 :     mooseError("This object does not support recovering");
     224             : 
     225             :   _dof_to_old_weighted_gap.clear();
     226             :   _dof_to_old_velocity.clear();
     227             :   _dof_to_nodal_old_wear_depth.clear();
     228             : 
     229        4300 :   for (auto & map_pr : _dof_to_weighted_gap)
     230        3660 :     _dof_to_old_weighted_gap.emplace(map_pr.first, std::move(map_pr.second.first));
     231             : 
     232        4300 :   for (auto & map_pr : _dof_to_velocity)
     233             :     _dof_to_old_velocity.emplace(map_pr);
     234             : 
     235        4300 :   for (auto & map_pr : _dof_to_nodal_wear_depth)
     236             :     _dof_to_nodal_old_wear_depth.emplace(map_pr);
     237         640 : }
     238             : 
     239             : void
     240       11137 : ComputeDynamicWeightedGapLMMechanicalContact::residualSetup()
     241             : {
     242             :   _dof_to_weighted_gap.clear();
     243             :   _dof_to_weighted_gap_dynamics.clear();
     244             :   _dof_to_velocity.clear();
     245             : 
     246             :   // Wear
     247             :   _dof_to_nodal_wear_depth.clear();
     248       11137 : }
     249             : 
     250             : void
     251         663 : ComputeDynamicWeightedGapLMMechanicalContact::post()
     252             : {
     253         663 :   Moose::Mortar::Contact::communicateGaps(
     254         663 :       _dof_to_weighted_gap, _mesh, _nodal, _normalize_c, _communicator, false);
     255             : 
     256         663 :   if (_has_wear)
     257           0 :     communicateWear();
     258             : 
     259             :   // There is a need for the dynamic constraint to uncouple the computation of the weighted gap from
     260             :   // the computation of the constraint itself since we are switching from gap constraint to
     261             :   // persistency constraint.
     262        4578 :   for (const auto & pr : _dof_to_weighted_gap)
     263             :   {
     264        3915 :     if (pr.first->processor_id() != this->processor_id())
     265             :       continue;
     266             : 
     267             :     //
     268        7830 :     _dof_to_weighted_gap[pr.first].first += _dof_to_nodal_wear_depth[pr.first];
     269             :     _dof_to_weighted_gap_dynamics[pr.first] +=
     270        7830 :         _newmark_gamma / _newmark_beta * _dof_to_nodal_wear_depth[pr.first] / _dt;
     271             :     //
     272             : 
     273             :     const auto is_dof_on_map = _dof_to_old_weighted_gap.find(pr.first);
     274             : 
     275             :     // If is_dof_on_map isn't on map, it means it's an initial step
     276        4851 :     if (is_dof_on_map == _dof_to_old_weighted_gap.end() ||
     277             :         _dof_to_old_weighted_gap[pr.first] > _capture_tolerance)
     278        2979 :       _weighted_gap_ptr = &pr.second.first;
     279             :     else
     280             :     {
     281        1872 :       ADReal term = -_newmark_gamma / _newmark_beta / _dt * _dof_to_old_weighted_gap[pr.first];
     282         936 :       term += _dof_to_old_velocity[pr.first];
     283         936 :       _dof_to_weighted_gap_dynamics[pr.first] += term;
     284         936 :       _weighted_gap_ptr = &_dof_to_weighted_gap_dynamics[pr.first];
     285             :     }
     286             : 
     287        3915 :     _normalization_ptr = &pr.second.second;
     288             : 
     289        3915 :     ComputeDynamicWeightedGapLMMechanicalContact::enforceConstraintOnDof(pr.first);
     290             :   }
     291         663 : }
     292             : 
     293             : void
     294       10474 : ComputeDynamicWeightedGapLMMechanicalContact::incorrectEdgeDroppingPost(
     295             :     const std::unordered_set<const Node *> & inactive_lm_nodes)
     296             : {
     297       10474 :   Moose::Mortar::Contact::communicateGaps(
     298       10474 :       _dof_to_weighted_gap, _mesh, _nodal, _normalize_c, _communicator, false);
     299             : 
     300       10474 :   if (_has_wear)
     301        1667 :     communicateWear();
     302             : 
     303       89589 :   for (const auto & pr : _dof_to_weighted_gap)
     304             :   {
     305       79115 :     if ((inactive_lm_nodes.find(static_cast<const Node *>(pr.first)) != inactive_lm_nodes.end()) ||
     306       78693 :         (pr.first->processor_id() != this->processor_id()))
     307             :       continue;
     308             : 
     309             :     //
     310      155802 :     _dof_to_weighted_gap[pr.first].first += _dof_to_nodal_wear_depth[pr.first];
     311             :     _dof_to_weighted_gap_dynamics[pr.first] +=
     312      155802 :         _newmark_gamma / _newmark_beta * _dof_to_nodal_wear_depth[pr.first] / _dt;
     313             :     const auto is_dof_on_map = _dof_to_old_weighted_gap.find(pr.first);
     314             : 
     315             :     // If is_dof_on_map isn't on map, it means it's an initial step
     316      137503 :     if (is_dof_on_map == _dof_to_old_weighted_gap.end() ||
     317             :         _dof_to_old_weighted_gap[pr.first] > _capture_tolerance)
     318             :     {
     319             :       // If this is the first step or the previous step gap is not identified as in contact, apply
     320             :       // regular conditions
     321       60190 :       _weighted_gap_ptr = &pr.second.first;
     322             :     }
     323             :     else
     324             :     {
     325       17711 :       ADReal term = _dof_to_weighted_gap[pr.first].first * _newmark_gamma / (_newmark_beta * _dt);
     326       35422 :       term -= _dof_to_old_weighted_gap[pr.first] * _newmark_gamma / (_newmark_beta * _dt);
     327       17711 :       term -= _dof_to_old_velocity[pr.first];
     328       17711 :       _dof_to_weighted_gap_dynamics[pr.first] = term;
     329             :       // Enable the application of persistency condition
     330       17711 :       _weighted_gap_ptr = &_dof_to_weighted_gap_dynamics[pr.first];
     331             :     }
     332             : 
     333       77901 :     _normalization_ptr = &pr.second.second;
     334             : 
     335       77901 :     ComputeDynamicWeightedGapLMMechanicalContact::enforceConstraintOnDof(pr.first);
     336             :   }
     337       10474 : }
     338             : 
     339             : void
     340        1667 : ComputeDynamicWeightedGapLMMechanicalContact::communicateWear()
     341             : {
     342             :   // We may have wear depth information that should go to other processes that own the dofs
     343             :   using Datum = std::pair<dof_id_type, ADReal>;
     344             :   std::unordered_map<processor_id_type, std::vector<Datum>> push_data;
     345             : 
     346       14031 :   for (auto & pr : _dof_to_nodal_wear_depth)
     347             :   {
     348       12364 :     const auto * const dof_object = pr.first;
     349       12364 :     const auto proc_id = dof_object->processor_id();
     350       12364 :     if (proc_id == this->processor_id())
     351       12364 :       continue;
     352             : 
     353           0 :     push_data[proc_id].push_back(std::make_pair(dof_object->id(), std::move(pr.second)));
     354             :   }
     355             : 
     356        1667 :   const auto & lm_mesh = _mesh.getMesh();
     357             : 
     358           0 :   auto action_functor = [this, &lm_mesh](const processor_id_type libmesh_dbg_var(pid),
     359             :                                          const std::vector<Datum> & sent_data)
     360             :   {
     361             :     mooseAssert(pid != this->processor_id(), "We do not send messages to ourself here");
     362           0 :     for (auto & pr : sent_data)
     363             :     {
     364           0 :       const auto dof_id = pr.first;
     365             :       const auto * const dof_object =
     366           0 :           _nodal ? static_cast<const DofObject *>(lm_mesh.node_ptr(dof_id))
     367           0 :                  : static_cast<const DofObject *>(lm_mesh.elem_ptr(dof_id));
     368             :       mooseAssert(dof_object, "This should be non-null");
     369           0 :       _dof_to_nodal_wear_depth[dof_object] += std::move(pr.second);
     370             :     }
     371        1667 :   };
     372             : 
     373        1667 :   TIMPI::push_parallel_vector_data(_communicator, push_data, action_functor);
     374        1667 : }
     375             : 
     376             : void
     377       81816 : ComputeDynamicWeightedGapLMMechanicalContact::enforceConstraintOnDof(const DofObject * const dof)
     378             : {
     379       81816 :   const auto & weighted_gap = *_weighted_gap_ptr;
     380       81816 :   const Real c = _normalize_c ? _c / *_normalization_ptr : _c;
     381             : 
     382       81816 :   const auto dof_index = dof->dof_number(_sys.number(), _var->number(), 0);
     383       81816 :   ADReal lm_value = (*_sys.currentSolution())(dof_index);
     384             :   Moose::derivInsert(lm_value.derivatives(), dof_index, 1.);
     385             : 
     386       81816 :   const ADReal dof_residual = std::min(lm_value, weighted_gap * c);
     387             : 
     388       81816 :   addResidualsAndJacobian(_assembly,
     389      163632 :                           std::array<ADReal, 1>{{dof_residual}},
     390       81816 :                           std::array<dof_id_type, 1>{{dof_index}},
     391       81816 :                           _var->scalingFactor());
     392       81816 : }
     393             : 
     394             : void
     395      744040 : ComputeDynamicWeightedGapLMMechanicalContact::computeResidual(const Moose::MortarType mortar_type)
     396             : {
     397      744040 :   if (mortar_type != Moose::MortarType::Lower)
     398             :     return;
     399             : 
     400             :   mooseAssert(_var, "LM variable is null");
     401             : 
     402     1064680 :   for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
     403             :   {
     404      801520 :     computeQpProperties();
     405     3505360 :     for (_i = 0; _i < _test.size(); ++_i)
     406     2703840 :       computeQpIProperties();
     407             :   }
     408             : }
     409             : 
     410             : void
     411        3083 : ComputeDynamicWeightedGapLMMechanicalContact::jacobianSetup()
     412             : {
     413        3083 :   residualSetup();
     414        3083 : }
     415             : 
     416             : void
     417      221024 : ComputeDynamicWeightedGapLMMechanicalContact::computeJacobian(const Moose::MortarType mortar_type)
     418             : {
     419             :   // During "computeResidual" and "computeJacobian" we are actually just computing properties on the
     420             :   // mortar segment element mesh. We are *not* actually assembling into the residual/Jacobian. For
     421             :   // the zero-penetration constraint, the property of interest is the map from node to weighted gap.
     422             :   // Computation of the properties proceeds identically for residual and Jacobian evaluation hence
     423             :   // why we simply call computeResidual here. We will assemble into the residual/Jacobian later from
     424             :   // the post() method
     425      221024 :   computeResidual(mortar_type);
     426      221024 : }

Generated by: LCOV version 1.14