-
Notifications
You must be signed in to change notification settings - Fork 5
/
Copy pathlas2.c
1804 lines (1516 loc) · 57.7 KB
/
las2.c
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
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*
Copyright © 2002, University of Tennessee Research Foundation.
All rights reserved.
Redistribution and use in source and binary forms, with or without
modification, are permitted provided that the following conditions are met:
* Redistributions of source code must retain the above copyright notice, this
list of conditions and the following disclaimer.
Redistributions in binary form must reproduce the above copyright notice,
this list of conditions and the following disclaimer in the documentation
and/or other materials provided with the distribution.
* Neither the name of the University of Tennessee nor the names of its
contributors may be used to endorse or promote products derived from this
software without specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
POSSIBILITY OF SUCH DAMAGE.
*/
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <errno.h>
#include <math.h>
#include <fcntl.h>
#include "svdlib.h"
#include "svdutil.h"
#define MAXLL 2
#define LMTNW 100000000 /* max. size of working area allowed */
enum storeVals {STORQ = 1, RETRQ, STORP, RETRP};
static char *error_msg[] = { /* error messages used by function *
* check_parameters */
NULL,
"",
"ENDL MUST BE LESS THAN ENDR",
"REQUESTED DIMENSIONS CANNOT EXCEED NUM ITERATIONS",
"ONE OF YOUR DIMENSIONS IS LESS THAN OR EQUAL TO ZERO",
"NUM ITERATIONS (NUMBER OF LANCZOS STEPS) IS INVALID",
"REQUESTED DIMENSIONS (NUMBER OF EIGENPAIRS DESIRED) IS INVALID",
"6*N+4*ITERATIONS+1 + ITERATIONS*ITERATIONS CANNOT EXCEED NW",
"6*N+4*ITERATIONS+1 CANNOT EXCEED NW", NULL};
double **LanStore, *OPBTemp;
double eps, eps1, reps, eps34;
long ierr;
/*
double rnm, anorm, tol;
FILE *fp_out1, *fp_out2;
*/
void purge(long n, long ll, double *r, double *q, double *ra,
double *qa, double *wrk, double *eta, double *oldeta, long step,
double *rnmp, double tol);
void ortbnd(double *alf, double *eta, double *oldeta, double *bet, long step,
double rnm);
double startv(SMat A, double *wptr[], long step, long n);
void store(long, long, long, double *);
void imtql2(long, long, double *, double *, double *);
void imtqlb(long n, double d[], double e[], double bnd[]);
void write_header(long, long, double, double, long, double, long, long,
long);
long check_parameters(SMat A, long dimensions, long iterations,
double endl, double endr, long vectors);
int lanso(SMat A, long iterations, long dimensions, double endl,
double endr, double *ritz, double *bnd, double *wptr[],
long *neigp, long n);
long ritvec(long n, SMat A, SVDRec R, double kappa, double *ritz,
double *bnd, double *alf, double *bet, double *w2,
long steps, long neig);
long lanczos_step(SMat A, long first, long last, double *wptr[],
double *alf, double *eta, double *oldeta,
double *bet, long *ll, long *enough, double *rnmp,
double *tolp, long n);
void stpone(SMat A, double *wrkptr[], double *rnmp, double *tolp, long n);
long error_bound(long *, double, double, double *, double *, long step,
double tol);
void machar(long *ibeta, long *it, long *irnd, long *machep, long *negep);
/***********************************************************************
* *
* main() *
* Sparse SVD(A) via Eigensystem of A'A symmetric Matrix *
* (double precision) *
* *
***********************************************************************/
/***********************************************************************
Description
-----------
This sample program uses landr to compute singular triplets of A via
the equivalent symmetric eigenvalue problem
B x = lambda x, where x' = (u',v'), lambda = sigma**2,
where sigma is a singular value of A,
B = A'A , and A is m (nrow) by n (ncol) (nrow >> ncol),
so that {u,sqrt(lambda),v} is a singular triplet of A.
(A' = transpose of A)
User supplied routines: svd_opa, opb, store, timer
svd_opa( x,y) takes an n-vector x and returns A*x in y.
svd_opb(ncol,x,y) takes an n-vector x and returns B*x in y.
Based on operation flag isw, store(n,isw,j,s) stores/retrieves
to/from storage a vector of length n in s.
User should edit timer() with an appropriate call to an intrinsic
timing routine that returns elapsed user time.
External parameters
-------------------
Defined and documented in las2.h
Local parameters
----------------
(input)
endl left end of interval containing unwanted eigenvalues of B
endr right end of interval containing unwanted eigenvalues of B
kappa relative accuracy of ritz values acceptable as eigenvalues
of B
vectors is not equal to 1
r work array
n dimension of the eigenproblem for matrix B (ncol)
dimensions upper limit of desired number of singular triplets of A
iterations upper limit of desired number of Lanczos steps
nnzero number of nonzeros in A
vectors 1 indicates both singular values and singular vectors are
wanted and they can be found in output file lav2;
0 indicates only singular values are wanted
(output)
ritz array of ritz values
bnd array of error bounds
d array of singular values
memory total memory allocated in bytes to solve the B-eigenproblem
Functions used
--------------
BLAS svd_daxpy, svd_dscal, svd_ddot
USER svd_opa, svd_opb, timer
MISC write_header, check_parameters
LAS2 landr
Precision
---------
All floating-point calculations are done in double precision;
variables are declared as long and double.
LAS2 development
----------------
LAS2 is a C translation of the Fortran-77 LAS2 from the SVDPACK
library written by Michael W. Berry, University of Tennessee,
Dept. of Computer Science, 107 Ayres Hall, Knoxville, TN, 37996-1301
31 Jan 1992: Date written
Theresa H. Do
University of Tennessee
Dept. of Computer Science
107 Ayres Hall
Knoxville, TN, 37996-1301
internet: [email protected]
***********************************************************************/
/***********************************************************************
* *
* check_parameters() *
* *
***********************************************************************/
/***********************************************************************
Description
-----------
Function validates input parameters and returns error code (long)
Parameters
----------
(input)
dimensions upper limit of desired number of eigenpairs of B
iterations upper limit of desired number of lanczos steps
n dimension of the eigenproblem for matrix B
endl left end of interval containing unwanted eigenvalues of B
endr right end of interval containing unwanted eigenvalues of B
vectors 1 indicates both eigenvalues and eigenvectors are wanted
and they can be found in lav2; 0 indicates eigenvalues only
nnzero number of nonzero elements in input matrix (matrix A)
***********************************************************************/
long check_parameters(SMat A, long dimensions, long iterations,
double endl, double endr, long vectors) {
long error_index;
error_index = 0;
if (endl >/*=*/ endr) error_index = 2;
else if (dimensions > iterations) error_index = 3;
else if (A->cols <= 0 || A->rows <= 0) error_index = 4;
/*else if (n > A->cols || n > A->rows) error_index = 1;*/
else if (iterations <= 0 || iterations > A->cols || iterations > A->rows)
error_index = 5;
else if (dimensions <= 0 || dimensions > iterations) error_index = 6;
if (error_index)
svd_error("svdLAS2 parameter error: %s\n", error_msg[error_index]);
return(error_index);
}
/***********************************************************************
* *
* write_header() *
* Function writes out header of output file containing ritz values *
* *
***********************************************************************/
void write_header(long iterations, long dimensions, double endl, double endr,
long vectors, double kappa, long nrow, long ncol,
long vals) {
printf("SOLVING THE [A^TA] EIGENPROBLEM\n");
printf("NO. OF ROWS = %6ld\n", nrow);
printf("NO. OF COLUMNS = %6ld\n", ncol);
printf("NO. OF NON-ZERO VALUES = %6ld\n", vals);
printf("MATRIX DENSITY = %6.2f%%\n",
((float) vals / nrow) * 100 / ncol);
/* printf("ORDER OF MATRIX A = %5ld\n", n); */
printf("MAX. NO. OF LANCZOS STEPS = %6ld\n", iterations);
printf("MAX. NO. OF EIGENPAIRS = %6ld\n", dimensions);
printf("LEFT END OF THE INTERVAL = %9.2E\n", endl);
printf("RIGHT END OF THE INTERVAL = %9.2E\n", endr);
printf("KAPPA = %9.2E\n", kappa);
/* printf("WANT S-VECTORS? [T/F] = %c\n", (vectors) ? 'T' : 'F'); */
printf("\n");
return;
}
/***********************************************************************
* *
* landr() *
* Lanczos algorithm with selective orthogonalization *
* Using Simon's Recurrence *
* (double precision) *
* *
***********************************************************************/
/***********************************************************************
Description
-----------
landr() is the LAS2 driver routine that, upon entry,
(1) checks for the validity of input parameters of the
B-eigenproblem
(2) determines several machine constants
(3) makes a Lanczos run
(4) calculates B-eigenvectors (singular vectors of A) if requested
by user
arguments
---------
(input)
n dimension of the eigenproblem for A'A
iterations upper limit of desired number of Lanczos steps
dimensions upper limit of desired number of eigenpairs
nnzero number of nonzeros in matrix A
endl left end of interval containing unwanted eigenvalues of B
endr right end of interval containing unwanted eigenvalues of B
vectors 1 indicates both eigenvalues and eigenvectors are wanted
and they can be found in output file lav2;
0 indicates only eigenvalues are wanted
kappa relative accuracy of ritz values acceptable as eigenvalues
of B (singular values of A)
r work array
(output)
j number of Lanczos steps actually taken
neig number of ritz values stabilized
ritz array to hold the ritz values
bnd array to hold the error bounds
External parameters
-------------------
Defined and documented in las2.h
local parameters
-------------------
ibeta radix for the floating-point representation
it number of base ibeta digits in the floating-point significand
irnd floating-point addition rounded or chopped
machep machine relative precision or round-off error
negeps largest negative integer
wptr array of pointers each pointing to a work space
Functions used
--------------
MISC svd_dmax, machar, check_parameters
LAS2 ritvec, lanso
***********************************************************************/
SVDRec svdLAS2A(SMat A, long dimensions) {
double end[2] = {-1.0e-30, 1.0e-30};
double kappa = 1e-6;
if (!A) {
svd_error("svdLAS2A called with NULL array\n");
return NULL;
}
return svdLAS2(A, dimensions, 0, end, kappa);
}
SVDRec svdLAS2(SMat A, long dimensions, long iterations, double end[2],
double kappa) {
char transpose = FALSE;
long ibeta, it, irnd, machep, negep, n, i, steps, nsig, neig, m;
double *wptr[10], *ritz, *bnd;
SVDRec R = NULL;
ierr = 0; // reset the global error flag
svdResetCounters();
m = svd_imin(A->rows, A->cols);
if (dimensions <= 0 || dimensions > m)
dimensions = m;
if (iterations <= 0 || iterations > m)
iterations = m;
if (iterations < dimensions) iterations = dimensions;
/* Write output header */
if (SVDVerbosity > 0)
write_header(iterations, dimensions, end[0], end[1], TRUE, kappa, A->rows,
A->cols, A->vals);
/* Check parameters */
if (check_parameters(A, dimensions, iterations, end[0], end[1], TRUE))
return NULL;
/* If A is wide, the SVD is computed on its transpose for speed. */
if (A->cols >= A->rows * 1.2) {
if (SVDVerbosity > 0) printf("TRANSPOSING THE MATRIX FOR SPEED\n");
transpose = TRUE;
A = svdTransposeS(A);
}
n = A->cols;
/* Compute machine precision */
machar(&ibeta, &it, &irnd, &machep, &negep);
eps1 = eps * sqrt((double) n);
reps = sqrt(eps);
eps34 = reps * sqrt(reps);
/* Allocate temporary space. */
if (!(wptr[0] = svd_doubleArray(n, TRUE, "las2: wptr[0]"))) goto abort;
if (!(wptr[1] = svd_doubleArray(n, FALSE, "las2: wptr[1]"))) goto abort;
if (!(wptr[2] = svd_doubleArray(n, FALSE, "las2: wptr[2]"))) goto abort;
if (!(wptr[3] = svd_doubleArray(n, FALSE, "las2: wptr[3]"))) goto abort;
if (!(wptr[4] = svd_doubleArray(n, FALSE, "las2: wptr[4]"))) goto abort;
if (!(wptr[5] = svd_doubleArray(n, FALSE, "las2: wptr[5]"))) goto abort;
if (!(wptr[6] = svd_doubleArray(iterations, FALSE, "las2: wptr[6]")))
goto abort;
if (!(wptr[7] = svd_doubleArray(iterations, FALSE, "las2: wptr[7]")))
goto abort;
if (!(wptr[8] = svd_doubleArray(iterations, FALSE, "las2: wptr[8]")))
goto abort;
if (!(wptr[9] = svd_doubleArray(iterations + 1, FALSE, "las2: wptr[9]")))
goto abort;
/* Calloc may be unnecessary: */
if (!(ritz = svd_doubleArray(iterations + 1, TRUE, "las2: ritz")))
goto abort;
/* Calloc may be unnecessary: */
if (!(bnd = svd_doubleArray(iterations + 1, TRUE, "las2: bnd")))
goto abort;
memset(bnd, 127, (iterations + 1) * sizeof(double));
if (!(LanStore = (double **) calloc(iterations + MAXLL, sizeof(double *))))
goto abort;
if (!(OPBTemp = svd_doubleArray(A->rows, FALSE, "las2: OPBTemp")))
goto abort;
/* Actually run the lanczos thing: */
steps = lanso(A, iterations, dimensions, end[0], end[1], ritz, bnd, wptr,
&neig, n);
/* Print some stuff. */
if (SVDVerbosity > 0) {
printf("NUMBER OF LANCZOS STEPS = %6ld\n"
"RITZ VALUES STABILIZED = %6ld\n", steps + 1, neig);
}
if (SVDVerbosity > 2) {
printf("\nCOMPUTED RITZ VALUES (ERROR BNDS)\n");
for (i = 0; i <= steps; i++)
printf("%3ld %22.14E (%11.2E)\n", i + 1, ritz[i], bnd[i]);
}
SAFE_FREE(wptr[0]);
SAFE_FREE(wptr[1]);
SAFE_FREE(wptr[2]);
SAFE_FREE(wptr[3]);
SAFE_FREE(wptr[4]);
SAFE_FREE(wptr[7]);
SAFE_FREE(wptr[8]);
/* Compute eigenvectors */
kappa = svd_dmax(fabs(kappa), eps34);
R = svdNewSVDRec();
if (!R) {
svd_error("svdLAS2: allocation of R failed");
goto cleanup;
}
R->d = /*svd_imin(nsig, dimensions)*/dimensions;
R->Ut = svdNewDMat(R->d, A->rows);
R->S = svd_doubleArray(R->d, TRUE, "las2: R->s");
R->Vt = svdNewDMat(R->d, A->cols);
if (!R->Ut || !R->S || !R->Vt) {
svd_error("svdLAS2: allocation of R failed");
goto cleanup;
}
nsig = ritvec(n, A, R, kappa, ritz, bnd, wptr[6], wptr[9], wptr[5], steps,
neig);
if (SVDVerbosity > 1) {
printf("\nSINGULAR VALUES: ");
svdWriteDenseArray(R->S, R->d, "-", FALSE);
if (SVDVerbosity > 2) {
printf("\nLEFT SINGULAR VECTORS (transpose of U): ");
svdWriteDenseMatrix(R->Ut, "-", SVD_F_DT);
printf("\nRIGHT SINGULAR VECTORS (transpose of V): ");
svdWriteDenseMatrix(R->Vt, "-", SVD_F_DT);
}
}
if (SVDVerbosity > 0) {
printf("SINGULAR VALUES FOUND = %6d\n"
"SIGNIFICANT VALUES = %6ld\n", R->d, nsig);
}
cleanup:
for (i = 0; i <= 9; i++)
SAFE_FREE(wptr[i]);
SAFE_FREE(ritz);
SAFE_FREE(bnd);
if (LanStore) {
for (i = 0; i < iterations + MAXLL; i++)
SAFE_FREE(LanStore[i]);
SAFE_FREE(LanStore);
}
SAFE_FREE(OPBTemp);
/* This swaps and transposes the singular matrices if A was transposed. */
if (R && transpose) {
DMat T;
svdFreeSMat(A);
T = R->Ut;
R->Ut = R->Vt;
R->Vt = T;
}
return R;
abort:
svd_error("svdLAS2: fatal error, aborting");
return NULL;
}
/***********************************************************************
* *
* ritvec() *
* Function computes the singular vectors of matrix A *
* *
***********************************************************************/
/***********************************************************************
Description
-----------
This function is invoked by landr() only if eigenvectors of the A'A
eigenproblem are desired. When called, ritvec() computes the
singular vectors of A and writes the result to an unformatted file.
Parameters
----------
(input)
nrow number of rows of A
steps number of Lanczos iterations performed
fp_out2 pointer to unformatted output file
n dimension of matrix A
kappa relative accuracy of ritz values acceptable as
eigenvalues of A'A
ritz array of ritz values
bnd array of error bounds
alf array of diagonal elements of the tridiagonal matrix T
bet array of off-diagonal elements of T
w1, w2 work space
(output)
xv1 array of eigenvectors of A'A (right singular vectors of A)
ierr error code
0 for normal return from imtql2()
k if convergence did not occur for k-th eigenvalue in
imtql2()
nsig number of accepted ritz values based on kappa
(local)
s work array which is initialized to the identity matrix
of order (j + 1) upon calling imtql2(). After the call,
s contains the orthonormal eigenvectors of the symmetric
tridiagonal matrix T
Functions used
--------------
BLAS svd_dscal, svd_dcopy, svd_daxpy
USER store
imtql2
***********************************************************************/
void rotateArray(double *a, int size, int x) {
int i, j, n, start;
double t1, t2;
if (x == 0) return;
j = start = 0;
t1 = a[0];
for (i = 0; i < size; i++) {
n = (j >= x) ? j - x : j + size - x;
t2 = a[n];
a[n] = t1;
t1 = t2;
j = n;
if (j == start) {
start = ++j;
t1 = a[j];
}
}
}
long ritvec(long n, SMat A, SVDRec R, double kappa, double *ritz, double *bnd,
double *alf, double *bet, double *w2, long steps, long neig) {
long js, jsq, i, k, /*size,*/ id2, tmp, nsig, x;
double *s, *xv2, tmp0, tmp1, xnorm, *w1 = R->Vt->value[0];
js = steps + 1;
jsq = js * js;
/*size = sizeof(double) * n;*/
s = svd_doubleArray(jsq, TRUE, "ritvec: s");
xv2 = svd_doubleArray(n, FALSE, "ritvec: xv2");
/* initialize s to an identity matrix */
for (i = 0; i < jsq; i+= (js+1)) s[i] = 1.0;
svd_dcopy(js, alf, 1, w1, -1);
svd_dcopy(steps, &bet[1], 1, &w2[1], -1);
/* on return from imtql2(), w1 contains eigenvalues in ascending
* order and s contains the corresponding eigenvectors */
imtql2(js, js, w1, w2, s);
/*fwrite((char *)&n, sizeof(n), 1, fp_out2);
fwrite((char *)&js, sizeof(js), 1, fp_out2);
fwrite((char *)&kappa, sizeof(kappa), 1, fp_out2);*/
/*id = 0;*/
nsig = 0;
if (ierr) {
R->d = 0;
} else {
x = 0;
id2 = jsq - js;
for (k = 0; k < js; k++) {
tmp = id2;
if (bnd[k] <= kappa * fabs(ritz[k]) && k > js-neig-1) {
if (--x < 0) x = R->d - 1;
w1 = R->Vt->value[x];
for (i = 0; i < n; i++) w1[i] = 0.0;
for (i = 0; i < js; i++) {
store(n, RETRQ, i, w2);
svd_daxpy(n, s[tmp], w2, 1, w1, 1);
tmp -= js;
}
/*fwrite((char *)w1, size, 1, fp_out2);*/
/* store the w1 vector row-wise in array xv1;
* size of xv1 is (steps+1) * (nrow+ncol) elements
* and each vector, even though only ncol long,
* will have (nrow+ncol) elements in xv1.
* It is as if xv1 is a 2-d array (steps+1) by
* (nrow+ncol) and each vector occupies a row */
/* j is the index in the R arrays, which are sorted by high to low
singular values. */
/*for (i = 0; i < n; i++) R->Vt->value[x]xv1[id++] = w1[i];*/
/*id += nrow;*/
nsig++;
}
id2++;
}
SAFE_FREE(s);
/* Rotate the singular vectors and values. */
/* x is now the location of the highest singular value. */
rotateArray(R->Vt->value[0], R->Vt->rows * R->Vt->cols,
x * R->Vt->cols);
R->d = svd_imin(R->d, nsig);
for (x = 0; x < R->d; x++) {
/* multiply by matrix B first */
svd_opb(A, R->Vt->value[x], xv2, OPBTemp);
tmp0 = svd_ddot(n, R->Vt->value[x], 1, xv2, 1);
svd_daxpy(n, -tmp0, R->Vt->value[x], 1, xv2, 1);
tmp0 = sqrt(tmp0);
xnorm = sqrt(svd_ddot(n, xv2, 1, xv2, 1));
/* multiply by matrix A to get (scaled) left s-vector */
svd_opa(A, R->Vt->value[x], R->Ut->value[x]);
tmp1 = 1.0 / tmp0;
svd_dscal(A->rows, tmp1, R->Ut->value[x], 1);
xnorm *= tmp1;
bnd[i] = xnorm;
R->S[x] = tmp0;
}
}
SAFE_FREE(s);
SAFE_FREE(xv2);
return nsig;
}
/***********************************************************************
* *
* lanso() *
* *
***********************************************************************/
/***********************************************************************
Description
-----------
Function determines when the restart of the Lanczos algorithm should
occur and when it should terminate.
Arguments
---------
(input)
n dimension of the eigenproblem for matrix B
iterations upper limit of desired number of lanczos steps
dimensions upper limit of desired number of eigenpairs
endl left end of interval containing unwanted eigenvalues
endr right end of interval containing unwanted eigenvalues
ritz array to hold the ritz values
bnd array to hold the error bounds
wptr array of pointers that point to work space:
wptr[0]-wptr[5] six vectors of length n
wptr[6] array to hold diagonal of the tridiagonal matrix T
wptr[9] array to hold off-diagonal of T
wptr[7] orthogonality estimate of Lanczos vectors at
step j
wptr[8] orthogonality estimate of Lanczos vectors at
step j-1
(output)
j number of Lanczos steps actually taken
neig number of ritz values stabilized
ritz array to hold the ritz values
bnd array to hold the error bounds
ierr (globally declared) error flag
ierr = 8192 if stpone() fails to find a starting vector
ierr = k if convergence did not occur for k-th eigenvalue
in imtqlb()
ierr = 0 otherwise
Functions used
--------------
LAS stpone, error_bound, lanczos_step
MISC svd_dsort2
UTILITY svd_imin, svd_imax
***********************************************************************/
int lanso(SMat A, long iterations, long dimensions, double endl,
double endr, double *ritz, double *bnd, double *wptr[],
long *neigp, long n) {
double *alf, *eta, *oldeta, *bet, *wrk, rnm, tol;
long ll, first, last, ENOUGH, id2, id3, i, l, neig, j = 0, intro = 0;
alf = wptr[6];
eta = wptr[7];
oldeta = wptr[8];
bet = wptr[9];
wrk = wptr[5];
/* take the first step */
stpone(A, wptr, &rnm, &tol, n);
if (!rnm || ierr) return 0;
eta[0] = eps1;
oldeta[0] = eps1;
ll = 0;
first = 1;
last = svd_imin(dimensions + svd_imax(8, dimensions), iterations);
ENOUGH = FALSE;
/*id1 = 0;*/
while (/*id1 < dimensions && */!ENOUGH) {
if (rnm <= tol) rnm = 0.0;
/* the actual lanczos loop */
j = lanczos_step(A, first, last, wptr, alf, eta, oldeta, bet, &ll,
&ENOUGH, &rnm, &tol, n);
if (ENOUGH) j = j - 1;
else j = last - 1;
first = j + 1;
bet[j+1] = rnm;
/* analyze T */
l = 0;
for (id2 = 0; id2 < j; id2++) {
if (l > j) break;
for (i = l; i <= j; i++) if (!bet[i+1]) break;
if (i > j) i = j;
/* now i is at the end of an unreduced submatrix */
svd_dcopy(i-l+1, &alf[l], 1, &ritz[l], -1);
svd_dcopy(i-l, &bet[l+1], 1, &wrk[l+1], -1);
imtqlb(i-l+1, &ritz[l], &wrk[l], &bnd[l]);
if (ierr) {
svd_error("svdLAS2: imtqlb failed to converge (ierr = %ld)\n", ierr);
svd_error(" l = %ld i = %ld\n", l, i);
for (id3 = l; id3 <= i; id3++)
svd_error(" %ld %lg %lg %lg\n",
id3, ritz[id3], wrk[id3], bnd[id3]);
}
for (id3 = l; id3 <= i; id3++)
bnd[id3] = rnm * fabs(bnd[id3]);
l = i + 1;
}
/* sort eigenvalues into increasing order */
svd_dsort2((j+1) / 2, j + 1, ritz, bnd);
/* for (i = 0; i < iterations; i++)
printf("%f ", ritz[i]);
printf("\n"); */
/* massage error bounds for very close ritz values */
neig = error_bound(&ENOUGH, endl, endr, ritz, bnd, j, tol);
*neigp = neig;
/* should we stop? */
if (neig < dimensions) {
if (!neig) {
last = first + 9;
intro = first;
} else last = first + svd_imax(3, 1 + ((j - intro) * (dimensions-neig)) /
neig);
last = svd_imin(last, iterations);
} else ENOUGH = TRUE;
ENOUGH = ENOUGH || first >= iterations;
/* id1++; */
/* printf("id1=%d dimen=%d first=%d\n", id1, dimensions, first); */
}
store(n, STORQ, j, wptr[1]);
return j;
}
/***********************************************************************
* *
* lanczos_step() *
* *
***********************************************************************/
/***********************************************************************
Description
-----------
Function embodies a single Lanczos step
Arguments
---------
(input)
n dimension of the eigenproblem for matrix B
first start of index through loop
last end of index through loop
wptr array of pointers pointing to work space
alf array to hold diagonal of the tridiagonal matrix T
eta orthogonality estimate of Lanczos vectors at step j
oldeta orthogonality estimate of Lanczos vectors at step j-1
bet array to hold off-diagonal of T
ll number of intitial Lanczos vectors in local orthog.
(has value of 0, 1 or 2)
enough stop flag
Functions used
--------------
BLAS svd_ddot, svd_dscal, svd_daxpy, svd_datx, svd_dcopy
USER store
LAS purge, ortbnd, startv
UTILITY svd_imin, svd_imax
***********************************************************************/
long lanczos_step(SMat A, long first, long last, double *wptr[],
double *alf, double *eta, double *oldeta,
double *bet, long *ll, long *enough, double *rnmp,
double *tolp, long n) {
double t, *mid, rnm = *rnmp, tol = *tolp, anorm;
long i, j;
for (j=first; j<last; j++) {
mid = wptr[2];
wptr[2] = wptr[1];
wptr[1] = mid;
mid = wptr[3];
wptr[3] = wptr[4];
wptr[4] = mid;
store(n, STORQ, j-1, wptr[2]);
if (j-1 < MAXLL) store(n, STORP, j-1, wptr[4]);
bet[j] = rnm;
/* restart if invariant subspace is found */
if (!bet[j]) {
rnm = startv(A, wptr, j, n);
if (ierr) return j;
if (!rnm) *enough = TRUE;
}
if (*enough) {
/* added by Doug... */
/* These lines fix a bug that occurs with low-rank matrices */
mid = wptr[2];
wptr[2] = wptr[1];
wptr[1] = mid;
/* ...added by Doug */
break;
}
/* take a lanczos step */
t = 1.0 / rnm;
svd_datx(n, t, wptr[0], 1, wptr[1], 1);
svd_dscal(n, t, wptr[3], 1);
svd_opb(A, wptr[3], wptr[0], OPBTemp);
svd_daxpy(n, -rnm, wptr[2], 1, wptr[0], 1);
alf[j] = svd_ddot(n, wptr[0], 1, wptr[3], 1);
svd_daxpy(n, -alf[j], wptr[1], 1, wptr[0], 1);
/* orthogonalize against initial lanczos vectors */
if (j <= MAXLL && (fabs(alf[j-1]) > 4.0 * fabs(alf[j])))
*ll = j;
for (i=0; i < svd_imin(*ll, j-1); i++) {
store(n, RETRP, i, wptr[5]);
t = svd_ddot(n, wptr[5], 1, wptr[0], 1);
store(n, RETRQ, i, wptr[5]);
svd_daxpy(n, -t, wptr[5], 1, wptr[0], 1);
eta[i] = eps1;
oldeta[i] = eps1;
}
/* extended local reorthogonalization */
t = svd_ddot(n, wptr[0], 1, wptr[4], 1);
svd_daxpy(n, -t, wptr[2], 1, wptr[0], 1);
if (bet[j] > 0.0) bet[j] = bet[j] + t;
t = svd_ddot(n, wptr[0], 1, wptr[3], 1);
svd_daxpy(n, -t, wptr[1], 1, wptr[0], 1);
alf[j] = alf[j] + t;
svd_dcopy(n, wptr[0], 1, wptr[4], 1);
rnm = sqrt(svd_ddot(n, wptr[0], 1, wptr[4], 1));
anorm = bet[j] + fabs(alf[j]) + rnm;
tol = reps * anorm;
/* update the orthogonality bounds */
ortbnd(alf, eta, oldeta, bet, j, rnm);
/* restore the orthogonality state when needed */
purge(n, *ll, wptr[0], wptr[1], wptr[4], wptr[3], wptr[5], eta, oldeta,
j, &rnm, tol);
if (rnm <= tol) rnm = 0.0;
}
*rnmp = rnm;
*tolp = tol;
return j;
}
/***********************************************************************
* *
* ortbnd() *
* *
***********************************************************************/
/***********************************************************************
Description
-----------
Funtion updates the eta recurrence
Arguments
---------
(input)
alf array to hold diagonal of the tridiagonal matrix T
eta orthogonality estimate of Lanczos vectors at step j
oldeta orthogonality estimate of Lanczos vectors at step j-1
bet array to hold off-diagonal of T
n dimension of the eigenproblem for matrix B
j dimension of T
rnm norm of the next residual vector
eps1 roundoff estimate for dot product of two unit vectors
(output)
eta orthogonality estimate of Lanczos vectors at step j+1
oldeta orthogonality estimate of Lanczos vectors at step j
Functions used
--------------
BLAS svd_dswap
***********************************************************************/
void ortbnd(double *alf, double *eta, double *oldeta, double *bet, long step,
double rnm) {
long i;
if (step < 1) return;
if (rnm) {
if (step > 1) {
oldeta[0] = (bet[1] * eta[1] + (alf[0]-alf[step]) * eta[0] -
bet[step] * oldeta[0]) / rnm + eps1;
}
for (i=1; i<=step-2; i++)
oldeta[i] = (bet[i+1] * eta[i+1] + (alf[i]-alf[step]) * eta[i] +
bet[i] * eta[i-1] - bet[step] * oldeta[i])/rnm + eps1;
}
oldeta[step-1] = eps1;
svd_dswap(step, oldeta, 1, eta, 1);
eta[step] = eps1;
return;
}
/***********************************************************************
* *
* purge() *
* *
***********************************************************************/
/***********************************************************************
Description
-----------
Function examines the state of orthogonality between the new Lanczos
vector and the previous ones to decide whether re-orthogonalization
should be performed
Arguments
---------