ScaLAPACK 2.1  2.1
ScaLAPACK: Scalable Linear Algebra PACKage
psormr2.f
Go to the documentation of this file.
1  SUBROUTINE psormr2( SIDE, TRANS, M, N, K, A, IA, JA, DESCA, TAU,
2  $ C, IC, JC, DESCC, WORK, LWORK, INFO )
3 *
4 * -- ScaLAPACK routine (version 1.7) --
5 * University of Tennessee, Knoxville, Oak Ridge National Laboratory,
6 * and University of California, Berkeley.
7 * May 25, 2001
8 *
9 * .. Scalar Arguments ..
10  CHARACTER SIDE, TRANS
11  INTEGER IA, IC, INFO, JA, JC, K, LWORK, M, N
12 * ..
13 * .. Array Arguments ..
14  INTEGER DESCA( * ), DESCC( * )
15  REAL A( * ), C( * ), TAU( * ), WORK( * )
16 * ..
17 *
18 * Purpose
19 * =======
20 *
21 * PSORMR2 overwrites the general real M-by-N distributed matrix
22 * sub( C ) = C(IC:IC+M-1,JC:JC+N-1) with
23 *
24 * SIDE = 'L' SIDE = 'R'
25 * TRANS = 'N': Q * sub( C ) sub( C ) * Q
26 * TRANS = 'T': Q**T * sub( C ) sub( C ) * Q**T
27 *
28 * where Q is a real orthogonal distributed matrix defined as the
29 * product of K elementary reflectors
30 *
31 * Q = H(1) H(2) . . . H(k)
32 *
33 * as returned by PSGERQF. Q is of order M if SIDE = 'L' and of order N
34 * if SIDE = 'R'.
35 *
36 * Notes
37 * =====
38 *
39 * Each global data object is described by an associated description
40 * vector. This vector stores the information required to establish
41 * the mapping between an object element and its corresponding process
42 * and memory location.
43 *
44 * Let A be a generic term for any 2D block cyclicly distributed array.
45 * Such a global array has an associated description vector DESCA.
46 * In the following comments, the character _ should be read as
47 * "of the global array".
48 *
49 * NOTATION STORED IN EXPLANATION
50 * --------------- -------------- --------------------------------------
51 * DTYPE_A(global) DESCA( DTYPE_ )The descriptor type. In this case,
52 * DTYPE_A = 1.
53 * CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
54 * the BLACS process grid A is distribu-
55 * ted over. The context itself is glo-
56 * bal, but the handle (the integer
57 * value) may vary.
58 * M_A (global) DESCA( M_ ) The number of rows in the global
59 * array A.
60 * N_A (global) DESCA( N_ ) The number of columns in the global
61 * array A.
62 * MB_A (global) DESCA( MB_ ) The blocking factor used to distribute
63 * the rows of the array.
64 * NB_A (global) DESCA( NB_ ) The blocking factor used to distribute
65 * the columns of the array.
66 * RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
67 * row of the array A is distributed.
68 * CSRC_A (global) DESCA( CSRC_ ) The process column over which the
69 * first column of the array A is
70 * distributed.
71 * LLD_A (local) DESCA( LLD_ ) The leading dimension of the local
72 * array. LLD_A >= MAX(1,LOCr(M_A)).
73 *
74 * Let K be the number of rows or columns of a distributed matrix,
75 * and assume that its process grid has dimension p x q.
76 * LOCr( K ) denotes the number of elements of K that a process
77 * would receive if K were distributed over the p processes of its
78 * process column.
79 * Similarly, LOCc( K ) denotes the number of elements of K that a
80 * process would receive if K were distributed over the q processes of
81 * its process row.
82 * The values of LOCr() and LOCc() may be determined via a call to the
83 * ScaLAPACK tool function, NUMROC:
84 * LOCr( M ) = NUMROC( M, MB_A, MYROW, RSRC_A, NPROW ),
85 * LOCc( N ) = NUMROC( N, NB_A, MYCOL, CSRC_A, NPCOL ).
86 * An upper bound for these quantities may be computed by:
87 * LOCr( M ) <= ceil( ceil(M/MB_A)/NPROW )*MB_A
88 * LOCc( N ) <= ceil( ceil(N/NB_A)/NPCOL )*NB_A
89 *
90 * Arguments
91 * =========
92 *
93 * SIDE (global input) CHARACTER
94 * = 'L': apply Q or Q**T from the Left;
95 * = 'R': apply Q or Q**T from the Right.
96 *
97 * TRANS (global input) CHARACTER
98 * = 'N': No transpose, apply Q;
99 * = 'T': Transpose, apply Q**T.
100 *
101 * M (global input) INTEGER
102 * The number of rows to be operated on i.e the number of rows
103 * of the distributed submatrix sub( C ). M >= 0.
104 *
105 * N (global input) INTEGER
106 * The number of columns to be operated on i.e the number of
107 * columns of the distributed submatrix sub( C ). N >= 0.
108 *
109 * K (global input) INTEGER
110 * The number of elementary reflectors whose product defines the
111 * matrix Q. If SIDE = 'L', M >= K >= 0, if SIDE = 'R',
112 * N >= K >= 0.
113 *
114 * A (local input) REAL pointer into the local memory
115 * to an array of dimension (LLD_A,LOCc(JA+M-1)) if SIDE='L',
116 * and (LLD_A,LOCc(JA+N-1)) if SIDE='R', where
117 * LLD_A >= MAX(1,LOCr(IA+K-1)); On entry, the i-th row must
118 * contain the vector which defines the elementary reflector
119 * H(i), IA <= i <= IA+K-1, as returned by PSGERQF in the
120 * K rows of its distributed matrix argument A(IA:IA+K-1,JA:*).
121 * A(IA:IA+K-1,JA:*) is modified by the routine but restored on
122 * exit.
123 *
124 * IA (global input) INTEGER
125 * The row index in the global array A indicating the first
126 * row of sub( A ).
127 *
128 * JA (global input) INTEGER
129 * The column index in the global array A indicating the
130 * first column of sub( A ).
131 *
132 * DESCA (global and local input) INTEGER array of dimension DLEN_.
133 * The array descriptor for the distributed matrix A.
134 *
135 * TAU (local input) REAL, array, dimension LOCc(IA+K-1).
136 * This array contains the scalar factors TAU(i) of the
137 * elementary reflectors H(i) as returned by PSGERQF.
138 * TAU is tied to the distributed matrix A.
139 *
140 * C (local input/local output) REAL pointer into the
141 * local memory to an array of dimension (LLD_C,LOCc(JC+N-1)).
142 * On entry, the local pieces of the distributed matrix sub(C).
143 * On exit, sub( C ) is overwritten by Q*sub( C ) or Q'*sub( C )
144 * or sub( C )*Q' or sub( C )*Q.
145 *
146 * IC (global input) INTEGER
147 * The row index in the global array C indicating the first
148 * row of sub( C ).
149 *
150 * JC (global input) INTEGER
151 * The column index in the global array C indicating the
152 * first column of sub( C ).
153 *
154 * DESCC (global and local input) INTEGER array of dimension DLEN_.
155 * The array descriptor for the distributed matrix C.
156 *
157 * WORK (local workspace/local output) REAL array,
158 * dimension (LWORK)
159 * On exit, WORK(1) returns the minimal and optimal LWORK.
160 *
161 * LWORK (local or global input) INTEGER
162 * The dimension of the array WORK.
163 * LWORK is local input and must be at least
164 * If SIDE = 'L', LWORK >= MpC0 + MAX( MAX( 1, NqC0 ), NUMROC(
165 * NUMROC( M+IROFFC,MB_A,0,0,NPROW ),MB_A,0,0,LCMP ) );
166 * if SIDE = 'R', LWORK >= NqC0 + MAX( 1, MpC0 );
167 *
168 * where LCMP = LCM / NPROW with LCM = ICLM( NPROW, NPCOL ),
169 *
170 * IROFFC = MOD( IC-1, MB_C ), ICOFFC = MOD( JC-1, NB_C ),
171 * ICROW = INDXG2P( IC, MB_C, MYROW, RSRC_C, NPROW ),
172 * ICCOL = INDXG2P( JC, NB_C, MYCOL, CSRC_C, NPCOL ),
173 * MpC0 = NUMROC( M+IROFFC, MB_C, MYROW, ICROW, NPROW ),
174 * NqC0 = NUMROC( N+ICOFFC, NB_C, MYCOL, ICCOL, NPCOL ),
175 *
176 * ILCM, INDXG2P and NUMROC are ScaLAPACK tool functions;
177 * MYROW, MYCOL, NPROW and NPCOL can be determined by calling
178 * the subroutine BLACS_GRIDINFO.
179 *
180 * If LWORK = -1, then LWORK is global input and a workspace
181 * query is assumed; the routine only calculates the minimum
182 * and optimal size for all work arrays. Each of these
183 * values is returned in the first entry of the corresponding
184 * work array, and no error message is issued by PXERBLA.
185 *
186 *
187 * INFO (local output) INTEGER
188 * = 0: successful exit
189 * < 0: If the i-th argument is an array and the j-entry had
190 * an illegal value, then INFO = -(i*100+j), if the i-th
191 * argument is a scalar and had an illegal value, then
192 * INFO = -i.
193 *
194 * Alignment requirements
195 * ======================
196 *
197 * The distributed submatrices A(IA:*, JA:*) and C(IC:IC+M-1,JC:JC+N-1)
198 * must verify some alignment properties, namely the following
199 * expressions should be true:
200 *
201 * If SIDE = 'L',
202 * ( NB_A.EQ.MB_C .AND. ICOFFA.EQ.IROFFC )
203 * If SIDE = 'R',
204 * ( NB_A.EQ.NB_C .AND. ICOFFA.EQ.ICOFFC .AND. IACOL.EQ.ICCOL )
205 *
206 * =====================================================================
207 *
208 * .. Parameters ..
209  INTEGER BLOCK_CYCLIC_2D, CSRC_, CTXT_, DLEN_, DTYPE_,
210  $ lld_, mb_, m_, nb_, n_, rsrc_
211  parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
212  $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
213  $ rsrc_ = 7, csrc_ = 8, lld_ = 9 )
214  REAL ONE
215  parameter( one = 1.0e+0 )
216 * ..
217 * .. Local Scalars ..
218  LOGICAL LEFT, LQUERY, NOTRAN
219  CHARACTER COLBTOP, ROWBTOP
220  INTEGER I, I1, I2, I3, IACOL, ICCOL, ICOFFA, ICOFFC,
221  $ icrow, ictxt, iroffc, lcm, lcmp, lwmin, mi,
222  $ mpc0, mycol, myrow, ni, npcol, nprow, nq, nqc0
223  REAL AII
224 * ..
225 * .. External Subroutines ..
226  EXTERNAL blacs_abort, blacs_gridinfo, chk1mat, pselset,
227  $ pselset2, pslarf, pb_topget, pb_topset, pxerbla
228 * ..
229 * .. External Functions ..
230  LOGICAL LSAME
231  INTEGER ILCM, INDXG2P, NUMROC
232  EXTERNAL ilcm, indxg2p, lsame, numroc
233 * ..
234 * .. Intrinsic Functions ..
235  INTRINSIC max, mod, real
236 * ..
237 * .. Executable Statements ..
238 *
239 * Get grid parameters
240 *
241  ictxt = desca( ctxt_ )
242  CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
243 *
244 * Test the input parameters
245 *
246  info = 0
247  IF( nprow.EQ.-1 ) THEN
248  info = -(900+ctxt_)
249  ELSE
250  left = lsame( side, 'L' )
251  notran = lsame( trans, 'N' )
252 *
253 * NQ is the order of Q
254 *
255  IF( left ) THEN
256  nq = m
257  CALL chk1mat( k, 5, m, 3, ia, ja, desca, 9, info )
258  ELSE
259  nq = n
260  CALL chk1mat( k, 5, n, 4, ia, ja, desca, 9, info )
261  END IF
262  CALL chk1mat( m, 3, n, 4, ic, jc, descc, 14, info )
263  IF( info.EQ.0 ) THEN
264  icoffa = mod( ja-1, desca( nb_ ) )
265  iroffc = mod( ic-1, descc( mb_ ) )
266  icoffc = mod( jc-1, descc( nb_ ) )
267  iacol = indxg2p( ja, desca( nb_ ), mycol, desca( csrc_ ),
268  $ npcol )
269  icrow = indxg2p( ic, descc( mb_ ), myrow, descc( rsrc_ ),
270  $ nprow )
271  iccol = indxg2p( jc, descc( nb_ ), mycol, descc( csrc_ ),
272  $ npcol )
273  mpc0 = numroc( m+iroffc, descc( mb_ ), myrow, icrow, nprow )
274  nqc0 = numroc( n+icoffc, descc( nb_ ), mycol, iccol, npcol )
275 *
276  IF( left ) THEN
277  lcm = ilcm( nprow, npcol )
278  lcmp = lcm / nprow
279  lwmin = mpc0 + max( max( 1, nqc0 ), numroc( numroc(
280  $ m+iroffc, desca( mb_ ), 0, 0, nprow ),
281  $ desca( mb_ ), 0, 0, lcmp ) )
282  ELSE
283  lwmin = nqc0 + max( 1, mpc0 )
284  END IF
285 *
286  work( 1 ) = real( lwmin )
287  lquery = ( lwork.EQ.-1 )
288  IF( .NOT.left .AND. .NOT.lsame( side, 'R' ) ) THEN
289  info = -1
290  ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) ) THEN
291  info = -2
292  ELSE IF( k.LT.0 .OR. k.GT.nq ) THEN
293  info = -5
294  ELSE IF( left .AND. desca( nb_ ).NE.descc( mb_ ) ) THEN
295  info = -(900+nb_)
296  ELSE IF( left .AND. icoffa.NE.iroffc ) THEN
297  info = -12
298  ELSE IF( .NOT.left .AND. icoffa.NE.icoffc ) THEN
299  info = -13
300  ELSE IF( .NOT.left .AND. iacol.NE.iccol ) THEN
301  info = -13
302  ELSE IF( .NOT.left .AND. desca( nb_ ).NE.descc( nb_ ) ) THEN
303  info = -(1400+nb_)
304  ELSE IF( ictxt.NE.descc( ctxt_ ) ) THEN
305  info = -(1400+ctxt_)
306  ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
307  info = -16
308  END IF
309  END IF
310  END IF
311 *
312  IF( info.NE.0 ) THEN
313  CALL pxerbla( ictxt, 'PSORMR2', -info )
314  CALL blacs_abort( ictxt, 1 )
315  RETURN
316  ELSE IF( lquery ) THEN
317  RETURN
318  END IF
319 *
320 * Quick return if possible
321 *
322  IF( m.EQ.0 .OR. n.EQ.0 .OR. k.EQ.0 )
323  $ RETURN
324 *
325  CALL pb_topget( ictxt, 'Broadcast', 'Rowwise', rowbtop )
326  CALL pb_topget( ictxt, 'Broadcast', 'Columnwise', colbtop )
327 *
328  IF( ( left .AND. .NOT.notran .OR. .NOT.left .AND. notran ) ) THEN
329  i1 = ia
330  i2 = ia + k - 1
331  i3 = 1
332  ELSE
333  i1 = ia + k - 1
334  i2 = ia
335  i3 = -1
336  END IF
337 *
338  IF( left ) THEN
339  ni = n
340  ELSE
341  mi = m
342  CALL pb_topset( ictxt, 'Broadcast', 'Rowwise', ' ' )
343  IF( notran ) THEN
344  CALL pb_topset( ictxt, 'Broadcast', 'Columnwise', 'I-ring' )
345  ELSE
346  CALL pb_topset( ictxt, 'Broadcast', 'Columnwise', 'D-ring' )
347  END IF
348  END IF
349 *
350  DO 10 i = i1, i2, i3
351  IF( left ) THEN
352 *
353 * H(i) or H(i)' is applied to C(ic:ic+m-k+i-ia,jc:jc+n-1)
354 *
355  mi = m - k + i - ia + 1
356  ELSE
357 *
358 * H(i) or H(i)' is applied to C(ic:ic+m-1,jc:jc+n-k+i-ia+1)
359 *
360  ni = n - k + i - ia + 1
361  END IF
362 *
363 * Apply H(i) or H(i)'
364 *
365  CALL pselset2( aii, a, i, ja+nq-k+i-ia, desca, one )
366  CALL pslarf( side, mi, ni, a, i, ja, desca, desca( m_ ),
367  $ tau, c, ic, jc, descc, work )
368  CALL pselset( a, i, ja+nq-k+i-ia, desca, aii )
369 *
370  10 CONTINUE
371 *
372  CALL pb_topset( ictxt, 'Broadcast', 'Rowwise', rowbtop )
373  CALL pb_topset( ictxt, 'Broadcast', 'Columnwise', colbtop )
374 *
375  work( 1 ) = real( lwmin )
376 *
377  RETURN
378 *
379 * End of PSORMR2
380 *
381  END
max
#define max(A, B)
Definition: pcgemr.c:180
psormr2
subroutine psormr2(SIDE, TRANS, M, N, K, A, IA, JA, DESCA, TAU, C, IC, JC, DESCC, WORK, LWORK, INFO)
Definition: psormr2.f:3
pselset
subroutine pselset(A, IA, JA, DESCA, ALPHA)
Definition: pselset.f:2
pselset2
subroutine pselset2(ALPHA, A, IA, JA, DESCA, BETA)
Definition: pselset2.f:2
pslarf
subroutine pslarf(SIDE, M, N, V, IV, JV, DESCV, INCV, TAU, C, IC, JC, DESCC, WORK)
Definition: pslarf.f:3
chk1mat
subroutine chk1mat(MA, MAPOS0, NA, NAPOS0, IA, JA, DESCA, DESCAPOS0, INFO)
Definition: chk1mat.f:3
pxerbla
subroutine pxerbla(ICTXT, SRNAME, INFO)
Definition: pxerbla.f:2