LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zptts2.f
Go to the documentation of this file.
1*> \brief \b ZPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZPTTS2 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zptts2.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zptts2.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zptts2.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZPTTS2( IUPLO, N, NRHS, D, E, B, LDB )
20*
21* .. Scalar Arguments ..
22* INTEGER IUPLO, LDB, N, NRHS
23* ..
24* .. Array Arguments ..
25* DOUBLE PRECISION D( * )
26* COMPLEX*16 B( LDB, * ), E( * )
27* ..
28*
29*
30*> \par Purpose:
31* =============
32*>
33*> \verbatim
34*>
35*> ZPTTS2 solves a tridiagonal system of the form
36*> A * X = B
37*> using the factorization A = U**H *D*U or A = L*D*L**H computed by ZPTTRF.
38*> D is a diagonal matrix specified in the vector D, U (or L) is a unit
39*> bidiagonal matrix whose superdiagonal (subdiagonal) is specified in
40*> the vector E, and X and B are N by NRHS matrices.
41*> \endverbatim
42*
43* Arguments:
44* ==========
45*
46*> \param[in] IUPLO
47*> \verbatim
48*> IUPLO is INTEGER
49*> Specifies the form of the factorization and whether the
50*> vector E is the superdiagonal of the upper bidiagonal factor
51*> U or the subdiagonal of the lower bidiagonal factor L.
52*> = 1: A = U**H *D*U, E is the superdiagonal of U
53*> = 0: A = L*D*L**H, E is the subdiagonal of L
54*> \endverbatim
55*>
56*> \param[in] N
57*> \verbatim
58*> N is INTEGER
59*> The order of the tridiagonal 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] D
70*> \verbatim
71*> D is DOUBLE PRECISION array, dimension (N)
72*> The n diagonal elements of the diagonal matrix D from the
73*> factorization A = U**H *D*U or A = L*D*L**H.
74*> \endverbatim
75*>
76*> \param[in] E
77*> \verbatim
78*> E is COMPLEX*16 array, dimension (N-1)
79*> If IUPLO = 1, the (n-1) superdiagonal elements of the unit
80*> bidiagonal factor U from the factorization A = U**H*D*U.
81*> If IUPLO = 0, the (n-1) subdiagonal elements of the unit
82*> bidiagonal factor L from the factorization A = L*D*L**H.
83*> \endverbatim
84*>
85*> \param[in,out] B
86*> \verbatim
87*> B is COMPLEX*16 array, dimension (LDB,NRHS)
88*> On entry, the right hand side vectors B for the system of
89*> linear equations.
90*> On exit, the solution vectors, X.
91*> \endverbatim
92*>
93*> \param[in] LDB
94*> \verbatim
95*> LDB is INTEGER
96*> The leading dimension of the array B. LDB >= max(1,N).
97*> \endverbatim
98*
99* Authors:
100* ========
101*
102*> \author Univ. of Tennessee
103*> \author Univ. of California Berkeley
104*> \author Univ. of Colorado Denver
105*> \author NAG Ltd.
106*
107*> \ingroup ptts2
108*
109* =====================================================================
110 SUBROUTINE zptts2( IUPLO, N, NRHS, D, E, B, LDB )
111*
112* -- LAPACK computational routine --
113* -- LAPACK is a software package provided by Univ. of Tennessee, --
114* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
115*
116* .. Scalar Arguments ..
117 INTEGER IUPLO, LDB, N, NRHS
118* ..
119* .. Array Arguments ..
120 DOUBLE PRECISION D( * )
121 COMPLEX*16 B( LDB, * ), E( * )
122* ..
123*
124* =====================================================================
125*
126* .. Local Scalars ..
127 INTEGER I, J
128* ..
129* .. External Subroutines ..
130 EXTERNAL zdscal
131* ..
132* .. Intrinsic Functions ..
133 INTRINSIC dconjg
134* ..
135* .. Executable Statements ..
136*
137* Quick return if possible
138*
139 IF( n.LE.1 ) THEN
140 IF( n.EQ.1 )
141 $ CALL zdscal( nrhs, 1.d0 / d( 1 ), b, ldb )
142 RETURN
143 END IF
144*
145 IF( iuplo.EQ.1 ) THEN
146*
147* Solve A * X = B using the factorization A = U**H *D*U,
148* overwriting each right hand side vector with its solution.
149*
150 IF( nrhs.LE.2 ) THEN
151 j = 1
152 10 CONTINUE
153*
154* Solve U**H * x = b.
155*
156 DO 20 i = 2, n
157 b( i, j ) = b( i, j ) - b( i-1, j )*dconjg( e( i-1 ) )
158 20 CONTINUE
159*
160* Solve D * U * x = b.
161*
162 DO 30 i = 1, n
163 b( i, j ) = b( i, j ) / d( i )
164 30 CONTINUE
165 DO 40 i = n - 1, 1, -1
166 b( i, j ) = b( i, j ) - b( i+1, j )*e( i )
167 40 CONTINUE
168 IF( j.LT.nrhs ) THEN
169 j = j + 1
170 GO TO 10
171 END IF
172 ELSE
173 DO 70 j = 1, nrhs
174*
175* Solve U**H * x = b.
176*
177 DO 50 i = 2, n
178 b( i, j ) = b( i, j ) - b( i-1, j )*dconjg( e( i-1 ) )
179 50 CONTINUE
180*
181* Solve D * U * x = b.
182*
183 b( n, j ) = b( n, j ) / d( n )
184 DO 60 i = n - 1, 1, -1
185 b( i, j ) = b( i, j ) / d( i ) - b( i+1, j )*e( i )
186 60 CONTINUE
187 70 CONTINUE
188 END IF
189 ELSE
190*
191* Solve A * X = B using the factorization A = L*D*L**H,
192* overwriting each right hand side vector with its solution.
193*
194 IF( nrhs.LE.2 ) THEN
195 j = 1
196 80 CONTINUE
197*
198* Solve L * x = b.
199*
200 DO 90 i = 2, n
201 b( i, j ) = b( i, j ) - b( i-1, j )*e( i-1 )
202 90 CONTINUE
203*
204* Solve D * L**H * x = b.
205*
206 DO 100 i = 1, n
207 b( i, j ) = b( i, j ) / d( i )
208 100 CONTINUE
209 DO 110 i = n - 1, 1, -1
210 b( i, j ) = b( i, j ) - b( i+1, j )*dconjg( e( i ) )
211 110 CONTINUE
212 IF( j.LT.nrhs ) THEN
213 j = j + 1
214 GO TO 80
215 END IF
216 ELSE
217 DO 140 j = 1, nrhs
218*
219* Solve L * x = b.
220*
221 DO 120 i = 2, n
222 b( i, j ) = b( i, j ) - b( i-1, j )*e( i-1 )
223 120 CONTINUE
224*
225* Solve D * L**H * x = b.
226*
227 b( n, j ) = b( n, j ) / d( n )
228 DO 130 i = n - 1, 1, -1
229 b( i, j ) = b( i, j ) / d( i ) -
230 $ b( i+1, j )*dconjg( e( i ) )
231 130 CONTINUE
232 140 CONTINUE
233 END IF
234 END IF
235*
236 RETURN
237*
238* End of ZPTTS2
239*
240 END
subroutine zptts2(iuplo, n, nrhs, d, e, b, ldb)
ZPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf...
Definition zptts2.f:111
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78