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.

sgghd3.f 32 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900
  1. *> \brief \b SGGHD3
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download SGGHD3 + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgghd3.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgghd3.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgghd3.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE SGGHD3( COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q,
  22. * LDQ, Z, LDZ, WORK, LWORK, INFO )
  23. *
  24. * .. Scalar Arguments ..
  25. * CHARACTER COMPQ, COMPZ
  26. * INTEGER IHI, ILO, INFO, LDA, LDB, LDQ, LDZ, N, LWORK
  27. * ..
  28. * .. Array Arguments ..
  29. * REAL A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
  30. * $ Z( LDZ, * ), WORK( * )
  31. * ..
  32. *
  33. *
  34. *> \par Purpose:
  35. * =============
  36. *>
  37. *> \verbatim
  38. *>
  39. *> SGGHD3 reduces a pair of real matrices (A,B) to generalized upper
  40. *> Hessenberg form using orthogonal transformations, where A is a
  41. *> general matrix and B is upper triangular. The form of the
  42. *> generalized eigenvalue problem is
  43. *> A*x = lambda*B*x,
  44. *> and B is typically made upper triangular by computing its QR
  45. *> factorization and moving the orthogonal matrix Q to the left side
  46. *> of the equation.
  47. *>
  48. *> This subroutine simultaneously reduces A to a Hessenberg matrix H:
  49. *> Q**T*A*Z = H
  50. *> and transforms B to another upper triangular matrix T:
  51. *> Q**T*B*Z = T
  52. *> in order to reduce the problem to its standard form
  53. *> H*y = lambda*T*y
  54. *> where y = Z**T*x.
  55. *>
  56. *> The orthogonal matrices Q and Z are determined as products of Givens
  57. *> rotations. They may either be formed explicitly, or they may be
  58. *> postmultiplied into input matrices Q1 and Z1, so that
  59. *>
  60. *> Q1 * A * Z1**T = (Q1*Q) * H * (Z1*Z)**T
  61. *>
  62. *> Q1 * B * Z1**T = (Q1*Q) * T * (Z1*Z)**T
  63. *>
  64. *> If Q1 is the orthogonal matrix from the QR factorization of B in the
  65. *> original equation A*x = lambda*B*x, then SGGHD3 reduces the original
  66. *> problem to generalized Hessenberg form.
  67. *>
  68. *> This is a blocked variant of SGGHRD, using matrix-matrix
  69. *> multiplications for parts of the computation to enhance performance.
  70. *> \endverbatim
  71. *
  72. * Arguments:
  73. * ==========
  74. *
  75. *> \param[in] COMPQ
  76. *> \verbatim
  77. *> COMPQ is CHARACTER*1
  78. *> = 'N': do not compute Q;
  79. *> = 'I': Q is initialized to the unit matrix, and the
  80. *> orthogonal matrix Q is returned;
  81. *> = 'V': Q must contain an orthogonal matrix Q1 on entry,
  82. *> and the product Q1*Q is returned.
  83. *> \endverbatim
  84. *>
  85. *> \param[in] COMPZ
  86. *> \verbatim
  87. *> COMPZ is CHARACTER*1
  88. *> = 'N': do not compute Z;
  89. *> = 'I': Z is initialized to the unit matrix, and the
  90. *> orthogonal matrix Z is returned;
  91. *> = 'V': Z must contain an orthogonal matrix Z1 on entry,
  92. *> and the product Z1*Z is returned.
  93. *> \endverbatim
  94. *>
  95. *> \param[in] N
  96. *> \verbatim
  97. *> N is INTEGER
  98. *> The order of the matrices A and B. N >= 0.
  99. *> \endverbatim
  100. *>
  101. *> \param[in] ILO
  102. *> \verbatim
  103. *> ILO is INTEGER
  104. *> \endverbatim
  105. *>
  106. *> \param[in] IHI
  107. *> \verbatim
  108. *> IHI is INTEGER
  109. *>
  110. *> ILO and IHI mark the rows and columns of A which are to be
  111. *> reduced. It is assumed that A is already upper triangular
  112. *> in rows and columns 1:ILO-1 and IHI+1:N. ILO and IHI are
  113. *> normally set by a previous call to SGGBAL; otherwise they
  114. *> should be set to 1 and N respectively.
  115. *> 1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0.
  116. *> \endverbatim
  117. *>
  118. *> \param[in,out] A
  119. *> \verbatim
  120. *> A is REAL array, dimension (LDA, N)
  121. *> On entry, the N-by-N general matrix to be reduced.
  122. *> On exit, the upper triangle and the first subdiagonal of A
  123. *> are overwritten with the upper Hessenberg matrix H, and the
  124. *> rest is set to zero.
  125. *> \endverbatim
  126. *>
  127. *> \param[in] LDA
  128. *> \verbatim
  129. *> LDA is INTEGER
  130. *> The leading dimension of the array A. LDA >= max(1,N).
  131. *> \endverbatim
  132. *>
  133. *> \param[in,out] B
  134. *> \verbatim
  135. *> B is REAL array, dimension (LDB, N)
  136. *> On entry, the N-by-N upper triangular matrix B.
  137. *> On exit, the upper triangular matrix T = Q**T B Z. The
  138. *> elements below the diagonal are set to zero.
  139. *> \endverbatim
  140. *>
  141. *> \param[in] LDB
  142. *> \verbatim
  143. *> LDB is INTEGER
  144. *> The leading dimension of the array B. LDB >= max(1,N).
  145. *> \endverbatim
  146. *>
  147. *> \param[in,out] Q
  148. *> \verbatim
  149. *> Q is REAL array, dimension (LDQ, N)
  150. *> On entry, if COMPQ = 'V', the orthogonal matrix Q1,
  151. *> typically from the QR factorization of B.
  152. *> On exit, if COMPQ='I', the orthogonal matrix Q, and if
  153. *> COMPQ = 'V', the product Q1*Q.
  154. *> Not referenced if COMPQ='N'.
  155. *> \endverbatim
  156. *>
  157. *> \param[in] LDQ
  158. *> \verbatim
  159. *> LDQ is INTEGER
  160. *> The leading dimension of the array Q.
  161. *> LDQ >= N if COMPQ='V' or 'I'; LDQ >= 1 otherwise.
  162. *> \endverbatim
  163. *>
  164. *> \param[in,out] Z
  165. *> \verbatim
  166. *> Z is REAL array, dimension (LDZ, N)
  167. *> On entry, if COMPZ = 'V', the orthogonal matrix Z1.
  168. *> On exit, if COMPZ='I', the orthogonal matrix Z, and if
  169. *> COMPZ = 'V', the product Z1*Z.
  170. *> Not referenced if COMPZ='N'.
  171. *> \endverbatim
  172. *>
  173. *> \param[in] LDZ
  174. *> \verbatim
  175. *> LDZ is INTEGER
  176. *> The leading dimension of the array Z.
  177. *> LDZ >= N if COMPZ='V' or 'I'; LDZ >= 1 otherwise.
  178. *> \endverbatim
  179. *>
  180. *> \param[out] WORK
  181. *> \verbatim
  182. *> WORK is REAL array, dimension (MAX(1,LWORK))
  183. *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
  184. *> \endverbatim
  185. *>
  186. *> \param[in] LWORK
  187. *> \verbatim
  188. *> LWORK is INTEGER
  189. *> The length of the array WORK. LWORK >= 1.
  190. *> For optimum performance LWORK >= 6*N*NB, where NB is the
  191. *> optimal blocksize.
  192. *>
  193. *> If LWORK = -1, then a workspace query is assumed; the routine
  194. *> only calculates the optimal size of the WORK array, returns
  195. *> this value as the first entry of the WORK array, and no error
  196. *> message related to LWORK is issued by XERBLA.
  197. *> \endverbatim
  198. *>
  199. *> \param[out] INFO
  200. *> \verbatim
  201. *> INFO is INTEGER
  202. *> = 0: successful exit.
  203. *> < 0: if INFO = -i, the i-th argument had an illegal value.
  204. *> \endverbatim
  205. *
  206. * Authors:
  207. * ========
  208. *
  209. *> \author Univ. of Tennessee
  210. *> \author Univ. of California Berkeley
  211. *> \author Univ. of Colorado Denver
  212. *> \author NAG Ltd.
  213. *
  214. *> \ingroup gghd3
  215. *
  216. *> \par Further Details:
  217. * =====================
  218. *>
  219. *> \verbatim
  220. *>
  221. *> This routine reduces A to Hessenberg form and maintains B in triangular form
  222. *> using a blocked variant of Moler and Stewart's original algorithm,
  223. *> as described by Kagstrom, Kressner, Quintana-Orti, and Quintana-Orti
  224. *> (BIT 2008).
  225. *> \endverbatim
  226. *>
  227. * =====================================================================
  228. SUBROUTINE SGGHD3( COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q,
  229. $ LDQ, Z, LDZ, WORK, LWORK, INFO )
  230. *
  231. * -- LAPACK computational routine --
  232. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  233. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  234. *
  235. IMPLICIT NONE
  236. *
  237. * .. Scalar Arguments ..
  238. CHARACTER COMPQ, COMPZ
  239. INTEGER IHI, ILO, INFO, LDA, LDB, LDQ, LDZ, N, LWORK
  240. * ..
  241. * .. Array Arguments ..
  242. REAL A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
  243. $ Z( LDZ, * ), WORK( * )
  244. * ..
  245. *
  246. * =====================================================================
  247. *
  248. * .. Parameters ..
  249. REAL ZERO, ONE
  250. PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 )
  251. * ..
  252. * .. Local Scalars ..
  253. LOGICAL BLK22, INITQ, INITZ, LQUERY, WANTQ, WANTZ
  254. CHARACTER*1 COMPQ2, COMPZ2
  255. INTEGER COLA, I, IERR, J, J0, JCOL, JJ, JROW, K,
  256. $ KACC22, LEN, LWKOPT, N2NB, NB, NBLST, NBMIN,
  257. $ NH, NNB, NX, PPW, PPWO, PW, TOP, TOPQ
  258. REAL C, C1, C2, S, S1, S2, TEMP, TEMP1, TEMP2, TEMP3
  259. * ..
  260. * .. External Functions ..
  261. LOGICAL LSAME
  262. INTEGER ILAENV
  263. REAL SROUNDUP_LWORK
  264. EXTERNAL ILAENV, LSAME, SROUNDUP_LWORK
  265. * ..
  266. * .. External Subroutines ..
  267. EXTERNAL SGGHRD, SLARTG, SLASET, SORM22, SROT, SGEMM,
  268. $ SGEMV, STRMV, SLACPY, XERBLA
  269. * ..
  270. * .. Intrinsic Functions ..
  271. INTRINSIC MAX
  272. * ..
  273. * .. Executable Statements ..
  274. *
  275. * Decode and test the input parameters.
  276. *
  277. INFO = 0
  278. NB = ILAENV( 1, 'SGGHD3', ' ', N, ILO, IHI, -1 )
  279. NH = IHI - ILO + 1
  280. IF( NH.LE.1 ) THEN
  281. LWKOPT = 1
  282. ELSE
  283. LWKOPT = 6*N*NB
  284. END IF
  285. WORK( 1 ) = SROUNDUP_LWORK( LWKOPT )
  286. INITQ = LSAME( COMPQ, 'I' )
  287. WANTQ = INITQ .OR. LSAME( COMPQ, 'V' )
  288. INITZ = LSAME( COMPZ, 'I' )
  289. WANTZ = INITZ .OR. LSAME( COMPZ, 'V' )
  290. LQUERY = ( LWORK.EQ.-1 )
  291. *
  292. IF( .NOT.LSAME( COMPQ, 'N' ) .AND. .NOT.WANTQ ) THEN
  293. INFO = -1
  294. ELSE IF( .NOT.LSAME( COMPZ, 'N' ) .AND. .NOT.WANTZ ) THEN
  295. INFO = -2
  296. ELSE IF( N.LT.0 ) THEN
  297. INFO = -3
  298. ELSE IF( ILO.LT.1 ) THEN
  299. INFO = -4
  300. ELSE IF( IHI.GT.N .OR. IHI.LT.ILO-1 ) THEN
  301. INFO = -5
  302. ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
  303. INFO = -7
  304. ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
  305. INFO = -9
  306. ELSE IF( ( WANTQ .AND. LDQ.LT.N ) .OR. LDQ.LT.1 ) THEN
  307. INFO = -11
  308. ELSE IF( ( WANTZ .AND. LDZ.LT.N ) .OR. LDZ.LT.1 ) THEN
  309. INFO = -13
  310. ELSE IF( LWORK.LT.1 .AND. .NOT.LQUERY ) THEN
  311. INFO = -15
  312. END IF
  313. IF( INFO.NE.0 ) THEN
  314. CALL XERBLA( 'SGGHD3', -INFO )
  315. RETURN
  316. ELSE IF( LQUERY ) THEN
  317. RETURN
  318. END IF
  319. *
  320. * Initialize Q and Z if desired.
  321. *
  322. IF( INITQ )
  323. $ CALL SLASET( 'All', N, N, ZERO, ONE, Q, LDQ )
  324. IF( INITZ )
  325. $ CALL SLASET( 'All', N, N, ZERO, ONE, Z, LDZ )
  326. *
  327. * Zero out lower triangle of B.
  328. *
  329. IF( N.GT.1 )
  330. $ CALL SLASET( 'Lower', N-1, N-1, ZERO, ZERO, B(2, 1), LDB )
  331. *
  332. * Quick return if possible
  333. *
  334. IF( NH.LE.1 ) THEN
  335. WORK( 1 ) = ONE
  336. RETURN
  337. END IF
  338. *
  339. * Determine the blocksize.
  340. *
  341. NBMIN = ILAENV( 2, 'SGGHD3', ' ', N, ILO, IHI, -1 )
  342. IF( NB.GT.1 .AND. NB.LT.NH ) THEN
  343. *
  344. * Determine when to use unblocked instead of blocked code.
  345. *
  346. NX = MAX( NB, ILAENV( 3, 'SGGHD3', ' ', N, ILO, IHI, -1 ) )
  347. IF( NX.LT.NH ) THEN
  348. *
  349. * Determine if workspace is large enough for blocked code.
  350. *
  351. IF( LWORK.LT.LWKOPT ) THEN
  352. *
  353. * Not enough workspace to use optimal NB: determine the
  354. * minimum value of NB, and reduce NB or force use of
  355. * unblocked code.
  356. *
  357. NBMIN = MAX( 2, ILAENV( 2, 'SGGHD3', ' ', N, ILO, IHI,
  358. $ -1 ) )
  359. IF( LWORK.GE.6*N*NBMIN ) THEN
  360. NB = LWORK / ( 6*N )
  361. ELSE
  362. NB = 1
  363. END IF
  364. END IF
  365. END IF
  366. END IF
  367. *
  368. IF( NB.LT.NBMIN .OR. NB.GE.NH ) THEN
  369. *
  370. * Use unblocked code below
  371. *
  372. JCOL = ILO
  373. *
  374. ELSE
  375. *
  376. * Use blocked code
  377. *
  378. KACC22 = ILAENV( 16, 'SGGHD3', ' ', N, ILO, IHI, -1 )
  379. BLK22 = KACC22.EQ.2
  380. DO JCOL = ILO, IHI-2, NB
  381. NNB = MIN( NB, IHI-JCOL-1 )
  382. *
  383. * Initialize small orthogonal factors that will hold the
  384. * accumulated Givens rotations in workspace.
  385. * N2NB denotes the number of 2*NNB-by-2*NNB factors
  386. * NBLST denotes the (possibly smaller) order of the last
  387. * factor.
  388. *
  389. N2NB = ( IHI-JCOL-1 ) / NNB - 1
  390. NBLST = IHI - JCOL - N2NB*NNB
  391. CALL SLASET( 'All', NBLST, NBLST, ZERO, ONE, WORK, NBLST )
  392. PW = NBLST * NBLST + 1
  393. DO I = 1, N2NB
  394. CALL SLASET( 'All', 2*NNB, 2*NNB, ZERO, ONE,
  395. $ WORK( PW ), 2*NNB )
  396. PW = PW + 4*NNB*NNB
  397. END DO
  398. *
  399. * Reduce columns JCOL:JCOL+NNB-1 of A to Hessenberg form.
  400. *
  401. DO J = JCOL, JCOL+NNB-1
  402. *
  403. * Reduce Jth column of A. Store cosines and sines in Jth
  404. * column of A and B, respectively.
  405. *
  406. DO I = IHI, J+2, -1
  407. TEMP = A( I-1, J )
  408. CALL SLARTG( TEMP, A( I, J ), C, S, A( I-1, J ) )
  409. A( I, J ) = C
  410. B( I, J ) = S
  411. END DO
  412. *
  413. * Accumulate Givens rotations into workspace array.
  414. *
  415. PPW = ( NBLST + 1 )*( NBLST - 2 ) - J + JCOL + 1
  416. LEN = 2 + J - JCOL
  417. JROW = J + N2NB*NNB + 2
  418. DO I = IHI, JROW, -1
  419. C = A( I, J )
  420. S = B( I, J )
  421. DO JJ = PPW, PPW+LEN-1
  422. TEMP = WORK( JJ + NBLST )
  423. WORK( JJ + NBLST ) = C*TEMP - S*WORK( JJ )
  424. WORK( JJ ) = S*TEMP + C*WORK( JJ )
  425. END DO
  426. LEN = LEN + 1
  427. PPW = PPW - NBLST - 1
  428. END DO
  429. *
  430. PPWO = NBLST*NBLST + ( NNB+J-JCOL-1 )*2*NNB + NNB
  431. J0 = JROW - NNB
  432. DO JROW = J0, J+2, -NNB
  433. PPW = PPWO
  434. LEN = 2 + J - JCOL
  435. DO I = JROW+NNB-1, JROW, -1
  436. C = A( I, J )
  437. S = B( I, J )
  438. DO JJ = PPW, PPW+LEN-1
  439. TEMP = WORK( JJ + 2*NNB )
  440. WORK( JJ + 2*NNB ) = C*TEMP - S*WORK( JJ )
  441. WORK( JJ ) = S*TEMP + C*WORK( JJ )
  442. END DO
  443. LEN = LEN + 1
  444. PPW = PPW - 2*NNB - 1
  445. END DO
  446. PPWO = PPWO + 4*NNB*NNB
  447. END DO
  448. *
  449. * TOP denotes the number of top rows in A and B that will
  450. * not be updated during the next steps.
  451. *
  452. IF( JCOL.LE.2 ) THEN
  453. TOP = 0
  454. ELSE
  455. TOP = JCOL
  456. END IF
  457. *
  458. * Propagate transformations through B and replace stored
  459. * left sines/cosines by right sines/cosines.
  460. *
  461. DO JJ = N, J+1, -1
  462. *
  463. * Update JJth column of B.
  464. *
  465. DO I = MIN( JJ+1, IHI ), J+2, -1
  466. C = A( I, J )
  467. S = B( I, J )
  468. TEMP = B( I, JJ )
  469. B( I, JJ ) = C*TEMP - S*B( I-1, JJ )
  470. B( I-1, JJ ) = S*TEMP + C*B( I-1, JJ )
  471. END DO
  472. *
  473. * Annihilate B( JJ+1, JJ ).
  474. *
  475. IF( JJ.LT.IHI ) THEN
  476. TEMP = B( JJ+1, JJ+1 )
  477. CALL SLARTG( TEMP, B( JJ+1, JJ ), C, S,
  478. $ B( JJ+1, JJ+1 ) )
  479. B( JJ+1, JJ ) = ZERO
  480. CALL SROT( JJ-TOP, B( TOP+1, JJ+1 ), 1,
  481. $ B( TOP+1, JJ ), 1, C, S )
  482. A( JJ+1, J ) = C
  483. B( JJ+1, J ) = -S
  484. END IF
  485. END DO
  486. *
  487. * Update A by transformations from right.
  488. * Explicit loop unrolling provides better performance
  489. * compared to SLASR.
  490. * CALL SLASR( 'Right', 'Variable', 'Backward', IHI-TOP,
  491. * $ IHI-J, A( J+2, J ), B( J+2, J ),
  492. * $ A( TOP+1, J+1 ), LDA )
  493. *
  494. JJ = MOD( IHI-J-1, 3 )
  495. DO I = IHI-J-3, JJ+1, -3
  496. C = A( J+1+I, J )
  497. S = -B( J+1+I, J )
  498. C1 = A( J+2+I, J )
  499. S1 = -B( J+2+I, J )
  500. C2 = A( J+3+I, J )
  501. S2 = -B( J+3+I, J )
  502. *
  503. DO K = TOP+1, IHI
  504. TEMP = A( K, J+I )
  505. TEMP1 = A( K, J+I+1 )
  506. TEMP2 = A( K, J+I+2 )
  507. TEMP3 = A( K, J+I+3 )
  508. A( K, J+I+3 ) = C2*TEMP3 + S2*TEMP2
  509. TEMP2 = -S2*TEMP3 + C2*TEMP2
  510. A( K, J+I+2 ) = C1*TEMP2 + S1*TEMP1
  511. TEMP1 = -S1*TEMP2 + C1*TEMP1
  512. A( K, J+I+1 ) = C*TEMP1 + S*TEMP
  513. A( K, J+I ) = -S*TEMP1 + C*TEMP
  514. END DO
  515. END DO
  516. *
  517. IF( JJ.GT.0 ) THEN
  518. DO I = JJ, 1, -1
  519. CALL SROT( IHI-TOP, A( TOP+1, J+I+1 ), 1,
  520. $ A( TOP+1, J+I ), 1, A( J+1+I, J ),
  521. $ -B( J+1+I, J ) )
  522. END DO
  523. END IF
  524. *
  525. * Update (J+1)th column of A by transformations from left.
  526. *
  527. IF ( J .LT. JCOL + NNB - 1 ) THEN
  528. LEN = 1 + J - JCOL
  529. *
  530. * Multiply with the trailing accumulated orthogonal
  531. * matrix, which takes the form
  532. *
  533. * [ U11 U12 ]
  534. * U = [ ],
  535. * [ U21 U22 ]
  536. *
  537. * where U21 is a LEN-by-LEN matrix and U12 is lower
  538. * triangular.
  539. *
  540. JROW = IHI - NBLST + 1
  541. CALL SGEMV( 'Transpose', NBLST, LEN, ONE, WORK,
  542. $ NBLST, A( JROW, J+1 ), 1, ZERO,
  543. $ WORK( PW ), 1 )
  544. PPW = PW + LEN
  545. DO I = JROW, JROW+NBLST-LEN-1
  546. WORK( PPW ) = A( I, J+1 )
  547. PPW = PPW + 1
  548. END DO
  549. CALL STRMV( 'Lower', 'Transpose', 'Non-unit',
  550. $ NBLST-LEN, WORK( LEN*NBLST + 1 ), NBLST,
  551. $ WORK( PW+LEN ), 1 )
  552. CALL SGEMV( 'Transpose', LEN, NBLST-LEN, ONE,
  553. $ WORK( (LEN+1)*NBLST - LEN + 1 ), NBLST,
  554. $ A( JROW+NBLST-LEN, J+1 ), 1, ONE,
  555. $ WORK( PW+LEN ), 1 )
  556. PPW = PW
  557. DO I = JROW, JROW+NBLST-1
  558. A( I, J+1 ) = WORK( PPW )
  559. PPW = PPW + 1
  560. END DO
  561. *
  562. * Multiply with the other accumulated orthogonal
  563. * matrices, which take the form
  564. *
  565. * [ U11 U12 0 ]
  566. * [ ]
  567. * U = [ U21 U22 0 ],
  568. * [ ]
  569. * [ 0 0 I ]
  570. *
  571. * where I denotes the (NNB-LEN)-by-(NNB-LEN) identity
  572. * matrix, U21 is a LEN-by-LEN upper triangular matrix
  573. * and U12 is an NNB-by-NNB lower triangular matrix.
  574. *
  575. PPWO = 1 + NBLST*NBLST
  576. J0 = JROW - NNB
  577. DO JROW = J0, JCOL+1, -NNB
  578. PPW = PW + LEN
  579. DO I = JROW, JROW+NNB-1
  580. WORK( PPW ) = A( I, J+1 )
  581. PPW = PPW + 1
  582. END DO
  583. PPW = PW
  584. DO I = JROW+NNB, JROW+NNB+LEN-1
  585. WORK( PPW ) = A( I, J+1 )
  586. PPW = PPW + 1
  587. END DO
  588. CALL STRMV( 'Upper', 'Transpose', 'Non-unit', LEN,
  589. $ WORK( PPWO + NNB ), 2*NNB, WORK( PW ),
  590. $ 1 )
  591. CALL STRMV( 'Lower', 'Transpose', 'Non-unit', NNB,
  592. $ WORK( PPWO + 2*LEN*NNB ),
  593. $ 2*NNB, WORK( PW + LEN ), 1 )
  594. CALL SGEMV( 'Transpose', NNB, LEN, ONE,
  595. $ WORK( PPWO ), 2*NNB, A( JROW, J+1 ), 1,
  596. $ ONE, WORK( PW ), 1 )
  597. CALL SGEMV( 'Transpose', LEN, NNB, ONE,
  598. $ WORK( PPWO + 2*LEN*NNB + NNB ), 2*NNB,
  599. $ A( JROW+NNB, J+1 ), 1, ONE,
  600. $ WORK( PW+LEN ), 1 )
  601. PPW = PW
  602. DO I = JROW, JROW+LEN+NNB-1
  603. A( I, J+1 ) = WORK( PPW )
  604. PPW = PPW + 1
  605. END DO
  606. PPWO = PPWO + 4*NNB*NNB
  607. END DO
  608. END IF
  609. END DO
  610. *
  611. * Apply accumulated orthogonal matrices to A.
  612. *
  613. COLA = N - JCOL - NNB + 1
  614. J = IHI - NBLST + 1
  615. CALL SGEMM( 'Transpose', 'No Transpose', NBLST,
  616. $ COLA, NBLST, ONE, WORK, NBLST,
  617. $ A( J, JCOL+NNB ), LDA, ZERO, WORK( PW ),
  618. $ NBLST )
  619. CALL SLACPY( 'All', NBLST, COLA, WORK( PW ), NBLST,
  620. $ A( J, JCOL+NNB ), LDA )
  621. PPWO = NBLST*NBLST + 1
  622. J0 = J - NNB
  623. DO J = J0, JCOL+1, -NNB
  624. IF ( BLK22 ) THEN
  625. *
  626. * Exploit the structure of
  627. *
  628. * [ U11 U12 ]
  629. * U = [ ]
  630. * [ U21 U22 ],
  631. *
  632. * where all blocks are NNB-by-NNB, U21 is upper
  633. * triangular and U12 is lower triangular.
  634. *
  635. CALL SORM22( 'Left', 'Transpose', 2*NNB, COLA, NNB,
  636. $ NNB, WORK( PPWO ), 2*NNB,
  637. $ A( J, JCOL+NNB ), LDA, WORK( PW ),
  638. $ LWORK-PW+1, IERR )
  639. ELSE
  640. *
  641. * Ignore the structure of U.
  642. *
  643. CALL SGEMM( 'Transpose', 'No Transpose', 2*NNB,
  644. $ COLA, 2*NNB, ONE, WORK( PPWO ), 2*NNB,
  645. $ A( J, JCOL+NNB ), LDA, ZERO, WORK( PW ),
  646. $ 2*NNB )
  647. CALL SLACPY( 'All', 2*NNB, COLA, WORK( PW ), 2*NNB,
  648. $ A( J, JCOL+NNB ), LDA )
  649. END IF
  650. PPWO = PPWO + 4*NNB*NNB
  651. END DO
  652. *
  653. * Apply accumulated orthogonal matrices to Q.
  654. *
  655. IF( WANTQ ) THEN
  656. J = IHI - NBLST + 1
  657. IF ( INITQ ) THEN
  658. TOPQ = MAX( 2, J - JCOL + 1 )
  659. NH = IHI - TOPQ + 1
  660. ELSE
  661. TOPQ = 1
  662. NH = N
  663. END IF
  664. CALL SGEMM( 'No Transpose', 'No Transpose', NH,
  665. $ NBLST, NBLST, ONE, Q( TOPQ, J ), LDQ,
  666. $ WORK, NBLST, ZERO, WORK( PW ), NH )
  667. CALL SLACPY( 'All', NH, NBLST, WORK( PW ), NH,
  668. $ Q( TOPQ, J ), LDQ )
  669. PPWO = NBLST*NBLST + 1
  670. J0 = J - NNB
  671. DO J = J0, JCOL+1, -NNB
  672. IF ( INITQ ) THEN
  673. TOPQ = MAX( 2, J - JCOL + 1 )
  674. NH = IHI - TOPQ + 1
  675. END IF
  676. IF ( BLK22 ) THEN
  677. *
  678. * Exploit the structure of U.
  679. *
  680. CALL SORM22( 'Right', 'No Transpose', NH, 2*NNB,
  681. $ NNB, NNB, WORK( PPWO ), 2*NNB,
  682. $ Q( TOPQ, J ), LDQ, WORK( PW ),
  683. $ LWORK-PW+1, IERR )
  684. ELSE
  685. *
  686. * Ignore the structure of U.
  687. *
  688. CALL SGEMM( 'No Transpose', 'No Transpose', NH,
  689. $ 2*NNB, 2*NNB, ONE, Q( TOPQ, J ), LDQ,
  690. $ WORK( PPWO ), 2*NNB, ZERO, WORK( PW ),
  691. $ NH )
  692. CALL SLACPY( 'All', NH, 2*NNB, WORK( PW ), NH,
  693. $ Q( TOPQ, J ), LDQ )
  694. END IF
  695. PPWO = PPWO + 4*NNB*NNB
  696. END DO
  697. END IF
  698. *
  699. * Accumulate right Givens rotations if required.
  700. *
  701. IF ( WANTZ .OR. TOP.GT.0 ) THEN
  702. *
  703. * Initialize small orthogonal factors that will hold the
  704. * accumulated Givens rotations in workspace.
  705. *
  706. CALL SLASET( 'All', NBLST, NBLST, ZERO, ONE, WORK,
  707. $ NBLST )
  708. PW = NBLST * NBLST + 1
  709. DO I = 1, N2NB
  710. CALL SLASET( 'All', 2*NNB, 2*NNB, ZERO, ONE,
  711. $ WORK( PW ), 2*NNB )
  712. PW = PW + 4*NNB*NNB
  713. END DO
  714. *
  715. * Accumulate Givens rotations into workspace array.
  716. *
  717. DO J = JCOL, JCOL+NNB-1
  718. PPW = ( NBLST + 1 )*( NBLST - 2 ) - J + JCOL + 1
  719. LEN = 2 + J - JCOL
  720. JROW = J + N2NB*NNB + 2
  721. DO I = IHI, JROW, -1
  722. C = A( I, J )
  723. A( I, J ) = ZERO
  724. S = B( I, J )
  725. B( I, J ) = ZERO
  726. DO JJ = PPW, PPW+LEN-1
  727. TEMP = WORK( JJ + NBLST )
  728. WORK( JJ + NBLST ) = C*TEMP - S*WORK( JJ )
  729. WORK( JJ ) = S*TEMP + C*WORK( JJ )
  730. END DO
  731. LEN = LEN + 1
  732. PPW = PPW - NBLST - 1
  733. END DO
  734. *
  735. PPWO = NBLST*NBLST + ( NNB+J-JCOL-1 )*2*NNB + NNB
  736. J0 = JROW - NNB
  737. DO JROW = J0, J+2, -NNB
  738. PPW = PPWO
  739. LEN = 2 + J - JCOL
  740. DO I = JROW+NNB-1, JROW, -1
  741. C = A( I, J )
  742. A( I, J ) = ZERO
  743. S = B( I, J )
  744. B( I, J ) = ZERO
  745. DO JJ = PPW, PPW+LEN-1
  746. TEMP = WORK( JJ + 2*NNB )
  747. WORK( JJ + 2*NNB ) = C*TEMP - S*WORK( JJ )
  748. WORK( JJ ) = S*TEMP + C*WORK( JJ )
  749. END DO
  750. LEN = LEN + 1
  751. PPW = PPW - 2*NNB - 1
  752. END DO
  753. PPWO = PPWO + 4*NNB*NNB
  754. END DO
  755. END DO
  756. ELSE
  757. *
  758. CALL SLASET( 'Lower', IHI - JCOL - 1, NNB, ZERO, ZERO,
  759. $ A( JCOL + 2, JCOL ), LDA )
  760. CALL SLASET( 'Lower', IHI - JCOL - 1, NNB, ZERO, ZERO,
  761. $ B( JCOL + 2, JCOL ), LDB )
  762. END IF
  763. *
  764. * Apply accumulated orthogonal matrices to A and B.
  765. *
  766. IF ( TOP.GT.0 ) THEN
  767. J = IHI - NBLST + 1
  768. CALL SGEMM( 'No Transpose', 'No Transpose', TOP,
  769. $ NBLST, NBLST, ONE, A( 1, J ), LDA,
  770. $ WORK, NBLST, ZERO, WORK( PW ), TOP )
  771. CALL SLACPY( 'All', TOP, NBLST, WORK( PW ), TOP,
  772. $ A( 1, J ), LDA )
  773. PPWO = NBLST*NBLST + 1
  774. J0 = J - NNB
  775. DO J = J0, JCOL+1, -NNB
  776. IF ( BLK22 ) THEN
  777. *
  778. * Exploit the structure of U.
  779. *
  780. CALL SORM22( 'Right', 'No Transpose', TOP, 2*NNB,
  781. $ NNB, NNB, WORK( PPWO ), 2*NNB,
  782. $ A( 1, J ), LDA, WORK( PW ),
  783. $ LWORK-PW+1, IERR )
  784. ELSE
  785. *
  786. * Ignore the structure of U.
  787. *
  788. CALL SGEMM( 'No Transpose', 'No Transpose', TOP,
  789. $ 2*NNB, 2*NNB, ONE, A( 1, J ), LDA,
  790. $ WORK( PPWO ), 2*NNB, ZERO,
  791. $ WORK( PW ), TOP )
  792. CALL SLACPY( 'All', TOP, 2*NNB, WORK( PW ), TOP,
  793. $ A( 1, J ), LDA )
  794. END IF
  795. PPWO = PPWO + 4*NNB*NNB
  796. END DO
  797. *
  798. J = IHI - NBLST + 1
  799. CALL SGEMM( 'No Transpose', 'No Transpose', TOP,
  800. $ NBLST, NBLST, ONE, B( 1, J ), LDB,
  801. $ WORK, NBLST, ZERO, WORK( PW ), TOP )
  802. CALL SLACPY( 'All', TOP, NBLST, WORK( PW ), TOP,
  803. $ B( 1, J ), LDB )
  804. PPWO = NBLST*NBLST + 1
  805. J0 = J - NNB
  806. DO J = J0, JCOL+1, -NNB
  807. IF ( BLK22 ) THEN
  808. *
  809. * Exploit the structure of U.
  810. *
  811. CALL SORM22( 'Right', 'No Transpose', TOP, 2*NNB,
  812. $ NNB, NNB, WORK( PPWO ), 2*NNB,
  813. $ B( 1, J ), LDB, WORK( PW ),
  814. $ LWORK-PW+1, IERR )
  815. ELSE
  816. *
  817. * Ignore the structure of U.
  818. *
  819. CALL SGEMM( 'No Transpose', 'No Transpose', TOP,
  820. $ 2*NNB, 2*NNB, ONE, B( 1, J ), LDB,
  821. $ WORK( PPWO ), 2*NNB, ZERO,
  822. $ WORK( PW ), TOP )
  823. CALL SLACPY( 'All', TOP, 2*NNB, WORK( PW ), TOP,
  824. $ B( 1, J ), LDB )
  825. END IF
  826. PPWO = PPWO + 4*NNB*NNB
  827. END DO
  828. END IF
  829. *
  830. * Apply accumulated orthogonal matrices to Z.
  831. *
  832. IF( WANTZ ) THEN
  833. J = IHI - NBLST + 1
  834. IF ( INITQ ) THEN
  835. TOPQ = MAX( 2, J - JCOL + 1 )
  836. NH = IHI - TOPQ + 1
  837. ELSE
  838. TOPQ = 1
  839. NH = N
  840. END IF
  841. CALL SGEMM( 'No Transpose', 'No Transpose', NH,
  842. $ NBLST, NBLST, ONE, Z( TOPQ, J ), LDZ,
  843. $ WORK, NBLST, ZERO, WORK( PW ), NH )
  844. CALL SLACPY( 'All', NH, NBLST, WORK( PW ), NH,
  845. $ Z( TOPQ, J ), LDZ )
  846. PPWO = NBLST*NBLST + 1
  847. J0 = J - NNB
  848. DO J = J0, JCOL+1, -NNB
  849. IF ( INITQ ) THEN
  850. TOPQ = MAX( 2, J - JCOL + 1 )
  851. NH = IHI - TOPQ + 1
  852. END IF
  853. IF ( BLK22 ) THEN
  854. *
  855. * Exploit the structure of U.
  856. *
  857. CALL SORM22( 'Right', 'No Transpose', NH, 2*NNB,
  858. $ NNB, NNB, WORK( PPWO ), 2*NNB,
  859. $ Z( TOPQ, J ), LDZ, WORK( PW ),
  860. $ LWORK-PW+1, IERR )
  861. ELSE
  862. *
  863. * Ignore the structure of U.
  864. *
  865. CALL SGEMM( 'No Transpose', 'No Transpose', NH,
  866. $ 2*NNB, 2*NNB, ONE, Z( TOPQ, J ), LDZ,
  867. $ WORK( PPWO ), 2*NNB, ZERO, WORK( PW ),
  868. $ NH )
  869. CALL SLACPY( 'All', NH, 2*NNB, WORK( PW ), NH,
  870. $ Z( TOPQ, J ), LDZ )
  871. END IF
  872. PPWO = PPWO + 4*NNB*NNB
  873. END DO
  874. END IF
  875. END DO
  876. END IF
  877. *
  878. * Use unblocked code to reduce the rest of the matrix
  879. * Avoid re-initialization of modified Q and Z.
  880. *
  881. COMPQ2 = COMPQ
  882. COMPZ2 = COMPZ
  883. IF ( JCOL.NE.ILO ) THEN
  884. IF ( WANTQ )
  885. $ COMPQ2 = 'V'
  886. IF ( WANTZ )
  887. $ COMPZ2 = 'V'
  888. END IF
  889. *
  890. IF ( JCOL.LT.IHI )
  891. $ CALL SGGHRD( COMPQ2, COMPZ2, N, JCOL, IHI, A, LDA, B, LDB, Q,
  892. $ LDQ, Z, LDZ, IERR )
  893. *
  894. WORK( 1 ) = SROUNDUP_LWORK( LWKOPT )
  895. *
  896. RETURN
  897. *
  898. * End of SGGHD3
  899. *
  900. END