LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zhesv_rook.f
Go to the documentation of this file.
1*> \brief \b ZHESV_ROOK computes the solution to a system of linear equations A * X = B for HE matrices using the bounded Bunch-Kaufman ("rook") diagonal pivoting method
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download ZHESV_ROOK + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zhesv_rook.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zhesv_rook.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zhesv_rook.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE ZHESV_ROOK( UPLO, N, NRHS, A, LDA, IPIV, B, LDB, WORK,
22* LWORK, INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER UPLO
26* INTEGER INFO, LDA, LDB, LWORK, N, NRHS
27* ..
28* .. Array Arguments ..
29* INTEGER IPIV( * )
30* COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> ZHESV_ROOK computes the solution to a complex system of linear equations
40*> A * X = B,
41*> where A is an N-by-N Hermitian matrix and X and B are N-by-NRHS
42*> matrices.
43*>
44*> The bounded Bunch-Kaufman ("rook") diagonal pivoting method is used
45*> to factor A as
46*> A = U * D * U**T, if UPLO = 'U', or
47*> A = L * D * L**T, if UPLO = 'L',
48*> where U (or L) is a product of permutation and unit upper (lower)
49*> triangular matrices, and D is Hermitian and block diagonal with
50*> 1-by-1 and 2-by-2 diagonal blocks.
51*>
52*> ZHETRF_ROOK is called to compute the factorization of a complex
53*> Hermition matrix A using the bounded Bunch-Kaufman ("rook") diagonal
54*> pivoting method.
55*>
56*> The factored form of A is then used to solve the system
57*> of equations A * X = B by calling ZHETRS_ROOK (uses BLAS 2).
58*> \endverbatim
59*
60* Arguments:
61* ==========
62*
63*> \param[in] UPLO
64*> \verbatim
65*> UPLO is CHARACTER*1
66*> = 'U': Upper triangle of A is stored;
67*> = 'L': Lower triangle of A is stored.
68*> \endverbatim
69*>
70*> \param[in] N
71*> \verbatim
72*> N is INTEGER
73*> The number of linear equations, i.e., the order of the
74*> matrix A. N >= 0.
75*> \endverbatim
76*>
77*> \param[in] NRHS
78*> \verbatim
79*> NRHS is INTEGER
80*> The number of right hand sides, i.e., the number of columns
81*> of the matrix B. NRHS >= 0.
82*> \endverbatim
83*>
84*> \param[in,out] A
85*> \verbatim
86*> A is COMPLEX*16 array, dimension (LDA,N)
87*> On entry, the Hermitian matrix A. If UPLO = 'U', the leading
88*> N-by-N upper triangular part of A contains the upper
89*> triangular part of the matrix A, and the strictly lower
90*> triangular part of A is not referenced. If UPLO = 'L', the
91*> leading N-by-N lower triangular part of A contains the lower
92*> triangular part of the matrix A, and the strictly upper
93*> triangular part of A is not referenced.
94*>
95*> On exit, if INFO = 0, the block diagonal matrix D and the
96*> multipliers used to obtain the factor U or L from the
97*> factorization A = U*D*U**H or A = L*D*L**H as computed by
98*> ZHETRF_ROOK.
99*> \endverbatim
100*>
101*> \param[in] LDA
102*> \verbatim
103*> LDA is INTEGER
104*> The leading dimension of the array A. LDA >= max(1,N).
105*> \endverbatim
106*>
107*> \param[out] IPIV
108*> \verbatim
109*> IPIV is INTEGER array, dimension (N)
110*> Details of the interchanges and the block structure of D.
111*>
112*> If UPLO = 'U':
113*> Only the last KB elements of IPIV are set.
114*>
115*> If IPIV(k) > 0, then rows and columns k and IPIV(k) were
116*> interchanged and D(k,k) is a 1-by-1 diagonal block.
117*>
118*> If IPIV(k) < 0 and IPIV(k-1) < 0, then rows and
119*> columns k and -IPIV(k) were interchanged and rows and
120*> columns k-1 and -IPIV(k-1) were inerchaged,
121*> D(k-1:k,k-1:k) is a 2-by-2 diagonal block.
122*>
123*> If UPLO = 'L':
124*> Only the first KB elements of IPIV are set.
125*>
126*> If IPIV(k) > 0, then rows and columns k and IPIV(k)
127*> were interchanged and D(k,k) is a 1-by-1 diagonal block.
128*>
129*> If IPIV(k) < 0 and IPIV(k+1) < 0, then rows and
130*> columns k and -IPIV(k) were interchanged and rows and
131*> columns k+1 and -IPIV(k+1) were inerchaged,
132*> D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
133*> \endverbatim
134*>
135*> \param[in,out] B
136*> \verbatim
137*> B is COMPLEX*16 array, dimension (LDB,NRHS)
138*> On entry, the N-by-NRHS right hand side matrix B.
139*> On exit, if INFO = 0, the N-by-NRHS solution matrix X.
140*> \endverbatim
141*>
142*> \param[in] LDB
143*> \verbatim
144*> LDB is INTEGER
145*> The leading dimension of the array B. LDB >= max(1,N).
146*> \endverbatim
147*>
148*> \param[out] WORK
149*> \verbatim
150*> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
151*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
152*> \endverbatim
153*>
154*> \param[in] LWORK
155*> \verbatim
156*> LWORK is INTEGER
157*> The length of WORK. LWORK >= 1, and for best performance
158*> LWORK >= max(1,N*NB), where NB is the optimal blocksize for
159*> ZHETRF_ROOK.
160*> for LWORK < N, TRS will be done with Level BLAS 2
161*> for LWORK >= N, TRS will be done with Level BLAS 3
162*>
163*> If LWORK = -1, then a workspace query is assumed; the routine
164*> only calculates the optimal size of the WORK array, returns
165*> this value as the first entry of the WORK array, and no error
166*> message related to LWORK is issued by XERBLA.
167*> \endverbatim
168*>
169*> \param[out] INFO
170*> \verbatim
171*> INFO is INTEGER
172*> = 0: successful exit
173*> < 0: if INFO = -i, the i-th argument had an illegal value
174*> > 0: if INFO = i, D(i,i) is exactly zero. The factorization
175*> has been completed, but the block diagonal matrix D is
176*> exactly singular, so the solution could not be computed.
177*> \endverbatim
178*
179* Authors:
180* ========
181*
182*> \author Univ. of Tennessee
183*> \author Univ. of California Berkeley
184*> \author Univ. of Colorado Denver
185*> \author NAG Ltd.
186*
187*> \ingroup hesv_rook
188*>
189*> \verbatim
190*>
191*> November 2013, Igor Kozachenko,
192*> Computer Science Division,
193*> University of California, Berkeley
194*>
195*> September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
196*> School of Mathematics,
197*> University of Manchester
198*>
199*> \endverbatim
200*
201*
202* =====================================================================
203 SUBROUTINE zhesv_rook( UPLO, N, NRHS, A, LDA, IPIV, B, LDB, WORK,
204 $ LWORK, INFO )
205*
206* -- LAPACK driver routine --
207* -- LAPACK is a software package provided by Univ. of Tennessee, --
208* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
209*
210* .. Scalar Arguments ..
211 CHARACTER UPLO
212 INTEGER INFO, LDA, LDB, LWORK, N, NRHS
213* ..
214* .. Array Arguments ..
215 INTEGER IPIV( * )
216 COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * )
217* ..
218*
219* =====================================================================
220*
221* .. Local Scalars ..
222 LOGICAL LQUERY
223 INTEGER LWKOPT, NB
224* ..
225* .. External Functions ..
226 LOGICAL LSAME
227 INTEGER ILAENV
228 EXTERNAL lsame, ilaenv
229* ..
230* .. External Subroutines ..
232* ..
233* .. Intrinsic Functions ..
234 INTRINSIC max
235* ..
236* .. Executable Statements ..
237*
238* Test the input parameters.
239*
240 info = 0
241 lquery = ( lwork.EQ.-1 )
242 IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
243 info = -1
244 ELSE IF( n.LT.0 ) THEN
245 info = -2
246 ELSE IF( nrhs.LT.0 ) THEN
247 info = -3
248 ELSE IF( lda.LT.max( 1, n ) ) THEN
249 info = -5
250 ELSE IF( ldb.LT.max( 1, n ) ) THEN
251 info = -8
252 ELSE IF( lwork.LT.1 .AND. .NOT.lquery ) THEN
253 info = -10
254 END IF
255*
256 IF( info.EQ.0 ) THEN
257 IF( n.EQ.0 ) THEN
258 lwkopt = 1
259 ELSE
260 nb = ilaenv( 1, 'ZHETRF_ROOK', uplo, n, -1, -1, -1 )
261 lwkopt = n*nb
262 END IF
263 work( 1 ) = lwkopt
264 END IF
265*
266 IF( info.NE.0 ) THEN
267 CALL xerbla( 'ZHESV_ROOK ', -info )
268 RETURN
269 ELSE IF( lquery ) THEN
270 RETURN
271 END IF
272*
273* Compute the factorization A = U*D*U**H or A = L*D*L**H.
274*
275 CALL zhetrf_rook( uplo, n, a, lda, ipiv, work, lwork, info )
276 IF( info.EQ.0 ) THEN
277*
278* Solve the system A*X = B, overwriting B with X.
279*
280* Solve with TRS ( Use Level BLAS 2)
281*
282 CALL zhetrs_rook( uplo, n, nrhs, a, lda, ipiv, b, ldb, info )
283*
284 END IF
285*
286 work( 1 ) = lwkopt
287*
288 RETURN
289*
290* End of ZHESV_ROOK
291*
292 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhesv_rook(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
ZHESV_ROOK computes the solution to a system of linear equations A * X = B for HE matrices using the ...
Definition zhesv_rook.f:205
subroutine zhetrf_rook(uplo, n, a, lda, ipiv, work, lwork, info)
ZHETRF_ROOK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bun...
subroutine zhetrs_rook(uplo, n, nrhs, a, lda, ipiv, b, ldb, info)
ZHETRS_ROOK computes the solution to a system of linear equations A * X = B for HE matrices using fac...