10 use oops_variables_mod
19 type(oops_variables),
public :: obsvars
20 type(oops_variables),
public :: geovars
21 character(len=MAXVARLEN),
public :: v_coord
36 use fckit_configuration_module,
only: fckit_configuration
40 type(fckit_configuration),
intent(in) :: yaml_conf
42 character(kind=c_char,len=:),
allocatable :: coord_name
46 if( yaml_conf%has(
"VertCoord") )
then
47 call yaml_conf%get_or_die(
"VertCoord",coord_name)
48 self%v_coord = coord_name
49 if( trim(self%v_coord) .ne.
var_z )
then
50 call abor1_ftn(
"ufo_radarradialvelocity: incorrect vertical coordinate specified")
56 call self%geovars%push_back(self%v_coord)
69 integer,
intent(in) :: nvars, nlocs
71 real(c_double),
intent(inout) :: hofx(nvars, nlocs)
72 type(c_ptr),
value,
intent(in) :: obss
75 integer :: iobs, ivar, nvars_geovars
76 real(kind_real),
dimension(:),
allocatable :: obsvcoord
77 real(kind_real),
dimension(:),
allocatable :: cosazm_costilt, sinazm_costilt, sintilt, vterminal
78 type(
ufo_geoval),
pointer :: vcoordprofile, profile
79 real(kind_real),
allocatable :: wf(:)
80 integer,
allocatable :: wi(:)
82 character(len=MAXVARLEN) :: geovar
84 real(kind_real),
allocatable :: tmp(:)
85 real(kind_real) :: tmp2
87 real(kind_real),
allocatable :: vfields(:,:)
95 allocate(obsvcoord(nlocs))
96 allocate(cosazm_costilt(nlocs))
97 allocate(sinazm_costilt(nlocs))
98 allocate(sintilt(nlocs))
99 allocate(vterminal(nlocs))
101 call obsspace_get_db(obss,
"MetaData",
"geometric_height", obsvcoord)
102 call obsspace_get_db(obss,
"MetaData",
"cosazm_costilt", cosazm_costilt)
103 call obsspace_get_db(obss,
"MetaData",
"sinazm_costilt", sinazm_costilt)
104 call obsspace_get_db(obss,
"MetaData",
"sintilt", sintilt)
114 allocate(tmp(vcoordprofile%nval))
116 tmp = vcoordprofile%vals(:,iobs)
117 tmp2 = obsvcoord(iobs)
122 nvars_geovars = self%geovars%nvars() - 1
123 allocate(vfields(nvars_geovars,nlocs))
125 do ivar = 1, nvars_geovars
127 geovar = self%geovars%variable(ivar)
135 & vfields(ivar,iobs), wi(iobs), wf(iobs))
142 hofx(ivar,iobs) = vfields(1,iobs)*cosazm_costilt(iobs) &
143 + vfields(2,iobs)*sinazm_costilt(iobs) &
144 + (vfields(3,iobs)-vterminal(iobs))*sintilt(iobs)
149 deallocate(obsvcoord)
150 deallocate(cosazm_costilt)
151 deallocate(sinazm_costilt)
153 deallocate(vterminal)
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Fortran module for radarradialvelocity observation operator.
character(len=maxvarlen), dimension(3), parameter geovars_default
subroutine ufo_radarradialvelocity_setup(self, yaml_conf)
subroutine ufo_radarradialvelocity_simobs(self, geovals, obss, nvars, nlocs, hofx)
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_weights(nlev, obl, vec, wi, wf)
subroutine vert_interp_apply(nlev, fvec, f, wi, wf)
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators
Fortran derived type for the observation type.