v0.14.0
Loading...
Searching...
No Matches
FlatPrismElementForcesAndSourcesCore.cpp
Go to the documentation of this file.
1/** \file FlatPrismElementForcesAndSourcesCore.cpp
2
3\brief Implementation of flat prism element
4
5*/
6
7
8
9namespace MoFEM {
10
12 Interface &m_field)
13 : ForcesAndSourcesCore(m_field),
14 meshPositionsFieldName("MESH_NODE_POSITIONS"),
15 opHOCoordsAndNormals(hoCoordsAtGaussPtsF3, nOrmals_at_GaussPtF3,
16 tAngent1_at_GaussPtF3, tAngent2_at_GaussPtF3,
17 hoCoordsAtGaussPtsF4, nOrmals_at_GaussPtF4,
18 tAngent1_at_GaussPtF4, tAngent2_at_GaussPtF4) {
20 boost::shared_ptr<BaseFunction>(new FlatPrismPolynomialBase());
22 "Problem with creation data on element")
23}
24
27
28 if (numeredEntFiniteElementPtr->getEntType() != MBPRISM)
30
34
36 int num_nodes;
37 const EntityHandle *conn;
38 CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
39 {
40 coords.resize(num_nodes * 3, false);
41 CHKERR mField.get_moab().get_coords(conn, num_nodes,
42 &*coords.data().begin());
43
44 double diff_n[6];
45 CHKERR ShapeDiffMBTRI(diff_n);
46 normal.resize(6, false);
47 CHKERR ShapeFaceNormalMBTRI(diff_n, &coords[0], &normal[0]);
48 CHKERR ShapeFaceNormalMBTRI(diff_n, &coords[9], &normal[3]);
49 aRea[0] = cblas_dnrm2(3, &normal[0], 1) * 0.5;
50 aRea[1] = cblas_dnrm2(3, &normal[3], 1) * 0.5;
51 }
52
54
55 // H1
56 if ((dataH1.spacesOnEntities[MBVERTEX]).test(H1)) {
57 CHKERR getEntitySense<MBEDGE>(dataH1);
58 CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
59 CHKERR getEntitySense<MBTRI>(dataH1);
60 CHKERR getEntityDataOrder<MBTRI>(dataH1, H1);
61 }
62
63 // H1
64 if ((dataH1.spacesOnEntities[MBEDGE]).test(HCURL)) {
65 CHKERR getEntitySense<MBEDGE>(data_curl);
66 CHKERR getEntityDataOrder<MBEDGE>(data_curl, HCURL);
67 CHKERR getEntitySense<MBTRI>(data_curl);
68 CHKERR getEntityDataOrder<MBTRI>(data_curl, HCURL);
69 }
70
71 // Hdiv
72 if ((dataH1.spacesOnEntities[MBTRI]).test(HDIV)) {
73 CHKERR getEntitySense<MBTRI>(data_div);
74 CHKERR getEntityDataOrder<MBTRI>(data_div, HDIV);
75 }
76
77 // L2
78 if ((dataH1.spacesOnEntities[MBTRI]).test(L2)) {
79 CHKERR getEntitySense<MBTRI>(data_l2);
80 CHKERR getEntityDataOrder<MBTRI>(data_l2, L2);
81 }
82
83 int order_data = getMaxDataOrder();
84 int order_row = getMaxRowOrder();
85 int order_col = getMaxColOrder();
86 int rule = getRule(order_row, order_col, order_data);
87 int nb_gauss_pts;
88 // int rule = getRule(order);
89 if (rule >= 0) {
90 if (rule < QUAD_2D_TABLE_SIZE) {
91 if (QUAD_2D_TABLE[rule]->dim != 2) {
92 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
93 }
94 if (QUAD_2D_TABLE[rule]->order < rule) {
95 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
96 "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
97 }
98 nb_gauss_pts = QUAD_2D_TABLE[rule]->npoints;
99 gaussPts.resize(3, nb_gauss_pts, false);
100 cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[1], 3,
101 &gaussPts(0, 0), 1);
102 cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[2], 3,
103 &gaussPts(1, 0), 1);
104 cblas_dcopy(nb_gauss_pts, QUAD_2D_TABLE[rule]->weights, 1,
105 &gaussPts(2, 0), 1);
106 dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 3,
107 false);
108 double *shape_ptr =
109 &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
110 cblas_dcopy(3 * nb_gauss_pts, QUAD_2D_TABLE[rule]->points, 1, shape_ptr,
111 1);
112 } else {
113 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
114 "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
115 nb_gauss_pts = 0;
116 }
117 } else {
118 CHKERR setGaussPts(order_row, order_col, order_data);
119 nb_gauss_pts = gaussPts.size2();
120 dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 3,
121 false);
122 if (nb_gauss_pts) {
124 &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin(),
125 &gaussPts(0, 0), &gaussPts(1, 0), nb_gauss_pts);
126 }
127 }
128
129 if (nb_gauss_pts == 0)
131
132 {
133 coordsAtGaussPts.resize(nb_gauss_pts, 6, false);
134 for (int gg = 0; gg < nb_gauss_pts; gg++) {
135 for (int dd = 0; dd < 3; dd++) {
136 coordsAtGaussPts(gg, dd) = cblas_ddot(
137 3, &dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0), 1,
138 &coords[dd], 3);
139 coordsAtGaussPts(gg, 3 + dd) = cblas_ddot(
140 3, &dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0), 1,
141 &coords[9 + dd], 3);
142 }
143 }
144 }
145
146 for (int space = HCURL; space != LASTSPACE; ++space)
147 if (dataOnElement[space]) {
148 dataH1.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
149 dataOnElement[H1]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
150 }
151
152 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
153 if (dataH1.bAse.test(b)) {
154 switch (static_cast<FieldApproximationBase>(b)) {
157 if (dataH1.spacesOnEntities[MBVERTEX].test(H1)) {
159 gaussPts,
160 boost::shared_ptr<BaseFunctionCtx>(new FlatPrismPolynomialBaseCtx(
162 H1, static_cast<FieldApproximationBase>(b), NOBASE)));
163 }
164 #ifndef NDEBUG
165 if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
166 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
167 "Not yet implemented");
168 }
169 if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
170 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
171 "Not yet implemented");
172 }
173 if (dataH1.spacesOnEntities[MBTET].test(L2)) {
174 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
175 "Not yet implemented");
176 }
177#endif
178 break;
179 default:
180 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
181 "Not yet implemented");
182 }
183 }
184 }
185
186 auto check_field = [&]() {
187 auto field_it =
189 if (field_it != fieldsPtr->get<FieldName_mi_tag>().end())
190 if ((numeredEntFiniteElementPtr->getBitFieldIdData() &
191 (*field_it)->getId())
192 .any())
193 return true;
194 return false;
195 };
196
197 // Check if field meshPositionsFieldName exist
198 if (check_field()) {
199
200 hoCoordsAtGaussPtsF3.resize(nb_gauss_pts, 3, false);
201 nOrmals_at_GaussPtF3.resize(nb_gauss_pts, 3, false);
202 tAngent1_at_GaussPtF3.resize(nb_gauss_pts, 3, false);
203 tAngent2_at_GaussPtF3.resize(nb_gauss_pts, 3, false);
204 hoCoordsAtGaussPtsF4.resize(nb_gauss_pts, 3, false);
205 nOrmals_at_GaussPtF4.resize(nb_gauss_pts, 3, false);
206 tAngent1_at_GaussPtF4.resize(nb_gauss_pts, 3, false);
207 tAngent2_at_GaussPtF4.resize(nb_gauss_pts, 3, false);
208 const auto bit_number = mField.get_field_bit_number(meshPositionsFieldName);
209 CHKERR getNodesFieldData(dataH1, bit_number);
210 CHKERR getEntityFieldData(dataH1, bit_number, MBEDGE);
211 CHKERR getEntityFieldData(dataH1, bit_number, MBEDGE);
214 } else {
215 hoCoordsAtGaussPtsF3.resize(0, 0, false);
216 nOrmals_at_GaussPtF3.resize(0, 0, false);
217 tAngent1_at_GaussPtF3.resize(0, 0, false);
218 tAngent2_at_GaussPtF3.resize(0, 0, false);
219 hoCoordsAtGaussPtsF4.resize(0, 0, false);
220 nOrmals_at_GaussPtF4.resize(0, 0, false);
221 tAngent1_at_GaussPtF4.resize(0, 0, false);
222 tAngent2_at_GaussPtF4.resize(0, 0, false);
223 }
224
225#ifndef NDEBUG
226 if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
227 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented yet");
228 }
229#endif
230
231 // Iterate over operators
233
235}
236
237MoFEMErrorCode FlatPrismElementForcesAndSourcesCore::UserDataOperator::setPtrFE(
240 if (!(ptrFE = dynamic_cast<FlatPrismElementForcesAndSourcesCore *>(ptr)))
241 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
242 "User operator and finite element do not work together");
244}
245
246} // namespace MoFEM
FieldApproximationBase
approximation base
Definition: definitions.h:58
@ LASTBASE
Definition: definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:62
@ NOBASE
Definition: definitions.h:59
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:595
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
@ L2
field with C-1 continuity
Definition: definitions.h:88
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:89
@ 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
#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
constexpr int order
const int dim
PetscErrorCode ShapeFaceNormalMBTRI(double *diffN, const double *coords, double *normal)
Definition: fem_tools.c:229
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:194
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
#define QUAD_2D_TABLE_SIZE
Definition: quad.h:174
static QUAD *const QUAD_2D_TABLE[]
Definition: quad.h:175
const Field_multiIndex * fieldsPtr
raw pointer to fields container
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode opRhs(EntitiesFieldData &data, const bool error_if_no_base=false)
Deprecated interface functions.
data structure for finite element entity
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
std::bitset< LASTBASE > bAse
bases on element
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
MultiIndex Tag for field name.
MoFEMErrorCode operator()()
function is run for every finite element
Class used to pass element data to calculate base functions on flat prism.
Calculate base functions on tetrahedralFIXME: Need moab and mofem finite element structure to work (t...
structure to get information form mofem into EntitiesFieldData
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
int getMaxRowOrder() const
Get max order of approximation for field in rows.
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
MoFEMErrorCode getEntityFieldData(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
MatrixDouble gaussPts
Matrix of integration points.
int getMaxColOrder() const
Get max order of approximation for field in columns.
MoFEMErrorCode getNodesFieldData(EntitiesFieldData &data, const int bit_number) const
Get data on nodes.
int getMaxDataOrder() const
Get max order of approximation for data fields.
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.
int npoints
Definition: quad.h:29