Reference documentation for deal.II version 8.5.1
data_out_dof_data.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2016 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 #include <deal.II/base/quadrature_lib.h>
17 #include <deal.II/base/work_stream.h>
18 #include <deal.II/base/memory_consumption.h>
19 #include <deal.II/lac/vector.h>
20 #include <deal.II/lac/block_vector.h>
21 #include <deal.II/lac/la_vector.h>
22 #include <deal.II/lac/la_parallel_vector.h>
23 #include <deal.II/lac/la_parallel_block_vector.h>
24 #include <deal.II/lac/petsc_vector.h>
25 #include <deal.II/lac/petsc_block_vector.h>
26 #include <deal.II/lac/trilinos_vector.h>
27 #include <deal.II/lac/trilinos_block_vector.h>
28 #include <deal.II/numerics/data_out.h>
29 #include <deal.II/numerics/data_out_dof_data.h>
30 #include <deal.II/dofs/dof_handler.h>
31 #include <deal.II/dofs/dof_accessor.h>
32 #include <deal.II/grid/tria.h>
33 #include <deal.II/grid/tria_iterator.h>
34 #include <deal.II/fe/fe.h>
35 #include <deal.II/fe/fe_dgq.h>
36 #include <deal.II/fe/fe_values.h>
37 #include <deal.II/hp/dof_handler.h>
38 #include <deal.II/hp/fe_values.h>
39 #include <deal.II/fe/mapping_q1.h>
40 
41 #include <sstream>
42 
43 DEAL_II_NAMESPACE_OPEN
44 
45 
46 namespace internal
47 {
48  namespace DataOut
49  {
50  template <int dim, int spacedim>
51  ParallelDataBase<dim,spacedim>::
52  ParallelDataBase (const unsigned int n_datasets,
53  const unsigned int n_subdivisions,
54  const std::vector<unsigned int> &n_postprocessor_outputs,
55  const Mapping<dim,spacedim> &mapping,
56  const std::vector<std_cxx11::shared_ptr<::hp::FECollection<dim,spacedim> > > &finite_elements,
57  const UpdateFlags update_flags,
58  const bool use_face_values)
59  :
60  n_datasets (n_datasets),
61  n_subdivisions (n_subdivisions),
62  postprocessed_values (n_postprocessor_outputs.size()),
63  mapping_collection (mapping),
64  finite_elements (finite_elements),
65  update_flags (update_flags)
66  {
67  unsigned int n_q_points = 0;
68  if (use_face_values == false)
69  {
71  quadrature(QIterated<dim>(QTrapez<1>(), n_subdivisions));
72  n_q_points = quadrature[0].size();
73  x_fe_values.resize(this->finite_elements.size());
74  for (unsigned int i=0; i<this->finite_elements.size(); ++i)
75  {
76  // check if there is a finite element that is equal to the
77  // present one, then we can re-use the FEValues object
78  for (unsigned int j=0; j<i; ++j)
79  if (this->finite_elements[i].get() ==
80  this->finite_elements[j].get())
81  {
82  x_fe_values[i] = x_fe_values[j];
83  break;
84  }
85  if (x_fe_values[i].get() == 0)
86  x_fe_values[i].reset(new ::hp::FEValues<dim,spacedim>
87  (this->mapping_collection,
88  *this->finite_elements[i],
89  quadrature,
90  this->update_flags));
91  }
92  }
93  else
94  {
95  ::hp::QCollection<dim-1>
96  quadrature(QIterated<dim-1>(QTrapez<1>(), n_subdivisions));
97  n_q_points = quadrature[0].size();
98  x_fe_face_values.resize(this->finite_elements.size());
99  for (unsigned int i=0; i<this->finite_elements.size(); ++i)
100  {
101  // check if there is a finite element that is equal to the
102  // present one, then we can re-use the FEValues object
103  for (unsigned int j=0; j<i; ++j)
104  if (this->finite_elements[i].get() ==
105  this->finite_elements[j].get())
106  {
107  x_fe_face_values[i] = x_fe_face_values[j];
108  break;
109  }
110  if (x_fe_face_values[i].get() == 0)
111  x_fe_face_values[i].reset(new ::hp::FEFaceValues<dim,spacedim>
112  (this->mapping_collection,
113  *this->finite_elements[i],
114  quadrature,
115  this->update_flags));
116  }
117  }
118 
119  patch_values_scalar.solution_values.resize (n_q_points);
120  patch_values_scalar.solution_gradients.resize (n_q_points);
121  patch_values_scalar.solution_hessians.resize (n_q_points);
122  patch_values_system.solution_values.resize (n_q_points);
123  patch_values_system.solution_gradients.resize (n_q_points);
124  patch_values_system.solution_hessians.resize (n_q_points);
125 
126  for (unsigned int dataset=0; dataset<n_postprocessor_outputs.size(); ++dataset)
127  if (n_postprocessor_outputs[dataset] != 0)
128  postprocessed_values[dataset]
129  .resize(n_q_points,
130  ::Vector<double>(n_postprocessor_outputs[dataset]));
131  }
132 
133 
134 
135 
136 
137  // implement copy constructor to create a thread's own version of
138  // x_fe_values
139  template <int dim, int spacedim>
140  ParallelDataBase<dim,spacedim>::
141  ParallelDataBase (const ParallelDataBase<dim,spacedim> &data)
142  :
143  n_datasets (data.n_datasets),
144  n_subdivisions (data.n_subdivisions),
145  patch_values_scalar (data.patch_values_scalar),
146  patch_values_system (data.patch_values_system),
147  postprocessed_values (data.postprocessed_values),
148  mapping_collection (data.mapping_collection),
149  finite_elements (data.finite_elements),
150  update_flags (data.update_flags)
151  {
152  if (data.x_fe_values.empty() == false)
153  {
154  Assert(data.x_fe_face_values.empty() == true, ExcInternalError());
156  quadrature(QIterated<dim>(QTrapez<1>(), n_subdivisions));
157  x_fe_values.resize(this->finite_elements.size());
158  for (unsigned int i=0; i<this->finite_elements.size(); ++i)
159  {
160  // check if there is a finite element that is equal to the
161  // present one, then we can re-use the FEValues object
162  for (unsigned int j=0; j<i; ++j)
163  if (this->finite_elements[i].get() ==
164  this->finite_elements[j].get())
165  {
166  x_fe_values[i] = x_fe_values[j];
167  break;
168  }
169  if (x_fe_values[i].get() == 0)
170  x_fe_values[i].reset(new ::hp::FEValues<dim,spacedim>
171  (this->mapping_collection,
172  *this->finite_elements[i],
173  quadrature,
174  this->update_flags));
175  }
176  }
177  else
178  {
179  ::hp::QCollection<dim-1>
180  quadrature(QIterated<dim-1>(QTrapez<1>(), n_subdivisions));
181  x_fe_face_values.resize(this->finite_elements.size());
182  for (unsigned int i=0; i<this->finite_elements.size(); ++i)
183  {
184  // check if there is a finite element that is equal to the
185  // present one, then we can re-use the FEValues object
186  for (unsigned int j=0; j<i; ++j)
187  if (this->finite_elements[i].get() ==
188  this->finite_elements[j].get())
189  {
190  x_fe_face_values[i] = x_fe_face_values[j];
191  break;
192  }
193  if (x_fe_face_values[i].get() == 0)
194  x_fe_face_values[i].reset(new ::hp::FEFaceValues<dim,spacedim>
195  (this->mapping_collection,
196  *this->finite_elements[i],
197  quadrature,
198  this->update_flags));
199  }
200  }
201  }
202 
203 
204 
205  template <int dim, int spacedim>
206  template <typename DoFHandlerType>
207  void
208  ParallelDataBase<dim,spacedim>::
209  reinit_all_fe_values(std::vector<std_cxx11::shared_ptr<DataEntryBase<DoFHandlerType> > > &dof_data,
210  const typename ::Triangulation<dim,spacedim>::cell_iterator &cell,
211  const unsigned int face)
212  {
213  for (unsigned int dataset=0; dataset<dof_data.size(); ++dataset)
214  {
215  bool duplicate = false;
216  for (unsigned int j=0; j<dataset; ++j)
217  if (finite_elements[dataset].get() == finite_elements[j].get())
218  duplicate = true;
219  if (duplicate == false)
220  {
221  typename DoFHandlerType::active_cell_iterator dh_cell(&cell->get_triangulation(),
222  cell->level(),
223  cell->index(),
224  dof_data[dataset]->dof_handler);
225  if (x_fe_values.empty())
226  {
227  AssertIndexRange(face,
229  x_fe_face_values[dataset]->reinit(dh_cell, face);
230  }
231  else
232  x_fe_values[dataset]->reinit (dh_cell);
233  }
234  }
235  if (dof_data.empty())
236  {
237  if (x_fe_values.empty())
238  {
239  AssertIndexRange(face,
241  x_fe_face_values[0]->reinit(cell, face);
242  }
243  else
244  x_fe_values[0]->reinit (cell);
245  }
246  }
247 
248 
249 
250  template <int dim, int spacedim>
252  ParallelDataBase<dim,spacedim>::
253  get_present_fe_values(const unsigned int dataset) const
254  {
255  AssertIndexRange(dataset, finite_elements.size());
256  if (x_fe_values.empty())
257  return x_fe_face_values[dataset]->get_present_fe_values();
258  else
259  return x_fe_values[dataset]->get_present_fe_values();
260  }
261 
262 
263 
264  template <int dim, int spacedim>
265  void
266  ParallelDataBase<dim,spacedim>::
267  resize_system_vectors(const unsigned int n_components)
268  {
269  Assert(patch_values_system.solution_values.size() > 0, ExcInternalError());
270  AssertDimension(patch_values_system.solution_values.size(),
271  patch_values_system.solution_gradients.size());
272  AssertDimension(patch_values_system.solution_values.size(),
273  patch_values_system.solution_hessians.size());
274  if (patch_values_system.solution_values[0].size() == n_components)
275  return;
276  for (unsigned int k=0; k<patch_values_system.solution_values.size(); ++k)
277  {
278  patch_values_system.solution_values[k].reinit(n_components);
279  patch_values_system.solution_gradients[k].resize(n_components);
280  patch_values_system.solution_hessians[k].resize(n_components);
281  }
282  }
283 
284 
285 
286 
291  template <int dim, int spacedim>
292  void
293  append_patch_to_list (const DataOutBase::Patch<dim,spacedim> &patch,
294  std::vector<DataOutBase::Patch<dim,spacedim> > &patches)
295  {
296  patches.push_back (patch);
297  patches.back().patch_index = patches.size()-1;
298  }
299  }
300 }
301 
302 namespace internal
303 {
304  namespace DataOut
305  {
306  template <typename DoFHandlerType>
308  (const DoFHandlerType *dofs,
309  const std::vector<std::string> &names_in,
310  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation)
311  :
313  names(names_in),
314  data_component_interpretation (data_component_interpretation),
315  postprocessor(0, typeid(*this).name()),
316  n_output_variables(names.size())
317  {
318  Assert (names.size() == data_component_interpretation.size(),
319  ExcDimensionMismatch(data_component_interpretation.size(),
320  names.size()));
321 
322  // check that the names use only allowed characters
323  for (unsigned int i=0; i<names.size(); ++i)
324  Assert (names[i].find_first_not_of("abcdefghijklmnopqrstuvwxyz"
325  "ABCDEFGHIJKLMNOPQRSTUVWXYZ"
326  "0123456789_<>()") == std::string::npos,
328  names[i].find_first_not_of("abcdefghijklmnopqrstuvwxyz"
329  "ABCDEFGHIJKLMNOPQRSTUVWXYZ"
330  "0123456789_<>()")));
331  }
332 
333 
334 
335  template <typename DoFHandlerType>
337  (const DoFHandlerType *dofs,
338  const DataPostprocessor<DoFHandlerType::space_dimension> *data_postprocessor)
339  :
341  names(data_postprocessor->get_names()),
342  data_component_interpretation (data_postprocessor->get_data_component_interpretation()),
343  postprocessor(data_postprocessor, typeid(*this).name()),
344  n_output_variables(names.size())
345  {
346  Assert (data_postprocessor->get_names().size()
347  ==
348  data_postprocessor->get_data_component_interpretation().size(),
349  ExcDimensionMismatch (data_postprocessor->get_names().size(),
350  data_postprocessor->get_data_component_interpretation().size()));
351 
352  // check that the names use only allowed characters
353  for (unsigned int i=0; i<names.size(); ++i)
354  Assert (names[i].find_first_not_of("abcdefghijklmnopqrstuvwxyz"
355  "ABCDEFGHIJKLMNOPQRSTUVWXYZ"
356  "0123456789_<>()") == std::string::npos,
358  names[i].find_first_not_of("abcdefghijklmnopqrstuvwxyz"
359  "ABCDEFGHIJKLMNOPQRSTUVWXYZ"
360  "0123456789_<>()")));
361  }
362 
363 
364 
365  template <typename DoFHandlerType>
367  {}
368 
369 
370 
377  template <typename DoFHandlerType, typename VectorType>
378  class DataEntry : public DataEntryBase<DoFHandlerType>
379  {
380  public:
386  DataEntry
387  (const DoFHandlerType *dofs,
388  const VectorType *data,
389  const std::vector<std::string> &names,
390  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation);
391 
397  DataEntry (const DoFHandlerType *dofs,
398  const VectorType *data,
399  const DataPostprocessor<DoFHandlerType::space_dimension> *data_postprocessor);
400 
405  virtual
406  double
407  get_cell_data_value (const unsigned int cell_number) const;
408 
413  virtual
414  void
417  std::vector<double> &patch_values) const;
418 
424  virtual
425  void
428  std::vector<::Vector<double> > &patch_values_system) const;
429 
434  virtual
435  void
438  std::vector<Tensor<1,DoFHandlerType::space_dimension> > &patch_gradients) const;
439 
445  virtual
446  void
449  std::vector<std::vector<Tensor<1,DoFHandlerType::space_dimension> > > &patch_gradients_system) const;
450 
455  virtual
456  void
459  std::vector<Tensor<2,DoFHandlerType::space_dimension> > &patch_hessians) const;
460 
466  virtual
467  void
470  std::vector<std::vector< Tensor<2,DoFHandlerType::space_dimension> > > &patch_hessians_system) const;
471 
475  virtual void clear ();
476 
481  virtual std::size_t memory_consumption () const;
482 
483  private:
488  const VectorType *vector;
489  };
490 
491 
492 
493  template <typename DoFHandlerType, typename VectorType>
495  DataEntry (const DoFHandlerType *dofs,
496  const VectorType *data,
497  const std::vector<std::string> &names,
498  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation)
499  :
500  DataEntryBase<DoFHandlerType> (dofs, names, data_component_interpretation),
501  vector (data)
502  {}
503 
504 
505 
506  template <typename DoFHandlerType, typename VectorType>
508  DataEntry (const DoFHandlerType *dofs,
509  const VectorType *data,
510  const DataPostprocessor<DoFHandlerType::space_dimension> *data_postprocessor)
511  :
512  DataEntryBase<DoFHandlerType> (dofs, data_postprocessor),
513  vector (data)
514  {}
515 
516 
517  namespace
518  {
519  template <typename VectorType>
520  double
521  get_vector_element (const VectorType &vector,
522  const unsigned int cell_number)
523  {
524  return vector[cell_number];
525  }
526 
527 
528  double
529  get_vector_element (const IndexSet &is,
530  const unsigned int cell_number)
531  {
532  return (is.is_element(cell_number) ? 1 : 0);
533  }
534  }
535 
536 
537 
538  template <typename DoFHandlerType, typename VectorType>
539  double
541  get_cell_data_value (const unsigned int cell_number) const
542  {
543  return get_vector_element(*vector, cell_number);
544  }
545 
546 
547 
548  template <typename DoFHandlerType, typename VectorType>
549  void
552  std::vector<::Vector<double> > &patch_values_system) const
553  {
554  // FIXME: FEValuesBase gives us data in types that match that of
555  // the solution vector. but this function needs to pass it back
556  // up as 'double' vectors. this requires the use of a temporary
557  // variable here if the data we get is not a 'double' vector.
558  // (of course, in reality, this also means that we may lose
559  // information to begin with.)
560  //
561  // the correct thing would be to also use the correct data type
562  // upstream somewhere, but this is complicated because we hide
563  // the actual data type from upstream. rather, we should at
564  // least make sure we can deal with complex numbers
565  if (typeid(typename VectorType::value_type) == typeid(double))
566  {
567  fe_patch_values.get_function_values (*vector,
568  // reinterpret output argument type; because of
569  // the 'if' statement above, this is the
570  // identity cast whenever the code is
571  // executed, but the cast is necessary
572  // to allow compilation even if we don't get here
573  reinterpret_cast<std::vector<::Vector<typename VectorType::value_type> >&>
574  (patch_values_system));
575  }
576  else
577  {
578  std::vector<::Vector<typename VectorType::value_type> > tmp(patch_values_system.size());
579  for (unsigned int i = 0; i < patch_values_system.size(); i++)
580  tmp[i].reinit(patch_values_system[i]);
581 
582  fe_patch_values.get_function_values (*vector, tmp);
583 
584  for (unsigned int i = 0; i < patch_values_system.size(); i++)
585  patch_values_system[i] = tmp[i];
586  }
587  }
588 
589 
590 
591  template <typename DoFHandlerType, typename VectorType>
592  void
595  std::vector<double> &patch_values) const
596  {
597  // FIXME: FEValuesBase gives us data in types that match that of
598  // the solution vector. but this function needs to pass it back
599  // up as 'double' vectors. this requires the use of a temporary
600  // variable here if the data we get is not a 'double' vector.
601  // (of course, in reality, this also means that we may lose
602  // information to begin with.)
603  //
604  // the correct thing would be to also use the correct data type
605  // upstream somewhere, but this is complicated because we hide
606  // the actual data type from upstream. rather, we should at
607  // least make sure we can deal with complex numbers
608  if (typeid(typename VectorType::value_type) == typeid(double))
609  {
610  fe_patch_values.get_function_values (*vector,
611  // reinterpret output argument type; because of
612  // the 'if' statement above, this is the
613  // identity cast whenever the code is
614  // executed, but the cast is necessary
615  // to allow compilation even if we don't get here
616  reinterpret_cast<std::vector<typename VectorType::value_type>&>
617  (patch_values));
618  }
619  else
620  {
621  std::vector<typename VectorType::value_type> tmp (patch_values.size());
622 
623  fe_patch_values.get_function_values (*vector, tmp);
624 
625  for (unsigned int i = 0; i < tmp.size(); i++)
626  patch_values[i] = tmp[i];
627  }
628  }
629 
630 
631 
632  template <typename DoFHandlerType, typename VectorType>
633  void
636  std::vector<std::vector<Tensor<1,DoFHandlerType::space_dimension> > > &patch_gradients_system) const
637  {
638  // FIXME: FEValuesBase gives us data in types that match that of
639  // the solution vector. but this function needs to pass it back
640  // up as 'double' vectors. this requires the use of a temporary
641  // variable here if the data we get is not a 'double' vector.
642  // (of course, in reality, this also means that we may lose
643  // information to begin with.)
644  //
645  // the correct thing would be to also use the correct data type
646  // upstream somewhere, but this is complicated because we hide
647  // the actual data type from upstream. rather, we should at
648  // least make sure we can deal with complex numbers
649  if (typeid(typename VectorType::value_type) == typeid(double))
650  {
651  fe_patch_values.get_function_gradients (*vector,
652  // reinterpret output argument type; because of
653  // the 'if' statement above, this is the
654  // identity cast whenever the code is
655  // executed, but the cast is necessary
656  // to allow compilation even if we don't get here
657  reinterpret_cast<std::vector<std::vector<Tensor<1,DoFHandlerType::space_dimension,typename VectorType::value_type> > >&>
658  (patch_gradients_system));
659  }
660  else
661  {
662  std::vector<std::vector<Tensor<1,DoFHandlerType::space_dimension,
663  typename VectorType::value_type> > >
664  tmp(patch_gradients_system.size());
665  for (unsigned int i = 0; i < tmp.size(); i++)
666  tmp[i].resize(patch_gradients_system[i].size());
667 
668  fe_patch_values.get_function_gradients (*vector, tmp);
669 
670  for (unsigned int i = 0; i < tmp.size(); i++)
671  for (unsigned int j = 0; j < tmp[i].size(); j++)
672  patch_gradients_system[i][j] = tmp[i][j];
673  }
674  }
675 
676 
677 
678  template <typename DoFHandlerType, typename VectorType>
679  void
682  std::vector<Tensor<1,DoFHandlerType::space_dimension> > &patch_gradients) const
683  {
684  // FIXME: FEValuesBase gives us data in types that match that of
685  // the solution vector. but this function needs to pass it back
686  // up as 'double' vectors. this requires the use of a temporary
687  // variable here if the data we get is not a 'double' vector.
688  // (of course, in reality, this also means that we may lose
689  // information to begin with.)
690  //
691  // the correct thing would be to also use the correct data type
692  // upstream somewhere, but this is complicated because we hide
693  // the actual data type from upstream. rather, we should at
694  // least make sure we can deal with complex numbers
695  if (typeid(typename VectorType::value_type) == typeid(double))
696  {
697  fe_patch_values.get_function_gradients (*vector,
698  // reinterpret output argument type; because of
699  // the 'if' statement above, this is the
700  // identity cast whenever the code is
701  // executed, but the cast is necessary
702  // to allow compilation even if we don't get here
703  reinterpret_cast<std::vector<Tensor<1,DoFHandlerType::space_dimension,typename VectorType::value_type> >&>
704  (patch_gradients));
705  }
706  else
707  {
708  std::vector<Tensor<1,DoFHandlerType::space_dimension,typename VectorType::value_type> > tmp;
709  tmp.resize(patch_gradients.size());
710 
711  fe_patch_values.get_function_gradients (*vector, tmp);
712 
713  for (unsigned int i = 0; i < tmp.size(); i++)
714  patch_gradients[i] = tmp[i];
715  }
716  }
717 
718 
719 
720  template <typename DoFHandlerType, typename VectorType>
721  void
724  std::vector<std::vector<Tensor<2,DoFHandlerType::space_dimension> > > &patch_hessians_system) const
725  {
726  // FIXME: FEValuesBase gives us data in types that match that of
727  // the solution vector. but this function needs to pass it back
728  // up as 'double' vectors. this requires the use of a temporary
729  // variable here if the data we get is not a 'double' vector.
730  // (of course, in reality, this also means that we may lose
731  // information to begin with.)
732  //
733  // the correct thing would be to also use the correct data type
734  // upstream somewhere, but this is complicated because we hide
735  // the actual data type from upstream. rather, we should at
736  // least make sure we can deal with complex numbers
737  if (typeid(typename VectorType::value_type) == typeid(double))
738  {
739  fe_patch_values.get_function_hessians (*vector,
740  // reinterpret output argument type; because of
741  // the 'if' statement above, this is the
742  // identity cast whenever the code is
743  // executed, but the cast is necessary
744  // to allow compilation even if we don't get here
745  reinterpret_cast<std::vector<std::vector<Tensor<2,DoFHandlerType::space_dimension,typename VectorType::value_type> > >&>
746  (patch_hessians_system));
747  }
748  else
749  {
750  std::vector<std::vector<Tensor<2,DoFHandlerType::space_dimension,
751  typename VectorType::value_type> > >
752  tmp(patch_hessians_system.size());
753  for (unsigned int i = 0; i < tmp.size(); i++)
754  tmp[i].resize(patch_hessians_system[i].size());
755 
756  fe_patch_values.get_function_hessians (*vector, tmp);
757 
758  for (unsigned int i = 0; i < tmp.size(); i++)
759  for (unsigned int j = 0; j < tmp[i].size(); j++)
760  patch_hessians_system[i][j] = tmp[i][j];
761  }
762  }
763 
764 
765 
766  template <typename DoFHandlerType, typename VectorType>
767  void
770  std::vector<Tensor<2,DoFHandlerType::space_dimension> > &patch_hessians) const
771  {
772  // FIXME: FEValuesBase gives us data in types that match that of
773  // the solution vector. but this function needs to pass it back
774  // up as 'double' vectors. this requires the use of a temporary
775  // variable here if the data we get is not a 'double' vector.
776  // (of course, in reality, this also means that we may lose
777  // information to begin with.)
778  //
779  // the correct thing would be to also use the correct data type
780  // upstream somewhere, but this is complicated because we hide
781  // the actual data type from upstream. rather, we should at
782  // least make sure we can deal with complex numbers
783  if (typeid(typename VectorType::value_type) == typeid(double))
784  {
785  fe_patch_values.get_function_hessians (*vector,
786  // reinterpret output argument type; because of
787  // the 'if' statement above, this is the
788  // identity cast whenever the code is
789  // executed, but the cast is necessary
790  // to allow compilation even if we don't get here
791  reinterpret_cast<std::vector<Tensor<2,DoFHandlerType
792  ::space_dimension,typename VectorType::value_type> >&>
793  (patch_hessians));
794  }
795  else
796  {
797  std::vector<Tensor<2,DoFHandlerType::space_dimension,typename VectorType::value_type> >
798  tmp(patch_hessians.size());
799 
800  fe_patch_values.get_function_hessians (*vector, tmp);
801 
802  for (unsigned int i = 0; i < tmp.size(); i++)
803  patch_hessians[i] = tmp[i];
804  }
805  }
806 
807 
808 
809  template <typename DoFHandlerType, typename VectorType>
810  std::size_t
812  {
813  return (sizeof (vector) +
815  }
816 
817 
818 
819  template <typename DoFHandlerType, typename VectorType>
820  void
822  {
823  vector = 0;
824  this->dof_handler = 0;
825  }
826  }
827 }
828 
829 
830 
831 template <typename DoFHandlerType,
832  int patch_dim, int patch_space_dim>
834  :
835  triangulation(0,typeid(*this).name()),
836  dofs(0,typeid(*this).name())
837 {}
838 
839 
840 
841 template <typename DoFHandlerType, int patch_dim, int patch_space_dim>
843 {
844  clear ();
845 }
846 
847 
848 
849 template <typename DoFHandlerType, int patch_dim, int patch_space_dim>
850 void
852 attach_dof_handler (const DoFHandlerType &d)
853 {
854  Assert (dof_data.size() == 0,
856  Assert (cell_data.size() == 0,
858 
859  triangulation = SmartPointer<const Triangulation<DoFHandlerType::dimension,
860  DoFHandlerType::space_dimension> >
861  (&d.get_triangulation(), typeid(*this).name());
862  dofs = SmartPointer<const DoFHandlerType>(&d, typeid(*this).name());
863 }
864 
865 
866 
867 template <typename DoFHandlerType, int patch_dim, int patch_space_dim>
868 void
871 {
872  Assert (dof_data.size() == 0,
874  Assert (cell_data.size() == 0,
876 
877  triangulation = SmartPointer<const Triangulation<DoFHandlerType::dimension,
878  DoFHandlerType::space_dimension> >
879  (&tria, typeid(*this).name());
880 }
881 
882 
883 
884 
885 template <typename DoFHandlerType,
886  int patch_dim, int patch_space_dim>
887 template <typename VectorType>
888 void
890 add_data_vector (const VectorType &vec,
891  const std::string &name,
892  const DataVectorType type,
893  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation)
894 {
895  Assert (triangulation != 0,
897  const unsigned int n_components =
898  dofs != 0 ? dofs->get_fe().n_components () : 1;
899 
900  std::vector<std::string> names;
901  // if only one component or vector is cell vector: we only need one name
902  if ((n_components == 1) ||
903  (vec.size() == triangulation->n_active_cells()))
904  {
905  names.resize (1, name);
906  }
907  else
908  // otherwise append _i to the given name
909  {
910  names.resize (n_components);
911  for (unsigned int i=0; i<n_components; ++i)
912  {
913  std::ostringstream namebuf;
914  namebuf << '_' << i;
915  names[i] = name + namebuf.str();
916  }
917  }
918 
919  add_data_vector (vec, names, type, data_component_interpretation);
920 }
921 
922 
923 
924 template <typename DoFHandlerType,
925  int patch_dim, int patch_space_dim>
926 template <typename VectorType>
927 void
929 add_data_vector (const VectorType &vec,
930  const std::vector<std::string> &names,
931  const DataVectorType type,
932  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation_)
933 {
934  Assert (triangulation != 0,
936 
937  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &
938  data_component_interpretation
939  = (data_component_interpretation_.size() != 0
940  ?
941  data_component_interpretation_
942  :
943  std::vector<DataComponentInterpretation::DataComponentInterpretation>
945 
946  // either cell data and one name,
947  // or dof data and n_components names
948  DataVectorType actual_type = type;
949  if (type == type_automatic)
950  {
951  // in the rare case that someone has a DGP(0) attached, we can not decide what she wants here:
952  Assert((dofs == 0) || (triangulation->n_active_cells() != dofs->n_dofs()),
953  ExcMessage("Unable to determine the type of vector automatically because the number of DoFs "
954  "is equal to the number of cells. Please specify DataVectorType."));
955 
956  if (vec.size() == triangulation->n_active_cells())
957  actual_type = type_cell_data;
958  else
959  actual_type = type_dof_data;
960  }
961 
962  switch (actual_type)
963  {
964  case type_cell_data:
965  Assert (vec.size() == triangulation->n_active_cells(),
966  ExcDimensionMismatch (vec.size(),
967  triangulation->n_active_cells()));
968  Assert (names.size() == 1,
970  break;
971 
972  case type_dof_data:
973  Assert (dofs != 0,
975  Assert (vec.size() == dofs->n_dofs(),
977  dofs->n_dofs(),
978  triangulation->n_active_cells()));
979  Assert (names.size() == dofs->get_fe().n_components(),
981  dofs->get_fe().n_components()));
982  break;
983 
984  case type_automatic:
985  // this case should have been handled above...
986  Assert (false, ExcInternalError());
987  }
988 
991  data_component_interpretation);
992  if (actual_type == type_dof_data)
993  dof_data.push_back (std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> >(new_entry));
994  else
995  cell_data.push_back (std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> >(new_entry));
996 }
997 
998 
999 
1000 template <typename DoFHandlerType,
1001  int patch_dim, int patch_space_dim>
1002 template <typename VectorType>
1003 void
1005 add_data_vector (const VectorType &vec,
1006  const DataPostprocessor<DoFHandlerType::space_dimension> &data_postprocessor)
1007 {
1008  // this is a specialized version of the other function where we have a
1009  // postprocessor. if we do, we know that we have type_dof_data, which makes
1010  // things a bit simpler, we also don't need to deal with some of the other
1011  // stuff and use a different constructor of DataEntry
1012 
1013  Assert (dofs != 0,
1015 
1016  Assert (vec.size() == dofs->n_dofs(),
1018  dofs->n_dofs(),
1019  dofs->get_triangulation().n_active_cells()));
1020 
1022  = new internal::DataOut::DataEntry<DoFHandlerType,VectorType>(dofs, &vec, &data_postprocessor);
1023  dof_data.push_back (std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> >(new_entry));
1024 }
1025 
1026 
1027 
1028 template <typename DoFHandlerType,
1029  int patch_dim, int patch_space_dim>
1030 template <typename VectorType>
1031 void
1033 add_data_vector (const DoFHandlerType &dof_handler,
1034  const VectorType &vec,
1035  const DataPostprocessor<DoFHandlerType::space_dimension> &data_postprocessor)
1036 {
1037  // this is a specialized version of the other function where we have a
1038  // postprocessor. if we do, we know that we have type_dof_data, which makes
1039  // things a bit simpler, we also don't need to deal with some of the other
1040  // stuff and use a different constructor of DataEntry
1041 
1042  AssertDimension (vec.size(), dof_handler.n_dofs());
1043 
1045  = new internal::DataOut::DataEntry<DoFHandlerType,VectorType>(&dof_handler, &vec, &data_postprocessor);
1046  dof_data.push_back (std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> >(new_entry));
1047 }
1048 
1049 
1050 
1051 template <typename DoFHandlerType,
1052  int patch_dim, int patch_space_dim>
1053 template <typename VectorType>
1054 void
1057 (const DoFHandlerType &dof_handler,
1058  const VectorType &data,
1059  const std::string &name,
1060  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation)
1061 {
1062  const unsigned int n_components = dof_handler.get_fe().n_components ();
1063 
1064  std::vector<std::string> names;
1065  // if only one component: we only need one name
1066  if (n_components == 1)
1067  names.resize (1, name);
1068  else
1069  // otherwise append _i to the given name
1070  {
1071  names.resize (n_components);
1072  for (unsigned int i=0; i<n_components; ++i)
1073  {
1074  std::ostringstream namebuf;
1075  namebuf << '_' << i;
1076  names[i] = name + namebuf.str();
1077  }
1078  }
1079 
1080  add_data_vector (dof_handler, data, names, data_component_interpretation);
1081 }
1082 
1083 
1084 
1085 template <typename DoFHandlerType,
1086  int patch_dim, int patch_space_dim>
1087 template <typename VectorType>
1088 void
1091 (const DoFHandlerType &dof_handler,
1092  const VectorType &data,
1093  const std::vector<std::string> &names,
1094  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation_)
1095 {
1096  // this is an extended version of the other functions where we pass a vector
1097  // together with its DoFHandler. if we do, we know that we have
1098  // type_dof_data, which makes things a bit simpler
1099  if (triangulation == 0)
1100  triangulation = SmartPointer<const Triangulation<DoFHandlerType::dimension,DoFHandlerType::space_dimension> >(&dof_handler.get_triangulation(), typeid(*this).name());
1101 
1102  Assert (&dof_handler.get_triangulation() == triangulation,
1103  ExcMessage("The triangulation attached to the DoFHandler does not "
1104  "match with the one set previously"));
1105 
1106  Assert (data.size() == dof_handler.n_dofs(),
1107  ExcDimensionMismatch (data.size(), dof_handler.n_dofs()));
1108  Assert (names.size() == dof_handler.get_fe().n_components(),
1110  dof_handler.get_fe().n_components()));
1111 
1112  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &
1113  data_component_interpretation
1114  = (data_component_interpretation_.size() != 0
1115  ?
1116  data_component_interpretation_
1117  :
1118  std::vector<DataComponentInterpretation::DataComponentInterpretation>
1120 
1122  = new internal::DataOut::DataEntry<DoFHandlerType,VectorType>(&dof_handler, &data, names,
1123  data_component_interpretation);
1124  dof_data.push_back (std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> >(new_entry));
1125 }
1126 
1127 
1128 
1129 template <typename DoFHandlerType,
1130  int patch_dim, int patch_space_dim>
1132 {
1133  dof_data.erase (dof_data.begin(), dof_data.end());
1134  cell_data.erase (cell_data.begin(), cell_data.end());
1135 
1136  // delete patches
1137  std::vector<Patch> dummy;
1138  patches.swap (dummy);
1139 }
1140 
1141 
1142 
1143 template <typename DoFHandlerType,
1144  int patch_dim, int patch_space_dim>
1145 void
1148 {
1149  for (unsigned int i=0; i<dof_data.size(); ++i)
1150  dof_data[i]->clear ();
1151 
1152  for (unsigned int i=0; i<cell_data.size(); ++i)
1153  cell_data[i]->clear ();
1154 
1155  if (dofs != 0)
1156  dofs = 0;
1157 }
1158 
1159 
1160 
1161 template <typename DoFHandlerType,
1162  int patch_dim, int patch_space_dim>
1163 void
1165 {
1166  dof_data.erase (dof_data.begin(), dof_data.end());
1167  cell_data.erase (cell_data.begin(), cell_data.end());
1168 
1169  if (dofs != 0)
1170  dofs = 0;
1171 
1172  // delete patches
1173  std::vector<Patch> dummy;
1174  patches.swap (dummy);
1175 }
1176 
1177 
1178 
1179 template <typename DoFHandlerType,
1180  int patch_dim, int patch_space_dim>
1181 std::vector<std::string>
1184 {
1185  std::vector<std::string> names;
1186  // collect the names of dof
1187  // and cell data
1188  typedef
1189  typename std::vector<std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> > >::const_iterator
1190  data_iterator;
1191 
1192  for (data_iterator d=dof_data.begin();
1193  d!=dof_data.end(); ++d)
1194  for (unsigned int i=0; i<(*d)->names.size(); ++i)
1195  names.push_back ((*d)->names[i]);
1196  for (data_iterator d=cell_data.begin(); d!=cell_data.end(); ++d)
1197  {
1198  Assert ((*d)->names.size() == 1, ExcInternalError());
1199  names.push_back ((*d)->names[0]);
1200  }
1201 
1202  return names;
1203 }
1204 
1205 
1206 
1207 template <typename DoFHandlerType,
1208  int patch_dim, int patch_space_dim>
1209 std::vector<std_cxx11::tuple<unsigned int, unsigned int, std::string> >
1211 {
1212  std::vector<std_cxx11::tuple<unsigned int, unsigned int, std::string> >
1213  ranges;
1214 
1215  // collect the ranges of dof
1216  // and cell data
1217  typedef
1218  typename std::vector<std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> > >::const_iterator
1219  data_iterator;
1220 
1221  unsigned int output_component = 0;
1222  for (data_iterator d=dof_data.begin();
1223  d!=dof_data.end(); ++d)
1224  for (unsigned int i=0; i<(*d)->n_output_variables;
1225  ++i, ++output_component)
1226  // see what kind of data we have
1227  // here. note that for the purpose of
1228  // the current function all we care
1229  // about is vector data
1230  if ((*d)->data_component_interpretation[i] ==
1232  {
1233  // ensure that there is a
1234  // continuous number of next
1235  // space_dim components that all
1236  // deal with vectors
1237  Assert (i+patch_space_dim <=
1238  (*d)->n_output_variables,
1240  (*d)->names[i]));
1241  for (unsigned int dd=1; dd<patch_space_dim; ++dd)
1242  Assert ((*d)->data_component_interpretation[i+dd]
1243  ==
1246  (*d)->names[i]));
1247 
1248  // all seems alright, so figure out
1249  // whether there is a common name
1250  // to these components. if not,
1251  // leave the name empty and let the
1252  // output format writer decide what
1253  // to do here
1254  std::string name = (*d)->names[i];
1255  for (unsigned int dd=1; dd<patch_space_dim; ++dd)
1256  if (name != (*d)->names[i+dd])
1257  {
1258  name = "";
1259  break;
1260  }
1261 
1262  // finally add a corresponding
1263  // range
1264  std_cxx11::tuple<unsigned int, unsigned int, std::string>
1265  range (output_component,
1266  output_component+patch_space_dim-1,
1267  name);
1268 
1269  ranges.push_back (range);
1270 
1271  // increase the 'component' counter
1272  // by the appropriate amount, same
1273  // for 'i', since we have already
1274  // dealt with all these components
1275  output_component += patch_space_dim-1;
1276  i += patch_space_dim-1;
1277  }
1278 
1279  // note that we do not have to traverse the
1280  // list of cell data here because cell data
1281  // is one value per (logical) cell and
1282  // therefore cannot be a vector
1283 
1284  // as a final check, the 'component'
1285  // counter should be at the total number of
1286  // components added up now
1287 #ifdef DEBUG
1288  unsigned int n_output_components = 0;
1289  for (data_iterator d=dof_data.begin();
1290  d!=dof_data.end(); ++d)
1291  n_output_components += (*d)->n_output_variables;
1292  Assert (output_component == n_output_components,
1293  ExcInternalError());
1294 #endif
1295 
1296  return ranges;
1297 }
1298 
1299 
1300 
1301 template <typename DoFHandlerType,
1302  int patch_dim, int patch_space_dim>
1303 const std::vector< ::DataOutBase::Patch<patch_dim, patch_space_dim> > &
1305 {
1306  return patches;
1307 }
1308 
1309 
1310 
1311 template <typename DoFHandlerType,
1312  int patch_dim, int patch_space_dim>
1313 std::vector<std_cxx11::shared_ptr<::hp::FECollection<DoFHandlerType::dimension,
1314  DoFHandlerType::space_dimension> > >
1316 {
1317  const unsigned int dhdim = DoFHandlerType::dimension;
1318  const unsigned int dhspacedim = DoFHandlerType::space_dimension;
1319  std::vector<std_cxx11::shared_ptr<::hp::FECollection<dhdim,dhspacedim> > >
1320  finite_elements(this->dof_data.size());
1321  for (unsigned int i=0; i<this->dof_data.size(); ++i)
1322  {
1323  Assert (dof_data[i]->dof_handler != 0,
1325 
1326  // avoid creating too many finite elements and doing a lot of work on
1327  // initializing FEValues downstream: if two DoFHandlers are the same
1328  // (checked by pointer comparison), we can re-use the shared_ptr object
1329  // for the second one. We cannot check for finite element equalities
1330  // because we need different FEValues objects for different dof
1331  // handlers.
1332  bool duplicate = false;
1333  for (unsigned int j=0; j<i; ++j)
1334  if (dof_data[i]->dof_handler == dof_data[j]->dof_handler)
1335  {
1336  finite_elements[i] = finite_elements[j];
1337  duplicate = true;
1338  }
1339  if (duplicate == false)
1340  finite_elements[i].reset(new ::hp::FECollection<dhdim,dhspacedim>
1341  (this->dof_data[i]->dof_handler->get_fe()));
1342  }
1343  if (this->dof_data.empty())
1344  {
1345  finite_elements.resize(1);
1346  finite_elements[0].reset(new ::hp::FECollection<dhdim,dhspacedim>
1348  }
1349  return finite_elements;
1350 }
1351 
1352 
1353 
1354 template <typename DoFHandlerType,
1355  int patch_dim, int patch_space_dim>
1356 std::size_t
1358 {
1362 }
1363 
1364 
1365 
1366 // explicit instantiations
1367 #include "data_out_dof_data.inst"
1368 
1369 DEAL_II_NAMESPACE_CLOSE
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
Definition: fe_values.cc:2850
static ::ExceptionBase & ExcInvalidVectorSize(int arg1, int arg2, int arg3)
void clear_input_data_references()
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
static ::ExceptionBase & ExcInvalidVectorDeclaration(int arg1, std::string arg2)
void attach_triangulation(const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &)
virtual void clear()
const std::vector< DataComponentInterpretation::DataComponentInterpretation > data_component_interpretation
Definition: fe_dgq.h:105
#define AssertIndexRange(index, range)
Definition: exceptions.h:1170
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
Definition: fe_values.cc:2689
void attach_dof_handler(const DoFHandlerType &)
static ::ExceptionBase & ExcOldDataStillPresent()
virtual std::vector< std_cxx11::tuple< unsigned int, unsigned int, std::string > > get_vector_data_ranges() const
DataEntryBase(const DoFHandlerType *dofs, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation)
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
Definition: fe_values.cc:2973
virtual std::size_t memory_consumption() const
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
virtual ~DataOut_DoFData()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: dof_tools.h:46
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
DataEntry(const DoFHandlerType *dofs, const VectorType *data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation)
virtual std::vector< std::string > get_dataset_names() const
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
virtual double get_cell_data_value(const unsigned int cell_number) const
std::size_t memory_consumption() const
virtual void get_function_hessians(const FEValuesBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &fe_patch_values, std::vector< Tensor< 2, DoFHandlerType::space_dimension > > &patch_hessians) const
std::vector< std_cxx11::shared_ptr<::hp::FECollection< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > get_finite_elements() const
Definition: mpi.h:41
virtual void get_function_gradients(const FEValuesBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &fe_patch_values, std::vector< Tensor< 1, DoFHandlerType::space_dimension > > &patch_gradients) const
virtual std::vector< DataComponentInterpretation::DataComponentInterpretation > get_data_component_interpretation() const
bool is_element(const size_type index) const
Definition: index_set.h:1489
static ::ExceptionBase & ExcInvalidCharacter(std::string arg1, size_t arg2)
virtual void get_function_values(const FEValuesBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &fe_patch_values, std::vector< double > &patch_values) const
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcNoTriangulationSelected()
const std::vector< std::string > names
virtual std::vector< std::string > get_names() const =0
static ::ExceptionBase & ExcNoDoFHandlerSelected()
static ::ExceptionBase & ExcInvalidNumberOfNames(int arg1, int arg2)
static ::ExceptionBase & ExcInternalError()
virtual const std::vector< Patch > & get_patches() const