Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+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\}}\)
Loading...
Searching...
No Matches
block_info.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2009 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
18
19#include <deal.II/fe/fe.h>
20#include <deal.II/fe/fe_tools.h>
21
23
25
26
27template <int dim, int spacedim>
28void
30 bool levels_only,
31 bool active_only)
32{
33 Assert(dof.has_hp_capabilities() == false,
35
36 if (!levels_only && dof.has_active_dofs())
37 {
38 const std::vector<types::global_dof_index> sizes =
40 bi_global.reinit(sizes);
41 }
42
43 if (!active_only && dof.has_level_dofs())
44 {
45 std::vector<std::vector<types::global_dof_index>> sizes(
47 std::vector<types::global_dof_index>(dof.get_fe().n_blocks()));
48
50 levels.resize(sizes.size());
51
52 for (unsigned int i = 0; i < sizes.size(); ++i)
53 levels[i].reinit(sizes[i]);
54 }
55}
56
57
58
59template <int dim, int spacedim>
60void
62{
63 Assert(dof.has_hp_capabilities() == false,
65
66 const FiniteElement<dim, spacedim> &fe = dof.get_fe();
67 std::vector<types::global_dof_index> sizes(fe.n_blocks());
68
69 base_elements.resize(fe.n_blocks());
70
71 for (unsigned int i = 0; i < base_elements.size(); ++i)
72 base_elements[i] = fe.block_to_base_index(i).first;
73
76 bi_local.reinit(sizes);
77}
78
79// explicit instantiations
80#include "block_info.inst"
81
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:213
std::vector< BlockIndices > levels
The multilevel block structure.
Definition block_info.h:197
BlockIndices bi_global
The block structure of the global system.
Definition block_info.h:193
BlockIndices bi_local
The block structure of the cell systems.
Definition block_info.h:202
std::vector< unsigned int > base_elements
Definition block_info.h:207
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:61
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:29
bool has_active_dofs() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
bool has_hp_capabilities() 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:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define Assert(cond, exc)
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 >())
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:1156