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

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

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

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

Public Member Functions

MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 Compute part of the left-hand side. More...
 
 OpContactMaterialSlaveOnFaceLhs_dX_dX (const string mesh_nodes_field_row, const string mesh_nodes_field_col, 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 normal vector from the expression for material traction contribution with respect to material coordinates on slave side.

Definition at line 2802 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpContactMaterialSlaveOnFaceLhs_dX_dX()

SimpleContactProblem::OpContactMaterialSlaveOnFaceLhs_dX_dX::OpContactMaterialSlaveOnFaceLhs_dX_dX ( const string  mesh_nodes_field_row,
const string  mesh_nodes_field_col,
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
mesh_nodes_field_colString of field name for material positions 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 3)

Definition at line 2844 of file SimpleContact.hpp.

2848 : OpContactMaterialLhs(mesh_nodes_field_row, mesh_nodes_field_col,
2849 common_data_contact, ContactOp::FACESLAVESLAVE,
2850 row_rank, col_rank) {
2851 sYmm = false; // This will make sure to loop over all intities (e.g.
2852 // for order=2 it will make doWork to loop 16 time)
2853 }
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::OpContactMaterialSlaveOnFaceLhs_dX_dX::iNtegrate ( EntData row_data,
EntData col_data 
)
virtual

Compute part of the left-hand side.

Computes the linearisation of the material component with respect to a variation of material coordinates \((\Delta{\mathbf{X}}^{(1)})\):

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

Here superscript \((1)\) denotes slave side coordinates and surfaces. Moreover, \(\lambda\) is the lagrange multiplier.

Reimplemented from SimpleContactProblem::OpContactMaterialLhs.

Definition at line 4173 of file SimpleContact.cpp.

4174 {
4175
4177
4181
4182 auto get_tensor2 = [](MatrixDouble &m, const int r, const int c) {
4184 &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 1, c + 0),
4185 &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 2, c + 0), &m(r + 2, c + 1),
4186 &m(r + 2, c + 2));
4187 };
4188
4189 auto get_tensor_vec = [](VectorDouble &n) {
4190 return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
4191 };
4192
4193 auto make_vec_der = [&](auto t_N, auto t_1, auto t_2) {
4195 t_n(i, j) = 0;
4196 t_n(i, j) += FTensor::levi_civita(i, j, k) * t_2(k) * t_N(0);
4197 t_n(i, j) -= FTensor::levi_civita(i, j, k) * t_1(k) * t_N(1);
4198 return t_n;
4199 };
4200
4201 // commonDataSimpleContact->faceRowData = nullptr;
4202
4203 auto t_F = getFTensor2FromMat<3, 3>(*commonDataSimpleContact->FMat);
4204
4205 auto t_w = getFTensor0IntegrationWeightMaster();
4206 auto lagrange_slave =
4207 getFTensor0FromVec(*commonDataSimpleContact->lagMultAtGaussPtsPtr);
4208 auto t_1 = get_tensor_vec(*commonDataSimpleContact->tangentOneVectorSlavePtr);
4209 auto t_2 = get_tensor_vec(*commonDataSimpleContact->tangentTwoVectorSlavePtr);
4210 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
4211
4212 auto t_N = col_data.getFTensor1DiffN<2>(gg, 0);
4213 const double val = 0.5 * t_w * lagrange_slave;
4214
4215 int bbc = 0;
4216 for (; bbc != nb_base_fun_col; ++bbc) {
4217
4218 FTensor::Tensor0<double *> t_base(&row_data.getN()(gg, 0));
4219
4220 int bbr = 0;
4221 for (; bbr != nb_base_fun_row; ++bbr) {
4222
4223 auto t_assemble = get_tensor2(matLhs, 3 * bbr, 3 * bbc);
4224 // TODO: handle hoGeometry
4225 auto t_d_n = make_vec_der(t_N, t_1, t_2);
4226 t_assemble(i, k) -= val * t_base * t_F(j, i) * t_d_n(j, k);
4227
4228 ++t_base;
4229 }
4230 ++t_N;
4231 }
4232 ++t_F;
4233 ++t_w;
4234 ++lagrange_slave;
4235 }
4236
4238}
#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
FTensor::Index< 'k', 3 > k
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
int r
Definition: sdf.py:5
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
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: