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

◆ claed0()

subroutine claed0 ( integer  qsiz,
integer  n,
real, dimension( * )  d,
real, dimension( * )  e,
complex, dimension( ldq, * )  q,
integer  ldq,
complex, dimension( ldqs, * )  qstore,
integer  ldqs,
real, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  info 
)

CLAED0 used by CSTEDC. Computes all eigenvalues and corresponding eigenvectors of an unreduced symmetric tridiagonal matrix using the divide and conquer method.

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

Purpose:
 Using the divide and conquer method, CLAED0 computes all eigenvalues
 of a symmetric tridiagonal matrix which is one diagonal block of
 those from reducing a dense or band Hermitian matrix and
 corresponding eigenvectors of the dense or band matrix.
Parameters
[in]QSIZ
          QSIZ is INTEGER
         The dimension of the unitary matrix used to reduce
         the full matrix to tridiagonal form.  QSIZ >= N if ICOMPQ = 1.
[in]N
          N is INTEGER
         The dimension of the symmetric tridiagonal matrix.  N >= 0.
[in,out]D
          D is REAL array, dimension (N)
         On entry, the diagonal elements of the tridiagonal matrix.
         On exit, the eigenvalues in ascending order.
[in,out]E
          E is REAL array, dimension (N-1)
         On entry, the off-diagonal elements of the tridiagonal matrix.
         On exit, E has been destroyed.
[in,out]Q
          Q is COMPLEX array, dimension (LDQ,N)
         On entry, Q must contain an QSIZ x N matrix whose columns
         unitarily orthonormal. It is a part of the unitary matrix
         that reduces the full dense Hermitian matrix to a
         (reducible) symmetric tridiagonal matrix.
[in]LDQ
          LDQ is INTEGER
         The leading dimension of the array Q.  LDQ >= max(1,N).
[out]IWORK
          IWORK is INTEGER array,
         the dimension of IWORK must be at least
                      6 + 6*N + 5*N*lg N
                      ( lg( N ) = smallest integer k
                                  such that 2^k >= N )
[out]RWORK
          RWORK is REAL array,
                               dimension (1 + 3*N + 2*N*lg N + 3*N**2)
                        ( lg( N ) = smallest integer k
                                    such that 2^k >= N )
[out]QSTORE
          QSTORE is COMPLEX array, dimension (LDQS, N)
         Used to store parts of
         the eigenvector matrix when the updating matrix multiplies
         take place.
[in]LDQS
          LDQS is INTEGER
         The leading dimension of the array QSTORE.
         LDQS >= max(1,N).
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
          > 0:  The algorithm failed to compute an eigenvalue while
                working on the submatrix lying in rows and columns
                INFO/(N+1) through mod(INFO,N+1).
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 143 of file claed0.f.

145*
146* -- LAPACK computational routine --
147* -- LAPACK is a software package provided by Univ. of Tennessee, --
148* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
149*
150* .. Scalar Arguments ..
151 INTEGER INFO, LDQ, LDQS, N, QSIZ
152* ..
153* .. Array Arguments ..
154 INTEGER IWORK( * )
155 REAL D( * ), E( * ), RWORK( * )
156 COMPLEX Q( LDQ, * ), QSTORE( LDQS, * )
157* ..
158*
159* =====================================================================
160*
161* Warning: N could be as big as QSIZ!
162*
163* .. Parameters ..
164 REAL TWO
165 parameter( two = 2.e+0 )
166* ..
167* .. Local Scalars ..
168 INTEGER CURLVL, CURPRB, CURR, I, IGIVCL, IGIVNM,
169 $ IGIVPT, INDXQ, IPERM, IPRMPT, IQ, IQPTR, IWREM,
170 $ J, K, LGN, LL, MATSIZ, MSD2, SMLSIZ, SMM1,
171 $ SPM1, SPM2, SUBMAT, SUBPBS, TLVLS
172 REAL TEMP
173* ..
174* .. External Subroutines ..
175 EXTERNAL ccopy, clacrm, claed7, scopy, ssteqr, xerbla
176* ..
177* .. External Functions ..
178 INTEGER ILAENV
179 EXTERNAL ilaenv
180* ..
181* .. Intrinsic Functions ..
182 INTRINSIC abs, int, log, max, real
183* ..
184* .. Executable Statements ..
185*
186* Test the input parameters.
187*
188 info = 0
189*
190* IF( ICOMPQ .LT. 0 .OR. ICOMPQ .GT. 2 ) THEN
191* INFO = -1
192* ELSE IF( ( ICOMPQ .EQ. 1 ) .AND. ( QSIZ .LT. MAX( 0, N ) ) )
193* $ THEN
194 IF( qsiz.LT.max( 0, n ) ) THEN
195 info = -1
196 ELSE IF( n.LT.0 ) THEN
197 info = -2
198 ELSE IF( ldq.LT.max( 1, n ) ) THEN
199 info = -6
200 ELSE IF( ldqs.LT.max( 1, n ) ) THEN
201 info = -8
202 END IF
203 IF( info.NE.0 ) THEN
204 CALL xerbla( 'CLAED0', -info )
205 RETURN
206 END IF
207*
208* Quick return if possible
209*
210 IF( n.EQ.0 )
211 $ RETURN
212*
213 smlsiz = ilaenv( 9, 'CLAED0', ' ', 0, 0, 0, 0 )
214*
215* Determine the size and placement of the submatrices, and save in
216* the leading elements of IWORK.
217*
218 iwork( 1 ) = n
219 subpbs = 1
220 tlvls = 0
221 10 CONTINUE
222 IF( iwork( subpbs ).GT.smlsiz ) THEN
223 DO 20 j = subpbs, 1, -1
224 iwork( 2*j ) = ( iwork( j )+1 ) / 2
225 iwork( 2*j-1 ) = iwork( j ) / 2
226 20 CONTINUE
227 tlvls = tlvls + 1
228 subpbs = 2*subpbs
229 GO TO 10
230 END IF
231 DO 30 j = 2, subpbs
232 iwork( j ) = iwork( j ) + iwork( j-1 )
233 30 CONTINUE
234*
235* Divide the matrix into SUBPBS submatrices of size at most SMLSIZ+1
236* using rank-1 modifications (cuts).
237*
238 spm1 = subpbs - 1
239 DO 40 i = 1, spm1
240 submat = iwork( i ) + 1
241 smm1 = submat - 1
242 d( smm1 ) = d( smm1 ) - abs( e( smm1 ) )
243 d( submat ) = d( submat ) - abs( e( smm1 ) )
244 40 CONTINUE
245*
246 indxq = 4*n + 3
247*
248* Set up workspaces for eigenvalues only/accumulate new vectors
249* routine
250*
251 temp = log( real( n ) ) / log( two )
252 lgn = int( temp )
253 IF( 2**lgn.LT.n )
254 $ lgn = lgn + 1
255 IF( 2**lgn.LT.n )
256 $ lgn = lgn + 1
257 iprmpt = indxq + n + 1
258 iperm = iprmpt + n*lgn
259 iqptr = iperm + n*lgn
260 igivpt = iqptr + n + 2
261 igivcl = igivpt + n*lgn
262*
263 igivnm = 1
264 iq = igivnm + 2*n*lgn
265 iwrem = iq + n**2 + 1
266* Initialize pointers
267 DO 50 i = 0, subpbs
268 iwork( iprmpt+i ) = 1
269 iwork( igivpt+i ) = 1
270 50 CONTINUE
271 iwork( iqptr ) = 1
272*
273* Solve each submatrix eigenproblem at the bottom of the divide and
274* conquer tree.
275*
276 curr = 0
277 DO 70 i = 0, spm1
278 IF( i.EQ.0 ) THEN
279 submat = 1
280 matsiz = iwork( 1 )
281 ELSE
282 submat = iwork( i ) + 1
283 matsiz = iwork( i+1 ) - iwork( i )
284 END IF
285 ll = iq - 1 + iwork( iqptr+curr )
286 CALL ssteqr( 'I', matsiz, d( submat ), e( submat ),
287 $ rwork( ll ), matsiz, rwork, info )
288 CALL clacrm( qsiz, matsiz, q( 1, submat ), ldq, rwork( ll ),
289 $ matsiz, qstore( 1, submat ), ldqs,
290 $ rwork( iwrem ) )
291 iwork( iqptr+curr+1 ) = iwork( iqptr+curr ) + matsiz**2
292 curr = curr + 1
293 IF( info.GT.0 ) THEN
294 info = submat*( n+1 ) + submat + matsiz - 1
295 RETURN
296 END IF
297 k = 1
298 DO 60 j = submat, iwork( i+1 )
299 iwork( indxq+j ) = k
300 k = k + 1
301 60 CONTINUE
302 70 CONTINUE
303*
304* Successively merge eigensystems of adjacent submatrices
305* into eigensystem for the corresponding larger matrix.
306*
307* while ( SUBPBS > 1 )
308*
309 curlvl = 1
310 80 CONTINUE
311 IF( subpbs.GT.1 ) THEN
312 spm2 = subpbs - 2
313 DO 90 i = 0, spm2, 2
314 IF( i.EQ.0 ) THEN
315 submat = 1
316 matsiz = iwork( 2 )
317 msd2 = iwork( 1 )
318 curprb = 0
319 ELSE
320 submat = iwork( i ) + 1
321 matsiz = iwork( i+2 ) - iwork( i )
322 msd2 = matsiz / 2
323 curprb = curprb + 1
324 END IF
325*
326* Merge lower order eigensystems (of size MSD2 and MATSIZ - MSD2)
327* into an eigensystem of size MATSIZ. CLAED7 handles the case
328* when the eigenvectors of a full or band Hermitian matrix (which
329* was reduced to tridiagonal form) are desired.
330*
331* I am free to use Q as a valuable working space until Loop 150.
332*
333 CALL claed7( matsiz, msd2, qsiz, tlvls, curlvl, curprb,
334 $ d( submat ), qstore( 1, submat ), ldqs,
335 $ e( submat+msd2-1 ), iwork( indxq+submat ),
336 $ rwork( iq ), iwork( iqptr ), iwork( iprmpt ),
337 $ iwork( iperm ), iwork( igivpt ),
338 $ iwork( igivcl ), rwork( igivnm ),
339 $ q( 1, submat ), rwork( iwrem ),
340 $ iwork( subpbs+1 ), info )
341 IF( info.GT.0 ) THEN
342 info = submat*( n+1 ) + submat + matsiz - 1
343 RETURN
344 END IF
345 iwork( i / 2+1 ) = iwork( i+2 )
346 90 CONTINUE
347 subpbs = subpbs / 2
348 curlvl = curlvl + 1
349 GO TO 80
350 END IF
351*
352* end while
353*
354* Re-merge the eigenvalues/vectors which were deflated at the final
355* merge step.
356*
357 DO 100 i = 1, n
358 j = iwork( indxq+i )
359 rwork( i ) = d( j )
360 CALL ccopy( qsiz, qstore( 1, j ), 1, q( 1, i ), 1 )
361 100 CONTINUE
362 CALL scopy( n, rwork, 1, d, 1 )
363*
364 RETURN
365*
366* End of CLAED0
367*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine clacrm(m, n, a, lda, b, ldb, c, ldc, rwork)
CLACRM multiplies a complex matrix by a square real matrix.
Definition clacrm.f:114
subroutine claed7(n, cutpnt, qsiz, tlvls, curlvl, curpbm, d, q, ldq, rho, indxq, qstore, qptr, prmptr, perm, givptr, givcol, givnum, work, rwork, iwork, info)
CLAED7 used by CSTEDC. Computes the updated eigensystem of a diagonal matrix after modification by a ...
Definition claed7.f:249
subroutine ssteqr(compz, n, d, e, z, ldz, work, info)
SSTEQR
Definition ssteqr.f:131
Here is the call graph for this function:
Here is the caller graph for this function: