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

Projection of edge entities with one mid-node on hierarchical basis. More...

#include "src/approximation/Projection10NodeCoordsOnField.hpp"

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

Public Member Functions

 Projection10NodeCoordsOnField (Interface &m_field, std::string field_name, int verb=0)
 
MoFEMErrorCode preProcess ()
 Pre-processing function executed at loop initialization.
 
MoFEMErrorCode operator() ()
 Main operator function executed for each loop iteration.
 
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.
 

Protected Attributes

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.
 
- 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.
 
- 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

Projection of edge entities with one mid-node on hierarchical basis.

Examples
mofem/atom_tests/forces_and_sources_testing_edge_element.cpp, mofem/atom_tests/forces_and_sources_testing_flat_prism_element.cpp, mofem/atom_tests/hcurl_curl_operator.cpp, mofem/atom_tests/hdiv_divergence_operator.cpp, mofem/atom_tests/simple_interface.cpp, mofem/atom_tests/tensor_divergence_operator.cpp, mofem/tutorials/adv-0/plastic.cpp, mofem/tutorials/adv-4/dynamic_first_order_con_law.cpp, mofem/tutorials/cor-10/navier_stokes.cpp, mofem/tutorials/cor-7/elasticity_mixed_formulation.cpp, mofem/tutorials/vec-3/nonlinear_dynamic_elastic.cpp, mofem/tutorials/vec-4/approx_sphere.cpp, mofem/tutorials/vec-7/adjoint.cpp, mofem/users_modules/adolc-plasticity/adolc_plasticity.cpp, mofem/users_modules/basic_finite_elements/atom_tests/testing_jacobian_of_hook_element.cpp, mofem/users_modules/basic_finite_elements/atom_tests/testing_jacobian_of_hook_scaled_with_density_element.cpp, mofem/users_modules/basic_finite_elements/elasticity/elasticity.cpp, mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp, mofem/users_modules/bone_remodelling/src/impl/Remodeling.cpp, and mofem/users_modules/eshelbian_plasticity/src/impl/EshelbianPlasticity.cpp.

Definition at line 24 of file Projection10NodeCoordsOnField.hpp.

Constructor & Destructor Documentation

◆ Projection10NodeCoordsOnField()

MoFEM::Projection10NodeCoordsOnField::Projection10NodeCoordsOnField ( Interface m_field,
std::string  field_name,
int  verb = 0 
)

Member Function Documentation

◆ operator()()

MoFEMErrorCode MoFEM::Projection10NodeCoordsOnField::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::BasicMethod.

Reimplemented in MoFEM::ProjectionFieldOn10NodeTet.

Definition at line 25 of file Projection10NodeCoordsOnField.cpp.

25 {
27
28 if (dofPtr == NULL) {
29 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
30 }
31 if (dofPtr->getName() != fieldName)
33 if (dofPtr->getEntType() == MBVERTEX) {
34 EntityHandle node = dofPtr->getEnt();
35 coords.resize(3);
36 CHKERR mField.get_moab().get_coords(&node, 1, &*coords.data().begin());
37 dofPtr->getFieldData() = coords[dofPtr->getDofCoeffIdx()];
38 if (vErbose > 0) {
39 MOFEM_TAG_AND_LOG_C("SELF", Sev::noisy, "Projection10NodeCoordsOnField",
40 "val = %6.7e\n", dofPtr->getFieldData());
41 MOFEM_LOG_CHANNEL("SELF");
42 }
43
45 }
46 if (dofPtr->getEntType() != MBEDGE) {
48 }
49 if (dofPtr->getEntDofIdx() != dofPtr->getDofCoeffIdx()) {
51 }
52 EntityHandle edge = dofPtr->getEnt();
53 if (type_from_handle(edge) != MBEDGE) {
54 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
55 "this method works only elements which are type of MBEDGE");
56 }
57 // coords
58 int num_nodes;
59 const EntityHandle *conn;
60 CHKERR mField.get_moab().get_connectivity(edge, conn, num_nodes, false);
61 if ((num_nodes != 2) && (num_nodes != 3)) {
62 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
63 "this method works only 4 node and 10 node tets");
64 }
65 coords.resize(num_nodes * 3);
66 CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
67 aveMidCoord.resize(3);
68 midNodeCoord.resize(3);
69 for (int dd = 0; dd < 3; dd++) {
70 aveMidCoord[dd] = (coords[0 * 3 + dd] + coords[1 * 3 + dd]) * 0.5;
71 if (num_nodes == 3) {
72 midNodeCoord[dd] = coords[2 * 3 + dd];
73 } else {
75 }
76 }
77
78 const auto base = dofPtr->getApproxBase();
79 double edge_shape_function_val;
80 switch (base) {
82 edge_shape_function_val = 0.25;
83 break;
85 edge_shape_function_val = 0.25 * LOBATTO_PHI0(0);
86 break;
88 double L[3];
89 CHKERR Legendre_polynomials(2, 0, NULL, L, NULL, 1);
90 edge_shape_function_val = 0.125 * LOBATTO_PHI0(0);
91 }; break;
92 default:
93 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not yet implemented");
94 }
95
96 diffNodeCoord.resize(3);
97 ublas::noalias(diffNodeCoord) = midNodeCoord - aveMidCoord;
98 dOf.resize(3);
99 ublas::noalias(dOf) = diffNodeCoord / edge_shape_function_val;
100 if (dofPtr->getNbOfCoeffs() > 3) {
101 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
102 "this method works only fields which are rank 3 or lower");
103 }
104 dofPtr->getFieldData() = dOf[dofPtr->getDofCoeffIdx()];
105
107}
#define MOFEM_TAG_AND_LOG_C(channel, severity, tag, format,...)
Tag and log in channel.
#define LOBATTO_PHI0(x)
Definitions taken from Hermes2d code.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ AINSWORTH_LOBATTO_BASE
Definition definitions.h:62
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#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 ...
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, int dim)
Calculate Legendre approximation basis.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition ddTensor0.hpp:33
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.

◆ postProcess()

MoFEMErrorCode MoFEM::Projection10NodeCoordsOnField::postProcess ( )
virtual

Post-processing function executed at loop completion.

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

  • Assembling matrices and vectors to global system
  • Computing global variables (e.g., total internal energy)
  • Finalizing computation results
  • Cleanup operations

Example of iterating over DOFs:

for(_IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP_(this,"DISPLACEMENT",it)) {
// Process each DOF
}
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Definition at line 109 of file Projection10NodeCoordsOnField.cpp.

◆ preProcess()

MoFEMErrorCode MoFEM::Projection10NodeCoordsOnField::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::BasicMethod.

Reimplemented in MoFEM::ProjectionFieldOn10NodeTet.

Definition at line 14 of file Projection10NodeCoordsOnField.cpp.

14 {
16 auto field_ptr = mField.get_field_structure(fieldName);
17 if (field_ptr->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
18 MOFEM_TAG_AND_LOG("WORLD", Sev::warning, "Projection10NodeCoordsOnField")
19 << "Only working well for first order AINSWORTH_BERNSTEIN_BEZIER_BASE!";
20 MOFEM_LOG_CHANNEL("WORLD");
21 }
23}
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition definitions.h:64
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure

Member Data Documentation

◆ aveMidCoord

VectorDouble3 MoFEM::Projection10NodeCoordsOnField::aveMidCoord
protected

Definition at line 42 of file Projection10NodeCoordsOnField.hpp.

◆ coords

VectorDouble MoFEM::Projection10NodeCoordsOnField::coords
protected

Definition at line 41 of file Projection10NodeCoordsOnField.hpp.

◆ diffNodeCoord

VectorDouble3 MoFEM::Projection10NodeCoordsOnField::diffNodeCoord
protected

Definition at line 44 of file Projection10NodeCoordsOnField.hpp.

◆ dOf

VectorDouble3 MoFEM::Projection10NodeCoordsOnField::dOf
protected

Definition at line 45 of file Projection10NodeCoordsOnField.hpp.

◆ fieldName

std::string MoFEM::Projection10NodeCoordsOnField::fieldName
protected

Definition at line 38 of file Projection10NodeCoordsOnField.hpp.

◆ mField

Interface& MoFEM::Projection10NodeCoordsOnField::mField
protected

Definition at line 37 of file Projection10NodeCoordsOnField.hpp.

◆ midNodeCoord

VectorDouble3 MoFEM::Projection10NodeCoordsOnField::midNodeCoord
protected

Definition at line 43 of file Projection10NodeCoordsOnField.hpp.

◆ vErbose

int MoFEM::Projection10NodeCoordsOnField::vErbose
protected

Definition at line 39 of file Projection10NodeCoordsOnField.hpp.


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