24 template <
int dim,
int spacedim>
34 , cell_quadrature(quadrature)
35 , face_quadrature(face_quadrature)
36 , cell_update_flags(update_flags)
37 , neighbor_cell_update_flags(update_flags)
38 , face_update_flags(face_update_flags)
39 , neighbor_face_update_flags(face_update_flags)
40 , local_dof_indices(fe.n_dofs_per_cell())
41 , neighbor_dof_indices(fe.n_dofs_per_cell())
46 template <
int dim,
int spacedim>
58 , cell_quadrature(quadrature)
59 , face_quadrature(face_quadrature)
60 , cell_update_flags(update_flags)
61 , neighbor_cell_update_flags(neighbor_update_flags)
62 , face_update_flags(face_update_flags)
63 , neighbor_face_update_flags(neighbor_face_update_flags)
64 , local_dof_indices(fe.n_dofs_per_cell())
65 , neighbor_dof_indices(fe.n_dofs_per_cell())
70 template <
int dim,
int spacedim>
88 template <
int dim,
int spacedim>
102 neighbor_update_flags,
105 neighbor_face_update_flags)
110 template <
int dim,
int spacedim>
113 : mapping(scratch.mapping)
115 , cell_quadrature(scratch.cell_quadrature)
116 , face_quadrature(scratch.face_quadrature)
117 , cell_update_flags(scratch.cell_update_flags)
118 , neighbor_cell_update_flags(scratch.neighbor_cell_update_flags)
119 , face_update_flags(scratch.face_update_flags)
120 , neighbor_face_update_flags(scratch.neighbor_face_update_flags)
121 , local_dof_indices(scratch.local_dof_indices)
122 , neighbor_dof_indices(scratch.neighbor_dof_indices)
123 , user_data_storage(scratch.user_data_storage)
124 , internal_data_storage(scratch.internal_data_storage)
129 template <
int dim,
int spacedim>
135 fe_values = std::make_unique<FEValues<dim, spacedim>>(*mapping,
140 fe_values->reinit(cell);
141 local_dof_indices.resize(fe_values->dofs_per_cell);
142 cell->get_dof_indices(local_dof_indices);
143 current_fe_values = fe_values.get();
149 template <
int dim,
int spacedim>
153 const unsigned int face_no)
156 fe_face_values = std::make_unique<FEFaceValues<dim, spacedim>>(
157 *mapping, *fe, face_quadrature, face_update_flags);
159 fe_face_values->reinit(cell, face_no);
160 local_dof_indices.resize(fe->n_dofs_per_cell());
161 cell->get_dof_indices(local_dof_indices);
162 current_fe_values = fe_face_values.get();
163 return *fe_face_values;
168 template <
int dim,
int spacedim>
172 const unsigned int face_no,
173 const unsigned int subface_no)
177 if (!fe_subface_values)
178 fe_subface_values = std::make_unique<FESubfaceValues<dim, spacedim>>(
179 *mapping, *fe, face_quadrature, face_update_flags);
180 fe_subface_values->reinit(cell, face_no, subface_no);
181 local_dof_indices.resize(fe->n_dofs_per_cell());
182 cell->get_dof_indices(local_dof_indices);
184 current_fe_values = fe_subface_values.get();
185 return *fe_subface_values;
188 return reinit(cell, face_no);
193 template <
int dim,
int spacedim>
197 const unsigned int face_no,
198 const unsigned int sub_face_no,
201 const unsigned int face_no_neighbor,
202 const unsigned int sub_face_no_neighbor)
204 if (!interface_fe_values)
205 interface_fe_values = std::make_unique<FEInterfaceValues<dim, spacedim>>(
206 *mapping, *fe, face_quadrature, face_update_flags);
207 interface_fe_values->reinit(cell,
212 sub_face_no_neighbor);
214 current_fe_values = &interface_fe_values->get_fe_face_values(0);
215 current_neighbor_fe_values = &interface_fe_values->get_fe_face_values(1);
217 cell_neighbor->get_dof_indices(neighbor_dof_indices);
218 local_dof_indices = interface_fe_values->get_interface_dof_indices();
219 return *interface_fe_values;
224 template <
int dim,
int spacedim>
229 if (!neighbor_fe_values)
230 neighbor_fe_values = std::make_unique<FEValues<dim, spacedim>>(
231 *mapping, *fe, cell_quadrature, neighbor_cell_update_flags);
233 neighbor_fe_values->reinit(cell);
234 cell->get_dof_indices(neighbor_dof_indices);
235 current_neighbor_fe_values = neighbor_fe_values.get();
236 return *neighbor_fe_values;
241 template <
int dim,
int spacedim>
245 const unsigned int face_no)
247 if (!neighbor_fe_face_values)
248 neighbor_fe_face_values = std::make_unique<FEFaceValues<dim, spacedim>>(
249 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
250 neighbor_fe_face_values->reinit(cell, face_no);
251 cell->get_dof_indices(neighbor_dof_indices);
252 current_neighbor_fe_values = neighbor_fe_face_values.get();
253 return *neighbor_fe_face_values;
258 template <
int dim,
int spacedim>
262 const unsigned int face_no,
263 const unsigned int subface_no)
267 if (!neighbor_fe_subface_values)
268 neighbor_fe_subface_values =
269 std::make_unique<FESubfaceValues<dim, spacedim>>(
270 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
271 neighbor_fe_subface_values->reinit(cell, face_no, subface_no);
272 cell->get_dof_indices(neighbor_dof_indices);
273 current_neighbor_fe_values = neighbor_fe_subface_values.get();
274 return *neighbor_fe_subface_values;
277 return reinit_neighbor(cell, face_no);
282 template <
int dim,
int spacedim>
286 Assert(current_fe_values !=
nullptr,
287 ExcMessage(
"You have to initialize the cache using one of the "
288 "reinit functions first!"));
289 return *current_fe_values;
294 template <
int dim,
int spacedim>
298 Assert(current_neighbor_fe_values !=
nullptr,
299 ExcMessage(
"You have to initialize the cache using one of the "
300 "reinit functions first!"));
301 return *current_neighbor_fe_values;
306 template <
int dim,
int spacedim>
307 const std::vector<Point<spacedim>> &
310 return get_current_fe_values().get_quadrature_points();
315 template <
int dim,
int spacedim>
316 const std::vector<double> &
319 return get_current_fe_values().get_JxW_values();
324 template <
int dim,
int spacedim>
325 const std::vector<double> &
328 return get_current_neighbor_fe_values().get_JxW_values();
333 template <
int dim,
int spacedim>
334 const std::vector<Tensor<1, spacedim>> &
337 return get_current_fe_values().get_normal_vectors();
342 template <
int dim,
int spacedim>
343 const std::vector<Tensor<1, spacedim>> &
346 return get_current_neighbor_fe_values().get_normal_vectors();
351 template <
int dim,
int spacedim>
352 const std::vector<types::global_dof_index> &
355 return local_dof_indices;
360 template <
int dim,
int spacedim>
361 const std::vector<types::global_dof_index> &
364 return neighbor_dof_indices;
369 template <
int dim,
int spacedim>
373 return user_data_storage;
378 template <
int dim,
int spacedim>
382 return user_data_storage;
387 template <
int dim,
int spacedim>
401#include "scratch_data.inst"
const FEValuesBase< dim, spacedim > & get_current_fe_values() const
GeneralDataStorage & get_general_data_storage()
const std::vector< Tensor< 1, spacedim > > & get_neighbor_normal_vectors()
const std::vector< double > & get_JxW_values() const
const Mapping< dim, spacedim > & get_mapping() const
const std::vector< double > & get_neighbor_JxW_values() const
const std::vector< types::global_dof_index > & get_neighbor_dof_indices() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
const FEValues< dim, spacedim > & reinit_neighbor(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
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
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)
const FEValues< dim, spacedim > & reinit(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
static const unsigned int invalid_unsigned_int