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

◆ ztrsyl()

subroutine ztrsyl ( character  trana,
character  tranb,
integer  isgn,
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,
double precision  scale,
integer  info 
)

ZTRSYL

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

Purpose:
 ZTRSYL solves the complex Sylvester matrix equation:

    op(A)*X + X*op(B) = scale*C or
    op(A)*X - X*op(B) = scale*C,

 where op(A) = A or A**H, and A and B are both upper triangular. A is
 M-by-M and B is N-by-N; the right hand side C and the solution X are
 M-by-N; and scale is an output scale factor, set <= 1 to avoid
 overflow in X.
Parameters
[in]TRANA
          TRANA is CHARACTER*1
          Specifies the option op(A):
          = 'N': op(A) = A    (No transpose)
          = 'C': op(A) = A**H (Conjugate transpose)
[in]TRANB
          TRANB is CHARACTER*1
          Specifies the option op(B):
          = 'N': op(B) = B    (No transpose)
          = 'C': op(B) = B**H (Conjugate transpose)
[in]ISGN
          ISGN is INTEGER
          Specifies the sign in the equation:
          = +1: solve op(A)*X + X*op(B) = scale*C
          = -1: solve op(A)*X - X*op(B) = scale*C
[in]M
          M is INTEGER
          The order of the matrix A, and the number of rows in the
          matrices X and C. M >= 0.
[in]N
          N is INTEGER
          The order of the matrix B, and the number of columns in the
          matrices X and C. N >= 0.
[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, the M-by-N right hand side matrix C.
          On exit, C is overwritten by the solution matrix X.
[in]LDC
          LDC is INTEGER
          The leading dimension of the array C. LDC >= max(1,M)
[out]SCALE
          SCALE is DOUBLE PRECISION
          The scale factor, scale, set <= 1 to avoid overflow in X.
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -i, the i-th argument had an illegal value
          = 1: A and B have common or very close eigenvalues; perturbed
               values were used to solve the equation (but the matrices
               A and B are unchanged).
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 155 of file ztrsyl.f.

157*
158* -- LAPACK computational routine --
159* -- LAPACK is a software package provided by Univ. of Tennessee, --
160* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
161*
162* .. Scalar Arguments ..
163 CHARACTER TRANA, TRANB
164 INTEGER INFO, ISGN, LDA, LDB, LDC, M, N
165 DOUBLE PRECISION SCALE
166* ..
167* .. Array Arguments ..
168 COMPLEX*16 A( LDA, * ), B( LDB, * ), C( LDC, * )
169* ..
170*
171* =====================================================================
172*
173* .. Parameters ..
174 DOUBLE PRECISION ONE
175 parameter( one = 1.0d+0 )
176* ..
177* .. Local Scalars ..
178 LOGICAL NOTRNA, NOTRNB
179 INTEGER J, K, L
180 DOUBLE PRECISION BIGNUM, DA11, DB, EPS, SCALOC, SGN, SMIN,
181 $ SMLNUM
182 COMPLEX*16 A11, SUML, SUMR, VEC, X11
183* ..
184* .. Local Arrays ..
185 DOUBLE PRECISION DUM( 1 )
186* ..
187* .. External Functions ..
188 LOGICAL LSAME
189 DOUBLE PRECISION DLAMCH, ZLANGE
190 COMPLEX*16 ZDOTC, ZDOTU, ZLADIV
191 EXTERNAL lsame, dlamch, zlange, zdotc, zdotu, zladiv
192* ..
193* .. External Subroutines ..
194 EXTERNAL xerbla, zdscal
195* ..
196* .. Intrinsic Functions ..
197 INTRINSIC abs, dble, dcmplx, dconjg, dimag, max, min
198* ..
199* .. Executable Statements ..
200*
201* Decode and Test input parameters
202*
203 notrna = lsame( trana, 'N' )
204 notrnb = lsame( tranb, 'N' )
205*
206 info = 0
207 IF( .NOT.notrna .AND. .NOT.lsame( trana, 'C' ) ) THEN
208 info = -1
209 ELSE IF( .NOT.notrnb .AND. .NOT.lsame( tranb, 'C' ) ) THEN
210 info = -2
211 ELSE IF( isgn.NE.1 .AND. isgn.NE.-1 ) THEN
212 info = -3
213 ELSE IF( m.LT.0 ) THEN
214 info = -4
215 ELSE IF( n.LT.0 ) THEN
216 info = -5
217 ELSE IF( lda.LT.max( 1, m ) ) THEN
218 info = -7
219 ELSE IF( ldb.LT.max( 1, n ) ) THEN
220 info = -9
221 ELSE IF( ldc.LT.max( 1, m ) ) THEN
222 info = -11
223 END IF
224 IF( info.NE.0 ) THEN
225 CALL xerbla( 'ZTRSYL', -info )
226 RETURN
227 END IF
228*
229* Quick return if possible
230*
231 scale = one
232 IF( m.EQ.0 .OR. n.EQ.0 )
233 $ RETURN
234*
235* Set constants to control overflow
236*
237 eps = dlamch( 'P' )
238 smlnum = dlamch( 'S' )
239 bignum = one / smlnum
240 smlnum = smlnum*dble( m*n ) / eps
241 bignum = one / smlnum
242 smin = max( smlnum, eps*zlange( 'M', m, m, a, lda, dum ),
243 $ eps*zlange( 'M', n, n, b, ldb, dum ) )
244 sgn = isgn
245*
246 IF( notrna .AND. notrnb ) THEN
247*
248* Solve A*X + ISGN*X*B = scale*C.
249*
250* The (K,L)th block of X is determined starting from
251* bottom-left corner column by column by
252*
253* A(K,K)*X(K,L) + ISGN*X(K,L)*B(L,L) = C(K,L) - R(K,L)
254*
255* Where
256* M L-1
257* R(K,L) = SUM [A(K,I)*X(I,L)] +ISGN*SUM [X(K,J)*B(J,L)].
258* I=K+1 J=1
259*
260 DO 30 l = 1, n
261 DO 20 k = m, 1, -1
262*
263 suml = zdotu( m-k, a( k, min( k+1, m ) ), lda,
264 $ c( min( k+1, m ), l ), 1 )
265 sumr = zdotu( l-1, c( k, 1 ), ldc, b( 1, l ), 1 )
266 vec = c( k, l ) - ( suml+sgn*sumr )
267*
268 scaloc = one
269 a11 = a( k, k ) + sgn*b( l, l )
270 da11 = abs( dble( a11 ) ) + abs( dimag( a11 ) )
271 IF( da11.LE.smin ) THEN
272 a11 = smin
273 da11 = smin
274 info = 1
275 END IF
276 db = abs( dble( vec ) ) + abs( dimag( vec ) )
277 IF( da11.LT.one .AND. db.GT.one ) THEN
278 IF( db.GT.bignum*da11 )
279 $ scaloc = one / db
280 END IF
281 x11 = zladiv( vec*dcmplx( scaloc ), a11 )
282*
283 IF( scaloc.NE.one ) THEN
284 DO 10 j = 1, n
285 CALL zdscal( m, scaloc, c( 1, j ), 1 )
286 10 CONTINUE
287 scale = scale*scaloc
288 END IF
289 c( k, l ) = x11
290*
291 20 CONTINUE
292 30 CONTINUE
293*
294 ELSE IF( .NOT.notrna .AND. notrnb ) THEN
295*
296* Solve A**H *X + ISGN*X*B = scale*C.
297*
298* The (K,L)th block of X is determined starting from
299* upper-left corner column by column by
300*
301* A**H(K,K)*X(K,L) + ISGN*X(K,L)*B(L,L) = C(K,L) - R(K,L)
302*
303* Where
304* K-1 L-1
305* R(K,L) = SUM [A**H(I,K)*X(I,L)] + ISGN*SUM [X(K,J)*B(J,L)]
306* I=1 J=1
307*
308 DO 60 l = 1, n
309 DO 50 k = 1, m
310*
311 suml = zdotc( k-1, a( 1, k ), 1, c( 1, l ), 1 )
312 sumr = zdotu( l-1, c( k, 1 ), ldc, b( 1, l ), 1 )
313 vec = c( k, l ) - ( suml+sgn*sumr )
314*
315 scaloc = one
316 a11 = dconjg( a( k, k ) ) + sgn*b( l, l )
317 da11 = abs( dble( a11 ) ) + abs( dimag( a11 ) )
318 IF( da11.LE.smin ) THEN
319 a11 = smin
320 da11 = smin
321 info = 1
322 END IF
323 db = abs( dble( vec ) ) + abs( dimag( vec ) )
324 IF( da11.LT.one .AND. db.GT.one ) THEN
325 IF( db.GT.bignum*da11 )
326 $ scaloc = one / db
327 END IF
328*
329 x11 = zladiv( vec*dcmplx( scaloc ), a11 )
330*
331 IF( scaloc.NE.one ) THEN
332 DO 40 j = 1, n
333 CALL zdscal( m, scaloc, c( 1, j ), 1 )
334 40 CONTINUE
335 scale = scale*scaloc
336 END IF
337 c( k, l ) = x11
338*
339 50 CONTINUE
340 60 CONTINUE
341*
342 ELSE IF( .NOT.notrna .AND. .NOT.notrnb ) THEN
343*
344* Solve A**H*X + ISGN*X*B**H = C.
345*
346* The (K,L)th block of X is determined starting from
347* upper-right corner column by column by
348*
349* A**H(K,K)*X(K,L) + ISGN*X(K,L)*B**H(L,L) = C(K,L) - R(K,L)
350*
351* Where
352* K-1
353* R(K,L) = SUM [A**H(I,K)*X(I,L)] +
354* I=1
355* N
356* ISGN*SUM [X(K,J)*B**H(L,J)].
357* J=L+1
358*
359 DO 90 l = n, 1, -1
360 DO 80 k = 1, m
361*
362 suml = zdotc( k-1, a( 1, k ), 1, c( 1, l ), 1 )
363 sumr = zdotc( n-l, c( k, min( l+1, n ) ), ldc,
364 $ b( l, min( l+1, n ) ), ldb )
365 vec = c( k, l ) - ( suml+sgn*dconjg( sumr ) )
366*
367 scaloc = one
368 a11 = dconjg( a( k, k )+sgn*b( l, l ) )
369 da11 = abs( dble( a11 ) ) + abs( dimag( a11 ) )
370 IF( da11.LE.smin ) THEN
371 a11 = smin
372 da11 = smin
373 info = 1
374 END IF
375 db = abs( dble( vec ) ) + abs( dimag( vec ) )
376 IF( da11.LT.one .AND. db.GT.one ) THEN
377 IF( db.GT.bignum*da11 )
378 $ scaloc = one / db
379 END IF
380*
381 x11 = zladiv( vec*dcmplx( scaloc ), a11 )
382*
383 IF( scaloc.NE.one ) THEN
384 DO 70 j = 1, n
385 CALL zdscal( m, scaloc, c( 1, j ), 1 )
386 70 CONTINUE
387 scale = scale*scaloc
388 END IF
389 c( k, l ) = x11
390*
391 80 CONTINUE
392 90 CONTINUE
393*
394 ELSE IF( notrna .AND. .NOT.notrnb ) THEN
395*
396* Solve A*X + ISGN*X*B**H = C.
397*
398* The (K,L)th block of X is determined starting from
399* bottom-left corner column by column by
400*
401* A(K,K)*X(K,L) + ISGN*X(K,L)*B**H(L,L) = C(K,L) - R(K,L)
402*
403* Where
404* M N
405* R(K,L) = SUM [A(K,I)*X(I,L)] + ISGN*SUM [X(K,J)*B**H(L,J)]
406* I=K+1 J=L+1
407*
408 DO 120 l = n, 1, -1
409 DO 110 k = m, 1, -1
410*
411 suml = zdotu( m-k, a( k, min( k+1, m ) ), lda,
412 $ c( min( k+1, m ), l ), 1 )
413 sumr = zdotc( n-l, c( k, min( l+1, n ) ), ldc,
414 $ b( l, min( l+1, n ) ), ldb )
415 vec = c( k, l ) - ( suml+sgn*dconjg( sumr ) )
416*
417 scaloc = one
418 a11 = a( k, k ) + sgn*dconjg( b( l, l ) )
419 da11 = abs( dble( a11 ) ) + abs( dimag( a11 ) )
420 IF( da11.LE.smin ) THEN
421 a11 = smin
422 da11 = smin
423 info = 1
424 END IF
425 db = abs( dble( vec ) ) + abs( dimag( vec ) )
426 IF( da11.LT.one .AND. db.GT.one ) THEN
427 IF( db.GT.bignum*da11 )
428 $ scaloc = one / db
429 END IF
430*
431 x11 = zladiv( vec*dcmplx( scaloc ), a11 )
432*
433 IF( scaloc.NE.one ) THEN
434 DO 100 j = 1, n
435 CALL zdscal( m, scaloc, c( 1, j ), 1 )
436 100 CONTINUE
437 scale = scale*scaloc
438 END IF
439 c( k, l ) = x11
440*
441 110 CONTINUE
442 120 CONTINUE
443*
444 END IF
445*
446 RETURN
447*
448* End of ZTRSYL
449*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
complex *16 function zdotc(n, zx, incx, zy, incy)
ZDOTC
Definition zdotc.f:83
complex *16 function zdotu(n, zx, incx, zy, incy)
ZDOTU
Definition zdotu.f:83
complex *16 function zladiv(x, y)
ZLADIV performs complex division in real arithmetic, avoiding unnecessary overflow.
Definition zladiv.f:64
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:115
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
Here is the call graph for this function:
Here is the caller graph for this function: