41 namespace DataOutImplementation
43 template <
int dim,
int spacedim>
45 const unsigned int n_datasets,
46 const unsigned int n_subdivisions,
47 const std::vector<unsigned int> &n_postprocessor_outputs,
48 const ::hp::MappingCollection<dim, spacedim> &mapping,
53 const std::vector<std::vector<unsigned int>> &cell_to_patch_index_map)
56 n_postprocessor_outputs,
61 , cell_to_patch_index_map(&cell_to_patch_index_map)
68template <
int dim,
int spacedim>
78 while ((cell != this->
triangulation->end()) && !cell->is_locally_owned())
87 while ((cell != this->
triangulation->end()) && !cell->is_locally_owned())
95template <
int dim,
int spacedim>
98 const std::pair<cell_iterator, unsigned int> * cell_and_index,
100 const unsigned int n_subdivisions,
128 for (
unsigned int dataset = 0;
144 if (curved_cell_region == curved_inner_cells ||
145 (curved_cell_region == curved_boundary &&
146 (cell_and_index->first->at_boundary() || (dim != spacedim))) ||
147 (cell_and_index->first->reference_cell() !=
148 ReferenceCells::get_hypercube<dim>()))
158 const std::vector<Point<spacedim>> &q_points =
161 patch.
data.reinit(scratch_data.
n_datasets + spacedim, n_q_points);
162 for (
unsigned int i = 0; i < spacedim; ++i)
163 for (
unsigned int q = 0; q < n_q_points; ++q)
164 patch.
data(patch.
data.size(0) - spacedim + i, q) = q_points[q][i];
177 unsigned int offset = 0;
180 unsigned int dataset_number = 0;
181 for (
const auto &dataset : this->dof_data)
185 const unsigned int n_components =
189 dataset->postprocessor;
191 if (postprocessor !=
nullptr)
198 if ((n_components == 1) &&
199 (dataset->is_complex_valued() ==
false))
208 dataset->get_function_values(
209 this_fe_patch_values,
215 dataset->get_function_gradients(
216 this_fe_patch_values,
222 dataset->get_function_hessians(
223 this_fe_patch_values,
235 dh_cell(&cell_and_index->first->get_triangulation(),
236 cell_and_index->first->level(),
237 cell_and_index->first->index(),
238 dataset->dof_handler);
259 if (dataset->is_complex_valued() ==
false)
264 dataset->get_function_values(
265 this_fe_patch_values,
271 dataset->get_function_gradients(
272 this_fe_patch_values,
278 dataset->get_function_hessians(
279 this_fe_patch_values,
290 if (n_components == 1)
299 dataset->get_function_values(
300 this_fe_patch_values,
302 ComponentExtractor::real_part,
306 for (
unsigned int i = 0;
308 .solution_values.size();
317 .solution_values[i][0] =
325 dataset->get_function_gradients(
326 this_fe_patch_values,
328 ComponentExtractor::real_part,
330 .solution_gradients);
332 for (
unsigned int i = 0;
334 .solution_gradients.size();
343 .solution_gradients[i][0] =
345 .solution_gradients[i];
351 dataset->get_function_hessians(
352 this_fe_patch_values,
354 ComponentExtractor::real_part,
358 for (
unsigned int i = 0;
360 .solution_hessians.size();
365 .solution_hessians[i]
369 .solution_hessians[i][0] =
371 .solution_hessians[i];
382 dataset->get_function_values(
383 this_fe_patch_values,
385 ComponentExtractor::imaginary_part,
389 for (
unsigned int i = 0;
391 .solution_values.size();
400 .solution_values[i][1] =
408 dataset->get_function_gradients(
409 this_fe_patch_values,
411 ComponentExtractor::imaginary_part,
413 .solution_gradients);
415 for (
unsigned int i = 0;
417 .solution_gradients.size();
426 .solution_gradients[i][1] =
428 .solution_gradients[i];
434 dataset->get_function_hessians(
435 this_fe_patch_values,
437 ComponentExtractor::imaginary_part,
441 for (
unsigned int i = 0;
443 .solution_hessians.size();
448 .solution_hessians[i]
452 .solution_hessians[i][1] =
454 .solution_hessians[i];
489 std::vector<Vector<double>> tmp(
495 dataset->get_function_values(
496 this_fe_patch_values,
498 ComponentExtractor::real_part,
503 for (
unsigned int i = 0;
505 .solution_values.size();
513 for (
unsigned int j = 0; j < n_components;
516 .solution_values[i][j] = tmp[i][j];
522 dataset->get_function_values(
523 this_fe_patch_values,
525 ComponentExtractor::imaginary_part,
528 for (
unsigned int i = 0;
530 .solution_values.size();
533 for (
unsigned int j = 0; j < n_components;
536 .solution_values[i][j + n_components] =
544 std::vector<std::vector<Tensor<1, spacedim>>> tmp(
546 .solution_gradients.size(),
550 dataset->get_function_gradients(
551 this_fe_patch_values,
553 ComponentExtractor::real_part,
556 for (
unsigned int i = 0;
558 .solution_gradients.size();
563 .solution_gradients[i]
566 for (
unsigned int j = 0; j < n_components;
569 .solution_gradients[i][j] = tmp[i][j];
573 dataset->get_function_gradients(
574 this_fe_patch_values,
576 ComponentExtractor::imaginary_part,
579 for (
unsigned int i = 0;
581 .solution_gradients.size();
584 for (
unsigned int j = 0; j < n_components;
587 .solution_gradients[i][j + n_components] =
595 std::vector<std::vector<Tensor<2, spacedim>>> tmp(
597 .solution_gradients.size(),
601 dataset->get_function_hessians(
602 this_fe_patch_values,
604 ComponentExtractor::real_part,
607 for (
unsigned int i = 0;
609 .solution_hessians.size();
614 .solution_hessians[i]
617 for (
unsigned int j = 0; j < n_components;
620 .solution_hessians[i][j] = tmp[i][j];
624 dataset->get_function_hessians(
625 this_fe_patch_values,
627 ComponentExtractor::imaginary_part,
630 for (
unsigned int i = 0;
632 .solution_hessians.size();
635 for (
unsigned int j = 0; j < n_components;
638 .solution_hessians[i][j + n_components] =
651 dh_cell(&cell_and_index->first->get_triangulation(),
652 cell_and_index->first->level(),
653 cell_and_index->first->index(),
654 dataset->dof_handler);
670 for (
unsigned int q = 0; q < n_q_points; ++q)
671 for (
unsigned int component = 0;
672 component < dataset->n_output_variables;
674 patch.
data(offset + component, q) =
680 offset += dataset->n_output_variables;
687 if (n_components == 1)
692 dataset->get_function_values(
693 this_fe_patch_values,
697 for (
unsigned int q = 0; q < n_q_points; ++q)
698 patch.
data(offset, q) =
707 if (dataset->is_complex_valued() ==
true)
709 dataset->get_function_values(
710 this_fe_patch_values,
714 for (
unsigned int q = 0; q < n_q_points; ++q)
715 patch.
data(offset, q) =
729 if (dataset->is_complex_valued() ==
false)
731 dataset->get_function_values(
732 this_fe_patch_values,
736 for (
unsigned int component = 0; component < n_components;
738 for (
unsigned int q = 0; q < n_q_points; ++q)
739 patch.
data(offset + component, q) =
744 offset += dataset->n_output_variables;
776 dataset->get_function_values(
777 this_fe_patch_values,
786 Assert(dataset->data_component_interpretation.size() ==
790 unsigned int destination = offset;
791 for (
unsigned int component = 0;
792 component < n_components;
796 dataset->data_component_interpretation[component])
809 for (
unsigned int q = 0; q < n_q_points;
811 patch.
data(destination, q) =
813 .solution_values[q](component);
835 const unsigned int size = spacedim;
836 for (
unsigned int c = 0; c < size; ++c)
837 for (
unsigned int q = 0; q < n_q_points;
839 patch.
data(destination + c, q) =
841 .solution_values[q](component + c);
844 destination += 2 * size;
853 const unsigned int size =
855 for (
unsigned int c = 0; c < size; ++c)
856 for (
unsigned int q = 0; q < n_q_points;
858 patch.
data(destination + c, q) =
860 .solution_values[q](component + c);
863 destination += 2 * size;
877 dataset->get_function_values(
878 this_fe_patch_values,
883 unsigned int destination = offset;
884 for (
unsigned int component = 0;
885 component < n_components;
889 dataset->data_component_interpretation[component])
898 for (
unsigned int q = 0; q < n_q_points;
900 patch.
data(destination + 1, q) =
902 .solution_values[q](component);
918 const unsigned int size = spacedim;
919 for (
unsigned int c = 0; c < size; ++c)
920 for (
unsigned int q = 0; q < n_q_points;
922 patch.
data(destination + size + c, q) =
924 .solution_values[q](component + c);
927 destination += 2 * size;
936 const unsigned int size =
938 for (
unsigned int c = 0; c < size; ++c)
939 for (
unsigned int q = 0; q < n_q_points;
941 patch.
data(destination + size + c, q) =
943 .solution_values[q](component + c);
946 destination += 2 * size;
961 offset += dataset->n_output_variables * 2;
973 if (this->cell_data.size() != 0)
977 for (
const auto &dataset : this->cell_data)
982 dataset->get_cell_data_value(cell_and_index->second,
984 ComponentExtractor::real_part);
985 for (
unsigned int q = 0; q < n_q_points; ++q)
990 if (dataset->is_complex_valued() ==
true)
992 const double value = dataset->get_cell_data_value(
993 cell_and_index->second,
996 for (
unsigned int q = 0; q < n_q_points; ++q)
1000 offset += (dataset->is_complex_valued() ? 2 : 1);
1006 for (
const unsigned int f : cell_and_index->first->face_indices())
1018 if (cell_and_index->first->at_boundary(f) ||
1019 (cell_and_index->first->neighbor(f)->level() !=
1020 cell_and_index->first->level()))
1026 const cell_iterator neighbor = cell_and_index->first->neighbor(f);
1027 Assert(
static_cast<unsigned int>(neighbor->level()) <
1030 if ((
static_cast<unsigned int>(neighbor->index()) >=
1033 [neighbor->index()] ==
1047 const unsigned int patch_idx =
1049 [cell_and_index->first->index()];
1057 this->patches[patch_idx].swap(patch);
1062template <
int dim,
int spacedim>
1069 .template get_default_linear_mapping<dim, spacedim>(),
1076template <
int dim,
int spacedim>
1079 const unsigned int n_subdivisions_,
1084 build_patches(mapping_collection, n_subdivisions_, curved_region);
1089template <
int dim,
int spacedim>
1093 const unsigned int n_subdivisions_,
1102 const unsigned int n_subdivisions =
1103 (n_subdivisions_ != 0) ? n_subdivisions_ : this->default_subdivisions;
1104 Assert(n_subdivisions >= 1,
1108 this->validate_dataset_names();
1126 std::vector<std::vector<unsigned int>> cell_to_patch_index_map;
1127 cell_to_patch_index_map.resize(this->
triangulation->n_levels());
1128 for (
unsigned int l = 0; l < this->
triangulation->n_levels(); ++l)
1131 unsigned int max_index = 0;
1135 if (
static_cast<unsigned int>(cell->level()) == l)
1137 std::max(max_index,
static_cast<unsigned int>(cell->index()));
1139 cell_to_patch_index_map[l].resize(
1144 std::vector<std::pair<cell_iterator, unsigned int>> all_cells;
1152 unsigned int active_index = 0;
1159 while (active_cell != this->
triangulation->end() && cell->is_active() &&
1160 decltype(active_cell)(cell) != active_cell)
1166 Assert(
static_cast<unsigned int>(cell->level()) <
1167 cell_to_patch_index_map.size(),
1169 Assert(
static_cast<unsigned int>(cell->index()) <
1170 cell_to_patch_index_map[cell->level()].size(),
1174 cell_to_patch_index_map[cell->level()][cell->index()] =
1177 all_cells.emplace_back(cell, active_index);
1181 this->patches.clear();
1182 this->patches.resize(all_cells.size());
1190 unsigned int n_datasets = 0;
1191 for (
unsigned int i = 0; i < this->cell_data.size(); ++i)
1192 n_datasets += (this->cell_data[i]->is_complex_valued() &&
1193 (this->cell_data[i]->postprocessor ==
nullptr) ?
1196 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
1197 n_datasets += (this->dof_data[i]->n_output_variables *
1198 (this->dof_data[i]->is_complex_valued() &&
1199 (this->dof_data[i]->postprocessor ==
nullptr) ?
1203 std::vector<unsigned int> n_postprocessor_outputs(this->dof_data.size());
1204 for (
unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
1205 if (this->dof_data[dataset]->postprocessor)
1206 n_postprocessor_outputs[dataset] =
1207 this->dof_data[dataset]->n_output_variables;
1209 n_postprocessor_outputs[dataset] = 0;
1212 (n_subdivisions < 2 ? no_curved_cells : curved_region);
1215 if (curved_cell_region != no_curved_cells)
1218 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
1219 if (this->dof_data[i]->postprocessor)
1221 this->dof_data[i]->postprocessor->get_needed_update_flags();
1227 "The update of normal vectors may not be requested for evaluation of "
1228 "data on cells via DataPostprocessor."));
1233 n_postprocessor_outputs,
1237 cell_to_patch_index_map);
1239 auto worker = [
this, n_subdivisions, curved_cell_region](
1240 const std::pair<cell_iterator, unsigned int> *cell_and_index,
1246 this->build_one_patch(cell_and_index,
1249 curved_cell_region);
1253 if (all_cells.size() > 0)
1255 all_cells.data() + all_cells.size(),
1258 std::function<
void(
const int)>(),
1273template <
int dim,
int spacedim>
1281 first_cell_function = first_cell;
1282 next_cell_function = next_cell;
1287template <
int dim,
int spacedim>
1292 const auto first_cell =
1303 const auto next_cell =
1323 set_cell_selection(first_cell, next_cell);
1328template <
int dim,
int spacedim>
1329const std::pair<typename DataOut<dim, spacedim>::FirstCellFunctionType,
1333 return std::make_pair(first_cell_function, next_cell_function);
1339#include "data_out.inst"
virtual void build_patches(const unsigned int n_subdivisions=0)
typename std::function< cell_iterator(const Triangulation< dim, spacedim > &, const cell_iterator &)> NextCellFunctionType
void build_one_patch(const std::pair< cell_iterator, unsigned int > *cell_and_index, internal::DataOutImplementation::ParallelData< dim, spacedim > &scratch_data, const unsigned int n_subdivisions, const CurvedCellRegion curved_cell_region)
const std::pair< FirstCellFunctionType, NextCellFunctionType > get_cell_selection() const
typename DataOut_DoFData< dim, dim, spacedim, spacedim >::cell_iterator cell_iterator
void set_cell_selection(const std::function< cell_iterator(const Triangulation< dim, spacedim > &)> &first_cell, const std::function< cell_iterator(const Triangulation< dim, spacedim > &, const cell_iterator &)> &next_cell)
virtual UpdateFlags get_needed_update_flags() const =0
virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
virtual void evaluate_scalar_field(const DataPostprocessorInputs::Scalar< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const Mapping< dim, spacedim > & get_mapping() const
const unsigned int n_quadrature_points
const FiniteElement< dim, spacedim > & get_fe() const
FilteredIterator & set_to_next_positive(const BaseIterator &bi)
unsigned int n_components() const
Abstract base class for mapping classes.
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
static unsigned int n_threads()
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
ReferenceCell reference_cell
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Point< spacedim > vertices[GeometryInfo< dim >::vertices_per_cell]
#define AssertDimension(dim1, dim2)
static const unsigned int space_dim
static ::ExceptionBase & ExcInternalError()
unsigned int n_subdivisions
static ::ExceptionBase & ExcNoTriangulationSelected()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
bool points_are_available
std::array< unsigned int, GeometryInfo< dim >::faces_per_cell > neighbors
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::cell_iterator cell_iterator
@ component_is_part_of_tensor
@ component_is_part_of_vector
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 const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
const FEValuesBase< dim, spacedim > & get_present_fe_values(const unsigned int dataset) const
const unsigned int n_datasets
DataPostprocessorInputs::Scalar< spacedim > patch_values_scalar
std::vector< std::vector<::Vector< double > > > postprocessed_values
DataPostprocessorInputs::Vector< spacedim > patch_values_system
void reinit_all_fe_values(std::vector< std::shared_ptr< DataEntryBase< dim, spacedim > > > &dof_data, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face=numbers::invalid_unsigned_int)
void resize_system_vectors(const unsigned int n_components)
const std::vector< std::vector< unsigned int > > * cell_to_patch_index_map
ParallelData(const unsigned int n_datasets, const unsigned int n_subdivisions, const std::vector< unsigned int > &n_postprocessor_outputs, const ::hp::MappingCollection< dim, spacedim > &mapping, const std::vector< std::shared_ptr<::hp::FECollection< dim, spacedim > > > &finite_elements, const UpdateFlags update_flags, const std::vector< std::vector< unsigned int > > &cell_to_patch_index_map)