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