40 namespace DataOutImplementation
42 template <
int dim,
int spacedim>
44 const unsigned int n_datasets,
45 const unsigned int n_subdivisions,
46 const std::vector<unsigned int> &n_postprocessor_outputs,
47 const ::hp::MappingCollection<dim, spacedim> &mapping,
52 const std::vector<std::vector<unsigned int>> &cell_to_patch_index_map)
55 n_postprocessor_outputs,
60 , cell_to_patch_index_map(&cell_to_patch_index_map)
67template <
int dim,
int spacedim>
77 while ((cell != this->
triangulation->end()) && !cell->is_locally_owned())
86 while ((cell != this->
triangulation->end()) && !cell->is_locally_owned())
94template <
int dim,
int spacedim>
97 const std::pair<cell_iterator, unsigned int> * cell_and_index,
99 const unsigned int n_subdivisions,
127 for (
unsigned int dataset = 0;
143 if (curved_cell_region == curved_inner_cells ||
144 (curved_cell_region == curved_boundary &&
145 (cell_and_index->first->at_boundary() || (dim != spacedim))) ||
146 (cell_and_index->first->reference_cell() !=
147 ReferenceCells::get_hypercube<dim>()))
157 const std::vector<Point<spacedim>> &q_points =
160 patch.
data.reinit(scratch_data.
n_datasets + spacedim, n_q_points);
161 for (
unsigned int i = 0; i < spacedim; ++i)
162 for (
unsigned int q = 0; q < n_q_points; ++q)
163 patch.
data(patch.
data.size(0) - spacedim + i, q) = q_points[q][i];
176 unsigned int offset = 0;
179 unsigned int dataset_number = 0;
180 for (
const auto &dataset : this->dof_data)
184 const unsigned int n_components =
188 dataset->postprocessor;
190 if (postprocessor !=
nullptr)
197 if ((n_components == 1) &&
198 (dataset->is_complex_valued() ==
false))
207 dataset->get_function_values(
208 this_fe_patch_values,
214 dataset->get_function_gradients(
215 this_fe_patch_values,
221 dataset->get_function_hessians(
222 this_fe_patch_values,
225 scratch_data.patch_values_scalar.solution_hessians);
234 dh_cell(&cell_and_index->first->get_triangulation(),
235 cell_and_index->first->level(),
236 cell_and_index->first->index(),
237 dataset->dof_handler);
258 if (dataset->is_complex_valued() ==
false)
263 dataset->get_function_values(
264 this_fe_patch_values,
267 scratch_data.patch_values_system.solution_values);
270 dataset->get_function_gradients(
271 this_fe_patch_values,
274 scratch_data.patch_values_system.solution_gradients);
277 dataset->get_function_hessians(
278 this_fe_patch_values,
281 scratch_data.patch_values_system.solution_hessians);
289 if (n_components == 1)
298 dataset->get_function_values(
299 this_fe_patch_values,
301 ComponentExtractor::real_part,
305 for (
unsigned int i = 0;
307 .solution_values.size();
316 .solution_values[i][0] =
324 dataset->get_function_gradients(
325 this_fe_patch_values,
327 ComponentExtractor::real_part,
329 .solution_gradients);
331 for (
unsigned int i = 0;
333 .solution_gradients.size();
342 .solution_gradients[i][0] =
344 .solution_gradients[i];
350 dataset->get_function_hessians(
351 this_fe_patch_values,
354 scratch_data.patch_values_scalar
357 for (
unsigned int i = 0;
359 .solution_hessians.size();
364 .solution_hessians[i]
368 .solution_hessians[i][0] =
370 .solution_hessians[i];
381 dataset->get_function_values(
382 this_fe_patch_values,
384 ComponentExtractor::imaginary_part,
388 for (
unsigned int i = 0;
390 .solution_values.size();
399 .solution_values[i][1] =
407 dataset->get_function_gradients(
408 this_fe_patch_values,
410 ComponentExtractor::imaginary_part,
411 scratch_data.patch_values_scalar
412 .solution_gradients);
414 for (
unsigned int i = 0;
416 .solution_gradients.size();
425 .solution_gradients[i][1] =
427 .solution_gradients[i];
433 dataset->get_function_hessians(
434 this_fe_patch_values,
437 scratch_data.patch_values_scalar
440 for (
unsigned int i = 0;
442 .solution_hessians.size();
447 .solution_hessians[i]
451 .solution_hessians[i][1] =
453 .solution_hessians[i];
488 std::vector<Vector<double>> tmp(
494 dataset->get_function_values(
495 this_fe_patch_values,
502 for (
unsigned int i = 0;
504 .solution_values.size();
512 for (
unsigned int j = 0; j < n_components;
515 .solution_values[i][j] = tmp[i][j];
521 dataset->get_function_values(
522 this_fe_patch_values,
527 for (
unsigned int i = 0;
529 .solution_values.size();
532 for (
unsigned int j = 0; j < n_components;
535 .solution_values[i][j + n_components] =
543 std::vector<std::vector<Tensor<1, spacedim>>> tmp(
545 .solution_gradients.size(),
549 dataset->get_function_gradients(
550 this_fe_patch_values,
555 for (
unsigned int i = 0;
557 .solution_gradients.size();
562 .solution_gradients[i]
565 for (
unsigned int j = 0; j < n_components;
568 .solution_gradients[i][j] = tmp[i][j];
572 dataset->get_function_gradients(
573 this_fe_patch_values,
578 for (
unsigned int i = 0;
580 .solution_gradients.size();
583 for (
unsigned int j = 0; j < n_components;
586 .solution_gradients[i][j + n_components] =
594 std::vector<std::vector<Tensor<2, spacedim>>> tmp(
596 .solution_gradients.size(),
600 dataset->get_function_hessians(
601 this_fe_patch_values,
606 for (
unsigned int i = 0;
608 .solution_hessians.size();
613 .solution_hessians[i]
616 for (
unsigned int j = 0; j < n_components;
619 .solution_hessians[i][j] = tmp[i][j];
623 dataset->get_function_hessians(
624 this_fe_patch_values,
629 for (
unsigned int i = 0;
631 .solution_hessians.size();
634 for (
unsigned int j = 0; j < n_components;
637 .solution_hessians[i][j + n_components] =
650 dh_cell(&cell_and_index->first->get_triangulation(),
651 cell_and_index->first->level(),
652 cell_and_index->first->index(),
653 dataset->dof_handler);
669 for (
unsigned int q = 0; q < n_q_points; ++q)
670 for (
unsigned int component = 0;
671 component < dataset->n_output_variables;
673 patch.
data(offset + component, q) =
679 offset += dataset->n_output_variables;
686 if (n_components == 1)
691 dataset->get_function_values(
692 this_fe_patch_values,
695 scratch_data.patch_values_scalar.solution_values);
696 for (
unsigned int q = 0; q < n_q_points; ++q)
697 patch.
data(offset, q) =
706 if (dataset->is_complex_valued() ==
true)
708 dataset->get_function_values(
709 this_fe_patch_values,
712 scratch_data.patch_values_scalar.solution_values);
713 for (
unsigned int q = 0; q < n_q_points; ++q)
714 patch.
data(offset, q) =
728 if (dataset->is_complex_valued() ==
false)
730 dataset->get_function_values(
731 this_fe_patch_values,
734 scratch_data.patch_values_system.solution_values);
735 for (
unsigned int component = 0; component < n_components;
737 for (
unsigned int q = 0; q < n_q_points; ++q)
738 patch.
data(offset + component, q) =
743 offset += dataset->n_output_variables;
775 dataset->get_function_values(
776 this_fe_patch_values,
779 scratch_data.patch_values_system.solution_values);
785 Assert(dataset->data_component_interpretation.size() ==
789 unsigned int destination = offset;
790 for (
unsigned int component = 0;
791 component < n_components;
795 dataset->data_component_interpretation[component])
808 for (
unsigned int q = 0; q < n_q_points;
810 patch.
data(destination, q) =
812 .solution_values[q](component);
834 const unsigned int size = spacedim;
835 for (
unsigned int c = 0; c < size; ++c)
836 for (
unsigned int q = 0; q < n_q_points;
838 patch.
data(destination + c, q) =
840 .solution_values[q](component + c);
843 destination += 2 * size;
852 const unsigned int size =
854 for (
unsigned int c = 0; c < size; ++c)
855 for (
unsigned int q = 0; q < n_q_points;
857 patch.
data(destination + c, q) =
859 .solution_values[q](component + c);
862 destination += 2 * size;
876 dataset->get_function_values(
877 this_fe_patch_values,
880 scratch_data.patch_values_system.solution_values);
882 unsigned int destination = offset;
883 for (
unsigned int component = 0;
884 component < n_components;
888 dataset->data_component_interpretation[component])
897 for (
unsigned int q = 0; q < n_q_points;
899 patch.
data(destination + 1, q) =
901 .solution_values[q](component);
917 const unsigned int size = spacedim;
918 for (
unsigned int c = 0; c < size; ++c)
919 for (
unsigned int q = 0; q < n_q_points;
921 patch.
data(destination + size + c, q) =
923 .solution_values[q](component + c);
926 destination += 2 * size;
935 const unsigned int size =
937 for (
unsigned int c = 0; c < size; ++c)
938 for (
unsigned int q = 0; q < n_q_points;
940 patch.
data(destination + size + c, q) =
942 .solution_values[q](component + c);
945 destination += 2 * size;
960 offset += dataset->n_output_variables * 2;
972 if (this->cell_data.size() != 0)
976 for (
const auto &dataset : this->cell_data)
981 dataset->get_cell_data_value(cell_and_index->second,
984 for (
unsigned int q = 0; q < n_q_points; ++q)
989 if (dataset->is_complex_valued() ==
true)
991 const double value = dataset->get_cell_data_value(
992 cell_and_index->second,
995 for (
unsigned int q = 0; q < n_q_points; ++q)
999 offset += (dataset->is_complex_valued() ? 2 : 1);
1005 for (
const unsigned int f : cell_and_index->
first->face_indices())
1017 if (cell_and_index->first->at_boundary(f) ||
1018 (cell_and_index->first->neighbor(f)->level() !=
1019 cell_and_index->first->level()))
1025 const cell_iterator neighbor = cell_and_index->first->neighbor(f);
1026 Assert(
static_cast<unsigned int>(neighbor->level()) <
1029 if ((
static_cast<unsigned int>(neighbor->index()) >=
1032 [neighbor->index()] ==
1046 const unsigned int patch_idx =
1048 [cell_and_index->first->index()];
1056 this->patches[patch_idx].swap(patch);
1061template <
int dim,
int spacedim>
1068 .template get_default_linear_mapping<dim, spacedim>(),
1075template <
int dim,
int spacedim>
1078 const unsigned int n_subdivisions_,
1083 build_patches(mapping_collection, n_subdivisions_, curved_region);
1088template <
int dim,
int spacedim>
1092 const unsigned int n_subdivisions_,
1101 const unsigned int n_subdivisions =
1102 (n_subdivisions_ != 0) ? n_subdivisions_ : this->default_subdivisions;
1103 Assert(n_subdivisions >= 1,
1107 this->validate_dataset_names();
1125 std::vector<std::vector<unsigned int>> cell_to_patch_index_map;
1126 cell_to_patch_index_map.resize(this->
triangulation->n_levels());
1127 for (
unsigned int l = 0; l < this->
triangulation->n_levels(); ++l)
1130 unsigned int max_index = 0;
1134 if (
static_cast<unsigned int>(cell->level()) == l)
1136 std::max(max_index,
static_cast<unsigned int>(cell->index()));
1138 cell_to_patch_index_map[l].resize(
1143 std::vector<std::pair<cell_iterator, unsigned int>> all_cells;
1151 unsigned int active_index = 0;
1158 while (active_cell != this->
triangulation->end() && cell->is_active() &&
1159 decltype(active_cell)(cell) != active_cell)
1165 Assert(
static_cast<unsigned int>(cell->level()) <
1166 cell_to_patch_index_map.size(),
1168 Assert(
static_cast<unsigned int>(cell->index()) <
1169 cell_to_patch_index_map[cell->level()].size(),
1173 cell_to_patch_index_map[cell->level()][cell->index()] =
1176 all_cells.emplace_back(cell, active_index);
1180 this->patches.clear();
1181 this->patches.resize(all_cells.size());
1189 unsigned int n_datasets = 0;
1190 for (
unsigned int i = 0; i < this->cell_data.size(); ++i)
1191 n_datasets += (this->cell_data[i]->is_complex_valued() &&
1192 (this->cell_data[i]->postprocessor ==
nullptr) ?
1195 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
1196 n_datasets += (this->dof_data[i]->n_output_variables *
1197 (this->dof_data[i]->is_complex_valued() &&
1198 (this->dof_data[i]->postprocessor ==
nullptr) ?
1202 std::vector<unsigned int> n_postprocessor_outputs(this->dof_data.size());
1203 for (
unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
1204 if (this->dof_data[dataset]->postprocessor)
1205 n_postprocessor_outputs[dataset] =
1206 this->dof_data[dataset]->n_output_variables;
1208 n_postprocessor_outputs[dataset] = 0;
1211 (n_subdivisions < 2 ? no_curved_cells : curved_region);
1214 if (curved_cell_region != no_curved_cells)
1217 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
1218 if (this->dof_data[i]->postprocessor)
1220 this->dof_data[i]->postprocessor->get_needed_update_flags();
1226 "The update of normal vectors may not be requested for evaluation of "
1227 "data on cells via DataPostprocessor."));
1232 n_postprocessor_outputs,
1236 cell_to_patch_index_map);
1238 auto worker = [
this, n_subdivisions, curved_cell_region](
1239 const std::pair<cell_iterator, unsigned int> *cell_and_index,
1245 this->build_one_patch(cell_and_index,
1248 curved_cell_region);
1252 if (all_cells.size() > 0)
1254 all_cells.data() + all_cells.size(),
1257 std::function<
void(
const int)>(),
1272template <
int dim,
int spacedim>
1280 first_cell_function = first_cell;
1281 next_cell_function = next_cell;
1286template <
int dim,
int spacedim>
1291 const auto first_cell =
1302 const auto next_cell =
1322 set_cell_selection(first_cell, next_cell);
1327template <
int dim,
int spacedim>
1328const std::pair<typename DataOut<dim, spacedim>::FirstCellFunctionType,
1332 return std::make_pair(first_cell_function, next_cell_function);
1338#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
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNoTriangulationSelected()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::cell_iterator cell_iterator
@ 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.
@ 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
ReferenceCell reference_cell
static const unsigned int space_dim
unsigned int n_subdivisions
std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > vertices
bool points_are_available
std::array< unsigned int, GeometryInfo< dim >::faces_per_cell > neighbors
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)