|
| v0.14.0
|
Calculate base functions on triangle.
More...
#include <src/approximation/QuadPolynomialBase.hpp>
|
using | DofsSideMap = multi_index_container< DofsSideMapData, indexed_by< ordered_non_unique< tag< TypeSide_mi_tag >, composite_key< DofsSideMapData, member< DofsSideMapData, EntityType, &DofsSideMapData::type >, member< DofsSideMapData, int, &DofsSideMapData::side > >>, ordered_unique< tag< EntDofIdx_mi_tag >, member< DofsSideMapData, int, &DofsSideMapData::dof > > > > |
| Map entity stype and side to element/entity dof index. More...
|
|
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...
|
|
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 430 of file QuadPolynomialBase.cpp.
436 int nb_gauss_pts = pts.size2();
443 "Wrong dimension of pts, should be at least 3 rows with coordinates");
455 if (base_shape.size1() != pts.size2())
457 "Number of base functions integration points not equal number of "
458 "set integration point");
459 if (base_shape.size2() != 4)
461 "Number of shape functions should be four");
462 if (diff_base.size1() != pts.size2())
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)
470 "Number of shape functions should be four");
◆ getValueH1()
◆ getValueH1AinsworthBase()
Definition at line 41 of file QuadPolynomialBase.cpp.
50 "Polynomial type not set");
51 PetscErrorCode (*base_polynomials)(
int p,
double s,
double *diff_s,
double *
L,
52 double *diffL,
const int dim) =
55 auto ©_base_fun = data.
dataOnEntities[MBVERTEX][0].getN(copy_base);
56 auto ©_diff_base_fun =
59 int nb_gauss_pts = pts.size2();
66 "should be four edges on quad");
68 int sense[4],
order[4];
69 double *H1edgeN[4], *diffH1edgeN[4];
70 for (
int ee = 0; ee != 4; ++ee) {
72 if (ent_dat.getSense() == 0)
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();
84 sense,
order, &*copy_base_fun.data().begin(),
85 &*copy_diff_base_fun.data().begin(), H1edgeN, diffH1edgeN, nb_gauss_pts,
94 "should be one quad to store bubble base on quad");
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,
◆ getValueH1DemkowiczBase()
Definition at line 112 of file QuadPolynomialBase.cpp.
119 int nb_gauss_pts = pts.size2();
120 auto ©_base_fun = data.
dataOnEntities[MBVERTEX][0].getN(copy_base);
121 auto ©_diff_base_fun =
128 "should be four edges on quad");
130 int sense[4],
order[4];
131 double *H1edgeN[4], *diffH1edgeN[4];
132 for (
int ee = 0; ee != 4; ++ee) {
134 if (ent_dat.getSense() == 0)
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();
146 sense,
order, &*copy_base_fun.data().begin(),
147 &*copy_diff_base_fun.data().begin(), H1edgeN, diffH1edgeN,
156 "should be one quad to store bubble base on quad");
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);
165 int face_nodes[] = {0, 1, 2, 3};
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);
◆ getValueHcurl()
◆ getValueHcurlDemkowiczBase()
Definition at line 248 of file QuadPolynomialBase.cpp.
255 int nb_gauss_pts = pts.size2();
256 auto ©_base_fun = data.
dataOnEntities[MBVERTEX][0].getN(copy_base);
257 auto ©_diff_base_fun =
265 "wrong number of edges on quad, should be 4 but is %d",
268 int sense[4],
order[4];
269 double *hcurl_edge_n[4];
270 double *diff_hcurl_edge_n[4];
272 for (
int ee = 0; ee != 4; ++ee) {
276 "orientation (sense) of edge is not set");
286 nb_gauss_pts, 3 * 2 * nb_dofs,
false);
290 diff_hcurl_edge_n[ee] =
294 sense,
order, &*copy_base_fun.data().begin(),
295 &*copy_diff_base_fun.data().begin(), hcurl_edge_n, diff_hcurl_edge_n,
304 "No data struture to keep base functions on face");
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);
312 int order[2] = {p, p};
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);
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);
334 double *ptr = &face_n(0, 0);
336 for (
int j = 0;
j != 3; ++
j) {
341 for (
int j = 0;
j != 3; ++
j) {
350 double *diff_ptr = &diff_face_n(0, 0);
352 for (
int j = 0;
j != 6; ++
j) {
353 *diff_ptr = *diff_ptr_f0;
357 for (
int j = 0;
j != 6; ++
j) {
358 *diff_ptr = *diff_ptr_f1;
◆ getValueHdiv()
◆ getValueHdivDemkowiczBase()
Definition at line 389 of file QuadPolynomialBase.cpp.
395 int nb_gauss_pts = pts.size2();
397 auto ©_base_fun = data.
dataOnEntities[MBVERTEX][0].getN(copy_base);
398 auto ©_diff_base_fun =
406 "No data struture to keep base functions on face");
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);
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);
◆ getValueL2()
◆ getValueL2DemkowiczBase()
Definition at line 200 of file QuadPolynomialBase.cpp.
207 int nb_gauss_pts = pts.size2();
208 auto ©_base_fun = data.
dataOnEntities[MBVERTEX][0].getN(copy_base);
209 auto ©_diff_base_fun =
213 int p = ent_dat.getOrder();
214 int order[2] = {p, 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);
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);
◆ 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 reference to pointer of interface.
#define NBEDGE_H1(P)
Number 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
const FieldContinuity spaceContinuity
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
@ CONTINUOUS
Regular field.
#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)