Reference documentation for deal.II version 9.3.2
\(\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.cc
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 
17 
18 #include <memory>
19 
21 
22 namespace MeshWorker
23 {
24  template <int dim, int spacedim>
26  const Mapping<dim, spacedim> & mapping,
28  const Quadrature<dim> & quadrature,
29  const UpdateFlags & update_flags,
30  const Quadrature<dim - 1> & face_quadrature,
31  const UpdateFlags & face_update_flags)
32  : mapping(&mapping)
33  , fe(&fe)
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())
42  {}
43 
44 
45 
46  template <int dim, int spacedim>
48  const Mapping<dim, spacedim> & mapping,
50  const Quadrature<dim> & quadrature,
51  const UpdateFlags & update_flags,
52  const UpdateFlags & neighbor_update_flags,
53  const Quadrature<dim - 1> & face_quadrature,
54  const UpdateFlags & face_update_flags,
55  const UpdateFlags & neighbor_face_update_flags)
56  : mapping(&mapping)
57  , fe(&fe)
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())
66  {}
67 
68 
69 
70  template <int dim, int spacedim>
73  const Quadrature<dim> & quadrature,
74  const UpdateFlags & update_flags,
75  const Quadrature<dim - 1> & face_quadrature,
76  const UpdateFlags & face_update_flags)
78  .template get_default_linear_mapping<dim, spacedim>(),
79  fe,
80  quadrature,
81  update_flags,
82  face_quadrature,
83  face_update_flags)
84  {}
85 
86 
87 
88  template <int dim, int spacedim>
91  const Quadrature<dim> & quadrature,
92  const UpdateFlags & update_flags,
93  const UpdateFlags & neighbor_update_flags,
94  const Quadrature<dim - 1> & face_quadrature,
95  const UpdateFlags & face_update_flags,
96  const UpdateFlags & neighbor_face_update_flags)
98  .template get_default_linear_mapping<dim, spacedim>(),
99  fe,
100  quadrature,
101  update_flags,
102  neighbor_update_flags,
103  face_quadrature,
104  face_update_flags,
105  neighbor_face_update_flags)
106  {}
107 
108 
109 
110  template <int dim, int spacedim>
112  const ScratchData<dim, spacedim> &scratch)
113  : mapping(scratch.mapping)
114  , fe(scratch.fe)
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)
125  {}
126 
127 
128 
129  template <int dim, int spacedim>
133  {
134  if (!fe_values)
135  fe_values = std::make_unique<FEValues<dim, spacedim>>(*mapping,
136  *fe,
137  cell_quadrature,
138  cell_update_flags);
139 
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();
144  return *fe_values;
145  }
146 
147 
148 
149  template <int dim, int spacedim>
153  const unsigned int face_no)
154  {
155  if (!fe_face_values)
156  fe_face_values = std::make_unique<FEFaceValues<dim, spacedim>>(
157  *mapping, *fe, face_quadrature, face_update_flags);
158 
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;
164  }
165 
166 
167 
168  template <int dim, int spacedim>
172  const unsigned int face_no,
173  const unsigned int subface_no)
174  {
175  if (subface_no != numbers::invalid_unsigned_int)
176  {
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);
183 
184  current_fe_values = fe_subface_values.get();
185  return *fe_subface_values;
186  }
187  else
188  return reinit(cell, face_no);
189  }
190 
191 
192 
193  template <int dim, int spacedim>
197  const unsigned int face_no,
198  const unsigned int sub_face_no,
200  & cell_neighbor,
201  const unsigned int face_no_neighbor,
202  const unsigned int sub_face_no_neighbor)
203  {
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,
208  face_no,
209  sub_face_no,
210  cell_neighbor,
211  face_no_neighbor,
212  sub_face_no_neighbor);
213 
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);
216 
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;
220  }
221 
222 
223 
224  template <int dim, int spacedim>
228  {
229  if (!neighbor_fe_values)
230  neighbor_fe_values = std::make_unique<FEValues<dim, spacedim>>(
231  *mapping, *fe, cell_quadrature, neighbor_cell_update_flags);
232 
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;
237  }
238 
239 
240 
241  template <int dim, int spacedim>
245  const unsigned int face_no)
246  {
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;
254  }
255 
256 
257 
258  template <int dim, int spacedim>
262  const unsigned int face_no,
263  const unsigned int subface_no)
264  {
265  if (subface_no != numbers::invalid_unsigned_int)
266  {
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;
275  }
276  else
277  return reinit_neighbor(cell, face_no);
278  }
279 
280 
281 
282  template <int dim, int spacedim>
285  {
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;
290  }
291 
292 
293 
294  template <int dim, int spacedim>
297  {
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;
302  }
303 
304 
305 
306  template <int dim, int spacedim>
307  const std::vector<Point<spacedim>> &
309  {
310  return get_current_fe_values().get_quadrature_points();
311  }
312 
313 
314 
315  template <int dim, int spacedim>
316  const std::vector<double> &
318  {
319  return get_current_fe_values().get_JxW_values();
320  }
321 
322 
323 
324  template <int dim, int spacedim>
325  const std::vector<double> &
327  {
328  return get_current_neighbor_fe_values().get_JxW_values();
329  }
330 
331 
332 
333  template <int dim, int spacedim>
334  const std::vector<Tensor<1, spacedim>> &
336  {
337  return get_current_fe_values().get_normal_vectors();
338  }
339 
340 
341 
342  template <int dim, int spacedim>
343  const std::vector<Tensor<1, spacedim>> &
345  {
346  return get_current_neighbor_fe_values().get_normal_vectors();
347  }
348 
349 
350 
351  template <int dim, int spacedim>
352  const std::vector<types::global_dof_index> &
354  {
355  return local_dof_indices;
356  }
357 
358 
359 
360  template <int dim, int spacedim>
361  const std::vector<types::global_dof_index> &
363  {
364  return neighbor_dof_indices;
365  }
366 
367 
368 
369  template <int dim, int spacedim>
372  {
373  return user_data_storage;
374  }
375 
376 
377 
378  template <int dim, int spacedim>
379  const GeneralDataStorage &
381  {
382  return user_data_storage;
383  }
384 
385 
386 
387  template <int dim, int spacedim>
388  const Mapping<dim, spacedim> &
390  {
391  return *mapping;
392  }
393 
394 } // namespace MeshWorker
396 
397 // Explicit instantiations
399 namespace MeshWorker
400 {
401 #include "scratch_data.inst"
402 }
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access >> &cell)
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)
Definition: scratch_data.cc:25
const FEValues< dim, spacedim > & reinit(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:396
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:397
UpdateFlags
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcMessage(std::string arg1)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition: mapping.cc:240
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618
static const unsigned int invalid_unsigned_int
Definition: types.h:196