Reference documentation for deal.II version 9.0.0
|
#include <deal.II/matrix_free/dof_info.h>
Public Types | |
enum | IndexStorageVariants : unsigned char { IndexStorageVariants::full, IndexStorageVariants::interleaved, IndexStorageVariants::contiguous } |
enum | DoFAccessIndex : unsigned char { dof_access_face_interior = 0, dof_access_face_exterior = 1, dof_access_cell = 2 } |
Public Member Functions | |
DoFInfo () | |
DoFInfo (const DoFInfo &)=default | |
void | clear () |
unsigned int | fe_index_from_degree (const unsigned int first_selected_component, const unsigned int fe_degree) const |
void | read_dof_indices (const std::vector< types::global_dof_index > &local_indices, const std::vector< unsigned int > &lexicographic_inv, const ConstraintMatrix &constraints, const unsigned int cell_number, ConstraintValues< double > &constraint_values, bool &cell_at_boundary) |
void | assign_ghosts (const std::vector< unsigned int > &boundary_cells) |
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) |
void | compute_cell_index_compression (const std::vector< unsigned char > &irregular_cells) |
template<int length> | |
void | compute_face_index_compression (const std::vector< FaceToCellTopology< length > > &faces) |
void | make_connectivity_graph (const TaskInfo &task_info, const std::vector< unsigned int > &renumbering, DynamicSparsityPattern &connectivity) const |
void | compute_dof_renumbering (std::vector< types::global_dof_index > &renumbering) |
template<int length> | |
void | compute_vector_zero_access_pattern (const TaskInfo &task_info, const std::vector< FaceToCellTopology< length > > &faces) |
std::size_t | memory_consumption () const |
template<typename StreamType > | |
void | print_memory_consumption (StreamType &out, const TaskInfo &size_info) const |
template<typename Number > | |
void | print (const std::vector< Number > &constraint_pool_data, const std::vector< unsigned int > &constraint_pool_row_index, std::ostream &out) const |
Public Attributes | |
unsigned int | dimension |
unsigned int | vectorization_length |
std::vector< IndexStorageVariants > | index_storage_variants [3] |
std::vector< std::pair< unsigned int, unsigned int > > | row_starts |
std::vector< unsigned int > | dof_indices |
std::vector< std::pair< unsigned short, unsigned short > > | constraint_indicator |
std::vector< unsigned int > | dof_indices_interleaved |
std::vector< unsigned int > | dof_indices_contiguous [3] |
std::vector< unsigned char > | n_vectorization_lanes_filled [3] |
std::shared_ptr< const Utilities::MPI::Partitioner > | vector_partitioner |
std::array< std::shared_ptr< const Utilities::MPI::Partitioner >, 3 > | vector_partitioner_face_variants |
std::vector< unsigned int > | constrained_dofs |
std::vector< unsigned int > | row_starts_plain_indices |
std::vector< unsigned int > | plain_dof_indices |
unsigned int | global_base_element_offset |
unsigned int | n_base_elements |
std::vector< unsigned int > | n_components |
std::vector< unsigned int > | start_components |
std::vector< unsigned int > | component_to_base_index |
std::vector< std::vector< unsigned int > > | component_dof_indices_offset |
std::vector< unsigned int > | dofs_per_cell |
std::vector< unsigned int > | dofs_per_face |
bool | store_plain_indices |
std::vector< unsigned int > | cell_active_fe_index |
unsigned int | max_fe_index |
std::vector< std::vector< unsigned int > > | fe_index_conversion |
std::vector< types::global_dof_index > | ghost_dofs |
std::vector< unsigned int > | vector_zero_range_list_index |
std::vector< unsigned int > | vector_zero_range_list |
Static Public Attributes | |
static const unsigned int | chunk_size_zero_vector = 8192 |
The class that stores the indices of the degrees of freedom for all the cells. Essentially, this is a smart number cache in the style of a DoFHandler that also embeds the description of constraints directly on the cell level without the need to refer to the external ConstraintMatrix.
This class only stores index relations. The weights for hanging node constraints are stored in a different field. This is because a different field allows for the same compressed weight data on different DoFHandlers for vector-valued problems. There, the indices might be constrained differently on different components (e.g. Dirichlet conditions only on selected components), whereas the weights from hanging nodes are the same and need to be stored only once. The combination will be handled in the MatrixFree class.
Definition at line 62 of file dof_info.h.
|
strong |
Enum for various storage variants of the indices. This storage format is used to implement more efficient indexing schemes in case the underlying data structures allow for them, and to inform the access functions in FEEvaluationBase::read_write_operation() on which array to get the data from. One example of more efficient storage is the enum value IndexStorageVariants::contiguous, which means that one can get the indices to all degrees of freedom of a cell by reading only the first index for each cell, whereas all subsequent indices are merely an offset from the first index.
Enumerator | |
---|---|
full | This value indicates that no index compression was found and the only valid storage is to access all indices present on the cell, possibly including constraints. For a cell/face of this index type, the data access in FEEvaluationBase is directed to the array |
interleaved | This value indicates that the indices are interleaved for access with vectorized gather and scatter operation. This storage variant is possible in case there are no constraints on the cell and the indices in the batch of cells are not pointing to the same global index in different slots of a vectorized array (in order to support scatter operations). For a cell/face of this index type, the data access in FEEvaluationBase is directed to the array |
contiguous | This value indicates that the indices within a cell are all contiguous, and one can get the index to the cell by reading that single value for each of the cells in the cell batch. For a cell/face of this index type, the data access in FEEvaluationBase is directed to the array |
Definition at line 218 of file dof_info.h.
enum internal::MatrixFreeFunctions::DoFInfo::DoFAccessIndex : unsigned char |
Enum used to distinguish the data arrays for the vectorization type in cells and faces.
Definition at line 256 of file dof_info.h.
internal::MatrixFreeFunctions::DoFInfo::DoFInfo | ( | ) |
Default empty constructor.
|
default |
Copy constructor.
void internal::MatrixFreeFunctions::DoFInfo::clear | ( | ) |
Clear all data fields in this class.
unsigned int internal::MatrixFreeFunctions::DoFInfo::fe_index_from_degree | ( | const unsigned int | first_selected_component, |
const unsigned int | fe_degree | ||
) | const |
Return the FE index for a given finite element degree. If not in hp mode, this function always returns index 0. If an index is not found in hp mode, it returns numbers::invalid_unsigned_int.
void internal::MatrixFreeFunctions::DoFInfo::read_dof_indices | ( | const std::vector< types::global_dof_index > & | local_indices, |
const std::vector< unsigned int > & | lexicographic_inv, | ||
const ConstraintMatrix & | constraints, | ||
const unsigned int | cell_number, | ||
ConstraintValues< double > & | constraint_values, | ||
bool & | cell_at_boundary | ||
) |
This internal method takes the local indices on a cell and fills them into this class. It resolves the constraints and distributes the results. Ghost indices, i.e., indices that are located on another processor, get a temporary number by this function, and will later be assigned the correct index after all the ghost indices have been collected by the call to assign_ghosts
.
void internal::MatrixFreeFunctions::DoFInfo::assign_ghosts | ( | const std::vector< unsigned int > & | boundary_cells | ) |
This method assigns the correct indices to ghost indices from the temporary numbering employed by the read_dof_indices
function. The numbers are localized with respect to the MPI process, and ghosts start at the end of the locally owned range. This way, we get direct access to all vector entries.
void internal::MatrixFreeFunctions::DoFInfo::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 | ||
) |
This method reorders the way cells are gone through based on a given renumbering of the cells. It also takes vectorization_length
cells together and interprets them as one cell only, as is needed for vectorization.
void internal::MatrixFreeFunctions::DoFInfo::compute_cell_index_compression | ( | const std::vector< unsigned char > & | irregular_cells | ) |
Finds possible compression for the cell indices that we can apply for increased efficiency. Run at the end of reorder_cells.
void internal::MatrixFreeFunctions::DoFInfo::compute_face_index_compression | ( | const std::vector< FaceToCellTopology< length > > & | faces | ) |
Finds possible compression for the face indices that we can apply for increased efficiency. Run at the end of reorder_cells.
void internal::MatrixFreeFunctions::DoFInfo::make_connectivity_graph | ( | const TaskInfo & | task_info, |
const std::vector< unsigned int > & | renumbering, | ||
DynamicSparsityPattern & | connectivity | ||
) | const |
This function computes the connectivity of the currently stored indices in terms of connections between the individual cells and fills the structure into a sparsity pattern.
void internal::MatrixFreeFunctions::DoFInfo::compute_dof_renumbering | ( | std::vector< types::global_dof_index > & | renumbering | ) |
Compute a renumbering of the degrees of freedom to improve the data access patterns for this class that can be utilized by the categories in the IndexStorageVariants enum. For example, the index ordering can be improved for typical DG elements by interleaving the degrees of freedom from batches of cells, which avoids the explicit data transposition in IndexStorageVariants::contiguous. Currently, these more advanced features are not implemented, so there is only limited value of this function.
void internal::MatrixFreeFunctions::DoFInfo::compute_vector_zero_access_pattern | ( | const TaskInfo & | task_info, |
const std::vector< FaceToCellTopology< length > > & | faces | ||
) |
Fills the array that defines how to zero selected ranges in the result vector within the cell loop, filling the two member variables vector_zero_range_list_index
and vector_zero_range_list
.
The intent of this pattern is to zero the vector entries in close temporal proximity to the first access and thus keeping the vector entries in cache.
std::size_t internal::MatrixFreeFunctions::DoFInfo::memory_consumption | ( | ) | const |
Return the memory consumption in bytes of this class.
void internal::MatrixFreeFunctions::DoFInfo::print_memory_consumption | ( | StreamType & | out, |
const TaskInfo & | size_info | ||
) | const |
Prints a detailed summary of memory consumption in the different structures of this class to the given output stream.
void internal::MatrixFreeFunctions::DoFInfo::print | ( | const std::vector< Number > & | constraint_pool_data, |
const std::vector< unsigned int > & | constraint_pool_row_index, | ||
std::ostream & | out | ||
) | const |
Prints a representation of the indices in the class to the given output stream.
|
static |
This value is used to define subranges in the vectors which we can zero inside the MatrixFree::loop() call. The goal is to only clear a part of the vector at a time to keep the values that are zeroed in caches, saving one global vector access for the case where this is applied rather than vector = 0.;
.
Size of the chunk is set to 64 kByte which generally fits to current caches.
Definition at line 74 of file dof_info.h.
unsigned int internal::MatrixFreeFunctions::DoFInfo::dimension |
Stores the dimension of the underlying DoFHandler. Since the indices are not templated, this is the variable that makes the dimension accessible in the (rare) cases it is needed inside this class.
Definition at line 277 of file dof_info.h.
unsigned int internal::MatrixFreeFunctions::DoFInfo::vectorization_length |
For efficiency reasons, always keep a fixed number of cells with similar properties together. This variable controls the number of cells batched together. As opposed to the other classes which are templated on the number type, this class as a pure index container is not templated, so we need to keep the information otherwise contained in VectorizedArray<Number>::n_array_elements.
Definition at line 287 of file dof_info.h.
std::vector<IndexStorageVariants> internal::MatrixFreeFunctions::DoFInfo::index_storage_variants[3] |
Stores the index storage variant of all cell and face batches.
The three arrays given here address the types for the faces decorated as interior (0), the faces decorated with as exterior (1), and the cells (2) according to CellOrFaceAccess.
Definition at line 296 of file dof_info.h.
std::vector<std::pair<unsigned int, unsigned int> > internal::MatrixFreeFunctions::DoFInfo::row_starts |
Stores the rowstart indices of the compressed row storage in the dof_indices
and constraint_indicator
fields. These two fields are always accessed together, so it is simpler to keep just one variable for them. This also obviates keeping two rowstart vectors in sync.
Definition at line 304 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::dof_indices |
Stores the indices of the degrees of freedom for each cell. These indices are computed in MPI-local index space, i.e., each processor stores the locally owned indices as numbers between 0
and n_locally_owned_dofs-1
and ghost indices in the range n_locally_owned_dofs
to n_locally_owned_dofs+n_ghost_dofs
. The translation between this MPI-local index space and the global numbering of degrees of freedom is stored in the vector_partitioner
data structure. This array also includes the indirect contributions from constraints, which are described by the constraint_indicator
field. Because of variable lengths of rows, this would be a vector of a vector. However, we use one contiguous memory region and store the rowstart in the variable row_starts
.
Definition at line 321 of file dof_info.h.
std::vector<std::pair<unsigned short,unsigned short> > internal::MatrixFreeFunctions::DoFInfo::constraint_indicator |
This variable describes the position of constraints in terms of the local numbering of degrees of freedom on a cell. The first number stores the distance from one constrained degree of freedom to the next. This allows to identify the position of constrained DoFs as we loop through the local degrees of freedom of the cell when reading from or writing to a vector. The second number stores the index of the constraint weights, stored in the variable constraint_pool_data.
Definition at line 332 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::dof_indices_interleaved |
Reordered index storage for IndexStorageVariants::interleaved
.
Definition at line 337 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::dof_indices_contiguous[3] |
Compressed index storage for faster access than through dof_indices
used according to the description in IndexStorageVariants.
The three arrays given here address the types for the faces decorated as interior (0), the faces decorated with as exterior (1), and the cells (2) according to CellOrFaceAccess.
Definition at line 347 of file dof_info.h.
std::vector<unsigned char> internal::MatrixFreeFunctions::DoFInfo::n_vectorization_lanes_filled[3] |
Caches the number of indices filled when vectorizing. This information can implicitly deduced from the row_starts data fields, but this field allows for faster access.
The three arrays given here address the types for the faces decorated as interior (0), the faces decorated with as exterior (1), and the cells (2) according to CellOrFaceAccess.
Definition at line 358 of file dof_info.h.
std::shared_ptr<const Utilities::MPI::Partitioner> internal::MatrixFreeFunctions::DoFInfo::vector_partitioner |
This stores the parallel partitioning that can be used to set up vectors. The partitioner includes the description of the local range in the vector, and also includes how the ghosts look like. This enables initialization of vectors based on the DoFInfo field.
Definition at line 366 of file dof_info.h.
std::array<std::shared_ptr<const Utilities::MPI::Partitioner>, 3> internal::MatrixFreeFunctions::DoFInfo::vector_partitioner_face_variants |
This partitioning selects a subset of ghost indices to the full vector partitioner stored in vector_partitioner
. These partitioners are used in specialized loops that only import parts of the ghosted region for reducing the amount of communication. There are three variants of the partitioner initialized, one that queries only the cell values, one that additionally describes the indices for evaluating the function values on the faces, and one that describes the indices for evaluation both the function values and the gradients on the faces adjacent to the locally owned cells.
Definition at line 379 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::constrained_dofs |
This stores a (sorted) list of all locally owned degrees of freedom that are constrained.
Definition at line 385 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::row_starts_plain_indices |
Stores the rowstart indices of the compressed row storage in the plain_dof_indices
fields.
Definition at line 391 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::plain_dof_indices |
Stores the indices of the degrees of freedom for each cell. This array does not include the indirect contributions from constraints, which are included in dof_indices
. Because of variable lengths of rows, this would be a vector of a vector. However, we use one contiguous memory region and store the rowstart in the variable row_starts_plain_indices
.
Definition at line 401 of file dof_info.h.
unsigned int internal::MatrixFreeFunctions::DoFInfo::global_base_element_offset |
Stores the offset in terms of the number of base elements over all DoFInfo objects.
Definition at line 407 of file dof_info.h.
unsigned int internal::MatrixFreeFunctions::DoFInfo::n_base_elements |
Stores the number of base elements in the DoFHandler where the indices have been read from.
Definition at line 413 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::n_components |
Stores the number of components of each base element in the finite element where the indices have been read from.
Definition at line 419 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::start_components |
The ith entry of this vector stores the component number of the given base element.
Definition at line 425 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::component_to_base_index |
For a given component in an FESystem, this variable tells which base element the index belongs to.
Definition at line 431 of file dof_info.h.
std::vector<std::vector<unsigned int> > internal::MatrixFreeFunctions::DoFInfo::component_dof_indices_offset |
For a vector-valued element, this gives the constant offset in the number of degrees of freedom starting at the given component, as the degrees are numbered by degrees of freedom. This data structure does not take possible constraints and thus, shorter or longer lists, into account. This information is encoded in the row_starts variables directly.
The outer vector goes through the various fe indices in the hp case, similarly to the dofs_per_cell
variable.
Definition at line 444 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::dofs_per_cell |
Stores the number of degrees of freedom per cell.
Definition at line 449 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::dofs_per_face |
Stores the number of degrees of freedom per face.
Definition at line 454 of file dof_info.h.
bool internal::MatrixFreeFunctions::DoFInfo::store_plain_indices |
Informs on whether plain indices are cached.
Definition at line 459 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::cell_active_fe_index |
Stores the index of the active finite element in the hp case.
Definition at line 464 of file dof_info.h.
unsigned int internal::MatrixFreeFunctions::DoFInfo::max_fe_index |
Stores the maximum degree of different finite elements for the hp case.
Definition at line 470 of file dof_info.h.
std::vector<std::vector<unsigned int> > internal::MatrixFreeFunctions::DoFInfo::fe_index_conversion |
To each of the slots in an hp adaptive case, the inner vector stores the corresponding element degree. This is used by the constructor of FEEvaluationBase to identify the correct data slot in the hp case.
Definition at line 477 of file dof_info.h.
std::vector<types::global_dof_index> internal::MatrixFreeFunctions::DoFInfo::ghost_dofs |
Temporarily stores the numbers of ghosts during setup. Cleared when calling assign_ghosts
. Then, all information is collected by the partitioner.
Definition at line 484 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::vector_zero_range_list_index |
Stores an integer to each partition in TaskInfo that indicates whether to clear certain parts in the result vector if the user requested it with the respective argument in the MatrixFree::loop.
Definition at line 491 of file dof_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::DoFInfo::vector_zero_range_list |
Stores the actual ranges in the vector to be cleared.
Definition at line 496 of file dof_info.h.