LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages

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