SCALAPACK 2.2.2
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ samntest()

subroutine samntest ( integer  outnum,
integer  verb,
integer  topsrepeat,
integer  topscohrnt,
integer  nscope,
character*1, dimension(nscope)  scope0,
integer  ntop,
character*1, dimension(ntop)  top0,
integer  nmat,
integer, dimension(nmat)  m0,
integer, dimension(nmat)  n0,
integer, dimension(nmat)  ldas0,
integer, dimension(nmat)  ldad0,
integer, dimension(nmat)  ldi0,
integer  ndest,
integer, dimension(ndest)  rdest0,
integer, dimension(ndest)  cdest0,
integer  ngrid,
integer, dimension(ngrid)  context0,
integer, dimension(ngrid)  p0,
integer, dimension(ngrid)  q0,
integer, dimension(*)  iseed,
integer, dimension(rclen)  rmem,
integer, dimension(rclen)  cmem,
integer  rclen,
real, dimension(memlen)  mem,
integer  memlen 
)

Definition at line 19136 of file blacstest.f.

19141*
19142* -- BLACS tester (version 1.0) --
19143* University of Tennessee
19144* December 15, 1994
19145*
19146*
19147* .. Scalar Arguments ..
19148 INTEGER MEMLEN, NDEST, NGRID, NMAT, NSCOPE, NTOP, OUTNUM, RCLEN,
19149 $ TOPSCOHRNT, TOPSREPEAT, VERB
19150* ..
19151* .. Array Arguments ..
19152 CHARACTER*1 SCOPE0(NSCOPE), TOP0(NTOP)
19153 INTEGER M0(NMAT), N0(NMAT), LDAS0(NMAT), LDAD0(NMAT), LDI0(NMAT)
19154 INTEGER RDEST0(NDEST), CDEST0(NDEST), CONTEXT0(NGRID)
19155 INTEGER P0(NGRID), Q0(NGRID), ISEED(*), RMEM(RCLEN), CMEM(RCLEN)
19156 REAL MEM(MEMLEN)
19157* ..
19158*
19159* Purpose
19160* =======
19161* STESTAMN: Test real AMN COMBINE
19162*
19163* Arguments
19164* =========
19165* OUTNUM (input) INTEGER
19166* The device number to write output to.
19167*
19168* VERB (input) INTEGER
19169* The level of verbosity (how much printing to do).
19170*
19171* NSCOPE (input) INTEGER
19172* The number of scopes to be tested.
19173*
19174* SCOPE0 (input) CHARACTER*1 array of dimension (NSCOPE)
19175* Values of the scopes to be tested.
19176*
19177* NTOP (input) INTEGER
19178* The number of topologies to be tested.
19179*
19180* TOP0 (input) CHARACTER*1 array of dimension (NTOP)
19181* Values of the topologies to be tested.
19182*
19183* NMAT (input) INTEGER
19184* The number of matrices to be tested.
19185*
19186* M0 (input) INTEGER array of dimension (NMAT)
19187* Values of M to be tested.
19188*
19189* M0 (input) INTEGER array of dimension (NMAT)
19190* Values of M to be tested.
19191*
19192* N0 (input) INTEGER array of dimension (NMAT)
19193* Values of N to be tested.
19194*
19195* LDAS0 (input) INTEGER array of dimension (NMAT)
19196* Values of LDAS (leading dimension of A on source process)
19197* to be tested.
19198*
19199* LDAD0 (input) INTEGER array of dimension (NMAT)
19200* Values of LDAD (leading dimension of A on destination
19201* process) to be tested.
19202* LDI0 (input) INTEGER array of dimension (NMAT)
19203* Values of LDI (leading dimension of RA/CA) to be tested.
19204* If LDI == -1, these RA/CA should not be accessed.
19205*
19206* NDEST (input) INTEGER
19207* The number of destinations to be tested.
19208*
19209* RDEST0 (input) INTEGER array of dimension (NNDEST)
19210* Values of RDEST (row coordinate of destination) to be
19211* tested.
19212*
19213* CDEST0 (input) INTEGER array of dimension (NNDEST)
19214* Values of CDEST (column coordinate of destination) to be
19215* tested.
19216*
19217* NGRID (input) INTEGER
19218* The number of process grids to be tested.
19219*
19220* CONTEXT0 (input) INTEGER array of dimension (NGRID)
19221* The BLACS context handles corresponding to the grids.
19222*
19223* P0 (input) INTEGER array of dimension (NGRID)
19224* Values of P (number of process rows, NPROW).
19225*
19226* Q0 (input) INTEGER array of dimension (NGRID)
19227* Values of Q (number of process columns, NPCOL).
19228*
19229* ISEED (workspace) INTEGER array of dimension ( MAX(NPROCS, NTESTS) )
19230* Workspace used to hold each process's random number SEED.
19231* This requires NPROCS (number of processor) elements.
19232* If VERB < 2, this workspace also serves to indicate which
19233* tests fail. This requires workspace of NTESTS
19234* (number of tests performed).
19235*
19236* RMEM (workspace) INTEGER array of dimension (RCLEN)
19237* Used for all RA arrays, and their pre and post padding.
19238*
19239* CMEM (workspace) INTEGER array of dimension (RCLEN)
19240* Used for all CA arrays, and their pre and post padding.
19241*
19242* RCLEN (input) INTEGER
19243* The length, in elements, of RMEM and CMEM.
19244*
19245* MEM (workspace) REAL array of dimension (MEMLEN)
19246* Used for all other workspaces, including the matrix A,
19247* and its pre and post padding.
19248*
19249* MEMLEN (input) INTEGER
19250* The length, in elements, of MEM.
19251*
19252* =====================================================================
19253*
19254* .. External Functions ..
19255 LOGICAL ALLPASS, LSAME
19256 INTEGER IBTMYPROC, IBTNPROCS, IBTSIZEOF
19258* ..
19259* .. External Subroutines ..
19260 EXTERNAL blacs_gridinfo, sgamn2d
19261 EXTERNAL sinitmat, schkpad, sbtcheckin
19262* ..
19263* .. Local Scalars ..
19264 CHARACTER*1 SCOPE, TOP
19265 LOGICAL INGRID, TESTOK, ALLRCV
19266 INTEGER APTR, CAPTR, CDEST, CDEST2, CONTEXT, ERRDPTR, ERRIPTR, I,
19267 $ IAM, ICHECKVAL, IDE, IGR, IMA, IPAD, IPOST, IPRE, ISC,
19268 $ ISIZE, ISTART, ISTOP, ITC, ITC1, ITC2, ITO, ITR, ITR1,
19269 $ ITR2, J, K, LDA, LDADST, LDASRC, LDI, M, MAXERR, MYCOL,
19270 $ MYROW, N, NERR, NFAIL, NPCOL, NPROW, NSKIP, PREAPTR,
19271 $ RAPTR, RDEST, RDEST2, SETWHAT, SSIZE, TESTNUM, VALPTR
19272 REAL CHECKVAL
19273* ..
19274* .. Executable Statements ..
19275*
19276* Choose padding value, and make it unique
19277*
19278 checkval = -0.61e0
19279 iam = ibtmyproc()
19280 checkval = iam * checkval
19281 isize = ibtsizeof('I')
19282 ssize = ibtsizeof('S')
19283 icheckval = -iam
19284*
19285* Verify file parameters
19286*
19287 IF( iam .EQ. 0 ) THEN
19288 WRITE(outnum, *) ' '
19289 WRITE(outnum, *) ' '
19290 WRITE(outnum, 1000 )
19291 IF( verb .GT. 0 ) THEN
19292 WRITE(outnum,*) ' '
19293 WRITE(outnum, 2000) 'NSCOPE:', nscope
19294 WRITE(outnum, 3000) ' SCOPE:', ( scope0(i), i = 1, nscope )
19295 WRITE(outnum, 2000) 'TReps :', topsrepeat
19296 WRITE(outnum, 2000) 'TCohr :', topscohrnt
19297 WRITE(outnum, 2000) 'NTOP :', ntop
19298 WRITE(outnum, 3000) ' TOP :', ( top0(i), i = 1, ntop )
19299 WRITE(outnum, 2000) 'NMAT :', nmat
19300 WRITE(outnum, 2000) ' M :', ( m0(i), i = 1, nmat )
19301 WRITE(outnum, 2000) ' N :', ( n0(i), i = 1, nmat )
19302 WRITE(outnum, 2000) ' LDAS :', ( ldas0(i), i = 1, nmat )
19303 WRITE(outnum, 2000) ' LDAD :', ( ldad0(i), i = 1, nmat )
19304 WRITE(outnum, 2000) ' LDI :', ( ldi0(i), i = 1, nmat )
19305 WRITE(outnum, 2000) 'NDEST :', ndest
19306 WRITE(outnum, 2000) ' RDEST:',( rdest0(i), i = 1, ndest )
19307 WRITE(outnum, 2000) ' CDEST:',( cdest0(i), i = 1, ndest )
19308 WRITE(outnum, 2000) 'NGRIDS:', ngrid
19309 WRITE(outnum, 2000) ' P :', ( p0(i), i = 1, ngrid )
19310 WRITE(outnum, 2000) ' Q :', ( q0(i), i = 1, ngrid )
19311 WRITE(outnum, 2000) 'VERB :', verb
19312 WRITE(outnum,*) ' '
19313 END IF
19314 IF( verb .GT. 1 ) THEN
19315 WRITE(outnum,4000)
19316 WRITE(outnum,5000)
19317 END IF
19318 END IF
19319 IF (topsrepeat.EQ.0) THEN
19320 itr1 = 0
19321 itr2 = 0
19322 ELSE IF (topsrepeat.EQ.1) THEN
19323 itr1 = 1
19324 itr2 = 1
19325 ELSE
19326 itr1 = 0
19327 itr2 = 1
19328 END IF
19329*
19330* Find biggest matrix, so we know where to stick error info
19331*
19332 i = 0
19333 DO 10 ima = 1, nmat
19334 ipad = 4 * m0(ima)
19335 k = n0(ima) * max0( ldas0(ima), ldad0(ima) ) + ipad
19336 IF( k .GT. i ) i = k
19337 10 CONTINUE
19338 i = i + ibtnprocs()
19339 maxerr = ( ssize * (memlen-i) ) / ( ssize*2 + isize*6 )
19340 IF( maxerr .LT. 1 ) THEN
19341 WRITE(outnum,*) 'ERROR: Not enough memory to run MIN tests.'
19342 CALL blacs_abort(-1, 1)
19343 END IF
19344 errdptr = i + 1
19345 erriptr = errdptr + maxerr
19346 nerr = 0
19347 testnum = 0
19348 nfail = 0
19349 nskip = 0
19350*
19351* Loop over grids of matrix
19352*
19353 DO 90 igr = 1, ngrid
19354*
19355* allocate process grid for the next batch of tests
19356*
19357 context = context0(igr)
19358 CALL blacs_gridinfo( context, nprow, npcol, myrow, mycol )
19359 ingrid = ( (myrow.LT.nprow) .AND. (mycol.LT.npcol) )
19360*
19361 DO 80 isc = 1, nscope
19362 scope = scope0(isc)
19363 DO 70 ito = 1, ntop
19364 top = top0(ito)
19365*
19366* If testing multiring ('M') or general tree ('T'), need to
19367* loop over calls to BLACS_SET to do full test
19368*
19369 IF( lsame(top, 'M') ) THEN
19370 setwhat = 13
19371 IF( scope .EQ. 'R' ) THEN
19372 istart = -(npcol - 1)
19373 istop = -istart
19374 ELSE IF (scope .EQ. 'C') THEN
19375 istart = -(nprow - 1)
19376 istop = -istart
19377 ELSE
19378 istart = -(nprow*npcol - 1)
19379 istop = -istart
19380 ENDIF
19381 ELSE IF( lsame(top, 'T') ) THEN
19382 setwhat = 14
19383 istart = 1
19384 IF( scope .EQ. 'R' ) THEN
19385 istop = npcol - 1
19386 ELSE IF (scope .EQ. 'C') THEN
19387 istop = nprow - 1
19388 ELSE
19389 istop = nprow*npcol - 1
19390 ENDIF
19391 ELSE
19392 setwhat = 0
19393 istart = 1
19394 istop = 1
19395 ENDIF
19396 DO 60 ima = 1, nmat
19397 m = m0(ima)
19398 n = n0(ima)
19399 ldasrc = ldas0(ima)
19400 ldadst = ldad0(ima)
19401 ldi = ldi0(ima)
19402 ipre = 2 * m
19403 ipost = ipre
19404 preaptr = 1
19405 aptr = preaptr + ipre
19406*
19407 DO 50 ide = 1, ndest
19408 testnum = testnum + 1
19409 rdest2 = rdest0(ide)
19410 cdest2 = cdest0(ide)
19411*
19412* If everyone gets the answer, create some bogus rdest/cdest
19413* so IF's are easier
19414*
19415 allrcv = ( (rdest2.EQ.-1) .OR. (cdest2.EQ.-1) )
19416 IF( allrcv ) THEN
19417 rdest = nprow - 1
19418 cdest = npcol - 1
19419 IF (topscohrnt.EQ.0) THEN
19420 itr1 = 0
19421 itr2 = 0
19422 ELSE IF (topscohrnt.EQ.1) THEN
19423 itr1 = 1
19424 itr2 = 1
19425 ELSE
19426 itr1 = 0
19427 itr2 = 1
19428 END IF
19429 ELSE
19430 rdest = rdest2
19431 cdest = cdest2
19432 itc1 = 0
19433 itc2 = 0
19434 END IF
19435 IF( rdest.GE.p0(igr) .OR. cdest.GE.q0(igr) ) THEN
19436 nskip = nskip + 1
19437 GOTO 50
19438 END IF
19439*
19440 IF( myrow.EQ.rdest .AND. mycol.EQ.cdest ) THEN
19441 lda = ldadst
19442 ELSE
19443 lda = ldasrc
19444 END IF
19445 valptr = aptr + ipost + n * lda
19446 IF( verb .GT. 1 ) THEN
19447 IF( iam .EQ. 0 ) THEN
19448 WRITE(outnum, 6000)
19449 $ testnum, 'RUNNING', scope, top, m, n,
19450 $ ldasrc, ldadst, ldi, rdest2, cdest2,
19451 $ nprow, npcol
19452 END IF
19453 END IF
19454*
19455* If I am in scope
19456*
19457 testok = .true.
19458 IF( ingrid ) THEN
19459 IF( (myrow.EQ.rdest .AND. scope.EQ.'R') .OR.
19460 $ (mycol.EQ.cdest .AND. scope.EQ.'C') .OR.
19461 $ (scope .EQ. 'A') ) THEN
19462*
19463 k = nerr
19464 DO 40 itr = itr1, itr2
19465 CALL blacs_set(context, 15, itr)
19466 DO 35 itc = itc1, itc2
19467 CALL blacs_set(context, 16, itc)
19468 DO 30 j = istart, istop
19469 IF( j.EQ.0) GOTO 30
19470 IF( setwhat.NE.0 )
19471 $ CALL blacs_set(context, setwhat, j)
19472*
19473*
19474* generate and pad matrix A
19475*
19476 CALL sinitmat('G','-', m, n, mem(preaptr),
19477 $ lda, ipre, ipost,
19478 $ checkval, testnum,
19479 $ myrow, mycol )
19480*
19481* If they exist, pad RA and CA arrays
19482*
19483 IF( ldi .NE. -1 ) THEN
19484 DO 15 i = 1, n*ldi + ipre + ipost
19485 rmem(i) = icheckval
19486 cmem(i) = icheckval
19487 15 CONTINUE
19488 raptr = 1 + ipre
19489 captr = 1 + ipre
19490 ELSE
19491 DO 20 i = 1, ipre+ipost
19492 rmem(i) = icheckval
19493 cmem(i) = icheckval
19494 20 CONTINUE
19495 raptr = 1
19496 captr = 1
19497 END IF
19498*
19499 CALL sgamn2d(context, scope, top, m, n,
19500 $ mem(aptr), lda, rmem(raptr),
19501 $ cmem(captr), ldi,
19502 $ rdest2, cdest2)
19503*
19504* If I've got the answer, check for errors in
19505* matrix or padding
19506*
19507 IF( (myrow.EQ.rdest .AND. mycol.EQ.cdest)
19508 $ .OR. allrcv ) THEN
19509 CALL schkpad('G','-', m, n,
19510 $ mem(preaptr), lda, rdest,
19511 $ cdest, myrow, mycol,
19512 $ ipre, ipost, checkval,
19513 $ testnum, maxerr, nerr,
19514 $ mem(erriptr),mem(errdptr))
19515 CALL schkamn(scope, context, m, n,
19516 $ mem(aptr), lda,
19517 $ rmem(raptr), cmem(captr),
19518 $ ldi, testnum, maxerr,nerr,
19519 $ mem(erriptr),mem(errdptr),
19520 $ iseed, mem(valptr))
19521 CALL srcchk(ipre, ipost, icheckval,
19522 $ m, n, rmem, cmem, ldi,
19523 $ myrow, mycol, testnum,
19524 $ maxerr, nerr,
19525 $ mem(erriptr), mem(errdptr))
19526 END IF
19527 30 CONTINUE
19528 CALL blacs_set(context, 16, 0)
19529 35 CONTINUE
19530 CALL blacs_set(context, 15, 0)
19531 40 CONTINUE
19532 testok = ( k .EQ. nerr )
19533 END IF
19534 END IF
19535*
19536 IF( verb .GT. 1 ) THEN
19537 i = nerr
19538 CALL sbtcheckin(0, outnum, maxerr, nerr,
19539 $ mem(erriptr), mem(errdptr), iseed)
19540 IF( iam .EQ. 0 ) THEN
19541 IF( testok .AND. nerr.EQ.i ) THEN
19542 WRITE(outnum,6000)testnum,'PASSED ',
19543 $ scope, top, m, n, ldasrc,
19544 $ ldadst, ldi, rdest2, cdest2,
19545 $ nprow, npcol
19546 ELSE
19547 nfail = nfail + 1
19548 WRITE(outnum,6000)testnum,'FAILED ',
19549 $ scope, top, m, n, ldasrc,
19550 $ ldadst, ldi, rdest2, cdest2,
19551 $ nprow, npcol
19552 END IF
19553 END IF
19554*
19555* Once we've printed out errors, can re-use buf space
19556*
19557 nerr = 0
19558 END IF
19559 50 CONTINUE
19560 60 CONTINUE
19561 70 CONTINUE
19562 80 CONTINUE
19563 90 CONTINUE
19564*
19565 IF( verb .LT. 2 ) THEN
19566 nfail = testnum
19567 CALL sbtcheckin( nfail, outnum, maxerr, nerr, mem(erriptr),
19568 $ mem(errdptr), iseed )
19569 END IF
19570 IF( iam .EQ. 0 ) THEN
19571 IF( verb .GT. 1 ) WRITE(outnum,*) ' '
19572 IF( nfail+nskip .EQ. 0 ) THEN
19573 WRITE(outnum, 7000 ) testnum
19574 ELSE
19575 WRITE(outnum, 8000 ) testnum, testnum-nskip-nfail,
19576 $ nskip, nfail
19577 END IF
19578 END IF
19579*
19580* Log whether their were any failures
19581*
19582 testok = allpass( (nfail.EQ.0) )
19583*
19584 1000 FORMAT('REAL AMN TESTS: BEGIN.' )
19585 2000 FORMAT(1x,a7,3x,10i6)
19586 3000 FORMAT(1x,a7,3x,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,
19587 $ 5x,a1,5x,a1)
19588 4000 FORMAT(' TEST# STATUS SCOPE TOP M N LDAS LDAD LDI ',
19589 $ 'RDEST CDEST P Q')
19590 5000 FORMAT(' ----- ------- ----- --- ----- ----- ----- ----- ----- ',
19591 $ '----- ----- ---- ----')
19592 6000 FORMAT(i6,1x,a7,5x,a1,3x,a1,7i6,2i5)
19593 7000 FORMAT('REAL AMN TESTS: PASSED ALL',
19594 $ i5, ' TESTS.')
19595 8000 FORMAT('REAL AMN TESTS:',i5,' TESTS;',i5,' PASSED,',
19596 $ i5,' SKIPPED,',i5,' FAILED.')
19597*
19598 RETURN
19599*
19600* End of STESTAMN.
19601*
subroutine sbtcheckin(nftests, outnum, maxerr, nerr, ierr, sval, tfailed)
Definition blacstest.f:7341
subroutine schkpad(uplo, diag, m, n, mem, lda, rsrc, csrc, myrow, mycol, ipre, ipost, checkval, testnum, maxerr, nerr, erribuf, errdbuf)
Definition blacstest.f:7746
subroutine schkamn(scope, ictxt, m, n, a, lda, ra, ca, ldi, testnum, maxerr, nerr, erribuf, errdbuf, iseed, vals)
logical function allpass(thistest)
Definition blacstest.f:1881
subroutine sinitmat(uplo, diag, m, n, mem, lda, ipre, ipost, checkval, testnum, myrow, mycol)
Definition blacstest.f:7463
subroutine srcchk(ipre, ipost, padval, m, n, ra, ca, ldi, myrow, mycol, testnum, maxerr, nerr, erribuf, errdbuf)
integer function ibtnprocs()
Definition btprim.f:81
integer function ibtmyproc()
Definition btprim.f:47
integer function ibtsizeof(type)
Definition btprim.f:286
logical function lsame(ca, cb)
Definition tools.f:1724
Here is the call graph for this function:
Here is the caller graph for this function: