32template <
typename VALUE,
typename TYPE>
35 bool const result = value <= std::numeric_limits<TYPE>::max();
38 ERR(
"The value {:d} is too large for conversion.", value);
39 ERR(
"Maximum available size is {:d}.",
40 std::numeric_limits<TYPE>::max());
62 int blocks[count] = {1, 3};
63 MPI_Datatype types[count] = {MPI_UNSIGNED_LONG, MPI_DOUBLE};
66 MPI_Type_create_struct(count, blocks, displacements, types,
72 const std::string& file_name_base)
78 std::string
const fname_new = file_name_base +
"_partitioned_msh_cfg" +
83 std::string
const fname_ascii = file_name_base +
84 "_partitioned_msh_cfg" +
88 ERR(
"Reading of ASCII mesh file {:s} is not supported since OGS "
92 OGS_FATAL(
"Required binary file {:s} does not exist.\n", fname_new);
95 INFO(
"Reading corresponding part of mesh data from binary file {:s} ...",
100 INFO(
"[time] Reading the mesh took {:f} s.", timer.
elapsed());
107template <
typename DATA>
116 ERR(
"The container size is too large for MPI_File_read() call.");
123 char* filename_char =
const_cast<char*
>(filename.data());
124 int const file_status =
126 MPI_INFO_NULL, &file);
128 if (file_status != 0)
130 ERR(
"Error opening file {:s}. MPI error code {:d}", filename,
136 char file_mode[] =
"native";
137 MPI_File_set_view(file, offset, type, type, file_mode, MPI_INFO_NULL);
139 MPI_File_read(file, data.data(),
static_cast<int>(data.size()), type,
141 MPI_File_close(&file);
147 const std::string& file_name_base)
151 const std::string fname_header = file_name_base +
"_partitioned_msh_";
152 const std::string fname_num_p_ext = std::to_string(
mpi_.
size) +
".bin";
157 fname_header +
"cfg" + fname_num_p_ext,
158 static_cast<MPI_Offset
>(
static_cast<unsigned>(
mpi_.
rank) *
172 std::vector<MeshLib::Node*> mesh_nodes;
173 std::vector<unsigned long> glb_node_ids;
174 setNodes(nodes, mesh_nodes, glb_node_ids);
182 MPI_LONG, elem_data))
185 std::vector<MeshLib::Element*> mesh_elems(
192 std::vector<unsigned long> ghost_elem_data(
197 MPI_LONG, ghost_elem_data))
200 const bool process_ghost =
true;
201 setElements(mesh_nodes, ghost_elem_data, mesh_elems, process_ghost);
208 glb_node_ids, mesh_elems, p);
212 const std::string& file_name_base)
const
227 const std::string fname_cfg = file_name_base +
"_partitioned_" + item_type +
230 std::ifstream is(fname_cfg.c_str(), std::ios::binary | std::ios::in);
234 "Could not open file '{:s}'.\n"
235 "\tYou can ignore this warning if the mesh does not contain {:s}-"
236 "wise property data.",
237 fname_cfg, item_type.data());
240 std::size_t number_of_properties = 0;
241 is.read(
reinterpret_cast<char*
>(&number_of_properties),
242 sizeof(std::size_t));
243 std::vector<std::optional<MeshLib::IO::PropertyVectorMetaData>> vec_pvmd(
244 number_of_properties);
245 for (std::size_t i(0); i < number_of_properties; ++i)
251 "Error in NodePartitionedMeshReader::readProperties: "
252 "Could not read the meta data for the PropertyVector {:d}",
256 for (std::size_t i(0); i < number_of_properties; ++i)
260 auto pos = is.tellg();
262 static_cast<long>(pos) +
266 std::optional<MeshLib::IO::PropertyVectorPartitionMetaData> pvpmd(
268 bool const all_pvpmd_read_ok =
270 if (!all_pvpmd_read_ok)
273 "Error in NodePartitionedMeshReader::readProperties: "
274 "Could not read the partition meta data for the mpi process {:d}",
277 DBUG(
"offset in the PropertyVector: {:d}", pvpmd->offset);
278 DBUG(
"{:d} tuples in partition.", pvpmd->number_of_tuples);
281 const std::string fname_val = file_name_base +
"_partitioned_" + item_type +
284 is.open(fname_val.c_str(), std::ios::binary | std::ios::in);
287 ERR(
"Could not open file '{:s}'\n."
288 "\tYou can ignore this warning if the mesh does not contain {:s}-"
289 "wise property data.",
290 fname_val, item_type.data());
297 std::vector<std::optional<MeshLib::IO::PropertyVectorMetaData>>
const&
304 unsigned long global_offset = 0;
305 std::size_t
const number_of_properties = vec_pvmd.size();
306 for (std::size_t i(0); i < number_of_properties; ++i)
308 DBUG(
"global offset: {:d}, offset within the PropertyVector: {:d}.",
310 global_offset + pvpmd.
offset * vec_pvmd[i]->number_of_components *
311 vec_pvmd[i]->data_type_size_in_bytes);
316 if (vec_pvmd[i]->property_name.find(
"_ip") == std::string::npos &&
322 global_offset += vec_pvmd[i]->data_type_size_in_bytes *
323 vec_pvmd[i]->number_of_tuples;
327 if (vec_pvmd[i]->is_int_type)
329 if (vec_pvmd[i]->is_data_type_signed)
331 if (vec_pvmd[i]->data_type_size_in_bytes ==
sizeof(
char))
336 else if (vec_pvmd[i]->data_type_size_in_bytes ==
sizeof(
int))
341 else if (vec_pvmd[i]->data_type_size_in_bytes ==
sizeof(
long))
347 "Implementation for reading signed integer property "
348 "vector '{:s}' is not available.",
349 vec_pvmd[i]->property_name);
354 if (vec_pvmd[i]->data_type_size_in_bytes ==
355 sizeof(
unsigned char))
358 is, *vec_pvmd[i], pvpmd, t, global_offset, p);
360 else if (vec_pvmd[i]->data_type_size_in_bytes ==
361 sizeof(
unsigned int))
363 is, *vec_pvmd[i], pvpmd, t, global_offset, p);
364 else if (vec_pvmd[i]->data_type_size_in_bytes ==
365 sizeof(
unsigned long))
367 is, *vec_pvmd[i], pvpmd, t, global_offset, p);
371 "Implementation for reading unsigned property vector "
372 "'{:s}' is not available.",
373 vec_pvmd[i]->property_name);
379 if (vec_pvmd[i]->data_type_size_in_bytes ==
sizeof(
float))
382 else if (vec_pvmd[i]->data_type_size_in_bytes ==
sizeof(
double))
388 "Implementation for reading floating point property vector "
389 "'{:s}' is not available.",
390 vec_pvmd[i]->property_name);
393 global_offset += vec_pvmd[i]->data_type_size_in_bytes *
394 vec_pvmd[i]->number_of_tuples *
395 vec_pvmd[i]->number_of_components;
400 std::string
const& mesh_name,
401 std::vector<MeshLib::Node*>
const& mesh_nodes,
402 std::vector<unsigned long>
const& glb_node_ids,
403 std::vector<MeshLib::Element*>
const& mesh_elems,
406 std::vector<std::size_t>
const gathered_n_regular_base_nodes =
409 std::vector<std::size_t> n_regular_base_nodes_at_rank =
412 std::size_t
const n_regular_high_order_nodes =
415 std::vector<std::size_t>
const gathered_n_regular_high_order_nodes =
418 std::vector<std::size_t> n_regular_high_order_nodes_at_rank =
422 mesh_name, mesh_nodes, glb_node_ids, mesh_elems, properties,
425 std::move(n_regular_base_nodes_at_rank),
426 std::move(n_regular_high_order_nodes_at_rank));
430 const std::vector<NodeData>& node_data,
431 std::vector<MeshLib::Node*>& mesh_node,
432 std::vector<unsigned long>& glb_node_ids)
const
437 for (std::size_t i = 0; i < mesh_node.size(); i++)
440 glb_node_ids[i] = nd.
index;
446 const std::vector<MeshLib::Node*>& mesh_nodes,
447 const std::vector<unsigned long>& elem_data,
448 std::vector<MeshLib::Element*>& mesh_elems,
const bool ghost)
const
453 unsigned long const id_offset_ghost =
456 for (
unsigned long i = 0; i < ne; i++)
458 unsigned long id_offset_elem = elem_data[i];
462 const unsigned mat_idx =
463 static_cast<unsigned>(elem_data[id_offset_elem++]);
466 const unsigned long e_type = elem_data[id_offset_elem++];
467 unsigned long const nnodes = elem_data[id_offset_elem++];
470 for (
unsigned long k = 0; k < nnodes; k++)
471 elem_nodes[k] = mesh_nodes[elem_data[id_offset_elem++]];
474 switch (
static_cast<CellType>(e_type))
477 mesh_elems[i + id_offset_ghost] =
481 mesh_elems[i + id_offset_ghost] =
new MeshLib::Line(elem_nodes);
484 mesh_elems[i + id_offset_ghost] =
488 mesh_elems[i + id_offset_ghost] =
new MeshLib::Quad(elem_nodes);
491 mesh_elems[i + id_offset_ghost] =
495 mesh_elems[i + id_offset_ghost] =
499 mesh_elems[i + id_offset_ghost] =
new MeshLib::Hex(elem_nodes);
502 mesh_elems[i + id_offset_ghost] =
507 "NodePartitionedMeshReader: construction of HEX27 element "
508 "with id {:d} is not implemented.",
512 mesh_elems[i + id_offset_ghost] =
new MeshLib::Tri(elem_nodes);
515 mesh_elems[i + id_offset_ghost] =
new MeshLib::Tri6(elem_nodes);
518 mesh_elems[i + id_offset_ghost] =
new MeshLib::Tet(elem_nodes);
521 mesh_elems[i + id_offset_ghost] =
525 mesh_elems[i + id_offset_ghost] =
529 mesh_elems[i + id_offset_ghost] =
533 mesh_elems[i + id_offset_ghost] =
537 mesh_elems[i + id_offset_ghost] =
542 "NodePartitionedMeshReader: construction of INVALID "
543 "element type with id {:d} is not possible.",
548 "NodePartitionedMeshReader: construction of element type "
549 "{:d} is not implemented.",
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition of mesh-related Enumerations.
Definition of the class Properties that implements a container of properties.
bool is_safely_convertable(VALUE const &value)
Declare a class to read node-wise partitioned mesh with MPI functions.
Definition of the RunTime class.
double elapsed() const
Get the elapsed time in seconds.
void start()
Start the timer.
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
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
MeshLib::NodePartitionedMesh * readMesh(const std::string &file_name_base)
Create a NodePartitionedMesh object, read binary mesh data in the manner of parallel,...
MeshLib::NodePartitionedMesh * 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 b...
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.
NodePartitionedMeshReader(MPI_Comm comm)
void createSpecificPropertyVectorPart(std::istream &is, MeshLib::IO::PropertyVectorMetaData const &pvmd, MeshLib::MeshItemType t, unsigned long global_offset, MeshLib::Properties &p) const
MPI_Datatype _mpi_node_type
MPI data type for struct NodeData.
void registerNodeDataMpiType()
Define MPI data type for NodeData struct.
~NodePartitionedMeshReader()
struct MeshLib::IO::NodePartitionedMeshReader::PartitionedMeshInfo _mesh_info
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.
BaseLib::MPI::Mpi mpi_
Pointer to MPI communicator, the rank and the size.
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
Property manager on mesh items. Class Properties manages scalar, vector or matrix properties....
static T allreduce(T const &value, MPI_Op const &mpi_op, Mpi const &mpi)
static std::vector< T > allgather(T const &value, Mpi const &mpi)
bool IsFileExisting(const std::string &strFilename)
Returns true if given file exists.
std::vector< ranges::range_value_t< R > > sizesToOffsets(R const &sizes)
std::string extractBaseName(std::string const &pathname)
void writePropertyVectorMetaData(std::ostream &os, PropertyVectorMetaData const &pvmd)
std::optional< PropertyVectorMetaData > readPropertyVectorMetaData(std::istream &is)
std::optional< PropertyVectorPartitionMetaData > readPropertyVectorPartitionMetaData(std::istream &is)
TemplateElement< MeshLib::QuadRule9 > Quad9
TemplateElement< MeshLib::TetRule10 > Tet10
CellType
Types of mesh elements supported by OpenGeoSys.
TemplateElement< MeshLib::HexRule20 > Hex20
TemplateElement< MeshLib::TetRule4 > Tet
TemplateElement< MeshLib::LineRule2 > Line
TemplateElement< MeshLib::QuadRule8 > Quad8
TemplateElement< MeshLib::PyramidRule13 > Pyramid13
TemplateElement< MeshLib::QuadRule4 > Quad
TemplateElement< MeshLib::TriRule3 > Tri
TemplateElement< MeshLib::PyramidRule5 > Pyramid
TemplateElement< MeshLib::PrismRule6 > Prism
TemplateElement< PointRule1 > Point
TemplateElement< MeshLib::PrismRule15 > Prism15
TemplateElement< MeshLib::TriRule6 > Tri6
static constexpr char const * toString(const MeshItemType t)
Returns a char array for a specific MeshItemType.
TemplateElement< MeshLib::LineRule3 > Line3
TemplateElement< MeshLib::HexRule8 > Hex
struct NodeData used for parallel reading and also partitioning
std::size_t index
Global node index.
unsigned long number_of_regular_base_nodes
unsigned long number_of_regular_nodes
unsigned long number_of_ghost_elements
unsigned long number_of_global_base_nodes
unsigned long number_of_regular_elements
unsigned long number_of_nodes
0: Number of all nodes of a partition,
unsigned long number_of_global_nodes
7: Number of all nodes of global mesh,