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

◆ cckgsv()

subroutine cckgsv ( integer nm,
integer, dimension( * ) mval,
integer, dimension( * ) pval,
integer, dimension( * ) nval,
integer nmats,
integer, dimension( 4 ) iseed,
real thresh,
integer nmax,
complex, dimension( * ) a,
complex, dimension( * ) af,
complex, dimension( * ) b,
complex, dimension( * ) bf,
complex, dimension( * ) u,
complex, dimension( * ) v,
complex, dimension( * ) q,
real, dimension( * ) alpha,
real, dimension( * ) beta,
complex, dimension( * ) r,
integer, dimension( * ) iwork,
complex, dimension( * ) work,
real, dimension( * ) rwork,
integer nin,
integer nout,
integer info )

CCKGSV

Purpose:
!>
!> CCKGSV tests CGGSVD:
!>        the GSVD for M-by-N matrix A and P-by-N matrix B.
!> 
Parameters
[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]PVAL
!>          PVAL is INTEGER array, dimension (NP)
!>          The values of the matrix row dimension P.
!> 
[in]NVAL
!>          NVAL is INTEGER array, dimension (NN)
!>          The values of the matrix column dimension N.
!> 
[in]NMATS
!>          NMATS is INTEGER
!>          The number of matrix types to be tested for each combination
!>          of matrix dimensions.  If NMATS >= NTYPES (the maximum
!>          number of matrix types), then all the different types are
!>          generated for testing.  If NMATS < NTYPES, another input line
!>          is read to get the numbers of the matrix types to be used.
!> 
[in,out]ISEED
!>          ISEED is INTEGER array, dimension (4)
!>          On entry, the seed of the random number generator.  The array
!>          elements should be between 0 and 4095, otherwise they will be
!>          reduced mod 4096, and ISEED(4) must be odd.
!>          On exit, the next seed in the random number sequence after
!>          all the test matrices have been generated.
!> 
[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]NMAX
!>          NMAX is INTEGER
!>          The maximum value permitted for M or N, used in dimensioning
!>          the work arrays.
!> 
[out]A
!>          A is COMPLEX array, dimension (NMAX*NMAX)
!> 
[out]AF
!>          AF is COMPLEX array, dimension (NMAX*NMAX)
!> 
[out]B
!>          B is COMPLEX array, dimension (NMAX*NMAX)
!> 
[out]BF
!>          BF is COMPLEX array, dimension (NMAX*NMAX)
!> 
[out]U
!>          U is COMPLEX array, dimension (NMAX*NMAX)
!> 
[out]V
!>          V is COMPLEX array, dimension (NMAX*NMAX)
!> 
[out]Q
!>          Q is COMPLEX array, dimension (NMAX*NMAX)
!> 
[out]ALPHA
!>          ALPHA is REAL array, dimension (NMAX)
!> 
[out]BETA
!>          BETA is REAL array, dimension (NMAX)
!> 
[out]R
!>          R is COMPLEX array, dimension (NMAX*NMAX)
!> 
[out]IWORK
!>          IWORK is INTEGER array, dimension (NMAX)
!> 
[out]WORK
!>          WORK is COMPLEX array, dimension (NMAX*NMAX)
!> 
[out]RWORK
!>          RWORK is REAL array, dimension (NMAX)
!> 
[in]NIN
!>          NIN is INTEGER
!>          The unit number for input.
!> 
[in]NOUT
!>          NOUT is INTEGER
!>          The unit number for output.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0 :  successful exit
!>          > 0 :  If CLATMS returns an error code, the absolute value
!>                 of it is returned.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 195 of file cckgsv.f.

198*
199* -- LAPACK test routine --
200* -- LAPACK is a software package provided by Univ. of Tennessee, --
201* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
202*
203* .. Scalar Arguments ..
204 INTEGER INFO, NIN, NM, NMATS, NMAX, NOUT
205 REAL THRESH
206* ..
207* .. Array Arguments ..
208 INTEGER ISEED( 4 ), IWORK( * ), MVAL( * ), NVAL( * ),
209 $ PVAL( * )
210 REAL ALPHA( * ), BETA( * ), RWORK( * )
211 COMPLEX A( * ), AF( * ), B( * ), BF( * ), Q( * ),
212 $ R( * ), U( * ), V( * ), WORK( * )
213* ..
214*
215* =====================================================================
216*
217* .. Parameters ..
218 INTEGER NTESTS
219 parameter( ntests = 12 )
220 INTEGER NTYPES
221 parameter( ntypes = 8 )
222* ..
223* .. Local Scalars ..
224 LOGICAL FIRSTT
225 CHARACTER DISTA, DISTB, TYPE
226 CHARACTER*3 PATH
227 INTEGER I, IINFO, IM, IMAT, KLA, KLB, KUA, KUB, LDA,
228 $ LDB, LDQ, LDR, LDU, LDV, LWORK, M, MODEA,
229 $ MODEB, N, NFAIL, NRUN, NT, P, K, L
230 REAL ANORM, BNORM, CNDNMA, CNDNMB
231* ..
232* .. Local Arrays ..
233 LOGICAL DOTYPE( NTYPES )
234 REAL RESULT( NTESTS )
235* ..
236* .. External Subroutines ..
238* ..
239* .. Intrinsic Functions ..
240 INTRINSIC abs
241* ..
242* .. Executable Statements ..
243*
244* Initialize constants and the random number seed.
245*
246 path( 1: 3 ) = 'GSV'
247 info = 0
248 nrun = 0
249 nfail = 0
250 firstt = .true.
251 CALL alareq( path, nmats, dotype, ntypes, nin, nout )
252 lda = nmax
253 ldb = nmax
254 ldu = nmax
255 ldv = nmax
256 ldq = nmax
257 ldr = nmax
258 lwork = nmax*nmax
259*
260* Specific cases
261*
262* Test: https://github.com/Reference-LAPACK/lapack/issues/411#issue-608776973
263*
264 m = 6
265 p = 6
266 n = 6
267 a(1:m*n) = cmplx(1.e0, 0.e0)
268 b(1:m*n) = cmplx(0.e0, 0.e0)
269 b(1+0*m) = cmplx(9.e19, 0.e0)
270 b(2+1*m) = cmplx(9.e18, 0.e0)
271 b(3+2*m) = cmplx(9.e17, 0.e0)
272 b(4+3*m) = cmplx(9.e16, 0.e0)
273 b(5+4*m) = cmplx(9.e15, 0.e0)
274 b(6+5*m) = cmplx(9.e14, 0.e0)
275 CALL cggsvd3('N','N','N', m, p, n, k, l, a, m, b, m,
276 $ alpha, beta, u, 1, v, 1, q, 1,
277 $ work, m*n, rwork, iwork, info)
278*
279* Print information there is a NAN in BETA
280 DO 40 i = 1, l
281 IF( beta(i).NE.beta(i) ) THEN
282 info = -i
283 EXIT
284 END IF
285 40 CONTINUE
286 IF( info.LT.0 ) THEN
287 IF( nfail.EQ.0 .AND. firstt ) THEN
288 firstt = .false.
289 CALL alahdg( nout, path )
290 END IF
291 WRITE( nout, fmt = 9997 ) -info
292 nfail = nfail + 1
293 END IF
294 nrun = nrun + 1
295 info = 0
296*
297* Do for each value of M in MVAL.
298*
299 DO 30 im = 1, nm
300 m = mval( im )
301 p = pval( im )
302 n = nval( im )
303*
304 DO 20 imat = 1, ntypes
305*
306* Do the tests only if DOTYPE( IMAT ) is true.
307*
308 IF( .NOT.dotype( imat ) )
309 $ GO TO 20
310*
311* Set up parameters with SLATB9 and generate test
312* matrices A and B with CLATMS.
313*
314 CALL slatb9( path, imat, m, p, n, TYPE, KLA, KUA, KLB, KUB,
315 $ ANORM, BNORM, MODEA, MODEB, CNDNMA, CNDNMB,
316 $ DISTA, DISTB )
317*
318* Generate M by N matrix A
319*
320 CALL clatms( m, n, dista, iseed, TYPE, RWORK, MODEA, CNDNMA,
321 $ ANORM, KLA, KUA, 'No packing', A, LDA, WORK,
322 $ IINFO )
323 IF( iinfo.NE.0 ) THEN
324 WRITE( nout, fmt = 9999 )iinfo
325 info = abs( iinfo )
326 GO TO 20
327 END IF
328*
329* Generate P by N matrix B
330*
331 CALL clatms( p, n, distb, iseed, TYPE, RWORK, MODEB, CNDNMB,
332 $ BNORM, KLB, KUB, 'No packing', B, LDB, WORK,
333 $ IINFO )
334 IF( iinfo.NE.0 ) THEN
335 WRITE( nout, fmt = 9999 )iinfo
336 info = abs( iinfo )
337 GO TO 20
338 END IF
339*
340 nt = 6
341*
342 CALL cgsvts3( m, p, n, a, af, lda, b, bf, ldb, u, ldu, v,
343 $ ldv, q, ldq, alpha, beta, r, ldr, iwork, work,
344 $ lwork, rwork, result )
345*
346* Print information about the tests that did not
347* pass the threshold.
348*
349 DO 10 i = 1, nt
350 IF( result( i ).GE.thresh ) THEN
351 IF( nfail.EQ.0 .AND. firstt ) THEN
352 firstt = .false.
353 CALL alahdg( nout, path )
354 END IF
355 WRITE( nout, fmt = 9998 )m, p, n, imat, i,
356 $ result( i )
357 nfail = nfail + 1
358 END IF
359 10 CONTINUE
360 nrun = nrun + nt
361*
362 20 CONTINUE
363 30 CONTINUE
364*
365* Print a summary of the results.
366*
367 CALL alasum( path, nout, nfail, nrun, 0 )
368*
369 9999 FORMAT( ' CLATMS in CCKGSV INFO = ', i5 )
370 9998 FORMAT( ' M=', i4, ' P=', i4, ', N=', i4, ', type ', i2,
371 $ ', test ', i2, ', ratio=', g13.6 )
372 9997 FORMAT( ' FOUND NaN in BETA(', i4,')' )
373 RETURN
374*
375* End of CCKGSV
376*
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 cgsvts3(m, p, n, a, af, lda, b, bf, ldb, u, ldu, v, ldv, q, ldq, alpha, beta, r, ldr, iwork, work, lwork, rwork, result)
CGSVTS3
Definition cgsvts3.f:209
subroutine clatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
CLATMS
Definition clatms.f:332
subroutine cggsvd3(jobu, jobv, jobq, m, n, p, k, l, a, lda, b, ldb, alpha, beta, u, ldu, v, ldv, q, ldq, work, lwork, rwork, iwork, info)
CGGSVD3 computes the singular value decomposition (SVD) for OTHER matrices
Definition cggsvd3.f:352
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
Here is the call graph for this function:
Here is the caller graph for this function: