v0.14.0
Public Member Functions | Private Attributes | List of all members
SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterSlave 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::OpCalContactAugmentedTractionOverSpatialMasterSlave:
[legend]
Collaboration diagram for SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterSlave:
[legend]

Public Member Functions

 OpCalContactAugmentedTractionOverSpatialMasterSlave (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 on slave 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 1776 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpCalContactAugmentedTractionOverSpatialMasterSlave()

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

Definition at line 1779 of file SimpleContact.hpp.

1784  field_name, field_name_2, UserDataOperator::OPROWCOL,
1786  FACEMASTERSLAVE),
1787  cN(cn_value), commonDataSimpleContact(common_data_contact) {
1788  sYmm = false; // This will make sure to loop over all intities (e.g.
1789  // for order=2 it will make doWork to loop 16 time)
1790  }

Member Function Documentation

◆ doWork()

MoFEMErrorCode SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterSlave::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 on slave side and assembles its components to LHS global matrix.

Computes linearisation of virtual work on master side integrated on the slave side and assembles the components of its derivative over spatial positions on slave side

\[ {\text D} {\delta W^{(2)}_{\text c}(\lambda, \delta \mathbf{x}^{(2)}})[\Delta {\mathbf{x}^{(1)}}] \,\, = \left\{ \begin{array}{ll} \int_{{\gamma}^{(1)}_{\text c}} -c_{\textrm n}\Delta {\mathbf{x}^{(1)}} \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)}\) and \(\mathbf{x}^{(1)}\) are the coordinates of the overlapping gauss points at master and slave triangles, respectively. Also, \( 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 1253 of file SimpleContact.cpp.

1255  {
1257 
1258  // Both sides are needed since both sides contribute their shape
1259  // function to the stiffness matrix
1260  const int nb_row = row_data.getIndices().size();
1261  if (!nb_row)
1263  const int nb_col = col_data.getIndices().size();
1264  if (!nb_col)
1266  const int nb_gauss_pts = row_data.getN().size1();
1267 
1268  int nb_base_fun_row = row_data.getFieldData().size() / 3;
1269  int nb_base_fun_col = col_data.getFieldData().size() / 3;
1270 
1271  auto get_tensor_from_mat = [](MatrixDouble &m, const int r, const int c) {
1273  &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 1, c + 0),
1274  &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 2, c + 0), &m(r + 2, c + 1),
1275  &m(r + 2, c + 2));
1276  };
1277 
1278  auto get_tensor_vec = [](VectorDouble &n) {
1279  return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
1280  };
1281 
1284 
1285  NN.resize(3 * nb_base_fun_row, 3 * nb_base_fun_col, false);
1286  NN.clear();
1287 
1288  auto t_w = getFTensor0IntegrationWeightSlave();
1289 
1290  const double area_master = commonDataSimpleContact->areaSlave;
1291 
1292  auto normal =
1293  get_tensor_vec(commonDataSimpleContact->normalVectorSlavePtr.get()[0]);
1294 
1295  auto t_aug_lambda_ptr =
1296  getFTensor0FromVec(*commonDataSimpleContact->augmentedLambdasPtr);
1297 
1298  for (int gg = 0; gg != nb_gauss_pts; gg++) {
1299 
1300  if (t_aug_lambda_ptr > 0 && std::abs(t_aug_lambda_ptr) > ALM_TOL) {
1301  ++t_w;
1302  ++t_aug_lambda_ptr;
1303  continue;
1304  }
1305 
1306  const double val_m = t_w * area_master * cN;
1307 
1308  FTensor::Tensor0<double *> t_base_slave_col(&col_data.getN()(gg, 0));
1309 
1310  for (int bbc = 0; bbc != nb_base_fun_col; bbc++) {
1311 
1312  FTensor::Tensor0<double *> t_base_master_row(&row_data.getN()(gg, 0));
1313  const double m = val_m * t_base_slave_col;
1314 
1315  for (int bbr = 0; bbr != nb_base_fun_row; bbr++) {
1316 
1317  auto t_assemble_s = get_tensor_from_mat(NN, 3 * bbr, 3 * bbc);
1318 
1319  t_assemble_s(i, j) -= m * normal(i) * normal(j) * t_base_master_row;
1320 
1321  ++t_base_master_row; // update rows
1322  }
1323  ++t_base_slave_col; // update cols slave
1324  }
1325  ++t_w;
1326  ++t_aug_lambda_ptr;
1327  }
1328 
1329  CHKERR MatSetValues(getSNESB(), row_data, col_data, &*NN.data().begin(),
1330  ADD_VALUES);
1331 
1333 }

Member Data Documentation

◆ cN

const double SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterSlave::cN
private

Definition at line 1833 of file SimpleContact.hpp.

◆ commonDataSimpleContact

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

Definition at line 1832 of file SimpleContact.hpp.

◆ NN

MatrixDouble SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterSlave::NN
private

Definition at line 1834 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
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::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1254
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: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
SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterSlave::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:1832
SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialMasterSlave::cN
const double cN
Definition: SimpleContact.hpp:1833
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:1318
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
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::OpCalContactAugmentedTractionOverSpatialMasterSlave::NN
MatrixDouble NN
Definition: SimpleContact.hpp:1834