ScaLAPACK 2.1  2.1
ScaLAPACK: Scalable Linear Algebra PACKage
pdsygst.f
Go to the documentation of this file.
1 *
2 *
3  SUBROUTINE pdsygst( IBTYPE, UPLO, N, A, IA, JA, DESCA, B, IB, JB,
4  $ DESCB, SCALE, INFO )
5 *
6 * -- ScaLAPACK routine (version 1.7) --
7 * University of Tennessee, Knoxville, Oak Ridge National Laboratory,
8 * and University of California, Berkeley.
9 * May 1, 1997
10 *
11 * .. Scalar Arguments ..
12  CHARACTER UPLO
13  INTEGER IA, IB, IBTYPE, INFO, JA, JB, N
14  DOUBLE PRECISION SCALE
15 * ..
16 * .. Array Arguments ..
17  INTEGER DESCA( * ), DESCB( * )
18  DOUBLE PRECISION A( * ), B( * )
19 * ..
20 *
21 * Purpose
22 * =======
23 *
24 * PDSYGST reduces a real symmetric-definite generalized eigenproblem
25 * to standard form.
26 *
27 * In the following sub( A ) denotes A( IA:IA+N-1, JA:JA+N-1 ) and
28 * sub( B ) denotes B( IB:IB+N-1, JB:JB+N-1 ).
29 *
30 * If IBTYPE = 1, the problem is sub( A )*x = lambda*sub( B )*x,
31 * and sub( A ) is overwritten by inv(U**T)*sub( A )*inv(U) or
32 * inv(L)*sub( A )*inv(L**T)
33 *
34 * If IBTYPE = 2 or 3, the problem is sub( A )*sub( B )*x = lambda*x or
35 * sub( B )*sub( A )*x = lambda*x, and sub( A ) is overwritten by
36 * U*sub( A )*U**T or L**T*sub( A )*L.
37 *
38 * sub( B ) must have been previously factorized as U**T*U or L*L**T by
39 * PDPOTRF.
40 *
41 * Notes
42 * =====
43 *
44 * Each global data object is described by an associated description
45 * vector. This vector stores the information required to establish
46 * the mapping between an object element and its corresponding process
47 * and memory location.
48 *
49 * Let A be a generic term for any 2D block cyclicly distributed array.
50 * Such a global array has an associated description vector DESCA.
51 * In the following comments, the character _ should be read as
52 * "of the global array".
53 *
54 * NOTATION STORED IN EXPLANATION
55 * --------------- -------------- --------------------------------------
56 * DTYPE_A(global) DESCA( DTYPE_ )The descriptor type. In this case,
57 * DTYPE_A = 1.
58 * CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
59 * the BLACS process grid A is distribu-
60 * ted over. The context itself is glo-
61 * bal, but the handle (the integer
62 * value) may vary.
63 * M_A (global) DESCA( M_ ) The number of rows in the global
64 * array A.
65 * N_A (global) DESCA( N_ ) The number of columns in the global
66 * array A.
67 * MB_A (global) DESCA( MB_ ) The blocking factor used to distribute
68 * the rows of the array.
69 * NB_A (global) DESCA( NB_ ) The blocking factor used to distribute
70 * the columns of the array.
71 * RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
72 * row of the array A is distributed.
73 * CSRC_A (global) DESCA( CSRC_ ) The process column over which the
74 * first column of the array A is
75 * distributed.
76 * LLD_A (local) DESCA( LLD_ ) The leading dimension of the local
77 * array. LLD_A >= MAX(1,LOCr(M_A)).
78 *
79 * Let K be the number of rows or columns of a distributed matrix,
80 * and assume that its process grid has dimension p x q.
81 * LOCr( K ) denotes the number of elements of K that a process
82 * would receive if K were distributed over the p processes of its
83 * process column.
84 * Similarly, LOCc( K ) denotes the number of elements of K that a
85 * process would receive if K were distributed over the q processes of
86 * its process row.
87 * The values of LOCr() and LOCc() may be determined via a call to the
88 * ScaLAPACK tool function, NUMROC:
89 * LOCr( M ) = NUMROC( M, MB_A, MYROW, RSRC_A, NPROW ),
90 * LOCc( N ) = NUMROC( N, NB_A, MYCOL, CSRC_A, NPCOL ).
91 * An upper bound for these quantities may be computed by:
92 * LOCr( M ) <= ceil( ceil(M/MB_A)/NPROW )*MB_A
93 * LOCc( N ) <= ceil( ceil(N/NB_A)/NPCOL )*NB_A
94 *
95 * Arguments
96 * =========
97 *
98 * IBTYPE (global input) INTEGER
99 * = 1: compute inv(U**T)*sub( A )*inv(U) or
100 * inv(L)*sub( A )*inv(L**T);
101 * = 2 or 3: compute U*sub( A )*U**T or L**T*sub( A )*L.
102 *
103 * UPLO (global input) CHARACTER
104 * = 'U': Upper triangle of sub( A ) is stored and sub( B ) is
105 * factored as U**T*U;
106 * = 'L': Lower triangle of sub( A ) is stored and sub( B ) is
107 * factored as L*L**T.
108 *
109 * N (global input) INTEGER
110 * The order of the matrices sub( A ) and sub( B ). N >= 0.
111 *
112 * A (local input/local output) DOUBLE PRECISION pointer into the
113 * local memory to an array of dimension (LLD_A, LOCc(JA+N-1)).
114 * On entry, this array contains the local pieces of the
115 * N-by-N symmetric distributed matrix sub( A ). If UPLO = 'U',
116 * the leading N-by-N upper triangular part of sub( A ) contains
117 * the upper triangular part of the matrix, and its strictly
118 * lower triangular part is not referenced. If UPLO = 'L', the
119 * leading N-by-N lower triangular part of sub( A ) contains
120 * the lower triangular part of the matrix, and its strictly
121 * upper triangular part is not referenced.
122 *
123 * On exit, if INFO = 0, the transformed matrix, stored in the
124 * same format as sub( A ).
125 *
126 * IA (global input) INTEGER
127 * A's global row index, which points to the beginning of the
128 * submatrix which is to be operated on.
129 *
130 * JA (global input) INTEGER
131 * A's global column index, which points to the beginning of
132 * the submatrix which is to be operated on.
133 *
134 * DESCA (global and local input) INTEGER array of dimension DLEN_.
135 * The array descriptor for the distributed matrix A.
136 *
137 * B (local input) DOUBLE PRECISION pointer into the local memory
138 * to an array of dimension (LLD_B, LOCc(JB+N-1)). On entry,
139 * this array contains the local pieces of the triangular factor
140 * from the Cholesky factorization of sub( B ), as returned by
141 * PDPOTRF.
142 *
143 * IB (global input) INTEGER
144 * B's global row index, which points to the beginning of the
145 * submatrix which is to be operated on.
146 *
147 * JB (global input) INTEGER
148 * B's global column index, which points to the beginning of
149 * the submatrix which is to be operated on.
150 *
151 * DESCB (global and local input) INTEGER array of dimension DLEN_.
152 * The array descriptor for the distributed matrix B.
153 *
154 * SCALE (global output) DOUBLE PRECISION
155 * Amount by which the eigenvalues should be scaled to
156 * compensate for the scaling performed in this routine.
157 * At present, SCALE is always returned as 1.0, it is
158 * returned here to allow for future enhancement.
159 *
160 * INFO (global output) INTEGER
161 * = 0: successful exit
162 * < 0: If the i-th argument is an array and the j-entry had
163 * an illegal value, then INFO = -(i*100+j), if the i-th
164 * argument is a scalar and had an illegal value, then
165 * INFO = -i.
166 *
167 * =====================================================================
168 *
169 * .. Parameters ..
170  INTEGER BLOCK_CYCLIC_2D, DLEN_, DTYPE_, CTXT_, M_, N_,
171  $ mb_, nb_, rsrc_, csrc_, lld_
172  parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
173  $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
174  $ rsrc_ = 7, csrc_ = 8, lld_ = 9 )
175  DOUBLE PRECISION ONE, HALF
176  parameter( one = 1.0d+0, half = 0.5d+0 )
177 * ..
178 * .. Local Scalars ..
179  LOGICAL UPPER
180  INTEGER IACOL, IAROW, IBCOL, IBROW, ICOFFA, ICOFFB,
181  $ ictxt, iroffa, iroffb, k, kb, mycol, myrow, nb,
182  $ npcol, nprow
183 * ..
184 * .. Local Arrays ..
185  INTEGER IDUM1( 2 ), IDUM2( 2 )
186 * ..
187 * .. External Subroutines ..
188  EXTERNAL blacs_gridinfo, chk1mat, pchk2mat, pdsygs2,
189  $ pdsymm, pdsyr2k, pdtrmm, pdtrsm, pxerbla
190 * ..
191 * .. Intrinsic Functions ..
192  INTRINSIC ichar, min, mod
193 * ..
194 * .. External Functions ..
195  LOGICAL LSAME
196  INTEGER ICEIL, INDXG2P
197  EXTERNAL lsame, iceil, indxg2p
198 * ..
199 * .. Executable Statements ..
200 * This is just to keep ftnchek happy
201  IF( block_cyclic_2d*csrc_*ctxt_*dlen_*dtype_*lld_*mb_*m_*nb_*n_*
202  $ rsrc_.LT.0 )RETURN
203 *
204 * Get grid parameters
205 *
206  scale = one
207 *
208  ictxt = desca( ctxt_ )
209  CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
210 *
211 * Test the input parameters
212 *
213  info = 0
214  IF( nprow.EQ.-1 ) THEN
215  info = -( 700+ctxt_ )
216  ELSE
217  upper = lsame( uplo, 'U' )
218  CALL chk1mat( n, 3, n, 3, ia, ja, desca, 7, info )
219  CALL chk1mat( n, 3, n, 3, ib, jb, descb, 11, info )
220  IF( info.EQ.0 ) THEN
221  iarow = indxg2p( ia, desca( mb_ ), myrow, desca( rsrc_ ),
222  $ nprow )
223  ibrow = indxg2p( ib, descb( mb_ ), myrow, descb( rsrc_ ),
224  $ nprow )
225  iacol = indxg2p( ja, desca( nb_ ), mycol, desca( csrc_ ),
226  $ npcol )
227  ibcol = indxg2p( jb, descb( nb_ ), mycol, descb( csrc_ ),
228  $ npcol )
229  iroffa = mod( ia-1, desca( mb_ ) )
230  icoffa = mod( ja-1, desca( nb_ ) )
231  iroffb = mod( ib-1, descb( mb_ ) )
232  icoffb = mod( jb-1, descb( nb_ ) )
233  IF( ibtype.LT.1 .OR. ibtype.GT.3 ) THEN
234  info = -1
235  ELSE IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
236  info = -2
237  ELSE IF( n.LT.0 ) THEN
238  info = -3
239  ELSE IF( iroffa.NE.0 ) THEN
240  info = -5
241  ELSE IF( icoffa.NE.0 ) THEN
242  info = -6
243  ELSE IF( desca( mb_ ).NE.desca( nb_ ) ) THEN
244  info = -( 700+nb_ )
245  ELSE IF( iroffb.NE.0 .OR. ibrow.NE.iarow ) THEN
246  info = -9
247  ELSE IF( icoffb.NE.0 .OR. ibcol.NE.iacol ) THEN
248  info = -10
249  ELSE IF( descb( mb_ ).NE.desca( mb_ ) ) THEN
250  info = -( 1100+mb_ )
251  ELSE IF( descb( nb_ ).NE.desca( nb_ ) ) THEN
252  info = -( 1100+nb_ )
253  ELSE IF( ictxt.NE.descb( ctxt_ ) ) THEN
254  info = -( 1100+ctxt_ )
255  END IF
256  END IF
257  idum1( 1 ) = ibtype
258  idum2( 1 ) = 1
259  IF( upper ) THEN
260  idum1( 2 ) = ichar( 'U' )
261  ELSE
262  idum1( 2 ) = ichar( 'L' )
263  END IF
264  idum2( 2 ) = 2
265  CALL pchk2mat( n, 3, n, 3, ia, ja, desca, 7, n, 3, n, 3, ib,
266  $ jb, descb, 11, 2, idum1, idum2, info )
267  END IF
268 *
269  IF( info.NE.0 ) THEN
270  CALL pxerbla( ictxt, 'PDSYGST', -info )
271  RETURN
272  END IF
273 *
274 * Quick return if possible
275 *
276  IF( n.EQ.0 )
277  $ RETURN
278 *
279  IF( ibtype.EQ.1 ) THEN
280  IF( upper ) THEN
281 *
282 * Compute inv(U')*sub( A )*inv(U)
283 *
284  k = 1
285  nb = desca( nb_ )
286  kb = min( iceil( ja, nb )*nb, ja+n-1 ) - ja + 1
287 *
288  10 CONTINUE
289 *
290 * Update the upper triangle of A(ia+k-1:ia+n-1,ja+k-1:ja+n-1)
291 *
292  CALL pdsygs2( ibtype, uplo, kb, a, ia+k-1, ja+k-1, desca, b,
293  $ ib+k-1, ib+k-1, descb, info )
294  IF( k+kb.LE.n ) THEN
295  CALL pdtrsm( 'Left', uplo, 'Transpose', 'Non-unit', kb,
296  $ n-k-kb+1, one, b, ib+k-1, jb+k-1, descb, a,
297  $ ia+k-1, ja+k+kb-1, desca )
298  CALL pdsymm( 'Left', uplo, kb, n-k-kb+1, -half, a,
299  $ ia+k-1, ja+k-1, desca, b, ib+k-1, jb+k+kb-1,
300  $ descb, one, a, ia+k-1, ja+k+kb-1, desca )
301  CALL pdsyr2k( uplo, 'Transpose', n-k-kb+1, kb, -one, a,
302  $ ia+k-1, ja+k+kb-1, desca, b, ib+k-1,
303  $ jb+k+kb-1, descb, one, a, ia+k+kb-1,
304  $ ja+k+kb-1, desca )
305  CALL pdsymm( 'Left', uplo, kb, n-k-kb+1, -half, a,
306  $ ia+k-1, ja+k-1, desca, b, ib+k-1, jb+k+kb-1,
307  $ descb, one, a, ia+k-1, ja+k+kb-1, desca )
308  CALL pdtrsm( 'Right', uplo, 'No transpose', 'Non-unit',
309  $ kb, n-k-kb+1, one, b, ib+k+kb-1, jb+k+kb-1,
310  $ descb, a, ia+k-1, ja+k+kb-1, desca )
311  END IF
312  k = k + kb
313  kb = min( n-k+1, nb )
314 *
315  IF( k.LE.n )
316  $ GO TO 10
317 *
318  ELSE
319 *
320 * Compute inv(L)*sub( A )*inv(L')
321 *
322  k = 1
323  nb = desca( mb_ )
324  kb = min( iceil( ia, nb )*nb, ia+n-1 ) - ia + 1
325 *
326  20 CONTINUE
327 *
328 * Update the lower triangle of A(ia+k-1:ia+n-1,ja+k-1:ja+n-1)
329 *
330  CALL pdsygs2( ibtype, uplo, kb, a, ia+k-1, ja+k-1, desca, b,
331  $ ib+k-1, jb+k-1, descb, info )
332  IF( k+kb.LE.n ) THEN
333  CALL pdtrsm( 'Right', uplo, 'Transpose', 'Non-unit',
334  $ n-k-kb+1, kb, one, b, ib+k-1, jb+k-1, descb,
335  $ a, ia+k+kb-1, ja+k-1, desca )
336  CALL pdsymm( 'Right', uplo, n-k-kb+1, kb, -half, a,
337  $ ia+k-1, ja+k-1, desca, b, ib+k+kb-1, jb+k-1,
338  $ descb, one, a, ia+k+kb-1, ja+k-1, desca )
339  CALL pdsyr2k( uplo, 'No transpose', n-k-kb+1, kb, -one,
340  $ a, ia+k+kb-1, ja+k-1, desca, b, ib+k+kb-1,
341  $ jb+k-1, descb, one, a, ia+k+kb-1,
342  $ ja+k+kb-1, desca )
343  CALL pdsymm( 'Right', uplo, n-k-kb+1, kb, -half, a,
344  $ ia+k-1, ja+k-1, desca, b, ib+k+kb-1, jb+k-1,
345  $ descb, one, a, ia+k+kb-1, ja+k-1, desca )
346  CALL pdtrsm( 'Left', uplo, 'No transpose', 'Non-unit',
347  $ n-k-kb+1, kb, one, b, ib+k+kb-1, jb+k+kb-1,
348  $ descb, a, ia+k+kb-1, ja+k-1, desca )
349  END IF
350  k = k + kb
351  kb = min( n-k+1, nb )
352 *
353  IF( k.LE.n )
354  $ GO TO 20
355 *
356  END IF
357 *
358  ELSE
359 *
360  IF( upper ) THEN
361 *
362 * Compute U*sub( A )*U'
363 *
364  k = 1
365  nb = desca( nb_ )
366  kb = min( iceil( ja, nb )*nb, ja+n-1 ) - ja + 1
367 *
368  30 CONTINUE
369 *
370 * Update the upper triangle of A(ia:ia+k+kb-2,ja:ja+k+kb-2)
371 *
372  CALL pdtrmm( 'Left', uplo, 'No transpose', 'Non-unit', k-1,
373  $ kb, one, b, ib, jb, descb, a, ia, ja+k-1,
374  $ desca )
375  CALL pdsymm( 'Right', uplo, k-1, kb, half, a, ia+k-1,
376  $ ja+k-1, desca, b, ib, jb+k-1, descb, one, a,
377  $ ia, ja+k-1, desca )
378  CALL pdsyr2k( uplo, 'No transpose', k-1, kb, one, a, ia,
379  $ ja+k-1, desca, b, ib, jb+k-1, descb, one, a,
380  $ ia, ja, desca )
381  CALL pdsymm( 'Right', uplo, k-1, kb, half, a, ia+k-1,
382  $ ja+k-1, desca, b, ib, jb+k-1, descb, one, a,
383  $ ia, ja+k-1, desca )
384  CALL pdtrmm( 'Right', uplo, 'Transpose', 'Non-unit', k-1,
385  $ kb, one, b, ib+k-1, jb+k-1, descb, a, ia,
386  $ ja+k-1, desca )
387  CALL pdsygs2( ibtype, uplo, kb, a, ia+k-1, ja+k-1, desca, b,
388  $ ib+k-1, jb+k-1, descb, info )
389 *
390  k = k + kb
391  kb = min( n-k+1, nb )
392 *
393  IF( k.LE.n )
394  $ GO TO 30
395 *
396  ELSE
397 *
398 * Compute L'*sub( A )*L
399 *
400  k = 1
401  nb = desca( mb_ )
402  kb = min( iceil( ia, nb )*nb, ia+n-1 ) - ia + 1
403 *
404  40 CONTINUE
405 *
406 * Update the lower triangle of A(ia:ia+k+kb-2,ja:ja+k+kb-2)
407 *
408  CALL pdtrmm( 'Right', uplo, 'No transpose', 'Non-unit', kb,
409  $ k-1, one, b, ib, jb, descb, a, ia+k-1, ja,
410  $ desca )
411  CALL pdsymm( 'Left', uplo, kb, k-1, half, a, ia+k-1, ja+k-1,
412  $ desca, b, ib+k-1, jb, descb, one, a, ia+k-1,
413  $ ja, desca )
414  CALL pdsyr2k( uplo, 'Transpose', k-1, kb, one, a, ia+k-1,
415  $ ja, desca, b, ib+k-1, jb, descb, one, a, ia,
416  $ ja, desca )
417  CALL pdsymm( 'Left', uplo, kb, k-1, half, a, ia+k-1, ja+k-1,
418  $ desca, b, ib+k-1, jb, descb, one, a, ia+k-1,
419  $ ja, desca )
420  CALL pdtrmm( 'Left', uplo, 'Transpose', 'Non-unit', kb, k-1,
421  $ one, b, ib+k-1, jb+k-1, descb, a, ia+k-1, ja,
422  $ desca )
423  CALL pdsygs2( ibtype, uplo, kb, a, ia+k-1, ja+k-1, desca, b,
424  $ ib+k-1, jb+k-1, descb, info )
425 *
426  k = k + kb
427  kb = min( n-k+1, nb )
428 *
429  IF( k.LE.n )
430  $ GO TO 40
431 *
432  END IF
433 *
434  END IF
435 *
436  RETURN
437 *
438 * End of PDSYGST
439 *
440  END
pchk2mat
subroutine pchk2mat(MA, MAPOS0, NA, NAPOS0, IA, JA, DESCA, DESCAPOS0, MB, MBPOS0, NB, NBPOS0, IB, JB, DESCB, DESCBPOS0, NEXTRA, EX, EXPOS, INFO)
Definition: pchkxmat.f:175
chk1mat
subroutine chk1mat(MA, MAPOS0, NA, NAPOS0, IA, JA, DESCA, DESCAPOS0, INFO)
Definition: chk1mat.f:3
pdsygst
subroutine pdsygst(IBTYPE, UPLO, N, A, IA, JA, DESCA, B, IB, JB, DESCB, SCALE, INFO)
Definition: pdsygst.f:5
pxerbla
subroutine pxerbla(ICTXT, SRNAME, INFO)
Definition: pxerbla.f:2
pdsygs2
subroutine pdsygs2(IBTYPE, UPLO, N, A, IA, JA, DESCA, B, IB, JB, DESCB, INFO)
Definition: pdsygs2.f:5
min
#define min(A, B)
Definition: pcgemr.c:181