Reference documentation for deal.II version 9.6.0
|
#include <deal.II/matrix_free/task_info.h>
Public Types | |
enum | TasksParallelScheme { none , partition_partition , partition_color , color } |
Public Member Functions | |
TaskInfo () | |
void | clear () |
void | loop (MFWorkerInterface &worker) const |
void | make_boundary_cells_divisible (std::vector< unsigned int > &boundary_cells) |
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) |
void | initial_setup_blocks_tasks (const std::vector< unsigned int > &boundary_cells, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &incompletely_filled_vectorization) |
void | guess_block_size (const unsigned int dofs_per_cell) |
void | make_thread_graph_partition_color (DynamicSparsityPattern &connectivity, std::vector< unsigned int > &renumbering, std::vector< unsigned char > &irregular_cells, const bool hp_bool) |
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) |
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) |
void | make_connectivity_cells_to_blocks (const std::vector< unsigned char > &irregular_cells, const DynamicSparsityPattern &connectivity_cells, DynamicSparsityPattern &connectivity_blocks) const |
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) |
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) |
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 |
void | update_task_info (const unsigned int partition) |
void | create_flow_graph () |
std::size_t | memory_consumption () const |
template<typename StreamType > | |
void | print_memory_statistics (StreamType &out, std::size_t data_length) const |
A struct that collects all information related to parallelization with threads: The work is subdivided into tasks that can be done independently.
Definition at line 106 of file task_info.h.
Enumerator | |
---|---|
none | |
partition_partition | |
partition_color | |
color |
Definition at line 112 of file task_info.h.
internal::MatrixFreeFunctions::TaskInfo::TaskInfo | ( | ) |
Constructor.
Definition at line 652 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::clear | ( | ) |
Clears all the data fields and resets them to zero.
Definition at line 660 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::loop | ( | MFWorkerInterface & | worker | ) | const |
Runs the matrix-free loop.
Definition at line 350 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::make_boundary_cells_divisible | ( | std::vector< unsigned int > & | boundary_cells | ) |
Make the number of cells which can only be treated in the communication overlap divisible by the vectorization length.
Definition at line 722 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::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 ) |
Sets up the blocks for running the cell loop based on the options controlled by the input arguments.
cells_with_comm | A list of cells that need to exchange data prior to performing computations. These will be given a certain id in the partitioning to make sure cell loops that overlap communication with communication have the ghost data ready. |
dofs_per_cell | Gives an expected value for the number of degrees of freedom on a cell, which is used to determine the block size for interleaving cell and face integrals. |
categories_are_hp | Defines whether cell_vectorization_categories is originating from a hp-adaptive computation with variable polynomial degree or a user-defined variant. |
cell_vectorization_categories | This set of categories defines the cells that should be grouped together inside the lanes of a vectorized array. This can be the polynomial degree in an hp-element or a user-provided grouping. |
cell_vectorization_categories_strict | Defines whether the categories defined by the previous variables should be separated strictly or whether it is allowed to insert lower categories into the next high one(s). |
parent_relation | This data field is used to specify which cells have the same parent cell. Cells with the same ancestor are grouped together into the same batch(es) with vectorization across cells. |
renumbering | When leaving this function, the vector contains a new numbering of the cells that aligns with the grouping stored in this class. |
incompletely_filled_vectorization | Given the vectorized layout of this class, some cell batches might have components in the vectorized array (SIMD lanes) that are not used and do not carray valid data. This array indicates the cell batches where this occurs according to the renumbering returned by this function. |
Definition at line 794 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::initial_setup_blocks_tasks | ( | const std::vector< unsigned int > & | boundary_cells, |
std::vector< unsigned int > & | renumbering, | ||
std::vector< unsigned char > & | incompletely_filled_vectorization ) |
First step in the block creation for the task-parallel blocking setup.
boundary_cells | A list of cells that need to exchange data prior to performing computations. These will be given a certain id in the partitioning. |
renumbering | When leaving this function, the vector contains a new numbering of the cells that aligns with the grouping stored in this class (before actually creating the tasks). |
incompletely_filled_vectorization | Given the vectorized layout of this class, some cell batches might have components in the vectorized array (SIMD lanes) that are not used and do not carray valid data. This array indicates the cell batches where this occurs according to the renumbering returned by this function. |
Definition at line 1150 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::guess_block_size | ( | const unsigned int | dofs_per_cell | ) |
This helper function determines a block size if the user decided not to force a block size through MatrixFree::AdditionalData. This is computed based on the number of hardware threads on the system and the number of cell batches that we should work on.
Definition at line 1216 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::make_thread_graph_partition_color | ( | DynamicSparsityPattern & | connectivity, |
std::vector< unsigned int > & | renumbering, | ||
std::vector< unsigned char > & | irregular_cells, | ||
const bool | hp_bool ) |
This method goes through all cells that have been filled into dof_indices
and finds out which cells can be worked on independently and which ones are neighboring and need to be done at different times when used in parallel.
The strategy is based on a two-level approach. The outer level is subdivided into partitions similar to the type of neighbors in Cuthill-McKee, and the inner level is subdivided via colors (for chunks within the same color, can work independently). One task is represented by a chunk of cells. The cell chunks are formed before subdivision into partitions and colors.
connectivity | (in/out) Determines whether cells i and j are conflicting, expressed by an entry in position (i,j). |
renumbering | (in/out) At output, the element j of this variable gives the original number of the cell that is reordered to place j by the ordering due to the thread graph. |
irregular_cells | (in/out) Informs the current function whether some SIMD lanes in VectorizedArray would not be filled for a given cell batch index. |
hp_bool | Defines whether we are in hp-mode or not |
Definition at line 1244 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::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 ) |
This function goes through all cells that have been filled into dof_indices
and finds out which cells can be worked on independently and which ones are neighboring and need to be done at different times when used in parallel.
The strategy is based on a two-level approach. The outer level is subdivided into partitions similar to the type of neighbors in Cuthill-McKee, and the inner level is again subdivided into Cuthill-McKee-like partitions (partitions whose level differs by more than 2 can be worked on independently). One task is represented by a chunk of cells. The cell chunks are formed after subdivision into the two levels of partitions.
cell_active_fe_index | The active FE index corresponding to the individual indices in the list of all cell indices, in order to be able to not place cells with different indices into the same cell batch with vectorization. |
connectivity | (in/out) Determines whether cells i and j are conflicting, expressed by an entry in position (i,j). |
renumbering | (in/out) At output, the element j of this variable gives the original number of the cell that is reordered to place j by the ordering due to the thread graph. |
irregular_cells | (in/out) Informs the current function whether some SIMD lanes in VectorizedArray would not be filled for a given cell batch index. |
hp_bool | Defines whether we are in hp-mode or not |
Definition at line 1586 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::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 ) |
Either calls make_thread_graph_partition_color() or make_thread_graph_partition_partition() accessible from the outside, depending on the setting in the data structure.
cell_active_fe_index | The active FE index corresponding to the individual indices in the list of all cell indices, in order to be able to not place cells with different indices into the same cell batch with vectorization. |
connectivity | (in/out) Determines whether cells i and j are conflicting, expressed by an entry in position (i,j). |
renumbering | (in/out) At output, the element j of this variable gives the original number of the cell that is reordered to place j by the ordering due to the thread graph. |
irregular_cells | (in/out) Informs the current function whether some SIMD lanes in VectorizedArray would not be filled for a given cell batch index. |
hp_bool | Defines whether we are in hp-mode or not |
Definition at line 1391 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::make_connectivity_cells_to_blocks | ( | const std::vector< unsigned char > & | irregular_cells, |
const DynamicSparsityPattern & | connectivity_cells, | ||
DynamicSparsityPattern & | connectivity_blocks ) const |
This function computes the connectivity between blocks of cells from the connectivity between the individual cells.
Definition at line 1654 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::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 ) |
Function to create coloring on the second layer within each partition.
Definition at line 2019 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::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 ) |
Function to create partitioning on the second layer within each partition.
Definition at line 1698 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::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 |
This function creates partitions according to the provided connectivity graph.
connectivity | Connectivity between (blocks of cells) |
cluster_size | The number of cells in each partition should be a multiple of cluster_size (for blocking later on) |
cell_partition | Saves of each (block of cells) to which partition the block belongs |
partition_list | partition_list[j] gives the old number of the block that should be renumbered to j due to the partitioning |
partition_size | Vector pointing to start of each partition (on output) |
partition | number of partitions created |
Definition at line 2095 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::update_task_info | ( | const unsigned int | partition | ) |
Update fields of task info for task graph set up in make_thread_graph.
Definition at line 2324 of file task_info.cc.
void internal::MatrixFreeFunctions::TaskInfo::create_flow_graph | ( | ) |
Creates a task graph from a connectivity structure.
std::size_t internal::MatrixFreeFunctions::TaskInfo::memory_consumption | ( | ) | const |
Returns the memory consumption of the class.
Definition at line 705 of file task_info.cc.
template void internal::MatrixFreeFunctions::TaskInfo::print_memory_statistics< ConditionalOStream > | ( | StreamType & | out, |
std::size_t | data_length ) const |
Prints minimum, average, and maximal memory consumption over the MPI processes.
Definition at line 690 of file task_info.cc.
unsigned int internal::MatrixFreeFunctions::TaskInfo::n_active_cells |
Number of physical cells in the mesh, not cell batches after vectorization
Definition at line 432 of file task_info.h.
unsigned int internal::MatrixFreeFunctions::TaskInfo::n_ghost_cells |
Number of physical ghost cells in the mesh which are subject to special treatment and should not be included in loops
Definition at line 438 of file task_info.h.
unsigned int internal::MatrixFreeFunctions::TaskInfo::vectorization_length |
Number of lanes in the SIMD array that are used for vectorization
Definition at line 443 of file task_info.h.
unsigned int internal::MatrixFreeFunctions::TaskInfo::block_size |
Block size information for multithreading
Definition at line 448 of file task_info.h.
unsigned int internal::MatrixFreeFunctions::TaskInfo::n_blocks |
Number of blocks for multithreading
Definition at line 453 of file task_info.h.
TasksParallelScheme internal::MatrixFreeFunctions::TaskInfo::scheme |
Parallel scheme applied by multithreading
Definition at line 458 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::partition_row_index |
The blocks are organized by a vector-of-vector concept, and this data field partition_row_index
stores the distance from one 'vector' to the next within the linear storage of all data to the two-level partitioning.
Definition at line 466 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::cell_partition_data |
This is a linear storage of all partitions, building a range of indices of the form cell_partition_data[idx] to cell_partition_data[idx+1] within the integer list of all cells in MatrixFree, subdivided into chunks by partition_row_index
.
Definition at line 474 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::cell_partition_data_hp |
Like cell_partition_data but with precomputed subranges for each active fe index. The start and end point of a partition is given by cell_partition_data_hp_ptr.
Definition at line 481 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::cell_partition_data_hp_ptr |
Pointers within cell_partition_data_hp, indicating the start and end of a partition.
Definition at line 487 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::face_partition_data |
This is a linear storage of all partitions of inner faces, building a range of indices of the form face_partition_data[idx] to face_partition_data[idx+1] within the integer list of all interior faces in MatrixFree, subdivided into chunks by partition_row_index
.
Definition at line 496 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::face_partition_data_hp |
Like face_partition_data but with precomputed subranges for each active fe index pair. The start and end point of a partition is given by face_partition_data_hp_ptr.
Definition at line 503 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::face_partition_data_hp_ptr |
Pointers within face_partition_data_hp, indicating the start and end of a partition.
Definition at line 509 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::boundary_partition_data |
This is a linear storage of all partitions of boundary faces, building a range of indices of the form boundary_partition_data[idx] to boundary_partition_data[idx+1] within the integer list of all boundary faces in MatrixFree, subdivided into chunks by partition_row_index
.
Definition at line 518 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::boundary_partition_data_hp |
Like boundary_partition_data but with precomputed subranges for each active fe index. The start and end point of a partition is given by boundary_partition_data_hp_ptr.
Definition at line 525 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::boundary_partition_data_hp_ptr |
Pointers within boundary_partition_data_hp, indicating the start and end of a partition.
Definition at line 531 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::ghost_face_partition_data |
This is a linear storage of all partitions of interior faces on boundaries to other processors that are not locally used, building a range of indices of the form ghost_face_partition_data[idx] to ghost_face_partition_data[idx+1] within the integer list of all such faces in MatrixFree, subdivided into chunks by partition_row_index
.
Definition at line 541 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::refinement_edge_face_partition_data |
This is a linear storage of all partitions of faces for multigrid levels that have a coarser neighbor and are only included in certain residual computations but not in smoothing, building a range of indices of the form refinement_edge_face_partition_data[idx] to refinement_edge_face_partition_data[idx+1] within the integer list of all such faces in MatrixFree, subdivided into chunks by partition_row_index
.
Definition at line 552 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::partition_evens |
Thread information (which chunk to start 'even' partitions from) to be handed to the dynamic task scheduler
Definition at line 558 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::partition_odds |
Thread information (which chunk to start 'odd' partitions from) to be handed to the dynamic task scheduler
Definition at line 564 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::partition_n_blocked_workers |
Thread information regarding the dependencies for partitions handed to the dynamic task scheduler
Definition at line 570 of file task_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::TaskInfo::partition_n_workers |
Thread information regarding the dependencies for partitions handed to the dynamic task scheduler
Definition at line 576 of file task_info.h.
unsigned int internal::MatrixFreeFunctions::TaskInfo::evens |
Number of even partitions accumulated over the field partition_evens
Definition at line 582 of file task_info.h.
unsigned int internal::MatrixFreeFunctions::TaskInfo::odds |
Number of odd partitions accumulated over the field partition_odds
Definition at line 588 of file task_info.h.
unsigned int internal::MatrixFreeFunctions::TaskInfo::n_blocked_workers |
Number of blocked workers accumulated over the field partition_n_blocked_workers
Definition at line 594 of file task_info.h.
unsigned int internal::MatrixFreeFunctions::TaskInfo::n_workers |
Number of workers accumulated over the field partition_n_workers
Definition at line 599 of file task_info.h.
std::vector<unsigned char> internal::MatrixFreeFunctions::TaskInfo::task_at_mpi_boundary |
Stores whether a particular task is at an MPI boundary and needs data exchange
Definition at line 605 of file task_info.h.
MPI_Comm internal::MatrixFreeFunctions::TaskInfo::communicator |
MPI communicator
Definition at line 610 of file task_info.h.
MPI_Comm internal::MatrixFreeFunctions::TaskInfo::communicator_sm |
Shared-memory MPI communicator
Definition at line 615 of file task_info.h.
bool internal::MatrixFreeFunctions::TaskInfo::allow_ghosted_vectors_in_loops |
Assert that vectors passed to the MatrixFree loops are not ghosted.
Definition at line 620 of file task_info.h.
unsigned int internal::MatrixFreeFunctions::TaskInfo::my_pid |
Rank of MPI process
Definition at line 625 of file task_info.h.
unsigned int internal::MatrixFreeFunctions::TaskInfo::n_procs |
Number of MPI rank for the current communicator
Definition at line 630 of file task_info.h.