LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ssytrd_sb2st.F
Go to the documentation of this file.
1*> \brief \b SSYTRD_SB2ST reduces a real symmetric band matrix A to real symmetric tridiagonal form T
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SSYTRD_SB2ST + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssytrd_sb2st.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssytrd_sb2st.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssytrd_sb2st.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SSYTRD_SB2ST( STAGE1, VECT, UPLO, N, KD, AB, LDAB,
20* D, E, HOUS, LHOUS, WORK, LWORK, INFO )
21*
22* #if defined(_OPENMP)
23* use omp_lib
24* #endif
25*
26* IMPLICIT NONE
27*
28* .. Scalar Arguments ..
29* CHARACTER STAGE1, UPLO, VECT
30* INTEGER N, KD, IB, LDAB, LHOUS, LWORK, INFO
31* ..
32* .. Array Arguments ..
33* REAL D( * ), E( * )
34* REAL AB( LDAB, * ), HOUS( * ), WORK( * )
35* ..
36*
37*
38*> \par Purpose:
39* =============
40*>
41*> \verbatim
42*>
43*> SSYTRD_SB2ST reduces a real symmetric band matrix A to real symmetric
44*> tridiagonal form T by a orthogonal similarity transformation:
45*> Q**T * A * Q = T.
46*> \endverbatim
47*
48* Arguments:
49* ==========
50*
51*> \param[in] STAGE1
52*> \verbatim
53*> STAGE1 is CHARACTER*1
54*> = 'N': "No": to mention that the stage 1 of the reduction
55*> from dense to band using the ssytrd_sy2sb routine
56*> was not called before this routine to reproduce AB.
57*> In other term this routine is called as standalone.
58*> = 'Y': "Yes": to mention that the stage 1 of the
59*> reduction from dense to band using the ssytrd_sy2sb
60*> routine has been called to produce AB (e.g., AB is
61*> the output of ssytrd_sy2sb.
62*> \endverbatim
63*>
64*> \param[in] VECT
65*> \verbatim
66*> VECT is CHARACTER*1
67*> = 'N': No need for the Housholder representation,
68*> and thus LHOUS is of size max(1, 4*N);
69*> = 'V': the Householder representation is needed to
70*> either generate or to apply Q later on,
71*> then LHOUS is to be queried and computed.
72*> (NOT AVAILABLE IN THIS RELEASE).
73*> \endverbatim
74*>
75*> \param[in] UPLO
76*> \verbatim
77*> UPLO is CHARACTER*1
78*> = 'U': Upper triangle of A is stored;
79*> = 'L': Lower triangle of A is stored.
80*> \endverbatim
81*>
82*> \param[in] N
83*> \verbatim
84*> N is INTEGER
85*> The order of the matrix A. N >= 0.
86*> \endverbatim
87*>
88*> \param[in] KD
89*> \verbatim
90*> KD is INTEGER
91*> The number of superdiagonals of the matrix A if UPLO = 'U',
92*> or the number of subdiagonals if UPLO = 'L'. KD >= 0.
93*> \endverbatim
94*>
95*> \param[in,out] AB
96*> \verbatim
97*> AB is REAL array, dimension (LDAB,N)
98*> On entry, the upper or lower triangle of the symmetric band
99*> matrix A, stored in the first KD+1 rows of the array. The
100*> j-th column of A is stored in the j-th column of the array AB
101*> as follows:
102*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
103*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
104*> On exit, the diagonal elements of AB are overwritten by the
105*> diagonal elements of the tridiagonal matrix T; if KD > 0, the
106*> elements on the first superdiagonal (if UPLO = 'U') or the
107*> first subdiagonal (if UPLO = 'L') are overwritten by the
108*> off-diagonal elements of T; the rest of AB is overwritten by
109*> values generated during the reduction.
110*> \endverbatim
111*>
112*> \param[in] LDAB
113*> \verbatim
114*> LDAB is INTEGER
115*> The leading dimension of the array AB. LDAB >= KD+1.
116*> \endverbatim
117*>
118*> \param[out] D
119*> \verbatim
120*> D is REAL array, dimension (N)
121*> The diagonal elements of the tridiagonal matrix T.
122*> \endverbatim
123*>
124*> \param[out] E
125*> \verbatim
126*> E is REAL array, dimension (N-1)
127*> The off-diagonal elements of the tridiagonal matrix T:
128*> E(i) = T(i,i+1) if UPLO = 'U'; E(i) = T(i+1,i) if UPLO = 'L'.
129*> \endverbatim
130*>
131*> \param[out] HOUS
132*> \verbatim
133*> HOUS is REAL array, dimension (MAX(1,LHOUS))
134*> Stores the Householder representation.
135*> \endverbatim
136*>
137*> \param[in] LHOUS
138*> \verbatim
139*> LHOUS is INTEGER
140*> The dimension of the array HOUS.
141*> If N = 0 or KD <= 1, LHOUS >= 1, else LHOUS = MAX(1, dimension)
142*>
143*> If LWORK = -1, or LHOUS = -1,
144*> then a query is assumed; the routine
145*> only calculates the optimal size of the HOUS array, returns
146*> this value as the first entry of the HOUS array, and no error
147*> message related to LHOUS is issued by XERBLA.
148*> LHOUS = MAX(1, dimension) where
149*> dimension = 4*N if VECT='N'
150*> not available now if VECT='H'
151*> \endverbatim
152*>
153*> \param[out] WORK
154*> \verbatim
155*> WORK is REAL array, dimension (MAX(1,LWORK))
156*> On exit, if INFO = 0, WORK(1) returns optimal LWORK.
157*> \endverbatim
158*>
159*> \param[in] LWORK
160*> \verbatim
161*> LWORK is INTEGER
162*> The dimension of the array WORK.
163*> IF N = 0 or KD <= 1, LWORK >= 1, else LWORK = MAX(1, dimension)
164*>
165*> If LWORK = -1, or LHOUS = -1,
166*> then a workspace query is assumed; the routine
167*> only calculates the optimal size of the WORK array, returns
168*> this value as the first entry of the WORK array, and no error
169*> message related to LWORK is issued by XERBLA.
170*> LWORK = MAX(1, dimension) where
171*> dimension = (2KD+1)*N + KD*NTHREADS
172*> where KD is the blocking size of the reduction,
173*> FACTOPTNB is the blocking used by the QR or LQ
174*> algorithm, usually FACTOPTNB=128 is a good choice
175*> NTHREADS is the number of threads used when
176*> openMP compilation is enabled, otherwise =1.
177*> \endverbatim
178*>
179*> \param[out] INFO
180*> \verbatim
181*> INFO is INTEGER
182*> = 0: successful exit
183*> < 0: if INFO = -i, the i-th argument had an illegal value
184*> \endverbatim
185*
186* Authors:
187* ========
188*
189*> \author Univ. of Tennessee
190*> \author Univ. of California Berkeley
191*> \author Univ. of Colorado Denver
192*> \author NAG Ltd.
193*
194*> \ingroup hetrd_hb2st
195*
196*> \par Further Details:
197* =====================
198*>
199*> \verbatim
200*>
201*> Implemented by Azzam Haidar.
202*>
203*> All details are available on technical report, SC11, SC13 papers.
204*>
205*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
206*> Parallel reduction to condensed forms for symmetric eigenvalue problems
207*> using aggregated fine-grained and memory-aware kernels. In Proceedings
208*> of 2011 International Conference for High Performance Computing,
209*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
210*> Article 8 , 11 pages.
211*> http://doi.acm.org/10.1145/2063384.2063394
212*>
213*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
214*> An improved parallel singular value algorithm and its implementation
215*> for multicore hardware, In Proceedings of 2013 International Conference
216*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
217*> Denver, Colorado, USA, 2013.
218*> Article 90, 12 pages.
219*> http://doi.acm.org/10.1145/2503210.2503292
220*>
221*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
222*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
223*> calculations based on fine-grained memory aware tasks.
224*> International Journal of High Performance Computing Applications.
225*> Volume 28 Issue 2, Pages 196-209, May 2014.
226*> http://hpc.sagepub.com/content/28/2/196
227*>
228*> \endverbatim
229*>
230* =====================================================================
231 SUBROUTINE ssytrd_sb2st( STAGE1, VECT, UPLO, N, KD, AB, LDAB,
232 $ D, E, HOUS, LHOUS, WORK, LWORK, INFO )
233*
234#if defined(_OPENMP)
235 use omp_lib
236#endif
237*
238 IMPLICIT NONE
239*
240* -- LAPACK computational routine --
241* -- LAPACK is a software package provided by Univ. of Tennessee, --
242* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
243*
244* .. Scalar Arguments ..
245 CHARACTER STAGE1, UPLO, VECT
246 INTEGER N, KD, LDAB, LHOUS, LWORK, INFO
247* ..
248* .. Array Arguments ..
249 REAL D( * ), E( * )
250 REAL AB( LDAB, * ), HOUS( * ), WORK( * )
251* ..
252*
253* =====================================================================
254*
255* .. Parameters ..
256 REAL RZERO
257 REAL ZERO, ONE
258 parameter( rzero = 0.0e+0,
259 $ zero = 0.0e+0,
260 $ one = 1.0e+0 )
261* ..
262* .. Local Scalars ..
263 LOGICAL LQUERY, WANTQ, UPPER, AFTERS1
264 INTEGER I, M, K, IB, SWEEPID, MYID, SHIFT, STT, ST,
265 $ ed, stind, edind, blklastind, colpt, thed,
266 $ stepercol, grsiz, thgrsiz, thgrnb, thgrid,
267 $ nbtiles, ttype, tid, nthreads,
268 $ abdpos, abofdpos, dpos, ofdpos, awpos,
269 $ inda, indw, apos, sizea, lda, indv, indtau,
270 $ sisev, sizetau, ldv, lhmin, lwmin
271* ..
272* .. External Subroutines ..
273 EXTERNAL ssb2st_kernels, slacpy,
274 $ slaset, xerbla
275* ..
276* .. Intrinsic Functions ..
277 INTRINSIC min, max, ceiling, real
278* ..
279* .. External Functions ..
280 LOGICAL LSAME
281 INTEGER ILAENV2STAGE
282 REAL SROUNDUP_LWORK
283 EXTERNAL lsame, ilaenv2stage, sroundup_lwork
284* ..
285* .. Executable Statements ..
286*
287* Determine the minimal workspace size required.
288* Test the input parameters
289*
290 info = 0
291 afters1 = lsame( stage1, 'Y' )
292 wantq = lsame( vect, 'V' )
293 upper = lsame( uplo, 'U' )
294 lquery = ( lwork.EQ.-1 ) .OR. ( lhous.EQ.-1 )
295*
296* Determine the block size, the workspace size and the hous size.
297*
298 ib = ilaenv2stage( 2, 'SSYTRD_SB2ST', vect, n, kd,
299 $ -1, -1 )
300 IF( n.EQ.0 .OR. kd.LE.1 ) THEN
301 lhmin = 1
302 lwmin = 1
303 ELSE
304 lhmin = ilaenv2stage( 3, 'SSYTRD_SB2ST', vect, n, kd, ib,
305 $ -1 )
306 lwmin = ilaenv2stage( 4, 'SSYTRD_SB2ST', vect, n, kd, ib,
307 $ -1 )
308 END IF
309*
310 IF( .NOT.afters1 .AND. .NOT.lsame( stage1, 'N' ) ) THEN
311 info = -1
312 ELSE IF( .NOT.lsame( vect, 'N' ) ) THEN
313 info = -2
314 ELSE IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
315 info = -3
316 ELSE IF( n.LT.0 ) THEN
317 info = -4
318 ELSE IF( kd.LT.0 ) THEN
319 info = -5
320 ELSE IF( ldab.LT.(kd+1) ) THEN
321 info = -7
322 ELSE IF( lhous.LT.lhmin .AND. .NOT.lquery ) THEN
323 info = -11
324 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
325 info = -13
326 END IF
327*
328 IF( info.EQ.0 ) THEN
329 hous( 1 ) = sroundup_lwork( lhmin )
330 work( 1 ) = sroundup_lwork( lwmin )
331 END IF
332*
333 IF( info.NE.0 ) THEN
334 CALL xerbla( 'SSYTRD_SB2ST', -info )
335 RETURN
336 ELSE IF( lquery ) THEN
337 RETURN
338 END IF
339*
340* Quick return if possible
341*
342 IF( n.EQ.0 ) THEN
343 hous( 1 ) = 1
344 work( 1 ) = 1
345 RETURN
346 END IF
347*
348* Determine pointer position
349*
350 ldv = kd + ib
351 sizetau = 2 * n
352 sisev = 2 * n
353 indtau = 1
354 indv = indtau + sizetau
355 lda = 2 * kd + 1
356 sizea = lda * n
357 inda = 1
358 indw = inda + sizea
359 nthreads = 1
360 tid = 0
361*
362 IF( upper ) THEN
363 apos = inda + kd
364 awpos = inda
365 dpos = apos + kd
366 ofdpos = dpos - 1
367 abdpos = kd + 1
368 abofdpos = kd
369 ELSE
370 apos = inda
371 awpos = inda + kd + 1
372 dpos = apos
373 ofdpos = dpos + 1
374 abdpos = 1
375 abofdpos = 2
376
377 ENDIF
378*
379* Case KD=0:
380* The matrix is diagonal. We just copy it (convert to "real" for
381* real because D is double and the imaginary part should be 0)
382* and store it in D. A sequential code here is better or
383* in a parallel environment it might need two cores for D and E
384*
385 IF( kd.EQ.0 ) THEN
386 DO 30 i = 1, n
387 d( i ) = ( ab( abdpos, i ) )
388 30 CONTINUE
389 DO 40 i = 1, n-1
390 e( i ) = rzero
391 40 CONTINUE
392*
393 hous( 1 ) = 1
394 work( 1 ) = 1
395 RETURN
396 END IF
397*
398* Case KD=1:
399* The matrix is already Tridiagonal. We have to make diagonal
400* and offdiagonal elements real, and store them in D and E.
401* For that, for real precision just copy the diag and offdiag
402* to D and E while for the COMPLEX case the bulge chasing is
403* performed to convert the hermetian tridiagonal to symmetric
404* tridiagonal. A simpler conversion formula might be used, but then
405* updating the Q matrix will be required and based if Q is generated
406* or not this might complicate the story.
407*
408 IF( kd.EQ.1 ) THEN
409 DO 50 i = 1, n
410 d( i ) = ( ab( abdpos, i ) )
411 50 CONTINUE
412*
413 IF( upper ) THEN
414 DO 60 i = 1, n-1
415 e( i ) = ( ab( abofdpos, i+1 ) )
416 60 CONTINUE
417 ELSE
418 DO 70 i = 1, n-1
419 e( i ) = ( ab( abofdpos, i ) )
420 70 CONTINUE
421 ENDIF
422*
423 hous( 1 ) = 1
424 work( 1 ) = 1
425 RETURN
426 END IF
427*
428* Main code start here.
429* Reduce the symmetric band of A to a tridiagonal matrix.
430*
431 thgrsiz = n
432 grsiz = 1
433 shift = 3
434 nbtiles = ceiling( real(n)/real(kd) )
435 stepercol = ceiling( real(shift)/real(grsiz) )
436 thgrnb = ceiling( real(n-1)/real(thgrsiz) )
437*
438 CALL slacpy( "A", kd+1, n, ab, ldab, work( apos ), lda )
439 CALL slaset( "A", kd, n, zero, zero, work( awpos ), lda )
440*
441*
442* openMP parallelisation start here
443*
444#if defined(_OPENMP)
445!$OMP PARALLEL PRIVATE( TID, THGRID, BLKLASTIND )
446!$OMP$ PRIVATE( THED, I, M, K, ST, ED, STT, SWEEPID )
447!$OMP$ PRIVATE( MYID, TTYPE, COLPT, STIND, EDIND )
448!$OMP$ SHARED ( UPLO, WANTQ, INDV, INDTAU, HOUS, WORK)
449!$OMP$ SHARED ( N, KD, IB, NBTILES, LDA, LDV, INDA )
450!$OMP$ SHARED ( STEPERCOL, THGRNB, THGRSIZ, GRSIZ, SHIFT )
451!$OMP MASTER
452#endif
453*
454* main bulge chasing loop
455*
456 DO 100 thgrid = 1, thgrnb
457 stt = (thgrid-1)*thgrsiz+1
458 thed = min( (stt + thgrsiz -1), (n-1))
459 DO 110 i = stt, n-1
460 ed = min( i, thed )
461 IF( stt.GT.ed ) EXIT
462 DO 120 m = 1, stepercol
463 st = stt
464 DO 130 sweepid = st, ed
465 DO 140 k = 1, grsiz
466 myid = (i-sweepid)*(stepercol*grsiz)
467 $ + (m-1)*grsiz + k
468 IF ( myid.EQ.1 ) THEN
469 ttype = 1
470 ELSE
471 ttype = mod( myid, 2 ) + 2
472 ENDIF
473
474 IF( ttype.EQ.2 ) THEN
475 colpt = (myid/2)*kd + sweepid
476 stind = colpt-kd+1
477 edind = min(colpt,n)
478 blklastind = colpt
479 ELSE
480 colpt = ((myid+1)/2)*kd + sweepid
481 stind = colpt-kd+1
482 edind = min(colpt,n)
483 IF( ( stind.GE.edind-1 ).AND.
484 $ ( edind.EQ.n ) ) THEN
485 blklastind = n
486 ELSE
487 blklastind = 0
488 ENDIF
489 ENDIF
490*
491* Call the kernel
492*
493#if defined(_OPENMP) && _OPENMP >= 201307
494 IF( ttype.NE.1 ) THEN
495!$OMP TASK DEPEND(in:WORK(MYID+SHIFT-1))
496!$OMP$ DEPEND(in:WORK(MYID-1))
497!$OMP$ DEPEND(out:WORK(MYID))
498 tid = omp_get_thread_num()
499 CALL ssb2st_kernels(
500 $ uplo, wantq, ttype,
501 $ stind, edind, sweepid, n, kd, ib,
502 $ work( inda ), lda,
503 $ hous( indv ), hous( indtau ), ldv,
504 $ work( indw + tid*kd ) )
505!$OMP END TASK
506 ELSE
507!$OMP TASK DEPEND(in:WORK(MYID+SHIFT-1))
508!$OMP$ DEPEND(out:WORK(MYID))
509 tid = omp_get_thread_num()
510 CALL ssb2st_kernels(
511 $ uplo, wantq, ttype,
512 $ stind, edind, sweepid, n, kd, ib,
513 $ work( inda ), lda,
514 $ hous( indv ), hous( indtau ), ldv,
515 $ work( indw + tid*kd ) )
516!$OMP END TASK
517 ENDIF
518#else
519 CALL ssb2st_kernels(
520 $ uplo, wantq, ttype,
521 $ stind, edind, sweepid, n, kd, ib,
522 $ work( inda ), lda,
523 $ hous( indv ), hous( indtau ), ldv,
524 $ work( indw ) )
525#endif
526 IF ( blklastind.GE.(n-1) ) THEN
527 stt = stt + 1
528 EXIT
529 ENDIF
530 140 CONTINUE
531 130 CONTINUE
532 120 CONTINUE
533 110 CONTINUE
534 100 CONTINUE
535*
536#if defined(_OPENMP)
537!$OMP END MASTER
538!$OMP END PARALLEL
539#endif
540*
541* Copy the diagonal from A to D. Note that D is REAL thus only
542* the Real part is needed, the imaginary part should be zero.
543*
544 DO 150 i = 1, n
545 d( i ) = ( work( dpos+(i-1)*lda ) )
546 150 CONTINUE
547*
548* Copy the off diagonal from A to E. Note that E is REAL thus only
549* the Real part is needed, the imaginary part should be zero.
550*
551 IF( upper ) THEN
552 DO 160 i = 1, n-1
553 e( i ) = ( work( ofdpos+i*lda ) )
554 160 CONTINUE
555 ELSE
556 DO 170 i = 1, n-1
557 e( i ) = ( work( ofdpos+(i-1)*lda ) )
558 170 CONTINUE
559 ENDIF
560*
561 work( 1 ) = sroundup_lwork( lwmin )
562 RETURN
563*
564* End of SSYTRD_SB2ST
565*
566 END
567
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssb2st_kernels(uplo, wantz, ttype, st, ed, sweep, n, nb, ib, a, lda, v, tau, ldvt, work)
SSB2ST_KERNELS
subroutine ssytrd_sb2st(stage1, vect, uplo, n, kd, ab, ldab, d, e, hous, lhous, work, lwork, info)
SSYTRD_SB2ST reduces a real symmetric band matrix A to real symmetric tridiagonal form T
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:101
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