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

◆ sspgst()

subroutine sspgst ( integer  itype,
character  uplo,
integer  n,
real, dimension( * )  ap,
real, dimension( * )  bp,
integer  info 
)

SSPGST

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

Purpose:
 SSPGST reduces a real symmetric-definite generalized eigenproblem
 to standard form, using packed storage.

 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 SPPTRF.
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]AP
          AP is REAL array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the symmetric matrix
          A, packed columnwise in a linear array.  The j-th column of A
          is stored in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.

          On exit, if INFO = 0, the transformed matrix, stored in the
          same format as A.
[in]BP
          BP is REAL array, dimension (N*(N+1)/2)
          The triangular factor from the Cholesky factorization of B,
          stored in the same format as A, as returned by SPPTRF.
[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 112 of file sspgst.f.

113*
114* -- LAPACK computational routine --
115* -- LAPACK is a software package provided by Univ. of Tennessee, --
116* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
117*
118* .. Scalar Arguments ..
119 CHARACTER UPLO
120 INTEGER INFO, ITYPE, N
121* ..
122* .. Array Arguments ..
123 REAL AP( * ), BP( * )
124* ..
125*
126* =====================================================================
127*
128* .. Parameters ..
129 REAL ONE, HALF
130 parameter( one = 1.0, half = 0.5 )
131* ..
132* .. Local Scalars ..
133 LOGICAL UPPER
134 INTEGER J, J1, J1J1, JJ, K, K1, K1K1, KK
135 REAL AJJ, AKK, BJJ, BKK, CT
136* ..
137* .. External Subroutines ..
138 EXTERNAL saxpy, sscal, sspmv, sspr2, stpmv, stpsv,
139 $ xerbla
140* ..
141* .. External Functions ..
142 LOGICAL LSAME
143 REAL SDOT
144 EXTERNAL lsame, sdot
145* ..
146* .. Executable Statements ..
147*
148* Test the input parameters.
149*
150 info = 0
151 upper = lsame( uplo, 'U' )
152 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
153 info = -1
154 ELSE IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
155 info = -2
156 ELSE IF( n.LT.0 ) THEN
157 info = -3
158 END IF
159 IF( info.NE.0 ) THEN
160 CALL xerbla( 'SSPGST', -info )
161 RETURN
162 END IF
163*
164 IF( itype.EQ.1 ) THEN
165 IF( upper ) THEN
166*
167* Compute inv(U**T)*A*inv(U)
168*
169* J1 and JJ are the indices of A(1,j) and A(j,j)
170*
171 jj = 0
172 DO 10 j = 1, n
173 j1 = jj + 1
174 jj = jj + j
175*
176* Compute the j-th column of the upper triangle of A
177*
178 bjj = bp( jj )
179 CALL stpsv( uplo, 'Transpose', 'Nonunit', j, bp,
180 $ ap( j1 ), 1 )
181 CALL sspmv( uplo, j-1, -one, ap, bp( j1 ), 1, one,
182 $ ap( j1 ), 1 )
183 CALL sscal( j-1, one / bjj, ap( j1 ), 1 )
184 ap( jj ) = ( ap( jj )-sdot( j-1, ap( j1 ), 1, bp( j1 ),
185 $ 1 ) ) / bjj
186 10 CONTINUE
187 ELSE
188*
189* Compute inv(L)*A*inv(L**T)
190*
191* KK and K1K1 are the indices of A(k,k) and A(k+1,k+1)
192*
193 kk = 1
194 DO 20 k = 1, n
195 k1k1 = kk + n - k + 1
196*
197* Update the lower triangle of A(k:n,k:n)
198*
199 akk = ap( kk )
200 bkk = bp( kk )
201 akk = akk / bkk**2
202 ap( kk ) = akk
203 IF( k.LT.n ) THEN
204 CALL sscal( n-k, one / bkk, ap( kk+1 ), 1 )
205 ct = -half*akk
206 CALL saxpy( n-k, ct, bp( kk+1 ), 1, ap( kk+1 ), 1 )
207 CALL sspr2( uplo, n-k, -one, ap( kk+1 ), 1,
208 $ bp( kk+1 ), 1, ap( k1k1 ) )
209 CALL saxpy( n-k, ct, bp( kk+1 ), 1, ap( kk+1 ), 1 )
210 CALL stpsv( uplo, 'No transpose', 'Non-unit', n-k,
211 $ bp( k1k1 ), ap( kk+1 ), 1 )
212 END IF
213 kk = k1k1
214 20 CONTINUE
215 END IF
216 ELSE
217 IF( upper ) THEN
218*
219* Compute U*A*U**T
220*
221* K1 and KK are the indices of A(1,k) and A(k,k)
222*
223 kk = 0
224 DO 30 k = 1, n
225 k1 = kk + 1
226 kk = kk + k
227*
228* Update the upper triangle of A(1:k,1:k)
229*
230 akk = ap( kk )
231 bkk = bp( kk )
232 CALL stpmv( uplo, 'No transpose', 'Non-unit', k-1, bp,
233 $ ap( k1 ), 1 )
234 ct = half*akk
235 CALL saxpy( k-1, ct, bp( k1 ), 1, ap( k1 ), 1 )
236 CALL sspr2( uplo, k-1, one, ap( k1 ), 1, bp( k1 ), 1,
237 $ ap )
238 CALL saxpy( k-1, ct, bp( k1 ), 1, ap( k1 ), 1 )
239 CALL sscal( k-1, bkk, ap( k1 ), 1 )
240 ap( kk ) = akk*bkk**2
241 30 CONTINUE
242 ELSE
243*
244* Compute L**T *A*L
245*
246* JJ and J1J1 are the indices of A(j,j) and A(j+1,j+1)
247*
248 jj = 1
249 DO 40 j = 1, n
250 j1j1 = jj + n - j + 1
251*
252* Compute the j-th column of the lower triangle of A
253*
254 ajj = ap( jj )
255 bjj = bp( jj )
256 ap( jj ) = ajj*bjj + sdot( n-j, ap( jj+1 ), 1,
257 $ bp( jj+1 ), 1 )
258 CALL sscal( n-j, bjj, ap( jj+1 ), 1 )
259 CALL sspmv( uplo, n-j, one, ap( j1j1 ), bp( jj+1 ), 1,
260 $ one, ap( jj+1 ), 1 )
261 CALL stpmv( uplo, 'Transpose', 'Non-unit', n-j+1,
262 $ bp( jj ), ap( jj ), 1 )
263 jj = j1j1
264 40 CONTINUE
265 END IF
266 END IF
267 RETURN
268*
269* End of SSPGST
270*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine saxpy(n, sa, sx, incx, sy, incy)
SAXPY
Definition saxpy.f:89
real function sdot(n, sx, incx, sy, incy)
SDOT
Definition sdot.f:82
subroutine sspmv(uplo, n, alpha, ap, x, incx, beta, y, incy)
SSPMV
Definition sspmv.f:147
subroutine sspr2(uplo, n, alpha, x, incx, y, incy, ap)
SSPR2
Definition sspr2.f:142
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine stpmv(uplo, trans, diag, n, ap, x, incx)
STPMV
Definition stpmv.f:142
subroutine stpsv(uplo, trans, diag, n, ap, x, incx)
STPSV
Definition stpsv.f:144
Here is the call graph for this function:
Here is the caller graph for this function: