v0.8.23
Public Member Functions | Private Member Functions | Private Attributes | List of all members
MoFEM::FlatPrismPolynomialBase Struct Reference

Calculate base functions on tetrahedralFIXME: Need moab and mofem finite element structure to work (that not perfect) More...

#include <src/approximation/FlatPrismPolynomialBase.hpp>

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

Public Member Functions

MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FlatPrismPolynomialBase ()
 
 ~FlatPrismPolynomialBase ()
 
MoFEMErrorCode getValue (MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 
- Public Member Functions inherited from MoFEM::BaseFunction
 BaseFunction ()
 
 ~BaseFunction ()
 
virtual MoFEMErrorCode getValue (MatrixDouble &pts_x, MatrixDouble &pts_t, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 
- 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 ()
 
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
 

Private Member Functions

MoFEMErrorCode getValueH1 (MatrixDouble &pts)
 
MoFEMErrorCode getValueL2 (MatrixDouble &pts)
 
MoFEMErrorCode getValueHdiv (MatrixDouble &pts)
 
MoFEMErrorCode getValueHcurl (MatrixDouble &pts)
 

Private Attributes

FlatPrismPolynomialBaseCtxcTx
 
int numNodes
 
const EntityHandleconnPrism
 
const EntityHandleconnFace3
 
const EntityHandleconnFace4
 
int faceNodes [2][3]
 
MatrixDouble N
 
MatrixDouble diffN
 

Additional Inherited Members

- 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

Calculate base functions on tetrahedral

FIXME: Need moab and mofem finite element structure to work (that not perfect)

Definition at line 59 of file FlatPrismPolynomialBase.hpp.

Constructor & Destructor Documentation

◆ FlatPrismPolynomialBase()

FlatPrismPolynomialBase::FlatPrismPolynomialBase ( )

Definition at line 65 of file FlatPrismPolynomialBase.cpp.

65 {}

◆ ~FlatPrismPolynomialBase()

FlatPrismPolynomialBase::~FlatPrismPolynomialBase ( )

Definition at line 64 of file FlatPrismPolynomialBase.cpp.

64 {}

Member Function Documentation

◆ getValue()

MoFEMErrorCode FlatPrismPolynomialBase::getValue ( MatrixDouble pts,
boost::shared_ptr< BaseFunctionCtx ctx_ptr 
)
virtual

Reimplemented from MoFEM::BaseFunction.

Definition at line 68 of file FlatPrismPolynomialBase.cpp.

69  {
70 
72 
74  CHKERR ctx_ptr->query_interface(IDD_FLATPRISM_BASE_FUNCTION, &iface);
75  cTx = reinterpret_cast<FlatPrismPolynomialBaseCtx *>(iface);
76  if (!cTx->fePtr) {
77  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
78  "Pointer to element should be given "
79  "when EntPolynomialBaseCtx is constructed "
80  "(use different constructor)");
81  }
82 
83  int nb_gauss_pts = pts.size2();
84  if (!nb_gauss_pts) {
86  }
87 
88  if (pts.size1() < 1) {
89  SETERRQ(
90  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
91  "Wrong dimension of pts, should be at least 3 rows with coordinates");
92  }
93 
94  const FieldApproximationBase base = cTx->bAse;
96 
97  if (cTx->copyNodeBase == LASTBASE) {
98  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
99  } else {
100  data.dataOnEntities[MBVERTEX][0].getNSharedPtr(base) =
101  data.dataOnEntities[MBVERTEX][0].getNSharedPtr(cTx->copyNodeBase);
102  }
103  data.dataOnEntities[MBVERTEX][0].getN(base).resize(nb_gauss_pts, 6, false);
104  data.dataOnEntities[MBVERTEX][0].getDiffN(base).resize(nb_gauss_pts, 12,
105  false);
106  if (data.dataOnEntities[MBVERTEX][0].getN(base).size1() !=
107  (unsigned int)nb_gauss_pts) {
108  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
109  "Base functions or nodes has wrong number of integration points "
110  "for base %s",
111  ApproximationBaseNames[base]);
112  }
113 
114  data.dataOnEntities[MBVERTEX][0].getN(base).resize(nb_gauss_pts, 6, false);
115  data.dataOnEntities[MBVERTEX][0].getDiffN(base).resize(nb_gauss_pts, 12,
116  false);
117  N.resize(nb_gauss_pts, 3, false);
118  diffN.resize(3, 2, false);
119  CHKERR ShapeMBTRI(&*N.data().begin(), &pts(0, 0), &pts(1, 0), nb_gauss_pts);
120  CHKERR ShapeDiffMBTRI(&*diffN.data().begin());
121 
122  // This is needed to have proper order of nodes on faces
123 
124  CHKERR cTx->mOab.get_connectivity(cTx->fePtr->getEnt(), connPrism, numNodes,
125  true);
126  SideNumber_multiIndex &side_table =
127  const_cast<SideNumber_multiIndex &>(cTx->fePtr->getSideNumberTable());
128  SideNumber_multiIndex::nth_index<1>::type::iterator siit3 =
129  side_table.get<1>().find(boost::make_tuple(MBTRI, 3));
130  SideNumber_multiIndex::nth_index<1>::type::iterator siit4 =
131  side_table.get<1>().find(boost::make_tuple(MBTRI, 4));
132  if (siit3 == side_table.get<1>().end())
133  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
134  if (siit4 == side_table.get<1>().end())
135  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
136  CHKERR cTx->mOab.get_connectivity(siit3->get()->ent, connFace3, numNodes,
137  true);
138  CHKERR cTx->mOab.get_connectivity(siit4->get()->ent, connFace4, numNodes,
139  true);
140 
141  for (int nn = 0; nn < 3; nn++) {
142  faceNodes[0][nn] = std::distance(
143  connPrism, std::find(connPrism, connPrism + 3, connFace3[nn]));
144  faceNodes[1][nn] = std::distance(
145  connPrism + 3, std::find(connPrism + 3, connPrism + 6, connFace4[nn]));
146  for (int gg = 0; gg < nb_gauss_pts; gg++) {
147  double val = N(gg, nn);
148  double val_x = diffN(nn, 0);
149  double val_y = diffN(nn, 1);
150  data.dataOnEntities[MBVERTEX][0].getN(base)(gg, nn) = val;
151  data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 2 * nn + 0) = val_x;
152  data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 2 * nn + 1) = val_y;
153  data.dataOnEntities[MBVERTEX][0].getN(base)(gg, 3 + nn) = val;
154  data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 6 + 2 * nn + 0) =
155  val_x;
156  data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 6 + 2 * nn + 1) =
157  val_y;
158  }
159  }
160  // for(int nn = 0;nn<3;nn++) {
161  // if(faceNodes[0][nn]!=faceNodes[1][nn]) {
162  // SETERRQ(
163  // PETSC_COMM_SELF,
164  // MOFEM_DATA_INCONSISTENCY,
165  // "Node order different on both faces"
166  // );
167  // }
168  // }
169 
170  switch (cTx->sPace) {
171  case H1:
172  CHKERR getValueH1(pts);
173  break;
174  case HDIV:
175  CHKERR getValueHdiv(pts);
176  break;
177  case HCURL:
178  CHKERR getValueHcurl(pts);
179  break;
180  case L2:
181  CHKERR getValueL2(pts);
182  break;
183  default:
184  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
185  }
186 
188 }
field with continuous normal traction
Definition: definitions.h:173
data structure for finite element entityIt keeps that about indices of degrees of freedom,...
static const MOFEMuuid IDD_FLATPRISM_BASE_FUNCTION
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
const FieldApproximationBase bAse
base class for all interface classes
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
const FieldApproximationBase copyNodeBase
MoFEMErrorCode getValueHcurl(MatrixDouble &pts)
static const char *const ApproximationBaseNames[]
Definition: definitions.h:158
MoFEMErrorCode getValueHdiv(MatrixDouble &pts)
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:206
FieldApproximationBase
approximation base
Definition: definitions.h:143
DataForcesAndSourcesCore & dAta
field with continuous tangents
Definition: definitions.h:172
const NumeredEntFiniteElement * fePtr
SideNumber_multiIndex & getSideNumberTable() const
#define CHKERR
Inline error check.
Definition: definitions.h:596
FlatPrismPolynomialBaseCtx * cTx
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:194
continuous field
Definition: definitions.h:171
multi_index_container< boost::shared_ptr< SideNumber >, indexed_by< hashed_unique< member< SideNumber, EntityHandle, &SideNumber::ent > >, ordered_non_unique< composite_key< SideNumber, const_mem_fun< SideNumber, EntityType, &SideNumber::getEntType >, member< SideNumber, char, &SideNumber::side_number > > >, ordered_non_unique< const_mem_fun< SideNumber, EntityType, &SideNumber::getEntType > > > > SideNumber_multiIndex
SideNumber_multiIndex for SideNumber.
EntityHandle getEnt() const
Get element entity.
MoFEMErrorCode getValueL2(MatrixDouble &pts)
MoFEMErrorCode getValueH1(MatrixDouble &pts)
field with C-1 continuity
Definition: definitions.h:174

◆ getValueH1()

MoFEMErrorCode FlatPrismPolynomialBase::getValueH1 ( MatrixDouble pts)
private

Definition at line 190 of file FlatPrismPolynomialBase.cpp.

190  {
191 
193 
195  const FieldApproximationBase base = cTx->bAse;
196  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s, double *L,
197  double *diffL, const int dim) =
199 
200  int nb_gauss_pts = pts.size2();
201 
202  // edges
203  if (data.dataOnEntities[MBEDGE].size() != 9) {
204  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
205  }
206  int valid_edges[] = {1, 1, 1, 0, 0, 0, 1, 1, 1};
207  int sense[9], order[9];
208  double *H1edgeN[9], *diffH1edgeN[9];
209  if ((data.spacesOnEntities[MBEDGE]).test(H1)) {
210  for (int ee = 0; ee < 9; ee++) {
211  if (!valid_edges[ee])
212  continue;
213  if (data.dataOnEntities[MBEDGE][ee].getSense() == 0) {
214  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
215  "data inconsistency");
216  }
217  sense[ee] = data.dataOnEntities[MBEDGE][ee].getSense();
218  order[ee] = data.dataOnEntities[MBEDGE][ee].getDataOrder();
219  int nb_dofs = NBEDGE_H1(data.dataOnEntities[MBEDGE][ee].getDataOrder());
220  data.dataOnEntities[MBEDGE][ee].getN(base).resize(nb_gauss_pts, nb_dofs,
221  false);
222  data.dataOnEntities[MBEDGE][ee].getDiffN(base).resize(nb_gauss_pts,
223  2 * nb_dofs, false);
224  H1edgeN[ee] = &*data.dataOnEntities[MBEDGE][ee].getN(base).data().begin();
225  diffH1edgeN[ee] =
226  &*data.dataOnEntities[MBEDGE][ee].getDiffN(base).data().begin();
227  }
228  // shape functions on face 3
230  &sense[0], &order[0], &*N.data().begin(), &*diffN.data().begin(),
231  &H1edgeN[0], &diffH1edgeN[0], nb_gauss_pts, base_polynomials);
232  // shape functions on face 4
234  &sense[6], &order[6], &*N.data().begin(), &*diffN.data().begin(),
235  &H1edgeN[6], &diffH1edgeN[6], nb_gauss_pts, base_polynomials);
236  }
237 
238  // face
239  if (data.dataOnEntities[MBTRI].size() != 5) {
240  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
241  }
242  if ((data.spacesOnEntities[MBTRI]).test(H1)) {
243  for (int ff = 3; ff <= 4; ff++) {
244  int nb_dofs = NBFACETRI_H1(data.dataOnEntities[MBTRI][ff].getDataOrder());
245  data.dataOnEntities[MBTRI][ff].getN(base).resize(nb_gauss_pts, nb_dofs,
246  false);
247  data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(nb_gauss_pts,
248  2 * nb_dofs, false);
250  faceNodes[ff - 3], data.dataOnEntities[MBTRI][ff].getDataOrder(),
251  &*N.data().begin(), &*diffN.data().begin(),
252  &*data.dataOnEntities[MBTRI][ff].getN(base).data().begin(),
253  &*data.dataOnEntities[MBTRI][ff].getDiffN(base).data().begin(),
254  nb_gauss_pts, base_polynomials);
255  }
256  }
257 
259 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
#define NBEDGE_H1(P)
Numer of base function on edge for H1 space.
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
data structure for finite element entityIt keeps that about indices of degrees of freedom,...
PetscErrorCode H1_EdgeShapeFunctions_MBTRI(int *sense, int *p, double *N, double *diffN, double *edgeN[3], double *diff_edgeN[3], int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
H1_EdgeShapeFunctions_MBTRI.
Definition: h1.c:27
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
const FieldApproximationBase bAse
PetscErrorCode H1_FaceShapeFunctions_MBTRI(const int *face_nodes, int p, double *N, double *diffN, double *faceN, double *diff_faceN, int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition: h1.c:162
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
FieldApproximationBase
approximation base
Definition: definitions.h:143
DataForcesAndSourcesCore & dAta
#define NBFACETRI_H1(P)
Number of base function on triangle for H1 space.
#define CHKERR
Inline error check.
Definition: definitions.h:596
FlatPrismPolynomialBaseCtx * cTx
continuous field
Definition: definitions.h:171
const int order
approximation order

◆ getValueHcurl()

MoFEMErrorCode FlatPrismPolynomialBase::getValueHcurl ( MatrixDouble pts)
private

Definition at line 271 of file FlatPrismPolynomialBase.cpp.

271  {
272  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
273 }

◆ getValueHdiv()

MoFEMErrorCode FlatPrismPolynomialBase::getValueHdiv ( MatrixDouble pts)
private

Definition at line 267 of file FlatPrismPolynomialBase.cpp.

267  {
268  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
269 }

◆ getValueL2()

MoFEMErrorCode FlatPrismPolynomialBase::getValueL2 ( MatrixDouble pts)
private

Definition at line 261 of file FlatPrismPolynomialBase.cpp.

261  {
263  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
265 }
#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

◆ query_interface()

MoFEMErrorCode FlatPrismPolynomialBase::query_interface ( const MOFEMuuid uuid,
MoFEM::UnknownInterface **  iface 
) const
virtual

Reimplemented from MoFEM::BaseFunction.

Definition at line 49 of file FlatPrismPolynomialBase.cpp.

50  {
51 
53  *iface = NULL;
54  if (uuid == IDD_FLATPRISM_BASE_FUNCTION) {
55  *iface = const_cast<FlatPrismPolynomialBase *>(this);
57  } else {
58  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "wrong interference");
59  }
62 }
static const MOFEMuuid IDD_FLATPRISM_BASE_FUNCTION
#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 CHKERR
Inline error check.
Definition: definitions.h:596
MoFEMErrorCode query_interface(const MOFEMuuid &uuid, MoFEM::UnknownInterface **iface) const

Member Data Documentation

◆ connFace3

const EntityHandle* MoFEM::FlatPrismPolynomialBase::connFace3
private

Definition at line 83 of file FlatPrismPolynomialBase.hpp.

◆ connFace4

const EntityHandle* MoFEM::FlatPrismPolynomialBase::connFace4
private

Definition at line 84 of file FlatPrismPolynomialBase.hpp.

◆ connPrism

const EntityHandle* MoFEM::FlatPrismPolynomialBase::connPrism
private

Definition at line 82 of file FlatPrismPolynomialBase.hpp.

◆ cTx

FlatPrismPolynomialBaseCtx* MoFEM::FlatPrismPolynomialBase::cTx
private

Definition at line 71 of file FlatPrismPolynomialBase.hpp.

◆ diffN

MatrixDouble MoFEM::FlatPrismPolynomialBase::diffN
private

Definition at line 87 of file FlatPrismPolynomialBase.hpp.

◆ faceNodes

int MoFEM::FlatPrismPolynomialBase::faceNodes[2][3]
private

Definition at line 85 of file FlatPrismPolynomialBase.hpp.

◆ N

MatrixDouble MoFEM::FlatPrismPolynomialBase::N
private

Definition at line 86 of file FlatPrismPolynomialBase.hpp.

◆ numNodes

int MoFEM::FlatPrismPolynomialBase::numNodes
private

Definition at line 81 of file FlatPrismPolynomialBase.hpp.


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