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

◆ zerrge()

subroutine zerrge ( character*3  path,
integer  nunit 
)

ZERRGE

Purpose:
 ZERRGE tests the error exits for the COMPLEX*16 routines
 for general matrices.
Parameters
[in]PATH
          PATH is CHARACTER*3
          The LAPACK path name for the routines to be tested.
[in]NUNIT
          NUNIT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 54 of file zerrge.f.

55*
56* -- LAPACK test routine --
57* -- LAPACK is a software package provided by Univ. of Tennessee, --
58* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
59*
60* .. Scalar Arguments ..
61 CHARACTER*3 PATH
62 INTEGER NUNIT
63* ..
64*
65* =====================================================================
66*
67* .. Parameters ..
68 INTEGER NMAX
69 parameter( nmax = 4 )
70* ..
71* .. Local Scalars ..
72 CHARACTER*2 C2
73 INTEGER I, INFO, J
74 DOUBLE PRECISION ANRM, CCOND, RCOND
75* ..
76* .. Local Arrays ..
77 INTEGER IP( NMAX )
78 DOUBLE PRECISION R( NMAX ), R1( NMAX ), R2( NMAX )
79 COMPLEX*16 A( NMAX, NMAX ), AF( NMAX, NMAX ), B( NMAX ),
80 $ W( 2*NMAX ), X( NMAX )
81* ..
82* .. External Functions ..
83 LOGICAL LSAMEN
84 EXTERNAL lsamen
85* ..
86* .. External Subroutines ..
87 EXTERNAL alaesm, chkxer, zgbcon, zgbequ, zgbrfs, zgbtf2,
90* ..
91* .. Scalars in Common ..
92 LOGICAL LERR, OK
93 CHARACTER*32 SRNAMT
94 INTEGER INFOT, NOUT
95* ..
96* .. Common blocks ..
97 COMMON / infoc / infot, nout, ok, lerr
98 COMMON / srnamc / srnamt
99* ..
100* .. Intrinsic Functions ..
101 INTRINSIC dble, dcmplx
102* ..
103* .. Executable Statements ..
104*
105 nout = nunit
106 WRITE( nout, fmt = * )
107 c2 = path( 2: 3 )
108*
109* Set the variables to innocuous values.
110*
111 DO 20 j = 1, nmax
112 DO 10 i = 1, nmax
113 a( i, j ) = dcmplx( 1.d0 / dble( i+j ),
114 $ -1.d0 / dble( i+j ) )
115 af( i, j ) = dcmplx( 1.d0 / dble( i+j ),
116 $ -1.d0 / dble( i+j ) )
117 10 CONTINUE
118 b( j ) = 0.d0
119 r1( j ) = 0.d0
120 r2( j ) = 0.d0
121 w( j ) = 0.d0
122 x( j ) = 0.d0
123 ip( j ) = j
124 20 CONTINUE
125 ok = .true.
126*
127* Test error exits of the routines that use the LU decomposition
128* of a general matrix.
129*
130 IF( lsamen( 2, c2, 'GE' ) ) THEN
131*
132* ZGETRF
133*
134 srnamt = 'ZGETRF'
135 infot = 1
136 CALL zgetrf( -1, 0, a, 1, ip, info )
137 CALL chkxer( 'ZGETRF', infot, nout, lerr, ok )
138 infot = 2
139 CALL zgetrf( 0, -1, a, 1, ip, info )
140 CALL chkxer( 'ZGETRF', infot, nout, lerr, ok )
141 infot = 4
142 CALL zgetrf( 2, 1, a, 1, ip, info )
143 CALL chkxer( 'ZGETRF', infot, nout, lerr, ok )
144*
145* ZGETF2
146*
147 srnamt = 'ZGETF2'
148 infot = 1
149 CALL zgetf2( -1, 0, a, 1, ip, info )
150 CALL chkxer( 'ZGETF2', infot, nout, lerr, ok )
151 infot = 2
152 CALL zgetf2( 0, -1, a, 1, ip, info )
153 CALL chkxer( 'ZGETF2', infot, nout, lerr, ok )
154 infot = 4
155 CALL zgetf2( 2, 1, a, 1, ip, info )
156 CALL chkxer( 'ZGETF2', infot, nout, lerr, ok )
157*
158* ZGETRI
159*
160 srnamt = 'ZGETRI'
161 infot = 1
162 CALL zgetri( -1, a, 1, ip, w, 1, info )
163 CALL chkxer( 'ZGETRI', infot, nout, lerr, ok )
164 infot = 3
165 CALL zgetri( 2, a, 1, ip, w, 2, info )
166 CALL chkxer( 'ZGETRI', infot, nout, lerr, ok )
167 infot = 6
168 CALL zgetri( 2, a, 2, ip, w, 1, info )
169 CALL chkxer( 'ZGETRI', infot, nout, lerr, ok )
170*
171* ZGETRS
172*
173 srnamt = 'ZGETRS'
174 infot = 1
175 CALL zgetrs( '/', 0, 0, a, 1, ip, b, 1, info )
176 CALL chkxer( 'ZGETRS', infot, nout, lerr, ok )
177 infot = 2
178 CALL zgetrs( 'N', -1, 0, a, 1, ip, b, 1, info )
179 CALL chkxer( 'ZGETRS', infot, nout, lerr, ok )
180 infot = 3
181 CALL zgetrs( 'N', 0, -1, a, 1, ip, b, 1, info )
182 CALL chkxer( 'ZGETRS', infot, nout, lerr, ok )
183 infot = 5
184 CALL zgetrs( 'N', 2, 1, a, 1, ip, b, 2, info )
185 CALL chkxer( 'ZGETRS', infot, nout, lerr, ok )
186 infot = 8
187 CALL zgetrs( 'N', 2, 1, a, 2, ip, b, 1, info )
188 CALL chkxer( 'ZGETRS', infot, nout, lerr, ok )
189*
190* ZGERFS
191*
192 srnamt = 'ZGERFS'
193 infot = 1
194 CALL zgerfs( '/', 0, 0, a, 1, af, 1, ip, b, 1, x, 1, r1, r2, w,
195 $ r, info )
196 CALL chkxer( 'ZGERFS', infot, nout, lerr, ok )
197 infot = 2
198 CALL zgerfs( 'N', -1, 0, a, 1, af, 1, ip, b, 1, x, 1, r1, r2,
199 $ w, r, info )
200 CALL chkxer( 'ZGERFS', infot, nout, lerr, ok )
201 infot = 3
202 CALL zgerfs( 'N', 0, -1, a, 1, af, 1, ip, b, 1, x, 1, r1, r2,
203 $ w, r, info )
204 CALL chkxer( 'ZGERFS', infot, nout, lerr, ok )
205 infot = 5
206 CALL zgerfs( 'N', 2, 1, a, 1, af, 2, ip, b, 2, x, 2, r1, r2, w,
207 $ r, info )
208 CALL chkxer( 'ZGERFS', infot, nout, lerr, ok )
209 infot = 7
210 CALL zgerfs( 'N', 2, 1, a, 2, af, 1, ip, b, 2, x, 2, r1, r2, w,
211 $ r, info )
212 CALL chkxer( 'ZGERFS', infot, nout, lerr, ok )
213 infot = 10
214 CALL zgerfs( 'N', 2, 1, a, 2, af, 2, ip, b, 1, x, 2, r1, r2, w,
215 $ r, info )
216 CALL chkxer( 'ZGERFS', infot, nout, lerr, ok )
217 infot = 12
218 CALL zgerfs( 'N', 2, 1, a, 2, af, 2, ip, b, 2, x, 1, r1, r2, w,
219 $ r, info )
220 CALL chkxer( 'ZGERFS', infot, nout, lerr, ok )
221*
222* ZGECON
223*
224 srnamt = 'ZGECON'
225 infot = 1
226 CALL zgecon( '/', 0, a, 1, anrm, rcond, w, r, info )
227 CALL chkxer( 'ZGECON', infot, nout, lerr, ok )
228 infot = 2
229 CALL zgecon( '1', -1, a, 1, anrm, rcond, w, r, info )
230 CALL chkxer( 'ZGECON', infot, nout, lerr, ok )
231 infot = 4
232 CALL zgecon( '1', 2, a, 1, anrm, rcond, w, r, info )
233 CALL chkxer( 'ZGECON', infot, nout, lerr, ok )
234*
235* ZGEEQU
236*
237 srnamt = 'ZGEEQU'
238 infot = 1
239 CALL zgeequ( -1, 0, a, 1, r1, r2, rcond, ccond, anrm, info )
240 CALL chkxer( 'ZGEEQU', infot, nout, lerr, ok )
241 infot = 2
242 CALL zgeequ( 0, -1, a, 1, r1, r2, rcond, ccond, anrm, info )
243 CALL chkxer( 'ZGEEQU', infot, nout, lerr, ok )
244 infot = 4
245 CALL zgeequ( 2, 2, a, 1, r1, r2, rcond, ccond, anrm, info )
246 CALL chkxer( 'ZGEEQU', infot, nout, lerr, ok )
247*
248* Test error exits of the routines that use the LU decomposition
249* of a general band matrix.
250*
251 ELSE IF( lsamen( 2, c2, 'GB' ) ) THEN
252*
253* ZGBTRF
254*
255 srnamt = 'ZGBTRF'
256 infot = 1
257 CALL zgbtrf( -1, 0, 0, 0, a, 1, ip, info )
258 CALL chkxer( 'ZGBTRF', infot, nout, lerr, ok )
259 infot = 2
260 CALL zgbtrf( 0, -1, 0, 0, a, 1, ip, info )
261 CALL chkxer( 'ZGBTRF', infot, nout, lerr, ok )
262 infot = 3
263 CALL zgbtrf( 1, 1, -1, 0, a, 1, ip, info )
264 CALL chkxer( 'ZGBTRF', infot, nout, lerr, ok )
265 infot = 4
266 CALL zgbtrf( 1, 1, 0, -1, a, 1, ip, info )
267 CALL chkxer( 'ZGBTRF', infot, nout, lerr, ok )
268 infot = 6
269 CALL zgbtrf( 2, 2, 1, 1, a, 3, ip, info )
270 CALL chkxer( 'ZGBTRF', infot, nout, lerr, ok )
271*
272* ZGBTF2
273*
274 srnamt = 'ZGBTF2'
275 infot = 1
276 CALL zgbtf2( -1, 0, 0, 0, a, 1, ip, info )
277 CALL chkxer( 'ZGBTF2', infot, nout, lerr, ok )
278 infot = 2
279 CALL zgbtf2( 0, -1, 0, 0, a, 1, ip, info )
280 CALL chkxer( 'ZGBTF2', infot, nout, lerr, ok )
281 infot = 3
282 CALL zgbtf2( 1, 1, -1, 0, a, 1, ip, info )
283 CALL chkxer( 'ZGBTF2', infot, nout, lerr, ok )
284 infot = 4
285 CALL zgbtf2( 1, 1, 0, -1, a, 1, ip, info )
286 CALL chkxer( 'ZGBTF2', infot, nout, lerr, ok )
287 infot = 6
288 CALL zgbtf2( 2, 2, 1, 1, a, 3, ip, info )
289 CALL chkxer( 'ZGBTF2', infot, nout, lerr, ok )
290*
291* ZGBTRS
292*
293 srnamt = 'ZGBTRS'
294 infot = 1
295 CALL zgbtrs( '/', 0, 0, 0, 1, a, 1, ip, b, 1, info )
296 CALL chkxer( 'ZGBTRS', infot, nout, lerr, ok )
297 infot = 2
298 CALL zgbtrs( 'N', -1, 0, 0, 1, a, 1, ip, b, 1, info )
299 CALL chkxer( 'ZGBTRS', infot, nout, lerr, ok )
300 infot = 3
301 CALL zgbtrs( 'N', 1, -1, 0, 1, a, 1, ip, b, 1, info )
302 CALL chkxer( 'ZGBTRS', infot, nout, lerr, ok )
303 infot = 4
304 CALL zgbtrs( 'N', 1, 0, -1, 1, a, 1, ip, b, 1, info )
305 CALL chkxer( 'ZGBTRS', infot, nout, lerr, ok )
306 infot = 5
307 CALL zgbtrs( 'N', 1, 0, 0, -1, a, 1, ip, b, 1, info )
308 CALL chkxer( 'ZGBTRS', infot, nout, lerr, ok )
309 infot = 7
310 CALL zgbtrs( 'N', 2, 1, 1, 1, a, 3, ip, b, 2, info )
311 CALL chkxer( 'ZGBTRS', infot, nout, lerr, ok )
312 infot = 10
313 CALL zgbtrs( 'N', 2, 0, 0, 1, a, 1, ip, b, 1, info )
314 CALL chkxer( 'ZGBTRS', infot, nout, lerr, ok )
315*
316* ZGBRFS
317*
318 srnamt = 'ZGBRFS'
319 infot = 1
320 CALL zgbrfs( '/', 0, 0, 0, 0, a, 1, af, 1, ip, b, 1, x, 1, r1,
321 $ r2, w, r, info )
322 CALL chkxer( 'ZGBRFS', infot, nout, lerr, ok )
323 infot = 2
324 CALL zgbrfs( 'N', -1, 0, 0, 0, a, 1, af, 1, ip, b, 1, x, 1, r1,
325 $ r2, w, r, info )
326 CALL chkxer( 'ZGBRFS', infot, nout, lerr, ok )
327 infot = 3
328 CALL zgbrfs( 'N', 1, -1, 0, 0, a, 1, af, 1, ip, b, 1, x, 1, r1,
329 $ r2, w, r, info )
330 CALL chkxer( 'ZGBRFS', infot, nout, lerr, ok )
331 infot = 4
332 CALL zgbrfs( 'N', 1, 0, -1, 0, a, 1, af, 1, ip, b, 1, x, 1, r1,
333 $ r2, w, r, info )
334 CALL chkxer( 'ZGBRFS', infot, nout, lerr, ok )
335 infot = 5
336 CALL zgbrfs( 'N', 1, 0, 0, -1, a, 1, af, 1, ip, b, 1, x, 1, r1,
337 $ r2, w, r, info )
338 CALL chkxer( 'ZGBRFS', infot, nout, lerr, ok )
339 infot = 7
340 CALL zgbrfs( 'N', 2, 1, 1, 1, a, 2, af, 4, ip, b, 2, x, 2, r1,
341 $ r2, w, r, info )
342 CALL chkxer( 'ZGBRFS', infot, nout, lerr, ok )
343 infot = 9
344 CALL zgbrfs( 'N', 2, 1, 1, 1, a, 3, af, 3, ip, b, 2, x, 2, r1,
345 $ r2, w, r, info )
346 CALL chkxer( 'ZGBRFS', infot, nout, lerr, ok )
347 infot = 12
348 CALL zgbrfs( 'N', 2, 0, 0, 1, a, 1, af, 1, ip, b, 1, x, 2, r1,
349 $ r2, w, r, info )
350 CALL chkxer( 'ZGBRFS', infot, nout, lerr, ok )
351 infot = 14
352 CALL zgbrfs( 'N', 2, 0, 0, 1, a, 1, af, 1, ip, b, 2, x, 1, r1,
353 $ r2, w, r, info )
354 CALL chkxer( 'ZGBRFS', infot, nout, lerr, ok )
355*
356* ZGBCON
357*
358 srnamt = 'ZGBCON'
359 infot = 1
360 CALL zgbcon( '/', 0, 0, 0, a, 1, ip, anrm, rcond, w, r, info )
361 CALL chkxer( 'ZGBCON', infot, nout, lerr, ok )
362 infot = 2
363 CALL zgbcon( '1', -1, 0, 0, a, 1, ip, anrm, rcond, w, r, info )
364 CALL chkxer( 'ZGBCON', infot, nout, lerr, ok )
365 infot = 3
366 CALL zgbcon( '1', 1, -1, 0, a, 1, ip, anrm, rcond, w, r, info )
367 CALL chkxer( 'ZGBCON', infot, nout, lerr, ok )
368 infot = 4
369 CALL zgbcon( '1', 1, 0, -1, a, 1, ip, anrm, rcond, w, r, info )
370 CALL chkxer( 'ZGBCON', infot, nout, lerr, ok )
371 infot = 6
372 CALL zgbcon( '1', 2, 1, 1, a, 3, ip, anrm, rcond, w, r, info )
373 CALL chkxer( 'ZGBCON', infot, nout, lerr, ok )
374*
375* ZGBEQU
376*
377 srnamt = 'ZGBEQU'
378 infot = 1
379 CALL zgbequ( -1, 0, 0, 0, a, 1, r1, r2, rcond, ccond, anrm,
380 $ info )
381 CALL chkxer( 'ZGBEQU', infot, nout, lerr, ok )
382 infot = 2
383 CALL zgbequ( 0, -1, 0, 0, a, 1, r1, r2, rcond, ccond, anrm,
384 $ info )
385 CALL chkxer( 'ZGBEQU', infot, nout, lerr, ok )
386 infot = 3
387 CALL zgbequ( 1, 1, -1, 0, a, 1, r1, r2, rcond, ccond, anrm,
388 $ info )
389 CALL chkxer( 'ZGBEQU', infot, nout, lerr, ok )
390 infot = 4
391 CALL zgbequ( 1, 1, 0, -1, a, 1, r1, r2, rcond, ccond, anrm,
392 $ info )
393 CALL chkxer( 'ZGBEQU', infot, nout, lerr, ok )
394 infot = 6
395 CALL zgbequ( 2, 2, 1, 1, a, 2, r1, r2, rcond, ccond, anrm,
396 $ info )
397 CALL chkxer( 'ZGBEQU', infot, nout, lerr, ok )
398 END IF
399*
400* Print a summary line.
401*
402 CALL alaesm( path, ok, nout )
403*
404 RETURN
405*
406* End of ZERRGE
407*
subroutine alaesm(path, ok, nout)
ALAESM
Definition alaesm.f:63
subroutine chkxer(srnamt, infot, nout, lerr, ok)
Definition cblat2.f:3224
subroutine zgbcon(norm, n, kl, ku, ab, ldab, ipiv, anorm, rcond, work, rwork, info)
ZGBCON
Definition zgbcon.f:147
subroutine zgbequ(m, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd, amax, info)
ZGBEQU
Definition zgbequ.f:154
subroutine zgbrfs(trans, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZGBRFS
Definition zgbrfs.f:206
subroutine zgbtf2(m, n, kl, ku, ab, ldab, ipiv, info)
ZGBTF2 computes the LU factorization of a general band matrix using the unblocked version of the algo...
Definition zgbtf2.f:145
subroutine zgbtrf(m, n, kl, ku, ab, ldab, ipiv, info)
ZGBTRF
Definition zgbtrf.f:144
subroutine zgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
ZGBTRS
Definition zgbtrs.f:138
subroutine zgecon(norm, n, a, lda, anorm, rcond, work, rwork, info)
ZGECON
Definition zgecon.f:132
subroutine zgeequ(m, n, a, lda, r, c, rowcnd, colcnd, amax, info)
ZGEEQU
Definition zgeequ.f:140
subroutine zgerfs(trans, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZGERFS
Definition zgerfs.f:186
subroutine zgetf2(m, n, a, lda, ipiv, info)
ZGETF2 computes the LU factorization of a general m-by-n matrix using partial pivoting with row inter...
Definition zgetf2.f:108
subroutine zgetrf(m, n, a, lda, ipiv, info)
ZGETRF
Definition zgetrf.f:108
subroutine zgetri(n, a, lda, ipiv, work, lwork, info)
ZGETRI
Definition zgetri.f:114
subroutine zgetrs(trans, n, nrhs, a, lda, ipiv, b, ldb, info)
ZGETRS
Definition zgetrs.f:121
logical function lsamen(n, ca, cb)
LSAMEN
Definition lsamen.f:74
Here is the call graph for this function:
Here is the caller graph for this function: