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

◆ dtrsen()

subroutine dtrsen ( character  job,
character  compq,
logical, dimension( * )  select,
integer  n,
double precision, dimension( ldt, * )  t,
integer  ldt,
double precision, dimension( ldq, * )  q,
integer  ldq,
double precision, dimension( * )  wr,
double precision, dimension( * )  wi,
integer  m,
double precision  s,
double precision  sep,
double precision, dimension( * )  work,
integer  lwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

DTRSEN

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

Purpose:
 DTRSEN reorders the real Schur factorization of a real matrix
 A = Q*T*Q**T, so that a selected cluster of eigenvalues appears in
 the leading diagonal blocks of the upper quasi-triangular matrix T,
 and the leading columns of Q form an orthonormal basis of the
 corresponding right invariant subspace.

 Optionally the routine computes the reciprocal condition numbers of
 the cluster of eigenvalues and/or the invariant subspace.

 T must be in Schur canonical form (as returned by DHSEQR), that is,
 block upper triangular with 1-by-1 and 2-by-2 diagonal blocks; each
 2-by-2 diagonal block has its diagonal elements equal and its
 off-diagonal elements of opposite sign.
Parameters
[in]JOB
          JOB is CHARACTER*1
          Specifies whether condition numbers are required for the
          cluster of eigenvalues (S) or the invariant subspace (SEP):
          = 'N': none;
          = 'E': for eigenvalues only (S);
          = 'V': for invariant subspace only (SEP);
          = 'B': for both eigenvalues and invariant subspace (S and
                 SEP).
[in]COMPQ
          COMPQ is CHARACTER*1
          = 'V': update the matrix Q of Schur vectors;
          = 'N': do not update Q.
[in]SELECT
          SELECT is LOGICAL array, dimension (N)
          SELECT specifies the eigenvalues in the selected cluster. To
          select a real eigenvalue w(j), SELECT(j) must be set to
          .TRUE.. To select a complex conjugate pair of eigenvalues
          w(j) and w(j+1), corresponding to a 2-by-2 diagonal block,
          either SELECT(j) or SELECT(j+1) or both must be set to
          .TRUE.; a complex conjugate pair of eigenvalues must be
          either both included in the cluster or both excluded.
[in]N
          N is INTEGER
          The order of the matrix T. N >= 0.
[in,out]T
          T is DOUBLE PRECISION array, dimension (LDT,N)
          On entry, the upper quasi-triangular matrix T, in Schur
          canonical form.
          On exit, T is overwritten by the reordered matrix T, again in
          Schur canonical form, with the selected eigenvalues in the
          leading diagonal blocks.
[in]LDT
          LDT is INTEGER
          The leading dimension of the array T. LDT >= max(1,N).
[in,out]Q
          Q is DOUBLE PRECISION array, dimension (LDQ,N)
          On entry, if COMPQ = 'V', the matrix Q of Schur vectors.
          On exit, if COMPQ = 'V', Q has been postmultiplied by the
          orthogonal transformation matrix which reorders T; the
          leading M columns of Q form an orthonormal basis for the
          specified invariant subspace.
          If COMPQ = 'N', Q is not referenced.
[in]LDQ
          LDQ is INTEGER
          The leading dimension of the array Q.
          LDQ >= 1; and if COMPQ = 'V', LDQ >= N.
[out]WR
          WR is DOUBLE PRECISION array, dimension (N)
[out]WI
          WI is DOUBLE PRECISION array, dimension (N)

          The real and imaginary parts, respectively, of the reordered
          eigenvalues of T. The eigenvalues are stored in the same
          order as on the diagonal of T, with WR(i) = T(i,i) and, if
          T(i:i+1,i:i+1) is a 2-by-2 diagonal block, WI(i) > 0 and
          WI(i+1) = -WI(i). Note that if a complex eigenvalue is
          sufficiently ill-conditioned, then its value may differ
          significantly from its value before reordering.
[out]M
          M is INTEGER
          The dimension of the specified invariant subspace.
          0 < = M <= N.
[out]S
          S is DOUBLE PRECISION
          If JOB = 'E' or 'B', S is a lower bound on the reciprocal
          condition number for the selected cluster of eigenvalues.
          S cannot underestimate the true reciprocal condition number
          by more than a factor of sqrt(N). If M = 0 or N, S = 1.
          If JOB = 'N' or 'V', S is not referenced.
[out]SEP
          SEP is DOUBLE PRECISION
          If JOB = 'V' or 'B', SEP is the estimated reciprocal
          condition number of the specified invariant subspace. If
          M = 0 or N, SEP = norm(T).
          If JOB = 'N' or 'E', SEP is not referenced.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.
          If JOB = 'N', LWORK >= max(1,N);
          if JOB = 'E', LWORK >= max(1,M*(N-M));
          if JOB = 'V' or 'B', LWORK >= max(1,2*M*(N-M)).

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the array IWORK.
          If JOB = 'N' or 'E', LIWORK >= 1;
          if JOB = 'V' or 'B', LIWORK >= max(1,M*(N-M)).

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal size of the IWORK array,
          returns this value as the first entry of the IWORK array, and
          no error message related to LIWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -i, the i-th argument had an illegal value
          = 1: reordering of T failed because some eigenvalues are too
               close to separate (the problem is very ill-conditioned);
               T may have been partially reordered, and WR and WI
               contain the eigenvalues in the same order as in T; S and
               SEP (if requested) are set to zero.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
  DTRSEN first collects the selected eigenvalues by computing an
  orthogonal transformation Z to move them to the top left corner of T.
  In other words, the selected eigenvalues are the eigenvalues of T11
  in:

          Z**T * T * Z = ( T11 T12 ) n1
                         (  0  T22 ) n2
                            n1  n2

  where N = n1+n2 and Z**T means the transpose of Z. The first n1 columns
  of Z span the specified invariant subspace of T.

  If T has been obtained from the real Schur factorization of a matrix
  A = Q*T*Q**T, then the reordered real Schur factorization of A is given
  by A = (Q*Z)*(Z**T*T*Z)*(Q*Z)**T, and the first n1 columns of Q*Z span
  the corresponding invariant subspace of A.

  The reciprocal condition number of the average of the eigenvalues of
  T11 may be returned in S. S lies between 0 (very badly conditioned)
  and 1 (very well conditioned). It is computed as follows. First we
  compute R so that

                         P = ( I  R ) n1
                             ( 0  0 ) n2
                               n1 n2

  is the projector on the invariant subspace associated with T11.
  R is the solution of the Sylvester equation:

                        T11*R - R*T22 = T12.

  Let F-norm(M) denote the Frobenius-norm of M and 2-norm(M) denote
  the two-norm of M. Then S is computed as the lower bound

                      (1 + F-norm(R)**2)**(-1/2)

  on the reciprocal of 2-norm(P), the true reciprocal condition number.
  S cannot underestimate 1 / 2-norm(P) by more than a factor of
  sqrt(N).

  An approximate error bound for the computed average of the
  eigenvalues of T11 is

                         EPS * norm(T) / S

  where EPS is the machine precision.

  The reciprocal condition number of the right invariant subspace
  spanned by the first n1 columns of Z (or of Q*Z) is returned in SEP.
  SEP is defined as the separation of T11 and T22:

                     sep( T11, T22 ) = sigma-min( C )

  where sigma-min(C) is the smallest singular value of the
  n1*n2-by-n1*n2 matrix

     C  = kprod( I(n2), T11 ) - kprod( transpose(T22), I(n1) )

  I(m) is an m by m identity matrix, and kprod denotes the Kronecker
  product. We estimate sigma-min(C) by the reciprocal of an estimate of
  the 1-norm of inverse(C). The true reciprocal 1-norm of inverse(C)
  cannot differ from sigma-min(C) by more than a factor of sqrt(n1*n2).

  When SEP is small, small changes in T can cause large changes in
  the invariant subspace. An approximate bound on the maximum angular
  error in the computed right invariant subspace is

                      EPS * norm(T) / SEP

Definition at line 311 of file dtrsen.f.

313*
314* -- LAPACK computational routine --
315* -- LAPACK is a software package provided by Univ. of Tennessee, --
316* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
317*
318* .. Scalar Arguments ..
319 CHARACTER COMPQ, JOB
320 INTEGER INFO, LDQ, LDT, LIWORK, LWORK, M, N
321 DOUBLE PRECISION S, SEP
322* ..
323* .. Array Arguments ..
324 LOGICAL SELECT( * )
325 INTEGER IWORK( * )
326 DOUBLE PRECISION Q( LDQ, * ), T( LDT, * ), WI( * ), WORK( * ),
327 $ WR( * )
328* ..
329*
330* =====================================================================
331*
332* .. Parameters ..
333 DOUBLE PRECISION ZERO, ONE
334 parameter( zero = 0.0d+0, one = 1.0d+0 )
335* ..
336* .. Local Scalars ..
337 LOGICAL LQUERY, PAIR, SWAP, WANTBH, WANTQ, WANTS,
338 $ WANTSP
339 INTEGER IERR, K, KASE, KK, KS, LIWMIN, LWMIN, N1, N2,
340 $ NN
341 DOUBLE PRECISION EST, RNORM, SCALE
342* ..
343* .. Local Arrays ..
344 INTEGER ISAVE( 3 )
345* ..
346* .. External Functions ..
347 LOGICAL LSAME
348 DOUBLE PRECISION DLANGE
349 EXTERNAL lsame, dlange
350* ..
351* .. External Subroutines ..
352 EXTERNAL dlacn2, dlacpy, dtrexc, dtrsyl, xerbla
353* ..
354* .. Intrinsic Functions ..
355 INTRINSIC abs, max, sqrt
356* ..
357* .. Executable Statements ..
358*
359* Decode and test the input parameters
360*
361 wantbh = lsame( job, 'B' )
362 wants = lsame( job, 'E' ) .OR. wantbh
363 wantsp = lsame( job, 'V' ) .OR. wantbh
364 wantq = lsame( compq, 'V' )
365*
366 info = 0
367 lquery = ( lwork.EQ.-1 )
368 IF( .NOT.lsame( job, 'N' ) .AND. .NOT.wants .AND. .NOT.wantsp )
369 $ THEN
370 info = -1
371 ELSE IF( .NOT.lsame( compq, 'N' ) .AND. .NOT.wantq ) THEN
372 info = -2
373 ELSE IF( n.LT.0 ) THEN
374 info = -4
375 ELSE IF( ldt.LT.max( 1, n ) ) THEN
376 info = -6
377 ELSE IF( ldq.LT.1 .OR. ( wantq .AND. ldq.LT.n ) ) THEN
378 info = -8
379 ELSE
380*
381* Set M to the dimension of the specified invariant subspace,
382* and test LWORK and LIWORK.
383*
384 m = 0
385 pair = .false.
386 DO 10 k = 1, n
387 IF( pair ) THEN
388 pair = .false.
389 ELSE
390 IF( k.LT.n ) THEN
391 IF( t( k+1, k ).EQ.zero ) THEN
392 IF( SELECT( k ) )
393 $ m = m + 1
394 ELSE
395 pair = .true.
396 IF( SELECT( k ) .OR. SELECT( k+1 ) )
397 $ m = m + 2
398 END IF
399 ELSE
400 IF( SELECT( n ) )
401 $ m = m + 1
402 END IF
403 END IF
404 10 CONTINUE
405*
406 n1 = m
407 n2 = n - m
408 nn = n1*n2
409*
410 IF( wantsp ) THEN
411 lwmin = max( 1, 2*nn )
412 liwmin = max( 1, nn )
413 ELSE IF( lsame( job, 'N' ) ) THEN
414 lwmin = max( 1, n )
415 liwmin = 1
416 ELSE IF( lsame( job, 'E' ) ) THEN
417 lwmin = max( 1, nn )
418 liwmin = 1
419 END IF
420*
421 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
422 info = -15
423 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
424 info = -17
425 END IF
426 END IF
427*
428 IF( info.EQ.0 ) THEN
429 work( 1 ) = lwmin
430 iwork( 1 ) = liwmin
431 END IF
432*
433 IF( info.NE.0 ) THEN
434 CALL xerbla( 'DTRSEN', -info )
435 RETURN
436 ELSE IF( lquery ) THEN
437 RETURN
438 END IF
439*
440* Quick return if possible.
441*
442 IF( m.EQ.n .OR. m.EQ.0 ) THEN
443 IF( wants )
444 $ s = one
445 IF( wantsp )
446 $ sep = dlange( '1', n, n, t, ldt, work )
447 GO TO 40
448 END IF
449*
450* Collect the selected blocks at the top-left corner of T.
451*
452 ks = 0
453 pair = .false.
454 DO 20 k = 1, n
455 IF( pair ) THEN
456 pair = .false.
457 ELSE
458 swap = SELECT( k )
459 IF( k.LT.n ) THEN
460 IF( t( k+1, k ).NE.zero ) THEN
461 pair = .true.
462 swap = swap .OR. SELECT( k+1 )
463 END IF
464 END IF
465 IF( swap ) THEN
466 ks = ks + 1
467*
468* Swap the K-th block to position KS.
469*
470 ierr = 0
471 kk = k
472 IF( k.NE.ks )
473 $ CALL dtrexc( compq, n, t, ldt, q, ldq, kk, ks, work,
474 $ ierr )
475 IF( ierr.EQ.1 .OR. ierr.EQ.2 ) THEN
476*
477* Blocks too close to swap: exit.
478*
479 info = 1
480 IF( wants )
481 $ s = zero
482 IF( wantsp )
483 $ sep = zero
484 GO TO 40
485 END IF
486 IF( pair )
487 $ ks = ks + 1
488 END IF
489 END IF
490 20 CONTINUE
491*
492 IF( wants ) THEN
493*
494* Solve Sylvester equation for R:
495*
496* T11*R - R*T22 = scale*T12
497*
498 CALL dlacpy( 'F', n1, n2, t( 1, n1+1 ), ldt, work, n1 )
499 CALL dtrsyl( 'N', 'N', -1, n1, n2, t, ldt, t( n1+1, n1+1 ),
500 $ ldt, work, n1, scale, ierr )
501*
502* Estimate the reciprocal of the condition number of the cluster
503* of eigenvalues.
504*
505 rnorm = dlange( 'F', n1, n2, work, n1, work )
506 IF( rnorm.EQ.zero ) THEN
507 s = one
508 ELSE
509 s = scale / ( sqrt( scale*scale / rnorm+rnorm )*
510 $ sqrt( rnorm ) )
511 END IF
512 END IF
513*
514 IF( wantsp ) THEN
515*
516* Estimate sep(T11,T22).
517*
518 est = zero
519 kase = 0
520 30 CONTINUE
521 CALL dlacn2( nn, work( nn+1 ), work, iwork, est, kase, isave )
522 IF( kase.NE.0 ) THEN
523 IF( kase.EQ.1 ) THEN
524*
525* Solve T11*R - R*T22 = scale*X.
526*
527 CALL dtrsyl( 'N', 'N', -1, n1, n2, t, ldt,
528 $ t( n1+1, n1+1 ), ldt, work, n1, scale,
529 $ ierr )
530 ELSE
531*
532* Solve T11**T*R - R*T22**T = scale*X.
533*
534 CALL dtrsyl( 'T', 'T', -1, n1, n2, t, ldt,
535 $ t( n1+1, n1+1 ), ldt, work, n1, scale,
536 $ ierr )
537 END IF
538 GO TO 30
539 END IF
540*
541 sep = scale / est
542 END IF
543*
544 40 CONTINUE
545*
546* Store the output eigenvalues in WR and WI.
547*
548 DO 50 k = 1, n
549 wr( k ) = t( k, k )
550 wi( k ) = zero
551 50 CONTINUE
552 DO 60 k = 1, n - 1
553 IF( t( k+1, k ).NE.zero ) THEN
554 wi( k ) = sqrt( abs( t( k, k+1 ) ) )*
555 $ sqrt( abs( t( k+1, k ) ) )
556 wi( k+1 ) = -wi( k )
557 END IF
558 60 CONTINUE
559*
560 work( 1 ) = lwmin
561 iwork( 1 ) = liwmin
562*
563 RETURN
564*
565* End of DTRSEN
566*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dlacn2(n, v, x, isgn, est, kase, isave)
DLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition dlacn2.f:136
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
double precision function dlange(norm, m, n, a, lda, work)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlange.f:114
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dtrexc(compq, n, t, ldt, q, ldq, ifst, ilst, work, info)
DTREXC
Definition dtrexc.f:148
subroutine dtrsyl(trana, tranb, isgn, m, n, a, lda, b, ldb, c, ldc, scale, info)
DTRSYL
Definition dtrsyl.f:164
Here is the call graph for this function:
Here is the caller graph for this function: