20 #include "libmesh/dense_matrix.h"
21 #include "libmesh/dense_vector.h"
22 #include "libmesh/int_range.h"
24 #if (LIBMESH_HAVE_PETSC)
25 # include "libmesh/petsc_macro.h"
26 # include "libmesh/ignore_warnings.h"
27 # include <petscblaslapack.h>
28 # include "libmesh/restore_warnings.h"
34 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
48 result_size = other.
m() * this->n();
49 if (other.
n() == this->m())
52 libmesh_fallthrough();
55 result_size = other.
n() * this->m();
56 if (other.
m() == this->n())
59 libmesh_fallthrough();
60 case LEFT_MULTIPLY_TRANSPOSE:
62 result_size = other.
n() * this->n();
63 if (other.
m() == this->m())
66 libmesh_fallthrough();
67 case RIGHT_MULTIPLY_TRANSPOSE:
69 result_size = other.
m() * this->m();
70 if (other.
n() == this->n())
73 libmesh_fallthrough();
75 libmesh_error_msg(
"Unknown flag selected or matrices are incompatible for multiplication.");
79 const DenseMatrix<T> * const_that = cast_ptr<const DenseMatrix<T> *>(&other);
96 if (flag==RIGHT_MULTIPLY || flag==RIGHT_MULTIPLY_TRANSPOSE)
113 PetscBLASInt M = static_cast<PetscBLASInt>(
A->n() );
122 PetscBLASInt N = static_cast<PetscBLASInt>(
B->m() );
132 PetscBLASInt K = static_cast<PetscBLASInt>(
A->m() );
136 PetscBLASInt LDA = static_cast<PetscBLASInt>(
A->n() );
140 PetscBLASInt LDB = static_cast<PetscBLASInt>(
B->n() );
142 if (flag == LEFT_MULTIPLY_TRANSPOSE)
145 N = static_cast<PetscBLASInt>(
B->n() );
148 else if (flag == RIGHT_MULTIPLY_TRANSPOSE)
156 PetscBLASInt LDC = M;
167 std::vector<T> result (result_size);
170 BLASgemm_(transa, transb, &M, &N, &K,
pPS(&alpha),
171 pPS(
A->get_values().data()), &LDA,
172 pPS(
B->get_values().data()), &LDB,
pPS(&beta),
173 pPS(result.data()), &LDC);
178 case LEFT_MULTIPLY: { this->_m = other.
m();
break; }
179 case RIGHT_MULTIPLY: { this->_n = other.
n();
break; }
180 case LEFT_MULTIPLY_TRANSPOSE: { this->_m = other.
n();
break; }
181 case RIGHT_MULTIPLY_TRANSPOSE: { this->_n = other.
m();
break; }
183 libmesh_error_msg(
"Unknown flag selected.");
187 this->_val.swap(result);
196 libmesh_error_msg(
"No PETSc-provided BLAS/LAPACK available!");
207 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
214 libmesh_assert_equal_to (this->_decomposition_type, NONE);
222 PetscBLASInt M = this->n();
227 PetscBLASInt N = this->m();
236 PetscBLASInt LDA = M;
242 this->_pivots.resize( std::min(M,N) );
251 PetscBLASInt INFO = 0;
254 LAPACKgetrf_(&M, &N,
pPS(this->_val.data()), &LDA, _pivots.data(),
259 libmesh_error_msg(
"INFO=" << INFO <<
", Error during Lapack LU factorization!");
262 this->_decomposition_type = LU_BLAS_LAPACK;
270 libmesh_error_msg(
"No PETSc-provided BLAS/LAPACK available!");
306 std::vector<Real> sigma_val;
307 std::vector<Number> U_val;
308 std::vector<Number> VT_val;
310 _svd_helper(JOBU, JOBVT, sigma_val, U_val, VT_val);
313 sigma.
resize(cast_int<unsigned int>(sigma_val.size()));
315 sigma(i) = sigma_val[i];
355 std::vector<Real> sigma_val;
358 int min_MN = (M < N) ? M : N;
377 sigma.
resize(cast_int<unsigned int>(sigma_val.size()));
379 sigma(i) = sigma_val[i];
382 #if (LIBMESH_HAVE_PETSC)
387 std::vector<Real> & sigma_val,
388 std::vector<Number> & U_val,
389 std::vector<Number> & VT_val)
395 PetscBLASInt M = this->n();
400 PetscBLASInt N = this->m();
402 PetscBLASInt min_MN = (M < N) ? M : N;
403 PetscBLASInt max_MN = (M > N) ? M : N;
418 PetscBLASInt LDA = M;
422 sigma_val.resize( min_MN );
427 PetscBLASInt LDU = M;
436 U_val.resize( LDU*min_MN );
438 U_val.resize( LDU*M );
443 PetscBLASInt LDVT = N;
453 VT_val.resize( LDVT*N );
464 PetscBLASInt larger = (3*min_MN+max_MN > 5*min_MN) ? 3*min_MN+max_MN : 5*min_MN;
465 PetscBLASInt LWORK = (larger > 1) ? larger : 1;
475 std::vector<Number> WORK( LWORK );
484 PetscBLASInt INFO = 0;
487 #ifdef LIBMESH_USE_REAL_NUMBERS
489 LAPACKgesvd_(&JOBU, &JOBVT, &M, &N,
pPS(_val.data()), &LDA,
490 pPS(sigma_val.data()),
pPS(U_val.data()), &LDU,
491 pPS(VT_val.data()), &LDVT,
pPS(WORK.data()), &LWORK,
497 std::vector<Number> val_copy(_val.size());
499 val_copy[i] = _val[i];
501 std::vector<Real> RWORK(5 * min_MN);
502 LAPACKgesvd_(&JOBU, &JOBVT, &M, &N, val_copy.data(), &LDA, sigma_val.data(), U_val.data(),
503 &LDU, VT_val.data(), &LDVT, WORK.data(), &LWORK, RWORK.data(), &INFO);
508 libmesh_error_msg(
"INFO=" << INFO <<
", Error during Lapack SVD calculation!");
518 std::vector<Number> &,
519 std::vector<Number> &)
521 libmesh_error_msg(
"No PETSc-provided BLAS/LAPACK available!");
528 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
529 #if !PETSC_VERSION_LESS_THAN(3,1,0)
545 this->get_transpose(A_trans);
550 PetscBLASInt M = A_trans.
n();
555 PetscBLASInt N = A_trans.
m();
558 PetscBLASInt max_MN = std::max(M,N);
559 PetscBLASInt min_MN = std::min(M,N);
566 PetscBLASInt NRHS = 1;
574 std::vector<T> & A_trans_vals = A_trans.
get_values();
578 PetscBLASInt LDA = M;
600 PetscBLASInt LDB = x.
size();
605 std::vector<T> S(min_MN);
611 PetscScalar RCOND = rcond;
616 PetscBLASInt RANK = 0;
630 PetscBLASInt LWORK = (3*min_MN + std::max(2*min_MN, std::max(max_MN, NRHS))) * 3/2;
634 std::vector<T> WORK(LWORK);
642 PetscBLASInt INFO = 0;
657 LAPACKgelss_(&M, &N, &NRHS,
pPS(A_trans_vals.data()), &LDA,
658 pPS(
B.data()), &LDB,
pPS(S.data()), &RCOND, &RANK,
659 pPS(WORK.data()), &LWORK, &INFO);
663 libmesh_error_msg(
"Error, argument " << -INFO <<
" to LAPACKgelss_ had an illegal value.");
665 libmesh_error_msg(
"The algorithm for computing the SVD failed to converge!");
691 libmesh_error_msg(
"svd_solve() requires PETSc >= 3.1!");
694 #endif // !PETSC_VERSION_LESS_THAN(3,1,0)
702 libmesh_not_implemented();
704 #endif // (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
708 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
717 if (this->m() != this->n())
718 libmesh_error_msg(
"Can only compute eigen-decompositions for square matrices.");
753 char JOBVL = VL ?
'V' :
'N';
758 char JOBVR = VR ?
'V' :
'N';
762 PetscBLASInt N = this->m();
770 PetscBLASInt LDA = N;
797 PetscBLASInt LDVL = VL ? N : 1;
814 PetscBLASInt LDVR = VR ? N : 1;
828 PetscBLASInt LWORK = (VR || VL) ? 4*N : 3*N;
829 std::vector<T> WORK(LWORK);
838 PetscBLASInt INFO = 0;
841 std::vector<T> & lambda_real_val = lambda_real.
get_values();
842 std::vector<T> & lambda_imag_val = lambda_imag.
get_values();
845 T * VR_ptr =
nullptr;
852 T * VL_ptr =
nullptr;
865 pPS(lambda_real_val.data()),
866 pPS(lambda_imag_val.data()),
877 libmesh_error_msg(
"INFO=" << INFO <<
", Error during Lapack eigenvalue calculation!");
907 libmesh_error_msg(
"_evd_lapack is currently only available when LIBMESH_USE_REAL_NUMBERS is defined!");
916 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
931 PetscBLASInt N = this->m();
937 PetscBLASInt NRHS = 1;
945 PetscBLASInt LDA = N;
967 PetscBLASInt LDB = N;
972 PetscBLASInt INFO = 0;
975 LAPACKgetrs_(TRANS, &N, &NRHS,
pPS(_val.data()), &LDA,
976 _pivots.data(),
pPS(x_vec.data()), &LDB, &INFO);
980 libmesh_error_msg(
"INFO=" << INFO <<
", Error during Lapack LU solve!");
997 libmesh_error_msg(
"No PETSc-provided BLAS/LAPACK available!");
1006 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
1008 template<
typename T>
1020 if ((dest.
size() != this->m()) || (arg.
size() != this->n()))
1021 libmesh_error_msg(
"Improper input argument sizes!");
1029 if ((dest.
size() != this->n()) || (arg.
size() != this->m()))
1030 libmesh_error_msg(
"Improper input argument sizes!");
1048 PetscBLASInt M = this->n();
1053 PetscBLASInt N = this->m();
1072 PetscBLASInt LDA = M;
1088 PetscBLASInt INCX = 1;
1108 PetscBLASInt INCY = 1;
1111 BLASgemv_(TRANS, &M, &N,
pPS(&alpha),
pPS(a.data()), &LDA,
1112 pPS(x.data()), &INCX,
pPS(&beta),
pPS(y.data()), &INCY);
1119 template<
typename T>
1126 libmesh_error_msg(
"No PETSc-provided BLAS/LAPACK available!");
1138 DenseVector<Real> &);
1141 DenseVector<Real> &,
1142 const DenseVector<Real> &,
1146 DenseMatrix<Number> &,
1147 DenseMatrix<Number> &);
1150 std::vector<Real> &,
1151 std::vector<Number> &,
1152 std::vector<Number> &);
1155 DenseVector<Real> &,
1156 DenseMatrix<Real> *,
1157 DenseMatrix<Real> *);
1159 #if !(LIBMESH_USE_REAL_NUMBERS)
1163 DenseVector<Number> &);
1166 DenseVector<Number> &,
1167 const DenseVector<Number> &,
1171 DenseMatrix<Number> &,
1172 DenseMatrix<Number> &);
1175 std::vector<Real> &,
1176 std::vector<Number> &,
1177 std::vector<Number> &);
1180 DenseVector<Number> &,
1181 DenseMatrix<Number> *,
1182 DenseMatrix<Number> *);