12 use missing_values_mod
19 type(oops_variables),
public :: obsvars
20 integer,
allocatable,
public :: obsvarindices(:)
22 type(oops_variables),
public :: geovars
23 integer :: nval, nlocs
24 real(kind_real),
allocatable :: wf(:)
25 integer,
allocatable :: wi(:)
26 character(len=MAXVARLEN),
public :: v_coord
27 character(len=MAXVARLEN),
public :: o_v_coord
28 logical,
public :: use_ln
43 use fckit_configuration_module,
only: fckit_configuration
46 type(fckit_configuration),
intent(in) :: grid_conf
48 character(kind=c_char,len=:),
allocatable :: coord_name
49 integer :: ivar, nvars
52 nvars = self%obsvars%nvars()
54 call self%geovars%push_back(self%obsvars%variable(ivar))
59 if( grid_conf%has(
"vertical coordinate") )
then
60 call grid_conf%get_or_die(
"vertical coordinate",coord_name)
61 self%v_coord = coord_name
62 if( (self%v_coord .eq.
var_prs) .or. (self%v_coord .eq.
var_prsi) )
then
73 if ( grid_conf%has(
"observation vertical coordinate") )
then
74 call grid_conf%get_or_die(
"observation vertical coordinate",coord_name)
75 self%o_v_coord = coord_name
77 self%o_v_coord = self%v_coord
89 type(c_ptr),
value,
intent(in) :: obss
91 real(kind_real),
allocatable :: obsvcoord(:)
94 real(kind_real),
allocatable :: tmp(:)
95 real(kind_real) :: tmp2
102 self%nval = vcoordprofile%nval
105 self%nlocs = obsspace_get_nlocs(obss)
106 allocate(obsvcoord(self%nlocs))
107 call obsspace_get_db(obss,
"MetaData", self%o_v_coord, obsvcoord)
110 allocate(self%wi(self%nlocs))
111 allocate(self%wf(self%nlocs))
114 allocate(tmp(vcoordprofile%nval))
115 do iobs = 1, self%nlocs
116 if (self%use_ln)
then
117 tmp = log(vcoordprofile%vals(:,iobs))
118 tmp2 = log(obsvcoord(iobs))
120 tmp = vcoordprofile%vals(:,iobs)
121 tmp2 = obsvcoord(iobs)
127 deallocate(obsvcoord)
138 integer,
intent(in) :: nvars, nlocs
139 real(c_double),
intent(inout) :: hofx(nvars, nlocs)
140 type(c_ptr),
value,
intent(in) :: obss
142 integer :: iobs, iobsvar, ivar
144 character(len=MAXVARLEN) :: geovar
146 do iobsvar = 1,
size(self%obsvarindices)
148 ivar = self%obsvarindices(iobsvar)
151 geovar = self%geovars%variable(iobsvar)
159 & hofx(ivar,iobs), self%wi(iobs), self%wf(iobs))
170 integer,
intent(in) :: nvars, nlocs
171 real(c_double),
intent(in) :: hofx(nvars, nlocs)
172 type(c_ptr),
value,
intent(in) :: obss
174 integer :: iobs, iobsvar, ivar
176 character(len=MAXVARLEN) :: geovar
177 real(c_double) :: missing
179 missing = missing_value(missing)
181 do iobsvar = 1,
size(self%obsvarindices)
183 ivar = self%obsvarindices(iobsvar)
186 geovar = self%geovars%variable(iobsvar)
192 do iobs = 1, self%nlocs
193 if (hofx(ivar,iobs) /= missing)
then
195 & hofx(ivar,iobs), self%wi(iobs), self%wf(iobs))
208 if (
allocated(self%wi))
deallocate(self%wi)
209 if (
allocated(self%wf))
deallocate(self%wf)
subroutine atmvertinterp_simobs_tl_(self, geovals, obss, nvars, nlocs, hofx)
subroutine atmvertinterp_tlad_setup_(self, grid_conf)
subroutine atmvertinterp_simobs_ad_(self, geovals, obss, nvars, nlocs, hofx)
subroutine atmvertinterp_tlad_cleanup_(self)
subroutine destructor(self)
subroutine atmvertinterp_tlad_settraj_(self, geovals, obss)
subroutine, public ufo_geovals_get_var(self, varname, geoval)
character(len=maxvarlen), parameter, public var_prsi
character(len=maxvarlen), parameter, public var_prs
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