da_to_zk_new.inc

References to this file elsewhere.
1 subroutine da_to_zk_new(obs_v, mdl_v, xp, v_interp_optn, zk,num,nlevels)
2 
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
6 
7    implicit none
8 
9    integer,          intent(in)  :: v_interp_optn
10    type(xpose_type), intent(in)  :: xp!Domain decomposition vars.
11    real,             intent(in)  :: obs_v(nlevels)
12    real,             intent(in)  :: mdl_v(xp%kms:xp%kme,num)
13    real,             intent(out) :: zk(nlevels,num)
14    integer,          intent(in)  :: num
15    integer,          intent(in)  :: nlevels
16 
17    integer                :: k,n,kk
18    real    :: r_ktsplus, r_kteminus
19 
20    if (trace_use) call da_trace_entry("da_to_zk_new")
21 
22    zk(:,:) = missing_r
23 
24    r_ktsplus  = real(xp%kts+1)
25    r_kteminus = real(xp%kte-1)
26 
27    if (v_interp_optn == v_interp_p) then
28       if (anal_type_verify) then
29          do n=1,num
30             do k=1,nlevels
31                if (obs_v(k) > mdl_v(xp%kts,n)) then
32                   ! below the lowest level:
33                   zk(k,n) = r_ktsplus - &
34                      (obs_v(k)-mdl_v(xp%kts+1,n))/(mdl_v(xp%kts,n)-mdl_v(xp%kts+1,n))
35                else if (obs_v(k) < mdl_v(xp%kts,n)) then
36                   ! above the highest level:
37                   zk(k,n) = r_kteminus + &
38                      (obs_v(k)-mdl_v(xp%kte-1,n))/(mdl_v(xp%kte,n)-mdl_v(xp%kte-1,n))
39                else
40                   do kk = xp%kts,xp%kte-1
41                      if (obs_v(k) <= mdl_v(kk,n) .and. obs_v(k) >= mdl_v(kk+1,n)) then
42                         zk(k,n) = real(kk) + (mdl_v(kk,n) - obs_v(k))/(mdl_v(kk,n) - mdl_v(kk+1,n))
43                         exit
44                      end if
45                   end do
46                end if
47             end do
48          end do
49       else
50          do n=1,num
51             do k=1,nlevels
52                do kk = xp%kts,xp%kte-1
53                   if (obs_v(k) <= mdl_v(kk,n) .and. obs_v(k) >= mdl_v(kk+1,n)) then
54                      zk(k,n) = real(kk) + (mdl_v(kk,n) - obs_v(k))/(mdl_v(kk,n) - mdl_v(kk+1,n))
55                      exit
56                   end if
57                end do
58             end do
59          end do
60       end if
61    else if(v_interp_optn == v_interp_h) then
62       if (anal_type_verify) then
63          do n=1,num
64             do k=1,nlevels
65                if (obs_v(k) < mdl_v(xp%kts,n)) then
66                    ! below the lowest level:
67                    zk(k,n) = r_ktsplus - &
68                      (obs_v(k)-mdl_v(xp%kts+1,n))/(mdl_v(xp%kts,n)-mdl_v(xp%kts+1,n))
69                else if (obs_v(k) > mdl_v(xp%kts,n)) then
70                   ! above the highest level:
71                   zk(k,n) = r_kteminus + &
72                      (obs_v(k)-mdl_v(xp%kte-1,n))/(mdl_v(xp%kte,n)-mdl_v(xp%kte-1,n))
73                else
74                   do kk = xp%kts,xp%kte-1
75                      if (obs_v(k) >= mdl_v(kk,n) .and. obs_v(k) <= mdl_v(kk+1,n)) then
76                         zk(k,n) = real(kk) + (mdl_v(kk,n) - obs_v(k))/(mdl_v(kk,n) - mdl_v(kk+1,n))
77                         exit
78                      end if
79                   end do
80                end if
81             end do
82          end do
83       else
84          do n=1,num
85             do k=1,nlevels
86                do kk = xp%kts,xp%kte-1
87                   if (obs_v(k) >= mdl_v(kk,n) .and. obs_v(k) <= mdl_v(kk+1,n)) then
88                      zk(k,n) = real(kk) + (mdl_v(kk,n) - obs_v(k))/(mdl_v(kk,n) - mdl_v(kk+1,n))
89                      exit
90                   end if
91                end do
92             end do
93          end do
94       end if
95    end if
96 
97    if (trace_use) call da_trace_exit("da_to_zk_new")
98 
99 end subroutine da_to_zk_new
100 
101