|
| v0.14.0
|
Calculate base functions on triangle.
More...
#include <src/approximation/TriPolynomialBase.hpp>
Calculate base functions on triangle.
Definition at line 16 of file TriPolynomialBase.hpp.
◆ TriPolynomialBase()
MoFEM::TriPolynomialBase::TriPolynomialBase |
( |
| ) |
|
|
default |
◆ ~TriPolynomialBase()
virtual MoFEM::TriPolynomialBase::~TriPolynomialBase |
( |
| ) |
|
|
virtualdefault |
◆ getValue()
Reimplemented from MoFEM::BaseFunction.
Definition at line 844 of file TriPolynomialBase.cpp.
850 int nb_gauss_pts = pts.size2();
858 "Wrong dimension of pts, should be at least 3 rows with coordinates");
865 data.
dataOnEntities[MBVERTEX][0].getN(base).resize(nb_gauss_pts, 3,
869 &pts(0, 0), &pts(1, 0), nb_gauss_pts);
870 data.
dataOnEntities[MBVERTEX][0].getDiffN(base).resize(3, 2,
false);
872 &*data.
dataOnEntities[MBVERTEX][0].getDiffN(base).data().begin());
880 (
unsigned int)nb_gauss_pts) {
882 "Base functions or nodes has wrong number of integration points "
887 auto set_node_derivative_for_all_gauss_pts = [&]() {
892 data.
dataOnEntities[MBVERTEX][0].getDiffN(base).resize(nb_gauss_pts, 6,
894 for (
int gg = 0; gg != nb_gauss_pts; ++gg)
901 CHKERR set_node_derivative_for_all_gauss_pts();
◆ getValueH1()
◆ getValueH1AinsworthBase()
Definition at line 248 of file TriPolynomialBase.cpp.
256 "Polynomial type not set");
258 PetscErrorCode (*base_polynomials)(
int p,
double s,
double *diff_s,
double *
L,
259 double *diffL,
const int dim) =
262 int nb_gauss_pts = pts.size2();
268 "expected size of data.dataOnEntities[MBEDGE] is 3");
270 int sense[3],
order[3];
271 double *H1edgeN[3], *diffH1edgeN[3];
272 for (
int ee = 0; ee < 3; ee++) {
276 "sense of the edge unknown");
281 data.
dataOnEntities[MBEDGE][ee].getN(base).resize(nb_gauss_pts, nb_dofs,
283 data.
dataOnEntities[MBEDGE][ee].getDiffN(base).resize(nb_gauss_pts,
285 H1edgeN[ee] = &*data.
dataOnEntities[MBEDGE][ee].getN(base).data().begin();
293 H1edgeN, diffH1edgeN, nb_gauss_pts, base_polynomials);
300 "expected that size data.dataOnEntities[MBTRI] is one");
304 data.
dataOnEntities[MBTRI][0].getN(base).resize(nb_gauss_pts, nb_dofs,
308 const int face_nodes[] = {0, 1, 2};
315 nb_gauss_pts, base_polynomials);
◆ getValueH1BernsteinBezierBase()
Definition at line 35 of file TriPolynomialBase.cpp.
39 int nb_gauss_pts = pts.size2();
42 auto &ptr = data.getBBAlphaIndicesSharedPtr(
field_name);
56 auto &ptr = data.getBBDiffNSharedPtr(
field_name);
62 auto get_alpha_by_name_ptr =
64 const std::string &
field_name) -> boost::shared_ptr<MatrixInt> & {
65 return data.getBBAlphaIndicesSharedPtr(
field_name);
68 auto get_base_by_name_ptr =
70 const std::string &
field_name) -> boost::shared_ptr<MatrixDouble> & {
74 auto get_diff_base_by_name_ptr =
76 const std::string &
field_name) -> boost::shared_ptr<MatrixDouble> & {
80 auto get_alpha_by_order_ptr =
81 [](
auto &data,
const size_t o) -> boost::shared_ptr<MatrixInt> & {
82 return data.getBBAlphaIndicesByOrderSharedPtr(o);
85 auto get_base_by_order_ptr =
86 [](
auto &data,
const size_t o) -> boost::shared_ptr<MatrixDouble> & {
87 return data.getBBNByOrderSharedPtr(o);
90 auto get_diff_base_by_order_ptr =
91 [](
auto &data,
const size_t o) -> boost::shared_ptr<MatrixDouble> & {
92 return data.getBBDiffNByOrderSharedPtr(o);
96 auto &vert_get_diff_n = get_diff_base(data.
dataOnEntities[MBVERTEX][0]);
97 vert_get_n.resize(nb_gauss_pts, 3,
false);
98 vert_get_diff_n.resize(nb_gauss_pts, 6,
false);
101 (
unsigned int)nb_gauss_pts)
103 "Base functions or nodes has wrong number of integration points "
109 vertex_alpha.resize(3, 3,
false);
110 vertex_alpha.clear();
111 for (
int n = 0;
n != 3; ++
n)
115 1,
lambda.size1(), vertex_alpha.size1(), &vertex_alpha(0, 0),
117 &vert_get_diff_n(0, 0));
118 for (
int n = 0;
n != 3; ++
n) {
119 const int f = boost::math::factorial<double>(
121 for (
int g = 0;
g != nb_gauss_pts; ++
g) {
122 vert_get_n(
g,
n) *=
f;
123 for (
int d = 0;
d != 2; ++
d)
124 vert_get_diff_n(
g, 2 *
n +
d) *=
f;
132 "Wrong size of ent data");
134 constexpr
int edges_nodes[3][2] = {{0, 1}, {1, 2}, {2, 0}};
135 for (
int ee = 0; ee != 3; ++ee) {
138 if (ent_data.getSense() == 0)
140 "Sense of the edge unknown");
142 const int sense = ent_data.getSense();
143 const int order = ent_data.getOrder();
144 const int nb_dofs =
NBEDGE_H1(ent_data.getOrder());
147 if (get_alpha_by_order_ptr(ent_data,
order)) {
149 get_alpha_by_order_ptr(ent_data,
order);
151 get_base_by_order_ptr(ent_data,
order);
152 get_diff_base_by_name_ptr(ent_data,
field_name) =
153 get_diff_base_by_order_ptr(ent_data,
order);
155 auto &get_n = get_base(ent_data);
156 auto &get_diff_n = get_diff_base(ent_data);
157 get_n.resize(nb_gauss_pts, nb_dofs,
false);
158 get_diff_n.resize(nb_gauss_pts, 2 * nb_dofs,
false);
161 edge_alpha.resize(nb_dofs, 3,
false);
165 for (
int i = 0;
i != edge_alpha.size1(); ++
i) {
166 int a = edge_alpha(
i, edges_nodes[ee][0]);
167 edge_alpha(
i, edges_nodes[ee][0]) =
168 edge_alpha(
i, edges_nodes[ee][1]);
169 edge_alpha(
i, edges_nodes[ee][1]) =
a;
173 order,
lambda.size1(), edge_alpha.size1(), &edge_alpha(0, 0),
177 get_alpha_by_order_ptr(ent_data,
order) =
179 get_base_by_order_ptr(ent_data,
order) =
181 get_diff_base_by_order_ptr(ent_data,
order) =
182 get_diff_base_by_name_ptr(ent_data,
field_name);
187 for (
int ee = 0; ee != 3; ++ee) {
189 ent_data.getBBAlphaIndicesSharedPtr(
field_name).reset();
190 auto &get_n = get_base(ent_data);
191 auto &get_diff_n = get_diff_base(ent_data);
192 get_n.resize(nb_gauss_pts, 0,
false);
193 get_diff_n.resize(nb_gauss_pts, 0,
false);
200 "Wrong size ent of ent data");
203 const int order = ent_data.getOrder();
205 if (get_alpha_by_order_ptr(ent_data,
order)) {
207 get_alpha_by_order_ptr(ent_data,
order);
209 get_base_by_order_ptr(ent_data,
order);
210 get_diff_base_by_name_ptr(ent_data,
field_name) =
211 get_diff_base_by_order_ptr(ent_data,
order);
214 auto &get_n = get_base(ent_data);
215 auto &get_diff_n = get_diff_base(ent_data);
216 get_n.resize(nb_gauss_pts, nb_dofs,
false);
217 get_diff_n.resize(nb_gauss_pts, 2 * nb_dofs,
false);
219 auto &tri_alpha = get_alpha(ent_data);
220 tri_alpha.resize(nb_dofs, 3,
false);
224 order,
lambda.size1(), tri_alpha.size1(), &tri_alpha(0, 0),
228 get_alpha_by_order_ptr(ent_data,
order) =
230 get_base_by_order_ptr(ent_data,
order) =
232 get_diff_base_by_order_ptr(ent_data,
order) =
233 get_diff_base_by_name_ptr(ent_data,
field_name);
238 ent_data.getBBAlphaIndicesSharedPtr(
field_name).reset();
239 auto &get_n = get_base(ent_data);
240 auto &get_diff_n = get_diff_base(ent_data);
241 get_n.resize(nb_gauss_pts, 0,
false);
242 get_diff_n.resize(nb_gauss_pts, 0,
false);
◆ getValueHcurl()
◆ getValueHcurlAinsworthBase()
Definition at line 645 of file TriPolynomialBase.cpp.
653 PetscErrorCode (*base_polynomials)(
int p,
double s,
double *diff_s,
double *
L,
654 double *diffL,
const int dim) =
657 int nb_gauss_pts = pts.size2();
664 int sense[3],
order[3];
665 double *hcurl_edge_n[3];
666 double *diff_hcurl_edge_n[3];
667 for (
int ee = 0; ee < 3; ee++) {
670 "data inconsistency");
679 nb_gauss_pts, 2 * 3 * nb_dofs,
false);
682 diff_hcurl_edge_n[ee] =
689 hcurl_edge_n, diff_hcurl_edge_n, nb_gauss_pts, base_polynomials);
691 for (
int ee = 0; ee < 3; ee++) {
692 data.
dataOnEntities[MBEDGE][ee].getN(base).resize(nb_gauss_pts, 0,
false);
693 data.
dataOnEntities[MBEDGE][ee].getDiffN(base).resize(nb_gauss_pts, 0,
709 data.
dataOnEntities[MBTRI][0].getN(base).resize(nb_gauss_pts, 3 * nb_dofs,
712 3 * 2 * nb_dofs,
false);
714 int face_nodes[] = {0, 1, 2};
721 nb_gauss_pts, base_polynomials);
724 data.
dataOnEntities[MBTRI][0].getN(base).resize(nb_gauss_pts, 0,
false);
725 data.
dataOnEntities[MBTRI][0].getDiffN(base).resize(nb_gauss_pts, 0,
false);
◆ getValueHcurlDemkowiczBase()
Definition at line 732 of file TriPolynomialBase.cpp.
738 int nb_gauss_pts = pts.size2();
746 "wrong number of data structures on edges, should be three but is %d",
749 int sense[3],
order[3];
750 double *hcurl_edge_n[3];
751 double *diff_hcurl_edge_n[3];
753 for (
int ee = 0; ee != 3; ++ee) {
757 "orientation (sense) of edge is not set");
766 nb_gauss_pts, 2 * 3 * nb_dofs,
false);
770 diff_hcurl_edge_n[ee] =
778 hcurl_edge_n, diff_hcurl_edge_n, nb_gauss_pts);
784 for (
int ee = 0; ee != 3; ++ee) {
785 data.
dataOnEntities[MBEDGE][ee].getN(base).resize(nb_gauss_pts, 0,
false);
786 data.
dataOnEntities[MBEDGE][ee].getDiffN(base).resize(nb_gauss_pts, 0,
796 "No data struture to keep base functions on face");
800 data.
dataOnEntities[MBTRI][0].getN(base).resize(nb_gauss_pts, 3 * nb_dofs,
803 3 * 2 * nb_dofs,
false);
805 int face_nodes[] = {0, 1, 2};
818 data.
dataOnEntities[MBTRI][0].getN(base).resize(nb_gauss_pts, 0,
false);
819 data.
dataOnEntities[MBTRI][0].getDiffN(base).resize(nb_gauss_pts, 0,
false);
◆ getValueHdiv()
◆ getValueHdivAinsworthBase()
Definition at line 525 of file TriPolynomialBase.cpp.
532 "Polynomial type not set");
533 PetscErrorCode (*base_polynomials)(
int p,
double s,
double *diff_s,
double *
L,
534 double *diffL,
const int dim) =
537 int nb_gauss_pts = pts.size2();
545 if (face_order > 0) {
547 for (
int ee = 0; ee < 3; ee++) {
556 int face_nodes[3] = {0, 1, 2};
558 face_nodes, face_order,
559 &data.
dataOnEntities[MBVERTEX][0].getN(base)(0, 0), NULL, PHI_f_e, NULL,
560 nb_gauss_pts, 3, base_polynomials);
562 face_nodes, face_order,
563 &data.
dataOnEntities[MBVERTEX][0].getN(base)(0, 0), NULL, PHI_f, NULL,
564 nb_gauss_pts, 3, base_polynomials);
573 for (
int oo = 0; oo < face_order; oo++) {
574 for (
int ee = 0; ee < 3; ee++) {
577 for (
int gg = 0; gg < nb_gauss_pts; gg++) {
585 for (
int gg = 0; gg < nb_gauss_pts; gg++) {
◆ getValueHdivDemkowiczBase()
Definition at line 596 of file TriPolynomialBase.cpp.
607 "This should be used only with DEMKOWICZ_JACOBI_BASE "
612 int nb_gauss_pts = pts.size2();
615 double *phi_f = &*data.
dataOnEntities[MBTRI][0].getN(base).data().begin();
618 int face_nodes[3] = {0, 1, 2};
622 &*data.
dataOnEntities[MBVERTEX][0].getDiffN(base).data().begin(), phi_f,
623 NULL, nb_gauss_pts, 3);
◆ getValueL2()
◆ getValueL2AinsworthBase()
Definition at line 339 of file TriPolynomialBase.cpp.
346 "Polynomial type not set");
347 PetscErrorCode (*base_polynomials)(
int p,
double s,
double *diff_s,
double *
L,
348 double *diffL,
const int dim) =
351 int nb_gauss_pts = pts.size2();
366 nb_gauss_pts, base_polynomials);
◆ getValueL2BernsteinBezierBase()
Definition at line 372 of file TriPolynomialBase.cpp.
377 int nb_gauss_pts = pts.size2();
380 auto &ptr = data.getBBAlphaIndicesSharedPtr(
field_name);
394 auto &ptr = data.getBBDiffNSharedPtr(
field_name);
400 auto get_alpha_by_name_ptr =
402 const std::string &
field_name) -> boost::shared_ptr<MatrixInt> & {
403 return data.getBBAlphaIndicesSharedPtr(
field_name);
406 auto get_base_by_name_ptr =
408 const std::string &
field_name) -> boost::shared_ptr<MatrixDouble> & {
412 auto get_diff_base_by_name_ptr =
414 const std::string &
field_name) -> boost::shared_ptr<MatrixDouble> & {
418 auto get_alpha_by_order_ptr =
419 [](
auto &data,
const size_t o) -> boost::shared_ptr<MatrixInt> & {
420 return data.getBBAlphaIndicesByOrderSharedPtr(o);
423 auto get_base_by_order_ptr =
424 [](
auto &data,
const size_t o) -> boost::shared_ptr<MatrixDouble> & {
425 return data.getBBNByOrderSharedPtr(o);
428 auto get_diff_base_by_order_ptr =
429 [](
auto &data,
const size_t o) -> boost::shared_ptr<MatrixDouble> & {
430 return data.getBBDiffNByOrderSharedPtr(o);
436 "Wrong size ent of ent data");
439 (
unsigned int)nb_gauss_pts)
441 "Base functions or nodes has wrong number of integration points "
447 const int order = ent_data.getOrder();
450 if (get_alpha_by_order_ptr(ent_data,
order)) {
452 get_alpha_by_order_ptr(ent_data,
order);
454 get_base_by_order_ptr(ent_data,
order);
455 get_diff_base_by_name_ptr(ent_data,
field_name) =
456 get_diff_base_by_order_ptr(ent_data,
order);
459 auto &get_n = get_base(ent_data);
460 auto &get_diff_n = get_diff_base(ent_data);
461 get_n.resize(nb_gauss_pts, nb_dofs,
false);
462 get_diff_n.resize(nb_gauss_pts, 2 * nb_dofs,
false);
469 "Inconsistent number of DOFs");
471 auto &tri_alpha = get_alpha(ent_data);
480 "Inconsistent number of DOFs");
482 auto &tri_alpha = get_alpha(ent_data);
483 tri_alpha.resize(nb_dofs, 3,
false);
490 std::array<int *, 3> face_edge_ptr{
494 face_n.data(), face_edge_ptr.data());
500 order,
lambda.size1(), tri_alpha.size1(), &tri_alpha(0, 0),
505 get_alpha_by_order_ptr(ent_data,
order) =
507 get_base_by_order_ptr(ent_data,
order) =
509 get_diff_base_by_order_ptr(ent_data,
order) =
510 get_diff_base_by_name_ptr(ent_data,
field_name);
515 ent_data.getBBAlphaIndicesSharedPtr(
field_name).reset();
516 auto &get_n = get_base(ent_data);
517 auto &get_diff_n = get_diff_base(ent_data);
518 get_n.resize(nb_gauss_pts, 0,
false);
519 get_diff_n.resize(nb_gauss_pts, 0,
false);
◆ query_interface()
◆ cTx
◆ diffN_face_bubble
ublas::vector<MatrixDouble> MoFEM::TriPolynomialBase::diffN_face_bubble |
|
private |
◆ diffN_face_edge
ublas::matrix<MatrixDouble> MoFEM::TriPolynomialBase::diffN_face_edge |
|
private |
◆ N_face_bubble
ublas::vector<MatrixDouble> MoFEM::TriPolynomialBase::N_face_bubble |
|
private |
◆ N_face_edge
ublas::matrix<MatrixDouble> MoFEM::TriPolynomialBase::N_face_edge |
|
private |
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.
EntPolynomialBaseCtx * cTx
#define NBEDGE_H1(P)
Numer of base function on edge for H1 space.
MoFEMErrorCode Hdiv_Demkowicz_Face_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int gdim, int nb)
MoFEMErrorCode getValueHdivDemkowiczBase(MatrixDouble &pts)
PetscErrorCode L2_Ainsworth_ShapeFunctions_MBTRI(int p, double *N, double *diffN, double *L2N, double *diff_L2N, int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Get base functions on triangle for L2 space.
Class used to pass element data to calculate base functions on tet,triangle,edge.
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.
ublas::matrix< MatrixDouble > N_face_edge
const static char *const ApproximationBaseNames[]
@ L2
field with C-1 continuity
#define NBFACETRI_AINSWORTH_EDGE_HDIV(P)
MoFEMErrorCode Hdiv_Ainsworth_EdgeFaceShapeFunctions_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f_e[3], double *diff_phi_f_e[3], int gdim, int nb, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Hdiv base functions, Edge-based face functions by Ainsworth .
UBlasMatrix< double > MatrixDouble
MoFEMErrorCode getValueHdiv(MatrixDouble &pts)
static MoFEMErrorCode generateIndicesTriTri(const int N, int *alpha)
MoFEMErrorCode Hcurl_Ainsworth_FaceFunctions_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face H-curl functions.
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
MoFEMErrorCode Hcurl_Demkowicz_EdgeBaseFunctions_MBTRI(int *sense, int *p, double *n, double *diff_n, double *phi[], double *diff_phi[], int nb_integration_pts)
Edge based H-curl base functions on teriangle.
const FieldApproximationBase copyNodeBase
MoFEMErrorCode getValueHdivAinsworthBase(MatrixDouble &pts)
#define CHKERR
Inline error check.
MoFEMErrorCode Hcurl_Ainsworth_EdgeBaseFunctions_MBTET_ON_FACE(int *sense, int *p, double *N, double *diffN, double *edgeN[], double *diff_edgeN[], int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Edge based H-curl base functions on face.
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))
MoFEMErrorCode getValueH1(MatrixDouble &pts)
MoFEMErrorCode getValueHcurlDemkowiczBase(MatrixDouble &pts)
Calculate base functions on triangle.
#define NBFACETRI_DEMKOWICZ_HDIV(P)
UBlasMatrix< int > MatrixInt
MoFEMErrorCode getValueL2BernsteinBezierBase(MatrixDouble &pts)
static MoFEMErrorCode generateIndicesEdgeTri(const int N[], int *alpha[])
#define NBEDGE_DEMKOWICZ_HCURL(P)
MoFEMErrorCode getValueL2AinsworthBase(MatrixDouble &pts)
MoFEMErrorCode Hdiv_Ainsworth_FaceBubbleShapeFunctions_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int gdim, int nb, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face bubble functions by Ainsworth .
MoFEMErrorCode getValueHcurl(MatrixDouble &pts)
FTensor::Index< 'i', SPACE_DIM > i
const FieldApproximationBase bAse
MoFEMErrorCode getValueH1AinsworthBase(MatrixDouble &pts)
constexpr auto field_name
#define NBFACETRI_DEMKOWICZ_HCURL(P)
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
#define NBEDGE_AINSWORTH_HCURL(P)
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
#define NBFACETRI_H1(P)
Number of base function on triangle for H1 space.
MoFEMErrorCode getValueH1BernsteinBezierBase(MatrixDouble &pts)
#define NBFACETRI_L2(P)
Number of base functions on triangle for L2 space.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
@ HCURL
field with continuous tangents
@ MOFEM_DATA_INCONSISTENCY
static MoFEMErrorCode baseFunctionsTri(const int N, const int gdim, const int n_alpha, const int *alpha, const double *lambda, const double *grad_lambda, double *base, double *grad_base)
const std::string fieldName
FieldApproximationBase
approximation base
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
#define NBFACETRI_AINSWORTH_FACE_HDIV(P)
ublas::vector< MatrixDouble > N_face_bubble
#define NBFACETRI_AINSWORTH_HCURL(P)
MoFEMErrorCode Hcurl_Demkowicz_FaceBaseFunctions_MBTRI(int *faces_nodes, int p, double *n, double *diff_n, double *phi, double *diff_phi, int nb_integration_pts)
Face base interior function.
data structure for finite element entity
static MoFEMErrorCode generateIndicesVertexTri(const int N, int *alpha)
MoFEMErrorCode getValueL2(MatrixDouble &pts)
#define NBFACETRI_AINSWORTH_HDIV(P)
MoFEMErrorCode getValueHcurlAinsworthBase(MatrixDouble &pts)
#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)