v0.15.0
Loading...
Searching...
No Matches
Functions | Variables
forces_and_sources_h1_continuity_check.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 
static const double face_coords [4][9]
 
static const double edge_coords [6][6]
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 22 of file forces_and_sources_h1_continuity_check.cpp.

22 {
23
24 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
25
26 try {
27
28 moab::Core mb_instance;
29 moab::Interface &moab = mb_instance;
30 int rank;
31 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
32
33 PetscBool flg = PETSC_TRUE;
34 char mesh_file_name[255];
35#if PETSC_VERSION_GE(3, 6, 4)
36 ierr = PetscOptionsGetString(PETSC_NULLPTR, "", "-my_file", mesh_file_name,
37 255, &flg);
39#else
40 ierr = PetscOptionsGetString(PETSC_NULLPTR, PETSC_NULLPTR, "-my_file",
41 mesh_file_name, 255, &flg);
43#endif
44 if (flg != PETSC_TRUE) {
45 SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
46 }
47 const char *option;
48 option = "";
49 CHKERR moab.load_file(mesh_file_name, 0, option);
50
51 // Create MoFEM (Joseph) database
52 MoFEM::Core core(moab);
53 MoFEM::Interface &m_field = core;
54
55 // set entitities bit level
56 BitRefLevel bit_level0;
57 bit_level0.set(0);
58 EntityHandle meshset_level0;
59 CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
60
61 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
62 0, 3, bit_level0);
63
64
65 // Fields
66 CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
67 3);
68
69 CHKERR m_field.add_field("H1", H1, AINSWORTH_LEGENDRE_BASE, 1);
70
71
72 // FE
73 CHKERR m_field.add_finite_element("TET_FE");
74 CHKERR m_field.add_finite_element("TRI_FE");
75 CHKERR m_field.add_finite_element("SKIN_FE");
76 CHKERR m_field.add_finite_element("EDGE_FE");
77
78
79 // Define rows/cols and element data
80 CHKERR m_field.modify_finite_element_add_field_row("TET_FE", "H1");
81 CHKERR m_field.modify_finite_element_add_field_col("TET_FE", "H1");
82 CHKERR m_field.modify_finite_element_add_field_data("TET_FE", "H1");
84 "MESH_NODE_POSITIONS");
85
86
87 CHKERR m_field.modify_finite_element_add_field_row("SKIN_FE", "H1");
88 CHKERR m_field.modify_finite_element_add_field_col("SKIN_FE", "H1");
89 CHKERR m_field.modify_finite_element_add_field_data("SKIN_FE", "H1");
91 "MESH_NODE_POSITIONS");
92
93
94 CHKERR m_field.modify_finite_element_add_field_row("TRI_FE", "H1");
95 CHKERR m_field.modify_finite_element_add_field_col("TRI_FE", "H1");
96 CHKERR m_field.modify_finite_element_add_field_data("TRI_FE", "H1");
98 "MESH_NODE_POSITIONS");
99
100
101 CHKERR m_field.modify_finite_element_add_field_row("EDGE_FE", "H1");
102 CHKERR m_field.modify_finite_element_add_field_col("EDGE_FE", "H1");
103 CHKERR m_field.modify_finite_element_add_field_data("EDGE_FE", "H1");
105 "MESH_NODE_POSITIONS");
106
107 // Problem
108 CHKERR m_field.add_problem("TEST_PROBLEM");
109
110
111 // set finite elements for problem
112 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TET_FE");
113 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "SKIN_FE");
114 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TRI_FE");
115 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "EDGE_FE");
116
117 // set refinement level for problem
118 CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
119
120
121 // meshset consisting all entities in mesh
122 EntityHandle root_set = moab.get_root_set();
123 // add entities to field
124 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "H1");
125
126
127 // add entities to finite element
128 CHKERR
129 m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "TET_FE");
130
131
132 Range tets;
133 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
134 BitRefLevel().set(0), BitRefLevel().set(), MBTET, tets);
135
136 Skinner skin(&moab);
137 Range skin_faces; // skin faces from 3d ents
138 CHKERR skin.find_skin(0, tets, false, skin_faces);
139 CHKERR m_field.add_ents_to_finite_element_by_type(skin_faces, MBTRI,
140 "SKIN_FE");
141
142 Range faces;
143 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
144 BitRefLevel().set(0), BitRefLevel().set(), MBTRI, faces);
145
146 faces = subtract(faces, skin_faces);
147 CHKERR m_field.add_ents_to_finite_element_by_type(faces, MBTRI, "TRI_FE");
148
149 Range edges;
150 CHKERR moab.get_adjacencies(faces, 1, false, edges, moab::Interface::UNION);
151 CHKERR m_field.add_ents_to_finite_element_by_type(edges, MBEDGE, "EDGE_FE");
152
153 Range skin_edges;
154 CHKERR moab.get_adjacencies(skin_faces, 1, false, skin_edges,
155 moab::Interface::UNION);
156 CHKERR m_field.add_ents_to_finite_element_by_type(edges, MBEDGE, "EDGE_FE");
157
158
159 // set app. order
160 int order = 4;
161 CHKERR m_field.set_field_order(root_set, MBTET, "H1", order);
162 CHKERR m_field.set_field_order(root_set, MBTRI, "H1", order);
163 CHKERR m_field.set_field_order(root_set, MBEDGE, "H1", order);
164
165 CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "MESH_NODE_POSITIONS");
166 CHKERR m_field.set_field_order(0, MBVERTEX, "MESH_NODE_POSITIONS", 1);
167 CHKERR m_field.set_field_order(0, MBEDGE, "MESH_NODE_POSITIONS", 2);
168 CHKERR m_field.set_field_order(0, MBTRI, "MESH_NODE_POSITIONS", 2);
169 CHKERR m_field.set_field_order(0, MBTET, "MESH_NODE_POSITIONS", 2);
170
171 /****/
172 // build database
173 // build field
174 CHKERR m_field.build_fields();
175
176 // build finite elemnts
178
179 // build adjacencies
180 CHKERR m_field.build_adjacencies(bit_level0);
181
182
183 ProblemsManager *prb_mng_ptr;
184 CHKERR m_field.getInterface(prb_mng_ptr);
185
186 // build problem
187 CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
188
189
190 // project geometry form 10 node tets on higher order approx. functions
191 Projection10NodeCoordsOnField ent_method(m_field, "MESH_NODE_POSITIONS");
192 CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method);
193
194 /****/
195 // mesh partitioning
196 // partition
197 CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
198 CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
199
200 // what are ghost nodes, see Petsc Manual
201 CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
202
203
204 Vec v;
205 CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
206 ROW, &v);
207 CHKERR VecSetRandom(v, PETSC_NULLPTR);
208
209 CHKERR m_field.getInterface<VecManager>()->setLocalGhostVector(
210 "TEST_PROBLEM", ROW, v, INSERT_VALUES, SCATTER_REVERSE);
211
212 CHKERR VecDestroy(&v);
213
214
215 typedef tee_device<std::ostream, std::ofstream> TeeDevice;
216 typedef stream<TeeDevice> TeeStream;
217 std::ofstream ofs("forces_and_sources_hdiv_continuity_check.txt");
218 TeeDevice my_tee(std::cout, ofs);
219 TeeStream my_split(my_tee);
220
221 struct OpTetFluxes
222 : public VolumeElementForcesAndSourcesCore::UserDataOperator {
223
224 MoFEM::Interface &mField;
225 Tag tH;
226
227 OpTetFluxes(MoFEM::Interface &m_field, Tag th)
229 "H1", UserDataOperator::OPROW),
230 mField(m_field), tH(th) {}
231
232 MoFEMErrorCode doWork(int side, EntityType type,
235
236 if (data.getFieldData().size() == 0)
238
239 if (type == MBTRI) {
240
241 boost::shared_ptr<const NumeredEntFiniteElement> mofem_fe =
242 getNumeredEntFiniteElementPtr();
243 SideNumber_multiIndex &side_table = mofem_fe->getSideNumberTable();
244 EntityHandle face = side_table.get<1>()
245 .find(boost::make_tuple(type, side))
246 ->get()
247 ->ent;
248 int sense = side_table.get<1>()
249 .find(boost::make_tuple(type, side))
250 ->get()
251 ->sense;
252
253 double t = 0;
254 int nb_dofs = data.getN().size2();
255 for (int dd = 0; dd < nb_dofs; dd++) {
256 t += data.getN(side)(dd) * data.getFieldData()[dd];
257 }
258
259 double *t_ptr;
260 CHKERR mField.get_moab().tag_get_by_ptr(tH, &face, 1,
261 (const void **)&t_ptr);
262
263 *t_ptr += sense * t;
264 }
265
266 if (type == MBEDGE) {
267
268 boost::shared_ptr<const NumeredEntFiniteElement> mofem_fe =
269 getNumeredEntFiniteElementPtr();
270 SideNumber_multiIndex &side_table = mofem_fe->getSideNumberTable();
271 EntityHandle edge = side_table.get<1>()
272 .find(boost::make_tuple(type, side))
273 ->get()
274 ->ent;
275 Range adj_tets;
276 CHKERR mField.get_moab().get_adjacencies(&edge, 1, 3, false, adj_tets,
277 moab::Interface::UNION);
278 const int nb_adj_tets = adj_tets.size();
279
280 double t = 0;
281 int nb_dofs = data.getN().size2();
282 for (int dd = 0; dd < nb_dofs; dd++) {
283 t += data.getN(4 + side)(dd) * data.getFieldData()[dd];
284 }
285
286 double *t_ptr;
287 CHKERR mField.get_moab().tag_get_by_ptr(tH, &edge, 1,
288 (const void **)&t_ptr);
289
290 *t_ptr += t / nb_adj_tets;
291 }
292
294 }
295 };
296
297 struct MyTetFE : public VolumeElementForcesAndSourcesCore {
298
299 MyTetFE(MoFEM::Interface &m_field)
301 int getRule(int order) { return -1; };
302
303 MatrixDouble N_tri;
306
307 N_tri.resize(1, 3);
308 CHKERR ShapeMBTRI(&N_tri(0, 0), G_TRI_X1, G_TRI_Y1, 1);
309
310 gaussPts.resize(4, 4 + 6);
311 int ff = 0;
312 for (; ff < 4; ff++) {
313 int dd = 0;
314 for (; dd < 3; dd++) {
315 gaussPts(dd, ff) =
316 cblas_ddot(3, &N_tri(0, 0), 1, &face_coords[ff][dd], 3);
317 }
318 gaussPts(3, ff) = G_TRI_W1[0];
319 }
320
321 int ee = 0;
322 for (; ee < 6; ee++) {
323 int dd = 0;
324 for (; dd < 3; dd++) {
325 gaussPts(dd, 4 + ee) =
326 (edge_coords[ee][0 + dd] + edge_coords[ee][3 + dd]) * 0.5;
327 }
328 gaussPts(3, 4 + ee) = 1;
329 }
330
332 }
333 };
334
335 struct OpFacesSkinFluxes
336 : public FaceElementForcesAndSourcesCore::UserDataOperator {
337
338 MoFEM::Interface &mField;
339 Tag tH1, tH2;
340 TeeStream &mySplit;
341
342 OpFacesSkinFluxes(MoFEM::Interface &m_field, Tag th1, Tag th2,
343 TeeStream &my_split)
345 "H1", UserDataOperator::OPROW),
346 mField(m_field), tH1(th1), tH2(th2), mySplit(my_split) {}
347
348 MoFEMErrorCode doWork(int side, EntityType type,
351
352 if (type != MBTRI)
354 EntityHandle face = getNumeredEntFiniteElementPtr()->getEnt();
355
356 double *t_ptr;
357 CHKERR mField.get_moab().tag_get_by_ptr(tH1, &face, 1,
358 (const void **)&t_ptr);
359
360 double *tn_ptr;
361 CHKERR mField.get_moab().tag_get_by_ptr(tH2, &face, 1,
362 (const void **)&tn_ptr);
363
364
365 *tn_ptr = *t_ptr;
366
367 int dd = 0;
368 int nb_dofs = data.getN().size2();
369 for (; dd < nb_dofs; dd++) {
370 double val = data.getFieldData()[dd];
371 *tn_ptr += -data.getN()(0, dd) * val;
372 }
373
374 const double eps = 1e-8;
375 if (fabs(*tn_ptr) > eps) {
376 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
377 "H1 continuity failed %6.4e", *tn_ptr);
378 }
379
380 mySplit.precision(5);
381 mySplit << face << " " << /*std::fixed <<*/ fabs(*tn_ptr) << std::endl;
382
384 }
385 };
386
387 struct OpFacesFluxes
388 : public FaceElementForcesAndSourcesCore::UserDataOperator {
389
390 MoFEM::Interface &mField;
391 Tag tH1, tH2;
392 TeeStream &mySplit;
393
394 OpFacesFluxes(MoFEM::Interface &m_field, Tag _th1, Tag _th2,
395 TeeStream &my_split)
397 "H1", UserDataOperator::OPROW),
398 mField(m_field), tH1(_th1), tH2(_th2), mySplit(my_split) {}
399
400 MoFEMErrorCode doWork(int side, EntityType type,
403
404 if (type != MBTRI)
406
407 EntityHandle face = getNumeredEntFiniteElementPtr()->getEnt();
408
409 double *t_ptr;
410 CHKERR mField.get_moab().tag_get_by_ptr(tH1, &face, 1,
411 (const void **)&t_ptr);
412
413 double *tn_ptr;
414 CHKERR mField.get_moab().tag_get_by_ptr(tH2, &face, 1,
415 (const void **)&tn_ptr);
416
417
418 *tn_ptr = *t_ptr;
419
420 const double eps = 1e-8;
421 if (fabs(*tn_ptr) > eps) {
422 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
423 "H1 continuity failed %6.4e", *tn_ptr);
424 }
425
426 mySplit.precision(5);
427
428 mySplit << face << " " << /*std::fixed <<*/ fabs(*tn_ptr) << std::endl;
429
431 }
432 };
433
434 struct MyTriFE : public FaceElementForcesAndSourcesCore {
435
436 MyTriFE(MoFEM::Interface &m_field)
438 int getRule(int order) { return -1; };
439
442
443 gaussPts.resize(3, 1);
444 gaussPts(0, 0) = G_TRI_X1[0];
445 gaussPts(1, 0) = G_TRI_Y1[0];
446 gaussPts(2, 0) = G_TRI_W1[0];
447
449 }
450 };
451
452 struct OpEdgesFluxes
453 : public EdgeElementForcesAndSourcesCore::UserDataOperator {
454
455 MoFEM::Interface &mField;
456 Tag tH1, tH2;
457 TeeStream &mySplit;
458
459 OpEdgesFluxes(MoFEM::Interface &m_field, Tag _th1, Tag _th2,
460 TeeStream &my_split)
462 "H1", UserDataOperator::OPROW),
463 mField(m_field), tH1(_th1), tH2(_th2), mySplit(my_split) {}
464
465 MoFEMErrorCode doWork(int side, EntityType type,
468
469 if (type != MBEDGE)
471
472 EntityHandle edge = getNumeredEntFiniteElementPtr()->getEnt();
473
474 double *t_ptr;
475 CHKERR mField.get_moab().tag_get_by_ptr(tH1, &edge, 1,
476 (const void **)&t_ptr);
477
478 double *tn_ptr;
479 CHKERR mField.get_moab().tag_get_by_ptr(tH2, &edge, 1,
480 (const void **)&tn_ptr);
481
482
483 *tn_ptr = *t_ptr;
484
485 double tn = 0;
486 int nb_dofs = data.getN().size2();
487 int dd = 0;
488 for (; dd < nb_dofs; dd++) {
489 double val = data.getFieldData()[dd];
490 tn += data.getN()(0, dd) * val;
491 }
492
493 *tn_ptr -= tn;
494
495 mySplit.precision(5);
496
497 mySplit << edge << " " << fabs(*tn_ptr) << std::endl;
498
500 }
501 };
502
503 struct MyEdgeFE : public EdgeElementForcesAndSourcesCore {
504
505 MyEdgeFE(MoFEM::Interface &m_field)
507 int getRule(int order) { return -1; };
508
511
512 gaussPts.resize(2, 1);
513 gaussPts(0, 0) = 0.5;
514 gaussPts(1, 0) = 1;
515
517 }
518 };
519
520 MyTetFE tet_fe(m_field);
521 MyTriFE tri_fe(m_field);
522 MyTriFE skin_fe(m_field);
523 MyEdgeFE edge_fe(m_field);
524
525 Tag th1;
526 double def_val[] = {0, 0, 0};
527 CHKERR moab.tag_get_handle("T", 3, MB_TYPE_DOUBLE, th1,
528 MB_TAG_CREAT | MB_TAG_SPARSE, &def_val);
529
530 tet_fe.getOpPtrVector().push_back(
531 new OpCalculateHOCoords("MESH_NODE_POSITIONS"));
532 tet_fe.getOpPtrVector().push_back(new OpTetFluxes(m_field, th1));
533
534 Tag th2;
535 CHKERR moab.tag_get_handle("TN", 1, MB_TYPE_DOUBLE, th2,
536 MB_TAG_CREAT | MB_TAG_SPARSE, &def_val);
537
538 tri_fe.getOpPtrVector().push_back(
539 new OpCalculateHOCoords("MESH_NODE_POSITIONS"));
540 tri_fe.getOpPtrVector().push_back(
541 new OpGetHONormalsOnFace("MESH_NODE_POSITIONS"));
542 tri_fe.getOpPtrVector().push_back(
543 new OpFacesFluxes(m_field, th1, th2, my_split));
544 skin_fe.getOpPtrVector().push_back(
545 new OpFacesSkinFluxes(m_field, th1, th2, my_split));
546
547 edge_fe.getOpPtrVector().push_back(
548 new OpGetHOTangentsOnEdge("MESH_NODE_POSITIONS"));
549 edge_fe.getOpPtrVector().push_back(
550 new OpCalculateHOCoords("MESH_NODE_POSITIONS"));
551 edge_fe.getOpPtrVector().push_back(
552 new OpEdgesFluxes(m_field, th1, th2, my_split));
553
554 for (Range::iterator fit = faces.begin(); fit != faces.end(); fit++) {
555 CHKERR moab.tag_set_data(th1, &*fit, 1, &def_val);
556
557 CHKERR moab.tag_set_data(th2, &*fit, 1, &def_val);
558
559 }
560
561 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TET_FE", tet_fe);
562 my_split << "internal\n";
563 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TRI_FE", tri_fe);
564 my_split << "skin\n";
565 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "SKIN_FE", skin_fe);
566 my_split << "edges\n";
567 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "EDGE_FE", edge_fe);
568
569
570 EntityHandle meshset;
571 CHKERR moab.create_meshset(MESHSET_SET, meshset);
572 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
573 BitRefLevel().set(0), BitRefLevel().set(), MBTRI, meshset);
574
575 CHKERR moab.write_file("out.vtk", "VTK", "", &meshset, 1);
576
577 }
579
581}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static const double eps
@ ROW
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ H1
continuous field
Definition definitions.h:85
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr int order
static const double G_TRI_Y1[]
Definition fem_tools.h:313
static const double G_TRI_W1[]
Definition fem_tools.h:314
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition fem_tools.c:182
static const double G_TRI_X1[]
Definition fem_tools.h:312
static const double face_coords[4][9]
static const double edge_coords[6][6]
tee_device< std::ostream, std::ofstream > TeeDevice
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.
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
const double v
phase velocity of light in medium (cm/ns)
char mesh_file_name[255]
const FTensor::Tensor2< T, Dim, Dim > Vec
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)
Definition ddTensor0.hpp:33
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
constexpr double t
plate stiffness
Definition plate.cpp:58
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MoFEMErrorCode add_field(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
Core (interface) class.
Definition Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorDouble & getFieldData() const
get dofs values
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
MatrixDouble gaussPts
Matrix of integration points.
Calculate HO coordinates at gauss points.
Calculate normals at Gauss points of triangle element.
Calculate tangent vector on edge form HO geometry approximation.
Problem manager is used to build and partition problems.
Projection of edge entities with one mid-node on hierarchical basis.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)

Variable Documentation

◆ edge_coords

const double edge_coords[6][6]
static
Initial value:
= {
{0, 0, 0, 1, 0, 0}, {1, 0, 0, 0, 1, 0}, {0, 1, 0, 0, 0, 0},
{0, 0, 0, 0, 0, 1}, {1, 0, 0, 0, 0, 1}, {0, 1, 0, 0, 0, 1}}
Examples
prism_elements_from_surface.cpp.

Definition at line 18 of file forces_and_sources_h1_continuity_check.cpp.

18 {
19 {0, 0, 0, 1, 0, 0}, {1, 0, 0, 0, 1, 0}, {0, 1, 0, 0, 0, 0},
20 {0, 0, 0, 0, 0, 1}, {1, 0, 0, 0, 0, 1}, {0, 1, 0, 0, 0, 1}};

◆ face_coords

const double face_coords[4][9]
static
Initial value:
= {{0, 0, 0, 1, 0, 0, 0, 0, 1},
{1, 0, 0, 0, 1, 0, 0, 0, 1},
{0, 0, 0, 0, 1, 0, 0, 0, 1},
{0, 0, 0, 1, 0, 0, 0, 1, 0}}

Definition at line 13 of file forces_and_sources_h1_continuity_check.cpp.

13 {{0, 0, 0, 1, 0, 0, 0, 0, 1},
14 {1, 0, 0, 0, 1, 0, 0, 0, 1},
15 {0, 0, 0, 0, 1, 0, 0, 0, 1},
16 {0, 0, 0, 1, 0, 0, 0, 1, 0}};

◆ help

char help[] = "...\n\n"
static

Definition at line 11 of file forces_and_sources_h1_continuity_check.cpp.