1 subroutine da_effang_tl(ifreq,theta,gx2,sigma,effangv,effangh, &
2 TGL_gx2,TGL_sigma,TGL_effangv,TGL_effangh )
4 !--------------------------------------------------------------------------
5 ! Purpose : Calculate the effective zenith angle of reflected microwave
6 ! radiation at SSM/I frequencies for vertical and horizontal polarization
8 ! Input :: TGL_gx2, TGL_sigma
9 ! Output :: TGL_effangv,TGL_effangh,effangv,effangh
10 !--------------------------------------------------------------------------
14 integer, intent(in) :: ifreq
15 real, intent(in) :: theta,gx2,sigma,TGL_gx2, TGL_sigma
16 real, intent(out) :: TGL_effangv,TGL_effangh,effangv,effangh
18 real c19v,c19h,c22v,c22h,c37v,c37h,c85v,c85h
19 real s19v(6),s19h(6),s22v(6),s22h(6), &
20 s37v(6),s37h(6),s85v(6),s85h(6)
22 real :: alnsig,gg,ggg,xd,xx
23 real :: z1,z2,z3,z4,z5,z6,alnthv
24 real :: y,dth,angh,angv,alnthh
25 real :: TGL_alnsig,TGL_gg,TGL_ggg,TGL_xd
26 real :: TGL_z1,TGL_z2,TGL_z3,TGL_z4,TGL_z5,TGL_z6,TGL_alnthv
27 real :: TGL_y,TGL_dth,TGL_angh,TGL_angv,TGL_xx,TGL_alnthh
29 data c19v,c19h,c22v,c22h,c37v,c37h,c85v,c85h &
30 /-.5108,.5306,-.5108,.5306,-.6931,.1823,-.9163,.3000/
31 data s19v /.225E+2,.698E+2,-.238E+2,-.648E+1,.402E+0,.262E+1/
32 data s19h /.743E+1,.507E+2,-.206E+2,-.191E+1,.648E-1,.291E+1/
33 data s22v /.249E+2,.701E+2,-.240E+2,-.714E+1,.405E+0,.256E+1/
34 data s22h /.498E+1,.442E+2,-.190E+2,-.129E+1,.803E-2,.277E+1/
35 data s37v /.215E+2,.573E+2,-.211E+2,-.670E+1,.443E+0,.253E+1/
36 data s37h /.869E+1,.571E+2,-.257E+2,-.302E+1,.237E+0,.386E+1/
37 data s85v /.116E+2,.263E+2,-.101E+2,-.358E+1,.270E+0,.175E+1/
38 data s85h /.736E+1,.568E+2,-.254E+2,-.248E+1,.196E+0,.387E+1/
40 if (trace_use) call da_trace_entry("da_effang_tl")
42 if (gx2 .le. 0.0 .or. sigma .le. 0.0) then
50 if (abs(sigma) .gt. 0.0) then
51 TGL_alnsig = TGL_sigma/sigma
56 TGL_gg = 2.0*gx2*TGL_gx2
58 TGL_ggg = TGL_gg*gx2 + gg*TGL_gx2
60 if (ifreq .eq. 1) then
65 TGL_xx = 2.0*xd*TGL_xd
67 TGL_z1 = TGL_xx*ggg + xx*TGL_ggg
69 TGL_z2 = TGL_xd*ggg + xd*TGL_ggg
71 TGL_z3 = TGL_xd*gg + xd*TGL_gg
73 TGL_z4 = TGL_xx*gg + xx*TGL_gg
75 TGL_z5 = TGL_xx*gx2 + xx*TGL_gx2
77 TGL_z6 = TGL_xd*gx2 + xd*TGL_gx2
78 alnthv = s19v(1)*z1 + s19v(2)*z2 + s19v(3)*z3 + &
79 s19v(4)*z4 + s19v(5)*z5 + s19v(6)*z6
80 TGL_alnthv = s19v(1)*TGL_z1 + s19v(2)*TGL_z2 + s19v(3)*TGL_z3 + &
81 s19v(4)*TGL_z4 + s19v(5)*TGL_z5 + s19v(6)*TGL_z6
82 alnthv = alnthv + 3.611
87 TGL_xx = 2.0*xd*TGL_xd
89 TGL_z1 = TGL_xx*ggg + xx*TGL_ggg
91 TGL_z2 = TGL_xd*ggg + xd*TGL_ggg
93 TGL_z3 = TGL_xd*gg + xd*TGL_gg
95 TGL_z4 = TGL_xx*gg + xx*TGL_gg
97 TGL_z5 = TGL_xx*gx2 + xx*TGL_gx2
99 TGL_z6 = TGL_xd*gx2 + xd*TGL_gx2
101 alnthh = s19h(1)*z1 + s19h(2)*z2 + s19h(3)*z3 + &
102 s19h(4)*z4 + s19h(5)*z5 + s19h(6)*z6
103 TGL_alnthh = s19h(1)*TGL_z1 + s19h(2)*TGL_z2 + s19h(3)*TGL_z3 + &
104 s19h(4)*TGL_z4 + s19h(5)*TGL_z5 + s19h(6)*TGL_z6
106 alnthh = alnthh + 3.611
108 else if (ifreq .eq. 2) then
112 TGL_xx = 2.0*xd*TGL_xd
114 TGL_z1 = TGL_xx*ggg + xx*TGL_ggg
116 TGL_z2 = TGL_xd*ggg + xd*TGL_ggg
118 TGL_z3 = TGL_xd*gg + xd*TGL_gg
120 TGL_z4 = TGL_xx*gg + xx*TGL_gg
122 TGL_z5 = TGL_xx*gx2 + xx*TGL_gx2
124 TGL_z6 = TGL_xd*gx2 + xd*TGL_gx2
125 alnthv = s22v(1)*z1 + s22v(2)*z2 + s22v(3)*z3 + &
126 s22v(4)*z4 + s22v(5)*z5 + s22v(6)*z6
127 TGL_alnthv = s22v(1)*TGL_z1 + s22v(2)*TGL_z2 + s22v(3)*TGL_z3 + &
128 s22v(4)*TGL_z4 + s22v(5)*TGL_z5 + s22v(6)*TGL_z6
129 alnthv = alnthv + 3.611
130 ! TGL_alnthv = TGL_alnthv
135 TGL_xx = 2.0*xd*TGL_xd
137 TGL_z1 = TGL_xx*ggg + xx*TGL_ggg
139 TGL_z2 = TGL_xd*ggg + xd*TGL_ggg
141 TGL_z3 = TGL_xd*gg + xd*TGL_gg
143 TGL_z4 = TGL_xx*gg + xx*TGL_gg
145 TGL_z5 = TGL_xx*gx2 + xx*TGL_gx2
147 TGL_z6 = TGL_xd*gx2 + xd*TGL_gx2
148 alnthh = s22h(1)*z1 + s22h(2)*z2 + s22h(3)*z3 + &
149 s22h(4)*z4 + s22h(5)*z5 + s22h(6)*z6
150 TGL_alnthh = s22h(1)*TGL_z1 + s22h(2)*TGL_z2 + s22h(3)*TGL_z3 + &
151 s22h(4)*TGL_z4 + s22h(5)*TGL_z5 + s22h(6)*TGL_z6
152 alnthh = alnthh + 3.611
153 ! TGL_alnthh = TGL_alnthh
154 else if (ifreq .eq. 3) then
158 TGL_xx = 2.0*xd*TGL_xd
160 TGL_z1 = TGL_xx*ggg + xx*TGL_ggg
162 TGL_z2 = TGL_xd*ggg + xd*TGL_ggg
164 TGL_z3 = TGL_xd*gg + xd*TGL_gg
166 TGL_z4 = TGL_xx*gg + xx*TGL_gg
168 TGL_z5 = TGL_xx*gx2 + xx*TGL_gx2
170 TGL_z6 = TGL_xd*gx2 + xd*TGL_gx2
171 alnthv = s37v(1)*z1 + s37v(2)*z2 + s37v(3)*z3 + &
172 s37v(4)*z4 + s37v(5)*z5 + s37v(6)*z6
173 TGL_alnthv = s37v(1)*TGL_z1 + s37v(2)*TGL_z2 + s37v(3)*TGL_z3 + &
174 s37v(4)*TGL_z4 + s37v(5)*TGL_z5 + s37v(6)*TGL_z6
175 alnthv = alnthv + 3.611
176 ! TGL_alnthv = TGL_alnthv
181 TGL_xx = 2.0*xd*TGL_xd
183 TGL_z1 = TGL_xx*ggg + xx*TGL_ggg
185 TGL_z2 = TGL_xd*ggg + xd*TGL_ggg
187 TGL_z3 = TGL_xd*gg + xd*TGL_gg
189 TGL_z4 = TGL_xx*gg + xx*TGL_gg
191 TGL_z5 = TGL_xx*gx2 + xx*TGL_gx2
193 TGL_z6 = TGL_xd*gx2 + xd*TGL_gx2
194 alnthh = s37h(1)*z1 + s37h(2)*z2 + s37h(3)*z3 + &
195 s37h(4)*z4 + s37h(5)*z5 + s37h(6)*z6
196 TGL_alnthh = s37h(1)*TGL_z1 + s37h(2)*TGL_z2 + s37h(3)*TGL_z3 + &
197 s37h(4)*TGL_z4 + s37h(5)*TGL_z5 + s37h(6)*TGL_z6
198 alnthh = alnthh + 3.611
199 ! TGL_alnthh = TGL_alnthh
200 else if (ifreq .eq. 4) then
204 TGL_xx = 2.0*xd*TGL_xd
206 TGL_z1 = TGL_xx*ggg + xx*TGL_ggg
208 TGL_z2 = TGL_xd*ggg + xd*TGL_ggg
210 TGL_z3 = TGL_xd*gg + xd*TGL_gg
212 TGL_z4 = TGL_xx*gg + xx*TGL_gg
214 TGL_z5 = TGL_xx*gx2 + xx*TGL_gx2
216 TGL_z6 = TGL_xd*gx2 + xd*TGL_gx2
217 alnthv = s85v(1)*z1 + s85v(2)*z2 + s85v(3)*z3 + &
218 s85v(4)*z4 + s85v(5)*z5 + s85v(6)*z6
219 TGL_alnthv = s85v(1)*TGL_z1 + s85v(2)*TGL_z2 + s85v(3)*TGL_z3 + &
220 s85v(4)*TGL_z4 + s85v(5)*TGL_z5 + s85v(6)*TGL_z6
221 alnthv = alnthv + 3.611
222 ! TGL_alnthv = TGL_alnthv
227 TGL_xx = 2.0*xd*TGL_xd
229 TGL_z1 = TGL_xx*ggg + xx*TGL_ggg
231 TGL_z2 = TGL_xd*ggg + xd*TGL_ggg
233 TGL_z3 = TGL_xd*gg + xd*TGL_gg
235 TGL_z4 = TGL_xx*gg + xx*TGL_gg
237 TGL_z5 = TGL_xx*gx2 + xx*TGL_gx2
239 TGL_z6 = TGL_xd*gx2 + xd*TGL_gx2
240 alnthh = s85h(1)*z1 + s85h(2)*z2 + s85h(3)*z3 + &
241 s85h(4)*z4 + s85h(5)*z5 + s85h(6)*z6
242 TGL_alnthh = s85h(1)*TGL_z1 + s85h(2)*TGL_z2 + s85h(3)*TGL_z3 + &
243 s85h(4)*TGL_z4 + s85h(5)*TGL_z5 + s85h(6)*TGL_z6
244 alnthh = alnthh + 3.611
245 ! TGL_alnthh = TGL_alnthh
247 angv = 90.0 - exp(alnthv)
248 TGL_angv = - TGL_alnthv*exp(alnthv)
249 angh = 90.0 - exp(alnthh)
250 TGL_angh = - TGL_alnthh*exp(alnthh)
252 TGL_y = - 28.0*TGL_gx2
257 dth = (theta - 53.0)*y
258 TGL_dth = (theta - 53.0)*TGL_y
260 TGL_effangv = TGL_angv + TGL_dth
262 TGL_effangh = TGL_angh + TGL_dth
264 if (trace_use) call da_trace_exit("da_effang_tl")
266 end subroutine da_effang_tl