LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ssytrd_sy2sb.f
Go to the documentation of this file.
1*> \brief \b SSYTRD_SY2SB
2*
3* @generated from zhetrd_he2hb.f, fortran z -> s, Wed Dec 7 08:22:40 2016
4*
5* =========== DOCUMENTATION ===========
6*
7* Online html documentation available at
8* http://www.netlib.org/lapack/explore-html/
9*
10*> Download SSYTRD_SY2SB + dependencies
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssytrd_sy2sb.f">
12*> [TGZ]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssytrd_sy2sb.f">
14*> [ZIP]</a>
15*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssytrd_sy2sb.f">
16*> [TXT]</a>
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE SSYTRD_SY2SB( UPLO, N, KD, A, LDA, AB, LDAB, TAU,
22* WORK, LWORK, INFO )
23*
24* IMPLICIT NONE
25*
26* .. Scalar Arguments ..
27* CHARACTER UPLO
28* INTEGER INFO, LDA, LDAB, LWORK, N, KD
29* ..
30* .. Array Arguments ..
31* REAL A( LDA, * ), AB( LDAB, * ),
32* TAU( * ), WORK( * )
33* ..
34*
35*
36*> \par Purpose:
37* =============
38*>
39*> \verbatim
40*>
41*> SSYTRD_SY2SB reduces a real symmetric matrix A to real symmetric
42*> band-diagonal form AB by a orthogonal similarity transformation:
43*> Q**T * A * Q = AB.
44*> \endverbatim
45*
46* Arguments:
47* ==========
48*
49*> \param[in] UPLO
50*> \verbatim
51*> UPLO is CHARACTER*1
52*> = 'U': Upper triangle of A is stored;
53*> = 'L': Lower triangle of A is stored.
54*> \endverbatim
55*>
56*> \param[in] N
57*> \verbatim
58*> N is INTEGER
59*> The order of the matrix A. N >= 0.
60*> \endverbatim
61*>
62*> \param[in] KD
63*> \verbatim
64*> KD is INTEGER
65*> The number of superdiagonals of the reduced matrix if UPLO = 'U',
66*> or the number of subdiagonals if UPLO = 'L'. KD >= 0.
67*> The reduced matrix is stored in the array AB.
68*> \endverbatim
69*>
70*> \param[in,out] A
71*> \verbatim
72*> A is REAL array, dimension (LDA,N)
73*> On entry, the symmetric matrix A. If UPLO = 'U', the leading
74*> N-by-N upper triangular part of A contains the upper
75*> triangular part of the matrix A, and the strictly lower
76*> triangular part of A is not referenced. If UPLO = 'L', the
77*> leading N-by-N lower triangular part of A contains the lower
78*> triangular part of the matrix A, and the strictly upper
79*> triangular part of A is not referenced.
80*> On exit, if UPLO = 'U', the diagonal and first superdiagonal
81*> of A are overwritten by the corresponding elements of the
82*> tridiagonal matrix T, and the elements above the first
83*> superdiagonal, with the array TAU, represent the orthogonal
84*> matrix Q as a product of elementary reflectors; if UPLO
85*> = 'L', the diagonal and first subdiagonal of A are over-
86*> written by the corresponding elements of the tridiagonal
87*> matrix T, and the elements below the first subdiagonal, with
88*> the array TAU, represent the orthogonal matrix Q as a product
89*> of elementary reflectors. See Further Details.
90*> \endverbatim
91*>
92*> \param[in] LDA
93*> \verbatim
94*> LDA is INTEGER
95*> The leading dimension of the array A. LDA >= max(1,N).
96*> \endverbatim
97*>
98*> \param[out] AB
99*> \verbatim
100*> AB is REAL array, dimension (LDAB,N)
101*> On exit, the upper or lower triangle of the symmetric band
102*> matrix A, stored in the first KD+1 rows of the array. The
103*> j-th column of A is stored in the j-th column of the array AB
104*> as follows:
105*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
106*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
107*> \endverbatim
108*>
109*> \param[in] LDAB
110*> \verbatim
111*> LDAB is INTEGER
112*> The leading dimension of the array AB. LDAB >= KD+1.
113*> \endverbatim
114*>
115*> \param[out] TAU
116*> \verbatim
117*> TAU is REAL array, dimension (N-KD)
118*> The scalar factors of the elementary reflectors (see Further
119*> Details).
120*> \endverbatim
121*>
122*> \param[out] WORK
123*> \verbatim
124*> WORK is REAL array, dimension (LWORK)
125*> On exit, if INFO = 0, or if LWORK = -1,
126*> WORK(1) returns the size of LWORK.
127*> \endverbatim
128*>
129*> \param[in] LWORK
130*> \verbatim
131*> LWORK is INTEGER
132*> The dimension of the array WORK which should be calculated
133*> by a workspace query.
134*> If N <= KD+1, LWORK >= 1, else LWORK = MAX(1, LWORK_QUERY)
135*>
136*> If LWORK = -1, then a workspace query is assumed; the routine
137*> only calculates the optimal size of the WORK array, returns
138*> this value as the first entry of the WORK array, and no error
139*> message related to LWORK is issued by XERBLA.
140*> LWORK_QUERY = N*KD + N*max(KD,FACTOPTNB) + 2*KD*KD
141*> where FACTOPTNB is the blocking used by the QR or LQ
142*> algorithm, usually FACTOPTNB=128 is a good choice otherwise
143*> putting LWORK=-1 will provide the size of WORK.
144*> \endverbatim
145*>
146*> \param[out] INFO
147*> \verbatim
148*> INFO is INTEGER
149*> = 0: successful exit
150*> < 0: if INFO = -i, the i-th argument had an illegal value
151*> \endverbatim
152*
153* Authors:
154* ========
155*
156*> \author Univ. of Tennessee
157*> \author Univ. of California Berkeley
158*> \author Univ. of Colorado Denver
159*> \author NAG Ltd.
160*
161*> \ingroup hetrd_he2hb
162*
163*> \par Further Details:
164* =====================
165*>
166*> \verbatim
167*>
168*> Implemented by Azzam Haidar.
169*>
170*> All details are available on technical report, SC11, SC13 papers.
171*>
172*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
173*> Parallel reduction to condensed forms for symmetric eigenvalue problems
174*> using aggregated fine-grained and memory-aware kernels. In Proceedings
175*> of 2011 International Conference for High Performance Computing,
176*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
177*> Article 8 , 11 pages.
178*> http://doi.acm.org/10.1145/2063384.2063394
179*>
180*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
181*> An improved parallel singular value algorithm and its implementation
182*> for multicore hardware, In Proceedings of 2013 International Conference
183*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
184*> Denver, Colorado, USA, 2013.
185*> Article 90, 12 pages.
186*> http://doi.acm.org/10.1145/2503210.2503292
187*>
188*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
189*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
190*> calculations based on fine-grained memory aware tasks.
191*> International Journal of High Performance Computing Applications.
192*> Volume 28 Issue 2, Pages 196-209, May 2014.
193*> http://hpc.sagepub.com/content/28/2/196
194*>
195*> \endverbatim
196*>
197*> \verbatim
198*>
199*> If UPLO = 'U', the matrix Q is represented as a product of elementary
200*> reflectors
201*>
202*> Q = H(k)**T . . . H(2)**T H(1)**T, where k = n-kd.
203*>
204*> Each H(i) has the form
205*>
206*> H(i) = I - tau * v * v**T
207*>
208*> where tau is a real scalar, and v is a real vector with
209*> v(1:i+kd-1) = 0 and v(i+kd) = 1; conjg(v(i+kd+1:n)) is stored on exit in
210*> A(i,i+kd+1:n), and tau in TAU(i).
211*>
212*> If UPLO = 'L', the matrix Q is represented as a product of elementary
213*> reflectors
214*>
215*> Q = H(1) H(2) . . . H(k), where k = n-kd.
216*>
217*> Each H(i) has the form
218*>
219*> H(i) = I - tau * v * v**T
220*>
221*> where tau is a real scalar, and v is a real vector with
222*> v(kd+1:i) = 0 and v(i+kd+1) = 1; v(i+kd+2:n) is stored on exit in
223*> A(i+kd+2:n,i), and tau in TAU(i).
224*>
225*> The contents of A on exit are illustrated by the following examples
226*> with n = 5:
227*>
228*> if UPLO = 'U': if UPLO = 'L':
229*>
230*> ( ab ab/v1 v1 v1 v1 ) ( ab )
231*> ( ab ab/v2 v2 v2 ) ( ab/v1 ab )
232*> ( ab ab/v3 v3 ) ( v1 ab/v2 ab )
233*> ( ab ab/v4 ) ( v1 v2 ab/v3 ab )
234*> ( ab ) ( v1 v2 v3 ab/v4 ab )
235*>
236*> where d and e denote diagonal and off-diagonal elements of T, and vi
237*> denotes an element of the vector defining H(i).
238*> \endverbatim
239*>
240* =====================================================================
241 SUBROUTINE ssytrd_sy2sb( UPLO, N, KD, A, LDA, AB, LDAB, TAU,
242 $ WORK, LWORK, INFO )
243*
244 IMPLICIT NONE
245*
246* -- LAPACK computational routine --
247* -- LAPACK is a software package provided by Univ. of Tennessee, --
248* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
249*
250* .. Scalar Arguments ..
251 CHARACTER UPLO
252 INTEGER INFO, LDA, LDAB, LWORK, N, KD
253* ..
254* .. Array Arguments ..
255 REAL A( LDA, * ), AB( LDAB, * ),
256 $ tau( * ), work( * )
257* ..
258*
259* =====================================================================
260*
261* .. Parameters ..
262 REAL RONE
263 REAL ZERO, ONE, HALF
264 parameter( rone = 1.0e+0,
265 $ zero = 0.0e+0,
266 $ one = 1.0e+0,
267 $ half = 0.5e+0 )
268* ..
269* .. Local Scalars ..
270 LOGICAL LQUERY, UPPER
271 INTEGER I, J, IINFO, LWMIN, PN, PK, LK,
272 $ ldt, ldw, lds2, lds1,
273 $ ls2, ls1, lw, lt,
274 $ tpos, wpos, s2pos, s1pos
275* ..
276* .. External Subroutines ..
277 EXTERNAL xerbla, ssyr2k, ssymm, sgemm,
278 $ scopy,
280* ..
281* .. Intrinsic Functions ..
282 INTRINSIC min, max
283* ..
284* .. External Functions ..
285 LOGICAL LSAME
286 INTEGER ILAENV2STAGE
287 REAL SROUNDUP_LWORK
288 EXTERNAL lsame, ilaenv2stage, sroundup_lwork
289* ..
290* .. Executable Statements ..
291*
292* Determine the minimal workspace size required
293* and test the input parameters
294*
295 info = 0
296 upper = lsame( uplo, 'U' )
297 lquery = ( lwork.EQ.-1 )
298 IF( n.LE.kd+1 ) THEN
299 lwmin = 1
300 ELSE
301 lwmin = ilaenv2stage( 4, 'SSYTRD_SY2SB', '', n, kd, -1, -1 )
302 END IF
303*
304 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
305 info = -1
306 ELSE IF( n.LT.0 ) THEN
307 info = -2
308 ELSE IF( kd.LT.0 ) THEN
309 info = -3
310 ELSE IF( lda.LT.max( 1, n ) ) THEN
311 info = -5
312 ELSE IF( ldab.LT.max( 1, kd+1 ) ) THEN
313 info = -7
314 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
315 info = -10
316 END IF
317*
318 IF( info.NE.0 ) THEN
319 CALL xerbla( 'SSYTRD_SY2SB', -info )
320 RETURN
321 ELSE IF( lquery ) THEN
322 work( 1 ) = sroundup_lwork( lwmin )
323 RETURN
324 END IF
325*
326* Quick return if possible
327* Copy the upper/lower portion of A into AB
328*
329 IF( n.LE.kd+1 ) THEN
330 IF( upper ) THEN
331 DO 100 i = 1, n
332 lk = min( kd+1, i )
333 CALL scopy( lk, a( i-lk+1, i ), 1,
334 $ ab( kd+1-lk+1, i ), 1 )
335 100 CONTINUE
336 ELSE
337 DO 110 i = 1, n
338 lk = min( kd+1, n-i+1 )
339 CALL scopy( lk, a( i, i ), 1, ab( 1, i ), 1 )
340 110 CONTINUE
341 ENDIF
342 work( 1 ) = 1
343 RETURN
344 END IF
345*
346* Determine the pointer position for the workspace
347*
348 ldt = kd
349 lds1 = kd
350 lt = ldt*kd
351 lw = n*kd
352 ls1 = lds1*kd
353 ls2 = lwmin - lt - lw - ls1
354* LS2 = N*MAX(KD,FACTOPTNB)
355 tpos = 1
356 wpos = tpos + lt
357 s1pos = wpos + lw
358 s2pos = s1pos + ls1
359 IF( upper ) THEN
360 ldw = kd
361 lds2 = kd
362 ELSE
363 ldw = n
364 lds2 = n
365 ENDIF
366*
367*
368* Set the workspace of the triangular matrix T to zero once such a
369* way every time T is generated the upper/lower portion will be always zero
370*
371 CALL slaset( "A", ldt, kd, zero, zero, work( tpos ), ldt )
372*
373 IF( upper ) THEN
374 DO 10 i = 1, n - kd, kd
375 pn = n-i-kd+1
376 pk = min( n-i-kd+1, kd )
377*
378* Compute the LQ factorization of the current block
379*
380 CALL sgelqf( kd, pn, a( i, i+kd ), lda,
381 $ tau( i ), work( s2pos ), ls2, iinfo )
382*
383* Copy the upper portion of A into AB
384*
385 DO 20 j = i, i+pk-1
386 lk = min( kd, n-j ) + 1
387 CALL scopy( lk, a( j, j ), lda, ab( kd+1, j ),
388 $ ldab-1 )
389 20 CONTINUE
390*
391 CALL slaset( 'Lower', pk, pk, zero, one,
392 $ a( i, i+kd ), lda )
393*
394* Form the matrix T
395*
396 CALL slarft( 'Forward', 'Rowwise', pn, pk,
397 $ a( i, i+kd ), lda, tau( i ),
398 $ work( tpos ), ldt )
399*
400* Compute W:
401*
402 CALL sgemm( 'Conjugate', 'No transpose', pk, pn, pk,
403 $ one, work( tpos ), ldt,
404 $ a( i, i+kd ), lda,
405 $ zero, work( s2pos ), lds2 )
406*
407 CALL ssymm( 'Right', uplo, pk, pn,
408 $ one, a( i+kd, i+kd ), lda,
409 $ work( s2pos ), lds2,
410 $ zero, work( wpos ), ldw )
411*
412 CALL sgemm( 'No transpose', 'Conjugate', pk, pk, pn,
413 $ one, work( wpos ), ldw,
414 $ work( s2pos ), lds2,
415 $ zero, work( s1pos ), lds1 )
416*
417 CALL sgemm( 'No transpose', 'No transpose', pk, pn, pk,
418 $ -half, work( s1pos ), lds1,
419 $ a( i, i+kd ), lda,
420 $ one, work( wpos ), ldw )
421*
422*
423* Update the unreduced submatrix A(i+kd:n,i+kd:n), using
424* an update of the form: A := A - V'*W - W'*V
425*
426 CALL ssyr2k( uplo, 'Conjugate', pn, pk,
427 $ -one, a( i, i+kd ), lda,
428 $ work( wpos ), ldw,
429 $ rone, a( i+kd, i+kd ), lda )
430 10 CONTINUE
431*
432* Copy the upper band to AB which is the band storage matrix
433*
434 DO 30 j = n-kd+1, n
435 lk = min(kd, n-j) + 1
436 CALL scopy( lk, a( j, j ), lda, ab( kd+1, j ), ldab-1 )
437 30 CONTINUE
438*
439 ELSE
440*
441* Reduce the lower triangle of A to lower band matrix
442*
443 DO 40 i = 1, n - kd, kd
444 pn = n-i-kd+1
445 pk = min( n-i-kd+1, kd )
446*
447* Compute the QR factorization of the current block
448*
449 CALL sgeqrf( pn, kd, a( i+kd, i ), lda,
450 $ tau( i ), work( s2pos ), ls2, iinfo )
451*
452* Copy the upper portion of A into AB
453*
454 DO 50 j = i, i+pk-1
455 lk = min( kd, n-j ) + 1
456 CALL scopy( lk, a( j, j ), 1, ab( 1, j ), 1 )
457 50 CONTINUE
458*
459 CALL slaset( 'Upper', pk, pk, zero, one,
460 $ a( i+kd, i ), lda )
461*
462* Form the matrix T
463*
464 CALL slarft( 'Forward', 'Columnwise', pn, pk,
465 $ a( i+kd, i ), lda, tau( i ),
466 $ work( tpos ), ldt )
467*
468* Compute W:
469*
470 CALL sgemm( 'No transpose', 'No transpose', pn, pk, pk,
471 $ one, a( i+kd, i ), lda,
472 $ work( tpos ), ldt,
473 $ zero, work( s2pos ), lds2 )
474*
475 CALL ssymm( 'Left', uplo, pn, pk,
476 $ one, a( i+kd, i+kd ), lda,
477 $ work( s2pos ), lds2,
478 $ zero, work( wpos ), ldw )
479*
480 CALL sgemm( 'Conjugate', 'No transpose', pk, pk, pn,
481 $ one, work( s2pos ), lds2,
482 $ work( wpos ), ldw,
483 $ zero, work( s1pos ), lds1 )
484*
485 CALL sgemm( 'No transpose', 'No transpose', pn, pk, pk,
486 $ -half, a( i+kd, i ), lda,
487 $ work( s1pos ), lds1,
488 $ one, work( wpos ), ldw )
489*
490*
491* Update the unreduced submatrix A(i+kd:n,i+kd:n), using
492* an update of the form: A := A - V*W' - W*V'
493*
494 CALL ssyr2k( uplo, 'No transpose', pn, pk,
495 $ -one, a( i+kd, i ), lda,
496 $ work( wpos ), ldw,
497 $ rone, a( i+kd, i+kd ), lda )
498* ==================================================================
499* RESTORE A FOR COMPARISON AND CHECKING TO BE REMOVED
500* DO 45 J = I, I+PK-1
501* LK = MIN( KD, N-J ) + 1
502* CALL SCOPY( LK, AB( 1, J ), 1, A( J, J ), 1 )
503* 45 CONTINUE
504* ==================================================================
505 40 CONTINUE
506*
507* Copy the lower band to AB which is the band storage matrix
508*
509 DO 60 j = n-kd+1, n
510 lk = min(kd, n-j) + 1
511 CALL scopy( lk, a( j, j ), 1, ab( 1, j ), 1 )
512 60 CONTINUE
513
514 END IF
515*
516 work( 1 ) = sroundup_lwork( lwmin )
517 RETURN
518*
519* End of SSYTRD_SY2SB
520*
521 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sgelqf(m, n, a, lda, tau, work, lwork, info)
SGELQF
Definition sgelqf.f:142
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:188
subroutine sgeqrf(m, n, a, lda, tau, work, lwork, info)
SGEQRF
Definition sgeqrf.f:144
subroutine ssymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
SSYMM
Definition ssymm.f:189
subroutine ssyr2k(uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SSYR2K
Definition ssyr2k.f:192
subroutine ssytrd_sy2sb(uplo, n, kd, a, lda, ab, ldab, tau, work, lwork, info)
SSYTRD_SY2SB
recursive subroutine slarft(direct, storev, n, k, v, ldv, tau, t, ldt)
SLARFT forms the triangular factor T of a block reflector H = I - vtvH
Definition slarft.f:162
subroutine slaset(uplo, m, n, alpha, beta, a, lda)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition slaset.f:108