updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_pseudo / da_get_innov_vector_pseudo.inc
blob39f1a4356d5a3ceea0a6cff3c3ccf35d54be6fac
1 subroutine da_get_innov_vector_pseudo(it, grid, ob, iv)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
7    implicit none
9    integer,          intent(in)    :: it          ! External iteration
10    type(domain),     intent(in)    :: grid        ! Background structure 
11    type(y_type),     intent(inout) :: ob          ! Observation structure.
12    type(iv_type),    intent(inout) :: iv          ! O-B structure.
14    integer :: n        ! Loop counter.
16    real, allocatable :: model_u(:,:)
17    real, allocatable :: model_v(:,:)
18    real, allocatable :: model_q(:,:)
19    real, allocatable :: model_p(:,:)
20    real, allocatable :: model_t(:,:)
22    real, allocatable :: model_qcw(:,:)
23    real, allocatable :: model_qci(:,:)
24    real, allocatable :: model_qrn(:,:)
25    real, allocatable :: model_qsn(:,:)
26    real, allocatable :: model_qgr(:,:)
28    if (trace_use_dull) call da_trace_entry("da_get_innov_vector_pseudo")
30    allocate (model_u(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
31    allocate (model_v(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
32    allocate (model_q(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
33    allocate (model_p(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
34    allocate (model_t(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
36    allocate (model_qcw(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
37    allocate (model_qci(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
38    allocate (model_qrn(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
39    allocate (model_qsn(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
40    allocate (model_qgr(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
42    call da_convert_zk (iv%info(pseudo))
44 #ifdef A2C
45    call da_interp_lin_3d (grid%xb%u, iv%info(pseudo), model_u,'u')
46    call da_interp_lin_3d (grid%xb%v, iv%info(pseudo), model_v,'v')
47 #else
48    call da_interp_lin_3d (grid%xb%u, iv%info(pseudo), model_u)
49    call da_interp_lin_3d (grid%xb%v, iv%info(pseudo), model_v)
50 #endif
51    call da_interp_lin_3d (grid%xb%t, iv%info(pseudo), model_t)
52    call da_interp_lin_3d (grid%xb%p, iv%info(pseudo), model_p)
53    call da_interp_lin_3d (grid%xb%q, iv%info(pseudo), model_q)
55    call da_interp_lin_3d (grid%xb%qcw, iv%info(pseudo), model_qcw)
56    call da_interp_lin_3d (grid%xb%qci, iv%info(pseudo), model_qci)
57    call da_interp_lin_3d (grid%xb%qrn, iv%info(pseudo), model_qrn)
58    call da_interp_lin_3d (grid%xb%qsn, iv%info(pseudo), model_qsn)
59    call da_interp_lin_3d (grid%xb%qgr, iv%info(pseudo), model_qgr)
61    do n=iv%info(pseudo)%n1,iv%info(pseudo)%n2
62       !---------------------------------------------------------------
63       ! [3.0] Calculate observation O = B +(O-B):
64       !---------------------------------------------------------------
66       ! inv is from namelist for the first outer-loop
67       select case(pseudo_var(1:1))
68       case ('u', 'U')
69          if ( it > 1 ) then
70              iv % pseudo(n) % u % inv = ob%pseudo(n)%u - model_u(1,n)
71          else
72              ob % pseudo(n) % u = model_u(1,n) + iv % pseudo(n) % u % inv
73          end if
74       case ('v', 'V')
75          if ( it > 1 ) then
76              iv % pseudo(n) % v % inv = ob%pseudo(n)%v - model_v(1,n)
77          else
78              ob % pseudo(n) % v = model_v(1,n) + iv % pseudo(n) % v % inv
79          end if
80       case ('t', 'T')
81          if ( it > 1 ) then
82              iv % pseudo(n) % t % inv = ob%pseudo(n)%t - model_t(1,n)
83          else
84              ob % pseudo(n) % t = model_t(1,n) + iv % pseudo(n) % t % inv
85          end if
86       case ('p', 'P')
87          if ( it > 1 ) then
88              iv % pseudo(n) % p % inv = ob%pseudo(n)%p - model_p(1,n)
89          else
90              ob % pseudo(n) % p = model_p(1,n) + iv % pseudo(n) % p % inv
91          end if
92       case ('q', 'Q')
93          if ( len_trim(adjustl(pseudo_var)) == 1 ) then
94             if ( it > 1 ) then
95                 iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_q(1,n)
96             else
97                 ob % pseudo(n) % q = model_q(1,n) + iv % pseudo(n) % q % inv
98             end if
99          end if
100       end select
102       select case(pseudo_var(1:3))
103       case ('qcw', 'QCW')
104          if ( it > 1 ) then
105              iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_qcw(1,n)
106          else
107              ob % pseudo(n) % q = model_qcw(1,n) + iv % pseudo(n) % q % inv
108          endif
109       case ('qci', 'QCI')
110          if ( it > 1 ) then
111              iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_qci(1,n)
112          else
113              ob % pseudo(n) % q = model_qci(1,n) + iv % pseudo(n) % q % inv
114          endif
115       case ('qrn', 'QRN')
116          if ( it > 1 ) then
117              iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_qrn(1,n)
118          else
119              ob % pseudo(n) % q = model_qrn(1,n) + iv % pseudo(n) % q % inv
120          endif
121       case ('qsn', 'QSN')
122          if ( it > 1 ) then
123              iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_qsn(1,n)
124          else
125              ob % pseudo(n) % q = model_qsn(1,n) + iv % pseudo(n) % q % inv
126          endif
127       case ('qgr', 'QGR')
128          if ( it > 1 ) then
129              iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_qgr(1,n)
130          else
131              ob % pseudo(n) % q = model_qgr(1,n) + iv % pseudo(n) % q % inv
132          endif
133       end select
135    end do
137    deallocate (model_u)
138    deallocate (model_v)
139    deallocate (model_q)
140    deallocate (model_p)
141    deallocate (model_t)
143    deallocate (model_qcw)
144    deallocate (model_qci)
145    deallocate (model_qrn)
146    deallocate (model_qsn)
147    deallocate (model_qgr)
149    if (trace_use_dull) call da_trace_exit("da_get_innov_vector_pseudo")
151 end subroutine da_get_innov_vector_pseudo