1c Copyright (C) 2008, 2009  VZLU Prague, a.s., Czech Republic
2c
3c Author: Jaroslav Hajek <highegg@gmail.com>
4c
5c This file is part of qrupdate.
6c
7c qrupdate is free software; you can redistribute it and/or modify
8c it under the terms of the GNU General Public License as published by
9c the Free Software Foundation; either version 3 of the License, or
10c (at your option) any later version.
11c
12c This program is distributed in the hope that it will be useful,
13c but WITHOUT ANY WARRANTY; without even the implied warranty of
14c MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
15c GNU General Public License for more details.
16c
17c You should have received a copy of the GNU General Public License
18c along with this software; see the file COPYING.  If not, see
19c <http://www.gnu.org/licenses/>.
20c
21      subroutine zch1up(n,R,ldr,u,w)
22c purpose:      given an upper triangular matrix R that is a Cholesky
23c               factor of a hermitian positive definite matrix A, i.e.
24c               A = R'*R, this subroutine updates R -> R1 so that
25c               R1'*R1 = A + u*u'
26c               (complex version)
27c arguments:
28c n (in)        the order of matrix R
29c R (io)        on entry, the upper triangular matrix R
30c               on exit, the updated matrix R1
31c ldr (in)      leading dimension of R. ldr >= n.
32c u (io)        the vector determining the rank-1 update
33c               on exit, u contains the rotation sines
34c               used to transform R to R1.
35c w (out)       cosine parts of rotations.
36c
37      integer n,ldr
38      double complex R(ldr,*),u(*)
39      double precision w(*)
40      external zlartg
41      double complex rr,ui,t
42      integer i,j
43
44      do i = 1,n
45c apply stored rotations, column-wise
46        ui = conjg(u(i))
47        do j = 1,i-1
48          t = w(j)*R(j,i) + u(j)*ui
49          ui = w(j)*ui - conjg(u(j))*R(j,i)
50          R(j,i) = t
51        end do
52c generate next rotation
53        call zlartg(R(i,i),ui,w(i),u(i),rr)
54        R(i,i) = rr
55      end do
56      end subroutine
57
58