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),
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(),
258 libmesh_error_msg_if(INFO != 0,
"INFO=" << INFO <<
", Error during Lapack LU factorization!");
261 this->_decomposition_type = LU_BLAS_LAPACK;
269 libmesh_error_msg(
"No PETSc-provided BLAS/LAPACK available!");
305 std::vector<Real> sigma_val;
306 std::vector<Number> U_val;
307 std::vector<Number> VT_val;
309 _svd_helper(JOBU, JOBVT, sigma_val, U_val, VT_val);
312 sigma.
resize(cast_int<unsigned int>(sigma_val.size()));
313 for (
auto i : make_range(sigma.
size()))
314 sigma(i) = sigma_val[i];
354 std::vector<Real> sigma_val;
357 int min_MN = (M < N) ? M : N;
376 sigma.
resize(cast_int<unsigned int>(sigma_val.size()));
377 for (
auto i : make_range(sigma.
size()))
378 sigma(i) = sigma_val[i];
381 #if (LIBMESH_HAVE_PETSC)
386 std::vector<Real> & sigma_val,
387 std::vector<Number> & U_val,
388 std::vector<Number> & VT_val)
394 PetscBLASInt M = this->n();
399 PetscBLASInt N = this->m();
401 PetscBLASInt min_MN = (M < N) ? M : N;
402 PetscBLASInt max_MN = (M > N) ? M : N;
417 PetscBLASInt LDA = M;
421 sigma_val.resize( min_MN );
426 PetscBLASInt LDU = M;
435 U_val.resize( LDU*min_MN );
437 U_val.resize( LDU*M );
442 PetscBLASInt LDVT = N;
452 VT_val.resize( LDVT*N );
463 PetscBLASInt larger = (3*min_MN+max_MN > 5*min_MN) ? 3*min_MN+max_MN : 5*min_MN;
464 PetscBLASInt LWORK = (larger > 1) ? larger : 1;
474 std::vector<Number> WORK( LWORK );
483 PetscBLASInt INFO = 0;
486 #ifdef LIBMESH_USE_REAL_NUMBERS
488 LAPACKgesvd_(&JOBU, &JOBVT, &M, &N,
pPS(_val.data()), &LDA,
489 pPS(sigma_val.data()),
pPS(U_val.data()), &LDU,
490 pPS(VT_val.data()), &LDVT,
pPS(WORK.data()), &LWORK,
496 std::vector<Number> val_copy(_val.size());
497 for (
auto i : make_range(_val.size()))
498 val_copy[i] = _val[i];
500 std::vector<Real> RWORK(5 * min_MN);
501 LAPACKgesvd_(&JOBU, &JOBVT, &M, &N, val_copy.data(), &LDA, sigma_val.data(), U_val.data(),
502 &LDU, VT_val.data(), &LDVT, WORK.data(), &LWORK, RWORK.data(), &INFO);
506 libmesh_error_msg_if(INFO != 0,
"INFO=" << INFO <<
", Error during Lapack SVD calculation!");
516 std::vector<Number> &,
517 std::vector<Number> &)
519 libmesh_error_msg(
"No PETSc-provided BLAS/LAPACK available!");
526 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
542 this->get_transpose(A_trans);
547 PetscBLASInt M = A_trans.
n();
552 PetscBLASInt N = A_trans.
m();
555 PetscBLASInt max_MN = std::max(M,N);
556 PetscBLASInt min_MN = std::min(M,N);
563 PetscBLASInt NRHS = 1;
571 std::vector<T> & A_trans_vals = A_trans.
get_values();
575 PetscBLASInt LDA = M;
589 for (
auto i : make_range(rhs.
size()))
597 PetscBLASInt LDB = x.
size();
602 std::vector<T> S(min_MN);
608 PetscScalar RCOND = rcond;
613 PetscBLASInt RANK = 0;
627 PetscBLASInt LWORK = (3*min_MN + std::max(2*min_MN, std::max(max_MN, NRHS))) * 3/2;
631 std::vector<T> WORK(LWORK);
639 PetscBLASInt INFO = 0;
654 LAPACKgelss_(&M, &N, &NRHS,
pPS(A_trans_vals.data()), &LDA,
655 pPS(B.data()), &LDB,
pPS(S.data()), &RCOND, &RANK,
656 pPS(WORK.data()), &LWORK, &INFO);
659 libmesh_error_msg_if(INFO < 0,
"Error, argument " << -INFO <<
" to LAPACKgelss_ had an illegal value.");
660 libmesh_error_msg_if(INFO > 0,
"The algorithm for computing the SVD failed to converge!");
686 libmesh_not_implemented();
688 #endif // (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
692 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
701 libmesh_error_msg_if(this->m() != this->n(),
"Can only compute eigen-decompositions for square matrices.");
712 for (
auto i : make_range(this->_m))
713 for (
auto j : make_range(i))
714 std::swap((*
this)(i,j), (*
this)(j,i));
736 char JOBVL = VL ?
'V' :
'N';
741 char JOBVR = VR ?
'V' :
'N';
745 PetscBLASInt N = this->m();
753 PetscBLASInt LDA = N;
780 PetscBLASInt LDVL = VL ? N : 1;
797 PetscBLASInt LDVR = VR ? N : 1;
811 PetscBLASInt LWORK = (VR || VL) ? 4*N : 3*N;
812 std::vector<T> WORK(LWORK);
821 PetscBLASInt INFO = 0;
824 std::vector<T> & lambda_real_val = lambda_real.
get_values();
825 std::vector<T> & lambda_imag_val = lambda_imag.
get_values();
828 T * VR_ptr =
nullptr;
835 T * VL_ptr =
nullptr;
848 pPS(lambda_real_val.data()),
849 pPS(lambda_imag_val.data()),
859 libmesh_error_msg_if(INFO != 0,
"INFO=" << INFO <<
", Error during Lapack eigenvalue calculation!");
868 for (
auto i : make_range(N))
869 for (
auto j : make_range(i))
870 std::swap((*VR)(i,j), (*VR)(j,i));
875 for (
auto i : make_range(N))
876 for (
auto j : make_range(i))
877 std::swap((*VL)(i,j), (*VL)(j,i));
889 libmesh_error_msg(
"_evd_lapack is currently only available when LIBMESH_USE_REAL_NUMBERS is defined!");
898 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
913 PetscBLASInt N = this->m();
919 PetscBLASInt NRHS = 1;
927 PetscBLASInt LDA = N;
949 PetscBLASInt LDB = N;
954 PetscBLASInt INFO = 0;
957 LAPACKgetrs_(TRANS, &N, &NRHS,
pPS(_val.data()), &LDA,
958 _pivots.data(),
pPS(x_vec.data()), &LDB, &INFO);
961 libmesh_error_msg_if(INFO != 0,
"INFO=" << INFO <<
", Error during Lapack LU solve!");
978 libmesh_error_msg(
"No PETSc-provided BLAS/LAPACK available!");
987 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
1001 libmesh_error_msg_if((dest.
size() != this->m()) || (arg.
size() != this->n()),
1002 "Improper input argument sizes!");
1010 libmesh_error_msg_if((dest.
size() != this->n()) || (arg.
size() != this->m()),
1011 "Improper input argument sizes!");
1029 PetscBLASInt M = this->n();
1034 PetscBLASInt N = this->m();
1053 PetscBLASInt LDA = M;
1069 PetscBLASInt INCX = 1;
1089 PetscBLASInt INCY = 1;
1092 BLASgemv_(TRANS, &M, &N,
pPS(&alpha),
pPS(a.data()), &LDA,
1093 pPS(x.data()), &INCX,
pPS(&beta),
pPS(y.data()), &INCY);
1100 template<
typename T>
1107 libmesh_error_msg(
"No PETSc-provided BLAS/LAPACK available!");
1127 DenseMatrix<Number> &,
1128 DenseMatrix<Number> &);
1131 std::vector<Real> &,
1132 std::vector<Number> &,
1133 std::vector<Number> &);
1140 #if !(LIBMESH_USE_REAL_NUMBERS)
1144 DenseVector<Number> &);
1147 DenseVector<Number> &,
1148 const DenseVector<Number> &,
1152 DenseMatrix<Number> &,
1153 DenseMatrix<Number> &);
1156 std::vector<Real> &,
1157 std::vector<Number> &,
1158 std::vector<Number> &);
1161 DenseVector<Number> &,
1162 DenseMatrix<Number> *,
1163 DenseMatrix<Number> *);
void _matvec_blas(T alpha, T beta, DenseVector< T > &dest, const DenseVector< T > &arg, bool trans=false) const
使用BLAS GEMV函数(通过PETSc)计算
unsigned int n() const
返回矩阵的列维度。
void _lu_decompose_lapack()
使用Lapack例程“getrf”计算矩阵的LU分解。此例程只能由lu_solve()函数的“use_blas_lapack”分支使用。 调用此函数后,矩阵将被其因式分解版本替换,并且Decomposi...
_BLAS_Multiply_Flag
用于确定_multiply_blas函数行为的枚举。
template class LIBMESH_EXPORT DenseVector< Real >
void resize(const unsigned int n)
调整向量的大小。将所有元素设置为0。
std::vector< T > & get_values()
PetscScalar * pPS(T *ptr)
void _multiply_blas(const DenseMatrixBase< T > &other, _BLAS_Multiply_Flag flag)
_multiply_blas函数使用BLAS gemm函数计算A <- op(A) * op(B)。 用于right_multiply()、left_multiply()、right_multiply_...
unsigned int m() const
返回矩阵的行维度。
void _lu_back_substitute_lapack(const DenseVector< T > &b, DenseVector< T > &x)
与_lu_decompose_lapack()相对应的伴随函数。不要直接使用,通过公共lu_solve()接口调用。 由于我们只是调用LAPACK例程,该函数在逻辑上是const的,因为它不修改矩阵, ...
template class LIBMESH_EXPORT DenseMatrix< Real >
void _svd_lapack(DenseVector< Real > &sigma)
使用Lapack例程“getsvd”计算矩阵的奇异值分解。 [ 在dense_matrix_blas_lapack.C中实现 ]
std::vector< T > & get_values()
返回对应于存储向量的引用的底层数据存储矢量。
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void _evd_lapack(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > *VL=nullptr, DenseMatrix< T > *VR=nullptr)
使用Lapack例程“DGEEV”计算矩阵的特征值。如果VR和/或VL不为nullptr, 则此函数还计算并返回右和/或左特征向量的矩阵。 [ 在dense_matrix_blas_lapack.C中实现 ]
void resize(const unsigned int new_m, const unsigned int new_n)
调整矩阵的大小,并调用 zero() 方法。
定义用于有限元计算的稠密向量类。该类基本上是为了补充 DenseMatrix 类而设计的。 它相对于 std::vector 具有额外的功能,使其在有限元中特别有用,特别是对于方程组。 所有重写的虚拟函...
void _svd_solve_lapack(const DenseVector< T > &rhs, DenseVector< T > &x, Real rcond) const
由svd_solve(rhs)调用。
为有限元类型的计算定义了一个抽象的稠密矩阵基类。例如 DenseSubMatrices 可以从这个类派生出来。
void _svd_helper(char JOBU, char JOBVT, std::vector< Real > &sigma_val, std::vector< Number > &U_val, std::vector< Number > &VT_val)
实际执行SVD的辅助函数。 [ 在dense_matrix_blas_lapack.C中实现 ]
定义用于有限元类型计算的密集矩阵。 用于在求和成全局矩阵之前存储单元刚度矩阵。所有被覆盖的虚函数都记录在dense_matrix_base.h中。
virtual unsigned int size() const overridefinal
template class LIBMESH_EXPORT DenseMatrixBase< Real >