v0.14.0
QuadPolynomialBase.cpp
Go to the documentation of this file.
1 /** \file QuadPolynomialBase.cpp
2 \brief Implementation of bases on a quad face
3 
4 */
5 
6 using namespace MoFEM;
7 
9 QuadPolynomialBase::query_interface(boost::typeindex::type_index type_index,
10  UnknownInterface **iface) const {
11  *iface = const_cast<QuadPolynomialBase *>(this);
12  return 0;
13 }
14 
17 
18  const FieldApproximationBase base = cTx->bAse;
19  EntitiesFieldData &data = cTx->dAta;
20  data.dataOnEntities[MBVERTEX][0].getNSharedPtr(base) =
21  data.dataOnEntities[MBVERTEX][0].getNSharedPtr(cTx->copyNodeBase);
22  data.dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(base) =
23  data.dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(cTx->copyNodeBase);
24 
25  switch (cTx->bAse) {
29  break;
32  break;
34  default:
35  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
36  }
37 
39 }
40 
43 
44  EntitiesFieldData &data = cTx->dAta;
45  const FieldApproximationBase base = cTx->bAse;
46  const auto copy_base = cTx->copyNodeBase;
47 
48  if (cTx->basePolynomialsType0 == NULL)
49  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
50  "Polynomial type not set");
51  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s, double *L,
52  double *diffL, const int dim) =
54 
55  auto &copy_base_fun = data.dataOnEntities[MBVERTEX][0].getN(copy_base);
56  auto &copy_diff_base_fun =
57  data.dataOnEntities[MBVERTEX][0].getDiffN(copy_base);
58 
59  int nb_gauss_pts = pts.size2();
60  auto &vert_dat = data.dataOnEntities[MBVERTEX][0];
61 
62  if (data.spacesOnEntities[MBEDGE].test(H1)) {
63  // edges
64  if (data.dataOnEntities[MBEDGE].size() != 4)
65  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
66  "should be four edges on quad");
67 
68  int sense[4], order[4];
69  double *H1edgeN[4], *diffH1edgeN[4];
70  for (int ee = 0; ee != 4; ++ee) {
71  auto &ent_dat = data.dataOnEntities[MBEDGE][ee];
72  if (ent_dat.getSense() == 0)
73  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "sense not set");
74 
75  sense[ee] = ent_dat.getSense();
76  order[ee] = ent_dat.getOrder();
77  int nb_dofs = NBEDGE_H1(ent_dat.getOrder());
78  ent_dat.getN(base).resize(nb_gauss_pts, nb_dofs, false);
79  ent_dat.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
80  H1edgeN[ee] = &*ent_dat.getN(base).data().begin();
81  diffH1edgeN[ee] = &*ent_dat.getDiffN(base).data().begin();
82  }
84  sense, order, &*copy_base_fun.data().begin(),
85  &*copy_diff_base_fun.data().begin(), H1edgeN, diffH1edgeN, nb_gauss_pts,
86  base_polynomials);
87  }
88 
89  if (data.spacesOnEntities[MBQUAD].test(H1)) {
90 
91  // face
92  if (data.dataOnEntities[MBQUAD].size() != 1)
93  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
94  "should be one quad to store bubble base on quad");
95 
96  auto &ent_dat = data.dataOnEntities[MBQUAD][0];
97  int nb_dofs = NBFACEQUAD_H1(ent_dat.getOrder());
98  ent_dat.getN(base).resize(nb_gauss_pts, nb_dofs, false);
99  ent_dat.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
100  int face_nodes[] = {0, 1, 2, 3};
102  face_nodes, ent_dat.getOrder(), &*vert_dat.getN(base).data().begin(),
103  &*vert_dat.getDiffN(base).data().begin(),
104  &*ent_dat.getN(base).data().begin(),
105  &*ent_dat.getDiffN(base).data().begin(), nb_gauss_pts,
106  base_polynomials);
107  }
108 
110 }
111 
114 
115  EntitiesFieldData &data = cTx->dAta;
116  const FieldApproximationBase base = cTx->bAse;
117  const auto copy_base = cTx->copyNodeBase;
118 
119  int nb_gauss_pts = pts.size2();
120  auto &copy_base_fun = data.dataOnEntities[MBVERTEX][0].getN(copy_base);
121  auto &copy_diff_base_fun =
122  data.dataOnEntities[MBVERTEX][0].getDiffN(copy_base);
123 
124  if (data.spacesOnEntities[MBEDGE].test(H1)) {
125  // edges
126  if (data.dataOnEntities[MBEDGE].size() != 4)
127  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
128  "should be four edges on quad");
129 
130  int sense[4], order[4];
131  double *H1edgeN[4], *diffH1edgeN[4];
132  for (int ee = 0; ee != 4; ++ee) {
133  auto &ent_dat = data.dataOnEntities[MBEDGE][ee];
134  if (ent_dat.getSense() == 0)
135  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "sense not set");
136 
137  sense[ee] = ent_dat.getSense();
138  order[ee] = ent_dat.getOrder();
139  int nb_dofs = NBEDGE_H1(ent_dat.getOrder());
140  ent_dat.getN(base).resize(nb_gauss_pts, nb_dofs, false);
141  ent_dat.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
142  H1edgeN[ee] = &*ent_dat.getN(base).data().begin();
143  diffH1edgeN[ee] = &*ent_dat.getDiffN(base).data().begin();
144  }
146  sense, order, &*copy_base_fun.data().begin(),
147  &*copy_diff_base_fun.data().begin(), H1edgeN, diffH1edgeN,
148  nb_gauss_pts);
149  }
150 
151  if (data.spacesOnEntities[MBQUAD].test(H1)) {
152 
153  // face
154  if (data.dataOnEntities[MBQUAD].size() != 1)
155  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
156  "should be one quad to store bubble base on quad");
157 
158  auto &ent_dat = data.dataOnEntities[MBQUAD][0];
159  int nb_dofs = NBFACEQUAD_H1(ent_dat.getOrder());
160  int p = ent_dat.getOrder();
161  int order[2] = {p, p};
162  ent_dat.getN(base).resize(nb_gauss_pts, nb_dofs, false);
163  ent_dat.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
164 
165  int face_nodes[] = {0, 1, 2, 3};
166  if (nb_dofs) {
168  face_nodes, order, &*copy_base_fun.data().begin(),
169  &*copy_diff_base_fun.data().begin(),
170  &*ent_dat.getN(base).data().begin(),
171  &*ent_dat.getDiffN(base).data().begin(), nb_gauss_pts);
172  }
173  }
174 
176 }
177 
180 
181  switch (cTx->bAse) {
184  break;
186  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
187  "Ainsworth not implemented on quad for L2 base");
188  break;
191  break;
193  default:
194  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
195  }
196 
198 }
199 
202 
203  EntitiesFieldData &data = cTx->dAta;
204  const FieldApproximationBase base = cTx->bAse;
205  const auto copy_base = cTx->copyNodeBase;
206 
207  int nb_gauss_pts = pts.size2();
208  auto &copy_base_fun = data.dataOnEntities[MBVERTEX][0].getN(copy_base);
209  auto &copy_diff_base_fun =
210  data.dataOnEntities[MBVERTEX][0].getDiffN(copy_base);
211 
212  auto &ent_dat = data.dataOnEntities[MBQUAD][0];
213  int p = ent_dat.getOrder();
214  int order[2] = {p, p};
215  int nb_dofs = NBFACEQUAD_L2(p);
216  ent_dat.getN(base).resize(nb_gauss_pts, nb_dofs, false);
217  ent_dat.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
218 
220  order, &*copy_base_fun.data().begin(),
221  &*copy_diff_base_fun.data().begin(), &*ent_dat.getN(base).data().begin(),
222  &*ent_dat.getDiffN(base).data().begin(), nb_gauss_pts);
223 
225 }
226 
229 
230  switch (cTx->bAse) {
233  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
234  "Ainsworth not implemented on quad for Hcurl base");
235  break;
238  break;
240  default:
241  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
242  }
243 
245 }
246 
250 
251  EntitiesFieldData &data = cTx->dAta;
252  const FieldApproximationBase base = cTx->bAse;
253  const auto copy_base = cTx->copyNodeBase;
254 
255  int nb_gauss_pts = pts.size2();
256  auto &copy_base_fun = data.dataOnEntities[MBVERTEX][0].getN(copy_base);
257  auto &copy_diff_base_fun =
258  data.dataOnEntities[MBVERTEX][0].getDiffN(copy_base);
259 
260  // Calculation H-curl on quad edges
261  if (data.spacesOnEntities[MBEDGE].test(HCURL)) {
262 
263  if (data.dataOnEntities[MBEDGE].size() != 4)
264  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
265  "wrong number of edges on quad, should be 4 but is %d",
266  data.dataOnEntities[MBEDGE].size());
267 
268  int sense[4], order[4];
269  double *hcurl_edge_n[4];
270  double *diff_hcurl_edge_n[4];
271 
272  for (int ee = 0; ee != 4; ++ee) {
273 
274  if (data.dataOnEntities[MBEDGE][ee].getSense() == 0)
275  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
276  "orientation (sense) of edge is not set");
277 
278  sense[ee] = data.dataOnEntities[MBEDGE][ee].getSense();
279  order[ee] = data.dataOnEntities[MBEDGE][ee].getOrder();
280  int nb_dofs =
281  NBEDGE_DEMKOWICZ_HCURL(data.dataOnEntities[MBEDGE][ee].getOrder());
282 
283  data.dataOnEntities[MBEDGE][ee].getN(base).resize(nb_gauss_pts,
284  3 * nb_dofs, false);
285  data.dataOnEntities[MBEDGE][ee].getDiffN(base).resize(
286  nb_gauss_pts, 3 * 2 * nb_dofs, false);
287 
288  hcurl_edge_n[ee] =
289  &*data.dataOnEntities[MBEDGE][ee].getN(base).data().begin();
290  diff_hcurl_edge_n[ee] =
291  &*data.dataOnEntities[MBEDGE][ee].getDiffN(base).data().begin();
292  }
294  sense, order, &*copy_base_fun.data().begin(),
295  &*copy_diff_base_fun.data().begin(), hcurl_edge_n, diff_hcurl_edge_n,
296  nb_gauss_pts);
297  }
298 
299  if (data.spacesOnEntities[MBQUAD].test(HCURL)) {
300 
301  // face
302  if (data.dataOnEntities[MBQUAD].size() != 1)
303  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
304  "No data struture to keep base functions on face");
305 
306  int p = data.dataOnEntities[MBQUAD][0].getOrder();
307  const int nb_dofs_family = NBFACEQUAD_DEMKOWICZ_FAMILY_HCURL(p, p);
308  if (nb_dofs_family) {
309  faceFamily.resize(2, 3 * nb_dofs_family * nb_gauss_pts, false);
310  diffFaceFamily.resize(2, 6 * nb_dofs_family * nb_gauss_pts, false);
311 
312  int order[2] = {p, p};
313  double *face_family_ptr[] = {&faceFamily(0, 0), &faceFamily(1, 0)};
314  double *diff_face_family_ptr[] = {&diffFaceFamily(0, 0),
315  &diffFaceFamily(1, 0)};
316  int face_nodes[] = {0, 1, 2, 3};
318  face_nodes, order, &*copy_base_fun.data().begin(),
319  &*copy_diff_base_fun.data().begin(), face_family_ptr,
320  diff_face_family_ptr, nb_gauss_pts);
321  }
322 
323  // put family back
324 
325  const int nb_dofs = NBFACEQUAD_DEMKOWICZ_HCURL(p);
326  auto &face_n = data.dataOnEntities[MBQUAD][0].getN(base);
327  auto &diff_face_n = data.dataOnEntities[MBQUAD][0].getDiffN(base);
328  face_n.resize(nb_gauss_pts, 3 * nb_dofs, false);
329  diff_face_n.resize(nb_gauss_pts, 3 * 2 * nb_dofs, false);
330 
331  if (nb_dofs) {
332  double *ptr_f0 = &faceFamily(0, 0);
333  double *ptr_f1 = &faceFamily(1, 0);
334  double *ptr = &face_n(0, 0);
335  for (int n = 0; n != faceFamily.size2() / 3; ++n) {
336  for (int j = 0; j != 3; ++j) {
337  *ptr = *ptr_f0;
338  ++ptr;
339  ++ptr_f0;
340  }
341  for (int j = 0; j != 3; ++j) {
342  *ptr = *ptr_f1;
343  ++ptr;
344  ++ptr_f1;
345  }
346  }
347 
348  double *diff_ptr_f0 = &diffFaceFamily(0, 0);
349  double *diff_ptr_f1 = &diffFaceFamily(1, 0);
350  double *diff_ptr = &diff_face_n(0, 0);
351  for (int n = 0; n != diffFaceFamily.size2() / 6; ++n) {
352  for (int j = 0; j != 6; ++j) {
353  *diff_ptr = *diff_ptr_f0;
354  ++diff_ptr;
355  ++diff_ptr_f0;
356  }
357  for (int j = 0; j != 6; ++j) {
358  *diff_ptr = *diff_ptr_f1;
359  ++diff_ptr;
360  ++diff_ptr_f1;
361  }
362  }
363  }
364  }
366 }
367 
370 
371  switch (cTx->bAse) {
374  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
375  "Ainsworth not implemented on quad for Hdiv base");
376  break;
379  break;
381  default:
382  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
383  }
384 
386 }
387 
391 
392  EntitiesFieldData &data = cTx->dAta;
393  const FieldApproximationBase base = cTx->bAse;
394  const auto copy_base = cTx->copyNodeBase;
395  int nb_gauss_pts = pts.size2();
396 
397  auto &copy_base_fun = data.dataOnEntities[MBVERTEX][0].getN(copy_base);
398  auto &copy_diff_base_fun =
399  data.dataOnEntities[MBVERTEX][0].getDiffN(copy_base);
400 
401  if (data.spacesOnEntities[MBQUAD].test(HDIV)) {
402 
403  // face
404  if (data.dataOnEntities[MBQUAD].size() != 1)
405  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
406  "No data struture to keep base functions on face");
407 
408  int p = data.dataOnEntities[MBQUAD][0].getOrder();
409  const int nb_dofs = NBFACEQUAD_DEMKOWICZ_HDIV(p);
410  auto &face_n = data.dataOnEntities[MBQUAD][0].getN(base);
411  auto &diff_face_n = data.dataOnEntities[MBQUAD][0].getDiffN(base);
412  face_n.resize(nb_gauss_pts, 3 * nb_dofs, false);
413  diff_face_n.resize(nb_gauss_pts, 6 * nb_dofs, false);
414 
415  if (nb_dofs) {
416 
417  std::array<int, 2> order = {p, p};
418  std::array<int, 6> face_nodes = {0, 1, 2, 3};
420  face_nodes.data(), order.data(), &*copy_base_fun.data().begin(),
421  &*copy_diff_base_fun.data().begin(), &*face_n.data().begin(),
422  &*diff_face_n.data().begin(), nb_gauss_pts);
423  }
424  }
425 
427 }
428 
431  boost::shared_ptr<BaseFunctionCtx> ctx_ptr) {
433 
434  cTx = ctx_ptr->getInterface<EntPolynomialBaseCtx>();
435 
436  int nb_gauss_pts = pts.size2();
437  if (!nb_gauss_pts)
439 
440  if (pts.size1() < 2)
441  SETERRQ(
442  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
443  "Wrong dimension of pts, should be at least 3 rows with coordinates");
444 
445  const FieldApproximationBase base = cTx->bAse;
446  EntitiesFieldData &data = cTx->dAta;
447  if (cTx->copyNodeBase != NOBASE)
448  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
449  "Shape base has to be on NOBASE", ApproximationBaseNames[base]);
450 
451  auto &base_shape = data.dataOnEntities[MBVERTEX][0].getN(cTx->copyNodeBase);
452  auto &diff_base =
453  data.dataOnEntities[MBVERTEX][0].getDiffN(cTx->copyNodeBase);
454 
455  if (base_shape.size1() != pts.size2())
456  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
457  "Number of base functions integration points not equal number of "
458  "set integration point");
459  if (base_shape.size2() != 4)
460  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
461  "Number of shape functions should be four");
462  if (diff_base.size1() != pts.size2())
463  SETERRQ2(
464  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
465  "Number of diff base functions integration points not equal number of "
466  "set integration point %d != %d",
467  diff_base.size1(), pts.size2());
468  if (diff_base.size2() != 8)
469  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
470  "Number of shape functions should be four");
471 
472  switch (cTx->spaceContinuity) {
473  case CONTINUOUS:
474  switch (cTx->sPace) {
475  case H1:
476  CHKERR getValueH1(pts);
477  break;
478  case L2:
479  CHKERR getValueL2(pts);
480  break;
481  case HCURL:
482  CHKERR getValueHcurl(pts);
483  break;
484  case HDIV:
485  CHKERR getValueHdiv(pts);
486  break;
487  default:
488  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
489  }
490 
491  break;
492 
493  default:
494  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Unknown continuity");
495  }
496 
498 }
UBlasMatrix< double >
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
H1
@ H1
continuous field
Definition: definitions.h:85
NBEDGE_H1
#define NBEDGE_H1(P)
Number of base function on edge for H1 space.
Definition: h1_hdiv_hcurl_l2.h:55
MoFEM::QuadPolynomialBase::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: QuadPolynomialBase.cpp:9
H1_EdgeShapeFunctions_MBQUAD
PetscErrorCode H1_EdgeShapeFunctions_MBQUAD(int *sense, int *p, double *N, double *diffN, double *edgeN[4], double *diff_edgeN[4], int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition: h1.c:1091
MoFEM::EntPolynomialBaseCtx
Class used to pass element data to calculate base functions on tet,triangle,edge.
Definition: EntPolynomialBaseCtx.hpp:22
MoFEM::QuadPolynomialBase::faceFamily
MatrixDouble faceFamily
Definition: QuadPolynomialBase.hpp:45
NOBASE
@ NOBASE
Definition: definitions.h:59
MoFEM::QuadPolynomialBase::getValueL2DemkowiczBase
MoFEMErrorCode getValueL2DemkowiczBase(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:200
MoFEM::DemkowiczHexAndQuad::L2_FaceShapeFunctions_ONQUAD
MoFEMErrorCode L2_FaceShapeFunctions_ONQUAD(int *p, double *N, double *diffN, double *face_buble, double *diff_face_bubble, int nb_integration_pts)
L2 Face base functions on Quad.
Definition: EdgeQuadHexPolynomials.cpp:335
ApproximationBaseNames
const static char *const ApproximationBaseNames[]
Definition: definitions.h:72
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
H1_QuadShapeFunctions_MBQUAD
PetscErrorCode H1_QuadShapeFunctions_MBQUAD(int *faces_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:959
MoFEM::EntPolynomialBaseCtx::dAta
EntitiesFieldData & dAta
Definition: EntPolynomialBaseCtx.hpp:36
MoFEM::QuadPolynomialBase::diffFaceFamily
MatrixDouble diffFaceFamily
Definition: QuadPolynomialBase.hpp:46
MoFEM::QuadPolynomialBase::getValueHdiv
MoFEMErrorCode getValueHdiv(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:368
NBFACEQUAD_DEMKOWICZ_HCURL
#define NBFACEQUAD_DEMKOWICZ_HCURL(P)
Definition: h1_hdiv_hcurl_l2.h:118
MoFEM::EntitiesFieldData::spacesOnEntities
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
Definition: EntitiesFieldData.hpp:49
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::QuadPolynomialBase::getValueH1AinsworthBase
MoFEMErrorCode getValueH1AinsworthBase(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:41
MoFEM::DemkowiczHexAndQuad::H1_EdgeShapeFunctions_ONQUAD
MoFEMErrorCode H1_EdgeShapeFunctions_ONQUAD(int *sense, int *p, double *N, double *diffN, double *edgeN[4], double *edgeDiffN[4], int nb_integration_pts)
H1 Edge base functions on Quad.
Definition: EdgeQuadHexPolynomials.cpp:202
MoFEM::QuadPolynomialBase::getValueHcurlDemkowiczBase
MoFEMErrorCode getValueHcurlDemkowiczBase(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:248
MoFEM::EntPolynomialBaseCtx::copyNodeBase
const FieldApproximationBase copyNodeBase
Definition: EntPolynomialBaseCtx.hpp:41
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::QuadPolynomialBase::getValueH1
MoFEMErrorCode getValueH1(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:15
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
NBFACEQUAD_L2
#define NBFACEQUAD_L2(P)
Number of base functions on quad for L2 space.
Definition: h1_hdiv_hcurl_l2.h:70
MoFEM::QuadPolynomialBase::getValueH1DemkowiczBase
MoFEMErrorCode getValueH1DemkowiczBase(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:112
MoFEM::DemkowiczHexAndQuad::Hdiv_FaceShapeFunctions_ONQUAD
MoFEMErrorCode Hdiv_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N, double *diffN, double *faceN, double *diff_faceN, int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:590
MoFEM::EntPolynomialBaseCtx::sPace
const FieldSpace sPace
Definition: EntPolynomialBaseCtx.hpp:37
NBEDGE_DEMKOWICZ_HCURL
#define NBEDGE_DEMKOWICZ_HCURL(P)
Definition: h1_hdiv_hcurl_l2.h:108
NBFACEQUAD_DEMKOWICZ_FAMILY_HCURL
#define NBFACEQUAD_DEMKOWICZ_FAMILY_HCURL(P, Q)
Number of base functions on quad for Hcurl space.
Definition: h1_hdiv_hcurl_l2.h:116
MoFEM::QuadPolynomialBase::getValue
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Definition: QuadPolynomialBase.cpp:430
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
AINSWORTH_LOBATTO_BASE
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:62
MoFEM::QuadPolynomialBase::getValueHcurl
MoFEMErrorCode getValueHcurl(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:227
NBFACEQUAD_H1
#define NBFACEQUAD_H1(P)
Number of base functions on quad for H1 space.
Definition: h1_hdiv_hcurl_l2.h:65
MoFEM::DemkowiczHexAndQuad::H1_FaceShapeFunctions_ONQUAD
MoFEMErrorCode H1_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N, double *diffN, double *faceN, double *diff_faceN, int nb_integration_pts)
H1 Face bubble functions on Quad.
Definition: EdgeQuadHexPolynomials.cpp:275
MoFEM::QuadPolynomialBase
Calculate base functions on triangle.
Definition: QuadPolynomialBase.hpp:21
MoFEM::EntPolynomialBaseCtx::bAse
const FieldApproximationBase bAse
Definition: EntPolynomialBaseCtx.hpp:39
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
convert.n
n
Definition: convert.py:82
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
MoFEM::EntPolynomialBaseCtx::spaceContinuity
const FieldContinuity spaceContinuity
Definition: EntPolynomialBaseCtx.hpp:38
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::QuadPolynomialBase::getValueL2
MoFEMErrorCode getValueL2(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:178
NBFACEQUAD_DEMKOWICZ_HDIV
#define NBFACEQUAD_DEMKOWICZ_HDIV(P)
Definition: h1_hdiv_hcurl_l2.h:145
MoFEM::DemkowiczHexAndQuad::Hcurl_FaceShapeFunctions_ONQUAD
MoFEMErrorCode Hcurl_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N, double *diffN, double *faceN[], double *diff_faceN[], int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:484
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::QuadPolynomialBase::getValueHdivDemkowiczBase
MoFEMErrorCode getValueHdivDemkowiczBase(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:389
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::QuadPolynomialBase::cTx
EntPolynomialBaseCtx * cTx
Definition: QuadPolynomialBase.hpp:32
MoFEM::EntitiesFieldData::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: EntitiesFieldData.hpp:57
MoFEM::DemkowiczHexAndQuad::Hcurl_EdgeShapeFunctions_ONQUAD
MoFEMErrorCode Hcurl_EdgeShapeFunctions_ONQUAD(int *sense, int *p, double *N, double *diffN, double *edgeN[4], double *curl_edgeN[4], int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:393
MoFEM::EntitiesFieldData
data structure for finite element entity
Definition: EntitiesFieldData.hpp:40
CONTINUOUS
@ CONTINUOUS
Regular field.
Definition: definitions.h:100
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::EntPolynomialBaseCtx::basePolynomialsType0
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Definition: EntPolynomialBaseCtx.hpp:27