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

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

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

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

Public Member Functions

MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 Integrates Lagrange multipliers virtual work on master side, \( \delta W_{\text c}\) derivative with respect to material positions on master side and assembles components to LHS global matrix. More...
 
 OpContactTractionMasterMaster_dX (const string field_name, const string mesh_nodes_field, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, const int row_rank, const int col_rank)
 
- Public Member Functions inherited from SimpleContactProblem::OpContactALELhs
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)
 
 OpContactALELhs (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)
 

Additional Inherited Members

- Public Attributes inherited from SimpleContactProblem::OpContactALELhs
boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 
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 simple contact element.

Integrates Lagrange multipliers virtual work on master side, \( \delta W_{\text c}\) derivative with respect to material positions on master side and assembles components of the RHS vector.

Definition at line 3168 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpContactTractionMasterMaster_dX()

SimpleContactProblem::OpContactTractionMasterMaster_dX::OpContactTractionMasterMaster_dX ( const string  field_name,
const string  mesh_nodes_field,
boost::shared_ptr< CommonDataSimpleContact common_data_contact,
const int  row_rank,
const int  col_rank 
)
inline
Parameters
field_nameString of field name for spatial positions for rows
mesh_nodes_fieldString 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 3216 of file SimpleContact.hpp.

3220 : OpContactALELhs(field_name, mesh_nodes_field, common_data_contact,
3221 ContactOp::FACEMASTERMASTER, row_rank, col_rank) {
3222 sYmm = false; // This will make sure to loop over all intities (e.g.
3223 // for order=2 it will make doWork to loop 16 time)
3224 }
constexpr auto field_name
OpContactALELhs(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)

Member Function Documentation

◆ iNtegrate()

MoFEMErrorCode SimpleContactProblem::OpContactTractionMasterMaster_dX::iNtegrate ( EntData row_data,
EntData col_data 
)
virtual

Integrates Lagrange multipliers virtual work on master side, \( \delta W_{\text c}\) derivative with respect to material positions on master side and assembles components to LHS global matrix.

\[ \textrm{D} \delta W_{\rm{c}}({\mathbf{x}}^{(2)}, {\mathbf{X}}^{(2)}, \lambda, \delta{\mathbf{x}}^{(2)}) [\Delta{\mathbf{X}}^{(2)}] = -\int\limits_{\mathcal{T}^{(2)}_{\xi}} \lambda \, \frac{{\mathbf N}^{(1)}(\xi, \eta)}{||{{\mathbf N}^{(1)}(\xi, \eta)}||} \left[ \frac{\partial{\mathbf{X}}^{(2)}} {\partial\xi} \cdot \left(\frac{\partial\Delta {\mathbf{X}}^{(2)}}{\partial\eta}\times\delta{\mathbf{x}}^{(2)}\right) -\frac{\partial{\mathbf{X}}^{(2)}} {\partial\eta} \cdot \left(\frac{\partial\Delta {\mathbf{X}}^{(2)}}{\partial\xi}\times \delta{\mathbf{x}}^{(2)}\right)\right] \cdot \frac{{\mathbf N}^{(2)}(\xi, \eta)}{||{{\mathbf N}^{(2)}(\xi, \eta)}||} \textrm{d}\xi\textrm{d}\eta \]

Here superscript \((1)\) and \((2)\) denote slave and master side coordinates and surfaces, respectively. Moreover, \(\lambda\) is the lagrange multiplier.

Reimplemented from SimpleContactProblem::OpContactALELhs.

Definition at line 3621 of file SimpleContact.cpp.

3622 {
3624
3625 // Both sides are needed since both sides contribute their shape
3626 // function to the stiffness matrix
3627 const int nb_row = row_data.getIndices().size();
3628 if (!nb_row)
3630 const int nb_col = col_data.getIndices().size();
3631 if (!nb_col)
3633 const int nb_gauss_pts = row_data.getN().size1();
3634
3635 int nb_base_fun_row = row_data.getFieldData().size() / 3;
3636 int nb_base_fun_col = col_data.getFieldData().size() / 3;
3637
3638 auto get_tensor_from_mat = [](MatrixDouble &m, const int r, const int c) {
3640 &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 1, c + 0),
3641 &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 2, c + 0), &m(r + 2, c + 1),
3642 &m(r + 2, c + 2));
3643 };
3644
3645 auto get_tensor_vec = [](VectorDouble &n) {
3646 return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
3647 };
3648
3649 auto get_tensor_vec_3 = [&](VectorDouble3 &n) {
3650 return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
3651 };
3652
3656
3657 auto make_vec_der = [&](auto t_N, auto t_1, auto t_2) {
3659 t_n(i, j) = 0;
3660 t_n(i, j) += FTensor::levi_civita(i, j, k) * t_2(k) * t_N(0);
3661 t_n(i, j) -= FTensor::levi_civita(i, j, k) * t_1(k) * t_N(1);
3662 return t_n;
3663 };
3664
3665 matLhs.resize(3 * nb_base_fun_row, 3 * nb_base_fun_col, false);
3666 matLhs.clear();
3667
3668 auto lagrange_slave =
3669 getFTensor0FromVec(*commonDataSimpleContact->lagMultAtGaussPtsPtr);
3670
3671 auto t_1 =
3672 get_tensor_vec(*commonDataSimpleContact->tangentOneVectorMasterPtr);
3673 auto t_2 =
3674 get_tensor_vec(*commonDataSimpleContact->tangentTwoVectorMasterPtr);
3675
3676 auto t_w = getFTensor0IntegrationWeightMaster();
3677
3678 auto t_const_unit_master =
3679 get_tensor_vec(*(commonDataSimpleContact->normalVectorMasterPtr));
3680
3681 auto t_const_unit_slave =
3682 get_tensor_vec(*(commonDataSimpleContact->normalVectorSlavePtr));
3683
3684 const double area_m = commonDataSimpleContact->areaMaster;
3685 const double area_s = commonDataSimpleContact->areaSlave;
3686
3687 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
3688
3689 auto t_N = col_data.getFTensor1DiffN<2>(gg, 0);
3690
3691 for (int bbc = 0; bbc != nb_base_fun_col; ++bbc) {
3692 const double mult_m = 0.5 * t_w * lagrange_slave;
3693 FTensor::Tensor0<double *> t_base_master(&row_data.getN()(gg, 0));
3694
3695 for (int bbr = 0; bbr != nb_base_fun_row; ++bbr) {
3696 const double s = mult_m * t_base_master;
3697
3698 auto t_d_n = make_vec_der(t_N, t_1, t_2);
3699
3700 auto t_assemble_s = get_tensor_from_mat(matLhs, 3 * bbr, 3 * bbc);
3701
3702 t_assemble_s(i, j) -=
3703 s * t_d_n(k, j) * t_const_unit_master(k) * t_const_unit_slave(i);
3704
3705 ++t_base_master; // update rows
3706 }
3707 ++t_N;
3708 }
3709 ++lagrange_slave;
3710 ++t_w;
3711 }
3712
3714}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#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
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:92
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....
const VectorDouble & getFieldData() const
get dofs values
const VectorInt & getIndices() const
Get global indices of dofs on entity.
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact

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