NLSL
stvect.f90
Go to the documentation of this file.
1 c NLSL Version 1.4 10/10/94
2 c**********************************************************************
3 c ===================
4 c SUBROUTINE : STVECT
5 c ===================
6 c
7 c subroutine for calculating the starting vector for eprll
8 c
9 c This differs from the standalone version in the addition of
10 c the parameter ierr. It is used to report non-convergence of
11 c Bessel function evaluations instead of halting when such
12 c errors occur.
13 c
14 c Notes:
15 c 1) The orientational integral is evaluated numerically
16 c using the Curtis-Clenshaw-Romberg extrapolation algorithm.
17 c For details, see Bruno's thesis.
18 c
19 c written by DJS 10-SEP-87
20 c
21 c Includes:
22 c nlsdim.inc
23 c eprprm.inc
24 c rndoff.inc
25 c
26 c Uses:
27 c ipar.f
28 c ccrint.f
29 c
30 c**********************************************************************
31 c
32  subroutine stvect(v,ierr)
33 c
34  use nlsdim
35  use eprprm
36  use errmsg
37  use rnddbl
38 c
39  implicit none
40  double precision v(2,mxdim)
41  integer ierr
42 c
43  logical evenk,fmpi0
44  integer i,iparkr,iparlr,jkr,jmr,nup,id,nrow,krmn,krmx,krsgn,
45  # ipnr,ipnrmx,ipnrmn,ipnrsg,iqnr,iqnrmx,iqnrmn,mr,mrmn,mrmx,
46  # mrsgn
47 c
48  double precision cnl,stv,stvlk,stvm,factor,dsq2,vnorm
49 c
50  double precision ACCRCY,SMALL
51  parameter(accrcy=1.0d-8,small=1.0d-10)
52 c
53  double precision ONE,ZERO
54  parameter(one=1.0d0,zero=0.0d0)
55 c
56  integer lr,kr,iberr
57  common /ifzdat/lr,kr,iberr
58 c
59  integer ipar
60  double precision fz
61  external fz,ipar
62 c
63 c######################################################################
64 c
65 c.... Debugging purposes only!
66 c open (unit=20,file='nlvec.tst',status='unknown',
67 c # access='sequential',form='formatted')
68 c............................
69 c
70  do i=1,ndim
71  v(1,i)=zero
72  v(2,i)=zero
73  end do
74 c
75  dsq2=sqrt(2.0d0)
76 c
77  nrow=0
78  nelv=0
79  vnorm=zero
80 c
81 c --------------------
82 c *** loop over lr ***
83 c --------------------
84  do 100 lr=0,lemx,ldelta
85  iparlr=ipar(lr)
86  if((lr.gt.lomx).and.(iparlr.ne.1)) go to 100
87 c
88  if (iparlr.eq.1) then
89  cnl=dsqrt(dble(2*lr+1))
90  else
91  cnl=zero
92  end if
93 c
94 c --------------------
95 c *** loop over kr ***
96 c --------------------
97  krmx=min(kmx,lr)
98  krmn=max(kmn,-lr)
99  do 110 krsgn=krmn,krmx,kdelta
100  kr=abs(krsgn)
101  iparkr=ipar(kr)
102  jkr=isign(1,krsgn)
103  if(kr.eq.0) jkr=iparlr
104  if(jkr.lt.jkmn) goto 110
105 c
106 c --- Only vector elements are for even L, K and jK=1
107 c
108  if(iparlr.eq.1 .and. iparkr.eq.1 .and. jkr.eq.1) then
109 c
110 c --- No potential: only elements are for L,K=0
111 c
112  if(lptmx.eq.0) then
113  if((lr.eq.0).and.(kr.eq.0)) then
114  stvlk=one
115  else
116  stvlk=zero
117  end if
118 c
119 c --- Axial potential: Only elements are for K=0
120 c
121  else if((kptmx.eq.0).and.(kr.ne.0)) then
122  stvlk=zero
123 c
124 c --- Integration to find vector element
125 c
126  else
127  call ccrint(zero,one,accrcy,small,stvlk,nup,fz,id)
128  if (iberr.ne.0) ierr=badbess
129  if(kr.ne.0) then
130  factor=one
131  do i=lr-kr+1,lr+kr
132  factor=factor*dble(i)
133  end do
134  factor=one/dsqrt(factor)
135  else
136  factor=one/dsq2
137  end if
138  stvlk=factor*stvlk*cnl
139  end if
140  else
141  stvlk=zero
142  end if
143 c
144 c ----------------------
145 c *** loop over mr ***
146 c ----------------------
147  mrmx=min0(mmx,lr)
148  mrmn=max0(mmn,-lr)
149  do 120 mrsgn=mrmn,mrmx
150  mr=abs(mrsgn)
151  jmr=isign(1,mrsgn)
152 c
153  if (mr.ne.0) then
154  stvm=zero
155  else
156  stvm=stvlk
157  end if
158 c
159  if(ipsi0.eq.0) then
160  ipnrmn=mr
161  ipnrmx=mr
162  else
163  ipnrmx=min0(in2,ipnmx)
164  if(mr.eq.0 .and. jmmn.eq.1) then
165  ipnrmn=0
166  else
167  ipnrmn=-ipnrmx
168  end if
169  end if
170 c
171 c ----------------------
172 c *** loop over ipnr ***
173 c ----------------------
174  do 130 ipnrsg=ipnrmn,ipnrmx
175  ipnr=abs(ipnrsg)
176  if (mr.eq.0) then
177  jmr=isign(1,ipnrsg)
178  if (ipnrsg.eq.0) jmr=iparlr
179  end if
180  if(jmr.lt.jmmn) goto 130
181 c
182  if (ipnr.ne.0 .or. jmr.ne.1) then
183  stv=zero
184  else
185  stv=stvm
186  end if
187 c
188  iqnrmx=in2-iabs(ipnr)
189  iqnrmn=-iqnrmx
190 c
191 c ----------------------
192 c *** loop over iqnr ***
193 c ----------------------
194  do 140 iqnr=iqnrmn,iqnrmx,2
195 c
196 c======================================================================
197 c center of loops
198 c======================================================================
199 c
200  nrow=nrow+1
201  v(1,nrow)=stv
202  vnorm=vnorm+stv*stv
203 c
204 c.................... Debugging purposes only!
205 c if (abs(v(1,nrow)).gt.RNDOFF) then
206 c write(20,7000) lr,krsgn,mrsgn,ipnrsg,iqnr,v(1,nrow)
207 c 7000 format('Re <v|',4(i3,','),i3,'> = ',2g14.7)
208 c end if
209 c.............................................
210 c
211 c======================================================================
212 c
213 c
214  140 continue
215  130 continue
216  120 continue
217  110 continue
218  100 continue
219 c
220 c----------------------------------------------------------------------
221 c normalize starting vector and zero out imaginary part
222 c----------------------------------------------------------------------
223 c
224  nelv=0
225  vnorm=one/dsqrt(vnorm)
226  do i=1,ndim
227  v(1,i)=v(1,i)*vnorm
228  if(abs(v(1,i)).gt.rndoff) then
229  nelv=nelv+1
230  else
231  v(1,i)=zero
232  end if
233  v(2,i)=zero
234  end do
235 c
236 c----------------------------------------------------------------------
237 c zero out remainder of vector
238 c----------------------------------------------------------------------
239 c
240  do i=ndim+1,mxdim
241  v(1,i)=zero
242  v(2,i)=zero
243  end do
244 c
245 c..........Debugging purposes only!
246 c close(20)
247 c...................................
248 c
249  return
250  end
integer, pointer, save ndim
Definition: eprprm.f90:69
integer, parameter badbess
Definition: errmsg.f90:51
integer, pointer, save mmx
Definition: eprprm.f90:69
integer, pointer, save jkmn
Definition: eprprm.f90:69
integer, save lptmx
Definition: eprprm.f90:82
integer, pointer, save mmn
Definition: eprprm.f90:69
integer, pointer, save ipnmx
Definition: eprprm.f90:69
integer, save nelv
Definition: eprprm.f90:82
integer, save kptmx
Definition: eprprm.f90:82
subroutine ccrint(bndlow, bndhi, epsiln, small, sum, neval, f, id)
Definition: ccrints.f90:32
integer, save ldelta
Definition: eprprm.f90:82
integer, pointer, save lomx
Definition: eprprm.f90:69
integer, pointer, save kmn
Definition: eprprm.f90:69
integer, pointer, save in2
Definition: eprprm.f90:69
integer, pointer, save kmx
Definition: eprprm.f90:69
integer, pointer, save jmmn
Definition: eprprm.f90:69
integer, pointer, save lemx
Definition: eprprm.f90:69
subroutine stvect(v, ierr)
Definition: stvect.f90:33
integer, parameter mxdim
Definition: nlsdim.f90:39
integer, save ipsi0
Definition: eprprm.f90:82
double precision, parameter rndoff
Definition: rnddbl.f90:86
integer, save kdelta
Definition: eprprm.f90:82