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

Public Member Functions

 OpCalContactAugmentedTractionOverLambdaMasterSlave (const string field_name, const string lagrange_field_name, boost::shared_ptr< CommonDataSimpleContact > common_data_contact)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 Integrates virtual work on master side , \( \delta W_{\text c}\), derivative with respect to Lagrange multipliers on slave side and assembles its components to LHS global matrix. More...
 

Private Attributes

boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 
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 Lagrange multipliers on master side and assembles components of the LHS global matrix.

Definition at line 1577 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpCalContactAugmentedTractionOverLambdaMasterSlave()

SimpleContactProblem::OpCalContactAugmentedTractionOverLambdaMasterSlave::OpCalContactAugmentedTractionOverLambdaMasterSlave ( const string  field_name,
const string  lagrange_field_name,
boost::shared_ptr< CommonDataSimpleContact common_data_contact 
)
inline

Definition at line 1580 of file SimpleContact.hpp.

1583  : ContactOp(field_name, lagrange_field_name, UserDataOperator::OPROWCOL,
1584  ContactOp::FACEMASTERSLAVE),
1585  commonDataSimpleContact(common_data_contact) {
1586  sYmm = false; // This will make sure to loop over all intities (e.g.
1587  // for order=2 it will make doWork to loop 16 time)
1588  }

Member Function Documentation

◆ doWork()

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

Integrates virtual work on master side , \( \delta W_{\text c}\), derivative with respect to Lagrange multipliers 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 Lagrange multipliers.

\[ {\text D} {\delta W^{(2)}_{\text c}(\lambda, \delta \mathbf{x}^{(2)}})[\Delta \lambda] \,\, = \left\{ \begin{array}{ll} \int_{{\gamma}^{(1)}_{\text c}} \Delta \lambda \mathbf{n}(\mathbf{x}^{(1)}) \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 804 of file SimpleContact.cpp.

805  {
807 
808  // Both sides are needed since both sides contribute their shape
809  // function to the stiffness matrix
810  const int nb_row = row_data.getIndices().size();
811  const int nb_col = col_data.getIndices().size();
812 
813  if (nb_row && nb_col) {
814  const int nb_gauss_pts = row_data.getN().size1();
815 
816  int nb_base_fun_row = row_data.getFieldData().size() / 3;
817  int nb_base_fun_col = col_data.getFieldData().size();
818 
819  const double area_master =
820  commonDataSimpleContact->areaSlave; // same area in master and slave
821 
822  auto get_tensor_from_mat = [](MatrixDouble &m, const int r, const int c) {
823  return FTensor::Tensor1<double *, 3>(&m(r + 0, c), &m(r + 1, c),
824  &m(r + 2, c));
825  };
826 
827  auto get_tensor_vec = [](VectorDouble &n) {
828  return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
829  };
830 
832 
833  NN.resize(3 * nb_base_fun_row, nb_base_fun_col, false);
834  NN.clear();
835 
836  auto const_unit_n =
837  get_tensor_vec(*(commonDataSimpleContact->normalVectorSlavePtr.get()));
838 
839  auto t_w = getFTensor0IntegrationWeightSlave();
840 
841  auto t_aug_lambda_ptr =
842  getFTensor0FromVec(*commonDataSimpleContact->augmentedLambdasPtr);
843 
844  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
845 
846  if (t_aug_lambda_ptr > 0 && std::abs(t_aug_lambda_ptr) > ALM_TOL) {
847  ++t_w;
848  ++t_aug_lambda_ptr;
849  continue;
850  }
851  double val_m = t_w * area_master;
852  auto t_base_master = row_data.getFTensor0N(gg, 0);
853 
854  for (int bbr = 0; bbr != nb_base_fun_row; ++bbr) {
855  auto t_assemble_m = get_tensor_from_mat(NN, 3 * bbr, 0);
856  auto t_base_lambda = col_data.getFTensor0N(gg, 0);
857  const double m = val_m * t_base_master;
858  for (int bbc = 0; bbc != nb_base_fun_col; ++bbc) {
859  const double n = m * t_base_lambda;
860  t_assemble_m(i) += n * const_unit_n(i);
861  ++t_assemble_m;
862  ++t_base_lambda; // update cols slave
863  }
864  ++t_base_master; // update rows master
865  }
866  ++t_w;
867  ++t_aug_lambda_ptr;
868  }
869 
870  CHKERR MatSetValues(getSNESB(), row_data, col_data, &*NN.data().begin(),
871  ADD_VALUES);
872  }
873 
875 }

Member Data Documentation

◆ commonDataSimpleContact

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

Definition at line 1629 of file SimpleContact.hpp.

◆ NN

MatrixDouble SimpleContactProblem::OpCalContactAugmentedTractionOverLambdaMasterSlave::NN
private

Definition at line 1630 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:1631
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
SimpleContactProblem::OpCalContactAugmentedTractionOverLambdaMasterSlave::NN
MatrixDouble NN
Definition: SimpleContact.hpp:1630
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
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
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1489
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
SimpleContactProblem::OpCalContactAugmentedTractionOverLambdaMasterSlave::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:1629
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: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
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:1305
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:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346