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,
206 std::vector<double> &patch_values)
const = 0;
264 virtual void clear () = 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;
350 std::vector<std::vector<::Vector<double> > > postprocessed_values;
352 const ::hp::MappingCollection<dim,spacedim> mapping_collection;
353 const std::vector<std_cxx11::shared_ptr<::hp::FECollection<dim,spacedim> > > finite_elements;
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;
506 template <
typename DoFHandlerType,
int patch_dim,
int patch_space_dim=patch_dim>
576 DoFHandlerType::space_dimension> &);
643 template <
class VectorType>
645 const std::vector<std::string> &names,
647 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
648 = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
666 template <
class VectorType>
668 const std::string &name,
670 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
671 = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
687 template <
class VectorType>
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>());
699 template <
class VectorType>
701 const VectorType &data,
702 const std::string &name,
703 const std::vector<DataComponentInterpretation::DataComponentInterpretation> &data_component_interpretation
704 = std::vector<DataComponentInterpretation::DataComponentInterpretation>());
734 template <
class VectorType>
744 template <
class VectorType>
746 const VectorType &data,
792 template <
typename DoFHandlerType2>
802 virtual void clear ();
814 typedef ::DataOutBase::Patch<patch_dim,patch_space_dim>
Patch;
829 std::vector<std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> > >
dof_data;
834 std::vector<std_cxx11::shared_ptr<internal::DataOut::DataEntryBase<DoFHandlerType> > >
cell_data;
861 std::vector<std_cxx11::shared_ptr<::hp::FECollection<DoFHandlerType::dimension,DoFHandlerType::space_dimension> > >
869 std::vector<std_cxx11::tuple<unsigned int, unsigned int, std::string> >
876 template <
class,
int,
int>
885 template <
typename DoFHandlerType,
int patch_dim,
int patch_space_dim>
886 template <
typename DoFHandlerType2>
892 const std::vector<Patch> source_patches = source.
get_patches ();
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."));
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(),
920 ExcMessage (
"Both sources need to declare the same components " 926 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 " 942 const unsigned int old_n_patches =
patches.size();
944 source_patches.begin(),
945 source_patches.end());
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;
954 for (
unsigned int i=old_n_patches; i<
patches.size(); ++i)
955 patches[i].patch_index += old_n_patches;
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)
961 patches[i].neighbors[n] += old_n_patches;
965 DEAL_II_NAMESPACE_CLOSE
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
::DataOutBase::Patch< patch_dim, patch_space_dim > Patch
#define DeclException2(Exception2, type1, type2, outsequence)
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 > &)
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)
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)
#define Assert(cond, exc)
virtual ~DataOut_DoFData()
Abstract base class for mapping classes.
#define DeclExceptionMsg(Exception, defaulttext)
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()
unsigned int n_output_variables
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
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)
static ::ExceptionBase & ExcInvalidCharacter(std::string arg1, size_t arg2)
void clear_data_vectors()
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