001:       SUBROUTINE SGSVJ1( JOBV, M, N, N1, A, LDA, D, SVA, MV, V, LDV,
002:      &            EPS, SFMIN, TOL, NSWEEP, WORK, LWORK, INFO )
003: *
004: *  -- LAPACK routine (version 3.2)                                    --
005: *
006: *  -- Contributed by Zlatko Drmac of the University of Zagreb and     --
007: *  -- Kresimir Veselic of the Fernuniversitaet Hagen                  --
008: *  -- November 2008                                                   --
009: *
010: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
011: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
012: *
013: * This routine is also part of SIGMA (version 1.23, October 23. 2008.)
014: * SIGMA is a library of algorithms for highly accurate algorithms for
015: * computation of SVD, PSVD, QSVD, (H,K)-SVD, and for solution of the
016: * eigenvalue problems Hx = lambda M x, H M x = lambda x with H, M > 0.
017: *
018: *     -#- Scalar Arguments -#-
019: *
020:       IMPLICIT    NONE
021:       REAL        EPS, SFMIN, TOL
022:       INTEGER     INFO, LDA, LDV, LWORK, M, MV, N, N1, NSWEEP
023:       CHARACTER*1 JOBV
024: *
025: *     -#- Array Arguments -#-
026: *
027:       REAL        A( LDA, * ), D( N ), SVA( N ), V( LDV, * ),
028:      &            WORK( LWORK )
029: *     ..
030: *
031: *  Purpose
032: *  ~~~~~~~
033: *  SGSVJ1 is called from SGESVJ as a pre-processor and that is its main
034: *  purpose. It applies Jacobi rotations in the same way as SGESVJ does, but
035: *  it targets only particular pivots and it does not check convergence
036: *  (stopping criterion). Few tunning parameters (marked by [TP]) are
037: *  available for the implementer.
038: *
039: *  Further details
040: *  ~~~~~~~~~~~~~~~
041: *  SGSVJ1 applies few sweeps of Jacobi rotations in the column space of
042: *  the input M-by-N matrix A. The pivot pairs are taken from the (1,2)
043: *  off-diagonal block in the corresponding N-by-N Gram matrix A^T * A. The
044: *  block-entries (tiles) of the (1,2) off-diagonal block are marked by the
045: *  [x]'s in the following scheme:
046: *
047: *     | *   *   * [x] [x] [x]|
048: *     | *   *   * [x] [x] [x]|    Row-cycling in the nblr-by-nblc [x] blocks.
049: *     | *   *   * [x] [x] [x]|    Row-cyclic pivoting inside each [x] block.
050: *     |[x] [x] [x] *   *   * |
051: *     |[x] [x] [x] *   *   * |
052: *     |[x] [x] [x] *   *   * |
053: *
054: *  In terms of the columns of A, the first N1 columns are rotated 'against'
055: *  the remaining N-N1 columns, trying to increase the angle between the
056: *  corresponding subspaces. The off-diagonal block is N1-by(N-N1) and it is
057: *  tiled using quadratic tiles of side KBL. Here, KBL is a tunning parmeter.
058: *  The number of sweeps is given in NSWEEP and the orthogonality threshold
059: *  is given in TOL.
060: *
061: *  Contributors
062: *  ~~~~~~~~~~~~
063: *  Zlatko Drmac (Zagreb, Croatia) and Kresimir Veselic (Hagen, Germany)
064: *
065: *  Arguments
066: *  ~~~~~~~~~
067: *
068: *  JOBV    (input) CHARACTER*1
069: *          Specifies whether the output from this procedure is used
070: *          to compute the matrix V:
071: *          = 'V': the product of the Jacobi rotations is accumulated
072: *                 by postmulyiplying the N-by-N array V.
073: *                (See the description of V.)
074: *          = 'A': the product of the Jacobi rotations is accumulated
075: *                 by postmulyiplying the MV-by-N array V.
076: *                (See the descriptions of MV and V.)
077: *          = 'N': the Jacobi rotations are not accumulated.
078: *
079: *  M       (input) INTEGER
080: *          The number of rows of the input matrix A.  M >= 0.
081: *
082: *  N       (input) INTEGER
083: *          The number of columns of the input matrix A.
084: *          M >= N >= 0.
085: *
086: *  N1      (input) INTEGER
087: *          N1 specifies the 2 x 2 block partition, the first N1 columns are
088: *          rotated 'against' the remaining N-N1 columns of A.
089: *
090: *  A       (input/output) REAL array, dimension (LDA,N)
091: *          On entry, M-by-N matrix A, such that A*diag(D) represents
092: *          the input matrix.
093: *          On exit,
094: *          A_onexit * D_onexit represents the input matrix A*diag(D)
095: *          post-multiplied by a sequence of Jacobi rotations, where the
096: *          rotation threshold and the total number of sweeps are given in
097: *          TOL and NSWEEP, respectively.
098: *          (See the descriptions of N1, D, TOL and NSWEEP.)
099: *
100: *  LDA     (input) INTEGER
101: *          The leading dimension of the array A.  LDA >= max(1,M).
102: *
103: *  D       (input/workspace/output) REAL array, dimension (N)
104: *          The array D accumulates the scaling factors from the fast scaled
105: *          Jacobi rotations.
106: *          On entry, A*diag(D) represents the input matrix.
107: *          On exit, A_onexit*diag(D_onexit) represents the input matrix
108: *          post-multiplied by a sequence of Jacobi rotations, where the
109: *          rotation threshold and the total number of sweeps are given in
110: *          TOL and NSWEEP, respectively.
111: *          (See the descriptions of N1, A, TOL and NSWEEP.)
112: *
113: *  SVA     (input/workspace/output) REAL array, dimension (N)
114: *          On entry, SVA contains the Euclidean norms of the columns of
115: *          the matrix A*diag(D).
116: *          On exit, SVA contains the Euclidean norms of the columns of
117: *          the matrix onexit*diag(D_onexit).
118: *
119: *  MV      (input) INTEGER
120: *          If JOBV .EQ. 'A', then MV rows of V are post-multipled by a
121: *                           sequence of Jacobi rotations.
122: *          If JOBV = 'N',   then MV is not referenced.
123: *
124: *  V       (input/output) REAL array, dimension (LDV,N)
125: *          If JOBV .EQ. 'V' then N rows of V are post-multipled by a
126: *                           sequence of Jacobi rotations.
127: *          If JOBV .EQ. 'A' then MV rows of V are post-multipled by a
128: *                           sequence of Jacobi rotations.
129: *          If JOBV = 'N',   then V is not referenced.
130: *
131: *  LDV     (input) INTEGER
132: *          The leading dimension of the array V,  LDV >= 1.
133: *          If JOBV = 'V', LDV .GE. N.
134: *          If JOBV = 'A', LDV .GE. MV.
135: *
136: *  EPS     (input) INTEGER
137: *          EPS = SLAMCH('Epsilon')
138: *
139: *  SFMIN   (input) INTEGER
140: *          SFMIN = SLAMCH('Safe Minimum')
141: *
142: *  TOL     (input) REAL
143: *          TOL is the threshold for Jacobi rotations. For a pair
144: *          A(:,p), A(:,q) of pivot columns, the Jacobi rotation is
145: *          applied only if ABS(COS(angle(A(:,p),A(:,q)))) .GT. TOL.
146: *
147: *  NSWEEP  (input) INTEGER
148: *          NSWEEP is the number of sweeps of Jacobi rotations to be
149: *          performed.
150: *
151: *  WORK    (workspace) REAL array, dimension LWORK.
152: *
153: *  LWORK   (input) INTEGER
154: *          LWORK is the dimension of WORK. LWORK .GE. M.
155: *
156: *  INFO    (output) INTEGER
157: *          = 0 : successful exit.
158: *          < 0 : if INFO = -i, then the i-th argument had an illegal value
159: *
160: *     -#- Local Parameters -#-
161: *
162:       REAL        ZERO,         HALF,         ONE,         TWO
163:       PARAMETER ( ZERO = 0.0E0, HALF = 0.5E0, ONE = 1.0E0, TWO = 2.0E0 )
164: 
165: *     -#- Local Scalars -#-
166: *
167:       REAL      AAPP, AAPP0, AAPQ, AAQQ, APOAQ, AQOAP, BIG, BIGTHETA,
168:      &          CS, LARGE, MXAAPQ, MXSINJ, ROOTBIG,ROOTEPS, ROOTSFMIN,
169:      &          ROOTTOL, SMALL, SN, T, TEMP1, THETA, THSIGN
170:       INTEGER   BLSKIP, EMPTSW, i, ibr, igl, IERR, IJBLSK, ISWROT, jbc,
171:      &          jgl, KBL, MVL, NOTROT, nblc, nblr, p, PSKIPPED, q,
172:      &          ROWSKIP, SWBAND
173:       LOGICAL   APPLV, ROTOK, RSVEC
174: *
175: *     Local Arrays
176:       REAL      FASTR(5)
177: *
178: *     Intrinsic Functions
179:       INTRINSIC ABS, AMAX1, FLOAT, MIN0, SIGN, SQRT
180: *
181: *     External Functions
182:       REAL             SDOT, SNRM2
183:       INTEGER          ISAMAX
184:       LOGICAL          LSAME
185:       EXTERNAL         ISAMAX, LSAME, SDOT, SNRM2
186: *
187: *     External Subroutines
188:       EXTERNAL  SAXPY, SCOPY, SLASCL, SLASSQ, SROTM, SSWAP
189: *
190: *
191:       APPLV = LSAME(JOBV,'A')
192:       RSVEC = LSAME(JOBV,'V')
193:       IF ( .NOT.( RSVEC .OR. APPLV .OR. LSAME(JOBV,'N'))) THEN
194:          INFO = -1
195:       ELSE IF ( M .LT. 0 ) THEN
196:          INFO = -2
197:       ELSE IF ( ( N .LT. 0 ) .OR. ( N .GT. M )) THEN
198:          INFO = -3
199:       ELSE IF ( N1 .LT. 0 ) THEN
200:          INFO = -4
201:       ELSE IF ( LDA .LT. M ) THEN
202:          INFO = -6
203:       ELSE IF ( MV .LT. 0 ) THEN
204:          INFO = -9
205:       ELSE IF ( LDV .LT. M ) THEN
206:          INFO = -11
207:       ELSE IF ( TOL .LE. EPS ) THEN
208:          INFO = -14
209:       ELSE IF ( NSWEEP .LT. 0 ) THEN
210:          INFO = -15
211:       ELSE IF ( LWORK .LT. M ) THEN
212:          INFO = -17
213:       ELSE
214:          INFO = 0
215:       END IF
216: *
217: *     #:(
218:       IF ( INFO .NE. 0 ) THEN
219:          CALL XERBLA( 'SGSVJ1', -INFO )
220:          RETURN
221:       END IF
222: *
223:       IF ( RSVEC ) THEN
224:          MVL = N
225:       ELSE IF ( APPLV ) THEN
226:          MVL = MV
227:       END IF
228:       RSVEC = RSVEC .OR. APPLV
229: 
230:          ROOTEPS     = SQRT(EPS)
231:          ROOTSFMIN   = SQRT(SFMIN)
232:          SMALL       = SFMIN  / EPS
233:          BIG         = ONE   / SFMIN
234:          ROOTBIG     = ONE  / ROOTSFMIN
235:          LARGE       = BIG / SQRT(FLOAT(M*N))
236:          BIGTHETA    = ONE  / ROOTEPS
237:          ROOTTOL = SQRT(TOL)
238: *
239: *     -#- Initialize the right singular vector matrix -#-
240: *
241: *     RSVEC = LSAME( JOBV, 'Y' )
242: *
243:       EMPTSW = N1 * ( N - N1 )
244:       NOTROT     = 0
245:       FASTR(1)   = ZERO
246: *
247: *     -#- Row-cyclic pivot strategy with de Rijk's pivoting -#-
248: *
249:       KBL = MIN0(8,N)
250:       NBLR = N1 / KBL
251:       IF ( ( NBLR * KBL ) .NE. N1 ) NBLR = NBLR + 1
252: 
253: *     .. the tiling is nblr-by-nblc [tiles]
254: 
255:       NBLC = ( N - N1 ) / KBL
256:       IF ( ( NBLC * KBL ) .NE. ( N - N1 ) ) NBLC = NBLC + 1
257:       BLSKIP = ( KBL**2 ) + 1
258: *[TP] BLKSKIP is a tuning parameter that depends on SWBAND and KBL.
259: 
260:       ROWSKIP = MIN0( 5, KBL )
261: *[TP] ROWSKIP is a tuning parameter.
262:       SWBAND = 0
263: *[TP] SWBAND is a tuning parameter. It is meaningful and effective
264: *     if SGESVJ is used as a computational routine in the preconditioned
265: *     Jacobi SVD algorithm SGESVJ.
266: *
267: *
268: *     | *   *   * [x] [x] [x]|
269: *     | *   *   * [x] [x] [x]|    Row-cycling in the nblr-by-nblc [x] blocks.
270: *     | *   *   * [x] [x] [x]|    Row-cyclic pivoting inside each [x] block.
271: *     |[x] [x] [x] *   *   * |
272: *     |[x] [x] [x] *   *   * |
273: *     |[x] [x] [x] *   *   * |
274: *
275: *
276:       DO 1993 i = 1, NSWEEP
277: *     .. go go go ...
278: *
279:       MXAAPQ = ZERO
280:       MXSINJ = ZERO
281:       ISWROT = 0
282: *
283:       NOTROT = 0
284:       PSKIPPED = 0
285: *
286:       DO 2000 ibr = 1, NBLR
287: 
288:       igl = ( ibr - 1 ) * KBL + 1
289: *
290: *
291: *........................................................
292: * ... go to the off diagonal blocks
293: 
294:       igl = ( ibr - 1 ) * KBL + 1
295: 
296:       DO 2010 jbc = 1, NBLC
297: 
298:          jgl = N1 + ( jbc - 1 ) * KBL + 1
299: 
300: *        doing the block at ( ibr, jbc )
301: 
302:          IJBLSK = 0
303:          DO 2100 p = igl, MIN0( igl + KBL - 1, N1 )
304: 
305:          AAPP = SVA(p)
306: 
307:          IF ( AAPP .GT. ZERO ) THEN
308: 
309:          PSKIPPED = 0
310: 
311:          DO 2200 q = jgl, MIN0( jgl + KBL - 1, N )
312: *
313:          AAQQ = SVA(q)
314: 
315:          IF ( AAQQ .GT. ZERO ) THEN
316:             AAPP0 = AAPP
317: *
318: *     -#- M x 2 Jacobi SVD -#-
319: *
320: *        -#- Safe Gram matrix computation -#-
321: *
322:          IF ( AAQQ .GE. ONE ) THEN
323:             IF ( AAPP .GE. AAQQ ) THEN
324:                ROTOK = ( SMALL*AAPP ) .LE. AAQQ
325:             ELSE
326:                ROTOK = ( SMALL*AAQQ ) .LE. AAPP
327:             END IF
328:             IF ( AAPP .LT. ( BIG / AAQQ ) ) THEN
329:                AAPQ = ( SDOT(M, A(1,p), 1, A(1,q), 1 ) *
330:      &                  D(p) * D(q) / AAQQ ) / AAPP
331:             ELSE
332:                CALL SCOPY( M, A(1,p), 1, WORK, 1 )
333:                CALL SLASCL( 'G', 0, 0, AAPP, D(p), M,
334:      &              1, WORK, LDA, IERR )
335:                AAPQ = SDOT( M, WORK, 1, A(1,q), 1 ) *
336:      &                D(q) / AAQQ
337:             END IF
338:          ELSE
339:             IF ( AAPP .GE. AAQQ ) THEN
340:                ROTOK = AAPP .LE. ( AAQQ / SMALL )
341:             ELSE
342:                ROTOK = AAQQ .LE. ( AAPP / SMALL )
343:             END IF
344:             IF ( AAPP .GT.  ( SMALL / AAQQ ) ) THEN
345:                AAPQ = ( SDOT( M, A(1,p), 1, A(1,q), 1 ) *
346:      &               D(p) * D(q) / AAQQ ) / AAPP
347:             ELSE
348:                CALL SCOPY( M, A(1,q), 1, WORK, 1 )
349:                CALL SLASCL( 'G', 0, 0, AAQQ, D(q), M, 1,
350:      &              WORK, LDA, IERR )
351:                AAPQ = SDOT(M,WORK,1,A(1,p),1) * D(p) / AAPP
352:             END IF
353:          END IF
354: 
355:          MXAAPQ = AMAX1( MXAAPQ, ABS(AAPQ) )
356: 
357: *        TO rotate or NOT to rotate, THAT is the question ...
358: *
359:          IF ( ABS( AAPQ ) .GT. TOL ) THEN
360:             NOTROT   = 0
361: *           ROTATED  = ROTATED + 1
362:             PSKIPPED = 0
363:             ISWROT   = ISWROT  + 1
364: *
365:             IF ( ROTOK ) THEN
366: *
367:                AQOAP = AAQQ / AAPP
368:                APOAQ = AAPP / AAQQ
369:                THETA = - HALF * ABS( AQOAP - APOAQ ) / AAPQ
370:                IF ( AAQQ .GT. AAPP0 ) THETA = - THETA
371: 
372:                IF ( ABS( THETA ) .GT. BIGTHETA ) THEN
373:                   T = HALF / THETA
374:                   FASTR(3) =  T * D(p) / D(q)
375:                   FASTR(4) = -T * D(q) / D(p)
376:                   CALL SROTM( M,   A(1,p), 1, A(1,q), 1, FASTR )
377:                   IF ( RSVEC )
378:      &            CALL SROTM( MVL, V(1,p), 1, V(1,q), 1, FASTR )
379:                   SVA(q) = AAQQ*SQRT( AMAX1(ZERO,ONE + T*APOAQ*AAPQ) )
380:                   AAPP   = AAPP*SQRT( AMAX1(ZERO,ONE - T*AQOAP*AAPQ) )
381:                   MXSINJ = AMAX1( MXSINJ, ABS(T) )
382:                ELSE
383: *
384: *                 .. choose correct signum for THETA and rotate
385: *
386:                   THSIGN = - SIGN(ONE,AAPQ)
387:                   IF ( AAQQ .GT. AAPP0 ) THSIGN = - THSIGN
388:                   T  = ONE / ( THETA + THSIGN*SQRT(ONE+THETA*THETA) )
389:                   CS = SQRT( ONE / ( ONE + T*T ) )
390:                   SN = T * CS
391:                   MXSINJ = AMAX1( MXSINJ, ABS(SN) )
392:                   SVA(q) = AAQQ*SQRT( AMAX1(ZERO, ONE+T*APOAQ*AAPQ) )
393:                   AAPP   = AAPP*SQRT( ONE - T*AQOAP*AAPQ)
394: 
395:                   APOAQ = D(p) / D(q)
396:                   AQOAP = D(q) / D(p)
397:                   IF ( D(p) .GE. ONE ) THEN
398: *
399:                      IF ( D(q) .GE.  ONE ) THEN
400:                         FASTR(3) =   T * APOAQ
401:                         FASTR(4) = - T * AQOAP
402:                         D(p)  = D(p) * CS
403:                         D(q)  = D(q) * CS
404:                         CALL SROTM( M,   A(1,p),1, A(1,q),1, FASTR )
405:                         IF ( RSVEC )
406:      &                  CALL SROTM( MVL, V(1,p),1, V(1,q),1, FASTR )
407:                      ELSE
408:                         CALL SAXPY( M,    -T*AQOAP, A(1,q),1, A(1,p),1 )
409:                         CALL SAXPY( M, CS*SN*APOAQ, A(1,p),1, A(1,q),1 )
410:                         IF ( RSVEC ) THEN
411:                         CALL SAXPY( MVL, -T*AQOAP,  V(1,q),1, V(1,p),1 )
412:                         CALL SAXPY( MVL,CS*SN*APOAQ,V(1,p),1, V(1,q),1 )
413:                         END IF
414:                         D(p) = D(p) * CS
415:                         D(q) = D(q) / CS
416:                      END IF
417:                   ELSE
418:                      IF ( D(q) .GE. ONE ) THEN
419:                         CALL SAXPY( M,     T*APOAQ, A(1,p),1, A(1,q),1 )
420:                         CALL SAXPY( M,-CS*SN*AQOAP, A(1,q),1, A(1,p),1 )
421:                         IF ( RSVEC ) THEN
422:                         CALL SAXPY(MVL,T*APOAQ,     V(1,p),1, V(1,q),1 )
423:                         CALL SAXPY(MVL,-CS*SN*AQOAP,V(1,q),1, V(1,p),1 )
424:                         END IF
425:                         D(p) = D(p) / CS
426:                         D(q) = D(q) * CS
427:                      ELSE
428:                         IF ( D(p) .GE. D(q) ) THEN
429:                            CALL SAXPY( M,-T*AQOAP,   A(1,q),1,A(1,p),1 )
430:                            CALL SAXPY( M,CS*SN*APOAQ,A(1,p),1,A(1,q),1 )
431:                            D(p) = D(p) * CS
432:                            D(q) = D(q) / CS
433:                            IF ( RSVEC ) THEN
434:                            CALL SAXPY( MVL, -T*AQOAP, V(1,q),1,V(1,p),1)
435:                            CALL SAXPY(MVL,CS*SN*APOAQ,V(1,p),1,V(1,q),1)
436:                            END IF
437:                         ELSE
438:                            CALL SAXPY(M, T*APOAQ,    A(1,p),1,A(1,q),1)
439:                            CALL SAXPY(M,-CS*SN*AQOAP,A(1,q),1,A(1,p),1)
440:                            D(p) = D(p) / CS
441:                            D(q) = D(q) * CS
442:                           IF ( RSVEC ) THEN
443:                           CALL SAXPY(MVL, T*APOAQ,    V(1,p),1,V(1,q),1)
444:                           CALL SAXPY(MVL,-CS*SN*AQOAP,V(1,q),1,V(1,p),1)
445:                           END IF
446:                         END IF
447:                      END IF
448:                   ENDIF
449:                END IF
450: 
451:             ELSE
452:                IF ( AAPP .GT. AAQQ ) THEN
453:                   CALL SCOPY( M, A(1,p), 1, WORK, 1 )
454:                   CALL SLASCL('G',0,0,AAPP,ONE,M,1,WORK,LDA,IERR)
455:                   CALL SLASCL('G',0,0,AAQQ,ONE,M,1,   A(1,q),LDA,IERR)
456:                   TEMP1 = -AAPQ * D(p) / D(q)
457:                   CALL SAXPY(M,TEMP1,WORK,1,A(1,q),1)
458:                   CALL SLASCL('G',0,0,ONE,AAQQ,M,1,A(1,q),LDA,IERR)
459:                   SVA(q) = AAQQ*SQRT(AMAX1(ZERO, ONE - AAPQ*AAPQ))
460:                   MXSINJ = AMAX1( MXSINJ, SFMIN )
461:                ELSE
462:                   CALL SCOPY( M, A(1,q), 1, WORK, 1 )
463:                   CALL SLASCL('G',0,0,AAQQ,ONE,M,1,WORK,LDA,IERR)
464:                   CALL SLASCL('G',0,0,AAPP,ONE,M,1,   A(1,p),LDA,IERR)
465:                   TEMP1 = -AAPQ * D(q) / D(p)
466:                   CALL SAXPY(M,TEMP1,WORK,1,A(1,p),1)
467:                   CALL SLASCL('G',0,0,ONE,AAPP,M,1,A(1,p),LDA,IERR)
468:                   SVA(p) = AAPP*SQRT(AMAX1(ZERO, ONE - AAPQ*AAPQ))
469:                   MXSINJ = AMAX1( MXSINJ, SFMIN )
470:                END IF
471:             END IF
472: *           END IF ROTOK THEN ... ELSE
473: *
474: *           In the case of cancellation in updating SVA(q)
475: *           .. recompute SVA(q)
476:             IF ( (SVA(q) / AAQQ )**2 .LE. ROOTEPS  ) THEN
477:                IF ((AAQQ .LT. ROOTBIG).AND.(AAQQ .GT. ROOTSFMIN)) THEN
478:                   SVA(q) = SNRM2( M, A(1,q), 1 ) * D(q)
479:                ELSE
480:                   T    = ZERO
481:                   AAQQ = ZERO
482:                   CALL SLASSQ( M, A(1,q), 1, T, AAQQ )
483:                   SVA(q) = T * SQRT(AAQQ) * D(q)
484:                END IF
485:             END IF
486:             IF ( (AAPP / AAPP0 )**2 .LE. ROOTEPS  ) THEN
487:                IF ((AAPP .LT. ROOTBIG).AND.(AAPP .GT. ROOTSFMIN)) THEN
488:                   AAPP = SNRM2( M, A(1,p), 1 ) * D(p)
489:                ELSE
490:                   T    = ZERO
491:                   AAPP = ZERO
492:                   CALL SLASSQ( M, A(1,p), 1, T, AAPP )
493:                   AAPP = T * SQRT(AAPP) * D(p)
494:                END IF
495:                SVA(p) = AAPP
496:             END IF
497: *              end of OK rotation
498:          ELSE
499:             NOTROT   = NOTROT   + 1
500: *           SKIPPED  = SKIPPED  + 1
501:             PSKIPPED = PSKIPPED + 1
502:             IJBLSK   = IJBLSK   + 1
503:          END IF
504:       ELSE
505:          NOTROT   = NOTROT   + 1
506:          PSKIPPED = PSKIPPED + 1
507:          IJBLSK   = IJBLSK   + 1
508:       END IF
509: 
510: *      IF ( NOTROT .GE. EMPTSW )  GO TO 2011
511:       IF ( ( i .LE. SWBAND ) .AND. ( IJBLSK .GE. BLSKIP ) ) THEN
512:          SVA(p) = AAPP
513:          NOTROT = 0
514:          GO TO 2011
515:       END IF
516:       IF ( ( i .LE. SWBAND ) .AND. ( PSKIPPED .GT. ROWSKIP ) ) THEN
517:          AAPP = -AAPP
518:          NOTROT = 0
519:          GO TO 2203
520:       END IF
521: 
522: *
523:  2200    CONTINUE
524: *        end of the q-loop
525:  2203    CONTINUE
526: 
527:          SVA(p) = AAPP
528: *
529:       ELSE
530:          IF ( AAPP .EQ. ZERO ) NOTROT=NOTROT+MIN0(jgl+KBL-1,N)-jgl+1
531:          IF ( AAPP .LT. ZERO ) NOTROT = 0
532: ***      IF ( NOTROT .GE. EMPTSW )  GO TO 2011
533:       END IF
534: 
535:  2100 CONTINUE
536: *     end of the p-loop
537:  2010 CONTINUE
538: *     end of the jbc-loop
539:  2011 CONTINUE
540: *2011 bailed out of the jbc-loop
541:       DO 2012 p = igl, MIN0( igl + KBL - 1, N )
542:          SVA(p) = ABS(SVA(p))
543:  2012 CONTINUE
544: ***   IF ( NOTROT .GE. EMPTSW ) GO TO 1994
545:  2000 CONTINUE
546: *2000 :: end of the ibr-loop
547: *
548: *     .. update SVA(N)
549:       IF ((SVA(N) .LT. ROOTBIG).AND.(SVA(N) .GT. ROOTSFMIN)) THEN
550:          SVA(N) = SNRM2( M, A(1,N), 1 ) * D(N)
551:       ELSE
552:          T    = ZERO
553:          AAPP = ZERO
554:          CALL SLASSQ( M, A(1,N), 1, T, AAPP )
555:          SVA(N) = T * SQRT(AAPP) * D(N)
556:       END IF
557: *
558: *     Additional steering devices
559: *
560:       IF ( ( i.LT.SWBAND ) .AND. ( ( MXAAPQ.LE.ROOTTOL ) .OR.
561:      &     ( ISWROT .LE. N ) ) )
562:      &   SWBAND = i
563: 
564:       IF ((i.GT.SWBAND+1).AND. (MXAAPQ.LT.FLOAT(N)*TOL).AND.
565:      &   (FLOAT(N)*MXAAPQ*MXSINJ.LT.TOL))THEN
566:         GO TO 1994
567:       END IF
568: 
569: *
570:       IF ( NOTROT .GE. EMPTSW ) GO TO 1994
571: 
572:  1993 CONTINUE
573: *     end i=1:NSWEEP loop
574: * #:) Reaching this point means that the procedure has completed the given
575: *     number of sweeps.
576:       INFO = NSWEEP - 1
577:       GO TO 1995
578:  1994 CONTINUE
579: * #:) Reaching this point means that during the i-th sweep all pivots were
580: *     below the given threshold, causing early exit.
581: 
582:       INFO = 0
583: * #:) INFO = 0 confirms successful iterations.
584:  1995 CONTINUE
585: *
586: *     Sort the vector D
587: *
588:       DO 5991 p = 1, N - 1
589:          q = ISAMAX( N-p+1, SVA(p), 1 ) + p - 1
590:          IF ( p .NE. q ) THEN
591:             TEMP1  = SVA(p)
592:             SVA(p) = SVA(q)
593:             SVA(q) = TEMP1
594:             TEMP1   = D(p)
595:             D(p) = D(q)
596:             D(q) = TEMP1
597:             CALL SSWAP( M, A(1,p), 1, A(1,q), 1 )
598:             IF ( RSVEC ) CALL SSWAP( MVL, V(1,p), 1, V(1,q), 1 )
599:          END IF
600:  5991 CONTINUE
601: *
602:       RETURN
603: *     ..
604: *     .. END OF SGSVJ1
605: *     ..
606:       END
607: *
608: