16 #include <deal.II/base/work_stream.h> 18 #include <deal.II/dofs/dof_accessor.h> 19 #include <deal.II/dofs/dof_handler.h> 21 #include <deal.II/fe/fe.h> 22 #include <deal.II/fe/fe_dgq.h> 23 #include <deal.II/fe/fe_values.h> 24 #include <deal.II/fe/mapping_q1.h> 26 #include <deal.II/grid/tria.h> 27 #include <deal.II/grid/tria_iterator.h> 29 #include <deal.II/hp/dof_handler.h> 30 #include <deal.II/hp/fe_values.h> 32 #include <deal.II/numerics/data_out.h> 36 DEAL_II_NAMESPACE_OPEN
41 namespace DataOutImplementation
43 template <
int dim,
int spacedim>
44 ParallelData<dim, spacedim>::ParallelData(
45 const unsigned int n_datasets,
46 const unsigned int n_subdivisions,
47 const std::vector<unsigned int> &n_postprocessor_outputs,
53 const std::vector<std::vector<unsigned int>> &cell_to_patch_index_map)
54 : ParallelDataBase<dim, spacedim>(n_datasets,
56 n_postprocessor_outputs,
61 , cell_to_patch_index_map(&cell_to_patch_index_map)
68 template <
int dim,
typename DoFHandlerType>
71 const std::pair<cell_iterator, unsigned int> *cell_and_index,
73 DoFHandlerType::space_dimension>
75 const unsigned int n_subdivisions,
80 DoFHandlerType::space_dimension>
87 for (
unsigned int vertex = 0;
88 vertex < GeometryInfo<DoFHandlerType::dimension>::vertices_per_cell;
90 if (scratch_data.mapping_collection[0].preserves_vertex_locations())
91 patch.vertices[vertex] = cell_and_index->first->vertex(vertex);
93 patch.vertices[vertex] =
94 scratch_data.mapping_collection[0].transform_unit_to_real_cell(
95 cell_and_index->first,
99 scratch_data.reinit_all_fe_values(this->dof_data, cell_and_index->first);
102 &fe_patch_values = scratch_data.get_present_fe_values(0);
113 if (curved_cell_region == curved_inner_cells ||
114 (curved_cell_region == curved_boundary &&
115 (cell_and_index->first->at_boundary() ||
116 (DoFHandlerType::dimension != DoFHandlerType::space_dimension))))
118 Assert(patch.space_dim == DoFHandlerType::space_dimension,
123 patch.points_are_available =
true;
127 const std::vector<Point<DoFHandlerType::space_dimension>> &q_points =
130 patch.data.reinit(scratch_data.n_datasets +
131 DoFHandlerType::space_dimension,
133 for (
unsigned int i = 0; i < DoFHandlerType::space_dimension; ++i)
134 for (
unsigned int q = 0; q < n_q_points; ++q)
135 patch.data(patch.data.size(0) - DoFHandlerType::space_dimension + i,
140 patch.data.reinit(scratch_data.n_datasets, n_q_points);
141 patch.points_are_available =
false;
145 if (scratch_data.n_datasets > 0)
148 unsigned int offset = 0;
151 for (
unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
154 DoFHandlerType::space_dimension>
155 &this_fe_patch_values = scratch_data.get_present_fe_values(dataset);
157 this_fe_patch_values.get_fe().n_components();
160 *postprocessor = this->dof_data[dataset]->postprocessor;
162 if (postprocessor !=
nullptr)
173 this->dof_data[dataset]->get_function_values(
174 this_fe_patch_values,
175 internal::DataOutImplementation::ComponentExtractor::
177 scratch_data.patch_values_scalar.solution_values);
179 this->dof_data[dataset]->get_function_gradients(
180 this_fe_patch_values,
181 internal::DataOutImplementation::ComponentExtractor::
183 scratch_data.patch_values_scalar.solution_gradients);
185 this->dof_data[dataset]->get_function_hessians(
186 this_fe_patch_values,
187 internal::DataOutImplementation::ComponentExtractor::
189 scratch_data.patch_values_scalar.solution_hessians);
192 scratch_data.patch_values_scalar.evaluation_points =
193 this_fe_patch_values.get_quadrature_points();
195 const typename DoFHandlerType::active_cell_iterator dh_cell(
196 &cell_and_index->first->get_triangulation(),
197 cell_and_index->first->level(),
198 cell_and_index->first->index(),
199 this->dof_data[dataset]->dof_handler);
200 scratch_data.patch_values_scalar
201 .template set_cell<DoFHandlerType>(dh_cell);
204 scratch_data.patch_values_scalar,
205 scratch_data.postprocessed_values[dataset]);
214 this->dof_data[dataset]->get_function_values(
215 this_fe_patch_values,
216 internal::DataOutImplementation::ComponentExtractor::
218 scratch_data.patch_values_system.solution_values);
220 this->dof_data[dataset]->get_function_gradients(
221 this_fe_patch_values,
222 internal::DataOutImplementation::ComponentExtractor::
224 scratch_data.patch_values_system.solution_gradients);
226 this->dof_data[dataset]->get_function_hessians(
227 this_fe_patch_values,
228 internal::DataOutImplementation::ComponentExtractor::
230 scratch_data.patch_values_system.solution_hessians);
233 scratch_data.patch_values_system.evaluation_points =
234 this_fe_patch_values.get_quadrature_points();
236 const typename DoFHandlerType::active_cell_iterator dh_cell(
237 &cell_and_index->first->get_triangulation(),
238 cell_and_index->first->level(),
239 cell_and_index->first->index(),
240 this->dof_data[dataset]->dof_handler);
241 scratch_data.patch_values_system
242 .template set_cell<DoFHandlerType>(dh_cell);
245 scratch_data.patch_values_system,
246 scratch_data.postprocessed_values[dataset]);
249 for (
unsigned int q = 0; q < n_q_points; ++q)
250 for (
unsigned int component = 0;
251 component < this->dof_data[dataset]->n_output_variables;
253 patch.data(offset + component, q) =
254 scratch_data.postprocessed_values[dataset][q](component);
263 this->dof_data[dataset]->get_function_values(
264 this_fe_patch_values,
265 internal::DataOutImplementation::ComponentExtractor::real_part,
266 scratch_data.patch_values_scalar.solution_values);
267 for (
unsigned int q = 0; q < n_q_points; ++q)
268 patch.data(offset, q) =
269 scratch_data.patch_values_scalar.solution_values[q];
272 if (this->dof_data[dataset]->is_complex_valued() ==
true)
274 this->dof_data[dataset]->get_function_values(
275 this_fe_patch_values,
276 internal::DataOutImplementation::ComponentExtractor::
278 scratch_data.patch_values_scalar.solution_values);
279 for (
unsigned int q = 0; q < n_q_points; ++q)
280 patch.data(offset + 1, q) =
281 scratch_data.patch_values_scalar.solution_values[q];
289 const unsigned int stride =
290 (this->dof_data[dataset]->is_complex_valued() ? 2 : 1);
291 this->dof_data[dataset]->get_function_values(
292 this_fe_patch_values,
293 internal::DataOutImplementation::ComponentExtractor::real_part,
294 scratch_data.patch_values_system.solution_values);
295 for (
unsigned int component = 0; component <
n_components;
297 for (
unsigned int q = 0; q < n_q_points; ++q)
298 patch.data(offset + component * stride, q) =
299 scratch_data.patch_values_system.solution_values[q](
303 if (this->dof_data[dataset]->is_complex_valued() ==
true)
305 this->dof_data[dataset]->get_function_values(
306 this_fe_patch_values,
307 internal::DataOutImplementation::ComponentExtractor::
309 scratch_data.patch_values_system.solution_values);
310 for (
unsigned int component = 0; component <
n_components;
312 for (
unsigned int q = 0; q < n_q_points; ++q)
313 patch.data(offset + component * stride + 1, q) =
314 scratch_data.patch_values_system.solution_values[q](
320 offset += this->dof_data[dataset]->n_output_variables *
321 (this->dof_data[dataset]->is_complex_valued() ? 2 : 1);
327 if (this->cell_data.size() != 0)
331 for (
unsigned int dataset = 0; dataset < this->cell_data.size();
337 this->cell_data[dataset]->get_cell_data_value(
338 cell_and_index->second,
339 internal::DataOutImplementation::ComponentExtractor::
341 for (
unsigned int q = 0; q < n_q_points; ++q)
342 patch.data(offset, q) = value;
346 if (this->cell_data[dataset]->is_complex_valued() ==
true)
349 this->cell_data[dataset]->get_cell_data_value(
350 cell_and_index->second,
351 internal::DataOutImplementation::ComponentExtractor::
353 for (
unsigned int q = 0; q < n_q_points; ++q)
354 patch.data(offset + 1, q) = value;
357 offset += (this->cell_data[dataset]->is_complex_valued() ? 2 : 1);
363 for (
unsigned int f = 0;
364 f < GeometryInfo<DoFHandlerType::dimension>::faces_per_cell;
377 if (cell_and_index->first->at_boundary(f) ||
378 (cell_and_index->first->neighbor(f)->level() !=
379 cell_and_index->first->level()))
385 const cell_iterator neighbor = cell_and_index->first->neighbor(f);
386 Assert(static_cast<unsigned int>(neighbor->level()) <
387 scratch_data.cell_to_patch_index_map->size(),
389 if ((static_cast<unsigned int>(neighbor->index()) >=
390 (*scratch_data.cell_to_patch_index_map)[neighbor->level()].size()) ||
391 ((*scratch_data.cell_to_patch_index_map)[neighbor->level()]
392 [neighbor->index()] ==
403 .cell_to_patch_index_map)[neighbor->level()][neighbor->index()];
406 const unsigned int patch_idx =
407 (*scratch_data.cell_to_patch_index_map)[cell_and_index->first->level()]
408 [cell_and_index->first->index()];
411 patch.patch_index = patch_idx;
416 this->patches[patch_idx].swap(patch);
421 template <
int dim,
typename DoFHandlerType>
426 DoFHandlerType::space_dimension>::mapping,
433 template <
int dim,
typename DoFHandlerType>
438 const unsigned int n_subdivisions_,
442 Assert(dim == DoFHandlerType::dimension,
445 Assert(this->triangulation !=
nullptr,
448 const unsigned int n_subdivisions =
449 (n_subdivisions_ != 0) ? n_subdivisions_ : this->default_subdivisions;
450 Assert(n_subdivisions >= 1,
454 this->validate_dataset_names();
472 std::vector<std::vector<unsigned int>> cell_to_patch_index_map;
473 cell_to_patch_index_map.resize(this->triangulation->n_levels());
474 for (
unsigned int l = 0; l < this->triangulation->n_levels(); ++l)
477 unsigned int max_index = 0;
479 cell != this->triangulation->end();
480 cell = next_locally_owned_cell(cell))
481 if (static_cast<unsigned int>(cell->level()) == l)
483 std::max(max_index, static_cast<unsigned int>(cell->index()));
485 cell_to_patch_index_map[l].resize(
488 DoFHandlerType::dimension,
489 DoFHandlerType::space_dimension>::no_neighbor);
493 std::vector<std::pair<cell_iterator, unsigned int>> all_cells;
500 active_cell_iterator active_cell = this->triangulation->begin_active();
501 unsigned int active_index = 0;
503 for (; cell != this->triangulation->end();
504 cell = next_locally_owned_cell(cell))
508 while (active_cell != this->triangulation->end() && cell->active() &&
509 active_cell_iterator(cell) != active_cell)
515 Assert(static_cast<unsigned int>(cell->level()) <
516 cell_to_patch_index_map.size(),
518 Assert(static_cast<unsigned int>(cell->index()) <
519 cell_to_patch_index_map[cell->level()].size(),
521 Assert(active_index < this->triangulation->n_active_cells(),
523 cell_to_patch_index_map[cell->level()][cell->index()] =
526 all_cells.emplace_back(cell, active_index);
530 this->patches.clear();
531 this->patches.resize(all_cells.size());
534 unsigned int n_datasets = 0;
535 for (
unsigned int i = 0; i < this->cell_data.size(); ++i)
536 n_datasets += (this->cell_data[i]->is_complex_valued() ? 2 : 1);
537 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
538 n_datasets += (this->dof_data[i]->n_output_variables *
539 (this->dof_data[i]->is_complex_valued() ? 2 : 1));
541 std::vector<unsigned int> n_postprocessor_outputs(this->dof_data.size());
542 for (
unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
543 if (this->dof_data[dataset]->postprocessor)
544 n_postprocessor_outputs[dataset] =
545 this->dof_data[dataset]->n_output_variables;
547 n_postprocessor_outputs[dataset] = 0;
550 (n_subdivisions < 2 ? no_curved_cells : curved_region);
553 if (curved_cell_region != no_curved_cells)
556 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
557 if (this->dof_data[i]->postprocessor)
559 this->dof_data[i]->postprocessor->get_needed_update_flags();
565 "The update of normal vectors may not be requested for evaluation of " 566 "data on cells via DataPostprocessor."));
569 DoFHandlerType::space_dimension>
570 thread_data(n_datasets,
572 n_postprocessor_outputs,
576 cell_to_patch_index_map);
579 if (all_cells.size() > 0)
582 all_cells.data() + all_cells.size(),
585 std::placeholders::_1,
586 std::placeholders::_2,
593 std::function<void(const int)>(),
608 template <
int dim,
typename DoFHandlerType>
612 return this->triangulation->begin_active();
617 template <
int dim,
typename DoFHandlerType>
625 DoFHandlerType::space_dimension>::active_cell_iterator
633 template <
int dim,
typename DoFHandlerType>
641 while ((cell != this->triangulation->end()) &&
642 (cell->has_children() ==
false) && !cell->is_locally_owned())
643 cell = next_cell(cell);
650 template <
int dim,
typename DoFHandlerType>
657 while ((cell != this->triangulation->end()) &&
658 (cell->has_children() ==
false) && !cell->is_locally_owned())
659 cell = next_cell(cell);
665 #include "data_out.inst" 667 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
static ::ExceptionBase & ExcNoTriangulationSelected()
void build_one_patch(const std::pair< cell_iterator, unsigned int > *cell_and_index, internal::DataOutImplementation::ParallelData< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &scratch_data, const unsigned int n_subdivisions, const CurvedCellRegion curved_cell_region)
virtual cell_iterator next_cell(const cell_iterator &cell)
typename DataOut_DoFData< DoFHandler< dim >, DoFHandler< dim > ::dimension, DoFHandler< dim > ::space_dimension >::cell_iterator cell_iterator
virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double >> &computed_quantities) const
Transformed quadrature points.
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
virtual void build_patches(const unsigned int n_subdivisions=0)
const std::vector< Point< spacedim > > & get_quadrature_points() const
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual cell_iterator first_locally_owned_cell()
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
virtual cell_iterator first_cell()
unsigned int n_subdivisions
Second derivatives of shape functions.
const unsigned int n_quadrature_points
Shape function gradients.
virtual void evaluate_scalar_field(const DataPostprocessorInputs::Scalar< dim > &input_data, std::vector< Vector< double >> &computed_quantities) const
static ::ExceptionBase & ExcNotImplemented()
void run(const std::vector< std::vector< Iterator >> &colored_iterators, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
static unsigned int n_threads()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual cell_iterator next_locally_owned_cell(const cell_iterator &cell)
static ::ExceptionBase & ExcInternalError()
virtual UpdateFlags get_needed_update_flags() const =0