NLSL
w3j.f90
Go to the documentation of this file.
1 c VERSION 1.1 1/12/93
2 c**********************************************************************
3 c
4 c ===================
5 c function W3J
6 c ===================
7 c
8 c This double precision function will calculate the values of
9 c the Wigner 3-J symbols used in the Stochastic Liouville matrix
10 c formulation of the slow-motional EPR calculation, i.e.
11 c
12 c ( J1 J2 J3 )
13 c ( M1 M2 M3 )
14 c
15 c For J2 <= 2 and arbitrary J1 and J3, this routine explicitly
16 c evaluates formulae given in Table 2 of "Angular Momentum in
17 c Quantum Mechanics" by A. R. Edmonds, revised 2nd printing,
18 c Princeton Univ. Press, 1968.
19 c
20 c For J2 > 2 and J1,J3 < Lmax (mxlval defined in 'maxl.inc'),
21 c this functions calls a modified version of code from
22 c Prof. C.C.J. Roothan, function wig3j, which appears in this file.
23 c
24 c NOTE: Before any calls to function w3j are made, it is
25 c necessary to set up a list of binomial coefficients in
26 c common block /bincom/ by calling subroutine bincof. This
27 c is now handled automatically by function wig3j (see below).
28 c
29 c Written by D.J. Schneider, 5/7/90
30 c Updated by D.E.Budil 1/10/93 to fix bug for J2=1
31 c
32 c Coded in this file:
33 c w3j(n1,n2,n3,n4,n5,n6)
34 c bincof()
35 c wig3j(n1,n2,n3,n4,n5,n6)
36 c
37 c Includes:
38 c maxl.inc
39 c
40 c**********************************************************************
41 c
42  function w3j(n1,n2,n3,n4,n5,n6)
43 c
44  use maxl
45 c
46  double precision w3j
47  integer n1,n2,n3,n4,n5,n6
48 c
49  integer j1,j2,j3,m1,m2,m3,jdelta,k
50  double precision phase,parity,x,y,z,ztemp1,ztemp2
51 c
52  double precision wig3j
53 c
54 c######################################################################
55 c
56 c----------------------------------------------------------------------
57 c check triangle conditions, etc.
58 c----------------------------------------------------------------------
59 c
60  if ((n1.lt.0).or.(n2.lt.0).or.(n3.lt.0).or.
61  # ((n2.gt.2).and.(n1+n2+n3+1.gt.2*(mxlval+8)+1))) then
62  write (*,*) '*** quantum nubers too large in w3j ***'
63  stop
64  end if
65 c
66  if ((abs(n4).gt.n1).or.(abs(n5).gt.n2).or.(abs(n6).gt.n3).or.
67  # ((n4+n5+n6).ne.0).or.
68  # ((n1+n2).lt.n3).or.((n1+n3).lt.n2).or.((n2+n3).lt.n1)) then
69  w3j=0.0d0
70  return
71  end if
72 c
73  j1=n1
74  j2=n2
75  j3=n3
76  m1=n4
77  m2=n5
78  m3=n6
79 c
80 c----------------------------------------------------------------------
81 c use wig3j to calculate if j2 > 2
82 c----------------------------------------------------------------------
83 c
84  if (n2.gt.2) then
85  w3j=wig3j(j1,j2,j3,m1,m2,m3)
86  return
87  end if
88 c
89 c----------------------------------------------------------------------
90 c permute variables if necessary to get m2 => 0 and j1 <= j3,
91 c keep track of phases, and store variables
92 c----------------------------------------------------------------------
93 c
94  j1=n1
95  j2=n2
96  j3=n3
97  m1=n4
98  m2=n5
99  m3=n6
100 c
101  phase=1.0d0
102 c
103  if (mod(j1+j2+j3,2).eq.0) then
104  parity=1.d0
105  else
106  parity=-1.d0
107  end if
108 c
109  if (m2.lt.0) then
110  m1=-m1
111  m2=-m2
112  m3=-m3
113  phase=parity
114  else
115  phase=1.d0
116  end if
117 c
118  if(j1.gt.j3) then
119  k=j1
120  j1=j3
121  j3=k
122  k=m1
123  m1=m3
124  m3=k
125  phase=phase*parity
126  end if
127 c
128  if (mod(j1-m3,2).ne.0) then
129  phase=-phase
130  end if
131 c
132 c----------------------------------------------------------------------
133 c calculate wigner 3-j symbols
134 c----------------------------------------------------------------------
135 c
136  jdelta=j3-j1
137  x=dble(j1+j3-1)
138  y=dble(j1-m3)
139  z=dble(j1+m3)
140 c
141  if (j2.eq.0) then
142  w3j=phase/sqrt(dble(2*j1+1))
143  else if (j2.eq.2) then
144  ztemp2=x*(x+1.0d0)*(x+2.0d0)*(x+3.0d0)*(x+4.0d0)
145  if (m2.eq.0) then
146  if (jdelta.eq.0) then
147  ztemp1=2.0d0*dble(3*m3*m3-j1*(j1+1))
148  w3j=phase*ztemp1/sqrt(ztemp2)
149  else if (jdelta.eq.1) then
150  ztemp1=6.0d0*(z+1.0d0)*(y+1.0d0)
151  w3j=-phase*2.0d0*dble(m3)*sqrt(ztemp1/ztemp2)
152  else
153  ztemp1=6.0d0*(z+2.0d0)*(z+1.0d0)
154  # *(y+2.0d0)*(y+1.0d0)
155  w3j=phase*sqrt(ztemp1/ztemp2)
156  end if
157  else if (m2.eq.1) then
158  if (jdelta.eq.0) then
159  ztemp1=6.0d0*(z+1.0d0)*y
160  w3j=phase*dble(2*m3+1)
161  # *sqrt(ztemp1/ztemp2)
162  else if (jdelta.eq.1) then
163  ztemp1=y*(y+1.0d0)
164  w3j=-phase*dble(2*j1+4*m3+4)*
165  # sqrt(ztemp1/ztemp2)
166  else
167  ztemp1=(z+2.0d0)*y*(y+1.0d0)*(y+2.0d0)
168  w3j=phase*2.0d0*sqrt(ztemp1/ztemp2)
169  end if
170  else
171  if (jdelta.eq.0) then
172  ztemp1=6.0d0*(y-1.0d0)*y
173  # *(z+1.0d0)*(z+2.0d0)
174  w3j=phase*sqrt(ztemp1/ztemp2)
175  else if (jdelta.eq.1) then
176  ztemp1=(y-1.0d0)*y*(y+1.0d0)*(z+2.0d0)
177  w3j=-phase*2.d0*sqrt(ztemp1/ztemp2)
178  else
179  ztemp1=(y-1.0d0)*y*(y+1.0d0)*(y+2.0d0)
180  w3j=phase*sqrt(ztemp1/ztemp2)
181  end if
182  end if
183 c
184  else
185  ztemp2=(x+1.0d0)*(x+2.0d0)*(x+3.0d0)
186  if (m2.eq.0) then
187  if (jdelta.eq.0) then
188  w3j=-phase*2.0*dble(m3)/sqrt(ztemp2)
189  else
190  ztemp1=2.0d0*(y+1.0d0)*(z+1.0d0)
191  w3j=-phase*sqrt(ztemp1/ztemp2)
192  end if
193  else
194  if (jdelta.eq.0) then
195  ztemp1=2.0d0*y*(z+1.0d0)
196  w3j=-phase*sqrt(ztemp1/ztemp2)
197  else
198  ztemp1=y*(y+1.0d0)
199  w3j=-phase*sqrt(ztemp1/ztemp2)
200  end if
201  end if
202  end if
203 c
204  return
205  end
206 
207 
208 
209 c----------------------------------------------------------------------
210 c
211 c =========================
212 c BINCOF
213 c =========================
214 c
215 c Generates an array of binomial coefficients for use with function
216 c This must be called before wig3j can calculate a 3-J symbol,
217 c and is automatically called by wig3j if the array has not been
218 c initialized.
219 c
220 c Designed and coded by
221 c Clemens C. J. Roothaan
222 c Department of Chemistry
223 c University of Chicago
224 c 5735 South Ellis Avenue
225 c Chicago, Illinois 60637
226 c December 9, 1988
227 c
228 c modified by DJS 5/7/90, DEB 1/10/93
229 c----------------------------------------------------------------------
230 
231  subroutine bincof()
232 c
233  use maxl
234  use bincom
235 c
236  integer i,j,ij
237  double precision bncf0,temp
238 c
239  ij=1
240  do 20 i=0,nb-1
241  bncfx(i+1)=ij
242  if (i.ne.0) then
243  bncf0=0.0d0
244  do 10 j=1,i
245  temp=bncf(ij-i)
246  bncf(ij)=bncf0+temp
247  bncf0=temp
248  ij=ij+1
249  10 continue
250  end if
251  bncf(ij)=1.0d0
252  ij=ij+1
253  20 continue
254 c
255  return
256  end
257 
258 c----------------------------------------------------------------------
259 c =========================
260 c function WIG3J
261 c =========================
262 c
263 c Calculates arbitrary 3-J symbol using binomial coefficients
264 c for J values up to an limit determined by the machine precision
265 c (approx. 125 for IEEE format real*8 floating-point numbers -DEB)
266 c
267 c This replaces the function of the same name written by G. Moro
268 c for use with EPR spectral simulations.
269 c
270 c Designed and coded by
271 c Clemens C. J. Roothaan
272 c Department of Chemistry
273 c University of Chicago
274 c 5735 South Ellis Avenue
275 c Chicago, Illinois 60637
276 c December 9, 1988
277 c
278 c modified by DJS 5/7/90, DEB 1/10/93
279 c----------------------------------------------------------------------
280 
281  function wig3j(j1,j2,j3,m1,m2,m3)
282 c
283  use maxl
284  use bincom
285 c
286  double precision wig3j
287  integer j1,j2,j3,m1,m2,m3
288  integer i,j,k,l,m,n,p,q,z,zmin,zmax,bp,bnj,bmk
289  double precision sum
290  logical notset
291  data notset / .true. /
292 c
293 c......................................................................
294 c
295 c
296 c
297  if (notset) then
298  call bincof
299  notset = .false.
300  end if
301 c
302  i=j1+m1
303  j=j1-m1
304  k=j2+m2
305  l=j2-m2
306  m=j2+j3-j1
307  n=j3+j1-j2
308  p=j1+j2-j3
309  q=j1+j2+j3 + 1
310 c
311  if(q+1.gt.nb) then
312  write (*,1000) j1,j2,j3,m1,m2,m3
313  return
314  endif
315 c
316  bp=bncfx(p+1)
317  bnj=bncfx(n+1)+n-j
318  bmk=bncfx(m+1)+m-k
319  zmin=max(0,j-n,k-m)
320  zmax=min(p,j,k)
321 c
322  sum=0.0d0
323  do 30 z=zmin,zmax
324  sum=-sum+bncf(bp+z)*bncf(bnj+z)*bncf(bmk+z)
325  30 continue
326 c
327  if (sum.ne.0.0d0) then
328  if (mod(i+l+zmax,2).ne.0) sum=-sum
329  sum=sum*sqrt(bncf(bncfx(q)+m)/
330  & bncf(bncfx(q)+i))/
331  & sqrt(bncf(bncfx(q-i)+k))*
332  & sqrt(bncf(bncfx(q-m)+n)/
333  & bncf(bncfx(q)+j))/
334  & sqrt(bncf(bncfx(q-j)+l))/
335  & sqrt(bncf(bncfx(q+1)+1))
336  endif
337 c
338  wig3j=sum
339  return
340 c
341 c######################################################################
342 c
343  1000 format(' wig3j called with (',6(i3,','),'):'/
344  # ' sum of L values exceeds limit')
345  end
double precision function wig3j(j1, j2, j3, m1, m2, m3)
Definition: w3j.f90:282
Definition: maxl.f90:22
double precision, dimension(nbncf), save bncf
Definition: bincom.f90:21
double precision function w3j(n1, n2, n3, n4, n5, n6)
Definition: w3j.f90:43
subroutine bincof()
Definition: w3j.f90:232
integer, dimension(nb), save bncfx
Definition: bincom.f90:20
integer, parameter mxlval
Definition: maxl.f90:25
integer, parameter nb
Definition: bincom.f90:17