Reference documentation for deal.II version 8.4.2
dof_info.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 2015 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/base/std_cxx11/shared_ptr.h>
23 #include <deal.II/dofs/block_info.h>
24 #include <deal.II/fe/fe_values.h>
25 #include <deal.II/meshworker/local_results.h>
26 #include <deal.II/meshworker/vector_selector.h>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 namespace MeshWorker
31 {
32  template <int dim, class DOFINFO> class DoFInfoBox;
33 
34 
67  template<int dim, int spacedim = dim, typename number = double>
68  class DoFInfo : public LocalResults<number>
69  {
70  public:
73 
76 
83  unsigned int face_number;
84 
91  unsigned int sub_number;
92 
93  /*
94  * The DoF indices of the
95  * current cell
96  */
97  std::vector<types::global_dof_index> indices;
98 
103  std::vector<std::vector<types::global_dof_index> > indices_by_block;
104 
108  DoFInfo(const BlockInfo &block_info);
109 
114  DoFInfo (const DoFHandler<dim, spacedim> &dof_handler);
115 
119  template <class DHCellIterator>
120  void reinit(const DHCellIterator &c);
121 
125  template <class DHCellIterator, class DHFaceIterator>
126  void reinit(const DHCellIterator &c,
127  const DHFaceIterator &f,
128  const unsigned int face_no);
129 
133  template <class DHCellIterator, class DHFaceIterator>
134  void reinit(const DHCellIterator &c,
135  const DHFaceIterator &f,
136  const unsigned int face_no,
137  const unsigned int subface_no);
138 
143  template <class DHFaceIterator>
144  void set_face (const DHFaceIterator &f,
145  const unsigned int face_no);
146 
151  template <class DHFaceIterator>
152  void set_subface (const DHFaceIterator &f,
153  const unsigned int face_no,
154  const unsigned int subface_no);
155 
156  const BlockIndices &local_indices() const;
157 
158 
161 
166 
167  private:
173  DoFInfo ();
174 
176  void set_block_indices ();
177 
179  template <class DHCellIterator>
180  void get_indices(const DHCellIterator &c);
181 
183  std::vector<types::global_dof_index> indices_org;
184 
191 
192  friend class DoFInfoBox<dim, DoFInfo<dim, spacedim, number> >;
193  };
194 
195 
207  template <int dim, class DOFINFO>
208  class DoFInfoBox
209  {
210  public:
214  DoFInfoBox(const DOFINFO &seed);
215 
221 
225  void reset();
226 
232  template <class ASSEMBLER>
233  void assemble(ASSEMBLER &ass) const;
234 
238  std::size_t memory_consumption () const;
239 
240 
244  DOFINFO cell;
253 
258  bool interior_face_available[GeometryInfo<dim>::faces_per_cell];
259 
264  bool exterior_face_available[GeometryInfo<dim>::faces_per_cell];
265 
270  };
271 
272 //----------------------------------------------------------------------//
273 
274  template <int dim, int spacedim, typename number>
276  :
277  level_cell (false)
278  {
279  std::vector<types::global_dof_index> aux(1);
280  aux[0] = dof_handler.get_fe().dofs_per_cell;
282  }
283 
284 
285  template <int dim, int spacedim, typename number>
286  template <class DHCellIterator>
287  inline void
289  {
290  indices.resize(c->get_fe().dofs_per_cell);
291  if (block_info == 0 || block_info->local().size() == 0)
292  c->get_active_or_mg_dof_indices(indices);
293  else
294  {
295  indices_org.resize(c->get_fe().dofs_per_cell);
296  c->get_active_or_mg_dof_indices(indices_org);
298  }
299  }
300 
301 
302  template <int dim, int spacedim, typename number>
303  template <class DHCellIterator>
304  inline void
305  DoFInfo<dim,spacedim,number>::reinit(const DHCellIterator &c)
306  {
307  get_indices(c);
308  level_cell = c->is_level_cell();
309 
313  if (block_info)
315  else
317  }
318 
319 
320  template<int dim, int spacedim, typename number>
321  template <class DHFaceIterator>
322  inline void
324  const DHFaceIterator &f,
325  const unsigned int face_no)
326  {
327  face = static_cast<typename Triangulation<dim>::face_iterator> (f);
328  face_number = face_no;
330  }
331 
332 
333  template<int dim, int spacedim, typename number>
334  template <class DHCellIterator, class DHFaceIterator>
335  inline void
337  const DHCellIterator &c,
338  const DHFaceIterator &f,
339  const unsigned int face_no)
340  {
341  if ((cell.state() != IteratorState::valid)
342  || cell != typename Triangulation<dim>::cell_iterator(*c))
343  get_indices(c);
344  level_cell = c->is_level_cell();
345 
347  set_face(f,face_no);
348 
349  if (block_info)
351  else
353  }
354 
355 
356  template<int dim, int spacedim, typename number>
357  template <class DHFaceIterator>
358  inline void
360  const DHFaceIterator &f,
361  const unsigned int face_no,
362  const unsigned int subface_no)
363  {
364  face = static_cast<typename Triangulation<dim>::face_iterator> (f);
365  face_number = face_no;
366  sub_number = subface_no;
367  }
368 
369 
370  template<int dim, int spacedim, typename number>
371  template <class DHCellIterator, class DHFaceIterator>
372  inline void
374  const DHCellIterator &c,
375  const DHFaceIterator &f,
376  const unsigned int face_no,
377  const unsigned int subface_no)
378  {
379  if (cell.state() != IteratorState::valid
380  || cell != static_cast<typename Triangulation<dim>::cell_iterator> (c))
381  get_indices(c);
382  level_cell = c->is_level_cell();
383 
384  cell = static_cast<typename Triangulation<dim>::cell_iterator> (c);
385  set_subface(f, face_no, subface_no);
386 
387  if (block_info)
389  else
391  }
392 
393 
394  template<int dim, int spacedim, typename number>
395  inline const BlockIndices &
397  {
398  if (block_info)
399  return block_info->local();
400  return aux_local_indices;
401  }
402 
403 //----------------------------------------------------------------------//
404 
405  template <int dim, class DOFINFO>
406  inline
408  :
409  cell(seed), cell_valid(true)
410  {
411  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
412  {
413  exterior[i] = seed;
414  interior[i] = seed;
415  interior_face_available[i] = false;
416  exterior_face_available[i] = false;
417  }
418  }
419 
420 
421  template <int dim, class DOFINFO>
422  inline
424  :
425  cell(other.cell), cell_valid(other.cell_valid)
426  {
427  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
428  {
429  exterior[i] = other.exterior[i];
430  interior[i] = other.interior[i];
431  interior_face_available[i] = false;
432  exterior_face_available[i] = false;
433  }
434  }
435 
436 
437  template <int dim, class DOFINFO>
438  inline void
440  {
441  cell_valid = false;
442  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
443  {
444  interior_face_available[i] = false;
445  exterior_face_available[i] = false;
446  }
447  }
448 
449 
450  template <int dim, class DOFINFO>
451  template <class ASSEMBLER>
452  inline void
453  DoFInfoBox<dim, DOFINFO>::assemble (ASSEMBLER &assembler) const
454  {
455  if (!cell_valid)
456  return;
457 
458  assembler.assemble(cell);
459  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
460  {
461  // Only do something if data available
463  {
464  // If both data
465  // available, it is an
466  // interior face
468  assembler.assemble(interior[i], exterior[i]);
469  else
470  assembler.assemble(interior[i]);
471  }
472  }
473  }
474 }
475 
476 DEAL_II_NAMESPACE_CLOSE
477 
478 #endif
void get_indices(const DHCellIterator &c)
Fill index vector with active indices.
Definition: dof_info.h:288
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition: dof_info.h:72
static const unsigned int invalid_unsigned_int
Definition: types.h:164
std::vector< types::global_dof_index > indices_org
Auxiliary vector.
Definition: dof_info.h:183
Auxiliary class aiding in the handling of block structures like in BlockVector or FESystem...
Definition: block_indices.h:54
void set_block_indices()
Set up local block indices.
unsigned int sub_number
Definition: dof_info.h:91
const FiniteElement< dim, spacedim > & get_fe() const
bool exterior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:264
BlockIndices aux_local_indices
Definition: dof_info.h:190
void set_subface(const DHFaceIterator &f, const unsigned int face_no, const unsigned int subface_no)
Definition: dof_info.h:359
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition: tria.h:1392
bool interior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:258
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
SmartPointer< const BlockInfo, DoFInfo< dim, spacedim > > block_info
The block structure of the system.
Definition: dof_info.h:160
DoFInfoBox(const DOFINFO &seed)
Definition: dof_info.h:407
void assemble(ASSEMBLER &ass) const
Definition: dof_info.h:453
std::vector< std::vector< types::global_dof_index > > indices_by_block
Definition: dof_info.h:103
void reinit(const BlockIndices &local_sizes)
Definition: mesh_worker.cc:27
std::size_t memory_consumption() const
Definition: mesh_worker.cc:45
DOFINFO interior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:248
unsigned int face_number
Definition: dof_info.h:83
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:252
Triangulation< dim, spacedim >::face_iterator face
The current face.
Definition: dof_info.h:75
void reinit(const DHCellIterator &c)
Definition: dof_info.h:305
void set_face(const DHFaceIterator &f, const unsigned int face_no)
Definition: dof_info.h:323