Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_ssmi / da_roughem_adj.inc
blob0ad3d3a4359013abe7d0b54893a9c12145fb1f4e
1 subroutine da_roughem_adj (ifreq,gx2,tk,theta,remv,remh,ADJ_gx2,ADJ_tk,ADJ_remv,ADJ_remh)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
7    implicit none
9    !-----------------------------------------------------------------
10    ! Input  :: ADJ_tk, ADJ_gx2
11    ! Output :: ADJ_remv,ADJ_remh, remv,remh
12    !
13    ! Calculates rough-surface emissivity of ocean surface at SSM/I
14    ! frequencies.
16    integer, intent(in)    :: ifreq
17    real   , intent(in)    :: tk, theta, gx2
18    real   , intent(in)    :: ADJ_remv,ADJ_remh
19    real   , intent(inout) :: ADJ_tk,ADJ_gx2
20    real,    intent(out)   :: remv,remh
22    real :: ssw
24    real :: a19v(4),a22v(4),a37v(4),a85v(4)
25    real :: a19h(4),a22h(4),a37h(4),a85h(4)
26    real :: f(4)
27    real :: semv,semh,ADJ_semv,ADJ_semh,remv_save,remh_save
28    real :: tp,g,x1,x2,x3,x4,dtheta
29    real :: ADJ_tp,ADJ_g,ADJ_x1,ADJ_x2,ADJ_x3,ADJ_x4
31    data a19v/  -0.111E+01,   0.713E+00,  -0.624E-01,   0.212E-01 /
32    data a19h/   0.812E+00,  -0.215E+00,   0.255E-01,   0.305E-02 /
33    data a22v/  -0.134E+01,   0.911E+00,  -0.893E-01,   0.463E-01 /
34    data a22h/   0.958E+00,  -0.350E+00,   0.566E-01,  -0.262E-01 /
35    data a37v/  -0.162E+01,   0.110E+01,  -0.730E-01,   0.298E-01 /
36    data a37h/   0.947E+00,  -0.320E+00,   0.624E-01,  -0.300E-01 /
37    data a85v/  -0.145E+01,   0.808E+00,  -0.147E-01,  -0.252E-01 /
38    data a85h/   0.717E+00,  -0.702E-01,   0.617E-01,  -0.243E-01 /
40    data f/ 19.35, 22.235, 37.0, 85.5 /
42    if (trace_use) call da_trace_entry("da_roughem_adj")
44    semv      = 0.0
45    semh      = 0.0
46    ADJ_semv  = 0.0
47    ADJ_semh  = 0.0
48    remv_save = 0.0
49    remh_save = 0.0
50    tp        = 0.0
51    g         = 0.0
52    x1        = 0.0
53    x2        = 0.0
54    x3        = 0.0
55    x4        = 0.0
56    dtheta    = 0.0
57    ADJ_tp    = 0.0
58    ADJ_g     = 0.0
59    ADJ_x1    = 0.0
60    ADJ_x2    = 0.0
61    ADJ_x3    = 0.0
62    ADJ_x4    = 0.0
64    tp     = tk/t_roughem
65    dtheta = theta-53.0
66    g      = 0.5* gx2
67    x1     = g
68    x2     = tp*g
69    x3     = dtheta* g
70    x4     = tp*x3
72    if (ifreq == 1) then
73       remv = x1*a19v(1) +     x2*a19v(2) +     x3*a19v(3) +     x4*a19v(4)
74       remh = x1*a19h(1) +     x2*a19h(2) +     x3*a19h(3) +     x4*a19h(4)
75    else if (ifreq == 2) then
76       remv = x1*a22v(1) +     x2*a22v(2) +     x3*a22v(3) +     x4*a22v(4)
77       remh = x1*a22h(1) +     x2*a22h(2) +     x3*a22h(3) +     x4*a22h(4)
78    else if (ifreq == 3) then
79       remv = x1*a37v(1) +     x2*a37v(2) +     x3*a37v(3) +     x4*a37v(4)
80       remh = x1*a37h(1) +     x2*a37h(2) +     x3*a37h(3) +     x4*a37h(4)
81    else if (ifreq == 4) then
82       remv = x1*a85v(1) +     x2*a85v(2) +     x3*a85v(3) +     x4*a85v(4)
83       remh = x1*a85h(1) +     x2*a85h(2) +     x3*a85h(3) +     x4*a85h(4)
84    end if
86    ssw=36.5
87    call spemiss(f(ifreq),tk,theta,ssw,semv,semh)
89    remv_save = remv
90    remh_save = remh
92    ! start
94    ADJ_semh = ADJ_remh
95      
96    ADJ_semv = ADJ_remv
98    remv = remv_save
99    remh = remh_save
101    call da_spemiss_adj(f(ifreq),tk,theta,ssw,semv,semh, ADJ_tk,ADJ_semv,ADJ_semh          )
104    if (ifreq == 1) then
105       ADJ_x1 = ADJ_remh*a19h(1)
106       ADJ_x2 = ADJ_remh*a19h(2)
107       ADJ_x3 = ADJ_remh*a19h(3)
108       ADJ_x4 = ADJ_remh*a19h(4)
110       ADJ_x1 = ADJ_remv*a19v(1) + ADJ_x1
111       ADJ_x2 = ADJ_remv*a19v(2) + ADJ_x2
112       ADJ_x3 = ADJ_remv*a19v(3) + ADJ_x3
113       ADJ_x4 = ADJ_remv*a19v(4) + ADJ_x4
114    else if (ifreq == 2) then
115       ADJ_x1 = ADJ_remh*a22h(1) 
116       ADJ_x2 = ADJ_remh*a22h(2)
117       ADJ_x3 = ADJ_remh*a22h(3)
118       ADJ_x4 = ADJ_remh*a22h(4)
120       ADJ_x1 = ADJ_remv*a22v(1) + ADJ_x1
121       ADJ_x2 = ADJ_remv*a22v(2) + ADJ_x2
122       ADJ_x3 = ADJ_remv*a22v(3) + ADJ_x3
123       ADJ_x4 = ADJ_remv*a22v(4) + ADJ_x4
124    else if (ifreq == 3) then
125       ADJ_x1 = ADJ_remh*a37h(1)
126       ADJ_x2 = ADJ_remh*a37h(2)
127       ADJ_x3 = ADJ_remh*a37h(3)
128       ADJ_x4 = ADJ_remh*a37h(4)         
130       ADJ_x1 = ADJ_remv*a37v(1) + ADJ_x1
131       ADJ_x2 = ADJ_remv*a37v(2) + ADJ_x2
132       ADJ_x3 = ADJ_remv*a37v(3) + ADJ_x3
133       ADJ_x4 = ADJ_remv*a37v(4) + ADJ_x4
134    else if (ifreq == 4) then
135       ADJ_x1 = ADJ_remh*a85h(1) 
136       ADJ_x2 = ADJ_remh*a85h(2) 
137       ADJ_x3 = ADJ_remh*a85h(3)
138       ADJ_x4 = ADJ_remh*a85h(4) 
140       ADJ_x1 = ADJ_remv*a85v(1) + ADJ_x1
141       ADJ_x2 = ADJ_remv*a85v(2) + ADJ_x2
142       ADJ_x3 = ADJ_remv*a85v(3) + ADJ_x3
143       ADJ_x4 = ADJ_remv*a85v(4) + ADJ_x4
144      end if
146    ADJ_tp  = ADJ_x4*x3
147    ADJ_x3  = tp*ADJ_x4     + ADJ_x3
148    ADJ_g   = dtheta*ADJ_x3
149    ADJ_tp  = ADJ_x2*g      + ADJ_tp
150    ADJ_g   = tp*ADJ_x2     + ADJ_g
151    ADJ_g   = ADJ_x1        + ADJ_g
152    ADJ_gx2 = 0.5*ADJ_g     + ADJ_gx2
153    ADJ_tk  = ADJ_tp/t_roughem  + ADJ_tk
155    if (trace_use) call da_trace_exit("da_roughem_adj")
157 end subroutine da_roughem_adj