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\}}\)
task_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_task_info_h
18 #define dealii_matrix_free_task_info_h
19 
20 
21 #include <deal.II/base/config.h>
22 
24 #include <deal.II/base/index_set.h>
26 #include <deal.II/base/tensor.h>
28 #include <deal.II/base/utilities.h>
30 
32 
33 
35 
36 
37 
38 namespace internal
39 {
47  {
48  public:
49  virtual ~MFWorkerInterface() = default;
50 
52  virtual void
54 
56  virtual void
58 
60  virtual void
62 
64  virtual void
66 
69  virtual void
70  zero_dst_vector_range(const unsigned int range_index) = 0;
71 
72  virtual void
73  cell_loop_pre_range(const unsigned int range_index) = 0;
74 
75  virtual void
76  cell_loop_post_range(const unsigned int range_index) = 0;
77 
80  virtual void
81  cell(const std::pair<unsigned int, unsigned int> &cell_range) = 0;
82 
85  virtual void
86  face(const std::pair<unsigned int, unsigned int> &face_range) = 0;
87 
90  virtual void
91  boundary(const std::pair<unsigned int, unsigned int> &face_range) = 0;
92  };
93 
94 
95 
96  namespace MatrixFreeFunctions
97  {
105  struct TaskInfo
106  {
107  // enum for choice of how to build the task graph. Odd add versions with
108  // preblocking and even versions with postblocking. partition_partition
109  // and partition_color are deprecated but kept for backward
110  // compatibility.
112  {
117  };
118 
122  TaskInfo();
123 
128  void
129  clear();
130 
134  void
135  loop(MFWorkerInterface &worker) const;
136 
141  void
142  make_boundary_cells_divisible(std::vector<unsigned int> &boundary_cells);
143 
186  void
188  const std::vector<unsigned int> &cells_with_comm,
189  const unsigned int dofs_per_cell,
190  const bool categories_are_hp,
191  const std::vector<unsigned int> &cell_vectorization_categories,
192  const bool cell_vectorization_categories_strict,
193  const std::vector<unsigned int> &parent_relation,
194  std::vector<unsigned int> & renumbering,
195  std::vector<unsigned char> & incompletely_filled_vectorization);
196 
214  void
216  const std::vector<unsigned int> &boundary_cells,
217  std::vector<unsigned int> & renumbering,
218  std::vector<unsigned char> & incompletely_filled_vectorization);
219 
226  void
227  guess_block_size(const unsigned int dofs_per_cell);
228 
255  void
257  DynamicSparsityPattern & connectivity,
258  std::vector<unsigned int> & renumbering,
259  std::vector<unsigned char> &irregular_cells,
260  const bool hp_bool);
261 
294  void
296  const std::vector<unsigned int> &cell_active_fe_index,
297  DynamicSparsityPattern & connectivity,
298  std::vector<unsigned int> & renumbering,
299  std::vector<unsigned char> & irregular_cells,
300  const bool hp_bool);
301 
325  void
326  make_thread_graph(const std::vector<unsigned int> &cell_active_fe_index,
327  DynamicSparsityPattern & connectivity,
328  std::vector<unsigned int> & renumbering,
329  std::vector<unsigned char> & irregular_cells,
330  const bool hp_bool);
331 
336  void
338  const std::vector<unsigned char> &irregular_cells,
339  const DynamicSparsityPattern & connectivity_cells,
340  DynamicSparsityPattern & connectivity_blocks) const;
341 
346  void
348  const DynamicSparsityPattern & connectivity,
349  const unsigned int partition,
350  const std::vector<unsigned int> &cell_partition,
351  const std::vector<unsigned int> &partition_list,
352  const std::vector<unsigned int> &partition_size,
353  std::vector<unsigned int> & partition_color_list);
354 
359  void
361  const DynamicSparsityPattern & connectivity,
362  const std::vector<unsigned int> &cell_active_fe_index,
363  const unsigned int partition,
364  const unsigned int cluster_size,
365  const bool hp_bool,
366  const std::vector<unsigned int> &cell_partition,
367  const std::vector<unsigned int> &partition_list,
368  const std::vector<unsigned int> &partition_size,
369  std::vector<unsigned int> & partition_partition_list,
370  std::vector<unsigned char> & irregular_cells);
371 
392  void
393  make_partitioning(const DynamicSparsityPattern &connectivity,
394  const unsigned int cluster_size,
395  std::vector<unsigned int> & cell_partition,
396  std::vector<unsigned int> & partition_list,
397  std::vector<unsigned int> & partition_size,
398  unsigned int & partition) const;
399 
404  void
405  update_task_info(const unsigned int partition);
406 
410  void
412 
416  std::size_t
417  memory_consumption() const;
418 
423  template <typename StreamType>
424  void
425  print_memory_statistics(StreamType &out, std::size_t data_length) const;
426 
431  unsigned int n_active_cells;
432 
437  unsigned int n_ghost_cells;
438 
442  unsigned int vectorization_length;
443 
447  unsigned int block_size;
448 
452  unsigned int n_blocks;
453 
458 
465  std::vector<unsigned int> partition_row_index;
466 
473  std::vector<unsigned int> cell_partition_data;
474 
482  std::vector<unsigned int> face_partition_data;
483 
491  std::vector<unsigned int> boundary_partition_data;
492 
501  std::vector<unsigned int> ghost_face_partition_data;
502 
512  std::vector<unsigned int> refinement_edge_face_partition_data;
513 
518  std::vector<unsigned int> partition_evens;
519 
524  std::vector<unsigned int> partition_odds;
525 
530  std::vector<unsigned int> partition_n_blocked_workers;
531 
536  std::vector<unsigned int> partition_n_workers;
537 
542  unsigned int evens;
543 
548  unsigned int odds;
549 
554  unsigned int n_blocked_workers;
555 
559  unsigned int n_workers;
560 
565  std::vector<unsigned char> task_at_mpi_boundary;
566 
571 
575  unsigned int my_pid;
576 
580  unsigned int n_procs;
581  };
582 
587 
588  } // end of namespace MatrixFreeFunctions
589 } // end of namespace internal
590 
592 
593 #endif
internal::MFWorkerInterface::face
virtual void face(const std::pair< unsigned int, unsigned int > &face_range)=0
internal::MatrixFreeFunctions::TaskInfo::make_thread_graph_partition_color
void make_thread_graph_partition_color(DynamicSparsityPattern &connectivity, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &irregular_cells, const bool hp_bool)
Definition: task_info.cc:1161
dynamic_sparsity_pattern.h
internal::MFWorkerInterface::vector_compress_start
virtual void vector_compress_start()=0
Starts the communication for the vector compress operation.
internal::MatrixFreeFunctions::TaskInfo::TaskInfo
TaskInfo()
Definition: task_info.cc:637
internal::MatrixFreeFunctions::TaskInfo::loop
void loop(MFWorkerInterface &worker) const
Definition: task_info.cc:339
internal::MatrixFreeFunctions::TaskInfo::make_partitioning
void make_partitioning(const DynamicSparsityPattern &connectivity, const unsigned int cluster_size, std::vector< unsigned int > &cell_partition, std::vector< unsigned int > &partition_list, std::vector< unsigned int > &partition_size, unsigned int &partition) const
Definition: task_info.cc:2012
internal::MatrixFreeFunctions::TaskInfo::boundary_partition_data
std::vector< unsigned int > boundary_partition_data
Definition: task_info.h:491
vectorization.h
internal::MFWorkerInterface::~MFWorkerInterface
virtual ~MFWorkerInterface()=default
internal::MatrixFreeFunctions::TaskInfo::partition_n_blocked_workers
std::vector< unsigned int > partition_n_blocked_workers
Definition: task_info.h:530
internal::MatrixFreeFunctions::TaskInfo::cell_partition_data
std::vector< unsigned int > cell_partition_data
Definition: task_info.h:473
internal::MatrixFreeFunctions::TaskInfo::communicator
MPI_Comm communicator
Definition: task_info.h:570
memory_consumption.h
utilities.h
internal::MatrixFreeFunctions::TaskInfo::none
@ none
Definition: task_info.h:113
internal::MatrixFreeFunctions::TaskInfo::n_blocked_workers
unsigned int n_blocked_workers
Definition: task_info.h:554
internal::MatrixFreeFunctions::TaskInfo::n_procs
unsigned int n_procs
Definition: task_info.h:580
internal::MatrixFreeFunctions::TaskInfo::make_partitioning_within_partitions_post_blocked
void make_partitioning_within_partitions_post_blocked(const DynamicSparsityPattern &connectivity, const std::vector< unsigned int > &cell_active_fe_index, const unsigned int partition, const unsigned int cluster_size, const bool hp_bool, const std::vector< unsigned int > &cell_partition, const std::vector< unsigned int > &partition_list, const std::vector< unsigned int > &partition_size, std::vector< unsigned int > &partition_partition_list, std::vector< unsigned char > &irregular_cells)
Definition: task_info.cc:1615
MPI_Comm
internal::MatrixFreeFunctions::TaskInfo::my_pid
unsigned int my_pid
Definition: task_info.h:575
internal::MFWorkerInterface
Definition: task_info.h:46
SparsityTools::partition
void partition(const SparsityPattern &sparsity_pattern, const unsigned int n_partitions, std::vector< unsigned int > &partition_indices, const Partitioner partitioner=Partitioner::metis)
Definition: sparsity_tools.cc:400
internal::MatrixFreeFunctions::TaskInfo
Definition: task_info.h:105
internal::MatrixFreeFunctions::TaskInfo::partition_odds
std::vector< unsigned int > partition_odds
Definition: task_info.h:524
thread_management.h
internal::MatrixFreeFunctions::TaskInfo::guess_block_size
void guess_block_size(const unsigned int dofs_per_cell)
Definition: task_info.cc:1133
internal::MatrixFreeFunctions::TaskInfo::n_ghost_cells
unsigned int n_ghost_cells
Definition: task_info.h:437
internal::MatrixFreeFunctions::TaskInfo::initial_setup_blocks_tasks
void initial_setup_blocks_tasks(const std::vector< unsigned int > &boundary_cells, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &incompletely_filled_vectorization)
Definition: task_info.cc:1067
tensor.h
internal::MatrixFreeFunctions::TaskInfo::evens
unsigned int evens
Definition: task_info.h:542
internal::MatrixFreeFunctions::TaskInfo::partition_partition
@ partition_partition
Definition: task_info.h:114
internal::MatrixFreeFunctions::TaskInfo::update_task_info
void update_task_info(const unsigned int partition)
Definition: task_info.cc:2241
internal::MatrixFreeFunctions::TaskInfo::task_at_mpi_boundary
std::vector< unsigned char > task_at_mpi_boundary
Definition: task_info.h:565
index_set.h
internal::MatrixFreeFunctions::TaskInfo::make_thread_graph
void make_thread_graph(const std::vector< unsigned int > &cell_active_fe_index, DynamicSparsityPattern &connectivity, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &irregular_cells, const bool hp_bool)
Definition: task_info.cc:1308
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
internal::MatrixFreeFunctions::TaskInfo::n_active_cells
unsigned int n_active_cells
Definition: task_info.h:431
internal::MatrixFreeFunctions::TaskInfo::make_thread_graph_partition_partition
void make_thread_graph_partition_partition(const std::vector< unsigned int > &cell_active_fe_index, DynamicSparsityPattern &connectivity, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &irregular_cells, const bool hp_bool)
Definition: task_info.cc:1503
internal::MatrixFreeFunctions::TaskInfo::memory_consumption
std::size_t memory_consumption() const
Definition: task_info.cc:690
internal::MatrixFreeFunctions::TaskInfo::print_memory_statistics
void print_memory_statistics(StreamType &out, std::size_t data_length) const
Definition: task_info.cc:675
DynamicSparsityPattern
Definition: dynamic_sparsity_pattern.h:323
internal::MatrixFreeFunctions::TaskInfo::n_blocks
unsigned int n_blocks
Definition: task_info.h:452
DEAL_II_DEPRECATED
#define DEAL_II_DEPRECATED
Definition: config.h:98
internal::MFWorkerInterface::boundary
virtual void boundary(const std::pair< unsigned int, unsigned int > &face_range)=0
internal::MatrixFreeFunctions::TaskInfo::ghost_face_partition_data
std::vector< unsigned int > ghost_face_partition_data
Definition: task_info.h:501
internal::MatrixFreeFunctions::TaskInfo::n_workers
unsigned int n_workers
Definition: task_info.h:559
internal::MFWorkerInterface::vector_update_ghosts_finish
virtual void vector_update_ghosts_finish()=0
Finishes the communication for the update ghost values operation.
exceptions.h
internal::MFWorkerInterface::cell_loop_pre_range
virtual void cell_loop_pre_range(const unsigned int range_index)=0
internal::MatrixFreeFunctions::TaskInfo::refinement_edge_face_partition_data
std::vector< unsigned int > refinement_edge_face_partition_data
Definition: task_info.h:512
internal::MFWorkerInterface::zero_dst_vector_range
virtual void zero_dst_vector_range(const unsigned int range_index)=0
internal::MFWorkerInterface::cell_loop_post_range
virtual void cell_loop_post_range(const unsigned int range_index)=0
internal::MFWorkerInterface::cell
virtual void cell(const std::pair< unsigned int, unsigned int > &cell_range)=0
internal::MatrixFreeFunctions::TaskInfo::face_partition_data
std::vector< unsigned int > face_partition_data
Definition: task_info.h:482
internal::MatrixFreeFunctions::TaskInfo::odds
unsigned int odds
Definition: task_info.h:548
internal::MatrixFreeFunctions::TaskInfo::partition_row_index
std::vector< unsigned int > partition_row_index
Definition: task_info.h:465
internal::MFWorkerInterface::vector_compress_finish
virtual void vector_compress_finish()=0
Finishes the communication for the vector compress operation.
internal::MatrixFreeFunctions::TaskInfo::TasksParallelScheme
TasksParallelScheme
Definition: task_info.h:111
config.h
internal::MatrixFreeFunctions::TaskInfo::vectorization_length
unsigned int vectorization_length
Definition: task_info.h:442
internal::MatrixFreeFunctions::TaskInfo::make_coloring_within_partitions_pre_blocked
void make_coloring_within_partitions_pre_blocked(const DynamicSparsityPattern &connectivity, const unsigned int partition, const std::vector< unsigned int > &cell_partition, const std::vector< unsigned int > &partition_list, const std::vector< unsigned int > &partition_size, std::vector< unsigned int > &partition_color_list)
Definition: task_info.cc:1936
internal::MatrixFreeFunctions::TaskInfo::partition_color
@ partition_color
Definition: task_info.h:115
internal
Definition: aligned_vector.h:369
internal::MatrixFreeFunctions::TaskInfo::color
@ color
Definition: task_info.h:116
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
internal::MatrixFreeFunctions::TaskInfo::block_size
unsigned int block_size
Definition: task_info.h:447
internal::MatrixFreeFunctions::TaskInfo::clear
void clear()
Definition: task_info.cc:645
internal::MatrixFreeFunctions::TaskInfo::create_blocks_serial
void create_blocks_serial(const std::vector< unsigned int > &cells_with_comm, const unsigned int dofs_per_cell, const bool categories_are_hp, const std::vector< unsigned int > &cell_vectorization_categories, const bool cell_vectorization_categories_strict, const std::vector< unsigned int > &parent_relation, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &incompletely_filled_vectorization)
Definition: task_info.cc:779
internal::MatrixFreeFunctions::TaskInfo::make_boundary_cells_divisible
void make_boundary_cells_divisible(std::vector< unsigned int > &boundary_cells)
Definition: task_info.cc:707
internal::MatrixFreeFunctions::TaskInfo::partition_evens
std::vector< unsigned int > partition_evens
Definition: task_info.h:518
internal::MatrixFreeFunctions::TaskInfo::scheme
TasksParallelScheme scheme
Definition: task_info.h:457
internal::MatrixFreeFunctions::TaskInfo::partition_n_workers
std::vector< unsigned int > partition_n_workers
Definition: task_info.h:536
internal::MFWorkerInterface::vector_update_ghosts_start
virtual void vector_update_ghosts_start()=0
Starts the communication for the update ghost values operation.
internal::MatrixFreeFunctions::TaskInfo::create_flow_graph
void create_flow_graph()
internal::MatrixFreeFunctions::TaskInfo::make_connectivity_cells_to_blocks
void make_connectivity_cells_to_blocks(const std::vector< unsigned char > &irregular_cells, const DynamicSparsityPattern &connectivity_cells, DynamicSparsityPattern &connectivity_blocks) const
Definition: task_info.cc:1571