v0.9.0
Classes | Public Member Functions | Public Attributes | Private Attributes | List of all members
PostProcFaceOnRefinedMesh Struct Reference

Postprocess on face. More...

#include <users_modules/basic_finite_elements/src/PostProcOnRefMesh.hpp>

Inheritance diagram for PostProcFaceOnRefinedMesh:
[legend]
Collaboration diagram for PostProcFaceOnRefinedMesh:
[legend]

Classes

struct  CommonData
 
struct  OpGetFieldGradientValuesOnSkin
 

Public Member Functions

 PostProcFaceOnRefinedMesh (MoFEM::Interface &m_field, bool six_node_post_proc_tris=true)
 
virtual ~PostProcFaceOnRefinedMesh ()
 
int getRule (int order)
 
MoFEMErrorCode generateReferenceElementMesh ()
 
MoFEMErrorCode setGaussPts (int order)
 
MoFEMErrorCode preProcess ()
 
MoFEMErrorCode postProcess ()
 
virtual PostProcCommonOnRefMesh::CommonDatagetCommonData ()
 
MoFEMErrorCode addFieldValuesGradientPostProcOnSkin (const std::string field_name, const std::string vol_fe_name, boost::shared_ptr< MatrixDouble > grad_mat_ptr=nullptr, bool save_on_tag=true)
 
- Public Member Functions inherited from PostProcTemplateOnRefineMesh< MoFEM::FaceElementForcesAndSourcesCore >
 PostProcTemplateOnRefineMesh (MoFEM::Interface &m_field)
 
MoFEMErrorCode addFieldValuesPostProc (const std::string field_name, Vec v=PETSC_NULL)
 Add operator to post-process L2, H1, Hdiv, Hcurl field value. More...
 
MoFEMErrorCode addFieldValuesPostProc (const std::string field_name, const std::string tag_name, Vec v=PETSC_NULL)
 Add operator to post-process L2 or H1 field value. More...
 
MoFEMErrorCode addFieldValuesGradientPostProc (const std::string field_name, Vec v=PETSC_NULL)
 Add operator to post-process L2 or H1 field gradient. More...
 
MoFEMErrorCode addFieldValuesGradientPostProc (const std::string field_name, const std::string tag_name, Vec v=PETSC_NULL)
 Add operator to post-process L2 or H1 field gradient. More...
 
MoFEMErrorCode writeFile (const std::string file_name)
 wrote results in (MOAB) format, use "file_name.h5m" More...
 

Public Attributes

bool sixNodePostProcTris
 
std::map< EntityHandle, EntityHandleelementsMap
 
CommonData commonData
 
- Public Attributes inherited from PostProcTemplateOnRefineMesh< MoFEM::FaceElementForcesAndSourcesCore >
moab::Core coreMesh
 
moab::Interface & postProcMesh
 
std::vector< EntityHandlemapGaussPts
 

Private Attributes

MatrixDouble gaussPtsTri
 
MatrixDouble gaussPtsQuad
 

Detailed Description

Postprocess on face.

Examples
cell_forces.cpp, minimal_surface_area.cpp, and reaction_diffusion_equation.cpp.

Definition at line 734 of file PostProcOnRefMesh.hpp.

Constructor & Destructor Documentation

◆ PostProcFaceOnRefinedMesh()

PostProcFaceOnRefinedMesh::PostProcFaceOnRefinedMesh ( MoFEM::Interface m_field,
bool  six_node_post_proc_tris = true 
)

◆ ~PostProcFaceOnRefinedMesh()

virtual PostProcFaceOnRefinedMesh::~PostProcFaceOnRefinedMesh ( )
virtual

Definition at line 745 of file PostProcOnRefMesh.hpp.

745  {
746  ParallelComm *pcomm_post_proc_mesh =
747  ParallelComm::get_pcomm(&postProcMesh, MYPCOMM_INDEX);
748  if (pcomm_post_proc_mesh != NULL) {
749  delete pcomm_post_proc_mesh;
750  }
751  }
moab::Interface & postProcMesh
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:285

Member Function Documentation

◆ addFieldValuesGradientPostProcOnSkin()

MoFEMErrorCode PostProcFaceOnRefinedMesh::addFieldValuesGradientPostProcOnSkin ( const std::string  field_name,
const std::string  vol_fe_name,
boost::shared_ptr< MatrixDouble >  grad_mat_ptr = nullptr,
bool  save_on_tag = true 
)

Definition at line 799 of file PostProcOnRefMesh.hpp.

802  {
804 
805  if (!grad_mat_ptr)
806  grad_mat_ptr = boost::make_shared<MatrixDouble>();
807 
808  boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> my_side_fe =
809  boost::make_shared<VolumeElementForcesAndSourcesCoreOnSide>(mField);
810 
811  // check number of coefficients
812  auto field_ptr = mField.get_field_structure(field_name);
813  const int nb_coefficients = field_ptr->getNbOfCoeffs();
814 
815  switch (nb_coefficients) {
816  case 1:
817  my_side_fe->getOpPtrVector().push_back(
818  new OpCalculateScalarFieldGradient<3>(field_name, grad_mat_ptr));
819  break;
820  case 3:
821  my_side_fe->getOpPtrVector().push_back(
822  new OpCalculateVectorFieldGradient<3, 3>(field_name, grad_mat_ptr));
823  break;
824  default:
825  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
826  "field with that number of coefficients is not implemented");
827  }
828 
829  FaceElementForcesAndSourcesCore::getOpPtrVector().push_back(
830  new OpGetFieldGradientValuesOnSkin(
831  postProcMesh, mapGaussPts, field_name, field_name + "_GRAD",
832  my_side_fe, vol_fe_name, grad_mat_ptr, save_on_tag));
833 
835  }
moab::Interface & postProcMesh
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
std::vector< EntityHandle > mapGaussPts
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ generateReferenceElementMesh()

MoFEMErrorCode PostProcFaceOnRefinedMesh::generateReferenceElementMesh ( )
Examples
minimal_surface_area.cpp.

Definition at line 523 of file PostProcOnRefMesh.cpp.

523  {
525 
526  auto generate_tri = [&](auto &gauss_pts) {
528  gauss_pts.resize(3, 3, false);
529  gauss_pts.clear();
530  gauss_pts(0, 0) = 0;
531  gauss_pts(1, 0) = 0;
532  gauss_pts(0, 1) = 1;
533  gauss_pts(1, 1) = 0;
534  gauss_pts(0, 2) = 0;
535  gauss_pts(1, 2) = 1;
536 
537  moab::Core core_ref;
538  moab::Interface &moab_ref = core_ref;
539  const EntityHandle *conn;
540  int num_nodes;
541  EntityHandle tri_conn[3];
542  MatrixDouble coords(6, 3);
543  for (int gg = 0; gg != 3; gg++) {
544  coords(gg, 0) = gauss_pts(0, gg);
545  coords(gg, 1) = gauss_pts(1, gg);
546  coords(gg, 2) = 0;
547  CHKERR moab_ref.create_vertex(&coords(gg, 0), tri_conn[gg]);
548  }
549 
550  EntityHandle tri;
551  CHKERR moab_ref.create_element(MBTRI, tri_conn, 3, tri);
552 
553  if (sixNodePostProcTris) {
554  Range edges;
555  CHKERR moab_ref.get_adjacencies(&tri, 1, 1, true, edges);
556  EntityHandle meshset;
557  CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
558  CHKERR moab_ref.add_entities(meshset, &tri, 1);
559  CHKERR moab_ref.add_entities(meshset, edges);
560  CHKERR moab_ref.convert_entities(meshset, true, false, false);
561  CHKERR moab_ref.get_connectivity(tri, conn, num_nodes, false);
562  CHKERR moab_ref.get_coords(conn, num_nodes, &coords(0, 0));
563  gauss_pts.resize(3, num_nodes, false);
564  gauss_pts.clear();
565  for (int nn = 0; nn < 3; nn++) {
566  gauss_pts(0, nn) = coords(nn, 0);
567  gauss_pts(1, nn) = coords(nn, 1);
568  gauss_pts(0, 3 + nn) = coords(3 + nn, 0);
569  gauss_pts(1, 3 + nn) = coords(3 + nn, 1);
570  }
571  } else {
572  CHKERR moab_ref.get_connectivity(tri, conn, num_nodes, false);
573  CHKERR moab_ref.get_coords(conn, num_nodes, &coords(0, 0));
574  gauss_pts.resize(3, num_nodes, false);
575  gauss_pts.clear();
576  for (int nn = 0; nn < 3; nn++) {
577  gauss_pts(0, nn) = coords(nn, 0);
578  gauss_pts(1, nn) = coords(nn, 1);
579  }
580  }
582  };
583 
584  auto generate_quad = [&](auto &gauss_pts) {
586 
587  gauss_pts.resize(3, 4, false);
588  gauss_pts.clear();
589 
590  gauss_pts(0, 0) = 0;
591  gauss_pts(1, 0) = 0;
592  gauss_pts(0, 1) = 1;
593  gauss_pts(1, 1) = 0;
594  gauss_pts(0, 2) = 1;
595  gauss_pts(1, 2) = 1;
596  gauss_pts(0, 3) = 0;
597  gauss_pts(1, 3) = 1;
598 
599  moab::Core core_ref;
600  moab::Interface &moab_ref = core_ref;
601 
602  const EntityHandle *conn;
603  int num_nodes;
604  std::array<EntityHandle, 4> quad_conn;
605  MatrixDouble coords(8, 3);
606  for (int gg = 0; gg != 4; gg++) {
607  coords(gg, 0) = gauss_pts(0, gg);
608  coords(gg, 1) = gauss_pts(1, gg);
609  coords(gg, 2) = 0;
610  CHKERR moab_ref.create_vertex(&coords(gg, 0), quad_conn[gg]);
611  }
612 
613  EntityHandle quad;
614  CHKERR moab_ref.create_element(MBQUAD, quad_conn.data(), 4, quad);
615 
616  if (sixNodePostProcTris) {
617  Range edges;
618  CHKERR moab_ref.get_adjacencies(&quad, 1, 1, true, edges);
619  EntityHandle meshset;
620  CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
621  CHKERR moab_ref.add_entities(meshset, &quad, 1);
622  CHKERR moab_ref.add_entities(meshset, edges);
623  CHKERR moab_ref.convert_entities(meshset, true, false, false);
624  CHKERR moab_ref.get_connectivity(quad, conn, num_nodes, false);
625  CHKERR moab_ref.get_coords(conn, num_nodes, &coords(0, 0));
626  gauss_pts.resize(3, num_nodes, false);
627  gauss_pts.clear();
628  for (int nn = 0; nn != 4; nn++) {
629  gauss_pts(0, nn) = coords(nn, 0);
630  gauss_pts(1, nn) = coords(nn, 1);
631  gauss_pts(0, 4 + nn) = coords(4 + nn, 0);
632  gauss_pts(1, 4 + nn) = coords(4 + nn, 1);
633  }
634  } else {
635  CHKERR moab_ref.get_connectivity(quad, conn, num_nodes, false);
636  CHKERR moab_ref.get_coords(conn, num_nodes, &coords(0, 0));
637  gauss_pts.resize(3, num_nodes, false);
638  gauss_pts.clear();
639  for (int nn = 0; nn != 4; nn++) {
640  gauss_pts(0, nn) = coords(nn, 0);
641  gauss_pts(1, nn) = coords(nn, 1);
642  }
643  }
644 
646  };
647 
648  CHKERR generate_tri(gaussPtsTri);
649  CHKERR generate_quad(gaussPtsQuad);
650 
652 }
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:74
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
MatrixDouble gaussPtsTri
MatrixDouble gaussPtsQuad
bool sixNodePostProcTris
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ getCommonData()

virtual PostProcCommonOnRefMesh::CommonData& PostProcFaceOnRefinedMesh::getCommonData ( )
virtual

Reimplemented from PostProcTemplateOnRefineMesh< MoFEM::FaceElementForcesAndSourcesCore >.

Definition at line 767 of file PostProcOnRefMesh.hpp.

767  {
768  return commonData;
769  }
CommonData commonData

◆ getRule()

int PostProcFaceOnRefinedMesh::getRule ( int  order)

Definition at line 754 of file PostProcOnRefMesh.hpp.

754 { return -1; };

◆ postProcess()

MoFEMErrorCode PostProcFaceOnRefinedMesh::postProcess ( )

Definition at line 794 of file PostProcOnRefMesh.cpp.

794  {
796  ParallelComm *pcomm =
797  ParallelComm::get_pcomm(&mField.get_moab(), MYPCOMM_INDEX);
798  ParallelComm *pcomm_post_proc_mesh =
799  ParallelComm::get_pcomm(&postProcMesh, MYPCOMM_INDEX);
800  if (pcomm_post_proc_mesh == NULL) {
801  pcomm_post_proc_mesh = new ParallelComm(&postProcMesh, mField.get_comm());
802  }
803  Range tris;
804  CHKERR postProcMesh.get_entities_by_type(0, MBTRI, tris, false);
805  int rank = pcomm->rank();
806  Range::iterator pit = tris.begin();
807  for (; pit != tris.end(); pit++) {
808  CHKERR postProcMesh.tag_set_data(pcomm_post_proc_mesh->part_tag(), &*pit, 1,
809  &rank);
810  }
811  CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
813 }
moab::Interface & postProcMesh
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:285
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ preProcess()

MoFEMErrorCode PostProcFaceOnRefinedMesh::preProcess ( )

Definition at line 784 of file PostProcOnRefMesh.cpp.

784  {
786  ParallelComm *pcomm_post_proc_mesh =
787  ParallelComm::get_pcomm(&postProcMesh, MYPCOMM_INDEX);
788  if (pcomm_post_proc_mesh != NULL) {
789  delete pcomm_post_proc_mesh;
790  }
792 }
moab::Interface & postProcMesh
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:285

◆ setGaussPts()

MoFEMErrorCode PostProcFaceOnRefinedMesh::setGaussPts ( int  order)

Definition at line 654 of file PostProcOnRefMesh.cpp.

654  {
656  if (gaussPtsTri.size1() == 0 && gaussPtsQuad.size1() == 0)
657  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
658  "post-process mesh not generated");
659 
660  auto create_tri = [&](auto &gauss_pts) {
661  EntityHandle tri;
662  std::array<EntityHandle, 3> tri_conn;
663  MatrixDouble3by3 coords_tri(3, 3);
664  VectorDouble3 coords(3);
665  CHKERR mField.get_moab().get_connectivity(
666  numeredEntFiniteElementPtr->getEnt(), conn, num_nodes, true);
667  CHKERR mField.get_moab().get_coords(conn, num_nodes, &coords_tri(0, 0));
668  for (int gg = 0; gg != 3; gg++) {
669  double ksi = gauss_pts(0, gg);
670  double eta = gauss_pts(1, gg);
671  double n0 = N_MBTRI0(ksi, eta);
672  double n1 = N_MBTRI1(ksi, eta);
673  double n2 = N_MBTRI2(ksi, eta);
674  double x =
675  n0 * coords_tri(0, 0) + n1 * coords_tri(1, 0) + n2 * coords_tri(2, 0);
676  double y =
677  n0 * coords_tri(0, 1) + n1 * coords_tri(1, 1) + n2 * coords_tri(2, 1);
678  double z =
679  n0 * coords_tri(0, 2) + n1 * coords_tri(1, 2) + n2 * coords_tri(2, 2);
680  coords[0] = x;
681  coords[1] = y;
682  coords[2] = z;
683  CHKERR postProcMesh.create_vertex(&coords[0], tri_conn[gg]);
684  }
685  CHKERR postProcMesh.create_element(MBTRI, &tri_conn[0], 3, tri);
686  return tri;
687  };
688 
689  auto create_quad = [&](auto &gauss_pts) {
690  EntityHandle quad;
691  std::array<EntityHandle, 4> quad_conn;
692  MatrixDouble coords_quad(4, 3);
693  VectorDouble coords(4);
694  CHKERR mField.get_moab().get_connectivity(
695  numeredEntFiniteElementPtr->getEnt(), conn, num_nodes, true);
696  CHKERR mField.get_moab().get_coords(conn, num_nodes, &coords_quad(0, 0));
697  for (int gg = 0; gg != 4; gg++) {
698  double ksi = gauss_pts(0, gg);
699  double eta = gauss_pts(1, gg);
700  double n0 = N_MBQUAD0(ksi, eta);
701  double n1 = N_MBQUAD1(ksi, eta);
702  double n2 = N_MBQUAD2(ksi, eta);
703  double n3 = N_MBQUAD3(ksi, eta);
704  double x = n0 * coords_quad(0, 0) + n1 * coords_quad(1, 0) +
705  n2 * coords_quad(2, 0) + n3 * coords_quad(3, 0);
706  double y = n0 * coords_quad(0, 1) + n1 * coords_quad(1, 1) +
707  n2 * coords_quad(2, 1) + n3 * coords_quad(3, 1);
708  double z = n0 * coords_quad(0, 2) + n1 * coords_quad(1, 2) +
709  n2 * coords_quad(2, 2) + n3 * coords_quad(3, 2);
710  coords[0] = x;
711  coords[1] = y;
712  coords[2] = z;
713  CHKERR postProcMesh.create_vertex(&coords[0], quad_conn[gg]);
714  }
715  CHKERR postProcMesh.create_element(MBQUAD, &quad_conn[0], 4, quad);
716  return quad;
717  };
718 
719  auto add_mid_nodes = [&](auto tri) {
721  Range edges;
722  CHKERR postProcMesh.get_adjacencies(&tri, 1, 1, true, edges);
723  EntityHandle meshset;
724  CHKERR postProcMesh.create_meshset(MESHSET_SET, meshset);
725  CHKERR postProcMesh.add_entities(meshset, &tri, 1);
726  CHKERR postProcMesh.add_entities(meshset, edges);
727  CHKERR postProcMesh.convert_entities(meshset, true, false, false);
728  CHKERR postProcMesh.delete_entities(&meshset, 1);
729  CHKERR postProcMesh.delete_entities(edges);
731  };
732 
733  EntityHandle tri;
734 
735  if (elementsMap.find(numeredEntFiniteElementPtr->getEnt()) !=
736  elementsMap.end()) {
737  tri = elementsMap[numeredEntFiniteElementPtr->getEnt()];
738  switch (numeredEntFiniteElementPtr->getEntType()) {
739  case MBTRI:
740  gaussPts.resize(gaussPtsTri.size1(), gaussPtsTri.size2(), false);
741  noalias(gaussPts) = gaussPtsTri;
742  break;
743  case MBQUAD:
744  gaussPts.resize(gaussPtsQuad.size1(), gaussPtsQuad.size2(), false);
745  noalias(gaussPts) = gaussPtsQuad;
746  break;
747  default:
748  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
749  "Not implemented element type");
750  }
751  } else {
752  switch (numeredEntFiniteElementPtr->getEntType()) {
753  case MBTRI:
754  gaussPts.resize(gaussPtsTri.size1(), gaussPtsTri.size2(), false);
755  noalias(gaussPts) = gaussPtsTri;
756  tri = create_tri(gaussPtsTri);
757  break;
758  case MBQUAD:
759  gaussPts.resize(gaussPtsQuad.size1(), gaussPtsQuad.size2(), false);
760  noalias(gaussPts) = gaussPtsQuad;
761  tri = create_quad(gaussPtsQuad);
762  break;
763  default:
764  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
765  "Not implemented element type");
766  }
768  CHKERR add_mid_nodes(tri);
769  elementsMap[numeredEntFiniteElementPtr->getEnt()] = tri;
770  }
771 
772  // Set values which map nodes with integration points on the prism
773  const EntityHandle *conn;
774  int num_nodes;
775 
776  CHKERR postProcMesh.get_connectivity(tri, conn, num_nodes, false);
777  mapGaussPts.resize(num_nodes);
778  for (int nn = 0; nn != num_nodes; nn++)
779  mapGaussPts[nn] = conn[nn];
780 
782 }
#define N_MBQUAD2(x, y)
quad shape function
Definition: fem_tools.h:67
moab::Interface & postProcMesh
std::map< EntityHandle, EntityHandle > elementsMap
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:74
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
MatrixDouble gaussPtsTri
#define N_MBQUAD3(x, y)
quad shape function
Definition: fem_tools.h:68
#define N_MBQUAD0(x, y)
quad shape function
Definition: fem_tools.h:65
MatrixDouble gaussPtsQuad
MatrixBoundedArray< double, 9 > MatrixDouble3by3
Definition: Types.hpp:97
bool sixNodePostProcTris
#define CHKERR
Inline error check.
Definition: definitions.h:596
std::vector< EntityHandle > mapGaussPts
#define N_MBTRI2(x, y)
triangle shape function
Definition: fem_tools.h:56
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:87
#define N_MBQUAD1(x, y)
quad shape function
Definition: fem_tools.h:66
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:72
#define N_MBTRI1(x, y)
triangle shape function
Definition: fem_tools.h:55
#define N_MBTRI0(x, y)
triangle shape function
Definition: fem_tools.h:54

Member Data Documentation

◆ commonData

CommonData PostProcFaceOnRefinedMesh::commonData

Definition at line 765 of file PostProcOnRefMesh.hpp.

◆ elementsMap

std::map<EntityHandle, EntityHandle> PostProcFaceOnRefinedMesh::elementsMap

Definition at line 759 of file PostProcOnRefMesh.hpp.

◆ gaussPtsQuad

MatrixDouble PostProcFaceOnRefinedMesh::gaussPtsQuad
private

Definition at line 839 of file PostProcOnRefMesh.hpp.

◆ gaussPtsTri

MatrixDouble PostProcFaceOnRefinedMesh::gaussPtsTri
private

Definition at line 838 of file PostProcOnRefMesh.hpp.

◆ sixNodePostProcTris

bool PostProcFaceOnRefinedMesh::sixNodePostProcTris

Definition at line 737 of file PostProcOnRefMesh.hpp.


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