36#ifndef OPENRS_BLAS_LAPACK_HEADER
37#define OPENRS_BLAS_LAPACK_HEADER
39#if ! __has_include(<FCMacros.h>)
40#error "Need FCMacros.h to be generated through FortranCInterface"
43#include <opm/common/ErrorMacros.hpp>
44#include <opm/porsol/common/fortran.hpp>
52#define F77_CHARACTER_TYPE const char*
55void FC_GLOBAL(dgemv,DGEMV)(F77_CHARACTER_TYPE,
56 const int* m ,
const int* n,
57 const double* a1 ,
const double* A,
const int* ldA ,
58 const double* x,
const int* incX,
59 const double* a2 ,
double* y,
const int* incY);
63void FC_GLOBAL(dgemm,DGEMM)(F77_CHARACTER_TYPE, F77_CHARACTER_TYPE,
64 const int* m ,
const int* n ,
const int* k ,
65 const double* a1 ,
const double* A ,
const int* ldA,
66 const double* B ,
const int* ldB,
67 const double* a2 ,
double* C ,
const int* ldC);
70void FC_GLOBAL(dsyrk,DSYRK)(F77_CHARACTER_TYPE, F77_CHARACTER_TYPE,
71 const int* n ,
const int* k ,
72 const double* a1 ,
const double* A ,
const int* ldA,
73 const double* a2 ,
double* C ,
const int* ldC);
76void FC_GLOBAL(dtrmm,DTRMM)(F77_CHARACTER_TYPE, F77_CHARACTER_TYPE,
77 F77_CHARACTER_TYPE, F77_CHARACTER_TYPE,
78 const int* m ,
const int* n ,
80 const double* A ,
const int* ldA ,
81 double* B ,
const int* ldB);
84void FC_GLOBAL(dgeqrf,DGEQRF)(
const int* m ,
const int* n ,
85 double* A ,
const int* ld ,
86 double* tau ,
double* work,
87 const int* lwork,
int* info);
89void FC_GLOBAL(dorgqr,DORGQR)(
const int* m ,
const int* n ,
const int* k ,
90 double* A ,
const int* ld ,
const double* tau,
91 double* work,
const int* lwork,
int* info);
93void FC_GLOBAL(dgetrf,DGETRF)(
const int* m ,
const int* n ,
94 double* A ,
const int* ld,
95 int* ipiv,
int* info);
97void FC_GLOBAL(dgetri,DGETRI)(
const int* n ,
98 double* A ,
const int* ld,
100 double* work,
int* lwork,
int* info);
107 namespace BLAS_LAPACK {
152 void GEMV(
const char* transA,
153 const int m ,
const int n,
154 const T& a1 ,
const T* A,
const int ldA ,
155 const T* x,
const int incX,
156 const T& a2 , T* y,
const int incY);
160 void GEMV<double>(
const char* transA,
161 const int m ,
const int n,
162 const double& a1 ,
const double* A,
const int ldA,
163 const double* x,
const int incX,
164 const double& a2 ,
double* y,
const int incY);
227 void GEMM(
const char* transA,
const char* transB,
228 const int m ,
const int n ,
const int k ,
229 const T& a1 ,
const T* A ,
const int ldA,
230 const T* B ,
const int ldB,
231 const T& a2 , T* C ,
const int ldC);
235 void GEMM<double>(
const char* transA,
const char* transB,
236 const int m ,
const int n ,
const int k ,
237 const double& a1 ,
const double* A ,
const int ldA,
238 const double* B ,
const int ldB,
239 const double& a2 ,
double* C ,
const int ldC);
247 void SYRK(
const char* uplo,
const char* trans,
248 const int n ,
const int k ,
249 const T& a1 ,
const T* A ,
const int ldA,
250 const T& a2 , T* C ,
const int ldC);
254 void SYRK<double>(
const char* uplo,
const char* trans,
255 const int n ,
const int k ,
256 const double& a1 ,
const double* A ,
const int ldA,
257 const double& a2 ,
double* C ,
const int ldC);
264 void TRMM(
const char* side ,
const char* uplo,
265 const char* transA,
const char* diag,
266 const int m ,
const int n ,
const T& a,
267 const T* A ,
const int ldA ,
268 T* B ,
const int ldB);
272 void TRMM<double>(
const char* side ,
const char* uplo,
273 const char* transA,
const char* diag,
274 const int m ,
const int n ,
const double& a,
275 const double* A ,
const int ldA ,
276 double* B ,
const int ldB);
283 void GEQRF(
const int m ,
const int n ,
284 T* A ,
const int ld ,
286 const int lwork,
int& info);
290 void GEQRF<double>(
const int m ,
const int n ,
291 double* A ,
const int ld ,
292 double* tau ,
double* work,
293 const int lwork,
int& info);
300 void ORGQR(
const int m ,
const int n ,
const int k ,
301 T* A ,
const int ld ,
const T* tau,
302 T* work,
const int lwork,
int& info);
308 void ORGQR<double>(
const int m ,
const int n ,
const int k ,
309 double* A ,
const int ld ,
const double* tau,
310 double* work,
const int lwork,
int& info);
317 void GETRF(
const int m,
const int n, T* A,
318 const int ld,
int* ipiv,
int& info);
323 void GETRF<double>(
const int m,
const int n ,
double* A,
324 const int ld,
int* ipiv,
int& info);
331 void GETRI(
const int n , T* A ,
const int ld,
332 const int* ipiv, T* work,
int lwork,
int& info);
336 void GETRI(
const int n ,
double* A ,
const int ld,
337 const int* ipiv,
double* work,
int lwork,
int& info);
Inverting small matrices.
Definition ImplicitAssembly.hpp:43