Chi-Tech
chi_meshcontinuum.h
Go to the documentation of this file.
1#ifndef CHI_MESHCONTINUUM_H_
2#define CHI_MESHCONTINUUM_H_
3
4#include <memory>
5#include <array>
6
7#include "../chi_mesh.h"
11
12#include "chi_mpi.h"
13
14namespace chi_data_types
15{
16template <typename T>
17class NDArray;
18} // namespace chi_data_types
19
20namespace chi
21{
22class ChiMPICommunicatorSet;
23}
24
25namespace chi_mesh
26{
27class GridFaceHistogram;
28class MeshGenerator;
29}
30
31// ######################################################### Class Definition
32/**Stores the relevant information for completely defining a computational
33 * domain. */
35{
36private:
37 typedef std::shared_ptr<chi::ChiMPICommunicatorSet> MPILocalCommSetPtr;
38
39private:
40 std::vector<std::unique_ptr<chi_mesh::Cell>>
41 local_cells_; ///< Actual local cells
42 std::vector<std::unique_ptr<chi_mesh::Cell>>
43 ghost_cells_; ///< Locally stored ghosts
44
45 std::map<uint64_t, uint64_t> global_cell_id_to_local_id_map_;
46 std::map<uint64_t, uint64_t> global_cell_id_to_nonlocal_id_map_;
47
49
50public:
54
55private:
57
58 struct
59 {
60 size_t Nx = 0;
61 size_t Ny = 0;
62 size_t Nz = 0;
64
65 std::map<uint64_t, std::string> boundary_id_map_;
66
67public:
74 {
75 }
76
77 void SetGlobalVertexCount(const uint64_t count)
78 {
80 }
81 uint64_t GetGlobalVertexCount() const { return global_vertex_count_; }
82
83 std::map<uint64_t, std::string>& GetBoundaryIDMap()
84 {
85 return boundary_id_map_;
86 }
87
88 const std::map<uint64_t, std::string>& GetBoundaryIDMap() const
89 {
90 return boundary_id_map_;
91 }
92
93 uint64_t MakeBoundaryID(const std::string& boundary_name) const;
94
95 static std::shared_ptr<MeshContinuum> New()
96 {
97 return std::make_shared<MeshContinuum>();
98 }
99
100 /**Method to be called if cells and nodes have been transferred
101 * to another grid.*/
103 {
104 local_cells_.clear();
105 ghost_cells_.clear();
108 vertices.Clear();
109 }
110
111 void ExportCellsToObj(const char* fileName,
112 bool per_material = false,
113 int options = 0) const;
114 void ExportCellsToVTK(const std::string& file_base_name) const;
115 void ExportCellsToExodus(const std::string& file_base_name,
116 bool suppress_node_sets = false,
117 bool suppress_side_sets = false) const;
118
119 std::shared_ptr<GridFaceHistogram>
120 MakeGridFaceHistogram(double master_tolerance = 100.0,
121 double slave_tolerance = 1.1) const;
122
123 bool IsCellLocal(uint64_t cell_global_index) const;
124 static int GetCellDimension(const chi_mesh::Cell& cell);
125
126 /**Creates a mapping of the current face local-ids to the
127 * adjacent face's local ids.*/
128 void FindAssociatedVertices(const chi_mesh::CellFace& cur_face,
129 std::vector<short>& dof_mapping) const;
130 /**Creates a mapping of the current face local-ids to the
131 * adjacent cell's local ids.*/
133 std::vector<short>& dof_mapping) const;
134 static size_t MapCellFace(const chi_mesh::Cell& cur_cell,
135 const chi_mesh::Cell& adj_cell,
136 unsigned int f);
137
138 /**Given a global-id of a cell, will return the local-id if the
139 * cell is local, otherwise will throw out_of_range.*/
140 size_t MapCellGlobalID2LocalID(uint64_t global_id) const;
141
143 ComputeCentroidFromListOfNodes(const std::vector<uint64_t>& list) const;
144
146
147 size_t GetGlobalNumberOfCells() const;
148
149 std::vector<uint64_t> GetDomainUniqueBoundaryIDs() const;
150
151 size_t
153 bool CheckPointInsideCell(const chi_mesh::Cell& cell,
154 const chi_mesh::Vector3& point) const;
155
157
158 std::array<size_t, 3> GetIJKInfo() const;
160 std::vector<chi_mesh::Vector3> MakeCellOrthoSizes() const;
161
162 std::pair<chi_mesh::Vector3, chi_mesh::Vector3> GetLocalBoundingBox() const;
163
164private:
168 std::array<size_t, 3> ortho_Nis = {0, 0, 0})
169 {
170 attributes = attributes | new_attribs;
171 ortho_attributes = {ortho_Nis[0], ortho_Nis[1], ortho_Nis[2]};
172 }
173};
174
175#endif // CHI_MESHCONTINUUM_H_
bool CheckPointInsideCell(const chi_mesh::Cell &cell, const chi_mesh::Vector3 &point) const
chi_mesh::Vector3 ComputeCentroidFromListOfNodes(const std::vector< uint64_t > &list) const
std::vector< std::unique_ptr< chi_mesh::Cell > > ghost_cells_
Locally stored ghosts.
void ExportCellsToExodus(const std::string &file_base_name, bool suppress_node_sets=false, bool suppress_side_sets=false) const
MeshAttributes Attributes() const
std::shared_ptr< GridFaceHistogram > MakeGridFaceHistogram(double master_tolerance=100.0, double slave_tolerance=1.1) const
LocalCellHandler local_cells
void FindAssociatedCellVertices(const chi_mesh::CellFace &cur_face, std::vector< short > &dof_mapping) const
MPILocalCommSetPtr MakeMPILocalCommunicatorSet() const
std::map< uint64_t, std::string > & GetBoundaryIDMap()
std::map< uint64_t, uint64_t > global_cell_id_to_local_id_map_
std::array< size_t, 3 > GetIJKInfo() const
std::vector< chi_mesh::Vector3 > MakeCellOrthoSizes() const
uint64_t GetGlobalVertexCount() const
chi_data_types::NDArray< uint64_t > MakeIJKToGlobalIDMapping() const
static std::shared_ptr< MeshContinuum > New()
std::map< uint64_t, std::string > boundary_id_map_
GlobalCellHandler cells
static int GetCellDimension(const chi_mesh::Cell &cell)
bool IsCellLocal(uint64_t cell_global_index) const
std::pair< chi_mesh::Vector3, chi_mesh::Vector3 > GetLocalBoundingBox() const
void SetAttributes(MeshAttributes new_attribs, std::array< size_t, 3 > ortho_Nis={0, 0, 0})
uint64_t MakeBoundaryID(const std::string &boundary_name) const
size_t MapCellGlobalID2LocalID(uint64_t global_id) const
void ExportCellsToVTK(const std::string &file_base_name) const
std::vector< std::unique_ptr< chi_mesh::Cell > > local_cells_
Actual local cells.
static size_t MapCellFace(const chi_mesh::Cell &cur_cell, const chi_mesh::Cell &adj_cell, unsigned int f)
const std::map< uint64_t, std::string > & GetBoundaryIDMap() const
std::shared_ptr< chi::ChiMPICommunicatorSet > MPILocalCommSetPtr
size_t CountCellsInLogicalVolume(const chi_mesh::LogicalVolume &log_vol) const
void ExportCellsToObj(const char *fileName, bool per_material=false, int options=0) const
void FindAssociatedVertices(const chi_mesh::CellFace &cur_face, std::vector< short > &dof_mapping) const
void SetGlobalVertexCount(const uint64_t count)
std::vector< uint64_t > GetDomainUniqueBoundaryIDs() const
std::map< uint64_t, uint64_t > global_cell_id_to_nonlocal_id_map_
struct chi_mesh::MeshContinuum::@0 ortho_attributes
std::shared_ptr< chi::ChiMPICommunicatorSet > MPILocalCommSetPtr
Definition: lbs_solver.h:31
MeshAttributes
Definition: chi_mesh.h:70
@ NONE
Definition: chi_mesh.h:71