A GluedContactConstraint forces the value of a variable to be the same on both sides of an interface.
More...
#include <GluedContactConstraint.h>
A GluedContactConstraint forces the value of a variable to be the same on both sides of an interface.
Definition at line 27 of file GluedContactConstraint.h.
◆ GluedContactConstraint()
GluedContactConstraint::GluedContactConstraint |
( |
const InputParameters & |
parameters | ) |
|
Definition at line 62 of file GluedContactConstraint.C.
64 _component(getParam<unsigned int>(
"component")),
65 _model(getParam<MooseEnum>(
"model").getEnum<ContactModel>()),
66 _formulation(getParam<MooseEnum>(
"formulation").getEnum<ContactFormulation>()),
72 _vars(3, libMesh::invalid_uint),
77 if (parameters.isParamValid(
"tangential_tolerance"))
78 _penetration_locator.setTangentialTolerance(getParam<Real>(
"tangential_tolerance"));
80 if (parameters.isParamValid(
"normal_smoothing_distance"))
81 _penetration_locator.setNormalSmoothingDistance(getParam<Real>(
"normal_smoothing_distance"));
83 if (parameters.isParamValid(
"normal_smoothing_method"))
84 _penetration_locator.setNormalSmoothingMethod(
85 parameters.get<MooseEnum>(
"normal_smoothing_method"));
87 if (isParamValid(
"displacements"))
90 for (
unsigned int i = 0; i < coupledComponents(
"displacements"); ++i)
91 _vars[i] = coupled(
"displacements", i);
96 if (isParamValid(
"disp_x"))
97 _vars[0] = coupled(
"disp_x");
98 if (isParamValid(
"disp_y"))
99 _vars[1] = coupled(
"disp_y");
100 if (isParamValid(
"disp_z"))
101 _vars[2] = coupled(
"disp_z");
103 mooseDeprecated(
"use the `displacements` parameter rather than the `disp_*` parameters (those "
104 "will go away with the deprecation of the Solid Mechanics module).");
107 _penetration_locator.setUpdate(
false);
◆ ~GluedContactConstraint()
virtual GluedContactConstraint::~GluedContactConstraint |
( |
| ) |
|
|
inlinevirtual |
◆ computeQpJacobian()
Real GluedContactConstraint::computeQpJacobian |
( |
Moose::ConstraintJacobianType |
type | ) |
|
|
virtual |
Reimplemented from SparsityBasedContactConstraint.
Definition at line 210 of file GluedContactConstraint.C.
214 case Moose::SlaveSlave:
215 return _penalty * _phi_slave[_j][_qp] * _test_slave[_i][_qp];
217 case Moose::SlaveMaster:
218 return _penalty * -_phi_master[_j][_qp] * _test_slave[_i][_qp];
220 case Moose::MasterSlave:
222 const double slave_jac = (*_jacobian)(_current_node->dof_number(0,
_vars[
_component], 0),
223 _connected_dof_indices[_j]);
224 return slave_jac * _test_master[_i][_qp];
227 case Moose::MasterMaster:
231 mooseError(
"Unhandled ConstraintJacobianType");
◆ computeQpOffDiagJacobian()
Real GluedContactConstraint::computeQpOffDiagJacobian |
( |
Moose::ConstraintJacobianType |
type, |
|
|
unsigned int |
jvar |
|
) |
| |
|
virtual |
Compute off-diagonal Jacobian entries.
- Parameters
-
type | The type of coupling |
jvar | The index of the coupled variable |
Definition at line 238 of file GluedContactConstraint.C.
245 case Moose::SlaveSlave:
246 case Moose::SlaveMaster:
247 case Moose::MasterMaster:
250 case Moose::MasterSlave:
252 const double slave_jac = (*_jacobian)(_current_node->dof_number(0,
_vars[
_component], 0),
253 _connected_dof_indices[_j]);
254 retVal = slave_jac * _test_master[_i][_qp];
259 mooseError(
"Unhandled ConstraintJacobianType");
◆ computeQpResidual()
Real GluedContactConstraint::computeQpResidual |
( |
Moose::ConstraintType |
type | ) |
|
|
virtual |
Reimplemented from SparsityBasedContactConstraint.
Definition at line 175 of file GluedContactConstraint.C.
181 PenetrationInfo * pinfo = _penetration_locator._penetration_info[_current_node->id()];
183 Real pen_force =
_penalty * distance;
185 Real resid = pen_force;
187 pinfo->_mech_status = PenetrationInfo::MS_STICKING;
189 return _test_slave[_i][_qp] * resid;
194 PenetrationInfo * pinfo = _penetration_locator._penetration_info[_current_node->id()];
196 dof_id_type dof_number = _current_node->dof_number(0,
_vars[
_component], 0);
200 pinfo->_mech_status = PenetrationInfo::MS_STICKING;
202 return _test_master[_i][_qp] * resid;
◆ computeQpSlaveValue()
Real GluedContactConstraint::computeQpSlaveValue |
( |
| ) |
|
|
virtual |
◆ getConnectedDofIndices()
void SparsityBasedContactConstraint::getConnectedDofIndices |
( |
| ) |
|
|
virtualinherited |
Gets the indices for all dofs conected to the constraint Get indices for all the slave Jacobian columns, not only those based on the mesh connectivity.
Definition at line 29 of file SparsityBasedContactConstraint.C.
31 #if defined(LIBMESH_HAVE_PETSC) && !PETSC_VERSION_LESS_THAN(3, 3, 0)
32 _connected_dof_indices.clear();
40 PetscMatrix<Number> * petsc_jacobian =
dynamic_cast<PetscMatrix<Number> *
>(_jacobian);
41 mooseAssert(petsc_jacobian,
"Expected a PETSc matrix");
42 Mat jac = petsc_jacobian->mat();
45 const PetscInt * cols;
46 ierr = MatGetRow(jac, static_cast<PetscInt>(_var.nodalDofIndex()), &ncols, &cols, PETSC_NULL);
47 CHKERRABORT(_communicator.get(), ierr);
51 libMesh::out <<
"_connected_dof_indices: adding " << ncols <<
" dofs from Jacobian row["
52 << _var.nodalDofIndex() <<
"] = [";
54 for (PetscInt i = 0; i < ncols; ++i)
58 libMesh::out << cols[i] <<
" ";
60 _connected_dof_indices.push_back(cols[i]);
64 libMesh::out <<
"]\n";
66 ierr = MatRestoreRow(jac, static_cast<PetscInt>(_var.nodalDofIndex()), &ncols, &cols, PETSC_NULL);
67 CHKERRABORT(_communicator.get(), ierr);
69 NodeFaceConstraint::getConnectedDofIndices();
◆ jacobianSetup()
void GluedContactConstraint::jacobianSetup |
( |
| ) |
|
|
virtual |
◆ shouldApply()
bool GluedContactConstraint::shouldApply |
( |
| ) |
|
Definition at line 162 of file GluedContactConstraint.C.
164 auto hpit = _penetration_locator._has_penetrated.find(_current_node->id());
165 return hpit != _penetration_locator._has_penetrated.end();
◆ timestepSetup()
void GluedContactConstraint::timestepSetup |
( |
| ) |
|
|
virtual |
◆ updateContactSet()
void GluedContactConstraint::updateContactSet |
( |
bool |
beginning_of_step = false | ) |
|
|
virtual |
Definition at line 133 of file GluedContactConstraint.C.
135 auto & has_penetrated = _penetration_locator._has_penetrated;
137 for (
auto & pinfo_pair : _penetration_locator._penetration_info)
139 PenetrationInfo * pinfo = pinfo_pair.second;
142 if (!pinfo || pinfo->_node->n_comp(_sys.number(),
_vars[
_component]) < 1)
145 const dof_id_type slave_node_num = pinfo_pair.first;
146 auto hpit = has_penetrated.find(slave_node_num);
148 if (beginning_of_step)
150 pinfo->_starting_elem = pinfo->_elem;
151 pinfo->_starting_side_num = pinfo->_side_num;
152 pinfo->_starting_closest_point_ref = pinfo->_closest_point_ref;
155 if (pinfo->_distance >= 0)
156 if (hpit == has_penetrated.end())
157 has_penetrated.insert(slave_node_num);
Referenced by jacobianSetup(), and timestepSetup().
◆ _aux_solution
const NumericVector<Number>* GluedContactConstraint::_aux_solution |
|
protected |
◆ _aux_system
SystemBase& GluedContactConstraint::_aux_system |
|
protected |
◆ _component
const unsigned int GluedContactConstraint::_component |
|
protected |
◆ _formulation
◆ _friction_coefficient
const Real GluedContactConstraint::_friction_coefficient |
|
protected |
◆ _model
◆ _nodal_area_var
MooseVariable* GluedContactConstraint::_nodal_area_var |
|
protected |
◆ _penalty
const Real GluedContactConstraint::_penalty |
|
protected |
◆ _residual_copy
NumericVector<Number>& GluedContactConstraint::_residual_copy |
|
protected |
◆ _tension_release
const Real GluedContactConstraint::_tension_release |
|
protected |
◆ _updateContactSet
bool GluedContactConstraint::_updateContactSet |
|
protected |
◆ _vars
std::vector<unsigned int> GluedContactConstraint::_vars |
|
protected |
◆ _x_var
unsigned int GluedContactConstraint::_x_var |
|
protected |
◆ _y_var
unsigned int GluedContactConstraint::_y_var |
|
protected |
◆ _z_var
unsigned int GluedContactConstraint::_z_var |
|
protected |
The documentation for this class was generated from the following files: