Reference documentation for deal.II version 8.4.2
point_value_history.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2015 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/parallel_vector.h>
20 #include <deal.II/lac/parallel_block_vector.h>
21 #include <deal.II/lac/petsc_vector.h>
22 #include <deal.II/lac/petsc_block_vector.h>
23 #include <deal.II/lac/trilinos_vector.h>
24 #include <deal.II/lac/trilinos_block_vector.h>
25 
26 #include <deal.II/numerics/vector_tools.h>
27 
28 #include <deal.II/numerics/point_value_history.h>
29 
30 #include <algorithm>
31 
32 
33 DEAL_II_NAMESPACE_OPEN
34 
35 
36 namespace internal
37 {
38  namespace PointValueHistory
39  {
41  template <int dim>
43  ::PointGeometryData (const Point <dim> &new_requested_location,
44  const std::vector <Point <dim> > &new_locations,
45  const std::vector <types::global_dof_index> &new_sol_indices)
46  {
47  requested_location = new_requested_location;
48  support_point_locations = new_locations;
49  solution_indices = new_sol_indices;
50  }
51  }
52 }
53 
54 
55 
56 template <int dim>
58 ::PointValueHistory (const unsigned int n_independent_variables) :
59  n_indep (n_independent_variables)
60 {
61  closed = false;
62  cleared = false;
63  triangulation_changed = false;
64  have_dof_handler = false;
65 
66  // make a vector for keys
67  dataset_key = std::vector <double> (); // initialize the std::vector
68 
69  // make a vector of independent values
70  independent_values
71  = std::vector<std::vector <double> > (n_indep, std::vector <double> (0));
72  indep_names = std::vector <std::string> ();
73 }
74 
75 
76 
77 template <int dim>
79  const unsigned int n_independent_variables) :
80  dof_handler (&dof_handler),
81  n_indep (n_independent_variables)
82 {
83  closed = false;
84  cleared = false;
85  triangulation_changed = false;
86  have_dof_handler = true;
87 
88  // make a vector to store keys
89  dataset_key = std::vector <double> (); // initialize the std::vector
90 
91  // make a vector for the independent values
93  = std::vector<std::vector <double> > (n_indep, std::vector <double> (0));
94  indep_names = std::vector <std::string> ();
95 
96  tria_listener = dof_handler.get_triangulation().signals.any_change.connect (std_cxx11::bind (&PointValueHistory<dim>::tria_change_listener,
97  std_cxx11::ref(*this)));
98 }
99 
100 
101 
102 template <int dim>
104 {
105  dataset_key = point_value_history.dataset_key;
106  independent_values = point_value_history.independent_values;
107  indep_names = point_value_history.indep_names;
108  data_store = point_value_history.data_store;
109  component_mask = point_value_history.component_mask;
110  component_names_map = point_value_history.component_names_map;
111  point_geometry_data = point_value_history.point_geometry_data;
112 
113  closed = point_value_history.closed;
114  cleared = point_value_history.cleared;
115 
116  dof_handler = point_value_history.dof_handler;
117 
118  triangulation_changed = point_value_history.triangulation_changed;
119  have_dof_handler = point_value_history.have_dof_handler;
120  n_indep = point_value_history.n_indep;
121 
122  // What to do with tria_listener?
123  // Presume subscribe new instance?
124  if (have_dof_handler)
125  {
127  std_cxx11::ref(*this)));
128  }
129 }
130 
131 
132 
133 template <int dim>
136 {
137  dataset_key = point_value_history.dataset_key;
138  independent_values = point_value_history.independent_values;
139  indep_names = point_value_history.indep_names;
140  data_store = point_value_history.data_store;
141  component_mask = point_value_history.component_mask;
142  component_names_map = point_value_history.component_names_map;
143  point_geometry_data = point_value_history.point_geometry_data;
144 
145  closed = point_value_history.closed;
146  cleared = point_value_history.cleared;
147 
148  dof_handler = point_value_history.dof_handler;
149 
150  triangulation_changed = point_value_history.triangulation_changed;
151  have_dof_handler = point_value_history.have_dof_handler;
152  n_indep = point_value_history.n_indep;
153 
154  // What to do with tria_listener?
155  // Presume subscribe new instance?
156  if (have_dof_handler)
157  {
159  std_cxx11::ref(*this)));
160  }
161 
162  return * this;
163 }
164 
165 
166 
167 template <int dim>
170 {
171  if (have_dof_handler)
172  {
173  tria_listener.disconnect ();
174  }
175 }
176 
177 
178 
179 template <int dim>
181 ::add_point (const Point <dim> &location)
182 {
183  // can't be closed to add additional points
184  // or vectors
185  AssertThrow (!closed, ExcInvalidState ());
186  AssertThrow (!cleared, ExcInvalidState ());
187  AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
188  AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
189 
190  // Implementation assumes that support
191  // points locations are dofs locations
192  AssertThrow (dof_handler->get_fe ().has_support_points (), ExcNotImplemented ());
193 
194  // FEValues object to extract quadrature
195  // points from
196  std::vector <Point <dim> >
197  unit_support_points = dof_handler->get_fe ().get_unit_support_points ();
198 
199  // While in general quadrature points seems
200  // to refer to Gauss quadrature points, in
201  // this case the quadrature points are
202  // forced to be the support points of the
203  // FE.
205  support_point_quadrature (dof_handler->get_fe ().get_unit_support_points ());
206  FEValues<dim> fe_values (dof_handler->get_fe (),
207  support_point_quadrature,
209  unsigned int n_support_points
210  = dof_handler->get_fe ().get_unit_support_points ().size ();
211  unsigned int n_components
213 
214  // set up a loop over all the cells in the
215  // DoFHandler
217  cell = dof_handler->begin_active ();
219  endc = dof_handler->end ();
220 
221  // default values to be replaced as closer
222  // points are found however they need to be
223  // consistent in case they are actually
224  // chosen
225  typename DoFHandler<dim>::active_cell_iterator current_cell = cell;
226  std::vector <unsigned int> current_fe_index (n_components, 0); // need one index per component
227  fe_values.reinit (cell);
228  std::vector <Point <dim> > current_points (n_components, Point <dim > ());
229  for (unsigned int support_point = 0;
230  support_point < n_support_points; support_point++)
231  {
232  // setup valid data in the empty
233  // vectors
234  unsigned int component
235  = dof_handler->get_fe ().system_to_component_index (support_point).first;
236  current_points [component] = fe_values.quadrature_point (support_point);
237  current_fe_index [component] = support_point;
238  }
239 
240  // check each cell to find a suitable
241  // support points
242  // GridTools::find_active_cell_around_point
243  // is an alternative. That method is not
244  // used here mostly because of the history
245  // of the class. The algorithm used in
246  // add_points below may be slightly more
247  // efficient than find_active_cell_around_point
248  // because it operates on a set of points.
249 
250  for (; cell != endc; cell++)
251  {
252  fe_values.reinit (cell);
253 
254  for (unsigned int support_point = 0;
255  support_point < n_support_points; support_point++)
256  {
257  unsigned int component
258  = dof_handler->get_fe ().system_to_component_index (support_point).first;
259  Point<dim> test_point
260  = fe_values.quadrature_point (support_point);
261 
262  if (location.distance (test_point) <
263  location.distance (current_points [component]))
264  {
265  // save the data
266  current_points [component] = test_point;
267  current_cell = cell;
268  current_fe_index [component] = support_point;
269  }
270  }
271  }
272 
273 
274  std::vector<types::global_dof_index>
275  local_dof_indices (dof_handler->get_fe ().dofs_per_cell);
276  std::vector <types::global_dof_index> new_solution_indices;
277  current_cell->get_dof_indices (local_dof_indices);
278  // there is an implicit assumption here
279  // that all the closest support point to
280  // the requested point for all finite
281  // element components lie in the same cell.
282  // this could possibly be violated if
283  // components use different fe orders,
284  // requested points are on the edge or
285  // vertex of a cell and we are unlucky with
286  // floating point rounding. Worst case
287  // scenario however is that the point
288  // selected isn't the closest possible, it
289  // will still lie within one cell distance.
290  // calling
291  // GridTools::find_active_cell_around_point
292  // to obtain a cell to search is an
293  // option for these methods, but currently
294  // the GridTools function does not cater for
295  // a vector of points, and does not seem to
296  // be intrinsicly faster than this method.
297  for (unsigned int component = 0;
298  component < dof_handler->get_fe ().n_components (); component++)
299  {
300  new_solution_indices
301  .push_back (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  std::map <std::string, std::vector <std::vector <double> > >::iterator
309  data_store_begin = data_store.begin ();
310  for (; data_store_begin != data_store.end (); data_store_begin++)
311  {
312  // add an extra row to each vector
313  // entry
314  const ComponentMask &current_mask = (component_mask.find (data_store_begin->first))->second;
315  unsigned int n_stored = current_mask.n_selected_components();
316  for (unsigned int component = 0; component < n_stored; component++)
317  {
318  data_store_begin->second.push_back (std::vector<double> (0));
319  }
320  }
321 }
322 
323 
324 
325 template <int dim>
327 ::add_points (const std::vector <Point <dim> > &locations)
328 {
329  // This algorithm adds points in the same
330  // order as they appear in the vector
331  // locations and users may depend on this
332  // so do not change order added!
333 
334  // can't be closed to add additional points or vectors
335  AssertThrow (!closed, ExcInvalidState ());
336  AssertThrow (!cleared, ExcInvalidState ());
337  AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
338  AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
339 
340 
341  // Implementation assumes that support
342  // points locations are dofs locations
343  AssertThrow (dof_handler->get_fe ().has_support_points (), ExcNotImplemented ());
344 
345  // FEValues object to extract quadrature
346  // points from
347  std::vector <Point <dim> > unit_support_points = dof_handler->get_fe ().get_unit_support_points ();
348 
349  // While in general quadrature points seems
350  // to refer to Gauss quadrature points, in
351  // this case the quadrature points are
352  // forced to be the support points of the
353  // FE.
354  Quadrature<dim> support_point_quadrature (dof_handler->get_fe ().get_unit_support_points ());
355  FEValues<dim> fe_values (dof_handler->get_fe (), support_point_quadrature, update_quadrature_points);
356  unsigned int n_support_points = dof_handler->get_fe ().get_unit_support_points ().size ();
357  unsigned int n_components = dof_handler->get_fe ().n_components ();
358 
359  // set up a loop over all the cells in the
360  // DoFHandler
363 
364  // default values to be replaced as closer
365  // points are found however they need to be
366  // consistent in case they are actually
367  // chosen vector <vector>s defined where
368  // previously single vectors were used
369 
370  // need to store one value per point per component
371  std::vector <typename DoFHandler<dim>::active_cell_iterator > current_cell (locations.size (), cell);
372 
373  fe_values.reinit (cell);
374  std::vector <Point <dim> > temp_points (n_components, Point <dim > ());
375  std::vector <unsigned int> temp_fe_index (n_components, 0);
376  for (unsigned int support_point = 0; support_point < n_support_points; support_point++)
377  {
378  // setup valid data in the empty
379  // vectors
380  unsigned int component = dof_handler->get_fe ().system_to_component_index (support_point).first;
381  temp_points [component] = fe_values.quadrature_point (support_point);
382  temp_fe_index [component] = support_point;
383  }
384  std::vector <std::vector <Point <dim> > > current_points (locations.size (), temp_points); // give a valid start point
385  std::vector <std::vector <unsigned int> > current_fe_index (locations.size (), temp_fe_index);
386 
387  // check each cell to find suitable support
388  // points
389  // GridTools::find_active_cell_around_point
390  // is an alternative. That method is not
391  // used here mostly because of the history
392  // of the class. The algorithm used here
393  // may be slightly more
394  // efficient than find_active_cell_around_point
395  // because it operates on a set of points.
396  for (; cell != endc; cell++)
397  {
398  fe_values.reinit (cell);
399  for (unsigned int support_point = 0; support_point < n_support_points; support_point++)
400  {
401  unsigned int component = dof_handler->get_fe ().system_to_component_index (support_point).first;
402  Point<dim> test_point = fe_values.quadrature_point (support_point);
403 
404  for (unsigned int point = 0; point < locations.size (); point++)
405  {
406  if (locations[point].distance (test_point) < locations[point].distance (current_points[point][component]))
407  {
408  // save the data
409  current_points[point][component] = test_point;
410  current_cell[point] = cell;
411  current_fe_index[point][component] = support_point;
412  }
413  }
414  }
415  }
416 
417  std::vector<types::global_dof_index> local_dof_indices (dof_handler->get_fe ().dofs_per_cell);
418  for (unsigned int point = 0; point < locations.size (); point++)
419  {
420  current_cell[point]->get_dof_indices (local_dof_indices);
421  std::vector<types::global_dof_index> new_solution_indices;
422 
423  for (unsigned int component = 0; component < dof_handler->get_fe ().n_components (); component++)
424  {
425  new_solution_indices.push_back (local_dof_indices[current_fe_index[point][component]]);
426  }
427 
428  internal::PointValueHistory::PointGeometryData<dim> new_point_geometry_data (locations[point], current_points[point], new_solution_indices);
429 
430  point_geometry_data.push_back (new_point_geometry_data);
431 
432  std::map <std::string, std::vector <std::vector <double> > >::iterator
433  data_store_begin = data_store.begin ();
434  for (; data_store_begin != data_store.end (); data_store_begin++)
435  {
436  // add an extra row to each vector
437  // entry
438  const ComponentMask current_mask = (component_mask.find (data_store_begin->first))->second;
439  unsigned int n_stored = current_mask.n_selected_components();
440  for (unsigned int component = 0; component < n_stored; component++)
441  {
442  data_store_begin->second.push_back (std::vector<double> (0));
443  }
444  }
445  }
446 }
447 
448 
449 
450 
451 
452 
453 template <int dim>
455 ::add_field_name (const std::string &vector_name,
456  const ComponentMask &mask)
457 {
458  // can't be closed to add additional points
459  // or vectors
460  AssertThrow (!closed, ExcInvalidState ());
461  AssertThrow (!cleared, ExcInvalidState ());
462  AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
463  AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
464 
465  // insert a component mask that is always of the right size
466  if (mask.represents_the_all_selected_mask() == false)
467  component_mask.insert (std::make_pair (vector_name, mask));
468  else
469  component_mask.insert (std::make_pair (vector_name,
470  ComponentMask(std::vector<bool>(dof_handler->get_fe().n_components(), true))));
471 
472  // insert an empty vector of strings
473  // to ensure each field has an entry
474  // in the map
475  std::pair <std::string, std::vector <std::string> >
476  empty_names (vector_name, std::vector <std::string> ());
477  component_names_map.insert (empty_names);
478 
479  // make and add a new vector
480  // point_geometry_data.size() long
481  std::pair<std::string, std::vector <std::vector <double> > > pair_data;
482  pair_data.first = vector_name;
483  const unsigned int n_stored = (mask.represents_the_all_selected_mask() == false
484  ?
485  mask.n_selected_components()
486  :
488 
489  int n_datastreams = point_geometry_data.size () * n_stored; // each point has n_stored sub parts
490  std::vector < std::vector <double> > vector_size (n_datastreams,
491  std::vector <double> (0));
492  pair_data.second = vector_size;
493  data_store.insert (pair_data);
494 }
495 
496 
497 template <int dim>
499 ::add_field_name(const std::string &vector_name, const unsigned int n_components)
500 {
501  std::vector <bool> temp_mask (n_components, true);
502  add_field_name (vector_name, temp_mask);
503 }
504 
505 
506 template <int dim>
508 ::add_component_names(const std::string &vector_name,
509  const std::vector <std::string> &component_names)
510 {
511  typename std::map <std::string, std::vector <std::string> >::iterator names = component_names_map.find(vector_name);
512  Assert (names != component_names_map.end(), ExcMessage("vector_name not in class"));
513 
514  typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(vector_name);
515  Assert (mask != component_mask.end(), ExcMessage("vector_name not in class"));
516  unsigned int n_stored = mask->second.n_selected_components();
517  (void)n_stored;
518  Assert (component_names.size() == n_stored, ExcDimensionMismatch (component_names.size(), n_stored));
519 
520  names->second = component_names;
521 }
522 
523 
524 template <int dim>
526 ::add_independent_names(const std::vector <std::string> &independent_names)
527 {
528  Assert (independent_names.size() == n_indep, ExcDimensionMismatch (independent_names.size(), n_indep));
529 
530  indep_names = independent_names;
531 }
532 
533 
534 template <int dim>
537 {
538  closed = true;
539 }
540 
541 
542 
543 template <int dim>
546 {
547  cleared = true;
548  dof_handler = 0;
549  have_dof_handler = false;
550 }
551 
552 // Need to test that the internal data has a full and complete dataset for
553 // each key. That is that the data has not got 'out of sync'. Testing that
554 // dataset_key is within 1 of independent_values is cheap and is done in all
555 // three methods. Evaluate_field will check that its vector_name is within 1
556 // of dataset_key. However this leaves the possibility that the user has
557 // neglected to call evaluate_field on one vector_name consistently. To catch
558 // this case start_new_dataset will call bool deap_check () which will test
559 // all vector_names and return a bool. This can be called from an Assert
560 // statement.
561 
562 
563 
564 template <int dim>
565 template <typename VectorType>
567 ::evaluate_field (const std::string &vector_name, const VectorType &solution)
568 {
569  // must be closed to add data to internal
570  // members.
571  Assert (closed, ExcInvalidState ());
572  Assert (!cleared, ExcInvalidState ());
573  AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
574  AssertThrow (!triangulation_changed, ExcDoFHandlerChanged ());
575 
576  if (n_indep != 0) // hopefully this will get optimized, can't test independent_values[0] unless n_indep > 0
577  {
578  Assert (std::abs ((int) dataset_key.size () - (int) independent_values[0].size ()) < 2, ExcDataLostSync ());
579  }
580  // Look up the field name and get an
581  // iterator for the map. Doing this
582  // up front means that it only needs
583  // to be done once and also allows us
584  // to check vector_name is in the map.
585  typename std::map <std::string, std::vector <std::vector <double> > >::iterator data_store_field = data_store.find(vector_name);
586  Assert (data_store_field != data_store.end(), ExcMessage("vector_name not in class"));
587  // Repeat for component_mask
588  typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(vector_name);
589  Assert (mask != component_mask.end(), ExcMessage("vector_name not in class"));
590 
591  unsigned int n_stored = mask->second.n_selected_components(dof_handler->get_fe ().n_components ());
592 
593  typename std::vector <internal::PointValueHistory::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
594  for (unsigned int data_store_index = 0; point != point_geometry_data.end (); point++, data_store_index++)
595  {
596  // Look up the components to add
597  // in the component_mask, and
598  // access the data associated with
599  // those components
600 
601  for (unsigned int store_index = 0, comp = 0; comp < dof_handler->get_fe ().n_components (); comp++)
602  {
603  if (mask->second[comp])
604  {
605  unsigned int solution_index = point->solution_indices[comp];
606  data_store_field->second[data_store_index * n_stored + store_index].push_back (solution (solution_index));
607  store_index++;
608  }
609  }
610  }
611 }
612 
613 
614 
615 
616 
617 template <int dim>
618 template <typename VectorType>
620 ::evaluate_field (const std::vector <std::string> &vector_names,
621  const VectorType &solution,
622  const DataPostprocessor< dim> &data_postprocessor,
623  const Quadrature<dim> &quadrature)
624 {
625  // must be closed to add data to internal
626  // members.
627  Assert (closed, ExcInvalidState ());
628  Assert (!cleared, ExcInvalidState ());
629  AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
630  if (n_indep != 0) // hopefully this will get optimized, can't test independent_values[0] unless n_indep > 0
631  {
632  Assert (std::abs ((int) dataset_key.size () - (int) independent_values[0].size ()) < 2, ExcDataLostSync ());
633  }
634 
635  // Make an FEValues object
636  const UpdateFlags update_flags = data_postprocessor.get_needed_update_flags() | update_quadrature_points;
637  Assert (!(update_flags & update_normal_vectors),
638  ExcMessage("The update of normal vectors may not be requested for evaluation of "
639  "data on cells via DataPostprocessor."));
640  FEValues<dim> fe_values (dof_handler->get_fe (), quadrature, update_flags);
641  unsigned int n_components = dof_handler->get_fe ().n_components ();
642  unsigned int n_quadrature_points = quadrature.size();
643 
644  unsigned int n_output_variables = data_postprocessor.get_names().size();
645 
646  // Loop over points and find correct cell
647  typename std::vector <internal::PointValueHistory::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
648  for (unsigned int data_store_index = 0; point != point_geometry_data.end (); point++, data_store_index++)
649  {
650  // we now have a point to query,
651  // need to know what cell it is in
652  Point <dim> requested_location = point->requested_location;
654 
655 
656  fe_values.reinit (cell);
657  std::vector< Vector< double > > computed_quantities (1, Vector <double> (n_output_variables)); // just one point needed
658 
659  // The case of a scalar FE
660  if (n_components == 1)
661  {
662  // Extract data for the
663  // PostProcessor object
664  std::vector< typename VectorType::value_type > uh (n_quadrature_points, 0.0);
665  std::vector< Tensor< 1, dim, typename VectorType::value_type > > duh (n_quadrature_points, Tensor <1, dim, typename VectorType::value_type> ());
666  std::vector< Tensor< 2, dim, typename VectorType::value_type > > dduh (n_quadrature_points, Tensor <2, dim, typename VectorType::value_type> ());
667  std::vector<Point<dim> > dummy_normals (1, Point<dim> ());
668  std::vector<Point<dim> > evaluation_points;
669  // at each point there is
670  // only one component of
671  // value, gradient etc.
672  if (update_flags & update_values)
673  fe_values.get_function_values (solution,
674  uh);
675  if (update_flags & update_gradients)
676  fe_values.get_function_gradients (solution,
677  duh);
678  if (update_flags & update_hessians)
679  fe_values.get_function_hessians (solution,
680  dduh);
681 
682  // find the closest quadrature point
683  evaluation_points = fe_values.get_quadrature_points();
684  double distance = cell->diameter ();
685  unsigned int selected_point = 0;
686  for (unsigned int q_point = 0; q_point < n_quadrature_points; q_point++)
687  {
688  if (requested_location.distance (evaluation_points[q_point]) < distance)
689  {
690  selected_point = q_point;
691  distance = requested_location.distance (evaluation_points[q_point]);
692  }
693  }
694 
695  // Call compute_derived_quantities_vector
696  // or compute_derived_quantities_scalar
697  // TODO this function should also operate with typename VectorType::value_type
698  data_postprocessor.
699  compute_derived_quantities_scalar(std::vector< double > (1, uh[selected_point]),
700  std::vector< Tensor< 1, dim > > (1, Tensor< 1, dim >(duh[selected_point]) ),
701  std::vector< Tensor< 2, dim > > (1, Tensor< 2, dim >(dduh[selected_point]) ),
702  dummy_normals,
703  std::vector<Point<dim> > (1, evaluation_points[selected_point]),
704  computed_quantities);
705 
706  }
707  else // The case of a vector FE
708  {
709  // Extract data for the PostProcessor object
710  std::vector< Vector< typename VectorType::value_type > > uh (n_quadrature_points, Vector <typename VectorType::value_type> (n_components));
711  std::vector< std::vector< Tensor< 1, dim, typename VectorType::value_type > > > duh (n_quadrature_points, std::vector< Tensor< 1, dim, typename VectorType::value_type > > (n_components, Tensor< 1, dim, typename VectorType::value_type >()));
712  std::vector< std::vector< Tensor< 2, dim, typename VectorType::value_type > > > dduh (n_quadrature_points, std::vector< Tensor< 2, dim, typename VectorType::value_type > > (n_components, Tensor< 2, dim, typename VectorType::value_type >()));
713  std::vector<Point<dim> > dummy_normals (1, Point<dim> ());
714  std::vector<Point<dim> > evaluation_points;
715  // at each point there is
716  // a vector valued
717  // function and its
718  // derivative...
719  if (update_flags & update_values)
720  fe_values.get_function_values (solution,
721  uh);
722  if (update_flags & update_gradients)
723  fe_values.get_function_gradients (solution,
724  duh);
725  if (update_flags & update_hessians)
726  fe_values.get_function_hessians (solution,
727  dduh);
728 
729  // find the closest quadrature point
730  evaluation_points = 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 (evaluation_points[q_point]) < distance)
736  {
737  selected_point = q_point;
738  distance = requested_location.distance (evaluation_points[q_point]);
739  }
740  }
741 
742  // FIXME: We need tmp vectors below because the data
743  // postprocessors are not equipped to deal with anything but
744  // doubles (scalars and tensors).
745  const Vector< typename VectorType::value_type > &uh_s = uh[selected_point];
746  const std::vector< Tensor< 1, dim, typename VectorType::value_type > > &duh_s = duh[selected_point];
747  const std::vector< Tensor< 2, dim, typename VectorType::value_type > > &dduh_s = dduh[selected_point];
748  std::vector< Tensor< 1, dim > > tmp_d (duh_s.size());
749  for (unsigned int i = 0; i < duh_s.size(); i++)
750  tmp_d[i] = duh_s[i];
751 
752  std::vector< Tensor< 2, dim > > tmp_dd (dduh_s.size());
753  for (unsigned int i = 0; i < dduh_s.size(); i++)
754  tmp_dd[i] = dduh_s[i];
755 
756  Vector< double > tmp(uh_s.size());
757  for (unsigned int i = 0; i < uh_s.size(); i++)
758  tmp[i] = uh_s[i];
759  // Call compute_derived_quantities_vector
760  // or compute_derived_quantities_scalar
761  data_postprocessor.
762  compute_derived_quantities_vector(std::vector< Vector< double > > (1, tmp),
763  std::vector< std::vector< Tensor< 1, dim > > > (1, tmp_d),
764  std::vector< std::vector< Tensor< 2, dim > > > (1, tmp_dd),
765  dummy_normals,
766  std::vector<Point<dim> > (1, evaluation_points[selected_point]),
767  computed_quantities);
768  }
769 
770 
771  // we now have the data and need to save it
772  // loop over data names
773  typename std::vector<std::string>::const_iterator name = vector_names.begin();
774  for (; name != vector_names.end(); name++)
775  {
776  typename std::map <std::string, std::vector <std::vector <double> > >::iterator data_store_field = data_store.find(*name);
777  Assert (data_store_field != data_store.end(), ExcMessage("vector_name not in class"));
778  // Repeat for component_mask
779  typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(*name);
780  Assert (mask != component_mask.end(), ExcMessage("vector_name not in class"));
781 
782  unsigned int n_stored = mask->second.n_selected_components(n_output_variables);
783 
784  // Push back computed quantities according
785  // to the component_mask.
786  for (unsigned int store_index = 0, comp = 0; comp < n_output_variables; comp++)
787  {
788  if (mask->second[comp])
789  {
790  data_store_field->second[data_store_index * n_stored + store_index].push_back (computed_quantities[0](comp));
791  store_index++;
792  }
793  }
794  }
795  } // end of loop over points
796 }
797 
798 
799 template <int dim>
800 template <typename VectorType>
802 ::evaluate_field (const std::string &vector_name,
803  const VectorType &solution,
804  const DataPostprocessor<dim> &data_postprocessor,
805  const Quadrature<dim> &quadrature)
806 {
807  std::vector <std::string> vector_names;
808  vector_names.push_back (vector_name);
809  evaluate_field (vector_names, solution, data_postprocessor, quadrature);
810 }
811 
812 
813 
814 template <int dim>
815 template <typename VectorType>
817 ::evaluate_field_at_requested_location (const std::string &vector_name,
818  const VectorType &solution)
819 {
820  // must be closed to add data to internal
821  // members.
822  Assert (closed, ExcInvalidState ());
823  Assert (!cleared, ExcInvalidState ());
824  AssertThrow (have_dof_handler, ExcDoFHandlerRequired ());
825 
826  if (n_indep != 0) // hopefully this will get optimized, can't test independent_values[0] unless n_indep > 0
827  {
828  Assert (std::abs ((int) dataset_key.size () - (int) independent_values[0].size ()) < 2, ExcDataLostSync ());
829  }
830  // Look up the field name and get an
831  // iterator for the map. Doing this
832  // up front means that it only needs
833  // to be done once and also allows us
834  // to check vector_name is in the map.
835  typename std::map <std::string, std::vector <std::vector <double> > >::iterator data_store_field = data_store.find(vector_name);
836  Assert (data_store_field != data_store.end(), ExcMessage("vector_name not in class"));
837  // Repeat for component_mask
838  typename std::map <std::string, ComponentMask>::iterator mask = component_mask.find(vector_name);
839  Assert (mask != component_mask.end(), ExcMessage("vector_name not in class"));
840 
841  unsigned int n_stored = mask->second.n_selected_components(dof_handler->get_fe ().n_components ());
842 
843  typename std::vector <internal::PointValueHistory::PointGeometryData <dim> >::iterator point = point_geometry_data.begin ();
845  for (unsigned int data_store_index = 0; point != point_geometry_data.end (); point++, data_store_index++)
846  {
847  // Make a Vector <double> for the value
848  // at the point. It will have as many
849  // components as there are in the fe.
850  VectorTools::point_value (*dof_handler, solution, point->requested_location, value);
851 
852  // Look up the component_mask and add
853  // in components according to that mask
854  for (unsigned int store_index = 0, comp = 0; comp < mask->second.size(); comp++)
855  {
856  if (mask->second[comp])
857  {
858  data_store_field->second[data_store_index * n_stored + store_index].push_back (value (comp));
859  store_index++;
860  }
861  }
862  }
863 }
864 
865 
866 template <int dim>
869 {
870  // must be closed to add data to internal
871  // members.
872  Assert (closed, ExcInvalidState ());
873  Assert (!cleared, ExcInvalidState ());
874  Assert (deep_check (false), ExcDataLostSync ());
875 
876  dataset_key.push_back (key);
877 }
878 
879 
880 
881 template <int dim>
883 ::push_back_independent (const std::vector <double> &indep_values)
884 {
885  // must be closed to add data to internal
886  // members.
887  Assert (closed, ExcInvalidState ());
888  Assert (!cleared, ExcInvalidState ());
889  Assert (indep_values.size () == n_indep, ExcDimensionMismatch (indep_values.size (), n_indep));
890  Assert (n_indep != 0, ExcNoIndependent ());
891  Assert (std::abs ((int) dataset_key.size () - (int) independent_values[0].size ()) < 2, ExcDataLostSync ());
892 
893  for (unsigned int component = 0; component < n_indep; component++)
894  independent_values[component].push_back (indep_values[component]);
895 }
896 
897 
898 
899 template <int dim>
901 ::write_gnuplot (const std::string &base_name, 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::PointValueHistory::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 ().n_components (); component++)
988  {
989  to_gnuplot << "# " << point->solution_indices[component] << " : " << point->support_point_locations [component] << "\n";
990  }
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];
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::PointValueHistory::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 ().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::PointValueHistory::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::PointValueHistory::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::PointValueHistory::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 ().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
void get_postprocessor_locations(const Quadrature< dim > &quadrature, std::vector< Point< dim > > &locations)
PointValueHistory(const unsigned int n_independent_variables=0)
::ExceptionBase & ExcMessage(std::string arg1)
cell_iterator end() const
Definition: dof_handler.cc:861
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)
std::vector< double > dataset_key
bool represents_the_all_selected_mask() const
Transformed quadrature points.
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
std::vector< std::string > indep_names
const FiniteElement< dim, spacedim > & get_fe() const
void add_points(const std::vector< Point< dim > > &locations)
void write_gnuplot(const std::string &base_name, const std::vector< Point< dim > > postprocessor_locations=std::vector< Point< dim > >())
bool has_support_points() const
Definition: fe.cc:972
std::vector< internal::PointValueHistory::PointGeometryData< dim > > point_geometry_data
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: dof_handler.cc:845
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
#define Assert(cond, exc)
Definition: exceptions.h:294
Signals signals
Definition: tria.h:2078
UpdateFlags
void evaluate_field_at_requested_location(const std::string &name, const VectorType &solution)
std::vector< std::vector< double > > independent_values
std::size_t size() const
void evaluate_field(const std::string &name, const VectorType &solution)
types::global_dof_index n_dofs() const
void start_new_dataset(const double key)
SmartPointer< const DoFHandler< dim >, PointValueHistory< dim > > dof_handler
PointValueHistory & operator=(const PointValueHistory &point_value_history)
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:82
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
void point_value(const DoFHandler< dim, spacedim > &dof, const VectorType &fe_function, const Point< spacedim > &point, Vector< double > &value)
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:283
const std::vector< Point< dim > > & get_unit_support_points() const
Definition: fe.cc:956
void add_independent_names(const std::vector< std::string > &independent_names)
Vector< double > mark_support_locations()
bool deep_check(const bool strict)
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
Definition: fe.h:2742
unsigned int n_components() const
Definition: mpi.h:48
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)
const Triangulation< dim, spacedim > & get_triangulation() const
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)
MeshType< dim, spacedim >::active_cell_iterator find_active_cell_around_point(const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p)
Definition: grid_tools.cc:1205
virtual std::vector< std::string > get_names() const =0
virtual UpdateFlags get_needed_update_flags() const =0
void push_back_independent(const std::vector< double > &independent_values)