ScaLAPACK 2.1  2.1
ScaLAPACK: Scalable Linear Algebra PACKage
pdorm2r.f
Go to the documentation of this file.
1  SUBROUTINE pdorm2r( 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  DOUBLE PRECISION A( * ), C( * ), TAU( * ), WORK( * )
16 * ..
17 *
18 * Purpose
19 * =======
20 *
21 * PDORM2R 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 PDGEQRF. 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) DOUBLE PRECISION pointer into the local memory
115 * to an array of dimension (LLD_A,LOCc(JA+K-1)). On entry, the
116 * j-th column must contain the vector which defines the elemen-
117 * tary reflector H(j), JA <= j <= JA+K-1, as returned by
118 * PDGEQRF in the K columns of its distributed matrix
119 * argument A(IA:*,JA:JA+K-1). A(IA:*,JA:JA+K-1) is modified by
120 * the routine but restored on exit.
121 * If SIDE = 'L', LLD_A >= MAX( 1, LOCr(IA+M-1) );
122 * if SIDE = 'R', LLD_A >= MAX( 1, LOCr(IA+N-1) ).
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) DOUBLE PRECISION array, dimension LOCc(JA+K-1).
136 * This array contains the scalar factors TAU(j) of the
137 * elementary reflectors H(j) as returned by PDGEQRF.
138 * TAU is tied to the distributed matrix A.
139 *
140 * C (local input/local output) DOUBLE PRECISION 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) DOUBLE PRECISION 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( 1, NqC0 );
165 * if SIDE = 'R', LWORK >= NqC0 + MAX( MAX( 1, MpC0 ), NUMROC(
166 * NUMROC( N+ICOFFC,NB_A,0,0,NPCOL ),NB_A,0,0,LCMQ ) );
167 *
168 * where LCMQ = LCM / NPCOL 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 * ( MB_A.EQ.MB_C .AND. IROFFA.EQ.IROFFC .AND. IAROW.EQ.ICROW )
203 * If SIDE = 'R',
204 * ( MB_A.EQ.NB_C .AND. IROFFA.EQ.ICOFFC )
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  DOUBLE PRECISION ONE
215  parameter( one = 1.0d+0 )
216 * ..
217 * .. Local Scalars ..
218  LOGICAL LEFT, LQUERY, NOTRAN
219  CHARACTER COLBTOP, ROWBTOP
220  INTEGER IACOL, IAROW, ICCOL, ICOFFC, ICROW, ICTXT, ICC,
221  $ ii, iroffa, iroffc, j, j1, j2, j3, jcc, jj,
222  $ lcm, lcmq, lwmin, mi, mp, mpc0, mycol, myrow,
223  $ ni, npcol, nprow, nq, nqc0
224  DOUBLE PRECISION AJJ
225 * ..
226 * .. External Subroutines ..
227  EXTERNAL blacs_abort, blacs_gridinfo, chk1mat, dgebr2d,
228  $ dgebs2d, dgerv2d, dgesd2d, dscal,
230  $ pb_topget, pb_topset, pxerbla
231 * ..
232 * .. External Functions ..
233  LOGICAL LSAME
234  INTEGER ILCM, INDXG2P, NUMROC
235  EXTERNAL ilcm, indxg2p, lsame, numroc
236 * ..
237 * .. Intrinsic Functions ..
238  INTRINSIC dble, max, mod
239 * ..
240 * .. Executable Statements ..
241 *
242 * Get grid parameters
243 *
244  ictxt = desca( ctxt_ )
245  CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
246 *
247 * Test the input parameters
248 *
249  info = 0
250  IF( nprow.EQ.-1 ) THEN
251  info = -(900+ctxt_)
252  ELSE
253  left = lsame( side, 'L' )
254  notran = lsame( trans, 'N' )
255 *
256 * NQ is the order of Q
257 *
258  IF( left ) THEN
259  nq = m
260  CALL chk1mat( m, 3, k, 5, ia, ja, desca, 9, info )
261  ELSE
262  nq = n
263  CALL chk1mat( n, 4, k, 5, ia, ja, desca, 9, info )
264  END IF
265  CALL chk1mat( m, 3, n, 4, ic, jc, descc, 14, info )
266  IF( info.EQ.0 ) THEN
267  iroffa = mod( ia-1, desca( mb_ ) )
268  iroffc = mod( ic-1, descc( mb_ ) )
269  icoffc = mod( jc-1, descc( nb_ ) )
270  iarow = indxg2p( ia, desca( mb_ ), myrow, desca( rsrc_ ),
271  $ nprow )
272  icrow = indxg2p( ic, descc( mb_ ), myrow, descc( rsrc_ ),
273  $ nprow )
274  iccol = indxg2p( jc, descc( nb_ ), mycol, descc( csrc_ ),
275  $ npcol )
276  mpc0 = numroc( m+iroffc, descc( mb_ ), myrow, icrow, nprow )
277  nqc0 = numroc( n+icoffc, descc( nb_ ), mycol, iccol, npcol )
278 *
279  IF( left ) THEN
280  lwmin = mpc0 + max( 1, nqc0 )
281  ELSE
282  lcm = ilcm( nprow, npcol )
283  lcmq = lcm / npcol
284  lwmin = nqc0 + max( max( 1, mpc0 ), numroc( numroc(
285  $ n+icoffc, desca( nb_ ), 0, 0, npcol ),
286  $ desca( nb_ ), 0, 0, lcmq ) )
287  END IF
288 *
289  work( 1 ) = dble( lwmin )
290  lquery = ( lwork.EQ.-1 )
291  IF( .NOT.left .AND. .NOT.lsame( side, 'R' ) ) THEN
292  info = -1
293  ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) ) THEN
294  info = -2
295  ELSE IF( k.LT.0 .OR. k.GT.nq ) THEN
296  info = -5
297  ELSE IF( .NOT.left .AND. desca( mb_ ).NE.descc( nb_ ) ) THEN
298  info = -(900+nb_)
299  ELSE IF( left .AND. iroffa.NE.iroffc ) THEN
300  info = -12
301  ELSE IF( left .AND. iarow.NE.icrow ) THEN
302  info = -12
303  ELSE IF( .NOT.left .AND. iroffa.NE.icoffc ) THEN
304  info = -13
305  ELSE IF( left .AND. desca( mb_ ).NE.descc( mb_ ) ) THEN
306  info = -(1400+mb_)
307  ELSE IF( ictxt.NE.descc( ctxt_ ) ) THEN
308  info = -(1400+ctxt_)
309  ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
310  info = -16
311  END IF
312  END IF
313  END IF
314 *
315  IF( info.NE.0 ) THEN
316  CALL pxerbla( ictxt, 'PDORM2R', -info )
317  CALL blacs_abort( ictxt, 1 )
318  RETURN
319  ELSE IF( lquery ) THEN
320  RETURN
321  END IF
322 *
323 * Quick return if possible
324 *
325  IF( m.EQ.0 .OR. n.EQ.0 .OR. k.EQ.0 )
326  $ RETURN
327 *
328  IF( desca( m_ ).EQ.1 ) THEN
329  CALL infog2l( ia, ja, desca, nprow, npcol, myrow, mycol, ii,
330  $ jj, iarow, iacol )
331  CALL infog2l( ic, jc, descc, nprow, npcol, myrow, mycol, icc,
332  $ jcc, icrow, iccol )
333  IF( left ) THEN
334  IF( myrow.EQ.iarow ) THEN
335  nq = numroc( jc+n-1, descc( nb_ ), mycol, descc( csrc_ ),
336  $ npcol )
337  IF( mycol.EQ.iacol ) THEN
338  ajj = one - tau( jj )
339  CALL dgebs2d( ictxt, 'Rowwise', ' ', 1, 1, ajj, 1 )
340  CALL dscal( nq-jcc+1, ajj,
341  $ c( icc+(jcc-1)*descc( lld_ ) ),
342  $ descc( lld_ ) )
343  ELSE
344  CALL dgebr2d( ictxt, 'Rowwise', ' ', 1, 1, ajj, 1,
345  $ iarow, iacol )
346  CALL dscal( nq-jcc+1, ajj,
347  $ c( icc+(jcc-1)*descc( lld_ ) ),
348  $ descc( lld_ ) )
349  END IF
350  END IF
351  ELSE
352  IF( mycol.EQ.iacol ) THEN
353  ajj = one - tau( jj )
354  END IF
355 *
356  IF( iacol.NE.iccol ) THEN
357  IF( mycol.EQ.iacol )
358  $ CALL dgesd2d( ictxt, 1, 1, ajj, 1, myrow, iccol )
359  IF( mycol.EQ.iccol )
360  $ CALL dgerv2d( ictxt, 1, 1, ajj, 1, myrow, iacol )
361  END IF
362 *
363  IF( mycol.EQ.iccol ) THEN
364  mp = numroc( ic+m-1, descc( mb_ ), myrow, descc( rsrc_ ),
365  $ nprow )
366  CALL dscal( mp-icc+1, ajj, c( icc+(jcc-1)*
367  $ descc( lld_ ) ), 1 )
368  END IF
369 *
370  END IF
371 *
372  ELSE
373 *
374  CALL pb_topget( ictxt, 'Broadcast', 'Rowwise', rowbtop )
375  CALL pb_topget( ictxt, 'Broadcast', 'Columnwise', colbtop )
376 *
377  IF( left .AND. .NOT.notran .OR. .NOT.left .AND. notran ) THEN
378  j1 = ja
379  j2 = ja+k-1
380  j3 = 1
381  ELSE
382  j1 = ja+k-1
383  j2 = ja
384  j3 = -1
385  END IF
386 *
387  IF( left ) THEN
388  ni = n
389  jcc = jc
390  IF( notran ) THEN
391  CALL pb_topset( ictxt, 'Broadcast', 'Rowwise', 'D-ring' )
392  ELSE
393  CALL pb_topset( ictxt, 'Broadcast', 'Rowwise', 'I-ring' )
394  END IF
395  CALL pb_topset( ictxt, 'Broadcast', 'Columnwise', ' ' )
396  ELSE
397  mi = m
398  icc = ic
399  END IF
400 *
401  DO 10 j = j1, j2, j3
402  IF( left ) THEN
403 *
404 * H(j) or H(j)' is applied to C(ic+j-ja:ic+m-1,jc:jc+n-1)
405 *
406  mi = m - j + ja
407  icc = ic + j - ja
408  ELSE
409 *
410 * H(j) or H(j)' is applied to C(ic:ic+m-1,jc+j-ja:jc+n-1)
411 *
412  ni = n - j + ja
413  jcc = jc + j - ja
414  END IF
415 *
416 * Apply H(j) or H(j)'
417 *
418  CALL pdelset2( ajj, a, ia+j-ja, j, desca, one )
419  CALL pdlarf( side, mi, ni, a, ia+j-ja, j, desca, 1, tau, c,
420  $ icc, jcc, descc, work )
421  CALL pdelset( a, ia+j-ja, j, desca, ajj )
422 *
423  10 CONTINUE
424 *
425  CALL pb_topset( ictxt, 'Broadcast', 'Rowwise', rowbtop )
426  CALL pb_topset( ictxt, 'Broadcast', 'Columnwise', colbtop )
427 *
428  END IF
429 *
430  work( 1 ) = dble( lwmin )
431 *
432  RETURN
433 *
434 * End of PDORM2R
435 *
436  END
pdlarf
subroutine pdlarf(SIDE, M, N, V, IV, JV, DESCV, INCV, TAU, C, IC, JC, DESCC, WORK)
Definition: pdlarf.f:3
max
#define max(A, B)
Definition: pcgemr.c:180
pdorm2r
subroutine pdorm2r(SIDE, TRANS, M, N, K, A, IA, JA, DESCA, TAU, C, IC, JC, DESCC, WORK, LWORK, INFO)
Definition: pdorm2r.f:3
infog2l
subroutine infog2l(GRINDX, GCINDX, DESC, NPROW, NPCOL, MYROW, MYCOL, LRINDX, LCINDX, RSRC, CSRC)
Definition: infog2l.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
pdelset
subroutine pdelset(A, IA, JA, DESCA, ALPHA)
Definition: pdelset.f:2
pdelset2
subroutine pdelset2(ALPHA, A, IA, JA, DESCA, BETA)
Definition: pdelset2.f:2