SCALAPACK 2.2.2
LAPACK: Linear Algebra PACKage
All Classes Files Functions Variables Typedefs Macros
infog2l.f
Go to the documentation of this file.
1 SUBROUTINE infog2l( GRINDX, GCINDX, DESC, NPROW, NPCOL, MYROW,
2 $ MYCOL, LRINDX, LCINDX, RSRC, CSRC )
3*
4* -- ScaLAPACK tools routine (version 1.7) --
5* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
6* and University of California, Berkeley.
7* May 1, 1997
8*
9* .. Scalar Arguments ..
10 INTEGER CSRC, GCINDX, GRINDX, LRINDX, LCINDX, MYCOL,
11 $ myrow, npcol, nprow, rsrc
12* ..
13* .. Array Arguments ..
14 INTEGER DESC( * )
15* ..
16*
17* Purpose
18* =======
19*
20* INFOG2L computes the starting local indexes LRINDX, LCINDX corres-
21* ponding to the distributed submatrix starting globally at the entry
22* pointed by GRINDX, GCINDX. This routine returns the coordinates in
23* the grid of the process owning the matrix entry of global indexes
24* GRINDX, GCINDX, namely RSRC and CSRC.
25*
26* Notes
27* =====
28*
29* Each global data object is described by an associated description
30* vector. This vector stores the information required to establish
31* the mapping between an object element and its corresponding process
32* and memory location.
33*
34* Let A be a generic term for any 2D block cyclicly distributed array.
35* Such a global array has an associated description vector DESCA.
36* In the following comments, the character _ should be read as
37* "of the global array".
38*
39* NOTATION STORED IN EXPLANATION
40* --------------- -------------- --------------------------------------
41* DTYPE_A(global) DESCA( DTYPE_ )The descriptor type. In this case,
42* DTYPE_A = 1.
43* CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
44* the BLACS process grid A is distribu-
45* ted over. The context itself is glo-
46* bal, but the handle (the integer
47* value) may vary.
48* M_A (global) DESCA( M_ ) The number of rows in the global
49* array A.
50* N_A (global) DESCA( N_ ) The number of columns in the global
51* array A.
52* MB_A (global) DESCA( MB_ ) The blocking factor used to distribute
53* the rows of the array.
54* NB_A (global) DESCA( NB_ ) The blocking factor used to distribute
55* the columns of the array.
56* RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
57* row of the array A is distributed.
58* CSRC_A (global) DESCA( CSRC_ ) The process column over which the
59* first column of the array A is
60* distributed.
61* LLD_A (local) DESCA( LLD_ ) The leading dimension of the local
62* array. LLD_A >= MAX(1,LOCr(M_A)).
63*
64* Let K be the number of rows or columns of a distributed matrix,
65* and assume that its process grid has dimension p x q.
66* LOCr( K ) denotes the number of elements of K that a process
67* would receive if K were distributed over the p processes of its
68* process column.
69* Similarly, LOCc( K ) denotes the number of elements of K that a
70* process would receive if K were distributed over the q processes of
71* its process row.
72* The values of LOCr() and LOCc() may be determined via a call to the
73* ScaLAPACK tool function, NUMROC:
74* LOCr( M ) = NUMROC( M, MB_A, MYROW, RSRC_A, NPROW ),
75* LOCc( N ) = NUMROC( N, NB_A, MYCOL, CSRC_A, NPCOL ).
76* An upper bound for these quantities may be computed by:
77* LOCr( M ) <= ceil( ceil(M/MB_A)/NPROW )*MB_A
78* LOCc( N ) <= ceil( ceil(N/NB_A)/NPCOL )*NB_A
79*
80* Arguments
81* =========
82*
83* GRINDX (global input) INTEGER
84* The global row starting index of the submatrix.
85*
86* GCINDX (global input) INTEGER
87* The global column starting index of the submatrix.
88*
89* DESC (input) INTEGER array of dimension DLEN_.
90* The array descriptor for the underlying distributed matrix.
91*
92* NPROW (global input) INTEGER
93* The total number of process rows over which the distributed
94* matrix is distributed.
95*
96* NPCOL (global input) INTEGER
97* The total number of process columns over which the
98* distributed matrix is distributed.
99*
100* MYROW (local input) INTEGER
101* The row coordinate of the process calling this routine.
102*
103* MYCOL (local input) INTEGER
104* The column coordinate of the process calling this routine.
105*
106* LRINDX (local output) INTEGER
107* The local rows starting index of the submatrix.
108*
109* LCINDX (local output) INTEGER
110* The local columns starting index of the submatrix.
111*
112* RSRC (global output) INTEGER
113* The row coordinate of the process that possesses the first
114* row and column of the submatrix.
115*
116* CSRC (global output) INTEGER
117* The column coordinate of the process that possesses the
118* first row and column of the submatrix.
119*
120* =====================================================================
121*
122* .. Parameters ..
123 INTEGER BLOCK_CYCLIC_2D, CSRC_, CTXT_, DLEN_, DTYPE_,
124 $ lld_, mb_, m_, nb_, n_, rsrc_
125 parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
126 $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
127 $ rsrc_ = 7, csrc_ = 8, lld_ = 9 )
128* ..
129* .. Local Scalars ..
130 INTEGER CBLK, GCCPY, GRCPY, RBLK
131* ..
132* .. Intrinsic Functions ..
133 INTRINSIC mod
134* ..
135* .. Executable Statements ..
136*
137 grcpy = grindx-1
138 gccpy = gcindx-1
139*
140 rblk = grcpy / desc(mb_)
141 cblk = gccpy / desc(nb_)
142 rsrc = mod( rblk + desc(rsrc_), nprow )
143 csrc = mod( cblk + desc(csrc_), npcol )
144*
145 lrindx = ( rblk / nprow + 1 ) * desc(mb_) + 1
146 lcindx = ( cblk / npcol + 1 ) * desc(nb_) + 1
147*
148 IF( mod( myrow+nprow-desc(rsrc_), nprow ) .GE.
149 $ mod( rblk, nprow ) ) THEN
150 IF( myrow.EQ.rsrc )
151 $ lrindx = lrindx + mod( grcpy, desc(mb_) )
152 lrindx = lrindx - desc(mb_)
153 END IF
154*
155 IF( mod( mycol+npcol-desc(csrc_), npcol ) .GE.
156 $ mod( cblk, npcol ) ) THEN
157 IF( mycol.EQ.csrc )
158 $ lcindx = lcindx + mod( gccpy, desc(nb_) )
159 lcindx = lcindx - desc(nb_)
160 END IF
161*
162 RETURN
163*
164* End of INFOG2L
165*
166 END
subroutine infog2l(grindx, gcindx, desc, nprow, npcol, myrow, mycol, lrindx, lcindx, rsrc, csrc)
Definition infog2l.f:3