20 use missing_values_mod
21 use fckit_log_module,
only : fckit_log
23 integer,
parameter :: max_string=800
28 integer :: nval, nlocs
29 real(kind_real),
allocatable :: prs(:,:), t(:,:), q(:,:), gph(:,:), gph_sfc(:,:)
46 type(c_ptr),
value,
intent(in) :: obss
47 character(len=*),
parameter :: myname_=
"ufo_gnssro_bndropp1d_tlad_settraj"
48 character(max_string) :: err_msg
49 type(
ufo_geoval),
pointer :: t, q, prs, gph, gph_sfc
52 write(err_msg,*)
"TRACE: ufo_gnssro_bndropp1d_tlad_settraj: begin"
53 call fckit_log%info(err_msg)
66 self%nlocs = obsspace_get_nlocs(obss)
68 allocate(self%t(self%nval,self%nlocs))
69 allocate(self%q(self%nval,self%nlocs))
70 allocate(self%prs(self%nval,self%nlocs))
71 allocate(self%gph(self%nval,self%nlocs))
72 allocate(self%gph_sfc(1,self%nlocs))
79 self%gph_sfc = gph_sfc%vals
91 real(kind_real),
intent(inout) :: hofx(:)
92 type(c_ptr),
value,
intent(in) :: obss
94 integer :: iobs,nlev, nlocs
96 character(len=*),
parameter :: myname_=
"ufo_gnssro_bndropp1d_simobs_tl"
97 character(max_string) :: err_msg
101 real(kind_real),
allocatable :: gph_d_zero(:)
102 real(kind_real) :: gph_sfc_d_zero
103 real(kind_real),
allocatable :: obslat(:), obslon(:), obsimpp(:), obslocr(:), obsgeoid(:)
106 write(err_msg,*)
"TRACE: ufo_gnssro_bndropp1d_simobs_tl: begin"
107 call fckit_log%info(err_msg)
110 if (.not. self%ltraj)
then
111 write(err_msg,*) myname_,
' trajectory wasnt set!'
112 call abor1_ftn(err_msg)
116 if (geovals%nlocs /=
size(hofx))
then
117 write(err_msg,*) myname_,
' error: nlocs inconsistent!'
118 call abor1_ftn(err_msg)
129 allocate(gph_d_zero(nlev))
134 allocate(obslon(nlocs))
135 allocate(obslat(nlocs))
136 allocate(obsimpp(nlocs))
137 allocate(obslocr(nlocs))
138 allocate(obsgeoid(nlocs))
139 call obsspace_get_db(obss,
"MetaData",
"longitude", obslon)
140 call obsspace_get_db(obss,
"MetaData",
"latitude", obslat)
141 call obsspace_get_db(obss,
"MetaData",
"impact_parameter", obsimpp)
142 call obsspace_get_db(obss,
"MetaData",
"earth_radius_of_curvature", obslocr)
143 call obsspace_get_db(obss,
"MetaData",
"geoid_height_above_reference_ellipsoid", obsgeoid)
152 write(err_msg,*)
"TRACE: ufo_gnssro_bndropp1d_simobs_tl: complete"
153 call fckit_log%info(err_msg)
165 real(kind_real),
intent(in) :: hofx(:)
166 type(c_ptr),
value,
intent(in) :: obss
167 real(c_double) :: missing
171 real(kind_real),
parameter :: gph_sfc_d_zero = 0.0
172 real(kind_real),
allocatable :: gph_d_zero(:)
174 real(kind_real),
allocatable :: obslat(:), obslon(:), obsimpp(:), obslocr(:), obsgeoid(:)
175 integer :: iobs,nlev, nlocs
176 character(len=*),
parameter :: myname_=
"ufo_gnssro_bndropp1d_simobs_ad"
177 character(max_string) :: err_msg
179 write(err_msg,*)
"TRACE: ufo_gnssro_bndropp1d_simobs_ad: begin"
180 call fckit_log%info(err_msg)
183 if (.not. self%ltraj)
then
184 write(err_msg,*) myname_,
' trajectory wasnt set!'
185 call abor1_ftn(err_msg)
188 if (geovals%nlocs /=
size(hofx))
then
189 write(err_msg,*) myname_,
' error: nlocs inconsistent!'
190 call abor1_ftn(err_msg)
201 allocate(gph_d_zero(nlev))
205 allocate(obslon(nlocs))
206 allocate(obslat(nlocs))
207 allocate(obsimpp(nlocs))
208 allocate(obslocr(nlocs))
209 allocate(obsgeoid(nlocs))
211 call obsspace_get_db(obss,
"MetaData",
"longitude", obslon)
212 call obsspace_get_db(obss,
"MetaData",
"latitude", obslat)
213 call obsspace_get_db(obss,
"MetaData",
"impact_parameter", obsimpp)
214 call obsspace_get_db(obss,
"MetaData",
"earth_radius_of_curvature", obslocr)
215 call obsspace_get_db(obss,
"MetaData",
"geoid_height_above_reference_ellipsoid", obsgeoid)
217 missing = missing_value(missing)
225 deallocate(gph_d_zero)
227 write(err_msg,*)
"TRACE: ufo_gnssro_bndropp1d_simobs_ad: complete"
228 call fckit_log%info(err_msg)
240 character(len=*),
parameter :: myname_=
"ufo_gnssro_bndropp_tlad_delete"
243 if (
allocated(self%prs))
deallocate(self%prs)
244 if (
allocated(self%t))
deallocate(self%t)
245 if (
allocated(self%q))
deallocate(self%q)
246 if (
allocated(self%gph))
deallocate(self%gph)
247 if (
allocated(self%gph_sfc))
deallocate(self%gph_sfc)
type(registry_t), public ufo_geovals_registry
Linked list interface - defines registry_t type.
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Fortran module for gnssro bending angle ropp1d tangent linear and adjoint following the ROPP (2018 Au...
subroutine ufo_gnssro_bndropp1d_simobs_ad(self, geovals, hofx, obss)
subroutine ufo_gnssro_bndropp1d_tlad_settraj(self, geovals, obss)
subroutine ufo_gnssro_bndropp1d_tlad_delete(self)
subroutine ufo_gnssro_bndropp1d_simobs_tl(self, geovals, hofx, obss)
character(len=maxvarlen), parameter, public var_prs
character(len=maxvarlen), parameter, public var_q
character(len=maxvarlen), parameter, public var_sfc_geomz
character(len=maxvarlen), parameter, public var_z
character(len=maxvarlen), parameter, public var_ts
Fortran module to perform linear interpolation.
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators
Fortran derived type for gnssro trajectory.