SCALAPACK 2.2.2
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
pdhseqr.f
Go to the documentation of this file.
1 SUBROUTINE pdhseqr( JOB, COMPZ, N, ILO, IHI, H, DESCH, WR, WI, Z,
2 $ DESCZ, WORK, LWORK, IWORK, LIWORK, INFO )
3*
4* Contribution from the Department of Computing Science and HPC2N,
5* Umea University, Sweden
6*
7* -- ScaLAPACK driver routine (version 2.0.1) --
8* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
9* Univ. of Colorado Denver and University of California, Berkeley.
10* January, 2012
11*
12 IMPLICIT NONE
13*
14* .. Scalar Arguments ..
15 INTEGER IHI, ILO, INFO, LWORK, LIWORK, N
16 CHARACTER COMPZ, JOB
17* ..
18* .. Array Arguments ..
19 INTEGER DESCH( * ) , DESCZ( * ), IWORK( * )
20 DOUBLE PRECISION H( * ), WI( N ), WORK( * ), WR( N ), Z( * )
21* ..
22* Purpose
23* =======
24*
25* PDHSEQR computes the eigenvalues of an upper Hessenberg matrix H
26* and, optionally, the matrices T and Z from the Schur decomposition
27* H = Z*T*Z**T, where T is an upper quasi-triangular matrix (the
28* Schur form), and Z is the orthogonal matrix of Schur vectors.
29*
30* Optionally Z may be postmultiplied into an input orthogonal
31* matrix Q so that this routine can give the Schur factorization
32* of a matrix A which has been reduced to the Hessenberg form H
33* by the orthogonal matrix Q: A = Q*H*Q**T = (QZ)*T*(QZ)**T.
34*
35* Notes
36* =====
37*
38* Each global data object is described by an associated description
39* vector. This vector stores the information required to establish
40* the mapping between an object element and its corresponding process
41* and memory location.
42*
43* Let A be a generic term for any 2D block cyclicly distributed array.
44* Such a global array has an associated description vector DESCA.
45* In the following comments, the character _ should be read as
46* "of the global array".
47*
48* NOTATION STORED IN EXPLANATION
49* --------------- -------------- --------------------------------------
50* DTYPE_A(global) DESCA( DTYPE_ )The descriptor type. In this case,
51* DTYPE_A = 1.
52* CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
53* the BLACS process grid A is distribu-
54* ted over. The context itself is glo-
55* bal, but the handle (the integer
56* value) may vary.
57* M_A (global) DESCA( M_ ) The number of rows in the global
58* array A.
59* N_A (global) DESCA( N_ ) The number of columns in the global
60* array A.
61* MB_A (global) DESCA( MB_ ) The blocking factor used to distribute
62* the rows of the array.
63* NB_A (global) DESCA( NB_ ) The blocking factor used to distribute
64* the columns of the array.
65* RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
66* row of the array A is distributed.
67* CSRC_A (global) DESCA( CSRC_ ) The process column over which the
68* first column of the array A is
69* distributed.
70* LLD_A (local) DESCA( LLD_ ) The leading dimension of the local
71* array. LLD_A >= MAX(1,LOCr(M_A)).
72*
73* Let K be the number of rows or columns of a distributed matrix,
74* and assume that its process grid has dimension p x q.
75* LOCr( K ) denotes the number of elements of K that a process
76* would receive if K were distributed over the p processes of its
77* process column.
78* Similarly, LOCc( K ) denotes the number of elements of K that a
79* process would receive if K were distributed over the q processes of
80* its process row.
81* The values of LOCr() and LOCc() may be determined via a call to the
82* ScaLAPACK tool function, NUMROC:
83* LOCr( M ) = NUMROC( M, MB_A, MYROW, RSRC_A, NPROW ),
84* LOCc( N ) = NUMROC( N, NB_A, MYCOL, CSRC_A, NPCOL ).
85* An upper bound for these quantities may be computed by:
86* LOCr( M ) <= ceil( ceil(M/MB_A)/NPROW )*MB_A
87* LOCc( N ) <= ceil( ceil(N/NB_A)/NPCOL )*NB_A
88*
89* Arguments
90* =========
91*
92* JOB (global input) CHARACTER*1
93* = 'E': compute eigenvalues only;
94* = 'S': compute eigenvalues and the Schur form T.
95*
96* COMPZ (global input) CHARACTER*1
97* = 'N': no Schur vectors are computed;
98* = 'I': Z is initialized to the unit matrix and the matrix Z
99* of Schur vectors of H is returned;
100* = 'V': Z must contain an orthogonal matrix Q on entry, and
101* the product Q*Z is returned.
102*
103* N (global input) INTEGER
104* The order of the Hessenberg matrix H (and Z if WANTZ).
105* N >= 0.
106*
107* ILO (global input) INTEGER
108* IHI (global input) INTEGER
109* It is assumed that H is already upper triangular in rows
110* and columns 1:ILO-1 and IHI+1:N. ILO and IHI are normally
111* set by a previous call to PDGEBAL, and then passed to PDGEHRD
112* when the matrix output by PDGEBAL is reduced to Hessenberg
113* form. Otherwise ILO and IHI should be set to 1 and N
114* respectively. If N.GT.0, then 1.LE.ILO.LE.IHI.LE.N.
115* If N = 0, then ILO = 1 and IHI = 0.
116*
117* H (global input/output) DOUBLE PRECISION array, dimension
118* (DESCH(LLD_),*)
119* On entry, the upper Hessenberg matrix H.
120* On exit, if JOB = 'S', H is upper quasi-triangular in
121* rows and columns ILO:IHI, with 1-by-1 and 2-by-2 blocks on
122* the main diagonal. The 2-by-2 diagonal blocks (corresponding
123* to complex conjugate pairs of eigenvalues) are returned in
124* standard form, with H(i,i) = H(i+1,i+1) and
125* H(i+1,i)*H(i,i+1).LT.0. If INFO = 0 and JOB = 'E', the
126* contents of H are unspecified on exit.
127*
128* DESCH (global and local input) INTEGER array of dimension DLEN_.
129* The array descriptor for the distributed matrix H.
130*
131* WR (global output) DOUBLE PRECISION array, dimension (N)
132* WI (global output) DOUBLE PRECISION array, dimension (N)
133* The real and imaginary parts, respectively, of the computed
134* eigenvalues ILO to IHI are stored in the corresponding
135* elements of WR and WI. If two eigenvalues are computed as a
136* complex conjugate pair, they are stored in consecutive
137* elements of WR and WI, say the i-th and (i+1)th, with
138* WI(i) > 0 and WI(i+1) < 0. If JOB = 'S', the
139* eigenvalues are stored in the same order as on the diagonal
140* of the Schur form returned in H.
141*
142* Z (global input/output) DOUBLE PRECISION array.
143* If COMPZ = 'V', on entry Z must contain the current
144* matrix Z of accumulated transformations from, e.g., PDGEHRD,
145* and on exit Z has been updated; transformations are applied
146* only to the submatrix Z(ILO:IHI,ILO:IHI).
147* If COMPZ = 'N', Z is not referenced.
148* If COMPZ = 'I', on entry Z need not be set and on exit,
149* if INFO = 0, Z contains the orthogonal matrix Z of the Schur
150* vectors of H.
151*
152* DESCZ (global and local input) INTEGER array of dimension DLEN_.
153* The array descriptor for the distributed matrix Z.
154*
155* WORK (local workspace) DOUBLE PRECISION array, dimension(LWORK)
156*
157* LWORK (local input) INTEGER
158* The length of the workspace array WORK.
159*
160* IWORK (local workspace) INTEGER array, dimension (LIWORK)
161*
162* LIWORK (local input) INTEGER
163* The length of the workspace array IWORK.
164*
165* INFO (output) INTEGER
166* = 0: successful exit
167* .LT. 0: if INFO = -i, the i-th argument had an illegal
168* value (see also below for -7777 and -8888).
169* .GT. 0: if INFO = i, PDHSEQR failed to compute all of
170* the eigenvalues. Elements 1:ilo-1 and i+1:n of WR
171* and WI contain those eigenvalues which have been
172* successfully computed. (Failures are rare.)
173*
174* If INFO .GT. 0 and JOB = 'E', then on exit, the
175* remaining unconverged eigenvalues are the eigen-
176* values of the upper Hessenberg matrix rows and
177* columns ILO through INFO of the final, output
178* value of H.
179*
180* If INFO .GT. 0 and JOB = 'S', then on exit
181*
182* (*) (initial value of H)*U = U*(final value of H)
183*
184* where U is an orthogonal matrix. The final
185* value of H is upper Hessenberg and quasi-triangular
186* in rows and columns INFO+1 through IHI.
187*
188* If INFO .GT. 0 and COMPZ = 'V', then on exit
189*
190* (final value of Z) = (initial value of Z)*U
191*
192* where U is the orthogonal matrix in (*) (regard-
193* less of the value of JOB.)
194*
195* If INFO .GT. 0 and COMPZ = 'I', then on exit
196* (final value of Z) = U
197* where U is the orthogonal matrix in (*) (regard-
198* less of the value of JOB.)
199*
200* If INFO .GT. 0 and COMPZ = 'N', then Z is not
201* accessed.
202*
203* = -7777: PDLAQR0 failed to converge and PDLAQR1 was called
204* instead. This could happen. Mostly due to a bug.
205* Please, send a bug report to the authors.
206* = -8888: PDLAQR1 failed to converge and PDLAQR0 was called
207* instead. This should not happen.
208*
209* ================================================================
210* Based on contributions by
211* Robert Granat, Department of Computing Science and HPC2N,
212* Umea University, Sweden.
213* ================================================================
214*
215* Restrictions: The block size in H and Z must be square and larger
216* than or equal to six (6) due to restrictions in PDLAQR1, PDLAQR5
217* and DLAQR6. Moreover, H and Z need to be distributed identically
218* with the same context.
219*
220* ================================================================
221* References:
222* K. Braman, R. Byers, and R. Mathias,
223* The Multi-Shift QR Algorithm Part I: Maintaining Well Focused
224* Shifts, and Level 3 Performance.
225* SIAM J. Matrix Anal. Appl., 23(4):929--947, 2002.
226*
227* K. Braman, R. Byers, and R. Mathias,
228* The Multi-Shift QR Algorithm Part II: Aggressive Early
229* Deflation.
230* SIAM J. Matrix Anal. Appl., 23(4):948--973, 2002.
231*
232* R. Granat, B. Kagstrom, and D. Kressner,
233* A Novel Parallel QR Algorithm for Hybrid Distributed Momory HPC
234* Systems.
235* SIAM J. Sci. Comput., 32(4):2345--2378, 2010.
236*
237* ================================================================
238* .. Parameters ..
239 INTEGER BLOCK_CYCLIC_2D, CSRC_, CTXT_, DLEN_, DTYPE_,
240 $ lld_, mb_, m_, nb_, n_, rsrc_
241 LOGICAL CRSOVER
242 parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
243 $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
244 $ rsrc_ = 7, csrc_ = 8, lld_ = 9,
245 $ crsover = .true. )
246 INTEGER NTINY
247 parameter( ntiny = 11 )
248 INTEGER NL
249 parameter( nl = 49 )
250 DOUBLE PRECISION ZERO, ONE
251 parameter( zero = 0.0d0, one = 1.0d0 )
252* ..
253* .. Local Scalars ..
254 INTEGER I, KBOT, NMIN, LLDH, LLDZ, ICTXT, NPROW, NPCOL,
255 $ myrow, mycol, hrows, hcols, ipw, nh, nb,
256 $ ii, jj, hrsrc, hcsrc, nprocs, iloc1, jloc1,
257 $ hrsrc1, hcsrc1, k, iloc2, jloc2, iloc3, jloc3,
258 $ iloc4, jloc4, hrsrc2, hcsrc2, hrsrc3, hcsrc3,
259 $ hrsrc4, hcsrc4, liwkopt
260 LOGICAL INITZ, LQUERY, WANTT, WANTZ, PAIR, BORDER
261 DOUBLE PRECISION TMP1, TMP2, TMP3, TMP4, DUM1, DUM2, DUM3,
262 $ dum4, elem1, elem4,
263 $ cs, sn, elem5, tmp, lwkopt
264* ..
265* .. Local Arrays ..
266 INTEGER DESCH2( DLEN_ )
267 DOUBLE PRECISION ELEM2( 1 ), ELEM3( 1 )
268* ..
269* .. External Functions ..
270 INTEGER PILAENVX, NUMROC, ICEIL
271 LOGICAL LSAME
272 EXTERNAL pilaenvx, lsame, numroc, iceil
273* ..
274* .. External Subroutines ..
276* ..
277* .. Intrinsic Functions ..
278 INTRINSIC dble, max, min
279* ..
280* .. Executable Statements ..
281*
282* Decode and check the input parameters.
283*
284 info = 0
285 ictxt = desch( ctxt_ )
286 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
287 nprocs = nprow*npcol
288 IF( nprow.EQ.-1 ) info = -(600+ctxt_)
289 IF( info.EQ.0 ) THEN
290 wantt = lsame( job, 'S' )
291 initz = lsame( compz, 'I' )
292 wantz = initz .OR. lsame( compz, 'V' )
293 lldh = desch( lld_ )
294 lldz = descz( lld_ )
295 nb = desch( mb_ )
296 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
297*
298 IF( .NOT.lsame( job, 'E' ) .AND. .NOT.wantt ) THEN
299 info = -1
300 ELSE IF( .NOT.lsame( compz, 'N' ) .AND. .NOT.wantz ) THEN
301 info = -2
302 ELSE IF( n.LT.0 ) THEN
303 info = -3
304 ELSE IF( ilo.LT.1 .OR. ilo.GT.max( 1, n ) ) THEN
305 info = -4
306 ELSE IF( ihi.LT.min( ilo, n ) .OR. ihi.GT.n ) THEN
307 info = -5
308 ELSEIF( descz( ctxt_ ).NE.desch( ctxt_ ) ) THEN
309 info = -( 1000+ctxt_ )
310 ELSEIF( desch( mb_ ).NE.desch( nb_ ) ) THEN
311 info = -( 700+nb_ )
312 ELSEIF( descz( mb_ ).NE.descz( nb_ ) ) THEN
313 info = -( 1000+nb_ )
314 ELSEIF( desch( mb_ ).NE.descz( mb_ ) ) THEN
315 info = -( 1000+mb_ )
316 ELSEIF( desch( mb_ ).LT.6 ) THEN
317 info = -( 700+nb_ )
318 ELSEIF( descz( mb_ ).LT.6 ) THEN
319 info = -( 1000+mb_ )
320 ELSE
321 CALL chk1mat( n, 3, n, 3, 1, 1, desch, 7, info )
322 IF( info.EQ.0 )
323 $ CALL chk1mat( n, 3, n, 3, 1, 1, descz, 11, info )
324 IF( info.EQ.0 )
325 $ CALL pchk2mat( n, 3, n, 3, 1, 1, desch, 7, n, 3, n, 3,
326 $ 1, 1, descz, 11, 0, iwork, iwork, info )
327 END IF
328 END IF
329*
330* Compute required workspace.
331*
332 CALL pdlaqr1( wantt, wantz, n, ilo, ihi, h, desch, wr, wi,
333 $ ilo, ihi, z, descz, work, -1, iwork, -1, info )
334 lwkopt = work(1)
335 liwkopt = iwork(1)
336 CALL pdlaqr0( wantt, wantz, n, ilo, ihi, h, desch, wr, wi,
337 $ ilo, ihi, z, descz, work, -1, iwork, -1, info, 0 )
338 IF( n.LT.nl ) THEN
339 hrows = numroc( nl, nb, myrow, desch(rsrc_), nprow )
340 hcols = numroc( nl, nb, mycol, desch(csrc_), npcol )
341 work(1) = work(1) + dble(2*hrows*hcols)
342 END IF
343 lwkopt = max( lwkopt, work(1) )
344 liwkopt = max( liwkopt, iwork(1) )
345 work(1) = lwkopt
346 iwork(1) = liwkopt
347*
348 IF( .NOT.lquery .AND. lwork.LT.int(lwkopt) ) THEN
349 info = -13
350 ELSEIF( .NOT.lquery .AND. liwork.LT.liwkopt ) THEN
351 info = -15
352 END IF
353*
354 IF( info.NE.0 ) THEN
355*
356* Quick return in case of invalid argument.
357*
358 CALL pxerbla( ictxt, 'PDHSEQR', -info )
359 RETURN
360*
361 ELSE IF( n.EQ.0 ) THEN
362*
363* Quick return in case N = 0; nothing to do.
364*
365 RETURN
366*
367 ELSE IF( lquery ) THEN
368*
369* Quick return in case of a workspace query.
370*
371 RETURN
372*
373 ELSE
374*
375* Copy eigenvalues isolated by PDGEBAL.
376*
377 DO 10 i = 1, ilo - 1
378 CALL infog2l( i, i, desch, nprow, npcol, myrow, mycol, ii,
379 $ jj, hrsrc, hcsrc )
380 IF( myrow.EQ.hrsrc .AND. mycol.EQ.hcsrc ) THEN
381 wr( i ) = h( (jj-1)*lldh + ii )
382 ELSE
383 wr( i ) = zero
384 END IF
385 wi( i ) = zero
386 10 CONTINUE
387 IF( ilo.GT.1 )
388 $ CALL dgsum2d( ictxt, 'All', '1-Tree', ilo-1, 1, wr, n, -1,
389 $ -1 )
390 DO 20 i = ihi + 1, n
391 CALL infog2l( i, i, desch, nprow, npcol, myrow, mycol, ii,
392 $ jj, hrsrc, hcsrc )
393 IF( myrow.EQ.hrsrc .AND. mycol.EQ.hcsrc ) THEN
394 wr( i ) = h( (jj-1)*lldh + ii )
395 ELSE
396 wr( i ) = zero
397 END IF
398 wi( i ) = zero
399 20 CONTINUE
400 IF( ihi.LT.n )
401 $ CALL dgsum2d( ictxt, 'All', '1-Tree', n-ihi, 1, wr(ihi+1),
402 $ n, -1, -1 )
403*
404* Initialize Z, if requested.
405*
406 IF( initz )
407 $ CALL pdlaset( 'A', n, n, zero, one, z, 1, 1, descz )
408*
409* Quick return if possible.
410*
411 nprocs = nprow*npcol
412 IF( ilo.EQ.ihi ) THEN
413 CALL infog2l( ilo, ilo, desch, nprow, npcol, myrow,
414 $ mycol, ii, jj, hrsrc, hcsrc )
415 IF( myrow.EQ.hrsrc .AND. mycol.EQ.hcsrc ) THEN
416 wr( ilo ) = h( (jj-1)*lldh + ii )
417 IF( nprocs.GT.1 )
418 $ CALL dgebs2d( ictxt, 'All', '1-Tree', 1, 1, wr(ilo),
419 $ 1 )
420 ELSE
421 CALL dgebr2d( ictxt, 'All', '1-Tree', 1, 1, wr(ilo),
422 $ 1, hrsrc, hcsrc )
423 END IF
424 wi( ilo ) = zero
425 RETURN
426 END IF
427*
428* PDLAQR1/PDLAQR0 crossover point.
429*
430 nh = ihi-ilo+1
431 nmin = pilaenvx( ictxt, 12, 'PDHSEQR',
432 $ job( : 1 ) // compz( : 1 ), n, ilo, ihi, lwork )
433 nmin = max( ntiny, nmin )
434*
435* PDLAQR0 for big matrices; PDLAQR1 for small ones.
436*
437 IF( (.NOT. crsover .AND. nh.GT.ntiny) .OR. nh.GT.nmin .OR.
438 $ desch(rsrc_).NE.0 .OR. desch(csrc_).NE.0 ) THEN
439 CALL pdlaqr0( wantt, wantz, n, ilo, ihi, h, desch, wr, wi,
440 $ ilo, ihi, z, descz, work, lwork, iwork, liwork, info,
441 $ 0 )
442 IF( info.GT.0 .AND. ( desch(rsrc_).NE.0 .OR.
443 $ desch(csrc_).NE.0 ) ) THEN
444*
445* A rare PDLAQR0 failure! PDLAQR1 sometimes succeeds
446* when PDLAQR0 fails.
447*
448 kbot = info
449 CALL pdlaqr1( wantt, wantz, n, ilo, ihi, h, desch, wr,
450 $ wi, ilo, ihi, z, descz, work, lwork, iwork,
451 $ liwork, info )
452 info = -7777
453 END IF
454 ELSE
455*
456* Small matrix.
457*
458 CALL pdlaqr1( wantt, wantz, n, ilo, ihi, h, desch, wr, wi,
459 $ ilo, ihi, z, descz, work, lwork, iwork, liwork, info )
460*
461 IF( info.GT.0 ) THEN
462*
463* A rare PDLAQR1 failure! PDLAQR0 sometimes succeeds
464* when PDLAQR1 fails.
465*
466 kbot = info
467*
468 IF( n.GE.nl ) THEN
469*
470* Larger matrices have enough subdiagonal scratch
471* space to call PDLAQR0 directly.
472*
473 CALL pdlaqr0( wantt, wantz, n, ilo, kbot, h, desch,
474 $ wr, wi, ilo, ihi, z, descz, work, lwork,
475 $ iwork, liwork, info, 0 )
476 ELSE
477*
478* Tiny matrices don't have enough subdiagonal
479* scratch space to benefit from PDLAQR0. Hence,
480* tiny matrices must be copied into a larger
481* array before calling PDLAQR0.
482*
483 hrows = numroc( nl, nb, myrow, desch(rsrc_), nprow )
484 hcols = numroc( nl, nb, mycol, desch(csrc_), npcol )
485 CALL descinit( desch2, nl, nl, nb, nb, desch(rsrc_),
486 $ desch(csrc_), ictxt, max(1, hrows), info )
487 CALL pdlacpy( 'All', n, n, h, 1, 1, desch, work, 1,
488 $ 1, desch2 )
489 CALL pdelset( work, n+1, n, desch2, zero )
490 CALL pdlaset( 'All', nl, nl-n, zero, zero, work, 1,
491 $ n+1, desch2 )
492 ipw = 1 + desch2(lld_)*hcols
493 CALL pdlaqr0( wantt, wantz, nl, ilo, kbot, work,
494 $ desch2, wr, wi, ilo, ihi, z, descz,
495 $ work(ipw), lwork-ipw+1, iwork,
496 $ liwork, info, 0 )
497 IF( wantt .OR. info.NE.0 )
498 $ CALL pdlacpy( 'All', n, n, work, 1, 1, desch2,
499 $ h, 1, 1, desch )
500 END IF
501 info = -8888
502 END IF
503 END IF
504*
505* Clear out the trash, if necessary.
506*
507 IF( ( wantt .OR. info.NE.0 ) .AND. n.GT.2 )
508 $ CALL pdlaset( 'L', n-2, n-2, zero, zero, h, 3, 1, desch )
509*
510* Force any 2-by-2 blocks to be complex conjugate pairs of
511* eigenvalues by removing false such blocks.
512*
513 DO 30 i = ilo, ihi-1
514 CALL pdelget( 'All', ' ', tmp3, h, i+1, i, desch )
515 IF( tmp3.NE.0.0d+00 ) THEN
516 CALL pdelget( 'All', ' ', tmp1, h, i, i, desch )
517 CALL pdelget( 'All', ' ', tmp2, h, i, i+1, desch )
518 CALL pdelget( 'All', ' ', tmp4, h, i+1, i+1, desch )
519 CALL dlanv2( tmp1, tmp2, tmp3, tmp4, dum1, dum2, dum3,
520 $ dum4, cs, sn )
521 IF( tmp3.EQ.0.0d+00 ) THEN
522 IF( wantt ) THEN
523 IF( i+2.LE.n )
524 $ CALL pdrot( n-i-1, h, i, i+2, desch,
525 $ desch(m_), h, i+1, i+2, desch, desch(m_),
526 $ cs, sn, work, lwork, info )
527 CALL pdrot( i-1, h, 1, i, desch, 1, h, 1, i+1,
528 $ desch, 1, cs, sn, work, lwork, info )
529 END IF
530 IF( wantz ) THEN
531 CALL pdrot( n, z, 1, i, descz, 1, z, 1, i+1, descz,
532 $ 1, cs, sn, work, lwork, info )
533 END IF
534 CALL pdelset( h, i, i, desch, tmp1 )
535 CALL pdelset( h, i, i+1, desch, tmp2 )
536 CALL pdelset( h, i+1, i, desch, tmp3 )
537 CALL pdelset( h, i+1, i+1, desch, tmp4 )
538 END IF
539 END IF
540 30 CONTINUE
541*
542* Read out eigenvalues: first let all the processes compute the
543* eigenvalue inside their diagonal blocks in parallel, except for
544* the eigenvalue located next to a block border. After that,
545* compute all eigenvalues located next to the block borders.
546* Finally, do a global summation over WR and WI so that all
547* processors receive the result.
548*
549 DO 40 k = ilo, ihi
550 wr( k ) = zero
551 wi( k ) = zero
552 40 CONTINUE
553 nb = desch( mb_ )
554*
555* Loop 50: extract eigenvalues from the blocks which are not laid
556* out across a border of the processor mesh, except for those 1x1
557* blocks on the border.
558*
559 pair = .false.
560 DO 50 k = ilo, ihi
561 IF( .NOT. pair ) THEN
562 border = mod( k, nb ).EQ.0 .OR. ( k.NE.1 .AND.
563 $ mod( k, nb ).EQ.1 )
564 IF( .NOT. border ) THEN
565 CALL infog2l( k, k, desch, nprow, npcol, myrow,
566 $ mycol, iloc1, jloc1, hrsrc1, hcsrc1 )
567 IF( myrow.EQ.hrsrc1 .AND. mycol.EQ.hcsrc1 ) THEN
568 elem1 = h((jloc1-1)*lldh+iloc1)
569 IF( k.LT.n ) THEN
570 elem3( 1 ) = h((jloc1-1)*lldh+iloc1+1)
571 ELSE
572 elem3( 1 ) = zero
573 END IF
574 IF( elem3( 1 ).NE.zero ) THEN
575 elem2( 1 ) = h((jloc1)*lldh+iloc1)
576 elem4 = h((jloc1)*lldh+iloc1+1)
577 CALL dlanv2( elem1, elem2( 1 ), elem3( 1 ),
578 $ elem4, wr( k ), wi( k ), wr( k+1 ),
579 $ wi( k+1 ), sn, cs )
580 pair = .true.
581 ELSE
582 IF( k.GT.1 ) THEN
583 tmp = h((jloc1-2)*lldh+iloc1)
584 IF( tmp.NE.zero ) THEN
585 elem1 = h((jloc1-2)*lldh+iloc1-1)
586 elem2( 1 ) = h((jloc1-1)*lldh+iloc1-1)
587 elem3( 1 ) = h((jloc1-2)*lldh+iloc1)
588 elem4 = h((jloc1-1)*lldh+iloc1)
589 CALL dlanv2( elem1, elem2( 1 ),
590 $ elem3( 1 ), elem4, wr( k-1 ),
591 $ wi( k-1 ), wr( k ), wi( k ), sn, cs )
592 ELSE
593 wr( k ) = elem1
594 END IF
595 ELSE
596 wr( k ) = elem1
597 END IF
598 END IF
599 END IF
600 END IF
601 ELSE
602 pair = .false.
603 END IF
604 50 CONTINUE
605*
606* Loop 60: extract eigenvalues from the blocks which are laid
607* out across a border of the processor mesh. The processors are
608* numbered as below:
609*
610* 1 | 2
611* --+--
612* 3 | 4
613*
614 DO 60 k = iceil(ilo,nb)*nb, ihi-1, nb
615 CALL infog2l( k, k, desch, nprow, npcol, myrow, mycol,
616 $ iloc1, jloc1, hrsrc1, hcsrc1 )
617 CALL infog2l( k, k+1, desch, nprow, npcol, myrow, mycol,
618 $ iloc2, jloc2, hrsrc2, hcsrc2 )
619 CALL infog2l( k+1, k, desch, nprow, npcol, myrow, mycol,
620 $ iloc3, jloc3, hrsrc3, hcsrc3 )
621 CALL infog2l( k+1, k+1, desch, nprow, npcol, myrow, mycol,
622 $ iloc4, jloc4, hrsrc4, hcsrc4 )
623 IF( myrow.EQ.hrsrc2 .AND. mycol.EQ.hcsrc2 ) THEN
624 elem2( 1 ) = h((jloc2-1)*lldh+iloc2)
625 IF( hrsrc1.NE.hrsrc2 .OR. hcsrc1.NE.hcsrc2 )
626 $ CALL dgesd2d( ictxt, 1, 1, elem2, 1, hrsrc1, hcsrc1)
627 END IF
628 IF( myrow.EQ.hrsrc3 .AND. mycol.EQ.hcsrc3 ) THEN
629 elem3( 1 ) = h((jloc3-1)*lldh+iloc3)
630 IF( hrsrc1.NE.hrsrc3 .OR. hcsrc1.NE.hcsrc3 )
631 $ CALL dgesd2d( ictxt, 1, 1, elem3, 1, hrsrc1, hcsrc1)
632 END IF
633 IF( myrow.EQ.hrsrc4 .AND. mycol.EQ.hcsrc4 ) THEN
634 work(1) = h((jloc4-1)*lldh+iloc4)
635 IF( k+1.LT.n ) THEN
636 work(2) = h((jloc4-1)*lldh+iloc4+1)
637 ELSE
638 work(2) = zero
639 END IF
640 IF( hrsrc1.NE.hrsrc4 .OR. hcsrc1.NE.hcsrc4 )
641 $ CALL dgesd2d( ictxt, 2, 1, work, 2, hrsrc1, hcsrc1 )
642 END IF
643 IF( myrow.EQ.hrsrc1 .AND. mycol.EQ.hcsrc1 ) THEN
644 elem1 = h((jloc1-1)*lldh+iloc1)
645 IF( hrsrc1.NE.hrsrc2 .OR. hcsrc1.NE.hcsrc2 )
646 $ CALL dgerv2d( ictxt, 1, 1, elem2, 1, hrsrc2, hcsrc2)
647 IF( hrsrc1.NE.hrsrc3 .OR. hcsrc1.NE.hcsrc3 )
648 $ CALL dgerv2d( ictxt, 1, 1, elem3, 1, hrsrc3, hcsrc3)
649 IF( hrsrc1.NE.hrsrc4 .OR. hcsrc1.NE.hcsrc4 )
650 $ CALL dgerv2d( ictxt, 2, 1, work, 2, hrsrc4, hcsrc4 )
651 elem4 = work(1)
652 elem5 = work(2)
653 IF( elem5.EQ.zero ) THEN
654 IF( wr( k ).EQ.zero .AND. wi( k ).EQ.zero ) THEN
655 CALL dlanv2( elem1, elem2( 1 ), elem3( 1 ), elem4,
656 $ wr( k ), wi( k ), wr( k+1 ), wi( k+1 ), sn,
657 $ cs )
658 ELSEIF( wr( k+1 ).EQ.zero .AND. wi( k+1 ).EQ.zero )
659 $ THEN
660 wr( k+1 ) = elem4
661 END IF
662 ELSEIF( wr( k ).EQ.zero .AND. wi( k ).EQ.zero )
663 $ THEN
664 wr( k ) = elem1
665 END IF
666 END IF
667 60 CONTINUE
668*
669 IF( nprocs.GT.1 ) THEN
670 CALL dgsum2d( ictxt, 'All', ' ', ihi-ilo+1, 1, wr(ilo), n,
671 $ -1, -1 )
672 CALL dgsum2d( ictxt, 'All', ' ', ihi-ilo+1, 1, wi(ilo), n,
673 $ -1, -1 )
674 END IF
675*
676 END IF
677*
678 work(1) = lwkopt
679 iwork(1) = liwkopt
680 RETURN
681*
682* End of PDHSEQR
683*
684 END
subroutine chk1mat(ma, mapos0, na, napos0, ia, ja, desca, descapos0, info)
Definition chk1mat.f:3
subroutine descinit(desc, m, n, mb, nb, irsrc, icsrc, ictxt, lld, info)
Definition descinit.f:3
subroutine infog2l(grindx, gcindx, desc, nprow, npcol, myrow, mycol, lrindx, lcindx, rsrc, csrc)
Definition infog2l.f:3
#define max(A, B)
Definition pcgemr.c:180
#define min(A, B)
Definition pcgemr.c:181
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
subroutine pdlaset(uplo, m, n, alpha, beta, a, ia, ja, desca)
Definition pdblastst.f:6862
subroutine pdelget(scope, top, alpha, a, ia, ja, desca)
Definition pdelget.f:2
subroutine pdelset(a, ia, ja, desca, alpha)
Definition pdelset.f:2
subroutine pdhseqr(job, compz, n, ilo, ihi, h, desch, wr, wi, z, descz, work, lwork, iwork, liwork, info)
Definition pdhseqr.f:3
subroutine pdlacpy(uplo, m, n, a, ia, ja, desca, b, ib, jb, descb)
Definition pdlacpy.f:3
recursive subroutine pdlaqr0(wantt, wantz, n, ilo, ihi, h, desch, wr, wi, iloz, ihiz, z, descz, work, lwork, iwork, liwork, info, reclevel)
Definition pdlaqr0.f:4
recursive subroutine pdlaqr1(wantt, wantz, n, ilo, ihi, a, desca, wr, wi, iloz, ihiz, z, descz, work, lwork, iwork, ilwork, info)
Definition pdlaqr1.f:5
subroutine pdrot(n, x, ix, jx, descx, incx, y, iy, jy, descy, incy, cs, sn, work, lwork, info)
Definition pdrot.f:3
subroutine pxerbla(ictxt, srname, info)
Definition pxerbla.f:2