YetAnotherCoupler 3.2.0_a
Loading...
Searching...
No Matches
yac_lapack_interface.c
Go to the documentation of this file.
1// Copyright (c) 2024 The YAC Authors
2//
3// SPDX-License-Identifier: BSD-3-Clause
4
6
7#if YAC_LAPACK_INTERFACE_ID == 3 // ATLAS CLAPACK
8
9#include <assert.h>
10
11#if 0
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 )
16{
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);
22}
23#endif
24
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 )
28{
29 assert(matrix_layout == LAPACK_COL_MAJOR);
30 if (sizeof(int) == sizeof(lapack_int))
31 {
32 return (lapack_int) clapack_dgesv(LAPACK_COL_MAJOR, (int) n, (int) nrhs,
33 a, (int) lda, (int*) ipiv, b, (int) ldb);
34 }
35 else
36 {
37 int i, result, ipiv_size = (int) n;
38 int ipiv_[ipiv_size];
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;
43 }
44}
45
46lapack_int LAPACKE_dgetrf( int matrix_layout, lapack_int m, lapack_int n,
47 double* a, lapack_int lda, lapack_int* ipiv )
48{
49 assert(matrix_layout == LAPACK_COL_MAJOR);
50 if (sizeof(int) == sizeof(lapack_int))
51 {
52 return (lapack_int) clapack_dgetrf(LAPACK_COL_MAJOR, (int) m, (int) n,
53 a, (int) lda, (int*) ipiv);
54 }
55 else
56 {
57 int i, result, ipiv_size = (int) (n < m ? n : m);
58 int ipiv_[ipiv_size];
59 result = clapack_dgetrf(LAPACK_COL_MAJOR, (int) m, (int) n,
60 a, (int) lda, ipiv_);
61 for (i = 0; i != ipiv_size; ++i) { ipiv[i] = (lapack_int) ipiv_[i]; }
62 return (lapack_int) result;
63 }
64}
65
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 )
69{
70 assert(matrix_layout == LAPACK_COL_MAJOR);
71 if (sizeof(int) == sizeof(lapack_int))
72 {
73 return (lapack_int) clapack_dgetri(LAPACK_COL_MAJOR, (int) n, a,
74 (int) lda, (int*) ipiv);
75 }
76 else
77 {
78 int i, ipiv_size = (int) n;
79 int ipiv_[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,
82 (int) lda, ipiv_);
83 }
84}
85
86#elif YAC_LAPACK_INTERFACE_ID == 4 // Netlib CLAPACK
87
88#include <assert.h>
89#include <clapack.h>
90
91#if 0
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 )
96{
97 integer info = 0;
98 assert(matrix_layout == LAPACK_COL_MAJOR &&
99 sizeof(doublereal) == sizeof(double));
100 if (sizeof(integer) == sizeof(lapack_int))
101 {
102 dgels_(&trans, (integer*) &m, (integer*) &n, (integer*) &nrhs,
103 (doublereal*) a, (integer*) &lda,
104 (doublereal*) b, (integer*) &ldb,
105 (doublereal*) work, (integer*) &lwork,
106 &info);
107 }
108 else
109 {
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_,
117 &info);
118 }
119 return (lapack_int) info;
120}
121#endif
122
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 )
126{
127 integer info = 0;
128 assert(matrix_layout == LAPACK_COL_MAJOR &&
129 sizeof(doublereal) == sizeof(double));
130 if (sizeof(integer) == sizeof(lapack_int))
131 {
132 dgesv_((integer*) &n, (integer*) &nrhs,
133 (doublereal*) a, (integer*) &lda , (integer*) ipiv,
134 (doublereal*) b, (integer*) &ldb, &info);
135 }
136 else
137 {
138 int i, ipiv_size = (int) n;
139 integer n_ = (integer) n, nrhs_ = (integer) nrhs,
140 lda_ = (integer) lda, ldb_ = (integer) ldb,
141 ipiv_[n];
142 dgesv_(&n_, &nrhs_,
143 (doublereal*) a, &lda_, ipiv_,
144 (doublereal*) b, &ldb_, &info);
145 for (i = 0; i != ipiv_size; ++i) { ipiv[i] = (lapack_int) ipiv_[i]; }
146 }
147 return (lapack_int) info;
148}
149
150lapack_int LAPACKE_dgetrf( int matrix_layout, lapack_int m, lapack_int n,
151 double* a, lapack_int lda, lapack_int* ipiv )
152{
153 integer info = 0;
154 assert(matrix_layout == LAPACK_COL_MAJOR &&
155 sizeof(doublereal) == sizeof(double));
156 if (sizeof(integer) == sizeof(lapack_int))
157 {
158 dgetrf_((integer*) &m, (integer*) &n,
159 (doublereal*) a, (integer*) &lda , (integer*) ipiv,
160 &info);
161 }
162 else
163 {
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];
167 dgetrf_(&m_, &n_,
168 (doublereal*) a, &lda_, ipiv_,
169 &info);
170 for (i = 0; i != ipiv_size; ++i) {ipiv[i] = (lapack_int) ipiv_[i]; }
171 }
172 return (lapack_int) info;
173}
174
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 )
178{
179 integer info = 0;
180 assert(matrix_layout == LAPACK_COL_MAJOR &&
181 sizeof(doublereal) == sizeof(double));
182 if (sizeof(integer) == sizeof(lapack_int))
183 {
184 dgetri_((integer*) &n, (doublereal*) a,
185 (integer*) &lda , (integer*) ipiv,
186 (doublereal*) work, (integer*) &lwork,
187 &info);
188 }
189 else
190 {
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,
196 &lda_, ipiv_,
197 (doublereal*) work, &lwork_,
198 &info);
199 }
200 return (lapack_int) info;
201}
202
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 )
206{
207 integer info = 0;
208 assert(matrix_layout == LAPACK_COL_MAJOR &&
209 sizeof(doublereal) == sizeof(double));
210 if (sizeof(integer) == sizeof(lapack_int))
211 {
212 dsytrf_(&uplo, (integer*) &n,
213 (doublereal*) a, (integer*) &lda, (integer*) ipiv,
214 (doublereal*) work, (integer*) &lwork, &info);
215 }
216 else
217 {
218 int i, size_ipiv = (int) n;
219 integer n_ = (integer) n, lda_ = (integer) lda,
220 lwork_ = (integer) lwork, ipiv_[size_ipiv];
221 dsytrf_(&uplo, &n_,
222 (doublereal*) a, &lda_, ipiv_,
223 (doublereal*) work, &lwork_, &info);
224 for (i = 0; i != size_ipiv; ++i) { ipiv[i] = (lapack_int) ipiv_[i]; }
225 }
226 return (lapack_int) info;
227}
228
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 )
232{
233 integer info = 0;
234 assert(matrix_layout == LAPACK_COL_MAJOR &&
235 sizeof(doublereal) == sizeof(double));
236 if (sizeof(integer) == sizeof(lapack_int))
237 {
238 dsytri_(&uplo, (integer*) &n,
239 (doublereal*) a, (integer*) &lda, (integer*) ipiv,
240 (doublereal*) work, &info);
241 }
242 else
243 {
244 int i, size_ipiv = (int) n;
245 integer n_ = (integer) n, lda_ = (integer) lda,
246 ipiv_[size_ipiv];
247 for (i = 0; i != size_ipiv; ++i) { ipiv_[i] = (integer) ipiv[i]; }
248 dsytri_(&uplo, &n_,
249 (doublereal*) a, &lda_, ipiv_,
250 (doublereal*) work, &info);
251 }
252 return (lapack_int) info;
253}
254
255#elif YAC_LAPACK_INTERFACE_ID == 5 // Fortran LAPACK
256
257#include <assert.h>
258
259#ifdef __cplusplus
260extern "C" {
261#endif
262
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)
269
270#if 0
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 );
274#endif
275void LAPACK_dgesv( lapack_int* n, lapack_int* nrhs, double* a, lapack_int* lda,
276 lapack_int* ipiv, double* b, lapack_int* ldb,
277 lapack_int *info );
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,
282 lapack_int *info );
283void LAPACK_dsytrf( char* uplo, lapack_int* n, double* a, lapack_int* lda,
284 lapack_int* ipiv, double* work, lapack_int* lwork,
285 lapack_int *info );
286void LAPACK_dsytri( char* uplo, lapack_int* n, double* a, lapack_int* lda,
287 const lapack_int* ipiv, double* work, lapack_int *info );
288
289#ifdef __cplusplus
290}
291#endif
292
293
294#if 0
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 )
299{
300 lapack_int info = 0;
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; }
304 return info;
305}
306#endif
307
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 )
311{
312 lapack_int info = 0;
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; }
316 return info;
317}
318
319lapack_int LAPACKE_dgetrf( int matrix_layout, lapack_int m, lapack_int n,
320 double* a, lapack_int lda, lapack_int* ipiv )
321{
322 lapack_int info = 0;
323 assert(matrix_layout == LAPACK_COL_MAJOR);
324 LAPACK_dgetrf( &m, &n, a, &lda, ipiv, &info );
325 if( info < 0 ) { info = info - 1; }
326 return info;
327}
328
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 )
332{
333 lapack_int info = 0;
334 assert(matrix_layout == LAPACK_COL_MAJOR);
335 LAPACK_dgetri( &n, a, &lda, ipiv, work, &lwork, &info );
336 if( info < 0 ) { info = info - 1; }
337 return info;
338}
339
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 )
343{
344 lapack_int info = 0;
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; }
348 return info;
349}
350
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 )
354{
355 lapack_int info = 0;
356 assert(matrix_layout == LAPACK_COL_MAJOR);
357 LAPACK_dsytri( &uplo, &n, a, &lda, ipiv, work, &info );
358 if( info < 0 ) { info = info - 1; }
359 return info;
360}
361
362#endif