Reference documentation for deal.II version 9.3.3
\(\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
37
38#include <array>
39#include <memory>
40
41
43
44namespace internal
45{
46 namespace MatrixFreeFunctions
47 {
52 template <typename Number>
54 {
56
64 template <typename number2>
65 unsigned short
67 const std::vector<std::pair<types::global_dof_index, number2>>
68 &entries);
69
70 std::vector<std::pair<types::global_dof_index, double>>
72 std::vector<types::global_dof_index> constraint_indices;
73
74 std::pair<std::vector<Number>, types::global_dof_index> next_constraint;
75 std::map<std::vector<Number>,
79 };
80
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 DoFInfo(DoFInfo &&) noexcept = default;
129
133 ~DoFInfo() = default;
134
138 DoFInfo &
139 operator=(const DoFInfo &) = default;
140
144 DoFInfo &
145 operator=(DoFInfo &&) noexcept = default;
146
150 void
151 clear();
152
158 unsigned int
159 fe_index_from_degree(const unsigned int first_selected_component,
160 const unsigned int fe_degree) const;
161
181 void
182 get_dof_indices_on_cell_batch(std::vector<unsigned int> &locall_indices,
183 const unsigned int cell,
184 const bool with_constraints = true) const;
185
194 template <typename number>
195 void
197 const std::vector<types::global_dof_index> &local_indices,
198 const std::vector<unsigned int> & lexicographic_inv,
199 const ::AffineConstraints<number> & constraints,
200 const unsigned int cell_number,
201 ConstraintValues<double> & constraint_values,
202 bool & cell_at_boundary);
203
211 void
212 assign_ghosts(const std::vector<unsigned int> &boundary_cells,
213 const MPI_Comm & communicator_sm,
214 const bool use_vector_data_exchanger_full);
215
222 void
223 reorder_cells(const TaskInfo & task_info,
224 const std::vector<unsigned int> & renumbering,
225 const std::vector<unsigned int> & constraint_pool_row_index,
226 const std::vector<unsigned char> &irregular_cells);
227
232 void
234 const std::vector<unsigned char> &irregular_cells);
235
240 template <int length>
241 void
243 const std::vector<FaceToCellTopology<length>> &faces);
244
250 void
251 make_connectivity_graph(const TaskInfo & task_info,
252 const std::vector<unsigned int> &renumbering,
253 DynamicSparsityPattern &connectivity) const;
254
260 void
262 const Table<2, ShapeInfo<double>> & shape_info,
263 const unsigned int n_owned_cells,
264 const unsigned int n_lanes,
265 const std::vector<FaceToCellTopology<1>> &inner_faces,
266 const std::vector<FaceToCellTopology<1>> &ghosted_faces,
267 const bool fill_cell_centric,
268 const MPI_Comm & communicator_sm,
269 const bool use_vector_data_exchanger_full);
270
276 void
278 std::array<std::vector<std::pair<unsigned int, unsigned int>>, 3>
279 &cell_indices_contiguous_sm);
280
291 void
293 std::vector<types::global_dof_index> &renumbering);
294
304 template <int length>
305 void
307 const TaskInfo & task_info,
308 const std::vector<FaceToCellTopology<length>> &faces);
309
313 std::size_t
314 memory_consumption() const;
315
320 template <typename StreamType>
321 void
322 print_memory_consumption(StreamType & out,
323 const TaskInfo &size_info) const;
324
329 template <typename Number>
330 void
331 print(const std::vector<Number> & constraint_pool_data,
332 const std::vector<unsigned int> &constraint_pool_row_index,
333 std::ostream & out) const;
334
346 enum class IndexStorageVariants : unsigned char
347 {
356 full,
368 interleaved,
377 contiguous,
401 interleaved_contiguous,
413 interleaved_contiguous_strided,
426 interleaved_contiguous_mixed_strides
427 };
428
433 enum DoFAccessIndex : unsigned char
434 {
447 };
448
454 unsigned int dimension;
455
465
473 std::array<std::vector<IndexStorageVariants>, 3> index_storage_variants;
474
481 std::vector<std::pair<unsigned int, unsigned int>> row_starts;
482
498 std::vector<unsigned int> dof_indices;
499
509 std::vector<std::pair<unsigned short, unsigned short>>
511
515 std::vector<unsigned int> dof_indices_interleaved;
516
525 std::array<std::vector<unsigned int>, 3> dof_indices_contiguous;
526
535 std::array<std::vector<std::pair<unsigned int, unsigned int>>, 3>
537
546 std::array<std::vector<unsigned int>, 3> dof_indices_interleave_strides;
547
557 std::array<std::vector<unsigned char>, 3> n_vectorization_lanes_filled;
558
565 std::shared_ptr<const Utilities::MPI::Partitioner> vector_partitioner;
566
570 std::shared_ptr<
573
593 std::array<
594 std::shared_ptr<
596 5>
598
603 std::vector<unsigned int> constrained_dofs;
604
609 std::vector<unsigned int> row_starts_plain_indices;
610
619 std::vector<unsigned int> plain_dof_indices;
620
626
631 unsigned int n_base_elements;
632
637 std::vector<unsigned int> n_components;
638
643 std::vector<unsigned int> start_components;
644
649 std::vector<unsigned int> component_to_base_index;
650
662 std::vector<std::vector<unsigned int>> component_dof_indices_offset;
663
667 std::vector<unsigned int> dofs_per_cell;
668
672 std::vector<unsigned int> dofs_per_face;
673
678
682 std::vector<unsigned int> cell_active_fe_index;
683
688 unsigned int max_fe_index;
689
695 std::vector<std::vector<unsigned int>> fe_index_conversion;
696
702 std::vector<types::global_dof_index> ghost_dofs;
703
709 std::vector<unsigned int> vector_zero_range_list_index;
710
714 std::vector<std::pair<unsigned int, unsigned int>> vector_zero_range_list;
715
721 std::vector<unsigned int> cell_loop_pre_list_index;
722
727 std::vector<std::pair<unsigned int, unsigned int>> cell_loop_pre_list;
728
734 std::vector<unsigned int> cell_loop_post_list_index;
735
740 std::vector<std::pair<unsigned int, unsigned int>> cell_loop_post_list;
741 };
742
743
744 /*-------------------------- Inline functions ---------------------------*/
745
746#ifndef DOXYGEN
747
748
749 inline unsigned int
750 DoFInfo::fe_index_from_degree(const unsigned int first_selected_component,
751 const unsigned int fe_degree) const
752 {
753 const unsigned int n_indices = fe_index_conversion.size();
754 if (n_indices <= 1)
755 return 0;
756 for (unsigned int i = 0; i < n_indices; ++i)
757 if (fe_index_conversion[i][first_selected_component] == fe_degree)
758 return i;
760 }
761
762#endif // ifndef DOXYGEN
763
764 } // end of namespace MatrixFreeFunctions
765} // end of namespace internal
766
768
769#endif
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
static const unsigned int invalid_unsigned_int
Definition: types.h:196
STL namespace.
Definition: types.h:32
unsigned int global_dof_index
Definition: types.h:76
std::vector< types::global_dof_index > constraint_indices
Definition: dof_info.h:72
unsigned short insert_entries(const std::vector< std::pair< types::global_dof_index, number2 > > &entries)
std::vector< std::pair< types::global_dof_index, double > > constraint_entries
Definition: dof_info.h:71
std::map< std::vector< Number >, types::global_dof_index, FPArrayComparator< Number > > constraints
Definition: dof_info.h:78
std::pair< std::vector< Number >, types::global_dof_index > next_constraint
Definition: dof_info.h:74
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:271
std::vector< std::pair< unsigned short, unsigned short > > constraint_indicator
Definition: dof_info.h:510
std::vector< unsigned int > cell_loop_pre_list_index
Definition: dof_info.h:721
std::vector< std::pair< unsigned int, unsigned int > > row_starts
Definition: dof_info.h:481
std::vector< std::vector< unsigned int > > component_dof_indices_offset
Definition: dof_info.h:662
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:473
void compute_face_index_compression(const std::vector< FaceToCellTopology< length > > &faces)
std::vector< unsigned int > cell_loop_post_list_index
Definition: dof_info.h:734
std::vector< unsigned int > dofs_per_cell
Definition: dof_info.h:667
void compute_dof_renumbering(std::vector< types::global_dof_index > &renumbering)
Definition: dof_info.cc:1379
static const unsigned int chunk_size_zero_vector
Definition: dof_info.h:113
std::shared_ptr< const internal::MatrixFreeFunctions::VectorDataExchange::Base > vector_exchanger
Definition: dof_info.h:572
std::vector< std::vector< unsigned int > > fe_index_conversion
Definition: dof_info.h:695
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:803
std::vector< unsigned int > dof_indices
Definition: dof_info.h:498
std::vector< std::pair< unsigned int, unsigned int > > vector_zero_range_list
Definition: dof_info.h:714
std::shared_ptr< const Utilities::MPI::Partitioner > vector_partitioner
Definition: dof_info.h:565
std::array< std::vector< unsigned int >, 3 > dof_indices_interleave_strides
Definition: dof_info.h:546
std::array< std::vector< std::pair< unsigned int, unsigned int > >, 3 > dof_indices_contiguous_sm
Definition: dof_info.h:536
std::vector< unsigned int > row_starts_plain_indices
Definition: dof_info.h:609
std::vector< unsigned int > dofs_per_face
Definition: dof_info.h:672
std::size_t memory_consumption() const
Definition: dof_info.cc:1437
std::vector< unsigned int > component_to_base_index
Definition: dof_info.h:649
void compute_vector_zero_access_pattern(const TaskInfo &task_info, const std::vector< FaceToCellTopology< length > > &faces)
std::vector< unsigned int > n_components
Definition: dof_info.h:637
void print_memory_consumption(StreamType &out, const TaskInfo &size_info) const
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)
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:141
std::vector< std::pair< unsigned int, unsigned int > > cell_loop_pre_list
Definition: dof_info.h:727
std::vector< types::global_dof_index > ghost_dofs
Definition: dof_info.h:702
std::array< std::vector< unsigned int >, 3 > dof_indices_contiguous
Definition: dof_info.h:525
DoFInfo(const DoFInfo &)=default
void get_dof_indices_on_cell_batch(std::vector< unsigned int > &locall_indices, const unsigned int cell, const bool with_constraints=true) const
Definition: dof_info.cc:73
std::vector< unsigned int > cell_active_fe_index
Definition: dof_info.h:682
std::array< std::shared_ptr< const internal::MatrixFreeFunctions::VectorDataExchange::Base >, 5 > vector_exchanger_face_variants
Definition: dof_info.h:597
std::vector< unsigned int > plain_dof_indices
Definition: dof_info.h:619
std::array< std::vector< unsigned char >, 3 > n_vectorization_lanes_filled
Definition: dof_info.h:557
std::vector< unsigned int > start_components
Definition: dof_info.h:643
std::vector< unsigned int > dof_indices_interleaved
Definition: dof_info.h:515
std::vector< unsigned int > constrained_dofs
Definition: dof_info.h:603
std::array< std::vector< IndexStorageVariants >, 3 > index_storage_variants
Definition: dof_info.h:473
void make_connectivity_graph(const TaskInfo &task_info, const std::vector< unsigned int > &renumbering, DynamicSparsityPattern &connectivity) const
Definition: dof_info.cc:1302
std::vector< unsigned int > vector_zero_range_list_index
Definition: dof_info.h:709
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:1153
std::vector< std::pair< unsigned int, unsigned int > > cell_loop_post_list
Definition: dof_info.h:740
DoFInfo(DoFInfo &&) noexcept=default