v0.14.0
extrude_prisms.cpp
Go to the documentation of this file.
1 /** \file extrude_prisms.cpp
2  * \brief Extrude prisms from surface elements block
3  * \example mesh_cut.cpp
4  *
5  * \ingroup mesh_cut
6  */
7 
8 
9 
10 #include <MoFEM.hpp>
11 
12 using namespace MoFEM;
13 
14 static char help[] = "mesh cutting\n\n";
15 
16 int main(int argc, char *argv[]) {
17 
18  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
19 
20  try {
21 
22  PetscBool flg_myfile = PETSC_TRUE;
23  char mesh_file_name[255];
24  char mesh_out_file[255] = "out.h5m";
25 
26  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "Prism surface", "none");
27  CHKERR PetscOptionsString("-my_file", "mesh file name", "", "mesh.h5m",
28  mesh_file_name, 255, &flg_myfile);
29  CHKERR PetscOptionsString("-output_file", "output mesh file name", "",
30  "out.h5m", mesh_out_file, 255, PETSC_NULL);
31  ierr = PetscOptionsEnd();
32  CHKERRG(ierr);
33 
34  if (flg_myfile != PETSC_TRUE)
35  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
36  "error -my_file (mesh file needed)");
37 
38  moab::Core mb_instance;
39  moab::Interface &moab = mb_instance;
40  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
41  if (pcomm == NULL)
42  pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
43 
44  const char *option;
45  option = "";
46  CHKERR moab.load_file(mesh_file_name, 0, option);
47 
48  MoFEM::Core core(moab);
49  MoFEM::CoreInterface &m_field =
51 
52  PrismsFromSurfaceInterface *prisms_from_surface_interface;
53  CHKERR m_field.getInterface(prisms_from_surface_interface);
54 
55  const std::string extrude_block_name = "EXTRUDE_PRISMS";
57  if (it->getName().compare(0, extrude_block_name.length(),
58  extrude_block_name) == 0) {
59  std::vector<double> thickness;
60  CHKERR it->getAttributes(thickness);
61  if (thickness.size() != 2)
62  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
63  "Data inconsistency");
64  Range tris;
65  CHKERR m_field.get_moab().get_entities_by_type(it->meshset, MBTRI, tris,
66  true);
67  Range block_prisms;
68  CHKERR prisms_from_surface_interface->createPrisms(
69  tris, PrismsFromSurfaceInterface::NO_SWAP, block_prisms);
70  CHKERR prisms_from_surface_interface->setNormalThickness(
71  block_prisms, thickness[0], thickness[1]);
72  CHKERR prisms_from_surface_interface->updateMeshestByEdgeBlock(
73  block_prisms);
74  CHKERR prisms_from_surface_interface->updateMeshestByTriBlock(
75  block_prisms);
76 
77  std::cout << "Extrude block " << it->getMeshsetId() << " set prisms "
78  << block_prisms.size() << endl;
79  }
80  }
81 
82  CHKERR moab.write_file(mesh_out_file);
83  }
85 
87 
88  return 0;
89 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::PrismsFromSurfaceInterface
merge node from two bit levels
Definition: PrismsFromSurfaceInterface.hpp:15
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
flg_myfile
PetscBool flg_myfile
Definition: mesh_smoothing.cpp:22
MoFEM.hpp
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::PrismsFromSurfaceInterface::updateMeshestByTriBlock
MoFEMErrorCode updateMeshestByTriBlock(const Range &prisms)
Add prism to bockset.
Definition: PrismsFromSurfaceInterface.cpp:287
help
static char help[]
Definition: extrude_prisms.cpp:14
main
int main(int argc, char *argv[])
Definition: extrude_prisms.cpp:16
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
Range
MoFEM::CoreTmp< 0 >::Initialize
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
Definition: MeshsetsManager.hpp:71
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
BLOCKSET
@ BLOCKSET
Definition: definitions.h:148
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEM::PrismsFromSurfaceInterface::NO_SWAP
@ NO_SWAP
Definition: PrismsFromSurfaceInterface.hpp:38
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::PrismsFromSurfaceInterface::updateMeshestByEdgeBlock
MoFEMErrorCode updateMeshestByEdgeBlock(const Range &prisms)
Add quads to bockset.
Definition: PrismsFromSurfaceInterface.cpp:260
MoFEM::PrismsFromSurfaceInterface::createPrisms
MoFEMErrorCode createPrisms(const Range &ents, const SwapType swap_type, Range &prisms, int verb=-1)
Make prisms from triangles.
Definition: PrismsFromSurfaceInterface.cpp:22
MoFEM::CoreInterface
Interface.
Definition: Interface.hpp:27
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
MoFEM::PrismsFromSurfaceInterface::setNormalThickness
MoFEMErrorCode setNormalThickness(const Range &prisms, double thickness3, double thickness4)
Definition: PrismsFromSurfaceInterface.cpp:195
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36