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\}}\)
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  if (!levels_only && dof.has_active_dofs())
35  {
36  const std::vector<types::global_dof_index> sizes =
38  bi_global.reinit(sizes);
39  }
40 
41  if (!active_only && dof.has_level_dofs())
42  {
43  std::vector<std::vector<types::global_dof_index>> sizes(
44  dof.get_triangulation().n_levels(),
45  std::vector<types::global_dof_index>(dof.get_fe().n_blocks()));
46 
48  levels.resize(sizes.size());
49 
50  for (unsigned int i = 0; i < sizes.size(); ++i)
51  levels[i].reinit(sizes[i]);
52  }
53 }
54 
55 
56 
57 template <int dim, int spacedim>
58 void
60  bool levels_only,
61  bool active_only)
62 {
64  (void)dof;
65  (void)levels_only;
66  (void)active_only;
67 }
68 
69 
70 
71 template <int dim, int spacedim>
72 void
74 {
75  const FiniteElement<dim, spacedim> & fe = dof.get_fe();
76  std::vector<types::global_dof_index> sizes(fe.n_blocks());
77 
78  base_elements.resize(fe.n_blocks());
79 
80  for (unsigned int i = 0; i < base_elements.size(); ++i)
81  base_elements[i] = fe.block_to_base_index(i).first;
82 
85  bi_local.reinit(sizes);
86 }
87 
88 
89 
90 template <int dim, int spacedim>
91 void
93 {
95  (void)dof;
96 }
97 
98 // explicit instantiations
99 #include "block_info.inst"
100 
BlockInfo::initialize_local
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:73
MGTools::count_dofs_per_block
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:1158
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
BlockInfo::bi_local
BlockIndices bi_local
The block structure of the cell systems.
Definition: block_info.h:226
FETools::compute_block_renumbering
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)
hp::DoFHandler
Definition: dof_handler.h:203
DoFHandler
Definition: dof_handler.h:205
FiniteElement::block_to_base_index
std::pair< unsigned int, unsigned int > block_to_base_index(const unsigned int block) const
Definition: fe.h:3215
FiniteElementData::n_blocks
unsigned int n_blocks() const
BlockInfo::levels
std::vector< BlockIndices > levels
The multilevel block structure.
Definition: block_info.h:221
BlockInfo::initialize
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
mg_tools.h
DoFHandler::has_level_dofs
bool has_level_dofs() const
DoFHandler::get_triangulation
const Triangulation< dim, spacedim > & get_triangulation() const
internal::reinit
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:621
BlockInfo::local_renumbering
std::vector< types::global_dof_index > local_renumbering
Definition: block_info.h:237
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
fe.h
FiniteElement
Definition: fe.h:648
DoFTools::count_dofs_per_fe_block
std::vector< types::global_dof_index > count_dofs_per_fe_block(const DoFHandlerType &dof, const std::vector< unsigned int > &target_block=std::vector< unsigned int >())
Definition: dof_tools.cc:2012
fe_tools.h
dof_tools.h
BlockInfo::base_elements
std::vector< unsigned int > base_elements
Definition: block_info.h:231
block_info.h
BlockIndices::reinit
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
Definition: block_indices.h:260
dof_handler.h
DoFHandler::has_active_dofs
bool has_active_dofs() const
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
FiniteElementData::n_dofs_per_cell
unsigned int n_dofs_per_cell() const
DoFHandler::get_fe
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
BlockInfo::bi_global
BlockIndices bi_global
The block structure of the global system.
Definition: block_info.h:217