1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
|
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <math.h>
#include <time.h>
#include "lapack-aux.h"
#define MACRO(B) do {B} while (0)
#define ERROR(CODE) MACRO(return CODE;)
#define REQUIRES(COND, CODE) MACRO(if(!(COND)) {ERROR(CODE);})
#define MIN(A,B) ((A)<(B)?(A):(B))
#define MAX(A,B) ((A)>(B)?(A):(B))
#ifdef DBGL
#define DEBUGMSG(M) printf("LAPACK Wrapper "M"\n: "); size_t t0 = time(NULL);
#define OK MACRO(printf("%ld s\n",time(0)-t0); return 0;);
#else
#define DEBUGMSG(M)
#define OK return 0;
#endif
#define CHECK(RES,CODE) MACRO(if(RES) return CODE;)
#define BAD_SIZE 2000
#define BAD_CODE 2001
#define MEM 2002
#define BAD_FILE 2003
#define SINGULAR 2004
#define NOCONVER 2005
#define NODEFPOS 2006
#define NOSPRTD 2007
//////////////////// real svd ////////////////////////////////////
int svd_l_R(KDMAT(a),DMAT(u), DVEC(s),DMAT(v)) {
integer m = ar;
integer n = ac;
integer q = MIN(m,n);
REQUIRES(ur==m && uc==m && sn==q && vr==n && vc==n,BAD_SIZE);
DEBUGMSG("svd_l_R");
double *B = (double*)malloc(m*n*sizeof(double));
CHECK(!B,MEM);
memcpy(B,ap,m*n*sizeof(double));
integer lwork = -1;
integer res;
// ask for optimal lwork
double ans;
//printf("ask zgesvd\n");
char* job = "A";
dgesvd_ (job,job,
&m,&n,B,&m,
sp,
up,&m,
vp,&n,
&ans, &lwork,
&res);
lwork = ceil(ans);
//printf("ans = %d\n",lwork);
double * work = (double*)malloc(lwork*sizeof(double));
CHECK(!work,MEM);
//printf("dgesdd\n");
dgesvd_ (job,job,
&m,&n,B,&m,
sp,
up,&m,
vp,&n,
work, &lwork,
&res);
CHECK(res,res);
free(work);
free(B);
OK
}
// (alternative version)
int svd_l_Rdd(KDMAT(a),DMAT(u), DVEC(s),DMAT(v)) {
integer m = ar;
integer n = ac;
integer q = MIN(m,n);
REQUIRES(ur==m && uc==m && sn==q && vr==n && vc==n,BAD_SIZE);
DEBUGMSG("svd_l_Rdd");
double *B = (double*)malloc(m*n*sizeof(double));
CHECK(!B,MEM);
memcpy(B,ap,m*n*sizeof(double));
integer* iwk = (integer*) malloc(8*q*sizeof(int));
CHECK(!iwk,MEM);
integer lwk = -1;
integer res;
// ask for optimal lwk
double ans;
//printf("ask dgesdd\n");
dgesdd_ ("A",&m,&n,B,&m,sp,up,&m,vp,&n,&ans,&lwk,iwk,&res);
lwk = 2*ceil(ans); // ????? otherwise 50x100 rejects lwk
//printf("lwk = %d\n",lwk);
double * workv = (double*)malloc(lwk*sizeof(double));
CHECK(!workv,MEM);
//printf("dgesdd\n");
dgesdd_ ("A",&m,&n,B,&m,sp,up,&m,vp,&n,workv,&lwk,iwk,&res);
CHECK(res,res);
free(iwk);
free(workv);
free(B);
OK
}
//////////////////// complex svd ////////////////////////////////////
// not in clapack.h
int zgesvd_(char *jobu, char *jobvt, integer *m, integer *n,
doublecomplex *a, integer *lda, doublereal *s, doublecomplex *u,
integer *ldu, doublecomplex *vt, integer *ldvt, doublecomplex *work,
integer *lwork, doublereal *rwork, integer *info);
int svd_l_C(KCMAT(a),CMAT(u), DVEC(s),CMAT(v)) {
integer m = ar;
integer n = ac;
integer q = MIN(m,n);
REQUIRES(ur==m && uc==m && sn==q && vr==n && vc==n,BAD_SIZE);
DEBUGMSG("svd_l_C");
double *B = (double*)malloc(2*m*n*sizeof(double));
CHECK(!B,MEM);
memcpy(B,ap,m*n*2*sizeof(double));
double *rwork = (double*) malloc(5*q*sizeof(double));
CHECK(!rwork,MEM);
integer lwork = -1;
integer res;
// ask for optimal lwork
doublecomplex ans;
//printf("ask zgesvd\n");
char* job = "A";
zgesvd_ (job,job,
&m,&n,(doublecomplex*)B,&m,
sp,
(doublecomplex*)up,&m,
(doublecomplex*)vp,&n,
&ans, &lwork,
rwork,
&res);
lwork = ceil(ans.r);
//printf("ans = %d\n",lwork);
doublecomplex * work = (doublecomplex*)malloc(lwork*2*sizeof(double));
CHECK(!work,MEM);
//printf("zgesvd\n");
zgesvd_ (job,job,
&m,&n,(doublecomplex*)B,&m,
sp,
(doublecomplex*)up,&m,
(doublecomplex*)vp,&n,
work, &lwork,
rwork,
&res);
CHECK(res,res);
free(work);
free(rwork);
free(B);
OK
}
//////////////////// general complex eigensystem ////////////
int eig_l_C(KCMAT(a),CMAT(u), CVEC(s),CMAT(v)) {
integer n = ar;
REQUIRES(n>=2 && ac==n && (ur==1 || (ur==n && uc==n)) && sn==n && (vr==1 || (vr==n && vc==n)),BAD_SIZE);
DEBUGMSG("eig_l_C");
double *B = (double*)malloc(2*n*n*sizeof(double));
CHECK(!B,MEM);
memcpy(B,ap,n*n*2*sizeof(double));
double *rwork = (double*) malloc(2*n*sizeof(double));
CHECK(!rwork,MEM);
integer lwork = -1;
char jobvl = ur==1?'N':'V';
char jobvr = vr==1?'N':'V';
integer res;
// ask for optimal lwork
doublecomplex ans;
//printf("ask zgeev\n");
zgeev_ (&jobvl,&jobvr,
&n,(doublecomplex*)B,&n,
(doublecomplex*)sp,
(doublecomplex*)up,&n,
(doublecomplex*)vp,&n,
&ans, &lwork,
rwork,
&res);
lwork = ceil(ans.r);
//printf("ans = %d\n",lwork);
doublecomplex * work = (doublecomplex*)malloc(lwork*2*sizeof(double));
CHECK(!work,MEM);
//printf("zgeev\n");
zgeev_ (&jobvl,&jobvr,
&n,(doublecomplex*)B,&n,
(doublecomplex*)sp,
(doublecomplex*)up,&n,
(doublecomplex*)vp,&n,
work, &lwork,
rwork,
&res);
CHECK(res,res);
free(work);
free(rwork);
free(B);
OK
}
//////////////////// general real eigensystem ////////////
int eig_l_R(KDMAT(a),DMAT(u), CVEC(s),DMAT(v)) {
integer n = ar;
REQUIRES(n>=2 && ac == n && (ur==1 || (ur==n && uc==n)) && sn==n && (vr==1 || (vr==n && vc==n)),BAD_SIZE);
DEBUGMSG("eig_l_R");
double *B = (double*)malloc(n*n*sizeof(double));
CHECK(!B,MEM);
memcpy(B,ap,n*n*sizeof(double));
integer lwork = -1;
char jobvl = ur==1?'N':'V';
char jobvr = vr==1?'N':'V';
integer res;
// ask for optimal lwork
double ans;
//printf("ask dgeev\n");
dgeev_ (&jobvl,&jobvr,
&n,B,&n,
sp, sp+n,
up,&n,
vp,&n,
&ans, &lwork,
&res);
lwork = ceil(ans);
//printf("ans = %d\n",lwork);
double * work = (double*)malloc(lwork*sizeof(double));
CHECK(!work,MEM);
//printf("dgeev\n");
dgeev_ (&jobvl,&jobvr,
&n,B,&n,
sp, sp+n,
up,&n,
vp,&n,
work, &lwork,
&res);
CHECK(res,res);
free(work);
free(B);
OK
}
//////////////////// symmetric real eigensystem ////////////
int eig_l_S(KDMAT(a),DVEC(s),DMAT(v)) {
integer n = ar;
REQUIRES(n>=2 && ac == n && sn==n && (vr==1 || (vr==n && vc==n)),BAD_SIZE);
DEBUGMSG("eig_l_S");
memcpy(vp,ap,n*n*sizeof(double));
integer lwork = -1;
char jobz = vr==1?'N':'V';
char uplo = 'U';
integer res;
// ask for optimal lwork
double ans;
//printf("ask dsyev\n");
dsyev_ (&jobz,&uplo,
&n,vp,&n,
sp,
&ans, &lwork,
&res);
lwork = ceil(ans);
//printf("ans = %d\n",lwork);
double * work = (double*)malloc(lwork*sizeof(double));
CHECK(!work,MEM);
dsyev_ (&jobz,&uplo,
&n,vp,&n,
sp,
work, &lwork,
&res);
CHECK(res,res);
free(work);
OK
}
//////////////////// hermitian complex eigensystem ////////////
int eig_l_H(KCMAT(a),DVEC(s),CMAT(v)) {
integer n = ar;
REQUIRES(n>=2 && ac==n && sn==n && (vr==1 || (vr==n && vc==n)),BAD_SIZE);
DEBUGMSG("eig_l_H");
memcpy(vp,ap,2*n*n*sizeof(double));
double *rwork = (double*) malloc((3*n-2)*sizeof(double));
CHECK(!rwork,MEM);
integer lwork = -1;
char jobz = vr==1?'N':'V';
char uplo = 'U';
integer res;
// ask for optimal lwork
doublecomplex ans;
//printf("ask zheev\n");
zheev_ (&jobz,&uplo,
&n,(doublecomplex*)vp,&n,
sp,
&ans, &lwork,
rwork,
&res);
lwork = ceil(ans.r);
//printf("ans = %d\n",lwork);
doublecomplex * work = (doublecomplex*)malloc(lwork*2*sizeof(double));
CHECK(!work,MEM);
zheev_ (&jobz,&uplo,
&n,(doublecomplex*)vp,&n,
sp,
work, &lwork,
rwork,
&res);
CHECK(res,res);
free(work);
free(rwork);
OK
}
//////////////////// general real linear system ////////////
int linearSolveR_l(KDMAT(a),KDMAT(b),DMAT(x)) {
integer n = ar;
integer nhrs = bc;
REQUIRES(n>=1 && ar==ac && ar==br,BAD_SIZE);
DEBUGMSG("linearSolveR_l");
double*AC = (double*)malloc(n*n*sizeof(double));
memcpy(AC,ap,n*n*sizeof(double));
memcpy(xp,bp,n*nhrs*sizeof(double));
integer * ipiv = (integer*)malloc(n*sizeof(integer));
integer res;
dgesv_ (&n,&nhrs,
AC, &n,
ipiv,
xp, &n,
&res);
if(res>0) {
return SINGULAR;
}
CHECK(res,res);
free(ipiv);
free(AC);
OK
}
//////////////////// general complex linear system ////////////
int linearSolveC_l(KCMAT(a),KCMAT(b),CMAT(x)) {
integer n = ar;
integer nhrs = bc;
REQUIRES(n>=1 && ar==ac && ar==br,BAD_SIZE);
DEBUGMSG("linearSolveC_l");
double*AC = (double*)malloc(2*n*n*sizeof(double));
memcpy(AC,ap,2*n*n*sizeof(double));
memcpy(xp,bp,2*n*nhrs*sizeof(double));
integer * ipiv = (integer*)malloc(n*sizeof(integer));
integer res;
zgesv_ (&n,&nhrs,
(doublecomplex*)AC, &n,
ipiv,
(doublecomplex*)xp, &n,
&res);
if(res>0) {
return SINGULAR;
}
CHECK(res,res);
free(ipiv);
free(AC);
OK
}
//////////////////// least squares real linear system ////////////
int linearSolveLSR_l(KDMAT(a),KDMAT(b),DMAT(x)) {
integer m = ar;
integer n = ac;
integer nrhs = bc;
integer ldb = xr;
REQUIRES(m>=1 && n>=1 && ar==br && xr==MAX(m,n) && xc == bc, BAD_SIZE);
DEBUGMSG("linearSolveLSR_l");
double*AC = (double*)malloc(m*n*sizeof(double));
memcpy(AC,ap,m*n*sizeof(double));
if (m>=n) {
memcpy(xp,bp,m*nrhs*sizeof(double));
} else {
int k;
for(k = 0; k<nrhs; k++) {
memcpy(xp+ldb*k,bp+m*k,m*sizeof(double));
}
}
integer res;
integer lwork = -1;
double ans;
dgels_ ("N",&m,&n,&nrhs,
AC,&m,
xp,&ldb,
&ans,&lwork,
&res);
lwork = ceil(ans);
//printf("ans = %d\n",lwork);
double * work = (double*)malloc(lwork*sizeof(double));
dgels_ ("N",&m,&n,&nrhs,
AC,&m,
xp,&ldb,
work,&lwork,
&res);
if(res>0) {
return SINGULAR;
}
CHECK(res,res);
free(work);
free(AC);
OK
}
//////////////////// least squares complex linear system ////////////
int linearSolveLSC_l(KCMAT(a),KCMAT(b),CMAT(x)) {
#ifdef _WIN32
return NOSPRTD;
#else
integer m = ar;
integer n = ac;
integer nrhs = bc;
integer ldb = xr;
REQUIRES(m>=1 && n>=1 && ar==br && xr==MAX(m,n) && xc == bc, BAD_SIZE);
DEBUGMSG("linearSolveLSC_l");
double*AC = (double*)malloc(2*m*n*sizeof(double));
memcpy(AC,ap,2*m*n*sizeof(double));
memcpy(AC,ap,2*m*n*sizeof(double));
if (m>=n) {
memcpy(xp,bp,2*m*nrhs*sizeof(double));
} else {
int k;
for(k = 0; k<nrhs; k++) {
memcpy(xp+2*ldb*k,bp+2*m*k,m*2*sizeof(double));
}
}
integer res;
integer lwork = -1;
doublecomplex ans;
zgels_ ("N",&m,&n,&nrhs,
(doublecomplex*)AC,&m,
(doublecomplex*)xp,&ldb,
&ans,&lwork,
&res);
lwork = ceil(ans.r);
//printf("ans = %d\n",lwork);
doublecomplex * work = (doublecomplex*)malloc(lwork*sizeof(doublecomplex));
zgels_ ("N",&m,&n,&nrhs,
(doublecomplex*)AC,&m,
(doublecomplex*)xp,&ldb,
work,&lwork,
&res);
if(res>0) {
return SINGULAR;
}
CHECK(res,res);
free(work);
free(AC);
OK
#endif
}
//////////////////// least squares real linear system using SVD ////////////
int linearSolveSVDR_l(double rcond,KDMAT(a),KDMAT(b),DMAT(x)) {
integer m = ar;
integer n = ac;
integer nrhs = bc;
integer ldb = xr;
REQUIRES(m>=1 && n>=1 && ar==br && xr==MAX(m,n) && xc == bc, BAD_SIZE);
DEBUGMSG("linearSolveSVDR_l");
double*AC = (double*)malloc(m*n*sizeof(double));
double*S = (double*)malloc(MIN(m,n)*sizeof(double));
memcpy(AC,ap,m*n*sizeof(double));
if (m>=n) {
memcpy(xp,bp,m*nrhs*sizeof(double));
} else {
int k;
for(k = 0; k<nrhs; k++) {
memcpy(xp+ldb*k,bp+m*k,m*sizeof(double));
}
}
integer res;
integer lwork = -1;
integer rank;
double ans;
dgelss_ (&m,&n,&nrhs,
AC,&m,
xp,&ldb,
S,
&rcond,&rank,
&ans,&lwork,
&res);
lwork = ceil(ans);
//printf("ans = %d\n",lwork);
double * work = (double*)malloc(lwork*sizeof(double));
dgelss_ (&m,&n,&nrhs,
AC,&m,
xp,&ldb,
S,
&rcond,&rank,
work,&lwork,
&res);
if(res>0) {
return NOCONVER;
}
CHECK(res,res);
free(work);
free(S);
free(AC);
OK
}
//////////////////// least squares complex linear system using SVD ////////////
// not in clapack.h
int zgelss_(integer *m, integer *n, integer *nhrs,
doublecomplex *a, integer *lda, doublecomplex *b, integer *ldb, doublereal *s,
doublereal *rcond, integer* rank,
doublecomplex *work, integer* lwork, doublereal* rwork,
integer *info);
int linearSolveSVDC_l(double rcond, KCMAT(a),KCMAT(b),CMAT(x)) {
#ifdef _WIN32
return NOSPRTD;
#else
integer m = ar;
integer n = ac;
integer nrhs = bc;
integer ldb = xr;
REQUIRES(m>=1 && n>=1 && ar==br && xr==MAX(m,n) && xc == bc, BAD_SIZE);
DEBUGMSG("linearSolveSVDC_l");
double*AC = (double*)malloc(2*m*n*sizeof(double));
double*S = (double*)malloc(MIN(m,n)*sizeof(double));
double*RWORK = (double*)malloc(5*MIN(m,n)*sizeof(double));
memcpy(AC,ap,2*m*n*sizeof(double));
if (m>=n) {
memcpy(xp,bp,2*m*nrhs*sizeof(double));
} else {
int k;
for(k = 0; k<nrhs; k++) {
memcpy(xp+2*ldb*k,bp+2*m*k,m*2*sizeof(double));
}
}
integer res;
integer lwork = -1;
integer rank;
doublecomplex ans;
zgelss_ (&m,&n,&nrhs,
(doublecomplex*)AC,&m,
(doublecomplex*)xp,&ldb,
S,
&rcond,&rank,
&ans,&lwork,
RWORK,
&res);
lwork = ceil(ans.r);
//printf("ans = %d\n",lwork);
doublecomplex * work = (doublecomplex*)malloc(lwork*sizeof(doublecomplex));
zgelss_ (&m,&n,&nrhs,
(doublecomplex*)AC,&m,
(doublecomplex*)xp,&ldb,
S,
&rcond,&rank,
work,&lwork,
RWORK,
&res);
if(res>0) {
return NOCONVER;
}
CHECK(res,res);
free(work);
free(RWORK);
free(S);
free(AC);
OK
#endif
}
//////////////////// Cholesky factorization /////////////////////////
int chol_l_H(KCMAT(a),CMAT(l)) {
integer n = ar;
REQUIRES(n>=1 && ac == n && lr==n && lc==n,BAD_SIZE);
DEBUGMSG("chol_l_H");
memcpy(lp,ap,n*n*sizeof(doublecomplex));
char uplo = 'U';
integer res;
zpotrf_ (&uplo,&n,(doublecomplex*)lp,&n,&res);
CHECK(res>0,NODEFPOS);
CHECK(res,res);
doublecomplex zero = {0.,0.};
int r,c;
for (r=0; r<lr-1; r++) {
for(c=r+1; c<lc; c++) {
((doublecomplex*)lp)[r*lc+c] = zero;
}
}
OK
}
int chol_l_S(KDMAT(a),DMAT(l)) {
integer n = ar;
REQUIRES(n>=1 && ac == n && lr==n && lc==n,BAD_SIZE);
DEBUGMSG("chol_l_S");
memcpy(lp,ap,n*n*sizeof(double));
char uplo = 'U';
integer res;
dpotrf_ (&uplo,&n,lp,&n,&res);
CHECK(res>0,NODEFPOS);
CHECK(res,res);
int r,c;
for (r=0; r<lr-1; r++) {
for(c=r+1; c<lc; c++) {
lp[r*lc+c] = 0.;
}
}
OK
}
//////////////////// QR factorization /////////////////////////
int qr_l_R(KDMAT(a), DVEC(tau), DMAT(r)) {
integer m = ar;
integer n = ac;
integer mn = MIN(m,n);
REQUIRES(m>=1 && n >=1 && rr== m && rc == n && taun == mn, BAD_SIZE);
DEBUGMSG("qr_l_R");
double *WORK = (double*)malloc(n*sizeof(double));
CHECK(!WORK,MEM);
memcpy(rp,ap,m*n*sizeof(double));
integer res;
dgeqr2_ (&m,&n,rp,&m,taup,WORK,&res);
CHECK(res,res);
free(WORK);
OK
}
int qr_l_C(KCMAT(a), CVEC(tau), CMAT(r)) {
integer m = ar;
integer n = ac;
integer mn = MIN(m,n);
REQUIRES(m>=1 && n >=1 && rr== m && rc == n && taun == mn, BAD_SIZE);
DEBUGMSG("qr_l_C");
doublecomplex *WORK = (doublecomplex*)malloc(n*sizeof(doublecomplex));
CHECK(!WORK,MEM);
memcpy(rp,ap,m*n*sizeof(doublecomplex));
integer res;
zgeqr2_ (&m,&n,(doublecomplex*)rp,&m,(doublecomplex*)taup,WORK,&res);
CHECK(res,res);
free(WORK);
OK
}
//////////////////// Hessenberg factorization /////////////////////////
int hess_l_R(KDMAT(a), DVEC(tau), DMAT(r)) {
integer m = ar;
integer n = ac;
integer mn = MIN(m,n);
REQUIRES(m>=1 && n == m && rr== m && rc == n && taun == mn-1, BAD_SIZE);
DEBUGMSG("hess_l_R");
integer lwork = 5*n; // fixme
double *WORK = (double*)malloc(lwork*sizeof(double));
CHECK(!WORK,MEM);
memcpy(rp,ap,m*n*sizeof(double));
integer res;
integer one = 1;
dgehrd_ (&n,&one,&n,rp,&n,taup,WORK,&lwork,&res);
CHECK(res,res);
free(WORK);
OK
}
int hess_l_C(KCMAT(a), CVEC(tau), CMAT(r)) {
integer m = ar;
integer n = ac;
integer mn = MIN(m,n);
REQUIRES(m>=1 && n == m && rr== m && rc == n && taun == mn-1, BAD_SIZE);
DEBUGMSG("hess_l_C");
integer lwork = 5*n; // fixme
doublecomplex *WORK = (doublecomplex*)malloc(lwork*sizeof(doublecomplex));
CHECK(!WORK,MEM);
memcpy(rp,ap,m*n*sizeof(doublecomplex));
integer res;
integer one = 1;
zgehrd_ (&n,&one,&n,(doublecomplex*)rp,&n,(doublecomplex*)taup,WORK,&lwork,&res);
CHECK(res,res);
free(WORK);
OK
}
//////////////////// Schur factorization /////////////////////////
int schur_l_R(KDMAT(a), DMAT(u), DMAT(s)) {
integer m = ar;
integer n = ac;
REQUIRES(m>=1 && n==m && ur==n && uc==n && sr==n && sc==n, BAD_SIZE);
DEBUGMSG("schur_l_R");
int k;
//printf("---------------------------\n");
//printf("%p: ",ap); for(k=0;k<n*n;k++) printf("%f ",ap[k]); printf("\n");
//printf("%p: ",up); for(k=0;k<n*n;k++) printf("%f ",up[k]); printf("\n");
//printf("%p: ",sp); for(k=0;k<n*n;k++) printf("%f ",sp[k]); printf("\n");
memcpy(sp,ap,n*n*sizeof(double));
integer lwork = 6*n; // fixme
double *WORK = (double*)malloc(lwork*sizeof(double));
double *WR = (double*)malloc(n*sizeof(double));
double *WI = (double*)malloc(n*sizeof(double));
// WR and WI not really required in this call
logical *BWORK = (logical*)malloc(n*sizeof(logical));
integer res;
integer sdim;
dgees_ ("V","N",NULL,&n,sp,&n,&sdim,WR,WI,up,&n,WORK,&lwork,BWORK,&res);
//printf("%p: ",ap); for(k=0;k<n*n;k++) printf("%f ",ap[k]); printf("\n");
//printf("%p: ",up); for(k=0;k<n*n;k++) printf("%f ",up[k]); printf("\n");
//printf("%p: ",sp); for(k=0;k<n*n;k++) printf("%f ",sp[k]); printf("\n");
if(res>0) {
return NOCONVER;
}
CHECK(res,res);
free(WR);
free(WI);
free(BWORK);
free(WORK);
OK
}
int schur_l_C(KCMAT(a), CMAT(u), CMAT(s)) {
#ifdef _WIN32
return NOSPRTD;
#else
integer m = ar;
integer n = ac;
REQUIRES(m>=1 && n==m && ur==n && uc==n && sr==n && sc==n, BAD_SIZE);
DEBUGMSG("schur_l_C");
memcpy(sp,ap,n*n*sizeof(doublecomplex));
integer lwork = 6*n; // fixme
doublecomplex *WORK = (doublecomplex*)malloc(lwork*sizeof(doublecomplex));
doublecomplex *W = (doublecomplex*)malloc(n*sizeof(doublecomplex));
// W not really required in this call
logical *BWORK = (logical*)malloc(n*sizeof(logical));
double *RWORK = (double*)malloc(n*sizeof(double));
integer res;
integer sdim;
zgees_ ("V","N",NULL,&n,(doublecomplex*)sp,&n,&sdim,W,
(doublecomplex*)up,&n,
WORK,&lwork,RWORK,BWORK,&res);
if(res>0) {
return NOCONVER;
}
CHECK(res,res);
free(W);
free(BWORK);
free(WORK);
OK
#endif
}
//////////////////// LU factorization /////////////////////////
int lu_l_R(KDMAT(a), DVEC(ipiv), DMAT(r)) {
integer m = ar;
integer n = ac;
integer mn = MIN(m,n);
REQUIRES(m>=1 && n >=1 && ipivn == mn, BAD_SIZE);
DEBUGMSG("lu_l_R");
integer* auxipiv = (integer*)malloc(mn*sizeof(integer));
memcpy(rp,ap,m*n*sizeof(double));
integer res;
dgetrf_ (&m,&n,rp,&m,auxipiv,&res);
if(res>0) {
res = 0; // fixme
}
CHECK(res,res);
int k;
for (k=0; k<mn; k++) {
ipivp[k] = auxipiv[k];
}
free(auxipiv);
OK
}
int lu_l_C(KCMAT(a), DVEC(ipiv), CMAT(r)) {
integer m = ar;
integer n = ac;
integer mn = MIN(m,n);
REQUIRES(m>=1 && n >=1 && ipivn == mn, BAD_SIZE);
DEBUGMSG("lu_l_C");
integer* auxipiv = (integer*)malloc(mn*sizeof(integer));
memcpy(rp,ap,m*n*sizeof(doublecomplex));
integer res;
zgetrf_ (&m,&n,(doublecomplex*)rp,&m,auxipiv,&res);
if(res>0) {
res = 0; // fixme
}
CHECK(res,res);
int k;
for (k=0; k<mn; k++) {
ipivp[k] = auxipiv[k];
}
free(auxipiv);
OK
}
|