hyper.deal
|
#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 ¯o_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::DoFInfo & | get_dof_info () const |
const internal::MatrixFreeFunctions::FaceInfo & | get_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 |
A matrix-free class for phase-space.
|
strong |
Type of access on faces.
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()).
bool hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::are_ghost_faces_supported | ( | ) | const |
Return if ghost faces or ghost cells are supported.
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.
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.
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.
const MPI_Comm & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_communicator | ( | ) | const |
Return global communicator.
const internal::MatrixFreeFunctions::DoFInfo & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_dof_info | ( | ) | const |
Return dof info.
const internal::MatrixFreeFunctions::FaceInfo & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_face_info | ( | ) | const |
Return face info.
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.
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.
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.
const internal::MatrixFreeFunctions::ShapeInfo< Number > & hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::get_shape_info | ( | ) | const |
Return shape info.
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.
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.
bool hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::is_ecl_supported | ( | ) | const |
Is ECL supported?
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.
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.
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.
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.
MemoryConsumption hyperdeal::MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::memory_consumption | ( | ) | const |
Return an estimate for the memory consumption, in bytes, of this object.
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.