You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

ssytrd_sb2st.F 19 kB

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