v0.14.0
Public Member Functions | Private Attributes | List of all members
SimpleContactProblem::OpCalIntCompFunSlave Struct Reference

RHS-operator for the simple contact element. More...

#include <users_modules/mortar_contact/src/SimpleContact.hpp>

Inheritance diagram for SimpleContactProblem::OpCalIntCompFunSlave:
[legend]
Collaboration diagram for SimpleContactProblem::OpCalIntCompFunSlave:
[legend]

Public Member Functions

 OpCalIntCompFunSlave (const string lagrange_field_name, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, boost::shared_ptr< double > cn)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 Integrates the complementarity function at slave face gauss points and assembles components to the RHS global vector. More...
 

Private Attributes

boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 
boost::shared_ptr< doublecNPtr
 
VectorDouble vecR
 

Detailed Description

RHS-operator for the simple contact element.

Integrates complementarity function that fulfills KKT conditions over slave contact area and assembles components of the RHS vector.

Definition at line 1160 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpCalIntCompFunSlave()

SimpleContactProblem::OpCalIntCompFunSlave::OpCalIntCompFunSlave ( const string  lagrange_field_name,
boost::shared_ptr< CommonDataSimpleContact common_data_contact,
boost::shared_ptr< double cn 
)
inline

Definition at line 1162 of file SimpleContact.hpp.

1166  : ContactOp(lagrange_field_name, UserDataOperator::OPCOL,
1167  ContactOp::FACESLAVE),
1168  commonDataSimpleContact(common_data_contact), cNPtr(cn) {}

Member Function Documentation

◆ doWork()

MoFEMErrorCode SimpleContactProblem::OpCalIntCompFunSlave::doWork ( int  side,
EntityType  type,
EntData data 
)

Integrates the complementarity function at slave face gauss points and assembles components to the RHS global vector.

Integrates the complementarity function to fulfil KKT conditions in the integral sense and assembles components to the RHS global vector

\[ {\overline C(\lambda, \mathbf{x}^{(i)}, \delta \lambda)} = \int_{{\gamma}^{(1)}_{\text c}} \left( \lambda + c_{\text n} g_{\textrm{n}} - \dfrac{1}{r}{\left| \lambda - c_{\text n} g_{\textrm{n}}\right|}^{r}\right) \delta{{\lambda}} \,\,{ {\text d} {\gamma}} \]

where \({\gamma}^{(1)}_{\text c}\) is the surface integration domain of the slave surface, \( \lambda\) is the Lagrange multiplier, \(\mathbf{x}^{(i)}\) are the coordinates of the overlapping gauss points at slave and master triangles for \(i = 1\) and \(i = 2\), respectively. Furthermore, \( c_{\text n}\) works as an augmentation parameter and affects convergence, \(r\) is regularisation parameter that can be chosen in \([1, 1.1]\) ( \(r = 1\)) is the default value) and \( g_{\textrm{n}}\) is the gap function evaluated at the slave triangle gauss points as:

\[ g_{\textrm{n}} = - \mathbf{n}(\mathbf{x}^{(1)}) \cdot \left( \mathbf{x}^{(1)} - \mathbf{x}^{(2)} \right) \]

Definition at line 1643 of file SimpleContact.cpp.

1644  {
1646 
1647  if (data.getIndices().size() == 0)
1649 
1650  const int nb_gauss_pts = data.getN().size1();
1651 
1652  int nb_base_fun_col = data.getFieldData().size();
1653  const double area_s =
1654  commonDataSimpleContact->areaSlave; // same area in master and slave
1655 
1656  vecR.resize(nb_base_fun_col, false);
1657  vecR.clear();
1658 
1659  auto t_lagrange_slave =
1660  getFTensor0FromVec(*commonDataSimpleContact->lagMultAtGaussPtsPtr);
1661  auto t_gap_gp = getFTensor0FromVec(*commonDataSimpleContact->gapPtr);
1662  auto t_w = getFTensor0IntegrationWeightSlave();
1663  const double cn_value = *cNPtr.get();
1664  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1665  const double val_s = t_w * area_s *
1667  cn_value, t_gap_gp, t_lagrange_slave);
1668  auto t_base_lambda = data.getFTensor0N(gg, 0);
1669  for (int bbr = 0; bbr != nb_base_fun_col; ++bbr) {
1670  vecR[bbr] += val_s * t_base_lambda;
1671 
1672  ++t_base_lambda; // update rows
1673  }
1674 
1675  ++t_lagrange_slave;
1676  ++t_gap_gp;
1677  ++t_w;
1678  } // for gauss points
1679 
1680  CHKERR VecSetValues(getSNESf(), data, &vecR[0], ADD_VALUES);
1681 
1683 }

Member Data Documentation

◆ cNPtr

boost::shared_ptr<double> SimpleContactProblem::OpCalIntCompFunSlave::cNPtr
private

Definition at line 1201 of file SimpleContact.hpp.

◆ commonDataSimpleContact

boost::shared_ptr<CommonDataSimpleContact> SimpleContactProblem::OpCalIntCompFunSlave::commonDataSimpleContact
private

Definition at line 1200 of file SimpleContact.hpp.

◆ vecR

VectorDouble SimpleContactProblem::OpCalIntCompFunSlave::vecR
private

Definition at line 1202 of file SimpleContact.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
SimpleContactProblem::OpCalIntCompFunSlave::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:1200
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1254
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1589
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1502
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
SimpleContactProblem::OpCalIntCompFunSlave::vecR
VectorDouble vecR
Definition: SimpleContact.hpp:1202
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPCOL
@ OPCOL
operator doWork function is executed on FE columns
Definition: ForcesAndSourcesCore.hpp:568
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1214
SimpleContactProblem::OpCalIntCompFunSlave::cNPtr
boost::shared_ptr< double > cNPtr
Definition: SimpleContact.hpp:1201
SimpleContactProblem::ContactOp
ContactPrismElementForcesAndSourcesCore::UserDataOperator ContactOp
Definition: SimpleContact.hpp:30
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1318
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
SimpleContactProblem::ConstrainFunction
static double ConstrainFunction(const double cn, const double g, const double l)
Definition: SimpleContact.hpp:3587
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359