Loading [MathJax]/extensions/tex2jax.js
OGS
MeshLib::IO Namespace Reference

Namespaces

namespace  Legacy
 

Classes

struct  FileCommunicator
 
struct  HdfData
 
class  HdfWriter
 
struct  MeshHdfData
 
struct  NodeData
 struct NodeData used for parallel reading and also partitioning More...
 
class  NodePartitionedMeshReader
 
struct  PartitionInfo
 
struct  PropertyVectorMetaData
 
struct  PropertyVectorPartitionMetaData
 
class  PVDFile
 
struct  TransformedMeshData
 
class  VtuInterface
 Reads and writes VtkXMLUnstructuredGrid-files (vtu) to and from OGS data structures. This class is currently not inherited from Writer because VTK will implement writing to a string from 6.2 onwards. More...
 
struct  XdmfData
 
struct  XdmfHdfData
 
struct  XdmfHdfMesh
 
class  XdmfHdfWriter
 
struct  XdmfTopology
 
class  XdmfWriter
 

Typedefs

using Hdf5DimType = hsize_t
 
using HDFAttributes = std::vector<HdfData>
 

Functions

void writePropertyVectorMetaData (std::ostream &os, PropertyVectorMetaData const &pvmd)
 
void writePropertyVectorMetaData (PropertyVectorMetaData const &pvmd)
 
std::optional< PropertyVectorMetaDatareadPropertyVectorMetaData (std::istream &is)
 
void writePropertyVectorPartitionMetaData (std::ostream &os, PropertyVectorPartitionMetaData const &pvpmd)
 
std::optional< PropertyVectorPartitionMetaDatareadPropertyVectorPartitionMetaData (std::istream &is)
 
MeshLib::MeshreadMeshFromFile (const std::string &file_name, bool const compute_element_neighbors)
 
std::string getVtuFileNameForPetscOutputWithoutExtension (std::string const &file_name)
 
int writeVtu (MeshLib::Mesh const &mesh, std::string const &file_name, int const data_mode)
 
int writeMeshToFile (const MeshLib::Mesh &mesh, std::filesystem::path const &file_path, std::set< std::string > variable_output_names)
 
int64_t createFile (std::filesystem::path const &filepath, unsigned int n_files)
 
int64_t openHDF5File (std::filesystem::path const &filepath, unsigned int n_files)
 
int64_t createHDF5TransferPolicy ()
 
static hid_t meshPropertyType2HdfType (MeshPropertyDataType const ogs_data_type)
 
std::filesystem::path partitionFilename (std::filesystem::path const &basic_filepath, int const file_group)
 
int getGroupIndex (int const input_index, int const input_size, int const new_group_size)
 
FileCommunicator getCommunicator (unsigned int const n_files)
 
bool isFileManager ()
 
PartitionInfo getPartitionInfo (std::size_t const size, unsigned int const n_files)
 
int64_t openHDF5File (std::filesystem::path const &filepath)
 
static constexpr auto elemOGSTypeToXDMFType ()
 
constexpr auto cellTypeOGS2XDMF (MeshLib::CellType const &cell_type)
 
std::optional< XdmfHdfDatatransformAttribute (std::pair< std::string, PropertyVectorBase * > const &property_pair, unsigned int const n_files, unsigned int const chunk_size_bytes)
 
std::vector< XdmfHdfDatatransformAttributes (MeshLib::Mesh const &mesh, unsigned int n_files, unsigned int chunk_size_bytes)
 Create meta data for attributes used for hdf5 and xdmf.
 
std::vector< double > transformToXDMFGeometry (MeshLib::Mesh const &mesh)
 Copies all node points into a new vector. Contiguous data used for writing. Conform with XDMF standard in https://xdmf.org/index.php/XDMF_Model_and_Format.
 
XdmfHdfData transformGeometry (MeshLib::Mesh const &mesh, double const *data_ptr, unsigned int n_files, unsigned int chunk_size_bytes)
 Create meta data for geometry used for hdf5 and xdmf.
 
ParentDataType getLocalTopologyType (MeshLib::Mesh const &mesh)
 
ParentDataType getTopologyType (MeshLib::Mesh const &mesh)
 
std::pair< std::vector< std::size_t >, ParentDataTypetransformToXDMFTopology (MeshLib::Mesh const &mesh, std::size_t const offset)
 Copies all cells into a new vector. Contiguous data used for writing. The topology is specific to xdmf because it contains the xdmf cell types!! See section topology in https://xdmf.org/index.php/XDMF_Model_and_Format.
 
XdmfHdfData transformTopology (std::vector< std::size_t > const &values, ParentDataType const parent_data_type, unsigned int n_files, unsigned int chunk_size_bytes)
 Create meta data for topology used for HDF5 and XDMF.
 
static std::string meshItemTypeString (std::optional< MeshItemType > const &item_type)
 
static auto meshPropertyDatatypeString ()
 
static std::string getPropertyDataTypeString (MeshPropertyDataType const &ogs_data_type)
 
static constexpr auto meshPropertyDatatypeSize ()
 
static std::string getPropertyDataTypeSize (MeshPropertyDataType const &ogs_data_type)
 
std::function< std::string(std::vector< double >)> write_xdmf (XdmfData const &geometry, XdmfData const &topology, std::vector< XdmfData > const &constant_attributes, std::vector< XdmfData > const &variable_attributes, std::string const &h5filename, std::string const &ogs_version, std::string const &mesh_name)
 Generator function that creates a function capturing the spatial data of a mesh Temporal data can later be passed as argument.
 
template<typename Data >
std::function< bool(Data)> isVariableAttribute (std::set< std::string > const &variable_output_names)
 

Variables

constexpr auto elem_type_ogs2xdmf = elemOGSTypeToXDMFType()
 
constexpr char const * mesh_item_type_strings []
 
auto ogs_to_xdmf_type_fn = meshPropertyDatatypeSize()
 

Typedef Documentation

◆ Hdf5DimType

using MeshLib::IO::Hdf5DimType = hsize_t

Definition at line 24 of file HdfData.h.

◆ HDFAttributes

using MeshLib::IO::HDFAttributes = std::vector<HdfData>

Definition at line 25 of file HdfWriter.h.

Function Documentation

◆ cellTypeOGS2XDMF()

auto MeshLib::IO::cellTypeOGS2XDMF ( MeshLib::CellType const & cell_type)
constexpr

Definition at line 87 of file transformData.cpp.

88{
89 return elem_type_ogs2xdmf[to_underlying(cell_type)];
90}
constexpr auto to_underlying(E e) noexcept
Converts an enumeration to its underlying type.
Definition cpp23.h:29
constexpr auto elem_type_ogs2xdmf

References elem_type_ogs2xdmf, and BaseLib::to_underlying().

Referenced by getLocalTopologyType(), and transformToXDMFTopology().

◆ createFile()

int64_t MeshLib::IO::createFile ( std::filesystem::path const & filepath,
unsigned int n_files )

Definition at line 38 of file fileIO.cpp.

40{
41 auto const communicator = getCommunicator(n_files);
42 MPI_Comm const comm = communicator.mpi_communicator;
43 MPI_Info const info = MPI_INFO_NULL;
44 hid_t const plist_id = H5Pcreate(H5P_FILE_ACCESS);
45
46 H5Pset_fapl_mpio(plist_id, comm, info);
47 H5Pset_coll_metadata_write(plist_id, true);
48
49 std::filesystem::path const partition_filename =
50 partitionFilename(filepath, communicator.color);
51 hid_t file = H5Fcreate(partition_filename.string().c_str(), H5F_ACC_TRUNC,
52 H5P_DEFAULT, plist_id);
53 H5Pclose(plist_id);
54
55 return file;
56}
std::filesystem::path partitionFilename(std::filesystem::path const &basic_filepath, int const file_group)
Definition fileIO.cpp:24
FileCommunicator getCommunicator(unsigned int const n_files)

References getCommunicator(), and partitionFilename().

◆ createHDF5TransferPolicy()

int64_t MeshLib::IO::createHDF5TransferPolicy ( )

Definition at line 70 of file fileIO.cpp.

71{
72 // property list for collective dataset write
73 hid_t io_transfer_property = H5Pcreate(H5P_DATASET_XFER);
74 H5Pset_dxpl_mpio(io_transfer_property, H5FD_MPIO_COLLECTIVE);
75 return io_transfer_property;
76}

Referenced by writeDataSet().

◆ elemOGSTypeToXDMFType()

static constexpr auto MeshLib::IO::elemOGSTypeToXDMFType ( )
staticconstexpr

Definition at line 42 of file transformData.cpp.

43{
45 elem_type{};
73 ParentDataType::WEDGE, 6}; // parallel triangle wedge
82 return elem_type;
83}

References EDGE_3, MeshLib::enum_length, MeshLib::HEX20, MeshLib::HEX27, MeshLib::HEX8, HEXAHEDRON, HEXAHEDRON_20, HEXAHEDRON_27, MeshLib::LINE2, MeshLib::LINE3, MeshLib::POINT1, POLYLINE, POLYVERTEX, MeshLib::PRISM15, MeshLib::PRISM18, MeshLib::PRISM6, PYRAMID, MeshLib::PYRAMID13, MeshLib::PYRAMID5, PYRAMID_13, MeshLib::QUAD4, MeshLib::QUAD8, MeshLib::QUAD9, QUADRILATERAL, QUADRILATERAL_8, QUADRILATERAL_9, MeshLib::TET10, MeshLib::TET4, TETRAHEDRON, TETRAHEDRON_10, BaseLib::to_underlying(), MeshLib::TRI3, MeshLib::TRI6, TRIANGLE, TRIANGLE_6, WEDGE, WEDGE_15, and WEDGE_18.

◆ getCommunicator()

FileCommunicator MeshLib::IO::getCommunicator ( unsigned int const n_files)

Definition at line 41 of file getCommunicator.cpp.

42{
43 int num_procs;
44 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
45 int rank_id;
46 MPI_Comm_rank(MPI_COMM_WORLD, &rank_id);
47 int const file_group_id = getGroupIndex(rank_id, num_procs, n_files);
48 MPI_Comm new_communicator;
49 MPI_Comm_split(MPI_COMM_WORLD, file_group_id, rank_id, &new_communicator);
50 return FileCommunicator{std::move(new_communicator),
51 std::move(file_group_id), ""};
52}
int getGroupIndex(int const input_index, int const input_size, int const new_group_size)

References getGroupIndex().

Referenced by createFile(), getPartitionInfo(), and openHDF5File().

◆ getGroupIndex()

int MeshLib::IO::getGroupIndex ( int const input_index,
int const input_size,
int const new_group_size )

Definition at line 27 of file getCommunicator.cpp.

29{
30 // A grouping algorithm that determines the number of groups and return the
31 // group idx of the specified input_index
32 assert(input_size >= new_group_size);
33 int const minimum_output_group_size =
34 std::lround(input_size / new_group_size);
35 int const maximum_output_group_size = (input_size % new_group_size)
36 ? minimum_output_group_size + 1
37 : minimum_output_group_size;
38 return std::lround(input_index / maximum_output_group_size);
39};

Referenced by getCommunicator().

◆ getLocalTopologyType()

ParentDataType MeshLib::IO::getLocalTopologyType ( MeshLib::Mesh const & mesh)

Definition at line 325 of file transformData.cpp.

326{
327 auto const& elements = mesh.getElements();
328
329 if (elements.empty())
330 {
332 }
333 auto const ogs_cell_type = elements[0]->getCellType();
334
335 if (std::any_of(elements.begin(), elements.end(),
336 [&](auto const& cell)
337 { return cell->getCellType() != ogs_cell_type; }))
338 {
340 }
341
342 return cellTypeOGS2XDMF(ogs_cell_type).id;
343}

References cellTypeOGS2XDMF(), MeshLib::Mesh::getElements(), and MIXED.

Referenced by getTopologyType().

◆ getPartitionInfo()

PartitionInfo MeshLib::IO::getPartitionInfo ( std::size_t const size,
unsigned int const n_files )

Definition at line 35 of file partition.cpp.

37{
39
40 std::vector<std::size_t> const partition_sizes =
41 BaseLib::MPI::allgather(size, mpi);
42
43 // the first partition's offset is zero, offsets for subsequent
44 // partitions are the accumulated sum of all preceding size.
45 std::vector<std::size_t> const partition_offsets =
46 BaseLib::sizesToOffsets(partition_sizes);
47
48 // chunked
49 std::size_t longest_partition =
50 *max_element(partition_sizes.begin(), partition_sizes.end());
51
52 // local_offset, local_length, longest_local_length, global_length
53 return {partition_offsets[mpi.rank], size, longest_partition,
54 partition_offsets.back()};
55}
static std::vector< T > allgather(T const &value, Mpi const &mpi)
Definition MPI.h:100
std::vector< ranges::range_value_t< R > > sizesToOffsets(R const &sizes)
Definition Algorithm.h:283

References BaseLib::MPI::allgather(), getCommunicator(), MeshLib::IO::FileCommunicator::mpi_communicator, and BaseLib::sizesToOffsets().

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

◆ getPropertyDataTypeSize()

static std::string MeshLib::IO::getPropertyDataTypeSize ( MeshPropertyDataType const & ogs_data_type)
static

Definition at line 99 of file writeXdmf.cpp.

101{
102 return fmt::format("{}", ogs_to_xdmf_type_fn[to_underlying(ogs_data_type)]);
103}

References ogs_to_xdmf_type_fn, and BaseLib::to_underlying().

Referenced by write_xdmf().

◆ getPropertyDataTypeString()

static std::string MeshLib::IO::getPropertyDataTypeString ( MeshPropertyDataType const & ogs_data_type)
static

Definition at line 74 of file writeXdmf.cpp.

76{
77 return meshPropertyDatatypeString()[to_underlying(ogs_data_type)];
78}
static auto meshPropertyDatatypeString()
Definition writeXdmf.cpp:56

References meshPropertyDatatypeString(), and BaseLib::to_underlying().

Referenced by write_xdmf().

◆ getTopologyType()

ParentDataType MeshLib::IO::getTopologyType ( MeshLib::Mesh const & mesh)

Definition at line 351 of file transformData.cpp.

352{
353 auto const local_topology_type = getLocalTopologyType(mesh);
354 // get the topology_type from other partitions
355 BaseLib::MPI::Mpi const mpi;
356
357 std::vector<int> const topology_types =
358 BaseLib::MPI::allgather(static_cast<int>(local_topology_type), mpi);
359
360 auto const common_topology_type =
361 std::all_of(topology_types.begin(), topology_types.end(),
362 [&local_topology_type](int other)
363 { return other == static_cast<int>(local_topology_type); })
364 ? local_topology_type
366 return common_topology_type;
367}
ParentDataType getLocalTopologyType(MeshLib::Mesh const &mesh)

References BaseLib::MPI::allgather(), getLocalTopologyType(), and MIXED.

Referenced by transformToXDMFTopology().

◆ getVtuFileNameForPetscOutputWithoutExtension()

std::string MeshLib::IO::getVtuFileNameForPetscOutputWithoutExtension ( std::string const & file_name)

Definition at line 122 of file VtuInterface.cpp.

124{
125 auto const file_name_extension = BaseLib::getFileExtension(file_name);
126 if (file_name_extension != ".vtu")
127 {
128 OGS_FATAL("Expected a .vtu file for petsc output.");
129 }
130
131 auto const file_name_base = boost::erase_last_copy(file_name, ".vtu");
132 auto basename = BaseLib::extractBaseName(file_name_base);
133
134 // Replace dots to underscores since the pvtu writing function drops all
135 // characters starting from a dot.
136 std::replace(basename.begin(), basename.end(), '.', '_');
137
138 // Restore the dirname if any.
139 auto const dirname = BaseLib::extractPath(file_name_base);
140 return BaseLib::joinPaths(dirname, basename);
141}
#define OGS_FATAL(...)
Definition Error.h:26
std::string getFileExtension(const std::string &path)
std::string extractPath(std::string const &pathname)
std::string joinPaths(std::string const &pathA, std::string const &pathB)
std::string extractBaseName(std::string const &pathname)

References BaseLib::extractBaseName(), BaseLib::extractPath(), BaseLib::getFileExtension(), BaseLib::joinPaths(), and OGS_FATAL.

Referenced by ApplicationsLib::TestDefinition::TestDefinition(), MeshLib::IO::PVDFile::addVTUFile(), and MeshLib::IO::VtuInterface::writeToFile().

◆ isFileManager()

bool MeshLib::IO::isFileManager ( )

Definition at line 28 of file partition.cpp.

29{
30 int mpi_rank;
31 MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
32 return mpi_rank == 0;
33}

Referenced by MeshLib::IO::XdmfHdfWriter::XdmfHdfWriter(), and MeshLib::IO::XdmfHdfWriter::writeStep().

◆ isVariableAttribute()

template<typename Data >
std::function< bool(Data)> MeshLib::IO::isVariableAttribute ( std::set< std::string > const & variable_output_names)

Definition at line 43 of file XdmfHdfWriter.cpp.

45{
46 if (variable_output_names.empty())
47 {
48 return [](Data const& data) -> bool
49 {
50 constexpr std::array constant_output_names{
51 "MaterialIDs"sv,
52 "topology"sv,
53 "geometry"sv,
54 "OGS_VERSION"sv,
59 return !ranges::contains(constant_output_names, data.name);
60 };
61 }
62 return [&variable_output_names](Data const& data) -> bool
63 { return variable_output_names.contains(data.name); };
64}
constexpr std::string_view getBulkIDString(MeshItemType mesh_item_type)
Definition Properties.h:185

References MeshLib::Cell, MeshLib::Edge, MeshLib::Face, MeshLib::getBulkIDString(), and MeshLib::Node.

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

◆ meshItemTypeString()

static std::string MeshLib::IO::meshItemTypeString ( std::optional< MeshItemType > const & item_type)
static

Definition at line 44 of file writeXdmf.cpp.

46{
47 if (item_type)
48 {
49 return mesh_item_type_strings[static_cast<int>(item_type.value())];
50 }
51 OGS_FATAL("Cannot convert an empty optional mesh item type.");
52}
constexpr char const * mesh_item_type_strings[]
Definition writeXdmf.cpp:40

References mesh_item_type_strings, and OGS_FATAL.

Referenced by write_xdmf().

◆ meshPropertyDatatypeSize()

static constexpr auto MeshLib::IO::meshPropertyDatatypeSize ( )
staticconstexpr

◆ meshPropertyDatatypeString()

static auto MeshLib::IO::meshPropertyDatatypeString ( )
static

Definition at line 56 of file writeXdmf.cpp.

57{
58 std::array<std::string, to_underlying(MeshPropertyDataType::enum_length)>
59 ogs_to_xdmf_type = {};
60 ogs_to_xdmf_type[to_underlying(MeshPropertyDataType::float64)] = "Float";
61 ogs_to_xdmf_type[to_underlying(MeshPropertyDataType::float32)] = "Float";
62 ogs_to_xdmf_type[to_underlying(MeshPropertyDataType::int32)] = "Int";
63 ogs_to_xdmf_type[to_underlying(MeshPropertyDataType::int64)] = "Int";
64 ogs_to_xdmf_type[to_underlying(MeshPropertyDataType::uint32)] = "UInt";
65 ogs_to_xdmf_type[to_underlying(MeshPropertyDataType::uint64)] = "UInt";
66 ogs_to_xdmf_type[to_underlying(MeshPropertyDataType::int8)] = "Int";
67 ogs_to_xdmf_type[to_underlying(MeshPropertyDataType::uint8)] = "UInt";
68 ogs_to_xdmf_type[to_underlying(MeshPropertyDataType::char_native)] = "Char";
69 ogs_to_xdmf_type[to_underlying(MeshPropertyDataType::uchar)] = "UChar";
70 return ogs_to_xdmf_type;
71}

References char_native, enum_length, float32, float64, int32, int64, int8, BaseLib::to_underlying(), uchar, uint32, uint64, and uint8.

Referenced by getPropertyDataTypeString().

◆ meshPropertyType2HdfType()

static hid_t MeshLib::IO::meshPropertyType2HdfType ( MeshPropertyDataType const ogs_data_type)
static

Definition at line 22 of file HdfData.cpp.

23{
24 std::map<MeshPropertyDataType const, hid_t> ogs_to_hdf_type = {
25 {MeshPropertyDataType::float64, H5T_NATIVE_DOUBLE},
26 {MeshPropertyDataType::float32, H5T_NATIVE_FLOAT},
27 {MeshPropertyDataType::int32, H5T_NATIVE_INT32},
28 {MeshPropertyDataType::int64, H5T_NATIVE_INT64},
29 {MeshPropertyDataType::uint32, H5T_NATIVE_UINT32},
30 {MeshPropertyDataType::uint64, H5T_NATIVE_UINT64},
31 {MeshPropertyDataType::int8, H5T_NATIVE_INT8},
32 {MeshPropertyDataType::uint8, H5T_NATIVE_UINT8},
33 {MeshPropertyDataType::char_native, H5T_NATIVE_CHAR},
34 {MeshPropertyDataType::uchar, H5T_NATIVE_UCHAR}};
35 try
36 {
37 return ogs_to_hdf_type.at(ogs_data_type);
38 }
39 catch (std::exception const& e)
40 {
41 OGS_FATAL("No known HDF5 type for OGS type. {:s}", e.what());
42 }
43}

References char_native, float32, float64, int32, int64, int8, OGS_FATAL, uchar, uint32, uint64, and uint8.

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

◆ openHDF5File() [1/2]

int64_t MeshLib::IO::openHDF5File ( std::filesystem::path const & filepath)

Definition at line 25 of file fileIO.cpp.

26{
27 return H5Fopen(filepath.string().c_str(), H5F_ACC_RDWR, H5P_DEFAULT);
28}

◆ openHDF5File() [2/2]

hid_t MeshLib::IO::openHDF5File ( std::filesystem::path const & filepath,
unsigned int n_files )

Definition at line 58 of file fileIO.cpp.

60{
61 MPI_Comm const comm = getCommunicator(n_files).mpi_communicator;
62 MPI_Info info = MPI_INFO_NULL;
63 hid_t const plist_id = H5Pcreate(H5P_FILE_ACCESS);
64 H5Pset_fapl_mpio(plist_id, comm, info);
65 hid_t file = H5Fopen(filepath.string().c_str(), H5F_ACC_RDWR, plist_id);
66 H5Pclose(plist_id);
67 return file;
68}

References getCommunicator(), and MeshLib::IO::FileCommunicator::mpi_communicator.

◆ partitionFilename()

std::filesystem::path MeshLib::IO::partitionFilename ( std::filesystem::path const & basic_filepath,
int const file_group )

Definition at line 24 of file fileIO.cpp.

26{
27 std::string const filename = (file_group > 0)
28 ? basic_filepath.stem().string() + '_' +
29 std::to_string(file_group) +
30 basic_filepath.extension().string()
31 : basic_filepath.filename().string();
32 std::filesystem::path const filepathwithextension =
33 basic_filepath.parent_path() / filename;
34 DBUG("HDF Filepath: {:s}.", filepathwithextension.string());
35 return filepathwithextension;
36};
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30

References DBUG().

Referenced by createFile().

◆ readMeshFromFile()

MeshLib::Mesh * MeshLib::IO::readMeshFromFile ( const std::string & file_name,
bool const compute_element_neighbors )

Definition at line 68 of file readMeshFromFile.cpp.

70{
71#ifdef USE_PETSC
73 if (mpi.size > 1)
74 {
76 const std::string file_name_base =
78 return read_pmesh.read(file_name_base);
79 }
80 if (mpi.size == 1)
81 {
82 std::unique_ptr<Mesh> mesh{
83 readMeshFromFileSerial(file_name, compute_element_neighbors)};
84
85 if (!mesh)
86 {
87 return nullptr;
88 }
89
90 return new MeshLib::NodePartitionedMesh(*mesh);
91 }
92 return nullptr;
93#endif
94 return readMeshFromFileSerial(file_name, compute_element_neighbors);
95}
std::string dropFileExtension(std::string const &filename)
MeshLib::Mesh * readMeshFromFileSerial(const std::string &file_name, bool const compute_element_neighbors)
MPI_Comm communicator
Definition MPI.h:60

References BaseLib::MPI::Mpi::communicator, BaseLib::dropFileExtension(), MeshLib::IO::NodePartitionedMeshReader::read(), and BaseLib::MPI::Mpi::size.

Referenced by createGeometries(), FileIO::createSurface(), MainWindow::loadFile(), main(), main(), MainWindow::mapGeometry(), anonymous_namespace{postLIE.cpp}::postVTU(), FileIO::XmlPrjInterface::readFile(), readMeshes(), and anonymous_namespace{ProjectData.cpp}::readSingleMesh().

◆ readPropertyVectorMetaData()

std::optional< PropertyVectorMetaData > MeshLib::IO::readPropertyVectorMetaData ( std::istream & is)
inline

Definition at line 85 of file PropertyVectorMetaData.h.

87{
88 // read the size of the name of the PropertyVector
89 std::string::size_type s = 0;
90 if (!is.read(reinterpret_cast<char*>(&s), sizeof(std::string::size_type)))
91 {
92 return std::optional<PropertyVectorMetaData>();
93 }
94
95 PropertyVectorMetaData pvmd;
96 char *dummy = new char[s];
97 if (!is.read(dummy, s))
98 {
99 return std::nullopt;
100 }
101 pvmd.property_name = std::string(dummy, s);
102 delete [] dummy;
103
104 if(!is.read(reinterpret_cast<char*>(&pvmd.is_int_type), sizeof(bool)))
105 return std::nullopt;
106 if(!is.read(reinterpret_cast<char*>(&pvmd.is_data_type_signed), sizeof(bool)))
107 return std::nullopt;
108 if(!is.read(reinterpret_cast<char*>(&pvmd.data_type_size_in_bytes),
109 sizeof(unsigned long)))
110 return std::nullopt;
111 if(!is.read(reinterpret_cast<char*>(&pvmd.number_of_components),
112 sizeof(unsigned long)))
113 return std::nullopt;
114 if(!is.read(reinterpret_cast<char*>(&pvmd.number_of_tuples),
115 sizeof(unsigned long)))
116 return std::nullopt;
117 return std::optional<PropertyVectorMetaData>(pvmd);
118}

References MeshLib::IO::PropertyVectorMetaData::data_type_size_in_bytes, MeshLib::IO::PropertyVectorMetaData::is_data_type_signed, MeshLib::IO::PropertyVectorMetaData::is_int_type, MeshLib::IO::PropertyVectorMetaData::number_of_components, MeshLib::IO::PropertyVectorMetaData::number_of_tuples, and MeshLib::IO::PropertyVectorMetaData::property_name.

Referenced by MeshLib::IO::NodePartitionedMeshReader::readProperties().

◆ readPropertyVectorPartitionMetaData()

std::optional< PropertyVectorPartitionMetaData > MeshLib::IO::readPropertyVectorPartitionMetaData ( std::istream & is)
inline

Definition at line 140 of file PropertyVectorMetaData.h.

141{
143 if (!is.read(reinterpret_cast<char*>(&pvpmd.offset), sizeof(unsigned long)))
144 {
145 return std::optional<PropertyVectorPartitionMetaData>();
146 }
147 if (!is.read(reinterpret_cast<char*>(&pvpmd.number_of_tuples),
148 sizeof(unsigned long)))
149 {
150 return std::optional<PropertyVectorPartitionMetaData>();
151 }
152 return std::optional<PropertyVectorPartitionMetaData>(pvpmd);
153}

References MeshLib::IO::PropertyVectorPartitionMetaData::number_of_tuples, and MeshLib::IO::PropertyVectorPartitionMetaData::offset.

Referenced by MeshLib::IO::NodePartitionedMeshReader::readProperties().

◆ transformAttribute()

std::optional< XdmfHdfData > MeshLib::IO::transformAttribute ( std::pair< std::string, PropertyVectorBase * > const & property_pair,
unsigned int const n_files,
unsigned int const chunk_size_bytes )

Definition at line 92 of file transformData.cpp.

95{
96 // 3 data that will be captured and written by lambda f below
98 std::size_t num_of_tuples = 0;
99 void const* data_ptr = 0;
100
101 // lambda f : Collects properties from the PropertyVectorBase. It captures
102 // (and overwrites) data that can only be collected via the typed property.
103 // It has boolean return type to allow kind of pipe using || operator.
104 auto f = [&data_type, &num_of_tuples, &data_ptr,
105 &property_pair](auto basic_type) -> bool
106 {
107 auto const property_base = property_pair.second;
108 auto const typed_property =
109 dynamic_cast<PropertyVector<decltype(basic_type)> const*>(
110 property_base);
111 if (typed_property == nullptr)
112 {
113 return false;
114 }
115 // overwrite captured data
116 num_of_tuples = typed_property->getNumberOfTuples();
117 data_ptr = typed_property->data();
118
119 if constexpr (std::is_same_v<double, decltype(basic_type)>)
120 {
121 // The standard 64-bit IEEE 754 floating-point type
122 // (double-precision) has a 53 bit fractional part (52 bits written,
123 // one implied)
124 static_assert((std::numeric_limits<double>::digits == 53),
125 "Double has 52 bits fractional part");
127 }
128 else if constexpr (std::is_same_v<float, decltype(basic_type)>)
129 {
130 // The standard 32-bit IEEE 754 floating-point type
131 // (single-precision) has a 24 bit fractional part (23 bits written,
132 // one implied)
133 static_assert((std::numeric_limits<float>::digits == 24),
134 "Float has 23 bits fractional part");
136 }
137 else if constexpr (std::is_same_v<int32_t, decltype(basic_type)>)
138 {
139 data_type = MeshPropertyDataType::int32;
140 }
141 else if constexpr (std::is_same_v<uint32_t, decltype(basic_type)>)
142 {
144 }
145 else if constexpr (std::is_same_v<int64_t, decltype(basic_type)>)
146 {
147 data_type = MeshPropertyDataType::int64;
148 }
149 else if constexpr (std::is_same_v<uint64_t, decltype(basic_type)>)
150 {
152 }
153 else if constexpr (std::is_same_v<int8_t, decltype(basic_type)>)
154 {
155 data_type = MeshPropertyDataType::int8;
156 }
157 else if constexpr (std::is_same_v<uint8_t, decltype(basic_type)>)
158 {
159 data_type = MeshPropertyDataType::uint8;
160 }
161 else if constexpr (std::is_same_v<char, decltype(basic_type)>)
162 {
164 }
165 else if constexpr (std::is_same_v<unsigned char, decltype(basic_type)>)
166 {
167 static_assert((std::numeric_limits<unsigned char>::digits == 8),
168 "Unsigned char has 8 bits");
169 data_type = MeshPropertyDataType::uchar;
170 }
171 else if constexpr (std::is_same_v<unsigned long, decltype(basic_type)>)
172 {
173 if (sizeof(unsigned long) == 8 &&
174 std::numeric_limits<unsigned long>::digits == 64)
175 {
177 }
178 else if (sizeof(unsigned long) == 4 &&
179 std::numeric_limits<unsigned long>::digits == 32)
180 {
182 }
183 else
184 {
185 return false;
186 }
187 }
188 else if constexpr (std::is_same_v<std::size_t, decltype(basic_type)>)
189 {
190 if (sizeof(std::size_t) == 8 &&
191 std::numeric_limits<std::size_t>::digits == 64)
192 {
194 }
195 else if (sizeof(std::size_t) == 4 &&
196 std::numeric_limits<std::size_t>::digits == 32)
197 {
199 }
200 else
201 {
202 return false;
203 }
204 }
205 else
206 {
207 return false;
208 }
209 return true;
210 };
211
212 f(double{}) || f(float{}) || f(int{}) || f(long{}) || f(unsigned{}) ||
213 f(long{}) || f(static_cast<unsigned long>(0)) || f(std::size_t{}) ||
214 f(char{}) || f(static_cast<unsigned char>(0));
215
216 if (data_type == MeshPropertyDataType::unknown)
217 {
218 return std::nullopt;
219 }
220
221 auto const& property_base = property_pair.second;
222 auto const& global_components =
223 property_base->getNumberOfGlobalComponents();
224 // TODO (tm) property_pair vector::getNumberOfGlobalComponents should return
225 // unsigned value. Then explicit cast from signed to unsigned int and
226 // assert can be removed here. Implicit cast to long long is fine and
227 // can be kept
228 assert(global_components >= 0);
229 auto const ui_global_components =
230 static_cast<unsigned int>(global_components);
231
232 MeshLib::MeshItemType const mesh_item_type =
233 property_base->getMeshItemType();
234
235 std::string const& name = property_base->getPropertyName();
236
237 HdfData hdf = {data_ptr, num_of_tuples, ui_global_components, name,
238 data_type, n_files, chunk_size_bytes};
239
240 XdmfData xdmf{num_of_tuples, ui_global_components, data_type,
241 name, mesh_item_type, 0,
242 n_files, std::nullopt};
243
244 return XdmfHdfData{std::move(hdf), std::move(xdmf)};
245}
MeshPropertyDataType
std::size_t getNumberOfTuples() const

References char_native, float32, float64, MeshLib::PropertyVector< T >::getNumberOfTuples(), int32, int64, int8, uchar, uint32, uint64, uint8, and unknown.

Referenced by transformAttributes().

◆ transformAttributes()

std::vector< XdmfHdfData > MeshLib::IO::transformAttributes ( MeshLib::Mesh const & mesh,
unsigned int n_files,
unsigned int chunk_size_bytes )

Create meta data for attributes used for hdf5 and xdmf.

Parameters
meshOGS mesh can be mesh or partitionedMesh
n_filesspecifies the number of files. If greater than 1 it groups the data of each process to n_files
chunk_size_bytesData will be split into chunks. The parameter specifies the size (in bytes) of the largest chunk.
Returns
vector of meta data

Definition at line 247 of file transformData.cpp.

250{
251 MeshLib::Properties const& properties = mesh.getProperties();
252
253 // \TODO (tm) use c++20 ranges
254 // a = p | filter (first!=OGS_VERSION) | filter null_opt | transformAttr |
255 std::vector<XdmfHdfData> attributes;
256 for (auto const& [name, property_base] : properties)
257 {
259 {
260 continue;
261 }
262
263 if (!property_base->is_for_output)
264 {
265 continue;
266 }
267
268 if (auto const attribute = transformAttribute(
269
270 std::pair(std::string(name), property_base), n_files,
271 chunk_size_bytes))
272
273 {
274 attributes.push_back(attribute.value());
275 }
276 else
277 {
278 WARN("Could not create attribute meta of {:s}.", name);
279 }
280 }
281 return attributes;
282}
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
Property manager on mesh items. Class Properties manages scalar, vector or matrix properties....
Definition Properties.h:33
const std::string OGS_VERSION
Definition GitInfo.cpp:20

References MeshLib::Mesh::getProperties(), GitInfoLib::GitInfo::OGS_VERSION, transformAttribute(), and WARN().

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

◆ transformGeometry()

XdmfHdfData MeshLib::IO::transformGeometry ( MeshLib::Mesh const & mesh,
double const * data_ptr,
unsigned int n_files,
unsigned int chunk_size_bytes )

Create meta data for geometry used for hdf5 and xdmf.

Parameters
meshOGS mesh can be mesh or partitionedMesh
data_ptrMemory location of geometry values.
n_filesspecifies the number of files. If greater than 1 it groups the data of each process to n_files
chunk_size_bytesData will be split into chunks. The parameter specifies the size (in bytes) of the largest chunk.
Returns
Geometry meta data

Definition at line 300 of file transformData.cpp.

303{
304 std::string const name = "geometry";
305 std::vector<MeshLib::Node*> const& nodes = mesh.getNodes();
306
307 int const point_size = 3;
308 auto const& partition_dim = nodes.size();
309
310 HdfData const hdf = {data_ptr,
311 partition_dim,
312 point_size,
313 name,
315 n_files,
316 chunk_size_bytes};
317 XdmfData const xdmf = {
318 partition_dim, point_size, MeshPropertyDataType::float64,
319 name, std::nullopt, 2,
320 n_files, std::nullopt};
321
322 return XdmfHdfData{std::move(hdf), std::move(xdmf)};
323}

References float64, and MeshLib::Mesh::getNodes().

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

◆ transformTopology()

XdmfHdfData MeshLib::IO::transformTopology ( std::vector< std::size_t > const & values,
ParentDataType const parent_data_type,
unsigned int n_files,
unsigned int chunk_size_bytes )

Create meta data for topology used for HDF5 and XDMF.

Parameters
valuesactual topology values to get size and memory location
parent_data_typeXDMF topological element data types as listed in the enum ParentDataTypei
n_filesspecifies the number of files. If greater than 1 it groups the data of each process to n_files
chunk_size_bytesData will be split into chunks. The parameter specifies the size (in bytes) of the largest chunk.
Returns
Topology meta data

Definition at line 415 of file transformData.cpp.

419{
420 std::string const name = "topology";
421 auto const tuple_size = ParentDataType2String(parent_data_type).second;
422 HdfData const hdf = {values.data(),
423 values.size() / tuple_size,
424 tuple_size,
425 name,
427 n_files,
428 chunk_size_bytes};
429 XdmfData const xdmf{values.size() / tuple_size,
430 tuple_size,
432 name,
433 std::nullopt,
434 3,
435 n_files,
436 parent_data_type};
437
438 return XdmfHdfData{std::move(hdf), std::move(xdmf)};
439}
std::pair< std::string, std::size_t > ParentDataType2String(ParentDataType p)

References ParentDataType2String(), and uint64.

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

◆ transformToXDMFGeometry()

std::vector< double > MeshLib::IO::transformToXDMFGeometry ( MeshLib::Mesh const & mesh)

Copies all node points into a new vector. Contiguous data used for writing. Conform with XDMF standard in https://xdmf.org/index.php/XDMF_Model_and_Format.

Parameters
meshOGS mesh can be mesh or partitionedMesh
Returns
vector containing a copy of the data

Definition at line 284 of file transformData.cpp.

285{
286 std::vector<MeshLib::Node*> const& nodes = mesh.getNodes();
287
288 int const point_size = 3;
289 std::vector<double> values;
290 values.reserve(nodes.size() * point_size);
291 for (auto const& n : nodes)
292 {
293 const double* x = n->data();
294 values.insert(values.cend(), x, x + point_size);
295 }
296
297 return values;
298}

References MeshLib::Mesh::getNodes().

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

◆ transformToXDMFTopology()

std::pair< std::vector< std::size_t >, ParentDataType > MeshLib::IO::transformToXDMFTopology ( MeshLib::Mesh const & mesh,
std::size_t const offset )

Copies all cells into a new vector. Contiguous data used for writing. The topology is specific to xdmf because it contains the xdmf cell types!! See section topology in https://xdmf.org/index.php/XDMF_Model_and_Format.

Parameters
meshOGS mesh can be mesh or partitionedMesh
offsetLocal offset to transform local to global cell ID. Offset must be zero in serial and must be defined for each process in parallel execution.
Returns
vector containing a copy of the data and the computed ParentDataType

Definition at line 370 of file transformData.cpp.

372{
373 std::vector<MeshLib::Element*> const& elements = mesh.getElements();
374 std::vector<std::size_t> values;
375
376 auto const push_cellnode_ids_to_vector =
377 [&values, &offset](auto const& cell)
378 {
379 values |= ranges::actions::push_back(
380 cell->nodes() | MeshLib::views::ids |
381 ranges::views::transform([&offset](auto const node_id)
382 { return node_id + offset; }));
383 };
384
385 auto const topology_type = getTopologyType(mesh);
386 if (topology_type == ParentDataType::MIXED)
387 {
388 values.reserve(elements.size() * 2); // each cell has at least two
389 // numbers
390 for (auto const& cell : elements)
391 {
392 auto const ogs_cell_type = cell->getCellType();
393 auto const xdmf_cell_id = cellTypeOGS2XDMF(ogs_cell_type).id;
394 values.push_back(to_underlying(xdmf_cell_id));
395 if (ogs_cell_type == MeshLib::CellType::LINE2)
396 {
397 values.push_back(2);
398 }
399 // LINE3 maps to ParentDataType::EDGE_3 -> no special care required
400 push_cellnode_ids_to_vector(cell);
401 }
402 }
403 else
404 {
405 values.reserve(elements.size() * elements[0]->getNumberOfNodes());
406 for (auto const& cell : elements)
407 {
408 push_cellnode_ids_to_vector(cell);
409 }
410 }
411
412 return {values, topology_type};
413}
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:227

References cellTypeOGS2XDMF(), MeshLib::Mesh::getElements(), getTopologyType(), MeshLib::views::ids, MeshLib::LINE2, MIXED, and BaseLib::to_underlying().

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

◆ write_xdmf()

std::function< std::string(std::vector< double >)> MeshLib::IO::write_xdmf ( XdmfData const & geometry,
XdmfData const & topology,
std::vector< XdmfData > const & constant_attributes,
std::vector< XdmfData > const & variable_attributes,
std::string const & h5filename,
std::string const & ogs_version,
std::string const & mesh_name )

Generator function that creates a function capturing the spatial data of a mesh Temporal data can later be passed as argument.

Parameters
geometryMetadata for the geometry (points) of the mesh
topologyMetadata for the topology of the mesh
variable_attributesMeta data for attributes changing over time
constant_attributesMeta data for attributes NOT changing over time
h5filenameName of the file where the actual data was written
ogs_versionOGS Version to be added to XdmfInformation tag
mesh_nameName of the output mesh
Returns
unary function with vector of time step values, returning XDMF string

Definition at line 107 of file writeXdmf.cpp.

113{
114 // Generates function that writes <DataItem>. Late binding needed because
115 // maximum number of steps unknown. Time step and h5filename are captured to
116 // return a unary function
117 // _a suffix is a user defined literal for fmt named arguments
118 auto const time_dataitem_genfn =
119 [](unsigned long long const time_step, int const max_step,
120 std::string const& h5filename, std::string const& mesh_name)
121 {
122 return
123 [time_step, max_step, h5filename, mesh_name](auto const& xdmfdata)
124 {
125 return fmt::format(
126 fmt::runtime("\n\t\t<DataItem DataType=\"{datatype}\" "
127 "Dimensions=\"{local_dimensions}\" "
128 "Format=\"HDF\" "
129 "Precision=\"{precision}\">"
130 "{filename}:/meshes/{meshname}/{datasetname}|"
131 "{time_step} {starts}:1 {strides}:1 "
132 "{local_dimensions}:{max_step} "
133 "{global_dimensions}</"
134 "DataItem>"),
135 "datatype"_a = getPropertyDataTypeString(xdmfdata.data_type),
136 "local_dimensions"_a =
137 fmt::join(xdmfdata.global_block_dims, " "),
138 "precision"_a = getPropertyDataTypeSize(xdmfdata.data_type),
139 "filename"_a = h5filename,
140 "meshname"_a = mesh_name,
141 "datasetname"_a = xdmfdata.name,
142 "starts"_a = fmt::join(xdmfdata.starts, " "),
143 "strides"_a = fmt::join(xdmfdata.strides, " "),
144 "global_dimensions"_a =
145 fmt::join(xdmfdata.global_block_dims, " "),
146 "time_step"_a = fmt::format("{}", time_step),
147 "max_step"_a = fmt::format("{}", max_step));
148 };
149 };
150
151 // string_join could be moved to ogs lib if used more
152 auto const string_join_fn = [](auto const& collection)
153 { return fmt::to_string(fmt::join(collection, "")); };
154
155 // m_bind could be moved to ogs lib if used more
156 auto const m_bind_fn = [](auto const& transform, auto const& join)
157 {
158 return [join, transform](auto const& collection)
159 {
160 std::vector<std::string> temp;
161 temp.reserve(collection.size());
162 std::transform(collection.begin(), collection.end(),
163 std::back_inserter(temp), transform);
164 return join(temp);
165 };
166 };
167
168 // XDMF if part of the data that is already written in any previous step
169 // subsequent steps can refer to this data collection of elements navigates
170 // the xml tree (take first child, then in child take 2nd child ...)
171 auto const pointer_transfrom = [](auto const& elements)
172 {
173 return fmt::format(
174 fmt::runtime("\n\t<xi:include xpointer=\"element(/{elements})\"/>"),
175 "elements"_a = fmt::join(elements, "/"));
176 };
177
178 // Defines content of <Attribute> in XDMF, child of Attribute is a single
179 // <DataItem>
180 auto const attribute_transform =
181 [](XdmfData const& attribute, auto const& dataitem_transform)
182 {
183 return fmt::format(
184 "\n\t<Attribute Center=\"{center}\" ElementCell=\"\" "
185 "ElementDegree=\"0\" "
186 "ElementFamily=\"\" ItemType=\"\" Name=\"{name}\" "
187 "Type=\"None\">{dataitem}\n\t</Attribute>",
188 "center"_a = meshItemTypeString(attribute.attribute_center),
189 "name"_a = attribute.name,
190 "dataitem"_a = dataitem_transform(attribute));
191 };
192
193 // Define content of <Geometry> in XDMF, same as attribute_transform
194 auto const geometry_transform =
195 [](XdmfData const& geometry, auto const& dataitem_transform)
196 {
197 return fmt::format(
198 fmt::runtime("\n\t<Geometry Origin=\"\" "
199 "Type=\"XYZ\">{dataitem}\n\t</Geometry>"),
200 "dataitem"_a = dataitem_transform(geometry));
201 };
202
203 auto tag_string = [](auto const& topology, int nodes_per_element,
204 auto const& dataitem_transform)
205 {
206 return fmt::format(
207 fmt::runtime("\n\t<Topology Type=\"{topology_type}\" "
208 "NodesPerElement=\"{nodes_per_element}\" "
209 "NumberOfElements=\"{number_of_elements}\">{dataitem}"
210 "\n\t</Topology>"),
211 "topology_type"_a =
212 ParentDataType2String(*topology.parent_data_type).first,
213 "dataitem"_a = dataitem_transform(topology),
214 "nodes_per_element"_a = nodes_per_element,
215 "number_of_elements"_a = topology.size_partitioned_dim);
216 };
217
218 // Define content of <Topology> in XDMF, same as attribute_transform
219 auto const topology_transform =
220 [&tag_string](XdmfData const& topology, auto const& dataitem_transform)
221 {
222 switch (*topology.parent_data_type)
223 {
225 return tag_string(topology, 1, dataitem_transform);
227 return tag_string(topology, 2, dataitem_transform);
245 return fmt::format(
246 fmt::runtime(
247 "\n\t<Topology "
248 "Type=\"{topology_type}\">{dataitem}\n\t</Topology>"),
249 "topology_type"_a =
250 ParentDataType2String(*topology.parent_data_type).first,
251 "dataitem"_a = dataitem_transform(topology));
252 }
253 OGS_FATAL("Could not transform unknown XDMF topology type");
254 return std::string{};
255 };
256
257 // Defines content of <Grid> of a single time step, takes specific transform
258 // functions for all of its child elements
259 auto const grid_transform =
260 [](double const& time_value, auto const& geometry, auto const& topology,
261 auto const& constant_attributes, auto const& variable_attributes)
262 {
263 return fmt::format(
264 R"(
265<Grid Name="Grid" GridType="Uniform">
266 <Time Value="{time_value:.{precision}g}"/>
267{geometry}
268{topology}
269{fix_attributes}
270{variable_attributes}
271</Grid>)",
272 "time_value"_a = time_value,
273 // Output of "Time Value" with sufficient precision.
274 "precision"_a = std::numeric_limits<double>::max_digits10,
275 "geometry"_a = geometry,
276 "topology"_a = topology,
277 "fix_attributes"_a = constant_attributes,
278 "variable_attributes"_a = variable_attributes);
279 };
280
281 // An attribute may change over time (variable) or stay constant
282 enum class time_attribute
283 {
284 constant,
285 variable
286 };
287
288 // Generates a function that either writes the data or points to existing
289 // data
290 auto const time_step_fn = [time_dataitem_genfn, pointer_transfrom,
291 h5filename,
292 mesh_name](auto const& component_transform,
293 unsigned long long const time_step,
294 int const max_step,
295 time_attribute const variable)
296 {
297 return [component_transform, time_dataitem_genfn, pointer_transfrom,
298 variable, time_step, max_step, h5filename,
299 mesh_name](XdmfData const& attr)
300 {
301 // new data arrived
302 bool const changed =
303 ((time_step > 0 && (variable == time_attribute::variable)) ||
304 time_step == 0);
305 if (changed)
306 {
307 auto dataitem = time_dataitem_genfn(time_step, max_step,
308 h5filename, mesh_name);
309 return component_transform(attr, dataitem);
310 }
311 else
312 {
313 std::array<unsigned int, 5> position = {1, 1, 2, 1, attr.index};
314 return pointer_transfrom(position);
315 };
316 };
317 };
318
319 // Top-Level transform function (take spatial and temporal transform
320 // functions) and return the time depended grid transform function
321 // ToDo (tm) Remove capturing m_bind and string_join as helper function
322
323 auto const time_grid_transform =
324 [time_step_fn, m_bind_fn, string_join_fn, grid_transform,
325 geometry_transform, topology_transform, attribute_transform](
326 double const time, int const step, int const max_step,
327 auto const& geometry, auto const& topology,
328 auto const& constant_attributes, auto const& variable_attributes)
329 {
330 auto const time_step_geometry_transform = time_step_fn(
331 geometry_transform, step, max_step, time_attribute::constant);
332 auto const time_step_topology_transform = time_step_fn(
333 topology_transform, step, max_step, time_attribute::constant);
334 auto const time_step_const_attribute_fn = time_step_fn(
335 attribute_transform, step, max_step, time_attribute::constant);
336 auto const time_step_variable_attribute_fn = time_step_fn(
337 attribute_transform, step, max_step, time_attribute::variable);
338
339 // lift both functions from handling single attributes to work with
340 // collection of attributes
341 auto const variable_attributes_transform =
342 m_bind_fn(time_step_variable_attribute_fn, string_join_fn);
343 auto const constant_attributes_transform =
344 m_bind_fn(time_step_const_attribute_fn, string_join_fn);
345
346 return grid_transform(
347 time, time_step_geometry_transform(geometry),
348 time_step_topology_transform(topology),
349 constant_attributes_transform(constant_attributes),
350 variable_attributes_transform(variable_attributes));
351 };
352
353 // Function that combines all the data from spatial (geometry, topology,
354 // attribute) and temporal domain (time) with the time domain
355 // (time_step_fn). And writes <Grid CollectionType="Temporal">
356 auto const temporal_grid_collection_transform =
357 [time_grid_transform](
358 auto const& times, auto const& geometry, auto const& topology,
359 auto const& constant_attributes, auto const& variable_attributes)
360 {
361 // c++20 ranges zip missing
362 auto const max_step = times.size();
363 std::vector<std::string> grids;
364 grids.reserve(max_step);
365 for (size_t time_step = 0; time_step < max_step; ++time_step)
366 {
367 grids.push_back(time_grid_transform(
368 times[time_step], time_step, max_step, geometry, topology,
369 constant_attributes, variable_attributes));
370 }
371 return fmt::format(
372 fmt::runtime(
373 "\n<Grid CollectionType=\"Temporal\" GridType=\"Collection\" "
374 "Name=\"Collection\">{grids}\n</Grid>\n"),
375 "grids"_a = fmt::join(grids, ""));
376 };
377
378 // Generator function that return a function that represents complete xdmf
379 // structure
380 auto const xdmf_writer_fn = [temporal_grid_collection_transform](
381 auto ogs_version, auto geometry,
382 auto topology, auto constant_attributes,
383 auto variable_attributes)
384 {
385 // This function will be the return of the surrounding named function
386 // capture by copy - it is the function that will survive this scope!!
387 // Use generalized lambda capture to move for optimization
388 return [temporal_grid_collection_transform,
389 ogs_version = std::move(ogs_version),
390 geometry = std::move(geometry), topology = std::move(topology),
391 constant_attributes = std::move(constant_attributes),
392 variable_attributes = std::move(variable_attributes)](
393 std::vector<double> const& times)
394 {
395 return fmt::format(
396 "<?xml version=\"1.0\" encoding=\"utf-8\"?>\n<Xdmf "
397 "xmlns:xi=\"http://www.w3.org/2001/XInclude\" "
398 "Version=\"3.0\">\n<Domain>\n<Information Name=\"OGS_VERSION\" "
399 "Value=\"{ogs_version}\"/>{grid_collection}</Domain>\n</Xdmf>",
400 "ogs_version"_a = ogs_version,
401 "grid_collection"_a = temporal_grid_collection_transform(
402 times, geometry, topology, constant_attributes,
GITINFOLIB_EXPORT const std::string ogs_version
static std::string getPropertyDataTypeString(MeshPropertyDataType const &ogs_data_type)
Definition writeXdmf.cpp:74
static std::string meshItemTypeString(std::optional< MeshItemType > const &item_type)
Definition writeXdmf.cpp:44
static std::string getPropertyDataTypeSize(MeshPropertyDataType const &ogs_data_type)
Definition writeXdmf.cpp:99

References MeshLib::IO::XdmfData::attribute_center, EDGE_3, getPropertyDataTypeSize(), getPropertyDataTypeString(), HEXAHEDRON, HEXAHEDRON_20, HEXAHEDRON_27, meshItemTypeString(), MIXED, MeshLib::IO::XdmfData::name, OGS_FATAL, MeshLib::IO::XdmfData::parent_data_type, ParentDataType2String(), POLYLINE, POLYVERTEX, PYRAMID, PYRAMID_13, QUADRILATERAL, QUADRILATERAL_8, QUADRILATERAL_9, MeshLib::IO::XdmfData::size_partitioned_dim, TETRAHEDRON, TETRAHEDRON_10, TRIANGLE, TRIANGLE_6, WEDGE, WEDGE_15, and WEDGE_18.

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

◆ writeMeshToFile()

int MeshLib::IO::writeMeshToFile ( const MeshLib::Mesh & mesh,
std::filesystem::path const & file_path,
std::set< std::string > variable_output_names )

Definition at line 24 of file writeMeshToFile.cpp.

28{
29 if (file_path.extension().string() == ".msh")
30 {
32 meshIO.setMesh(&mesh);
34 return 0;
35 }
36 if (file_path.extension().string() == ".vtu")
37 {
38 MeshLib::IO::VtuInterface writer(&mesh);
39 auto const result = writer.writeToFile(file_path);
40 if (!result)
41 {
42 ERR("writeMeshToFile(): Could not write mesh to '{:s}'.",
43 file_path.string());
44 return -1;
45 }
46 return 0;
47 }
48 if (file_path.extension().string() == ".xdmf")
49 {
50 std::vector<std::reference_wrapper<const MeshLib::Mesh>> meshes;
51 const std::reference_wrapper<const MeshLib::Mesh> mr = mesh;
52 meshes.push_back(mr);
53 MeshLib::IO::XdmfHdfWriter(std::move(meshes), file_path, 0, 0.0,
54 variable_output_names, true, 1, 1048576);
55 return 0;
56 }
57 ERR("writeMeshToFile(): Unknown file extension '{:s}'. Can not write file "
58 "'{:s}'.",
59 file_path.extension().string(), file_path.string());
60 return 0;
61}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
std::string writeToString()
Writes the object to a string.
Definition Writer.cpp:31
Interface for handling mesh files from OGS-5 and below. (*.msh files)
Definition MeshIO.h:37
void setMesh(const MeshLib::Mesh *mesh)
Set mesh for writing.
Definition MeshIO.cpp:437
Reads and writes VtkXMLUnstructuredGrid-files (vtu) to and from OGS data structures....
int writeStringToFile(std::string_view content, std::filesystem::path const &file_path)
Definition Writer.cpp:45

References ERR(), MeshLib::IO::Legacy::MeshIO::setMesh(), BaseLib::IO::writeStringToFile(), MeshLib::IO::VtuInterface::writeToFile(), and BaseLib::IO::Writer::writeToString().

Referenced by FileIO::Gocad::generateFaceSets(), identifyAndWriteBoundaryMeshes(), main(), main(), anonymous_namespace{postLIE.cpp}::postVTU(), and FileIO::XmlPrjInterface::write().

◆ writePropertyVectorMetaData() [1/2]

void MeshLib::IO::writePropertyVectorMetaData ( PropertyVectorMetaData const & pvmd)
inline

Definition at line 74 of file PropertyVectorMetaData.h.

75{
76 DBUG(
77 "name: '{:s}':\t is_int_data_type: {:d}, is_data_type_signed: "
78 "{:d}, data_type_size_in_bytes: {:d}, number of components / "
79 "tuples: {:d} / {:d}",
80 pvmd.property_name, pvmd.is_int_type, pvmd.is_data_type_signed,
81 pvmd.data_type_size_in_bytes, pvmd.number_of_components,
82 pvmd.number_of_tuples);
83}

References MeshLib::IO::PropertyVectorMetaData::data_type_size_in_bytes, DBUG(), MeshLib::IO::PropertyVectorMetaData::is_data_type_signed, MeshLib::IO::PropertyVectorMetaData::is_int_type, MeshLib::IO::PropertyVectorMetaData::number_of_components, MeshLib::IO::PropertyVectorMetaData::number_of_tuples, and MeshLib::IO::PropertyVectorMetaData::property_name.

◆ writePropertyVectorMetaData() [2/2]

void MeshLib::IO::writePropertyVectorMetaData ( std::ostream & os,
PropertyVectorMetaData const & pvmd )
inline

Definition at line 47 of file PropertyVectorMetaData.h.

49{
50 std::string::size_type s(pvmd.property_name.length());
51 os.write(reinterpret_cast<char*>(&s), sizeof(std::string::size_type));
52
53 os.write(
54 const_cast<char*>(
55 const_cast<PropertyVectorMetaData&>(pvmd).property_name.data()),
56 s);
57 os.write(reinterpret_cast<char*>(
58 &const_cast<PropertyVectorMetaData&>(pvmd).is_int_type),
59 sizeof(bool));
60 os.write(reinterpret_cast<char*>(&const_cast<PropertyVectorMetaData&>(
61 pvmd).is_data_type_signed),
62 sizeof(bool));
63 os.write(reinterpret_cast<char*>(&const_cast<PropertyVectorMetaData&>(
64 pvmd).data_type_size_in_bytes),
65 sizeof(unsigned long));
66 os.write(reinterpret_cast<char*>(&const_cast<PropertyVectorMetaData&>(
67 pvmd).number_of_components),
68 sizeof(unsigned long));
69 os.write(reinterpret_cast<char*>(
70 &const_cast<PropertyVectorMetaData&>(pvmd).number_of_tuples),
71 sizeof(unsigned long));
72}

References MeshLib::IO::PropertyVectorMetaData::property_name.

Referenced by MeshLib::IO::NodePartitionedMeshReader::readProperties(), and ApplicationUtils::writePropertyVector().

◆ writePropertyVectorPartitionMetaData()

void MeshLib::IO::writePropertyVectorPartitionMetaData ( std::ostream & os,
PropertyVectorPartitionMetaData const & pvpmd )
inline

Definition at line 126 of file PropertyVectorMetaData.h.

128{
129 os.write(reinterpret_cast<char*>(
130 &const_cast<PropertyVectorPartitionMetaData&>(pvpmd)
131 .offset),
132 sizeof(unsigned long));
133 os.write(reinterpret_cast<char*>(
134 &const_cast<PropertyVectorPartitionMetaData&>(pvpmd)
135 .number_of_tuples),
136 sizeof(unsigned long));
137}

Referenced by ApplicationUtils::writeProperties().

◆ writeVtu()

int MeshLib::IO::writeVtu ( MeshLib::Mesh const & mesh,
std::string const & file_name,
int const data_mode )

Definition at line 160 of file VtuInterface.cpp.

162{
163 MeshLib::IO::VtuInterface writer(&mesh, data_mode);
164 auto const result = writer.writeToFile(file_name);
165 if (!result)
166 {
167 ERR("writeMeshToFile(): Could not write mesh to '{:s}'.", file_name);
168 return -1;
169 }
170 return 0;
171}

References ERR(), and MeshLib::IO::VtuInterface::writeToFile().

Referenced by main().

Variable Documentation

◆ elem_type_ogs2xdmf

auto MeshLib::IO::elem_type_ogs2xdmf = elemOGSTypeToXDMFType()
constexpr

Definition at line 85 of file transformData.cpp.

Referenced by cellTypeOGS2XDMF().

◆ mesh_item_type_strings

char const* MeshLib::IO::mesh_item_type_strings[]
constexpr
Initial value:
= {"Node", "Edge", "Face",
"Cell", "Other"}

Definition at line 40 of file writeXdmf.cpp.

40 {"Node", "Edge", "Face",
41 "Cell", "Other"};

Referenced by meshItemTypeString().

◆ ogs_to_xdmf_type_fn

auto MeshLib::IO::ogs_to_xdmf_type_fn = meshPropertyDatatypeSize()
inline

Definition at line 96 of file writeXdmf.cpp.

Referenced by getPropertyDataTypeSize().