wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_ships / da_get_innov_vector_ships.inc
blobbf50ce6c4e9b7076e1231ff01f8f9628969eb819
1 subroutine da_get_innov_vector_ships( it,num_qcstat_conv, grid, ob, iv)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD     
5    !    Updated for Analysis on Arakawa-C grid
6    !    Author: Syed RH Rizvi,  MMM/ESSL/NCAR,  Date: 10/22/2008
7    !-----------------------------------------------------------------------
9    implicit none
11    integer,          intent(in)    :: it      ! External iteration.
12    type(domain),     intent(in)    :: grid    ! first guess state.
13    type(y_type),     intent(inout) :: ob      ! Observation structure.
14    type(iv_type),    intent(inout) :: iv      ! O-B structure.
15    integer,          intent(inout) :: num_qcstat_conv(:,:,:,:)
17    integer :: n        ! Loop counter.
18    integer :: i, j, k  ! Index dimension.
19    real    :: dx, dxm  ! Interpolation weights.
20    real    :: dy, dym  ! Interpolation weights.
21    real, allocatable :: model_u(:,:)  ! Model value u at oblocation.
22    real, allocatable :: model_v(:,:)  ! Model value v at oblocation.
23    real, allocatable :: model_t(:,:)  ! Model value t at oblocation.
24    real, allocatable :: model_p(:,:)  ! Model value p at oblocation.
25    real, allocatable :: model_q(:,:)  ! Model value q at oblocation.
26    real, allocatable :: model_hsm(:,:)
28    real    :: v_h(kms:kme)   ! Model value h at ob hor. location.
29    real    :: v_p(kms:kme)   ! Model value p at ob hor. location.
31    real    :: hd, psfcm
32    real    :: ho, to, qo
33    
34    if (trace_use_dull) call da_trace_entry("da_get_innov_vector_ships")
36    allocate (model_u(1,iv%info(ships)%n1:iv%info(ships)%n2))
37    allocate (model_v(1,iv%info(ships)%n1:iv%info(ships)%n2))
38    allocate (model_t(1,iv%info(ships)%n1:iv%info(ships)%n2))
39    allocate (model_p(1,iv%info(ships)%n1:iv%info(ships)%n2))
40    allocate (model_q(1,iv%info(ships)%n1:iv%info(ships)%n2))
41    allocate (model_hsm(1,iv%info(ships)%n1:iv%info(ships)%n2))
43    if ( it > 1 ) then
44       do n=iv%info(ships)%n1,iv%info(ships)%n2
45          if (iv%ships(n)%u%qc == fails_error_max) iv%ships(n)%u%qc = 0
46          if (iv%ships(n)%v%qc == fails_error_max) iv%ships(n)%v%qc = 0
47          if (iv%ships(n)%t%qc == fails_error_max) iv%ships(n)%t%qc = 0
48          if (iv%ships(n)%p%qc == fails_error_max) iv%ships(n)%p%qc = 0
49          if (iv%ships(n)%q%qc == fails_error_max) iv%ships(n)%q%qc = 0
50       end do
51    end if
53    if (sfc_assi_options == sfc_assi_options_1) then
54       do n=iv%info(ships)%n1,iv%info(ships)%n2
56          ! [1.1] Get horizontal interpolation weights:
58          i   = iv%info(ships)%i(1,n)
59          j   = iv%info(ships)%j(1,n)
60          dx  = iv%info(ships)%dx(1,n)
61          dy  = iv%info(ships)%dy(1,n)
62          dxm = iv%info(ships)%dxm(1,n)
63          dym = iv%info(ships)%dym(1,n)
65          ! Surface correction
67          iv%ships(n)%p%inv = ob%ships(n)%p
68          iv%ships(n)%t%inv = ob%ships(n)%t
69          iv%ships(n)%q%inv = ob%ships(n)%q
70          iv%ships(n)%u%inv = ob%ships(n)%u
71          iv%ships(n)%v%inv = ob%ships(n)%v
73          if (iv % ships(n) % h > missing_r) then
74             do k=kts,kte
75               v_h(k) = dym*(dxm*grid%xb%h(i,j  ,k) + dx*grid%xb%h(i+1,j  ,k)) &
76                  + dy *(dxm*grid%xb%h(i,j+1,k) + dx*grid%xb%h(i+1,j+1,k))
77             end do
79             hd = v_h(kts) - iv % ships(n) % h
81             if (abs(hd) <= Max_StHeight_Diff .or. anal_type_verify) then
82                if (iv % ships(n) % h < v_h(kts)) then
83                   iv%info(ships)%zk(:,n) = 1.0+1.0e-6
84                   call da_obs_sfc_correction(iv%info(ships), iv%ships(n), n, grid%xb)
85                else
86                   call da_to_zk(iv % ships(n) % h, v_h, v_interp_h, iv%info(ships)%zk(1,n))
87                end if
88             end if
89          else if (ob % ships(n) % p > 1.0) then
90             do k=kts,kte
91               v_p(k) = dym*(dxm*grid%xb%p(i,j  ,k) + dx*grid%xb%p(i+1,j  ,k)) &
92                        + dy *(dxm*grid%xb%p(i,j+1,k) + dx*grid%xb%p(i+1,j+1,k))
93             end do
95             call da_to_zk(ob % ships(n) % p, v_p, v_interp_p, iv%info(ships)%zk(1,n))
97             if (iv%info(ships)%zk(1,n) < 0.0 .and.  .not.anal_type_verify) then
98                iv % ships(n) % p % inv = missing_r
99                iv % ships(n) % p % qc  = missing_data
100                iv%info(ships)%zk(:,n) = 1.0+1.0e-6
101             end if
102          end if
103       end do
105       call da_convert_zk (iv%info(ships))
107       if (.not.anal_type_verify) then
108          do n=iv%info(ships)%n1,iv%info(ships)%n2
109             if (iv%info(ships)%zk(1,n) < 0.0) then
110                iv % ships(n) % u % qc = missing_data
111                iv % ships(n) % v % qc = missing_data
112                iv % ships(n) % t % qc = missing_data
113                iv % ships(n) % q % qc = missing_data
114                iv % ships(n) % p % qc = missing_data
115             end if
116          end do
117       end if
119       ! Interpolate horizontally:
120 #ifdef A2C
121       call da_interp_lin_3d (grid%xb%u, iv%info(ships), model_u,'u')
122       call da_interp_lin_3d (grid%xb%v, iv%info(ships), model_v,'v')
123 #else
124       call da_interp_lin_3d (grid%xb%u, iv%info(ships), model_u)
125       call da_interp_lin_3d (grid%xb%v, iv%info(ships), model_v)
126 #endif
127       call da_interp_lin_3d (grid%xb%t, iv%info(ships), model_t)
128       call da_interp_lin_3d (grid%xb%q, iv%info(ships), model_q)
129       call da_interp_lin_3d (grid%xb%p, iv%info(ships), model_p)
131    else if (sfc_assi_options == sfc_assi_options_2) then
132       ! Surface data assimilation approach 2
134       ! 1.2.1 Surface assmiilation approach 2(10-m u, v, 2-m t, q, and 
135       ! sfc_p)
137       call da_interp_lin_2d (grid%xb%u10,  iv%info(ships), 1, model_u)
138       call da_interp_lin_2d (grid%xb%v10,  iv%info(ships), 1, model_v)
139       call da_interp_lin_2d (grid%xb%t2,   iv%info(ships), 1, model_t)
140       call da_interp_lin_2d (grid%xb%q2,   iv%info(ships), 1, model_q)
141       call da_interp_lin_2d (grid%xb%psfc, iv%info(ships), 1, model_p)
143       do n=iv%info(ships)%n1,iv%info(ships)%n2
144          iv%ships(n)%p%inv = ob%ships(n)%p
145          iv%ships(n)%t%inv = ob%ships(n)%t
146          iv%ships(n)%q%inv = ob%ships(n)%q
147          iv%ships(n)%u%inv = ob%ships(n)%u
148          iv%ships(n)%v%inv = ob%ships(n)%v
150          if (iv%ships(n)%p%qc >= 0) then
152             ! model surface p, t, q, h at observed site:
154             call da_interp_lin_2d_partial (grid%xb%terr, iv%info(ships), 1, n, n, model_hsm(:,n))
156             ho = iv%ships(n)%h
157             to = -888888.0
158             qo = -888888.0
160             if (iv%ships(n)%t%qc >= 0 .and. iv%ships(n)%q%qc >= 0) then
161                to = ob%ships(n)%t
162                qo = ob%ships(n)%q
163                call da_sfc_pre(psfcm, model_p(1,n), model_t(1,n), model_q(1,n), model_hsm(1,n), ho, to, qo)
164             else if (iv%ships(n)%t%qc >= 0 .and. iv%ships(n)%q%qc < 0) then
165                 to = ob%ships(n)%t
166                 call da_sfc_pre(psfcm, model_p(1,n), model_t(1,n), model_q(1,n), model_hsm(1,n), ho, to)
167             else
168                 call da_sfc_pre(psfcm, model_p(1,n), model_t(1,n), model_q(1,n), model_hsm(1,n), ho)
169             end if
171             ! Pressure at the observed height:
172             model_p(1,n) = psfcm
173          end if
174       end do
175    end if
177    do n=iv%info(ships)%n1,iv%info(ships)%n2
179       !-----------------------------------------------------------------------
180       ! [3.0] Fast interpolation:
181       !-----------------------------------------------------------------------
183       if (ob % ships(n) % u > missing_r .AND. iv % ships(n) % u % qc >= obs_qc_pointer) then
184          iv % ships(n) % u % inv = iv%ships(n)%u%inv - model_u(1,n)
185       else
186          iv % ships(n) % u % inv = 0.0
187       end if
189       if (ob % ships(n) % v > missing_r .AND. iv % ships(n) % v % qc >= obs_qc_pointer) then
190          iv % ships(n) % v % inv = iv%ships(n)%v%inv - model_v(1,n)
191       else
192          iv % ships(n) % v % inv = 0.0
193       end if
195       if (ob % ships(n) % p > 0.0 .AND. iv % ships(n) % p % qc >= obs_qc_pointer) then
196          iv % ships(n) % p % inv = iv%ships(n)%p%inv - model_p(1,n)
197       else
198          iv % ships(n) % p % inv = 0.0
199       end if
201       if (ob % ships(n) % t > 0.0 .AND. iv % ships(n) % t % qc >= obs_qc_pointer) then
202          iv % ships(n) % t % inv = iv%ships(n)%t%inv - model_t(1,n)
203       else
204          iv % ships(n) % t % inv = 0.0
205       end if
207       if (ob % ships(n) % q > 0.0 .AND. iv % ships(n) % q % qc >= obs_qc_pointer) then
208          iv % ships(n) % q % inv = iv%ships(n)%q%inv - model_q(1,n)
209       else
210          iv % ships(n) % q % inv = 0.0
211       end if
212    end do
214    !---------------------------------------------------------------------
215    ! [5.0] Perform optional maximum error check:
216    !---------------------------------------------------------------------
218    if ( check_max_iv ) &
219       call da_check_max_iv_ships(iv,ob, it, num_qcstat_conv)
220    
221    deallocate (model_u)
222    deallocate (model_v)
223    deallocate (model_t)
224    deallocate (model_p)
225    deallocate (model_q)
226    deallocate (model_hsm)
227    
228    if (trace_use_dull) call da_trace_exit("da_get_innov_vector_ships")
230 end subroutine da_get_innov_vector_ships