7#if YAC_LAPACK_INTERFACE_ID == 3
12lapack_int LAPACKE_dgels_work(
int matrix_layout,
char trans, lapack_int m,
13 lapack_int n, lapack_int nrhs,
double* a,
14 lapack_int lda,
double* b, lapack_int ldb,
15 double* work, lapack_int lwork )
17 assert(matrix_layout == LAPACK_COL_MAJOR);
18 return (lapack_int) clapack_dgels(LAPACK_COL_MAJOR,
19 trans ==
'N' ? CblasNoTrans : CblasTrans,
20 (ATL_INT) m, (ATL_INT) n, (ATL_INT) nrhs,
21 a, (ATL_INT) lda, b, (int) ldb);
25lapack_int LAPACKE_dgesv(
int matrix_layout, lapack_int n, lapack_int nrhs,
26 double* a, lapack_int lda, lapack_int* ipiv,
27 double* b, lapack_int ldb )
29 assert(matrix_layout == LAPACK_COL_MAJOR);
30 if (
sizeof(
int) ==
sizeof(lapack_int))
32 return (lapack_int) clapack_dgesv(LAPACK_COL_MAJOR, (
int) n, (
int) nrhs,
33 a, (
int) lda, (
int*) ipiv, b, (
int) ldb);
37 int i, result, ipiv_size = (int) n;
39 result = clapack_dgesv(LAPACK_COL_MAJOR, (
int) n, (
int) nrhs,
40 a, (
int) lda, ipiv_, b, (
int) ldb);
41 for (i = 0; i != ipiv_size; ++i) { ipiv[i] = (lapack_int) ipiv_[i]; }
42 return (lapack_int) result;
46lapack_int LAPACKE_dgetrf(
int matrix_layout, lapack_int m, lapack_int n,
47 double* a, lapack_int lda, lapack_int* ipiv )
49 assert(matrix_layout == LAPACK_COL_MAJOR);
50 if (
sizeof(
int) ==
sizeof(lapack_int))
52 return (lapack_int) clapack_dgetrf(LAPACK_COL_MAJOR, (
int) m, (
int) n,
53 a, (
int) lda, (
int*) ipiv);
57 int i, result, ipiv_size = (int) (n < m ? n : m);
59 result = clapack_dgetrf(LAPACK_COL_MAJOR, (
int) m, (
int) n,
61 for (i = 0; i != ipiv_size; ++i) { ipiv[i] = (lapack_int) ipiv_[i]; }
62 return (lapack_int) result;
66lapack_int LAPACKE_dgetri_work(
int matrix_layout, lapack_int n,
double* a,
67 lapack_int lda,
const lapack_int* ipiv,
68 double* work, lapack_int lwork )
70 assert(matrix_layout == LAPACK_COL_MAJOR);
71 if (
sizeof(
int) ==
sizeof(lapack_int))
73 return (lapack_int) clapack_dgetri(LAPACK_COL_MAJOR, (
int) n, a,
74 (
int) lda, (
int*) ipiv);
78 int i, ipiv_size = (int) n;
80 for (i = 0; i != ipiv_size; ++i) { ipiv_[i] = (int) ipiv[i]; }
81 return (lapack_int) clapack_dgetri(LAPACK_COL_MAJOR, (
int) n, a,
86#elif YAC_LAPACK_INTERFACE_ID == 4
92lapack_int LAPACKE_dgels_work(
int matrix_layout,
char trans, lapack_int m,
93 lapack_int n, lapack_int nrhs,
double* a,
94 lapack_int lda,
double* b, lapack_int ldb,
95 double* work, lapack_int lwork )
98 assert(matrix_layout == LAPACK_COL_MAJOR &&
99 sizeof(doublereal) ==
sizeof(
double));
100 if (
sizeof(integer) ==
sizeof(lapack_int))
102 dgels_(&trans, (integer*) &m, (integer*) &n, (integer*) &nrhs,
103 (doublereal*) a, (integer*) &lda,
104 (doublereal*) b, (integer*) &ldb,
105 (doublereal*) work, (integer*) &lwork,
110 integer m_ = (integer) m, n_ = (integer) n,
111 nrhs_ = (integer) nrhs, lda_ = (integer) lda,
112 ldb_ = (integer) ldb, lwork_ = (integer) lwork;
113 dgels_(&trans, &m_, &n_, &nrhs_,
114 (doublereal*) a, &lda_,
115 (doublereal*) b, &ldb_,
116 (doublereal*) work, &lwork_,
119 return (lapack_int) info;
123lapack_int LAPACKE_dgesv(
int matrix_layout, lapack_int n, lapack_int nrhs,
124 double* a, lapack_int lda, lapack_int* ipiv,
125 double* b, lapack_int ldb )
128 assert(matrix_layout == LAPACK_COL_MAJOR &&
129 sizeof(doublereal) ==
sizeof(
double));
130 if (
sizeof(integer) ==
sizeof(lapack_int))
132 dgesv_((integer*) &n, (integer*) &nrhs,
133 (doublereal*) a, (integer*) &lda , (integer*) ipiv,
134 (doublereal*) b, (integer*) &ldb, &info);
138 int i, ipiv_size = (int) n;
139 integer n_ = (integer) n, nrhs_ = (integer) nrhs,
140 lda_ = (integer) lda, ldb_ = (integer) ldb,
143 (doublereal*) a, &lda_, ipiv_,
144 (doublereal*) b, &ldb_, &info);
145 for (i = 0; i != ipiv_size; ++i) { ipiv[i] = (lapack_int) ipiv_[i]; }
147 return (lapack_int) info;
150lapack_int LAPACKE_dgetrf(
int matrix_layout, lapack_int m, lapack_int n,
151 double* a, lapack_int lda, lapack_int* ipiv )
154 assert(matrix_layout == LAPACK_COL_MAJOR &&
155 sizeof(doublereal) ==
sizeof(
double));
156 if (
sizeof(integer) ==
sizeof(lapack_int))
158 dgetrf_((integer*) &m, (integer*) &n,
159 (doublereal*) a, (integer*) &lda , (integer*) ipiv,
164 int i, ipiv_size = (int) (n < m ? n : m);
165 integer m_ = (integer) m, n_ = (integer) n,
166 lda_ = (integer) lda, ipiv_[ipiv_size];
168 (doublereal*) a, &lda_, ipiv_,
170 for (i = 0; i != ipiv_size; ++i) {ipiv[i] = (lapack_int) ipiv_[i]; }
172 return (lapack_int) info;
175lapack_int LAPACKE_dgetri_work(
int matrix_layout, lapack_int n,
double* a,
176 lapack_int lda,
const lapack_int* ipiv,
177 double* work, lapack_int lwork )
180 assert(matrix_layout == LAPACK_COL_MAJOR &&
181 sizeof(doublereal) ==
sizeof(
double));
182 if (
sizeof(integer) ==
sizeof(lapack_int))
184 dgetri_((integer*) &n, (doublereal*) a,
185 (integer*) &lda , (integer*) ipiv,
186 (doublereal*) work, (integer*) &lwork,
191 int i, size_ipiv = (int) n;
192 integer n_ = (integer) n, lda_ = (integer) lda,
193 lwork_ = (integer) lwork, ipiv_[size_ipiv];
194 for (i = 0; i != size_ipiv; ++i) { ipiv_[i] = (integer) ipiv[i]; }
195 dgetri_(&n_, (doublereal*) a,
197 (doublereal*) work, &lwork_,
200 return (lapack_int) info;
203lapack_int LAPACKE_dsytrf_work(
int matrix_layout,
char uplo, lapack_int n,
204 double* a, lapack_int lda, lapack_int* ipiv,
205 double* work, lapack_int lwork )
208 assert(matrix_layout == LAPACK_COL_MAJOR &&
209 sizeof(doublereal) ==
sizeof(
double));
210 if (
sizeof(integer) ==
sizeof(lapack_int))
212 dsytrf_(&uplo, (integer*) &n,
213 (doublereal*) a, (integer*) &lda, (integer*) ipiv,
214 (doublereal*) work, (integer*) &lwork, &info);
218 int i, size_ipiv = (int) n;
219 integer n_ = (integer) n, lda_ = (integer) lda,
220 lwork_ = (integer) lwork, ipiv_[size_ipiv];
222 (doublereal*) a, &lda_, ipiv_,
223 (doublereal*) work, &lwork_, &info);
224 for (i = 0; i != size_ipiv; ++i) { ipiv[i] = (lapack_int) ipiv_[i]; }
226 return (lapack_int) info;
229lapack_int LAPACKE_dsytri_work(
int matrix_layout,
char uplo, lapack_int n,
230 double* a, lapack_int lda,
231 const lapack_int* ipiv,
double* work )
234 assert(matrix_layout == LAPACK_COL_MAJOR &&
235 sizeof(doublereal) ==
sizeof(
double));
236 if (
sizeof(integer) ==
sizeof(lapack_int))
238 dsytri_(&uplo, (integer*) &n,
239 (doublereal*) a, (integer*) &lda, (integer*) ipiv,
240 (doublereal*) work, &info);
244 int i, size_ipiv = (int) n;
245 integer n_ = (integer) n, lda_ = (integer) lda,
247 for (i = 0; i != size_ipiv; ++i) { ipiv_[i] = (integer) ipiv[i]; }
249 (doublereal*) a, &lda_, ipiv_,
250 (doublereal*) work, &info);
252 return (lapack_int) info;
255#elif YAC_LAPACK_INTERFACE_ID == 5
263#define LAPACK_dgels YAC_FC_GLOBAL(dgels,DGELS)
264#define LAPACK_dgesv YAC_FC_GLOBAL(dgesv,DGESV)
265#define LAPACK_dgetrf YAC_FC_GLOBAL(dgetrf,DGETRF)
266#define LAPACK_dgetri YAC_FC_GLOBAL(dgetri,DGETRI)
267#define LAPACK_dsytrf YAC_FC_GLOBAL(dsytrf,DSYTRF)
268#define LAPACK_dsytri YAC_FC_GLOBAL(dsytri,DSYTRI)
271void LAPACK_dgels(
char* trans, lapack_int* m, lapack_int* n, lapack_int* nrhs,
272 double* a, lapack_int* lda,
double* b, lapack_int* ldb,
273 double* work, lapack_int* lwork, lapack_int *info );
275void LAPACK_dgesv( lapack_int* n, lapack_int* nrhs,
double* a, lapack_int* lda,
276 lapack_int* ipiv,
double* b, lapack_int* ldb,
278void LAPACK_dgetrf( lapack_int* m, lapack_int* n,
double* a, lapack_int* lda,
279 lapack_int* ipiv, lapack_int *info );
280void LAPACK_dgetri( lapack_int* n,
double* a, lapack_int* lda,
281 const lapack_int* ipiv,
double* work, lapack_int* lwork,
283void LAPACK_dsytrf(
char* uplo, lapack_int* n,
double* a, lapack_int* lda,
284 lapack_int* ipiv,
double* work, lapack_int* lwork,
286void LAPACK_dsytri(
char* uplo, lapack_int* n,
double* a, lapack_int* lda,
287 const lapack_int* ipiv,
double* work, lapack_int *info );
295lapack_int LAPACKE_dgels_work(
int matrix_layout,
char trans, lapack_int m,
296 lapack_int n, lapack_int nrhs,
double* a,
297 lapack_int lda,
double* b, lapack_int ldb,
298 double* work, lapack_int lwork )
301 assert(matrix_layout == LAPACK_COL_MAJOR);
302 LAPACK_dgels(&trans, &m, &n, &nrhs, a, &lda, b, &ldb, work, &lwork, &info);
303 if( info < 0 ) { info = info - 1; }
308lapack_int LAPACKE_dgesv(
int matrix_layout, lapack_int n, lapack_int nrhs,
309 double* a, lapack_int lda, lapack_int* ipiv,
310 double* b, lapack_int ldb )
313 assert(matrix_layout == LAPACK_COL_MAJOR);
314 LAPACK_dgesv( &n, &nrhs, a, &lda, ipiv, b, &ldb, &info );
315 if( info < 0 ) { info = info - 1; }
319lapack_int LAPACKE_dgetrf(
int matrix_layout, lapack_int m, lapack_int n,
320 double* a, lapack_int lda, lapack_int* ipiv )
323 assert(matrix_layout == LAPACK_COL_MAJOR);
324 LAPACK_dgetrf( &m, &n, a, &lda, ipiv, &info );
325 if( info < 0 ) { info = info - 1; }
329lapack_int LAPACKE_dgetri_work(
int matrix_layout, lapack_int n,
double* a,
330 lapack_int lda,
const lapack_int* ipiv,
331 double* work, lapack_int lwork )
334 assert(matrix_layout == LAPACK_COL_MAJOR);
335 LAPACK_dgetri( &n, a, &lda, ipiv, work, &lwork, &info );
336 if( info < 0 ) { info = info - 1; }
340lapack_int LAPACKE_dsytrf_work(
int matrix_layout,
char uplo, lapack_int n,
341 double* a, lapack_int lda, lapack_int* ipiv,
342 double* work, lapack_int lwork )
345 assert(matrix_layout == LAPACK_COL_MAJOR);
346 LAPACK_dsytrf( &uplo, &n, a, &lda, ipiv, work, &lwork, &info );
347 if( info < 0 ) { info = info - 1; }
351lapack_int LAPACKE_dsytri_work(
int matrix_layout,
char uplo, lapack_int n,
352 double* a, lapack_int lda,
353 const lapack_int* ipiv,
double* work )
356 assert(matrix_layout == LAPACK_COL_MAJOR);
357 LAPACK_dsytri( &uplo, &n, a, &lda, ipiv, work, &info );
358 if( info < 0 ) { info = info - 1; }