v0.14.0
Public Member Functions | Private Member Functions | Private Attributes | List of all members
VtkInterface Struct Reference
Collaboration diagram for VtkInterface:
[legend]

Public Member Functions

 VtkInterface (MoFEM::Interface &m_field, moab::Interface &moab)
 
MoFEMErrorCode readVtk ()
 

Private Member Functions

MoFEMErrorCode readMesh ()
 
MoFEMErrorCode getBoundaryConditions ()
 
MoFEMErrorCode setBoundaryConditions ()
 
MoFEMErrorCode getMaterialProperties ()
 
MoFEMErrorCode getBlockSetBCFromMesh (std::string)
 
MoFEMErrorCode writeOutput ()
 
MoFEMErrorCode checkResults ()
 

Private Attributes

MoFEM::InterfacemField
 
moab::Interface & mOab
 
SimplesimpleInterface
 

Detailed Description

Definition at line 11 of file read_vtk.cpp.

Constructor & Destructor Documentation

◆ VtkInterface()

VtkInterface::VtkInterface ( MoFEM::Interface m_field,
moab::Interface &  moab 
)
inline

Definition at line 13 of file read_vtk.cpp.

14  : mField(m_field), mOab(moab) {}

Member Function Documentation

◆ checkResults()

MoFEMErrorCode VtkInterface::checkResults ( )
private

Definition at line 237 of file read_vtk.cpp.

237  {
239  PetscBool atom_test = PETSC_FALSE;
240  PetscOptionsGetBool(PETSC_NULL, "", "-atom_test", &atom_test, PETSC_NULL);
241  if (atom_test) {
242  // check mesh has been read correctly
243  Range ents;
244  CHKERR mOab.get_entities_by_dimension(0, 3, ents);
245  if (ents.size() != 125) {
246  SETERRQ3(
247  PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
248  "atom test %d failed: Wrong number of elements %d, should be 125",
249  atom_test, ents.size(), 125);
250  }
251  // check boundary conditions have been set correctly
252  MeshsetsManager *meshsets_manager_ptr;
253  CHKERR mField.getInterface(meshsets_manager_ptr);
254  // iterate meshsets
256  Range ents;
257  CHKERR mField.get_moab().get_entities_by_dimension(it->meshset, 2, ents,
258  true);
259  if (it->getMeshsetId() == 1) {
260  if (it->getName() != "FIX_ALL") {
261  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
262  "atom test %d failed: Wrong name of meshset %d, should be "
263  "FIX_ALL is %s",
264  atom_test, it->getMeshsetId(), it->getName());
265  }
266  if (ents.size() != 25) {
267  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
268  "atom test %d failed: Wrong number of entities in meshset "
269  "%d with %d, should be 25",
270  atom_test, it->getMeshsetId(), ents.size());
271  }
272  } else if (it->getMeshsetId() == 2) {
273  if (it->getName() != "FIX_X") {
274  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
275  "atom test %d failed: Wrong name of meshset %d, should be "
276  "FIX_X is %s",
277  atom_test, it->getMeshsetId(), it->getName());
278  }
279  if (ents.size() != 25) {
280  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
281  "atom test %d failed: Wrong number of entities in meshset "
282  "%d with %d, should be 25",
283  atom_test, it->getMeshsetId(), ents.size());
284  }
285  } else if (it->getMeshsetId() == 100) {
286  if (it->getName() != "ADOLCMAT") {
287  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
288  "atom test %d failed: Wrong name of meshset %d, should be "
289  "ADOLCMAT is %s",
290  atom_test, it->getMeshsetId(), it->getName());
291  }
292  }
293  }
294  }
296 }

◆ getBlockSetBCFromMesh()

MoFEMErrorCode VtkInterface::getBlockSetBCFromMesh ( std::string  block_name)
private

Definition at line 53 of file read_vtk.cpp.

53  {
55 
56  MeshsetsManager *meshsets_mng;
57  CHKERR mField.getInterface(meshsets_mng);
58 
59  Tag mat_tag;
60  rval = mOab.tag_get_handle(block_name.c_str(), 1, MB_TYPE_INTEGER, mat_tag,
61  MB_TAG_EXCL | MB_TAG_SPARSE);
62  if (rval != MB_ALREADY_ALLOCATED) {
63  // std::cerr << "Error creating tag " << block_name << std::endl;
65  }
66  Range block_ents;
67  for (auto mbtype : {MBTET, MBTRI, MBHEX, MBQUAD, MBEDGE, MBVERTEX}) {
68  CHKERR mOab.get_entities_by_type_and_tag(0, mbtype, &mat_tag, NULL, 1,
69  block_ents, moab::Interface::UNION);
70  }
71 
72  if (!block_ents.empty()) {
73  std::map<int, Range> mat_map;
74  for (auto &ent : block_ents) {
75  int mat_id;
76  CHKERR mOab.tag_get_data(mat_tag, &ent, 1, &mat_id);
77  mat_map[mat_id].insert(ent);
78  }
79 
80  for (auto &[mat_id, mat_ents] : mat_map) {
81  std::string meshset_name = block_name + "_" + std::to_string(mat_id);
82  CHKERR meshsets_mng->addMeshset(BLOCKSET, mat_id, meshset_name);
83  CHKERR meshsets_mng->addEntitiesToMeshset(BLOCKSET, mat_id, mat_ents);
84  MOFEM_LOG("WORLD", Sev::inform)
85  << "Blockset " << meshset_name << " set added.";
86  }
87  }
88 
90 }

◆ getBoundaryConditions()

MoFEMErrorCode VtkInterface::getBoundaryConditions ( )
private

Definition at line 92 of file read_vtk.cpp.

92  {
94  MeshsetsManager *meshsets_manager_ptr;
95  CHKERR mField.getInterface(meshsets_manager_ptr);
96 
97  Tag bc_tag;
98  CHKERR mOab.tag_get_handle("BOUNDARY_CONDITIONS", 1, MB_TYPE_INTEGER, bc_tag,
99  MB_TAG_CREAT | MB_TAG_SPARSE);
100 
101  std::vector<std::tuple<int, std::string, std::string>> conditions = {
102  {1, "FIX_ALL", "FIX_ALL"},
103  {2, "FIX_X", "FIX_X"},
104  {3, "FIX_Y", "FIX_Y"},
105  {4, "FIX_Z", "FIX_Z"},
106  {5, "FIX_X_2", "DISP_X"},
107  {6, "FIX_Y_2", "DISP_Y"},
108  {7, "FIX_Z_2", "DISP_Z"},
109  {8, "FORCE_X", "FORCE_X"},
110  {9, "FORCE_Y", "FORCE_Y"},
111  {10, "FORCE_Z", "FORCE_Z"},
112  {11, "VEL_X", "VEL_X"},
113  {12, "VEL_Y", "VEL_Y"},
114  {13, "VEL_Z", "VEL_Z"},
115  {14, "ACCL_X", "ACCL_X"},
116  {15, "ACCL_Y", "ACCL_Y"},
117  {16, "ACCL_Z", "ACCL_Z"},
118  {17, "TEMP", "TEMP"},
119  {18, "PRESSURE", "PRESSURE"},
120  {19, "HEAT_FLUX", "HEAT_FLUX"},
121  {20, "CONTACT", "CONTACT"},
122  {21, "DISPLACEMENT", "DISPLACEMENT"},
123  {22, "ROTATE_ALL", "ROTATE_ALL"},
124  {23, "ROTATE_X", "ROTATE_X"},
125  {24, "ROTATE_Y", "ROTATE_Y"},
126  {25, "ROTATE_Z", "ROTATE_Z"},
127  {26, "TEMPERATURE", "TEMPERATURE"},
128  {50, "TIE_MATRIX", "TIE_MATRIX"}};
129 
130  for (auto &condition : conditions) {
131  int desired_val = std::get<0>(condition);
132  const void *desired_val_ptr = &desired_val;
133  Range bc;
134  CHKERR mOab.get_entities_by_type_and_tag(0, MBVERTEX, &bc_tag,
135  &desired_val_ptr, 1, bc);
136 
137  if (!bc.empty()) {
138  Range faces;
139  CHKERR mOab.get_adjacencies(bc, 2, true, faces, moab::Interface::UNION);
140  Range edges;
141  CHKERR mOab.get_adjacencies(bc, 1, true, edges, moab::Interface::UNION);
142 
143  Range adj_nodes;
144  CHKERR mOab.get_adjacencies(faces, 0, false, adj_nodes,
145  moab::Interface::UNION);
146 
147  Range non_BC_nodes;
148  non_BC_nodes = subtract(adj_nodes, bc);
149  Range non_BC_faces;
150  CHKERR mOab.get_adjacencies(non_BC_nodes, 2, false, non_BC_faces,
151  moab::Interface::UNION);
152  Range non_BC_edges;
153  CHKERR mOab.get_adjacencies(non_BC_nodes, 1, false, non_BC_edges,
154  moab::Interface::UNION);
155 
156  faces = subtract(faces, non_BC_faces);
157  edges = subtract(edges, non_BC_edges);
158 
159  MOFEM_LOG_C("WORLD", Sev::inform, "Number of nodes assigned to %s = %d",
160  std::get<2>(condition).c_str(), bc.size());
161  MOFEM_LOG_C("WORLD", Sev::inform, "Number of edges assigned to %s = %d",
162  std::get<2>(condition).c_str(), edges.size());
163  MOFEM_LOG_C("WORLD", Sev::inform, "Number of faces assigned to %s = %d",
164  std::get<2>(condition).c_str(), faces.size());
165 
166  CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, desired_val,
167  std::get<1>(condition));
168  CHKERR meshsets_manager_ptr->addEntitiesToMeshset(BLOCKSET, desired_val,
169  bc);
170  CHKERR meshsets_manager_ptr->addEntitiesToMeshset(BLOCKSET, desired_val,
171  edges);
172  CHKERR meshsets_manager_ptr->addEntitiesToMeshset(BLOCKSET, desired_val,
173  faces);
174  } else {
175  CHKERR getBlockSetBCFromMesh(std::get<1>(condition));
176  }
177  }
178 
180 }

◆ getMaterialProperties()

MoFEMErrorCode VtkInterface::getMaterialProperties ( )
private

Definition at line 201 of file read_vtk.cpp.

201  {
203  MeshsetsManager *meshsets_manager_ptr;
204  CHKERR mField.getInterface(meshsets_manager_ptr);
205 
206  Range ents;
207  rval = mOab.get_entities_by_dimension(0, 3, ents, true);
208  if (rval != MB_SUCCESS) {
209  MOFEM_LOG("WORLD", Sev::warning)
210  << "No hexes/tets in the mesh, no material block set. Not Implemented";
211  } else {
212  CHKERR getBlockSetBCFromMesh("MAT_ELASTIC");
213  // for backward compatibility
214  CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 100, "ADOLCMAT");
215  CHKERR meshsets_manager_ptr->addEntitiesToMeshset(BLOCKSET, 100, ents);
216  MOFEM_LOG("WORLD", Sev::inform) << "Material block ADOLCMAT set added.";
217  }
219 }

◆ readMesh()

MoFEMErrorCode VtkInterface::readMesh ( )
private

Definition at line 43 of file read_vtk.cpp.

◆ readVtk()

MoFEMErrorCode VtkInterface::readVtk ( )

Definition at line 32 of file read_vtk.cpp.

◆ setBoundaryConditions()

MoFEMErrorCode VtkInterface::setBoundaryConditions ( )
private

Definition at line 182 of file read_vtk.cpp.

182  {
184  // Add meshsets if config file provided
185  MeshsetsManager *meshsets_interface_ptr;
186  CHKERR mField.getInterface(meshsets_interface_ptr);
187  CHKERR meshsets_interface_ptr->setMeshsetFromFile();
188 
189  MOFEM_LOG_CHANNEL("WORLD");
190  MOFEM_LOG_TAG("WORLD", "read_vtk")
191  MOFEM_LOG("WORLD", Sev::inform)
192  << "Print all meshsets (old and added from meshsets "
193  "configurational file)";
194  for (auto cit = meshsets_interface_ptr->getBegin();
195  cit != meshsets_interface_ptr->getEnd(); cit++)
196  MOFEM_LOG("WORLD", Sev::inform) << *cit;
197 
199 }

◆ writeOutput()

MoFEMErrorCode VtkInterface::writeOutput ( )
private

Definition at line 221 of file read_vtk.cpp.

221  {
223 
224  char mesh_out_file[255] = "out_vtk.h5m";
225 
226  CHKERR PetscOptionsBegin(mField.get_comm(), "", "Read vtk tool", "none");
227  CHKERR PetscOptionsString("-output_file", "output mesh file name", "",
228  "out.h5m", mesh_out_file, 255, PETSC_NULL);
229  ierr = PetscOptionsEnd();
230  CHKERRQ(ierr);
231 
232  CHKERR mOab.write_file(mesh_out_file);
233 
235 }

Member Data Documentation

◆ mField

MoFEM::Interface& VtkInterface::mField
private

Definition at line 19 of file read_vtk.cpp.

◆ mOab

moab::Interface& VtkInterface::mOab
private

Definition at line 20 of file read_vtk.cpp.

◆ simpleInterface

Simple* VtkInterface::simpleInterface
private

Definition at line 21 of file read_vtk.cpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
VtkInterface::getBlockSetBCFromMesh
MoFEMErrorCode getBlockSetBCFromMesh(std::string)
Definition: read_vtk.cpp:53
VtkInterface::getBoundaryConditions
MoFEMErrorCode getBoundaryConditions()
Definition: read_vtk.cpp:92
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:221
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
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
MoFEM::Simple::loadFile
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
VtkInterface::checkResults
MoFEMErrorCode checkResults()
Definition: read_vtk.cpp:237
VtkInterface::getMaterialProperties
MoFEMErrorCode getMaterialProperties()
Definition: read_vtk.cpp:201
MoFEM::Simple::getOptions
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
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::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
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
VtkInterface::mField
MoFEM::Interface & mField
Definition: read_vtk.cpp:19
MoFEM::LogManager::SeverityLevel
SeverityLevel
Severity levels.
Definition: LogManager.hpp:33
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
VtkInterface::readMesh
MoFEMErrorCode readMesh()
Definition: read_vtk.cpp:43
VtkInterface::simpleInterface
Simple * simpleInterface
Definition: read_vtk.cpp:21
Range
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
_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
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
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
VtkInterface::setBoundaryConditions
MoFEMErrorCode setBoundaryConditions()
Definition: read_vtk.cpp:182
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