FV3 Bundle
vert_interp.F90
Go to the documentation of this file.
1 ! (C) Copyright 2018 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
6 !> Fortran module to perform linear interpolation
7 
9 
10 use kinds, only: kind_real
11 
12 implicit none
13 public
14 
15 contains
16 
17 ! ------------------------------------------------------------------------------
18 
19 subroutine vert_interp_weights(nlev,obl,vec,wi,wf)
20 
21 implicit none
22 integer, intent(in ) :: nlev !Number of model levels
23 real(kind_real), intent(in ) :: obl !Observation location
24 real(kind_real), intent(in ) :: vec(nlev) !Structured vector of grid points
25 integer, intent(out) :: wi !Index for interpolation
26 real(kind_real), intent(out) :: wf !Weight for interpolation
27 
28 integer :: k
29 
30 if (vec(1) < vec(nlev)) then !Pressure increases with index
31 
32  if (obl < vec(1)) then
33  wi = 1
34  elseif (obl > vec(nlev)) then
35  wi = nlev - 1
36  else
37  do k = 1,nlev-1
38  if (obl >= vec(k) .and. obl <= vec(k+1)) then
39  wi = k
40  endif
41  enddo
42  endif
43 
44 else !Pressure decreases with index
45 
46  if (obl > vec(1)) then
47  wi = 1
48  elseif (obl < vec(nlev)) then
49  wi = nlev - 1
50  else
51  do k = 1,nlev-1
52  if (obl >= vec(k+1) .and. obl <= vec(k)) then
53  wi = k
54  endif
55  enddo
56  endif
57 
58 endif
59 
60 wf = (vec(wi+1) - obl)/(vec(wi+1) - vec(wi))
61 
62 end subroutine vert_interp_weights
63 
64 ! ------------------------------------------------------------------------------
65 
66 subroutine vert_interp_apply(nlev, fvec, f, wi, wf)
67 
68 implicit none
69 integer, intent(in ) :: nlev !Number of model levels
70 real(kind_real), intent(in ) :: fvec(nlev) !Field at grid points
71 integer, intent(in ) :: wi !Index for interpolation
72 real(kind_real), intent(in ) :: wf !Weight for interpolation
73 real(kind_real), intent(out) :: f !Output at obs location using linear interp
74 
75 f = fvec(wi)*wf + fvec(wi+1)*(1.0-wf)
76 
77 end subroutine vert_interp_apply
78 
79 ! ------------------------------------------------------------------------------
80 
81 subroutine vert_interp_apply_tl(nlev, fvec_tl, f_tl, wi, wf)
82 
83 implicit none
84 integer, intent(in) :: nlev
85 real(kind_real), intent(in) :: fvec_tl(nlev)
86 integer, intent(in) :: wi
87 real(kind_real), intent(in) :: wf
88 real(kind_real), intent(out) :: f_tl
89 
90 f_tl = fvec_tl(wi)*wf + fvec_tl(wi+1)*(1.0_kind_real-wf)
91 
92 end subroutine vert_interp_apply_tl
93 
94 ! ------------------------------------------------------------------------------
95 
96 subroutine vert_interp_apply_ad(nlev, fvec_ad, f_ad, wi, wf)
97 
98 implicit none
99 integer, intent(in) :: nlev
100 real(kind_real), intent(inout) :: fvec_ad(nlev)
101 integer, intent(in) :: wi
102 real(kind_real), intent(in) :: wf
103 real(kind_real), intent(in) :: f_ad
104 
105 fvec_ad(wi ) = fvec_ad(wi ) + f_ad*wf
106 fvec_ad(wi+1) = fvec_ad(wi+1) + f_ad*(1.0_kind_real-wf)
107 
108 end subroutine vert_interp_apply_ad
109 
110 ! ------------------------------------------------------------------------------
111 
112 end module vert_interp_mod
subroutine vert_interp_apply(nlev, fvec, f, wi, wf)
Definition: vert_interp.F90:67
subroutine vert_interp_apply_ad(nlev, fvec_ad, f_ad, wi, wf)
Definition: vert_interp.F90:97
subroutine vert_interp_apply_tl(nlev, fvec_tl, f_tl, wi, wf)
Definition: vert_interp.F90:82
Fortran module to perform linear interpolation.
Definition: vert_interp.F90:8
integer, parameter, public kind_real
subroutine vert_interp_weights(nlev, obl, vec, wi, wf)
Definition: vert_interp.F90:20