LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ cpbstf()

subroutine cpbstf ( character  uplo,
integer  n,
integer  kd,
complex, dimension( ldab, * )  ab,
integer  ldab,
integer  info 
)

CPBSTF

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

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

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

 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 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.
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 152 of file cpbstf.f.

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