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

◆ cgsvj1()

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 tuning 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 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 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 = '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 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 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)))) > 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 >= 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 232 of file cgsvj1.f.

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