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

◆ dchklq()

subroutine dchklq ( 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( * ) al,
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 nout )

DCHKLQ

Purpose:
!>
!> DCHKLQ tests DGELQF, DORGLQ and DORMLQ.
!> 
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]AL
!>          AL 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)
!> 
[in]NOUT
!>          NOUT is INTEGER
!>          The unit number for output.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 193 of file dchklq.f.

196*
197* -- LAPACK test routine --
198* -- LAPACK is a software package provided by Univ. of Tennessee, --
199* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
200*
201* .. Scalar Arguments ..
202 LOGICAL TSTERR
203 INTEGER NM, NMAX, NN, NNB, NOUT, NRHS
204 DOUBLE PRECISION THRESH
205* ..
206* .. Array Arguments ..
207 LOGICAL DOTYPE( * )
208 INTEGER MVAL( * ), NBVAL( * ), NVAL( * ),
209 $ NXVAL( * )
210 DOUBLE PRECISION A( * ), AC( * ), AF( * ), AL( * ), AQ( * ),
211 $ B( * ), RWORK( * ), TAU( * ), WORK( * ),
212 $ X( * ), XACT( * )
213* ..
214*
215* =====================================================================
216*
217* .. Parameters ..
218 INTEGER NTESTS
219 parameter( ntests = 7 )
220 INTEGER NTYPES
221 parameter( ntypes = 8 )
222 DOUBLE PRECISION ZERO
223 parameter( zero = 0.0d0 )
224* ..
225* .. Local Scalars ..
226 CHARACTER DIST, TYPE
227 CHARACTER*3 PATH
228 INTEGER I, IK, IM, IMAT, IN, INB, INFO, K, KL, KU, LDA,
229 $ LWORK, M, MINMN, MODE, N, NB, NERRS, NFAIL, NK,
230 $ NRUN, NT, NX
231 DOUBLE PRECISION ANORM, CNDNUM
232* ..
233* .. Local Arrays ..
234 INTEGER ISEED( 4 ), ISEEDY( 4 ), KVAL( 4 )
235 DOUBLE PRECISION RESULT( NTESTS )
236* ..
237* .. External Subroutines ..
238 EXTERNAL alaerh, alahd, alasum, derrlq, dgels, dget02,
240 $ dlqt03, xlaenv
241* ..
242* .. Intrinsic Functions ..
243 INTRINSIC max, min
244* ..
245* .. Scalars in Common ..
246 LOGICAL LERR, OK
247 CHARACTER*32 SRNAMT
248 INTEGER INFOT, NUNIT
249* ..
250* .. Common blocks ..
251 COMMON / infoc / infot, nunit, ok, lerr
252 COMMON / srnamc / srnamt
253* ..
254* .. Data statements ..
255 DATA iseedy / 1988, 1989, 1990, 1991 /
256* ..
257* .. Executable Statements ..
258*
259* Initialize constants and the random number seed.
260*
261 path( 1: 1 ) = 'Double precision'
262 path( 2: 3 ) = 'LQ'
263 nrun = 0
264 nfail = 0
265 nerrs = 0
266 DO 10 i = 1, 4
267 iseed( i ) = iseedy( i )
268 10 CONTINUE
269*
270* Test the error exits
271*
272 IF( tsterr )
273 $ CALL derrlq( path, nout )
274 infot = 0
275 CALL xlaenv( 2, 2 )
276*
277 lda = nmax
278 lwork = nmax*max( nmax, nrhs )
279*
280* Do for each value of M in MVAL.
281*
282 DO 70 im = 1, nm
283 m = mval( im )
284*
285* Do for each value of N in NVAL.
286*
287 DO 60 in = 1, nn
288 n = nval( in )
289 minmn = min( m, n )
290 DO 50 imat = 1, ntypes
291*
292* Do the tests only if DOTYPE( IMAT ) is true.
293*
294 IF( .NOT.dotype( imat ) )
295 $ GO TO 50
296*
297* Set up parameters with DLATB4 and generate a test matrix
298* with DLATMS.
299*
300 CALL dlatb4( path, imat, m, n, TYPE, KL, KU, ANORM, MODE,
301 $ CNDNUM, DIST )
302*
303 srnamt = 'DLATMS'
304 CALL dlatms( m, n, dist, iseed, TYPE, RWORK, MODE,
305 $ CNDNUM, ANORM, KL, KU, 'No packing', A, LDA,
306 $ WORK, INFO )
307*
308* Check error code from DLATMS.
309*
310 IF( info.NE.0 ) THEN
311 CALL alaerh( path, 'DLATMS', info, 0, ' ', m, n, -1,
312 $ -1, -1, imat, nfail, nerrs, nout )
313 GO TO 50
314 END IF
315*
316* Set some values for K: the first value must be MINMN,
317* corresponding to the call of DLQT01; other values are
318* used in the calls of DLQT02, and must not exceed MINMN.
319*
320 kval( 1 ) = minmn
321 kval( 2 ) = 0
322 kval( 3 ) = 1
323 kval( 4 ) = minmn / 2
324 IF( minmn.EQ.0 ) THEN
325 nk = 1
326 ELSE IF( minmn.EQ.1 ) THEN
327 nk = 2
328 ELSE IF( minmn.LE.3 ) THEN
329 nk = 3
330 ELSE
331 nk = 4
332 END IF
333*
334* Do for each value of K in KVAL
335*
336 DO 40 ik = 1, nk
337 k = kval( ik )
338*
339* Do for each pair of values (NB,NX) in NBVAL and NXVAL.
340*
341 DO 30 inb = 1, nnb
342 nb = nbval( inb )
343 CALL xlaenv( 1, nb )
344 nx = nxval( inb )
345 CALL xlaenv( 3, nx )
346 DO i = 1, ntests
347 result( i ) = zero
348 END DO
349 nt = 2
350 IF( ik.EQ.1 ) THEN
351*
352* Test DGELQF
353*
354 CALL dlqt01( m, n, a, af, aq, al, lda, tau,
355 $ work, lwork, rwork, result( 1 ) )
356 ELSE IF( m.LE.n ) THEN
357*
358* Test DORGLQ, using factorization
359* returned by DLQT01
360*
361 CALL dlqt02( m, n, k, a, af, aq, al, lda, tau,
362 $ work, lwork, rwork, result( 1 ) )
363 ELSE
364 result( 1 ) = zero
365 result( 2 ) = zero
366 END IF
367 IF( m.GE.k ) THEN
368*
369* Test DORMLQ, using factorization returned
370* by DLQT01
371*
372 CALL dlqt03( m, n, k, af, ac, al, aq, lda, tau,
373 $ work, lwork, rwork, result( 3 ) )
374 nt = nt + 4
375*
376* If M<=N and K=M, call DGELS to solve a system
377* with NRHS right hand sides and compute the
378* residual.
379*
380 IF( k.EQ.m .AND. inb.EQ.1 ) THEN
381*
382* Generate a solution and set the right
383* hand side.
384*
385 srnamt = 'DLARHS'
386 CALL dlarhs( path, 'New', 'Full',
387 $ 'No transpose', m, n, 0, 0,
388 $ nrhs, a, lda, xact, lda, b, lda,
389 $ iseed, info )
390*
391 CALL dlacpy( 'Full', m, nrhs, b, lda, x,
392 $ lda )
393*
394* Reset AF to the original matrix. DGELS
395* factors the matrix before solving the system.
396*
397 CALL dlacpy( 'Full', m, n, a, lda, af, lda )
398*
399 srnamt = 'DGELS'
400 CALL dgels( 'No transpose', m, n, nrhs, af,
401 $ lda, x, lda, work, lwork, info )
402*
403* Check error code from DGELS.
404*
405 IF( info.NE.0 )
406 $ CALL alaerh( path, 'DGELS', info, 0, 'N',
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 ELSE
415 result( 7 ) = zero
416 END IF
417 ELSE
418 result( 3 ) = zero
419 result( 4 ) = zero
420 result( 5 ) = zero
421 result( 6 ) = zero
422 END IF
423*
424* Print information about the tests that did not
425* pass the threshold.
426*
427 DO 20 i = 1, nt
428 IF( result( i ).GE.thresh ) THEN
429 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
430 $ CALL alahd( nout, path )
431 WRITE( nout, fmt = 9999 )m, n, k, nb, nx,
432 $ imat, i, result( i )
433 nfail = nfail + 1
434 END IF
435 20 CONTINUE
436 nrun = nrun + nt
437 30 CONTINUE
438 40 CONTINUE
439 50 CONTINUE
440 60 CONTINUE
441 70 CONTINUE
442*
443* Print a summary of the results.
444*
445 CALL alasum( path, nout, nfail, nrun, nerrs )
446*
447 9999 FORMAT( ' M=', i5, ', N=', i5, ', K=', i5, ', NB=', i4, ', NX=',
448 $ i5, ', type ', i2, ', test(', i2, ')=', g12.5 )
449 RETURN
450*
451* End of DCHKLQ
452*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine dget02(trans, m, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
DGET02
Definition dget02.f:135
subroutine dlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
DLARHS
Definition dlarhs.f:205
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
subroutine derrlq(path, nunit)
DERRLQ
Definition derrlq.f:55
subroutine dlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
DLATB4
Definition dlatb4.f:120
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
Definition dlatms.f:321
subroutine dlqt01(m, n, a, af, q, l, lda, tau, work, lwork, rwork, result)
DLQT01
Definition dlqt01.f:126
subroutine dlqt02(m, n, k, a, af, q, l, lda, tau, work, lwork, rwork, result)
DLQT02
Definition dlqt02.f:135
subroutine dlqt03(m, n, k, af, c, cc, q, lda, tau, work, lwork, rwork, result)
DLQT03
Definition dlqt03.f:136
subroutine dgels(trans, m, n, nrhs, a, lda, b, ldb, work, lwork, info)
DGELS solves overdetermined or underdetermined systems for GE matrices
Definition dgels.f:192
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
Here is the call graph for this function:
Here is the caller graph for this function: