Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
scalapack.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_scalapack_h
17 #define dealii_scalapack_h
18 
19 #include <deal.II/base/config.h>
20 
21 #ifdef DEAL_II_WITH_SCALAPACK
22 
23 # include <deal.II/base/exceptions.h>
24 # include <deal.II/base/mpi.h>
27 
28 # include <deal.II/lac/full_matrix.h>
31 
32 # include <mpi.h>
33 
34 # include <memory>
35 
37 
76 template <typename NumberType>
77 class ScaLAPACKMatrix : protected TransposeTable<NumberType>
78 {
79 public:
83  using size_type = unsigned int;
84 
94  const size_type n_rows,
95  const size_type n_columns,
96  const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
97  const size_type row_block_size = 32,
98  const size_type column_block_size = 32,
100 
111  const size_type size,
112  const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
113  const size_type block_size = 32,
116 
130  const std::string & filename,
131  const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
132  const size_type row_block_size = 32,
133  const size_type column_block_size = 32);
134 
138  ~ScaLAPACKMatrix() override = default;
139 
148  void
149  reinit(
150  const size_type n_rows,
151  const size_type n_columns,
152  const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
153  const size_type row_block_size = 32,
154  const size_type column_block_size = 32,
156 
164  void
165  reinit(const size_type size,
166  const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
167  const size_type block_size = 32,
170 
174  void
176 
181  get_property() const;
182 
187  get_state() const;
188 
197 
207  void
209  const unsigned int rank);
210 
217  void
219 
229  void
230  copy_to(LAPACKFullMatrix<NumberType> &matrix, const unsigned int rank) const;
231 
237  void
239 
264  void
266  const std::pair<unsigned int, unsigned int> &offset_A,
267  const std::pair<unsigned int, unsigned int> &offset_B,
268  const std::pair<unsigned int, unsigned int> &submatrix_size) const;
269 
278  void
280 
292  void
294  const NumberType a = 0.,
295  const NumberType b = 1.,
296  const bool transpose_B = false);
297 
307  void
308  add(const NumberType b, const ScaLAPACKMatrix<NumberType> &B);
309 
319  void
320  Tadd(const NumberType b, const ScaLAPACKMatrix<NumberType> &B);
321 
342  void
343  mult(const NumberType b,
345  const NumberType c,
347  const bool transpose_A = false,
348  const bool transpose_B = false) const;
349 
367  void
370  const bool adding = false) const;
371 
389  void
392  const bool adding = false) const;
393 
411  void
414  const bool adding = false) const;
415 
434  void
437  const bool adding = false) const;
438 
459  void
460  save(const std::string & filename,
461  const std::pair<unsigned int, unsigned int> &chunk_size =
462  std::make_pair(numbers::invalid_unsigned_int,
464 
476  void
477  load(const std::string &filename);
478 
484  void
486 
492  void
494 
507  void
508  invert();
509 
525  std::vector<NumberType>
527  const std::pair<unsigned int, unsigned int> &index_limits,
528  const bool compute_eigenvectors);
529 
538  std::vector<NumberType>
540  const std::pair<NumberType, NumberType> &value_limits,
541  const bool compute_eigenvectors);
542 
559  std::vector<NumberType>
561  const std::pair<unsigned int, unsigned int> &index_limits,
562  const bool compute_eigenvectors);
563 
574  std::vector<NumberType>
576  const std::pair<NumberType, NumberType> &value_limits,
577  const bool compute_eigenvectors);
578 
607  std::vector<NumberType>
609  ScaLAPACKMatrix<NumberType> *VT = nullptr);
610 
651  void
652  least_squares(ScaLAPACKMatrix<NumberType> &B, const bool transpose = false);
653 
677  unsigned int
678  pseudoinverse(const NumberType ratio);
679 
694  NumberType
695  reciprocal_condition_number(const NumberType a_norm) const;
696 
700  NumberType
701  l1_norm() const;
702 
706  NumberType
707  linfty_norm() const;
708 
712  NumberType
713  frobenius_norm() const;
714 
718  size_type
719  m() const;
720 
724  size_type
725  n() const;
726 
730  unsigned int
731  local_m() const;
732 
736  unsigned int
737  local_n() const;
738 
742  unsigned int
743  global_row(const unsigned int loc_row) const;
744 
748  unsigned int
749  global_column(const unsigned int loc_column) const;
750 
754  NumberType
755  local_el(const unsigned int loc_row, const unsigned int loc_column) const;
756 
760  NumberType &
761  local_el(const unsigned int loc_row, const unsigned int loc_column);
762 
773  template <class InputVector>
774  void
775  scale_columns(const InputVector &factors);
776 
787  template <class InputVector>
788  void
789  scale_rows(const InputVector &factors);
790 
791 private:
796  NumberType
797  norm_symmetric(const char type) const;
798 
803  NumberType
804  norm_general(const char type) const;
805 
816  std::vector<NumberType>
818  const bool compute_eigenvectors,
819  const std::pair<unsigned int, unsigned int> &index_limits =
820  std::make_pair(numbers::invalid_unsigned_int,
822  const std::pair<NumberType, NumberType> &value_limits =
823  std::make_pair(std::numeric_limits<NumberType>::quiet_NaN(),
824  std::numeric_limits<NumberType>::quiet_NaN()));
825 
845  std::vector<NumberType>
847  const bool compute_eigenvectors,
848  const std::pair<unsigned int, unsigned int> &index_limits =
849  std::make_pair(numbers::invalid_unsigned_int,
851  const std::pair<NumberType, NumberType> &value_limits =
852  std::make_pair(std::numeric_limits<NumberType>::quiet_NaN(),
853  std::numeric_limits<NumberType>::quiet_NaN()));
854 
855  /*
856  * Stores the distributed matrix in @p filename
857  * using serial routines
858  */
859  void
860  save_serial(const std::string & filename,
861  const std::pair<unsigned int, unsigned int> &chunk_size) const;
862 
863  /*
864  * Loads the distributed matrix from file @p filename
865  * using serial routines
866  */
867  void
868  load_serial(const std::string &filename);
869 
870  /*
871  * Stores the distributed matrix in @p filename
872  * using parallel routines
873  */
874  void
875  save_parallel(const std::string & filename,
876  const std::pair<unsigned int, unsigned int> &chunk_size) const;
877 
878  /*
879  * Loads the distributed matrix from file @p filename
880  * using parallel routines
881  */
882  void
883  load_parallel(const std::string &filename);
884 
890 
896 
902  std::shared_ptr<const Utilities::MPI::ProcessGrid> grid;
903 
907  int n_rows;
908 
913 
918 
923 
928 
933 
937  int descriptor[9];
938 
942  mutable std::vector<NumberType> work;
943 
947  mutable std::vector<int> iwork;
948 
953  std::vector<int> ipiv;
954 
959  const char uplo;
960 
965  const int first_process_row;
966 
972 
977  const int submatrix_row;
978 
983  const int submatrix_column;
984 
989 };
990 
991 // ----------------------- inline functions ----------------------------
992 
993 # ifndef DOXYGEN
994 
995 template <typename NumberType>
996 inline NumberType
997 ScaLAPACKMatrix<NumberType>::local_el(const unsigned int loc_row,
998  const unsigned int loc_column) const
999 {
1000  return (*this)(loc_row, loc_column);
1001 }
1002 
1003 
1004 
1005 template <typename NumberType>
1006 inline NumberType &
1007 ScaLAPACKMatrix<NumberType>::local_el(const unsigned int loc_row,
1008  const unsigned int loc_column)
1009 {
1010  return (*this)(loc_row, loc_column);
1011 }
1012 
1013 
1014 template <typename NumberType>
1015 inline unsigned int
1017 {
1018  return n_rows;
1019 }
1020 
1021 
1022 
1023 template <typename NumberType>
1024 inline unsigned int
1026 {
1027  return n_columns;
1028 }
1029 
1030 
1031 
1032 template <typename NumberType>
1033 unsigned int
1035 {
1036  return n_local_rows;
1037 }
1038 
1039 
1040 
1041 template <typename NumberType>
1042 unsigned int
1044 {
1045  return n_local_columns;
1046 }
1047 
1048 
1049 # endif // DOXYGEN
1050 
1052 
1053 #endif // DEAL_II_WITH_SCALAPACK
1054 
1055 #endif
ScaLAPACKMatrix::reciprocal_condition_number
NumberType reciprocal_condition_number(const NumberType a_norm) const
Definition: scalapack.cc:2312
ScaLAPACKMatrix::eigenpairs_symmetric_by_index
std::vector< NumberType > eigenpairs_symmetric_by_index(const std::pair< unsigned int, unsigned int > &index_limits, const bool compute_eigenvectors)
Definition: scalapack.cc:1429
ScaLAPACKMatrix::compute_lu_factorization
void compute_lu_factorization()
Definition: scalapack.cc:1272
ScaLAPACKMatrix::~ScaLAPACKMatrix
~ScaLAPACKMatrix() override=default
ScaLAPACKMatrix::descriptor
int descriptor[9]
Definition: scalapack.h:937
ScaLAPACKMatrix::scale_rows
void scale_rows(const InputVector &factors)
Definition: scalapack.cc:3506
ScaLAPACKMatrix::l1_norm
NumberType l1_norm() const
Definition: scalapack.cc:2374
ScaLAPACKMatrix::compute_cholesky_factorization
void compute_cholesky_factorization()
Definition: scalapack.cc:1239
lapack_support.h
ScaLAPACKMatrix::operator=
ScaLAPACKMatrix< NumberType > & operator=(const FullMatrix< NumberType > &)
Definition: scalapack.cc:330
CUDAWrappers::block_size
constexpr int block_size
Definition: cuda_size.h:29
ScaLAPACKMatrix::grid
std::shared_ptr< const Utilities::MPI::ProcessGrid > grid
Definition: scalapack.h:902
ScaLAPACKMatrix::state
LAPACKSupport::State state
Definition: scalapack.h:889
ScaLAPACKMatrix::first_process_row
const int first_process_row
Definition: scalapack.h:965
Threads::Mutex
Definition: thread_management.h:91
ScaLAPACKMatrix::copy_from
void copy_from(const LAPACKFullMatrix< NumberType > &matrix, const unsigned int rank)
Definition: scalapack.cc:361
ScaLAPACKMatrix::eigenpairs_symmetric_by_value
std::vector< NumberType > eigenpairs_symmetric_by_value(const std::pair< NumberType, NumberType > &value_limits, const bool compute_eigenvectors)
Definition: scalapack.cc:1452
ScaLAPACKMatrix::norm_general
NumberType norm_general(const char type) const
Definition: scalapack.cc:2416
ScaLAPACKMatrix::ScaLAPACKMatrix
ScaLAPACKMatrix(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)
Definition: scalapack.cc:79
ScaLAPACKMatrix::column_block_size
int column_block_size
Definition: scalapack.h:922
ScaLAPACKMatrix::save_parallel
void save_parallel(const std::string &filename, const std::pair< unsigned int, unsigned int > &chunk_size) const
Definition: scalapack.cc:2806
ScaLAPACKMatrix::mult
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
Definition: scalapack.cc:1066
ScaLAPACKMatrix::submatrix_row
const int submatrix_row
Definition: scalapack.h:977
LAPACKSupport::U
static const char U
Definition: lapack_support.h:167
ScaLAPACKMatrix::save_serial
void save_serial(const std::string &filename, const std::pair< unsigned int, unsigned int > &chunk_size) const
Definition: scalapack.cc:2649
ScaLAPACKMatrix::copy_to
void copy_to(FullMatrix< NumberType > &matrix) const
Definition: scalapack.cc:663
ScaLAPACKMatrix::row_block_size
int row_block_size
Definition: scalapack.h:917
Physics::Elasticity::Kinematics::C
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
ScaLAPACKMatrix::reinit
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)
Definition: scalapack.cc:215
ScaLAPACKMatrix::uplo
const char uplo
Definition: scalapack.h:959
LAPACKFullMatrix
Definition: lapack_full_matrix.h:60
thread_management.h
ScaLAPACKMatrix
Definition: scalapack.h:77
ScaLAPACKMatrix::norm_symmetric
NumberType norm_symmetric(const char type) const
Definition: scalapack.cc:2475
ScaLAPACKMatrix::eigenpairs_symmetric_by_value_MRRR
std::vector< NumberType > eigenpairs_symmetric_by_value_MRRR(const std::pair< NumberType, NumberType > &value_limits, const bool compute_eigenvectors)
Definition: scalapack.cc:1785
ScaLAPACKMatrix::TmTmult
void TmTmult(ScaLAPACKMatrix< NumberType > &C, const ScaLAPACKMatrix< NumberType > &B, const bool adding=false) const
Definition: scalapack.cc:1225
process_grid.h
ScaLAPACKMatrix::copy_transposed
void copy_transposed(const ScaLAPACKMatrix< NumberType > &B)
Definition: scalapack.cc:980
ScaLAPACKMatrix::least_squares
void least_squares(ScaLAPACKMatrix< NumberType > &B, const bool transpose=false)
Definition: scalapack.cc:2133
LAPACKSupport::symmetric
@ symmetric
Matrix is symmetric.
Definition: lapack_support.h:115
ScaLAPACKMatrix::save
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
Definition: scalapack.cc:2607
ScaLAPACKMatrix::set_property
void set_property(const LAPACKSupport::Property property)
Definition: scalapack.cc:302
ScaLAPACKMatrix::mmult
void mmult(ScaLAPACKMatrix< NumberType > &C, const ScaLAPACKMatrix< NumberType > &B, const bool adding=false) const
Definition: scalapack.cc:1183
ScaLAPACKMatrix::property
LAPACKSupport::Property property
Definition: scalapack.h:895
ScaLAPACKMatrix::iwork
std::vector< int > iwork
Definition: scalapack.h:947
ScaLAPACKMatrix::n
size_type n() const
ScaLAPACKMatrix::local_el
NumberType local_el(const unsigned int loc_row, const unsigned int loc_column) const
ScaLAPACKMatrix::n_columns
int n_columns
Definition: scalapack.h:912
LAPACKSupport::State
State
Definition: lapack_support.h:57
mpi.h
ScaLAPACKMatrix::global_row
unsigned int global_row(const unsigned int loc_row) const
Definition: scalapack.cc:495
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
transpose
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
Definition: derivative_form.h:470
Physics::Elasticity::Kinematics::b
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
TransposeTable
Definition: table.h:1931
ScaLAPACKMatrix::load_serial
void load_serial(const std::string &filename)
Definition: scalapack.cc:3063
LAPACKSupport::matrix
@ matrix
Contents is actually a matrix.
Definition: lapack_support.h:60
ScaLAPACKMatrix::mTmult
void mTmult(ScaLAPACKMatrix< NumberType > &C, const ScaLAPACKMatrix< NumberType > &B, const bool adding=false) const
Definition: scalapack.cc:1211
ScaLAPACKMatrix::mutex
Threads::Mutex mutex
Definition: scalapack.h:988
ScaLAPACKMatrix::global_column
unsigned int global_column(const unsigned int loc_column) const
Definition: scalapack.cc:512
LAPACKSupport::general
@ general
No special properties.
Definition: lapack_support.h:113
lapack_full_matrix.h
ScaLAPACKMatrix::eigenpairs_symmetric_MRRR
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()))
Definition: scalapack.cc:1803
ScaLAPACKMatrix::ipiv
std::vector< int > ipiv
Definition: scalapack.h:953
exceptions.h
ScaLAPACKMatrix::get_state
LAPACKSupport::State get_state() const
Definition: scalapack.cc:321
ScaLAPACKMatrix::scale_columns
void scale_columns(const InputVector &factors)
Definition: scalapack.cc:3495
unsigned int
ScaLAPACKMatrix::load
void load(const std::string &filename)
Definition: scalapack.cc:3042
ScaLAPACKMatrix::local_m
unsigned int local_m() const
ScaLAPACKMatrix::invert
void invert()
Definition: scalapack.cc:1313
ScaLAPACKMatrix::Tmmult
void Tmmult(ScaLAPACKMatrix< NumberType > &C, const ScaLAPACKMatrix< NumberType > &B, const bool adding=false) const
Definition: scalapack.cc:1197
ScaLAPACKMatrix::add
void add(const ScaLAPACKMatrix< NumberType > &B, const NumberType a=0., const NumberType b=1., const bool transpose_B=false)
Definition: scalapack.cc:990
ScaLAPACKMatrix::m
size_type m() const
CUDAWrappers::chunk_size
constexpr int chunk_size
Definition: cuda_size.h:35
ScaLAPACKMatrix::pseudoinverse
unsigned int pseudoinverse(const NumberType ratio)
Definition: scalapack.cc:2225
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
ScaLAPACKMatrix::local_n
unsigned int local_n() const
ScaLAPACKMatrix::eigenpairs_symmetric
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()))
Definition: scalapack.cc:1472
ScaLAPACKMatrix::first_process_column
const int first_process_column
Definition: scalapack.h:971
ScaLAPACKMatrix::Tadd
void Tadd(const NumberType b, const ScaLAPACKMatrix< NumberType > &B)
Definition: scalapack.cc:1056
ScaLAPACKMatrix::frobenius_norm
NumberType frobenius_norm() const
Definition: scalapack.cc:2402
config.h
ScaLAPACKMatrix::compute_SVD
std::vector< NumberType > compute_SVD(ScaLAPACKMatrix< NumberType > *U=nullptr, ScaLAPACKMatrix< NumberType > *VT=nullptr)
Definition: scalapack.cc:2013
FullMatrix
Definition: full_matrix.h:71
ScaLAPACKMatrix::work
std::vector< NumberType > work
Definition: scalapack.h:942
ScaLAPACKMatrix::linfty_norm
NumberType linfty_norm() const
Definition: scalapack.cc:2388
LAPACKSupport::Property
Property
Definition: lapack_support.h:110
ScaLAPACKMatrix::eigenpairs_symmetric_by_index_MRRR
std::vector< NumberType > eigenpairs_symmetric_by_index_MRRR(const std::pair< unsigned int, unsigned int > &index_limits, const bool compute_eigenvectors)
Definition: scalapack.cc:1762
TableBase< 2, NumberType >::size
const TableIndices< N > & size() const
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
ScaLAPACKMatrix::n_rows
int n_rows
Definition: scalapack.h:907
ScaLAPACKMatrix::submatrix_column
const int submatrix_column
Definition: scalapack.h:983
ScaLAPACKMatrix::load_parallel
void load_parallel(const std::string &filename)
Definition: scalapack.cc:3241
ScaLAPACKMatrix::get_property
LAPACKSupport::Property get_property() const
Definition: scalapack.cc:312
full_matrix.h
ScaLAPACKMatrix::n_local_rows
int n_local_rows
Definition: scalapack.h:927
ScaLAPACKMatrix::n_local_columns
int n_local_columns
Definition: scalapack.h:932
int