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

◆ zdrvpt()

subroutine zdrvpt ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nrhs,
double precision  thresh,
logical  tsterr,
complex*16, dimension( * )  a,
double precision, dimension( * )  d,
complex*16, dimension( * )  e,
complex*16, dimension( * )  b,
complex*16, dimension( * )  x,
complex*16, dimension( * )  xact,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork,
integer  nout 
)

ZDRVPT

Purpose:
 ZDRVPT tests ZPTSV and -SVX.
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]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 dimension N.
[in]NRHS
          NRHS is INTEGER
          The number of right hand side vectors to be generated for
          each linear system.
[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]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[out]A
          A is COMPLEX*16 array, dimension (NMAX*2)
[out]D
          D is DOUBLE PRECISION array, dimension (NMAX*2)
[out]E
          E is COMPLEX*16 array, dimension (NMAX*2)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]XACT
          XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
[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 138 of file zdrvpt.f.

140*
141* -- LAPACK test routine --
142* -- LAPACK is a software package provided by Univ. of Tennessee, --
143* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
144*
145* .. Scalar Arguments ..
146 LOGICAL TSTERR
147 INTEGER NN, NOUT, NRHS
148 DOUBLE PRECISION THRESH
149* ..
150* .. Array Arguments ..
151 LOGICAL DOTYPE( * )
152 INTEGER NVAL( * )
153 DOUBLE PRECISION D( * ), RWORK( * )
154 COMPLEX*16 A( * ), B( * ), E( * ), WORK( * ), X( * ),
155 $ XACT( * )
156* ..
157*
158* =====================================================================
159*
160* .. Parameters ..
161 DOUBLE PRECISION ONE, ZERO
162 parameter( one = 1.0d+0, zero = 0.0d+0 )
163 INTEGER NTYPES
164 parameter( ntypes = 12 )
165 INTEGER NTESTS
166 parameter( ntests = 6 )
167* ..
168* .. Local Scalars ..
169 LOGICAL ZEROT
170 CHARACTER DIST, FACT, TYPE
171 CHARACTER*3 PATH
172 INTEGER I, IA, IFACT, IMAT, IN, INFO, IX, IZERO, J, K,
173 $ K1, KL, KU, LDA, MODE, N, NERRS, NFAIL, NIMAT,
174 $ NRUN, NT
175 DOUBLE PRECISION AINVNM, ANORM, COND, DMAX, RCOND, RCONDC
176* ..
177* .. Local Arrays ..
178 INTEGER ISEED( 4 ), ISEEDY( 4 )
179 DOUBLE PRECISION RESULT( NTESTS ), Z( 3 )
180* ..
181* .. External Functions ..
182 INTEGER IDAMAX
183 DOUBLE PRECISION DGET06, DZASUM, ZLANHT
184 EXTERNAL idamax, dget06, dzasum, zlanht
185* ..
186* .. External Subroutines ..
187 EXTERNAL aladhd, alaerh, alasvm, dcopy, dlarnv, dscal,
191* ..
192* .. Intrinsic Functions ..
193 INTRINSIC abs, dcmplx, max
194* ..
195* .. Scalars in Common ..
196 LOGICAL LERR, OK
197 CHARACTER*32 SRNAMT
198 INTEGER INFOT, NUNIT
199* ..
200* .. Common blocks ..
201 COMMON / infoc / infot, nunit, ok, lerr
202 COMMON / srnamc / srnamt
203* ..
204* .. Data statements ..
205 DATA iseedy / 0, 0, 0, 1 /
206* ..
207* .. Executable Statements ..
208*
209 path( 1: 1 ) = 'Zomplex precision'
210 path( 2: 3 ) = 'PT'
211 nrun = 0
212 nfail = 0
213 nerrs = 0
214 DO 10 i = 1, 4
215 iseed( i ) = iseedy( i )
216 10 CONTINUE
217*
218* Test the error exits
219*
220 IF( tsterr )
221 $ CALL zerrvx( path, nout )
222 infot = 0
223*
224 DO 120 in = 1, nn
225*
226* Do for each value of N in NVAL.
227*
228 n = nval( in )
229 lda = max( 1, n )
230 nimat = ntypes
231 IF( n.LE.0 )
232 $ nimat = 1
233*
234 DO 110 imat = 1, nimat
235*
236* Do the tests only if DOTYPE( IMAT ) is true.
237*
238 IF( n.GT.0 .AND. .NOT.dotype( imat ) )
239 $ GO TO 110
240*
241* Set up parameters with ZLATB4.
242*
243 CALL zlatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
244 $ COND, DIST )
245*
246 zerot = imat.GE.8 .AND. imat.LE.10
247 IF( imat.LE.6 ) THEN
248*
249* Type 1-6: generate a symmetric tridiagonal matrix of
250* known condition number in lower triangular band storage.
251*
252 srnamt = 'ZLATMS'
253 CALL zlatms( n, n, dist, iseed, TYPE, RWORK, MODE, COND,
254 $ ANORM, KL, KU, 'B', A, 2, WORK, INFO )
255*
256* Check the error code from ZLATMS.
257*
258 IF( info.NE.0 ) THEN
259 CALL alaerh( path, 'ZLATMS', info, 0, ' ', n, n, kl,
260 $ ku, -1, imat, nfail, nerrs, nout )
261 GO TO 110
262 END IF
263 izero = 0
264*
265* Copy the matrix to D and E.
266*
267 ia = 1
268 DO 20 i = 1, n - 1
269 d( i ) = dble( a( ia ) )
270 e( i ) = a( ia+1 )
271 ia = ia + 2
272 20 CONTINUE
273 IF( n.GT.0 )
274 $ d( n ) = dble( a( ia ) )
275 ELSE
276*
277* Type 7-12: generate a diagonally dominant matrix with
278* unknown condition number in the vectors D and E.
279*
280 IF( .NOT.zerot .OR. .NOT.dotype( 7 ) ) THEN
281*
282* Let D and E have values from [-1,1].
283*
284 CALL dlarnv( 2, iseed, n, d )
285 CALL zlarnv( 2, iseed, n-1, e )
286*
287* Make the tridiagonal matrix diagonally dominant.
288*
289 IF( n.EQ.1 ) THEN
290 d( 1 ) = abs( d( 1 ) )
291 ELSE
292 d( 1 ) = abs( d( 1 ) ) + abs( e( 1 ) )
293 d( n ) = abs( d( n ) ) + abs( e( n-1 ) )
294 DO 30 i = 2, n - 1
295 d( i ) = abs( d( i ) ) + abs( e( i ) ) +
296 $ abs( e( i-1 ) )
297 30 CONTINUE
298 END IF
299*
300* Scale D and E so the maximum element is ANORM.
301*
302 ix = idamax( n, d, 1 )
303 dmax = d( ix )
304 CALL dscal( n, anorm / dmax, d, 1 )
305 IF( n.GT.1 )
306 $ CALL zdscal( n-1, anorm / dmax, e, 1 )
307*
308 ELSE IF( izero.GT.0 ) THEN
309*
310* Reuse the last matrix by copying back the zeroed out
311* elements.
312*
313 IF( izero.EQ.1 ) THEN
314 d( 1 ) = z( 2 )
315 IF( n.GT.1 )
316 $ e( 1 ) = z( 3 )
317 ELSE IF( izero.EQ.n ) THEN
318 e( n-1 ) = z( 1 )
319 d( n ) = z( 2 )
320 ELSE
321 e( izero-1 ) = z( 1 )
322 d( izero ) = z( 2 )
323 e( izero ) = z( 3 )
324 END IF
325 END IF
326*
327* For types 8-10, set one row and column of the matrix to
328* zero.
329*
330 izero = 0
331 IF( imat.EQ.8 ) THEN
332 izero = 1
333 z( 2 ) = d( 1 )
334 d( 1 ) = zero
335 IF( n.GT.1 ) THEN
336 z( 3 ) = dble( e( 1 ) )
337 e( 1 ) = zero
338 END IF
339 ELSE IF( imat.EQ.9 ) THEN
340 izero = n
341 IF( n.GT.1 ) THEN
342 z( 1 ) = dble( e( n-1 ) )
343 e( n-1 ) = zero
344 END IF
345 z( 2 ) = d( n )
346 d( n ) = zero
347 ELSE IF( imat.EQ.10 ) THEN
348 izero = ( n+1 ) / 2
349 IF( izero.GT.1 ) THEN
350 z( 1 ) = dble( e( izero-1 ) )
351 e( izero-1 ) = zero
352 z( 3 ) = dble( e( izero ) )
353 e( izero ) = zero
354 END IF
355 z( 2 ) = d( izero )
356 d( izero ) = zero
357 END IF
358 END IF
359*
360* Generate NRHS random solution vectors.
361*
362 ix = 1
363 DO 40 j = 1, nrhs
364 CALL zlarnv( 2, iseed, n, xact( ix ) )
365 ix = ix + lda
366 40 CONTINUE
367*
368* Set the right hand side.
369*
370 CALL zlaptm( 'Lower', n, nrhs, one, d, e, xact, lda, zero,
371 $ b, lda )
372*
373 DO 100 ifact = 1, 2
374 IF( ifact.EQ.1 ) THEN
375 fact = 'F'
376 ELSE
377 fact = 'N'
378 END IF
379*
380* Compute the condition number for comparison with
381* the value returned by ZPTSVX.
382*
383 IF( zerot ) THEN
384 IF( ifact.EQ.1 )
385 $ GO TO 100
386 rcondc = zero
387*
388 ELSE IF( ifact.EQ.1 ) THEN
389*
390* Compute the 1-norm of A.
391*
392 anorm = zlanht( '1', n, d, e )
393*
394 CALL dcopy( n, d, 1, d( n+1 ), 1 )
395 IF( n.GT.1 )
396 $ CALL zcopy( n-1, e, 1, e( n+1 ), 1 )
397*
398* Factor the matrix A.
399*
400 CALL zpttrf( n, d( n+1 ), e( n+1 ), info )
401*
402* Use ZPTTRS to solve for one column at a time of
403* inv(A), computing the maximum column sum as we go.
404*
405 ainvnm = zero
406 DO 60 i = 1, n
407 DO 50 j = 1, n
408 x( j ) = zero
409 50 CONTINUE
410 x( i ) = one
411 CALL zpttrs( 'Lower', n, 1, d( n+1 ), e( n+1 ), x,
412 $ lda, info )
413 ainvnm = max( ainvnm, dzasum( n, x, 1 ) )
414 60 CONTINUE
415*
416* Compute the 1-norm condition number of A.
417*
418 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
419 rcondc = one
420 ELSE
421 rcondc = ( one / anorm ) / ainvnm
422 END IF
423 END IF
424*
425 IF( ifact.EQ.2 ) THEN
426*
427* --- Test ZPTSV --
428*
429 CALL dcopy( n, d, 1, d( n+1 ), 1 )
430 IF( n.GT.1 )
431 $ CALL zcopy( n-1, e, 1, e( n+1 ), 1 )
432 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
433*
434* Factor A as L*D*L' and solve the system A*X = B.
435*
436 srnamt = 'ZPTSV '
437 CALL zptsv( n, nrhs, d( n+1 ), e( n+1 ), x, lda,
438 $ info )
439*
440* Check error code from ZPTSV .
441*
442 IF( info.NE.izero )
443 $ CALL alaerh( path, 'ZPTSV ', info, izero, ' ', n,
444 $ n, 1, 1, nrhs, imat, nfail, nerrs,
445 $ nout )
446 nt = 0
447 IF( izero.EQ.0 ) THEN
448*
449* Check the factorization by computing the ratio
450* norm(L*D*L' - A) / (n * norm(A) * EPS )
451*
452 CALL zptt01( n, d, e, d( n+1 ), e( n+1 ), work,
453 $ result( 1 ) )
454*
455* Compute the residual in the solution.
456*
457 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
458 CALL zptt02( 'Lower', n, nrhs, d, e, x, lda, work,
459 $ lda, result( 2 ) )
460*
461* Check solution from generated exact solution.
462*
463 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
464 $ result( 3 ) )
465 nt = 3
466 END IF
467*
468* Print information about the tests that did not pass
469* the threshold.
470*
471 DO 70 k = 1, nt
472 IF( result( k ).GE.thresh ) THEN
473 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
474 $ CALL aladhd( nout, path )
475 WRITE( nout, fmt = 9999 )'ZPTSV ', n, imat, k,
476 $ result( k )
477 nfail = nfail + 1
478 END IF
479 70 CONTINUE
480 nrun = nrun + nt
481 END IF
482*
483* --- Test ZPTSVX ---
484*
485 IF( ifact.GT.1 ) THEN
486*
487* Initialize D( N+1:2*N ) and E( N+1:2*N ) to zero.
488*
489 DO 80 i = 1, n - 1
490 d( n+i ) = zero
491 e( n+i ) = zero
492 80 CONTINUE
493 IF( n.GT.0 )
494 $ d( n+n ) = zero
495 END IF
496*
497 CALL zlaset( 'Full', n, nrhs, dcmplx( zero ),
498 $ dcmplx( zero ), x, lda )
499*
500* Solve the system and compute the condition number and
501* error bounds using ZPTSVX.
502*
503 srnamt = 'ZPTSVX'
504 CALL zptsvx( fact, n, nrhs, d, e, d( n+1 ), e( n+1 ), b,
505 $ lda, x, lda, rcond, rwork, rwork( nrhs+1 ),
506 $ work, rwork( 2*nrhs+1 ), info )
507*
508* Check the error code from ZPTSVX.
509*
510 IF( info.NE.izero )
511 $ CALL alaerh( path, 'ZPTSVX', info, izero, fact, n, n,
512 $ 1, 1, nrhs, imat, nfail, nerrs, nout )
513 IF( izero.EQ.0 ) THEN
514 IF( ifact.EQ.2 ) THEN
515*
516* Check the factorization by computing the ratio
517* norm(L*D*L' - A) / (n * norm(A) * EPS )
518*
519 k1 = 1
520 CALL zptt01( n, d, e, d( n+1 ), e( n+1 ), work,
521 $ result( 1 ) )
522 ELSE
523 k1 = 2
524 END IF
525*
526* Compute the residual in the solution.
527*
528 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
529 CALL zptt02( 'Lower', n, nrhs, d, e, x, lda, work,
530 $ lda, result( 2 ) )
531*
532* Check solution from generated exact solution.
533*
534 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
535 $ result( 3 ) )
536*
537* Check error bounds from iterative refinement.
538*
539 CALL zptt05( n, nrhs, d, e, b, lda, x, lda, xact, lda,
540 $ rwork, rwork( nrhs+1 ), result( 4 ) )
541 ELSE
542 k1 = 6
543 END IF
544*
545* Check the reciprocal of the condition number.
546*
547 result( 6 ) = dget06( rcond, rcondc )
548*
549* Print information about the tests that did not pass
550* the threshold.
551*
552 DO 90 k = k1, 6
553 IF( result( k ).GE.thresh ) THEN
554 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
555 $ CALL aladhd( nout, path )
556 WRITE( nout, fmt = 9998 )'ZPTSVX', fact, n, imat,
557 $ k, result( k )
558 nfail = nfail + 1
559 END IF
560 90 CONTINUE
561 nrun = nrun + 7 - k1
562 100 CONTINUE
563 110 CONTINUE
564 120 CONTINUE
565*
566* Print a summary of the results.
567*
568 CALL alasvm( path, nout, nfail, nrun, nerrs )
569*
570 9999 FORMAT( 1x, a, ', N =', i5, ', type ', i2, ', test ', i2,
571 $ ', ratio = ', g12.5 )
572 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', N =', i5, ', type ', i2,
573 $ ', test ', i2, ', ratio = ', g12.5 )
574 RETURN
575*
576* End of ZDRVPT
577*
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine aladhd(iounit, path)
ALADHD
Definition aladhd.f:90
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
double precision function dget06(rcond, rcondc)
DGET06
Definition dget06.f:55
double precision function dzasum(n, zx, incx)
DZASUM
Definition dzasum.f:72
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
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
double precision function zlanht(norm, n, d, e)
ZLANHT returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanht.f:101
subroutine zlarnv(idist, iseed, n, x)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition zlarnv.f:99
subroutine dlarnv(idist, iseed, n, x)
DLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition dlarnv.f:97
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 zptsv(n, nrhs, d, e, b, ldb, info)
ZPTSV computes the solution to system of linear equations A * X = B for PT matrices
Definition zptsv.f:115
subroutine zptsvx(fact, n, nrhs, d, e, df, ef, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
ZPTSVX computes the solution to system of linear equations A * X = B for PT matrices
Definition zptsvx.f:234
subroutine zpttrf(n, d, e, info)
ZPTTRF
Definition zpttrf.f:92
subroutine zpttrs(uplo, n, nrhs, d, e, b, ldb, info)
ZPTTRS
Definition zpttrs.f:121
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
subroutine zerrvx(path, nunit)
ZERRVX
Definition zerrvx.f:55
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zlaptm(uplo, n, nrhs, alpha, d, e, x, ldx, beta, b, ldb)
ZLAPTM
Definition zlaptm.f:129
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 zptt01(n, d, e, df, ef, work, resid)
ZPTT01
Definition zptt01.f:92
subroutine zptt02(uplo, n, nrhs, d, e, x, ldx, b, ldb, resid)
ZPTT02
Definition zptt02.f:115
subroutine zptt05(n, nrhs, d, e, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZPTT05
Definition zptt05.f:150
Here is the call graph for this function:
Here is the caller graph for this function: