Reference documentation for deal.II version 9.6.0
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1999 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_data_out_dof_data_h
16#define dealii_data_out_dof_data_h
17
18
19
20#include <deal.II/base/config.h>
21
25
27
28#include <deal.II/fe/mapping.h>
29
30#include <deal.II/grid/tria.h>
31
36
39
40#include <memory>
41
43
44namespace Exceptions
45{
50 namespace DataOutImplementation
51 {
56 int,
57 << "The number of subdivisions per patch, " << arg1
58 << ", is not valid. It needs to be greater or equal to "
59 "one, or zero if you want it to be determined "
60 "automatically.");
61
66 "For the operation you are attempting, you first need to "
67 "tell the DataOut or related object which DoFHandler or "
68 "triangulation you would like to work on. This is "
69 "generally done using the 'attach_dof_handler()' or "
70 "'attach_triangulation()' member functions.");
71
76 "For the operation you are attempting, you first need to "
77 "tell the DataOut or related object which DoFHandler "
78 "you would like to work on. This is "
79 "generally done using the 'attach_dof_handler()' "
80 "member function.");
81
86 int,
87 int,
88 int,
89 << "The vector has size " << arg1
90 << " but the DoFHandler object says that there are " << arg2
91 << " degrees of freedom and there are " << arg3
92 << " active cells. The size of your vector needs to be"
93 << " either equal to the number of degrees of freedom (when"
94 << " the data is of type type_dof_data), or equal to the"
95 << " number of active cells (when the data is of type "
96 << " type_cell_data).");
102 std::string,
103 size_t,
104 << "Please use only the characters [a-zA-Z0-9_<>()] for" << std::endl
105 << "description strings since some graphics formats will only accept these."
106 << std::endl
107 << "The string you gave was <" << arg1
108 << ">, within which the invalid character is <" << arg1[arg2] << ">."
109 << std::endl);
115 "When attaching a triangulation or DoFHandler object, it is "
116 "not allowed if old data vectors are still referenced. If "
117 "you want to reuse an object of the current type, you first "
118 "need to call the 'clear_data_vector()' function.");
123 int,
124 int,
125 << "You have to give one name per component in your "
126 << "data vector. The number you gave was " << arg1
127 << ", but the number of components is " << arg2 << '.');
132 "While merging sets of patches, the two sets to be merged "
133 "need to refer to data that agrees on the names of the "
134 "various variables represented. In other words, you "
135 "cannot merge sets of patches that originate from "
136 "entirely unrelated simulations.");
141 "While merging sets of patches, the two sets to be merged "
142 "need to refer to data that agrees on the number of "
143 "subdivisions and other properties. In other words, you "
144 "cannot merge sets of patches that originate from "
145 "entirely unrelated simulations.");
146
148 int,
149 std::string,
150 << "When declaring that a number of components in a data "
151 << "set to be output logically form a vector instead of "
152 << "simply a set of scalar fields, you need to specify "
153 << "this for all relevant components. Furthermore, "
154 << "vectors must always consist of exactly <dim> "
155 << "components. However, the vector component at "
156 << "position " << arg1 << " with name <" << arg2
157 << "> does not satisfy these conditions.");
158
160 int,
161 std::string,
162 << "When declaring that a number of components in a data "
163 << "set to be output logically form a tensor instead of "
164 << "simply a set of scalar fields, you need to specify "
165 << "this for all relevant components. Furthermore, "
166 << "tensors must always consist of exactly <dim*dim> "
167 << "components. However, the tensor component at "
168 << "position " << arg1 << " with name <" << arg2
169 << "> does not satisfy these conditions.");
170
171 } // namespace DataOutImplementation
172} // namespace Exceptions
173
174
175namespace internal
176{
177 namespace DataOutImplementation
178 {
203 {
204 real_part,
206 };
207
208
224 template <int dim, int spacedim>
226 {
227 public:
234 const std::vector<std::string> &names,
235 const std::vector<
238
245 const DataPostprocessor<spacedim> *data_postprocessor);
246
250 virtual ~DataEntryBase() = default;
251
256 virtual double
257 get_cell_data_value(const unsigned int cell_number,
258 const ComponentExtractor extract_component) const = 0;
259
264 virtual void
266 const ComponentExtractor extract_component,
267 std::vector<double> &patch_values) const = 0;
268
274 virtual void
276 const FEValuesBase<dim, spacedim> &fe_patch_values,
277 const ComponentExtractor extract_component,
278 std::vector<::Vector<double>> &patch_values_system) const = 0;
279
284 virtual void
286 const FEValuesBase<dim, spacedim> &fe_patch_values,
287 const ComponentExtractor extract_component,
288 std::vector<Tensor<1, spacedim>> &patch_gradients) const = 0;
289
295 virtual void
297 const ComponentExtractor extract_component,
298 std::vector<std::vector<Tensor<1, spacedim>>>
299 &patch_gradients_system) const = 0;
300
305 virtual void
307 const FEValuesBase<dim, spacedim> &fe_patch_values,
308 const ComponentExtractor extract_component,
309 std::vector<Tensor<2, spacedim>> &patch_hessians) const = 0;
310
316 virtual void
318 const ComponentExtractor extract_component,
319 std::vector<std::vector<Tensor<2, spacedim>>>
320 &patch_hessians_system) const = 0;
321
326 virtual bool
327 is_complex_valued() const = 0;
328
332 virtual void
333 clear() = 0;
334
339 virtual std::size_t
341
346
350 const std::vector<std::string> names;
351
357 const std::vector<
360
366
373 unsigned int n_output_variables;
374 };
375
376
393 template <int dim, int spacedim>
395 {
397 const unsigned int n_datasets,
398 const unsigned int n_subdivisions,
399 const std::vector<unsigned int> &n_postprocessor_outputs,
400 const Mapping<dim, spacedim> &mapping,
401 const std::vector<
402 std::shared_ptr<::hp::FECollection<dim, spacedim>>>
405 const bool use_face_values);
406
408 const unsigned int n_datasets,
409 const unsigned int n_subdivisions,
410 const std::vector<unsigned int> &n_postprocessor_outputs,
411 const ::hp::MappingCollection<dim, spacedim> &mapping,
412 const std::vector<
413 std::shared_ptr<::hp::FECollection<dim, spacedim>>>
416 const bool use_face_values);
417
419
420 void
422 std::vector<std::shared_ptr<DataEntryBase<dim, spacedim>>> &dof_data,
423 const typename ::Triangulation<dim, spacedim>::cell_iterator
424 &cell,
425 const unsigned int face = numbers::invalid_unsigned_int);
426
428 get_present_fe_values(const unsigned int dataset) const;
429
430 void
431 resize_system_vectors(const unsigned int n_components);
432
433 const unsigned int n_datasets;
434 const unsigned int n_subdivisions;
435
438 std::vector<std::vector<::Vector<double>>> postprocessed_values;
439
440 const ::hp::MappingCollection<dim, spacedim> mapping_collection;
441 const std::vector<
442 std::shared_ptr<::hp::FECollection<dim, spacedim>>>
445
446 std::vector<std::shared_ptr<::hp::FEValues<dim, spacedim>>>
448 std::vector<std::shared_ptr<::hp::FEFaceValues<dim, spacedim>>>
450 };
451 } // namespace DataOutImplementation
452} // namespace internal
453
454
455// TODO: Most of the documentation of DataOut_DoFData applies to DataOut.
456
592template <int dim,
593 int patch_dim,
594 int spacedim = dim,
595 int patch_spacedim = patch_dim>
596class DataOut_DoFData : public DataOutInterface<patch_dim, patch_spacedim>
597{
598public:
604
605public:
631
636
640 virtual ~DataOut_DoFData() override;
641
650 void
652
662 void
664
728 template <typename VectorType>
729 void
731 const VectorType &data,
732 const std::vector<std::string> &names,
733 const DataVectorType type = type_automatic,
734 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
735 &data_component_interpretation = {});
736
753 template <typename VectorType>
754 void
756 const VectorType &data,
757 const std::string &name,
758 const DataVectorType type = type_automatic,
759 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
760 &data_component_interpretation = {});
761
777 template <typename VectorType>
778 void
780 const DoFHandler<dim, spacedim> &dof_handler,
781 const VectorType &data,
782 const std::vector<std::string> &names,
783 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
784 &data_component_interpretation = {});
785
786
791 template <typename VectorType>
792 void
794 const DoFHandler<dim, spacedim> &dof_handler,
795 const VectorType &data,
796 const std::string &name,
797 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
798 &data_component_interpretation = {});
799
826 template <typename VectorType>
827 void
828 add_data_vector(const VectorType &data,
829 const DataPostprocessor<spacedim> &data_postprocessor);
830
837 template <typename VectorType>
838 void
840 const VectorType &data,
841 const DataPostprocessor<spacedim> &data_postprocessor);
842
860 template <typename VectorType>
861 void
863 const DoFHandler<dim, spacedim> &dof_handler,
864 const MGLevelObject<VectorType> &data,
865 const std::vector<std::string> &names,
866 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
867 &data_component_interpretation = std::vector<
869
873 template <typename VectorType>
874 void
876 const MGLevelObject<VectorType> &data,
877 const std::string &name);
878
885 void
887
898 void
900
924 template <int dim2, int spacedim2>
925 void
929
936 virtual void
938
943 std::size_t
945
946protected:
951
956
961
965 std::vector<std::shared_ptr<
968
972 std::vector<std::shared_ptr<
975
981 std::vector<Patch> patches;
982
987public:
988 virtual const std::vector<Patch> &
989 get_patches() const override;
990
991protected:
996 virtual std::vector<std::string>
997 get_dataset_names() const override;
998
1003 virtual std::vector<
1004 std::tuple<unsigned int,
1005 unsigned int,
1006 std::string,
1009
1010protected:
1015 std::vector<std::shared_ptr<::hp::FECollection<dim, spacedim>>>
1016 get_fes() const;
1017
1018 // Make all template siblings friends. Needed for the merge_patches()
1019 // function.
1020 template <int, int, int, int>
1021 friend class DataOut_DoFData;
1022
1025 template <int, class>
1026 friend class MGDataOut;
1027
1028private:
1032 template <typename VectorType>
1033 void
1035 const DoFHandler<dim, spacedim> *dof_handler,
1036 const VectorType &data,
1037 const std::vector<std::string> &names,
1038 const DataVectorType type,
1039 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1040 &data_component_interpretation,
1041 const bool deduce_output_names);
1042};
1043
1044
1045
1046// -------------------- template and inline functions ------------------------
1047template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1048template <typename VectorType>
1049void
1051 const VectorType &vec,
1052 const std::string &name,
1053 const DataVectorType type,
1054 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1055 &data_component_interpretation)
1056{
1057 Assert(triangulation != nullptr,
1059 std::vector<std::string> names(1, name);
1060 add_data_vector_internal(
1061 dofs, vec, names, type, data_component_interpretation, true);
1062}
1063
1064
1065
1066template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1067template <typename VectorType>
1068void
1070 const VectorType &vec,
1071 const std::vector<std::string> &names,
1072 const DataVectorType type,
1073 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1074 &data_component_interpretation)
1075{
1076 Assert(triangulation != nullptr,
1078 add_data_vector_internal(
1079 dofs, vec, names, type, data_component_interpretation, false);
1080}
1081
1082
1083
1084template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1085template <typename VectorType>
1086void
1088 const DoFHandler<dim, spacedim> &dof_handler,
1089 const VectorType &data,
1090 const std::string &name,
1091 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1092 &data_component_interpretation)
1093{
1094 std::vector<std::string> names(1, name);
1095 add_data_vector_internal(&dof_handler,
1096 data,
1097 names,
1098 type_dof_data,
1099 data_component_interpretation,
1100 true);
1101}
1102
1103
1104
1105template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1106template <typename VectorType>
1107void
1109 const DoFHandler<dim, spacedim> &dof_handler,
1110 const VectorType &data,
1111 const std::vector<std::string> &names,
1112 const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1113 &data_component_interpretation)
1114{
1115 add_data_vector_internal(&dof_handler,
1116 data,
1117 names,
1118 type_dof_data,
1119 data_component_interpretation,
1120 false);
1121}
1122
1123
1124
1125template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1126template <typename VectorType>
1127void
1129 const VectorType &vec,
1130 const DataPostprocessor<spacedim> &data_postprocessor)
1131{
1132 Assert(dofs != nullptr,
1134 add_data_vector(*dofs, vec, data_postprocessor);
1135}
1136
1137
1138
1139template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1140template <int dim2, int spacedim2>
1141void
1144 const Point<patch_spacedim> &shift)
1145{
1146 const std::vector<Patch> &source_patches = source.get_patches();
1147 Assert((patches.size() != 0) && (source_patches.size() != 0),
1148 ExcMessage("When calling this function, both the current "
1149 "object and the one being merged need to have a "
1150 "nonzero number of patches associated with it. "
1151 "Either you called this function on objects that "
1152 "are empty, or you may have forgotten to call "
1153 "the 'build_patches()' function."));
1154 // Check equality of component names
1155 Assert(get_dataset_names() == source.get_dataset_names(),
1157
1158 // Make sure patches are compatible. Ideally, we would check that all input
1159 // patches from both collections are all compatible, but we'll be content
1160 // with checking that just the first ones from both sources are.
1161 //
1162 // We check compatibility by testing that both sets of patches result
1163 // from the same number of subdivisions, and that they have the same
1164 // number of source vectors (they really should, since we already checked
1165 // that there are the same number of source components above, but you
1166 // never know). This implies that the data should have the same number of
1167 // columns. They should really have the same number of rows as well,
1168 // but depending on whether a patch has points included or not, the
1169 // number of rows may or may not include coordinates for the points,
1170 // and the comparison has to account for that because in each source
1171 // stream, the patches may include some that have points included.
1172 Assert(patches[0].n_subdivisions == source_patches[0].n_subdivisions,
1174 Assert(patches[0].data.n_cols() == source_patches[0].data.n_cols(),
1176 Assert((patches[0].data.n_rows() +
1177 (patches[0].points_are_available ? 0 : patch_spacedim)) ==
1178 (source_patches[0].data.n_rows() +
1179 (source_patches[0].points_are_available ? 0 : patch_spacedim)),
1181
1182 // check equality of the vector data
1183 // specifications
1184 Assert(get_nonscalar_data_ranges().size() ==
1185 source.get_nonscalar_data_ranges().size(),
1186 ExcMessage("Both sources need to declare the same components "
1187 "as vectors."));
1188 for (unsigned int i = 0; i < get_nonscalar_data_ranges().size(); ++i)
1189 {
1190 Assert(std::get<0>(get_nonscalar_data_ranges()[i]) ==
1191 std::get<0>(source.get_nonscalar_data_ranges()[i]),
1192 ExcMessage("Both sources need to declare the same components "
1193 "as vectors."));
1194 Assert(std::get<1>(get_nonscalar_data_ranges()[i]) ==
1195 std::get<1>(source.get_nonscalar_data_ranges()[i]),
1196 ExcMessage("Both sources need to declare the same components "
1197 "as vectors."));
1198 Assert(std::get<2>(get_nonscalar_data_ranges()[i]) ==
1199 std::get<2>(source.get_nonscalar_data_ranges()[i]),
1200 ExcMessage("Both sources need to declare the same components "
1201 "as vectors."));
1202 }
1203
1204 // merge patches. store old number
1205 // of elements, since we need to
1206 // adjust patch numbers, etc
1207 // afterwards
1208 const unsigned int old_n_patches = patches.size();
1209 patches.insert(patches.end(), source_patches.begin(), source_patches.end());
1210
1211 // perform shift, if so desired
1212 if (shift != Point<patch_spacedim>())
1213 for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1214 for (const unsigned int v : GeometryInfo<patch_dim>::vertex_indices())
1215 patches[i].vertices[v] += shift;
1216
1217 // adjust patch numbers
1218 for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1219 patches[i].patch_index += old_n_patches;
1220
1221 // adjust patch neighbors
1222 for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1223 for (const unsigned int n : GeometryInfo<patch_dim>::face_indices())
1224 if (patches[i].neighbors[n] != Patch::no_neighbor)
1225 patches[i].neighbors[n] += old_n_patches;
1226}
1227
1228
1230
1231#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 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={})
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
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:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
Point< 3 > vertices[4]
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)
static ::ExceptionBase & ExcNoDoFHandlerSelected()
static ::ExceptionBase & ExcOldDataStillPresent()
static ::ExceptionBase & ExcInvalidTensorDeclaration(int arg1, std::string arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:539
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
static ::ExceptionBase & ExcInvalidVectorDeclaration(int arg1, std::string arg2)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition exceptions.h:562
static ::ExceptionBase & ExcNoTriangulationSelected()
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:516
static ::ExceptionBase & ExcInvalidVectorSize(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
UpdateFlags
static const unsigned int invalid_unsigned_int
Definition types.h:220
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices()
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)