LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dlavsy ( character  UPLO,
character  TRANS,
character  DIAG,
integer  N,
integer  NRHS,
double precision, dimension( lda, * )  A,
integer  LDA,
integer, dimension( * )  IPIV,
double precision, dimension( ldb, * )  B,
integer  LDB,
integer  INFO 
)

DLAVSY

Purpose:
 DLAVSY  performs one of the matrix-vector operations
    x := A*x  or  x := A'*x,
 where x is an N element vector and A is one of the factors
 from the block U*D*U' or L*D*L' factorization computed by DSYTRF.

 If TRANS = 'N', multiplies by U  or U * D  (or L  or L * D)
 If TRANS = 'T', multiplies by U' or D * U' (or L' or D * L')
 If TRANS = 'C', multiplies by U' or D * U' (or L' or D * L')
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the factor stored in A is upper or lower
          triangular.
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[in]TRANS
          TRANS is CHARACTER*1
          Specifies the operation to be performed:
          = 'N':  x := A*x
          = 'T':  x := A'*x
          = 'C':  x := A'*x
[in]DIAG
          DIAG is CHARACTER*1
          Specifies whether or not the diagonal blocks are unit
          matrices.  If the diagonal blocks are assumed to be unit,
          then A = U or A = L, otherwise A = U*D or A = L*D.
          = 'U':  Diagonal blocks are assumed to be unit matrices.
          = 'N':  Diagonal blocks are assumed to be non-unit matrices.
[in]N
          N is INTEGER
          The number of rows and columns of the matrix A.  N >= 0.
[in]NRHS
          NRHS is INTEGER
          The number of right hand sides, i.e., the number of vectors
          x to be multiplied by A.  NRHS >= 0.
[in]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          The block diagonal matrix D and the multipliers used to
          obtain the factor U or L as computed by DSYTRF.
          Stored as a 2-D triangular matrix.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in]IPIV
          IPIV is INTEGER array, dimension (N)
          Details of the interchanges and the block structure of D,
          as determined by DSYTRF.

          If UPLO = 'U':
               If IPIV(k) > 0, then rows and columns k and IPIV(k)
               were interchanged and D(k,k) is a 1-by-1 diagonal block.
               (If IPIV( k ) = k, no interchange was done).

               If IPIV(k) = IPIV(k-1) < 0, then rows and
               columns k-1 and -IPIV(k) were interchanged,
               D(k-1:k,k-1:k) is a 2-by-2 diagonal block.

          If UPLO = 'L':
               If IPIV(k) > 0, then rows and columns k and IPIV(k)
               were interchanged and D(k,k) is a 1-by-1 diagonal block.
               (If IPIV( k ) = k, no interchange was done).

               If IPIV(k) = IPIV(k+1) < 0, then rows and
               columns k+1 and -IPIV(k) were interchanged,
               D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
[in,out]B
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
          On entry, B contains NRHS vectors of length N.
          On exit, B is overwritten with the product A * B.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -k, the k-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2013

Definition at line 157 of file dlavsy.f.

157 *
158 * -- LAPACK test routine (version 3.5.0) --
159 * -- LAPACK is a software package provided by Univ. of Tennessee, --
160 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
161 * November 2013
162 *
163 * .. Scalar Arguments ..
164  CHARACTER diag, trans, uplo
165  INTEGER info, lda, ldb, n, nrhs
166 * ..
167 * .. Array Arguments ..
168  INTEGER ipiv( * )
169  DOUBLE PRECISION a( lda, * ), b( ldb, * )
170 * ..
171 *
172 * =====================================================================
173 *
174 * .. Parameters ..
175  DOUBLE PRECISION one
176  parameter ( one = 1.0d+0 )
177 * ..
178 * .. Local Scalars ..
179  LOGICAL nounit
180  INTEGER j, k, kp
181  DOUBLE PRECISION d11, d12, d21, d22, t1, t2
182 * ..
183 * .. External Functions ..
184  LOGICAL lsame
185  EXTERNAL lsame
186 * ..
187 * .. External Subroutines ..
188  EXTERNAL dgemv, dger, dscal, dswap, xerbla
189 * ..
190 * .. Intrinsic Functions ..
191  INTRINSIC abs, max
192 * ..
193 * .. Executable Statements ..
194 *
195 * Test the input parameters.
196 *
197  info = 0
198  IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
199  info = -1
200  ELSE IF( .NOT.lsame( trans, 'N' ) .AND. .NOT.
201  $ lsame( trans, 'T' ) .AND. .NOT.lsame( trans, 'C' ) ) THEN
202  info = -2
203  ELSE IF( .NOT.lsame( diag, 'U' ) .AND. .NOT.lsame( diag, 'N' ) )
204  $ THEN
205  info = -3
206  ELSE IF( n.LT.0 ) THEN
207  info = -4
208  ELSE IF( lda.LT.max( 1, n ) ) THEN
209  info = -6
210  ELSE IF( ldb.LT.max( 1, n ) ) THEN
211  info = -9
212  END IF
213  IF( info.NE.0 ) THEN
214  CALL xerbla( 'DLAVSY ', -info )
215  RETURN
216  END IF
217 *
218 * Quick return if possible.
219 *
220  IF( n.EQ.0 )
221  $ RETURN
222 *
223  nounit = lsame( diag, 'N' )
224 *------------------------------------------
225 *
226 * Compute B := A * B (No transpose)
227 *
228 *------------------------------------------
229  IF( lsame( trans, 'N' ) ) THEN
230 *
231 * Compute B := U*B
232 * where U = P(m)*inv(U(m))* ... *P(1)*inv(U(1))
233 *
234  IF( lsame( uplo, 'U' ) ) THEN
235 *
236 * Loop forward applying the transformations.
237 *
238  k = 1
239  10 CONTINUE
240  IF( k.GT.n )
241  $ GO TO 30
242  IF( ipiv( k ).GT.0 ) THEN
243 *
244 * 1 x 1 pivot block
245 *
246 * Multiply by the diagonal element if forming U * D.
247 *
248  IF( nounit )
249  $ CALL dscal( nrhs, a( k, k ), b( k, 1 ), ldb )
250 *
251 * Multiply by P(K) * inv(U(K)) if K > 1.
252 *
253  IF( k.GT.1 ) THEN
254 *
255 * Apply the transformation.
256 *
257  CALL dger( k-1, nrhs, one, a( 1, k ), 1, b( k, 1 ),
258  $ ldb, b( 1, 1 ), ldb )
259 *
260 * Interchange if P(K) .ne. I.
261 *
262  kp = ipiv( k )
263  IF( kp.NE.k )
264  $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
265  END IF
266  k = k + 1
267  ELSE
268 *
269 * 2 x 2 pivot block
270 *
271 * Multiply by the diagonal block if forming U * D.
272 *
273  IF( nounit ) THEN
274  d11 = a( k, k )
275  d22 = a( k+1, k+1 )
276  d12 = a( k, k+1 )
277  d21 = d12
278  DO 20 j = 1, nrhs
279  t1 = b( k, j )
280  t2 = b( k+1, j )
281  b( k, j ) = d11*t1 + d12*t2
282  b( k+1, j ) = d21*t1 + d22*t2
283  20 CONTINUE
284  END IF
285 *
286 * Multiply by P(K) * inv(U(K)) if K > 1.
287 *
288  IF( k.GT.1 ) THEN
289 *
290 * Apply the transformations.
291 *
292  CALL dger( k-1, nrhs, one, a( 1, k ), 1, b( k, 1 ),
293  $ ldb, b( 1, 1 ), ldb )
294  CALL dger( k-1, nrhs, one, a( 1, k+1 ), 1,
295  $ b( k+1, 1 ), ldb, b( 1, 1 ), ldb )
296 *
297 * Interchange if P(K) .ne. I.
298 *
299  kp = abs( ipiv( k ) )
300  IF( kp.NE.k )
301  $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
302  END IF
303  k = k + 2
304  END IF
305  GO TO 10
306  30 CONTINUE
307 *
308 * Compute B := L*B
309 * where L = P(1)*inv(L(1))* ... *P(m)*inv(L(m)) .
310 *
311  ELSE
312 *
313 * Loop backward applying the transformations to B.
314 *
315  k = n
316  40 CONTINUE
317  IF( k.LT.1 )
318  $ GO TO 60
319 *
320 * Test the pivot index. If greater than zero, a 1 x 1
321 * pivot was used, otherwise a 2 x 2 pivot was used.
322 *
323  IF( ipiv( k ).GT.0 ) THEN
324 *
325 * 1 x 1 pivot block:
326 *
327 * Multiply by the diagonal element if forming L * D.
328 *
329  IF( nounit )
330  $ CALL dscal( nrhs, a( k, k ), b( k, 1 ), ldb )
331 *
332 * Multiply by P(K) * inv(L(K)) if K < N.
333 *
334  IF( k.NE.n ) THEN
335  kp = ipiv( k )
336 *
337 * Apply the transformation.
338 *
339  CALL dger( n-k, nrhs, one, a( k+1, k ), 1, b( k, 1 ),
340  $ ldb, b( k+1, 1 ), ldb )
341 *
342 * Interchange if a permutation was applied at the
343 * K-th step of the factorization.
344 *
345  IF( kp.NE.k )
346  $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
347  END IF
348  k = k - 1
349 *
350  ELSE
351 *
352 * 2 x 2 pivot block:
353 *
354 * Multiply by the diagonal block if forming L * D.
355 *
356  IF( nounit ) THEN
357  d11 = a( k-1, k-1 )
358  d22 = a( k, k )
359  d21 = a( k, k-1 )
360  d12 = d21
361  DO 50 j = 1, nrhs
362  t1 = b( k-1, j )
363  t2 = b( k, j )
364  b( k-1, j ) = d11*t1 + d12*t2
365  b( k, j ) = d21*t1 + d22*t2
366  50 CONTINUE
367  END IF
368 *
369 * Multiply by P(K) * inv(L(K)) if K < N.
370 *
371  IF( k.NE.n ) THEN
372 *
373 * Apply the transformation.
374 *
375  CALL dger( n-k, nrhs, one, a( k+1, k ), 1, b( k, 1 ),
376  $ ldb, b( k+1, 1 ), ldb )
377  CALL dger( n-k, nrhs, one, a( k+1, k-1 ), 1,
378  $ b( k-1, 1 ), ldb, b( k+1, 1 ), ldb )
379 *
380 * Interchange if a permutation was applied at the
381 * K-th step of the factorization.
382 *
383  kp = abs( ipiv( k ) )
384  IF( kp.NE.k )
385  $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
386  END IF
387  k = k - 2
388  END IF
389  GO TO 40
390  60 CONTINUE
391  END IF
392 *----------------------------------------
393 *
394 * Compute B := A' * B (transpose)
395 *
396 *----------------------------------------
397  ELSE
398 *
399 * Form B := U'*B
400 * where U = P(m)*inv(U(m))* ... *P(1)*inv(U(1))
401 * and U' = inv(U'(1))*P(1)* ... *inv(U'(m))*P(m)
402 *
403  IF( lsame( uplo, 'U' ) ) THEN
404 *
405 * Loop backward applying the transformations.
406 *
407  k = n
408  70 CONTINUE
409  IF( k.LT.1 )
410  $ GO TO 90
411 *
412 * 1 x 1 pivot block.
413 *
414  IF( ipiv( k ).GT.0 ) THEN
415  IF( k.GT.1 ) THEN
416 *
417 * Interchange if P(K) .ne. I.
418 *
419  kp = ipiv( k )
420  IF( kp.NE.k )
421  $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
422 *
423 * Apply the transformation
424 *
425  CALL dgemv( 'Transpose', k-1, nrhs, one, b, ldb,
426  $ a( 1, k ), 1, one, b( k, 1 ), ldb )
427  END IF
428  IF( nounit )
429  $ CALL dscal( nrhs, a( k, k ), b( k, 1 ), ldb )
430  k = k - 1
431 *
432 * 2 x 2 pivot block.
433 *
434  ELSE
435  IF( k.GT.2 ) THEN
436 *
437 * Interchange if P(K) .ne. I.
438 *
439  kp = abs( ipiv( k ) )
440  IF( kp.NE.k-1 )
441  $ CALL dswap( nrhs, b( k-1, 1 ), ldb, b( kp, 1 ),
442  $ ldb )
443 *
444 * Apply the transformations
445 *
446  CALL dgemv( 'Transpose', k-2, nrhs, one, b, ldb,
447  $ a( 1, k ), 1, one, b( k, 1 ), ldb )
448  CALL dgemv( 'Transpose', k-2, nrhs, one, b, ldb,
449  $ a( 1, k-1 ), 1, one, b( k-1, 1 ), ldb )
450  END IF
451 *
452 * Multiply by the diagonal block if non-unit.
453 *
454  IF( nounit ) THEN
455  d11 = a( k-1, k-1 )
456  d22 = a( k, k )
457  d12 = a( k-1, k )
458  d21 = d12
459  DO 80 j = 1, nrhs
460  t1 = b( k-1, j )
461  t2 = b( k, j )
462  b( k-1, j ) = d11*t1 + d12*t2
463  b( k, j ) = d21*t1 + d22*t2
464  80 CONTINUE
465  END IF
466  k = k - 2
467  END IF
468  GO TO 70
469  90 CONTINUE
470 *
471 * Form B := L'*B
472 * where L = P(1)*inv(L(1))* ... *P(m)*inv(L(m))
473 * and L' = inv(L'(m))*P(m)* ... *inv(L'(1))*P(1)
474 *
475  ELSE
476 *
477 * Loop forward applying the L-transformations.
478 *
479  k = 1
480  100 CONTINUE
481  IF( k.GT.n )
482  $ GO TO 120
483 *
484 * 1 x 1 pivot block
485 *
486  IF( ipiv( k ).GT.0 ) THEN
487  IF( k.LT.n ) THEN
488 *
489 * Interchange if P(K) .ne. I.
490 *
491  kp = ipiv( k )
492  IF( kp.NE.k )
493  $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
494 *
495 * Apply the transformation
496 *
497  CALL dgemv( 'Transpose', n-k, nrhs, one, b( k+1, 1 ),
498  $ ldb, a( k+1, k ), 1, one, b( k, 1 ), ldb )
499  END IF
500  IF( nounit )
501  $ CALL dscal( nrhs, a( k, k ), b( k, 1 ), ldb )
502  k = k + 1
503 *
504 * 2 x 2 pivot block.
505 *
506  ELSE
507  IF( k.LT.n-1 ) THEN
508 *
509 * Interchange if P(K) .ne. I.
510 *
511  kp = abs( ipiv( k ) )
512  IF( kp.NE.k+1 )
513  $ CALL dswap( nrhs, b( k+1, 1 ), ldb, b( kp, 1 ),
514  $ ldb )
515 *
516 * Apply the transformation
517 *
518  CALL dgemv( 'Transpose', n-k-1, nrhs, one,
519  $ b( k+2, 1 ), ldb, a( k+2, k+1 ), 1, one,
520  $ b( k+1, 1 ), ldb )
521  CALL dgemv( 'Transpose', n-k-1, nrhs, one,
522  $ b( k+2, 1 ), ldb, a( k+2, k ), 1, one,
523  $ b( k, 1 ), ldb )
524  END IF
525 *
526 * Multiply by the diagonal block if non-unit.
527 *
528  IF( nounit ) THEN
529  d11 = a( k, k )
530  d22 = a( k+1, k+1 )
531  d21 = a( k+1, k )
532  d12 = d21
533  DO 110 j = 1, nrhs
534  t1 = b( k, j )
535  t2 = b( k+1, j )
536  b( k, j ) = d11*t1 + d12*t2
537  b( k+1, j ) = d21*t1 + d22*t2
538  110 CONTINUE
539  END IF
540  k = k + 2
541  END IF
542  GO TO 100
543  120 CONTINUE
544  END IF
545 *
546  END IF
547  RETURN
548 *
549 * End of DLAVSY
550 *
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine dger(M, N, ALPHA, X, INCX, Y, INCY, A, LDA)
DGER
Definition: dger.f:132
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: