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

SCHKRQ

Purpose:
 SCHKRQ tests SGERQF, SORGRQ and SORMRQ.
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 REAL
          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 REAL array, dimension (NMAX*NMAX)
[out]AF
          AF is REAL array, dimension (NMAX*NMAX)
[out]AQ
          AQ is REAL array, dimension (NMAX*NMAX)
[out]AR
          AR is REAL array, dimension (NMAX*NMAX)
[out]AC
          AC is REAL array, dimension (NMAX*NMAX)
[out]B
          B is REAL array, dimension (NMAX*NRHS)
[out]X
          X is REAL array, dimension (NMAX*NRHS)
[out]XACT
          XACT is REAL array, dimension (NMAX*NRHS)
[out]TAU
          TAU is REAL array, dimension (NMAX)
[out]WORK
          WORK is REAL array, dimension (NMAX*NMAX)
[out]RWORK
          RWORK is REAL 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 2015

Definition at line 203 of file schkrq.f.

203 *
204 * -- LAPACK test routine (version 3.6.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 2015
208 *
209 * .. Scalar Arguments ..
210  LOGICAL tsterr
211  INTEGER nm, nmax, nn, nnb, nout, nrhs
212  REAL thresh
213 * ..
214 * .. Array Arguments ..
215  LOGICAL dotype( * )
216  INTEGER iwork( * ), mval( * ), nbval( * ), nval( * ),
217  $ nxval( * )
218  REAL 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  REAL zero
231  parameter ( zero = 0.0e0 )
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  REAL anorm, cndnum
240 * ..
241 * .. Local Arrays ..
242  INTEGER iseed( 4 ), iseedy( 4 ), kval( 4 )
243  REAL result( ntests )
244 * ..
245 * .. External Subroutines ..
246  EXTERNAL alaerh, alahd, alasum, serrrq, sgerqs, sget02,
248  $ srqt03, 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 ) = 'Single 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 serrrq( 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 SLATB4 and generate a test matrix
306 * with SLATMS.
307 *
308  CALL slatb4( path, imat, m, n, TYPE, kl, ku, anorm, mode,
309  $ cndnum, dist )
310 *
311  srnamt = 'SLATMS'
312  CALL slatms( 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 SLATMS.
317 *
318  IF( info.NE.0 ) THEN
319  CALL alaerh( path, 'SLATMS', 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 SRQT01; other values are
326 * used in the calls of SRQT02, 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 SGERQF
361 *
362  CALL srqt01( m, n, a, af, aq, ar, lda, tau,
363  $ work, lwork, rwork, result( 1 ) )
364  ELSE IF( m.LE.n ) THEN
365 *
366 * Test SORGRQ, using factorization
367 * returned by SRQT01
368 *
369  CALL srqt02( m, n, k, a, af, aq, ar, lda, tau,
370  $ work, lwork, rwork, result( 1 ) )
371  END IF
372  IF( m.GE.k ) THEN
373 *
374 * Test SORMRQ, using factorization returned
375 * by SRQT01
376 *
377  CALL srqt03( m, n, k, af, ac, ar, aq, lda, tau,
378  $ work, lwork, rwork, result( 3 ) )
379  nt = nt + 4
380 *
381 * If M>=N and K=N, call SGERQS to solve a system
382 * with NRHS right hand sides and compute the
383 * residual.
384 *
385  IF( k.EQ.m .AND. inb.EQ.1 ) THEN
386 *
387 * Generate a solution and set the right
388 * hand side.
389 *
390  srnamt = 'SLARHS'
391  CALL slarhs( path, 'New', 'Full',
392  $ 'No transpose', m, n, 0, 0,
393  $ nrhs, a, lda, xact, lda, b, lda,
394  $ iseed, info )
395 *
396  CALL slacpy( 'Full', m, nrhs, b, lda,
397  $ x( n-m+1 ), lda )
398  srnamt = 'SGERQS'
399  CALL sgerqs( m, n, nrhs, af, lda, tau, x,
400  $ lda, work, lwork, info )
401 *
402 * Check error code from SGERQS.
403 *
404  IF( info.NE.0 )
405  $ CALL alaerh( path, 'SGERQS', info, 0, ' ',
406  $ m, n, nrhs, -1, nb, imat,
407  $ nfail, nerrs, nout )
408 *
409  CALL sget02( 'No transpose', m, n, nrhs, a,
410  $ lda, x, lda, b, lda, rwork,
411  $ result( 7 ) )
412  nt = nt + 1
413  END IF
414  END IF
415 *
416 * Print information about the tests that did not
417 * pass the threshold.
418 *
419  DO 20 i = 1, nt
420  IF( result( i ).GE.thresh ) THEN
421  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
422  $ CALL alahd( nout, path )
423  WRITE( nout, fmt = 9999 )m, n, k, nb, nx,
424  $ imat, i, result( i )
425  nfail = nfail + 1
426  END IF
427  20 CONTINUE
428  nrun = nrun + nt
429  30 CONTINUE
430  40 CONTINUE
431  50 CONTINUE
432  60 CONTINUE
433  70 CONTINUE
434 *
435 * Print a summary of the results.
436 *
437  CALL alasum( path, nout, nfail, nrun, nerrs )
438 *
439  9999 FORMAT( ' M=', i5, ', N=', i5, ', K=', i5, ', NB=', i4, ', NX=',
440  $ i5, ', type ', i2, ', test(', i2, ')=', g12.5 )
441  RETURN
442 *
443 * End of SCHKRQ
444 *
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 sget02(TRANS, M, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
SGET02
Definition: sget02.f:135
subroutine slatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
SLATB4
Definition: slatb4.f:122
subroutine srqt02(M, N, K, A, AF, Q, R, LDA, TAU, WORK, LWORK, RWORK, RESULT)
SRQT02
Definition: srqt02.f:138
subroutine slarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
SLARHS
Definition: slarhs.f:206
subroutine sgerqs(M, N, NRHS, A, LDA, TAU, B, LDB, WORK, LWORK, INFO)
SGERQS
Definition: sgerqs.f:124
subroutine srqt03(M, N, K, AF, C, CC, Q, LDA, TAU, WORK, LWORK, RWORK, RESULT)
SRQT03
Definition: srqt03.f:138
subroutine serrrq(PATH, NUNIT)
SERRRQ
Definition: serrrq.f:57
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
subroutine slatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
SLATMS
Definition: slatms.f:323
subroutine srqt01(M, N, A, AF, Q, R, LDA, TAU, WORK, LWORK, RWORK, RESULT)
SRQT01
Definition: srqt01.f:128
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: