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