https://mooseframework.inl.gov
WeightedGapVelAux.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 
10 #include "WeightedGapVelAux.h"
11 #include "MooseVariableFE.h"
12 #include "FEProblemBase.h"
13 #include "Assembly.h"
14 
16 
19 {
21  params.addClassDescription(
22  "Returns the weighted gap velocity at a node. This quantity is useful for mortar contact, "
23  "particularly when dual basis functions are used in contact mechanics");
24  params.addCoupledVar("v",
25  "Optional variable to take the value of. If omitted the value of the "
26  "`variable` itself is returned.");
27  params.addRequiredCoupledVar("disp_x", "The x displacement variable");
28  params.addRequiredCoupledVar("disp_y", "The y displacement variable");
29  params.addCoupledVar("disp_z", "The z displacement variable");
30  params.set<bool>("interpolate_normals") = false;
31  params.addParam<bool>("use_displaced_mesh",
32  true,
33  "Whether to use the displaced mesh to compute the auxiliary kernel value.");
34  return params;
35 }
36 
38  : MortarNodalAuxKernel(parameters),
39  _has_disp_z(isCoupled("disp_z")),
40  _disp_x(*getVar("disp_x", 0)),
41  _disp_y(*getVar("disp_y", 0)),
42  _disp_z(getVar("disp_z", 0)),
43  _secondary_x_dot(_disp_x.adUDot()),
44  _primary_x_dot(_disp_x.adUDotNeighbor()),
45  _secondary_y_dot(_disp_y.adUDot()),
46  _primary_y_dot(_disp_y.adUDotNeighbor()),
47  _secondary_z_dot(_has_disp_z ? &_disp_z->adUDot() : nullptr),
48  _primary_z_dot(_has_disp_z ? &_disp_z->adUDotNeighbor() : nullptr),
49  _weighted_gap_velocity(0),
50  _qp_gap_velocity(0),
51  _qp_gap_velocity_nodal(0)
52 {
53  if (!_displaced)
55  "use_displaced_mesh",
56  "This auxiliary kernel typically requires the use of displaced meshes to compute the "
57  "weighted gap velocity.");
58 }
59 
60 Real
62 {
64  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
65  {
67  for (_i = 0; _i < _test_lower.size(); ++_i)
69  }
70 
72 }
73 
74 void
76 {
77  RealVectorValue gap_velocity_vec;
80 
81  if (_has_disp_z)
82  gap_velocity_vec(2) = MetaPhysicL::raw_value((*_secondary_z_dot)[_qp] - (*_primary_z_dot)[_qp]);
83 
84  _qp_gap_velocity_nodal = gap_velocity_vec * (_JxW_msm[_qp] * _coord_msm[_qp]);
85 
87 }
88 
89 void
91 {
93 }
registerMooseObject("ContactApp", WeightedGapVelAux)
static InputParameters validParams()
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
RealVectorValue _qp_gap_velocity_nodal
The gap velocity vector at the current quadrature point, used when we are not interpolating the norma...
static InputParameters validParams()
T & set(const std::string &name, bool quiet_mode=false)
auto raw_value(const Eigen::Map< T > &in)
const ADVariableValue *const _secondary_z_dot
z-velocity on the secondary face
const bool _has_disp_z
For 2D mortar contact no displacement will be specified, so const pointers used.
const std::vector< Real > & _JxW_msm
const libMesh::QBase *const & _qrule_msm
const OutputTools< ComputeValueType >::VariableTestValue & _test_lower
const ADVariableValue & _primary_x_dot
x-velocity on the primary face
const ADVariableValue *const _primary_z_dot
z-velocity on the primary face
std::vector< Point > _normals
const MooseArray< Real > & _coord_msm
unsigned int n_points() const
Real _weighted_gap_velocity
The weighted gap velocity.
WeightedGapVelAux(const InputParameters &parameters)
Factory constructor, takes parameters so that all derived classes can be built using the same constru...
void addCoupledVar(const std::string &name, const std::string &doc_string)
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const ADVariableValue & _secondary_x_dot
x-velocity on the secondary face
void addClassDescription(const std::string &doc_string)
const ADVariableValue & _primary_y_dot
y-velocity on the primary face
void paramWarning(const std::string &param, Args... args) const
Compute nodal weighted gap velocity based on a mortar discretization.
Real computeValue() override
const ADVariableValue & _secondary_y_dot
y-velocity on the secondary face