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

◆ dchkgb()

subroutine dchkgb ( logical, dimension( * )  dotype,
integer  nm,
integer, dimension( * )  mval,
integer  nn,
integer, dimension( * )  nval,
integer  nnb,
integer, dimension( * )  nbval,
integer  nns,
integer, dimension( * )  nsval,
double precision  thresh,
logical  tsterr,
double precision, dimension( * )  a,
integer  la,
double precision, dimension( * )  afac,
integer  lafac,
double precision, dimension( * )  b,
double precision, dimension( * )  x,
double precision, dimension( * )  xact,
double precision, dimension( * )  work,
double precision, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  nout 
)

DCHKGB

Purpose:
 DCHKGB tests DGBTRF, -TRS, -RFS, and -CON
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]NM
          NM is INTEGER
          The number of values of M contained in the vector MVAL.
[in]MVAL
          MVAL is INTEGER array, dimension (NM)
          The values of the matrix row dimension M.
[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 column 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 blocksize NB.
[in]NNS
          NNS is INTEGER
          The number of values of NRHS contained in the vector NSVAL.
[in]NSVAL
          NSVAL is INTEGER array, dimension (NNS)
          The values of the number of right hand sides NRHS.
[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 DOUBLE PRECISION array, dimension (LA)
[in]LA
          LA is INTEGER
          The length of the array A.  LA >= (KLMAX+KUMAX+1)*NMAX
          where KLMAX is the largest entry in the local array KLVAL,
                KUMAX is the largest entry in the local array KUVAL and
                NMAX is the largest entry in the input array NVAL.
[out]AFAC
          AFAC is DOUBLE PRECISION array, dimension (LAFAC)
[in]LAFAC
          LAFAC is INTEGER
          The length of the array AFAC. LAFAC >= (2*KLMAX+KUMAX+1)*NMAX
          where KLMAX is the largest entry in the local array KLVAL,
                KUMAX is the largest entry in the local array KUVAL and
                NMAX is the largest entry in the input array NVAL.
[out]B
          B is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension
                      (NMAX*max(3,NSMAX,NMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      (NMAX+2*NSMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (2*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 188 of file dchkgb.f.

191*
192* -- LAPACK test routine --
193* -- LAPACK is a software package provided by Univ. of Tennessee, --
194* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
195*
196* .. Scalar Arguments ..
197 LOGICAL TSTERR
198 INTEGER LA, LAFAC, NM, NN, NNB, NNS, NOUT
199 DOUBLE PRECISION THRESH
200* ..
201* .. Array Arguments ..
202 LOGICAL DOTYPE( * )
203 INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
204 $ NVAL( * )
205 DOUBLE PRECISION A( * ), AFAC( * ), B( * ), RWORK( * ),
206 $ WORK( * ), X( * ), XACT( * )
207* ..
208*
209* =====================================================================
210*
211* .. Parameters ..
212 DOUBLE PRECISION ONE, ZERO
213 parameter( one = 1.0d+0, zero = 0.0d+0 )
214 INTEGER NTYPES, NTESTS
215 parameter( ntypes = 8, ntests = 7 )
216 INTEGER NBW, NTRAN
217 parameter( nbw = 4, ntran = 3 )
218* ..
219* .. Local Scalars ..
220 LOGICAL TRFCON, ZEROT
221 CHARACTER DIST, NORM, TRANS, TYPE, XTYPE
222 CHARACTER*3 PATH
223 INTEGER I, I1, I2, IKL, IKU, IM, IMAT, IN, INB, INFO,
224 $ IOFF, IRHS, ITRAN, IZERO, J, K, KL, KOFF, KU,
225 $ LDA, LDAFAC, LDB, M, MODE, N, NB, NERRS, NFAIL,
226 $ NIMAT, NKL, NKU, NRHS, NRUN
227 DOUBLE PRECISION AINVNM, ANORM, ANORMI, ANORMO, CNDNUM, RCOND,
228 $ RCONDC, RCONDI, RCONDO
229* ..
230* .. Local Arrays ..
231 CHARACTER TRANSS( NTRAN )
232 INTEGER ISEED( 4 ), ISEEDY( 4 ), KLVAL( NBW ),
233 $ KUVAL( NBW )
234 DOUBLE PRECISION RESULT( NTESTS )
235* ..
236* .. External Functions ..
237 DOUBLE PRECISION DGET06, DLANGB, DLANGE
238 EXTERNAL dget06, dlangb, dlange
239* ..
240* .. External Subroutines ..
241 EXTERNAL alaerh, alahd, alasum, dcopy, derrge, dgbcon,
244 $ xlaenv
245* ..
246* .. Intrinsic Functions ..
247 INTRINSIC max, min
248* ..
249* .. Scalars in Common ..
250 LOGICAL LERR, OK
251 CHARACTER*32 SRNAMT
252 INTEGER INFOT, NUNIT
253* ..
254* .. Common blocks ..
255 COMMON / infoc / infot, nunit, ok, lerr
256 COMMON / srnamc / srnamt
257* ..
258* .. Data statements ..
259 DATA iseedy / 1988, 1989, 1990, 1991 / ,
260 $ transs / 'N', 'T', 'C' /
261* ..
262* .. Executable Statements ..
263*
264* Initialize constants and the random number seed.
265*
266 path( 1: 1 ) = 'Double precision'
267 path( 2: 3 ) = 'GB'
268 nrun = 0
269 nfail = 0
270 nerrs = 0
271 DO 10 i = 1, 4
272 iseed( i ) = iseedy( i )
273 10 CONTINUE
274*
275* Test the error exits
276*
277 IF( tsterr )
278 $ CALL derrge( path, nout )
279 infot = 0
280 CALL xlaenv( 2, 2 )
281*
282* Initialize the first value for the lower and upper bandwidths.
283*
284 klval( 1 ) = 0
285 kuval( 1 ) = 0
286*
287* Do for each value of M in MVAL
288*
289 DO 160 im = 1, nm
290 m = mval( im )
291*
292* Set values to use for the lower bandwidth.
293*
294 klval( 2 ) = m + ( m+1 ) / 4
295*
296* KLVAL( 2 ) = MAX( M-1, 0 )
297*
298 klval( 3 ) = ( 3*m-1 ) / 4
299 klval( 4 ) = ( m+1 ) / 4
300*
301* Do for each value of N in NVAL
302*
303 DO 150 in = 1, nn
304 n = nval( in )
305 xtype = 'N'
306*
307* Set values to use for the upper bandwidth.
308*
309 kuval( 2 ) = n + ( n+1 ) / 4
310*
311* KUVAL( 2 ) = MAX( N-1, 0 )
312*
313 kuval( 3 ) = ( 3*n-1 ) / 4
314 kuval( 4 ) = ( n+1 ) / 4
315*
316* Set limits on the number of loop iterations.
317*
318 nkl = min( m+1, 4 )
319 IF( n.EQ.0 )
320 $ nkl = 2
321 nku = min( n+1, 4 )
322 IF( m.EQ.0 )
323 $ nku = 2
324 nimat = ntypes
325 IF( m.LE.0 .OR. n.LE.0 )
326 $ nimat = 1
327*
328 DO 140 ikl = 1, nkl
329*
330* Do for KL = 0, (5*M+1)/4, (3M-1)/4, and (M+1)/4. This
331* order makes it easier to skip redundant values for small
332* values of M.
333*
334 kl = klval( ikl )
335 DO 130 iku = 1, nku
336*
337* Do for KU = 0, (5*N+1)/4, (3N-1)/4, and (N+1)/4. This
338* order makes it easier to skip redundant values for
339* small values of N.
340*
341 ku = kuval( iku )
342*
343* Check that A and AFAC are big enough to generate this
344* matrix.
345*
346 lda = kl + ku + 1
347 ldafac = 2*kl + ku + 1
348 IF( ( lda*n ).GT.la .OR. ( ldafac*n ).GT.lafac ) THEN
349 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
350 $ CALL alahd( nout, path )
351 IF( n*( kl+ku+1 ).GT.la ) THEN
352 WRITE( nout, fmt = 9999 )la, m, n, kl, ku,
353 $ n*( kl+ku+1 )
354 nerrs = nerrs + 1
355 END IF
356 IF( n*( 2*kl+ku+1 ).GT.lafac ) THEN
357 WRITE( nout, fmt = 9998 )lafac, m, n, kl, ku,
358 $ n*( 2*kl+ku+1 )
359 nerrs = nerrs + 1
360 END IF
361 GO TO 130
362 END IF
363*
364 DO 120 imat = 1, nimat
365*
366* Do the tests only if DOTYPE( IMAT ) is true.
367*
368 IF( .NOT.dotype( imat ) )
369 $ GO TO 120
370*
371* Skip types 2, 3, or 4 if the matrix size is too
372* small.
373*
374 zerot = imat.GE.2 .AND. imat.LE.4
375 IF( zerot .AND. n.LT.imat-1 )
376 $ GO TO 120
377*
378 IF( .NOT.zerot .OR. .NOT.dotype( 1 ) ) THEN
379*
380* Set up parameters with DLATB4 and generate a
381* test matrix with DLATMS.
382*
383 CALL dlatb4( path, imat, m, n, TYPE, KL, KU,
384 $ ANORM, MODE, CNDNUM, DIST )
385*
386 koff = max( 1, ku+2-n )
387 DO 20 i = 1, koff - 1
388 a( i ) = zero
389 20 CONTINUE
390 srnamt = 'DLATMS'
391 CALL dlatms( m, n, dist, iseed, TYPE, RWORK,
392 $ MODE, CNDNUM, ANORM, KL, KU, 'Z',
393 $ A( KOFF ), LDA, WORK, INFO )
394*
395* Check the error code from DLATMS.
396*
397 IF( info.NE.0 ) THEN
398 CALL alaerh( path, 'DLATMS', info, 0, ' ', m,
399 $ n, kl, ku, -1, imat, nfail,
400 $ nerrs, nout )
401 GO TO 120
402 END IF
403 ELSE IF( izero.GT.0 ) THEN
404*
405* Use the same matrix for types 3 and 4 as for
406* type 2 by copying back the zeroed out column.
407*
408 CALL dcopy( i2-i1+1, b, 1, a( ioff+i1 ), 1 )
409 END IF
410*
411* For types 2, 3, and 4, zero one or more columns of
412* the matrix to test that INFO is returned correctly.
413*
414 izero = 0
415 IF( zerot ) THEN
416 IF( imat.EQ.2 ) THEN
417 izero = 1
418 ELSE IF( imat.EQ.3 ) THEN
419 izero = min( m, n )
420 ELSE
421 izero = min( m, n ) / 2 + 1
422 END IF
423 ioff = ( izero-1 )*lda
424 IF( imat.LT.4 ) THEN
425*
426* Store the column to be zeroed out in B.
427*
428 i1 = max( 1, ku+2-izero )
429 i2 = min( kl+ku+1, ku+1+( m-izero ) )
430 CALL dcopy( i2-i1+1, a( ioff+i1 ), 1, b, 1 )
431*
432 DO 30 i = i1, i2
433 a( ioff+i ) = zero
434 30 CONTINUE
435 ELSE
436 DO 50 j = izero, n
437 DO 40 i = max( 1, ku+2-j ),
438 $ min( kl+ku+1, ku+1+( m-j ) )
439 a( ioff+i ) = zero
440 40 CONTINUE
441 ioff = ioff + lda
442 50 CONTINUE
443 END IF
444 END IF
445*
446* These lines, if used in place of the calls in the
447* loop over INB, cause the code to bomb on a Sun
448* SPARCstation.
449*
450* ANORMO = DLANGB( 'O', N, KL, KU, A, LDA, RWORK )
451* ANORMI = DLANGB( 'I', N, KL, KU, A, LDA, RWORK )
452*
453* Do for each blocksize in NBVAL
454*
455 DO 110 inb = 1, nnb
456 nb = nbval( inb )
457 CALL xlaenv( 1, nb )
458*
459* Compute the LU factorization of the band matrix.
460*
461 IF( m.GT.0 .AND. n.GT.0 )
462 $ CALL dlacpy( 'Full', kl+ku+1, n, a, lda,
463 $ afac( kl+1 ), ldafac )
464 srnamt = 'DGBTRF'
465 CALL dgbtrf( m, n, kl, ku, afac, ldafac, iwork,
466 $ info )
467*
468* Check error code from DGBTRF.
469*
470 IF( info.NE.izero )
471 $ CALL alaerh( path, 'DGBTRF', info, izero,
472 $ ' ', m, n, kl, ku, nb, imat,
473 $ nfail, nerrs, nout )
474 trfcon = .false.
475*
476*+ TEST 1
477* Reconstruct matrix from factors and compute
478* residual.
479*
480 CALL dgbt01( m, n, kl, ku, a, lda, afac, ldafac,
481 $ iwork, work, result( 1 ) )
482*
483* Print information about the tests so far that
484* did not pass the threshold.
485*
486 IF( result( 1 ).GE.thresh ) THEN
487 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
488 $ CALL alahd( nout, path )
489 WRITE( nout, fmt = 9997 )m, n, kl, ku, nb,
490 $ imat, 1, result( 1 )
491 nfail = nfail + 1
492 END IF
493 nrun = nrun + 1
494*
495* Skip the remaining tests if this is not the
496* first block size or if M .ne. N.
497*
498 IF( inb.GT.1 .OR. m.NE.n )
499 $ GO TO 110
500*
501 anormo = dlangb( 'O', n, kl, ku, a, lda, rwork )
502 anormi = dlangb( 'I', n, kl, ku, a, lda, rwork )
503*
504 IF( info.EQ.0 ) THEN
505*
506* Form the inverse of A so we can get a good
507* estimate of CNDNUM = norm(A) * norm(inv(A)).
508*
509 ldb = max( 1, n )
510 CALL dlaset( 'Full', n, n, zero, one, work,
511 $ ldb )
512 srnamt = 'DGBTRS'
513 CALL dgbtrs( 'No transpose', n, kl, ku, n,
514 $ afac, ldafac, iwork, work, ldb,
515 $ info )
516*
517* Compute the 1-norm condition number of A.
518*
519 ainvnm = dlange( 'O', n, n, work, ldb,
520 $ rwork )
521 IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
522 rcondo = one
523 ELSE
524 rcondo = ( one / anormo ) / ainvnm
525 END IF
526*
527* Compute the infinity-norm condition number of
528* A.
529*
530 ainvnm = dlange( 'I', n, n, work, ldb,
531 $ rwork )
532 IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
533 rcondi = one
534 ELSE
535 rcondi = ( one / anormi ) / ainvnm
536 END IF
537 ELSE
538*
539* Do only the condition estimate if INFO.NE.0.
540*
541 trfcon = .true.
542 rcondo = zero
543 rcondi = zero
544 END IF
545*
546* Skip the solve tests if the matrix is singular.
547*
548 IF( trfcon )
549 $ GO TO 90
550*
551 DO 80 irhs = 1, nns
552 nrhs = nsval( irhs )
553 xtype = 'N'
554*
555 DO 70 itran = 1, ntran
556 trans = transs( itran )
557 IF( itran.EQ.1 ) THEN
558 rcondc = rcondo
559 norm = 'O'
560 ELSE
561 rcondc = rcondi
562 norm = 'I'
563 END IF
564*
565*+ TEST 2:
566* Solve and compute residual for op(A) * X = B.
567*
568 srnamt = 'DLARHS'
569 CALL dlarhs( path, xtype, ' ', trans, n,
570 $ n, kl, ku, nrhs, a, lda,
571 $ xact, ldb, b, ldb, iseed,
572 $ info )
573 xtype = 'C'
574 CALL dlacpy( 'Full', n, nrhs, b, ldb, x,
575 $ ldb )
576*
577 srnamt = 'DGBTRS'
578 CALL dgbtrs( trans, n, kl, ku, nrhs, afac,
579 $ ldafac, iwork, x, ldb, info )
580*
581* Check error code from DGBTRS.
582*
583 IF( info.NE.0 )
584 $ CALL alaerh( path, 'DGBTRS', info, 0,
585 $ trans, n, n, kl, ku, -1,
586 $ imat, nfail, nerrs, nout )
587*
588 CALL dlacpy( 'Full', n, nrhs, b, ldb,
589 $ work, ldb )
590 CALL dgbt02( trans, m, n, kl, ku, nrhs, a,
591 $ lda, x, ldb, work, ldb,
592 $ rwork, result( 2 ) )
593*
594*+ TEST 3:
595* Check solution from generated exact
596* solution.
597*
598 CALL dget04( n, nrhs, x, ldb, xact, ldb,
599 $ rcondc, result( 3 ) )
600*
601*+ TESTS 4, 5, 6:
602* Use iterative refinement to improve the
603* solution.
604*
605 srnamt = 'DGBRFS'
606 CALL dgbrfs( trans, n, kl, ku, nrhs, a,
607 $ lda, afac, ldafac, iwork, b,
608 $ ldb, x, ldb, rwork,
609 $ rwork( nrhs+1 ), work,
610 $ iwork( n+1 ), info )
611*
612* Check error code from DGBRFS.
613*
614 IF( info.NE.0 )
615 $ CALL alaerh( path, 'DGBRFS', info, 0,
616 $ trans, n, n, kl, ku, nrhs,
617 $ imat, nfail, nerrs, nout )
618*
619 CALL dget04( n, nrhs, x, ldb, xact, ldb,
620 $ rcondc, result( 4 ) )
621 CALL dgbt05( trans, n, kl, ku, nrhs, a,
622 $ lda, b, ldb, x, ldb, xact,
623 $ ldb, rwork, rwork( nrhs+1 ),
624 $ result( 5 ) )
625 DO 60 k = 2, 6
626 IF( result( k ).GE.thresh ) THEN
627 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
628 $ CALL alahd( nout, path )
629 WRITE( nout, fmt = 9996 )trans, n,
630 $ kl, ku, nrhs, imat, k,
631 $ result( k )
632 nfail = nfail + 1
633 END IF
634 60 CONTINUE
635 nrun = nrun + 5
636 70 CONTINUE
637 80 CONTINUE
638*
639*+ TEST 7:
640* Get an estimate of RCOND = 1/CNDNUM.
641*
642 90 CONTINUE
643 DO 100 itran = 1, 2
644 IF( itran.EQ.1 ) THEN
645 anorm = anormo
646 rcondc = rcondo
647 norm = 'O'
648 ELSE
649 anorm = anormi
650 rcondc = rcondi
651 norm = 'I'
652 END IF
653 srnamt = 'DGBCON'
654 CALL dgbcon( norm, n, kl, ku, afac, ldafac,
655 $ iwork, anorm, rcond, work,
656 $ iwork( n+1 ), info )
657*
658* Check error code from DGBCON.
659*
660 IF( info.NE.0 )
661 $ CALL alaerh( path, 'DGBCON', info, 0,
662 $ norm, n, n, kl, ku, -1, imat,
663 $ nfail, nerrs, nout )
664*
665 result( 7 ) = dget06( rcond, rcondc )
666*
667* Print information about the tests that did
668* not pass the threshold.
669*
670 IF( result( 7 ).GE.thresh ) THEN
671 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
672 $ CALL alahd( nout, path )
673 WRITE( nout, fmt = 9995 )norm, n, kl, ku,
674 $ imat, 7, result( 7 )
675 nfail = nfail + 1
676 END IF
677 nrun = nrun + 1
678 100 CONTINUE
679*
680 110 CONTINUE
681 120 CONTINUE
682 130 CONTINUE
683 140 CONTINUE
684 150 CONTINUE
685 160 CONTINUE
686*
687* Print a summary of the results.
688*
689 CALL alasum( path, nout, nfail, nrun, nerrs )
690*
691 9999 FORMAT( ' *** In DCHKGB, LA=', i5, ' is too small for M=', i5,
692 $ ', N=', i5, ', KL=', i4, ', KU=', i4,
693 $ / ' ==> Increase LA to at least ', i5 )
694 9998 FORMAT( ' *** In DCHKGB, LAFAC=', i5, ' is too small for M=', i5,
695 $ ', N=', i5, ', KL=', i4, ', KU=', i4,
696 $ / ' ==> Increase LAFAC to at least ', i5 )
697 9997 FORMAT( ' M =', i5, ', N =', i5, ', KL=', i5, ', KU=', i5,
698 $ ', NB =', i4, ', type ', i1, ', test(', i1, ')=', g12.5 )
699 9996 FORMAT( ' TRANS=''', a1, ''', N=', i5, ', KL=', i5, ', KU=', i5,
700 $ ', NRHS=', i3, ', type ', i1, ', test(', i1, ')=', g12.5 )
701 9995 FORMAT( ' NORM =''', a1, ''', N=', i5, ', KL=', i5, ', KU=', i5,
702 $ ',', 10x, ' type ', i1, ', test(', i1, ')=', g12.5 )
703*
704 RETURN
705*
706* End of DCHKGB
707*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
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 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 derrge(path, nunit)
DERRGE
Definition derrge.f:55
subroutine dgbt01(m, n, kl, ku, a, lda, afac, ldafac, ipiv, work, resid)
DGBT01
Definition dgbt01.f:126
subroutine dgbt02(trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
DGBT02
Definition dgbt02.f:149
subroutine dgbt05(trans, n, kl, ku, nrhs, ab, ldab, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
DGBT05
Definition dgbt05.f:176
subroutine dget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
DGET04
Definition dget04.f:102
double precision function dget06(rcond, rcondc)
DGET06
Definition dget06.f:55
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 dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgbcon(norm, n, kl, ku, ab, ldab, ipiv, anorm, rcond, work, iwork, info)
DGBCON
Definition dgbcon.f:146
subroutine dgbrfs(trans, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv, b, ldb, x, ldx, ferr, berr, work, iwork, info)
DGBRFS
Definition dgbrfs.f:205
subroutine dgbtrf(m, n, kl, ku, ab, ldab, ipiv, info)
DGBTRF
Definition dgbtrf.f:144
subroutine dgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
DGBTRS
Definition dgbtrs.f:138
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
double precision function dlangb(norm, n, kl, ku, ab, ldab, work)
DLANGB returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlangb.f:124
double precision function dlange(norm, m, n, a, lda, work)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlange.f:114
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
Here is the call graph for this function:
Here is the caller graph for this function: