Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_ssmi / effang.inc
blobb610b5ebe76ef7bc978121fc2e6895f9d6286a8a
1       subroutine effang(ifreq,theta,gx2,sigma,effangv,effangh)
3 ! Calculated the effective zenith angle of reflected microwave radiation
4 ! at SSM/I frequencies for vertical and horizontal polarization
6       real theta, gx2,sigma,effangv,effangh
7       integer ifreq
8       real gg,ggg,xx,xd,z1,z2,z3,z4,z5,z6,alnthv
9       real alnthh,angv,angh,y,alnsig,dth
11       real c19v,c19h,c22v,c22h,c37v,c37h,c85v,c85h
12       real s19v(6),s19h(6),s22v(6),s22h(6), &
13            s37v(6),s37h(6),s85v(6),s85h(6)
15       data c19v,c19h,c22v,c22h,c37v,c37h,c85v,c85h &
16         /-.5108,.5306,-.5108,.5306,-.6931,.1823,-.9163,.3000/
17       data s19v /.225E+2,.698E+2,-.238E+2,-.648E+1,.402E+0,.262E+1/
18       data s19h /.743E+1,.507E+2,-.206E+2,-.191E+1,.648E-1,.291E+1/
19       data s22v /.249E+2,.701E+2,-.240E+2,-.714E+1,.405E+0,.256E+1/
20       data s22h /.498E+1,.442E+2,-.190E+2,-.129E+1,.803E-2,.277E+1/
21       data s37v /.215E+2,.573E+2,-.211E+2,-.670E+1,.443E+0,.253E+1/
22       data s37h /.869E+1,.571E+2,-.257E+2,-.302E+1,.237E+0,.386E+1/
23       data s85v /.116E+2,.263E+2,-.101E+2,-.358E+1,.270E+0,.175E+1/
24       data s85h /.736E+1,.568E+2,-.254E+2,-.248E+1,.196E+0,.387E+1/
26       if (gx2 .le. 0.0 .or. sigma .le. 0.0) then
27         effangv = theta
28         effangh = theta
29         return
30       end if
31       alnsig = alog(sigma)
32       gg = gx2*gx2
33       ggg = gg*gx2
34       if (ifreq .eq. 1) then 
35          xd = alnsig - c19v
36          xx = xd*xd
37          z1 =  xx*ggg
38          z2 =  xd*ggg
39          z3 =  xd*gg
40          z4 =  xx*gg
41          z5 =  xx*gx2
42          z6 =  xd*gx2
43          alnthv = s19v(1)*z1 + s19v(2)*z2 +s19v(3)*z3 + &
44            s19v(4)*z4 + s19v(5)*z5 + s19v(6)*z6
45          alnthv = alnthv + 3.611
47          xd = alnsig - c19h
48          xx = xd*xd
49          z1 =  xx*ggg
50          z2 =  xd*ggg
51          z3 =  xd*gg
52          z4 =  xx*gg
53          z5 =  xx*gx2
54          z6 =  xd*gx2
55          alnthh = s19h(1)*z1 + s19h(2)*z2 +s19h(3)*z3 + &
56            s19h(4)*z4 + s19h(5)*z5 + s19h(6)*z6
57          alnthh = alnthh + 3.611
58       else if (ifreq .eq. 2) then 
59          xd = alnsig - c22v
60          xx = xd*xd
61          z1 =  xx*ggg
62          z2 =  xd*ggg
63          z3 =  xd*gg
64          z4 =  xx*gg
65          z5 =  xx*gx2
66          z6 =  xd*gx2
67          alnthv = s22v(1)*z1 + s22v(2)*z2 +s22v(3)*z3 + &
68            s22v(4)*z4 + s22v(5)*z5 + s22v(6)*z6
69          alnthv = alnthv + 3.611
71          xd = alnsig - c22h
72          xx = xd*xd
73          z1 =  xx*ggg
74          z2 =  xd*ggg
75          z3 =  xd*gg
76          z4 =  xx*gg
77          z5 =  xx*gx2
78          z6 =  xd*gx2
79          alnthh = s22h(1)*z1 + s22h(2)*z2 +s22h(3)*z3 + &
80            s22h(4)*z4 + s22h(5)*z5 + s22h(6)*z6
81          alnthh = alnthh + 3.611
82       else if (ifreq .eq. 3) then 
83          xd = alnsig - c37v
84          xx = xd*xd
85          z1 =  xx*ggg
86          z2 =  xd*ggg
87          z3 =  xd*gg
88          z4 =  xx*gg
89          z5 =  xx*gx2
90          z6 =  xd*gx2
91          alnthv = s37v(1)*z1 + s37v(2)*z2 +s37v(3)*z3 + &
92            s37v(4)*z4 + s37v(5)*z5 + s37v(6)*z6
93          alnthv = alnthv + 3.611
95          xd = alnsig - c37h
96          xx = xd*xd
97          z1 =  xx*ggg
98          z2 =  xd*ggg
99          z3 =  xd*gg
100          z4 =  xx*gg
101          z5 =  xx*gx2
102          z6 =  xd*gx2
103          alnthh = s37h(1)*z1 + s37h(2)*z2 +s37h(3)*z3 + &
104            s37h(4)*z4 + s37h(5)*z5 + s37h(6)*z6
105          alnthh = alnthh + 3.611
106       else if (ifreq .eq. 4) then 
107          xd = alnsig - c85v
108          xx = xd*xd
109          z1 =  xx*ggg
110          z2 =  xd*ggg
111          z3 =  xd*gg
112          z4 =  xx*gg
113          z5 =  xx*gx2
114          z6 =  xd*gx2
115          alnthv = s85v(1)*z1 + s85v(2)*z2 +s85v(3)*z3 + &
116            s85v(4)*z4 + s85v(5)*z5 + s85v(6)*z6
117          alnthv = alnthv + 3.611
119          xd = alnsig - c85h
120          xx = xd*xd
121          z1 =  xx*ggg
122          z2 =  xd*ggg
123          z3 =  xd*gg
124          z4 =  xx*gg
125          z5 =  xx*gx2
126          z6 =  xd*gx2
127          alnthh = s85h(1)*z1 + s85h(2)*z2 +s85h(3)*z3 + &
128            s85h(4)*z4 + s85h(5)*z5 + s85h(6)*z6
129          alnthh = alnthh + 3.611
130       end if
131       angv = 90.0 - exp(alnthv)
132       angh = 90.0 - exp(alnthh)
133       y  = 1.0 - 28.0*gx2
134       if (y .lt. 0.0) y = 0.0
135       dth = (theta - 53.0)*y
136       effangv = angv + dth
137       effangh = angh + dth
139       end subroutine effang