Reference documentation for deal.II version 9.0.0
dof_tools.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2018 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 #ifndef dealii_dof_tools_h
17 #define dealii_dof_tools_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/index_set.h>
23 #include <deal.II/base/point.h>
24 #include <deal.II/lac/constraint_matrix.h>
25 #include <deal.II/dofs/dof_handler.h>
26 #include <deal.II/fe/component_mask.h>
27 #include <deal.II/hp/dof_handler.h>
28 
29 #include <vector>
30 #include <set>
31 #include <map>
32 #include <ostream>
33 
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 class BlockMask;
38 template <int dim, typename RangeNumberType> class Function;
39 template <int dim, int spacedim> class FiniteElement;
40 namespace hp
41 {
42  template <int dim, int spacedim> class MappingCollection;
43  template <int dim, int spacedim> class FECollection;
44 }
45 template <class MeshType> class InterGridMap;
46 template <int dim, int spacedim> class Mapping;
47 class SparsityPattern;
48 template <int dim, class T> class Table;
49 template <typename Number> class Vector;
50 
51 namespace GridTools
52 {
53  template <typename CellIterator> struct PeriodicFacePair;
54 }
55 
56 
176 namespace DoFTools
177 {
178 
190  enum Coupling
191  {
206  };
207 
222  template <int dim, int spacedim>
223  void
225  const Table<2, Coupling> &table_by_component,
226  std::vector<Table<2,Coupling> > &tables_by_block);
227 
237  template <int dim, int spacedim>
238  void
240  const Table<2, Coupling> &table_by_component,
241  std::vector<Table<2,Coupling> > &tables_by_block);
242 
248  template <int dim, int spacedim>
251  const Table<2,Coupling> &component_couplings);
252 
260  template <int dim, int spacedim>
261  std::vector<Table<2,Coupling> >
263  const Table<2,Coupling> &component_couplings);
279  template <int dim, int spacedim>
280  DEAL_II_DEPRECATED
281  unsigned int
283 
291  template <int dim, int spacedim>
292  DEAL_II_DEPRECATED
293  unsigned int
295 
296 
307  template <int dim, int spacedim>
308  DEAL_II_DEPRECATED
309  unsigned int
311 
322  template <int dim, int spacedim>
323  DEAL_II_DEPRECATED
324  unsigned int
326 
337  template <int dim, int spacedim>
338  DEAL_II_DEPRECATED
339  unsigned int
341 
352  template <int dim, int spacedim>
353  DEAL_II_DEPRECATED
354  unsigned int
356 
368  template <int dim, int spacedim>
369  DEAL_II_DEPRECATED
370  unsigned int
372 
384  template <int dim, int spacedim>
385  DEAL_II_DEPRECATED
386  unsigned int
388 
400  template <int dim, int spacedim>
401  DEAL_II_DEPRECATED
402  bool
404 
416  template <int dim, int spacedim>
417  DEAL_II_DEPRECATED
418  bool
420 
545  template <typename DoFHandlerType, typename SparsityPatternType>
546  void
547  make_sparsity_pattern (const DoFHandlerType &dof_handler,
548  SparsityPatternType &sparsity_pattern,
549  const ConstraintMatrix &constraints = ConstraintMatrix(),
550  const bool keep_constrained_dofs = true,
552 
618  template <typename DoFHandlerType, typename SparsityPatternType>
619  void
620  make_sparsity_pattern (const DoFHandlerType &dof_handler,
621  const Table<2, Coupling> &coupling,
622  SparsityPatternType &sparsity_pattern,
623  const ConstraintMatrix &constraints = ConstraintMatrix(),
624  const bool keep_constrained_dofs = true,
626 
647  template <typename DoFHandlerType, typename SparsityPatternType>
648  void
649  make_sparsity_pattern (const DoFHandlerType &dof_row,
650  const DoFHandlerType &dof_col,
651  SparsityPatternType &sparsity);
652 
698  template <typename DoFHandlerType, typename SparsityPatternType>
699  void
700  make_flux_sparsity_pattern (const DoFHandlerType &dof_handler,
701  SparsityPatternType &sparsity_pattern);
702 
711  template <typename DoFHandlerType, typename SparsityPatternType>
712  void
713  make_flux_sparsity_pattern (const DoFHandlerType &dof_handler,
714  SparsityPatternType &sparsity_pattern,
715  const ConstraintMatrix &constraints,
716  const bool keep_constrained_dofs = true,
718 
738  template <typename DoFHandlerType, typename SparsityPatternType>
739  void
740  make_flux_sparsity_pattern (const DoFHandlerType &dof,
741  SparsityPatternType &sparsity,
742  const Table<2,Coupling> &cell_integrals_mask,
743  const Table<2,Coupling> &face_integrals_mask,
753  template <typename DoFHandlerType, typename SparsityPatternType>
754  void
755  make_flux_sparsity_pattern (const DoFHandlerType &dof,
756  SparsityPatternType &sparsity,
757  const ConstraintMatrix &constraints,
758  const bool keep_constrained_dofs,
759  const Table<2,Coupling> &couplings,
760  const Table<2,Coupling> &face_couplings,
761  const types::subdomain_id subdomain_id);
762 
772  template <typename DoFHandlerType, typename SparsityPatternType>
773  void
774  make_boundary_sparsity_pattern (const DoFHandlerType &dof,
775  const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
776  SparsityPatternType &sparsity_pattern);
777 
794  template <typename DoFHandlerType, typename SparsityPatternType, typename number>
795  void
797  (const DoFHandlerType &dof,
798  const std::map<types::boundary_id, const Function<DoFHandlerType::space_dimension,number>*> &boundary_ids,
799  const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
800  SparsityPatternType &sparsity);
801 
852  template <typename DoFHandlerType>
853  void
854  make_hanging_node_constraints (const DoFHandlerType &dof_handler,
855  ConstraintMatrix &constraints);
856 
924  template <int dim, int spacedim>
925  void
927  const unsigned int coarse_component,
928  const DoFHandler<dim,spacedim> &fine_grid,
929  const unsigned int fine_component,
930  const InterGridMap<DoFHandler<dim,spacedim> > &coarse_to_fine_grid_map,
931  ConstraintMatrix &constraints);
932 
933 
950  template <int dim, int spacedim>
951  void
953  const unsigned int coarse_component,
954  const DoFHandler<dim,spacedim> &fine_grid,
955  const unsigned int fine_component,
956  const InterGridMap<DoFHandler<dim,spacedim> > &coarse_to_fine_grid_map,
957  std::vector<std::map<types::global_dof_index, float> > &transfer_representation);
958 
1121  template <typename FaceIterator>
1122  void
1124  (const FaceIterator &face_1,
1125  const typename identity<FaceIterator>::type &face_2,
1126  ::ConstraintMatrix &constraint_matrix,
1127  const ComponentMask &component_mask = ComponentMask(),
1128  const bool face_orientation = true,
1129  const bool face_flip = false,
1130  const bool face_rotation = false,
1131  const FullMatrix<double> &matrix = FullMatrix<double>(),
1132  const std::vector<unsigned int> &first_vector_components = std::vector<unsigned int>());
1133 
1134 
1135 
1157  template <typename DoFHandlerType>
1158  void
1161  &periodic_faces,
1162  ::ConstraintMatrix &constraint_matrix,
1163  const ComponentMask &component_mask = ComponentMask(),
1164  const std::vector<unsigned int> &first_vector_components = std::vector<unsigned int>());
1165 
1166 
1167 
1198  template <typename DoFHandlerType>
1199  void
1201  (const DoFHandlerType &dof_handler,
1202  const types::boundary_id b_id1,
1203  const types::boundary_id b_id2,
1204  const int direction,
1205  ::ConstraintMatrix &constraint_matrix,
1206  const ComponentMask &component_mask = ComponentMask());
1207 
1208 
1209 
1234  template <typename DoFHandlerType>
1235  void
1237  (const DoFHandlerType &dof_handler,
1238  const types::boundary_id b_id,
1239  const int direction,
1240  ::ConstraintMatrix &constraint_matrix,
1241  const ComponentMask &component_mask = ComponentMask());
1242 
1267  template <int dim, int spacedim>
1268  void
1270  std::vector<bool> &selected_dofs);
1271 
1276  template <int dim, int spacedim>
1277  IndexSet
1279 
1313  template <int dim, int spacedim>
1314  void
1315  extract_dofs (const DoFHandler<dim,spacedim> &dof_handler,
1316  const ComponentMask &component_mask,
1317  std::vector<bool> &selected_dofs);
1318 
1322  template <int dim, int spacedim>
1323  void
1324  extract_dofs (const hp::DoFHandler<dim,spacedim> &dof_handler,
1325  const ComponentMask &component_mask,
1326  std::vector<bool> &selected_dofs);
1327 
1352  template <int dim, int spacedim>
1353  void
1354  extract_dofs (const DoFHandler<dim,spacedim> &dof_handler,
1355  const BlockMask &block_mask,
1356  std::vector<bool> &selected_dofs);
1357 
1361  template <int dim, int spacedim>
1362  void
1363  extract_dofs (const hp::DoFHandler<dim,spacedim> &dof_handler,
1364  const BlockMask &block_mask,
1365  std::vector<bool> &selected_dofs);
1366 
1371  template <typename DoFHandlerType>
1372  void
1373  extract_level_dofs (const unsigned int level,
1374  const DoFHandlerType &dof,
1375  const ComponentMask &component_mask,
1376  std::vector<bool> &selected_dofs);
1377 
1382  template <typename DoFHandlerType>
1383  void
1384  extract_level_dofs (const unsigned int level,
1385  const DoFHandlerType &dof,
1386  const BlockMask &component_mask,
1387  std::vector<bool> &selected_dofs);
1388 
1440  template <typename DoFHandlerType>
1441  void
1442  extract_boundary_dofs (const DoFHandlerType &dof_handler,
1443  const ComponentMask &component_mask,
1444  std::vector<bool> &selected_dofs,
1445  const std::set<types::boundary_id> &boundary_ids = std::set<types::boundary_id>());
1446 
1477  template <typename DoFHandlerType>
1478  void
1479  extract_boundary_dofs (const DoFHandlerType &dof_handler,
1480  const ComponentMask &component_mask,
1481  IndexSet &selected_dofs,
1482  const std::set<types::boundary_id> &boundary_ids = std::set<types::boundary_id>());
1483 
1500  template <typename DoFHandlerType>
1501  void
1502  extract_dofs_with_support_on_boundary (const DoFHandlerType &dof_handler,
1503  const ComponentMask &component_mask,
1504  std::vector<bool> &selected_dofs,
1505  const std::set<types::boundary_id> &boundary_ids = std::set<types::boundary_id>());
1506 
1534  template <typename DoFHandlerType>
1535  IndexSet
1536  extract_dofs_with_support_contained_within(const DoFHandlerType &dof_handler,
1537  const std::function< bool(const typename DoFHandlerType::active_cell_iterator &)> &predicate,
1538  const ConstraintMatrix &cm = ConstraintMatrix());
1539 
1573  template <typename DoFHandlerType>
1574  void
1575  extract_constant_modes (const DoFHandlerType &dof_handler,
1576  const ComponentMask &component_mask,
1577  std::vector<std::vector<bool> > &constant_modes);
1579 
1593  template <typename DoFHandlerType>
1594  void
1595  extract_subdomain_dofs (const DoFHandlerType &dof_handler,
1596  const types::subdomain_id subdomain_id,
1597  std::vector<bool> &selected_dofs);
1598 
1599 
1609  template <typename DoFHandlerType>
1610  DEAL_II_DEPRECATED
1611  void
1612  extract_locally_owned_dofs (const DoFHandlerType &dof_handler,
1613  IndexSet &dof_set);
1614 
1615 
1630  template <typename DoFHandlerType>
1631  void
1632  extract_locally_active_dofs (const DoFHandlerType &dof_handler,
1633  IndexSet &dof_set);
1634 
1644  template <typename DoFHandlerType>
1645  void
1646  extract_locally_relevant_dofs (const DoFHandlerType &dof_handler,
1647  IndexSet &dof_set);
1648 
1665  template <typename DoFHandlerType>
1666  std::vector<IndexSet>
1667  locally_owned_dofs_per_subdomain (const DoFHandlerType &dof_handler);
1668 
1685  template <typename DoFHandlerType>
1686  std::vector<IndexSet>
1687  locally_relevant_dofs_per_subdomain (const DoFHandlerType &dof_handler);
1688 
1689 
1694  template <typename DoFHandlerType>
1695  void
1696  extract_locally_relevant_level_dofs (const DoFHandlerType &dof_handler,
1697  const unsigned int level,
1698  IndexSet &dof_set);
1699 
1700 
1731  template <typename DoFHandlerType>
1732  void
1733  get_subdomain_association (const DoFHandlerType &dof_handler,
1734  std::vector<types::subdomain_id> &subdomain);
1735 
1761  template <typename DoFHandlerType>
1762  unsigned int
1763  count_dofs_with_subdomain_association (const DoFHandlerType &dof_handler,
1764  const types::subdomain_id subdomain);
1765 
1785  template <typename DoFHandlerType>
1786  void
1787  count_dofs_with_subdomain_association (const DoFHandlerType &dof_handler,
1788  const types::subdomain_id subdomain,
1789  std::vector<unsigned int> &n_dofs_on_subdomain);
1790 
1812  template <typename DoFHandlerType>
1813  IndexSet
1814  dof_indices_with_subdomain_association (const DoFHandlerType &dof_handler,
1815  const types::subdomain_id subdomain);
1816  // @}
1826 
1878  template <typename DoFHandlerType>
1879  std::vector<types::global_dof_index>
1880  get_dofs_on_patch (const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
1881 
1902  template <int dim, int spacedim>
1903  void make_cell_patches(SparsityPattern &block_list,
1904  const DoFHandler<dim,spacedim> &dof_handler,
1905  const unsigned int level,
1906  const std::vector<bool> &selected_dofs = std::vector<bool>(),
1907  const types::global_dof_index offset = 0);
1908 
1959  template <typename DoFHandlerType>
1960  std::vector<unsigned int>
1962  const DoFHandlerType &dof_handler,
1963  const unsigned int level,
1964  const bool interior_dofs_only,
1965  const bool boundary_patches = false,
1966  const bool level_boundary_patches = false,
1967  const bool single_cell_patches = false,
1968  const bool invert_vertex_mapping = false);
1969 
1983  template <typename DoFHandlerType>
1984  std::vector<unsigned int>
1986  const DoFHandlerType &dof_handler,
1987  const unsigned int level,
1988  const BlockMask &exclude_boundary_dofs = BlockMask(),
1989  const bool boundary_patches = false,
1990  const bool level_boundary_patches = false,
1991  const bool single_cell_patches = false,
1992  const bool invert_vertex_mapping = false);
1993 
2031  template <typename DoFHandlerType>
2032  void make_child_patches(SparsityPattern &block_list,
2033  const DoFHandlerType &dof_handler,
2034  const unsigned int level,
2035  const bool interior_dofs_only,
2036  const bool boundary_dofs = false);
2037 
2057  template <typename DoFHandlerType>
2058  void make_single_patch(SparsityPattern &block_list,
2059  const DoFHandlerType &dof_handler,
2060  const unsigned int level,
2061  const bool interior_dofs_only = false);
2062 
2105  template <typename DoFHandlerType>
2106  void
2107  count_dofs_per_component (const DoFHandlerType &dof_handler,
2108  std::vector<types::global_dof_index> &dofs_per_component,
2109  const bool vector_valued_once = false,
2110  std::vector<unsigned int> target_component
2111  = std::vector<unsigned int>());
2112 
2129  template <typename DoFHandlerType>
2130  void
2131  count_dofs_per_block (const DoFHandlerType &dof,
2132  std::vector<types::global_dof_index> &dofs_per_block,
2133  const std::vector<unsigned int> &target_block
2134  = std::vector<unsigned int>());
2135 
2146  template <typename DoFHandlerType>
2147  void
2148  get_active_fe_indices (const DoFHandlerType &dof_handler,
2149  std::vector<unsigned int> &active_fe_indices);
2150 
2185  template <typename DoFHandlerType>
2186  unsigned int
2187  count_dofs_on_patch (const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
2188 
2211  template <typename DoFHandlerType>
2212  void
2213  map_dof_to_boundary_indices (const DoFHandlerType &dof_handler,
2214  std::vector<types::global_dof_index> &mapping);
2215 
2226  template <typename DoFHandlerType>
2227  void
2228  map_dof_to_boundary_indices (const DoFHandlerType &dof_handler,
2229  const std::set<types::boundary_id> &boundary_ids,
2230  std::vector<types::global_dof_index> &mapping);
2231 
2249  template <int dim, int spacedim>
2250  void
2252  const DoFHandler<dim,spacedim> &dof_handler,
2253  std::vector<Point<spacedim> > &support_points);
2254 
2259  template <int dim, int spacedim>
2260  void
2261  map_dofs_to_support_points (const ::hp::MappingCollection<dim,spacedim> &mapping,
2262  const hp::DoFHandler<dim,spacedim> &dof_handler,
2263  std::vector<Point<spacedim> > &support_points);
2264 
2292  template <int dim, int spacedim>
2293  void
2295  const DoFHandler<dim,spacedim> &dof_handler,
2296  std::map<types::global_dof_index, Point<spacedim> > &support_points);
2297 
2301  template <int dim, int spacedim>
2302  void
2303  map_dofs_to_support_points (const ::hp::MappingCollection<dim,spacedim> &mapping,
2304  const hp::DoFHandler<dim,spacedim> &dof_handler,
2305  std::map<types::global_dof_index, Point<spacedim> > &support_points);
2306 
2307 
2328  template <typename DoFHandlerType, class Comp>
2329  void
2332  const DoFHandlerType &dof_handler,
2333  std::map<Point<DoFHandlerType::space_dimension>, types::global_dof_index, Comp> &point_to_index_map);
2371  template <typename DoFHandlerType, typename Number>
2372  void
2373  distribute_cell_to_dof_vector (const DoFHandlerType &dof_handler,
2374  const Vector<Number> &cell_data,
2375  Vector<double> &dof_data,
2376  const unsigned int component = 0);
2377 
2378 
2451  template <int spacedim>
2452  void
2453  write_gnuplot_dof_support_point_info(std::ostream &out,
2454  const std::map<types::global_dof_index, Point<spacedim> > &support_points);
2455 
2456 
2497  template <int dim, int spacedim, template <int, int> class DoFHandlerType>
2498  void
2499  make_zero_boundary_constraints (const DoFHandlerType<dim,spacedim> &dof,
2500  const types::boundary_id boundary_id,
2501  ConstraintMatrix &zero_boundary_constraints,
2502  const ComponentMask &component_mask = ComponentMask());
2503 
2514  template <int dim, int spacedim, template <int, int> class DoFHandlerType>
2515  void
2516  make_zero_boundary_constraints (const DoFHandlerType<dim,spacedim> &dof,
2517  ConstraintMatrix &zero_boundary_constraints,
2518  const ComponentMask &component_mask = ComponentMask());
2519 
2564 }
2565 
2566 
2567 
2568 /* ------------------------- inline functions -------------- */
2569 
2570 #ifndef DOXYGEN
2571 
2572 namespace DoFTools
2573 {
2579  inline
2581  const Coupling c2)
2582  {
2583  if (c2 == always)
2584  c1 = always;
2585  else if (c1 != always && c2 == nonzero)
2586  return c1 = nonzero;
2587  return c1;
2588  }
2589 
2590 
2596  inline
2597  Coupling operator | (const Coupling c1,
2598  const Coupling c2)
2599  {
2600  if (c1 == always || c2 == always)
2601  return always;
2602  if (c1 == nonzero || c2 == nonzero)
2603  return nonzero;
2604  return none;
2605  }
2606 
2607 
2608 // ---------------------- inline and template functions --------------------
2609 
2610  template <int dim, int spacedim>
2611  inline
2612  unsigned int
2614  {
2615  return dh.get_fe().dofs_per_cell;
2616  }
2617 
2618 
2619  template <int dim, int spacedim>
2620  inline
2621  unsigned int
2623  {
2624  return dh.get_fe().dofs_per_face;
2625  }
2626 
2627 
2628  template <int dim, int spacedim>
2629  inline
2630  unsigned int
2632  {
2633  return dh.get_fe().dofs_per_vertex;
2634  }
2635 
2636 
2637  template <int dim, int spacedim>
2638  inline
2639  unsigned int
2641  {
2642  return dh.get_fe().n_components();
2643  }
2644 
2645 
2646 
2647  template <int dim, int spacedim>
2648  inline
2649  bool
2651  {
2652  return dh.get_fe().is_primitive();
2653  }
2654 
2655 
2656  template <int dim, int spacedim>
2657  inline
2658  unsigned int
2660  {
2661  return dh.get_fe_collection().max_dofs_per_cell ();
2662  }
2663 
2664 
2665  template <int dim, int spacedim>
2666  inline
2667  unsigned int
2669  {
2670  return dh.get_fe_collection().max_dofs_per_face ();
2671  }
2672 
2673 
2674  template <int dim, int spacedim>
2675  inline
2676  unsigned int
2678  {
2679  return dh.get_fe_collection().max_dofs_per_vertex ();
2680  }
2681 
2682 
2683  template <int dim, int spacedim>
2684  inline
2685  unsigned int
2687  {
2688  return dh.get_fe(0).n_components();
2689  }
2690 
2691 
2692  template <int dim, int spacedim>
2693  inline
2694  bool
2696  {
2697  return dh.get_fe(0).is_primitive();
2698  }
2699 
2700 
2701  template <typename DoFHandlerType, class Comp>
2702  void
2704  (
2706  const DoFHandlerType &dof_handler,
2707  std::map<Point<DoFHandlerType::space_dimension>, types::global_dof_index, Comp> &point_to_index_map)
2708  {
2709  // let the checking of arguments be
2710  // done by the function first
2711  // called
2712  std::vector<Point<DoFHandlerType::space_dimension> > support_points (dof_handler.n_dofs());
2713  map_dofs_to_support_points (mapping, dof_handler, support_points);
2714  // now copy over the results of the
2715  // previous function into the
2716  // output arg
2717  point_to_index_map.clear ();
2718  for (types::global_dof_index i=0; i<dof_handler.n_dofs(); ++i)
2719  point_to_index_map[support_points[i]] = i;
2720  }
2721 }
2722 
2723 #endif
2724 
2725 DEAL_II_NAMESPACE_CLOSE
2726 
2727 #endif
void make_child_patches(SparsityPattern &block_list, const DoFHandlerType &dof_handler, const unsigned int level, const bool interior_dofs_only, const bool boundary_dofs=false)
Definition: dof_tools.cc:2341
unsigned int max_dofs_per_vertex(const DoFHandler< dim, spacedim > &dh)
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)
std::vector< IndexSet > locally_owned_dofs_per_subdomain(const DoFHandlerType &dof_handler)
Definition: dof_tools.cc:1276
const types::subdomain_id invalid_subdomain_id
Definition: types.h:248
static ::ExceptionBase & ExcGridsDontMatch()
void make_boundary_sparsity_pattern(const DoFHandlerType &dof, const std::vector< types::global_dof_index > &dof_to_boundary_mapping, SparsityPatternType &sparsity_pattern)
void count_dofs_per_component(const DoFHandlerType &dof_handler, std::vector< types::global_dof_index > &dofs_per_component, const bool vector_valued_once=false, std::vector< unsigned int > target_component=std::vector< unsigned int >())
Definition: dof_tools.cc:1749
void map_dofs_to_support_points(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, std::vector< Point< spacedim > > &support_points)
Definition: dof_tools.cc:2077
void get_subdomain_association(const DoFHandlerType &dof_handler, std::vector< types::subdomain_id > &subdomain)
Definition: dof_tools.cc:1430
std::vector< types::global_dof_index > get_dofs_on_patch(const std::vector< typename DoFHandlerType::active_cell_iterator > &patch)
Definition: dof_tools.cc:2594
unsigned int boundary_id
Definition: types.h:110
void map_support_points_to_dofs(const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const DoFHandlerType &dof_handler, std::map< Point< DoFHandlerType::space_dimension >, types::global_dof_index, Comp > &point_to_index_map)
bool fe_is_primitive(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcFiniteElementsDontMatch()
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:2156
unsigned int count_dofs_with_subdomain_association(const DoFHandlerType &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1517
void map_dof_to_boundary_indices(const DoFHandlerType &dof_handler, std::vector< types::global_dof_index > &mapping)
Definition: dof_tools.cc:1916
Table< 2, Coupling > dof_couplings_from_component_couplings(const FiniteElement< dim, spacedim > &fe, const Table< 2, Coupling > &component_couplings)
void extract_locally_relevant_level_dofs(const DoFHandlerType &dof_handler, const unsigned int level, IndexSet &dof_set)
Definition: dof_tools.cc:1127
IndexSet dof_indices_with_subdomain_association(const DoFHandlerType &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1532
IndexSet extract_dofs_with_support_contained_within(const DoFHandlerType &dof_handler, const std::function< bool(const typename DoFHandlerType::active_cell_iterator &)> &predicate, const ConstraintMatrix &cm=ConstraintMatrix())
Definition: dof_tools.cc:788
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
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, ConstraintMatrix &constraints)
const hp::FECollection< dim, spacedim > & get_fe() const
Definition: tria.h:76
void make_hanging_node_constraints(const DoFHandlerType &dof_handler, ConstraintMatrix &constraints)
unsigned int max_dofs_per_face(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcGridNotCoarser()
unsigned int global_dof_index
Definition: types.h:88
void extract_locally_owned_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:1039
Abstract base class for mapping classes.
Definition: dof_tools.h:46
void make_cell_patches(SparsityPattern &block_list, const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, const std::vector< bool > &selected_dofs=std::vector< bool >(), const types::global_dof_index offset=0)
Definition: dof_tools.cc:2256
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
unsigned int max_dofs_per_cell(const DoFHandler< dim, spacedim > &dh)
#define DeclException0(Exception0)
Definition: exceptions.h:323
void make_flux_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern)
void extract_locally_relevant_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:1087
void extract_boundary_dofs(const DoFHandlerType &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:573
void extract_hanging_node_dofs(const DoFHandler< dim, spacedim > &dof_handler, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:981
void extract_dofs_with_support_on_boundary(const DoFHandlerType &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:707
unsigned int subdomain_id
Definition: types.h:42
std::vector< unsigned int > make_vertex_patches(SparsityPattern &block_list, const DoFHandlerType &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:2405
const hp::FECollection< dim, spacedim > & get_fe_collection() const
Definition: hp.h:102
void convert_couplings_to_blocks(const hp::DoFHandler< dim, spacedim > &dof_handler, const Table< 2, Coupling > &table_by_component, std::vector< Table< 2, Coupling > > &tables_by_block)
Definition: dof_tools.cc:2227
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
std::vector< IndexSet > locally_relevant_dofs_per_subdomain(const DoFHandlerType &dof_handler)
Definition: dof_tools.cc:1367
void extract_locally_active_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:1051
void extract_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:390
void extract_subdomain_dofs(const DoFHandlerType &dof_handler, const types::subdomain_id subdomain_id, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:1006
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
void make_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern, const ConstraintMatrix &constraints=ConstraintMatrix(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
unsigned int count_dofs_on_patch(const std::vector< typename DoFHandlerType::active_cell_iterator > &patch)
Definition: dof_tools.cc:2566
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
Definition: table.h:34
void extract_level_dofs(const unsigned int level, const DoFHandlerType &dof, const ComponentMask &component_mask, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:507
void make_zero_boundary_constraints(const DoFHandlerType< dim, spacedim > &dof, const types::boundary_id boundary_id, ConstraintMatrix &zero_boundary_constraints, const ComponentMask &component_mask=ComponentMask())
void make_periodicity_constraints(const FaceIterator &face_1, const typename identity< FaceIterator >::type &face_2, ::ConstraintMatrix &constraint_matrix, 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 >())
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
void get_active_fe_indices(const DoFHandlerType &dof_handler, std::vector< unsigned int > &active_fe_indices)
Definition: dof_tools.cc:1262
void distribute_cell_to_dof_vector(const DoFHandlerType &dof_handler, const Vector< Number > &cell_data, Vector< double > &dof_data, const unsigned int component=0)
Definition: dof_tools.cc:308
void extract_constant_modes(const DoFHandlerType &dof_handler, const ComponentMask &component_mask, std::vector< std::vector< bool > > &constant_modes)
Definition: dof_tools.cc:1175
static ::ExceptionBase & ExcNoFESelected()
void make_single_patch(SparsityPattern &block_list, const DoFHandlerType &dof_handler, const unsigned int level, const bool interior_dofs_only=false)
Definition: dof_tools.cc:2298