LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ddrvac.f
Go to the documentation of this file.
1*> \brief \b DDRVAC
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 DDRVAC( DOTYPE, NM, MVAL, NNS, NSVAL, THRESH, NMAX,
12* A, AFAC, B, X, WORK,
13* RWORK, SWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* INTEGER NMAX, NM, NNS, NOUT
17* DOUBLE PRECISION THRESH
18* ..
19* .. Array Arguments ..
20* LOGICAL DOTYPE( * )
21* INTEGER MVAL( * ), NSVAL( * )
22* REAL SWORK(*)
23* DOUBLE PRECISION A( * ), AFAC( * ), B( * ),
24* $ RWORK( * ), WORK( * ), X( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*>
33*> DDRVAC tests DSPOSV.
34*> \endverbatim
35*
36* Arguments:
37* ==========
38*
39*> \param[in] DOTYPE
40*> \verbatim
41*> DOTYPE is LOGICAL array, dimension (NTYPES)
42*> The matrix types to be used for testing. Matrices of type j
43*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
44*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
45*> \endverbatim
46*>
47*> \param[in] NM
48*> \verbatim
49*> NM is INTEGER
50*> The number of values of N contained in the vector MVAL.
51*> \endverbatim
52*>
53*> \param[in] MVAL
54*> \verbatim
55*> MVAL is INTEGER array, dimension (NM)
56*> The values of the matrix dimension N.
57*> \endverbatim
58*>
59*> \param[in] NNS
60*> \verbatim
61*> NNS is INTEGER
62*> The number of values of NRHS contained in the vector NSVAL.
63*> \endverbatim
64*>
65*> \param[in] NSVAL
66*> \verbatim
67*> NSVAL is INTEGER array, dimension (NNS)
68*> The values of the number of right hand sides NRHS.
69*> \endverbatim
70*>
71*> \param[in] THRESH
72*> \verbatim
73*> THRESH is DOUBLE PRECISION
74*> The threshold value for the test ratios. A result is
75*> included in the output file if RESULT >= THRESH. To have
76*> every test ratio printed, use THRESH = 0.
77*> \endverbatim
78*>
79*> \param[in] NMAX
80*> \verbatim
81*> NMAX is INTEGER
82*> The maximum value permitted for N, used in dimensioning the
83*> work arrays.
84*> \endverbatim
85*>
86*> \param[out] A
87*> \verbatim
88*> A is DOUBLE PRECISION array, dimension (NMAX*NMAX)
89*> \endverbatim
90*>
91*> \param[out] AFAC
92*> \verbatim
93*> AFAC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
94*> \endverbatim
95*>
96*> \param[out] B
97*> \verbatim
98*> B is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
99*> \endverbatim
100*>
101*> \param[out] X
102*> \verbatim
103*> X is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
104*> \endverbatim
105*>
106*> \param[out] WORK
107*> \verbatim
108*> WORK is DOUBLE PRECISION array, dimension
109*> (NMAX*max(3,NSMAX))
110*> \endverbatim
111*>
112*> \param[out] RWORK
113*> \verbatim
114*> RWORK is DOUBLE PRECISION array, dimension
115*> (max(2*NMAX,2*NSMAX+NWORK))
116*> \endverbatim
117*>
118*> \param[out] SWORK
119*> \verbatim
120*> SWORK is REAL array, dimension
121*> (NMAX*(NSMAX+NMAX))
122*> \endverbatim
123*>
124*> \param[in] NOUT
125*> \verbatim
126*> NOUT is INTEGER
127*> The unit number for output.
128*> \endverbatim
129*
130* Authors:
131* ========
132*
133*> \author Univ. of Tennessee
134*> \author Univ. of California Berkeley
135*> \author Univ. of Colorado Denver
136*> \author NAG Ltd.
137*
138*> \ingroup double_lin
139*
140* =====================================================================
141 SUBROUTINE ddrvac( DOTYPE, NM, MVAL, NNS, NSVAL, THRESH, NMAX,
142 $ A, AFAC, B, X, WORK,
143 $ RWORK, SWORK, NOUT )
144*
145* -- LAPACK test routine --
146* -- LAPACK is a software package provided by Univ. of Tennessee, --
147* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
148*
149* .. Scalar Arguments ..
150 INTEGER NMAX, NM, NNS, NOUT
151 DOUBLE PRECISION THRESH
152* ..
153* .. Array Arguments ..
154 LOGICAL DOTYPE( * )
155 INTEGER MVAL( * ), NSVAL( * )
156 REAL SWORK(*)
157 DOUBLE PRECISION A( * ), AFAC( * ), B( * ),
158 $ rwork( * ), work( * ), x( * )
159* ..
160*
161* =====================================================================
162*
163* .. Parameters ..
164 DOUBLE PRECISION ZERO
165 PARAMETER ( ZERO = 0.0d+0 )
166 INTEGER NTYPES
167 parameter( ntypes = 9 )
168 INTEGER NTESTS
169 parameter( ntests = 1 )
170* ..
171* .. Local Scalars ..
172 LOGICAL ZEROT
173 CHARACTER DIST, TYPE, UPLO, XTYPE
174 CHARACTER*3 PATH
175 INTEGER I, IM, IMAT, INFO, IOFF, IRHS, IUPLO,
176 $ izero, kl, ku, lda, mode, n,
177 $ nerrs, nfail, nimat, nrhs, nrun
178 DOUBLE PRECISION ANORM, CNDNUM
179* ..
180* .. Local Arrays ..
181 CHARACTER UPLOS( 2 )
182 INTEGER ISEED( 4 ), ISEEDY( 4 )
183 DOUBLE PRECISION RESULT( NTESTS )
184* ..
185* .. Local Variables ..
186 INTEGER ITER, KASE
187* ..
188* .. External Functions ..
189 LOGICAL LSAME
190 EXTERNAL LSAME
191* ..
192* .. External Subroutines ..
193 EXTERNAL alaerh, dlacpy,
195 $ dpot06, dsposv
196* ..
197* .. Intrinsic Functions ..
198 INTRINSIC dble, max, sqrt
199* ..
200* .. Scalars in Common ..
201 LOGICAL LERR, OK
202 CHARACTER*32 SRNAMT
203 INTEGER INFOT, NUNIT
204* ..
205* .. Common blocks ..
206 COMMON / infoc / infot, nunit, ok, lerr
207 COMMON / srnamc / srnamt
208* ..
209* .. Data statements ..
210 DATA iseedy / 1988, 1989, 1990, 1991 /
211 DATA uplos / 'U', 'L' /
212* ..
213* .. Executable Statements ..
214*
215* Initialize constants and the random number seed.
216*
217 kase = 0
218 path( 1: 1 ) = 'Double precision'
219 path( 2: 3 ) = 'PO'
220 nrun = 0
221 nfail = 0
222 nerrs = 0
223 DO 10 i = 1, 4
224 iseed( i ) = iseedy( i )
225 10 CONTINUE
226*
227 infot = 0
228*
229* Do for each value of N in MVAL
230*
231 DO 120 im = 1, nm
232 n = mval( im )
233 lda = max( n, 1 )
234 nimat = ntypes
235 IF( n.LE.0 )
236 $ nimat = 1
237*
238 DO 110 imat = 1, nimat
239*
240* Do the tests only if DOTYPE( IMAT ) is true.
241*
242 IF( .NOT.dotype( imat ) )
243 $ GO TO 110
244*
245* Skip types 3, 4, or 5 if the matrix size is too small.
246*
247 zerot = imat.GE.3 .AND. imat.LE.5
248 IF( zerot .AND. n.LT.imat-2 )
249 $ GO TO 110
250*
251* Do first for UPLO = 'U', then for UPLO = 'L'
252*
253 DO 100 iuplo = 1, 2
254 uplo = uplos( iuplo )
255*
256* Set up parameters with DLATB4 and generate a test matrix
257* with DLATMS.
258*
259 CALL dlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
260 $ cndnum, dist )
261*
262 srnamt = 'DLATMS'
263 CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
264 $ cndnum, anorm, kl, ku, uplo, a, lda, work,
265 $ info )
266*
267* Check error code from DLATMS.
268*
269 IF( info.NE.0 ) THEN
270 CALL alaerh( path, 'DLATMS', info, 0, uplo, n, n, -1,
271 $ -1, -1, imat, nfail, nerrs, nout )
272 GO TO 100
273 END IF
274*
275* For types 3-5, zero one row and column of the matrix to
276* test that INFO is returned correctly.
277*
278 IF( zerot ) THEN
279 IF( imat.EQ.3 ) THEN
280 izero = 1
281 ELSE IF( imat.EQ.4 ) THEN
282 izero = n
283 ELSE
284 izero = n / 2 + 1
285 END IF
286 ioff = ( izero-1 )*lda
287*
288* Set row and column IZERO of A to 0.
289*
290 IF( iuplo.EQ.1 ) THEN
291 DO 20 i = 1, izero - 1
292 a( ioff+i ) = zero
293 20 CONTINUE
294 ioff = ioff + izero
295 DO 30 i = izero, n
296 a( ioff ) = zero
297 ioff = ioff + lda
298 30 CONTINUE
299 ELSE
300 ioff = izero
301 DO 40 i = 1, izero - 1
302 a( ioff ) = zero
303 ioff = ioff + lda
304 40 CONTINUE
305 ioff = ioff - izero
306 DO 50 i = izero, n
307 a( ioff+i ) = zero
308 50 CONTINUE
309 END IF
310 ELSE
311 izero = 0
312 END IF
313*
314 DO 60 irhs = 1, nns
315 nrhs = nsval( irhs )
316 xtype = 'N'
317*
318* Form an exact solution and set the right hand side.
319*
320 srnamt = 'DLARHS'
321 CALL dlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
322 $ nrhs, a, lda, x, lda, b, lda,
323 $ iseed, info )
324*
325* Compute the L*L' or U'*U factorization of the
326* matrix and solve the system.
327*
328 srnamt = 'DSPOSV '
329 kase = kase + 1
330*
331 CALL dlacpy( 'All', n, n, a, lda, afac, lda)
332*
333 CALL dsposv( uplo, n, nrhs, afac, lda, b, lda, x, lda,
334 $ work, swork, iter, info )
335
336 IF (iter.LT.0) THEN
337 CALL dlacpy( 'All', n, n, a, lda, afac, lda )
338 ENDIF
339*
340* Check error code from DSPOSV .
341*
342 IF( info.NE.izero ) THEN
343*
344 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
345 $ CALL alahd( nout, path )
346 nerrs = nerrs + 1
347*
348 IF( info.NE.izero .AND. izero.NE.0 ) THEN
349 WRITE( nout, fmt = 9988 )'DSPOSV',info,izero,n,
350 $ imat
351 ELSE
352 WRITE( nout, fmt = 9975 )'DSPOSV',info,n,imat
353 END IF
354 END IF
355*
356* Skip the remaining test if the matrix is singular.
357*
358 IF( info.NE.0 )
359 $ GO TO 110
360*
361* Check the quality of the solution
362*
363 CALL dlacpy( 'All', n, nrhs, b, lda, work, lda )
364*
365 CALL dpot06( uplo, n, nrhs, a, lda, x, lda, work,
366 $ lda, rwork, result( 1 ) )
367*
368* Check if the test passes the testing.
369* Print information about the tests that did not
370* pass the testing.
371*
372* If iterative refinement has been used and claimed to
373* be successful (ITER>0), we want
374* NORM1(B - A*X)/(NORM1(A)*NORM1(X)*EPS*SRQT(N)) < 1
375*
376* If double precision has been used (ITER<0), we want
377* NORM1(B - A*X)/(NORM1(A)*NORM1(X)*EPS) < THRES
378* (Cf. the linear solver testing routines)
379*
380 IF ((thresh.LE.0.0e+00)
381 $ .OR.((iter.GE.0).AND.(n.GT.0)
382 $ .AND.(result(1).GE.sqrt(dble(n))))
383 $ .OR.((iter.LT.0).AND.(result(1).GE.thresh))) THEN
384*
385 IF( nfail.EQ.0 .AND. nerrs.EQ.0 ) THEN
386 WRITE( nout, fmt = 8999 )'DPO'
387 WRITE( nout, fmt = '( '' Matrix types:'' )' )
388 WRITE( nout, fmt = 8979 )
389 WRITE( nout, fmt = '( '' Test ratios:'' )' )
390 WRITE( nout, fmt = 8960 )1
391 WRITE( nout, fmt = '( '' Messages:'' )' )
392 END IF
393*
394 WRITE( nout, fmt = 9998 )uplo, n, nrhs, imat, 1,
395 $ result( 1 )
396*
397 nfail = nfail + 1
398*
399 END IF
400*
401 nrun = nrun + 1
402*
403 60 CONTINUE
404 100 CONTINUE
405 110 CONTINUE
406 120 CONTINUE
407*
408* Print a summary of the results.
409*
410 IF( nfail.GT.0 ) THEN
411 WRITE( nout, fmt = 9996 )'DSPOSV', nfail, nrun
412 ELSE
413 WRITE( nout, fmt = 9995 )'DSPOSV', nrun
414 END IF
415 IF( nerrs.GT.0 ) THEN
416 WRITE( nout, fmt = 9994 )nerrs
417 END IF
418*
419 9998 FORMAT( ' UPLO=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
420 $ i2, ', test(', i2, ') =', g12.5 )
421 9996 FORMAT( 1x, a6, ': ', i6, ' out of ', i6,
422 $ ' tests failed to pass the threshold' )
423 9995 FORMAT( /1x, 'All tests for ', a6,
424 $ ' routines passed the threshold ( ', i6, ' tests run)' )
425 9994 FORMAT( 6x, i6, ' error messages recorded' )
426*
427* SUBNAM, INFO, INFOE, N, IMAT
428*
429 9988 FORMAT( ' *** ', a6, ' returned with INFO =', i5, ' instead of ',
430 $ i5, / ' ==> N =', i5, ', type ',
431 $ i2 )
432*
433* SUBNAM, INFO, N, IMAT
434*
435 9975 FORMAT( ' *** Error code from ', a6, '=', i5, ' for M=', i5,
436 $ ', type ', i2 )
437 8999 FORMAT( / 1x, a3, ': positive definite dense matrices' )
438 8979 FORMAT( 4x, '1. Diagonal', 24x, '7. Last n/2 columns zero', / 4x,
439 $ '2. Upper triangular', 16x,
440 $ '8. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
441 $ '3. Lower triangular', 16x, '9. Random, CNDNUM = 0.1/EPS',
442 $ / 4x, '4. Random, CNDNUM = 2', 13x,
443 $ '10. Scaled near underflow', / 4x, '5. First column zero',
444 $ 14x, '11. Scaled near overflow', / 4x,
445 $ '6. Last column zero' )
446 8960 FORMAT( 3x, i2, ': norm_1( B - A * X ) / ',
447 $ '( norm_1(A) * norm_1(X) * EPS * SQRT(N) ) > 1 if ITERREF',
448 $ / 4x, 'or norm_1( B - A * X ) / ',
449 $ '( norm_1(A) * norm_1(X) * EPS ) > THRES if DPOTRF' )
450
451 RETURN
452*
453* End of DDRVAC
454*
455 END
subroutine dlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
DLARHS
Definition dlarhs.f:205
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 ddrvac(dotype, nm, mval, nns, nsval, thresh, nmax, a, afac, b, x, work, rwork, swork, nout)
DDRVAC
Definition ddrvac.f:144
subroutine dlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
DLATB4
Definition dlatb4.f:120
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
Definition dlatms.f:321
subroutine dpot06(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
DPOT06
Definition dpot06.f:127
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
subroutine dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition dlaset.f:110
subroutine dsposv(uplo, n, nrhs, a, lda, b, ldb, x, ldx, work, swork, iter, info)
DSPOSV computes the solution to system of linear equations A * X = B for PO matrices
Definition dsposv.f:199