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

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

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

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

Public Member Functions

 OpGapConstraintAugmentedOverLambda (const string lagrange_field_name, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, const double cn)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 Integrates the conditions that fulfil KKT conditions at slave face gauss points and assembles components to LHS global matrix. More...
 

Private Attributes

boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 
const double cN
 
MatrixDouble NN
 

Detailed Description

LHS-operator for the simple contact element.

Integrates variation of the conditions that fulfil KKT conditions with respect to Lagrange multipliers on slave side and assembles components to LHS global matrix.

Definition at line 1987 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpGapConstraintAugmentedOverLambda()

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

Definition at line 1989 of file SimpleContact.hpp.

1993  : ContactOp(lagrange_field_name, UserDataOperator::OPROWCOL,
1994  ContactOp::FACESLAVESLAVE),
1995  commonDataSimpleContact(common_data_contact), cN(cn) {
1996  sYmm = false; // This will make sure to loop over all entities (e.g.
1997  // for order=2 it will make doWork to loop 16 time)
1998  }

Member Function Documentation

◆ doWork()

MoFEMErrorCode SimpleContactProblem::OpGapConstraintAugmentedOverLambda::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntData row_data,
EntData col_data 
)

Integrates the conditions that fulfil KKT conditions at slave face gauss points and assembles components to LHS global matrix.

Integrates variation of the expresion that fulfils KKT conditions with respect to Lagrange multipliers and assembles components to LHS global matrix.

\[ {\text D}{\overline C(\lambda, \mathbf{x}^{(1)}, \delta \lambda)}[\Delta \lambda] = \left\{ \begin{array}{ll} 0 & \lambda + c_{\text n} g_{\textrm{n}}\leq 0 \\ \int_{{\gamma}^{(1)}_{\text c}} -\dfrac{1}{c_{\text n}} \Delta \lambda \delta{\lambda}\,\,{ {\text d} {\gamma}} & \lambda + c_{\text n} g_{\textrm{n}}> 0 \\ \end{array} \right. \]

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, 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 953 of file SimpleContact.cpp.

955  {
957 
958  const int nb_row = row_data.getIndices().size();
959  const int nb_col = col_data.getIndices().size();
960 
961  if (nb_row && nb_col) {
962  const int nb_gauss_pts = row_data.getN().size1();
963  const double area_slave =
964  commonDataSimpleContact->areaSlave; // same area in master and slave
965 
966  NN.resize(nb_row, nb_col, false);
967  NN.clear();
968 
969  auto t_lagrange_slave =
970  getFTensor0FromVec(*commonDataSimpleContact->lagMultAtGaussPtsPtr);
971  auto t_w = getFTensor0IntegrationWeightSlave();
972 
973  auto t_aug_lambda_ptr =
974  getFTensor0FromVec(*commonDataSimpleContact->augmentedLambdasPtr);
975 
976  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
977 
978  if (t_aug_lambda_ptr <= 0) {
979  ++t_w;
980  ++t_aug_lambda_ptr;
981  ++t_lagrange_slave;
982  continue;
983  }
984 
985  const double val_s = -t_w * area_slave / cN;
986 
988  &*NN.data().begin());
989 
990  auto t_base_lambda_row = row_data.getFTensor0N(gg, 0);
991  for (int bbr = 0; bbr != nb_row; ++bbr) {
992  auto t_base_lambda_col = col_data.getFTensor0N(gg, 0);
993  const double s = val_s * t_base_lambda_row;
994  for (int bbc = 0; bbc != nb_col; ++bbc) {
995 
996  t_mat += s * t_base_lambda_col;
997 
998  ++t_mat;
999  ++t_base_lambda_col; // update cols
1000  }
1001  ++t_base_lambda_row; // update rows
1002  }
1003 
1004  ++t_lagrange_slave;
1005  ++t_w;
1006  ++t_aug_lambda_ptr;
1007  }
1008  CHKERR MatSetValues(getSNESB(), row_data, col_data, &*NN.data().begin(),
1009  ADD_VALUES);
1010  }
1011 
1013 }

Member Data Documentation

◆ cN

const double SimpleContactProblem::OpGapConstraintAugmentedOverLambda::cN
private

Definition at line 2039 of file SimpleContact.hpp.

◆ commonDataSimpleContact

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

Definition at line 2038 of file SimpleContact.hpp.

◆ NN

MatrixDouble SimpleContactProblem::OpGapConstraintAugmentedOverLambda::NN
private

Definition at line 2040 of file SimpleContact.hpp.


The documentation for this struct was generated from the following files:
MoFEM::MatSetValues
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
Definition: EntitiesFieldData.hpp:1644
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
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
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::ContactOp
ContactPrismElementForcesAndSourcesCore::UserDataOperator ContactOp
Definition: SimpleContact.hpp:30
FTensor::Tensor0
Definition: Tensor0.hpp:16
SimpleContactProblem::OpGapConstraintAugmentedOverLambda::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:2038
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
SimpleContactProblem::OpGapConstraintAugmentedOverLambda::cN
const double cN
Definition: SimpleContact.hpp:2039
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
SimpleContactProblem::OpGapConstraintAugmentedOverLambda::NN
MatrixDouble NN
Definition: SimpleContact.hpp:2040