da_get_innov_vector_pseudo.inc
References to this file elsewhere.
1 subroutine da_get_innov_vector_pseudo( xb, xp, ob, iv)
2
3 !-----------------------------------------------------------------------
4 ! Purpose: TBD
5 !-----------------------------------------------------------------------
6
7 implicit none
8
9 type(xb_type), intent(in) :: xb ! Background structure
10 type(xpose_type), intent(in):: xp ! Domain decomposition vars.
11 type( y_type), intent(inout):: ob ! Observation structure.
12 type(ob_type), intent(inout):: iv ! O-B structure.
13
14 integer :: n ! Loop counter.
15 integer :: i, j ! Index dimension.
16 real :: dx, dxm ! Interpolation weights.
17 real :: dy, dym ! Interpolation weights.
18 real :: model_u ! Model value u at oblocation.
19 real :: model_v ! Model value v at oblocation.
20 real :: model_t ! Model value t at oblocation.
21 real :: model_p ! Model value p at oblocation.
22 real :: model_q ! Model value q at oblocation.
23
24 if (iv % num_pseudo > 0) then
25 do n=1, iv % num_pseudo
26 ! [1.1] Get horizontal interpolation weights:
27
28 i = iv%pseudo(n)%loc%i
29 dy = iv%pseudo(n)%loc%dy
30 dym = iv%pseudo(n)%loc%dym
31 j = iv%pseudo(n)%loc%j
32 dx = iv%pseudo(n)%loc%dx
33 dxm = iv%pseudo(n)%loc%dxm
34
35 write(unit=stdout, fmt='(/a,2i5/a,4f15.5/)') &
36 'pseudo: i,j=', i,j, &
37 'pseudo: dy, dym, dx, dxm=', dy, dym, dx, dxm
38
39 ! [1.2] Interpolate horizontally:
40 call da_interp_obs_lin_2d( xb % u, xp, i, j, dx, dy, dxm, dym, &
41 model_u, iv%pseudo(n)%zk)
42 call da_interp_obs_lin_2d( xb % v, xp, i, j, dx, dy, dxm, dym, &
43 model_v, iv%pseudo(n)%zk)
44 call da_interp_obs_lin_2d( xb % t, xp, i, j, dx, dy, dxm, dym, &
45 model_t, iv%pseudo(n)%zk)
46 call da_interp_obs_lin_2d( xb % p, xp, i, j, dx, dy, dxm, dym, &
47 model_p, iv%pseudo(n)%zk)
48 call da_interp_obs_lin_2d( xb % q, xp, i, j, dx, dy, dxm, dym, &
49 model_q, iv%pseudo(n)%zk)
50
51 !---------------------------------------------------------------
52 ! [3.0] Calculate observation O = B +(O-B):
53 !---------------------------------------------------------------
54
55 if (pseudo_var(1:1) == 'u' .or. pseudo_var(1:1) == 'U') then
56 ob % pseudo(n) % u = model_u + iv % pseudo(n) % u % inv
57 else if (pseudo_var(1:1) == 'v' .or. pseudo_var(1:1) == 'V') then
58 ob % pseudo(n) % v = model_v + iv % pseudo(n) % v % inv
59 else if (pseudo_var(1:1) == 't' .or. pseudo_var(1:1) == 'T') then
60 ob % pseudo(n) % t = model_t + iv % pseudo(n) % t % inv
61 else if (pseudo_var(1:1) == 'p' .or. pseudo_var(1:1) == 'P') then
62 ob % pseudo(n) % p = model_p + iv % pseudo(n) % p % inv
63 else if (pseudo_var(1:1) == 'q' .or. pseudo_var(1:1) == 'Q') then
64 ob % pseudo(n) % q = model_q + iv % pseudo(n) % q % inv
65 end if
66 end do
67 end if
68
69 end subroutine da_get_innov_vector_pseudo
70
71