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

◆ dlaeda()

subroutine dlaeda ( integer  n,
integer  tlvls,
integer  curlvl,
integer  curpbm,
integer, dimension( * )  prmptr,
integer, dimension( * )  perm,
integer, dimension( * )  givptr,
integer, dimension( 2, * )  givcol,
double precision, dimension( 2, * )  givnum,
double precision, dimension( * )  q,
integer, dimension( * )  qptr,
double precision, dimension( * )  z,
double precision, dimension( * )  ztemp,
integer  info 
)

DLAEDA used by DSTEDC. Computes the Z vector determining the rank-one modification of the diagonal matrix. Used when the original matrix is dense.

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

Purpose:
 DLAEDA computes the Z vector corresponding to the merge step in the
 CURLVLth step of the merge process with TLVLS steps for the CURPBMth
 problem.
Parameters
[in]N
          N is INTEGER
         The dimension of the symmetric tridiagonal matrix.  N >= 0.
[in]TLVLS
          TLVLS is INTEGER
         The total number of merging levels in the overall divide and
         conquer tree.
[in]CURLVL
          CURLVL is INTEGER
         The current level in the overall merge routine,
         0 <= curlvl <= tlvls.
[in]CURPBM
          CURPBM is INTEGER
         The current problem in the current level in the overall
         merge routine (counting from upper left to lower right).
[in]PRMPTR
          PRMPTR is INTEGER array, dimension (N lg N)
         Contains a list of pointers which indicate where in PERM a
         level's permutation is stored.  PRMPTR(i+1) - PRMPTR(i)
         indicates the size of the permutation and incidentally the
         size of the full, non-deflated problem.
[in]PERM
          PERM is INTEGER array, dimension (N lg N)
         Contains the permutations (from deflation and sorting) to be
         applied to each eigenblock.
[in]GIVPTR
          GIVPTR is INTEGER array, dimension (N lg N)
         Contains a list of pointers which indicate where in GIVCOL a
         level's Givens rotations are stored.  GIVPTR(i+1) - GIVPTR(i)
         indicates the number of Givens rotations.
[in]GIVCOL
          GIVCOL is INTEGER array, dimension (2, N lg N)
         Each pair of numbers indicates a pair of columns to take place
         in a Givens rotation.
[in]GIVNUM
          GIVNUM is DOUBLE PRECISION array, dimension (2, N lg N)
         Each number indicates the S value to be used in the
         corresponding Givens rotation.
[in]Q
          Q is DOUBLE PRECISION array, dimension (N**2)
         Contains the square eigenblocks from previous levels, the
         starting positions for blocks are given by QPTR.
[in]QPTR
          QPTR is INTEGER array, dimension (N+2)
         Contains a list of pointers which indicate where in Q an
         eigenblock is stored.  SQRT( QPTR(i+1) - QPTR(i) ) indicates
         the size of the block.
[out]Z
          Z is DOUBLE PRECISION array, dimension (N)
         On output this vector contains the updating vector (the last
         row of the first sub-eigenvector matrix and the first row of
         the second sub-eigenvector matrix).
[out]ZTEMP
          ZTEMP is DOUBLE PRECISION array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Definition at line 164 of file dlaeda.f.

166*
167* -- LAPACK computational routine --
168* -- LAPACK is a software package provided by Univ. of Tennessee, --
169* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
170*
171* .. Scalar Arguments ..
172 INTEGER CURLVL, CURPBM, INFO, N, TLVLS
173* ..
174* .. Array Arguments ..
175 INTEGER GIVCOL( 2, * ), GIVPTR( * ), PERM( * ),
176 $ PRMPTR( * ), QPTR( * )
177 DOUBLE PRECISION GIVNUM( 2, * ), Q( * ), Z( * ), ZTEMP( * )
178* ..
179*
180* =====================================================================
181*
182* .. Parameters ..
183 DOUBLE PRECISION ZERO, HALF, ONE
184 parameter( zero = 0.0d0, half = 0.5d0, one = 1.0d0 )
185* ..
186* .. Local Scalars ..
187 INTEGER BSIZ1, BSIZ2, CURR, I, K, MID, PSIZ1, PSIZ2,
188 $ PTR, ZPTR1
189* ..
190* .. External Subroutines ..
191 EXTERNAL dcopy, dgemv, drot, xerbla
192* ..
193* .. Intrinsic Functions ..
194 INTRINSIC dble, int, sqrt
195* ..
196* .. Executable Statements ..
197*
198* Test the input parameters.
199*
200 info = 0
201*
202 IF( n.LT.0 ) THEN
203 info = -1
204 END IF
205 IF( info.NE.0 ) THEN
206 CALL xerbla( 'DLAEDA', -info )
207 RETURN
208 END IF
209*
210* Quick return if possible
211*
212 IF( n.EQ.0 )
213 $ RETURN
214*
215* Determine location of first number in second half.
216*
217 mid = n / 2 + 1
218*
219* Gather last/first rows of appropriate eigenblocks into center of Z
220*
221 ptr = 1
222*
223* Determine location of lowest level subproblem in the full storage
224* scheme
225*
226 curr = ptr + curpbm*2**curlvl + 2**( curlvl-1 ) - 1
227*
228* Determine size of these matrices. We add HALF to the value of
229* the SQRT in case the machine underestimates one of these square
230* roots.
231*
232 bsiz1 = int( half+sqrt( dble( qptr( curr+1 )-qptr( curr ) ) ) )
233 bsiz2 = int( half+sqrt( dble( qptr( curr+2 )-qptr( curr+1 ) ) ) )
234 DO 10 k = 1, mid - bsiz1 - 1
235 z( k ) = zero
236 10 CONTINUE
237 CALL dcopy( bsiz1, q( qptr( curr )+bsiz1-1 ), bsiz1,
238 $ z( mid-bsiz1 ), 1 )
239 CALL dcopy( bsiz2, q( qptr( curr+1 ) ), bsiz2, z( mid ), 1 )
240 DO 20 k = mid + bsiz2, n
241 z( k ) = zero
242 20 CONTINUE
243*
244* Loop through remaining levels 1 -> CURLVL applying the Givens
245* rotations and permutation and then multiplying the center matrices
246* against the current Z.
247*
248 ptr = 2**tlvls + 1
249 DO 70 k = 1, curlvl - 1
250 curr = ptr + curpbm*2**( curlvl-k ) + 2**( curlvl-k-1 ) - 1
251 psiz1 = prmptr( curr+1 ) - prmptr( curr )
252 psiz2 = prmptr( curr+2 ) - prmptr( curr+1 )
253 zptr1 = mid - psiz1
254*
255* Apply Givens at CURR and CURR+1
256*
257 DO 30 i = givptr( curr ), givptr( curr+1 ) - 1
258 CALL drot( 1, z( zptr1+givcol( 1, i )-1 ), 1,
259 $ z( zptr1+givcol( 2, i )-1 ), 1, givnum( 1, i ),
260 $ givnum( 2, i ) )
261 30 CONTINUE
262 DO 40 i = givptr( curr+1 ), givptr( curr+2 ) - 1
263 CALL drot( 1, z( mid-1+givcol( 1, i ) ), 1,
264 $ z( mid-1+givcol( 2, i ) ), 1, givnum( 1, i ),
265 $ givnum( 2, i ) )
266 40 CONTINUE
267 psiz1 = prmptr( curr+1 ) - prmptr( curr )
268 psiz2 = prmptr( curr+2 ) - prmptr( curr+1 )
269 DO 50 i = 0, psiz1 - 1
270 ztemp( i+1 ) = z( zptr1+perm( prmptr( curr )+i )-1 )
271 50 CONTINUE
272 DO 60 i = 0, psiz2 - 1
273 ztemp( psiz1+i+1 ) = z( mid+perm( prmptr( curr+1 )+i )-1 )
274 60 CONTINUE
275*
276* Multiply Blocks at CURR and CURR+1
277*
278* Determine size of these matrices. We add HALF to the value of
279* the SQRT in case the machine underestimates one of these
280* square roots.
281*
282 bsiz1 = int( half+sqrt( dble( qptr( curr+1 )-qptr( curr ) ) ) )
283 bsiz2 = int( half+sqrt( dble( qptr( curr+2 )-qptr( curr+
284 $ 1 ) ) ) )
285 IF( bsiz1.GT.0 ) THEN
286 CALL dgemv( 'T', bsiz1, bsiz1, one, q( qptr( curr ) ),
287 $ bsiz1, ztemp( 1 ), 1, zero, z( zptr1 ), 1 )
288 END IF
289 CALL dcopy( psiz1-bsiz1, ztemp( bsiz1+1 ), 1, z( zptr1+bsiz1 ),
290 $ 1 )
291 IF( bsiz2.GT.0 ) THEN
292 CALL dgemv( 'T', bsiz2, bsiz2, one, q( qptr( curr+1 ) ),
293 $ bsiz2, ztemp( psiz1+1 ), 1, zero, z( mid ), 1 )
294 END IF
295 CALL dcopy( psiz2-bsiz2, ztemp( psiz1+bsiz2+1 ), 1,
296 $ z( mid+bsiz2 ), 1 )
297*
298 ptr = ptr + 2**( tlvls-k )
299 70 CONTINUE
300*
301 RETURN
302*
303* End of DLAEDA
304*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
Definition dgemv.f:158
subroutine drot(n, dx, incx, dy, incy, c, s)
DROT
Definition drot.f:92
Here is the call graph for this function:
Here is the caller graph for this function: