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