v0.14.0
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 ()
 function is run at the beginning of loop More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
- Public Member Functions inherited from MoFEM::DofMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 DofMethod ()=default
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities. More...
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements. More...
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements. More...
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 

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 }
 pass information about context of KSP/DM for with finite element is computed 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_X_T = 1 << 4, CTX_SET_X_TT = 1 << 6, CTX_SET_TIME = 1 << 7
}
 
using Switches = std::bitset< 8 >
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION, CTX_SNESSETJACOBIAN, CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION, CTX_TSSETRHSJACOBIAN, CTX_TSSETIFUNCTION, CTX_TSSETIJACOBIAN,
  CTX_TSTSMONITORSET, CTX_TSNONE
}
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 
- Public Attributes inherited from MoFEM::DofMethod
boost::shared_ptr< FieldfieldPtr
 
boost::shared_ptr< DofEntitydofPtr
 
boost::shared_ptr< NumeredDofEntitydofNumeredPtr
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities. More...
 
int rAnk
 processor rank More...
 
int sIze
 number of processors in communicator More...
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities More...
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities More...
 
const ProblemproblemPtr
 raw pointer to problem More...
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container More...
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities More...
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs More...
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements More...
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing. More...
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing. More...
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator. More...
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec & ksp_f
 
Mat & ksp_A
 
Mat & ksp_B
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 
Vec f
 
Mat A
 
Mat B
 
Vec x
 
Vec x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 snes solver More...
 
Vec & snes_x
 state vector More...
 
Vec & snes_f
 residual More...
 
Mat & snes_A
 jacobian matrix More...
 
Mat & snes_B
 preconditioner of jacobian matrix More...
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 time solver More...
 
TSContext ts_ctx
 
PetscInt ts_step
 time step number More...
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver) More...
 
PetscReal ts_aa
 shift for U_tt shift for U_tt More...
 
PetscReal ts_t
 time More...
 
PetscReal ts_dt
 time step size More...
 
Vec & ts_u
 state vector More...
 
Vec & ts_u_t
 time derivative of state vector More...
 
Vec & ts_u_tt
 second time derivative of state vector More...
 
Vec & ts_F
 residual vector More...
 
Mat & ts_A
 
Mat & ts_B
 Preconditioner for ts_A. More...
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 

Detailed Description

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

Examples
approx_sphere.cpp, dynamic_first_order_con_law.cpp, elasticity.cpp, elasticity_mixed_formulation.cpp, forces_and_sources_testing_edge_element.cpp, forces_and_sources_testing_flat_prism_element.cpp, hcurl_curl_operator.cpp, hdiv_divergence_operator.cpp, navier_stokes.cpp, nonlinear_dynamics.cpp, plastic.cpp, Remodeling.cpp, simple_contact.cpp, simple_contact_thermal.cpp, simple_interface.cpp, tensor_divergence_operator.cpp, test_jacobian_of_simple_contact_element.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.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 
)

Definition at line 10 of file Projection10NodeCoordsOnField.cpp.

12  : mField(m_field), fieldName(field_name), vErbose(verb) {}

Member Function Documentation

◆ operator()()

MoFEMErrorCode MoFEM::Projection10NodeCoordsOnField::operator() ( )
virtual

function is run for every finite element

It is used to calculate element local matrices and assembly. It can be used for post-processing.

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;
87  case DEMKOWICZ_JACOBI_BASE: {
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 }

◆ postProcess()

MoFEMErrorCode MoFEM::Projection10NodeCoordsOnField::postProcess ( )
virtual

function is run at the end of loop

It is used to assembly matrices and vectors, calculating global variables, f.e. total internal energy, ect.

Iterating over dofs: Example1 iterating over dofs in row by name of the field for(IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP(this,"DISPLACEMENT",it)) { ... }

Reimplemented from MoFEM::BasicMethod.

Definition at line 109 of file Projection10NodeCoordsOnField.cpp.

109  {
112 }

◆ preProcess()

MoFEMErrorCode MoFEM::Projection10NodeCoordsOnField::preProcess ( )
virtual

function is run at the beginning of loop

It is used to zeroing matrices and vectors, calculation of shape functions on reference element, preprocessing boundary conditions, etc.

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 }

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:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
EntityHandle
MoFEM::Projection10NodeCoordsOnField::fieldName
std::string fieldName
Definition: Projection10NodeCoordsOnField.hpp:38
MoFEM::CoreInterface::get_field_structure
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
MoFEM::Projection10NodeCoordsOnField::mField
Interface & mField
Definition: Projection10NodeCoordsOnField.hpp:37
MoFEM::DofMethod::dofPtr
boost::shared_ptr< DofEntity > dofPtr
Definition: LoopMethods.hpp:505
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM::Projection10NodeCoordsOnField::vErbose
int vErbose
Definition: Projection10NodeCoordsOnField.hpp:39
MoFEM::Projection10NodeCoordsOnField::aveMidCoord
VectorDouble3 aveMidCoord
Definition: Projection10NodeCoordsOnField.hpp:42
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
AINSWORTH_LOBATTO_BASE
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:62
MoFEM::type_from_handle
auto type_from_handle(const EntityHandle h)
get type from entity handle
Definition: Templates.hpp:1869
Legendre_polynomials
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Legendre approximation basis.
Definition: base_functions.c:15
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
LOBATTO_PHI0
#define LOBATTO_PHI0(x)
Definitions taken from Hermes2d code.
Definition: base_functions.h:126
FTensor::dd
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
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::Projection10NodeCoordsOnField::dOf
VectorDouble3 dOf
Definition: Projection10NodeCoordsOnField.hpp:45
MOFEM_TAG_AND_LOG_C
#define MOFEM_TAG_AND_LOG_C(channel, severity, tag, format,...)
Tag and log in channel.
Definition: LogManager.hpp:370
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::Projection10NodeCoordsOnField::coords
VectorDouble coords
Definition: Projection10NodeCoordsOnField.hpp:41
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEM::Projection10NodeCoordsOnField::midNodeCoord
VectorDouble3 midNodeCoord
Definition: Projection10NodeCoordsOnField.hpp:43
MoFEM::Projection10NodeCoordsOnField::diffNodeCoord
VectorDouble3 diffNodeCoord
Definition: Projection10NodeCoordsOnField.hpp:44