OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
Variants Computational routines

Functions

subroutine cpotrf (uplo, n, a, lda, info)
 CPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.
subroutine dpotrf (uplo, n, a, lda, info)
 DPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.
subroutine spotrf (uplo, n, a, lda, info)
 SPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.
subroutine zpotrf (uplo, n, a, lda, info)
 ZPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.

Detailed Description

This is the group of Variants Computational routines

Function Documentation

◆ cpotrf()

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.

CPOTRF VARIANT: top-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
December 2016

Purpose:

!>
!> CPOTRF 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 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 cpotrf.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 A( LDA, * )
113* ..
114*
115* =====================================================================
116*
117* .. Parameters ..
118 REAL ONE
119 COMPLEX CONE
120 parameter( one = 1.0e+0, cone = ( 1.0e+0, 0.0e+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 cgemm, cpotf2, cherk, ctrsm, 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( 'CPOTRF', -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, 'CPOTRF', uplo, n, -1, -1, -1 )
163 IF( nb.LE.1 .OR. nb.GE.n ) THEN
164*
165* Use unblocked code.
166*
167 CALL cpotf2( 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 cpotf2( '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 ctrsm( '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 cherk( '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 cpotf2( '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 ctrsm( '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 cherk( '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 CPOTRF
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 cpotf2(uplo, n, a, lda, info)
CPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblock...
Definition cpotf2.f:109
subroutine ctrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
CTRSM
Definition ctrsm.f:180
subroutine cherk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
CHERK
Definition cherk.f:173
subroutine cgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CGEMM
Definition cgemm.f:187
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21

◆ dpotrf()

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.

DPOTRF VARIANT: top-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
December 2016

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

Definition at line 101 of file dpotrf.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 DOUBLE PRECISION A( LDA, * )
113* ..
114*
115* =====================================================================
116*
117* .. Parameters ..
118 DOUBLE PRECISION ONE
119 parameter( one = 1.0d+0 )
120* ..
121* .. Local Scalars ..
122 LOGICAL UPPER
123 INTEGER J, JB, NB
124* ..
125* .. External Functions ..
126 LOGICAL LSAME
127 INTEGER ILAENV
128 EXTERNAL lsame, ilaenv
129* ..
130* .. External Subroutines ..
131 EXTERNAL dgemm, dpotf2, dsyrk, dtrsm, xerbla
132* ..
133* .. Intrinsic Functions ..
134 INTRINSIC max, min
135* ..
136* .. Executable Statements ..
137*
138* Test the input parameters.
139*
140 info = 0
141 upper = lsame( uplo, 'U' )
142 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
143 info = -1
144 ELSE IF( n.LT.0 ) THEN
145 info = -2
146 ELSE IF( lda.LT.max( 1, n ) ) THEN
147 info = -4
148 END IF
149 IF( info.NE.0 ) THEN
150 CALL xerbla( 'DPOTRF', -info )
151 RETURN
152 END IF
153*
154* Quick return if possible
155*
156 IF( n.EQ.0 )
157 $ RETURN
158*
159* Determine the block size for this environment.
160*
161 nb = ilaenv( 1, 'DPOTRF', uplo, n, -1, -1, -1 )
162 IF( nb.LE.1 .OR. nb.GE.n ) THEN
163*
164* Use unblocked code.
165*
166 CALL dpotf2( uplo, n, a, lda, info )
167 ELSE
168*
169* Use blocked code.
170*
171 IF( upper ) THEN
172*
173* Compute the Cholesky factorization A = U'*U.
174*
175 DO 10 j = 1, n, nb
176*
177* Update and factorize the current diagonal block and test
178* for non-positive-definiteness.
179*
180 jb = min( nb, n-j+1 )
181
182 CALL dpotf2( 'Upper', jb, a( j, j ), lda, info )
183
184 IF( info.NE.0 )
185 $ GO TO 30
186
187 IF( j+jb.LE.n ) THEN
188*
189* Updating the trailing submatrix.
190*
191 CALL dtrsm( 'Left', 'Upper', 'Transpose', 'Non-unit',
192 $ jb, n-j-jb+1, one, a( j, j ), lda,
193 $ a( j, j+jb ), lda )
194 CALL dsyrk( 'Upper', 'Transpose', n-j-jb+1, jb, -one,
195 $ a( j, j+jb ), lda,
196 $ one, a( j+jb, j+jb ), lda )
197 END IF
198 10 CONTINUE
199*
200 ELSE
201*
202* Compute the Cholesky factorization A = L*L'.
203*
204 DO 20 j = 1, n, nb
205*
206* Update and factorize the current diagonal block and test
207* for non-positive-definiteness.
208*
209 jb = min( nb, n-j+1 )
210
211 CALL dpotf2( 'Lower', jb, a( j, j ), lda, info )
212
213 IF( info.NE.0 )
214 $ GO TO 30
215
216 IF( j+jb.LE.n ) THEN
217*
218* Updating the trailing submatrix.
219*
220 CALL dtrsm( 'Right', 'Lower', 'Transpose', 'Non-unit',
221 $ n-j-jb+1, jb, one, a( j, j ), lda,
222 $ a( j+jb, j ), lda )
223
224 CALL dsyrk( 'Lower', 'No Transpose', n-j-jb+1, jb,
225 $ -one, a( j+jb, j ), lda,
226 $ one, a( j+jb, j+jb ), lda )
227 END IF
228 20 CONTINUE
229 END IF
230 END IF
231 GO TO 40
232*
233 30 CONTINUE
234 info = info + j - 1
235*
236 40 CONTINUE
237 RETURN
238*
239* End of DPOTRF
240*
subroutine dpotf2(uplo, n, a, lda, info)
DPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblock...
Definition dpotf2.f:109
subroutine dsyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
DSYRK
Definition dsyrk.f:169
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:187
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181

◆ spotrf()

subroutine spotrf ( character uplo,
integer n,
real, dimension( lda, * ) a,
integer lda,
integer info )

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

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

Purpose:

!>
!> SPOTRF 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 REAL 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
December 2016

Purpose:

!>
!> SPOTRF 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 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 REAL 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
December 2016

Definition at line 101 of file spotrf.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 REAL A( LDA, * )
113* ..
114*
115* =====================================================================
116*
117* .. Parameters ..
118 REAL ONE
119 parameter( one = 1.0e+0 )
120* ..
121* .. Local Scalars ..
122 LOGICAL UPPER
123 INTEGER J, JB, NB
124* ..
125* .. External Functions ..
126 LOGICAL LSAME
127 INTEGER ILAENV
128 EXTERNAL lsame, ilaenv
129* ..
130* .. External Subroutines ..
131 EXTERNAL sgemm, spotf2, ssyrk, strsm, xerbla
132* ..
133* .. Intrinsic Functions ..
134 INTRINSIC max, min
135* ..
136* .. Executable Statements ..
137*
138* Test the input parameters.
139*
140 info = 0
141 upper = lsame( uplo, 'U' )
142 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
143 info = -1
144 ELSE IF( n.LT.0 ) THEN
145 info = -2
146 ELSE IF( lda.LT.max( 1, n ) ) THEN
147 info = -4
148 END IF
149 IF( info.NE.0 ) THEN
150 CALL xerbla( 'SPOTRF', -info )
151 RETURN
152 END IF
153*
154* Quick return if possible
155*
156 IF( n.EQ.0 )
157 $ RETURN
158*
159* Determine the block size for this environment.
160*
161 nb = ilaenv( 1, 'SPOTRF', uplo, n, -1, -1, -1 )
162 IF( nb.LE.1 .OR. nb.GE.n ) THEN
163*
164* Use unblocked code.
165*
166 CALL spotf2( uplo, n, a, lda, info )
167 ELSE
168*
169* Use blocked code.
170*
171 IF( upper ) THEN
172*
173* Compute the Cholesky factorization A = U'*U.
174*
175 DO 10 j = 1, n, nb
176*
177* Update and factorize the current diagonal block and test
178* for non-positive-definiteness.
179*
180 jb = min( nb, n-j+1 )
181
182 CALL spotf2( 'Upper', jb, a( j, j ), lda, info )
183
184 IF( info.NE.0 )
185 $ GO TO 30
186
187 IF( j+jb.LE.n ) THEN
188*
189* Updating the trailing submatrix.
190*
191 CALL strsm( 'Left', 'Upper', 'Transpose', 'Non-unit',
192 $ jb, n-j-jb+1, one, a( j, j ), lda,
193 $ a( j, j+jb ), lda )
194 CALL ssyrk( 'Upper', 'Transpose', n-j-jb+1, jb, -one,
195 $ a( j, j+jb ), lda,
196 $ one, a( j+jb, j+jb ), lda )
197 END IF
198 10 CONTINUE
199*
200 ELSE
201*
202* Compute the Cholesky factorization A = L*L'.
203*
204 DO 20 j = 1, n, nb
205*
206* Update and factorize the current diagonal block and test
207* for non-positive-definiteness.
208*
209 jb = min( nb, n-j+1 )
210
211 CALL spotf2( 'Lower', jb, a( j, j ), lda, info )
212
213 IF( info.NE.0 )
214 $ GO TO 30
215
216 IF( j+jb.LE.n ) THEN
217*
218* Updating the trailing submatrix.
219*
220 CALL strsm( 'Right', 'Lower', 'Transpose', 'Non-unit',
221 $ n-j-jb+1, jb, one, a( j, j ), lda,
222 $ a( j+jb, j ), lda )
223
224 CALL ssyrk( 'Lower', 'No Transpose', n-j-jb+1, jb,
225 $ -one, a( j+jb, j ), lda,
226 $ one, a( j+jb, j+jb ), lda )
227 END IF
228 20 CONTINUE
229 END IF
230 END IF
231 GO TO 40
232*
233 30 CONTINUE
234 info = info + j - 1
235*
236 40 CONTINUE
237 RETURN
238*
239* End of SPOTRF
240*
subroutine spotf2(uplo, n, a, lda, info)
SPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblock...
Definition spotf2.f:109
subroutine ssyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
SSYRK
Definition ssyrk.f:169
subroutine strsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRSM
Definition strsm.f:181
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:187

◆ 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*
subroutine zpotf2(uplo, n, a, lda, info)
ZPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblock...
Definition zpotf2.f:109
subroutine zherk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
ZHERK
Definition zherk.f:173
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:187
subroutine ztrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRSM
Definition ztrsm.f:180