OGS
BinaryToPVTU.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#include <tclap/CmdLine.h>
5#include <vtkMPIController.h>
6#include <vtkSmartPointer.h>
7
8#include "BaseLib/CPUTime.h"
9#include "BaseLib/FileTools.h"
10#include "BaseLib/Logging.h"
11#include "BaseLib/MPI.h"
12#include "BaseLib/RunTime.h"
14#include "InfoLib/GitInfo.h"
18
19int main(int argc, char* argv[])
20{
21 TCLAP::CmdLine cmd(
22 "Reads the binary mesh format and writes the data as vtus/pvtu.\n"
23 "OpenGeoSys-6 software, version " +
25 ".\n"
26 "Copyright (c) 2012-2026, OpenGeoSys Community "
27 "(http://www.opengeosys.org)",
29 TCLAP::ValueArg<std::string> mesh_input(
30 "i", "mesh-input-file-base",
31 "Input (.bin). The base name of the files containing the input mesh, "
32 "i.e., the file "
33 "name without the string beginning with '_partitioned' and ending "
34 "with ",
35 true, "", "BASE_FILENAME_INPUT_MESH");
36 cmd.add(mesh_input);
37
38 TCLAP::ValueArg<std::string> output_directory_arg(
39 "o", "output",
40 "output directory name and output base file name without extensions",
41 true, "", "directory/base_file_name_without_extensions");
42 cmd.add(output_directory_arg);
43
44 auto log_level_arg = BaseLib::makeLogLevelArg();
45 cmd.add(log_level_arg);
46
47 cmd.parse(argc, argv);
48
49 BaseLib::MPI::Setup mpi_setup(argc, argv);
50 BaseLib::initOGSLogger(log_level_arg.getValue());
51
52 // start the timer
53 BaseLib::RunTime run_timer;
54 run_timer.start();
55 BaseLib::CPUTime CPU_timer;
56 CPU_timer.start();
57
58 // init vtkMPI
59 vtkSmartPointer<vtkMPIController> controller =
60 vtkSmartPointer<vtkMPIController>::New();
61 controller->Initialize(&argc, &argv, 1);
62 vtkMPIController::SetGlobalController(controller);
63
65 dynamic_cast<MeshLib::NodePartitionedMesh*>(
66 MeshLib::IO::readMeshFromFile(mesh_input.getValue()));
67 INFO("Mesh '{:s}' read: {:d} nodes, {:d} elements.",
68 mesh->getName(),
69 mesh->getNumberOfNodes(),
70 mesh->getNumberOfElements());
71 INFO("Mesh read runtime: {:g} s.", run_timer.elapsed());
72 INFO("Mesh read CPU time: {:g} s.", CPU_timer.elapsed());
73
74 // Write output file
75 auto const file_name = output_directory_arg.getValue() + ".vtu";
76 DBUG("Writing output to '{:s}'.", file_name);
77
78 MeshLib::IO::VtuInterface vtu_interface(mesh);
79 vtu_interface.writeToFile(file_name);
80
81 INFO("Total runtime: {:g} s.", run_timer.elapsed());
82 INFO("Total CPU time: {:g} s.", CPU_timer.elapsed());
83 return EXIT_SUCCESS;
84}
int main(int argc, char *argv[])
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
Count CPU time.
Definition CPUTime.h:12
void start()
Start the timer.
Definition CPUTime.h:15
double elapsed() const
Get the elapsed time after started.
Definition CPUTime.h:18
Count the running time.
Definition RunTime.h:18
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:31
void start()
Start the timer.
Definition RunTime.h:21
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:94
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:91
std::size_t getNumberOfElements() const
Get the number of elements.
Definition Mesh.h:88
TCLAP::ValueArg< std::string > makeLogLevelArg()
void initOGSLogger(std::string const &log_level)
Definition Logging.cpp:56
GITINFOLIB_EXPORT const std::string ogs_version
MeshLib::Mesh * readMeshFromFile(const std::string &file_name, bool const compute_element_neighbors)