OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
bhmie.py
Go to the documentation of this file.
1 from numpy import *
2 
3 def bhmie(x,refrel,nang):
4 
5 # This file is converted from mie.m, see http://atol.ucsd.edu/scatlib/index.htm
6 # Bohren and Huffman originally published the code in their book on light scattering
7 
8 # Calculation based on Mie scattering theory
9 # input:
10 # x - size parameter = k*radius = 2pi/lambda * radius
11 # (lambda is the wavelength in the medium around the scatterers)
12 # refrel - refraction index (n in complex form for example: 1.5+0.02*i;
13 # nang - number of angles for S1 and S2 function in range from 0 to pi/2
14 # output:
15 # S1, S2 - funtion which correspond to the (complex) phase functions
16 # Qext - extinction efficiency
17 # Qsca - scattering efficiency
18 # Qback - backscatter efficiency
19 # gsca - asymmetry parameter
20 
21 
22  nmxx=150000
23 
24  # Require NANG>1 in order to calculate scattering intensities
25  if (nang < 2):
26  nang = 2
27 
28  s1_1=zeros(nang,dtype=complex128)
29  s1_2=zeros(nang,dtype=complex128)
30  s2_1=zeros(nang,dtype=complex128)
31  s2_2=zeros(nang,dtype=complex128)
32  pi=zeros(nang,dtype=complex128)
33  tau=zeros(nang,dtype=complex128)
34 
35  if (nang > 1000):
36  print ('error: nang > mxnang=1000 in bhmie')
37  return
38 
39  pii = 4.*arctan(1.)
40  dx = x
41 
42  drefrl = refrel
43  y = x*drefrl
44  ymod = abs(y)
45 
46 
47  # Series expansion terminated after NSTOP terms
48  # Logarithmic derivatives calculated from NMX on down
49 
50  xstop = x + 4.*x**0.3333 + 2.0
51  #xstop = x + 4.*x**0.3333 + 10.0
52  nmx = max(xstop,ymod) + 15.0
53  nmx=fix(nmx)
54 
55  # BTD experiment 91/1/15: add one more term to series and compare resu<s
56  # NMX=AMAX1(XSTOP,YMOD)+16
57  # test: compute 7001 wavelen>hs between .0001 and 1000 micron
58  # for a=1.0micron SiC grain. When NMX increased by 1, only a single
59  # computed number changed (out of 4*7001) and it only changed by 1/8387
60  # conclusion: we are indeed retaining enough terms in series!
61 
62  nstop = int(xstop)
63 
64  if (nmx > nmxx):
65  print ( "error: nmx > nmxx=%f for |m|x=%f" % ( nmxx, ymod) )
66  return
67 
68  dang = .5*pii/ (nang-1)
69 
70 
71  amu=arange(0.0,nang,1)
72  amu=cos(amu*dang)
73 
74  pi0=zeros(nang,dtype=complex128)
75  pi1=ones(nang,dtype=complex128)
76 
77  # Logarithmic derivative D(J) calculated by downward recurrence
78  # beginning with initial value (0.,0.) at J=NMX
79 
80  nn = int(nmx)-1
81  d=zeros(nn+1,dtype=complex128)
82  for n in range(0,nn):
83  en = nmx - n
84  d[nn-n-1] = (en/y) - (1./ (d[nn-n]+en/y))
85 
86 
87  #*** Riccati-Bessel functions with real argument X
88  # calculated by upward recurrence
89 
90  psi0 = cos(dx)
91  psi1 = sin(dx)
92  chi0 = -sin(dx)
93  chi1 = cos(dx)
94  xi1 = psi1-chi1*1j
95  qsca = 0.
96  gsca = 0.
97  p = -1
98 
99  for n in range(0,nstop):
100  en = n+1.0
101  fn = (2.*en+1.)/(en* (en+1.))
102 
103  # for given N, PSI = psi_n CHI = chi_n
104  # PSI1 = psi_{n-1} CHI1 = chi_{n-1}
105  # PSI0 = psi_{n-2} CHI0 = chi_{n-2}
106  # Calculate psi_n and chi_n
107  psi = (2.*en-1.)*psi1/dx - psi0
108  chi = (2.*en-1.)*chi1/dx - chi0
109  xi = psi-chi*1j
110 
111  #*** Store previous values of AN and BN for use
112  # in computation of g=<cos(theta)>
113  if (n > 0):
114  an1 = an
115  bn1 = bn
116 
117  #*** Compute AN and BN:
118  an = (d[n]/drefrl+en/dx)*psi - psi1
119  an = an/ ((d[n]/drefrl+en/dx)*xi-xi1)
120  bn = (drefrl*d[n]+en/dx)*psi - psi1
121  bn = bn/ ((drefrl*d[n]+en/dx)*xi-xi1)
122 
123  #*** Augment sums for Qsca and g=<cos(theta)>
124  qsca += (2.*en+1.)* (abs(an)**2+abs(bn)**2)
125  gsca += ((2.*en+1.)/ (en* (en+1.)))*( real(an)* real(bn)+imag(an)*imag(bn))
126 
127  if (n > 0):
128  gsca += ((en-1.)* (en+1.)/en)*( real(an1)* real(an)+imag(an1)*imag(an)+real(bn1)* real(bn)+imag(bn1)*imag(bn))
129 
130 
131  #*** Now calculate scattering intensity pattern
132  # First do angles from 0 to 90
133  pi=0+pi1 # 0+pi1 because we want a hard copy of the values
134  tau=en*amu*pi-(en+1.)*pi0
135  s1_1 += fn* (an*pi+bn*tau)
136  s2_1 += fn* (an*tau+bn*pi)
137 
138  #*** Now do angles greater than 90 using PI and TAU from
139  # angles less than 90.
140  # P=1 for N=1,3,...% P=-1 for N=2,4,...
141  # remember that we have to reverse the order of the elements
142  # of the second part of s1 and s2 after the calculation
143  p = -p
144  s1_2+= fn*p* (an*pi-bn*tau)
145  s2_2+= fn*p* (bn*pi-an*tau)
146 
147  psi0 = psi1
148  psi1 = psi
149  chi0 = chi1
150  chi1 = chi
151  xi1 = psi1-chi1*1j
152 
153  #*** Compute pi_n for next value of n
154  # For each angle J, compute pi_n+1
155  # from PI = pi_n , PI0 = pi_n-1
156  pi1 = ((2.*en+1.)*amu*pi- (en+1.)*pi0)/ en
157  pi0 = 0+pi # 0+pi because we want a hard copy of the values
158 
159  #*** Have summed sufficient terms.
160  # Now compute QSCA,QEXT,QBACK,and GSCA
161 
162  # we have to reverse the order of the elements of the second part of s1 and s2
163  s1=concatenate((s1_1,s1_2[-2::-1]))
164  s2=concatenate((s2_1,s2_2[-2::-1]))
165  gsca = 2.*gsca/qsca
166  qsca = (2./ (dx*dx))*qsca
167  qext = (4./ (dx*dx))* real(s1[0])
168 
169  # more common definition of the backscattering efficiency,
170  # so that the backscattering cross section really
171  # has dimension of length squared
172  qback = 4*(abs(s1[2*nang-2])/dx)**2
173  #qback = ((abs(s1[2*nang-2])/dx)**2 )/pii #old form
174 
175  return s1,s2,qext,qsca,qback,gsca
176 
def bhmie(x, refrel, nang)
Definition: bhmie.py:3
#define real
Definition: DbAlgOcean.cpp:26
#define abs(a)
Definition: misc.h:90