LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dchkq3.f
Go to the documentation of this file.
1*> \brief \b DCHKQ3
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE DCHKQ3( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL,
12* THRESH, A, COPYA, S, TAU, WORK, IWORK,
13* NOUT )
14*
15* .. Scalar Arguments ..
16* INTEGER NM, NN, NNB, NOUT
17* DOUBLE PRECISION THRESH
18* ..
19* .. Array Arguments ..
20* LOGICAL DOTYPE( * )
21* INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NVAL( * ),
22* $ NXVAL( * )
23* DOUBLE PRECISION A( * ), COPYA( * ), S( * ),
24* $ TAU( * ), WORK( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*>
33*> DCHKQ3 tests DGEQP3.
34*> \endverbatim
35*
36* Arguments:
37* ==========
38*
39*> \param[in] DOTYPE
40*> \verbatim
41*> DOTYPE is LOGICAL array, dimension (NTYPES)
42*> The matrix types to be used for testing. Matrices of type j
43*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
44*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
45*> \endverbatim
46*>
47*> \param[in] NM
48*> \verbatim
49*> NM is INTEGER
50*> The number of values of M contained in the vector MVAL.
51*> \endverbatim
52*>
53*> \param[in] MVAL
54*> \verbatim
55*> MVAL is INTEGER array, dimension (NM)
56*> The values of the matrix row dimension M.
57*> \endverbatim
58*>
59*> \param[in] NN
60*> \verbatim
61*> NN is INTEGER
62*> The number of values of N contained in the vector NVAL.
63*> \endverbatim
64*>
65*> \param[in] NVAL
66*> \verbatim
67*> NVAL is INTEGER array, dimension (NN)
68*> The values of the matrix column dimension N.
69*> \endverbatim
70*>
71*> \param[in] NNB
72*> \verbatim
73*> NNB is INTEGER
74*> The number of values of NB and NX contained in the
75*> vectors NBVAL and NXVAL. The blocking parameters are used
76*> in pairs (NB,NX).
77*> \endverbatim
78*>
79*> \param[in] NBVAL
80*> \verbatim
81*> NBVAL is INTEGER array, dimension (NNB)
82*> The values of the blocksize NB.
83*> \endverbatim
84*>
85*> \param[in] NXVAL
86*> \verbatim
87*> NXVAL is INTEGER array, dimension (NNB)
88*> The values of the crossover point NX.
89*> \endverbatim
90*>
91*> \param[in] THRESH
92*> \verbatim
93*> THRESH is DOUBLE PRECISION
94*> The threshold value for the test ratios. A result is
95*> included in the output file if RESULT >= THRESH. To have
96*> every test ratio printed, use THRESH = 0.
97*> \endverbatim
98*>
99*> \param[out] A
100*> \verbatim
101*> A is DOUBLE PRECISION array, dimension (MMAX*NMAX)
102*> where MMAX is the maximum value of M in MVAL and NMAX is the
103*> maximum value of N in NVAL.
104*> \endverbatim
105*>
106*> \param[out] COPYA
107*> \verbatim
108*> COPYA is DOUBLE PRECISION array, dimension (MMAX*NMAX)
109*> \endverbatim
110*>
111*> \param[out] S
112*> \verbatim
113*> S is DOUBLE PRECISION array, dimension
114*> (min(MMAX,NMAX))
115*> \endverbatim
116*>
117*> \param[out] TAU
118*> \verbatim
119*> TAU is DOUBLE PRECISION array, dimension (MMAX)
120*> \endverbatim
121*>
122*> \param[out] WORK
123*> \verbatim
124*> WORK is DOUBLE PRECISION array, dimension
125*> (MMAX*NMAX + 4*NMAX + MMAX)
126*> \endverbatim
127*>
128*> \param[out] IWORK
129*> \verbatim
130*> IWORK is INTEGER array, dimension (2*NMAX)
131*> \endverbatim
132*>
133*> \param[in] NOUT
134*> \verbatim
135*> NOUT is INTEGER
136*> The unit number for output.
137*> \endverbatim
138*
139* Authors:
140* ========
141*
142*> \author Univ. of Tennessee
143*> \author Univ. of California Berkeley
144*> \author Univ. of Colorado Denver
145*> \author NAG Ltd.
146*
147*> \ingroup double_lin
148*
149* =====================================================================
150 SUBROUTINE dchkq3( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL,
151 $ THRESH, A, COPYA, S, TAU, WORK, IWORK,
152 $ NOUT )
153*
154* -- LAPACK test routine --
155* -- LAPACK is a software package provided by Univ. of Tennessee, --
156* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
157*
158* .. Scalar Arguments ..
159 INTEGER NM, NN, NNB, NOUT
160 DOUBLE PRECISION THRESH
161* ..
162* .. Array Arguments ..
163 LOGICAL DOTYPE( * )
164 INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NVAL( * ),
165 $ nxval( * )
166 DOUBLE PRECISION A( * ), COPYA( * ), S( * ),
167 $ TAU( * ), WORK( * )
168* ..
169*
170* =====================================================================
171*
172* .. Parameters ..
173 INTEGER NTYPES
174 PARAMETER ( NTYPES = 6 )
175 INTEGER NTESTS
176 parameter( ntests = 3 )
177 DOUBLE PRECISION ONE, ZERO
178 parameter( one = 1.0d0, zero = 0.0d0 )
179* ..
180* .. Local Scalars ..
181 CHARACTER*3 PATH
182 INTEGER I, IHIGH, ILOW, IM, IMODE, IN, INB, INFO,
183 $ istep, k, lda, lw, lwork, m, mnmin, mode, n,
184 $ nb, nerrs, nfail, nrun, nx
185 DOUBLE PRECISION EPS
186* ..
187* .. Local Arrays ..
188 INTEGER ISEED( 4 ), ISEEDY( 4 )
189 DOUBLE PRECISION RESULT( NTESTS )
190* ..
191* .. External Functions ..
192 DOUBLE PRECISION DLAMCH, DQPT01, DQRT11, DQRT12
193 EXTERNAL DLAMCH, DQPT01, DQRT11, DQRT12
194* ..
195* .. External Subroutines ..
196 EXTERNAL alahd, alasum, dgeqp3, dlacpy, dlaord, dlaset,
198* ..
199* .. Intrinsic Functions ..
200 INTRINSIC max, min
201* ..
202* .. Scalars in Common ..
203 LOGICAL LERR, OK
204 CHARACTER*32 SRNAMT
205 INTEGER INFOT, IOUNIT
206* ..
207* .. Common blocks ..
208 COMMON / infoc / infot, iounit, ok, lerr
209 COMMON / srnamc / srnamt
210* ..
211* .. Data statements ..
212 DATA iseedy / 1988, 1989, 1990, 1991 /
213* ..
214* .. Executable Statements ..
215*
216* Initialize constants and the random number seed.
217*
218 path( 1: 1 ) = 'Double precision'
219 path( 2: 3 ) = 'Q3'
220 nrun = 0
221 nfail = 0
222 nerrs = 0
223 DO 10 i = 1, 4
224 iseed( i ) = iseedy( i )
225 10 CONTINUE
226 eps = dlamch( 'Epsilon' )
227 infot = 0
228*
229 DO 90 im = 1, nm
230*
231* Do for each value of M in MVAL.
232*
233 m = mval( im )
234 lda = max( 1, m )
235*
236 DO 80 in = 1, nn
237*
238* Do for each value of N in NVAL.
239*
240 n = nval( in )
241 mnmin = min( m, n )
242 lwork = max( 1, m*max( m, n )+4*mnmin+max( m, n ),
243 $ m*n + 2*mnmin + 4*n )
244*
245 DO 70 imode = 1, ntypes
246 IF( .NOT.dotype( imode ) )
247 $ GO TO 70
248*
249* Do for each type of matrix
250* 1: zero matrix
251* 2: one small singular value
252* 3: geometric distribution of singular values
253* 4: first n/2 columns fixed
254* 5: last n/2 columns fixed
255* 6: every second column fixed
256*
257 mode = imode
258 IF( imode.GT.3 )
259 $ mode = 1
260*
261* Generate test matrix of size m by n using
262* singular value distribution indicated by `mode'.
263*
264 DO 20 i = 1, n
265 iwork( i ) = 0
266 20 CONTINUE
267 IF( imode.EQ.1 ) THEN
268 CALL dlaset( 'Full', m, n, zero, zero, copya, lda )
269 DO 30 i = 1, mnmin
270 s( i ) = zero
271 30 CONTINUE
272 ELSE
273 CALL dlatms( m, n, 'Uniform', iseed, 'Nonsymm', s,
274 $ mode, one / eps, one, m, n, 'No packing',
275 $ copya, lda, work, info )
276 IF( imode.GE.4 ) THEN
277 IF( imode.EQ.4 ) THEN
278 ilow = 1
279 istep = 1
280 ihigh = max( 1, n / 2 )
281 ELSE IF( imode.EQ.5 ) THEN
282 ilow = max( 1, n / 2 )
283 istep = 1
284 ihigh = n
285 ELSE IF( imode.EQ.6 ) THEN
286 ilow = 1
287 istep = 2
288 ihigh = n
289 END IF
290 DO 40 i = ilow, ihigh, istep
291 iwork( i ) = 1
292 40 CONTINUE
293 END IF
294 CALL dlaord( 'Decreasing', mnmin, s, 1 )
295 END IF
296*
297 DO 60 inb = 1, nnb
298*
299* Do for each pair of values (NB,NX) in NBVAL and NXVAL.
300*
301 nb = nbval( inb )
302 CALL xlaenv( 1, nb )
303 nx = nxval( inb )
304 CALL xlaenv( 3, nx )
305*
306* Get a working copy of COPYA into A and a copy of
307* vector IWORK.
308*
309 CALL dlacpy( 'All', m, n, copya, lda, a, lda )
310 CALL icopy( n, iwork( 1 ), 1, iwork( n+1 ), 1 )
311*
312* Compute the QR factorization with pivoting of A
313*
314 lw = max( 1, 2*n+nb*( n+1 ) )
315*
316* Compute the QP3 factorization of A
317*
318 srnamt = 'DGEQP3'
319 CALL dgeqp3( m, n, a, lda, iwork( n+1 ), tau, work,
320 $ lw, info )
321*
322* Compute norm(svd(a) - svd(r))
323*
324 result( 1 ) = dqrt12( m, n, a, lda, s, work,
325 $ lwork )
326*
327* Compute norm( A*P - Q*R )
328*
329 result( 2 ) = dqpt01( m, n, mnmin, copya, a, lda, tau,
330 $ iwork( n+1 ), work, lwork )
331*
332* Compute Q'*Q
333*
334 result( 3 ) = dqrt11( m, mnmin, a, lda, tau, work,
335 $ lwork )
336*
337* Print information about the tests that did not pass
338* the threshold.
339*
340 DO 50 k = 1, ntests
341 IF( result( k ).GE.thresh ) THEN
342 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
343 $ CALL alahd( nout, path )
344 WRITE( nout, fmt = 9999 )'DGEQP3', m, n, nb,
345 $ imode, k, result( k )
346 nfail = nfail + 1
347 END IF
348 50 CONTINUE
349 nrun = nrun + ntests
350*
351 60 CONTINUE
352 70 CONTINUE
353 80 CONTINUE
354 90 CONTINUE
355*
356* Print a summary of the results.
357*
358 CALL alasum( path, nout, nfail, nrun, nerrs )
359*
360 9999 FORMAT( 1x, a, ' M =', i5, ', N =', i5, ', NB =', i4, ', type ',
361 $ i2, ', test ', i2, ', ratio =', g12.5 )
362*
363* End of DCHKQ3
364*
365 END
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
subroutine dchkq3(dotype, nm, mval, nn, nval, nnb, nbval, nxval, thresh, a, copya, s, tau, work, iwork, nout)
DCHKQ3
Definition dchkq3.f:153
subroutine dlaord(job, n, x, incx)
DLAORD
Definition dlaord.f:73
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
Definition dlatms.f:321
subroutine dgeqp3(m, n, a, lda, jpvt, tau, work, lwork, info)
DGEQP3
Definition dgeqp3.f:151
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
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:110
subroutine icopy(n, sx, incx, sy, incy)
ICOPY
Definition icopy.f:75