LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dsytrs_aa.f
Go to the documentation of this file.
1*> \brief \b DSYTRS_AA
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download DSYTRS_AA + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsytrs_aa.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsytrs_aa.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsytrs_aa.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DSYTRS_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* DOUBLE PRECISION A( LDA, * ), B( LDB, * ), WORK( * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> DSYTRS_AA solves a system of linear equations A*X = B with a real
40*> symmetric matrix A using the factorization A = U**T*T*U or
41*> A = L*T*L**T computed by DSYTRF_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 DOUBLE PRECISION array, dimension (LDA,N)
72*> Details of factors computed by DSYTRF_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 DSYTRF_AA.
85*> \endverbatim
86*>
87*> \param[in,out] B
88*> \verbatim
89*> B is DOUBLE PRECISION 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 DOUBLE PRECISION 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 dsytrs_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 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), WORK( * )
145* ..
146*
147* =====================================================================
148*
149 DOUBLE PRECISION ONE
150 parameter( one = 1.0d+0 )
151* ..
152* .. Local Scalars ..
153 LOGICAL LQUERY, UPPER
154 INTEGER K, KP, LWKOPT
155* ..
156* .. External Functions ..
157 LOGICAL LSAME
158 EXTERNAL lsame
159* ..
160* .. External Subroutines ..
161 EXTERNAL dlacpy, dgtsv, dswap, dtrsm, xerbla
162* ..
163* .. Intrinsic Functions ..
164 INTRINSIC max
165* ..
166* .. Executable Statements ..
167*
168 info = 0
169 upper = lsame( uplo, 'U' )
170 lquery = ( lwork.EQ.-1 )
171 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
172 info = -1
173 ELSE IF( n.LT.0 ) THEN
174 info = -2
175 ELSE IF( nrhs.LT.0 ) THEN
176 info = -3
177 ELSE IF( lda.LT.max( 1, n ) ) THEN
178 info = -5
179 ELSE IF( ldb.LT.max( 1, n ) ) THEN
180 info = -8
181 ELSE IF( lwork.LT.max( 1, 3*n-2 ) .AND. .NOT.lquery ) THEN
182 info = -10
183 END IF
184 IF( info.NE.0 ) THEN
185 CALL xerbla( 'DSYTRS_AA', -info )
186 RETURN
187 ELSE IF( lquery ) THEN
188 lwkopt = (3*n-2)
189 work( 1 ) = lwkopt
190 RETURN
191 END IF
192*
193* Quick return if possible
194*
195 IF( n.EQ.0 .OR. nrhs.EQ.0 )
196 $ RETURN
197*
198 IF( upper ) THEN
199*
200* Solve A*X = B, where A = U**T*T*U.
201*
202* 1) Forward substitution with U**T
203*
204 IF( n.GT.1 ) THEN
205*
206* Pivot, P**T * B -> B
207*
208 DO k = 1, n
209 kp = ipiv( k )
210 IF( kp.NE.k )
211 $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
212 END DO
213*
214* Compute U**T \ B -> B [ (U**T \P**T * B) ]
215*
216 CALL dtrsm('L', 'U', 'T', 'U', n-1, nrhs, one, a( 1, 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 dlacpy( 'F', 1, n, a( 1, 1 ), lda+1, work( n ), 1)
225 IF( n.GT.1 ) THEN
226 CALL dlacpy( 'F', 1, n-1, a( 1, 2 ), lda+1, work( 1 ), 1 )
227 CALL dlacpy( 'F', 1, n-1, a( 1, 2 ), lda+1, work( 2*n ), 1 )
228 END IF
229 CALL dgtsv( n, nrhs, work( 1 ), work( n ), work( 2*n ), b, ldb,
230 $ info )
231*
232* 3) Backward substitution with U
233*
234 IF( n.GT.1 ) THEN
235*
236* Compute U \ B -> B [ U \ (T \ (U**T \P**T * B) ) ]
237*
238 CALL dtrsm( 'L', 'U', 'N', 'U', n-1, nrhs, one, a( 1, 2 ),
239 $ lda, b( 2, 1 ), ldb)
240*
241* Pivot, P * B -> B [ P * (U \ (T \ (U**T \P**T * B) )) ]
242*
243 DO k = n, 1, -1
244 kp = ipiv( k )
245 IF( kp.NE.k )
246 $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
247 END DO
248 END IF
249*
250 ELSE
251*
252* Solve A*X = B, where A = L*T*L**T.
253*
254* 1) Forward substitution with L
255*
256 IF( n.GT.1 ) THEN
257*
258* Pivot, P**T * B -> B
259*
260 DO k = 1, n
261 kp = ipiv( k )
262 IF( kp.NE.k )
263 $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
264 END DO
265*
266* Compute L \ B -> B [ (L \P**T * B) ]
267*
268 CALL dtrsm( 'L', 'L', 'N', 'U', n-1, nrhs, one, a( 2, 1 ),
269 $ lda, b( 2, 1 ), ldb)
270 END IF
271*
272* 2) Solve with triangular matrix T
273*
274* Compute T \ B -> B [ T \ (L \P**T * B) ]
275*
276 CALL dlacpy( 'F', 1, n, a(1, 1), lda+1, work(n), 1)
277 IF( n.GT.1 ) THEN
278 CALL dlacpy( 'F', 1, n-1, a( 2, 1 ), lda+1, work( 1 ), 1 )
279 CALL dlacpy( 'F', 1, n-1, a( 2, 1 ), lda+1, work( 2*n ), 1 )
280 END IF
281 CALL dgtsv( n, nrhs, work( 1 ), work(n), work( 2*n ), b, ldb,
282 $ info)
283*
284* 3) Backward substitution with L**T
285*
286 IF( n.GT.1 ) THEN
287*
288* Compute (L**T \ B) -> B [ L**T \ (T \ (L \P**T * B) ) ]
289*
290 CALL dtrsm( 'L', 'L', 'T', 'U', n-1, nrhs, one, a( 2, 1 ),
291 $ lda, b( 2, 1 ), ldb)
292*
293* Pivot, P * B -> B [ P * (L**T \ (T \ (L \P**T * B) )) ]
294*
295 DO k = n, 1, -1
296 kp = ipiv( k )
297 IF( kp.NE.k )
298 $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
299 END DO
300 END IF
301*
302 END IF
303*
304 RETURN
305*
306* End of DSYTRS_AA
307*
308 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgtsv(n, nrhs, dl, d, du, b, ldb, info)
DGTSV computes the solution to system of linear equations A * X = B for GT matrices
Definition dgtsv.f:127
subroutine dsytrs_aa(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
DSYTRS_AA
Definition dsytrs_aa.f:131
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181