LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
c_cblas3.c
Go to the documentation of this file.
1/*
2 * Written by D.P. Manley, Digital Equipment Corporation.
3 * Prefixed "C_" to BLAS routines and their declarations.
4 *
5 * Modified by T. H. Do, 4/15/98, SGI/CRAY Research.
6 */
7#include <stdlib.h>
8#include "cblas.h"
9#include "cblas_test.h"
10#define TEST_COL_MJR 0
11#define TEST_ROW_MJR 1
12#define UNDEFINED -1
13
14void F77_cgemm(CBLAS_INT *layout, char *transpa, char *transpb, CBLAS_INT *m, CBLAS_INT *n,
19 , FORTRAN_STRLEN transpa_len, FORTRAN_STRLEN transpb_len
20#endif
21) {
22
23 CBLAS_TEST_COMPLEX *A, *B, *C;
24 CBLAS_INT i,j,LDA, LDB, LDC;
25 CBLAS_TRANSPOSE transa, transb;
26
27 get_transpose_type(transpa, &transa);
28 get_transpose_type(transpb, &transb);
29
30 if (*layout == TEST_ROW_MJR) {
31 if (transa == CblasNoTrans) {
32 LDA = *k+1;
33 A=(CBLAS_TEST_COMPLEX*)malloc((*m)*LDA*sizeof(CBLAS_TEST_COMPLEX));
34 for( i=0; i<*m; i++ )
35 for( j=0; j<*k; j++ ) {
36 A[i*LDA+j].real=a[j*(*lda)+i].real;
37 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
38 }
39 }
40 else {
41 LDA = *m+1;
42 A=(CBLAS_TEST_COMPLEX* )malloc(LDA*(*k)*sizeof(CBLAS_TEST_COMPLEX));
43 for( i=0; i<*k; i++ )
44 for( j=0; j<*m; j++ ) {
45 A[i*LDA+j].real=a[j*(*lda)+i].real;
46 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
47 }
48 }
49
50 if (transb == CblasNoTrans) {
51 LDB = *n+1;
52 B=(CBLAS_TEST_COMPLEX* )malloc((*k)*LDB*sizeof(CBLAS_TEST_COMPLEX) );
53 for( i=0; i<*k; i++ )
54 for( j=0; j<*n; j++ ) {
55 B[i*LDB+j].real=b[j*(*ldb)+i].real;
56 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
57 }
58 }
59 else {
60 LDB = *k+1;
61 B=(CBLAS_TEST_COMPLEX* )malloc(LDB*(*n)*sizeof(CBLAS_TEST_COMPLEX));
62 for( i=0; i<*n; i++ )
63 for( j=0; j<*k; j++ ) {
64 B[i*LDB+j].real=b[j*(*ldb)+i].real;
65 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
66 }
67 }
68
69 LDC = *n+1;
70 C=(CBLAS_TEST_COMPLEX* )malloc((*m)*LDC*sizeof(CBLAS_TEST_COMPLEX));
71 for( j=0; j<*n; j++ )
72 for( i=0; i<*m; i++ ) {
73 C[i*LDC+j].real=c[j*(*ldc)+i].real;
74 C[i*LDC+j].imag=c[j*(*ldc)+i].imag;
75 }
76 cblas_cgemm( CblasRowMajor, transa, transb, *m, *n, *k, alpha, A, LDA,
77 B, LDB, beta, C, LDC );
78 for( j=0; j<*n; j++ )
79 for( i=0; i<*m; i++ ) {
80 c[j*(*ldc)+i].real=C[i*LDC+j].real;
81 c[j*(*ldc)+i].imag=C[i*LDC+j].imag;
82 }
83 free(A);
84 free(B);
85 free(C);
86 }
87 else if (*layout == TEST_COL_MJR)
88 cblas_cgemm( CblasColMajor, transa, transb, *m, *n, *k, alpha, a, *lda,
89 b, *ldb, beta, c, *ldc );
90 else
91 cblas_cgemm( UNDEFINED, transa, transb, *m, *n, *k, alpha, a, *lda,
92 b, *ldb, beta, c, *ldc );
93}
94
95void F77_cgemmtr(CBLAS_INT *layout, char *uplop, char *transpa, char *transpb, CBLAS_INT *n,
98 CBLAS_TEST_COMPLEX *c, CBLAS_INT *ldc ) {
99
100 CBLAS_TEST_COMPLEX *A, *B, *C;
101 CBLAS_INT i,j,LDA, LDB, LDC;
102 CBLAS_TRANSPOSE transa, transb;
103 CBLAS_UPLO uplo;
104
105 get_transpose_type(transpa, &transa);
106 get_transpose_type(transpb, &transb);
107 get_uplo_type(uplop, &uplo);
108
109 if (*layout == TEST_ROW_MJR) {
110 if (transa == CblasNoTrans) {
111 LDA = *k+1;
112 A=(CBLAS_TEST_COMPLEX*)malloc((*n)*LDA*sizeof(CBLAS_TEST_COMPLEX));
113 for( i=0; i<*n; i++ )
114 for( j=0; j<*k; j++ ) {
115 A[i*LDA+j].real=a[j*(*lda)+i].real;
116 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
117 }
118 }
119 else {
120 LDA = *n+1;
121 A=(CBLAS_TEST_COMPLEX* )malloc(LDA*(*k)*sizeof(CBLAS_TEST_COMPLEX));
122 for( i=0; i<*k; i++ )
123 for( j=0; j<*n; j++ ) {
124 A[i*LDA+j].real=a[j*(*lda)+i].real;
125 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
126 }
127 }
128
129 if (transb == CblasNoTrans) {
130 LDB = *n+1;
131 B=(CBLAS_TEST_COMPLEX* )malloc((*k)*LDB*sizeof(CBLAS_TEST_COMPLEX) );
132 for( i=0; i<*k; i++ )
133 for( j=0; j<*n; j++ ) {
134 B[i*LDB+j].real=b[j*(*ldb)+i].real;
135 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
136 }
137 }
138 else {
139 LDB = *k+1;
140 B=(CBLAS_TEST_COMPLEX* )malloc(LDB*(*n)*sizeof(CBLAS_TEST_COMPLEX));
141 for( i=0; i<*n; i++ )
142 for( j=0; j<*k; j++ ) {
143 B[i*LDB+j].real=b[j*(*ldb)+i].real;
144 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
145 }
146 }
147
148 LDC = *n+1;
149 C=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDC*sizeof(CBLAS_TEST_COMPLEX));
150 for( j=0; j<*n; j++ )
151 for( i=0; i<*n; i++ ) {
152 C[i*LDC+j].real=c[j*(*ldc)+i].real;
153 C[i*LDC+j].imag=c[j*(*ldc)+i].imag;
154 }
155 cblas_cgemmtr( CblasRowMajor, uplo, transa, transb, *n, *k, alpha, A, LDA,
156 B, LDB, beta, C, LDC );
157 for( j=0; j<*n; j++ )
158 for( i=0; i<*n; i++ ) {
159 c[j*(*ldc)+i].real=C[i*LDC+j].real;
160 c[j*(*ldc)+i].imag=C[i*LDC+j].imag;
161 }
162 free(A);
163 free(B);
164 free(C);
165 }
166 else if (*layout == TEST_COL_MJR)
167 cblas_cgemmtr( CblasColMajor, uplo, transa, transb, *n, *k, alpha, a, *lda,
168 b, *ldb, beta, c, *ldc );
169 else
170 cblas_cgemmtr( UNDEFINED, uplo, transa, transb, *n, *k, alpha, a, *lda,
171 b, *ldb, beta, c, *ldc );
172}
173
174
175void F77_chemm(CBLAS_INT *layout, char *rtlf, char *uplow, CBLAS_INT *m, CBLAS_INT *n,
180 , FORTRAN_STRLEN rtlf_len, FORTRAN_STRLEN uplow_len
181#endif
182) {
183
184 CBLAS_TEST_COMPLEX *A, *B, *C;
185 CBLAS_INT i,j,LDA, LDB, LDC;
186 CBLAS_UPLO uplo;
187 CBLAS_SIDE side;
188
189 get_uplo_type(uplow,&uplo);
190 get_side_type(rtlf,&side);
191
192 if (*layout == TEST_ROW_MJR) {
193 if (side == CblasLeft) {
194 LDA = *m+1;
195 A= (CBLAS_TEST_COMPLEX* )malloc((*m)*LDA*sizeof(CBLAS_TEST_COMPLEX));
196 for( i=0; i<*m; i++ )
197 for( j=0; j<*m; j++ ) {
198 A[i*LDA+j].real=a[j*(*lda)+i].real;
199 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
200 }
201 }
202 else{
203 LDA = *n+1;
204 A=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDA*sizeof(CBLAS_TEST_COMPLEX ) );
205 for( i=0; i<*n; i++ )
206 for( j=0; j<*n; j++ ) {
207 A[i*LDA+j].real=a[j*(*lda)+i].real;
208 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
209 }
210 }
211 LDB = *n+1;
212 B=(CBLAS_TEST_COMPLEX* )malloc( (*m)*LDB*sizeof(CBLAS_TEST_COMPLEX ) );
213 for( i=0; i<*m; i++ )
214 for( j=0; j<*n; j++ ) {
215 B[i*LDB+j].real=b[j*(*ldb)+i].real;
216 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
217 }
218 LDC = *n+1;
219 C=(CBLAS_TEST_COMPLEX* )malloc((*m)*LDC*sizeof(CBLAS_TEST_COMPLEX ) );
220 for( j=0; j<*n; j++ )
221 for( i=0; i<*m; i++ ) {
222 C[i*LDC+j].real=c[j*(*ldc)+i].real;
223 C[i*LDC+j].imag=c[j*(*ldc)+i].imag;
224 }
225 cblas_chemm( CblasRowMajor, side, uplo, *m, *n, alpha, A, LDA, B, LDB,
226 beta, C, LDC );
227 for( j=0; j<*n; j++ )
228 for( i=0; i<*m; i++ ) {
229 c[j*(*ldc)+i].real=C[i*LDC+j].real;
230 c[j*(*ldc)+i].imag=C[i*LDC+j].imag;
231 }
232 free(A);
233 free(B);
234 free(C);
235 }
236 else if (*layout == TEST_COL_MJR)
237 cblas_chemm( CblasColMajor, side, uplo, *m, *n, alpha, a, *lda, b, *ldb,
238 beta, c, *ldc );
239 else
240 cblas_chemm( UNDEFINED, side, uplo, *m, *n, alpha, a, *lda, b, *ldb,
241 beta, c, *ldc );
242}
243void F77_csymm(CBLAS_INT *layout, char *rtlf, char *uplow, CBLAS_INT *m, CBLAS_INT *n,
248 , FORTRAN_STRLEN rtlf_len, FORTRAN_STRLEN uplow_len
249#endif
250) {
251
252 CBLAS_TEST_COMPLEX *A, *B, *C;
253 CBLAS_INT i,j,LDA, LDB, LDC;
254 CBLAS_UPLO uplo;
255 CBLAS_SIDE side;
256
257 get_uplo_type(uplow,&uplo);
258 get_side_type(rtlf,&side);
259
260 if (*layout == TEST_ROW_MJR) {
261 if (side == CblasLeft) {
262 LDA = *m+1;
263 A=(CBLAS_TEST_COMPLEX* )malloc((*m)*LDA*sizeof(CBLAS_TEST_COMPLEX));
264 for( i=0; i<*m; i++ )
265 for( j=0; j<*m; j++ )
266 A[i*LDA+j]=a[j*(*lda)+i];
267 }
268 else{
269 LDA = *n+1;
270 A=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDA*sizeof(CBLAS_TEST_COMPLEX ) );
271 for( i=0; i<*n; i++ )
272 for( j=0; j<*n; j++ )
273 A[i*LDA+j]=a[j*(*lda)+i];
274 }
275 LDB = *n+1;
276 B=(CBLAS_TEST_COMPLEX* )malloc((*m)*LDB*sizeof(CBLAS_TEST_COMPLEX ));
277 for( i=0; i<*m; i++ )
278 for( j=0; j<*n; j++ )
279 B[i*LDB+j]=b[j*(*ldb)+i];
280 LDC = *n+1;
281 C=(CBLAS_TEST_COMPLEX* )malloc((*m)*LDC*sizeof(CBLAS_TEST_COMPLEX));
282 for( j=0; j<*n; j++ )
283 for( i=0; i<*m; i++ )
284 C[i*LDC+j]=c[j*(*ldc)+i];
285 cblas_csymm( CblasRowMajor, side, uplo, *m, *n, alpha, A, LDA, B, LDB,
286 beta, C, LDC );
287 for( j=0; j<*n; j++ )
288 for( i=0; i<*m; i++ )
289 c[j*(*ldc)+i]=C[i*LDC+j];
290 free(A);
291 free(B);
292 free(C);
293 }
294 else if (*layout == TEST_COL_MJR)
295 cblas_csymm( CblasColMajor, side, uplo, *m, *n, alpha, a, *lda, b, *ldb,
296 beta, c, *ldc );
297 else
298 cblas_csymm( UNDEFINED, side, uplo, *m, *n, alpha, a, *lda, b, *ldb,
299 beta, c, *ldc );
300}
301
302void F77_cherk(CBLAS_INT *layout, char *uplow, char *transp, CBLAS_INT *n, CBLAS_INT *k,
303 float *alpha, CBLAS_TEST_COMPLEX *a, CBLAS_INT *lda,
304 float *beta, CBLAS_TEST_COMPLEX *c, CBLAS_INT *ldc
306 , FORTRAN_STRLEN uplow_len, FORTRAN_STRLEN transp_len
307#endif
308) {
309
310 CBLAS_INT i,j,LDA,LDC;
311 CBLAS_TEST_COMPLEX *A, *C;
312 CBLAS_UPLO uplo;
313 CBLAS_TRANSPOSE trans;
314
315 get_uplo_type(uplow,&uplo);
316 get_transpose_type(transp,&trans);
317
318 if (*layout == TEST_ROW_MJR) {
319 if (trans == CblasNoTrans) {
320 LDA = *k+1;
321 A=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDA*sizeof(CBLAS_TEST_COMPLEX ) );
322 for( i=0; i<*n; i++ )
323 for( j=0; j<*k; j++ ) {
324 A[i*LDA+j].real=a[j*(*lda)+i].real;
325 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
326 }
327 }
328 else{
329 LDA = *n+1;
330 A=(CBLAS_TEST_COMPLEX* )malloc((*k)*LDA*sizeof(CBLAS_TEST_COMPLEX ) );
331 for( i=0; i<*k; i++ )
332 for( j=0; j<*n; j++ ) {
333 A[i*LDA+j].real=a[j*(*lda)+i].real;
334 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
335 }
336 }
337 LDC = *n+1;
338 C=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDC*sizeof(CBLAS_TEST_COMPLEX ) );
339 for( i=0; i<*n; i++ )
340 for( j=0; j<*n; j++ ) {
341 C[i*LDC+j].real=c[j*(*ldc)+i].real;
342 C[i*LDC+j].imag=c[j*(*ldc)+i].imag;
343 }
344 cblas_cherk(CblasRowMajor, uplo, trans, *n, *k, *alpha, A, LDA, *beta,
345 C, LDC );
346 for( j=0; j<*n; j++ )
347 for( i=0; i<*n; i++ ) {
348 c[j*(*ldc)+i].real=C[i*LDC+j].real;
349 c[j*(*ldc)+i].imag=C[i*LDC+j].imag;
350 }
351 free(A);
352 free(C);
353 }
354 else if (*layout == TEST_COL_MJR)
355 cblas_cherk(CblasColMajor, uplo, trans, *n, *k, *alpha, a, *lda, *beta,
356 c, *ldc );
357 else
358 cblas_cherk(UNDEFINED, uplo, trans, *n, *k, *alpha, a, *lda, *beta,
359 c, *ldc );
360}
361
362void F77_csyrk(CBLAS_INT *layout, char *uplow, char *transp, CBLAS_INT *n, CBLAS_INT *k,
366 , FORTRAN_STRLEN uplow_len, FORTRAN_STRLEN transp_len
367#endif
368) {
369
370 CBLAS_INT i,j,LDA,LDC;
371 CBLAS_TEST_COMPLEX *A, *C;
372 CBLAS_UPLO uplo;
373 CBLAS_TRANSPOSE trans;
374
375 get_uplo_type(uplow,&uplo);
376 get_transpose_type(transp,&trans);
377
378 if (*layout == TEST_ROW_MJR) {
379 if (trans == CblasNoTrans) {
380 LDA = *k+1;
381 A=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDA*sizeof(CBLAS_TEST_COMPLEX));
382 for( i=0; i<*n; i++ )
383 for( j=0; j<*k; j++ ) {
384 A[i*LDA+j].real=a[j*(*lda)+i].real;
385 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
386 }
387 }
388 else{
389 LDA = *n+1;
390 A=(CBLAS_TEST_COMPLEX* )malloc((*k)*LDA*sizeof(CBLAS_TEST_COMPLEX ) );
391 for( i=0; i<*k; i++ )
392 for( j=0; j<*n; j++ ) {
393 A[i*LDA+j].real=a[j*(*lda)+i].real;
394 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
395 }
396 }
397 LDC = *n+1;
398 C=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDC*sizeof(CBLAS_TEST_COMPLEX ) );
399 for( i=0; i<*n; i++ )
400 for( j=0; j<*n; j++ ) {
401 C[i*LDC+j].real=c[j*(*ldc)+i].real;
402 C[i*LDC+j].imag=c[j*(*ldc)+i].imag;
403 }
404 cblas_csyrk(CblasRowMajor, uplo, trans, *n, *k, alpha, A, LDA, beta,
405 C, LDC );
406 for( j=0; j<*n; j++ )
407 for( i=0; i<*n; i++ ) {
408 c[j*(*ldc)+i].real=C[i*LDC+j].real;
409 c[j*(*ldc)+i].imag=C[i*LDC+j].imag;
410 }
411 free(A);
412 free(C);
413 }
414 else if (*layout == TEST_COL_MJR)
415 cblas_csyrk(CblasColMajor, uplo, trans, *n, *k, alpha, a, *lda, beta,
416 c, *ldc );
417 else
418 cblas_csyrk(UNDEFINED, uplo, trans, *n, *k, alpha, a, *lda, beta,
419 c, *ldc );
420}
421void F77_cher2k(CBLAS_INT *layout, char *uplow, char *transp, CBLAS_INT *n, CBLAS_INT *k,
423 CBLAS_TEST_COMPLEX *b, CBLAS_INT *ldb, float *beta,
426 , FORTRAN_STRLEN uplow_len, FORTRAN_STRLEN transp_len
427#endif
428) {
429 CBLAS_INT i,j,LDA,LDB,LDC;
430 CBLAS_TEST_COMPLEX *A, *B, *C;
431 CBLAS_UPLO uplo;
432 CBLAS_TRANSPOSE trans;
433
434 get_uplo_type(uplow,&uplo);
435 get_transpose_type(transp,&trans);
436
437 if (*layout == TEST_ROW_MJR) {
438 if (trans == CblasNoTrans) {
439 LDA = *k+1;
440 LDB = *k+1;
441 A=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDA*sizeof(CBLAS_TEST_COMPLEX ));
442 B=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDB*sizeof(CBLAS_TEST_COMPLEX ));
443 for( i=0; i<*n; i++ )
444 for( j=0; j<*k; j++ ) {
445 A[i*LDA+j].real=a[j*(*lda)+i].real;
446 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
447 B[i*LDB+j].real=b[j*(*ldb)+i].real;
448 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
449 }
450 }
451 else {
452 LDA = *n+1;
453 LDB = *n+1;
454 A=(CBLAS_TEST_COMPLEX* )malloc( LDA*(*k)*sizeof(CBLAS_TEST_COMPLEX ) );
455 B=(CBLAS_TEST_COMPLEX* )malloc( LDB*(*k)*sizeof(CBLAS_TEST_COMPLEX ) );
456 for( i=0; i<*k; i++ )
457 for( j=0; j<*n; j++ ){
458 A[i*LDA+j].real=a[j*(*lda)+i].real;
459 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
460 B[i*LDB+j].real=b[j*(*ldb)+i].real;
461 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
462 }
463 }
464 LDC = *n+1;
465 C=(CBLAS_TEST_COMPLEX* )malloc( (*n)*LDC*sizeof(CBLAS_TEST_COMPLEX ) );
466 for( i=0; i<*n; i++ )
467 for( j=0; j<*n; j++ ) {
468 C[i*LDC+j].real=c[j*(*ldc)+i].real;
469 C[i*LDC+j].imag=c[j*(*ldc)+i].imag;
470 }
471 cblas_cher2k(CblasRowMajor, uplo, trans, *n, *k, alpha, A, LDA,
472 B, LDB, *beta, C, LDC );
473 for( j=0; j<*n; j++ )
474 for( i=0; i<*n; i++ ) {
475 c[j*(*ldc)+i].real=C[i*LDC+j].real;
476 c[j*(*ldc)+i].imag=C[i*LDC+j].imag;
477 }
478 free(A);
479 free(B);
480 free(C);
481 }
482 else if (*layout == TEST_COL_MJR)
483 cblas_cher2k(CblasColMajor, uplo, trans, *n, *k, alpha, a, *lda,
484 b, *ldb, *beta, c, *ldc );
485 else
486 cblas_cher2k(UNDEFINED, uplo, trans, *n, *k, alpha, a, *lda,
487 b, *ldb, *beta, c, *ldc );
488}
489void F77_csyr2k(CBLAS_INT *layout, char *uplow, char *transp, CBLAS_INT *n, CBLAS_INT *k,
494 , FORTRAN_STRLEN uplow_len, FORTRAN_STRLEN transp_len
495#endif
496) {
497 CBLAS_INT i,j,LDA,LDB,LDC;
498 CBLAS_TEST_COMPLEX *A, *B, *C;
499 CBLAS_UPLO uplo;
500 CBLAS_TRANSPOSE trans;
501
502 get_uplo_type(uplow,&uplo);
503 get_transpose_type(transp,&trans);
504
505 if (*layout == TEST_ROW_MJR) {
506 if (trans == CblasNoTrans) {
507 LDA = *k+1;
508 LDB = *k+1;
509 A=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDA*sizeof(CBLAS_TEST_COMPLEX));
510 B=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDB*sizeof(CBLAS_TEST_COMPLEX));
511 for( i=0; i<*n; i++ )
512 for( j=0; j<*k; j++ ) {
513 A[i*LDA+j].real=a[j*(*lda)+i].real;
514 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
515 B[i*LDB+j].real=b[j*(*ldb)+i].real;
516 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
517 }
518 }
519 else {
520 LDA = *n+1;
521 LDB = *n+1;
522 A=(CBLAS_TEST_COMPLEX* )malloc(LDA*(*k)*sizeof(CBLAS_TEST_COMPLEX));
523 B=(CBLAS_TEST_COMPLEX* )malloc(LDB*(*k)*sizeof(CBLAS_TEST_COMPLEX));
524 for( i=0; i<*k; i++ )
525 for( j=0; j<*n; j++ ){
526 A[i*LDA+j].real=a[j*(*lda)+i].real;
527 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
528 B[i*LDB+j].real=b[j*(*ldb)+i].real;
529 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
530 }
531 }
532 LDC = *n+1;
533 C=(CBLAS_TEST_COMPLEX* )malloc( (*n)*LDC*sizeof(CBLAS_TEST_COMPLEX));
534 for( i=0; i<*n; i++ )
535 for( j=0; j<*n; j++ ) {
536 C[i*LDC+j].real=c[j*(*ldc)+i].real;
537 C[i*LDC+j].imag=c[j*(*ldc)+i].imag;
538 }
539 cblas_csyr2k(CblasRowMajor, uplo, trans, *n, *k, alpha, A, LDA,
540 B, LDB, beta, C, LDC );
541 for( j=0; j<*n; j++ )
542 for( i=0; i<*n; i++ ) {
543 c[j*(*ldc)+i].real=C[i*LDC+j].real;
544 c[j*(*ldc)+i].imag=C[i*LDC+j].imag;
545 }
546 free(A);
547 free(B);
548 free(C);
549 }
550 else if (*layout == TEST_COL_MJR)
551 cblas_csyr2k(CblasColMajor, uplo, trans, *n, *k, alpha, a, *lda,
552 b, *ldb, beta, c, *ldc );
553 else
554 cblas_csyr2k(UNDEFINED, uplo, trans, *n, *k, alpha, a, *lda,
555 b, *ldb, beta, c, *ldc );
556}
557void F77_ctrmm(CBLAS_INT *layout, char *rtlf, char *uplow, char *transp, char *diagn,
561 , FORTRAN_STRLEN rtlf_len, FORTRAN_STRLEN uplow_len, FORTRAN_STRLEN transp_len, FORTRAN_STRLEN diagn_len
562#endif
563) {
564 CBLAS_INT i,j,LDA,LDB;
565 CBLAS_TEST_COMPLEX *A, *B;
566 CBLAS_SIDE side;
567 CBLAS_DIAG diag;
568 CBLAS_UPLO uplo;
569 CBLAS_TRANSPOSE trans;
570
571 get_uplo_type(uplow,&uplo);
572 get_transpose_type(transp,&trans);
573 get_diag_type(diagn,&diag);
574 get_side_type(rtlf,&side);
575
576 if (*layout == TEST_ROW_MJR) {
577 if (side == CblasLeft) {
578 LDA = *m+1;
579 A=(CBLAS_TEST_COMPLEX* )malloc((*m)*LDA*sizeof(CBLAS_TEST_COMPLEX));
580 for( i=0; i<*m; i++ )
581 for( j=0; j<*m; j++ ) {
582 A[i*LDA+j].real=a[j*(*lda)+i].real;
583 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
584 }
585 }
586 else{
587 LDA = *n+1;
588 A=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDA*sizeof(CBLAS_TEST_COMPLEX));
589 for( i=0; i<*n; i++ )
590 for( j=0; j<*n; j++ ) {
591 A[i*LDA+j].real=a[j*(*lda)+i].real;
592 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
593 }
594 }
595 LDB = *n+1;
596 B=(CBLAS_TEST_COMPLEX* )malloc((*m)*LDB*sizeof(CBLAS_TEST_COMPLEX));
597 for( i=0; i<*m; i++ )
598 for( j=0; j<*n; j++ ) {
599 B[i*LDB+j].real=b[j*(*ldb)+i].real;
600 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
601 }
602 cblas_ctrmm(CblasRowMajor, side, uplo, trans, diag, *m, *n, alpha,
603 A, LDA, B, LDB );
604 for( j=0; j<*n; j++ )
605 for( i=0; i<*m; i++ ) {
606 b[j*(*ldb)+i].real=B[i*LDB+j].real;
607 b[j*(*ldb)+i].imag=B[i*LDB+j].imag;
608 }
609 free(A);
610 free(B);
611 }
612 else if (*layout == TEST_COL_MJR)
613 cblas_ctrmm(CblasColMajor, side, uplo, trans, diag, *m, *n, alpha,
614 a, *lda, b, *ldb);
615 else
616 cblas_ctrmm(UNDEFINED, side, uplo, trans, diag, *m, *n, alpha,
617 a, *lda, b, *ldb);
618}
619
620void F77_ctrsm(CBLAS_INT *layout, char *rtlf, char *uplow, char *transp, char *diagn,
624 , FORTRAN_STRLEN rtlf_len, FORTRAN_STRLEN uplow_len, FORTRAN_STRLEN transp_len, FORTRAN_STRLEN diagn_len
625#endif
626) {
627 CBLAS_INT i,j,LDA,LDB;
628 CBLAS_TEST_COMPLEX *A, *B;
629 CBLAS_SIDE side;
630 CBLAS_DIAG diag;
631 CBLAS_UPLO uplo;
632 CBLAS_TRANSPOSE trans;
633
634 get_uplo_type(uplow,&uplo);
635 get_transpose_type(transp,&trans);
636 get_diag_type(diagn,&diag);
637 get_side_type(rtlf,&side);
638
639 if (*layout == TEST_ROW_MJR) {
640 if (side == CblasLeft) {
641 LDA = *m+1;
642 A=(CBLAS_TEST_COMPLEX* )malloc( (*m)*LDA*sizeof(CBLAS_TEST_COMPLEX ) );
643 for( i=0; i<*m; i++ )
644 for( j=0; j<*m; j++ ) {
645 A[i*LDA+j].real=a[j*(*lda)+i].real;
646 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
647 }
648 }
649 else{
650 LDA = *n+1;
651 A=(CBLAS_TEST_COMPLEX* )malloc((*n)*LDA*sizeof(CBLAS_TEST_COMPLEX));
652 for( i=0; i<*n; i++ )
653 for( j=0; j<*n; j++ ) {
654 A[i*LDA+j].real=a[j*(*lda)+i].real;
655 A[i*LDA+j].imag=a[j*(*lda)+i].imag;
656 }
657 }
658 LDB = *n+1;
659 B=(CBLAS_TEST_COMPLEX* )malloc((*m)*LDB*sizeof(CBLAS_TEST_COMPLEX));
660 for( i=0; i<*m; i++ )
661 for( j=0; j<*n; j++ ) {
662 B[i*LDB+j].real=b[j*(*ldb)+i].real;
663 B[i*LDB+j].imag=b[j*(*ldb)+i].imag;
664 }
665 cblas_ctrsm(CblasRowMajor, side, uplo, trans, diag, *m, *n, alpha,
666 A, LDA, B, LDB );
667 for( j=0; j<*n; j++ )
668 for( i=0; i<*m; i++ ) {
669 b[j*(*ldb)+i].real=B[i*LDB+j].real;
670 b[j*(*ldb)+i].imag=B[i*LDB+j].imag;
671 }
672 free(A);
673 free(B);
674 }
675 else if (*layout == TEST_COL_MJR)
676 cblas_ctrsm(CblasColMajor, side, uplo, trans, diag, *m, *n, alpha,
677 a, *lda, b, *ldb);
678 else
679 cblas_ctrsm(UNDEFINED, side, uplo, trans, diag, *m, *n, alpha,
680 a, *lda, b, *ldb);
681}
#define UNDEFINED
Definition c_cblas3.c:12
#define TEST_ROW_MJR
Definition c_cblas3.c:11
#define TEST_COL_MJR
Definition c_cblas3.c:10
void cblas_cherk(CBLAS_LAYOUT layout, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE Trans, const CBLAS_INT N, const CBLAS_INT K, const float alpha, const void *A, const CBLAS_INT lda, const float beta, void *C, const CBLAS_INT ldc)
Definition cblas_cherk.c:12
void cblas_cher2k(CBLAS_LAYOUT layout, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE Trans, const CBLAS_INT N, const CBLAS_INT K, const void *alpha, const void *A, const CBLAS_INT lda, const void *B, const CBLAS_INT ldb, const float beta, void *C, const CBLAS_INT ldc)
CBLAS_UPLO
Definition cblas.h:41
void cblas_csyr2k(CBLAS_LAYOUT layout, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE Trans, const CBLAS_INT N, const CBLAS_INT K, const void *alpha, const void *A, const CBLAS_INT lda, const void *B, const CBLAS_INT ldb, const void *beta, void *C, const CBLAS_INT ldc)
void cblas_csymm(CBLAS_LAYOUT layout, CBLAS_SIDE Side, CBLAS_UPLO Uplo, const CBLAS_INT M, const CBLAS_INT N, const void *alpha, const void *A, const CBLAS_INT lda, const void *B, const CBLAS_INT ldb, const void *beta, void *C, const CBLAS_INT ldc)
Definition cblas_csymm.c:12
CBLAS_TRANSPOSE
Definition cblas.h:40
@ CblasNoTrans
Definition cblas.h:40
CBLAS_SIDE
Definition cblas.h:43
@ CblasLeft
Definition cblas.h:43
void cblas_ctrmm(CBLAS_LAYOUT layout, CBLAS_SIDE Side, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE TransA, CBLAS_DIAG Diag, const CBLAS_INT M, const CBLAS_INT N, const void *alpha, const void *A, const CBLAS_INT lda, void *B, const CBLAS_INT ldb)
Definition cblas_ctrmm.c:12
@ CblasColMajor
Definition cblas.h:39
@ CblasRowMajor
Definition cblas.h:39
CBLAS_DIAG
Definition cblas.h:42
void cblas_csyrk(CBLAS_LAYOUT layout, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE Trans, const CBLAS_INT N, const CBLAS_INT K, const void *alpha, const void *A, const CBLAS_INT lda, const void *beta, void *C, const CBLAS_INT ldc)
Definition cblas_csyrk.c:12
void cblas_chemm(CBLAS_LAYOUT layout, CBLAS_SIDE Side, CBLAS_UPLO Uplo, const CBLAS_INT M, const CBLAS_INT N, const void *alpha, const void *A, const CBLAS_INT lda, const void *B, const CBLAS_INT ldb, const void *beta, void *C, const CBLAS_INT ldc)
Definition cblas_chemm.c:12
void cblas_cgemm(CBLAS_LAYOUT layout, CBLAS_TRANSPOSE TransA, CBLAS_TRANSPOSE TransB, const CBLAS_INT M, const CBLAS_INT N, const CBLAS_INT K, const void *alpha, const void *A, const CBLAS_INT lda, const void *B, const CBLAS_INT ldb, const void *beta, void *C, const CBLAS_INT ldc)
Definition cblas_cgemm.c:12
void cblas_ctrsm(CBLAS_LAYOUT layout, CBLAS_SIDE Side, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE TransA, CBLAS_DIAG Diag, const CBLAS_INT M, const CBLAS_INT N, const void *alpha, const void *A, const CBLAS_INT lda, void *B, const CBLAS_INT ldb)
Definition cblas_ctrsm.c:12
#define CBLAS_INT
Definition cblas.h:24
void cblas_cgemmtr(CBLAS_LAYOUT layout, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE TransA, CBLAS_TRANSPOSE TransB, const CBLAS_INT N, const CBLAS_INT K, const void *alpha, const void *A, const CBLAS_INT lda, const void *B, const CBLAS_INT ldb, const void *beta, void *C, const CBLAS_INT ldc)
#define F77_csyrk(...)
Definition cblas_f77.h:419
#define F77_cgemm(...)
Definition cblas_f77.h:415
#define F77_cgemmtr(...)
Definition cblas_f77.h:416
#define F77_csyr2k(...)
Definition cblas_f77.h:421
#define BLAS_FORTRAN_STRLEN_END
Definition cblas_f77.h:18
#define FORTRAN_STRLEN
Definition cblas_f77.h:21
#define F77_cher2k(...)
Definition cblas_f77.h:422
#define F77_ctrsm(...)
Definition cblas_f77.h:424
#define F77_ctrmm(...)
Definition cblas_f77.h:423
#define F77_chemm(...)
Definition cblas_f77.h:418
#define F77_cherk(...)
Definition cblas_f77.h:420
#define F77_csymm(...)
Definition cblas_f77.h:417
void get_diag_type(char *type, CBLAS_DIAG *diag)
Definition auxiliary.c:25
void get_side_type(char *type, CBLAS_SIDE *side)
Definition auxiliary.c:32
void get_uplo_type(char *type, CBLAS_UPLO *uplo)
Definition auxiliary.c:18
void get_transpose_type(char *type, CBLAS_TRANSPOSE *trans)
Definition auxiliary.c:8