1*> \brief \b ZGSVJ1 pre-processor for the routine zgesvj, applies Jacobi rotations targeting only particular pivots. 2* 3* =========== DOCUMENTATION =========== 4* 5* Online html documentation available at 6* http://www.netlib.org/lapack/explore-html/ 7* 8*> \htmlonly 9*> Download ZGSVJ1 + dependencies 10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgsvj1.f"> 11*> [TGZ]</a> 12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgsvj1.f"> 13*> [ZIP]</a> 14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgsvj1.f"> 15*> [TXT]</a> 16*> \endhtmlonly 17* 18* Definition: 19* =========== 20* 21* SUBROUTINE ZGSVJ1( JOBV, M, N, N1, A, LDA, D, SVA, MV, V, LDV, 22* EPS, SFMIN, TOL, NSWEEP, WORK, LWORK, INFO ) 23* 24* .. Scalar Arguments .. 25* DOUBLE PRECISION EPS, SFMIN, TOL 26* INTEGER INFO, LDA, LDV, LWORK, M, MV, N, N1, NSWEEP 27* CHARACTER*1 JOBV 28* .. 29* .. Array Arguments .. 30* COMPLEX*16 A( LDA, * ), D( N ), V( LDV, * ), WORK( LWORK ) 31* DOUBLE PRECISION SVA( N ) 32* .. 33* 34* 35*> \par Purpose: 36* ============= 37*> 38*> \verbatim 39*> 40*> ZGSVJ1 is called from ZGESVJ as a pre-processor and that is its main 41*> purpose. It applies Jacobi rotations in the same way as ZGESVJ does, but 42*> it targets only particular pivots and it does not check convergence 43*> (stopping criterion). Few tuning parameters (marked by [TP]) are 44*> available for the implementer. 45*> 46*> Further Details 47*> ~~~~~~~~~~~~~~~ 48*> ZGSVJ1 applies few sweeps of Jacobi rotations in the column space of 49*> the input M-by-N matrix A. The pivot pairs are taken from the (1,2) 50*> off-diagonal block in the corresponding N-by-N Gram matrix A^T * A. The 51*> block-entries (tiles) of the (1,2) off-diagonal block are marked by the 52*> [x]'s in the following scheme: 53*> 54*> | * * * [x] [x] [x]| 55*> | * * * [x] [x] [x]| Row-cycling in the nblr-by-nblc [x] blocks. 56*> | * * * [x] [x] [x]| Row-cyclic pivoting inside each [x] block. 57*> |[x] [x] [x] * * * | 58*> |[x] [x] [x] * * * | 59*> |[x] [x] [x] * * * | 60*> 61*> In terms of the columns of A, the first N1 columns are rotated 'against' 62*> the remaining N-N1 columns, trying to increase the angle between the 63*> corresponding subspaces. The off-diagonal block is N1-by(N-N1) and it is 64*> tiled using quadratic tiles of side KBL. Here, KBL is a tuning parameter. 65*> The number of sweeps is given in NSWEEP and the orthogonality threshold 66*> is given in TOL. 67*> \endverbatim 68* 69* Arguments: 70* ========== 71* 72*> \param[in] JOBV 73*> \verbatim 74*> JOBV is CHARACTER*1 75*> Specifies whether the output from this procedure is used 76*> to compute the matrix V: 77*> = 'V': the product of the Jacobi rotations is accumulated 78*> by postmulyiplying the N-by-N array V. 79*> (See the description of V.) 80*> = 'A': the product of the Jacobi rotations is accumulated 81*> by postmulyiplying the MV-by-N array V. 82*> (See the descriptions of MV and V.) 83*> = 'N': the Jacobi rotations are not accumulated. 84*> \endverbatim 85*> 86*> \param[in] M 87*> \verbatim 88*> M is INTEGER 89*> The number of rows of the input matrix A. M >= 0. 90*> \endverbatim 91*> 92*> \param[in] N 93*> \verbatim 94*> N is INTEGER 95*> The number of columns of the input matrix A. 96*> M >= N >= 0. 97*> \endverbatim 98*> 99*> \param[in] N1 100*> \verbatim 101*> N1 is INTEGER 102*> N1 specifies the 2 x 2 block partition, the first N1 columns are 103*> rotated 'against' the remaining N-N1 columns of A. 104*> \endverbatim 105*> 106*> \param[in,out] A 107*> \verbatim 108*> A is COMPLEX*16 array, dimension (LDA,N) 109*> On entry, M-by-N matrix A, such that A*diag(D) represents 110*> the input matrix. 111*> On exit, 112*> A_onexit * D_onexit represents the input matrix A*diag(D) 113*> post-multiplied by a sequence of Jacobi rotations, where the 114*> rotation threshold and the total number of sweeps are given in 115*> TOL and NSWEEP, respectively. 116*> (See the descriptions of N1, D, TOL and NSWEEP.) 117*> \endverbatim 118*> 119*> \param[in] LDA 120*> \verbatim 121*> LDA is INTEGER 122*> The leading dimension of the array A. LDA >= max(1,M). 123*> \endverbatim 124*> 125*> \param[in,out] D 126*> \verbatim 127*> D is COMPLEX*16 array, dimension (N) 128*> The array D accumulates the scaling factors from the fast scaled 129*> Jacobi rotations. 130*> On entry, A*diag(D) represents the input matrix. 131*> On exit, A_onexit*diag(D_onexit) represents the input matrix 132*> post-multiplied by a sequence of Jacobi rotations, where the 133*> rotation threshold and the total number of sweeps are given in 134*> TOL and NSWEEP, respectively. 135*> (See the descriptions of N1, A, TOL and NSWEEP.) 136*> \endverbatim 137*> 138*> \param[in,out] SVA 139*> \verbatim 140*> SVA is DOUBLE PRECISION array, dimension (N) 141*> On entry, SVA contains the Euclidean norms of the columns of 142*> the matrix A*diag(D). 143*> On exit, SVA contains the Euclidean norms of the columns of 144*> the matrix onexit*diag(D_onexit). 145*> \endverbatim 146*> 147*> \param[in] MV 148*> \verbatim 149*> MV is INTEGER 150*> If JOBV = 'A', then MV rows of V are post-multipled by a 151*> sequence of Jacobi rotations. 152*> If JOBV = 'N', then MV is not referenced. 153*> \endverbatim 154*> 155*> \param[in,out] V 156*> \verbatim 157*> V is COMPLEX*16 array, dimension (LDV,N) 158*> If JOBV = 'V' then N rows of V are post-multipled by a 159*> sequence of Jacobi rotations. 160*> If JOBV = 'A' then MV rows of V are post-multipled by a 161*> sequence of Jacobi rotations. 162*> If JOBV = 'N', then V is not referenced. 163*> \endverbatim 164*> 165*> \param[in] LDV 166*> \verbatim 167*> LDV is INTEGER 168*> The leading dimension of the array V, LDV >= 1. 169*> If JOBV = 'V', LDV >= N. 170*> If JOBV = 'A', LDV >= MV. 171*> \endverbatim 172*> 173*> \param[in] EPS 174*> \verbatim 175*> EPS is DOUBLE PRECISION 176*> EPS = DLAMCH('Epsilon') 177*> \endverbatim 178*> 179*> \param[in] SFMIN 180*> \verbatim 181*> SFMIN is DOUBLE PRECISION 182*> SFMIN = DLAMCH('Safe Minimum') 183*> \endverbatim 184*> 185*> \param[in] TOL 186*> \verbatim 187*> TOL is DOUBLE PRECISION 188*> TOL is the threshold for Jacobi rotations. For a pair 189*> A(:,p), A(:,q) of pivot columns, the Jacobi rotation is 190*> applied only if ABS(COS(angle(A(:,p),A(:,q)))) > TOL. 191*> \endverbatim 192*> 193*> \param[in] NSWEEP 194*> \verbatim 195*> NSWEEP is INTEGER 196*> NSWEEP is the number of sweeps of Jacobi rotations to be 197*> performed. 198*> \endverbatim 199*> 200*> \param[out] WORK 201*> \verbatim 202*> WORK is COMPLEX*16 array, dimension (LWORK) 203*> \endverbatim 204*> 205*> \param[in] LWORK 206*> \verbatim 207*> LWORK is INTEGER 208*> LWORK is the dimension of WORK. LWORK >= M. 209*> \endverbatim 210*> 211*> \param[out] INFO 212*> \verbatim 213*> INFO is INTEGER 214*> = 0: successful exit. 215*> < 0: if INFO = -i, then the i-th argument had an illegal value 216*> \endverbatim 217* 218* Authors: 219* ======== 220* 221*> \author Univ. of Tennessee 222*> \author Univ. of California Berkeley 223*> \author Univ. of Colorado Denver 224*> \author NAG Ltd. 225* 226*> \ingroup complex16OTHERcomputational 227* 228*> \par Contributor: 229* ================== 230*> 231*> Zlatko Drmac (Zagreb, Croatia) 232* 233* ===================================================================== 234 SUBROUTINE ZGSVJ1( JOBV, M, N, N1, A, LDA, D, SVA, MV, V, LDV, 235 $ EPS, SFMIN, TOL, NSWEEP, WORK, LWORK, INFO ) 236* 237* -- LAPACK computational routine -- 238* -- LAPACK is a software package provided by Univ. of Tennessee, -- 239* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 240* 241 IMPLICIT NONE 242* .. Scalar Arguments .. 243 DOUBLE PRECISION EPS, SFMIN, TOL 244 INTEGER INFO, LDA, LDV, LWORK, M, MV, N, N1, NSWEEP 245 CHARACTER*1 JOBV 246* .. 247* .. Array Arguments .. 248 COMPLEX*16 A( LDA, * ), D( N ), V( LDV, * ), WORK( LWORK ) 249 DOUBLE PRECISION SVA( N ) 250* .. 251* 252* ===================================================================== 253* 254* .. Local Parameters .. 255 DOUBLE PRECISION ZERO, HALF, ONE 256 PARAMETER ( ZERO = 0.0D0, HALF = 0.5D0, ONE = 1.0D0) 257* .. 258* .. Local Scalars .. 259 COMPLEX*16 AAPQ, OMPQ 260 DOUBLE PRECISION AAPP, AAPP0, AAPQ1, AAQQ, APOAQ, AQOAP, BIG, 261 $ BIGTHETA, CS, MXAAPQ, MXSINJ, ROOTBIG, 262 $ ROOTEPS, ROOTSFMIN, ROOTTOL, SMALL, SN, T, 263 $ TEMP1, THETA, THSIGN 264 INTEGER BLSKIP, EMPTSW, i, ibr, igl, IERR, IJBLSK, 265 $ ISWROT, jbc, jgl, KBL, MVL, NOTROT, nblc, nblr, 266 $ p, PSKIPPED, q, ROWSKIP, SWBAND 267 LOGICAL APPLV, ROTOK, RSVEC 268* .. 269* .. 270* .. Intrinsic Functions .. 271 INTRINSIC ABS, CONJG, MAX, DBLE, MIN, SIGN, SQRT 272* .. 273* .. External Functions .. 274 DOUBLE PRECISION DZNRM2 275 COMPLEX*16 ZDOTC 276 INTEGER IDAMAX 277 LOGICAL LSAME 278 EXTERNAL IDAMAX, LSAME, ZDOTC, DZNRM2 279* .. 280* .. External Subroutines .. 281* .. from BLAS 282 EXTERNAL ZCOPY, ZROT, ZSWAP, ZAXPY 283* .. from LAPACK 284 EXTERNAL ZLASCL, ZLASSQ, XERBLA 285* .. 286* .. Executable Statements .. 287* 288* Test the input parameters. 289* 290 APPLV = LSAME( JOBV, 'A' ) 291 RSVEC = LSAME( JOBV, 'V' ) 292 IF( .NOT.( RSVEC .OR. APPLV .OR. LSAME( JOBV, 'N' ) ) ) THEN 293 INFO = -1 294 ELSE IF( M.LT.0 ) THEN 295 INFO = -2 296 ELSE IF( ( N.LT.0 ) .OR. ( N.GT.M ) ) THEN 297 INFO = -3 298 ELSE IF( N1.LT.0 ) THEN 299 INFO = -4 300 ELSE IF( LDA.LT.M ) THEN 301 INFO = -6 302 ELSE IF( ( RSVEC.OR.APPLV ) .AND. ( MV.LT.0 ) ) THEN 303 INFO = -9 304 ELSE IF( ( RSVEC.AND.( LDV.LT.N ) ).OR. 305 $ ( APPLV.AND.( LDV.LT.MV ) ) ) THEN 306 INFO = -11 307 ELSE IF( TOL.LE.EPS ) THEN 308 INFO = -14 309 ELSE IF( NSWEEP.LT.0 ) THEN 310 INFO = -15 311 ELSE IF( LWORK.LT.M ) THEN 312 INFO = -17 313 ELSE 314 INFO = 0 315 END IF 316* 317* #:( 318 IF( INFO.NE.0 ) THEN 319 CALL XERBLA( 'ZGSVJ1', -INFO ) 320 RETURN 321 END IF 322* 323 IF( RSVEC ) THEN 324 MVL = N 325 ELSE IF( APPLV ) THEN 326 MVL = MV 327 END IF 328 RSVEC = RSVEC .OR. APPLV 329 330 ROOTEPS = SQRT( EPS ) 331 ROOTSFMIN = SQRT( SFMIN ) 332 SMALL = SFMIN / EPS 333 BIG = ONE / SFMIN 334 ROOTBIG = ONE / ROOTSFMIN 335* LARGE = BIG / SQRT( DBLE( M*N ) ) 336 BIGTHETA = ONE / ROOTEPS 337 ROOTTOL = SQRT( TOL ) 338* 339* .. Initialize the right singular vector matrix .. 340* 341* RSVEC = LSAME( JOBV, 'Y' ) 342* 343 EMPTSW = N1*( N-N1 ) 344 NOTROT = 0 345* 346* .. Row-cyclic pivot strategy with de Rijk's pivoting .. 347* 348 KBL = MIN( 8, N ) 349 NBLR = N1 / KBL 350 IF( ( NBLR*KBL ).NE.N1 )NBLR = NBLR + 1 351 352* .. the tiling is nblr-by-nblc [tiles] 353 354 NBLC = ( N-N1 ) / KBL 355 IF( ( NBLC*KBL ).NE.( N-N1 ) )NBLC = NBLC + 1 356 BLSKIP = ( KBL**2 ) + 1 357*[TP] BLKSKIP is a tuning parameter that depends on SWBAND and KBL. 358 359 ROWSKIP = MIN( 5, KBL ) 360*[TP] ROWSKIP is a tuning parameter. 361 SWBAND = 0 362*[TP] SWBAND is a tuning parameter. It is meaningful and effective 363* if ZGESVJ is used as a computational routine in the preconditioned 364* Jacobi SVD algorithm ZGEJSV. 365* 366* 367* | * * * [x] [x] [x]| 368* | * * * [x] [x] [x]| Row-cycling in the nblr-by-nblc [x] blocks. 369* | * * * [x] [x] [x]| Row-cyclic pivoting inside each [x] block. 370* |[x] [x] [x] * * * | 371* |[x] [x] [x] * * * | 372* |[x] [x] [x] * * * | 373* 374* 375 DO 1993 i = 1, NSWEEP 376* 377* .. go go go ... 378* 379 MXAAPQ = ZERO 380 MXSINJ = ZERO 381 ISWROT = 0 382* 383 NOTROT = 0 384 PSKIPPED = 0 385* 386* Each sweep is unrolled using KBL-by-KBL tiles over the pivot pairs 387* 1 <= p < q <= N. This is the first step toward a blocked implementation 388* of the rotations. New implementation, based on block transformations, 389* is under development. 390* 391 DO 2000 ibr = 1, NBLR 392* 393 igl = ( ibr-1 )*KBL + 1 394* 395 396* 397* ... go to the off diagonal blocks 398* 399 igl = ( ibr-1 )*KBL + 1 400* 401* DO 2010 jbc = ibr + 1, NBL 402 DO 2010 jbc = 1, NBLC 403* 404 jgl = ( jbc-1 )*KBL + N1 + 1 405* 406* doing the block at ( ibr, jbc ) 407* 408 IJBLSK = 0 409 DO 2100 p = igl, MIN( igl+KBL-1, N1 ) 410* 411 AAPP = SVA( p ) 412 IF( AAPP.GT.ZERO ) THEN 413* 414 PSKIPPED = 0 415* 416 DO 2200 q = jgl, MIN( jgl+KBL-1, N ) 417* 418 AAQQ = SVA( q ) 419 IF( AAQQ.GT.ZERO ) THEN 420 AAPP0 = AAPP 421* 422* .. M x 2 Jacobi SVD .. 423* 424* Safe Gram matrix computation 425* 426 IF( AAQQ.GE.ONE ) THEN 427 IF( AAPP.GE.AAQQ ) THEN 428 ROTOK = ( SMALL*AAPP ).LE.AAQQ 429 ELSE 430 ROTOK = ( SMALL*AAQQ ).LE.AAPP 431 END IF 432 IF( AAPP.LT.( BIG / AAQQ ) ) THEN 433 AAPQ = ( ZDOTC( M, A( 1, p ), 1, 434 $ A( 1, q ), 1 ) / AAQQ ) / AAPP 435 ELSE 436 CALL ZCOPY( M, A( 1, p ), 1, 437 $ WORK, 1 ) 438 CALL ZLASCL( 'G', 0, 0, AAPP, 439 $ ONE, M, 1, 440 $ WORK, LDA, IERR ) 441 AAPQ = ZDOTC( M, WORK, 1, 442 $ A( 1, q ), 1 ) / AAQQ 443 END IF 444 ELSE 445 IF( AAPP.GE.AAQQ ) THEN 446 ROTOK = AAPP.LE.( AAQQ / SMALL ) 447 ELSE 448 ROTOK = AAQQ.LE.( AAPP / SMALL ) 449 END IF 450 IF( AAPP.GT.( SMALL / AAQQ ) ) THEN 451 AAPQ = ( ZDOTC( M, A( 1, p ), 1, 452 $ A( 1, q ), 1 ) / MAX(AAQQ,AAPP) ) 453 $ / MIN(AAQQ,AAPP) 454 ELSE 455 CALL ZCOPY( M, A( 1, q ), 1, 456 $ WORK, 1 ) 457 CALL ZLASCL( 'G', 0, 0, AAQQ, 458 $ ONE, M, 1, 459 $ WORK, LDA, IERR ) 460 AAPQ = ZDOTC( M, A( 1, p ), 1, 461 $ WORK, 1 ) / AAPP 462 END IF 463 END IF 464* 465* AAPQ = AAPQ * CONJG(CWORK(p))*CWORK(q) 466 AAPQ1 = -ABS(AAPQ) 467 MXAAPQ = MAX( MXAAPQ, -AAPQ1 ) 468* 469* TO rotate or NOT to rotate, THAT is the question ... 470* 471 IF( ABS( AAPQ1 ).GT.TOL ) THEN 472 OMPQ = AAPQ / ABS(AAPQ) 473 NOTROT = 0 474*[RTD] ROTATED = ROTATED + 1 475 PSKIPPED = 0 476 ISWROT = ISWROT + 1 477* 478 IF( ROTOK ) THEN 479* 480 AQOAP = AAQQ / AAPP 481 APOAQ = AAPP / AAQQ 482 THETA = -HALF*ABS( AQOAP-APOAQ )/ AAPQ1 483 IF( AAQQ.GT.AAPP0 )THETA = -THETA 484* 485 IF( ABS( THETA ).GT.BIGTHETA ) THEN 486 T = HALF / THETA 487 CS = ONE 488 CALL ZROT( M, A(1,p), 1, A(1,q), 1, 489 $ CS, CONJG(OMPQ)*T ) 490 IF( RSVEC ) THEN 491 CALL ZROT( MVL, V(1,p), 1, 492 $ V(1,q), 1, CS, CONJG(OMPQ)*T ) 493 END IF 494 SVA( q ) = AAQQ*SQRT( MAX( ZERO, 495 $ ONE+T*APOAQ*AAPQ1 ) ) 496 AAPP = AAPP*SQRT( MAX( ZERO, 497 $ ONE-T*AQOAP*AAPQ1 ) ) 498 MXSINJ = MAX( MXSINJ, ABS( T ) ) 499 ELSE 500* 501* .. choose correct signum for THETA and rotate 502* 503 THSIGN = -SIGN( ONE, AAPQ1 ) 504 IF( AAQQ.GT.AAPP0 )THSIGN = -THSIGN 505 T = ONE / ( THETA+THSIGN* 506 $ SQRT( ONE+THETA*THETA ) ) 507 CS = SQRT( ONE / ( ONE+T*T ) ) 508 SN = T*CS 509 MXSINJ = MAX( MXSINJ, ABS( SN ) ) 510 SVA( q ) = AAQQ*SQRT( MAX( ZERO, 511 $ ONE+T*APOAQ*AAPQ1 ) ) 512 AAPP = AAPP*SQRT( MAX( ZERO, 513 $ ONE-T*AQOAP*AAPQ1 ) ) 514* 515 CALL ZROT( M, A(1,p), 1, A(1,q), 1, 516 $ CS, CONJG(OMPQ)*SN ) 517 IF( RSVEC ) THEN 518 CALL ZROT( MVL, V(1,p), 1, 519 $ V(1,q), 1, CS, CONJG(OMPQ)*SN ) 520 END IF 521 END IF 522 D(p) = -D(q) * OMPQ 523* 524 ELSE 525* .. have to use modified Gram-Schmidt like transformation 526 IF( AAPP.GT.AAQQ ) THEN 527 CALL ZCOPY( M, A( 1, p ), 1, 528 $ WORK, 1 ) 529 CALL ZLASCL( 'G', 0, 0, AAPP, ONE, 530 $ M, 1, WORK,LDA, 531 $ IERR ) 532 CALL ZLASCL( 'G', 0, 0, AAQQ, ONE, 533 $ M, 1, A( 1, q ), LDA, 534 $ IERR ) 535 CALL ZAXPY( M, -AAPQ, WORK, 536 $ 1, A( 1, q ), 1 ) 537 CALL ZLASCL( 'G', 0, 0, ONE, AAQQ, 538 $ M, 1, A( 1, q ), LDA, 539 $ IERR ) 540 SVA( q ) = AAQQ*SQRT( MAX( ZERO, 541 $ ONE-AAPQ1*AAPQ1 ) ) 542 MXSINJ = MAX( MXSINJ, SFMIN ) 543 ELSE 544 CALL ZCOPY( M, A( 1, q ), 1, 545 $ WORK, 1 ) 546 CALL ZLASCL( 'G', 0, 0, AAQQ, ONE, 547 $ M, 1, WORK,LDA, 548 $ IERR ) 549 CALL ZLASCL( 'G', 0, 0, AAPP, ONE, 550 $ M, 1, A( 1, p ), LDA, 551 $ IERR ) 552 CALL ZAXPY( M, -CONJG(AAPQ), 553 $ WORK, 1, A( 1, p ), 1 ) 554 CALL ZLASCL( 'G', 0, 0, ONE, AAPP, 555 $ M, 1, A( 1, p ), LDA, 556 $ IERR ) 557 SVA( p ) = AAPP*SQRT( MAX( ZERO, 558 $ ONE-AAPQ1*AAPQ1 ) ) 559 MXSINJ = MAX( MXSINJ, SFMIN ) 560 END IF 561 END IF 562* END IF ROTOK THEN ... ELSE 563* 564* In the case of cancellation in updating SVA(q), SVA(p) 565* .. recompute SVA(q), SVA(p) 566 IF( ( SVA( q ) / AAQQ )**2.LE.ROOTEPS ) 567 $ THEN 568 IF( ( AAQQ.LT.ROOTBIG ) .AND. 569 $ ( AAQQ.GT.ROOTSFMIN ) ) THEN 570 SVA( q ) = DZNRM2( M, A( 1, q ), 1) 571 ELSE 572 T = ZERO 573 AAQQ = ONE 574 CALL ZLASSQ( M, A( 1, q ), 1, T, 575 $ AAQQ ) 576 SVA( q ) = T*SQRT( AAQQ ) 577 END IF 578 END IF 579 IF( ( AAPP / AAPP0 )**2.LE.ROOTEPS ) THEN 580 IF( ( AAPP.LT.ROOTBIG ) .AND. 581 $ ( AAPP.GT.ROOTSFMIN ) ) THEN 582 AAPP = DZNRM2( M, A( 1, p ), 1 ) 583 ELSE 584 T = ZERO 585 AAPP = ONE 586 CALL ZLASSQ( M, A( 1, p ), 1, T, 587 $ AAPP ) 588 AAPP = T*SQRT( AAPP ) 589 END IF 590 SVA( p ) = AAPP 591 END IF 592* end of OK rotation 593 ELSE 594 NOTROT = NOTROT + 1 595*[RTD] SKIPPED = SKIPPED + 1 596 PSKIPPED = PSKIPPED + 1 597 IJBLSK = IJBLSK + 1 598 END IF 599 ELSE 600 NOTROT = NOTROT + 1 601 PSKIPPED = PSKIPPED + 1 602 IJBLSK = IJBLSK + 1 603 END IF 604* 605 IF( ( i.LE.SWBAND ) .AND. ( IJBLSK.GE.BLSKIP ) ) 606 $ THEN 607 SVA( p ) = AAPP 608 NOTROT = 0 609 GO TO 2011 610 END IF 611 IF( ( i.LE.SWBAND ) .AND. 612 $ ( PSKIPPED.GT.ROWSKIP ) ) THEN 613 AAPP = -AAPP 614 NOTROT = 0 615 GO TO 2203 616 END IF 617* 618 2200 CONTINUE 619* end of the q-loop 620 2203 CONTINUE 621* 622 SVA( p ) = AAPP 623* 624 ELSE 625* 626 IF( AAPP.EQ.ZERO )NOTROT = NOTROT + 627 $ MIN( jgl+KBL-1, N ) - jgl + 1 628 IF( AAPP.LT.ZERO )NOTROT = 0 629* 630 END IF 631* 632 2100 CONTINUE 633* end of the p-loop 634 2010 CONTINUE 635* end of the jbc-loop 636 2011 CONTINUE 637*2011 bailed out of the jbc-loop 638 DO 2012 p = igl, MIN( igl+KBL-1, N ) 639 SVA( p ) = ABS( SVA( p ) ) 640 2012 CONTINUE 641*** 642 2000 CONTINUE 643*2000 :: end of the ibr-loop 644* 645* .. update SVA(N) 646 IF( ( SVA( N ).LT.ROOTBIG ) .AND. ( SVA( N ).GT.ROOTSFMIN ) ) 647 $ THEN 648 SVA( N ) = DZNRM2( M, A( 1, N ), 1 ) 649 ELSE 650 T = ZERO 651 AAPP = ONE 652 CALL ZLASSQ( M, A( 1, N ), 1, T, AAPP ) 653 SVA( N ) = T*SQRT( AAPP ) 654 END IF 655* 656* Additional steering devices 657* 658 IF( ( i.LT.SWBAND ) .AND. ( ( MXAAPQ.LE.ROOTTOL ) .OR. 659 $ ( ISWROT.LE.N ) ) )SWBAND = i 660* 661 IF( ( i.GT.SWBAND+1 ) .AND. ( MXAAPQ.LT.SQRT( DBLE( N ) )* 662 $ TOL ) .AND. ( DBLE( N )*MXAAPQ*MXSINJ.LT.TOL ) ) THEN 663 GO TO 1994 664 END IF 665* 666 IF( NOTROT.GE.EMPTSW )GO TO 1994 667* 668 1993 CONTINUE 669* end i=1:NSWEEP loop 670* 671* #:( Reaching this point means that the procedure has not converged. 672 INFO = NSWEEP - 1 673 GO TO 1995 674* 675 1994 CONTINUE 676* #:) Reaching this point means numerical convergence after the i-th 677* sweep. 678* 679 INFO = 0 680* #:) INFO = 0 confirms successful iterations. 681 1995 CONTINUE 682* 683* Sort the vector SVA() of column norms. 684 DO 5991 p = 1, N - 1 685 q = IDAMAX( N-p+1, SVA( p ), 1 ) + p - 1 686 IF( p.NE.q ) THEN 687 TEMP1 = SVA( p ) 688 SVA( p ) = SVA( q ) 689 SVA( q ) = TEMP1 690 AAPQ = D( p ) 691 D( p ) = D( q ) 692 D( q ) = AAPQ 693 CALL ZSWAP( M, A( 1, p ), 1, A( 1, q ), 1 ) 694 IF( RSVEC )CALL ZSWAP( MVL, V( 1, p ), 1, V( 1, q ), 1 ) 695 END IF 696 5991 CONTINUE 697* 698* 699 RETURN 700* .. 701* .. END OF ZGSVJ1 702* .. 703 END 704