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,
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)
68 template <
int dim,
typename DoFHandlerType>
77 return this->first_locally_owned_cell();
80 return this->next_locally_owned_cell(cell);
86 template <
int dim,
typename DoFHandlerType>
89 const std::pair<cell_iterator, unsigned int> *cell_and_index,
91 DoFHandlerType::space_dimension>
93 const unsigned int n_subdivisions,
98 DoFHandlerType::space_dimension>
105 for (
const unsigned int vertex :
107 if (scratch_data.mapping_collection[0].preserves_vertex_locations())
108 patch.vertices[vertex] = cell_and_index->first->vertex(vertex);
110 patch.vertices[vertex] =
111 scratch_data.mapping_collection[0].transform_unit_to_real_cell(
112 cell_and_index->first,
116 scratch_data.reinit_all_fe_values(this->dof_data, cell_and_index->first);
119 &fe_patch_values = scratch_data.get_present_fe_values(0);
133 if (curved_cell_region == curved_inner_cells ||
134 (curved_cell_region == curved_boundary &&
135 (cell_and_index->first->at_boundary() ||
136 (DoFHandlerType::dimension != DoFHandlerType::space_dimension))))
138 Assert(patch.space_dim == DoFHandlerType::space_dimension,
143 patch.points_are_available =
true;
147 const std::vector<Point<DoFHandlerType::space_dimension>> &q_points =
150 patch.data.reinit(scratch_data.n_datasets +
151 DoFHandlerType::space_dimension,
153 for (
unsigned int i = 0; i < DoFHandlerType::space_dimension; ++i)
154 for (
unsigned int q = 0; q < n_q_points; ++q)
155 patch.data(patch.data.size(0) - DoFHandlerType::space_dimension + i,
160 patch.data.reinit(scratch_data.n_datasets, n_q_points);
161 patch.points_are_available =
false;
166 if (scratch_data.n_datasets > 0)
169 unsigned int offset = 0;
172 unsigned int dataset_number = 0;
173 for (
const auto &dataset : this->dof_data)
176 DoFHandlerType::space_dimension>
177 &this_fe_patch_values =
178 scratch_data.get_present_fe_values(dataset_number);
180 this_fe_patch_values.get_fe().n_components();
183 *postprocessor = dataset->postprocessor;
185 if (postprocessor !=
nullptr)
193 (dataset->is_complex_valued() ==
false))
202 dataset->get_function_values(
203 this_fe_patch_values,
206 scratch_data.patch_values_scalar.solution_values);
209 dataset->get_function_gradients(
210 this_fe_patch_values,
213 scratch_data.patch_values_scalar.solution_gradients);
216 dataset->get_function_hessians(
217 this_fe_patch_values,
220 scratch_data.patch_values_scalar.solution_hessians);
225 scratch_data.patch_values_scalar.evaluation_points =
226 this_fe_patch_values.get_quadrature_points();
228 const typename DoFHandlerType::cell_iterator dh_cell(
229 &cell_and_index->first->get_triangulation(),
230 cell_and_index->first->level(),
231 cell_and_index->first->index(),
232 dataset->dof_handler);
233 scratch_data.patch_values_scalar
234 .template set_cell<DoFHandlerType>(dh_cell);
239 scratch_data.patch_values_scalar,
240 scratch_data.postprocessed_values[dataset_number]);
253 if (dataset->is_complex_valued() ==
false)
258 dataset->get_function_values(
259 this_fe_patch_values,
262 scratch_data.patch_values_system.solution_values);
265 dataset->get_function_gradients(
266 this_fe_patch_values,
269 scratch_data.patch_values_system.solution_gradients);
272 dataset->get_function_hessians(
273 this_fe_patch_values,
276 scratch_data.patch_values_system.solution_hessians);
286 scratch_data.resize_system_vectors(2);
293 dataset->get_function_values(
294 this_fe_patch_values,
296 ComponentExtractor::real_part,
297 scratch_data.patch_values_scalar
300 for (
unsigned int i = 0;
301 i < scratch_data.patch_values_scalar
302 .solution_values.size();
306 scratch_data.patch_values_system
310 scratch_data.patch_values_system
311 .solution_values[i][0] =
312 scratch_data.patch_values_scalar
319 dataset->get_function_gradients(
320 this_fe_patch_values,
322 ComponentExtractor::real_part,
323 scratch_data.patch_values_scalar
324 .solution_gradients);
326 for (
unsigned int i = 0;
327 i < scratch_data.patch_values_scalar
328 .solution_gradients.size();
332 scratch_data.patch_values_system
336 scratch_data.patch_values_system
337 .solution_gradients[i][0] =
338 scratch_data.patch_values_scalar
339 .solution_gradients[i];
345 dataset->get_function_hessians(
346 this_fe_patch_values,
348 ComponentExtractor::real_part,
349 scratch_data.patch_values_scalar
352 for (
unsigned int i = 0;
353 i < scratch_data.patch_values_scalar
354 .solution_hessians.size();
358 scratch_data.patch_values_system
359 .solution_hessians[i]
362 scratch_data.patch_values_system
363 .solution_hessians[i][0] =
364 scratch_data.patch_values_scalar
365 .solution_hessians[i];
376 dataset->get_function_values(
377 this_fe_patch_values,
379 ComponentExtractor::imaginary_part,
380 scratch_data.patch_values_scalar
383 for (
unsigned int i = 0;
384 i < scratch_data.patch_values_scalar
385 .solution_values.size();
389 scratch_data.patch_values_system
393 scratch_data.patch_values_system
394 .solution_values[i][1] =
395 scratch_data.patch_values_scalar
402 dataset->get_function_gradients(
403 this_fe_patch_values,
405 ComponentExtractor::imaginary_part,
406 scratch_data.patch_values_scalar
407 .solution_gradients);
409 for (
unsigned int i = 0;
410 i < scratch_data.patch_values_scalar
411 .solution_gradients.size();
415 scratch_data.patch_values_system
419 scratch_data.patch_values_system
420 .solution_gradients[i][1] =
421 scratch_data.patch_values_scalar
422 .solution_gradients[i];
428 dataset->get_function_hessians(
429 this_fe_patch_values,
431 ComponentExtractor::imaginary_part,
432 scratch_data.patch_values_scalar
435 for (
unsigned int i = 0;
436 i < scratch_data.patch_values_scalar
437 .solution_hessians.size();
441 scratch_data.patch_values_system
442 .solution_hessians[i]
445 scratch_data.patch_values_system
446 .solution_hessians[i][1] =
447 scratch_data.patch_values_scalar
448 .solution_hessians[i];
483 std::vector<Vector<double>> tmp(
484 scratch_data.patch_values_system.solution_values
489 dataset->get_function_values(
490 this_fe_patch_values,
492 ComponentExtractor::real_part,
497 for (
unsigned int i = 0;
498 i < scratch_data.patch_values_system
499 .solution_values.size();
503 scratch_data.patch_values_system
509 scratch_data.patch_values_system
510 .solution_values[i][j] = tmp[i][j];
516 dataset->get_function_values(
517 this_fe_patch_values,
519 ComponentExtractor::imaginary_part,
522 for (
unsigned int i = 0;
523 i < scratch_data.patch_values_system
524 .solution_values.size();
529 scratch_data.patch_values_system
538 std::vector<std::vector<
541 scratch_data.patch_values_system
542 .solution_gradients.size(),
548 dataset->get_function_gradients(
549 this_fe_patch_values,
551 ComponentExtractor::real_part,
554 for (
unsigned int i = 0;
555 i < scratch_data.patch_values_system
556 .solution_gradients.size();
560 scratch_data.patch_values_system
561 .solution_gradients[i]
566 scratch_data.patch_values_system
567 .solution_gradients[i][j] = tmp[i][j];
571 dataset->get_function_gradients(
572 this_fe_patch_values,
574 ComponentExtractor::imaginary_part,
577 for (
unsigned int i = 0;
578 i < scratch_data.patch_values_system
579 .solution_gradients.size();
584 scratch_data.patch_values_system
593 std::vector<std::vector<
596 scratch_data.patch_values_system
597 .solution_gradients.size(),
603 dataset->get_function_hessians(
604 this_fe_patch_values,
606 ComponentExtractor::real_part,
609 for (
unsigned int i = 0;
610 i < scratch_data.patch_values_system
611 .solution_hessians.size();
615 scratch_data.patch_values_system
616 .solution_hessians[i]
621 scratch_data.patch_values_system
622 .solution_hessians[i][j] = tmp[i][j];
626 dataset->get_function_hessians(
627 this_fe_patch_values,
629 ComponentExtractor::imaginary_part,
632 for (
unsigned int i = 0;
633 i < scratch_data.patch_values_system
634 .solution_hessians.size();
639 scratch_data.patch_values_system
649 scratch_data.patch_values_system.evaluation_points =
650 this_fe_patch_values.get_quadrature_points();
652 const typename DoFHandlerType::cell_iterator dh_cell(
653 &cell_and_index->first->get_triangulation(),
654 cell_and_index->first->level(),
655 cell_and_index->first->index(),
656 dataset->dof_handler);
657 scratch_data.patch_values_system
658 .template set_cell<DoFHandlerType>(dh_cell);
665 scratch_data.patch_values_system,
666 scratch_data.postprocessed_values[dataset_number]);
672 for (
unsigned int q = 0; q < n_q_points; ++q)
673 for (
unsigned int component = 0;
674 component < dataset->n_output_variables;
676 patch.data(offset + component, q) =
677 scratch_data.postprocessed_values[dataset_number][q](
682 offset += dataset->n_output_variables;
694 dataset->get_function_values(
695 this_fe_patch_values,
698 scratch_data.patch_values_scalar.solution_values);
699 for (
unsigned int q = 0; q < n_q_points; ++q)
700 patch.data(offset, q) =
701 scratch_data.patch_values_scalar.solution_values[q];
709 if (dataset->is_complex_valued() ==
true)
711 dataset->get_function_values(
712 this_fe_patch_values,
715 scratch_data.patch_values_scalar.solution_values);
716 for (
unsigned int q = 0; q < n_q_points; ++q)
717 patch.data(offset, q) =
718 scratch_data.patch_values_scalar.solution_values[q];
731 if (dataset->is_complex_valued() ==
false)
733 dataset->get_function_values(
734 this_fe_patch_values,
737 scratch_data.patch_values_system.solution_values);
738 for (
unsigned int component = 0; component <
n_components;
740 for (
unsigned int q = 0; q < n_q_points; ++q)
741 patch.data(offset + component, q) =
742 scratch_data.patch_values_system.solution_values[q](
746 offset += dataset->n_output_variables;
778 dataset->get_function_values(
779 this_fe_patch_values,
782 scratch_data.patch_values_system.solution_values);
788 Assert(dataset->data_component_interpretation.size() ==
792 unsigned int destination = offset;
793 for (
unsigned int component = 0;
798 dataset->data_component_interpretation[component])
811 for (
unsigned int q = 0; q < n_q_points;
813 patch.data(destination, q) =
814 scratch_data.patch_values_system
815 .solution_values[q](component);
837 const unsigned int size =
838 DoFHandlerType::space_dimension;
839 for (
unsigned int c = 0; c < size; ++c)
840 for (
unsigned int q = 0; q < n_q_points;
842 patch.data(destination + c, q) =
843 scratch_data.patch_values_system
844 .solution_values[q](component + c);
847 destination += 2 * size;
856 const unsigned int size =
857 DoFHandlerType::space_dimension *
858 DoFHandlerType::space_dimension;
859 for (
unsigned int c = 0; c < size; ++c)
860 for (
unsigned int q = 0; q < n_q_points;
862 patch.data(destination + c, q) =
863 scratch_data.patch_values_system
864 .solution_values[q](component + c);
867 destination += 2 * size;
881 dataset->get_function_values(
882 this_fe_patch_values,
885 scratch_data.patch_values_system.solution_values);
887 unsigned int destination = offset;
888 for (
unsigned int component = 0;
893 dataset->data_component_interpretation[component])
902 for (
unsigned int q = 0; q < n_q_points;
904 patch.data(destination + 1, q) =
905 scratch_data.patch_values_system
906 .solution_values[q](component);
922 const unsigned int size =
923 DoFHandlerType::space_dimension;
924 for (
unsigned int c = 0; c < size; ++c)
925 for (
unsigned int q = 0; q < n_q_points;
927 patch.data(destination + size + c, q) =
928 scratch_data.patch_values_system
929 .solution_values[q](component + c);
932 destination += 2 * size;
941 const unsigned int size =
942 DoFHandlerType::space_dimension *
943 DoFHandlerType::space_dimension;
944 for (
unsigned int c = 0; c < size; ++c)
945 for (
unsigned int q = 0; q < n_q_points;
947 patch.data(destination + size + c, q) =
948 scratch_data.patch_values_system
949 .solution_values[q](component + c);
952 destination += 2 * size;
967 offset += dataset->n_output_variables * 2;
979 if (this->cell_data.size() != 0)
983 for (
const auto &dataset : this->cell_data)
988 dataset->get_cell_data_value(cell_and_index->second,
990 ComponentExtractor::real_part);
991 for (
unsigned int q = 0; q < n_q_points; ++q)
992 patch.data(offset, q) =
value;
996 if (dataset->is_complex_valued() ==
true)
998 const double value = dataset->get_cell_data_value(
999 cell_and_index->second,
1002 for (
unsigned int q = 0; q < n_q_points; ++q)
1003 patch.data(offset + 1, q) =
value;
1006 offset += (dataset->is_complex_valued() ? 2 : 1);
1012 for (
const unsigned int f :
1025 if (cell_and_index->first->at_boundary(f) ||
1026 (cell_and_index->first->neighbor(f)->level() !=
1027 cell_and_index->first->level()))
1033 const cell_iterator neighbor = cell_and_index->first->neighbor(f);
1034 Assert(
static_cast<unsigned int>(neighbor->level()) <
1035 scratch_data.cell_to_patch_index_map->size(),
1037 if ((
static_cast<unsigned int>(neighbor->index()) >=
1038 (*scratch_data.cell_to_patch_index_map)[neighbor->level()].size()) ||
1039 ((*scratch_data.cell_to_patch_index_map)[neighbor->level()]
1040 [neighbor->index()] ==
1049 patch.neighbors[f] =
1051 .cell_to_patch_index_map)[neighbor->level()][neighbor->index()];
1054 const unsigned int patch_idx =
1055 (*scratch_data.cell_to_patch_index_map)[cell_and_index->first->level()]
1056 [cell_and_index->first->index()];
1059 patch.patch_index = patch_idx;
1064 this->patches[patch_idx].swap(patch);
1069 template <
int dim,
typename DoFHandlerType>
1074 DoFHandlerType::space_dimension>::mapping,
1081 template <
int dim,
typename DoFHandlerType>
1086 const unsigned int n_subdivisions_,
1090 Assert(dim == DoFHandlerType::dimension,
1096 const unsigned int n_subdivisions =
1097 (n_subdivisions_ != 0) ? n_subdivisions_ : this->default_subdivisions;
1098 Assert(n_subdivisions >= 1,
1102 this->validate_dataset_names();
1120 std::vector<std::vector<unsigned int>> cell_to_patch_index_map;
1121 cell_to_patch_index_map.resize(this->
triangulation->n_levels());
1125 unsigned int max_index = 0;
1129 if (
static_cast<unsigned int>(cell->level()) ==
l)
1131 std::max(max_index,
static_cast<unsigned int>(cell->index()));
1133 cell_to_patch_index_map[
l].resize(
1136 DoFHandlerType::dimension,
1137 DoFHandlerType::space_dimension>::no_neighbor);
1141 std::vector<std::pair<cell_iterator, unsigned int>> all_cells;
1149 unsigned int active_index = 0;
1156 while (active_cell != this->
triangulation->end() && cell->is_active() &&
1157 decltype(active_cell)(cell) != active_cell)
1163 Assert(
static_cast<unsigned int>(cell->level()) <
1164 cell_to_patch_index_map.size(),
1166 Assert(
static_cast<unsigned int>(cell->index()) <
1167 cell_to_patch_index_map[cell->level()].size(),
1171 cell_to_patch_index_map[cell->level()][cell->index()] =
1174 all_cells.emplace_back(cell, active_index);
1178 this->patches.clear();
1179 this->patches.resize(all_cells.size());
1187 unsigned int n_datasets = 0;
1188 for (
unsigned int i = 0; i < this->cell_data.size(); ++i)
1189 n_datasets += (this->cell_data[i]->is_complex_valued() &&
1190 (this->cell_data[i]->postprocessor ==
nullptr) ?
1193 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
1194 n_datasets += (this->dof_data[i]->n_output_variables *
1195 (this->dof_data[i]->is_complex_valued() &&
1196 (this->dof_data[i]->postprocessor ==
nullptr) ?
1200 std::vector<unsigned int> n_postprocessor_outputs(this->dof_data.size());
1201 for (
unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
1202 if (this->dof_data[dataset]->postprocessor)
1203 n_postprocessor_outputs[dataset] =
1204 this->dof_data[dataset]->n_output_variables;
1206 n_postprocessor_outputs[dataset] = 0;
1209 (n_subdivisions < 2 ? no_curved_cells : curved_region);
1212 if (curved_cell_region != no_curved_cells)
1215 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
1216 if (this->dof_data[i]->postprocessor)
1218 this->dof_data[i]->postprocessor->get_needed_update_flags();
1224 "The update of normal vectors may not be requested for evaluation of "
1225 "data on cells via DataPostprocessor."));
1228 DoFHandlerType::space_dimension>
1229 thread_data(n_datasets,
1231 n_postprocessor_outputs,
1235 cell_to_patch_index_map);
1237 auto worker = [
this, n_subdivisions, curved_cell_region](
1238 const std::pair<cell_iterator, unsigned int> *cell_and_index,
1240 DoFHandlerType::dimension,
1241 DoFHandlerType::space_dimension> &scratch_data,
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)>(),
1272 template <
int dim,
typename DoFHandlerType>
1280 first_cell_function = first_cell;
1281 next_cell_function = next_cell;
1286 template <
int dim,
typename DoFHandlerType>
1291 const auto first_cell =
1302 const auto next_cell =
1322 set_cell_selection(first_cell, next_cell);
1327 template <
int dim,
typename DoFHandlerType>
1328 const std::pair<typename DataOut<dim, DoFHandlerType>::FirstCellFunctionType,
1332 return std::make_pair(first_cell_function, next_cell_function);
1337 template <
int dim,
typename DoFHandlerType>
1346 template <
int dim,
typename DoFHandlerType>
1354 DoFHandlerType::space_dimension>::active_cell_iterator
1362 template <
int dim,
typename DoFHandlerType>
1370 while ((cell != this->
triangulation->end()) && cell->is_active() &&
1371 !cell->is_locally_owned())
1372 cell = next_cell(cell);
1379 template <
int dim,
typename DoFHandlerType>
1385 next_cell(old_cell);
1386 while ((cell != this->
triangulation->end()) && cell->is_active() &&
1387 !cell->is_locally_owned())
1388 cell = next_cell(cell);
1394 #include "data_out.inst"