Reference documentation for deal.II version 9.0.0
block_info.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2017 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 #include <deal.II/dofs/block_info.h>
17 #include <deal.II/dofs/dof_handler.h>
18 #include <deal.II/dofs/dof_tools.h>
19 #include <deal.II/multigrid/mg_tools.h>
20 #include <deal.II/fe/fe.h>
21 #include <deal.II/fe/fe_tools.h>
22 
23 DEAL_II_NAMESPACE_OPEN
24 
25 
26 template <int dim, int spacedim>
27 void
28 BlockInfo::initialize(const DoFHandler<dim, spacedim> &dof, bool levels_only, bool active_only)
29 {
30  if (!levels_only && dof.has_active_dofs())
31  {
32  const FiniteElement<dim, spacedim> &fe = dof.get_fe();
33  std::vector<types::global_dof_index> sizes(fe.n_blocks());
35  bi_global.reinit(sizes);
36  }
37 
38  if (!active_only && dof.has_level_dofs())
39  {
40  std::vector<std::vector<types::global_dof_index> > sizes (dof.get_triangulation().n_levels ());
41 
42  for (unsigned int i = 0; i < sizes.size (); ++i)
43  sizes[i].resize (dof.get_fe ().n_blocks ());
44 
45  MGTools::count_dofs_per_block (dof, sizes);
46  levels.resize (sizes.size ());
47 
48  for (unsigned int i = 0; i < sizes.size (); ++i)
49  levels[i].reinit (sizes[i]);
50  }
51 }
52 
53 
54 template <int dim, int spacedim>
55 void
57 {
58  const FiniteElement<dim, spacedim> &fe = dof.get_fe();
59  std::vector<types::global_dof_index> sizes(fe.n_blocks());
60 
61  base_elements.resize(fe.n_blocks());
62 
63  for (unsigned int i=0; i<base_elements.size(); ++i)
64  base_elements[i] = fe.block_to_base_index(i).first;
65 
69  sizes, false);
70  bi_local.reinit(sizes);
71 }
72 
73 
74 // explicit instantiations
75 #include "block_info.inst"
76 
77 DEAL_II_NAMESPACE_CLOSE
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:200
unsigned int n_blocks() const
bool has_active_dofs() const
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
std::pair< unsigned int, unsigned int > block_to_base_index(const unsigned int block) const
Definition: fe.h:3166
BlockIndices bi_global
The block structure of the global system.
Definition: block_info.h:180
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:28
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:56
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:1836
std::vector< BlockIndices > levels
The multilevel block structure.
Definition: block_info.h:184
unsigned int n_dofs_per_cell() const
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=std::vector< unsigned int >())
Definition: mg_tools.cc:1125
const Triangulation< dim, spacedim > & get_triangulation() const
BlockIndices bi_local
The block structure of the cell systems.
Definition: block_info.h:189
bool has_level_dofs() const
std::vector< unsigned int > base_elements
Definition: block_info.h:194