v0.14.0
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)
 

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 410 of file PostProcBrokenMeshInMoabBase.cpp.

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

◆ PostProcGenerateRefMeshBase()

MoFEM::PostProcGenerateRefMeshBase::PostProcGenerateRefMeshBase

Definition at line 13 of file PostProcBrokenMeshInMoabBase.cpp.

14  : hoNodes(PETSC_TRUE), defMaxLevel(0), countEle(0), countVertEle(0),
15  nbVertices(0) {}

The documentation for this struct was generated from the following files:
EntityHandle
MoFEM::PostProcGenerateRefMeshBase::nbVertices
int nbVertices
Definition: PostProcBrokenMeshInMoabBase.hpp:40
N_MBQUAD2
#define N_MBQUAD2(x, y)
quad shape function
Definition: fem_tools.h:59
zeta
double zeta
Viscous hardening.
Definition: plastic.cpp:177
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
MoFEM::PostProcGenerateRefMeshBase::levelShapeFunctions
std::vector< MatrixDouble > levelShapeFunctions
Definition: PostProcBrokenMeshInMoabBase.hpp:26
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::PostProcGenerateRefMeshBase::hoNodes
PetscBool hoNodes
Definition: PostProcBrokenMeshInMoabBase.hpp:49
N_MBQUAD1
#define N_MBQUAD1(x, y)
quad shape function
Definition: fem_tools.h:58
MoFEM::PostProcGenerateRefMeshBase::levelGaussPtsOnRefMesh
std::vector< MatrixDouble > levelGaussPtsOnRefMesh
Definition: PostProcBrokenMeshInMoabBase.hpp:29
MoFEM::PostProcGenerateRefMeshBase::countVertEle
int countVertEle
Definition: PostProcBrokenMeshInMoabBase.hpp:38
MoFEM::PostProcGenerateRefMeshBase::defMaxLevel
int defMaxLevel
Definition: PostProcBrokenMeshInMoabBase.hpp:50
Range
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::PostProcGenerateRefMeshBase::levelRef
std::vector< ublas::matrix< int > > levelRef
Definition: PostProcBrokenMeshInMoabBase.hpp:30
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
MoFEM::PostProcGenerateRefMeshBase::countEle
int countEle
Definition: PostProcBrokenMeshInMoabBase.hpp:37
N_MBQUAD0
#define N_MBQUAD0(x, y)
quad shape function
Definition: fem_tools.h:57
N_MBQUAD3
#define N_MBQUAD3(x, y)
quad shape function
Definition: fem_tools.h:60
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346