Reference documentation for deal.II version GIT 92344eb39a 2023-03-23 02:40:02+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\}}\)
dof_info.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2011 - 2022 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii_matrix_free_dof_info_h
18 #define dealii_matrix_free_dof_info_h
19 
20 
21 #include <deal.II/base/config.h>
22 
25 
29 
30 #include <array>
31 #include <memory>
32 
33 
35 
36 #ifndef DOXYGEN
37 
38 // forward declarations
39 
40 namespace internal
41 {
42  namespace MatrixFreeFunctions
43  {
44  template <int dim>
45  class HangingNodes;
46 
47  struct TaskInfo;
48 
49  template <typename Number>
50  struct ConstraintValues;
51  } // namespace MatrixFreeFunctions
52 } // namespace internal
53 
54 template <typename>
55 class AffineConstraints;
56 
58 
59 template <typename>
60 class TriaIterator;
61 
62 template <int, int, bool>
63 class DoFCellAccessor;
64 
65 namespace Utilities
66 {
67  namespace MPI
68  {
69  class Partitioner;
70  }
71 } // namespace Utilities
72 
73 #endif
74 
75 namespace internal
76 {
77  namespace MatrixFreeFunctions
78  {
83  using compressed_constraint_kind = std::uint8_t;
84 
103  struct DoFInfo
104  {
117  static const unsigned int chunk_size_zero_vector = 64;
118 
122  DoFInfo();
123 
127  DoFInfo(const DoFInfo &) = default;
128 
132  DoFInfo(DoFInfo &&) noexcept = default;
133 
137  ~DoFInfo() = default;
138 
142  DoFInfo &
143  operator=(const DoFInfo &) = default;
144 
148  DoFInfo &
149  operator=(DoFInfo &&) noexcept = default;
150 
154  void
155  clear();
156 
162  unsigned int
163  fe_index_from_degree(const unsigned int first_selected_component,
164  const unsigned int fe_degree) const;
165 
185  void
186  get_dof_indices_on_cell_batch(std::vector<unsigned int> &local_indices,
187  const unsigned int cell_batch,
188  const bool with_constraints = true) const;
189 
199  template <typename number>
200  void
202  const std::vector<types::global_dof_index> &local_indices_resolved,
203  const std::vector<types::global_dof_index> &local_indices,
204  const bool cell_has_hanging_nodes,
205  const ::AffineConstraints<number> & constraints,
206  const unsigned int cell_number,
207  ConstraintValues<double> & constraint_values,
208  bool & cell_at_boundary);
209 
214  template <int dim>
215  bool
217  const HangingNodes<dim> & hanging_nodes,
218  const std::vector<std::vector<unsigned int>> &lexicographic_mapping,
219  const unsigned int cell_number,
220  const TriaIterator<DoFCellAccessor<dim, dim, false>> &cell,
221  std::vector<types::global_dof_index> & dof_indices);
222 
230  void
231  assign_ghosts(const std::vector<unsigned int> &boundary_cells,
232  const MPI_Comm & communicator_sm,
233  const bool use_vector_data_exchanger_full);
234 
241  void
242  reorder_cells(const TaskInfo & task_info,
243  const std::vector<unsigned int> & renumbering,
244  const std::vector<unsigned int> & constraint_pool_row_index,
245  const std::vector<unsigned char> &irregular_cells);
246 
251  void
253  const std::vector<unsigned char> &irregular_cells);
254 
259  template <int length>
260  void
262  const std::vector<FaceToCellTopology<length>> &faces);
263 
269  void
270  make_connectivity_graph(const TaskInfo & task_info,
271  const std::vector<unsigned int> &renumbering,
272  DynamicSparsityPattern &connectivity) const;
273 
279  void
281  const Table<2, ShapeInfo<double>> & shape_info,
282  const unsigned int n_owned_cells,
283  const unsigned int n_lanes,
284  const std::vector<FaceToCellTopology<1>> &inner_faces,
285  const std::vector<FaceToCellTopology<1>> &ghosted_faces,
286  const bool fill_cell_centric,
287  const MPI_Comm & communicator_sm,
288  const bool use_vector_data_exchanger_full);
289 
295  void
297  std::array<std::vector<std::pair<unsigned int, unsigned int>>, 3>
298  &cell_indices_contiguous_sm);
299 
310  void
312  std::vector<types::global_dof_index> &renumbering);
313 
323  template <int length>
324  void
326  const TaskInfo & task_info,
327  const std::vector<FaceToCellTopology<length>> &faces);
328 
332  std::size_t
333  memory_consumption() const;
334 
339  template <typename StreamType>
340  void
341  print_memory_consumption(StreamType & out,
342  const TaskInfo &size_info) const;
343 
348  template <typename Number>
349  void
350  print(const std::vector<Number> & constraint_pool_data,
351  const std::vector<unsigned int> &constraint_pool_row_index,
352  std::ostream & out) const;
353 
365  enum class IndexStorageVariants : unsigned char
366  {
375  full,
387  interleaved,
396  contiguous,
420  interleaved_contiguous,
432  interleaved_contiguous_strided,
445  interleaved_contiguous_mixed_strides
446  };
447 
452  enum DoFAccessIndex : unsigned char
453  {
465  dof_access_cell = 2
466  };
467 
476  unsigned int vectorization_length;
477 
485  std::array<std::vector<IndexStorageVariants>, 3> index_storage_variants;
486 
493  std::vector<std::pair<unsigned int, unsigned int>> row_starts;
494 
510  std::vector<unsigned int> dof_indices;
511 
516  std::vector<std::vector<bool>> hanging_node_constraint_masks_comp;
517 
522  std::vector<compressed_constraint_kind> hanging_node_constraint_masks;
523 
533  std::vector<std::pair<unsigned short, unsigned short>>
535 
539  std::vector<unsigned int> dof_indices_interleaved;
540 
549  std::array<std::vector<unsigned int>, 3> dof_indices_contiguous;
550 
559  std::array<std::vector<std::pair<unsigned int, unsigned int>>, 3>
561 
570  std::array<std::vector<unsigned int>, 3> dof_indices_interleave_strides;
571 
581  std::array<std::vector<unsigned char>, 3> n_vectorization_lanes_filled;
582 
589  std::shared_ptr<const Utilities::MPI::Partitioner> vector_partitioner;
590 
594  std::shared_ptr<
597 
617  std::array<
618  std::shared_ptr<
620  5>
622 
627  std::vector<unsigned int> constrained_dofs;
628 
633  std::vector<unsigned int> row_starts_plain_indices;
634 
643  std::vector<unsigned int> plain_dof_indices;
644 
650 
655  unsigned int n_base_elements;
656 
661  std::vector<unsigned int> n_components;
662 
667  std::vector<unsigned int> start_components;
668 
673  std::vector<unsigned int> component_to_base_index;
674 
686  std::vector<std::vector<unsigned int>> component_dof_indices_offset;
687 
691  std::vector<unsigned int> dofs_per_cell;
692 
696  std::vector<unsigned int> dofs_per_face;
697 
702 
706  std::vector<unsigned int> cell_active_fe_index;
707 
712  unsigned int max_fe_index;
713 
719  std::vector<std::vector<unsigned int>> fe_index_conversion;
720 
726  std::vector<types::global_dof_index> ghost_dofs;
727 
733  std::vector<unsigned int> vector_zero_range_list_index;
734 
738  std::vector<std::pair<unsigned int, unsigned int>> vector_zero_range_list;
739 
745  std::vector<unsigned int> cell_loop_pre_list_index;
746 
751  std::vector<std::pair<unsigned int, unsigned int>> cell_loop_pre_list;
752 
758  std::vector<unsigned int> cell_loop_post_list_index;
759 
764  std::vector<std::pair<unsigned int, unsigned int>> cell_loop_post_list;
765  };
766 
767 
768  /*-------------------------- Inline functions ---------------------------*/
769 
770 #ifndef DOXYGEN
771 
772 
773  inline unsigned int
774  DoFInfo::fe_index_from_degree(const unsigned int first_selected_component,
775  const unsigned int fe_degree) const
776  {
777  const unsigned int n_indices = fe_index_conversion.size();
778  if (n_indices <= 1)
779  return 0;
780  for (unsigned int i = 0; i < n_indices; ++i)
781  if (fe_index_conversion[i][first_selected_component] == fe_degree)
782  return i;
784  }
785 
786 #endif // ifndef DOXYGEN
787 
788  } // end of namespace MatrixFreeFunctions
789 } // end of namespace internal
790 
792 
793 #endif
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
std::uint8_t compressed_constraint_kind
Definition: dof_info.h:83
static const unsigned int invalid_unsigned_int
Definition: types.h:213
Definition: types.h:33
unsigned int global_dof_index
Definition: types.h:82
void print(const std::vector< Number > &constraint_pool_data, const std::vector< unsigned int > &constraint_pool_row_index, std::ostream &out) const
void reorder_cells(const TaskInfo &task_info, const std::vector< unsigned int > &renumbering, const std::vector< unsigned int > &constraint_pool_row_index, const std::vector< unsigned char > &irregular_cells)
Definition: dof_info.cc:298
std::vector< std::pair< unsigned short, unsigned short > > constraint_indicator
Definition: dof_info.h:534
std::vector< unsigned int > cell_loop_pre_list_index
Definition: dof_info.h:745
std::vector< std::pair< unsigned int, unsigned int > > row_starts
Definition: dof_info.h:493
std::vector< std::vector< unsigned int > > component_dof_indices_offset
Definition: dof_info.h:686
std::vector< std::vector< bool > > hanging_node_constraint_masks_comp
Definition: dof_info.h:516
void get_dof_indices_on_cell_batch(std::vector< unsigned int > &local_indices, const unsigned int cell_batch, const bool with_constraints=true) const
Definition: dof_info.cc:82
unsigned int fe_index_from_degree(const unsigned int first_selected_component, const unsigned int fe_degree) const
void compute_cell_index_compression(const std::vector< unsigned char > &irregular_cells)
Definition: dof_info.cc:530
std::vector< unsigned int > cell_loop_post_list_index
Definition: dof_info.h:758
std::vector< unsigned int > dofs_per_cell
Definition: dof_info.h:691
void compute_shared_memory_contiguous_indices(std::array< std::vector< std::pair< unsigned int, unsigned int >>, 3 > &cell_indices_contiguous_sm)
Definition: dof_info.cc:1243
void compute_dof_renumbering(std::vector< types::global_dof_index > &renumbering)
Definition: dof_info.cc:1469
static const unsigned int chunk_size_zero_vector
Definition: dof_info.h:117
std::shared_ptr< const internal::MatrixFreeFunctions::VectorDataExchange::Base > vector_exchanger
Definition: dof_info.h:596
std::vector< std::vector< unsigned int > > fe_index_conversion
Definition: dof_info.h:719
std::vector< unsigned int > dof_indices
Definition: dof_info.h:510
std::vector< std::pair< unsigned int, unsigned int > > vector_zero_range_list
Definition: dof_info.h:738
std::shared_ptr< const Utilities::MPI::Partitioner > vector_partitioner
Definition: dof_info.h:589
std::vector< compressed_constraint_kind > hanging_node_constraint_masks
Definition: dof_info.h:522
std::array< std::vector< unsigned int >, 3 > dof_indices_interleave_strides
Definition: dof_info.h:570
std::array< std::vector< std::pair< unsigned int, unsigned int > >, 3 > dof_indices_contiguous_sm
Definition: dof_info.h:560
std::vector< unsigned int > row_starts_plain_indices
Definition: dof_info.h:633
std::vector< unsigned int > dofs_per_face
Definition: dof_info.h:696
std::size_t memory_consumption() const
Definition: dof_info.cc:1527
std::vector< unsigned int > component_to_base_index
Definition: dof_info.h:673
std::vector< unsigned int > n_components
Definition: dof_info.h:661
void print_memory_consumption(StreamType &out, const TaskInfo &size_info) const
void assign_ghosts(const std::vector< unsigned int > &boundary_cells, const MPI_Comm &communicator_sm, const bool use_vector_data_exchanger_full)
Definition: dof_info.cc:155
std::vector< std::pair< unsigned int, unsigned int > > cell_loop_pre_list
Definition: dof_info.h:751
std::vector< types::global_dof_index > ghost_dofs
Definition: dof_info.h:726
void read_dof_indices(const std::vector< types::global_dof_index > &local_indices_resolved, const std::vector< types::global_dof_index > &local_indices, const bool cell_has_hanging_nodes, const ::AffineConstraints< number > &constraints, const unsigned int cell_number, ConstraintValues< double > &constraint_values, bool &cell_at_boundary)
std::array< std::vector< unsigned int >, 3 > dof_indices_contiguous
Definition: dof_info.h:549
void compute_tight_partitioners(const Table< 2, ShapeInfo< double >> &shape_info, const unsigned int n_owned_cells, const unsigned int n_lanes, const std::vector< FaceToCellTopology< 1 >> &inner_faces, const std::vector< FaceToCellTopology< 1 >> &ghosted_faces, const bool fill_cell_centric, const MPI_Comm &communicator_sm, const bool use_vector_data_exchanger_full)
Definition: dof_info.cc:888
DoFInfo(const DoFInfo &)=default
std::vector< unsigned int > cell_active_fe_index
Definition: dof_info.h:706
std::array< std::shared_ptr< const internal::MatrixFreeFunctions::VectorDataExchange::Base >, 5 > vector_exchanger_face_variants
Definition: dof_info.h:621
std::vector< unsigned int > plain_dof_indices
Definition: dof_info.h:643
std::array< std::vector< unsigned char >, 3 > n_vectorization_lanes_filled
Definition: dof_info.h:581
std::vector< unsigned int > start_components
Definition: dof_info.h:667
std::vector< unsigned int > dof_indices_interleaved
Definition: dof_info.h:539
std::vector< unsigned int > constrained_dofs
Definition: dof_info.h:627
std::array< std::vector< IndexStorageVariants >, 3 > index_storage_variants
Definition: dof_info.h:485
bool process_hanging_node_constraints(const HangingNodes< dim > &hanging_nodes, const std::vector< std::vector< unsigned int >> &lexicographic_mapping, const unsigned int cell_number, const TriaIterator< DoFCellAccessor< dim, dim, false >> &cell, std::vector< types::global_dof_index > &dof_indices)
void make_connectivity_graph(const TaskInfo &task_info, const std::vector< unsigned int > &renumbering, DynamicSparsityPattern &connectivity) const
Definition: dof_info.cc:1392
std::vector< unsigned int > vector_zero_range_list_index
Definition: dof_info.h:733
void compute_face_index_compression(const std::vector< FaceToCellTopology< length >> &faces)
void compute_vector_zero_access_pattern(const TaskInfo &task_info, const std::vector< FaceToCellTopology< length >> &faces)
std::vector< std::pair< unsigned int, unsigned int > > cell_loop_post_list
Definition: dof_info.h:764
DoFInfo(DoFInfo &&) noexcept=default