OGS
Layers2Grid.cpp File Reference

Detailed Description

Definition in file Layers2Grid.cpp.

#include <algorithm>
#include <memory>
#include <string>
#include <vector>
#include <tclap/CmdLine.h>
#include "BaseLib/IO/readStringListFromFile.h"
#include "BaseLib/MPI.h"
#include "GeoLib/AABB.h"
#include "InfoLib/GitInfo.h"
#include "MathLib/Point3d.h"
#include "MeshLib/Elements/Element.h"
#include "MeshLib/IO/VtkIO/VtuInterface.h"
#include "MeshLib/IO/readMeshFromFile.h"
#include "MeshToolsLib/MeshGenerators/VoxelGridFromLayeredMeshes.h"
Include dependency graph for Layers2Grid.cpp:

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Function Documentation

◆ main()

int main ( int argc,
char * argv[] )

Definition at line 28 of file Layers2Grid.cpp.

29{
30 TCLAP::CmdLine cmd(
31 "Reads a list of 2D unstructured mesh layers and samples them onto a "
32 "structured grid of the same extent. Note, that a large cube size may "
33 "result in an undersampling of the original structure.\nCube sizes are "
34 "defines by x/y/z-parameters. For equilateral cubes, only the "
35 "x-parameter needs to be set.\n\n"
36 "OpenGeoSys-6 software, version " +
38 ".\n"
39 "Copyright (c) 2012-2024, OpenGeoSys Community "
40 "(http://www.opengeosys.org)",
42 TCLAP::SwitchArg dilate_arg(
43 "d", "dilate",
44 "assign mat IDs based on single nodes instead of a majority of nodes, "
45 "which can result in a slightly increased voxel grid extent",
46 false);
47 cmd.add(dilate_arg);
48
49 TCLAP::ValueArg<double> z_arg("z", "cellsize-z",
50 "edge length of cubes in z-direction (depth)",
51 false, 1000, "floating point number");
52 cmd.add(z_arg);
53
54 TCLAP::ValueArg<double> y_arg(
55 "y", "cellsize-y", "edge length of cubes in y-direction (latitude)",
56 false, 1000, "floating point number");
57 cmd.add(y_arg);
58
59 TCLAP::ValueArg<double> x_arg(
60 "x", "cellsize-x",
61 "edge length of cubes in x-direction (longitude) or all directions, if "
62 "y and z are not set",
63 true, 1000, "floating point number");
64 cmd.add(x_arg);
65
66 TCLAP::ValueArg<std::string> output_arg(
67 "o", "output", "name of output mesh (*.vtu)", true, "", "string");
68 cmd.add(output_arg);
69
70 TCLAP::ValueArg<std::string> input_arg(
71 "i", "input",
72 "name of the input file list containing the paths the all input layers "
73 "in correct order from top to bottom",
74 true, "", "string");
75 cmd.add(input_arg);
76 cmd.parse(argc, argv);
77
78 BaseLib::MPI::Setup mpi_setup(argc, argv);
79
80 if ((y_arg.isSet() && !z_arg.isSet()) ||
81 ((!y_arg.isSet() && z_arg.isSet())))
82 {
83 ERR("For equilateral cubes, only x needs to be set. For unequal "
84 "cuboids, all three edge lengths (x/y/z) need to be specified.");
85 return EXIT_FAILURE;
86 }
87
88 double const x_size = x_arg.getValue();
89 double const y_size = (y_arg.isSet()) ? y_arg.getValue() : x_arg.getValue();
90 double const z_size = (z_arg.isSet()) ? z_arg.getValue() : x_arg.getValue();
91 std::array<double, 3> const cellsize = {x_size, y_size, z_size};
92
93 std::string const input_name = input_arg.getValue();
94 std::string const output_name = output_arg.getValue();
95 auto const layer_names = BaseLib::IO::readStringListFromFile(input_name);
96 if (layer_names.size() < 2)
97 {
98 ERR("At least two layers are required to create a 3D Mesh");
99 return EXIT_FAILURE;
100 }
101
102 std::vector<std::unique_ptr<MeshLib::Mesh>> layers;
103 layers.reserve(layer_names.size());
104 constexpr double minval = std::numeric_limits<double>::max();
105 constexpr double maxval = std::numeric_limits<double>::lowest();
106 std::pair<MathLib::Point3d, MathLib::Point3d> extent(
107 MathLib::Point3d{{minval, minval, minval}},
108 MathLib::Point3d{{maxval, maxval, maxval}});
109
110 for (auto const& layer : layer_names)
111 {
112 auto mesh(MeshLib::IO::readMeshFromFile(layer));
113 if (mesh == nullptr)
114 {
115 ERR("Input layer '{:s}' not found. Aborting...", layer);
116 return EXIT_FAILURE;
117 }
118 layers.emplace_back(mesh);
119 }
120 std::vector<MeshLib::Mesh const*> layers_ptr;
121 std::transform(std::begin(layers), std::end(layers),
122 std::back_inserter(layers_ptr),
123 [](auto const& layer) { return layer.get(); });
124 bool const dilate = dilate_arg.getValue();
126 createVoxelFromLayeredMesh(extent, layers_ptr, cellsize, dilate);
127 if (mesh == nullptr)
128 {
129 ERR("The VoxelGrid could not be created.");
130 return EXIT_FAILURE;
131 }
132 MeshLib::IO::VtuInterface vtu(mesh.get());
133 vtu.writeToFile(output_name);
134 return EXIT_SUCCESS;
135}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
Reads and writes VtkXMLUnstructuredGrid-files (vtu) to and from OGS data structures....
std::vector< std::string > readStringListFromFile(std::string const &filename)
Reads non-empty lines from a list of strings from a file into a vector.
GITINFOLIB_EXPORT const std::string ogs_version
MeshLib::Mesh * readMeshFromFile(const std::string &file_name, bool const compute_element_neighbors)
std::unique_ptr< MeshLib::Mesh > createVoxelFromLayeredMesh(std::pair< MathLib::Point3d, MathLib::Point3d > &extent, std::vector< MeshLib::Mesh const * > const &layers, std::array< double, 3 > const cellsize, bool const dilate)

References MeshToolsLib::MeshGenerators::VoxelFromLayeredMeshes::createVoxelFromLayeredMesh(), ERR(), GitInfoLib::GitInfo::ogs_version, MeshLib::IO::readMeshFromFile(), BaseLib::IO::readStringListFromFile(), and MeshLib::IO::VtuInterface::writeToFile().