v0.15.5
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Protected Attributes | List of all members
MoFEM::ProjectionFieldOn10NodeTet Struct Reference

#include "src/approximation/Projection10NodeCoordsOnField.hpp"

Inheritance diagram for MoFEM::ProjectionFieldOn10NodeTet:
[legend]
Collaboration diagram for MoFEM::ProjectionFieldOn10NodeTet:
[legend]

Public Member Functions

 ProjectionFieldOn10NodeTet (Interface &m_field, std::string _fieldName, bool set_nodes, bool on_coords, std::string on_tag="NoNE")
 
MoFEMErrorCode preProcess ()
 Pre-processing function executed at loop initialization.
 
MoFEMErrorCode operator() ()
 Main operator function executed for each loop iteration.
 
- Public Member Functions inherited from MoFEM::Projection10NodeCoordsOnField
 Projection10NodeCoordsOnField (Interface &m_field, std::string field_name, int verb=0)
 
MoFEMErrorCode postProcess ()
 Post-processing function executed at loop completion.
 
- Public Member Functions inherited from MoFEM::DofMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 DofMethod ()=default
 Default constructor.
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 Default constructor.
 
virtual ~BasicMethod ()=default
 Virtual destructor.
 
int getNinTheLoop () const
 Get current loop iteration index.
 
int getLoopSize () const
 Get total loop size.
 
auto getLoHiFERank () const
 Get processor rank range for finite element iteration.
 
auto getLoFERank () const
 Get lower processor rank for finite element iteration.
 
auto getHiFERank () const
 Get upper processor rank for finite element iteration.
 
unsigned int getFieldBitNumber (std::string field_name) const
 Get bit number for a specific field by name.
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from another BasicMethod instance.
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak pointer object.
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 KspMethod ()
 Default constructor.
 
virtual ~KspMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 Copy data from another KSP method.
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 Default constructor.
 
virtual ~PetscData ()=default
 Virtual destructor.
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 Copy PETSc data from another instance.
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 SnesMethod ()
 Default constructor.
 
virtual ~SnesMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy SNES data from another instance.
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TSMethod ()
 Default constructor.
 
virtual ~TSMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data from another instance.
 
- Public Member Functions inherited from MoFEM::TaoMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TaoMethod ()
 Default constructor.
 
virtual ~TaoMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTao (const TaoMethod &tao)
 Copy TAO data from another instance.
 

Public Attributes

bool setNodes
 
bool onCoords
 
std::string onTag
 
const int maxApproximationOrder
 
Tag th
 
- Public Attributes inherited from MoFEM::DofMethod
boost::shared_ptr< FieldfieldPtr
 Shared pointer to field information.
 
boost::shared_ptr< DofEntitydofPtr
 Shared pointer to DOF entity data.
 
boost::shared_ptr< NumeredDofEntitydofNumeredPtr
 Shared pointer to numbered DOF entity data.
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 Current index of processed method in the loop.
 
int loopSize
 Total number of methods to process in the loop.
 
std::pair< int, int > loHiFERank
 Processor rank range for distributed finite element iteration.
 
int rAnk
 Current processor rank in MPI communicator.
 
int sIze
 Total number of processors in MPI communicator.
 
const RefEntity_multiIndexrefinedEntitiesPtr
 Pointer to container of refined MoFEM DOF entities.
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 Pointer to container of refined finite element entities.
 
const ProblemproblemPtr
 Raw pointer to current MoFEM problem instance.
 
const Field_multiIndexfieldsPtr
 Raw pointer to fields multi-index container.
 
const FieldEntity_multiIndexentitiesPtr
 Raw pointer to container of field entities.
 
const DofEntity_multiIndexdofsPtr
 Raw pointer to container of degree of freedom entities.
 
const FiniteElement_multiIndexfiniteElementsPtr
 Raw pointer to container of finite elements.
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 Raw pointer to container of finite element entities.
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 Raw pointer to container of adjacencies between DOFs and finite elements.
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing operations.
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing operations.
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for main operator execution.
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 Switch for vector assembly operations.
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 Switch for matrix assembly operations.
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 Weak pointer to cached entity data.
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Current KSP computation context.
 
KSP ksp
 PETSc KSP linear solver object.
 
Vec & ksp_f
 Reference to residual vector in KSP context.
 
Mat & ksp_A
 Reference to system matrix in KSP context.
 
Mat & ksp_B
 Reference to preconditioner matrix in KSP context.
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 Current data context switches.
 
Vec f
 PETSc residual vector.
 
Mat A
 PETSc Jacobian matrix.
 
Mat B
 PETSc preconditioner matrix.
 
Vec x
 PETSc solution vector.
 
Vec dx
 PETSc solution increment vector.
 
Vec x_t
 PETSc first time derivative vector.
 
Vec x_tt
 PETSc second time derivative vector.
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 Current SNES computation context.
 
SNES snes
 PETSc SNES nonlinear solver object.
 
Vec & snes_x
 Reference to current solution state vector.
 
Vec & snes_dx
 Reference to solution update/increment vector.
 
Vec & snes_f
 Reference to residual vector.
 
Mat & snes_A
 Reference to Jacobian matrix.
 
Mat & snes_B
 Reference to preconditioner of Jacobian matrix.
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 PETSc time stepping solver object.
 
TSContext ts_ctx
 Current TS computation context.
 
PetscInt ts_step
 Current time step number.
 
PetscReal ts_a
 Shift parameter for U_t (see PETSc Time Solver documentation)
 
PetscReal ts_aa
 Shift parameter for U_tt (second time derivative)
 
PetscReal ts_t
 Current time value.
 
PetscReal ts_dt
 Current time step size.
 
Vec & ts_u
 Reference to current state vector U(t)
 
Vec & ts_u_t
 Reference to first time derivative of state vector dU/dt.
 
Vec & ts_u_tt
 Reference to second time derivative of state vector d²U/dt²
 
Vec & ts_F
 Reference to residual vector F(t,U,U_t,U_tt)
 
Mat & ts_A
 Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt.
 
Mat & ts_B
 Reference to preconditioner matrix for ts_A.
 
- Public Attributes inherited from MoFEM::TaoMethod
TAOContext tao_ctx
 Current TAO computation context.
 
Tao tao
 PETSc TAO optimization solver object.
 
Vec & tao_x
 Reference to optimization variables vector.
 
Vec & tao_f
 Reference to gradient vector.
 
Mat & tao_A
 Reference to Hessian matrix.
 
Mat & tao_B
 Reference to preconditioner matrix for Hessian.
 

Protected Attributes

Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
 
VectorDouble L
 
VectorDouble K
 
- Protected Attributes inherited from MoFEM::Projection10NodeCoordsOnField
InterfacemField
 
std::string fieldName
 
int vErbose
 
VectorDouble coords
 
VectorDouble3 aveMidCoord
 
VectorDouble3 midNodeCoord
 
VectorDouble3 diffNodeCoord
 
VectorDouble3 dOf
 

Additional Inherited Members

- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 Context enumeration for KSP solver phases. More...
 
- Public Types inherited from MoFEM::PetscData
enum  DataContext {
  CTX_SET_NONE = 0 , CTX_SET_F = 1 << 0 , CTX_SET_A = 1 << 1 , CTX_SET_B = 1 << 2 ,
  CTX_SET_X = 1 << 3 , CTX_SET_DX = 1 << 4 , CTX_SET_X_T = 1 << 5 , CTX_SET_X_TT = 1 << 6 ,
  CTX_SET_TIME = 1 << 7
}
 Enumeration for data context flags. More...
 
using Switches = std::bitset< 8 >
 Bitset type for context switches.
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 Context enumeration for SNES solver phases. More...
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 Context enumeration for TS solver phases. More...
 
- Public Types inherited from MoFEM::TaoMethod
enum  TAOContext { CTX_TAO_OBJECTIVE , CTX_TAO_GRADIENT , CTX_TAO_HESSIAN , CTX_TAO_NONE }
 Context enumeration for TAO solver phases. More...
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 No data switch.
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 Residual vector switch.
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 Jacobian matrix switch.
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 Preconditioner matrix switch.
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 Solution vector switch.
 
static constexpr Switches CtxSetDX = PetscData::Switches(CTX_SET_DX)
 Solution increment switch.
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 First time derivative switch.
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 Second time derivative switch.
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 Time value switch.
 

Detailed Description

Examples
mofem/tutorials/cor-10_navier_stokes/navier_stokes.cpp.

Definition at line 95 of file Projection10NodeCoordsOnField.hpp.

Constructor & Destructor Documentation

◆ ProjectionFieldOn10NodeTet()

MoFEM::ProjectionFieldOn10NodeTet::ProjectionFieldOn10NodeTet ( Interface m_field,
std::string  _fieldName,
bool  set_nodes,
bool  on_coords,
std::string  on_tag = "NoNE" 
)

Member Function Documentation

◆ operator()()

MoFEMErrorCode MoFEM::ProjectionFieldOn10NodeTet::operator() ( )
virtual

Main operator function executed for each loop iteration.

This virtual function is called for every item (finite element, entity, etc.) in the loop. It is the core computation function typically used for:

  • Calculating element local matrices and vectors
  • Assembling contributions to global system
  • Element-level post-processing operations
Returns
Error code indicating success or failure

Reimplemented from MoFEM::Projection10NodeCoordsOnField.

Definition at line 502 of file Projection10NodeCoordsOnField.cpp.

502 {
504 if (dofPtr->getName() != fieldName)
506 if (setNodes) {
507 if (dofPtr->getEntType() == MBVERTEX) {
508 EntityHandle node = dofPtr->getEnt();
509 if (onCoords) {
510 coords.resize(3);
511 CHKERR mField.get_moab().get_coords(&node, 1, &*coords.data().begin());
512 coords[dofPtr->getDofCoeffIdx()] = dofPtr->getFieldData();
513 CHKERR mField.get_moab().set_coords(&node, 1, &*coords.data().begin());
514 } else {
515 int field_rank = (*field_it)->getNbOfCoeffs();
516 if (field_rank != dofPtr->getNbOfCoeffs()) {
517 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
518 "data inconsistency");
519 }
520 double *tag_value;
521 int tag_size;
522 CHKERR mField.get_moab().tag_get_by_ptr(
523 th, &node, 1, (const void **)&tag_value, &tag_size);
524 if (tag_size != dofPtr->getNbOfCoeffs()) {
525 SETERRQ(PETSC_COMM_SELF, 1, "data inconsistency");
526 }
527 tag_value[dofPtr->getDofCoeffIdx()] = dofPtr->getFieldData();
528 }
529 }
531 }
532 if (dofPtr->getEntType() != MBEDGE) {
534 }
535 EntityHandle edge = dofPtr->getEnt();
536 if (type_from_handle(edge) != MBEDGE) {
537 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
538 "this method works only elements which are type of MBEDGE");
539 }
540
541 int num_nodes;
542 const EntityHandle *conn;
543 CHKERR mField.get_moab().get_connectivity(edge, conn, num_nodes, false);
544 if ((num_nodes != 2) && (num_nodes != 3)) {
545 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
546 "this method works only 4 node and 10 node tets");
547 }
548 if (num_nodes == 2) {
550 }
551
552 if (dofPtr->getDofOrder() >= maxApproximationOrder) {
553 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
554 "too big approximation order, increase constant "
555 "max_ApproximationOrder");
556 }
557 double approx_val = 0;
558 FieldApproximationBase base = dofPtr->getApproxBase();
559 switch (base) {
561 approx_val = 0.25 * L[dofPtr->getDofOrder() - 2] * dofPtr->getFieldData();
562 break;
564 approx_val = 0.25 * K[dofPtr->getDofOrder() - 2] * dofPtr->getFieldData();
565 break;
566 default:
567 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not yet implemented");
568 }
569
570 if (onCoords) {
571 coords.resize(num_nodes * 3);
572 CHKERR mField.get_moab().get_coords(conn, num_nodes,
573 &*coords.data().begin());
574 if (dofPtr->getEntDofIdx() == dofPtr->getDofCoeffIdx()) {
575 // add only one when higher order terms present
576 double ave_mid = (coords[3 * 0 + dofPtr->getDofCoeffIdx()] +
577 coords[3 * 1 + dofPtr->getDofCoeffIdx()]) *
578 0.5;
579 coords[2 * 3 + dofPtr->getDofCoeffIdx()] = ave_mid;
580 }
581 coords[2 * 3 + dofPtr->getDofCoeffIdx()] += approx_val;
582 CHKERR mField.get_moab().set_coords(&conn[2], 1, &coords[3 * 2]);
583 } else {
584 int tag_size;
585 double *tag_value[num_nodes];
586 CHKERR mField.get_moab().tag_get_by_ptr(
587 th, conn, num_nodes, (const void **)tag_value, &tag_size);
588 if (tag_size != dofPtr->getNbOfCoeffs()) {
589 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
590 }
591 if (dofPtr->getEntDofIdx() == dofPtr->getDofCoeffIdx()) {
592 // add only one when higher order terms present
593 double ave_mid = (tag_value[0][dofPtr->getDofCoeffIdx()] +
594 tag_value[1][dofPtr->getDofCoeffIdx()]) *
595 0.5;
596 tag_value[2][dofPtr->getDofCoeffIdx()] = ave_mid;
597 }
598 tag_value[2][dofPtr->getDofCoeffIdx()] += approx_val;
599 }
601}
FieldApproximationBase
approximation base
Definition definitions.h:58
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ AINSWORTH_LOBATTO_BASE
Definition definitions.h:62
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
auto type_from_handle(const EntityHandle h)
get type from entity handle
virtual moab::Interface & get_moab()=0
boost::shared_ptr< DofEntity > dofPtr
Shared pointer to DOF entity data.

◆ preProcess()

MoFEMErrorCode MoFEM::ProjectionFieldOn10NodeTet::preProcess ( )
virtual

Pre-processing function executed at loop initialization.

This virtual function is called once at the beginning of the loop. It is typically used for:

  • Zeroing matrices and vectors
  • Computing shape functions on reference elements
  • Preprocessing boundary conditions
  • Setting up initial computation state
Returns
Error code indicating success or failure

Reimplemented from MoFEM::Projection10NodeCoordsOnField.

Definition at line 473 of file Projection10NodeCoordsOnField.cpp.

473 {
475 if (!onCoords) {
476 if (onTag == "NoNE") {
477 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
478 "tag name not specified");
479 }
480 field_it = fieldsPtr->get<FieldName_mi_tag>().find(fieldName);
481 if (field_it == fieldsPtr->get<FieldName_mi_tag>().end()) {
482 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "field not found %s",
483 fieldName.c_str());
484 }
485 int field_rank = (*field_it)->getNbOfCoeffs();
486 VectorDouble def_VAL = ublas::zero_vector<double>(field_rank);
487 CHKERR mField.get_moab().tag_get_handle(
488 onTag.c_str(), field_rank, MB_TYPE_DOUBLE, th,
489 MB_TAG_CREAT | MB_TAG_SPARSE, &*def_VAL.data().begin());
491 }
492
493 L.resize(maxApproximationOrder + 1);
495 &*L.data().begin(), NULL, 3);
496 K.resize(10);
497 CHKERR LobattoKernel_polynomials(9, 0., NULL, &*K.data().begin(), NULL, 3);
498
500}
#define MOAB_THROW(err)
Check error code of MoAB function and throw MoFEM exception.
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, int dim)
Calculate Legendre approximation basis.
PetscErrorCode LobattoKernel_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Kernel Lobatto base functions.
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
UBlasVector< double > VectorDouble
Definition Types.hpp:68
const Field_multiIndex * fieldsPtr
Raw pointer to fields multi-index container.
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it

Member Data Documentation

◆ field_it

Field_multiIndex::index<FieldName_mi_tag>::type::iterator MoFEM::ProjectionFieldOn10NodeTet::field_it
protected

Definition at line 117 of file Projection10NodeCoordsOnField.hpp.

◆ K

VectorDouble MoFEM::ProjectionFieldOn10NodeTet::K
protected

Definition at line 119 of file Projection10NodeCoordsOnField.hpp.

◆ L

VectorDouble MoFEM::ProjectionFieldOn10NodeTet::L
protected

Definition at line 118 of file Projection10NodeCoordsOnField.hpp.

◆ maxApproximationOrder

const int MoFEM::ProjectionFieldOn10NodeTet::maxApproximationOrder

Definition at line 111 of file Projection10NodeCoordsOnField.hpp.

◆ onCoords

bool MoFEM::ProjectionFieldOn10NodeTet::onCoords

Definition at line 108 of file Projection10NodeCoordsOnField.hpp.

◆ onTag

std::string MoFEM::ProjectionFieldOn10NodeTet::onTag

Definition at line 109 of file Projection10NodeCoordsOnField.hpp.

◆ setNodes

bool MoFEM::ProjectionFieldOn10NodeTet::setNodes

Definition at line 107 of file Projection10NodeCoordsOnField.hpp.

◆ th

Tag MoFEM::ProjectionFieldOn10NodeTet::th

Definition at line 113 of file Projection10NodeCoordsOnField.hpp.


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