12 use missing_values_mod
19 type(oops_variables),
public :: obsvars
20 type(oops_variables),
public :: geovars
21 integer :: nval, nlocs
22 character(len=MAXVARLEN),
public :: v_coord
23 real(kind_real),
allocatable :: wf(:)
24 integer,
allocatable :: wi(:)
25 real(kind_real),
dimension(:),
allocatable :: cosazm_costilt, sinazm_costilt, sintilt, vterminal
44 use fckit_configuration_module,
only: fckit_configuration
47 type(fckit_configuration),
intent(in) :: yaml_conf
49 character(kind=c_char,len=:),
allocatable :: coord_name
51 integer :: ivar, nvars
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")
75 type(c_ptr),
value,
intent(in) :: obss
78 integer :: iobs, ivar, nvars_geovars
79 real(kind_real),
dimension(:),
allocatable :: obsvcoord
81 type(
ufo_geoval),
pointer :: vcoordprofile, profile
83 real(kind_real),
allocatable :: tmp(:)
84 real(kind_real) :: tmp2
87 real(kind_real),
allocatable :: vfields(:,:)
94 self%nval = vcoordprofile%nval
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))
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)
111 allocate(self%wi(self%nlocs))
112 allocate(self%wf(self%nlocs))
115 allocate(tmp(vcoordprofile%nval))
116 do iobs = 1, self%nlocs
117 tmp = vcoordprofile%vals(:,iobs)
118 tmp2 = obsvcoord(iobs)
123 deallocate(obsvcoord)
134 integer,
intent(in) :: nvars, nlocs
135 real(c_double),
intent(inout) :: hofx(nvars, nlocs)
136 type(c_ptr),
value,
intent(in) :: obss
138 integer :: iobs, ivar, nvars_geovars
140 character(len=MAXVARLEN) :: geovar
142 real(kind_real),
allocatable :: vfields(:,:)
146 nvars_geovars = self%geovars%nvars() - 1
147 allocate(vfields(nvars_geovars,nlocs))
150 do ivar = 1, nvars_geovars
151 geovar = self%geovars%variable(ivar)
159 & vfields(ivar,iobs), self%wi(iobs), self%wf(iobs))
165 hofx(ivar,iobs) = vfields(1,iobs)*self%cosazm_costilt(iobs) &
166 + vfields(2,iobs)*self%sinazm_costilt(iobs)
180 integer,
intent(in) :: nvars, nlocs
181 real(c_double),
intent(in) :: hofx(nvars, nlocs)
182 type(c_ptr),
value,
intent(in) :: obss
184 integer :: iobs, ivar, nvars_geovars
186 character(len=MAXVARLEN) :: geovar
187 real(c_double) :: missing
188 real(kind_real),
allocatable :: vfields(:,:)
190 missing = missing_value(missing)
193 nvars_geovars = self%geovars%nvars() - 1
194 allocate(vfields(nvars_geovars,nlocs))
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)
206 do ivar = 1, nvars_geovars
208 geovar = self%geovars%variable(ivar)
216 & vfields(ivar,iobs), self%wi(iobs), self%wf(iobs))
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)
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
subroutine radarradialvelocity_tlad_cleanup_(self)
subroutine destructor(self)
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.
subroutine vert_interp_apply_ad(nlev, fvec_ad, f_ad, wi, wf)
subroutine vert_interp_weights(nlev, obl, vec, wi, wf)
subroutine vert_interp_apply_tl(nlev, fvec_tl, f_tl, wi, wf)
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators