31template <
typename VALUE,
typename TYPE>
34 bool const result = value <= std::numeric_limits<TYPE>::max();
37 ERR(
"The value {:d} is too large for conversion.", value);
38 ERR(
"Maximum available size is {:d}.",
39 std::numeric_limits<TYPE>::max());
65 int blocks[count] = {1, 3};
66 MPI_Datatype types[count] = {MPI_UNSIGNED_LONG, MPI_DOUBLE};
69 MPI_Type_create_struct(count, blocks, displacements, types,
75 const std::string& file_name_base)
81 std::string
const fname_new = file_name_base +
"_partitioned_msh_cfg" +
86 std::string
const fname_ascii = file_name_base +
87 "_partitioned_msh_cfg" +
91 ERR(
"Reading of ASCII mesh file {:s} is not supported since OGS "
95 OGS_FATAL(
"Required binary file {:s} does not exist.\n", fname_new);
98 INFO(
"Reading corresponding part of mesh data from binary file {:s} ...",
103 INFO(
"[time] Reading the mesh took {:f} s.", timer.
elapsed());
110template <
typename DATA>
119 ERR(
"The container size is too large for MPI_File_read() call.");
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);
130 if (file_status != 0)
132 ERR(
"Error opening file {:s}. MPI error code {:d}", filename,
138 char file_mode[] =
"native";
139 MPI_File_set_view(file, offset, type, type, file_mode, MPI_INFO_NULL);
141 MPI_File_read(file, data.data(),
static_cast<int>(data.size()), type,
143 MPI_File_close(&file);
149 const std::string& file_name_base)
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";
159 fname_header +
"cfg" + fname_num_p_ext,
160 static_cast<MPI_Offset
>(
static_cast<unsigned>(
_mpi_rank) *
174 std::vector<MeshLib::Node*> mesh_nodes;
175 std::vector<unsigned long> glb_node_ids;
176 setNodes(nodes, mesh_nodes, glb_node_ids);
184 MPI_LONG, elem_data))
187 std::vector<MeshLib::Element*> mesh_elems(
194 std::vector<unsigned long> ghost_elem_data(
199 MPI_LONG, ghost_elem_data))
202 const bool process_ghost =
true;
203 setElements(mesh_nodes, ghost_elem_data, mesh_elems, process_ghost);
210 glb_node_ids, mesh_elems, p);
214 const std::string& file_name_base)
const
229 const std::string fname_cfg = file_name_base +
"_partitioned_" + item_type +
232 std::ifstream is(fname_cfg.c_str(), std::ios::binary | std::ios::in);
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());
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)
253 "Error in NodePartitionedMeshReader::readProperties: "
254 "Could not read the meta data for the PropertyVector {:d}",
258 for (std::size_t i(0); i < number_of_properties; ++i)
262 auto pos = is.tellg();
264 static_cast<long>(pos) +
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,
274 if (!all_pvpmd_read_ok)
277 "Error in NodePartitionedMeshReader::readProperties: "
278 "Could not read the partition meta data for the mpi process {:d}",
281 DBUG(
"offset in the PropertyVector: {:d}", pvpmd->offset);
282 DBUG(
"{:d} tuples in partition.", pvpmd->number_of_tuples);
285 const std::string fname_val = file_name_base +
"_partitioned_" + item_type +
288 is.open(fname_val.c_str(), std::ios::binary | std::ios::in);
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());
301 std::vector<std::optional<MeshLib::IO::PropertyVectorMetaData>>
const&
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)
312 DBUG(
"global offset: {:d}, offset within the PropertyVector: {:d}.",
314 global_offset + pvpmd.
offset * vec_pvmd[i]->number_of_components *
315 vec_pvmd[i]->data_type_size_in_bytes);
320 if (vec_pvmd[i]->property_name.find(
"_ip") == std::string::npos &&
326 global_offset += vec_pvmd[i]->data_type_size_in_bytes *
327 vec_pvmd[i]->number_of_tuples;
331 if (vec_pvmd[i]->is_int_type)
333 if (vec_pvmd[i]->is_data_type_signed)
335 if (vec_pvmd[i]->data_type_size_in_bytes ==
sizeof(
char))
340 else if (vec_pvmd[i]->data_type_size_in_bytes ==
sizeof(
int))
345 else if (vec_pvmd[i]->data_type_size_in_bytes ==
sizeof(
long))
351 "Implementation for reading signed integer property "
352 "vector '{:s}' is not available.",
353 vec_pvmd[i]->property_name);
358 if (vec_pvmd[i]->data_type_size_in_bytes ==
359 sizeof(
unsigned char))
362 is, *vec_pvmd[i], pvpmd, t, global_offset, p);
364 else if (vec_pvmd[i]->data_type_size_in_bytes ==
365 sizeof(
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))
371 is, *vec_pvmd[i], pvpmd, t, global_offset, p);
375 "Implementation for reading unsigned property vector "
376 "'{:s}' is not available.",
377 vec_pvmd[i]->property_name);
383 if (vec_pvmd[i]->data_type_size_in_bytes ==
sizeof(
float))
386 else if (vec_pvmd[i]->data_type_size_in_bytes ==
sizeof(
double))
392 "Implementation for reading floating point property vector "
393 "'{:s}' is not available.",
394 vec_pvmd[i]->property_name);
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;
404 std::string
const& mesh_name,
405 std::vector<MeshLib::Node*>
const& mesh_nodes,
406 std::vector<unsigned long>
const& glb_node_ids,
407 std::vector<MeshLib::Element*>
const& mesh_elems,
410 std::vector<std::size_t> gathered_n_regular_base_nodes(
_mpi_comm_size);
415 gathered_n_regular_base_nodes.data(),
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));
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,
434 gathered_n_regular_high_order_nodes.data(),
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));
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));
454 const std::vector<NodeData>& node_data,
455 std::vector<MeshLib::Node*>& mesh_node,
456 std::vector<unsigned long>& glb_node_ids)
const
461 for (std::size_t i = 0; i < mesh_node.size(); i++)
464 glb_node_ids[i] = nd.
index;
470 const std::vector<MeshLib::Node*>& mesh_nodes,
471 const std::vector<unsigned long>& elem_data,
472 std::vector<MeshLib::Element*>& mesh_elems,
const bool ghost)
const
477 unsigned long const id_offset_ghost =
480 for (
unsigned long i = 0; i < ne; i++)
482 unsigned long id_offset_elem = elem_data[i];
486 const unsigned mat_idx =
487 static_cast<unsigned>(elem_data[id_offset_elem++]);
490 const unsigned long e_type = elem_data[id_offset_elem++];
491 unsigned long const nnodes = elem_data[id_offset_elem++];
494 for (
unsigned long k = 0; k < nnodes; k++)
495 elem_nodes[k] = mesh_nodes[elem_data[id_offset_elem++]];
498 switch (
static_cast<CellType>(e_type))
501 mesh_elems[i + id_offset_ghost] =
505 mesh_elems[i + id_offset_ghost] =
new MeshLib::Line(elem_nodes);
508 mesh_elems[i + id_offset_ghost] =
512 mesh_elems[i + id_offset_ghost] =
new MeshLib::Quad(elem_nodes);
515 mesh_elems[i + id_offset_ghost] =
519 mesh_elems[i + id_offset_ghost] =
523 mesh_elems[i + id_offset_ghost] =
new MeshLib::Hex(elem_nodes);
526 mesh_elems[i + id_offset_ghost] =
531 "NodePartitionedMeshReader: construction of HEX27 element "
532 "with id {:d} is not implemented.",
536 mesh_elems[i + id_offset_ghost] =
new MeshLib::Tri(elem_nodes);
539 mesh_elems[i + id_offset_ghost] =
new MeshLib::Tri6(elem_nodes);
542 mesh_elems[i + id_offset_ghost] =
new MeshLib::Tet(elem_nodes);
545 mesh_elems[i + id_offset_ghost] =
549 mesh_elems[i + id_offset_ghost] =
553 mesh_elems[i + id_offset_ghost] =
557 mesh_elems[i + id_offset_ghost] =
561 mesh_elems[i + id_offset_ghost] =
566 "NodePartitionedMeshReader: construction of INVALID "
567 "element type with id {:d} is not possible.",
572 "NodePartitionedMeshReader: construction of element type "
573 "{: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.
int _mpi_comm_size
Number of processes in the communicator: _mpi_comm.
NodePartitionedMeshReader(MPI_Comm comm)
MPI_Comm _mpi_comm
Pointer to MPI communicator.
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
int _mpi_rank
Rank of compute core.
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.
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....
bool IsFileExisting(const std::string &strFilename)
Returns true if given file exists.
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,