v0.14.0
read_vtk.cpp
Go to the documentation of this file.
1 /** \file read_vtk.cpp
2 
3  \brief Tool to read mesh files as vtk files assuming Boundary conditions
4  are set with a field named "BOUNDARY_CONDITIONS"
5 
6 */
7 
8 #include <MoFEM.hpp>
9 using namespace MoFEM;
10 
11 struct VtkInterface {
12 
14  : mField(m_field), mOab(moab) {}
15 
16  MoFEMErrorCode readVtk();
17 
18 private:
22 
23  MoFEMErrorCode readMesh();
24  MoFEMErrorCode getBoundaryConditions();
25  MoFEMErrorCode setBoundaryConditions();
26  MoFEMErrorCode getMaterialProperties();
27  MoFEMErrorCode writeOutput();
28  MoFEMErrorCode checkResults();
29 };
30 
33  CHKERR readMesh();
34  CHKERR getBoundaryConditions();
35  CHKERR setBoundaryConditions();
36  CHKERR getMaterialProperties();
37  CHKERR writeOutput();
38  CHKERR checkResults();
40 }
41 
44 
45  CHKERR mField.getInterface(simpleInterface);
46  CHKERR simpleInterface->getOptions();
47  CHKERR simpleInterface->loadFile();
48 
50 }
51 
54  MeshsetsManager *meshsets_manager_ptr;
55  CHKERR mField.getInterface(meshsets_manager_ptr);
56 
57  Tag bc_tag;
58  CHKERR mOab.tag_get_handle("BOUNDARY_CONDITIONS", 1, MB_TYPE_INTEGER, bc_tag,
59  MB_TAG_CREAT | MB_TAG_SPARSE);
60 
61  std::vector<std::tuple<int, std::string, std::string>> conditions = {
62  {1, "FIX_ALL", "FIX_ALL"},
63  {2, "FIX_X", "FIX_X"},
64  {3, "FIX_Y", "FIX_Y"},
65  {4, "FIX_Z", "FIX_Z"},
66  {5, "FIX_X_2", "DISP_X"},
67  {6, "FIX_Y_2", "DISP_Y"},
68  {7, "FIX_Z_2", "DISP_Z"},
69  {8, "FORCE_X", "FORCE_X"},
70  {9, "FORCE_Y", "FORCE_Y"},
71  {10, "FORCE_Z", "FORCE_Z"},
72  {11, "VEL_X", "VEL_X"},
73  {12, "VEL_Y", "VEL_Y"},
74  {13, "VEL_Z", "VEL_Z"},
75  {14, "ACCL_X", "ACCL_X"},
76  {15, "ACCL_Y", "ACCL_Y"},
77  {16, "ACCL_Z", "ACCL_Z"},
78  {17, "TEMP", "TEMP"},
79  {18, "PRESSURE", "PRESSURE"},
80  {19, "HEAT_FLUX", "HEAT_FLUX"},
81  {20, "CONTACT", "CONTACT"},
82  {50, "TIE_MATRIX", "TIE_MATRIX"}};
83 
84  for (auto &condition : conditions) {
85  int desired_val = std::get<0>(condition);
86  const void *desired_val_ptr = &desired_val;
87  Range bc;
88  CHKERR mOab.get_entities_by_type_and_tag(0, MBVERTEX, &bc_tag,
89  &desired_val_ptr, 1, bc);
90 
91  if (!bc.empty()) {
92  Range faces;
93  CHKERR mOab.get_adjacencies(bc, 2, true, faces, moab::Interface::UNION);
94  Range edges;
95  CHKERR mOab.get_adjacencies(bc, 1, true, edges, moab::Interface::UNION);
96 
97  Range adj_nodes;
98  CHKERR mOab.get_adjacencies(faces, 0, false, adj_nodes,
99  moab::Interface::UNION);
100 
101  Range non_BC_nodes;
102  non_BC_nodes = subtract(adj_nodes, bc);
103  Range non_BC_faces;
104  CHKERR mOab.get_adjacencies(non_BC_nodes, 2, false, non_BC_faces,
105  moab::Interface::UNION);
106  Range non_BC_edges;
107  CHKERR mOab.get_adjacencies(non_BC_nodes, 1, false, non_BC_edges,
108  moab::Interface::UNION);
109 
110  faces = subtract(faces, non_BC_faces);
111  edges = subtract(edges, non_BC_edges);
112 
113  MOFEM_LOG_C("WORLD", Sev::inform, "Number of nodes assigned to %s = %d",
114  std::get<2>(condition).c_str(), bc.size());
115  MOFEM_LOG_C("WORLD", Sev::inform, "Number of edges assigned to %s = %d",
116  std::get<2>(condition).c_str(), edges.size());
117  MOFEM_LOG_C("WORLD", Sev::inform, "Number of faces assigned to %s = %d",
118  std::get<2>(condition).c_str(), faces.size());
119 
120  CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, desired_val,
121  std::get<1>(condition));
122  CHKERR meshsets_manager_ptr->addEntitiesToMeshset(BLOCKSET, desired_val,
123  bc);
124  CHKERR meshsets_manager_ptr->addEntitiesToMeshset(BLOCKSET, desired_val,
125  edges);
126  CHKERR meshsets_manager_ptr->addEntitiesToMeshset(BLOCKSET, desired_val,
127  faces);
128  }
129  }
130 
132 }
133 
136  // Add meshsets if config file provided
137  MeshsetsManager *meshsets_interface_ptr;
138  CHKERR mField.getInterface(meshsets_interface_ptr);
139  CHKERR meshsets_interface_ptr->setMeshsetFromFile();
140 
141  MOFEM_LOG_CHANNEL("WORLD");
142  MOFEM_LOG_TAG("WORLD", "read_vtk")
143  MOFEM_LOG("WORLD", Sev::inform)
144  << "Print all meshsets (old and added from meshsets "
145  "configurational file)";
146  for (auto cit = meshsets_interface_ptr->getBegin();
147  cit != meshsets_interface_ptr->getEnd(); cit++)
148  MOFEM_LOG("WORLD", Sev::inform) << *cit;
149 
151 }
152 
155  MeshsetsManager *meshsets_manager_ptr;
156  CHKERR mField.getInterface(meshsets_manager_ptr);
157 
158  Range ents;
159  rval = mOab.get_entities_by_type(0, MBHEX, ents, true);
160  if (rval != MB_SUCCESS) {
161  MOFEM_LOG("WORLD", Sev::warning)
162  << "No hexes in the mesh, no material block set. Not Implemented";
163  } else {
164  CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 100, "ADOLCMAT");
165  CHKERR meshsets_manager_ptr->addEntitiesToMeshset(BLOCKSET, 100, ents);
166  MOFEM_LOG("WORLD", Sev::inform) << "Material block ADOLCMAT set added";
167  MOFEM_LOG("WORLD", Sev::warning)
168  << "Other material block sets not implemented";
169  }
171 }
172 
175 
176  char mesh_out_file[255] = "out_vtk.h5m";
177 
178  CHKERR PetscOptionsBegin(mField.get_comm(), "", "Read vtk tool", "none");
179  CHKERR PetscOptionsString("-output_file", "output mesh file name", "",
180  "out.h5m", mesh_out_file, 255, PETSC_NULL);
181  ierr = PetscOptionsEnd();
182  CHKERRQ(ierr);
183 
184  CHKERR mOab.write_file(mesh_out_file);
185 
187 }
188 
191  PetscBool atom_test = PETSC_FALSE;
192  PetscOptionsGetBool(PETSC_NULL, "", "-atom_test", &atom_test, PETSC_NULL);
193  if (atom_test) {
194  // check mesh has been read correctly
195  Range ents;
196  CHKERR mOab.get_entities_by_dimension(0, 3, ents);
197  if (ents.size() != 125) {
198  SETERRQ3(
199  PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
200  "atom test %d failed: Wrong number of elements %d, should be 125",
201  atom_test, ents.size(), 125);
202  }
203  // check boundary conditions have been set correctly
204  MeshsetsManager *meshsets_manager_ptr;
205  CHKERR mField.getInterface(meshsets_manager_ptr);
206  // iterate meshsets
208  Range ents;
209  CHKERR mField.get_moab().get_entities_by_dimension(it->meshset, 2, ents,
210  true);
211  if (it->getMeshsetId() == 1) {
212  if (it->getName() != "FIX_ALL") {
213  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
214  "atom test %d failed: Wrong name of meshset %d, should be "
215  "FIX_ALL is %s",
216  atom_test, it->getMeshsetId(), it->getName());
217  }
218  if (ents.size() != 25) {
219  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
220  "atom test %d failed: Wrong number of entities in meshset "
221  "%d with %d, should be 25",
222  atom_test, it->getMeshsetId(), ents.size());
223  }
224  } else if (it->getMeshsetId() == 2) {
225  if (it->getName() != "FIX_X") {
226  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
227  "atom test %d failed: Wrong name of meshset %d, should be "
228  "FIX_X is %s",
229  atom_test, it->getMeshsetId(), it->getName());
230  }
231  if (ents.size() != 25) {
232  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
233  "atom test %d failed: Wrong number of entities in meshset "
234  "%d with %d, should be 25",
235  atom_test, it->getMeshsetId(), ents.size());
236  }
237  } else if (it->getMeshsetId() == 100) {
238  if (it->getName() != "ADOLCMAT") {
239  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
240  "atom test %d failed: Wrong name of meshset %d, should be "
241  "ADOLCMAT is %s",
242  atom_test, it->getMeshsetId(), it->getName());
243  }
244  }
245  }
246  }
248 }
249 
250 static char help[] = "...\n\n";
251 
252 int main(int argc, char *argv[]) {
253  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
254 
255  try {
256 
257  moab::Core mb_instance;
258  moab::Interface &moab = mb_instance;
259  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
260  if (pcomm == NULL)
261  pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
262 
263  // Create MoFEM database
264  MoFEM::Core core(moab);
265  MoFEM::Interface &m_field = core;
266 
267  VtkInterface read_vtk(m_field, moab);
268  CHKERR read_vtk.readVtk();
269  }
270  CATCH_ERRORS;
271 
273 
274  return 0;
275 }
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
VtkInterface::getBoundaryConditions
MoFEMErrorCode getBoundaryConditions()
Definition: read_vtk.cpp:52
MoFEM::MeshsetsManager::addEntitiesToMeshset
MoFEMErrorCode addEntitiesToMeshset(const CubitBCType cubit_bc_type, const int ms_id, const Range &ents)
add entities to cubit meshset
Definition: MeshsetsManager.cpp:418
VtkInterface::writeOutput
MoFEMErrorCode writeOutput()
Definition: read_vtk.cpp:173
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
VtkInterface::mOab
moab::Interface & mOab
Definition: read_vtk.cpp:20
VtkInterface::readVtk
MoFEMErrorCode readVtk()
Definition: read_vtk.cpp:31
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM.hpp
VtkInterface::VtkInterface
VtkInterface(MoFEM::Interface &m_field, moab::Interface &moab)
Definition: read_vtk.cpp:13
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
VtkInterface::checkResults
MoFEMErrorCode checkResults()
Definition: read_vtk.cpp:189
VtkInterface::getMaterialProperties
MoFEMErrorCode getMaterialProperties()
Definition: read_vtk.cpp:153
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
atom_test
int atom_test
Definition: contact.cpp:97
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::MeshsetsManager::getEnd
CubitMeshSet_multiIndex::iterator getEnd() const
get begin iterator of cubit mehset of given type (instead you can use IT_CUBITMESHSETS_TYPE_FOR_LOOP(...
Definition: MeshsetsManager.hpp:257
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::MeshsetsManager::getBegin
CubitMeshSet_multiIndex::iterator getBegin() const
get begin iterator of cubit mehset of given type (instead you can use IT_CUBITMESHSETS_TYPE_FOR_LOOP(...
Definition: MeshsetsManager.hpp:243
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
MoFEM::MeshsetsManager::setMeshsetFromFile
MoFEMErrorCode setMeshsetFromFile(const string file_name, const bool clean_file_options=true)
add blocksets reading config file
Definition: MeshsetsManager.cpp:791
MoFEM::MeshsetsManager::addMeshset
MoFEMErrorCode addMeshset(const CubitBCType cubit_bc_type, const int ms_id, const std::string name="")
add cubit meshset
Definition: MeshsetsManager.cpp:385
main
int main(int argc, char *argv[])
Definition: read_vtk.cpp:252
VtkInterface::mField
MoFEM::Interface & mField
Definition: read_vtk.cpp:19
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
VtkInterface::readMesh
MoFEMErrorCode readMesh()
Definition: read_vtk.cpp:42
VtkInterface::simpleInterface
Simple * simpleInterface
Definition: read_vtk.cpp:21
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
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
_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
VtkInterface
Definition: read_vtk.cpp:11
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
help
static char help[]
Definition: read_vtk.cpp:250
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
VtkInterface::setBoundaryConditions
MoFEMErrorCode setBoundaryConditions()
Definition: read_vtk.cpp:134
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182