1 *> \brief \b DSYTRF_AASEN
3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download DSYTRF_AASEN + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsytrf_aasen.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsytrf_aasen.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsytrf_aasen.f">
21 * SUBROUTINE DSYTRF_AASEN( UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO )
23 * .. Scalar Arguments ..
25 * INTEGER N, LDA, LWORK, INFO
27 * .. Array Arguments ..
29 * DOUBLE PRECISION A( LDA, * ), WORK( * )
37 *> DSYTRF_AASEN computes the factorization of a real symmetric matrix A
38 *> using the Aasen's algorithm. The form of the factorization is
40 *> A = U*T*U**T or A = L*T*L**T
42 *> where U (or L) is a product of permutation and unit upper (lower)
43 *> triangular matrices, and T is a symmetric tridiagonal matrix.
45 *> This is the blocked version of the algorithm, calling Level 3 BLAS.
53 *> UPLO is CHARACTER*1
54 *> = 'U': Upper triangle of A is stored;
55 *> = 'L': Lower triangle of A is stored.
61 *> The order of the matrix A. N >= 0.
66 *> A is DOUBLE PRECISION array, dimension (LDA,N)
67 *> On entry, the symmetric matrix A. If UPLO = 'U', the leading
68 *> N-by-N upper triangular part of A contains the upper
69 *> triangular part of the matrix A, and the strictly lower
70 *> triangular part of A is not referenced. If UPLO = 'L', the
71 *> leading N-by-N lower triangular part of A contains the lower
72 *> triangular part of the matrix A, and the strictly upper
73 *> triangular part of A is not referenced.
75 *> On exit, the tridiagonal matrix is stored in the diagonals
76 *> and the subdiagonals of A just below (or above) the diagonals,
77 *> and L is stored below (or above) the subdiaonals, when UPLO
84 *> The leading dimension of the array A. LDA >= max(1,N).
89 *> IPIV is INTEGER array, dimension (N)
90 *> On exit, it contains the details of the interchanges, i.e.,
91 *> the row and column k of A were interchanged with the
92 *> row and column IPIV(k).
97 *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
98 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
104 *> The length of WORK. LWORK >=2*N. For optimum performance
105 *> LWORK >= N*(1+NB), where NB is the optimal blocksize.
107 *> If LWORK = -1, then a workspace query is assumed; the routine
108 *> only calculates the optimal size of the WORK array, returns
109 *> this value as the first entry of the WORK array, and no error
110 *> message related to LWORK is issued by XERBLA.
116 *> = 0: successful exit
117 *> < 0: if INFO = -i, the i-th argument had an illegal value
118 *> > 0: if INFO = i, D(i,i) is exactly zero. The factorization
119 *> has been completed, but the block diagonal matrix D is
120 *> exactly singular, and division by zero will occur if it
121 *> is used to solve a system of equations.
127 *> \author Univ. of Tennessee
128 *> \author Univ. of California Berkeley
129 *> \author Univ. of Colorado Denver
132 *> \date November 2016
134 *> \ingroup doubleSYcomputational
136 * @precisions fortran d -> s
138 * =====================================================================
139 SUBROUTINE DSYTRF_AASEN( UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO)
141 * -- LAPACK computational routine (version 3.4.0) --
142 * -- LAPACK is a software package provided by Univ. of Tennessee, --
143 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
148 * .. Scalar Arguments ..
150 INTEGER N, LDA, LWORK, INFO
152 * .. Array Arguments ..
154 DOUBLE PRECISION A( LDA, * ), WORK( * )
157 * =====================================================================
159 DOUBLE PRECISION ZERO, ONE
160 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
162 * .. Local Scalars ..
163 LOGICAL LQUERY, UPPER
164 INTEGER J, LWKOPT, IINFO
165 INTEGER NB, MJ, NJ, K1, K2, J1, J2, J3, JB
166 DOUBLE PRECISION ALPHA
168 * .. External Functions ..
171 EXTERNAL LSAME, ILAENV
173 * .. External Subroutines ..
176 * .. Intrinsic Functions ..
179 * .. Executable Statements ..
181 * Determine the block size
183 NB = ILAENV( 1, 'DSYTRF', UPLO, N, -1, -1, -1 )
185 * Test the input parameters.
188 UPPER = LSAME( UPLO, 'U' )
189 LQUERY = ( LWORK.EQ.-1 )
190 IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
192 ELSE IF( N.LT.0 ) THEN
194 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
196 ELSE IF( LWORK.LT.( 2*N ) .AND. .NOT.LQUERY ) THEN
206 CALL XERBLA( 'DSYTRF_AASEN', -INFO )
208 ELSE IF( LQUERY ) THEN
219 IF ( A( 1, 1 ).EQ.ZERO ) THEN
225 * Adjubst block size based on the workspace size
227 IF( LWORK.LT.((1+NB)*N) ) THEN
233 * .....................................................
234 * Factorize A as L*D*L**T using the upper triangle of A
235 * .....................................................
237 * Copy first row A(1, 1:N) into H(1:n) (stored in WORK(1:N))
239 CALL DCOPY( N, A( 1, 1 ), LDA, WORK( 1 ), 1 )
241 * J is the main loop index, increasing from 1 to N in steps of
242 * JB, where JB is the number of columns factorized by DLASYF;
243 * JB is either NB, or N-J+1 for the last block
250 * each step of the main loop
251 * J is the last column of the previous panel
252 * J1 is the first column of the current panel
253 * K1 identifies if the previous column of the panel has been
254 * explicitly stored, e.g., K1=1 for the first panel, and
258 JB = MIN( N-J1+1, NB )
261 * Panel factorization
263 CALL DLASYF_AASEN( UPLO, 2-K1, N-J, JB,
264 $ A( MAX(1, J), J+1 ), LDA,
265 $ IPIV( J+1 ), WORK, N, WORK( N*NB+1 ),
267 IF( (IINFO.GT.0) .AND. (INFO.EQ.0) ) THEN
271 * Ajust IPIV and apply it back (J-th step picks (J+1)-th pivot)
273 DO J2 = J+2, MIN(N, J+JB+1)
274 IPIV( J2 ) = IPIV( J2 ) + J
275 IF( (J2.NE.IPIV(J2)) .AND. ((J1-K1).GT.2) ) THEN
276 CALL DSWAP( J1-K1-2, A( 1, J2 ), 1,
277 $ A( 1, IPIV(J2) ), 1 )
282 * Trailing submatrix update, where
283 * the row A(J1-1, J2-1:N) stores U(J1, J2+1:N) and
284 * WORK stores the current block of the auxiriarly matrix H
288 * If first panel and JB=1 (NB=1), then nothing to do
290 IF( J1.GT.1 .OR. JB.GT.1 ) THEN
292 * Merge rank-1 update with BLAS-3 update
296 CALL DCOPY( N-J, A( J-1, J+1 ), LDA,
297 $ WORK( (J+1-J1+1)+JB*N ), 1 )
298 CALL DSCAL( N-J, ALPHA, WORK( (J+1-J1+1)+JB*N ), 1 )
300 * K1 identifies if the previous column of the panel has been
301 * explicitly stored, e.g., K1=1 and K2= 0 for the first panel,
302 * while K1=0 and K2=1 for the rest
315 * First update skips the first column
321 NJ = MIN( NB, N-J2+1 )
323 * Update (J2, J2) diagonal block with DGEMV
327 CALL DGEMV( 'No transpose', MJ, JB+1,
328 $ -ONE, WORK( J3-J1+1+K1*N ), N,
330 $ ONE, A( J3, J3 ), LDA )
334 * Update off-diagonal block of J2-th block row with DGEMM
336 CALL DGEMM( 'Transpose', 'Transpose',
338 $ -ONE, A( J1-K2, J2 ), LDA,
339 $ WORK( J3-J1+1+K1*N ), N,
340 $ ONE, A( J2, J3 ), LDA )
343 * Recover T( J, J+1 )
348 * WORK(J+1, 1) stores H(J+1, 1)
350 CALL DCOPY( N-J, A( J+1, J+1 ), LDA, WORK( 1 ), 1 )
355 * .....................................................
356 * Factorize A as L*D*L**T using the lower triangle of A
357 * .....................................................
359 * copy first column A(1:N, 1) into H(1:N, 1)
360 * (stored in WORK(1:N))
362 CALL DCOPY( N, A( 1, 1 ), 1, WORK( 1 ), 1 )
364 * J is the main loop index, increasing from 1 to N in steps of
365 * JB, where JB is the number of columns factorized by DLASYF;
366 * JB is either NB, or N-J+1 for the last block
373 * each step of the main loop
374 * J is the last column of the previous panel
375 * J1 is the first column of the current panel
376 * K1 identifies if the previous column of the panel has been
377 * explicitly stored, e.g., K1=1 for the first panel, and
381 JB = MIN( N-J1+1, NB )
384 * Panel factorization
386 CALL DLASYF_AASEN( UPLO, 2-K1, N-J, JB,
387 $ A( J+1, MAX(1, J) ), LDA,
388 $ IPIV( J+1 ), WORK, N, WORK( N*NB+1 ), IINFO)
389 IF( (IINFO.GT.0) .AND. (INFO.EQ.0) ) THEN
393 * Ajust IPIV and apply it back (J-th step picks (J+1)-th pivot)
395 DO J2 = J+2, MIN(N, J+JB+1)
396 IPIV( J2 ) = IPIV( J2 ) + J
397 IF( (J2.NE.IPIV(J2)) .AND. ((J1-K1).GT.2) ) THEN
398 CALL DSWAP( J1-K1-2, A( J2, 1 ), LDA,
399 $ A( IPIV(J2), 1 ), LDA )
404 * Trailing submatrix update, where
405 * A(J2+1, J1-1) stores L(J2+1, J1) and
406 * WORK(J2+1, 1) stores H(J2+1, 1)
410 * if first panel and JB=1 (NB=1), then nothing to do
412 IF( J1.GT.1 .OR. JB.GT.1 ) THEN
414 * Merge rank-1 update with BLAS-3 update
418 CALL DCOPY( N-J, A( J+1, J-1 ), 1,
419 $ WORK( (J+1-J1+1)+JB*N ), 1 )
420 CALL DSCAL( N-J, ALPHA, WORK( (J+1-J1+1)+JB*N ), 1 )
422 * K1 identifies if the previous column of the panel has been
423 * explicitly stored, e.g., K1=1 and K2= 0 for the first panel,
424 * while K1=0 and K2=1 for the rest
437 * First update skips the first column
443 NJ = MIN( NB, N-J2+1 )
445 * Update (J2, J2) diagonal block with DGEMV
449 CALL DGEMV( 'No transpose', MJ, JB+1,
450 $ -ONE, WORK( J3-J1+1+K1*N ), N,
451 $ A( J3, J1-K2 ), LDA,
452 $ ONE, A( J3, J3 ), 1 )
456 * Update off-diagonal block in J2-th block column with DGEMM
458 CALL DGEMM( 'No transpose', 'Transpose',
460 $ -ONE, WORK( J3-J1+1+K1*N ), N,
461 $ A( J2, J1-K2 ), LDA,
462 $ ONE, A( J3, J2 ), LDA )
465 * Recover T( J+1, J )
470 * WORK(J+1, 1) stores H(J+1, 1)
472 CALL DCOPY( N-J, A( J+1, J+1 ), 1, WORK( 1 ), 1 )
480 * End of DSYTRF_AASEN