LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dlaeda.f
Go to the documentation of this file.
1*> \brief \b DLAEDA used by DSTEDC. Computes the Z vector determining the rank-one modification of the diagonal matrix. Used when the original matrix is dense.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DLAEDA + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlaeda.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlaeda.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlaeda.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DLAEDA( N, TLVLS, CURLVL, CURPBM, PRMPTR, PERM, GIVPTR,
20* GIVCOL, GIVNUM, Q, QPTR, Z, ZTEMP, INFO )
21*
22* .. Scalar Arguments ..
23* INTEGER CURLVL, CURPBM, INFO, N, TLVLS
24* ..
25* .. Array Arguments ..
26* INTEGER GIVCOL( 2, * ), GIVPTR( * ), PERM( * ),
27* $ PRMPTR( * ), QPTR( * )
28* DOUBLE PRECISION GIVNUM( 2, * ), Q( * ), Z( * ), ZTEMP( * )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> DLAEDA computes the Z vector corresponding to the merge step in the
38*> CURLVLth step of the merge process with TLVLS steps for the CURPBMth
39*> problem.
40*> \endverbatim
41*
42* Arguments:
43* ==========
44*
45*> \param[in] N
46*> \verbatim
47*> N is INTEGER
48*> The dimension of the symmetric tridiagonal matrix. N >= 0.
49*> \endverbatim
50*>
51*> \param[in] TLVLS
52*> \verbatim
53*> TLVLS is INTEGER
54*> The total number of merging levels in the overall divide and
55*> conquer tree.
56*> \endverbatim
57*>
58*> \param[in] CURLVL
59*> \verbatim
60*> CURLVL is INTEGER
61*> The current level in the overall merge routine,
62*> 0 <= curlvl <= tlvls.
63*> \endverbatim
64*>
65*> \param[in] CURPBM
66*> \verbatim
67*> CURPBM is INTEGER
68*> The current problem in the current level in the overall
69*> merge routine (counting from upper left to lower right).
70*> \endverbatim
71*>
72*> \param[in] PRMPTR
73*> \verbatim
74*> PRMPTR is INTEGER array, dimension (N lg N)
75*> Contains a list of pointers which indicate where in PERM a
76*> level's permutation is stored. PRMPTR(i+1) - PRMPTR(i)
77*> indicates the size of the permutation and incidentally the
78*> size of the full, non-deflated problem.
79*> \endverbatim
80*>
81*> \param[in] PERM
82*> \verbatim
83*> PERM is INTEGER array, dimension (N lg N)
84*> Contains the permutations (from deflation and sorting) to be
85*> applied to each eigenblock.
86*> \endverbatim
87*>
88*> \param[in] GIVPTR
89*> \verbatim
90*> GIVPTR is INTEGER array, dimension (N lg N)
91*> Contains a list of pointers which indicate where in GIVCOL a
92*> level's Givens rotations are stored. GIVPTR(i+1) - GIVPTR(i)
93*> indicates the number of Givens rotations.
94*> \endverbatim
95*>
96*> \param[in] GIVCOL
97*> \verbatim
98*> GIVCOL is INTEGER array, dimension (2, N lg N)
99*> Each pair of numbers indicates a pair of columns to take place
100*> in a Givens rotation.
101*> \endverbatim
102*>
103*> \param[in] GIVNUM
104*> \verbatim
105*> GIVNUM is DOUBLE PRECISION array, dimension (2, N lg N)
106*> Each number indicates the S value to be used in the
107*> corresponding Givens rotation.
108*> \endverbatim
109*>
110*> \param[in] Q
111*> \verbatim
112*> Q is DOUBLE PRECISION array, dimension (N**2)
113*> Contains the square eigenblocks from previous levels, the
114*> starting positions for blocks are given by QPTR.
115*> \endverbatim
116*>
117*> \param[in] QPTR
118*> \verbatim
119*> QPTR is INTEGER array, dimension (N+2)
120*> Contains a list of pointers which indicate where in Q an
121*> eigenblock is stored. SQRT( QPTR(i+1) - QPTR(i) ) indicates
122*> the size of the block.
123*> \endverbatim
124*>
125*> \param[out] Z
126*> \verbatim
127*> Z is DOUBLE PRECISION array, dimension (N)
128*> On output this vector contains the updating vector (the last
129*> row of the first sub-eigenvector matrix and the first row of
130*> the second sub-eigenvector matrix).
131*> \endverbatim
132*>
133*> \param[out] ZTEMP
134*> \verbatim
135*> ZTEMP is DOUBLE PRECISION array, dimension (N)
136*> \endverbatim
137*>
138*> \param[out] INFO
139*> \verbatim
140*> INFO is INTEGER
141*> = 0: successful exit.
142*> < 0: if INFO = -i, the i-th argument had an illegal value.
143*> \endverbatim
144*
145* Authors:
146* ========
147*
148*> \author Univ. of Tennessee
149*> \author Univ. of California Berkeley
150*> \author Univ. of Colorado Denver
151*> \author NAG Ltd.
152*
153*> \ingroup laeda
154*
155*> \par Contributors:
156* ==================
157*>
158*> Jeff Rutter, Computer Science Division, University of California
159*> at Berkeley, USA
160*
161* =====================================================================
162 SUBROUTINE dlaeda( N, TLVLS, CURLVL, CURPBM, PRMPTR, PERM,
163 $ GIVPTR,
164 $ GIVCOL, GIVNUM, Q, QPTR, Z, ZTEMP, INFO )
165*
166* -- LAPACK computational routine --
167* -- LAPACK is a software package provided by Univ. of Tennessee, --
168* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
169*
170* .. Scalar Arguments ..
171 INTEGER CURLVL, CURPBM, INFO, N, TLVLS
172* ..
173* .. Array Arguments ..
174 INTEGER GIVCOL( 2, * ), GIVPTR( * ), PERM( * ),
175 $ PRMPTR( * ), QPTR( * )
176 DOUBLE PRECISION GIVNUM( 2, * ), Q( * ), Z( * ), ZTEMP( * )
177* ..
178*
179* =====================================================================
180*
181* .. Parameters ..
182 DOUBLE PRECISION ZERO, HALF, ONE
183 PARAMETER ( ZERO = 0.0d0, half = 0.5d0, one = 1.0d0 )
184* ..
185* .. Local Scalars ..
186 INTEGER BSIZ1, BSIZ2, CURR, I, K, MID, PSIZ1, PSIZ2,
187 $ PTR, ZPTR1
188* ..
189* .. External Subroutines ..
190 EXTERNAL dcopy, dgemv, drot, xerbla
191* ..
192* .. Intrinsic Functions ..
193 INTRINSIC dble, int, sqrt
194* ..
195* .. Executable Statements ..
196*
197* Test the input parameters.
198*
199 info = 0
200*
201 IF( n.LT.0 ) THEN
202 info = -1
203 END IF
204 IF( info.NE.0 ) THEN
205 CALL xerbla( 'DLAEDA', -info )
206 RETURN
207 END IF
208*
209* Quick return if possible
210*
211 IF( n.EQ.0 )
212 $ RETURN
213*
214* Determine location of first number in second half.
215*
216 mid = n / 2 + 1
217*
218* Gather last/first rows of appropriate eigenblocks into center of Z
219*
220 ptr = 1
221*
222* Determine location of lowest level subproblem in the full storage
223* scheme
224*
225 curr = ptr + curpbm*2**curlvl + 2**( curlvl-1 ) - 1
226*
227* Determine size of these matrices. We add HALF to the value of
228* the SQRT in case the machine underestimates one of these square
229* roots.
230*
231 bsiz1 = int( half+sqrt( dble( qptr( curr+1 )-qptr( curr ) ) ) )
232 bsiz2 = int( half+sqrt( dble( qptr( curr+2 )-qptr( curr+1 ) ) ) )
233 DO 10 k = 1, mid - bsiz1 - 1
234 z( k ) = zero
235 10 CONTINUE
236 CALL dcopy( bsiz1, q( qptr( curr )+bsiz1-1 ), bsiz1,
237 $ z( mid-bsiz1 ), 1 )
238 CALL dcopy( bsiz2, q( qptr( curr+1 ) ), bsiz2, z( mid ), 1 )
239 DO 20 k = mid + bsiz2, n
240 z( k ) = zero
241 20 CONTINUE
242*
243* Loop through remaining levels 1 -> CURLVL applying the Givens
244* rotations and permutation and then multiplying the center matrices
245* against the current Z.
246*
247 ptr = 2**tlvls + 1
248 DO 70 k = 1, curlvl - 1
249 curr = ptr + curpbm*2**( curlvl-k ) + 2**( curlvl-k-1 ) - 1
250 psiz1 = prmptr( curr+1 ) - prmptr( curr )
251 psiz2 = prmptr( curr+2 ) - prmptr( curr+1 )
252 zptr1 = mid - psiz1
253*
254* Apply Givens at CURR and CURR+1
255*
256 DO 30 i = givptr( curr ), givptr( curr+1 ) - 1
257 CALL drot( 1, z( zptr1+givcol( 1, i )-1 ), 1,
258 $ z( zptr1+givcol( 2, i )-1 ), 1, givnum( 1, i ),
259 $ givnum( 2, i ) )
260 30 CONTINUE
261 DO 40 i = givptr( curr+1 ), givptr( curr+2 ) - 1
262 CALL drot( 1, z( mid-1+givcol( 1, i ) ), 1,
263 $ z( mid-1+givcol( 2, i ) ), 1, givnum( 1, i ),
264 $ givnum( 2, i ) )
265 40 CONTINUE
266 psiz1 = prmptr( curr+1 ) - prmptr( curr )
267 psiz2 = prmptr( curr+2 ) - prmptr( curr+1 )
268 DO 50 i = 0, psiz1 - 1
269 ztemp( i+1 ) = z( zptr1+perm( prmptr( curr )+i )-1 )
270 50 CONTINUE
271 DO 60 i = 0, psiz2 - 1
272 ztemp( psiz1+i+1 ) = z( mid+perm( prmptr( curr+1 )+i )-1 )
273 60 CONTINUE
274*
275* Multiply Blocks at CURR and CURR+1
276*
277* Determine size of these matrices. We add HALF to the value of
278* the SQRT in case the machine underestimates one of these
279* square roots.
280*
281 bsiz1 = int( half+sqrt( dble( qptr( curr+1 )-qptr( curr ) ) ) )
282 bsiz2 = int( half+sqrt( dble( qptr( curr+2 )-qptr( curr+
283 $ 1 ) ) ) )
284 IF( bsiz1.GT.0 ) THEN
285 CALL dgemv( 'T', bsiz1, bsiz1, one, q( qptr( curr ) ),
286 $ bsiz1, ztemp( 1 ), 1, zero, z( zptr1 ), 1 )
287 END IF
288 CALL dcopy( psiz1-bsiz1, ztemp( bsiz1+1 ), 1,
289 $ 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*
305 END
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 dlaeda(n, tlvls, curlvl, curpbm, prmptr, perm, givptr, givcol, givnum, q, qptr, z, ztemp, info)
DLAEDA used by DSTEDC. Computes the Z vector determining the rank-one modification of the diagonal ma...
Definition dlaeda.f:165
subroutine drot(n, dx, incx, dy, incy, c, s)
DROT
Definition drot.f:92