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
|
SUBROUTINE DGSVJ0( JOBV, M, N, A, LDA, D, SVA, MV, V, LDV, EPS,
+ SFMIN, TOL, NSWEEP, WORK, LWORK, INFO )
*
* -- LAPACK routine (version 3.2) --
*
* -- Contributed by Zlatko Drmac of the University of Zagreb and --
* -- Kresimir Veselic of the Fernuniversitaet Hagen --
* -- November 2008 --
*
* -- LAPACK is a software package provided by Univ. of Tennessee, --
* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
*
* This routine is also part of SIGMA (version 1.23, October 23. 2008.)
* SIGMA is a library of algorithms for highly accurate algorithms for
* computation of SVD, PSVD, QSVD, (H,K)-SVD, and for solution of the
* eigenvalue problems Hx = lambda M x, H M x = lambda x with H, M > 0.
*
IMPLICIT NONE
* .. Scalar Arguments ..
INTEGER INFO, LDA, LDV, LWORK, M, MV, N, NSWEEP
DOUBLE PRECISION EPS, SFMIN, TOL
CHARACTER*1 JOBV
* ..
* .. Array Arguments ..
DOUBLE PRECISION A( LDA, * ), SVA( N ), D( N ), V( LDV, * ),
+ WORK( LWORK )
* ..
*
* Purpose
* =======
*
* DGSVJ0 is called from DGESVJ as a pre-processor and that is its main
* purpose. It applies Jacobi rotations in the same way as DGESVJ does, but
* it does not check convergence (stopping criterion). Few tuning
* parameters (marked by [TP]) are available for the implementer.
*
* Further details
* ~~~~~~~~~~~~~~~
* DGSVJ0 is used just to enable SGESVJ to call a simplified version of
* itself to work on a submatrix of the original matrix.
*
* Contributors
* ~~~~~~~~~~~~
* Zlatko Drmac (Zagreb, Croatia) and Kresimir Veselic (Hagen, Germany)
*
* Bugs, Examples and Comments
* ~~~~~~~~~~~~~~~~~~~~~~~~~~~
* Please report all bugs and send interesting test examples and comments to
* drmac@math.hr. Thank you.
*
* Arguments
* =========
*
* JOBV (input) CHARACTER*1
* Specifies whether the output from this procedure is used
* to compute the matrix V:
* = 'V': the product of the Jacobi rotations is accumulated
* by postmulyiplying the N-by-N array V.
* (See the description of V.)
* = 'A': the product of the Jacobi rotations is accumulated
* by postmulyiplying the MV-by-N array V.
* (See the descriptions of MV and V.)
* = 'N': the Jacobi rotations are not accumulated.
*
* M (input) INTEGER
* The number of rows of the input matrix A. M >= 0.
*
* N (input) INTEGER
* The number of columns of the input matrix A.
* M >= N >= 0.
*
* A (input/output) REAL array, dimension (LDA,N)
* On entry, M-by-N matrix A, such that A*diag(D) represents
* the input matrix.
* On exit,
* A_onexit * D_onexit represents the input matrix A*diag(D)
* post-multiplied by a sequence of Jacobi rotations, where the
* rotation threshold and the total number of sweeps are given in
* TOL and NSWEEP, respectively.
* (See the descriptions of D, TOL and NSWEEP.)
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,M).
*
* D (input/workspace/output) REAL array, dimension (N)
* The array D accumulates the scaling factors from the fast scaled
* Jacobi rotations.
* On entry, A*diag(D) represents the input matrix.
* On exit, A_onexit*diag(D_onexit) represents the input matrix
* post-multiplied by a sequence of Jacobi rotations, where the
* rotation threshold and the total number of sweeps are given in
* TOL and NSWEEP, respectively.
* (See the descriptions of A, TOL and NSWEEP.)
*
* SVA (input/workspace/output) REAL array, dimension (N)
* On entry, SVA contains the Euclidean norms of the columns of
* the matrix A*diag(D).
* On exit, SVA contains the Euclidean norms of the columns of
* the matrix onexit*diag(D_onexit).
*
* MV (input) INTEGER
* If JOBV .EQ. 'A', then MV rows of V are post-multipled by a
* sequence of Jacobi rotations.
* If JOBV = 'N', then MV is not referenced.
*
* V (input/output) REAL array, dimension (LDV,N)
* If JOBV .EQ. 'V' then N rows of V are post-multipled by a
* sequence of Jacobi rotations.
* If JOBV .EQ. 'A' then MV rows of V are post-multipled by a
* sequence of Jacobi rotations.
* If JOBV = 'N', then V is not referenced.
*
* LDV (input) INTEGER
* The leading dimension of the array V, LDV >= 1.
* If JOBV = 'V', LDV .GE. N.
* If JOBV = 'A', LDV .GE. MV.
*
* EPS (input) INTEGER
* EPS = SLAMCH('Epsilon')
*
* SFMIN (input) INTEGER
* SFMIN = SLAMCH('Safe Minimum')
*
* TOL (input) REAL
* TOL is the threshold for Jacobi rotations. For a pair
* A(:,p), A(:,q) of pivot columns, the Jacobi rotation is
* applied only if DABS(COS(angle(A(:,p),A(:,q)))) .GT. TOL.
*
* NSWEEP (input) INTEGER
* NSWEEP is the number of sweeps of Jacobi rotations to be
* performed.
*
* WORK (workspace) REAL array, dimension LWORK.
*
* LWORK (input) INTEGER
* LWORK is the dimension of WORK. LWORK .GE. M.
*
* INFO (output) INTEGER
* = 0 : successful exit.
* < 0 : if INFO = -i, then the i-th argument had an illegal value
*
* =====================================================================
*
* .. Local Parameters ..
DOUBLE PRECISION ZERO, HALF, ONE, TWO
PARAMETER ( ZERO = 0.0D0, HALF = 0.5D0, ONE = 1.0D0,
+ TWO = 2.0D0 )
* ..
* .. Local Scalars ..
DOUBLE PRECISION AAPP, AAPP0, AAPQ, AAQQ, APOAQ, AQOAP, BIG,
+ BIGTHETA, CS, MXAAPQ, MXSINJ, ROOTBIG, ROOTEPS,
+ ROOTSFMIN, ROOTTOL, SMALL, SN, T, TEMP1, THETA,
+ THSIGN
INTEGER BLSKIP, EMPTSW, i, ibr, IERR, igl, IJBLSK, ir1,
+ ISWROT, jbc, jgl, KBL, LKAHEAD, MVL, NBL,
+ NOTROT, p, PSKIPPED, q, ROWSKIP, SWBAND
LOGICAL APPLV, ROTOK, RSVEC
* ..
* .. Local Arrays ..
DOUBLE PRECISION FASTR( 5 )
* ..
* .. Intrinsic Functions ..
INTRINSIC DABS, DMAX1, DBLE, MIN0, DSIGN, DSQRT
* ..
* .. External Functions ..
DOUBLE PRECISION DDOT, DNRM2
INTEGER IDAMAX
LOGICAL LSAME
EXTERNAL IDAMAX, LSAME, DDOT, DNRM2
* ..
* .. External Subroutines ..
EXTERNAL DAXPY, DCOPY, DLASCL, DLASSQ, DROTM, DSWAP
* ..
* .. Executable Statements ..
*
APPLV = LSAME( JOBV, 'A' )
RSVEC = LSAME( JOBV, 'V' )
IF( .NOT.( RSVEC .OR. APPLV .OR. LSAME( JOBV, 'N' ) ) ) THEN
INFO = -1
ELSE IF( M.LT.0 ) THEN
INFO = -2
ELSE IF( ( N.LT.0 ) .OR. ( N.GT.M ) ) THEN
INFO = -3
ELSE IF( LDA.LT.M ) THEN
INFO = -5
ELSE IF( MV.LT.0 ) THEN
INFO = -8
ELSE IF( LDV.LT.M ) THEN
INFO = -10
ELSE IF( TOL.LE.EPS ) THEN
INFO = -13
ELSE IF( NSWEEP.LT.0 ) THEN
INFO = -14
ELSE IF( LWORK.LT.M ) THEN
INFO = -16
ELSE
INFO = 0
END IF
*
* #:(
IF( INFO.NE.0 ) THEN
CALL XERBLA( 'DGSVJ0', -INFO )
RETURN
END IF
*
IF( RSVEC ) THEN
MVL = N
ELSE IF( APPLV ) THEN
MVL = MV
END IF
RSVEC = RSVEC .OR. APPLV
ROOTEPS = DSQRT( EPS )
ROOTSFMIN = DSQRT( SFMIN )
SMALL = SFMIN / EPS
BIG = ONE / SFMIN
ROOTBIG = ONE / ROOTSFMIN
BIGTHETA = ONE / ROOTEPS
ROOTTOL = DSQRT( TOL )
*
*
* -#- Row-cyclic Jacobi SVD algorithm with column pivoting -#-
*
EMPTSW = ( N*( N-1 ) ) / 2
NOTROT = 0
FASTR( 1 ) = ZERO
*
* -#- Row-cyclic pivot strategy with de Rijk's pivoting -#-
*
SWBAND = 0
*[TP] SWBAND is a tuning parameter. It is meaningful and effective
* if SGESVJ is used as a computational routine in the preconditioned
* Jacobi SVD algorithm SGESVJ. For sweeps i=1:SWBAND the procedure
* ......
KBL = MIN0( 8, N )
*[TP] KBL is a tuning parameter that defines the tile size in the
* tiling of the p-q loops of pivot pairs. In general, an optimal
* value of KBL depends on the matrix dimensions and on the
* parameters of the computer's memory.
*
NBL = N / KBL
IF( ( NBL*KBL ).NE.N )NBL = NBL + 1
BLSKIP = ( KBL**2 ) + 1
*[TP] BLKSKIP is a tuning parameter that depends on SWBAND and KBL.
ROWSKIP = MIN0( 5, KBL )
*[TP] ROWSKIP is a tuning parameter.
LKAHEAD = 1
*[TP] LKAHEAD is a tuning parameter.
SWBAND = 0
PSKIPPED = 0
*
DO 1993 i = 1, NSWEEP
* .. go go go ...
*
MXAAPQ = ZERO
MXSINJ = ZERO
ISWROT = 0
*
NOTROT = 0
PSKIPPED = 0
*
DO 2000 ibr = 1, NBL
igl = ( ibr-1 )*KBL + 1
*
DO 1002 ir1 = 0, MIN0( LKAHEAD, NBL-ibr )
*
igl = igl + ir1*KBL
*
DO 2001 p = igl, MIN0( igl+KBL-1, N-1 )
* .. de Rijk's pivoting
q = IDAMAX( N-p+1, SVA( p ), 1 ) + p - 1
IF( p.NE.q ) THEN
CALL DSWAP( M, A( 1, p ), 1, A( 1, q ), 1 )
IF( RSVEC )CALL DSWAP( MVL, V( 1, p ), 1,
+ V( 1, q ), 1 )
TEMP1 = SVA( p )
SVA( p ) = SVA( q )
SVA( q ) = TEMP1
TEMP1 = D( p )
D( p ) = D( q )
D( q ) = TEMP1
END IF
*
IF( ir1.EQ.0 ) THEN
*
* Column norms are periodically updated by explicit
* norm computation.
* Caveat:
* Some BLAS implementations compute DNRM2(M,A(1,p),1)
* as DSQRT(DDOT(M,A(1,p),1,A(1,p),1)), which may result in
* overflow for ||A(:,p)||_2 > DSQRT(overflow_threshold), and
* undeflow for ||A(:,p)||_2 < DSQRT(underflow_threshold).
* Hence, DNRM2 cannot be trusted, not even in the case when
* the true norm is far from the under(over)flow boundaries.
* If properly implemented DNRM2 is available, the IF-THEN-ELSE
* below should read "AAPP = DNRM2( M, A(1,p), 1 ) * D(p)".
*
IF( ( SVA( p ).LT.ROOTBIG ) .AND.
+ ( SVA( p ).GT.ROOTSFMIN ) ) THEN
SVA( p ) = DNRM2( M, A( 1, p ), 1 )*D( p )
ELSE
TEMP1 = ZERO
AAPP = ZERO
CALL DLASSQ( M, A( 1, p ), 1, TEMP1, AAPP )
SVA( p ) = TEMP1*DSQRT( AAPP )*D( p )
END IF
AAPP = SVA( p )
ELSE
AAPP = SVA( p )
END IF
*
IF( AAPP.GT.ZERO ) THEN
*
PSKIPPED = 0
*
DO 2002 q = p + 1, MIN0( igl+KBL-1, N )
*
AAQQ = SVA( q )
IF( AAQQ.GT.ZERO ) THEN
*
AAPP0 = AAPP
IF( AAQQ.GE.ONE ) THEN
ROTOK = ( SMALL*AAPP ).LE.AAQQ
IF( AAPP.LT.( BIG / AAQQ ) ) THEN
AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
+ q ), 1 )*D( p )*D( q ) / AAQQ )
+ / AAPP
ELSE
CALL DCOPY( M, A( 1, p ), 1, WORK, 1 )
CALL DLASCL( 'G', 0, 0, AAPP, D( p ),
+ M, 1, WORK, LDA, IERR )
AAPQ = DDOT( M, WORK, 1, A( 1, q ),
+ 1 )*D( q ) / AAQQ
END IF
ELSE
ROTOK = AAPP.LE.( AAQQ / SMALL )
IF( AAPP.GT.( SMALL / AAQQ ) ) THEN
AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
+ q ), 1 )*D( p )*D( q ) / AAQQ )
+ / AAPP
ELSE
CALL DCOPY( M, A( 1, q ), 1, WORK, 1 )
CALL DLASCL( 'G', 0, 0, AAQQ, D( q ),
+ M, 1, WORK, LDA, IERR )
AAPQ = DDOT( M, WORK, 1, A( 1, p ),
+ 1 )*D( p ) / AAPP
END IF
END IF
*
MXAAPQ = DMAX1( MXAAPQ, DABS( AAPQ ) )
*
* TO rotate or NOT to rotate, THAT is the question ...
*
IF( DABS( AAPQ ).GT.TOL ) THEN
*
* .. rotate
* ROTATED = ROTATED + ONE
*
IF( ir1.EQ.0 ) THEN
NOTROT = 0
PSKIPPED = 0
ISWROT = ISWROT + 1
END IF
*
IF( ROTOK ) THEN
*
AQOAP = AAQQ / AAPP
APOAQ = AAPP / AAQQ
THETA = -HALF*DABS( AQOAP-APOAQ ) /
+ AAPQ
*
IF( DABS( THETA ).GT.BIGTHETA ) THEN
*
T = HALF / THETA
FASTR( 3 ) = T*D( p ) / D( q )
FASTR( 4 ) = -T*D( q ) / D( p )
CALL DROTM( M, A( 1, p ), 1,
+ A( 1, q ), 1, FASTR )
IF( RSVEC )CALL DROTM( MVL,
+ V( 1, p ), 1,
+ V( 1, q ), 1,
+ FASTR )
SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
+ ONE+T*APOAQ*AAPQ ) )
AAPP = AAPP*DSQRT( ONE-T*AQOAP*
+ AAPQ )
MXSINJ = DMAX1( MXSINJ, DABS( T ) )
*
ELSE
*
* .. choose correct signum for THETA and rotate
*
THSIGN = -DSIGN( ONE, AAPQ )
T = ONE / ( THETA+THSIGN*
+ DSQRT( ONE+THETA*THETA ) )
CS = DSQRT( ONE / ( ONE+T*T ) )
SN = T*CS
*
MXSINJ = DMAX1( MXSINJ, DABS( SN ) )
SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
+ ONE+T*APOAQ*AAPQ ) )
AAPP = AAPP*DSQRT( DMAX1( ZERO,
+ ONE-T*AQOAP*AAPQ ) )
*
APOAQ = D( p ) / D( q )
AQOAP = D( q ) / D( p )
IF( D( p ).GE.ONE ) THEN
IF( D( q ).GE.ONE ) THEN
FASTR( 3 ) = T*APOAQ
FASTR( 4 ) = -T*AQOAP
D( p ) = D( p )*CS
D( q ) = D( q )*CS
CALL DROTM( M, A( 1, p ), 1,
+ A( 1, q ), 1,
+ FASTR )
IF( RSVEC )CALL DROTM( MVL,
+ V( 1, p ), 1, V( 1, q ),
+ 1, FASTR )
ELSE
CALL DAXPY( M, -T*AQOAP,
+ A( 1, q ), 1,
+ A( 1, p ), 1 )
CALL DAXPY( M, CS*SN*APOAQ,
+ A( 1, p ), 1,
+ A( 1, q ), 1 )
D( p ) = D( p )*CS
D( q ) = D( q ) / CS
IF( RSVEC ) THEN
CALL DAXPY( MVL, -T*AQOAP,
+ V( 1, q ), 1,
+ V( 1, p ), 1 )
CALL DAXPY( MVL,
+ CS*SN*APOAQ,
+ V( 1, p ), 1,
+ V( 1, q ), 1 )
END IF
END IF
ELSE
IF( D( q ).GE.ONE ) THEN
CALL DAXPY( M, T*APOAQ,
+ A( 1, p ), 1,
+ A( 1, q ), 1 )
CALL DAXPY( M, -CS*SN*AQOAP,
+ A( 1, q ), 1,
+ A( 1, p ), 1 )
D( p ) = D( p ) / CS
D( q ) = D( q )*CS
IF( RSVEC ) THEN
CALL DAXPY( MVL, T*APOAQ,
+ V( 1, p ), 1,
+ V( 1, q ), 1 )
CALL DAXPY( MVL,
+ -CS*SN*AQOAP,
+ V( 1, q ), 1,
+ V( 1, p ), 1 )
END IF
ELSE
IF( D( p ).GE.D( q ) ) THEN
CALL DAXPY( M, -T*AQOAP,
+ A( 1, q ), 1,
+ A( 1, p ), 1 )
CALL DAXPY( M, CS*SN*APOAQ,
+ A( 1, p ), 1,
+ A( 1, q ), 1 )
D( p ) = D( p )*CS
D( q ) = D( q ) / CS
IF( RSVEC ) THEN
CALL DAXPY( MVL,
+ -T*AQOAP,
+ V( 1, q ), 1,
+ V( 1, p ), 1 )
CALL DAXPY( MVL,
+ CS*SN*APOAQ,
+ V( 1, p ), 1,
+ V( 1, q ), 1 )
END IF
ELSE
CALL DAXPY( M, T*APOAQ,
+ A( 1, p ), 1,
+ A( 1, q ), 1 )
CALL DAXPY( M,
+ -CS*SN*AQOAP,
+ A( 1, q ), 1,
+ A( 1, p ), 1 )
D( p ) = D( p ) / CS
D( q ) = D( q )*CS
IF( RSVEC ) THEN
CALL DAXPY( MVL,
+ T*APOAQ, V( 1, p ),
+ 1, V( 1, q ), 1 )
CALL DAXPY( MVL,
+ -CS*SN*AQOAP,
+ V( 1, q ), 1,
+ V( 1, p ), 1 )
END IF
END IF
END IF
END IF
END IF
*
ELSE
* .. have to use modified Gram-Schmidt like transformation
CALL DCOPY( M, A( 1, p ), 1, WORK, 1 )
CALL DLASCL( 'G', 0, 0, AAPP, ONE, M,
+ 1, WORK, LDA, IERR )
CALL DLASCL( 'G', 0, 0, AAQQ, ONE, M,
+ 1, A( 1, q ), LDA, IERR )
TEMP1 = -AAPQ*D( p ) / D( q )
CALL DAXPY( M, TEMP1, WORK, 1,
+ A( 1, q ), 1 )
CALL DLASCL( 'G', 0, 0, ONE, AAQQ, M,
+ 1, A( 1, q ), LDA, IERR )
SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
+ ONE-AAPQ*AAPQ ) )
MXSINJ = DMAX1( MXSINJ, SFMIN )
END IF
* END IF ROTOK THEN ... ELSE
*
* In the case of cancellation in updating SVA(q), SVA(p)
* recompute SVA(q), SVA(p).
IF( ( SVA( q ) / AAQQ )**2.LE.ROOTEPS )
+ THEN
IF( ( AAQQ.LT.ROOTBIG ) .AND.
+ ( AAQQ.GT.ROOTSFMIN ) ) THEN
SVA( q ) = DNRM2( M, A( 1, q ), 1 )*
+ D( q )
ELSE
T = ZERO
AAQQ = ZERO
CALL DLASSQ( M, A( 1, q ), 1, T,
+ AAQQ )
SVA( q ) = T*DSQRT( AAQQ )*D( q )
END IF
END IF
IF( ( AAPP / AAPP0 ).LE.ROOTEPS ) THEN
IF( ( AAPP.LT.ROOTBIG ) .AND.
+ ( AAPP.GT.ROOTSFMIN ) ) THEN
AAPP = DNRM2( M, A( 1, p ), 1 )*
+ D( p )
ELSE
T = ZERO
AAPP = ZERO
CALL DLASSQ( M, A( 1, p ), 1, T,
+ AAPP )
AAPP = T*DSQRT( AAPP )*D( p )
END IF
SVA( p ) = AAPP
END IF
*
ELSE
* A(:,p) and A(:,q) already numerically orthogonal
IF( ir1.EQ.0 )NOTROT = NOTROT + 1
PSKIPPED = PSKIPPED + 1
END IF
ELSE
* A(:,q) is zero column
IF( ir1.EQ.0 )NOTROT = NOTROT + 1
PSKIPPED = PSKIPPED + 1
END IF
*
IF( ( i.LE.SWBAND ) .AND.
+ ( PSKIPPED.GT.ROWSKIP ) ) THEN
IF( ir1.EQ.0 )AAPP = -AAPP
NOTROT = 0
GO TO 2103
END IF
*
2002 CONTINUE
* END q-LOOP
*
2103 CONTINUE
* bailed out of q-loop
SVA( p ) = AAPP
ELSE
SVA( p ) = AAPP
IF( ( ir1.EQ.0 ) .AND. ( AAPP.EQ.ZERO ) )
+ NOTROT = NOTROT + MIN0( igl+KBL-1, N ) - p
END IF
*
2001 CONTINUE
* end of the p-loop
* end of doing the block ( ibr, ibr )
1002 CONTINUE
* end of ir1-loop
*
*........................................................
* ... go to the off diagonal blocks
*
igl = ( ibr-1 )*KBL + 1
*
DO 2010 jbc = ibr + 1, NBL
*
jgl = ( jbc-1 )*KBL + 1
*
* doing the block at ( ibr, jbc )
*
IJBLSK = 0
DO 2100 p = igl, MIN0( igl+KBL-1, N )
*
AAPP = SVA( p )
*
IF( AAPP.GT.ZERO ) THEN
*
PSKIPPED = 0
*
DO 2200 q = jgl, MIN0( jgl+KBL-1, N )
*
AAQQ = SVA( q )
*
IF( AAQQ.GT.ZERO ) THEN
AAPP0 = AAPP
*
* -#- M x 2 Jacobi SVD -#-
*
* -#- Safe Gram matrix computation -#-
*
IF( AAQQ.GE.ONE ) THEN
IF( AAPP.GE.AAQQ ) THEN
ROTOK = ( SMALL*AAPP ).LE.AAQQ
ELSE
ROTOK = ( SMALL*AAQQ ).LE.AAPP
END IF
IF( AAPP.LT.( BIG / AAQQ ) ) THEN
AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
+ q ), 1 )*D( p )*D( q ) / AAQQ )
+ / AAPP
ELSE
CALL DCOPY( M, A( 1, p ), 1, WORK, 1 )
CALL DLASCL( 'G', 0, 0, AAPP, D( p ),
+ M, 1, WORK, LDA, IERR )
AAPQ = DDOT( M, WORK, 1, A( 1, q ),
+ 1 )*D( q ) / AAQQ
END IF
ELSE
IF( AAPP.GE.AAQQ ) THEN
ROTOK = AAPP.LE.( AAQQ / SMALL )
ELSE
ROTOK = AAQQ.LE.( AAPP / SMALL )
END IF
IF( AAPP.GT.( SMALL / AAQQ ) ) THEN
AAPQ = ( DDOT( M, A( 1, p ), 1, A( 1,
+ q ), 1 )*D( p )*D( q ) / AAQQ )
+ / AAPP
ELSE
CALL DCOPY( M, A( 1, q ), 1, WORK, 1 )
CALL DLASCL( 'G', 0, 0, AAQQ, D( q ),
+ M, 1, WORK, LDA, IERR )
AAPQ = DDOT( M, WORK, 1, A( 1, p ),
+ 1 )*D( p ) / AAPP
END IF
END IF
*
MXAAPQ = DMAX1( MXAAPQ, DABS( AAPQ ) )
*
* TO rotate or NOT to rotate, THAT is the question ...
*
IF( DABS( AAPQ ).GT.TOL ) THEN
NOTROT = 0
* ROTATED = ROTATED + 1
PSKIPPED = 0
ISWROT = ISWROT + 1
*
IF( ROTOK ) THEN
*
AQOAP = AAQQ / AAPP
APOAQ = AAPP / AAQQ
THETA = -HALF*DABS( AQOAP-APOAQ ) /
+ AAPQ
IF( AAQQ.GT.AAPP0 )THETA = -THETA
*
IF( DABS( THETA ).GT.BIGTHETA ) THEN
T = HALF / THETA
FASTR( 3 ) = T*D( p ) / D( q )
FASTR( 4 ) = -T*D( q ) / D( p )
CALL DROTM( M, A( 1, p ), 1,
+ A( 1, q ), 1, FASTR )
IF( RSVEC )CALL DROTM( MVL,
+ V( 1, p ), 1,
+ V( 1, q ), 1,
+ FASTR )
SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
+ ONE+T*APOAQ*AAPQ ) )
AAPP = AAPP*DSQRT( DMAX1( ZERO,
+ ONE-T*AQOAP*AAPQ ) )
MXSINJ = DMAX1( MXSINJ, DABS( T ) )
ELSE
*
* .. choose correct signum for THETA and rotate
*
THSIGN = -DSIGN( ONE, AAPQ )
IF( AAQQ.GT.AAPP0 )THSIGN = -THSIGN
T = ONE / ( THETA+THSIGN*
+ DSQRT( ONE+THETA*THETA ) )
CS = DSQRT( ONE / ( ONE+T*T ) )
SN = T*CS
MXSINJ = DMAX1( MXSINJ, DABS( SN ) )
SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
+ ONE+T*APOAQ*AAPQ ) )
AAPP = AAPP*DSQRT( ONE-T*AQOAP*
+ AAPQ )
*
APOAQ = D( p ) / D( q )
AQOAP = D( q ) / D( p )
IF( D( p ).GE.ONE ) THEN
*
IF( D( q ).GE.ONE ) THEN
FASTR( 3 ) = T*APOAQ
FASTR( 4 ) = -T*AQOAP
D( p ) = D( p )*CS
D( q ) = D( q )*CS
CALL DROTM( M, A( 1, p ), 1,
+ A( 1, q ), 1,
+ FASTR )
IF( RSVEC )CALL DROTM( MVL,
+ V( 1, p ), 1, V( 1, q ),
+ 1, FASTR )
ELSE
CALL DAXPY( M, -T*AQOAP,
+ A( 1, q ), 1,
+ A( 1, p ), 1 )
CALL DAXPY( M, CS*SN*APOAQ,
+ A( 1, p ), 1,
+ A( 1, q ), 1 )
IF( RSVEC ) THEN
CALL DAXPY( MVL, -T*AQOAP,
+ V( 1, q ), 1,
+ V( 1, p ), 1 )
CALL DAXPY( MVL,
+ CS*SN*APOAQ,
+ V( 1, p ), 1,
+ V( 1, q ), 1 )
END IF
D( p ) = D( p )*CS
D( q ) = D( q ) / CS
END IF
ELSE
IF( D( q ).GE.ONE ) THEN
CALL DAXPY( M, T*APOAQ,
+ A( 1, p ), 1,
+ A( 1, q ), 1 )
CALL DAXPY( M, -CS*SN*AQOAP,
+ A( 1, q ), 1,
+ A( 1, p ), 1 )
IF( RSVEC ) THEN
CALL DAXPY( MVL, T*APOAQ,
+ V( 1, p ), 1,
+ V( 1, q ), 1 )
CALL DAXPY( MVL,
+ -CS*SN*AQOAP,
+ V( 1, q ), 1,
+ V( 1, p ), 1 )
END IF
D( p ) = D( p ) / CS
D( q ) = D( q )*CS
ELSE
IF( D( p ).GE.D( q ) ) THEN
CALL DAXPY( M, -T*AQOAP,
+ A( 1, q ), 1,
+ A( 1, p ), 1 )
CALL DAXPY( M, CS*SN*APOAQ,
+ A( 1, p ), 1,
+ A( 1, q ), 1 )
D( p ) = D( p )*CS
D( q ) = D( q ) / CS
IF( RSVEC ) THEN
CALL DAXPY( MVL,
+ -T*AQOAP,
+ V( 1, q ), 1,
+ V( 1, p ), 1 )
CALL DAXPY( MVL,
+ CS*SN*APOAQ,
+ V( 1, p ), 1,
+ V( 1, q ), 1 )
END IF
ELSE
CALL DAXPY( M, T*APOAQ,
+ A( 1, p ), 1,
+ A( 1, q ), 1 )
CALL DAXPY( M,
+ -CS*SN*AQOAP,
+ A( 1, q ), 1,
+ A( 1, p ), 1 )
D( p ) = D( p ) / CS
D( q ) = D( q )*CS
IF( RSVEC ) THEN
CALL DAXPY( MVL,
+ T*APOAQ, V( 1, p ),
+ 1, V( 1, q ), 1 )
CALL DAXPY( MVL,
+ -CS*SN*AQOAP,
+ V( 1, q ), 1,
+ V( 1, p ), 1 )
END IF
END IF
END IF
END IF
END IF
*
ELSE
IF( AAPP.GT.AAQQ ) THEN
CALL DCOPY( M, A( 1, p ), 1, WORK,
+ 1 )
CALL DLASCL( 'G', 0, 0, AAPP, ONE,
+ M, 1, WORK, LDA, IERR )
CALL DLASCL( 'G', 0, 0, AAQQ, ONE,
+ M, 1, A( 1, q ), LDA,
+ IERR )
TEMP1 = -AAPQ*D( p ) / D( q )
CALL DAXPY( M, TEMP1, WORK, 1,
+ A( 1, q ), 1 )
CALL DLASCL( 'G', 0, 0, ONE, AAQQ,
+ M, 1, A( 1, q ), LDA,
+ IERR )
SVA( q ) = AAQQ*DSQRT( DMAX1( ZERO,
+ ONE-AAPQ*AAPQ ) )
MXSINJ = DMAX1( MXSINJ, SFMIN )
ELSE
CALL DCOPY( M, A( 1, q ), 1, WORK,
+ 1 )
CALL DLASCL( 'G', 0, 0, AAQQ, ONE,
+ M, 1, WORK, LDA, IERR )
CALL DLASCL( 'G', 0, 0, AAPP, ONE,
+ M, 1, A( 1, p ), LDA,
+ IERR )
TEMP1 = -AAPQ*D( q ) / D( p )
CALL DAXPY( M, TEMP1, WORK, 1,
+ A( 1, p ), 1 )
CALL DLASCL( 'G', 0, 0, ONE, AAPP,
+ M, 1, A( 1, p ), LDA,
+ IERR )
SVA( p ) = AAPP*DSQRT( DMAX1( ZERO,
+ ONE-AAPQ*AAPQ ) )
MXSINJ = DMAX1( MXSINJ, SFMIN )
END IF
END IF
* END IF ROTOK THEN ... ELSE
*
* In the case of cancellation in updating SVA(q)
* .. recompute SVA(q)
IF( ( SVA( q ) / AAQQ )**2.LE.ROOTEPS )
+ THEN
IF( ( AAQQ.LT.ROOTBIG ) .AND.
+ ( AAQQ.GT.ROOTSFMIN ) ) THEN
SVA( q ) = DNRM2( M, A( 1, q ), 1 )*
+ D( q )
ELSE
T = ZERO
AAQQ = ZERO
CALL DLASSQ( M, A( 1, q ), 1, T,
+ AAQQ )
SVA( q ) = T*DSQRT( AAQQ )*D( q )
END IF
END IF
IF( ( AAPP / AAPP0 )**2.LE.ROOTEPS ) THEN
IF( ( AAPP.LT.ROOTBIG ) .AND.
+ ( AAPP.GT.ROOTSFMIN ) ) THEN
AAPP = DNRM2( M, A( 1, p ), 1 )*
+ D( p )
ELSE
T = ZERO
AAPP = ZERO
CALL DLASSQ( M, A( 1, p ), 1, T,
+ AAPP )
AAPP = T*DSQRT( AAPP )*D( p )
END IF
SVA( p ) = AAPP
END IF
* end of OK rotation
ELSE
NOTROT = NOTROT + 1
PSKIPPED = PSKIPPED + 1
IJBLSK = IJBLSK + 1
END IF
ELSE
NOTROT = NOTROT + 1
PSKIPPED = PSKIPPED + 1
IJBLSK = IJBLSK + 1
END IF
*
IF( ( i.LE.SWBAND ) .AND. ( IJBLSK.GE.BLSKIP ) )
+ THEN
SVA( p ) = AAPP
NOTROT = 0
GO TO 2011
END IF
IF( ( i.LE.SWBAND ) .AND.
+ ( PSKIPPED.GT.ROWSKIP ) ) THEN
AAPP = -AAPP
NOTROT = 0
GO TO 2203
END IF
*
2200 CONTINUE
* end of the q-loop
2203 CONTINUE
*
SVA( p ) = AAPP
*
ELSE
IF( AAPP.EQ.ZERO )NOTROT = NOTROT +
+ MIN0( jgl+KBL-1, N ) - jgl + 1
IF( AAPP.LT.ZERO )NOTROT = 0
END IF
2100 CONTINUE
* end of the p-loop
2010 CONTINUE
* end of the jbc-loop
2011 CONTINUE
*2011 bailed out of the jbc-loop
DO 2012 p = igl, MIN0( igl+KBL-1, N )
SVA( p ) = DABS( SVA( p ) )
2012 CONTINUE
*
2000 CONTINUE
*2000 :: end of the ibr-loop
*
* .. update SVA(N)
IF( ( SVA( N ).LT.ROOTBIG ) .AND. ( SVA( N ).GT.ROOTSFMIN ) )
+ THEN
SVA( N ) = DNRM2( M, A( 1, N ), 1 )*D( N )
ELSE
T = ZERO
AAPP = ZERO
CALL DLASSQ( M, A( 1, N ), 1, T, AAPP )
SVA( N ) = T*DSQRT( AAPP )*D( N )
END IF
*
* Additional steering devices
*
IF( ( i.LT.SWBAND ) .AND. ( ( MXAAPQ.LE.ROOTTOL ) .OR.
+ ( ISWROT.LE.N ) ) )SWBAND = i
*
IF( ( i.GT.SWBAND+1 ) .AND. ( MXAAPQ.LT.DBLE( N )*TOL ) .AND.
+ ( DBLE( N )*MXAAPQ*MXSINJ.LT.TOL ) ) THEN
GO TO 1994
END IF
*
IF( NOTROT.GE.EMPTSW )GO TO 1994
1993 CONTINUE
* end i=1:NSWEEP loop
* #:) Reaching this point means that the procedure has comleted the given
* number of iterations.
INFO = NSWEEP - 1
GO TO 1995
1994 CONTINUE
* #:) Reaching this point means that during the i-th sweep all pivots were
* below the given tolerance, causing early exit.
*
INFO = 0
* #:) INFO = 0 confirms successful iterations.
1995 CONTINUE
*
* Sort the vector D.
DO 5991 p = 1, N - 1
q = IDAMAX( N-p+1, SVA( p ), 1 ) + p - 1
IF( p.NE.q ) THEN
TEMP1 = SVA( p )
SVA( p ) = SVA( q )
SVA( q ) = TEMP1
TEMP1 = D( p )
D( p ) = D( q )
D( q ) = TEMP1
CALL DSWAP( M, A( 1, p ), 1, A( 1, q ), 1 )
IF( RSVEC )CALL DSWAP( MVL, V( 1, p ), 1, V( 1, q ), 1 )
END IF
5991 CONTINUE
*
RETURN
* ..
* .. END OF DGSVJ0
* ..
END
|