LAPACK  3.7.1
LAPACK: Linear Algebra PACKage

◆ dsytrf_aa()

subroutine dsytrf_aa ( character  UPLO,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
integer, dimension( * )  IPIV,
double precision, dimension( * )  WORK,
integer  LWORK,
integer  INFO 
)

DSYTRF_AA

Download DSYTRF_AA + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 DSYTRF_AA computes the factorization of a real symmetric matrix A
 using the Aasen's algorithm.  The form of the factorization is

    A = U*T*U**T  or  A = L*T*L**T

 where U (or L) is a product of permutation and unit upper (lower)
 triangular matrices, and T is a symmetric tridiagonal matrix.

 This is the blocked 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, the tridiagonal matrix is stored in the diagonals
          and the subdiagonals of A just below (or above) the diagonals,
          and L is stored below (or above) the subdiaonals, when UPLO
          is 'L' (or 'U').
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]IPIV
          IPIV is INTEGER array, dimension (N)
          On exit, it contains the details of the interchanges, i.e.,
          the row and column k of A were interchanged with the
          row and column IPIV(k).
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The length of WORK.  LWORK >= MAX(1,2*N). For optimum performance
          LWORK >= N*(1+NB), where NB is the optimal blocksize.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
December 2016

Definition at line 134 of file dsytrf_aa.f.

134 *
135 * -- LAPACK computational routine (version 3.7.0) --
136 * -- LAPACK is a software package provided by Univ. of Tennessee, --
137 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
138 * December 2016
139 *
140  IMPLICIT NONE
141 *
142 * .. Scalar Arguments ..
143  CHARACTER uplo
144  INTEGER n, lda, lwork, info
145 * ..
146 * .. Array Arguments ..
147  INTEGER ipiv( * )
148  DOUBLE PRECISION a( lda, * ), work( * )
149 * ..
150 *
151 * =====================================================================
152 * .. Parameters ..
153  DOUBLE PRECISION zero, one
154  parameter( zero = 0.0d+0, one = 1.0d+0 )
155 *
156 * .. Local Scalars ..
157  LOGICAL lquery, upper
158  INTEGER j, lwkopt
159  INTEGER nb, mj, nj, k1, k2, j1, j2, j3, jb
160  DOUBLE PRECISION alpha
161 * ..
162 * .. External Functions ..
163  LOGICAL lsame
164  INTEGER ilaenv
165  EXTERNAL lsame, ilaenv
166 * ..
167 * .. External Subroutines ..
168  EXTERNAL xerbla
169 * ..
170 * .. Intrinsic Functions ..
171  INTRINSIC max
172 * ..
173 * .. Executable Statements ..
174 *
175 * Determine the block size
176 *
177  nb = ilaenv( 1, 'DSYTRF', uplo, n, -1, -1, -1 )
178 *
179 * Test the input parameters.
180 *
181  info = 0
182  upper = lsame( uplo, 'U' )
183  lquery = ( lwork.EQ.-1 )
184  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
185  info = -1
186  ELSE IF( n.LT.0 ) THEN
187  info = -2
188  ELSE IF( lda.LT.max( 1, n ) ) THEN
189  info = -4
190  ELSE IF( lwork.LT.max( 1, 2*n ) .AND. .NOT.lquery ) THEN
191  info = -7
192  END IF
193 *
194  IF( info.EQ.0 ) THEN
195  lwkopt = (nb+1)*n
196  work( 1 ) = lwkopt
197  END IF
198 *
199  IF( info.NE.0 ) THEN
200  CALL xerbla( 'DSYTRF_AA', -info )
201  RETURN
202  ELSE IF( lquery ) THEN
203  RETURN
204  END IF
205 *
206 * Quick return
207 *
208  IF ( n.EQ.0 ) THEN
209  RETURN
210  ENDIF
211  ipiv( 1 ) = 1
212  IF ( n.EQ.1 ) THEN
213  RETURN
214  END IF
215 *
216 * Adjust block size based on the workspace size
217 *
218  IF( lwork.LT.((1+nb)*n) ) THEN
219  nb = ( lwork-n ) / n
220  END IF
221 *
222  IF( upper ) THEN
223 *
224 * .....................................................
225 * Factorize A as L*D*L**T using the upper triangle of A
226 * .....................................................
227 *
228 * Copy first row A(1, 1:N) into H(1:n) (stored in WORK(1:N))
229 *
230  CALL dcopy( n, a( 1, 1 ), lda, work( 1 ), 1 )
231 *
232 * J is the main loop index, increasing from 1 to N in steps of
233 * JB, where JB is the number of columns factorized by DLASYF;
234 * JB is either NB, or N-J+1 for the last block
235 *
236  j = 0
237  10 CONTINUE
238  IF( j.GE.n )
239  $ GO TO 20
240 *
241 * each step of the main loop
242 * J is the last column of the previous panel
243 * J1 is the first column of the current panel
244 * K1 identifies if the previous column of the panel has been
245 * explicitly stored, e.g., K1=1 for the first panel, and
246 * K1=0 for the rest
247 *
248  j1 = j + 1
249  jb = min( n-j1+1, nb )
250  k1 = max(1, j)-j
251 *
252 * Panel factorization
253 *
254  CALL dlasyf_aa( uplo, 2-k1, n-j, jb,
255  $ a( max(1, j), j+1 ), lda,
256  $ ipiv( j+1 ), work, n, work( n*nb+1 ) )
257 *
258 * Ajust IPIV and apply it back (J-th step picks (J+1)-th pivot)
259 *
260  DO j2 = j+2, min(n, j+jb+1)
261  ipiv( j2 ) = ipiv( j2 ) + j
262  IF( (j2.NE.ipiv(j2)) .AND. ((j1-k1).GT.2) ) THEN
263  CALL dswap( j1-k1-2, a( 1, j2 ), 1,
264  $ a( 1, ipiv(j2) ), 1 )
265  END IF
266  END DO
267  j = j + jb
268 *
269 * Trailing submatrix update, where
270 * the row A(J1-1, J2-1:N) stores U(J1, J2+1:N) and
271 * WORK stores the current block of the auxiriarly matrix H
272 *
273  IF( j.LT.n ) THEN
274 *
275 * If first panel and JB=1 (NB=1), then nothing to do
276 *
277  IF( j1.GT.1 .OR. jb.GT.1 ) THEN
278 *
279 * Merge rank-1 update with BLAS-3 update
280 *
281  alpha = a( j, j+1 )
282  a( j, j+1 ) = one
283  CALL dcopy( n-j, a( j-1, j+1 ), lda,
284  $ work( (j+1-j1+1)+jb*n ), 1 )
285  CALL dscal( n-j, alpha, work( (j+1-j1+1)+jb*n ), 1 )
286 *
287 * K1 identifies if the previous column of the panel has been
288 * explicitly stored, e.g., K1=1 and K2= 0 for the first panel,
289 * while K1=0 and K2=1 for the rest
290 *
291  IF( j1.GT.1 ) THEN
292 *
293 * Not first panel
294 *
295  k2 = 1
296  ELSE
297 *
298 * First panel
299 *
300  k2 = 0
301 *
302 * First update skips the first column
303 *
304  jb = jb - 1
305  END IF
306 *
307  DO j2 = j+1, n, nb
308  nj = min( nb, n-j2+1 )
309 *
310 * Update (J2, J2) diagonal block with DGEMV
311 *
312  j3 = j2
313  DO mj = nj-1, 1, -1
314  CALL dgemv( 'No transpose', mj, jb+1,
315  $ -one, work( j3-j1+1+k1*n ), n,
316  $ a( j1-k2, j3 ), 1,
317  $ one, a( j3, j3 ), lda )
318  j3 = j3 + 1
319  END DO
320 *
321 * Update off-diagonal block of J2-th block row with DGEMM
322 *
323  CALL dgemm( 'Transpose', 'Transpose',
324  $ nj, n-j3+1, jb+1,
325  $ -one, a( j1-k2, j2 ), lda,
326  $ work( j3-j1+1+k1*n ), n,
327  $ one, a( j2, j3 ), lda )
328  END DO
329 *
330 * Recover T( J, J+1 )
331 *
332  a( j, j+1 ) = alpha
333  END IF
334 *
335 * WORK(J+1, 1) stores H(J+1, 1)
336 *
337  CALL dcopy( n-j, a( j+1, j+1 ), lda, work( 1 ), 1 )
338  END IF
339  GO TO 10
340  ELSE
341 *
342 * .....................................................
343 * Factorize A as L*D*L**T using the lower triangle of A
344 * .....................................................
345 *
346 * copy first column A(1:N, 1) into H(1:N, 1)
347 * (stored in WORK(1:N))
348 *
349  CALL dcopy( n, a( 1, 1 ), 1, work( 1 ), 1 )
350 *
351 * J is the main loop index, increasing from 1 to N in steps of
352 * JB, where JB is the number of columns factorized by DLASYF;
353 * JB is either NB, or N-J+1 for the last block
354 *
355  j = 0
356  11 CONTINUE
357  IF( j.GE.n )
358  $ GO TO 20
359 *
360 * each step of the main loop
361 * J is the last column of the previous panel
362 * J1 is the first column of the current panel
363 * K1 identifies if the previous column of the panel has been
364 * explicitly stored, e.g., K1=1 for the first panel, and
365 * K1=0 for the rest
366 *
367  j1 = j+1
368  jb = min( n-j1+1, nb )
369  k1 = max(1, j)-j
370 *
371 * Panel factorization
372 *
373  CALL dlasyf_aa( uplo, 2-k1, n-j, jb,
374  $ a( j+1, max(1, j) ), lda,
375  $ ipiv( j+1 ), work, n, work( n*nb+1 ) )
376 *
377 * Ajust IPIV and apply it back (J-th step picks (J+1)-th pivot)
378 *
379  DO j2 = j+2, min(n, j+jb+1)
380  ipiv( j2 ) = ipiv( j2 ) + j
381  IF( (j2.NE.ipiv(j2)) .AND. ((j1-k1).GT.2) ) THEN
382  CALL dswap( j1-k1-2, a( j2, 1 ), lda,
383  $ a( ipiv(j2), 1 ), lda )
384  END IF
385  END DO
386  j = j + jb
387 *
388 * Trailing submatrix update, where
389 * A(J2+1, J1-1) stores L(J2+1, J1) and
390 * WORK(J2+1, 1) stores H(J2+1, 1)
391 *
392  IF( j.LT.n ) THEN
393 *
394 * if first panel and JB=1 (NB=1), then nothing to do
395 *
396  IF( j1.GT.1 .OR. jb.GT.1 ) THEN
397 *
398 * Merge rank-1 update with BLAS-3 update
399 *
400  alpha = a( j+1, j )
401  a( j+1, j ) = one
402  CALL dcopy( n-j, a( j+1, j-1 ), 1,
403  $ work( (j+1-j1+1)+jb*n ), 1 )
404  CALL dscal( n-j, alpha, work( (j+1-j1+1)+jb*n ), 1 )
405 *
406 * K1 identifies if the previous column of the panel has been
407 * explicitly stored, e.g., K1=1 and K2= 0 for the first panel,
408 * while K1=0 and K2=1 for the rest
409 *
410  IF( j1.GT.1 ) THEN
411 *
412 * Not first panel
413 *
414  k2 = 1
415  ELSE
416 *
417 * First panel
418 *
419  k2 = 0
420 *
421 * First update skips the first column
422 *
423  jb = jb - 1
424  END IF
425 *
426  DO j2 = j+1, n, nb
427  nj = min( nb, n-j2+1 )
428 *
429 * Update (J2, J2) diagonal block with DGEMV
430 *
431  j3 = j2
432  DO mj = nj-1, 1, -1
433  CALL dgemv( 'No transpose', mj, jb+1,
434  $ -one, work( j3-j1+1+k1*n ), n,
435  $ a( j3, j1-k2 ), lda,
436  $ one, a( j3, j3 ), 1 )
437  j3 = j3 + 1
438  END DO
439 *
440 * Update off-diagonal block in J2-th block column with DGEMM
441 *
442  CALL dgemm( 'No transpose', 'Transpose',
443  $ n-j3+1, nj, jb+1,
444  $ -one, work( j3-j1+1+k1*n ), n,
445  $ a( j2, j1-k2 ), lda,
446  $ one, a( j3, j2 ), lda )
447  END DO
448 *
449 * Recover T( J+1, J )
450 *
451  a( j+1, j ) = alpha
452  END IF
453 *
454 * WORK(J+1, 1) stores H(J+1, 1)
455 *
456  CALL dcopy( n-j, a( j+1, j+1 ), 1, work( 1 ), 1 )
457  END IF
458  GO TO 11
459  END IF
460 *
461  20 CONTINUE
462  RETURN
463 *
464 * End of DSYTRF_AA
465 *
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:84
subroutine dlasyf_aa(UPLO, J1, M, NB, A, LDA, IPIV, H, LDH, WORK)
DLASYF_AA
Definition: dlasyf_aa.f:146
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: tstiee.f:83
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:84
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:81
Here is the call graph for this function:
Here is the caller graph for this function: