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

◆ ctgsyl()

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

CTGSYL

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

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