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

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

Purpose:

 CPOTRF computes the Cholesky factorization of a real Hermitian
 positive definite matrix A.

 The factorization has the form
    A = U**H * U,  if UPLO = 'U', or
    A = L  * L**H,  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 COMPLEX array, dimension (LDA,N)
          On entry, the Hermitian 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**H*U or A = L*L**H.
[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 cpotrf.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  COMPLEX a( lda, * )
114 * ..
115 *
116 * =====================================================================
117 *
118 * .. Parameters ..
119  REAL one
120  COMPLEX cone
121  parameter ( one = 1.0e+0, cone = ( 1.0e+0, 0.0e+0 ) )
122 * ..
123 * .. Local Scalars ..
124  LOGICAL upper
125  INTEGER j, jb, nb
126 * ..
127 * .. External Functions ..
128  LOGICAL lsame
129  INTEGER ilaenv
130  EXTERNAL lsame, ilaenv
131 * ..
132 * .. External Subroutines ..
133  EXTERNAL cgemm, cpotf2, cherk, ctrsm, xerbla
134 * ..
135 * .. Intrinsic Functions ..
136  INTRINSIC max, min
137 * ..
138 * .. Executable Statements ..
139 *
140 * Test the input parameters.
141 *
142  info = 0
143  upper = lsame( uplo, 'U' )
144  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
145  info = -1
146  ELSE IF( n.LT.0 ) THEN
147  info = -2
148  ELSE IF( lda.LT.max( 1, n ) ) THEN
149  info = -4
150  END IF
151  IF( info.NE.0 ) THEN
152  CALL xerbla( 'CPOTRF', -info )
153  RETURN
154  END IF
155 *
156 * Quick return if possible
157 *
158  IF( n.EQ.0 )
159  $ RETURN
160 *
161 * Determine the block size for this environment.
162 *
163  nb = ilaenv( 1, 'CPOTRF', uplo, n, -1, -1, -1 )
164  IF( nb.LE.1 .OR. nb.GE.n ) THEN
165 *
166 * Use unblocked code.
167 *
168  CALL cpotf2( uplo, n, a, lda, info )
169  ELSE
170 *
171 * Use blocked code.
172 *
173  IF( upper ) THEN
174 *
175 * Compute the Cholesky factorization A = U'*U.
176 *
177  DO 10 j = 1, n, nb
178 *
179 * Update and factorize the current diagonal block and test
180 * for non-positive-definiteness.
181 *
182  jb = min( nb, n-j+1 )
183 
184  CALL cpotf2( 'Upper', jb, a( j, j ), lda, info )
185 
186  IF( info.NE.0 )
187  $ GO TO 30
188 
189  IF( j+jb.LE.n ) THEN
190 *
191 * Updating the trailing submatrix.
192 *
193  CALL ctrsm( 'Left', 'Upper', 'Conjugate Transpose',
194  $ 'Non-unit', jb, n-j-jb+1, cone, a( j, j ),
195  $ lda, a( j, j+jb ), lda )
196  CALL cherk( 'Upper', 'Conjugate transpose', n-j-jb+1,
197  $ jb, -one, a( j, j+jb ), lda,
198  $ one, a( j+jb, j+jb ), lda )
199  END IF
200  10 CONTINUE
201 *
202  ELSE
203 *
204 * Compute the Cholesky factorization A = L*L'.
205 *
206  DO 20 j = 1, n, nb
207 *
208 * Update and factorize the current diagonal block and test
209 * for non-positive-definiteness.
210 *
211  jb = min( nb, n-j+1 )
212 
213  CALL cpotf2( 'Lower', jb, a( j, j ), lda, info )
214 
215  IF( info.NE.0 )
216  $ GO TO 30
217 
218  IF( j+jb.LE.n ) THEN
219 *
220 * Updating the trailing submatrix.
221 *
222  CALL ctrsm( 'Right', 'Lower', 'Conjugate Transpose',
223  $ 'Non-unit', n-j-jb+1, jb, cone, a( j, j ),
224  $ lda, a( j+jb, j ), lda )
225 
226  CALL cherk( 'Lower', 'No Transpose', n-j-jb+1, jb,
227  $ -one, a( j+jb, j ), lda,
228  $ one, a( j+jb, j+jb ), lda )
229  END IF
230  20 CONTINUE
231  END IF
232  END IF
233  GO TO 40
234 *
235  30 CONTINUE
236  info = info + j - 1
237 *
238  40 CONTINUE
239  RETURN
240 *
241 * End of CPOTRF
242 *
subroutine cherk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
CHERK
Definition: cherk.f:175
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ctrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRSM
Definition: ctrsm.f:182
subroutine cpotf2(UPLO, N, A, LDA, INFO)
CPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblock...
Definition: cpotf2.f:111
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM
Definition: cgemm.f:189
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function: