OGS
Vtu2Grid.cpp
Go to the documentation of this file.
1
10#include <tclap/CmdLine.h>
11
12#include "BaseLib/Logging.h"
13#include "BaseLib/MPI.h"
15#include "InfoLib/GitInfo.h"
18#include "MeshLib/Mesh.h"
23
24int main(int argc, char* argv[])
25{
26 TCLAP::CmdLine cmd(
27 "Reads a 3D unstructured mesh and samples it onto a structured grid of "
28 "the same extent. Cell properties are mapped onto the grid (sampled at "
29 "the centre-points of each cube), node properties are ignored. Note, "
30 "that a large cube size may result in an undersampling of the original "
31 "mesh structure.\nCube sizes are defines by x/y/z-parameters. For "
32 "equilateral cubes, only the x-parameter needs to be set.\n\n"
33 "OpenGeoSys-6 software, version " +
35 ".\n"
36 "Copyright (c) 2012-2025, OpenGeoSys Community "
37 "(http://www.opengeosys.org)",
39
40 TCLAP::ValueArg<double> z_arg("z", "cellsize-z",
41 "edge length of cubes in z-direction "
42 "(depth), (min = 0)",
43 false, 1000, "CELLSIZE_Z");
44 cmd.add(z_arg);
45
46 TCLAP::ValueArg<double> y_arg("y", "cellsize-y",
47 "edge length of cubes in y-direction "
48 "(latitude), (min = 0)",
49 false, 1000, "CELLSIZE_Y");
50 cmd.add(y_arg);
51
52 TCLAP::ValueArg<double> x_arg(
53 "x", "cellsize-x",
54 "edge length of cubes in x-direction (longitude) or all directions, if "
55 "y and z are not set, (min = 0)",
56 true, 1000, "CELLSIZE_X");
57 cmd.add(x_arg);
58
59 TCLAP::ValueArg<std::string> output_arg(
60 "o", "output", "Output (.vtu). The output grid file", true, "",
61 "OUTPUT_FILE");
62 cmd.add(output_arg);
63
64 TCLAP::ValueArg<std::string> input_arg(
65 "i", "input", "Input (.vtu | .msh). The 3D input mesh file", true, "",
66 "INPUT_FILE");
67 cmd.add(input_arg);
68 auto log_level_arg = BaseLib::makeLogLevelArg();
69 cmd.add(log_level_arg);
70 cmd.parse(argc, argv);
71
72 BaseLib::MPI::Setup mpi_setup(argc, argv);
73 BaseLib::initOGSLogger(log_level_arg.getValue());
74
75 if ((y_arg.isSet() && !z_arg.isSet()) ||
76 ((!y_arg.isSet() && z_arg.isSet())))
77 {
78 ERR("For equilateral cubes, only x needs to be set. For unequal "
79 "cuboids, all three edge lengths (x/y/z) need to be specified.");
80 return -1;
81 }
82 using namespace MeshToolsLib::MeshGenerator;
83
84 double const x_size = x_arg.getValue();
85 double const y_size = (y_arg.isSet()) ? y_arg.getValue() : x_arg.getValue();
86 double const z_size = (z_arg.isSet()) ? z_arg.getValue() : x_arg.getValue();
87
88 if (x_size <= 0 || y_size <= 0 || z_size <= 0)
89 {
90 ERR("A cellsize ({},{},{}) is not allowed to be <= 0", x_size, y_size,
91 z_size);
92 return -1;
93 }
94
95 std::array<double, 3> const cellsize = {x_size, y_size, z_size};
96
97 vtkSmartPointer<vtkUnstructuredGrid> mesh =
99 input_arg.getValue());
100 if (mesh == nullptr)
101 {
102 return EXIT_FAILURE;
103 }
104
105 double* const bounds = mesh->GetBounds();
106 MathLib::Point3d const min(
107 std::array<double, 3>{bounds[0], bounds[2], bounds[4]});
108 MathLib::Point3d const max(
109 std::array<double, 3>{bounds[1], bounds[3], bounds[5]});
110 std::array<double, 3> ranges = {max[0] - min[0], max[1] - min[1],
111 max[2] - min[2]};
112 if (ranges[0] < 0 || ranges[1] < 0 || ranges[2] < 0)
113 {
114 ERR("The range ({},{},{}) is not allowed to be < 0", ranges[0],
115 ranges[1], ranges[2]);
116 return -1;
117 }
118 std::array<std::size_t, 3> const dims =
119 VoxelGridFromMesh::getNumberOfVoxelPerDimension(ranges, cellsize);
120 std::unique_ptr<MeshLib::Mesh> grid(
122 dims[0], dims[1], dims[2], cellsize[0], cellsize[1], cellsize[2],
123 min, "grid"));
124
125 std::vector<int> const tmp_ids =
126 VoxelGridFromMesh::assignCellIds(mesh, min, dims, cellsize);
127 auto* const cell_ids = grid->getProperties().createNewPropertyVector<int>(
128 VoxelGridFromMesh::cell_id_name, MeshLib::MeshItemType::Cell, 1);
129 assert(cell_ids);
130 cell_ids->assign(tmp_ids);
131
132 if (!VoxelGridFromMesh::removeUnusedGridCells(mesh, grid))
133 {
134 return EXIT_FAILURE;
135 }
136
137 VoxelGridFromMesh::mapMeshArraysOntoGrid(mesh, grid);
138
139 if (MeshLib::IO::writeMeshToFile(*grid, output_arg.getValue()) != 0)
140 {
141 return EXIT_FAILURE;
142 }
143
144 return EXIT_SUCCESS;
145}
Git information.
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:48
Definition of the Mesh class.
int main(int argc, char *argv[])
Definition Vtu2Grid.cpp:24
Implementation of the VtuInterface class.
static vtkSmartPointer< vtkUnstructuredGrid > readVtuFileToVtkUnstructuredGrid(std::string const &file_name)
TCLAP::ValueArg< std::string > makeLogLevelArg()
void initOGSLogger(std::string const &log_level)
Definition Logging.cpp:64
GITINFOLIB_EXPORT const std::string ogs_version
int writeMeshToFile(const MeshLib::Mesh &mesh, std::filesystem::path const &file_path, std::set< std::string > variable_output_names)
MeshLib::Mesh * generateRegularHexMesh(const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, const BaseLib::ISubdivision &div_z, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")