Reference documentation for deal.II version 9.4.1
\(\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
dof_tools.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1999 - 2022 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#ifndef dealii_dof_tools_h
17#define dealii_dof_tools_h
18
19
20#include <deal.II/base/config.h>
21
24#include <deal.II/base/point.h>
25
27
29
31
33
34#include <map>
35#include <ostream>
36#include <set>
37#include <vector>
38
39
41
42// Forward declarations
43#ifndef DOXYGEN
44class BlockMask;
45template <int dim, typename RangeNumberType>
46class Function;
47template <int dim, int spacedim>
48class FiniteElement;
49namespace hp
50{
51 template <int dim, int spacedim>
52 class MappingCollection;
53 template <int dim, int spacedim>
54 class FECollection;
55} // namespace hp
56template <class MeshType>
57class InterGridMap;
58template <int dim, int spacedim>
59class Mapping;
60class SparsityPattern;
61template <int dim, class T>
62class Table;
63template <typename Number>
64class Vector;
65
66namespace GridTools
67{
68 template <typename CellIterator>
69 struct PeriodicFacePair;
70}
71
72namespace DoFTools
73{
74 namespace internal
75 {
76 /*
77 * Default value of the face_has_flux_coupling parameter of
78 * make_flux_sparsity_pattern. Defined here (instead of using a default
79 * lambda in the parameter list) to avoid a bug in gcc where the same lambda
80 * gets defined multiple times.
81 */
82 template <int dim, int spacedim>
83 inline bool
84 always_couple_on_faces(
86 const unsigned int)
87 {
88 return true;
89 }
90 } // namespace internal
91} // namespace DoFTools
92
93#endif
94
213namespace DoFTools
214{
227 {
243 nonzero
244 };
245
258 template <int dim, int spacedim>
259 void
261 const Table<2, Coupling> &table_by_component,
262 std::vector<Table<2, Coupling>> &tables_by_block);
263
269 template <int dim, int spacedim>
273 const Table<2, Coupling> & component_couplings);
274
282 template <int dim, int spacedim>
283 std::vector<Table<2, Coupling>>
286 const Table<2, Coupling> & component_couplings);
412 template <int dim,
413 int spacedim,
414 typename SparsityPatternType,
415 typename number = double>
416 void
418 const DoFHandler<dim, spacedim> &dof_handler,
419 SparsityPatternType & sparsity_pattern,
421 const bool keep_constrained_dofs = true,
423
489 template <int dim,
490 int spacedim,
491 typename SparsityPatternType,
492 typename number = double>
493 void
495 const DoFHandler<dim, spacedim> &dof_handler,
496 const Table<2, Coupling> & coupling,
497 SparsityPatternType & sparsity_pattern,
499 const bool keep_constrained_dofs = true,
501
522 template <int dim, int spacedim, typename SparsityPatternType>
523 void
525 const DoFHandler<dim, spacedim> &dof_col,
526 SparsityPatternType & sparsity);
527
573 template <int dim, int spacedim, typename SparsityPatternType>
574 void
576 SparsityPatternType & sparsity_pattern);
577
586 template <int dim,
587 int spacedim,
588 typename SparsityPatternType,
589 typename number>
590 void
592 const DoFHandler<dim, spacedim> &dof_handler,
593 SparsityPatternType & sparsity_pattern,
594 const AffineConstraints<number> &constraints,
595 const bool keep_constrained_dofs = true,
597
598
618 template <int dim, int spacedim, typename SparsityPatternType>
619 void
621 const DoFHandler<dim, spacedim> &dof,
622 SparsityPatternType & sparsity,
623 const Table<2, Coupling> & cell_integrals_mask,
624 const Table<2, Coupling> & face_integrals_mask,
626
627
654 template <int dim,
655 int spacedim,
656 typename SparsityPatternType,
657 typename number>
658 void
660 const DoFHandler<dim, spacedim> &dof,
661 SparsityPatternType & sparsity,
662 const AffineConstraints<number> &constraints,
663 const bool keep_constrained_dofs,
664 const Table<2, Coupling> & couplings,
665 const Table<2, Coupling> & face_couplings,
666 const types::subdomain_id subdomain_id,
667 const std::function<
669 const unsigned int)> &face_has_flux_coupling =
670 &internal::always_couple_on_faces<dim, spacedim>);
671
681 template <int dim, int spacedim, typename SparsityPatternType>
682 void
684 const DoFHandler<dim, spacedim> & dof,
685 const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
686 SparsityPatternType & sparsity_pattern);
687
705 template <int dim,
706 int spacedim,
707 typename SparsityPatternType,
708 typename number>
709 void
711 const DoFHandler<dim, spacedim> &dof,
712 const std::map<types::boundary_id, const Function<spacedim, number> *>
713 & boundary_ids,
714 const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
715 SparsityPatternType & sparsity);
716
764 template <int dim, int spacedim, typename number>
765 void
767 AffineConstraints<number> & constraints);
768
836 template <int dim, int spacedim>
837 void
839 const DoFHandler<dim, spacedim> & coarse_grid,
840 const unsigned int coarse_component,
841 const DoFHandler<dim, spacedim> & fine_grid,
842 const unsigned int fine_component,
843 const InterGridMap<DoFHandler<dim, spacedim>> &coarse_to_fine_grid_map,
844 AffineConstraints<double> & constraints);
845
846
863 template <int dim, int spacedim>
864 void
866 const DoFHandler<dim, spacedim> & coarse_grid,
867 const unsigned int coarse_component,
868 const DoFHandler<dim, spacedim> & fine_grid,
869 const unsigned int fine_component,
870 const InterGridMap<DoFHandler<dim, spacedim>> &coarse_to_fine_grid_map,
871 std::vector<std::map<types::global_dof_index, float>>
872 &transfer_representation);
873
1048 template <typename FaceIterator, typename number>
1049 void
1051 const FaceIterator & face_1,
1052 const typename identity<FaceIterator>::type &face_2,
1053 AffineConstraints<number> & constraints,
1054 const ComponentMask & component_mask = ComponentMask(),
1055 const bool face_orientation = true,
1056 const bool face_flip = false,
1057 const bool face_rotation = false,
1058 const FullMatrix<double> & matrix = FullMatrix<double>(),
1059 const std::vector<unsigned int> &first_vector_components =
1060 std::vector<unsigned int>(),
1061 const number periodicity_factor = 1.);
1062
1063
1064
1084 template <int dim, int spacedim, typename number>
1085 void
1087 const std::vector<GridTools::PeriodicFacePair<
1088 typename DoFHandler<dim, spacedim>::cell_iterator>> &periodic_faces,
1089 AffineConstraints<number> & constraints,
1090 const ComponentMask & component_mask = ComponentMask(),
1091 const std::vector<unsigned int> &first_vector_components =
1092 std::vector<unsigned int>(),
1093 const number periodicity_factor = 1.);
1094
1101 template <typename DoFHandlerType, typename number>
1104 const std::vector<
1106 & periodic_faces,
1107 AffineConstraints<number> & constraints,
1108 const ComponentMask & component_mask = ComponentMask(),
1109 const std::vector<unsigned int> &first_vector_components =
1110 std::vector<unsigned int>(),
1111 const number periodicity_factor = 1.);
1112
1113
1114
1144 template <int dim, int spacedim, typename number>
1145 void
1147 const DoFHandler<dim, spacedim> &dof_handler,
1148 const types::boundary_id b_id1,
1149 const types::boundary_id b_id2,
1150 const unsigned int direction,
1151 AffineConstraints<number> & constraints,
1152 const ComponentMask & component_mask = ComponentMask(),
1153 const number periodicity_factor = 1.);
1154
1155
1156
1182 template <int dim, int spacedim, typename number>
1183 void
1185 const DoFHandler<dim, spacedim> &dof_handler,
1186 const types::boundary_id b_id,
1187 const unsigned int direction,
1188 AffineConstraints<number> & constraints,
1189 const ComponentMask & component_mask = ComponentMask(),
1190 const number periodicity_factor = 1.);
1191
1209 template <int dim, int spacedim>
1210 IndexSet
1212
1245 template <int dim, int spacedim>
1246 IndexSet
1247 extract_dofs(const DoFHandler<dim, spacedim> &dof_handler,
1248 const ComponentMask & component_mask);
1249
1273 template <int dim, int spacedim>
1274 IndexSet
1275 extract_dofs(const DoFHandler<dim, spacedim> &dof_handler,
1276 const BlockMask & block_mask);
1277
1282 template <int dim, int spacedim>
1283 void
1284 extract_level_dofs(const unsigned int level,
1285 const DoFHandler<dim, spacedim> &dof,
1286 const ComponentMask & component_mask,
1287 std::vector<bool> & selected_dofs);
1288
1293 template <int dim, int spacedim>
1294 void
1295 extract_level_dofs(const unsigned int level,
1296 const DoFHandler<dim, spacedim> &dof,
1297 const BlockMask & component_mask,
1298 std::vector<bool> & selected_dofs);
1299
1354 template <int dim, int spacedim>
1357 const ComponentMask & component_mask,
1358 std::vector<bool> & selected_dofs,
1359 const std::set<types::boundary_id> &boundary_ids = {});
1360
1403 template <int dim, int spacedim>
1404 IndexSet
1406 const ComponentMask &component_mask = ComponentMask(),
1407 const std::set<types::boundary_id> &boundary_ids = {});
1408
1415 template <int dim, int spacedim>
1418 const ComponentMask & component_mask,
1419 IndexSet & selected_dofs,
1420 const std::set<types::boundary_id> &boundary_ids = {});
1421
1438 template <int dim, int spacedim>
1439 void
1441 const DoFHandler<dim, spacedim> & dof_handler,
1442 const ComponentMask & component_mask,
1443 std::vector<bool> & selected_dofs,
1444 const std::set<types::boundary_id> &boundary_ids =
1445 std::set<types::boundary_id>());
1446
1474 template <int dim, int spacedim, typename number = double>
1475 IndexSet
1477 const DoFHandler<dim, spacedim> &dof_handler,
1478 const std::function<
1480 & predicate,
1482
1518 template <int dim, int spacedim>
1519 void
1521 const ComponentMask & component_mask,
1522 std::vector<std::vector<bool>> & constant_modes);
1524
1538 template <int dim, int spacedim>
1539 void
1541 const types::subdomain_id subdomain_id,
1542 std::vector<bool> & selected_dofs);
1543
1558 template <int dim, int spacedim>
1559 IndexSet
1561
1578 template <int dim, int spacedim>
1579 void
1581 IndexSet & dof_set);
1582
1589 template <int dim, int spacedim>
1590 IndexSet
1592 const DoFHandler<dim, spacedim> &dof_handler,
1593 const unsigned int level);
1594
1603 template <int dim, int spacedim>
1604 void
1606 const DoFHandler<dim, spacedim> &dof_handler,
1607 IndexSet & dof_set,
1608 const unsigned int level);
1609
1619 template <int dim, int spacedim>
1620 IndexSet
1622
1634 template <int dim, int spacedim>
1635 void
1637 IndexSet & dof_set);
1638
1648 template <int dim, int spacedim>
1649 std::vector<IndexSet>
1651 const DoFHandler<dim, spacedim> &dof_handler,
1652 const ComponentMask & components = ComponentMask());
1653
1667 template <int dim, int spacedim>
1668 std::vector<IndexSet>
1670 const DoFHandler<dim, spacedim> &dof_handler);
1671
1686 template <int dim, int spacedim>
1687 std::vector<IndexSet>
1689 const DoFHandler<dim, spacedim> &dof_handler);
1690
1695 template <int dim, int spacedim>
1696 IndexSet
1698 const DoFHandler<dim, spacedim> &dof_handler,
1699 const unsigned int level);
1700
1707 template <int dim, int spacedim>
1708 void
1710 const DoFHandler<dim, spacedim> &dof_handler,
1711 const unsigned int level,
1712 IndexSet & dof_set);
1713
1714
1745 template <int dim, int spacedim>
1746 void
1748 std::vector<types::subdomain_id> &subdomain);
1749
1775 template <int dim, int spacedim>
1776 unsigned int
1778 const DoFHandler<dim, spacedim> &dof_handler,
1779 const types::subdomain_id subdomain);
1780
1800 template <int dim, int spacedim>
1801 void
1803 const DoFHandler<dim, spacedim> &dof_handler,
1804 const types::subdomain_id subdomain,
1805 std::vector<unsigned int> & n_dofs_on_subdomain);
1806
1828 template <int dim, int spacedim>
1829 IndexSet
1831 const DoFHandler<dim, spacedim> &dof_handler,
1832 const types::subdomain_id subdomain);
1833 // @}
1843
1888 template <int dim, int spacedim>
1889 std::vector<types::global_dof_index>
1891 const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
1892 &patch);
1893
1900 template <typename DoFHandlerType>
1901 DEAL_II_DEPRECATED std::vector<types::global_dof_index>
1903 const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
1904
1925 template <int dim, int spacedim>
1926 void
1928 const DoFHandler<dim, spacedim> &dof_handler,
1929 const unsigned int level,
1930 const std::vector<bool> & selected_dofs = {},
1931 const types::global_dof_index offset = 0);
1932
1984 template <int dim, int spacedim>
1985 std::vector<unsigned int>
1987 const DoFHandler<dim, spacedim> &dof_handler,
1988 const unsigned int level,
1989 const bool interior_dofs_only,
1990 const bool boundary_patches = false,
1991 const bool level_boundary_patches = false,
1992 const bool single_cell_patches = false,
1993 const bool invert_vertex_mapping = false);
1994
2009 template <int dim, int spacedim>
2010 std::vector<unsigned int>
2012 const DoFHandler<dim, spacedim> &dof_handler,
2013 const unsigned int level,
2014 const BlockMask &exclude_boundary_dofs = BlockMask(),
2015 const bool boundary_patches = false,
2016 const bool level_boundary_patches = false,
2017 const bool single_cell_patches = false,
2018 const bool invert_vertex_mapping = false);
2019
2057 template <int dim, int spacedim>
2058 void
2060 const DoFHandler<dim, spacedim> &dof_handler,
2061 const unsigned int level,
2062 const bool interior_dofs_only,
2063 const bool boundary_dofs = false);
2064
2084 template <int dim, int spacedim>
2085 void
2087 const DoFHandler<dim, spacedim> &dof_handler,
2088 const unsigned int level,
2089 const bool interior_dofs_only = false);
2090
2133 template <int dim, int spacedim>
2134 std::vector<types::global_dof_index>
2136 const DoFHandler<dim, spacedim> &dof_handler,
2137 const bool vector_valued_once = false,
2138 const std::vector<unsigned int> &target_component = {});
2139
2156 template <int dim, int spacedim>
2157 std::vector<types::global_dof_index>
2159 const std::vector<unsigned int> &target_block =
2160 std::vector<unsigned int>());
2161
2171 template <int dim, int spacedim>
2172 void
2174 std::vector<unsigned int> & active_fe_indices);
2175
2203 template <int dim, int spacedim>
2204 unsigned int
2206 const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
2207 &patch);
2208
2215 template <typename DoFHandlerType>
2216 DEAL_II_DEPRECATED unsigned int
2218 const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
2219
2242 template <int dim, int spacedim>
2243 void
2245 std::vector<types::global_dof_index> &mapping);
2246
2257 template <int dim, int spacedim>
2258 void
2260 const std::set<types::boundary_id> &boundary_ids,
2261 std::vector<types::global_dof_index> &mapping);
2262
2291 template <int dim, int spacedim>
2292 void
2294 const DoFHandler<dim, spacedim> &dof_handler,
2295 std::vector<Point<spacedim>> & support_points,
2296 const ComponentMask &mask = ComponentMask());
2297
2301 template <int dim, int spacedim>
2302 void
2304 const ::hp::MappingCollection<dim, spacedim> &mapping,
2305 const DoFHandler<dim, spacedim> & dof_handler,
2306 std::vector<Point<spacedim>> & support_points,
2307 const ComponentMask & mask = ComponentMask());
2308
2339 template <int dim, int spacedim>
2340 void
2342 const Mapping<dim, spacedim> & mapping,
2343 const DoFHandler<dim, spacedim> & dof_handler,
2344 std::map<types::global_dof_index, Point<spacedim>> &support_points,
2345 const ComponentMask & mask = ComponentMask());
2346
2350 template <int dim, int spacedim>
2351 void
2353 const ::hp::MappingCollection<dim, spacedim> &mapping,
2354 const DoFHandler<dim, spacedim> & dof_handler,
2355 std::map<types::global_dof_index, Point<spacedim>> &support_points,
2356 const ComponentMask & mask = ComponentMask());
2357
2358
2379 template <int dim, int spacedim, class Comp>
2380 void
2382 const Mapping<dim, spacedim> & mapping,
2383 const DoFHandler<dim, spacedim> &dof_handler,
2385 &point_to_index_map);
2423 template <int dim, int spacedim, typename Number>
2424 void
2426 const Vector<Number> & cell_data,
2427 Vector<double> & dof_data,
2428 const unsigned int component = 0);
2429
2430
2507 template <int spacedim>
2508 void
2510 std::ostream & out,
2511 const std::map<types::global_dof_index, Point<spacedim>> &support_points);
2512
2513
2554 template <int dim, int spacedim, typename number>
2555 void
2557 const DoFHandler<dim, spacedim> &dof,
2558 const types::boundary_id boundary_id,
2559 AffineConstraints<number> & zero_boundary_constraints,
2560 const ComponentMask & component_mask = ComponentMask());
2561
2572 template <int dim, int spacedim, typename number>
2573 void
2575 const DoFHandler<dim, spacedim> &dof,
2576 AffineConstraints<number> & zero_boundary_constraints,
2577 const ComponentMask & component_mask = ComponentMask());
2578
2623} // namespace DoFTools
2624
2625
2626
2627/* ------------------------- inline functions -------------- */
2628
2629#ifndef DOXYGEN
2630
2631namespace DoFTools
2632{
2638 inline Coupling
2639 operator|=(Coupling &c1, const Coupling c2)
2640 {
2641 if (c2 == always)
2642 c1 = always;
2643 else if (c1 != always && c2 == nonzero)
2644 return c1 = nonzero;
2645 return c1;
2646 }
2647
2648
2654 inline Coupling
2655 operator|(const Coupling c1, const Coupling c2)
2656 {
2657 if (c1 == always || c2 == always)
2658 return always;
2659 if (c1 == nonzero || c2 == nonzero)
2660 return nonzero;
2661 return none;
2662 }
2663
2664
2665 // ---------------------- inline and template functions --------------------
2666
2667 template <int dim, int spacedim, class Comp>
2668 void
2670 const Mapping<dim, spacedim> & mapping,
2671 const DoFHandler<dim, spacedim> &dof_handler,
2673 &point_to_index_map)
2674 {
2675 // let the checking of arguments be
2676 // done by the function first
2677 // called
2678 std::vector<Point<spacedim>> support_points(dof_handler.n_dofs());
2679 map_dofs_to_support_points(mapping, dof_handler, support_points);
2680 // now copy over the results of the
2681 // previous function into the
2682 // output arg
2683 point_to_index_map.clear();
2684 for (types::global_dof_index i = 0; i < dof_handler.n_dofs(); ++i)
2685 point_to_index_map[support_points[i]] = i;
2686 }
2687
2688
2689
2690 template <typename DoFHandlerType, typename number>
2691 inline void
2693 const std::vector<
2695 & periodic_faces,
2696 AffineConstraints<number> & constraints,
2697 const ComponentMask & component_mask,
2698 const std::vector<unsigned int> &first_vector_components,
2699 const number periodicity_factor)
2700 {
2701 make_periodicity_constraints<DoFHandlerType::dimension,
2702 DoFHandlerType::space_dimension>(
2703 periodic_faces,
2704 constraints,
2705 component_mask,
2706 first_vector_components,
2707 periodicity_factor);
2708 }
2709
2710
2711
2712 template <typename DoFHandlerType>
2713 inline std::vector<types::global_dof_index>
2715 const std::vector<typename DoFHandlerType::active_cell_iterator> &patch)
2716 {
2717 return get_dofs_on_patch<DoFHandlerType::dimension,
2718 DoFHandlerType::space_dimension>(patch);
2719 }
2720
2721
2722
2723 template <typename DoFHandlerType>
2724 inline unsigned int
2726 const std::vector<typename DoFHandlerType::active_cell_iterator> &patch)
2727 {
2728 return count_dofs_on_patch<DoFHandlerType::dimension,
2729 DoFHandlerType::space_dimension>(patch);
2730 }
2731} // namespace DoFTools
2732
2733#endif
2734
2736
2737#endif
types::global_dof_index n_dofs() const
Abstract base class for mapping classes.
Definition: mapping.h:311
Definition: point.h:111
Definition: vector.h:109
#define DEAL_II_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
unsigned int level
Definition: grid_out.cc:4606
#define DeclException0(Exception0)
Definition: exceptions.h:464
static ::ExceptionBase & ExcGridNotCoarser()
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
static ::ExceptionBase & ExcFiniteElementsDontMatch()
static ::ExceptionBase & ExcNoFESelected()
static ::ExceptionBase & ExcGridsDontMatch()
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:466
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:438
void compute_intergrid_transfer_representation(const DoFHandler< dim, spacedim > &coarse_grid, const unsigned int coarse_component, const DoFHandler< dim, spacedim > &fine_grid, const unsigned int fine_component, const InterGridMap< DoFHandler< dim, spacedim > > &coarse_to_fine_grid_map, std::vector< std::map< types::global_dof_index, float > > &transfer_representation)
void make_zero_boundary_constraints(const DoFHandler< dim, spacedim > &dof, const types::boundary_id boundary_id, AffineConstraints< number > &zero_boundary_constraints, const ComponentMask &component_mask=ComponentMask())
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void compute_intergrid_constraints(const DoFHandler< dim, spacedim > &coarse_grid, const unsigned int coarse_component, const DoFHandler< dim, spacedim > &fine_grid, const unsigned int fine_component, const InterGridMap< DoFHandler< dim, spacedim > > &coarse_to_fine_grid_map, AffineConstraints< double > &constraints)
void make_flux_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity_pattern)
void make_boundary_sparsity_pattern(const DoFHandler< dim, spacedim > &dof, const std::vector< types::global_dof_index > &dof_to_boundary_mapping, SparsityPatternType &sparsity_pattern)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
void get_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, std::vector< types::subdomain_id > &subdomain)
Definition: dof_tools.cc:1554
void extract_boundary_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< bool > &selected_dofs, const std::set< types::boundary_id > &boundary_ids={})
Definition: dof_tools.cc:549
IndexSet extract_dofs_with_support_contained_within(const DoFHandler< dim, spacedim > &dof_handler, const std::function< bool(const typename DoFHandler< dim, spacedim >::active_cell_iterator &)> &predicate, const AffineConstraints< number > &constraints=AffineConstraints< number >())
Definition: dof_tools.cc:800
IndexSet dof_indices_with_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1661
std::vector< IndexSet > locally_owned_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1383
void map_dofs_to_support_points(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, std::vector< Point< spacedim > > &support_points, const ComponentMask &mask=ComponentMask())
Definition: dof_tools.cc:2251
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1144
void map_support_points_to_dofs(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, std::map< Point< spacedim >, types::global_dof_index, Comp > &point_to_index_map)
void get_active_fe_indices(const DoFHandler< dim, spacedim > &dof_handler, std::vector< unsigned int > &active_fe_indices)
Definition: dof_tools.cc:1371
void make_cell_patches(SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const std::vector< bool > &selected_dofs={}, const types::global_dof_index offset=0)
Definition: dof_tools.cc:2433
IndexSet extract_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask)
Definition: dof_tools.cc:390
IndexSet extract_locally_active_dofs(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1052
IndexSet extract_locally_relevant_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
Definition: dof_tools.cc:1194
std::vector< IndexSet > locally_owned_dofs_per_component(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &components=ComponentMask())
Definition: dof_tools.cc:443
void write_gnuplot_dof_support_point_info(std::ostream &out, const std::map< types::global_dof_index, Point< spacedim > > &support_points)
Definition: dof_tools.cc:2341
void extract_subdomain_dofs(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain_id, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:1021
IndexSet extract_locally_active_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
Definition: dof_tools.cc:1095
std::vector< unsigned int > make_vertex_patches(SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const bool interior_dofs_only, const bool boundary_patches=false, const bool level_boundary_patches=false, const bool single_cell_patches=false, const bool invert_vertex_mapping=false)
Definition: dof_tools.cc:2581
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:1980
void extract_level_dofs(const unsigned int level, const DoFHandler< dim, spacedim > &dof, const ComponentMask &component_mask, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:477
void extract_dofs_with_support_on_boundary(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< bool > &selected_dofs, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
Definition: dof_tools.cc:718
void make_periodicity_constraints(const FaceIterator &face_1, const typename identity< FaceIterator >::type &face_2, AffineConstraints< number > &constraints, const ComponentMask &component_mask=ComponentMask(), const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const FullMatrix< double > &matrix=FullMatrix< double >(), const std::vector< unsigned int > &first_vector_components=std::vector< unsigned int >(), const number periodicity_factor=1.)
void distribute_cell_to_dof_vector(const DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &cell_data, Vector< double > &dof_data, const unsigned int component=0)
Definition: dof_tools.cc:306
std::vector< types::global_dof_index > count_dofs_per_fe_component(const DoFHandler< dim, spacedim > &dof_handler, const bool vector_valued_once=false, const std::vector< unsigned int > &target_component={})
Definition: dof_tools.cc:1888
void make_child_patches(SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const bool interior_dofs_only, const bool boundary_dofs=false)
Definition: dof_tools.cc:2518
void map_dof_to_boundary_indices(const DoFHandler< dim, spacedim > &dof_handler, std::vector< types::global_dof_index > &mapping)
Definition: dof_tools.cc:2069
std::vector< IndexSet > locally_relevant_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1478
void make_single_patch(SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const bool interior_dofs_only=false)
Definition: dof_tools.cc:2476
std::vector< types::global_dof_index > get_dofs_on_patch(const std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > &patch)
Definition: dof_tools.cc:2782
unsigned int count_dofs_with_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1644
unsigned int count_dofs_on_patch(const std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > &patch)
Definition: dof_tools.cc:2751
void convert_couplings_to_blocks(const DoFHandler< dim, spacedim > &dof_handler, const Table< 2, Coupling > &table_by_component, std::vector< Table< 2, Coupling > > &tables_by_block)
Definition: dof_tools.cc:2380
void extract_constant_modes(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< std::vector< bool > > &constant_modes)
Definition: dof_tools.cc:1255
Table< 2, Coupling > dof_couplings_from_component_couplings(const FiniteElement< dim, spacedim > &fe, const Table< 2, Coupling > &component_couplings)
IndexSet extract_hanging_node_dofs(const DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_tools.cc:1012
Definition: hp.h:118
const types::subdomain_id invalid_subdomain_id
Definition: types.h:281
unsigned int global_dof_index
Definition: types.h:76
unsigned int boundary_id
Definition: types.h:129
ParameterHandler::OutputStyle operator|(const ParameterHandler::OutputStyle f1, const ParameterHandler::OutputStyle f2)