9 use fckit_configuration_module,
only: fckit_configuration
19 use missing_values_mod
25 integer :: nval, nlocs
26 real(kind_real),
allocatable :: wf(:)
27 integer,
allocatable :: wi(:)
28 real(kind_real),
allocatable :: jac_t(:), jac_q(:), jac_prs(:)
43 type(fckit_configuration),
intent(in) :: f_conf
56 type(c_ptr),
value,
intent(in) :: obss
58 character(len=*),
parameter :: myname_=
"ufo_gnssro_ref_tlad_settraj"
61 real(kind_real),
allocatable :: obsZ(:), obsLat(:)
62 real(kind_real) :: obsH,gesT,gesQ,gesP
63 real(kind_real) :: Tv, Tv0
77 self%nlocs = obsspace_get_nlocs(obss)
79 allocate(self%wi(self%nlocs))
80 allocate(self%wf(self%nlocs))
81 allocate(self%jac_t(self%nlocs))
82 allocate(self%jac_q(self%nlocs))
83 allocate(self%jac_prs(self%nlocs))
84 allocate(obsz(self%nlocs))
85 allocate(obslat(self%nlocs))
88 call obsspace_get_db(obss,
"MetaData",
"altitude", obsz)
89 call obsspace_get_db(obss,
"MetaData",
"latitude", obslat)
92 do iobs = 1, self%nlocs
104 tv = gest*( one + (rv_over_rd-one)*gesq/(1.0-gesq) )
105 tv0 = t%vals(wi0,iobs)*(one + (rv_over_rd-one)*q%vals(wi0,iobs)/(1.0-q%vals(wi0,iobs) ))
106 gesp = prs%vals(wi0,iobs)/exp(two*grav*(obsh-gph%vals(wi0,iobs))/(rd*(tv+tv0)))
108 self%jac_t(iobs) = -
n_a*gesp/gest**2 &
109 -
n_b*two*gesp*gesq/( ((1-rd_over_rv)*gesq+rd_over_rv)*gest**3 ) &
110 -
n_c*gesp*gesq/( ((1-rd_over_rv)*gesq+rd_over_rv)*gest**2 )
112 self%jac_q(iobs) =
n_b*gesp/( gest**2*( (1-rd_over_rv)*gesq+rd_over_rv)**2 ) &
114 +
n_c*gesp/( gest *( (1-rd_over_rv)*gesq+rd_over_rv)**2 ) &
116 self%jac_prs(iobs) =
n_a/gest &
117 +
n_b*gesq/ ( ((1-rd_over_rv)*gesq+rd_over_rv)*gest**2 ) &
118 +
n_c*gesq/ ( ((1-rd_over_rv)*gesq+rd_over_rv)*gest )
136 real(kind_real),
intent(inout) :: hofx(:)
137 type(c_ptr),
value,
intent(in) :: obss
139 character(len=*),
parameter :: myname=
"ufo_gnssro_ref_tlad_tl"
140 character(max_string) :: err_msg
142 type(
ufo_geoval),
pointer :: t_tl, q_tl, prs_tl
143 real(kind_real) :: gesT_tl, gesQ_tl, gesP_tl
147 if (.not. self%ltraj)
then
148 write(err_msg,*) myname,
' trajectory wasnt set!'
149 call abor1_ftn(err_msg)
153 if (geovals%nlocs /=
size(hofx))
then
154 write(err_msg,*) myname,
' error: nlocs inconsistent!'
155 call abor1_ftn(err_msg)
164 do iobs = 1, geovals%nlocs
167 call vert_interp_apply_tl(prs_tl%nval,prs_tl%vals(:,iobs), gesp_tl, self%wi(iobs),self%wf(iobs))
168 hofx(iobs) = self%jac_t(iobs)*gest_tl + self%jac_q(iobs)*gesq_tl + self%jac_prs(iobs)*gesp_tl
178 real(kind_real),
intent(in) :: hofx(:)
179 type(c_ptr),
value,
intent(in) :: obss
181 character(len=*),
parameter :: myname=
"ufo_gnssro_ref_tlad_ad"
182 character(max_string) :: err_msg
184 real(kind_real) :: gesT_d, gesQ_d, gesP_d
185 real(c_double) :: missing
189 if (.not. self%ltraj)
then
190 write(err_msg,*) myname,
' trajectory wasnt set!'
191 call abor1_ftn(err_msg)
195 if (geovals%nlocs /=
size(hofx))
then
196 write(err_msg,*) myname,
' error: nlocs inconsistent!'
197 call abor1_ftn(err_msg)
205 missing = missing_value(missing)
207 do iobs = 1, geovals%nlocs
209 if (hofx(iobs) .ne. missing)
then
210 gest_d = 0.0_kind_real
211 gesq_d = 0.0_kind_real
212 gesp_d = 0.0_kind_real
213 gest_d = gest_d + hofx(iobs)*self%jac_t(iobs)
214 gesq_d = gesq_d + hofx(iobs)*self%jac_q(iobs)
215 gesp_d = gesp_d + hofx(iobs)*self%jac_prs(iobs)
229 character(len=*),
parameter :: myname_=
"ufo_gnssro_ref_tlad_delete"
232 if (
allocated(self%wi))
deallocate(self%wi)
233 if (
allocated(self%wf))
deallocate(self%wf)
234 if (
allocated(self%jac_t))
deallocate(self%jac_t)
235 if (
allocated(self%jac_q))
deallocate(self%jac_q)
236 if (
allocated(self%jac_prs))
deallocate(self%jac_prs)
subroutine, public gnssro_conf_setup(roconf, f_conf)
subroutine, public gnssro_ref_constants(use_compress)
real(kind_real), public n_a
real(kind_real), public n_b
real(kind_real), public n_c
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 refractivity tangent linear and adjoint.
subroutine ufo_gnssro_ref_tlad_settraj(self, geovals, obss)
subroutine ufo_gnssro_ref_simobs_ad(self, geovals, hofx, obss)
subroutine ufo_gnssro_ref_tlad_setup(self, f_conf)
subroutine ufo_gnssro_ref_tlad_delete(self)
subroutine ufo_gnssro_ref_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_z
character(len=maxvarlen), parameter, public var_ts
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(nlev, fvec, f, 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
Fortran derived type for gnssro trajectory.