18#ifdef DEAL_II_WITH_SCALAPACK
22# include <deal.II/base/mpi.templates.h>
24# include <deal.II/lac/scalapack.templates.h>
26# 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,
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,
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_;
244 if (
grid->mpi_process_is_active)
249 &(
grid->this_process_row),
251 &(
grid->n_process_rows));
254 &(
grid->this_process_column),
256 &(
grid->n_process_columns));
270 &(
grid->blacs_context),
288template <
typename NumberType>
292 const std::shared_ptr<const Utilities::MPI::ProcessGrid> &process_grid,
301template <
typename NumberType>
306 property = property_;
311template <
typename NumberType>
320template <
typename NumberType>
329template <
typename NumberType>
342 if (
grid->mpi_process_is_active)
350 local_el(i, j) = matrix(glob_i, glob_j);
360template <
typename NumberType>
363 const unsigned int rank)
368 const unsigned int this_mpi_process(
375 "All processes have to call routine with identical rank"));
378 "All processes have to call routine with identical rank"));
382 if (this_mpi_process == rank)
393 MPI_Comm_group(this->
grid->mpi_communicator, &group_A);
395 const std::vector<int> ranks(
n, rank);
397 MPI_Group_incl(group_A,
n, ranks.data(), &group_B);
401 const int ierr = MPI_Comm_create_group(this->
grid->mpi_communicator,
406 int n_proc_rows_B = 1, n_proc_cols_B = 1;
407 int this_process_row_B = -1, this_process_column_B = -1;
408 int blacs_context_B = -1;
409 if (MPI_COMM_NULL != communicator_B)
412 blacs_context_B = Csys2blacs_handle(communicator_B);
413 const char *order =
"Col";
414 Cblacs_gridinit(&blacs_context_B, order, n_proc_rows_B, n_proc_cols_B);
415 Cblacs_gridinfo(blacs_context_B,
419 &this_process_column_B);
425 const bool mpi_process_is_active_B =
426 (this_process_row_B >= 0 && this_process_column_B >= 0);
429 std::vector<int> descriptor_B(9, -1);
430 const int first_process_row_B = 0, first_process_col_B = 0;
432 if (mpi_process_is_active_B)
435 int n_local_rows_B = numroc_(&
n_rows,
438 &first_process_row_B,
442 &this_process_column_B,
443 &first_process_col_B,
448 int lda =
std::max(1, n_local_rows_B);
450 descinit_(descriptor_B.data(),
455 &first_process_row_B,
456 &first_process_col_B,
462 if (this->
grid->mpi_process_is_active)
465 NumberType *loc_vals_A =
466 this->
values.size() > 0 ? this->
values.data() :
nullptr;
467 const NumberType *loc_vals_B =
468 mpi_process_is_active_B ? &(B(0, 0)) :
nullptr;
482 &(this->grid->blacs_context));
484 if (mpi_process_is_active_B)
485 Cblacs_gridexit(blacs_context_B);
487 MPI_Group_free(&group_A);
488 MPI_Group_free(&group_B);
489 if (MPI_COMM_NULL != communicator_B)
497template <
typename NumberType>
503 const int i = loc_row + 1;
506 &(
grid->this_process_row),
508 &(
grid->n_process_rows)) -
514template <
typename NumberType>
521 const int j = loc_column + 1;
524 &(
grid->this_process_column),
526 &(
grid->n_process_columns)) -
532template <
typename NumberType>
535 const unsigned int rank)
const
540 const unsigned int this_mpi_process(
547 "All processes have to call routine with identical rank"));
550 "All processes have to call routine with identical rank"));
553 if (this_mpi_process == rank)
566 MPI_Comm_group(this->
grid->mpi_communicator, &group_A);
568 const std::vector<int> ranks(
n, rank);
570 MPI_Group_incl(group_A,
n, ranks.data(), &group_B);
574 const int ierr = MPI_Comm_create_group(this->
grid->mpi_communicator,
579 int n_proc_rows_B = 1, n_proc_cols_B = 1;
580 int this_process_row_B = -1, this_process_column_B = -1;
581 int blacs_context_B = -1;
582 if (MPI_COMM_NULL != communicator_B)
585 blacs_context_B = Csys2blacs_handle(communicator_B);
586 const char *order =
"Col";
587 Cblacs_gridinit(&blacs_context_B, order, n_proc_rows_B, n_proc_cols_B);
588 Cblacs_gridinfo(blacs_context_B,
592 &this_process_column_B);
598 const bool mpi_process_is_active_B =
599 (this_process_row_B >= 0 && this_process_column_B >= 0);
602 std::vector<int> descriptor_B(9, -1);
603 const int first_process_row_B = 0, first_process_col_B = 0;
605 if (mpi_process_is_active_B)
608 int n_local_rows_B = numroc_(&
n_rows,
611 &first_process_row_B,
615 &this_process_column_B,
616 &first_process_col_B,
621 int lda =
std::max(1, n_local_rows_B);
624 descinit_(descriptor_B.data(),
629 &first_process_row_B,
630 &first_process_col_B,
638 if (this->
grid->mpi_process_is_active)
641 const NumberType *loc_vals_A =
642 this->
values.size() > 0 ? this->
values.data() :
nullptr;
643 NumberType *loc_vals_B = mpi_process_is_active_B ? &(B(0, 0)) :
nullptr;
655 &(this->grid->blacs_context));
657 if (mpi_process_is_active_B)
658 Cblacs_gridexit(blacs_context_B);
660 MPI_Group_free(&group_A);
661 MPI_Group_free(&group_B);
662 if (MPI_COMM_NULL != communicator_B)
668template <
typename NumberType>
681 if (
grid->mpi_process_is_active)
689 matrix(glob_i, glob_j) =
local_el(i, j);
701 for (
unsigned int i = 0; i < matrix.n(); ++i)
702 for (
unsigned int j = i + 1; j < matrix.m(); ++j)
705 for (
unsigned int i = 0; i < matrix.n(); ++i)
706 for (
unsigned int j = 0; j < i; ++j)
713 for (
unsigned int i = 0; i < matrix.n(); ++i)
714 for (
unsigned int j = i + 1; j < matrix.m(); ++j)
715 matrix(i, j) = matrix(j, i);
716 else if (
uplo ==
'U')
717 for (
unsigned int i = 0; i < matrix.n(); ++i)
718 for (
unsigned int j = 0; j < i; ++j)
719 matrix(i, j) = matrix(j, i);
725template <
typename NumberType>
729 const std::pair<unsigned int, unsigned int> &offset_A,
730 const std::pair<unsigned int, unsigned int> &offset_B,
731 const std::pair<unsigned int, unsigned int> &submatrix_size)
const
734 if (submatrix_size.first == 0 || submatrix_size.second == 0)
747 int ierr, comparison;
748 ierr = MPI_Comm_compare(
grid->mpi_communicator,
752 Assert(comparison == MPI_IDENT,
753 ExcMessage(
"Matrix A and B must have a common MPI Communicator"));
761 int union_blacs_context = Csys2blacs_handle(this->
grid->mpi_communicator);
762 const char *order =
"Col";
763 int union_n_process_rows =
765 int union_n_process_columns = 1;
766 Cblacs_gridinit(&union_blacs_context,
768 union_n_process_rows,
769 union_n_process_columns);
771 int n_grid_rows_A, n_grid_columns_A, my_row_A, my_column_A;
772 Cblacs_gridinfo(this->
grid->blacs_context,
779 const bool in_context_A =
780 (my_row_A >= 0 && my_row_A < n_grid_rows_A) &&
781 (my_column_A >= 0 && my_column_A < n_grid_columns_A);
783 int n_grid_rows_B, n_grid_columns_B, my_row_B, my_column_B;
791 const bool in_context_B =
792 (my_row_B >= 0 && my_row_B < n_grid_rows_B) &&
793 (my_column_B >= 0 && my_column_B < n_grid_columns_B);
795 const int n_rows_submatrix = submatrix_size.first;
796 const int n_columns_submatrix = submatrix_size.second;
799 int ia = offset_A.first + 1, ja = offset_A.second + 1;
800 int ib = offset_B.first + 1, jb = offset_B.second + 1;
802 std::array<int, 9> desc_A, desc_B;
804 const NumberType *loc_vals_A =
nullptr;
805 NumberType *loc_vals_B =
nullptr;
814 if (this->
values.size() != 0)
815 loc_vals_A = this->
values.data();
817 for (
unsigned int i = 0; i < desc_A.size(); ++i)
828 for (
unsigned int i = 0; i < desc_B.size(); ++i)
834 pgemr2d(&n_rows_submatrix,
835 &n_columns_submatrix,
844 &union_blacs_context);
849 Cblacs_gridexit(union_blacs_context);
854template <
typename NumberType>
862 if (this->
grid->mpi_process_is_active)
866 "Copying of ScaLAPACK matrices only implemented for dense matrices"));
871 "Copying of ScaLAPACK matrices only implemented for dense matrices"));
887 MPI_Group group_source, group_dest, group_union;
888 ierr = MPI_Comm_group(this->
grid->mpi_communicator, &group_source);
892 ierr = MPI_Group_union(group_source, group_dest, &group_union);
909 ierr = MPI_Comm_create_group(MPI_COMM_WORLD,
912 &mpi_communicator_union);
920 int union_blacs_context = Csys2blacs_handle(mpi_communicator_union);
921 const char *order =
"Col";
922 int union_n_process_rows =
924 int union_n_process_columns = 1;
925 Cblacs_gridinit(&union_blacs_context,
927 union_n_process_rows,
928 union_n_process_columns);
930 const NumberType *loc_vals_source =
nullptr;
931 NumberType *loc_vals_dest =
nullptr;
933 if (this->
grid->mpi_process_is_active && (this->values.size() > 0))
937 "source: process is active but local matrix empty"));
938 loc_vals_source = this->
values.data();
945 "destination: process is active but local matrix empty"));
958 &union_blacs_context);
960 Cblacs_gridexit(union_blacs_context);
962 if (mpi_communicator_union != MPI_COMM_NULL)
964 ierr = MPI_Group_free(&group_source);
966 ierr = MPI_Group_free(&group_dest);
968 ierr = MPI_Group_free(&group_union);
973 if (this->
grid->mpi_process_is_active)
982template <
typename NumberType>
992template <
typename NumberType>
995 const NumberType alpha,
996 const NumberType beta,
997 const bool transpose_B)
1019 ExcMessage(
"The matrices A and B need to have the same process grid"));
1021 if (this->
grid->mpi_process_is_active)
1023 char trans_b = transpose_B ?
'T' :
'N';
1025 (this->
values.size() > 0) ? this->
values.data() :
nullptr;
1026 const NumberType *B_loc =
1048template <
typename NumberType>
1053 add(B, 1, a,
false);
1058template <
typename NumberType>
1068template <
typename NumberType>
1074 const bool transpose_A,
1075 const bool transpose_B)
const
1078 ExcMessage(
"The matrices A and B need to have the same process grid"));
1080 ExcMessage(
"The matrices B and C need to have the same process grid"));
1084 if (!transpose_A && !transpose_B)
1099 else if (transpose_A && !transpose_B)
1114 else if (!transpose_A && transpose_B)
1146 if (this->
grid->mpi_process_is_active)
1148 char trans_a = transpose_A ?
'T' :
'N';
1149 char trans_b = transpose_B ?
'T' :
'N';
1151 const NumberType *A_loc =
1152 (this->
values.size() > 0) ? this->
values.data() :
nullptr;
1153 const NumberType *B_loc =
1155 NumberType *C_loc = (C.values.size() > 0) ? C.values.data() :
nullptr;
1157 int n = C.n_columns;
1177 &C.submatrix_column,
1185template <
typename NumberType>
1189 const bool adding)
const
1192 mult(1., B, 1., C,
false,
false);
1194 mult(1., B, 0, C,
false,
false);
1199template <
typename NumberType>
1203 const bool adding)
const
1206 mult(1., B, 1., C,
true,
false);
1208 mult(1., B, 0, C,
true,
false);
1213template <
typename NumberType>
1217 const bool adding)
const
1220 mult(1., B, 1., C,
false,
true);
1222 mult(1., B, 0, C,
false,
true);
1227template <
typename NumberType>
1231 const bool adding)
const
1234 mult(1., B, 1., C,
true,
true);
1236 mult(1., B, 0, C,
true,
true);
1241template <
typename NumberType>
1248 "Cholesky factorization can be applied to symmetric matrices only."));
1251 "Matrix has to be in Matrix state before calling this function."));
1253 if (
grid->mpi_process_is_active)
1256 NumberType *A_loc = this->
values.data();
1274template <
typename NumberType>
1280 "Matrix has to be in Matrix state before calling this function."));
1282 if (
grid->mpi_process_is_active)
1285 NumberType *A_loc = this->
values.data();
1289 &(
grid->this_process_row),
1291 &(
grid->n_process_rows));
1292 const int mp = numroc_(&
n_rows,
1294 &(
grid->this_process_row),
1296 &(
grid->n_process_rows));
1315template <
typename NumberType>
1333 if (
grid->mpi_process_is_active)
1335 const char uploTriangular =
1337 const char diag =
'N';
1339 NumberType *A_loc = this->
values.data();
1340 ptrtri(&uploTriangular,
1366 if (
grid->mpi_process_is_active)
1369 NumberType *A_loc = this->
values.data();
1386 int lwork = -1, liwork = -1;
1404 lwork =
static_cast<int>(
work[0]);
1407 iwork.resize(liwork);
1431template <
typename NumberType>
1432std::vector<NumberType>
1434 const std::pair<unsigned int, unsigned int> &index_limits,
1435 const bool compute_eigenvectors)
1441 std::pair<unsigned int, unsigned int> idx =
1442 std::make_pair(
std::min(index_limits.first, index_limits.second),
1443 std::max(index_limits.first, index_limits.second));
1446 if (idx.first == 0 && idx.second ==
static_cast<unsigned int>(
n_rows - 1))
1454template <
typename NumberType>
1455std::vector<NumberType>
1457 const std::pair<NumberType, NumberType> &value_limits,
1458 const bool compute_eigenvectors)
1460 Assert(!std::isnan(value_limits.first),
1462 Assert(!std::isnan(value_limits.second),
1465 std::pair<unsigned int, unsigned int> indices =
1474template <
typename NumberType>
1475std::vector<NumberType>
1477 const bool compute_eigenvectors,
1478 const std::pair<unsigned int, unsigned int> &eigenvalue_idx,
1479 const std::pair<NumberType, NumberType> &eigenvalue_limits)
1483 "Matrix has to be in Matrix state before calling this function."));
1485 ExcMessage(
"Matrix has to be symmetric for this operation."));
1487 std::lock_guard<std::mutex> lock(
mutex);
1489 const bool use_values = (std::isnan(eigenvalue_limits.first) ||
1490 std::isnan(eigenvalue_limits.second)) ?
1493 const bool use_indices =
1500 !(use_values && use_indices),
1502 "Prescribing both the index and value range for the eigenvalues is ambiguous"));
1506 std::unique_ptr<ScaLAPACKMatrix<NumberType>>
eigenvectors =
1507 compute_eigenvectors ?
1508 std::make_unique<ScaLAPACKMatrix<NumberType>>(
n_rows,
1511 std::make_unique<ScaLAPACKMatrix<NumberType>>(
1512 grid->n_process_rows,
grid->n_process_columns,
grid, 1, 1);
1519 std::vector<NumberType> ev(
n_rows);
1521 if (
grid->mpi_process_is_active)
1528 char jobz = compute_eigenvectors ?
'V' :
'N';
1531 bool all_eigenpairs =
true;
1532 NumberType vl = NumberType(), vu = NumberType();
1538 NumberType abstol = NumberType();
1545 NumberType orfac = 0;
1547 std::vector<int> ifail;
1554 std::vector<int> iclustr;
1570 all_eigenpairs =
true;
1575 all_eigenpairs =
false;
1576 vl =
std::min(eigenvalue_limits.first, eigenvalue_limits.second);
1577 vu =
std::max(eigenvalue_limits.first, eigenvalue_limits.second);
1583 all_eigenpairs =
false;
1586 il =
std::min(eigenvalue_idx.first, eigenvalue_idx.second) + 1;
1587 iu =
std::max(eigenvalue_idx.first, eigenvalue_idx.second) + 1;
1589 NumberType *A_loc = this->
values.data();
1596 NumberType *eigenvectors_loc =
1597 (compute_eigenvectors ?
eigenvectors->values.data() :
nullptr);
1632 char cmach = compute_eigenvectors ?
'U' :
'S';
1633 plamch(&(this->grid->blacs_context), &cmach, abstol);
1636 iclustr.resize(2 *
grid->n_process_rows *
grid->n_process_columns);
1637 gap.resize(
grid->n_process_rows *
grid->n_process_columns);
1670 lwork =
static_cast<int>(
work[0]);
1697 iwork.resize(liwork);
1734 if (compute_eigenvectors)
1738 while (ev.size() >
static_cast<size_type>(
m))
1744 grid->send_to_inactive(&
m, 1);
1749 if (!
grid->mpi_process_is_active)
1754 grid->send_to_inactive(ev.data(), ev.size());
1761 if (compute_eigenvectors)
1774template <
typename NumberType>
1775std::vector<NumberType>
1777 const std::pair<unsigned int, unsigned int> &index_limits,
1778 const bool compute_eigenvectors)
1784 const std::pair<unsigned int, unsigned int> idx =
1785 std::make_pair(
std::min(index_limits.first, index_limits.second),
1786 std::max(index_limits.first, index_limits.second));
1789 if (idx.first == 0 && idx.second ==
static_cast<unsigned int>(
n_rows - 1))
1797template <
typename NumberType>
1798std::vector<NumberType>
1800 const std::pair<NumberType, NumberType> &value_limits,
1801 const bool compute_eigenvectors)
1806 const std::pair<unsigned int, unsigned int> indices =
1815template <
typename NumberType>
1816std::vector<NumberType>
1818 const bool compute_eigenvectors,
1819 const std::pair<unsigned int, unsigned int> &eigenvalue_idx,
1820 const std::pair<NumberType, NumberType> &eigenvalue_limits)
1824 "Matrix has to be in Matrix state before calling this function."));
1826 ExcMessage(
"Matrix has to be symmetric for this operation."));
1828 std::lock_guard<std::mutex> lock(
mutex);
1830 const bool use_values = (std::isnan(eigenvalue_limits.first) ||
1831 std::isnan(eigenvalue_limits.second)) ?
1834 const bool use_indices =
1841 !(use_values && use_indices),
1843 "Prescribing both the index and value range for the eigenvalues is ambiguous"));
1847 std::unique_ptr<ScaLAPACKMatrix<NumberType>>
eigenvectors =
1848 compute_eigenvectors ?
1849 std::make_unique<ScaLAPACKMatrix<NumberType>>(
n_rows,
1852 std::make_unique<ScaLAPACKMatrix<NumberType>>(
1853 grid->n_process_rows,
grid->n_process_columns,
grid, 1, 1);
1859 std::vector<NumberType> ev(
n_rows);
1866 if (
grid->mpi_process_is_active)
1873 char jobz = compute_eigenvectors ?
'V' :
'N';
1877 NumberType vl = NumberType(), vu = NumberType();
1892 vl =
std::min(eigenvalue_limits.first, eigenvalue_limits.second);
1893 vu =
std::max(eigenvalue_limits.first, eigenvalue_limits.second);
1901 il =
std::min(eigenvalue_idx.first, eigenvalue_idx.second) + 1;
1902 iu =
std::max(eigenvalue_idx.first, eigenvalue_idx.second) + 1;
1904 NumberType *A_loc = this->
values.data();
1912 NumberType *eigenvectors_loc =
1913 (compute_eigenvectors ?
eigenvectors->values.data() :
nullptr);
1944 lwork =
static_cast<int>(
work[0]);
1947 iwork.resize(liwork);
1976 if (compute_eigenvectors)
1980 "psyevr failed to compute all eigenvectors for the selected eigenvalues"));
1985 if (compute_eigenvectors)
1989 while (ev.size() >
static_cast<size_type>(
m))
1995 grid->send_to_inactive(&
m, 1);
2000 if (!
grid->mpi_process_is_active)
2005 grid->send_to_inactive(ev.data(), ev.size());
2012 if (compute_eigenvectors)
2025template <
typename NumberType>
2026std::vector<NumberType>
2032 "Matrix has to be in Matrix state before calling this function."));
2036 const bool left_singular_vectors = (U !=
nullptr) ?
true :
false;
2037 const bool right_singular_vectors = (VT !=
nullptr) ?
true :
false;
2039 if (left_singular_vectors)
2042 Assert(U->n_rows == U->n_columns,
2048 Assert(
grid->blacs_context == U->grid->blacs_context,
2051 if (right_singular_vectors)
2065 std::lock_guard<std::mutex> lock(
mutex);
2069 if (
grid->mpi_process_is_active)
2071 char jobu = left_singular_vectors ?
'V' :
'N';
2072 char jobvt = right_singular_vectors ?
'V' :
'N';
2073 NumberType *A_loc = this->
values.data();
2074 NumberType *U_loc = left_singular_vectors ? U->values.data() :
nullptr;
2075 NumberType *VT_loc = right_singular_vectors ? VT->
values.
data() :
nullptr;
2095 &U->submatrix_column,
2106 lwork =
static_cast<int>(
work[0]);
2120 &U->submatrix_column,
2135 grid->send_to_inactive(sv.data(), sv.size());
2145template <
typename NumberType>
2151 ExcMessage(
"The matrices A and B need to have the same process grid"));
2154 "Matrix has to be in Matrix state before calling this function."));
2157 "Matrix B has to be in Matrix state before calling this function."));
2172 "Use identical block sizes for rows and columns of matrix A"));
2175 "Use identical block sizes for rows and columns of matrix B"));
2178 "Use identical block-cyclic distribution for matrices A and B"));
2180 std::lock_guard<std::mutex> lock(
mutex);
2182 if (
grid->mpi_process_is_active)
2185 NumberType *A_loc = this->
values.data();
2212 lwork =
static_cast<int>(
work[0]);
2237template <
typename NumberType>
2243 "Matrix has to be in Matrix state before calling this function."));
2246 "Use identical block sizes for rows and columns of matrix A"));
2248 ratio > 0. && ratio < 1.,
2250 "input parameter ratio has to be larger than zero and smaller than 1"));
2264 std::vector<NumberType> sv = this->
compute_SVD(&U, &VT);
2265 AssertThrow(sv[0] > std::numeric_limits<NumberType>::min(),
2272 unsigned int n_sv = 1;
2273 std::vector<NumberType> inv_sigma;
2274 inv_sigma.push_back(1 / sv[0]);
2276 for (
unsigned int i = 1; i < sv.size(); ++i)
2277 if (sv[i] > sv[0] * ratio)
2280 inv_sigma.push_back(1 / sv[i]);
2302 std::make_pair(0, 0),
2303 std::make_pair(0, 0),
2304 std::make_pair(
n_rows, n_sv));
2306 std::make_pair(0, 0),
2307 std::make_pair(0, 0),
2317 VT_R.
mult(1, U_R, 0, *
this,
true,
true);
2324template <
typename NumberType>
2327 const NumberType a_norm)
const
2331 "Matrix has to be in Cholesky state before calling this function."));
2332 std::lock_guard<std::mutex> lock(
mutex);
2333 NumberType rcond = 0.;
2335 if (
grid->mpi_process_is_active)
2338 iwork.resize(liwork);
2341 const NumberType *A_loc = this->
values.data();
2361 lwork =
static_cast<int>(std::ceil(
work[0]));
2380 grid->send_to_inactive(&rcond);
2386template <
typename NumberType>
2390 const char type(
'O');
2400template <
typename NumberType>
2404 const char type(
'I');
2414template <
typename NumberType>
2418 const char type(
'F');
2428template <
typename NumberType>
2434 ExcMessage(
"norms can be called in matrix state only."));
2435 std::lock_guard<std::mutex> lock(
mutex);
2436 NumberType res = 0.;
2438 if (
grid->mpi_process_is_active)
2442 &(
grid->this_process_row),
2444 &(
grid->n_process_rows));
2447 &(
grid->this_process_column),
2449 &(
grid->n_process_columns));
2450 const int mp0 = numroc_(&
n_rows,
2452 &(
grid->this_process_row),
2454 &(
grid->n_process_rows));
2457 &(
grid->this_process_column),
2459 &(
grid->n_process_columns));
2465 if (type ==
'O' || type ==
'1')
2467 else if (type ==
'I')
2471 const NumberType *A_loc = this->
values.begin();
2481 grid->send_to_inactive(&res);
2487template <
typename NumberType>
2493 ExcMessage(
"norms can be called in matrix state only."));
2495 ExcMessage(
"Matrix has to be symmetric for this operation."));
2496 std::lock_guard<std::mutex> lock(
mutex);
2497 NumberType res = 0.;
2499 if (
grid->mpi_process_is_active)
2504 ilcm_(&(
grid->n_process_rows), &(
grid->n_process_columns));
2505 const int v2 = lcm / (
grid->n_process_rows);
2509 &(
grid->this_process_row),
2511 &(
grid->n_process_rows));
2514 &(
grid->this_process_column),
2516 &(
grid->n_process_columns));
2519 &(
grid->this_process_row),
2521 &(
grid->n_process_rows));
2524 &(
grid->this_process_column),
2526 &(
grid->n_process_columns));
2534 (type ==
'M' || type ==
'F' || type ==
'E') ? 0 : 2 * Nq0 + Np0 + ldw;
2536 const NumberType *A_loc = this->
values.begin();
2546 grid->send_to_inactive(&res);
2552# ifdef DEAL_II_WITH_HDF5
2558 create_HDF5_state_enum_id(hid_t &state_enum_id)
2564 herr_t status = H5Tenum_insert(state_enum_id,
"cholesky", &val);
2567 status = H5Tenum_insert(state_enum_id,
"eigenvalues", &val);
2570 status = H5Tenum_insert(state_enum_id,
"inverse_matrix", &val);
2573 status = H5Tenum_insert(state_enum_id,
"inverse_svd", &val);
2576 status = H5Tenum_insert(state_enum_id,
"lu", &val);
2579 status = H5Tenum_insert(state_enum_id,
"matrix", &val);
2582 status = H5Tenum_insert(state_enum_id,
"svd", &val);
2585 status = H5Tenum_insert(state_enum_id,
"unusable", &val);
2590 create_HDF5_property_enum_id(hid_t &property_enum_id)
2595 herr_t status = H5Tenum_insert(property_enum_id,
"diagonal", &prop);
2598 status = H5Tenum_insert(property_enum_id,
"general", &prop);
2601 status = H5Tenum_insert(property_enum_id,
"hessenberg", &prop);
2604 status = H5Tenum_insert(property_enum_id,
"lower_triangular", &prop);
2607 status = H5Tenum_insert(property_enum_id,
"symmetric", &prop);
2610 status = H5Tenum_insert(property_enum_id,
"upper_triangular", &prop);
2619template <
typename NumberType>
2622 const std::string &filename,
2623 const std::pair<unsigned int, unsigned int> &chunk_size)
const
2625# ifndef DEAL_II_WITH_HDF5
2631 std::pair<unsigned int, unsigned int> chunks_size_ = chunk_size;
2637 chunks_size_.first =
n_rows;
2638 chunks_size_.second = 1;
2640 Assert(chunks_size_.first > 0,
2641 ExcMessage(
"The row chunk size must be larger than 0."));
2643 Assert(chunks_size_.second > 0,
2644 ExcMessage(
"The column chunk size must be larger than 0."));
2647# ifdef H5_HAVE_PARALLEL
2661template <
typename NumberType>
2664 const std::string &filename,
2665 const std::pair<unsigned int, unsigned int> &chunk_size)
const
2667# ifndef DEAL_II_WITH_HDF5
2682 const auto column_grid =
2683 std::make_shared<Utilities::MPI::ProcessGrid>(this->
grid->mpi_communicator,
2699 H5Fcreate(filename.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
2702 hsize_t chunk_dims[2];
2705 chunk_dims[0] = chunk_size.second;
2706 chunk_dims[1] = chunk_size.first;
2707 hid_t data_property = H5Pcreate(H5P_DATASET_CREATE);
2708 status = H5Pset_chunk(data_property, 2, chunk_dims);
2717 hid_t dataspace_id = H5Screate_simple(2, dims,
nullptr);
2721 hid_t dataset_id = H5Dcreate2(file_id,
2731 dataset_id, type_id, H5S_ALL, H5S_ALL, H5P_DEFAULT, tmp.
values.
data());
2736 hid_t state_enum_id, property_enum_id;
2737 internal::create_HDF5_state_enum_id(state_enum_id);
2738 internal::create_HDF5_property_enum_id(property_enum_id);
2741 hsize_t dims_state[1];
2743 hid_t state_enum_dataspace = H5Screate_simple(1, dims_state,
nullptr);
2745 hid_t state_enum_dataset = H5Dcreate2(file_id,
2748 state_enum_dataspace,
2753 status = H5Dwrite(state_enum_dataset,
2762 hsize_t dims_property[1];
2763 dims_property[0] = 1;
2764 hid_t property_enum_dataspace =
2765 H5Screate_simple(1, dims_property,
nullptr);
2767 hid_t property_enum_dataset = H5Dcreate2(file_id,
2770 property_enum_dataspace,
2775 status = H5Dwrite(property_enum_dataset,
2784 status = H5Dclose(dataset_id);
2786 status = H5Dclose(state_enum_dataset);
2788 status = H5Dclose(property_enum_dataset);
2792 status = H5Sclose(dataspace_id);
2794 status = H5Sclose(state_enum_dataspace);
2796 status = H5Sclose(property_enum_dataspace);
2800 status = H5Tclose(state_enum_id);
2802 status = H5Tclose(property_enum_id);
2806 status = H5Pclose(data_property);
2810 status = H5Fclose(file_id);
2818template <
typename NumberType>
2821 const std::string &filename,
2822 const std::pair<unsigned int, unsigned int> &chunk_size)
const
2824# ifndef DEAL_II_WITH_HDF5
2830 const unsigned int n_mpi_processes(
2832 MPI_Info info = MPI_INFO_NULL;
2840 const auto column_grid =
2841 std::make_shared<Utilities::MPI::ProcessGrid>(this->
grid->mpi_communicator,
2858 const int NB =
std::max(
static_cast<int>(std::ceil(
2859 static_cast<double>(
n_columns) / n_mpi_processes)),
2873 hid_t plist_id = H5Pcreate(H5P_FILE_ACCESS);
2879 H5Fcreate(filename.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, plist_id);
2880 status = H5Pclose(plist_id);
2889 hid_t filespace = H5Screate_simple(2, dims,
nullptr);
2892 hsize_t chunk_dims[2];
2894 chunk_dims[0] = chunk_size.second;
2895 chunk_dims[1] = chunk_size.first;
2896 plist_id = H5Pcreate(H5P_DATASET_CREATE);
2897 H5Pset_chunk(plist_id, 2, chunk_dims);
2899 hid_t dset_id = H5Dcreate2(
2900 file_id,
"/matrix", type_id, filespace, H5P_DEFAULT, plist_id, H5P_DEFAULT);
2902 status = H5Sclose(filespace);
2905 status = H5Pclose(plist_id);
2909 std::vector<int> proc_n_local_rows(n_mpi_processes),
2910 proc_n_local_columns(n_mpi_processes);
2914 proc_n_local_rows.data(),
2922 proc_n_local_columns.data(),
2928 const unsigned int my_rank(
2937 hid_t memspace = H5Screate_simple(2, count,
nullptr);
2939 hsize_t offset[2] = {0};
2940 for (
unsigned int i = 0; i < my_rank; ++i)
2941 offset[0] += proc_n_local_columns[i];
2944 filespace = H5Dget_space(dset_id);
2945 status = H5Sselect_hyperslab(
2946 filespace, H5S_SELECT_SET, offset,
nullptr, count,
nullptr);
2950 plist_id = H5Pcreate(H5P_DATASET_XFER);
2951 status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT);
2957 status = H5Dwrite(dset_id, type_id, memspace, filespace, plist_id, data);
2961 status = H5Dclose(dset_id);
2963 status = H5Sclose(filespace);
2965 status = H5Sclose(memspace);
2967 status = H5Pclose(plist_id);
2969 status = H5Fclose(file_id);
2981 hid_t file_id_reopen =
2982 H5Fopen(filename.c_str(), H5F_ACC_RDWR, H5P_DEFAULT);
2986 hid_t state_enum_id, property_enum_id;
2987 internal::create_HDF5_state_enum_id(state_enum_id);
2988 internal::create_HDF5_property_enum_id(property_enum_id);
2991 hsize_t dims_state[1];
2993 hid_t state_enum_dataspace = H5Screate_simple(1, dims_state,
nullptr);
2995 hid_t state_enum_dataset = H5Dcreate2(file_id_reopen,
2998 state_enum_dataspace,
3003 status = H5Dwrite(state_enum_dataset,
3012 hsize_t dims_property[1];
3013 dims_property[0] = 1;
3014 hid_t property_enum_dataspace =
3015 H5Screate_simple(1, dims_property,
nullptr);
3017 hid_t property_enum_dataset = H5Dcreate2(file_id_reopen,
3020 property_enum_dataspace,
3025 status = H5Dwrite(property_enum_dataset,
3033 status = H5Dclose(state_enum_dataset);
3035 status = H5Dclose(property_enum_dataset);
3037 status = H5Sclose(state_enum_dataspace);
3039 status = H5Sclose(property_enum_dataspace);
3041 status = H5Tclose(state_enum_id);
3043 status = H5Tclose(property_enum_id);
3045 status = H5Fclose(file_id_reopen);
3054template <
typename NumberType>
3058# ifndef DEAL_II_WITH_HDF5
3062# ifdef H5_HAVE_PARALLEL
3075template <
typename NumberType>
3079# ifndef DEAL_II_WITH_HDF5
3090 const auto one_grid =
3091 std::make_shared<Utilities::MPI::ProcessGrid>(this->
grid->mpi_communicator,
3099 int property_int = -1;
3108 hid_t file_id = H5Fopen(filename.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
3111 hid_t dataset_id = H5Dopen2(file_id,
"/matrix", H5P_DEFAULT);
3117 hid_t datatype = H5Dget_type(dataset_id);
3118 H5T_class_t t_class_in = H5Tget_class(datatype);
3121 t_class_in == t_class,
3123 "The data type of the matrix to be read does not match the archive"));
3126 hid_t dataspace_id = H5Dget_space(dataset_id);
3128 const int ndims = H5Sget_simple_extent_ndims(dataspace_id);
3132 H5Sget_simple_extent_dims(dataspace_id, dims,
nullptr);
3136 "The number of columns of the matrix does not match the content of the archive"));
3138 static_cast<int>(dims[1]) ==
n_rows,
3140 "The number of rows of the matrix does not match the content of the archive"));
3143 status = H5Dread(dataset_id,
3153 hid_t state_enum_id, property_enum_id;
3154 internal::create_HDF5_state_enum_id(state_enum_id);
3155 internal::create_HDF5_property_enum_id(property_enum_id);
3158 hid_t dataset_state_id = H5Dopen2(file_id,
"/state", H5P_DEFAULT);
3159 hid_t datatype_state = H5Dget_type(dataset_state_id);
3160 H5T_class_t t_class_state = H5Tget_class(datatype_state);
3163 hid_t dataset_property_id = H5Dopen2(file_id,
"/property", H5P_DEFAULT);
3164 hid_t datatype_property = H5Dget_type(dataset_property_id);
3165 H5T_class_t t_class_property = H5Tget_class(datatype_property);
3169 hid_t dataspace_state = H5Dget_space(dataset_state_id);
3170 hid_t dataspace_property = H5Dget_space(dataset_property_id);
3172 const int ndims_state = H5Sget_simple_extent_ndims(dataspace_state);
3174 const int ndims_property = H5Sget_simple_extent_ndims(dataspace_property);
3177 hsize_t dims_state[1];
3178 H5Sget_simple_extent_dims(dataspace_state, dims_state,
nullptr);
3180 hsize_t dims_property[1];
3181 H5Sget_simple_extent_dims(dataspace_property, dims_property,
nullptr);
3185 status = H5Dread(dataset_state_id,
3195 state_int =
static_cast<int>(tmp.
state);
3197 status = H5Dread(dataset_property_id,
3207 property_int =
static_cast<int>(tmp.
property);
3210 status = H5Sclose(dataspace_id);
3212 status = H5Sclose(dataspace_state);
3214 status = H5Sclose(dataspace_property);
3218 status = H5Tclose(datatype);
3220 status = H5Tclose(state_enum_id);
3222 status = H5Tclose(property_enum_id);
3226 status = H5Dclose(dataset_state_id);
3228 status = H5Dclose(dataset_id);
3230 status = H5Dclose(dataset_property_id);
3234 status = H5Fclose(file_id);
3253template <
typename NumberType>
3257# ifndef DEAL_II_WITH_HDF5
3261# ifndef H5_HAVE_PARALLEL
3266 const unsigned int n_mpi_processes(
3268 MPI_Info info = MPI_INFO_NULL;
3275 const auto column_grid =
3276 std::make_shared<Utilities::MPI::ProcessGrid>(this->
grid->mpi_communicator,
3282 const int NB =
std::max(
static_cast<int>(std::ceil(
3283 static_cast<double>(
n_columns) / n_mpi_processes)),
3294 hid_t plist_id = H5Pcreate(H5P_FILE_ACCESS);
3300 hid_t file_id = H5Fopen(filename.c_str(), H5F_ACC_RDONLY, plist_id);
3301 status = H5Pclose(plist_id);
3305 hid_t dataset_id = H5Dopen2(file_id,
"/matrix", H5P_DEFAULT);
3312 hid_t datatype_inp = H5Dget_type(dataset_id);
3313 H5T_class_t t_class_inp = H5Tget_class(datatype_inp);
3314 H5T_class_t t_class = H5Tget_class(datatype);
3316 t_class_inp == t_class,
3318 "The data type of the matrix to be read does not match the archive"));
3322 hid_t dataspace_id = H5Dget_space(dataset_id);
3324 const int ndims = H5Sget_simple_extent_ndims(dataspace_id);
3328 status = H5Sget_simple_extent_dims(dataspace_id, dims,
nullptr);
3333 "The number of columns of the matrix does not match the content of the archive"));
3335 static_cast<int>(dims[1]) ==
n_rows,
3337 "The number of rows of the matrix does not match the content of the archive"));
3340 std::vector<int> proc_n_local_rows(n_mpi_processes),
3341 proc_n_local_columns(n_mpi_processes);
3345 proc_n_local_rows.data(),
3353 proc_n_local_columns.data(),
3359 const unsigned int my_rank(
3369 hsize_t offset[2] = {0};
3370 for (
unsigned int i = 0; i < my_rank; ++i)
3371 offset[0] += proc_n_local_columns[i];
3374 status = H5Sselect_hyperslab(
3375 dataspace_id, H5S_SELECT_SET, offset,
nullptr, count,
nullptr);
3379 hid_t memspace = H5Screate_simple(2, count,
nullptr);
3383 H5Dread(dataset_id, datatype, memspace, dataspace_id, H5P_DEFAULT, data);
3387 hid_t state_enum_id, property_enum_id;
3388 internal::create_HDF5_state_enum_id(state_enum_id);
3389 internal::create_HDF5_property_enum_id(property_enum_id);
3392 hid_t dataset_state_id = H5Dopen2(file_id,
"/state", H5P_DEFAULT);
3393 hid_t datatype_state = H5Dget_type(dataset_state_id);
3394 H5T_class_t t_class_state = H5Tget_class(datatype_state);
3397 hid_t dataset_property_id = H5Dopen2(file_id,
"/property", H5P_DEFAULT);
3398 hid_t datatype_property = H5Dget_type(dataset_property_id);
3399 H5T_class_t t_class_property = H5Tget_class(datatype_property);
3403 hid_t dataspace_state = H5Dget_space(dataset_state_id);
3404 hid_t dataspace_property = H5Dget_space(dataset_property_id);
3406 const int ndims_state = H5Sget_simple_extent_ndims(dataspace_state);
3408 const int ndims_property = H5Sget_simple_extent_ndims(dataspace_property);
3411 hsize_t dims_state[1];
3412 H5Sget_simple_extent_dims(dataspace_state, dims_state,
nullptr);
3414 hsize_t dims_property[1];
3415 H5Sget_simple_extent_dims(dataspace_property, dims_property,
nullptr);
3420 dataset_state_id, state_enum_id, H5S_ALL, H5S_ALL, H5P_DEFAULT, &tmp.
state);
3423 status = H5Dread(dataset_property_id,
3432 status = H5Sclose(memspace);
3434 status = H5Dclose(dataset_id);
3436 status = H5Dclose(dataset_state_id);
3438 status = H5Dclose(dataset_property_id);
3440 status = H5Sclose(dataspace_id);
3442 status = H5Sclose(dataspace_state);
3444 status = H5Sclose(dataspace_property);
3448 status = H5Tclose(state_enum_id);
3450 status = H5Tclose(property_enum_id);
3452 status = H5Fclose(file_id);
3468 template <
typename NumberType>
3476 for (
unsigned int i = 0; i < matrix.local_n(); ++i)
3478 const NumberType s = factors[matrix.global_column(i)];
3480 for (
unsigned int j = 0; j < matrix.local_m(); ++j)
3481 matrix.local_el(j, i) *= s;
3485 template <
typename NumberType>
3487 scale_rows(ScaLAPACKMatrix<NumberType> &matrix,
3488 const ArrayView<const NumberType> &factors)
3493 for (
unsigned int i = 0; i <
matrix.local_m(); ++i)
3497 for (
unsigned int j = 0; j <
matrix.local_n(); ++j)
3498 matrix.local_el(i, j) *= s;
3507template <
typename NumberType>
3508template <
class InputVector>
3512 if (this->
grid->mpi_process_is_active)
3518template <
typename NumberType>
3519template <
class InputVector>
3523 if (this->
grid->mpi_process_is_active)
3530# include "lac/scalapack.inst"
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, 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)
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
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
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
void compute_cholesky_factorization()
void copy_transposed(const ScaLAPACKMatrix< NumberType > &B)
NumberType linfty_norm() const
void scale_columns(const InputVector &factors)
AlignedVector< NumberType > values
size_type size(const unsigned int i) const
void reinit(const size_type size1, const size_type size2, const bool omit_default_initialization=false)
MPI_Comm mpi_communicator
bool mpi_process_is_active
void send_to_inactive(NumberType *value, const int count=1) const
const unsigned int this_mpi_process
#define DEAL_II_NAMESPACE_OPEN
constexpr bool running_in_debug_mode()
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcIO()
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 & ExcNeedsHDF5()
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.
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
T max(const T &t, const MPI_Comm mpi_communicator)
T min(const T &t, const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
const MPI_Datatype mpi_type_id_for_type
void free_communicator(MPI_Comm mpi_communicator)
constexpr 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 > >, dim > eigenvectors(const SymmetricTensor< 2, dim, Number > &T, const SymmetricTensorEigenvectorMethod method=SymmetricTensorEigenvectorMethod::ql_implicit_shifts)