UFO
ufo_gnssro_bndropp2d_tlad_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2017-2018 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
6 !> Fortran module for gnssro bending angle ropp2d tangent linear and adjoint
7 !> following the ROPP (2018 Aug) implementation
8 
10 
11 use fckit_configuration_module, only: fckit_configuration
12 !use iso_c_binding
13 use kinds
14 use ufo_vars_mod
19 use obsspace_mod
21 use missing_values_mod
24 use fckit_log_module, only : fckit_log
25 
26 integer, parameter :: max_string=800
27 
28 !> Fortran derived type for gnssro trajectory
30  private
31  integer :: nval, nlocs
32  real(kind_real), allocatable :: prs(:,:), t(:,:), q(:,:), gph(:,:), gph_sfc(:,:)
33  integer :: iflip ! geoval ascending order flag
34  type(gnssro_conf) :: roconf ! ro configuration
35  real(kind_real), allocatable :: obslon2d(:), obslat2d(:) !2d locations - nlocs*n_horiz
36  contains
37  procedure :: setup => ufo_gnssro_bndropp2d_tlad_setup
38  procedure :: delete => ufo_gnssro_bndropp2d_tlad_delete
39  procedure :: settraj => ufo_gnssro_bndropp2d_tlad_settraj
40  procedure :: simobs_tl => ufo_gnssro_bndropp2d_simobs_tl
41  procedure :: simobs_ad => ufo_gnssro_bndropp2d_simobs_ad
43 
44 contains
45 
46 ! ------------------------------------------------------------------------------
47 subroutine ufo_gnssro_bndropp2d_tlad_setup(self, f_conf)
48  implicit none
49  class(ufo_gnssro_bndropp2d_tlad), intent(inout) :: self
50  type(fckit_configuration), intent(in) :: f_conf
51 
52  call gnssro_conf_setup(self%roconf,f_conf)
53 
55 
56 ! ------------------------------------------------------------------------------
57 subroutine ufo_gnssro_bndropp2d_tlad_settraj(self, geovals, obss)
58 
59  implicit none
60  class(ufo_gnssro_bndropp2d_tlad), intent(inout) :: self
61  type(ufo_geovals), intent(in) :: geovals
62  type(c_ptr), value, intent(in) :: obss
63  character(len=*), parameter :: myname_="ufo_gnssro_bndropp2d_tlad_settraj"
64  character(max_string) :: err_msg
65  type(ufo_geoval), pointer :: t, q, prs, gph, gph_sfc
66  integer :: i, kerror
67  real(kind_real), allocatable :: obsAzim(:) ! nlocs
68  real(kind_real), allocatable :: obsLat(:), obsLon(:) ! nlocs
69  real(kind_real), allocatable :: obsLonnh(:),obsLatnh(:) ! n_horiz
70  integer :: n_horiz
71  real(kind_real) :: dtheta
72 
73  write(err_msg,*) "TRACE: ufo_gnssro_bndropp2d_tlad_settraj: begin"
74  call fckit_log%info(err_msg)
75 
76 ! get model state variables from geovals
77  call ufo_geovals_get_var(geovals, var_ts, t) ! temperature
78  call ufo_geovals_get_var(geovals, var_q, q) ! specific humidity
79  call ufo_geovals_get_var(geovals, var_prs, prs) ! pressure
80  call ufo_geovals_get_var(geovals, var_z, gph) ! geopotential height
81  call ufo_geovals_get_var(geovals, var_sfc_geomz, gph_sfc) ! surface geopotential height
82  call self%delete()
83 
84  self%nval = prs%nval
85  self%nlocs = obsspace_get_nlocs(obss)
86  self%iflip = 0
87 
88  n_horiz = self%roconf%n_horiz
89  dtheta = self%roconf%dtheta
90 
91  if (prs%vals(1,1) .lt. prs%vals(prs%nval,1) ) then
92  self%iflip = 1
93  write(err_msg,'(a)') ' ufo_gnssro_bndropp2d_tlad_settraj:'//new_line('a')// &
94  ' Model vertical height profile is in descending order,'//new_line('a')// &
95  ' but ROPP requires it to be ascending order, need flip'
96  call fckit_log%info(err_msg)
97  end if
98 
99  allocate(self%obsLat2d(self%nlocs*n_horiz))
100  allocate(self%obsLon2d(self%nlocs*n_horiz))
101 
102  allocate(obslon(self%nlocs))
103  allocate(obslat(self%nlocs))
104  allocate(obsazim(self%nlocs))
105 
106  call obsspace_get_db(obss, "MetaData", "longitude", obslon)
107  call obsspace_get_db(obss, "MetaData", "latitude", obslat)
108  call obsspace_get_db(obss, "MetaData", "sensor_azimuth_angle", obsazim)
109 
110  allocate(obslatnh(n_horiz))
111  allocate(obslonnh(n_horiz))
112 
113  do i = 1, self%nlocs
114  call ropp_fm_2d_plane(obslat(i),obslon(i),obsazim(i),dtheta,n_horiz,obslatnh,obslonnh,kerror)
115  self%obsLon2d((i-1)*n_horiz+1:i*n_horiz) = obslonnh
116  self%obsLat2d((i-1)*n_horiz+1:i*n_horiz) = obslatnh
117  end do
118 
119  deallocate(obslat)
120  deallocate(obslon)
121  deallocate(obslonnh)
122  deallocate(obslatnh)
123  deallocate(obsazim)
124 
125  allocate(self%t(self%nval,self%nlocs*n_horiz))
126  allocate(self%q(self%nval,self%nlocs*n_horiz))
127  allocate(self%prs(self%nval,self%nlocs*n_horiz))
128  allocate(self%gph(self%nval,self%nlocs*n_horiz))
129  allocate(self%gph_sfc(1,self%nlocs*n_horiz))
130 
131 ! allocate
132  self%gph = gph%vals
133  self%t = t%vals
134  self%q = q%vals
135  self%prs = prs%vals
136  self%gph_sfc = gph_sfc%vals
137 
138  self%ltraj = .true.
139 
141 
142 ! ------------------------------------------------------------------------------
143 ! ------------------------------------------------------------------------------
144 subroutine ufo_gnssro_bndropp2d_simobs_tl(self, geovals, hofx, obss)
145 
146  use ropp_fm_types, only: state2dfm, state1dfm
147  use ropp_fm_types, only: obs1dbangle
148  use datetimetypes, only: dp
149  implicit none
150  class(ufo_gnssro_bndropp2d_tlad), intent(in) :: self
151  type(ufo_geovals), intent(in) :: geovals ! perturbed quantities
152  real(kind_real), intent(inout) :: hofx(:)
153  type(c_ptr), value, intent(in) :: obss
154  real(c_double) :: missing
155 
156  type(state2dfm) :: x,x_tl
157  type(state1dfm) :: x1d,x1d_tl
158  type(obs1dbangle) :: y,y_tl
159 
160  integer :: iobs,nlev, nlocs,nvprof
161 
162  character(len=*), parameter :: myname_="ufo_gnssro_bndropp2d_simobs_tl"
163  character(max_string) :: err_msg
164  type(ufo_geoval), pointer :: t_d, q_d, prs_d
165 
166 ! hack - set local geopotential height to zero for ropp routines
167  real(kind_real), allocatable :: gph_d_zero(:,:)
168  real(kind_real) :: gph_sfc_d_zero
169 ! hack - set local geopotential height to zero for ropp routines
170  real(kind_real), allocatable :: obsImpP(:),obsLocR(:),obsGeoid(:),obsAzim(:) !nlocs
171  real(kind_real), allocatable :: obsLat(:),obsLon(:) !nlocs
172  integer :: n_horiz
173  real(kind_real) :: dtheta
174  real(kind_real) :: ob_time
175 
176  missing = missing_value(missing)
177 
178  n_horiz = self%roconf%n_horiz
179  dtheta = self%roconf%dtheta
180 
181  write(err_msg,*) "TRACE: ufo_gnssro_bndropp2d_simobs_tl: begin"
182  call fckit_log%info(err_msg)
183 
184 ! check if trajectory was set
185  if (.not. self%ltraj) then
186  write(err_msg,*) myname_, ' trajectory wasnt set!'
187  call abor1_ftn(err_msg)
188  endif
189 
190 ! check if nlocs is consistent in geovals & hofx
191  if (geovals%nlocs /= size(hofx)*n_horiz ) then
192  write(err_msg,*) myname_, ' error: 2d nlocs inconsistent! geovals%nlocs, size(hofx), &
193  and n_horiz are', geovals%nlocs, size(hofx), n_horiz
194  call abor1_ftn(err_msg)
195  endif
196 
197 ! get variables from geovals
198  call ufo_geovals_get_var(geovals, var_ts, t_d) ! temperature
199  call ufo_geovals_get_var(geovals, var_q, q_d) ! specific humidity
200  call ufo_geovals_get_var(geovals, var_prs, prs_d) ! pressure
201 
202  nlev = self%nval
203  nlocs = self%nlocs
204 
205  allocate(gph_d_zero(nlev,nlocs*n_horiz))
206  gph_d_zero = 0.0
207  gph_sfc_d_zero = 0.0
208 
209 ! set obs space struture
210  allocate(obslon(nlocs))
211  allocate(obslat(nlocs))
212  allocate(obsimpp(nlocs))
213  allocate(obslocr(nlocs))
214  allocate(obsgeoid(nlocs))
215  allocate(obsazim(nlocs))
216  call obsspace_get_db(obss, "MetaData", "longitude", obslon)
217  call obsspace_get_db(obss, "MetaData", "latitude", obslat)
218  call obsspace_get_db(obss, "MetaData", "impact_parameter", obsimpp)
219  call obsspace_get_db(obss, "MetaData", "earth_radius_of_curvature", obslocr)
220  call obsspace_get_db(obss, "MetaData", "geoid_height_above_reference_ellipsoid", obsgeoid)
221  call obsspace_get_db(obss, "MetaData", "sensor_azimuth_angle", obsazim)
222 
223  nvprof = 1 ! no. of bending angles in profile
224  ob_time = 0.0
225 
226 ! loop through the obs
227  obs_loop: do iobs = 1, nlocs ! order of loop doesn't matter
228 
229  if ( ( obsimpp(iobs)-obslocr(iobs)-obsgeoid(iobs) ) <= self%roconf%top_2d .and. &
230  obsazim(iobs) /= missing ) then
231 
232 ! map the trajectory to ROPP 2D structure x
233  call init_ropp_2d_statevec(self%obsLon2d( (iobs-1)*n_horiz+1:iobs*n_horiz ), &
234  self%obsLat2d( (iobs-1)*n_horiz+1:iobs*n_horiz ), &
235  self%t(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
236  self%q(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
237  self%prs(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
238  self%gph(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
239  nlev, x, n_horiz, dtheta, self%iflip)
240 
241 ! hack -- make non zero humidity to avoid zero denominator in tangent linear
242 ! see ropp_fm/bangle_1d/ropp_fm_bangle_1d_tl.f90
243  where(x%shum .le. 1e-8) x%shum = 1e-8
244 ! hack -- make non zero humidity to avoid zero denominator in tangent linear
245 
246  call init_ropp_2d_statevec(self%obsLon2d( (iobs-1)*n_horiz+1:iobs*n_horiz ), &
247  self%obsLat2d( (iobs-1)*n_horiz+1:iobs*n_horiz ), &
248  t_d%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
249  q_d%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
250  prs_d%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
251  gph_d_zero(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
252  nlev, x_tl, n_horiz, dtheta, self%iflip)
253 
254 ! set both y and y_tl structures
255  call init_ropp_2d_obvec_tlad(iobs, nvprof, &
256  obsimpp(iobs), &
257  obslat(iobs), &
258  obslon(iobs), &
259  obslocr(iobs), &
260  obsgeoid(iobs), &
261  y,y_tl)
262 
263 ! now call TL of forward model
264  call ropp_fm_bangle_2d_tl(x,x_tl,y, y_tl)
265  hofx(iobs) = y_tl%bangle(nvprof) ! this will need to change if profile is passed
266 
267 ! tidy up -deallocate ropp structures
268  call ropp_tidy_up_tlad_2d(x,x_tl,y,y_tl)
269 
270  else ! apply ropp1d above top_2d or when azimuth angle is missing
271 
272 ! map the trajectory to ROPP 1D structure x1d
273  call init_ropp_1d_statevec(ob_time, &
274  obslon(iobs), &
275  obslat(iobs), &
276  self%t(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
277  self%q(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
278  self%prs(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
279  self%gph(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
280  nlev, &
281  self%gph_sfc(1,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
282  x1d, self%iflip)
283 
284  where(x1d%shum .le. 1e-8) x1d%shum = 1e-8
285 
286  call init_ropp_1d_statevec( ob_time, &
287  obslon(iobs), &
288  obslat(iobs), &
289  t_d%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
290  q_d%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
291  prs_d%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
292  gph_d_zero(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
293  nlev, &
294  gph_sfc_d_zero, &
295  x1d_tl, self%iflip)
296 
297 ! y and y_tl structures
298  call init_ropp_1d_obvec_tlad(iobs, nvprof, &
299  obsimpp(iobs), &
300  obslat(iobs), &
301  obslon(iobs), &
302  obslocr(iobs), &
303  obsgeoid(iobs), &
304  y,y_tl)
305 
306 ! TL
307  call ropp_fm_bangle_1d_tl(x1d,x1d_tl,y,y_tl%bangle(nvprof))
308  hofx(iobs) = y_tl%bangle(nvprof)
309 
310 ! tidy up
311  call ropp_tidy_up_tlad_1d(x1d,x1d_tl,y,y_tl)
312  end if
313  end do obs_loop
314 
315 ! tidy up - deallocate obsspace structures
316  deallocate(obslat)
317  deallocate(obslon)
318  deallocate(obsimpp)
319  deallocate(obslocr)
320  deallocate(obsgeoid)
321  deallocate(obsazim)
322  deallocate(gph_d_zero)
323 
324  write(err_msg,*) "TRACE: ufo_gnssro_bndropp2d_simobs_tl: complete"
325  call fckit_log%info(err_msg)
326 
327  return
328 
329 end subroutine ufo_gnssro_bndropp2d_simobs_tl
330 
331 ! ------------------------------------------------------------------------------
332 ! ------------------------------------------------------------------------------
333 subroutine ufo_gnssro_bndropp2d_simobs_ad(self, geovals, hofx, obss)
334 
335  use ropp_fm_types, only: state2dfm, state1dfm
336  use ropp_fm_types, only: obs1dbangle
337  use typesizes, only: wp => eightbytereal
338  use datetimetypes, only: dp
339 
340  implicit none
341  class(ufo_gnssro_bndropp2d_tlad), intent(in) :: self
342  type(ufo_geovals), intent(inout) :: geovals ! perturbed quantities
343  real(kind_real), intent(in) :: hofx(:)
344  type(c_ptr), value, intent(in) :: obss
345  real(c_double) :: missing
346 
347  type(ufo_geoval), pointer :: t_d, q_d, prs_d
348 
349 ! set local geopotential height to zero for ropp routines
350  real(kind_real), allocatable :: gph_d_zero(:,:)
351  real(kind_real) :: gph_sfc_d_zero
352 
353  real(kind_real), allocatable :: obsLat(:), obsLon(:), obsImpP(:), obsLocR(:), obsGeoid(:)
354  real(kind_real), allocatable :: obsAzim(:)
355  type(state2dfm) :: x,x_ad
356  type(state1dfm) :: x1d,x1d_ad
357  type(obs1dbangle) :: y,y_ad
358  integer :: iobs,nlev,nlocs,nvprof
359  character(len=*), parameter :: myname_="ufo_gnssro_bndropp2d_simobs_ad"
360  character(max_string) :: err_msg
361  integer :: n_horiz
362  real(kind_real) :: dtheta
363  real(kind_real) :: ob_time
364 
365  n_horiz = self%roconf%n_horiz
366  dtheta = self%roconf%dtheta
367 
368  write(err_msg,*) "TRACE: ufo_gnssro_bndropp2d_simobs_ad: begin"
369  call fckit_log%info(err_msg)
370 
371 ! check if trajectory was set
372  if (.not. self%ltraj) then
373  write(err_msg,*) myname_, ' trajectory wasnt set!'
374  call abor1_ftn(err_msg)
375  endif
376 ! check if nlocs is consistent in geovals & hofx
377  if (geovals%nlocs /= size(hofx)*n_horiz) then
378  write(err_msg,*) myname_, ' error: 2d nlocs inconsistent!'
379  call abor1_ftn(err_msg)
380  endif
381 
382 ! get variables from geovals
383  call ufo_geovals_get_var(geovals, var_ts, t_d) ! temperature
384  call ufo_geovals_get_var(geovals, var_q, q_d) ! specific humidity
385  call ufo_geovals_get_var(geovals, var_prs, prs_d) ! pressure
386 
387 ! allocate if not yet allocated
388  if (.not. allocated(t_d%vals)) then
389  t_d%nlocs = self%nlocs*n_horiz
390  t_d%nval = self%nval
391  allocate(t_d%vals(t_d%nval,t_d%nlocs))
392  t_d%vals = 0.0_kind_real
393  endif
394 
395  if (.not. allocated(prs_d%vals)) then
396  prs_d%nlocs = self%nlocs*n_horiz
397  prs_d%nval = self%nval
398  allocate(prs_d%vals(prs_d%nval,prs_d%nlocs))
399  prs_d%vals = 0.0_kind_real
400  endif
401 
402  if (.not. allocated(q_d%vals)) then
403  q_d%nlocs = self%nlocs*n_horiz
404  q_d%nval = self%nval
405  allocate(q_d%vals(q_d%nval,q_d%nlocs))
406  q_d%vals = 0.0_kind_real
407  endif
408 
409  if (.not. geovals%linit ) geovals%linit=.true.
410 
411  nlev = self%nval
412  nlocs = self%nlocs
413 
414  allocate(gph_d_zero(nlev,nlocs*n_horiz))
415  gph_d_zero = 0.0
416  gph_sfc_d_zero = 0.0
417 
418 ! set obs space struture
419  allocate(obslon(nlocs))
420  allocate(obslat(nlocs))
421  allocate(obsimpp(nlocs))
422  allocate(obslocr(nlocs))
423  allocate(obsgeoid(nlocs))
424  allocate(obsazim(nlocs))
425 
426  call obsspace_get_db(obss, "MetaData", "longitude", obslon)
427  call obsspace_get_db(obss, "MetaData", "latitude", obslat)
428  call obsspace_get_db(obss, "MetaData", "impact_parameter", obsimpp)
429  call obsspace_get_db(obss, "MetaData", "earth_radius_of_curvature", obslocr)
430  call obsspace_get_db(obss, "MetaData", "geoid_height_above_reference_ellipsoid", obsgeoid)
431  call obsspace_get_db(obss, "MetaData", "sensor_azimuth_angle", obsazim)
432 
433  missing = missing_value(missing)
434 
435 ! loop through the obs
436  nvprof = 1 ! no. of bending angles in profile
437  ob_time = 0.0
438 
439  obs_loop: do iobs = 1, nlocs
440 
441  if (hofx(iobs) .gt. missing) then
442  if ( ( obsimpp(iobs)-obslocr(iobs)-obsgeoid(iobs) ) <= self%roconf%top_2d .and. &
443  obsazim(iobs) /= missing ) then
444 
445 ! map the trajectory to ROPP structure x
446  call init_ropp_2d_statevec(self%obsLon2d((iobs-1)*n_horiz+1:iobs*n_horiz), &
447  self%obsLat2d((iobs-1)*n_horiz+1:iobs*n_horiz), &
448  self%t(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
449  self%q(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
450  self%prs(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
451  self%gph(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
452  nlev, x, n_horiz, dtheta, self%iflip)
453 
454  call init_ropp_2d_statevec(self%obsLon2d( (iobs-1)*n_horiz+1:iobs*n_horiz ), &
455  self%obsLat2d( (iobs-1)*n_horiz+1:iobs*n_horiz ), &
456  t_d%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
457  q_d%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
458  prs_d%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
459  gph_d_zero(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
460  nlev, x_ad, n_horiz, dtheta, self%iflip)
461 
462  ! x_ad is local so initialise to 0.0
463  x_ad%temp(:,:) = 0.0_wp
464  x_ad%pres(:,:) = 0.0_wp
465  x_ad%shum(:,:) = 0.0_wp
466  x_ad%geop(:,:) = 0.0_wp
467 
468  ! set both y and y_ad structures
469  call init_ropp_2d_obvec_tlad(iobs, nvprof, &
470  obsimpp(iobs), &
471  obslat(iobs), &
472  obslon(iobs), &
473  obslocr(iobs), &
474  obsgeoid(iobs), &
475  y,y_ad)
476 
477 
478 ! local variable initialise
479  y_ad%bangle(:) = 0.0_wp
480 
481 ! now call AD of forward model
482  y_ad%bangle(nvprof) = y_ad%bangle(nvprof) + hofx(iobs)
483  call ropp_fm_bangle_2d_ad(x,x_ad,y,y_ad)
484 
486  t_d%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
487  q_d%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
488  prs_d%vals(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
489  gph_d_zero(:,(iobs-1)*n_horiz+1:iobs*n_horiz), &
490 
491  nlev, x_ad, n_horiz,self%iflip)
492 
493 ! tidy up - deallocate ropp structures
494  call ropp_tidy_up_tlad_2d(x,x_ad,y,y_ad)
495 
496  else ! apply ropp1d above top_2d or when azimuth angle is missing
497 
498 ! map the trajectory to ROPP 1D structure x1d
499  call init_ropp_1d_statevec( ob_time, &
500  obslon(iobs), &
501  obslat(iobs), &
502  self%t(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
503  self%q(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
504  self%prs(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
505  self%gph(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
506  nlev, &
507  self%gph_sfc(1,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
508  x1d, self%iflip)
509 
510  call init_ropp_1d_statevec( ob_time, &
511  obslon(iobs), &
512  obslat(iobs), &
513  t_d%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
514  q_d%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
515  prs_d%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
516  gph_d_zero(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
517  nlev, &
518  gph_sfc_d_zero, &
519  x1d_ad, self%iflip)
520 
521 
522  ! x_ad is local so initialise to 0.0
523  x1d_ad%temp(:) = 0.0_wp
524  x1d_ad%pres(:) = 0.0_wp
525  x1d_ad%shum(:) = 0.0_wp
526  x1d_ad%geop(:) = 0.0_wp
527 
528  ! set both y and y_ad structures
529  call init_ropp_1d_obvec_tlad(iobs, nvprof, &
530  obsimpp(iobs), &
531  obslat(iobs), &
532  obslon(iobs), &
533  obslocr(iobs), &
534  obsgeoid(iobs), &
535  y,y_ad)
536 
537 
538 ! local variable initialise
539  y_ad%bangle(:) = 0.0_wp
540 
541 ! now call AD of forward model
542  y_ad%bangle(nvprof) = y_ad%bangle(nvprof) + hofx(iobs)
543  call ropp_fm_bangle_1d_ad(x1d,x1d_ad,y,y_ad)
545  t_d%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
546  q_d%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
547  prs_d%vals(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
548  gph_d_zero(:,(iobs-1)*n_horiz+1+(n_horiz-1)/2), &
549  nlev, x1d_ad, self%iflip)
550 
551 ! tidy up
552  call ropp_tidy_up_tlad_1d(x1d,x1d_ad,y,y_ad)
553 
554  end if ! end top_2d check
555 
556  end if ! end missing value check
557 
558  end do obs_loop
559 
560 ! tidy up - deallocate obsspace structures
561  deallocate(obslat)
562  deallocate(obslon)
563  deallocate(obsimpp)
564  deallocate(obslocr)
565  deallocate(obsgeoid)
566  deallocate(obsazim)
567  deallocate(gph_d_zero)
568 
569  write(err_msg,*) "TRACE: ufo_gnssro_bndropp2d_simobs_ad: complete"
570  call fckit_log%info(err_msg)
571 
572  return
573 
574 end subroutine ufo_gnssro_bndropp2d_simobs_ad
575 
576 !-------------------------------------------------------------------------
577 !-------------------------------------------------------------------------
579 
580  implicit none
581  class(ufo_gnssro_bndropp2d_tlad), intent(inout) :: self
582  character(len=*), parameter :: myname_="ufo_gnssro_bndropp_tlad_delete"
583 
584  self%nval = 0
585  if (allocated(self%prs)) deallocate(self%prs)
586  if (allocated(self%t)) deallocate(self%t)
587  if (allocated(self%q)) deallocate(self%q)
588  if (allocated(self%gph)) deallocate(self%gph)
589  if (allocated(self%gph_sfc)) deallocate(self%gph_sfc)
590  if (allocated(self%obsLat2d)) deallocate(self%obsLat2d)
591  if (allocated(self%obsLon2d)) deallocate(self%obsLon2d)
592 
593  self%ltraj = .false.
594 
596 
597 !-------------------------------------------------------------------------
598 
ufo_avgkernel_mod::max_string
integer, parameter max_string
Definition: ufo_avgkernel_mod.F90:17
ufo_gnssro_ropp2d_utils_mod::init_ropp_2d_obvec_tlad
subroutine, public init_ropp_2d_obvec_tlad(iloop, nvprof, obs_impact, rlat, rlon, roc, undulat, y, y_p)
Definition: ufo_gnssro_ropp2d_utils_mod.F90:271
ufo_basis_tlad_mod
Definition: ufo_basis_tlad_mod.F90:6
ufo_gnssro_ropp2d_utils_mod::ropp_tidy_up_tlad_2d
subroutine, public ropp_tidy_up_tlad_2d(x, x_p, y, y_p)
Definition: ufo_gnssro_ropp2d_utils_mod.F90:374
ufo_gnssro_ropp1d_utils_mod::init_ropp_1d_statevec
subroutine, public init_ropp_1d_statevec(step_time, rlon, rlat, temp, shum, pres, phi, lm, phi_sfc, x, iflip)
Definition: ufo_gnssro_ropp1d_utils_mod.F90:37
ufo_gnssro_bndropp2d_tlad_mod::ufo_gnssro_bndropp2d_tlad_setup
subroutine ufo_gnssro_bndropp2d_tlad_setup(self, f_conf)
Definition: ufo_gnssro_bndropp2d_tlad_mod.F90:48
ufo_gnssro_ropp1d_utils_mod::init_ropp_1d_obvec_tlad
subroutine, public init_ropp_1d_obvec_tlad(iloop, nvprof, obs_impact, rlat, rlon, roc, undulat, y, y_p)
Definition: ufo_gnssro_ropp1d_utils_mod.F90:326
ufo_gnssro_ropp1d_utils_mod::ropp_tidy_up_tlad_1d
subroutine, public ropp_tidy_up_tlad_1d(x, x_p, y, y_p)
Definition: ufo_gnssro_ropp1d_utils_mod.F90:399
ufo_gnssro_ropp2d_utils_mod::init_ropp_2d_statevec_ad
subroutine, public init_ropp_2d_statevec_ad(temp_d, shum_d, pres_d, phi_d, lm, x_ad, n_horiz, iflip)
Definition: ufo_gnssro_ropp2d_utils_mod.F90:125
ufo_gnssro_bndropp2d_tlad_mod::ufo_gnssro_bndropp2d_tlad_delete
subroutine ufo_gnssro_bndropp2d_tlad_delete(self)
Definition: ufo_gnssro_bndropp2d_tlad_mod.F90:579
ufo_gnssro_bndropp2d_tlad_mod::ufo_gnssro_bndropp2d_tlad_settraj
subroutine ufo_gnssro_bndropp2d_tlad_settraj(self, geovals, obss)
Definition: ufo_gnssro_bndropp2d_tlad_mod.F90:58
ufo_gnssro_ropp2d_utils_mod::init_ropp_2d_statevec
subroutine, public init_ropp_2d_statevec(rlon, rlat, temp, shum, pres, phi, lm, x, n_horiz, dtheta, iflip)
Definition: ufo_gnssro_ropp2d_utils_mod.F90:36
ufo_gnssro_ropp1d_utils_mod::init_ropp_1d_statevec_ad
subroutine, public init_ropp_1d_statevec_ad(temp_d, shum_d, pres_d, phi_d, lm, x_ad, iflip)
Definition: ufo_gnssro_ropp1d_utils_mod.F90:166
ufo_gnssro_bndropp2d_tlad_mod::ufo_gnssro_bndropp2d_simobs_tl
subroutine ufo_gnssro_bndropp2d_simobs_tl(self, geovals, hofx, obss)
Definition: ufo_gnssro_bndropp2d_tlad_mod.F90:145
gnssro_mod_conf::gnssro_conf
Definition: gnssro_mod_conf.F90:14
ufo_gnssro_ropp2d_utils_mod
Fortran module to handle gnssro bending angle observations following the ROPP (2018 Aug) implementati...
Definition: ufo_gnssro_ropp2d_utils_mod.F90:9
vert_interp_mod
Fortran module to perform linear interpolation.
Definition: vert_interp.F90:8
ufo_basis_tlad_mod::ufo_basis_tlad
Definition: ufo_basis_tlad_mod.F90:12
ufo_geovals_mod
Definition: ufo_geovals_mod.F90:7
ufo_geovals_mod_c
Definition: GeoVaLs.interface.F90:7
ufo_gnssro_bndropp2d_tlad_mod
Fortran module for gnssro bending angle ropp2d tangent linear and adjoint following the ROPP (2018 Au...
Definition: ufo_gnssro_bndropp2d_tlad_mod.F90:9
ufo_gnssro_bndropp2d_tlad_mod::ufo_gnssro_bndropp2d_tlad
Fortran derived type for gnssro trajectory.
Definition: ufo_gnssro_bndropp2d_tlad_mod.F90:29
ufo_vars_mod::var_sfc_geomz
character(len=maxvarlen), parameter, public var_sfc_geomz
Definition: ufo_variables_mod.F90:70
gnssro_mod_conf::gnssro_conf_setup
subroutine, public gnssro_conf_setup(roconf, f_conf)
Definition: gnssro_mod_conf.F90:35
ufo_vars_mod
Definition: ufo_variables_mod.F90:8
ufo_vars_mod::var_z
character(len=maxvarlen), parameter, public var_z
Definition: ufo_variables_mod.F90:29
ufo_geovals_mod::ufo_geovals_get_var
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Definition: ufo_geovals_mod.F90:128
ufo_vars_mod::var_q
character(len=maxvarlen), parameter, public var_q
Definition: ufo_variables_mod.F90:22
ufo_vars_mod::var_ts
character(len=maxvarlen), parameter, public var_ts
Definition: ufo_variables_mod.F90:19
ufo_geovals_mod::ufo_geovals
type to hold interpolated fields required by the obs operators
Definition: ufo_geovals_mod.F90:47
gnssro_mod_conf
Definition: gnssro_mod_conf.F90:2
ufo_gnssro_ropp1d_utils_mod
Fortran module to handle gnssro bending angle observations following the ROPP (2018 Aug) implementati...
Definition: ufo_gnssro_ropp1d_utils_mod.F90:9
ufo_gnssro_bndropp2d_tlad_mod::ufo_gnssro_bndropp2d_simobs_ad
subroutine ufo_gnssro_bndropp2d_simobs_ad(self, geovals, hofx, obss)
Definition: ufo_gnssro_bndropp2d_tlad_mod.F90:334
ufo_geovals_mod::ufo_geoval
type to hold interpolated field for one variable, one observation
Definition: ufo_geovals_mod.F90:40
ufo_vars_mod::var_prs
character(len=maxvarlen), parameter, public var_prs
Definition: ufo_variables_mod.F90:25
ufo_geovals_mod_c::ufo_geovals_registry
type(registry_t), public ufo_geovals_registry
Linked list interface - defines registry_t type.
Definition: GeoVaLs.interface.F90:30