LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
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**H*T*U or
42*> A = L*T*L**H 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**H*T*U;
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] 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[out] WORK
102*> \verbatim
103*> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
104*> \endverbatim
105*>
106*> \param[in] LWORK
107*> \verbatim
108*> LWORK is INTEGER
109*> The dimension of the array WORK. LWORK >= max(1,3*N-2).
110*> \endverbatim
111*>
112*> \param[out] INFO
113*> \verbatim
114*> INFO is INTEGER
115*> = 0: successful exit
116*> < 0: if INFO = -i, the i-th argument had an illegal value
117*> \endverbatim
118*
119* Authors:
120* ========
121*
122*> \author Univ. of Tennessee
123*> \author Univ. of California Berkeley
124*> \author Univ. of Colorado Denver
125*> \author NAG Ltd.
126*
127*> \ingroup hetrs_aa
128*
129* =====================================================================
130 SUBROUTINE zhetrs_aa( UPLO, N, NRHS, A, LDA, IPIV, B, LDB,
131 $ WORK, LWORK, INFO )
132*
133* -- LAPACK computational routine --
134* -- LAPACK is a software package provided by Univ. of Tennessee, --
135* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
136*
137 IMPLICIT NONE
138*
139* .. Scalar Arguments ..
140 CHARACTER UPLO
141 INTEGER N, NRHS, LDA, LDB, LWORK, INFO
142* ..
143* .. Array Arguments ..
144 INTEGER IPIV( * )
145 COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * )
146* ..
147*
148* =====================================================================
149*
150 COMPLEX*16 ONE
151 parameter( one = 1.0d+0 )
152* ..
153* .. Local Scalars ..
154 LOGICAL LQUERY, UPPER
155 INTEGER K, KP, LWKOPT
156* ..
157* .. External Functions ..
158 LOGICAL LSAME
159 EXTERNAL lsame
160* ..
161* .. External Subroutines ..
162 EXTERNAL zgtsv, zswap, ztrsm, zlacgv, zlacpy, xerbla
163* ..
164* .. Intrinsic Functions ..
165 INTRINSIC max
166* ..
167* .. Executable Statements ..
168*
169 info = 0
170 upper = lsame( uplo, 'U' )
171 lquery = ( lwork.EQ.-1 )
172 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
173 info = -1
174 ELSE IF( n.LT.0 ) THEN
175 info = -2
176 ELSE IF( nrhs.LT.0 ) THEN
177 info = -3
178 ELSE IF( lda.LT.max( 1, n ) ) THEN
179 info = -5
180 ELSE IF( ldb.LT.max( 1, n ) ) THEN
181 info = -8
182 ELSE IF( lwork.LT.max( 1, 3*n-2 ) .AND. .NOT.lquery ) THEN
183 info = -10
184 END IF
185 IF( info.NE.0 ) THEN
186 CALL xerbla( 'ZHETRS_AA', -info )
187 RETURN
188 ELSE IF( lquery ) THEN
189 lwkopt = (3*n-2)
190 work( 1 ) = lwkopt
191 RETURN
192 END IF
193*
194* Quick return if possible
195*
196 IF( n.EQ.0 .OR. nrhs.EQ.0 )
197 $ RETURN
198*
199 IF( upper ) THEN
200*
201* Solve A*X = B, where A = U**H*T*U.
202*
203* 1) Forward substitution with U**H
204*
205 IF( n.GT.1 ) THEN
206*
207* Pivot, P**T * B -> B
208*
209 DO k = 1, n
210 kp = ipiv( k )
211 IF( kp.NE.k )
212 $ CALL zswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
213 END DO
214*
215* Compute U**H \ B -> B [ (U**H \P**T * B) ]
216*
217 CALL ztrsm( 'L', 'U', 'C', 'U', n-1, nrhs, one, a( 1, 2 ),
218 $ lda, b( 2, 1 ), ldb )
219 END IF
220*
221* 2) Solve with triangular matrix T
222*
223* Compute T \ B -> B [ T \ (U**H \P**T * B) ]
224*
225 CALL zlacpy( 'F', 1, n, a(1, 1), lda+1, work(n), 1 )
226 IF( n.GT.1 ) THEN
227 CALL zlacpy( 'F', 1, n-1, a( 1, 2 ), lda+1, work( 2*n ), 1)
228 CALL zlacpy( 'F', 1, n-1, a( 1, 2 ), lda+1, work( 1 ), 1 )
229 CALL zlacgv( n-1, work( 1 ), 1 )
230 END IF
231 CALL zgtsv( n, nrhs, work(1), work(n), work(2*n), b, ldb,
232 $ info )
233*
234* 3) Backward substitution with U
235*
236 IF( n.GT.1 ) THEN
237*
238* Compute U \ B -> B [ U \ (T \ (U**H \P**T * B) ) ]
239*
240 CALL ztrsm( 'L', 'U', 'N', 'U', n-1, nrhs, one, a( 1, 2 ),
241 $ lda, b(2, 1), ldb)
242*
243* Pivot, P * B [ P * (U**H \ (T \ (U \P**T * B) )) ]
244*
245 DO k = n, 1, -1
246 kp = ipiv( k )
247 IF( kp.NE.k )
248 $ CALL zswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
249 END DO
250 END IF
251*
252 ELSE
253*
254* Solve A*X = B, where A = L*T*L**H.
255*
256* 1) Forward substitution with L
257*
258 IF( n.GT.1 ) THEN
259*
260* Pivot, P**T * B -> B
261*
262 DO k = 1, n
263 kp = ipiv( k )
264 IF( kp.NE.k )
265 $ CALL zswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
266 END DO
267*
268* Compute L \ B -> B [ (L \P**T * B) ]
269*
270 CALL ztrsm( 'L', 'L', 'N', 'U', n-1, nrhs, one, a( 2, 1 ),
271 $ lda, b(2, 1), ldb)
272 END IF
273*
274* 2) Solve with triangular matrix T
275*
276* Compute T \ B -> B [ T \ (L \P**T * B) ]
277*
278 CALL zlacpy( 'F', 1, n, a(1, 1), lda+1, work(n), 1)
279 IF( n.GT.1 ) THEN
280 CALL zlacpy( 'F', 1, n-1, a( 2, 1 ), lda+1, work( 1 ), 1)
281 CALL zlacpy( 'F', 1, n-1, a( 2, 1 ), lda+1, work( 2*n ), 1)
282 CALL zlacgv( n-1, work( 2*n ), 1 )
283 END IF
284 CALL zgtsv(n, nrhs, work(1), work(n), work(2*n), b, ldb,
285 $ info)
286*
287* 3) Backward substitution with L**H
288*
289 IF( n.GT.1 ) THEN
290*
291* Compute L**H \ B -> B [ L**H \ (T \ (L \P**T * B) ) ]
292*
293 CALL ztrsm( 'L', 'L', 'C', 'U', n-1, nrhs, one, a( 2, 1 ),
294 $ lda, b( 2, 1 ), ldb)
295*
296* Pivot, P * B [ P * (L**H \ (T \ (L \P**T * B) )) ]
297*
298 DO k = n, 1, -1
299 kp = ipiv( k )
300 IF( kp.NE.k )
301 $ CALL zswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
302 END DO
303 END IF
304*
305 END IF
306*
307 RETURN
308*
309* End of ZHETRS_AA
310*
311 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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:124
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:74
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
subroutine ztrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRSM
Definition ztrsm.f:180