https://mooseframework.inl.gov
WeightedVelocitiesUserObject.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 
11 #include "MooseVariableField.h"
12 #include "SubProblem.h"
13 #include "SystemBase.h"
14 #include "MortarUtils.h"
15 #include "MooseUtils.h"
16 #include "MortarContactUtils.h"
17 #include "libmesh/quadrature.h"
18 
21 {
23  params.addRequiredParam<VariableName>("secondary_variable",
24  "Primal variable on secondary surface.");
25  params.addParam<VariableName>(
26  "primary_variable",
27  "Primal variable on primary surface. If this parameter is not provided then the primary "
28  "variable will be initialized to the secondary variable");
29  return params;
30 }
31 
33  : WeightedGapUserObject(parameters),
34  _sys(*getCheckedPointerParam<SystemBase *>("_sys")),
35  _secondary_var(
36  isParamValid("secondary_variable")
37  ? _sys.getActualFieldVariable<Real>(_tid, parameters.getMooseType("secondary_variable"))
38  : _sys.getActualFieldVariable<Real>(_tid, parameters.getMooseType("primary_variable"))),
39  _primary_var(
40  isParamValid("primary_variable")
41  ? _sys.getActualFieldVariable<Real>(_tid, parameters.getMooseType("primary_variable"))
42  : _secondary_var),
43  _secondary_x_dot(_secondary_var.adUDot()),
44  _primary_x_dot(_primary_var.adUDotNeighbor()),
45  _secondary_y_dot(adCoupledDot("disp_y")),
46  _primary_y_dot(adCoupledNeighborValueDot("disp_y")),
47  _secondary_z_dot(_has_disp_z ? &adCoupledDot("disp_z") : nullptr),
48  _primary_z_dot(_has_disp_z ? &adCoupledNeighborValueDot("disp_z") : nullptr),
49  _3d(_has_disp_z)
50 {
51  if (!getParam<bool>("use_displaced_mesh"))
52  paramError("use_displaced_mesh",
53  "'use_displaced_mesh' must be true for the WeightedVelocitiesUserObject object");
54 }
55 
56 void
58 {
59  // Compute the value of _qp_gap
61 
62  // Trim derivatives
63  const auto & primary_ip_lowerd_map = amg().getPrimaryIpToLowerElementMap(
65  const auto & secondary_ip_lowerd_map =
67 
68  std::array<const MooseVariable *, 3> var_array{{_disp_x_var, _disp_y_var, _disp_z_var}};
69  std::array<ADReal, 3> primary_disp_dot{
70  {_primary_x_dot[_qp], _primary_y_dot[_qp], _has_disp_z ? (*_primary_z_dot)[_qp] : 0}};
71  std::array<ADReal, 3> secondary_disp_dot{
72  {_secondary_x_dot[_qp], _secondary_y_dot[_qp], _has_disp_z ? (*_secondary_z_dot)[_qp] : 0}};
73 
74  trimInteriorNodeDerivatives(primary_ip_lowerd_map, var_array, primary_disp_dot, false);
75  trimInteriorNodeDerivatives(secondary_ip_lowerd_map, var_array, secondary_disp_dot, true);
76 
77  const ADReal & prim_x_dot = primary_disp_dot[0];
78  const ADReal & prim_y_dot = primary_disp_dot[1];
79  const ADReal * prim_z_dot = nullptr;
80  if (_has_disp_z)
81  prim_z_dot = &primary_disp_dot[2];
82 
83  const ADReal & sec_x_dot = secondary_disp_dot[0];
84  const ADReal & sec_y_dot = secondary_disp_dot[1];
85  const ADReal * sec_z_dot = nullptr;
86  if (_has_disp_z)
87  sec_z_dot = &secondary_disp_dot[2];
88 
89  // Build relative velocity vector
90  ADRealVectorValue relative_velocity;
91 
92  if (_3d)
93  relative_velocity = {sec_x_dot - prim_x_dot, sec_y_dot - prim_y_dot, *sec_z_dot - *prim_z_dot};
94  else
95  relative_velocity = {sec_x_dot - prim_x_dot, sec_y_dot - prim_y_dot, 0.0};
96 
97  // Geometry is averaged and used at the nodes for constraint enforcement.
98  _qp_real_tangential_velocity_nodal = relative_velocity;
99  _qp_tangential_velocity_nodal = relative_velocity * (_JxW_msm[_qp] * _coord[_qp]);
100 }
101 
102 void
104 {
106 
107  const auto & nodal_tangents = amg().getNodalTangents(*_lower_secondary_elem);
108  // Get the _dof_to_weighted_tangential_velocity map
109  const DofObject * const dof =
110  _is_weighted_gap_nodal ? static_cast<const DofObject *>(_lower_secondary_elem->node_ptr(_i))
111  : static_cast<const DofObject *>(_lower_secondary_elem);
112 
114  (*_test)[_i][_qp] * _qp_tangential_velocity_nodal * nodal_tangents[0][_i];
116  (*_test)[_i][_qp] * _qp_real_tangential_velocity_nodal * nodal_tangents[0][_i];
117 
118  // Get the _dof_to_weighted_tangential_velocity map for a second direction
119  if (_3d)
120  {
122  (*_test)[_i][_qp] * _qp_tangential_velocity_nodal * nodal_tangents[1][_i];
123 
125  (*_test)[_i][_qp] * _qp_real_tangential_velocity_nodal * nodal_tangents[1][_i];
126  }
127 }
128 
129 void
131 {
134 }
135 
136 void
138 {
139  // Clear weighted gaps
141  selfInitialize();
142 }
143 
144 void
146 {
148 
149  // If the constraint is performed by the owner, then we don't need any data sent back; the owner
150  // will take care of it. But if the constraint is not performed by the owner and we might have to
151  // do some of the constraining ourselves, then we need data sent back to us
152  const bool send_data_back = !constrainedByOwner();
153 
155  _subproblem.mesh(),
156  _nodal,
158  send_data_back);
161 }
virtual MooseMesh & mesh()=0
ADRealVectorValue _qp_tangential_velocity_nodal
The value of the tangential velocity vectors at the current node.
void communicateVelocities(std::unordered_map< const DofObject *, T > &dof_map, const MooseMesh &mesh, const bool nodal, const Parallel::Communicator &communicator, const bool send_data_back)
This function is used to communicate velocities across processes.
unsigned int _qp
Quadrature point index for the mortar segments.
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
const ADVariableValue & _secondary_y_dot
y-velocity on the secondary face
std::map< unsigned int, unsigned int > getPrimaryIpToLowerElementMap(const Elem &primary_elem, const Elem &primary_elem_ip, const Elem &lower_secondary_elem) const
std::unordered_map< const DofObject *, std::array< ADReal, 2 > > _dof_to_real_tangential_velocity
A map from node to two interpolated, physical tangential velocities.
virtual void initialize() override
const MooseVariable *const _disp_z_var
The z displacement variable.
bool _3d
Automatic flag to determine whether we are doing three-dimensional work.
static InputParameters validParams()
Elem const *const & _lower_primary_elem
const bool _nodal
Whether the dof objects are nodal; if they&#39;re not, then they&#39;re elemental.
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 Parallel::Communicator & _communicator
DualNumber< Real, DNDerivativeType, true > ADReal
void addRequiredParam(const std::string &name, const std::string &doc_string)
ADRealVectorValue _qp_real_tangential_velocity_nodal
The value of the real tangential velocity vectors at the current node.
SubProblem & _subproblem
Elem const *const & _lower_secondary_elem
const ADVariableValue & _primary_y_dot
y-velocity on the primary face
const AutomaticMortarGeneration & amg() const
static void trimInteriorNodeDerivatives(const std::map< unsigned int, unsigned int > &primary_ip_lowerd_map, const Variables &moose_var, DualNumbers &ad_vars, const bool is_secondary)
std::array< MooseUtils::SemidynamicVector< Point, 9 >, 2 > getNodalTangents(const Elem &secondary_elem) const
const ADVariableValue & _primary_x_dot
x-velocity on the primary face
virtual void computeQpProperties()
Computes properties that are functions only of the current quadrature point (_qp), e.g.
std::map< unsigned int, unsigned int > getSecondaryIpToLowerElementMap(const Elem &lower_secondary_elem) const
virtual void computeQpProperties() override
Computes properties that are functions only of the current quadrature point (_qp), e.g.
void paramError(const std::string &param, Args... args) const
const MooseArray< Real > & _coord
Member for handling change of coordinate systems (xyz, rz, spherical)
Creates dof object to weighted gap map.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool _is_weighted_gap_nodal
Whether the weighted gap is associated with nodes or elements (like for a CONSTANT MONOMIAL Lagrange ...
const ADVariableValue & _secondary_x_dot
x-velocity on the secondary face
const MooseVariable *const _disp_y_var
The y displacement variable.
WeightedVelocitiesUserObject(const InputParameters &parameters)
virtual bool constrainedByOwner() const =0
unsigned int _i
Test function index.
std::unordered_map< const DofObject *, std::array< ADReal, 2 > > _dof_to_weighted_tangential_velocity
A map from node to two weighted tangential velocities.
const MooseVariable *const _disp_x_var
The x displacement variable.
virtual void computeQpIProperties()
Computes properties that are functions both of _qp and _i, for example the weighted gap...
virtual void computeQpIProperties() override
Computes properties that are functions both of _qp and _i, for example the weighted gap...
virtual void finalize() override