updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_ssmi / da_tb_tl.inc
blob6adba7100db3c8c5606a30529774b50f0135da13
1 subroutine da_tb_tl(ifreq,theta,p0,ta,gamma,sst,wv,hwv,u,alw,zcld,            &
2 !  tbv,tbh,                                                  &
3   TGL_p0,TGL_ta,TGL_gamma,TGL_sst,TGL_wv,                   &
4   TGL_hwv,TGL_u,TGL_alw,TGL_zcld,TGL_tbv,TGL_tbh           )\
6    !---------------------------------------------------------------------------
7    ! Purpose: TBD
8    ! Input  : TGL_p0,  TGL_ta,   TGL_gamma, TGL_sst, TGL_wv, TGL_hwv, TGL_u
9    !          TGL_alw, TGL_zcld
10    ! Output : TGL_tbv, TGL_tbh,  tbv,  tbh
11    ! ---------------------------------------------------------------------------
13    implicit none
15    integer, intent(in) :: ifreq
16    real,    intent(in) :: theta,p0,ta,gamma,sst,wv,hwv,u,alw,zcld
17    real,    intent(in) :: TGL_p0,TGL_ta,TGL_gamma,TGL_sst,TGL_wv, TGL_hwv,TGL_u,TGL_alw,TGL_zcld
18 !   real   , intent(in) :: tbv,tbh
19    real,    intent(out) :: TGL_tbv,TGL_tbh 
21    real :: freq(4),ebiasv(4),ebiash(4),cf1(4),cf2(4),cg(4)
23    real :: f,costhet,gx2,tbup,tbdn,tauatm,sigma,remv,remh
24    real :: effangv,effangh,tbdnv,dum,foam,emissv,emissh 
25    real :: refv,refh,semv,semh,scatv,scath,tbdnh
26    real :: TGL_gx2,TGL_tbup,TGL_tbdn,TGL_tauatm,TGL_sigma
27    real :: TGL_remh,TGL_effangv,TGL_effangh
28    ! real :: TGL_tremv
29    real :: TGL_tbdnh,TGL_dum,TGL_foam,TGL_emissv
30    real :: TGL_emissh,TGL_refv,TGL_refh,TGL_semv
31    real :: TGL_semh,TGL_scatv,TGL_scath,TGL_remv,TGL_tbdnv
32    real :: TGL_theta
34    real :: fem
36    data freq/19.35,22.235,37.0,85.5/
38    ! empirical bias corrections for surface emissivity
40    data ebiasv/0.0095, 0.005,-0.014, -0.0010/
41    data ebiash/0.004,   0.0,-0.023, 0.043/
43    data cf1/0.0015,0.004,0.0027,0.005/
44    data cg/4.50e-3, 5.200e-3, 5.5e-3, 7.2e-3 /
46    data cf2/0.0023,0.000,0.0002,-0.006/
48    ! 'foam' emissivity
49    data fem /1.0/
51    if (trace_use) call da_trace_entry("da_tb_tl")
53    f = freq(ifreq)
54    costhet = cos(theta*0.017453)
56    ! effective surface slope variance
58    gx2 = cg(ifreq)*    u
59    TGL_gx2 = cg(ifreq)*TGL_u
61    ! get upwelling atmospheric brightness temperature
63    TGL_theta=0.0
64    call da_tbatmos_tl(ifreq,theta,p0,wv,hwv,ta,gamma,alw,zcld, &
65                        tbup,tbdn,tauatm,                        &
66                        TGL_theta,TGL_p0,TGL_wv,TGL_hwv,TGL_ta,TGL_gamma,  &
67                        TGL_alw,TGL_zcld,TGL_tbup,TGL_tbdn,      &
68                        TGL_tauatm                              )
70    ! convert transmittance to optical depth
72    sigma = -alog(tauatm)*costhet
73    TGL_sigma = -costhet*TGL_tauatm/tauatm
75    ! get rough surface emissivity
77    call da_roughem_tl(ifreq,gx2,sst,theta,remv,remh,         &
78                        TGL_gx2,TGL_sst,TGL_remv,TGL_remh     )
80    ! get effective zenith angles for scattered radiation at surface
82    call da_effang_tl(ifreq,theta,gx2,sigma,effangv,effangh,    &
83                       TGL_gx2,TGL_sigma,TGL_effangv,TGL_effangh)
85    ! get effective sky brightness temperatures for scattered radiation
87    call da_tbatmos_tl(ifreq,effangv,p0,wv,hwv,ta,gamma,alw,    &
88                        zcld,dum,tbdnv,dum,                      &
89                        TGL_effangv,TGL_p0,TGL_wv,TGL_hwv,       &
90                        TGL_ta,TGL_gamma,TGL_alw,TGL_zcld,       & 
91                        TGL_dum,TGL_tbdnv,TGL_dum               )
93    call da_tbatmos_tl(ifreq,effangh,p0,wv,hwv,ta,gamma,alw,    &
94                        zcld,dum,tbdnh,dum,                      &
95                        TGL_effangh,TGL_p0,TGL_wv,TGL_hwv,       &
96                        TGL_ta,TGL_gamma,TGL_alw,TGL_zcld,       &
97                        TGL_dum,TGL_tbdnh,TGL_dum               )
99    ! compute 'foam' coverage
101    foam = cf1(ifreq)*    u
102    TGL_foam = cf1(ifreq)*TGL_u
104    if (u .gt. 5.0) then
105       TGL_foam = TGL_foam + cf2(ifreq)*TGL_u
106       foam =     foam + cf2(ifreq)*(  u-5.0)
107    end if
109    ! compute surface emissivities and reflectivity
111    emissv     =     foam*fem + (1.0 - foam)*(remv + ebiasv(ifreq))
112    TGL_emissv = TGL_foam*fem - TGL_foam*(remv + ebiasv(ifreq)) &
113                                 + (1.0 - foam)*TGL_remv
114    emissh     =     foam*fem + (1.0 - foam)*(remh + ebiash(ifreq))
115    TGL_emissh = TGL_foam*fem - TGL_foam*(remh + ebiash(ifreq)) &
116                                 + (1.0 - foam)*TGL_remh
117    refv     =   1.0 - emissv
118    TGL_refv = - TGL_emissv
119    refh     =   1.0 - emissh
120    TGL_refh = - TGL_emissh
122    ! compute surface emission term
124    semv     = sst*emissv
125    TGL_semv = TGL_sst*emissv + sst*TGL_emissv
126    semh     = sst*emissh
127    TGL_semh = TGL_sst*emissh + sst*TGL_emissh
129    ! compute surface scattering term
131    scatv     = refv*tbdnv
132    TGL_scatv = TGL_refv*tbdnv + refv*TGL_tbdnv
133    scath     = refh*tbdnh
134    TGL_scath = TGL_refh*tbdnh + refh*TGL_tbdnh
136    ! combine to get space-observed brightness temperature
138    ! tbv     =     tbup + tauatm*(semv + scatv)
139    TGL_tbv = TGL_tbup + TGL_tauatm*(semv + scatv)   &
140                          + tauatm*(TGL_semv + TGL_scatv)
141    ! tbh     =     tbup + tauatm*(semh + scath)
142    TGL_tbh = TGL_tbup + TGL_tauatm*(semh + scath)   &
143                          + tauatm*(TGL_semh + TGL_scath)
145    if (trace_use) call da_trace_exit("da_tb_tl")
147 end subroutine da_tb_tl