LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dlarrb ( integer  N,
double precision, dimension( * )  D,
double precision, dimension( * )  LLD,
integer  IFIRST,
integer  ILAST,
double precision  RTOL1,
double precision  RTOL2,
integer  OFFSET,
double precision, dimension( * )  W,
double precision, dimension( * )  WGAP,
double precision, dimension( * )  WERR,
double precision, dimension( * )  WORK,
integer, dimension( * )  IWORK,
double precision  PIVMIN,
double precision  SPDIAM,
integer  TWIST,
integer  INFO 
)

DLARRB provides limited bisection to locate eigenvalues for more accuracy.

Download DLARRB + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 Given the relatively robust representation(RRR) L D L^T, DLARRB
 does "limited" bisection to refine the eigenvalues of L D L^T,
 W( IFIRST-OFFSET ) through W( ILAST-OFFSET ), to more accuracy. Initial
 guesses for these eigenvalues are input in W, the corresponding estimate
 of the error in these guesses and their gaps are input in WERR
 and WGAP, respectively. During bisection, intervals
 [left, right] are maintained by storing their mid-points and
 semi-widths in the arrays W and WERR respectively.
Parameters
[in]N
          N is INTEGER
          The order of the matrix.
[in]D
          D is DOUBLE PRECISION array, dimension (N)
          The N diagonal elements of the diagonal matrix D.
[in]LLD
          LLD is DOUBLE PRECISION array, dimension (N-1)
          The (N-1) elements L(i)*L(i)*D(i).
[in]IFIRST
          IFIRST is INTEGER
          The index of the first eigenvalue to be computed.
[in]ILAST
          ILAST is INTEGER
          The index of the last eigenvalue to be computed.
[in]RTOL1
          RTOL1 is DOUBLE PRECISION
[in]RTOL2
          RTOL2 is DOUBLE PRECISION
          Tolerance for the convergence of the bisection intervals.
          An interval [LEFT,RIGHT] has converged if
          RIGHT-LEFT.LT.MAX( RTOL1*GAP, RTOL2*MAX(|LEFT|,|RIGHT|) )
          where GAP is the (estimated) distance to the nearest
          eigenvalue.
[in]OFFSET
          OFFSET is INTEGER
          Offset for the arrays W, WGAP and WERR, i.e., the IFIRST-OFFSET
          through ILAST-OFFSET elements of these arrays are to be used.
[in,out]W
          W is DOUBLE PRECISION array, dimension (N)
          On input, W( IFIRST-OFFSET ) through W( ILAST-OFFSET ) are
          estimates of the eigenvalues of L D L^T indexed IFIRST throug
          ILAST.
          On output, these estimates are refined.
[in,out]WGAP
          WGAP is DOUBLE PRECISION array, dimension (N-1)
          On input, the (estimated) gaps between consecutive
          eigenvalues of L D L^T, i.e., WGAP(I-OFFSET) is the gap between
          eigenvalues I and I+1. Note that if IFIRST.EQ.ILAST
          then WGAP(IFIRST-OFFSET) must be set to ZERO.
          On output, these gaps are refined.
[in,out]WERR
          WERR is DOUBLE PRECISION array, dimension (N)
          On input, WERR( IFIRST-OFFSET ) through WERR( ILAST-OFFSET ) are
          the errors in the estimates of the corresponding elements in W.
          On output, these errors are refined.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (2*N)
          Workspace.
[out]IWORK
          IWORK is INTEGER array, dimension (2*N)
          Workspace.
[in]PIVMIN
          PIVMIN is DOUBLE PRECISION
          The minimum pivot in the Sturm sequence.
[in]SPDIAM
          SPDIAM is DOUBLE PRECISION
          The spectral diameter of the matrix.
[in]TWIST
          TWIST is INTEGER
          The twist index for the twisted factorization that is used
          for the negcount.
          TWIST = N: Compute negcount from L D L^T - LAMBDA I = L+ D+ L+^T
          TWIST = 1: Compute negcount from L D L^T - LAMBDA I = U- D- U-^T
          TWIST = R: Compute negcount from L D L^T - LAMBDA I = N(r) D(r) N(r)
[out]INFO
          INFO is INTEGER
          Error flag.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012
Contributors:
Beresford Parlett, University of California, Berkeley, USA
Jim Demmel, University of California, Berkeley, USA
Inderjit Dhillon, University of Texas, Austin, USA
Osni Marques, LBNL/NERSC, USA
Christof Voemel, University of California, Berkeley, USA

Definition at line 198 of file dlarrb.f.

198 *
199 * -- LAPACK auxiliary routine (version 3.4.2) --
200 * -- LAPACK is a software package provided by Univ. of Tennessee, --
201 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
202 * September 2012
203 *
204 * .. Scalar Arguments ..
205  INTEGER ifirst, ilast, info, n, offset, twist
206  DOUBLE PRECISION pivmin, rtol1, rtol2, spdiam
207 * ..
208 * .. Array Arguments ..
209  INTEGER iwork( * )
210  DOUBLE PRECISION d( * ), lld( * ), w( * ),
211  $ werr( * ), wgap( * ), work( * )
212 * ..
213 *
214 * =====================================================================
215 *
216 * .. Parameters ..
217  DOUBLE PRECISION zero, two, half
218  parameter ( zero = 0.0d0, two = 2.0d0,
219  $ half = 0.5d0 )
220  INTEGER maxitr
221 * ..
222 * .. Local Scalars ..
223  INTEGER i, i1, ii, ip, iter, k, negcnt, next, nint,
224  $ olnint, prev, r
225  DOUBLE PRECISION back, cvrgd, gap, left, lgap, mid, mnwdth,
226  $ rgap, right, tmp, width
227 * ..
228 * .. External Functions ..
229  INTEGER dlaneg
230  EXTERNAL dlaneg
231 *
232 * ..
233 * .. Intrinsic Functions ..
234  INTRINSIC abs, max, min
235 * ..
236 * .. Executable Statements ..
237 *
238  info = 0
239 *
240  maxitr = int( ( log( spdiam+pivmin )-log( pivmin ) ) /
241  $ log( two ) ) + 2
242  mnwdth = two * pivmin
243 *
244  r = twist
245  IF((r.LT.1).OR.(r.GT.n)) r = n
246 *
247 * Initialize unconverged intervals in [ WORK(2*I-1), WORK(2*I) ].
248 * The Sturm Count, Count( WORK(2*I-1) ) is arranged to be I-1, while
249 * Count( WORK(2*I) ) is stored in IWORK( 2*I ). The integer IWORK( 2*I-1 )
250 * for an unconverged interval is set to the index of the next unconverged
251 * interval, and is -1 or 0 for a converged interval. Thus a linked
252 * list of unconverged intervals is set up.
253 *
254  i1 = ifirst
255 * The number of unconverged intervals
256  nint = 0
257 * The last unconverged interval found
258  prev = 0
259 
260  rgap = wgap( i1-offset )
261  DO 75 i = i1, ilast
262  k = 2*i
263  ii = i - offset
264  left = w( ii ) - werr( ii )
265  right = w( ii ) + werr( ii )
266  lgap = rgap
267  rgap = wgap( ii )
268  gap = min( lgap, rgap )
269 
270 * Make sure that [LEFT,RIGHT] contains the desired eigenvalue
271 * Compute negcount from dstqds facto L+D+L+^T = L D L^T - LEFT
272 *
273 * Do while( NEGCNT(LEFT).GT.I-1 )
274 *
275  back = werr( ii )
276  20 CONTINUE
277  negcnt = dlaneg( n, d, lld, left, pivmin, r )
278  IF( negcnt.GT.i-1 ) THEN
279  left = left - back
280  back = two*back
281  GO TO 20
282  END IF
283 *
284 * Do while( NEGCNT(RIGHT).LT.I )
285 * Compute negcount from dstqds facto L+D+L+^T = L D L^T - RIGHT
286 *
287  back = werr( ii )
288  50 CONTINUE
289 
290  negcnt = dlaneg( n, d, lld, right, pivmin, r )
291  IF( negcnt.LT.i ) THEN
292  right = right + back
293  back = two*back
294  GO TO 50
295  END IF
296  width = half*abs( left - right )
297  tmp = max( abs( left ), abs( right ) )
298  cvrgd = max(rtol1*gap,rtol2*tmp)
299  IF( width.LE.cvrgd .OR. width.LE.mnwdth ) THEN
300 * This interval has already converged and does not need refinement.
301 * (Note that the gaps might change through refining the
302 * eigenvalues, however, they can only get bigger.)
303 * Remove it from the list.
304  iwork( k-1 ) = -1
305 * Make sure that I1 always points to the first unconverged interval
306  IF((i.EQ.i1).AND.(i.LT.ilast)) i1 = i + 1
307  IF((prev.GE.i1).AND.(i.LE.ilast)) iwork( 2*prev-1 ) = i + 1
308  ELSE
309 * unconverged interval found
310  prev = i
311  nint = nint + 1
312  iwork( k-1 ) = i + 1
313  iwork( k ) = negcnt
314  END IF
315  work( k-1 ) = left
316  work( k ) = right
317  75 CONTINUE
318 
319 *
320 * Do while( NINT.GT.0 ), i.e. there are still unconverged intervals
321 * and while (ITER.LT.MAXITR)
322 *
323  iter = 0
324  80 CONTINUE
325  prev = i1 - 1
326  i = i1
327  olnint = nint
328 
329  DO 100 ip = 1, olnint
330  k = 2*i
331  ii = i - offset
332  rgap = wgap( ii )
333  lgap = rgap
334  IF(ii.GT.1) lgap = wgap( ii-1 )
335  gap = min( lgap, rgap )
336  next = iwork( k-1 )
337  left = work( k-1 )
338  right = work( k )
339  mid = half*( left + right )
340 
341 * semiwidth of interval
342  width = right - mid
343  tmp = max( abs( left ), abs( right ) )
344  cvrgd = max(rtol1*gap,rtol2*tmp)
345  IF( ( width.LE.cvrgd ) .OR. ( width.LE.mnwdth ).OR.
346  $ ( iter.EQ.maxitr ) )THEN
347 * reduce number of unconverged intervals
348  nint = nint - 1
349 * Mark interval as converged.
350  iwork( k-1 ) = 0
351  IF( i1.EQ.i ) THEN
352  i1 = next
353  ELSE
354 * Prev holds the last unconverged interval previously examined
355  IF(prev.GE.i1) iwork( 2*prev-1 ) = next
356  END IF
357  i = next
358  GO TO 100
359  END IF
360  prev = i
361 *
362 * Perform one bisection step
363 *
364  negcnt = dlaneg( n, d, lld, mid, pivmin, r )
365  IF( negcnt.LE.i-1 ) THEN
366  work( k-1 ) = mid
367  ELSE
368  work( k ) = mid
369  END IF
370  i = next
371  100 CONTINUE
372  iter = iter + 1
373 * do another loop if there are still unconverged intervals
374 * However, in the last iteration, all intervals are accepted
375 * since this is the best we can do.
376  IF( ( nint.GT.0 ).AND.(iter.LE.maxitr) ) GO TO 80
377 *
378 *
379 * At this point, all the intervals have converged
380  DO 110 i = ifirst, ilast
381  k = 2*i
382  ii = i - offset
383 * All intervals marked by '0' have been refined.
384  IF( iwork( k-1 ).EQ.0 ) THEN
385  w( ii ) = half*( work( k-1 )+work( k ) )
386  werr( ii ) = work( k ) - w( ii )
387  END IF
388  110 CONTINUE
389 *
390  DO 111 i = ifirst+1, ilast
391  k = 2*i
392  ii = i - offset
393  wgap( ii-1 ) = max( zero,
394  $ w(ii) - werr(ii) - w( ii-1 ) - werr( ii-1 ))
395  111 CONTINUE
396 
397  RETURN
398 *
399 * End of DLARRB
400 *
integer function dlaneg(N, D, LLD, SIGMA, PIVMIN, R)
DLANEG computes the Sturm count.
Definition: dlaneg.f:120

Here is the caller graph for this function: