LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
slasq2.f
Go to the documentation of this file.
1*> \brief \b SLASQ2 computes all the eigenvalues of the symmetric positive definite tridiagonal matrix associated with the qd Array Z to high relative accuracy. Used by sbdsqr and sstegr.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SLASQ2 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/slasq2.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/slasq2.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/slasq2.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SLASQ2( N, Z, INFO )
20*
21* .. Scalar Arguments ..
22* INTEGER INFO, N
23* ..
24* .. Array Arguments ..
25* REAL Z( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> SLASQ2 computes all the eigenvalues of the symmetric positive
35*> definite tridiagonal matrix associated with the qd array Z to high
36*> relative accuracy are computed to high relative accuracy, in the
37*> absence of denormalization, underflow and overflow.
38*>
39*> To see the relation of Z to the tridiagonal matrix, let L be a
40*> unit lower bidiagonal matrix with subdiagonals Z(2,4,6,,..) and
41*> let U be an upper bidiagonal matrix with 1's above and diagonal
42*> Z(1,3,5,,..). The tridiagonal is L*U or, if you prefer, the
43*> symmetric tridiagonal to which it is similar.
44*>
45*> Note : SLASQ2 defines a logical variable, IEEE, which is true
46*> on machines which follow ieee-754 floating-point standard in their
47*> handling of infinities and NaNs, and false otherwise. This variable
48*> is passed to SLASQ3.
49*> \endverbatim
50*
51* Arguments:
52* ==========
53*
54*> \param[in] N
55*> \verbatim
56*> N is INTEGER
57*> The number of rows and columns in the matrix. N >= 0.
58*> \endverbatim
59*>
60*> \param[in,out] Z
61*> \verbatim
62*> Z is REAL array, dimension ( 4*N )
63*> On entry Z holds the qd array. On exit, entries 1 to N hold
64*> the eigenvalues in decreasing order, Z( 2*N+1 ) holds the
65*> trace, and Z( 2*N+2 ) holds the sum of the eigenvalues. If
66*> N > 2, then Z( 2*N+3 ) holds the iteration count, Z( 2*N+4 )
67*> holds NDIVS/NIN^2, and Z( 2*N+5 ) holds the percentage of
68*> shifts that failed.
69*> \endverbatim
70*>
71*> \param[out] INFO
72*> \verbatim
73*> INFO is INTEGER
74*> = 0: successful exit
75*> < 0: if the i-th argument is a scalar and had an illegal
76*> value, then INFO = -i, if the i-th argument is an
77*> array and the j-entry had an illegal value, then
78*> INFO = -(i*100+j)
79*> > 0: the algorithm failed
80*> = 1, a split was marked by a positive value in E
81*> = 2, current block of Z not diagonalized after 100*N
82*> iterations (in inner while loop). On exit Z holds
83*> a qd array with the same eigenvalues as the given Z.
84*> = 3, termination criterion of outer while loop not met
85*> (program created more than N unreduced blocks)
86*> \endverbatim
87*
88* Authors:
89* ========
90*
91*> \author Univ. of Tennessee
92*> \author Univ. of California Berkeley
93*> \author Univ. of Colorado Denver
94*> \author NAG Ltd.
95*
96*> \ingroup lasq2
97*
98*> \par Further Details:
99* =====================
100*>
101*> \verbatim
102*>
103*> Local Variables: I0:N0 defines a current unreduced segment of Z.
104*> The shifts are accumulated in SIGMA. Iteration count is in ITER.
105*> Ping-pong is controlled by PP (alternates between 0 and 1).
106*> \endverbatim
107*>
108* =====================================================================
109 SUBROUTINE slasq2( N, Z, INFO )
110*
111* -- LAPACK computational routine --
112* -- LAPACK is a software package provided by Univ. of Tennessee, --
113* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
114*
115* .. Scalar Arguments ..
116 INTEGER INFO, N
117* ..
118* .. Array Arguments ..
119 REAL Z( * )
120* ..
121*
122* =====================================================================
123*
124* .. Parameters ..
125 REAL CBIAS
126 parameter( cbias = 1.50e0 )
127 REAL ZERO, HALF, ONE, TWO, FOUR, HUNDRD
128 parameter( zero = 0.0e0, half = 0.5e0, one = 1.0e0,
129 $ two = 2.0e0, four = 4.0e0, hundrd = 100.0e0 )
130* ..
131* .. Local Scalars ..
132 LOGICAL IEEE
133 INTEGER I0, I4, IINFO, IPN4, ITER, IWHILA, IWHILB, K,
134 $ KMIN, N0, NBIG, NDIV, NFAIL, PP, SPLT, TTYPE,
135 $ I1, N1
136 REAL D, DEE, DEEMIN, DESIG, DMIN, DMIN1, DMIN2, DN,
137 $ DN1, DN2, E, EMAX, EMIN, EPS, G, OLDEMN, QMAX,
138 $ QMIN, S, SAFMIN, SIGMA, T, TAU, TEMP, TOL,
139 $ TOL2, TRACE, ZMAX, TEMPE, TEMPQ
140* ..
141* .. External Subroutines ..
142 EXTERNAL slasq3, slasrt, xerbla
143* ..
144* .. External Functions ..
145 REAL SLAMCH
146 EXTERNAL slamch
147* ..
148* .. Intrinsic Functions ..
149 INTRINSIC abs, max, min, real, sqrt
150* ..
151* .. Executable Statements ..
152*
153* Test the input arguments.
154* (in case SLASQ2 is not called by SLASQ1)
155*
156 info = 0
157 eps = slamch( 'Precision' )
158 safmin = slamch( 'Safe minimum' )
159 tol = eps*hundrd
160 tol2 = tol**2
161*
162 IF( n.LT.0 ) THEN
163 info = -1
164 CALL xerbla( 'SLASQ2', 1 )
165 RETURN
166 ELSE IF( n.EQ.0 ) THEN
167 RETURN
168 ELSE IF( n.EQ.1 ) THEN
169*
170* 1-by-1 case.
171*
172 IF( z( 1 ).LT.zero ) THEN
173 info = -201
174 CALL xerbla( 'SLASQ2', 2 )
175 END IF
176 RETURN
177 ELSE IF( n.EQ.2 ) THEN
178*
179* 2-by-2 case.
180*
181 IF( z( 1 ).LT.zero ) THEN
182 info = -201
183 CALL xerbla( 'SLASQ2', 2 )
184 RETURN
185 ELSE IF( z( 2 ).LT.zero ) THEN
186 info = -202
187 CALL xerbla( 'SLASQ2', 2 )
188 RETURN
189 ELSE IF( z( 3 ).LT.zero ) THEN
190 info = -203
191 CALL xerbla( 'SLASQ2', 2 )
192 RETURN
193 ELSE IF( z( 3 ).GT.z( 1 ) ) THEN
194 d = z( 3 )
195 z( 3 ) = z( 1 )
196 z( 1 ) = d
197 END IF
198 z( 5 ) = z( 1 ) + z( 2 ) + z( 3 )
199 IF( z( 2 ).GT.z( 3 )*tol2 ) THEN
200 t = half*( ( z( 1 )-z( 3 ) )+z( 2 ) )
201 s = z( 3 )*( z( 2 ) / t )
202 IF( s.LE.t ) THEN
203 s = z( 3 )*( z( 2 ) / ( t*( one+sqrt( one+s / t ) ) ) )
204 ELSE
205 s = z( 3 )*( z( 2 ) / ( t+sqrt( t )*sqrt( t+s ) ) )
206 END IF
207 t = z( 1 ) + ( s+z( 2 ) )
208 z( 3 ) = z( 3 )*( z( 1 ) / t )
209 z( 1 ) = t
210 END IF
211 z( 2 ) = z( 3 )
212 z( 6 ) = z( 2 ) + z( 1 )
213 RETURN
214 END IF
215*
216* Check for negative data and compute sums of q's and e's.
217*
218 z( 2*n ) = zero
219 emin = z( 2 )
220 qmax = zero
221 zmax = zero
222 d = zero
223 e = zero
224*
225 DO 10 k = 1, 2*( n-1 ), 2
226 IF( z( k ).LT.zero ) THEN
227 info = -( 200+k )
228 CALL xerbla( 'SLASQ2', 2 )
229 RETURN
230 ELSE IF( z( k+1 ).LT.zero ) THEN
231 info = -( 200+k+1 )
232 CALL xerbla( 'SLASQ2', 2 )
233 RETURN
234 END IF
235 d = d + z( k )
236 e = e + z( k+1 )
237 qmax = max( qmax, z( k ) )
238 emin = min( emin, z( k+1 ) )
239 zmax = max( qmax, zmax, z( k+1 ) )
240 10 CONTINUE
241 IF( z( 2*n-1 ).LT.zero ) THEN
242 info = -( 200+2*n-1 )
243 CALL xerbla( 'SLASQ2', 2 )
244 RETURN
245 END IF
246 d = d + z( 2*n-1 )
247 qmax = max( qmax, z( 2*n-1 ) )
248 zmax = max( qmax, zmax )
249*
250* Check for diagonality.
251*
252 IF( e.EQ.zero ) THEN
253 DO 20 k = 2, n
254 z( k ) = z( 2*k-1 )
255 20 CONTINUE
256 CALL slasrt( 'D', n, z, iinfo )
257 z( 2*n-1 ) = d
258 RETURN
259 END IF
260*
261 trace = d + e
262*
263* Check for zero data.
264*
265 IF( trace.EQ.zero ) THEN
266 z( 2*n-1 ) = zero
267 RETURN
268 END IF
269*
270* Check whether the machine is IEEE conformable.
271*
272* IEEE = ( ILAENV( 10, 'SLASQ2', 'N', 1, 2, 3, 4 ).EQ.1 )
273*
274* [11/15/2008] The case IEEE=.TRUE. has a problem in single precision with
275* some the test matrices of type 16. The double precision code is fine.
276*
277 ieee = .false.
278*
279* Rearrange data for locality: Z=(q1,qq1,e1,ee1,q2,qq2,e2,ee2,...).
280*
281 DO 30 k = 2*n, 2, -2
282 z( 2*k ) = zero
283 z( 2*k-1 ) = z( k )
284 z( 2*k-2 ) = zero
285 z( 2*k-3 ) = z( k-1 )
286 30 CONTINUE
287*
288 i0 = 1
289 n0 = n
290*
291* Reverse the qd-array, if warranted.
292*
293 IF( cbias*z( 4*i0-3 ).LT.z( 4*n0-3 ) ) THEN
294 ipn4 = 4*( i0+n0 )
295 DO 40 i4 = 4*i0, 2*( i0+n0-1 ), 4
296 temp = z( i4-3 )
297 z( i4-3 ) = z( ipn4-i4-3 )
298 z( ipn4-i4-3 ) = temp
299 temp = z( i4-1 )
300 z( i4-1 ) = z( ipn4-i4-5 )
301 z( ipn4-i4-5 ) = temp
302 40 CONTINUE
303 END IF
304*
305* Initial split checking via dqd and Li's test.
306*
307 pp = 0
308*
309 DO 80 k = 1, 2
310*
311 d = z( 4*n0+pp-3 )
312 DO 50 i4 = 4*( n0-1 ) + pp, 4*i0 + pp, -4
313 IF( z( i4-1 ).LE.tol2*d ) THEN
314 z( i4-1 ) = -zero
315 d = z( i4-3 )
316 ELSE
317 d = z( i4-3 )*( d / ( d+z( i4-1 ) ) )
318 END IF
319 50 CONTINUE
320*
321* dqd maps Z to ZZ plus Li's test.
322*
323 emin = z( 4*i0+pp+1 )
324 d = z( 4*i0+pp-3 )
325 DO 60 i4 = 4*i0 + pp, 4*( n0-1 ) + pp, 4
326 z( i4-2*pp-2 ) = d + z( i4-1 )
327 IF( z( i4-1 ).LE.tol2*d ) THEN
328 z( i4-1 ) = -zero
329 z( i4-2*pp-2 ) = d
330 z( i4-2*pp ) = zero
331 d = z( i4+1 )
332 ELSE IF( safmin*z( i4+1 ).LT.z( i4-2*pp-2 ) .AND.
333 $ safmin*z( i4-2*pp-2 ).LT.z( i4+1 ) ) THEN
334 temp = z( i4+1 ) / z( i4-2*pp-2 )
335 z( i4-2*pp ) = z( i4-1 )*temp
336 d = d*temp
337 ELSE
338 z( i4-2*pp ) = z( i4+1 )*( z( i4-1 ) / z( i4-2*pp-2 ) )
339 d = z( i4+1 )*( d / z( i4-2*pp-2 ) )
340 END IF
341 emin = min( emin, z( i4-2*pp ) )
342 60 CONTINUE
343 z( 4*n0-pp-2 ) = d
344*
345* Now find qmax.
346*
347 qmax = z( 4*i0-pp-2 )
348 DO 70 i4 = 4*i0 - pp + 2, 4*n0 - pp - 2, 4
349 qmax = max( qmax, z( i4 ) )
350 70 CONTINUE
351*
352* Prepare for the next iteration on K.
353*
354 pp = 1 - pp
355 80 CONTINUE
356*
357* Initialise variables to pass to SLASQ3.
358*
359 ttype = 0
360 dmin1 = zero
361 dmin2 = zero
362 dn = zero
363 dn1 = zero
364 dn2 = zero
365 g = zero
366 tau = zero
367*
368 iter = 2
369 nfail = 0
370 ndiv = 2*( n0-i0 )
371*
372 DO 160 iwhila = 1, n + 1
373 IF( n0.LT.1 )
374 $ GO TO 170
375*
376* While array unfinished do
377*
378* E(N0) holds the value of SIGMA when submatrix in I0:N0
379* splits from the rest of the array, but is negated.
380*
381 desig = zero
382 IF( n0.EQ.n ) THEN
383 sigma = zero
384 ELSE
385 sigma = -z( 4*n0-1 )
386 END IF
387 IF( sigma.LT.zero ) THEN
388 info = 1
389 RETURN
390 END IF
391*
392* Find last unreduced submatrix's top index I0, find QMAX and
393* EMIN. Find Gershgorin-type bound if Q's much greater than E's.
394*
395 emax = zero
396 IF( n0.GT.i0 ) THEN
397 emin = abs( z( 4*n0-5 ) )
398 ELSE
399 emin = zero
400 END IF
401 qmin = z( 4*n0-3 )
402 qmax = qmin
403 DO 90 i4 = 4*n0, 8, -4
404 IF( z( i4-5 ).LE.zero )
405 $ GO TO 100
406 IF( qmin.GE.four*emax ) THEN
407 qmin = min( qmin, z( i4-3 ) )
408 emax = max( emax, z( i4-5 ) )
409 END IF
410 qmax = max( qmax, z( i4-7 )+z( i4-5 ) )
411 emin = min( emin, z( i4-5 ) )
412 90 CONTINUE
413 i4 = 4
414*
415 100 CONTINUE
416 i0 = i4 / 4
417 pp = 0
418*
419 IF( n0-i0.GT.1 ) THEN
420 dee = z( 4*i0-3 )
421 deemin = dee
422 kmin = i0
423 DO 110 i4 = 4*i0+1, 4*n0-3, 4
424 dee = z( i4 )*( dee /( dee+z( i4-2 ) ) )
425 IF( dee.LE.deemin ) THEN
426 deemin = dee
427 kmin = ( i4+3 )/4
428 END IF
429 110 CONTINUE
430 IF( (kmin-i0)*2.LT.n0-kmin .AND.
431 $ deemin.LE.half*z(4*n0-3) ) THEN
432 ipn4 = 4*( i0+n0 )
433 pp = 2
434 DO 120 i4 = 4*i0, 2*( i0+n0-1 ), 4
435 temp = z( i4-3 )
436 z( i4-3 ) = z( ipn4-i4-3 )
437 z( ipn4-i4-3 ) = temp
438 temp = z( i4-2 )
439 z( i4-2 ) = z( ipn4-i4-2 )
440 z( ipn4-i4-2 ) = temp
441 temp = z( i4-1 )
442 z( i4-1 ) = z( ipn4-i4-5 )
443 z( ipn4-i4-5 ) = temp
444 temp = z( i4 )
445 z( i4 ) = z( ipn4-i4-4 )
446 z( ipn4-i4-4 ) = temp
447 120 CONTINUE
448 END IF
449 END IF
450*
451* Put -(initial shift) into DMIN.
452*
453 dmin = -max( zero, qmin-two*sqrt( qmin )*sqrt( emax ) )
454*
455* Now I0:N0 is unreduced.
456* PP = 0 for ping, PP = 1 for pong.
457* PP = 2 indicates that flipping was applied to the Z array and
458* and that the tests for deflation upon entry in SLASQ3
459* should not be performed.
460*
461 nbig = 100*( n0-i0+1 )
462 DO 140 iwhilb = 1, nbig
463 IF( i0.GT.n0 )
464 $ GO TO 150
465*
466* While submatrix unfinished take a good dqds step.
467*
468 CALL slasq3( i0, n0, z, pp, dmin, sigma, desig, qmax,
469 $ nfail,
470 $ iter, ndiv, ieee, ttype, dmin1, dmin2, dn, dn1,
471 $ dn2, g, tau )
472*
473 pp = 1 - pp
474*
475* When EMIN is very small check for splits.
476*
477 IF( pp.EQ.0 .AND. n0-i0.GE.3 ) THEN
478 IF( z( 4*n0 ).LE.tol2*qmax .OR.
479 $ z( 4*n0-1 ).LE.tol2*sigma ) THEN
480 splt = i0 - 1
481 qmax = z( 4*i0-3 )
482 emin = z( 4*i0-1 )
483 oldemn = z( 4*i0 )
484 DO 130 i4 = 4*i0, 4*( n0-3 ), 4
485 IF( z( i4 ).LE.tol2*z( i4-3 ) .OR.
486 $ z( i4-1 ).LE.tol2*sigma ) THEN
487 z( i4-1 ) = -sigma
488 splt = i4 / 4
489 qmax = zero
490 emin = z( i4+3 )
491 oldemn = z( i4+4 )
492 ELSE
493 qmax = max( qmax, z( i4+1 ) )
494 emin = min( emin, z( i4-1 ) )
495 oldemn = min( oldemn, z( i4 ) )
496 END IF
497 130 CONTINUE
498 z( 4*n0-1 ) = emin
499 z( 4*n0 ) = oldemn
500 i0 = splt + 1
501 END IF
502 END IF
503*
504 140 CONTINUE
505*
506 info = 2
507*
508* Maximum number of iterations exceeded, restore the shift
509* SIGMA and place the new d's and e's in a qd array.
510* This might need to be done for several blocks
511*
512 i1 = i0
513 n1 = n0
514 145 CONTINUE
515 tempq = z( 4*i0-3 )
516 z( 4*i0-3 ) = z( 4*i0-3 ) + sigma
517 DO k = i0+1, n0
518 tempe = z( 4*k-5 )
519 z( 4*k-5 ) = z( 4*k-5 ) * (tempq / z( 4*k-7 ))
520 tempq = z( 4*k-3 )
521 z( 4*k-3 ) = z( 4*k-3 ) + sigma + tempe - z( 4*k-5 )
522 END DO
523*
524* Prepare to do this on the previous block if there is one
525*
526 IF( i1.GT.1 ) THEN
527 n1 = i1-1
528 DO WHILE( ( i1.GE.2 ) .AND. ( z(4*i1-5).GE.zero ) )
529 i1 = i1 - 1
530 END DO
531 IF( i1.GE.1 ) THEN
532 sigma = -z(4*n1-1)
533 GO TO 145
534 END IF
535 END IF
536
537 DO k = 1, n
538 z( 2*k-1 ) = z( 4*k-3 )
539*
540* Only the block 1..N0 is unfinished. The rest of the e's
541* must be essentially zero, although sometimes other data
542* has been stored in them.
543*
544 IF( k.LT.n0 ) THEN
545 z( 2*k ) = z( 4*k-1 )
546 ELSE
547 z( 2*k ) = 0
548 END IF
549 END DO
550 RETURN
551*
552* end IWHILB
553*
554 150 CONTINUE
555*
556 160 CONTINUE
557*
558 info = 3
559 RETURN
560*
561* end IWHILA
562*
563 170 CONTINUE
564*
565* Move q's to the front.
566*
567 DO 180 k = 2, n
568 z( k ) = z( 4*k-3 )
569 180 CONTINUE
570*
571* Sort and compute sum of eigenvalues.
572*
573 CALL slasrt( 'D', n, z, iinfo )
574*
575 e = zero
576 DO 190 k = n, 1, -1
577 e = e + z( k )
578 190 CONTINUE
579*
580* Store trace, sum(eigenvalues) and information on performance.
581*
582 z( 2*n+1 ) = trace
583 z( 2*n+2 ) = e
584 z( 2*n+3 ) = real( iter )
585 z( 2*n+4 ) = real( ndiv ) / real( n**2 )
586 z( 2*n+5 ) = hundrd*real( nfail / iter )
587 RETURN
588*
589* End of SLASQ2
590*
591 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine slasq2(n, z, info)
SLASQ2 computes all the eigenvalues of the symmetric positive definite tridiagonal matrix associated ...
Definition slasq2.f:110
subroutine slasq3(i0, n0, z, pp, dmin, sigma, desig, qmax, nfail, iter, ndiv, ieee, ttype, dmin1, dmin2, dn, dn1, dn2, g, tau)
SLASQ3 checks for deflation, computes a shift and calls dqds. Used by sbdsqr.
Definition slasq3.f:181
subroutine slasrt(id, n, d, info)
SLASRT sorts numbers in increasing or decreasing order.
Definition slasrt.f:86