LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dchkrq ( logical, dimension( * )  DOTYPE,
integer  NM,
integer, dimension( * )  MVAL,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNB,
integer, dimension( * )  NBVAL,
integer, dimension( * )  NXVAL,
integer  NRHS,
double precision  THRESH,
logical  TSTERR,
integer  NMAX,
double precision, dimension( * )  A,
double precision, dimension( * )  AF,
double precision, dimension( * )  AQ,
double precision, dimension( * )  AR,
double precision, dimension( * )  AC,
double precision, dimension( * )  B,
double precision, dimension( * )  X,
double precision, dimension( * )  XACT,
double precision, dimension( * )  TAU,
double precision, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

DCHKRQ

Purpose:
 DCHKRQ tests DGERQF, DORGRQ and DORMRQ.
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]NRHS
          NRHS is INTEGER
          The number of right hand side vectors to be generated for
          each linear system.
[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.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[in]NMAX
          NMAX is INTEGER
          The maximum value permitted for M or N, used in dimensioning
          the work arrays.
[out]A
          A is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AF
          AF is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AQ
          AQ is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AR
          AR is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AC
          AC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]B
          B is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]X
          X is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]XACT
          XACT is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]TAU
          TAU is DOUBLE PRECISION array, dimension (NMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (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 203 of file dchkrq.f.

203 *
204 * -- LAPACK test routine (version 3.4.0) --
205 * -- LAPACK is a software package provided by Univ. of Tennessee, --
206 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
207 * November 2011
208 *
209 * .. Scalar Arguments ..
210  LOGICAL tsterr
211  INTEGER nm, nmax, nn, nnb, nout, nrhs
212  DOUBLE PRECISION thresh
213 * ..
214 * .. Array Arguments ..
215  LOGICAL dotype( * )
216  INTEGER iwork( * ), mval( * ), nbval( * ), nval( * ),
217  $ nxval( * )
218  DOUBLE PRECISION a( * ), ac( * ), af( * ), aq( * ), ar( * ),
219  $ b( * ), rwork( * ), tau( * ), work( * ),
220  $ x( * ), xact( * )
221 * ..
222 *
223 * =====================================================================
224 *
225 * .. Parameters ..
226  INTEGER ntests
227  parameter ( ntests = 7 )
228  INTEGER ntypes
229  parameter ( ntypes = 8 )
230  DOUBLE PRECISION zero
231  parameter ( zero = 0.0d0 )
232 * ..
233 * .. Local Scalars ..
234  CHARACTER dist, type
235  CHARACTER*3 path
236  INTEGER i, ik, im, imat, in, inb, info, k, kl, ku, lda,
237  $ lwork, m, minmn, mode, n, nb, nerrs, nfail, nk,
238  $ nrun, nt, nx
239  DOUBLE PRECISION anorm, cndnum
240 * ..
241 * .. Local Arrays ..
242  INTEGER iseed( 4 ), iseedy( 4 ), kval( 4 )
243  DOUBLE PRECISION result( ntests )
244 * ..
245 * .. External Subroutines ..
246  EXTERNAL alaerh, alahd, alasum, derrrq, dgerqs, dget02,
248  $ drqt03, xlaenv
249 * ..
250 * .. Intrinsic Functions ..
251  INTRINSIC max, min
252 * ..
253 * .. Scalars in Common ..
254  LOGICAL lerr, ok
255  CHARACTER*32 srnamt
256  INTEGER infot, nunit
257 * ..
258 * .. Common blocks ..
259  COMMON / infoc / infot, nunit, ok, lerr
260  COMMON / srnamc / srnamt
261 * ..
262 * .. Data statements ..
263  DATA iseedy / 1988, 1989, 1990, 1991 /
264 * ..
265 * .. Executable Statements ..
266 *
267 * Initialize constants and the random number seed.
268 *
269  path( 1: 1 ) = 'Double precision'
270  path( 2: 3 ) = 'RQ'
271  nrun = 0
272  nfail = 0
273  nerrs = 0
274  DO 10 i = 1, 4
275  iseed( i ) = iseedy( i )
276  10 CONTINUE
277 *
278 * Test the error exits
279 *
280  IF( tsterr )
281  $ CALL derrrq( path, nout )
282  infot = 0
283  CALL xlaenv( 2, 2 )
284 *
285  lda = nmax
286  lwork = nmax*max( nmax, nrhs )
287 *
288 * Do for each value of M in MVAL.
289 *
290  DO 70 im = 1, nm
291  m = mval( im )
292 *
293 * Do for each value of N in NVAL.
294 *
295  DO 60 in = 1, nn
296  n = nval( in )
297  minmn = min( m, n )
298  DO 50 imat = 1, ntypes
299 *
300 * Do the tests only if DOTYPE( IMAT ) is true.
301 *
302  IF( .NOT.dotype( imat ) )
303  $ GO TO 50
304 *
305 * Set up parameters with DLATB4 and generate a test matrix
306 * with DLATMS.
307 *
308  CALL dlatb4( path, imat, m, n, TYPE, kl, ku, anorm, mode,
309  $ cndnum, dist )
310 *
311  srnamt = 'DLATMS'
312  CALL dlatms( m, n, dist, iseed, TYPE, rwork, mode,
313  $ cndnum, anorm, kl, ku, 'No packing', a, lda,
314  $ work, info )
315 *
316 * Check error code from DLATMS.
317 *
318  IF( info.NE.0 ) THEN
319  CALL alaerh( path, 'DLATMS', info, 0, ' ', m, n, -1,
320  $ -1, -1, imat, nfail, nerrs, nout )
321  GO TO 50
322  END IF
323 *
324 * Set some values for K: the first value must be MINMN,
325 * corresponding to the call of DRQT01; other values are
326 * used in the calls of DRQT02, and must not exceed MINMN.
327 *
328  kval( 1 ) = minmn
329  kval( 2 ) = 0
330  kval( 3 ) = 1
331  kval( 4 ) = minmn / 2
332  IF( minmn.EQ.0 ) THEN
333  nk = 1
334  ELSE IF( minmn.EQ.1 ) THEN
335  nk = 2
336  ELSE IF( minmn.LE.3 ) THEN
337  nk = 3
338  ELSE
339  nk = 4
340  END IF
341 *
342 * Do for each value of K in KVAL
343 *
344  DO 40 ik = 1, nk
345  k = kval( ik )
346 *
347 * Do for each pair of values (NB,NX) in NBVAL and NXVAL.
348 *
349  DO 30 inb = 1, nnb
350  nb = nbval( inb )
351  CALL xlaenv( 1, nb )
352  nx = nxval( inb )
353  CALL xlaenv( 3, nx )
354  DO i = 1, ntests
355  result( i ) = zero
356  END DO
357  nt = 2
358  IF( ik.EQ.1 ) THEN
359 *
360 * Test DGERQF
361 *
362  CALL drqt01( m, n, a, af, aq, ar, lda, tau,
363  $ work, lwork, rwork, result( 1 ) )
364  ELSE IF( m.LE.n ) THEN
365 *
366 * Test DORGRQ, using factorization
367 * returned by DRQT01
368 *
369  CALL drqt02( m, n, k, a, af, aq, ar, lda, tau,
370  $ work, lwork, rwork, result( 1 ) )
371 
372  END IF
373  IF( m.GE.k ) THEN
374 *
375 * Test DORMRQ, using factorization returned
376 * by DRQT01
377 *
378  CALL drqt03( m, n, k, af, ac, ar, aq, lda, tau,
379  $ work, lwork, rwork, result( 3 ) )
380  nt = nt + 4
381 *
382 * If M>=N and K=N, call DGERQS to solve a system
383 * with NRHS right hand sides and compute the
384 * residual.
385 *
386  IF( k.EQ.m .AND. inb.EQ.1 ) THEN
387 *
388 * Generate a solution and set the right
389 * hand side.
390 *
391  srnamt = 'DLARHS'
392  CALL dlarhs( path, 'New', 'Full',
393  $ 'No transpose', m, n, 0, 0,
394  $ nrhs, a, lda, xact, lda, b, lda,
395  $ iseed, info )
396 *
397  CALL dlacpy( 'Full', m, nrhs, b, lda,
398  $ x( n-m+1 ), lda )
399  srnamt = 'DGERQS'
400  CALL dgerqs( m, n, nrhs, af, lda, tau, x,
401  $ lda, work, lwork, info )
402 *
403 * Check error code from DGERQS.
404 *
405  IF( info.NE.0 )
406  $ CALL alaerh( path, 'DGERQS', info, 0, ' ',
407  $ m, n, nrhs, -1, nb, imat,
408  $ nfail, nerrs, nout )
409 *
410  CALL dget02( 'No transpose', m, n, nrhs, a,
411  $ lda, x, lda, b, lda, rwork,
412  $ result( 7 ) )
413  nt = nt + 1
414  END IF
415  END IF
416 *
417 * Print information about the tests that did not
418 * pass the threshold.
419 *
420  DO 20 i = 1, nt
421  IF( result( i ).GE.thresh ) THEN
422  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
423  $ CALL alahd( nout, path )
424  WRITE( nout, fmt = 9999 )m, n, k, nb, nx,
425  $ imat, i, result( i )
426  nfail = nfail + 1
427  END IF
428  20 CONTINUE
429  nrun = nrun + nt
430  30 CONTINUE
431  40 CONTINUE
432  50 CONTINUE
433  60 CONTINUE
434  70 CONTINUE
435 *
436 * Print a summary of the results.
437 *
438  CALL alasum( path, nout, nfail, nrun, nerrs )
439 *
440  9999 FORMAT( ' M=', i5, ', N=', i5, ', K=', i5, ', NB=', i4, ', NX=',
441  $ i5, ', type ', i2, ', test(', i2, ')=', g12.5 )
442  RETURN
443 *
444 * End of DCHKRQ
445 *
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine dlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
DLARHS
Definition: dlarhs.f:206
subroutine drqt01(M, N, A, AF, Q, R, LDA, TAU, WORK, LWORK, RWORK, RESULT)
DRQT01
Definition: drqt01.f:128
subroutine drqt02(M, N, K, A, AF, Q, R, LDA, TAU, WORK, LWORK, RWORK, RESULT)
DRQT02
Definition: drqt02.f:138
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 xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine dgerqs(M, N, NRHS, A, LDA, TAU, B, LDB, WORK, LWORK, INFO)
DGERQS
Definition: dgerqs.f:124
subroutine dget02(TRANS, M, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
DGET02
Definition: dget02.f:135
subroutine dlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
DLATB4
Definition: dlatb4.f:122
subroutine derrrq(PATH, NUNIT)
DERRRQ
Definition: derrrq.f:57
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine drqt03(M, N, K, AF, C, CC, Q, LDA, TAU, WORK, LWORK, RWORK, RESULT)
DRQT03
Definition: drqt03.f:138
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: