v0.15.5
Loading...
Searching...
No Matches
Classes | Public Member Functions | Public Attributes | List of all members
SurfaceSlidingConstrains Struct Reference

Shape preserving constrains, i.e. nodes sliding on body surface. More...

#include "tools/src/SurfaceSlidingConstrains.hpp"

Inheritance diagram for SurfaceSlidingConstrains:
[legend]
Collaboration diagram for SurfaceSlidingConstrains:
[legend]

Classes

struct  DriverElementOrientation
 Class implemented by user to detect face orientation. More...
 
struct  MyTriangleFE
 
struct  OpJacobian
 

Public Member Functions

MyTriangleFEgetLoopFeRhs ()
 
MyTriangleFEgetLoopFeLhs ()
 
MoFEMErrorCode getOptions ()
 
 SurfaceSlidingConstrains (MoFEM::Interface &m_field, DriverElementOrientation &orientation)
 
virtual ~SurfaceSlidingConstrains ()=default
 
MoFEMErrorCode setOperators (int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name, const double *alpha=nullptr)
 
MoFEMErrorCode setOperatorsConstrainOnly (int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name)
 
MyTriangleFEgetLoopFeRhs ()
 
MyTriangleFEgetLoopFeLhs ()
 
MoFEMErrorCode getOptions ()
 
 SurfaceSlidingConstrains (MoFEM::Interface &m_field, DriverElementOrientation &orientation)
 
virtual ~SurfaceSlidingConstrains ()=default
 
MoFEMErrorCode setOperators (int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name, const double *alpha=nullptr)
 
MoFEMErrorCode setOperatorsConstrainOnly (int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name)
 
- Public Member Functions inherited from GenericSliding
 GenericSliding ()=default
 
virtual ~GenericSliding ()=default
 
 GenericSliding ()=default
 
virtual ~GenericSliding ()=default
 

Public Attributes

MoFEM::InterfacemField
 
boost::shared_ptr< MyTriangleFEfeRhsPtr
 
boost::shared_ptr< MyTriangleFEfeLhsPtr
 
MyTriangleFEfeRhs
 
MyTriangleFEfeLhs
 
DriverElementOrientationcrackFrontOrientation
 
double aLpha
 

Detailed Description

Shape preserving constrains, i.e. nodes sliding on body surface.

Derivation and implementation of constrains preserving body surface, i.e. body shape and volume.

The idea starts form observation that body shape can be globally characterized by constant calculated as volume over its area

\[ \frac{V}{A} = C \]

Above equation expressed in integral form is

\[ \int_\Omega \textrm{d}V = C \int_\Gamma \textrm{d}S \]

where notting that,

\[ \frac{1}{3} \int_\Omega \textrm{div}[\mathbf{X}] \textrm{d}\Omega = C \int_\Gamma \textrm{d}S \]

and applying Gauss theorem we get

\[ \int_\Gamma \mathbf{X}\cdot \frac{\mathbf{N}}{\|\mathbf{N}\|} \textrm{d}\Gamma = 3C \int_\Gamma \textrm{d}S. \]

Drooping integrals on both sides, and linearizing equation, we get

\[ \frac{\mathbf{N}}{\|\mathbf{N}\|} \cdot \delta \mathbf{X} = 3C - \frac{\mathbf{N}}{\|\mathbf{N}\|}\cdot \mathbf{X} \]

where \(\delta \mathbf{X}\) is displacement sub-inctrement. Above equation is a constrain if satisfied in body shape and volume is conserved. Final form of constrain equation is

\[ \mathcal{r} = \frac{\mathbf{N}}{\|\mathbf{N}\|}\cdot \mathbf{X} - \frac{\mathbf{N_0}}{\|\mathbf{N_0}\|}\cdot \mathbf{X}_0 = \frac{\mathbf{N}}{\|\mathbf{N}\|}\cdot (\mathbf{X}-\mathbf{X}_0) \]

In the spirit of finite element method the constrain equation is multiplied by shape functions and enforce using Lagrange multiplier method

\[ \int_\Gamma \mathbf{N}^\mathsf{T}_\lambda \left( \frac{\mathbf{N}}{\|\mathbf{N}\|}\mathbf{N}_\mathbf{X}\cdot (\overline{\mathbf{X}}-\overline{\mathbf{X}}_0) \right) \|\mathbf{N}\| \textrm{d}\Gamma = \mathbf{0}. \]

Above equation is nonlinear, applying to it Taylor expansion, we can get form which can be used with Newton interactive method

\[ \begin{split} &\int_\Gamma \mathbf{N}^\mathsf{T}_\lambda \left\{ \mathbf{N}\mathbf{N}_\mathbf{X} + \left(\mathbf{X}-\mathbf{X}_0\right) \cdot \left( \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta - \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi \right) \right\} \textrm{d}\Gamma \cdot \delta \overline{\mathbf{X}}\\ = &\int_\Gamma \mathbf{N}^\mathsf{T}_\lambda \mathbf{N}\cdot(\mathbf{X}-\mathbf{X}_0) \textrm{d}\Gamma \end{split}. \]

Equation expressing forces on shape as result of constrains, as result Lagrange multiplier method have following form

\[ \begin{split} &\int_\Gamma \mathbf{N}^\mathsf{T}_\mathbf{X} \cdot \mathbf{N} \mathbf{N}_\lambda \textrm{d}\Gamma \cdot \delta\overline{\lambda}\\ + &\int_\Gamma \lambda \mathbf{N}^\mathsf{T}_\mathbf{X} \left( \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta - \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi \right) \textrm{d}\Gamma \delta\overline{\mathbf{X}}\\ = &\int_\Gamma \lambda \mathbf{N}^\mathsf{T}_\mathbf{X} \cdot \mathbf{N} \textrm{d}\Gamma \end{split} \]

Above equations are assembled into global system of equations as following

\[ \left[ \begin{array}{cc} \mathbf{K} + \mathbf{B} & \mathbf{C}^\mathsf{T} \\ \mathbf{C} + \mathbf{A} & 0 \end{array} \right] \left\{ \begin{array}{c} \delta \overline{\mathbf{X}} \\ \delta \overline{\lambda} \end{array} \right\}= \left[ \begin{array}{c} \mathbf{f} - \mathbf{C}^\mathsf{T}\overline{\lambda} \\ \overline{\mathbf{r}} \end{array} \right] \]

where

\[ \mathbf{C}= \int_\Gamma \mathbf{N}_\lambda^\mathsf{T} \mathbf{N} \cdot \mathbf{N}_\mathbf{X} \textrm{d}\Gamma, \]

\[ \mathbf{B}= \int_\Gamma \lambda (\mathbf{X}-\mathbf{X}_0)\cdot \left( \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta - \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi \right) \textrm{d}\Gamma \]

and

\[ \mathbf{A}= \int_\Gamma \mathbf{N}^\mathsf{T}_\lambda \left(\mathbf{X}-\mathbf{X}_0\right) \cdot \left( \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta - \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi \right). \]

Derivation and implementation of constrains preserving body surface, i.e. body shape and volume.

The idea starts form observation that body shape can be globally characterized by constant calculated as volume over its area

\[ \frac{V}{A} = C \]

Above equation expressed in integral form is

\[ \int_\Omega \textrm{d}V = C \int_\Gamma \textrm{d}S \]

where notting that,

\[ \frac{1}{3} \int_\Omega \textrm{div}[\mathbf{X}] \textrm{d}\Omega = C \int_\Gamma \textrm{d}S \]

and applying Gauss theorem we get

\[ \int_\Gamma \mathbf{X}\cdot \frac{\mathbf{N}}{\|\mathbf{N}\|} \textrm{d}\Gamma = 3C \int_\Gamma \textrm{d}S. \]

Drooping integrals on both sides, and linearizing equation, we get

\[ \frac{\mathbf{N}}{\|\mathbf{N}\|} \cdot \delta \mathbf{X} = 3C - \frac{\mathbf{N}}{\|\mathbf{N}\|}\cdot \mathbf{X} \]

where \(\delta \mathbf{X}\) is displacement sub-inctrement. Above equation is a constrain if satisfied in body shape and volume is conserved. Final form of constrain equation is

\[ \mathcal{r} = \frac{\mathbf{N}}{\|\mathbf{N}\|}\cdot \mathbf{X} - \frac{\mathbf{N_0}}{\|\mathbf{N_0}\|}\cdot \mathbf{X}_0 = \frac{\mathbf{N}}{\|\mathbf{N}\|}\cdot (\mathbf{X}-\mathbf{X}_0) \]

In the spirit of finite element method the constrain equation is multiplied by shape functions and enforce using Lagrange multiplier method

\[ \int_\Gamma \mathbf{N}^\mathsf{T}_\lambda \left( \frac{\mathbf{N}}{\|\mathbf{N}\|}\mathbf{N}_\mathbf{X}\cdot (\overline{\mathbf{X}}-\overline{\mathbf{X}}_0) \right) \|\mathbf{N}\| \textrm{d}\Gamma = \mathbf{0}. \]

Above equation is nonlinear, applying to it Taylor expansion, we can get form which can be used with Newton interactive method

\[ \begin{split} &\int_\Gamma \mathbf{N}^\mathsf{T}_\lambda \left\{ \mathbf{N}\mathbf{N}_\mathbf{X} + \left(\mathbf{X}-\mathbf{X}_0\right) \cdot \left( \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta - \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi \right) \right\} \textrm{d}\Gamma \cdot \delta \overline{\mathbf{X}}\\ = &\int_\Gamma \mathbf{N}^\mathsf{T}_\lambda \mathbf{N}\cdot(\mathbf{X}-\mathbf{X}_0) \textrm{d}\Gamma \end{split}. \]

Equation expressing forces on shape as result of constrains, as result Lagrange multiplier method have following form

\[ \begin{split} &\int_\Gamma \mathbf{N}^\mathsf{T}_\mathbf{X} \cdot \mathbf{N} \mathbf{N}_\lambda \textrm{d}\Gamma \cdot \delta\overline{\lambda}\\ + &\int_\Gamma \lambda \mathbf{N}^\mathsf{T}_\mathbf{X} \left( \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta - \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi \right) \textrm{d}\Gamma \delta\overline{\mathbf{X}}\\ = &\int_\Gamma \lambda \mathbf{N}^\mathsf{T}_\mathbf{X} \cdot \mathbf{N} \textrm{d}\Gamma \end{split} \]

Above equations are assembled into global system of equations as following

\[ \left[ \begin{array}{cc} \mathbf{K} + \mathbf{B} & \mathbf{C}^\mathsf{T} \\ \mathbf{C} + \mathbf{A} & 0 \end{array} \right] \left\{ \begin{array}{c} \delta \overline{\mathbf{X}} \\ \delta \overline{\lambda} \end{array} \right\}= \left[ \begin{array}{c} \mathbf{f} - \mathbf{C}^\mathsf{T}\overline{\lambda} \\ \overline{\mathbf{r}} \end{array} \right] \]

where

\[ \mathbf{C}= \int_\Gamma \mathbf{N}_\lambda^\mathsf{T} \mathbf{N} \cdot \mathbf{N}_\mathbf{X} \textrm{d}\Gamma, \]

\[ \mathbf{B}= \int_\Gamma \lambda (\mathbf{X}-\mathbf{X}_0)\cdot \left( \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta - \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi \right) \textrm{d}\Gamma \]

and

\[ \mathbf{A}= \int_\Gamma \mathbf{N}^\mathsf{T}_\lambda \left(\mathbf{X}-\mathbf{X}_0\right) \cdot \left( \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta - \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi \right). \]

Definition at line 321 of file SurfaceSlidingConstrains.hpp.

Constructor & Destructor Documentation

◆ SurfaceSlidingConstrains() [1/2]

SurfaceSlidingConstrains::SurfaceSlidingConstrains ( MoFEM::Interface m_field,
DriverElementOrientation orientation 
)
inline

Definition at line 401 of file SurfaceSlidingConstrains.hpp.

403 : mField(m_field), feRhsPtr(new MyTriangleFE(m_field)),
404 feLhsPtr(new MyTriangleFE(m_field)), feRhs(*feRhsPtr), feLhs(*feLhsPtr),
405 crackFrontOrientation(orientation), aLpha(1) {
406
407 ierr = getOptions();
408 CHKERRABORT(PETSC_COMM_WORLD, ierr);
409 }
static PetscErrorCode ierr
boost::shared_ptr< MyTriangleFE > feRhsPtr
DriverElementOrientation & crackFrontOrientation
boost::shared_ptr< MyTriangleFE > feLhsPtr

◆ ~SurfaceSlidingConstrains() [1/2]

virtual SurfaceSlidingConstrains::~SurfaceSlidingConstrains ( )
virtualdefault

◆ SurfaceSlidingConstrains() [2/2]

SurfaceSlidingConstrains::SurfaceSlidingConstrains ( MoFEM::Interface m_field,
DriverElementOrientation orientation 
)
inline

Definition at line 403 of file SurfaceSlidingConstrains.hpp.

405 : mField(m_field), feRhsPtr(new MyTriangleFE(m_field)),
406 feLhsPtr(new MyTriangleFE(m_field)), feRhs(*feRhsPtr), feLhs(*feLhsPtr),
407 crackFrontOrientation(orientation), aLpha(1) {
408
409 ierr = getOptions();
410 CHKERRABORT(PETSC_COMM_WORLD, ierr);
411 }

◆ ~SurfaceSlidingConstrains() [2/2]

virtual SurfaceSlidingConstrains::~SurfaceSlidingConstrains ( )
virtualdefault

Member Function Documentation

◆ getLoopFeLhs() [1/2]

MyTriangleFE & SurfaceSlidingConstrains::getLoopFeLhs ( )
inline

Definition at line 385 of file SurfaceSlidingConstrains.hpp.

385{ return feLhs; }

◆ getLoopFeLhs() [2/2]

MyTriangleFE & SurfaceSlidingConstrains::getLoopFeLhs ( )
inline

Definition at line 387 of file SurfaceSlidingConstrains.hpp.

387{ return feLhs; }

◆ getLoopFeRhs() [1/2]

MyTriangleFE & SurfaceSlidingConstrains::getLoopFeRhs ( )
inline

Definition at line 383 of file SurfaceSlidingConstrains.hpp.

383{ return feRhs; }

◆ getLoopFeRhs() [2/2]

MyTriangleFE & SurfaceSlidingConstrains::getLoopFeRhs ( )
inline

Definition at line 385 of file SurfaceSlidingConstrains.hpp.

385{ return feRhs; }

◆ getOptions() [1/2]

MoFEMErrorCode SurfaceSlidingConstrains::getOptions ( )
inline

Definition at line 390 of file SurfaceSlidingConstrains.hpp.

390 {
392 PetscOptionsBegin(PETSC_COMM_WORLD, "",
393 "Get surface sliding constrains element scaling",
394 "none");
395 CHKERR PetscOptionsScalar("-surface_sliding_alpha", "scaling parameter", "",
396 aLpha, &aLpha, PETSC_NULLPTR);
397 PetscOptionsEnd();
399 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.

◆ getOptions() [2/2]

MoFEMErrorCode SurfaceSlidingConstrains::getOptions ( )
inline

Definition at line 392 of file SurfaceSlidingConstrains.hpp.

392 {
394 PetscOptionsBegin(PETSC_COMM_WORLD, "",
395 "Get surface sliding constrains element scaling",
396 "none");
397 CHKERR PetscOptionsScalar("-surface_sliding_alpha", "scaling parameter", "",
398 aLpha, &aLpha, PETSC_NULLPTR);
399 PetscOptionsEnd();
401 }

◆ setOperators() [1/2]

MoFEMErrorCode SurfaceSlidingConstrains::setOperators ( int  tag,
const std::string  lagrange_multipliers_field_name,
const std::string  material_field_name,
const double alpha = nullptr 
)
inline

Definition at line 554 of file SurfaceSlidingConstrains.hpp.

557 {
559
560 if (alpha != nullptr) {
561 aLpha = *alpha;
562 }
563
564 boost::shared_ptr<VectorDouble> active_variables_ptr(
565 new VectorDouble(3 + 9));
566 boost::shared_ptr<VectorDouble> results_ptr(new VectorDouble(3 + 9));
567 boost::shared_ptr<MatrixDouble> jacobian_ptr(
568 new MatrixDouble(3 + 9, 3 + 9));
569
570 feRhs.getOpPtrVector().clear();
571 feRhs.getOpPtrVector().push_back(new OpGetActiveDofsLambda(
572 lagrange_multipliers_field_name, active_variables_ptr));
573 feRhs.getOpPtrVector().push_back(new OpGetActiveDofsPositions<3>(
574 material_field_name, active_variables_ptr));
575 feRhs.getOpPtrVector().push_back(new OpJacobian(
576 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
577 jacobian_ptr, crackFrontOrientation, false, aLpha));
578 feRhs.getOpPtrVector().push_back(
579 new OpAssembleRhs<3, 9>(lagrange_multipliers_field_name, results_ptr));
580 feRhs.getOpPtrVector().push_back(
581 new OpAssembleRhs<3, 9>(material_field_name, results_ptr));
582
583 // Adding operators to calculate the left hand side
584 feLhs.getOpPtrVector().clear();
585 feLhs.getOpPtrVector().push_back(new OpGetActiveDofsLambda(
586 lagrange_multipliers_field_name, active_variables_ptr));
587 feLhs.getOpPtrVector().push_back(new OpGetActiveDofsPositions<3>(
588 material_field_name, active_variables_ptr));
589 feLhs.getOpPtrVector().push_back(new OpJacobian(
590 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
591 jacobian_ptr, crackFrontOrientation, true, aLpha));
592 feLhs.getOpPtrVector().push_back(new OpAssembleLhs<3, 9>(
593 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
594 feLhs.getOpPtrVector().push_back(new OpAssembleLhs<3, 9>(
595 material_field_name, lagrange_multipliers_field_name, jacobian_ptr));
596 feLhs.getOpPtrVector().push_back(new OpAssembleLhs<3, 9>(
597 material_field_name, material_field_name, jacobian_ptr));
598
600 }
UBlasMatrix< double > MatrixDouble
Definition Types.hpp:77
UBlasVector< double > VectorDouble
Definition Types.hpp:68
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.

◆ setOperators() [2/2]

MoFEMErrorCode SurfaceSlidingConstrains::setOperators ( int  tag,
const std::string  lagrange_multipliers_field_name,
const std::string  material_field_name,
const double alpha = nullptr 
)
inline

Definition at line 558 of file SurfaceSlidingConstrains.hpp.

561 {
563
564 if (alpha != nullptr) {
565 aLpha = *alpha;
566 }
567
568 boost::shared_ptr<VectorDouble> active_variables_ptr(
569 new VectorDouble(3 + 9));
570 boost::shared_ptr<VectorDouble> results_ptr(new VectorDouble(3 + 9));
571 boost::shared_ptr<MatrixDouble> jacobian_ptr(
572 new MatrixDouble(3 + 9, 3 + 9));
573
574 feRhs.getOpPtrVector().clear();
575 feRhs.getOpPtrVector().push_back(new OpGetActiveDofsLambda(
576 lagrange_multipliers_field_name, active_variables_ptr));
577 feRhs.getOpPtrVector().push_back(new OpGetActiveDofsPositions<3>(
578 material_field_name, active_variables_ptr));
579 feRhs.getOpPtrVector().push_back(new OpJacobian(
580 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
581 jacobian_ptr, crackFrontOrientation, false, aLpha));
582 feRhs.getOpPtrVector().push_back(
583 new OpAssembleRhs<3, 9>(lagrange_multipliers_field_name, results_ptr));
584 feRhs.getOpPtrVector().push_back(
585 new OpAssembleRhs<3, 9>(material_field_name, results_ptr));
586
587 // Adding operators to calculate the left hand side
588 feLhs.getOpPtrVector().clear();
589 feLhs.getOpPtrVector().push_back(new OpGetActiveDofsLambda(
590 lagrange_multipliers_field_name, active_variables_ptr));
591 feLhs.getOpPtrVector().push_back(new OpGetActiveDofsPositions<3>(
592 material_field_name, active_variables_ptr));
593 feLhs.getOpPtrVector().push_back(new OpJacobian(
594 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
595 jacobian_ptr, crackFrontOrientation, true, aLpha));
596 feLhs.getOpPtrVector().push_back(new OpAssembleLhs<3, 9>(
597 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
598 feLhs.getOpPtrVector().push_back(new OpAssembleLhs<3, 9>(
599 material_field_name, lagrange_multipliers_field_name, jacobian_ptr));
600 feLhs.getOpPtrVector().push_back(new OpAssembleLhs<3, 9>(
601 material_field_name, material_field_name, jacobian_ptr));
602
604 }

◆ setOperatorsConstrainOnly() [1/2]

MoFEMErrorCode SurfaceSlidingConstrains::setOperatorsConstrainOnly ( int  tag,
const std::string  lagrange_multipliers_field_name,
const std::string  material_field_name 
)
inline

Definition at line 603 of file SurfaceSlidingConstrains.hpp.

605 {
607
608 boost::shared_ptr<VectorDouble> active_variables_ptr(
609 new VectorDouble(3 + 9));
610 boost::shared_ptr<VectorDouble> results_ptr(new VectorDouble(3 + 9));
611 boost::shared_ptr<MatrixDouble> jacobian_ptr(
612 new MatrixDouble(3 + 9, 3 + 9));
613
614 // Adding operators to calculate the left hand side
615 feLhs.getOpPtrVector().clear();
616 feLhs.getOpPtrVector().push_back(new OpGetActiveDofsLambda(
617 lagrange_multipliers_field_name, active_variables_ptr));
618 feLhs.getOpPtrVector().push_back(new OpGetActiveDofsPositions<3>(
619 material_field_name, active_variables_ptr));
620 feLhs.getOpPtrVector().push_back(new OpJacobian(
621 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
622 jacobian_ptr, crackFrontOrientation, true, aLpha));
623 feLhs.getOpPtrVector().push_back(new OpAssembleLhs<3, 9>(
624 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
625
627 }

◆ setOperatorsConstrainOnly() [2/2]

MoFEMErrorCode SurfaceSlidingConstrains::setOperatorsConstrainOnly ( int  tag,
const std::string  lagrange_multipliers_field_name,
const std::string  material_field_name 
)
inline

Definition at line 607 of file SurfaceSlidingConstrains.hpp.

609 {
611
612 boost::shared_ptr<VectorDouble> active_variables_ptr(
613 new VectorDouble(3 + 9));
614 boost::shared_ptr<VectorDouble> results_ptr(new VectorDouble(3 + 9));
615 boost::shared_ptr<MatrixDouble> jacobian_ptr(
616 new MatrixDouble(3 + 9, 3 + 9));
617
618 // Adding operators to calculate the left hand side
619 feLhs.getOpPtrVector().clear();
620 feLhs.getOpPtrVector().push_back(new OpGetActiveDofsLambda(
621 lagrange_multipliers_field_name, active_variables_ptr));
622 feLhs.getOpPtrVector().push_back(new OpGetActiveDofsPositions<3>(
623 material_field_name, active_variables_ptr));
624 feLhs.getOpPtrVector().push_back(new OpJacobian(
625 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
626 jacobian_ptr, crackFrontOrientation, true, aLpha));
627 feLhs.getOpPtrVector().push_back(new OpAssembleLhs<3, 9>(
628 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
629
631 }

Member Data Documentation

◆ aLpha

double SurfaceSlidingConstrains::aLpha

Definition at line 389 of file SurfaceSlidingConstrains.hpp.

◆ crackFrontOrientation

DriverElementOrientation & SurfaceSlidingConstrains::crackFrontOrientation

Definition at line 387 of file SurfaceSlidingConstrains.hpp.

◆ feLhs

MyTriangleFE & SurfaceSlidingConstrains::feLhs

Definition at line 384 of file SurfaceSlidingConstrains.hpp.

◆ feLhsPtr

boost::shared_ptr< MyTriangleFE > SurfaceSlidingConstrains::feLhsPtr

Definition at line 380 of file SurfaceSlidingConstrains.hpp.

◆ feRhs

MyTriangleFE & SurfaceSlidingConstrains::feRhs

Definition at line 382 of file SurfaceSlidingConstrains.hpp.

◆ feRhsPtr

boost::shared_ptr< MyTriangleFE > SurfaceSlidingConstrains::feRhsPtr

Definition at line 380 of file SurfaceSlidingConstrains.hpp.

◆ mField

MoFEM::Interface & SurfaceSlidingConstrains::mField

Definition at line 351 of file SurfaceSlidingConstrains.hpp.


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