UFO
ufo_radiancerttov_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 radiancerttov tl/ad observation operator
7 
9 
10  use fckit_configuration_module, only: fckit_configuration
11  use fckit_log_module, only : fckit_log
12  use iso_c_binding
13  use kinds
14  use missing_values_mod
15 
16  use obsspace_mod
17 
19  use ufo_vars_mod
21 
22  use rttov_types
23  use rttov_const ! errorstatus and gas_id
24  use rttov_unix_env
25 
26  implicit none
27  private
28 
29  !> Fortran derived type for radiancerttov trajectory
30  type, public :: ufo_radiancerttov_tlad
31  private
32  character(len=MAXVARLEN), public, allocatable :: varin(:) ! variables requested from the model
33  integer, allocatable :: channels(:)
34  type(rttov_conf) :: conf
35  type(rttov_conf) :: conf_traj
36  type(ufo_rttov_io) :: rtprof_k
37 
38  integer :: nprofiles
39  integer :: nchan_total
40  integer :: nlevels
41 
42  logical :: ltraj
43 
44  contains
45  procedure :: setup => ufo_radiancerttov_tlad_setup
46  procedure :: delete => ufo_radiancerttov_tlad_delete
47  procedure :: settraj => ufo_radiancerttov_tlad_settraj
48  procedure :: simobs_tl => ufo_radiancerttov_simobs_tl
49  procedure :: simobs_ad => ufo_radiancerttov_simobs_ad
50  end type ufo_radiancerttov_tlad
51 
52  character(len=maxvarlen), dimension(1), parameter :: varin_default_tlad = &
53  (/var_ts/)
54 
55 contains
56 
57  ! ------------------------------------------------------------------------------
58  subroutine ufo_radiancerttov_tlad_setup(self, f_confOper, channels)
59  implicit none
60 
61  class(ufo_radiancerttov_tlad), intent(inout) :: self
62  type(fckit_configuration), intent(in) :: f_confOper
63  integer(c_int), intent(in) :: channels(:) !List of channels to use
64 
65  type(fckit_configuration) :: f_confOpts ! RTcontrol
66  type(fckit_configuration) :: f_confLinOper
67  integer :: nvars_in
68  integer :: ind, jspec
69 
70  call f_confoper % get_or_die("obs options",f_confopts)
71 
72  call rttov_conf_setup(self % conf_traj, f_confopts,f_confoper)
73 
74  if ( f_confoper%has("linear obs operator") ) then
75  call f_confoper%get_or_die("linear obs operator",f_conflinoper)
76  call rttov_conf_setup(self%conf, f_confopts, f_conflinoper)
77  else
78  call rttov_conf_setup(self%conf, f_confopts, f_confoper)
79  end if
80 
81  !DAR what is the RTTOV equivalant of making sure that humidity and ozone data are present
82  if ( ufo_vars_getindex(self%conf%Absorbers, var_mixr) < 1 .and. &
83  ufo_vars_getindex(self%conf%Absorbers, var_q) < 1 ) then
84  write(message,*) 'ufo_radiancerttov_setup error: H2O must be included in RTTOV Absorbers'
85  call abor1_ftn(message)
86  end if
87 
88  nvars_in = size(varin_default_tlad) + self%conf%ngas + 5 ! 5 near-surface parameters
89 
90  allocate(self%varin(nvars_in))
91  self%varin(1:size(varin_default_tlad)) = varin_default_tlad
92  ind = size(varin_default_tlad) + 1
93 
94  do jspec = 1, self%conf%ngas
95  self%varin(ind) = self%conf%Absorbers(jspec)
96  ind = ind + 1
97  end do
98 
99  self%varin(ind) = var_sfc_t2m
100  ind = ind + 1
101  self%varin(ind) = var_sfc_q2m
102  ind = ind + 1
103  self%varin(ind) = var_sfc_tskin
104  ind = ind + 1
105  self%varin(ind) = var_sfc_u10
106  ind = ind + 1
107  self%varin(ind) = var_sfc_v10
108 
109  ! save channels
110  allocate(self%channels(size(channels)))
111  self%channels(:) = channels(:)
112 
113  end subroutine ufo_radiancerttov_tlad_setup
114 
115  ! ------------------------------------------------------------------------------
117  implicit none
118 
119  class(ufo_radiancerttov_tlad), intent(inout) :: self
120  integer(kind=jpim) :: errorstatus
121 
122  self % ltraj = .false.
123 
124  call rttov_conf_delete(self % conf)
125  call rttov_conf_delete(self % conf_traj)
126 
127  end subroutine ufo_radiancerttov_tlad_delete
128 
129  ! ------------------------------------------------------------------------------
130  subroutine ufo_radiancerttov_tlad_settraj(self, geovals, obss, hofxdiags)
131 
132  use fckit_mpi_module, only: fckit_mpi_comm
133 
134  implicit none
135 
136  class(ufo_radiancerttov_tlad), intent(inout) :: self
137  type(ufo_geovals), intent(in) :: geovals
138  type(c_ptr), value, intent(in) :: obss
139  type(ufo_geovals), intent(inout) :: hofxdiags !non-h(x) diagnostics
140 
141  real(c_double) :: missing
142  type(fckit_mpi_comm) :: f_comm
143 
144  ! Local Variables
145  character(*), parameter :: routine_name = 'ufo_radiancerttov_tlad_settraj'
146  type(rttov_chanprof), allocatable :: chanprof(:)
147  type(ufo_geoval), pointer :: geoval_temp
148 
149  integer(kind=jpim) :: errorstatus ! Return error status of RTTOV subroutine calls
150 
151  integer :: i_inst, nchan_total, ichan, iprof, prof, iprof_rttov
152  integer :: nprof_sim, nprof_max_sim, ichan_sim
153  integer :: prof_start, prof_end
154 
155  logical :: jacobian_needed
156 
157  include 'rttov_k.interface'
158 
159  !Initialisations
160  missing = missing_value(missing)
161 
162  !Return the name and name length of obsspace communicator (from ioda)
163  call obsspace_get_comm(obss, f_comm)
164 
165  !! Parse hofxdiags%variables into independent/dependent variables and channel assumed formats:
166  ! Note this sets the jacobian_needed flag
167  !! jacobian var --> <ystr>_jacobian_<xstr>_<chstr>
168  !! non-jacobian var --> <ystr>_<chstr>
169  call parse_hofxdiags(hofxdiags, jacobian_needed)
170 
171  ! Get number of profiles and levels from geovals
172  self % nprofiles = geovals % nlocs
173  call ufo_geovals_get_var(geovals, var_ts, geoval_temp)
174  self % nlevels = geoval_temp % nval
175  nullify(geoval_temp)
176 
177  ! Allocate RTTOV profiles for ALL geovals for the direct calculation
178  write(message,'(A, A, I0, A, I0, A)') &
179  trim(routine_name), ': Allocating ', self % nprofiles, ' profiles with ', self % nlevels, ' levels'
180  call fckit_log%debug(message)
181  call self % RTprof_K % alloc_profs(errorstatus, self % conf, self % nprofiles, self % nlevels, init=.true., asw=1)
182 
183  !Assign the atmospheric and surface data from the GeoVaLs
184  write(message,'(A, A, I0, A, I0, A)') trim(routine_name), ': Creating RTTOV profiles from geovals'
185  call fckit_log%debug(message)
186  call self % RTprof_K % setup(geovals,obss,self % conf)
187 
188  !DAR: Removing sensor_loop until it's demonstrated to be needed and properly tested
189  ! at the moment self % channels is a single 1D array so cannot adequately contain more than one set of channels
190 ! Sensor_Loop:do i_inst = 1, self % conf % nSensors
191  i_inst = 1
192 
193  ! Number of channels to be simulated for this instrument (from the configuration, not necessarily the full instrument complement)
194  nchan_inst = size(self % channels)
195 
196  ! Allocate memory for *ALL* RTTOV_K channels
197  write(message,'(A,A,I0,A)') &
198  trim(routine_name), ': Allocating Trajectory resources for RTTOV K: ', self % nprofiles * nchan_inst, ' total channels'
199  call self % RTprof_K % alloc_profs_K(errorstatus, self % conf, self % nprofiles * nchan_inst, self % nlevels, init=.true., asw=1)
200 
201  ! Used for keeping track of profiles for setting emissivity
202  allocate(self % RTprof_K % chanprof ( self % nprofiles * nchan_inst ))
203 
204  ! Maximum number of profiles to be processed by RTTOV per pass
205  if(self % conf % prof_by_prof) then
206  nprof_max_sim = 1
207  else
208  nprof_max_sim = max(1,self % conf % nchan_max_sim / nchan_inst)
209  endif
210  nprof_sim = min(nprof_max_sim, self % nprofiles)
211 
212  ! Determine the total number of radiances to simulate (nchan_sim).
213  nchan_sim = nprof_sim * size(self%channels)
214 
215  ! Allocate structures for RTTOV direct and K code
216  write(message,'(A,A,I0,A,I0,A)') &
217  trim(routine_name), ': Allocating resources for RTTOV direct (K): ', nprof_sim, ' and ', nchan_sim, ' channels'
218  call fckit_log%debug(message)
219  call self % RTprof_K % alloc_direct(errorstatus, self % conf, nprof_sim, nchan_sim, self % nlevels, init=.true., asw=1)
220 
221  write(message,'(A,A,I0,A,I0,A)') &
222  trim(routine_name), ': Allocating resources for RTTOV K code: ', nprof_sim, ' and ', nchan_sim, ' channels'
223  call fckit_log%debug(message)
224  call self % RTprof_K % alloc_k(errorstatus, self % conf, nprof_sim, nchan_sim, self % nlevels, init=.true., asw=1)
225 
226  prof_start = 1; prof_end = self % nprofiles
227  nchan_total = 0
228 
229  rttov_loop : do while (prof_start <= prof_end)
230 
231  ! Reduce number of simulated profiles/channel if at end of the of profiles to be processed
232  nprof_sim = min(nprof_sim, prof_end - prof_start + 1)
233  nchan_sim = nprof_sim * size(self%channels)
234 
235  ! allocate and initialise local chanprof structure
236  allocate(chanprof( nchan_sim ))
237  chanprof(:) % prof = 0
238  chanprof(:) % chan = 0
239 
240  ! index for simulated channel
241  ichan_sim = 0_jpim
242 
243  ! Build the list of profile/channel indices in chanprof
244  do iprof_rttov = 1, nprof_sim
245  errorstatus = errorstatus_success
246 
247  ! iprof is the index for the full set of RTTOV profiles
248  iprof = prof_start + iprof_rttov - 1
249 
250  ! print profile information if requested
251  if(any(self % conf % inspect == iprof)) call self % RTprof_K % print(self % conf, iprof, i_inst)
252 
253  ! check RTTOV profile and flag it if it fails the check
254  if(self % conf % RTTOV_profile_checkinput) call self % RTprof_K % check(self % conf, iprof, i_inst, errorstatus)
255 
256  if (errorstatus == errorstatus_success) then
257  do ichan = 1, nchan_inst
258  ichan_sim = ichan_sim + 1_jpim
259  chanprof(ichan_sim) % prof = iprof_rttov ! this refers to the slice of the RTprofile array passed to RTTOV
260  chanprof(ichan_sim) % chan = self % channels(ichan)
261  self % RTprof_K % chanprof(nchan_total + ichan_sim) % prof = iprof ! this refers to the index of the profile from the geoval
262  self % RTprof_K % chanprof(nchan_total + ichan_sim) % chan = self % channels(ichan)
263  end do
264  nchan_sim = ichan_sim
265  endif
266 
267  end do
268 
269  ! Set surface emissivity
270  call self % RTProf_K % init_emissivity(self % conf, prof_start)
271 
272  ! --------------------------------------------------------------------------
273  ! Call RTTOV K model
274  ! --------------------------------------------------------------------------
275 
276  call rttov_k( &
277  errorstatus, &! out error flag
278  chanprof(1:nchan_sim), &! in channel and profile index structure
279  self % conf % rttov_opts, &! in options structure
280  self % RTprof_K % profiles(prof_start:prof_start + nprof_sim - 1), &! in profile array
281  self % RTprof_K % profiles_k(nchan_total + 1 : nchan_total + nchan_sim), &! in profile array
282  self % conf % rttov_coef_array(i_inst), &! in coefficients structure
283  self % RTprof_K % transmission, &! inout computed transmittances
284  self % RTprof_K % transmission_k, &! inout computed transmittances
285  self % RTprof_K % radiance, &! inout computed radiances
286  self % RTprof_K % radiance_k, &! inout computed radiances
287  calcemis = self % RTprof_K % calcemis(1:nchan_sim), &! in flag for internal emissivity calcs
288  emissivity = self % RTprof_K % emissivity(1:nchan_sim), &!, &! inout input/output emissivities per channel
289  emissivity_k = self % RTprof_K % emissivity_k(1:nchan_sim))!, &! inout input/output emissivities per channel
290 
291  if ( errorstatus /= errorstatus_success ) then
292  write(message,'(A, A, 2I6)') trim(routine_name), 'after rttov_k: error ', errorstatus, i_inst, &
293  ' skipping profiles ', prof_start, ' -- ', prof_start + nprof_sim - 1
294  call fckit_log%info(message)
295  else
296  ! Put simulated diagnostics into hofxdiags
297  ! ----------------------------------------------
298  if(hofxdiags%nvar > 0) call populate_hofxdiags(self % RTprof_K, chanprof, self % conf, prof_start, hofxdiags)
299  end if
300 
301  ! increment profile and channel counters
302  nchan_total = nchan_total + nchan_sim
303  prof_start = prof_start + nprof_sim
304 
305  self % nchan_total = nchan_total
306  deallocate (chanprof)
307  end do rttov_loop
308 
309  ! end do Sensor_Loop
310  ! Deallocate structures for rttov_direct
311  call self % RTprof_K % alloc_k(errorstatus, self % conf, -1, -1, -1, asw=0)
312  call self % RTprof_K % alloc_direct(errorstatus, self % conf, -1, -1, -1, asw=0)
313  call self % RTprof_K % alloc_profs(errorstatus, self % conf, -1, -1, asw=0)
314 
315 
316  ! Set flag that the tracectory was set
317  ! ------------------------------------
318  self % ltraj = .true.
319 
320 end subroutine ufo_radiancerttov_tlad_settraj
321 
322  ! ------------------------------------------------------------------------------
323  subroutine ufo_radiancerttov_simobs_tl(self, geovals, obss, nvars, nlocs, hofx)
324 
325  use ufo_constants_mod, only : zero, g_to_kg
326 
327  implicit none
328 
329  class(ufo_radiancerttov_tlad), intent(in) :: self
330  type(ufo_geovals), intent(in) :: geovals
331  type(c_ptr), value, intent(in) :: obss
332  integer, intent(in) :: nvars, nlocs
333  real(c_double), intent(inout) :: hofx(nvars, nlocs)
334 
335  character(len=*), parameter :: myname_="ufo_radiancerttov_simobs_tl"
336  integer :: ichan, jchan, prof, jspec
337 
338  type(ufo_geoval), pointer :: geoval_d, geoval_d2
339 
340  ! Initial checks
341  ! --------------
342 
343  ! Check if trajectory was set
344  if (.not. self % ltraj) then
345  write(message,*) myname_, ' trajectory wasnt set!'
346  call abor1_ftn(message)
347  end if
348 
349  ! Check if nlocs is consistent in geovals & hofx
350  if (geovals % nlocs /= self % nprofiles) then
351  write(message,*) myname_, ' error: nlocs inconsistent!'
352  call abor1_ftn(message)
353  end if
354 
355  ! Initialize hofx
356  ! ---------------
357  hofx(:,:) = zero
358 
359  ! Temperature
360  ! -----------
361  call ufo_geovals_get_var(geovals, var_ts, geoval_d) ! var_ts = air_temperature
362 
363  ! Check model levels is consistent in geovals
364  if (geoval_d % nval /= self % nlevels) then
365  write(message,*) myname_, ' error: layers inconsistent!'
366  call abor1_ftn(message)
367 
368  end if
369 
370  do ichan = 1, self % nchan_total, size(self%channels)
371  prof = self % RTprof_K % chanprof(ichan) % prof
372  do jchan = 1, size(self%channels)
373  hofx(jchan,prof) = hofx(jchan,prof) + &
374  sum(self % RTprof_K % profiles_k(ichan+jchan-1) % t(self % nlevels:1:-1) * geoval_d % vals(1:geoval_d % nval,prof))
375  enddo
376  end do
377 
378  do jspec = 1, self%conf%ngas
379  call ufo_geovals_get_var(geovals, self%conf%Absorbers(jspec), geoval_d)
380 
381  ! Check model levels is consistent in geovals
382  if (geoval_d % nval /= self % nlevels) then
383  write(message,*) myname_, ' error: layers inconsistent!'
384  call abor1_ftn(message)
385  end if
386 
387  ! Absorbers
388  ! ---------
389  ! This is where CO2 and friends will live as well as CLW
390  do ichan = 1, self % nchan_total, size(self%channels)
391  prof = self % RTprof_K % chanprof(ichan) % prof
392  do jchan = 1, size(self%channels)
393  if(self%conf%Absorbers(jspec) == var_q) then
394  hofx(jchan,prof) = hofx(jchan,prof) + &
395  sum(self % RTprof_K % profiles_k(ichan+jchan-1) % q(self % nlevels:1:-1) * &
396  geoval_d % vals(1:geoval_d % nval,prof)) * self%conf%scale_fac(gas_id_watervapour)
397  elseif(self%conf%Absorbers(jspec) == var_mixr) then
398  hofx(jchan,prof) = hofx(jchan,prof) + &
399  sum(self % RTprof_K % profiles_k(ichan+jchan-1) % q(self % nlevels:1:-1) * &
400  geoval_d % vals(1:geoval_d % nval,prof)) * self%conf%scale_fac(gas_id_watervapour) / g_to_kg
401  elseif(self%conf%Absorbers(jspec) == var_clw) then
402  hofx(jchan,prof) = hofx(jchan,prof) + &
403  sum(self % RTprof_K % profiles_k(ichan+jchan-1) % clw(self % nlevels:1:-1) * &
404  geoval_d % vals(1:geoval_d % nval,prof))
405  endif
406  enddo
407  end do
408  enddo
409 
410  ! Cloud
411  ! --------------------------
412  !IR
413 
414 
415  ! Surface + Single-valued Variables
416  ! --------------------------
417  !CTP/cloudfrac
418  !O3total
419  !LWP
420 
421  !T2m
422  call ufo_geovals_get_var(geovals, var_sfc_t2m, geoval_d)
423 
424  do ichan = 1, self % nchan_total, size(self%channels)
425  prof = self % RTprof_K % chanprof(ichan) % prof
426  do jchan = 1, size(self%channels)
427  hofx(jchan,prof) = hofx(jchan,prof) + &
428  self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % t * geoval_d % vals(1,prof)
429  enddo
430  end do
431 
432  !q2m
433  call ufo_geovals_get_var(geovals, var_sfc_q2m, geoval_d)
434 
435  do ichan = 1, self % nchan_total, size(self%channels)
436  prof = self % RTprof_K % chanprof(ichan) % prof
437  do jchan = 1, size(self%channels)
438  hofx(jchan,prof) = hofx(jchan,prof) + &
439  self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % q * &
440  geoval_d % vals(1,prof) * self%conf%scale_fac(gas_id_watervapour)
441  enddo
442  end do
443 
444  !windspeed
445  call ufo_geovals_get_var(geovals, var_sfc_u10, geoval_d)
446  call ufo_geovals_get_var(geovals, var_sfc_v10, geoval_d2)
447 
448  do ichan = 1, self % nchan_total, size(self%channels)
449  prof = self % RTprof_K % chanprof(ichan) % prof
450  do jchan = 1, size(self%channels)
451  hofx(jchan,prof) = hofx(jchan,prof) + &
452  self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % u * geoval_d % vals(1,prof) + &
453  self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % v * geoval_d2 % vals(1,prof)
454  enddo
455  end do
456 
457  !Tskin
458  call ufo_geovals_get_var(geovals, var_sfc_tskin, geoval_d)
459 
460  do ichan = 1, self % nchan_total, size(self%channels)
461  prof = self % RTprof_K % chanprof(ichan) % prof
462  do jchan = 1, size(self%channels)
463  hofx(jchan,prof) = hofx(jchan,prof) + &
464  self % RTprof_K % profiles_k(ichan+jchan-1) % skin % t * geoval_d % vals(1,prof)
465  enddo
466  end do
467 
468  end subroutine ufo_radiancerttov_simobs_tl
469 
470  ! ------------------------------------------------------------------------------
471  subroutine ufo_radiancerttov_simobs_ad(self, geovals, obss, nvars, nlocs, hofx)
472 
473  use ufo_constants_mod, only : zero, g_to_kg
474 
475  implicit none
476 
477  class(ufo_radiancerttov_tlad), intent(in) :: self
478  type(ufo_geovals), intent(inout) :: geovals
479  type(c_ptr), value, intent(in) :: obss
480  integer, intent(in) :: nvars, nlocs
481  real(c_double), intent(in) :: hofx(nvars, nlocs)
482 
483  type(ufo_geoval), pointer :: geoval_d, geoval_d2
484 
485  real(c_double) :: missing
486  integer :: ichan, jchan, prof, jspec
487 
488  character(len=*), parameter :: myname_ = "ufo_radiancerttov_simobs_ad"
489 
490  ! Set missing value
491  missing = missing_value(missing)
492 
493  ! Initial checks
494  ! --------------
495 
496  ! Check if trajectory was set
497  if (.not. self % ltraj) then
498  write(message,*) myname_, ' trajectory wasnt set!'
499  call abor1_ftn(message)
500  end if
501 
502  ! Check if nlocs is consistent in geovals & hofx
503  if (geovals % nlocs /= self % nprofiles) then
504  write(message,*) myname_, ' error: nlocs inconsistent!'
505  call abor1_ftn(message)
506  end if
507 
508  ! Temperature
509  ! -----------
510  call ufo_geovals_get_var(geovals, var_ts, geoval_d) ! var_ts = air_temperature
511 
512  do ichan = 1, self % nchan_total, size(self%channels)
513  prof = self % RTprof_K % chanprof(ichan) % prof
514  do jchan = 1, size(self%channels)
515  if (hofx(jchan, prof) /= missing) then
516  geoval_d % vals(:,prof) = geoval_d % vals(:,prof) + &
517  self % RTprof_K % profiles_k(ichan+jchan-1) % t(self % nlevels:1:-1) * hofx(jchan,prof)
518  endif
519  enddo
520  end do
521 
522  ! Absorbers
523  ! ---------
524  ! This is where CO2 and friends will live as well as CLW
525 
526  do jspec = 1, self%conf%ngas
527  call ufo_geovals_get_var(geovals, self%conf%Absorbers(jspec), geoval_d)
528 
529  do ichan = 1, self % nchan_total, size(self%channels)
530  prof = self % RTprof_K % chanprof(ichan) % prof
531  do jchan = 1, size(self%channels)
532  if (hofx(jchan, prof) /= missing) then
533 
534  if(self%conf%Absorbers(jspec) == var_q) then
535  geoval_d % vals(:,prof) = geoval_d % vals(:,prof) + &
536  self % RTprof_K % profiles_k(ichan+jchan-1) % q(self % nlevels:1:-1) * hofx(jchan,prof) * &
537  self%conf%scale_fac(gas_id_watervapour)
538  elseif(self%conf%Absorbers(jspec) == var_mixr) then
539  geoval_d % vals(:,prof) = geoval_d % vals(:,prof) + &
540  (self % RTprof_K % profiles_k(ichan+jchan-1) % q(self % nlevels:1:-1) * hofx(jchan,prof)) * &
541  self%conf%scale_fac(gas_id_watervapour) / g_to_kg
542  elseif(self%conf%Absorbers(jspec) == var_clw) then
543  geoval_d % vals(:,prof) = geoval_d % vals(:,prof) + &
544  self % RTprof_K % profiles_k(ichan+jchan-1) % clw(self % nlevels:1:-1) * hofx(jchan,prof)
545  endif
546  endif
547  enddo
548  enddo
549  enddo
550 
551  ! Cloud
552  ! --------------------------
553  !IR
554 
555  ! Surface + Single-valued Variables
556  ! --------------------------
557  !CTP/cloudfrac
558  !O3total
559  !LWP
560  !T2m
561 
562  call ufo_geovals_get_var(geovals, var_sfc_t2m, geoval_d)
563 
564  do ichan = 1, self % nchan_total, size(self%channels)
565  prof = self % RTprof_K % chanprof(ichan) % prof
566  do jchan = 1, size(self%channels)
567  if (hofx(jchan, prof) /= missing) then
568  geoval_d % vals(1,prof) = geoval_d % vals(1,prof) + &
569  self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % t * hofx(jchan,prof)
570  endif
571  enddo
572  end do
573 
574  !q2m
575  call ufo_geovals_get_var(geovals, var_sfc_q2m, geoval_d)
576 
577  do ichan = 1, self % nchan_total, size(self%channels)
578  prof = self % RTprof_K % chanprof(ichan) % prof
579  do jchan = 1, size(self%channels)
580  if (hofx(jchan, prof) /= missing) then
581  geoval_d % vals(1,prof) = geoval_d % vals(1,prof) + &
582  self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % q * &
583  hofx(jchan,prof) * self%conf%scale_fac(gas_id_watervapour)
584  endif
585  enddo
586  end do
587 
588  !windspeed
589  call ufo_geovals_get_var(geovals, var_sfc_u10, geoval_d)
590  call ufo_geovals_get_var(geovals, var_sfc_v10, geoval_d2)
591 
592  do ichan = 1, self % nchan_total, size(self%channels)
593  prof = self % RTprof_K % chanprof(ichan) % prof
594  do jchan = 1, size(self%channels)
595  if (hofx(jchan, prof) /= missing) then
596  geoval_d % vals(1,prof) = geoval_d % vals(1,prof) + &
597  self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % u * hofx(jchan,prof)
598 
599  geoval_d2 % vals(1,prof) = geoval_d2 % vals(1,prof) + &
600  self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % v * hofx(jchan,prof)
601  endif
602  enddo
603  end do
604 
605  !Tskin
606  call ufo_geovals_get_var(geovals, var_sfc_tskin, geoval_d)
607 
608  do ichan = 1, self % nchan_total, size(self%channels)
609  prof = self % RTprof_K % chanprof(ichan) % prof
610  do jchan = 1, size(self%channels)
611  if (hofx(jchan, prof) /= missing) then
612  geoval_d % vals(1,prof) = geoval_d % vals(1,prof) + &
613  self % RTprof_K % profiles_k(ichan+jchan-1) % skin % t * hofx(jchan,prof)
614  endif
615  enddo
616  end do
617 
618  ! Once all geovals set replace flag
619  ! ---------------------------------
620  if (.not. geovals % linit ) geovals % linit=.true.
621 
622  end subroutine ufo_radiancerttov_simobs_ad
623 
624  ! ------------------------------------------------------------------------------
625 
real(kind_real), parameter, public g_to_kg
real(kind_real), parameter, public zero
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Fortran module for radiancerttov tl/ad observation operator.
subroutine ufo_radiancerttov_tlad_setup(self, f_confOper, channels)
subroutine ufo_radiancerttov_tlad_settraj(self, geovals, obss, hofxdiags)
subroutine ufo_radiancerttov_simobs_ad(self, geovals, obss, nvars, nlocs, hofx)
character(len=maxvarlen), dimension(1), parameter varin_default_tlad
subroutine ufo_radiancerttov_simobs_tl(self, geovals, obss, nvars, nlocs, hofx)
Fortran module to provide code shared between nonlinear and tlm/adm radiance calculations.
subroutine, public parse_hofxdiags(hofxdiags, jacobian_needed)
subroutine, public rttov_conf_setup(conf, f_confOpts, f_confOper)
character(len=max_string), public message
subroutine, public rttov_conf_delete(conf)
subroutine, public populate_hofxdiags(RTProf, chanprof, conf, prof_start, hofxdiags)
integer function, public ufo_vars_getindex(vars, varname)
character(len=maxvarlen), parameter, public var_clw
character(len=maxvarlen), parameter, public var_sfc_v10
character(len=maxvarlen), parameter, public var_q
character(len=maxvarlen), parameter, public var_sfc_u10
character(len=maxvarlen), parameter, public var_sfc_q2m
character(len=maxvarlen), parameter, public var_mixr
character(len=maxvarlen), parameter, public var_sfc_tskin
character(len=maxvarlen), parameter, public var_ts
character(len=maxvarlen), parameter, public var_sfc_t2m
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators
Fortran derived type for radiancerttov trajectory.