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

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

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

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

Public Member Functions

 OpCalAugmentedTractionRhsMaster (const string field_name, boost::shared_ptr< CommonDataSimpleContact > common_data_contact)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Integrates Lagrange multipliers virtual work on slave surface and assembles components to the RHS global vector. More...
 

Private Attributes

boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 
VectorDouble vecF
 

Detailed Description

RHS-operator for the simple contact element.

Integrates Augmented Lagrange multipliers virtual work on master surface and assemble components to RHS global vector.

Definition at line 1051 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpCalAugmentedTractionRhsMaster()

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

Member Function Documentation

◆ doWork()

MoFEMErrorCode SimpleContactProblem::OpCalAugmentedTractionRhsMaster::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)

Integrates Lagrange multipliers virtual work on slave surface and assembles components to the RHS global vector.

Integrates Lagrange multipliers virtual work \( \delta W_{\text c}\) on slave surface and assembles components to the RHS global vector

\[ {\delta W^{(2)}_{\text c}(\lambda, \delta \mathbf{x}^{(2)}}) \,\, = \left\{ \begin{array}{ll} \int_{{\gamma}^{(1)}_{\text c}} (\lambda + c_{\textrm n} g_{\textrm{n}}) \mathbf{n}(\mathbf{x}^{(1)}) \cdot \delta{\mathbf{x}^{(2)}} \,\,{ {\text d} {\gamma}} & \lambda + c_{\text n} g_{\textrm{n}}\leq 0 \\ 0 & \lambda + c_{\text n} g_{\textrm{n}} > 0 \\ \end{array} \right. \]

where \({\gamma}^{(1)}_{\text c}\) is the surface integration domain of the slave surface, \( \lambda\) is the Lagrange multiplier, \( c_{\textrm n}\) is the regularisation/augmentation parameter of stress dimensions, \( g_{\textrm{n}}\) is the value of gap at the associated gauss point, \(\mathbf{x}^{(2)}\) are the coordinates of the overlapping gauss points at master triangles.

Definition at line 1862 of file SimpleContact.cpp.

1863  {
1865 
1866  if (data.getIndices().size() == 0)
1868 
1869  const int nb_dofs = data.getIndices().size();
1870  if (nb_dofs) {
1871 
1872  const int nb_gauss_pts = data.getN().size1();
1873  const int nb_base_fun_col = data.getFieldData().size() / 3;
1874 
1875  vecF.resize(nb_dofs,
1876  false); // the last false in ublas
1877  // resize will destroy (not
1878  // preserved) the old
1879  // values
1880  vecF.clear();
1881 
1882  const double area_m =
1883  commonDataSimpleContact->areaSlave; // same area in master and slave
1884 
1885  auto get_tensor_vec = [](VectorDouble &n, const int r) {
1886  return FTensor::Tensor1<double *, 3>(&n(r + 0), &n(r + 1), &n(r + 2));
1887  };
1888 
1890 
1891  auto const_unit_n = get_tensor_vec(
1892  commonDataSimpleContact->normalVectorSlavePtr.get()[0], 0);
1893 
1894  auto t_w = getFTensor0IntegrationWeightSlave();
1895  auto t_aug_lambda_ptr =
1896  getFTensor0FromVec(*commonDataSimpleContact->augmentedLambdasPtr);
1897 
1898  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1899  if (t_aug_lambda_ptr > 0 && std::abs(t_aug_lambda_ptr) > ALM_TOL) {
1900  ++t_aug_lambda_ptr;
1901  ++t_w;
1902  continue;
1903  }
1904  const double val_m = t_w * area_m;
1905 
1906  FTensor::Tensor0<double *> t_base_master(&data.getN()(gg, 0));
1907 
1908  for (int bbc = 0; bbc != nb_base_fun_col; ++bbc) {
1909 
1910  auto t_assemble_m = get_tensor_vec(vecF, 3 * bbc);
1911 
1912  t_assemble_m(i) +=
1913  val_m * const_unit_n(i) * t_aug_lambda_ptr * t_base_master;
1914 
1915  ++t_base_master;
1916  }
1917  ++t_aug_lambda_ptr;
1918  ++t_w;
1919  } // for gauss points
1920 
1921  CHKERR VecSetValues(getSNESf(), data, &*vecF.begin(), ADD_VALUES);
1922  }
1924 }

Member Data Documentation

◆ commonDataSimpleContact

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

Definition at line 1093 of file SimpleContact.hpp.

◆ vecF

VectorDouble SimpleContactProblem::OpCalAugmentedTractionRhsMaster::vecF
private

Definition at line 1094 of file SimpleContact.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1241
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1576
sdf.r
int r
Definition: sdf.py:8
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPCOL
@ OPCOL
operator doWork function is executed on FE columns
Definition: ForcesAndSourcesCore.hpp:568
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
SimpleContactProblem::ALM_TOL
static constexpr double ALM_TOL
Definition: SimpleContact.hpp:65
SimpleContactProblem::OpCalAugmentedTractionRhsMaster::commonDataSimpleContact
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Definition: SimpleContact.hpp:1093
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
FTensor::Tensor0
Definition: Tensor0.hpp:16
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
SimpleContactProblem::OpCalAugmentedTractionRhsMaster::vecF
VectorDouble vecF
Definition: SimpleContact.hpp:1094
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
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