Reference documentation for deal.II version 9.4.1
\(\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\}}\)
Loading...
Searching...
No Matches
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
26
31
32#include <array>
33#include <memory>
34
35
37
38#ifndef DOXYGEN
39
40// forward declarations
41
42namespace internal
43{
44 namespace MatrixFreeFunctions
45 {
46 template <int dim>
47 class HangingNodes;
48
49 template <typename, typename>
50 struct FPArrayComparator;
51 } // namespace MatrixFreeFunctions
52} // namespace internal
53
54template <typename>
56
58
59template <typename>
60class TriaIterator;
61
62template <int, int, bool>
63class DoFCellAccessor;
64
65#endif
66
67namespace internal
68{
69 namespace MatrixFreeFunctions
70 {
75 using compressed_constraint_kind = std::uint8_t;
76
81 template <typename Number>
83 {
85
93 template <typename number2>
94 unsigned short
96 const std::vector<std::pair<types::global_dof_index, number2>>
97 &entries);
98
99 std::vector<std::pair<types::global_dof_index, double>>
101 std::vector<types::global_dof_index> constraint_indices;
102
103 std::pair<std::vector<Number>, types::global_dof_index> next_constraint;
104 std::map<std::vector<Number>,
108 };
109
128 struct DoFInfo
129 {
142 static const unsigned int chunk_size_zero_vector = 64;
143
147 DoFInfo();
148
152 DoFInfo(const DoFInfo &) = default;
153
157 DoFInfo(DoFInfo &&) noexcept = default;
158
162 ~DoFInfo() = default;
163
167 DoFInfo &
168 operator=(const DoFInfo &) = default;
169
173 DoFInfo &
174 operator=(DoFInfo &&) noexcept = default;
175
179 void
180 clear();
181
187 unsigned int
188 fe_index_from_degree(const unsigned int first_selected_component,
189 const unsigned int fe_degree) const;
190
210 void
211 get_dof_indices_on_cell_batch(std::vector<unsigned int> &local_indices,
212 const unsigned int cell_batch,
213 const bool with_constraints = true) const;
214
224 template <typename number>
225 void
227 const std::vector<types::global_dof_index> &local_indices_resolved,
228 const std::vector<types::global_dof_index> &local_indices,
229 const bool cell_has_hanging_nodes,
230 const ::AffineConstraints<number> & constraints,
231 const unsigned int cell_number,
232 ConstraintValues<double> & constraint_values,
233 bool & cell_at_boundary);
234
239 template <int dim>
240 bool
242 const HangingNodes<dim> & hanging_nodes,
243 const std::vector<std::vector<unsigned int>> &lexicographic_mapping,
244 const unsigned int cell_number,
245 const TriaIterator<DoFCellAccessor<dim, dim, false>> &cell,
246 std::vector<types::global_dof_index> & dof_indices);
247
255 void
256 assign_ghosts(const std::vector<unsigned int> &boundary_cells,
257 const MPI_Comm & communicator_sm,
258 const bool use_vector_data_exchanger_full);
259
266 void
267 reorder_cells(const TaskInfo & task_info,
268 const std::vector<unsigned int> & renumbering,
269 const std::vector<unsigned int> & constraint_pool_row_index,
270 const std::vector<unsigned char> &irregular_cells);
271
276 void
278 const std::vector<unsigned char> &irregular_cells);
279
284 template <int length>
285 void
287 const std::vector<FaceToCellTopology<length>> &faces);
288
294 void
295 make_connectivity_graph(const TaskInfo & task_info,
296 const std::vector<unsigned int> &renumbering,
297 DynamicSparsityPattern &connectivity) const;
298
304 void
306 const Table<2, ShapeInfo<double>> & shape_info,
307 const unsigned int n_owned_cells,
308 const unsigned int n_lanes,
309 const std::vector<FaceToCellTopology<1>> &inner_faces,
310 const std::vector<FaceToCellTopology<1>> &ghosted_faces,
311 const bool fill_cell_centric,
312 const MPI_Comm & communicator_sm,
313 const bool use_vector_data_exchanger_full);
314
320 void
322 std::array<std::vector<std::pair<unsigned int, unsigned int>>, 3>
323 &cell_indices_contiguous_sm);
324
335 void
337 std::vector<types::global_dof_index> &renumbering);
338
348 template <int length>
349 void
351 const TaskInfo & task_info,
352 const std::vector<FaceToCellTopology<length>> &faces);
353
357 std::size_t
358 memory_consumption() const;
359
364 template <typename StreamType>
365 void
366 print_memory_consumption(StreamType & out,
367 const TaskInfo &size_info) const;
368
373 template <typename Number>
374 void
375 print(const std::vector<Number> & constraint_pool_data,
376 const std::vector<unsigned int> &constraint_pool_row_index,
377 std::ostream & out) const;
378
390 enum class IndexStorageVariants : unsigned char
391 {
400 full,
471 };
472
477 enum DoFAccessIndex : unsigned char
478 {
491 };
492
498 unsigned int dimension;
499
509
517 std::array<std::vector<IndexStorageVariants>, 3> index_storage_variants;
518
525 std::vector<std::pair<unsigned int, unsigned int>> row_starts;
526
542 std::vector<unsigned int> dof_indices;
543
548 std::vector<std::vector<bool>> hanging_node_constraint_masks_comp;
549
554 std::vector<compressed_constraint_kind> hanging_node_constraint_masks;
555
565 std::vector<std::pair<unsigned short, unsigned short>>
567
571 std::vector<unsigned int> dof_indices_interleaved;
572
581 std::array<std::vector<unsigned int>, 3> dof_indices_contiguous;
582
591 std::array<std::vector<std::pair<unsigned int, unsigned int>>, 3>
593
602 std::array<std::vector<unsigned int>, 3> dof_indices_interleave_strides;
603
613 std::array<std::vector<unsigned char>, 3> n_vectorization_lanes_filled;
614
621 std::shared_ptr<const Utilities::MPI::Partitioner> vector_partitioner;
622
626 std::shared_ptr<
629
649 std::array<
650 std::shared_ptr<
652 5>
654
659 std::vector<unsigned int> constrained_dofs;
660
665 std::vector<unsigned int> row_starts_plain_indices;
666
675 std::vector<unsigned int> plain_dof_indices;
676
682
687 unsigned int n_base_elements;
688
693 std::vector<unsigned int> n_components;
694
699 std::vector<unsigned int> start_components;
700
705 std::vector<unsigned int> component_to_base_index;
706
718 std::vector<std::vector<unsigned int>> component_dof_indices_offset;
719
723 std::vector<unsigned int> dofs_per_cell;
724
728 std::vector<unsigned int> dofs_per_face;
729
734
738 std::vector<unsigned int> cell_active_fe_index;
739
744 unsigned int max_fe_index;
745
751 std::vector<std::vector<unsigned int>> fe_index_conversion;
752
758 std::vector<types::global_dof_index> ghost_dofs;
759
765 std::vector<unsigned int> vector_zero_range_list_index;
766
770 std::vector<std::pair<unsigned int, unsigned int>> vector_zero_range_list;
771
777 std::vector<unsigned int> cell_loop_pre_list_index;
778
783 std::vector<std::pair<unsigned int, unsigned int>> cell_loop_pre_list;
784
790 std::vector<unsigned int> cell_loop_post_list_index;
791
796 std::vector<std::pair<unsigned int, unsigned int>> cell_loop_post_list;
797 };
798
799
800 /*-------------------------- Inline functions ---------------------------*/
801
802#ifndef DOXYGEN
803
804
805 inline unsigned int
806 DoFInfo::fe_index_from_degree(const unsigned int first_selected_component,
807 const unsigned int fe_degree) const
808 {
809 const unsigned int n_indices = fe_index_conversion.size();
810 if (n_indices <= 1)
811 return 0;
812 for (unsigned int i = 0; i < n_indices; ++i)
813 if (fe_index_conversion[i][first_selected_component] == fe_degree)
814 return i;
816 }
817
818#endif // ifndef DOXYGEN
819
820 } // end of namespace MatrixFreeFunctions
821} // end of namespace internal
822
824
825#endif
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
std::uint8_t compressed_constraint_kind
Definition: dof_info.h:75
static const unsigned int invalid_unsigned_int
Definition: types.h:201
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:101
std::map< std::vector< Number >, types::global_dof_index, FPArrayComparator< Number, VectorizedArray< Number > > > constraints
Definition: dof_info.h:107
std::vector< std::pair< types::global_dof_index, double > > constraint_entries
Definition: dof_info.h:100
std::pair< std::vector< Number >, types::global_dof_index > next_constraint
Definition: dof_info.h:103
unsigned short insert_entries(const std::vector< std::pair< types::global_dof_index, number2 > > &entries)
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:296
std::vector< std::pair< unsigned short, unsigned short > > constraint_indicator
Definition: dof_info.h:566
std::vector< unsigned int > cell_loop_pre_list_index
Definition: dof_info.h:777
std::vector< std::pair< unsigned int, unsigned int > > row_starts
Definition: dof_info.h:525
std::vector< std::vector< unsigned int > > component_dof_indices_offset
Definition: dof_info.h:718
std::vector< std::vector< bool > > hanging_node_constraint_masks_comp
Definition: dof_info.h:548
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:80
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:528
void compute_face_index_compression(const std::vector< FaceToCellTopology< length > > &faces)
std::vector< unsigned int > cell_loop_post_list_index
Definition: dof_info.h:790
std::vector< unsigned int > dofs_per_cell
Definition: dof_info.h:723
void compute_dof_renumbering(std::vector< types::global_dof_index > &renumbering)
Definition: dof_info.cc:1436
static const unsigned int chunk_size_zero_vector
Definition: dof_info.h:142
std::shared_ptr< const internal::MatrixFreeFunctions::VectorDataExchange::Base > vector_exchanger
Definition: dof_info.h:628
std::vector< std::vector< unsigned int > > fe_index_conversion
Definition: dof_info.h:751
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:860
std::vector< unsigned int > dof_indices
Definition: dof_info.h:542
std::vector< std::pair< unsigned int, unsigned int > > vector_zero_range_list
Definition: dof_info.h:770
std::shared_ptr< const Utilities::MPI::Partitioner > vector_partitioner
Definition: dof_info.h:621
std::vector< compressed_constraint_kind > hanging_node_constraint_masks
Definition: dof_info.h:554
std::array< std::vector< unsigned int >, 3 > dof_indices_interleave_strides
Definition: dof_info.h:602
std::array< std::vector< std::pair< unsigned int, unsigned int > >, 3 > dof_indices_contiguous_sm
Definition: dof_info.h:592
std::vector< unsigned int > row_starts_plain_indices
Definition: dof_info.h:665
std::vector< unsigned int > dofs_per_face
Definition: dof_info.h:728
std::size_t memory_consumption() const
Definition: dof_info.cc:1494
std::vector< unsigned int > component_to_base_index
Definition: dof_info.h:705
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:693
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:153
std::vector< std::pair< unsigned int, unsigned int > > cell_loop_pre_list
Definition: dof_info.h:783
std::vector< types::global_dof_index > ghost_dofs
Definition: dof_info.h:758
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:581
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)
DoFInfo(const DoFInfo &)=default
std::vector< unsigned int > cell_active_fe_index
Definition: dof_info.h:738
std::array< std::shared_ptr< const internal::MatrixFreeFunctions::VectorDataExchange::Base >, 5 > vector_exchanger_face_variants
Definition: dof_info.h:653
std::vector< unsigned int > plain_dof_indices
Definition: dof_info.h:675
std::array< std::vector< unsigned char >, 3 > n_vectorization_lanes_filled
Definition: dof_info.h:613
std::vector< unsigned int > start_components
Definition: dof_info.h:699
std::vector< unsigned int > dof_indices_interleaved
Definition: dof_info.h:571
std::vector< unsigned int > constrained_dofs
Definition: dof_info.h:659
std::array< std::vector< IndexStorageVariants >, 3 > index_storage_variants
Definition: dof_info.h:517
void make_connectivity_graph(const TaskInfo &task_info, const std::vector< unsigned int > &renumbering, DynamicSparsityPattern &connectivity) const
Definition: dof_info.cc:1359
std::vector< unsigned int > vector_zero_range_list_index
Definition: dof_info.h:765
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:1210
std::vector< std::pair< unsigned int, unsigned int > > cell_loop_post_list
Definition: dof_info.h:796
DoFInfo(DoFInfo &&) noexcept=default