updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_ssmi / da_transform_xtoy_ssmi_tb.inc
blob967bd26b93d3a00a711bc00025cf84bedaa5d7e6
1 subroutine da_transform_xtoy_ssmi_tb(grid, iv, y)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
7    implicit none
9    type (domain),  intent(in)    :: grid
10    type (iv_type), intent(in)    :: iv       ! obs. increment vector (o-b).
11    type (y_type),  intent(inout) :: y        ! y = h (grid%xa)
13    integer                      :: n        ! loop counter.
15    ! real, dimension(mix,mjy)     :: TGL_dotspeed
16    ! real, dimension(mix,mjy)     :: psfc,ta,gamma,sst,htpw,speed,alw,zcld,tpw
17    ! real, dimension(mix,mjy)     :: TGL_psfc,TGL_ta,TGL_gamma,TGL_sst,TGL_tpw
18    ! real, dimension(mix,mjy)     :: TGL_htpw,TGL_speed,TGL_alw,TGL_zcld         
20    ! real, dimension(mix,mjy)     :: tb19v,tb19h, &
21    !                                  tb22v,       &
22    !                                  tb37v,tb37h, &
23    !                                  tb85v,tb85h
25    ! real, dimension(mix,mjy)     :: TGL_tb19v,TGL_tb19h, &
26    !                                  TGL_tb22v,           &
27    !                                  TGL_tb37v,TGL_tb37h, &
28    !                                  TGL_tb85v,TGL_tb85h
30    ! real, dimension(mkz)         :: zh
31    ! real, dimension(mkz+1)       :: zf
32    !----------------------------------------------------------------------------
34    ! Wind at 1st level at dot pivnts
36    ! speed, surface wind speed (m/sec)    (0 - 30) , take 10 m wind
38    ! total precipitable water in cm
40    ! call da_transform_xtotpw(grid%xa, grid%xb)
42    ! do j=1,mjy-1
43    !    do i=1,mix-1
44    !       zf(1) = grid%xb%ztop
45    !       do k=1,mkz
46    !          zh(k)=grid%xb%h(i,j,k)
47    !          zf(k+1)=grid%xb%hf(i,j,k)
48    !       end do
50    !       ! surface pressure (mb) (940 -1030)
52    !       psfc(i,j)     = (grid%xb%p(i,j,mkz)+grid%xb%Psac(i,j)*(1.0-grid%xb % sigmah(mkz)))*0.001
53    !       TGL_psfc(i,j) =  grid%xa%p  (i,j,mkz)*0.01
55    !       ! sea surface temperature (k) (273 - 303) (doesnot change) 
57    !       sst(i,j)      = grid%xb%tgrn(i,j)
58    !       TGL_sst(i,j)  = grid%xa%tgrn(1,1)
59    !       TGL_sst(i,j)  = 0.0
61    !       ! effective surface air temperature (263 - 303)
63    !       ta(i,j)     = grid%xb%tgrn(i,j) + &
64    !                   (grid%xb%t(i,j,mkz)-grid%xb%tgrn(i,j))*log(2.0/0.0001)/ &
65    !                   log(zh(mkz)/0.0001)
67    !       TGL_ta(i,j) = grid%xa%tgrn(1,1) + &
68    !                   (grid%xa%t(i,j,mkz)-grid%xa%tgrn(1,1))*log(2.0/0.0001)/ &
69    !                   log(zh(mkz)/0.0001)
71    !       TGL_ta(i,j) = (grid%xa%t(i,j,mkz)-0.)*log(2.0/0.0001)/ &
72    !                   log(zh(mkz)/0.0001)
74    !       ! gamma is an emperical formula and zcld is given for now
76    !       gamma(i,j) = (ta(i,j)-270)*0.023 + 5.03  ! effective lapse rate   (km) (4.0 - 6.5)
77    !       zcld(i,j)  = 1                           ! effective cloud height (km)
78                                                      ! = 1 if no cloud infomation
79    !       TGL_gamma(i,j) = TGL_ta(i,j)*0.023
80    !       TGL_zcld(i,j)  = 0.0
82    !       ! total precipitable water in (kg/m**2) (0 - 70)
84    !       tpw(i,j)     = grid%xb%tpw(i,j)*10.0
85    !       TGL_tpw(i,j) = grid%xa%tpw(i,j)*10.0
87    !       ! Column liquid water (kg/m**2)   (0 - 0.5) (no data now. So, do it later.)
89    !       alw(i,j)     = 0.0
90    !       TGL_alw(i,j) = 0.0
92    !       ! Column height weighted mivsture density on the grid locally 
94    !       call da_transform_xtozrhoq(grid%xb, i, j, zh, zf, zrhom)
95    !       call da_transform_xtozrhoq_lin(grid%xb, grid%xa, i, j, zh, zf, TGL_zrhom)
97    !       ! Column mivsture density on the grid locally
99    !       htpw(i,j)     = zrhom/tpw(i,j)/1000.0
100    !       TGL_htpw(i,j) = TGL_zrhom/tpw(i,j)/1000.0 &
101    !                   - TGL_tpw(i,j)/tpw(i,j)*htpw(i,j)
103    !    end do
104    ! end do
106    ! do j=1,mjy-1
107    !    do i=1,mix-1
109    !       call tgl_tb(1,53.,psfc(i,j),ta(i,j),gamma(i,j),sst(i,j),tpw(i,j),  &
110    !          htpw(i,j),speed(i,j),alw(i,j),zcld(i,j),               &
111    !          grid%xb%tb19v(i,j),grid%xb%tb19h(i,j),                           &
112    !          TGL_psfc(i,j),TGL_ta(i,j),TGL_gamma(i,j),TGL_sst(i,j), &
113    !          TGL_tpw(i,j),TGL_htpw(i,j),TGL_speed(i,j),TGL_alw(i,j),&
114    !          TGL_zcld(i,j),TGL_tb19v(i,j),TGL_tb19h(i,j)           )
116    !       call TGL_tb(2,53.,psfc(i,j),ta(i,j),gamma(i,j),sst(i,j),tpw(i,j),  &
117    !          htpw(i,j),speed(i,j),alw(i,j),zcld(i,j),               &
118    !          grid%xb%tb22v(i,j),dum1,                                    &
119    !          TGL_psfc(i,j),TGL_ta(i,j),TGL_gamma(i,j),TGL_sst(i,j), &
120    !          TGL_tpw(i,j),TGL_htpw(i,j),TGL_speed(i,j),TGL_alw(i,j),&
121    !          TGL_zcld(i,j),TGL_tb22v(i,j),dum2                     )
123    !       call TGL_tb(3,53.,psfc(i,j),ta(i,j),gamma(i,j),sst(i,j),tpw(i,j),  &
124    !          htpw(i,j),speed(i,j),alw(i,j),zcld(i,j),               &
125    !          grid%xb%tb37v(i,j),grid%xb%tb37h(i,j),                           &
126    !          TGL_psfc(i,j),TGL_ta(i,j),TGL_gamma(i,j),TGL_sst(i,j), &
127    !          TGL_tpw(i,j),TGL_htpw(i,j),TGL_speed(i,j),TGL_alw(i,j),&
128    !          TGL_zcld(i,j),TGL_tb37v(i,j),TGL_tb37h(i,j)           )
130    !       call TGL_tb(4,53.,psfc(i,j),ta(i,j),gamma(i,j),sst(i,j),tpw(i,j),  &
131    !          htpw(i,j),speed(i,j),alw(i,j),zcld(i,j),               &
132    !          grid%xb%tb85v(i,j),grid%xb%tb85h(i,j),                           &
133    !          TGL_psfc(i,j),TGL_ta(i,j),TGL_gamma(i,j),TGL_sst(i,j), &
134    !          TGL_tpw(i,j),TGL_htpw(i,j),TGL_speed(i,j),TGL_alw(i,j),&
135    !          TGL_zcld(i,j),TGL_tb85v(i,j),TGL_tb85h(i,j)           )
136    !    end do
137    ! end do
139    real, allocatable :: tb19v(:)
140    real, allocatable :: tb19h(:)
141    real, allocatable :: tb22v(:)
142    real, allocatable :: tb37v(:)
143    real, allocatable :: tb37h(:)
144    real, allocatable :: tb85v(:)
145    real, allocatable :: tb85h(:)
147    if (trace_use) call da_trace_entry("da_transform_xtoy_ssmi_tb")
149    allocate(tb19v(iv%info(ssmi_tb)%n1:iv%info(ssmi_tb)%n2))
150    allocate(tb19h(iv%info(ssmi_tb)%n1:iv%info(ssmi_tb)%n2))
151    allocate(tb22v(iv%info(ssmi_tb)%n1:iv%info(ssmi_tb)%n2))
152    allocate(tb37v(iv%info(ssmi_tb)%n1:iv%info(ssmi_tb)%n2))
153    allocate(tb37h(iv%info(ssmi_tb)%n1:iv%info(ssmi_tb)%n2))
154    allocate(tb85v(iv%info(ssmi_tb)%n1:iv%info(ssmi_tb)%n2))
155    allocate(tb85h(iv%info(ssmi_tb)%n1:iv%info(ssmi_tb)%n2))
157    call da_interp_lin_2d (grid%xa%tb19v, iv%info(ssmi_tb), 1, tb19v)
158    call da_interp_lin_2d (grid%xa%tb19h, iv%info(ssmi_tb), 1, tb19h)
159    call da_interp_lin_2d (grid%xa%tb22v, iv%info(ssmi_tb), 1, tb22v)
160    call da_interp_lin_2d (grid%xa%tb37v, iv%info(ssmi_tb), 1, tb37v)
161    call da_interp_lin_2d (grid%xa%tb37h, iv%info(ssmi_tb), 1, tb37h)
162    call da_interp_lin_2d (grid%xa%tb85v, iv%info(ssmi_tb), 1, tb85v)
163    call da_interp_lin_2d (grid%xa%tb85h, iv%info(ssmi_tb), 1, tb85h)
165    do n=iv%info(ssmi_tb)%n1,iv%info(ssmi_tb)%n2
166       y%ssmi_tb(n)%tb19v = tb19v(n)
167       y%ssmi_tb(n)%tb19h = tb19h(n)
168       y%ssmi_tb(n)%tb22v = tb22v(n)
169       y%ssmi_tb(n)%tb37v = tb37v(n)
170       y%ssmi_tb(n)%tb37h = tb37h(n)
171       y%ssmi_tb(n)%tb85v = tb85v(n)
172       y%ssmi_tb(n)%tb85h = tb85h(n)
173    end do
175    deallocate(tb19v)
176    deallocate(tb19h)
177    deallocate(tb22v)
178    deallocate(tb37v)
179    deallocate(tb37h)
180    deallocate(tb85v)
181    deallocate(tb85h)
182     
183    ! WHY?
184    !          y%ssmi_tb(n)%tb19v = hor_interp(dxm,dx,dym,dy, &
185    !              TGL_tb19v(i,j ),TGL_tb19v(i+1,j ), &
186    !              TGL_tb19v(i,j+1),TGL_tb19v(i+1,j+1) )
188    !          y%ssmi_tb(n)%tb19h = hor_interp(dxm,dx,dym,dy, &
189    !              TGL_tb19h(i,j ),TGL_tb19h(i+1,j ), &
190    !              TGL_tb19h(i,j+1),TGL_tb19h(i+1,j+1) )
192    !          y%ssmi_tb(n)%tb22v = hor_interp(dxm,dx,dym,dy, &
193    !              TGL_tb22v(i,j ),TGL_tb22v(i+1,j ), &
194    !              TGL_tb22v(i,j+1),TGL_tb22v(i+1,j+1) )
196    !          y%ssmi_tb(n)%tb37v = hor_interp(dxm,dx,dym,dy, &
197    !              TGL_tb37v(i,j ),TGL_tb37v(i+1,j ), &
198    !              TGL_tb37v(i,j+1),TGL_tb37v(i+1,j+1) )
200    !          y%ssmi_tb(n)%tb37h = hor_interp(dxm,dx,dym,dy, &
201    !              TGL_tb37h(i,j ),TGL_tb37h(i+1,j ), &
202    !              TGL_tb37h(i,j+1),TGL_tb37h(i+1,j+1) )
204    !          y%ssmi_tb(n)%tb85v = hor_interp(dxm,dx,dym,dy, &
205    !              TGL_tb85v(i,j ),TGL_tb85v(i+1,j ), &
206    !              TGL_tb85v(i,j+1),TGL_tb85v(i+1,j+1) )
208    !          y%ssmi_tb(n)%tb85h = hor_interp(dxm,dx,dym,dy, &
209    !              TGL_tb85h(i,j ),TGL_tb85h(i+1,j ), &
210    !              TGL_tb85h(i,j+1),TGL_tb85h(i+1,j+1) )
211    !       end if
212    !    end if
214    ! end do
216    if (trace_use) call da_trace_exit("da_transform_xtoy_ssmi_tb")
218 end subroutine da_transform_xtoy_ssmi_tb