v0.13.0
Public Member Functions | Public Attributes | List of all members
NitscheMethod::OpRhsNormal Struct Reference

Calculate Nitsche method terms on right hand side. More...

#include <users_modules/homogenisation/src/NitscheMethod.hpp>

Inheritance diagram for NitscheMethod::OpRhsNormal:
[legend]
Collaboration diagram for NitscheMethod::OpRhsNormal:
[legend]

Public Member Functions

 OpRhsNormal (const std::string field_name, BlockData &nitsche_block_data, CommonData &nitsche_common_data, NonlinearElasticElement::BlockData &data, NonlinearElasticElement::CommonData &common_data, bool field_disp)
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, DataForcesAndSourcesCore::EntData &row_data)
 
- Public Member Functions inherited from NitscheMethod::OpCommon
 OpCommon (const std::string field_name, BlockData &nitsche_block_data, CommonData &nitsche_common_data, NonlinearElasticElement::BlockData &data, NonlinearElasticElement::CommonData &common_data, bool field_disp, const char type)
 
virtual MoFEMErrorCode calculateP (int gg, int fgg, int ff)
 
MoFEMErrorCode getJac (DataForcesAndSourcesCore::EntData &data, int gg, MatrixDouble &jac)
 
MoFEMErrorCode getTractionVariance (int gg, int fgg, int ff, MatrixDouble &jac, MatrixDouble &trac)
 
- Public Member Functions inherited from NitscheMethod::OpBasicCommon
 OpBasicCommon (const std::string field_name, BlockData &nitsche_block_data, CommonData &nitsche_common_data, bool field_disp, const char type)
 
virtual MoFEMErrorCode getFaceRadius (int ff)
 
virtual MoFEMErrorCode getGammaH (double gamma, int gg)
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
double & getVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoreBasegetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string &field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
const std::string & getFEName () const
 Get name of the element. More...
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
double & getMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Public Attributes

VectorDouble nF
 
- Public Attributes inherited from NitscheMethod::OpCommon
NonlinearElasticElement::BlockDatadAta
 
NonlinearElasticElement::CommonDatacommonData
 
VectorDouble dIsp
 
VectorDouble tRaction
 
MatrixDouble jAc_row
 
MatrixDouble jAc_col
 
MatrixDouble tRac_v
 
MatrixDouble tRac_u
 
- Public Attributes inherited from NitscheMethod::OpBasicCommon
BlockDatanitscheBlockData
 
CommonDatanitscheCommonData
 
bool fieldDisp
 
double faceRadius
 
double gammaH
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType { OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 }
 Controls loop over entities on element. More...
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Calculate Nitsche method terms on right hand side.

Definition at line 686 of file NitscheMethod.hpp.

Constructor & Destructor Documentation

◆ OpRhsNormal()

NitscheMethod::OpRhsNormal::OpRhsNormal ( const std::string  field_name,
BlockData nitsche_block_data,
CommonData nitsche_common_data,
NonlinearElasticElement::BlockData data,
NonlinearElasticElement::CommonData common_data,
bool  field_disp 
)

Definition at line 688 of file NitscheMethod.hpp.

693  : OpCommon(field_name, nitsche_block_data, nitsche_common_data, data,
694  common_data, field_disp, UserDataOperator::OPROW) {}
OpCommon(const std::string field_name, BlockData &nitsche_block_data, CommonData &nitsche_common_data, NonlinearElasticElement::BlockData &data, NonlinearElasticElement::CommonData &common_data, bool field_disp, const char type)

Member Function Documentation

◆ doWork()

MoFEMErrorCode NitscheMethod::OpRhsNormal::doWork ( int  row_side,
EntityType  row_type,
DataForcesAndSourcesCore::EntData row_data 
)

Definition at line 698 of file NitscheMethod.hpp.

699  {
701 
702  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
703  dAta.tEts.end()) {
705  }
706  if (row_data.getIndices().size() == 0)
708  int nb_dofs_row = row_data.getIndices().size();
709  double gamma = nitscheBlockData.gamma;
710  double phi = nitscheBlockData.phi;
711 
712  try {
713 
714  nF.resize(nb_dofs_row, false);
715  nF.clear();
716 
717  int gg = 0;
718  for (int ff = 0; ff < 4; ff++) {
719  if (!nitscheCommonData.facesFePtr[ff])
720  continue;
721  int nb_face_gauss_pts = nitscheCommonData.faceGaussPts[ff].size2();
722  ierr = getFaceRadius(ff);
723  CHKERRG(ierr);
724  for (int fgg = 0; fgg < nb_face_gauss_pts; fgg++, gg++) {
725 
726  ierr = getGammaH(gamma, gg);
727  CHKERRG(ierr);
728  double val = getGaussPts()(3, gg);
729  ierr = getJac(row_data, gg, jAc_row);
730  CHKERRG(ierr);
731  ierr = getTractionVariance(gg, fgg, ff, jAc_row, tRac_v);
732  CHKERRG(ierr);
734  3, ublas::shallow_array_adaptor<double>(
735  3, &nitscheCommonData.faceNormals[ff](fgg, 0)));
736  double area = cblas_dnrm2(3, &normal[0], 1);
737 
738  ierr = calculateP(gg, fgg, ff);
739  CHKERRG(ierr);
740  MatrixDouble &P = nitscheCommonData.P[ff][fgg];
741 
742  VectorDouble &u =
744  VectorDouble t = prod(trans(commonData.sTress[gg]),
745  normal); // this is Piola-Kirchhoff I stress
746 
747  for (int dd1 = 0; dd1 < nb_dofs_row / 3; dd1++) {
748  double n_row = row_data.getN()(gg, dd1);
749  for (int dd2 = 0; dd2 < 3; dd2++) {
750  nF[3 * dd1 + dd2] +=
751  gammaH *
752  (val * area * n_row *
753  (P(dd2, 0) * u[0] + P(dd2, 1) * u[1] + P(dd2, 2) * u[2]));
754  nF[3 * dd1 + dd2] +=
755  val * n_row *
756  (P(dd2, 0) * t[0] + P(dd2, 1) * t[1] + P(dd2, 2) * t[2]);
757  nF[3 * dd1 + dd2] +=
758  phi * tRac_v(dd2, 3 * dd1 + dd2) *
759  (P(dd2, 0) * u[0] + P(dd2, 1) * u[1] + P(dd2, 2) * u[2]);
760  }
761  }
762  }
763  }
764 
765  ierr = VecSetValues(getFEMethod()->snes_f, nb_dofs_row,
766  &row_data.getIndices()[0], &nF[0], ADD_VALUES);
767  CHKERRG(ierr);
768 
769  } catch (const std::exception &ex) {
770  std::ostringstream ss;
771  ss << "throw in method: " << ex.what() << std::endl;
772  SETERRQ(PETSC_COMM_SELF, 1, ss.str().c_str());
773  }
774 
776  }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
FTensor::Tensor1< double, SPACE_DIM > normal(FTensor::Tensor1< T1, 3 > &t_coords, FTensor::Tensor1< T2, SPACE_DIM > &t_disp)
Definition: ContactOps.hpp:77
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:87
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:126
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:88
UBlasVector< double > VectorDouble
Definition: Types.hpp:79
constexpr double t
plate stiffness
Definition: plate.cpp:76
static double phi
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
FTensor::Tensor2< double *, 3, 3 > & getJac()
get element Jacobian
double gamma
Penalty term, see .
double phi
Nitsche method parameter, see .
std::vector< boost::shared_ptr< const NumeredEntFiniteElement > > facesFePtr
std::vector< MatrixDouble > faceNormals
std::vector< std::vector< MatrixDouble > > P
projection matrix
std::vector< MatrixDouble > faceGaussPts
virtual MoFEMErrorCode getFaceRadius(int ff)
virtual MoFEMErrorCode getGammaH(double gamma, int gg)
NonlinearElasticElement::CommonData & commonData
virtual MoFEMErrorCode calculateP(int gg, int fgg, int ff)
NonlinearElasticElement::BlockData & dAta
MoFEMErrorCode getTractionVariance(int gg, int fgg, int ff, MatrixDouble &jac, MatrixDouble &trac)
Range tEts
constrains elements in block set
std::map< std::string, std::vector< VectorDouble > > dataAtGaussPts
std::vector< MatrixDouble3by3 > sTress

Member Data Documentation

◆ nF

VectorDouble NitscheMethod::OpRhsNormal::nF

Definition at line 696 of file NitscheMethod.hpp.


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