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\}}\)
scratch_data.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2019 - 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_meshworker_scratch_data_h
17#define dealii_meshworker_scratch_data_h
18
19#include <deal.II/base/config.h>
20
22
24
28
29#include <boost/any.hpp>
30
31#include <algorithm>
32#include <map>
33#include <vector>
34
36
37namespace MeshWorker
38{
219 template <int dim, int spacedim = dim>
221 {
222 public:
245 const Quadrature<dim> & quadrature,
246 const UpdateFlags & update_flags,
249
269 const Quadrature<dim> & quadrature,
270 const UpdateFlags & update_flags,
271 const UpdateFlags & neighbor_update_flags,
275
290 const Quadrature<dim> & quadrature,
291 const UpdateFlags & update_flags,
294
311 const Quadrature<dim> & quadrature,
312 const UpdateFlags & update_flags,
313 const UpdateFlags & neighbor_update_flags,
317
322 // CurrentCellMethods
327
336 reinit(
338
348 const unsigned int face_no);
349
362 const unsigned int face_no,
363 const unsigned int subface_no);
364
381 const unsigned int face_no,
382 const unsigned int sub_face_no,
384 & cell_neighbor,
385 const unsigned int face_no_neighbor,
386 const unsigned int sub_face_no_neighbor);
387
399 get_current_fe_values() const;
400
404 const std::vector<Point<spacedim>> &
405 get_quadrature_points() const;
406
410 const std::vector<double> &
411 get_JxW_values() const;
412
416 const std::vector<Tensor<1, spacedim>> &
417 get_normal_vectors() const;
418
423 const std::vector<types::global_dof_index> &
424 get_local_dof_indices() const;
425 // CurrentCellMethods
427 // NeighborCellMethods
432
443
454 const unsigned int face_no);
455
469 const unsigned int face_no,
470 const unsigned int subface_no);
471
485
489 const std::vector<double> &
491
495 const std::vector<Tensor<1, spacedim>> &
497
502 const std::vector<types::global_dof_index> &
504 // NeighborCellMethods
506
514
520 const GeneralDataStorage &
522 // CurrentCellEvaluation
527
557 template <typename VectorType, typename Number = double>
558 void
559 extract_local_dof_values(const std::string &global_vector_name,
560 const VectorType & input_vector,
561 const Number dummy = Number(0));
562
571 template <typename Number = double>
572 const std::vector<Number> &
573 get_local_dof_values(const std::string &global_vector_name,
574 Number dummy = Number(0)) const;
575
597 template <typename Extractor, typename Number = double>
598 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
599 template solution_value_type<Number>> &
600 get_values(const std::string &global_vector_name,
601 const Extractor & variable,
602 const Number dummy = Number(0));
603
625 template <typename Extractor, typename Number = double>
626 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
627 template solution_gradient_type<Number>> &
628 get_gradients(const std::string &global_vector_name,
629 const Extractor & variable,
630 const Number dummy = Number(0));
631
654 template <typename Extractor, typename Number = double>
655 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
656 template solution_symmetric_gradient_type<Number>> &
657 get_symmetric_gradients(const std::string &global_vector_name,
658 const Extractor & variable,
659 const Number dummy = Number(0));
660
682 template <typename Extractor, typename Number = double>
683 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
684 template solution_divergence_type<Number>> &
685 get_divergences(const std::string &global_vector_name,
686 const Extractor & variable,
687 const Number dummy = Number(0));
688
710 template <typename Extractor, typename Number = double>
711 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
712 template solution_curl_type<Number>> &
713 get_curls(const std::string &global_vector_name,
714 const Extractor & variable,
715 const Number dummy = Number(0));
716
738 template <typename Extractor, typename Number = double>
739 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
740 template solution_hessian_type<Number>> &
741 get_hessians(const std::string &global_vector_name,
742 const Extractor & variable,
743 const Number dummy = Number(0));
744
766 template <typename Extractor, typename Number = double>
767 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
768 template solution_laplacian_type<Number>> &
769 get_laplacians(const std::string &global_vector_name,
770 const Extractor & variable,
771 const Number dummy = Number(0));
772
794 template <typename Extractor, typename Number = double>
795 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
796 template solution_third_derivative_type<Number>> &
797 get_third_derivatives(const std::string &global_vector_name,
798 const Extractor & variable,
799 const Number dummy = Number(0));
800 // CurrentCellEvaluation
802
807 get_mapping() const;
808
809 private:
814 template <typename Extractor, typename Number = double>
815 std::string
816 get_unique_name(const std::string &global_vector_name,
817 const Extractor & variable,
818 const std::string &object_type,
819 const unsigned int size,
820 const Number & exemplar_number) const;
821
825 template <typename Number = double>
826 std::string
827 get_unique_dofs_name(const std::string &global_vector_name,
828 const unsigned int size,
829 const Number & exemplar_number) const;
830
836
842
848
854
859
864
870
876
880 std::unique_ptr<FEValues<dim, spacedim>> fe_values;
881
885 std::unique_ptr<FEFaceValues<dim, spacedim>> fe_face_values;
886
890 std::unique_ptr<FESubfaceValues<dim, spacedim>> fe_subface_values;
891
895 std::unique_ptr<FEValues<dim, spacedim>> neighbor_fe_values;
896
900 std::unique_ptr<FEFaceValues<dim, spacedim>> neighbor_fe_face_values;
901
905 std::unique_ptr<FESubfaceValues<dim, spacedim>> neighbor_fe_subface_values;
906
910 std::unique_ptr<FEInterfaceValues<dim, spacedim>> interface_fe_values;
911
915 std::vector<types::global_dof_index> local_dof_indices;
916
920 std::vector<types::global_dof_index> neighbor_dof_indices;
921
926
931
937
943 };
944
945#ifndef DOXYGEN
946 template <int dim, int spacedim>
947 template <typename Extractor, typename Number>
948 std::string
950 const std::string &global_vector_name,
951 const Extractor & variable,
952 const std::string &object_type,
953 const unsigned int size,
954 const Number & exemplar_number) const
955 {
956 return global_vector_name + "_" + variable.get_name() + "_" + object_type +
957 "_" + Utilities::int_to_string(size) + "_" +
958 Utilities::type_to_string(exemplar_number);
959 }
960
961
962
963 template <int dim, int spacedim>
964 template <typename Number>
965 std::string
967 const std::string &global_vector_name,
968 const unsigned int size,
969 const Number & exemplar_number) const
970 {
971 return global_vector_name + "_independent_local_dofs_" +
972 Utilities::int_to_string(size) + "_" +
973 Utilities::type_to_string(exemplar_number);
974 }
975
976
977
978 template <int dim, int spacedim>
979 template <typename VectorType, typename Number>
980 void
982 const std::string &global_vector_name,
983 const VectorType & input_vector,
984 const Number dummy)
985 {
986 const unsigned int n_dofs = local_dof_indices.size();
987
988 const std::string name =
989 get_unique_dofs_name(global_vector_name, n_dofs, dummy);
990
991 auto &independent_local_dofs =
992 internal_data_storage
993 .template get_or_add_object_with_name<std::vector<Number>>(name,
994 n_dofs);
995
996 AssertDimension(independent_local_dofs.size(), n_dofs);
997
999 for (unsigned int i = 0; i < n_dofs; ++i)
1001 input_vector(local_dof_indices[i]),
1002 i,
1003 n_dofs,
1004 independent_local_dofs[i]);
1005 else
1006 for (unsigned int i = 0; i < n_dofs; ++i)
1007 independent_local_dofs[i] = input_vector(local_dof_indices[i]);
1008 }
1009
1010
1011
1012 template <int dim, int spacedim>
1013 template <typename Number>
1014 const std::vector<Number> &
1016 const std::string &global_vector_name,
1017 Number dummy) const
1018 {
1019 const unsigned int n_dofs =
1020 get_current_fe_values().get_fe().n_dofs_per_cell();
1021
1022 const std::string dofs_name =
1023 get_unique_dofs_name(global_vector_name, n_dofs, dummy);
1024
1025 Assert(
1026 internal_data_storage.stores_object_with_name(dofs_name),
1027 ExcMessage(
1028 "You did not call yet extract_local_dof_values with the right types!"));
1029
1030 return internal_data_storage
1031 .template get_object_with_name<std::vector<Number>>(dofs_name);
1032 }
1033
1034
1035
1036 template <int dim, int spacedim>
1037 template <typename Extractor, typename Number>
1038 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1039 template solution_value_type<Number>> &
1040 ScratchData<dim, spacedim>::get_values(const std::string &global_vector_name,
1041 const Extractor & variable,
1042 const Number dummy)
1043 {
1044 const std::vector<Number> &independent_local_dofs =
1045 get_local_dof_values(global_vector_name, dummy);
1046
1047 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1048
1049 const unsigned int n_q_points = fev.n_quadrature_points;
1050
1051 const std::string name = get_unique_name(
1052 global_vector_name, variable, "_values_q", n_q_points, dummy);
1053
1054 // Now build the return type
1055 using RetType =
1056 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1057 template solution_value_type<Number>>;
1058
1059 RetType &ret =
1060 internal_data_storage.template get_or_add_object_with_name<RetType>(
1061 name, n_q_points);
1062
1063 AssertDimension(ret.size(), n_q_points);
1064
1065 fev[variable].get_function_values_from_local_dof_values(
1066 independent_local_dofs, ret);
1067 return ret;
1068 }
1069
1070
1071
1072 template <int dim, int spacedim>
1073 template <typename Extractor, typename Number>
1074 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1075 template solution_gradient_type<Number>> &
1077 const std::string &global_vector_name,
1078 const Extractor & variable,
1079 const Number dummy)
1080 {
1081 const std::vector<Number> &independent_local_dofs =
1082 get_local_dof_values(global_vector_name, dummy);
1083
1084 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1085
1086 const unsigned int n_q_points = fev.n_quadrature_points;
1087
1088 const std::string name = get_unique_name(
1089 global_vector_name, variable, "_gradients_q", n_q_points, dummy);
1090
1091 // Now build the return type
1092 using RetType =
1093 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1094 template solution_gradient_type<Number>>;
1095
1096 RetType &ret =
1097 internal_data_storage.template get_or_add_object_with_name<RetType>(
1098 name, n_q_points);
1099
1100 AssertDimension(ret.size(), n_q_points);
1101
1102 fev[variable].get_function_gradients_from_local_dof_values(
1103 independent_local_dofs, ret);
1104 return ret;
1105 }
1106
1107
1108
1109 template <int dim, int spacedim>
1110 template <typename Extractor, typename Number>
1111 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1112 template solution_hessian_type<Number>> &
1114 const std::string &global_vector_name,
1115 const Extractor & variable,
1116 const Number dummy)
1117 {
1118 const std::vector<Number> &independent_local_dofs =
1119 get_local_dof_values(global_vector_name, dummy);
1120
1121 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1122
1123 const unsigned int n_q_points = fev.n_quadrature_points;
1124
1125 const std::string name = get_unique_name(
1126 global_vector_name, variable, "_hessians_q", n_q_points, dummy);
1127
1128 // Now build the return type
1129 using RetType =
1130 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1131 template solution_hessian_type<Number>>;
1132
1133 RetType &ret =
1134 internal_data_storage.template get_or_add_object_with_name<RetType>(
1135 name, n_q_points);
1136
1137
1138 AssertDimension(ret.size(), n_q_points);
1139
1140 fev[variable].get_function_hessians_from_local_dof_values(
1141 independent_local_dofs, ret);
1142 return ret;
1143 }
1144
1145
1146
1147 template <int dim, int spacedim>
1148 template <typename Extractor, typename Number>
1149 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1150 template solution_laplacian_type<Number>> &
1152 const std::string &global_vector_name,
1153 const Extractor & variable,
1154 const Number dummy)
1155 {
1156 const std::vector<Number> &independent_local_dofs =
1157 get_local_dof_values(global_vector_name, dummy);
1158
1159 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1160
1161 const unsigned int n_q_points = fev.n_quadrature_points;
1162
1163 const std::string name = get_unique_name(
1164 global_vector_name, variable, "_laplacians_q", n_q_points, dummy);
1165
1166 // Now build the return type
1167 using RetType =
1168 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1169 template solution_laplacian_type<Number>>;
1170
1171 RetType &ret =
1172 internal_data_storage.template get_or_add_object_with_name<RetType>(
1173 name, n_q_points);
1174
1175
1176 AssertDimension(ret.size(), n_q_points);
1177
1178 fev[variable].get_function_laplacians_from_local_dof_values(
1179 independent_local_dofs, ret);
1180 return ret;
1181 }
1182
1183
1184
1185 template <int dim, int spacedim>
1186 template <typename Extractor, typename Number>
1187 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1188 template solution_third_derivative_type<Number>> &
1190 const std::string &global_vector_name,
1191 const Extractor & variable,
1192 const Number dummy)
1193 {
1194 const std::vector<Number> &independent_local_dofs =
1195 get_local_dof_values(global_vector_name, dummy);
1196
1197 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1198
1199 const unsigned int n_q_points = fev.n_quadrature_points;
1200
1201 const std::string name = get_unique_name(
1202 global_vector_name, variable, "_third_derivatives_q", n_q_points, dummy);
1203
1204 // Now build the return type
1205 using RetType =
1206 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1207 template solution_third_derivative_type<Number>>;
1208
1209 RetType &ret =
1210 internal_data_storage.template get_or_add_object_with_name<RetType>(
1211 name, n_q_points);
1212
1213
1214 AssertDimension(ret.size(), n_q_points);
1215
1216 fev[variable].get_function_third_derivatives_from_local_dof_values(
1217 independent_local_dofs, ret);
1218 return ret;
1219 }
1220
1221
1222
1223 template <int dim, int spacedim>
1224 template <typename Extractor, typename Number>
1225 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1226 template solution_symmetric_gradient_type<Number>> &
1228 const std::string &global_vector_name,
1229 const Extractor & variable,
1230 const Number dummy)
1231 {
1232 const std::vector<Number> &independent_local_dofs =
1233 get_local_dof_values(global_vector_name, dummy);
1234
1235 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1236
1237 const unsigned int n_q_points = fev.n_quadrature_points;
1238
1239 const std::string name = get_unique_name(
1240 global_vector_name, variable, "_symmetric_gradient_q", n_q_points, dummy);
1241
1242
1243 // Now build the return type
1244 using RetType =
1245 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1246 template solution_symmetric_gradient_type<Number>>;
1247
1248 RetType &ret =
1249 internal_data_storage.template get_or_add_object_with_name<RetType>(
1250 name, n_q_points);
1251
1252
1253 AssertDimension(ret.size(), n_q_points);
1254
1255 fev[variable].get_function_symmetric_gradients_from_local_dof_values(
1256 independent_local_dofs, ret);
1257 return ret;
1258 }
1259
1260
1261 template <int dim, int spacedim>
1262 template <typename Extractor, typename Number>
1263 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1264 template solution_divergence_type<Number>> &
1266 const std::string &global_vector_name,
1267 const Extractor & variable,
1268 const Number dummy)
1269 {
1270 const std::vector<Number> &independent_local_dofs =
1271 get_local_dof_values(global_vector_name, dummy);
1272
1273 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1274
1275 const unsigned int n_q_points = fev.n_quadrature_points;
1276
1277 const std::string name = get_unique_name(
1278 global_vector_name, variable, "_divergence_q", n_q_points, dummy);
1279
1280 // Now build the return type
1281 using RetType =
1282 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1283 template solution_divergence_type<Number>>;
1284
1285 RetType &ret =
1286 internal_data_storage.template get_or_add_object_with_name<RetType>(
1287 name, n_q_points);
1288
1289
1290 AssertDimension(ret.size(), n_q_points);
1291
1292 fev[variable].get_function_divergences_from_local_dof_values(
1293 independent_local_dofs, ret);
1294 return ret;
1295 }
1296
1297
1298
1299 template <int dim, int spacedim>
1300 template <typename Extractor, typename Number>
1301 const std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1302 template solution_curl_type<Number>> &
1303 ScratchData<dim, spacedim>::get_curls(const std::string &global_vector_name,
1304 const Extractor & variable,
1305 const Number dummy)
1306 {
1307 const std::vector<Number> &independent_local_dofs =
1308 get_local_dof_values(global_vector_name, dummy);
1309
1310 const FEValuesBase<dim, spacedim> &fev = get_current_fe_values();
1311
1312 const unsigned int n_q_points = fev.n_quadrature_points;
1313
1314 const std::string name = get_unique_name(
1315 global_vector_name, variable, "_curl_q", n_q_points, dummy);
1316
1317 // Now build the return type
1318 using RetType =
1319 std::vector<typename FEValuesViews::View<dim, spacedim, Extractor>::
1320 template solution_curl_type<Number>>;
1321
1322 RetType &ret =
1323 internal_data_storage.template get_or_add_object_with_name<RetType>(
1324 name, n_q_points);
1325
1326 AssertDimension(ret.size(), n_q_points);
1327
1328 fev[variable].get_function_curls_from_local_dof_values(
1329 independent_local_dofs, ret);
1330 return ret;
1331 }
1332
1333#endif
1334
1335} // namespace MeshWorker
1336
1338
1339#endif
const unsigned int n_quadrature_points
Definition: fe_values.h:2432
const FiniteElement< dim, spacedim > & get_fe() const
unsigned int n_dofs_per_cell() const
GeneralDataStorage user_data_storage
Definition: scratch_data.h:925
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_value_type< Number > > & get_values(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const FEValuesBase< dim, spacedim > & get_current_fe_values() const
GeneralDataStorage & get_general_data_storage()
std::string get_unique_name(const std::string &global_vector_name, const Extractor &variable, const std::string &object_type, const unsigned int size, const Number &exemplar_number) const
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_divergence_type< Number > > & get_divergences(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const std::vector< Tensor< 1, spacedim > > & get_neighbor_normal_vectors()
std::unique_ptr< FEValues< dim, spacedim > > neighbor_fe_values
Definition: scratch_data.h:895
std::unique_ptr< FEFaceValues< dim, spacedim > > fe_face_values
Definition: scratch_data.h:885
std::string get_unique_dofs_name(const std::string &global_vector_name, const unsigned int size, const Number &exemplar_number) const
void extract_local_dof_values(const std::string &global_vector_name, const VectorType &input_vector, const Number dummy=Number(0))
std::unique_ptr< FEValues< dim, spacedim > > fe_values
Definition: scratch_data.h:880
UpdateFlags cell_update_flags
Definition: scratch_data.h:858
const std::vector< double > & get_JxW_values() const
const Mapping< dim, spacedim > & get_mapping() const
std::unique_ptr< FEInterfaceValues< dim, spacedim > > interface_fe_values
Definition: scratch_data.h:910
const std::vector< double > & get_neighbor_JxW_values() const
Quadrature< dim > cell_quadrature
Definition: scratch_data.h:847
std::unique_ptr< FEFaceValues< dim, spacedim > > neighbor_fe_face_values
Definition: scratch_data.h:900
std::vector< types::global_dof_index > neighbor_dof_indices
Definition: scratch_data.h:920
SmartPointer< const FiniteElement< dim, spacedim > > fe
Definition: scratch_data.h:841
UpdateFlags face_update_flags
Definition: scratch_data.h:869
const std::vector< types::global_dof_index > & get_neighbor_dof_indices() const
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_laplacian_type< Number > > & get_laplacians(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
UpdateFlags neighbor_cell_update_flags
Definition: scratch_data.h:863
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_third_derivative_type< Number > > & get_third_derivatives(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_curl_type< Number > > & get_curls(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
const FEValues< dim, spacedim > & reinit_neighbor(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_hessian_type< Number > > & get_hessians(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
std::vector< types::global_dof_index > local_dof_indices
Definition: scratch_data.h:915
SmartPointer< const Mapping< dim, spacedim > > mapping
Definition: scratch_data.h:835
SmartPointer< const FEValuesBase< dim, spacedim > > current_fe_values
Definition: scratch_data.h:936
const std::vector< Number > & get_local_dof_values(const std::string &global_vector_name, Number dummy=Number(0)) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const FEValuesBase< dim, spacedim > & get_current_neighbor_fe_values() const
const std::vector< types::global_dof_index > & get_local_dof_indices() const
SmartPointer< const FEValuesBase< dim, spacedim > > current_neighbor_fe_values
Definition: scratch_data.h:942
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_symmetric_gradient_type< Number > > & get_symmetric_gradients(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
Quadrature< dim - 1 > face_quadrature
Definition: scratch_data.h:853
const std::vector< typename FEValuesViews::View< dim, spacedim, Extractor >::template solution_gradient_type< Number > > & get_gradients(const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
UpdateFlags neighbor_face_update_flags
Definition: scratch_data.h:875
ScratchData(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags &update_flags, const Quadrature< dim - 1 > &face_quadrature=Quadrature< dim - 1 >(), const UpdateFlags &face_update_flags=update_default)
Definition: scratch_data.cc:25
std::unique_ptr< FESubfaceValues< dim, spacedim > > neighbor_fe_subface_values
Definition: scratch_data.h:905
const FEValues< dim, spacedim > & reinit(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
std::unique_ptr< FESubfaceValues< dim, spacedim > > fe_subface_values
Definition: scratch_data.h:890
GeneralDataStorage internal_data_storage
Definition: scratch_data.h:930
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
UpdateFlags
@ update_default
No update.
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
static ::ExceptionBase & ExcMessage(std::string arg1)
std::string type_to_string(const T &t)
Definition: utilities.h:1102
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:473