da_to_zk_new.inc

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