1 subroutine da_get_innov_vector_pseudo(it, grid, ob, iv)
3 !-----------------------------------------------------------------------
5 !-----------------------------------------------------------------------
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))
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')
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)
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))
70 iv % pseudo(n) % u % inv = ob%pseudo(n)%u - model_u(1,n)
72 ob % pseudo(n) % u = model_u(1,n) + iv % pseudo(n) % u % inv
76 iv % pseudo(n) % v % inv = ob%pseudo(n)%v - model_v(1,n)
78 ob % pseudo(n) % v = model_v(1,n) + iv % pseudo(n) % v % inv
82 iv % pseudo(n) % t % inv = ob%pseudo(n)%t - model_t(1,n)
84 ob % pseudo(n) % t = model_t(1,n) + iv % pseudo(n) % t % inv
88 iv % pseudo(n) % p % inv = ob%pseudo(n)%p - model_p(1,n)
90 ob % pseudo(n) % p = model_p(1,n) + iv % pseudo(n) % p % inv
93 if ( len_trim(adjustl(pseudo_var)) == 1 ) then
95 iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_q(1,n)
97 ob % pseudo(n) % q = model_q(1,n) + iv % pseudo(n) % q % inv
102 select case(pseudo_var(1:3))
105 iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_qcw(1,n)
107 ob % pseudo(n) % q = model_qcw(1,n) + iv % pseudo(n) % q % inv
111 iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_qci(1,n)
113 ob % pseudo(n) % q = model_qci(1,n) + iv % pseudo(n) % q % inv
117 iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_qrn(1,n)
119 ob % pseudo(n) % q = model_qrn(1,n) + iv % pseudo(n) % q % inv
123 iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_qsn(1,n)
125 ob % pseudo(n) % q = model_qsn(1,n) + iv % pseudo(n) % q % inv
129 iv % pseudo(n) % q % inv = ob%pseudo(n)%q - model_qgr(1,n)
131 ob % pseudo(n) % q = model_qgr(1,n) + iv % pseudo(n) % q % inv
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