v0.14.0
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 ()
 function is run at the beginning of loop More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::Projection10NodeCoordsOnField
 Projection10NodeCoordsOnField (Interface &m_field, std::string field_name, int verb=0)
 
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 reference 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...
 

Public Attributes

bool setNodes
 
bool onCoords
 
std::string onTag
 
const int maxApproximationOrder
 
Tag th
 
- 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...
 

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

Examples
navier_stokes.cpp.

Definition at line 49 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" 
)

Definition at line 114 of file Projection10NodeCoordsOnField.cpp.

119  : Projection10NodeCoordsOnField(m_field, _fieldName), setNodes(set_nodes),
120  onCoords(on_coords), onTag(on_tag), maxApproximationOrder(20) {}

Member Function Documentation

◆ operator()()

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

Definition at line 156 of file Projection10NodeCoordsOnField.cpp.

156  {
158  if (dofPtr->getName() != fieldName)
160  if (setNodes) {
161  if (dofPtr->getEntType() == MBVERTEX) {
162  EntityHandle node = dofPtr->getEnt();
163  if (onCoords) {
164  coords.resize(3);
165  CHKERR mField.get_moab().get_coords(&node, 1, &*coords.data().begin());
166  coords[dofPtr->getDofCoeffIdx()] = dofPtr->getFieldData();
167  CHKERR mField.get_moab().set_coords(&node, 1, &*coords.data().begin());
168  } else {
169  int field_rank = (*field_it)->getNbOfCoeffs();
170  if (field_rank != dofPtr->getNbOfCoeffs()) {
171  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
172  "data inconsistency");
173  }
174  double *tag_value;
175  int tag_size;
176  CHKERR mField.get_moab().tag_get_by_ptr(
177  th, &node, 1, (const void **)&tag_value, &tag_size);
178  if (tag_size != dofPtr->getNbOfCoeffs()) {
179  SETERRQ(PETSC_COMM_SELF, 1, "data inconsistency");
180  }
181  tag_value[dofPtr->getDofCoeffIdx()] = dofPtr->getFieldData();
182  }
183  }
185  }
186  if (dofPtr->getEntType() != MBEDGE) {
188  }
189  EntityHandle edge = dofPtr->getEnt();
190  if (type_from_handle(edge) != MBEDGE) {
191  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
192  "this method works only elements which are type of MBEDGE");
193  }
194 
195  int num_nodes;
196  const EntityHandle *conn;
197  CHKERR mField.get_moab().get_connectivity(edge, conn, num_nodes, false);
198  if ((num_nodes != 2) && (num_nodes != 3)) {
199  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
200  "this method works only 4 node and 10 node tets");
201  }
202  if (num_nodes == 2) {
204  }
205 
206  if (dofPtr->getDofOrder() >= maxApproximationOrder) {
207  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
208  "too big approximation order, increase constant "
209  "max_ApproximationOrder");
210  }
211  double approx_val = 0;
212  FieldApproximationBase base = dofPtr->getApproxBase();
213  switch (base) {
215  approx_val = 0.25 * L[dofPtr->getDofOrder() - 2] * dofPtr->getFieldData();
216  break;
218  approx_val = 0.25 * K[dofPtr->getDofOrder() - 2] * dofPtr->getFieldData();
219  break;
220  default:
221  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not yet implemented");
222  }
223 
224  if (onCoords) {
225  coords.resize(num_nodes * 3);
226  CHKERR mField.get_moab().get_coords(conn, num_nodes,
227  &*coords.data().begin());
228  if (dofPtr->getEntDofIdx() == dofPtr->getDofCoeffIdx()) {
229  // add only one when higher order terms present
230  double ave_mid = (coords[3 * 0 + dofPtr->getDofCoeffIdx()] +
231  coords[3 * 1 + dofPtr->getDofCoeffIdx()]) *
232  0.5;
233  coords[2 * 3 + dofPtr->getDofCoeffIdx()] = ave_mid;
234  }
235  coords[2 * 3 + dofPtr->getDofCoeffIdx()] += approx_val;
236  CHKERR mField.get_moab().set_coords(&conn[2], 1, &coords[3 * 2]);
237  } else {
238  int tag_size;
239  double *tag_value[num_nodes];
240  CHKERR mField.get_moab().tag_get_by_ptr(
241  th, conn, num_nodes, (const void **)tag_value, &tag_size);
242  if (tag_size != dofPtr->getNbOfCoeffs()) {
243  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
244  }
245  if (dofPtr->getEntDofIdx() == dofPtr->getDofCoeffIdx()) {
246  // add only one when higher order terms present
247  double ave_mid = (tag_value[0][dofPtr->getDofCoeffIdx()] +
248  tag_value[1][dofPtr->getDofCoeffIdx()]) *
249  0.5;
250  tag_value[2][dofPtr->getDofCoeffIdx()] = ave_mid;
251  }
252  tag_value[2][dofPtr->getDofCoeffIdx()] += approx_val;
253  }
255 }

◆ preProcess()

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

Definition at line 127 of file Projection10NodeCoordsOnField.cpp.

127  {
129  if (!onCoords) {
130  if (onTag == "NoNE") {
131  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
132  "tag name not specified");
133  }
134  field_it = fieldsPtr->get<FieldName_mi_tag>().find(fieldName);
135  if (field_it == fieldsPtr->get<FieldName_mi_tag>().end()) {
136  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "field not found %s",
137  fieldName.c_str());
138  }
139  int field_rank = (*field_it)->getNbOfCoeffs();
140  VectorDouble def_VAL = ublas::zero_vector<double>(field_rank);
141  CHKERR mField.get_moab().tag_get_handle(
142  onTag.c_str(), field_rank, MB_TYPE_DOUBLE, th,
143  MB_TAG_CREAT | MB_TAG_SPARSE, &*def_VAL.data().begin());
144  MOAB_THROW(rval);
145  }
146 
147  L.resize(maxApproximationOrder + 1);
149  &*L.data().begin(), NULL, 3);
150  K.resize(10);
151  CHKERR LobattoKernel_polynomials(9, 0., NULL, &*K.data().begin(), NULL, 3);
152 
154 }

Member Data Documentation

◆ field_it

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

Definition at line 71 of file Projection10NodeCoordsOnField.hpp.

◆ K

VectorDouble MoFEM::ProjectionFieldOn10NodeTet::K
protected

Definition at line 73 of file Projection10NodeCoordsOnField.hpp.

◆ L

VectorDouble MoFEM::ProjectionFieldOn10NodeTet::L
protected

Definition at line 72 of file Projection10NodeCoordsOnField.hpp.

◆ maxApproximationOrder

const int MoFEM::ProjectionFieldOn10NodeTet::maxApproximationOrder

Definition at line 65 of file Projection10NodeCoordsOnField.hpp.

◆ onCoords

bool MoFEM::ProjectionFieldOn10NodeTet::onCoords

Definition at line 62 of file Projection10NodeCoordsOnField.hpp.

◆ onTag

std::string MoFEM::ProjectionFieldOn10NodeTet::onTag

Definition at line 63 of file Projection10NodeCoordsOnField.hpp.

◆ setNodes

bool MoFEM::ProjectionFieldOn10NodeTet::setNodes

Definition at line 61 of file Projection10NodeCoordsOnField.hpp.

◆ th

Tag MoFEM::ProjectionFieldOn10NodeTet::th

Definition at line 67 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:460
EntityHandle
MoFEM::Projection10NodeCoordsOnField::fieldName
std::string fieldName
Definition: Projection10NodeCoordsOnField.hpp:38
MoFEM::ProjectionFieldOn10NodeTet::maxApproximationOrder
const int maxApproximationOrder
Definition: Projection10NodeCoordsOnField.hpp:65
MOAB_THROW
#define MOAB_THROW(err)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:554
MoFEM::Projection10NodeCoordsOnField::mField
Interface & mField
Definition: Projection10NodeCoordsOnField.hpp:37
MoFEM::DofMethod::dofPtr
boost::shared_ptr< DofEntity > dofPtr
Definition: LoopMethods.hpp:505
MoFEM::ProjectionFieldOn10NodeTet::field_it
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
Definition: Projection10NodeCoordsOnField.hpp:71
MoFEM::BasicMethod::fieldsPtr
const Field_multiIndex * fieldsPtr
raw pointer to fields container
Definition: LoopMethods.hpp:252
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
MoFEM::ProjectionFieldOn10NodeTet::K
VectorDouble K
Definition: Projection10NodeCoordsOnField.hpp:73
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM::ProjectionFieldOn10NodeTet::onCoords
bool onCoords
Definition: Projection10NodeCoordsOnField.hpp:62
MoFEM::ProjectionFieldOn10NodeTet::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.hpp:72
MoFEM::ProjectionFieldOn10NodeTet::setNodes
bool setNodes
Definition: Projection10NodeCoordsOnField.hpp:61
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:1898
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
LobattoKernel_polynomials
PetscErrorCode LobattoKernel_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Kernel Lobatto base functions.
Definition: base_functions.c:328
MoFEM::ProjectionFieldOn10NodeTet::onTag
std::string onTag
Definition: Projection10NodeCoordsOnField.hpp:63
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
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
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:453
MoFEM::ProjectionFieldOn10NodeTet::th
Tag th
Definition: Projection10NodeCoordsOnField.hpp:67
MoFEM::Projection10NodeCoordsOnField::Projection10NodeCoordsOnField
Projection10NodeCoordsOnField(Interface &m_field, std::string field_name, int verb=0)
Definition: Projection10NodeCoordsOnField.cpp:10
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32