www.mooseframework.org
Public Member Functions | Protected Attributes | List of all members
GluedContactConstraint Class Reference

A GluedContactConstraint forces the value of a variable to be the same on both sides of an interface. More...

#include <GluedContactConstraint.h>

Inheritance diagram for GluedContactConstraint:
[legend]

Public Member Functions

 GluedContactConstraint (const InputParameters &parameters)
 
virtual ~GluedContactConstraint ()
 
virtual void timestepSetup ()
 
virtual void jacobianSetup ()
 
virtual void updateContactSet (bool beginning_of_step=false)
 
virtual Real computeQpSlaveValue ()
 
virtual Real computeQpResidual (Moose::ConstraintType type)
 
virtual Real computeQpJacobian (Moose::ConstraintJacobianType type)
 
virtual Real computeQpOffDiagJacobian (Moose::ConstraintJacobianType type, unsigned int jvar)
 Compute off-diagonal Jacobian entries. More...
 
bool shouldApply ()
 
virtual void getConnectedDofIndices ()
 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. More...
 

Protected Attributes

const unsigned int _component
 
const ContactModel _model
 
const ContactFormulation _formulation
 
const Real _penalty
 
const Real _friction_coefficient
 
const Real _tension_release
 
bool _updateContactSet
 
NumericVector< Number > & _residual_copy
 
unsigned int _x_var
 
unsigned int _y_var
 
unsigned int _z_var
 
std::vector< unsigned int > _vars
 
MooseVariable * _nodal_area_var
 
SystemBase & _aux_system
 
const NumericVector< Number > * _aux_solution
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ GluedContactConstraint()

GluedContactConstraint::GluedContactConstraint ( const InputParameters &  parameters)

Definition at line 62 of file GluedContactConstraint.C.

63  : SparsityBasedContactConstraint(parameters),
64  _component(getParam<unsigned int>("component")),
65  _model(getParam<MooseEnum>("model").getEnum<ContactModel>()),
66  _formulation(getParam<MooseEnum>("formulation").getEnum<ContactFormulation>()),
67  _penalty(getParam<Real>("penalty")),
68  _friction_coefficient(getParam<Real>("friction_coefficient")),
69  _tension_release(getParam<Real>("tension_release")),
70  _updateContactSet(true),
71  _residual_copy(_sys.residualGhosted()),
72  _vars(3, libMesh::invalid_uint),
73  _nodal_area_var(getVar("nodal_area", 0)),
75  _aux_solution(_aux_system.currentSolution())
76 {
77  if (parameters.isParamValid("tangential_tolerance"))
78  _penetration_locator.setTangentialTolerance(getParam<Real>("tangential_tolerance"));
79 
80  if (parameters.isParamValid("normal_smoothing_distance"))
81  _penetration_locator.setNormalSmoothingDistance(getParam<Real>("normal_smoothing_distance"));
82 
83  if (parameters.isParamValid("normal_smoothing_method"))
84  _penetration_locator.setNormalSmoothingMethod(
85  parameters.get<MooseEnum>("normal_smoothing_method"));
86 
87  if (isParamValid("displacements"))
88  {
89  // modern parameter scheme for displacements
90  for (unsigned int i = 0; i < coupledComponents("displacements"); ++i)
91  _vars[i] = coupled("displacements", i);
92  }
93  else
94  {
95  // Legacy parameter scheme for displacements
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");
102 
103  mooseDeprecated("use the `displacements` parameter rather than the `disp_*` parameters (those "
104  "will go away with the deprecation of the Solid Mechanics module).");
105  }
106 
107  _penetration_locator.setUpdate(false);
108 }

◆ ~GluedContactConstraint()

virtual GluedContactConstraint::~GluedContactConstraint ( )
inlinevirtual

Definition at line 31 of file GluedContactConstraint.h.

31 {}

Member Function Documentation

◆ computeQpJacobian()

Real GluedContactConstraint::computeQpJacobian ( Moose::ConstraintJacobianType  type)
virtual

Reimplemented from SparsityBasedContactConstraint.

Definition at line 210 of file GluedContactConstraint.C.

211 {
212  switch (type)
213  {
214  case Moose::SlaveSlave:
215  return _penalty * _phi_slave[_j][_qp] * _test_slave[_i][_qp];
216 
217  case Moose::SlaveMaster:
218  return _penalty * -_phi_master[_j][_qp] * _test_slave[_i][_qp];
219 
220  case Moose::MasterSlave:
221  {
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];
225  }
226 
227  case Moose::MasterMaster:
228  return 0.0;
229 
230  default:
231  mooseError("Unhandled ConstraintJacobianType");
232  }
233 
234  return 0.0;
235 }

◆ computeQpOffDiagJacobian()

Real GluedContactConstraint::computeQpOffDiagJacobian ( Moose::ConstraintJacobianType  type,
unsigned int  jvar 
)
virtual

Compute off-diagonal Jacobian entries.

Parameters
typeThe type of coupling
jvarThe index of the coupled variable

Definition at line 238 of file GluedContactConstraint.C.

240 {
241  Real retVal = 0.0;
242 
243  switch (type)
244  {
245  case Moose::SlaveSlave:
246  case Moose::SlaveMaster:
247  case Moose::MasterMaster:
248  break;
249 
250  case Moose::MasterSlave:
251  {
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];
255  break;
256  }
257 
258  default:
259  mooseError("Unhandled ConstraintJacobianType");
260  }
261 
262  return retVal;
263 }

◆ computeQpResidual()

Real GluedContactConstraint::computeQpResidual ( Moose::ConstraintType  type)
virtual

Reimplemented from SparsityBasedContactConstraint.

Definition at line 175 of file GluedContactConstraint.C.

176 {
177  switch (type)
178  {
179  case Moose::Slave:
180  {
181  PenetrationInfo * pinfo = _penetration_locator._penetration_info[_current_node->id()];
182  Real distance = (*_current_node)(_component)-pinfo->_closest_point(_component);
183  Real pen_force = _penalty * distance;
184 
185  Real resid = pen_force;
186  pinfo->_contact_force(_component) = resid;
187  pinfo->_mech_status = PenetrationInfo::MS_STICKING;
188 
189  return _test_slave[_i][_qp] * resid;
190  }
191 
192  case Moose::Master:
193  {
194  PenetrationInfo * pinfo = _penetration_locator._penetration_info[_current_node->id()];
195 
196  dof_id_type dof_number = _current_node->dof_number(0, _vars[_component], 0);
197  Real resid = _residual_copy(dof_number);
198 
199  pinfo->_contact_force(_component) = -resid;
200  pinfo->_mech_status = PenetrationInfo::MS_STICKING;
201 
202  return _test_master[_i][_qp] * resid;
203  }
204  }
205 
206  return 0.0;
207 }

◆ computeQpSlaveValue()

Real GluedContactConstraint::computeQpSlaveValue ( )
virtual

Reimplemented from SparsityBasedContactConstraint.

Definition at line 169 of file GluedContactConstraint.C.

170 {
171  return _u_slave[_qp];
172 }

◆ 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.

30 {
31 #if defined(LIBMESH_HAVE_PETSC) && !PETSC_VERSION_LESS_THAN(3, 3, 0)
32  _connected_dof_indices.clear();
33 
34  // An ugly hack: have to extract the row and look at it's sparsity structure, since otherwise I
35  // won't get the dofs connected to this row by virtue of intervariable coupling
36  // This is easier than sifting through all coupled variables, selecting those active on the
37  // current subdomain, dealing with the scalar variables, etc.
38  // Also, importantly, this will miss coupling to variables that might have introduced by prior
39  // constraints similar to this one!
40  PetscMatrix<Number> * petsc_jacobian = dynamic_cast<PetscMatrix<Number> *>(_jacobian);
41  mooseAssert(petsc_jacobian, "Expected a PETSc matrix");
42  Mat jac = petsc_jacobian->mat();
43  PetscErrorCode ierr;
44  PetscInt ncols;
45  const PetscInt * cols;
46  ierr = MatGetRow(jac, static_cast<PetscInt>(_var.nodalDofIndex()), &ncols, &cols, PETSC_NULL);
47  CHKERRABORT(_communicator.get(), ierr);
48  bool debug = false;
49  if (debug)
50  {
51  libMesh::out << "_connected_dof_indices: adding " << ncols << " dofs from Jacobian row["
52  << _var.nodalDofIndex() << "] = [";
53  }
54  for (PetscInt i = 0; i < ncols; ++i)
55  {
56  if (debug)
57  {
58  libMesh::out << cols[i] << " ";
59  }
60  _connected_dof_indices.push_back(cols[i]);
61  }
62  if (debug)
63  {
64  libMesh::out << "]\n";
65  }
66  ierr = MatRestoreRow(jac, static_cast<PetscInt>(_var.nodalDofIndex()), &ncols, &cols, PETSC_NULL);
67  CHKERRABORT(_communicator.get(), ierr);
68 #else
69  NodeFaceConstraint::getConnectedDofIndices();
70 #endif
71 }

◆ jacobianSetup()

void GluedContactConstraint::jacobianSetup ( )
virtual

Definition at line 121 of file GluedContactConstraint.C.

122 {
123  if (_component == 0)
124  {
125  if (_updateContactSet)
127 
128  _updateContactSet = true;
129  }
130 }

◆ shouldApply()

bool GluedContactConstraint::shouldApply ( )

Definition at line 162 of file GluedContactConstraint.C.

163 {
164  auto hpit = _penetration_locator._has_penetrated.find(_current_node->id());
165  return hpit != _penetration_locator._has_penetrated.end();
166 }

◆ timestepSetup()

void GluedContactConstraint::timestepSetup ( )
virtual

Definition at line 111 of file GluedContactConstraint.C.

112 {
113  if (_component == 0)
114  {
115  updateContactSet(true);
116  _updateContactSet = false;
117  }
118 }

◆ updateContactSet()

void GluedContactConstraint::updateContactSet ( bool  beginning_of_step = false)
virtual

Definition at line 133 of file GluedContactConstraint.C.

134 {
135  auto & has_penetrated = _penetration_locator._has_penetrated;
136 
137  for (auto & pinfo_pair : _penetration_locator._penetration_info)
138  {
139  PenetrationInfo * pinfo = pinfo_pair.second;
140 
141  // Skip this pinfo if there are no DOFs on this node.
142  if (!pinfo || pinfo->_node->n_comp(_sys.number(), _vars[_component]) < 1)
143  continue;
144 
145  const dof_id_type slave_node_num = pinfo_pair.first;
146  auto hpit = has_penetrated.find(slave_node_num);
147 
148  if (beginning_of_step)
149  {
150  pinfo->_starting_elem = pinfo->_elem;
151  pinfo->_starting_side_num = pinfo->_side_num;
152  pinfo->_starting_closest_point_ref = pinfo->_closest_point_ref;
153  }
154 
155  if (pinfo->_distance >= 0)
156  if (hpit == has_penetrated.end())
157  has_penetrated.insert(slave_node_num);
158  }
159 }

Referenced by jacobianSetup(), and timestepSetup().

Member Data Documentation

◆ _aux_solution

const NumericVector<Number>* GluedContactConstraint::_aux_solution
protected

Definition at line 73 of file GluedContactConstraint.h.

◆ _aux_system

SystemBase& GluedContactConstraint::_aux_system
protected

Definition at line 72 of file GluedContactConstraint.h.

◆ _component

const unsigned int GluedContactConstraint::_component
protected

◆ _formulation

const ContactFormulation GluedContactConstraint::_formulation
protected

Definition at line 56 of file GluedContactConstraint.h.

◆ _friction_coefficient

const Real GluedContactConstraint::_friction_coefficient
protected

Definition at line 59 of file GluedContactConstraint.h.

◆ _model

const ContactModel GluedContactConstraint::_model
protected

Definition at line 55 of file GluedContactConstraint.h.

◆ _nodal_area_var

MooseVariable* GluedContactConstraint::_nodal_area_var
protected

Definition at line 71 of file GluedContactConstraint.h.

◆ _penalty

const Real GluedContactConstraint::_penalty
protected

Definition at line 58 of file GluedContactConstraint.h.

Referenced by computeQpJacobian(), and computeQpResidual().

◆ _residual_copy

NumericVector<Number>& GluedContactConstraint::_residual_copy
protected

Definition at line 63 of file GluedContactConstraint.h.

Referenced by computeQpResidual().

◆ _tension_release

const Real GluedContactConstraint::_tension_release
protected

Definition at line 60 of file GluedContactConstraint.h.

◆ _updateContactSet

bool GluedContactConstraint::_updateContactSet
protected

Definition at line 61 of file GluedContactConstraint.h.

Referenced by jacobianSetup(), and timestepSetup().

◆ _vars

std::vector<unsigned int> GluedContactConstraint::_vars
protected

◆ _x_var

unsigned int GluedContactConstraint::_x_var
protected

Definition at line 65 of file GluedContactConstraint.h.

◆ _y_var

unsigned int GluedContactConstraint::_y_var
protected

Definition at line 66 of file GluedContactConstraint.h.

◆ _z_var

unsigned int GluedContactConstraint::_z_var
protected

Definition at line 67 of file GluedContactConstraint.h.


The documentation for this class was generated from the following files:
GluedContactConstraint::_formulation
const ContactFormulation _formulation
Definition: GluedContactConstraint.h:56
GluedContactConstraint::_residual_copy
NumericVector< Number > & _residual_copy
Definition: GluedContactConstraint.h:63
GluedContactConstraint::_friction_coefficient
const Real _friction_coefficient
Definition: GluedContactConstraint.h:59
GluedContactConstraint::_tension_release
const Real _tension_release
Definition: GluedContactConstraint.h:60
GluedContactConstraint::_model
const ContactModel _model
Definition: GluedContactConstraint.h:55
GluedContactConstraint::_penalty
const Real _penalty
Definition: GluedContactConstraint.h:58
GluedContactConstraint::_aux_solution
const NumericVector< Number > * _aux_solution
Definition: GluedContactConstraint.h:73
GluedContactConstraint::_nodal_area_var
MooseVariable * _nodal_area_var
Definition: GluedContactConstraint.h:71
SparsityBasedContactConstraint::SparsityBasedContactConstraint
SparsityBasedContactConstraint(const InputParameters &parameters)
Definition: SparsityBasedContactConstraint.h:24
GluedContactConstraint::_component
const unsigned int _component
Definition: GluedContactConstraint.h:54
GluedContactConstraint::_vars
std::vector< unsigned int > _vars
Definition: GluedContactConstraint.h:69
GluedContactConstraint::_aux_system
SystemBase & _aux_system
Definition: GluedContactConstraint.h:72
GluedContactConstraint::updateContactSet
virtual void updateContactSet(bool beginning_of_step=false)
Definition: GluedContactConstraint.C:133
GluedContactConstraint::_updateContactSet
bool _updateContactSet
Definition: GluedContactConstraint.h:61