Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
data_out_dof_data.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1999 - 2021 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.md at
12// the top level directory of deal.II.
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
26
28
29#include <deal.II/fe/mapping.h>
30
31#include <deal.II/grid/tria.h>
32
37
40
41#include <memory>
42
44
45namespace Exceptions
46{
51 namespace DataOutImplementation
52 {
57 int,
58 << "The number of subdivisions per patch, " << arg1
59 << ", is not valid. It needs to be greater or equal to "
60 "one, or zero if you want it to be determined "
61 "automatically.");
62
67 "For the operation you are attempting, you first need to "
68 "tell the DataOut or related object which DoFHandler or "
69 "triangulation you would like to work on.");
70
75 "For the operation you are attempting, you first need to "
76 "tell the DataOut or related object which DoFHandler "
77 "you would like to work on.");
78
83 int,
84 int,
85 int,
86 << "The vector has size " << arg1
87 << " but the DoFHandler object says that there are " << arg2
88 << " degrees of freedom and there are " << arg3
89 << " active cells. The size of your vector needs to be"
90 << " either equal to the number of degrees of freedom (when"
91 << " the data is of type type_dof_data), or equal to the"
92 << " number of active cells (when the data is of type "
93 << " type_cell_data).");
99 std::string,
100 size_t,
101 << "Please use only the characters [a-zA-Z0-9_<>()] for" << std::endl
102 << "description strings since some graphics formats will only accept these."
103 << std::endl
104 << "The string you gave was <" << arg1
105 << ">, within which the invalid character is <" << arg1[arg2] << ">."
106 << std::endl);
112 "When attaching a triangulation or DoFHandler object, it is "
113 "not allowed if old data vectors are still referenced. If "
114 "you want to reuse an object of the current type, you first "
115 "need to call the 'clear_data_vector()' function.");
120 int,
121 int,
122 << "You have to give one name per component in your "
123 << "data vector. The number you gave was " << arg1
124 << ", but the number of components is " << arg2 << ".");
129 "While merging sets of patches, the two sets to be merged "
130 "need to refer to data that agrees on the names of the "
131 "various variables represented. In other words, you "
132 "cannot merge sets of patches that originate from "
133 "entirely unrelated simulations.");
138 "While merging sets of patches, the two sets to be merged "
139 "need to refer to data that agrees on the number of "
140 "subdivisions and other properties. In other words, you "
141 "cannot merge sets of patches that originate from "
142 "entirely unrelated simulations.");
143
145 int,
146 std::string,
147 << "When declaring that a number of components in a data "
148 << "set to be output logically form a vector instead of "
149 << "simply a set of scalar fields, you need to specify "
150 << "this for all relevant components. Furthermore, "
151 << "vectors must always consist of exactly <dim> "
152 << "components. However, the vector component at "
153 << "position " << arg1 << " with name <" << arg2
154 << "> does not satisfy these conditions.");
155
157 int,
158 std::string,
159 << "When declaring that a number of components in a data "
160 << "set to be output logically form a tensor instead of "
161 << "simply a set of scalar fields, you need to specify "
162 << "this for all relevant components. Furthermore, "
163 << "tensors must always consist of exactly <dim*dim> "
164 << "components. However, the tensor component at "
165 << "position " << arg1 << " with name <" << arg2
166 << "> does not satisfy these conditions.");
167
168 } // namespace DataOutImplementation
169} // namespace Exceptions
170
171
172namespace internal
173{
174 namespace DataOutImplementation
175 {
200 {
201 real_part,
203 };
204
205
221 template <int dim, int spacedim>
223 {
224 public:
231 const std::vector<std::string> & names,
232 const std::vector<
235
242 const DataPostprocessor<spacedim> *data_postprocessor);
243
247 virtual ~DataEntryBase() = default;
248
253 virtual double
254 get_cell_data_value(const unsigned int cell_number,
255 const ComponentExtractor extract_component) const = 0;
256
261 virtual void
263 const ComponentExtractor extract_component,
264 std::vector<double> &patch_values) const = 0;
265
271 virtual void
273 const FEValuesBase<dim, spacedim> & fe_patch_values,
274 const ComponentExtractor extract_component,
275 std::vector<::Vector<double>> &patch_values_system) const = 0;
276
281 virtual void
283 const FEValuesBase<dim, spacedim> &fe_patch_values,
284 const ComponentExtractor extract_component,
285 std::vector<Tensor<1, spacedim>> & patch_gradients) const = 0;
286
292 virtual void
294 const ComponentExtractor extract_component,
295 std::vector<std::vector<Tensor<1, spacedim>>>
296 &patch_gradients_system) const = 0;
297
302 virtual void
304 const FEValuesBase<dim, spacedim> &fe_patch_values,
305 const ComponentExtractor extract_component,
306 std::vector<Tensor<2, spacedim>> & patch_hessians) const = 0;
307
313 virtual void
315 const ComponentExtractor extract_component,
316 std::vector<std::vector<Tensor<2, spacedim>>>
317 &patch_hessians_system) const = 0;
318
323 virtual bool
324 is_complex_valued() const = 0;
325
329 virtual void
330 clear() = 0;
331
336 virtual std::size_t
338
343
347 const std::vector<std::string> names;
348
354 const std::vector<
357
363
370 unsigned int n_output_variables;
371 };
372
373
390 template <int dim, int spacedim>
392 {
394 const unsigned int n_datasets,
395 const unsigned int n_subdivisions,
396 const std::vector<unsigned int> &n_postprocessor_outputs,
397 const Mapping<dim, spacedim> & mapping,
398 const std::vector<
399 std::shared_ptr<::hp::FECollection<dim, spacedim>>>
402 const bool use_face_values);
403
405 const unsigned int n_datasets,
406 const unsigned int n_subdivisions,
407 const std::vector<unsigned int> &n_postprocessor_outputs,
408 const ::hp::MappingCollection<dim, spacedim> &mapping,
409 const std::vector<
410 std::shared_ptr<::hp::FECollection<dim, spacedim>>>
413 const bool use_face_values);
414
416
417 void
419 std::vector<std::shared_ptr<DataEntryBase<dim, spacedim>>> &dof_data,
420 const typename ::Triangulation<dim, spacedim>::cell_iterator
421 & cell,
422 const unsigned int face = numbers::invalid_unsigned_int);
423
425 get_present_fe_values(const unsigned int dataset) const;
426
427 void
428 resize_system_vectors(const unsigned int n_components);
429
430 const unsigned int n_datasets;
431 const unsigned int n_subdivisions;
432
435 std::vector<std::vector<::Vector<double>>> postprocessed_values;
436
437 const ::hp::MappingCollection<dim, spacedim> mapping_collection;
438 const std::vector<
439 std::shared_ptr<::hp::FECollection<dim, spacedim>>>
442
443 std::vector<std::shared_ptr<::hp::FEValues<dim, spacedim>>>
445 std::vector<std::shared_ptr<::hp::FEFaceValues<dim, spacedim>>>
447 };
448 } // namespace DataOutImplementation
449} // namespace internal
450
451
452// TODO: Most of the documentation of DataOut_DoFData applies to DataOut.
453
594template <typename DoFHandlerType,
595 int patch_dim,
596 int patch_space_dim = patch_dim>
597class DataOut_DoFData : public DataOutInterface<patch_dim, patch_space_dim>
598{
599public:
605 typename Triangulation<DoFHandlerType::dimension,
606 DoFHandlerType::space_dimension>::cell_iterator;
607
608public:
618 {
623
628
633 };
634
639
643 virtual ~DataOut_DoFData() override;
644
653 void
654 attach_dof_handler(const DoFHandlerType &);
655
665 void
666 attach_triangulation(const Triangulation<DoFHandlerType::dimension,
667 DoFHandlerType::space_dimension> &);
668
737 template <class VectorType>
738 void
740 const VectorType & data,
741 const std::vector<std::string> &names,
742 const DataVectorType type = type_automatic,
743 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
744 &data_component_interpretation = std::vector<
746
763 template <class VectorType>
764 void
766 const VectorType & data,
767 const std::string & name,
768 const DataVectorType type = type_automatic,
769 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
770 &data_component_interpretation = std::vector<
772
788 template <class VectorType>
789 void
791 const DoFHandlerType & dof_handler,
792 const VectorType & data,
793 const std::vector<std::string> &names,
794 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
795 &data_component_interpretation = std::vector<
797
798
803 template <class VectorType>
804 void
806 const DoFHandlerType &dof_handler,
807 const VectorType & data,
808 const std::string & name,
809 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
810 &data_component_interpretation = std::vector<
812
839 template <class VectorType>
840 void
841 add_data_vector(const VectorType &data,
843 &data_postprocessor);
844
851 template <class VectorType>
852 void
853 add_data_vector(const DoFHandlerType &dof_handler,
854 const VectorType & data,
856 &data_postprocessor);
857
875 template <class VectorType>
876 void
878 const DoFHandlerType & dof_handler,
879 const MGLevelObject<VectorType> &data,
880 const std::vector<std::string> & names,
881 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
882 &data_component_interpretation = std::vector<
884
888 template <class VectorType>
889 void
890 add_mg_data_vector(const DoFHandlerType & dof_handler,
891 const MGLevelObject<VectorType> &data,
892 const std::string & name);
893
900 void
902
913 void
915
939 template <typename DoFHandlerType2>
940 void
944
951 virtual void
953
958 std::size_t
960
961protected:
966
970 SmartPointer<const Triangulation<DoFHandlerType::dimension,
971 DoFHandlerType::space_dimension>>
973
978
982 std::vector<std::shared_ptr<internal::DataOutImplementation::DataEntryBase<
983 DoFHandlerType::dimension,
984 DoFHandlerType::space_dimension>>>
986
990 std::vector<std::shared_ptr<internal::DataOutImplementation::DataEntryBase<
991 DoFHandlerType::dimension,
992 DoFHandlerType::space_dimension>>>
994
1000 std::vector<Patch> patches;
1001
1006 virtual const std::vector<Patch> &
1007 get_patches() const override;
1008
1013 virtual std::vector<std::string>
1014 get_dataset_names() const override;
1015
1020 std::vector<
1021 std::shared_ptr<::hp::FECollection<DoFHandlerType::dimension,
1022 DoFHandlerType::space_dimension>>>
1023 get_fes() const;
1024
1029 virtual std::vector<
1030 std::tuple<unsigned int,
1031 unsigned int,
1032 std::string,
1035
1036 // Make all template siblings friends. Needed for the merge_patches()
1037 // function.
1038 template <class, int, int>
1039 friend class DataOut_DoFData;
1040
1043 template <int, class>
1044 friend class MGDataOut;
1045
1046private:
1050 template <class VectorType>
1051 void
1053 const DoFHandlerType * dof_handler,
1054 const VectorType & data,
1055 const std::vector<std::string> &names,
1056 const DataVectorType type,
1057 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1058 & data_component_interpretation,
1059 const bool deduce_output_names);
1060};
1061
1062
1063
1064// -------------------- template and inline functions ------------------------
1065template <typename DoFHandlerType, int patch_dim, int patch_space_dim>
1066template <typename VectorType>
1067void
1069 const VectorType & vec,
1070 const std::string & name,
1071 const DataVectorType type,
1072 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1073 &data_component_interpretation)
1074{
1075 Assert(triangulation != nullptr,
1077 std::vector<std::string> names(1, name);
1078 add_data_vector_internal(
1079 dofs, vec, names, type, data_component_interpretation, true);
1080}
1081
1082
1083
1084template <typename DoFHandlerType, int patch_dim, int patch_space_dim>
1085template <typename VectorType>
1086void
1088 const VectorType & vec,
1089 const std::vector<std::string> &names,
1090 const DataVectorType type,
1091 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1092 &data_component_interpretation)
1093{
1094 Assert(triangulation != nullptr,
1096 add_data_vector_internal(
1097 dofs, vec, names, type, data_component_interpretation, false);
1098}
1099
1100
1101
1102template <typename DoFHandlerType, int patch_dim, int patch_space_dim>
1103template <typename VectorType>
1104void
1106 const DoFHandlerType &dof_handler,
1107 const VectorType & data,
1108 const std::string & name,
1109 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1110 &data_component_interpretation)
1111{
1112 std::vector<std::string> names(1, name);
1113 add_data_vector_internal(&dof_handler,
1114 data,
1115 names,
1116 type_dof_data,
1117 data_component_interpretation,
1118 true);
1119}
1120
1121
1122
1123template <typename DoFHandlerType, int patch_dim, int patch_space_dim>
1124template <typename VectorType>
1125void
1127 const DoFHandlerType & dof_handler,
1128 const VectorType & data,
1129 const std::vector<std::string> &names,
1130 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1131 &data_component_interpretation)
1132{
1133 add_data_vector_internal(&dof_handler,
1134 data,
1135 names,
1136 type_dof_data,
1137 data_component_interpretation,
1138 false);
1139}
1140
1141
1142
1143template <typename DoFHandlerType, int patch_dim, int patch_space_dim>
1144template <typename VectorType>
1145void
1147 const VectorType & vec,
1149{
1150 Assert(dofs != nullptr,
1152 add_data_vector(*dofs, vec, data_postprocessor);
1153}
1154
1155
1156
1157template <typename DoFHandlerType, int patch_dim, int patch_space_dim>
1158template <typename DoFHandlerType2>
1159void
1163{
1164 const std::vector<Patch> &source_patches = source.get_patches();
1165 Assert((patches.size() != 0) && (source_patches.size() != 0),
1166 ExcMessage("When calling this function, both the current "
1167 "object and the one being merged need to have a "
1168 "nonzero number of patches associated with it. "
1169 "Either you called this function on objects that "
1170 "are empty, or you may have forgotten to call "
1171 "the 'build_patches()' function."));
1172 // Check equality of component names
1173 Assert(get_dataset_names() == source.get_dataset_names(),
1175
1176 // Make sure patches are compatible. Ideally, we would check that all input
1177 // patches from both collections are all compatible, but we'll be content
1178 // with checking that just the first ones from both sources are.
1179 //
1180 // We check compatibility by testing that both sets of patches result
1181 // from the same number of subdivisions, and that they have the same
1182 // number of source vectors (they really should, since we already checked
1183 // that there are the same number of source components above, but you
1184 // never know). This implies that the data should have the same number of
1185 // columns. They should really have the same number of rows as well,
1186 // but depending on whether a patch has points included or not, the
1187 // number of rows may or may not include coordinates for the points,
1188 // and the comparison has to account for that because in each source
1189 // stream, the patches may include some that have points included.
1190 Assert(patches[0].n_subdivisions == source_patches[0].n_subdivisions,
1192 Assert(patches[0].data.n_cols() == source_patches[0].data.n_cols(),
1194 Assert((patches[0].data.n_rows() +
1195 (patches[0].points_are_available ? 0 : patch_space_dim)) ==
1196 (source_patches[0].data.n_rows() +
1197 (source_patches[0].points_are_available ? 0 : patch_space_dim)),
1199
1200 // check equality of the vector data
1201 // specifications
1202 Assert(get_nonscalar_data_ranges().size() ==
1203 source.get_nonscalar_data_ranges().size(),
1204 ExcMessage("Both sources need to declare the same components "
1205 "as vectors."));
1206 for (unsigned int i = 0; i < get_nonscalar_data_ranges().size(); ++i)
1207 {
1208 Assert(std::get<0>(get_nonscalar_data_ranges()[i]) ==
1209 std::get<0>(source.get_nonscalar_data_ranges()[i]),
1210 ExcMessage("Both sources need to declare the same components "
1211 "as vectors."));
1212 Assert(std::get<1>(get_nonscalar_data_ranges()[i]) ==
1213 std::get<1>(source.get_nonscalar_data_ranges()[i]),
1214 ExcMessage("Both sources need to declare the same components "
1215 "as vectors."));
1216 Assert(std::get<2>(get_nonscalar_data_ranges()[i]) ==
1217 std::get<2>(source.get_nonscalar_data_ranges()[i]),
1218 ExcMessage("Both sources need to declare the same components "
1219 "as vectors."));
1220 }
1221
1222 // merge patches. store old number
1223 // of elements, since we need to
1224 // adjust patch numbers, etc
1225 // afterwards
1226 const unsigned int old_n_patches = patches.size();
1227 patches.insert(patches.end(), source_patches.begin(), source_patches.end());
1228
1229 // perform shift, if so desired
1231 for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1232 for (const unsigned int v : GeometryInfo<patch_dim>::vertex_indices())
1233 patches[i].vertices[v] += shift;
1234
1235 // adjust patch numbers
1236 for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1237 patches[i].patch_index += old_n_patches;
1238
1239 // adjust patch neighbors
1240 for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1241 for (const unsigned int n : GeometryInfo<patch_dim>::face_indices())
1242 if (patches[i].neighbors[n] != Patch::no_neighbor)
1243 patches[i].neighbors[n] += old_n_patches;
1244}
1245
1246namespace Legacy
1247{
1254 template <typename DoFHandlerType,
1255 int patch_dim,
1256 int patch_space_dim = patch_dim>
1259} // namespace Legacy
1260
1261
1263
1264#endif
typename Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension >::cell_iterator cell_iterator
void attach_triangulation(const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &)
virtual std::vector< std::string > get_dataset_names() const override
void clear_input_data_references()
void clear_data_vectors()
void add_data_vector(const DoFHandlerType &dof_handler, const VectorType &data, const DataPostprocessor< DoFHandlerType::space_dimension > &data_postprocessor)
virtual void clear()
void add_data_vector(const DoFHandlerType &dof_handler, const VectorType &data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
void merge_patches(const DataOut_DoFData< DoFHandlerType2, patch_dim, patch_space_dim > &source, const Point< patch_space_dim > &shift=Point< patch_space_dim >())
void add_data_vector(const VectorType &data, const DataPostprocessor< DoFHandlerType::space_dimension > &data_postprocessor)
void add_data_vector(const DoFHandlerType &dof_handler, const VectorType &data, const std::string &name, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
void add_data_vector(const VectorType &data, const std::string &name, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
SmartPointer< const DoFHandlerType > dofs
virtual const std::vector< Patch > & get_patches() const override
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > cell_data
void add_mg_data_vector(const DoFHandlerType &dof_handler, const MGLevelObject< VectorType > &data, const std::string &name)
void add_data_vector_internal(const DoFHandlerType *dof_handler, const VectorType &data, const std::vector< std::string > &names, const DataVectorType type, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation, const bool deduce_output_names)
void add_mg_data_vector(const DoFHandlerType &dof_handler, const MGLevelObject< VectorType > &data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
virtual ~DataOut_DoFData() override
SmartPointer< const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > triangulation
void attach_dof_handler(const DoFHandlerType &)
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > dof_data
std::vector< std::shared_ptr<::hp::FECollection< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > get_fes() const
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 >())
std::vector< Patch > patches
virtual std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > get_nonscalar_data_ranges() const override
friend class MGDataOut
std::size_t memory_consumption() const
Definition: point.h:111
Definition: vector.h:110
const std::vector< DataComponentInterpretation::DataComponentInterpretation > data_component_interpretation
virtual void get_function_gradients(const FEValuesBase< dim, spacedim > &fe_patch_values, const ComponentExtractor extract_component, std::vector< std::vector< Tensor< 1, spacedim > > > &patch_gradients_system) const =0
virtual void get_function_gradients(const FEValuesBase< dim, spacedim > &fe_patch_values, const ComponentExtractor extract_component, std::vector< Tensor< 1, spacedim > > &patch_gradients) const =0
virtual void get_function_values(const FEValuesBase< dim, spacedim > &fe_patch_values, const ComponentExtractor extract_component, std::vector< double > &patch_values) const =0
SmartPointer< const DoFHandler< dim, spacedim > > dof_handler
virtual double get_cell_data_value(const unsigned int cell_number, const ComponentExtractor extract_component) const =0
virtual void get_function_hessians(const FEValuesBase< dim, spacedim > &fe_patch_values, const ComponentExtractor extract_component, std::vector< Tensor< 2, spacedim > > &patch_hessians) const =0
virtual std::size_t memory_consumption() const =0
virtual void get_function_hessians(const FEValuesBase< dim, spacedim > &fe_patch_values, const ComponentExtractor extract_component, std::vector< std::vector< Tensor< 2, spacedim > > > &patch_hessians_system) const =0
SmartPointer< const ::DataPostprocessor< spacedim > > postprocessor
DataEntryBase(const DoFHandler< dim, spacedim > *dofs, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation)
virtual void get_function_values(const FEValuesBase< dim, spacedim > &fe_patch_values, const ComponentExtractor extract_component, std::vector<::Vector< double > > &patch_values_system) const =0
DataEntryBase(const DoFHandler< dim, spacedim > *dofs, const DataPostprocessor< spacedim > *data_postprocessor)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
UpdateFlags
static ::ExceptionBase & ExcIncompatiblePatchLists()
static ::ExceptionBase & ExcInvalidCharacter(std::string arg1, size_t arg2)
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
static ::ExceptionBase & ExcInvalidNumberOfNames(int arg1, int arg2)
static ::ExceptionBase & ExcIncompatibleDatasetNames()
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcNoDoFHandlerSelected()
static ::ExceptionBase & ExcOldDataStillPresent()
static ::ExceptionBase & ExcInvalidTensorDeclaration(int arg1, std::string arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:538
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcInvalidVectorDeclaration(int arg1, std::string arg2)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:561
static ::ExceptionBase & ExcNoTriangulationSelected()
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
static ::ExceptionBase & ExcInvalidVectorSize(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
void shift(const Tensor< 1, spacedim > &shift_vector, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2022
static const unsigned int invalid_unsigned_int
Definition: types.h:196
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
const FEValuesBase< dim, spacedim > & get_present_fe_values(const unsigned int dataset) const
std::vector< std::shared_ptr<::hp::FEFaceValues< dim, spacedim > > > x_fe_face_values
ParallelDataBase(const unsigned int n_datasets, const unsigned int n_subdivisions, const std::vector< unsigned int > &n_postprocessor_outputs, const ::hp::MappingCollection< dim, spacedim > &mapping, const std::vector< std::shared_ptr<::hp::FECollection< dim, spacedim > > > &finite_elements, const UpdateFlags update_flags, const bool use_face_values)
DataPostprocessorInputs::Scalar< spacedim > patch_values_scalar
ParallelDataBase(const unsigned int n_datasets, const unsigned int n_subdivisions, const std::vector< unsigned int > &n_postprocessor_outputs, const Mapping< dim, spacedim > &mapping, const std::vector< std::shared_ptr<::hp::FECollection< dim, spacedim > > > &finite_elements, const UpdateFlags update_flags, const bool use_face_values)
const std::vector< std::shared_ptr<::hp::FECollection< dim, spacedim > > > finite_elements
std::vector< std::vector<::Vector< double > > > postprocessed_values
const ::hp::MappingCollection< dim, spacedim > mapping_collection
std::vector< std::shared_ptr<::hp::FEValues< dim, spacedim > > > x_fe_values
DataPostprocessorInputs::Vector< spacedim > patch_values_system
ParallelDataBase(const ParallelDataBase &data)
void reinit_all_fe_values(std::vector< std::shared_ptr< DataEntryBase< dim, spacedim > > > &dof_data, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face=numbers::invalid_unsigned_int)
void resize_system_vectors(const unsigned int n_components)