Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | List of all members
internal::MatrixFreeFunctions::TaskInfo Struct Reference

#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
 

Public Attributes

unsigned int n_active_cells
 
unsigned int n_ghost_cells
 
unsigned int vectorization_length
 
unsigned int block_size
 
unsigned int n_blocks
 
TasksParallelScheme scheme
 
std::vector< unsigned intpartition_row_index
 
std::vector< unsigned intcell_partition_data
 
std::vector< unsigned intcell_partition_data_hp
 
std::vector< unsigned intcell_partition_data_hp_ptr
 
std::vector< unsigned intface_partition_data
 
std::vector< unsigned intface_partition_data_hp
 
std::vector< unsigned intface_partition_data_hp_ptr
 
std::vector< unsigned intboundary_partition_data
 
std::vector< unsigned intboundary_partition_data_hp
 
std::vector< unsigned intboundary_partition_data_hp_ptr
 
std::vector< unsigned intghost_face_partition_data
 
std::vector< unsigned intrefinement_edge_face_partition_data
 
std::vector< unsigned intpartition_evens
 
std::vector< unsigned intpartition_odds
 
std::vector< unsigned intpartition_n_blocked_workers
 
std::vector< unsigned intpartition_n_workers
 
unsigned int evens
 
unsigned int odds
 
unsigned int n_blocked_workers
 
unsigned int n_workers
 
std::vector< unsigned char > task_at_mpi_boundary
 
MPI_Comm communicator
 
MPI_Comm communicator_sm
 
bool allow_ghosted_vectors_in_loops
 
unsigned int my_pid
 
unsigned int n_procs
 

Detailed Description

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 107 of file task_info.h.

Member Enumeration Documentation

◆ TasksParallelScheme

Enumerator
none 
partition_partition 
partition_color 
color 

Definition at line 113 of file task_info.h.

Constructor & Destructor Documentation

◆ TaskInfo()

internal::MatrixFreeFunctions::TaskInfo::TaskInfo ( )

Constructor.

Definition at line 653 of file task_info.cc.

Member Function Documentation

◆ clear()

void internal::MatrixFreeFunctions::TaskInfo::clear ( )

Clears all the data fields and resets them to zero.

Definition at line 661 of file task_info.cc.

◆ loop()

void internal::MatrixFreeFunctions::TaskInfo::loop ( MFWorkerInterface worker) const

Runs the matrix-free loop.

Definition at line 351 of file task_info.cc.

◆ make_boundary_cells_divisible()

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 723 of file task_info.cc.

◆ create_blocks_serial()

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.

Parameters
cells_with_commA 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_cellGives 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_hpDefines whether cell_vectorization_categories is originating from a hp-adaptive computation with variable polynomial degree or a user-defined variant.
cell_vectorization_categoriesThis 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_strictDefines 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_relationThis 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.
renumberingWhen leaving this function, the vector contains a new numbering of the cells that aligns with the grouping stored in this class.
incompletely_filled_vectorizationGiven 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 795 of file task_info.cc.

◆ initial_setup_blocks_tasks()

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.

Parameters
boundary_cellsA list of cells that need to exchange data prior to performing computations. These will be given a certain id in the partitioning.
renumberingWhen 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_vectorizationGiven 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 1151 of file task_info.cc.

◆ guess_block_size()

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 1217 of file task_info.cc.

◆ make_thread_graph_partition_color()

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.

Parameters
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_boolDefines whether we are in hp-mode or not

Definition at line 1245 of file task_info.cc.

◆ make_thread_graph_partition_partition()

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.

Parameters
cell_active_fe_indexThe 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_boolDefines whether we are in hp-mode or not

Definition at line 1587 of file task_info.cc.

◆ make_thread_graph()

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.

Parameters
cell_active_fe_indexThe 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_boolDefines whether we are in hp-mode or not

Definition at line 1392 of file task_info.cc.

◆ make_connectivity_cells_to_blocks()

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 1655 of file task_info.cc.

◆ make_coloring_within_partitions_pre_blocked()

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 2020 of file task_info.cc.

◆ make_partitioning_within_partitions_post_blocked()

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 1699 of file task_info.cc.

◆ make_partitioning()

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.

Parameters
connectivityConnectivity between (blocks of cells)
cluster_sizeThe number of cells in each partition should be a multiple of cluster_size (for blocking later on)
cell_partitionSaves of each (block of cells) to which partition the block belongs
partition_listpartition_list[j] gives the old number of the block that should be renumbered to j due to the partitioning
partition_sizeVector pointing to start of each partition (on output)
partitionnumber of partitions created

Definition at line 2096 of file task_info.cc.

◆ update_task_info()

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 2325 of file task_info.cc.

◆ create_flow_graph()

void internal::MatrixFreeFunctions::TaskInfo::create_flow_graph ( )

Creates a task graph from a connectivity structure.

◆ memory_consumption()

std::size_t internal::MatrixFreeFunctions::TaskInfo::memory_consumption ( ) const

Returns the memory consumption of the class.

Definition at line 706 of file task_info.cc.

◆ print_memory_statistics()

template<typename StreamType >
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 691 of file task_info.cc.

Member Data Documentation

◆ n_active_cells

unsigned int internal::MatrixFreeFunctions::TaskInfo::n_active_cells

Number of physical cells in the mesh, not cell batches after vectorization

Definition at line 433 of file task_info.h.

◆ n_ghost_cells

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 439 of file task_info.h.

◆ vectorization_length

unsigned int internal::MatrixFreeFunctions::TaskInfo::vectorization_length

Number of lanes in the SIMD array that are used for vectorization

Definition at line 444 of file task_info.h.

◆ block_size

unsigned int internal::MatrixFreeFunctions::TaskInfo::block_size

Block size information for multithreading

Definition at line 449 of file task_info.h.

◆ n_blocks

unsigned int internal::MatrixFreeFunctions::TaskInfo::n_blocks

Number of blocks for multithreading

Definition at line 454 of file task_info.h.

◆ scheme

TasksParallelScheme internal::MatrixFreeFunctions::TaskInfo::scheme

Parallel scheme applied by multithreading

Definition at line 459 of file task_info.h.

◆ partition_row_index

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 467 of file task_info.h.

◆ cell_partition_data

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 475 of file task_info.h.

◆ cell_partition_data_hp

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 482 of file task_info.h.

◆ cell_partition_data_hp_ptr

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 488 of file task_info.h.

◆ face_partition_data

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 497 of file task_info.h.

◆ face_partition_data_hp

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 504 of file task_info.h.

◆ face_partition_data_hp_ptr

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 510 of file task_info.h.

◆ boundary_partition_data

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 519 of file task_info.h.

◆ boundary_partition_data_hp

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 526 of file task_info.h.

◆ boundary_partition_data_hp_ptr

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 532 of file task_info.h.

◆ ghost_face_partition_data

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 542 of file task_info.h.

◆ refinement_edge_face_partition_data

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 553 of file task_info.h.

◆ partition_evens

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 559 of file task_info.h.

◆ partition_odds

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 565 of file task_info.h.

◆ partition_n_blocked_workers

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 571 of file task_info.h.

◆ partition_n_workers

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 577 of file task_info.h.

◆ evens

unsigned int internal::MatrixFreeFunctions::TaskInfo::evens

Number of even partitions accumulated over the field partition_evens

Definition at line 583 of file task_info.h.

◆ odds

unsigned int internal::MatrixFreeFunctions::TaskInfo::odds

Number of odd partitions accumulated over the field partition_odds

Definition at line 589 of file task_info.h.

◆ n_blocked_workers

unsigned int internal::MatrixFreeFunctions::TaskInfo::n_blocked_workers

Number of blocked workers accumulated over the field partition_n_blocked_workers

Definition at line 595 of file task_info.h.

◆ n_workers

unsigned int internal::MatrixFreeFunctions::TaskInfo::n_workers

Number of workers accumulated over the field partition_n_workers

Definition at line 600 of file task_info.h.

◆ task_at_mpi_boundary

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 606 of file task_info.h.

◆ communicator

MPI_Comm internal::MatrixFreeFunctions::TaskInfo::communicator

MPI communicator

Definition at line 611 of file task_info.h.

◆ communicator_sm

MPI_Comm internal::MatrixFreeFunctions::TaskInfo::communicator_sm

Shared-memory MPI communicator

Definition at line 616 of file task_info.h.

◆ allow_ghosted_vectors_in_loops

bool internal::MatrixFreeFunctions::TaskInfo::allow_ghosted_vectors_in_loops

Assert that vectors passed to the MatrixFree loops are not ghosted.

Definition at line 621 of file task_info.h.

◆ my_pid

unsigned int internal::MatrixFreeFunctions::TaskInfo::my_pid

Rank of MPI process

Definition at line 626 of file task_info.h.

◆ n_procs

unsigned int internal::MatrixFreeFunctions::TaskInfo::n_procs

Number of MPI rank for the current communicator

Definition at line 631 of file task_info.h.


The documentation for this struct was generated from the following files: