Reference documentation for deal.II version 9.2.0
\(\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\}}\)
dof_info.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 2020 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 
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 
23 
25 
26 #include <deal.II/fe/fe_values.h>
27 
30 
31 #include <memory>
32 
34 
35 namespace MeshWorker
36 {
37  // Forward declaration
38 #ifndef DOXYGEN
39  template <int dim, class DOFINFO>
40  class DoFInfoBox;
41 #endif
42 
75  template <int dim, int spacedim = dim, typename number = double>
76  class DoFInfo : public LocalResults<number>
77  {
78  public:
81 
84 
91  unsigned int face_number;
92 
99  unsigned int sub_number;
100 
105  std::vector<types::global_dof_index> indices;
106 
111  std::vector<std::vector<types::global_dof_index>> indices_by_block;
112 
116  DoFInfo(const BlockInfo &block_info);
117 
122  DoFInfo(const DoFHandler<dim, spacedim> &dof_handler);
123 
127  template <class DHCellIterator>
128  void
129  reinit(const DHCellIterator &c);
130 
134  template <class DHCellIterator, class DHFaceIterator>
135  void
136  reinit(const DHCellIterator &c,
137  const DHFaceIterator &f,
138  const unsigned int face_no);
139 
143  template <class DHCellIterator, class DHFaceIterator>
144  void
145  reinit(const DHCellIterator &c,
146  const DHFaceIterator &f,
147  const unsigned int face_no,
148  const unsigned int subface_no);
149 
154  template <class DHFaceIterator>
155  void
156  set_face(const DHFaceIterator &f, const unsigned int face_no);
157 
162  template <class DHFaceIterator>
163  void
164  set_subface(const DHFaceIterator &f,
165  const unsigned int face_no,
166  const unsigned int subface_no);
167 
168  const BlockIndices &
169  local_indices() const;
170 
171 
174 
179 
180  private:
186  DoFInfo();
187 
189  void
191 
193  template <class DHCellIterator>
194  void
195  get_indices(const DHCellIterator &c);
196 
198  std::vector<types::global_dof_index> indices_org;
199 
206 
207  friend class DoFInfoBox<dim, DoFInfo<dim, spacedim, number>>;
208  };
209 
210 
222  template <int dim, class DOFINFO>
224  {
225  public:
229  DoFInfoBox(const DOFINFO &seed);
230 
236 
240  DoFInfoBox &
242 
246  void
247  reset();
248 
254  template <class ASSEMBLER>
255  void
256  assemble(ASSEMBLER &ass) const;
257 
258 
262  DOFINFO cell;
271 
277 
283 
288  };
289 
290  //----------------------------------------------------------------------//
291 
292  template <int dim, int spacedim, typename number>
294  : face_number(numbers::invalid_unsigned_int)
295  , sub_number(numbers::invalid_unsigned_int)
296  , level_cell(false)
297  {}
298 
299 
300 
301  template <int dim, int spacedim, typename number>
303  const DoFHandler<dim, spacedim> &dof_handler)
304  : face_number(numbers::invalid_unsigned_int)
305  , sub_number(numbers::invalid_unsigned_int)
306  , level_cell(false)
307  {
308  std::vector<types::global_dof_index> aux(1);
309  aux[0] = dof_handler.get_fe().dofs_per_cell;
311  }
312 
313 
314  template <int dim, int spacedim, typename number>
315  template <class DHCellIterator>
316  inline void
318  {
319  indices.resize(c->get_fe().dofs_per_cell);
320  if (block_info == nullptr || block_info->local().size() == 0)
321  c->get_active_or_mg_dof_indices(indices);
322  else
323  {
324  indices_org.resize(c->get_fe().dofs_per_cell);
325  c->get_active_or_mg_dof_indices(indices_org);
326  set_block_indices();
327  }
328  }
329 
330 
331  template <int dim, int spacedim, typename number>
332  template <class DHCellIterator>
333  inline void
334  DoFInfo<dim, spacedim, number>::reinit(const DHCellIterator &c)
335  {
336  get_indices(c);
337  level_cell = c->is_level_cell();
338 
340  face_number = numbers::invalid_unsigned_int;
341  sub_number = numbers::invalid_unsigned_int;
342  if (block_info)
343  LocalResults<number>::reinit(block_info->local());
344  else
345  LocalResults<number>::reinit(aux_local_indices);
346  }
347 
348 
349  template <int dim, int spacedim, typename number>
350  template <class DHFaceIterator>
351  inline void
353  const unsigned int face_no)
354  {
355  face = static_cast<typename Triangulation<dim, spacedim>::face_iterator>(f);
356  face_number = face_no;
357  sub_number = numbers::invalid_unsigned_int;
358  }
359 
360 
361  template <int dim, int spacedim, typename number>
362  template <class DHCellIterator, class DHFaceIterator>
363  inline void
364  DoFInfo<dim, spacedim, number>::reinit(const DHCellIterator &c,
365  const DHFaceIterator &f,
366  const unsigned int face_no)
367  {
368  if ((cell.state() != IteratorState::valid) ||
369  cell != typename Triangulation<dim, spacedim>::cell_iterator(*c))
370  get_indices(c);
371  level_cell = c->is_level_cell();
372 
374  set_face(f, face_no);
375 
376  if (block_info)
377  LocalResults<number>::reinit(block_info->local());
378  else
379  LocalResults<number>::reinit(aux_local_indices);
380  }
381 
382 
383  template <int dim, int spacedim, typename number>
384  template <class DHFaceIterator>
385  inline void
387  const unsigned int face_no,
388  const unsigned int subface_no)
389  {
390  face = static_cast<typename Triangulation<dim, spacedim>::face_iterator>(f);
391  face_number = face_no;
392  sub_number = subface_no;
393  }
394 
395 
396  template <int dim, int spacedim, typename number>
397  template <class DHCellIterator, class DHFaceIterator>
398  inline void
399  DoFInfo<dim, spacedim, number>::reinit(const DHCellIterator &c,
400  const DHFaceIterator &f,
401  const unsigned int face_no,
402  const unsigned int subface_no)
403  {
404  if (cell.state() != IteratorState::valid ||
405  cell !=
406  static_cast<typename Triangulation<dim, spacedim>::cell_iterator>(c))
407  get_indices(c);
408  level_cell = c->is_level_cell();
409 
410  cell = static_cast<typename Triangulation<dim, spacedim>::cell_iterator>(c);
411  set_subface(f, face_no, subface_no);
412 
413  if (block_info)
414  LocalResults<number>::reinit(block_info->local());
415  else
416  LocalResults<number>::reinit(aux_local_indices);
417  }
418 
419 
420  template <int dim, int spacedim, typename number>
421  inline const BlockIndices &
423  {
424  if (block_info)
425  return block_info->local();
426  return aux_local_indices;
427  }
428 
429  //----------------------------------------------------------------------//
430 
431  template <int dim, class DOFINFO>
432  inline DoFInfoBox<dim, DOFINFO>::DoFInfoBox(const DOFINFO &seed)
433  : cell(seed)
434  , cell_valid(true)
435  {
436  for (unsigned int i : GeometryInfo<dim>::face_indices())
437  {
438  exterior[i] = seed;
439  interior[i] = seed;
440  interior_face_available[i] = false;
441  exterior_face_available[i] = false;
442  }
443  }
444 
445 
446  template <int dim, class DOFINFO>
448  const DoFInfoBox<dim, DOFINFO> &other)
449  : cell(other.cell)
450  , cell_valid(other.cell_valid)
451  {
452  for (unsigned int i : GeometryInfo<dim>::face_indices())
453  {
454  exterior[i] = other.exterior[i];
455  interior[i] = other.interior[i];
456  interior_face_available[i] = false;
457  exterior_face_available[i] = false;
458  }
459  }
460 
461 
462  template <int dim, class DOFINFO>
463  inline DoFInfoBox<dim, DOFINFO> &
465  {
466  cell = other.cell;
467  cell_valid = other.cell_valid;
468  for (unsigned int i : GeometryInfo<dim>::face_indices())
469  {
470  exterior[i] = other.exterior[i];
471  interior[i] = other.interior[i];
472  interior_face_available[i] = false;
473  exterior_face_available[i] = false;
474  }
475  return *this;
476  }
477 
478 
479  template <int dim, class DOFINFO>
480  inline void
482  {
483  cell_valid = false;
484  for (unsigned int i : GeometryInfo<dim>::face_indices())
485  {
486  interior_face_available[i] = false;
487  exterior_face_available[i] = false;
488  }
489  }
490 
491 
492  template <int dim, class DOFINFO>
493  template <class ASSEMBLER>
494  inline void
495  DoFInfoBox<dim, DOFINFO>::assemble(ASSEMBLER &assembler) const
496  {
497  if (!cell_valid)
498  return;
499 
500  assembler.assemble(cell);
501  for (unsigned int i : GeometryInfo<dim>::face_indices())
502  {
503  // Only do something if data available
504  if (interior_face_available[i])
505  {
506  // If both data
507  // available, it is an
508  // interior face
509  if (exterior_face_available[i])
510  assembler.assemble(interior[i], exterior[i]);
511  else
512  assembler.assemble(interior[i]);
513  }
514  }
515  }
516 } // namespace MeshWorker
517 
519 
520 #endif
MeshWorker::DoFInfoBox::interior_face_available
bool interior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:276
MeshWorker::DoFInfo::set_block_indices
void set_block_indices()
Set up local block indices.
fe_values.h
MeshWorker::DoFInfoBox::DoFInfoBox
DoFInfoBox(const DOFINFO &seed)
Definition: dof_info.h:432
vector_selector.h
MeshWorker::DoFInfoBox::cell
DOFINFO cell
Definition: dof_info.h:262
MeshWorker
Definition: assemble_flags.h:30
MeshWorker::DoFInfo::reinit
void reinit(const DHCellIterator &c)
Definition: dof_info.h:334
MeshWorker::DoFInfo::set_subface
void set_subface(const DHFaceIterator &f, const unsigned int face_no, const unsigned int subface_no)
Definition: dof_info.h:386
MeshWorker::DoFInfoBox::interior
DOFINFO interior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:266
local_results.h
GeometryInfo
Definition: geometry_info.h:1224
MeshWorker::DoFInfo::block_info
SmartPointer< const BlockInfo, DoFInfo< dim, spacedim > > block_info
The block structure of the system.
Definition: dof_info.h:173
IteratorState::valid
@ valid
Iterator points to a valid object.
Definition: tria_iterator_base.h:38
MeshWorker::DoFInfoBox::exterior
DOFINFO exterior[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:270
MeshWorker::DoFInfo::sub_number
unsigned int sub_number
Definition: dof_info.h:99
MeshWorker::DoFInfo::aux_local_indices
BlockIndices aux_local_indices
Definition: dof_info.h:205
MeshWorker::DoFInfoBox::reset
void reset()
Definition: dof_info.h:481
DoFHandler
Definition: dof_handler.h:205
quadrature_lib.h
MeshWorker::DoFInfoBox::cell_valid
bool cell_valid
Definition: dof_info.h:287
MeshWorker::DoFInfo::local_indices
const BlockIndices & local_indices() const
Definition: dof_info.h:422
MeshWorker::LocalResults
Definition: local_results.h:226
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
MeshWorker::DoFInfoBox::operator=
DoFInfoBox & operator=(const DoFInfoBox< dim, DOFINFO > &)
Definition: dof_info.h:464
numbers
Definition: numbers.h:207
MeshWorker::DoFInfo::indices
std::vector< types::global_dof_index > indices
Definition: dof_info.h:105
MeshWorker::DoFInfo::level_cell
bool level_cell
Definition: dof_info.h:178
SmartPointer
Definition: smartpointer.h:68
Triangulation::cell_iterator
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition: tria.h:1343
MeshWorker::DoFInfoBox::assemble
void assemble(ASSEMBLER &ass) const
Definition: dof_info.h:495
MeshWorker::DoFInfo::set_face
void set_face(const DHFaceIterator &f, const unsigned int face_no)
Definition: dof_info.h:352
block_info.h
BlockIndices::reinit
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
Definition: block_indices.h:260
MeshWorker::DoFInfo::face_number
unsigned int face_number
Definition: dof_info.h:91
MeshWorker::LocalResults::reinit
void reinit(const BlockIndices &local_sizes)
Definition: mesh_worker.cc:28
MeshWorker::DoFInfo::get_indices
void get_indices(const DHCellIterator &c)
Fill index vector with active indices.
Definition: dof_info.h:317
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
MeshWorker::DoFInfo::indices_org
std::vector< types::global_dof_index > indices_org
Auxiliary vector.
Definition: dof_info.h:198
MeshWorker::DoFInfoBox::exterior_face_available
bool exterior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition: dof_info.h:282
config.h
BlockInfo
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
Definition: block_info.h:99
MeshWorker::DoFInfo::indices_by_block
std::vector< std::vector< types::global_dof_index > > indices_by_block
Definition: dof_info.h:111
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
TriaIterator
Definition: tria_iterator.h:578
DoFHandler::get_fe
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
MeshWorker::DoFInfo::face
Triangulation< dim, spacedim >::face_iterator face
The current face.
Definition: dof_info.h:83
MeshWorker::DoFInfo
Definition: dof_info.h:76
MeshWorker::DoFInfo::cell
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition: dof_info.h:80
BlockIndices
Definition: block_indices.h:60
MeshWorker::DoFInfoBox
Definition: dof_info.h:223
MeshWorker::DoFInfo::DoFInfo
DoFInfo()
Definition: dof_info.h:293