19#ifdef DEAL_II_WITH_SCALAPACK
23# include <deal.II/base/mpi.templates.h>
25# include <deal.II/lac/scalapack.templates.h>
27# ifdef DEAL_II_WITH_HDF5
35# ifdef DEAL_II_WITH_HDF5
37template <
typename number>
49 return H5T_NATIVE_DOUBLE;
55 return H5T_NATIVE_FLOAT;
61 return H5T_NATIVE_INT;
67 return H5T_NATIVE_UINT;
73 return H5T_NATIVE_CHAR;
79template <
typename NumberType>
83 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
90 , first_process_column(0)
104template <
typename NumberType>
107 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
120template <
typename NumberType>
122 const std::string & filename,
123 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
129 , first_process_column(0)
131 , submatrix_column(1)
133# ifndef DEAL_II_WITH_HDF5
141 "This function is only available when deal.II is configured with HDF5"));
144 const unsigned int this_mpi_process(
150 if (this_mpi_process == 0)
155 hid_t file = H5Fopen(filename.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
159 hid_t dataset = H5Dopen2(file,
"/matrix", H5P_DEFAULT);
163 hid_t filespace = H5Dget_space(dataset);
166 int rank = H5Sget_simple_extent_ndims(filespace);
169 status = H5Sget_simple_extent_dims(filespace, dims,
nullptr);
178 status = H5Sclose(filespace);
180 status = H5Dclose(dataset);
182 status = H5Fclose(file);
185 int ierr = MPI_Bcast(&
n_rows,
189 process_grid->mpi_communicator);
196 process_grid->mpi_communicator);
207 load(filename.c_str());
214template <
typename NumberType>
219 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
224 Assert(row_block_size_ > 0,
ExcMessage(
"Row block size has to be positive."));
225 Assert(column_block_size_ > 0,
226 ExcMessage(
"Column block size has to be positive."));
228 row_block_size_ <= n_rows_,
230 "Row block size can not be greater than the number of rows of the matrix"));
232 column_block_size_ <= n_columns_,
234 "Column block size can not be greater than the number of columns of the matrix"));
237 property = property_;
240 n_columns = n_columns_;
241 row_block_size = row_block_size_;
242 column_block_size = column_block_size_;
244 if (grid->mpi_process_is_active)
247 n_local_rows = numroc_(&n_rows,
249 &(grid->this_process_row),
251 &(grid->n_process_rows));
252 n_local_columns = numroc_(&n_columns,
254 &(grid->this_process_column),
255 &first_process_column,
256 &(grid->n_process_columns));
260 int lda =
std::max(1, n_local_rows);
263 descinit_(descriptor,
269 &first_process_column,
270 &(grid->blacs_context),
281 n_local_columns = -1;
282 std::fill(std::begin(descriptor), std::end(descriptor), -1);
288template <
typename NumberType>
292 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
296 reinit(size, size, process_grid, block_size, block_size, property);
301template <
typename NumberType>
306 property = property_;
311template <
typename NumberType>
320template <
typename NumberType>
329template <
typename NumberType>
339 Assert(n_columns ==
int(matrix.n()),
342 if (grid->mpi_process_is_active)
344 for (
int i = 0; i < n_local_rows; ++i)
346 const int glob_i = global_row(i);
347 for (
int j = 0; j < n_local_columns; ++j)
349 const int glob_j = global_column(j);
350 local_el(i, j) = matrix(glob_i, glob_j);
360template <
typename NumberType>
363 const unsigned int rank)
365 if (n_rows * n_columns == 0)
368 const unsigned int this_mpi_process(
373 ExcMessage(
"All processes have to call routine with identical rank"));
375 ExcMessage(
"All processes have to call routine with identical rank"));
379 if (this_mpi_process == rank)
390 MPI_Comm_group(this->grid->mpi_communicator, &group_A);
392 const std::vector<int> ranks(n, rank);
394 MPI_Group_incl(group_A, n, ranks.data(), &group_B);
398 const int ierr = MPI_Comm_create_group(this->grid->mpi_communicator,
403 int n_proc_rows_B = 1, n_proc_cols_B = 1;
404 int this_process_row_B = -1, this_process_column_B = -1;
405 int blacs_context_B = -1;
406 if (MPI_COMM_NULL != communicator_B)
409 blacs_context_B = Csys2blacs_handle(communicator_B);
410 const char *order =
"Col";
411 Cblacs_gridinit(&blacs_context_B, order, n_proc_rows_B, n_proc_cols_B);
412 Cblacs_gridinfo(blacs_context_B,
416 &this_process_column_B);
422 const bool mpi_process_is_active_B =
423 (this_process_row_B >= 0 && this_process_column_B >= 0);
426 std::vector<int> descriptor_B(9, -1);
427 const int first_process_row_B = 0, first_process_col_B = 0;
429 if (mpi_process_is_active_B)
432 int n_local_rows_B = numroc_(&n_rows,
435 &first_process_row_B,
437 int n_local_cols_B = numroc_(&n_columns,
439 &this_process_column_B,
440 &first_process_col_B,
444 (void)n_local_cols_B;
446 int lda =
std::max(1, n_local_rows_B);
448 descinit_(descriptor_B.data(),
453 &first_process_row_B,
454 &first_process_col_B,
460 if (this->grid->mpi_process_is_active)
463 NumberType *loc_vals_A =
464 this->values.size() > 0 ? this->values.data() :
nullptr;
465 const NumberType *loc_vals_B =
466 mpi_process_is_active_B ? &(B(0, 0)) :
nullptr;
480 &(this->grid->blacs_context));
482 if (mpi_process_is_active_B)
483 Cblacs_gridexit(blacs_context_B);
485 MPI_Group_free(&group_A);
486 MPI_Group_free(&group_B);
487 if (MPI_COMM_NULL != communicator_B)
495template <
typename NumberType>
499 Assert(n_local_rows >= 0 && loc_row <
static_cast<unsigned int>(n_local_rows),
501 const int i = loc_row + 1;
504 &(grid->this_process_row),
506 &(grid->n_process_rows)) -
512template <
typename NumberType>
516 Assert(n_local_columns >= 0 &&
517 loc_column <
static_cast<unsigned int>(n_local_columns),
519 const int j = loc_column + 1;
522 &(grid->this_process_column),
523 &first_process_column,
524 &(grid->n_process_columns)) -
530template <
typename NumberType>
533 const unsigned int rank)
const
535 if (n_rows * n_columns == 0)
538 const unsigned int this_mpi_process(
543 ExcMessage(
"All processes have to call routine with identical rank"));
545 ExcMessage(
"All processes have to call routine with identical rank"));
548 if (this_mpi_process == rank)
561 MPI_Comm_group(this->grid->mpi_communicator, &group_A);
563 const std::vector<int> ranks(n, rank);
565 MPI_Group_incl(group_A, n, ranks.data(), &group_B);
569 const int ierr = MPI_Comm_create_group(this->grid->mpi_communicator,
574 int n_proc_rows_B = 1, n_proc_cols_B = 1;
575 int this_process_row_B = -1, this_process_column_B = -1;
576 int blacs_context_B = -1;
577 if (MPI_COMM_NULL != communicator_B)
580 blacs_context_B = Csys2blacs_handle(communicator_B);
581 const char *order =
"Col";
582 Cblacs_gridinit(&blacs_context_B, order, n_proc_rows_B, n_proc_cols_B);
583 Cblacs_gridinfo(blacs_context_B,
587 &this_process_column_B);
593 const bool mpi_process_is_active_B =
594 (this_process_row_B >= 0 && this_process_column_B >= 0);
597 std::vector<int> descriptor_B(9, -1);
598 const int first_process_row_B = 0, first_process_col_B = 0;
600 if (mpi_process_is_active_B)
603 int n_local_rows_B = numroc_(&n_rows,
606 &first_process_row_B,
608 int n_local_cols_B = numroc_(&n_columns,
610 &this_process_column_B,
611 &first_process_col_B,
615 (void)n_local_cols_B;
617 int lda =
std::max(1, n_local_rows_B);
620 descinit_(descriptor_B.data(),
625 &first_process_row_B,
626 &first_process_col_B,
634 if (this->grid->mpi_process_is_active)
637 const NumberType *loc_vals_A =
638 this->values.size() > 0 ? this->values.data() :
nullptr;
639 NumberType *loc_vals_B = mpi_process_is_active_B ? &(B(0, 0)) :
nullptr;
651 &(this->grid->blacs_context));
653 if (mpi_process_is_active_B)
654 Cblacs_gridexit(blacs_context_B);
656 MPI_Group_free(&group_A);
657 MPI_Group_free(&group_B);
658 if (MPI_COMM_NULL != communicator_B)
664template <
typename NumberType>
673 Assert(n_columns ==
int(matrix.n()),
677 if (grid->mpi_process_is_active)
679 for (
int i = 0; i < n_local_rows; ++i)
681 const int glob_i = global_row(i);
682 for (
int j = 0; j < n_local_columns; ++j)
684 const int glob_j = global_column(j);
685 matrix(glob_i, glob_j) = local_el(i, j);
697 for (
unsigned int i = 0; i < matrix.n(); ++i)
698 for (
unsigned int j = i + 1; j < matrix.m(); ++j)
701 for (
unsigned int i = 0; i < matrix.n(); ++i)
702 for (
unsigned int j = 0; j < i; ++j)
709 for (
unsigned int i = 0; i < matrix.n(); ++i)
710 for (
unsigned int j = i + 1; j < matrix.m(); ++j)
711 matrix(i, j) = matrix(j, i);
712 else if (uplo ==
'U')
713 for (
unsigned int i = 0; i < matrix.n(); ++i)
714 for (
unsigned int j = 0; j < i; ++j)
715 matrix(i, j) = matrix(j, i);
721template <
typename NumberType>
725 const std::pair<unsigned int, unsigned int> &offset_A,
726 const std::pair<unsigned int, unsigned int> &offset_B,
727 const std::pair<unsigned int, unsigned int> &submatrix_size)
const
730 if (submatrix_size.first == 0 || submatrix_size.second == 0)
743 int ierr, comparison;
744 ierr = MPI_Comm_compare(grid->mpi_communicator,
745 B.
grid->mpi_communicator,
748 Assert(comparison == MPI_IDENT,
749 ExcMessage(
"Matrix A and B must have a common MPI Communicator"));
757 int union_blacs_context = Csys2blacs_handle(this->grid->mpi_communicator);
758 const char *order =
"Col";
759 int union_n_process_rows =
761 int union_n_process_columns = 1;
762 Cblacs_gridinit(&union_blacs_context,
764 union_n_process_rows,
765 union_n_process_columns);
767 int n_grid_rows_A, n_grid_columns_A, my_row_A, my_column_A;
768 Cblacs_gridinfo(this->grid->blacs_context,
775 const bool in_context_A =
776 (my_row_A >= 0 && my_row_A < n_grid_rows_A) &&
777 (my_column_A >= 0 && my_column_A < n_grid_columns_A);
779 int n_grid_rows_B, n_grid_columns_B, my_row_B, my_column_B;
780 Cblacs_gridinfo(B.
grid->blacs_context,
787 const bool in_context_B =
788 (my_row_B >= 0 && my_row_B < n_grid_rows_B) &&
789 (my_column_B >= 0 && my_column_B < n_grid_columns_B);
791 const int n_rows_submatrix = submatrix_size.first;
792 const int n_columns_submatrix = submatrix_size.second;
795 int ia = offset_A.first + 1, ja = offset_A.second + 1;
796 int ib = offset_B.first + 1, jb = offset_B.second + 1;
798 std::array<int, 9> desc_A, desc_B;
800 const NumberType *loc_vals_A =
nullptr;
801 NumberType * loc_vals_B =
nullptr;
810 if (this->values.size() != 0)
811 loc_vals_A = this->values.data();
813 for (
unsigned int i = 0; i < desc_A.size(); ++i)
814 desc_A[i] = this->descriptor[i];
822 loc_vals_B = B.
values.data();
824 for (
unsigned int i = 0; i < desc_B.size(); ++i)
830 pgemr2d(&n_rows_submatrix,
831 &n_columns_submatrix,
840 &union_blacs_context);
845 Cblacs_gridexit(union_blacs_context);
850template <
typename NumberType>
858 if (this->grid->mpi_process_is_active)
860 this->descriptor[0] == 1,
862 "Copying of ScaLAPACK matrices only implemented for dense matrices"));
863 if (dest.
grid->mpi_process_is_active)
867 "Copying of ScaLAPACK matrices only implemented for dense matrices"));
883 MPI_Group group_source, group_dest, group_union;
884 ierr = MPI_Comm_group(this->grid->mpi_communicator, &group_source);
886 ierr = MPI_Comm_group(dest.
grid->mpi_communicator, &group_dest);
888 ierr = MPI_Group_union(group_source, group_dest, &group_union);
905 ierr = MPI_Comm_create_group(MPI_COMM_WORLD,
908 &mpi_communicator_union);
916 int union_blacs_context = Csys2blacs_handle(mpi_communicator_union);
917 const char *order =
"Col";
918 int union_n_process_rows =
920 int union_n_process_columns = 1;
921 Cblacs_gridinit(&union_blacs_context,
923 union_n_process_rows,
924 union_n_process_columns);
926 const NumberType *loc_vals_source =
nullptr;
927 NumberType * loc_vals_dest =
nullptr;
929 if (this->grid->mpi_process_is_active && (this->values.size() > 0))
933 "source: process is active but local matrix empty"));
934 loc_vals_source = this->values.data();
936 if (dest.
grid->mpi_process_is_active && (dest.
values.size() > 0))
941 "destination: process is active but local matrix empty"));
942 loc_vals_dest = dest.
values.data();
954 &union_blacs_context);
956 Cblacs_gridexit(union_blacs_context);
958 if (mpi_communicator_union != MPI_COMM_NULL)
960 ierr = MPI_Group_free(&group_source);
962 ierr = MPI_Group_free(&group_dest);
964 ierr = MPI_Group_free(&group_union);
969 if (this->grid->mpi_process_is_active)
970 dest.
values = this->values;
978template <
typename NumberType>
988template <
typename NumberType>
991 const NumberType alpha,
992 const NumberType beta,
993 const bool transpose_B)
1015 ExcMessage(
"The matrices A and B need to have the same process grid"));
1017 if (this->grid->mpi_process_is_active)
1019 char trans_b = transpose_B ?
'T' :
'N';
1021 (this->values.size() > 0) ? this->values.data() :
nullptr;
1022 const NumberType *B_loc =
1044template <
typename NumberType>
1049 add(B, 1, a,
false);
1054template <
typename NumberType>
1064template <
typename NumberType>
1070 const bool transpose_A,
1071 const bool transpose_B)
const
1074 ExcMessage(
"The matrices A and B need to have the same process grid"));
1076 ExcMessage(
"The matrices B and C need to have the same process grid"));
1080 if (!transpose_A && !transpose_B)
1084 Assert(this->n_rows == C.n_rows,
1088 Assert(this->row_block_size == C.row_block_size,
1095 else if (transpose_A && !transpose_B)
1099 Assert(this->n_columns == C.n_rows,
1103 Assert(this->column_block_size == C.row_block_size,
1110 else if (!transpose_A && transpose_B)
1114 Assert(this->n_rows == C.n_rows,
1118 Assert(this->row_block_size == C.row_block_size,
1130 Assert(this->n_columns == C.n_rows,
1134 Assert(this->column_block_size == C.row_block_size,
1142 if (this->grid->mpi_process_is_active)
1144 char trans_a = transpose_A ?
'T' :
'N';
1145 char trans_b = transpose_B ?
'T' :
'N';
1147 const NumberType *A_loc =
1148 (this->values.size() > 0) ? this->values.data() :
nullptr;
1149 const NumberType *B_loc =
1151 NumberType *C_loc = (C.values.size() > 0) ? C.values.data() :
nullptr;
1153 int n = C.n_columns;
1154 int k = transpose_A ? this->n_rows : this->n_columns;
1163 &(this->submatrix_row),
1164 &(this->submatrix_column),
1173 &C.submatrix_column,
1181template <
typename NumberType>
1185 const bool adding)
const
1188 mult(1., B, 1., C,
false,
false);
1190 mult(1., B, 0, C,
false,
false);
1195template <
typename NumberType>
1199 const bool adding)
const
1202 mult(1., B, 1., C,
true,
false);
1204 mult(1., B, 0, C,
true,
false);
1209template <
typename NumberType>
1213 const bool adding)
const
1216 mult(1., B, 1., C,
false,
true);
1218 mult(1., B, 0, C,
false,
true);
1223template <
typename NumberType>
1227 const bool adding)
const
1230 mult(1., B, 1., C,
true,
true);
1232 mult(1., B, 0, C,
true,
true);
1237template <
typename NumberType>
1244 "Cholesky factorization can be applied to symmetric matrices only."));
1247 "Matrix has to be in Matrix state before calling this function."));
1249 if (grid->mpi_process_is_active)
1252 NumberType *A_loc = this->values.data();
1270template <
typename NumberType>
1276 "Matrix has to be in Matrix state before calling this function."));
1278 if (grid->mpi_process_is_active)
1281 NumberType *A_loc = this->values.data();
1283 const int iarow = indxg2p_(&submatrix_row,
1285 &(grid->this_process_row),
1287 &(grid->n_process_rows));
1288 const int mp = numroc_(&n_rows,
1290 &(grid->this_process_row),
1292 &(grid->n_process_rows));
1293 ipiv.resize(mp + row_block_size);
1311template <
typename NumberType>
1329 if (grid->mpi_process_is_active)
1331 const char uploTriangular =
1333 const char diag =
'N';
1335 NumberType *A_loc = this->values.data();
1336 ptrtri(&uploTriangular,
1358 compute_cholesky_factorization();
1360 compute_lu_factorization();
1362 if (grid->mpi_process_is_active)
1365 NumberType *A_loc = this->values.data();
1382 int lwork = -1, liwork = -1;
1400 lwork =
static_cast<int>(work[0]);
1403 iwork.resize(liwork);
1427template <
typename NumberType>
1428std::vector<NumberType>
1430 const std::pair<unsigned int, unsigned int> &index_limits,
1431 const bool compute_eigenvectors)
1437 std::pair<unsigned int, unsigned int> idx =
1438 std::make_pair(
std::min(index_limits.first, index_limits.second),
1439 std::max(index_limits.first, index_limits.second));
1442 if (idx.first == 0 && idx.second ==
static_cast<unsigned int>(n_rows - 1))
1443 return eigenpairs_symmetric(compute_eigenvectors);
1445 return eigenpairs_symmetric(compute_eigenvectors, idx);
1450template <
typename NumberType>
1451std::vector<NumberType>
1453 const std::pair<NumberType, NumberType> &value_limits,
1454 const bool compute_eigenvectors)
1456 Assert(!std::isnan(value_limits.first),
1458 Assert(!std::isnan(value_limits.second),
1461 std::pair<unsigned int, unsigned int> indices =
1465 return eigenpairs_symmetric(compute_eigenvectors, indices, value_limits);
1470template <
typename NumberType>
1471std::vector<NumberType>
1473 const bool compute_eigenvectors,
1474 const std::pair<unsigned int, unsigned int> &eigenvalue_idx,
1475 const std::pair<NumberType, NumberType> & eigenvalue_limits)
1479 "Matrix has to be in Matrix state before calling this function."));
1481 ExcMessage(
"Matrix has to be symmetric for this operation."));
1483 std::lock_guard<std::mutex> lock(mutex);
1485 const bool use_values = (std::isnan(eigenvalue_limits.first) ||
1486 std::isnan(eigenvalue_limits.second)) ?
1489 const bool use_indices =
1496 !(use_values && use_indices),
1498 "Prescribing both the index and value range for the eigenvalues is ambiguous"));
1502 std::unique_ptr<ScaLAPACKMatrix<NumberType>>
eigenvectors =
1503 compute_eigenvectors ?
1504 std::make_unique<ScaLAPACKMatrix<NumberType>>(n_rows,
1508 grid->n_process_rows, grid->n_process_columns, grid, 1, 1);
1515 std::vector<NumberType> ev(n_rows);
1517 if (grid->mpi_process_is_active)
1524 char jobz = compute_eigenvectors ?
'V' :
'N';
1527 bool all_eigenpairs =
true;
1528 NumberType vl = NumberType(), vu = NumberType();
1534 NumberType abstol = NumberType();
1541 NumberType orfac = 0;
1543 std::vector<int> ifail;
1550 std::vector<int> iclustr;
1556 std::vector<NumberType> gap(n_local_rows * n_local_columns);
1566 all_eigenpairs =
true;
1571 all_eigenpairs =
false;
1572 vl =
std::min(eigenvalue_limits.first, eigenvalue_limits.second);
1573 vu =
std::max(eigenvalue_limits.first, eigenvalue_limits.second);
1579 all_eigenpairs =
false;
1582 il =
std::min(eigenvalue_idx.first, eigenvalue_idx.second) + 1;
1583 iu =
std::max(eigenvalue_idx.first, eigenvalue_idx.second) + 1;
1585 NumberType *A_loc = this->values.data();
1592 NumberType *eigenvectors_loc =
1593 (compute_eigenvectors ?
eigenvectors->values.data() :
nullptr);
1618 char cmach = compute_eigenvectors ?
'U' :
'S';
1619 plamch(&(this->grid->blacs_context), &cmach, abstol);
1621 ifail.resize(n_rows);
1622 iclustr.resize(2 * grid->n_process_rows * grid->n_process_columns);
1623 gap.resize(grid->n_process_rows * grid->n_process_columns);
1656 lwork =
static_cast<int>(work[0]);
1683 iwork.resize(liwork);
1720 if (compute_eigenvectors)
1724 while (ev.size() >
static_cast<size_type>(m))
1730 grid->send_to_inactive(&m, 1);
1735 if (!grid->mpi_process_is_active)
1740 grid->send_to_inactive(ev.data(), ev.size());
1747 if (compute_eigenvectors)
1760template <
typename NumberType>
1761std::vector<NumberType>
1763 const std::pair<unsigned int, unsigned int> &index_limits,
1764 const bool compute_eigenvectors)
1770 const std::pair<unsigned int, unsigned int> idx =
1771 std::make_pair(
std::min(index_limits.first, index_limits.second),
1772 std::max(index_limits.first, index_limits.second));
1775 if (idx.first == 0 && idx.second ==
static_cast<unsigned int>(n_rows - 1))
1776 return eigenpairs_symmetric_MRRR(compute_eigenvectors);
1778 return eigenpairs_symmetric_MRRR(compute_eigenvectors, idx);
1783template <
typename NumberType>
1784std::vector<NumberType>
1786 const std::pair<NumberType, NumberType> &value_limits,
1787 const bool compute_eigenvectors)
1792 const std::pair<unsigned int, unsigned int> indices =
1796 return eigenpairs_symmetric_MRRR(compute_eigenvectors, indices, value_limits);
1801template <
typename NumberType>
1802std::vector<NumberType>
1804 const bool compute_eigenvectors,
1805 const std::pair<unsigned int, unsigned int> &eigenvalue_idx,
1806 const std::pair<NumberType, NumberType> & eigenvalue_limits)
1810 "Matrix has to be in Matrix state before calling this function."));
1812 ExcMessage(
"Matrix has to be symmetric for this operation."));
1814 std::lock_guard<std::mutex> lock(mutex);
1816 const bool use_values = (std::isnan(eigenvalue_limits.first) ||
1817 std::isnan(eigenvalue_limits.second)) ?
1820 const bool use_indices =
1827 !(use_values && use_indices),
1829 "Prescribing both the index and value range for the eigenvalues is ambiguous"));
1833 std::unique_ptr<ScaLAPACKMatrix<NumberType>>
eigenvectors =
1834 compute_eigenvectors ?
1835 std::make_unique<ScaLAPACKMatrix<NumberType>>(n_rows,
1839 grid->n_process_rows, grid->n_process_columns, grid, 1, 1);
1845 std::vector<NumberType> ev(n_rows);
1852 if (grid->mpi_process_is_active)
1859 char jobz = compute_eigenvectors ?
'V' :
'N';
1863 NumberType vl = NumberType(), vu = NumberType();
1878 vl =
std::min(eigenvalue_limits.first, eigenvalue_limits.second);
1879 vu =
std::max(eigenvalue_limits.first, eigenvalue_limits.second);
1887 il =
std::min(eigenvalue_idx.first, eigenvalue_idx.second) + 1;
1888 iu =
std::max(eigenvalue_idx.first, eigenvalue_idx.second) + 1;
1890 NumberType *A_loc = this->values.data();
1898 NumberType *eigenvectors_loc =
1899 (compute_eigenvectors ?
eigenvectors->values.data() :
nullptr);
1930 lwork =
static_cast<int>(work[0]);
1933 iwork.resize(liwork);
1962 if (compute_eigenvectors)
1966 "psyevr failed to compute all eigenvectors for the selected eigenvalues"));
1971 if (compute_eigenvectors)
1975 while (ev.size() >
static_cast<size_type>(m))
1981 grid->send_to_inactive(&m, 1);
1986 if (!grid->mpi_process_is_active)
1991 grid->send_to_inactive(ev.data(), ev.size());
1998 if (compute_eigenvectors)
2011template <
typename NumberType>
2012std::vector<NumberType>
2018 "Matrix has to be in Matrix state before calling this function."));
2019 Assert(row_block_size == column_block_size,
2022 const bool left_singluar_vectors = (U !=
nullptr) ?
true :
false;
2023 const bool right_singluar_vectors = (VT !=
nullptr) ?
true :
false;
2025 if (left_singluar_vectors)
2034 Assert(grid->blacs_context == U->
grid->blacs_context,
2037 if (right_singluar_vectors)
2047 Assert(grid->blacs_context == VT->
grid->blacs_context,
2049 VT->
grid->blacs_context));
2051 std::lock_guard<std::mutex> lock(mutex);
2053 std::vector<NumberType> sv(
std::min(n_rows, n_columns));
2055 if (grid->mpi_process_is_active)
2057 char jobu = left_singluar_vectors ?
'V' :
'N';
2058 char jobvt = right_singluar_vectors ?
'V' :
'N';
2059 NumberType *A_loc = this->values.data();
2060 NumberType *U_loc = left_singluar_vectors ? U->
values.data() :
nullptr;
2061 NumberType *VT_loc = right_singluar_vectors ? VT->
values.data() :
nullptr;
2092 lwork =
static_cast<int>(work[0]);
2121 grid->send_to_inactive(sv.data(), sv.size());
2131template <
typename NumberType>
2137 ExcMessage(
"The matrices A and B need to have the same process grid"));
2140 "Matrix has to be in Matrix state before calling this function."));
2143 "Matrix B has to be in Matrix state before calling this function."));
2156 Assert(row_block_size == column_block_size,
2158 "Use identical block sizes for rows and columns of matrix A"));
2161 "Use identical block sizes for rows and columns of matrix B"));
2164 "Use identical block-cyclic distribution for matrices A and B"));
2166 std::lock_guard<std::mutex> lock(mutex);
2168 if (grid->mpi_process_is_active)
2171 NumberType *A_loc = this->values.data();
2172 NumberType *B_loc = B.
values.data();
2198 lwork =
static_cast<int>(work[0]);
2223template <
typename NumberType>
2229 "Matrix has to be in Matrix state before calling this function."));
2230 Assert(row_block_size == column_block_size,
2232 "Use identical block sizes for rows and columns of matrix A"));
2234 ratio > 0. && ratio < 1.,
2236 "input parameter ratio has to be larger than zero and smaller than 1"));
2250 std::vector<NumberType> sv = this->compute_SVD(&U, &VT);
2251 AssertThrow(sv[0] > std::numeric_limits<NumberType>::min(),
2258 unsigned int n_sv = 1;
2259 std::vector<NumberType> inv_sigma;
2260 inv_sigma.push_back(1 / sv[0]);
2262 for (
unsigned int i = 1; i < sv.size(); ++i)
2263 if (sv[i] > sv[0] * ratio)
2266 inv_sigma.push_back(1 / sv[i]);
2288 std::make_pair(0, 0),
2289 std::make_pair(0, 0),
2290 std::make_pair(n_rows, n_sv));
2292 std::make_pair(0, 0),
2293 std::make_pair(0, 0),
2294 std::make_pair(n_sv, n_columns));
2297 this->reinit(n_columns,
2303 VT_R.
mult(1, U_R, 0, *
this,
true,
true);
2310template <
typename NumberType>
2313 const NumberType a_norm)
const
2317 "Matrix has to be in Cholesky state before calling this function."));
2318 std::lock_guard<std::mutex> lock(mutex);
2319 NumberType rcond = 0.;
2321 if (grid->mpi_process_is_active)
2323 int liwork = n_local_rows;
2324 iwork.resize(liwork);
2327 const NumberType *A_loc = this->values.data();
2347 lwork =
static_cast<int>(std::ceil(work[0]));
2366 grid->send_to_inactive(&rcond);
2372template <
typename NumberType>
2376 const char type(
'O');
2379 return norm_symmetric(type);
2381 return norm_general(type);
2386template <
typename NumberType>
2390 const char type(
'I');
2393 return norm_symmetric(type);
2395 return norm_general(type);
2400template <
typename NumberType>
2404 const char type(
'F');
2407 return norm_symmetric(type);
2409 return norm_general(type);
2414template <
typename NumberType>
2420 ExcMessage(
"norms can be called in matrix state only."));
2421 std::lock_guard<std::mutex> lock(mutex);
2422 NumberType res = 0.;
2424 if (grid->mpi_process_is_active)
2426 const int iarow = indxg2p_(&submatrix_row,
2428 &(grid->this_process_row),
2430 &(grid->n_process_rows));
2431 const int iacol = indxg2p_(&submatrix_column,
2433 &(grid->this_process_column),
2434 &first_process_column,
2435 &(grid->n_process_columns));
2436 const int mp0 = numroc_(&n_rows,
2438 &(grid->this_process_row),
2440 &(grid->n_process_rows));
2441 const int nq0 = numroc_(&n_columns,
2443 &(grid->this_process_column),
2445 &(grid->n_process_columns));
2451 if (type ==
'O' || type ==
'1')
2453 else if (type ==
'I')
2457 const NumberType *A_loc = this->values.begin();
2467 grid->send_to_inactive(&res);
2473template <
typename NumberType>
2479 ExcMessage(
"norms can be called in matrix state only."));
2481 ExcMessage(
"Matrix has to be symmetric for this operation."));
2482 std::lock_guard<std::mutex> lock(mutex);
2483 NumberType res = 0.;
2485 if (grid->mpi_process_is_active)
2490 ilcm_(&(grid->n_process_rows), &(grid->n_process_columns));
2491 const int v2 = lcm / (grid->n_process_rows);
2493 const int IAROW = indxg2p_(&submatrix_row,
2495 &(grid->this_process_row),
2497 &(grid->n_process_rows));
2498 const int IACOL = indxg2p_(&submatrix_column,
2500 &(grid->this_process_column),
2501 &first_process_column,
2502 &(grid->n_process_columns));
2503 const int Np0 = numroc_(&n_columns ,
2505 &(grid->this_process_row),
2507 &(grid->n_process_rows));
2508 const int Nq0 = numroc_(&n_columns ,
2510 &(grid->this_process_column),
2512 &(grid->n_process_columns));
2514 const int v1 = iceil_(&Np0, &row_block_size);
2515 const int ldw = (n_local_rows == n_local_columns) ?
2517 row_block_size * iceil_(&
v1, &v2);
2520 (type ==
'M' || type ==
'F' || type ==
'E') ? 0 : 2 * Nq0 + Np0 + ldw;
2522 const NumberType *A_loc = this->values.begin();
2532 grid->send_to_inactive(&res);
2538# ifdef DEAL_II_WITH_HDF5
2544 create_HDF5_state_enum_id(hid_t &state_enum_id)
2550 herr_t status = H5Tenum_insert(state_enum_id,
"cholesky", &val);
2553 status = H5Tenum_insert(state_enum_id,
"eigenvalues", &val);
2556 status = H5Tenum_insert(state_enum_id,
"inverse_matrix", &val);
2559 status = H5Tenum_insert(state_enum_id,
"inverse_svd", &val);
2562 status = H5Tenum_insert(state_enum_id,
"lu", &val);
2565 status = H5Tenum_insert(state_enum_id,
"matrix", &val);
2568 status = H5Tenum_insert(state_enum_id,
"svd", &val);
2571 status = H5Tenum_insert(state_enum_id,
"unusable", &val);
2576 create_HDF5_property_enum_id(hid_t &property_enum_id)
2581 herr_t status = H5Tenum_insert(property_enum_id,
"diagonal", &prop);
2584 status = H5Tenum_insert(property_enum_id,
"general", &prop);
2587 status = H5Tenum_insert(property_enum_id,
"hessenberg", &prop);
2590 status = H5Tenum_insert(property_enum_id,
"lower_triangular", &prop);
2593 status = H5Tenum_insert(property_enum_id,
"symmetric", &prop);
2596 status = H5Tenum_insert(property_enum_id,
"upper_triangular", &prop);
2605template <
typename NumberType>
2608 const std::string & filename,
2609 const std::pair<unsigned int, unsigned int> &chunk_size)
const
2611# ifndef DEAL_II_WITH_HDF5
2617 std::pair<unsigned int, unsigned int> chunks_size_ = chunk_size;
2623 chunks_size_.first = n_rows;
2624 chunks_size_.second = 1;
2626 Assert(chunks_size_.first > 0,
2627 ExcMessage(
"The row chunk size must be larger than 0."));
2629 Assert(chunks_size_.second > 0,
2630 ExcMessage(
"The column chunk size must be larger than 0."));
2633# ifdef H5_HAVE_PARALLEL
2635 save_parallel(filename, chunks_size_);
2639 save_serial(filename, chunks_size_);
2647template <
typename NumberType>
2650 const std::string & filename,
2651 const std::pair<unsigned int, unsigned int> &chunk_size)
const
2653# ifndef DEAL_II_WITH_HDF5
2668 const auto column_grid =
2669 std::make_shared<Utilities::MPI::ProcessGrid>(this->grid->mpi_communicator,
2673 const int MB = n_rows, NB = n_columns;
2679 if (tmp.
grid->mpi_process_is_active)
2685 H5Fcreate(filename.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
2688 hsize_t chunk_dims[2];
2691 chunk_dims[0] = chunk_size.second;
2692 chunk_dims[1] = chunk_size.first;
2693 hid_t data_property = H5Pcreate(H5P_DATASET_CREATE);
2694 status = H5Pset_chunk(data_property, 2, chunk_dims);
2701 dims[0] = n_columns;
2703 hid_t dataspace_id = H5Screate_simple(2, dims,
nullptr);
2707 hid_t dataset_id = H5Dcreate2(file_id,
2717 dataset_id, type_id, H5S_ALL, H5S_ALL, H5P_DEFAULT, tmp.
values.data());
2722 hid_t state_enum_id, property_enum_id;
2723 internal::create_HDF5_state_enum_id(state_enum_id);
2724 internal::create_HDF5_property_enum_id(property_enum_id);
2727 hsize_t dims_state[1];
2729 hid_t state_enum_dataspace = H5Screate_simple(1, dims_state,
nullptr);
2731 hid_t state_enum_dataset = H5Dcreate2(file_id,
2734 state_enum_dataspace,
2739 status = H5Dwrite(state_enum_dataset,
2748 hsize_t dims_property[1];
2749 dims_property[0] = 1;
2750 hid_t property_enum_dataspace =
2751 H5Screate_simple(1, dims_property,
nullptr);
2753 hid_t property_enum_dataset = H5Dcreate2(file_id,
2756 property_enum_dataspace,
2761 status = H5Dwrite(property_enum_dataset,
2770 status = H5Dclose(dataset_id);
2772 status = H5Dclose(state_enum_dataset);
2774 status = H5Dclose(property_enum_dataset);
2778 status = H5Sclose(dataspace_id);
2780 status = H5Sclose(state_enum_dataspace);
2782 status = H5Sclose(property_enum_dataspace);
2786 status = H5Tclose(state_enum_id);
2788 status = H5Tclose(property_enum_id);
2792 status = H5Pclose(data_property);
2796 status = H5Fclose(file_id);
2804template <
typename NumberType>
2807 const std::string & filename,
2808 const std::pair<unsigned int, unsigned int> &chunk_size)
const
2810# ifndef DEAL_II_WITH_HDF5
2816 const unsigned int n_mpi_processes(
2818 MPI_Info info = MPI_INFO_NULL;
2826 const auto column_grid =
2827 std::make_shared<Utilities::MPI::ProcessGrid>(this->grid->mpi_communicator,
2831 const int MB = n_rows;
2844 const int NB =
std::max(
static_cast<int>(std::ceil(
2845 static_cast<double>(n_columns) / n_mpi_processes)),
2852 NumberType *data = (tmp.
values.size() > 0) ? tmp.
values.data() :
nullptr;
2859 hid_t plist_id = H5Pcreate(H5P_FILE_ACCESS);
2860 status = H5Pset_fapl_mpio(plist_id, tmp.
grid->mpi_communicator, info);
2865 H5Fcreate(filename.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, plist_id);
2866 status = H5Pclose(plist_id);
2875 hid_t filespace = H5Screate_simple(2, dims,
nullptr);
2878 hsize_t chunk_dims[2];
2880 chunk_dims[0] = chunk_size.second;
2881 chunk_dims[1] = chunk_size.first;
2882 plist_id = H5Pcreate(H5P_DATASET_CREATE);
2883 H5Pset_chunk(plist_id, 2, chunk_dims);
2885 hid_t dset_id = H5Dcreate2(
2886 file_id,
"/matrix", type_id, filespace, H5P_DEFAULT, plist_id, H5P_DEFAULT);
2888 status = H5Sclose(filespace);
2891 status = H5Pclose(plist_id);
2895 std::vector<int> proc_n_local_rows(n_mpi_processes),
2896 proc_n_local_columns(n_mpi_processes);
2900 proc_n_local_rows.data(),
2903 tmp.
grid->mpi_communicator);
2908 proc_n_local_columns.data(),
2911 tmp.
grid->mpi_communicator);
2914 const unsigned int my_rank(
2923 hid_t memspace = H5Screate_simple(2, count,
nullptr);
2925 hsize_t offset[2] = {0};
2926 for (
unsigned int i = 0; i < my_rank; ++i)
2927 offset[0] += proc_n_local_columns[i];
2930 filespace = H5Dget_space(dset_id);
2931 status = H5Sselect_hyperslab(
2932 filespace, H5S_SELECT_SET, offset,
nullptr, count,
nullptr);
2936 plist_id = H5Pcreate(H5P_DATASET_XFER);
2937 status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT);
2941 if (tmp.
values.size() > 0)
2943 status = H5Dwrite(dset_id, type_id, memspace, filespace, plist_id, data);
2947 status = H5Dclose(dset_id);
2949 status = H5Sclose(filespace);
2951 status = H5Sclose(memspace);
2953 status = H5Pclose(plist_id);
2955 status = H5Fclose(file_id);
2960 ierr = MPI_Barrier(tmp.
grid->mpi_communicator);
2964 if (tmp.
grid->this_mpi_process == 0)
2967 hid_t file_id_reopen =
2968 H5Fopen(filename.c_str(), H5F_ACC_RDWR, H5P_DEFAULT);
2972 hid_t state_enum_id, property_enum_id;
2973 internal::create_HDF5_state_enum_id(state_enum_id);
2974 internal::create_HDF5_property_enum_id(property_enum_id);
2977 hsize_t dims_state[1];
2979 hid_t state_enum_dataspace = H5Screate_simple(1, dims_state,
nullptr);
2981 hid_t state_enum_dataset = H5Dcreate2(file_id_reopen,
2984 state_enum_dataspace,
2989 status = H5Dwrite(state_enum_dataset,
2998 hsize_t dims_property[1];
2999 dims_property[0] = 1;
3000 hid_t property_enum_dataspace =
3001 H5Screate_simple(1, dims_property,
nullptr);
3003 hid_t property_enum_dataset = H5Dcreate2(file_id_reopen,
3006 property_enum_dataspace,
3011 status = H5Dwrite(property_enum_dataset,
3019 status = H5Dclose(state_enum_dataset);
3021 status = H5Dclose(property_enum_dataset);
3023 status = H5Sclose(state_enum_dataspace);
3025 status = H5Sclose(property_enum_dataspace);
3027 status = H5Tclose(state_enum_id);
3029 status = H5Tclose(property_enum_id);
3031 status = H5Fclose(file_id_reopen);
3040template <
typename NumberType>
3044# ifndef DEAL_II_WITH_HDF5
3048# ifdef H5_HAVE_PARALLEL
3050 load_parallel(filename);
3054 load_serial(filename);
3061template <
typename NumberType>
3065# ifndef DEAL_II_WITH_HDF5
3076 const auto one_grid =
3077 std::make_shared<Utilities::MPI::ProcessGrid>(this->grid->mpi_communicator,
3081 const int MB = n_rows, NB = n_columns;
3085 int property_int = -1;
3089 if (tmp.
grid->mpi_process_is_active)
3094 hid_t file_id = H5Fopen(filename.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
3097 hid_t dataset_id = H5Dopen2(file_id,
"/matrix", H5P_DEFAULT);
3103 hid_t datatype = H5Dget_type(dataset_id);
3104 H5T_class_t t_class_in = H5Tget_class(datatype);
3107 t_class_in == t_class,
3109 "The data type of the matrix to be read does not match the archive"));
3112 hid_t dataspace_id = H5Dget_space(dataset_id);
3114 const int ndims = H5Sget_simple_extent_ndims(dataspace_id);
3118 H5Sget_simple_extent_dims(dataspace_id, dims,
nullptr);
3120 static_cast<int>(dims[0]) == n_columns,
3122 "The number of columns of the matrix does not match the content of the archive"));
3124 static_cast<int>(dims[1]) == n_rows,
3126 "The number of rows of the matrix does not match the content of the archive"));
3129 status = H5Dread(dataset_id,
3139 hid_t state_enum_id, property_enum_id;
3140 internal::create_HDF5_state_enum_id(state_enum_id);
3141 internal::create_HDF5_property_enum_id(property_enum_id);
3144 hid_t dataset_state_id = H5Dopen2(file_id,
"/state", H5P_DEFAULT);
3145 hid_t datatype_state = H5Dget_type(dataset_state_id);
3146 H5T_class_t t_class_state = H5Tget_class(datatype_state);
3149 hid_t dataset_property_id = H5Dopen2(file_id,
"/property", H5P_DEFAULT);
3150 hid_t datatype_property = H5Dget_type(dataset_property_id);
3151 H5T_class_t t_class_property = H5Tget_class(datatype_property);
3155 hid_t dataspace_state = H5Dget_space(dataset_state_id);
3156 hid_t dataspace_property = H5Dget_space(dataset_property_id);
3158 const int ndims_state = H5Sget_simple_extent_ndims(dataspace_state);
3160 const int ndims_property = H5Sget_simple_extent_ndims(dataspace_property);
3163 hsize_t dims_state[1];
3164 H5Sget_simple_extent_dims(dataspace_state, dims_state,
nullptr);
3166 hsize_t dims_property[1];
3167 H5Sget_simple_extent_dims(dataspace_property, dims_property,
nullptr);
3171 status = H5Dread(dataset_state_id,
3181 state_int =
static_cast<int>(tmp.
state);
3183 status = H5Dread(dataset_property_id,
3193 property_int =
static_cast<int>(tmp.
property);
3196 status = H5Sclose(dataspace_id);
3198 status = H5Sclose(dataspace_state);
3200 status = H5Sclose(dataspace_property);
3204 status = H5Tclose(datatype);
3206 status = H5Tclose(state_enum_id);
3208 status = H5Tclose(property_enum_id);
3212 status = H5Dclose(dataset_state_id);
3214 status = H5Dclose(dataset_id);
3216 status = H5Dclose(dataset_property_id);
3220 status = H5Fclose(file_id);
3224 tmp.
grid->send_to_inactive(&state_int, 1);
3227 tmp.
grid->send_to_inactive(&property_int, 1);
3239template <
typename NumberType>
3243# ifndef DEAL_II_WITH_HDF5
3247# ifndef H5_HAVE_PARALLEL
3251 const unsigned int n_mpi_processes(
3253 MPI_Info info = MPI_INFO_NULL;
3260 const auto column_grid =
3261 std::make_shared<Utilities::MPI::ProcessGrid>(this->grid->mpi_communicator,
3265 const int MB = n_rows;
3267 const int NB =
std::max(
static_cast<int>(std::ceil(
3268 static_cast<double>(n_columns) / n_mpi_processes)),
3274 NumberType *data = (tmp.
values.size() > 0) ? tmp.
values.data() :
nullptr;
3279 hid_t plist_id = H5Pcreate(H5P_FILE_ACCESS);
3280 status = H5Pset_fapl_mpio(plist_id, tmp.
grid->mpi_communicator, info);
3285 hid_t file_id = H5Fopen(filename.c_str(), H5F_ACC_RDONLY, plist_id);
3286 status = H5Pclose(plist_id);
3290 hid_t dataset_id = H5Dopen2(file_id,
"/matrix", H5P_DEFAULT);
3297 hid_t datatype_inp = H5Dget_type(dataset_id);
3298 H5T_class_t t_class_inp = H5Tget_class(datatype_inp);
3299 H5T_class_t t_class = H5Tget_class(datatype);
3301 t_class_inp == t_class,
3303 "The data type of the matrix to be read does not match the archive"));
3307 hid_t dataspace_id = H5Dget_space(dataset_id);
3309 const int ndims = H5Sget_simple_extent_ndims(dataspace_id);
3313 status = H5Sget_simple_extent_dims(dataspace_id, dims,
nullptr);
3316 static_cast<int>(dims[0]) == n_columns,
3318 "The number of columns of the matrix does not match the content of the archive"));
3320 static_cast<int>(dims[1]) == n_rows,
3322 "The number of rows of the matrix does not match the content of the archive"));
3325 std::vector<int> proc_n_local_rows(n_mpi_processes),
3326 proc_n_local_columns(n_mpi_processes);
3330 proc_n_local_rows.data(),
3333 tmp.
grid->mpi_communicator);
3338 proc_n_local_columns.data(),
3341 tmp.
grid->mpi_communicator);
3344 const unsigned int my_rank(
3354 hsize_t offset[2] = {0};
3355 for (
unsigned int i = 0; i < my_rank; ++i)
3356 offset[0] += proc_n_local_columns[i];
3359 status = H5Sselect_hyperslab(
3360 dataspace_id, H5S_SELECT_SET, offset,
nullptr, count,
nullptr);
3364 hid_t memspace = H5Screate_simple(2, count,
nullptr);
3368 H5Dread(dataset_id, datatype, memspace, dataspace_id, H5P_DEFAULT, data);
3372 hid_t state_enum_id, property_enum_id;
3373 internal::create_HDF5_state_enum_id(state_enum_id);
3374 internal::create_HDF5_property_enum_id(property_enum_id);
3377 hid_t dataset_state_id = H5Dopen2(file_id,
"/state", H5P_DEFAULT);
3378 hid_t datatype_state = H5Dget_type(dataset_state_id);
3379 H5T_class_t t_class_state = H5Tget_class(datatype_state);
3382 hid_t dataset_property_id = H5Dopen2(file_id,
"/property", H5P_DEFAULT);
3383 hid_t datatype_property = H5Dget_type(dataset_property_id);
3384 H5T_class_t t_class_property = H5Tget_class(datatype_property);
3388 hid_t dataspace_state = H5Dget_space(dataset_state_id);
3389 hid_t dataspace_property = H5Dget_space(dataset_property_id);
3391 const int ndims_state = H5Sget_simple_extent_ndims(dataspace_state);
3393 const int ndims_property = H5Sget_simple_extent_ndims(dataspace_property);
3396 hsize_t dims_state[1];
3397 H5Sget_simple_extent_dims(dataspace_state, dims_state,
nullptr);
3399 hsize_t dims_property[1];
3400 H5Sget_simple_extent_dims(dataspace_property, dims_property,
nullptr);
3405 dataset_state_id, state_enum_id, H5S_ALL, H5S_ALL, H5P_DEFAULT, &tmp.
state);
3408 status = H5Dread(dataset_property_id,
3417 status = H5Sclose(memspace);
3419 status = H5Dclose(dataset_id);
3421 status = H5Dclose(dataset_state_id);
3423 status = H5Dclose(dataset_property_id);
3425 status = H5Sclose(dataspace_id);
3427 status = H5Sclose(dataspace_state);
3429 status = H5Sclose(dataspace_property);
3433 status = H5Tclose(state_enum_id);
3435 status = H5Tclose(property_enum_id);
3437 status = H5Fclose(file_id);
3453 template <
typename NumberType>
3461 for (
unsigned int i = 0; i < matrix.local_n(); ++i)
3463 const NumberType s = factors[matrix.global_column(i)];
3465 for (
unsigned int j = 0; j < matrix.local_m(); ++j)
3466 matrix.local_el(j, i) *= s;
3470 template <
typename NumberType>
3478 for (
unsigned int i = 0; i <
matrix.local_m(); ++i)
3480 const NumberType s = factors[
matrix.global_row(i)];
3482 for (
unsigned int j = 0; j <
matrix.local_n(); ++j)
3483 matrix.local_el(i, j) *= s;
3492template <
typename NumberType>
3493template <
class InputVector>
3497 if (this->grid->mpi_process_is_active)
3503template <
typename NumberType>
3504template <
class InputVector>
3508 if (this->grid->mpi_process_is_active)
3515# include "scalapack.inst"
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
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)
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)
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)
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)
std::shared_ptr< const Utilities::MPI::ProcessGrid > grid
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)
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
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
void compute_cholesky_factorization()
void copy_transposed(const ScaLAPACKMatrix< NumberType > &B)
NumberType linfty_norm() const
void scale_columns(const InputVector &factors)
AlignedVector< T > values
void reinit(const size_type size1, const size_type size2, const bool omit_default_initialization=false)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcErrorCode(std::string arg1, types::blas_int arg2)
#define Assert(cond, exc)
#define AssertIsFinite(number)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ cholesky
Contents is a Cholesky decomposition.
@ lu
Contents is an LU decomposition.
@ matrix
Contents is actually a matrix.
@ unusable
Contents is something useless.
@ inverse_matrix
Contents is the inverse of a matrix.
@ svd
Matrix contains singular value decomposition,.
@ inverse_svd
Matrix is the inverse of a singular value decomposition.
@ eigenvalues
Eigenvalue vector is filled.
@ symmetric
Matrix is symmetric.
@ hessenberg
Matrix is in upper Hessenberg form.
@ diagonal
Matrix is diagonal.
@ upper_triangular
Matrix is upper triangular.
@ lower_triangular
Matrix is lower triangular.
@ general
No special properties.
void free_communicator(MPI_Comm &mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
T min(const T &t, const MPI_Comm &mpi_communicator)
T sum(const T &t, const MPI_Comm &mpi_communicator)
const MPI_Datatype mpi_type_id_for_type
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
T max(const T &t, const MPI_Comm &mpi_communicator)
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
hid_t hdf5_type_id(const number *)
std::array< std::pair< Number, Tensor< 1, dim, Number > >, std::integral_constant< int, dim >::value > eigenvectors(const SymmetricTensor< 2, dim, Number > &T, const SymmetricTensorEigenvectorMethod method=SymmetricTensorEigenvectorMethod::ql_implicit_shifts)