da_effang_adj.inc

References to this file elsewhere.
1 subroutine da_effang_adj(ifreq,theta,gx2,sigma,effangv,effangh,     &
2                       ADJ_gx2,ADJ_sigma,ADJ_effangv,ADJ_effangh  )
3 
4    implicit none
5 
6    !-----------------------------------------------------------------
7    ! Purpose: Calculate the effective zenith angle of reflected microwave
8    ! radiation at SSM/I frequencies for vertical and horizontal polarization
9    !
10    ! Output: ADJ_gx2, ADJ_sigma
11    ! Input: ADJ_effangv,ADJ_effangh,effangv,effangh
12    !-----------------------------------------------------------------
13  
14    integer, intent(in   ) :: ifreq
15    real   , intent(in   ) :: theta,gx2,sigma
16    real   , intent(inout) :: ADJ_gx2,ADJ_sigma
17    real   , intent(inout) :: ADJ_effangv,ADJ_effangh
18    real   , intent(out)   :: effangv,effangh
19 
20    real c19v,c19h,c22v,c22h,c37v,c37h,c85v,c85h
21    real s19v(6),s19h(6),s22v(6),s22h(6), &
22         s37v(6),s37h(6),s85v(6),s85h(6)
23 
24    real :: alnsig,gg,ggg,xd,xx,xd_save,xx_save
25    real :: z1,z2,z3,z4,z5,z6
26    real :: y,dth,angh,angv,alnthv,alnthh,alnthv_save,alnthh_save
27    real :: ADJ_alnsig,ADJ_gg,ADJ_ggg,ADJ_xd
28    real :: ADJ_z1,ADJ_z2,ADJ_z3,ADJ_z4,ADJ_z5,ADJ_z6,ADJ_alnthv
29    real :: ADJ_y,ADJ_dth,ADJ_angh,ADJ_angv,ADJ_xx,ADJ_alnthh
30 
31    data c19v,c19h,c22v,c22h,c37v,c37h,c85v,c85h &
32      /-.5108,.5306,-.5108,.5306,-.6931,.1823,-.9163,.3000/
33    data s19v /.225E+2,.698E+2,-.238E+2,-.648E+1,.402E+0,.262E+1/
34    data s19h /.743E+1,.507E+2,-.206E+2,-.191E+1,.648E-1,.291E+1/
35    data s22v /.249E+2,.701E+2,-.240E+2,-.714E+1,.405E+0,.256E+1/
36    data s22h /.498E+1,.442E+2,-.190E+2,-.129E+1,.803E-2,.277E+1/
37    data s37v /.215E+2,.573E+2,-.211E+2,-.670E+1,.443E+0,.253E+1/
38    data s37h /.869E+1,.571E+2,-.257E+2,-.302E+1,.237E+0,.386E+1/
39    data s85v /.116E+2,.263E+2,-.101E+2,-.358E+1,.270E+0,.175E+1/
40    data s85h /.736E+1,.568E+2,-.254E+2,-.248E+1,.196E+0,.387E+1/
41 
42    alnsig=0.;gg=0.;ggg=0.;xd=0.;xx=0.;xd_save=0.;xx_save=0.
43    z1=0.;z2=0.;z3=0.;z4=0.;z5=0.;z6=0.;y=0.;dth=0.;angh=0.
44    angv=0.;alnthv=0.;alnthh=0.;alnthv_save=0.;alnthh_save=0.
45    ADJ_alnsig=0.;ADJ_gg=0.;ADJ_ggg=0.;ADJ_xd=0.
46    ADJ_z1=0.;ADJ_z2=0.;ADJ_z3=0.;ADJ_z4=0.;ADJ_z5=0.
47    ADJ_z6=0.;ADJ_alnthv=0.;ADJ_y=0.;ADJ_dth=0.;ADJ_angh=0.
48    ADJ_angv=0.;ADJ_xx=0.;ADJ_alnthh=0.
49 
50 
51    if (gx2 .le. 0.0 .or. sigma .le. 0.0) then
52       effangv = theta
53       effangh = theta
54       return
55    end if
56    alnsig = alog(sigma)
57    gg  = gx2*gx2
58    ggg = gg*gx2
59 
60    if (ifreq .eq. 1) then 
61       xd =      alnsig - c19v
62       xx =  xd*xd
63       z1 =  xx*ggg
64       z2 =  xd*ggg
65       z3 =  xd*gg
66       z4 =  xx*gg
67       z5 =  xx*gx2
68       z6 =  xd*gx2
69       alnthv =  s19v(1)*z1 + s19v(2)*z2 + s19v(3)*z3 + &
70                 s19v(4)*z4 + s19v(5)*z5 + s19v(6)*z6
71       alnthv_save = alnthv
72       alnthv =  alnthv + 3.611
73 
74       xd_save = xd
75       xx_save = xx
76 
77       xd =  alnsig - c19h
78       xx =  xd*xd
79       z1 =  xx*ggg
80       z2 =  xd*ggg
81       z3 =  xd*gg
82       z4 =  xx*gg
83       z5 =  xx*gx2
84       z6 =  xd*gx2
85 
86       alnthh =  s19h(1)*z1 + s19h(2)*z2 + s19h(3)*z3 + &
87                 s19h(4)*z4 + s19h(5)*z5 + s19h(6)*z6
88       alnthh_save = alnthh
89       alnthh =  alnthh + 3.611
90 
91    else if (ifreq .eq. 2) then 
92       xd =      alnsig - c22v
93       xx =  xd*xd
94       z1 =  xx*ggg
95       z2 =  xd*ggg
96       z3 =  xd*gg
97       z4 =  xx*gg
98       z5 =  xx*gx2
99       z6 =  xd*gx2
100       alnthv =  s22v(1)*z1 + s22v(2)*z2 + s22v(3)*z3 + &
101                 s22v(4)*z4 + s22v(5)*z5 + s22v(6)*z6
102 
103       alnthv_save = alnthv
104       alnthv =  alnthv + 3.611
105  
106       xd_save = xd
107       xx_save = xx
108 
109       xd =      alnsig - c22h
110       xx =  xd*xd
111       z1 =  xx*ggg
112       z2 =  xd*ggg
113       z3 =  xd*gg
114       z4 =  xx*gg
115       z5 =  xx*gx2
116       z6 =  xd*gx2
117       alnthh =  s22h(1)*z1 + s22h(2)*z2 + s22h(3)*z3 + &
118                    s22h(4)*z4 + s22h(5)*z5 + s22h(6)*z6
119 
120       alnthh_save = alnthh
121       alnthh =  alnthh + 3.611
122 
123    else if (ifreq .eq. 3) then 
124       xd =      alnsig - c37v
125       xx =  xd*xd
126       z1 =  xx*ggg
127       z2 =  xd*ggg
128       z3 =  xd*gg
129       z4 =  xx*gg
130       z5 =  xx*gx2
131       z6 =  xd*gx2
132       alnthv =  s37v(1)*z1 + s37v(2)*z2 + s37v(3)*z3 + &
133                 s37v(4)*z4 + s37v(5)*z5 + s37v(6)*z6
134 
135       alnthv_save = alnthv
136       alnthv =  alnthv + 3.611
137  
138       xd_save = xd
139       xx_save = xx
140 
141       xd =      alnsig - c37h
142       xx =  xd*xd
143       z1 =  xx*ggg
144       z2 =  xd*ggg
145       z3 =  xd*gg
146       z4 =  xx*gg
147       z5 =  xx*gx2
148       z6 =  xd*gx2
149       alnthh =  s37h(1)*z1 + s37h(2)*z2 + s37h(3)*z3 + &
150                 s37h(4)*z4 + s37h(5)*z5 + s37h(6)*z6
151  
152       alnthh_save = alnthh
153       alnthh =  alnthh + 3.611
154 
155    else if (ifreq .eq. 4) then 
156       xd =      alnsig - c85v
157       xx =  xd*xd
158       z1 =  xx*ggg
159       z2 =  xd*ggg
160       z3 =  xd*gg
161       z4 =  xx*gg
162       z5 =  xx*gx2
163       z6 =  xd*gx2
164       alnthv =  s85v(1)*z1 + s85v(2)*z2 + s85v(3)*z3 + &
165                 s85v(4)*z4 + s85v(5)*z5 + s85v(6)*z6
166 
167       alnthv_save = alnthv
168       alnthv =  alnthv + 3.611
169 
170       xd_save = xd
171       xx_save = xx
172 
173       xd =      alnsig - c85h
174       xx =  xd*xd
175       z1 =  xx*ggg
176       z2 =  xd*ggg
177       z3 =  xd*gg
178       z4 =  xx*gg
179       z5 =  xx*gx2
180       z6 =  xd*gx2
181       alnthh =  s85h(1)*z1 + s85h(2)*z2 + s85h(3)*z3 + &
182                 s85h(4)*z4 + s85h(5)*z5 + s85h(6)*z6
183 
184       alnthh_save = alnthh
185       alnthh =  alnthh + 3.611
186    end if
187 
188    angv =   90.0 - exp(alnthv)
189    angh =   90.0 - exp(alnthh)
190        y    =   1.0 - 28.0*gx2
191    if (y .lt. 0.0) then
192        y = 0.0
193    end if
194 
195    dth     = (theta - 53.0)*y
196    effangv = angv + dth
197    effangh = angh + dth
198 
199    ! start
200 
201    if (gx2 .le. 0.0 .or. sigma .le. 0.0) then
202       ADJ_effangv = 0.0
203       ADJ_effangh = 0.0
204       return
205    end if
206 
207    ADJ_angh  = ADJ_effangh
208    ADJ_dth   = ADJ_effangh
209 
210    ADJ_angv  = ADJ_effangv
211    ADJ_dth   = ADJ_effangv  + ADJ_dth
212 
213    ADJ_y     = (theta - 53.0)*ADJ_dth  
214 
215    if (y .lt. 0.0) then
216       ADJ_y = 0.0
217    end if
218 
219    ADJ_gx2  = - 28.0*ADJ_y + ADJ_gx2
220 
221    ADJ_alnthh = - ADJ_angh*exp(alnthh)
222    ADJ_alnthv = - ADJ_angv*exp(alnthv)
223 
224    if (ifreq .eq. 1) then 
225 
226       alnthh = alnthh_save
227 
228       ADJ_z1  = s19h(1)*ADJ_alnthh
229       ADJ_z2  = s19h(2)*ADJ_alnthh
230       ADJ_z3  = s19h(3)*ADJ_alnthh
231       ADJ_z4  = s19h(4)*ADJ_alnthh
232       ADJ_z5  = s19h(5)*ADJ_alnthh
233       ADJ_z6  = s19h(6)*ADJ_alnthh
234 
235       ADJ_xd  = ADJ_z6*gx2
236       ADJ_gx2 = xd*ADJ_z6  + ADJ_gx2
237       ADJ_xx  = ADJ_z5*gx2
238       ADJ_gx2 = xx*ADJ_z5  + ADJ_gx2
239       ADJ_xx  = ADJ_z4*gg  + ADJ_xx
240       ADJ_gg  = xx*ADJ_z4
241       ADJ_xd  = ADJ_z3*gg  + ADJ_xd
242       ADJ_gg  = xd*ADJ_z3  + ADJ_gg
243       ADJ_xd  = ADJ_z2*ggg + ADJ_xd
244       ADJ_ggg = xd*ADJ_z2
245       ADJ_xx  = ADJ_z1*ggg + ADJ_xx
246       ADJ_ggg = xx*ADJ_z1  + ADJ_ggg
247 
248       ADJ_xd  = 2.*xd*ADJ_xx + ADJ_xd
249       ADJ_alnsig = ADJ_xd
250 
251       ! vertical
252       xd =  xd_save
253       xx =  xx_save
254 
255       alnthv = alnthv_save
256 
257       ADJ_z1   = s19v(1)*ADJ_alnthv  
258       ADJ_z2   = s19v(2)*ADJ_alnthv 
259       ADJ_z3   = s19v(3)*ADJ_alnthv 
260       ADJ_z4   = s19v(4)*ADJ_alnthv
261       ADJ_z5   = s19v(5)*ADJ_alnthv
262       ADJ_z6   = s19v(6)*ADJ_alnthv
263 
264       ADJ_xd  = ADJ_z6*gx2 !+ ADJ_xd
265       ADJ_gx2 = xd*ADJ_z6   + ADJ_gx2
266       ADJ_xx  = ADJ_z5*gx2 !+ ADJ_xx
267       ADJ_gx2 = xx*ADJ_z5   + ADJ_gx2
268       ADJ_xx  = ADJ_z4*gg   + ADJ_xx
269       ADJ_gg  = xx*ADJ_z4  + ADJ_gg
270       ADJ_xd  = ADJ_z3*gg  + ADJ_xd
271       ADJ_gg  = xd*ADJ_z3  + ADJ_gg
272       ADJ_xd  = ADJ_z2*ggg + ADJ_xd
273       ADJ_ggg = xd*ADJ_z2  + ADJ_ggg
274       ADJ_xx  = ADJ_z1*ggg + ADJ_xx
275       ADJ_ggg = xx*ADJ_z1  + ADJ_ggg
276 
277       ADJ_xd  = 2.*xd*ADJ_xx + ADJ_xd
278 
279       ADJ_alnsig = ADJ_xd + ADJ_alnsig
280 
281    else if (ifreq .eq. 2) then 
282 
283       ADJ_z1 = s22h(1)*ADJ_alnthh
284       ADJ_z2 = s22h(2)*ADJ_alnthh
285       ADJ_z3 = s22h(3)*ADJ_alnthh
286       ADJ_z4 = s22h(4)*ADJ_alnthh
287       ADJ_z5 = s22h(5)*ADJ_alnthh
288       ADJ_z6 = s22h(6)*ADJ_alnthh
289  
290       ADJ_xd  = ADJ_z6*gx2
291       ADJ_gx2 = xd*ADJ_z6  + ADJ_gx2
292       ADJ_xx  = ADJ_z5*gx2
293       ADJ_gx2 = xx*ADJ_z5  + ADJ_gx2
294       ADJ_xx  = ADJ_z4*gg  + ADJ_xx
295       ADJ_gg  = xx*ADJ_z4  + ADJ_gg
296       ADJ_xd  = ADJ_z3*gg  + ADJ_xd
297       ADJ_gg  = xd*ADJ_z3  + ADJ_gg
298       ADJ_xd  = ADJ_z2*ggg + ADJ_xd
299       ADJ_ggg = xd*ADJ_z2  + ADJ_ggg
300       ADJ_xx  = ADJ_z1*ggg + ADJ_xx 
301       ADJ_ggg = xx*ADJ_z1  + ADJ_ggg
302 
303       ADJ_xd  =  2.*xd*ADJ_xx + ADJ_xd
304       ADJ_alnsig = ADJ_xd
305 
306       ! vertical
307       xd =  xd_save
308       xx =  xx_save
309 
310       alnthv = alnthv_save
311 
312       ADJ_z1  = s22v(1)*ADJ_alnthv
313       ADJ_z2  = s22v(2)*ADJ_alnthv
314       ADJ_z3  = s22v(3)*ADJ_alnthv
315       ADJ_z4  = s22v(4)*ADJ_alnthv
316       ADJ_z5  = s22v(5)*ADJ_alnthv
317       ADJ_z6  = s22v(6)*ADJ_alnthv
318 
319       ADJ_xd  = ADJ_z6*gx2
320       ADJ_gx2 = xd*ADJ_z6  + ADJ_gx2
321       ADJ_xx  = ADJ_z5*gx2 
322       ADJ_gx2 = xx*ADJ_z5  + ADJ_gx2
323       ADJ_xx  = ADJ_z4*gg  + ADJ_xx
324       ADJ_gg  = xx*ADJ_z4  + ADJ_gg
325       ADJ_xd  = ADJ_z3*gg  + ADJ_xd
326       ADJ_gg  = xd*ADJ_z3  + ADJ_gg
327       ADJ_xd  = ADJ_z2*ggg + ADJ_xd
328       ADJ_ggg = xd*ADJ_z2  + ADJ_ggg
329       ADJ_xx  = ADJ_z1*ggg + ADJ_xx
330       ADJ_ggg = xx*ADJ_z1  + ADJ_ggg
331       ADJ_xd  =  2.*xd*ADJ_xx + ADJ_xd
332       ADJ_alnsig = ADJ_xd     + ADJ_alnsig
333 
334    else if (ifreq .eq. 3) then 
335 
336       ADJ_z1  = s37h(1)*ADJ_alnthh
337       ADJ_z2  = s37h(2)*ADJ_alnthh
338       ADJ_z3  = s37h(3)*ADJ_alnthh
339       ADJ_z4  = s37h(4)*ADJ_alnthh
340       ADJ_z5  = s37h(5)*ADJ_alnthh
341       ADJ_z6  = s37h(6)*ADJ_alnthh
342 
343       ADJ_xd  = ADJ_z6*gx2
344       ADJ_gx2 = xd*ADJ_z6  + ADJ_gx2
345       ADJ_xx  = ADJ_z5*gx2 
346       ADJ_gx2 = xx*ADJ_z5  + ADJ_gx2
347       ADJ_xx  = ADJ_z4*gg  + ADJ_xx
348       ADJ_gg  = xx*ADJ_z4  + ADJ_gg
349       ADJ_xd  = ADJ_z3*gg  + ADJ_xd
350       ADJ_gg  = xd*ADJ_z3  + ADJ_gg
351       ADJ_xd  = ADJ_z2*ggg + ADJ_xd
352       ADJ_ggg = xd*ADJ_z2  + ADJ_ggg
353       ADJ_xx  = ADJ_z1*ggg + ADJ_xx
354       ADJ_ggg = xx*ADJ_z1  + ADJ_ggg
355       ADJ_xd  = 2.*xd*ADJ_xx + ADJ_xd
356       ADJ_alnsig = ADJ_xd
357 
358       ! vertical
359       xd =  xd_save
360       xx =  xx_save
361 
362       alnthv = alnthv_save
363 
364       ADJ_z1  = s37v(1)*ADJ_alnthv
365       ADJ_z2  = s37v(2)*ADJ_alnthv
366       ADJ_z3  = s37v(3)*ADJ_alnthv
367       ADJ_z4  = s37v(4)*ADJ_alnthv
368       ADJ_z5  = s37v(5)*ADJ_alnthv
369       ADJ_z6  = s37v(6)*ADJ_alnthv
370   
371       ADJ_xd  = ADJ_z6*gx2
372       ADJ_gx2 = xd*ADJ_z6  + ADJ_gx2
373       ADJ_xx  = ADJ_z5*gx2
374       ADJ_gx2 = xx*ADJ_z5  + ADJ_gx2
375       ADJ_xx  = ADJ_z4*gg  + ADJ_xx
376       ADJ_gg  = xx*ADJ_z4  + ADJ_gg
377       ADJ_xd  =ADJ_z3*gg   + ADJ_xd
378       ADJ_gg  = xd*ADJ_z3  + ADJ_gg
379       ADJ_xd  = ADJ_z2*ggg + ADJ_xd 
380       ADJ_ggg = xd*ADJ_z2  + ADJ_ggg
381       ADJ_xx  = ADJ_z1*ggg + ADJ_xx
382       ADJ_ggg = xx*ADJ_z1  + ADJ_ggg
383       ADJ_xd  = 2.*xd*ADJ_xx + ADJ_xd
384       ADJ_alnsig = ADJ_xd  + ADJ_alnsig
385 
386    else if (ifreq .eq. 4) then 
387 
388       ADJ_z1  = s85h(1)*ADJ_alnthh
389       ADJ_z2  = s85h(2)*ADJ_alnthh
390       ADJ_z3  = s85h(3)*ADJ_alnthh
391       ADJ_z4  = s85h(4)*ADJ_alnthh
392       ADJ_z5  = s85h(5)*ADJ_alnthh
393       ADJ_z6  = s85h(6)*ADJ_alnthh
394 
395       ADJ_xd  = ADJ_z6*gx2
396       ADJ_gx2 = xd*ADJ_z6  + ADJ_gx2
397       ADJ_xx  = ADJ_z5*gx2
398       ADJ_gx2 = xx*ADJ_z5  + ADJ_gx2
399       ADJ_xx  = ADJ_z4*gg  + ADJ_xx
400       ADJ_gg  = xx*ADJ_z4  + ADJ_gg
401       ADJ_xd  = ADJ_z3*gg  + ADJ_xd
402       ADJ_gg  = xd*ADJ_z3  + ADJ_gg
403       ADJ_xd  = ADJ_z2*ggg + ADJ_xd
404       ADJ_ggg = xd*ADJ_z2  + ADJ_ggg
405       ADJ_xx  = ADJ_z1*ggg + ADJ_xx
406       ADJ_ggg = xx*ADJ_z1  + ADJ_ggg
407       ADJ_xd  = 2.*xd*ADJ_xx + ADJ_xd
408       ADJ_alnsig = ADJ_xd
409 
410       ! vertical
411       xd =  xd_save
412       xx =  xx_save
413 
414       alnthv = alnthv_save
415 
416       ADJ_z1  = s85v(1)*ADJ_alnthv
417       ADJ_z2  = s85v(2)*ADJ_alnthv
418       ADJ_z3  = s85v(3)*ADJ_alnthv
419       ADJ_z4  = s85v(4)*ADJ_alnthv
420       ADJ_z5  = s85v(5)*ADJ_alnthv
421       ADJ_z6  = s85v(6)*ADJ_alnthv
422 
423       ADJ_xd  = ADJ_z6*gx2
424       ADJ_gx2 = xd*ADJ_z6  + ADJ_gx2
425       ADJ_xx  = ADJ_z5*gx2
426       ADJ_gx2 = xx*ADJ_z5  + ADJ_gx2
427       ADJ_xx  = ADJ_z4*gg  + ADJ_xx
428       ADJ_gg  = xx*ADJ_z4  + ADJ_gg
429       ADJ_xd  = ADJ_z3*gg  + ADJ_xd
430       ADJ_gg  = xd*ADJ_z3  + ADJ_gg
431       ADJ_xd  = ADJ_z2*ggg + ADJ_xd
432       ADJ_ggg = xd*ADJ_z2  + ADJ_ggg
433       ADJ_xx  = ADJ_z1*ggg + ADJ_xx
434       ADJ_ggg = xx*ADJ_z1  + ADJ_ggg
435       ADJ_xd  = 2.*xd*ADJ_xx  + ADJ_xd
436       ADJ_alnsig = ADJ_xd  + ADJ_alnsig
437    end if
438 
439    ADJ_gg  = ADJ_ggg*gx2   + ADJ_gg
440    ADJ_gx2 = gg*ADJ_ggg    + ADJ_gx2
441    ADJ_gx2 = 2.*gx2*ADJ_gg + ADJ_gx2
442 
443    if (abs(sigma) .gt. 0.) then
444       ADJ_sigma = ADJ_alnsig/sigma + ADJ_sigma
445    ! else
446    !    ADJ_sigma = 0.
447    end if
448 
449 end subroutine da_effang_adj
450 
451