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

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

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

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

Public Member Functions

MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 Integrates linearisation of the complementarity function at slave face gauss points and assembles components to LHS global matrix. More...
 
 OpCalDerIntCompFunSlaveSlave_dX (const string lagrange_field_name, const string mesh_nodes_field, boost::shared_ptr< double > cn, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, int row_rank, const int col_rank)
 
- Public Member Functions inherited from SimpleContactProblem::OpContactALELhs
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 
MoFEMErrorCode aSsemble (EntData &row_data, EntData &col_data)
 
 OpContactALELhs (const string field_name_1, const string field_name_2, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, const ContactOp::FaceType face_type, const int rank_row, const int rank_col)
 

Private Attributes

boost::shared_ptr< doublecNPtr
 

Additional Inherited Members

- Public Attributes inherited from SimpleContactProblem::OpContactALELhs
boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 
MatrixDouble matLhs
 
VectorInt rowIndices
 
VectorInt colIndices
 
int row_nb_dofs
 
int col_nb_dofs
 
int nb_gauss_pts
 
int nb_base_fun_row
 
int nb_base_fun_col
 
int rankRow
 
int rankCol
 

Detailed Description

LHS-operator for the simple contact element.

Integrates the variation with respect to slave material positions of the complementarity function to fulfill KKT conditions in the integral sense and assembles components to LHS global matrix.

Definition at line 3236 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpCalDerIntCompFunSlaveSlave_dX()

SimpleContactProblem::OpCalDerIntCompFunSlaveSlave_dX::OpCalDerIntCompFunSlaveSlave_dX ( const string  lagrange_field_name,
const string  mesh_nodes_field,
boost::shared_ptr< double cn,
boost::shared_ptr< CommonDataSimpleContact common_data_contact,
int  row_rank,
const int  col_rank 
)
inline
Parameters
lagrange_field_nameString of field name for lagrange multipliers for rows
mesh_nodes_fieldString of field name for material positions for columns
cnregularisation/augmentation parameter affecting convergence
common_data_contactPointer to the common data for simple contact element
row_rankParameter setting the dimension of the associated field for rows (in this case is 1)
col_rankParameter setting the dimension of the associated field for cols (in this case is 3)

Definition at line 3315 of file SimpleContact.hpp.

3320  : OpContactALELhs(lagrange_field_name, mesh_nodes_field,
3321  common_data_contact, ContactOp::FACESLAVESLAVE,
3322  row_rank, col_rank),
3323  cNPtr(cn) {
3324  sYmm = false; // This will make sure to loop over all intities (e.g.
3325  // for order=2 it will make doWork to loop 16 time)
3326  }

Member Function Documentation

◆ iNtegrate()

MoFEMErrorCode SimpleContactProblem::OpCalDerIntCompFunSlaveSlave_dX::iNtegrate ( EntData row_data,
EntData col_data 
)
virtual

Integrates linearisation of the complementarity function at slave face gauss points and assembles components to LHS global matrix.

Integrates and assembles the variation with respect to slave spatial positions of the complementarity function to fulfill KKT conditions in the integral sense and assembles components to LHS global matrix.

\[ {\text D}{\overline C(\lambda, \mathbf{x}^{(i)}, \mathbf{X}^{(1)}, \delta \lambda)}[\Delta \mathbf{X}^{(1)}] = \int_{{{\mathcal{T}}^{(1)}_{\xi}}} -\left[ \textrm{D}\mathbf{N}^{(1)}(\mathbf{X}^{(1)})[\Delta\mathbf{X}^{(1)}] - \frac{\mathbf{N}^{(1)}(\mathbf{X}^{(1)})}{||{{\mathbf N}^{(1)}(\xi, \eta)}||} \textrm{D}\mathbf{N}^{(1)}(\mathbf{X}^{(1)})[\Delta\mathbf{X}^{(1)}] \cdot \frac{\mathbf{N}^{(1)}(\mathbf{X}^{(1)})}{||{{\mathbf N}^{(1)}(\xi, \eta)}||} \right] \cdot \left( {\mathbf{x}}^{(1)} - {\mathbf{x}}^{(2)}\right) c_{\text n} \left( 1 + {\text {sign}}\left( \lambda - c_{\text n} g_{\textrm{n}} \right) {\left| \lambda - c_{\text n} g_{\textrm{n}}\right|}^{r-1}\right) \delta{{\lambda}} \,\,{ {\text d} {\xi} {\text d} {\eta}} \\ + \int_{{{\mathcal{T}}^{(1)}_{\xi}}} \overline C(\lambda, \mathbf{x}^{(i)}, \mathbf{X}^{(1)}, \delta \lambda) \textrm{D}\mathbf{N}^{(1)}(\mathbf{X}^{(1)})[\Delta\mathbf{X}^{(1)}] \cdot \frac{\mathbf{N}^{(1)}(\mathbf{X}^{(1)})}{||{{\mathbf N}^{(1)}(\xi, \eta)}||} \delta{{\lambda}} \,\,{ {\text d} {\xi} {\text d} {\eta}} \]

where

\[ \textrm{D}\mathbf{N}^{(1)}(\mathbf{X}^{(1)})[\Delta\mathbf{X}^{(1)}] = \frac{\partial{\Delta\mathbf{X}}^{(1)}} {\partial\xi} \times \frac{\partial {\mathbf{X}}^{(1)}}{\partial\eta} + \frac{\partial{\mathbf{X}}^{(1)}} {\partial\xi} \times \frac{\partial {\Delta\mathbf{X}}^{(1)}}{\partial\eta} \]

where \(\xi, \eta\) are coordinates in the parent space \(({\mathcal{T}}^{(1)}_\xi)\) 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, \(r\) is regularisation parameter that here is chosen to be \(r = 1\) 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) \]

Reimplemented from SimpleContactProblem::OpContactALELhs.

Definition at line 3716 of file SimpleContact.cpp.

3717  {
3719 
3720  const int nb_row = row_data.getIndices().size();
3721  if (!nb_row)
3723  const int nb_col = col_data.getIndices().size();
3724  if (!nb_col)
3726  const int nb_gauss_pts = row_data.getN().size1();
3727  int nb_base_fun_row = row_data.getFieldData().size();
3728  int nb_base_fun_col = col_data.getFieldData().size() / 3;
3729 
3730  matLhs.resize(nb_base_fun_row, 3 * nb_base_fun_col, false);
3731  matLhs.clear();
3732 
3733  auto get_tensor_vec = [](VectorDouble &n) {
3734  return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
3735  };
3736 
3737  auto get_vec_from_mat = [](MatrixDouble &m, const int r, const int c) {
3738  return FTensor::Tensor1<double *, 3>(&m(r + 0, c + 0), &m(r + 0, c + 1),
3739  &m(r + 0, c + 2));
3740  };
3741 
3745 
3746  auto make_vec_der = [&](auto t_N, auto t_1, auto t_2) {
3748  t_n(i, j) = 0;
3749  t_n(i, j) += FTensor::levi_civita(i, j, k) * t_2(k) * t_N(0);
3750  t_n(i, j) -= FTensor::levi_civita(i, j, k) * t_1(k) * t_N(1);
3751  return t_n;
3752  };
3753 
3754  auto x_m = getFTensor1FromMat<3>(
3755  *commonDataSimpleContact->positionAtGaussPtsMasterPtr);
3756  auto x_s = getFTensor1FromMat<3>(
3757  *commonDataSimpleContact->positionAtGaussPtsSlavePtr);
3758  auto t_lagrange_slave =
3759  getFTensor0FromVec(*commonDataSimpleContact->lagMultAtGaussPtsPtr);
3760 
3761  const double length_normal = commonDataSimpleContact->areaSlave;
3762 
3763  auto normal_at_gp =
3764  get_tensor_vec(*commonDataSimpleContact->normalVectorSlavePtr);
3765 
3766  auto t_w = getFTensor0IntegrationWeightSlave();
3767  auto t_1 = get_tensor_vec(*commonDataSimpleContact->tangentOneVectorSlavePtr);
3768  auto t_2 = get_tensor_vec(*commonDataSimpleContact->tangentTwoVectorSlavePtr);
3769  auto t_gap_gp = getFTensor0FromVec(*commonDataSimpleContact->gapPtr);
3770  const double cn_value = *cNPtr.get();
3771  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
3772  double val_s = t_w * 0.5;
3773  auto t_N = col_data.getFTensor1DiffN<2>(gg, 0);
3774 
3775  for (int bbc = 0; bbc != nb_base_fun_col; ++bbc) {
3776  FTensor::Tensor0<double *> t_base_lambda(&row_data.getN()(gg, 0));
3777 
3778  for (int bbr = 0; bbr != nb_base_fun_row; ++bbr) {
3779  const double s = val_s * t_base_lambda;
3780 
3781  auto t_d_n = make_vec_der(t_N, t_1, t_2);
3782 
3783  auto assemble_mat = get_vec_from_mat(matLhs, bbr, 3 * bbc);
3784 
3785  assemble_mat(j) -=
3786  0.5 *
3787  (t_d_n(i, j) - normal_at_gp(i) * t_d_n(k, j) * normal_at_gp(k)) *
3788  (x_s(i) - x_m(i)) * s * cn_value *
3789  (1 + SimpleContactProblem::Sign(t_lagrange_slave -
3790  cn_value * t_gap_gp));
3791 
3792  assemble_mat(j) += t_d_n(i, j) * normal_at_gp(i) * s *
3794  cn_value, t_gap_gp, t_lagrange_slave);
3795 
3796  ++t_base_lambda; // update rows
3797  }
3798  ++t_N;
3799  }
3800 
3801  ++x_m;
3802  ++x_s;
3803  ++t_lagrange_slave;
3804  ++t_gap_gp;
3805  ++t_w;
3806  }
3807 
3809 }

Member Data Documentation

◆ cNPtr

boost::shared_ptr<double> SimpleContactProblem::OpCalDerIntCompFunSlaveSlave_dX::cNPtr
private

Definition at line 3329 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
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
SimpleContactProblem::OpContactALELhs::nb_base_fun_col
int nb_base_fun_col
Definition: SimpleContact.hpp:560
SimpleContactProblem::OpCalDerIntCompFunSlaveSlave_dX::cNPtr
boost::shared_ptr< double > cNPtr
Definition: SimpleContact.hpp:3329
SimpleContactProblem::Sign
static double Sign(double x)
Definition: SimpleContact.hpp:3568
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
FTensor::levi_civita
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
Definition: Levi_Civita.hpp:617
SimpleContactProblem::OpContactALELhs::matLhs
MatrixDouble matLhs
Definition: SimpleContact.hpp:551
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
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
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
FTensor::Index< 'i', 3 >
convert.n
n
Definition: convert.py:82
MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
Definition: EntitiesFieldData.cpp:526
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
SimpleContactProblem::OpContactALELhs::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:549
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::OpContactALELhs::nb_gauss_pts
int nb_gauss_pts
Definition: SimpleContact.hpp:557
SimpleContactProblem::OpContactALELhs::nb_base_fun_row
int nb_base_fun_row
Definition: SimpleContact.hpp:559
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
SimpleContactProblem::ConstrainFunction
static double ConstrainFunction(const double cn, const double g, const double l)
Definition: SimpleContact.hpp:3587
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
SimpleContactProblem::OpContactALELhs::OpContactALELhs
OpContactALELhs(const string field_name_1, const string field_name_2, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, const ContactOp::FaceType face_type, const int rank_row, const int rank_col)
Definition: SimpleContact.hpp:576