ScaLAPACK 2.1  2.1
ScaLAPACK: Scalable Linear Algebra PACKage
pdpotf2.f
Go to the documentation of this file.
1  SUBROUTINE pdpotf2( UPLO, N, A, IA, JA, DESCA, INFO )
2 *
3 * -- ScaLAPACK routine (version 1.7) --
4 * University of Tennessee, Knoxville, Oak Ridge National Laboratory,
5 * and University of California, Berkeley.
6 * May 1, 1997
7 *
8 * .. Scalar Arguments ..
9  CHARACTER UPLO
10  INTEGER IA, INFO, JA, N
11 * ..
12 * .. Array Arguments ..
13  INTEGER DESCA( * )
14  DOUBLE PRECISION A( * )
15 * ..
16 *
17 * Purpose
18 * =======
19 *
20 * PDPOTF2 computes the Cholesky factorization of a real symmetric
21 * positive definite distributed matrix sub( A )=A(IA:IA+N-1,JA:JA+N-1).
22 *
23 * The factorization has the form
24 *
25 * sub( A ) = U' * U , if UPLO = 'U', or
26 *
27 * sub( A ) = L * L', if UPLO = 'L',
28 *
29 * where U is an upper triangular matrix and L is lower triangular.
30 *
31 * Notes
32 * =====
33 *
34 * Each global data object is described by an associated description
35 * vector. This vector stores the information required to establish
36 * the mapping between an object element and its corresponding process
37 * and memory location.
38 *
39 * Let A be a generic term for any 2D block cyclicly distributed array.
40 * Such a global array has an associated description vector DESCA.
41 * In the following comments, the character _ should be read as
42 * "of the global array".
43 *
44 * NOTATION STORED IN EXPLANATION
45 * --------------- -------------- --------------------------------------
46 * DTYPE_A(global) DESCA( DTYPE_ )The descriptor type. In this case,
47 * DTYPE_A = 1.
48 * CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
49 * the BLACS process grid A is distribu-
50 * ted over. The context itself is glo-
51 * bal, but the handle (the integer
52 * value) may vary.
53 * M_A (global) DESCA( M_ ) The number of rows in the global
54 * array A.
55 * N_A (global) DESCA( N_ ) The number of columns in the global
56 * array A.
57 * MB_A (global) DESCA( MB_ ) The blocking factor used to distribute
58 * the rows of the array.
59 * NB_A (global) DESCA( NB_ ) The blocking factor used to distribute
60 * the columns of the array.
61 * RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
62 * row of the array A is distributed.
63 * CSRC_A (global) DESCA( CSRC_ ) The process column over which the
64 * first column of the array A is
65 * distributed.
66 * LLD_A (local) DESCA( LLD_ ) The leading dimension of the local
67 * array. LLD_A >= MAX(1,LOCr(M_A)).
68 *
69 * Let K be the number of rows or columns of a distributed matrix,
70 * and assume that its process grid has dimension p x q.
71 * LOCr( K ) denotes the number of elements of K that a process
72 * would receive if K were distributed over the p processes of its
73 * process column.
74 * Similarly, LOCc( K ) denotes the number of elements of K that a
75 * process would receive if K were distributed over the q processes of
76 * its process row.
77 * The values of LOCr() and LOCc() may be determined via a call to the
78 * ScaLAPACK tool function, NUMROC:
79 * LOCr( M ) = NUMROC( M, MB_A, MYROW, RSRC_A, NPROW ),
80 * LOCc( N ) = NUMROC( N, NB_A, MYCOL, CSRC_A, NPCOL ).
81 * An upper bound for these quantities may be computed by:
82 * LOCr( M ) <= ceil( ceil(M/MB_A)/NPROW )*MB_A
83 * LOCc( N ) <= ceil( ceil(N/NB_A)/NPCOL )*NB_A
84 *
85 * This routine requires N <= NB_A-MOD(JA-1, NB_A) and square block
86 * decomposition ( MB_A = NB_A ).
87 *
88 * Arguments
89 * =========
90 *
91 * UPLO (global input) CHARACTER
92 * = 'U': Upper triangle of sub( A ) is stored;
93 * = 'L': Lower triangle of sub( A ) is stored.
94 *
95 * N (global input) INTEGER
96 * The number of rows and columns to be operated on, i.e. the
97 * order of the distributed submatrix sub( A ). N >= 0.
98 *
99 * A (local input/local output) DOUBLE PRECISION pointer into the
100 * local memory to an array of dimension (LLD_A, LOCc(JA+N-1)).
101 * On entry, this array contains the local pieces of the
102 * N-by-N symmetric distributed matrix sub( A ) to be factored.
103 * If UPLO = 'U', the leading N-by-N upper triangular part of
104 * sub( A ) contains the upper triangular part of the matrix,
105 * and its strictly lower triangular part is not referenced.
106 * If UPLO = 'L', the leading N-by-N lower triangular part of
107 * sub( A ) contains the lower triangular part of the distribu-
108 * ted matrix, and its strictly upper triangular part is not
109 * referenced. On exit, if UPLO = 'U', the upper triangular
110 * part of the distributed matrix contains the Cholesky factor
111 * U, if UPLO = 'L', the lower triangular part of the distribu-
112 * ted matrix contains the Cholesky factor L.
113 *
114 * IA (global input) INTEGER
115 * The row index in the global array A indicating the first
116 * row of sub( A ).
117 *
118 * JA (global input) INTEGER
119 * The column index in the global array A indicating the
120 * first column of sub( A ).
121 *
122 * DESCA (global and local input) INTEGER array of dimension DLEN_.
123 * The array descriptor for the distributed matrix A.
124 *
125 * INFO (local output) INTEGER
126 * = 0: successful exit
127 * < 0: If the i-th argument is an array and the j-entry had
128 * an illegal value, then INFO = -(i*100+j), if the i-th
129 * argument is a scalar and had an illegal value, then
130 * INFO = -i.
131 * > 0: If INFO = K, the leading minor of order K,
132 * A(IA:IA+K-1,JA:JA+K-1) is not positive definite, and
133 * the factorization could not be completed.
134 *
135 * =====================================================================
136 *
137 * .. Parameters ..
138  INTEGER BLOCK_CYCLIC_2D, CSRC_, CTXT_, DLEN_, DTYPE_,
139  $ LLD_, MB_, M_, NB_, N_, RSRC_
140  parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
141  $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
142  $ rsrc_ = 7, csrc_ = 8, lld_ = 9 )
143  DOUBLE PRECISION ONE, ZERO
144  parameter( one = 1.0d+0, zero = 0.0d+0 )
145 * ..
146 * .. Local Scalars ..
147  LOGICAL UPPER
148  CHARACTER COLBTOP, ROWBTOP
149  INTEGER IACOL, IAROW, ICOFF, ICTXT, ICURR, IDIAG, IIA,
150  $ IOFFA, IROFF, J, JJA, LDA, MYCOL, MYROW,
151  $ NPCOL, NPROW
152  DOUBLE PRECISION AJJ
153 * ..
154 * .. External Subroutines ..
155  EXTERNAL blacs_abort, blacs_gridinfo, chk1mat, dgemv,
156  $ dscal, igebr2d, igebs2d, infog2l, pb_topget,
157  $ pxerbla
158 * ..
159 * .. Intrinsic Functions ..
160  INTRINSIC mod, sqrt
161 * ..
162 * .. External Functions ..
163  LOGICAL LSAME
164  DOUBLE PRECISION DDOT
165  EXTERNAL lsame, ddot
166 * ..
167 * .. Executable Statements ..
168 *
169 * Get grid parameters
170 *
171  ictxt = desca( ctxt_ )
172  CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
173 *
174 * Test the input parameters.
175 *
176  info = 0
177  IF( nprow.EQ.-1 ) THEN
178  info = -(600+ctxt_)
179  ELSE
180  CALL chk1mat( n, 2, n, 2, ia, ja, desca, 6, info )
181  IF( info.EQ.0 ) THEN
182  upper = lsame( uplo, 'U' )
183  iroff = mod( ia-1, desca( mb_ ) )
184  icoff = mod( ja-1, desca( nb_ ) )
185  IF ( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
186  info = -1
187  ELSE IF( n+icoff.GT.desca( nb_ ) ) THEN
188  info = -2
189  ELSE IF( iroff.NE.0 ) THEN
190  info = -4
191  ELSE IF( icoff.NE.0 ) THEN
192  info = -5
193  ELSE IF( desca( mb_ ).NE.desca( nb_ ) ) THEN
194  info = -(600+nb_)
195  END IF
196  END IF
197  END IF
198 *
199  IF( info.NE.0 ) THEN
200  CALL pxerbla( ictxt, 'PDPOTF2', -info )
201  CALL blacs_abort( ictxt, 1 )
202  RETURN
203  END IF
204 *
205 * Quick return if possible
206 *
207  IF( n.EQ.0 )
208  $ RETURN
209 *
210 * Compute local information
211 *
212  CALL infog2l( ia, ja, desca, nprow, npcol, myrow, mycol, iia, jja,
213  $ iarow, iacol )
214  CALL pb_topget( ictxt, 'Broadcast', 'Rowwise', rowbtop )
215  CALL pb_topget( ictxt, 'Broadcast', 'Columnwise', colbtop )
216 *
217  IF ( upper ) THEN
218 *
219 * Process (IAROW, IACOL) owns block to be factorized
220 *
221  IF( myrow.EQ.iarow ) THEN
222  IF( mycol.EQ.iacol ) THEN
223 *
224 * Compute the Cholesky factorization A = U'*U.
225 *
226  lda = desca( lld_ )
227  idiag = iia + ( jja - 1 ) * lda
228  ioffa = idiag
229 *
230  DO 10 j = ja, ja+n-1
231 *
232 * Compute U(J,J) and test for non-positive-definiteness.
233 *
234  ajj = a( idiag ) -
235  $ ddot( j-ja, a( ioffa ), 1, a( ioffa ), 1 )
236  IF( ajj.LE.zero ) THEN
237  a( idiag ) = ajj
238  info = j - ja + 1
239  GO TO 20
240  END IF
241  ajj = sqrt( ajj )
242  a( idiag ) = ajj
243 *
244 * Compute elements J+1:JA+N-1 of row J.
245 *
246  IF( j.LT.ja+n-1 ) THEN
247  icurr = idiag + lda
248  CALL dgemv( 'Transpose', j-ja, ja+n-j-1, -one,
249  $ a( ioffa+lda ), lda, a( ioffa ), 1,
250  $ one, a( icurr ), lda )
251  CALL dscal( n-j+ja-1, one / ajj, a( icurr ), lda )
252  END IF
253  idiag = idiag + lda + 1
254  ioffa = ioffa + lda
255  10 CONTINUE
256 *
257  20 CONTINUE
258 *
259 * Broadcast INFO to all processes in my IAROW.
260 *
261  CALL igebs2d( ictxt, 'Rowwise', rowbtop, 1, 1, info, 1 )
262 *
263  ELSE
264 *
265  CALL igebr2d( ictxt, 'Rowwise', rowbtop, 1, 1, info, 1,
266  $ myrow, iacol )
267  END IF
268 *
269 * IAROW bcasts along columns so that everyone has INFO
270 *
271  CALL igebs2d( ictxt, 'Columnwise', colbtop, 1, 1, info, 1 )
272 *
273  ELSE
274 *
275  CALL igebr2d( ictxt, 'Columnwise', colbtop, 1, 1, info, 1,
276  $ iarow, mycol )
277 *
278  END IF
279 *
280  ELSE
281 *
282 * Process (IAROW, IACOL) owns block to be factorized
283 *
284  IF( mycol.EQ.iacol ) THEN
285  IF( myrow.EQ.iarow ) THEN
286 *
287 * Compute the Cholesky factorization A = L*L'.
288 *
289  lda = desca( lld_ )
290  idiag = iia + ( jja - 1 ) * lda
291  ioffa = idiag
292 *
293  DO 30 j = ja, ja+n-1
294 *
295 * Compute L(J,J) and test for non-positive-definiteness.
296 *
297  ajj = a( idiag ) -
298  $ ddot( j-ja, a( ioffa ), lda, a( ioffa ), lda )
299  IF ( ajj.LE.zero ) THEN
300  a( idiag ) = ajj
301  info = j - ja + 1
302  GO TO 40
303  END IF
304  ajj = sqrt( ajj )
305  a( idiag ) = ajj
306 *
307 * Compute elements J+1:JA+N-1 of column J.
308 *
309  IF( j.LT.ja+n-1 ) THEN
310  icurr = idiag + 1
311  CALL dgemv( 'No transpose', ja+n-j-1, j-ja, -one,
312  $ a( ioffa+1 ), lda, a( ioffa ), lda,
313  $ one, a( icurr ), 1 )
314  CALL dscal( ja+n-j-1, one / ajj, a( icurr ), 1 )
315  END IF
316  idiag = idiag + lda + 1
317  ioffa = ioffa + 1
318  30 CONTINUE
319 *
320  40 CONTINUE
321 *
322 * Broadcast INFO to everyone in IACOL
323 *
324  CALL igebs2d( ictxt, 'Columnwise', colbtop, 1, 1, info,
325  $ 1 )
326 *
327  ELSE
328 *
329  CALL igebr2d( ictxt, 'Columnwise', colbtop, 1, 1, info,
330  $ 1, iarow, mycol )
331 *
332  END IF
333 *
334 * IACOL bcasts INFO along rows so that everyone has it
335 *
336  CALL igebs2d( ictxt, 'Rowwise', rowbtop, 1, 1, info, 1 )
337 *
338  ELSE
339 *
340  CALL igebr2d( ictxt, 'Rowwise', rowbtop, 1, 1, info, 1,
341  $ myrow, iacol )
342 *
343  END IF
344 *
345  END IF
346 *
347  RETURN
348 *
349 * End of PDPOTF2
350 *
351  END
infog2l
subroutine infog2l(GRINDX, GCINDX, DESC, NPROW, NPCOL, MYROW, MYCOL, LRINDX, LCINDX, RSRC, CSRC)
Definition: infog2l.f:3
pdpotf2
subroutine pdpotf2(UPLO, N, A, IA, JA, DESCA, INFO)
Definition: pdpotf2.f:2
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