LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zdrvac ( 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  NOUT 
)

ZDRVAC

Purpose:
 ZDRVAC tests ZCPOSV.
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 N contained in the vector MVAL.
[in]MVAL
          MVAL is INTEGER array, dimension (NM)
          The values of the matrix dimension N.
[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 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)
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      (max(2*NMAX,2*NSMAX+NWORK))
[out]SWORK
          SWORK is COMPLEX array, dimension
                      (NMAX*(NSMAX+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 147 of file zdrvac.f.

147 *
148 * -- LAPACK test routine (version 3.4.0) --
149 * -- LAPACK is a software package provided by Univ. of Tennessee, --
150 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
151 * November 2011
152 *
153 * .. Scalar Arguments ..
154  INTEGER nmax, nm, nns, nout
155  DOUBLE PRECISION thresh
156 * ..
157 * .. Array Arguments ..
158  LOGICAL dotype( * )
159  INTEGER mval( * ), nsval( * )
160  DOUBLE PRECISION rwork( * )
161  COMPLEX swork(*)
162  COMPLEX*16 a( * ), afac( * ), b( * ),
163  $ work( * ), x( * )
164 * ..
165 *
166 * =====================================================================
167 *
168 * .. Parameters ..
169  DOUBLE PRECISION zero
170  parameter ( zero = 0.0d+0 )
171  INTEGER ntypes
172  parameter ( ntypes = 9 )
173  INTEGER ntests
174  parameter ( ntests = 1 )
175 * ..
176 * .. Local Scalars ..
177  LOGICAL zerot
178  CHARACTER dist, TYPE, uplo, xtype
179  CHARACTER*3 path
180  INTEGER i, im, imat, info, ioff, irhs, iuplo,
181  $ izero, kl, ku, lda, mode, n,
182  $ nerrs, nfail, nimat, nrhs, nrun
183  DOUBLE PRECISION anorm, cndnum
184 * ..
185 * .. Local Arrays ..
186  CHARACTER uplos( 2 )
187  INTEGER iseed( 4 ), iseedy( 4 )
188  DOUBLE PRECISION result( ntests )
189 * ..
190 * .. Local Variables ..
191  INTEGER iter, kase
192 * ..
193 * .. External Subroutines ..
194  EXTERNAL alaerh, zlacpy, zlaipd,
195  $ zlarhs, zlatb4, zlatms,
196  $ zpot06, zcposv
197 * ..
198 * .. Intrinsic Functions ..
199  INTRINSIC dble, max, sqrt
200 * ..
201 * .. Scalars in Common ..
202  LOGICAL lerr, ok
203  CHARACTER*32 srnamt
204  INTEGER infot, nunit
205 * ..
206 * .. Common blocks ..
207  COMMON / infoc / infot, nunit, ok, lerr
208  COMMON / srnamc / srnamt
209 * ..
210 * .. Data statements ..
211  DATA iseedy / 1988, 1989, 1990, 1991 /
212  DATA uplos / 'U', 'L' /
213 * ..
214 * .. Executable Statements ..
215 *
216 * Initialize constants and the random number seed.
217 *
218  kase = 0
219  path( 1: 1 ) = 'Zomplex precision'
220  path( 2: 3 ) = 'PO'
221  nrun = 0
222  nfail = 0
223  nerrs = 0
224  DO 10 i = 1, 4
225  iseed( i ) = iseedy( i )
226  10 CONTINUE
227 *
228  infot = 0
229 *
230 * Do for each value of N in MVAL
231 *
232  DO 120 im = 1, nm
233  n = mval( im )
234  lda = max( n, 1 )
235  nimat = ntypes
236  IF( n.LE.0 )
237  $ nimat = 1
238 *
239  DO 110 imat = 1, nimat
240 *
241 * Do the tests only if DOTYPE( IMAT ) is true.
242 *
243  IF( .NOT.dotype( imat ) )
244  $ GO TO 110
245 *
246 * Skip types 3, 4, or 5 if the matrix size is too small.
247 *
248  zerot = imat.GE.3 .AND. imat.LE.5
249  IF( zerot .AND. n.LT.imat-2 )
250  $ GO TO 110
251 *
252 * Do first for UPLO = 'U', then for UPLO = 'L'
253 *
254  DO 100 iuplo = 1, 2
255  uplo = uplos( iuplo )
256 *
257 * Set up parameters with ZLATB4 and generate a test matrix
258 * with ZLATMS.
259 *
260  CALL zlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
261  $ cndnum, dist )
262 *
263  srnamt = 'ZLATMS'
264  CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
265  $ cndnum, anorm, kl, ku, uplo, a, lda, work,
266  $ info )
267 *
268 * Check error code from ZLATMS.
269 *
270  IF( info.NE.0 ) THEN
271  CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n, -1,
272  $ -1, -1, imat, nfail, nerrs, nout )
273  GO TO 100
274  END IF
275 *
276 * For types 3-5, zero one row and column of the matrix to
277 * test that INFO is returned correctly.
278 *
279  IF( zerot ) THEN
280  IF( imat.EQ.3 ) THEN
281  izero = 1
282  ELSE IF( imat.EQ.4 ) THEN
283  izero = n
284  ELSE
285  izero = n / 2 + 1
286  END IF
287  ioff = ( izero-1 )*lda
288 *
289 * Set row and column IZERO of A to 0.
290 *
291  IF( iuplo.EQ.1 ) THEN
292  DO 20 i = 1, izero - 1
293  a( ioff+i ) = zero
294  20 CONTINUE
295  ioff = ioff + izero
296  DO 30 i = izero, n
297  a( ioff ) = zero
298  ioff = ioff + lda
299  30 CONTINUE
300  ELSE
301  ioff = izero
302  DO 40 i = 1, izero - 1
303  a( ioff ) = zero
304  ioff = ioff + lda
305  40 CONTINUE
306  ioff = ioff - izero
307  DO 50 i = izero, n
308  a( ioff+i ) = zero
309  50 CONTINUE
310  END IF
311  ELSE
312  izero = 0
313  END IF
314 *
315 * Set the imaginary part of the diagonals.
316 *
317  CALL zlaipd( n, a, lda+1, 0 )
318 *
319  DO 60 irhs = 1, nns
320  nrhs = nsval( irhs )
321  xtype = 'N'
322 *
323 * Form an exact solution and set the right hand side.
324 *
325  srnamt = 'ZLARHS'
326  CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
327  $ nrhs, a, lda, x, lda, b, lda,
328  $ iseed, info )
329 *
330 * Compute the L*L' or U'*U factorization of the
331 * matrix and solve the system.
332 *
333  srnamt = 'ZCPOSV '
334  kase = kase + 1
335 *
336  CALL zlacpy( 'All', n, n, a, lda, afac, lda)
337 *
338  CALL zcposv( uplo, n, nrhs, afac, lda, b, lda, x, lda,
339  $ work, swork, rwork, iter, info )
340 *
341  IF (iter.LT.0) THEN
342  CALL zlacpy( 'All', n, n, a, lda, afac, lda )
343  ENDIF
344 *
345 * Check error code from ZCPOSV .
346 *
347  IF( info.NE.izero ) THEN
348 *
349  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
350  $ CALL alahd( nout, path )
351  nerrs = nerrs + 1
352 *
353  IF( info.NE.izero .AND. izero.NE.0 ) THEN
354  WRITE( nout, fmt = 9988 )'ZCPOSV',info,izero,n,
355  $ imat
356  ELSE
357  WRITE( nout, fmt = 9975 )'ZCPOSV',info,n,imat
358  END IF
359  END IF
360 *
361 * Skip the remaining test if the matrix is singular.
362 *
363  IF( info.NE.0 )
364  $ GO TO 110
365 *
366 * Check the quality of the solution
367 *
368  CALL zlacpy( 'All', n, nrhs, b, lda, work, lda )
369 *
370  CALL zpot06( uplo, n, nrhs, a, lda, x, lda, work,
371  $ lda, rwork, result( 1 ) )
372 *
373 * Check if the test passes the tesing.
374 * Print information about the tests that did not
375 * pass the testing.
376 *
377 * If iterative refinement has been used and claimed to
378 * be successful (ITER>0), we want
379 * NORM1(B - A*X)/(NORM1(A)*NORM1(X)*EPS*SRQT(N)) < 1
380 *
381 * If double precision has been used (ITER<0), we want
382 * NORM1(B - A*X)/(NORM1(A)*NORM1(X)*EPS) < THRES
383 * (Cf. the linear solver testing routines)
384 *
385  IF ((thresh.LE.0.0e+00)
386  $ .OR.((iter.GE.0).AND.(n.GT.0)
387  $ .AND.(result(1).GE.sqrt(dble(n))))
388  $ .OR.((iter.LT.0).AND.(result(1).GE.thresh))) THEN
389 *
390  IF( nfail.EQ.0 .AND. nerrs.EQ.0 ) THEN
391  WRITE( nout, fmt = 8999 )'ZPO'
392  WRITE( nout, fmt = '( '' Matrix types:'' )' )
393  WRITE( nout, fmt = 8979 )
394  WRITE( nout, fmt = '( '' Test ratios:'' )' )
395  WRITE( nout, fmt = 8960 )1
396  WRITE( nout, fmt = '( '' Messages:'' )' )
397  END IF
398 *
399  WRITE( nout, fmt = 9998 )uplo, n, nrhs, imat, 1,
400  $ result( 1 )
401 *
402  nfail = nfail + 1
403 *
404  END IF
405 *
406  nrun = nrun + 1
407 *
408  60 CONTINUE
409  100 CONTINUE
410  110 CONTINUE
411  120 CONTINUE
412 *
413 * Print a summary of the results.
414 *
415  IF( nfail.GT.0 ) THEN
416  WRITE( nout, fmt = 9996 )'ZCPOSV', nfail, nrun
417  ELSE
418  WRITE( nout, fmt = 9995 )'ZCPOSV', nrun
419  END IF
420  IF( nerrs.GT.0 ) THEN
421  WRITE( nout, fmt = 9994 )nerrs
422  END IF
423 *
424  9998 FORMAT( ' UPLO=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
425  $ i2, ', test(', i2, ') =', g12.5 )
426  9996 FORMAT( 1x, a6, ': ', i6, ' out of ', i6,
427  $ ' tests failed to pass the threshold' )
428  9995 FORMAT( /1x, 'All tests for ', a6,
429  $ ' routines passed the threshold ( ', i6, ' tests run)' )
430  9994 FORMAT( 6x, i6, ' error messages recorded' )
431 *
432 * SUBNAM, INFO, INFOE, N, IMAT
433 *
434  9988 FORMAT( ' *** ', a6, ' returned with INFO =', i5, ' instead of ',
435  $ i5, / ' ==> N =', i5, ', type ',
436  $ i2 )
437 *
438 * SUBNAM, INFO, N, IMAT
439 *
440  9975 FORMAT( ' *** Error code from ', a6, '=', i5, ' for M=', i5,
441  $ ', type ', i2 )
442  8999 FORMAT( / 1x, a3, ': positive definite dense matrices' )
443  8979 FORMAT( 4x, '1. Diagonal', 24x, '7. Last n/2 columns zero', / 4x,
444  $ '2. Upper triangular', 16x,
445  $ '8. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
446  $ '3. Lower triangular', 16x, '9. Random, CNDNUM = 0.1/EPS',
447  $ / 4x, '4. Random, CNDNUM = 2', 13x,
448  $ '10. Scaled near underflow', / 4x, '5. First column zero',
449  $ 14x, '11. Scaled near overflow', / 4x,
450  $ '6. Last column zero' )
451  8960 FORMAT( 3x, i2, ': norm_1( B - A * X ) / ',
452  $ '( norm_1(A) * norm_1(X) * EPS * SQRT(N) ) > 1 if ITERREF',
453  $ / 4x, 'or norm_1( B - A * X ) / ',
454  $ '( norm_1(A) * norm_1(X) * EPS ) > THRES if ZPOTRF' )
455 
456  RETURN
457 *
458 * End of ZDRVAC
459 *
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine zlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
ZLARHS
Definition: zlarhs.f:211
subroutine zlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
ZLATB4
Definition: zlatb4.f:123
subroutine zlaipd(N, A, INDA, VINDA)
ZLAIPD
Definition: zlaipd.f:85
subroutine zcposv(UPLO, N, NRHS, A, LDA, B, LDB, X, LDX, WORK, SWORK, RWORK, ITER, INFO)
ZCPOSV computes the solution to system of linear equations A * X = B for PO matrices ...
Definition: zcposv.f:211
subroutine zpot06(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
ZPOT06
Definition: zpot06.f:129
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334

Here is the call graph for this function:

Here is the caller graph for this function: