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