20 #define Subroutine static inline void
25 Subroutine dgecon (
char norm,
int n,
const Real* A,
int lda,
26 Real anorm,
Real* rcond,
Real* work,
int* iwork,
int& info)
27 {
dgecon_(&norm,&n,
const_cast<Real*
>(A),&lda,&anorm,rcond,work,iwork,&info); }
32 Subroutine dgesv (
int n,
int nrhs,
33 Real* A,
int lda,
int* ipiv,
34 Real* B,
int ldb,
int& info)
35 {
dgesv_(&n,&nrhs,A,&lda,ipiv,B,&ldb,&info); }
40 Subroutine dgetrf (
int m,
int n,
Real* A,
int lda,
42 {
dgetrf_(&m,&n,A,&lda,ipiv,&info); }
47 Subroutine dgetri (
int n,
Real* A,
int lda,
48 int* ipiv,
Real* work,
int lwork,
int& info)
49 {
dgetri_(&n,A,&lda,ipiv,work,&lwork,&info); }
54 Subroutine dgetrs (
char trans,
int n,
int nrhs,
55 Real* A,
int lda,
int* ipiv,
56 Real* B,
int ldb,
int& info)
57 {
dgetrs_(&trans,&n,&nrhs,A,&lda,ipiv,B,&ldb,&info); }
62 static inline double dlange (
char norm,
int m,
int n,
const Real* A,
64 {
return dlange_(&norm,&m,&n,
const_cast<Real*
>(A),&lda,work); }
69 Subroutine dposv (
char uplo,
int n,
int nrhs,
70 Real* A,
int lda,
Real* B,
int ldb,
int& info)
71 {
dposv_(&uplo,&n,&nrhs,A,&lda,B,&ldb,&info); }
76 Subroutine dpotrs (
char uplo,
int n,
int nrhs,
77 Real* A,
int lda,
Real* B,
int ldb,
int& info)
78 {
dpotrs_(&uplo,&n,&nrhs,A,&lda,B,&ldb,&info); }
83 Subroutine dsyev (
char jobz,
char uplo,
84 int n,
double* A,
int lda,
double* w,
85 double* work,
int lwork,
int& info)
86 {
dsyev_(&jobz,&uplo,&n,A,&lda,w,work,&lwork,&info); }
91 Subroutine dsyevx (
char jobz,
char range,
char uplo,
92 int n,
Real* A,
int lda,
95 Real* work,
int lwork,
int* iwork,
int* ifail,
int& info)
96 {
dsyevx_(&jobz,&range,&uplo,&n,A,&lda,&vl,&vu,&il,&iu,&abstol,
97 &m,w,z,&ldz,work,&lwork,iwork,ifail,&info); }
102 Subroutine dsygvx (
int itype,
char jobz,
char range,
103 char uplo,
int n,
Real* A,
int lda,
105 int il,
int iu,
Real abstol,
107 Real* work,
int lwork,
int* iwork,
int* ifail,
int& info)
108 {
dsygvx_(&itype,&jobz,&range,&uplo,&n,A,&lda,B,&ldb,&vl,&vu,&il,&iu,&abstol,
109 &m,w,z,&ldz,work,&lwork,iwork,ifail,&info); }
114 Subroutine dgeev (
char jobvl,
char jobvr,
115 int n,
Real* A,
int lda,
118 Real* work,
int lwork,
int& info)
119 {
dgeev_(&jobvl,&jobvr,&n,A,&lda,wr,wi,vl,&ldvl,vr,&ldvr,work,&lwork,&info); }
123 #if defined(_WIN32) && !defined(__MINGW32__) && !defined(__MINGW64__)
124 #define dgecon_ DGECON
126 #define dgetrf_ DGETRF
127 #define dgetri_ DGETRI
128 #define dgetrs_ DGETRS
129 #define dlange_ DLANGE
131 #define dpotrs_ DPOTRS
133 #define dsyevx_ DSYEVX
134 #define dsygvx_ DSYGVX
137 #define dgecon_ dgecon
139 #define dgetrf_ dgetrf
140 #define dgetri_ dgetri
141 #define dgetrs_ dgetrs
142 #define dlange_ dlange
144 #define dpotrs_ dpotrs
146 #define dsyevx_ dsyevx
147 #define dsygvx_ dsygvx
155 void dgecon_(
const char& norm,
const int& n,
const Real* A,
const int& lda,
156 const Real& anorm,
Real* rcond,
Real* work,
int* iwork,
int& info);
161 void dgesv_(
const int& n,
const int& nrhs,
162 Real* A,
const int& lda,
int* ipiv,
163 Real* B,
const int& ldb,
int& info);
169 int* ipiv,
int& info);
175 int* ipiv,
Real* work,
const int& lwork,
int& info);
180 void dgetrs_(
const char& trans,
const int& n,
const int& nrhs,
181 Real* A,
const int& lda,
int* ipiv,
182 Real* B,
const int& ldb,
int& info);
187 double dlange_(
const char& norm,
const int& m,
const int& n,
const Real* A,
188 const int& lda,
Real* work);
193 void dposv_(
const char& uplo,
const int& n,
const int& nrhs,
194 Real* A,
const int& lda,
Real* B,
const int& ldb,
int& info);
199 void dpotrs_(
const char& uplo,
const int& n,
const int& nrhs,
200 Real* A,
const int& lda,
Real* B,
const int& ldb,
int& info);
205 void dsyev_(
const char& jobz,
const char& uplo,
206 const int& n,
double* A,
const int& lda,
double* w,
207 double* work,
const int& lwork,
int& info);
212 void dsyevx_(
const char& jobz,
const char& range,
const char& uplo,
213 const int& n,
Real* a,
const int& lda,
214 const Real& vl,
const Real& vu,
const int& il,
const int& iu,
215 const Real& abstol,
int& m,
Real* w,
Real* z,
const int& ldz,
216 Real* work,
const int& lwork,
int* iwork,
int* ifail,
int& info);
221 void dsygvx_(
const int& itype,
const char& jobz,
const char& range,
222 const char& uplo,
const int& n,
Real* a,
const int& lda,
223 Real* b,
const int& ldb,
const Real& vl,
const Real& vu,
224 const int& il,
const int& iu,
const Real& abstol,
225 int& m,
Real* w,
Real* z,
const int& ldz,
226 Real* work,
const int& lwork,
int* iwork,
int* ifail,
int& info);
231 void dgeev_(
const char& jobvl,
const char& jobvr,
232 const int& n,
Real* a,
const int& lda,
234 Real* vr,
const int& ldvr,
235 Real* work,
const int& lwork,
int& info);
BLAS support for various platforms.
#define Real
The floating point type to use.
Definition: ImmersedBoundaries.h:18
void dgetri_(const int &n, Real *A, const int &lda, int *ipiv, Real *work, const int &lwork, int &info)
Computes the inverse of a matrix based on its LU factorization.
void dsyevx_(const char &jobz, const char &range, const char &uplo, const int &n, Real *a, const int &lda, const Real &vl, const Real &vu, const int &il, const int &iu, const Real &abstol, int &m, Real *w, Real *z, const int &ldz, Real *work, const int &lwork, int *iwork, int *ifail, int &info)
Solves the standard eigenproblem A*x=(lambda)*x.
void dgetrf_(const int &m, const int &n, Real *A, const int &lda, int *ipiv, int &info)
Computes an LU factorization of a general M-by-N matrix A.
void dgetrs_(const char &trans, const int &n, const int &nrhs, Real *A, const int &lda, int *ipiv, Real *B, const int &ldb, int &info)
Solves the equation system A*x=b when A is already factorized.
void dsygvx_(const int &itype, const char &jobz, const char &range, const char &uplo, const int &n, Real *a, const int &lda, Real *b, const int &ldb, const Real &vl, const Real &vu, const int &il, const int &iu, const Real &abstol, int &m, Real *w, Real *z, const int &ldz, Real *work, const int &lwork, int *iwork, int *ifail, int &info)
Solves the generalized eigenproblem A*x=(lambda)*B*x.
void dsyev_(const char &jobz, const char &uplo, const int &n, double *A, const int &lda, double *w, double *work, const int &lwork, int &info)
Solves the standard eigenproblem A*x=(lambda)*x.
void dposv_(const char &uplo, const int &n, const int &nrhs, Real *A, const int &lda, Real *B, const int &ldb, int &info)
Solves the symmetric linear equation system A*x=b.
void dgecon_(const char &norm, const int &n, const Real *A, const int &lda, const Real &anorm, Real *rcond, Real *work, int *iwork, int &info)
Estimates the reciprocal of the condition number of the matrix A.
void dgesv_(const int &n, const int &nrhs, Real *A, const int &lda, int *ipiv, Real *B, const int &ldb, int &info)
Solves the linear equation system A*x=b.
double dlange_(const char &norm, const int &m, const int &n, const Real *A, const int &lda, Real *work)
Returns a norm of the matrix A.
void dgeev_(const char &jobvl, const char &jobvr, const int &n, Real *a, const int &lda, Real *wr, Real *wi, Real *vl, const int &ldvl, Real *vr, const int &ldvr, Real *work, const int &lwork, int &info)
Solves the non-symmetric eigenproblem A*x=(lambda)*x.
void dpotrs_(const char &uplo, const int &n, const int &nrhs, Real *A, const int &lda, Real *B, const int &ldb, int &info)
Solves the symmetric equation system A*x=b for prefactored A.