v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | List of all members
MoFEM::PostProcGenerateRefMesh< MBQUAD > Struct Reference

#include <src/post_proc/PostProcBrokenMeshInMoabBase.hpp>

Inheritance diagram for MoFEM::PostProcGenerateRefMesh< MBQUAD >:
[legend]
Collaboration diagram for MoFEM::PostProcGenerateRefMesh< MBQUAD >:
[legend]

Public Member Functions

MoFEMErrorCode generateReferenceElementMesh ()
 
 PostProcGenerateRefMeshBase ()
 
- Public Member Functions inherited from MoFEM::PostProcGenerateRefMeshBase
 PostProcGenerateRefMeshBase ()
 
virtual ~PostProcGenerateRefMeshBase ()=default
 
MoFEMErrorCode getOptions ()
 
virtual MoFEMErrorCode generateReferenceElementMesh ()=0
 

Additional Inherited Members

- Public Attributes inherited from MoFEM::PostProcGenerateRefMeshBase
std::vector< MatrixDoublelevelShapeFunctions
 
std::vector< MatrixDoublelevelGaussPtsOnRefMesh
 
std::vector< ublas::matrix< int > > levelRef
 
EntityHandle startingVertEleHandle
 
std::vector< double * > verticesOnEleArrays
 
EntityHandle startingEleHandle
 
EntityHandleeleConn
 
int countEle
 
int countVertEle
 
int nbVertices
 
int nbEles
 
PetscBool hoNodes
 
int defMaxLevel
 
std::string optPrefix
 

Detailed Description

Definition at line 84 of file PostProcBrokenMeshInMoabBase.hpp.

Member Function Documentation

◆ generateReferenceElementMesh()

MoFEMErrorCode MoFEM::PostProcGenerateRefMesh< MBQUAD >::generateReferenceElementMesh ( )
virtual

Implements MoFEM::PostProcGenerateRefMeshBase.

Definition at line 411 of file PostProcBrokenMeshInMoabBase.cpp.

411 {
413
415#ifndef NDEBUG
416 if (defMaxLevel > 0)
417 MOFEM_LOG("WORLD", Sev::warning)
418 << "Refinement for quad is not implemented";
419#endif
420
421 moab::Core core_ref;
422 moab::Interface &moab_ref = core_ref;
423
424 auto create_reference_element = [&moab_ref]() {
426 constexpr double base_coords[] = {
427
428 0, 0,
429 0,
430
431 1, 0,
432 0,
433
434 1, 1,
435 0,
436
437 0, 1,
438 0
439
440 };
441 EntityHandle nodes[4];
442 for (int nn = 0; nn < 4; nn++)
443 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
444 EntityHandle quad;
445 CHKERR moab_ref.create_element(MBQUAD, nodes, 4, quad);
447 };
448
449 auto add_ho_nodes = [&]() {
451 Range quads;
452 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
453 EntityHandle meshset;
454 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
455 CHKERR moab_ref.add_entities(meshset, quads);
456 CHKERR moab_ref.convert_entities(meshset, true, true, true);
457 CHKERR moab_ref.delete_entities(&meshset, 1);
459 };
460
461 auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
463 Range quads;
464 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
465 Range quads_nodes;
466 CHKERR moab_ref.get_connectivity(quads, quads_nodes, false);
467 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
468 gauss_pts.resize(quads_nodes.size(), 4, false);
469 size_t gg = 0;
470 for (auto node : quads_nodes) {
471 CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
472 little_map[node] = gg;
473 ++gg;
474 }
475 gauss_pts = trans(gauss_pts);
477 };
478
479 auto set_ref_quads = [&](std::map<EntityHandle, int> &little_map) {
481 Range quads;
482 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
483 size_t hh = 0;
484 auto &ref_quads = levelRef[0];
485 for (auto quad : quads) {
486 const EntityHandle *conn;
487 int num_nodes;
488 CHKERR moab_ref.get_connectivity(quad, conn, num_nodes, false);
489 if (ref_quads.size2() != num_nodes) {
490 ref_quads.resize(quads.size(), num_nodes);
491 }
492 for (int nn = 0; nn != num_nodes; ++nn) {
493 ref_quads(hh, nn) = little_map[conn[nn]];
494 }
495 ++hh;
496 }
498 };
499
500 auto set_shape_functions = [&]() {
502 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
503 auto &shape_functions = levelShapeFunctions[0];
504 const auto nb_gauss_pts = gauss_pts.size2();
505 shape_functions.resize(nb_gauss_pts, 4);
506 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
507 const double ksi = gauss_pts(0, gg);
508 const double zeta = gauss_pts(1, gg);
509 shape_functions(gg, 0) = N_MBQUAD0(ksi, zeta);
510 shape_functions(gg, 1) = N_MBQUAD1(ksi, zeta);
511 shape_functions(gg, 2) = N_MBQUAD2(ksi, zeta);
512 shape_functions(gg, 3) = N_MBQUAD3(ksi, zeta);
513 }
515 };
516
517 levelRef.resize(1);
518 levelGaussPtsOnRefMesh.resize(1);
519 levelShapeFunctions.resize(1);
520
521 CHKERR create_reference_element();
522 if (hoNodes)
523 CHKERR add_ho_nodes();
524 std::map<EntityHandle, int> little_map;
525 CHKERR set_gauss_pts(little_map);
526 CHKERR set_ref_quads(little_map);
527 CHKERR set_shape_functions();
528
530}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#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 N_MBQUAD3(x, y)
quad shape function
Definition: fem_tools.h:60
#define N_MBQUAD0(x, y)
quad shape function
Definition: fem_tools.h:57
#define N_MBQUAD2(x, y)
quad shape function
Definition: fem_tools.h:59
#define N_MBQUAD1(x, y)
quad shape function
Definition: fem_tools.h:58
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
double zeta
Definition: plastic.cpp:107
std::vector< ublas::matrix< int > > levelRef
PetscBool hoNodes
int defMaxLevel
std::vector< MatrixDouble > levelShapeFunctions
std::vector< MatrixDouble > levelGaussPtsOnRefMesh
MoFEMErrorCode getOptions()

◆ PostProcGenerateRefMeshBase()

MoFEM::PostProcGenerateRefMeshBase::PostProcGenerateRefMeshBase ( )

The documentation for this struct was generated from the following files: