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

LHS-operator for the simple contact element with Augmented Lagrangian Method. More...

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

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

Public Member Functions

 OpCalContactAugmentedTractionOverSpatialMasterMaster (const string field_name, const string field_name_2, const double cn_value, boost::shared_ptr< CommonDataSimpleContact > common_data_contact)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Integrates virtual work on master side , \( \delta W_{\text c}\), derivative with respect to spatial positions of the master side and assembles its 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 with Augmented Lagrangian Method.

Integrates Lagrange multipliers virtual work, \( \delta W_{\text c}\) derivative with respect to spatial positions on master side and assembles components of the LHS matrix.

Definition at line 1706 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpCalContactAugmentedTractionOverSpatialMasterMaster()

SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterMaster::OpCalContactAugmentedTractionOverSpatialMasterMaster ( const string  field_name,
const string  field_name_2,
const double  cn_value,
boost::shared_ptr< CommonDataSimpleContact common_data_contact 
)
inline

Definition at line 1709 of file SimpleContact.hpp.

1714  field_name, field_name_2, UserDataOperator::OPROWCOL,
1716  FACEMASTERMASTER),
1717  cN(cn_value), commonDataSimpleContact(common_data_contact) {
1718  sYmm = false; // This will make sure to loop over all intities (e.g.
1719  // for order=2 it will make doWork to loop 16 time)
1720  }

Member Function Documentation

◆ doWork()

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

Integrates virtual work on master side , \( \delta W_{\text c}\), derivative with respect to spatial positions of the master side and assembles its components to LHS global matrix.

Computes linearisation of integrated on slave side complementarity function and assembles Lagrange multipliers virtual work \( \delta W_{\text c}\) on master side with respect to spatial positions of the master side and assembles components to LHS global matrix

\[ {\text D} {\delta W^{(2)}_{\text c}(\lambda, \delta \mathbf{x}^{(2)}})[\Delta {\mathbf{x}^{(2)}}] \,\, = \left\{ \begin{array}{ll} \int_{{\gamma}^{(1)}_{\text c}} c_{\textrm n}\Delta {\mathbf{x}^{(2)}} \cdot [{\mathbf{n}}_{\rm c} \otimes {\mathbf{n}}_{\rm c}] \cdot \delta{\mathbf{x}^{(2)}}\,\,{ {\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}^{(2)}\) are the coordinates of the overlapping gauss points at master triangles, \( c_{\textrm n}\) is the regularisation/augmentation parameter of stress dimensions and \( g_{\textrm{n}}\) is the gap evaluated on the corresponding slave side.

Definition at line 1169 of file SimpleContact.cpp.

1171  {
1173 
1174  // Both sides are needed since both sides contribute their shape
1175  // function to the stiffness matrix
1176  const int nb_row = row_data.getIndices().size();
1177  if (!nb_row)
1179  const int nb_col = col_data.getIndices().size();
1180  if (!nb_col)
1182  const int nb_gauss_pts = row_data.getN().size1();
1183 
1184  int nb_base_fun_row = row_data.getFieldData().size() / 3;
1185  int nb_base_fun_col = col_data.getFieldData().size() / 3;
1186 
1187  auto get_tensor_from_mat = [](MatrixDouble &m, const int r, const int c) {
1189  &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 1, c + 0),
1190  &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 2, c + 0), &m(r + 2, c + 1),
1191  &m(r + 2, c + 2));
1192  };
1193 
1194  auto get_tensor_vec = [](VectorDouble &n) {
1195  return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
1196  };
1197 
1200 
1201  NN.resize(3 * nb_base_fun_row, 3 * nb_base_fun_col, false);
1202  NN.clear();
1203 
1204  auto t_w = getFTensor0IntegrationWeightSlave();
1205 
1206  const double area_master = commonDataSimpleContact->areaSlave;
1207 
1208  auto normal =
1209  get_tensor_vec(commonDataSimpleContact->normalVectorSlavePtr.get()[0]);
1210 
1211  auto t_aug_lambda_ptr =
1212  getFTensor0FromVec(*commonDataSimpleContact->augmentedLambdasPtr);
1213 
1214  for (int gg = 0; gg != nb_gauss_pts; gg++) {
1215 
1216  if (t_aug_lambda_ptr > 0 && std::abs(t_aug_lambda_ptr) > ALM_TOL) {
1217  ++t_w;
1218  ++t_aug_lambda_ptr;
1219  continue;
1220  }
1221 
1222  const double val_m = t_w * area_master * cN;
1223 
1224  FTensor::Tensor0<double *> t_base_master_col(&col_data.getN()(gg, 0));
1225 
1226  for (int bbc = 0; bbc != nb_base_fun_col; bbc++) {
1227 
1228  FTensor::Tensor0<double *> t_base_master_row(&row_data.getN()(gg, 0));
1229  const double m = val_m * t_base_master_col;
1230 
1231  for (int bbr = 0; bbr != nb_base_fun_row; bbr++) {
1232 
1233  auto t_assemble_s = get_tensor_from_mat(NN, 3 * bbr, 3 * bbc);
1234 
1235  t_assemble_s(i, j) += m * normal(i) * normal(j) * t_base_master_row;
1236 
1237  ++t_base_master_row; // update rows
1238  }
1239  ++t_base_master_col; // update cols slave
1240  }
1241  ++t_w;
1242  ++t_aug_lambda_ptr;
1243  }
1244 
1245  CHKERR MatSetValues(getSNESB(), row_data, col_data, &*NN.data().begin(),
1246  ADD_VALUES);
1247 
1249 }

Member Data Documentation

◆ cN

const double SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterMaster::cN
private

Definition at line 1763 of file SimpleContact.hpp.

◆ commonDataSimpleContact

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

Definition at line 1762 of file SimpleContact.hpp.

◆ NN

MatrixDouble SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterMaster::NN
private

Definition at line 1764 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:447
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:1631
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterMaster::NN
MatrixDouble NN
Definition: SimpleContact.hpp:1764
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1241
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
sdf.r
int r
Definition: sdf.py:8
FTensor::Tensor2< double *, 3, 3 >
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterMaster::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:1762
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1201
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
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', 3 >
convert.n
n
Definition: convert.py:82
FTensor::Tensor0
Definition: Tensor0.hpp:16
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:1305
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterMaster::cN
const double cN
Definition: SimpleContact.hpp:1763
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346