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