Reference documentation for deal.II version 9.0.0
dof_info.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 2018 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii_mesh_worker_dof_info_h
18 #define dealii_mesh_worker_dof_info_h
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/quadrature_lib.h>
22 #include <deal.II/dofs/block_info.h>
23 #include <deal.II/fe/fe_values.h>
24 #include <deal.II/meshworker/local_results.h>
25 #include <deal.II/meshworker/vector_selector.h>
26 
27 #include <memory>
28 
29 DEAL_II_NAMESPACE_OPEN
30 
31 namespace MeshWorker
32 {
33  template <int dim, class DOFINFO> class DoFInfoBox;
34 
35 
68  template <int dim, int spacedim = dim, typename number = double>
69  class DoFInfo : public LocalResults<number>
70  {
71  public:
74 
77 
84  unsigned int face_number;
85 
92  unsigned int sub_number;
93 
94  /*
95  * The DoF indices of the
96  * current cell
97  */
98  std::vector<types::global_dof_index> indices;
99 
104  std::vector<std::vector<types::global_dof_index> > indices_by_block;
105 
109  DoFInfo(const BlockInfo &block_info);
110 
115  DoFInfo (const DoFHandler<dim, spacedim> &dof_handler);
116 
120  template <class DHCellIterator>
121  void reinit(const DHCellIterator &c);
122 
126  template <class DHCellIterator, class DHFaceIterator>
127  void reinit(const DHCellIterator &c,
128  const DHFaceIterator &f,
129  const unsigned int face_no);
130 
134  template <class DHCellIterator, class DHFaceIterator>
135  void reinit(const DHCellIterator &c,
136  const DHFaceIterator &f,
137  const unsigned int face_no,
138  const unsigned int subface_no);
139 
144  template <class DHFaceIterator>
145  void set_face (const DHFaceIterator &f,
146  const unsigned int face_no);
147 
152  template <class DHFaceIterator>
153  void set_subface (const DHFaceIterator &f,
154  const unsigned int face_no,
155  const unsigned int subface_no);
156 
157  const BlockIndices &local_indices() const;
158 
159 
162 
167 
168  private:
174  DoFInfo ();
175 
177  void set_block_indices ();
178 
180  template <class DHCellIterator>
181  void get_indices(const DHCellIterator &c);
182 
184  std::vector<types::global_dof_index> indices_org;
185 
192 
193  friend class DoFInfoBox<dim, DoFInfo<dim, spacedim, number> >;
194  };
195 
196 
208  template <int dim, class DOFINFO>
209  class DoFInfoBox
210  {
211  public:
215  DoFInfoBox(const DOFINFO &seed);
216 
222 
226  void reset();
227 
233  template <class ASSEMBLER>
234  void assemble(ASSEMBLER &ass) const;
235 
236 
240  DOFINFO cell;
249 
255 
261 
266  };
267 
268 //----------------------------------------------------------------------//
269 
270  template <int dim, int spacedim, typename number>
272  :
273  face_number(numbers::invalid_unsigned_int),
274  sub_number(numbers::invalid_unsigned_int),
275  level_cell(false)
276  {}
277 
278 
279 
280  template <int dim, int spacedim, typename number>
282  :
283  face_number (numbers::invalid_unsigned_int),
284  sub_number (numbers::invalid_unsigned_int),
285  level_cell (false)
286  {
287  std::vector<types::global_dof_index> aux(1);
288  aux[0] = dof_handler.get_fe().dofs_per_cell;
290  }
291 
292 
293  template <int dim, int spacedim, typename number>
294  template <class DHCellIterator>
295  inline void
297  {
298  indices.resize(c->get_fe().dofs_per_cell);
299  if (block_info == nullptr || block_info->local().size() == 0)
300  c->get_active_or_mg_dof_indices(indices);
301  else
302  {
303  indices_org.resize(c->get_fe().dofs_per_cell);
304  c->get_active_or_mg_dof_indices(indices_org);
305  set_block_indices();
306  }
307  }
308 
309 
310  template <int dim, int spacedim, typename number>
311  template <class DHCellIterator>
312  inline void
313  DoFInfo<dim,spacedim,number>::reinit(const DHCellIterator &c)
314  {
315  get_indices(c);
316  level_cell = c->is_level_cell();
317 
318  cell = typename Triangulation<dim,spacedim>::cell_iterator(*c);
319  face_number = numbers::invalid_unsigned_int;
320  sub_number = numbers::invalid_unsigned_int;
321  if (block_info)
322  LocalResults<number>::reinit(block_info->local());
323  else
324  LocalResults<number>::reinit(aux_local_indices);
325  }
326 
327 
328  template <int dim, int spacedim, typename number>
329  template <class DHFaceIterator>
330  inline void
332  const DHFaceIterator &f,
333  const unsigned int face_no)
334  {
335  face = static_cast<typename Triangulation<dim, spacedim>::face_iterator> (f);
336  face_number = face_no;
337  sub_number = numbers::invalid_unsigned_int;
338  }
339 
340 
341  template <int dim, int spacedim, typename number>
342  template <class DHCellIterator, class DHFaceIterator>
343  inline void
345  const DHCellIterator &c,
346  const DHFaceIterator &f,
347  const unsigned int face_no)
348  {
349  if ((cell.state() != IteratorState::valid)
350  || cell != typename Triangulation<dim, spacedim>::cell_iterator(*c))
351  get_indices(c);
352  level_cell = c->is_level_cell();
353 
355  set_face(f,face_no);
356 
357  if (block_info)
358  LocalResults<number>::reinit(block_info->local());
359  else
360  LocalResults<number>::reinit(aux_local_indices);
361  }
362 
363 
364  template <int dim, int spacedim, typename number>
365  template <class DHFaceIterator>
366  inline void
368  const DHFaceIterator &f,
369  const unsigned int face_no,
370  const unsigned int subface_no)
371  {
372  face = static_cast<typename Triangulation<dim, spacedim>::face_iterator> (f);
373  face_number = face_no;
374  sub_number = subface_no;
375  }
376 
377 
378  template <int dim, int spacedim, typename number>
379  template <class DHCellIterator, class DHFaceIterator>
380  inline void
382  const DHCellIterator &c,
383  const DHFaceIterator &f,
384  const unsigned int face_no,
385  const unsigned int subface_no)
386  {
387  if (cell.state() != IteratorState::valid
388  || cell != static_cast<typename Triangulation<dim, spacedim>::cell_iterator> (c))
389  get_indices(c);
390  level_cell = c->is_level_cell();
391 
392  cell = static_cast<typename Triangulation<dim, spacedim>::cell_iterator> (c);
393  set_subface(f, face_no, subface_no);
394 
395  if (block_info)
396  LocalResults<number>::reinit(block_info->local());
397  else
398  LocalResults<number>::reinit(aux_local_indices);
399  }
400 
401 
402  template <int dim, int spacedim, typename number>
403  inline const BlockIndices &
405  {
406  if (block_info)
407  return block_info->local();
408  return aux_local_indices;
409  }
410 
411 //----------------------------------------------------------------------//
412 
413  template <int dim, class DOFINFO>
414  inline
416  :
417  cell(seed), cell_valid(true)
418  {
419  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
420  {
421  exterior[i] = seed;
422  interior[i] = seed;
423  interior_face_available[i] = false;
424  exterior_face_available[i] = false;
425  }
426  }
427 
428 
429  template <int dim, class DOFINFO>
430  inline
432  :
433  cell(other.cell), cell_valid(other.cell_valid)
434  {
435  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
436  {
437  exterior[i] = other.exterior[i];
438  interior[i] = other.interior[i];
439  interior_face_available[i] = false;
440  exterior_face_available[i] = false;
441  }
442  }
443 
444 
445  template <int dim, class DOFINFO>
446  inline void
448  {
449  cell_valid = false;
450  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
451  {
452  interior_face_available[i] = false;
453  exterior_face_available[i] = false;
454  }
455  }
456 
457 
458  template <int dim, class DOFINFO>
459  template <class ASSEMBLER>
460  inline void
461  DoFInfoBox<dim, DOFINFO>::assemble (ASSEMBLER &assembler) const
462  {
463  if (!cell_valid)
464  return;
465 
466  assembler.assemble(cell);
467  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
468  {
469  // Only do something if data available
470  if (interior_face_available[i])
471  {
472  // If both data
473  // available, it is an
474  // interior face
475  if (exterior_face_available[i])
476  assembler.assemble(interior[i], exterior[i]);
477  else
478  assembler.assemble(interior[i]);
479  }
480  }
481  }
482 }
483 
484 DEAL_II_NAMESPACE_CLOSE
485 
486 #endif
void get_indices(const DHCellIterator &c)
Fill index vector with active indices.
Definition: dof_info.h:296
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition: dof_info.h:73
static const unsigned int invalid_unsigned_int
Definition: types.h:173
std::vector< types::global_dof_index > indices_org
Auxiliary vector.
Definition: dof_info.h:184
void set_block_indices()
Set up local block indices.
unsigned int sub_number
Definition: dof_info.h:92
bool exterior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:260
BlockIndices aux_local_indices
Definition: dof_info.h:191
void set_subface(const DHFaceIterator &f, const unsigned int face_no, const unsigned int subface_no)
Definition: dof_info.h:367
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition: tria.h:1509
bool interior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:254
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
SmartPointer< const BlockInfo, DoFInfo< dim, spacedim > > block_info
The block structure of the system.
Definition: dof_info.h:161
DoFInfoBox(const DOFINFO &seed)
Definition: dof_info.h:415
void assemble(ASSEMBLER &ass) const
Definition: dof_info.h:461
std::vector< std::vector< types::global_dof_index > > indices_by_block
Definition: dof_info.h:104
void reinit(const BlockIndices &local_sizes)
Definition: mesh_worker.cc:27
DOFINFO interior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:244
unsigned int face_number
Definition: dof_info.h:84
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
Definition: block_info.h:93
Iterator points to a valid object.
DOFINFO exterior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:248
Triangulation< dim, spacedim >::face_iterator face
The current face.
Definition: dof_info.h:76
void reinit(const DHCellIterator &c)
Definition: dof_info.h:313
void set_face(const DHFaceIterator &f, const unsigned int face_no)
Definition: dof_info.h:331