Reference documentation for deal.II version 8.5.1
Public Types | Public Member Functions | Public Attributes | List of all members
MatrixFree< dim, Number >::AdditionalData Struct Reference

#include <deal.II/matrix_free/matrix_free.h>

Public Types

enum  TasksParallelScheme { none, partition_partition, partition_color, color }
 

Public Member Functions

 AdditionalData (const TasksParallelScheme tasks_parallel_scheme=partition_partition, const unsigned int tasks_block_size=0, const UpdateFlags mapping_update_flags=update_gradients|update_JxW_values, const unsigned int level_mg_handler=numbers::invalid_unsigned_int, const bool store_plain_indices=true, const bool initialize_indices=true, const bool initialize_mapping=true)
 
 AdditionalData (const MPI_Comm mpi_communicator, const TasksParallelScheme tasks_parallel_scheme=partition_partition, const unsigned int tasks_block_size=0, const UpdateFlags mapping_update_flags=update_gradients|update_JxW_values, const unsigned int level_mg_handler=numbers::invalid_unsigned_int, const bool store_plain_indices=true, const bool initialize_indices=true, const bool initialize_mapping=true)
 

Public Attributes

MPI_Comm mpi_communicator
 
TasksParallelScheme tasks_parallel_scheme
 
unsigned int tasks_block_size
 
UpdateFlags mapping_update_flags
 
unsigned int level_mg_handler
 
bool store_plain_indices
 
bool initialize_indices
 
bool initialize_mapping
 

Detailed Description

template<int dim, typename Number = double>
struct MatrixFree< dim, Number >::AdditionalData

Collects the options for initialization of the MatrixFree class. The first parameter specifies the MPI communicator to be used, the second the parallelization options in shared memory (task-based parallelism, where one can choose between no parallelism and three schemes that avoid that cells with access to the same vector entries are accessed simultaneously), the third with the block size for task parallel scheduling, the fourth the update flags that should be stored by this class.

The fifth parameter specifies the level in the triangulation from which the indices are to be used. If the level is set to numbers::invalid_unsigned_int, the active cells are traversed, and otherwise the cells in the given level. This option has no effect in case a DoFHandler or hp::DoFHandler is given.

The parameter initialize_plain_indices indicates whether the DoFInfo class should also allow for access to vectors without resolving constraints.

The last two parameters allow the user to disable some of the initialization processes. For example, if only the scheduling that avoids touching the same vector/matrix indices simultaneously is to be found, the mapping needs not be initialized. Likewise, if the mapping has changed from one iteration to the next but the topology has not (like when using a deforming mesh with MappingQEulerian), it suffices to initialize the mapping only.

Definition at line 141 of file matrix_free.h.

Member Enumeration Documentation

◆ TasksParallelScheme

template<int dim, typename Number = double>
enum MatrixFree::AdditionalData::TasksParallelScheme

Collects options for task parallelism. See the documentation of the member variable MatrixFree::AdditionalData::tasks_parallel_scheme for a thorough description.

Enumerator
none 

Perform application in serial.

partition_partition 

Partition the cells into two levels and afterwards form chunks.

partition_color 

Partition on the global level and color cells within the partitions.

color 

Use the traditional coloring algorithm: this is like TasksParallelScheme::partition_color, but only uses one partition.

Definition at line 148 of file matrix_free.h.

Constructor & Destructor Documentation

◆ AdditionalData() [1/2]

template<int dim, typename Number = double>
MatrixFree< dim, Number >::AdditionalData::AdditionalData ( const TasksParallelScheme  tasks_parallel_scheme = partition_partition,
const unsigned int  tasks_block_size = 0,
const UpdateFlags  mapping_update_flags = update_gradients | update_JxW_values,
const unsigned int  level_mg_handler = numbers::invalid_unsigned_int,
const bool  store_plain_indices = true,
const bool  initialize_indices = true,
const bool  initialize_mapping = true 
)
inline

Constructor for AdditionalData.

Definition at line 175 of file matrix_free.h.

◆ AdditionalData() [2/2]

template<int dim, typename Number = double>
MatrixFree< dim, Number >::AdditionalData::AdditionalData ( const MPI_Comm  mpi_communicator,
const TasksParallelScheme  tasks_parallel_scheme = partition_partition,
const unsigned int  tasks_block_size = 0,
const UpdateFlags  mapping_update_flags = update_gradients | update_JxW_values,
const unsigned int  level_mg_handler = numbers::invalid_unsigned_int,
const bool  store_plain_indices = true,
const bool  initialize_indices = true,
const bool  initialize_mapping = true 
)
inline

Constructor for AdditionalData.

Deprecated:
mpi_communicator should not be specified here since it is not used in the MatrixFree class.

Definition at line 199 of file matrix_free.h.

Member Data Documentation

◆ mpi_communicator

template<int dim, typename Number = double>
MPI_Comm MatrixFree< dim, Number >::AdditionalData::mpi_communicator

Set the MPI communicator that the parallel layout of the operator should be based upon. Defaults to MPI_COMM_SELF, but should be set to a communicator similar to the one used for a distributed triangulation in order to inform this class over all cells that are present.

Deprecated:
This variable is not used anymore. The mpi_communicator is automatically set to the one used by the triangulation.

Definition at line 227 of file matrix_free.h.

◆ tasks_parallel_scheme

template<int dim, typename Number = double>
TasksParallelScheme MatrixFree< dim, Number >::AdditionalData::tasks_parallel_scheme

Set the scheme for task parallelism. There are four options available. If set to none, the operator application is done in serial without shared memory parallelism. If this class is used together with MPI and MPI is also used for parallelism within the nodes, this flag should be set to none. The default value is partition_partition, i.e. we actually use multithreading with the first option described below.

The first option partition_partition is to partition the cells on two levels in onion-skin-like partitions and forming chunks of tasks_block_size after the partitioning. The partitioning finds sets of independent cells that enable working in parallel without accessing the same vector entries at the same time.

The second option partition_color is to use a partition on the global level and color cells within the partitions (where all chunks within a color are independent). Here, the subdivision into chunks of cells is done before the partitioning, which might give worse partitions but better cache performance if degrees of freedom are not renumbered.

The third option color is to use a traditional algorithm of coloring on the global level. This scheme is a special case of the second option where only one partition is present. Note that for problems with hanging nodes, there are quite many colors (50 or more in 3D), which might degrade parallel performance (bad cache behavior, many synchronization points).

Definition at line 257 of file matrix_free.h.

◆ tasks_block_size

template<int dim, typename Number = double>
unsigned int MatrixFree< dim, Number >::AdditionalData::tasks_block_size

Set the number of so-called macro cells that should form one partition. If zero size is given, the class tries to find a good size for the blocks based on MultithreadInfo::n_threads() and the number of cells present. Otherwise, the given number is used. If the given number is larger than one third of the number of total cells, this means no parallelism. Note that in the case vectorization is used, a macro cell consists of more than one physical cell.

Definition at line 268 of file matrix_free.h.

◆ mapping_update_flags

template<int dim, typename Number = double>
UpdateFlags MatrixFree< dim, Number >::AdditionalData::mapping_update_flags

This flag is used to determine which quantities should be cached. This class can cache data needed for gradient computations (inverse Jacobians), Jacobian determinants (JxW), quadrature points as well as data for Hessians (derivative of Jacobians). By default, only data for gradients and Jacobian determinants times quadrature weights, JxW, are cached. If quadrature points or second derivatives are needed, they must be specified by this field (even though second derivatives might still be evaluated on Cartesian cells without this option set here, since there the Jacobian describes the mapping completely).

Definition at line 281 of file matrix_free.h.

◆ level_mg_handler

template<int dim, typename Number = double>
unsigned int MatrixFree< dim, Number >::AdditionalData::level_mg_handler

This option can be used to define whether we work on a certain level of the mesh, and not the active cells. If set to invalid_unsigned_int (which is the default value), the active cells are gone through, otherwise the level given by this parameter. Note that if you specify to work on a level, its dofs must be distributed by using dof_handler.distribute_mg_dofs(fe);.

Definition at line 291 of file matrix_free.h.

◆ store_plain_indices

template<int dim, typename Number = double>
bool MatrixFree< dim, Number >::AdditionalData::store_plain_indices

Controls whether to allow reading from vectors without resolving constraints, i.e., just read the local values of the vector. By default, this option is disabled, so if you want to use FEEvaluationBase::read_dof_values_plain, this flag needs to be set.

Definition at line 299 of file matrix_free.h.

◆ initialize_indices

template<int dim, typename Number = double>
bool MatrixFree< dim, Number >::AdditionalData::initialize_indices

Option to control whether the indices stored in the DoFHandler should be read and the pattern for task parallelism should be set up in the initialize method of MatrixFree. Defaults to true. Can be disabled in case the mapping should be recomputed (e.g. when using a deforming mesh described through MappingEulerian) but the topology of cells has remained the same.

Definition at line 309 of file matrix_free.h.

◆ initialize_mapping

template<int dim, typename Number = double>
bool MatrixFree< dim, Number >::AdditionalData::initialize_mapping

Option to control whether the mapping information should be computed in the initialize method of MatrixFree. Defaults to true. Can be disabled when only some indices should be set up (e.g. when only a set of independent cells should be computed).

Definition at line 317 of file matrix_free.h.


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