16 #ifndef dealii__data_out_dof_data_h 17 #define dealii__data_out_dof_data_h 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> 34 #include <deal.II/base/std_cxx11/shared_ptr.h> 36 DEAL_II_NAMESPACE_OPEN
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 " 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.");
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.");
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.");
91 <<
"Please use only the characters [a-zA-Z0-9_<>()] for" << std::endl
92 <<
"description strings since some graphics formats will only accept these." 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.");
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
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.");
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.");
164 template <
typename DoFHandlerType>
174 const std::vector<std::string> &names,
175 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation);
196 get_cell_data_value (
const unsigned int cell_number)
const = 0;
206 std::vector<double> &patch_values)
const = 0;
225 get_function_gradients
236 get_function_gradients
246 get_function_hessians
257 get_function_hessians
264 virtual void clear () = 0;
270 virtual std::size_t memory_consumption ()
const = 0;
280 const std::vector<std::string>
names;
287 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
322 template <
int dim,
int spacedim>
326 const unsigned int n_subdivisions,
327 const std::vector<unsigned int> &n_postprocessor_outputs,
331 const bool use_face_values);
335 template <
typename DoFHandlerType>
337 const typename ::Triangulation<dim,spacedim>::cell_iterator &cell,
341 get_present_fe_values(
const unsigned int dataset)
const;
343 void resize_system_vectors(
const unsigned int n_components);
345 const unsigned int n_datasets;
346 const unsigned int n_subdivisions;
348 std::vector<double> patch_values;
349 std::vector<::Vector<double> > patch_values_system;
350 std::vector<Tensor<1,spacedim> > patch_gradients;
351 std::vector<std::vector<Tensor<1,spacedim> > > patch_gradients_system;
352 std::vector<Tensor<2,spacedim> > patch_hessians;
353 std::vector<std::vector<Tensor<2,spacedim> > > patch_hessians_system;
354 std::vector<std::vector<::Vector<double> > > postprocessed_values;
356 const ::hp::MappingCollection<dim,spacedim> mapping_collection;
357 const std::vector<std_cxx11::shared_ptr<::hp::FECollection<dim,spacedim> > > finite_elements;
360 std::vector<std_cxx11::shared_ptr<::hp::FEValues<dim,spacedim> > > x_fe_values;
361 std::vector<std_cxx11::shared_ptr<::hp::FEFaceValues<dim,spacedim> > > x_fe_face_values;
510 template <
typename DoFHandlerType,
int patch_dim,
int patch_space_dim=patch_dim>
568 void attach_dof_handler (
const DoFHandlerType &);
579 void attach_triangulation (
const Triangulation<DoFHandlerType::dimension,
580 DoFHandlerType::space_dimension> &);
647 template <
class VectorType>
648 void add_data_vector (
const VectorType &data,
649 const std::vector<std::string> &names,
651 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
652 = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
670 template <
class VectorType>
671 void add_data_vector (
const VectorType &data,
672 const std::string &name,
674 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
675 = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
691 template <
class VectorType>
692 void add_data_vector (
const DoFHandlerType &dof_handler,
693 const VectorType &data,
694 const std::vector<std::string> &names,
695 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
696 = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
703 template <
class VectorType>
704 void add_data_vector (
const DoFHandlerType &dof_handler,
705 const VectorType &data,
706 const std::string &name,
707 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
708 = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
738 template <
class VectorType>
739 void add_data_vector (
const VectorType &data,
748 template <
class VectorType>
749 void add_data_vector (
const DoFHandlerType &dof_handler,
750 const VectorType &data,
759 void clear_data_vectors ();
771 void clear_input_data_references ();
796 template <
typename DoFHandlerType2>
806 virtual void clear ();
812 std::size_t memory_consumption ()
const;
818 typedef ::DataOutBase::Patch<patch_dim,patch_space_dim>
Patch;
833 std::vector<std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> > >
dof_data;
838 std::vector<std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> > >
cell_data;
852 const std::vector<Patch> &get_patches ()
const;
859 std::vector<std::string> get_dataset_names ()
const;
865 std::vector<std_cxx11::shared_ptr<::hp::FECollection<DoFHandlerType::dimension,DoFHandlerType::space_dimension> > >
866 get_finite_elements()
const;
873 std::vector<std_cxx11::tuple<unsigned int, unsigned int, std::string> >
874 get_vector_data_ranges ()
const;
880 template <
class,
int,
int>
889 template <
typename DoFHandlerType,
int patch_dim,
int patch_space_dim>
890 template <
typename DoFHandlerType2>
896 const std::vector<Patch> source_patches = source.
get_patches ();
898 (source_patches.size () != 0),
899 ExcMessage (
"When calling this function, both the current " 900 "object and the one being merged need to have a " 901 "nonzero number of patches associated with it. " 902 "Either you called this function on objects that " 903 "are empty, or you may have forgotten to call " 904 "the 'build_patches()' function."));
908 Exceptions::DataOut::ExcIncompatibleDatasetNames());
913 Assert (
patches[0].n_subdivisions == source_patches[0].n_subdivisions,
914 Exceptions::DataOut::ExcIncompatiblePatchLists());
915 Assert (
patches[0].data.n_rows() == source_patches[0].data.n_rows(),
916 Exceptions::DataOut::ExcIncompatiblePatchLists());
917 Assert (
patches[0].data.n_cols() == source_patches[0].data.n_cols(),
918 Exceptions::DataOut::ExcIncompatiblePatchLists());
924 ExcMessage (
"Both sources need to declare the same components " 930 ExcMessage (
"Both sources need to declare the same components " 934 ExcMessage (
"Both sources need to declare the same components " 938 ExcMessage (
"Both sources need to declare the same components " 946 const unsigned int old_n_patches =
patches.size();
948 source_patches.begin(),
949 source_patches.end());
953 for (
unsigned int i=old_n_patches; i<
patches.size(); ++i)
954 for (
unsigned int v=0; v<GeometryInfo<patch_dim>::vertices_per_cell; ++v)
955 patches[i].vertices[v] += shift;
958 for (
unsigned int i=old_n_patches; i<
patches.size(); ++i)
959 patches[i].patch_index += old_n_patches;
962 for (
unsigned int i=old_n_patches; i<
patches.size(); ++i)
963 for (
unsigned int n=0; n<GeometryInfo<patch_dim>::faces_per_cell; ++n)
965 patches[i].neighbors[n] += old_n_patches;
969 DEAL_II_NAMESPACE_CLOSE
std::vector< std_cxx11::shared_ptr< internal::DataOut::DataEntryBase< DoFHandlerType > > > cell_data
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
::DataOutBase::Patch< patch_dim, patch_space_dim > Patch
DeclException1(ExcInvalidNumberOfSubdivisions, int,<< "The number of subdivisions per patch, "<< arg1<< ", is not valid. It needs to be greater or equal to " "one, or zero if you want it to be determined " "automatically.")
const std::vector< DataComponentInterpretation::DataComponentInterpretation > data_component_interpretation
::ExceptionBase & ExcMessage(std::string arg1)
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
DeclExceptionMsg(ExcNoTriangulationSelected, "For the operation you are attempting, you first need to " "tell the DataOut or related object which DoFHandler or " "triangulation you would like to work on.")
#define Assert(cond, exc)
Abstract base class for mapping classes.
SmartPointer< const DoFHandlerType > dof_handler
std::vector< Patch > patches
unsigned int n_output_variables
SmartPointer< const ::DataPostprocessor< DoFHandlerType::space_dimension > > postprocessor
virtual std::vector< std::string > get_dataset_names() const
static const unsigned int no_neighbor
SmartPointer< const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > triangulation
DeclException3(ExcInvalidVectorSize, int, int, int,<< "The vector has size "<< arg1<< " but the DoFHandler object says that there are "<< arg2<< " degrees of freedom and there are "<< arg3<< " active cells. The size of your vector needs to be"<< " either equal to the number of degrees of freedom, or"<< " equal to the number of active cells.")
std::vector< std_cxx11::shared_ptr< internal::DataOut::DataEntryBase< DoFHandlerType > > > dof_data
SmartPointer< const DoFHandlerType > dofs
const std::vector< std::string > names
DeclException2(ExcInvalidCharacter, std::string, size_t,<< "Please use only the characters [a-zA-Z0-9_<>()] for"<< std::endl<< "description strings since some graphics formats will only accept these."<< std::endl<< "The string you gave was <"<< arg1<< ">, within which the invalid character is <"<< arg1[arg2]<< ">."<< std::endl)
virtual const std::vector< Patch > & get_patches() const