|
| v0.14.0
|
Calculate base functions on triangle.
More...
#include <src/approximation/QuadPolynomialBase.hpp>
Calculate base functions on triangle.
Definition at line 21 of file QuadPolynomialBase.hpp.
◆ QuadPolynomialBase()
MoFEM::QuadPolynomialBase::QuadPolynomialBase |
( |
| ) |
|
|
default |
◆ ~QuadPolynomialBase()
MoFEM::QuadPolynomialBase::~QuadPolynomialBase |
( |
| ) |
|
|
default |
◆ getValue()
Reimplemented from MoFEM::BaseFunction.
Definition at line 434 of file QuadPolynomialBase.cpp.
440 int nb_gauss_pts = pts.size2();
447 "Wrong dimension of pts, should be at least 3 rows with coordinates");
459 if (base_shape.size1() != pts.size2())
461 "Number of base functions integration points not equal number of "
462 "set integration point");
463 if (base_shape.size2() != 4)
465 "Number of shape functions should be four");
466 if (diff_base.size1() != pts.size2())
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)
474 "Number of shape functions should be four");
◆ getValueH1()
◆ getValueH1AinsworthBase()
Definition at line 43 of file QuadPolynomialBase.cpp.
52 "Polynomial type not set");
53 PetscErrorCode (*base_polynomials)(
int p,
double s,
double *diff_s,
double *
L,
54 double *diffL,
const int dim) =
57 auto ©_base_fun = data.
dataOnEntities[MBVERTEX][0].getN(copy_base);
58 auto ©_diff_base_fun =
61 int nb_gauss_pts = pts.size2();
68 "should be four edges on quad");
70 int sense[4],
order[4];
71 double *H1edgeN[4], *diffH1edgeN[4];
72 for (
int ee = 0; ee != 4; ++ee) {
74 if (ent_dat.getSense() == 0)
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();
86 sense,
order, &*copy_base_fun.data().begin(),
87 &*copy_diff_base_fun.data().begin(), H1edgeN, diffH1edgeN, nb_gauss_pts,
96 "should be one quad to store bubble base on quad");
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,
◆ getValueH1DemkowiczBase()
Definition at line 115 of file QuadPolynomialBase.cpp.
122 int nb_gauss_pts = pts.size2();
123 auto ©_base_fun = data.
dataOnEntities[MBVERTEX][0].getN(copy_base);
124 auto ©_diff_base_fun =
131 "should be four edges on quad");
133 int sense[4],
order[4];
134 double *H1edgeN[4], *diffH1edgeN[4];
135 for (
int ee = 0; ee != 4; ++ee) {
137 if (ent_dat.getSense() == 0)
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();
149 sense,
order, &*copy_base_fun.data().begin(),
150 &*copy_diff_base_fun.data().begin(), H1edgeN, diffH1edgeN,
159 "should be one quad to store bubble base on quad");
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);
168 int face_nodes[] = {0, 1, 2, 3};
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);
◆ getValueHcurl()
◆ getValueHcurlDemkowiczBase()
Definition at line 251 of file QuadPolynomialBase.cpp.
258 int nb_gauss_pts = pts.size2();
259 auto ©_base_fun = data.
dataOnEntities[MBVERTEX][0].getN(copy_base);
260 auto ©_diff_base_fun =
268 "wrong number of edges on quad, should be 4 but is %d",
271 int sense[4],
order[4];
272 double *hcurl_edge_n[4];
273 double *diff_hcurl_edge_n[4];
275 for (
int ee = 0; ee != 4; ++ee) {
279 "orientation (sense) of edge is not set");
289 nb_gauss_pts, 3 * 2 * nb_dofs,
false);
293 diff_hcurl_edge_n[ee] =
297 sense,
order, &*copy_base_fun.data().begin(),
298 &*copy_diff_base_fun.data().begin(), hcurl_edge_n, diff_hcurl_edge_n,
307 "No data struture to keep base functions on face");
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);
315 int order[2] = {p, p};
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);
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);
337 double *ptr = &face_n(0, 0);
339 for (
int j = 0;
j != 3; ++
j) {
344 for (
int j = 0;
j != 3; ++
j) {
353 double *diff_ptr = &diff_face_n(0, 0);
355 for (
int j = 0;
j != 6; ++
j) {
356 *diff_ptr = *diff_ptr_f0;
360 for (
int j = 0;
j != 6; ++
j) {
361 *diff_ptr = *diff_ptr_f1;
◆ getValueHdiv()
◆ getValueHdivDemkowiczBase()
Definition at line 392 of file QuadPolynomialBase.cpp.
398 int nb_gauss_pts = pts.size2();
400 auto ©_base_fun = data.
dataOnEntities[MBVERTEX][0].getN(copy_base);
401 auto ©_diff_base_fun =
409 "No data struture to keep base functions on face");
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);
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);
◆ getValueL2()
◆ getValueL2DemkowiczBase()
Definition at line 203 of file QuadPolynomialBase.cpp.
210 int nb_gauss_pts = pts.size2();
211 auto ©_base_fun = data.
dataOnEntities[MBVERTEX][0].getN(copy_base);
212 auto ©_diff_base_fun =
216 int p = ent_dat.getOrder();
217 int order[2] = {p, 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);
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);
◆ query_interface()
◆ cTx
◆ diffFaceFamily
◆ faceFamily
The documentation for this struct was generated from the following files:
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
#define NBEDGE_H1(P)
Numer of base function on edge for H1 space.
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))
Class used to pass element data to calculate base functions on tet,triangle,edge.
MoFEMErrorCode getValueL2DemkowiczBase(MatrixDouble &pts)
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.
const static char *const ApproximationBaseNames[]
@ L2
field with C-1 continuity
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))
MatrixDouble diffFaceFamily
MoFEMErrorCode getValueHdiv(MatrixDouble &pts)
#define NBFACEQUAD_DEMKOWICZ_HCURL(P)
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
MoFEMErrorCode getValueH1AinsworthBase(MatrixDouble &pts)
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.
MoFEMErrorCode getValueHcurlDemkowiczBase(MatrixDouble &pts)
const FieldApproximationBase copyNodeBase
#define CHKERR
Inline error check.
MoFEMErrorCode getValueH1(MatrixDouble &pts)
#define NBFACEQUAD_L2(P)
Number of base functions on quad for L2 space.
MoFEMErrorCode getValueH1DemkowiczBase(MatrixDouble &pts)
MoFEMErrorCode Hdiv_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N, double *diffN, double *faceN, double *diff_faceN, int nb_integration_pts)
#define NBEDGE_DEMKOWICZ_HCURL(P)
#define NBFACEQUAD_DEMKOWICZ_FAMILY_HCURL(P, Q)
Number of base functions on quad for Hcurl space.
MoFEMErrorCode getValueHcurl(MatrixDouble &pts)
#define NBFACEQUAD_H1(P)
Number of base functions on quad for H1 space.
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.
Calculate base functions on triangle.
const FieldApproximationBase bAse
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
MoFEMErrorCode getValueL2(MatrixDouble &pts)
#define NBFACEQUAD_DEMKOWICZ_HDIV(P)
MoFEMErrorCode Hcurl_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N, double *diffN, double *faceN[], double *diff_faceN[], int nb_integration_pts)
FTensor::Index< 'j', 3 > j
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
@ HCURL
field with continuous tangents
@ MOFEM_DATA_INCONSISTENCY
MoFEMErrorCode getValueHdivDemkowiczBase(MatrixDouble &pts)
FieldApproximationBase
approximation base
EntPolynomialBaseCtx * cTx
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
MoFEMErrorCode Hcurl_EdgeShapeFunctions_ONQUAD(int *sense, int *p, double *N, double *diffN, double *edgeN[4], double *curl_edgeN[4], int nb_integration_pts)
data structure for finite element entity
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HDIV
field with continuous normal traction
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)