LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
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*> \htmlonly
9*> Download CSYTRS_AA + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/csytrs_aa.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/csytrs_aa.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/csytrs_aa.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE CSYTRS_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 A( LDA, * ), B( LDB, * ), WORK( * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> CSYTRS_AA solves a system of linear equations A*X = B with a complex
40*> symmetric matrix A using the factorization A = U**T*T*U or
41*> A = L*T*L**T computed by CSYTRF_AA.
42*> \endverbatim
43*
44* Arguments:
45* ==========
46*
47*> \param[in] UPLO
48*> \verbatim
49*> UPLO is CHARACTER*1
50*> Specifies whether the details of the factorization are stored
51*> as an upper or lower triangular matrix.
52*> = 'U': Upper triangular, form is A = U**T*T*U;
53*> = 'L': Lower triangular, form is A = L*T*L**T.
54*> \endverbatim
55*>
56*> \param[in] N
57*> \verbatim
58*> N is INTEGER
59*> The order of the matrix A. N >= 0.
60*> \endverbatim
61*>
62*> \param[in] NRHS
63*> \verbatim
64*> NRHS is INTEGER
65*> The number of right hand sides, i.e., the number of columns
66*> of the matrix B. NRHS >= 0.
67*> \endverbatim
68*>
69*> \param[in] A
70*> \verbatim
71*> A is COMPLEX array, dimension (LDA,N)
72*> Details of factors computed by CSYTRF_AA.
73*> \endverbatim
74*>
75*> \param[in] LDA
76*> \verbatim
77*> LDA is INTEGER
78*> The leading dimension of the array A. LDA >= max(1,N).
79*> \endverbatim
80*>
81*> \param[in] IPIV
82*> \verbatim
83*> IPIV is INTEGER array, dimension (N)
84*> Details of the interchanges as computed by CSYTRF_AA.
85*> \endverbatim
86*>
87*> \param[in,out] B
88*> \verbatim
89*> B is COMPLEX array, dimension (LDB,NRHS)
90*> On entry, the right hand side matrix B.
91*> On exit, the solution matrix X.
92*> \endverbatim
93*>
94*> \param[in] LDB
95*> \verbatim
96*> LDB is INTEGER
97*> The leading dimension of the array B. LDB >= max(1,N).
98*> \endverbatim
99*>
100*> \param[out] WORK
101*> \verbatim
102*> WORK is COMPLEX array, dimension (MAX(1,LWORK))
103*> \endverbatim
104*>
105*> \param[in] LWORK
106*> \verbatim
107*> LWORK is INTEGER
108*> The dimension of the array WORK. LWORK >= max(1,3*N-2).
109*> \endverbatim
110*>
111*> \param[out] INFO
112*> \verbatim
113*> INFO is INTEGER
114*> = 0: successful exit
115*> < 0: if INFO = -i, the i-th argument had an illegal value
116*> \endverbatim
117*
118* Authors:
119* ========
120*
121*> \author Univ. of Tennessee
122*> \author Univ. of California Berkeley
123*> \author Univ. of Colorado Denver
124*> \author NAG Ltd.
125*
126*> \ingroup hetrs_aa
127*
128* =====================================================================
129 SUBROUTINE csytrs_aa( UPLO, N, NRHS, A, LDA, IPIV, B, LDB,
130 $ WORK, LWORK, INFO )
131*
132* -- LAPACK computational routine --
133* -- LAPACK is a software package provided by Univ. of Tennessee, --
134* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
135*
136 IMPLICIT NONE
137*
138* .. Scalar Arguments ..
139 CHARACTER UPLO
140 INTEGER N, NRHS, LDA, LDB, LWORK, INFO
141* ..
142* .. Array Arguments ..
143 INTEGER IPIV( * )
144 COMPLEX A( LDA, * ), B( LDB, * ), WORK( * )
145* ..
146*
147* =====================================================================
148*
149 COMPLEX ONE
150 parameter( one = 1.0e+0 )
151* ..
152* .. Local Scalars ..
153 LOGICAL LQUERY, UPPER
154 INTEGER K, KP, LWKOPT
155* ..
156* .. External Functions ..
157 LOGICAL LSAME
158 REAL SROUNDUP_LWORK
159 EXTERNAL lsame, sroundup_lwork
160* ..
161* .. External Subroutines ..
162 EXTERNAL clacpy, cgtsv, cswap, ctrsm, 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( 'CSYTRS_AA', -info )
187 RETURN
188 ELSE IF( lquery ) THEN
189 lwkopt = (3*n-2)
190 work( 1 ) = sroundup_lwork(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**T*T*U.
202*
203* 1) Forward substitution with U**T
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 cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
213 END DO
214*
215* Compute U**T \ B -> B [ (U**T \P**T * B) ]
216*
217 CALL ctrsm( 'L', 'U', 'T', '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**T \P**T * B) ]
224*
225 CALL clacpy( 'F', 1, n, a( 1, 1 ), lda+1, work( n ), 1)
226 IF( n.GT.1 ) THEN
227 CALL clacpy( 'F', 1, n-1, a( 1, 2 ), lda+1, work( 1 ), 1 )
228 CALL clacpy( 'F', 1, n-1, a( 1, 2 ), lda+1, work( 2*n ), 1 )
229 END IF
230 CALL cgtsv( n, nrhs, work( 1 ), work( n ), work( 2*n ), b, ldb,
231 $ info )
232*
233* 3) Backward substitution with U
234*
235 IF( n.GT.1 ) THEN
236*
237* Compute U \ B -> B [ U \ (T \ (U**T \P**T * B) ) ]
238*
239 CALL ctrsm( 'L', 'U', 'N', 'U', n-1, nrhs, one, a( 1, 2 ),
240 $ lda, b( 2, 1 ), ldb)
241*
242* Pivot, P * B -> B [ P * (U**T \ (T \ (U \P**T * B) )) ]
243*
244 DO k = n, 1, -1
245 kp = ipiv( k )
246 IF( kp.NE.k )
247 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
248 END DO
249 END IF
250*
251 ELSE
252*
253* Solve A*X = B, where A = L*T*L**T.
254*
255* 1) Forward substitution with L
256*
257 IF( n.GT.1 ) THEN
258*
259* Pivot, P**T * B -> B
260*
261 DO k = 1, n
262 kp = ipiv( k )
263 IF( kp.NE.k )
264 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
265 END DO
266*
267* Compute L \ B -> B [ (L \P**T * B) ]
268*
269 CALL ctrsm( 'L', 'L', 'N', 'U', n-1, nrhs, one, a( 2, 1 ),
270 $ lda, b( 2, 1 ), ldb)
271 END IF
272*
273* 2) Solve with triangular matrix T
274*
275*
276* Compute T \ B -> B [ T \ (L \P**T * B) ]
277*
278 CALL clacpy( 'F', 1, n, a(1, 1), lda+1, work(n), 1)
279 IF( n.GT.1 ) THEN
280 CALL clacpy( 'F', 1, n-1, a( 2, 1 ), lda+1, work( 1 ), 1 )
281 CALL clacpy( 'F', 1, n-1, a( 2, 1 ), lda+1, work( 2*n ), 1 )
282 END IF
283 CALL cgtsv( n, nrhs, work( 1 ), work(n), work( 2*n ), b, ldb,
284 $ info)
285*
286* 3) Backward substitution with L**T
287*
288 IF( n.GT.1 ) THEN
289*
290* Compute (L**T \ B) -> B [ L**T \ (T \ (L \P**T * B) ) ]
291*
292 CALL ctrsm( 'L', 'L', 'T', 'U', n-1, nrhs, one, a( 2, 1 ),
293 $ lda, b( 2, 1 ), ldb)
294*
295* Pivot, P * B -> B [ P * (L**T \ (T \ (L \P**T * B) )) ]
296*
297 DO k = n, 1, -1
298 kp = ipiv( k )
299 IF( kp.NE.k )
300 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
301 END DO
302 END IF
303*
304 END IF
305*
306 RETURN
307*
308* End of CSYTRS_AA
309*
310 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:124
subroutine csytrs_aa(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
CSYTRS_AA
Definition csytrs_aa.f:131
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
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