v0.14.0
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 <users_modules/basic_finite_elements/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)
 
- Public Member Functions inherited from GenericSliding
 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). \]

Examples
mesh_smoothing.cpp.

Definition at line 323 of file SurfaceSlidingConstrains.hpp.

Constructor & Destructor Documentation

◆ SurfaceSlidingConstrains()

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

Definition at line 405 of file SurfaceSlidingConstrains.hpp.

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

◆ ~SurfaceSlidingConstrains()

virtual SurfaceSlidingConstrains::~SurfaceSlidingConstrains ( )
virtualdefault

Member Function Documentation

◆ getLoopFeLhs()

MyTriangleFE & SurfaceSlidingConstrains::getLoopFeLhs ( )
inline

Definition at line 387 of file SurfaceSlidingConstrains.hpp.

387{ return feLhs; }

◆ getLoopFeRhs()

MyTriangleFE & SurfaceSlidingConstrains::getLoopFeRhs ( )
inline

Definition at line 385 of file SurfaceSlidingConstrains.hpp.

385{ return feRhs; }

◆ getOptions()

MoFEMErrorCode SurfaceSlidingConstrains::getOptions ( )
inline

Definition at line 392 of file SurfaceSlidingConstrains.hpp.

392 {
394 ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "",
395 "Get surface sliding constrains element scaling",
396 "none");
397 CHKERRQ(ierr);
398 CHKERR PetscOptionsScalar("-surface_sliding_alpha", "scaling parameter", "",
399 aLpha, &aLpha, PETSC_NULL);
400 ierr = PetscOptionsEnd();
401 CHKERRQ(ierr);
403 }
#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

◆ setOperators()

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 563 of file SurfaceSlidingConstrains.hpp.

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

◆ setOperatorsConstrainOnly()

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

Definition at line 612 of file SurfaceSlidingConstrains.hpp.

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

Member Data Documentation

◆ aLpha

double SurfaceSlidingConstrains::aLpha

Definition at line 391 of file SurfaceSlidingConstrains.hpp.

◆ crackFrontOrientation

DriverElementOrientation& SurfaceSlidingConstrains::crackFrontOrientation

Definition at line 389 of file SurfaceSlidingConstrains.hpp.

◆ feLhs

MyTriangleFE& SurfaceSlidingConstrains::feLhs

Definition at line 386 of file SurfaceSlidingConstrains.hpp.

◆ feLhsPtr

boost::shared_ptr<MyTriangleFE> SurfaceSlidingConstrains::feLhsPtr

Definition at line 382 of file SurfaceSlidingConstrains.hpp.

◆ feRhs

MyTriangleFE& SurfaceSlidingConstrains::feRhs

Definition at line 384 of file SurfaceSlidingConstrains.hpp.

◆ feRhsPtr

boost::shared_ptr<MyTriangleFE> SurfaceSlidingConstrains::feRhsPtr

Definition at line 382 of file SurfaceSlidingConstrains.hpp.

◆ mField

MoFEM::Interface& SurfaceSlidingConstrains::mField

Definition at line 353 of file SurfaceSlidingConstrains.hpp.


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