hyper.deal
Loading...
Searching...
No Matches
hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > Class Template Reference

#include <matrix_free.h>

Classes

struct  AdditionalData
 

Public Types

enum class  DataAccessOnFaces { none , values , gradients , unspecified }
 
using ID = TensorID
 
using VectorizedArrayTypeX = VectorizedArrayType
 
using VectorizedArrayTypeV = dealii::VectorizedArray<Number, 1>
 

Public Member Functions

 MatrixFree (const MPI_Comm comm, const MPI_Comm comm_sm, const dealii::MatrixFree< dim_x, Number, VectorizedArrayTypeX > &matrix_free_x, const dealii::MatrixFree< dim_v, Number, VectorizedArrayTypeV > &matrix_free_v)
 
void reinit (const AdditionalData &ad=AdditionalData())
 
template<typename OutVector , typename InVector >
void cell_loop (const std::function< void(const MatrixFree &, OutVector &, const InVector &, const ID)> &cell_operation, OutVector &dst, const InVector &src) const
 
template<typename CLASS , typename OutVector , typename InVector >
void cell_loop (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const ID), CLASS *owning_class, OutVector &dst, const InVector &src) const
 
template<typename OutVector , typename InVector >
void loop_cell_centric (const std::function< void(const MatrixFree &, OutVector &, const InVector &, const ID)> &cell_operation, OutVector &dst, const InVector &src, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified, Timers *timers=nullptr) const
 
template<typename CLASS , typename OutVector , typename InVector >
void loop_cell_centric (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const ID), CLASS *owning_class, OutVector &dst, const InVector &src, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified, Timers *timers=nullptr) const
 
template<typename OutVector , typename InVector >
void loop (const std::function< void(const MatrixFree &, OutVector &, const InVector &, const ID)> &cell_operation, const std::function< void(const MatrixFree &, OutVector &, const InVector &, const ID)> &face_operation, const std::function< void(const MatrixFree &, OutVector &, const InVector &, const ID)> &boundary_operation, OutVector &dst, const InVector &src, const DataAccessOnFaces dst_vector_face_access=DataAccessOnFaces::unspecified, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified, Timers *timers=nullptr) const
 
template<typename CLASS , typename OutVector , typename InVector >
void loop (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const ID), void(CLASS::*face_operation)(const MatrixFree &, OutVector &, const InVector &, const ID), void(CLASS::*boundary_operation)(const MatrixFree &, OutVector &, const InVector &, const ID), CLASS *owning_class, OutVector &dst, const InVector &src, const DataAccessOnFaces dst_vector_face_access=DataAccessOnFaces::unspecified, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified, Timers *timers=nullptr) const
 
const MPI_Comm & get_communicator () const
 
void initialize_dof_vector (dealii::LinearAlgebra::distributed::Vector< Number > &vec, const unsigned int dof_handler_index=0, const bool do_ghosts=true, const bool zero_out_values=true) const
 
dealii::types::boundary_id get_boundary_id (const ID macro_face) const
 
dealii::types::boundary_id get_faces_by_cells_boundary_id (const TensorID &macro_cell, const unsigned int face_number) const
 
bool is_ecl_supported () const
 
bool are_ghost_faces_supported () const
 
const dealii::MatrixFree< dim_x, Number, VectorizedArrayTypeX > & get_matrix_free_x () const
 
const dealii::MatrixFree< dim_v, Number, VectorizedArrayTypeV > & get_matrix_free_v () const
 
const internal::MatrixFreeFunctions::DoFInfoget_dof_info () const
 
const internal::MatrixFreeFunctions::FaceInfoget_face_info () const
 
const internal::MatrixFreeFunctions::ShapeInfo< Number > & get_shape_info () const
 
MemoryConsumption memory_consumption () const
 
const std::shared_ptr< const dealii::internal::MatrixFreeFunctions::VectorDataExchange::Base > & get_vector_partitioner () const
 

Detailed Description

template<int dim_x, int dim_v, typename Number, typename VectorizedArrayType>
class hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >

A matrix-free class for phase-space.

Member Enumeration Documentation

◆ DataAccessOnFaces

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
enum class hyperdeal::MatrixFree::DataAccessOnFaces
strong

Type of access on faces.

Note
Currently only none and values (ghost-face access) are supported.

Constructor & Destructor Documentation

◆ MatrixFree()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::MatrixFree ( const MPI_Comm comm,
const MPI_Comm comm_sm,
const dealii::MatrixFree< dim_x, Number, VectorizedArrayTypeX > & matrix_free_x,
const dealii::MatrixFree< dim_v, Number, VectorizedArrayTypeV > & matrix_free_v )

Constructor (does nothing - see reinit()).

Member Function Documentation

◆ are_ghost_faces_supported()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
bool hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::are_ghost_faces_supported ( ) const

Return if ghost faces or ghost cells are supported.

◆ cell_loop() [1/2]

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
template<typename OutVector , typename InVector >
void hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::cell_loop ( const std::function< void(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &, OutVector &, const InVector &, const ID)> & cell_operation,
OutVector & dst,
const InVector & src ) const

Loop over all cell pairs. No communication performed.

◆ cell_loop() [2/2]

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
template<typename CLASS , typename OutVector , typename InVector >
void hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::cell_loop ( void(CLASS::* cell_operation )(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &, OutVector &, const InVector &, const ID),
CLASS * owning_class,
OutVector & dst,
const InVector & src ) const

The same as above, but without std::function.

◆ get_boundary_id()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
dealii::types::boundary_id hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_boundary_id ( const ID macro_face) const

Return boundary id of face.

◆ get_communicator()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
const MPI_Comm & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_communicator ( ) const

Return global communicator.

◆ get_dof_info()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
const internal::MatrixFreeFunctions::DoFInfo & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_dof_info ( ) const

Return dof info.

◆ get_face_info()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
const internal::MatrixFreeFunctions::FaceInfo & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_face_info ( ) const

Return face info.

◆ get_faces_by_cells_boundary_id()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
dealii::types::boundary_id hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_faces_by_cells_boundary_id ( const TensorID & macro_cell,
const unsigned int face_number ) const

Return boundary id of face of cell.

Note
In constrast to dealii::MatrixFree::get_faces_by_cells_boundary_id() we only return a single number since we require that the cells of a macro cell all have the same type of boundary faces.

◆ get_matrix_free_v()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
const dealii::MatrixFree< dim_v, Number, typename MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::VectorizedArrayTypeV > & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_matrix_free_v ( ) const

Return dealii::MatrixFree for v-space.

◆ get_matrix_free_x()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
const dealii::MatrixFree< dim_x, Number, typename MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::VectorizedArrayTypeX > & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_matrix_free_x ( ) const

Return dealii::MatrixFree for x-space.

◆ get_shape_info()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
const internal::MatrixFreeFunctions::ShapeInfo< Number > & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_shape_info ( ) const

Return shape info.

◆ get_vector_partitioner()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
const std::shared_ptr< const dealii::internal::MatrixFreeFunctions::VectorDataExchange::Base > & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_vector_partitioner ( ) const

Return partitioner of the vectors.

◆ initialize_dof_vector()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
void hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::initialize_dof_vector ( dealii::LinearAlgebra::distributed::Vector< Number > & vec,
const unsigned int dof_handler_index = 0,
const bool do_ghosts = true,
const bool zero_out_values = true ) const

Allocate shared-memory vector. Additionally, all values are zeroed out so that out-of-memory errors are back-trackable.

Note
If this function has been called the first time, also the partitioner is set up.

◆ is_ecl_supported()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
bool hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::is_ecl_supported ( ) const

Is ECL supported?

◆ loop() [1/2]

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
template<typename OutVector , typename InVector >
void hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::loop ( const std::function< void(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &, OutVector &, const InVector &, const ID)> & cell_operation,
const std::function< void(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &, OutVector &, const InVector &, const ID)> & face_operation,
const std::function< void(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &, OutVector &, const InVector &, const ID)> & boundary_operation,
OutVector & dst,
const InVector & src,
const DataAccessOnFaces dst_vector_face_access = DataAccessOnFaces::unspecified,
const DataAccessOnFaces src_vector_face_access = DataAccessOnFaces::unspecified,
Timers * timers = nullptr ) const

Loop over all cell pairs in an face-centric fashion (FCL). It includes a ghost value update of the source vector and a compression of the destination vector.

◆ loop() [2/2]

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
template<typename CLASS , typename OutVector , typename InVector >
void hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::loop ( void(CLASS::* cell_operation )(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &, OutVector &, const InVector &, const ID),
void(CLASS::* face_operation )(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &, OutVector &, const InVector &, const ID),
void(CLASS::* boundary_operation )(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &, OutVector &, const InVector &, const ID),
CLASS * owning_class,
OutVector & dst,
const InVector & src,
const DataAccessOnFaces dst_vector_face_access = DataAccessOnFaces::unspecified,
const DataAccessOnFaces src_vector_face_access = DataAccessOnFaces::unspecified,
Timers * timers = nullptr ) const

Loop over all cell pairs in an face-centric fashion (FCL). It includes a ghost value update of the source vector and a compression of the destination vector.

◆ loop_cell_centric() [1/2]

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
template<typename OutVector , typename InVector >
void hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::loop_cell_centric ( const std::function< void(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &, OutVector &, const InVector &, const ID)> & cell_operation,
OutVector & dst,
const InVector & src,
const DataAccessOnFaces src_vector_face_access = DataAccessOnFaces::unspecified,
Timers * timers = nullptr ) const

Loop over all cell pairs in an element-centric fashion (ECL). It includes a ghost value update of the source vector.

◆ loop_cell_centric() [2/2]

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
template<typename CLASS , typename OutVector , typename InVector >
void hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::loop_cell_centric ( void(CLASS::* cell_operation )(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &, OutVector &, const InVector &, const ID),
CLASS * owning_class,
OutVector & dst,
const InVector & src,
const DataAccessOnFaces src_vector_face_access = DataAccessOnFaces::unspecified,
Timers * timers = nullptr ) const

The same as above, but without std::function.

◆ memory_consumption()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
MemoryConsumption hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::memory_consumption ( ) const

Return an estimate for the memory consumption, in bytes, of this object.

◆ reinit()

template<int dim_x, int dim_v, typename Number , typename VectorizedArrayType >
void hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::reinit ( const AdditionalData & ad = AdditionalData())

Actually setup internal data structures (except partitioner).

this->use_ecl &&

Global ID of cell each face/cell is related to.

Note
(0) FCL (interior); (1) FCL (exterior); (2) cell; (3) ECL (exterior)

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