PRISMS-PF Manual
Loading...
Searching...
No Matches
read_vtk.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/exceptions.h>
7#include <deal.II/base/point.h>
8#include <deal.II/lac/vector.h>
9
11
12#include <vtkCellLocator.h>
13#include <vtkDataArray.h>
14#include <vtkGenericCell.h>
15#include <vtkNew.h>
16#include <vtkPointData.h>
17#include <vtkUnstructuredGrid.h>
18#include <vtkUnstructuredGridReader.h>
19
21
22template <unsigned int dim, typename number>
23class ReadUnstructuredVTK : public ReadFieldBase<dim, number>
24{
25public:
31
36 get_output();
37
41 [[nodiscard]] dealii::types::global_dof_index
42 get_n_points() const;
43
47 [[nodiscard]] dealii::types::global_dof_index
48 get_n_cells() const;
49
53 void
54 print_file() override;
55
59 std::vector<std::string>
61
65 std::vector<std::string>
67
71 number
72 get_scalar_value(const dealii::Point<dim> &point,
73 const std::string &scalar_name) override;
74
78 dealii::Vector<number>
79 get_vector_value(const dealii::Point<dim> &point,
80 const std::string &vector_name) override;
81
82private:
87
91 dealii::types::global_dof_index n_points;
92
96 dealii::types::global_dof_index n_cells;
97
101 unsigned int n_scalars;
102
106 unsigned int n_vectors;
107
111 const unsigned int n_points_per_hex_cell = 8;
112
116 const unsigned int n_space_coordinates = 3;
117};
118
119template <unsigned int dim, typename number>
124{
125 // Create a reader for the vtk file and update it
126 // vtkNew is a smart pointer so we don't need to manage it with delete
128 reader->SetFileName(this->ic_file.filename.c_str());
129 reader->Update();
130
131 // Check that the file is an unstructured grid
132 AssertThrow(reader->IsFileUnstructuredGrid(),
133 dealii::ExcMessage("The vtk file must be an unstructured grid"));
134
135 // Check that we only have one cell type
136 auto *output = reader->GetOutput();
138 output->IsHomogeneous(),
139 dealii::ExcMessage(
140 "The vtk file must have homogeneous cells of type VTK_HEXAHEDRON or VTK_QUAD"));
141
142 // Check that the cells are hexahedra or quads
143 if constexpr (dim == 3)
144 {
145 AssertThrow(output->GetCellType(0) == VTK_HEXAHEDRON,
146 dealii::ExcMessage(
147 "For 3D meshes, the cells must be of type VTK_HEXAHEDRON "));
148 }
149 else if constexpr (dim == 2)
150 {
151 AssertThrow(output->GetCellType(0) == VTK_QUAD,
152 dealii::ExcMessage(
153 "For 2D meshes, the cells must be of type VTK_QUAD"));
154 }
155 else
156 {
157 AssertThrow(false,
158 dealii::ExcMessage("File read-in is not supported for 1D meshes"));
159 }
160
161 // Get the number of points and cells. We first fill the variables in the same type as
162 // the VTK return type so we can check for types mismatches with deal.II
163 const vtkIdType n_points_vtk = output->GetNumberOfPoints();
164 const vtkIdType n_cells_vtk = output->GetNumberOfCells();
165
166 // Check that the number of points and cells are not too large
167 AssertThrow(n_points_vtk < std::numeric_limits<dealii::types::global_dof_index>::max(),
168 dealii::ExcMessage(
169 "The number of points being read-in from the vtk file is too large. Try "
170 "recompiling deal.II with 64-bit indices."));
171 AssertThrow(n_cells_vtk < std::numeric_limits<dealii::types::global_dof_index>::max(),
172 dealii::ExcMessage(
173 "The number of cells being read-in from the vtk file is too large. Try "
174 "recompiling deal.II with 64-bit indices."));
175
176 // Convert the number of points and cells to the dealii type
179
180 // Get the number of scalars and vectors
181 n_scalars = reader->GetNumberOfScalarsInFile();
182 n_vectors = reader->GetNumberOfVectorsInFile();
183}
184
185template <unsigned int dim, typename number>
186inline vtkUnstructuredGrid *
188{
189 return reader->GetOutput();
190}
191
192template <unsigned int dim, typename number>
193inline dealii::types::global_dof_index
195{
196 return n_points;
197}
198
199template <unsigned int dim, typename number>
200inline dealii::types::global_dof_index
202{
203 return n_cells;
204}
205
206template <unsigned int dim, typename number>
207inline void
209{
210 // TODO (landinjm): Should we print only on rank 0?.
211 reader->GetOutput()->PrintSelf(std::cout, vtkIndent());
212}
213
214template <unsigned int dim, typename number>
215inline std::vector<std::string>
217{
218 std::vector<std::string> scalars_names(n_scalars);
219 for (unsigned int i = 0; i < n_scalars; ++i)
220 {
221 scalars_names[i] = reader->GetScalarsNameInFile(static_cast<int>(i));
222 }
223 return scalars_names;
224}
225
226template <unsigned int dim, typename number>
227inline std::vector<std::string>
229{
230 std::vector<std::string> vectors_names(n_vectors);
231 for (unsigned int i = 0; i < n_vectors; ++i)
232 {
233 vectors_names[i] = reader->GetVectorsNameInFile(static_cast<int>(i));
234 }
235 return vectors_names;
236}
237
238template <unsigned int dim, typename number>
239inline number
241 const std::string &scalar_name)
242{
243 // Check that the scalar name is in the vtk file
244 auto scalars_names = get_scalars_names();
245 AssertThrow(std::find(scalars_names.begin(), scalars_names.end(), scalar_name) !=
246 scalars_names.end(),
247 dealii::ExcMessage(
248 "The provided vtk dataset does not contain a field named " +
249 scalar_name));
250
251 // Convert the dealii point to a vector
253
254 // Set the active scalar and update the reader
255 reader->SetScalarsName(scalar_name.c_str());
256 reader->Update();
257
258 // Grab the output and point data
259 vtkUnstructuredGrid *output = reader->GetOutput();
260 vtkPointData *point_data = output->GetPointData();
261
262 // Find the point id in the vtk file
263 const vtkIdType point_id = output->FindPoint(point_vector.data());
264
265 // Check that point is inside the grid
266 AssertThrow(point_id >= 0, dealii::ExcMessage("No matching point found in VTK grid"));
267
268 // Check that the point is within some tolerance to know whether we have to interpolate
269 // or not
270 std::vector<double> point_in_dataset(n_space_coordinates);
271 output->GetPoint(point_id, point_in_dataset.data());
272 bool interpolate = false;
273 for (unsigned int i = 0; i < dim; i++)
274 {
275 // NOLINTBEGIN(cppcoreguidelines-pro-bounds-constant-array-index)
277 {
278 interpolate = true;
279 }
280 // NOLINTEND(cppcoreguidelines-pro-bounds-constant-array-index)
281 }
282
283 // Get the data array
284 vtkDataArray *data_array = point_data->GetArray(scalar_name.c_str());
285 AssertThrow(data_array != nullptr,
286 dealii::ExcMessage(std::string("Data array not found: ") + scalar_name));
287
288 // Get the value of the scalar at the point
289 if (interpolate)
290 {
292 cell_locator->SetDataSet(output);
293 cell_locator->BuildLocator();
294
295 std::vector<double> pcoords(n_space_coordinates);
296 std::vector<double> weights(n_points_per_hex_cell);
297
298 int sub_id = 0;
299
300 vtkGenericCell *cell = vtkGenericCell::New();
301 // NOLINTBEGIN(cppcoreguidelines-pro-bounds-array-to-pointer-decay,hicpp-no-array-decay)
302 const vtkIdType cell_id = cell_locator->FindCell(point_vector.data(),
304 cell,
305 sub_id,
306 pcoords.data(),
307 weights.data());
308 // NOLINTEND(cppcoreguidelines-pro-bounds-array-to-pointer-decay,hicpp-no-array-decay)
309
310 AssertThrow(cell_id >= 0,
311 dealii::ExcMessage("Point not inside any cell for interpolation"));
312
313 // Interpolate scalar value using weights and nodal values
314 vtkIdList *point_ids = output->GetCell(cell_id)->GetPointIds();
316 for (vtkIdType id = 0; id < point_ids->GetNumberOfIds(); ++id)
317 {
318 const vtkIdType pt_id = point_ids->GetId(id);
319 // NOLINTBEGIN(cppcoreguidelines-pro-bounds-constant-array-index)
320 interpolated_value += weights[id] * data_array->GetComponent(pt_id, 0);
321 // NOLINTEND(cppcoreguidelines-pro-bounds-constant-array-index)
322 }
323
324 return interpolated_value;
325 }
326
327 // If we are not interpolating, we can just get the value at the point
328 return data_array->GetComponent(point_id, 0);
329}
330
331template <unsigned int dim, typename number>
332inline dealii::Vector<number>
334 const std::string &vector_name)
335{
336 // Check that the scalar name is in the vtk file
337 auto vectors_names = get_vectors_names();
338 AssertThrow(std::find(vectors_names.begin(), vectors_names.end(), vector_name) !=
339 vectors_names.end(),
340 dealii::ExcMessage(
341 "The provided vtk dataset does not contain a field named " +
342 vector_name));
343
344 // Convert the dealii point to a vector
346
347 // Set the active vector and update the reader
348 reader->SetVectorsName(vector_name.c_str());
349 reader->Update();
350
351 // Grab the output and point data
352 vtkUnstructuredGrid *output = reader->GetOutput();
353 vtkPointData *point_data = output->GetPointData();
354
355 // Find the point id in the vtk file
356 vtkIdType point_id = output->FindPoint(point_vector.data());
357
358 // Check that point is inside the grid
359 AssertThrow(point_id >= 0, dealii::ExcMessage("No matching point found in VTK grid"));
360
361 // Check that the point is within some tolerance to know whether we have to interpolate
362 // or not
363 std::vector<double> point_in_dataset(n_space_coordinates);
364 output->GetPoint(point_id, point_in_dataset.data());
365 bool interpolate = false;
366 for (unsigned int i = 0; i < dim; i++)
367 {
368 // NOLINTBEGIN(cppcoreguidelines-pro-bounds-constant-array-index)
370 {
371 interpolate = true;
372 }
373 // NOLINTEND(cppcoreguidelines-pro-bounds-constant-array-index)
374 }
375
376 // Get the data array
377 vtkDataArray *data_array = point_data->GetArray(vector_name.c_str());
378 AssertThrow(data_array != nullptr,
379 dealii::ExcMessage(std::string("Data array not found: ") + vector_name));
380
381 // Get the value of the vector at the point
382 dealii::Vector<number> vector_value(dim);
383
384 for (unsigned int i = 0; i < dim; i++)
385 {
386 if (interpolate)
387 {
389 cell_locator->SetDataSet(output);
390 cell_locator->BuildLocator();
391
392 std::vector<double> pcoords(n_space_coordinates);
393 std::vector<double> weights(n_points_per_hex_cell);
394
395 int sub_id = 0;
396
397 vtkGenericCell *cell = vtkGenericCell::New();
398 // NOLINTBEGIN(cppcoreguidelines-pro-bounds-array-to-pointer-decay,hicpp-no-array-decay)
399 const vtkIdType cell_id = cell_locator->FindCell(point_vector.data(),
401 cell,
402 sub_id,
403 pcoords.data(),
404 weights.data());
405 // NOLINTEND(cppcoreguidelines-pro-bounds-array-to-pointer-decay,hicpp-no-array-decay)
406
407 AssertThrow(cell_id >= 0,
408 dealii::ExcMessage("Point not inside any cell for interpolation"));
409
410 // Interpolate scalar value using weights and nodal values
411 vtkIdList *point_ids = output->GetCell(cell_id)->GetPointIds();
413 for (vtkIdType id = 0; id < point_ids->GetNumberOfIds(); ++id)
414 {
416 dealii::ExcMessage("Index out of bounds for data array components"));
417 const vtkIdType pt_id = point_ids->GetId(id);
418 // NOLINTBEGIN(cppcoreguidelines-pro-bounds-constant-array-index)
420 weights[id] * data_array->GetComponent(pt_id, static_cast<int>(id));
421 // NOLINTEND(cppcoreguidelines-pro-bounds-constant-array-index)
422 }
423
424 vector_value[i] = interpolated_value;
425 }
426 else
427 {
428 vector_value[i] = data_array->GetComponent(point_id, i);
429 }
430 }
431 return vector_value;
432}
433
Definition read_field_base.h:23
const InitialConditionFile & ic_file
Initial condition file object.
Definition read_field_base.h:86
Definition read_vtk.h:24
dealii::types::global_dof_index get_n_points() const
Get the number of points.
Definition read_vtk.h:194
unsigned int n_scalars
Number of scalars in file.
Definition read_vtk.h:101
number get_scalar_value(const dealii::Point< dim > &point, const std::string &scalar_name) override
Get scalar value for a given point.
Definition read_vtk.h:240
std::vector< std::string > get_scalars_names()
Get the names of the scalars in the vtk file.
Definition read_vtk.h:216
std::vector< std::string > get_vectors_names()
Get the names of the vectors in the vtk file.
Definition read_vtk.h:228
ReadUnstructuredVTK(const InitialConditionFile &_ic_file, const SpatialDiscretization< dim > &_spatial_discretization)
Constructor.
Definition read_vtk.h:120
dealii::types::global_dof_index n_points
Number of points.
Definition read_vtk.h:91
vtkUnstructuredGrid * get_output()
Get the vtk output.
Definition read_vtk.h:187
dealii::types::global_dof_index get_n_cells() const
Get the number of cells.
Definition read_vtk.h:201
dealii::types::global_dof_index n_cells
Number of cells.
Definition read_vtk.h:96
void print_file() override
Print the vtk file for debugging.
Definition read_vtk.h:208
vtkNew< vtkUnstructuredGridReader > reader
Reader for the vtk file.
Definition read_vtk.h:86
const unsigned int n_points_per_hex_cell
Number of points in a hex cell.
Definition read_vtk.h:111
const unsigned int n_space_coordinates
Number of space coordinates in a point.
Definition read_vtk.h:116
dealii::Vector< number > get_vector_value(const dealii::Point< dim > &point, const std::string &vector_name) override
Get vector value for a given point.
Definition read_vtk.h:333
unsigned int n_vectors
Number of vectors in file.
Definition read_vtk.h:106
@ Value
Use value of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:31
static const double mesh_tolerance
Default mesh tolerance.
Definition types.h:62
Definition conditional_ostreams.cc:20
Struct that store the read-in information for a single file.
Definition load_initial_condition_parameters.h:25
std::string filename
Definition load_initial_condition_parameters.h:27
Struct that holds spatial discretization parameters.
Definition spatial_discretization.h:230