v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Protected Attributes | List of all members
MoFEM::OpSetContravariantPiolaTransformOnFace2DImpl< 2 > Struct Reference

#include <src/finite_elements/UserDataOperators.hpp>

Inheritance diagram for MoFEM::OpSetContravariantPiolaTransformOnFace2DImpl< 2 >:
[legend]
Collaboration diagram for MoFEM::OpSetContravariantPiolaTransformOnFace2DImpl< 2 >:
[legend]

Public Member Functions

 OpSetContravariantPiolaTransformOnFace2DImpl (boost::shared_ptr< MatrixDouble > jac_ptr)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
- Public Member Functions inherited from MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
double getArea ()
 get area of face More...
 
double getMeasure ()
 get measure of element More...
 
VectorDoublegetNormal ()
 get triangle normal More...
 
VectorDoublegetTangent1 ()
 get triangle tangent 1 More...
 
VectorDoublegetTangent2 ()
 get triangle tangent 2 More...
 
auto getFTensor1Normal ()
 get normal as tensor More...
 
auto getFTensor1Tangent1 ()
 get tangentOne as tensor More...
 
auto getFTensor1Tangent2 ()
 get tangentTwo as tensor More...
 
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
VectorDoublegetCoords ()
 get triangle coordinates More...
 
auto getFTensor1Coords ()
 get get coords at gauss points More...
 
MatrixDoublegetNormalsAtGaussPts ()
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
MatrixDoublegetTangent1AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
auto getFTensor1NormalsAtGaussPts ()
 get normal at integration points More...
 
auto getFTensor1Tangent1AtGaussPts ()
 get tangent 1 at integration points More...
 
auto getFTensor1Tangent2AtGaussPts ()
 get tangent 2 at integration points More...
 
FaceElementForcesAndSourcesCoregetFaceFE ()
 return pointer to Generic Triangle Finite Element object More...
 
MoFEMErrorCode loopSideVolumes (const string fe_name, VolumeElementForcesAndSourcesCoreOnSide &fe_method)
 
- 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...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
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 getTStimeStep () 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...
 
doublegetMeasure ()
 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...
 

Protected Attributes

boost::shared_ptr< MatrixDoublejacPtr
 
MatrixDouble piolaN
 
MatrixDouble piolaDiffN
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 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)>
 
- 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...
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 

Detailed Description

Examples
continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, hcurl_check_approx_in_2d.cpp, mixed_poisson.cpp, phase.cpp, plot_base.cpp, and thermo_elastic.cpp.

Definition at line 3035 of file UserDataOperators.hpp.

Constructor & Destructor Documentation

◆ OpSetContravariantPiolaTransformOnFace2DImpl()

MoFEM::OpSetContravariantPiolaTransformOnFace2DImpl< 2 >::OpSetContravariantPiolaTransformOnFace2DImpl ( boost::shared_ptr< MatrixDouble jac_ptr)
inline

Definition at line 3038 of file UserDataOperators.hpp.

3040 : FaceElementForcesAndSourcesCore::UserDataOperator(HCURL),
3041 jacPtr(jac_ptr) {}
@ HCURL
field with continuous tangents
Definition: definitions.h:86

Member Function Documentation

◆ doWork()

Operator for linear form, usually to calculate values on right hand side.

Reimplemented from MoFEM::DataOperator.

Reimplemented in MoFEM::OpSetContravariantPiolaTransformOnFace2DImpl< 3 >.

Definition at line 394 of file UserDataOperators.cpp.

395 {
396
398
399 if (type != MBEDGE && type != MBTRI && type != MBQUAD)
401
405
406 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
407
409
410 const size_t nb_base_functions = data.getN(base).size2() / 3;
411 if (nb_base_functions) {
412
413 const size_t nb_gauss_pts = data.getN(base).size1();
414 piolaN.resize(nb_gauss_pts, data.getN(base).size2(), false);
415 if (data.getN(base).size2() > 0) {
416 auto t_n = data.getFTensor1N<3>(base);
417 double *t_transformed_n_ptr = &*piolaN.data().begin();
419 t_transformed_n_ptr, // HVEC0
420 &t_transformed_n_ptr[HVEC1], &t_transformed_n_ptr[HVEC2]);
421 auto t_jac = getFTensor2FromMat<2, 2>(*jacPtr);
422 for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg, ++t_jac) {
423 double det;
424 CHKERR determinantTensor2by2(t_jac, det);
425 for (unsigned int bb = 0; bb != nb_base_functions; ++bb) {
426 t_transformed_n(i) = t_jac(i, k) * t_n(k) / det;
427 ++t_n;
428 ++t_transformed_n;
429 }
430 }
431 data.getN(base).swap(piolaN);
432 }
433
434 piolaDiffN.resize(nb_gauss_pts, data.getDiffN(base).size2(), false);
435 if (data.getDiffN(base).data().size() > 0) {
436 auto t_diff_n = data.getFTensor2DiffN<3, 2>(base);
437 double *t_transformed_diff_n_ptr = &*piolaDiffN.data().begin();
439 t_transformed_diff_n(t_transformed_diff_n_ptr,
440 &t_transformed_diff_n_ptr[HVEC0_1],
441
442 &t_transformed_diff_n_ptr[HVEC1_0],
443 &t_transformed_diff_n_ptr[HVEC1_1],
444
445 &t_transformed_diff_n_ptr[HVEC2_0],
446 &t_transformed_diff_n_ptr[HVEC2_1]);
447 auto t_jac = getFTensor2FromMat<2, 2>(*jacPtr);
448 for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg, ++t_jac) {
449 double det;
450 CHKERR determinantTensor2by2(t_jac, det);
451 for (unsigned int bb = 0; bb != nb_base_functions; ++bb) {
452 t_transformed_diff_n(i, k) = t_jac(i, j) * t_diff_n(j, k) / det;
453 ++t_diff_n;
454 ++t_transformed_diff_n;
455 }
456 }
457 data.getDiffN(base).swap(piolaDiffN);
458 }
459 }
460 }
461
463}
FieldApproximationBase
approximation base
Definition: definitions.h:58
@ LASTBASE
Definition: definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
#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
@ HVEC1
Definition: definitions.h:186
@ HVEC2
Definition: definitions.h:186
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
@ HVEC1_1
Definition: definitions.h:196
@ HVEC0_1
Definition: definitions.h:195
@ HVEC1_0
Definition: definitions.h:193
@ HVEC2_1
Definition: definitions.h:197
@ HVEC2_0
Definition: definitions.h:194
#define CHKERR
Inline error check.
Definition: definitions.h:535
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
MoFEMErrorCode determinantTensor2by2(T1 &t, T2 &det)
Calculate determinant 2 by 2.
Definition: Templates.hpp:1393

Member Data Documentation

◆ jacPtr

boost::shared_ptr<MatrixDouble> MoFEM::OpSetContravariantPiolaTransformOnFace2DImpl< 2 >::jacPtr
protected

Definition at line 3047 of file UserDataOperators.hpp.

◆ piolaDiffN

Definition at line 3049 of file UserDataOperators.hpp.

◆ piolaN

Definition at line 3048 of file UserDataOperators.hpp.


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