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

CLAVHE

Purpose:
 CLAVHE performs one of the matrix-vector operations
    x := A*x  or  x := A^H*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 CHETRF.

 If TRANS = 'N', multiplies by U  or U * D  (or L  or L * D)
 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
          = 'C':   x := A^H*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 COMPLEX array, dimension (LDA,N)
          The block diagonal matrix D and the multipliers used to
          obtain the factor U or L as computed by CHETRF_ROOK.
          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 CHETRF.

          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 COMPLEX 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 155 of file clavhe.f.

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