LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cchkq3 ( logical, dimension( * )  DOTYPE,
integer  NM,
integer, dimension( * )  MVAL,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNB,
integer, dimension( * )  NBVAL,
integer, dimension( * )  NXVAL,
real  THRESH,
complex, dimension( * )  A,
complex, dimension( * )  COPYA,
real, dimension( * )  S,
complex, dimension( * )  TAU,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

CCHKQ3

Purpose:
 CCHKQ3 tests CGEQP3.
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]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.
[out]A
          A is COMPLEX array, dimension (MMAX*NMAX)
          where MMAX is the maximum value of M in MVAL and NMAX is the
          maximum value of N in NVAL.
[out]COPYA
          COPYA is COMPLEX array, dimension (MMAX*NMAX)
[out]S
          S is REAL array, dimension
                      (min(MMAX,NMAX))
[out]TAU
          TAU is COMPLEX array, dimension (MMAX)
[out]WORK
          WORK is COMPLEX array, dimension
                      (max(M*max(M,N) + 4*min(M,N) + max(M,N)))
[out]RWORK
          RWORK is REAL array, dimension (4*NMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (2*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 160 of file cchkq3.f.

160 *
161 * -- LAPACK test routine (version 3.4.0) --
162 * -- LAPACK is a software package provided by Univ. of Tennessee, --
163 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
164 * November 2011
165 *
166 * .. Scalar Arguments ..
167  INTEGER nm, nn, nnb, nout
168  REAL thresh
169 * ..
170 * .. Array Arguments ..
171  LOGICAL dotype( * )
172  INTEGER iwork( * ), mval( * ), nbval( * ), nval( * ),
173  $ nxval( * )
174  REAL s( * ), rwork( * )
175  COMPLEX a( * ), copya( * ), tau( * ), work( * )
176 * ..
177 *
178 * =====================================================================
179 *
180 * .. Parameters ..
181  INTEGER ntypes
182  parameter ( ntypes = 6 )
183  INTEGER ntests
184  parameter ( ntests = 3 )
185  REAL one, zero
186  COMPLEX czero
187  parameter ( one = 1.0e+0, zero = 0.0e+0,
188  $ czero = ( 0.0e+0, 0.0e+0 ) )
189 * ..
190 * .. Local Scalars ..
191  CHARACTER*3 path
192  INTEGER i, ihigh, ilow, im, imode, in, inb, info,
193  $ istep, k, lda, lw, lwork, m, mnmin, mode, n,
194  $ nb, nerrs, nfail, nrun, nx
195  REAL eps
196 * ..
197 * .. Local Arrays ..
198  INTEGER iseed( 4 ), iseedy( 4 )
199  REAL result( ntests )
200 * ..
201 * .. External Functions ..
202  REAL cqpt01, cqrt11, cqrt12, slamch
203  EXTERNAL cqpt01, cqrt11, cqrt12, slamch
204 * ..
205 * .. External Subroutines ..
206  EXTERNAL alahd, alasum, cgeqp3, clacpy, claset, clatms,
207  $ icopy, slaord, xlaenv
208 * ..
209 * .. Intrinsic Functions ..
210  INTRINSIC max, min
211 * ..
212 * .. Scalars in Common ..
213  LOGICAL lerr, ok
214  CHARACTER*32 srnamt
215  INTEGER infot, iounit
216 * ..
217 * .. Common blocks ..
218  COMMON / infoc / infot, iounit, ok, lerr
219  COMMON / srnamc / srnamt
220 * ..
221 * .. Data statements ..
222  DATA iseedy / 1988, 1989, 1990, 1991 /
223 * ..
224 * .. Executable Statements ..
225 *
226 * Initialize constants and the random number seed.
227 *
228  path( 1: 1 ) = 'Complex precision'
229  path( 2: 3 ) = 'Q3'
230  nrun = 0
231  nfail = 0
232  nerrs = 0
233  DO 10 i = 1, 4
234  iseed( i ) = iseedy( i )
235  10 CONTINUE
236  eps = slamch( 'Epsilon' )
237  infot = 0
238 *
239  DO 90 im = 1, nm
240 *
241 * Do for each value of M in MVAL.
242 *
243  m = mval( im )
244  lda = max( 1, m )
245 *
246  DO 80 in = 1, nn
247 *
248 * Do for each value of N in NVAL.
249 *
250  n = nval( in )
251  mnmin = min( m, n )
252  lwork = max( 1, m*max( m, n )+4*mnmin+max( m, n ) )
253 *
254  DO 70 imode = 1, ntypes
255  IF( .NOT.dotype( imode ) )
256  $ GO TO 70
257 *
258 * Do for each type of matrix
259 * 1: zero matrix
260 * 2: one small singular value
261 * 3: geometric distribution of singular values
262 * 4: first n/2 columns fixed
263 * 5: last n/2 columns fixed
264 * 6: every second column fixed
265 *
266  mode = imode
267  IF( imode.GT.3 )
268  $ mode = 1
269 *
270 * Generate test matrix of size m by n using
271 * singular value distribution indicated by `mode'.
272 *
273  DO 20 i = 1, n
274  iwork( i ) = 0
275  20 CONTINUE
276  IF( imode.EQ.1 ) THEN
277  CALL claset( 'Full', m, n, czero, czero, copya, lda )
278  DO 30 i = 1, mnmin
279  s( i ) = zero
280  30 CONTINUE
281  ELSE
282  CALL clatms( m, n, 'Uniform', iseed, 'Nonsymm', s,
283  $ mode, one / eps, one, m, n, 'No packing',
284  $ copya, lda, work, info )
285  IF( imode.GE.4 ) THEN
286  IF( imode.EQ.4 ) THEN
287  ilow = 1
288  istep = 1
289  ihigh = max( 1, n / 2 )
290  ELSE IF( imode.EQ.5 ) THEN
291  ilow = max( 1, n / 2 )
292  istep = 1
293  ihigh = n
294  ELSE IF( imode.EQ.6 ) THEN
295  ilow = 1
296  istep = 2
297  ihigh = n
298  END IF
299  DO 40 i = ilow, ihigh, istep
300  iwork( i ) = 1
301  40 CONTINUE
302  END IF
303  CALL slaord( 'Decreasing', mnmin, s, 1 )
304  END IF
305 *
306  DO 60 inb = 1, nnb
307 *
308 * Do for each pair of values (NB,NX) in NBVAL and NXVAL.
309 *
310  nb = nbval( inb )
311  CALL xlaenv( 1, nb )
312  nx = nxval( inb )
313  CALL xlaenv( 3, nx )
314 *
315 * Save A and its singular values and a copy of
316 * vector IWORK.
317 *
318  CALL clacpy( 'All', m, n, copya, lda, a, lda )
319  CALL icopy( n, iwork( 1 ), 1, iwork( n+1 ), 1 )
320 *
321 * Workspace needed.
322 *
323  lw = nb*( n+1 )
324 *
325  srnamt = 'CGEQP3'
326  CALL cgeqp3( m, n, a, lda, iwork( n+1 ), tau, work,
327  $ lw, rwork, info )
328 *
329 * Compute norm(svd(a) - svd(r))
330 *
331  result( 1 ) = cqrt12( m, n, a, lda, s, work,
332  $ lwork, rwork )
333 *
334 * Compute norm( A*P - Q*R )
335 *
336  result( 2 ) = cqpt01( m, n, mnmin, copya, a, lda, tau,
337  $ iwork( n+1 ), work, lwork )
338 *
339 * Compute Q'*Q
340 *
341  result( 3 ) = cqrt11( m, mnmin, a, lda, tau, work,
342  $ lwork )
343 *
344 * Print information about the tests that did not pass
345 * the threshold.
346 *
347  DO 50 k = 1, ntests
348  IF( result( k ).GE.thresh ) THEN
349  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
350  $ CALL alahd( nout, path )
351  WRITE( nout, fmt = 9999 )'CGEQP3', m, n, nb,
352  $ imode, k, result( k )
353  nfail = nfail + 1
354  END IF
355  50 CONTINUE
356  nrun = nrun + ntests
357 *
358  60 CONTINUE
359  70 CONTINUE
360  80 CONTINUE
361  90 CONTINUE
362 *
363 * Print a summary of the results.
364 *
365  CALL alasum( path, nout, nfail, nrun, nerrs )
366 *
367  9999 FORMAT( 1x, a, ' M =', i5, ', N =', i5, ', NB =', i4, ', type ',
368  $ i2, ', test ', i2, ', ratio =', g12.5 )
369 *
370 * End of CCHKQ3
371 *
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine icopy(N, SX, INCX, SY, INCY)
ICOPY
Definition: icopy.f:77
real function cqpt01(M, N, K, A, AF, LDA, TAU, JPVT, WORK, LWORK)
CQPT01
Definition: cqpt01.f:122
real function cqrt11(M, K, A, LDA, TAU, WORK, LWORK)
CQRT11
Definition: cqrt11.f:100
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine claset(UPLO, M, N, ALPHA, BETA, A, LDA)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: claset.f:108
subroutine cgeqp3(M, N, A, LDA, JPVT, TAU, WORK, LWORK, RWORK, INFO)
CGEQP3
Definition: cgeqp3.f:161
subroutine clatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
CLATMS
Definition: clatms.f:334
subroutine slaord(JOB, N, X, INCX)
SLAORD
Definition: slaord.f:75
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
real function cqrt12(M, N, A, LDA, S, WORK, LWORK, RWORK)
CQRT12
Definition: cqrt12.f:99
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
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: