v0.14.0
Public Member Functions | Private Attributes | List of all members
SimpleContactProblem::OpLhsConvectIntegrationPtsConstrainMasterGap Struct Reference

Tangent opeerator for contrains equation for change of spatial positions on master and slave. More...

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

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

Public Member Functions

 OpLhsConvectIntegrationPtsConstrainMasterGap (const string lagrange_field_name, const string field_name, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, boost::shared_ptr< double > cn, const ContactOp::FaceType face_type, boost::shared_ptr< MatrixDouble > diff_convect)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 

Private Attributes

MatrixDouble matLhs
 
boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 
boost::shared_ptr< doublecNPtr
 
boost::shared_ptr< MatrixDouble > diffConvect
 

Detailed Description

Tangent opeerator for contrains equation for change of spatial positions on master and slave.

Definition at line 2427 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpLhsConvectIntegrationPtsConstrainMasterGap()

SimpleContactProblem::OpLhsConvectIntegrationPtsConstrainMasterGap::OpLhsConvectIntegrationPtsConstrainMasterGap ( const string  lagrange_field_name,
const string  field_name,
boost::shared_ptr< CommonDataSimpleContact common_data_contact,
boost::shared_ptr< double cn,
const ContactOp::FaceType  face_type,
boost::shared_ptr< MatrixDouble >  diff_convect 
)
inline

Definition at line 2429 of file SimpleContact.hpp.

2434  : ContactOp(lagrange_field_name, field_name, UserDataOperator::OPROWCOL,
2435  face_type),
2436  commonDataSimpleContact(common_data_contact), cNPtr(cn),
2437  diffConvect(diff_convect) {
2438  sYmm = false;
2439  }

Member Function Documentation

◆ doWork()

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

Definition at line 2987 of file SimpleContact.cpp.

2989  {
2991 
2992  const int nb_row = row_data.getIndices().size();
2993  const int nb_col = col_data.getIndices().size();
2994 
2995  if (nb_row && nb_col && col_type == MBVERTEX) {
2996 
2997  const int nb_gauss_pts = row_data.getN().size1();
2998  int nb_base_fun_row = nb_row;
2999  int nb_base_fun_col = nb_col / 3;
3000 
3001  const double area_slave =
3002  commonDataSimpleContact->areaSlave; // same area in master and slave
3003 
3004  matLhs.resize(nb_row, nb_col, false);
3005  matLhs.clear();
3006 
3007  auto get_diff_ksi = [](auto &m, auto gg) {
3009  &m(0, gg), &m(1, gg), &m(2, gg), &m(3, gg), &m(4, gg), &m(5, gg));
3010  };
3011 
3012  auto get_tensor_from_vec = [](VectorDouble &n) {
3013  return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
3014  };
3015 
3019 
3020  auto t_const_unit_n =
3021  get_tensor_from_vec(*(commonDataSimpleContact->normalVectorSlavePtr));
3022 
3023  auto t_lagrange_slave =
3024  getFTensor0FromVec(*commonDataSimpleContact->lagMultAtGaussPtsPtr);
3025  auto t_gap_gp = getFTensor0FromVec(*commonDataSimpleContact->gapPtr);
3026  auto &xi_grad_mat =
3027  *(commonDataSimpleContact->gradKsiPositionAtGaussPtsPtr);
3028  auto t_grad = getFTensor2FromMat<3, 2>(xi_grad_mat);
3029 
3030  auto t_w = getFTensor0IntegrationWeightSlave();
3031  const double cn_value = *cNPtr.get();
3032  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
3033  const double val_m = SimpleContactProblem::ConstrainFunction(
3034  cn_value, t_gap_gp, t_lagrange_slave) *
3035  t_w * area_slave;
3036  const double val_diff_m_l = SimpleContactProblem::ConstrainFunction_dl(
3037  cn_value, t_gap_gp, t_lagrange_slave) *
3038  t_w * area_slave;
3039  const double val_diff_m_g = SimpleContactProblem::ConstrainFunction_dg(
3040  cn_value, t_gap_gp, t_lagrange_slave) *
3041  t_w * area_slave;
3042 
3043  auto t_base_lambda = row_data.getFTensor0N(gg, 0);
3044  auto t_base_diff_lambda = row_data.getFTensor1DiffN<2>(gg, 0);
3045 
3047  &matLhs(0, 0), &matLhs(0, 1), &matLhs(0, 2)};
3048 
3049  for (int bbr = 0; bbr != nb_base_fun_row; ++bbr) {
3050 
3051  auto t_base_diff_disp = col_data.getFTensor1DiffN<2>(gg, 0);
3052  auto t_diff_convect = get_diff_ksi(*diffConvect, 3 * gg);
3053 
3054  for (int bbc = 0; bbc != nb_base_fun_col; ++bbc) {
3055 
3056  t_mat(i) += t_base_lambda * val_diff_m_g * t_const_unit_n(j) *
3057  t_grad(j, I) * t_diff_convect(I, i);
3058 
3059  ++t_base_diff_disp;
3060  ++t_diff_convect;
3061  ++t_mat;
3062  }
3063 
3064  ++t_base_lambda;
3065  ++t_base_diff_lambda;
3066  }
3067 
3068  ++t_gap_gp;
3069  ++t_lagrange_slave;
3070  ++t_grad;
3071  ++t_w;
3072  }
3073 
3074  CHKERR MatSetValues(getSNESB(), row_data, col_data, &*matLhs.data().begin(),
3075  ADD_VALUES);
3076  }
3078 }

Member Data Documentation

◆ cNPtr

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

Definition at line 2448 of file SimpleContact.hpp.

◆ commonDataSimpleContact

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

Definition at line 2447 of file SimpleContact.hpp.

◆ diffConvect

boost::shared_ptr<MatrixDouble> SimpleContactProblem::OpLhsConvectIntegrationPtsConstrainMasterGap::diffConvect
private

Definition at line 2449 of file SimpleContact.hpp.

◆ matLhs

MatrixDouble SimpleContactProblem::OpLhsConvectIntegrationPtsConstrainMasterGap::matLhs
private

Definition at line 2446 of file SimpleContact.hpp.


The documentation for this struct was generated from the following files:
SimpleContactProblem::ConstrainFunction_dg
static double ConstrainFunction_dg(const double cn, const double g, const double l)
Definition: SimpleContact.hpp:3595
MoFEM::MatSetValues
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
Definition: EntitiesFieldData.hpp:1631
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
SimpleContactProblem::OpLhsConvectIntegrationPtsConstrainMasterGap::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:2447
SimpleContactProblem::OpLhsConvectIntegrationPtsConstrainMasterGap::matLhs
MatrixDouble matLhs
Definition: SimpleContact.hpp:2446
SimpleContactProblem::ConstrainFunction_dl
static double ConstrainFunction_dl(const double cn, const double g, const double l)
Definition: SimpleContact.hpp:3601
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
FTensor::Tensor2
Definition: Tensor2_value.hpp:16
I
constexpr IntegrationType I
Definition: operators_tests.cpp:31
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1489
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1201
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
SimpleContactProblem::OpLhsConvectIntegrationPtsConstrainMasterGap::diffConvect
boost::shared_ptr< MatrixDouble > diffConvect
Definition: SimpleContact.hpp:2449
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', 3 >
convert.n
n
Definition: convert.py:82
SimpleContactProblem::OpLhsConvectIntegrationPtsConstrainMasterGap::cNPtr
boost::shared_ptr< double > cNPtr
Definition: SimpleContact.hpp:2448
SimpleContactProblem::ContactOp
ContactPrismElementForcesAndSourcesCore::UserDataOperator ContactOp
Definition: SimpleContact.hpp:30
MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
Definition: EntitiesFieldData.cpp:526
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1305
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
SimpleContactProblem::ConstrainFunction
static double ConstrainFunction(const double cn, const double g, const double l)
Definition: SimpleContact.hpp:3587
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346