LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zhetrd_2stage.f
Go to the documentation of this file.
1*> \brief \b ZHETRD_2STAGE
2*
3* @precisions fortran z -> s d c
4*
5* =========== DOCUMENTATION ===========
6*
7* Online html documentation available at
8* http://www.netlib.org/lapack/explore-html/
9*
10*> Download ZHETRD_2STAGE + dependencies
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zhetrd_2stage.f">
12*> [TGZ]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zhetrd_2stage.f">
14*> [ZIP]</a>
15*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zhetrd_2stage.f">
16*> [TXT]</a>
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE ZHETRD_2STAGE( VECT, UPLO, N, A, LDA, D, E, TAU,
22* HOUS2, LHOUS2, WORK, LWORK, INFO )
23*
24* IMPLICIT NONE
25*
26* .. Scalar Arguments ..
27* CHARACTER VECT, UPLO
28* INTEGER N, LDA, LWORK, LHOUS2, INFO
29* ..
30* .. Array Arguments ..
31* DOUBLE PRECISION D( * ), E( * )
32* COMPLEX*16 A( LDA, * ), TAU( * ),
33* HOUS2( * ), WORK( * )
34* ..
35*
36*
37*> \par Purpose:
38* =============
39*>
40*> \verbatim
41*>
42*> ZHETRD_2STAGE reduces a complex Hermitian matrix A to real symmetric
43*> tridiagonal form T by a unitary similarity transformation:
44*> Q1**H Q2**H* A * Q2 * Q1 = T.
45*> \endverbatim
46*
47* Arguments:
48* ==========
49*
50*> \param[in] VECT
51*> \verbatim
52*> VECT is CHARACTER*1
53*> = 'N': No need for the Housholder representation,
54*> in particular for the second stage (Band to
55*> tridiagonal) and thus LHOUS2 is of size max(1, 4*N);
56*> = 'V': the Householder representation is needed to
57*> either generate Q1 Q2 or to apply Q1 Q2,
58*> then LHOUS2 is to be queried and computed.
59*> (NOT AVAILABLE IN THIS RELEASE).
60*> \endverbatim
61*>
62*> \param[in] UPLO
63*> \verbatim
64*> UPLO is CHARACTER*1
65*> = 'U': Upper triangle of A is stored;
66*> = 'L': Lower triangle of A is stored.
67*> \endverbatim
68*>
69*> \param[in] N
70*> \verbatim
71*> N is INTEGER
72*> The order of the matrix A. N >= 0.
73*> \endverbatim
74*>
75*> \param[in,out] A
76*> \verbatim
77*> A is COMPLEX*16 array, dimension (LDA,N)
78*> On entry, the Hermitian matrix A. If UPLO = 'U', the leading
79*> N-by-N upper triangular part of A contains the upper
80*> triangular part of the matrix A, and the strictly lower
81*> triangular part of A is not referenced. If UPLO = 'L', the
82*> leading N-by-N lower triangular part of A contains the lower
83*> triangular part of the matrix A, and the strictly upper
84*> triangular part of A is not referenced.
85*> On exit, if UPLO = 'U', the band superdiagonal
86*> of A are overwritten by the corresponding elements of the
87*> internal band-diagonal matrix AB, and the elements above
88*> the KD superdiagonal, with the array TAU, represent the unitary
89*> matrix Q1 as a product of elementary reflectors; if UPLO
90*> = 'L', the diagonal and band subdiagonal of A are over-
91*> written by the corresponding elements of the internal band-diagonal
92*> matrix AB, and the elements below the KD subdiagonal, with
93*> the array TAU, represent the unitary matrix Q1 as a product
94*> of elementary reflectors. See Further Details.
95*> \endverbatim
96*>
97*> \param[in] LDA
98*> \verbatim
99*> LDA is INTEGER
100*> The leading dimension of the array A. LDA >= max(1,N).
101*> \endverbatim
102*>
103*> \param[out] D
104*> \verbatim
105*> D is DOUBLE PRECISION array, dimension (N)
106*> The diagonal elements of the tridiagonal matrix T.
107*> \endverbatim
108*>
109*> \param[out] E
110*> \verbatim
111*> E is DOUBLE PRECISION array, dimension (N-1)
112*> The off-diagonal elements of the tridiagonal matrix T.
113*> \endverbatim
114*>
115*> \param[out] TAU
116*> \verbatim
117*> TAU is COMPLEX*16 array, dimension (N-KD)
118*> The scalar factors of the elementary reflectors of
119*> the first stage (see Further Details).
120*> \endverbatim
121*>
122*> \param[out] HOUS2
123*> \verbatim
124*> HOUS2 is COMPLEX*16 array, dimension (MAX(1,LHOUS2))
125*> Stores the Householder representation of the stage2
126*> band to tridiagonal.
127*> \endverbatim
128*>
129*> \param[in] LHOUS2
130*> \verbatim
131*> LHOUS2 is INTEGER
132*> The dimension of the array HOUS2.
133*> LHOUS2 >= 1.
134*>
135*> If LWORK = -1, or LHOUS2 = -1,
136*> then a query is assumed; the routine
137*> only calculates the optimal size of the HOUS2 array, returns
138*> this value as the first entry of the HOUS2 array, and no error
139*> message related to LHOUS2 is issued by XERBLA.
140*> If VECT='N', LHOUS2 = max(1, 4*n);
141*> if VECT='V', option not yet available.
142*> \endverbatim
143*>
144*> \param[out] WORK
145*> \verbatim
146*> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
147*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
148*> \endverbatim
149*>
150*> \param[in] LWORK
151*> \verbatim
152*> LWORK is INTEGER
153*> The dimension of the array WORK.
154*> If N = 0, LWORK >= 1, else LWORK = MAX(1, dimension).
155*>
156*> If LWORK = -1, or LHOUS2 = -1,
157*> then a workspace query is assumed; the routine
158*> only calculates the optimal size of the WORK array, returns
159*> this value as the first entry of the WORK array, and no error
160*> message related to LWORK is issued by XERBLA.
161*> LWORK = MAX(1, dimension) where
162*> dimension = max(stage1,stage2) + (KD+1)*N
163*> = N*KD + N*max(KD+1,FACTOPTNB)
164*> + max(2*KD*KD, KD*NTHREADS)
165*> + (KD+1)*N
166*> where KD is the blocking size of the reduction,
167*> FACTOPTNB is the blocking used by the QR or LQ
168*> algorithm, usually FACTOPTNB=128 is a good choice
169*> NTHREADS is the number of threads used when
170*> openMP compilation is enabled, otherwise =1.
171*> \endverbatim
172*>
173*> \param[out] INFO
174*> \verbatim
175*> INFO is INTEGER
176*> = 0: successful exit
177*> < 0: if INFO = -i, the i-th argument had an illegal value
178*> \endverbatim
179*
180* Authors:
181* ========
182*
183*> \author Univ. of Tennessee
184*> \author Univ. of California Berkeley
185*> \author Univ. of Colorado Denver
186*> \author NAG Ltd.
187*
188*> \ingroup hetrd_2stage
189*
190*> \par Further Details:
191* =====================
192*>
193*> \verbatim
194*>
195*> Implemented by Azzam Haidar.
196*>
197*> All details are available on technical report, SC11, SC13 papers.
198*>
199*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
200*> Parallel reduction to condensed forms for symmetric eigenvalue problems
201*> using aggregated fine-grained and memory-aware kernels. In Proceedings
202*> of 2011 International Conference for High Performance Computing,
203*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
204*> Article 8 , 11 pages.
205*> http://doi.acm.org/10.1145/2063384.2063394
206*>
207*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
208*> An improved parallel singular value algorithm and its implementation
209*> for multicore hardware, In Proceedings of 2013 International Conference
210*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
211*> Denver, Colorado, USA, 2013.
212*> Article 90, 12 pages.
213*> http://doi.acm.org/10.1145/2503210.2503292
214*>
215*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
216*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
217*> calculations based on fine-grained memory aware tasks.
218*> International Journal of High Performance Computing Applications.
219*> Volume 28 Issue 2, Pages 196-209, May 2014.
220*> http://hpc.sagepub.com/content/28/2/196
221*>
222*> \endverbatim
223*>
224* =====================================================================
225 SUBROUTINE zhetrd_2stage( VECT, UPLO, N, A, LDA, D, E, TAU,
226 $ HOUS2, LHOUS2, WORK, LWORK, INFO )
227*
228 IMPLICIT NONE
229*
230* -- LAPACK computational routine --
231* -- LAPACK is a software package provided by Univ. of Tennessee, --
232* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
233*
234* .. Scalar Arguments ..
235 CHARACTER VECT, UPLO
236 INTEGER N, LDA, LWORK, LHOUS2, INFO
237* ..
238* .. Array Arguments ..
239 DOUBLE PRECISION D( * ), E( * )
240 COMPLEX*16 A( LDA, * ), TAU( * ),
241 $ hous2( * ), work( * )
242* ..
243*
244* =====================================================================
245* ..
246* .. Local Scalars ..
247 LOGICAL LQUERY, UPPER, WANTQ
248 INTEGER KD, IB, LWMIN, LHMIN, LWRK, LDAB, WPOS, ABPOS
249* ..
250* .. External Subroutines ..
252* ..
253* .. External Functions ..
254 LOGICAL LSAME
255 INTEGER ILAENV2STAGE
256 EXTERNAL lsame, ilaenv2stage
257* ..
258* .. Executable Statements ..
259*
260* Test the input parameters
261*
262 info = 0
263 wantq = lsame( vect, 'V' )
264 upper = lsame( uplo, 'U' )
265 lquery = ( lwork.EQ.-1 ) .OR. ( lhous2.EQ.-1 )
266*
267* Determine the block size, the workspace size and the hous size.
268*
269 kd = ilaenv2stage( 1, 'ZHETRD_2STAGE', vect, n, -1, -1,
270 $ -1 )
271 ib = ilaenv2stage( 2, 'ZHETRD_2STAGE', vect, n, kd, -1,
272 $ -1 )
273 IF( n.EQ.0 ) THEN
274 lhmin = 1
275 lwmin = 1
276 ELSE
277 lhmin = ilaenv2stage( 3, 'ZHETRD_2STAGE', vect, n, kd, ib,
278 $ -1 )
279 lwmin = ilaenv2stage( 4, 'ZHETRD_2STAGE', vect, n, kd, ib,
280 $ -1 )
281 END IF
282*
283 IF( .NOT.lsame( vect, 'N' ) ) THEN
284 info = -1
285 ELSE IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
286 info = -2
287 ELSE IF( n.LT.0 ) THEN
288 info = -3
289 ELSE IF( lda.LT.max( 1, n ) ) THEN
290 info = -5
291 ELSE IF( lhous2.LT.lhmin .AND. .NOT.lquery ) THEN
292 info = -10
293 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
294 info = -12
295 END IF
296*
297 IF( info.EQ.0 ) THEN
298 hous2( 1 ) = lhmin
299 work( 1 ) = lwmin
300 END IF
301*
302 IF( info.NE.0 ) THEN
303 CALL xerbla( 'ZHETRD_2STAGE', -info )
304 RETURN
305 ELSE IF( lquery ) THEN
306 RETURN
307 END IF
308*
309* Quick return if possible
310*
311 IF( n.EQ.0 ) THEN
312 work( 1 ) = 1
313 RETURN
314 END IF
315*
316* Determine pointer position
317*
318 ldab = kd+1
319 lwrk = lwork-ldab*n
320 abpos = 1
321 wpos = abpos + ldab*n
322 CALL zhetrd_he2hb( uplo, n, kd, a, lda, work( abpos ), ldab,
323 $ tau, work( wpos ), lwrk, info )
324 IF( info.NE.0 ) THEN
325 CALL xerbla( 'ZHETRD_HE2HB', -info )
326 RETURN
327 END IF
328 CALL zhetrd_hb2st( 'Y', vect, uplo, n, kd,
329 $ work( abpos ), ldab, d, e,
330 $ hous2, lhous2, work( wpos ), lwrk, info )
331 IF( info.NE.0 ) THEN
332 CALL xerbla( 'ZHETRD_HB2ST', -info )
333 RETURN
334 END IF
335*
336*
337 work( 1 ) = lwmin
338 RETURN
339*
340* End of ZHETRD_2STAGE
341*
342 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhetrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
ZHETRD_2STAGE
subroutine zhetrd_hb2st(stage1, vect, uplo, n, kd, ab, ldab, d, e, hous, lhous, work, lwork, info)
ZHETRD_HB2ST reduces a complex Hermitian band matrix A to real symmetric tridiagonal form T
subroutine zhetrd_he2hb(uplo, n, kd, a, lda, ab, ldab, tau, work, lwork, info)
ZHETRD_HE2HB