LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dchkq3 ( logical, dimension( * )  DOTYPE,
integer  NM,
integer, dimension( * )  MVAL,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNB,
integer, dimension( * )  NBVAL,
integer, dimension( * )  NXVAL,
double precision  THRESH,
double precision, dimension( * )  A,
double precision, dimension( * )  COPYA,
double precision, dimension( * )  S,
double precision, dimension( * )  TAU,
double precision, dimension( * )  WORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

DCHKQ3

Purpose:
 DCHKQ3 tests DGEQP3.
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
[in]NM
          NM is INTEGER
          The number of values of M contained in the vector MVAL.
[in]MVAL
          MVAL is INTEGER array, dimension (NM)
          The values of the matrix row dimension M.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix column dimension N.
[in]NNB
          NNB is INTEGER
          The number of values of NB and NX contained in the
          vectors NBVAL and NXVAL.  The blocking parameters are used
          in pairs (NB,NX).
[in]NBVAL
          NBVAL is INTEGER array, dimension (NNB)
          The values of the blocksize NB.
[in]NXVAL
          NXVAL is INTEGER array, dimension (NNB)
          The values of the crossover point NX.
[in]THRESH
          THRESH is DOUBLE PRECISION
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[out]A
          A is DOUBLE PRECISION array, dimension (MMAX*NMAX)
          where MMAX is the maximum value of M in MVAL and NMAX is the
          maximum value of N in NVAL.
[out]COPYA
          COPYA is DOUBLE PRECISION array, dimension (MMAX*NMAX)
[out]S
          S is DOUBLE PRECISION array, dimension
                      (min(MMAX,NMAX))
[out]TAU
          TAU is DOUBLE PRECISION array, dimension (MMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension
                      (MMAX*NMAX + 4*NMAX + MMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (2*NMAX)
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 155 of file dchkq3.f.

155 *
156 * -- LAPACK test routine (version 3.4.0) --
157 * -- LAPACK is a software package provided by Univ. of Tennessee, --
158 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
159 * November 2011
160 *
161 * .. Scalar Arguments ..
162  INTEGER nm, nn, nnb, nout
163  DOUBLE PRECISION thresh
164 * ..
165 * .. Array Arguments ..
166  LOGICAL dotype( * )
167  INTEGER iwork( * ), mval( * ), nbval( * ), nval( * ),
168  $ nxval( * )
169  DOUBLE PRECISION a( * ), copya( * ), s( * ),
170  $ tau( * ), work( * )
171 * ..
172 *
173 * =====================================================================
174 *
175 * .. Parameters ..
176  INTEGER ntypes
177  parameter ( ntypes = 6 )
178  INTEGER ntests
179  parameter ( ntests = 3 )
180  DOUBLE PRECISION one, zero
181  parameter ( one = 1.0d0, zero = 0.0d0 )
182 * ..
183 * .. Local Scalars ..
184  CHARACTER*3 path
185  INTEGER i, ihigh, ilow, im, imode, in, inb, info,
186  $ istep, k, lda, lw, lwork, m, mnmin, mode, n,
187  $ nb, nerrs, nfail, nrun, nx
188  DOUBLE PRECISION eps
189 * ..
190 * .. Local Arrays ..
191  INTEGER iseed( 4 ), iseedy( 4 )
192  DOUBLE PRECISION result( ntests )
193 * ..
194 * .. External Functions ..
195  DOUBLE PRECISION dlamch, dqpt01, dqrt11, dqrt12
196  EXTERNAL dlamch, dqpt01, dqrt11, dqrt12
197 * ..
198 * .. External Subroutines ..
199  EXTERNAL alahd, alasum, dgeqp3, dlacpy, dlaord, dlaset,
200  $ dlatms, icopy, xlaenv
201 * ..
202 * .. Intrinsic Functions ..
203  INTRINSIC max, min
204 * ..
205 * .. Scalars in Common ..
206  LOGICAL lerr, ok
207  CHARACTER*32 srnamt
208  INTEGER infot, iounit
209 * ..
210 * .. Common blocks ..
211  COMMON / infoc / infot, iounit, ok, lerr
212  COMMON / srnamc / srnamt
213 * ..
214 * .. Data statements ..
215  DATA iseedy / 1988, 1989, 1990, 1991 /
216 * ..
217 * .. Executable Statements ..
218 *
219 * Initialize constants and the random number seed.
220 *
221  path( 1: 1 ) = 'Double precision'
222  path( 2: 3 ) = 'Q3'
223  nrun = 0
224  nfail = 0
225  nerrs = 0
226  DO 10 i = 1, 4
227  iseed( i ) = iseedy( i )
228  10 CONTINUE
229  eps = dlamch( 'Epsilon' )
230  infot = 0
231 *
232  DO 90 im = 1, nm
233 *
234 * Do for each value of M in MVAL.
235 *
236  m = mval( im )
237  lda = max( 1, m )
238 *
239  DO 80 in = 1, nn
240 *
241 * Do for each value of N in NVAL.
242 *
243  n = nval( in )
244  mnmin = min( m, n )
245  lwork = max( 1, m*max( m, n )+4*mnmin+max( m, n ),
246  $ m*n + 2*mnmin + 4*n )
247 *
248  DO 70 imode = 1, ntypes
249  IF( .NOT.dotype( imode ) )
250  $ GO TO 70
251 *
252 * Do for each type of matrix
253 * 1: zero matrix
254 * 2: one small singular value
255 * 3: geometric distribution of singular values
256 * 4: first n/2 columns fixed
257 * 5: last n/2 columns fixed
258 * 6: every second column fixed
259 *
260  mode = imode
261  IF( imode.GT.3 )
262  $ mode = 1
263 *
264 * Generate test matrix of size m by n using
265 * singular value distribution indicated by `mode'.
266 *
267  DO 20 i = 1, n
268  iwork( i ) = 0
269  20 CONTINUE
270  IF( imode.EQ.1 ) THEN
271  CALL dlaset( 'Full', m, n, zero, zero, copya, lda )
272  DO 30 i = 1, mnmin
273  s( i ) = zero
274  30 CONTINUE
275  ELSE
276  CALL dlatms( m, n, 'Uniform', iseed, 'Nonsymm', s,
277  $ mode, one / eps, one, m, n, 'No packing',
278  $ copya, lda, work, info )
279  IF( imode.GE.4 ) THEN
280  IF( imode.EQ.4 ) THEN
281  ilow = 1
282  istep = 1
283  ihigh = max( 1, n / 2 )
284  ELSE IF( imode.EQ.5 ) THEN
285  ilow = max( 1, n / 2 )
286  istep = 1
287  ihigh = n
288  ELSE IF( imode.EQ.6 ) THEN
289  ilow = 1
290  istep = 2
291  ihigh = n
292  END IF
293  DO 40 i = ilow, ihigh, istep
294  iwork( i ) = 1
295  40 CONTINUE
296  END IF
297  CALL dlaord( 'Decreasing', mnmin, s, 1 )
298  END IF
299 *
300  DO 60 inb = 1, nnb
301 *
302 * Do for each pair of values (NB,NX) in NBVAL and NXVAL.
303 *
304  nb = nbval( inb )
305  CALL xlaenv( 1, nb )
306  nx = nxval( inb )
307  CALL xlaenv( 3, nx )
308 *
309 * Get a working copy of COPYA into A and a copy of
310 * vector IWORK.
311 *
312  CALL dlacpy( 'All', m, n, copya, lda, a, lda )
313  CALL icopy( n, iwork( 1 ), 1, iwork( n+1 ), 1 )
314 *
315 * Compute the QR factorization with pivoting of A
316 *
317  lw = max( 1, 2*n+nb*( n+1 ) )
318 *
319 * Compute the QP3 factorization of A
320 *
321  srnamt = 'DGEQP3'
322  CALL dgeqp3( m, n, a, lda, iwork( n+1 ), tau, work,
323  $ lw, info )
324 *
325 * Compute norm(svd(a) - svd(r))
326 *
327  result( 1 ) = dqrt12( m, n, a, lda, s, work,
328  $ lwork )
329 *
330 * Compute norm( A*P - Q*R )
331 *
332  result( 2 ) = dqpt01( m, n, mnmin, copya, a, lda, tau,
333  $ iwork( n+1 ), work, lwork )
334 *
335 * Compute Q'*Q
336 *
337  result( 3 ) = dqrt11( m, mnmin, a, lda, tau, work,
338  $ lwork )
339 *
340 * Print information about the tests that did not pass
341 * the threshold.
342 *
343  DO 50 k = 1, ntests
344  IF( result( k ).GE.thresh ) THEN
345  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
346  $ CALL alahd( nout, path )
347  WRITE( nout, fmt = 9999 )'DGEQP3', m, n, nb,
348  $ imode, k, result( k )
349  nfail = nfail + 1
350  END IF
351  50 CONTINUE
352  nrun = nrun + ntests
353 *
354  60 CONTINUE
355  70 CONTINUE
356  80 CONTINUE
357  90 CONTINUE
358 *
359 * Print a summary of the results.
360 *
361  CALL alasum( path, nout, nfail, nrun, nerrs )
362 *
363  9999 FORMAT( 1x, a, ' M =', i5, ', N =', i5, ', NB =', i4, ', type ',
364  $ i2, ', test ', i2, ', ratio =', g12.5 )
365 *
366 * End of DCHKQ3
367 *
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: dlaset.f:112
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine icopy(N, SX, INCX, SY, INCY)
ICOPY
Definition: icopy.f:77
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
double precision function dqrt12(M, N, A, LDA, S, WORK, LWORK)
DQRT12
Definition: dqrt12.f:91
subroutine dgeqp3(M, N, A, LDA, JPVT, TAU, WORK, LWORK, INFO)
DGEQP3
Definition: dgeqp3.f:153
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dlaord(JOB, N, X, INCX)
DLAORD
Definition: dlaord.f:75
double precision function dqrt11(M, K, A, LDA, TAU, WORK, LWORK)
DQRT11
Definition: dqrt11.f:100
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
double precision function dqpt01(M, N, K, A, AF, LDA, TAU, JPVT, WORK, LWORK)
DQPT01
Definition: dqpt01.f:122
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:75

Here is the call graph for this function:

Here is the caller graph for this function: