v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | List of all members
SimpleContactProblem::OpContactMaterialSlaveSlaveLhs_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::OpContactMaterialSlaveSlaveLhs_dX_dLagmult:
[legend]
Collaboration diagram for SimpleContactProblem::OpContactMaterialSlaveSlaveLhs_dX_dLagmult:
[legend]

Public Member Functions

MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 Compute part of the left-hand side. More...
 
 OpContactMaterialSlaveSlaveLhs_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 slave side.

Definition at line 2944 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpContactMaterialSlaveSlaveLhs_dX_dLagmult()

SimpleContactProblem::OpContactMaterialSlaveSlaveLhs_dX_dLagmult::OpContactMaterialSlaveSlaveLhs_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 multipliers 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 q)

Definition at line 3003 of file SimpleContact.hpp.

3007 : OpContactMaterialLhs(mesh_nodes_field_row, lagrange_field_name,
3008 common_data_contact, ContactOp::FACESLAVESLAVE,
3009 row_rank, col_rank) {
3010 sYmm = false; // This will make sure to loop over all intities (e.g.
3011 // for order=2 it will make doWork to loop 16 time)
3012 }
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::OpContactMaterialSlaveSlaveLhs_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 slave side with respect to a variation of lagrange multipliers \((\Delta\lambda)\):

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

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

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

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 4299 of file SimpleContact.cpp.

4301 {
4302
4304
4307
4308 auto get_tensor_vec = [](VectorDouble &n) {
4309 return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
4310 };
4311
4312 auto get_tensor_from_mat = [](MatrixDouble &m, const int r, const int c) {
4313 return FTensor::Tensor1<double *, 3>(&m(r + 0, c + 0), &m(r + 1, c + 0),
4314 &m(r + 2, c + 0));
4315 };
4316
4317 auto t_F = getFTensor2FromMat<3, 3>(*commonDataSimpleContact->FMat);
4318
4319 auto t_w = getFTensor0IntegrationWeightSlave();
4320
4321 auto normal_master_at_gp =
4322 get_tensor_vec(*commonDataSimpleContact->normalVectorSlavePtr);
4323
4324 const double area_m = commonDataSimpleContact->areaSlave;
4325
4326 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
4327
4328 FTensor::Tensor0<double *> t_col_base(&col_data.getN()(gg, 0));
4329
4330 const double val = t_w * area_m;
4331
4332 int bbc = 0;
4333 for (; bbc != nb_base_fun_col; ++bbc) {
4334
4335 FTensor::Tensor0<double *> t_row_base(&row_data.getN()(gg, 0));
4336
4337 int bbr = 0;
4338 for (; bbr != nb_base_fun_row; ++bbr) {
4339
4340 auto t_assemble = get_tensor_from_mat(matLhs, 3 * bbr, bbc);
4341
4342 t_assemble(i) -=
4343 val * t_row_base * t_F(j, i) * normal_master_at_gp(j) * t_col_base;
4344
4345 ++t_row_base;
4346 }
4347 ++t_col_base;
4348 }
4349 ++t_F;
4350 ++t_w;
4351 }
4352
4354}
#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:8
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: