LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zpbstf ( character  UPLO,
integer  N,
integer  KD,
complex*16, dimension( ldab, * )  AB,
integer  LDAB,
integer  INFO 
)

ZPBSTF

Download ZPBSTF + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 ZPBSTF computes a split Cholesky factorization of a complex
 Hermitian positive definite band matrix A.

 This routine is designed to be used in conjunction with ZHBGST.

 The factorization has the form  A = S**H*S  where S is a band matrix
 of the same bandwidth as A and the following structure:

   S = ( U    )
       ( M  L )

 where U is upper triangular of order m = (n+kd)/2, and L is lower
 triangular of order n-m.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]KD
          KD is INTEGER
          The number of superdiagonals of the matrix A if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'.  KD >= 0.
[in,out]AB
          AB is COMPLEX*16 array, dimension (LDAB,N)
          On entry, the upper or lower triangle of the Hermitian band
          matrix A, stored in the first kd+1 rows of the array.  The
          j-th column of A is stored in the j-th column of the array AB
          as follows:
          if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+kd).

          On exit, if INFO = 0, the factor S from the split Cholesky
          factorization A = S**H*S. See Further Details.
[in]LDAB
          LDAB is INTEGER
          The leading dimension of the array AB.  LDAB >= KD+1.
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -i, the i-th argument had an illegal value
          > 0: if INFO = i, the factorization could not be completed,
               because the updated element a(i,i) was negative; the
               matrix A is not positive definite.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
Further Details:
  The band storage scheme is illustrated by the following example, when
  N = 7, KD = 2:

  S = ( s11  s12  s13                     )
      (      s22  s23  s24                )
      (           s33  s34                )
      (                s44                )
      (           s53  s54  s55           )
      (                s64  s65  s66      )
      (                     s75  s76  s77 )

  If UPLO = 'U', the array AB holds:

  on entry:                          on exit:

   *    *   a13  a24  a35  a46  a57   *    *   s13  s24  s53**H s64**H s75**H
   *   a12  a23  a34  a45  a56  a67   *   s12  s23  s34  s54**H s65**H s76**H
  a11  a22  a33  a44  a55  a66  a77  s11  s22  s33  s44  s55    s66    s77

  If UPLO = 'L', the array AB holds:

  on entry:                          on exit:

  a11  a22  a33  a44  a55  a66  a77  s11    s22    s33    s44  s55  s66  s77
  a21  a32  a43  a54  a65  a76   *   s12**H s23**H s34**H s54  s65  s76   *
  a31  a42  a53  a64  a64   *    *   s13**H s24**H s53    s64  s75   *    *

  Array elements marked * are not used by the routine; s12**H denotes
  conjg(s12); the diagonal elements of S are real.

Definition at line 155 of file zpbstf.f.

155 *
156 * -- LAPACK computational routine (version 3.4.0) --
157 * -- LAPACK is a software package provided by Univ. of Tennessee, --
158 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
159 * November 2011
160 *
161 * .. Scalar Arguments ..
162  CHARACTER uplo
163  INTEGER info, kd, ldab, n
164 * ..
165 * .. Array Arguments ..
166  COMPLEX*16 ab( ldab, * )
167 * ..
168 *
169 * =====================================================================
170 *
171 * .. Parameters ..
172  DOUBLE PRECISION one, zero
173  parameter ( one = 1.0d+0, zero = 0.0d+0 )
174 * ..
175 * .. Local Scalars ..
176  LOGICAL upper
177  INTEGER j, kld, km, m
178  DOUBLE PRECISION ajj
179 * ..
180 * .. External Functions ..
181  LOGICAL lsame
182  EXTERNAL lsame
183 * ..
184 * .. External Subroutines ..
185  EXTERNAL xerbla, zdscal, zher, zlacgv
186 * ..
187 * .. Intrinsic Functions ..
188  INTRINSIC dble, max, min, sqrt
189 * ..
190 * .. Executable Statements ..
191 *
192 * Test the input parameters.
193 *
194  info = 0
195  upper = lsame( uplo, 'U' )
196  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
197  info = -1
198  ELSE IF( n.LT.0 ) THEN
199  info = -2
200  ELSE IF( kd.LT.0 ) THEN
201  info = -3
202  ELSE IF( ldab.LT.kd+1 ) THEN
203  info = -5
204  END IF
205  IF( info.NE.0 ) THEN
206  CALL xerbla( 'ZPBSTF', -info )
207  RETURN
208  END IF
209 *
210 * Quick return if possible
211 *
212  IF( n.EQ.0 )
213  $ RETURN
214 *
215  kld = max( 1, ldab-1 )
216 *
217 * Set the splitting point m.
218 *
219  m = ( n+kd ) / 2
220 *
221  IF( upper ) THEN
222 *
223 * Factorize A(m+1:n,m+1:n) as L**H*L, and update A(1:m,1:m).
224 *
225  DO 10 j = n, m + 1, -1
226 *
227 * Compute s(j,j) and test for non-positive-definiteness.
228 *
229  ajj = dble( ab( kd+1, j ) )
230  IF( ajj.LE.zero ) THEN
231  ab( kd+1, j ) = ajj
232  GO TO 50
233  END IF
234  ajj = sqrt( ajj )
235  ab( kd+1, j ) = ajj
236  km = min( j-1, kd )
237 *
238 * Compute elements j-km:j-1 of the j-th column and update the
239 * the leading submatrix within the band.
240 *
241  CALL zdscal( km, one / ajj, ab( kd+1-km, j ), 1 )
242  CALL zher( 'Upper', km, -one, ab( kd+1-km, j ), 1,
243  $ ab( kd+1, j-km ), kld )
244  10 CONTINUE
245 *
246 * Factorize the updated submatrix A(1:m,1:m) as U**H*U.
247 *
248  DO 20 j = 1, m
249 *
250 * Compute s(j,j) and test for non-positive-definiteness.
251 *
252  ajj = dble( ab( kd+1, j ) )
253  IF( ajj.LE.zero ) THEN
254  ab( kd+1, j ) = ajj
255  GO TO 50
256  END IF
257  ajj = sqrt( ajj )
258  ab( kd+1, j ) = ajj
259  km = min( kd, m-j )
260 *
261 * Compute elements j+1:j+km of the j-th row and update the
262 * trailing submatrix within the band.
263 *
264  IF( km.GT.0 ) THEN
265  CALL zdscal( km, one / ajj, ab( kd, j+1 ), kld )
266  CALL zlacgv( km, ab( kd, j+1 ), kld )
267  CALL zher( 'Upper', km, -one, ab( kd, j+1 ), kld,
268  $ ab( kd+1, j+1 ), kld )
269  CALL zlacgv( km, ab( kd, j+1 ), kld )
270  END IF
271  20 CONTINUE
272  ELSE
273 *
274 * Factorize A(m+1:n,m+1:n) as L**H*L, and update A(1:m,1:m).
275 *
276  DO 30 j = n, m + 1, -1
277 *
278 * Compute s(j,j) and test for non-positive-definiteness.
279 *
280  ajj = dble( ab( 1, j ) )
281  IF( ajj.LE.zero ) THEN
282  ab( 1, j ) = ajj
283  GO TO 50
284  END IF
285  ajj = sqrt( ajj )
286  ab( 1, j ) = ajj
287  km = min( j-1, kd )
288 *
289 * Compute elements j-km:j-1 of the j-th row and update the
290 * trailing submatrix within the band.
291 *
292  CALL zdscal( km, one / ajj, ab( km+1, j-km ), kld )
293  CALL zlacgv( km, ab( km+1, j-km ), kld )
294  CALL zher( 'Lower', km, -one, ab( km+1, j-km ), kld,
295  $ ab( 1, j-km ), kld )
296  CALL zlacgv( km, ab( km+1, j-km ), kld )
297  30 CONTINUE
298 *
299 * Factorize the updated submatrix A(1:m,1:m) as U**H*U.
300 *
301  DO 40 j = 1, m
302 *
303 * Compute s(j,j) and test for non-positive-definiteness.
304 *
305  ajj = dble( ab( 1, j ) )
306  IF( ajj.LE.zero ) THEN
307  ab( 1, j ) = ajj
308  GO TO 50
309  END IF
310  ajj = sqrt( ajj )
311  ab( 1, j ) = ajj
312  km = min( kd, m-j )
313 *
314 * Compute elements j+1:j+km of the j-th column and update the
315 * trailing submatrix within the band.
316 *
317  IF( km.GT.0 ) THEN
318  CALL zdscal( km, one / ajj, ab( 2, j ), 1 )
319  CALL zher( 'Lower', km, -one, ab( 2, j ), 1,
320  $ ab( 1, j+1 ), kld )
321  END IF
322  40 CONTINUE
323  END IF
324  RETURN
325 *
326  50 CONTINUE
327  info = j
328  RETURN
329 *
330 * End of ZPBSTF
331 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zher(UPLO, N, ALPHA, X, INCX, A, LDA)
ZHER
Definition: zher.f:137
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zlacgv(N, X, INCX)
ZLACGV conjugates a complex vector.
Definition: zlacgv.f:76

Here is the call graph for this function:

Here is the caller graph for this function: