v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
SimpleContactProblem::OpGapConstraintAugmentedOverSpatialMaster Struct Reference

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

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

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

Public Member Functions

 OpGapConstraintAugmentedOverSpatialMaster (const string field_name, const string lagrange_field_name, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, const double cn)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 Integrates the conditions that fulfil KKT conditions at master face gauss points and assembles 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.

Integrates variation on the slave sid the conditions that fulfil KKT conditions with respect to Spatial positions on the master side and assembles components to LHS global matrix.

Definition at line 2051 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpGapConstraintAugmentedOverSpatialMaster()

SimpleContactProblem::OpGapConstraintAugmentedOverSpatialMaster::OpGapConstraintAugmentedOverSpatialMaster ( const string  field_name,
const string  lagrange_field_name,
boost::shared_ptr< CommonDataSimpleContact common_data_contact,
const double  cn 
)
inline

Definition at line 2053 of file SimpleContact.hpp.

2057 : ContactOp(lagrange_field_name, field_name, UserDataOperator::OPROWCOL,
2058 ContactOp::FACESLAVEMASTER),
2059 commonDataSimpleContact(common_data_contact), cN(cn) {
2060 sYmm = false; // This will make sure to loop over all entities (e.g.
2061 // for order=2 it will make doWork to loop 16 time)
2062 }
constexpr auto field_name
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
ContactPrismElementForcesAndSourcesCore::UserDataOperator ContactOp

Member Function Documentation

◆ doWork()

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

Integrates the conditions that fulfil KKT conditions at master face gauss points and assembles components to LHS global matrix.

Integrates variation of the expresion that fulfils KKT conditions with respect to spatial positions in the integral sense and assembles components to LHS global matrix.

\[ {\text D}{\overline C(\lambda, \mathbf{x}^{(1)}, \delta \lambda)}[\Delta \mathbf{x}^{(2)}] = \left\{ \begin{array}{ll} \int_{{\gamma}^{(1)}_{\text c}} c_{\text n} \Delta \mathbf{x}^{(2)} \cdot {\mathbf{n}}_{\rm c} \delta{\lambda}\,\,{ {\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 and master triangles for \(i = 1\) and \(i = 2\), respectively. Furthermore, \( c_{\text n}\) works as an augmentation parameter and affects convergence, 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) \]

Definition at line 1016 of file SimpleContact.cpp.

1018 {
1020
1021 const int nb_row = row_data.getIndices().size();
1022 const int nb_col = col_data.getIndices().size();
1023
1024 if (nb_row && nb_col) {
1025
1026 const int nb_gauss_pts = row_data.getN().size1();
1027 int nb_base_fun_row = row_data.getFieldData().size();
1028 int nb_base_fun_col = col_data.getFieldData().size() / 3;
1029
1030 const double area_master =
1031 commonDataSimpleContact->areaSlave; // same area in master and slave
1032
1033 const double area_slave =
1034 commonDataSimpleContact->areaSlave; // same area in master and slave
1035
1036 NN.resize(nb_base_fun_row, 3 * nb_base_fun_col, false);
1037 NN.clear();
1038
1039 auto get_tensor_from_vec = [](VectorDouble &n) {
1040 return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
1041 };
1042
1044
1045 auto t_const_unit_n =
1046 get_tensor_from_vec(*(commonDataSimpleContact->normalVectorSlavePtr));
1047
1048 auto t_const_unit_master =
1049 get_tensor_from_vec(*(commonDataSimpleContact->normalVectorMasterPtr));
1050
1051 auto t_aug_lambda_ptr =
1052 getFTensor0FromVec(*commonDataSimpleContact->augmentedLambdasPtr);
1053
1054 auto t_w = getFTensor0IntegrationWeightSlave();
1055 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1056
1057 if (t_aug_lambda_ptr > 0 && std::abs(t_aug_lambda_ptr) > ALM_TOL) {
1058 ++t_w;
1059 ++t_aug_lambda_ptr;
1060 continue;
1061 }
1062
1063 const double val_m = t_w * area_slave;
1064
1065 auto t_base_lambda = row_data.getFTensor0N(gg, 0);
1066 for (int bbr = 0; bbr != nb_base_fun_row; ++bbr) {
1067
1068 auto t_base_master = col_data.getFTensor0N(gg, 0);
1070 &NN(bbr, 0), &NN(bbr, 1), &NN(bbr, 2)};
1071 const double m = val_m * t_base_lambda;
1072
1073 for (int bbc = 0; bbc != nb_base_fun_col; ++bbc) {
1074
1075 t_mat(i) += t_const_unit_n(i) * m * t_base_master;
1076
1077 ++t_base_master; // update rows
1078 ++t_mat;
1079 }
1080 ++t_base_lambda; // update cols master
1081 }
1082
1083 ++t_aug_lambda_ptr;
1084 ++t_w;
1085 }
1086
1087 CHKERR MatSetValues(getSNESB(), row_data, col_data, &*NN.data().begin(),
1088 ADD_VALUES);
1089 }
1090
1092}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
FTensor::Index< 'n', SPACE_DIM > n
FTensor::Index< 'm', SPACE_DIM > m
FTensor::Index< 'i', SPACE_DIM > i
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorDouble & getFieldData() const
get dofs values
const VectorInt & getIndices() const
Get global indices of dofs on entity.
static constexpr double ALM_TOL

Member Data Documentation

◆ cN

const double SimpleContactProblem::OpGapConstraintAugmentedOverSpatialMaster::cN
private

Definition at line 2104 of file SimpleContact.hpp.

◆ commonDataSimpleContact

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

Definition at line 2103 of file SimpleContact.hpp.

◆ NN

MatrixDouble SimpleContactProblem::OpGapConstraintAugmentedOverSpatialMaster::NN
private

Definition at line 2105 of file SimpleContact.hpp.


The documentation for this struct was generated from the following files: