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.

dsbgvd.f 11 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363
  1. *> \brief \b DSBGVD
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download DSBGVD + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsbgvd.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsbgvd.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsbgvd.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE DSBGVD( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W,
  22. * Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO )
  23. *
  24. * .. Scalar Arguments ..
  25. * CHARACTER JOBZ, UPLO
  26. * INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LWORK, N
  27. * ..
  28. * .. Array Arguments ..
  29. * INTEGER IWORK( * )
  30. * DOUBLE PRECISION AB( LDAB, * ), BB( LDBB, * ), W( * ),
  31. * $ WORK( * ), Z( LDZ, * )
  32. * ..
  33. *
  34. *
  35. *> \par Purpose:
  36. * =============
  37. *>
  38. *> \verbatim
  39. *>
  40. *> DSBGVD computes all the eigenvalues, and optionally, the eigenvectors
  41. *> of a real generalized symmetric-definite banded eigenproblem, of the
  42. *> form A*x=(lambda)*B*x. Here A and B are assumed to be symmetric and
  43. *> banded, and B is also positive definite. If eigenvectors are
  44. *> desired, it uses a divide and conquer algorithm.
  45. *>
  46. *> \endverbatim
  47. *
  48. * Arguments:
  49. * ==========
  50. *
  51. *> \param[in] JOBZ
  52. *> \verbatim
  53. *> JOBZ is CHARACTER*1
  54. *> = 'N': Compute eigenvalues only;
  55. *> = 'V': Compute eigenvalues and eigenvectors.
  56. *> \endverbatim
  57. *>
  58. *> \param[in] UPLO
  59. *> \verbatim
  60. *> UPLO is CHARACTER*1
  61. *> = 'U': Upper triangles of A and B are stored;
  62. *> = 'L': Lower triangles of A and B are stored.
  63. *> \endverbatim
  64. *>
  65. *> \param[in] N
  66. *> \verbatim
  67. *> N is INTEGER
  68. *> The order of the matrices A and B. N >= 0.
  69. *> \endverbatim
  70. *>
  71. *> \param[in] KA
  72. *> \verbatim
  73. *> KA is INTEGER
  74. *> The number of superdiagonals of the matrix A if UPLO = 'U',
  75. *> or the number of subdiagonals if UPLO = 'L'. KA >= 0.
  76. *> \endverbatim
  77. *>
  78. *> \param[in] KB
  79. *> \verbatim
  80. *> KB is INTEGER
  81. *> The number of superdiagonals of the matrix B if UPLO = 'U',
  82. *> or the number of subdiagonals if UPLO = 'L'. KB >= 0.
  83. *> \endverbatim
  84. *>
  85. *> \param[in,out] AB
  86. *> \verbatim
  87. *> AB is DOUBLE PRECISION array, dimension (LDAB, N)
  88. *> On entry, the upper or lower triangle of the symmetric band
  89. *> matrix A, stored in the first ka+1 rows of the array. The
  90. *> j-th column of A is stored in the j-th column of the array AB
  91. *> as follows:
  92. *> if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
  93. *> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
  94. *>
  95. *> On exit, the contents of AB are destroyed.
  96. *> \endverbatim
  97. *>
  98. *> \param[in] LDAB
  99. *> \verbatim
  100. *> LDAB is INTEGER
  101. *> The leading dimension of the array AB. LDAB >= KA+1.
  102. *> \endverbatim
  103. *>
  104. *> \param[in,out] BB
  105. *> \verbatim
  106. *> BB is DOUBLE PRECISION array, dimension (LDBB, N)
  107. *> On entry, the upper or lower triangle of the symmetric band
  108. *> matrix B, stored in the first kb+1 rows of the array. The
  109. *> j-th column of B is stored in the j-th column of the array BB
  110. *> as follows:
  111. *> if UPLO = 'U', BB(ka+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
  112. *> if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
  113. *>
  114. *> On exit, the factor S from the split Cholesky factorization
  115. *> B = S**T*S, as returned by DPBSTF.
  116. *> \endverbatim
  117. *>
  118. *> \param[in] LDBB
  119. *> \verbatim
  120. *> LDBB is INTEGER
  121. *> The leading dimension of the array BB. LDBB >= KB+1.
  122. *> \endverbatim
  123. *>
  124. *> \param[out] W
  125. *> \verbatim
  126. *> W is DOUBLE PRECISION array, dimension (N)
  127. *> If INFO = 0, the eigenvalues in ascending order.
  128. *> \endverbatim
  129. *>
  130. *> \param[out] Z
  131. *> \verbatim
  132. *> Z is DOUBLE PRECISION array, dimension (LDZ, N)
  133. *> If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
  134. *> eigenvectors, with the i-th column of Z holding the
  135. *> eigenvector associated with W(i). The eigenvectors are
  136. *> normalized so Z**T*B*Z = I.
  137. *> If JOBZ = 'N', then Z is not referenced.
  138. *> \endverbatim
  139. *>
  140. *> \param[in] LDZ
  141. *> \verbatim
  142. *> LDZ is INTEGER
  143. *> The leading dimension of the array Z. LDZ >= 1, and if
  144. *> JOBZ = 'V', LDZ >= max(1,N).
  145. *> \endverbatim
  146. *>
  147. *> \param[out] WORK
  148. *> \verbatim
  149. *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
  150. *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
  151. *> \endverbatim
  152. *>
  153. *> \param[in] LWORK
  154. *> \verbatim
  155. *> LWORK is INTEGER
  156. *> The dimension of the array WORK.
  157. *> If N <= 1, LWORK >= 1.
  158. *> If JOBZ = 'N' and N > 1, LWORK >= 2*N.
  159. *> If JOBZ = 'V' and N > 1, LWORK >= 1 + 5*N + 2*N**2.
  160. *>
  161. *> If LWORK = -1, then a workspace query is assumed; the routine
  162. *> only calculates the optimal sizes of the WORK and IWORK
  163. *> arrays, returns these values as the first entries of the WORK
  164. *> and IWORK arrays, and no error message related to LWORK or
  165. *> LIWORK is issued by XERBLA.
  166. *> \endverbatim
  167. *>
  168. *> \param[out] IWORK
  169. *> \verbatim
  170. *> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
  171. *> On exit, if LIWORK > 0, IWORK(1) returns the optimal LIWORK.
  172. *> \endverbatim
  173. *>
  174. *> \param[in] LIWORK
  175. *> \verbatim
  176. *> LIWORK is INTEGER
  177. *> The dimension of the array IWORK.
  178. *> If JOBZ = 'N' or N <= 1, LIWORK >= 1.
  179. *> If JOBZ = 'V' and N > 1, LIWORK >= 3 + 5*N.
  180. *>
  181. *> If LIWORK = -1, then a workspace query is assumed; the
  182. *> routine only calculates the optimal sizes of the WORK and
  183. *> IWORK arrays, returns these values as the first entries of
  184. *> the WORK and IWORK arrays, and no error message related to
  185. *> LWORK or LIWORK is issued by XERBLA.
  186. *> \endverbatim
  187. *>
  188. *> \param[out] INFO
  189. *> \verbatim
  190. *> INFO is INTEGER
  191. *> = 0: successful exit
  192. *> < 0: if INFO = -i, the i-th argument had an illegal value
  193. *> > 0: if INFO = i, and i is:
  194. *> <= N: the algorithm failed to converge:
  195. *> i off-diagonal elements of an intermediate
  196. *> tridiagonal form did not converge to zero;
  197. *> > N: if INFO = N + i, for 1 <= i <= N, then DPBSTF
  198. *> returned INFO = i: B is not positive definite.
  199. *> The factorization of B could not be completed and
  200. *> no eigenvalues or eigenvectors were computed.
  201. *> \endverbatim
  202. *
  203. * Authors:
  204. * ========
  205. *
  206. *> \author Univ. of Tennessee
  207. *> \author Univ. of California Berkeley
  208. *> \author Univ. of Colorado Denver
  209. *> \author NAG Ltd.
  210. *
  211. *> \ingroup doubleOTHEReigen
  212. *
  213. *> \par Contributors:
  214. * ==================
  215. *>
  216. *> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
  217. *
  218. * =====================================================================
  219. SUBROUTINE DSBGVD( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W,
  220. $ Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO )
  221. *
  222. * -- LAPACK driver routine --
  223. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  224. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  225. *
  226. * .. Scalar Arguments ..
  227. CHARACTER JOBZ, UPLO
  228. INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LWORK, N
  229. * ..
  230. * .. Array Arguments ..
  231. INTEGER IWORK( * )
  232. DOUBLE PRECISION AB( LDAB, * ), BB( LDBB, * ), W( * ),
  233. $ WORK( * ), Z( LDZ, * )
  234. * ..
  235. *
  236. * =====================================================================
  237. *
  238. * .. Parameters ..
  239. DOUBLE PRECISION ONE, ZERO
  240. PARAMETER ( ONE = 1.0D+0, ZERO = 0.0D+0 )
  241. * ..
  242. * .. Local Scalars ..
  243. LOGICAL LQUERY, UPPER, WANTZ
  244. CHARACTER VECT
  245. INTEGER IINFO, INDE, INDWK2, INDWRK, LIWMIN, LLWRK2,
  246. $ LWMIN
  247. * ..
  248. * .. External Functions ..
  249. LOGICAL LSAME
  250. EXTERNAL LSAME
  251. * ..
  252. * .. External Subroutines ..
  253. EXTERNAL DGEMM, DLACPY, DPBSTF, DSBGST, DSBTRD, DSTEDC,
  254. $ DSTERF, XERBLA
  255. * ..
  256. * .. Executable Statements ..
  257. *
  258. * Test the input parameters.
  259. *
  260. WANTZ = LSAME( JOBZ, 'V' )
  261. UPPER = LSAME( UPLO, 'U' )
  262. LQUERY = ( LWORK.EQ.-1 .OR. LIWORK.EQ.-1 )
  263. *
  264. INFO = 0
  265. IF( N.LE.1 ) THEN
  266. LIWMIN = 1
  267. LWMIN = 1
  268. ELSE IF( WANTZ ) THEN
  269. LIWMIN = 3 + 5*N
  270. LWMIN = 1 + 5*N + 2*N**2
  271. ELSE
  272. LIWMIN = 1
  273. LWMIN = 2*N
  274. END IF
  275. *
  276. IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
  277. INFO = -1
  278. ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
  279. INFO = -2
  280. ELSE IF( N.LT.0 ) THEN
  281. INFO = -3
  282. ELSE IF( KA.LT.0 ) THEN
  283. INFO = -4
  284. ELSE IF( KB.LT.0 .OR. KB.GT.KA ) THEN
  285. INFO = -5
  286. ELSE IF( LDAB.LT.KA+1 ) THEN
  287. INFO = -7
  288. ELSE IF( LDBB.LT.KB+1 ) THEN
  289. INFO = -9
  290. ELSE IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
  291. INFO = -12
  292. END IF
  293. *
  294. IF( INFO.EQ.0 ) THEN
  295. WORK( 1 ) = LWMIN
  296. IWORK( 1 ) = LIWMIN
  297. *
  298. IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY ) THEN
  299. INFO = -14
  300. ELSE IF( LIWORK.LT.LIWMIN .AND. .NOT.LQUERY ) THEN
  301. INFO = -16
  302. END IF
  303. END IF
  304. *
  305. IF( INFO.NE.0 ) THEN
  306. CALL XERBLA( 'DSBGVD', -INFO )
  307. RETURN
  308. ELSE IF( LQUERY ) THEN
  309. RETURN
  310. END IF
  311. *
  312. * Quick return if possible
  313. *
  314. IF( N.EQ.0 )
  315. $ RETURN
  316. *
  317. * Form a split Cholesky factorization of B.
  318. *
  319. CALL DPBSTF( UPLO, N, KB, BB, LDBB, INFO )
  320. IF( INFO.NE.0 ) THEN
  321. INFO = N + INFO
  322. RETURN
  323. END IF
  324. *
  325. * Transform problem to standard eigenvalue problem.
  326. *
  327. INDE = 1
  328. INDWRK = INDE + N
  329. INDWK2 = INDWRK + N*N
  330. LLWRK2 = LWORK - INDWK2 + 1
  331. CALL DSBGST( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Z, LDZ,
  332. $ WORK, IINFO )
  333. *
  334. * Reduce to tridiagonal form.
  335. *
  336. IF( WANTZ ) THEN
  337. VECT = 'U'
  338. ELSE
  339. VECT = 'N'
  340. END IF
  341. CALL DSBTRD( VECT, UPLO, N, KA, AB, LDAB, W, WORK( INDE ), Z, LDZ,
  342. $ WORK( INDWRK ), IINFO )
  343. *
  344. * For eigenvalues only, call DSTERF. For eigenvectors, call SSTEDC.
  345. *
  346. IF( .NOT.WANTZ ) THEN
  347. CALL DSTERF( N, W, WORK( INDE ), INFO )
  348. ELSE
  349. CALL DSTEDC( 'I', N, W, WORK( INDE ), WORK( INDWRK ), N,
  350. $ WORK( INDWK2 ), LLWRK2, IWORK, LIWORK, INFO )
  351. CALL DGEMM( 'N', 'N', N, N, N, ONE, Z, LDZ, WORK( INDWRK ), N,
  352. $ ZERO, WORK( INDWK2 ), N )
  353. CALL DLACPY( 'A', N, N, WORK( INDWK2 ), N, Z, LDZ )
  354. END IF
  355. *
  356. WORK( 1 ) = LWMIN
  357. IWORK( 1 ) = LIWMIN
  358. *
  359. RETURN
  360. *
  361. * End of DSBGVD
  362. *
  363. END