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

Public Member Functions

 OpCalContactAugmentedTractionOverSpatialSlaveSlave (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 slave side, \( \delta W_{\text c}\), derivative with respect to slave spatial positions 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 Spatial position on slave side 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 1846 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpCalContactAugmentedTractionOverSpatialSlaveSlave()

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

Definition at line 1849 of file SimpleContact.hpp.

1854  field_name, field_name_2, UserDataOperator::OPROWCOL,
1856  FACESLAVESLAVE),
1857  cN(cn_value), commonDataSimpleContact(common_data_contact) {
1858  sYmm = false; // This will make sure to loop over all intities (e.g.
1859  // for order=2 it will make doWork to loop 16 time)
1860  }

Member Function Documentation

◆ doWork()

MoFEMErrorCode SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialSlaveSlave::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 slave side, \( \delta W_{\text c}\), derivative with respect to slave spatial positions and assembles its components to LHS global matrix.

Computes linearisation of virtual work on slave side integrated on the slave side and assembles the components of its derivative over Lagrange multipliers.

\[ {\text D} {\delta W^{(1)}_{\text c}(\lambda, \delta \mathbf{x}^{(1)}})[\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}^{(1)}}\,\,{ {\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 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 1337 of file SimpleContact.cpp.

1339  {
1341 
1342  // Both sides are needed since both sides contribute their shape
1343  // function to the stiffness matrix
1344  const int nb_row = row_data.getIndices().size();
1345  if (!nb_row)
1347  const int nb_col = col_data.getIndices().size();
1348  if (!nb_col)
1350  const int nb_gauss_pts = row_data.getN().size1();
1351 
1352  int nb_base_fun_row = row_data.getFieldData().size() / 3;
1353  int nb_base_fun_col = col_data.getFieldData().size() / 3;
1354 
1355  auto get_tensor_from_mat = [](MatrixDouble &m, const int r, const int c) {
1357  &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 1, c + 0),
1358  &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 2, c + 0), &m(r + 2, c + 1),
1359  &m(r + 2, c + 2));
1360  };
1361 
1362  auto get_tensor_vec = [](VectorDouble &n) {
1363  return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
1364  };
1365 
1368 
1369  NN.resize(3 * nb_base_fun_row, 3 * nb_base_fun_col, false);
1370  NN.clear();
1371 
1372  auto t_w = getFTensor0IntegrationWeightSlave();
1373 
1374  const double area_slave = commonDataSimpleContact->areaSlave;
1375 
1376  auto normal =
1377  get_tensor_vec(commonDataSimpleContact->normalVectorSlavePtr.get()[0]);
1378 
1379  auto t_aug_lambda_ptr =
1380  getFTensor0FromVec(*commonDataSimpleContact->augmentedLambdasPtr);
1381 
1382  for (int gg = 0; gg != nb_gauss_pts; gg++) {
1383 
1384  if (t_aug_lambda_ptr > 0 && std::abs(t_aug_lambda_ptr) > ALM_TOL) {
1385  ++t_w;
1386  ++t_aug_lambda_ptr;
1387  continue;
1388  }
1389 
1390  const double val_s = t_w * area_slave * cN;
1391 
1392  FTensor::Tensor0<double *> t_base_slave_col(&col_data.getN()(gg, 0));
1393 
1394  for (int bbc = 0; bbc != nb_base_fun_col; bbc++) {
1395 
1396  FTensor::Tensor0<double *> t_base_slave_row(&row_data.getN()(gg, 0));
1397  const double s = val_s * t_base_slave_col;
1398 
1399  for (int bbr = 0; bbr != nb_base_fun_row; bbr++) {
1400 
1401  auto t_assemble_s = get_tensor_from_mat(NN, 3 * bbr, 3 * bbc);
1402 
1403  t_assemble_s(i, j) += s * normal(i) * normal(j) * t_base_slave_row;
1404 
1405  ++t_base_slave_row; // update rows
1406  }
1407  ++t_base_slave_col; // update cols slave
1408  }
1409  ++t_w;
1410  ++t_aug_lambda_ptr;
1411  }
1412 
1413  CHKERR MatSetValues(getSNESB(), row_data, col_data, &*NN.data().begin(),
1414  ADD_VALUES);
1415 
1417 }

Member Data Documentation

◆ cN

const double SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialSlaveSlave::cN
private

Definition at line 1903 of file SimpleContact.hpp.

◆ commonDataSimpleContact

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

Definition at line 1902 of file SimpleContact.hpp.

◆ NN

MatrixDouble SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialSlaveSlave::NN
private

Definition at line 1904 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::OpCalContactAugmentedTractionOverSpatialSlaveSlave::NN
MatrixDouble NN
Definition: SimpleContact.hpp:1904
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
SimpleContactProblem::OpCalContactAugmentedTractionOverSpatialSlaveSlave::cN
const double cN
Definition: SimpleContact.hpp:1903
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::OpCalContactAugmentedTractionOverSpatialSlaveSlave::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:1902
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: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