v0.9.0
Public Member Functions | Public 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 (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

bool setNodes
 
bool onCoords
 
std::string onTag
 
const int maxApproximationOrder
 
Tag th
 
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
 
VectorDouble L
 
VectorDouble K
 
- Public Attributes inherited from MoFEM::Projection10NodeCoordsOnField
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

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

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

184  {
186  if (dofPtr->getName() != fieldName)
188  if (setNodes) {
189  if (dofPtr->getEntType() == MBVERTEX) {
190  EntityHandle node = dofPtr->getEnt();
191  if (onCoords) {
192  coords.resize(3);
193  CHKERR mField.get_moab().get_coords(&node, 1,
194  &*coords.data().begin());
195  coords[dofPtr->getDofCoeffIdx()] = dofPtr->getFieldData();
196  CHKERR mField.get_moab().set_coords(&node, 1,
197  &*coords.data().begin());
198  } else {
199  int field_rank = (*field_it)->getNbOfCoeffs();
200  if (field_rank != dofPtr->getNbOfCoeffs()) {
201  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
202  "data inconsistency");
203  }
204  double *tag_value;
205  int tag_size;
206  CHKERR mField.get_moab().tag_get_by_ptr(
207  th, &node, 1, (const void **)&tag_value, &tag_size);
208  if (tag_size != dofPtr->getNbOfCoeffs()) {
209  SETERRQ(PETSC_COMM_SELF, 1, "data inconsistency");
210  }
211  tag_value[dofPtr->getDofCoeffIdx()] = dofPtr->getFieldData();
212  }
213  }
215  }
216  if (dofPtr->getEntType() != MBEDGE) {
218  }
219  EntityHandle edge = dofPtr->getEnt();
220  if (mField.get_moab().type_from_handle(edge) != MBEDGE) {
221  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
222  "this method works only elements which are type of MBEDGE");
223  }
224 
225  int num_nodes;
226  const EntityHandle *conn;
227  CHKERR mField.get_moab().get_connectivity(edge, conn, num_nodes, false);
228  if ((num_nodes != 2) && (num_nodes != 3)) {
229  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
230  "this method works only 4 node and 10 node tets");
231  }
232  if (num_nodes == 2) {
234  }
235 
236  if (dofPtr->getDofOrder() >= maxApproximationOrder) {
237  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
238  "too big approximation order, increase constant "
239  "max_ApproximationOrder");
240  }
241  double approx_val = 0;
242  FieldApproximationBase base = dofPtr->getApproxBase();
243  switch (base) {
245  approx_val = 0.25 * L[dofPtr->getDofOrder() - 2] * dofPtr->getFieldData();
246  break;
248  approx_val = 0.25 * K[dofPtr->getDofOrder() - 2] * dofPtr->getFieldData();
249  break;
250  default:
251  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not yet implemented");
252  }
253 
254  if (onCoords) {
255  coords.resize(num_nodes * 3);
256  CHKERR mField.get_moab().get_coords(conn, num_nodes,
257  &*coords.data().begin());
258  if (dofPtr->getEntDofIdx() == dofPtr->getDofCoeffIdx()) {
259  // add only one when higher order terms present
260  double ave_mid = (coords[3 * 0 + dofPtr->getDofCoeffIdx()] +
261  coords[3 * 1 + dofPtr->getDofCoeffIdx()]) *
262  0.5;
263  coords[2 * 3 + dofPtr->getDofCoeffIdx()] = ave_mid;
264  }
265  coords[2 * 3 + dofPtr->getDofCoeffIdx()] += approx_val;
266  CHKERR mField.get_moab().set_coords(&conn[2], 1, &coords[3 * 2]);
267  } else {
268  int tag_size;
269  double *tag_value[num_nodes];
270  CHKERR mField.get_moab().tag_get_by_ptr(
271  th, conn, num_nodes, (const void **)tag_value, &tag_size);
272  if (tag_size != dofPtr->getNbOfCoeffs()) {
273  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
274  "data inconsistency");
275  }
276  if (dofPtr->getEntDofIdx() == dofPtr->getDofCoeffIdx()) {
277  // add only one when higher order terms present
278  double ave_mid = (tag_value[0][dofPtr->getDofCoeffIdx()] +
279  tag_value[1][dofPtr->getDofCoeffIdx()]) *
280  0.5;
281  tag_value[2][dofPtr->getDofCoeffIdx()] = ave_mid;
282  }
283  tag_value[2][dofPtr->getDofCoeffIdx()] += approx_val;
284  }
286  }
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
FieldApproximationBase
approximation base
Definition: definitions.h:144
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:146
#define CHKERR
Inline error check.
Definition: definitions.h:596

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

155  {
157  if (!onCoords) {
158  if (onTag == "NoNE") {
159  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
160  "tag name not specified");
161  }
162  field_it = fieldsPtr->get<FieldName_mi_tag>().find(fieldName);
163  if (field_it == fieldsPtr->get<FieldName_mi_tag>().end()) {
164  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
165  "field not found %s", fieldName.c_str());
166  }
167  int field_rank = (*field_it)->getNbOfCoeffs();
168  VectorDouble def_VAL = ublas::zero_vector<double>(field_rank);
169  CHKERR mField.get_moab().tag_get_handle(
170  onTag.c_str(), field_rank, MB_TYPE_DOUBLE, th,
171  MB_TAG_CREAT | MB_TAG_SPARSE, &*def_VAL.data().begin());
172  MOAB_THROW(rval);
173  }
174 
175  L.resize(maxApproximationOrder + 1);
177  &*L.data().begin(), NULL, 3);
178  K.resize(10);
179  CHKERR LobattoKernel_polynomials(9, 0., NULL, &*K.data().begin(), NULL, 3);
180 
182  }
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Legendre approximation basis.
virtual moab::Interface & get_moab()=0
const Field_multiIndex * fieldsPtr
raw pointer to fields container
#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
#define MOAB_THROW(a)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:602
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:84
PetscErrorCode LobattoKernel_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Kernel Lobatto base functions.
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
#define CHKERR
Inline error check.
Definition: definitions.h:596
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:72

Member Data Documentation

◆ field_it

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

Definition at line 151 of file Projection10NodeCoordsOnField.hpp.

◆ K

VectorDouble MoFEM::ProjectionFieldOn10NodeTet::K

Definition at line 153 of file Projection10NodeCoordsOnField.hpp.

◆ L

VectorDouble MoFEM::ProjectionFieldOn10NodeTet::L

Definition at line 152 of file Projection10NodeCoordsOnField.hpp.

◆ maxApproximationOrder

const int MoFEM::ProjectionFieldOn10NodeTet::maxApproximationOrder

Definition at line 142 of file Projection10NodeCoordsOnField.hpp.

◆ onCoords

bool MoFEM::ProjectionFieldOn10NodeTet::onCoords

Definition at line 139 of file Projection10NodeCoordsOnField.hpp.

◆ onTag

std::string MoFEM::ProjectionFieldOn10NodeTet::onTag

Definition at line 140 of file Projection10NodeCoordsOnField.hpp.

◆ setNodes

bool MoFEM::ProjectionFieldOn10NodeTet::setNodes

Definition at line 138 of file Projection10NodeCoordsOnField.hpp.

◆ th

Tag MoFEM::ProjectionFieldOn10NodeTet::th

Definition at line 150 of file Projection10NodeCoordsOnField.hpp.


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