v0.14.0
Loading...
Searching...
No Matches
FlatPrismPolynomialBase.cpp
Go to the documentation of this file.
1/** \file FlatPrismPolynomialBase.cpp
2\brief Implementation of Ainsworth-Cole H1 base on edge
3*/
4
5
6
7using namespace MoFEM;
8
10 boost::typeindex::type_index type_index, UnknownInterface **iface) const {
11 *iface = const_cast<FlatPrismPolynomialBaseCtx *>(this);
12 return 0;
13}
14
16 EntitiesFieldData &data, moab::Interface &moab,
17 const NumeredEntFiniteElement *fe_ptr, const FieldSpace space,
18 const FieldApproximationBase base,
19 const FieldApproximationBase copy_node_base)
20 : EntPolynomialBaseCtx(data, space, base, copy_node_base), mOab(moab),
21 fePtr(fe_ptr) {
23 CHKERRABORT(PETSC_COMM_WORLD, ierr);
24}
26
28 boost::typeindex::type_index type_index, UnknownInterface **iface) const {
29 *iface = const_cast<FlatPrismPolynomialBase *>(this);
30 return 0;
31}
32
35
38 boost::shared_ptr<BaseFunctionCtx> ctx_ptr) {
40
42
43 if (!cTx->fePtr)
44 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
45 "Pointer to element should be given "
46 "when EntPolynomialBaseCtx is constructed "
47 "(use different constructor)");
48
49 int nb_gauss_pts = pts.size2();
50 if (!nb_gauss_pts)
52
53 if (pts.size1() < 1)
54 SETERRQ(
55 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
56 "Wrong dimension of pts, should be at least 3 rows with coordinates");
57
58 const FieldApproximationBase base = cTx->bAse;
59 EntitiesFieldData &data = cTx->dAta;
60
62 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
63 else
64 data.dataOnEntities[MBVERTEX][0].getNSharedPtr(base) =
65 data.dataOnEntities[MBVERTEX][0].getNSharedPtr(cTx->copyNodeBase);
66
67 data.dataOnEntities[MBVERTEX][0].getN(base).resize(nb_gauss_pts, 6, false);
68 data.dataOnEntities[MBVERTEX][0].getDiffN(base).resize(nb_gauss_pts, 12,
69 false);
70 if (data.dataOnEntities[MBVERTEX][0].getN(base).size1() !=
71 (unsigned int)nb_gauss_pts) {
72 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
73 "Base functions or nodes has wrong number of integration points "
74 "for base %s",
76 }
77
78 data.dataOnEntities[MBVERTEX][0].getN(base).resize(nb_gauss_pts, 6, false);
79 data.dataOnEntities[MBVERTEX][0].getDiffN(base).resize(nb_gauss_pts, 12,
80 false);
81 N.resize(nb_gauss_pts, 3, false);
82 diffN.resize(3, 2, false);
83 CHKERR ShapeMBTRI(&*N.data().begin(), &pts(0, 0), &pts(1, 0), nb_gauss_pts);
85 &*diffN.data().begin());
86
87 // This is needed to have proper order of nodes on faces
88 CHKERR cTx->mOab.get_connectivity(cTx->fePtr->getEnt(), connPrism, numNodes,
89 true);
90 SideNumber_multiIndex &side_table =
92 SideNumber_multiIndex::nth_index<1>::type::iterator siit3 =
93 side_table.get<1>().find(boost::make_tuple(MBTRI, 3));
94 SideNumber_multiIndex::nth_index<1>::type::iterator siit4 =
95 side_table.get<1>().find(boost::make_tuple(MBTRI, 4));
96 if (siit3 == side_table.get<1>().end())
97 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
98 if (siit4 == side_table.get<1>().end())
99 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
100 CHKERR cTx->mOab.get_connectivity(siit3->get()->ent, connFace3, numNodes,
101 true);
102 CHKERR cTx->mOab.get_connectivity(siit4->get()->ent, connFace4, numNodes,
103 true);
104
105 for (int nn = 0; nn < 3; nn++) {
106 faceNodes[0][nn] = std::distance(
107 connPrism, std::find(connPrism, connPrism + 3, connFace3[nn]));
108 faceNodes[1][nn] = std::distance(
109 connPrism + 3, std::find(connPrism + 3, connPrism + 6, connFace4[nn]));
110 for (int gg = 0; gg < nb_gauss_pts; gg++) {
111 double val = N(gg, nn);
112 double val_x = diffN(nn, 0);
113 double val_y = diffN(nn, 1);
114 data.dataOnEntities[MBVERTEX][0].getN(base)(gg, nn) = val;
115 data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 2 * nn + 0) = val_x;
116 data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 2 * nn + 1) = val_y;
117 data.dataOnEntities[MBVERTEX][0].getN(base)(gg, 3 + nn) = val;
118 data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 6 + 2 * nn + 0) =
119 val_x;
120 data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 6 + 2 * nn + 1) =
121 val_y;
122 }
123 }
124
125 switch (cTx->sPace) {
126 case H1:
127 CHKERR getValueH1(pts);
128 break;
129 case HDIV:
130 CHKERR getValueHdiv(pts);
131 break;
132 case HCURL:
134 break;
135 case L2:
136 CHKERR getValueL2(pts);
137 break;
138 default:
139 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
140 }
141
143}
144
147
148 EntitiesFieldData &data = cTx->dAta;
149 const FieldApproximationBase base = cTx->bAse;
150 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s, double *L,
151 double *diffL, const int dim) =
153
154 int nb_gauss_pts = pts.size2();
155
156 // edges
157 if (data.dataOnEntities[MBEDGE].size() != 9)
158 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
159
160 int sense[9], order[9];
161 double *H1edgeN[9], *diffH1edgeN[9];
162
163 auto set_edge_base_data = [&](const int ee) {
165 auto &ent_data = data.dataOnEntities[MBEDGE][ee];
166 if (ent_data.getSense() == 0)
167 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
168 sense[ee] = ent_data.getSense();
169 order[ee] = ent_data.getOrder();
170 int nb_dofs = NBEDGE_H1(ent_data.getOrder());
171 ent_data.getN(base).resize(nb_gauss_pts, nb_dofs, false);
172 ent_data.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
173 H1edgeN[ee] = &*ent_data.getN(base).data().begin();
174 diffH1edgeN[ee] = &*ent_data.getDiffN(base).data().begin();
176 };
177
178 if ((data.spacesOnEntities[MBEDGE]).test(H1)) {
179
180 for (int ee = 0; ee != 3; ++ee)
181 CHKERR set_edge_base_data(ee);
182 for (int ee = 6; ee != 9; ++ee)
183 CHKERR set_edge_base_data(ee);
184
185 // shape functions on face 3
187 &sense[0], &order[0], &*N.data().begin(), &*diffN.data().begin(),
188 &H1edgeN[0], &diffH1edgeN[0], nb_gauss_pts, base_polynomials);
189 // shape functions on face 4
191 &sense[6], &order[6], &*N.data().begin(), &*diffN.data().begin(),
192 &H1edgeN[6], &diffH1edgeN[6], nb_gauss_pts, base_polynomials);
193 }
194
195 // face
196 if (data.dataOnEntities[MBTRI].size() != 5)
197 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
198
199 if ((data.spacesOnEntities[MBTRI]).test(H1)) {
200 for (int ff = 3; ff <= 4; ff++) {
201 int nb_dofs = NBFACETRI_H1(data.dataOnEntities[MBTRI][ff].getOrder());
202 data.dataOnEntities[MBTRI][ff].getN(base).resize(nb_gauss_pts, nb_dofs,
203 false);
204 data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(nb_gauss_pts,
205 2 * nb_dofs, false);
207 faceNodes[ff - 3], data.dataOnEntities[MBTRI][ff].getOrder(),
208 &*N.data().begin(), &*diffN.data().begin(),
209 &*data.dataOnEntities[MBTRI][ff].getN(base).data().begin(),
210 &*data.dataOnEntities[MBTRI][ff].getDiffN(base).data().begin(),
211 nb_gauss_pts, base_polynomials);
212 }
213 }
214
216}
217
220 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
222}
223
225 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
226}
227
229 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
230}
static Index< 'p', 3 > p
FieldApproximationBase
approximation base
Definition: definitions.h:58
@ LASTBASE
Definition: definitions.h:69
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
FieldSpace
approximation spaces
Definition: definitions.h:82
@ L2
field with C-1 continuity
Definition: definitions.h:88
@ H1
continuous field
Definition: definitions.h:85
@ HCURL
field with continuous tangents
Definition: definitions.h:86
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
static const char *const ApproximationBaseNames[]
Definition: definitions.h:72
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
const int dim
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182
multi_index_container< boost::shared_ptr< SideNumber >, indexed_by< ordered_unique< member< SideNumber, EntityHandle, &SideNumber::ent > >, ordered_non_unique< composite_key< SideNumber, const_mem_fun< SideNumber, EntityType, &SideNumber::getEntType >, member< SideNumber, signed char, &SideNumber::side_number > > > > > SideNumber_multiIndex
SideNumber_multiIndex for SideNumber.
#define NBEDGE_H1(P)
Numer of base function on edge for H1 space.
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.
Definition: h1.c:17
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))
Definition: h1.c:191
#define NBFACETRI_H1(P)
Number of base function on triangle for H1 space.
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
Class used to pass element data to calculate base functions on tet,triangle,edge.
const FieldApproximationBase copyNodeBase
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
const FieldApproximationBase bAse
data structure for finite element entity
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Class used to pass element data to calculate base functions on flat prism.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
const NumeredEntFiniteElement * fePtr
FlatPrismPolynomialBaseCtx(EntitiesFieldData &data, moab::Interface &moab, const NumeredEntFiniteElement *fe_ptr, const FieldSpace space, const FieldApproximationBase base, const FieldApproximationBase copy_node_base=LASTBASE)
Calculate base functions on tetrahedralFIXME: Need moab and mofem finite element structure to work (t...
MoFEMErrorCode getValueHcurl(MatrixDouble &pts)
FlatPrismPolynomialBaseCtx * cTx
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
MoFEMErrorCode getValueHdiv(MatrixDouble &pts)
MoFEMErrorCode getValueH1(MatrixDouble &pts)
MoFEMErrorCode getValueL2(MatrixDouble &pts)
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Partitioned (Indexed) Finite Element in Problem.
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:104
base class for all interface classes
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
SideNumber_multiIndex & getSideNumberTable() const