1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
|
SUBROUTINE DPOTRF( UPLO, N, A, LDA, INFO )
*
* -- LAPACK routine (version 3.1) --
* Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd..
* November 2006
*
* .. Scalar Arguments ..
CHARACTER UPLO
INTEGER INFO, LDA, N
* ..
* .. Array Arguments ..
DOUBLE PRECISION A( LDA, * )
* ..
*
* 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 block version of the algorithm, calling Level 3 BLAS.
*
* Arguments
* =========
*
* UPLO (input) CHARACTER*1
* = 'U': Upper triangle of A is stored;
* = 'L': Lower triangle of A is stored.
*
* N (input) INTEGER
* The order of the matrix A. N >= 0.
*
* A (input/output) 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.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* INFO (output) 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.
*
* =====================================================================
*
* .. Parameters ..
DOUBLE PRECISION ONE
PARAMETER ( ONE = 1.0D+0 )
* ..
* .. Local Scalars ..
LOGICAL UPPER
INTEGER J, JB, NB
* ..
* .. External Functions ..
LOGICAL LSAME
INTEGER ILAENV
EXTERNAL LSAME, ILAENV
* ..
* .. External Subroutines ..
EXTERNAL DGEMM, DPOTF2, DSYRK, DTRSM, XERBLA
* ..
* .. Intrinsic Functions ..
INTRINSIC MAX, MIN
* ..
* .. Executable Statements ..
*
* Test the input parameters.
*
INFO = 0
UPPER = LSAME( UPLO, 'U' )
IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
INFO = -1
ELSE IF( N.LT.0 ) THEN
INFO = -2
ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
INFO = -4
END IF
IF( INFO.NE.0 ) THEN
CALL XERBLA( 'DPOTRF', -INFO )
RETURN
END IF
*
* Quick return if possible
*
IF( N.EQ.0 )
$ RETURN
*
* Determine the block size for this environment.
*
NB = ILAENV( 1, 'DPOTRF', UPLO, N, -1, -1, -1 )
IF( NB.LE.1 .OR. NB.GE.N ) THEN
*
* Use unblocked code.
*
CALL DPOTF2( UPLO, N, A, LDA, INFO )
ELSE
*
* Use blocked code.
*
IF( UPPER ) THEN
*
* Compute the Cholesky factorization A = U'*U.
*
DO 10 J = 1, N, NB
*
* Update and factorize the current diagonal block and test
* for non-positive-definiteness.
*
JB = MIN( NB, N-J+1 )
CALL DSYRK( 'Upper', 'Transpose', JB, J-1, -ONE,
$ A( 1, J ), LDA, ONE, A( J, J ), LDA )
CALL DPOTF2( 'Upper', JB, A( J, J ), LDA, INFO )
IF( INFO.NE.0 )
$ GO TO 30
IF( J+JB.LE.N ) THEN
*
* Compute the current block row.
*
CALL DGEMM( 'Transpose', 'No transpose', JB, N-J-JB+1,
$ J-1, -ONE, A( 1, J ), LDA, A( 1, J+JB ),
$ LDA, ONE, A( J, J+JB ), LDA )
CALL DTRSM( 'Left', 'Upper', 'Transpose', 'Non-unit',
$ JB, N-J-JB+1, ONE, A( J, J ), LDA,
$ A( J, J+JB ), LDA )
END IF
10 CONTINUE
*
ELSE
*
* Compute the Cholesky factorization A = L*L'.
*
DO 20 J = 1, N, NB
*
* Update and factorize the current diagonal block and test
* for non-positive-definiteness.
*
JB = MIN( NB, N-J+1 )
CALL DSYRK( 'Lower', 'No transpose', JB, J-1, -ONE,
$ A( J, 1 ), LDA, ONE, A( J, J ), LDA )
CALL DPOTF2( 'Lower', JB, A( J, J ), LDA, INFO )
IF( INFO.NE.0 )
$ GO TO 30
IF( J+JB.LE.N ) THEN
*
* Compute the current block column.
*
CALL DGEMM( 'No transpose', 'Transpose', N-J-JB+1, JB,
$ J-1, -ONE, A( J+JB, 1 ), LDA, A( J, 1 ),
$ LDA, ONE, A( J+JB, J ), LDA )
CALL DTRSM( 'Right', 'Lower', 'Transpose', 'Non-unit',
$ N-J-JB+1, JB, ONE, A( J, J ), LDA,
$ A( J+JB, J ), LDA )
END IF
20 CONTINUE
END IF
END IF
GO TO 40
*
30 CONTINUE
INFO = INFO + J - 1
*
40 CONTINUE
RETURN
*
* End of DPOTRF
*
END
|