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

ZPOTRF

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

Purpose:
 ZPOTRF computes the Cholesky factorization of a complex 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 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*16 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 2015

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

Here is the call graph for this function: