LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sckgqr.f
Go to the documentation of this file.
1*> \brief \b SCKGQR
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 SCKGQR( NM, MVAL, NP, PVAL, NN, NVAL, NMATS, ISEED,
12* THRESH, NMAX, A, AF, AQ, AR, TAUA, B, BF, BZ,
13* BT, BWK, TAUB, WORK, RWORK, NIN, NOUT, INFO )
14*
15* .. Scalar Arguments ..
16* INTEGER INFO, NIN, NM, NMATS, NMAX, NN, NOUT, NP
17* REAL THRESH
18* ..
19* .. Array Arguments ..
20* INTEGER ISEED( 4 ), MVAL( * ), NVAL( * ), PVAL( * )
21* REAL A( * ), AF( * ), AQ( * ), AR( * ), B( * ),
22* $ BF( * ), BT( * ), BWK( * ), BZ( * ),
23* $ RWORK( * ), TAUA( * ), TAUB( * ), WORK( * )
24* ..
25*
26*
27*> \par Purpose:
28* =============
29*>
30*> \verbatim
31*>
32*> SCKGQR tests
33*> SGGQRF: GQR factorization for N-by-M matrix A and N-by-P matrix B,
34*> SGGRQF: GRQ factorization for M-by-N matrix A and P-by-N matrix B.
35*> \endverbatim
36*
37* Arguments:
38* ==========
39*
40*> \param[in] NM
41*> \verbatim
42*> NM is INTEGER
43*> The number of values of M contained in the vector MVAL.
44*> \endverbatim
45*>
46*> \param[in] MVAL
47*> \verbatim
48*> MVAL is INTEGER array, dimension (NM)
49*> The values of the matrix row(column) dimension M.
50*> \endverbatim
51*>
52*> \param[in] NP
53*> \verbatim
54*> NP is INTEGER
55*> The number of values of P contained in the vector PVAL.
56*> \endverbatim
57*>
58*> \param[in] PVAL
59*> \verbatim
60*> PVAL is INTEGER array, dimension (NP)
61*> The values of the matrix row(column) dimension P.
62*> \endverbatim
63*>
64*> \param[in] NN
65*> \verbatim
66*> NN is INTEGER
67*> The number of values of N contained in the vector NVAL.
68*> \endverbatim
69*>
70*> \param[in] NVAL
71*> \verbatim
72*> NVAL is INTEGER array, dimension (NN)
73*> The values of the matrix column(row) dimension N.
74*> \endverbatim
75*>
76*> \param[in] NMATS
77*> \verbatim
78*> NMATS is INTEGER
79*> The number of matrix types to be tested for each combination
80*> of matrix dimensions. If NMATS >= NTYPES (the maximum
81*> number of matrix types), then all the different types are
82*> generated for testing. If NMATS < NTYPES, another input line
83*> is read to get the numbers of the matrix types to be used.
84*> \endverbatim
85*>
86*> \param[in,out] ISEED
87*> \verbatim
88*> ISEED is INTEGER array, dimension (4)
89*> On entry, the seed of the random number generator. The array
90*> elements should be between 0 and 4095, otherwise they will be
91*> reduced mod 4096, and ISEED(4) must be odd.
92*> On exit, the next seed in the random number sequence after
93*> all the test matrices have been generated.
94*> \endverbatim
95*>
96*> \param[in] THRESH
97*> \verbatim
98*> THRESH is REAL
99*> The threshold value for the test ratios. A result is
100*> included in the output file if RESULT >= THRESH. To have
101*> every test ratio printed, use THRESH = 0.
102*> \endverbatim
103*>
104*> \param[in] NMAX
105*> \verbatim
106*> NMAX is INTEGER
107*> The maximum value permitted for M or N, used in dimensioning
108*> the work arrays.
109*> \endverbatim
110*>
111*> \param[out] A
112*> \verbatim
113*> A is REAL array, dimension (NMAX*NMAX)
114*> \endverbatim
115*>
116*> \param[out] AF
117*> \verbatim
118*> AF is REAL array, dimension (NMAX*NMAX)
119*> \endverbatim
120*>
121*> \param[out] AQ
122*> \verbatim
123*> AQ is REAL array, dimension (NMAX*NMAX)
124*> \endverbatim
125*>
126*> \param[out] AR
127*> \verbatim
128*> AR is REAL array, dimension (NMAX*NMAX)
129*> \endverbatim
130*>
131*> \param[out] TAUA
132*> \verbatim
133*> TAUA is REAL array, dimension (NMAX)
134*> \endverbatim
135*>
136*> \param[out] B
137*> \verbatim
138*> B is REAL array, dimension (NMAX*NMAX)
139*> \endverbatim
140*>
141*> \param[out] BF
142*> \verbatim
143*> BF is REAL array, dimension (NMAX*NMAX)
144*> \endverbatim
145*>
146*> \param[out] BZ
147*> \verbatim
148*> BZ is REAL array, dimension (NMAX*NMAX)
149*> \endverbatim
150*>
151*> \param[out] BT
152*> \verbatim
153*> BT is REAL array, dimension (NMAX*NMAX)
154*> \endverbatim
155*>
156*> \param[out] BWK
157*> \verbatim
158*> BWK is REAL array, dimension (NMAX*NMAX)
159*> \endverbatim
160*>
161*> \param[out] TAUB
162*> \verbatim
163*> TAUB is REAL array, dimension (NMAX)
164*> \endverbatim
165*>
166*> \param[out] WORK
167*> \verbatim
168*> WORK is REAL array, dimension (NMAX*NMAX)
169*> \endverbatim
170*>
171*> \param[out] RWORK
172*> \verbatim
173*> RWORK is REAL array, dimension (NMAX)
174*> \endverbatim
175*>
176*> \param[in] NIN
177*> \verbatim
178*> NIN is INTEGER
179*> The unit number for input.
180*> \endverbatim
181*>
182*> \param[in] NOUT
183*> \verbatim
184*> NOUT is INTEGER
185*> The unit number for output.
186*> \endverbatim
187*>
188*> \param[out] INFO
189*> \verbatim
190*> INFO is INTEGER
191*> = 0 : successful exit
192*> > 0 : If SLATMS returns an error code, the absolute value
193*> of it is returned.
194*> \endverbatim
195*
196* Authors:
197* ========
198*
199*> \author Univ. of Tennessee
200*> \author Univ. of California Berkeley
201*> \author Univ. of Colorado Denver
202*> \author NAG Ltd.
203*
204*> \ingroup single_eig
205*
206* =====================================================================
207 SUBROUTINE sckgqr( NM, MVAL, NP, PVAL, NN, NVAL, NMATS, ISEED,
208 $ THRESH, NMAX, A, AF, AQ, AR, TAUA, B, BF, BZ,
209 $ BT, BWK, TAUB, WORK, RWORK, NIN, NOUT, INFO )
210*
211* -- LAPACK test routine --
212* -- LAPACK is a software package provided by Univ. of Tennessee, --
213* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
214*
215* .. Scalar Arguments ..
216 INTEGER INFO, NIN, NM, NMATS, NMAX, NN, NOUT, NP
217 REAL THRESH
218* ..
219* .. Array Arguments ..
220 INTEGER ISEED( 4 ), MVAL( * ), NVAL( * ), PVAL( * )
221 REAL A( * ), AF( * ), AQ( * ), AR( * ), B( * ),
222 $ bf( * ), bt( * ), bwk( * ), bz( * ),
223 $ rwork( * ), taua( * ), taub( * ), work( * )
224* ..
225*
226* =====================================================================
227*
228* .. Parameters ..
229 INTEGER NTESTS
230 PARAMETER ( NTESTS = 7 )
231 INTEGER NTYPES
232 parameter( ntypes = 8 )
233* ..
234* .. Local Scalars ..
235 LOGICAL FIRSTT
236 CHARACTER DISTA, DISTB, TYPE
237 CHARACTER*3 PATH
238 INTEGER I, IINFO, IM, IMAT, IN, IP, KLA, KLB, KUA, KUB,
239 $ lda, ldb, lwork, m, modea, modeb, n, nfail,
240 $ nrun, nt, p
241 REAL ANORM, BNORM, CNDNMA, CNDNMB
242* ..
243* .. Local Arrays ..
244 LOGICAL DOTYPE( NTYPES )
245 REAL RESULT( NTESTS )
246* ..
247* .. External Subroutines ..
248 EXTERNAL alahdg, alareq, alasum, sgqrts, sgrqts, slatb9,
249 $ slatms
250* ..
251* .. Intrinsic Functions ..
252 INTRINSIC abs
253* ..
254* .. Executable Statements ..
255*
256* Initialize constants.
257*
258 path( 1: 3 ) = 'GQR'
259 info = 0
260 nrun = 0
261 nfail = 0
262 firstt = .true.
263 CALL alareq( path, nmats, dotype, ntypes, nin, nout )
264 lda = nmax
265 ldb = nmax
266 lwork = nmax*nmax
267*
268* Do for each value of M in MVAL.
269*
270 DO 60 im = 1, nm
271 m = mval( im )
272*
273* Do for each value of P in PVAL.
274*
275 DO 50 ip = 1, np
276 p = pval( ip )
277*
278* Do for each value of N in NVAL.
279*
280 DO 40 in = 1, nn
281 n = nval( in )
282*
283 DO 30 imat = 1, ntypes
284*
285* Do the tests only if DOTYPE( IMAT ) is true.
286*
287 IF( .NOT.dotype( imat ) )
288 $ GO TO 30
289*
290* Test SGGRQF
291*
292* Set up parameters with SLATB9 and generate test
293* matrices A and B with SLATMS.
294*
295 CALL slatb9( 'GRQ', imat, m, p, n, TYPE, kla, kua,
296 $ klb, kub, anorm, bnorm, modea, modeb,
297 $ cndnma, cndnmb, dista, distb )
298*
299* Generate M by N matrix A
300*
301 CALL slatms( m, n, dista, iseed, TYPE, rwork, modea,
302 $ cndnma, anorm, kla, kua, 'No packing', a,
303 $ lda, work, iinfo )
304 IF( iinfo.NE.0 ) THEN
305 WRITE( nout, fmt = 9999 )iinfo
306 info = abs( iinfo )
307 GO TO 30
308 END IF
309*
310* Generate P by N matrix B
311*
312 CALL slatms( p, n, distb, iseed, TYPE, rwork, modeb,
313 $ cndnmb, bnorm, klb, kub, 'No packing', b,
314 $ ldb, work, iinfo )
315 IF( iinfo.NE.0 ) THEN
316 WRITE( nout, fmt = 9999 )iinfo
317 info = abs( iinfo )
318 GO TO 30
319 END IF
320*
321 nt = 4
322*
323 CALL sgrqts( m, p, n, a, af, aq, ar, lda, taua, b, bf,
324 $ bz, bt, bwk, ldb, taub, work, lwork,
325 $ rwork, result )
326*
327* Print information about the tests that did not
328* pass the threshold.
329*
330 DO 10 i = 1, nt
331 IF( result( i ).GE.thresh ) THEN
332 IF( nfail.EQ.0 .AND. firstt ) THEN
333 firstt = .false.
334 CALL alahdg( nout, 'GRQ' )
335 END IF
336 WRITE( nout, fmt = 9998 )m, p, n, imat, i,
337 $ result( i )
338 nfail = nfail + 1
339 END IF
340 10 CONTINUE
341 nrun = nrun + nt
342*
343* Test SGGQRF
344*
345* Set up parameters with SLATB9 and generate test
346* matrices A and B with SLATMS.
347*
348 CALL slatb9( 'GQR', imat, m, p, n, TYPE, kla, kua,
349 $ klb, kub, anorm, bnorm, modea, modeb,
350 $ cndnma, cndnmb, dista, distb )
351*
352* Generate N-by-M matrix A
353*
354 CALL slatms( n, m, dista, iseed, TYPE, rwork, modea,
355 $ cndnma, anorm, kla, kua, 'No packing', a,
356 $ lda, work, iinfo )
357 IF( iinfo.NE.0 ) THEN
358 WRITE( nout, fmt = 9999 )iinfo
359 info = abs( iinfo )
360 GO TO 30
361 END IF
362*
363* Generate N-by-P matrix B
364*
365 CALL slatms( n, p, distb, iseed, TYPE, rwork, modea,
366 $ cndnma, bnorm, klb, kub, 'No packing', b,
367 $ ldb, work, iinfo )
368 IF( iinfo.NE.0 ) THEN
369 WRITE( nout, fmt = 9999 )iinfo
370 info = abs( iinfo )
371 GO TO 30
372 END IF
373*
374 nt = 4
375*
376 CALL sgqrts( n, m, p, a, af, aq, ar, lda, taua, b, bf,
377 $ bz, bt, bwk, ldb, taub, work, lwork,
378 $ rwork, result )
379*
380* Print information about the tests that did not
381* pass the threshold.
382*
383 DO 20 i = 1, nt
384 IF( result( i ).GE.thresh ) THEN
385 IF( nfail.EQ.0 .AND. firstt ) THEN
386 firstt = .false.
387 CALL alahdg( nout, path )
388 END IF
389 WRITE( nout, fmt = 9997 )n, m, p, imat, i,
390 $ result( i )
391 nfail = nfail + 1
392 END IF
393 20 CONTINUE
394 nrun = nrun + nt
395*
396 30 CONTINUE
397 40 CONTINUE
398 50 CONTINUE
399 60 CONTINUE
400*
401* Print a summary of the results.
402*
403 CALL alasum( path, nout, nfail, nrun, 0 )
404*
405 9999 FORMAT( ' SLATMS in SCKGQR: INFO = ', i5 )
406 9998 FORMAT( ' M=', i4, ' P=', i4, ', N=', i4, ', type ', i2,
407 $ ', test ', i2, ', ratio=', g13.6 )
408 9997 FORMAT( ' N=', i4, ' M=', i4, ', P=', i4, ', type ', i2,
409 $ ', test ', i2, ', ratio=', g13.6 )
410 RETURN
411*
412* End of SCKGQR
413*
414 END
subroutine alareq(path, nmats, dotype, ntypes, nin, nout)
ALAREQ
Definition alareq.f:90
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine alahdg(iounit, path)
ALAHDG
Definition alahdg.f:62
subroutine sckgqr(nm, mval, np, pval, nn, nval, nmats, iseed, thresh, nmax, a, af, aq, ar, taua, b, bf, bz, bt, bwk, taub, work, rwork, nin, nout, info)
SCKGQR
Definition sckgqr.f:210
subroutine sgqrts(n, m, p, a, af, q, r, lda, taua, b, bf, z, t, bwk, ldb, taub, work, lwork, rwork, result)
SGQRTS
Definition sgqrts.f:176
subroutine sgrqts(m, p, n, a, af, q, r, lda, taua, b, bf, z, t, bwk, ldb, taub, work, lwork, rwork, result)
SGRQTS
Definition sgrqts.f:177
subroutine slatb9(path, imat, m, p, n, type, kla, kua, klb, kub, anorm, bnorm, modea, modeb, cndnma, cndnmb, dista, distb)
SLATB9
Definition slatb9.f:170
subroutine slatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
SLATMS
Definition slatms.f:321