OpenCAEPoro  v0.5.0
A simulator for multicomponent porous media flow
DenseMat.hpp
Go to the documentation of this file.
1 
12 #ifndef __DENSEMAT_HEADER__
13 #define __DENSEMAT_HEADER__
14 
15 // Standard header files
16 #include <algorithm>
17 #include <iomanip>
18 #include <iostream>
19 #include <string>
20 
21 using namespace std;
22 
23 extern "C" {
24 
26 
28 void dscal_(const int* n, const double* alpha, double* x, const int* incx);
29 
31 double ddot_(const int* n, double* a, const int* inca, double* b, const int* incb);
32 
34 int dcopy_(
35  const int* n, const double* src, const int* incx, double* dst, const int* incy);
36 
38 int daxpy_(const int* n,
39  const double* alpha,
40  const double* x,
41  const int* incx,
42  double* y,
43  const int* incy);
44 
46 double dnrm2_(const int* n, double* x, const int* incx);
47 
49 double dasum_(const int* n, double* x, const int* incx);
50 
52 int idamax_(const int* n, double* x, const int* incx);
53 
55 int dgemm_(const char* transa,
56  const char* transb,
57  const int* m,
58  const int* n,
59  const int* k,
60  const double* alpha,
61  const double* A,
62  const int* lda,
63  const double* B,
64  const int* ldb,
65  const double* beta,
66  double* C,
67  const int* ldc);
68 
70 
72 int dgesv_(const int* n,
73  const int* nrhs,
74  double* A,
75  const int* lda,
76  int* ipiv,
77  double* b,
78  const int* ldb,
79  int* info);
80 
82 int dsysv_(const char* uplo,
83  const int* n,
84  const int* nrhs,
85  double* A,
86  const int* lda,
87  int* ipiv,
88  double* b,
89  const int* ldb,
90  double* work,
91  const int* lwork,
92  int* info);
93 
96 int ssyevd_(char* jobz,
97  char* uplo,
98  const int* n,
99  float* A,
100  const int* lda,
101  float* w,
102  float* work,
103  const int* lwork,
104  int* iwork,
105  const int* liwork,
106  int* info);
107 }
108 
110 void CalEigenSY(const int& N, float* A, float* w, float* work, const int& lwork);
111 
113 // void MinEigenS(const int& N, float* a, float* w);
114 
116 void Dcopy(const int& N, double* dst, const double* src);
117 
119 double Ddot(int n, double* a, double* b);
120 
122 double Dnorm1(const int& N, double* x);
123 
125 double Dnorm2(const int& N, double* x);
126 
128 void Dscalar(const int& n, const double& alpha, double* x);
129 
131 void Daxpy(const int& n, const double& alpha, const double* x, double* y);
132 
134 void DaABpbC(const int& m,
135  const int& n,
136  const int& k,
137  const double& alpha,
138  const double* A,
139  const double* B,
140  const double& beta,
141  double* C);
142 
143 // test
144 void myDABpC(const int& m,
145  const int& n,
146  const int& k,
147  const double* A,
148  const double* B,
149  double* C);
150 void myDABpCp(const int& m,
151  const int& n,
152  const int& k,
153  const double* A,
154  const double* B,
155  double* C,
156  const int* flag,
157  const int N);
158 void myDABpCp1(const int& m,
159  const int& n,
160  const int& k,
161  const double* A,
162  const double* B,
163  double* C,
164  const int* flag,
165  const int N);
166 void myDABpCp2(const int& m,
167  const int& n,
168  const int& k,
169  const double* A,
170  const double* B,
171  double* C,
172  const int* flag,
173  const int N);
174 
176 void DaAxpby(const int& m,
177  const int& n,
178  const double& a,
179  const double* A,
180  const double* x,
181  const double& b,
182  double* y);
183 
185 void LUSolve(const int& nrhs, const int& N, double* A, double* b, int* pivot);
186 
188 int SYSSolve(const int& nrhs,
189  const char* uplo,
190  const int& N,
191  double* A,
192  double* b,
193  int* pivot,
194  double* work,
195  const int& lwork);
196 
198 template <typename T>
199 void PrintDX(const int& N, const T* x)
200 {
201  for (int i = 0; i < N; i++) {
202  cout << i << " " << setprecision(16) << x[i] << endl;
203  }
204  cout << endl;
205 }
206 
208 template <typename T>
209 bool CheckNan(const int& N, const T* x)
210 {
211  for (int i = 0; i < N; i++) {
212  if (!isfinite(x[i])) {
213  return false;
214  }
215  }
216  return true;
217 }
218 
219 #endif
220 
221 /*----------------------------------------------------------------------------*/
222 /* Brief Change History of This File */
223 /*----------------------------------------------------------------------------*/
224 /* Author Date Actions */
225 /*----------------------------------------------------------------------------*/
226 /* Shizhe Li Oct/24/2021 Create file */
227 /* Chensong Zhang Jan/16/2022 Update Doxygen */
228 /*----------------------------------------------------------------------------*/
double dasum_(const int *n, double *x, const int *incx)
Computes the sum of the absolute values of a vector.
int idamax_(const int *n, double *x, const int *incx)
Finds the index of element having max absolute value.
int dsysv_(const char *uplo, const int *n, const int *nrhs, double *A, const int *lda, int *ipiv, double *b, const int *ldb, double *work, const int *lwork, int *info)
Computes the solution to system of linear equations A * X = B for symm matrices.
int dgesv_(const int *n, const int *nrhs, double *A, const int *lda, int *ipiv, double *b, const int *ldb, int *info)
Computes the solution to system of linear equations A * X = B for general matrices.
double Ddot(int n, double *a, double *b)
Dot product of two double vectors stored as pointers.
Definition: DenseMat.cpp:42
void DaABpbC(const int &m, const int &n, const int &k, const double &alpha, const double *A, const double *B, const double &beta, double *C)
Computes C' = alpha B'A' + beta C', all matrices are column-major.
Definition: DenseMat.cpp:75
int ssyevd_(char *jobz, char *uplo, const int *n, float *A, const int *lda, float *w, float *work, const int *lwork, int *iwork, const int *liwork, int *info)
int daxpy_(const int *n, const double *alpha, const double *x, const int *incx, double *y, const int *incy)
Constant times a vector plus a vector.
void Dscalar(const int &n, const double &alpha, double *x)
Scales a vector by a constant.
Definition: DenseMat.cpp:61
void CalEigenSY(const int &N, float *A, float *w, float *work, const int &lwork)
Calculate the minimal eigenvalue for symmetric matrix with mkl lapack.
Definition: DenseMat.cpp:14
void DaAxpby(const int &m, const int &n, const double &a, const double *A, const double *x, const double &b, double *y)
Computes y = a A x + b y.
Definition: DenseMat.cpp:203
void LUSolve(const int &nrhs, const int &N, double *A, double *b, int *pivot)
Calls dgesv to solve the linear system for general matrices.
Definition: DenseMat.cpp:221
bool CheckNan(const int &N, const T *x)
check NaN
Definition: DenseMat.hpp:209
void Daxpy(const int &n, const double &alpha, const double *x, double *y)
Constant times a vector plus a vector.
Definition: DenseMat.cpp:68
int dcopy_(const int *n, const double *src, const int *incx, double *dst, const int *incy)
Copies a vector, src, to a vector, dst.
int dgemm_(const char *transa, const char *transb, const int *m, const int *n, const int *k, const double *alpha, const double *A, const int *lda, const double *B, const int *ldb, const double *beta, double *C, const int *ldc)
Performs matrix-matrix operations C : = alpha * op(A) * op(B) + beta * C.
double dnrm2_(const int *n, double *x, const int *incx)
Computes the Euclidean norm of a vector.
void dscal_(const int *n, const double *alpha, double *x, const int *incx)
Scales a vector by a constant.
double Dnorm1(const int &N, double *x)
Computes the L1-norm of a vector.
Definition: DenseMat.cpp:49
void PrintDX(const int &N, const T *x)
Prints a vector.
Definition: DenseMat.hpp:199
int SYSSolve(const int &nrhs, const char *uplo, const int &N, double *A, double *b, int *pivot, double *work, const int &lwork)
Calls dsysy to solve the linear system for symm matrices.
Definition: DenseMat.cpp:234
void Dcopy(const int &N, double *dst, const double *src)
Calculate the minimal eigenvalue for symmetric matrix with mkl lapack.
Definition: DenseMat.cpp:36
double ddot_(const int *n, double *a, const int *inca, double *b, const int *incb)
Forms the dot product of two vectors.
double Dnorm2(const int &N, double *x)
Computes the L2-norm of a vector.
Definition: DenseMat.cpp:55