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 sqrshc(m,n,k,Q,ldq,R,ldr,i,j,w)
22c purpose:      updates a QR factorization after circular shift of
23c               columns.
24c               i.e., given an m-by-k orthogonal matrix Q, an k-by-n
25c               upper trapezoidal matrix R and index j in the range
26c               1:n+1, this subroutine updates the matrix Q -> Q1 and
27c               R -> R1 so that Q1 is again orthogonal, R1 upper
28c               trapezoidal, and
29c               Q1*R1 = A(:,p), where A = Q*R and p is the permutation
30c               [1:i-1,shift(i:j,-1),j+1:n] if i < j  or
31c               [1:j-1,shift(j:i,+1),i+1:n] if j < i.
32c               (real version)
33c arguments:
34c m (in)        number of rows of the matrix Q.
35c n (in)        number of columns of the matrix R.
36c k (in)        number of columns of Q1, and rows of R1. Must be
37c               either k = m (full Q) or k = n <= m (economical form).
38c Q (io)        on entry, the unitary m-by-k matrix Q.
39c               on exit, the updated matrix Q1.
40c ldq (in)      leading dimension of Q. ldq >= m.
41c R (io)        on entry, the original matrix R.
42c               on exit, the updated matrix R1.
43c ldr (in)      leading dimension of R. ldr >= k.
44c i (in)        the first index determining the range (see above)
45c j (in)        the second index determining the range (see above)
46c w (o)         a workspace vector of size 2*k.
47c
48      integer m,n,k,ldq,ldr,i,j
49      real Q(ldq,*),R(ldr,*),w(*)
50      external xerbla,scopy,sqrtv1,sqrqh,sqhqr
51      integer info,jj,kk,l
52c quick return if possible.
53      if (m == 0 .or. n == 1) return
54      info = 0
55c check arguments.
56      if (m < 0) then
57        info = 1
58      else if (n < 0) then
59        info = 2
60      else if (k /= m .and. (k /= n .or. n > m)) then
61        info = 3
62      else if (i < 1 .or. i > n) then
63        info = 6
64      else if (j < 1 .or. j > n) then
65        info = 7
66      end if
67      if (info /= 0) then
68        call xerbla('SQRSHC',info)
69        return
70      end if
71
72      if (i < j) then
73c shift columns
74        call scopy(k,R(1,i),1,w,1)
75        do l = i,j-1
76          call scopy(k,R(1,l+1),1,R(1,l),1)
77        end do
78        call scopy(k,w,1,R(1,j),1)
79c retriangularize
80        if (i < k) then
81          kk = min(k,j)
82          call sqhqr(kk+1-i,n+1-i,R(i,i),ldr,w(k+1),w)
83c apply rotations to Q.
84          call sqrot('F',m,kk+1-i,Q(1,i),ldq,w(k+1),w)
85        end if
86      else if (j < i) then
87c shift columns
88        call scopy(k,R(1,i),1,w,1)
89        do l = i,j+1,-1
90          call scopy(k,R(1,l-1),1,R(1,l),1)
91        end do
92        call scopy(k,w,1,R(1,j),1)
93c retriangularize
94        if (j < k) then
95          jj = min(j+1,n)
96          kk = min(k,i)
97c eliminate the introduced spike.
98          call sqrtv1(kk+1-j,R(j,j),w(k+1))
99c apply rotations to R
100          call sqrqh(kk+1-j,n-j,R(j,jj),ldr,w(k+1),R(j+1,j))
101c apply rotations to Q
102          call sqrot('B',m,kk+1-j,Q(1,j),ldq,w(k+1),R(j+1,j))
103c zero spike.
104          do l = j+1,kk
105            R(l,j) = 0e0
106          end do
107        end if
108      end if
109      end subroutine
110