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

◆ dlaqr3()

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 <= NW <= (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 <= 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 <= ILOZ <= IHIZ <= 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 <= 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
!>          The leading dimension of V just as declared in the
!>          calling subroutine.  NW <= LDV
!> 
[in]NH
!>          NH is INTEGER
!>          The number of columns of T.  NH >= 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 <= LDT
!> 
[in]NV
!>          NV is INTEGER
!>          The number of rows of work array WV available for
!>          workspace.  NV >= 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 <= 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.
Contributors:
Karen Braman and Ralph Byers, Department of Mathematics, University of Kansas, USA

Definition at line 270 of file dlaqr3.f.

274*
275* -- LAPACK auxiliary routine --
276* -- LAPACK is a software package provided by Univ. of Tennessee, --
277* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
278*
279* .. Scalar Arguments ..
280 INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
281 $ LDZ, LWORK, N, ND, NH, NS, NV, NW
282 LOGICAL WANTT, WANTZ
283* ..
284* .. Array Arguments ..
285 DOUBLE PRECISION H( LDH, * ), SI( * ), SR( * ), T( LDT, * ),
286 $ V( LDV, * ), WORK( * ), WV( LDWV, * ),
287 $ Z( LDZ, * )
288* ..
289*
290* ================================================================
291* .. Parameters ..
292 DOUBLE PRECISION ZERO, ONE
293 parameter( zero = 0.0d0, one = 1.0d0 )
294* ..
295* .. Local Scalars ..
296 DOUBLE PRECISION AA, BB, BETA, CC, CS, DD, EVI, EVK, FOO, S,
297 $ SAFMAX, SAFMIN, SMLNUM, SN, TAU, ULP
298 INTEGER I, IFST, ILST, INFO, INFQR, J, JW, K, KCOL,
299 $ KEND, KLN, KROW, KWTOP, LTOP, LWK1, LWK2, LWK3,
300 $ LWKOPT, NMIN
301 LOGICAL BULGE, SORTED
302* ..
303* .. External Functions ..
304 DOUBLE PRECISION DLAMCH
305 INTEGER ILAENV
306 EXTERNAL dlamch, ilaenv
307* ..
308* .. External Subroutines ..
309 EXTERNAL dcopy, dgehrd, dgemm, dlacpy, dlahqr,
310 $ dlanv2,
312* ..
313* .. Intrinsic Functions ..
314 INTRINSIC abs, dble, int, max, min, sqrt
315* ..
316* .. Executable Statements ..
317*
318* ==== Estimate optimal workspace. ====
319*
320 jw = min( nw, kbot-ktop+1 )
321 IF( jw.LE.2 ) THEN
322 lwkopt = 1
323 ELSE
324*
325* ==== Workspace query call to DGEHRD ====
326*
327 CALL dgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
328 lwk1 = int( work( 1 ) )
329*
330* ==== Workspace query call to DORMHR ====
331*
332 CALL dormhr( 'R', 'N', jw, jw, 1, jw-1, t, ldt, work, v,
333 $ ldv,
334 $ work, -1, info )
335 lwk2 = int( work( 1 ) )
336*
337* ==== Workspace query call to DLAQR4 ====
338*
339 CALL dlaqr4( .true., .true., jw, 1, jw, t, ldt, sr, si, 1,
340 $ jw,
341 $ v, ldv, work, -1, infqr )
342 lwk3 = int( work( 1 ) )
343*
344* ==== Optimal workspace ====
345*
346 lwkopt = max( jw+max( lwk1, lwk2 ), lwk3 )
347 END IF
348*
349* ==== Quick return in case of workspace query. ====
350*
351 IF( lwork.EQ.-1 ) THEN
352 work( 1 ) = dble( lwkopt )
353 RETURN
354 END IF
355*
356* ==== Nothing to do ...
357* ... for an empty active block ... ====
358 ns = 0
359 nd = 0
360 work( 1 ) = one
361 IF( ktop.GT.kbot )
362 $ RETURN
363* ... nor for an empty deflation window. ====
364 IF( nw.LT.1 )
365 $ RETURN
366*
367* ==== Machine constants ====
368*
369 safmin = dlamch( 'SAFE MINIMUM' )
370 safmax = one / safmin
371 ulp = dlamch( 'PRECISION' )
372 smlnum = safmin*( dble( n ) / ulp )
373*
374* ==== Setup deflation window ====
375*
376 jw = min( nw, kbot-ktop+1 )
377 kwtop = kbot - jw + 1
378 IF( kwtop.EQ.ktop ) THEN
379 s = zero
380 ELSE
381 s = h( kwtop, kwtop-1 )
382 END IF
383*
384 IF( kbot.EQ.kwtop ) THEN
385*
386* ==== 1-by-1 deflation window: not much to do ====
387*
388 sr( kwtop ) = h( kwtop, kwtop )
389 si( kwtop ) = zero
390 ns = 1
391 nd = 0
392 IF( abs( s ).LE.max( smlnum, ulp*abs( h( kwtop, kwtop ) ) ) )
393 $ THEN
394 ns = 0
395 nd = 1
396 IF( kwtop.GT.ktop )
397 $ h( kwtop, kwtop-1 ) = zero
398 END IF
399 work( 1 ) = one
400 RETURN
401 END IF
402*
403* ==== Convert to spike-triangular form. (In case of a
404* . rare QR failure, this routine continues to do
405* . aggressive early deflation using that part of
406* . the deflation window that converged using INFQR
407* . here and there to keep track.) ====
408*
409 CALL dlacpy( 'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
410 CALL dcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ),
411 $ ldt+1 )
412*
413 CALL dlaset( 'A', jw, jw, zero, one, v, ldv )
414 nmin = ilaenv( 12, 'DLAQR3', 'SV', jw, 1, jw, lwork )
415 IF( jw.GT.nmin ) THEN
416 CALL dlaqr4( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
417 $ si( kwtop ), 1, jw, v, ldv, work, lwork, infqr )
418 ELSE
419 CALL dlahqr( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
420 $ si( kwtop ), 1, jw, v, ldv, infqr )
421 END IF
422*
423* ==== DTREXC needs a clean margin near the diagonal ====
424*
425 DO 10 j = 1, jw - 3
426 t( j+2, j ) = zero
427 t( j+3, j ) = zero
428 10 CONTINUE
429 IF( jw.GT.2 )
430 $ t( jw, jw-2 ) = zero
431*
432* ==== Deflation detection loop ====
433*
434 ns = jw
435 ilst = infqr + 1
436 20 CONTINUE
437 IF( ilst.LE.ns ) THEN
438 IF( ns.EQ.1 ) THEN
439 bulge = .false.
440 ELSE
441 bulge = t( ns, ns-1 ).NE.zero
442 END IF
443*
444* ==== Small spike tip test for deflation ====
445*
446 IF( .NOT. bulge ) THEN
447*
448* ==== Real eigenvalue ====
449*
450 foo = abs( t( ns, ns ) )
451 IF( foo.EQ.zero )
452 $ foo = abs( s )
453 IF( abs( s*v( 1, ns ) ).LE.max( smlnum, ulp*foo ) ) THEN
454*
455* ==== Deflatable ====
456*
457 ns = ns - 1
458 ELSE
459*
460* ==== Undeflatable. Move it up out of the way.
461* . (DTREXC can not fail in this case.) ====
462*
463 ifst = ns
464 CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst,
465 $ work,
466 $ info )
467 ilst = ilst + 1
468 END IF
469 ELSE
470*
471* ==== Complex conjugate pair ====
472*
473 foo = abs( t( ns, ns ) ) + sqrt( abs( t( ns, ns-1 ) ) )*
474 $ sqrt( abs( t( ns-1, ns ) ) )
475 IF( foo.EQ.zero )
476 $ foo = abs( s )
477 IF( max( abs( s*v( 1, ns ) ), abs( s*v( 1, ns-1 ) ) ).LE.
478 $ max( smlnum, ulp*foo ) ) THEN
479*
480* ==== Deflatable ====
481*
482 ns = ns - 2
483 ELSE
484*
485* ==== Undeflatable. Move them up out of the way.
486* . Fortunately, DTREXC does the right thing with
487* . ILST in case of a rare exchange failure. ====
488*
489 ifst = ns
490 CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst,
491 $ 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,
554 $ work,
555 $ info )
556 IF( info.EQ.0 ) THEN
557 i = ilst
558 ELSE
559 i = k
560 END IF
561 END IF
562 IF( i.EQ.kend ) THEN
563 k = i + 1
564 ELSE IF( t( i+1, i ).EQ.zero ) THEN
565 k = i + 1
566 ELSE
567 k = i + 2
568 END IF
569 GO TO 40
570 END IF
571 GO TO 30
572 50 CONTINUE
573 END IF
574*
575* ==== Restore shift/eigenvalue array from T ====
576*
577 i = jw
578 60 CONTINUE
579 IF( i.GE.infqr+1 ) THEN
580 IF( i.EQ.infqr+1 ) THEN
581 sr( kwtop+i-1 ) = t( i, i )
582 si( kwtop+i-1 ) = zero
583 i = i - 1
584 ELSE IF( t( i, i-1 ).EQ.zero ) THEN
585 sr( kwtop+i-1 ) = t( i, i )
586 si( kwtop+i-1 ) = zero
587 i = i - 1
588 ELSE
589 aa = t( i-1, i-1 )
590 cc = t( i, i-1 )
591 bb = t( i-1, i )
592 dd = t( i, i )
593 CALL dlanv2( aa, bb, cc, dd, sr( kwtop+i-2 ),
594 $ si( kwtop+i-2 ), sr( kwtop+i-1 ),
595 $ si( kwtop+i-1 ), cs, sn )
596 i = i - 2
597 END IF
598 GO TO 60
599 END IF
600*
601 IF( ns.LT.jw .OR. s.EQ.zero ) THEN
602 IF( ns.GT.1 .AND. s.NE.zero ) THEN
603*
604* ==== Reflect spike back into lower triangle ====
605*
606 CALL dcopy( ns, v, ldv, work, 1 )
607 beta = work( 1 )
608 CALL dlarfg( ns, beta, work( 2 ), 1, tau )
609*
610 CALL dlaset( 'L', jw-2, jw-2, zero, zero, t( 3, 1 ),
611 $ ldt )
612*
613 CALL dlarf1f( 'L', ns, jw, work, 1, tau, t, ldt,
614 $ work( jw+1 ) )
615 CALL dlarf1f( 'R', ns, ns, work, 1, tau, t, ldt,
616 $ work( jw+1 ) )
617 CALL dlarf1f( 'R', jw, ns, work, 1, tau, v, ldv,
618 $ work( jw+1 ) )
619*
620 CALL dgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
621 $ lwork-jw, info )
622 END IF
623*
624* ==== Copy updated reduced window into place ====
625*
626 IF( kwtop.GT.1 )
627 $ h( kwtop, kwtop-1 ) = s*v( 1, 1 )
628 CALL dlacpy( 'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
629 CALL dcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
630 $ ldh+1 )
631*
632* ==== Accumulate orthogonal matrix in order update
633* . H and Z, if requested. ====
634*
635 IF( ns.GT.1 .AND. s.NE.zero )
636 $ CALL dormhr( 'R', 'N', jw, ns, 1, ns, t, ldt, work, v,
637 $ ldv,
638 $ work( jw+1 ), lwork-jw, info )
639*
640* ==== Update vertical slab in H ====
641*
642 IF( wantt ) THEN
643 ltop = 1
644 ELSE
645 ltop = ktop
646 END IF
647 DO 70 krow = ltop, kwtop - 1, nv
648 kln = min( nv, kwtop-krow )
649 CALL dgemm( 'N', 'N', kln, jw, jw, one, h( krow, kwtop ),
650 $ ldh, v, ldv, zero, wv, ldwv )
651 CALL dlacpy( 'A', kln, jw, wv, ldwv, h( krow, kwtop ),
652 $ ldh )
653 70 CONTINUE
654*
655* ==== Update horizontal slab in H ====
656*
657 IF( wantt ) THEN
658 DO 80 kcol = kbot + 1, n, nh
659 kln = min( nh, n-kcol+1 )
660 CALL dgemm( 'C', 'N', jw, kln, jw, one, v, ldv,
661 $ h( kwtop, kcol ), ldh, zero, t, ldt )
662 CALL dlacpy( 'A', jw, kln, t, ldt, h( kwtop, kcol ),
663 $ ldh )
664 80 CONTINUE
665 END IF
666*
667* ==== Update vertical slab in Z ====
668*
669 IF( wantz ) THEN
670 DO 90 krow = iloz, ihiz, nv
671 kln = min( nv, ihiz-krow+1 )
672 CALL dgemm( 'N', 'N', kln, jw, jw, one, z( krow,
673 $ kwtop ),
674 $ ldz, v, ldv, zero, wv, ldwv )
675 CALL dlacpy( 'A', kln, jw, wv, ldwv, z( krow, kwtop ),
676 $ ldz )
677 90 CONTINUE
678 END IF
679 END IF
680*
681* ==== Return the number of deflations ... ====
682*
683 nd = jw - ns
684*
685* ==== ... and the number of shifts. (Subtracting
686* . INFQR from the spike length takes care
687* . of the case of a rare QR failure while
688* . calculating eigenvalues of the deflation
689* . window.) ====
690*
691 ns = ns - infqr
692*
693* ==== Return optimal workspace. ====
694*
695 work( 1 ) = dble( lwkopt )
696*
697* ==== End of DLAQR3 ====
698*
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
DGEHRD
Definition dgehrd.f:166
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
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,...
Definition dlahqr.f:205
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
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:125
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:261
subroutine dlarf1f(side, m, n, v, incv, tau, c, ldc, work)
DLARF1F applies an elementary reflector to a general rectangular
Definition dlarf1f.f:157
subroutine dlarfg(n, alpha, x, incx, tau)
DLARFG generates an elementary reflector (Householder matrix).
Definition dlarfg.f:104
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:108
subroutine dtrexc(compq, n, t, ldt, q, ldq, ifst, ilst, work, info)
DTREXC
Definition dtrexc.f:146
subroutine dormhr(side, trans, m, n, ilo, ihi, a, lda, tau, c, ldc, work, lwork, info)
DORMHR
Definition dormhr.f:176
Here is the call graph for this function:
Here is the caller graph for this function: