LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine spotf2 ( character  UPLO,
integer  N,
real, dimension( lda, * )  A,
integer  LDA,
integer  INFO 
)

SPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblocked algorithm).

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

Purpose:
 SPOTF2 computes the Cholesky factorization of a real symmetric
 positive definite matrix A.

 The factorization has the form
    A = U**T * U ,  if UPLO = 'U', or
    A = L  * L**T,  if UPLO = 'L',
 where U is an upper triangular matrix and L is lower triangular.

 This is the unblocked version of the algorithm, calling Level 2 BLAS.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the upper or lower triangular part of the
          symmetric matrix A is stored.
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[in]N
          N is INTEGER
          The order of the matrix A.  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 factor U or L from the Cholesky
          factorization A = U**T *U  or A = L*L**T.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -k, the k-th argument had an illegal value
          > 0: if INFO = k, the leading minor of order k is not
               positive definite, and the factorization could not be
               completed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 111 of file spotf2.f.

111 *
112 * -- LAPACK computational routine (version 3.4.2) --
113 * -- LAPACK is a software package provided by Univ. of Tennessee, --
114 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
115 * September 2012
116 *
117 * .. Scalar Arguments ..
118  CHARACTER uplo
119  INTEGER info, lda, n
120 * ..
121 * .. Array Arguments ..
122  REAL a( lda, * )
123 * ..
124 *
125 * =====================================================================
126 *
127 * .. Parameters ..
128  REAL one, zero
129  parameter ( one = 1.0e+0, zero = 0.0e+0 )
130 * ..
131 * .. Local Scalars ..
132  LOGICAL upper
133  INTEGER j
134  REAL ajj
135 * ..
136 * .. External Functions ..
137  LOGICAL lsame, sisnan
138  REAL sdot
139  EXTERNAL lsame, sdot, sisnan
140 * ..
141 * .. External Subroutines ..
142  EXTERNAL sgemv, sscal, xerbla
143 * ..
144 * .. Intrinsic Functions ..
145  INTRINSIC max, sqrt
146 * ..
147 * .. Executable Statements ..
148 *
149 * Test the input parameters.
150 *
151  info = 0
152  upper = lsame( uplo, 'U' )
153  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
154  info = -1
155  ELSE IF( n.LT.0 ) THEN
156  info = -2
157  ELSE IF( lda.LT.max( 1, n ) ) THEN
158  info = -4
159  END IF
160  IF( info.NE.0 ) THEN
161  CALL xerbla( 'SPOTF2', -info )
162  RETURN
163  END IF
164 *
165 * Quick return if possible
166 *
167  IF( n.EQ.0 )
168  $ RETURN
169 *
170  IF( upper ) THEN
171 *
172 * Compute the Cholesky factorization A = U**T *U.
173 *
174  DO 10 j = 1, n
175 *
176 * Compute U(J,J) and test for non-positive-definiteness.
177 *
178  ajj = a( j, j ) - sdot( j-1, a( 1, j ), 1, a( 1, j ), 1 )
179  IF( ajj.LE.zero.OR.sisnan( ajj ) ) THEN
180  a( j, j ) = ajj
181  GO TO 30
182  END IF
183  ajj = sqrt( ajj )
184  a( j, j ) = ajj
185 *
186 * Compute elements J+1:N of row J.
187 *
188  IF( j.LT.n ) THEN
189  CALL sgemv( 'Transpose', j-1, n-j, -one, a( 1, j+1 ),
190  $ lda, a( 1, j ), 1, one, a( j, j+1 ), lda )
191  CALL sscal( n-j, one / ajj, a( j, j+1 ), lda )
192  END IF
193  10 CONTINUE
194  ELSE
195 *
196 * Compute the Cholesky factorization A = L*L**T.
197 *
198  DO 20 j = 1, n
199 *
200 * Compute L(J,J) and test for non-positive-definiteness.
201 *
202  ajj = a( j, j ) - sdot( j-1, a( j, 1 ), lda, a( j, 1 ),
203  $ lda )
204  IF( ajj.LE.zero.OR.sisnan( ajj ) ) THEN
205  a( j, j ) = ajj
206  GO TO 30
207  END IF
208  ajj = sqrt( ajj )
209  a( j, j ) = ajj
210 *
211 * Compute elements J+1:N of column J.
212 *
213  IF( j.LT.n ) THEN
214  CALL sgemv( 'No transpose', n-j, j-1, -one, a( j+1, 1 ),
215  $ lda, a( j, 1 ), lda, one, a( j+1, j ), 1 )
216  CALL sscal( n-j, one / ajj, a( j+1, j ), 1 )
217  END IF
218  20 CONTINUE
219  END IF
220  GO TO 40
221 *
222  30 CONTINUE
223  info = j
224 *
225  40 CONTINUE
226  RETURN
227 *
228 * End of SPOTF2
229 *
real function sdot(N, SX, INCX, SY, INCY)
SDOT
Definition: sdot.f:53
logical function sisnan(SIN)
SISNAN tests input for NaN.
Definition: sisnan.f:61
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine sgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
SGEMV
Definition: sgemv.f:158
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: