A MultiDContactConstraint forces the value of a variable to be the same on both sides of an interface.
More...
#include <MultiDContactConstraint.h>
A MultiDContactConstraint forces the value of a variable to be the same on both sides of an interface.
Definition at line 27 of file MultiDContactConstraint.h.
◆ MultiDContactConstraint()
MultiDContactConstraint::MultiDContactConstraint |
( |
const InputParameters & |
parameters | ) |
|
Definition at line 58 of file MultiDContactConstraint.C.
59 : NodeFaceConstraint(parameters),
62 _component(getParam<unsigned int>(
"component")),
66 _vars(3, libMesh::invalid_uint)
68 _overwrite_slave_residual =
false;
70 if (isParamValid(
"displacements"))
73 for (
unsigned int i = 0; i < coupledComponents(
"displacements"); ++i)
74 _vars[i] = coupled(
"displacements", i);
79 if (isParamValid(
"disp_x"))
80 _vars[0] = coupled(
"disp_x");
81 if (isParamValid(
"disp_y"))
82 _vars[1] = coupled(
"disp_y");
83 if (isParamValid(
"disp_z"))
84 _vars[2] = coupled(
"disp_z");
86 mooseDeprecated(
"use the `displacements` parameter rather than the `disp_*` parameters (those "
87 "will go away with the deprecation of the Solid Mechanics module).");
◆ ~MultiDContactConstraint()
virtual MultiDContactConstraint::~MultiDContactConstraint |
( |
| ) |
|
|
inlinevirtual |
◆ computeQpJacobian()
Real MultiDContactConstraint::computeQpJacobian |
( |
Moose::ConstraintJacobianType |
type | ) |
|
|
virtual |
Definition at line 264 of file MultiDContactConstraint.C.
266 PenetrationInfo * pinfo = _penetration_locator._penetration_info[_current_node->id()];
268 Real slave_jac = 0.0;
271 case Moose::SlaveSlave:
278 (*_jacobian)(_current_node->dof_number(0, _var.number(), 0),
279 _connected_dof_indices[_j]));
287 mooseError(
"Invalid or unavailable contact model");
290 return _test_slave[_i][_qp] * slave_jac;
292 case Moose::SlaveMaster:
310 mooseError(
"Invalid or unavailable contact model");
313 return _test_slave[_i][_qp] * slave_jac;
315 case Moose::MasterSlave:
317 (*_jacobian)(_current_node->dof_number(0, _var.number(), 0), _connected_dof_indices[_j]);
318 return slave_jac * _test_master[_i][_qp];
320 case Moose::MasterMaster:
324 mooseError(
"Unhandled ConstraintJacobianType");
◆ computeQpResidual()
Real MultiDContactConstraint::computeQpResidual |
( |
Moose::ConstraintType |
type | ) |
|
|
virtual |
Definition at line 208 of file MultiDContactConstraint.C.
210 PenetrationInfo * pinfo = _penetration_locator._penetration_info[_current_node->id()];
211 const Node * node = pinfo->_node;
213 RealVectorValue res_vec;
217 dof_id_type dof_number = node->dof_number(0,
_vars[i], 0);
221 const RealVectorValue distance_vec(_mesh.nodeRef(node->id()) - pinfo->_closest_point);
222 const RealVectorValue pen_force(
_penalty * distance_vec);
231 resid = pinfo->_normal(
_component) * (pinfo->_normal * (pen_force - res_vec));
239 mooseError(
"Invalid or unavailable contact model");
242 return _test_slave[_i][_qp] * resid;
247 resid = pinfo->_normal(
_component) * (pinfo->_normal * res_vec);
255 mooseError(
"Invalid or unavailable contact model");
258 return _test_master[_i][_qp] * resid;
◆ computeQpSlaveValue()
Real MultiDContactConstraint::computeQpSlaveValue |
( |
| ) |
|
|
virtual |
Definition at line 184 of file MultiDContactConstraint.C.
186 PenetrationInfo * pinfo = _penetration_locator._penetration_info[_current_node->id()];
◆ jacobianSetup()
void MultiDContactConstraint::jacobianSetup |
( |
| ) |
|
|
virtual |
◆ shouldApply()
bool MultiDContactConstraint::shouldApply |
( |
| ) |
|
Definition at line 176 of file MultiDContactConstraint.C.
178 std::set<dof_id_type>::iterator hpit =
179 _penetration_locator._has_penetrated.find(_current_node->id());
180 return (hpit != _penetration_locator._has_penetrated.end());
◆ timestepSetup()
void MultiDContactConstraint::timestepSetup |
( |
| ) |
|
|
virtual |
◆ updateContactSet()
void MultiDContactConstraint::updateContactSet |
( |
| ) |
|
|
virtual |
Definition at line 106 of file MultiDContactConstraint.C.
108 auto & has_penetrated = _penetration_locator._has_penetrated;
110 for (
auto & pinfo_pair : _penetration_locator._penetration_info)
112 PenetrationInfo * pinfo = pinfo_pair.second;
115 if (!pinfo || pinfo->_node->n_comp(_sys.number(),
_vars[
_component]) < 1)
118 const Node * node = pinfo->_node;
120 dof_id_type slave_node_num = pinfo_pair.first;
121 auto hpit = has_penetrated.find(slave_node_num);
123 RealVectorValue res_vec;
127 dof_id_type dof_number = node->dof_number(0,
_vars[i], 0);
143 mooseError(
"Invalid or unavailable contact model");
162 if (pinfo->_distance > 0 &&
163 hpit == has_penetrated.end())
170 has_penetrated.insert(slave_node_num);
Referenced by jacobianSetup(), and timestepSetup().
◆ _component
const unsigned int MultiDContactConstraint::_component |
|
protected |
◆ _jacobian_update
bool MultiDContactConstraint::_jacobian_update |
|
protected |
◆ _mesh_dimension
const unsigned int MultiDContactConstraint::_mesh_dimension |
|
protected |
◆ _model
◆ _penalty
Real MultiDContactConstraint::_penalty |
|
protected |
◆ _residual_copy
NumericVector<Number>& MultiDContactConstraint::_residual_copy |
|
protected |
◆ _vars
std::vector<unsigned int> MultiDContactConstraint::_vars |
|
protected |
◆ _x_var
unsigned int MultiDContactConstraint::_x_var |
|
protected |
◆ _y_var
unsigned int MultiDContactConstraint::_y_var |
|
protected |
◆ _z_var
unsigned int MultiDContactConstraint::_z_var |
|
protected |
The documentation for this class was generated from the following files: