NLSL
cscg.f90
Go to the documentation of this file.
1 c NLSL Version 1.5.1 beta 11/25/95
2 c**********************************************************************
3 c
4 c SUBROUTINE : CSCG
5 c -----------------
6 c
7 c (C)omplex (S)ymmetric (C)onjugate (G)radient Algorithm
8 c
9 c This subroutine attempts a solution of an Ax=b type problem
10 c where A is a known sparse complex symmetric matrix, b is a
11 c known vector and the approximate solution vector is obtained
12 c via a simple complex symmetric conjugate gradient procedure.
13 c The algorithm used here is based on algorithm 10.3-1 in
14 c "Matrix Computations", first edition, by G. Golub and C.
15 c Van Loan, Johns Hopkins Univ. Press, 1983. The "Lanczos"
16 c tridiagonal matrix in constructed using formula 10.2-14 in
17 c the reference above.
18 c
19 c It is assumed in this routine that the right hand vector is
20 c normalized and the solution vector is initialized at the start
21 c of the routine. The matrix and index arrays are passed to
22 c this routine through the inclusion of the file eprmat.inc.
23 c
24 c arguments
25 c ---------
26 c
27 c b : right hand vector
28 c ndim : number of rows and columns in the matrix
29 c mxcgs : maximum number of conjugate gradient steps
30 c allowed
31 c cgtol : maximum residual allowed for termination
32 c shift : a complex shift value used to avoid extraneous
33 c divisions by zero
34 c x : approximate solution vector
35 c al,bl : quantities that can be used to construct the
36 c "Lanczos" tridiagonal matrix
37 c ndone : number of cg steps executed (positive if
38 c converged, negative otherwise)
39 c error : pseudonorm of the residual vector at return
40 c
41 c a : diagonal elements of tridiagonal matrix
42 c b : off-diagonal elements of tridiagonal matrix
43 c z : elements of input matrix
44 c izmat : column number and indicator for ends of rows and
45 c real and imaginary parts of matrix elements of
46 c stored in z
47 c if izmat(i)<0 then z(i) is pure real
48 c if izmat(i)=0 then end of row
49 c if izmat(i)>0 then z(i) is pure imaginary
50 c zdiag : the diagonal elements of the matrix
51 c
52 c
53 c Includes:
54 c nlsdim.inc
55 c rndoff.inc
56 c eprmat.inc
57 c
58 c Uses:
59 c zdotu.f
60 c zaypx.f
61 c scmvm.f
62 c zaxpy.f
63 c
64 c**********************************************************************
65 c
66  subroutine cscg(b,ndim,mxcgs,cgtol,shift,
67  # x,al,bl,ndone,error)
68 c
69  use nlsdim
70  use rnddbl
71  use eprmat
72  use stdio
73 c
74  integer ndim,mxcgs,ndone
75  double precision cgtol,error
76  complex*16 shift
77  complex*16 x,b,al,bl
78  dimension x(mxdim),b(mxdim),al(mxstep),bl(mxstep)
79  integer i,nstep
80 c
81  complex*16 rho1,rho2,CZERO
82  parameter(czero=(0.0d0,0.0d0))
83 c
84  complex*16 p,r,w,z,alpha,beta
85  common /cgdata/ p(mxdim),r(mxdim),w(mxdim),z(mxdim)
86 c
87  complex*16 zdotu
88  external zdotu
89 c
90 c######################################################################
91 c
92 c.......... Debugging purposes only ..........
93 c integer j
94 c open (unit=20,file='nlcg.tst',status='unknown',
95 c # access='sequential',form='formatted')
96 c write (20,7002)
97 c do i=1,ndim
98 c if (abs(b(i)).gt. RNDOFF) write (20,7003) i,real(b(i))
99 c end do
100 c
101 c write (20,7001) 'Imaginary'
102 c do i=1,ndim
103 c write (20,7000) i,i,zdiag(2,i)
104 c write (20,7000) (i,izmat(j),zmat(j),j=jzmat(i),jzmat(i+1)-1)
105 c end do
106 c
107 c write (20,7001) 'Real'
108 c do i=1,ndim
109 c write (20,7000) i,i,zdiag(1,i)
110 c do j=kzmat(i),kzmat(i+1)-1
111 c k=MXEL-j+1
112 c write (20,7000) i,izmat(k),zmat(k)
113 c end do
114 c end do
115 c close(20)
116 c 7000 format(i5,',',i5,2x,g14.7)
117 c 7001 format('c --- ',a,' matrix elements')
118 c 7002 format('c --- vector elements')
119 c 7003 format(i5,2x,g14.7)
120 c........................................
121 c
122 c----------------------------------------------------------------------
123 c setup r,x and p
124 c----------------------------------------------------------------------
125 c
126  do i=1,ndim
127  r(i)=b(i)
128  end do
129 c
130  do i=1,ndim
131  p(i)=b(i)
132  end do
133 c
134  do i=1,ndim
135  x(i)=czero
136  end do
137 c
138  rho1=zdotu(r,r,ndim)
139  rho2=czero
140 c
141 c======================================================================
142 c begin loop over CG steps
143 c======================================================================
144 c
145  nstep=0
146  100 nstep=nstep+1
147 c
148 c----------------------------------------------------------------------
149 c calculate beta
150 c----------------------------------------------------------------------
151 c
152  if (nstep.eq.1) then
153  beta=czero
154  else
155  rho2=rho1
156  rho1=zdotu(r,r,ndim)
157  beta=rho1/rho2
158  end if
159 c
160  bl(nstep)=beta
161 c
162 c----------------------------------------------------------------------
163 c check for convergence
164 c----------------------------------------------------------------------
165 c
166  error=sqrt(abs(rho1))
167 c
168  if (error.lt.cgtol) go to 200
169 c
170 c----------------------------------------------------------------------
171 c update p ( p <- r+beta*p )
172 c----------------------------------------------------------------------
173 c
174  if (nstep.ne.1) call zaypx(r,p,beta,ndim)
175 c
176 c----------------------------------------------------------------------
177 c calculate w ( w <- A*p )
178 c----------------------------------------------------------------------
179 c
180 c
181 c Check for window input and/or user halt before proceeding
182 c
183  call wpoll
184  if (hltcmd.ne.0 .or. hltfit.ne.0) return
185 
186  call scmvm(p,w,ndim)
187  call zaxpy(p,w,shift,ndim)
188 c
189 c----------------------------------------------------------------------
190 c calculate alpha
191 c----------------------------------------------------------------------
192 c
193  alpha=rho1/zdotu(p,w,ndim)
194 c
195  al(nstep)=1.0d0/alpha
196 c
197 c----------------------------------------------------------------------
198 c update x ( x <- x+alpha*p )
199 c----------------------------------------------------------------------
200 c
201  call zaxpy(p,x,alpha,ndim)
202 c
203 c----------------------------------------------------------------------
204 c update r ( r <- r-alpha*w )
205 c----------------------------------------------------------------------
206 c
207  call zaxpy(w,r,-alpha,ndim)
208 c
209 c======================================================================
210 c end of loop over CG steps
211 c======================================================================
212 c
213  if (nstep.lt.mxcgs) go to 100
214 c
215 c----------------------------------------------------------------------
216 c construct Lanczos tridiagonal matrix from CG quantities
217 c----------------------------------------------------------------------
218 c
219  200 continue
220 c
221  call cgltri(nstep,al,bl,dreal(shift),dimag(shift))
222 c
223 c----------------------------------------------------------------------
224 c return error code if not converged
225 c----------------------------------------------------------------------
226 c
227  if (error.le.cgtol) then
228  ndone=nstep
229  else
230  ndone=-nstep
231  end if
232 c
233 c----------------------------------------------------------------------
234 c return to calling routine
235 c----------------------------------------------------------------------
236 c
237  return
238  end
void FORTRAN() wpoll()
Definition: pltx.c:860
subroutine cgltri(ndone, a, b, shiftr, shifti)
Definition: cgltri.f90:29
integer, parameter mxstep
Definition: nlsdim.f90:39
integer, save hltcmd
Definition: stdio.f90:32
integer, save hltfit
Definition: stdio.f90:32
subroutine cscg(b, ndim, mxcgs, cgtol, shift, x, al, bl, ndone, error)
Definition: cscg.f90:68
Definition: stdio.f90:26
subroutine zaxpy(x, y, scale, ndim)
Definition: zaxpy.f90:24
subroutine zaypx(x, y, scale, ndim)
Definition: zaypx.f90:19
subroutine scmvm(x, y, ndim)
Definition: scmvm.f90:87
integer, parameter mxdim
Definition: nlsdim.f90:39