Chi-Tech
chi_unpartitioned_mesh.h
Go to the documentation of this file.
1#ifndef CHI_MESH_UNPARTITIONED_MESH_H
2#define CHI_MESH_UNPARTITIONED_MESH_H
3
4#include "mesh/chi_mesh.h"
5#include "mesh/Cell/cell.h"
6
7class vtkCell;
8class vtkUnstructuredGrid;
9template <class T>
10class vtkSmartPointer;
11
12#include <map>
13#include <array>
14
15// ###################################################################
16/**This object is intented for unpartitioned meshes that still require
17 * partitioning.*/
19{
20public:
22 {
23 std::vector<uint64_t> vertex_ids;
24 bool has_neighbor = false;
25 uint64_t neighbor = 0;
26
27 LightWeightFace() = default;
28 explicit LightWeightFace(std::vector<uint64_t> in_vertex_ids)
29 : vertex_ids(std::move(in_vertex_ids))
30 {
31 }
32 };
34 {
38 int material_id = -1;
39 std::vector<uint64_t> vertex_ids;
40 std::vector<LightWeightFace> faces;
41
43 chi_mesh::CellType in_sub_type)
44 : type(in_type), sub_type(in_sub_type)
45 {
46 }
47 };
48
49 struct Options
50 {
51 std::string file_name;
52 std::string material_id_fieldname = "BlockID";
54 double scale = 1.0;
55 size_t ortho_Nx = 0;
56 size_t ortho_Ny = 0;
57 size_t ortho_Nz = 0;
58
59 std::map<uint64_t, std::string> boundary_id_map;
60 };
61
62 struct BoundBox
63 {
64 double xmin = 0.0, xmax = 0.0, ymin = 0.0, ymax = 0.0, zmin = 0.0,
65 zmax = 0.0;
66 };
67
68protected:
69 std::vector<chi_mesh::Vertex> vertices_;
70 std::vector<LightWeightCell*> raw_cells_;
71 std::vector<LightWeightCell*> raw_boundary_cells_;
72 std::vector<std::set<uint64_t>> vertex_cell_subscriptions_;
73
76 std::shared_ptr<BoundBox> bound_box_ = nullptr;
77
78protected:
79 static LightWeightCell* CreateCellFromVTKPolyhedron(vtkCell* vtk_cell);
80 static LightWeightCell* CreateCellFromVTKPolygon(vtkCell* vtk_cell);
81 static LightWeightCell* CreateCellFromVTKLine(vtkCell* vtk_cell);
82 static LightWeightCell* CreateCellFromVTKVertex(vtkCell* vtk_cell);
83
85 typedef std::pair<vtkUGridPtr, std::string> vtkUGridPtrAndName;
86 void CopyUGridCellsAndPoints(vtkUnstructuredGrid& ugrid,
87 double scale,
88 int dimension_to_copy);
89
90 void SetMaterialIDsFromList(const std::vector<int>& material_ids);
91
92 void
93 SetBoundaryIDsFromBlocks(std::vector<vtkUGridPtrAndName>& bndry_grid_blocks);
94
95public:
96 const BoundBox& GetBoundBox() const { return *bound_box_; }
97
99 const Options& GetMeshOptions() const { return mesh_options_; }
100
103
104 const std::vector<std::set<uint64_t>>& GetVertextCellSubscriptions() const
105 {
107 }
108
109 void AddCell(LightWeightCell*& cell) { raw_cells_.push_back(cell); }
110 size_t GetNumberOfCells() const { return raw_cells_.size(); }
111 std::vector<LightWeightCell*>& GetRawCells() { return raw_cells_; }
112 const std::vector<LightWeightCell*>& GetRawCells() const
113 {
114 return raw_cells_;
115 }
116
117 const std::vector<chi_mesh::Vertex>& GetVertices() const { return vertices_; }
118 std::vector<chi_mesh::Vertex>& GetVertices() { return vertices_; }
119
122 /**Makes or gets a boundary that uniquely identifies the given name.*/
123 uint64_t MakeBoundaryID(const std::string& boundary_name);
124 void SetAttributes(MeshAttributes new_attribs,
125 std::array<size_t, 3> ortho_Nis = {0, 0, 0});
126
127 void ReadFromVTU(const Options& options);
128 void ReadFromPVTU(const Options& options);
129 void ReadFromEnsightGold(const Options& options);
130 void ReadFromWavefrontOBJ(const Options& options);
131
132 void ReadFromMsh(const Options& options);
133
134 void ReadFromExodus(const Options& options);
135
136 void PushProxyCell(const std::string& type_str,
137 const std::string& sub_type_str,
138 int cell_num_faces,
139 int cell_material_id,
140 const std::vector<std::vector<uint64_t>>& proxy_faces);
141
143
144 void CleanUp()
145 {
146 for (auto& cell : raw_cells_)
147 delete cell;
148 for (auto& cell : raw_boundary_cells_)
149 delete cell;
150 vertices_.clear();
151 vertices_.shrink_to_fit();
152 raw_cells_.clear();
153 raw_cells_.shrink_to_fit();
154 raw_boundary_cells_.clear();
155 raw_boundary_cells_.shrink_to_fit();
157 vertex_cell_subscriptions_.shrink_to_fit();
158 }
159};
160
161#endif // CHI_MESH_UNPARTITIONED_MESH_H
const MeshAttributes & GetMeshAttributes() const
std::pair< vtkUGridPtr, std::string > vtkUGridPtrAndName
static LightWeightCell * CreateCellFromVTKLine(vtkCell *vtk_cell)
const BoundBox & GetBoundBox() const
void SetAttributes(MeshAttributes new_attribs, std::array< size_t, 3 > ortho_Nis={0, 0, 0})
vtkSmartPointer< vtkUnstructuredGrid > vtkUGridPtr
void ReadFromExodus(const Options &options)
std::vector< chi_mesh::Vertex > & GetVertices()
const Options & GetMeshOptions() const
void ReadFromVTU(const Options &options)
const std::vector< chi_mesh::Vertex > & GetVertices() const
void CopyUGridCellsAndPoints(vtkUnstructuredGrid &ugrid, double scale, int dimension_to_copy)
void ReadFromEnsightGold(const Options &options)
std::vector< LightWeightCell * > & GetRawCells()
static LightWeightCell * CreateCellFromVTKVertex(vtkCell *vtk_cell)
void ReadFromPVTU(const Options &options)
void SetMaterialIDsFromList(const std::vector< int > &material_ids)
std::vector< chi_mesh::Vertex > vertices_
std::vector< std::set< uint64_t > > vertex_cell_subscriptions_
void SetBoundaryIDsFromBlocks(std::vector< vtkUGridPtrAndName > &bndry_grid_blocks)
static LightWeightCell * CreateCellFromVTKPolygon(vtkCell *vtk_cell)
const std::vector< LightWeightCell * > & GetRawCells() const
void ReadFromWavefrontOBJ(const Options &options)
void ReadFromMsh(const Options &options)
std::vector< LightWeightCell * > raw_cells_
std::vector< LightWeightCell * > raw_boundary_cells_
void PushProxyCell(const std::string &type_str, const std::string &sub_type_str, int cell_num_faces, int cell_material_id, const std::vector< std::vector< uint64_t > > &proxy_faces)
static LightWeightCell * CreateCellFromVTKPolyhedron(vtkCell *vtk_cell)
std::shared_ptr< BoundBox > bound_box_
uint64_t MakeBoundaryID(const std::string &boundary_name)
const std::vector< std::set< uint64_t > > & GetVertextCellSubscriptions() const
void AddCell(LightWeightCell *&cell)
CellType
Definition: cell.h:12
MeshAttributes
Definition: chi_mesh.h:70
@ NONE
Definition: chi_mesh.h:71
LightWeightCell(chi_mesh::CellType in_type, chi_mesh::CellType in_sub_type)
LightWeightFace(std::vector< uint64_t > in_vertex_ids)
std::map< uint64_t, std::string > boundary_id_map