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

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

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

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

Public Member Functions

 OpCalContactTractionOverLambdaSlaveSlave (const string field_name, const string lagrange_field_name, boost::shared_ptr< CommonDataSimpleContact > common_data_contact)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 Integrates and assembles Lagrange multipliers virtual work, \( \delta W_{\text c}\) derivative with respect to Lagrange multipliers with respect to Lagrange multipliers on slave side and assembles components to LHS global matrix. More...
 

Private Attributes

boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 
MatrixDouble NN
 

Detailed Description

LHS-operator for the simple contact element.

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

Definition at line 1332 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpCalContactTractionOverLambdaSlaveSlave()

SimpleContactProblem::OpCalContactTractionOverLambdaSlaveSlave::OpCalContactTractionOverLambdaSlaveSlave ( const string  field_name,
const string  lagrange_field_name,
boost::shared_ptr< CommonDataSimpleContact common_data_contact 
)
inline

Definition at line 1334 of file SimpleContact.hpp.

1337  : ContactOp(field_name, lagrange_field_name, UserDataOperator::OPROWCOL,
1338  ContactOp::FACESLAVESLAVE),
1339  commonDataSimpleContact(common_data_contact) {
1340  sYmm = false; // This will make sure to loop over all intities (e.g.
1341  // for order=2 it will make doWork to loop 16 time)
1342  }

Member Function Documentation

◆ doWork()

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

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

Computes linearisation of integrated on slave side complementarity function and assembles Lagrange multipliers virtual work, \( \delta W_{\text c}\) with respect to Lagrange multipliers

\[ {\text D} {\delta W^{(1)}_{\text c}(\lambda, \delta \mathbf{x}^{(1)}})[\Delta \lambda] \,\, = \int_{{\gamma}^{(1)}_{\text c}} -\Delta \lambda \delta{\mathbf{x}^{(1)}} \,\,{ {\text d} {\gamma}} \]

where \({\gamma}^{(1)}_{\text c}\) is the surface integration domain of the slave surface, \( \lambda\) is the Lagrange multiplier, \(\mathbf{x}^{(1)}\) are the coordinates of the overlapping gauss points at slave triangles.

Definition at line 1927 of file SimpleContact.cpp.

1929  {
1931 
1932  // Both sides are needed since both sides contribute their shape
1933  // function to the stiffness matrix
1934  const int nb_row = row_data.getIndices().size();
1935  const int nb_col = col_data.getIndices().size();
1936 
1937  if (nb_row && nb_col) {
1938  const int nb_gauss_pts = row_data.getN().size1();
1939 
1940  int nb_base_fun_row = row_data.getFieldData().size() / 3;
1941  int nb_base_fun_col = col_data.getFieldData().size();
1942 
1943  const double area_slave =
1944  commonDataSimpleContact->areaSlave; // same area in master and slave
1945 
1946  auto get_tensor_from_mat = [](MatrixDouble &m, const int r, const int c) {
1947  return FTensor::Tensor1<double *, 3>(&m(r + 0, c), &m(r + 1, c),
1948  &m(r + 2, c));
1949  };
1950 
1951  auto get_tensor_vec = [](VectorDouble &n) {
1952  return FTensor::Tensor1<double *, 3>(&n(0), &n(1), &n(2));
1953  };
1954 
1956 
1957  NN.resize(3 * nb_base_fun_row, nb_base_fun_col, false);
1958  NN.clear();
1959 
1960  auto const_unit_n =
1961  get_tensor_vec(*(commonDataSimpleContact->normalVectorSlavePtr.get()));
1962 
1963  auto t_w = getFTensor0IntegrationWeightSlave();
1964 
1965  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1966 
1967  double val_m = t_w * area_slave;
1968  auto t_base_master = row_data.getFTensor0N(gg, 0);
1969 
1970  for (int bbr = 0; bbr != nb_base_fun_row; ++bbr) {
1971  auto t_assemble_m = get_tensor_from_mat(NN, 3 * bbr, 0);
1972  auto t_base_lambda = col_data.getFTensor0N(gg, 0);
1973  const double m = val_m * t_base_master;
1974  for (int bbc = 0; bbc != nb_base_fun_col; ++bbc) {
1975  const double n = m * t_base_lambda;
1976  t_assemble_m(i) += n * const_unit_n(i);
1977  ++t_assemble_m;
1978  ++t_base_lambda; // update cols slave
1979  }
1980  ++t_base_master; // update rows master
1981  }
1982  ++t_w;
1983  }
1984 
1985  CHKERR MatSetValues(getSNESB(), row_data, col_data, &*NN.data().begin(),
1986  ADD_VALUES);
1987  }
1989 }

Member Data Documentation

◆ commonDataSimpleContact

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

Definition at line 1373 of file SimpleContact.hpp.

◆ NN

MatrixDouble SimpleContactProblem::OpCalContactTractionOverLambdaSlaveSlave::NN
private

Definition at line 1374 of file SimpleContact.hpp.


The documentation for this struct was generated from the following files:
SimpleContactProblem::OpCalContactTractionOverLambdaSlaveSlave::NN
MatrixDouble NN
Definition: SimpleContact.hpp:1374
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
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1241
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
sdf.r
int r
Definition: sdf.py:8
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1489
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
SimpleContactProblem::OpCalContactTractionOverLambdaSlaveSlave::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:1373
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
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', 3 >
convert.n
n
Definition: convert.py:82
SimpleContactProblem::ContactOp
ContactPrismElementForcesAndSourcesCore::UserDataOperator ContactOp
Definition: SimpleContact.hpp:30
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
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
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346