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

DPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.

Purpose:

 DPOTRF 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 right looking block version of the algorithm, calling Level 3 BLAS.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]A
          A is DOUBLE PRECISION 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 = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, the leading minor of order i 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
November 2011

Definition at line 102 of file dpotrf.f.

102 *
103 * -- LAPACK computational routine (version 3.1) --
104 * -- LAPACK is a software package provided by Univ. of Tennessee, --
105 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
106 * November 2011
107 *
108 * .. Scalar Arguments ..
109  CHARACTER uplo
110  INTEGER info, lda, n
111 * ..
112 * .. Array Arguments ..
113  DOUBLE PRECISION a( lda, * )
114 * ..
115 *
116 * =====================================================================
117 *
118 * .. Parameters ..
119  DOUBLE PRECISION one
120  parameter ( one = 1.0d+0 )
121 * ..
122 * .. Local Scalars ..
123  LOGICAL upper
124  INTEGER j, jb, nb
125 * ..
126 * .. External Functions ..
127  LOGICAL lsame
128  INTEGER ilaenv
129  EXTERNAL lsame, ilaenv
130 * ..
131 * .. External Subroutines ..
132  EXTERNAL dgemm, dpotf2, dsyrk, dtrsm, xerbla
133 * ..
134 * .. Intrinsic Functions ..
135  INTRINSIC max, min
136 * ..
137 * .. Executable Statements ..
138 *
139 * Test the input parameters.
140 *
141  info = 0
142  upper = lsame( uplo, 'U' )
143  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
144  info = -1
145  ELSE IF( n.LT.0 ) THEN
146  info = -2
147  ELSE IF( lda.LT.max( 1, n ) ) THEN
148  info = -4
149  END IF
150  IF( info.NE.0 ) THEN
151  CALL xerbla( 'DPOTRF', -info )
152  RETURN
153  END IF
154 *
155 * Quick return if possible
156 *
157  IF( n.EQ.0 )
158  $ RETURN
159 *
160 * Determine the block size for this environment.
161 *
162  nb = ilaenv( 1, 'DPOTRF', uplo, n, -1, -1, -1 )
163  IF( nb.LE.1 .OR. nb.GE.n ) THEN
164 *
165 * Use unblocked code.
166 *
167  CALL dpotf2( uplo, n, a, lda, info )
168  ELSE
169 *
170 * Use blocked code.
171 *
172  IF( upper ) THEN
173 *
174 * Compute the Cholesky factorization A = U'*U.
175 *
176  DO 10 j = 1, n, nb
177 *
178 * Update and factorize the current diagonal block and test
179 * for non-positive-definiteness.
180 *
181  jb = min( nb, n-j+1 )
182 
183  CALL dpotf2( 'Upper', jb, a( j, j ), lda, info )
184 
185  IF( info.NE.0 )
186  $ GO TO 30
187 
188  IF( j+jb.LE.n ) THEN
189 *
190 * Updating the trailing submatrix.
191 *
192  CALL dtrsm( 'Left', 'Upper', 'Transpose', 'Non-unit',
193  $ jb, n-j-jb+1, one, a( j, j ), lda,
194  $ a( j, j+jb ), lda )
195  CALL dsyrk( 'Upper', 'Transpose', n-j-jb+1, jb, -one,
196  $ a( j, j+jb ), lda,
197  $ one, a( j+jb, j+jb ), lda )
198  END IF
199  10 CONTINUE
200 *
201  ELSE
202 *
203 * Compute the Cholesky factorization A = L*L'.
204 *
205  DO 20 j = 1, n, nb
206 *
207 * Update and factorize the current diagonal block and test
208 * for non-positive-definiteness.
209 *
210  jb = min( nb, n-j+1 )
211 
212  CALL dpotf2( 'Lower', jb, a( j, j ), lda, info )
213 
214  IF( info.NE.0 )
215  $ GO TO 30
216 
217  IF( j+jb.LE.n ) THEN
218 *
219 * Updating the trailing submatrix.
220 *
221  CALL dtrsm( 'Right', 'Lower', 'Transpose', 'Non-unit',
222  $ n-j-jb+1, jb, one, a( j, j ), lda,
223  $ a( j+jb, j ), lda )
224 
225  CALL dsyrk( 'Lower', 'No Transpose', n-j-jb+1, jb,
226  $ -one, a( j+jb, j ), lda,
227  $ one, a( j+jb, j+jb ), lda )
228  END IF
229  20 CONTINUE
230  END IF
231  END IF
232  GO TO 40
233 *
234  30 CONTINUE
235  info = info + j - 1
236 *
237  40 CONTINUE
238  RETURN
239 *
240 * End of DPOTRF
241 *
subroutine dtrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRSM
Definition: dtrsm.f:183
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine dsyrk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
DSYRK
Definition: dsyrk.f:171
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dpotf2(UPLO, N, A, LDA, INFO)
DPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblock...
Definition: dpotf2.f:111
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function: