OGS
BinaryToPVTU.cpp
Go to the documentation of this file.
1
13#include <spdlog/spdlog.h>
14#include <tclap/CmdLine.h>
15#include <vtkMPIController.h>
16#include <vtkSmartPointer.h>
17
18#include "BaseLib/CPUTime.h"
19#include "BaseLib/FileTools.h"
20#include "BaseLib/MPI.h"
21#include "BaseLib/RunTime.h"
22#include "InfoLib/GitInfo.h"
26
27int main(int argc, char* argv[])
28{
29 TCLAP::CmdLine cmd(
30 "Reads the binary mesh format and writes the data as vtus/pvtu.\n"
31 "OpenGeoSys-6 software, version " +
33 ".\n"
34 "Copyright (c) 2012-2025, OpenGeoSys Community "
35 "(http://www.opengeosys.org)",
37 TCLAP::ValueArg<std::string> mesh_input(
38 "i", "mesh-input-file-base",
39 "the base name of the files containing the input mesh, i.e., the file "
40 "name without the string beginning with '_partitioned' and ending with "
41 "'.bin'",
42 true, "", "base_file_name_of_input_mesh");
43 cmd.add(mesh_input);
44
45 TCLAP::ValueArg<std::string> output_directory_arg(
46 "o", "output",
47 "output directory name and output base file name without extensions",
48 true, "", "directory/base_file_name_without_extensions");
49 cmd.add(output_directory_arg);
50
51 TCLAP::ValueArg<std::string> log_level_arg(
52 "l", "log-level",
53 "the verbosity of logging messages: none, error, warn, info, debug, "
54 "all",
55 false,
56#ifdef NDEBUG
57 "info",
58#else
59 "all",
60#endif
61 "LOG_LEVEL");
62 cmd.add(log_level_arg);
63
64 cmd.parse(argc, argv);
65
66 BaseLib::setConsoleLogLevel(log_level_arg.getValue());
67 spdlog::set_pattern("%^%l:%$ %v");
68 spdlog::set_error_handler(
69 [](const std::string& msg)
70 {
71 std::cerr << "spdlog error: " << msg << std::endl;
72 OGS_FATAL("spdlog logger error occurred.");
73 });
74
75 BaseLib::MPI::Setup mpi_setup(argc, argv);
76 // start the timer
77 BaseLib::RunTime run_timer;
78 run_timer.start();
79 BaseLib::CPUTime CPU_timer;
80 CPU_timer.start();
81
82 // add mpi_rank to logger output
83 int mpi_rank;
84 MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
85 spdlog::set_pattern(fmt::format("[{}] %^%l:%$ %v", mpi_rank));
86
87 // init vtkMPI
88 vtkSmartPointer<vtkMPIController> controller =
89 vtkSmartPointer<vtkMPIController>::New();
90 controller->Initialize(&argc, &argv, 1);
91 vtkMPIController::SetGlobalController(controller);
92
94 dynamic_cast<MeshLib::NodePartitionedMesh*>(
95 MeshLib::IO::readMeshFromFile(mesh_input.getValue()));
96 INFO("Mesh '{:s}' read: {:d} nodes, {:d} elements.",
97 mesh->getName(),
98 mesh->getNumberOfNodes(),
99 mesh->getNumberOfElements());
100 INFO("Mesh read runtime: {:g} s.", run_timer.elapsed());
101 INFO("Mesh read CPU time: {:g} s.", CPU_timer.elapsed());
102
103 // Write output file
104 auto const file_name = output_directory_arg.getValue() + ".vtu";
105 DBUG("Writing output to '{:s}'.", file_name);
106
107 MeshLib::IO::VtuInterface vtu_interface(mesh);
108 vtu_interface.writeToFile(file_name);
109
110 INFO("Total runtime: {:g} s.", run_timer.elapsed());
111 INFO("Total CPU time: {:g} s.", CPU_timer.elapsed());
112 return EXIT_SUCCESS;
113}
int main(int argc, char *argv[])
Definition of the CPUTime class.
#define OGS_FATAL(...)
Definition Error.h:26
Filename manipulation routines.
Git information.
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Definition of mesh class for partitioned mesh (by node) for parallel computing within the framework o...
Definition of the RunTime class.
Implementation of the VtuInterface class.
Count CPU time.
Definition CPUTime.h:23
void start()
Start the timer.
Definition CPUTime.h:26
double elapsed() const
Get the elapsed time after started.
Definition CPUTime.h:29
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
Reads and writes VtkXMLUnstructuredGrid-files (vtu) to and from OGS data structures....
bool writeToFile(std::filesystem::path const &file_path)
const std::string getName() const
Get name of the mesh.
Definition Mesh.h:103
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:100
std::size_t getNumberOfElements() const
Get the number of elements.
Definition Mesh.h:97
void setConsoleLogLevel(std::string const &level_string)
Definition Logging.cpp:37
GITINFOLIB_EXPORT const std::string ogs_version
MeshLib::Mesh * readMeshFromFile(const std::string &file_name, bool const compute_element_neighbors)
Definition of readMeshFromFile function.