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

◆ zbsbrtest()

subroutine zbsbrtest ( 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,
double complex, dimension(memlen)  mem,
integer  memlen 
)

Definition at line 5581 of file blacstest.f.

5585*
5586* -- BLACS tester (version 1.0) --
5587* University of Tennessee
5588* December 15, 1994
5589*
5590*
5591* .. Scalar Arguments ..
5592 INTEGER OUTNUM, VERB, NSCOPE, NTOP, NSHAPE, NMAT, NSRC, NGRID
5593 INTEGER MEMLEN
5594* ..
5595* .. Array Arguments ..
5596 CHARACTER*1 SCOPE0(NSCOPE), TOP0(NTOP)
5597 CHARACTER*1 UPLO0(NSHAPE), DIAG0(NSHAPE)
5598 INTEGER M0(NMAT), N0(NMAT), LDAS0(NMAT), LDAD0(NMAT)
5599 INTEGER RSRC0(NSRC), CSRC0(NSRC), CONTEXT0(NGRID)
5600 INTEGER P0(NGRID), Q0(NGRID), TFAIL(*)
5601 DOUBLE COMPLEX MEM(MEMLEN)
5602* ..
5603*
5604* Purpose
5605* =======
5606* ZTESTBSBR: Test double complex broadcast
5607*
5608* Arguments
5609* =========
5610* OUTNUM (input) INTEGER
5611* The device number to write output to.
5612*
5613* VERB (input) INTEGER
5614* The level of verbosity (how much printing to do).
5615*
5616* NSCOPE (input) INTEGER
5617* The number of scopes to be tested.
5618*
5619* SCOPE0 (input) CHARACTER*1 array of dimension (NSCOPE)
5620* Values of the scopes to be tested.
5621*
5622* NTOP (input) INTEGER
5623* The number of topologies to be tested.
5624*
5625* TOP0 (input) CHARACTER*1 array of dimension (NTOP)
5626* Values of the topologies to be tested.
5627*
5628* NSHAPE (input) INTEGER
5629* The number of matrix shapes to be tested.
5630*
5631* UPLO0 (input) CHARACTER*1 array of dimension (NSHAPE)
5632* Values of UPLO to be tested.
5633*
5634* DIAG0 (input) CHARACTER*1 array of dimension (NSHAPE)
5635* Values of DIAG to be tested.
5636*
5637* NMAT (input) INTEGER
5638* The number of matrices to be tested.
5639*
5640* M0 (input) INTEGER array of dimension (NMAT)
5641* Values of M to be tested.
5642*
5643* M0 (input) INTEGER array of dimension (NMAT)
5644* Values of M to be tested.
5645*
5646* N0 (input) INTEGER array of dimension (NMAT)
5647* Values of N to be tested.
5648*
5649* LDAS0 (input) INTEGER array of dimension (NMAT)
5650* Values of LDAS (leading dimension of A on source process)
5651* to be tested.
5652*
5653* LDAD0 (input) INTEGER array of dimension (NMAT)
5654* Values of LDAD (leading dimension of A on destination
5655* process) to be tested.
5656* NSRC (input) INTEGER
5657* The number of sources to be tested.
5658*
5659* RSRC0 (input) INTEGER array of dimension (NDEST)
5660* Values of RSRC (row coordinate of source) to be tested.
5661*
5662* CSRC0 (input) INTEGER array of dimension (NDEST)
5663* Values of CSRC (column coordinate of source) to be tested.
5664*
5665* NGRID (input) INTEGER
5666* The number of process grids to be tested.
5667*
5668* CONTEXT0 (input) INTEGER array of dimension (NGRID)
5669* The BLACS context handles corresponding to the grids.
5670*
5671* P0 (input) INTEGER array of dimension (NGRID)
5672* Values of P (number of process rows, NPROW).
5673*
5674* Q0 (input) INTEGER array of dimension (NGRID)
5675* Values of Q (number of process columns, NPCOL).
5676*
5677* TFAIL (workspace) INTEGER array of dimension (NTESTS)
5678* If VERB < 2, serves to indicate which tests fail. This
5679* requires workspace of NTESTS (number of tests performed).
5680*
5681* MEM (workspace) DOUBLE COMPLEX array of dimension (MEMLEN)
5682* Used for all other workspaces, including the matrix A,
5683* and its pre and post padding.
5684*
5685* MEMLEN (input) INTEGER
5686* The length, in elements, of MEM.
5687*
5688* =====================================================================
5689*
5690* .. External Functions ..
5691 LOGICAL ALLPASS, LSAME
5692 INTEGER IBTMYPROC, IBTSIZEOF
5693 EXTERNAL allpass, lsame, ibtmyproc, ibtsizeof
5694* ..
5695* .. External Subroutines ..
5696 EXTERNAL blacs_gridinfo
5697 EXTERNAL ztrbs2d, zgebs2d, ztrbr2d, zgebr2d
5699* ..
5700* .. Local Scalars ..
5701 CHARACTER*1 SCOPE, TOP, UPLO, DIAG
5702 LOGICAL TESTOK, INGRID
5703 INTEGER IAM, I, K, J, IGR, ISH, IMA, ISO, ISC, ITO
5704 INTEGER M, N, NPROW, NPCOL, MYROW, MYCOL, RSRC, CSRC
5705 INTEGER ISTART, ISTOP, IPRE, IPOST, SETWHAT
5706 INTEGER NERR, NSKIP, NFAIL, TESTNUM, CONTEXT, MAXERR, LDASRC
5707 INTEGER LDADST, ERRDPTR, APTR, ERRIPTR, ISIZE, ZSIZE
5708 DOUBLE COMPLEX SCHECKVAL, RCHECKVAL
5709* ..
5710* .. Executable Statements ..
5711*
5712 scheckval = dcmplx( -0.01d0, -0.01d0 )
5713 rcheckval = dcmplx( -0.02d0, -0.02d0 )
5714*
5715 iam = ibtmyproc()
5716 isize = ibtsizeof('I')
5717 zsize = ibtsizeof('Z')
5718*
5719* Verify file parameters
5720*
5721 IF( iam .EQ. 0 ) THEN
5722 WRITE(outnum, *) ' '
5723 WRITE(outnum, *) ' '
5724 WRITE(outnum, 1000 )
5725 IF( verb .GT. 0 ) THEN
5726 WRITE(outnum,*) ' '
5727 WRITE(outnum, 2000) 'NSCOPE:', nscope
5728 WRITE(outnum, 3000) ' SCOPE:', ( scope0(i), i = 1, nscope )
5729 WRITE(outnum, 2000) 'NTOP :', ntop
5730 WRITE(outnum, 3000) ' TOP :', ( top0(i), i = 1, ntop )
5731 WRITE(outnum, 2000) 'NSHAPE:', nshape
5732 WRITE(outnum, 3000) ' UPLO :', ( uplo0(i), i = 1, nshape )
5733 WRITE(outnum, 3000) ' DIAG :', ( diag0(i), i = 1, nshape )
5734 WRITE(outnum, 2000) 'NMAT :', nmat
5735 WRITE(outnum, 2000) ' M :', ( m0(i), i = 1, nmat )
5736 WRITE(outnum, 2000) ' N :', ( n0(i), i = 1, nmat )
5737 WRITE(outnum, 2000) ' LDAS :', ( ldas0(i), i = 1, nmat )
5738 WRITE(outnum, 2000) ' LDAD :', ( ldad0(i), i = 1, nmat )
5739 WRITE(outnum, 2000) 'NSRC :', nsrc
5740 WRITE(outnum, 2000) ' RSRC :',( rsrc0(i), i = 1, nsrc )
5741 WRITE(outnum, 2000) ' CSRC :',( csrc0(i), i = 1, nsrc )
5742 WRITE(outnum, 2000) 'NGRIDS:', ngrid
5743 WRITE(outnum, 2000) ' P :', ( p0(i), i = 1, ngrid )
5744 WRITE(outnum, 2000) ' Q :', ( q0(i), i = 1, ngrid )
5745 WRITE(outnum, 2000) 'VERB :', verb
5746 WRITE(outnum,*) ' '
5747 END IF
5748 IF( verb .GT. 1 ) THEN
5749 WRITE(outnum,5000)
5750 WRITE(outnum,6000)
5751 END IF
5752 END IF
5753*
5754* Find biggest matrix, so we know where to stick error info
5755*
5756 i = 0
5757 DO 10 ima = 1, nmat
5758 k = n0(ima) * max0( ldas0(ima), ldad0(ima) ) + 4 * m0(ima)
5759 IF( k .GT. i ) i = k
5760 10 CONTINUE
5761 maxerr = ( zsize * (memlen-i) ) / ( zsize*2 + isize*6 )
5762 IF( maxerr .LT. 1 ) THEN
5763 WRITE(outnum,*) 'ERROR: Not enough memory to run BSBR tests.'
5764 CALL blacs_abort(-1, 1)
5765 END IF
5766 errdptr = i + 1
5767 erriptr = errdptr + maxerr
5768 nerr = 0
5769 testnum = 0
5770 nfail = 0
5771 nskip = 0
5772*
5773* Loop over grids of matrix
5774*
5775 DO 110 igr = 1, ngrid
5776*
5777 context = context0(igr)
5778 CALL blacs_gridinfo( context, nprow, npcol, myrow, mycol )
5779*
5780 ingrid = ( nprow .GT. 0 )
5781*
5782 DO 100 isc = 1, nscope
5783 scope = scope0(isc)
5784 DO 90 ito = 1, ntop
5785 top = top0(ito)
5786*
5787* If testing multipath ('M') or general tree ('T'),
5788* need to loop over calls to BLACS_SET
5789*
5790 IF( lsame(top, 'M') ) THEN
5791 setwhat = 11
5792 IF( scope .EQ. 'R' ) THEN
5793 istart = -(npcol - 1)
5794 istop = -istart
5795 ELSE IF (scope .EQ. 'C') THEN
5796 istart = -(nprow - 1)
5797 istop = -istart
5798 ELSE
5799 istart = -(nprow*npcol - 1)
5800 istop = -istart
5801 ENDIF
5802 ELSE IF( lsame(top, 'T') ) THEN
5803 setwhat = 12
5804 istart = 1
5805 IF( scope .EQ. 'R' ) THEN
5806 istop = npcol - 1
5807 ELSE IF (scope .EQ. 'C') THEN
5808 istop = nprow - 1
5809 ELSE
5810 istop = nprow*npcol - 1
5811 ENDIF
5812 ELSE
5813 setwhat = 0
5814 istart = 1
5815 istop = 1
5816 ENDIF
5817 DO 80 ish = 1, nshape
5818 uplo = uplo0(ish)
5819 diag = diag0(ish)
5820*
5821 DO 70 ima = 1, nmat
5822 m = m0(ima)
5823 n = n0(ima)
5824 ldasrc = ldas0(ima)
5825 ldadst = ldad0(ima)
5826*
5827 DO 60 iso = 1, nsrc
5828 testnum = testnum + 1
5829 rsrc = rsrc0(iso)
5830 csrc = csrc0(iso)
5831 IF( rsrc.GE.p0(igr) .OR. csrc.GE.q0(igr) ) THEN
5832 nskip = nskip + 1
5833 GOTO 60
5834 END IF
5835 IF( verb .GT. 1 ) THEN
5836 IF( iam .EQ. 0 ) THEN
5837 WRITE(outnum, 7000)
5838 $ testnum, 'RUNNING',scope, top, uplo, diag,
5839 $ m, n, ldasrc, ldadst, rsrc, csrc,
5840 $ nprow, npcol
5841 END IF
5842 END IF
5843*
5844 testok = .true.
5845 ipre = 2 * m
5846 ipost = ipre
5847 aptr = ipre + 1
5848*
5849* If I am in scope
5850*
5851 IF( (myrow.EQ.rsrc .AND. scope.EQ.'R') .OR.
5852 $ (mycol.EQ.csrc .AND. scope.EQ.'C') .OR.
5853 $ (scope .EQ. 'A') ) THEN
5854*
5855* source process generates matrix and sends it
5856*
5857 IF( myrow.EQ.rsrc .AND. mycol.EQ.csrc ) THEN
5858 CALL zinitmat(uplo, diag, m, n, mem,
5859 $ ldasrc, ipre, ipost,
5860 $ scheckval, testnum,
5861 $ myrow, mycol )
5862*
5863 DO 20 j = istart, istop
5864 IF( j.EQ.0 ) GOTO 20
5865 IF( setwhat.NE.0 )
5866 $ CALL blacs_set(context, setwhat, j)
5867 IF( uplo.EQ.'U' .OR. uplo.EQ.'L' ) THEN
5868 CALL ztrbs2d(context, scope, top,
5869 $ uplo, diag, m, n,
5870 $ mem(aptr), ldasrc )
5871 ELSE
5872 CALL zgebs2d(context, scope, top,
5873 $ m, n, mem(aptr),
5874 $ ldasrc )
5875 END IF
5876 20 CONTINUE
5877*
5878* Destination processes
5879*
5880 ELSE IF( ingrid ) THEN
5881 DO 40 j = istart, istop
5882 IF( j.EQ.0 ) GOTO 40
5883 IF( setwhat.NE.0 )
5884 $ CALL blacs_set(context, setwhat, j)
5885*
5886* Pad entire matrix area
5887*
5888 DO 30 k = 1, ipre+ipost+ldadst*n
5889 mem(k) = rcheckval
5890 30 CONTINUE
5891*
5892* Receive matrix
5893*
5894 IF( uplo.EQ.'U' .OR. uplo.EQ.'L' ) THEN
5895 CALL ztrbr2d(context, scope, top,
5896 $ uplo, diag, m, n,
5897 $ mem(aptr), ldadst,
5898 $ rsrc, csrc)
5899 ELSE
5900 CALL zgebr2d(context, scope, top,
5901 $ m, n, mem(aptr),
5902 $ ldadst, rsrc, csrc)
5903 END IF
5904*
5905* Check for errors in matrix or padding
5906*
5907 i = nerr
5908 CALL zchkmat(uplo, diag, m, n,
5909 $ mem(aptr), ldadst, rsrc, csrc,
5910 $ myrow, mycol, testnum, maxerr,
5911 $ nerr, mem(erriptr),
5912 $ mem(errdptr))
5913*
5914 CALL zchkpad(uplo, diag, m, n, mem,
5915 $ ldadst, rsrc, csrc, myrow,
5916 $ mycol, ipre, ipost, rcheckval,
5917 $ testnum, maxerr, nerr,
5918 $ mem(erriptr), mem(errdptr))
5919 40 CONTINUE
5920 testok = ( i .EQ. nerr )
5921 END IF
5922 END IF
5923*
5924 IF( verb .GT. 1 ) THEN
5925 i = nerr
5926 CALL zbtcheckin(0, outnum, maxerr, nerr,
5927 $ mem(erriptr), mem(errdptr),
5928 $ tfail)
5929 IF( iam .EQ. 0 ) THEN
5930 testok = ( testok .AND. (i.EQ.nerr) )
5931 IF( testok ) THEN
5932 WRITE(outnum,7000)testnum,'PASSED ',
5933 $ scope, top, uplo, diag, m, n,
5934 $ ldasrc, ldadst, rsrc, csrc,
5935 $ nprow, npcol
5936 ELSE
5937 nfail = nfail + 1
5938 WRITE(outnum,7000)testnum,'FAILED ',
5939 $ scope, top, uplo, diag, m, n,
5940 $ ldasrc, ldadst, rsrc, csrc,
5941 $ nprow, npcol
5942 END IF
5943 END IF
5944*
5945* Once we've printed out errors, can re-use buf space
5946*
5947 nerr = 0
5948 END IF
5949 60 CONTINUE
5950 70 CONTINUE
5951 80 CONTINUE
5952 90 CONTINUE
5953 100 CONTINUE
5954 110 CONTINUE
5955*
5956 IF( verb .LT. 2 ) THEN
5957 nfail = testnum
5958 CALL zbtcheckin( nfail, outnum, maxerr, nerr, mem(erriptr),
5959 $ mem(errdptr), tfail )
5960 END IF
5961 IF( iam .EQ. 0 ) THEN
5962 IF( verb .GT. 1 ) WRITE(outnum,*) ' '
5963 IF( nfail+nskip .EQ. 0 ) THEN
5964 WRITE(outnum, 8000 ) testnum
5965 ELSE
5966 WRITE(outnum, 9000 ) testnum, testnum-nskip-nfail,
5967 $ nskip, nfail
5968 END IF
5969 END IF
5970*
5971* Log whether their were any failures
5972*
5973 testok = allpass( (nfail.EQ.0) )
5974*
5975 1000 FORMAT('DOUBLE COMPLEX BSBR TESTS: BEGIN.' )
5976 2000 FORMAT(1x,a7,3x,10i6)
5977 3000 FORMAT(1x,a7,3x,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,5x,a1,
5978 $ 5x,a1,5x,a1)
5979 5000 FORMAT(' TEST# STATUS SCOPE TOP UPLO DIAG M N LDAS ',
5980 $ ' LDAD RSRC CSRC P Q')
5981 6000 FORMAT(' ----- ------- ----- --- ---- ---- ----- ----- ----- ',
5982 $ '----- ---- ---- ---- ----')
5983 7000 FORMAT(i6,1x,a7,5x,a1,3x,a1,2(4x,a1), 4i6, 4i5)
5984 8000 FORMAT('DOUBLE COMPLEX BSBR TESTS: PASSED ALL',
5985 $ i5, ' TESTS.')
5986 9000 FORMAT('DOUBLE COMPLEX BSBR TESTS:',i5,' TESTS;',i5,' PASSED,',
5987 $ i5,' SKIPPED,',i5,' FAILED.')
5988*
5989 RETURN
5990*
5991* End of ZBSBRTEST.
5992*
logical function allpass(thistest)
Definition blacstest.f:1881
subroutine zchkmat(uplo, diag, m, n, a, lda, rsrc, csrc, myrow, mycol, testnum, maxerr, nerr, erribuf, errdbuf)
subroutine zinitmat(uplo, diag, m, n, mem, lda, ipre, ipost, checkval, testnum, myrow, mycol)
subroutine zchkpad(uplo, diag, m, n, mem, lda, rsrc, csrc, myrow, mycol, ipre, ipost, checkval, testnum, maxerr, nerr, erribuf, errdbuf)
subroutine zbtcheckin(nftests, outnum, maxerr, nerr, ierr, zval, tfailed)
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: