v0.14.0
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
 
virtual MoFEMErrorCode getOptions (std::string prefix)
 
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
 

Detailed Description

Definition at line 83 of file PostProcBrokenMeshInMoabBase.hpp.

Member Function Documentation

◆ generateReferenceElementMesh()

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

Implements MoFEM::PostProcGenerateRefMeshBase.

Definition at line 408 of file PostProcBrokenMeshInMoabBase.cpp.

408 {
410
411#ifndef NDEBUG
412 if (defMaxLevel > 0)
413 MOFEM_LOG("WORLD", Sev::warning)
414 << "Refinement for quad is not implemented";
415#endif
416
417 moab::Core core_ref;
418 moab::Interface &moab_ref = core_ref;
419
420 auto create_reference_element = [&moab_ref]() {
422 constexpr double base_coords[] = {
423
424 0, 0,
425 0,
426
427 1, 0,
428 0,
429
430 1, 1,
431 0,
432
433 0, 1,
434 0
435
436 };
437 EntityHandle nodes[4];
438 for (int nn = 0; nn < 4; nn++)
439 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
440 EntityHandle quad;
441 CHKERR moab_ref.create_element(MBQUAD, nodes, 4, quad);
443 };
444
445 auto add_ho_nodes = [&]() {
447 Range quads;
448 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
449 EntityHandle meshset;
450 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
451 CHKERR moab_ref.add_entities(meshset, quads);
452 CHKERR moab_ref.convert_entities(meshset, true, true, true);
453 CHKERR moab_ref.delete_entities(&meshset, 1);
455 };
456
457 auto set_gauss_pts = [&](std::map<EntityHandle, int> &little_map) {
459 Range quads;
460 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
461 Range quads_nodes;
462 CHKERR moab_ref.get_connectivity(quads, quads_nodes, false);
463 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
464 gauss_pts.resize(quads_nodes.size(), 4, false);
465 size_t gg = 0;
466 for (auto node : quads_nodes) {
467 CHKERR moab_ref.get_coords(&node, 1, &gauss_pts(gg, 0));
468 little_map[node] = gg;
469 ++gg;
470 }
471 gauss_pts = trans(gauss_pts);
473 };
474
475 auto set_ref_quads = [&](std::map<EntityHandle, int> &little_map) {
477 Range quads;
478 CHKERR moab_ref.get_entities_by_type(0, MBQUAD, quads, true);
479 size_t hh = 0;
480 auto &ref_quads = levelRef[0];
481 for (auto quad : quads) {
482 const EntityHandle *conn;
483 int num_nodes;
484 CHKERR moab_ref.get_connectivity(quad, conn, num_nodes, false);
485 if (ref_quads.size2() != num_nodes) {
486 ref_quads.resize(quads.size(), num_nodes);
487 }
488 for (int nn = 0; nn != num_nodes; ++nn) {
489 ref_quads(hh, nn) = little_map[conn[nn]];
490 }
491 ++hh;
492 }
494 };
495
496 auto set_shape_functions = [&]() {
498 auto &gauss_pts = levelGaussPtsOnRefMesh[0];
499 auto &shape_functions = levelShapeFunctions[0];
500 const auto nb_gauss_pts = gauss_pts.size2();
501 shape_functions.resize(nb_gauss_pts, 4);
502 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
503 const double ksi = gauss_pts(0, gg);
504 const double zeta = gauss_pts(1, gg);
505 shape_functions(gg, 0) = N_MBQUAD0(ksi, zeta);
506 shape_functions(gg, 1) = N_MBQUAD1(ksi, zeta);
507 shape_functions(gg, 2) = N_MBQUAD2(ksi, zeta);
508 shape_functions(gg, 3) = N_MBQUAD3(ksi, zeta);
509 }
511 };
512
513 levelRef.resize(1);
514 levelGaussPtsOnRefMesh.resize(1);
515 levelShapeFunctions.resize(1);
516
517 CHKERR create_reference_element();
518 if (hoNodes)
519 CHKERR add_ho_nodes();
520 std::map<EntityHandle, int> little_map;
521 CHKERR set_gauss_pts(little_map);
522 CHKERR set_ref_quads(little_map);
523 CHKERR set_shape_functions();
524
526}
#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
Viscous hardening.
Definition: plastic.cpp:173
std::vector< ublas::matrix< int > > levelRef
PetscBool hoNodes
int defMaxLevel
std::vector< MatrixDouble > levelShapeFunctions
std::vector< MatrixDouble > levelGaussPtsOnRefMesh

◆ PostProcGenerateRefMeshBase()

MoFEM::PostProcGenerateRefMeshBase::PostProcGenerateRefMeshBase ( )

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