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 < 1) return
25    
26    if (trace_use) call da_trace_entry("da_get_innov_vector_pseudo")
27 
28    do n=1, iv % num_pseudo
29       ! [1.1] Get horizontal interpolation weights:
30 
31       i = iv%pseudo(n)%loc%i
32       dy = iv%pseudo(n)%loc%dy
33       dym = iv%pseudo(n)%loc%dym
34       j = iv%pseudo(n)%loc%j
35       dx = iv%pseudo(n)%loc%dx
36       dxm = iv%pseudo(n)%loc%dxm
37 
38       write(unit=stdout, fmt='(/a,2i5/a,4f15.5/)') &
39            'pseudo: i,j=', i,j, &
40            'pseudo: dy, dym, dx, dxm=', dy, dym, dx, dxm
41 
42       ! [1.2] Interpolate horizontally:
43       call da_interp_obs_lin_2d( xb % u, xp, i, j, dx, dy, dxm, dym, &
44                               model_u, iv%pseudo(n)%zk)
45       call da_interp_obs_lin_2d( xb % v, xp, i, j, dx, dy, dxm, dym, &
46                               model_v, iv%pseudo(n)%zk)
47       call da_interp_obs_lin_2d( xb % t, xp, i, j, dx, dy, dxm, dym, &
48                               model_t, iv%pseudo(n)%zk)
49       call da_interp_obs_lin_2d( xb % p, xp, i, j, dx, dy, dxm, dym, &
50                               model_p, iv%pseudo(n)%zk)
51       call da_interp_obs_lin_2d( xb % q, xp, i, j, dx, dy, dxm, dym, &
52                               model_q, iv%pseudo(n)%zk)
53 
54       !---------------------------------------------------------------
55       ! [3.0] Calculate observation O = B +(O-B):
56       !---------------------------------------------------------------
57 
58       if (pseudo_var(1:1) == 'u' .or.  pseudo_var(1:1) == 'U') then
59          ob % pseudo(n) % u = model_u + iv % pseudo(n) % u % inv
60       else if (pseudo_var(1:1) == 'v' .or.  pseudo_var(1:1) == 'V') then
61          ob % pseudo(n) % v = model_v + iv % pseudo(n) % v % inv
62       else if (pseudo_var(1:1) == 't' .or.  pseudo_var(1:1) == 'T') then
63          ob % pseudo(n) % t = model_t + iv % pseudo(n) % t % inv
64       else if (pseudo_var(1:1) == 'p' .or.  pseudo_var(1:1) == 'P') then
65          ob % pseudo(n) % p = model_p + iv % pseudo(n) % p % inv
66       else if (pseudo_var(1:1) == 'q' .or.  pseudo_var(1:1) == 'Q') then
67          ob % pseudo(n) % q = model_q + iv % pseudo(n) % q % inv
68       end if 
69    end do
70    
71    if (trace_use) call da_trace_exit("da_get_innov_vector_pseudo")
72 
73 end subroutine da_get_innov_vector_pseudo
74 
75