OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
cpotrf.f
Go to the documentation of this file.
1C> \brief \b CPOTRF VARIANT: top-looking block version of the algorithm, calling Level 3 BLAS.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE CPOTRF ( UPLO, N, A, LDA, INFO )
12*
13* .. Scalar Arguments ..
14* CHARACTER UPLO
15* INTEGER INFO, LDA, N
16* ..
17* .. Array Arguments ..
18* COMPLEX A( LDA, * )
19* ..
20*
21* Purpose
22* =======
23*
24C>\details \b Purpose:
25C>\verbatim
26C>
27C> CPOTRF computes the Cholesky factorization of a real symmetric
28C> positive definite matrix A.
29C>
30C> The factorization has the form
31C> A = U**H * U, if UPLO = 'U', or
32C> A = L * L**H, if UPLO = 'L',
33C> where U is an upper triangular matrix and L is lower triangular.
34C>
35C> This is the top-looking block version of the algorithm, calling Level 3 BLAS.
36C>
37C>\endverbatim
38*
39* Arguments:
40* ==========
41*
42C> \param[in] UPLO
43C> \verbatim
44C> UPLO is CHARACTER*1
45C> = 'U': Upper triangle of A is stored;
46C> = 'L': Lower triangle of A is stored.
47C> \endverbatim
48C>
49C> \param[in] N
50C> \verbatim
51C> N is INTEGER
52C> The order of the matrix A. N >= 0.
53C> \endverbatim
54C>
55C> \param[in,out] A
56C> \verbatim
57C> A is COMPLEX array, dimension (LDA,N)
58C> On entry, the symmetric matrix A. If UPLO = 'U', the leading
59C> N-by-N upper triangular part of A contains the upper
60C> triangular part of the matrix A, and the strictly lower
61C> triangular part of A is not referenced. If UPLO = 'L', the
62C> leading N-by-N lower triangular part of A contains the lower
63C> triangular part of the matrix A, and the strictly upper
64C> triangular part of A is not referenced.
65C> \endverbatim
66C> \verbatim
67C> On exit, if INFO = 0, the factor U or L from the Cholesky
68C> factorization A = U**H*U or A = L*L**H.
69C> \endverbatim
70C>
71C> \param[in] LDA
72C> \verbatim
73C> LDA is INTEGER
74C> The leading dimension of the array A. LDA >= max(1,N).
75C> \endverbatim
76C>
77C> \param[out] INFO
78C> \verbatim
79C> INFO is INTEGER
80C> = 0: successful exit
81C> < 0: if INFO = -i, the i-th argument had an illegal value
82C> > 0: if INFO = i, the leading minor of order i is not
83C> positive definite, and the factorization could not be
84C> completed.
85C> \endverbatim
86C>
87*
88* Authors:
89* ========
90*
91C> \author Univ. of Tennessee
92C> \author Univ. of California Berkeley
93C> \author Univ. of Colorado Denver
94C> \author NAG Ltd.
95*
96C> \date December 2016
97*
98C> \ingroup variantsPOcomputational
99*
100* =====================================================================
101 SUBROUTINE cpotrf ( UPLO, N, A, LDA, INFO )
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.LT. ELSE IF( N0 ) THEN
146 INFO = -2
147.LT. ELSE IF( LDAMAX( 1, N ) ) THEN
148 INFO = -4
149 END IF
150.NE. IF( INFO0 ) THEN
151 CALL XERBLA( 'cpotrf', -INFO )
152 RETURN
153 END IF
154*
155* Quick return if possible
156*
157.EQ. IF( N0 )
158 $ RETURN
159*
160* Determine the block size for this environment.
161*
162 NB = ILAENV( 1, 'cpotrf', UPLO, N, -1, -1, -1 )
163.LE..OR..GE. IF( NB1 NBN ) 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 JB = MIN( NB, N-J+1 )
179*
180* Compute the current block.
181*
182 CALL CTRSM( 'left', 'upper', 'conjugate transpose',
183 $ 'non-unit', J-1, JB, CONE, A( 1, 1 ), LDA,
184 $ A( 1, J ), LDA )
185
186 CALL CHERK( 'upper', 'conjugate transpose', JB, J-1,
187 $ -ONE, A( 1, J ), LDA, ONE, A( J, J ), LDA )
188*
189* Update and factorize the current diagonal block and test
190* for non-positive-definiteness.
191*
192 CALL CPOTF2( 'upper', JB, A( J, J ), LDA, INFO )
193.NE. IF( INFO0 )
194 $ GO TO 30
195
196 10 CONTINUE
197*
198 ELSE
199*
200* Compute the Cholesky factorization A = L*L'.
201*
202 DO 20 J = 1, N, NB
203
204 JB = MIN( NB, N-J+1 )
205*
206* Compute the current block.
207*
208 CALL CTRSM( 'right', 'lower', 'conjugate transpose',
209 $ 'non-unit', JB, J-1, CONE, A( 1, 1 ), LDA,
210 $ A( J, 1 ), LDA )
211
212 CALL CHERK( 'lower', 'no transpose', JB, J-1,
213 $ -ONE, A( J, 1 ), LDA,
214 $ ONE, A( J, J ), LDA )
215*
216* Update and factorize the current diagonal block and test
217* for non-positive-definiteness.
218*
219 CALL CPOTF2( 'lower', JB, A( J, J ), LDA, INFO )
220.NE. IF( INFO0 )
221 $ GO TO 30
222
223 20 CONTINUE
224 END IF
225 END IF
226 GO TO 40
227*
228 30 CONTINUE
229 INFO = INFO + J - 1
230*
231 40 CONTINUE
232 RETURN
233*
234* End of CPOTRF
235*
236 END
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 cpotrf(uplo, n, a, lda, info)
CPOTRF
Definition cpotrf.f:107
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