File:  [local] / rpl / lapack / lapack / dtgsyl.f
Revision 1.19: download - view: text, annotated - select for diffs - revision graph
Mon Aug 7 08:39:12 2023 UTC (8 months, 3 weeks ago) by bertrand
Branches: MAIN
CVS tags: rpl-4_1_35, rpl-4_1_34, HEAD
Première mise à jour de lapack et blas.

    1: *> \brief \b DTGSYL
    2: *
    3: *  =========== DOCUMENTATION ===========
    4: *
    5: * Online html documentation available at
    6: *            http://www.netlib.org/lapack/explore-html/
    7: *
    8: *> \htmlonly
    9: *> Download DTGSYL + dependencies
   10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtgsyl.f">
   11: *> [TGZ]</a>
   12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtgsyl.f">
   13: *> [ZIP]</a>
   14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtgsyl.f">
   15: *> [TXT]</a>
   16: *> \endhtmlonly
   17: *
   18: *  Definition:
   19: *  ===========
   20: *
   21: *       SUBROUTINE DTGSYL( TRANS, IJOB, M, N, A, LDA, B, LDB, C, LDC, D,
   22: *                          LDD, E, LDE, F, LDF, SCALE, DIF, WORK, LWORK,
   23: *                          IWORK, INFO )
   24: *
   25: *       .. Scalar Arguments ..
   26: *       CHARACTER          TRANS
   27: *       INTEGER            IJOB, INFO, LDA, LDB, LDC, LDD, LDE, LDF,
   28: *      $                   LWORK, M, N
   29: *       DOUBLE PRECISION   DIF, SCALE
   30: *       ..
   31: *       .. Array Arguments ..
   32: *       INTEGER            IWORK( * )
   33: *       DOUBLE PRECISION   A( LDA, * ), B( LDB, * ), C( LDC, * ),
   34: *      $                   D( LDD, * ), E( LDE, * ), F( LDF, * ),
   35: *      $                   WORK( * )
   36: *       ..
   37: *
   38: *
   39: *> \par Purpose:
   40: *  =============
   41: *>
   42: *> \verbatim
   43: *>
   44: *> DTGSYL solves the generalized Sylvester equation:
   45: *>
   46: *>             A * R - L * B = scale * C                 (1)
   47: *>             D * R - L * E = scale * F
   48: *>
   49: *> where R and L are unknown m-by-n matrices, (A, D), (B, E) and
   50: *> (C, F) are given matrix pairs of size m-by-m, n-by-n and m-by-n,
   51: *> respectively, with real entries. (A, D) and (B, E) must be in
   52: *> generalized (real) Schur canonical form, i.e. A, B are upper quasi
   53: *> triangular and D, E are upper triangular.
   54: *>
   55: *> The solution (R, L) overwrites (C, F). 0 <= SCALE <= 1 is an output
   56: *> scaling factor chosen to avoid overflow.
   57: *>
   58: *> In matrix notation (1) is equivalent to solve  Zx = scale b, where
   59: *> Z is defined as
   60: *>
   61: *>            Z = [ kron(In, A)  -kron(B**T, Im) ]         (2)
   62: *>                [ kron(In, D)  -kron(E**T, Im) ].
   63: *>
   64: *> Here Ik is the identity matrix of size k and X**T is the transpose of
   65: *> X. kron(X, Y) is the Kronecker product between the matrices X and Y.
   66: *>
   67: *> If TRANS = 'T', DTGSYL solves the transposed system Z**T*y = scale*b,
   68: *> which is equivalent to solve for R and L in
   69: *>
   70: *>             A**T * R + D**T * L = scale * C           (3)
   71: *>             R * B**T + L * E**T = scale * -F
   72: *>
   73: *> This case (TRANS = 'T') is used to compute an one-norm-based estimate
   74: *> of Dif[(A,D), (B,E)], the separation between the matrix pairs (A,D)
   75: *> and (B,E), using DLACON.
   76: *>
   77: *> If IJOB >= 1, DTGSYL computes a Frobenius norm-based estimate
   78: *> of Dif[(A,D),(B,E)]. That is, the reciprocal of a lower bound on the
   79: *> reciprocal of the smallest singular value of Z. See [1-2] for more
   80: *> information.
   81: *>
   82: *> This is a level 3 BLAS algorithm.
   83: *> \endverbatim
   84: *
   85: *  Arguments:
   86: *  ==========
   87: *
   88: *> \param[in] TRANS
   89: *> \verbatim
   90: *>          TRANS is CHARACTER*1
   91: *>          = 'N': solve the generalized Sylvester equation (1).
   92: *>          = 'T': solve the 'transposed' system (3).
   93: *> \endverbatim
   94: *>
   95: *> \param[in] IJOB
   96: *> \verbatim
   97: *>          IJOB is INTEGER
   98: *>          Specifies what kind of functionality to be performed.
   99: *>          = 0: solve (1) only.
  100: *>          = 1: The functionality of 0 and 3.
  101: *>          = 2: The functionality of 0 and 4.
  102: *>          = 3: Only an estimate of Dif[(A,D), (B,E)] is computed.
  103: *>               (look ahead strategy IJOB  = 1 is used).
  104: *>          = 4: Only an estimate of Dif[(A,D), (B,E)] is computed.
  105: *>               ( DGECON on sub-systems is used ).
  106: *>          Not referenced if TRANS = 'T'.
  107: *> \endverbatim
  108: *>
  109: *> \param[in] M
  110: *> \verbatim
  111: *>          M is INTEGER
  112: *>          The order of the matrices A and D, and the row dimension of
  113: *>          the matrices C, F, R and L.
  114: *> \endverbatim
  115: *>
  116: *> \param[in] N
  117: *> \verbatim
  118: *>          N is INTEGER
  119: *>          The order of the matrices B and E, and the column dimension
  120: *>          of the matrices C, F, R and L.
  121: *> \endverbatim
  122: *>
  123: *> \param[in] A
  124: *> \verbatim
  125: *>          A is DOUBLE PRECISION array, dimension (LDA, M)
  126: *>          The upper quasi triangular matrix A.
  127: *> \endverbatim
  128: *>
  129: *> \param[in] LDA
  130: *> \verbatim
  131: *>          LDA is INTEGER
  132: *>          The leading dimension of the array A. LDA >= max(1, M).
  133: *> \endverbatim
  134: *>
  135: *> \param[in] B
  136: *> \verbatim
  137: *>          B is DOUBLE PRECISION array, dimension (LDB, N)
  138: *>          The upper quasi triangular matrix B.
  139: *> \endverbatim
  140: *>
  141: *> \param[in] LDB
  142: *> \verbatim
  143: *>          LDB is INTEGER
  144: *>          The leading dimension of the array B. LDB >= max(1, N).
  145: *> \endverbatim
  146: *>
  147: *> \param[in,out] C
  148: *> \verbatim
  149: *>          C is DOUBLE PRECISION array, dimension (LDC, N)
  150: *>          On entry, C contains the right-hand-side of the first matrix
  151: *>          equation in (1) or (3).
  152: *>          On exit, if IJOB = 0, 1 or 2, C has been overwritten by
  153: *>          the solution R. If IJOB = 3 or 4 and TRANS = 'N', C holds R,
  154: *>          the solution achieved during the computation of the
  155: *>          Dif-estimate.
  156: *> \endverbatim
  157: *>
  158: *> \param[in] LDC
  159: *> \verbatim
  160: *>          LDC is INTEGER
  161: *>          The leading dimension of the array C. LDC >= max(1, M).
  162: *> \endverbatim
  163: *>
  164: *> \param[in] D
  165: *> \verbatim
  166: *>          D is DOUBLE PRECISION array, dimension (LDD, M)
  167: *>          The upper triangular matrix D.
  168: *> \endverbatim
  169: *>
  170: *> \param[in] LDD
  171: *> \verbatim
  172: *>          LDD is INTEGER
  173: *>          The leading dimension of the array D. LDD >= max(1, M).
  174: *> \endverbatim
  175: *>
  176: *> \param[in] E
  177: *> \verbatim
  178: *>          E is DOUBLE PRECISION array, dimension (LDE, N)
  179: *>          The upper triangular matrix E.
  180: *> \endverbatim
  181: *>
  182: *> \param[in] LDE
  183: *> \verbatim
  184: *>          LDE is INTEGER
  185: *>          The leading dimension of the array E. LDE >= max(1, N).
  186: *> \endverbatim
  187: *>
  188: *> \param[in,out] F
  189: *> \verbatim
  190: *>          F is DOUBLE PRECISION array, dimension (LDF, N)
  191: *>          On entry, F contains the right-hand-side of the second matrix
  192: *>          equation in (1) or (3).
  193: *>          On exit, if IJOB = 0, 1 or 2, F has been overwritten by
  194: *>          the solution L. If IJOB = 3 or 4 and TRANS = 'N', F holds L,
  195: *>          the solution achieved during the computation of the
  196: *>          Dif-estimate.
  197: *> \endverbatim
  198: *>
  199: *> \param[in] LDF
  200: *> \verbatim
  201: *>          LDF is INTEGER
  202: *>          The leading dimension of the array F. LDF >= max(1, M).
  203: *> \endverbatim
  204: *>
  205: *> \param[out] DIF
  206: *> \verbatim
  207: *>          DIF is DOUBLE PRECISION
  208: *>          On exit DIF is the reciprocal of a lower bound of the
  209: *>          reciprocal of the Dif-function, i.e. DIF is an upper bound of
  210: *>          Dif[(A,D), (B,E)] = sigma_min(Z), where Z as in (2).
  211: *>          IF IJOB = 0 or TRANS = 'T', DIF is not touched.
  212: *> \endverbatim
  213: *>
  214: *> \param[out] SCALE
  215: *> \verbatim
  216: *>          SCALE is DOUBLE PRECISION
  217: *>          On exit SCALE is the scaling factor in (1) or (3).
  218: *>          If 0 < SCALE < 1, C and F hold the solutions R and L, resp.,
  219: *>          to a slightly perturbed system but the input matrices A, B, D
  220: *>          and E have not been changed. If SCALE = 0, C and F hold the
  221: *>          solutions R and L, respectively, to the homogeneous system
  222: *>          with C = F = 0. Normally, SCALE = 1.
  223: *> \endverbatim
  224: *>
  225: *> \param[out] WORK
  226: *> \verbatim
  227: *>          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
  228: *>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
  229: *> \endverbatim
  230: *>
  231: *> \param[in] LWORK
  232: *> \verbatim
  233: *>          LWORK is INTEGER
  234: *>          The dimension of the array WORK. LWORK > = 1.
  235: *>          If IJOB = 1 or 2 and TRANS = 'N', LWORK >= max(1,2*M*N).
  236: *>
  237: *>          If LWORK = -1, then a workspace query is assumed; the routine
  238: *>          only calculates the optimal size of the WORK array, returns
  239: *>          this value as the first entry of the WORK array, and no error
  240: *>          message related to LWORK is issued by XERBLA.
  241: *> \endverbatim
  242: *>
  243: *> \param[out] IWORK
  244: *> \verbatim
  245: *>          IWORK is INTEGER array, dimension (M+N+6)
  246: *> \endverbatim
  247: *>
  248: *> \param[out] INFO
  249: *> \verbatim
  250: *>          INFO is INTEGER
  251: *>            =0: successful exit
  252: *>            <0: If INFO = -i, the i-th argument had an illegal value.
  253: *>            >0: (A, D) and (B, E) have common or close eigenvalues.
  254: *> \endverbatim
  255: *
  256: *  Authors:
  257: *  ========
  258: *
  259: *> \author Univ. of Tennessee
  260: *> \author Univ. of California Berkeley
  261: *> \author Univ. of Colorado Denver
  262: *> \author NAG Ltd.
  263: *
  264: *> \ingroup doubleSYcomputational
  265: *
  266: *> \par Contributors:
  267: *  ==================
  268: *>
  269: *>     Bo Kagstrom and Peter Poromaa, Department of Computing Science,
  270: *>     Umea University, S-901 87 Umea, Sweden.
  271: *
  272: *> \par References:
  273: *  ================
  274: *>
  275: *> \verbatim
  276: *>
  277: *>  [1] B. Kagstrom and P. Poromaa, LAPACK-Style Algorithms and Software
  278: *>      for Solving the Generalized Sylvester Equation and Estimating the
  279: *>      Separation between Regular Matrix Pairs, Report UMINF - 93.23,
  280: *>      Department of Computing Science, Umea University, S-901 87 Umea,
  281: *>      Sweden, December 1993, Revised April 1994, Also as LAPACK Working
  282: *>      Note 75.  To appear in ACM Trans. on Math. Software, Vol 22,
  283: *>      No 1, 1996.
  284: *>
  285: *>  [2] B. Kagstrom, A Perturbation Analysis of the Generalized Sylvester
  286: *>      Equation (AR - LB, DR - LE ) = (C, F), SIAM J. Matrix Anal.
  287: *>      Appl., 15(4):1045-1060, 1994
  288: *>
  289: *>  [3] B. Kagstrom and L. Westin, Generalized Schur Methods with
  290: *>      Condition Estimators for Solving the Generalized Sylvester
  291: *>      Equation, IEEE Transactions on Automatic Control, Vol. 34, No. 7,
  292: *>      July 1989, pp 745-751.
  293: *> \endverbatim
  294: *>
  295: *  =====================================================================
  296:       SUBROUTINE DTGSYL( TRANS, IJOB, M, N, A, LDA, B, LDB, C, LDC, D,
  297:      $                   LDD, E, LDE, F, LDF, SCALE, DIF, WORK, LWORK,
  298:      $                   IWORK, INFO )
  299: *
  300: *  -- LAPACK computational routine --
  301: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
  302: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  303: *
  304: *     .. Scalar Arguments ..
  305:       CHARACTER          TRANS
  306:       INTEGER            IJOB, INFO, LDA, LDB, LDC, LDD, LDE, LDF,
  307:      $                   LWORK, M, N
  308:       DOUBLE PRECISION   DIF, SCALE
  309: *     ..
  310: *     .. Array Arguments ..
  311:       INTEGER            IWORK( * )
  312:       DOUBLE PRECISION   A( LDA, * ), B( LDB, * ), C( LDC, * ),
  313:      $                   D( LDD, * ), E( LDE, * ), F( LDF, * ),
  314:      $                   WORK( * )
  315: *     ..
  316: *
  317: *  =====================================================================
  318: *  Replaced various illegal calls to DCOPY by calls to DLASET.
  319: *  Sven Hammarling, 1/5/02.
  320: *
  321: *     .. Parameters ..
  322:       DOUBLE PRECISION   ZERO, ONE
  323:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
  324: *     ..
  325: *     .. Local Scalars ..
  326:       LOGICAL            LQUERY, NOTRAN
  327:       INTEGER            I, IE, IFUNC, IROUND, IS, ISOLVE, J, JE, JS, K,
  328:      $                   LINFO, LWMIN, MB, NB, P, PPQQ, PQ, Q
  329:       DOUBLE PRECISION   DSCALE, DSUM, SCALE2, SCALOC
  330: *     ..
  331: *     .. External Functions ..
  332:       LOGICAL            LSAME
  333:       INTEGER            ILAENV
  334:       EXTERNAL           LSAME, ILAENV
  335: *     ..
  336: *     .. External Subroutines ..
  337:       EXTERNAL           DGEMM, DLACPY, DLASET, DSCAL, DTGSY2, XERBLA
  338: *     ..
  339: *     .. Intrinsic Functions ..
  340:       INTRINSIC          DBLE, MAX, SQRT
  341: *     ..
  342: *     .. Executable Statements ..
  343: *
  344: *     Decode and test input parameters
  345: *
  346:       INFO = 0
  347:       NOTRAN = LSAME( TRANS, 'N' )
  348:       LQUERY = ( LWORK.EQ.-1 )
  349: *
  350:       IF( .NOT.NOTRAN .AND. .NOT.LSAME( TRANS, 'T' ) ) THEN
  351:          INFO = -1
  352:       ELSE IF( NOTRAN ) THEN
  353:          IF( ( IJOB.LT.0 ) .OR. ( IJOB.GT.4 ) ) THEN
  354:             INFO = -2
  355:          END IF
  356:       END IF
  357:       IF( INFO.EQ.0 ) THEN
  358:          IF( M.LE.0 ) THEN
  359:             INFO = -3
  360:          ELSE IF( N.LE.0 ) THEN
  361:             INFO = -4
  362:          ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
  363:             INFO = -6
  364:          ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
  365:             INFO = -8
  366:          ELSE IF( LDC.LT.MAX( 1, M ) ) THEN
  367:             INFO = -10
  368:          ELSE IF( LDD.LT.MAX( 1, M ) ) THEN
  369:             INFO = -12
  370:          ELSE IF( LDE.LT.MAX( 1, N ) ) THEN
  371:             INFO = -14
  372:          ELSE IF( LDF.LT.MAX( 1, M ) ) THEN
  373:             INFO = -16
  374:          END IF
  375:       END IF
  376: *
  377:       IF( INFO.EQ.0 ) THEN
  378:          IF( NOTRAN ) THEN
  379:             IF( IJOB.EQ.1 .OR. IJOB.EQ.2 ) THEN
  380:                LWMIN = MAX( 1, 2*M*N )
  381:             ELSE
  382:                LWMIN = 1
  383:             END IF
  384:          ELSE
  385:             LWMIN = 1
  386:          END IF
  387:          WORK( 1 ) = LWMIN
  388: *
  389:          IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY ) THEN
  390:             INFO = -20
  391:          END IF
  392:       END IF
  393: *
  394:       IF( INFO.NE.0 ) THEN
  395:          CALL XERBLA( 'DTGSYL', -INFO )
  396:          RETURN
  397:       ELSE IF( LQUERY ) THEN
  398:          RETURN
  399:       END IF
  400: *
  401: *     Quick return if possible
  402: *
  403:       IF( M.EQ.0 .OR. N.EQ.0 ) THEN
  404:          SCALE = 1
  405:          IF( NOTRAN ) THEN
  406:             IF( IJOB.NE.0 ) THEN
  407:                DIF = 0
  408:             END IF
  409:          END IF
  410:          RETURN
  411:       END IF
  412: *
  413: *     Determine optimal block sizes MB and NB
  414: *
  415:       MB = ILAENV( 2, 'DTGSYL', TRANS, M, N, -1, -1 )
  416:       NB = ILAENV( 5, 'DTGSYL', TRANS, M, N, -1, -1 )
  417: *
  418:       ISOLVE = 1
  419:       IFUNC = 0
  420:       IF( NOTRAN ) THEN
  421:          IF( IJOB.GE.3 ) THEN
  422:             IFUNC = IJOB - 2
  423:             CALL DLASET( 'F', M, N, ZERO, ZERO, C, LDC )
  424:             CALL DLASET( 'F', M, N, ZERO, ZERO, F, LDF )
  425:          ELSE IF( IJOB.GE.1 ) THEN
  426:             ISOLVE = 2
  427:          END IF
  428:       END IF
  429: *
  430:       IF( ( MB.LE.1 .AND. NB.LE.1 ) .OR. ( MB.GE.M .AND. NB.GE.N ) )
  431:      $     THEN
  432: *
  433:          DO 30 IROUND = 1, ISOLVE
  434: *
  435: *           Use unblocked Level 2 solver
  436: *
  437:             DSCALE = ZERO
  438:             DSUM = ONE
  439:             PQ = 0
  440:             CALL DTGSY2( TRANS, IFUNC, M, N, A, LDA, B, LDB, C, LDC, D,
  441:      $                   LDD, E, LDE, F, LDF, SCALE, DSUM, DSCALE,
  442:      $                   IWORK, PQ, 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: *
  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 DLACPY( 'F', M, N, C, LDC, WORK, M )
  457:                CALL DLACPY( 'F', M, N, F, LDF, WORK( M*N+1 ), M )
  458:                CALL DLASET( 'F', M, N, ZERO, ZERO, C, LDC )
  459:                CALL DLASET( 'F', M, N, ZERO, ZERO, F, LDF )
  460:             ELSE IF( ISOLVE.EQ.2 .AND. IROUND.EQ.2 ) THEN
  461:                CALL DLACPY( 'F', M, N, WORK, M, C, LDC )
  462:                CALL DLACPY( 'F', M, N, WORK( M*N+1 ), M, F, LDF )
  463:                SCALE = SCALE2
  464:             END IF
  465:    30    CONTINUE
  466: *
  467:          RETURN
  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:       IF( A( I, I-1 ).NE.ZERO )
  483:      $   I = I + 1
  484:       GO TO 40
  485:    50 CONTINUE
  486: *
  487:       IWORK( P+1 ) = M + 1
  488:       IF( IWORK( P ).EQ.IWORK( P+1 ) )
  489:      $   P = P - 1
  490: *
  491: *     Determine block structure of B
  492: *
  493:       Q = P + 1
  494:       J = 1
  495:    60 CONTINUE
  496:       IF( J.GT.N )
  497:      $   GO TO 70
  498:       Q = Q + 1
  499:       IWORK( Q ) = J
  500:       J = J + NB
  501:       IF( J.GE.N )
  502:      $   GO TO 70
  503:       IF( B( J, J-1 ).NE.ZERO )
  504:      $   J = J + 1
  505:       GO TO 60
  506:    70 CONTINUE
  507: *
  508:       IWORK( Q+1 ) = N + 1
  509:       IF( IWORK( Q ).EQ.IWORK( Q+1 ) )
  510:      $   Q = Q - 1
  511: *
  512:       IF( NOTRAN ) THEN
  513: *
  514:          DO 150 IROUND = 1, ISOLVE
  515: *
  516: *           Solve (I, J)-subsystem
  517: *               A(I, I) * R(I, J) - L(I, J) * B(J, J) = C(I, J)
  518: *               D(I, I) * R(I, J) - L(I, J) * E(J, J) = F(I, J)
  519: *           for I = P, P - 1,..., 1; J = 1, 2,..., Q
  520: *
  521:             DSCALE = ZERO
  522:             DSUM = ONE
  523:             PQ = 0
  524:             SCALE = ONE
  525:             DO 130 J = P + 2, Q
  526:                JS = IWORK( J )
  527:                JE = IWORK( J+1 ) - 1
  528:                NB = JE - JS + 1
  529:                DO 120 I = P, 1, -1
  530:                   IS = IWORK( I )
  531:                   IE = IWORK( I+1 ) - 1
  532:                   MB = IE - IS + 1
  533:                   PPQQ = 0
  534:                   CALL DTGSY2( TRANS, IFUNC, MB, NB, A( IS, IS ), LDA,
  535:      $                         B( JS, JS ), LDB, C( IS, JS ), LDC,
  536:      $                         D( IS, IS ), LDD, E( JS, JS ), LDE,
  537:      $                         F( IS, JS ), LDF, SCALOC, DSUM, DSCALE,
  538:      $                         IWORK( Q+2 ), PPQQ, LINFO )
  539:                   IF( LINFO.GT.0 )
  540:      $               INFO = LINFO
  541: *
  542:                   PQ = PQ + PPQQ
  543:                   IF( SCALOC.NE.ONE ) THEN
  544:                      DO 80 K = 1, JS - 1
  545:                         CALL DSCAL( M, SCALOC, C( 1, K ), 1 )
  546:                         CALL DSCAL( M, SCALOC, F( 1, K ), 1 )
  547:    80                CONTINUE
  548:                      DO 90 K = JS, JE
  549:                         CALL DSCAL( IS-1, SCALOC, C( 1, K ), 1 )
  550:                         CALL DSCAL( IS-1, SCALOC, F( 1, K ), 1 )
  551:    90                CONTINUE
  552:                      DO 100 K = JS, JE
  553:                         CALL DSCAL( M-IE, SCALOC, C( IE+1, K ), 1 )
  554:                         CALL DSCAL( M-IE, SCALOC, F( IE+1, K ), 1 )
  555:   100                CONTINUE
  556:                      DO 110 K = JE + 1, N
  557:                         CALL DSCAL( M, SCALOC, C( 1, K ), 1 )
  558:                         CALL DSCAL( M, SCALOC, F( 1, K ), 1 )
  559:   110                CONTINUE
  560:                      SCALE = SCALE*SCALOC
  561:                   END IF
  562: *
  563: *                 Substitute R(I, J) and L(I, J) into remaining
  564: *                 equation.
  565: *
  566:                   IF( I.GT.1 ) THEN
  567:                      CALL DGEMM( 'N', 'N', IS-1, NB, MB, -ONE,
  568:      $                           A( 1, IS ), LDA, C( IS, JS ), LDC, ONE,
  569:      $                           C( 1, JS ), LDC )
  570:                      CALL DGEMM( 'N', 'N', IS-1, NB, MB, -ONE,
  571:      $                           D( 1, IS ), LDD, C( IS, JS ), LDC, ONE,
  572:      $                           F( 1, JS ), LDF )
  573:                   END IF
  574:                   IF( J.LT.Q ) THEN
  575:                      CALL DGEMM( 'N', 'N', MB, N-JE, NB, ONE,
  576:      $                           F( IS, JS ), LDF, B( JS, JE+1 ), LDB,
  577:      $                           ONE, C( IS, JE+1 ), LDC )
  578:                      CALL DGEMM( 'N', 'N', MB, N-JE, NB, ONE,
  579:      $                           F( IS, JS ), LDF, E( JS, JE+1 ), LDE,
  580:      $                           ONE, F( IS, JE+1 ), LDF )
  581:                   END IF
  582:   120          CONTINUE
  583:   130       CONTINUE
  584:             IF( DSCALE.NE.ZERO ) THEN
  585:                IF( IJOB.EQ.1 .OR. IJOB.EQ.3 ) THEN
  586:                   DIF = SQRT( DBLE( 2*M*N ) ) / ( DSCALE*SQRT( DSUM ) )
  587:                ELSE
  588:                   DIF = SQRT( DBLE( PQ ) ) / ( DSCALE*SQRT( DSUM ) )
  589:                END IF
  590:             END IF
  591:             IF( ISOLVE.EQ.2 .AND. IROUND.EQ.1 ) THEN
  592:                IF( NOTRAN ) THEN
  593:                   IFUNC = IJOB
  594:                END IF
  595:                SCALE2 = SCALE
  596:                CALL DLACPY( 'F', M, N, C, LDC, WORK, M )
  597:                CALL DLACPY( 'F', M, N, F, LDF, WORK( M*N+1 ), M )
  598:                CALL DLASET( 'F', M, N, ZERO, ZERO, C, LDC )
  599:                CALL DLASET( 'F', M, N, ZERO, ZERO, F, LDF )
  600:             ELSE IF( ISOLVE.EQ.2 .AND. IROUND.EQ.2 ) THEN
  601:                CALL DLACPY( 'F', M, N, WORK, M, C, LDC )
  602:                CALL DLACPY( 'F', M, N, WORK( M*N+1 ), M, F, LDF )
  603:                SCALE = SCALE2
  604:             END IF
  605:   150    CONTINUE
  606: *
  607:       ELSE
  608: *
  609: *        Solve transposed (I, J)-subsystem
  610: *             A(I, I)**T * R(I, J)  + D(I, I)**T * L(I, J)  =  C(I, J)
  611: *             R(I, J)  * B(J, J)**T + L(I, J)  * E(J, J)**T = -F(I, J)
  612: *        for I = 1,2,..., P; J = Q, Q-1,..., 1
  613: *
  614:          SCALE = ONE
  615:          DO 210 I = 1, P
  616:             IS = IWORK( I )
  617:             IE = IWORK( I+1 ) - 1
  618:             MB = IE - IS + 1
  619:             DO 200 J = Q, P + 2, -1
  620:                JS = IWORK( J )
  621:                JE = IWORK( J+1 ) - 1
  622:                NB = JE - JS + 1
  623:                CALL DTGSY2( TRANS, IFUNC, MB, NB, A( IS, IS ), LDA,
  624:      $                      B( JS, JS ), LDB, C( IS, JS ), LDC,
  625:      $                      D( IS, IS ), LDD, E( JS, JS ), LDE,
  626:      $                      F( IS, JS ), LDF, SCALOC, DSUM, DSCALE,
  627:      $                      IWORK( Q+2 ), PPQQ, LINFO )
  628:                IF( LINFO.GT.0 )
  629:      $            INFO = LINFO
  630:                IF( SCALOC.NE.ONE ) THEN
  631:                   DO 160 K = 1, JS - 1
  632:                      CALL DSCAL( M, SCALOC, C( 1, K ), 1 )
  633:                      CALL DSCAL( M, SCALOC, F( 1, K ), 1 )
  634:   160             CONTINUE
  635:                   DO 170 K = JS, JE
  636:                      CALL DSCAL( IS-1, SCALOC, C( 1, K ), 1 )
  637:                      CALL DSCAL( IS-1, SCALOC, F( 1, K ), 1 )
  638:   170             CONTINUE
  639:                   DO 180 K = JS, JE
  640:                      CALL DSCAL( M-IE, SCALOC, C( IE+1, K ), 1 )
  641:                      CALL DSCAL( M-IE, SCALOC, F( IE+1, K ), 1 )
  642:   180             CONTINUE
  643:                   DO 190 K = JE + 1, N
  644:                      CALL DSCAL( M, SCALOC, C( 1, K ), 1 )
  645:                      CALL DSCAL( M, SCALOC, F( 1, K ), 1 )
  646:   190             CONTINUE
  647:                   SCALE = SCALE*SCALOC
  648:                END IF
  649: *
  650: *              Substitute R(I, J) and L(I, J) into remaining equation.
  651: *
  652:                IF( J.GT.P+2 ) THEN
  653:                   CALL DGEMM( 'N', 'T', MB, JS-1, NB, ONE, C( IS, JS ),
  654:      $                        LDC, B( 1, JS ), LDB, ONE, F( IS, 1 ),
  655:      $                        LDF )
  656:                   CALL DGEMM( 'N', 'T', MB, JS-1, NB, ONE, F( IS, JS ),
  657:      $                        LDF, E( 1, JS ), LDE, ONE, F( IS, 1 ),
  658:      $                        LDF )
  659:                END IF
  660:                IF( I.LT.P ) THEN
  661:                   CALL DGEMM( 'T', 'N', M-IE, NB, MB, -ONE,
  662:      $                        A( IS, IE+1 ), LDA, C( IS, JS ), LDC, ONE,
  663:      $                        C( IE+1, JS ), LDC )
  664:                   CALL DGEMM( 'T', 'N', M-IE, NB, MB, -ONE,
  665:      $                        D( IS, IE+1 ), LDD, F( IS, JS ), LDF, ONE,
  666:      $                        C( IE+1, JS ), LDC )
  667:                END IF
  668:   200       CONTINUE
  669:   210    CONTINUE
  670: *
  671:       END IF
  672: *
  673:       WORK( 1 ) = LWMIN
  674: *
  675:       RETURN
  676: *
  677: *     End of DTGSYL
  678: *
  679:       END

CVSweb interface <joel.bertrand@systella.fr>