v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | List of all members
SimpleContactProblem::OpContactMaterialMasterSlaveLhs_dX_dLagmult Struct Reference

LHS-operator for the contact element (material configuration) More...

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

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

Public Member Functions

MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 Compute part of the left-hand side. More...
 
 OpContactMaterialMasterSlaveLhs_dX_dLagmult (const string mesh_nodes_field_row, const string lagrange_field_name, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, const int row_rank, const int col_rank)
 
- Public Member Functions inherited from SimpleContactProblem::OpContactMaterialLhs
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 
virtual MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 
MoFEMErrorCode aSsemble (EntData &row_data, EntData &col_data)
 
 OpContactMaterialLhs (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, boost::shared_ptr< VolumeElementForcesAndSourcesCoreOnContactPrismSide > side_fe=NULL, const string side_fe_name="")
 

Additional Inherited Members

- Public Attributes inherited from SimpleContactProblem::OpContactMaterialLhs
boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 
boost::shared_ptr< VolumeElementForcesAndSourcesCoreOnContactPrismSide > sideFe
 
string sideFeName
 
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 contact element (material configuration)

Computes linearisation of the expression for material traction contribution with respect to material coordinates on master side.

Definition at line 2864 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpContactMaterialMasterSlaveLhs_dX_dLagmult()

SimpleContactProblem::OpContactMaterialMasterSlaveLhs_dX_dLagmult::OpContactMaterialMasterSlaveLhs_dX_dLagmult ( const string  mesh_nodes_field_row,
const string  lagrange_field_name,
boost::shared_ptr< CommonDataSimpleContact common_data_contact,
const int  row_rank,
const int  col_rank 
)
inline
Parameters
mesh_nodes_field_rowString of field name for material positions for rows
lagrange_field_nameString of field name for lagrange field name for columns
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 3)
col_rankParameter setting the dimension of the associated field for cols (in this case is 1)

Definition at line 2924 of file SimpleContact.hpp.

2928 : OpContactMaterialLhs(mesh_nodes_field_row, lagrange_field_name,
2929 common_data_contact, ContactOp::FACEMASTERSLAVE,
2930 row_rank, col_rank) {
2931 sYmm = false; // This will make sure to loop over all intities (e.g.
2932 // for order=2 it will make doWork to loop 16 time)
2933 }
OpContactMaterialLhs(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, boost::shared_ptr< VolumeElementForcesAndSourcesCoreOnContactPrismSide > side_fe=NULL, const string side_fe_name="")

Member Function Documentation

◆ iNtegrate()

MoFEMErrorCode SimpleContactProblem::OpContactMaterialMasterSlaveLhs_dX_dLagmult::iNtegrate ( EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)
virtual

Compute part of the left-hand side.

Computes the linearisation of the material component of contact tractions on master side with respect to a variation of lagrange multipliers \((\Delta\lambda)\):

\[ \delta W^\text{material}_p({\mathbf{x}}^{(2)}, {\mathbf{X}}^{(2)}, \delta{\mathbf{X}}^{(2)}, \lambda)[\Delta\lambda] = -\int\limits_{{\mathcal{T}}^{(2)}} \Delta\lambda \left\{{\left( \mathbf{F}^{(2)} \right)}^{\intercal}\cdot \mathbf{N}({\mathbf{X}}^{(2)}) \right\} \cdot \delta{\mathbf{X}}^{(2)}\, \textrm{d}{\mathcal{T}}^{(2)} = -\int\limits_{{\mathcal{T}}^{(2)}_{\xi}} \Delta\lambda \left\{{\left( \mathbf{F}^{(2)} \right)}^{\intercal}\cdot \left(\frac{\partial\mathbf{X}^{(2)}} {\partial\xi}\times\frac{\partial {\mathbf{X}}^{(2)}} {\partial\eta}\right) \right\} \cdot \delta{\mathbf{X}}^{(2)}\, \textrm{d}\xi\textrm{d}\eta \]

where \((2)\) denotes that the variable belongs to the master side, \( \lambda \) is contact traction on slave surface, \({\mathbf{N}}({\mathbf{X}}^{(2)})\) is a normal to the face in the material configuration, \(\xi, \eta\) are coordinates in the parent space \(({\mathcal{T}}^{(2)}_\xi)\) and \({\mathbf{F}}^{(2)}\) is the deformation gradient:

\[ {\mathbf{F}}^{(2)} = \mathbf{h}({\mathbf{x}}^{(2)})\,\mathbf{H}({\mathbf{X}}^{(2)})^{-1} = \frac{\partial{\mathbf{x}}^{(2)}}{\partial{\boldsymbol{\chi}}^{(2)}} \frac{\partial{\boldsymbol{\chi}}^{(2)}}{\partial{\mathbf{X}}^{(2)}} \]

where \(\mathbf{h}\) and \(\mathbf{H}\) are the gradients of the spatial and material maps, respectively, and \(\boldsymbol{\chi}\) are the reference coordinates.

Reimplemented from SimpleContactProblem::OpContactMaterialLhs.

Definition at line 4241 of file SimpleContact.cpp.

4243 {
4244
4246
4249
4250 auto get_tensor_vec = [](VectorDouble &n) {
4251 return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
4252 };
4253
4254 auto get_tensor_from_mat = [](MatrixDouble &m, const int r, const int c) {
4255 return FTensor::Tensor1<double *, 3>(&m(r + 0, c + 0), &m(r + 1, c + 0),
4256 &m(r + 2, c + 0));
4257 };
4258
4259 auto t_F = getFTensor2FromMat<3, 3>(*commonDataSimpleContact->FMat);
4260
4261 auto t_w = getFTensor0IntegrationWeightMaster();
4262
4263 auto normal_master_at_gp =
4264 get_tensor_vec(*commonDataSimpleContact->normalVectorMasterPtr);
4265
4266 const double area_m = commonDataSimpleContact->areaMaster;
4267
4268 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
4269
4270 FTensor::Tensor0<double *> t_col_base(&col_data.getN()(gg, 0));
4271
4272 const double val = t_w * area_m;
4273
4274 int bbc = 0;
4275 for (; bbc != nb_base_fun_col; ++bbc) {
4276
4277 FTensor::Tensor0<double *> t_row_base(&row_data.getN()(gg, 0));
4278
4279 int bbr = 0;
4280 for (; bbr != nb_base_fun_row; ++bbr) {
4281
4282 auto t_assemble = get_tensor_from_mat(matLhs, 3 * bbr, bbc);
4283
4284 t_assemble(i) -=
4285 val * t_row_base * t_F(j, i) * normal_master_at_gp(j) * t_col_base;
4286
4287 ++t_row_base;
4288 }
4289 ++t_col_base;
4290 }
4291 ++t_F;
4292 ++t_w;
4293 }
4294
4296}
#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
FTensor::Index< 'n', SPACE_DIM > n
FTensor::Index< 'm', SPACE_DIM > m
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
FTensor::Index< 'j', 3 > j
int r
Definition: sdf.py:5
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact

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