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

RHS-operator for the simple contact element for Augmented Lagrangian Method. More...

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

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

Public Member Functions

 OpGapConstraintAugmentedRhs (const string lagrange_field_name, boost::shared_ptr< CommonDataSimpleContact > common_data_contact, const double cn)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 

Private Attributes

boost::shared_ptr< CommonDataSimpleContactcommonDataSimpleContact
 Integrates KKT conditions for Augmented Lagrangian formulation at slave face gauss points and assembles components to the RHS global vector. More...
 
const double cN
 
VectorDouble vecR
 

Detailed Description

RHS-operator for the simple contact element for Augmented Lagrangian Method.

Integrates ALM constraints that fulfill KKT conditions over slave contact area and assembles components of the RHS vector.

Definition at line 1214 of file SimpleContact.hpp.

Constructor & Destructor Documentation

◆ OpGapConstraintAugmentedRhs()

SimpleContactProblem::OpGapConstraintAugmentedRhs::OpGapConstraintAugmentedRhs ( const string  lagrange_field_name,
boost::shared_ptr< CommonDataSimpleContact common_data_contact,
const double  cn 
)
inline

Definition at line 1217 of file SimpleContact.hpp.

1221 : ContactPrismElementForcesAndSourcesCore::UserDataOperator(
1222 lagrange_field_name, UserDataOperator::OPCOL,
1223 ContactOp::FACESLAVE),
1224 commonDataSimpleContact(common_data_contact), cN(cn) {}
double cn
Definition: contact.cpp:124
boost::shared_ptr< CommonDataSimpleContact > commonDataSimpleContact
Integrates KKT conditions for Augmented Lagrangian formulation at slave face gauss points and assembl...

Member Function Documentation

◆ doWork()

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

Definition at line 1685 of file SimpleContact.cpp.

1686 {
1688
1689 if (data.getIndices().size() == 0)
1691
1692 const int nb_gauss_pts = data.getN().size1();
1693
1694 int nb_base_fun_col = data.getFieldData().size();
1695 const double area_s =
1696 commonDataSimpleContact->areaSlave; // same area in master and slave
1697
1698 vecR.resize(nb_base_fun_col, false);
1699 vecR.clear();
1700
1701 auto t_lagrange_slave =
1702 getFTensor0FromVec(*commonDataSimpleContact->lagMultAtGaussPtsPtr);
1703 auto t_gap_gp = getFTensor0FromVec(*commonDataSimpleContact->gapPtr);
1704 auto t_w = getFTensor0IntegrationWeightSlave();
1705
1706 auto t_aug_lambda_ptr =
1707 getFTensor0FromVec(*commonDataSimpleContact->augmentedLambdasPtr);
1708
1709 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1710
1711 double branch_gg;
1712 if (t_aug_lambda_ptr > 0 && std::abs(t_aug_lambda_ptr) > ALM_TOL) {
1713 branch_gg = -t_lagrange_slave / cN;
1714 } else {
1715 branch_gg = t_gap_gp;
1716 }
1717
1718 const double val_s = t_w * area_s * branch_gg;
1719 auto t_base_lambda = data.getFTensor0N(gg, 0);
1720 for (int bbr = 0; bbr != nb_base_fun_col; ++bbr) {
1721 vecR[bbr] += val_s * t_base_lambda;
1722
1723 ++t_base_lambda; // update rows
1724 }
1725
1726 ++t_lagrange_slave;
1727 ++t_gap_gp;
1728 ++t_aug_lambda_ptr;
1729 ++t_w;
1730 } // for gauss points
1731
1732 CHKERR VecSetValues(getSNESf(), data, &vecR[0], ADD_VALUES);
1733
1735}
#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
#define CHKERR
Inline error check.
Definition: definitions.h:535
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
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.
static constexpr double ALM_TOL

Member Data Documentation

◆ cN

const double SimpleContactProblem::OpGapConstraintAugmentedRhs::cN
private

Definition at line 1265 of file SimpleContact.hpp.

◆ commonDataSimpleContact

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

Integrates KKT conditions for Augmented Lagrangian formulation at slave face gauss points and assembles components to the RHS global vector.

Integrates the Augmented Lagrangian multipliers formulation to fulfil KKT conditions in the integral sense and assembles components to the RHS global vector

\[ {\overline C(\lambda, \mathbf{x}^{(i)}, \delta \lambda)} = \left\{ \begin{array}{ll} \int_{{\gamma}^{(1)}_{\text c}} g_{\textrm{n}} \delta{{\lambda}} \,\,{ {\text d} {\gamma}} & \lambda + c_{\text n} g_{\textrm{n}}\leq 0 \\ \int_{{\gamma}^{(1)}_{\text c}} -\dfrac{1}{c_{\text n}} \lambda \delta{{\lambda}} \,\,{ {\text d} {\gamma}} & \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, \(\mathbf{x}^{(i)}\) are the coordinates of the overlapping gauss points at slave and master triangles for \(i = 1\) and \(i = 2\), respectively. Furthermore, \( c_{\text n}\) works as an augmentation parameter and affects convergence, \(r\) is regularisation parameter that can be chosen in \([1, 1.1]\) ( \(r = 1\)) is the default value) and \( g_{\textrm{n}}\) is the gap function evaluated at the slave triangle gauss points as:

\[ g_{\textrm{n}} = - \mathbf{n}(\mathbf{x}^{(1)}) \cdot \left( \mathbf{x}^{(1)} - \mathbf{x}^{(2)} \right) \]

Definition at line 1264 of file SimpleContact.hpp.

◆ vecR

VectorDouble SimpleContactProblem::OpGapConstraintAugmentedRhs::vecR
private

Definition at line 1266 of file SimpleContact.hpp.


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