Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
dof_tools.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1999 - 2023 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
32
33#include <map>
34#include <ostream>
35#include <set>
36#include <vector>
37
38
40
41// Forward declarations
42#ifndef DOXYGEN
43class BlockMask;
44template <int dim, typename RangeNumberType>
45class Function;
46template <int dim, int spacedim>
47class FiniteElement;
48namespace hp
49{
50 template <int dim, int spacedim>
51 class MappingCollection;
52 template <int dim, int spacedim>
53 class FECollection;
54} // namespace hp
55template <class MeshType>
57class InterGridMap;
58template <int dim, int spacedim>
59class Mapping;
60template <int dim, class T>
61class Table;
62template <typename Number>
63class Vector;
64
65namespace GridTools
66{
67 template <typename CellIterator>
68 struct PeriodicFacePair;
69}
70
71namespace DoFTools
72{
73 namespace internal
74 {
75 /*
76 * Default value of the face_has_flux_coupling parameter of
77 * make_flux_sparsity_pattern. Defined here (instead of using a default
78 * lambda in the parameter list) to avoid a bug in gcc where the same lambda
79 * gets defined multiple times.
80 */
81 template <int dim, int spacedim>
82 inline bool
83 always_couple_on_faces(
85 const unsigned int)
86 {
87 return true;
88 }
89 } // namespace internal
90} // namespace DoFTools
91
92#endif
93
212namespace DoFTools
213{
226 {
242 nonzero
243 };
244
257 template <int dim, int spacedim>
258 void
260 const Table<2, Coupling> &table_by_component,
261 std::vector<Table<2, Coupling>> &tables_by_block);
262
268 template <int dim, int spacedim>
272 const Table<2, Coupling> & component_couplings);
273
281 template <int dim, int spacedim>
282 std::vector<Table<2, Coupling>>
285 const Table<2, Coupling> & component_couplings);
409 template <int dim, int spacedim, typename number = double>
410 void
412 const DoFHandler<dim, spacedim> &dof_handler,
413 SparsityPatternBase & sparsity_pattern,
415 const bool keep_constrained_dofs = true,
417
483 template <int dim, int spacedim, typename number = double>
484 void
486 const DoFHandler<dim, spacedim> &dof_handler,
487 const Table<2, Coupling> & coupling,
488 SparsityPatternBase & sparsity_pattern,
490 const bool keep_constrained_dofs = true,
492
513 template <int dim, int spacedim>
514 void
516 const DoFHandler<dim, spacedim> &dof_col,
517 SparsityPatternBase & sparsity);
518
564 template <int dim, int spacedim>
565 void
567 SparsityPatternBase & sparsity_pattern);
568
577 template <int dim, int spacedim, typename number>
578 void
580 const DoFHandler<dim, spacedim> &dof_handler,
581 SparsityPatternBase & sparsity_pattern,
582 const AffineConstraints<number> &constraints,
583 const bool keep_constrained_dofs = true,
585
586
606 template <int dim, int spacedim>
607 void
609 const DoFHandler<dim, spacedim> &dof,
610 SparsityPatternBase & sparsity,
611 const Table<2, Coupling> & cell_integrals_mask,
612 const Table<2, Coupling> & face_integrals_mask,
614
615
642 template <int dim, int spacedim, typename number>
643 void
645 const DoFHandler<dim, spacedim> &dof,
646 SparsityPatternBase & sparsity,
647 const AffineConstraints<number> &constraints,
648 const bool keep_constrained_dofs,
649 const Table<2, Coupling> & couplings,
650 const Table<2, Coupling> & face_couplings,
651 const types::subdomain_id subdomain_id,
652 const std::function<
654 const unsigned int)> &face_has_flux_coupling =
655 &internal::always_couple_on_faces<dim, spacedim>);
656
666 template <int dim, int spacedim>
667 void
669 const DoFHandler<dim, spacedim> & dof,
670 const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
671 SparsityPatternBase & sparsity_pattern);
672
690 template <int dim, int spacedim, typename number>
691 void
693 const DoFHandler<dim, spacedim> &dof,
694 const std::map<types::boundary_id, const Function<spacedim, number> *>
695 & boundary_ids,
696 const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
697 SparsityPatternBase & sparsity);
698
746 template <int dim, int spacedim, typename number>
747 void
749 AffineConstraints<number> & constraints);
750
818 template <int dim, int spacedim>
819 void
821 const DoFHandler<dim, spacedim> & coarse_grid,
822 const unsigned int coarse_component,
823 const DoFHandler<dim, spacedim> & fine_grid,
824 const unsigned int fine_component,
825 const InterGridMap<DoFHandler<dim, spacedim>> &coarse_to_fine_grid_map,
826 AffineConstraints<double> & constraints);
827
828
845 template <int dim, int spacedim>
846 void
848 const DoFHandler<dim, spacedim> & coarse_grid,
849 const unsigned int coarse_component,
850 const DoFHandler<dim, spacedim> & fine_grid,
851 const unsigned int fine_component,
852 const InterGridMap<DoFHandler<dim, spacedim>> &coarse_to_fine_grid_map,
853 std::vector<std::map<types::global_dof_index, float>>
854 &transfer_representation);
855
1030 template <typename FaceIterator, typename number>
1031 void
1033 const FaceIterator & face_1,
1035 AffineConstraints<number> & constraints,
1036 const ComponentMask & component_mask = ComponentMask(),
1037 const bool face_orientation = true,
1038 const bool face_flip = false,
1039 const bool face_rotation = false,
1040 const FullMatrix<double> & matrix = FullMatrix<double>(),
1041 const std::vector<unsigned int> &first_vector_components =
1042 std::vector<unsigned int>(),
1043 const number periodicity_factor = 1.);
1044
1045
1065 template <int dim, int spacedim, typename number>
1066 void
1068 const std::vector<GridTools::PeriodicFacePair<
1069 typename DoFHandler<dim, spacedim>::cell_iterator>> &periodic_faces,
1070 AffineConstraints<number> & constraints,
1071 const ComponentMask & component_mask = ComponentMask(),
1072 const std::vector<unsigned int> &first_vector_components =
1073 std::vector<unsigned int>(),
1074 const number periodicity_factor = 1.);
1075
1076
1077
1107 template <int dim, int spacedim, typename number>
1108 void
1110 const DoFHandler<dim, spacedim> &dof_handler,
1111 const types::boundary_id b_id1,
1112 const types::boundary_id b_id2,
1113 const unsigned int direction,
1114 AffineConstraints<number> & constraints,
1115 const ComponentMask & component_mask = ComponentMask(),
1116 const number periodicity_factor = 1.);
1117
1118
1119
1145 template <int dim, int spacedim, typename number>
1146 void
1148 const DoFHandler<dim, spacedim> &dof_handler,
1149 const types::boundary_id b_id,
1150 const unsigned int direction,
1151 AffineConstraints<number> & constraints,
1152 const ComponentMask & component_mask = ComponentMask(),
1153 const number periodicity_factor = 1.);
1154
1172 template <int dim, int spacedim>
1173 IndexSet
1175
1208 template <int dim, int spacedim>
1209 IndexSet
1210 extract_dofs(const DoFHandler<dim, spacedim> &dof_handler,
1211 const ComponentMask & component_mask);
1212
1236 template <int dim, int spacedim>
1237 IndexSet
1238 extract_dofs(const DoFHandler<dim, spacedim> &dof_handler,
1239 const BlockMask & block_mask);
1240
1245 template <int dim, int spacedim>
1246 void
1247 extract_level_dofs(const unsigned int level,
1248 const DoFHandler<dim, spacedim> &dof,
1249 const ComponentMask & component_mask,
1250 std::vector<bool> & selected_dofs);
1251
1256 template <int dim, int spacedim>
1257 void
1258 extract_level_dofs(const unsigned int level,
1259 const DoFHandler<dim, spacedim> &dof,
1260 const BlockMask & component_mask,
1261 std::vector<bool> & selected_dofs);
1262
1317 template <int dim, int spacedim>
1320 const ComponentMask & component_mask,
1321 std::vector<bool> & selected_dofs,
1322 const std::set<types::boundary_id> &boundary_ids = {});
1323
1366 template <int dim, int spacedim>
1367 IndexSet
1369 const ComponentMask &component_mask = ComponentMask(),
1370 const std::set<types::boundary_id> &boundary_ids = {});
1371
1378 template <int dim, int spacedim>
1381 const ComponentMask & component_mask,
1382 IndexSet & selected_dofs,
1383 const std::set<types::boundary_id> &boundary_ids = {});
1384
1401 template <int dim, int spacedim>
1402 void
1404 const DoFHandler<dim, spacedim> & dof_handler,
1405 const ComponentMask & component_mask,
1406 std::vector<bool> & selected_dofs,
1407 const std::set<types::boundary_id> &boundary_ids =
1408 std::set<types::boundary_id>());
1409
1437 template <int dim, int spacedim, typename number = double>
1438 IndexSet
1440 const DoFHandler<dim, spacedim> &dof_handler,
1441 const std::function<
1443 & predicate,
1445
1481 template <int dim, int spacedim>
1482 void
1484 const ComponentMask & component_mask,
1485 std::vector<std::vector<bool>> & constant_modes);
1501 template <int dim, int spacedim>
1502 void
1504 const types::subdomain_id subdomain_id,
1505 std::vector<bool> & selected_dofs);
1506
1521 template <int dim, int spacedim>
1522 IndexSet
1524
1541 template <int dim, int spacedim>
1542 void
1544 IndexSet & dof_set);
1545
1552 template <int dim, int spacedim>
1553 IndexSet
1555 const DoFHandler<dim, spacedim> &dof_handler,
1556 const unsigned int level);
1557
1566 template <int dim, int spacedim>
1567 void
1569 const DoFHandler<dim, spacedim> &dof_handler,
1570 IndexSet & dof_set,
1571 const unsigned int level);
1572
1582 template <int dim, int spacedim>
1583 IndexSet
1585
1596 template <int dim, int spacedim>
1597 void
1599 IndexSet & dof_set);
1600
1610 template <int dim, int spacedim>
1611 std::vector<IndexSet>
1613 const DoFHandler<dim, spacedim> &dof_handler,
1614 const ComponentMask & components = ComponentMask());
1615
1629 template <int dim, int spacedim>
1630 std::vector<IndexSet>
1632 const DoFHandler<dim, spacedim> &dof_handler);
1633
1648 template <int dim, int spacedim>
1649 std::vector<IndexSet>
1651 const DoFHandler<dim, spacedim> &dof_handler);
1652
1657 template <int dim, int spacedim>
1658 IndexSet
1660 const DoFHandler<dim, spacedim> &dof_handler,
1661 const unsigned int level);
1662
1669 template <int dim, int spacedim>
1670 void
1672 const DoFHandler<dim, spacedim> &dof_handler,
1673 const unsigned int level,
1674 IndexSet & dof_set);
1675
1676
1707 template <int dim, int spacedim>
1708 void
1710 std::vector<types::subdomain_id> &subdomain);
1711
1737 template <int dim, int spacedim>
1738 unsigned int
1740 const DoFHandler<dim, spacedim> &dof_handler,
1741 const types::subdomain_id subdomain);
1742
1762 template <int dim, int spacedim>
1763 void
1765 const DoFHandler<dim, spacedim> &dof_handler,
1766 const types::subdomain_id subdomain,
1767 std::vector<unsigned int> & n_dofs_on_subdomain);
1768
1790 template <int dim, int spacedim>
1791 IndexSet
1793 const DoFHandler<dim, spacedim> &dof_handler,
1794 const types::subdomain_id subdomain);
1851 template <int dim, int spacedim>
1852 std::vector<types::global_dof_index>
1854 const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
1855 &patch);
1856
1877 template <int dim, int spacedim>
1878 void
1880 const DoFHandler<dim, spacedim> &dof_handler,
1881 const unsigned int level,
1882 const std::vector<bool> & selected_dofs = {},
1883 const types::global_dof_index offset = 0);
1884
1936 template <int dim, int spacedim>
1937 std::vector<unsigned int>
1939 const DoFHandler<dim, spacedim> &dof_handler,
1940 const unsigned int level,
1941 const bool interior_dofs_only,
1942 const bool boundary_patches = false,
1943 const bool level_boundary_patches = false,
1944 const bool single_cell_patches = false,
1945 const bool invert_vertex_mapping = false);
1946
1961 template <int dim, int spacedim>
1962 std::vector<unsigned int>
1964 const DoFHandler<dim, spacedim> &dof_handler,
1965 const unsigned int level,
1966 const BlockMask &exclude_boundary_dofs = BlockMask(),
1967 const bool boundary_patches = false,
1968 const bool level_boundary_patches = false,
1969 const bool single_cell_patches = false,
1970 const bool invert_vertex_mapping = false);
1971
2009 template <int dim, int spacedim>
2010 void
2012 const DoFHandler<dim, spacedim> &dof_handler,
2013 const unsigned int level,
2014 const bool interior_dofs_only,
2015 const bool boundary_dofs = false);
2016
2036 template <int dim, int spacedim>
2037 void
2039 const DoFHandler<dim, spacedim> &dof_handler,
2040 const unsigned int level,
2041 const bool interior_dofs_only = false);
2042
2085 template <int dim, int spacedim>
2086 std::vector<types::global_dof_index>
2088 const DoFHandler<dim, spacedim> &dof_handler,
2089 const bool vector_valued_once = false,
2090 const std::vector<unsigned int> &target_component = {});
2091
2108 template <int dim, int spacedim>
2109 std::vector<types::global_dof_index>
2111 const std::vector<unsigned int> &target_block =
2112 std::vector<unsigned int>());
2113
2129 template <int dim, int spacedim>
2130 DEAL_II_DEPRECATED_EARLY void
2132 std::vector<unsigned int> & active_fe_indices);
2133
2161 template <int dim, int spacedim>
2162 unsigned int
2164 const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
2165 &patch);
2166
2189 template <int dim, int spacedim>
2190 void
2192 std::vector<types::global_dof_index> &mapping);
2193
2204 template <int dim, int spacedim>
2205 void
2207 const std::set<types::boundary_id> &boundary_ids,
2208 std::vector<types::global_dof_index> &mapping);
2209
2238 template <int dim, int spacedim>
2239 void
2241 const DoFHandler<dim, spacedim> &dof_handler,
2242 std::vector<Point<spacedim>> & support_points,
2243 const ComponentMask &mask = ComponentMask());
2244
2248 template <int dim, int spacedim>
2249 void
2252 const DoFHandler<dim, spacedim> & dof_handler,
2253 std::vector<Point<spacedim>> & support_points,
2254 const ComponentMask & mask = ComponentMask());
2255
2285 template <int dim, int spacedim>
2286 std::map<types::global_dof_index, Point<spacedim>>
2288 const DoFHandler<dim, spacedim> &dof_handler,
2289 const ComponentMask &mask = ComponentMask());
2290
2294 template <int dim, int spacedim>
2295 std::map<types::global_dof_index, Point<spacedim>>
2298 const DoFHandler<dim, spacedim> & dof_handler,
2299 const ComponentMask & mask = ComponentMask());
2300
2306 template <int dim, int spacedim>
2307 DEAL_II_DEPRECATED_EARLY void
2309 const Mapping<dim, spacedim> & mapping,
2310 const DoFHandler<dim, spacedim> & dof_handler,
2311 std::map<types::global_dof_index, Point<spacedim>> &support_points,
2312 const ComponentMask & mask = ComponentMask());
2313
2319 template <int dim, int spacedim>
2320 DEAL_II_DEPRECATED_EARLY void
2323 const DoFHandler<dim, spacedim> & dof_handler,
2324 std::map<types::global_dof_index, Point<spacedim>> &support_points,
2325 const ComponentMask & mask = ComponentMask());
2326
2327
2348 template <int dim, int spacedim, class Comp>
2349 void
2351 const Mapping<dim, spacedim> & mapping,
2352 const DoFHandler<dim, spacedim> &dof_handler,
2354 &point_to_index_map);
2392 template <int dim, int spacedim, typename Number>
2393 void
2395 const Vector<Number> & cell_data,
2396 Vector<double> & dof_data,
2397 const unsigned int component = 0);
2398
2399
2477 template <int spacedim>
2478 void
2480 std::ostream & out,
2481 const std::map<types::global_dof_index, Point<spacedim>> &support_points);
2482
2483
2524 template <int dim, int spacedim, typename number>
2525 void
2527 const DoFHandler<dim, spacedim> &dof,
2528 const types::boundary_id boundary_id,
2529 AffineConstraints<number> & zero_boundary_constraints,
2530 const ComponentMask & component_mask = ComponentMask());
2531
2542 template <int dim, int spacedim, typename number>
2543 void
2545 const DoFHandler<dim, spacedim> &dof,
2546 AffineConstraints<number> & zero_boundary_constraints,
2547 const ComponentMask & component_mask = ComponentMask());
2548
2593} // namespace DoFTools
2594
2595
2596
2597/* ------------------------- inline functions -------------- */
2598
2599#ifndef DOXYGEN
2600
2601namespace DoFTools
2602{
2608 inline Coupling
2609 operator|=(Coupling &c1, const Coupling c2)
2610 {
2611 if (c2 == always)
2612 c1 = always;
2613 else if (c1 != always && c2 == nonzero)
2614 return c1 = nonzero;
2615 return c1;
2616 }
2617
2618
2624 inline Coupling
2625 operator|(const Coupling c1, const Coupling c2)
2626 {
2627 if (c1 == always || c2 == always)
2628 return always;
2629 if (c1 == nonzero || c2 == nonzero)
2630 return nonzero;
2631 return none;
2632 }
2633
2634
2635 // ---------------------- inline and template functions --------------------
2636
2637 template <int dim, int spacedim, class Comp>
2638 void
2640 const Mapping<dim, spacedim> & mapping,
2641 const DoFHandler<dim, spacedim> &dof_handler,
2643 &point_to_index_map)
2644 {
2645 // let the checking of arguments be
2646 // done by the function first
2647 // called
2648 std::vector<Point<spacedim>> support_points(dof_handler.n_dofs());
2649 map_dofs_to_support_points(mapping, dof_handler, support_points);
2650 // now copy over the results of the
2651 // previous function into the
2652 // output arg
2653 point_to_index_map.clear();
2654 for (types::global_dof_index i = 0; i < dof_handler.n_dofs(); ++i)
2655 point_to_index_map[support_points[i]] = i;
2656 }
2657} // namespace DoFTools
2658
2659#endif
2660
2662
2663#endif
types::global_dof_index n_dofs() const
Abstract base class for mapping classes.
Definition mapping.h:317
Definition point.h:112
#define DEAL_II_DEPRECATED
Definition config.h:172
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:160
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
unsigned int level
Definition grid_out.cc:4618
#define DeclException0(Exception0)
Definition exceptions.h:465
static ::ExceptionBase & ExcGridNotCoarser()
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
static ::ExceptionBase & ExcFiniteElementsDontMatch()
static ::ExceptionBase & ExcNoFESelected()
static ::ExceptionBase & ExcGridsDontMatch()
typename ActiveSelector::cell_iterator cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
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, SparsityPatternBase &sparsity_pattern)
void make_boundary_sparsity_pattern(const DoFHandler< dim, spacedim > &dof, const std::vector< types::global_dof_index > &dof_to_boundary_mapping, SparsityPatternBase &sparsity_pattern)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
void get_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, std::vector< types::subdomain_id > &subdomain)
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:545
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:796
IndexSet dof_indices_with_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain)
std::vector< IndexSet > locally_owned_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
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())
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
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)
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)
IndexSet extract_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask)
Definition dof_tools.cc:386
IndexSet extract_locally_active_dofs(const DoFHandler< dim, spacedim > &dof_handler)
IndexSet extract_locally_relevant_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
std::vector< IndexSet > locally_owned_dofs_per_component(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &components=ComponentMask())
Definition dof_tools.cc:439
void write_gnuplot_dof_support_point_info(std::ostream &out, const std::map< types::global_dof_index, Point< spacedim > > &support_points)
void extract_subdomain_dofs(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain_id, std::vector< bool > &selected_dofs)
IndexSet extract_locally_active_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
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)
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 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:473
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:714
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:302
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={})
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)
void map_dof_to_boundary_indices(const DoFHandler< dim, spacedim > &dof_handler, std::vector< types::global_dof_index > &mapping)
std::vector< IndexSet > locally_relevant_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
void make_single_patch(SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const bool interior_dofs_only=false)
void make_periodicity_constraints(const FaceIterator &face_1, const std_cxx20::type_identity_t< FaceIterator > &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.)
std::vector< types::global_dof_index > get_dofs_on_patch(const std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > &patch)
unsigned int count_dofs_with_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain)
unsigned int count_dofs_on_patch(const std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > &patch)
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)
void extract_constant_modes(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< std::vector< bool > > &constant_modes)
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 hp.h:118
const types::subdomain_id invalid_subdomain_id
Definition types.h:298
typename type_identity< T >::type type_identity_t
Definition type_traits.h:96
unsigned int global_dof_index
Definition types.h:82
unsigned int boundary_id
Definition types.h:141
ParameterHandler::OutputStyle operator|(const ParameterHandler::OutputStyle f1, const ParameterHandler::OutputStyle f2)