deal.II version GIT relicensing-2017-g7ae82fad8b 2024-10-22 19:20:00+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
dof_tools.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1999 - 2024 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
15#ifndef dealii_dof_tools_h
16#define dealii_dof_tools_h
17
18
19#include <deal.II/base/config.h>
20
23#include <deal.II/base/point.h>
24
26
28
31
32#include <map>
33#include <ostream>
34#include <set>
35#include <vector>
36
37
39
40// Forward declarations
41#ifndef DOXYGEN
42class BlockMask;
43template <int dim, typename RangeNumberType>
44class Function;
45template <int dim, int spacedim>
46class FiniteElement;
47namespace hp
48{
49 template <int dim, int spacedim>
50 class MappingCollection;
51 template <int dim, int spacedim>
52 class FECollection;
53} // namespace hp
54template <typename MeshType>
56class InterGridMap;
57template <int dim, int spacedim>
58class Mapping;
59template <int dim, class T>
60class Table;
61template <typename Number>
62class Vector;
63
64namespace GridTools
65{
66 template <typename CellIterator>
67 struct PeriodicFacePair;
68}
69
70namespace DoFTools
71{
72 namespace internal
73 {
74 /*
75 * Default value of the face_has_flux_coupling parameter of
76 * make_flux_sparsity_pattern. Defined here (instead of using a default
77 * lambda in the parameter list) to avoid a bug in gcc where the same lambda
78 * gets defined multiple times.
79 */
80 template <int dim, int spacedim>
81 inline bool
82 always_couple_on_faces(
84 const unsigned int)
85 {
86 return true;
87 }
88 } // namespace internal
89} // namespace DoFTools
90
91#endif
92
225namespace DoFTools
226{
257
270 template <int dim, int spacedim>
271 void
273 const Table<2, Coupling> &table_by_component,
274 std::vector<Table<2, Coupling>> &tables_by_block);
275
281 template <int dim, int spacedim>
285 const Table<2, Coupling> &component_couplings);
286
294 template <int dim, int spacedim>
295 std::vector<Table<2, Coupling>>
298 const Table<2, Coupling> &component_couplings);
422 template <int dim, int spacedim, typename number = double>
423 void
425 const DoFHandler<dim, spacedim> &dof_handler,
426 SparsityPatternBase &sparsity_pattern,
427 const AffineConstraints<number> &constraints = {},
428 const bool keep_constrained_dofs = true,
430
496 template <int dim, int spacedim, typename number = double>
497 void
499 const DoFHandler<dim, spacedim> &dof_handler,
500 const Table<2, Coupling> &coupling,
501 SparsityPatternBase &sparsity_pattern,
502 const AffineConstraints<number> &constraints = {},
503 const bool keep_constrained_dofs = true,
505
526 template <int dim, int spacedim>
527 void
529 const DoFHandler<dim, spacedim> &dof_col,
530 SparsityPatternBase &sparsity);
531
577 template <int dim, int spacedim>
578 void
580 SparsityPatternBase &sparsity_pattern);
581
590 template <int dim, int spacedim, typename number>
591 void
593 const DoFHandler<dim, spacedim> &dof_handler,
594 SparsityPatternBase &sparsity_pattern,
595 const AffineConstraints<number> &constraints,
596 const bool keep_constrained_dofs = true,
598
599
619 template <int dim, int spacedim>
620 void
622 const DoFHandler<dim, spacedim> &dof,
623 SparsityPatternBase &sparsity,
624 const Table<2, Coupling> &cell_integrals_mask,
625 const Table<2, Coupling> &face_integrals_mask,
627
628
655 template <int dim, int spacedim, typename number>
656 void
658 const DoFHandler<dim, spacedim> &dof,
659 SparsityPatternBase &sparsity,
660 const AffineConstraints<number> &constraints,
661 const bool keep_constrained_dofs,
662 const Table<2, Coupling> &couplings,
663 const Table<2, Coupling> &face_couplings,
664 const types::subdomain_id subdomain_id,
665 const std::function<
667 const unsigned int)> &face_has_flux_coupling =
668 &internal::always_couple_on_faces<dim, spacedim>);
669
679 template <int dim, int spacedim>
680 void
682 const DoFHandler<dim, spacedim> &dof,
683 const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
684 SparsityPatternBase &sparsity_pattern);
685
703 template <int dim, int spacedim, typename number>
704 void
706 const DoFHandler<dim, spacedim> &dof,
707 const std::map<types::boundary_id, const Function<spacedim, number> *>
708 &boundary_ids,
709 const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
710 SparsityPatternBase &sparsity);
711
759 template <int dim, int spacedim, typename number>
760 void
762 AffineConstraints<number> &constraints);
763
831 template <int dim, int spacedim>
832 void
834 const DoFHandler<dim, spacedim> &coarse_grid,
835 const unsigned int coarse_component,
836 const DoFHandler<dim, spacedim> &fine_grid,
837 const unsigned int fine_component,
838 const InterGridMap<DoFHandler<dim, spacedim>> &coarse_to_fine_grid_map,
839 AffineConstraints<double> &constraints);
840
841
858 template <int dim, int spacedim>
859 void
861 const DoFHandler<dim, spacedim> &coarse_grid,
862 const unsigned int coarse_component,
863 const DoFHandler<dim, spacedim> &fine_grid,
864 const unsigned int fine_component,
865 const InterGridMap<DoFHandler<dim, spacedim>> &coarse_to_fine_grid_map,
866 std::vector<std::map<types::global_dof_index, float>>
867 &transfer_representation);
868
955 template <typename FaceIterator, typename number>
956 void
958 const FaceIterator &face_1,
960 AffineConstraints<number> &constraints,
961 const ComponentMask &component_mask = {},
962 const unsigned char combined_orientation =
965 const std::vector<unsigned int> &first_vector_components =
966 std::vector<unsigned int>(),
967 const number periodicity_factor = 1.);
968
969
989 template <int dim, int spacedim, typename number>
990 void
992 const std::vector<GridTools::PeriodicFacePair<
993 typename DoFHandler<dim, spacedim>::cell_iterator>> &periodic_faces,
994 AffineConstraints<number> &constraints,
995 const ComponentMask &component_mask = {},
996 const std::vector<unsigned int> &first_vector_components =
997 std::vector<unsigned int>(),
998 const number periodicity_factor = 1.);
999
1000
1001
1031 template <int dim, int spacedim, typename number>
1032 void
1034 const types::boundary_id b_id1,
1035 const types::boundary_id b_id2,
1036 const unsigned int direction,
1037 AffineConstraints<number> &constraints,
1038 const ComponentMask &component_mask = {},
1039 const number periodicity_factor = 1.);
1040
1041
1042
1068 template <int dim, int spacedim, typename number>
1069 void
1071 const types::boundary_id b_id,
1072 const unsigned int direction,
1073 AffineConstraints<number> &constraints,
1074 const ComponentMask &component_mask = {},
1075 const number periodicity_factor = 1.);
1076
1077 namespace internal
1078 {
1105 template <typename FaceIterator, typename number>
1106 void
1108 const FaceIterator &face_1,
1110 const FullMatrix<double> &transformation,
1111 AffineConstraints<number> &affine_constraints,
1112 const ComponentMask &component_mask,
1113 const unsigned char combined_orientation,
1114 const number periodicity_factor,
1115 const unsigned int level = numbers::invalid_unsigned_int);
1116 } // namespace internal
1117
1135 template <int dim, int spacedim>
1136 IndexSet
1138
1174 template <int dim, int spacedim>
1175 IndexSet
1176 extract_dofs(const DoFHandler<dim, spacedim> &dof_handler,
1177 const ComponentMask &component_mask);
1178
1202 template <int dim, int spacedim>
1203 IndexSet
1204 extract_dofs(const DoFHandler<dim, spacedim> &dof_handler,
1205 const BlockMask &block_mask);
1206
1211 template <int dim, int spacedim>
1212 void
1213 extract_level_dofs(const unsigned int level,
1214 const DoFHandler<dim, spacedim> &dof,
1215 const ComponentMask &component_mask,
1216 std::vector<bool> &selected_dofs);
1217
1222 template <int dim, int spacedim>
1223 void
1224 extract_level_dofs(const unsigned int level,
1225 const DoFHandler<dim, spacedim> &dof,
1226 const BlockMask &component_mask,
1227 std::vector<bool> &selected_dofs);
1228
1271 template <int dim, int spacedim>
1272 IndexSet
1274 const ComponentMask &component_mask = {},
1275 const std::set<types::boundary_id> &boundary_ids = {});
1276
1293 template <int dim, int spacedim>
1294 void
1296 const DoFHandler<dim, spacedim> &dof_handler,
1297 const ComponentMask &component_mask,
1298 std::vector<bool> &selected_dofs,
1299 const std::set<types::boundary_id> &boundary_ids =
1300 std::set<types::boundary_id>());
1301
1329 template <int dim, int spacedim, typename number = double>
1330 IndexSet
1332 const DoFHandler<dim, spacedim> &dof_handler,
1333 const std::function<
1335 &predicate,
1336 const AffineConstraints<number> &constraints = {});
1337
1373 template <int dim, int spacedim>
1374 std::vector<std::vector<bool>>
1376 const ComponentMask &component_mask = {});
1377
1383 template <int dim, int spacedim>
1384 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
1385 "Use the other function that returns the constant modes by value, rather than via an argument.")
1386 void extract_constant_modes(const DoFHandler<dim, spacedim> &dof_handler,
1387 const ComponentMask &component_mask,
1388 std::vector<std::vector<bool>> &constant_modes);
1389
1393 template <int dim, int spacedim>
1394 std::vector<std::vector<bool>>
1395 extract_level_constant_modes(const unsigned int level,
1396 const DoFHandler<dim, spacedim> &dof_handler,
1397 const ComponentMask &component_mask = {});
1398
1404 template <int dim, int spacedim>
1405 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
1406 "Use the other function that returns the constant modes by value, rather than via an argument.")
1408 const unsigned int level,
1409 const DoFHandler<dim, spacedim> &dof_handler,
1410 const ComponentMask &component_mask,
1411 std::vector<std::vector<bool>> &constant_modes);
1412
1418 template <int dim, int spacedim>
1419 std::vector<std::vector<double>>
1420 extract_rigid_body_modes(const Mapping<dim, spacedim> &mapping,
1421 const DoFHandler<dim, spacedim> &dof_handler,
1422 const ComponentMask &component_mask = {});
1423
1427 template <int dim, int spacedim>
1428 std::vector<std::vector<double>>
1429 extract_level_rigid_body_modes(const unsigned int level,
1430 const Mapping<dim, spacedim> &mapping,
1431 const DoFHandler<dim, spacedim> &dof_handler,
1432 const ComponentMask &component_mask = {});
1519 template <int dim, int spacedim>
1520 std::map<typename DoFHandler<dim - 1, spacedim>::active_cell_iterator,
1521 std::pair<typename DoFHandler<dim, spacedim>::active_cell_iterator,
1522 unsigned int>>
1526 &c1_to_c0,
1527 const DoFHandler<dim, spacedim> &c0_dh,
1528 const DoFHandler<dim - 1, spacedim> &c1_dh);
1529
1547 template <int dim, int spacedim>
1548 void
1550 const types::subdomain_id subdomain_id,
1551 std::vector<bool> &selected_dofs);
1552
1567 template <int dim, int spacedim>
1568 IndexSet
1570
1587 template <int dim, int spacedim>
1590 IndexSet &dof_set);
1591
1598 template <int dim, int spacedim>
1599 IndexSet
1601 const DoFHandler<dim, spacedim> &dof_handler,
1602 const unsigned int level);
1603
1612 template <int dim, int spacedim>
1615 const DoFHandler<dim, spacedim> &dof_handler,
1616 IndexSet &dof_set,
1617 const unsigned int level);
1618
1628 template <int dim, int spacedim>
1629 IndexSet
1631
1642 template <int dim, int spacedim>
1645 IndexSet &dof_set);
1646
1656 template <int dim, int spacedim>
1657 std::vector<IndexSet>
1659 const ComponentMask &components = {});
1660
1674 template <int dim, int spacedim>
1675 std::vector<IndexSet>
1677 const DoFHandler<dim, spacedim> &dof_handler);
1678
1693 template <int dim, int spacedim>
1694 std::vector<IndexSet>
1696 const DoFHandler<dim, spacedim> &dof_handler);
1697
1702 template <int dim, int spacedim>
1703 IndexSet
1705 const DoFHandler<dim, spacedim> &dof_handler,
1706 const unsigned int level);
1707
1714 template <int dim, int spacedim>
1717 const DoFHandler<dim, spacedim> &dof_handler,
1718 const unsigned int level,
1719 IndexSet &dof_set);
1720
1721
1752 template <int dim, int spacedim>
1753 void
1755 std::vector<types::subdomain_id> &subdomain);
1756
1782 template <int dim, int spacedim>
1783 unsigned int
1785 const DoFHandler<dim, spacedim> &dof_handler,
1786 const types::subdomain_id subdomain);
1787
1807 template <int dim, int spacedim>
1808 void
1810 const DoFHandler<dim, spacedim> &dof_handler,
1811 const types::subdomain_id subdomain,
1812 std::vector<unsigned int> &n_dofs_on_subdomain);
1813
1835 template <int dim, int spacedim>
1836 IndexSet
1838 const DoFHandler<dim, spacedim> &dof_handler,
1839 const types::subdomain_id subdomain);
1896 template <int dim, int spacedim>
1897 std::vector<types::global_dof_index>
1899 const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
1900 &patch);
1901
1922 template <int dim, int spacedim>
1923 void
1925 const DoFHandler<dim, spacedim> &dof_handler,
1926 const unsigned int level,
1927 const std::vector<bool> &selected_dofs = {},
1928 const types::global_dof_index offset = 0);
1929
1981 template <int dim, int spacedim>
1982 std::vector<unsigned int>
1984 const DoFHandler<dim, spacedim> &dof_handler,
1985 const unsigned int level,
1986 const bool interior_dofs_only,
1987 const bool boundary_patches = false,
1988 const bool level_boundary_patches = false,
1989 const bool single_cell_patches = false,
1990 const bool invert_vertex_mapping = false);
1991
2006 template <int dim, int spacedim>
2007 std::vector<unsigned int>
2009 const DoFHandler<dim, spacedim> &dof_handler,
2010 const unsigned int level,
2011 const BlockMask &exclude_boundary_dofs = BlockMask(),
2012 const bool boundary_patches = false,
2013 const bool level_boundary_patches = false,
2014 const bool single_cell_patches = false,
2015 const bool invert_vertex_mapping = false);
2016
2054 template <int dim, int spacedim>
2055 void
2057 const DoFHandler<dim, spacedim> &dof_handler,
2058 const unsigned int level,
2059 const bool interior_dofs_only,
2060 const bool boundary_dofs = false);
2061
2081 template <int dim, int spacedim>
2082 void
2084 const DoFHandler<dim, spacedim> &dof_handler,
2085 const unsigned int level,
2086 const bool interior_dofs_only = false);
2087
2130 template <int dim, int spacedim>
2131 std::vector<types::global_dof_index>
2133 const DoFHandler<dim, spacedim> &dof_handler,
2134 const bool vector_valued_once = false,
2135 const std::vector<unsigned int> &target_component = {});
2136
2153 template <int dim, int spacedim>
2154 std::vector<types::global_dof_index>
2156 const std::vector<unsigned int> &target_block =
2157 std::vector<unsigned int>());
2158
2174 template <int dim, int spacedim>
2177 std::vector<unsigned int> &active_fe_indices);
2178
2206 template <int dim, int spacedim>
2207 unsigned int
2209 const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
2210 &patch);
2211
2234 template <int dim, int spacedim>
2235 void
2237 std::vector<types::global_dof_index> &mapping);
2238
2249 template <int dim, int spacedim>
2250 void
2252 const std::set<types::boundary_id> &boundary_ids,
2253 std::vector<types::global_dof_index> &mapping);
2254
2283 template <int dim, int spacedim>
2284 void
2286 const DoFHandler<dim, spacedim> &dof_handler,
2287 std::vector<Point<spacedim>> &support_points,
2288 const ComponentMask &mask = {});
2289
2293 template <int dim, int spacedim>
2294 void
2297 const DoFHandler<dim, spacedim> &dof_handler,
2298 std::vector<Point<spacedim>> &support_points,
2299 const ComponentMask &mask = {});
2300
2330 template <int dim, int spacedim>
2331 std::map<types::global_dof_index, Point<spacedim>>
2333 const DoFHandler<dim, spacedim> &dof_handler,
2334 const ComponentMask &mask = {});
2335
2339 template <int dim, int spacedim>
2340 std::map<types::global_dof_index, Point<spacedim>>
2343 const DoFHandler<dim, spacedim> &dof_handler,
2344 const ComponentMask &mask = {});
2345
2351 template <int dim, int spacedim>
2354 const Mapping<dim, spacedim> &mapping,
2355 const DoFHandler<dim, spacedim> &dof_handler,
2356 std::map<types::global_dof_index, Point<spacedim>> &support_points,
2357 const ComponentMask &mask = {});
2358
2364 template <int dim, int spacedim>
2368 const DoFHandler<dim, spacedim> &dof_handler,
2369 std::map<types::global_dof_index, Point<spacedim>> &support_points,
2370 const ComponentMask &mask = {});
2371
2372
2393 template <int dim, int spacedim, class Comp>
2394 void
2396 const Mapping<dim, spacedim> &mapping,
2397 const DoFHandler<dim, spacedim> &dof_handler,
2399 &point_to_index_map);
2437 template <int dim, int spacedim, typename Number>
2438 void
2440 const Vector<Number> &cell_data,
2441 Vector<double> &dof_data,
2442 const unsigned int component = 0);
2443
2444
2522 template <int spacedim>
2523 void
2525 std::ostream &out,
2526 const std::map<types::global_dof_index, Point<spacedim>> &support_points);
2527
2528
2569 template <int dim, int spacedim, typename number>
2570 void
2572 const DoFHandler<dim, spacedim> &dof,
2573 const types::boundary_id boundary_id,
2574 AffineConstraints<number> &zero_boundary_constraints,
2575 const ComponentMask &component_mask = {});
2576
2587 template <int dim, int spacedim, typename number>
2588 void
2590 const DoFHandler<dim, spacedim> &dof,
2591 AffineConstraints<number> &zero_boundary_constraints,
2592 const ComponentMask &component_mask = {});
2593
2638} // namespace DoFTools
2639
2640
2641
2642/* ------------------------- inline functions -------------- */
2643
2644#ifndef DOXYGEN
2645
2646namespace DoFTools
2647{
2653 inline Coupling
2654 operator|=(Coupling &c1, const Coupling c2)
2655 {
2656 if (c2 == always)
2657 c1 = always;
2658 else if (c1 != always && c2 == nonzero)
2659 return c1 = nonzero;
2660 return c1;
2661 }
2662
2663
2669 inline Coupling
2670 operator|(const Coupling c1, const Coupling c2)
2671 {
2672 if (c1 == always || c2 == always)
2673 return always;
2674 if (c1 == nonzero || c2 == nonzero)
2675 return nonzero;
2676 return none;
2677 }
2678
2679
2680 // ---------------------- inline and template functions --------------------
2681
2682 template <int dim, int spacedim, class Comp>
2683 void
2685 const Mapping<dim, spacedim> &mapping,
2686 const DoFHandler<dim, spacedim> &dof_handler,
2688 &point_to_index_map)
2689 {
2690 // let the checking of arguments be
2691 // done by the function first
2692 // called
2693 std::vector<Point<spacedim>> support_points(dof_handler.n_dofs());
2694 map_dofs_to_support_points(mapping, dof_handler, support_points);
2695 // now copy over the results of the
2696 // previous function into the
2697 // output arg
2698 point_to_index_map.clear();
2699 for (types::global_dof_index i = 0; i < dof_handler.n_dofs(); ++i)
2700 point_to_index_map[support_points[i]] = i;
2701 }
2702} // namespace DoFTools
2703
2704#endif
2705
2707
2708#endif
types::global_dof_index n_dofs() const
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
static constexpr unsigned char default_combined_face_orientation()
#define DEAL_II_DEPRECATED
Definition config.h:205
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:175
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
unsigned int level
Definition grid_out.cc:4626
#define DeclException0(Exception0)
Definition exceptions.h:466
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_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
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_zero_boundary_constraints(const DoFHandler< dim, spacedim > &dof, const types::boundary_id boundary_id, AffineConstraints< number > &zero_boundary_constraints, const ComponentMask &component_mask={})
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)
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
InitializeLibrary operator|(const InitializeLibrary f1, const InitializeLibrary f2)
void set_periodicity_constraints(const FaceIterator &face_1, const std_cxx20::type_identity_t< FaceIterator > &face_2, const FullMatrix< double > &transformation, AffineConstraints< number > &affine_constraints, const ComponentMask &component_mask, const unsigned char combined_orientation, const number periodicity_factor, const unsigned int level=numbers::invalid_unsigned_int)
void get_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, std::vector< types::subdomain_id > &subdomain)
IndexSet dof_indices_with_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain)
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={})
Definition dof_tools.cc:826
std::vector< IndexSet > locally_owned_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
IndexSet extract_boundary_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask={}, const std::set< types::boundary_id > &boundary_ids={})
Definition dof_tools.cc:616
std::vector< std::vector< bool > > extract_level_constant_modes(const unsigned int level, const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask={})
void make_periodicity_constraints(const FaceIterator &face_1, const std_cxx20::type_identity_t< FaceIterator > &face_2, AffineConstraints< number > &constraints, const ComponentMask &component_mask={}, const unsigned char combined_orientation=ReferenceCell::default_combined_face_orientation(), const FullMatrix< double > &matrix=FullMatrix< double >(), const std::vector< unsigned int > &first_vector_components=std::vector< unsigned int >(), const number periodicity_factor=1.)
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
std::vector< std::vector< bool > > extract_constant_modes(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask={})
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)
std::map< typename DoFHandler< dim - 1, spacedim >::active_cell_iterator, std::pair< typename DoFHandler< dim, spacedim >::active_cell_iterator, unsigned int > > map_boundary_to_bulk_dof_iterators(const std::map< typename Triangulation< dim - 1, spacedim >::cell_iterator, typename Triangulation< dim, spacedim >::face_iterator > &c1_to_c0, const DoFHandler< dim, spacedim > &c0_dh, const DoFHandler< dim - 1, spacedim > &c1_dh)
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:416
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< std::vector< double > > extract_rigid_body_modes(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask={})
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={})
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:503
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:744
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:332
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)
std::vector< types::global_dof_index > get_dofs_on_patch(const std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > &patch)
std::vector< IndexSet > locally_owned_dofs_per_component(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &components={})
Definition dof_tools.cc:469
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)
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)
std::vector< std::vector< double > > extract_level_rigid_body_modes(const unsigned int level, const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask={})
@ matrix
Contents is actually a matrix.
Definition hp.h:117
const types::subdomain_id invalid_subdomain_id
Definition types.h:341
static const unsigned int invalid_unsigned_int
Definition types.h:220
typename type_identity< T >::type type_identity_t
Definition type_traits.h:95
STL namespace.
unsigned int subdomain_id
Definition types.h:43