LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zgsvj1()

subroutine zgsvj1 ( character*1  jobv,
integer  m,
integer  n,
integer  n1,
complex*16, dimension( lda, * )  a,
integer  lda,
complex*16, dimension( n )  d,
double precision, dimension( n )  sva,
integer  mv,
complex*16, dimension( ldv, * )  v,
integer  ldv,
double precision  eps,
double precision  sfmin,
double precision  tol,
integer  nsweep,
complex*16, dimension( lwork )  work,
integer  lwork,
integer  info 
)

ZGSVJ1 pre-processor for the routine zgesvj, applies Jacobi rotations targeting only particular pivots.

Download ZGSVJ1 + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 ZGSVJ1 is called from ZGESVJ as a pre-processor and that is its main
 purpose. It applies Jacobi rotations in the same way as ZGESVJ does, but
 it targets only particular pivots and it does not check convergence
 (stopping criterion). Few tuning parameters (marked by [TP]) are
 available for the implementer.

 Further Details
 ~~~~~~~~~~~~~~~
 ZGSVJ1 applies few sweeps of Jacobi rotations in the column space of
 the input M-by-N matrix A. The pivot pairs are taken from the (1,2)
 off-diagonal block in the corresponding N-by-N Gram matrix A^T * A. The
 block-entries (tiles) of the (1,2) off-diagonal block are marked by the
 [x]'s in the following scheme:

    | *  *  * [x] [x] [x]|
    | *  *  * [x] [x] [x]|    Row-cycling in the nblr-by-nblc [x] blocks.
    | *  *  * [x] [x] [x]|    Row-cyclic pivoting inside each [x] block.
    |[x] [x] [x] *  *  * |
    |[x] [x] [x] *  *  * |
    |[x] [x] [x] *  *  * |

 In terms of the columns of A, the first N1 columns are rotated 'against'
 the remaining N-N1 columns, trying to increase the angle between the
 corresponding subspaces. The off-diagonal block is N1-by(N-N1) and it is
 tiled using quadratic tiles of side KBL. Here, KBL is a tuning parameter.
 The number of sweeps is given in NSWEEP and the orthogonality threshold
 is given in TOL.
Parameters
[in]JOBV
          JOBV is 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 postmultiplying the N-by-N array V.
                (See the description of V.)
          = 'A': the product of the Jacobi rotations is accumulated
                 by postmultiplying the MV-by-N array V.
                (See the descriptions of MV and V.)
          = 'N': the Jacobi rotations are not accumulated.
[in]M
          M is INTEGER
          The number of rows of the input matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the input matrix A.
          M >= N >= 0.
[in]N1
          N1 is INTEGER
          N1 specifies the 2 x 2 block partition, the first N1 columns are
          rotated 'against' the remaining N-N1 columns of A.
[in,out]A
          A is COMPLEX*16 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 N1, D, TOL and NSWEEP.)
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[in,out]D
          D is COMPLEX*16 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 N1, A, TOL and NSWEEP.)
[in,out]SVA
          SVA is DOUBLE PRECISION 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).
[in]MV
          MV is INTEGER
          If JOBV = 'A', then MV rows of V are post-multiplied by a
                           sequence of Jacobi rotations.
          If JOBV = 'N',   then MV is not referenced.
[in,out]V
          V is COMPLEX*16 array, dimension (LDV,N)
          If JOBV = 'V' then N rows of V are post-multiplied by a
                           sequence of Jacobi rotations.
          If JOBV = 'A' then MV rows of V are post-multiplied by a
                           sequence of Jacobi rotations.
          If JOBV = 'N',   then V is not referenced.
[in]LDV
          LDV is INTEGER
          The leading dimension of the array V,  LDV >= 1.
          If JOBV = 'V', LDV >= N.
          If JOBV = 'A', LDV >= MV.
[in]EPS
          EPS is DOUBLE PRECISION
          EPS = DLAMCH('Epsilon')
[in]SFMIN
          SFMIN is DOUBLE PRECISION
          SFMIN = DLAMCH('Safe Minimum')
[in]TOL
          TOL is DOUBLE PRECISION
          TOL is the threshold for Jacobi rotations. For a pair
          A(:,p), A(:,q) of pivot columns, the Jacobi rotation is
          applied only if ABS(COS(angle(A(:,p),A(:,q)))) > TOL.
[in]NSWEEP
          NSWEEP is INTEGER
          NSWEEP is the number of sweeps of Jacobi rotations to be
          performed.
[out]WORK
          WORK is COMPLEX*16 array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          LWORK is the dimension of WORK. LWORK >= M.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, then the i-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributor:
Zlatko Drmac (Zagreb, Croatia)

Definition at line 234 of file zgsvj1.f.

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* ..
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zaxpy(n, za, zx, incx, zy, incy)
ZAXPY
Definition zaxpy.f:88
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
complex *16 function zdotc(n, zx, incx, zy, incy)
ZDOTC
Definition zdotc.f:83
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
subroutine zlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition zlascl.f:143
subroutine zlassq(n, x, incx, scale, sumsq)
ZLASSQ updates a sum of squares represented in scaled form.
Definition zlassq.f90:124
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real(wp) function dznrm2(n, x, incx)
DZNRM2
Definition dznrm2.f90:90
subroutine zrot(n, cx, incx, cy, incy, c, s)
ZROT applies a plane rotation with real cosine and complex sine to a pair of complex vectors.
Definition zrot.f:103
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
Here is the call graph for this function:
Here is the caller graph for this function: