LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zpotrf()

subroutine zpotrf ( character  UPLO,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer  INFO 
)

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

ZPOTRF VARIANT: top-looking block version of the algorithm, calling Level 3 BLAS.

Purpose:

 ZPOTRF 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*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
December 2016

Purpose:

 ZPOTRF computes the Cholesky factorization of a real symmetric
 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 top-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*16 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**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
December 2016

Definition at line 101 of file zpotrf.f.

102*
103* -- LAPACK computational routine --
104* -- LAPACK is a software package provided by Univ. of Tennessee, --
105* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
106*
107* .. Scalar Arguments ..
108 CHARACTER UPLO
109 INTEGER INFO, LDA, N
110* ..
111* .. Array Arguments ..
112 COMPLEX*16 A( LDA, * )
113* ..
114*
115* =====================================================================
116*
117* .. Parameters ..
118 DOUBLE PRECISION ONE
119 COMPLEX*16 CONE
120 parameter( one = 1.0d+0, cone = ( 1.0d+0, 0.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 zgemm, zpotf2, zherk, ztrsm, 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( 'ZPOTRF', -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, 'ZPOTRF', uplo, n, -1, -1, -1 )
163 IF( nb.LE.1 .OR. nb.GE.n ) THEN
164*
165* Use unblocked code.
166*
167 CALL zpotf2( 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 zpotf2( '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 ztrsm( 'Left', 'Upper', 'Conjugate Transpose',
193 $ 'Non-unit', jb, n-j-jb+1, cone, a( j, j ),
194 $ lda, a( j, j+jb ), lda )
195 CALL zherk( 'Upper', 'Conjugate transpose', n-j-jb+1,
196 $ jb, -one, 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 zpotf2( '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 ztrsm( 'Right', 'Lower', 'Conjugate Transpose',
222 $ 'Non-unit', n-j-jb+1, jb, cone, a( j, j ),
223 $ lda, a( j+jb, j ), lda )
224
225 CALL zherk( '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 ZPOTRF
241*
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: ilaenv.f:162
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:187
subroutine zherk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
ZHERK
Definition: zherk.f:173
subroutine ztrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRSM
Definition: ztrsm.f:180
subroutine zpotf2(UPLO, N, A, LDA, INFO)
ZPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblock...
Definition: zpotf2.f:109
Here is the call graph for this function:
Here is the caller graph for this function: