v0.9.0
Public Member Functions | Public 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 (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 DofMethod ()
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()
 
MoFEMErrorCode setKspCtx (const KSPContext &ctx)
 set operator type More...
 
MoFEMErrorCode setKsp (KSP ksp)
 set solver More...
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of 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
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()
 
MoFEMErrorCode setSnesCtx (const SNESContext &ctx)
 Set SNES context. More...
 
MoFEMErrorCode setSnes (SNES snes)
 Set SNES instance. More...
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. More...
 

Public Attributes

InterfacemField
 
std::string fieldName
 
int vErbose
 
VectorDouble coords
 
VectorDouble3 aveMidCoord
 
VectorDouble3 midNodeCoord
 
VectorDouble3 diffNodeCoord
 
VectorDouble3 dOf
 
- 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...
 
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
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec ksp_f
 the right hand side vector More...
 
Mat ksp_A
 matrix More...
 
Mat ksp_B
 preconditioner matrix More...
 
- 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
TSContext ts_ctx
 
TS ts
 time solver 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...
 
PetscInt ts_step
 time step More...
 
PetscReal ts_a
 shift for U_tt (see PETSc Time Solver) More...
 
PetscReal ts_v
 shift for U_t shift for U_t More...
 
PetscReal ts_t
 time More...
 

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::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
}
 
- Protected Member Functions inherited from MoFEM::UnknownInterface
boost::typeindex::type_index getClassIdx (const MOFEMuuid &uid) const
 Get type name for interface Id. More...
 
MOFEMuuid getUId (const boost::typeindex::type_index &class_idx) const
 Get interface Id for class name. More...
 

Detailed Description

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

Examples
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, Remodeling.cpp, simple_interface.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.

Definition at line 36 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

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 57 of file Projection10NodeCoordsOnField.hpp.

57  {
59  if (dofPtr == NULL) {
60  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
61  }
62  if (dofPtr->getName() != fieldName)
64  if (dofPtr->getEntType() == MBVERTEX) {
65  EntityHandle node = dofPtr->getEnt();
66  coords.resize(3);
67  CHKERR mField.get_moab().get_coords(&node, 1, &*coords.data().begin());
68  dofPtr->getFieldData() = coords[dofPtr->getDofCoeffIdx()];
69  if (vErbose > 0) {
70  PetscPrintf(mField.get_comm(), "val = %6.7e\n", dofPtr->getFieldData());
71  }
73  }
74  if (dofPtr->getEntType() != MBEDGE) {
76  }
77  if (dofPtr->getEntDofIdx() != dofPtr->getDofCoeffIdx()) {
79  }
80  EntityHandle edge = dofPtr->getEnt();
81  if (mField.get_moab().type_from_handle(edge) != MBEDGE) {
82  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
83  "this method works only elements which are type of MBEDGE");
84  }
85  // coords
86  int num_nodes;
87  const EntityHandle *conn;
88  CHKERR mField.get_moab().get_connectivity(edge, conn, num_nodes, false);
89  if ((num_nodes != 2) && (num_nodes != 3)) {
90  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
91  "this method works only 4 node and 10 node tets");
92  }
93  coords.resize(num_nodes * 3);
94  CHKERR mField.get_moab().get_coords(conn, num_nodes,
95  &*coords.data().begin());
96  aveMidCoord.resize(3);
97  midNodeCoord.resize(3);
98  for (int dd = 0; dd < 3; dd++) {
99  aveMidCoord[dd] = (coords[0 * 3 + dd] + coords[1 * 3 + dd]) * 0.5;
100  if (num_nodes == 3) {
101  midNodeCoord[dd] = coords[2 * 3 + dd];
102  } else {
104  }
105  }
106  double edge_shape_function_val = 0.25;
107  FieldApproximationBase base = dofPtr->getApproxBase();
108  switch (base) {
110  break;
112  edge_shape_function_val *= LOBATTO_PHI0(0);
113  break;
114  default:
115  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not yet implemented");
116  }
117 
118  diffNodeCoord.resize(3);
119  ublas::noalias(diffNodeCoord) = midNodeCoord - aveMidCoord;
120  dOf.resize(3);
121  ublas::noalias(dOf) = diffNodeCoord / edge_shape_function_val;
122  if (dofPtr->getNbOfCoeffs() != 3) {
123  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
124  "this method works only fields which are rank 3");
125  }
126  dofPtr->getFieldData() = dOf[dofPtr->getDofCoeffIdx()];
128  }
boost::shared_ptr< DofEntity > dofPtr
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
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
FieldApproximationBase
approximation base
Definition: definitions.h:144
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:146
#define LOBATTO_PHI0(x)
Definitions taken from Hermes2d code.
#define CHKERR
Inline error check.
Definition: definitions.h:596
virtual MPI_Comm & get_comm() const =0

◆ 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 130 of file Projection10NodeCoordsOnField.hpp.

130  {
133  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508

◆ 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 46 of file Projection10NodeCoordsOnField.hpp.

46  {
49  }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508

Member Data Documentation

◆ aveMidCoord

VectorDouble3 MoFEM::Projection10NodeCoordsOnField::aveMidCoord

Definition at line 52 of file Projection10NodeCoordsOnField.hpp.

◆ coords

VectorDouble MoFEM::Projection10NodeCoordsOnField::coords

Definition at line 51 of file Projection10NodeCoordsOnField.hpp.

◆ diffNodeCoord

VectorDouble3 MoFEM::Projection10NodeCoordsOnField::diffNodeCoord

Definition at line 54 of file Projection10NodeCoordsOnField.hpp.

◆ dOf

VectorDouble3 MoFEM::Projection10NodeCoordsOnField::dOf

Definition at line 55 of file Projection10NodeCoordsOnField.hpp.

◆ fieldName

std::string MoFEM::Projection10NodeCoordsOnField::fieldName

Definition at line 39 of file Projection10NodeCoordsOnField.hpp.

◆ mField

Interface& MoFEM::Projection10NodeCoordsOnField::mField

Definition at line 38 of file Projection10NodeCoordsOnField.hpp.

◆ midNodeCoord

VectorDouble3 MoFEM::Projection10NodeCoordsOnField::midNodeCoord

Definition at line 53 of file Projection10NodeCoordsOnField.hpp.

◆ vErbose

int MoFEM::Projection10NodeCoordsOnField::vErbose

Definition at line 40 of file Projection10NodeCoordsOnField.hpp.


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