Reference documentation for deal.II version 9.4.1
\(\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\}}\)
Loading...
Searching...
No Matches
data_out_dof_data.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1999 - 2022 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. This is "
70 "generally done using the 'attach_dof_handler()' or "
71 "'attach_triangulation()' member functions.");
72
77 "For the operation you are attempting, you first need to "
78 "tell the DataOut or related object which DoFHandler "
79 "you would like to work on. This is "
80 "generally done using the 'attach_dof_handler()' "
81 "member function.");
82
87 int,
88 int,
89 int,
90 << "The vector has size " << arg1
91 << " but the DoFHandler object says that there are " << arg2
92 << " degrees of freedom and there are " << arg3
93 << " active cells. The size of your vector needs to be"
94 << " either equal to the number of degrees of freedom (when"
95 << " the data is of type type_dof_data), or equal to the"
96 << " number of active cells (when the data is of type "
97 << " type_cell_data).");
103 std::string,
104 size_t,
105 << "Please use only the characters [a-zA-Z0-9_<>()] for" << std::endl
106 << "description strings since some graphics formats will only accept these."
107 << std::endl
108 << "The string you gave was <" << arg1
109 << ">, within which the invalid character is <" << arg1[arg2] << ">."
110 << std::endl);
116 "When attaching a triangulation or DoFHandler object, it is "
117 "not allowed if old data vectors are still referenced. If "
118 "you want to reuse an object of the current type, you first "
119 "need to call the 'clear_data_vector()' function.");
124 int,
125 int,
126 << "You have to give one name per component in your "
127 << "data vector. The number you gave was " << arg1
128 << ", but the number of components is " << arg2 << '.');
133 "While merging sets of patches, the two sets to be merged "
134 "need to refer to data that agrees on the names of the "
135 "various variables represented. In other words, you "
136 "cannot merge sets of patches that originate from "
137 "entirely unrelated simulations.");
142 "While merging sets of patches, the two sets to be merged "
143 "need to refer to data that agrees on the number of "
144 "subdivisions and other properties. In other words, you "
145 "cannot merge sets of patches that originate from "
146 "entirely unrelated simulations.");
147
149 int,
150 std::string,
151 << "When declaring that a number of components in a data "
152 << "set to be output logically form a vector instead of "
153 << "simply a set of scalar fields, you need to specify "
154 << "this for all relevant components. Furthermore, "
155 << "vectors must always consist of exactly <dim> "
156 << "components. However, the vector component at "
157 << "position " << arg1 << " with name <" << arg2
158 << "> does not satisfy these conditions.");
159
161 int,
162 std::string,
163 << "When declaring that a number of components in a data "
164 << "set to be output logically form a tensor instead of "
165 << "simply a set of scalar fields, you need to specify "
166 << "this for all relevant components. Furthermore, "
167 << "tensors must always consist of exactly <dim*dim> "
168 << "components. However, the tensor component at "
169 << "position " << arg1 << " with name <" << arg2
170 << "> does not satisfy these conditions.");
171
172 } // namespace DataOutImplementation
173} // namespace Exceptions
174
175
176namespace internal
177{
178 namespace DataOutImplementation
179 {
204 {
205 real_part,
207 };
208
209
225 template <int dim, int spacedim>
227 {
228 public:
235 const std::vector<std::string> & names,
236 const std::vector<
239
246 const DataPostprocessor<spacedim> *data_postprocessor);
247
251 virtual ~DataEntryBase() = default;
252
257 virtual double
258 get_cell_data_value(const unsigned int cell_number,
259 const ComponentExtractor extract_component) const = 0;
260
265 virtual void
267 const ComponentExtractor extract_component,
268 std::vector<double> &patch_values) const = 0;
269
275 virtual void
277 const FEValuesBase<dim, spacedim> & fe_patch_values,
278 const ComponentExtractor extract_component,
279 std::vector<::Vector<double>> &patch_values_system) const = 0;
280
285 virtual void
287 const FEValuesBase<dim, spacedim> &fe_patch_values,
288 const ComponentExtractor extract_component,
289 std::vector<Tensor<1, spacedim>> & patch_gradients) const = 0;
290
296 virtual void
298 const ComponentExtractor extract_component,
299 std::vector<std::vector<Tensor<1, spacedim>>>
300 &patch_gradients_system) const = 0;
301
306 virtual void
308 const FEValuesBase<dim, spacedim> &fe_patch_values,
309 const ComponentExtractor extract_component,
310 std::vector<Tensor<2, spacedim>> & patch_hessians) const = 0;
311
317 virtual void
319 const ComponentExtractor extract_component,
320 std::vector<std::vector<Tensor<2, spacedim>>>
321 &patch_hessians_system) const = 0;
322
327 virtual bool
328 is_complex_valued() const = 0;
329
333 virtual void
334 clear() = 0;
335
340 virtual std::size_t
342
347
351 const std::vector<std::string> names;
352
358 const std::vector<
361
367
374 unsigned int n_output_variables;
375 };
376
377
394 template <int dim, int spacedim>
396 {
398 const unsigned int n_datasets,
399 const unsigned int n_subdivisions,
400 const std::vector<unsigned int> &n_postprocessor_outputs,
401 const Mapping<dim, spacedim> & mapping,
402 const std::vector<
403 std::shared_ptr<::hp::FECollection<dim, spacedim>>>
406 const bool use_face_values);
407
409 const unsigned int n_datasets,
410 const unsigned int n_subdivisions,
411 const std::vector<unsigned int> &n_postprocessor_outputs,
412 const ::hp::MappingCollection<dim, spacedim> &mapping,
413 const std::vector<
414 std::shared_ptr<::hp::FECollection<dim, spacedim>>>
417 const bool use_face_values);
418
420
421 void
423 std::vector<std::shared_ptr<DataEntryBase<dim, spacedim>>> &dof_data,
424 const typename ::Triangulation<dim, spacedim>::cell_iterator
425 & cell,
426 const unsigned int face = numbers::invalid_unsigned_int);
427
429 get_present_fe_values(const unsigned int dataset) const;
430
431 void
432 resize_system_vectors(const unsigned int n_components);
433
434 const unsigned int n_datasets;
435 const unsigned int n_subdivisions;
436
439 std::vector<std::vector<::Vector<double>>> postprocessed_values;
440
441 const ::hp::MappingCollection<dim, spacedim> mapping_collection;
442 const std::vector<
443 std::shared_ptr<::hp::FECollection<dim, spacedim>>>
446
447 std::vector<std::shared_ptr<::hp::FEValues<dim, spacedim>>>
449 std::vector<std::shared_ptr<::hp::FEFaceValues<dim, spacedim>>>
451 };
452 } // namespace DataOutImplementation
453} // namespace internal
454
455
456// TODO: Most of the documentation of DataOut_DoFData applies to DataOut.
457
593template <int dim,
594 int patch_dim,
595 int spacedim = dim,
596 int patch_spacedim = patch_dim>
597class DataOut_DoFData : public DataOutInterface<patch_dim, patch_spacedim>
598{
599public:
605
606public:
616 {
621
626
631 };
632
637
641 virtual ~DataOut_DoFData() override;
642
651 void
653
663 void
665
729 template <class VectorType>
730 void
732 const VectorType & data,
733 const std::vector<std::string> &names,
734 const DataVectorType type = type_automatic,
735 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
736 &data_component_interpretation = {});
737
754 template <class VectorType>
755 void
757 const VectorType & data,
758 const std::string & name,
759 const DataVectorType type = type_automatic,
760 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
761 &data_component_interpretation = {});
762
778 template <class VectorType>
779 void
781 const DoFHandler<dim, spacedim> &dof_handler,
782 const VectorType & data,
783 const std::vector<std::string> & names,
784 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
785 &data_component_interpretation = {});
786
787
792 template <class VectorType>
793 void
795 const DoFHandler<dim, spacedim> &dof_handler,
796 const VectorType & data,
797 const std::string & name,
798 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
799 &data_component_interpretation = {});
800
827 template <class VectorType>
828 void
829 add_data_vector(const VectorType & data,
830 const DataPostprocessor<spacedim> &data_postprocessor);
831
838 template <class VectorType>
839 void
841 const VectorType & data,
842 const DataPostprocessor<spacedim> &data_postprocessor);
843
861 template <class VectorType>
862 void
864 const DoFHandler<dim, spacedim> &dof_handler,
865 const MGLevelObject<VectorType> &data,
866 const std::vector<std::string> & names,
867 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
868 &data_component_interpretation = std::vector<
870
874 template <class VectorType>
875 void
877 const MGLevelObject<VectorType> &data,
878 const std::string & name);
879
886 void
888
899 void
901
925 template <int dim2, int spacedim2>
926 void
930
935 template <typename DoFHandlerType2>
937 merge_patches(const DataOut_DoFData<DoFHandlerType2::dimension,
938 patch_dim,
939 DoFHandlerType2::space_dimension,
940 patch_spacedim> &source,
942
949 virtual void
951
956 std::size_t
958
959protected:
964
969
974
978 std::vector<std::shared_ptr<
981
985 std::vector<std::shared_ptr<
988
994 std::vector<Patch> patches;
995
1000public:
1001 virtual const std::vector<Patch> &
1002 get_patches() const override;
1003
1004protected:
1009 virtual std::vector<std::string>
1010 get_dataset_names() const override;
1011
1016 virtual std::vector<
1017 std::tuple<unsigned int,
1018 unsigned int,
1019 std::string,
1022
1023protected:
1028 std::vector<std::shared_ptr<::hp::FECollection<dim, spacedim>>>
1029 get_fes() const;
1030
1031 // Make all template siblings friends. Needed for the merge_patches()
1032 // function.
1033 template <int, int, int, int>
1034 friend class DataOut_DoFData;
1035
1038 template <int, class>
1039 friend class MGDataOut;
1040
1041private:
1045 template <class VectorType>
1046 void
1048 const DoFHandler<dim, spacedim> *dof_handler,
1049 const VectorType & data,
1050 const std::vector<std::string> & names,
1051 const DataVectorType type,
1052 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1053 & data_component_interpretation,
1054 const bool deduce_output_names);
1055};
1056
1057
1058
1059// -------------------- template and inline functions ------------------------
1060template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1061template <typename VectorType>
1062void
1064 const VectorType & vec,
1065 const std::string & name,
1066 const DataVectorType type,
1067 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1068 &data_component_interpretation)
1069{
1070 Assert(triangulation != nullptr,
1072 std::vector<std::string> names(1, name);
1073 add_data_vector_internal(
1074 dofs, vec, names, type, data_component_interpretation, true);
1075}
1076
1077
1078
1079template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1080template <typename VectorType>
1081void
1083 const VectorType & vec,
1084 const std::vector<std::string> &names,
1085 const DataVectorType type,
1086 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1087 &data_component_interpretation)
1088{
1089 Assert(triangulation != nullptr,
1091 add_data_vector_internal(
1092 dofs, vec, names, type, data_component_interpretation, false);
1093}
1094
1095
1096
1097template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1098template <typename VectorType>
1099void
1101 const DoFHandler<dim, spacedim> &dof_handler,
1102 const VectorType & data,
1103 const std::string & name,
1104 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1105 &data_component_interpretation)
1106{
1107 std::vector<std::string> names(1, name);
1108 add_data_vector_internal(&dof_handler,
1109 data,
1110 names,
1111 type_dof_data,
1112 data_component_interpretation,
1113 true);
1114}
1115
1116
1117
1118template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1119template <typename VectorType>
1120void
1122 const DoFHandler<dim, spacedim> &dof_handler,
1123 const VectorType & data,
1124 const std::vector<std::string> & names,
1125 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1126 &data_component_interpretation)
1127{
1128 add_data_vector_internal(&dof_handler,
1129 data,
1130 names,
1131 type_dof_data,
1132 data_component_interpretation,
1133 false);
1134}
1135
1136
1137
1138template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1139template <typename VectorType>
1140void
1142 const VectorType & vec,
1143 const DataPostprocessor<spacedim> &data_postprocessor)
1144{
1145 Assert(dofs != nullptr,
1147 add_data_vector(*dofs, vec, data_postprocessor);
1148}
1149
1150
1151
1152template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1153template <int dim2, int spacedim2>
1154void
1157 const Point<patch_spacedim> & shift)
1158{
1159 const std::vector<Patch> &source_patches = source.get_patches();
1160 Assert((patches.size() != 0) && (source_patches.size() != 0),
1161 ExcMessage("When calling this function, both the current "
1162 "object and the one being merged need to have a "
1163 "nonzero number of patches associated with it. "
1164 "Either you called this function on objects that "
1165 "are empty, or you may have forgotten to call "
1166 "the 'build_patches()' function."));
1167 // Check equality of component names
1168 Assert(get_dataset_names() == source.get_dataset_names(),
1170
1171 // Make sure patches are compatible. Ideally, we would check that all input
1172 // patches from both collections are all compatible, but we'll be content
1173 // with checking that just the first ones from both sources are.
1174 //
1175 // We check compatibility by testing that both sets of patches result
1176 // from the same number of subdivisions, and that they have the same
1177 // number of source vectors (they really should, since we already checked
1178 // that there are the same number of source components above, but you
1179 // never know). This implies that the data should have the same number of
1180 // columns. They should really have the same number of rows as well,
1181 // but depending on whether a patch has points included or not, the
1182 // number of rows may or may not include coordinates for the points,
1183 // and the comparison has to account for that because in each source
1184 // stream, the patches may include some that have points included.
1185 Assert(patches[0].n_subdivisions == source_patches[0].n_subdivisions,
1187 Assert(patches[0].data.n_cols() == source_patches[0].data.n_cols(),
1189 Assert((patches[0].data.n_rows() +
1190 (patches[0].points_are_available ? 0 : patch_spacedim)) ==
1191 (source_patches[0].data.n_rows() +
1192 (source_patches[0].points_are_available ? 0 : patch_spacedim)),
1194
1195 // check equality of the vector data
1196 // specifications
1197 Assert(get_nonscalar_data_ranges().size() ==
1198 source.get_nonscalar_data_ranges().size(),
1199 ExcMessage("Both sources need to declare the same components "
1200 "as vectors."));
1201 for (unsigned int i = 0; i < get_nonscalar_data_ranges().size(); ++i)
1202 {
1203 Assert(std::get<0>(get_nonscalar_data_ranges()[i]) ==
1204 std::get<0>(source.get_nonscalar_data_ranges()[i]),
1205 ExcMessage("Both sources need to declare the same components "
1206 "as vectors."));
1207 Assert(std::get<1>(get_nonscalar_data_ranges()[i]) ==
1208 std::get<1>(source.get_nonscalar_data_ranges()[i]),
1209 ExcMessage("Both sources need to declare the same components "
1210 "as vectors."));
1211 Assert(std::get<2>(get_nonscalar_data_ranges()[i]) ==
1212 std::get<2>(source.get_nonscalar_data_ranges()[i]),
1213 ExcMessage("Both sources need to declare the same components "
1214 "as vectors."));
1215 }
1216
1217 // merge patches. store old number
1218 // of elements, since we need to
1219 // adjust patch numbers, etc
1220 // afterwards
1221 const unsigned int old_n_patches = patches.size();
1222 patches.insert(patches.end(), source_patches.begin(), source_patches.end());
1223
1224 // perform shift, if so desired
1225 if (shift != Point<patch_spacedim>())
1226 for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1227 for (const unsigned int v : GeometryInfo<patch_dim>::vertex_indices())
1228 patches[i].vertices[v] += shift;
1229
1230 // adjust patch numbers
1231 for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1232 patches[i].patch_index += old_n_patches;
1233
1234 // adjust patch neighbors
1235 for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1236 for (const unsigned int n : GeometryInfo<patch_dim>::face_indices())
1237 if (patches[i].neighbors[n] != Patch::no_neighbor)
1238 patches[i].neighbors[n] += old_n_patches;
1239}
1240
1241
1242
1243template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1244template <typename DoFHandlerType2>
1245void
1247 const DataOut_DoFData<DoFHandlerType2::dimension,
1248 patch_dim,
1249 DoFHandlerType2::space_dimension,
1250 patch_spacedim> &source,
1251 const Point<patch_spacedim> & shift)
1252{
1253 this->merge_patches<DoFHandlerType2::dimension,
1254 DoFHandlerType2::space_dimension>(source, shift);
1255}
1256
1257
1258
1259namespace Legacy
1260{
1265 template <typename DoFHandlerType,
1266 int patch_dim,
1267 int patch_space_dim = patch_dim>
1269 ::DataOut_DoFData<DoFHandlerType::dimension,
1270 patch_dim,
1271 DoFHandlerType::space_dimension,
1272 patch_space_dim>;
1273} // namespace Legacy
1274
1275
1277
1278#endif
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< std::shared_ptr<::hp::FECollection< dim, spacedim > > > get_fes() const
virtual std::vector< std::string > get_dataset_names() const override
void clear_data_vectors()
void add_data_vector(const DoFHandler< dim, spacedim > &dof_handler, const VectorType &data, const std::string &name, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< dim, spacedim > > > cell_data
std::vector< Patch > patches
void add_data_vector(const DoFHandler< dim, spacedim > &dof_handler, const VectorType &data, const DataPostprocessor< spacedim > &data_postprocessor)
void merge_patches(const DataOut_DoFData< dim2, patch_dim, spacedim2, patch_spacedim > &source, const Point< patch_spacedim > &shift=Point< patch_spacedim >())
void add_data_vector_internal(const DoFHandler< dim, spacedim > *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_data_vector(const VectorType &data, const DataPostprocessor< spacedim > &data_postprocessor)
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< dim, spacedim > > > dof_data
SmartPointer< const DoFHandler< dim, spacedim > > dofs
virtual void clear()
virtual std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > get_nonscalar_data_ranges() const override
virtual const std::vector< Patch > & get_patches() const override
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
void add_data_vector(const DoFHandler< dim, spacedim > &dof_handler, const VectorType &data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
std::size_t memory_consumption() const
void add_mg_data_vector(const DoFHandler< dim, spacedim > &dof_handler, const MGLevelObject< VectorType > &data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
void clear_input_data_references()
typename Triangulation< dim, spacedim >::cell_iterator cell_iterator
void add_mg_data_vector(const DoFHandler< dim, spacedim > &dof_handler, const MGLevelObject< VectorType > &data, const std::string &name)
void merge_patches(const DataOut_DoFData< DoFHandlerType2::dimension, patch_dim, DoFHandlerType2::space_dimension, patch_spacedim > &source, const Point< patch_spacedim > &shift=Point< patch_spacedim >())
void attach_triangulation(const Triangulation< dim, spacedim > &)
SmartPointer< const Triangulation< dim, spacedim > > triangulation
virtual ~DataOut_DoFData() override
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={})
friend class MGDataOut
Abstract base class for mapping classes.
Definition: mapping.h:311
Definition: point.h:111
Definition: tensor.h:503
Definition: vector.h:109
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_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
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:1473
static ::ExceptionBase & ExcNoDoFHandlerSelected()
static ::ExceptionBase & ExcOldDataStillPresent()
static ::ExceptionBase & ExcInvalidTensorDeclaration(int arg1, std::string arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:532
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:487
static ::ExceptionBase & ExcInvalidVectorDeclaration(int arg1, std::string arg2)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:555
static ::ExceptionBase & ExcNoTriangulationSelected()
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:509
static ::ExceptionBase & ExcInvalidVectorSize(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
static const unsigned int invalid_unsigned_int
Definition: types.h:201
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)