LAPACK  3.10.0
LAPACK: Linear Algebra PACKage
spotrf.f
Go to the documentation of this file.
1 C> \brief \b SPOTRF VARIANT: right 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 SPOTRF ( UPLO, N, A, LDA, INFO )
12 *
13 * .. Scalar Arguments ..
14 * CHARACTER UPLO
15 * INTEGER INFO, LDA, N
16 * ..
17 * .. Array Arguments ..
18 * REAL A( LDA, * )
19 * ..
20 *
21 * Purpose
22 * =======
23 *
24 C>\details \b Purpose:
25 C>\verbatim
26 C>
27 C> SPOTRF computes the Cholesky factorization of a real symmetric
28 C> positive definite matrix A.
29 C>
30 C> The factorization has the form
31 C> A = U**T * U, if UPLO = 'U', or
32 C> A = L * L**T, if UPLO = 'L',
33 C> where U is an upper triangular matrix and L is lower triangular.
34 C>
35 C> This is the right looking block version of the algorithm, calling Level 3 BLAS.
36 C>
37 C>\endverbatim
38 *
39 * Arguments:
40 * ==========
41 *
42 C> \param[in] UPLO
43 C> \verbatim
44 C> UPLO is CHARACTER*1
45 C> = 'U': Upper triangle of A is stored;
46 C> = 'L': Lower triangle of A is stored.
47 C> \endverbatim
48 C>
49 C> \param[in] N
50 C> \verbatim
51 C> N is INTEGER
52 C> The order of the matrix A. N >= 0.
53 C> \endverbatim
54 C>
55 C> \param[in,out] A
56 C> \verbatim
57 C> A is REAL array, dimension (LDA,N)
58 C> On entry, the symmetric matrix A. If UPLO = 'U', the leading
59 C> N-by-N upper triangular part of A contains the upper
60 C> triangular part of the matrix A, and the strictly lower
61 C> triangular part of A is not referenced. If UPLO = 'L', the
62 C> leading N-by-N lower triangular part of A contains the lower
63 C> triangular part of the matrix A, and the strictly upper
64 C> triangular part of A is not referenced.
65 C> \endverbatim
66 C> \verbatim
67 C> On exit, if INFO = 0, the factor U or L from the Cholesky
68 C> factorization A = U**T*U or A = L*L**T.
69 C> \endverbatim
70 C>
71 C> \param[in] LDA
72 C> \verbatim
73 C> LDA is INTEGER
74 C> The leading dimension of the array A. LDA >= max(1,N).
75 C> \endverbatim
76 C>
77 C> \param[out] INFO
78 C> \verbatim
79 C> INFO is INTEGER
80 C> = 0: successful exit
81 C> < 0: if INFO = -i, the i-th argument had an illegal value
82 C> > 0: if INFO = i, the leading minor of order i is not
83 C> positive definite, and the factorization could not be
84 C> completed.
85 C> \endverbatim
86 C>
87 *
88 * Authors:
89 * ========
90 *
91 C> \author Univ. of Tennessee
92 C> \author Univ. of California Berkeley
93 C> \author Univ. of Colorado Denver
94 C> \author NAG Ltd.
95 *
96 C> \date December 2016
97 *
98 C> \ingroup variantsPOcomputational
99 *
100 * =====================================================================
101  SUBROUTINE spotrf ( 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  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 *
241  END
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
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 spotrf(UPLO, N, A, LDA, INFO)
SPOTRF
Definition: spotrf.f:107
subroutine strsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
STRSM
Definition: strsm.f:181
subroutine ssyrk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
SSYRK
Definition: ssyrk.f:169
subroutine sgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
SGEMM
Definition: sgemm.f:187