LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
zhetrs_aa.f
Go to the documentation of this file.
1 *> \brief \b ZHETRS_AA
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZHETRS_AA + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zhetrs_aa.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zhetrs_aa.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zhetrs_aa.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZHETRS_AA( UPLO, N, NRHS, A, LDA, IPIV, B, LDB,
22 * WORK, LWORK, INFO )
23 *
24 * .. Scalar Arguments ..
25 * CHARACTER UPLO
26 * INTEGER N, NRHS, LDA, LDB, LWORK, INFO
27 * ..
28 * .. Array Arguments ..
29 * INTEGER IPIV( * )
30 * COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * )
31 * ..
32 *
33 *
34 *
35 *> \par Purpose:
36 * =============
37 *>
38 *> \verbatim
39 *>
40 *> ZHETRS_AA solves a system of linear equations A*X = B with a complex
41 *> hermitian matrix A using the factorization A = U*T*U**H or
42 *> A = L*T*L**T computed by ZHETRF_AA.
43 *> \endverbatim
44 *
45 * Arguments:
46 * ==========
47 *
48 *> \param[in] UPLO
49 *> \verbatim
50 *> UPLO is CHARACTER*1
51 *> Specifies whether the details of the factorization are stored
52 *> as an upper or lower triangular matrix.
53 *> = 'U': Upper triangular, form is A = U*T*U**H;
54 *> = 'L': Lower triangular, form is A = L*T*L**H.
55 *> \endverbatim
56 *>
57 *> \param[in] N
58 *> \verbatim
59 *> N is INTEGER
60 *> The order of the matrix A. N >= 0.
61 *> \endverbatim
62 *>
63 *> \param[in] NRHS
64 *> \verbatim
65 *> NRHS is INTEGER
66 *> The number of right hand sides, i.e., the number of columns
67 *> of the matrix B. NRHS >= 0.
68 *> \endverbatim
69 *>
70 *> \param[in,out] A
71 *> \verbatim
72 *> A is COMPLEX*16 array, dimension (LDA,N)
73 *> Details of factors computed by ZHETRF_AA.
74 *> \endverbatim
75 *>
76 *> \param[in] LDA
77 *> \verbatim
78 *> LDA is INTEGER
79 *> The leading dimension of the array A. LDA >= max(1,N).
80 *> \endverbatim
81 *>
82 *> \param[in] IPIV
83 *> \verbatim
84 *> IPIV is INTEGER array, dimension (N)
85 *> Details of the interchanges as computed by ZHETRF_AA.
86 *> \endverbatim
87 *>
88 *> \param[in,out] B
89 *> \verbatim
90 *> B is COMPLEX*16 array, dimension (LDB,NRHS)
91 *> On entry, the right hand side matrix B.
92 *> On exit, the solution matrix X.
93 *> \endverbatim
94 *>
95 *> \param[in] LDB
96 *> \verbatim
97 *> LDB is INTEGER
98 *> The leading dimension of the array B. LDB >= max(1,N).
99 *> \endverbatim
100 *>
101 *> \param[in] WORK
102 *> \verbatim
103 *> WORK is DOUBLE array, dimension (MAX(1,LWORK))
104 *> \endverbatim
105 *>
106 *> \param[in] LWORK
107 *> \verbatim
108 *> LWORK is INTEGER, LWORK >= MAX(1,3*N-2).
109 *>
110 *> \param[out] INFO
111 *> \verbatim
112 *> INFO is INTEGER
113 *> = 0: successful exit
114 *> < 0: if INFO = -i, the i-th argument had an illegal value
115 *> \endverbatim
116 *
117 * Authors:
118 * ========
119 *
120 *> \author Univ. of Tennessee
121 *> \author Univ. of California Berkeley
122 *> \author Univ. of Colorado Denver
123 *> \author NAG Ltd.
124 *
125 *> \date December 2016
126 *
127 *> \ingroup complex16HEcomputational
128 *
129 * =====================================================================
130  SUBROUTINE zhetrs_aa( UPLO, N, NRHS, A, LDA, IPIV, B, LDB,
131  $ work, lwork, info )
132 *
133 * -- LAPACK computational routine (version 3.7.0) --
134 * -- LAPACK is a software package provided by Univ. of Tennessee, --
135 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
136 * December 2016
137 *
138  IMPLICIT NONE
139 *
140 * .. Scalar Arguments ..
141  CHARACTER UPLO
142  INTEGER N, NRHS, LDA, LDB, LWORK, INFO
143 * ..
144 * .. Array Arguments ..
145  INTEGER IPIV( * )
146  COMPLEX*16 A( lda, * ), B( ldb, * ), WORK( * )
147 * ..
148 *
149 * =====================================================================
150 *
151  COMPLEX*16 ONE
152  parameter ( one = 1.0d+0 )
153 * ..
154 * .. Local Scalars ..
155  LOGICAL LQUERY, UPPER
156  INTEGER K, KP, LWKOPT
157 * ..
158 * .. External Functions ..
159  LOGICAL LSAME
160  EXTERNAL lsame
161 * ..
162 * .. External Subroutines ..
163  EXTERNAL zgtsv, zswap, ztrsm, xerbla
164 * ..
165 * .. Intrinsic Functions ..
166  INTRINSIC max
167 * ..
168 * .. Executable Statements ..
169 *
170  info = 0
171  upper = lsame( uplo, 'U' )
172  lquery = ( lwork.EQ.-1 )
173  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
174  info = -1
175  ELSE IF( n.LT.0 ) THEN
176  info = -2
177  ELSE IF( nrhs.LT.0 ) THEN
178  info = -3
179  ELSE IF( lda.LT.max( 1, n ) ) THEN
180  info = -5
181  ELSE IF( ldb.LT.max( 1, n ) ) THEN
182  info = -8
183  ELSE IF( lwork.LT.max( 1, 3*n-2 ) .AND. .NOT.lquery ) THEN
184  info = -10
185  END IF
186  IF( info.NE.0 ) THEN
187  CALL xerbla( 'ZHETRS_AA', -info )
188  RETURN
189  ELSE IF( lquery ) THEN
190  lwkopt = (3*n-2)
191  work( 1 ) = lwkopt
192  RETURN
193  END IF
194 *
195 * Quick return if possible
196 *
197  IF( n.EQ.0 .OR. nrhs.EQ.0 )
198  $ RETURN
199 *
200  IF( upper ) THEN
201 *
202 * Solve A*X = B, where A = U*T*U**T.
203 *
204 * Pivot, P**T * B
205 *
206  DO k = 1, n
207  kp = ipiv( k )
208  IF( kp.NE.k )
209  $ CALL zswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
210  END DO
211 *
212 * Compute (U \P**T * B) -> B [ (U \P**T * B) ]
213 *
214  CALL ztrsm('L', 'U', 'C', 'U', n-1, nrhs, one, a( 1, 2 ), lda,
215  $ b( 2, 1 ), ldb)
216 *
217 * Compute T \ B -> B [ T \ (U \P**T * B) ]
218 *
219  CALL zlacpy( 'F', 1, n, a(1, 1), lda+1, work(n), 1)
220  IF( n.GT.1 ) THEN
221  CALL zlacpy( 'F', 1, n-1, a( 1, 2 ), lda+1, work( 2*n ), 1)
222  CALL zlacpy( 'F', 1, n-1, a( 1, 2 ), lda+1, work( 1 ), 1)
223  CALL zlacgv( n-1, work( 1 ), 1 )
224  END IF
225  CALL zgtsv(n, nrhs, work(1), work(n), work(2*n), b, ldb,
226  $ info)
227 *
228 * Compute (U**T \ B) -> B [ U**T \ (T \ (U \P**T * B) ) ]
229 *
230  CALL ztrsm( 'L', 'U', 'N', 'U', n-1, nrhs, one, a( 1, 2 ), lda,
231  $ b(2, 1), ldb)
232 *
233 * Pivot, P * B [ P * (U**T \ (T \ (U \P**T * B) )) ]
234 *
235  DO k = n, 1, -1
236  kp = ipiv( k )
237  IF( kp.NE.k )
238  $ CALL zswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
239  END DO
240 *
241  ELSE
242 *
243 * Solve A*X = B, where A = L*T*L**T.
244 *
245 * Pivot, P**T * B
246 *
247  DO k = 1, n
248  kp = ipiv( k )
249  IF( kp.NE.k )
250  $ CALL zswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
251  END DO
252 *
253 * Compute (L \P**T * B) -> B [ (L \P**T * B) ]
254 *
255  CALL ztrsm( 'L', 'L', 'N', 'U', n-1, nrhs, one, a( 2, 1 ), lda,
256  $ b(2, 1), ldb)
257 *
258 * Compute T \ B -> B [ T \ (L \P**T * B) ]
259 *
260  CALL zlacpy( 'F', 1, n, a(1, 1), lda+1, work(n), 1)
261  IF( n.GT.1 ) THEN
262  CALL zlacpy( 'F', 1, n-1, a( 2, 1 ), lda+1, work( 1 ), 1)
263  CALL zlacpy( 'F', 1, n-1, a( 2, 1 ), lda+1, work( 2*n ), 1)
264  CALL zlacgv( n-1, work( 2*n ), 1 )
265  END IF
266  CALL zgtsv(n, nrhs, work(1), work(n), work(2*n), b, ldb,
267  $ info)
268 *
269 * Compute (L**T \ B) -> B [ L**T \ (T \ (L \P**T * B) ) ]
270 *
271  CALL ztrsm( 'L', 'L', 'C', 'U', n-1, nrhs, one, a( 2, 1 ), lda,
272  $ b( 2, 1 ), ldb)
273 *
274 * Pivot, P * B [ P * (L**T \ (T \ (L \P**T * B) )) ]
275 *
276  DO k = n, 1, -1
277  kp = ipiv( k )
278  IF( kp.NE.k )
279  $ CALL zswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
280  END DO
281 *
282  END IF
283 *
284  RETURN
285 *
286 * End of ZHETRS_AA
287 *
288  END
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zhetrs_aa(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, WORK, LWORK, INFO)
ZHETRS_AA
Definition: zhetrs_aa.f:132
subroutine zlacgv(N, X, INCX)
ZLACGV conjugates a complex vector.
Definition: zlacgv.f:76
subroutine ztrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRSM
Definition: ztrsm.f:182
subroutine zgtsv(N, NRHS, DL, D, DU, B, LDB, INFO)
ZGTSV computes the solution to system of linear equations A * X = B for GT matrices ...
Definition: zgtsv.f:126