LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
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 tunning 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 tunning parmeter.
 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 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.
[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 .EQ. 'A', then MV rows of V are post-multipled 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 .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.
[in]LDV
          LDV is INTEGER
          The leading dimension of the array V,  LDV >= 1.
          If JOBV = 'V', LDV .GE. N.
          If JOBV = 'A', LDV .GE. 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)))) .GT. 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 .GE. 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.
Date
June 2016
Contributor:
Zlatko Drmac (Zagreb, Croatia)

Definition at line 238 of file zgsvj1.f.

238 *
239 * -- LAPACK computational routine (version 3.7.0) --
240 * -- LAPACK is a software package provided by Univ. of Tennessee, --
241 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
242 * June 2016
243 *
244  IMPLICIT NONE
245 * .. Scalar Arguments ..
246  DOUBLE PRECISION eps, sfmin, tol
247  INTEGER info, lda, ldv, lwork, m, mv, n, n1, nsweep
248  CHARACTER*1 jobv
249 * ..
250 * .. Array Arguments ..
251  COMPLEX*16 a( lda, * ), d( n ), v( ldv, * ), work( lwork )
252  DOUBLE PRECISION sva( n )
253 * ..
254 *
255 * =====================================================================
256 *
257 * .. Local Parameters ..
258  DOUBLE PRECISION zero, half, one
259  parameter ( zero = 0.0d0, half = 0.5d0, one = 1.0d0)
260 * ..
261 * .. Local Scalars ..
262  COMPLEX*16 aapq, ompq
263  DOUBLE PRECISION aapp, aapp0, aapq1, aaqq, apoaq, aqoap, big,
264  $ bigtheta, cs, mxaapq, mxsinj, rootbig,
265  $ rooteps, rootsfmin, roottol, small, sn, t,
266  $ temp1, theta, thsign
267  INTEGER blskip, emptsw, i, ibr, igl, ierr, ijblsk,
268  $ iswrot, jbc, jgl, kbl, mvl, notrot, nblc, nblr,
269  $ p, pskipped, q, rowskip, swband
270  LOGICAL applv, rotok, rsvec
271 * ..
272 * ..
273 * .. Intrinsic Functions ..
274  INTRINSIC abs, conjg, max, dble, min, sign, sqrt
275 * ..
276 * .. External Functions ..
277  DOUBLE PRECISION dznrm2
278  COMPLEX*16 zdotc
279  INTEGER idamax
280  LOGICAL lsame
281  EXTERNAL idamax, lsame, zdotc, dznrm2
282 * ..
283 * .. External Subroutines ..
284 * .. from BLAS
285  EXTERNAL zcopy, zrot, zswap
286 * .. from LAPACK
287  EXTERNAL zlascl, zlassq, xerbla
288 * ..
289 * .. Executable Statements ..
290 *
291 * Test the input parameters.
292 *
293  applv = lsame( jobv, 'A' )
294  rsvec = lsame( jobv, 'V' )
295  IF( .NOT.( rsvec .OR. applv .OR. lsame( jobv, 'N' ) ) ) THEN
296  info = -1
297  ELSE IF( m.LT.0 ) THEN
298  info = -2
299  ELSE IF( ( n.LT.0 ) .OR. ( n.GT.m ) ) THEN
300  info = -3
301  ELSE IF( n1.LT.0 ) THEN
302  info = -4
303  ELSE IF( lda.LT.m ) THEN
304  info = -6
305  ELSE IF( ( rsvec.OR.applv ) .AND. ( mv.LT.0 ) ) THEN
306  info = -9
307  ELSE IF( ( rsvec.AND.( ldv.LT.n ) ).OR.
308  $ ( applv.AND.( ldv.LT.mv ) ) ) THEN
309  info = -11
310  ELSE IF( tol.LE.eps ) THEN
311  info = -14
312  ELSE IF( nsweep.LT.0 ) THEN
313  info = -15
314  ELSE IF( lwork.LT.m ) THEN
315  info = -17
316  ELSE
317  info = 0
318  END IF
319 *
320 * #:(
321  IF( info.NE.0 ) THEN
322  CALL xerbla( 'ZGSVJ1', -info )
323  RETURN
324  END IF
325 *
326  IF( rsvec ) THEN
327  mvl = n
328  ELSE IF( applv ) THEN
329  mvl = mv
330  END IF
331  rsvec = rsvec .OR. applv
332 
333  rooteps = sqrt( eps )
334  rootsfmin = sqrt( sfmin )
335  small = sfmin / eps
336  big = one / sfmin
337  rootbig = one / rootsfmin
338 * LARGE = BIG / SQRT( DBLE( M*N ) )
339  bigtheta = one / rooteps
340  roottol = sqrt( tol )
341 *
342 * .. Initialize the right singular vector matrix ..
343 *
344 * RSVEC = LSAME( JOBV, 'Y' )
345 *
346  emptsw = n1*( n-n1 )
347  notrot = 0
348 *
349 * .. Row-cyclic pivot strategy with de Rijk's pivoting ..
350 *
351  kbl = min( 8, n )
352  nblr = n1 / kbl
353  IF( ( nblr*kbl ).NE.n1 )nblr = nblr + 1
354 
355 * .. the tiling is nblr-by-nblc [tiles]
356 
357  nblc = ( n-n1 ) / kbl
358  IF( ( nblc*kbl ).NE.( n-n1 ) )nblc = nblc + 1
359  blskip = ( kbl**2 ) + 1
360 *[TP] BLKSKIP is a tuning parameter that depends on SWBAND and KBL.
361 
362  rowskip = min( 5, kbl )
363 *[TP] ROWSKIP is a tuning parameter.
364  swband = 0
365 *[TP] SWBAND is a tuning parameter. It is meaningful and effective
366 * if ZGESVJ is used as a computational routine in the preconditioned
367 * Jacobi SVD algorithm ZGEJSV.
368 *
369 *
370 * | * * * [x] [x] [x]|
371 * | * * * [x] [x] [x]| Row-cycling in the nblr-by-nblc [x] blocks.
372 * | * * * [x] [x] [x]| Row-cyclic pivoting inside each [x] block.
373 * |[x] [x] [x] * * * |
374 * |[x] [x] [x] * * * |
375 * |[x] [x] [x] * * * |
376 *
377 *
378  DO 1993 i = 1, nsweep
379 *
380 * .. go go go ...
381 *
382  mxaapq = zero
383  mxsinj = zero
384  iswrot = 0
385 *
386  notrot = 0
387  pskipped = 0
388 *
389 * Each sweep is unrolled using KBL-by-KBL tiles over the pivot pairs
390 * 1 <= p < q <= N. This is the first step toward a blocked implementation
391 * of the rotations. New implementation, based on block transformations,
392 * is under development.
393 *
394  DO 2000 ibr = 1, nblr
395 *
396  igl = ( ibr-1 )*kbl + 1
397 *
398 
399 *
400 * ... go to the off diagonal blocks
401 *
402  igl = ( ibr-1 )*kbl + 1
403 *
404 * DO 2010 jbc = ibr + 1, NBL
405  DO 2010 jbc = 1, nblc
406 *
407  jgl = ( jbc-1 )*kbl + n1 + 1
408 *
409 * doing the block at ( ibr, jbc )
410 *
411  ijblsk = 0
412  DO 2100 p = igl, min( igl+kbl-1, n1 )
413 *
414  aapp = sva( p )
415  IF( aapp.GT.zero ) THEN
416 *
417  pskipped = 0
418 *
419  DO 2200 q = jgl, min( jgl+kbl-1, n )
420 *
421  aaqq = sva( q )
422  IF( aaqq.GT.zero ) THEN
423  aapp0 = aapp
424 *
425 * .. M x 2 Jacobi SVD ..
426 *
427 * Safe Gram matrix computation
428 *
429  IF( aaqq.GE.one ) THEN
430  IF( aapp.GE.aaqq ) THEN
431  rotok = ( small*aapp ).LE.aaqq
432  ELSE
433  rotok = ( small*aaqq ).LE.aapp
434  END IF
435  IF( aapp.LT.( big / aaqq ) ) THEN
436  aapq = ( zdotc( m, a( 1, p ), 1,
437  $ a( 1, q ), 1 ) / aaqq ) / aapp
438  ELSE
439  CALL zcopy( m, a( 1, p ), 1,
440  $ work, 1 )
441  CALL zlascl( 'G', 0, 0, aapp,
442  $ one, m, 1,
443  $ work, lda, ierr )
444  aapq = zdotc( m, work, 1,
445  $ a( 1, q ), 1 ) / aaqq
446  END IF
447  ELSE
448  IF( aapp.GE.aaqq ) THEN
449  rotok = aapp.LE.( aaqq / small )
450  ELSE
451  rotok = aaqq.LE.( aapp / small )
452  END IF
453  IF( aapp.GT.( small / aaqq ) ) THEN
454  aapq = ( zdotc( m, a( 1, p ), 1,
455  $ a( 1, q ), 1 ) / max(aaqq,aapp) )
456  $ / min(aaqq,aapp)
457  ELSE
458  CALL zcopy( m, a( 1, q ), 1,
459  $ work, 1 )
460  CALL zlascl( 'G', 0, 0, aaqq,
461  $ one, m, 1,
462  $ work, lda, ierr )
463  aapq = zdotc( m, a( 1, p ), 1,
464  $ work, 1 ) / aapp
465  END IF
466  END IF
467 *
468 * AAPQ = AAPQ * CONJG(CWORK(p))*CWORK(q)
469  aapq1 = -abs(aapq)
470  mxaapq = max( mxaapq, -aapq1 )
471 *
472 * TO rotate or NOT to rotate, THAT is the question ...
473 *
474  IF( abs( aapq1 ).GT.tol ) THEN
475  ompq = aapq / abs(aapq)
476  notrot = 0
477 *[RTD] ROTATED = ROTATED + 1
478  pskipped = 0
479  iswrot = iswrot + 1
480 *
481  IF( rotok ) THEN
482 *
483  aqoap = aaqq / aapp
484  apoaq = aapp / aaqq
485  theta = -half*abs( aqoap-apoaq )/ aapq1
486  IF( aaqq.GT.aapp0 )theta = -theta
487 *
488  IF( abs( theta ).GT.bigtheta ) THEN
489  t = half / theta
490  cs = one
491  CALL zrot( m, a(1,p), 1, a(1,q), 1,
492  $ cs, conjg(ompq)*t )
493  IF( rsvec ) THEN
494  CALL zrot( mvl, v(1,p), 1,
495  $ v(1,q), 1, cs, conjg(ompq)*t )
496  END IF
497  sva( q ) = aaqq*sqrt( max( zero,
498  $ one+t*apoaq*aapq1 ) )
499  aapp = aapp*sqrt( max( zero,
500  $ one-t*aqoap*aapq1 ) )
501  mxsinj = max( mxsinj, abs( t ) )
502  ELSE
503 *
504 * .. choose correct signum for THETA and rotate
505 *
506  thsign = -sign( one, aapq1 )
507  IF( aaqq.GT.aapp0 )thsign = -thsign
508  t = one / ( theta+thsign*
509  $ sqrt( one+theta*theta ) )
510  cs = sqrt( one / ( one+t*t ) )
511  sn = t*cs
512  mxsinj = max( mxsinj, abs( sn ) )
513  sva( q ) = aaqq*sqrt( max( zero,
514  $ one+t*apoaq*aapq1 ) )
515  aapp = aapp*sqrt( max( zero,
516  $ one-t*aqoap*aapq1 ) )
517 *
518  CALL zrot( m, a(1,p), 1, a(1,q), 1,
519  $ cs, conjg(ompq)*sn )
520  IF( rsvec ) THEN
521  CALL zrot( mvl, v(1,p), 1,
522  $ v(1,q), 1, cs, conjg(ompq)*sn )
523  END IF
524  END IF
525  d(p) = -d(q) * ompq
526 *
527  ELSE
528 * .. have to use modified Gram-Schmidt like transformation
529  IF( aapp.GT.aaqq ) THEN
530  CALL zcopy( m, a( 1, p ), 1,
531  $ work, 1 )
532  CALL zlascl( 'G', 0, 0, aapp, one,
533  $ m, 1, work,lda,
534  $ ierr )
535  CALL zlascl( 'G', 0, 0, aaqq, one,
536  $ m, 1, a( 1, q ), lda,
537  $ ierr )
538  CALL zaxpy( m, -aapq, work,
539  $ 1, a( 1, q ), 1 )
540  CALL zlascl( 'G', 0, 0, one, aaqq,
541  $ m, 1, a( 1, q ), lda,
542  $ ierr )
543  sva( q ) = aaqq*sqrt( max( zero,
544  $ one-aapq1*aapq1 ) )
545  mxsinj = max( mxsinj, sfmin )
546  ELSE
547  CALL zcopy( m, a( 1, q ), 1,
548  $ work, 1 )
549  CALL zlascl( 'G', 0, 0, aaqq, one,
550  $ m, 1, work,lda,
551  $ ierr )
552  CALL zlascl( 'G', 0, 0, aapp, one,
553  $ m, 1, a( 1, p ), lda,
554  $ ierr )
555  CALL zaxpy( m, -conjg(aapq),
556  $ work, 1, a( 1, p ), 1 )
557  CALL zlascl( 'G', 0, 0, one, aapp,
558  $ m, 1, a( 1, p ), lda,
559  $ ierr )
560  sva( p ) = aapp*sqrt( max( zero,
561  $ one-aapq1*aapq1 ) )
562  mxsinj = max( mxsinj, sfmin )
563  END IF
564  END IF
565 * END IF ROTOK THEN ... ELSE
566 *
567 * In the case of cancellation in updating SVA(q), SVA(p)
568 * .. recompute SVA(q), SVA(p)
569  IF( ( sva( q ) / aaqq )**2.LE.rooteps )
570  $ THEN
571  IF( ( aaqq.LT.rootbig ) .AND.
572  $ ( aaqq.GT.rootsfmin ) ) THEN
573  sva( q ) = dznrm2( m, a( 1, q ), 1)
574  ELSE
575  t = zero
576  aaqq = one
577  CALL zlassq( m, a( 1, q ), 1, t,
578  $ aaqq )
579  sva( q ) = t*sqrt( aaqq )
580  END IF
581  END IF
582  IF( ( aapp / aapp0 )**2.LE.rooteps ) THEN
583  IF( ( aapp.LT.rootbig ) .AND.
584  $ ( aapp.GT.rootsfmin ) ) THEN
585  aapp = dznrm2( m, a( 1, p ), 1 )
586  ELSE
587  t = zero
588  aapp = one
589  CALL zlassq( m, a( 1, p ), 1, t,
590  $ aapp )
591  aapp = t*sqrt( aapp )
592  END IF
593  sva( p ) = aapp
594  END IF
595 * end of OK rotation
596  ELSE
597  notrot = notrot + 1
598 *[RTD] SKIPPED = SKIPPED + 1
599  pskipped = pskipped + 1
600  ijblsk = ijblsk + 1
601  END IF
602  ELSE
603  notrot = notrot + 1
604  pskipped = pskipped + 1
605  ijblsk = ijblsk + 1
606  END IF
607 *
608  IF( ( i.LE.swband ) .AND. ( ijblsk.GE.blskip ) )
609  $ THEN
610  sva( p ) = aapp
611  notrot = 0
612  GO TO 2011
613  END IF
614  IF( ( i.LE.swband ) .AND.
615  $ ( pskipped.GT.rowskip ) ) THEN
616  aapp = -aapp
617  notrot = 0
618  GO TO 2203
619  END IF
620 *
621  2200 CONTINUE
622 * end of the q-loop
623  2203 CONTINUE
624 *
625  sva( p ) = aapp
626 *
627  ELSE
628 *
629  IF( aapp.EQ.zero )notrot = notrot +
630  $ min( jgl+kbl-1, n ) - jgl + 1
631  IF( aapp.LT.zero )notrot = 0
632 *
633  END IF
634 *
635  2100 CONTINUE
636 * end of the p-loop
637  2010 CONTINUE
638 * end of the jbc-loop
639  2011 CONTINUE
640 *2011 bailed out of the jbc-loop
641  DO 2012 p = igl, min( igl+kbl-1, n )
642  sva( p ) = abs( sva( p ) )
643  2012 CONTINUE
644 ***
645  2000 CONTINUE
646 *2000 :: end of the ibr-loop
647 *
648 * .. update SVA(N)
649  IF( ( sva( n ).LT.rootbig ) .AND. ( sva( n ).GT.rootsfmin ) )
650  $ THEN
651  sva( n ) = dznrm2( m, a( 1, n ), 1 )
652  ELSE
653  t = zero
654  aapp = one
655  CALL zlassq( m, a( 1, n ), 1, t, aapp )
656  sva( n ) = t*sqrt( aapp )
657  END IF
658 *
659 * Additional steering devices
660 *
661  IF( ( i.LT.swband ) .AND. ( ( mxaapq.LE.roottol ) .OR.
662  $ ( iswrot.LE.n ) ) )swband = i
663 *
664  IF( ( i.GT.swband+1 ) .AND. ( mxaapq.LT.sqrt( dble( n ) )*
665  $ tol ) .AND. ( dble( n )*mxaapq*mxsinj.LT.tol ) ) THEN
666  GO TO 1994
667  END IF
668 *
669  IF( notrot.GE.emptsw )GO TO 1994
670 *
671  1993 CONTINUE
672 * end i=1:NSWEEP loop
673 *
674 * #:( Reaching this point means that the procedure has not converged.
675  info = nsweep - 1
676  GO TO 1995
677 *
678  1994 CONTINUE
679 * #:) Reaching this point means numerical convergence after the i-th
680 * sweep.
681 *
682  info = 0
683 * #:) INFO = 0 confirms successful iterations.
684  1995 CONTINUE
685 *
686 * Sort the vector SVA() of column norms.
687  DO 5991 p = 1, n - 1
688  q = idamax( n-p+1, sva( p ), 1 ) + p - 1
689  IF( p.NE.q ) THEN
690  temp1 = sva( p )
691  sva( p ) = sva( q )
692  sva( q ) = temp1
693  aapq = d( p )
694  d( p ) = d( q )
695  d( q ) = aapq
696  CALL zswap( m, a( 1, p ), 1, a( 1, q ), 1 )
697  IF( rsvec )CALL zswap( mvl, v( 1, p ), 1, v( 1, q ), 1 )
698  END IF
699  5991 CONTINUE
700 *
701 *
702  RETURN
703 * ..
704 * .. END OF ZGSVJ1
705 * ..
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:53
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
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:145
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:105
complex *16 function zdotc(N, ZX, INCX, ZY, INCY)
ZDOTC
Definition: zdotc.f:54
double precision function dznrm2(N, X, INCX)
DZNRM2
Definition: dznrm2.f:56
subroutine zlassq(N, X, INCX, SCALE, SUMSQ)
ZLASSQ updates a sum of squares represented in scaled form.
Definition: zlassq.f:108
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zaxpy(N, ZA, ZX, INCX, ZY, INCY)
ZAXPY
Definition: zaxpy.f:53

Here is the call graph for this function:

Here is the caller graph for this function: