deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
point_value_history.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2009 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
16
24#include <deal.II/lac/vector.h>
26
29
30#include <algorithm>
31
32
34
35
36namespace internal
37{
39 {
40 template <int dim>
42 const Point<dim> &new_requested_location,
43 const std::vector<Point<dim>> &new_locations,
44 const std::vector<types::global_dof_index> &new_sol_indices)
45 {
46 requested_location = new_requested_location;
47 support_point_locations = new_locations;
48 solution_indices = new_sol_indices;
49 }
50 } // namespace PointValueHistoryImplementation
51} // namespace internal
52
53
54
55template <int dim>
57 const unsigned int n_independent_variables)
58 : n_indep(n_independent_variables)
59{
60 closed = false;
61 cleared = false;
63 have_dof_handler = false;
64
65 // make a vector for keys
66 dataset_key = std::vector<double>(); // initialize the std::vector
67
68 // make a vector of independent values
70 std::vector<std::vector<double>>(n_indep, std::vector<double>(0));
71 indep_names = std::vector<std::string>();
72}
73
74
75
76template <int dim>
79 const unsigned int n_independent_variables)
81 , n_indep(n_independent_variables)
82{
83 closed = false;
84 cleared = false;
86 have_dof_handler = true;
87
88 // make a vector to store keys
89 dataset_key = std::vector<double>(); // initialize the std::vector
90
91 // make a vector for the independent values
93 std::vector<std::vector<double>>(n_indep, std::vector<double>(0));
94 indep_names = std::vector<std::string>();
95
96 tria_listener = dof_handler.get_triangulation().signals.any_change.connect(
97 [this]() { this->tria_change_listener(); });
98}
99
100
101
102template <int dim>
104 const PointValueHistory &point_value_history)
105{
106 dataset_key = point_value_history.dataset_key;
107 independent_values = point_value_history.independent_values;
108 indep_names = point_value_history.indep_names;
109 data_store = point_value_history.data_store;
110 component_mask = point_value_history.component_mask;
111 component_names_map = point_value_history.component_names_map;
112 point_geometry_data = point_value_history.point_geometry_data;
113
114 closed = point_value_history.closed;
115 cleared = point_value_history.cleared;
116
117 dof_handler = point_value_history.dof_handler;
118
119 triangulation_changed = point_value_history.triangulation_changed;
120 have_dof_handler = point_value_history.have_dof_handler;
121 n_indep = point_value_history.n_indep;
122
123 // What to do with tria_listener?
124 // Presume subscribe new instance?
126 {
128 dof_handler->get_triangulation().signals.any_change.connect(
129 [this]() { this->tria_change_listener(); });
130 }
131}
132
133
134
135template <int dim>
138{
139 dataset_key = point_value_history.dataset_key;
140 independent_values = point_value_history.independent_values;
141 indep_names = point_value_history.indep_names;
142 data_store = point_value_history.data_store;
143 component_mask = point_value_history.component_mask;
144 component_names_map = point_value_history.component_names_map;
145 point_geometry_data = point_value_history.point_geometry_data;
146
147 closed = point_value_history.closed;
148 cleared = point_value_history.cleared;
149
150 dof_handler = point_value_history.dof_handler;
151
152 triangulation_changed = point_value_history.triangulation_changed;
153 have_dof_handler = point_value_history.have_dof_handler;
154 n_indep = point_value_history.n_indep;
155
156 // What to do with tria_listener?
157 // Presume subscribe new instance?
159 {
161 dof_handler->get_triangulation().signals.any_change.connect(
162 [this]() { this->tria_change_listener(); });
163 }
164
165 return *this;
166}
167
168
169
170template <int dim>
172{
174 {
175 tria_listener.disconnect();
176 }
177}
178
179
180
181template <int dim>
182void
184{
185 // can't be closed to add additional points
186 // or vectors
191
192 // Implementation assumes that support
193 // points locations are dofs locations
194 AssertThrow(dof_handler->get_fe().has_support_points(), ExcNotImplemented());
195
196 // While in general quadrature points seems
197 // to refer to Gauss quadrature points, in
198 // this case the quadrature points are
199 // forced to be the support points of the
200 // FE.
201 Quadrature<dim> support_point_quadrature(
202 dof_handler->get_fe().get_unit_support_points());
203 FEValues<dim> fe_values(dof_handler->get_fe(),
204 support_point_quadrature,
206 unsigned int n_support_points =
207 dof_handler->get_fe().get_unit_support_points().size();
208 unsigned int n_components = dof_handler->get_fe(0).n_components();
209
210 // set up a loop over all the cells in the
211 // DoFHandler
213 dof_handler->begin_active();
215
216 // default values to be replaced as closer
217 // points are found however they need to be
218 // consistent in case they are actually
219 // chosen
220 typename DoFHandler<dim>::active_cell_iterator current_cell = cell;
221 std::vector<unsigned int> current_fe_index(n_components,
222 0); // need one index per component
223 fe_values.reinit(cell);
224 std::vector<Point<dim>> current_points(n_components, Point<dim>());
225 for (unsigned int support_point = 0; support_point < n_support_points;
226 support_point++)
227 {
228 // set up valid data in the empty vectors
229 unsigned int component =
230 dof_handler->get_fe().system_to_component_index(support_point).first;
231 current_points[component] = fe_values.quadrature_point(support_point);
232 current_fe_index[component] = support_point;
233 }
234
235 // check each cell to find a suitable
236 // support points
237 // GridTools::find_active_cell_around_point
238 // is an alternative. That method is not
239 // used here mostly because of the history
240 // of the class. The algorithm used in
241 // add_points below may be slightly more
242 // efficient than find_active_cell_around_point
243 // because it operates on a set of points.
244
245 for (; cell != endc; ++cell)
246 {
247 fe_values.reinit(cell);
248
249 for (unsigned int support_point = 0; support_point < n_support_points;
250 support_point++)
251 {
252 unsigned int component = dof_handler->get_fe()
253 .system_to_component_index(support_point)
254 .first;
255 const Point<dim> &test_point =
256 fe_values.quadrature_point(support_point);
257
258 if (location.distance(test_point) <
259 location.distance(current_points[component]))
260 {
261 // save the data
262 current_points[component] = test_point;
263 current_cell = cell;
264 current_fe_index[component] = support_point;
265 }
266 }
267 }
268
269
270 std::vector<types::global_dof_index> local_dof_indices(
271 dof_handler->get_fe().n_dofs_per_cell());
272 std::vector<types::global_dof_index> new_solution_indices;
273 current_cell->get_dof_indices(local_dof_indices);
274 // there is an implicit assumption here
275 // that all the closest support point to
276 // the requested point for all finite
277 // element components lie in the same cell.
278 // this could possibly be violated if
279 // components use different FE orders,
280 // requested points are on the edge or
281 // vertex of a cell and we are unlucky with
282 // floating point rounding. Worst case
283 // scenario however is that the point
284 // selected isn't the closest possible, it
285 // will still lie within one cell distance.
286 // calling
287 // GridTools::find_active_cell_around_point
288 // to obtain a cell to search is an
289 // option for these methods, but currently
290 // the GridTools function does not cater for
291 // a vector of points, and does not seem to
292 // be intrinsicly faster than this method.
293 new_solution_indices.reserve(dof_handler->get_fe(0).n_components());
294 for (unsigned int component = 0;
295 component < dof_handler->get_fe(0).n_components();
296 component++)
297 {
298 new_solution_indices.push_back(
299 local_dof_indices[current_fe_index[component]]);
300 }
301
303 new_point_geometry_data(location, current_points, new_solution_indices);
304 point_geometry_data.push_back(new_point_geometry_data);
305
306 for (auto &data_entry : data_store)
307 {
308 // add an extra row to each vector entry
309 const ComponentMask &current_mask =
310 (component_mask.find(data_entry.first))->second;
311 unsigned int n_stored = current_mask.n_selected_components();
312 data_entry.second.resize(data_entry.second.size() + n_stored);
313 }
314}
315
316
317
318template <int dim>
319void
320PointValueHistory<dim>::add_points(const std::vector<Point<dim>> &locations)
321{
322 // This algorithm adds points in the same
323 // order as they appear in the vector
324 // locations and users may depend on this
325 // so do not change order added!
326
327 // can't be closed to add additional points or vectors
332
333
334 // Implementation assumes that support
335 // points locations are dofs locations
336 AssertThrow(dof_handler->get_fe().has_support_points(), ExcNotImplemented());
337
338 // While in general quadrature points seems
339 // to refer to Gauss quadrature points, in
340 // this case the quadrature points are
341 // forced to be the support points of the
342 // FE.
343 Quadrature<dim> support_point_quadrature(
344 dof_handler->get_fe().get_unit_support_points());
345 FEValues<dim> fe_values(dof_handler->get_fe(),
346 support_point_quadrature,
348 unsigned int n_support_points =
349 dof_handler->get_fe().get_unit_support_points().size();
350 unsigned int n_components = dof_handler->get_fe(0).n_components();
351
352 // set up a loop over all the cells in the
353 // DoFHandler
355 dof_handler->begin_active();
357
358 // default values to be replaced as closer
359 // points are found however they need to be
360 // consistent in case they are actually
361 // chosen vector <vector>s defined where
362 // previously single vectors were used
363
364 // need to store one value per point per component
365 std::vector<typename DoFHandler<dim>::active_cell_iterator> current_cell(
366 locations.size(), cell);
367
368 fe_values.reinit(cell);
369 std::vector<Point<dim>> temp_points(n_components, Point<dim>());
370 std::vector<unsigned int> temp_fe_index(n_components, 0);
371 for (unsigned int support_point = 0; support_point < n_support_points;
372 support_point++)
373 {
374 // set up valid data in the empty vectors
375 unsigned int component =
376 dof_handler->get_fe().system_to_component_index(support_point).first;
377 temp_points[component] = fe_values.quadrature_point(support_point);
378 temp_fe_index[component] = support_point;
379 }
380 std::vector<std::vector<Point<dim>>> current_points(
381 locations.size(), temp_points); // give a valid start point
382 std::vector<std::vector<unsigned int>> current_fe_index(locations.size(),
383 temp_fe_index);
384
385 // check each cell to find suitable support
386 // points
387 // GridTools::find_active_cell_around_point
388 // is an alternative. That method is not
389 // used here mostly because of the history
390 // of the class. The algorithm used here
391 // may be slightly more
392 // efficient than find_active_cell_around_point
393 // because it operates on a set of points.
394 for (; cell != endc; ++cell)
395 {
396 fe_values.reinit(cell);
397 for (unsigned int support_point = 0; support_point < n_support_points;
398 support_point++)
399 {
400 unsigned int component = dof_handler->get_fe()
401 .system_to_component_index(support_point)
402 .first;
403 const Point<dim> &test_point =
404 fe_values.quadrature_point(support_point);
405
406 for (unsigned int point = 0; point < locations.size(); ++point)
407 {
408 if (locations[point].distance(test_point) <
409 locations[point].distance(current_points[point][component]))
410 {
411 // save the data
412 current_points[point][component] = test_point;
413 current_cell[point] = cell;
414 current_fe_index[point][component] = support_point;
415 }
416 }
417 }
418 }
419
420 std::vector<types::global_dof_index> local_dof_indices(
421 dof_handler->get_fe().n_dofs_per_cell());
422 for (unsigned int point = 0; point < locations.size(); ++point)
423 {
424 current_cell[point]->get_dof_indices(local_dof_indices);
425 std::vector<types::global_dof_index> new_solution_indices;
426
427 new_solution_indices.reserve(dof_handler->get_fe(0).n_components());
428 for (unsigned int component = 0;
429 component < dof_handler->get_fe(0).n_components();
430 component++)
431 {
432 new_solution_indices.push_back(
433 local_dof_indices[current_fe_index[point][component]]);
434 }
435
437 new_point_geometry_data(locations[point],
438 current_points[point],
439 new_solution_indices);
440
441 point_geometry_data.push_back(new_point_geometry_data);
442
443 for (auto &data_entry : data_store)
444 {
445 // add an extra row to each vector entry
446 const ComponentMask current_mask =
447 (component_mask.find(data_entry.first))->second;
448 unsigned int n_stored = current_mask.n_selected_components();
449 data_entry.second.resize(data_entry.second.size() + n_stored);
450 }
451 }
452}
453
454
455
456template <int dim>
457void
458PointValueHistory<dim>::add_field_name(const std::string &vector_name,
459 const ComponentMask &mask)
460{
461 // can't be closed to add additional points
462 // or vectors
467
468 // insert a component mask that is always of the right size
469 if (mask.represents_the_all_selected_mask() == false)
470 component_mask.insert(std::make_pair(vector_name, mask));
471 else
472 component_mask.insert(
473 std::make_pair(vector_name,
474 ComponentMask(std::vector<bool>(
475 dof_handler->get_fe(0).n_components(), true))));
476
477 // insert an empty vector of strings
478 // to ensure each field has an entry
479 // in the map
480 std::pair<std::string, std::vector<std::string>> empty_names(
481 vector_name, std::vector<std::string>());
482 component_names_map.insert(empty_names);
483
484 // make and add a new vector
485 // point_geometry_data.size() long
486 std::pair<std::string, std::vector<std::vector<double>>> pair_data;
487 pair_data.first = vector_name;
488 const unsigned int n_stored =
489 (mask.represents_the_all_selected_mask() == false ?
490 mask.n_selected_components() :
491 dof_handler->get_fe(0).n_components());
492
493 int n_datastreams =
494 point_geometry_data.size() * n_stored; // each point has n_stored sub parts
495 std::vector<std::vector<double>> vector_size(n_datastreams,
496 std::vector<double>(0));
497 pair_data.second = vector_size;
498 data_store.insert(pair_data);
499}
500
501
502template <int dim>
503void
504PointValueHistory<dim>::add_field_name(const std::string &vector_name,
505 const unsigned int n_components)
506{
507 ComponentMask temp_mask(std::vector<bool>(n_components, true));
508 add_field_name(vector_name, temp_mask);
509}
510
511
512template <int dim>
513void
515 const std::string &vector_name,
516 const std::vector<std::string> &component_names)
517{
518 typename std::map<std::string, std::vector<std::string>>::iterator names =
519 component_names_map.find(vector_name);
520 Assert(names != component_names_map.end(),
521 ExcMessage("vector_name not in class"));
522
523 typename std::map<std::string, ComponentMask>::iterator mask =
524 component_mask.find(vector_name);
525 Assert(mask != component_mask.end(), ExcMessage("vector_name not in class"));
526 unsigned int n_stored = mask->second.n_selected_components();
527 Assert(component_names.size() == n_stored,
528 ExcDimensionMismatch(component_names.size(), n_stored));
529
530 names->second = component_names;
531}
532
533
534template <int dim>
535void
537 const std::vector<std::string> &independent_names)
538{
539 Assert(independent_names.size() == n_indep,
540 ExcDimensionMismatch(independent_names.size(), n_indep));
541
542 indep_names = independent_names;
543}
544
545
546template <int dim>
547void
549{
550 closed = true;
551}
552
553
554
555template <int dim>
556void
558{
559 cleared = true;
560 dof_handler = nullptr;
561 have_dof_handler = false;
562}
563
564// Need to test that the internal data has a full and complete dataset for
565// each key. That is that the data has not got 'out of sync'. Testing that
566// dataset_key is within 1 of independent_values is cheap and is done in all
567// three methods. Evaluate_field will check that its vector_name is within 1
568// of dataset_key. However this leaves the possibility that the user has
569// neglected to call evaluate_field on one vector_name consistently. To catch
570// this case start_new_dataset will call bool deap_check () which will test
571// all vector_names and return a bool. This can be called from an Assert
572// statement.
573
574
575
576template <int dim>
577template <typename VectorType>
578void
579PointValueHistory<dim>::evaluate_field(const std::string &vector_name,
580 const VectorType &solution)
581{
582 // must be closed to add data to internal
583 // members.
588
589 if (n_indep != 0) // hopefully this will get optimized, can't test
590 // independent_values[0] unless n_indep > 0
591 {
592 Assert(std::abs(static_cast<int>(dataset_key.size()) -
593 static_cast<int>(independent_values[0].size())) < 2,
595 }
596 // Look up the field name and get an
597 // iterator for the map. Doing this
598 // up front means that it only needs
599 // to be done once and also allows us
600 // to check vector_name is in the map.
601 typename std::map<std::string, std::vector<std::vector<double>>>::iterator
602 data_store_field = data_store.find(vector_name);
603 Assert(data_store_field != data_store.end(),
604 ExcMessage("vector_name not in class"));
605 // Repeat for component_mask
606 typename std::map<std::string, ComponentMask>::iterator mask =
607 component_mask.find(vector_name);
608 Assert(mask != component_mask.end(), ExcMessage("vector_name not in class"));
609
610 unsigned int n_stored =
611 mask->second.n_selected_components(dof_handler->get_fe(0).n_components());
612
613 typename std::vector<
615 point = point_geometry_data.begin();
616 for (unsigned int data_store_index = 0; point != point_geometry_data.end();
617 ++point, ++data_store_index)
618 {
619 // Look up the components to add
620 // in the component_mask, and
621 // access the data associated with
622 // those components
623
624 for (unsigned int store_index = 0, comp = 0;
625 comp < dof_handler->get_fe(0).n_components();
626 comp++)
627 {
628 if (mask->second[comp])
629 {
630 unsigned int solution_index = point->solution_indices[comp];
631 data_store_field
632 ->second[data_store_index * n_stored + store_index]
633 .push_back(
635 solution_index));
636 ++store_index;
637 }
638 }
639 }
640}
641
642
643
644template <int dim>
645template <typename VectorType>
646void
648 const std::vector<std::string> &vector_names,
649 const VectorType &solution,
650 const DataPostprocessor<dim> &data_postprocessor,
651 const Quadrature<dim> &quadrature)
652{
653 // must be closed to add data to internal
654 // members.
658 if (n_indep != 0) // hopefully this will get optimized, can't test
659 // independent_values[0] unless n_indep > 0
660 {
661 Assert(std::abs(static_cast<int>(dataset_key.size()) -
662 static_cast<int>(independent_values[0].size())) < 2,
664 }
665
666 // Make an FEValues object
667 const UpdateFlags update_flags =
669 Assert(
670 !(update_flags & update_normal_vectors),
672 "The update of normal vectors may not be requested for evaluation of "
673 "data on cells via DataPostprocessor."));
674 FEValues<dim> fe_values(dof_handler->get_fe(), quadrature, update_flags);
675 unsigned int n_components = dof_handler->get_fe(0).n_components();
676 unsigned int n_quadrature_points = quadrature.size();
677
678 unsigned int n_output_variables = data_postprocessor.get_names().size();
679
680 // declare some temp objects for evaluating the solution at quadrature
681 // points. we will either need the scalar or vector version in the code
682 // below
683 std::vector<typename VectorType::value_type> scalar_solution_values(
684 n_quadrature_points);
685 std::vector<Tensor<1, dim, typename VectorType::value_type>>
686 scalar_solution_gradients(n_quadrature_points);
687 std::vector<Tensor<2, dim, typename VectorType::value_type>>
688 scalar_solution_hessians(n_quadrature_points);
689
690 std::vector<Vector<typename VectorType::value_type>> vector_solution_values(
691 n_quadrature_points, Vector<typename VectorType::value_type>(n_components));
692
693 std::vector<std::vector<Tensor<1, dim, typename VectorType::value_type>>>
694 vector_solution_gradients(
695 n_quadrature_points,
698
699 std::vector<std::vector<Tensor<2, dim, typename VectorType::value_type>>>
700 vector_solution_hessians(
701 n_quadrature_points,
704
705 // Loop over points and find correct cell
706 typename std::vector<
708 point = point_geometry_data.begin();
709 Assert(!dof_handler->get_triangulation().is_mixed_mesh(),
711 const auto reference_cell =
712 dof_handler->get_triangulation().get_reference_cells()[0];
713 for (unsigned int data_store_index = 0; point != point_geometry_data.end();
714 ++point, ++data_store_index)
715 {
716 // we now have a point to query, need to know what cell it is in
717 const Point<dim> requested_location = point->requested_location;
718 const typename DoFHandler<dim>::active_cell_iterator cell =
720 reference_cell.template get_default_linear_mapping<dim, dim>(),
722 requested_location)
723 .first;
724
725
726 fe_values.reinit(cell);
727 std::vector<Vector<double>> computed_quantities(
728 1, Vector<double>(n_output_variables)); // just one point needed
729
730 // find the closest quadrature point
731 std::vector<Point<dim>> quadrature_points =
732 fe_values.get_quadrature_points();
733 double distance = cell->diameter();
734 unsigned int selected_point = 0;
735 for (unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
736 {
737 if (requested_location.distance(quadrature_points[q_point]) <
738 distance)
739 {
740 selected_point = q_point;
741 distance =
742 requested_location.distance(quadrature_points[q_point]);
743 }
744 }
745
746
747 // The case of a scalar FE
748 if (n_components == 1)
749 {
750 // Extract data for the DataPostprocessor object
751 DataPostprocessorInputs::Scalar<dim> postprocessor_input;
752
753 // for each quantity that is requested (values, gradients, hessians),
754 // first get them at all quadrature points, then restrict to the
755 // one value on the quadrature point closest to the evaluation
756 // point in question
757 //
758 // we need temporary objects because the underlying scalar
759 // type of the solution vector may be different from 'double',
760 // but the DataPostprocessorInputs only allow for 'double'
761 // data types
762 if (update_flags & update_values)
763 {
764 fe_values.get_function_values(solution, scalar_solution_values);
765 postprocessor_input.solution_values =
766 std::vector<double>(1, scalar_solution_values[selected_point]);
767 }
768 if (update_flags & update_gradients)
769 {
770 fe_values.get_function_gradients(solution,
771 scalar_solution_gradients);
772 postprocessor_input.solution_gradients =
773 std::vector<Tensor<1, dim>>(
774 1, scalar_solution_gradients[selected_point]);
775 }
776 if (update_flags & update_hessians)
777 {
778 fe_values.get_function_hessians(solution,
779 scalar_solution_hessians);
780 postprocessor_input.solution_hessians =
781 std::vector<Tensor<2, dim>>(
782 1, scalar_solution_hessians[selected_point]);
783 }
784
785 // then also set the single evaluation point
786 postprocessor_input.evaluation_points =
787 std::vector<Point<dim>>(1, quadrature_points[selected_point]);
788
789 // and finally do the postprocessing
790 data_postprocessor.evaluate_scalar_field(postprocessor_input,
791 computed_quantities);
792 }
793 else // The case of a vector FE
794 {
795 // exact same idea as above
796 DataPostprocessorInputs::Vector<dim> postprocessor_input;
797
798 if (update_flags & update_values)
799 {
800 fe_values.get_function_values(solution, vector_solution_values);
801 postprocessor_input.solution_values.resize(
802 1, Vector<double>(n_components));
803 std::copy(vector_solution_values[selected_point].begin(),
804 vector_solution_values[selected_point].end(),
805 postprocessor_input.solution_values[0].begin());
806 }
807 if (update_flags & update_gradients)
808 {
809 fe_values.get_function_gradients(solution,
810 vector_solution_gradients);
811 postprocessor_input.solution_gradients.resize(
812 1, std::vector<Tensor<1, dim>>(n_components));
813 std::copy(vector_solution_gradients[selected_point].begin(),
814 vector_solution_gradients[selected_point].end(),
815 postprocessor_input.solution_gradients[0].begin());
816 }
817 if (update_flags & update_hessians)
818 {
819 fe_values.get_function_hessians(solution,
820 vector_solution_hessians);
821 postprocessor_input.solution_hessians.resize(
822 1, std::vector<Tensor<2, dim>>(n_components));
823 std::copy(vector_solution_hessians[selected_point].begin(),
824 vector_solution_hessians[selected_point].end(),
825 postprocessor_input.solution_hessians[0].begin());
826 }
827
828 postprocessor_input.evaluation_points =
829 std::vector<Point<dim>>(1, quadrature_points[selected_point]);
830
831 data_postprocessor.evaluate_vector_field(postprocessor_input,
832 computed_quantities);
833 }
834
835
836 // we now have the data and need to save it
837 // loop over data names
838 typename std::vector<std::string>::const_iterator name =
839 vector_names.begin();
840 for (; name != vector_names.end(); ++name)
841 {
842 typename std::map<std::string,
843 std::vector<std::vector<double>>>::iterator
844 data_store_field = data_store.find(*name);
845 Assert(data_store_field != data_store.end(),
846 ExcMessage("vector_name not in class"));
847 // Repeat for component_mask
848 typename std::map<std::string, ComponentMask>::iterator mask =
849 component_mask.find(*name);
850 Assert(mask != component_mask.end(),
851 ExcMessage("vector_name not in class"));
852
853 unsigned int n_stored =
854 mask->second.n_selected_components(n_output_variables);
855
856 // Push back computed quantities according
857 // to the component_mask.
858 for (unsigned int store_index = 0, comp = 0;
859 comp < n_output_variables;
860 comp++)
861 {
862 if (mask->second[comp])
863 {
864 data_store_field
865 ->second[data_store_index * n_stored + store_index]
866 .push_back(computed_quantities[0](comp));
867 ++store_index;
868 }
869 }
870 }
871 } // end of loop over points
872}
873
874
875
876template <int dim>
877template <typename VectorType>
878void
880 const std::string &vector_name,
881 const VectorType &solution,
882 const DataPostprocessor<dim> &data_postprocessor,
883 const Quadrature<dim> &quadrature)
884{
885 std::vector<std::string> vector_names;
886 vector_names.push_back(vector_name);
887 evaluate_field(vector_names, solution, data_postprocessor, quadrature);
888}
889
890
891
892template <int dim>
893template <typename VectorType>
894void
896 const std::string &vector_name,
897 const VectorType &solution)
898{
899 using number = typename VectorType::value_type;
900 // must be closed to add data to internal
901 // members.
905
906 if (n_indep != 0) // hopefully this will get optimized, can't test
907 // independent_values[0] unless n_indep > 0
908 {
909 Assert(std::abs(static_cast<int>(dataset_key.size()) -
910 static_cast<int>(independent_values[0].size())) < 2,
912 }
913 // Look up the field name and get an
914 // iterator for the map. Doing this
915 // up front means that it only needs
916 // to be done once and also allows us
917 // to check vector_name is in the map.
918 typename std::map<std::string, std::vector<std::vector<double>>>::iterator
919 data_store_field = data_store.find(vector_name);
920 Assert(data_store_field != data_store.end(),
921 ExcMessage("vector_name not in class"));
922 // Repeat for component_mask
923 typename std::map<std::string, ComponentMask>::iterator mask =
924 component_mask.find(vector_name);
925 Assert(mask != component_mask.end(), ExcMessage("vector_name not in class"));
926
927 unsigned int n_stored =
928 mask->second.n_selected_components(dof_handler->get_fe(0).n_components());
929
930 typename std::vector<
932 point = point_geometry_data.begin();
933 Vector<number> value(dof_handler->get_fe(0).n_components());
934 for (unsigned int data_store_index = 0; point != point_geometry_data.end();
935 ++point, ++data_store_index)
936 {
937 // Make a Vector <double> for the value
938 // at the point. It will have as many
939 // components as there are in the fe.
941 solution,
942 point->requested_location,
943 value);
944
945 // Look up the component_mask and add
946 // in components according to that mask
947 for (unsigned int store_index = 0, comp = 0; comp < mask->second.size();
948 comp++)
949 {
950 if (mask->second[comp])
951 {
952 data_store_field
953 ->second[data_store_index * n_stored + store_index]
954 .push_back(value(comp));
955 ++store_index;
956 }
957 }
958 }
959}
960
961
962template <int dim>
963void
965{
966 // must be closed to add data to internal
967 // members.
971
972 dataset_key.push_back(key);
973}
974
975
976
977template <int dim>
978void
980 const std::vector<double> &indep_values)
981{
982 // must be closed to add data to internal
983 // members.
986 Assert(indep_values.size() == n_indep,
987 ExcDimensionMismatch(indep_values.size(), n_indep));
989 Assert(std::abs(static_cast<int>(dataset_key.size()) -
990 static_cast<int>(independent_values[0].size())) < 2,
992
993 for (unsigned int component = 0; component < n_indep; ++component)
994 independent_values[component].push_back(indep_values[component]);
995}
996
997
998
999template <int dim>
1000void
1002 const std::string &base_name,
1003 const std::vector<Point<dim>> &postprocessor_locations)
1004{
1008
1009 // write inputs to a file
1010 if (n_indep != 0)
1011 {
1012 std::string filename = base_name + "_indep.gpl";
1013 std::ofstream to_gnuplot(filename);
1014
1015 to_gnuplot << "# Data independent of mesh location\n";
1016
1017 // write column headings
1018 to_gnuplot << "# <Key> ";
1019
1020 if (indep_names.size() > 0)
1021 {
1022 for (const auto &indep_name : indep_names)
1023 {
1024 to_gnuplot << "<" << indep_name << "> ";
1025 }
1026 to_gnuplot << '\n';
1027 }
1028 else
1029 {
1030 for (unsigned int component = 0; component < n_indep; ++component)
1031 {
1032 to_gnuplot << "<Indep_" << component << "> ";
1033 }
1034 to_gnuplot << '\n';
1035 }
1036 // write general data stored
1037 for (unsigned int key = 0; key < dataset_key.size(); ++key)
1038 {
1039 to_gnuplot << dataset_key[key];
1040
1041 for (unsigned int component = 0; component < n_indep; ++component)
1042 {
1043 to_gnuplot << " " << independent_values[component][key];
1044 }
1045 to_gnuplot << '\n';
1046 }
1047
1048 to_gnuplot.close();
1049 }
1050
1051
1052
1053 // write points to a file
1054 if (have_dof_handler)
1055 {
1057 AssertThrow(postprocessor_locations.empty() ||
1058 postprocessor_locations.size() ==
1059 point_geometry_data.size(),
1060 ExcDimensionMismatch(postprocessor_locations.size(),
1061 point_geometry_data.size()));
1062 // We previously required the
1063 // number of dofs to remain the
1064 // same to provide some sort of
1065 // test on the relevance of the
1066 // support point indices stored.
1067 // We now relax that to allow
1068 // adaptive refinement strategies
1069 // to make use of the
1070 // evaluate_field_requested_locations
1071 // method. Note that the support point
1072 // information is not meaningful if
1073 // the number of dofs has changed.
1074 // AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
1075
1076 typename std::vector<internal::PointValueHistoryImplementation::
1077 PointGeometryData<dim>>::iterator point =
1078 point_geometry_data.begin();
1079 for (unsigned int data_store_index = 0;
1080 point != point_geometry_data.end();
1081 ++point, ++data_store_index)
1082 {
1083 // for each point, open a file to
1084 // be written to
1085 std::string filename = base_name + "_" +
1086 Utilities::int_to_string(data_store_index, 2) +
1087 ".gpl"; // store by order pushed back
1088 // due to
1089 // Utilities::int_to_string(data_store_index,
1090 // 2) call, can handle up to 100
1091 // points
1092 std::ofstream to_gnuplot(filename);
1093
1094 // put helpful info about the
1095 // support point into the file as
1096 // comments
1097 to_gnuplot << "# Requested location: " << point->requested_location
1098 << '\n';
1099 to_gnuplot << "# DoF_index : Support location (for each component)\n";
1100 for (unsigned int component = 0;
1101 component < dof_handler->get_fe(0).n_components();
1102 component++)
1103 {
1104 to_gnuplot << "# " << point->solution_indices[component] << " : "
1105 << point->support_point_locations[component] << '\n';
1106 }
1108 to_gnuplot
1109 << "# (Original components and locations, may be invalidated by mesh change.)\n";
1110
1111 if (postprocessor_locations.size() != 0)
1112 {
1113 to_gnuplot << "# Postprocessor location: "
1114 << postprocessor_locations[data_store_index];
1116 to_gnuplot << " (may be approximate)\n";
1117 }
1118 to_gnuplot << "#\n";
1119
1120
1121 // write column headings
1122 to_gnuplot << "# <Key> ";
1123
1124 if (indep_names.size() > 0)
1125 {
1126 for (const auto &indep_name : indep_names)
1127 {
1128 to_gnuplot << "<" << indep_name << "> ";
1129 }
1130 }
1131 else
1132 {
1133 for (unsigned int component = 0; component < n_indep; ++component)
1134 {
1135 to_gnuplot << "<Indep_" << component << "> ";
1136 }
1137 }
1138
1139 for (const auto &data_entry : data_store)
1140 {
1141 typename std::map<std::string, ComponentMask>::iterator mask =
1142 component_mask.find(data_entry.first);
1143 unsigned int n_stored = mask->second.n_selected_components();
1144 std::vector<std::string> names =
1145 (component_names_map.find(data_entry.first))->second;
1146
1147 if (names.size() > 0)
1148 {
1149 AssertThrow(names.size() == n_stored,
1150 ExcDimensionMismatch(names.size(), n_stored));
1151 for (const auto &name : names)
1152 {
1153 to_gnuplot << "<" << name << "> ";
1154 }
1155 }
1156 else
1157 {
1158 for (unsigned int component = 0; component < n_stored;
1159 component++)
1160 {
1161 to_gnuplot << "<" << data_entry.first << "_" << component
1162 << "> ";
1163 }
1164 }
1165 }
1166 to_gnuplot << '\n';
1167
1168 // write data stored for the point
1169 for (unsigned int key = 0; key < dataset_key.size(); ++key)
1170 {
1171 to_gnuplot << dataset_key[key];
1172
1173 for (unsigned int component = 0; component < n_indep; ++component)
1174 {
1175 to_gnuplot << " " << independent_values[component][key];
1176 }
1177
1178 for (const auto &data_entry : data_store)
1179 {
1180 typename std::map<std::string, ComponentMask>::iterator mask =
1181 component_mask.find(data_entry.first);
1182 unsigned int n_stored = mask->second.n_selected_components();
1183
1184 for (unsigned int component = 0; component < n_stored;
1185 component++)
1186 {
1187 to_gnuplot
1188 << " "
1189 << (data_entry.second)[data_store_index * n_stored +
1190 component][key];
1191 }
1192 }
1193 to_gnuplot << '\n';
1194 }
1195
1196 to_gnuplot.close();
1197 }
1198 }
1199}
1200
1201
1202
1203template <int dim>
1206{
1207 // a method to put a one at each point on
1208 // the grid where a location is defined
1212
1213 Vector<double> dof_vector(dof_handler->n_dofs());
1214
1215 typename std::vector<
1217 point = point_geometry_data.begin();
1218 for (; point != point_geometry_data.end(); ++point)
1219 {
1220 for (unsigned int component = 0;
1221 component < dof_handler->get_fe(0).n_components();
1222 component++)
1223 {
1224 dof_vector(point->solution_indices[component]) = 1;
1225 }
1226 }
1227 return dof_vector;
1228}
1229
1230
1231template <int dim>
1232void
1234 std::vector<std::vector<Point<dim>>> &locations)
1235{
1239
1240 std::vector<std::vector<Point<dim>>> actual_points;
1241 typename std::vector<
1243 point = point_geometry_data.begin();
1244
1245 for (; point != point_geometry_data.end(); ++point)
1246 {
1247 actual_points.push_back(point->support_point_locations);
1248 }
1249 locations = actual_points;
1250}
1251
1252
1253
1254template <int dim>
1255void
1257 const Quadrature<dim> &quadrature,
1258 std::vector<Point<dim>> &locations)
1259{
1262
1263 locations = std::vector<Point<dim>>();
1264
1265 FEValues<dim> fe_values(dof_handler->get_fe(),
1266 quadrature,
1268 unsigned int n_quadrature_points = quadrature.size();
1269 std::vector<Point<dim>> evaluation_points;
1270
1271 // Loop over points and find correct cell
1272 typename std::vector<
1274 point = point_geometry_data.begin();
1275 Assert(!dof_handler->get_triangulation().is_mixed_mesh(),
1277 const auto reference_cell =
1278 dof_handler->get_triangulation().get_reference_cells()[0];
1279 for (unsigned int data_store_index = 0; point != point_geometry_data.end();
1280 ++point, ++data_store_index)
1281 {
1282 // we now have a point to query,
1283 // need to know what cell it is in
1284 Point<dim> requested_location = point->requested_location;
1287 reference_cell.template get_default_linear_mapping<dim, dim>(),
1288 *dof_handler,
1289 requested_location)
1290 .first;
1291 fe_values.reinit(cell);
1292
1293 evaluation_points = fe_values.get_quadrature_points();
1294 double distance = cell->diameter();
1295 unsigned int selected_point = 0;
1296
1297 for (unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
1298 {
1299 if (requested_location.distance(evaluation_points[q_point]) <
1300 distance)
1301 {
1302 selected_point = q_point;
1303 distance =
1304 requested_location.distance(evaluation_points[q_point]);
1305 }
1306 }
1307
1308 locations.push_back(evaluation_points[selected_point]);
1309 }
1310}
1311
1312
1313template <int dim>
1314void
1316{
1317 out << "***PointValueHistory status output***\n\n";
1318 out << "Closed: " << closed << '\n';
1319 out << "Cleared: " << cleared << '\n';
1320 out << "Triangulation_changed: " << triangulation_changed << '\n';
1321 out << "Have_dof_handler: " << have_dof_handler << '\n';
1322 out << "Geometric Data" << '\n';
1323
1324 typename std::vector<
1326 point = point_geometry_data.begin();
1327 if (point == point_geometry_data.end())
1328 {
1329 out << "No points stored currently\n";
1330 }
1331 else
1332 {
1333 if (!cleared)
1334 {
1335 for (; point != point_geometry_data.end(); ++point)
1336 {
1337 out << "# Requested location: " << point->requested_location
1338 << '\n';
1339 out << "# DoF_index : Support location (for each component)\n";
1340 for (unsigned int component = 0;
1341 component < dof_handler->get_fe(0).n_components();
1342 component++)
1343 {
1344 out << point->solution_indices[component] << " : "
1345 << point->support_point_locations[component] << '\n';
1346 }
1347 out << '\n';
1348 }
1349 }
1350 else
1351 {
1352 out << "#Cannot access DoF_indices once cleared\n";
1353 }
1354 }
1355 out << '\n';
1356
1357 if (independent_values.size() != 0)
1358 {
1359 out << "Independent value(s): " << independent_values.size() << " : "
1360 << independent_values[0].size() << '\n';
1361 if (indep_names.size() > 0)
1362 {
1363 out << "Names: ";
1364 for (const auto &indep_name : indep_names)
1365 {
1366 out << "<" << indep_name << "> ";
1367 }
1368 out << '\n';
1369 }
1370 }
1371 else
1372 {
1373 out << "No independent values stored\n";
1374 }
1375
1376 if (data_store.begin() != data_store.end())
1377 {
1378 out
1379 << "Mnemonic: data set size (mask size, n true components) : n data sets\n";
1380 }
1381 for (const auto &data_entry : data_store)
1382 {
1383 // Find field mnemonic
1384 std::string vector_name = data_entry.first;
1385 typename std::map<std::string, ComponentMask>::iterator mask =
1386 component_mask.find(vector_name);
1387 Assert(mask != component_mask.end(),
1388 ExcMessage("vector_name not in class"));
1389 typename std::map<std::string, std::vector<std::string>>::iterator
1390 component_names = component_names_map.find(vector_name);
1391 Assert(component_names != component_names_map.end(),
1392 ExcMessage("vector_name not in class"));
1393
1394 if (data_entry.second.size() != 0)
1395 {
1396 out << data_entry.first << ": " << data_entry.second.size() << " (";
1397 out << mask->second.size() << ", "
1398 << mask->second.n_selected_components() << ") : ";
1399 out << (data_entry.second)[0].size() << '\n';
1400 }
1401 else
1402 {
1403 out << data_entry.first << ": " << data_entry.second.size() << " (";
1404 out << mask->second.size() << ", "
1405 << mask->second.n_selected_components() << ") : ";
1406 out << "No points added" << '\n';
1407 }
1408 // add names, if available
1409 if (component_names->second.size() > 0)
1410 {
1411 for (const auto &name : component_names->second)
1412 {
1413 out << "<" << name << "> ";
1414 }
1415 out << '\n';
1416 }
1417 }
1418 out << '\n';
1419 out << "***end of status output***\n\n";
1420}
1421
1422
1423
1424template <int dim>
1425bool
1427{
1428 // test ways that it can fail, if control
1429 // reaches last statement return true
1430 if (strict)
1431 {
1432 if (n_indep != 0)
1433 {
1434 if (dataset_key.size() != independent_values[0].size())
1435 {
1436 return false;
1437 }
1438 }
1439 if (have_dof_handler)
1440 {
1441 for (const auto &data_entry : data_store)
1442 {
1443 Assert(data_entry.second.size() > 0, ExcInternalError());
1444 if ((data_entry.second)[0].size() != dataset_key.size())
1445 return false;
1446 // this loop only tests one
1447 // member for each name,
1448 // i.e. checks the user it will
1449 // not catch internal errors
1450 // which do not update all
1451 // fields for a name.
1452 }
1453 }
1454 return true;
1455 }
1456 if (n_indep != 0)
1457 {
1458 if (std::abs(static_cast<int>(dataset_key.size()) -
1459 static_cast<int>(independent_values[0].size())) >= 2)
1460 {
1461 return false;
1462 }
1463 }
1464
1465 if (have_dof_handler)
1466 {
1467 for (const auto &data_entry : data_store)
1468 {
1469 Assert(data_entry.second.size() > 0, ExcInternalError());
1470
1471 if (std::abs(static_cast<int>((data_entry.second)[0].size()) -
1472 static_cast<int>(dataset_key.size())) >= 2)
1473 return false;
1474 // this loop only tests one member
1475 // for each name, i.e. checks the
1476 // user it will not catch internal
1477 // errors which do not update all
1478 // fields for a name.
1479 }
1480 }
1481 return true;
1482}
1483
1484
1485
1486template <int dim>
1487void
1489{
1490 // this function is called by the
1491 // Triangulation whenever something
1492 // changes, by virtue of having
1493 // attached the function to the
1494 // signal handler in the
1495 // triangulation object
1496
1497 // we record the fact that the mesh
1498 // has changed. we need to take
1499 // this into account next time we
1500 // evaluate the solution
1501 triangulation_changed = true;
1502}
1503
1504
1505// explicit instantiations
1506#include "numerics/point_value_history.inst"
1507
1508
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
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
virtual std::vector< std::string > get_names() const =0
void get_function_values(const ReadVector< Number > &fe_function, std::vector< Number > &values) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
void get_function_hessians(const ReadVector< Number > &fe_function, std::vector< Tensor< 2, spacedim, Number > > &hessians) const
const Point< spacedim > & quadrature_point(const unsigned int q_point) const
void get_function_gradients(const ReadVector< Number > &fe_function, std::vector< Tensor< 1, spacedim, Number > > &gradients) const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
void evaluate_field_at_requested_location(const std::string &name, const VectorType &solution)
void add_field_name(const std::string &vector_name, const ComponentMask &component_mask={})
boost::signals2::connection tria_listener
void add_point(const Point< dim > &location)
std::map< std::string, ComponentMask > component_mask
std::vector< internal::PointValueHistoryImplementation::PointGeometryData< dim > > point_geometry_data
std::map< std::string, std::vector< std::string > > component_names_map
PointValueHistory & operator=(const PointValueHistory &point_value_history)
ObserverPointer< const DoFHandler< dim >, PointValueHistory< dim > > dof_handler
void evaluate_field(const std::string &name, const VectorType &solution)
Vector< double > mark_support_locations()
std::vector< std::string > indep_names
std::vector< double > dataset_key
void write_gnuplot(const std::string &base_name, const std::vector< Point< dim > > &postprocessor_locations=std::vector< Point< dim > >())
void status(std::ostream &out)
void add_independent_names(const std::vector< std::string > &independent_names)
PointValueHistory(const unsigned int n_independent_variables=0)
void get_support_locations(std::vector< std::vector< Point< dim > > > &locations)
std::map< std::string, std::vector< std::vector< double > > > data_store
void add_component_names(const std::string &vector_name, const std::vector< std::string > &component_names)
void start_new_dataset(const double key)
void add_points(const std::vector< Point< dim > > &locations)
std::vector< std::vector< double > > independent_values
void push_back_independent(const std::vector< double > &independent_values)
void get_postprocessor_locations(const Quadrature< dim > &quadrature, std::vector< Point< dim > > &locations)
bool deep_check(const bool strict)
Definition point.h:113
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
unsigned int size() const
PointGeometryData(const Point< dim > &new_requested_location, const std::vector< Point< dim > > &new_locations, const std::vector< types::global_dof_index > &new_sol_indices)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
static ::ExceptionBase & ExcNoIndependent()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcDataLostSync()
static ::ExceptionBase & ExcDoFHandlerRequired()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDoFHandlerChanged()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcInvalidState()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
const bool IsBlockVector< VectorType >::value
@ 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.
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:316
std::pair< typename MeshType< dim, spacedim >::active_cell_iterator, Point< dim > > find_active_cell_around_point(const Mapping< dim, spacedim > &mapping, const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p, const std::vector< bool > &marked_vertices={}, const double tolerance=1.e-10)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:466
void point_value(const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim, double > &point, Vector< typename VectorType::value_type > &value)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
std::vector< Point< spacedim > > evaluation_points
std::vector< double > solution_values
std::vector< Tensor< 2, spacedim > > solution_hessians
std::vector< Tensor< 1, spacedim > > solution_gradients
std::vector< std::vector< Tensor< 2, spacedim > > > solution_hessians
std::vector<::Vector< double > > solution_values
std::vector< std::vector< Tensor< 1, spacedim > > > solution_gradients
static VectorType::value_type get(const VectorType &V, const types::global_dof_index i)