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

◆ dqrt17()

double precision function dqrt17 ( character  trans,
integer  iresid,
integer  m,
integer  n,
integer  nrhs,
double precision, dimension( lda, * )  a,
integer  lda,
double precision, dimension( ldx, * )  x,
integer  ldx,
double precision, dimension( ldb, * )  b,
integer  ldb,
double precision, dimension( ldb, * )  c,
double precision, dimension( lwork )  work,
integer  lwork 
)

DQRT17

Purpose:
 DQRT17 computes the ratio

    norm(R**T * op(A)) / ( norm(A) * alpha * max(M,N,NRHS) * EPS ),

 where R = B - op(A)*X, op(A) is A or A**T, depending on TRANS, EPS
 is the machine epsilon, and

    alpha = norm(B) if IRESID = 1 (zero-residual problem)
    alpha = norm(R) if IRESID = 2 (otherwise).

 The norm used is the 1-norm.
Parameters
[in]TRANS
          TRANS is CHARACTER*1
          Specifies whether or not the transpose of A is used.
          = 'N':  No transpose, op(A) = A.
          = 'T':  Transpose, op(A) = A**T.
[in]IRESID
          IRESID is INTEGER
          IRESID = 1 indicates zero-residual problem.
          IRESID = 2 indicates non-zero residual.
[in]M
          M is INTEGER
          The number of rows of the matrix A.
          If TRANS = 'N', the number of rows of the matrix B.
          If TRANS = 'T', the number of rows of the matrix X.
[in]N
          N is INTEGER
          The number of columns of the matrix  A.
          If TRANS = 'N', the number of rows of the matrix X.
          If TRANS = 'T', the number of rows of the matrix B.
[in]NRHS
          NRHS is INTEGER
          The number of columns of the matrices X and B.
[in]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          The m-by-n matrix A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A. LDA >= M.
[in]X
          X is DOUBLE PRECISION array, dimension (LDX,NRHS)
          If TRANS = 'N', the n-by-nrhs matrix X.
          If TRANS = 'T', the m-by-nrhs matrix X.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.
          If TRANS = 'N', LDX >= N.
          If TRANS = 'T', LDX >= M.
[in]B
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
          If TRANS = 'N', the m-by-nrhs matrix B.
          If TRANS = 'T', the n-by-nrhs matrix B.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.
          If TRANS = 'N', LDB >= M.
          If TRANS = 'T', LDB >= N.
[out]C
          C is DOUBLE PRECISION array, dimension (LDB,NRHS)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The length of the array WORK.  LWORK >= NRHS*(M+N).
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 151 of file dqrt17.f.

153*
154* -- LAPACK test routine --
155* -- LAPACK is a software package provided by Univ. of Tennessee, --
156* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
157*
158* .. Scalar Arguments ..
159 CHARACTER TRANS
160 INTEGER IRESID, LDA, LDB, LDX, LWORK, M, N, NRHS
161* ..
162* .. Array Arguments ..
163 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), C( LDB, * ),
164 $ WORK( LWORK ), X( LDX, * )
165* ..
166*
167* =====================================================================
168*
169* .. Parameters ..
170 DOUBLE PRECISION ZERO, ONE
171 parameter( zero = 0.0d0, one = 1.0d0 )
172* ..
173* .. Local Scalars ..
174 INTEGER INFO, ISCL, NCOLS, NROWS
175 DOUBLE PRECISION ERR, NORMA, NORMB, NORMRS, SMLNUM
176* ..
177* .. Local Arrays ..
178 DOUBLE PRECISION RWORK( 1 )
179* ..
180* .. External Functions ..
181 LOGICAL LSAME
182 DOUBLE PRECISION DLAMCH, DLANGE
183 EXTERNAL lsame, dlamch, dlange
184* ..
185* .. External Subroutines ..
186 EXTERNAL dgemm, dlacpy, dlascl, xerbla
187* ..
188* .. Intrinsic Functions ..
189 INTRINSIC dble, max
190* ..
191* .. Executable Statements ..
192*
193 dqrt17 = zero
194*
195 IF( lsame( trans, 'N' ) ) THEN
196 nrows = m
197 ncols = n
198 ELSE IF( lsame( trans, 'T' ) ) THEN
199 nrows = n
200 ncols = m
201 ELSE
202 CALL xerbla( 'DQRT17', 1 )
203 RETURN
204 END IF
205*
206 IF( lwork.LT.ncols*nrhs ) THEN
207 CALL xerbla( 'DQRT17', 13 )
208 RETURN
209 END IF
210*
211 IF( m.LE.0 .OR. n.LE.0 .OR. nrhs.LE.0 ) THEN
212 RETURN
213 END IF
214*
215 norma = dlange( 'One-norm', m, n, a, lda, rwork )
216 smlnum = dlamch( 'Safe minimum' ) / dlamch( 'Precision' )
217 iscl = 0
218*
219* compute residual and scale it
220*
221 CALL dlacpy( 'All', nrows, nrhs, b, ldb, c, ldb )
222 CALL dgemm( trans, 'No transpose', nrows, nrhs, ncols, -one, a,
223 $ lda, x, ldx, one, c, ldb )
224 normrs = dlange( 'Max', nrows, nrhs, c, ldb, rwork )
225 IF( normrs.GT.smlnum ) THEN
226 iscl = 1
227 CALL dlascl( 'General', 0, 0, normrs, one, nrows, nrhs, c, ldb,
228 $ info )
229 END IF
230*
231* compute R**T * op(A)
232*
233 CALL dgemm( 'Transpose', trans, nrhs, ncols, nrows, one, c, ldb,
234 $ a, lda, zero, work, nrhs )
235*
236* compute and properly scale error
237*
238 err = dlange( 'One-norm', nrhs, ncols, work, nrhs, rwork )
239 IF( norma.NE.zero )
240 $ err = err / norma
241*
242 IF( iscl.EQ.1 )
243 $ err = err*normrs
244*
245 IF( iresid.EQ.1 ) THEN
246 normb = dlange( 'One-norm', nrows, nrhs, b, ldb, rwork )
247 IF( normb.NE.zero )
248 $ err = err / normb
249 ELSE
250 IF( normrs.NE.zero )
251 $ err = err / normrs
252 END IF
253*
254 dqrt17 = err / ( dlamch( 'Epsilon' )*dble( max( m, n, nrhs ) ) )
255 RETURN
256*
257* End of DQRT17
258*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
double precision function dqrt17(trans, iresid, m, n, nrhs, a, lda, x, ldx, b, ldb, c, work, lwork)
DQRT17
Definition dqrt17.f:153
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlange(norm, m, n, a, lda, work)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlange.f:114
subroutine dlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition dlascl.f:143
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: