25# include <vtkCellData.h>
26# include <vtkCleanUnstructuredGrid.h>
27# include <vtkDataArray.h>
28# include <vtkPointData.h>
29# include <vtkSmartPointer.h>
30# include <vtkUnstructuredGrid.h>
31# include <vtkUnstructuredGridReader.h>
37 template <
int dim,
int spacedim>
39 read_vtk(
const std::string &vtk_filename,
40 Triangulation<dim, spacedim> &tria,
43 auto reader = vtkSmartPointer<vtkUnstructuredGridReader>::New();
44 reader->SetFileName(vtk_filename.c_str());
46 vtkUnstructuredGrid *grid = reader->GetOutput();
47 AssertThrow(grid, ExcMessage(
"Failed to read VTK file: " + vtk_filename));
49 auto cleaner = vtkSmartPointer<vtkCleanUnstructuredGrid>::New();
54 cleaner->SetInputData(grid);
56 grid = cleaner->GetOutput();
58 ExcMessage(
"Failed to clean VTK file: " + vtk_filename));
62 vtkPoints *vtk_points = grid->GetPoints();
63 const vtkIdType n_points = vtk_points->GetNumberOfPoints();
64 std::vector<Point<spacedim>> points(n_points);
65 for (vtkIdType i = 0; i < n_points; ++i)
67 std::array<double, 3> coords = {{0, 0, 0}};
68 vtk_points->GetPoint(i, coords.data());
69 for (
unsigned int d = 0;
d < spacedim; ++
d)
70 points[i][d] = coords[d];
74 std::vector<CellData<dim>> cells;
75 const vtkIdType
n_cells = grid->GetNumberOfCells();
76 for (vtkIdType i = 0; i <
n_cells; ++i)
78 vtkCell *cell = grid->GetCell(i);
79 if constexpr (dim == 1)
81 if (cell->GetCellType() != VTK_LINE)
84 "Unsupported cell type in 1D VTK file: only "
85 "VTK_LINE is supported."));
88 "Only line cells with 2 points are supported."));
89 CellData<1> cell_data;
90 for (
unsigned int j = 0; j < 2; ++j)
91 cell_data.
vertices[j] = cell->GetPointId(j);
93 cells.push_back(cell_data);
95 else if constexpr (dim == 2)
97 if (cell->GetCellType() == VTK_QUAD)
101 "Only quad cells with 4 points are supported."));
102 CellData<2> cell_data;
103 for (
unsigned int j = 0; j < 4; ++j)
104 cell_data.
vertices[j] = cell->GetPointId(j);
106 cells.push_back(cell_data);
108 else if (cell->GetCellType() == VTK_TRIANGLE)
111 cell->GetNumberOfPoints() == 3,
113 "Only triangle cells with 3 points are supported."));
114 CellData<2> cell_data;
115 for (
unsigned int j = 0; j < 3; ++j)
116 cell_data.
vertices[j] = cell->GetPointId(j);
118 cells.push_back(cell_data);
123 "Unsupported cell type in 2D VTK file: only "
124 "VTK_QUAD and VTK_TRIANGLE are supported."));
126 else if constexpr (dim == 3)
128 if (cell->GetCellType() == VTK_HEXAHEDRON)
132 "Only hex cells with 8 points are supported."));
133 CellData<3> cell_data;
134 for (
unsigned int j = 0; j < 8; ++j)
135 cell_data.
vertices[j] = cell->GetPointId(j);
140 cells.push_back(cell_data);
142 else if (cell->GetCellType() == VTK_TETRA)
145 cell->GetNumberOfPoints() == 4,
147 "Only tetrahedron cells with 4 points are supported."));
148 CellData<3> cell_data;
149 for (
unsigned int j = 0; j < 4; ++j)
150 cell_data.
vertices[j] = cell->GetPointId(j);
152 cells.push_back(cell_data);
154 else if (cell->GetCellType() == VTK_WEDGE)
158 "Only prism cells with 6 points are supported."));
159 CellData<3> cell_data;
160 for (
unsigned int j = 0; j < 6; ++j)
161 cell_data.
vertices[j] = cell->GetPointId(j);
163 cells.push_back(cell_data);
165 else if (cell->GetCellType() == VTK_PYRAMID)
168 cell->GetNumberOfPoints() == 5,
170 "Only pyramid cells with 5 points are supported."));
171 CellData<3> cell_data;
172 for (
unsigned int j = 0; j < 5; ++j)
173 cell_data.
vertices[j] = cell->GetPointId(j);
175 cells.push_back(cell_data);
181 "Unsupported cell type in 3D VTK file: only "
182 "VTK_HEXAHEDRON, VTK_TETRA, VTK_WEDGE, and VTK_PYRAMID are supported."));
186 AssertThrow(
false, ExcMessage(
"Unsupported dimension."));
195 read_cell_data(
const std::string &vtk_filename,
196 const std::string &cell_data_name,
197 Vector<double> &output_vector)
199 auto reader = vtkSmartPointer<vtkUnstructuredGridReader>::New();
200 reader->SetFileName(vtk_filename.c_str());
202 vtkUnstructuredGrid *grid = reader->GetOutput();
203 AssertThrow(grid, ExcMessage(
"Failed to read VTK file: " + vtk_filename));
204 vtkDataArray *data_array =
205 grid->GetCellData()->GetArray(cell_data_name.c_str());
207 ExcMessage(
"Cell data array '" + cell_data_name +
208 "' not found in VTK file: " + vtk_filename));
209 vtkIdType n_tuples = data_array->GetNumberOfTuples();
210 int n_components = data_array->GetNumberOfComponents();
211 output_vector.
reinit(n_tuples * n_components);
212 for (vtkIdType i = 0; i < n_tuples; ++i)
213 for (
int j = 0; j < n_components; ++j)
214 output_vector[i * n_components + j] = data_array->GetComponent(i, j);
218 read_vertex_data(
const std::string &vtk_filename,
219 const std::string &point_data_name,
220 Vector<double> &output_vector)
222 auto reader = vtkSmartPointer<vtkUnstructuredGridReader>::New();
223 reader->SetFileName(vtk_filename.c_str());
225 vtkUnstructuredGrid *grid = reader->GetOutput();
226 AssertThrow(grid, ExcMessage(
"Failed to read VTK file: " + vtk_filename));
227 vtkDataArray *data_array =
228 grid->GetPointData()->GetArray(point_data_name.c_str());
230 ExcMessage(
"Point data array '" + point_data_name +
231 "' not found in VTK file: " + vtk_filename));
232 vtkIdType n_tuples = data_array->GetNumberOfTuples();
233 int n_components = data_array->GetNumberOfComponents();
234 output_vector.
reinit(n_tuples * n_components);
235 for (vtkIdType i = 0; i < n_tuples; ++i)
236 for (
int j = 0; j < n_components; ++j)
237 output_vector[i * n_components + j] = data_array->GetComponent(i, j);
241 read_data(
const std::string &vtk_filename, Vector<double> &output_vector)
243 auto reader = vtkSmartPointer<vtkUnstructuredGridReader>::New();
244 reader->SetFileName(vtk_filename.c_str());
246 vtkUnstructuredGrid *grid = reader->GetOutput();
247 AssertThrow(grid, ExcMessage(
"Failed to read VTK file: " + vtk_filename));
249 std::vector<double>
data;
251 vtkPointData *point_data = grid->GetPointData();
254 for (
int i = 0; i < point_data->GetNumberOfArrays(); ++i)
256 vtkDataArray *data_array = point_data->GetArray(i);
259 vtkIdType n_tuples = data_array->GetNumberOfTuples();
260 int n_components = data_array->GetNumberOfComponents();
261 unsigned int current_size =
data.size();
262 data.resize(current_size + n_tuples * n_components, 0.0);
263 for (vtkIdType tuple_idx = 0; tuple_idx < n_tuples; ++tuple_idx)
264 for (
int comp_idx = 0; comp_idx < n_components; ++comp_idx)
265 data[current_size + tuple_idx * n_components + comp_idx] =
266 data_array->GetComponent(tuple_idx, comp_idx);
270 vtkCellData *cell_data = grid->GetCellData();
273 for (
int i = 0; i < cell_data->GetNumberOfArrays(); ++i)
275 vtkDataArray *data_array = cell_data->GetArray(i);
278 vtkIdType n_tuples = data_array->GetNumberOfTuples();
279 int n_components = data_array->GetNumberOfComponents();
280 unsigned int current_size =
data.size();
281 data.resize(current_size + n_tuples * n_components,
true);
282 for (vtkIdType tuple_idx = 0; tuple_idx < n_tuples; ++tuple_idx)
283 for (
int comp_idx = 0; comp_idx < n_components; ++comp_idx)
284 data[current_size + tuple_idx * n_components + comp_idx] =
285 data_array->GetComponent(tuple_idx, comp_idx);
292 template <
int dim,
int spacedim>
293 std::pair<std::unique_ptr<FiniteElement<dim, spacedim>>,
294 std::vector<std::string>>
295 vtk_to_finite_element(
const std::string &vtk_filename)
297 std::vector<std::string> data_names;
298 auto reader = vtkSmartPointer<vtkUnstructuredGridReader>::New();
299 reader->SetFileName(vtk_filename.c_str());
301 vtkUnstructuredGrid *grid = reader->GetOutput();
302 AssertThrow(grid, ExcMessage(
"Failed to read VTK file: " + vtk_filename));
304 vtkCellData *cell_data = grid->GetCellData();
305 vtkPointData *point_data = grid->GetPointData();
307 std::vector<std::shared_ptr<FiniteElement<dim, spacedim>>> fe_collection;
308 std::vector<unsigned int> n_components_collection;
311 for (
int i = 0; i < point_data->GetNumberOfArrays(); ++i)
313 vtkDataArray *arr = point_data->GetArray(i);
316 std::string name = arr->GetName();
317 int n_comp = arr->GetNumberOfComponents();
319 fe_collection.push_back(std::make_shared<FE_Q<dim, spacedim>>(1));
322 fe_collection.push_back(
323 std::make_shared<FESystem<dim, spacedim>>(FE_Q<dim, spacedim>(1),
325 n_components_collection.push_back(n_comp);
326 data_names.push_back(name);
330 for (
int i = 0; i < cell_data->GetNumberOfArrays(); ++i)
332 vtkDataArray *arr = cell_data->GetArray(i);
335 std::string name = arr->GetName();
336 int n_comp = arr->GetNumberOfComponents();
338 fe_collection.push_back(std::make_shared<FE_DGQ<dim, spacedim>>(0));
341 fe_collection.push_back(
342 std::make_shared<FESystem<dim, spacedim>>(FE_DGQ<dim, spacedim>(0),
344 n_components_collection.push_back(n_comp);
345 data_names.push_back(name);
350 std::vector<const FiniteElement<dim, spacedim> *> fe_ptrs;
351 std::vector<unsigned int> multiplicities;
352 for (
const auto &fe : fe_collection)
354 fe_ptrs.push_back(fe.get());
355 multiplicities.push_back(1);
358 return std::make_pair(std::make_unique<FE_Nothing<dim, spacedim>>(),
359 std::vector<std::string>());
362 return std::make_pair(
363 std::make_unique<FESystem<dim, spacedim>>(fe_ptrs, multiplicities),
370 template <
int dim,
int spacedim>
372 read_vtk(
const std::string &vtk_filename,
373 DoFHandler<dim, spacedim> &dof_handler,
374 Vector<double> &output_vector,
375 std::vector<std::string> &data_names)
378 auto &tria =
const_cast<Triangulation<dim, spacedim> &
>(
383 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *
>(&tria);
386 "The input triangulation must be a serial triangulation."));
391 read_vtk(vtk_filename, tria,
true);
393 Vector<double> raw_data_vector;
394 read_data(vtk_filename, raw_data_vector);
396 auto [fe, data_names_from_fe] =
397 vtk_to_finite_element<dim, spacedim>(vtk_filename);
401 data_to_dealii_vector(tria, raw_data_vector, dof_handler, output_vector);
405 data_names = data_names_from_fe;
408 template <
int dim,
int spacedim>
410 serial_vector_to_distributed_vector(
411 const DoFHandler<dim, spacedim> &serial_dh,
412 const DoFHandler<dim, spacedim> ¶llel_dh,
413 const Vector<double> &serial_vec,
414 LinearAlgebra::distributed::Vector<double> ¶llel_vec)
422 ExcMessage(
"The finite element systems of the serial and "
423 "parallel DoFHandlers must be the same."));
425 std::vector<types::global_dof_index> serial_dof_indices(
427 std::vector<types::global_dof_index> parallel_dof_indices(
433 for (; parallel_cell != parallel_dh.
end(); ++parallel_cell)
434 if (parallel_cell->is_locally_owned())
438 while (serial_cell->id() < parallel_cell->id())
440 serial_cell->get_dof_indices(serial_dof_indices);
441 parallel_cell->get_dof_indices(parallel_dof_indices);
442 unsigned int serial_index = 0;
443 for (
const auto &i : parallel_dof_indices)
448 serial_vec[serial_dof_indices[serial_index]];
456 template <
int dim,
int spacedim>
457 std::vector<types::global_vertex_index>
458 distributed_to_serial_vertex_indices(
459 const Triangulation<dim, spacedim> &serial_tria,
460 const Triangulation<dim, spacedim> ¶llel_tria)
462 const auto locally_owned_indices =
464 std::vector<types::global_vertex_index>
465 distributed_to_serial_vertex_indices(parallel_tria.
n_vertices(),
471 for (; parallel_cell != parallel_tria.
end(); ++parallel_cell)
472 if (parallel_cell->is_locally_owned())
476 while (serial_cell->id() < parallel_cell->id())
478 for (
const unsigned int &v : serial_cell->vertex_indices())
480 const auto serial_index = serial_cell->vertex_index(v);
481 const auto parallel_index = parallel_cell->vertex_index(v);
482 if (locally_owned_indices[parallel_index])
483 distributed_to_serial_vertex_indices[parallel_index] =
487 return distributed_to_serial_vertex_indices;
507template std::pair<std::unique_ptr<FiniteElement<1, 1>>,
508 std::vector<std::string>>
509VTKUtils::vtk_to_finite_element(
const std::string &);
510template std::pair<std::unique_ptr<FiniteElement<1, 2>>,
511 std::vector<std::string>>
512VTKUtils::vtk_to_finite_element(
const std::string &);
513template std::pair<std::unique_ptr<FiniteElement<1, 3>>,
514 std::vector<std::string>>
515VTKUtils::vtk_to_finite_element(
const std::string &);
516template std::pair<std::unique_ptr<FiniteElement<2, 2>>,
517 std::vector<std::string>>
518VTKUtils::vtk_to_finite_element(
const std::string &);
519template std::pair<std::unique_ptr<FiniteElement<2, 3>>,
520 std::vector<std::string>>
521VTKUtils::vtk_to_finite_element(
const std::string &);
522template std::pair<std::unique_ptr<FiniteElement<3, 3>>,
523 std::vector<std::string>>
524VTKUtils::vtk_to_finite_element(
const std::string &);
527VTKUtils::read_vtk(
const std::string &,
530 std::vector<std::string> &);
532VTKUtils::read_vtk(
const std::string &,
535 std::vector<std::string> &);
537VTKUtils::read_vtk(
const std::string &,
540 std::vector<std::string> &);
542VTKUtils::read_vtk(
const std::string &,
545 std::vector<std::string> &);
547VTKUtils::read_vtk(
const std::string &,
550 std::vector<std::string> &);
552VTKUtils::read_vtk(
const std::string &,
555 std::vector<std::string> &);
558VTKUtils::serial_vector_to_distributed_vector(
564VTKUtils::serial_vector_to_distributed_vector(
570VTKUtils::serial_vector_to_distributed_vector(
576VTKUtils::serial_vector_to_distributed_vector(
582VTKUtils::serial_vector_to_distributed_vector(
588VTKUtils::serial_vector_to_distributed_vector(
594template std::vector<types::global_vertex_index>
597template std::vector<types::global_vertex_index>
600template std::vector<types::global_vertex_index>
603template std::vector<types::global_vertex_index>
606template std::vector<types::global_vertex_index>
609template std::vector<types::global_vertex_index>
cell_iterator end() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
void distribute_dofs(const FiniteElement< dim, spacedim > &fe)
const Triangulation< dim, spacedim > & get_triangulation() const
active_cell_iterator begin_active(const unsigned int level=0) const
types::global_dof_index n_dofs() const
unsigned int n_dofs_per_cell() const
unsigned int n_blocks() const
void compress(VectorOperation::values operation)
virtual size_type size() const override
bool in_local_range(const size_type global_index) const
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
cell_iterator end() const
unsigned int n_vertices() const
active_cell_iterator begin_active(const unsigned int level=0) const
virtual size_type size() const override
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
#define AssertDimension(dim1, dim2)
#define AssertThrow(cond, exc)
std::vector< index_type > data
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
constexpr unsigned int invalid_unsigned_int
std::vector< unsigned int > vertices
types::material_id material_id