wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_pseudo / da_get_innov_vector_pseudo.inc
blobd7503d957596076c836a70a285b687c4709adc12
1 subroutine da_get_innov_vector_pseudo(grid, ob, iv)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
7    implicit none
8    
9    type(domain),     intent(in)    :: grid        ! Background structure 
10    type(y_type),     intent(inout) :: ob          ! Observation structure.
11    type(iv_type),    intent(inout) :: iv          ! O-B structure.
13    integer :: n        ! Loop counter.
15    real, allocatable :: model_u(:,:)
16    real, allocatable :: model_v(:,:)
17    real, allocatable :: model_q(:,:)
18    real, allocatable :: model_p(:,:)
19    real, allocatable :: model_t(:,:)
20    
21    if (trace_use_dull) call da_trace_entry("da_get_innov_vector_pseudo")
24    allocate (model_u(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
25    allocate (model_v(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
26    allocate (model_q(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
27    allocate (model_p(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
28    allocate (model_t(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
30    call da_convert_zk (iv%info(pseudo))
32 #ifdef A2C
33    call da_interp_lin_3d (grid%xb%u, iv%info(pseudo), model_u,'u')
34    call da_interp_lin_3d (grid%xb%v, iv%info(pseudo), model_v,'v')
35 #else
36    call da_interp_lin_3d (grid%xb%u, iv%info(pseudo), model_u)
37    call da_interp_lin_3d (grid%xb%v, iv%info(pseudo), model_v)
38 #endif
39    call da_interp_lin_3d (grid%xb%t, iv%info(pseudo), model_t)
40    call da_interp_lin_3d (grid%xb%p, iv%info(pseudo), model_p)
41    call da_interp_lin_3d (grid%xb%q, iv%info(pseudo), model_q)
43    do n=iv%info(pseudo)%n1,iv%info(pseudo)%n2
44       !---------------------------------------------------------------
45       ! [3.0] Calculate observation O = B +(O-B):
46       !---------------------------------------------------------------
48       select case(pseudo_var(1:1))
49       case ('u', 'U')
50          if (ob % pseudo(n) % u > missing_r) then
51              iv % pseudo(n) % u % inv = ob%pseudo(n)%u - model_u(1,n)
52          else
53              ob % pseudo(n) % u = model_u(1,n) + iv % pseudo(n) % u % inv
54          endif   
55       case ('v', 'V')
56          if (ob % pseudo(n) % v > missing_r) then
57              iv % pseudo(n) % v % inv = ob%pseudo(n)%v - model_v(1,n)
58          else
59              ob % pseudo(n) % v = model_v(1,n) + iv % pseudo(n) % v % inv
60          endif   
61       case ('t', 'T')
62          if (ob % pseudo(n) % t > missing_r) then
63              iv % pseudo(n) % t % inv = ob%pseudo(n)%t - model_t(1,n)
64          else
65              ob % pseudo(n) % t = model_t(1,n) + iv % pseudo(n) % t % inv
66          endif   
67       case ('p', 'P')
68          if (ob % pseudo(n) % p > missing_r) then
69              iv % pseudo(n) % p % inv = ob%pseudo(n)%p - model_p(1,n)
70          else
71              ob % pseudo(n) % p = model_p(1,n) + iv % pseudo(n) % p % inv
72          endif   
73       case ('q', 'Q')
74          if (ob % pseudo(n) % q > missing_r) then
75              iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_q(1,n)
76          else
77              ob % pseudo(n) % q = model_q(1,n) + iv % pseudo(n) % q % inv
78          endif   
79       end select 
80    end do
82    deallocate (model_u)
83    deallocate (model_v)
84    deallocate (model_q)
85    deallocate (model_p)
86    deallocate (model_t)
87    
88    if (trace_use_dull) call da_trace_exit("da_get_innov_vector_pseudo")
90 end subroutine da_get_innov_vector_pseudo