LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ ssygst()

subroutine ssygst ( integer itype,
character uplo,
integer n,
real, dimension( lda, * ) a,
integer lda,
real, dimension( ldb, * ) b,
integer ldb,
integer info )

SSYGST

Download SSYGST + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> SSYGST reduces a real symmetric-definite generalized eigenproblem
!> to standard form.
!>
!> If ITYPE = 1, the problem is A*x = lambda*B*x,
!> and A is overwritten by inv(U**T)*A*inv(U) or inv(L)*A*inv(L**T)
!>
!> If ITYPE = 2 or 3, the problem is A*B*x = lambda*x or
!> B*A*x = lambda*x, and A is overwritten by U*A*U**T or L**T*A*L.
!>
!> B must have been previously factorized as U**T*U or L*L**T by SPOTRF.
!> 
Parameters
[in]ITYPE
!>          ITYPE is INTEGER
!>          = 1: compute inv(U**T)*A*inv(U) or inv(L)*A*inv(L**T);
!>          = 2 or 3: compute U*A*U**T or L**T*A*L.
!> 
[in]UPLO
!>          UPLO is CHARACTER*1
!>          = 'U':  Upper triangle of A is stored and B is factored as
!>                  U**T*U;
!>          = 'L':  Lower triangle of A is stored and B is factored as
!>                  L*L**T.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrices A and B.  N >= 0.
!> 
[in,out]A
!>          A is REAL array, dimension (LDA,N)
!>          On entry, the symmetric matrix A.  If UPLO = 'U', the leading
!>          N-by-N upper triangular part of A contains the upper
!>          triangular part of the matrix A, and the strictly lower
!>          triangular part of A is not referenced.  If UPLO = 'L', the
!>          leading N-by-N lower triangular part of A contains the lower
!>          triangular part of the matrix A, and the strictly upper
!>          triangular part of A is not referenced.
!>
!>          On exit, if INFO = 0, the transformed matrix, stored in the
!>          same format as A.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in]B
!>          B is REAL array, dimension (LDB,N)
!>          The triangular factor from the Cholesky factorization of B,
!>          as returned by SPOTRF.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= max(1,N).
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 124 of file ssygst.f.

125*
126* -- LAPACK computational routine --
127* -- LAPACK is a software package provided by Univ. of Tennessee, --
128* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
129*
130* .. Scalar Arguments ..
131 CHARACTER UPLO
132 INTEGER INFO, ITYPE, LDA, LDB, N
133* ..
134* .. Array Arguments ..
135 REAL A( LDA, * ), B( LDB, * )
136* ..
137*
138* =====================================================================
139*
140* .. Parameters ..
141 REAL ONE, HALF
142 parameter( one = 1.0, half = 0.5 )
143* ..
144* .. Local Scalars ..
145 LOGICAL UPPER
146 INTEGER K, KB, NB
147* ..
148* .. External Subroutines ..
149 EXTERNAL ssygs2, ssymm, ssyr2k, strmm, strsm,
150 $ xerbla
151* ..
152* .. Intrinsic Functions ..
153 INTRINSIC max, min
154* ..
155* .. External Functions ..
156 LOGICAL LSAME
157 INTEGER ILAENV
158 EXTERNAL lsame, ilaenv
159* ..
160* .. Executable Statements ..
161*
162* Test the input parameters.
163*
164 info = 0
165 upper = lsame( uplo, 'U' )
166 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
167 info = -1
168 ELSE IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
169 info = -2
170 ELSE IF( n.LT.0 ) THEN
171 info = -3
172 ELSE IF( lda.LT.max( 1, n ) ) THEN
173 info = -5
174 ELSE IF( ldb.LT.max( 1, n ) ) THEN
175 info = -7
176 END IF
177 IF( info.NE.0 ) THEN
178 CALL xerbla( 'SSYGST', -info )
179 RETURN
180 END IF
181*
182* Quick return if possible
183*
184 IF( n.EQ.0 )
185 $ RETURN
186*
187* Determine the block size for this environment.
188*
189 nb = ilaenv( 1, 'SSYGST', uplo, n, -1, -1, -1 )
190*
191 IF( nb.LE.1 .OR. nb.GE.n ) THEN
192*
193* Use unblocked code
194*
195 CALL ssygs2( itype, uplo, n, a, lda, b, ldb, info )
196 ELSE
197*
198* Use blocked code
199*
200 IF( itype.EQ.1 ) THEN
201 IF( upper ) THEN
202*
203* Compute inv(U**T)*A*inv(U)
204*
205 DO 10 k = 1, n, nb
206 kb = min( n-k+1, nb )
207*
208* Update the upper triangle of A(k:n,k:n)
209*
210 CALL ssygs2( itype, uplo, kb, a( k, k ), lda,
211 $ b( k, k ), ldb, info )
212 IF( k+kb.LE.n ) THEN
213 CALL strsm( 'Left', uplo, 'Transpose',
214 $ 'Non-unit',
215 $ kb, n-k-kb+1, one, b( k, k ), ldb,
216 $ a( k, k+kb ), lda )
217 CALL ssymm( 'Left', uplo, kb, n-k-kb+1, -half,
218 $ a( k, k ), lda, b( k, k+kb ), ldb, one,
219 $ a( k, k+kb ), lda )
220 CALL ssyr2k( uplo, 'Transpose', n-k-kb+1, kb,
221 $ -one,
222 $ a( k, k+kb ), lda, b( k, k+kb ), ldb,
223 $ one, a( k+kb, k+kb ), lda )
224 CALL ssymm( 'Left', uplo, kb, n-k-kb+1, -half,
225 $ a( k, k ), lda, b( k, k+kb ), ldb, one,
226 $ a( k, k+kb ), lda )
227 CALL strsm( 'Right', uplo, 'No transpose',
228 $ 'Non-unit', kb, n-k-kb+1, one,
229 $ b( k+kb, k+kb ), ldb, a( k, k+kb ),
230 $ lda )
231 END IF
232 10 CONTINUE
233 ELSE
234*
235* Compute inv(L)*A*inv(L**T)
236*
237 DO 20 k = 1, n, nb
238 kb = min( n-k+1, nb )
239*
240* Update the lower triangle of A(k:n,k:n)
241*
242 CALL ssygs2( itype, uplo, kb, a( k, k ), lda,
243 $ b( k, k ), ldb, info )
244 IF( k+kb.LE.n ) THEN
245 CALL strsm( 'Right', uplo, 'Transpose',
246 $ 'Non-unit',
247 $ n-k-kb+1, kb, one, b( k, k ), ldb,
248 $ a( k+kb, k ), lda )
249 CALL ssymm( 'Right', uplo, n-k-kb+1, kb, -half,
250 $ a( k, k ), lda, b( k+kb, k ), ldb, one,
251 $ a( k+kb, k ), lda )
252 CALL ssyr2k( uplo, 'No transpose', n-k-kb+1, kb,
253 $ -one, a( k+kb, k ), lda, b( k+kb, k ),
254 $ ldb, one, a( k+kb, k+kb ), lda )
255 CALL ssymm( 'Right', uplo, n-k-kb+1, kb, -half,
256 $ a( k, k ), lda, b( k+kb, k ), ldb, one,
257 $ a( k+kb, k ), lda )
258 CALL strsm( 'Left', uplo, 'No transpose',
259 $ 'Non-unit', n-k-kb+1, kb, one,
260 $ b( k+kb, k+kb ), ldb, a( k+kb, k ),
261 $ lda )
262 END IF
263 20 CONTINUE
264 END IF
265 ELSE
266 IF( upper ) THEN
267*
268* Compute U*A*U**T
269*
270 DO 30 k = 1, n, nb
271 kb = min( n-k+1, nb )
272*
273* Update the upper triangle of A(1:k+kb-1,1:k+kb-1)
274*
275 CALL strmm( 'Left', uplo, 'No transpose',
276 $ 'Non-unit',
277 $ k-1, kb, one, b, ldb, a( 1, k ), lda )
278 CALL ssymm( 'Right', uplo, k-1, kb, half, a( k,
279 $ k ),
280 $ lda, b( 1, k ), ldb, one, a( 1, k ), lda )
281 CALL ssyr2k( uplo, 'No transpose', k-1, kb, one,
282 $ a( 1, k ), lda, b( 1, k ), ldb, one, a,
283 $ lda )
284 CALL ssymm( 'Right', uplo, k-1, kb, half, a( k,
285 $ k ),
286 $ lda, b( 1, k ), ldb, one, a( 1, k ), lda )
287 CALL strmm( 'Right', uplo, 'Transpose', 'Non-unit',
288 $ k-1, kb, one, b( k, k ), ldb, a( 1, k ),
289 $ lda )
290 CALL ssygs2( itype, uplo, kb, a( k, k ), lda,
291 $ b( k, k ), ldb, info )
292 30 CONTINUE
293 ELSE
294*
295* Compute L**T*A*L
296*
297 DO 40 k = 1, n, nb
298 kb = min( n-k+1, nb )
299*
300* Update the lower triangle of A(1:k+kb-1,1:k+kb-1)
301*
302 CALL strmm( 'Right', uplo, 'No transpose',
303 $ 'Non-unit',
304 $ kb, k-1, one, b, ldb, a( k, 1 ), lda )
305 CALL ssymm( 'Left', uplo, kb, k-1, half, a( k, k ),
306 $ lda, b( k, 1 ), ldb, one, a( k, 1 ), lda )
307 CALL ssyr2k( uplo, 'Transpose', k-1, kb, one,
308 $ a( k, 1 ), lda, b( k, 1 ), ldb, one, a,
309 $ lda )
310 CALL ssymm( 'Left', uplo, kb, k-1, half, a( k, k ),
311 $ lda, b( k, 1 ), ldb, one, a( k, 1 ), lda )
312 CALL strmm( 'Left', uplo, 'Transpose', 'Non-unit',
313 $ kb,
314 $ k-1, one, b( k, k ), ldb, a( k, 1 ), lda )
315 CALL ssygs2( itype, uplo, kb, a( k, k ), lda,
316 $ b( k, k ), ldb, info )
317 40 CONTINUE
318 END IF
319 END IF
320 END IF
321 RETURN
322*
323* End of SSYGST
324*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssygs2(itype, uplo, n, a, lda, b, ldb, info)
SSYGS2 reduces a symmetric definite generalized eigenproblem to standard form, using the factorizatio...
Definition ssygs2.f:125
subroutine ssymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
SSYMM
Definition ssymm.f:189
subroutine ssyr2k(uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SSYR2K
Definition ssyr2k.f:192
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine strmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRMM
Definition strmm.f:177
subroutine strsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRSM
Definition strsm.f:181
Here is the call graph for this function:
Here is the caller graph for this function: