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

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

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

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

Public Member Functions

 OpGapConstraintAugmentedOverSpatialSlave (const string field_name, 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 on the slave side variation of the conditions that fulfil KKT conditions with respect to Spatial positions on the slave side and assembles components to LHS global matrix.

Definition at line 2116 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpGapConstraintAugmentedOverSpatialSlave()

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

Definition at line 2118 of file SimpleContact.hpp.

2122  : ContactOp(lagrange_field_name, field_name, UserDataOperator::OPROWCOL,
2123  ContactOp::FACESLAVESLAVE),
2124  cN(cn), commonDataSimpleContact(common_data_contact) {
2125  sYmm = false; // This will make sure to loop over all entities (e.g.
2126  // for order=2 it will make doWork to loop 16 time)
2127  }

Member Function Documentation

◆ doWork()

MoFEMErrorCode SimpleContactProblem::OpGapConstraintAugmentedOverSpatialSlave::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 spatial positions on slave side components to LHS global matrix.

\[ {\text D}{\overline C(\lambda, \mathbf{x}^{(1)}, \delta \lambda)}[\Delta \mathbf{x}^{(1)}] = \left\{ \begin{array}{ll} \int_{{\gamma}^{(1)}_{\text c}} -c_{\text n} \Delta \mathbf{x}^{(1)} \cdot {\mathbf{n}}_{\rm c} \delta{\lambda}\,\,{ {\text d} {\gamma}} & \lambda + c_{\text n} g_{\textrm{n}}\leq 0 \\ 0 & \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}^{(1)}\) 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 1095 of file SimpleContact.cpp.

1097  {
1099 
1100  const int nb_row = row_data.getIndices().size();
1101  const int nb_col = col_data.getIndices().size();
1102 
1103  if (nb_row && nb_col) {
1104 
1105  const int nb_gauss_pts = row_data.getN().size1();
1106  int nb_base_fun_row = row_data.getFieldData().size();
1107  int nb_base_fun_col = col_data.getFieldData().size() / 3;
1108 
1109  const double area_slave =
1110  commonDataSimpleContact->areaSlave; // same area in master and slave
1111 
1112  NN.resize(nb_base_fun_row, 3 * nb_base_fun_col, false);
1113  NN.clear();
1114 
1115  auto get_tensor_from_vec = [](VectorDouble &n) {
1116  return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
1117  };
1118 
1120 
1121  auto t_const_unit_n =
1122  get_tensor_from_vec(*(commonDataSimpleContact->normalVectorSlavePtr));
1123 
1124  auto t_aug_lambda_ptr =
1125  getFTensor0FromVec(*commonDataSimpleContact->augmentedLambdasPtr);
1126 
1127  auto t_w = getFTensor0IntegrationWeightSlave();
1128  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1129 
1130  if (t_aug_lambda_ptr > 0 && std::abs(t_aug_lambda_ptr) > ALM_TOL) {
1131  ++t_w;
1132  ++t_aug_lambda_ptr;
1133  continue;
1134  }
1135 
1136  const double val_m = t_w * area_slave;
1137 
1138  auto t_base_lambda = row_data.getFTensor0N(gg, 0);
1139  for (int bbr = 0; bbr != nb_base_fun_row; ++bbr) {
1140 
1141  auto t_base_slave = col_data.getFTensor0N(gg, 0);
1143  &NN(bbr, 0), &NN(bbr, 1), &NN(bbr, 2)};
1144  const double m = val_m * t_base_lambda;
1145 
1146  for (int bbc = 0; bbc != nb_base_fun_col; ++bbc) {
1147 
1148  t_mat(i) -= t_const_unit_n(i) * m * t_base_slave;
1149 
1150  ++t_base_slave; // update rows
1151  ++t_mat;
1152  }
1153  ++t_base_lambda; // update cols master
1154  }
1155 
1156  ++t_aug_lambda_ptr;
1157  ++t_w;
1158  }
1159 
1160  CHKERR MatSetValues(getSNESB(), row_data, col_data, &*NN.data().begin(),
1161  ADD_VALUES);
1162  }
1163 
1165 }

Member Data Documentation

◆ cN

const double SimpleContactProblem::OpGapConstraintAugmentedOverSpatialSlave::cN
private

Definition at line 2170 of file SimpleContact.hpp.

◆ commonDataSimpleContact

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

Definition at line 2169 of file SimpleContact.hpp.

◆ NN

MatrixDouble SimpleContactProblem::OpGapConstraintAugmentedOverSpatialSlave::NN
private

Definition at line 2171 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
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1254
SimpleContactProblem::OpGapConstraintAugmentedOverSpatialSlave::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:2169
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
SimpleContactProblem::OpGapConstraintAugmentedOverSpatialSlave::NN
MatrixDouble NN
Definition: SimpleContact.hpp:2171
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::ALM_TOL
static constexpr double ALM_TOL
Definition: SimpleContact.hpp:65
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
SimpleContactProblem::OpGapConstraintAugmentedOverSpatialSlave::cN
const double cN
Definition: SimpleContact.hpp:2170
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', 3 >
convert.n
n
Definition: convert.py:82
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
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
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