LAPACK 3.12.1
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*> Download DSYTRS_AA + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsytrs_aa.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsytrs_aa.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsytrs_aa.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DSYTRS_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* DOUBLE PRECISION A( LDA, * ), B( LDB, * ), WORK( * )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> DSYTRS_AA solves a system of linear equations A*X = B with a real
38*> symmetric matrix A using the factorization A = U**T*T*U or
39*> A = L*T*L**T computed by DSYTRF_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 DOUBLE PRECISION array, dimension (LDA,N)
70*> Details of factors computed by DSYTRF_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 DSYTRF_AA.
83*> \endverbatim
84*>
85*> \param[in,out] B
86*> \verbatim
87*> B is DOUBLE PRECISION 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 DOUBLE PRECISION 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.
107*> If MIN(N,NRHS) = 0, LWORK >= 1, else LWORK >= 3*N-2.
108*>
109*> If LWORK = -1, then a workspace query is assumed; the routine
110*> only calculates the minimal size of the WORK array, returns
111*> this value as the first entry of the WORK array, and no error
112*> message related to LWORK is issued by XERBLA.
113*> \endverbatim
114*>
115*> \param[out] INFO
116*> \verbatim
117*> INFO is INTEGER
118*> = 0: successful exit
119*> < 0: if INFO = -i, the i-th argument had an illegal value
120*> \endverbatim
121*
122* Authors:
123* ========
124*
125*> \author Univ. of Tennessee
126*> \author Univ. of California Berkeley
127*> \author Univ. of Colorado Denver
128*> \author NAG Ltd.
129*
130*> \ingroup hetrs_aa
131*
132* =====================================================================
133 SUBROUTINE dsytrs_aa( UPLO, N, NRHS, A, LDA, IPIV, B, LDB,
134 $ WORK, LWORK, INFO )
135*
136* -- LAPACK computational routine --
137* -- LAPACK is a software package provided by Univ. of Tennessee, --
138* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
139*
140 IMPLICIT NONE
141*
142* .. Scalar Arguments ..
143 CHARACTER UPLO
144 INTEGER N, NRHS, LDA, LDB, LWORK, INFO
145* ..
146* .. Array Arguments ..
147 INTEGER IPIV( * )
148 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), WORK( * )
149* ..
150*
151* =====================================================================
152*
153 DOUBLE PRECISION ONE
154 parameter( one = 1.0d+0 )
155* ..
156* .. Local Scalars ..
157 LOGICAL LQUERY, UPPER
158 INTEGER K, KP, LWKMIN
159* ..
160* .. External Functions ..
161 LOGICAL LSAME
162 EXTERNAL lsame
163* ..
164* .. External Subroutines ..
165 EXTERNAL dlacpy, dgtsv, dswap, dtrsm, xerbla
166* ..
167* .. Intrinsic Functions ..
168 INTRINSIC min, max
169* ..
170* .. Executable Statements ..
171*
172 info = 0
173 upper = lsame( uplo, 'U' )
174 lquery = ( lwork.EQ.-1 )
175 IF( min( n, nrhs ).EQ.0 ) THEN
176 lwkmin = 1
177 ELSE
178 lwkmin = 3*n-2
179 END IF
180*
181 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
182 info = -1
183 ELSE IF( n.LT.0 ) THEN
184 info = -2
185 ELSE IF( nrhs.LT.0 ) THEN
186 info = -3
187 ELSE IF( lda.LT.max( 1, n ) ) THEN
188 info = -5
189 ELSE IF( ldb.LT.max( 1, n ) ) THEN
190 info = -8
191 ELSE IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
192 info = -10
193 END IF
194 IF( info.NE.0 ) THEN
195 CALL xerbla( 'DSYTRS_AA', -info )
196 RETURN
197 ELSE IF( lquery ) THEN
198 work( 1 ) = lwkmin
199 RETURN
200 END IF
201*
202* Quick return if possible
203*
204 IF( min( n, nrhs ).EQ.0 )
205 $ RETURN
206*
207 IF( upper ) THEN
208*
209* Solve A*X = B, where A = U**T*T*U.
210*
211* 1) Forward substitution with U**T
212*
213 IF( n.GT.1 ) THEN
214*
215* Pivot, P**T * B -> B
216*
217 DO k = 1, n
218 kp = ipiv( k )
219 IF( kp.NE.k )
220 $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ),
221 $ ldb )
222 END DO
223*
224* Compute U**T \ B -> B [ (U**T \P**T * B) ]
225*
226 CALL dtrsm('L', 'U', 'T', 'U', n-1, nrhs, one, a( 1, 2 ),
227 $ lda, b( 2, 1 ), ldb)
228 END IF
229*
230* 2) Solve with triangular matrix T
231*
232* Compute T \ B -> B [ T \ (U**T \P**T * B) ]
233*
234 CALL dlacpy( 'F', 1, n, a( 1, 1 ), lda+1, work( n ), 1)
235 IF( n.GT.1 ) THEN
236 CALL dlacpy( 'F', 1, n-1, a( 1, 2 ), lda+1, work( 1 ),
237 $ 1 )
238 CALL dlacpy( 'F', 1, n-1, a( 1, 2 ), lda+1, work( 2*n ),
239 $ 1 )
240 END IF
241 CALL dgtsv( n, nrhs, work( 1 ), work( n ), work( 2*n ), b,
242 $ ldb,
243 $ info )
244*
245* 3) Backward substitution with U
246*
247 IF( n.GT.1 ) THEN
248*
249* Compute U \ B -> B [ U \ (T \ (U**T \P**T * B) ) ]
250*
251 CALL dtrsm( 'L', 'U', 'N', 'U', n-1, nrhs, one, a( 1,
252 $ 2 ),
253 $ lda, b( 2, 1 ), ldb)
254*
255* Pivot, P * B -> B [ P * (U \ (T \ (U**T \P**T * B) )) ]
256*
257 DO k = n, 1, -1
258 kp = ipiv( k )
259 IF( kp.NE.k )
260 $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
261 END DO
262 END IF
263*
264 ELSE
265*
266* Solve A*X = B, where A = L*T*L**T.
267*
268* 1) Forward substitution with L
269*
270 IF( n.GT.1 ) THEN
271*
272* Pivot, P**T * B -> B
273*
274 DO k = 1, n
275 kp = ipiv( k )
276 IF( kp.NE.k )
277 $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
278 END DO
279*
280* Compute L \ B -> B [ (L \P**T * B) ]
281*
282 CALL dtrsm( 'L', 'L', 'N', 'U', n-1, nrhs, one, a( 2,
283 $ 1 ),
284 $ lda, b( 2, 1 ), ldb)
285 END IF
286*
287* 2) Solve with triangular matrix T
288*
289* Compute T \ B -> B [ T \ (L \P**T * B) ]
290*
291 CALL dlacpy( 'F', 1, n, a(1, 1), lda+1, work(n), 1)
292 IF( n.GT.1 ) THEN
293 CALL dlacpy( 'F', 1, n-1, a( 2, 1 ), lda+1, work( 1 ),
294 $ 1 )
295 CALL dlacpy( 'F', 1, n-1, a( 2, 1 ), lda+1, work( 2*n ),
296 $ 1 )
297 END IF
298 CALL dgtsv( n, nrhs, work( 1 ), work(n), work( 2*n ), b,
299 $ ldb,
300 $ info)
301*
302* 3) Backward substitution with L**T
303*
304 IF( n.GT.1 ) THEN
305*
306* Compute (L**T \ B) -> B [ L**T \ (T \ (L \P**T * B) ) ]
307*
308 CALL dtrsm( 'L', 'L', 'T', 'U', n-1, nrhs, one, a( 2,
309 $ 1 ),
310 $ lda, b( 2, 1 ), ldb)
311*
312* Pivot, P * B -> B [ P * (L**T \ (T \ (L \P**T * B) )) ]
313*
314 DO k = n, 1, -1
315 kp = ipiv( k )
316 IF( kp.NE.k )
317 $ CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
318 END DO
319 END IF
320*
321 END IF
322*
323 RETURN
324*
325* End of DSYTRS_AA
326*
327 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:125
subroutine dsytrs_aa(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
DSYTRS_AA
Definition dsytrs_aa.f:135
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
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