LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zdrvpo.f
Go to the documentation of this file.
1*> \brief \b ZDRVPO
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE ZDRVPO( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
12* A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
13* RWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER NMAX, NN, NOUT, NRHS
18* DOUBLE PRECISION THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER NVAL( * )
23* DOUBLE PRECISION RWORK( * ), S( * )
24* COMPLEX*16 A( * ), AFAC( * ), ASAV( * ), B( * ),
25* $ BSAV( * ), WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> ZDRVPO tests the driver routines ZPOSV and -SVX.
35*> \endverbatim
36*
37* Arguments:
38* ==========
39*
40*> \param[in] DOTYPE
41*> \verbatim
42*> DOTYPE is LOGICAL array, dimension (NTYPES)
43*> The matrix types to be used for testing. Matrices of type j
44*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
45*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
46*> \endverbatim
47*>
48*> \param[in] NN
49*> \verbatim
50*> NN is INTEGER
51*> The number of values of N contained in the vector NVAL.
52*> \endverbatim
53*>
54*> \param[in] NVAL
55*> \verbatim
56*> NVAL is INTEGER array, dimension (NN)
57*> The values of the matrix dimension N.
58*> \endverbatim
59*>
60*> \param[in] NRHS
61*> \verbatim
62*> NRHS is INTEGER
63*> The number of right hand side vectors to be generated for
64*> each linear system.
65*> \endverbatim
66*>
67*> \param[in] THRESH
68*> \verbatim
69*> THRESH is DOUBLE PRECISION
70*> The threshold value for the test ratios. A result is
71*> included in the output file if RESULT >= THRESH. To have
72*> every test ratio printed, use THRESH = 0.
73*> \endverbatim
74*>
75*> \param[in] TSTERR
76*> \verbatim
77*> TSTERR is LOGICAL
78*> Flag that indicates whether error exits are to be tested.
79*> \endverbatim
80*>
81*> \param[in] NMAX
82*> \verbatim
83*> NMAX is INTEGER
84*> The maximum value permitted for N, used in dimensioning the
85*> work arrays.
86*> \endverbatim
87*>
88*> \param[out] A
89*> \verbatim
90*> A is COMPLEX*16 array, dimension (NMAX*NMAX)
91*> \endverbatim
92*>
93*> \param[out] AFAC
94*> \verbatim
95*> AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
96*> \endverbatim
97*>
98*> \param[out] ASAV
99*> \verbatim
100*> ASAV is COMPLEX*16 array, dimension (NMAX*NMAX)
101*> \endverbatim
102*>
103*> \param[out] B
104*> \verbatim
105*> B is COMPLEX*16 array, dimension (NMAX*NRHS)
106*> \endverbatim
107*>
108*> \param[out] BSAV
109*> \verbatim
110*> BSAV is COMPLEX*16 array, dimension (NMAX*NRHS)
111*> \endverbatim
112*>
113*> \param[out] X
114*> \verbatim
115*> X is COMPLEX*16 array, dimension (NMAX*NRHS)
116*> \endverbatim
117*>
118*> \param[out] XACT
119*> \verbatim
120*> XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
121*> \endverbatim
122*>
123*> \param[out] S
124*> \verbatim
125*> S is DOUBLE PRECISION array, dimension (NMAX)
126*> \endverbatim
127*>
128*> \param[out] WORK
129*> \verbatim
130*> WORK is COMPLEX*16 array, dimension
131*> (NMAX*max(3,NRHS))
132*> \endverbatim
133*>
134*> \param[out] RWORK
135*> \verbatim
136*> RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
137*> \endverbatim
138*>
139*> \param[in] NOUT
140*> \verbatim
141*> NOUT is INTEGER
142*> The unit number for output.
143*> \endverbatim
144*
145* Authors:
146* ========
147*
148*> \author Univ. of Tennessee
149*> \author Univ. of California Berkeley
150*> \author Univ. of Colorado Denver
151*> \author NAG Ltd.
152*
153*> \ingroup complex16_lin
154*
155* =====================================================================
156 SUBROUTINE zdrvpo( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
157 $ A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
158 $ RWORK, NOUT )
159*
160* -- LAPACK test routine --
161* -- LAPACK is a software package provided by Univ. of Tennessee, --
162* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
163*
164* .. Scalar Arguments ..
165 LOGICAL TSTERR
166 INTEGER NMAX, NN, NOUT, NRHS
167 DOUBLE PRECISION THRESH
168* ..
169* .. Array Arguments ..
170 LOGICAL DOTYPE( * )
171 INTEGER NVAL( * )
172 DOUBLE PRECISION RWORK( * ), S( * )
173 COMPLEX*16 A( * ), AFAC( * ), ASAV( * ), B( * ),
174 $ bsav( * ), work( * ), x( * ), xact( * )
175* ..
176*
177* =====================================================================
178*
179* .. Parameters ..
180 DOUBLE PRECISION ONE, ZERO
181 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
182 INTEGER NTYPES
183 parameter( ntypes = 9 )
184 INTEGER NTESTS
185 parameter( ntests = 6 )
186* ..
187* .. Local Scalars ..
188 LOGICAL EQUIL, NOFACT, PREFAC, ZEROT
189 CHARACTER DIST, EQUED, FACT, TYPE, UPLO, XTYPE
190 CHARACTER*3 PATH
191 INTEGER I, IEQUED, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
192 $ izero, k, k1, kl, ku, lda, mode, n, nb, nbmin,
193 $ nerrs, nfact, nfail, nimat, nrun, nt
194 DOUBLE PRECISION AINVNM, AMAX, ANORM, CNDNUM, RCOND, RCONDC,
195 $ ROLDC, SCOND
196* ..
197* .. Local Arrays ..
198 CHARACTER EQUEDS( 2 ), FACTS( 3 ), UPLOS( 2 )
199 INTEGER ISEED( 4 ), ISEEDY( 4 )
200 DOUBLE PRECISION RESULT( NTESTS )
201* ..
202* .. External Functions ..
203 LOGICAL LSAME
204 DOUBLE PRECISION DGET06, ZLANHE
205 EXTERNAL lsame, dget06, zlanhe
206* ..
207* .. External Subroutines ..
208 EXTERNAL aladhd, alaerh, alasvm, xlaenv, zerrvx, zget04,
212* ..
213* .. Scalars in Common ..
214 LOGICAL LERR, OK
215 CHARACTER*32 SRNAMT
216 INTEGER INFOT, NUNIT
217* ..
218* .. Common blocks ..
219 COMMON / infoc / infot, nunit, ok, lerr
220 COMMON / srnamc / srnamt
221* ..
222* .. Intrinsic Functions ..
223 INTRINSIC dcmplx, max
224* ..
225* .. Data statements ..
226 DATA iseedy / 1988, 1989, 1990, 1991 /
227 DATA uplos / 'U', 'L' /
228 DATA facts / 'F', 'N', 'E' /
229 DATA equeds / 'N', 'Y' /
230* ..
231* .. Executable Statements ..
232*
233* Initialize constants and the random number seed.
234*
235 path( 1: 1 ) = 'Zomplex precision'
236 path( 2: 3 ) = 'PO'
237 nrun = 0
238 nfail = 0
239 nerrs = 0
240 DO 10 i = 1, 4
241 iseed( i ) = iseedy( i )
242 10 CONTINUE
243*
244* Test the error exits
245*
246 IF( tsterr )
247 $ CALL zerrvx( path, nout )
248 infot = 0
249*
250* Set the block size and minimum block size for testing.
251*
252 nb = 1
253 nbmin = 2
254 CALL xlaenv( 1, nb )
255 CALL xlaenv( 2, nbmin )
256*
257* Do for each value of N in NVAL
258*
259 DO 130 in = 1, nn
260 n = nval( in )
261 lda = max( n, 1 )
262 xtype = 'N'
263 nimat = ntypes
264 IF( n.LE.0 )
265 $ nimat = 1
266*
267 DO 120 imat = 1, nimat
268*
269* Do the tests only if DOTYPE( IMAT ) is true.
270*
271 IF( .NOT.dotype( imat ) )
272 $ GO TO 120
273*
274* Skip types 3, 4, or 5 if the matrix size is too small.
275*
276 zerot = imat.GE.3 .AND. imat.LE.5
277 IF( zerot .AND. n.LT.imat-2 )
278 $ GO TO 120
279*
280* Do first for UPLO = 'U', then for UPLO = 'L'
281*
282 DO 110 iuplo = 1, 2
283 uplo = uplos( iuplo )
284*
285* Set up parameters with ZLATB4 and generate a test matrix
286* with ZLATMS.
287*
288 CALL zlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
289 $ cndnum, dist )
290*
291 srnamt = 'ZLATMS'
292 CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
293 $ cndnum, anorm, kl, ku, uplo, a, lda, work,
294 $ info )
295*
296* Check error code from ZLATMS.
297*
298 IF( info.NE.0 ) THEN
299 CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n, -1,
300 $ -1, -1, imat, nfail, nerrs, nout )
301 GO TO 110
302 END IF
303*
304* For types 3-5, zero one row and column of the matrix to
305* test that INFO is returned correctly.
306*
307 IF( zerot ) THEN
308 IF( imat.EQ.3 ) THEN
309 izero = 1
310 ELSE IF( imat.EQ.4 ) THEN
311 izero = n
312 ELSE
313 izero = n / 2 + 1
314 END IF
315 ioff = ( izero-1 )*lda
316*
317* Set row and column IZERO of A to 0.
318*
319 IF( iuplo.EQ.1 ) THEN
320 DO 20 i = 1, izero - 1
321 a( ioff+i ) = zero
322 20 CONTINUE
323 ioff = ioff + izero
324 DO 30 i = izero, n
325 a( ioff ) = zero
326 ioff = ioff + lda
327 30 CONTINUE
328 ELSE
329 ioff = izero
330 DO 40 i = 1, izero - 1
331 a( ioff ) = zero
332 ioff = ioff + lda
333 40 CONTINUE
334 ioff = ioff - izero
335 DO 50 i = izero, n
336 a( ioff+i ) = zero
337 50 CONTINUE
338 END IF
339 ELSE
340 izero = 0
341 END IF
342*
343* Set the imaginary part of the diagonals.
344*
345 CALL zlaipd( n, a, lda+1, 0 )
346*
347* Save a copy of the matrix A in ASAV.
348*
349 CALL zlacpy( uplo, n, n, a, lda, asav, lda )
350*
351 DO 100 iequed = 1, 2
352 equed = equeds( iequed )
353 IF( iequed.EQ.1 ) THEN
354 nfact = 3
355 ELSE
356 nfact = 1
357 END IF
358*
359 DO 90 ifact = 1, nfact
360 fact = facts( ifact )
361 prefac = lsame( fact, 'F' )
362 nofact = lsame( fact, 'N' )
363 equil = lsame( fact, 'E' )
364*
365 IF( zerot ) THEN
366 IF( prefac )
367 $ GO TO 90
368 rcondc = zero
369*
370 ELSE IF( .NOT.lsame( fact, 'N' ) ) THEN
371*
372* Compute the condition number for comparison with
373* the value returned by ZPOSVX (FACT = 'N' reuses
374* the condition number from the previous iteration
375* with FACT = 'F').
376*
377 CALL zlacpy( uplo, n, n, asav, lda, afac, lda )
378 IF( equil .OR. iequed.GT.1 ) THEN
379*
380* Compute row and column scale factors to
381* equilibrate the matrix A.
382*
383 CALL zpoequ( n, afac, lda, s, scond, amax,
384 $ info )
385 IF( info.EQ.0 .AND. n.GT.0 ) THEN
386 IF( iequed.GT.1 )
387 $ scond = zero
388*
389* Equilibrate the matrix.
390*
391 CALL zlaqhe( uplo, n, afac, lda, s, scond,
392 $ amax, equed )
393 END IF
394 END IF
395*
396* Save the condition number of the
397* non-equilibrated system for use in ZGET04.
398*
399 IF( equil )
400 $ roldc = rcondc
401*
402* Compute the 1-norm of A.
403*
404 anorm = zlanhe( '1', uplo, n, afac, lda, rwork )
405*
406* Factor the matrix A.
407*
408 CALL zpotrf( uplo, n, afac, lda, info )
409*
410* Form the inverse of A.
411*
412 CALL zlacpy( uplo, n, n, afac, lda, a, lda )
413 CALL zpotri( uplo, n, a, lda, info )
414*
415* Compute the 1-norm condition number of A.
416*
417 ainvnm = zlanhe( '1', uplo, n, a, lda, rwork )
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* Restore the matrix A.
426*
427 CALL zlacpy( uplo, n, n, asav, lda, a, lda )
428*
429* Form an exact solution and set the right hand side.
430*
431 srnamt = 'ZLARHS'
432 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
433 $ nrhs, a, lda, xact, lda, b, lda,
434 $ iseed, info )
435 xtype = 'C'
436 CALL zlacpy( 'Full', n, nrhs, b, lda, bsav, lda )
437*
438 IF( nofact ) THEN
439*
440* --- Test ZPOSV ---
441*
442* Compute the L*L' or U'*U factorization of the
443* matrix and solve the system.
444*
445 CALL zlacpy( uplo, n, n, a, lda, afac, lda )
446 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
447*
448 srnamt = 'ZPOSV '
449 CALL zposv( uplo, n, nrhs, afac, lda, x, lda,
450 $ info )
451*
452* Check error code from ZPOSV .
453*
454 IF( info.NE.izero ) THEN
455 CALL alaerh( path, 'ZPOSV ', info, izero,
456 $ uplo, n, n, -1, -1, nrhs, imat,
457 $ nfail, nerrs, nout )
458 GO TO 70
459 ELSE IF( info.NE.0 ) THEN
460 GO TO 70
461 END IF
462*
463* Reconstruct matrix from factors and compute
464* residual.
465*
466 CALL zpot01( uplo, n, a, lda, afac, lda, rwork,
467 $ result( 1 ) )
468*
469* Compute residual of the computed solution.
470*
471 CALL zlacpy( 'Full', n, nrhs, b, lda, work,
472 $ lda )
473 CALL zpot02( uplo, n, nrhs, a, lda, x, lda,
474 $ work, lda, rwork, result( 2 ) )
475*
476* Check solution from generated exact solution.
477*
478 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
479 $ result( 3 ) )
480 nt = 3
481*
482* Print information about the tests that did not
483* pass the threshold.
484*
485 DO 60 k = 1, nt
486 IF( result( k ).GE.thresh ) THEN
487 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
488 $ CALL aladhd( nout, path )
489 WRITE( nout, fmt = 9999 )'ZPOSV ', uplo,
490 $ n, imat, k, result( k )
491 nfail = nfail + 1
492 END IF
493 60 CONTINUE
494 nrun = nrun + nt
495 70 CONTINUE
496 END IF
497*
498* --- Test ZPOSVX ---
499*
500 IF( .NOT.prefac )
501 $ CALL zlaset( uplo, n, n, dcmplx( zero ),
502 $ dcmplx( zero ), afac, lda )
503 CALL zlaset( 'Full', n, nrhs, dcmplx( zero ),
504 $ dcmplx( zero ), x, lda )
505 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
506*
507* Equilibrate the matrix if FACT='F' and
508* EQUED='Y'.
509*
510 CALL zlaqhe( uplo, n, a, lda, s, scond, amax,
511 $ equed )
512 END IF
513*
514* Solve the system and compute the condition number
515* and error bounds using ZPOSVX.
516*
517 srnamt = 'ZPOSVX'
518 CALL zposvx( fact, uplo, n, nrhs, a, lda, afac,
519 $ lda, equed, s, b, lda, x, lda, rcond,
520 $ rwork, rwork( nrhs+1 ), work,
521 $ rwork( 2*nrhs+1 ), info )
522*
523* Check the error code from ZPOSVX.
524*
525 IF( info.NE.izero ) THEN
526 CALL alaerh( path, 'ZPOSVX', info, izero,
527 $ fact // uplo, n, n, -1, -1, nrhs,
528 $ imat, nfail, nerrs, nout )
529 GO TO 90
530 END IF
531*
532 IF( info.EQ.0 ) THEN
533 IF( .NOT.prefac ) THEN
534*
535* Reconstruct matrix from factors and compute
536* residual.
537*
538 CALL zpot01( uplo, n, a, lda, afac, lda,
539 $ rwork( 2*nrhs+1 ), result( 1 ) )
540 k1 = 1
541 ELSE
542 k1 = 2
543 END IF
544*
545* Compute residual of the computed solution.
546*
547 CALL zlacpy( 'Full', n, nrhs, bsav, lda, work,
548 $ lda )
549 CALL zpot02( uplo, n, nrhs, asav, lda, x, lda,
550 $ work, lda, rwork( 2*nrhs+1 ),
551 $ result( 2 ) )
552*
553* Check solution from generated exact solution.
554*
555 IF( nofact .OR. ( prefac .AND. lsame( equed,
556 $ 'N' ) ) ) THEN
557 CALL zget04( n, nrhs, x, lda, xact, lda,
558 $ rcondc, result( 3 ) )
559 ELSE
560 CALL zget04( n, nrhs, x, lda, xact, lda,
561 $ roldc, result( 3 ) )
562 END IF
563*
564* Check the error bounds from iterative
565* refinement.
566*
567 CALL zpot05( uplo, n, nrhs, asav, lda, b, lda,
568 $ x, lda, xact, lda, rwork,
569 $ rwork( nrhs+1 ), result( 4 ) )
570 ELSE
571 k1 = 6
572 END IF
573*
574* Compare RCOND from ZPOSVX with the computed value
575* in RCONDC.
576*
577 result( 6 ) = dget06( rcond, rcondc )
578*
579* Print information about the tests that did not pass
580* the threshold.
581*
582 DO 80 k = k1, 6
583 IF( result( k ).GE.thresh ) THEN
584 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
585 $ CALL aladhd( nout, path )
586 IF( prefac ) THEN
587 WRITE( nout, fmt = 9997 )'ZPOSVX', fact,
588 $ uplo, n, equed, imat, k, result( k )
589 ELSE
590 WRITE( nout, fmt = 9998 )'ZPOSVX', fact,
591 $ uplo, n, imat, k, result( k )
592 END IF
593 nfail = nfail + 1
594 END IF
595 80 CONTINUE
596 nrun = nrun + 7 - k1
597 90 CONTINUE
598 100 CONTINUE
599 110 CONTINUE
600 120 CONTINUE
601 130 CONTINUE
602*
603* Print a summary of the results.
604*
605 CALL alasvm( path, nout, nfail, nrun, nerrs )
606*
607 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i1,
608 $ ', test(', i1, ')=', g12.5 )
609 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
610 $ ', type ', i1, ', test(', i1, ')=', g12.5 )
611 9997 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
612 $ ', EQUED=''', a1, ''', type ', i1, ', test(', i1, ') =',
613 $ g12.5 )
614 RETURN
615*
616* End of ZDRVPO
617*
618 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
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 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
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 zlaqhe(uplo, n, a, lda, s, scond, amax, equed)
ZLAQHE scales a Hermitian matrix.
Definition zlaqhe.f:134
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 zpoequ(n, a, lda, s, scond, amax, info)
ZPOEQU
Definition zpoequ.f:113
subroutine zposv(uplo, n, nrhs, a, lda, b, ldb, info)
ZPOSV computes the solution to system of linear equations A * X = B for PO matrices
Definition zposv.f:130
subroutine zposvx(fact, uplo, n, nrhs, a, lda, af, ldaf, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
ZPOSVX computes the solution to system of linear equations A * X = B for PO matrices
Definition zposvx.f:306
subroutine zpotrf(uplo, n, a, lda, info)
ZPOTRF
Definition zpotrf.f:107
subroutine zpotri(uplo, n, a, lda, info)
ZPOTRI
Definition zpotri.f:95
subroutine zdrvpo(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, asav, b, bsav, x, xact, s, work, rwork, nout)
ZDRVPO
Definition zdrvpo.f:159
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 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 zpot01(uplo, n, a, lda, afac, ldafac, rwork, resid)
ZPOT01
Definition zpot01.f:106
subroutine zpot02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
ZPOT02
Definition zpot02.f:127
subroutine zpot05(uplo, n, nrhs, a, lda, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZPOT05
Definition zpot05.f:165