LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zlar1v()

subroutine zlar1v ( integer  n,
integer  b1,
integer  bn,
double precision  lambda,
double precision, dimension( * )  d,
double precision, dimension( * )  l,
double precision, dimension( * )  ld,
double precision, dimension( * )  lld,
double precision  pivmin,
double precision  gaptol,
complex*16, dimension( * )  z,
logical  wantnc,
integer  negcnt,
double precision  ztz,
double precision  mingma,
integer  r,
integer, dimension( * )  isuppz,
double precision  nrminv,
double precision  resid,
double precision  rqcorr,
double precision, dimension( * )  work 
)

ZLAR1V computes the (scaled) r-th column of the inverse of the submatrix in rows b1 through bn of the tridiagonal matrix LDLT - λI.

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

Purpose:
 ZLAR1V computes the (scaled) r-th column of the inverse of
 the sumbmatrix in rows B1 through BN of the tridiagonal matrix
 L D L**T - sigma I. When sigma is close to an eigenvalue, the
 computed vector is an accurate eigenvector. Usually, r corresponds
 to the index where the eigenvector is largest in magnitude.
 The following steps accomplish this computation :
 (a) Stationary qd transform,  L D L**T - sigma I = L(+) D(+) L(+)**T,
 (b) Progressive qd transform, L D L**T - sigma I = U(-) D(-) U(-)**T,
 (c) Computation of the diagonal elements of the inverse of
     L D L**T - sigma I by combining the above transforms, and choosing
     r as the index where the diagonal of the inverse is (one of the)
     largest in magnitude.
 (d) Computation of the (scaled) r-th column of the inverse using the
     twisted factorization obtained by combining the top part of the
     the stationary and the bottom part of the progressive transform.
Parameters
[in]N
          N is INTEGER
           The order of the matrix L D L**T.
[in]B1
          B1 is INTEGER
           First index of the submatrix of L D L**T.
[in]BN
          BN is INTEGER
           Last index of the submatrix of L D L**T.
[in]LAMBDA
          LAMBDA is DOUBLE PRECISION
           The shift. In order to compute an accurate eigenvector,
           LAMBDA should be a good approximation to an eigenvalue
           of L D L**T.
[in]L
          L is DOUBLE PRECISION array, dimension (N-1)
           The (n-1) subdiagonal elements of the unit bidiagonal matrix
           L, in elements 1 to N-1.
[in]D
          D is DOUBLE PRECISION array, dimension (N)
           The n diagonal elements of the diagonal matrix D.
[in]LD
          LD is DOUBLE PRECISION array, dimension (N-1)
           The n-1 elements L(i)*D(i).
[in]LLD
          LLD is DOUBLE PRECISION array, dimension (N-1)
           The n-1 elements L(i)*L(i)*D(i).
[in]PIVMIN
          PIVMIN is DOUBLE PRECISION
           The minimum pivot in the Sturm sequence.
[in]GAPTOL
          GAPTOL is DOUBLE PRECISION
           Tolerance that indicates when eigenvector entries are negligible
           w.r.t. their contribution to the residual.
[in,out]Z
          Z is COMPLEX*16 array, dimension (N)
           On input, all entries of Z must be set to 0.
           On output, Z contains the (scaled) r-th column of the
           inverse. The scaling is such that Z(R) equals 1.
[in]WANTNC
          WANTNC is LOGICAL
           Specifies whether NEGCNT has to be computed.
[out]NEGCNT
          NEGCNT is INTEGER
           If WANTNC is .TRUE. then NEGCNT = the number of pivots < pivmin
           in the  matrix factorization L D L**T, and NEGCNT = -1 otherwise.
[out]ZTZ
          ZTZ is DOUBLE PRECISION
           The square of the 2-norm of Z.
[out]MINGMA
          MINGMA is DOUBLE PRECISION
           The reciprocal of the largest (in magnitude) diagonal
           element of the inverse of L D L**T - sigma I.
[in,out]R
          R is INTEGER
           The twist index for the twisted factorization used to
           compute Z.
           On input, 0 <= R <= N. If R is input as 0, R is set to
           the index where (L D L**T - sigma I)^{-1} is largest
           in magnitude. If 1 <= R <= N, R is unchanged.
           On output, R contains the twist index used to compute Z.
           Ideally, R designates the position of the maximum entry in the
           eigenvector.
[out]ISUPPZ
          ISUPPZ is INTEGER array, dimension (2)
           The support of the vector in Z, i.e., the vector Z is
           nonzero only in elements ISUPPZ(1) through ISUPPZ( 2 ).
[out]NRMINV
          NRMINV is DOUBLE PRECISION
           NRMINV = 1/SQRT( ZTZ )
[out]RESID
          RESID is DOUBLE PRECISION
           The residual of the FP vector.
           RESID = ABS( MINGMA )/SQRT( ZTZ )
[out]RQCORR
          RQCORR is DOUBLE PRECISION
           The Rayleigh Quotient correction to LAMBDA.
           RQCORR = MINGMA*TMP
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (4*N)
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
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 227 of file zlar1v.f.

230*
231* -- LAPACK auxiliary routine --
232* -- LAPACK is a software package provided by Univ. of Tennessee, --
233* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
234*
235* .. Scalar Arguments ..
236 LOGICAL WANTNC
237 INTEGER B1, BN, N, NEGCNT, R
238 DOUBLE PRECISION GAPTOL, LAMBDA, MINGMA, NRMINV, PIVMIN, RESID,
239 $ RQCORR, ZTZ
240* ..
241* .. Array Arguments ..
242 INTEGER ISUPPZ( * )
243 DOUBLE PRECISION D( * ), L( * ), LD( * ), LLD( * ),
244 $ WORK( * )
245 COMPLEX*16 Z( * )
246* ..
247*
248* =====================================================================
249*
250* .. Parameters ..
251 DOUBLE PRECISION ZERO, ONE
252 parameter( zero = 0.0d0, one = 1.0d0 )
253 COMPLEX*16 CONE
254 parameter( cone = ( 1.0d0, 0.0d0 ) )
255
256* ..
257* .. Local Scalars ..
258 LOGICAL SAWNAN1, SAWNAN2
259 INTEGER I, INDLPL, INDP, INDS, INDUMN, NEG1, NEG2, R1,
260 $ R2
261 DOUBLE PRECISION DMINUS, DPLUS, EPS, S, TMP
262* ..
263* .. External Functions ..
264 LOGICAL DISNAN
265 DOUBLE PRECISION DLAMCH
266 EXTERNAL disnan, dlamch
267* ..
268* .. Intrinsic Functions ..
269 INTRINSIC abs, dble
270* ..
271* .. Executable Statements ..
272*
273 eps = dlamch( 'Precision' )
274
275
276 IF( r.EQ.0 ) THEN
277 r1 = b1
278 r2 = bn
279 ELSE
280 r1 = r
281 r2 = r
282 END IF
283
284* Storage for LPLUS
285 indlpl = 0
286* Storage for UMINUS
287 indumn = n
288 inds = 2*n + 1
289 indp = 3*n + 1
290
291 IF( b1.EQ.1 ) THEN
292 work( inds ) = zero
293 ELSE
294 work( inds+b1-1 ) = lld( b1-1 )
295 END IF
296
297*
298* Compute the stationary transform (using the differential form)
299* until the index R2.
300*
301 sawnan1 = .false.
302 neg1 = 0
303 s = work( inds+b1-1 ) - lambda
304 DO 50 i = b1, r1 - 1
305 dplus = d( i ) + s
306 work( indlpl+i ) = ld( i ) / dplus
307 IF(dplus.LT.zero) neg1 = neg1 + 1
308 work( inds+i ) = s*work( indlpl+i )*l( i )
309 s = work( inds+i ) - lambda
310 50 CONTINUE
311 sawnan1 = disnan( s )
312 IF( sawnan1 ) GOTO 60
313 DO 51 i = r1, r2 - 1
314 dplus = d( i ) + s
315 work( indlpl+i ) = ld( i ) / dplus
316 work( inds+i ) = s*work( indlpl+i )*l( i )
317 s = work( inds+i ) - lambda
318 51 CONTINUE
319 sawnan1 = disnan( s )
320*
321 60 CONTINUE
322 IF( sawnan1 ) THEN
323* Runs a slower version of the above loop if a NaN is detected
324 neg1 = 0
325 s = work( inds+b1-1 ) - lambda
326 DO 70 i = b1, r1 - 1
327 dplus = d( i ) + s
328 IF(abs(dplus).LT.pivmin) dplus = -pivmin
329 work( indlpl+i ) = ld( i ) / dplus
330 IF(dplus.LT.zero) neg1 = neg1 + 1
331 work( inds+i ) = s*work( indlpl+i )*l( i )
332 IF( work( indlpl+i ).EQ.zero )
333 $ work( inds+i ) = lld( i )
334 s = work( inds+i ) - lambda
335 70 CONTINUE
336 DO 71 i = r1, r2 - 1
337 dplus = d( i ) + s
338 IF(abs(dplus).LT.pivmin) dplus = -pivmin
339 work( indlpl+i ) = ld( i ) / dplus
340 work( inds+i ) = s*work( indlpl+i )*l( i )
341 IF( work( indlpl+i ).EQ.zero )
342 $ work( inds+i ) = lld( i )
343 s = work( inds+i ) - lambda
344 71 CONTINUE
345 END IF
346*
347* Compute the progressive transform (using the differential form)
348* until the index R1
349*
350 sawnan2 = .false.
351 neg2 = 0
352 work( indp+bn-1 ) = d( bn ) - lambda
353 DO 80 i = bn - 1, r1, -1
354 dminus = lld( i ) + work( indp+i )
355 tmp = d( i ) / dminus
356 IF(dminus.LT.zero) neg2 = neg2 + 1
357 work( indumn+i ) = l( i )*tmp
358 work( indp+i-1 ) = work( indp+i )*tmp - lambda
359 80 CONTINUE
360 tmp = work( indp+r1-1 )
361 sawnan2 = disnan( tmp )
362
363 IF( sawnan2 ) THEN
364* Runs a slower version of the above loop if a NaN is detected
365 neg2 = 0
366 DO 100 i = bn-1, r1, -1
367 dminus = lld( i ) + work( indp+i )
368 IF(abs(dminus).LT.pivmin) dminus = -pivmin
369 tmp = d( i ) / dminus
370 IF(dminus.LT.zero) neg2 = neg2 + 1
371 work( indumn+i ) = l( i )*tmp
372 work( indp+i-1 ) = work( indp+i )*tmp - lambda
373 IF( tmp.EQ.zero )
374 $ work( indp+i-1 ) = d( i ) - lambda
375 100 CONTINUE
376 END IF
377*
378* Find the index (from R1 to R2) of the largest (in magnitude)
379* diagonal element of the inverse
380*
381 mingma = work( inds+r1-1 ) + work( indp+r1-1 )
382 IF( mingma.LT.zero ) neg1 = neg1 + 1
383 IF( wantnc ) THEN
384 negcnt = neg1 + neg2
385 ELSE
386 negcnt = -1
387 ENDIF
388 IF( abs(mingma).EQ.zero )
389 $ mingma = eps*work( inds+r1-1 )
390 r = r1
391 DO 110 i = r1, r2 - 1
392 tmp = work( inds+i ) + work( indp+i )
393 IF( tmp.EQ.zero )
394 $ tmp = eps*work( inds+i )
395 IF( abs( tmp ).LE.abs( mingma ) ) THEN
396 mingma = tmp
397 r = i + 1
398 END IF
399 110 CONTINUE
400*
401* Compute the FP vector: solve N^T v = e_r
402*
403 isuppz( 1 ) = b1
404 isuppz( 2 ) = bn
405 z( r ) = cone
406 ztz = one
407*
408* Compute the FP vector upwards from R
409*
410 IF( .NOT.sawnan1 .AND. .NOT.sawnan2 ) THEN
411 DO 210 i = r-1, b1, -1
412 z( i ) = -( work( indlpl+i )*z( i+1 ) )
413 IF( (abs(z(i))+abs(z(i+1)))* abs(ld(i)).LT.gaptol )
414 $ THEN
415 z( i ) = zero
416 isuppz( 1 ) = i + 1
417 GOTO 220
418 ENDIF
419 ztz = ztz + dble( z( i )*z( i ) )
420 210 CONTINUE
421 220 CONTINUE
422 ELSE
423* Run slower loop if NaN occurred.
424 DO 230 i = r - 1, b1, -1
425 IF( z( i+1 ).EQ.zero ) THEN
426 z( i ) = -( ld( i+1 ) / ld( i ) )*z( i+2 )
427 ELSE
428 z( i ) = -( work( indlpl+i )*z( i+1 ) )
429 END IF
430 IF( (abs(z(i))+abs(z(i+1)))* abs(ld(i)).LT.gaptol )
431 $ THEN
432 z( i ) = zero
433 isuppz( 1 ) = i + 1
434 GO TO 240
435 END IF
436 ztz = ztz + dble( z( i )*z( i ) )
437 230 CONTINUE
438 240 CONTINUE
439 ENDIF
440
441* Compute the FP vector downwards from R in blocks of size BLKSIZ
442 IF( .NOT.sawnan1 .AND. .NOT.sawnan2 ) THEN
443 DO 250 i = r, bn-1
444 z( i+1 ) = -( work( indumn+i )*z( i ) )
445 IF( (abs(z(i))+abs(z(i+1)))* abs(ld(i)).LT.gaptol )
446 $ THEN
447 z( i+1 ) = zero
448 isuppz( 2 ) = i
449 GO TO 260
450 END IF
451 ztz = ztz + dble( z( i+1 )*z( i+1 ) )
452 250 CONTINUE
453 260 CONTINUE
454 ELSE
455* Run slower loop if NaN occurred.
456 DO 270 i = r, bn - 1
457 IF( z( i ).EQ.zero ) THEN
458 z( i+1 ) = -( ld( i-1 ) / ld( i ) )*z( i-1 )
459 ELSE
460 z( i+1 ) = -( work( indumn+i )*z( i ) )
461 END IF
462 IF( (abs(z(i))+abs(z(i+1)))* abs(ld(i)).LT.gaptol )
463 $ THEN
464 z( i+1 ) = zero
465 isuppz( 2 ) = i
466 GO TO 280
467 END IF
468 ztz = ztz + dble( z( i+1 )*z( i+1 ) )
469 270 CONTINUE
470 280 CONTINUE
471 END IF
472*
473* Compute quantities for convergence test
474*
475 tmp = one / ztz
476 nrminv = sqrt( tmp )
477 resid = abs( mingma )*nrminv
478 rqcorr = mingma*tmp
479*
480*
481 RETURN
482*
483* End of ZLAR1V
484*
logical function disnan(din)
DISNAN tests input for NaN.
Definition disnan.f:59
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
Here is the caller graph for this function: