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