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

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

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

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

Public Member Functions

 OpCalDerIntCompFunOverSpatPosSlaveMaster (const string lagrange_field_name, const string field_name, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, boost::shared_ptr< double > cn)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, 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...
 

Private Attributes

boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 
boost::shared_ptr< doublecNPtr
 
MatrixDouble NN
 

Detailed Description

LHS-operator for the simple contact element.

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

Definition at line 1448 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpCalDerIntCompFunOverSpatPosSlaveMaster()

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

Definition at line 1450 of file SimpleContact.hpp.

1454 : ContactOp(lagrange_field_name, field_name, UserDataOperator::OPROWCOL,
1455 ContactOp::FACESLAVEMASTER),
1456 commonDataSimpleContact(common_data_contact), cNPtr(cn) {
1457 sYmm = false; // This will make sure to loop over all entities (e.g.
1458 // for order=2 it will make doWork to loop 16 time)
1459 }
constexpr auto field_name
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
ContactPrismElementForcesAndSourcesCore::UserDataOperator ContactOp

Member Function Documentation

◆ doWork()

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

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

Integrates the variation with respect to master 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)}, \delta \lambda)}[\Delta \mathbf{x}^{(2)}] = \int_{{\gamma}^{(1)}_{\text c}} \Delta \mathbf{x}^{(2)} \cdot \mathbf{n}(\mathbf{x}^{(1)}) 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} {\gamma}} \]

where \({\gamma}^{(1)}_{\text c}\) is the surface integration domain 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 can be chosen in \([1, 1.1]\) ( \(r = 1\) is the default value) 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 2049 of file SimpleContact.cpp.

2051 {
2053
2054 const int nb_row = row_data.getIndices().size();
2055 const int nb_col = col_data.getIndices().size();
2056
2057 if (nb_row && nb_col) {
2058
2059 const int nb_gauss_pts = row_data.getN().size1();
2060 int nb_base_fun_row = row_data.getFieldData().size();
2061 int nb_base_fun_col = col_data.getFieldData().size() / 3;
2062
2063 const double area_master =
2064 commonDataSimpleContact->areaSlave; // same area in master and slave
2065
2066 const double area_slave =
2067 commonDataSimpleContact->areaSlave; // same area in master and slave
2068
2069 NN.resize(nb_base_fun_row, 3 * nb_base_fun_col, false);
2070 NN.clear();
2071
2072 auto get_tensor_from_vec = [](VectorDouble &n) {
2073 return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
2074 };
2075
2077
2078 auto t_const_unit_n =
2079 get_tensor_from_vec(*(commonDataSimpleContact->normalVectorSlavePtr));
2080
2081 auto t_const_unit_master =
2082 get_tensor_from_vec(*(commonDataSimpleContact->normalVectorMasterPtr));
2083
2084 auto t_lagrange_slave =
2085 getFTensor0FromVec(*commonDataSimpleContact->lagMultAtGaussPtsPtr);
2086 auto t_gap_gp = getFTensor0FromVec(*commonDataSimpleContact->gapPtr);
2087
2088 auto t_w = getFTensor0IntegrationWeightSlave();
2089 const double cn_value = *cNPtr.get();
2090
2091 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
2092 const double val_m = SimpleContactProblem::ConstrainFunction_dg(
2093 cn_value, t_gap_gp, t_lagrange_slave) *
2094 t_w * area_slave;
2095
2096 auto t_base_lambda = row_data.getFTensor0N(gg, 0);
2097 for (int bbr = 0; bbr != nb_base_fun_row; ++bbr) {
2098
2099 auto t_base_master = col_data.getFTensor0N(gg, 0);
2101 &NN(bbr, 0), &NN(bbr, 1), &NN(bbr, 2)};
2102 const double m = val_m * t_base_lambda;
2103
2104 for (int bbc = 0; bbc != nb_base_fun_col; ++bbc) {
2105
2106 t_mat(i) += t_const_unit_n(i) * m * t_base_master;
2107
2108 ++t_base_master; // update rows
2109 ++t_mat;
2110 }
2111 ++t_base_lambda; // update cols master
2112 }
2113
2114 ++t_gap_gp;
2115 ++t_lagrange_slave;
2116 ++t_w;
2117 }
2118
2119 CHKERR MatSetValues(getSNESB(), row_data, col_data, &*NN.data().begin(),
2120 ADD_VALUES);
2121 }
2122
2124}
#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 double ConstrainFunction_dg(const double cn, const double g, const double l)

Member Data Documentation

◆ cNPtr

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

Definition at line 1499 of file SimpleContact.hpp.

◆ commonDataSimpleContact

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

Definition at line 1498 of file SimpleContact.hpp.

◆ NN

MatrixDouble SimpleContactProblem::OpCalDerIntCompFunOverSpatPosSlaveMaster::NN
private

Definition at line 1500 of file SimpleContact.hpp.


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