LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine cgsvj1 ( character*1  JOBV,
integer  M,
integer  N,
integer  N1,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( n )  D,
real, dimension( n )  SVA,
integer  MV,
complex, dimension( ldv, * )  V,
integer  LDV,
real  EPS,
real  SFMIN,
real  TOL,
integer  NSWEEP,
complex, dimension( lwork )  WORK,
integer  LWORK,
integer  INFO 
)

CGSVJ1 pre-processor for the routine cgesvj, applies Jacobi rotations targeting only particular pivots.

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

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

Here is the call graph for this function:

Here is the caller graph for this function: