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.

dsytrd_sb2st.F 19 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553
  1. *> \brief \b DSYTRD_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 DSYTRD_SB2ST + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsytrd_sb2st.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsytrd_sb2st.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsytrd_sb2st.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE DSYTRD_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. * DOUBLE PRECISION D( * ), E( * )
  36. * DOUBLE PRECISION AB( LDAB, * ), HOUS( * ), WORK( * )
  37. * ..
  38. *
  39. *
  40. *> \par Purpose:
  41. * =============
  42. *>
  43. *> \verbatim
  44. *>
  45. *> DSYTRD_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] STAGE1
  54. *> \verbatim
  55. *> STAGE1 is CHARACTER*1
  56. *> = 'N': "No": to mention that the stage 1 of the reduction
  57. *> from dense to band using the dsytrd_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 dsytrd_sy2sb
  62. *> routine has been called to produce AB (e.g., AB is
  63. *> the output of dsytrd_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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (N)
  123. *> The diagonal elements of the tridiagonal matrix T.
  124. *> \endverbatim
  125. *>
  126. *> \param[out] E
  127. *> \verbatim
  128. *> E is DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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. *> \ingroup real16OTHERcomputational
  192. *
  193. *> \par Further Details:
  194. * =====================
  195. *>
  196. *> \verbatim
  197. *>
  198. *> Implemented by Azzam Haidar.
  199. *>
  200. *> All details are available on technical report, SC11, SC13 papers.
  201. *>
  202. *> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
  203. *> Parallel reduction to condensed forms for symmetric eigenvalue problems
  204. *> using aggregated fine-grained and memory-aware kernels. In Proceedings
  205. *> of 2011 International Conference for High Performance Computing,
  206. *> Networking, Storage and Analysis (SC '11), New York, NY, USA,
  207. *> Article 8 , 11 pages.
  208. *> http://doi.acm.org/10.1145/2063384.2063394
  209. *>
  210. *> A. Haidar, J. Kurzak, P. Luszczek, 2013.
  211. *> An improved parallel singular value algorithm and its implementation
  212. *> for multicore hardware, In Proceedings of 2013 International Conference
  213. *> for High Performance Computing, Networking, Storage and Analysis (SC '13).
  214. *> Denver, Colorado, USA, 2013.
  215. *> Article 90, 12 pages.
  216. *> http://doi.acm.org/10.1145/2503210.2503292
  217. *>
  218. *> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
  219. *> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
  220. *> calculations based on fine-grained memory aware tasks.
  221. *> International Journal of High Performance Computing Applications.
  222. *> Volume 28 Issue 2, Pages 196-209, May 2014.
  223. *> http://hpc.sagepub.com/content/28/2/196
  224. *>
  225. *> \endverbatim
  226. *>
  227. * =====================================================================
  228. SUBROUTINE DSYTRD_SB2ST( STAGE1, VECT, UPLO, N, KD, AB, LDAB,
  229. $ D, E, HOUS, LHOUS, WORK, LWORK, INFO )
  230. *
  231. #if defined(_OPENMP)
  232. use omp_lib
  233. #endif
  234. *
  235. IMPLICIT NONE
  236. *
  237. * -- LAPACK computational routine --
  238. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  239. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  240. *
  241. * .. Scalar Arguments ..
  242. CHARACTER STAGE1, UPLO, VECT
  243. INTEGER N, KD, LDAB, LHOUS, LWORK, INFO
  244. * ..
  245. * .. Array Arguments ..
  246. DOUBLE PRECISION D( * ), E( * )
  247. DOUBLE PRECISION AB( LDAB, * ), HOUS( * ), WORK( * )
  248. * ..
  249. *
  250. * =====================================================================
  251. *
  252. * .. Parameters ..
  253. DOUBLE PRECISION RZERO
  254. DOUBLE PRECISION ZERO, ONE
  255. PARAMETER ( RZERO = 0.0D+0,
  256. $ ZERO = 0.0D+0,
  257. $ ONE = 1.0D+0 )
  258. * ..
  259. * .. Local Scalars ..
  260. LOGICAL LQUERY, WANTQ, UPPER, AFTERS1
  261. INTEGER I, M, K, IB, SWEEPID, MYID, SHIFT, STT, ST,
  262. $ ED, STIND, EDIND, BLKLASTIND, COLPT, THED,
  263. $ STEPERCOL, GRSIZ, THGRSIZ, THGRNB, THGRID,
  264. $ NBTILES, TTYPE, TID, NTHREADS, DEBUG,
  265. $ ABDPOS, ABOFDPOS, DPOS, OFDPOS, AWPOS,
  266. $ INDA, INDW, APOS, SIZEA, LDA, INDV, INDTAU,
  267. $ SIDEV, SIZETAU, LDV, LHMIN, LWMIN
  268. * ..
  269. * .. External Subroutines ..
  270. EXTERNAL DSB2ST_KERNELS, DLACPY, DLASET, XERBLA
  271. * ..
  272. * .. Intrinsic Functions ..
  273. INTRINSIC MIN, MAX, CEILING, REAL
  274. * ..
  275. * .. External Functions ..
  276. LOGICAL LSAME
  277. INTEGER ILAENV2STAGE
  278. EXTERNAL LSAME, ILAENV2STAGE
  279. * ..
  280. * .. Executable Statements ..
  281. *
  282. * Determine the minimal workspace size required.
  283. * Test the input parameters
  284. *
  285. DEBUG = 0
  286. INFO = 0
  287. AFTERS1 = LSAME( STAGE1, 'Y' )
  288. WANTQ = LSAME( VECT, 'V' )
  289. UPPER = LSAME( UPLO, 'U' )
  290. LQUERY = ( LWORK.EQ.-1 ) .OR. ( LHOUS.EQ.-1 )
  291. *
  292. * Determine the block size, the workspace size and the hous size.
  293. *
  294. IB = ILAENV2STAGE( 2, 'DSYTRD_SB2ST', VECT, N, KD, -1, -1 )
  295. LHMIN = ILAENV2STAGE( 3, 'DSYTRD_SB2ST', VECT, N, KD, IB, -1 )
  296. LWMIN = ILAENV2STAGE( 4, 'DSYTRD_SB2ST', VECT, N, KD, IB, -1 )
  297. *
  298. IF( .NOT.AFTERS1 .AND. .NOT.LSAME( STAGE1, 'N' ) ) THEN
  299. INFO = -1
  300. ELSE IF( .NOT.LSAME( VECT, 'N' ) ) THEN
  301. INFO = -2
  302. ELSE IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
  303. INFO = -3
  304. ELSE IF( N.LT.0 ) THEN
  305. INFO = -4
  306. ELSE IF( KD.LT.0 ) THEN
  307. INFO = -5
  308. ELSE IF( LDAB.LT.(KD+1) ) THEN
  309. INFO = -7
  310. ELSE IF( LHOUS.LT.LHMIN .AND. .NOT.LQUERY ) THEN
  311. INFO = -11
  312. ELSE IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY ) THEN
  313. INFO = -13
  314. END IF
  315. *
  316. IF( INFO.EQ.0 ) THEN
  317. HOUS( 1 ) = LHMIN
  318. WORK( 1 ) = LWMIN
  319. END IF
  320. *
  321. IF( INFO.NE.0 ) THEN
  322. CALL XERBLA( 'DSYTRD_SB2ST', -INFO )
  323. RETURN
  324. ELSE IF( LQUERY ) THEN
  325. RETURN
  326. END IF
  327. *
  328. * Quick return if possible
  329. *
  330. IF( N.EQ.0 ) THEN
  331. HOUS( 1 ) = 1
  332. WORK( 1 ) = 1
  333. RETURN
  334. END IF
  335. *
  336. * Determine pointer position
  337. *
  338. LDV = KD + IB
  339. SIZETAU = 2 * N
  340. SIDEV = 2 * N
  341. INDTAU = 1
  342. INDV = INDTAU + SIZETAU
  343. LDA = 2 * KD + 1
  344. SIZEA = LDA * N
  345. INDA = 1
  346. INDW = INDA + SIZEA
  347. NTHREADS = 1
  348. TID = 0
  349. *
  350. IF( UPPER ) THEN
  351. APOS = INDA + KD
  352. AWPOS = INDA
  353. DPOS = APOS + KD
  354. OFDPOS = DPOS - 1
  355. ABDPOS = KD + 1
  356. ABOFDPOS = KD
  357. ELSE
  358. APOS = INDA
  359. AWPOS = INDA + KD + 1
  360. DPOS = APOS
  361. OFDPOS = DPOS + 1
  362. ABDPOS = 1
  363. ABOFDPOS = 2
  364. ENDIF
  365. *
  366. * Case KD=0:
  367. * The matrix is diagonal. We just copy it (convert to "real" for
  368. * real because D is double and the imaginary part should be 0)
  369. * and store it in D. A sequential code here is better or
  370. * in a parallel environment it might need two cores for D and E
  371. *
  372. IF( KD.EQ.0 ) THEN
  373. DO 30 I = 1, N
  374. D( I ) = ( AB( ABDPOS, I ) )
  375. 30 CONTINUE
  376. DO 40 I = 1, N-1
  377. E( I ) = RZERO
  378. 40 CONTINUE
  379. *
  380. HOUS( 1 ) = 1
  381. WORK( 1 ) = 1
  382. RETURN
  383. END IF
  384. *
  385. * Case KD=1:
  386. * The matrix is already Tridiagonal. We have to make diagonal
  387. * and offdiagonal elements real, and store them in D and E.
  388. * For that, for real precision just copy the diag and offdiag
  389. * to D and E while for the COMPLEX case the bulge chasing is
  390. * performed to convert the hermetian tridiagonal to symmetric
  391. * tridiagonal. A simpler conversion formula might be used, but then
  392. * updating the Q matrix will be required and based if Q is generated
  393. * or not this might complicate the story.
  394. *
  395. IF( KD.EQ.1 ) THEN
  396. DO 50 I = 1, N
  397. D( I ) = ( AB( ABDPOS, I ) )
  398. 50 CONTINUE
  399. *
  400. IF( UPPER ) THEN
  401. DO 60 I = 1, N-1
  402. E( I ) = ( AB( ABOFDPOS, I+1 ) )
  403. 60 CONTINUE
  404. ELSE
  405. DO 70 I = 1, N-1
  406. E( I ) = ( AB( ABOFDPOS, I ) )
  407. 70 CONTINUE
  408. ENDIF
  409. *
  410. HOUS( 1 ) = 1
  411. WORK( 1 ) = 1
  412. RETURN
  413. END IF
  414. *
  415. * Main code start here.
  416. * Reduce the symmetric band of A to a tridiagonal matrix.
  417. *
  418. THGRSIZ = N
  419. GRSIZ = 1
  420. SHIFT = 3
  421. NBTILES = CEILING( REAL(N)/REAL(KD) )
  422. STEPERCOL = CEILING( REAL(SHIFT)/REAL(GRSIZ) )
  423. THGRNB = CEILING( REAL(N-1)/REAL(THGRSIZ) )
  424. *
  425. CALL DLACPY( "A", KD+1, N, AB, LDAB, WORK( APOS ), LDA )
  426. CALL DLASET( "A", KD, N, ZERO, ZERO, WORK( AWPOS ), LDA )
  427. *
  428. *
  429. * openMP parallelisation start here
  430. *
  431. #if defined(_OPENMP)
  432. !$OMP PARALLEL PRIVATE( TID, THGRID, BLKLASTIND )
  433. !$OMP$ PRIVATE( THED, I, M, K, ST, ED, STT, SWEEPID )
  434. !$OMP$ PRIVATE( MYID, TTYPE, COLPT, STIND, EDIND )
  435. !$OMP$ SHARED ( UPLO, WANTQ, INDV, INDTAU, HOUS, WORK)
  436. !$OMP$ SHARED ( N, KD, IB, NBTILES, LDA, LDV, INDA )
  437. !$OMP$ SHARED ( STEPERCOL, THGRNB, THGRSIZ, GRSIZ, SHIFT )
  438. !$OMP MASTER
  439. #endif
  440. *
  441. * main bulge chasing loop
  442. *
  443. DO 100 THGRID = 1, THGRNB
  444. STT = (THGRID-1)*THGRSIZ+1
  445. THED = MIN( (STT + THGRSIZ -1), (N-1))
  446. DO 110 I = STT, N-1
  447. ED = MIN( I, THED )
  448. IF( STT.GT.ED ) EXIT
  449. DO 120 M = 1, STEPERCOL
  450. ST = STT
  451. DO 130 SWEEPID = ST, ED
  452. DO 140 K = 1, GRSIZ
  453. MYID = (I-SWEEPID)*(STEPERCOL*GRSIZ)
  454. $ + (M-1)*GRSIZ + K
  455. IF ( MYID.EQ.1 ) THEN
  456. TTYPE = 1
  457. ELSE
  458. TTYPE = MOD( MYID, 2 ) + 2
  459. ENDIF
  460. IF( TTYPE.EQ.2 ) THEN
  461. COLPT = (MYID/2)*KD + SWEEPID
  462. STIND = COLPT-KD+1
  463. EDIND = MIN(COLPT,N)
  464. BLKLASTIND = COLPT
  465. ELSE
  466. COLPT = ((MYID+1)/2)*KD + SWEEPID
  467. STIND = COLPT-KD+1
  468. EDIND = MIN(COLPT,N)
  469. IF( ( STIND.GE.EDIND-1 ).AND.
  470. $ ( EDIND.EQ.N ) ) THEN
  471. BLKLASTIND = N
  472. ELSE
  473. BLKLASTIND = 0
  474. ENDIF
  475. ENDIF
  476. *
  477. * Call the kernel
  478. *
  479. #if defined(_OPENMP) && _OPENMP >= 201307
  480. IF( TTYPE.NE.1 ) THEN
  481. !$OMP TASK DEPEND(in:WORK(MYID+SHIFT-1))
  482. !$OMP$ DEPEND(in:WORK(MYID-1))
  483. !$OMP$ DEPEND(out:WORK(MYID))
  484. TID = OMP_GET_THREAD_NUM()
  485. CALL DSB2ST_KERNELS( UPLO, WANTQ, TTYPE,
  486. $ STIND, EDIND, SWEEPID, N, KD, IB,
  487. $ WORK ( INDA ), LDA,
  488. $ HOUS( INDV ), HOUS( INDTAU ), LDV,
  489. $ WORK( INDW + TID*KD ) )
  490. !$OMP END TASK
  491. ELSE
  492. !$OMP TASK DEPEND(in:WORK(MYID+SHIFT-1))
  493. !$OMP$ DEPEND(out:WORK(MYID))
  494. TID = OMP_GET_THREAD_NUM()
  495. CALL DSB2ST_KERNELS( UPLO, WANTQ, TTYPE,
  496. $ STIND, EDIND, SWEEPID, N, KD, IB,
  497. $ WORK ( INDA ), LDA,
  498. $ HOUS( INDV ), HOUS( INDTAU ), LDV,
  499. $ WORK( INDW + TID*KD ) )
  500. !$OMP END TASK
  501. ENDIF
  502. #else
  503. CALL DSB2ST_KERNELS( UPLO, WANTQ, TTYPE,
  504. $ STIND, EDIND, SWEEPID, N, KD, IB,
  505. $ WORK ( INDA ), LDA,
  506. $ HOUS( INDV ), HOUS( INDTAU ), LDV,
  507. $ WORK( INDW + TID*KD ) )
  508. #endif
  509. IF ( BLKLASTIND.GE.(N-1) ) THEN
  510. STT = STT + 1
  511. EXIT
  512. ENDIF
  513. 140 CONTINUE
  514. 130 CONTINUE
  515. 120 CONTINUE
  516. 110 CONTINUE
  517. 100 CONTINUE
  518. *
  519. #if defined(_OPENMP)
  520. !$OMP END MASTER
  521. !$OMP END PARALLEL
  522. #endif
  523. *
  524. * Copy the diagonal from A to D. Note that D is REAL thus only
  525. * the Real part is needed, the imaginary part should be zero.
  526. *
  527. DO 150 I = 1, N
  528. D( I ) = ( WORK( DPOS+(I-1)*LDA ) )
  529. 150 CONTINUE
  530. *
  531. * Copy the off diagonal from A to E. Note that E is REAL thus only
  532. * the Real part is needed, the imaginary part should be zero.
  533. *
  534. IF( UPPER ) THEN
  535. DO 160 I = 1, N-1
  536. E( I ) = ( WORK( OFDPOS+I*LDA ) )
  537. 160 CONTINUE
  538. ELSE
  539. DO 170 I = 1, N-1
  540. E( I ) = ( WORK( OFDPOS+(I-1)*LDA ) )
  541. 170 CONTINUE
  542. ENDIF
  543. *
  544. HOUS( 1 ) = LHMIN
  545. WORK( 1 ) = LWMIN
  546. RETURN
  547. *
  548. * End of DSYTRD_SB2ST
  549. *
  550. END