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

◆ zdrvab()

subroutine zdrvab ( logical, dimension( * )  dotype,
integer  nm,
integer, dimension( * )  mval,
integer  nns,
integer, dimension( * )  nsval,
double precision  thresh,
integer  nmax,
complex*16, dimension( * )  a,
complex*16, dimension( * )  afac,
complex*16, dimension( * )  b,
complex*16, dimension( * )  x,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork,
complex, dimension( * )  swork,
integer, dimension( * )  iwork,
integer  nout 
)

ZDRVAB

Purpose:
 ZDRVAB tests ZCGESV
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]NNS
          NNS is INTEGER
          The number of values of NRHS contained in the vector NSVAL.
[in]NSVAL
          NSVAL is INTEGER array, dimension (NNS)
          The values of the number of right hand sides NRHS.
[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]NMAX
          NMAX is INTEGER
          The maximum value permitted for M or N, used in dimensioning
          the work arrays.
[out]A
          A is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(3,NSMAX*2))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      NMAX
[out]SWORK
          SWORK is COMPLEX array, dimension
                      (NMAX*(NSMAX+NMAX))
[out]IWORK
          IWORK is INTEGER 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 149 of file zdrvab.f.

152*
153* -- LAPACK test routine --
154* -- LAPACK is a software package provided by Univ. of Tennessee, --
155* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
156*
157* .. Scalar Arguments ..
158 INTEGER NM, NMAX, NNS, NOUT
159 DOUBLE PRECISION THRESH
160* ..
161* .. Array Arguments ..
162 LOGICAL DOTYPE( * )
163 INTEGER MVAL( * ), NSVAL( * ), IWORK( * )
164 DOUBLE PRECISION RWORK( * )
165 COMPLEX SWORK( * )
166 COMPLEX*16 A( * ), AFAC( * ), B( * ),
167 $ WORK( * ), X( * )
168* ..
169*
170* =====================================================================
171*
172* .. Parameters ..
173 DOUBLE PRECISION ZERO
174 parameter( zero = 0.0d+0 )
175 INTEGER NTYPES
176 parameter( ntypes = 11 )
177 INTEGER NTESTS
178 parameter( ntests = 1 )
179* ..
180* .. Local Scalars ..
181 LOGICAL ZEROT
182 CHARACTER DIST, TRANS, TYPE, XTYPE
183 CHARACTER*3 PATH
184 INTEGER I, IM, IMAT, INFO, IOFF, IRHS,
185 $ IZERO, KL, KU, LDA, M, MODE, N,
186 $ NERRS, NFAIL, NIMAT, NRHS, NRUN
187 DOUBLE PRECISION ANORM, CNDNUM
188* ..
189* .. Local Arrays ..
190 INTEGER ISEED( 4 ), ISEEDY( 4 )
191 DOUBLE PRECISION RESULT( NTESTS )
192* ..
193* .. Local Variables ..
194 INTEGER ITER, KASE
195* ..
196* .. External Subroutines ..
197 EXTERNAL alaerh, alahd, zget08, zlacpy, zlarhs, zlaset,
198 $ zlatb4, zlatms
199* ..
200* .. Intrinsic Functions ..
201 INTRINSIC dcmplx, dble, max, min, sqrt
202* ..
203* .. Scalars in Common ..
204 LOGICAL LERR, OK
205 CHARACTER*32 SRNAMT
206 INTEGER INFOT, NUNIT
207* ..
208* .. Common blocks ..
209 COMMON / infoc / infot, nunit, ok, lerr
210 COMMON / srnamc / srnamt
211* ..
212* .. Data statements ..
213 DATA iseedy / 2006, 2007, 2008, 2009 /
214* ..
215* .. Executable Statements ..
216*
217* Initialize constants and the random number seed.
218*
219 kase = 0
220 path( 1: 1 ) = 'Zomplex precision'
221 path( 2: 3 ) = 'GE'
222 nrun = 0
223 nfail = 0
224 nerrs = 0
225 DO 10 i = 1, 4
226 iseed( i ) = iseedy( i )
227 10 CONTINUE
228*
229 infot = 0
230*
231* Do for each value of M in MVAL
232*
233 DO 120 im = 1, nm
234 m = mval( im )
235 lda = max( 1, m )
236*
237 n = m
238 nimat = ntypes
239 IF( m.LE.0 .OR. n.LE.0 )
240 $ nimat = 1
241*
242 DO 100 imat = 1, nimat
243*
244* Do the tests only if DOTYPE( IMAT ) is true.
245*
246 IF( .NOT.dotype( imat ) )
247 $ GO TO 100
248*
249* Skip types 5, 6, or 7 if the matrix size is too small.
250*
251 zerot = imat.GE.5 .AND. imat.LE.7
252 IF( zerot .AND. n.LT.imat-4 )
253 $ GO TO 100
254*
255* Set up parameters with ZLATB4 and generate a test matrix
256* with ZLATMS.
257*
258 CALL zlatb4( path, imat, m, n, TYPE, KL, KU, ANORM, MODE,
259 $ CNDNUM, DIST )
260*
261 srnamt = 'ZLATMS'
262 CALL zlatms( m, n, dist, iseed, TYPE, RWORK, MODE,
263 $ CNDNUM, ANORM, KL, KU, 'No packing', A, LDA,
264 $ WORK, INFO )
265*
266* Check error code from ZLATMS.
267*
268 IF( info.NE.0 ) THEN
269 CALL alaerh( path, 'ZLATMS', info, 0, ' ', m, n, -1,
270 $ -1, -1, imat, nfail, nerrs, nout )
271 GO TO 100
272 END IF
273*
274* For types 5-7, zero one or more columns of the matrix to
275* test that INFO is returned correctly.
276*
277 IF( zerot ) THEN
278 IF( imat.EQ.5 ) THEN
279 izero = 1
280 ELSE IF( imat.EQ.6 ) THEN
281 izero = min( m, n )
282 ELSE
283 izero = min( m, n ) / 2 + 1
284 END IF
285 ioff = ( izero-1 )*lda
286 IF( imat.LT.7 ) THEN
287 DO 20 i = 1, m
288 a( ioff+i ) = zero
289 20 CONTINUE
290 ELSE
291 CALL zlaset( 'Full', m, n-izero+1, dcmplx(zero),
292 $ dcmplx(zero), a( ioff+1 ), lda )
293 END IF
294 ELSE
295 izero = 0
296 END IF
297*
298 DO 60 irhs = 1, nns
299 nrhs = nsval( irhs )
300 xtype = 'N'
301 trans = 'N'
302*
303 srnamt = 'ZLARHS'
304 CALL zlarhs( path, xtype, ' ', trans, n, n, kl,
305 $ ku, nrhs, a, lda, x, lda, b,
306 $ lda, iseed, info )
307*
308 srnamt = 'ZCGESV'
309*
310 kase = kase + 1
311*
312 CALL zlacpy( 'Full', m, n, a, lda, afac, lda )
313*
314 CALL zcgesv( n, nrhs, a, lda, iwork, b, lda, x, lda,
315 $ work, swork, rwork, iter, info)
316*
317 IF (iter.LT.0) THEN
318 CALL zlacpy( 'Full', m, n, afac, lda, a, lda )
319 ENDIF
320*
321* Check error code from ZCGESV. This should be the same as
322* the one of DGETRF.
323*
324 IF( info.NE.izero ) THEN
325*
326 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
327 $ CALL alahd( nout, path )
328 nerrs = nerrs + 1
329*
330 IF( info.NE.izero .AND. izero.NE.0 ) THEN
331 WRITE( nout, fmt = 9988 )'ZCGESV',info,
332 $ izero,m,imat
333 ELSE
334 WRITE( nout, fmt = 9975 )'ZCGESV',info,
335 $ m, imat
336 END IF
337 END IF
338*
339* Skip the remaining test if the matrix is singular.
340*
341 IF( info.NE.0 )
342 $ GO TO 100
343*
344* Check the quality of the solution
345*
346 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
347*
348 CALL zget08( trans, n, n, nrhs, a, lda, x, lda, work,
349 $ lda, rwork, result( 1 ) )
350*
351* Check if the test passes the testing.
352* Print information about the tests that did not
353* pass the testing.
354*
355* If iterative refinement has been used and claimed to
356* be successful (ITER>0), we want
357* NORMI(B - A*X)/(NORMI(A)*NORMI(X)*EPS*SRQT(N)) < 1
358*
359* If double precision has been used (ITER<0), we want
360* NORMI(B - A*X)/(NORMI(A)*NORMI(X)*EPS) < THRES
361* (Cf. the linear solver testing routines)
362*
363 IF ((thresh.LE.0.0e+00)
364 $ .OR.((iter.GE.0).AND.(n.GT.0)
365 $ .AND.(result(1).GE.sqrt(dble(n))))
366 $ .OR.((iter.LT.0).AND.(result(1).GE.thresh))) THEN
367*
368 IF( nfail.EQ.0 .AND. nerrs.EQ.0 ) THEN
369 WRITE( nout, fmt = 8999 )'DGE'
370 WRITE( nout, fmt = '( '' Matrix types:'' )' )
371 WRITE( nout, fmt = 8979 )
372 WRITE( nout, fmt = '( '' Test ratios:'' )' )
373 WRITE( nout, fmt = 8960 )1
374 WRITE( nout, fmt = '( '' Messages:'' )' )
375 END IF
376*
377 WRITE( nout, fmt = 9998 )trans, n, nrhs,
378 $ imat, 1, result( 1 )
379 nfail = nfail + 1
380 END IF
381 nrun = nrun + 1
382 60 CONTINUE
383 100 CONTINUE
384 120 CONTINUE
385*
386* Print a summary of the results.
387*
388 IF( nfail.GT.0 ) THEN
389 WRITE( nout, fmt = 9996 )'ZCGESV', nfail, nrun
390 ELSE
391 WRITE( nout, fmt = 9995 )'ZCGESV', nrun
392 END IF
393 IF( nerrs.GT.0 ) THEN
394 WRITE( nout, fmt = 9994 )nerrs
395 END IF
396*
397 9998 FORMAT( ' TRANS=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
398 $ i2, ', test(', i2, ') =', g12.5 )
399 9996 FORMAT( 1x, a6, ': ', i6, ' out of ', i6,
400 $ ' tests failed to pass the threshold' )
401 9995 FORMAT( /1x, 'All tests for ', a6,
402 $ ' routines passed the threshold ( ', i6, ' tests run)' )
403 9994 FORMAT( 6x, i6, ' error messages recorded' )
404*
405* SUBNAM, INFO, INFOE, M, IMAT
406*
407 9988 FORMAT( ' *** ', a6, ' returned with INFO =', i5, ' instead of ',
408 $ i5, / ' ==> M =', i5, ', type ',
409 $ i2 )
410*
411* SUBNAM, INFO, M, IMAT
412*
413 9975 FORMAT( ' *** Error code from ', a6, '=', i5, ' for M=', i5,
414 $ ', type ', i2 )
415 8999 FORMAT( / 1x, a3, ': General dense matrices' )
416 8979 FORMAT( 4x, '1. Diagonal', 24x, '7. Last n/2 columns zero', / 4x,
417 $ '2. Upper triangular', 16x,
418 $ '8. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
419 $ '3. Lower triangular', 16x, '9. Random, CNDNUM = 0.1/EPS',
420 $ / 4x, '4. Random, CNDNUM = 2', 13x,
421 $ '10. Scaled near underflow', / 4x, '5. First column zero',
422 $ 14x, '11. Scaled near overflow', / 4x,
423 $ '6. Last column zero' )
424 8960 FORMAT( 3x, i2, ': norm_1( B - A * X ) / ',
425 $ '( norm_1(A) * norm_1(X) * EPS * SQRT(N) ) > 1 if ITERREF',
426 $ / 4x, 'or norm_1( B - A * X ) / ',
427 $ '( norm_1(A) * norm_1(X) * EPS ) > THRES if DGETRF' )
428 RETURN
429*
430* End of ZDRVAB
431*
subroutine zlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
ZLARHS
Definition zlarhs.f:208
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 zcgesv(n, nrhs, a, lda, ipiv, b, ldb, x, ldx, work, swork, rwork, iter, info)
ZCGESV computes the solution to system of linear equations A * X = B for GE matrices (mixed precision...
Definition zcgesv.f:201
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine zget08(trans, m, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
ZGET08
Definition zget08.f:133
subroutine zlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
ZLATB4
Definition zlatb4.f:121
subroutine zlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
ZLATMS
Definition zlatms.f:332
Here is the call graph for this function:
Here is the caller graph for this function: