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

DPOTRF

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

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 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 2015

Definition at line 109 of file dpotrf.f.

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