Reference documentation for deal.II version 8.5.1
dof_info.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 2016 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 
259 
265 
270  };
271 
272 //----------------------------------------------------------------------//
273 
274  template <int dim, int spacedim, typename number>
276  :
277  face_number (numbers::invalid_unsigned_int),
278  sub_number (numbers::invalid_unsigned_int),
279  level_cell (false)
280  {
281  std::vector<types::global_dof_index> aux(1);
282  aux[0] = dof_handler.get_fe().dofs_per_cell;
284  }
285 
286 
287  template <int dim, int spacedim, typename number>
288  template <class DHCellIterator>
289  inline void
291  {
292  indices.resize(c->get_fe().dofs_per_cell);
293  if (block_info == 0 || block_info->local().size() == 0)
294  c->get_active_or_mg_dof_indices(indices);
295  else
296  {
297  indices_org.resize(c->get_fe().dofs_per_cell);
298  c->get_active_or_mg_dof_indices(indices_org);
299  set_block_indices();
300  }
301  }
302 
303 
304  template <int dim, int spacedim, typename number>
305  template <class DHCellIterator>
306  inline void
307  DoFInfo<dim,spacedim,number>::reinit(const DHCellIterator &c)
308  {
309  get_indices(c);
310  level_cell = c->is_level_cell();
311 
312  cell = typename Triangulation<dim,spacedim>::cell_iterator(*c);
313  face_number = numbers::invalid_unsigned_int;
314  sub_number = numbers::invalid_unsigned_int;
315  if (block_info)
316  LocalResults<number>::reinit(block_info->local());
317  else
318  LocalResults<number>::reinit(aux_local_indices);
319  }
320 
321 
322  template<int dim, int spacedim, typename number>
323  template <class DHFaceIterator>
324  inline void
326  const DHFaceIterator &f,
327  const unsigned int face_no)
328  {
329  face = static_cast<typename Triangulation<dim>::face_iterator> (f);
330  face_number = face_no;
331  sub_number = numbers::invalid_unsigned_int;
332  }
333 
334 
335  template<int dim, int spacedim, typename number>
336  template <class DHCellIterator, class DHFaceIterator>
337  inline void
339  const DHCellIterator &c,
340  const DHFaceIterator &f,
341  const unsigned int face_no)
342  {
343  if ((cell.state() != IteratorState::valid)
344  || cell != typename Triangulation<dim>::cell_iterator(*c))
345  get_indices(c);
346  level_cell = c->is_level_cell();
347 
348  cell = typename Triangulation<dim>::cell_iterator(*c);
349  set_face(f,face_no);
350 
351  if (block_info)
352  LocalResults<number>::reinit(block_info->local());
353  else
354  LocalResults<number>::reinit(aux_local_indices);
355  }
356 
357 
358  template<int dim, int spacedim, typename number>
359  template <class DHFaceIterator>
360  inline void
362  const DHFaceIterator &f,
363  const unsigned int face_no,
364  const unsigned int subface_no)
365  {
366  face = static_cast<typename Triangulation<dim>::face_iterator> (f);
367  face_number = face_no;
368  sub_number = subface_no;
369  }
370 
371 
372  template<int dim, int spacedim, typename number>
373  template <class DHCellIterator, class DHFaceIterator>
374  inline void
376  const DHCellIterator &c,
377  const DHFaceIterator &f,
378  const unsigned int face_no,
379  const unsigned int subface_no)
380  {
381  if (cell.state() != IteratorState::valid
382  || cell != static_cast<typename Triangulation<dim>::cell_iterator> (c))
383  get_indices(c);
384  level_cell = c->is_level_cell();
385 
386  cell = static_cast<typename Triangulation<dim>::cell_iterator> (c);
387  set_subface(f, face_no, subface_no);
388 
389  if (block_info)
390  LocalResults<number>::reinit(block_info->local());
391  else
392  LocalResults<number>::reinit(aux_local_indices);
393  }
394 
395 
396  template<int dim, int spacedim, typename number>
397  inline const BlockIndices &
399  {
400  if (block_info)
401  return block_info->local();
402  return aux_local_indices;
403  }
404 
405 //----------------------------------------------------------------------//
406 
407  template <int dim, class DOFINFO>
408  inline
410  :
411  cell(seed), cell_valid(true)
412  {
413  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
414  {
415  exterior[i] = seed;
416  interior[i] = seed;
417  interior_face_available[i] = false;
418  exterior_face_available[i] = false;
419  }
420  }
421 
422 
423  template <int dim, class DOFINFO>
424  inline
426  :
427  cell(other.cell), cell_valid(other.cell_valid)
428  {
429  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
430  {
431  exterior[i] = other.exterior[i];
432  interior[i] = other.interior[i];
433  interior_face_available[i] = false;
434  exterior_face_available[i] = false;
435  }
436  }
437 
438 
439  template <int dim, class DOFINFO>
440  inline void
442  {
443  cell_valid = false;
444  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
445  {
446  interior_face_available[i] = false;
447  exterior_face_available[i] = false;
448  }
449  }
450 
451 
452  template <int dim, class DOFINFO>
453  template <class ASSEMBLER>
454  inline void
455  DoFInfoBox<dim, DOFINFO>::assemble (ASSEMBLER &assembler) const
456  {
457  if (!cell_valid)
458  return;
459 
460  assembler.assemble(cell);
461  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
462  {
463  // Only do something if data available
464  if (interior_face_available[i])
465  {
466  // If both data
467  // available, it is an
468  // interior face
469  if (exterior_face_available[i])
470  assembler.assemble(interior[i], exterior[i]);
471  else
472  assembler.assemble(interior[i]);
473  }
474  }
475  }
476 }
477 
478 DEAL_II_NAMESPACE_CLOSE
479 
480 #endif
void get_indices(const DHCellIterator &c)
Fill index vector with active indices.
Definition: dof_info.h:290
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition: dof_info.h:72
static const unsigned int invalid_unsigned_int
Definition: types.h:170
std::vector< types::global_dof_index > indices_org
Auxiliary vector.
Definition: dof_info.h:183
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:361
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition: tria.h:1481
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:409
void assemble(ASSEMBLER &ass) const
Definition: dof_info.h:455
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
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
std::size_t memory_consumption() const
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:307
void set_face(const DHFaceIterator &f, const unsigned int face_no)
Definition: dof_info.h:325