Reference documentation for deal.II version 8.5.1
data_out_dof_data.h
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 #ifndef dealii__data_out_dof_data_h
17 #define dealii__data_out_dof_data_h
18 
19 
20 
21 #include <deal.II/base/config.h>
22 #include <deal.II/base/smartpointer.h>
23 #include <deal.II/base/data_out_base.h>
24 #include <deal.II/dofs/dof_handler.h>
25 #include <deal.II/grid/tria.h>
26 #include <deal.II/fe/mapping.h>
27 #include <deal.II/hp/q_collection.h>
28 #include <deal.II/hp/fe_collection.h>
29 #include <deal.II/hp/mapping_collection.h>
30 #include <deal.II/hp/fe_values.h>
31 #include <deal.II/numerics/data_postprocessor.h>
32 #include <deal.II/numerics/data_component_interpretation.h>
33 
34 #include <deal.II/base/std_cxx11/shared_ptr.h>
35 
36 DEAL_II_NAMESPACE_OPEN
37 
38 template <int, int> class FEValuesBase;
39 
40 
41 namespace Exceptions
42 {
47  namespace DataOut
48  {
53  int,
54  << "The number of subdivisions per patch, " << arg1
55  << ", is not valid. It needs to be greater or equal to "
56  "one, or zero if you want it to be determined "
57  "automatically.");
58 
63  "For the operation you are attempting, you first need to "
64  "tell the DataOut or related object which DoFHandler or "
65  "triangulation you would like to work on.");
66 
71  "For the operation you are attempting, you first need to "
72  "tell the DataOut or related object which DoFHandler "
73  "you would like to work on.");
74 
79  int, int, int,
80  << "The vector has size " << arg1
81  << " but the DoFHandler object says that there are " << arg2
82  << " degrees of freedom and there are " << arg3
83  << " active cells. The size of your vector needs to be"
84  << " either equal to the number of degrees of freedom, or"
85  << " equal to the number of active cells.");
90  std::string, size_t,
91  << "Please use only the characters [a-zA-Z0-9_<>()] for" << std::endl
92  << "description strings since some graphics formats will only accept these."
93  << std::endl
94  << "The string you gave was <" << arg1
95  << ">, within which the invalid character is <" << arg1[arg2]
96  << ">." << std::endl);
101  "When attaching a triangulation or DoFHandler object, it is "
102  "not allowed if old data vectors are still referenced. If "
103  "you want to reuse an object of the current type, you first "
104  "need to call the 'clear_data_vector()' function.");
109  int, int,
110  << "You have to give one name per component in your "
111  << "data vector. The number you gave was " << arg1
112  << ", but the number of components is " << arg2
113  << ".");
118  "While merging sets of patches, the two sets to be merged "
119  "need to refer to data that agrees on the names of the "
120  "various variables represented. In other words, you "
121  "cannot merge sets of patches that originate from "
122  "entirely unrelated simulations.");
127  "While merging sets of patches, the two sets to be merged "
128  "need to refer to data that agrees on the number of "
129  "subdivisions and other properties. In other words, you "
130  "cannot merge sets of patches that originate from "
131  "entirely unrelated simulations.");
132 
134  int, std::string,
135  << "When declaring that a number of components in a data "
136  << "set to be output logically form a vector instead of "
137  << "simply a set of scalar fields, you need to specify "
138  << "this for all relevant components. Furthermore, "
139  << "vectors must always consist of exactly <dim> "
140  << "components. However, the vector component at "
141  << "position " << arg1 << " with name <" << arg2
142  << "> does not satisfy these conditions.");
143  }
144 }
145 
146 
147 namespace internal
148 {
149  namespace DataOut
150  {
164  template <typename DoFHandlerType>
166  {
167  public:
173  DataEntryBase (const DoFHandlerType *dofs,
174  const std::vector<std::string> &names,
175  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation);
176 
182  DataEntryBase (const DoFHandlerType *dofs,
183  const DataPostprocessor<DoFHandlerType::space_dimension> *data_postprocessor);
184 
188  virtual ~DataEntryBase ();
189 
194  virtual
195  double
196  get_cell_data_value (const unsigned int cell_number) const = 0;
197 
202  virtual
203  void
206  std::vector<double> &patch_values) const = 0;
207 
213  virtual
214  void
217  std::vector<::Vector<double> > &patch_values_system) const = 0;
218 
223  virtual
224  void
227  std::vector<Tensor<1,DoFHandlerType::space_dimension> > &patch_gradients) const = 0;
228 
234  virtual
235  void
238  std::vector<std::vector<Tensor<1,DoFHandlerType::space_dimension> > > &patch_gradients_system) const = 0;
239 
244  virtual
245  void
248  std::vector<Tensor<2,DoFHandlerType::space_dimension> > &patch_hessians) const = 0;
249 
255  virtual
256  void
259  std::vector<std::vector< Tensor<2,DoFHandlerType::space_dimension> > > &patch_hessians_system) const = 0;
260 
264  virtual void clear () = 0;
265 
270  virtual std::size_t memory_consumption () const = 0;
271 
276 
280  const std::vector<std::string> names;
281 
287  const std::vector<DataComponentInterpretation::DataComponentInterpretation>
289 
295 
302  unsigned int n_output_variables;
303  };
304 
305 
322  template <int dim, int spacedim>
324  {
325  ParallelDataBase (const unsigned int n_datasets,
326  const unsigned int n_subdivisions,
327  const std::vector<unsigned int> &n_postprocessor_outputs,
328  const Mapping<dim,spacedim> &mapping,
329  const std::vector<std_cxx11::shared_ptr<::hp::FECollection<dim,spacedim> > > &finite_elements,
330  const UpdateFlags update_flags,
331  const bool use_face_values);
332 
333  ParallelDataBase (const ParallelDataBase &data);
334 
335  template <typename DoFHandlerType>
336  void reinit_all_fe_values(std::vector<std_cxx11::shared_ptr<DataEntryBase<DoFHandlerType> > > &dof_data,
337  const typename ::Triangulation<dim,spacedim>::cell_iterator &cell,
338  const unsigned int face = numbers::invalid_unsigned_int);
339 
341  get_present_fe_values(const unsigned int dataset) const;
342 
343  void resize_system_vectors(const unsigned int n_components);
344 
345  const unsigned int n_datasets;
346  const unsigned int n_subdivisions;
347 
348  DataPostprocessorInputs::Scalar<spacedim> patch_values_scalar;
349  DataPostprocessorInputs::Vector<spacedim> patch_values_system;
350  std::vector<std::vector<::Vector<double> > > postprocessed_values;
351 
352  const ::hp::MappingCollection<dim,spacedim> mapping_collection;
353  const std::vector<std_cxx11::shared_ptr<::hp::FECollection<dim,spacedim> > > finite_elements;
354  const UpdateFlags update_flags;
355 
356  std::vector<std_cxx11::shared_ptr<::hp::FEValues<dim,spacedim> > > x_fe_values;
357  std::vector<std_cxx11::shared_ptr<::hp::FEFaceValues<dim,spacedim> > > x_fe_face_values;
358  };
359  }
360 }
361 
362 
363 //TODO: Most of the documentation of DataOut_DoFData applies to DataOut.
364 
506 template <typename DoFHandlerType, int patch_dim, int patch_space_dim=patch_dim>
507 class DataOut_DoFData : public DataOutInterface<patch_dim,patch_space_dim>
508 {
509 public:
510 
517 
518 public:
519 
529  {
534 
539 
544  };
545 
549  DataOut_DoFData ();
550 
554  virtual ~DataOut_DoFData ();
555 
564  void attach_dof_handler (const DoFHandlerType &);
565 
575  void attach_triangulation (const Triangulation<DoFHandlerType::dimension,
576  DoFHandlerType::space_dimension> &);
577 
643  template <class VectorType>
644  void add_data_vector (const VectorType &data,
645  const std::vector<std::string> &names,
646  const DataVectorType type = type_automatic,
647  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
648  = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
649 
666  template <class VectorType>
667  void add_data_vector (const VectorType &data,
668  const std::string &name,
669  const DataVectorType type = type_automatic,
670  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
671  = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
672 
687  template <class VectorType>
688  void add_data_vector (const DoFHandlerType &dof_handler,
689  const VectorType &data,
690  const std::vector<std::string> &names,
691  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
692  = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
693 
694 
699  template <class VectorType>
700  void add_data_vector (const DoFHandlerType &dof_handler,
701  const VectorType &data,
702  const std::string &name,
703  const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
704  = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
705 
734  template <class VectorType>
735  void add_data_vector (const VectorType &data,
736  const DataPostprocessor<DoFHandlerType::space_dimension> &data_postprocessor);
737 
744  template <class VectorType>
745  void add_data_vector (const DoFHandlerType &dof_handler,
746  const VectorType &data,
747  const DataPostprocessor<DoFHandlerType::space_dimension> &data_postprocessor);
748 
755  void clear_data_vectors ();
756 
768 
792  template <typename DoFHandlerType2>
795 
802  virtual void clear ();
803 
808  std::size_t memory_consumption () const;
809 
810 protected:
814  typedef ::DataOutBase::Patch<patch_dim,patch_space_dim> Patch;
815 
820 
825 
829  std::vector<std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> > > dof_data;
830 
834  std::vector<std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> > > cell_data;
835 
841  std::vector<Patch> patches;
842 
847  virtual
848  const std::vector<Patch> &get_patches () const;
849 
854  virtual
855  std::vector<std::string> get_dataset_names () const;
856 
861  std::vector<std_cxx11::shared_ptr<::hp::FECollection<DoFHandlerType::dimension,DoFHandlerType::space_dimension> > >
862  get_finite_elements() const;
863 
868  virtual
869  std::vector<std_cxx11::tuple<unsigned int, unsigned int, std::string> >
870  get_vector_data_ranges () const;
871 
876  template <class, int, int>
877  friend class DataOut_DoFData;
878 };
879 
880 
881 
882 // -------------------- template and inline functions ------------------------
883 
884 
885 template <typename DoFHandlerType, int patch_dim, int patch_space_dim>
886 template <typename DoFHandlerType2>
887 void
890  const Point<patch_space_dim> &shift)
891 {
892  const std::vector<Patch> source_patches = source.get_patches ();
893  Assert ((patches.size () != 0) &&
894  (source_patches.size () != 0),
895  ExcMessage ("When calling this function, both the current "
896  "object and the one being merged need to have a "
897  "nonzero number of patches associated with it. "
898  "Either you called this function on objects that "
899  "are empty, or you may have forgotten to call "
900  "the 'build_patches()' function."));
901  // check equality of component
902  // names
905  // make sure patches are compatible. we'll
906  // assume that if the first respective
907  // patches are ok that all the other ones
908  // are ok as well
909  Assert (patches[0].n_subdivisions == source_patches[0].n_subdivisions,
911  Assert (patches[0].data.n_rows() == source_patches[0].data.n_rows(),
913  Assert (patches[0].data.n_cols() == source_patches[0].data.n_cols(),
915 
916  // check equality of the vector data
917  // specifications
918  Assert (get_vector_data_ranges().size() ==
919  source.get_vector_data_ranges().size(),
920  ExcMessage ("Both sources need to declare the same components "
921  "as vectors."));
922  for (unsigned int i=0; i<get_vector_data_ranges().size(); ++i)
923  {
924  Assert (std_cxx11::get<0>(get_vector_data_ranges()[i]) ==
925  std_cxx11::get<0>(source.get_vector_data_ranges()[i]),
926  ExcMessage ("Both sources need to declare the same components "
927  "as vectors."));
928  Assert (std_cxx11::get<1>(get_vector_data_ranges()[i]) ==
929  std_cxx11::get<1>(source.get_vector_data_ranges()[i]),
930  ExcMessage ("Both sources need to declare the same components "
931  "as vectors."));
932  Assert (std_cxx11::get<2>(get_vector_data_ranges()[i]) ==
933  std_cxx11::get<2>(source.get_vector_data_ranges()[i]),
934  ExcMessage ("Both sources need to declare the same components "
935  "as vectors."));
936  }
937 
938  // merge patches. store old number
939  // of elements, since we need to
940  // adjust patch numbers, etc
941  // afterwards
942  const unsigned int old_n_patches = patches.size();
943  patches.insert (patches.end(),
944  source_patches.begin(),
945  source_patches.end());
946 
947  // perform shift, if so desired
948  if (shift != Point<patch_space_dim>())
949  for (unsigned int i=old_n_patches; i<patches.size(); ++i)
950  for (unsigned int v=0; v<GeometryInfo<patch_dim>::vertices_per_cell; ++v)
951  patches[i].vertices[v] += shift;
952 
953  // adjust patch numbers
954  for (unsigned int i=old_n_patches; i<patches.size(); ++i)
955  patches[i].patch_index += old_n_patches;
956 
957  // adjust patch neighbors
958  for (unsigned int i=old_n_patches; i<patches.size(); ++i)
959  for (unsigned int n=0; n<GeometryInfo<patch_dim>::faces_per_cell; ++n)
960  if (patches[i].neighbors[n] != Patch::no_neighbor)
961  patches[i].neighbors[n] += old_n_patches;
962 }
963 
964 
965 DEAL_II_NAMESPACE_CLOSE
966 
967 #endif
std::vector< std_cxx11::shared_ptr< internal::DataOut::DataEntryBase< DoFHandlerType > > > cell_data
static ::ExceptionBase & ExcInvalidVectorSize(int arg1, int arg2, int arg3)
void clear_input_data_references()
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
void merge_patches(const DataOut_DoFData< DoFHandlerType2, patch_dim, patch_space_dim > &source, const Point< patch_space_dim > &shift=Point< patch_space_dim >())
static const unsigned int invalid_unsigned_int
Definition: types.h:170
::DataOutBase::Patch< patch_dim, patch_space_dim > Patch
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
static ::ExceptionBase & ExcIncompatiblePatchLists()
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 =0
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
virtual std::size_t memory_consumption() const =0
void attach_dof_handler(const DoFHandlerType &)
static ::ExceptionBase & ExcOldDataStillPresent()
Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension >::cell_iterator cell_iterator
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)
Definition: point.h:89
static ::ExceptionBase & ExcMessage(std::string arg1)
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 =0
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
virtual ~DataOut_DoFData()
Abstract base class for mapping classes.
Definition: dof_tools.h:46
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:553
SmartPointer< const DoFHandlerType > dof_handler
std::vector< Patch > patches
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 >())
static ::ExceptionBase & ExcIncompatibleDatasetNames()
SmartPointer< const ::DataPostprocessor< DoFHandlerType::space_dimension > > postprocessor
virtual void get_function_values(const FEValuesBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &fe_patch_values, std::vector< double > &patch_values) const =0
virtual std::vector< std::string > get_dataset_names() const
std::size_t memory_consumption() const
static const unsigned int no_neighbor
std::vector< std_cxx11::shared_ptr<::hp::FECollection< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > get_finite_elements() const
Definition: mpi.h:41
SmartPointer< const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > triangulation
virtual double get_cell_data_value(const unsigned int cell_number) const =0
std::vector< std_cxx11::shared_ptr< internal::DataOut::DataEntryBase< DoFHandlerType > > > dof_data
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:588
static ::ExceptionBase & ExcInvalidCharacter(std::string arg1, size_t arg2)
static ::ExceptionBase & ExcNoTriangulationSelected()
SmartPointer< const DoFHandlerType > dofs
const std::vector< std::string > names
static ::ExceptionBase & ExcNoDoFHandlerSelected()
static ::ExceptionBase & ExcInvalidNumberOfNames(int arg1, int arg2)
virtual const std::vector< Patch > & get_patches() const