chombo-discharge
CD_LaPackUtils.H
Go to the documentation of this file.
1 /* chombo-discharge
2  * Copyright © 2021 SINTEF Energy Research.
3  * Please refer to Copyright.txt and LICENSE in the chombo-discharge root directory.
4  */
5 
13 #ifndef CD_LaPackUtils_H
14 #define CD_LaPackUtils_H
15 
16 // Std includes
17 #include <vector>
18 
19 // Our includes
20 #include <CD_NamespaceHeader.H>
21 
26 extern "C" void
27 dgelss_(int* a_M,
28  int* a_N,
29  int* a_nRHS,
30  double* a_A,
31  int* a_LDA,
32  double* a_B,
33  int* a_LDB,
34  double* a_S,
35  double* a_rcond,
36  int* a_RANK,
37  double* a_WORK,
38  int* a_LWORK,
39  int* a_INFO);
40 
45 extern "C" void
46 dgesdd_(char* a_JOBZ,
47  int* a_M,
48  int* a_N,
49  double* a_A,
50  int* a_LDA,
51  double* a_S,
52  double* a_U,
53  int* a_LDU,
54  double* a_VT,
55  int* a_LDVT,
56  double* a_WORK,
57  int* a_LWORK,
58  int* a_IWORK,
59  int* a_INFO);
60 
65 extern "C" void
66 sgesdd_(char* a_JOBZ,
67  int* a_M,
68  int* a_N,
69  float* a_A,
70  int* a_LDA,
71  float* a_S,
72  float* a_U,
73  int* a_LDU,
74  float* a_VT,
75  int* a_LDVT,
76  float* a_WORK,
77  int* a_LWORK,
78  int* a_IWORK,
79  int* a_INFO);
80 
84 extern "C" void
85 dgemm_(char* a_TRANSA,
86  char* a_TRANSB,
87  int* a_M,
88  int* a_N,
89  int* a_K,
90  double* a_ALPHA,
91  double* a_A,
92  int* a_LDA,
93  double* a_B,
94  int* a_LDB,
95  double* a_BETA,
96  double* a_C,
97  int* a_LDC);
98 
103 extern "C" void
104 sgemm_(char* a_TRANSA,
105  char* a_TRANSB,
106  int* a_M,
107  int* a_N,
108  int* a_K,
109  float* a_ALPHA,
110  float* a_A,
111  int* a_LDA,
112  float* a_B,
113  int* a_LDB,
114  float* a_BETA,
115  float* a_C,
116  int* a_LDC);
117 
121 extern "C" void
122 dgesv_(int* N, int* NRHS, double* A, int* LDA, int* IPIV, double* B, int* LDB, int* INFO);
123 
127 namespace LaPackUtils {
128 
136  int
137  linearIndex(const int irow, const int jcol, const int M, const int N);
138 
153  bool
154  computeSVD(std::vector<double>& a_linU,
155  std::vector<double>& a_linSigma,
156  std::vector<double>& a_linVtran,
157  const std::vector<double>& a_linA,
158  const int& a_M,
159  const int& a_N);
160 
175  bool
176  computeSVD(std::vector<float>& a_linU,
177  std::vector<float>& a_linSigma,
178  std::vector<float>& a_linVtran,
179  const std::vector<float>& a_linA,
180  const int& a_M,
181  const int& a_N);
182 
190  bool
191  computePseudoInverse(std::vector<double>& a_linAplus,
192  const std::vector<double>& a_linA,
193  const int& a_M,
194  const int& a_N);
195 
203  bool
204  computePseudoInverse(std::vector<float>& a_linAplus,
205  const std::vector<float>& a_linA,
206  const int& a_M,
207  const int& a_N);
208 
218  void
219  linearizeMatrix(std::vector<double>& a_linA,
220  int& a_M,
221  int& a_N,
222  const std::vector<std::vector<double>>& a_A,
223  const char& a_format);
224 
232  void
233  linearizeColumnMajorMatrix(std::vector<double>& a_linA,
234  int& a_M,
235  int& a_N,
236  const std::vector<std::vector<double>>& a_A);
237 
245  void
246  linearizeRowMajorMatrix(std::vector<double>& a_linA, int& a_M, int& a_N, const std::vector<std::vector<double>>& a_A);
247 
257  void
258  deLinearizeMatrix(std::vector<std::vector<double>>& a_A,
259  const int& a_M,
260  const int& a_N,
261  const std::vector<double>& a_linA,
262  const char& a_format);
263 
271  void
272  deLinearizeColumnMajorMatrix(std::vector<std::vector<double>>& a_A,
273  const int& a_M,
274  const int& a_N,
275  const std::vector<double>& a_linA);
276 
284  void
285  deLinearizeRowMajorMatrix(std::vector<std::vector<double>>& a_A,
286  const int& a_M,
287  const int& a_N,
288  const std::vector<double>& a_linA);
289 } // namespace LaPackUtils
290 
291 #include <CD_NamespaceFooter.H>
292 
293 #endif
void dgesv_(int *N, int *NRHS, double *A, int *LDA, int *IPIV, double *B, int *LDB, int *INFO)
Interface to LaPack for computing solution to Ax=b.
void dgelss_(int *a_M, int *a_N, int *a_nRHS, double *a_A, int *a_LDA, double *a_B, int *a_LDB, double *a_S, double *a_rcond, int *a_RANK, double *a_WORK, int *a_LWORK, int *a_INFO)
This is an interface to LAPACK for solving least squares problems using singular-value decomposition.
void sgemm_(char *a_TRANSA, char *a_TRANSB, int *a_M, int *a_N, int *a_K, float *a_ALPHA, float *a_A, int *a_LDA, float *a_B, int *a_LDB, float *a_BETA, float *a_C, int *a_LDC)
Interface to LaPack for matrix multiplication.
void dgesdd_(char *a_JOBZ, int *a_M, int *a_N, double *a_A, int *a_LDA, double *a_S, double *a_U, int *a_LDU, double *a_VT, int *a_LDVT, double *a_WORK, int *a_LWORK, int *a_IWORK, int *a_INFO)
Interface to LaPack for computing the singular value decomposition of a matrix.
void sgesdd_(char *a_JOBZ, int *a_M, int *a_N, float *a_A, int *a_LDA, float *a_S, float *a_U, int *a_LDU, float *a_VT, int *a_LDVT, float *a_WORK, int *a_LWORK, int *a_IWORK, int *a_INFO)
Interface to LaPack for computing the singular value decomposition of a matrix.
void dgemm_(char *a_TRANSA, char *a_TRANSB, int *a_M, int *a_N, int *a_K, double *a_ALPHA, double *a_A, int *a_LDA, double *a_B, int *a_LDB, double *a_BETA, double *a_C, int *a_LDC)
Interface to LaPack for matrix multiplication.
Namespace containing various useful linear algebra routines using LaPACK.
Definition: CD_LaPackUtils.H:127
void deLinearizeColumnMajorMatrix(std::vector< std::vector< double >> &a_A, const int &a_M, const int &a_N, const std::vector< double > &a_linA)
Delinearize a linearized matrix from column major Fortran to column major form.
Definition: CD_LaPackUtils.cpp:417
void linearizeColumnMajorMatrix(std::vector< double > &a_linA, int &a_M, int &a_N, const std::vector< std::vector< double >> &a_A)
Linearize a matrix to column major Fortran form by column major format of the input matrix.
Definition: CD_LaPackUtils.cpp:359
void deLinearizeRowMajorMatrix(std::vector< std::vector< double >> &a_A, const int &a_M, const int &a_N, const std::vector< double > &a_linA)
Delinearize a linearized matrix from column major Fortran to row major form.
Definition: CD_LaPackUtils.cpp:438
bool computeSVD(std::vector< double > &a_linU, std::vector< double > &a_linSigma, std::vector< double > &a_linVtran, const std::vector< double > &a_linA, const int &a_M, const int &a_N)
Compute the singular value decomposition of a matrix.
Definition: CD_LaPackUtils.cpp:29
int linearIndex(const int irow, const int jcol, const int M, const int N)
Get the index in the linearized matrix of size MxN.
Definition: CD_LaPackUtils.cpp:23
void linearizeMatrix(std::vector< double > &a_linA, int &a_M, int &a_N, const std::vector< std::vector< double >> &a_A, const char &a_format)
Linearize a matrix to column major Fortran form by assuming row or major colum format of the input ma...
Definition: CD_LaPackUtils.cpp:399
bool computePseudoInverse(std::vector< double > &a_linAplus, const std::vector< double > &a_linA, const int &a_M, const int &a_N)
Compute the pseudoinverse of matrix through singular value decomposition.
Definition: CD_LaPackUtils.cpp:163
void linearizeRowMajorMatrix(std::vector< double > &a_linA, int &a_M, int &a_N, const std::vector< std::vector< double >> &a_A)
Linearize a matrix to column major Fortran form by row major format of the input matrix.
Definition: CD_LaPackUtils.cpp:379
void deLinearizeMatrix(std::vector< std::vector< double >> &a_A, const int &a_M, const int &a_N, const std::vector< double > &a_linA, const char &a_format)
Delinearize a linearized matrix from column major Fortran to row or major colum matrix format.
Definition: CD_LaPackUtils.cpp:459