UFO
ufo_radarradialvelocity_tlad_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2017-2019 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 
7 
8  use oops_variables_mod
9  use ufo_vars_mod
10  use ufo_geovals_mod
11  use vert_interp_mod
12  use missing_values_mod
13 
14 
15 ! ------------------------------------------------------------------------------
16 
18  private
19  type(oops_variables), public :: obsvars
20  type(oops_variables), public :: geovars
21  integer :: nval, nlocs
22  character(len=MAXVARLEN), public :: v_coord ! GeoVaL to use to interpolate in vertical
23  real(kind_real), allocatable :: wf(:)
24  integer, allocatable :: wi(:)
25  real(kind_real), dimension(:), allocatable :: cosazm_costilt, sinazm_costilt, sintilt, vterminal
26  contains
27  procedure :: setup => radarradialvelocity_tlad_setup_
28  procedure :: cleanup => radarradialvelocity_tlad_cleanup_
29  procedure :: settraj => radarradialvelocity_tlad_settraj_
30  procedure :: simobs_tl => radarradialvelocity_simobs_tl_
31  procedure :: simobs_ad => radarradialvelocity_simobs_ad_
32  final :: destructor
34 
35  character(len=maxvarlen), dimension(3), parameter :: geovars_default = (/var_u, &
36  var_v, &
37  var_w /)
38 
39 ! ------------------------------------------------------------------------------
40 contains
41 ! ------------------------------------------------------------------------------
42 
43 subroutine radarradialvelocity_tlad_setup_(self, yaml_conf)
44  use fckit_configuration_module, only: fckit_configuration
45  implicit none
46  class(ufo_radarradialvelocity_tlad), intent(inout) :: self
47  type(fckit_configuration), intent(in) :: yaml_conf
48 
49  character(kind=c_char,len=:), allocatable :: coord_name
50 
51  integer :: ivar, nvars
52  call self%geovars%push_back(geovars_default)
53 
54  if( yaml_conf%has("VertCoord") ) then
55  call yaml_conf%get_or_die("VertCoord",coord_name)
56  self%v_coord = coord_name
57  if( trim(self%v_coord) .ne. var_z ) then
58  call abor1_ftn("ufo_radarradialvelocity: incorrect vertical coordinate specified")
59  endif
60  else ! default
61  !self%v_coord = var_z
62  self%v_coord = var_zm
63  endif
64 
65 
67 
68 ! ------------------------------------------------------------------------------
69 
70 subroutine radarradialvelocity_tlad_settraj_(self, geovals, obss)
71  use obsspace_mod
72  implicit none
73  class(ufo_radarradialvelocity_tlad), intent(inout) :: self
74  type(ufo_geovals), intent(in) :: geovals
75  type(c_ptr), value, intent(in) :: obss
76 
77  !local variables
78  integer :: iobs, ivar, nvars_geovars
79  real(kind_real), dimension(:), allocatable :: obsvcoord
80  !real(kind_real), dimension(:), allocatable :: cosazm_costilt, sinazm_costilt, sintilt, vterminal
81  type(ufo_geoval), pointer :: vcoordprofile, profile
82 
83  real(kind_real), allocatable :: tmp(:)
84  real(kind_real) :: tmp2
85  integer:: i
86 
87  real(kind_real), allocatable :: vfields(:,:) ! background fields (u,v,w) interplated vertically to the observation height
88 
89  ! Make sure nothing already allocated
90  call self%cleanup()
91 
92  ! Get height profiles from geovals
93  call ufo_geovals_get_var(geovals, self%v_coord, vcoordprofile)
94  self%nval = vcoordprofile%nval
95 
96  ! Get the observation vertical coordinates
97  self%nlocs = obsspace_get_nlocs(obss)
98  allocate(obsvcoord(self%nlocs))
99  allocate(self%cosazm_costilt(self%nlocs))
100  allocate(self%sinazm_costilt(self%nlocs))
101  allocate(self%sintilt(self%nlocs))
102  allocate(self%vterminal(self%nlocs))
103 
104  call obsspace_get_db(obss, "MetaData", "geometric_height", obsvcoord)
105  call obsspace_get_db(obss, "MetaData", "cosazm_costilt", self%cosazm_costilt)
106  call obsspace_get_db(obss, "MetaData", "sinazm_costilt", self%sinazm_costilt)
107  call obsspace_get_db(obss, "MetaData", "sintilt", self%sintilt)
108 ! call obsspace_get_db(obss, "MetaData", "vterminal", self%vterminal)
109 
110  ! Allocate arrays for interpolation weights
111  allocate(self%wi(self%nlocs))
112  allocate(self%wf(self%nlocs))
113 
114  ! Calculate the interpolation weights
115  allocate(tmp(vcoordprofile%nval))
116  do iobs = 1, self%nlocs
117  tmp = vcoordprofile%vals(:,iobs)
118  tmp2 = obsvcoord(iobs)
119  call vert_interp_weights(vcoordprofile%nval, tmp2, tmp, self%wi(iobs), self%wf(iobs))
120  enddo
121 
122 ! Cleanup memory
123  deallocate(obsvcoord)
124  deallocate(tmp)
125 
127 
128 ! ------------------------------------------------------------------------------
129 
130 subroutine radarradialvelocity_simobs_tl_(self, geovals, obss, nvars, nlocs, hofx)
131  implicit none
132  class(ufo_radarradialvelocity_tlad), intent(in) :: self
133  type(ufo_geovals), intent(in) :: geovals
134  integer, intent(in) :: nvars, nlocs
135  real(c_double), intent(inout) :: hofx(nvars, nlocs)
136  type(c_ptr), value, intent(in) :: obss
137 
138  integer :: iobs, ivar, nvars_geovars
139  type(ufo_geoval), pointer :: profile
140  character(len=MAXVARLEN) :: geovar
141 
142  real(kind_real), allocatable :: vfields(:,:) ! background fields (u,v,w) interplated vertically to the observation height
143 
144 
145 ! Number of variables in geovars (without the vertical coordinate)
146  nvars_geovars = self%geovars%nvars() - 1
147  allocate(vfields(nvars_geovars,nlocs))
148  vfields=0.0
149 
150  do ivar = 1, nvars_geovars
151  geovar = self%geovars%variable(ivar)
152 
153 ! Get profile for this variable from geovals
154  call ufo_geovals_get_var(geovals, geovar, profile)
155 
156 ! Interpolate from geovals to observational location into hofx
157  do iobs = 1, nlocs
158  call vert_interp_apply_tl(profile%nval, profile%vals(:,iobs), &
159  & vfields(ivar,iobs), self%wi(iobs), self%wf(iobs))
160  enddo
161  enddo
162 
163  do ivar = 1, nvars
164  do iobs=1,nlocs
165  hofx(ivar,iobs) = vfields(1,iobs)*self%cosazm_costilt(iobs) &
166  + vfields(2,iobs)*self%sinazm_costilt(iobs)
167  enddo
168  end do
169 
170  deallocate(vfields)
171 
172 end subroutine radarradialvelocity_simobs_tl_
173 
174 ! ------------------------------------------------------------------------------
175 
176 subroutine radarradialvelocity_simobs_ad_(self, geovals, obss, nvars, nlocs, hofx)
177  implicit none
178  class(ufo_radarradialvelocity_tlad), intent(in) :: self
179  type(ufo_geovals), intent(inout) :: geovals
180  integer, intent(in) :: nvars, nlocs
181  real(c_double), intent(in) :: hofx(nvars, nlocs)
182  type(c_ptr), value, intent(in) :: obss
183 
184  integer :: iobs, ivar, nvars_geovars
185  type(ufo_geoval), pointer :: profile
186  character(len=MAXVARLEN) :: geovar
187  real(c_double) :: missing
188  real(kind_real), allocatable :: vfields(:,:) ! background fields (u,v,w) interplated vertically to the observation height
189 
190  missing = missing_value(missing)
191 
192 ! Number of variables in geovars (without the vertical coordinate)
193  nvars_geovars = self%geovars%nvars() - 1
194  allocate(vfields(nvars_geovars,nlocs))
195 
196  vfields=0.0
197  do ivar = 1, nvars
198  do iobs=1,nlocs
199  ! no vertical velocity and terminal velocity in GSI rw observer, it can add
200  ! in future after acceptance test
201  vfields(1,iobs) = vfields(1,iobs) + hofx(ivar,iobs)*self%cosazm_costilt(iobs)
202  vfields(2,iobs) = vfields(2,iobs) + hofx(ivar,iobs)*self%sinazm_costilt(iobs)
203  enddo
204  end do
205 
206  do ivar = 1, nvars_geovars
207 ! Get the name of input variable in geovals
208  geovar = self%geovars%variable(ivar)
209 
210 ! Get profile for this variable from geovals
211  call ufo_geovals_get_var(geovals, geovar, profile)
212 
213 ! Interpolate from geovals to observational location into hofx
214  do iobs = 1, nlocs
215  call vert_interp_apply_ad(profile%nval, profile%vals(:,iobs), &
216  & vfields(ivar,iobs), self%wi(iobs), self%wf(iobs))
217  enddo
218  enddo
219 
220 
221  deallocate(vfields)
222 end subroutine radarradialvelocity_simobs_ad_
223 
224 ! ------------------------------------------------------------------------------
225 
227  implicit none
228  class(ufo_radarradialvelocity_tlad), intent(inout) :: self
229  self%nval = 0
230  self%nlocs = 0
231  if (allocated(self%wi)) deallocate(self%wi)
232  if (allocated(self%wf)) deallocate(self%wf)
233  if (allocated(self%cosazm_costilt)) deallocate(self%cosazm_costilt)
234  if (allocated(self%sinazm_costilt)) deallocate(self%sinazm_costilt)
235  if (allocated(self%sintilt)) deallocate(self%sintilt)
236  if (allocated(self%vterminal)) deallocate(self%vterminal)
238 
239 ! ------------------------------------------------------------------------------
240 
241 subroutine destructor(self)
242  type(ufo_radarradialvelocity_tlad), intent(inout) :: self
243 
244  call self%cleanup()
245 
246 end subroutine destructor
247 
248 ! ------------------------------------------------------------------------------
249 
subroutine, public ufo_geovals_get_var(self, varname, geoval)
subroutine radarradialvelocity_tlad_settraj_(self, geovals, obss)
subroutine radarradialvelocity_simobs_tl_(self, geovals, obss, nvars, nlocs, hofx)
subroutine radarradialvelocity_tlad_setup_(self, yaml_conf)
subroutine radarradialvelocity_simobs_ad_(self, geovals, obss, nvars, nlocs, hofx)
character(len=maxvarlen), dimension(3), parameter geovars_default
character(len=maxvarlen), parameter, public var_zm
character(len=maxvarlen), parameter, public var_v
character(len=maxvarlen), parameter, public var_z
character(len=maxvarlen), parameter, public var_u
character(len=maxvarlen), parameter, public var_w
Fortran module to perform linear interpolation.
Definition: vert_interp.F90:8
subroutine vert_interp_apply_ad(nlev, fvec_ad, f_ad, wi, wf)
subroutine vert_interp_weights(nlev, obl, vec, wi, wf)
Definition: vert_interp.F90:22
subroutine vert_interp_apply_tl(nlev, fvec_tl, f_tl, wi, wf)
Definition: vert_interp.F90:92
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators