Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
dof_info.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2011 - 2020 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 
26 
28 
31 
35 
36 #include <array>
37 #include <memory>
38 
39 
41 
42 namespace internal
43 {
44  namespace MatrixFreeFunctions
45  {
50  template <typename Number>
52  {
54 
62  template <typename number2>
63  unsigned short
65  const std::vector<std::pair<types::global_dof_index, number2>>
66  &entries);
67 
68  std::vector<std::pair<types::global_dof_index, double>>
70  std::vector<types::global_dof_index> constraint_indices;
71 
72  std::pair<std::vector<Number>, types::global_dof_index> next_constraint;
73  std::map<std::vector<Number>,
77  };
78 
99  struct DoFInfo
100  {
113  static const unsigned int chunk_size_zero_vector = 64;
114 
118  DoFInfo();
119 
123  DoFInfo(const DoFInfo &) = default;
124 
128  void
129  clear();
130 
136  unsigned int
137  fe_index_from_degree(const unsigned int first_selected_component,
138  const unsigned int fe_degree) const;
139 
159  void
160  get_dof_indices_on_cell_batch(std::vector<unsigned int> &locall_indices,
161  const unsigned int cell,
162  const bool with_constraints = true) const;
163 
172  template <typename number>
173  void
175  const std::vector<types::global_dof_index> &local_indices,
176  const std::vector<unsigned int> & lexicographic_inv,
177  const AffineConstraints<number> & constraints,
178  const unsigned int cell_number,
179  ConstraintValues<double> & constraint_values,
180  bool & cell_at_boundary);
181 
189  void
190  assign_ghosts(const std::vector<unsigned int> &boundary_cells);
191 
198  void
199  reorder_cells(const TaskInfo & task_info,
200  const std::vector<unsigned int> & renumbering,
201  const std::vector<unsigned int> & constraint_pool_row_index,
202  const std::vector<unsigned char> &irregular_cells);
203 
208  void
210  const std::vector<unsigned char> &irregular_cells);
211 
216  template <int length>
217  void
219  const std::vector<FaceToCellTopology<length>> &faces);
220 
226  void
227  make_connectivity_graph(const TaskInfo & task_info,
228  const std::vector<unsigned int> &renumbering,
229  DynamicSparsityPattern &connectivity) const;
230 
241  void
243  std::vector<types::global_dof_index> &renumbering);
244 
254  template <int length>
255  void
257  const TaskInfo & task_info,
258  const std::vector<FaceToCellTopology<length>> &faces);
259 
263  std::size_t
264  memory_consumption() const;
265 
270  template <typename StreamType>
271  void
272  print_memory_consumption(StreamType & out,
273  const TaskInfo &size_info) const;
274 
279  template <typename Number>
280  void
281  print(const std::vector<Number> & constraint_pool_data,
282  const std::vector<unsigned int> &constraint_pool_row_index,
283  std::ostream & out) const;
284 
296  enum class IndexStorageVariants : unsigned char
297  {
306  full,
318  interleaved,
327  contiguous,
377  };
378 
383  enum DoFAccessIndex : unsigned char
384  {
397  };
398 
404  unsigned int dimension;
405 
414  unsigned int vectorization_length;
415 
423  std::vector<IndexStorageVariants> index_storage_variants[3];
424 
431  std::vector<std::pair<unsigned int, unsigned int>> row_starts;
432 
448  std::vector<unsigned int> dof_indices;
449 
459  std::vector<std::pair<unsigned short, unsigned short>>
461 
465  std::vector<unsigned int> dof_indices_interleaved;
466 
475  std::vector<unsigned int> dof_indices_contiguous[3];
476 
485  std::vector<unsigned int> dof_indices_interleave_strides[3];
486 
496  std::vector<unsigned char> n_vectorization_lanes_filled[3];
497 
504  std::shared_ptr<const Utilities::MPI::Partitioner> vector_partitioner;
505 
524  std::array<std::shared_ptr<const Utilities::MPI::Partitioner>, 5>
526 
531  std::vector<unsigned int> constrained_dofs;
532 
537  std::vector<unsigned int> row_starts_plain_indices;
538 
547  std::vector<unsigned int> plain_dof_indices;
548 
554 
559  unsigned int n_base_elements;
560 
565  std::vector<unsigned int> n_components;
566 
571  std::vector<unsigned int> start_components;
572 
577  std::vector<unsigned int> component_to_base_index;
578 
590  std::vector<std::vector<unsigned int>> component_dof_indices_offset;
591 
595  std::vector<unsigned int> dofs_per_cell;
596 
600  std::vector<unsigned int> dofs_per_face;
601 
606 
610  std::vector<unsigned int> cell_active_fe_index;
611 
616  unsigned int max_fe_index;
617 
623  std::vector<std::vector<unsigned int>> fe_index_conversion;
624 
630  std::vector<types::global_dof_index> ghost_dofs;
631 
637  std::vector<unsigned int> vector_zero_range_list_index;
638 
642  std::vector<std::pair<unsigned int, unsigned int>> vector_zero_range_list;
643 
649  std::vector<unsigned int> cell_loop_pre_list_index;
650 
655  std::vector<std::pair<unsigned int, unsigned int>> cell_loop_pre_list;
656 
662  std::vector<unsigned int> cell_loop_post_list_index;
663 
668  std::vector<std::pair<unsigned int, unsigned int>> cell_loop_post_list;
669  };
670 
671 
672  /*-------------------------- Inline functions ---------------------------*/
673 
674 #ifndef DOXYGEN
675 
676 
677  inline unsigned int
678  DoFInfo::fe_index_from_degree(const unsigned int first_selected_component,
679  const unsigned int fe_degree) const
680  {
681  const unsigned int n_indices = fe_index_conversion.size();
682  if (n_indices <= 1)
683  return 0;
684  for (unsigned int i = 0; i < n_indices; ++i)
685  if (fe_index_conversion[i][first_selected_component] == fe_degree)
686  return i;
688  }
689 
690 #endif // ifndef DOXYGEN
691 
692  } // end of namespace MatrixFreeFunctions
693 } // end of namespace internal
694 
696 
697 #endif
internal::MatrixFreeFunctions::DoFInfo::dof_indices_contiguous
std::vector< unsigned int > dof_indices_contiguous[3]
Definition: dof_info.h:475
internal::MatrixFreeFunctions::DoFInfo::IndexStorageVariants::interleaved_contiguous
@ interleaved_contiguous
internal::MatrixFreeFunctions::DoFInfo::IndexStorageVariants::interleaved_contiguous_mixed_strides
@ interleaved_contiguous_mixed_strides
dynamic_sparsity_pattern.h
internal::MatrixFreeFunctions::ConstraintValues::constraints
std::map< std::vector< Number >, types::global_dof_index, FPArrayComparator< Number > > constraints
Definition: dof_info.h:76
internal::MatrixFreeFunctions::DoFInfo::dof_indices_interleaved
std::vector< unsigned int > dof_indices_interleaved
Definition: dof_info.h:465
internal::MatrixFreeFunctions::DoFInfo::component_dof_indices_offset
std::vector< std::vector< unsigned int > > component_dof_indices_offset
Definition: dof_info.h:590
internal::MatrixFreeFunctions::DoFInfo::dof_access_cell
@ dof_access_cell
Definition: dof_info.h:396
internal::MatrixFreeFunctions::DoFInfo::constrained_dofs
std::vector< unsigned int > constrained_dofs
Definition: dof_info.h:531
internal::MatrixFreeFunctions::DoFInfo::dof_access_face_interior
@ dof_access_face_interior
Definition: dof_info.h:388
vectorization.h
internal::MatrixFreeFunctions::ConstraintValues::ConstraintValues
ConstraintValues()
internal::MatrixFreeFunctions::FaceToCellTopology
Definition: face_info.h:56
internal::MatrixFreeFunctions::DoFInfo::dof_indices_interleave_strides
std::vector< unsigned int > dof_indices_interleave_strides[3]
Definition: dof_info.h:485
internal::MatrixFreeFunctions::DoFInfo::memory_consumption
std::size_t memory_consumption() const
internal::MatrixFreeFunctions::DoFInfo::cell_active_fe_index
std::vector< unsigned int > cell_active_fe_index
Definition: dof_info.h:610
internal::MatrixFreeFunctions::DoFInfo::vector_partitioner_face_variants
std::array< std::shared_ptr< const Utilities::MPI::Partitioner >, 5 > vector_partitioner_face_variants
Definition: dof_info.h:525
partitioner.h
internal::MatrixFreeFunctions::DoFInfo::component_to_base_index
std::vector< unsigned int > component_to_base_index
Definition: dof_info.h:577
internal::MatrixFreeFunctions::DoFInfo::dof_indices
std::vector< unsigned int > dof_indices
Definition: dof_info.h:448
face_info.h
internal::MatrixFreeFunctions::DoFInfo::IndexStorageVariants::contiguous
@ contiguous
internal::MatrixFreeFunctions::ConstraintValues::constraint_indices
std::vector< types::global_dof_index > constraint_indices
Definition: dof_info.h:70
internal::MatrixFreeFunctions::DoFInfo::compute_dof_renumbering
void compute_dof_renumbering(std::vector< types::global_dof_index > &renumbering)
internal::MatrixFreeFunctions::DoFInfo::IndexStorageVariants::interleaved_contiguous_strided
@ interleaved_contiguous_strided
internal::MatrixFreeFunctions::DoFInfo::get_dof_indices_on_cell_batch
void get_dof_indices_on_cell_batch(std::vector< unsigned int > &locall_indices, const unsigned int cell, const bool with_constraints=true) const
internal::MatrixFreeFunctions::DoFInfo
Definition: dof_info.h:99
internal::MatrixFreeFunctions::DoFInfo::DoFInfo
DoFInfo()
internal::MatrixFreeFunctions::DoFInfo::n_vectorization_lanes_filled
std::vector< unsigned char > n_vectorization_lanes_filled[3]
Definition: dof_info.h:496
internal::MatrixFreeFunctions::ConstraintValues::constraint_entries
std::vector< std::pair< types::global_dof_index, double > > constraint_entries
Definition: dof_info.h:69
internal::MatrixFreeFunctions::DoFInfo::constraint_indicator
std::vector< std::pair< unsigned short, unsigned short > > constraint_indicator
Definition: dof_info.h:460
internal::MatrixFreeFunctions::TaskInfo
Definition: task_info.h:105
internal::MatrixFreeFunctions::DoFInfo::vector_partitioner
std::shared_ptr< const Utilities::MPI::Partitioner > vector_partitioner
Definition: dof_info.h:504
internal::MatrixFreeFunctions::DoFInfo::cell_loop_post_list_index
std::vector< unsigned int > cell_loop_post_list_index
Definition: dof_info.h:662
internal::MatrixFreeFunctions::DoFInfo::vector_zero_range_list_index
std::vector< unsigned int > vector_zero_range_list_index
Definition: dof_info.h:637
internal::MatrixFreeFunctions::DoFInfo::IndexStorageVariants::full
@ full
internal::MatrixFreeFunctions::DoFInfo::ghost_dofs
std::vector< types::global_dof_index > ghost_dofs
Definition: dof_info.h:630
internal::MatrixFreeFunctions::DoFInfo::start_components
std::vector< unsigned int > start_components
Definition: dof_info.h:571
internal::MatrixFreeFunctions::DoFInfo::global_base_element_offset
unsigned int global_base_element_offset
Definition: dof_info.h:553
types::global_dof_index
unsigned int global_dof_index
Definition: types.h:76
internal::MatrixFreeFunctions::DoFInfo::index_storage_variants
std::vector< IndexStorageVariants > index_storage_variants[3]
Definition: dof_info.h:423
internal::MatrixFreeFunctions::DoFInfo::cell_loop_pre_list_index
std::vector< unsigned int > cell_loop_pre_list_index
Definition: dof_info.h:649
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
internal::MatrixFreeFunctions::DoFInfo::n_components
std::vector< unsigned int > n_components
Definition: dof_info.h:565
internal::MatrixFreeFunctions::DoFInfo::compute_face_index_compression
void compute_face_index_compression(const std::vector< FaceToCellTopology< length >> &faces)
internal::MatrixFreeFunctions::DoFInfo::store_plain_indices
bool store_plain_indices
Definition: dof_info.h:605
DynamicSparsityPattern
Definition: dynamic_sparsity_pattern.h:323
internal::MatrixFreeFunctions::ConstraintValues::insert_entries
unsigned short insert_entries(const std::vector< std::pair< types::global_dof_index, number2 >> &entries)
internal::MatrixFreeFunctions::DoFInfo::plain_dof_indices
std::vector< unsigned int > plain_dof_indices
Definition: dof_info.h:547
task_info.h
internal::MatrixFreeFunctions::DoFInfo::dofs_per_cell
std::vector< unsigned int > dofs_per_cell
Definition: dof_info.h:595
exceptions.h
mapping_info.h
internal::MatrixFreeFunctions::ConstraintValues::next_constraint
std::pair< std::vector< Number >, types::global_dof_index > next_constraint
Definition: dof_info.h:72
unsigned int
internal::MatrixFreeFunctions::DoFInfo::dof_access_face_exterior
@ dof_access_face_exterior
Definition: dof_info.h:392
AffineConstraints
Definition: affine_constraints.h:180
internal::MatrixFreeFunctions::DoFInfo::cell_loop_pre_list
std::vector< std::pair< unsigned int, unsigned int > > cell_loop_pre_list
Definition: dof_info.h:655
internal::MatrixFreeFunctions::DoFInfo::print
void print(const std::vector< Number > &constraint_pool_data, const std::vector< unsigned int > &constraint_pool_row_index, std::ostream &out) const
dof_handler.h
internal::MatrixFreeFunctions::DoFInfo::make_connectivity_graph
void make_connectivity_graph(const TaskInfo &task_info, const std::vector< unsigned int > &renumbering, DynamicSparsityPattern &connectivity) const
internal::MatrixFreeFunctions::DoFInfo::IndexStorageVariants
IndexStorageVariants
Definition: dof_info.h:296
affine_constraints.h
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
internal::MatrixFreeFunctions::DoFInfo::vectorization_length
unsigned int vectorization_length
Definition: dof_info.h:414
internal::MatrixFreeFunctions::DoFInfo::IndexStorageVariants::interleaved
@ interleaved
config.h
internal::MatrixFreeFunctions::DoFInfo::dimension
unsigned int dimension
Definition: dof_info.h:404
internal::MatrixFreeFunctions::DoFInfo::compute_vector_zero_access_pattern
void compute_vector_zero_access_pattern(const TaskInfo &task_info, const std::vector< FaceToCellTopology< length >> &faces)
internal::MatrixFreeFunctions::DoFInfo::clear
void clear()
internal::MatrixFreeFunctions::DoFInfo::fe_index_from_degree
unsigned int fe_index_from_degree(const unsigned int first_selected_component, const unsigned int fe_degree) const
internal
Definition: aligned_vector.h:369
internal::MatrixFreeFunctions::DoFInfo::print_memory_consumption
void print_memory_consumption(StreamType &out, const TaskInfo &size_info) const
internal::MatrixFreeFunctions::DoFInfo::assign_ghosts
void assign_ghosts(const std::vector< unsigned int > &boundary_cells)
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
internal::MatrixFreeFunctions::DoFInfo::dofs_per_face
std::vector< unsigned int > dofs_per_face
Definition: dof_info.h:600
internal::MatrixFreeFunctions::FPArrayComparator
Definition: mapping_info.h:566
internal::MatrixFreeFunctions::DoFInfo::cell_loop_post_list
std::vector< std::pair< unsigned int, unsigned int > > cell_loop_post_list
Definition: dof_info.h:668
internal::MatrixFreeFunctions::DoFInfo::DoFAccessIndex
DoFAccessIndex
Definition: dof_info.h:383
internal::MatrixFreeFunctions::DoFInfo::row_starts_plain_indices
std::vector< unsigned int > row_starts_plain_indices
Definition: dof_info.h:537
internal::MatrixFreeFunctions::DoFInfo::max_fe_index
unsigned int max_fe_index
Definition: dof_info.h:616
internal::MatrixFreeFunctions::DoFInfo::fe_index_conversion
std::vector< std::vector< unsigned int > > fe_index_conversion
Definition: dof_info.h:623
internal::MatrixFreeFunctions::DoFInfo::compute_cell_index_compression
void compute_cell_index_compression(const std::vector< unsigned char > &irregular_cells)
internal::MatrixFreeFunctions::DoFInfo::read_dof_indices
void read_dof_indices(const std::vector< types::global_dof_index > &local_indices, const std::vector< unsigned int > &lexicographic_inv, const AffineConstraints< number > &constraints, const unsigned int cell_number, ConstraintValues< double > &constraint_values, bool &cell_at_boundary)
internal::MatrixFreeFunctions::DoFInfo::vector_zero_range_list
std::vector< std::pair< unsigned int, unsigned int > > vector_zero_range_list
Definition: dof_info.h:642
internal::MatrixFreeFunctions::DoFInfo::n_base_elements
unsigned int n_base_elements
Definition: dof_info.h:559
internal::MatrixFreeFunctions::DoFInfo::row_starts
std::vector< std::pair< unsigned int, unsigned int > > row_starts
Definition: dof_info.h:431
internal::MatrixFreeFunctions::ConstraintValues
Definition: dof_info.h:51
internal::MatrixFreeFunctions::DoFInfo::reorder_cells
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)
internal::MatrixFreeFunctions::DoFInfo::chunk_size_zero_vector
static const unsigned int chunk_size_zero_vector
Definition: dof_info.h:113