LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.
Date
November 2011

Definition at line 115 of file sspgst.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: