Reference documentation for deal.II version GIT 29f9da0a34 2023-12-07 10:00:01+00:00
\(\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\}}\)
block_info.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 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 
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 
24 
26 
27 
28 template <int dim, int spacedim>
29 void
31  bool levels_only,
32  bool active_only)
33 {
34  Assert(dof.has_hp_capabilities() == false,
36 
37  if (!levels_only && dof.has_active_dofs())
38  {
39  const std::vector<types::global_dof_index> sizes =
41  bi_global.reinit(sizes);
42  }
43 
44  if (!active_only && dof.has_level_dofs())
45  {
46  std::vector<std::vector<types::global_dof_index>> sizes(
48  std::vector<types::global_dof_index>(dof.get_fe().n_blocks()));
49 
51  levels.resize(sizes.size());
52 
53  for (unsigned int i = 0; i < sizes.size(); ++i)
54  levels[i].reinit(sizes[i]);
55  }
56 }
57 
58 
59 
60 template <int dim, int spacedim>
61 void
63 {
64  Assert(dof.has_hp_capabilities() == false,
66 
67  const FiniteElement<dim, spacedim> &fe = dof.get_fe();
68  std::vector<types::global_dof_index> sizes(fe.n_blocks());
69 
70  base_elements.resize(fe.n_blocks());
71 
72  for (unsigned int i = 0; i < base_elements.size(); ++i)
73  base_elements[i] = fe.block_to_base_index(i).first;
74 
77  bi_local.reinit(sizes);
78 }
79 
80 // explicit instantiations
81 #include "block_info.inst"
82 
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
std::vector< types::global_dof_index > local_renumbering
Definition: block_info.h:214
std::vector< BlockIndices > levels
The multilevel block structure.
Definition: block_info.h:198
BlockIndices bi_global
The block structure of the global system.
Definition: block_info.h:194
BlockIndices bi_local
The block structure of the cell systems.
Definition: block_info.h:203
std::vector< unsigned int > base_elements
Definition: block_info.h:208
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:62
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
bool has_active_dofs() const
const Triangulation< dim, spacedim > & get_triangulation() const
bool has_hp_capabilities() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
bool has_level_dofs() const
unsigned int n_dofs_per_cell() const
unsigned int n_blocks() const
std::pair< unsigned int, unsigned int > block_to_base_index(const unsigned int block) const
unsigned int n_levels() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:477
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:478
#define Assert(cond, exc)
Definition: exceptions.h:1631
std::vector< types::global_dof_index > count_dofs_per_fe_block(const DoFHandler< dim, spacedim > &dof, const std::vector< unsigned int > &target_block=std::vector< unsigned int >())
Definition: dof_tools.cc:2017
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)
void count_dofs_per_block(const DoFHandler< dim, spacedim > &dof_handler, std::vector< std::vector< types::global_dof_index >> &dofs_per_block, std::vector< unsigned int > target_block={})
Definition: mg_tools.cc:1157
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618