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

Calculate base functions on triangle. More...

#include <src/approximation/QuadPolynomialBase.hpp>

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

Public Member Functions

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

Private Member Functions

MoFEMErrorCode getValueH1 (MatrixDouble &pts)
 
MoFEMErrorCode getValueL2 (MatrixDouble &pts)
 
MoFEMErrorCode getValueHcurl (MatrixDouble &pts)
 
MoFEMErrorCode getValueHdiv (MatrixDouble &pts)
 
MoFEMErrorCode getValueH1AinsworthBase (MatrixDouble &pts)
 
MoFEMErrorCode getValueH1DemkowiczBase (MatrixDouble &pts)
 
MoFEMErrorCode getValueL2DemkowiczBase (MatrixDouble &pts)
 
MoFEMErrorCode getValueHcurlDemkowiczBase (MatrixDouble &pts)
 
MoFEMErrorCode getValueHdivDemkowiczBase (MatrixDouble &pts)
 

Private Attributes

EntPolynomialBaseCtxcTx
 
MatrixDouble faceFamily
 
MatrixDouble diffFaceFamily
 

Additional Inherited Members

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

Detailed Description

Calculate base functions on triangle.

Definition at line 21 of file QuadPolynomialBase.hpp.

Constructor & Destructor Documentation

◆ QuadPolynomialBase()

MoFEM::QuadPolynomialBase::QuadPolynomialBase ( )
default

◆ ~QuadPolynomialBase()

MoFEM::QuadPolynomialBase::~QuadPolynomialBase ( )
default

Member Function Documentation

◆ getValue()

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

Reimplemented from MoFEM::BaseFunction.

Definition at line 434 of file QuadPolynomialBase.cpp.

435  {
437 
438  cTx = ctx_ptr->getInterface<EntPolynomialBaseCtx>();
439 
440  int nb_gauss_pts = pts.size2();
441  if (!nb_gauss_pts)
443 
444  if (pts.size1() < 2)
445  SETERRQ(
446  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
447  "Wrong dimension of pts, should be at least 3 rows with coordinates");
448 
449  const FieldApproximationBase base = cTx->bAse;
450  EntitiesFieldData &data = cTx->dAta;
451  if (cTx->copyNodeBase != NOBASE)
452  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
453  "Shape base has to be on NOBASE", ApproximationBaseNames[base]);
454 
455  auto &base_shape = data.dataOnEntities[MBVERTEX][0].getN(cTx->copyNodeBase);
456  auto &diff_base =
457  data.dataOnEntities[MBVERTEX][0].getDiffN(cTx->copyNodeBase);
458 
459  if (base_shape.size1() != pts.size2())
460  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
461  "Number of base functions integration points not equal number of "
462  "set integration point");
463  if (base_shape.size2() != 4)
464  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
465  "Number of shape functions should be four");
466  if (diff_base.size1() != pts.size2())
467  SETERRQ2(
468  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
469  "Number of diff base functions integration points not equal number of "
470  "set integration point %d != %d",
471  diff_base.size1(), pts.size2());
472  if (diff_base.size2() != 8)
473  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
474  "Number of shape functions should be four");
475 
476  switch (cTx->sPace) {
477  case H1:
478  CHKERR getValueH1(pts);
479  break;
480  case L2:
481  CHKERR getValueL2(pts);
482  break;
483  case HCURL:
484  CHKERR getValueHcurl(pts);
485  break;
486  case HDIV:
487  CHKERR getValueHdiv(pts);
488  break;
489  default:
490  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
491  }
492 
494 }

◆ getValueH1()

MoFEMErrorCode QuadPolynomialBase::getValueH1 ( MatrixDouble pts)
private

Definition at line 17 of file QuadPolynomialBase.cpp.

17  {
19 
20  const FieldApproximationBase base = cTx->bAse;
21  EntitiesFieldData &data = cTx->dAta;
22  data.dataOnEntities[MBVERTEX][0].getNSharedPtr(base) =
23  data.dataOnEntities[MBVERTEX][0].getNSharedPtr(cTx->copyNodeBase);
24  data.dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(base) =
25  data.dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(cTx->copyNodeBase);
26 
27  switch (cTx->bAse) {
31  break;
34  break;
36  default:
37  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
38  }
39 
41 }

◆ getValueH1AinsworthBase()

MoFEMErrorCode QuadPolynomialBase::getValueH1AinsworthBase ( MatrixDouble pts)
private

Definition at line 43 of file QuadPolynomialBase.cpp.

43  {
45 
46  EntitiesFieldData &data = cTx->dAta;
47  const FieldApproximationBase base = cTx->bAse;
48  const auto copy_base = cTx->copyNodeBase;
49 
50  if (cTx->basePolynomialsType0 == NULL)
51  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
52  "Polynomial type not set");
53  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s, double *L,
54  double *diffL, const int dim) =
56 
57  auto &copy_base_fun = data.dataOnEntities[MBVERTEX][0].getN(copy_base);
58  auto &copy_diff_base_fun =
59  data.dataOnEntities[MBVERTEX][0].getDiffN(copy_base);
60 
61  int nb_gauss_pts = pts.size2();
62  auto &vert_dat = data.dataOnEntities[MBVERTEX][0];
63 
64  if (data.spacesOnEntities[MBEDGE].test(H1)) {
65  // edges
66  if (data.dataOnEntities[MBEDGE].size() != 4)
67  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
68  "should be four edges on quad");
69 
70  int sense[4], order[4];
71  double *H1edgeN[4], *diffH1edgeN[4];
72  for (int ee = 0; ee != 4; ++ee) {
73  auto &ent_dat = data.dataOnEntities[MBEDGE][ee];
74  if (ent_dat.getSense() == 0)
75  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "sense not set");
76 
77  sense[ee] = ent_dat.getSense();
78  order[ee] = ent_dat.getOrder();
79  int nb_dofs = NBEDGE_H1(ent_dat.getOrder());
80  ent_dat.getN(base).resize(nb_gauss_pts, nb_dofs, false);
81  ent_dat.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
82  H1edgeN[ee] = &*ent_dat.getN(base).data().begin();
83  diffH1edgeN[ee] = &*ent_dat.getDiffN(base).data().begin();
84  }
86  sense, order, &*copy_base_fun.data().begin(),
87  &*copy_diff_base_fun.data().begin(), H1edgeN, diffH1edgeN, nb_gauss_pts,
88  base_polynomials);
89  }
90 
91  if (data.spacesOnEntities[MBQUAD].test(H1)) {
92 
93  // face
94  if (data.dataOnEntities[MBQUAD].size() != 1)
95  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
96  "should be one quad to store bubble base on quad");
97 
98  auto &ent_dat = data.dataOnEntities[MBQUAD][0];
99  int nb_dofs = NBFACEQUAD_H1(ent_dat.getOrder());
100  ent_dat.getN(base).resize(nb_gauss_pts, nb_dofs, false);
101  ent_dat.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
102  int face_nodes[] = {0, 1, 2, 3};
104  face_nodes, ent_dat.getOrder(),
105  &*vert_dat.getN(base).data().begin(),
106  &*vert_dat.getDiffN(base).data().begin(),
107  &*ent_dat.getN(base).data().begin(),
108  &*ent_dat.getDiffN(base).data().begin(), nb_gauss_pts,
109  base_polynomials);
110  }
111 
113 }

◆ getValueH1DemkowiczBase()

MoFEMErrorCode QuadPolynomialBase::getValueH1DemkowiczBase ( MatrixDouble pts)
private

Definition at line 115 of file QuadPolynomialBase.cpp.

115  {
117 
118  EntitiesFieldData &data = cTx->dAta;
119  const FieldApproximationBase base = cTx->bAse;
120  const auto copy_base = cTx->copyNodeBase;
121 
122  int nb_gauss_pts = pts.size2();
123  auto &copy_base_fun = data.dataOnEntities[MBVERTEX][0].getN(copy_base);
124  auto &copy_diff_base_fun =
125  data.dataOnEntities[MBVERTEX][0].getDiffN(copy_base);
126 
127  if (data.spacesOnEntities[MBEDGE].test(H1)) {
128  // edges
129  if (data.dataOnEntities[MBEDGE].size() != 4)
130  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
131  "should be four edges on quad");
132 
133  int sense[4], order[4];
134  double *H1edgeN[4], *diffH1edgeN[4];
135  for (int ee = 0; ee != 4; ++ee) {
136  auto &ent_dat = data.dataOnEntities[MBEDGE][ee];
137  if (ent_dat.getSense() == 0)
138  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "sense not set");
139 
140  sense[ee] = ent_dat.getSense();
141  order[ee] = ent_dat.getOrder();
142  int nb_dofs = NBEDGE_H1(ent_dat.getOrder());
143  ent_dat.getN(base).resize(nb_gauss_pts, nb_dofs, false);
144  ent_dat.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
145  H1edgeN[ee] = &*ent_dat.getN(base).data().begin();
146  diffH1edgeN[ee] = &*ent_dat.getDiffN(base).data().begin();
147  }
149  sense, order, &*copy_base_fun.data().begin(),
150  &*copy_diff_base_fun.data().begin(), H1edgeN, diffH1edgeN,
151  nb_gauss_pts);
152  }
153 
154  if (data.spacesOnEntities[MBQUAD].test(H1)) {
155 
156  // face
157  if (data.dataOnEntities[MBQUAD].size() != 1)
158  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
159  "should be one quad to store bubble base on quad");
160 
161  auto &ent_dat = data.dataOnEntities[MBQUAD][0];
162  int nb_dofs = NBFACEQUAD_H1(ent_dat.getOrder());
163  int p = ent_dat.getOrder();
164  int order[2] = {p, p};
165  ent_dat.getN(base).resize(nb_gauss_pts, nb_dofs, false);
166  ent_dat.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
167 
168  int face_nodes[] = {0, 1, 2, 3};
169  if (nb_dofs) {
171  face_nodes, order, &*copy_base_fun.data().begin(),
172  &*copy_diff_base_fun.data().begin(),
173  &*ent_dat.getN(base).data().begin(),
174  &*ent_dat.getDiffN(base).data().begin(), nb_gauss_pts);
175  }
176  }
177 
179 }

◆ getValueHcurl()

MoFEMErrorCode QuadPolynomialBase::getValueHcurl ( MatrixDouble pts)
private

Definition at line 230 of file QuadPolynomialBase.cpp.

230  {
232 
233  switch (cTx->bAse) {
236  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
237  "Ainsworth not implemented on quad for Hcurl base");
238  break;
241  break;
243  default:
244  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
245  }
246 
248 }

◆ getValueHcurlDemkowiczBase()

MoFEMErrorCode QuadPolynomialBase::getValueHcurlDemkowiczBase ( MatrixDouble pts)
private

Definition at line 251 of file QuadPolynomialBase.cpp.

251  {
253 
254  EntitiesFieldData &data = cTx->dAta;
255  const FieldApproximationBase base = cTx->bAse;
256  const auto copy_base = cTx->copyNodeBase;
257 
258  int nb_gauss_pts = pts.size2();
259  auto &copy_base_fun = data.dataOnEntities[MBVERTEX][0].getN(copy_base);
260  auto &copy_diff_base_fun =
261  data.dataOnEntities[MBVERTEX][0].getDiffN(copy_base);
262 
263  // Calculation H-curl on quad edges
264  if (data.spacesOnEntities[MBEDGE].test(HCURL)) {
265 
266  if (data.dataOnEntities[MBEDGE].size() != 4)
267  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
268  "wrong number of edges on quad, should be 4 but is %d",
269  data.dataOnEntities[MBEDGE].size());
270 
271  int sense[4], order[4];
272  double *hcurl_edge_n[4];
273  double *diff_hcurl_edge_n[4];
274 
275  for (int ee = 0; ee != 4; ++ee) {
276 
277  if (data.dataOnEntities[MBEDGE][ee].getSense() == 0)
278  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
279  "orientation (sense) of edge is not set");
280 
281  sense[ee] = data.dataOnEntities[MBEDGE][ee].getSense();
282  order[ee] = data.dataOnEntities[MBEDGE][ee].getOrder();
283  int nb_dofs = NBEDGE_DEMKOWICZ_HCURL(
284  data.dataOnEntities[MBEDGE][ee].getOrder());
285 
286  data.dataOnEntities[MBEDGE][ee].getN(base).resize(nb_gauss_pts,
287  3 * nb_dofs, false);
288  data.dataOnEntities[MBEDGE][ee].getDiffN(base).resize(
289  nb_gauss_pts, 3 * 2 * nb_dofs, false);
290 
291  hcurl_edge_n[ee] =
292  &*data.dataOnEntities[MBEDGE][ee].getN(base).data().begin();
293  diff_hcurl_edge_n[ee] =
294  &*data.dataOnEntities[MBEDGE][ee].getDiffN(base).data().begin();
295  }
297  sense, order, &*copy_base_fun.data().begin(),
298  &*copy_diff_base_fun.data().begin(), hcurl_edge_n, diff_hcurl_edge_n,
299  nb_gauss_pts);
300  }
301 
302  if (data.spacesOnEntities[MBQUAD].test(HCURL)) {
303 
304  // face
305  if (data.dataOnEntities[MBQUAD].size() != 1)
306  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
307  "No data struture to keep base functions on face");
308 
309  int p = data.dataOnEntities[MBQUAD][0].getOrder();
310  const int nb_dofs_family = NBFACEQUAD_DEMKOWICZ_FAMILY_HCURL(p, p);
311  if (nb_dofs_family) {
312  faceFamily.resize(2, 3 * nb_dofs_family * nb_gauss_pts, false);
313  diffFaceFamily.resize(2, 6 * nb_dofs_family * nb_gauss_pts, false);
314 
315  int order[2] = {p, p};
316  double *face_family_ptr[] = {&faceFamily(0, 0), &faceFamily(1, 0)};
317  double *diff_face_family_ptr[] = {&diffFaceFamily(0, 0),
318  &diffFaceFamily(1, 0)};
319  int face_nodes[] = {0, 1, 2, 3};
321  face_nodes, order, &*copy_base_fun.data().begin(),
322  &*copy_diff_base_fun.data().begin(), face_family_ptr,
323  diff_face_family_ptr, nb_gauss_pts);
324  }
325 
326  // put family back
327 
328  const int nb_dofs = NBFACEQUAD_DEMKOWICZ_HCURL(p);
329  auto &face_n = data.dataOnEntities[MBQUAD][0].getN(base);
330  auto &diff_face_n = data.dataOnEntities[MBQUAD][0].getDiffN(base);
331  face_n.resize(nb_gauss_pts, 3 * nb_dofs, false);
332  diff_face_n.resize(nb_gauss_pts, 3 * 2 * nb_dofs, false);
333 
334  if (nb_dofs) {
335  double *ptr_f0 = &faceFamily(0, 0);
336  double *ptr_f1 = &faceFamily(1, 0);
337  double *ptr = &face_n(0, 0);
338  for (int n = 0; n != faceFamily.size2() / 3; ++n) {
339  for (int j = 0; j != 3; ++j) {
340  *ptr = *ptr_f0;
341  ++ptr;
342  ++ptr_f0;
343  }
344  for (int j = 0; j != 3; ++j) {
345  *ptr = *ptr_f1;
346  ++ptr;
347  ++ptr_f1;
348  }
349  }
350 
351  double *diff_ptr_f0 = &diffFaceFamily(0, 0);
352  double *diff_ptr_f1 = &diffFaceFamily(1, 0);
353  double *diff_ptr = &diff_face_n(0, 0);
354  for (int n = 0; n != diffFaceFamily.size2() / 6; ++n) {
355  for (int j = 0; j != 6; ++j) {
356  *diff_ptr = *diff_ptr_f0;
357  ++diff_ptr;
358  ++diff_ptr_f0;
359  }
360  for (int j = 0; j != 6; ++j) {
361  *diff_ptr = *diff_ptr_f1;
362  ++diff_ptr;
363  ++diff_ptr_f1;
364  }
365  }
366  }
367  }
369 }

◆ getValueHdiv()

MoFEMErrorCode QuadPolynomialBase::getValueHdiv ( MatrixDouble pts)
private

Definition at line 371 of file QuadPolynomialBase.cpp.

371  {
373 
374  switch (cTx->bAse) {
377  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
378  "Ainsworth not implemented on quad for Hdiv base");
379  break;
382  break;
384  default:
385  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
386  }
387 
389 }

◆ getValueHdivDemkowiczBase()

MoFEMErrorCode QuadPolynomialBase::getValueHdivDemkowiczBase ( MatrixDouble pts)
private

Definition at line 392 of file QuadPolynomialBase.cpp.

392  {
394 
395  EntitiesFieldData &data = cTx->dAta;
396  const FieldApproximationBase base = cTx->bAse;
397  const auto copy_base = cTx->copyNodeBase;
398  int nb_gauss_pts = pts.size2();
399 
400  auto &copy_base_fun = data.dataOnEntities[MBVERTEX][0].getN(copy_base);
401  auto &copy_diff_base_fun =
402  data.dataOnEntities[MBVERTEX][0].getDiffN(copy_base);
403 
404  if (data.spacesOnEntities[MBQUAD].test(HDIV)) {
405 
406  // face
407  if (data.dataOnEntities[MBQUAD].size() != 1)
408  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
409  "No data struture to keep base functions on face");
410 
411  int p = data.dataOnEntities[MBQUAD][0].getOrder();
412  const int nb_dofs = NBFACEQUAD_DEMKOWICZ_HDIV(p);
413  auto &face_n = data.dataOnEntities[MBQUAD][0].getN(base);
414  auto &diff_face_n = data.dataOnEntities[MBQUAD][0].getDiffN(base);
415  face_n.resize(nb_gauss_pts, 3 * nb_dofs, false);
416  diff_face_n.resize(nb_gauss_pts, 6 * nb_dofs, false);
417 
418  if (nb_dofs) {
419 
420  std::array<int, 2> order = {p, p};
421  std::array<int, 6> face_nodes = {0, 1, 2, 3};
423  face_nodes.data(), order.data(), &*copy_base_fun.data().begin(),
424  &*copy_diff_base_fun.data().begin(), &*face_n.data().begin(),
425  &*diff_face_n.data().begin(), nb_gauss_pts);
426 
427  }
428  }
429 
431 }

◆ getValueL2()

MoFEMErrorCode QuadPolynomialBase::getValueL2 ( MatrixDouble pts)
private

Definition at line 181 of file QuadPolynomialBase.cpp.

181  {
183 
184  switch (cTx->bAse) {
187  break;
189  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
190  "Ainsworth not implemented on quad for L2 base");
191  break;
194  break;
196  default:
197  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
198  }
199 
201 }

◆ getValueL2DemkowiczBase()

MoFEMErrorCode QuadPolynomialBase::getValueL2DemkowiczBase ( MatrixDouble pts)
private

Definition at line 203 of file QuadPolynomialBase.cpp.

203  {
205 
206  EntitiesFieldData &data = cTx->dAta;
207  const FieldApproximationBase base = cTx->bAse;
208  const auto copy_base = cTx->copyNodeBase;
209 
210  int nb_gauss_pts = pts.size2();
211  auto &copy_base_fun = data.dataOnEntities[MBVERTEX][0].getN(copy_base);
212  auto &copy_diff_base_fun =
213  data.dataOnEntities[MBVERTEX][0].getDiffN(copy_base);
214 
215  auto &ent_dat = data.dataOnEntities[MBQUAD][0];
216  int p = ent_dat.getOrder();
217  int order[2] = {p, p};
218  int nb_dofs = NBFACEQUAD_L2(p);
219  ent_dat.getN(base).resize(nb_gauss_pts, nb_dofs, false);
220  ent_dat.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
221 
223  order, &*copy_base_fun.data().begin(),
224  &*copy_diff_base_fun.data().begin(), &*ent_dat.getN(base).data().begin(),
225  &*ent_dat.getDiffN(base).data().begin(), nb_gauss_pts);
226 
228 }

◆ query_interface()

MoFEMErrorCode QuadPolynomialBase::query_interface ( boost::typeindex::type_index  type_index,
UnknownInterface **  iface 
) const
virtual

Reimplemented from MoFEM::BaseFunction.

Definition at line 11 of file QuadPolynomialBase.cpp.

12  {
13  *iface = const_cast<QuadPolynomialBase *>(this);
14  return 0;
15 }

Member Data Documentation

◆ cTx

EntPolynomialBaseCtx* MoFEM::QuadPolynomialBase::cTx
private

Definition at line 32 of file QuadPolynomialBase.hpp.

◆ diffFaceFamily

MatrixDouble MoFEM::QuadPolynomialBase::diffFaceFamily
private

Definition at line 46 of file QuadPolynomialBase.hpp.

◆ faceFamily

MatrixDouble MoFEM::QuadPolynomialBase::faceFamily
private

Definition at line 45 of file QuadPolynomialBase.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:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
H1
@ H1
continuous field
Definition: definitions.h:85
NBEDGE_H1
#define NBEDGE_H1(P)
Numer of base function on edge for H1 space.
Definition: h1_hdiv_hcurl_l2.h:55
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
MoFEM::QuadPolynomialBase::getValueL2DemkowiczBase
MoFEMErrorCode getValueL2DemkowiczBase(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:203
NOBASE
@ NOBASE
Definition: definitions.h:59
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
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:371
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:50
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::QuadPolynomialBase::getValueH1AinsworthBase
MoFEMErrorCode getValueH1AinsworthBase(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:43
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:251
MoFEM::EntPolynomialBaseCtx::copyNodeBase
const FieldApproximationBase copyNodeBase
Definition: EntPolynomialBaseCtx.hpp:40
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::QuadPolynomialBase::getValueH1
MoFEMErrorCode getValueH1(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:17
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:115
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::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:230
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:38
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
convert.n
n
Definition: convert.py:82
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::QuadPolynomialBase::getValueL2
MoFEMErrorCode getValueL2(MatrixDouble &pts)
Definition: QuadPolynomialBase.cpp:181
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:392
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:56
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
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
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:346
MoFEM::EntPolynomialBaseCtx::basePolynomialsType0
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Definition: EntPolynomialBaseCtx.hpp:27