15#ifndef dealii_scalapack_h
16#define dealii_scalapack_h
20#ifdef DEAL_II_WITH_SCALAPACK
74template <
typename NumberType>
94 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
109 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
127 const std::string &filename,
128 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
149 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
163 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
206 const unsigned int rank);
263 const std::pair<unsigned int, unsigned int> &offset_A,
264 const std::pair<unsigned int, unsigned int> &offset_B,
265 const std::pair<unsigned int, unsigned int> &submatrix_size)
const;
291 const NumberType a = 0.,
292 const NumberType b = 1.,
293 const bool transpose_B =
false);
340 mult(
const NumberType b,
344 const bool transpose_A =
false,
345 const bool transpose_B =
false)
const;
367 const bool adding =
false)
const;
389 const bool adding =
false)
const;
411 const bool adding =
false)
const;
434 const bool adding =
false)
const;
457 save(
const std::string &filename,
458 const std::pair<unsigned int, unsigned int> &chunk_size =
474 load(
const std::string &filename);
522 std::vector<NumberType>
524 const std::pair<unsigned int, unsigned int> &index_limits,
525 const bool compute_eigenvectors);
535 std::vector<NumberType>
537 const std::pair<NumberType, NumberType> &value_limits,
538 const bool compute_eigenvectors);
556 std::vector<NumberType>
558 const std::pair<unsigned int, unsigned int> &index_limits,
559 const bool compute_eigenvectors);
571 std::vector<NumberType>
573 const std::pair<NumberType, NumberType> &value_limits,
574 const bool compute_eigenvectors);
604 std::vector<NumberType>
752 local_el(
const unsigned int loc_row,
const unsigned int loc_column)
const;
758 local_el(
const unsigned int loc_row,
const unsigned int loc_column);
770 template <
class InputVector>
784 template <
class InputVector>
813 std::vector<NumberType>
815 const bool compute_eigenvectors,
816 const std::pair<unsigned int, unsigned int> &index_limits =
819 const std::pair<NumberType, NumberType> &value_limits =
820 std::make_pair(std::numeric_limits<NumberType>::quiet_NaN(),
821 std::numeric_limits<NumberType>::quiet_NaN()));
842 std::vector<NumberType>
844 const bool compute_eigenvectors,
845 const std::pair<unsigned int, unsigned int> &index_limits =
848 const std::pair<NumberType, NumberType> &value_limits =
849 std::make_pair(std::numeric_limits<NumberType>::quiet_NaN(),
850 std::numeric_limits<NumberType>::quiet_NaN()));
858 const std::pair<unsigned int, unsigned int> &chunk_size)
const;
873 const std::pair<unsigned int, unsigned int> &chunk_size)
const;
899 std::shared_ptr<const Utilities::MPI::ProcessGrid>
grid;
939 mutable std::vector<NumberType>
work;
992template <
typename NumberType>
995 const unsigned int loc_column)
const
997 return (*
this)(loc_row, loc_column);
1002template <
typename NumberType>
1005 const unsigned int loc_column)
1007 return (*
this)(loc_row, loc_column);
1011template <
typename NumberType>
1020template <
typename NumberType>
1029template <
typename NumberType>
1033 return n_local_rows;
1038template <
typename NumberType>
1042 return n_local_columns;
std::vector< NumberType > compute_SVD(ScaLAPACKMatrix< NumberType > *U=nullptr, ScaLAPACKMatrix< NumberType > *VT=nullptr)
std::vector< NumberType > eigenpairs_symmetric_by_value(const std::pair< NumberType, NumberType > &value_limits, const bool compute_eigenvectors)
unsigned int local_n() const
NumberType frobenius_norm() const
unsigned int pseudoinverse(const NumberType ratio)
std::vector< NumberType > eigenpairs_symmetric_by_value_MRRR(const std::pair< NumberType, NumberType > &value_limits, const bool compute_eigenvectors)
const int first_process_column
void copy_from(const LAPACKFullMatrix< NumberType > &matrix, const unsigned int rank)
void save_parallel(const std::string &filename, const std::pair< unsigned int, unsigned int > &chunk_size) const
void least_squares(ScaLAPACKMatrix< NumberType > &B, const bool transpose=false)
NumberType local_el(const unsigned int loc_row, const unsigned int loc_column) const
ScaLAPACKMatrix< NumberType > & operator=(const FullMatrix< NumberType > &)
void Tadd(const NumberType b, const ScaLAPACKMatrix< NumberType > &B)
void mTmult(ScaLAPACKMatrix< NumberType > &C, const ScaLAPACKMatrix< NumberType > &B, const bool adding=false) const
void add(const ScaLAPACKMatrix< NumberType > &B, const NumberType a=0., const NumberType b=1., const bool transpose_B=false)
LAPACKSupport::State get_state() const
LAPACKSupport::Property get_property() const
void Tmmult(ScaLAPACKMatrix< NumberType > &C, const ScaLAPACKMatrix< NumberType > &B, const bool adding=false) const
std::vector< NumberType > eigenpairs_symmetric_MRRR(const bool compute_eigenvectors, const std::pair< unsigned int, unsigned int > &index_limits=std::make_pair(numbers::invalid_unsigned_int, numbers::invalid_unsigned_int), const std::pair< NumberType, NumberType > &value_limits=std::make_pair(std::numeric_limits< NumberType >::quiet_NaN(), std::numeric_limits< NumberType >::quiet_NaN()))
void scale_rows(const InputVector &factors)
const int first_process_row
NumberType & local_el(const unsigned int loc_row, const unsigned int loc_column)
std::shared_ptr< const Utilities::MPI::ProcessGrid > grid
void load(const std::string &filename)
const int submatrix_column
void mmult(ScaLAPACKMatrix< NumberType > &C, const ScaLAPACKMatrix< NumberType > &B, const bool adding=false) const
std::vector< NumberType > eigenpairs_symmetric(const bool compute_eigenvectors, const std::pair< unsigned int, unsigned int > &index_limits=std::make_pair(numbers::invalid_unsigned_int, numbers::invalid_unsigned_int), const std::pair< NumberType, NumberType > &value_limits=std::make_pair(std::numeric_limits< NumberType >::quiet_NaN(), std::numeric_limits< NumberType >::quiet_NaN()))
void save_serial(const std::string &filename, const std::pair< unsigned int, unsigned int > &chunk_size) const
NumberType norm_general(const char type) const
~ScaLAPACKMatrix() override=default
std::vector< NumberType > work
void save(const std::string &filename, const std::pair< unsigned int, unsigned int > &chunk_size=std::make_pair(numbers::invalid_unsigned_int, numbers::invalid_unsigned_int)) const
void load_parallel(const std::string &filename)
NumberType l1_norm() const
void compute_lu_factorization()
NumberType norm_symmetric(const char type) const
void mult(const NumberType b, const ScaLAPACKMatrix< NumberType > &B, const NumberType c, ScaLAPACKMatrix< NumberType > &C, const bool transpose_A=false, const bool transpose_B=false) const
LAPACKSupport::Property property
void set_property(const LAPACKSupport::Property property)
void reinit(const size_type n_rows, const size_type n_columns, const std::shared_ptr< const Utilities::MPI::ProcessGrid > &process_grid, const size_type row_block_size=32, const size_type column_block_size=32, const LAPACKSupport::Property property=LAPACKSupport::Property::general)
void load_serial(const std::string &filename)
std::vector< NumberType > eigenpairs_symmetric_by_index_MRRR(const std::pair< unsigned int, unsigned int > &index_limits, const bool compute_eigenvectors)
NumberType reciprocal_condition_number(const NumberType a_norm) const
LAPACKSupport::State state
void TmTmult(ScaLAPACKMatrix< NumberType > &C, const ScaLAPACKMatrix< NumberType > &B, const bool adding=false) const
std::vector< NumberType > eigenpairs_symmetric_by_index(const std::pair< unsigned int, unsigned int > &index_limits, const bool compute_eigenvectors)
unsigned int global_column(const unsigned int loc_column) const
void copy_to(FullMatrix< NumberType > &matrix) const
unsigned int global_row(const unsigned int loc_row) const
unsigned int local_m() const
void compute_cholesky_factorization()
void copy_transposed(const ScaLAPACKMatrix< NumberType > &B)
NumberType linfty_norm() const
void scale_columns(const InputVector &factors)
const TableIndices< N > & size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ symmetric
Matrix is symmetric.
@ general
No special properties.
static const unsigned int invalid_unsigned_int