OGS
MeshLib::IO::NodePartitionedMeshReader Class Referencefinal

Detailed Description

Class for parallel reading of binary partitioned mesh files into a NodePartitionedMesh via MPI.

Definition at line 37 of file NodePartitionedMeshReader.h.

#include <NodePartitionedMeshReader.h>

Collaboration diagram for MeshLib::IO::NodePartitionedMeshReader:
[legend]

Classes

struct  PartitionedMeshInfo
 A collection of integers that configure the partitioned mesh data. More...
 

Public Member Functions

 NodePartitionedMeshReader (MPI_Comm comm)
 
 ~NodePartitionedMeshReader ()
 
MeshLib::NodePartitionedMeshread (const std::string &file_name_base)
 Create a NodePartitionedMesh object, read data to it, and return a pointer to it. Data files are in binary format.
 

Private Member Functions

void registerNodeDataMpiType ()
 Define MPI data type for NodeData struct.
 
MeshLib::NodePartitionedMeshnewMesh (std::string const &mesh_name, std::vector< MeshLib::Node * > const &mesh_nodes, std::vector< unsigned long > const &glb_node_ids, std::vector< MeshLib::Element * > const &mesh_elems, MeshLib::Properties const &properties) const
 Create a new mesh of NodePartitionedMesh after reading and processing the data.
 
template<typename DATA >
bool readDataFromFile (std::string const &filename, MPI_Offset offset, MPI_Datatype type, DATA &data) const
 
MeshLib::NodePartitionedMeshreadMesh (const std::string &file_name_base)
 Create a NodePartitionedMesh object, read binary mesh data in the manner of parallel, and return a pointer to it. Four binary files have to been read in this function named as:
 
MeshLib::Properties readProperties (const std::string &file_name_base) const
 
void readProperties (const std::string &file_name_base, MeshLib::MeshItemType t, MeshLib::Properties &p) const
 
void readDomainSpecificPartOfPropertyVectors (std::vector< std::optional< MeshLib::IO::PropertyVectorMetaData > > const &vec_pvmd, MeshLib::IO::PropertyVectorPartitionMetaData const &pvpmd, MeshLib::MeshItemType t, std::istream &is, MeshLib::Properties &p) const
 
template<typename T >
void createPropertyVectorPart (std::istream &is, MeshLib::IO::PropertyVectorMetaData const &pvmd, MeshLib::IO::PropertyVectorPartitionMetaData const &pvpmd, MeshLib::MeshItemType t, unsigned long global_offset, MeshLib::Properties &p) const
 
template<typename T >
void createSpecificPropertyVectorPart (std::istream &is, MeshLib::IO::PropertyVectorMetaData const &pvmd, MeshLib::MeshItemType t, unsigned long global_offset, MeshLib::Properties &p) const
 
void setNodes (const std::vector< NodeData > &node_data, std::vector< MeshLib::Node * > &mesh_node, std::vector< unsigned long > &glb_node_ids) const
 Set mesh nodes from a temporary array containing node data read from file.
 
void setElements (const std::vector< MeshLib::Node * > &mesh_nodes, const std::vector< unsigned long > &elem_data, std::vector< MeshLib::Element * > &mesh_elems, const bool ghost=false) const
 Set mesh elements from a temporary array containing node data read from file.
 

Private Attributes

MPI_Comm _mpi_comm
 Pointer to MPI communicator.
 
int _mpi_comm_size
 Number of processes in the communicator: _mpi_comm.
 
int _mpi_rank
 Rank of compute core.
 
MPI_Datatype _mpi_node_type
 MPI data type for struct NodeData.
 
struct MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo _mesh_info
 

Constructor & Destructor Documentation

◆ NodePartitionedMeshReader()

MeshLib::IO::NodePartitionedMeshReader::NodePartitionedMeshReader ( MPI_Comm comm)
explicit
Parameters
commMPI communicator.

Definition at line 48 of file NodePartitionedMeshReader.cpp.

49 : _mpi_comm(comm)
50{
51 MPI_Comm_size(_mpi_comm, &_mpi_comm_size);
52 MPI_Comm_rank(_mpi_comm, &_mpi_rank);
53
55}
int _mpi_comm_size
Number of processes in the communicator: _mpi_comm.
MPI_Comm _mpi_comm
Pointer to MPI communicator.
void registerNodeDataMpiType()
Define MPI data type for NodeData struct.

References _mpi_comm, _mpi_comm_size, _mpi_rank, and registerNodeDataMpiType().

◆ ~NodePartitionedMeshReader()

MeshLib::IO::NodePartitionedMeshReader::~NodePartitionedMeshReader ( )

Definition at line 57 of file NodePartitionedMeshReader.cpp.

58{
59 MPI_Type_free(&_mpi_node_type);
60}
MPI_Datatype _mpi_node_type
MPI data type for struct NodeData.

References _mpi_node_type.

Member Function Documentation

◆ createPropertyVectorPart()

template<typename T >
void MeshLib::IO::NodePartitionedMeshReader::createPropertyVectorPart ( std::istream & is,
MeshLib::IO::PropertyVectorMetaData const & pvmd,
MeshLib::IO::PropertyVectorPartitionMetaData const & pvpmd,
MeshLib::MeshItemType t,
unsigned long global_offset,
MeshLib::Properties & p ) const
inlineprivate

Definition at line 188 of file NodePartitionedMeshReader.h.

193 {
194 MeshLib::PropertyVector<T>* pv = p.createNewPropertyVector<T>(
195 pvmd.property_name, t, pvmd.number_of_components);
196 pv->resize(pvpmd.number_of_tuples * pvmd.number_of_components);
197
198 // Locate the start position of the data in the file for the current
199 // rank.
200 is.seekg(global_offset +
201 pvpmd.offset * pvmd.number_of_components * sizeof(T));
202 // read the values
203 if (!is.read(reinterpret_cast<char*>(pv->data()),
204 pv->size() * sizeof(T)))
205 OGS_FATAL(
206 "Error in NodePartitionedMeshReader::readProperties: "
207 "Could not read part {:d} of the PropertyVector.",
208 _mpi_rank);
209 }
#define OGS_FATAL(...)
Definition Error.h:26
std::size_t size() const

References _mpi_rank, MeshLib::IO::PropertyVectorMetaData::number_of_components, MeshLib::IO::PropertyVectorPartitionMetaData::number_of_tuples, MeshLib::IO::PropertyVectorPartitionMetaData::offset, OGS_FATAL, MeshLib::IO::PropertyVectorMetaData::property_name, and MeshLib::PropertyVector< PROP_VAL_TYPE >::size().

◆ createSpecificPropertyVectorPart()

template<typename T >
void MeshLib::IO::NodePartitionedMeshReader::createSpecificPropertyVectorPart ( std::istream & is,
MeshLib::IO::PropertyVectorMetaData const & pvmd,
MeshLib::MeshItemType t,
unsigned long global_offset,
MeshLib::Properties & p ) const
inlineprivate

Read data for property OGS_VERSION or IntegrationPointMetaData, and create property vector for it.

Definition at line 214 of file NodePartitionedMeshReader.h.

218 {
219 MeshLib::PropertyVector<T>* pv = p.createNewPropertyVector<T>(
220 pvmd.property_name, t, pvmd.number_of_components);
221
222 std::size_t const property_vector_size =
223 pvmd.number_of_tuples / _mpi_comm_size;
224 pv->resize(property_vector_size);
225
226 // Locate the start position of the data in the file for the current
227 // rank.
228 is.seekg(global_offset + property_vector_size * sizeof(T) * _mpi_rank);
229
230 // read the values
231 if (!is.read(reinterpret_cast<char*>(pv->data()),
232 pv->size() * sizeof(T)))
233 {
234 OGS_FATAL(
235 "Error in NodePartitionedMeshReader::readProperties: "
236 "Could not read part {:d} of the PropertyVector.",
237 _mpi_rank);
238 }
239 }

References _mpi_comm_size, _mpi_rank, MeshLib::IO::PropertyVectorMetaData::number_of_components, MeshLib::IO::PropertyVectorMetaData::number_of_tuples, OGS_FATAL, MeshLib::IO::PropertyVectorMetaData::property_name, and MeshLib::PropertyVector< PROP_VAL_TYPE >::size().

◆ newMesh()

MeshLib::NodePartitionedMesh * MeshLib::IO::NodePartitionedMeshReader::newMesh ( std::string const & mesh_name,
std::vector< MeshLib::Node * > const & mesh_nodes,
std::vector< unsigned long > const & glb_node_ids,
std::vector< MeshLib::Element * > const & mesh_elems,
MeshLib::Properties const & properties ) const
private

Create a new mesh of NodePartitionedMesh after reading and processing the data.

Parameters
mesh_nameName assigned to the new mesh.
mesh_nodesNode data.
glb_node_idsGlobal IDs of nodes.
mesh_elemsElement data.
propertiesCollection of PropertyVector's assigned to the mesh.
Returns
Returns a pointer to a NodePartitionedMesh

Definition at line 403 of file NodePartitionedMeshReader.cpp.

409{
410 std::vector<std::size_t> gathered_n_regular_base_nodes(_mpi_comm_size);
411
413 1,
414 MPI_UNSIGNED_LONG,
415 gathered_n_regular_base_nodes.data(),
416 1,
417 MPI_UNSIGNED_LONG,
418 _mpi_comm);
419
420 std::vector<std::size_t> n_regular_base_nodes_at_rank;
421 n_regular_base_nodes_at_rank.push_back(0);
422 std::partial_sum(begin(gathered_n_regular_base_nodes),
423 end(gathered_n_regular_base_nodes),
424 back_inserter(n_regular_base_nodes_at_rank));
425
426 std::vector<std::size_t> gathered_n_regular_high_order_nodes(
428 std::size_t const n_regular_high_order_nodes =
431 MPI_Allgather(&n_regular_high_order_nodes,
432 1,
433 MPI_UNSIGNED_LONG,
434 gathered_n_regular_high_order_nodes.data(),
435 1,
436 MPI_UNSIGNED_LONG,
437 _mpi_comm);
438
439 std::vector<std::size_t> n_regular_high_order_nodes_at_rank;
440 n_regular_high_order_nodes_at_rank.push_back(0);
441 std::partial_sum(begin(gathered_n_regular_high_order_nodes),
442 end(gathered_n_regular_high_order_nodes),
443 back_inserter(n_regular_high_order_nodes_at_rank));
444
446 mesh_name, mesh_nodes, glb_node_ids, mesh_elems, properties,
449 std::move(n_regular_base_nodes_at_rank),
450 std::move(n_regular_high_order_nodes_at_rank));
451}
struct MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo _mesh_info
unsigned long number_of_global_nodes
7: Number of all nodes of global mesh,

References _mesh_info, _mpi_comm, _mpi_comm_size, MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::number_of_global_base_nodes, MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::number_of_global_nodes, MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::number_of_regular_base_nodes, and MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::number_of_regular_nodes.

Referenced by readMesh().

◆ read()

MeshLib::NodePartitionedMesh * MeshLib::IO::NodePartitionedMeshReader::read ( const std::string & file_name_base)

Create a NodePartitionedMesh object, read data to it, and return a pointer to it. Data files are in binary format.

Parameters
file_name_baseName of file to be read, and it must be base name without name extension.
Returns
Pointer to Mesh object. If the creation of mesh object fails, it returns a null pointer.

Definition at line 74 of file NodePartitionedMeshReader.cpp.

76{
77 BaseLib::RunTime timer;
78 timer.start();
79
80 // Always try binary file first
81 std::string const fname_new = file_name_base + "_partitioned_msh_cfg" +
82 std::to_string(_mpi_comm_size) + ".bin";
83
84 if (!BaseLib::IsFileExisting(fname_new)) // binary file does not exist.
85 {
86 std::string const fname_ascii = file_name_base +
87 "_partitioned_msh_cfg" +
88 std::to_string(_mpi_comm_size) + ".msh";
89 if (BaseLib::IsFileExisting(fname_ascii))
90 {
91 ERR("Reading of ASCII mesh file {:s} is not supported since OGS "
92 "version 6.3.3.",
93 fname_ascii);
94 }
95 OGS_FATAL("Required binary file {:s} does not exist.\n", fname_new);
96 }
97
98 INFO("Reading corresponding part of mesh data from binary file {:s} ...",
99 file_name_base);
100
101 MeshLib::NodePartitionedMesh* mesh = readMesh(file_name_base);
102
103 INFO("[time] Reading the mesh took {:f} s.", timer.elapsed());
104
105 MPI_Barrier(_mpi_comm);
106
107 return mesh;
108}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
Count the running time.
Definition RunTime.h:29
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:42
void start()
Start the timer.
Definition RunTime.h:32
MeshLib::NodePartitionedMesh * readMesh(const std::string &file_name_base)
Create a NodePartitionedMesh object, read binary mesh data in the manner of parallel,...
bool IsFileExisting(const std::string &strFilename)
Returns true if given file exists.
Definition FileTools.cpp:47

References _mpi_comm, _mpi_comm_size, BaseLib::RunTime::elapsed(), ERR(), INFO(), BaseLib::IsFileExisting(), OGS_FATAL, readMesh(), and BaseLib::RunTime::start().

Referenced by MeshLib::IO::readMeshFromFile().

◆ readDataFromFile()

template<typename DATA >
bool MeshLib::IO::NodePartitionedMeshReader::readDataFromFile ( std::string const & filename,
MPI_Offset offset,
MPI_Datatype type,
DATA & data ) const
private

Parallel reading of a binary file via MPI_File_read, reading mesh data head, nodes, non-ghost elements and ghost elements.

Note
In case of failure during opening of the file, an error message is printed.
If the number of elements in container is larger than MPI_file_read() supports (maximum of current int type), an error is printed.
Parameters
filenameFile name containing data.
offsetDisplacement of the data accessible from the view. see MPI_File_set_view() documentation.
typeType of data.
dataA container to be filled with data. Its size is used to determine how many values should be read.
Template Parameters
DATAA homogeneous container type supporting data() and size().
Returns
True on success and false otherwise.

Definition at line 111 of file NodePartitionedMeshReader.cpp.

115{
116 // Check container size
117 if (!is_safely_convertable<std::size_t, int>(data.size()))
118 {
119 ERR("The container size is too large for MPI_File_read() call.");
120 return false;
121 }
122
123 // Open file
124 MPI_File file;
125
126 char* filename_char = const_cast<char*>(filename.data());
127 int const file_status = MPI_File_open(
128 _mpi_comm, filename_char, MPI_MODE_RDONLY, MPI_INFO_NULL, &file);
129
130 if (file_status != 0)
131 {
132 ERR("Error opening file {:s}. MPI error code {:d}", filename,
133 file_status);
134 return false;
135 }
136
137 // Read data
138 char file_mode[] = "native";
139 MPI_File_set_view(file, offset, type, type, file_mode, MPI_INFO_NULL);
140 // The static cast is checked above.
141 MPI_File_read(file, data.data(), static_cast<int>(data.size()), type,
142 MPI_STATUS_IGNORE);
143 MPI_File_close(&file);
144
145 return true;
146}

References _mpi_comm, and ERR().

Referenced by readMesh().

◆ readDomainSpecificPartOfPropertyVectors()

void MeshLib::IO::NodePartitionedMeshReader::readDomainSpecificPartOfPropertyVectors ( std::vector< std::optional< MeshLib::IO::PropertyVectorMetaData > > const & vec_pvmd,
MeshLib::IO::PropertyVectorPartitionMetaData const & pvpmd,
MeshLib::MeshItemType t,
std::istream & is,
MeshLib::Properties & p ) const
private

Definition at line 300 of file NodePartitionedMeshReader.cpp.

307{
308 unsigned long global_offset = 0;
309 std::size_t const number_of_properties = vec_pvmd.size();
310 for (std::size_t i(0); i < number_of_properties; ++i)
311 {
312 DBUG("global offset: {:d}, offset within the PropertyVector: {:d}.",
313 global_offset,
314 global_offset + pvpmd.offset * vec_pvmd[i]->number_of_components *
315 vec_pvmd[i]->data_type_size_in_bytes);
316
317 // Special field data such as OGS_VERSION, IntegrationPointMetaData,
318 // etc., which are not "real" integration points, are copied "as is"
319 // (i.e. fully) for every partition.
320 if (vec_pvmd[i]->property_name.find("_ip") == std::string::npos &&
322 {
323 createSpecificPropertyVectorPart<char>(is, *vec_pvmd[i], t,
324 global_offset, p);
325
326 global_offset += vec_pvmd[i]->data_type_size_in_bytes *
327 vec_pvmd[i]->number_of_tuples;
328 continue;
329 }
330
331 if (vec_pvmd[i]->is_int_type)
332 {
333 if (vec_pvmd[i]->is_data_type_signed)
334 {
335 if (vec_pvmd[i]->data_type_size_in_bytes == sizeof(char))
336 {
337 createPropertyVectorPart<char>(is, *vec_pvmd[i], pvpmd, t,
338 global_offset, p);
339 }
340 else if (vec_pvmd[i]->data_type_size_in_bytes == sizeof(int))
341 {
342 createPropertyVectorPart<int>(is, *vec_pvmd[i], pvpmd, t,
343 global_offset, p);
344 }
345 else if (vec_pvmd[i]->data_type_size_in_bytes == sizeof(long))
346 createPropertyVectorPart<long>(is, *vec_pvmd[i], pvpmd, t,
347 global_offset, p);
348 else
349 {
350 WARN(
351 "Implementation for reading signed integer property "
352 "vector '{:s}' is not available.",
353 vec_pvmd[i]->property_name);
354 }
355 }
356 else
357 {
358 if (vec_pvmd[i]->data_type_size_in_bytes ==
359 sizeof(unsigned char))
360 {
361 createPropertyVectorPart<unsigned char>(
362 is, *vec_pvmd[i], pvpmd, t, global_offset, p);
363 }
364 else if (vec_pvmd[i]->data_type_size_in_bytes ==
365 sizeof(unsigned int))
366 createPropertyVectorPart<unsigned int>(
367 is, *vec_pvmd[i], pvpmd, t, global_offset, p);
368 else if (vec_pvmd[i]->data_type_size_in_bytes ==
369 sizeof(unsigned long))
370 createPropertyVectorPart<unsigned long>(
371 is, *vec_pvmd[i], pvpmd, t, global_offset, p);
372 else
373 {
374 WARN(
375 "Implementation for reading unsigned property vector "
376 "'{:s}' is not available.",
377 vec_pvmd[i]->property_name);
378 }
379 }
380 }
381 else
382 {
383 if (vec_pvmd[i]->data_type_size_in_bytes == sizeof(float))
384 createPropertyVectorPart<float>(is, *vec_pvmd[i], pvpmd, t,
385 global_offset, p);
386 else if (vec_pvmd[i]->data_type_size_in_bytes == sizeof(double))
387 createPropertyVectorPart<double>(is, *vec_pvmd[i], pvpmd, t,
388 global_offset, p);
389 else
390 {
391 WARN(
392 "Implementation for reading floating point property vector "
393 "'{:s}' is not available.",
394 vec_pvmd[i]->property_name);
395 }
396 }
397 global_offset += vec_pvmd[i]->data_type_size_in_bytes *
398 vec_pvmd[i]->number_of_tuples *
399 vec_pvmd[i]->number_of_components;
400 }
401}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40

References DBUG(), MeshLib::IntegrationPoint, MeshLib::IO::PropertyVectorPartitionMetaData::offset, and WARN().

Referenced by readProperties().

◆ readMesh()

MeshLib::NodePartitionedMesh * MeshLib::IO::NodePartitionedMeshReader::readMesh ( const std::string & file_name_base)
private

Create a NodePartitionedMesh object, read binary mesh data in the manner of parallel, and return a pointer to it. Four binary files have to been read in this function named as:

  • file_name_base+_partitioned_msh_cfg[number of partitions].bin
  • file_name_base+_partitioned_msh_nod[number of partitions].bin
  • file_name_base+_partitioned_msh_ele[number of partitions].bin
  • file_name_base+_partitioned_msh_ele_g[number of partitions].bin

in which, the first file contains an array of integers for the PartitionMeshInfo for all partitions, the second file contains a struct type (long, double double double) array of nodes information of global IDs and coordinates of all partitions, the third file contains a long type integer array of element information of material ID, element type and node IDs of each non-ghost element of all partitions, and the forth file contains a long type integer array of element information of material ID, element type and node IDs of each ghost element of all partitions.

Parameters
file_name_baseName of file to be read, which must be a name with the path to the file and without file extension.
Returns
Pointer to Mesh object.

Definition at line 148 of file NodePartitionedMeshReader.cpp.

150{
151 //----------------------------------------------------------------------------------
152 // Read headers
153 const std::string fname_header = file_name_base + "_partitioned_msh_";
154 const std::string fname_num_p_ext = std::to_string(_mpi_comm_size) + ".bin";
155
156 // Read the config meta data from *cfg* file into struct PartitionedMeshInfo
157 // _mesh_info
158 if (!readDataFromFile(
159 fname_header + "cfg" + fname_num_p_ext,
160 static_cast<MPI_Offset>(static_cast<unsigned>(_mpi_rank) *
161 sizeof(_mesh_info)),
162 MPI_LONG, _mesh_info))
163 return nullptr;
164
165 //----------------------------------------------------------------------------------
166 // Read Nodes
167 std::vector<NodeData> nodes(_mesh_info.number_of_nodes);
168
169 if (!readDataFromFile(fname_header + "nod" + fname_num_p_ext,
170 static_cast<MPI_Offset>(_mesh_info.offset[2]),
171 _mpi_node_type, nodes))
172 return nullptr;
173
174 std::vector<MeshLib::Node*> mesh_nodes;
175 std::vector<unsigned long> glb_node_ids;
176 setNodes(nodes, mesh_nodes, glb_node_ids);
177
178 //----------------------------------------------------------------------------------
179 // Read non-ghost elements
180 std::vector<unsigned long> elem_data(_mesh_info.number_of_regular_elements +
181 _mesh_info.offset[0]);
182 if (!readDataFromFile(fname_header + "ele" + fname_num_p_ext,
183 static_cast<MPI_Offset>(_mesh_info.offset[3]),
184 MPI_LONG, elem_data))
185 return nullptr;
186
187 std::vector<MeshLib::Element*> mesh_elems(
190 setElements(mesh_nodes, elem_data, mesh_elems);
191
192 //----------------------------------------------------------------------------------
193 // Read ghost element
194 std::vector<unsigned long> ghost_elem_data(
196
197 if (!readDataFromFile(fname_header + "ele_g" + fname_num_p_ext,
198 static_cast<MPI_Offset>(_mesh_info.offset[4]),
199 MPI_LONG, ghost_elem_data))
200 return nullptr;
201
202 const bool process_ghost = true;
203 setElements(mesh_nodes, ghost_elem_data, mesh_elems, process_ghost);
204
205 //----------------------------------------------------------------------------------
206 // read the properties
207 MeshLib::Properties p(readProperties(file_name_base));
208
209 return newMesh(BaseLib::extractBaseName(file_name_base), mesh_nodes,
210 glb_node_ids, mesh_elems, p);
211}
void setElements(const std::vector< MeshLib::Node * > &mesh_nodes, const std::vector< unsigned long > &elem_data, std::vector< MeshLib::Element * > &mesh_elems, const bool ghost=false) const
Set mesh elements from a temporary array containing node data read from file.
bool readDataFromFile(std::string const &filename, MPI_Offset offset, MPI_Datatype type, DATA &data) const
void setNodes(const std::vector< NodeData > &node_data, std::vector< MeshLib::Node * > &mesh_node, std::vector< unsigned long > &glb_node_ids) const
Set mesh nodes from a temporary array containing node data read from file.
MeshLib::Properties readProperties(const std::string &file_name_base) const
MeshLib::NodePartitionedMesh * newMesh(std::string const &mesh_name, std::vector< MeshLib::Node * > const &mesh_nodes, std::vector< unsigned long > const &glb_node_ids, std::vector< MeshLib::Element * > const &mesh_elems, MeshLib::Properties const &properties) const
Create a new mesh of NodePartitionedMesh after reading and processing the data.
Property manager on mesh items. Class Properties manages scalar, vector or matrix properties....
Definition Properties.h:36
std::string extractBaseName(std::string const &pathname)
unsigned long number_of_nodes
0: Number of all nodes of a partition,

References _mesh_info, _mpi_comm_size, _mpi_node_type, _mpi_rank, BaseLib::extractBaseName(), newMesh(), MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::number_of_ghost_elements, MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::number_of_nodes, MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::number_of_regular_elements, MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::offset, readDataFromFile(), readProperties(), setElements(), and setNodes().

Referenced by read().

◆ readProperties() [1/2]

MeshLib::Properties MeshLib::IO::NodePartitionedMeshReader::readProperties ( const std::string & file_name_base) const
private

◆ readProperties() [2/2]

void MeshLib::IO::NodePartitionedMeshReader::readProperties ( const std::string & file_name_base,
MeshLib::MeshItemType t,
MeshLib::Properties & p ) const
private

Definition at line 223 of file NodePartitionedMeshReader.cpp.

226{
227 const std::string item_type = MeshLib::toString(t);
228
229 const std::string fname_cfg = file_name_base + "_partitioned_" + item_type +
230 "_properties_cfg" +
231 std::to_string(_mpi_comm_size) + ".bin";
232 std::ifstream is(fname_cfg.c_str(), std::ios::binary | std::ios::in);
233 if (!is)
234 {
235 WARN(
236 "Could not open file '{:s}'.\n"
237 "\tYou can ignore this warning if the mesh does not contain {:s}-"
238 "wise property data.",
239 fname_cfg, item_type.data());
240 return;
241 }
242 std::size_t number_of_properties = 0;
243 is.read(reinterpret_cast<char*>(&number_of_properties),
244 sizeof(std::size_t));
245 std::vector<std::optional<MeshLib::IO::PropertyVectorMetaData>> vec_pvmd(
246 number_of_properties);
247 for (std::size_t i(0); i < number_of_properties; ++i)
248 {
250 if (!vec_pvmd[i])
251 {
252 OGS_FATAL(
253 "Error in NodePartitionedMeshReader::readProperties: "
254 "Could not read the meta data for the PropertyVector {:d}",
255 i);
256 }
257 }
258 for (std::size_t i(0); i < number_of_properties; ++i)
259 {
261 }
262 auto pos = is.tellg();
263 auto offset =
264 static_cast<long>(pos) +
265 static_cast<long>(_mpi_rank *
267 is.seekg(offset);
268 std::optional<MeshLib::IO::PropertyVectorPartitionMetaData> pvpmd(
270 bool pvpmd_read_ok = static_cast<bool>(pvpmd);
271 bool all_pvpmd_read_ok;
272 MPI_Allreduce(&pvpmd_read_ok, &all_pvpmd_read_ok, 1, MPI_C_BOOL, MPI_LOR,
273 _mpi_comm);
274 if (!all_pvpmd_read_ok)
275 {
276 OGS_FATAL(
277 "Error in NodePartitionedMeshReader::readProperties: "
278 "Could not read the partition meta data for the mpi process {:d}",
279 _mpi_rank);
280 }
281 DBUG("offset in the PropertyVector: {:d}", pvpmd->offset);
282 DBUG("{:d} tuples in partition.", pvpmd->number_of_tuples);
283 is.close();
284
285 const std::string fname_val = file_name_base + "_partitioned_" + item_type +
286 "_properties_val" +
287 std::to_string(_mpi_comm_size) + ".bin";
288 is.open(fname_val.c_str(), std::ios::binary | std::ios::in);
289 if (!is)
290 {
291 ERR("Could not open file '{:s}'\n."
292 "\tYou can ignore this warning if the mesh does not contain {:s}-"
293 "wise property data.",
294 fname_val, item_type.data());
295 }
296
297 readDomainSpecificPartOfPropertyVectors(vec_pvmd, *pvpmd, t, is, p);
298}
void readDomainSpecificPartOfPropertyVectors(std::vector< std::optional< MeshLib::IO::PropertyVectorMetaData > > const &vec_pvmd, MeshLib::IO::PropertyVectorPartitionMetaData const &pvpmd, MeshLib::MeshItemType t, std::istream &is, MeshLib::Properties &p) const
void writePropertyVectorMetaData(std::ostream &os, PropertyVectorMetaData const &pvmd)
std::optional< PropertyVectorMetaData > readPropertyVectorMetaData(std::istream &is)
std::optional< PropertyVectorPartitionMetaData > readPropertyVectorPartitionMetaData(std::istream &is)
static constexpr char const * toString(const MeshItemType t)
Returns a char array for a specific MeshItemType.
Definition Location.h:28

References _mpi_comm, _mpi_comm_size, _mpi_rank, DBUG(), ERR(), OGS_FATAL, readDomainSpecificPartOfPropertyVectors(), MeshLib::IO::readPropertyVectorMetaData(), MeshLib::IO::readPropertyVectorPartitionMetaData(), MeshLib::toString(), WARN(), and MeshLib::IO::writePropertyVectorMetaData().

◆ registerNodeDataMpiType()

void MeshLib::IO::NodePartitionedMeshReader::registerNodeDataMpiType ( )
private

Define MPI data type for NodeData struct.

Definition at line 62 of file NodePartitionedMeshReader.cpp.

63{
64 int const count = 2;
65 int blocks[count] = {1, 3};
66 MPI_Datatype types[count] = {MPI_UNSIGNED_LONG, MPI_DOUBLE};
67 MPI_Aint displacements[count] = {0, sizeof(NodeData::index)};
68
69 MPI_Type_create_struct(count, blocks, displacements, types,
71 MPI_Type_commit(&_mpi_node_type);
72}
std::size_t index
Global node index.
Definition NodeData.h:26

References _mpi_node_type, and MeshLib::IO::NodeData::index.

Referenced by NodePartitionedMeshReader().

◆ setElements()

void MeshLib::IO::NodePartitionedMeshReader::setElements ( const std::vector< MeshLib::Node * > & mesh_nodes,
const std::vector< unsigned long > & elem_data,
std::vector< MeshLib::Element * > & mesh_elems,
const bool ghost = false ) const
private

Set mesh elements from a temporary array containing node data read from file.

Parameters
mesh_nodesVector of mesh nodes used to set element nodes.
elem_dataVector containing element data read from file.
mesh_elemsVector of mesh elements to be set.
ghostFlag of processing ghost elements.

Definition at line 469 of file NodePartitionedMeshReader.cpp.

473{
474 // Number of elements, either ghost or regular
475 unsigned long const ne = ghost ? _mesh_info.number_of_ghost_elements
477 unsigned long const id_offset_ghost =
479
480 for (unsigned long i = 0; i < ne; i++)
481 {
482 unsigned long id_offset_elem = elem_data[i];
483
484 // Unused for now, keep for elem_data documentation purpose here.
485 {
486 const unsigned mat_idx =
487 static_cast<unsigned>(elem_data[id_offset_elem++]);
488 (void)mat_idx;
489 }
490 const unsigned long e_type = elem_data[id_offset_elem++];
491 unsigned long const nnodes = elem_data[id_offset_elem++];
492
493 MeshLib::Node** elem_nodes = new MeshLib::Node*[nnodes];
494 for (unsigned long k = 0; k < nnodes; k++)
495 elem_nodes[k] = mesh_nodes[elem_data[id_offset_elem++]];
496
497 // The element types below are defined by the MeshLib::CellType.
498 switch (static_cast<CellType>(e_type))
499 {
500 case CellType::POINT1:
501 mesh_elems[i + id_offset_ghost] =
502 new MeshLib::Point(elem_nodes);
503 break;
504 case CellType::LINE2:
505 mesh_elems[i + id_offset_ghost] = new MeshLib::Line(elem_nodes);
506 break;
507 case CellType::LINE3:
508 mesh_elems[i + id_offset_ghost] =
509 new MeshLib::Line3(elem_nodes);
510 break;
511 case CellType::QUAD4:
512 mesh_elems[i + id_offset_ghost] = new MeshLib::Quad(elem_nodes);
513 break;
514 case CellType::QUAD8:
515 mesh_elems[i + id_offset_ghost] =
516 new MeshLib::Quad8(elem_nodes);
517 break;
518 case CellType::QUAD9:
519 mesh_elems[i + id_offset_ghost] =
520 new MeshLib::Quad9(elem_nodes);
521 break;
522 case CellType::HEX8:
523 mesh_elems[i + id_offset_ghost] = new MeshLib::Hex(elem_nodes);
524 break;
525 case CellType::HEX20:
526 mesh_elems[i + id_offset_ghost] =
527 new MeshLib::Hex20(elem_nodes);
528 break;
529 case CellType::HEX27:
530 OGS_FATAL(
531 "NodePartitionedMeshReader: construction of HEX27 element "
532 "with id {:d} is not implemented.",
533 i);
534 break;
535 case CellType::TRI3:
536 mesh_elems[i + id_offset_ghost] = new MeshLib::Tri(elem_nodes);
537 break;
538 case CellType::TRI6:
539 mesh_elems[i + id_offset_ghost] = new MeshLib::Tri6(elem_nodes);
540 break;
541 case CellType::TET4:
542 mesh_elems[i + id_offset_ghost] = new MeshLib::Tet(elem_nodes);
543 break;
544 case CellType::TET10:
545 mesh_elems[i + id_offset_ghost] =
546 new MeshLib::Tet10(elem_nodes);
547 break;
548 case CellType::PRISM6:
549 mesh_elems[i + id_offset_ghost] =
550 new MeshLib::Prism(elem_nodes);
551 break;
553 mesh_elems[i + id_offset_ghost] =
554 new MeshLib::Prism15(elem_nodes);
555 break;
557 mesh_elems[i + id_offset_ghost] =
558 new MeshLib::Pyramid(elem_nodes);
559 break;
561 mesh_elems[i + id_offset_ghost] =
562 new MeshLib::Pyramid13(elem_nodes);
563 break;
565 OGS_FATAL(
566 "NodePartitionedMeshReader: construction of INVALID "
567 "element type with id {:d} is not possible.",
568 i);
569 break;
570 default:
571 OGS_FATAL(
572 "NodePartitionedMeshReader: construction of element type "
573 "{:d} is not implemented.",
574 e_type);
575 }
576 }
577}
TemplateElement< MeshLib::QuadRule9 > Quad9
Definition Quad.h:30
TemplateElement< MeshLib::TetRule10 > Tet10
Definition Tet.h:26
CellType
Types of mesh elements supported by OpenGeoSys.
Definition MeshEnums.h:43
TemplateElement< MeshLib::HexRule20 > Hex20
Definition Hex.h:26
TemplateElement< MeshLib::TetRule4 > Tet
Definition Tet.h:25
TemplateElement< MeshLib::LineRule2 > Line
Definition Line.h:25
TemplateElement< MeshLib::QuadRule8 > Quad8
Definition Quad.h:29
TemplateElement< MeshLib::PyramidRule13 > Pyramid13
Definition Pyramid.h:26
TemplateElement< MeshLib::QuadRule4 > Quad
Definition Quad.h:28
TemplateElement< MeshLib::TriRule3 > Tri
Definition Tri.h:26
TemplateElement< MeshLib::PyramidRule5 > Pyramid
Definition Pyramid.h:25
TemplateElement< MeshLib::PrismRule6 > Prism
Definition Prism.h:25
TemplateElement< PointRule1 > Point
Definition Point.h:20
TemplateElement< MeshLib::PrismRule15 > Prism15
Definition Prism.h:26
TemplateElement< MeshLib::TriRule6 > Tri6
Definition Tri.h:27
TemplateElement< MeshLib::LineRule3 > Line3
Definition Line.h:26
TemplateElement< MeshLib::HexRule8 > Hex
Definition Hex.h:25

References _mesh_info, MeshLib::HEX20, MeshLib::HEX27, MeshLib::HEX8, MeshLib::INVALID, MeshLib::LINE2, MeshLib::LINE3, MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::number_of_ghost_elements, MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::number_of_regular_elements, OGS_FATAL, MeshLib::POINT1, MeshLib::PRISM15, MeshLib::PRISM6, MeshLib::PYRAMID13, MeshLib::PYRAMID5, MeshLib::QUAD4, MeshLib::QUAD8, MeshLib::QUAD9, MeshLib::TET10, MeshLib::TET4, MeshLib::TRI3, and MeshLib::TRI6.

Referenced by readMesh().

◆ setNodes()

void MeshLib::IO::NodePartitionedMeshReader::setNodes ( const std::vector< NodeData > & node_data,
std::vector< MeshLib::Node * > & mesh_node,
std::vector< unsigned long > & glb_node_ids ) const
private

Set mesh nodes from a temporary array containing node data read from file.

Parameters
node_dataVector containing node data read from file.
mesh_nodeVector of mesh nodes to be set.
glb_node_idsGlobal IDs of nodes of a partition.

Definition at line 453 of file NodePartitionedMeshReader.cpp.

457{
458 mesh_node.resize(_mesh_info.number_of_nodes);
459 glb_node_ids.resize(_mesh_info.number_of_nodes);
460
461 for (std::size_t i = 0; i < mesh_node.size(); i++)
462 {
463 NodeData const& nd = node_data[i];
464 glb_node_ids[i] = nd.index;
465 mesh_node[i] = new MeshLib::Node(nd.x, nd.y, nd.z, i);
466 }
467}

References _mesh_info, MeshLib::IO::NodeData::index, MeshLib::Node, MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo::number_of_nodes, MeshLib::IO::NodeData::x, MeshLib::IO::NodeData::y, and MeshLib::IO::NodeData::z.

Referenced by readMesh().

Member Data Documentation

◆ _mesh_info

struct MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo MeshLib::IO::NodePartitionedMeshReader::_mesh_info
private

◆ _mpi_comm

MPI_Comm MeshLib::IO::NodePartitionedMeshReader::_mpi_comm
private

Pointer to MPI communicator.

Definition at line 57 of file NodePartitionedMeshReader.h.

Referenced by NodePartitionedMeshReader(), newMesh(), read(), readDataFromFile(), and readProperties().

◆ _mpi_comm_size

int MeshLib::IO::NodePartitionedMeshReader::_mpi_comm_size
private

Number of processes in the communicator: _mpi_comm.

Definition at line 60 of file NodePartitionedMeshReader.h.

Referenced by NodePartitionedMeshReader(), createSpecificPropertyVectorPart(), newMesh(), read(), readMesh(), and readProperties().

◆ _mpi_node_type

MPI_Datatype MeshLib::IO::NodePartitionedMeshReader::_mpi_node_type
private

MPI data type for struct NodeData.

Definition at line 66 of file NodePartitionedMeshReader.h.

Referenced by ~NodePartitionedMeshReader(), readMesh(), and registerNodeDataMpiType().

◆ _mpi_rank

int MeshLib::IO::NodePartitionedMeshReader::_mpi_rank
private

The documentation for this class was generated from the following files: