6#include <deal.II/base/exceptions.h>
7#include <deal.II/base/point.h>
8#include <deal.II/lac/vector.h>
12#include <vtkCellLocator.h>
13#include <vtkDataArray.h>
14#include <vtkGenericCell.h>
16#include <vtkPointData.h>
17#include <vtkUnstructuredGrid.h>
18#include <vtkUnstructuredGridReader.h>
22template <
unsigned int dim,
typename number>
41 [[
nodiscard]] dealii::types::global_dof_index
47 [[
nodiscard]] dealii::types::global_dof_index
59 std::vector<std::string>
65 std::vector<std::string>
78 dealii::Vector<number>
119template <
unsigned int dim,
typename number>
133 dealii::ExcMessage(
"The vtk file must be an unstructured grid"));
140 "The vtk file must have homogeneous cells of type VTK_HEXAHEDRON or VTK_QUAD"));
143 if constexpr (
dim == 3)
147 "For 3D meshes, the cells must be of type VTK_HEXAHEDRON "));
149 else if constexpr (
dim == 2)
153 "For 2D meshes, the cells must be of type VTK_QUAD"));
158 dealii::ExcMessage(
"File read-in is not supported for 1D meshes"));
169 "The number of points being read-in from the vtk file is too large. Try "
170 "recompiling deal.II with 64-bit indices."));
173 "The number of cells being read-in from the vtk file is too large. Try "
174 "recompiling deal.II with 64-bit indices."));
185template <
unsigned int dim,
typename number>
189 return reader->GetOutput();
192template <
unsigned int dim,
typename number>
193inline dealii::types::global_dof_index
199template <
unsigned int dim,
typename number>
200inline dealii::types::global_dof_index
206template <
unsigned int dim,
typename number>
211 reader->GetOutput()->PrintSelf(std::cout,
vtkIndent());
214template <
unsigned int dim,
typename number>
215inline std::vector<std::string>
219 for (
unsigned int i = 0;
i < n_scalars; ++
i)
226template <
unsigned int dim,
typename number>
227inline std::vector<std::string>
231 for (
unsigned int i = 0;
i < n_vectors; ++
i)
238template <
unsigned int dim,
typename number>
248 "The provided vtk dataset does not contain a field named " +
272 bool interpolate =
false;
273 for (
unsigned int i = 0;
i <
dim;
i++)
286 dealii::ExcMessage(std::string(
"Data array not found: ") +
scalar_name));
295 std::vector<double>
pcoords(n_space_coordinates);
296 std::vector<double>
weights(n_points_per_hex_cell);
311 dealii::ExcMessage(
"Point not inside any cell for interpolation"));
331template <
unsigned int dim,
typename number>
332inline dealii::Vector<number>
341 "The provided vtk dataset does not contain a field named " +
365 bool interpolate =
false;
366 for (
unsigned int i = 0;
i <
dim;
i++)
379 dealii::ExcMessage(std::string(
"Data array not found: ") +
vector_name));
382 dealii::Vector<number> vector_value(
dim);
384 for (
unsigned int i = 0;
i <
dim;
i++)
392 std::vector<double>
pcoords(n_space_coordinates);
393 std::vector<double>
weights(n_points_per_hex_cell);
408 dealii::ExcMessage(
"Point not inside any cell for interpolation"));
416 dealii::ExcMessage(
"Index out of bounds for data array components"));
Definition read_field_base.h:23
const InitialConditionFile & ic_file
Initial condition file object.
Definition read_field_base.h:86
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