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

◆ zchkps()

subroutine zchkps ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nnb,
integer, dimension( * )  nbval,
integer  nrank,
integer, dimension( * )  rankval,
double precision  thresh,
logical  tsterr,
integer  nmax,
complex*16, dimension( * )  a,
complex*16, dimension( * )  afac,
complex*16, dimension( * )  perm,
integer, dimension( * )  piv,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork,
integer  nout 
)

ZCHKPS

Purpose:
 ZCHKPS tests ZPSTRF.
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]NNB
          NNB is INTEGER
          The number of values of NB contained in the vector NBVAL.
[in]NBVAL
          NBVAL is INTEGER array, dimension (NNB)
          The values of the block size NB.
[in]NRANK
          NRANK is INTEGER
          The number of values of RANK contained in the vector RANKVAL.
[in]RANKVAL
          RANKVAL is INTEGER array, dimension (NBVAL)
          The values of the block size NB.
[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.
[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]PERM
          PERM is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]PIV
          PIV is INTEGER array, dimension (NMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension (NMAX*3)
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (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 151 of file zchkps.f.

154*
155* -- LAPACK test routine --
156* -- LAPACK is a software package provided by Univ. of Tennessee, --
157* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
158*
159* .. Scalar Arguments ..
160 DOUBLE PRECISION THRESH
161 INTEGER NMAX, NN, NNB, NOUT, NRANK
162 LOGICAL TSTERR
163* ..
164* .. Array Arguments ..
165 COMPLEX*16 A( * ), AFAC( * ), PERM( * ), WORK( * )
166 DOUBLE PRECISION RWORK( * )
167 INTEGER NBVAL( * ), NVAL( * ), PIV( * ), RANKVAL( * )
168 LOGICAL DOTYPE( * )
169* ..
170*
171* =====================================================================
172*
173* .. Parameters ..
174 DOUBLE PRECISION ONE
175 parameter( one = 1.0e+0 )
176 INTEGER NTYPES
177 parameter( ntypes = 9 )
178* ..
179* .. Local Scalars ..
180 DOUBLE PRECISION ANORM, CNDNUM, RESULT, TOL
181 INTEGER COMPRANK, I, IMAT, IN, INB, INFO, IRANK, IUPLO,
182 $ IZERO, KL, KU, LDA, MODE, N, NB, NERRS, NFAIL,
183 $ NIMAT, NRUN, RANK, RANKDIFF
184 CHARACTER DIST, TYPE, UPLO
185 CHARACTER*3 PATH
186* ..
187* .. Local Arrays ..
188 INTEGER ISEED( 4 ), ISEEDY( 4 )
189 CHARACTER UPLOS( 2 )
190* ..
191* .. External Subroutines ..
192 EXTERNAL alaerh, alahd, alasum, xlaenv, zerrps, zlacpy,
194* ..
195* .. Scalars in Common ..
196 INTEGER INFOT, NUNIT
197 LOGICAL LERR, OK
198 CHARACTER*32 SRNAMT
199* ..
200* .. Common blocks ..
201 COMMON / infoc / infot, nunit, ok, lerr
202 COMMON / srnamc / srnamt
203* ..
204* .. Intrinsic Functions ..
205 INTRINSIC dble, max, ceiling
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 path( 1: 1 ) = 'Zomplex Precision'
216 path( 2: 3 ) = 'PS'
217 nrun = 0
218 nfail = 0
219 nerrs = 0
220 DO 100 i = 1, 4
221 iseed( i ) = iseedy( i )
222 100 CONTINUE
223*
224* Test the error exits
225*
226 IF( tsterr )
227 $ CALL zerrps( path, nout )
228 infot = 0
229*
230* Do for each value of N in NVAL
231*
232 DO 150 in = 1, nn
233 n = nval( in )
234 lda = max( n, 1 )
235 nimat = ntypes
236 IF( n.LE.0 )
237 $ nimat = 1
238*
239 izero = 0
240 DO 140 imat = 1, nimat
241*
242* Do the tests only if DOTYPE( IMAT ) is true.
243*
244 IF( .NOT.dotype( imat ) )
245 $ GO TO 140
246*
247* Do for each value of RANK in RANKVAL
248*
249 DO 130 irank = 1, nrank
250*
251* Only repeat test 3 to 5 for different ranks
252* Other tests use full rank
253*
254 IF( ( imat.LT.3 .OR. imat.GT.5 ) .AND. irank.GT.1 )
255 $ GO TO 130
256*
257 rank = ceiling( ( n * dble( rankval( irank ) ) )
258 $ / 100.e+0 )
259*
260*
261* Do first for UPLO = 'U', then for UPLO = 'L'
262*
263 DO 120 iuplo = 1, 2
264 uplo = uplos( iuplo )
265*
266* Set up parameters with ZLATB5 and generate a test matrix
267* with ZLATMT.
268*
269 CALL zlatb5( path, imat, n, TYPE, KL, KU, ANORM,
270 $ MODE, CNDNUM, DIST )
271*
272 srnamt = 'ZLATMT'
273 CALL zlatmt( n, n, dist, iseed, TYPE, RWORK, MODE,
274 $ CNDNUM, ANORM, RANK, KL, KU, UPLO, A,
275 $ LDA, WORK, INFO )
276*
277* Check error code from ZLATMT.
278*
279 IF( info.NE.0 ) THEN
280 CALL alaerh( path, 'ZLATMT', info, 0, uplo, n,
281 $ n, -1, -1, -1, imat, nfail, nerrs,
282 $ nout )
283 GO TO 120
284 END IF
285*
286* Do for each value of NB in NBVAL
287*
288 DO 110 inb = 1, nnb
289 nb = nbval( inb )
290 CALL xlaenv( 1, nb )
291*
292* Compute the pivoted L*L' or U'*U factorization
293* of the matrix.
294*
295 CALL zlacpy( uplo, n, n, a, lda, afac, lda )
296 srnamt = 'ZPSTRF'
297*
298* Use default tolerance
299*
300 tol = -one
301 CALL zpstrf( uplo, n, afac, lda, piv, comprank,
302 $ tol, rwork, info )
303*
304* Check error code from ZPSTRF.
305*
306 IF( (info.LT.izero)
307 $ .OR.(info.NE.izero.AND.rank.EQ.n)
308 $ .OR.(info.LE.izero.AND.rank.LT.n) ) THEN
309 CALL alaerh( path, 'ZPSTRF', info, izero,
310 $ uplo, n, n, -1, -1, nb, imat,
311 $ nfail, nerrs, nout )
312 GO TO 110
313 END IF
314*
315* Skip the test if INFO is not 0.
316*
317 IF( info.NE.0 )
318 $ GO TO 110
319*
320* Reconstruct matrix from factors and compute residual.
321*
322* PERM holds permuted L*L^T or U^T*U
323*
324 CALL zpst01( uplo, n, a, lda, afac, lda, perm, lda,
325 $ piv, rwork, result, comprank )
326*
327* Print information about the tests that did not pass
328* the threshold or where computed rank was not RANK.
329*
330 IF( n.EQ.0 )
331 $ comprank = 0
332 rankdiff = rank - comprank
333 IF( result.GE.thresh ) THEN
334 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
335 $ CALL alahd( nout, path )
336 WRITE( nout, fmt = 9999 )uplo, n, rank,
337 $ rankdiff, nb, imat, result
338 nfail = nfail + 1
339 END IF
340 nrun = nrun + 1
341 110 CONTINUE
342*
343 120 CONTINUE
344 130 CONTINUE
345 140 CONTINUE
346 150 CONTINUE
347*
348* Print a summary of the results.
349*
350 CALL alasum( path, nout, nfail, nrun, nerrs )
351*
352 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', RANK =', i3,
353 $ ', Diff =', i5, ', NB =', i4, ', type ', i2, ', Ratio =',
354 $ g12.5 )
355 RETURN
356*
357* End of ZCHKPS
358*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
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 zpstrf(uplo, n, a, lda, piv, rank, tol, work, info)
ZPSTRF computes the Cholesky factorization with complete pivoting of a complex Hermitian positive sem...
Definition zpstrf.f:142
subroutine zerrps(path, nunit)
ZERRPS
Definition zerrps.f:55
subroutine zlatb5(path, imat, n, type, kl, ku, anorm, mode, cndnum, dist)
ZLATB5
Definition zlatb5.f:114
subroutine zlatmt(m, n, dist, iseed, sym, d, mode, cond, dmax, rank, kl, ku, pack, a, lda, work, info)
ZLATMT
Definition zlatmt.f:340
subroutine zpst01(uplo, n, a, lda, afac, ldafac, perm, ldperm, piv, rwork, resid, rank)
ZPST01
Definition zpst01.f:136
Here is the call graph for this function:
Here is the caller graph for this function: