LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
chetrs_3.f
Go to the documentation of this file.
1*> \brief \b CHETRS_3
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download CHETRS_3 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/chetrs_3.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/chetrs_3.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/chetrs_3.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE CHETRS_3( UPLO, N, NRHS, A, LDA, E, IPIV, B, LDB,
20* INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER UPLO
24* INTEGER INFO, LDA, LDB, N, NRHS
25* ..
26* .. Array Arguments ..
27* INTEGER IPIV( * )
28* COMPLEX A( LDA, * ), B( LDB, * ), E( * )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*> CHETRS_3 solves a system of linear equations A * X = B with a complex
37*> Hermitian matrix A using the factorization computed
38*> by CHETRF_RK or CHETRF_BK:
39*>
40*> A = P*U*D*(U**H)*(P**T) or A = P*L*D*(L**H)*(P**T),
41*>
42*> where U (or L) is unit upper (or lower) triangular matrix,
43*> U**H (or L**H) is the conjugate of U (or L), P is a permutation
44*> matrix, P**T is the transpose of P, and D is Hermitian and block
45*> diagonal with 1-by-1 and 2-by-2 diagonal blocks.
46*>
47*> This algorithm is using Level 3 BLAS.
48*> \endverbatim
49*
50* Arguments:
51* ==========
52*
53*> \param[in] UPLO
54*> \verbatim
55*> UPLO is CHARACTER*1
56*> Specifies whether the details of the factorization are
57*> stored as an upper or lower triangular matrix:
58*> = 'U': Upper triangular, form is A = P*U*D*(U**H)*(P**T);
59*> = 'L': Lower triangular, form is A = P*L*D*(L**H)*(P**T).
60*> \endverbatim
61*>
62*> \param[in] N
63*> \verbatim
64*> N is INTEGER
65*> The order of the matrix A. N >= 0.
66*> \endverbatim
67*>
68*> \param[in] NRHS
69*> \verbatim
70*> NRHS is INTEGER
71*> The number of right hand sides, i.e., the number of columns
72*> of the matrix B. NRHS >= 0.
73*> \endverbatim
74*>
75*> \param[in] A
76*> \verbatim
77*> A is COMPLEX array, dimension (LDA,N)
78*> Diagonal of the block diagonal matrix D and factors U or L
79*> as computed by CHETRF_RK and CHETRF_BK:
80*> a) ONLY diagonal elements of the Hermitian block diagonal
81*> matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
82*> (superdiagonal (or subdiagonal) elements of D
83*> should be provided on entry in array E), and
84*> b) If UPLO = 'U': factor U in the superdiagonal part of A.
85*> If UPLO = 'L': factor L in the subdiagonal part of A.
86*> \endverbatim
87*>
88*> \param[in] LDA
89*> \verbatim
90*> LDA is INTEGER
91*> The leading dimension of the array A. LDA >= max(1,N).
92*> \endverbatim
93*>
94*> \param[in] E
95*> \verbatim
96*> E is COMPLEX array, dimension (N)
97*> On entry, contains the superdiagonal (or subdiagonal)
98*> elements of the Hermitian block diagonal matrix D
99*> with 1-by-1 or 2-by-2 diagonal blocks, where
100*> If UPLO = 'U': E(i) = D(i-1,i),i=2:N, E(1) not referenced;
101*> If UPLO = 'L': E(i) = D(i+1,i),i=1:N-1, E(N) not referenced.
102*>
103*> NOTE: For 1-by-1 diagonal block D(k), where
104*> 1 <= k <= N, the element E(k) is not referenced in both
105*> UPLO = 'U' or UPLO = 'L' cases.
106*> \endverbatim
107*>
108*> \param[in] IPIV
109*> \verbatim
110*> IPIV is INTEGER array, dimension (N)
111*> Details of the interchanges and the block structure of D
112*> as determined by CHETRF_RK or CHETRF_BK.
113*> \endverbatim
114*>
115*> \param[in,out] B
116*> \verbatim
117*> B is COMPLEX array, dimension (LDB,NRHS)
118*> On entry, the right hand side matrix B.
119*> On exit, the solution matrix X.
120*> \endverbatim
121*>
122*> \param[in] LDB
123*> \verbatim
124*> LDB is INTEGER
125*> The leading dimension of the array B. LDB >= max(1,N).
126*> \endverbatim
127*>
128*> \param[out] INFO
129*> \verbatim
130*> INFO is INTEGER
131*> = 0: successful exit
132*> < 0: if INFO = -i, the i-th argument had an illegal value
133*> \endverbatim
134*
135* Authors:
136* ========
137*
138*> \author Univ. of Tennessee
139*> \author Univ. of California Berkeley
140*> \author Univ. of Colorado Denver
141*> \author NAG Ltd.
142*
143*> \ingroup hetrs_3
144*
145*> \par Contributors:
146* ==================
147*>
148*> \verbatim
149*>
150*> June 2017, Igor Kozachenko,
151*> Computer Science Division,
152*> University of California, Berkeley
153*>
154*> September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
155*> School of Mathematics,
156*> University of Manchester
157*>
158*> \endverbatim
159*
160* =====================================================================
161 SUBROUTINE chetrs_3( UPLO, N, NRHS, A, LDA, E, IPIV, B, LDB,
162 $ INFO )
163*
164* -- LAPACK computational routine --
165* -- LAPACK is a software package provided by Univ. of Tennessee, --
166* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
167*
168* .. Scalar Arguments ..
169 CHARACTER UPLO
170 INTEGER INFO, LDA, LDB, N, NRHS
171* ..
172* .. Array Arguments ..
173 INTEGER IPIV( * )
174 COMPLEX A( LDA, * ), B( LDB, * ), E( * )
175* ..
176*
177* =====================================================================
178*
179* .. Parameters ..
180 COMPLEX ONE
181 parameter( one = ( 1.0e+0,0.0e+0 ) )
182* ..
183* .. Local Scalars ..
184 LOGICAL UPPER
185 INTEGER I, J, K, KP
186 REAL S
187 COMPLEX AK, AKM1, AKM1K, BK, BKM1, DENOM
188* ..
189* .. External Functions ..
190 LOGICAL LSAME
191 EXTERNAL lsame
192* ..
193* .. External Subroutines ..
194 EXTERNAL csscal, cswap, ctrsm, xerbla
195* ..
196* .. Intrinsic Functions ..
197 INTRINSIC abs, conjg, max, real
198* ..
199* .. Executable Statements ..
200*
201 info = 0
202 upper = lsame( uplo, 'U' )
203 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
204 info = -1
205 ELSE IF( n.LT.0 ) THEN
206 info = -2
207 ELSE IF( nrhs.LT.0 ) THEN
208 info = -3
209 ELSE IF( lda.LT.max( 1, n ) ) THEN
210 info = -5
211 ELSE IF( ldb.LT.max( 1, n ) ) THEN
212 info = -9
213 END IF
214 IF( info.NE.0 ) THEN
215 CALL xerbla( 'CHETRS_3', -info )
216 RETURN
217 END IF
218*
219* Quick return if possible
220*
221 IF( n.EQ.0 .OR. nrhs.EQ.0 )
222 $ RETURN
223*
224 IF( upper ) THEN
225*
226* Begin Upper
227*
228* Solve A*X = B, where A = U*D*U**H.
229*
230* P**T * B
231*
232* Interchange rows K and IPIV(K) of matrix B in the same order
233* that the formation order of IPIV(I) vector for Upper case.
234*
235* (We can do the simple loop over IPIV with decrement -1,
236* since the ABS value of IPIV(I) represents the row index
237* of the interchange with row i in both 1x1 and 2x2 pivot cases)
238*
239 DO k = n, 1, -1
240 kp = abs( ipiv( k ) )
241 IF( kp.NE.k ) THEN
242 CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
243 END IF
244 END DO
245*
246* Compute (U \P**T * B) -> B [ (U \P**T * B) ]
247*
248 CALL ctrsm( 'L', 'U', 'N', 'U', n, nrhs, one, a, lda, b,
249 $ ldb )
250*
251* Compute D \ B -> B [ D \ (U \P**T * B) ]
252*
253 i = n
254 DO WHILE ( i.GE.1 )
255 IF( ipiv( i ).GT.0 ) THEN
256 s = real( one ) / real( a( i, i ) )
257 CALL csscal( nrhs, s, b( i, 1 ), ldb )
258 ELSE IF ( i.GT.1 ) THEN
259 akm1k = e( i )
260 akm1 = a( i-1, i-1 ) / akm1k
261 ak = a( i, i ) / conjg( akm1k )
262 denom = akm1*ak - one
263 DO j = 1, nrhs
264 bkm1 = b( i-1, j ) / akm1k
265 bk = b( i, j ) / conjg( akm1k )
266 b( i-1, j ) = ( ak*bkm1-bk ) / denom
267 b( i, j ) = ( akm1*bk-bkm1 ) / denom
268 END DO
269 i = i - 1
270 END IF
271 i = i - 1
272 END DO
273*
274* Compute (U**H \ B) -> B [ U**H \ (D \ (U \P**T * B) ) ]
275*
276 CALL ctrsm( 'L', 'U', 'C', 'U', n, nrhs, one, a, lda, b,
277 $ ldb )
278*
279* P * B [ P * (U**H \ (D \ (U \P**T * B) )) ]
280*
281* Interchange rows K and IPIV(K) of matrix B in reverse order
282* from the formation order of IPIV(I) vector for Upper case.
283*
284* (We can do the simple loop over IPIV with increment 1,
285* since the ABS value of IPIV(I) represents the row index
286* of the interchange with row i in both 1x1 and 2x2 pivot cases)
287*
288 DO k = 1, n, 1
289 kp = abs( ipiv( k ) )
290 IF( kp.NE.k ) THEN
291 CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
292 END IF
293 END DO
294*
295 ELSE
296*
297* Begin Lower
298*
299* Solve A*X = B, where A = L*D*L**H.
300*
301* P**T * B
302* Interchange rows K and IPIV(K) of matrix B in the same order
303* that the formation order of IPIV(I) vector for Lower case.
304*
305* (We can do the simple loop over IPIV with increment 1,
306* since the ABS value of IPIV(I) represents the row index
307* of the interchange with row i in both 1x1 and 2x2 pivot cases)
308*
309 DO k = 1, n, 1
310 kp = abs( ipiv( k ) )
311 IF( kp.NE.k ) THEN
312 CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
313 END IF
314 END DO
315*
316* Compute (L \P**T * B) -> B [ (L \P**T * B) ]
317*
318 CALL ctrsm( 'L', 'L', 'N', 'U', n, nrhs, one, a, lda, b,
319 $ ldb )
320*
321* Compute D \ B -> B [ D \ (L \P**T * B) ]
322*
323 i = 1
324 DO WHILE ( i.LE.n )
325 IF( ipiv( i ).GT.0 ) THEN
326 s = real( one ) / real( a( i, i ) )
327 CALL csscal( nrhs, s, b( i, 1 ), ldb )
328 ELSE IF( i.LT.n ) THEN
329 akm1k = e( i )
330 akm1 = a( i, i ) / conjg( akm1k )
331 ak = a( i+1, i+1 ) / akm1k
332 denom = akm1*ak - one
333 DO j = 1, nrhs
334 bkm1 = b( i, j ) / conjg( akm1k )
335 bk = b( i+1, j ) / akm1k
336 b( i, j ) = ( ak*bkm1-bk ) / denom
337 b( i+1, j ) = ( akm1*bk-bkm1 ) / denom
338 END DO
339 i = i + 1
340 END IF
341 i = i + 1
342 END DO
343*
344* Compute (L**H \ B) -> B [ L**H \ (D \ (L \P**T * B) ) ]
345*
346 CALL ctrsm('L', 'L', 'C', 'U', n, nrhs, one, a, lda, b,
347 $ ldb )
348*
349* P * B [ P * (L**H \ (D \ (L \P**T * B) )) ]
350*
351* Interchange rows K and IPIV(K) of matrix B in reverse order
352* from the formation order of IPIV(I) vector for Lower case.
353*
354* (We can do the simple loop over IPIV with decrement -1,
355* since the ABS value of IPIV(I) represents the row index
356* of the interchange with row i in both 1x1 and 2x2 pivot cases)
357*
358 DO k = n, 1, -1
359 kp = abs( ipiv( k ) )
360 IF( kp.NE.k ) THEN
361 CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
362 END IF
363 END DO
364*
365* END Lower
366*
367 END IF
368*
369 RETURN
370*
371* End of CHETRS_3
372*
373 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine chetrs_3(uplo, n, nrhs, a, lda, e, ipiv, b, ldb, info)
CHETRS_3
Definition chetrs_3.f:163
subroutine csscal(n, sa, cx, incx)
CSSCAL
Definition csscal.f:78
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