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

◆ pclagen()

subroutine pclagen ( logical  inplace,
character*1  aform,
character*1  diag,
integer  offa,
integer  m,
integer  n,
integer  ia,
integer  ja,
integer, dimension( * )  desca,
integer  iaseed,
complex, dimension( lda, * )  a,
integer  lda 
)

Definition at line 8489 of file pcblastst.f.

8491*
8492* -- PBLAS test routine (version 2.0) --
8493* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
8494* and University of California, Berkeley.
8495* April 1, 1998
8496*
8497* .. Scalar Arguments ..
8498 LOGICAL INPLACE
8499 CHARACTER*1 AFORM, DIAG
8500 INTEGER IA, IASEED, JA, LDA, M, N, OFFA
8501* ..
8502* .. Array Arguments ..
8503 INTEGER DESCA( * )
8504 COMPLEX A( LDA, * )
8505* ..
8506*
8507* Purpose
8508* =======
8509*
8510* PCLAGEN generates (or regenerates) a submatrix sub( A ) denoting
8511* A(IA:IA+M-1,JA:JA+N-1).
8512*
8513* Notes
8514* =====
8515*
8516* A description vector is associated with each 2D block-cyclicly dis-
8517* tributed matrix. This vector stores the information required to
8518* establish the mapping between a matrix entry and its corresponding
8519* process and memory location.
8520*
8521* In the following comments, the character _ should be read as
8522* "of the distributed matrix". Let A be a generic term for any 2D
8523* block cyclicly distributed matrix. Its description vector is DESCA:
8524*
8525* NOTATION STORED IN EXPLANATION
8526* ---------------- --------------- ------------------------------------
8527* DTYPE_A (global) DESCA( DTYPE_ ) The descriptor type.
8528* CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
8529* the NPROW x NPCOL BLACS process grid
8530* A is distributed over. The context
8531* itself is global, but the handle
8532* (the integer value) may vary.
8533* M_A (global) DESCA( M_ ) The number of rows in the distribu-
8534* ted matrix A, M_A >= 0.
8535* N_A (global) DESCA( N_ ) The number of columns in the distri-
8536* buted matrix A, N_A >= 0.
8537* IMB_A (global) DESCA( IMB_ ) The number of rows of the upper left
8538* block of the matrix A, IMB_A > 0.
8539* INB_A (global) DESCA( INB_ ) The number of columns of the upper
8540* left block of the matrix A,
8541* INB_A > 0.
8542* MB_A (global) DESCA( MB_ ) The blocking factor used to distri-
8543* bute the last M_A-IMB_A rows of A,
8544* MB_A > 0.
8545* NB_A (global) DESCA( NB_ ) The blocking factor used to distri-
8546* bute the last N_A-INB_A columns of
8547* A, NB_A > 0.
8548* RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
8549* row of the matrix A is distributed,
8550* NPROW > RSRC_A >= 0.
8551* CSRC_A (global) DESCA( CSRC_ ) The process column over which the
8552* first column of A is distributed.
8553* NPCOL > CSRC_A >= 0.
8554* LLD_A (local) DESCA( LLD_ ) The leading dimension of the local
8555* array storing the local blocks of
8556* the distributed matrix A,
8557* IF( Lc( 1, N_A ) > 0 )
8558* LLD_A >= MAX( 1, Lr( 1, M_A ) )
8559* ELSE
8560* LLD_A >= 1.
8561*
8562* Let K be the number of rows of a matrix A starting at the global in-
8563* dex IA,i.e, A( IA:IA+K-1, : ). Lr( IA, K ) denotes the number of rows
8564* that the process of row coordinate MYROW ( 0 <= MYROW < NPROW ) would
8565* receive if these K rows were distributed over NPROW processes. If K
8566* is the number of columns of a matrix A starting at the global index
8567* JA, i.e, A( :, JA:JA+K-1, : ), Lc( JA, K ) denotes the number of co-
8568* lumns that the process MYCOL ( 0 <= MYCOL < NPCOL ) would receive if
8569* these K columns were distributed over NPCOL processes.
8570*
8571* The values of Lr() and Lc() may be determined via a call to the func-
8572* tion PB_NUMROC:
8573* Lr( IA, K ) = PB_NUMROC( K, IA, IMB_A, MB_A, MYROW, RSRC_A, NPROW )
8574* Lc( JA, K ) = PB_NUMROC( K, JA, INB_A, NB_A, MYCOL, CSRC_A, NPCOL )
8575*
8576* Arguments
8577* =========
8578*
8579* INPLACE (global input) LOGICAL
8580* On entry, INPLACE specifies if the matrix should be generated
8581* in place or not. If INPLACE is .TRUE., the local random array
8582* to be generated will start in memory at the local memory lo-
8583* cation A( 1, 1 ), otherwise it will start at the local posi-
8584* tion induced by IA and JA.
8585*
8586* AFORM (global input) CHARACTER*1
8587* On entry, AFORM specifies the type of submatrix to be genera-
8588* ted as follows:
8589* AFORM = 'S', sub( A ) is a symmetric matrix,
8590* AFORM = 'H', sub( A ) is a Hermitian matrix,
8591* AFORM = 'T', sub( A ) is overrwritten with the transpose
8592* of what would normally be generated,
8593* AFORM = 'C', sub( A ) is overwritten with the conjugate
8594* transpose of what would normally be genera-
8595* ted.
8596* AFORM = 'N', a random submatrix is generated.
8597*
8598* DIAG (global input) CHARACTER*1
8599* On entry, DIAG specifies if the generated submatrix is diago-
8600* nally dominant or not as follows:
8601* DIAG = 'D' : sub( A ) is diagonally dominant,
8602* DIAG = 'N' : sub( A ) is not diagonally dominant.
8603*
8604* OFFA (global input) INTEGER
8605* On entry, OFFA specifies the offdiagonal of the underlying
8606* matrix A(1:DESCA(M_),1:DESCA(N_)) of interest when the subma-
8607* trix is symmetric, Hermitian or diagonally dominant. OFFA = 0
8608* specifies the main diagonal, OFFA > 0 specifies a subdiago-
8609* nal, and OFFA < 0 specifies a superdiagonal (see further de-
8610* tails).
8611*
8612* M (global input) INTEGER
8613* On entry, M specifies the global number of matrix rows of the
8614* submatrix sub( A ) to be generated. M must be at least zero.
8615*
8616* N (global input) INTEGER
8617* On entry, N specifies the global number of matrix columns of
8618* the submatrix sub( A ) to be generated. N must be at least
8619* zero.
8620*
8621* IA (global input) INTEGER
8622* On entry, IA specifies A's global row index, which points to
8623* the beginning of the submatrix sub( A ).
8624*
8625* JA (global input) INTEGER
8626* On entry, JA specifies A's global column index, which points
8627* to the beginning of the submatrix sub( A ).
8628*
8629* DESCA (global and local input) INTEGER array
8630* On entry, DESCA is an integer array of dimension DLEN_. This
8631* is the array descriptor for the matrix A.
8632*
8633* IASEED (global input) INTEGER
8634* On entry, IASEED specifies the seed number to generate the
8635* matrix A. IASEED must be at least zero.
8636*
8637* A (local output) COMPLEX array
8638* On entry, A is an array of dimension (LLD_A, Ka), where Ka is
8639* at least Lc( 1, JA+N-1 ). On exit, this array contains the
8640* local entries of the randomly generated submatrix sub( A ).
8641*
8642* LDA (local input) INTEGER
8643* On entry, LDA specifies the local leading dimension of the
8644* array A. When INPLACE is .FALSE., LDA is usually DESCA(LLD_).
8645* This restriction is however not enforced, and this subroutine
8646* requires only that LDA >= MAX( 1, Mp ) where
8647*
8648* Mp = PB_NUMROC( M, IA, IMB_A, MB_A, MYROW, RSRC_A, NPROW ).
8649*
8650* PB_NUMROC is a ScaLAPACK tool function; MYROW, MYCOL, NPROW
8651* and NPCOL can be determined by calling the BLACS subroutine
8652* BLACS_GRIDINFO.
8653*
8654* Further Details
8655* ===============
8656*
8657* OFFD is tied to the matrix described by DESCA, as opposed to the
8658* piece that is currently (re)generated. This is a global information
8659* independent from the distribution parameters. Below are examples of
8660* the meaning of OFFD for a global 7 by 5 matrix:
8661*
8662* ---------------------------------------------------------------------
8663* OFFD | 0 -1 -2 -3 -4 0 -1 -2 -3 -4 0 -1 -2 -3 -4
8664* -------|-------------------------------------------------------------
8665* | | OFFD=-1 | OFFD=0 OFFD=2
8666* | V V
8667* 0 | . d . . . -> d . . . . . . . . .
8668* 1 | . . d . . . d . . . . . . . .
8669* 2 | . . . d . . . d . . -> d . . . .
8670* 3 | . . . . d . . . d . . d . . .
8671* 4 | . . . . . . . . . d . . d . .
8672* 5 | . . . . . . . . . . . . . d .
8673* 6 | . . . . . . . . . . . . . . d
8674* ---------------------------------------------------------------------
8675*
8676* -- Written on April 1, 1998 by
8677* Antoine Petitet, University of Tennessee, Knoxville 37996, USA.
8678*
8679* =====================================================================
8680*
8681* .. Parameters ..
8682 INTEGER BLOCK_CYCLIC_2D_INB, CSRC_, CTXT_, DLEN_,
8683 $ DTYPE_, IMB_, INB_, LLD_, MB_, M_, NB_, N_,
8684 $ RSRC_
8685 parameter( block_cyclic_2d_inb = 2, dlen_ = 11,
8686 $ dtype_ = 1, ctxt_ = 2, m_ = 3, n_ = 4,
8687 $ imb_ = 5, inb_ = 6, mb_ = 7, nb_ = 8,
8688 $ rsrc_ = 9, csrc_ = 10, lld_ = 11 )
8689 INTEGER JMP_1, JMP_COL, JMP_IMBV, JMP_INBV, JMP_LEN,
8690 $ JMP_MB, JMP_NB, JMP_NPIMBLOC, JMP_NPMB,
8691 $ JMP_NQINBLOC, JMP_NQNB, JMP_ROW
8692 parameter( jmp_1 = 1, jmp_row = 2, jmp_col = 3,
8693 $ jmp_mb = 4, jmp_imbv = 5, jmp_npmb = 6,
8694 $ jmp_npimbloc = 7, jmp_nb = 8, jmp_inbv = 9,
8695 $ jmp_nqnb = 10, jmp_nqinbloc = 11,
8696 $ jmp_len = 11 )
8697 REAL ZERO
8698 parameter( zero = 0.0e+0 )
8699* ..
8700* .. Local Scalars ..
8701 LOGICAL DIAGDO, SYMM, HERM, NOTRAN
8702 INTEGER CSRC, I, IACOL, IAROW, ICTXT, IIA, ILOCBLK,
8703 $ ILOCOFF, ILOW, IMB, IMB1, IMBLOC, IMBVIR, INB,
8704 $ INB1, INBLOC, INBVIR, INFO, IOFFDA, ITMP, IUPP,
8705 $ IVIR, JJA, JLOCBLK, JLOCOFF, JVIR, LCMT00,
8706 $ LMBLOC, LNBLOC, LOW, MAXMN, MB, MBLKS, MP,
8707 $ MRCOL, MRROW, MYCDIST, MYCOL, MYRDIST, MYROW,
8708 $ NB, NBLKS, NPCOL, NPROW, NQ, NVIR, RSRC, UPP
8709 COMPLEX ALPHA
8710* ..
8711* .. Local Arrays ..
8712 INTEGER DESCA2( DLEN_ ), IMULADD( 4, JMP_LEN ),
8713 $ IRAN( 2 ), JMP( JMP_LEN ), MULADD0( 4 )
8714* ..
8715* .. External Subroutines ..
8716 EXTERNAL blacs_gridinfo, pb_ainfog2l, pb_binfo,
8720* ..
8721* .. External Functions ..
8722 LOGICAL LSAME
8723 EXTERNAL lsame
8724* ..
8725* .. Intrinsic Functions ..
8726 INTRINSIC cmplx, max, min, real
8727* ..
8728* .. Data Statements ..
8729 DATA ( muladd0( i ), i = 1, 4 ) / 20077, 16838,
8730 $ 12345, 0 /
8731* ..
8732* .. Executable Statements ..
8733*
8734* Convert descriptor
8735*
8736 CALL pb_desctrans( desca, desca2 )
8737*
8738* Test the input arguments
8739*
8740 ictxt = desca2( ctxt_ )
8741 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
8742*
8743* Test the input parameters
8744*
8745 info = 0
8746 IF( nprow.EQ.-1 ) THEN
8747 info = -( 1000 + ctxt_ )
8748 ELSE
8749 symm = lsame( aform, 'S' )
8750 herm = lsame( aform, 'H' )
8751 notran = lsame( aform, 'N' )
8752 diagdo = lsame( diag, 'D' )
8753 IF( .NOT.( symm.OR.herm.OR.notran ) .AND.
8754 $ .NOT.( lsame( aform, 'T' ) ) .AND.
8755 $ .NOT.( lsame( aform, 'C' ) ) ) THEN
8756 info = -2
8757 ELSE IF( ( .NOT.diagdo ) .AND.
8758 $ ( .NOT.lsame( diag, 'N' ) ) ) THEN
8759 info = -3
8760 END IF
8761 CALL pb_chkmat( ictxt, m, 5, n, 6, ia, ja, desca2, 10, info )
8762 END IF
8763*
8764 IF( info.NE.0 ) THEN
8765 CALL pxerbla( ictxt, 'PCLAGEN', -info )
8766 RETURN
8767 END IF
8768*
8769* Quick return if possible
8770*
8771 IF( ( m.LE.0 ).OR.( n.LE.0 ) )
8772 $ RETURN
8773*
8774* Start the operations
8775*
8776 mb = desca2( mb_ )
8777 nb = desca2( nb_ )
8778 imb = desca2( imb_ )
8779 inb = desca2( inb_ )
8780 rsrc = desca2( rsrc_ )
8781 csrc = desca2( csrc_ )
8782*
8783* Figure out local information about the distributed matrix operand
8784*
8785 CALL pb_ainfog2l( m, n, ia, ja, desca2, nprow, npcol, myrow,
8786 $ mycol, imb1, inb1, mp, nq, iia, jja, iarow,
8787 $ iacol, mrrow, mrcol )
8788*
8789* Decide where the entries shall be stored in memory
8790*
8791 IF( inplace ) THEN
8792 iia = 1
8793 jja = 1
8794 END IF
8795*
8796* Initialize LCMT00, MBLKS, NBLKS, IMBLOC, INBLOC, LMBLOC, LNBLOC,
8797* ILOW, LOW, IUPP, and UPP.
8798*
8799 ioffda = ja + offa - ia
8800 CALL pb_binfo( ioffda, mp, nq, imb1, inb1, mb, nb, mrrow,
8801 $ mrcol, lcmt00, mblks, nblks, imbloc, inbloc,
8802 $ lmbloc, lnbloc, ilow, low, iupp, upp )
8803*
8804* Initialize ILOCBLK, ILOCOFF, MYRDIST, JLOCBLK, JLOCOFF, MYCDIST
8805* This values correspond to the square virtual underlying matrix
8806* of size MAX( M_ + MAX( 0, -OFFA ), N_ + MAX( 0, OFFA ) ) used
8807* to set up the random sequence. For practical purposes, the size
8808* of this virtual matrix is upper bounded by M_ + N_ - 1.
8809*
8810 itmp = max( 0, -offa )
8811 ivir = ia + itmp
8812 imbvir = imb + itmp
8813 nvir = desca2( m_ ) + itmp
8814*
8815 CALL pb_locinfo( ivir, imbvir, mb, myrow, rsrc, nprow, ilocblk,
8816 $ ilocoff, myrdist )
8817*
8818 itmp = max( 0, offa )
8819 jvir = ja + itmp
8820 inbvir = inb + itmp
8821 nvir = max( max( nvir, desca2( n_ ) + itmp ),
8822 $ desca2( m_ ) + desca2( n_ ) - 1 )
8823*
8824 CALL pb_locinfo( jvir, inbvir, nb, mycol, csrc, npcol, jlocblk,
8825 $ jlocoff, mycdist )
8826*
8827 IF( symm .OR. herm .OR. notran ) THEN
8828*
8829 CALL pb_initjmp( .true., nvir, imbvir, inbvir, imbloc, inbloc,
8830 $ mb, nb, rsrc, csrc, nprow, npcol, 2, jmp )
8831*
8832* Compute constants to jump JMP( * ) numbers in the sequence
8833*
8834 CALL pb_initmuladd( muladd0, jmp, imuladd )
8835*
8836* Compute and set the random value corresponding to A( IA, JA )
8837*
8838 CALL pb_setlocran( iaseed, ilocblk, jlocblk, ilocoff, jlocoff,
8839 $ myrdist, mycdist, nprow, npcol, jmp,
8840 $ imuladd, iran )
8841*
8842 CALL pb_clagen( 'Lower', aform, a( iia, jja ), lda, lcmt00,
8843 $ iran, mblks, imbloc, mb, lmbloc, nblks, inbloc,
8844 $ nb, lnbloc, jmp, imuladd )
8845*
8846 END IF
8847*
8848 IF( symm .OR. herm .OR. ( .NOT. notran ) ) THEN
8849*
8850 CALL pb_initjmp( .false., nvir, imbvir, inbvir, imbloc, inbloc,
8851 $ mb, nb, rsrc, csrc, nprow, npcol, 2, jmp )
8852*
8853* Compute constants to jump JMP( * ) numbers in the sequence
8854*
8855 CALL pb_initmuladd( muladd0, jmp, imuladd )
8856*
8857* Compute and set the random value corresponding to A( IA, JA )
8858*
8859 CALL pb_setlocran( iaseed, ilocblk, jlocblk, ilocoff, jlocoff,
8860 $ myrdist, mycdist, nprow, npcol, jmp,
8861 $ imuladd, iran )
8862*
8863 CALL pb_clagen( 'Upper', aform, a( iia, jja ), lda, lcmt00,
8864 $ iran, mblks, imbloc, mb, lmbloc, nblks, inbloc,
8865 $ nb, lnbloc, jmp, imuladd )
8866*
8867 END IF
8868*
8869 IF( diagdo ) THEN
8870*
8871 maxmn = max( desca2( m_ ), desca2( n_ ) )
8872 IF( herm ) THEN
8873 alpha = cmplx( real( 2 * maxmn ), zero )
8874 ELSE
8875 alpha = cmplx( real( maxmn ), real( maxmn ) )
8876 END IF
8877*
8878 IF( ioffda.GE.0 ) THEN
8879 CALL pcladom( inplace, min( max( 0, m-ioffda ), n ), alpha,
8880 $ a, min( ia+ioffda, ia+m-1 ), ja, desca )
8881 ELSE
8882 CALL pcladom( inplace, min( m, max( 0, n+ioffda ) ), alpha,
8883 $ a, ia, min( ja-ioffda, ja+n-1 ), desca )
8884 END IF
8885*
8886 END IF
8887*
8888 RETURN
8889*
8890* End of PCLAGEN
8891*
float cmplx[2]
Definition pblas.h:136
subroutine pb_ainfog2l(m, n, i, j, desc, nprow, npcol, myrow, mycol, imb1, inb1, mp, nq, ii, jj, prow, pcol, rprow, rpcol)
Definition pblastst.f:2023
subroutine pb_binfo(offd, m, n, imb1, inb1, mb, nb, mrrow, mrcol, lcmt00, mblks, nblks, imbloc, inbloc, lmbloc, lnbloc, ilow, low, iupp, upp)
Definition pblastst.f:3577
subroutine pb_setran(iran, iac)
Definition pblastst.f:4759
subroutine pb_locinfo(i, inb, nb, myroc, srcproc, nprocs, ilocblk, ilocoff, mydist)
Definition pblastst.f:3910
subroutine pb_chkmat(ictxt, m, mpos0, n, npos0, ia, ja, desca, dpos0, info)
Definition pblastst.f:2742
subroutine pb_jump(k, muladd, irann, iranm, ima)
Definition pblastst.f:4648
subroutine pb_setlocran(seed, ilocblk, jlocblk, ilocoff, jlocoff, myrdist, mycdist, nprow, npcol, jmp, imuladd, iran)
Definition pblastst.f:4302
subroutine pb_initmuladd(muladd0, jmp, imuladd)
Definition pblastst.f:4196
subroutine pb_desctrans(descin, descout)
Definition pblastst.f:2964
subroutine pb_initjmp(colmaj, nvir, imbvir, inbvir, imbloc, inbloc, mb, nb, rsrc, csrc, nprow, npcol, stride, jmp)
Definition pblastst.f:4045
subroutine pb_jumpit(muladd, irann, iranm)
Definition pblastst.f:4822
subroutine pcladom(inplace, n, alpha, a, ia, ja, desca)
Definition pcblastst.f:8894
subroutine pb_clagen(uplo, aform, a, lda, lcmt00, iran, mblks, imbloc, mb, lmbloc, nblks, inbloc, nb, lnbloc, jmp, imuladd)
#define max(A, B)
Definition pcgemr.c:180
#define min(A, B)
Definition pcgemr.c:181
subroutine pxerbla(ictxt, srname, info)
Definition pxerbla.f:2
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: