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

◆ ztgsyl()

subroutine ztgsyl ( character trans,
integer ijob,
integer m,
integer n,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( ldc, * ) c,
integer ldc,
complex*16, dimension( ldd, * ) d,
integer ldd,
complex*16, dimension( lde, * ) e,
integer lde,
complex*16, dimension( ldf, * ) f,
integer ldf,
double precision scale,
double precision dif,
complex*16, dimension( * ) work,
integer lwork,
integer, dimension( * ) iwork,
integer info )

ZTGSYL

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

Purpose:
!>
!> ZTGSYL solves the generalized Sylvester equation:
!>
!>             A * R - L * B = scale * C            (1)
!>             D * R - L * E = scale * F
!>
!> where R and L are unknown m-by-n matrices, (A, D), (B, E) and
!> (C, F) are given matrix pairs of size m-by-m, n-by-n and m-by-n,
!> respectively, with complex entries. A, B, D and E are upper
!> triangular (i.e., (A,D) and (B,E) in generalized Schur form).
!>
!> The solution (R, L) overwrites (C, F). 0 <= SCALE <= 1
!> is an output scaling factor chosen to avoid overflow.
!>
!> In matrix notation (1) is equivalent to solve Zx = scale*b, where Z
!> is defined as
!>
!>        Z = [ kron(In, A)  -kron(B**H, Im) ]        (2)
!>            [ kron(In, D)  -kron(E**H, Im) ],
!>
!> Here Ix is the identity matrix of size x and X**H is the conjugate
!> transpose of X. Kron(X, Y) is the Kronecker product between the
!> matrices X and Y.
!>
!> If TRANS = 'C', y in the conjugate transposed system Z**H *y = scale*b
!> is solved for, which is equivalent to solve for R and L in
!>
!>             A**H * R + D**H * L = scale * C           (3)
!>             R * B**H + L * E**H = scale * -F
!>
!> This case (TRANS = 'C') is used to compute an one-norm-based estimate
!> of Dif[(A,D), (B,E)], the separation between the matrix pairs (A,D)
!> and (B,E), using ZLACON.
!>
!> If IJOB >= 1, ZTGSYL computes a Frobenius norm-based estimate of
!> Dif[(A,D),(B,E)]. That is, the reciprocal of a lower bound on the
!> reciprocal of the smallest singular value of Z.
!>
!> This is a level-3 BLAS algorithm.
!> 
Parameters
[in]TRANS
!>          TRANS is CHARACTER*1
!>          = 'N': solve the generalized sylvester equation (1).
!>          = 'C': solve the  system (3).
!> 
[in]IJOB
!>          IJOB is INTEGER
!>          Specifies what kind of functionality to be performed.
!>          =0: solve (1) only.
!>          =1: The functionality of 0 and 3.
!>          =2: The functionality of 0 and 4.
!>          =3: Only an estimate of Dif[(A,D), (B,E)] is computed.
!>              (look ahead strategy is used).
!>          =4: Only an estimate of Dif[(A,D), (B,E)] is computed.
!>              (ZGECON on sub-systems is used).
!>          Not referenced if TRANS = 'C'.
!> 
[in]M
!>          M is INTEGER
!>          The order of the matrices A and D, and the row dimension of
!>          the matrices C, F, R and L.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrices B and E, and the column dimension
!>          of the matrices C, F, R and L.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension (LDA, M)
!>          The upper triangular matrix A.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A. LDA >= max(1, M).
!> 
[in]B
!>          B is COMPLEX*16 array, dimension (LDB, N)
!>          The upper triangular matrix B.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B. LDB >= max(1, N).
!> 
[in,out]C
!>          C is COMPLEX*16 array, dimension (LDC, N)
!>          On entry, C contains the right-hand-side of the first matrix
!>          equation in (1) or (3).
!>          On exit, if IJOB = 0, 1 or 2, C has been overwritten by
!>          the solution R. If IJOB = 3 or 4 and TRANS = 'N', C holds R,
!>          the solution achieved during the computation of the
!>          Dif-estimate.
!> 
[in]LDC
!>          LDC is INTEGER
!>          The leading dimension of the array C. LDC >= max(1, M).
!> 
[in]D
!>          D is COMPLEX*16 array, dimension (LDD, M)
!>          The upper triangular matrix D.
!> 
[in]LDD
!>          LDD is INTEGER
!>          The leading dimension of the array D. LDD >= max(1, M).
!> 
[in]E
!>          E is COMPLEX*16 array, dimension (LDE, N)
!>          The upper triangular matrix E.
!> 
[in]LDE
!>          LDE is INTEGER
!>          The leading dimension of the array E. LDE >= max(1, N).
!> 
[in,out]F
!>          F is COMPLEX*16 array, dimension (LDF, N)
!>          On entry, F contains the right-hand-side of the second matrix
!>          equation in (1) or (3).
!>          On exit, if IJOB = 0, 1 or 2, F has been overwritten by
!>          the solution L. If IJOB = 3 or 4 and TRANS = 'N', F holds L,
!>          the solution achieved during the computation of the
!>          Dif-estimate.
!> 
[in]LDF
!>          LDF is INTEGER
!>          The leading dimension of the array F. LDF >= max(1, M).
!> 
[out]DIF
!>          DIF is DOUBLE PRECISION
!>          On exit DIF is the reciprocal of a lower bound of the
!>          reciprocal of the Dif-function, i.e. DIF is an upper bound of
!>          Dif[(A,D), (B,E)] = sigma-min(Z), where Z as in (2).
!>          IF IJOB = 0 or TRANS = 'C', DIF is not referenced.
!> 
[out]SCALE
!>          SCALE is DOUBLE PRECISION
!>          On exit SCALE is the scaling factor in (1) or (3).
!>          If 0 < SCALE < 1, C and F hold the solutions R and L, resp.,
!>          to a slightly perturbed system but the input matrices A, B,
!>          D and E have not been changed. If SCALE = 0, R and L will
!>          hold the solutions to the homogeneous system with C = F = 0.
!> 
[out]WORK
!>          WORK is COMPLEX*16 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. LWORK > = 1.
!>          If IJOB = 1 or 2 and TRANS = 'N', LWORK >= max(1,2*M*N).
!>
!>          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 (M+N+2)
!> 
[out]INFO
!>          INFO is INTEGER
!>            =0: successful exit
!>            <0: If INFO = -i, the i-th argument had an illegal value.
!>            >0: (A, D) and (B, E) have common or very close
!>                eigenvalues.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Bo Kagstrom and Peter Poromaa, Department of Computing Science, Umea University, S-901 87 Umea, Sweden.
References:
[1] B. Kagstrom and P. Poromaa, LAPACK-Style Algorithms and Software for Solving the Generalized Sylvester Equation and Estimating the Separation between Regular Matrix Pairs, Report UMINF - 93.23, Department of Computing Science, Umea University, S-901 87 Umea, Sweden, December 1993, Revised April 1994, Also as LAPACK Working Note 75. To appear in ACM Trans. on Math. Software, Vol 22, No 1, 1996.
[2] B. Kagstrom, A Perturbation Analysis of the Generalized Sylvester Equation (AR - LB, DR - LE ) = (C, F), SIAM J. Matrix Anal. Appl., 15(4):1045-1060, 1994.
[3] B. Kagstrom and L. Westin, Generalized Schur Methods with Condition Estimators for Solving the Generalized Sylvester Equation, IEEE Transactions on Automatic Control, Vol. 34, No. 7, July 1989, pp 745-751.

Definition at line 290 of file ztgsyl.f.

294*
295* -- LAPACK computational routine --
296* -- LAPACK is a software package provided by Univ. of Tennessee, --
297* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
298*
299* .. Scalar Arguments ..
300 CHARACTER TRANS
301 INTEGER IJOB, INFO, LDA, LDB, LDC, LDD, LDE, LDF,
302 $ LWORK, M, N
303 DOUBLE PRECISION DIF, SCALE
304* ..
305* .. Array Arguments ..
306 INTEGER IWORK( * )
307 COMPLEX*16 A( LDA, * ), B( LDB, * ), C( LDC, * ),
308 $ D( LDD, * ), E( LDE, * ), F( LDF, * ),
309 $ WORK( * )
310* ..
311*
312* =====================================================================
313* Replaced various illegal calls to CCOPY by calls to CLASET.
314* Sven Hammarling, 1/5/02.
315*
316* .. Parameters ..
317 DOUBLE PRECISION ZERO, ONE
318 parameter( zero = 0.0d+0, one = 1.0d+0 )
319 COMPLEX*16 CZERO
320 parameter( czero = (0.0d+0, 0.0d+0) )
321* ..
322* .. Local Scalars ..
323 LOGICAL LQUERY, NOTRAN
324 INTEGER I, IE, IFUNC, IROUND, IS, ISOLVE, J, JE, JS, K,
325 $ LINFO, LWMIN, MB, NB, P, PQ, Q
326 DOUBLE PRECISION DSCALE, DSUM, SCALE2, SCALOC
327* ..
328* .. External Functions ..
329 LOGICAL LSAME
330 INTEGER ILAENV
331 EXTERNAL lsame, ilaenv
332* ..
333* .. External Subroutines ..
334 EXTERNAL xerbla, zgemm, zlacpy, zlaset, zscal,
335 $ ztgsy2
336* ..
337* .. Intrinsic Functions ..
338 INTRINSIC dble, dcmplx, max, sqrt
339* ..
340* .. Executable Statements ..
341*
342* Decode and test input parameters
343*
344 info = 0
345 notran = lsame( trans, 'N' )
346 lquery = ( lwork.EQ.-1 )
347*
348 IF( .NOT.notran .AND. .NOT.lsame( trans, 'C' ) ) THEN
349 info = -1
350 ELSE IF( notran ) THEN
351 IF( ( ijob.LT.0 ) .OR. ( ijob.GT.4 ) ) THEN
352 info = -2
353 END IF
354 END IF
355 IF( info.EQ.0 ) THEN
356 IF( m.LE.0 ) THEN
357 info = -3
358 ELSE IF( n.LE.0 ) THEN
359 info = -4
360 ELSE IF( lda.LT.max( 1, m ) ) THEN
361 info = -6
362 ELSE IF( ldb.LT.max( 1, n ) ) THEN
363 info = -8
364 ELSE IF( ldc.LT.max( 1, m ) ) THEN
365 info = -10
366 ELSE IF( ldd.LT.max( 1, m ) ) THEN
367 info = -12
368 ELSE IF( lde.LT.max( 1, n ) ) THEN
369 info = -14
370 ELSE IF( ldf.LT.max( 1, m ) ) THEN
371 info = -16
372 END IF
373 END IF
374*
375 IF( info.EQ.0 ) THEN
376 IF( notran ) THEN
377 IF( ijob.EQ.1 .OR. ijob.EQ.2 ) THEN
378 lwmin = max( 1, 2*m*n )
379 ELSE
380 lwmin = 1
381 END IF
382 ELSE
383 lwmin = 1
384 END IF
385 work( 1 ) = lwmin
386*
387 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
388 info = -20
389 END IF
390 END IF
391*
392 IF( info.NE.0 ) THEN
393 CALL xerbla( 'ZTGSYL', -info )
394 RETURN
395 ELSE IF( lquery ) THEN
396 RETURN
397 END IF
398*
399* Quick return if possible
400*
401 IF( m.EQ.0 .OR. n.EQ.0 ) THEN
402 scale = 1
403 IF( notran ) THEN
404 IF( ijob.NE.0 ) THEN
405 dif = 0
406 END IF
407 END IF
408 RETURN
409 END IF
410*
411* Determine optimal block sizes MB and NB
412*
413 mb = ilaenv( 2, 'ZTGSYL', trans, m, n, -1, -1 )
414 nb = ilaenv( 5, 'ZTGSYL', trans, m, n, -1, -1 )
415*
416 isolve = 1
417 ifunc = 0
418 IF( notran ) THEN
419 IF( ijob.GE.3 ) THEN
420 ifunc = ijob - 2
421 CALL zlaset( 'F', m, n, czero, czero, c, ldc )
422 CALL zlaset( 'F', m, n, czero, czero, f, ldf )
423 ELSE IF( ijob.GE.1 .AND. notran ) THEN
424 isolve = 2
425 END IF
426 END IF
427*
428 IF( ( mb.LE.1 .AND. nb.LE.1 ) .OR. ( mb.GE.m .AND. nb.GE.n ) )
429 $ THEN
430*
431* Use unblocked Level 2 solver
432*
433 DO 30 iround = 1, isolve
434*
435 scale = one
436 dscale = zero
437 dsum = one
438 pq = m*n
439 CALL ztgsy2( trans, ifunc, m, n, a, lda, b, ldb, c, ldc,
440 $ d,
441 $ ldd, e, lde, f, ldf, scale, dsum, dscale,
442 $ info )
443 IF( dscale.NE.zero ) THEN
444 IF( ijob.EQ.1 .OR. ijob.EQ.3 ) THEN
445 dif = sqrt( dble( 2*m*n ) ) / ( dscale*sqrt( dsum ) )
446 ELSE
447 dif = sqrt( dble( pq ) ) / ( dscale*sqrt( dsum ) )
448 END IF
449 END IF
450 IF( isolve.EQ.2 .AND. iround.EQ.1 ) THEN
451 IF( notran ) THEN
452 ifunc = ijob
453 END IF
454 scale2 = scale
455 CALL zlacpy( 'F', m, n, c, ldc, work, m )
456 CALL zlacpy( 'F', m, n, f, ldf, work( m*n+1 ), m )
457 CALL zlaset( 'F', m, n, czero, czero, c, ldc )
458 CALL zlaset( 'F', m, n, czero, czero, f, ldf )
459 ELSE IF( isolve.EQ.2 .AND. iround.EQ.2 ) THEN
460 CALL zlacpy( 'F', m, n, work, m, c, ldc )
461 CALL zlacpy( 'F', m, n, work( m*n+1 ), m, f, ldf )
462 scale = scale2
463 END IF
464 30 CONTINUE
465*
466 RETURN
467*
468 END IF
469*
470* Determine block structure of A
471*
472 p = 0
473 i = 1
474 40 CONTINUE
475 IF( i.GT.m )
476 $ GO TO 50
477 p = p + 1
478 iwork( p ) = i
479 i = i + mb
480 IF( i.GE.m )
481 $ GO TO 50
482 GO TO 40
483 50 CONTINUE
484 iwork( p+1 ) = m + 1
485 IF( iwork( p ).EQ.iwork( p+1 ) )
486 $ p = p - 1
487*
488* Determine block structure of B
489*
490 q = p + 1
491 j = 1
492 60 CONTINUE
493 IF( j.GT.n )
494 $ GO TO 70
495*
496 q = q + 1
497 iwork( q ) = j
498 j = j + nb
499 IF( j.GE.n )
500 $ GO TO 70
501 GO TO 60
502*
503 70 CONTINUE
504 iwork( q+1 ) = n + 1
505 IF( iwork( q ).EQ.iwork( q+1 ) )
506 $ q = q - 1
507*
508 IF( notran ) THEN
509 DO 150 iround = 1, isolve
510*
511* Solve (I, J) - subsystem
512* A(I, I) * R(I, J) - L(I, J) * B(J, J) = C(I, J)
513* D(I, I) * R(I, J) - L(I, J) * E(J, J) = F(I, J)
514* for I = P, P - 1, ..., 1; J = 1, 2, ..., Q
515*
516 pq = 0
517 scale = one
518 dscale = zero
519 dsum = one
520 DO 130 j = p + 2, q
521 js = iwork( j )
522 je = iwork( j+1 ) - 1
523 nb = je - js + 1
524 DO 120 i = p, 1, -1
525 is = iwork( i )
526 ie = iwork( i+1 ) - 1
527 mb = ie - is + 1
528 CALL ztgsy2( trans, ifunc, mb, nb, a( is, is ),
529 $ lda,
530 $ b( js, js ), ldb, c( is, js ), ldc,
531 $ d( is, is ), ldd, e( js, js ), lde,
532 $ f( is, js ), ldf, scaloc, dsum, dscale,
533 $ linfo )
534 IF( linfo.GT.0 )
535 $ info = linfo
536 pq = pq + mb*nb
537 IF( scaloc.NE.one ) THEN
538 DO 80 k = 1, js - 1
539 CALL zscal( m, dcmplx( scaloc, zero ),
540 $ c( 1, k ), 1 )
541 CALL zscal( m, dcmplx( scaloc, zero ),
542 $ f( 1, k ), 1 )
543 80 CONTINUE
544 DO 90 k = js, je
545 CALL zscal( is-1, dcmplx( scaloc, zero ),
546 $ c( 1, k ), 1 )
547 CALL zscal( is-1, dcmplx( scaloc, zero ),
548 $ f( 1, k ), 1 )
549 90 CONTINUE
550 DO 100 k = js, je
551 CALL zscal( m-ie, dcmplx( scaloc, zero ),
552 $ c( ie+1, k ), 1 )
553 CALL zscal( m-ie, dcmplx( scaloc, zero ),
554 $ f( ie+1, k ), 1 )
555 100 CONTINUE
556 DO 110 k = je + 1, n
557 CALL zscal( m, dcmplx( scaloc, zero ),
558 $ c( 1, k ), 1 )
559 CALL zscal( m, dcmplx( scaloc, zero ),
560 $ f( 1, k ), 1 )
561 110 CONTINUE
562 scale = scale*scaloc
563 END IF
564*
565* Substitute R(I,J) and L(I,J) into remaining equation.
566*
567 IF( i.GT.1 ) THEN
568 CALL zgemm( 'N', 'N', is-1, nb, mb,
569 $ dcmplx( -one, zero ), a( 1, is ), lda,
570 $ c( is, js ), ldc, dcmplx( one, zero ),
571 $ c( 1, js ), ldc )
572 CALL zgemm( 'N', 'N', is-1, nb, mb,
573 $ dcmplx( -one, zero ), d( 1, is ), ldd,
574 $ c( is, js ), ldc, dcmplx( one, zero ),
575 $ f( 1, js ), ldf )
576 END IF
577 IF( j.LT.q ) THEN
578 CALL zgemm( 'N', 'N', mb, n-je, nb,
579 $ dcmplx( one, zero ), f( is, js ), ldf,
580 $ b( js, je+1 ), ldb,
581 $ dcmplx( one, zero ), c( is, je+1 ),
582 $ ldc )
583 CALL zgemm( 'N', 'N', mb, n-je, nb,
584 $ dcmplx( one, zero ), f( is, js ), ldf,
585 $ e( js, je+1 ), lde,
586 $ dcmplx( one, zero ), f( is, je+1 ),
587 $ ldf )
588 END IF
589 120 CONTINUE
590 130 CONTINUE
591 IF( dscale.NE.zero ) THEN
592 IF( ijob.EQ.1 .OR. ijob.EQ.3 ) THEN
593 dif = sqrt( dble( 2*m*n ) ) / ( dscale*sqrt( dsum ) )
594 ELSE
595 dif = sqrt( dble( pq ) ) / ( dscale*sqrt( dsum ) )
596 END IF
597 END IF
598 IF( isolve.EQ.2 .AND. iround.EQ.1 ) THEN
599 IF( notran ) THEN
600 ifunc = ijob
601 END IF
602 scale2 = scale
603 CALL zlacpy( 'F', m, n, c, ldc, work, m )
604 CALL zlacpy( 'F', m, n, f, ldf, work( m*n+1 ), m )
605 CALL zlaset( 'F', m, n, czero, czero, c, ldc )
606 CALL zlaset( 'F', m, n, czero, czero, f, ldf )
607 ELSE IF( isolve.EQ.2 .AND. iround.EQ.2 ) THEN
608 CALL zlacpy( 'F', m, n, work, m, c, ldc )
609 CALL zlacpy( 'F', m, n, work( m*n+1 ), m, f, ldf )
610 scale = scale2
611 END IF
612 150 CONTINUE
613 ELSE
614*
615* Solve transposed (I, J)-subsystem
616* A(I, I)**H * R(I, J) + D(I, I)**H * L(I, J) = C(I, J)
617* R(I, J) * B(J, J) + L(I, J) * E(J, J) = -F(I, J)
618* for I = 1,2,..., P; J = Q, Q-1,..., 1
619*
620 scale = one
621 DO 210 i = 1, p
622 is = iwork( i )
623 ie = iwork( i+1 ) - 1
624 mb = ie - is + 1
625 DO 200 j = q, p + 2, -1
626 js = iwork( j )
627 je = iwork( j+1 ) - 1
628 nb = je - js + 1
629 CALL ztgsy2( trans, ifunc, mb, nb, a( is, is ), lda,
630 $ b( js, js ), ldb, c( is, js ), ldc,
631 $ d( is, is ), ldd, e( js, js ), lde,
632 $ f( is, js ), ldf, scaloc, dsum, dscale,
633 $ linfo )
634 IF( linfo.GT.0 )
635 $ info = linfo
636 IF( scaloc.NE.one ) THEN
637 DO 160 k = 1, js - 1
638 CALL zscal( m, dcmplx( scaloc, zero ), c( 1,
639 $ k ),
640 $ 1 )
641 CALL zscal( m, dcmplx( scaloc, zero ), f( 1,
642 $ k ),
643 $ 1 )
644 160 CONTINUE
645 DO 170 k = js, je
646 CALL zscal( is-1, dcmplx( scaloc, zero ),
647 $ c( 1, k ), 1 )
648 CALL zscal( is-1, dcmplx( scaloc, zero ),
649 $ f( 1, k ), 1 )
650 170 CONTINUE
651 DO 180 k = js, je
652 CALL zscal( m-ie, dcmplx( scaloc, zero ),
653 $ c( ie+1, k ), 1 )
654 CALL zscal( m-ie, dcmplx( scaloc, zero ),
655 $ f( ie+1, k ), 1 )
656 180 CONTINUE
657 DO 190 k = je + 1, n
658 CALL zscal( m, dcmplx( scaloc, zero ), c( 1,
659 $ k ),
660 $ 1 )
661 CALL zscal( m, dcmplx( scaloc, zero ), f( 1,
662 $ k ),
663 $ 1 )
664 190 CONTINUE
665 scale = scale*scaloc
666 END IF
667*
668* Substitute R(I,J) and L(I,J) into remaining equation.
669*
670 IF( j.GT.p+2 ) THEN
671 CALL zgemm( 'N', 'C', mb, js-1, nb,
672 $ dcmplx( one, zero ), c( is, js ), ldc,
673 $ b( 1, js ), ldb, dcmplx( one, zero ),
674 $ f( is, 1 ), ldf )
675 CALL zgemm( 'N', 'C', mb, js-1, nb,
676 $ dcmplx( one, zero ), f( is, js ), ldf,
677 $ e( 1, js ), lde, dcmplx( one, zero ),
678 $ f( is, 1 ), ldf )
679 END IF
680 IF( i.LT.p ) THEN
681 CALL zgemm( 'C', 'N', m-ie, nb, mb,
682 $ dcmplx( -one, zero ), a( is, ie+1 ), lda,
683 $ c( is, js ), ldc, dcmplx( one, zero ),
684 $ c( ie+1, js ), ldc )
685 CALL zgemm( 'C', 'N', m-ie, nb, mb,
686 $ dcmplx( -one, zero ), d( is, ie+1 ), ldd,
687 $ f( is, js ), ldf, dcmplx( one, zero ),
688 $ c( ie+1, js ), ldc )
689 END IF
690 200 CONTINUE
691 210 CONTINUE
692 END IF
693*
694 work( 1 ) = lwmin
695*
696 RETURN
697*
698* End of ZTGSYL
699*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
logical function lde(ri, rj, lr)
Definition dblat2.f:2970
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:104
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zscal(n, za, zx, incx)
ZSCAL
Definition zscal.f:78
subroutine ztgsy2(trans, ijob, m, n, a, lda, b, ldb, c, ldc, d, ldd, e, lde, f, ldf, scale, rdsum, rdscal, info)
ZTGSY2 solves the generalized Sylvester equation (unblocked algorithm).
Definition ztgsy2.f:258
Here is the call graph for this function:
Here is the caller graph for this function: