OGS
Vtu2GridDialog.cpp
Go to the documentation of this file.
1
14#include "Vtu2GridDialog.h"
15
16#include <vtkXMLUnstructuredGridWriter.h>
17
18#include <QStringList>
19#include <QStringListModel>
20#include <optional>
21#include <string>
22
23#include "Base/OGSError.h"
25#include "GeoLib/AABB.h"
27#include "MeshLib/Mesh.h"
28#include "MeshLib/Node.h"
30#include "MeshModel.h"
33
34Vtu2GridDialog::Vtu2GridDialog(MeshModel& mesh_model, QDialog* parent)
35 : QDialog(parent), _mesh_model(mesh_model)
36{
37 setupUi(this);
38 QStringList MeshList;
39
40 for (int model_index = 0; model_index < mesh_model.rowCount();
41 ++model_index)
42 {
43 auto const* mesh = mesh_model.getMesh(mesh_model.index(model_index, 0));
44 MeshList.append(QString::fromStdString(mesh->getName()));
45 }
46
47 if (MeshList.empty())
48 {
49 MeshList.append("[No Mesh available.]");
50 this->expectedVoxelLabel->setText(
51 "Expected number of voxels: undefined");
52 }
53
54 _allMeshes.setStringList(MeshList);
55 this->meshListBox->addItems(_allMeshes.stringList());
56 this->xlineEdit->setFocus();
57}
58
59std::optional<std::array<double, 3>> fillXYZ(QString xin, QString yin,
60 QString zin)
61{
62 bool ok;
63 if (!xin.toDouble(&ok))
64 {
65 return std::nullopt;
66 }
67 double const xinput = xin.toDouble();
68 double const yinput = (yin.toDouble(&ok)) ? yin.toDouble() : xinput;
69 double const zinput = (zin.toDouble(&ok)) ? zin.toDouble() : xinput;
70
71 if (xinput <= 0 || yinput <= 0 || zinput <= 0)
72 {
73 return std::nullopt;
74 }
75
76 return std::optional<std::array<double, 3>>{{xinput, yinput, zinput}};
77}
78
79Eigen::Vector3d getMeshExtent(MeshLib::Mesh const* _mesh)
80{
81 auto const& nodes = _mesh->getNodes();
82 GeoLib::AABB const aabb(nodes.cbegin(), nodes.cend());
83 auto const& min = aabb.getMinPoint();
84 auto const& max = aabb.getMaxPoint();
85 return max - min;
86}
87
89{
90 if (_allMeshes.stringList()[0] == "[No Mesh available.]")
91 {
92 this->expectedVoxelLabel->setText("approximated Voxel: undefined");
93 return;
94 }
95
96 auto const opt_xyz =
97 fillXYZ(this->xlineEdit->text(), this->ylineEdit->text(),
98 this->zlineEdit->text());
99
100 if (!opt_xyz)
101 {
102 this->expectedVoxelLabel->setText("approximated Voxel: undefined");
103 return;
104 }
105
106 auto const delta = getMeshExtent(
107 _mesh_model.getMesh(this->meshListBox->currentText().toStdString()));
108 double const expected_voxel = (delta[0]) * (delta[1]) * (delta[2]) /
109 (*opt_xyz)[0] / (*opt_xyz)[1] / (*opt_xyz)[2];
110
111 int const exponent = std::floor(std::log10(std::abs(expected_voxel)));
112 this->expectedVoxelLabel->setText(
113 "approximated Voxel = " +
114 QString::number(std::round(expected_voxel / std::pow(10, exponent))) +
115 " x 10^" + QString::number(exponent));
116}
117
122
127
132
134{
135 using namespace MeshToolsLib::MeshGenerator;
136 if (this->meshListBox->currentText().toStdString() ==
137 "[No Mesh available.]")
138 {
140 "Please specify the input meshes. It has to be a 3D mesh.");
141 return;
142 }
143
144 auto cellsize = fillXYZ(this->xlineEdit->text(), this->ylineEdit->text(),
145 this->zlineEdit->text());
146
147 if (!cellsize)
148 {
150 "At least the x-length of a voxel must be specified and > 0.\n If "
151 "y-/z-input "
152 "are not specified, equal to 0, or not a real number, they are "
153 "treated as "
154 "the x-input.");
155 }
156 auto _mesh(
157 _mesh_model.getMesh(this->meshListBox->currentText().toStdString()));
158
159 if (_mesh->MeshLib::Mesh::getDimension() < 3)
160 {
161 OGSError::box("The dimension of the mesh has to be 3.");
162 return;
163 }
164
165 vtkNew<MeshLib::VtkMappedMeshSource> vtkSource;
166 vtkSource->SetMesh(_mesh);
167 vtkSource->Update();
168 vtkSmartPointer<vtkUnstructuredGrid> mesh = vtkSource->GetOutput();
169
170 double* const bounds = mesh->GetBounds();
171 MathLib::Point3d const min(
172 std::array<double, 3>{bounds[0], bounds[2], bounds[4]});
173 MathLib::Point3d const max(
174 std::array<double, 3>{bounds[1], bounds[3], bounds[5]});
175 std::array<double, 3> ranges = {max[0] - min[0], max[1] - min[1],
176 max[2] - min[2]};
177 if (ranges[0] < 0 || ranges[1] < 0 || ranges[2] < 0)
178 {
180 "The range (max-min of the bounding box) is not allowed to be < 0");
181 }
182 std::array<std::size_t, 3> const dims =
183 VoxelGridFromMesh::getNumberOfVoxelPerDimension(ranges, *cellsize);
184 std::unique_ptr<MeshLib::Mesh> grid(
185 generateRegularHexMesh(dims[0], dims[1], dims[2], (*cellsize)[0],
186 (*cellsize)[1], (*cellsize)[2], min, "grid"));
187 if (grid == nullptr)
188 {
189 OGSError::box(QString::fromStdString(
190 fmt::format("Could not generate regular hex mesh. With "
191 "parameters dims={} {} {}, cellsize={} {} {}",
192 dims[0], dims[1], dims[2], (*cellsize)[0],
193 (*cellsize)[1], (*cellsize)[2])));
194 }
195
196 std::vector<int>* cell_ids =
197 grid->getProperties().createNewPropertyVector<int>(
198 VoxelGridFromMesh::cell_id_name, MeshLib::MeshItemType::Cell, 1);
199 if (cell_ids == nullptr)
200 {
201 OGSError::box("Could not create cell ids.");
202 }
203 *cell_ids = VoxelGridFromMesh::assignCellIds(mesh, min, dims, *cellsize);
204 if (!VoxelGridFromMesh::removeUnusedGridCells(mesh, grid))
205 {
206 return;
207 }
208
209 VoxelGridFromMesh::mapMeshArraysOntoGrid(mesh, grid);
210
211 if (grid == nullptr)
212 {
213 OGSError::box("No voxelgrid could be created from the mesh.");
214 return;
215 }
216
217 _mesh_model.addMesh(grid.release());
218 this->done(QDialog::Accepted);
219}
Definition of the AABB class.
Definition of the MeshModel class.
Definition of the Mesh class.
Definition of the Node class.
Definition of the OGSError class.
Implementation of the StrictDoubleValidator class.
VtkMappedMeshSource is a source class to transform OGS meshes into complete vtkUnstructuredGrids....
std::optional< std::array< double, 3 > > fillXYZ(QString xin, QString yin, QString zin)
Eigen::Vector3d getMeshExtent(MeshLib::Mesh const *_mesh)
Definition of the Vtu2GridDialog class.
Class AABB is an axis aligned bounding box around a given set of geometric points of (template) type ...
Definition AABB.h:56
Eigen::Vector3d const & getMaxPoint() const
Definition AABB.h:187
Eigen::Vector3d const & getMinPoint() const
Definition AABB.h:180
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
void addMesh(std::unique_ptr< MeshLib::Mesh > mesh)
Adds a new mesh.
Definition MeshModel.cpp:52
const MeshLib::Mesh * getMesh(const QModelIndex &idx) const
Returns the mesh with the given index.
Definition MeshModel.cpp:96
static void box(const QString &e)
Definition OGSError.cpp:23
int rowCount(const QModelIndex &parent=QModelIndex()) const override
QModelIndex index(int row, int column, const QModelIndex &parent=QModelIndex()) const override
Definition TreeModel.cpp:50
void on_xlineEdit_textChanged()
Vtu2GridDialog(MeshModel &mesh_model, QDialog *parent=nullptr)
void on_zlineEdit_textChanged()
void updateExpectedVoxel()
As the x/y/z input changes an estimation of the expected Voxel is given.
MeshModel & _mesh_model
void on_ylineEdit_textChanged()
void accept() override
Instructions if the OK-Button has been pressed.
QStringListModel _allMeshes