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

◆ zdrvac()

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.

Definition at line 142 of file zdrvac.f.

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