Reference documentation for deal.II version 8.5.1
dof_tools.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2017 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 Number> 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);
277  template <int dim, int spacedim>
278  unsigned int
279  max_dofs_per_cell (const DoFHandler<dim,spacedim> &dh);
280 
286  template <int dim, int spacedim>
287  unsigned int
288  max_dofs_per_cell (const hp::DoFHandler<dim,spacedim> &dh);
289 
290 
299  template <int dim, int spacedim>
300  unsigned int
301  max_dofs_per_face (const DoFHandler<dim,spacedim> &dh);
302 
311  template <int dim, int spacedim>
312  unsigned int
313  max_dofs_per_face (const hp::DoFHandler<dim,spacedim> &dh);
314 
323  template <int dim, int spacedim>
324  unsigned int
325  max_dofs_per_vertex (const DoFHandler<dim,spacedim> &dh);
326 
335  template <int dim, int spacedim>
336  unsigned int
337  max_dofs_per_vertex (const hp::DoFHandler<dim,spacedim> &dh);
338 
348  template <int dim, int spacedim>
349  unsigned int
350  n_components (const DoFHandler<dim,spacedim> &dh);
351 
361  template <int dim, int spacedim>
362  unsigned int
363  n_components (const hp::DoFHandler<dim,spacedim> &dh);
364 
374  template <int dim, int spacedim>
375  bool
376  fe_is_primitive (const DoFHandler<dim,spacedim> &dh);
377 
387  template <int dim, int spacedim>
388  bool
389  fe_is_primitive (const hp::DoFHandler<dim,spacedim> &dh);
390 
515  template <typename DoFHandlerType, typename SparsityPatternType>
516  void
517  make_sparsity_pattern (const DoFHandlerType &dof_handler,
518  SparsityPatternType &sparsity_pattern,
519  const ConstraintMatrix &constraints = ConstraintMatrix(),
520  const bool keep_constrained_dofs = true,
522 
588  template <typename DoFHandlerType, typename SparsityPatternType>
589  void
590  make_sparsity_pattern (const DoFHandlerType &dof_handler,
591  const Table<2, Coupling> &coupling,
592  SparsityPatternType &sparsity_pattern,
593  const ConstraintMatrix &constraints = ConstraintMatrix(),
594  const bool keep_constrained_dofs = true,
596 
617  template <typename DoFHandlerType, typename SparsityPatternType>
618  void
619  make_sparsity_pattern (const DoFHandlerType &dof_row,
620  const DoFHandlerType &dof_col,
621  SparsityPatternType &sparsity);
622 
668  template<typename DoFHandlerType, typename SparsityPatternType>
669  void
670  make_flux_sparsity_pattern (const DoFHandlerType &dof_handler,
671  SparsityPatternType &sparsity_pattern);
672 
681  template<typename DoFHandlerType, typename SparsityPatternType>
682  void
683  make_flux_sparsity_pattern (const DoFHandlerType &dof_handler,
684  SparsityPatternType &sparsity_pattern,
685  const ConstraintMatrix &constraints,
686  const bool keep_constrained_dofs = true,
688 
708  template <typename DoFHandlerType, typename SparsityPatternType>
709  void
710  make_flux_sparsity_pattern (const DoFHandlerType &dof,
711  SparsityPatternType &sparsity,
712  const Table<2,Coupling> &cell_integrals_mask,
713  const Table<2,Coupling> &face_integrals_mask,
723  template <typename DoFHandlerType, typename SparsityPatternType>
724  void
725  make_flux_sparsity_pattern (const DoFHandlerType &dof,
726  SparsityPatternType &sparsity,
727  const ConstraintMatrix &constraints,
728  const bool keep_constrained_dofs,
729  const Table<2,Coupling> &couplings,
730  const Table<2,Coupling> &face_couplings,
731  const types::subdomain_id subdomain_id);
732 
742  template <typename DoFHandlerType, typename SparsityPatternType>
743  void
744  make_boundary_sparsity_pattern (const DoFHandlerType &dof,
745  const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
746  SparsityPatternType &sparsity_pattern);
747 
764  template <typename DoFHandlerType, typename SparsityPatternType, typename number>
765  void
767  (const DoFHandlerType &dof,
768  const std::map<types::boundary_id, const Function<DoFHandlerType::space_dimension,number>*> &boundary_ids,
769  const std::vector<types::global_dof_index> &dof_to_boundary_mapping,
770  SparsityPatternType &sparsity);
771 
818  template <typename DoFHandlerType>
819  void
820  make_hanging_node_constraints (const DoFHandlerType &dof_handler,
821  ConstraintMatrix &constraints);
822 
890  template <int dim, int spacedim>
891  void
893  const unsigned int coarse_component,
894  const DoFHandler<dim,spacedim> &fine_grid,
895  const unsigned int fine_component,
896  const InterGridMap<DoFHandler<dim,spacedim> > &coarse_to_fine_grid_map,
897  ConstraintMatrix &constraints);
898 
899 
916  template <int dim, int spacedim>
917  void
919  const unsigned int coarse_component,
920  const DoFHandler<dim,spacedim> &fine_grid,
921  const unsigned int fine_component,
922  const InterGridMap<DoFHandler<dim,spacedim> > &coarse_to_fine_grid_map,
923  std::vector<std::map<types::global_dof_index, float> > &transfer_representation);
924 
1084  template<typename FaceIterator>
1085  void
1087  (const FaceIterator &face_1,
1088  const typename identity<FaceIterator>::type &face_2,
1089  ::ConstraintMatrix &constraint_matrix,
1090  const ComponentMask &component_mask = ComponentMask(),
1091  const bool face_orientation = true,
1092  const bool face_flip = false,
1093  const bool face_rotation = false,
1094  const FullMatrix<double> &matrix = FullMatrix<double>(),
1095  const std::vector<unsigned int> &first_vector_components = std::vector<unsigned int>());
1096 
1097 
1098 
1120  template<typename DoFHandlerType>
1121  void
1124  &periodic_faces,
1125  ::ConstraintMatrix &constraint_matrix,
1126  const ComponentMask &component_mask = ComponentMask(),
1127  const std::vector<unsigned int> &first_vector_components = std::vector<unsigned int>());
1128 
1129 
1130 
1161  template<typename DoFHandlerType>
1162  void
1164  (const DoFHandlerType &dof_handler,
1165  const types::boundary_id b_id1,
1166  const types::boundary_id b_id2,
1167  const int direction,
1168  ::ConstraintMatrix &constraint_matrix,
1169  const ComponentMask &component_mask = ComponentMask());
1170 
1171 
1172 
1197  template<typename DoFHandlerType>
1198  void
1200  (const DoFHandlerType &dof_handler,
1201  const types::boundary_id b_id,
1202  const int direction,
1203  ::ConstraintMatrix &constraint_matrix,
1204  const ComponentMask &component_mask = ComponentMask());
1205 
1222  template <int dim, int spacedim>
1223  void
1225  std::vector<bool> &selected_dofs);
1226 
1251  template <int dim, int spacedim>
1252  void
1253  extract_dofs (const DoFHandler<dim,spacedim> &dof_handler,
1254  const ComponentMask &component_mask,
1255  std::vector<bool> &selected_dofs);
1256 
1260  template <int dim, int spacedim>
1261  void
1262  extract_dofs (const hp::DoFHandler<dim,spacedim> &dof_handler,
1263  const ComponentMask &component_mask,
1264  std::vector<bool> &selected_dofs);
1265 
1284  template <int dim, int spacedim>
1285  void
1286  extract_dofs (const DoFHandler<dim,spacedim> &dof_handler,
1287  const BlockMask &block_mask,
1288  std::vector<bool> &selected_dofs);
1289 
1293  template <int dim, int spacedim>
1294  void
1295  extract_dofs (const hp::DoFHandler<dim,spacedim> &dof_handler,
1296  const BlockMask &block_mask,
1297  std::vector<bool> &selected_dofs);
1298 
1303  template <typename DoFHandlerType>
1304  void
1305  extract_level_dofs (const unsigned int level,
1306  const DoFHandlerType &dof,
1307  const ComponentMask &component_mask,
1308  std::vector<bool> &selected_dofs);
1309 
1314  template <typename DoFHandlerType>
1315  void
1316  extract_level_dofs (const unsigned int level,
1317  const DoFHandlerType &dof,
1318  const BlockMask &component_mask,
1319  std::vector<bool> &selected_dofs);
1320 
1372  template <typename DoFHandlerType>
1373  void
1374  extract_boundary_dofs (const DoFHandlerType &dof_handler,
1375  const ComponentMask &component_mask,
1376  std::vector<bool> &selected_dofs,
1377  const std::set<types::boundary_id> &boundary_ids = std::set<types::boundary_id>());
1378 
1409  template <typename DoFHandlerType>
1410  void
1411  extract_boundary_dofs (const DoFHandlerType &dof_handler,
1412  const ComponentMask &component_mask,
1413  IndexSet &selected_dofs,
1414  const std::set<types::boundary_id> &boundary_ids = std::set<types::boundary_id>());
1415 
1432  template <typename DoFHandlerType>
1433  void
1434  extract_dofs_with_support_on_boundary (const DoFHandlerType &dof_handler,
1435  const ComponentMask &component_mask,
1436  std::vector<bool> &selected_dofs,
1437  const std::set<types::boundary_id> &boundary_ids = std::set<types::boundary_id>());
1438 
1472  template <typename DoFHandlerType>
1473  void
1474  extract_constant_modes (const DoFHandlerType &dof_handler,
1475  const ComponentMask &component_mask,
1476  std::vector<std::vector<bool> > &constant_modes);
1478 
1492  template <typename DoFHandlerType>
1493  void
1494  extract_subdomain_dofs (const DoFHandlerType &dof_handler,
1495  const types::subdomain_id subdomain_id,
1496  std::vector<bool> &selected_dofs);
1497 
1498 
1505  template <typename DoFHandlerType>
1506  void
1507  extract_locally_owned_dofs (const DoFHandlerType &dof_handler,
1508  IndexSet &dof_set);
1509 
1510 
1525  template <typename DoFHandlerType>
1526  void
1527  extract_locally_active_dofs (const DoFHandlerType &dof_handler,
1528  IndexSet &dof_set);
1529 
1539  template <typename DoFHandlerType>
1540  void
1541  extract_locally_relevant_dofs (const DoFHandlerType &dof_handler,
1542  IndexSet &dof_set);
1543 
1560  template <typename DoFHandlerType>
1561  std::vector<IndexSet>
1562  locally_owned_dofs_per_subdomain (const DoFHandlerType &dof_handler);
1563 
1580  template <typename DoFHandlerType>
1581  std::vector<IndexSet>
1582  locally_relevant_dofs_per_subdomain (const DoFHandlerType &dof_handler);
1583 
1584 
1589  template <typename DoFHandlerType>
1590  void
1591  extract_locally_relevant_level_dofs (const DoFHandlerType &dof_handler,
1592  const unsigned int level,
1593  IndexSet &dof_set);
1594 
1595 
1626  template <typename DoFHandlerType>
1627  void
1628  get_subdomain_association (const DoFHandlerType &dof_handler,
1629  std::vector<types::subdomain_id> &subdomain);
1630 
1656  template <typename DoFHandlerType>
1657  unsigned int
1658  count_dofs_with_subdomain_association (const DoFHandlerType &dof_handler,
1659  const types::subdomain_id subdomain);
1660 
1680  template <typename DoFHandlerType>
1681  void
1682  count_dofs_with_subdomain_association (const DoFHandlerType &dof_handler,
1683  const types::subdomain_id subdomain,
1684  std::vector<unsigned int> &n_dofs_on_subdomain);
1685 
1706  template <typename DoFHandlerType>
1707  IndexSet
1708  dof_indices_with_subdomain_association (const DoFHandlerType &dof_handler,
1709  const types::subdomain_id subdomain);
1710  // @}
1720 
1772  template <typename DoFHandlerType>
1773  std::vector<types::global_dof_index>
1774  get_dofs_on_patch (const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
1775 
1796  template <int dim, int spacedim>
1797  void make_cell_patches(SparsityPattern &block_list,
1798  const DoFHandler<dim,spacedim> &dof_handler,
1799  const unsigned int level,
1800  const std::vector<bool> &selected_dofs = std::vector<bool>(),
1801  const types::global_dof_index offset = 0);
1802 
1853  template <typename DoFHandlerType>
1854  std::vector<unsigned int>
1856  const DoFHandlerType &dof_handler,
1857  const unsigned int level,
1858  const bool interior_dofs_only,
1859  const bool boundary_patches = false,
1860  const bool level_boundary_patches = false,
1861  const bool single_cell_patches = false,
1862  const bool invert_vertex_mapping = false);
1863 
1877  template <typename DoFHandlerType>
1878  std::vector<unsigned int>
1880  const DoFHandlerType &dof_handler,
1881  const unsigned int level,
1882  const BlockMask &exclude_boundary_dofs = BlockMask(),
1883  const bool boundary_patches = false,
1884  const bool level_boundary_patches = false,
1885  const bool single_cell_patches = false,
1886  const bool invert_vertex_mapping = false);
1887 
1925  template <typename DoFHandlerType>
1926  void make_child_patches(SparsityPattern &block_list,
1927  const DoFHandlerType &dof_handler,
1928  const unsigned int level,
1929  const bool interior_dofs_only,
1930  const bool boundary_dofs = false);
1931 
1951  template <typename DoFHandlerType>
1952  void make_single_patch(SparsityPattern &block_list,
1953  const DoFHandlerType &dof_handler,
1954  const unsigned int level,
1955  const bool interior_dofs_only = false);
1956 
1999  template <typename DoFHandlerType>
2000  void
2001  count_dofs_per_component (const DoFHandlerType &dof_handler,
2002  std::vector<types::global_dof_index> &dofs_per_component,
2003  const bool vector_valued_once = false,
2004  std::vector<unsigned int> target_component
2005  = std::vector<unsigned int>());
2006 
2023  template <typename DoFHandlerType>
2024  void
2025  count_dofs_per_block (const DoFHandlerType &dof,
2026  std::vector<types::global_dof_index> &dofs_per_block,
2027  const std::vector<unsigned int> &target_block
2028  = std::vector<unsigned int>());
2029 
2040  template <typename DoFHandlerType>
2041  void
2042  get_active_fe_indices (const DoFHandlerType &dof_handler,
2043  std::vector<unsigned int> &active_fe_indices);
2044 
2079  template <typename DoFHandlerType>
2080  unsigned int
2081  count_dofs_on_patch (const std::vector<typename DoFHandlerType::active_cell_iterator> &patch);
2082 
2105  template <typename DoFHandlerType>
2106  void
2107  map_dof_to_boundary_indices (const DoFHandlerType &dof_handler,
2108  std::vector<types::global_dof_index> &mapping);
2109 
2120  template <typename DoFHandlerType>
2121  void
2122  map_dof_to_boundary_indices (const DoFHandlerType &dof_handler,
2123  const std::set<types::boundary_id> &boundary_ids,
2124  std::vector<types::global_dof_index> &mapping);
2125 
2143  template <int dim, int spacedim>
2144  void
2146  const DoFHandler<dim,spacedim> &dof_handler,
2147  std::vector<Point<spacedim> > &support_points);
2148 
2153  template <int dim, int spacedim>
2154  void
2155  map_dofs_to_support_points (const ::hp::MappingCollection<dim,spacedim> &mapping,
2156  const hp::DoFHandler<dim,spacedim> &dof_handler,
2157  std::vector<Point<spacedim> > &support_points);
2158 
2186  template <int dim, int spacedim>
2187  void
2189  const DoFHandler<dim,spacedim> &dof_handler,
2190  std::map<types::global_dof_index, Point<spacedim> > &support_points);
2191 
2195  template <int dim, int spacedim>
2196  void
2197  map_dofs_to_support_points (const ::hp::MappingCollection<dim,spacedim> &mapping,
2198  const hp::DoFHandler<dim,spacedim> &dof_handler,
2199  std::map<types::global_dof_index, Point<spacedim> > &support_points);
2200 
2201 
2222  template <typename DoFHandlerType, class Comp>
2223  void
2226  const DoFHandlerType &dof_handler,
2227  std::map<Point<DoFHandlerType::space_dimension>, types::global_dof_index, Comp> &point_to_index_map);
2265  template <typename DoFHandlerType, typename Number>
2266  void
2267  distribute_cell_to_dof_vector (const DoFHandlerType &dof_handler,
2268  const Vector<Number> &cell_data,
2269  Vector<double> &dof_data,
2270  const unsigned int component = 0);
2271 
2272 
2298  template <int spacedim>
2299  void
2300  write_gnuplot_dof_support_point_info(std::ostream &out,
2301  const std::map<types::global_dof_index, Point<spacedim> > &support_points);
2302 
2303 
2344  template <int dim, int spacedim, template <int, int> class DoFHandlerType>
2345  void
2346  make_zero_boundary_constraints (const DoFHandlerType<dim,spacedim> &dof,
2347  const types::boundary_id boundary_id,
2348  ConstraintMatrix &zero_boundary_constraints,
2349  const ComponentMask &component_mask = ComponentMask());
2350 
2361  template <int dim, int spacedim, template <int, int> class DoFHandlerType>
2362  void
2363  make_zero_boundary_constraints (const DoFHandlerType<dim,spacedim> &dof,
2364  ConstraintMatrix &zero_boundary_constraints,
2365  const ComponentMask &component_mask = ComponentMask());
2366 
2411 }
2412 
2413 
2414 
2415 /* ------------------------- inline functions -------------- */
2416 
2417 #ifndef DOXYGEN
2418 
2419 namespace DoFTools
2420 {
2426  inline
2428  const Coupling c2)
2429  {
2430  if (c2 == always)
2431  c1 = always;
2432  else if (c1 != always && c2 == nonzero)
2433  return c1 = nonzero;
2434  return c1;
2435  }
2436 
2437 
2443  inline
2444  Coupling operator | (const Coupling c1,
2445  const Coupling c2)
2446  {
2447  if (c1 == always || c2 == always)
2448  return always;
2449  if (c1 == nonzero || c2 == nonzero)
2450  return nonzero;
2451  return none;
2452  }
2453 
2454 
2455 // ---------------------- inline and template functions --------------------
2456 
2457  template <int dim, int spacedim>
2458  inline
2459  unsigned int
2460  max_dofs_per_cell (const DoFHandler<dim,spacedim> &dh)
2461  {
2462  return dh.get_fe().dofs_per_cell;
2463  }
2464 
2465 
2466  template <int dim, int spacedim>
2467  inline
2468  unsigned int
2469  max_dofs_per_face (const DoFHandler<dim,spacedim> &dh)
2470  {
2471  return dh.get_fe().dofs_per_face;
2472  }
2473 
2474 
2475  template <int dim, int spacedim>
2476  inline
2477  unsigned int
2478  max_dofs_per_vertex (const DoFHandler<dim,spacedim> &dh)
2479  {
2480  return dh.get_fe().dofs_per_vertex;
2481  }
2482 
2483 
2484  template <int dim, int spacedim>
2485  inline
2486  unsigned int
2487  n_components (const DoFHandler<dim,spacedim> &dh)
2488  {
2489  return dh.get_fe().n_components();
2490  }
2491 
2492 
2493 
2494  template <int dim, int spacedim>
2495  inline
2496  bool
2497  fe_is_primitive (const DoFHandler<dim,spacedim> &dh)
2498  {
2499  return dh.get_fe().is_primitive();
2500  }
2501 
2502 
2503  template <int dim, int spacedim>
2504  inline
2505  unsigned int
2506  max_dofs_per_cell (const hp::DoFHandler<dim,spacedim> &dh)
2507  {
2508  return dh.get_fe().max_dofs_per_cell ();
2509  }
2510 
2511 
2512  template <int dim, int spacedim>
2513  inline
2514  unsigned int
2515  max_dofs_per_face (const hp::DoFHandler<dim,spacedim> &dh)
2516  {
2517  return dh.get_fe().max_dofs_per_face ();
2518  }
2519 
2520 
2521  template <int dim, int spacedim>
2522  inline
2523  unsigned int
2524  max_dofs_per_vertex (const hp::DoFHandler<dim,spacedim> &dh)
2525  {
2526  return dh.get_fe().max_dofs_per_vertex ();
2527  }
2528 
2529 
2530  template <int dim, int spacedim>
2531  inline
2532  unsigned int
2533  n_components (const hp::DoFHandler<dim,spacedim> &dh)
2534  {
2535  return dh.get_fe()[0].n_components();
2536  }
2537 
2538 
2539  template <int dim, int spacedim>
2540  inline
2541  bool
2542  fe_is_primitive (const hp::DoFHandler<dim,spacedim> &dh)
2543  {
2544  return dh.get_fe()[0].is_primitive();
2545  }
2546 
2547 
2548  template <typename DoFHandlerType, class Comp>
2549  void
2551  (
2553  const DoFHandlerType &dof_handler,
2554  std::map<Point<DoFHandlerType::space_dimension>, types::global_dof_index, Comp> &point_to_index_map)
2555  {
2556  // let the checking of arguments be
2557  // done by the function first
2558  // called
2559  std::vector<Point<DoFHandlerType::space_dimension> > support_points (dof_handler.n_dofs());
2560  map_dofs_to_support_points (mapping, dof_handler, support_points);
2561  // now copy over the results of the
2562  // previous function into the
2563  // output arg
2564  point_to_index_map.clear ();
2565  for (types::global_dof_index i=0; i<dof_handler.n_dofs(); ++i)
2566  point_to_index_map[support_points[i]] = i;
2567  }
2568 }
2569 
2570 #endif
2571 
2572 DEAL_II_NAMESPACE_CLOSE
2573 
2574 #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:2224
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:1166
const types::subdomain_id invalid_subdomain_id
Definition: types.h:245
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:1632
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:1960
void get_subdomain_association(const DoFHandlerType &dof_handler, std::vector< types::subdomain_id > &subdomain)
Definition: dof_tools.cc:1309
std::vector< types::global_dof_index > get_dofs_on_patch(const std::vector< typename DoFHandlerType::active_cell_iterator > &patch)
Definition: dof_tools.cc:2477
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)
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:2039
unsigned int count_dofs_with_subdomain_association(const DoFHandlerType &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1400
const FiniteElement< dim, spacedim > & get_fe() const
void map_dof_to_boundary_indices(const DoFHandlerType &dof_handler, std::vector< types::global_dof_index > &mapping)
Definition: dof_tools.cc:1799
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:1019
IndexSet dof_indices_with_subdomain_association(const DoFHandlerType &dof_handler, const types::subdomain_id subdomain)
Definition: dof_tools.cc:1415
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)
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:932
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:2139
#define DeclException0(Exception0)
Definition: exceptions.h:541
void make_flux_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern)
UpdateFlags operator|(UpdateFlags f1, UpdateFlags f2)
void extract_locally_relevant_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:980
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:555
void extract_hanging_node_dofs(const DoFHandler< dim, spacedim > &dof_handler, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:888
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:689
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:2288
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:2110
std::vector< IndexSet > locally_relevant_dofs_per_subdomain(const DoFHandlerType &dof_handler)
Definition: dof_tools.cc:1246
void extract_locally_active_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:944
void extract_dofs(const DoFHandler< dim, spacedim > &dof_handler, const ComponentMask &component_mask, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:372
void extract_subdomain_dofs(const DoFHandlerType &dof_handler, const types::subdomain_id subdomain_id, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:899
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:1719
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:2449
UpdateFlags & operator|=(UpdateFlags &f1, UpdateFlags f2)
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
Definition: table.h:33
unsigned char boundary_id
Definition: types.h:110
void extract_level_dofs(const unsigned int level, const DoFHandlerType &dof, const ComponentMask &component_mask, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:489
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 >())
void get_active_fe_indices(const DoFHandlerType &dof_handler, std::vector< unsigned int > &active_fe_indices)
Definition: dof_tools.cc:1152
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:290
void extract_constant_modes(const DoFHandlerType &dof_handler, const ComponentMask &component_mask, std::vector< std::vector< bool > > &constant_modes)
Definition: dof_tools.cc:1065
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:2181