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

◆ cbsbrtest()

subroutine cbsbrtest ( integer  outnum,
integer  verb,
integer  nscope,
character*1, dimension(nscope)  scope0,
integer  ntop,
character*1, dimension(ntop)  top0,
integer  nshape,
character*1, dimension(nshape)  uplo0,
character*1, dimension(nshape)  diag0,
integer  nmat,
integer, dimension(nmat)  m0,
integer, dimension(nmat)  n0,
integer, dimension(nmat)  ldas0,
integer, dimension(nmat)  ldad0,
integer  nsrc,
integer, dimension(nsrc)  rsrc0,
integer, dimension(nsrc)  csrc0,
integer  ngrid,
integer, dimension(ngrid)  context0,
integer, dimension(ngrid)  p0,
integer, dimension(ngrid)  q0,
integer, dimension(*)  tfail,
complex, dimension(memlen)  mem,
integer  memlen 
)

Definition at line 5166 of file blacstest.f.

5170*
5171* -- BLACS tester (version 1.0) --
5172* University of Tennessee
5173* December 15, 1994
5174*
5175*
5176* .. Scalar Arguments ..
5177 INTEGER OUTNUM, VERB, NSCOPE, NTOP, NSHAPE, NMAT, NSRC, NGRID
5178 INTEGER MEMLEN
5179* ..
5180* .. Array Arguments ..
5181 CHARACTER*1 SCOPE0(NSCOPE), TOP0(NTOP)
5182 CHARACTER*1 UPLO0(NSHAPE), DIAG0(NSHAPE)
5183 INTEGER M0(NMAT), N0(NMAT), LDAS0(NMAT), LDAD0(NMAT)
5184 INTEGER RSRC0(NSRC), CSRC0(NSRC), CONTEXT0(NGRID)
5185 INTEGER P0(NGRID), Q0(NGRID), TFAIL(*)
5186 COMPLEX MEM(MEMLEN)
5187* ..
5188*
5189* Purpose
5190* =======
5191* CTESTBSBR: Test complex broadcast
5192*
5193* Arguments
5194* =========
5195* OUTNUM (input) INTEGER
5196* The device number to write output to.
5197*
5198* VERB (input) INTEGER
5199* The level of verbosity (how much printing to do).
5200*
5201* NSCOPE (input) INTEGER
5202* The number of scopes to be tested.
5203*
5204* SCOPE0 (input) CHARACTER*1 array of dimension (NSCOPE)
5205* Values of the scopes to be tested.
5206*
5207* NTOP (input) INTEGER
5208* The number of topologies to be tested.
5209*
5210* TOP0 (input) CHARACTER*1 array of dimension (NTOP)
5211* Values of the topologies to be tested.
5212*
5213* NSHAPE (input) INTEGER
5214* The number of matrix shapes to be tested.
5215*
5216* UPLO0 (input) CHARACTER*1 array of dimension (NSHAPE)
5217* Values of UPLO to be tested.
5218*
5219* DIAG0 (input) CHARACTER*1 array of dimension (NSHAPE)
5220* Values of DIAG to be tested.
5221*
5222* NMAT (input) INTEGER
5223* The number of matrices to be tested.
5224*
5225* M0 (input) INTEGER array of dimension (NMAT)
5226* Values of M to be tested.
5227*
5228* M0 (input) INTEGER array of dimension (NMAT)
5229* Values of M to be tested.
5230*
5231* N0 (input) INTEGER array of dimension (NMAT)
5232* Values of N to be tested.
5233*
5234* LDAS0 (input) INTEGER array of dimension (NMAT)
5235* Values of LDAS (leading dimension of A on source process)
5236* to be tested.
5237*
5238* LDAD0 (input) INTEGER array of dimension (NMAT)
5239* Values of LDAD (leading dimension of A on destination
5240* process) to be tested.
5241* NSRC (input) INTEGER
5242* The number of sources to be tested.
5243*
5244* RSRC0 (input) INTEGER array of dimension (NDEST)
5245* Values of RSRC (row coordinate of source) to be tested.
5246*
5247* CSRC0 (input) INTEGER array of dimension (NDEST)
5248* Values of CSRC (column coordinate of source) to be tested.
5249*
5250* NGRID (input) INTEGER
5251* The number of process grids to be tested.
5252*
5253* CONTEXT0 (input) INTEGER array of dimension (NGRID)
5254* The BLACS context handles corresponding to the grids.
5255*
5256* P0 (input) INTEGER array of dimension (NGRID)
5257* Values of P (number of process rows, NPROW).
5258*
5259* Q0 (input) INTEGER array of dimension (NGRID)
5260* Values of Q (number of process columns, NPCOL).
5261*
5262* TFAIL (workspace) INTEGER array of dimension (NTESTS)
5263* If VERB < 2, serves to indicate which tests fail. This
5264* requires workspace of NTESTS (number of tests performed).
5265*
5266* MEM (workspace) COMPLEX array of dimension (MEMLEN)
5267* Used for all other workspaces, including the matrix A,
5268* and its pre and post padding.
5269*
5270* MEMLEN (input) INTEGER
5271* The length, in elements, of MEM.
5272*
5273* =====================================================================
5274*
5275* .. External Functions ..
5276 LOGICAL ALLPASS, LSAME
5277 INTEGER IBTMYPROC, IBTSIZEOF
5278 EXTERNAL allpass, lsame, ibtmyproc, ibtsizeof
5279* ..
5280* .. External Subroutines ..
5281 EXTERNAL blacs_gridinfo
5282 EXTERNAL ctrbs2d, cgebs2d, ctrbr2d, cgebr2d
5284* ..
5285* .. Local Scalars ..
5286 CHARACTER*1 SCOPE, TOP, UPLO, DIAG
5287 LOGICAL TESTOK, INGRID
5288 INTEGER IAM, I, K, J, IGR, ISH, IMA, ISO, ISC, ITO
5289 INTEGER M, N, NPROW, NPCOL, MYROW, MYCOL, RSRC, CSRC
5290 INTEGER ISTART, ISTOP, IPRE, IPOST, SETWHAT
5291 INTEGER NERR, NSKIP, NFAIL, TESTNUM, CONTEXT, MAXERR, LDASRC
5292 INTEGER LDADST, ERRDPTR, APTR, ERRIPTR, ISIZE, CSIZE
5293 COMPLEX SCHECKVAL, RCHECKVAL
5294* ..
5295* .. Executable Statements ..
5296*
5297 scheckval = cmplx( -0.01, -0.01 )
5298 rcheckval = cmplx( -0.02, -0.02 )
5299*
5300 iam = ibtmyproc()
5301 isize = ibtsizeof('I')
5302 csize = ibtsizeof('C')
5303*
5304* Verify file parameters
5305*
5306 IF( iam .EQ. 0 ) THEN
5307 WRITE(outnum, *) ' '
5308 WRITE(outnum, *) ' '
5309 WRITE(outnum, 1000 )
5310 IF( verb .GT. 0 ) THEN
5311 WRITE(outnum,*) ' '
5312 WRITE(outnum, 2000) 'NSCOPE:', nscope
5313 WRITE(outnum, 3000) ' SCOPE:', ( scope0(i), i = 1, nscope )
5314 WRITE(outnum, 2000) 'NTOP :', ntop
5315 WRITE(outnum, 3000) ' TOP :', ( top0(i), i = 1, ntop )
5316 WRITE(outnum, 2000) 'NSHAPE:', nshape
5317 WRITE(outnum, 3000) ' UPLO :', ( uplo0(i), i = 1, nshape )
5318 WRITE(outnum, 3000) ' DIAG :', ( diag0(i), i = 1, nshape )
5319 WRITE(outnum, 2000) 'NMAT :', nmat
5320 WRITE(outnum, 2000) ' M :', ( m0(i), i = 1, nmat )
5321 WRITE(outnum, 2000) ' N :', ( n0(i), i = 1, nmat )
5322 WRITE(outnum, 2000) ' LDAS :', ( ldas0(i), i = 1, nmat )
5323 WRITE(outnum, 2000) ' LDAD :', ( ldad0(i), i = 1, nmat )
5324 WRITE(outnum, 2000) 'NSRC :', nsrc
5325 WRITE(outnum, 2000) ' RSRC :',( rsrc0(i), i = 1, nsrc )
5326 WRITE(outnum, 2000) ' CSRC :',( csrc0(i), i = 1, nsrc )
5327 WRITE(outnum, 2000) 'NGRIDS:', ngrid
5328 WRITE(outnum, 2000) ' P :', ( p0(i), i = 1, ngrid )
5329 WRITE(outnum, 2000) ' Q :', ( q0(i), i = 1, ngrid )
5330 WRITE(outnum, 2000) 'VERB :', verb
5331 WRITE(outnum,*) ' '
5332 END IF
5333 IF( verb .GT. 1 ) THEN
5334 WRITE(outnum,5000)
5335 WRITE(outnum,6000)
5336 END IF
5337 END IF
5338*
5339* Find biggest matrix, so we know where to stick error info
5340*
5341 i = 0
5342 DO 10 ima = 1, nmat
5343 k = n0(ima) * max0( ldas0(ima), ldad0(ima) ) + 4 * m0(ima)
5344 IF( k .GT. i ) i = k
5345 10 CONTINUE
5346 maxerr = ( csize * (memlen-i) ) / ( csize*2 + isize*6 )
5347 IF( maxerr .LT. 1 ) THEN
5348 WRITE(outnum,*) 'ERROR: Not enough memory to run BSBR tests.'
5349 CALL blacs_abort(-1, 1)
5350 END IF
5351 errdptr = i + 1
5352 erriptr = errdptr + maxerr
5353 nerr = 0
5354 testnum = 0
5355 nfail = 0
5356 nskip = 0
5357*
5358* Loop over grids of matrix
5359*
5360 DO 110 igr = 1, ngrid
5361*
5362 context = context0(igr)
5363 CALL blacs_gridinfo( context, nprow, npcol, myrow, mycol )
5364*
5365 ingrid = ( nprow .GT. 0 )
5366*
5367 DO 100 isc = 1, nscope
5368 scope = scope0(isc)
5369 DO 90 ito = 1, ntop
5370 top = top0(ito)
5371*
5372* If testing multipath ('M') or general tree ('T'),
5373* need to loop over calls to BLACS_SET
5374*
5375 IF( lsame(top, 'M') ) THEN
5376 setwhat = 11
5377 IF( scope .EQ. 'R' ) THEN
5378 istart = -(npcol - 1)
5379 istop = -istart
5380 ELSE IF (scope .EQ. 'C') THEN
5381 istart = -(nprow - 1)
5382 istop = -istart
5383 ELSE
5384 istart = -(nprow*npcol - 1)
5385 istop = -istart
5386 ENDIF
5387 ELSE IF( lsame(top, 'T') ) THEN
5388 setwhat = 12
5389 istart = 1
5390 IF( scope .EQ. 'R' ) THEN
5391 istop = npcol - 1
5392 ELSE IF (scope .EQ. 'C') THEN
5393 istop = nprow - 1
5394 ELSE
5395 istop = nprow*npcol - 1
5396 ENDIF
5397 ELSE
5398 setwhat = 0
5399 istart = 1
5400 istop = 1
5401 ENDIF
5402 DO 80 ish = 1, nshape
5403 uplo = uplo0(ish)
5404 diag = diag0(ish)
5405*
5406 DO 70 ima = 1, nmat
5407 m = m0(ima)
5408 n = n0(ima)
5409 ldasrc = ldas0(ima)
5410 ldadst = ldad0(ima)
5411*
5412 DO 60 iso = 1, nsrc
5413 testnum = testnum + 1
5414 rsrc = rsrc0(iso)
5415 csrc = csrc0(iso)
5416 IF( rsrc.GE.p0(igr) .OR. csrc.GE.q0(igr) ) THEN
5417 nskip = nskip + 1
5418 GOTO 60
5419 END IF
5420 IF( verb .GT. 1 ) THEN
5421 IF( iam .EQ. 0 ) THEN
5422 WRITE(outnum, 7000)
5423 $ testnum, 'RUNNING',scope, top, uplo, diag,
5424 $ m, n, ldasrc, ldadst, rsrc, csrc,
5425 $ nprow, npcol
5426 END IF
5427 END IF
5428*
5429 testok = .true.
5430 ipre = 2 * m
5431 ipost = ipre
5432 aptr = ipre + 1
5433*
5434* If I am in scope
5435*
5436 IF( (myrow.EQ.rsrc .AND. scope.EQ.'R') .OR.
5437 $ (mycol.EQ.csrc .AND. scope.EQ.'C') .OR.
5438 $ (scope .EQ. 'A') ) THEN
5439*
5440* source process generates matrix and sends it
5441*
5442 IF( myrow.EQ.rsrc .AND. mycol.EQ.csrc ) THEN
5443 CALL cinitmat(uplo, diag, m, n, mem,
5444 $ ldasrc, ipre, ipost,
5445 $ scheckval, testnum,
5446 $ myrow, mycol )
5447*
5448 DO 20 j = istart, istop
5449 IF( j.EQ.0 ) GOTO 20
5450 IF( setwhat.NE.0 )
5451 $ CALL blacs_set(context, setwhat, j)
5452 IF( uplo.EQ.'U' .OR. uplo.EQ.'L' ) THEN
5453 CALL ctrbs2d(context, scope, top,
5454 $ uplo, diag, m, n,
5455 $ mem(aptr), ldasrc )
5456 ELSE
5457 CALL cgebs2d(context, scope, top,
5458 $ m, n, mem(aptr),
5459 $ ldasrc )
5460 END IF
5461 20 CONTINUE
5462*
5463* Destination processes
5464*
5465 ELSE IF( ingrid ) THEN
5466 DO 40 j = istart, istop
5467 IF( j.EQ.0 ) GOTO 40
5468 IF( setwhat.NE.0 )
5469 $ CALL blacs_set(context, setwhat, j)
5470*
5471* Pad entire matrix area
5472*
5473 DO 30 k = 1, ipre+ipost+ldadst*n
5474 mem(k) = rcheckval
5475 30 CONTINUE
5476*
5477* Receive matrix
5478*
5479 IF( uplo.EQ.'U' .OR. uplo.EQ.'L' ) THEN
5480 CALL ctrbr2d(context, scope, top,
5481 $ uplo, diag, m, n,
5482 $ mem(aptr), ldadst,
5483 $ rsrc, csrc)
5484 ELSE
5485 CALL cgebr2d(context, scope, top,
5486 $ m, n, mem(aptr),
5487 $ ldadst, rsrc, csrc)
5488 END IF
5489*
5490* Check for errors in matrix or padding
5491*
5492 i = nerr
5493 CALL cchkmat(uplo, diag, m, n,
5494 $ mem(aptr), ldadst, rsrc, csrc,
5495 $ myrow, mycol, testnum, maxerr,
5496 $ nerr, mem(erriptr),
5497 $ mem(errdptr))
5498*
5499 CALL cchkpad(uplo, diag, m, n, mem,
5500 $ ldadst, rsrc, csrc, myrow,
5501 $ mycol, ipre, ipost, rcheckval,
5502 $ testnum, maxerr, nerr,
5503 $ mem(erriptr), mem(errdptr))
5504 40 CONTINUE
5505 testok = ( i .EQ. nerr )
5506 END IF
5507 END IF
5508*
5509 IF( verb .GT. 1 ) THEN
5510 i = nerr
5511 CALL cbtcheckin(0, outnum, maxerr, nerr,
5512 $ mem(erriptr), mem(errdptr),
5513 $ tfail)
5514 IF( iam .EQ. 0 ) THEN
5515 testok = ( testok .AND. (i.EQ.nerr) )
5516 IF( testok ) THEN
5517 WRITE(outnum,7000)testnum,'PASSED ',
5518 $ scope, top, uplo, diag, m, n,
5519 $ ldasrc, ldadst, rsrc, csrc,
5520 $ nprow, npcol
5521 ELSE
5522 nfail = nfail + 1
5523 WRITE(outnum,7000)testnum,'FAILED ',
5524 $ scope, top, uplo, diag, m, n,
5525 $ ldasrc, ldadst, rsrc, csrc,
5526 $ nprow, npcol
5527 END IF
5528 END IF
5529*
5530* Once we've printed out errors, can re-use buf space
5531*
5532 nerr = 0
5533 END IF
5534 60 CONTINUE
5535 70 CONTINUE
5536 80 CONTINUE
5537 90 CONTINUE
5538 100 CONTINUE
5539 110 CONTINUE
5540*
5541 IF( verb .LT. 2 ) THEN
5542 nfail = testnum
5543 CALL cbtcheckin( nfail, outnum, maxerr, nerr, mem(erriptr),
5544 $ mem(errdptr), tfail )
5545 END IF
5546 IF( iam .EQ. 0 ) THEN
5547 IF( verb .GT. 1 ) WRITE(outnum,*) ' '
5548 IF( nfail+nskip .EQ. 0 ) THEN
5549 WRITE(outnum, 8000 ) testnum
5550 ELSE
5551 WRITE(outnum, 9000 ) testnum, testnum-nskip-nfail,
5552 $ nskip, nfail
5553 END IF
5554 END IF
5555*
5556* Log whether their were any failures
5557*
5558 testok = allpass( (nfail.EQ.0) )
5559*
5560 1000 FORMAT('COMPLEX BSBR TESTS: BEGIN.' )
5561 2000 FORMAT(1x,a7,3x,10i6)
5562 3000 FORMAT(1x,a7,3x,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,
5563 $ 5x,a1,5x,a1)
5564 5000 FORMAT(' TEST# STATUS SCOPE TOP UPLO DIAG M N LDAS ',
5565 $ ' LDAD RSRC CSRC P Q')
5566 6000 FORMAT(' ----- ------- ----- --- ---- ---- ----- ----- ----- ',
5567 $ '----- ---- ---- ---- ----')
5568 7000 FORMAT(i6,1x,a7,5x,a1,3x,a1,2(4x,a1), 4i6, 4i5)
5569 8000 FORMAT('COMPLEX BSBR TESTS: PASSED ALL',
5570 $ i5, ' TESTS.')
5571 9000 FORMAT('COMPLEX BSBR TESTS:',i5,' TESTS;',i5,' PASSED,',
5572 $ i5,' SKIPPED,',i5,' FAILED.')
5573*
5574 RETURN
5575*
5576* End of CBSBRTEST.
5577*
float cmplx[2]
Definition pblas.h:136
subroutine cchkpad(uplo, diag, m, n, mem, lda, rsrc, csrc, myrow, mycol, ipre, ipost, checkval, testnum, maxerr, nerr, erribuf, errdbuf)
Definition blacstest.f:9872
logical function allpass(thistest)
Definition blacstest.f:1881
subroutine cinitmat(uplo, diag, m, n, mem, lda, ipre, ipost, checkval, testnum, myrow, mycol)
Definition blacstest.f:9591
subroutine cbtcheckin(nftests, outnum, maxerr, nerr, ierr, cval, tfailed)
Definition blacstest.f:9469
subroutine cchkmat(uplo, diag, m, n, a, lda, rsrc, csrc, myrow, mycol, testnum, maxerr, nerr, erribuf, errdbuf)
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 caller graph for this function: