da_convertor_v_interp.inc
References to this file elsewhere.
1 subroutine da_convertor_v_interp( v_in, z_in, n_in, &
2 v_out, z_out, n_out )
3
4 implicit none
5
6 integer, intent(in) :: n_in, n_out
7 real, dimension(n_in), intent(in) :: v_in, z_in
8 real, dimension(n_out), intent(in) :: z_out
9 real, dimension(n_out), intent(out) :: v_out
10
11 integer :: k, kabv, kblw, i
12 real :: w1, w2
13
14 logical :: increase_in_vertical
15
16 ! does vertical coordinate increase or decrease with increasing k?
17 ! set offset appropriately
18
19 if (z_in(n_in) > z_in(1)) then
20 increase_in_vertical = .true.
21 else
22 increase_in_vertical = .false.
23 end if
24
25 if (increase_in_vertical) then
26 do k=1, n_out
27 if (z_out(k) <= z_in(1)) then
28 kabv = 2
29 else if (z_out(k) >= z_in(n_in)) then
30 kabv = n_in
31 else
32 i_loop: do i=2, n_in
33 if (z_out(k) <= z_in(i)) then
34 kabv = i
35 exit i_loop
36 end if
37 end do i_loop
38 end if
39
40 kblw = kabv - 1
41 w2 = (z_in(kabv)-z_out(k))/(z_in(kabv)-z_in(kblw))
42 w1 = 1.0-w2
43 v_out(k) = w1*v_in(kabv) + w2*v_in(kblw)
44 end do
45 else
46 do k=1, n_out
47 if (z_out(k) >= z_in(1)) then
48 kabv = 2
49 else if (z_out(k) <= z_in(n_in)) then
50 kabv = n_in
51 else
52 d_loop: do i=2, n_in
53 if (z_out(k) >= z_in(i)) then
54 kabv = i
55 exit d_loop
56 end if
57 end do d_loop
58 end if
59
60 kblw = kabv - 1
61 w2 = (z_in(kabv)-z_out(k))/(z_in(kabv)-z_in(kblw))
62 w1 = 1.0-w2
63 v_out(k) = w1*v_in(kabv) + w2*v_in(kblw)
64 end do
65 end if
66
67 end subroutine da_convertor_v_interp
68