LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
clatrs.f
Go to the documentation of this file.
1 *> \brief \b CLATRS solves a triangular system of equations with the scale factor set to prevent overflow.
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download CLATRS + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/clatrs.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/clatrs.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/clatrs.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE CLATRS( UPLO, TRANS, DIAG, NORMIN, N, A, LDA, X, SCALE,
22 * CNORM, INFO )
23 *
24 * .. Scalar Arguments ..
25 * CHARACTER DIAG, NORMIN, TRANS, UPLO
26 * INTEGER INFO, LDA, N
27 * REAL SCALE
28 * ..
29 * .. Array Arguments ..
30 * REAL CNORM( * )
31 * COMPLEX A( LDA, * ), X( * )
32 * ..
33 *
34 *
35 *> \par Purpose:
36 * =============
37 *>
38 *> \verbatim
39 *>
40 *> CLATRS solves one of the triangular systems
41 *>
42 *> A * x = s*b, A**T * x = s*b, or A**H * x = s*b,
43 *>
44 *> with scaling to prevent overflow. Here A is an upper or lower
45 *> triangular matrix, A**T denotes the transpose of A, A**H denotes the
46 *> conjugate transpose of A, x and b are n-element vectors, and s is a
47 *> scaling factor, usually less than or equal to 1, chosen so that the
48 *> components of x will be less than the overflow threshold. If the
49 *> unscaled problem will not cause overflow, the Level 2 BLAS routine
50 *> CTRSV is called. If the matrix A is singular (A(j,j) = 0 for some j),
51 *> then s is set to 0 and a non-trivial solution to A*x = 0 is returned.
52 *> \endverbatim
53 *
54 * Arguments:
55 * ==========
56 *
57 *> \param[in] UPLO
58 *> \verbatim
59 *> UPLO is CHARACTER*1
60 *> Specifies whether the matrix A is upper or lower triangular.
61 *> = 'U': Upper triangular
62 *> = 'L': Lower triangular
63 *> \endverbatim
64 *>
65 *> \param[in] TRANS
66 *> \verbatim
67 *> TRANS is CHARACTER*1
68 *> Specifies the operation applied to A.
69 *> = 'N': Solve A * x = s*b (No transpose)
70 *> = 'T': Solve A**T * x = s*b (Transpose)
71 *> = 'C': Solve A**H * x = s*b (Conjugate transpose)
72 *> \endverbatim
73 *>
74 *> \param[in] DIAG
75 *> \verbatim
76 *> DIAG is CHARACTER*1
77 *> Specifies whether or not the matrix A is unit triangular.
78 *> = 'N': Non-unit triangular
79 *> = 'U': Unit triangular
80 *> \endverbatim
81 *>
82 *> \param[in] NORMIN
83 *> \verbatim
84 *> NORMIN is CHARACTER*1
85 *> Specifies whether CNORM has been set or not.
86 *> = 'Y': CNORM contains the column norms on entry
87 *> = 'N': CNORM is not set on entry. On exit, the norms will
88 *> be computed and stored in CNORM.
89 *> \endverbatim
90 *>
91 *> \param[in] N
92 *> \verbatim
93 *> N is INTEGER
94 *> The order of the matrix A. N >= 0.
95 *> \endverbatim
96 *>
97 *> \param[in] A
98 *> \verbatim
99 *> A is COMPLEX array, dimension (LDA,N)
100 *> The triangular matrix A. If UPLO = 'U', the leading n by n
101 *> upper triangular part of the array A contains the upper
102 *> triangular matrix, and the strictly lower triangular part of
103 *> A is not referenced. If UPLO = 'L', the leading n by n lower
104 *> triangular part of the array A contains the lower triangular
105 *> matrix, and the strictly upper triangular part of A is not
106 *> referenced. If DIAG = 'U', the diagonal elements of A are
107 *> also not referenced and are assumed to be 1.
108 *> \endverbatim
109 *>
110 *> \param[in] LDA
111 *> \verbatim
112 *> LDA is INTEGER
113 *> The leading dimension of the array A. LDA >= max (1,N).
114 *> \endverbatim
115 *>
116 *> \param[in,out] X
117 *> \verbatim
118 *> X is COMPLEX array, dimension (N)
119 *> On entry, the right hand side b of the triangular system.
120 *> On exit, X is overwritten by the solution vector x.
121 *> \endverbatim
122 *>
123 *> \param[out] SCALE
124 *> \verbatim
125 *> SCALE is REAL
126 *> The scaling factor s for the triangular system
127 *> A * x = s*b, A**T * x = s*b, or A**H * x = s*b.
128 *> If SCALE = 0, the matrix A is singular or badly scaled, and
129 *> the vector x is an exact or approximate solution to A*x = 0.
130 *> \endverbatim
131 *>
132 *> \param[in,out] CNORM
133 *> \verbatim
134 *> CNORM is REAL array, dimension (N)
135 *>
136 *> If NORMIN = 'Y', CNORM is an input argument and CNORM(j)
137 *> contains the norm of the off-diagonal part of the j-th column
138 *> of A. If TRANS = 'N', CNORM(j) must be greater than or equal
139 *> to the infinity-norm, and if TRANS = 'T' or 'C', CNORM(j)
140 *> must be greater than or equal to the 1-norm.
141 *>
142 *> If NORMIN = 'N', CNORM is an output argument and CNORM(j)
143 *> returns the 1-norm of the offdiagonal part of the j-th column
144 *> of A.
145 *> \endverbatim
146 *>
147 *> \param[out] INFO
148 *> \verbatim
149 *> INFO is INTEGER
150 *> = 0: successful exit
151 *> < 0: if INFO = -k, the k-th argument had an illegal value
152 *> \endverbatim
153 *
154 * Authors:
155 * ========
156 *
157 *> \author Univ. of Tennessee
158 *> \author Univ. of California Berkeley
159 *> \author Univ. of Colorado Denver
160 *> \author NAG Ltd.
161 *
162 *> \date September 2012
163 *
164 *> \ingroup complexOTHERauxiliary
165 *
166 *> \par Further Details:
167 * =====================
168 *>
169 *> \verbatim
170 *>
171 *> A rough bound on x is computed; if that is less than overflow, CTRSV
172 *> is called, otherwise, specific code is used which checks for possible
173 *> overflow or divide-by-zero at every operation.
174 *>
175 *> A columnwise scheme is used for solving A*x = b. The basic algorithm
176 *> if A is lower triangular is
177 *>
178 *> x[1:n] := b[1:n]
179 *> for j = 1, ..., n
180 *> x(j) := x(j) / A(j,j)
181 *> x[j+1:n] := x[j+1:n] - x(j) * A[j+1:n,j]
182 *> end
183 *>
184 *> Define bounds on the components of x after j iterations of the loop:
185 *> M(j) = bound on x[1:j]
186 *> G(j) = bound on x[j+1:n]
187 *> Initially, let M(0) = 0 and G(0) = max{x(i), i=1,...,n}.
188 *>
189 *> Then for iteration j+1 we have
190 *> M(j+1) <= G(j) / | A(j+1,j+1) |
191 *> G(j+1) <= G(j) + M(j+1) * | A[j+2:n,j+1] |
192 *> <= G(j) ( 1 + CNORM(j+1) / | A(j+1,j+1) | )
193 *>
194 *> where CNORM(j+1) is greater than or equal to the infinity-norm of
195 *> column j+1 of A, not counting the diagonal. Hence
196 *>
197 *> G(j) <= G(0) product ( 1 + CNORM(i) / | A(i,i) | )
198 *> 1<=i<=j
199 *> and
200 *>
201 *> |x(j)| <= ( G(0) / |A(j,j)| ) product ( 1 + CNORM(i) / |A(i,i)| )
202 *> 1<=i< j
203 *>
204 *> Since |x(j)| <= M(j), we use the Level 2 BLAS routine CTRSV if the
205 *> reciprocal of the largest M(j), j=1,..,n, is larger than
206 *> max(underflow, 1/overflow).
207 *>
208 *> The bound on x(j) is also used to determine when a step in the
209 *> columnwise method can be performed without fear of overflow. If
210 *> the computed bound is greater than a large constant, x is scaled to
211 *> prevent overflow, but if the bound overflows, x is set to 0, x(j) to
212 *> 1, and scale to 0, and a non-trivial solution to A*x = 0 is found.
213 *>
214 *> Similarly, a row-wise scheme is used to solve A**T *x = b or
215 *> A**H *x = b. The basic algorithm for A upper triangular is
216 *>
217 *> for j = 1, ..., n
218 *> x(j) := ( b(j) - A[1:j-1,j]' * x[1:j-1] ) / A(j,j)
219 *> end
220 *>
221 *> We simultaneously compute two bounds
222 *> G(j) = bound on ( b(i) - A[1:i-1,i]' * x[1:i-1] ), 1<=i<=j
223 *> M(j) = bound on x(i), 1<=i<=j
224 *>
225 *> The initial values are G(0) = 0, M(0) = max{b(i), i=1,..,n}, and we
226 *> add the constraint G(j) >= G(j-1) and M(j) >= M(j-1) for j >= 1.
227 *> Then the bound on x(j) is
228 *>
229 *> M(j) <= M(j-1) * ( 1 + CNORM(j) ) / | A(j,j) |
230 *>
231 *> <= M(0) * product ( ( 1 + CNORM(i) ) / |A(i,i)| )
232 *> 1<=i<=j
233 *>
234 *> and we can safely call CTRSV if 1/M(n) and 1/G(n) are both greater
235 *> than max(underflow, 1/overflow).
236 *> \endverbatim
237 *>
238 * =====================================================================
239  SUBROUTINE clatrs( UPLO, TRANS, DIAG, NORMIN, N, A, LDA, X, SCALE,
240  $ cnorm, info )
241 *
242 * -- LAPACK auxiliary routine (version 3.4.2) --
243 * -- LAPACK is a software package provided by Univ. of Tennessee, --
244 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
245 * September 2012
246 *
247 * .. Scalar Arguments ..
248  CHARACTER diag, normin, trans, uplo
249  INTEGER info, lda, n
250  REAL scale
251 * ..
252 * .. Array Arguments ..
253  REAL cnorm( * )
254  COMPLEX a( lda, * ), x( * )
255 * ..
256 *
257 * =====================================================================
258 *
259 * .. Parameters ..
260  REAL zero, half, one, two
261  parameter( zero = 0.0e+0, half = 0.5e+0, one = 1.0e+0,
262  $ two = 2.0e+0 )
263 * ..
264 * .. Local Scalars ..
265  LOGICAL notran, nounit, upper
266  INTEGER i, imax, j, jfirst, jinc, jlast
267  REAL bignum, grow, rec, smlnum, tjj, tmax, tscal,
268  $ xbnd, xj, xmax
269  COMPLEX csumj, tjjs, uscal, zdum
270 * ..
271 * .. External Functions ..
272  LOGICAL lsame
273  INTEGER icamax, isamax
274  REAL scasum, slamch
275  COMPLEX cdotc, cdotu, cladiv
276  EXTERNAL lsame, icamax, isamax, scasum, slamch, cdotc,
277  $ cdotu, cladiv
278 * ..
279 * .. External Subroutines ..
280  EXTERNAL caxpy, csscal, ctrsv, slabad, sscal, xerbla
281 * ..
282 * .. Intrinsic Functions ..
283  INTRINSIC abs, aimag, cmplx, conjg, max, min, real
284 * ..
285 * .. Statement Functions ..
286  REAL cabs1, cabs2
287 * ..
288 * .. Statement Function definitions ..
289  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
290  cabs2( zdum ) = abs( REAL( ZDUM ) / 2. ) +
291  $ abs( aimag( zdum ) / 2. )
292 * ..
293 * .. Executable Statements ..
294 *
295  info = 0
296  upper = lsame( uplo, 'U' )
297  notran = lsame( trans, 'N' )
298  nounit = lsame( diag, 'N' )
299 *
300 * Test the input parameters.
301 *
302  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
303  info = -1
304  ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
305  $ lsame( trans, 'C' ) ) THEN
306  info = -2
307  ELSE IF( .NOT.nounit .AND. .NOT.lsame( diag, 'U' ) ) THEN
308  info = -3
309  ELSE IF( .NOT.lsame( normin, 'Y' ) .AND. .NOT.
310  $ lsame( normin, 'N' ) ) THEN
311  info = -4
312  ELSE IF( n.LT.0 ) THEN
313  info = -5
314  ELSE IF( lda.LT.max( 1, n ) ) THEN
315  info = -7
316  END IF
317  IF( info.NE.0 ) THEN
318  CALL xerbla( 'CLATRS', -info )
319  return
320  END IF
321 *
322 * Quick return if possible
323 *
324  IF( n.EQ.0 )
325  $ return
326 *
327 * Determine machine dependent parameters to control overflow.
328 *
329  smlnum = slamch( 'Safe minimum' )
330  bignum = one / smlnum
331  CALL slabad( smlnum, bignum )
332  smlnum = smlnum / slamch( 'Precision' )
333  bignum = one / smlnum
334  scale = one
335 *
336  IF( lsame( normin, 'N' ) ) THEN
337 *
338 * Compute the 1-norm of each column, not including the diagonal.
339 *
340  IF( upper ) THEN
341 *
342 * A is upper triangular.
343 *
344  DO 10 j = 1, n
345  cnorm( j ) = scasum( j-1, a( 1, j ), 1 )
346  10 continue
347  ELSE
348 *
349 * A is lower triangular.
350 *
351  DO 20 j = 1, n - 1
352  cnorm( j ) = scasum( n-j, a( j+1, j ), 1 )
353  20 continue
354  cnorm( n ) = zero
355  END IF
356  END IF
357 *
358 * Scale the column norms by TSCAL if the maximum element in CNORM is
359 * greater than BIGNUM/2.
360 *
361  imax = isamax( n, cnorm, 1 )
362  tmax = cnorm( imax )
363  IF( tmax.LE.bignum*half ) THEN
364  tscal = one
365  ELSE
366  tscal = half / ( smlnum*tmax )
367  CALL sscal( n, tscal, cnorm, 1 )
368  END IF
369 *
370 * Compute a bound on the computed solution vector to see if the
371 * Level 2 BLAS routine CTRSV can be used.
372 *
373  xmax = zero
374  DO 30 j = 1, n
375  xmax = max( xmax, cabs2( x( j ) ) )
376  30 continue
377  xbnd = xmax
378 *
379  IF( notran ) THEN
380 *
381 * Compute the growth in A * x = b.
382 *
383  IF( upper ) THEN
384  jfirst = n
385  jlast = 1
386  jinc = -1
387  ELSE
388  jfirst = 1
389  jlast = n
390  jinc = 1
391  END IF
392 *
393  IF( tscal.NE.one ) THEN
394  grow = zero
395  go to 60
396  END IF
397 *
398  IF( nounit ) THEN
399 *
400 * A is non-unit triangular.
401 *
402 * Compute GROW = 1/G(j) and XBND = 1/M(j).
403 * Initially, G(0) = max{x(i), i=1,...,n}.
404 *
405  grow = half / max( xbnd, smlnum )
406  xbnd = grow
407  DO 40 j = jfirst, jlast, jinc
408 *
409 * Exit the loop if the growth factor is too small.
410 *
411  IF( grow.LE.smlnum )
412  $ go to 60
413 *
414  tjjs = a( j, j )
415  tjj = cabs1( tjjs )
416 *
417  IF( tjj.GE.smlnum ) THEN
418 *
419 * M(j) = G(j-1) / abs(A(j,j))
420 *
421  xbnd = min( xbnd, min( one, tjj )*grow )
422  ELSE
423 *
424 * M(j) could overflow, set XBND to 0.
425 *
426  xbnd = zero
427  END IF
428 *
429  IF( tjj+cnorm( j ).GE.smlnum ) THEN
430 *
431 * G(j) = G(j-1)*( 1 + CNORM(j) / abs(A(j,j)) )
432 *
433  grow = grow*( tjj / ( tjj+cnorm( j ) ) )
434  ELSE
435 *
436 * G(j) could overflow, set GROW to 0.
437 *
438  grow = zero
439  END IF
440  40 continue
441  grow = xbnd
442  ELSE
443 *
444 * A is unit triangular.
445 *
446 * Compute GROW = 1/G(j), where G(0) = max{x(i), i=1,...,n}.
447 *
448  grow = min( one, half / max( xbnd, smlnum ) )
449  DO 50 j = jfirst, jlast, jinc
450 *
451 * Exit the loop if the growth factor is too small.
452 *
453  IF( grow.LE.smlnum )
454  $ go to 60
455 *
456 * G(j) = G(j-1)*( 1 + CNORM(j) )
457 *
458  grow = grow*( one / ( one+cnorm( j ) ) )
459  50 continue
460  END IF
461  60 continue
462 *
463  ELSE
464 *
465 * Compute the growth in A**T * x = b or A**H * x = b.
466 *
467  IF( upper ) THEN
468  jfirst = 1
469  jlast = n
470  jinc = 1
471  ELSE
472  jfirst = n
473  jlast = 1
474  jinc = -1
475  END IF
476 *
477  IF( tscal.NE.one ) THEN
478  grow = zero
479  go to 90
480  END IF
481 *
482  IF( nounit ) THEN
483 *
484 * A is non-unit triangular.
485 *
486 * Compute GROW = 1/G(j) and XBND = 1/M(j).
487 * Initially, M(0) = max{x(i), i=1,...,n}.
488 *
489  grow = half / max( xbnd, smlnum )
490  xbnd = grow
491  DO 70 j = jfirst, jlast, jinc
492 *
493 * Exit the loop if the growth factor is too small.
494 *
495  IF( grow.LE.smlnum )
496  $ go to 90
497 *
498 * G(j) = max( G(j-1), M(j-1)*( 1 + CNORM(j) ) )
499 *
500  xj = one + cnorm( j )
501  grow = min( grow, xbnd / xj )
502 *
503  tjjs = a( j, j )
504  tjj = cabs1( tjjs )
505 *
506  IF( tjj.GE.smlnum ) THEN
507 *
508 * M(j) = M(j-1)*( 1 + CNORM(j) ) / abs(A(j,j))
509 *
510  IF( xj.GT.tjj )
511  $ xbnd = xbnd*( tjj / xj )
512  ELSE
513 *
514 * M(j) could overflow, set XBND to 0.
515 *
516  xbnd = zero
517  END IF
518  70 continue
519  grow = min( grow, xbnd )
520  ELSE
521 *
522 * A is unit triangular.
523 *
524 * Compute GROW = 1/G(j), where G(0) = max{x(i), i=1,...,n}.
525 *
526  grow = min( one, half / max( xbnd, smlnum ) )
527  DO 80 j = jfirst, jlast, jinc
528 *
529 * Exit the loop if the growth factor is too small.
530 *
531  IF( grow.LE.smlnum )
532  $ go to 90
533 *
534 * G(j) = ( 1 + CNORM(j) )*G(j-1)
535 *
536  xj = one + cnorm( j )
537  grow = grow / xj
538  80 continue
539  END IF
540  90 continue
541  END IF
542 *
543  IF( ( grow*tscal ).GT.smlnum ) THEN
544 *
545 * Use the Level 2 BLAS solve if the reciprocal of the bound on
546 * elements of X is not too small.
547 *
548  CALL ctrsv( uplo, trans, diag, n, a, lda, x, 1 )
549  ELSE
550 *
551 * Use a Level 1 BLAS solve, scaling intermediate results.
552 *
553  IF( xmax.GT.bignum*half ) THEN
554 *
555 * Scale X so that its components are less than or equal to
556 * BIGNUM in absolute value.
557 *
558  scale = ( bignum*half ) / xmax
559  CALL csscal( n, scale, x, 1 )
560  xmax = bignum
561  ELSE
562  xmax = xmax*two
563  END IF
564 *
565  IF( notran ) THEN
566 *
567 * Solve A * x = b
568 *
569  DO 110 j = jfirst, jlast, jinc
570 *
571 * Compute x(j) = b(j) / A(j,j), scaling x if necessary.
572 *
573  xj = cabs1( x( j ) )
574  IF( nounit ) THEN
575  tjjs = a( j, j )*tscal
576  ELSE
577  tjjs = tscal
578  IF( tscal.EQ.one )
579  $ go to 105
580  END IF
581  tjj = cabs1( tjjs )
582  IF( tjj.GT.smlnum ) THEN
583 *
584 * abs(A(j,j)) > SMLNUM:
585 *
586  IF( tjj.LT.one ) THEN
587  IF( xj.GT.tjj*bignum ) THEN
588 *
589 * Scale x by 1/b(j).
590 *
591  rec = one / xj
592  CALL csscal( n, rec, x, 1 )
593  scale = scale*rec
594  xmax = xmax*rec
595  END IF
596  END IF
597  x( j ) = cladiv( x( j ), tjjs )
598  xj = cabs1( x( j ) )
599  ELSE IF( tjj.GT.zero ) THEN
600 *
601 * 0 < abs(A(j,j)) <= SMLNUM:
602 *
603  IF( xj.GT.tjj*bignum ) THEN
604 *
605 * Scale x by (1/abs(x(j)))*abs(A(j,j))*BIGNUM
606 * to avoid overflow when dividing by A(j,j).
607 *
608  rec = ( tjj*bignum ) / xj
609  IF( cnorm( j ).GT.one ) THEN
610 *
611 * Scale by 1/CNORM(j) to avoid overflow when
612 * multiplying x(j) times column j.
613 *
614  rec = rec / cnorm( j )
615  END IF
616  CALL csscal( n, rec, x, 1 )
617  scale = scale*rec
618  xmax = xmax*rec
619  END IF
620  x( j ) = cladiv( x( j ), tjjs )
621  xj = cabs1( x( j ) )
622  ELSE
623 *
624 * A(j,j) = 0: Set x(1:n) = 0, x(j) = 1, and
625 * scale = 0, and compute a solution to A*x = 0.
626 *
627  DO 100 i = 1, n
628  x( i ) = zero
629  100 continue
630  x( j ) = one
631  xj = one
632  scale = zero
633  xmax = zero
634  END IF
635  105 continue
636 *
637 * Scale x if necessary to avoid overflow when adding a
638 * multiple of column j of A.
639 *
640  IF( xj.GT.one ) THEN
641  rec = one / xj
642  IF( cnorm( j ).GT.( bignum-xmax )*rec ) THEN
643 *
644 * Scale x by 1/(2*abs(x(j))).
645 *
646  rec = rec*half
647  CALL csscal( n, rec, x, 1 )
648  scale = scale*rec
649  END IF
650  ELSE IF( xj*cnorm( j ).GT.( bignum-xmax ) ) THEN
651 *
652 * Scale x by 1/2.
653 *
654  CALL csscal( n, half, x, 1 )
655  scale = scale*half
656  END IF
657 *
658  IF( upper ) THEN
659  IF( j.GT.1 ) THEN
660 *
661 * Compute the update
662 * x(1:j-1) := x(1:j-1) - x(j) * A(1:j-1,j)
663 *
664  CALL caxpy( j-1, -x( j )*tscal, a( 1, j ), 1, x,
665  $ 1 )
666  i = icamax( j-1, x, 1 )
667  xmax = cabs1( x( i ) )
668  END IF
669  ELSE
670  IF( j.LT.n ) THEN
671 *
672 * Compute the update
673 * x(j+1:n) := x(j+1:n) - x(j) * A(j+1:n,j)
674 *
675  CALL caxpy( n-j, -x( j )*tscal, a( j+1, j ), 1,
676  $ x( j+1 ), 1 )
677  i = j + icamax( n-j, x( j+1 ), 1 )
678  xmax = cabs1( x( i ) )
679  END IF
680  END IF
681  110 continue
682 *
683  ELSE IF( lsame( trans, 'T' ) ) THEN
684 *
685 * Solve A**T * x = b
686 *
687  DO 150 j = jfirst, jlast, jinc
688 *
689 * Compute x(j) = b(j) - sum A(k,j)*x(k).
690 * k<>j
691 *
692  xj = cabs1( x( j ) )
693  uscal = tscal
694  rec = one / max( xmax, one )
695  IF( cnorm( j ).GT.( bignum-xj )*rec ) THEN
696 *
697 * If x(j) could overflow, scale x by 1/(2*XMAX).
698 *
699  rec = rec*half
700  IF( nounit ) THEN
701  tjjs = a( j, j )*tscal
702  ELSE
703  tjjs = tscal
704  END IF
705  tjj = cabs1( tjjs )
706  IF( tjj.GT.one ) THEN
707 *
708 * Divide by A(j,j) when scaling x if A(j,j) > 1.
709 *
710  rec = min( one, rec*tjj )
711  uscal = cladiv( uscal, tjjs )
712  END IF
713  IF( rec.LT.one ) THEN
714  CALL csscal( n, rec, x, 1 )
715  scale = scale*rec
716  xmax = xmax*rec
717  END IF
718  END IF
719 *
720  csumj = zero
721  IF( uscal.EQ.cmplx( one ) ) THEN
722 *
723 * If the scaling needed for A in the dot product is 1,
724 * call CDOTU to perform the dot product.
725 *
726  IF( upper ) THEN
727  csumj = cdotu( j-1, a( 1, j ), 1, x, 1 )
728  ELSE IF( j.LT.n ) THEN
729  csumj = cdotu( n-j, a( j+1, j ), 1, x( j+1 ), 1 )
730  END IF
731  ELSE
732 *
733 * Otherwise, use in-line code for the dot product.
734 *
735  IF( upper ) THEN
736  DO 120 i = 1, j - 1
737  csumj = csumj + ( a( i, j )*uscal )*x( i )
738  120 continue
739  ELSE IF( j.LT.n ) THEN
740  DO 130 i = j + 1, n
741  csumj = csumj + ( a( i, j )*uscal )*x( i )
742  130 continue
743  END IF
744  END IF
745 *
746  IF( uscal.EQ.cmplx( tscal ) ) THEN
747 *
748 * Compute x(j) := ( x(j) - CSUMJ ) / A(j,j) if 1/A(j,j)
749 * was not used to scale the dotproduct.
750 *
751  x( j ) = x( j ) - csumj
752  xj = cabs1( x( j ) )
753  IF( nounit ) THEN
754  tjjs = a( j, j )*tscal
755  ELSE
756  tjjs = tscal
757  IF( tscal.EQ.one )
758  $ go to 145
759  END IF
760 *
761 * Compute x(j) = x(j) / A(j,j), scaling if necessary.
762 *
763  tjj = cabs1( tjjs )
764  IF( tjj.GT.smlnum ) THEN
765 *
766 * abs(A(j,j)) > SMLNUM:
767 *
768  IF( tjj.LT.one ) THEN
769  IF( xj.GT.tjj*bignum ) THEN
770 *
771 * Scale X by 1/abs(x(j)).
772 *
773  rec = one / xj
774  CALL csscal( n, rec, x, 1 )
775  scale = scale*rec
776  xmax = xmax*rec
777  END IF
778  END IF
779  x( j ) = cladiv( x( j ), tjjs )
780  ELSE IF( tjj.GT.zero ) THEN
781 *
782 * 0 < abs(A(j,j)) <= SMLNUM:
783 *
784  IF( xj.GT.tjj*bignum ) THEN
785 *
786 * Scale x by (1/abs(x(j)))*abs(A(j,j))*BIGNUM.
787 *
788  rec = ( tjj*bignum ) / xj
789  CALL csscal( n, rec, x, 1 )
790  scale = scale*rec
791  xmax = xmax*rec
792  END IF
793  x( j ) = cladiv( x( j ), tjjs )
794  ELSE
795 *
796 * A(j,j) = 0: Set x(1:n) = 0, x(j) = 1, and
797 * scale = 0 and compute a solution to A**T *x = 0.
798 *
799  DO 140 i = 1, n
800  x( i ) = zero
801  140 continue
802  x( j ) = one
803  scale = zero
804  xmax = zero
805  END IF
806  145 continue
807  ELSE
808 *
809 * Compute x(j) := x(j) / A(j,j) - CSUMJ if the dot
810 * product has already been divided by 1/A(j,j).
811 *
812  x( j ) = cladiv( x( j ), tjjs ) - csumj
813  END IF
814  xmax = max( xmax, cabs1( x( j ) ) )
815  150 continue
816 *
817  ELSE
818 *
819 * Solve A**H * x = b
820 *
821  DO 190 j = jfirst, jlast, jinc
822 *
823 * Compute x(j) = b(j) - sum A(k,j)*x(k).
824 * k<>j
825 *
826  xj = cabs1( x( j ) )
827  uscal = tscal
828  rec = one / max( xmax, one )
829  IF( cnorm( j ).GT.( bignum-xj )*rec ) THEN
830 *
831 * If x(j) could overflow, scale x by 1/(2*XMAX).
832 *
833  rec = rec*half
834  IF( nounit ) THEN
835  tjjs = conjg( a( j, j ) )*tscal
836  ELSE
837  tjjs = tscal
838  END IF
839  tjj = cabs1( tjjs )
840  IF( tjj.GT.one ) THEN
841 *
842 * Divide by A(j,j) when scaling x if A(j,j) > 1.
843 *
844  rec = min( one, rec*tjj )
845  uscal = cladiv( uscal, tjjs )
846  END IF
847  IF( rec.LT.one ) THEN
848  CALL csscal( n, rec, x, 1 )
849  scale = scale*rec
850  xmax = xmax*rec
851  END IF
852  END IF
853 *
854  csumj = zero
855  IF( uscal.EQ.cmplx( one ) ) THEN
856 *
857 * If the scaling needed for A in the dot product is 1,
858 * call CDOTC to perform the dot product.
859 *
860  IF( upper ) THEN
861  csumj = cdotc( j-1, a( 1, j ), 1, x, 1 )
862  ELSE IF( j.LT.n ) THEN
863  csumj = cdotc( n-j, a( j+1, j ), 1, x( j+1 ), 1 )
864  END IF
865  ELSE
866 *
867 * Otherwise, use in-line code for the dot product.
868 *
869  IF( upper ) THEN
870  DO 160 i = 1, j - 1
871  csumj = csumj + ( conjg( a( i, j ) )*uscal )*
872  $ x( i )
873  160 continue
874  ELSE IF( j.LT.n ) THEN
875  DO 170 i = j + 1, n
876  csumj = csumj + ( conjg( a( i, j ) )*uscal )*
877  $ x( i )
878  170 continue
879  END IF
880  END IF
881 *
882  IF( uscal.EQ.cmplx( tscal ) ) THEN
883 *
884 * Compute x(j) := ( x(j) - CSUMJ ) / A(j,j) if 1/A(j,j)
885 * was not used to scale the dotproduct.
886 *
887  x( j ) = x( j ) - csumj
888  xj = cabs1( x( j ) )
889  IF( nounit ) THEN
890  tjjs = conjg( a( j, j ) )*tscal
891  ELSE
892  tjjs = tscal
893  IF( tscal.EQ.one )
894  $ go to 185
895  END IF
896 *
897 * Compute x(j) = x(j) / A(j,j), scaling if necessary.
898 *
899  tjj = cabs1( tjjs )
900  IF( tjj.GT.smlnum ) THEN
901 *
902 * abs(A(j,j)) > SMLNUM:
903 *
904  IF( tjj.LT.one ) THEN
905  IF( xj.GT.tjj*bignum ) THEN
906 *
907 * Scale X by 1/abs(x(j)).
908 *
909  rec = one / xj
910  CALL csscal( n, rec, x, 1 )
911  scale = scale*rec
912  xmax = xmax*rec
913  END IF
914  END IF
915  x( j ) = cladiv( x( j ), tjjs )
916  ELSE IF( tjj.GT.zero ) THEN
917 *
918 * 0 < abs(A(j,j)) <= SMLNUM:
919 *
920  IF( xj.GT.tjj*bignum ) THEN
921 *
922 * Scale x by (1/abs(x(j)))*abs(A(j,j))*BIGNUM.
923 *
924  rec = ( tjj*bignum ) / xj
925  CALL csscal( n, rec, x, 1 )
926  scale = scale*rec
927  xmax = xmax*rec
928  END IF
929  x( j ) = cladiv( x( j ), tjjs )
930  ELSE
931 *
932 * A(j,j) = 0: Set x(1:n) = 0, x(j) = 1, and
933 * scale = 0 and compute a solution to A**H *x = 0.
934 *
935  DO 180 i = 1, n
936  x( i ) = zero
937  180 continue
938  x( j ) = one
939  scale = zero
940  xmax = zero
941  END IF
942  185 continue
943  ELSE
944 *
945 * Compute x(j) := x(j) / A(j,j) - CSUMJ if the dot
946 * product has already been divided by 1/A(j,j).
947 *
948  x( j ) = cladiv( x( j ), tjjs ) - csumj
949  END IF
950  xmax = max( xmax, cabs1( x( j ) ) )
951  190 continue
952  END IF
953  scale = scale / tscal
954  END IF
955 *
956 * Scale the column norms by 1/TSCAL for return.
957 *
958  IF( tscal.NE.one ) THEN
959  CALL sscal( n, one / tscal, cnorm, 1 )
960  END IF
961 *
962  return
963 *
964 * End of CLATRS
965 *
966  END