Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
block_info.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2019 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 #include <deal.II/dofs/block_info.h>
17 #include <deal.II/dofs/dof_handler.h>
18 #include <deal.II/dofs/dof_tools.h>
19 
20 #include <deal.II/fe/fe.h>
21 #include <deal.II/fe/fe_tools.h>
22 
23 #include <deal.II/multigrid/mg_tools.h>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 
28 template <int dim, int spacedim>
29 void
31  bool levels_only,
32  bool active_only)
33 {
34  if (!levels_only && dof.has_active_dofs())
35  {
36  const FiniteElement<dim, spacedim> & fe = dof.get_fe();
37  std::vector<types::global_dof_index> sizes(fe.n_blocks());
39  bi_global.reinit(sizes);
40  }
41 
42  if (!active_only && dof.has_level_dofs())
43  {
44  std::vector<std::vector<types::global_dof_index>> sizes(
45  dof.get_triangulation().n_levels(),
46  std::vector<types::global_dof_index>(dof.get_fe().n_blocks()));
47 
49  levels.resize(sizes.size());
50 
51  for (unsigned int i = 0; i < sizes.size(); ++i)
52  levels[i].reinit(sizes[i]);
53  }
54 }
55 
56 
57 template <int dim, int spacedim>
58 void
60 {
61  const FiniteElement<dim, spacedim> & fe = dof.get_fe();
62  std::vector<types::global_dof_index> sizes(fe.n_blocks());
63 
64  base_elements.resize(fe.n_blocks());
65 
66  for (unsigned int i = 0; i < base_elements.size(); ++i)
67  base_elements[i] = fe.block_to_base_index(i).first;
68 
71  bi_local.reinit(sizes);
72 }
73 
74 
75 // explicit instantiations
76 #include "block_info.inst"
77 
78 DEAL_II_NAMESPACE_CLOSE
const Triangulation< dim, spacedim > & get_triangulation() const
void compute_block_renumbering(const FiniteElement< dim, spacedim > &fe, std::vector< types::global_dof_index > &renumbering, std::vector< types::global_dof_index > &block_data, bool return_start_indices=true)
std::vector< types::global_dof_index > local_renumbering
Definition: block_info.h:214
void count_dofs_per_block(const DoFHandlerType &dof_handler, std::vector< std::vector< types::global_dof_index >> &dofs_per_block, std::vector< unsigned int > target_block={})
Definition: mg_tools.cc:1170
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
unsigned int n_blocks() const
bool has_active_dofs() const
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
std::pair< unsigned int, unsigned int > block_to_base_index(const unsigned int block) const
Definition: fe.h:3248
BlockIndices bi_global
The block structure of the global system.
Definition: block_info.h:194
void initialize(const DoFHandler< dim, spacedim > &, bool levels_only=false, bool active_only=false)
Fill the object with values describing block structure of the DoFHandler.
Definition: block_info.cc:30
void initialize_local(const DoFHandler< dim, spacedim > &)
Initialize block structure on cells and compute renumbering between cell dofs and block cell dofs...
Definition: block_info.cc:59
void count_dofs_per_block(const DoFHandlerType &dof, std::vector< types::global_dof_index > &dofs_per_block, const std::vector< unsigned int > &target_block=std::vector< unsigned int >())
Definition: dof_tools.cc:1953
std::vector< BlockIndices > levels
The multilevel block structure.
Definition: block_info.h:198
unsigned int n_dofs_per_cell() const
BlockIndices bi_local
The block structure of the cell systems.
Definition: block_info.h:203
bool has_level_dofs() const
std::vector< unsigned int > base_elements
Definition: block_info.h:208