Reference documentation for deal.II version Git 7a0e96d111 2021-06-21 21:20:26 -0400
\(\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 
34 #include <deal.II/hp/fe_values.h>
37 
40 
41 #include <memory>
42 
44 
45 namespace 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 
172 namespace 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<
234  &data_component_interpretation);
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
262  get_function_values(const FEValuesBase<dim, spacedim> &fe_patch_values,
263  const ComponentExtractor extract_component,
264  std::vector<double> &patch_values) const = 0;
265 
271  virtual void
272  get_function_values(
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
282  get_function_gradients(
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
293  get_function_gradients(const FEValuesBase<dim, spacedim> &fe_patch_values,
294  const ComponentExtractor extract_component,
295  std::vector<std::vector<Tensor<1, spacedim>>>
296  &patch_gradients_system) const = 0;
297 
302  virtual void
303  get_function_hessians(
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
314  get_function_hessians(const FEValuesBase<dim, spacedim> &fe_patch_values,
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
337  memory_consumption() const = 0;
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>>>
400  & finite_elements,
401  const UpdateFlags update_flags,
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>>>
411  & finite_elements,
412  const UpdateFlags update_flags,
413  const bool use_face_values);
414 
415  ParallelDataBase(const ParallelDataBase &data);
416 
417  void
418  reinit_all_fe_values(
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 
589 template <int dim,
590  int patch_dim,
591  int spacedim = dim,
592  int patch_spacedim = patch_dim>
593 class DataOut_DoFData : public DataOutInterface<patch_dim, patch_spacedim>
594 {
595 public:
601 
602 public:
612  {
617 
622 
626  type_automatic
627  };
628 
632  DataOut_DoFData();
633 
637  virtual ~DataOut_DoFData() override;
638 
647  void
648  attach_dof_handler(const DoFHandler<dim, spacedim> &);
649 
659  void
660  attach_triangulation(const Triangulation<dim, spacedim> &);
661 
725  template <class VectorType>
726  void
727  add_data_vector(
728  const VectorType & data,
729  const std::vector<std::string> &names,
730  const DataVectorType type = type_automatic,
731  const std::vector<DataComponentInterpretation::DataComponentInterpretation>
732  &data_component_interpretation = std::vector<
734 
751  template <class VectorType>
752  void
753  add_data_vector(
754  const VectorType & data,
755  const std::string & name,
756  const DataVectorType type = type_automatic,
757  const std::vector<DataComponentInterpretation::DataComponentInterpretation>
758  &data_component_interpretation = std::vector<
760 
776  template <class VectorType>
777  void
778  add_data_vector(
779  const DoFHandler<dim, spacedim> &dof_handler,
780  const VectorType & data,
781  const std::vector<std::string> & names,
782  const std::vector<DataComponentInterpretation::DataComponentInterpretation>
783  &data_component_interpretation = std::vector<
785 
786 
791  template <class VectorType>
792  void
793  add_data_vector(
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 = std::vector<
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
840  add_data_vector(const DoFHandler<dim, spacedim> & dof_handler,
841  const VectorType & data,
842  const DataPostprocessor<spacedim> &data_postprocessor);
843 
861  template <class VectorType>
862  void
863  add_mg_data_vector(
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
876  add_mg_data_vector(const DoFHandler<dim, spacedim> &dof_handler,
877  const MGLevelObject<VectorType> &data,
878  const std::string & name);
879 
886  void
887  clear_data_vectors();
888 
899  void
900  clear_input_data_references();
901 
925  template <int dim2, int spacedim2>
926  void
927  merge_patches(
930 
935  template <typename DoFHandlerType2>
936  DEAL_II_DEPRECATED void
937  merge_patches(const DataOut_DoFData<DoFHandlerType2::dimension,
938  patch_dim,
939  DoFHandlerType2::space_dimension,
940  patch_spacedim> &source,
942 
949  virtual void
950  clear();
951 
956  std::size_t
957  memory_consumption() const;
958 
959 protected:
964 
969 
974 
978  std::vector<std::shared_ptr<
981 
985  std::vector<std::shared_ptr<
988 
994  std::vector<Patch> patches;
995 
1000 public:
1001  virtual const std::vector<Patch> &
1002  get_patches() const override;
1003 
1008  virtual std::vector<std::string>
1009  get_dataset_names() const override;
1010 
1015  virtual std::vector<
1016  std::tuple<unsigned int,
1017  unsigned int,
1018  std::string,
1020  get_nonscalar_data_ranges() const override;
1021 
1022 protected:
1027  std::vector<std::shared_ptr<::hp::FECollection<dim, spacedim>>>
1028  get_fes() const;
1029 
1030  // Make all template siblings friends. Needed for the merge_patches()
1031  // function.
1032  template <int, int, int, int>
1033  friend class DataOut_DoFData;
1034 
1037  template <int, class>
1038  friend class MGDataOut;
1039 
1040 private:
1044  template <class VectorType>
1045  void
1046  add_data_vector_internal(
1047  const DoFHandler<dim, spacedim> *dof_handler,
1048  const VectorType & data,
1049  const std::vector<std::string> & names,
1050  const DataVectorType type,
1051  const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1052  & data_component_interpretation,
1053  const bool deduce_output_names);
1054 };
1055 
1056 
1057 
1058 // -------------------- template and inline functions ------------------------
1059 template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1060 template <typename VectorType>
1061 void
1063  const VectorType & vec,
1064  const std::string & name,
1065  const DataVectorType type,
1066  const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1067  &data_component_interpretation)
1068 {
1069  Assert(triangulation != nullptr,
1071  std::vector<std::string> names(1, name);
1072  add_data_vector_internal(
1073  dofs, vec, names, type, data_component_interpretation, true);
1074 }
1075 
1076 
1077 
1078 template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1079 template <typename VectorType>
1080 void
1082  const VectorType & vec,
1083  const std::vector<std::string> &names,
1084  const DataVectorType type,
1085  const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1086  &data_component_interpretation)
1087 {
1088  Assert(triangulation != nullptr,
1090  add_data_vector_internal(
1091  dofs, vec, names, type, data_component_interpretation, false);
1092 }
1093 
1094 
1095 
1096 template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1097 template <typename VectorType>
1098 void
1100  const DoFHandler<dim, spacedim> &dof_handler,
1101  const VectorType & data,
1102  const std::string & name,
1103  const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1104  &data_component_interpretation)
1105 {
1106  std::vector<std::string> names(1, name);
1107  add_data_vector_internal(&dof_handler,
1108  data,
1109  names,
1110  type_dof_data,
1111  data_component_interpretation,
1112  true);
1113 }
1114 
1115 
1116 
1117 template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1118 template <typename VectorType>
1119 void
1121  const DoFHandler<dim, spacedim> &dof_handler,
1122  const VectorType & data,
1123  const std::vector<std::string> & names,
1124  const std::vector<DataComponentInterpretation::DataComponentInterpretation>
1125  &data_component_interpretation)
1126 {
1127  add_data_vector_internal(&dof_handler,
1128  data,
1129  names,
1130  type_dof_data,
1131  data_component_interpretation,
1132  false);
1133 }
1134 
1135 
1136 
1137 template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1138 template <typename VectorType>
1139 void
1141  const VectorType & vec,
1142  const DataPostprocessor<spacedim> &data_postprocessor)
1143 {
1144  Assert(dofs != nullptr,
1146  add_data_vector(*dofs, vec, data_postprocessor);
1147 }
1148 
1149 
1150 
1151 template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1152 template <int dim2, int spacedim2>
1153 void
1156  const Point<patch_spacedim> & shift)
1157 {
1158  const std::vector<Patch> &source_patches = source.get_patches();
1159  Assert((patches.size() != 0) && (source_patches.size() != 0),
1160  ExcMessage("When calling this function, both the current "
1161  "object and the one being merged need to have a "
1162  "nonzero number of patches associated with it. "
1163  "Either you called this function on objects that "
1164  "are empty, or you may have forgotten to call "
1165  "the 'build_patches()' function."));
1166  // Check equality of component names
1167  Assert(get_dataset_names() == source.get_dataset_names(),
1169 
1170  // Make sure patches are compatible. Ideally, we would check that all input
1171  // patches from both collections are all compatible, but we'll be content
1172  // with checking that just the first ones from both sources are.
1173  //
1174  // We check compatibility by testing that both sets of patches result
1175  // from the same number of subdivisions, and that they have the same
1176  // number of source vectors (they really should, since we already checked
1177  // that there are the same number of source components above, but you
1178  // never know). This implies that the data should have the same number of
1179  // columns. They should really have the same number of rows as well,
1180  // but depending on whether a patch has points included or not, the
1181  // number of rows may or may not include coordinates for the points,
1182  // and the comparison has to account for that because in each source
1183  // stream, the patches may include some that have points included.
1184  Assert(patches[0].n_subdivisions == source_patches[0].n_subdivisions,
1186  Assert(patches[0].data.n_cols() == source_patches[0].data.n_cols(),
1188  Assert((patches[0].data.n_rows() +
1189  (patches[0].points_are_available ? 0 : patch_spacedim)) ==
1190  (source_patches[0].data.n_rows() +
1191  (source_patches[0].points_are_available ? 0 : patch_spacedim)),
1193 
1194  // check equality of the vector data
1195  // specifications
1196  Assert(get_nonscalar_data_ranges().size() ==
1197  source.get_nonscalar_data_ranges().size(),
1198  ExcMessage("Both sources need to declare the same components "
1199  "as vectors."));
1200  for (unsigned int i = 0; i < get_nonscalar_data_ranges().size(); ++i)
1201  {
1202  Assert(std::get<0>(get_nonscalar_data_ranges()[i]) ==
1203  std::get<0>(source.get_nonscalar_data_ranges()[i]),
1204  ExcMessage("Both sources need to declare the same components "
1205  "as vectors."));
1206  Assert(std::get<1>(get_nonscalar_data_ranges()[i]) ==
1207  std::get<1>(source.get_nonscalar_data_ranges()[i]),
1208  ExcMessage("Both sources need to declare the same components "
1209  "as vectors."));
1210  Assert(std::get<2>(get_nonscalar_data_ranges()[i]) ==
1211  std::get<2>(source.get_nonscalar_data_ranges()[i]),
1212  ExcMessage("Both sources need to declare the same components "
1213  "as vectors."));
1214  }
1215 
1216  // merge patches. store old number
1217  // of elements, since we need to
1218  // adjust patch numbers, etc
1219  // afterwards
1220  const unsigned int old_n_patches = patches.size();
1221  patches.insert(patches.end(), source_patches.begin(), source_patches.end());
1222 
1223  // perform shift, if so desired
1224  if (shift != Point<patch_spacedim>())
1225  for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1226  for (const unsigned int v : GeometryInfo<patch_dim>::vertex_indices())
1227  patches[i].vertices[v] += shift;
1228 
1229  // adjust patch numbers
1230  for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1231  patches[i].patch_index += old_n_patches;
1232 
1233  // adjust patch neighbors
1234  for (unsigned int i = old_n_patches; i < patches.size(); ++i)
1235  for (const unsigned int n : GeometryInfo<patch_dim>::face_indices())
1236  if (patches[i].neighbors[n] != Patch::no_neighbor)
1237  patches[i].neighbors[n] += old_n_patches;
1238 }
1239 
1240 
1241 
1242 template <int dim, int patch_dim, int spacedim, int patch_spacedim>
1243 template <typename DoFHandlerType2>
1244 void
1246  const DataOut_DoFData<DoFHandlerType2::dimension,
1247  patch_dim,
1248  DoFHandlerType2::space_dimension,
1249  patch_spacedim> &source,
1250  const Point<patch_spacedim> & shift)
1251 {
1252  this->merge_patches<DoFHandlerType2::dimension,
1253  DoFHandlerType2::space_dimension>(source, shift);
1254 }
1255 
1256 
1257 
1258 namespace Legacy
1259 {
1264  template <typename DoFHandlerType,
1265  int patch_dim,
1266  int patch_space_dim = patch_dim>
1268  ::DataOut_DoFData<DoFHandlerType::dimension,
1269  patch_dim,
1270  DoFHandlerType::space_dimension,
1271  patch_space_dim>;
1272 } // namespace Legacy
1273 
1274 
1276 
1277 #endif
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< dim, spacedim > > > cell_data
static const unsigned int invalid_unsigned_int
Definition: types.h:196
static ::ExceptionBase & ExcNoTriangulationSelected()
static ::ExceptionBase & ExcIncompatiblePatchLists()
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:538
const std::vector< std::shared_ptr<::hp::FECollection< dim, spacedim > > > finite_elements
::DataOut_DoFData< DoFHandlerType::dimension, patch_dim, DoFHandlerType::space_dimension, patch_space_dim > DataOut_DoFData
std::vector< Patch > patches
static ::ExceptionBase & ExcOldDataStillPresent()
static ::ExceptionBase & ExcInvalidVectorSize(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
Definition: point.h:110
static ::ExceptionBase & ExcIncompatibleDatasetNames()
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual const std::vector< Patch > & get_patches() const override
static ::ExceptionBase & ExcNoDoFHandlerSelected()
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
SmartPointer< const Triangulation< dim, spacedim > > triangulation
#define Assert(cond, exc)
Definition: exceptions.h:1465
UpdateFlags
Abstract base class for mapping classes.
Definition: mapping.h:303
static ::ExceptionBase & ExcInvalidNumberOfNames(int arg1, int arg2)
typename Triangulation< dim, spacedim >::cell_iterator cell_iterator
static ::ExceptionBase & ExcInvalidVectorDeclaration(int arg1, std::string arg2)
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
virtual std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > get_nonscalar_data_ranges() const override
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< dim, spacedim > > > dof_data
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:399
static ::ExceptionBase & ExcInvalidCharacter(std::string arg1, size_t arg2)
Point< 3 > vertices[4]
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 >())
DataPostprocessorInputs::Scalar< spacedim > patch_values_scalar
const ::hp::MappingCollection< dim, spacedim > mapping_collection
static ::ExceptionBase & ExcInvalidTensorDeclaration(int arg1, std::string arg2)
const std::vector< DataComponentInterpretation::DataComponentInterpretation > data_component_interpretation
DataPostprocessorInputs::Vector< spacedim > patch_values_system
SmartPointer< const DoFHandler< dim, spacedim > > dof_handler
SmartPointer< const ::DataPostprocessor< spacedim > > postprocessor
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:398
SmartPointer< const DoFHandler< dim, spacedim > > dofs
void merge_patches(const DataOut_DoFData< dim2, patch_dim, spacedim2, patch_spacedim > &source, const Point< patch_spacedim > &shift=Point< patch_spacedim >())
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:561
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::vector< std::shared_ptr<::hp::FEValues< dim, spacedim > > > x_fe_values
virtual std::vector< std::string > get_dataset_names() const override
#define DEAL_II_DEPRECATED
Definition: config.h:158
std::vector< std::vector<::Vector< double > > > postprocessed_values
std::vector< std::shared_ptr<::hp::FEFaceValues< dim, spacedim > > > x_fe_face_values
void shift(const Tensor< 1, spacedim > &shift_vector, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2022
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)