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

◆ dqrt14()

double precision function dqrt14 ( character  trans,
integer  m,
integer  n,
integer  nrhs,
double precision, dimension( lda, * )  a,
integer  lda,
double precision, dimension( ldx, * )  x,
integer  ldx,
double precision, dimension( lwork )  work,
integer  lwork 
)

DQRT14

Purpose:
 DQRT14 checks whether X is in the row space of A or A'.  It does so
 by scaling both X and A such that their norms are in the range
 [sqrt(eps), 1/sqrt(eps)], then computing a QR factorization of [A,X]
 (if TRANS = 'T') or an LQ factorization of [A',X]' (if TRANS = 'N'),
 and returning the norm of the trailing triangle, scaled by
 MAX(M,N,NRHS)*eps.
Parameters
[in]TRANS
          TRANS is CHARACTER*1
          = 'N':  No transpose, check for X in the row space of A
          = 'T':  Transpose, check for X in the row space of A'.
[in]M
          M is INTEGER
          The number of rows of the matrix A.
[in]N
          N is INTEGER
          The number of columns of the matrix A.
[in]NRHS
          NRHS is INTEGER
          The number of right hand sides, i.e., the number of columns
          of X.
[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.
[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.
[out]WORK
          WORK is DOUBLE PRECISION array dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          length of workspace array required
          If TRANS = 'N', LWORK >= (M+NRHS)*(N+2);
          if TRANS = 'T', LWORK >= (N+NRHS)*(M+2).
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 114 of file dqrt14.f.

116*
117* -- LAPACK test routine --
118* -- LAPACK is a software package provided by Univ. of Tennessee, --
119* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
120*
121* .. Scalar Arguments ..
122 CHARACTER TRANS
123 INTEGER LDA, LDX, LWORK, M, N, NRHS
124* ..
125* .. Array Arguments ..
126 DOUBLE PRECISION A( LDA, * ), WORK( LWORK ), X( LDX, * )
127* ..
128*
129* =====================================================================
130*
131* .. Parameters ..
132 DOUBLE PRECISION ZERO, ONE
133 parameter( zero = 0.0d0, one = 1.0d0 )
134* ..
135* .. Local Scalars ..
136 LOGICAL TPSD
137 INTEGER I, INFO, J, LDWORK
138 DOUBLE PRECISION ANRM, ERR, XNRM
139* ..
140* .. Local Arrays ..
141 DOUBLE PRECISION RWORK( 1 )
142* ..
143* .. External Functions ..
144 LOGICAL LSAME
145 DOUBLE PRECISION DLAMCH, DLANGE
146 EXTERNAL lsame, dlamch, dlange
147* ..
148* .. External Subroutines ..
149 EXTERNAL dgelq2, dgeqr2, dlacpy, dlascl, xerbla
150* ..
151* .. Intrinsic Functions ..
152 INTRINSIC abs, dble, max, min
153* ..
154* .. Executable Statements ..
155*
156 dqrt14 = zero
157 IF( lsame( trans, 'N' ) ) THEN
158 ldwork = m + nrhs
159 tpsd = .false.
160 IF( lwork.LT.( m+nrhs )*( n+2 ) ) THEN
161 CALL xerbla( 'DQRT14', 10 )
162 RETURN
163 ELSE IF( n.LE.0 .OR. nrhs.LE.0 ) THEN
164 RETURN
165 END IF
166 ELSE IF( lsame( trans, 'T' ) ) THEN
167 ldwork = m
168 tpsd = .true.
169 IF( lwork.LT.( n+nrhs )*( m+2 ) ) THEN
170 CALL xerbla( 'DQRT14', 10 )
171 RETURN
172 ELSE IF( m.LE.0 .OR. nrhs.LE.0 ) THEN
173 RETURN
174 END IF
175 ELSE
176 CALL xerbla( 'DQRT14', 1 )
177 RETURN
178 END IF
179*
180* Copy and scale A
181*
182 CALL dlacpy( 'All', m, n, a, lda, work, ldwork )
183 anrm = dlange( 'M', m, n, work, ldwork, rwork )
184 IF( anrm.NE.zero )
185 $ CALL dlascl( 'G', 0, 0, anrm, one, m, n, work, ldwork, info )
186*
187* Copy X or X' into the right place and scale it
188*
189 IF( tpsd ) THEN
190*
191* Copy X into columns n+1:n+nrhs of work
192*
193 CALL dlacpy( 'All', m, nrhs, x, ldx, work( n*ldwork+1 ),
194 $ ldwork )
195 xnrm = dlange( 'M', m, nrhs, work( n*ldwork+1 ), ldwork,
196 $ rwork )
197 IF( xnrm.NE.zero )
198 $ CALL dlascl( 'G', 0, 0, xnrm, one, m, nrhs,
199 $ work( n*ldwork+1 ), ldwork, info )
200*
201* Compute QR factorization of X
202*
203 CALL dgeqr2( m, n+nrhs, work, ldwork,
204 $ work( ldwork*( n+nrhs )+1 ),
205 $ work( ldwork*( n+nrhs )+min( m, n+nrhs )+1 ),
206 $ info )
207*
208* Compute largest entry in upper triangle of
209* work(n+1:m,n+1:n+nrhs)
210*
211 err = zero
212 DO 20 j = n + 1, n + nrhs
213 DO 10 i = n + 1, min( m, j )
214 err = max( err, abs( work( i+( j-1 )*m ) ) )
215 10 CONTINUE
216 20 CONTINUE
217*
218 ELSE
219*
220* Copy X' into rows m+1:m+nrhs of work
221*
222 DO 40 i = 1, n
223 DO 30 j = 1, nrhs
224 work( m+j+( i-1 )*ldwork ) = x( i, j )
225 30 CONTINUE
226 40 CONTINUE
227*
228 xnrm = dlange( 'M', nrhs, n, work( m+1 ), ldwork, rwork )
229 IF( xnrm.NE.zero )
230 $ CALL dlascl( 'G', 0, 0, xnrm, one, nrhs, n, work( m+1 ),
231 $ ldwork, info )
232*
233* Compute LQ factorization of work
234*
235 CALL dgelq2( ldwork, n, work, ldwork, work( ldwork*n+1 ),
236 $ work( ldwork*( n+1 )+1 ), info )
237*
238* Compute largest entry in lower triangle in
239* work(m+1:m+nrhs,m+1:n)
240*
241 err = zero
242 DO 60 j = m + 1, n
243 DO 50 i = j, ldwork
244 err = max( err, abs( work( i+( j-1 )*ldwork ) ) )
245 50 CONTINUE
246 60 CONTINUE
247*
248 END IF
249*
250 dqrt14 = err / ( dble( max( m, n, nrhs ) )*dlamch( 'Epsilon' ) )
251*
252 RETURN
253*
254* End of DQRT14
255*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
double precision function dqrt14(trans, m, n, nrhs, a, lda, x, ldx, work, lwork)
DQRT14
Definition dqrt14.f:116
subroutine dgelq2(m, n, a, lda, tau, work, info)
DGELQ2 computes the LQ factorization of a general rectangular matrix using an unblocked algorithm.
Definition dgelq2.f:129
subroutine dgeqr2(m, n, a, lda, tau, work, info)
DGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.
Definition dgeqr2.f:130
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: