LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dlaqr3 ( logical  WANTT,
logical  WANTZ,
integer  N,
integer  KTOP,
integer  KBOT,
integer  NW,
double precision, dimension( ldh, * )  H,
integer  LDH,
integer  ILOZ,
integer  IHIZ,
double precision, dimension( ldz, * )  Z,
integer  LDZ,
integer  NS,
integer  ND,
double precision, dimension( * )  SR,
double precision, dimension( * )  SI,
double precision, dimension( ldv, * )  V,
integer  LDV,
integer  NH,
double precision, dimension( ldt, * )  T,
integer  LDT,
integer  NV,
double precision, dimension( ldwv, * )  WV,
integer  LDWV,
double precision, dimension( * )  WORK,
integer  LWORK 
)

DLAQR3 performs the orthogonal similarity transformation of a Hessenberg matrix to detect and deflate fully converged eigenvalues from a trailing principal submatrix (aggressive early deflation).

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

Purpose:
    Aggressive early deflation:

    DLAQR3 accepts as input an upper Hessenberg matrix
    H and performs an orthogonal similarity transformation
    designed to detect and deflate fully converged eigenvalues from
    a trailing principal submatrix.  On output H has been over-
    written by a new Hessenberg matrix that is a perturbation of
    an orthogonal similarity transformation of H.  It is to be
    hoped that the final version of H has many zero subdiagonal
    entries.
Parameters
[in]WANTT
          WANTT is LOGICAL
          If .TRUE., then the Hessenberg matrix H is fully updated
          so that the quasi-triangular Schur factor may be
          computed (in cooperation with the calling subroutine).
          If .FALSE., then only enough of H is updated to preserve
          the eigenvalues.
[in]WANTZ
          WANTZ is LOGICAL
          If .TRUE., then the orthogonal matrix Z is updated so
          so that the orthogonal Schur factor may be computed
          (in cooperation with the calling subroutine).
          If .FALSE., then Z is not referenced.
[in]N
          N is INTEGER
          The order of the matrix H and (if WANTZ is .TRUE.) the
          order of the orthogonal matrix Z.
[in]KTOP
          KTOP is INTEGER
          It is assumed that either KTOP = 1 or H(KTOP,KTOP-1)=0.
          KBOT and KTOP together determine an isolated block
          along the diagonal of the Hessenberg matrix.
[in]KBOT
          KBOT is INTEGER
          It is assumed without a check that either
          KBOT = N or H(KBOT+1,KBOT)=0.  KBOT and KTOP together
          determine an isolated block along the diagonal of the
          Hessenberg matrix.
[in]NW
          NW is INTEGER
          Deflation window size.  1 .LE. NW .LE. (KBOT-KTOP+1).
[in,out]H
          H is DOUBLE PRECISION array, dimension (LDH,N)
          On input the initial N-by-N section of H stores the
          Hessenberg matrix undergoing aggressive early deflation.
          On output H has been transformed by an orthogonal
          similarity transformation, perturbed, and the returned
          to Hessenberg form that (it is to be hoped) has some
          zero subdiagonal entries.
[in]LDH
          LDH is integer
          Leading dimension of H just as declared in the calling
          subroutine.  N .LE. LDH
[in]ILOZ
          ILOZ is INTEGER
[in]IHIZ
          IHIZ is INTEGER
          Specify the rows of Z to which transformations must be
          applied if WANTZ is .TRUE.. 1 .LE. ILOZ .LE. IHIZ .LE. N.
[in,out]Z
          Z is DOUBLE PRECISION array, dimension (LDZ,N)
          IF WANTZ is .TRUE., then on output, the orthogonal
          similarity transformation mentioned above has been
          accumulated into Z(ILOZ:IHIZ,ILOZ:IHIZ) from the right.
          If WANTZ is .FALSE., then Z is unreferenced.
[in]LDZ
          LDZ is integer
          The leading dimension of Z just as declared in the
          calling subroutine.  1 .LE. LDZ.
[out]NS
          NS is integer
          The number of unconverged (ie approximate) eigenvalues
          returned in SR and SI that may be used as shifts by the
          calling subroutine.
[out]ND
          ND is integer
          The number of converged eigenvalues uncovered by this
          subroutine.
[out]SR
          SR is DOUBLE PRECISION array, dimension (KBOT)
[out]SI
          SI is DOUBLE PRECISION array, dimension (KBOT)
          On output, the real and imaginary parts of approximate
          eigenvalues that may be used for shifts are stored in
          SR(KBOT-ND-NS+1) through SR(KBOT-ND) and
          SI(KBOT-ND-NS+1) through SI(KBOT-ND), respectively.
          The real and imaginary parts of converged eigenvalues
          are stored in SR(KBOT-ND+1) through SR(KBOT) and
          SI(KBOT-ND+1) through SI(KBOT), respectively.
[out]V
          V is DOUBLE PRECISION array, dimension (LDV,NW)
          An NW-by-NW work array.
[in]LDV
          LDV is integer scalar
          The leading dimension of V just as declared in the
          calling subroutine.  NW .LE. LDV
[in]NH
          NH is integer scalar
          The number of columns of T.  NH.GE.NW.
[out]T
          T is DOUBLE PRECISION array, dimension (LDT,NW)
[in]LDT
          LDT is integer
          The leading dimension of T just as declared in the
          calling subroutine.  NW .LE. LDT
[in]NV
          NV is integer
          The number of rows of work array WV available for
          workspace.  NV.GE.NW.
[out]WV
          WV is DOUBLE PRECISION array, dimension (LDWV,NW)
[in]LDWV
          LDWV is integer
          The leading dimension of W just as declared in the
          calling subroutine.  NW .LE. LDV
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (LWORK)
          On exit, WORK(1) is set to an estimate of the optimal value
          of LWORK for the given values of N, NW, KTOP and KBOT.
[in]LWORK
          LWORK is integer
          The dimension of the work array WORK.  LWORK = 2*NW
          suffices, but greater efficiency may result from larger
          values of LWORK.

          If LWORK = -1, then a workspace query is assumed; DLAQR3
          only estimates the optimal workspace size for the given
          values of N, NW, KTOP and KBOT.  The estimate is returned
          in WORK(1).  No error message related to LWORK is issued
          by XERBLA.  Neither H nor Z are accessed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
June 2016
Contributors:
Karen Braman and Ralph Byers, Department of Mathematics, University of Kansas, USA

Definition at line 277 of file dlaqr3.f.

277 *
278 * -- LAPACK auxiliary routine (version 3.6.1) --
279 * -- LAPACK is a software package provided by Univ. of Tennessee, --
280 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
281 * June 2016
282 *
283 * .. Scalar Arguments ..
284  INTEGER ihiz, iloz, kbot, ktop, ldh, ldt, ldv, ldwv,
285  $ ldz, lwork, n, nd, nh, ns, nv, nw
286  LOGICAL wantt, wantz
287 * ..
288 * .. Array Arguments ..
289  DOUBLE PRECISION h( ldh, * ), si( * ), sr( * ), t( ldt, * ),
290  $ v( ldv, * ), work( * ), wv( ldwv, * ),
291  $ z( ldz, * )
292 * ..
293 *
294 * ================================================================
295 * .. Parameters ..
296  DOUBLE PRECISION zero, one
297  parameter ( zero = 0.0d0, one = 1.0d0 )
298 * ..
299 * .. Local Scalars ..
300  DOUBLE PRECISION aa, bb, beta, cc, cs, dd, evi, evk, foo, s,
301  $ safmax, safmin, smlnum, sn, tau, ulp
302  INTEGER i, ifst, ilst, info, infqr, j, jw, k, kcol,
303  $ kend, kln, krow, kwtop, ltop, lwk1, lwk2, lwk3,
304  $ lwkopt, nmin
305  LOGICAL bulge, sorted
306 * ..
307 * .. External Functions ..
308  DOUBLE PRECISION dlamch
309  INTEGER ilaenv
310  EXTERNAL dlamch, ilaenv
311 * ..
312 * .. External Subroutines ..
313  EXTERNAL dcopy, dgehrd, dgemm, dlabad, dlacpy, dlahqr,
315  $ dtrexc
316 * ..
317 * .. Intrinsic Functions ..
318  INTRINSIC abs, dble, int, max, min, sqrt
319 * ..
320 * .. Executable Statements ..
321 *
322 * ==== Estimate optimal workspace. ====
323 *
324  jw = min( nw, kbot-ktop+1 )
325  IF( jw.LE.2 ) THEN
326  lwkopt = 1
327  ELSE
328 *
329 * ==== Workspace query call to DGEHRD ====
330 *
331  CALL dgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
332  lwk1 = int( work( 1 ) )
333 *
334 * ==== Workspace query call to DORMHR ====
335 *
336  CALL dormhr( 'R', 'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
337  $ work, -1, info )
338  lwk2 = int( work( 1 ) )
339 *
340 * ==== Workspace query call to DLAQR4 ====
341 *
342  CALL dlaqr4( .true., .true., jw, 1, jw, t, ldt, sr, si, 1, jw,
343  $ v, ldv, work, -1, infqr )
344  lwk3 = int( work( 1 ) )
345 *
346 * ==== Optimal workspace ====
347 *
348  lwkopt = max( jw+max( lwk1, lwk2 ), lwk3 )
349  END IF
350 *
351 * ==== Quick return in case of workspace query. ====
352 *
353  IF( lwork.EQ.-1 ) THEN
354  work( 1 ) = dble( lwkopt )
355  RETURN
356  END IF
357 *
358 * ==== Nothing to do ...
359 * ... for an empty active block ... ====
360  ns = 0
361  nd = 0
362  work( 1 ) = one
363  IF( ktop.GT.kbot )
364  $ RETURN
365 * ... nor for an empty deflation window. ====
366  IF( nw.LT.1 )
367  $ RETURN
368 *
369 * ==== Machine constants ====
370 *
371  safmin = dlamch( 'SAFE MINIMUM' )
372  safmax = one / safmin
373  CALL dlabad( safmin, safmax )
374  ulp = dlamch( 'PRECISION' )
375  smlnum = safmin*( dble( n ) / ulp )
376 *
377 * ==== Setup deflation window ====
378 *
379  jw = min( nw, kbot-ktop+1 )
380  kwtop = kbot - jw + 1
381  IF( kwtop.EQ.ktop ) THEN
382  s = zero
383  ELSE
384  s = h( kwtop, kwtop-1 )
385  END IF
386 *
387  IF( kbot.EQ.kwtop ) THEN
388 *
389 * ==== 1-by-1 deflation window: not much to do ====
390 *
391  sr( kwtop ) = h( kwtop, kwtop )
392  si( kwtop ) = zero
393  ns = 1
394  nd = 0
395  IF( abs( s ).LE.max( smlnum, ulp*abs( h( kwtop, kwtop ) ) ) )
396  $ THEN
397  ns = 0
398  nd = 1
399  IF( kwtop.GT.ktop )
400  $ h( kwtop, kwtop-1 ) = zero
401  END IF
402  work( 1 ) = one
403  RETURN
404  END IF
405 *
406 * ==== Convert to spike-triangular form. (In case of a
407 * . rare QR failure, this routine continues to do
408 * . aggressive early deflation using that part of
409 * . the deflation window that converged using INFQR
410 * . here and there to keep track.) ====
411 *
412  CALL dlacpy( 'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
413  CALL dcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
414 *
415  CALL dlaset( 'A', jw, jw, zero, one, v, ldv )
416  nmin = ilaenv( 12, 'DLAQR3', 'SV', jw, 1, jw, lwork )
417  IF( jw.GT.nmin ) THEN
418  CALL dlaqr4( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
419  $ si( kwtop ), 1, jw, v, ldv, work, lwork, infqr )
420  ELSE
421  CALL dlahqr( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
422  $ si( kwtop ), 1, jw, v, ldv, infqr )
423  END IF
424 *
425 * ==== DTREXC needs a clean margin near the diagonal ====
426 *
427  DO 10 j = 1, jw - 3
428  t( j+2, j ) = zero
429  t( j+3, j ) = zero
430  10 CONTINUE
431  IF( jw.GT.2 )
432  $ t( jw, jw-2 ) = zero
433 *
434 * ==== Deflation detection loop ====
435 *
436  ns = jw
437  ilst = infqr + 1
438  20 CONTINUE
439  IF( ilst.LE.ns ) THEN
440  IF( ns.EQ.1 ) THEN
441  bulge = .false.
442  ELSE
443  bulge = t( ns, ns-1 ).NE.zero
444  END IF
445 *
446 * ==== Small spike tip test for deflation ====
447 *
448  IF( .NOT. bulge ) THEN
449 *
450 * ==== Real eigenvalue ====
451 *
452  foo = abs( t( ns, ns ) )
453  IF( foo.EQ.zero )
454  $ foo = abs( s )
455  IF( abs( s*v( 1, ns ) ).LE.max( smlnum, ulp*foo ) ) THEN
456 *
457 * ==== Deflatable ====
458 *
459  ns = ns - 1
460  ELSE
461 *
462 * ==== Undeflatable. Move it up out of the way.
463 * . (DTREXC can not fail in this case.) ====
464 *
465  ifst = ns
466  CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, work,
467  $ info )
468  ilst = ilst + 1
469  END IF
470  ELSE
471 *
472 * ==== Complex conjugate pair ====
473 *
474  foo = abs( t( ns, ns ) ) + sqrt( abs( t( ns, ns-1 ) ) )*
475  $ sqrt( abs( t( ns-1, ns ) ) )
476  IF( foo.EQ.zero )
477  $ foo = abs( s )
478  IF( max( abs( s*v( 1, ns ) ), abs( s*v( 1, ns-1 ) ) ).LE.
479  $ max( smlnum, ulp*foo ) ) THEN
480 *
481 * ==== Deflatable ====
482 *
483  ns = ns - 2
484  ELSE
485 *
486 * ==== Undeflatable. Move them up out of the way.
487 * . Fortunately, DTREXC does the right thing with
488 * . ILST in case of a rare exchange failure. ====
489 *
490  ifst = ns
491  CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, work,
492  $ info )
493  ilst = ilst + 2
494  END IF
495  END IF
496 *
497 * ==== End deflation detection loop ====
498 *
499  GO TO 20
500  END IF
501 *
502 * ==== Return to Hessenberg form ====
503 *
504  IF( ns.EQ.0 )
505  $ s = zero
506 *
507  IF( ns.LT.jw ) THEN
508 *
509 * ==== sorting diagonal blocks of T improves accuracy for
510 * . graded matrices. Bubble sort deals well with
511 * . exchange failures. ====
512 *
513  sorted = .false.
514  i = ns + 1
515  30 CONTINUE
516  IF( sorted )
517  $ GO TO 50
518  sorted = .true.
519 *
520  kend = i - 1
521  i = infqr + 1
522  IF( i.EQ.ns ) THEN
523  k = i + 1
524  ELSE IF( t( i+1, i ).EQ.zero ) THEN
525  k = i + 1
526  ELSE
527  k = i + 2
528  END IF
529  40 CONTINUE
530  IF( k.LE.kend ) THEN
531  IF( k.EQ.i+1 ) THEN
532  evi = abs( t( i, i ) )
533  ELSE
534  evi = abs( t( i, i ) ) + sqrt( abs( t( i+1, i ) ) )*
535  $ sqrt( abs( t( i, i+1 ) ) )
536  END IF
537 *
538  IF( k.EQ.kend ) THEN
539  evk = abs( t( k, k ) )
540  ELSE IF( t( k+1, k ).EQ.zero ) THEN
541  evk = abs( t( k, k ) )
542  ELSE
543  evk = abs( t( k, k ) ) + sqrt( abs( t( k+1, k ) ) )*
544  $ sqrt( abs( t( k, k+1 ) ) )
545  END IF
546 *
547  IF( evi.GE.evk ) THEN
548  i = k
549  ELSE
550  sorted = .false.
551  ifst = i
552  ilst = k
553  CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, work,
554  $ info )
555  IF( info.EQ.0 ) THEN
556  i = ilst
557  ELSE
558  i = k
559  END IF
560  END IF
561  IF( i.EQ.kend ) THEN
562  k = i + 1
563  ELSE IF( t( i+1, i ).EQ.zero ) THEN
564  k = i + 1
565  ELSE
566  k = i + 2
567  END IF
568  GO TO 40
569  END IF
570  GO TO 30
571  50 CONTINUE
572  END IF
573 *
574 * ==== Restore shift/eigenvalue array from T ====
575 *
576  i = jw
577  60 CONTINUE
578  IF( i.GE.infqr+1 ) THEN
579  IF( i.EQ.infqr+1 ) THEN
580  sr( kwtop+i-1 ) = t( i, i )
581  si( kwtop+i-1 ) = zero
582  i = i - 1
583  ELSE IF( t( i, i-1 ).EQ.zero ) THEN
584  sr( kwtop+i-1 ) = t( i, i )
585  si( kwtop+i-1 ) = zero
586  i = i - 1
587  ELSE
588  aa = t( i-1, i-1 )
589  cc = t( i, i-1 )
590  bb = t( i-1, i )
591  dd = t( i, i )
592  CALL dlanv2( aa, bb, cc, dd, sr( kwtop+i-2 ),
593  $ si( kwtop+i-2 ), sr( kwtop+i-1 ),
594  $ si( kwtop+i-1 ), cs, sn )
595  i = i - 2
596  END IF
597  GO TO 60
598  END IF
599 *
600  IF( ns.LT.jw .OR. s.EQ.zero ) THEN
601  IF( ns.GT.1 .AND. s.NE.zero ) THEN
602 *
603 * ==== Reflect spike back into lower triangle ====
604 *
605  CALL dcopy( ns, v, ldv, work, 1 )
606  beta = work( 1 )
607  CALL dlarfg( ns, beta, work( 2 ), 1, tau )
608  work( 1 ) = one
609 *
610  CALL dlaset( 'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
611 *
612  CALL dlarf( 'L', ns, jw, work, 1, tau, t, ldt,
613  $ work( jw+1 ) )
614  CALL dlarf( 'R', ns, ns, work, 1, tau, t, ldt,
615  $ work( jw+1 ) )
616  CALL dlarf( 'R', jw, ns, work, 1, tau, v, ldv,
617  $ work( jw+1 ) )
618 *
619  CALL dgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
620  $ lwork-jw, info )
621  END IF
622 *
623 * ==== Copy updated reduced window into place ====
624 *
625  IF( kwtop.GT.1 )
626  $ h( kwtop, kwtop-1 ) = s*v( 1, 1 )
627  CALL dlacpy( 'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
628  CALL dcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
629  $ ldh+1 )
630 *
631 * ==== Accumulate orthogonal matrix in order update
632 * . H and Z, if requested. ====
633 *
634  IF( ns.GT.1 .AND. s.NE.zero )
635  $ CALL dormhr( 'R', 'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
636  $ work( jw+1 ), lwork-jw, info )
637 *
638 * ==== Update vertical slab in H ====
639 *
640  IF( wantt ) THEN
641  ltop = 1
642  ELSE
643  ltop = ktop
644  END IF
645  DO 70 krow = ltop, kwtop - 1, nv
646  kln = min( nv, kwtop-krow )
647  CALL dgemm( 'N', 'N', kln, jw, jw, one, h( krow, kwtop ),
648  $ ldh, v, ldv, zero, wv, ldwv )
649  CALL dlacpy( 'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
650  70 CONTINUE
651 *
652 * ==== Update horizontal slab in H ====
653 *
654  IF( wantt ) THEN
655  DO 80 kcol = kbot + 1, n, nh
656  kln = min( nh, n-kcol+1 )
657  CALL dgemm( 'C', 'N', jw, kln, jw, one, v, ldv,
658  $ h( kwtop, kcol ), ldh, zero, t, ldt )
659  CALL dlacpy( 'A', jw, kln, t, ldt, h( kwtop, kcol ),
660  $ ldh )
661  80 CONTINUE
662  END IF
663 *
664 * ==== Update vertical slab in Z ====
665 *
666  IF( wantz ) THEN
667  DO 90 krow = iloz, ihiz, nv
668  kln = min( nv, ihiz-krow+1 )
669  CALL dgemm( 'N', 'N', kln, jw, jw, one, z( krow, kwtop ),
670  $ ldz, v, ldv, zero, wv, ldwv )
671  CALL dlacpy( 'A', kln, jw, wv, ldwv, z( krow, kwtop ),
672  $ ldz )
673  90 CONTINUE
674  END IF
675  END IF
676 *
677 * ==== Return the number of deflations ... ====
678 *
679  nd = jw - ns
680 *
681 * ==== ... and the number of shifts. (Subtracting
682 * . INFQR from the spike length takes care
683 * . of the case of a rare QR failure while
684 * . calculating eigenvalues of the deflation
685 * . window.) ====
686 *
687  ns = ns - infqr
688 *
689 * ==== Return optimal workspace. ====
690 *
691  work( 1 ) = dble( lwkopt )
692 *
693 * ==== End of DLAQR3 ====
694 *
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: dlaset.f:112
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dtrexc(COMPQ, N, T, LDT, Q, LDQ, IFST, ILST, WORK, INFO)
DTREXC
Definition: dtrexc.f:148
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
DGEHRD
Definition: dgehrd.f:169
subroutine dlarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
DLARF applies an elementary reflector to a general rectangular matrix.
Definition: dlarf.f:126
subroutine dlahqr(WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI, ILOZ, IHIZ, Z, LDZ, INFO)
DLAHQR computes the eigenvalues and Schur factorization of an upper Hessenberg matrix, using the double-shift/single-shift QR algorithm.
Definition: dlahqr.f:209
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:76
subroutine dlarfg(N, ALPHA, X, INCX, TAU)
DLARFG generates an elementary reflector (Householder matrix).
Definition: dlarfg.f:108
subroutine dlaqr4(WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI, ILOZ, IHIZ, Z, LDZ, WORK, LWORK, INFO)
DLAQR4 computes the eigenvalues of a Hessenberg matrix, and optionally the matrices from the Schur de...
Definition: dlaqr4.f:265
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine dormhr(SIDE, TRANS, M, N, ILO, IHI, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMHR
Definition: dormhr.f:180
subroutine dlanv2(A, B, C, D, RT1R, RT1I, RT2R, RT2I, CS, SN)
DLANV2 computes the Schur factorization of a real 2-by-2 nonsymmetric matrix in standard form...
Definition: dlanv2.f:129

Here is the call graph for this function:

Here is the caller graph for this function: