16 #include <deal.II/meshworker/scratch_data.h> 18 DEAL_II_NAMESPACE_OPEN
22 template <
int dim,
int spacedim>
32 , cell_quadrature(quadrature)
33 , face_quadrature(face_quadrature)
34 , cell_update_flags(update_flags)
35 , neighbor_cell_update_flags(update_flags)
36 , face_update_flags(face_update_flags)
37 , neighbor_face_update_flags(face_update_flags)
38 , local_dof_indices(fe.dofs_per_cell)
39 , neighbor_dof_indices(fe.dofs_per_cell)
44 template <
int dim,
int spacedim>
56 , cell_quadrature(quadrature)
57 , face_quadrature(face_quadrature)
58 , cell_update_flags(update_flags)
59 , neighbor_cell_update_flags(neighbor_update_flags)
60 , face_update_flags(face_update_flags)
61 , neighbor_face_update_flags(neighbor_face_update_flags)
62 , local_dof_indices(fe.dofs_per_cell)
63 , neighbor_dof_indices(fe.dofs_per_cell)
68 template <
int dim,
int spacedim>
85 template <
int dim,
int spacedim>
98 neighbor_update_flags,
101 neighbor_face_update_flags)
106 template <
int dim,
int spacedim>
109 : mapping(scratch.mapping)
111 , cell_quadrature(scratch.cell_quadrature)
112 , face_quadrature(scratch.face_quadrature)
113 , cell_update_flags(scratch.cell_update_flags)
114 , neighbor_cell_update_flags(scratch.neighbor_cell_update_flags)
115 , face_update_flags(scratch.face_update_flags)
116 , neighbor_face_update_flags(scratch.neighbor_face_update_flags)
117 , local_dof_indices(scratch.local_dof_indices)
118 , neighbor_dof_indices(scratch.neighbor_dof_indices)
119 , user_data_storage(scratch.user_data_storage)
120 , internal_data_storage(scratch.internal_data_storage)
125 template <
int dim,
int spacedim>
131 fe_values = std_cxx14::make_unique<FEValues<dim, spacedim>>(
132 *mapping, *fe, cell_quadrature, cell_update_flags);
134 fe_values->reinit(cell);
135 cell->get_dof_indices(local_dof_indices);
136 current_fe_values = fe_values.get();
142 template <
int dim,
int spacedim>
146 const unsigned int face_no)
149 fe_face_values = std_cxx14::make_unique<FEFaceValues<dim, spacedim>>(
150 *mapping, *fe, face_quadrature, face_update_flags);
152 fe_face_values->reinit(cell, face_no);
153 cell->get_dof_indices(local_dof_indices);
154 current_fe_values = fe_face_values.get();
155 return *fe_face_values;
160 template <
int dim,
int spacedim>
164 const unsigned int face_no,
165 const unsigned int subface_no)
169 if (!fe_subface_values)
171 std_cxx14::make_unique<FESubfaceValues<dim, spacedim>>(
172 *mapping, *fe, face_quadrature, face_update_flags);
173 fe_subface_values->reinit(cell, face_no, subface_no);
174 cell->get_dof_indices(local_dof_indices);
176 current_fe_values = fe_subface_values.get();
177 return *fe_subface_values;
180 return reinit(cell, face_no);
185 template <
int dim,
int spacedim>
190 if (!neighbor_fe_values)
191 neighbor_fe_values = std_cxx14::make_unique<FEValues<dim, spacedim>>(
192 *mapping, *fe, cell_quadrature, neighbor_cell_update_flags);
194 neighbor_fe_values->reinit(cell);
195 cell->get_dof_indices(neighbor_dof_indices);
196 current_neighbor_fe_values = neighbor_fe_values.get();
197 return *neighbor_fe_values;
202 template <
int dim,
int spacedim>
206 const unsigned int face_no)
208 if (!neighbor_fe_face_values)
209 neighbor_fe_face_values =
210 std_cxx14::make_unique<FEFaceValues<dim, spacedim>>(
211 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
212 neighbor_fe_face_values->reinit(cell, face_no);
213 cell->get_dof_indices(neighbor_dof_indices);
214 current_neighbor_fe_values = neighbor_fe_face_values.get();
215 return *neighbor_fe_face_values;
220 template <
int dim,
int spacedim>
224 const unsigned int face_no,
225 const unsigned int subface_no)
229 if (!neighbor_fe_subface_values)
230 neighbor_fe_subface_values =
231 std_cxx14::make_unique<FESubfaceValues<dim, spacedim>>(
232 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
233 neighbor_fe_subface_values->reinit(cell, face_no, subface_no);
234 cell->get_dof_indices(neighbor_dof_indices);
235 current_neighbor_fe_values = neighbor_fe_subface_values.get();
236 return *neighbor_fe_subface_values;
239 return reinit_neighbor(cell, face_no);
244 template <
int dim,
int spacedim>
248 Assert(current_fe_values !=
nullptr,
249 ExcMessage(
"You have to initialize the cache using one of the " 250 "reinit functions first!"));
251 return *current_fe_values;
256 template <
int dim,
int spacedim>
260 Assert(current_neighbor_fe_values !=
nullptr,
261 ExcMessage(
"You have to initialize the cache using one of the " 262 "reinit functions first!"));
263 return *current_neighbor_fe_values;
268 template <
int dim,
int spacedim>
269 const std::vector<Point<spacedim>> &
272 return get_current_fe_values().get_quadrature_points();
277 template <
int dim,
int spacedim>
278 const std::vector<double> &
281 return get_current_fe_values().get_JxW_values();
286 template <
int dim,
int spacedim>
287 const std::vector<double> &
290 return get_current_neighbor_fe_values().get_JxW_values();
295 template <
int dim,
int spacedim>
296 const std::vector<Tensor<1, spacedim>> &
299 return get_current_fe_values().get_normal_vectors();
304 template <
int dim,
int spacedim>
305 const std::vector<Tensor<1, spacedim>> &
308 return get_current_neighbor_fe_values().get_normal_vectors();
313 template <
int dim,
int spacedim>
314 const std::vector<types::global_dof_index> &
317 return local_dof_indices;
322 template <
int dim,
int spacedim>
323 const std::vector<types::global_dof_index> &
326 return neighbor_dof_indices;
331 template <
int dim,
int spacedim>
335 return user_data_storage;
340 template <
int dim,
int spacedim>
348 DEAL_II_NAMESPACE_CLOSE
351 DEAL_II_NAMESPACE_OPEN
354 #include "scratch_data.inst" 356 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
const Mapping< dim, spacedim > & get_mapping() const
const std::vector< types::global_dof_index > & get_local_dof_indices() const
const FEValuesBase< dim, spacedim > & get_current_neighbor_fe_values() const
const FEValues< dim, spacedim > & reinit(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
const std::vector< Tensor< 1, spacedim > > & get_neighbor_normal_vectors()
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)
GeneralDataStorage & get_general_data_storage()
static ::ExceptionBase & ExcMessage(std::string arg1)
const FEValues< dim, spacedim > & reinit_neighbor(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
#define Assert(cond, exc)
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
Abstract base class for mapping classes.
const std::vector< Point< spacedim > > & get_quadrature_points() const
const std::vector< types::global_dof_index > & get_neighbor_dof_indices() const
const FEValuesBase< dim, spacedim > & get_current_fe_values() const
const std::vector< double > & get_neighbor_JxW_values() const
const std::vector< double > & get_JxW_values() const
typename ActiveSelector::active_cell_iterator active_cell_iterator