1*> \brief \b DTPT03
2*
3*  =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6*            http://www.netlib.org/lapack/explore-html/
7*
8*  Definition:
9*  ===========
10*
11*       SUBROUTINE DTPT03( UPLO, TRANS, DIAG, N, NRHS, AP, SCALE, CNORM,
12*                          TSCAL, X, LDX, B, LDB, WORK, RESID )
13*
14*       .. Scalar Arguments ..
15*       CHARACTER          DIAG, TRANS, UPLO
16*       INTEGER            LDB, LDX, N, NRHS
17*       DOUBLE PRECISION   RESID, SCALE, TSCAL
18*       ..
19*       .. Array Arguments ..
20*       DOUBLE PRECISION   AP( * ), B( LDB, * ), CNORM( * ), WORK( * ),
21*      $                   X( LDX, * )
22*       ..
23*
24*
25*> \par Purpose:
26*  =============
27*>
28*> \verbatim
29*>
30*> DTPT03 computes the residual for the solution to a scaled triangular
31*> system of equations A*x = s*b  or  A'*x = s*b  when the triangular
32*> matrix A is stored in packed format.  Here A' is the transpose of A,
33*> s is a scalar, and x and b are N by NRHS matrices.  The test ratio is
34*> the maximum over the number of right hand sides of
35*>    norm(s*b - op(A)*x) / ( norm(op(A)) * norm(x) * EPS ),
36*> where op(A) denotes A or A' and EPS is the machine epsilon.
37*> \endverbatim
38*
39*  Arguments:
40*  ==========
41*
42*> \param[in] UPLO
43*> \verbatim
44*>          UPLO is CHARACTER*1
45*>          Specifies whether the matrix A is upper or lower triangular.
46*>          = 'U':  Upper triangular
47*>          = 'L':  Lower triangular
48*> \endverbatim
49*>
50*> \param[in] TRANS
51*> \verbatim
52*>          TRANS is CHARACTER*1
53*>          Specifies the operation applied to A.
54*>          = 'N':  A *x = s*b  (No transpose)
55*>          = 'T':  A'*x = s*b  (Transpose)
56*>          = 'C':  A'*x = s*b  (Conjugate transpose = Transpose)
57*> \endverbatim
58*>
59*> \param[in] DIAG
60*> \verbatim
61*>          DIAG is CHARACTER*1
62*>          Specifies whether or not the matrix A is unit triangular.
63*>          = 'N':  Non-unit triangular
64*>          = 'U':  Unit triangular
65*> \endverbatim
66*>
67*> \param[in] N
68*> \verbatim
69*>          N is INTEGER
70*>          The order of the matrix A.  N >= 0.
71*> \endverbatim
72*>
73*> \param[in] NRHS
74*> \verbatim
75*>          NRHS is INTEGER
76*>          The number of right hand sides, i.e., the number of columns
77*>          of the matrices X and B.  NRHS >= 0.
78*> \endverbatim
79*>
80*> \param[in] AP
81*> \verbatim
82*>          AP is DOUBLE PRECISION array, dimension (N*(N+1)/2)
83*>          The upper or lower triangular matrix A, packed columnwise in
84*>          a linear array.  The j-th column of A is stored in the array
85*>          AP as follows:
86*>          if UPLO = 'U', AP((j-1)*j/2 + i) = A(i,j) for 1<=i<=j;
87*>          if UPLO = 'L',
88*>             AP((j-1)*(n-j) + j*(j+1)/2 + i-j) = A(i,j) for j<=i<=n.
89*> \endverbatim
90*>
91*> \param[in] SCALE
92*> \verbatim
93*>          SCALE is DOUBLE PRECISION
94*>          The scaling factor s used in solving the triangular system.
95*> \endverbatim
96*>
97*> \param[in] CNORM
98*> \verbatim
99*>          CNORM is DOUBLE PRECISION array, dimension (N)
100*>          The 1-norms of the columns of A, not counting the diagonal.
101*> \endverbatim
102*>
103*> \param[in] TSCAL
104*> \verbatim
105*>          TSCAL is DOUBLE PRECISION
106*>          The scaling factor used in computing the 1-norms in CNORM.
107*>          CNORM actually contains the column norms of TSCAL*A.
108*> \endverbatim
109*>
110*> \param[in] X
111*> \verbatim
112*>          X is DOUBLE PRECISION array, dimension (LDX,NRHS)
113*>          The computed solution vectors for the system of linear
114*>          equations.
115*> \endverbatim
116*>
117*> \param[in] LDX
118*> \verbatim
119*>          LDX is INTEGER
120*>          The leading dimension of the array X.  LDX >= max(1,N).
121*> \endverbatim
122*>
123*> \param[in] B
124*> \verbatim
125*>          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
126*>          The right hand side vectors for the system of linear
127*>          equations.
128*> \endverbatim
129*>
130*> \param[in] LDB
131*> \verbatim
132*>          LDB is INTEGER
133*>          The leading dimension of the array B.  LDB >= max(1,N).
134*> \endverbatim
135*>
136*> \param[out] WORK
137*> \verbatim
138*>          WORK is DOUBLE PRECISION array, dimension (N)
139*> \endverbatim
140*>
141*> \param[out] RESID
142*> \verbatim
143*>          RESID is DOUBLE PRECISION
144*>          The maximum over the number of right hand sides of
145*>          norm(op(A)*x - s*b) / ( norm(op(A)) * norm(x) * EPS ).
146*> \endverbatim
147*
148*  Authors:
149*  ========
150*
151*> \author Univ. of Tennessee
152*> \author Univ. of California Berkeley
153*> \author Univ. of Colorado Denver
154*> \author NAG Ltd.
155*
156*> \ingroup double_lin
157*
158*  =====================================================================
159      SUBROUTINE DTPT03( UPLO, TRANS, DIAG, N, NRHS, AP, SCALE, CNORM,
160     $                   TSCAL, X, LDX, B, LDB, WORK, RESID )
161*
162*  -- LAPACK test routine --
163*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
164*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
165*
166*     .. Scalar Arguments ..
167      CHARACTER          DIAG, TRANS, UPLO
168      INTEGER            LDB, LDX, N, NRHS
169      DOUBLE PRECISION   RESID, SCALE, TSCAL
170*     ..
171*     .. Array Arguments ..
172      DOUBLE PRECISION   AP( * ), B( LDB, * ), CNORM( * ), WORK( * ),
173     $                   X( LDX, * )
174*     ..
175*
176*  =====================================================================
177*
178*     .. Parameters ..
179      DOUBLE PRECISION   ONE, ZERO
180      PARAMETER          ( ONE = 1.0D+0, ZERO = 0.0D+0 )
181*     ..
182*     .. Local Scalars ..
183      INTEGER            IX, J, JJ
184      DOUBLE PRECISION   BIGNUM, EPS, ERR, SMLNUM, TNORM, XNORM, XSCAL
185*     ..
186*     .. External Functions ..
187      LOGICAL            LSAME
188      INTEGER            IDAMAX
189      DOUBLE PRECISION   DLAMCH
190      EXTERNAL           LSAME, IDAMAX, DLAMCH
191*     ..
192*     .. External Subroutines ..
193      EXTERNAL           DAXPY, DCOPY, DLABAD, DSCAL, DTPMV
194*     ..
195*     .. Intrinsic Functions ..
196      INTRINSIC          ABS, DBLE, MAX
197*     ..
198*     .. Executable Statements ..
199*
200*     Quick exit if N = 0.
201*
202      IF( N.LE.0 .OR. NRHS.LE.0 ) THEN
203         RESID = ZERO
204         RETURN
205      END IF
206      EPS = DLAMCH( 'Epsilon' )
207      SMLNUM = DLAMCH( 'Safe minimum' )
208      BIGNUM = ONE / SMLNUM
209      CALL DLABAD( SMLNUM, BIGNUM )
210*
211*     Compute the norm of the triangular matrix A using the column
212*     norms already computed by DLATPS.
213*
214      TNORM = ZERO
215      IF( LSAME( DIAG, 'N' ) ) THEN
216         IF( LSAME( UPLO, 'U' ) ) THEN
217            JJ = 1
218            DO 10 J = 1, N
219               TNORM = MAX( TNORM, TSCAL*ABS( AP( JJ ) )+CNORM( J ) )
220               JJ = JJ + J + 1
221   10       CONTINUE
222         ELSE
223            JJ = 1
224            DO 20 J = 1, N
225               TNORM = MAX( TNORM, TSCAL*ABS( AP( JJ ) )+CNORM( J ) )
226               JJ = JJ + N - J + 1
227   20       CONTINUE
228         END IF
229      ELSE
230         DO 30 J = 1, N
231            TNORM = MAX( TNORM, TSCAL+CNORM( J ) )
232   30    CONTINUE
233      END IF
234*
235*     Compute the maximum over the number of right hand sides of
236*        norm(op(A)*x - s*b) / ( norm(op(A)) * norm(x) * EPS ).
237*
238      RESID = ZERO
239      DO 40 J = 1, NRHS
240         CALL DCOPY( N, X( 1, J ), 1, WORK, 1 )
241         IX = IDAMAX( N, WORK, 1 )
242         XNORM = MAX( ONE, ABS( X( IX, J ) ) )
243         XSCAL = ( ONE / XNORM ) / DBLE( N )
244         CALL DSCAL( N, XSCAL, WORK, 1 )
245         CALL DTPMV( UPLO, TRANS, DIAG, N, AP, WORK, 1 )
246         CALL DAXPY( N, -SCALE*XSCAL, B( 1, J ), 1, WORK, 1 )
247         IX = IDAMAX( N, WORK, 1 )
248         ERR = TSCAL*ABS( WORK( IX ) )
249         IX = IDAMAX( N, X( 1, J ), 1 )
250         XNORM = ABS( X( IX, J ) )
251         IF( ERR*SMLNUM.LE.XNORM ) THEN
252            IF( XNORM.GT.ZERO )
253     $         ERR = ERR / XNORM
254         ELSE
255            IF( ERR.GT.ZERO )
256     $         ERR = ONE / EPS
257         END IF
258         IF( ERR*SMLNUM.LE.TNORM ) THEN
259            IF( TNORM.GT.ZERO )
260     $         ERR = ERR / TNORM
261         ELSE
262            IF( ERR.GT.ZERO )
263     $         ERR = ONE / EPS
264         END IF
265         RESID = MAX( RESID, ERR )
266   40 CONTINUE
267*
268      RETURN
269*
270*     End of DTPT03
271*
272      END
273