UFO
ufo_radiancerttov_utils_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 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 to provide code shared between nonlinear and tlm/adm radiance calculations
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 rttov_types, only : rttov_options, rttov_profile, rttov_coefs, &
17  rttov_radiance, rttov_transmission, rttov_emissivity, &
18  rttov_chanprof
19 
20  use rttov_const ! gas_ids
21 
22  use ufo_vars_mod
24  use ufo_basis_mod, only: ufo_basis
25  use obsspace_mod
26 
27  implicit none
28  private
29 
30  public rttov_conf
31  public rttov_conf_setup
32  public rttov_conf_delete
33  public load_atm_data_rttov
35  public parse_hofxdiags
36  public populate_hofxdiags
38 
39  integer, parameter, public :: max_string=800
40  integer, parameter, public :: maxvarin = 50
41 
42  !DARFIX this should go somewhere else
43 
44  character(len=maxvarlen), public :: varin_temp(maxvarin)
45  character(len=max_string), public :: message
46 
47  integer, public :: nvars_in
48  integer, public :: rttov_errorstatus
49 
50  character(len=MAXVARLEN), allocatable :: ystr_diags(:), xstr_diags(:)
51  integer, allocatable :: ch_diags(:)
52 
53  real(c_double) :: missing
54 
55  character(len=maxvarlen), dimension(21), target :: varin_default_crtm = &
56  (/var_ts, var_prs, var_prsi, &
62 
63  !var_ps
64  character(len=maxvarlen), dimension(11), target :: varin_default_satrad = &
68 
69  character(len=maxvarlen), pointer, public :: varin_default(:)
70 
71  ! copy of ABSORBER_ID_NAME defined in rttov_const
72  character(len=*), parameter :: &
73  rttov_absorbers(ngases_max+1) = &
74  [gas_name(1:ngases_max),'CLW']
75 
76  integer, parameter :: &
77  rttov_absorber_id(ngases_max+1) = &
78  [gas_id_mixed, &
79  gas_id_watervapour, &
80  gas_id_ozone, &
81  gas_id_wvcont, &
82  gas_id_co2, &
83  gas_id_n2o, &
84  gas_id_co, &
85  gas_id_ch4, &
86  gas_id_so2, 99]
87 
88  character(len=MAXVARLEN), parameter :: null_str = ''
89 
90  !DARFIX: need to get correct names (correct units for RTTOV) in ufo_vars_mod
91  character(len=MAXVARLEN), parameter :: &
92  ufo_absorbers(ngases_max+1) = &
93  [ null_str, var_mixr, var_oz, null_str, var_co2, 'mole_fraction_of_nitrous_oxide_in_air', &
94  'mole_fraction_of_carbon_monoxide_in_air', 'mole_fraction_of_methane_in_air', &
95  'mole_fraction_of_sulfur_dioxide_in_air', var_clw]
96 
97  integer, public :: nchan_inst ! number of channels being simulated (may be less than full instrument)
98  integer, public :: nchan_sim ! total number of 'obs' = nprofiles * nchannels
99  integer, public :: nlocs_total ! nprofiles (including skipped)
100  logical, public :: debug
101 !Common counters
102  integer :: iprof
103 
104  type, public :: ufo_rttov_io
105  logical, pointer :: calcemis(:) ! Flag to indicate calculation of emissivity within RTTOV
106 
107  type(rttov_emissivity), pointer :: emissivity(:) ! Input/output surface emissivity
108  type(rttov_profile), pointer :: profiles(:) ! Input profiles
109  type(rttov_profile), pointer :: profiles_k(:) ! Input profiles
110  type(rttov_chanprof), pointer :: chanprof(:) ! Input profiles
111  type(rttov_transmission) :: transmission ! Output transmittances
112  type(rttov_radiance) :: radiance ! Output radiances
113 
114  type(rttov_emissivity), pointer :: emissivity_k(:) !Input/output surface emissivity
115  type(rttov_transmission) :: transmission_k ! Output transmittances
116  type(rttov_radiance) :: radiance_k ! Output radiances
117 
118  contains
119 
120  procedure :: alloc => alloc_rttov
121  procedure :: init_emissivity => ufo_rttov_init_emissivity
122 
123  end type ufo_rttov_io
124 
125  !Type for general config
127  integer :: nsensors
128  integer :: ngas
129 
130  character(len=MAXVARLEN), allocatable :: absorbers(:)
131  integer, allocatable :: absorber_id(:)
132 
133  character(len=255), allocatable :: sensor_id(:)
134  character(len=255) :: coefficient_path
135 
136  type(rttov_coefs), allocatable :: rttov_coef_array(:)
137  character(len=8) :: rttov_default_opts = 'RTTOV'
138  type(rttov_options) :: rttov_opts
139  logical :: rttov_is_setup = .false.
140 
141  logical :: satrad_compatibility = .true.
142  logical :: useqtsplitrain, splitqtotal = .false. ! true for MW, false otherwise
143  logical :: rttov_profile_checkinput = .false.
144 
145  integer :: inspect
146  integer :: nchan_max_sim
147 
148  contains
149 
150  procedure :: set_options => set_options_rttov
151  procedure :: setup => setup_rttov
152  procedure :: set_defaults => set_defaults_rttov
153 
154  end type rttov_conf
155 
156 contains
157 
158  ! ------------------------------------------------------------------------------
159 
160  subroutine rttov_conf_setup(conf, f_confOpts, f_confOper)
161  implicit none
162 
163  type(rttov_conf), intent(inout) :: conf
164  type(fckit_configuration), intent(in) :: f_confopts ! RTcontrol
165  type(fckit_configuration), intent(in) :: f_confoper ! what is this
166 
167  character(*), parameter :: routine_name = 'rttov_conf_setup'
168  integer :: ivar, jspec
169  character(len=:), allocatable :: str
170  character(len=:), allocatable :: str_array(:)
171  logical :: varin_satrad = .false.
172 
173  !Number of sensors, each call to RTTOV will be for a single sensor
174  !type (zenith/scan angle will be different)
175  conf % nSensors = 1
176 
177  if (f_confoper%has("Debug")) then
178  call f_confoper % get_or_die("Debug",debug)
179  else
180  debug = .false. ! default
181  endif
182 
183  if (f_confoper%has("GeoVal_type")) then
184  call f_confoper%get_or_die("GeoVal_type",str)
185  if (trim(str) == 'MetO' .or. trim(str) == 'SatRad') then
187  varin_satrad = .true.
188  elseif(trim(str) == 'CRTM') then
190  varin_satrad = .false.
191  else
192  write(message,*) trim(routine_name),' error: ',trim(str),' is not a supported GeoVal type'
193  call abor1_ftn(message)
194  endif
195  else
197  endif
198 
199  ! Absorbers
200  !----------
201  conf%ngas = 0
202  if (f_confoper%has("Absorbers")) &
203  conf%ngas = conf%ngas + f_confoper%get_size("Absorbers")
204 
205  allocate(conf%Absorbers( conf%ngas ), &
206  conf%Absorber_Id( conf%ngas ))
207 
208  if (conf%ngas > 0) then
209  call f_confoper%get_or_die("Absorbers",str_array)
210  conf%Absorbers(1:conf%ngas) = str_array
211 
212  end if
213 
214  ! check for duplications
215  do jspec = 2, conf%ngas
216  if ( any(conf%Absorbers(jspec-1) == conf%Absorbers(jspec:conf%ngas)) ) then
217  write(message,*) trim(routine_name),' error: ',trim(conf%Absorbers(jspec)),' is duplicated in Absorbers'
218  call abor1_ftn(message)
219  end if
220  end do
221 
222  ! convert from CRTM names to UFO CF names and define Id and Units
223  do jspec = 1, conf%ngas
224  ivar = ufo_vars_getindex(rttov_absorbers, conf%Absorbers(jspec))
225 
226  if (ivar < 1 .or. ivar > size(ufo_absorbers)) then
227  write(message,*) trim(routine_name),' error: ',trim(conf%Absorbers(jspec)),' not supported by UFO_Absorbers'
228  call abor1_ftn(message)
229  end if
230 
231  conf%Absorbers(jspec) = ufo_absorbers(ivar)
232 
233  ! DARFIX replace humidity_mixing_ratio with specific_humidity
234  ! this is starting to get a bit messy and dangerous
235  if(conf%Absorbers(jspec) == var_mixr .and. varin_satrad) then
236  conf%Absorbers(jspec) = var_q
237  endif
238  conf%Absorber_Id(jspec) = rttov_absorber_id(ivar)
239  end do
240 
241 
242  ! Allocate SENSOR_ID
243  allocate(conf % SENSOR_ID(conf % nSensors))
244 
245  ! Get sensor ID from config
246  call f_confopts % get_or_die("Sensor_ID",str)
247  conf % SENSOR_ID(conf%nSensors) = str
248 
249  ! Path to coefficient files
250  call f_confopts % get_or_die("CoefficientPath",str)
251  conf % COEFFICIENT_PATH = str
252 
253  if(f_confopts % has("RTTOV_default_opts")) then
254  call f_confopts % get_or_die("RTTOV_default_opts",str) !test, OPS, RTTOV
255  conf % RTTOV_default_opts = str
256  endif
257 
258  if(f_confopts % has("SatRad_compatibility")) then
259  call f_confopts % get_or_die("SatRad_compatibility",conf % SatRad_compatibility)
260  endif
261 
262  if(f_confopts % has("max_channels_per_batch")) then
263  call f_confopts % get_or_die("max_channels_per_batch",conf % nchan_max_sim)
264  else
265  conf % nchan_max_sim = 10000
266  endif
267 
268  if( .not. conf % rttov_is_setup) then
269  call conf % setup(f_confopts, asw=1)
270  end if
271 
272  !DARFIX THIS ONLY WORKS FOR ONE INSTRUMENT
273  !
274  if (conf % rttov_coef_array(1) % coef % id_sensor == sensor_id_mw) then
275  if(conf % rttov_opts % rt_mw % clw_data .and. conf % SatRad_compatibility) then
276  conf % UseQtsplitRain = .true.
277  conf % splitQtotal = .true.
278  endif
279 
280  conf % rttov_opts % rt_ir % ozone_data = .false.
281  conf % rttov_opts % rt_ir % co2_data = .false.
282  conf % rttov_opts % rt_ir % n2o_data = .false.
283  conf % rttov_opts % rt_ir % ch4_data = .false.
284  conf % rttov_opts % rt_ir % so2_data = .false.
285  else
286  conf % UseQtsplitRain = .false.
287  conf % splitQtotal = .false.
288  endif
289 
290  if(f_confopts % has("RTTOV_profile_checkinput")) then
291  call f_confopts % get_or_die("RTTOV_profile_checkinput",conf % RTTOV_profile_checkinput) !test, OPS, RTTOV
292  endif
293 
294  conf % inspect = 0
295  if (f_confopts%has("InspectProfileNumber")) then
296  call f_confopts % get_or_die("InspectProfileNumber",conf % inspect)
297  conf % RTTOV_profile_checkinput = .true.
298  endif
299 
300  end subroutine rttov_conf_setup
301 
302  ! -----------------------------------------------------------------------------
303 
304  subroutine rttov_conf_delete(conf)
305 
306  implicit none
307  type(rttov_conf), intent(inout) :: conf
308 
309  deallocate(conf%SENSOR_ID)
310  deallocate(conf%Absorbers)
311  deallocate(conf%Absorber_Id)
312 
313  end subroutine rttov_conf_delete
314 
315  ! -----------------------------------------------------------------------------
316 
317  subroutine load_atm_data_rttov(geovals,obss,profiles,prof_start,conf,layer_quantities,obs_info)
318 
320 
321  implicit none
322 
323  type(ufo_geovals), intent(in) :: geovals
324  type(c_ptr), value, intent(in) :: obss
325  type(rttov_profile), intent(inout) :: profiles(:)
326  integer, intent(in) :: prof_start
327  type(rttov_conf), intent(in) :: conf
328  logical, intent(inout) :: layer_quantities
329  logical, optional, intent(in) :: obs_info
330 
331  ! Local variables
332  integer :: jspec
333  integer :: nlevels
334  integer :: nprofiles
335 
336  type(ufo_geoval), pointer :: geoval
337  character(MAXVARLEN) :: varname
338 
339  real(kind_real) :: ifrac, sfrac, lfrac, wfrac
340  real(kind_real) :: itmp, stmp, ltmp
341 
342  real(kind_real), allocatable :: tmpvar(:), windsp(:), p(:)
343  logical :: variable_present
344 
345  integer :: top_level, bottom_level, stride
346  real(kind_real) :: tstar, newt
347  integer :: level_1000hpa, level_950hpa
348 
349  real(kind_real), allocatable :: q_temp(:), clw_temp(:), ciw_temp(:), qtotal(:)
350 
351  if(present(obs_info)) then
352  nlocs_total = 1
353  else
354  nlocs_total = obsspace_get_nlocs(obss)
355  end if
356 
357  nprofiles = min(size(profiles), geovals%nlocs - prof_start + 1)
358 
359  nlevels = size(profiles(1)%p)
360 
361  ! Assume that the pressure profile coming from the geovals is increasing in pressure (ToA->surface)...
362  top_level = 1
363  bottom_level = nlevels
364  stride=1
365 
366  if (ufo_vars_getindex(geovals%variables, var_prs) > 0) then
367  call ufo_geovals_get_var(geovals, var_prs, geoval)
368 
369  allocate(p(nlevels))
370  do iprof = 1, nprofiles
371  p = geoval%vals(geoval%nval-nlevels+1:geoval%nval,prof_start + iprof - 1) * pa_to_hpa ! for RTTOV
372  if (iprof == 1) then
373  if (p(1) > p(2)) then ! ...but be ready to switch. Assume if one profile is 'upside-down' then they all are
374  top_level = nlevels
375  bottom_level = 1
376  stride = -1
377  endif
378  endif
379  profiles(iprof)%p(top_level:bottom_level:stride) = p(:)
380  end do
381  deallocate(p)
382 
383  ! DARFIX - it might be that all of the temperatures are on layers so it may be that we need to do interpolation
384  ! but it might be better to force this using an option if we know where the geovals are coming from
385  else
386  call ufo_geovals_get_var(geovals, var_prsi, geoval)
387  allocate(p(nlevels-1))
388  do iprof = 1, nprofiles
389  p = geoval%vals(geoval%nval-(nlevels-1)+1:geoval%nval,prof_start + iprof - 1) * pa_to_hpa
390  if (iprof == 1) then
391  if (p(1) > p(2)) then !upside-down
392  top_level = nlevels - 1
393  bottom_level = 1
394  stride = -1
395  else
396  top_level = 1
397  bottom_level = nlevels - 1
398  stride = 1
399  endif
400  endif
401  profiles(iprof)%p(top_level+stride:bottom_level:stride) = half * (p(top_level:bottom_level-stride:stride) + p(top_level+stride:bottom_level:stride))
402  profiles(iprof)%p(1) = max( profiles(iprof)%p(2) - half * (profiles(iprof)%p(3) - profiles(iprof)%p(2)),half * profiles(iprof)%p(2))
403  profiles(iprof)%p(nlevels) = profiles(iprof)%p(nlevels-1) - half * (profiles(iprof)%p(nlevels-2) - profiles(iprof)%p(nlevels-1))
404  end do
405  deallocate(p)
406  endif
407 
408 ! Get temperature
409  varname = var_ts
410  call ufo_geovals_get_var(geovals, varname, geoval)
411 
412 ! Check if temperatures are provided on levels as required for RTTOV, otherwise assume that temperatures are layer quantities
413 ! (and all future atmospheric variables) and do some interpolation to prepare for RTTOV.
414 ! TODO: Put a warning in here that this is happening
415  if( size(geoval%vals(:,1)) == nlevels) then
416  do iprof = 1, nprofiles
417  profiles(iprof)%t(top_level:bottom_level:stride) = geoval%vals(:,prof_start + iprof - 1) ! K
418  end do
419  else
420  layer_quantities = .true.
421  bottom_level = nlevels-1
422 
423  do iprof = 1, nprofiles
424  profiles(iprof)%t(top_level+stride:bottom_level:stride) = half * (profiles(iprof)%t(top_level:bottom_level-stride:stride) + &
425  profiles(iprof)%t(top_level+stride:bottom_level:stride))
426  profiles(iprof)%t(1) = profiles(iprof)%t(2)
427  profiles(iprof)%t(nlevels) = profiles(iprof)%t(nlevels-1) - half * (profiles(iprof)%t(nlevels-2) - profiles(iprof)%t(nlevels-1))
428  end do
429  endif
430 
431 ! Get absorbers. Assume mass mixing ratio (moist air). The distinction has been made that q will be in kg/kg as per OPS convention
432 ! And mixr will be in g/kg.
433 ! For now all other gases will be in kg/kg and this needs to be handled carefully.
434 
435  profiles(1:nprofiles)%gas_units = 1
436 
437  do jspec = 1, conf%ngas
438  call ufo_geovals_get_var(geovals,conf%Absorbers(jspec) , geoval)
439 
440  select case (conf%Absorbers(jspec))
441  case (var_mixr)
442  do iprof = 1, nprofiles
443  profiles(iprof)%q(top_level:bottom_level:stride) = geoval%vals(:,prof_start + iprof - 1) * g_to_kg
444  end do
445  case (var_q)
446  do iprof = 1, nprofiles
447  profiles(iprof)%q(top_level:bottom_level:stride) = geoval%vals(:,prof_start + iprof - 1) ! kg/kg
448  end do
449  case (var_oz)
450  if (associated(profiles(1)%o3)) then
451  do iprof = 1, nprofiles
452  profiles(iprof)%o3(top_level:bottom_level:stride) = geoval%vals(:,prof_start + iprof - 1) ! kg/kg
453  end do
454  endif
455  case (var_co2)
456  if (associated(profiles(1)%co2)) then
457  do iprof = 1, nprofiles
458  profiles(iprof)%co2(top_level:bottom_level:stride) = geoval%vals(:,prof_start + iprof - 1) ! kg/kg
459  end do
460  endif
461  case (var_clw)
462  if (associated(profiles(1)%clw)) then
463  do iprof = 1, nprofiles
464  profiles(iprof)%clw(top_level:bottom_level:stride) = geoval%vals(:,prof_start + iprof - 1) ! kg/kg
465  end do
466  endif
467  case default
468 
469  end select
470 
471  end do
472 
473 ! Get near-surface variables (s2m)
474 
475  varname = var_sfc_p2m
476  if (ufo_vars_getindex(geovals%variables, varname) > 0) then
477  call ufo_geovals_get_var(geovals, varname, geoval)
478 
479  profiles(1:nprofiles)%s2m%p = geoval%vals(1,prof_start:prof_start + nprofiles - 1) * pa_to_hpa
480  else
481  write(message,'(A)') 'No near-surface pressure. Using bottom pressure level'
482  call fckit_log%info(message)
483 
484  do iprof = 1, nprofiles
485  profiles(iprof)%s2m%p = profiles(iprof)%p(nlevels)
486  enddo
487  endif
488 
489  varname = var_sfc_t2m ! 2m temperature
490  if (ufo_vars_getindex(geovals%variables, varname) > 0) then
491  call ufo_geovals_get_var(geovals, varname, geoval)
492  profiles(1:nprofiles)%s2m%t = geoval%vals(1,prof_start:prof_start + nprofiles - 1)
493  else
494  write(message,'(A)') 'No near-surface temperature. Using bottom temperature level'
495  call fckit_log%info(message)
496 
497  do iprof = 1, nprofiles
498  profiles(iprof)%s2m%t = profiles(iprof)%t(nlevels)
499  enddo
500  endif
501 
502  varname = var_sfc_q2m ! 2m specific humidity (kg/kg)
503  if (ufo_vars_getindex(geovals%variables, varname) > 0) then
504  call ufo_geovals_get_var(geovals, varname, geoval) ! lfric
505 
506  profiles(1:nprofiles)%s2m%q = geoval%vals(1,prof_start:prof_start + nprofiles - 1)
507  else
508  write(message,'(A)') 'No near-surface specific humidity. Using bottom q level'
509  call fckit_log%info(message)
510 
511  do iprof = 1, nprofiles
512  profiles(iprof)%s2m%q = profiles(iprof)%q(nlevels)
513  enddo
514  endif
515 
516  varname = var_u ! Eastward-wind in m/s
517  if (ufo_vars_getindex(geovals%variables, varname) > 0) then
518  call ufo_geovals_get_var(geovals, varname, geoval)
519 
520  profiles(1:nprofiles)%s2m%u = geoval%vals(1,prof_start:prof_start + nprofiles - 1)
521  !assume if eastward then northward too
522 
523  varname = var_v ! Northward-wind in m/s
524  call ufo_geovals_get_var(geovals, varname, geoval)
525 
526  profiles(1:nprofiles)%s2m%v = geoval%vals(1,prof_start:prof_start + nprofiles - 1)
527  else !! use windspeed and direction instead
528  allocate(windsp(nprofiles))
529  call ufo_geovals_get_var(geovals, var_sfc_wspeed, geoval)
530 
531  windsp(1:nprofiles) = geoval%vals(1,prof_start:prof_start + nprofiles - 1)
532 
533  call ufo_geovals_get_var(geovals, var_sfc_wdir, geoval)
534 
535  do iprof = 1, nprofiles
536  profiles(iprof)%s2m%u = windsp(iprof) * cos(geoval%vals(1,prof_start + iprof - 1) * deg2rad)
537  profiles(iprof)%s2m%v = windsp(iprof) * sin(geoval%vals(1,prof_start + iprof - 1) * deg2rad)
538  enddo
539  deallocate(windsp)
540  endif
541 
542 ! Get surface type from geoval if it exists else diagnose surface type from water fraction.
543 ! If RTTOV surface type exists then it is expected that the related variables will also be available
544 ! e.g. water type, Tskin
545 ! DARFIX : This will be moved out to a separate subroutine in the next release to support all instruments
546 ! particularly where we have additional surface metadata for certain instruments
547 
548  varname = var_surf_type_rttov ! RTTOV surface type: 0 (land), 1 (water), 2 (sea-ice)
549  if (ufo_vars_getindex(geovals%variables, varname) > 0) then
550  call ufo_geovals_get_var(geovals, varname, geoval)
551  profiles(1:nprofiles)%skin%surftype = int(geoval%vals(1,prof_start:prof_start + nprofiles - 1),kind=jpim)
552 
553  varname = var_water_type_rttov ! RTTOV water type: 0 (fresh), 1 (sea)
554  call ufo_geovals_get_var(geovals, varname, geoval)
555  profiles(1:nprofiles)%skin%watertype = int(geoval%vals(1,prof_start:prof_start + nprofiles - 1),kind=jpim)
556 
557  varname = var_sfc_tskin !Skin (surface) temperature (K)
558  call ufo_geovals_get_var(geovals, varname, geoval)
559  profiles(1:nprofiles)%skin%t = geoval%vals(1,prof_start:prof_start + nprofiles - 1)
560 
561 
562  else
563 
564  ! Try to diagnose RTTOV surface from water fraction
565  !DARFIX: this is not consistent with CRTM in any way. Need to find out how they select surface type.
566 
567  varname = var_sfc_wfrac
568  if (ufo_vars_getindex(geovals%variables, varname) > 0) then
569  call ufo_geovals_get_var(geovals, varname, geoval)
570 
571  do iprof = 1, nprofiles
572  !Land point or sea point
573  wfrac = geoval%vals(1,iprof)
574  if (wfrac > half) then
575  profiles(iprof)%skin%surftype = surftype_sea
576  call ufo_geovals_get_var(geovals, var_sfc_wtmp, geoval)
577  profiles(iprof)%skin%t = geoval%vals(1,prof_start + iprof - 1)
578  else
579  !maybe it's predominantly land or ice
580  ! !determine land, snow and ice fractions and temperatures to determine average temperature
581  profiles(iprof)%skin%surftype = surftype_land ! land
582 
583  call ufo_geovals_get_var(geovals, var_sfc_lfrac, geoval)
584  lfrac = geoval%vals(1,prof_start + iprof - 1)
585 
586  call ufo_geovals_get_var(geovals, var_sfc_sfrac, geoval)
587  sfrac = geoval%vals(1,prof_start + iprof - 1)
588 
589  call ufo_geovals_get_var(geovals, var_sfc_ifrac, geoval)
590  ifrac = geoval%vals(1,prof_start + iprof - 1)
591 
592  call ufo_geovals_get_var(geovals, var_sfc_ltmp, geoval)
593  ltmp = geoval%vals(1,prof_start + iprof - 1)
594 
595  call ufo_geovals_get_var(geovals, var_sfc_stmp, geoval)
596  stmp = geoval%vals(1,prof_start + iprof - 1)
597 
598  call ufo_geovals_get_var(geovals, var_sfc_itmp, geoval)
599  itmp = geoval%vals(1,prof_start + iprof - 1)
600 
601  !Skin temperature is a combination of (i)ce temp, (l)and temp and (s)now temp
602  profiles(iprof)%skin%t = (lfrac * ltmp + sfrac * stmp + ifrac * itmp) / (lfrac + sfrac + ifrac)
603  endif
604  end do
605  endif
606  endif
607 
608  !DAR: Salinity fixed for now too
609  profiles(1:nprofiles)%skin%salinity = 35.0_kind_real
610 
611  !DAR: Default fastem parameters. We are not using FASTEM over land so these are unused
612  do iprof = 1,nprofiles
613 ! profiles(iprof)%skin%fastem = [3.0, 5.0, 15.0, 0.1, 0.3]
614  profiles(iprof)%skin%fastem = [0,0,0,0,0]
615  end do
616 
617 ! ---------------------------
618 ! SatRad profile manipulation
619 ! ---------------------------
620 
621  do iprof = 1, nprofiles
622  if(conf % SatRad_compatibility) then
623 
624  !----
625  ! Reset low level temperatures over seaice and cold, low land as per Ops_SatRad_SetUpRTprofBg.F90
626  ! N.B. I think this should be flagged so it's clear that the background has been modified
627  !----
628 
629  if(profiles(iprof)%skin%surftype /= surftype_sea) then
630  if(profiles(iprof)%skin%t < 271.4_kind_real .and. &
631  profiles(iprof)%s2m%p > 950.0_kind_real) then
632 
633  level_1000hpa = minloc(abs(profiles(iprof)%p - 1000.0_kind_real),dim=1)
634  level_950hpa = minloc(abs(profiles(iprof)%p - 950.0_kind_real),dim=1)
635 
636  newt = profiles(iprof)%t(level_950hpa)
637  if(profiles(iprof)%s2m%p > 1000.0_kind_real) &
638  newt = max(newt,profiles(iprof)%t(level_1000hpa))
639  newt = min(newt, 271.4_kind_real)
640 
641  profiles(iprof)%t(level_1000hpa) = max(profiles(iprof)%t(level_1000hpa), newt)
642  profiles(iprof)%s2m%t = max(profiles(iprof)%s2m%t, newt)
643  tstar = profiles(iprof)%skin%t
644  profiles(iprof)%skin%t = max(profiles(iprof)%skin%t, newt)
645  endif
646  endif
647 
648  !min_q fix
649  where(profiles(iprof)%q < min_q) profiles(iprof)%q = min_q
650  if(profiles(iprof)%s2m%q < min_q) profiles(iprof)%s2m%q = min_q
651 
652  endif
653 
654  enddo
655 
656 ! ----------------------------
657 ! Interpolate layer quantities
658 ! ----------------------------
659 ! Continuation of the profile interpolation if absorber profiles are layer quantities (RTTOV expects levels).
660 ! This is done once any profile manipulation (for q primarily) has been done.
661 
662  if(layer_quantities) then
663  do iprof=1,nprofiles
664  profiles(iprof)%q(top_level+stride:bottom_level:stride) = half * (profiles(iprof)%q(top_level:bottom_level-stride:stride) + &
665  profiles(iprof)%q(top_level+stride:bottom_level:stride))
666  profiles(iprof)%q(1) = profiles(iprof)%q(2)
667  profiles(iprof)%q(nlevels) = profiles(iprof)%q(nlevels-1) - half * (profiles(iprof)%q(nlevels-2) - profiles(iprof)%q(nlevels-1))
668 
669  if (associated(profiles(1)%o3)) then
670  profiles(iprof)%o3(top_level+stride:bottom_level:stride) = half * &
671  (profiles(iprof)%o3(top_level:bottom_level-stride:stride) + profiles(iprof)%o3(top_level+stride:bottom_level:stride))
672  profiles(iprof)%o3(1) = profiles(iprof)%o3(2)
673  profiles(iprof)%o3(nlevels) = profiles(iprof)%o3(nlevels-1) - &
674  half * (profiles(iprof)%o3(nlevels-2) - profiles(iprof)%o3(nlevels-1))
675  endif
676 
677  if (associated(profiles(1)%clw)) then
678  profiles(iprof)%clw(top_level+stride:bottom_level:stride) = half * &
679  (profiles(iprof)%clw(top_level:bottom_level-stride:stride) + profiles(iprof)%clw(top_level+stride:bottom_level:stride))
680  profiles(iprof)%clw(1) = profiles(iprof)%clw(2)
681  profiles(iprof)%clw(nlevels) = profiles(iprof)%clw(nlevels-1) - &
682  half * (profiles(iprof)%clw(nlevels-2) - profiles(iprof)%clw(nlevels-1))
683  endif
684  enddo
685  endif
686 
687 
688  !non mwscatt only at the moment
689  if(conf % SplitQtotal) then
690 
691  allocate(qtotal(nlevels), q_temp(nlevels), clw_temp(nlevels), ciw_temp(nlevels))
692 
693  do iprof = 1, nprofiles
694  ! compute bg qtotal using q and clw only
695  ! currently ice is ignored
696  qtotal(:) = profiles(iprof) % q(:) + profiles(iprof) % clw(:)
697 
698  ! generate first guess cloud and q based on the qtotal physics
699  call ops_satrad_qsplit (1, & ! in
700  profiles(iprof) % p(:), & ! in
701  profiles(iprof) % t(:), & ! in
702  qtotal(:), & ! in
703  q_temp(:), & ! out
704  clw_temp(:), & ! out
705  ciw_temp(:), & ! out
706  useqtsplitrain = conf % UseQtSplitRain)
707 
708  ! store in the profile
709  profiles(iprof) % clw(:) = clw_temp(:)
710  profiles(iprof) % q(:) = q_temp(:)
711 
712  ! !store non active variable
713  ! IF (ALLOCATED(CloudIce)) THEN
714 
715  ! CloudIce(toplevel_q_1dvar:) = ciw_temp(:)
716 
717  ! ENDIF
718  enddo
719 
720  deallocate(qtotal, q_temp, clw_temp, ciw_temp)
721 
722  end if
723 
724  if(present(obs_info)) then
725  ! profiles(1)%elevation = obs_info%elevation / 1000.0_kind_real ! m -> km
726  ! profiles(1)%latitude = obs_info%latitude
727  ! profiles(1)%longitude = obs_info%longitude
728 
729  else
730 
731  allocate(tmpvar(nlocs_total))
732  variable_present = obsspace_has(obss, "MetaData", "elevation")
733  if (variable_present) then
734  call obsspace_get_db(obss, "MetaData", "elevation", tmpvar)
735  profiles(1:nprofiles)%elevation = tmpvar(prof_start:prof_start + nprofiles - 1) * m_to_km !for RTTOV
736  else
737  variable_present = obsspace_has(obss, "MetaData", "surface_height")
738  if (variable_present) then
739  call obsspace_get_db(obss, "MetaData", "surface_height", tmpvar)
740  profiles(1:nprofiles)%elevation = tmpvar(prof_start:prof_start + nprofiles - 1) * m_to_km !for RTTOV
741  else ! from geoval
742  if (ufo_vars_getindex(geovals%variables, 'surface_altitude') > 0) then
743  call ufo_geovals_get_var(geovals, 'surface_altitude', geoval)
744  profiles(1:nprofiles)%elevation = geoval%vals(1,prof_start:prof_start + nprofiles - 1) * m_to_km
745  else
746  write(message,'(A)') &
747  'MetaData elevation not in database: check implicit filtering'
748  call fckit_log%info(message)
749  endif
750  endif
751  end if
752 
753  variable_present = obsspace_has(obss, "MetaData", "latitude")
754  if (variable_present) then
755  call obsspace_get_db(obss, "MetaData", "latitude", tmpvar )
756  profiles(1:nprofiles)%latitude = tmpvar(prof_start:prof_start + nprofiles - 1)
757  else
758  write(message,'(A)') &
759  'MetaData latitude not in database: check implicit filtering'
760  call fckit_log%info(message)
761  end if
762 
763  variable_present = obsspace_has(obss, "MetaData", "longitude")
764  if (variable_present) then
765  call obsspace_get_db(obss, "MetaData", "longitude", tmpvar )
766  profiles(1:nprofiles)%longitude = tmpvar(prof_start:prof_start + nprofiles - 1)
767  else
768  write(message,'(A)') &
769  'MetaData longitude not in database: check implicit filtering'
770  call fckit_log%info(message)
771  end if
772 
773  end if
774 
775  end subroutine load_atm_data_rttov
776 
777  ! Internal subprogam to load some test geometry data
778  subroutine load_geom_data_rttov(obss,profiles,prof_start1,obs_info)
779 
780  use obsspace_mod, only : obsspace_get_nlocs, obsspace_get_db
781 
782  implicit none
783 
784  type(c_ptr), value, intent(in) :: obss
785  type(rttov_profile), intent(inout) :: profiles(:)
786  integer, optional, intent(in) :: prof_start1
787  logical, optional, intent(in) :: obs_info ! DAR temporary
788 
789  real(kind_real), allocatable :: tmpvar(:)
790 
791  integer :: prof_start
792  integer :: nprofiles
793  integer :: nlevels
794 
795  if(present(prof_start1)) then
796  prof_start = prof_start1
797  else
798  prof_start = 1
799  end if
800 
801  if(present(obs_info)) then
802 
803  ! nprofiles = 1
804  ! nlevels = SIZE(profiles(1)%p)
805 
806  ! profiles(1)%zenangle = obs_info%sensor_zenith_angle
807  ! profiles(1)%azangle = obs_info%sensor_azimuth_angle
808  ! profiles(1)%sunzenangle = obs_info%solar_zenith_angle
809  ! profiles(1)%sunazangle = obs_info%solar_azimuth_angle
810 
811  else
812 
813  nprofiles = min(size(profiles), nlocs_total - prof_start + 1) !nlocs_total set in read_atm
814 
815  allocate(tmpvar(nlocs_total))
816 
817  nlevels = size(profiles(1)%p)
818 
819  !RTTOV convention 0-max (coef dependent). Nadir = 0 deg
820  call obsspace_get_db(obss, "MetaData", "sensor_zenith_angle", tmpvar)
821  if (any(tmpvar(prof_start:prof_start + nprofiles - 1) < 0.)) then
822  profiles(1:nprofiles)%zenangle = abs(tmpvar(prof_start:prof_start + nprofiles - 1))
823  ! DAR need to make a modification to azangle here then?
824  else
825  profiles(1:nprofiles)%zenangle = abs(tmpvar(prof_start:prof_start + nprofiles - 1))
826  endif
827 
828  !RTTOV convention is 0-360 deg. E=+90
829  call obsspace_get_db(obss, "MetaData", "sensor_azimuth_angle", tmpvar)
830  profiles(1:nprofiles)%azangle = tmpvar(prof_start:prof_start + nprofiles - 1)
831 
832  call obsspace_get_db(obss, "MetaData", "solar_zenith_angle", tmpvar)
833  if (any(tmpvar(prof_start:prof_start + nprofiles - 1) < 0.)) then
834  profiles(1:nprofiles)%sunzenangle = abs(tmpvar(prof_start:prof_start + nprofiles - 1))
835  ! DAR need to make a modification to sunazangle here then?
836  else
837  profiles(1:nprofiles)%sunzenangle = abs(tmpvar(prof_start:prof_start + nprofiles - 1))
838  endif
839 
840  call obsspace_get_db(obss, "MetaData", "solar_azimuth_angle", tmpvar)
841  profiles(1:nprofiles)%sunazangle = tmpvar(prof_start:prof_start + nprofiles - 1)
842 
843  deallocate(tmpvar)
844 
845  end if
846 
847  end subroutine load_geom_data_rttov
848 
849  ! ------------------------------------------------------------------------------
850 
851  subroutine set_options_rttov(self, f_confOpts)
852  implicit none
853 
854  class(rttov_conf), intent(inout) :: self
855  type(fckit_configuration), intent(in) :: f_confOpts ! RTcontrol
856 
857  include 'rttov_print_opts.interface'
858 
859  call self % set_defaults(self%RTTOV_default_opts) ! test, OPS, RTTOV
860 
861  !< Switch to enable atmospheric refraction
862  if (f_confopts % has("RTTOV_addrefrac")) then
863  call f_confopts % get_or_die("RTTOV_addrefrac", self % rttov_opts % rt_all % addrefrac)
864  end if
865 
866  !< Switch for input units in AD/K models
867  if (f_confopts % has("RTTOV_switchrad")) then
868  call f_confopts % get_or_die("RTTOV_switchrad", self % rttov_opts % rt_all % switchrad)
869  end if
870 
871  !< Switch to enable use of 2m q variable
872  if (f_confopts % has("RTTOV_use_q2m")) then
873  call f_confopts % get_or_die("RTTOV_use_q2m", self % rttov_opts % rt_all % use_q2m)
874  end if
875 
876  !< Switch for setting Lambertian reflection (IR and MW)
877  if (f_confopts % has("RTTOV_do_lambertian")) then
878  call f_confopts % get_or_die("RTTOV_do_lambertian", self % rttov_opts % rt_all % do_lambertian)
879  end if
880 
881  !< Switch for fixed/parameterised effective angle for Lambertian option
882  if (f_confopts % has("RTTOV_lambertian_fixed_angle")) then
883  call f_confopts % get_or_die("RTTOV_lambertian_fixed_angle", self % rttov_opts % rt_all % lambertian_fixed_angle)
884  end if
885 
886  !< Switch to ignore atmospheric curvature
887  if (f_confopts % has("RTTOV_plane_parallel")) then
888  call f_confopts % get_or_die("RTTOV_plane_parallel", self % rttov_opts % rt_all % plane_parallel)
889  end if
890 
891  !< Linear-in-tau or layer-mean for downwelling radiances
892  if (f_confopts % has("RTTOV_rad_down_lin_tau")) then
893  call f_confopts % get_or_die("RTTOV_rad_down_lin_tau", self % rttov_opts % rt_all % rad_down_lin_tau)
894  end if
895 
896  !< Switch to apply dtau test in transmit/integrate calculations
897  if (f_confopts % has("RTTOV_dtau_test")) then
898  call f_confopts % get_or_die("RTTOV_dtau_test", self % rttov_opts % rt_all % dtau_test)
899  end if
900 
901  !< FASTEM version (0-6); 0 => TESSEM2
902  if (f_confopts % has("RTTOV_fastem_version")) then
903  call f_confopts % get_or_die("RTTOV_fastem_version", self % rttov_opts % rt_mw % fastem_version)
904  end if
905 
906  !< Supply a foam fraction to FASTEM
907  if (f_confopts % has("RTTOV_supply_foam_fraction")) then
908  call f_confopts % get_or_die("RTTOV_supply_foam_fraction", self % rttov_opts % rt_mw % supply_foam_fraction)
909  end if
910 
911  !< Switch to enable input of cloud liquid water profile
912  if (f_confopts % has("RTTOV_clw_data")) then
913  call f_confopts % get_or_die("RTTOV_clw_data", self % rttov_opts % rt_mw % clw_data)
914  end if
915 
916  !< MW CLW scheme: 1 => Liebe, 2 => Rosenkranz, 3 => TKC
917  if (f_confopts % has("RTTOV_clw_scheme")) then
918  call f_confopts % get_or_die("RTTOV_clw_scheme", self % rttov_opts % rt_mw % clw_scheme)
919  end if
920 
921  !< Apply MW CLW calculations on coef/user levels (true/false resp.)
922  if (f_confopts % has("RTTOV_clw_calc_on_coef_lev")) then
923  call f_confopts % get_or_die("RTTOV_clw_calc_on_coef_lev", self % rttov_opts % rt_mw % clw_calc_on_coef_lev)
924  end if
925 
926  !< Lower pressure limit for MW CLW calculations (hPa)
927  if (f_confopts % has("RTTOV_clw_cloud_top")) then
928  call f_confopts % get_or_die("RTTOV_clw_cloud_top", self % rttov_opts % rt_mw % clw_cloud_top)
929  end if
930 
931  !< Apply band-correction for Planck radiance and BT calculations
932  if (f_confopts % has("RTTOV_apply_band_correction")) then
933  call f_confopts % get_or_die("RTTOV_apply_band_correction", self % rttov_opts % rt_mw % apply_band_correction)
934  end if
935 
936  !< Switch to enable RTTOV interpolator
937  if (f_confopts % has("RTTOV_addinterp")) then
938  call f_confopts % get_or_die("RTTOV_addinterp", self % rttov_opts % interpolation % addinterp)
939  end if
940 
941  !< Interpolation mode (1-5, see user guide)
942  if (f_confopts % has("RTTOV_interp_mode")) then
943  call f_confopts % get_or_die("RTTOV_interp_mode", self % rttov_opts % interpolation % interp_mode)
944  end if
945 
946  !< Switch to make pressure an active variable in TL/AD/K models
947  if (f_confopts % has("RTTOV_lgradp")) then
948  call f_confopts % get_or_die("RTTOV_lgradp", self % rttov_opts % interpolation % lgradp)
949  end if
950 
951  !< Switch to assume space boundary at top-most input pressure l
952  if (f_confopts % has("RTTOV_spacetop")) then
953  call f_confopts % get_or_die("RTTOV_spacetop", self % rttov_opts % interpolation % spacetop)
954  end if
955 
956  !< Switch to extrapolate input profiles using regression limits
957  if (f_confopts % has("RTTOV_reg_limit_extrap")) then
958  call f_confopts % get_or_die("RTTOV_reg_limit_extrap", self % rttov_opts % interpolation % reg_limit_extrap)
959  end if
960 
961  if (f_confopts % has("RTTOV_fix_hgpl")) then
962  call f_confopts % get_or_die("RTTOV_fix_hgpl", self % rttov_opts % config % fix_hgpl)
963  end if
964 
965  if (f_confopts % has("RTTOV_verbose")) then
966  call f_confopts % get_or_die("RTTOV_verbose", self % rttov_opts % config % verbose)
967  end if
968 
969  if (f_confopts % has("RTTOV_do_checkinput")) then
970  call f_confopts % get_or_die("RTTOV_do_checkinput", self % rttov_opts % config % do_checkinput)
971  end if
972 
973  if (f_confopts % has("RTTOV_apply_reg_limits")) then
974  call f_confopts % get_or_die("RTTOV_apply_reg_limits", self % rttov_opts % config % apply_reg_limits)
975  end if
976 
977  !< Solar sea BRDF model (1-2)
978  if (f_confopts % has("RTTOV_solar_sea_brdf_model")) then
979  call f_confopts % get_or_die("RTTOV_solar_sea_brdf_model", self % rttov_opts % rt_ir % solar_sea_brdf_model)
980  end if
981 
982  !< IR sea emissivity model (1-2)
983  if (f_confopts % has("RTTOV_ir_sea_emis_model")) then
984  call f_confopts % get_or_die("RTTOV_ir_sea_emis_model", self % rttov_opts % rt_ir % ir_sea_emis_model)
985  end if
986 
987  !< Switch to enable solar simulations
988  if (f_confopts % has("RTTOV_addsolar")) then
989  call f_confopts % get_or_die("RTTOV_addsolar", self % rttov_opts % rt_ir % addsolar)
990  end if
991 
992  !< Switch to enable Rayleigh single-scattering for VIS/NIR channel
993  if (f_confopts % has("RTTOV_rayleigh_single_scatt")) then
994  call f_confopts % get_or_die("RTTOV_rayleigh_single_scatt", self % rttov_opts % rt_ir % rayleigh_single_scatt)
995  end if
996 
997  !< Switch to enable NLTE bias correction
998  if (f_confopts % has("RTTOV_do_nlte_correction")) then
999  call f_confopts % get_or_die("RTTOV_do_nlte_correction", self % rttov_opts % rt_ir % do_nlte_correction)
1000  end if
1001 
1002  !< Switch to enable IR aerosol calculations
1003  if (f_confopts % has("RTTOV_addaerosl")) then
1004  call f_confopts % get_or_die("RTTOV_addaerosl", self % rttov_opts % rt_ir % addaerosl)
1005  end if
1006 
1007  !< Switch to supply aerosol optical properties explicitly per channel
1008  if (f_confopts % has("RTTOV_user_aer_opt_param")) then
1009  call f_confopts % get_or_die("RTTOV_user_aer_opt_param", self % rttov_opts % rt_ir % user_aer_opt_param)
1010  end if
1011 
1012  !< Switch to enable IR cloudy calculations
1013  if (f_confopts % has("RTTOV_addclouds")) then
1014  call f_confopts % get_or_die("RTTOV_addclouds", self % rttov_opts % rt_ir % addclouds)
1015  end if
1016 
1017  !< Switch to supply cloud optical properties explicitly per channel
1018  if (f_confopts % has("RTTOV_user_cld_opt_param")) then
1019  call f_confopts % get_or_die("RTTOV_", self % rttov_opts % rt_ir % user_cld_opt_param)
1020  end if
1021 
1022  !< Switch to supply grid-box average cloud concentration or cloud
1023  if (f_confopts % has("RTTOV_grid_box_avg_cloud")) then
1024  call f_confopts % get_or_die("RTTOV_grid_box_avg_cloud", self % rttov_opts % rt_ir % grid_box_avg_cloud)
1025  end if
1026 
1027  !! concentration in cloudy fraction of each layer
1028  !< Ignore cloud streams with weights lower than this
1029  if (f_confopts % has("RTTOV_cldstr_threshold")) then
1030  call f_confopts % get_or_die("RTTOV_cldstr_threshold", self % rttov_opts % rt_ir % cldstr_threshold)
1031  end if
1032 
1033  !< Switch for simplified cloud stream option - USE WITH CAUTION
1034  if (f_confopts % has("RTTOV_cldstr_simple")) then
1035  call f_confopts % get_or_die("RTTOV_cldstr_simple", self % rttov_opts % rt_ir % cldstr_simple)
1036  end if
1037 
1038  !< Upper pressure limit for cldstr_simple option (hPa)
1039  if (f_confopts % has("RTTOV_cldstr_low_cloud_top")) then
1040  call f_confopts % get_or_die("RTTOV_cldstr_low_cloud_top", self % rttov_opts % rt_ir % cldstr_low_cloud_top)
1041  end if
1042 
1043  !< IR scattering model to use
1044  if (f_confopts % has("RTTOV_ir_scatt_model")) then
1045  call f_confopts % get_or_die("RTTOV_ir_scatt_model", self % rttov_opts % rt_ir % ir_scatt_model)
1046  end if
1047 
1048  !< VIS/NIR scattering model to use
1049  if (f_confopts % has("RTTOV_vis_scatt_model")) then
1050  call f_confopts % get_or_die("RTTOV_vis_scatt_model", self % rttov_opts % rt_ir % vis_scatt_model)
1051  end if
1052 
1053  !< Number of DOM streams, must be even and not less than 2
1054  if (f_confopts % has("RTTOV_dom_nstreams")) then
1055  call f_confopts % get_or_die("RTTOV_dom_nstreams", self % rttov_opts % rt_ir % dom_nstreams)
1056  end if
1057 
1058  !< Convergence criterion for termination of DOM azimuthal loop
1059  if (f_confopts % has("RTTOV_dom_accuracy")) then
1060  call f_confopts % get_or_die("RTTOV_dom_accuracy", self % rttov_opts % rt_ir % dom_accuracy)
1061  end if
1062 
1063  !< DOM ignores levels below this optical depth:
1064  if (f_confopts % has("RTTOV_dom_opdep_threshold")) then
1065  call f_confopts % get_or_die("RTTOV_dom_opdep_threshold", self % rttov_opts % rt_ir % dom_opdep_threshold)
1066  end if
1067 
1068  !< Switch to enable input of O3 profile
1069  if (f_confopts % has("RTTOV_ozone_data")) then
1070  call f_confopts % get_or_die("RTTOV_ozone_data", self % rttov_opts % rt_ir % ozone_data)
1071  end if
1072  !< Switch to enable input of CO2 profile
1073  if (f_confopts % has("RTTOV_co2_data")) then
1074  call f_confopts % get_or_die("RTTOV_co2_data", self % rttov_opts % rt_ir % co2_data)
1075  end if
1076 
1077  !< Switch to enable input of N2O profile
1078  if (f_confopts % has("RTTOV_n2o_data")) then
1079  call f_confopts % get_or_die("RTTOV_n2o_data", self % rttov_opts % rt_ir % n2o_data)
1080  end if
1081 
1082  !< Switch to enable input of CO profile
1083  if (f_confopts % has("RTTOV_co_data")) then
1084  call f_confopts % get_or_die("RTTOV_co_data", self % rttov_opts % rt_ir % co_data)
1085  end if
1086 
1087  !< Switch to enable input of CH4 profile
1088  if (f_confopts % has("RTTOV_ch4_data")) then
1089  call f_confopts % get_or_die("RTTOV_ch4_data", self % rttov_opts % rt_ir % ch4_data)
1090  end if
1091 
1092  !< Switch to enable input of SO2 profile
1093  if (f_confopts % has("RTTOV_so2_data")) then
1094  call f_confopts % get_or_die("RTTOV_so2_data", self % rttov_opts % rt_ir % so2_data)
1095  end if
1096 
1097  ! OPEN(666,file='RTTOV_options.log')
1098  ! CALL rttov_print_opts(self % rttov_opts,lu=666)
1099  ! CLOSE(666)
1100  end subroutine set_options_rttov
1101 
1102  ! ------------------------------------------------------------------------------
1103 
1104  subroutine setup_rttov(self, f_confOpts, asw)
1105  class(rttov_conf), intent(inout) :: self
1106  type(fckit_configuration), intent(in) :: f_confOpts ! RTcontrol
1107  integer, intent(in) :: asw !allocate switch
1108 
1109  character(len=255) :: coef_filename
1110  character(len=4) :: coef_ext
1111  integer :: i_inst
1112 
1113  include 'rttov_read_coefs.interface'
1114 
1115  rttov_errorstatus = 0
1116  coef_ext = '.dat'
1117 
1118  if(asw == 1) then
1119 
1120  ! --------------------------------------------------------------------------
1121  ! 1. Setup rttov options
1122  ! --------------------------------------------------------------------------
1123  call self % set_options(f_confopts)
1124 
1125  ! --------------------------------------------------------------------------
1126  ! 2. Read coefficients
1127  ! --------------------------------------------------------------------------
1128  allocate(self % rttov_coef_array(self % nSensors))
1129 
1130  do i_inst = 1, self%nSensors
1131  coef_filename = &
1132  trim(self % COEFFICIENT_PATH) // 'rtcoef_' // trim(self%SENSOR_ID(i_inst)) // trim(coef_ext)
1133 
1134  call rttov_read_coefs(rttov_errorstatus, & !out
1135  self % rttov_coef_array(i_inst), & !inout
1136  self % rttov_opts, & !in
1137  file_coef = coef_filename) !in
1138 
1139  if (rttov_errorstatus /= errorstatus_success) then
1140  write(*,*) 'fatal error reading coefficients'
1141  else
1142  write(*,*) 'successfully read' // coef_filename
1143  end if
1144 
1145  end do
1146 
1147  self % rttov_is_setup =.true.
1148  else !asw == 0
1149  deallocate(self % rttov_coef_array)
1150  self%rttov_is_setup =.false.
1151  end if
1152 
1153  end subroutine setup_rttov
1154 
1155  ! ------------------------------------------------------------------------------
1156 
1157  subroutine get_var_name(n,varname)
1158 
1159  integer, intent(in) :: n
1160  character(len=*), intent(out) :: varname
1161 
1162  character(len=6) :: chan
1163 
1164  write(chan, '(I0)') n
1165  varname = 'brightness_temperature_' // trim(chan)
1166 
1167  end subroutine get_var_name
1168 
1169  !ufo_rttov_alloc is a wrapper for RTTOV12/13 allocation
1170  subroutine alloc_rttov(self, errorstatus, conf, nprofiles, nchannels, nlevels, init, asw)
1171  implicit none
1172 
1173  class(ufo_rttov_io), intent(inout) :: self
1174  type(rttov_conf), intent(in) :: conf
1175 
1176  integer, intent(out) :: errorstatus ! Return error status of RTTOV subroutine calls
1177  integer, intent(in) :: asw ! allocate switch
1178  logical, optional, intent(in) :: init ! initialise, default yes
1179  integer , intent(in) :: nchannels
1180  integer , intent(in) :: nprofiles
1181  integer , intent(in) :: nlevels
1182 
1183  logical :: init1
1184  integer :: asw_direct, asw_k
1185 
1186  integer, save :: called_from
1187 
1188  include 'rttov_alloc_direct.interface'
1189  include 'rttov_alloc_k.interface'
1190 
1191  asw_direct = -1
1192  asw_k = -1
1193 
1194  if(asw == 1) then
1195  asw_direct = 1
1196  called_from = 1
1197  endif
1198  if(asw == 2) then
1199  asw_k = 1
1200  called_from = 2
1201  endif
1202 
1203  if (present(init)) then
1204  init1 = init
1205  else
1206  init1 = .true.
1207  endif
1208 
1209  if(asw_direct /= -1 .and. called_from == 1) then
1210  ! Allocate structures for rttov_direct
1211  call rttov_alloc_direct( &
1212  errorstatus, &
1213  asw_direct, & ! 1 => allocate
1214  nprofiles, &
1215  nchannels, &
1216  nlevels, &
1217  self % chanprof, &
1218  conf % rttov_opts, &
1219  self % profiles, &
1220  conf % rttov_coef_array(1), &
1221  self % transmission, &
1222  self % radiance, &
1223  calcemis = self % calcemis, &
1224  emissivity = self % emissivity, &
1225  init=init1)
1226 
1227  if (errorstatus /= errorstatus_success) then
1228  write(message,'(A, I6)') 'after rttov_alloc_direct error = ', errorstatus
1229  call abor1_ftn(message)
1230  end if
1231 
1232  endif
1233 
1234  if (asw_k /= -1 .and. called_from == 2) then
1235  ! Allocate structures for rttov_k
1236  call rttov_alloc_k( &
1237  errorstatus, &
1238  asw_k, & ! 1 => allocate
1239  nprofiles, &
1240  nchannels, &
1241  nlevels, &
1242  self % chanprof, &
1243  conf % rttov_opts, &
1244  self % profiles, &
1245  self % profiles_k, &
1246  conf % rttov_coef_array(1), &
1247  self % transmission, &
1248  self % transmission_k, &
1249  self % radiance, &
1250  self % radiance_k, &
1251  calcemis = self % calcemis, &
1252  emissivity = self % emissivity, &
1253  emissivity_k = self % emissivity_k, &
1254  init=init1)
1255 
1256  if (errorstatus /= errorstatus_success) then
1257  write(message,'(A, I6)') 'after rttov_alloc_k error = ', errorstatus
1258  call abor1_ftn(message)
1259  end if
1260  endif
1261 
1262  if ( asw /= 0 ) then
1263  self % calcemis = .false.
1264  self % emissivity % emis_in = -1.0_kind_real
1265  self % emissivity % emis_out = -1.0_kind_real
1266  self % profiles(:) % skin % surftype = -1_jpim
1267  endif
1268 
1269  end subroutine alloc_rttov
1270 
1271  subroutine ufo_rttov_skip_profiles(nprofiles,nchannels,channels,obss,Skip_Profiles)
1272  ! Profiles are skipped when the ObsValue of all channels is missing.
1273  ! TODO: Use complete QC information
1274  ! It would be more comprehensive to use EffectiveQC or EffectiveError. That
1275  ! would require those ObsSpace values to be initialized before calls to
1276  ! this subroutine within ufo_radiancerttov_simobs+ufo_radiancerttov_tlad_settraj.
1277  implicit none
1278 
1279  integer, intent(in) :: nprofiles, nchannels
1280  type(c_ptr), value, intent(in) :: obss
1281  integer(c_int), intent(in) :: channels(:)
1282  logical, intent(inout) :: skip_profiles(:)
1283 
1284  integer :: jprofile, jchannel
1285  character(len=MAXVARLEN) :: varname
1286  real(kind_real) :: obsval(nprofiles,nchannels)
1287  !real(kind_real) :: EffObsErr(n_Profiles,n_Channels)
1288  !integer :: EffQC(n_Profiles,n_Channels)
1289 
1290  ! Set missing value
1291  missing = missing_value(missing)
1292 
1293  obsval = missing
1294  ! EffObsErr = missing
1295  ! EffQC = 0
1296 
1297  do jchannel = 1, nchannels
1298  call get_var_name(channels(jchannel),varname)
1299  call obsspace_get_db(obss, "ObsValue", varname, obsval(:,jchannel))
1300  ! call obsspace_get_db(obss, "EffectiveError", varname, EffObsErr(:,jchannel))
1301  ! call obsspace_get_db(obss, "EffectiveQC{iter}", varname, EffQC(:,jchannel))
1302  enddo
1303 
1304  !Loop over all n_Profiles, i.e. number of locations
1305  do jprofile = 1, nprofiles
1306  skip_profiles(jprofile) = all(obsval(jprofile,:) == missing) .or. all(obsval(jprofile,:) > 400)
1307  ! .OR. all(EffObsErr(jprofile,:) == missing) &
1308  ! .OR. all(EffQC(jprofile,:) /= 0)
1309  end do
1310 
1311  end subroutine ufo_rttov_skip_profiles
1312 
1313  subroutine ufo_rttov_init_emissivity(self, conf)
1314  class(ufo_rttov_io), intent(inout) :: self
1315  type(rttov_conf), intent(in) :: conf
1316 
1317  integer :: prof, ichan
1318 
1319  self % emissivity(:) % emis_in = -1.0_kind_real ! RMDI?
1320  self % emissivity(:) % emis_out = -1.0_kind_real ! RMDI?
1321  self % calcemis(:) = .false.
1322 
1323 !Emissivity and calcemis are only set for used channels. So if a profile is skipped then you must not set emis data for the channels that are skipped
1324 
1325  if ( conf % rttov_coef_array(1) % coef % id_sensor == sensor_id_mw) then
1326  do ichan = 1, nchan_sim, nchan_inst ! all channels initialised equally
1327  prof = self % chanprof(ichan)%prof
1328  if (self % profiles(prof) % skin % surftype == surftype_sea) then
1329 
1330  self % emissivity(ichan:ichan + nchan_inst - 1) % emis_in = 0.0_kind_real
1331  self % calcemis(ichan:ichan + nchan_inst - 1) = .true.
1332  else
1333  !IF ATLAS
1334 
1335  !ELSE
1336  if (self % profiles(prof) % skin % surftype == surftype_land) then
1337 
1338  !IF FASTEM (not implemented at MetO) ELSE
1339  self % emissivity(ichan:ichan + nchan_inst - 1) % emis_in = 0.95_kind_real
1340  elseif (self % profiles(prof) % skin % surftype == surftype_seaice) then
1341 
1342  !IF FASTEM (not implemented at MetO) ELSE
1343  self % emissivity(ichan:ichan + nchan_inst - 1) % emis_in = 0.92_kind_real
1344  endif
1345  !ENDIF !ATLAS
1346  endif
1347  enddo
1348  elseif ( conf % rttov_coef_array(1) % coef % id_sensor == sensor_id_ir .or. &
1349  conf % rttov_coef_array(1) % coef % id_sensor == sensor_id_hi) then
1350 
1351  do ichan = 1, size (self % chanprof(:)), nchan_inst ! all channels initialised equally
1352  prof = self % chanprof(ichan)%prof
1353  if (self % profiles(prof) % skin % surftype == surftype_sea) then
1354  ! Calculate by SSIREM or IREMIS
1355  self % emissivity(ichan:ichan + nchan_inst - 1) % emis_in = 0.0_kind_real
1356  self % calcemis(ichan:ichan + nchan_inst - 1) = .true.
1357  else
1358  !IF ATLAS ! CAMEL
1359  !ELSE
1360  if (self % profiles(prof) % skin % surftype == surftype_land) then
1361  self % emissivity(ichan:ichan + nchan_inst - 1) % emis_in = 0.95_kind_real
1362  elseif (self % profiles(prof) % skin % surftype == surftype_seaice) then
1363  !IF FASTEM (not implemented at MetO) ELSE
1364  self % emissivity(ichan:ichan + nchan_inst - 1) % emis_in = 0.92_kind_real
1365  endif
1366  !ENDIF !ATLAS
1367  endif
1368  enddo
1369  endif
1370 
1371  end subroutine ufo_rttov_init_emissivity
1372 
1373  subroutine set_defaults_rttov(self, default_opts_set)
1374 
1375  class(rttov_conf), intent(inout) :: self
1376  character(8), intent(in) :: default_opts_set
1377 
1378  write(message,'(A, A)') 'Setting RTTOV default options to ', default_opts_set
1379  call fckit_log%info(message)
1380 
1381  select case (trim(default_opts_set))
1382  case ('RTTOV') ! RTTOV defaults - needn't set this...
1383  self % rttov_opts % rt_all % addrefrac = .false. !< Switch to enable atmospheric refraction
1384  self % rttov_opts % rt_all % switchrad = .false. !< Switch for input units in AD/K models
1385  self % rttov_opts % rt_all % use_q2m = .true. !< Switch to enable use of 2m q variable
1386  self % rttov_opts % rt_all % do_lambertian = .false. !< Switch for setting Lambertian reflection (IR and MW)
1387  self % rttov_opts % rt_all % lambertian_fixed_angle = .true. !< Switch for fixed/parameterised effective angle for Lambertian option
1388  self % rttov_opts % rt_all % plane_parallel = .false. !< Switch to ignore atmospheric curvature
1389  self % rttov_opts % rt_all % rad_down_lin_tau = .true. !< Linear-in-tau or layer-mean for downwelling radiances
1390  self % rttov_opts % rt_all % dtau_test = .true. !< Switch to apply dtau test in transmit/integrate calculations
1391 
1392  self % rttov_opts % rt_ir % solar_sea_brdf_model = 1 !< Solar sea BRDF model (1-2)
1393  self % rttov_opts % rt_ir % ir_sea_emis_model = 2 !< IR sea emissivity model (1-2)
1394  self % rttov_opts % rt_ir % addsolar = .false. !< Switch to enable solar simulations
1395  self % rttov_opts % rt_ir % rayleigh_single_scatt = .true. !< Switch to enable Rayleigh single-scattering for VIS/NIR channels
1396  self % rttov_opts % rt_ir % do_nlte_correction = .false. !< Switch to enable NLTE bias correction
1397  self % rttov_opts % rt_ir % addaerosl = .false. !< Switch to enable IR aerosol calculations
1398  self % rttov_opts % rt_ir % user_aer_opt_param = .false. !< Switch to supply aerosol optical properties explicitly per channel
1399  self % rttov_opts % rt_ir % addclouds = .false. !< Switch to enable IR cloudy calculations
1400  self % rttov_opts % rt_ir % user_cld_opt_param = .false. !< Switch to supply cloud optical properties explicitly per channel
1401  self % rttov_opts % rt_ir % grid_box_avg_cloud = .false. !< Switch to supply grid-box average cloud concentration or cloud
1402  !! concentration in cloudy fraction of each layer
1403  self % rttov_opts % rt_ir % cldstr_threshold = -1.0_kind_real !< Ignore cloud streams with weights lower than this
1404  self % rttov_opts % rt_ir % cldstr_simple = .false. !< Switch for simplified cloud stream option - USE WITH CAUTION
1405  self % rttov_opts % rt_ir % cldstr_low_cloud_top = 750._kind_real !< Upper pressure limit for cldstr_simple option (hPa)
1406  self % rttov_opts % rt_ir % ir_scatt_model = ir_scatt_chou !< IR scattering model to use
1407  self % rttov_opts % rt_ir % vis_scatt_model = vis_scatt_dom !< VIS/NIR scattering model to use
1408  self % rttov_opts % rt_ir % dom_nstreams = 8 !< Number of DOM streams, must be even and not less than 2
1409  self % rttov_opts % rt_ir % dom_accuracy = 0._kind_real !< Convergence criterion for termination of DOM azimuthal loop
1410  self % rttov_opts % rt_ir % dom_opdep_threshold = 0._kind_real !< DOM ignores levels below this optical depth:
1411  !! 10. is reasonable, not applied if < = 0
1412  self % rttov_opts % rt_ir % ozone_data = .false. !< Switch to enable input of O3 profile
1413  self % rttov_opts % rt_ir % co2_data = .false. !< Switch to enable input of CO2 profile
1414  self % rttov_opts % rt_ir % n2o_data = .false. !< Switch to enable input of N2O profile
1415  self % rttov_opts % rt_ir % co_data = .false. !< Switch to enable input of CO profile
1416  self % rttov_opts % rt_ir % ch4_data = .false. !< Switch to enable input of CH4 profile
1417  self % rttov_opts % rt_ir % so2_data = .false. !< Switch to enable input of SO2 profile
1418 
1419  self % rttov_opts % rt_ir % pc % addpc = .false. !< Switch to enable PC-RTTOV
1420  self % rttov_opts % rt_ir % pc % ipcbnd = -1 !< PC spectral band
1421  self % rttov_opts % rt_ir % pc % ipcreg = -1 !< PC predictor channel set
1422  self % rttov_opts % rt_ir % pc % npcscores = -1 !< Number of PC scores to compute, if less than 1 npcscores is derived
1423  !! from the size of the pccomp%pcscores array
1424  self % rttov_opts % rt_ir % pc % addradrec = .false. !< Switch for calculation of reconstructed radiances
1425 
1426  !> MW-only radiative transfer options
1427  self % rttov_opts % rt_mw % fastem_version = 6 !< FASTEM version (0-6); 0 => TESSEM2
1428  self % rttov_opts % rt_mw % supply_foam_fraction = .false. !< Supply a foam fraction to FASTEM
1429  self % rttov_opts % rt_mw % clw_data = .false. !< Switch to enable input of cloud liquid water profile
1430  self % rttov_opts % rt_mw % clw_scheme = 1 !< MW CLW scheme: 1 => Liebe, 2 => Rosenkranz, 3 => TKC
1431  self % rttov_opts % rt_mw % clw_calc_on_coef_lev = .true. !< Apply MW CLW calculations on coef/user levels (true/false resp.)
1432  self % rttov_opts % rt_mw % clw_cloud_top = 322 !< Lower pressure limit for MW CLW calculations (hPa)
1433  self % rttov_opts % rt_mw % apply_band_correction = .true. !< Apply band-correction for Planck radiance and BT calculations
1434 
1435  self % rttov_opts % interpolation % addinterp = .false. !< Switch to enable RTTOV interpolator
1436  self % rttov_opts % interpolation % interp_mode = interp_rochon !< Interpolation mode (1-5, see user guide)
1437  self % rttov_opts % interpolation % lgradp = .false. !< Switch to make pressure an active variable in TL/AD/K models
1438  self % rttov_opts % interpolation % spacetop = .true. !< Switch to assume space boundary at top-most input pressure level
1439  self % rttov_opts % interpolation % reg_limit_extrap = .false. !< Switch to extrapolate input profiles using regression limits
1440 
1441  !> HTFRTC options structure
1442  self % rttov_opts % htfrtc_opts % htfrtc = .false. !< Switch to use htfrtc
1443  self % rttov_opts % htfrtc_opts % n_pc_in = -1 !< Number of principal components to be used
1444  self % rttov_opts % htfrtc_opts % reconstruct = .false. !< Switch to select reconstructed radiances
1445  self % rttov_opts % htfrtc_opts % simple_cloud = .false. !< Calculate simple cloud
1446  self % rttov_opts % htfrtc_opts % overcast = .false. !< Calculate overcast cloud on all levels
1447 
1448  case ('OPS','PS45')
1449  self % rttov_opts % config % apply_reg_limits = .true. ! set as true for all instruments
1450  self % rttov_opts % config % fix_hgpl = .false. ! true @ PS44
1451  self % rttov_opts % config % verbose = .false. ! true if (ProcessMode > VerboseMode .OR. RTTOV_Verbosity > 0)
1452  self % rttov_opts % config % do_checkinput = .false. ! we will use the more thorough and verbose
1453  ! rttov_user_profile_checkinput instead to save ourselves an
1454  ! RTTOV call where the profile is no good
1455 
1456  self % rttov_opts % rt_all % rad_down_lin_tau = .true. !FALSE @PS44
1457  self % rttov_opts % rt_all % dtau_test = .true. !FALSE @PS44
1458 
1459  if(trim(default_opts_set) == 'OPS') then
1460  self % rttov_opts % rt_all % addrefrac = .false. ! At PS44 .TRUE. @ PS45?
1461  else
1462  self % rttov_opts % rt_all % addrefrac = .true. ! At PS44 .TRUE. @ PS45?
1463  endif
1464 
1465  self % rttov_opts % rt_all % switchrad = .true.
1466  if(trim(default_opts_set) == 'OPS') then
1467  self % rttov_opts % rt_all % use_q2m = .false. ! At PS44 .TRUE. @ PS45?
1468  endif
1469  self % rttov_opts % rt_all % do_lambertian = .false. !
1470  ! self % rttov_opts % rt_all % plane_parallel = .FALSE.
1471 
1472  self % rttov_opts % interpolation % addinterp = .true. ! Allow interpolation of input profile
1473  self % rttov_opts % interpolation % interp_mode = 4 ! Set interpolation method (4 for all insts at PS44)
1474  ! self % rttov_opts % interpolation % lgradp = .FALSE.
1475  if(trim(default_opts_set) == 'OPS') then
1476  self % rttov_opts % interpolation % spacetop = .false. ! At PS44 .TRUE. @ PS45?
1477  self % rttov_opts % interpolation % reg_limit_extrap = .false. ! At PS44 .TRUE. @ PS45?
1478  else
1479  self % rttov_opts % interpolation % reg_limit_extrap = .true. ! At PS44 .TRUE. @ PS45?
1480  endif
1481 
1482  ! self % rttov_opts % rt_ir % addsolar = .FALSE.
1483  ! self % rttov_opts % rt_ir % do_nlte_correction = .FALSE.
1484  self % rttov_opts % rt_ir % ir_sea_emis_model = 1 ! At PS44 2 @ PS45?
1485 
1486  ! self % rttov_opts % rt_ir % user_aer_opt_param = .FALSE. ! user specifies aerosol scattering optical parameters
1487  ! self % rttov_opts % rt_ir % user_cld_opt_param = .FALSE. ! user specifies cloud scattering optical parameters
1488  ! self % rttov_opts % rt_ir % ir_scatt_model = 2 ! Chou-scaling
1489  self % rttov_opts % rt_ir % grid_box_avg_cloud = .true. ! Assume grid-box average for cloud concentrations
1490 
1491  self % rttov_opts % rt_ir % ozone_data = .true. ! There WILL be an ozone profile
1492  self % rttov_opts % rt_ir % co2_data = .false. ! We do not have profiles
1493  self % rttov_opts % rt_ir % n2o_data = .false. ! for any other constituents
1494  self % rttov_opts % rt_ir % ch4_data = .false.
1495  self % rttov_opts % rt_ir % co_data = .false.
1496  self % rttov_opts % rt_ir % so2_data = .false.
1497 
1498  ! RTTOV v12 Discrete ordinates method (DOM) not implemented for now
1499 
1500  self % rttov_opts % rt_mw % clw_calc_on_coef_lev = .true. ! .FALSE. @PS44
1501  self % rttov_opts % rt_mw % clw_data = .true. ! Set to true for allocation
1502  if(trim(default_opts_set) == 'OPS') then
1503  self % rttov_opts % rt_mw % fastem_version = 2 ! at PS44 6 at PS45?
1504  endif
1505  !self % rttov_opts % rt_mw % supply_foam_fraction = RTTOV_Supply_Foam_Fraction
1506 
1507  ! IF (ANY (MwscattSwitch)) THEN
1508  ! RTTOV12_opts_scatt % config = self % rttov_opts % config
1509  ! RTTOV12_opts_scatt % use_q2m = self % rttov_opts % rt_all % use_q2m
1510  ! RTTOV12_opts_scatt % fastem_version = self % rttov_opts % rt_mw % fastem_version
1511  ! RTTOV12_opts_scatt % supply_foam_fraction = self % rttov_opts % rt_mw % supply_foam_fraction
1512  ! RTTOV12_opts_scatt % interp_mode = self % rttov_opts % interpolation % interp_mode
1513  ! RTTOV12_opts_scatt % reg_limit_extrap = self % rttov_opts % interpolation % reg_limit_extrap
1514  ! RTTOV12_opts_scatt % rttov9_fastem_bug = self % rttov_opts % compatibility % rttov9_fastem_bug
1515  ! END IF
1516 
1517  case ('test')
1518  self % rttov_opts % rt_ir % addsolar = .false. ! Do not include solar radiation
1519  self % rttov_opts % interpolation % addinterp = .true. ! Allow interpolation of input profile
1520  self % rttov_opts % interpolation % interp_mode = 4 ! Set interpolation method
1521  self % rttov_opts % interpolation % reg_limit_extrap = .true. ! reg_limit_extrap
1522  self % rttov_opts % rt_all % addrefrac = .true. ! Include refraction in path calc
1523  self % rttov_opts % rt_all % switchrad = .true. ! Include refraction in path calc
1524  self % rttov_opts % rt_all % dtau_test = .false.
1525  self % rttov_opts % rt_ir % addclouds = .false. ! Don't include cloud effects
1526  self % rttov_opts % rt_ir % addaerosl = .false. ! Don't include aerosol effects
1527 
1528  self % rttov_opts % rt_ir % ozone_data = .false. ! Set the relevant flag to .TRUE.
1529  self % rttov_opts % rt_ir % co2_data = .false. ! when supplying a profile of the
1530  self % rttov_opts % rt_ir % n2o_data = .false. ! given trace gas (ensure the
1531  self % rttov_opts % rt_ir % ch4_data = .false. ! coef file supports the gas)
1532  self % rttov_opts % rt_ir % co_data = .false. !
1533  self % rttov_opts % rt_ir % so2_data = .false. !
1534  self % rttov_opts % rt_mw % clw_data = .false. !
1535 
1536  self % rttov_opts % config % verbose = .true. ! Enable printing of warnings
1537  self % rttov_opts % config % apply_reg_limits = .true.
1538  self % rttov_opts % config % do_checkinput = .false.
1539  self % rttov_opts % config % fix_hgpl = .false.
1540  case('default')
1541  ! RTTOV defaults
1542  end select
1543 
1544  end subroutine set_defaults_rttov
1545 
1546  subroutine populate_hofxdiags(RTProf, chanprof, hofxdiags)
1547  type(ufo_rttov_io), intent(in) :: rtprof
1548  type(rttov_chanprof), intent(in) :: chanprof(:)
1549  type(ufo_geovals), intent(inout) :: hofxdiags !non-h(x) diagnostics
1550 
1551  integer :: jvar, chan, prof, ichan
1552  integer :: nlayers, nchanprof, nlevels, nprofiles
1553  real(kind_real), allocatable :: od_level(:), wfunc(:)
1554 
1555  include 'rttov_calc_weighting_fn.interface'
1556 
1557  allocate(od_level(size(rtprof % transmission%tau_levels(:,1))))
1558  allocate(wfunc(size(rtprof % transmission%tau_levels(:,1))))
1559 
1560  missing = missing_value(missing)
1561 
1562  nchanprof = size(chanprof)
1563  nlevels = size(rtprof % profiles(1) % p)
1564  nprofiles = maxval(chanprof(:)%prof) !SIZE(RTProf % profiles)
1565 
1566  do jvar = 1, hofxdiags%nvar
1567  if (len(trim(hofxdiags%variables(jvar))) < 1) cycle
1568 
1569  !============================================
1570  ! Diagnostics used for QC and bias correction
1571  !============================================
1572 
1573  if (trim(xstr_diags(jvar)) == "") then
1574  ! forward h(x) diags
1575  select case(trim(ystr_diags(jvar)))
1576 
1577  ! variable: optical_thickness_of_atmosphere_layer_CH
1578  ! variable: transmittances_of_atmosphere_layer_CH
1579  ! variable: weightingfunction_of_atmosphere_layer_CH
1581 
1582  nlayers = nlevels - 1
1583  hofxdiags%geovals(jvar)%nval = nlayers
1584  allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,nprofiles))
1585  hofxdiags%geovals(jvar)%vals = missing
1586  ! get channel/profile
1587  do ichan = 1, nchanprof
1588  chan = chanprof(ichan)%chan
1589  prof = chanprof(ichan)%prof
1590 
1591  if(chan == ch_diags(jvar)) then
1592  ! if profile not skipped
1593  if(ystr_diags(jvar) == var_opt_depth) then
1594  od_level(:) = log(rtprof % transmission%tau_levels(:,chan)) !level->TOA transmittances -> od
1595  hofxdiags%geovals(jvar)%vals(:,prof) = od_level(1:nlevels-1) - od_level(2:nlevels) ! defined +ve
1596  else if (ystr_diags(jvar) == var_lvl_transmit) then
1597  hofxdiags%geovals(jvar)%vals(:,prof) = rtprof % transmission % tau_levels(1:nlevels-1,chan) - rtprof % transmission%tau_levels(2:,chan)
1598  else if (ystr_diags(jvar) == var_lvl_weightfunc) then
1599  od_level(:) = log(rtprof % transmission%tau_levels(:,chan)) !level->TOA transmittances -> od
1600  call rttov_calc_weighting_fn(rttov_errorstatus, rtprof % profiles(prof)%p, od_level(:), &
1601  hofxdiags%geovals(jvar)%vals(:,prof))
1602 
1603  endif
1604  !endif
1605  endif
1606  enddo
1607 
1608  ! variable: toa_outgoing_radiance_per_unit_wavenumber_CH [mW / (m^2 sr cm^-1)] (nval=1)
1609  ! variable: brightness_temperature_assuming_clear_sky_CH
1610  ! variable: pressure_level_at_peak_of_weightingfunction_CH
1612  ! always returned
1613  hofxdiags%geovals(jvar)%nval = 1
1614  allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,nprofiles))
1615  hofxdiags%geovals(jvar)%vals = missing
1616 
1617  do ichan = 1, nchanprof
1618  chan = chanprof(ichan)%chan
1619  prof = chanprof(ichan)%prof
1620 
1621  if(chan == ch_diags(jvar)) then
1622  if(ystr_diags(jvar) == var_radiance) then
1623  hofxdiags%geovals(jvar)%vals(1,prof) = rtprof % radiance % total(ichan)
1624  else if(ystr_diags(jvar) == var_tb_clr) then
1625  hofxdiags%geovals(jvar)%vals(1,prof) = rtprof % radiance % bt_clear(ichan)
1626  else if(ystr_diags(jvar) == var_tb) then
1627  hofxdiags%geovals(jvar)%vals(1,prof) = rtprof % radiance % bt(ichan)
1628  else if(ystr_diags(jvar) == var_pmaxlev_weightfunc ) then
1629  call rttov_calc_weighting_fn(rttov_errorstatus, rtprof % profiles(prof)%p, od_level(:), &
1630  wfunc(:))
1631  hofxdiags%geovals(jvar)%vals(1,prof) = maxloc(wfunc(:), dim=1) ! scalar not array(1)
1632  end if
1633  endif
1634  end do
1635 
1636  case default
1637  write(message,*) 'ufo_radiancerttov_simobs: //&
1638  & ObsDiagnostic is unsupported, ', &
1639  & hofxdiags%variables(jvar)
1640  call fckit_log%info(message)
1641  end select
1642 
1643  else if (trim(ystr_diags(jvar)) == var_tb) then
1644  ! var_tb jacobians
1645  select case (trim(xstr_diags(jvar)))
1646 
1647  case (var_ts,var_mixr,var_clw)
1648  nlayers = nlevels - 1
1649  hofxdiags%geovals(jvar)%nval = nlayers
1650  allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,nprofiles))
1651  hofxdiags%geovals(jvar)%vals = missing
1652 
1653  do ichan = 1, nchanprof
1654  chan = chanprof(ichan)%chan
1655  prof = chanprof(ichan)%prof
1656 
1657  if(chan == ch_diags(jvar)) then
1658  if(ystr_diags(jvar) == var_ts) then
1659  hofxdiags%geovals(jvar)%vals(:,prof) = &
1660  rtprof % profiles_k(ichan) % t(:)
1661  else if(ystr_diags(jvar) == var_mixr) then
1662  hofxdiags%geovals(jvar)%vals(:,prof) = &
1663  rtprof % profiles_k(ichan) % q(:)
1664  else if(ystr_diags(jvar) == var_clw) then
1665  hofxdiags%geovals(jvar)%vals(:,prof) = &
1666  rtprof % profiles_k(ichan) % clw(:)
1667  endif
1668  endif
1669  enddo
1670 
1672  hofxdiags%geovals(jvar)%nval = 1
1673  allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,nprofiles))
1674  hofxdiags%geovals(jvar)%vals = missing
1675 
1676  do ichan = 1, nchanprof
1677  chan = chanprof(ichan)%chan
1678  prof = chanprof(ichan)%prof
1679 
1680  if(chan == ch_diags(jvar)) then
1681  if(ystr_diags(jvar) == var_sfc_tskin) then
1682  hofxdiags%geovals(jvar)%vals(1,prof) = &
1683  rtprof % profiles_k(ichan) % skin % t
1684  else if (ystr_diags(jvar) == var_sfc_t2m) then
1685  hofxdiags%geovals(jvar)%vals(1,prof) = &
1686  rtprof % profiles_k(ichan) % s2m % t
1687  else if (ystr_diags(jvar) == var_sfc_p2m) then
1688  hofxdiags%geovals(jvar)%vals(1,prof) = &
1689  rtprof % profiles_k(ichan) % s2m % p
1690  else if (ystr_diags(jvar) == var_sfc_q2m) then
1691  hofxdiags%geovals(jvar)%vals(1,prof) = &
1692  rtprof % profiles_k(ichan) % s2m % q
1693  else if (ystr_diags(jvar) == var_u) then
1694  hofxdiags%geovals(jvar)%vals(1,prof) = &
1695  rtprof % profiles_k(ichan) % s2m % u
1696  else if (ystr_diags(jvar) == var_v) then
1697  hofxdiags%geovals(jvar)%vals(1,prof) = &
1698  rtprof % profiles_k(ichan) % s2m % v
1699  else if (ystr_diags(jvar) == var_sfc_emiss) then
1700  hofxdiags%geovals(jvar)%vals(1,prof) = &
1701  rtprof % emissivity_k(ichan) % emis_in
1702  end if
1703  end if
1704  end do
1705 
1706  case default
1707  write(message,*) 'ufo_radiancerttov_simobs: //&
1708  & ObsDiagnostic is unsupported, ', &
1709  & hofxdiags%variables(jvar)
1710  call fckit_log%info(message)
1711  end select
1712  else
1713  write(message,*) 'ufo_radiancerttov_simobs: //&
1714  & ObsDiagnostic is unsupported, ', &
1715  & hofxdiags%variables(jvar)
1716  call fckit_log%info(message)
1717  end if
1718 
1719  enddo
1720 
1721  deallocate(od_level,wfunc)
1722 
1723  end subroutine populate_hofxdiags
1724 
1725  subroutine parse_hofxdiags(hofxdiags, jacobian_needed)
1726  implicit none
1727 
1728  type(ufo_geovals), intent(inout) :: hofxdiags !non-h(x) diagnostics
1729  logical, intent(out) :: jacobian_needed
1730  character(10), parameter :: jacobianstr = "_jacobian_"
1731 
1732  integer :: str_pos(4)
1733  character(len=MAXVARLEN) :: varstr
1734  integer :: jvar
1735  character(len=max_string) :: err_msg
1736 
1737  jacobian_needed = .false.
1738 
1739  if(hofxdiags%nvar > 0) then
1740  allocate (ystr_diags(hofxdiags%nvar), &
1741  xstr_diags(hofxdiags%nvar), &
1742  ch_diags(hofxdiags%nvar))
1743 
1744  ch_diags = -9999
1745 
1746  do jvar = 1, hofxdiags%nvar
1747  varstr = hofxdiags%variables(jvar)
1748 
1749  str_pos(4) = len_trim(varstr)
1750  if (str_pos(4) < 1) cycle
1751  str_pos(3) = index(varstr,"_",back=.true.) !final "_" before channel
1752  read(varstr(str_pos(3)+1:str_pos(4)),*, err=999) ch_diags(jvar)
1753 999 str_pos(1) = index(varstr,jacobianstr) - 1 !position before jacobianstr
1754  if (str_pos(1) == 0) then
1755  write(err_msg,*) 'parse_hofxdiags: _jacobian_ must be // &
1756  & preceded by dependent variable in config: ', &
1757  & hofxdiags%variables(jvar)
1758  call abor1_ftn(err_msg)
1759  else if (str_pos(1) > 0) then
1760  !Diagnostic is a Jacobian member (dy/dx)
1761  ystr_diags(jvar) = varstr(1:str_pos(1))
1762  str_pos(2) = str_pos(1) + len(jacobianstr) + 1 !begin xstr_diags
1763  jacobian_needed = .true.
1764  str_pos(4) = str_pos(3) - str_pos(2)
1765  xstr_diags(jvar)(1:str_pos(4)) = varstr(str_pos(2):str_pos(3)-1)
1766  xstr_diags(jvar)(str_pos(4)+1:) = ""
1767  else !null
1768  !Diagnostic is a dependent variable (y)
1769 
1770  xstr_diags(jvar) = ""
1771  ystr_diags(jvar)(1:str_pos(3)-1) = varstr(1:str_pos(3)-1)
1772  ystr_diags(jvar)(str_pos(3):) = ""
1773  if (ch_diags(jvar) < 0) ystr_diags(jvar) = varstr
1774  end if
1775  end do
1776  end if
1777 
1778  end subroutine parse_hofxdiags
1779 
1780  !-------------------------------------------------------------------------------
1781  ! (C) Crown copyright Met Office. All rights reserved.
1782  ! Refer to COPYRIGHT.txt of this distribution for details.
1783  !-------------------------------------------------------------------------------
1784  ! IF output_type=1 : Split total water content (qtotal) into
1785  ! water vapor content (q) and
1786  ! cloud liquid water content (ql) and
1787  ! cloud ice water content (qi)
1788  ! If output_type ne 1 : Compute derivatives: (q) =dq/dqtotal
1789  ! (ql)=dql/dqtotal
1790  ! (qi)=dqi/dqtotal
1791  !
1792  ! WARNING: The derivatives are not valid if LtemperatureVar=.TRUE. since
1793  ! qsaturated depends on temperature.
1794  !
1795  ! The partitioning of the excess moisture between ice and clw uses a temperature
1796  ! based parametrization based on aircraft data Ref: Jones DC Reading phdthesis
1797  ! p126
1798  !-------------------------------------------------------------------------------
1799 
1800  subroutine ops_satrad_qsplit ( &
1801  output_type, &
1802  p, &
1803  t, &
1804  qtotal, &
1805  q, &
1806  ql, &
1807  qi, &
1808  UseQtSplitRain)
1809 
1810  use ufo_utils_mod, only : ops_qsat
1811  use ufo_constants_mod, only : zero, zerodegc, min_q
1812 
1813  implicit none
1814 
1815  ! Subroutine arguments:
1816  integer, intent(in) :: output_type
1817  real(kind_real), intent(in) :: p(:)
1818  real(kind_real), intent(in) :: t(:)
1819  real(kind_real), intent(in) :: qtotal(:)
1820  real(kind_real), intent(out) :: q(size(qtotal)) ! humidity component q
1821  real(kind_real), intent(out) :: ql(size(qtotal)) ! liquid component ql
1822  real(kind_real), intent(out) :: qi(size(qtotal)) ! ice component qi
1823  logical, intent(in) :: UseQtSplitRain
1824 
1825  ! Local declarations:
1826  integer :: i
1827  real(kind_real), parameter :: lower_rh = 0.95
1828  real(kind_real), parameter :: upper_rh = 1.05
1829  real(kind_real), parameter :: Split_Factor = 0.5
1830  real(kind_real), parameter :: MinTempQl = 233.15 ! temperature (K) below which all cloud is ice
1831  real(kind_real) :: qsaturated(size(qtotal))
1832  real(kind_real) :: RH_qtotal(size(qtotal))
1833  real(kind_real) :: qnv(size(qtotal)) ! non vapour component
1834  real(kind_real) :: qc(size(qtotal)) ! cloud component
1835  real(kind_real) :: V1(size(qtotal))
1836  real(kind_real) :: V2(size(qtotal))
1837  real(kind_real) :: V1zero
1838  real(kind_real) :: V2zero
1839  real(kind_real) :: W(size(qtotal))
1840  real(kind_real) :: Y1,Y2,Y3,Y4
1841  real(kind_real) :: IntConst
1842  real(kind_real) :: Aconst
1843  real(kind_real) :: Bconst
1844  real(kind_real) :: Cconst
1845  real(kind_real) :: Dconst
1846  real(kind_real) :: Denom
1847  real(kind_real) :: SmallValue
1848  real(kind_real) :: LF(size(qtotal)) ! fraction of ql to ql+qi
1849  character(len=*), parameter :: RoutineName = "Ops_SatRad_Qsplit"
1850 
1851  real(kind_real) :: QsplitRainParamA !Parameters used to define proportion of
1852  real(kind_real) :: QsplitRainParamB !qt that is partitioned into a rain compnenent
1853  real(kind_real) :: QsplitRainParamC !
1854 
1855  !Qsplit_MixPhaseParam = 1 !Jones' method as default
1856  qsplitrainparama = 0.15_kind_real
1857  qsplitrainparamb = 0.09_kind_real
1858  qsplitrainparamc = 50.0_kind_real
1859 
1860  ! Compute saturated water vapor profile for nlevels_q only
1861 
1862  call ops_qsat (qsaturated(:), & ! out
1863  t(:), & ! in
1864  p(:), & ! in
1865  size(qtotal)) ! in
1866 
1867  smallvalue = 1.0_kind_real / 8.5_kind_real
1868  denom = smallvalue * (upper_rh - lower_rh)
1869 
1870  ! don't let rh exceed 2.0 to avoid cosh function blowing up
1871  rh_qtotal(:) = min(qtotal(:) / qsaturated(:), 2.0_kind_real)
1872 
1873  v1(:) = (rh_qtotal(:) - lower_rh) / denom
1874  v2(:) = (rh_qtotal(:) - upper_rh) / denom
1875 
1876  y1 = 1.0_kind_real
1877  y2 = split_factor
1878  y3 = 0.0_kind_real
1879  y4 = split_factor
1880 
1881  aconst = (y2 - y1) / 2.0_kind_real
1882  bconst = -(y4 - y3) / 2.0_kind_real
1883  cconst = (y2 + y1) / 2.0_kind_real
1884  dconst = -(y4 + y3) / 2.0_kind_real
1885 
1886  ! Compute fraction of ql to ql+qi based on temperature profile
1887 
1888  where (t(:) - zerodegc >= -0.01_kind_real) ! -0.01degc and above
1889 
1890  ! all ql
1891  lf(:) = 1.0_kind_real
1892 
1893  end where
1894 
1895  where (t(:) <= mintempql)
1896 
1897  ! all qi
1898  lf(:) = 0.0_kind_real
1899 
1900  end where
1901 
1902  where (t(:) > mintempql .and. t(:) - zerodegc < -0.01_kind_real)
1903 
1904  ! Jones' parametrization
1905  lf(:) = sqrt(-1.0_kind_real * log(-0.025_kind_real * (t(:) - zerodegc)) / 70.0_kind_real)
1906 
1907  end where
1908 
1909  ! finally set LF to 0.0_kind_real for the rttov levels on which clw jacobians are not
1910  ! calculated since nlevels_mwclw < nlevels_q
1911 
1912  v1zero = -1.0_kind_real * lower_rh / denom
1913  v2zero = -1.0_kind_real * upper_rh / denom
1914  intconst = -(aconst * denom * log(cosh(v1zero)) + bconst * denom * log(cosh(v2zero)))
1915  w(:) = aconst * denom * log(cosh(v1(:))) + bconst * denom * log(cosh(v2(:))) + &
1916  (cconst + dconst) * rh_qtotal(:) + intconst
1917 
1918  ! store the components of qtotal
1919  ! ensuring that they are above lower limits
1920 
1921  if (useqtsplitrain) then
1922 
1923  ! Split qtotal into q and qnv (non-vapour part - includes
1924  ! ql, qi, qr)
1925 
1926  ! Split qtotal into q, ql ,qi
1927 
1928  do i = 1, size(qtotal)
1929 
1930  q(i) = max(w(i) * qsaturated(i), min_q)
1931  qnv(i) = max(qtotal(i) - q(i), 0.0_kind_real)
1932 
1933  ! Split qnv into a cloud and precipitation part
1934 
1935  qc(i) = max(qsplitrainparama * (qsplitrainparamb - (qsplitrainparamb / ((qsplitrainparamc * qnv(i)) + 1))), 0.0_kind_real)
1936 
1937  ! Finally split non-precip part into liquid and ice
1938 
1939  ql(i) = max(lf(i) * qc(i), 0.0_kind_real)
1940  qi(i) = max((1.0_kind_real - lf(i)) * (qc(i)), 0.0_kind_real)
1941 
1942  end do
1943 
1944  else
1945  do i = 1, size(qtotal)
1946 
1947  q(i) = max(w(i) * qsaturated(i), min_q)
1948  ql(i) = max(lf(i) * (qtotal(i) - q(i)), 0.0_kind_real)
1949  qi(i) = max((1.0_kind_real - lf(i)) * (qtotal(i) - q(i)), 0.0_kind_real)
1950 
1951  end do
1952 
1953  end if
1954 
1955  ! Values of q, ql and qi are overwritten if output_type /= 1
1956  ! and replaced with the derivatives
1957 
1958  if (output_type /= 1) then
1959 
1960  ! Compute derivates
1961  ! q = dq/dqtotal, ql = dql/dqtotal, qi=dqi/dqtotal
1962 
1963  q(:) = aconst * tanh(v1(:)) + cconst + bconst * tanh(v2(:)) + dconst
1964 
1965  if (useqtsplitrain) then
1966 
1967  ql(:) = lf(:) * qsplitrainparama * qsplitrainparamb * qsplitrainparamc * (1.0_kind_real - q(:)) / &
1968  ((qsplitrainparamc * qnv(:)) + 1.0_kind_real) ** 2
1969  qi(:) = (1.0_kind_real - lf(:)) * qsplitrainparama * qsplitrainparamb * qsplitrainparamc * (1.0_kind_real - q(:)) / &
1970  ((qsplitrainparamc * qnv(:)) + 1.0_kind_real) ** 2
1971 
1972  else
1973 
1974  ql(:) = lf(:) * (1.0_kind_real - q(:))
1975  qi(:) = (1.0_kind_real - lf(:)) * (1.0_kind_real - q(:))
1976 
1977  end if
1978 
1979  end if
1980 
1981  end subroutine ops_satrad_qsplit
1982 
1983 
1984 end module ufo_radiancerttov_utils_mod
ufo_avgkernel_mod::max_string
integer, parameter max_string
Definition: ufo_avgkernel_mod.F90:17
ufo_constants_mod::m_to_km
real(kind_real), parameter, public m_to_km
Definition: ufo_constants_mod.F90:65
ufo_vars_mod::var_sfc_vegtyp
character(len=maxvarlen), parameter, public var_sfc_vegtyp
Definition: ufo_variables_mod.F90:67
ufo_utils_mod::ops_qsat
subroutine, public ops_qsat(QS, T, P, NPNTS)
Definition: ufo_utils_mod.F90:36
ufo_radiancerttov_utils_mod::rttov_absorbers
character(len= *), dimension(ngases_max+1), parameter rttov_absorbers
Definition: ufo_radiancerttov_utils_mod.F90:72
ufo_vars_mod::var_v
character(len=maxvarlen), parameter, public var_v
Definition: ufo_variables_mod.F90:24
ufo_vars_mod::var_sfc_vegfrac
character(len=maxvarlen), parameter, public var_sfc_vegfrac
Definition: ufo_variables_mod.F90:60
ufo_utils_mod
Definition: ufo_utils_mod.F90:1
ufo_vars_mod::var_radiance
character(len=maxvarlen), parameter, public var_radiance
Definition: ufo_variables_mod.F90:77
ufo_basis_mod
Definition: ufo_basis_mod.F90:6
ufo_radiancerttov_utils_mod::nchan_sim
integer, public nchan_sim
Definition: ufo_radiancerttov_utils_mod.F90:98
ufo_vars_mod::var_sfc_soilm
character(len=maxvarlen), parameter, public var_sfc_soilm
Definition: ufo_variables_mod.F90:64
ufo_radiancerttov_utils_mod::ufo_absorbers
character(len=maxvarlen), dimension(ngases_max+1), parameter ufo_absorbers
Definition: ufo_radiancerttov_utils_mod.F90:91
ufo_vars_mod::var_surf_type_rttov
character(len=maxvarlen), parameter, public var_surf_type_rttov
Definition: ufo_variables_mod.F90:92
ufo_vars_mod::var_sfc_wdir
character(len=maxvarlen), parameter, public var_sfc_wdir
Definition: ufo_variables_mod.F90:62
ufo_radiancerttov_utils_mod::debug
logical, public debug
Definition: ufo_radiancerttov_utils_mod.F90:100
ufo_radiancerttov_utils_mod::load_atm_data_rttov
subroutine, public load_atm_data_rttov(geovals, obss, profiles, prof_start, conf, layer_quantities, obs_info)
Definition: ufo_radiancerttov_utils_mod.F90:318
ufo_gnssro_bndnbam_tlad_mod::missing
real(c_double) missing
Definition: ufo_gnssro_bndnbam_tlad_mod.F90:23
ufo_vars_mod::var_lvl_weightfunc
character(len=maxvarlen), parameter, public var_lvl_weightfunc
Definition: ufo_variables_mod.F90:81
ufo_vars_mod::var_sfc_ifrac
character(len=maxvarlen), parameter, public var_sfc_ifrac
Definition: ufo_variables_mod.F90:53
ufo_radiancerttov_utils_mod::rttov_conf
Definition: ufo_radiancerttov_utils_mod.F90:126
ufo_constants_mod::half
real(kind_real), parameter, public half
Definition: ufo_constants_mod.F90:26
ufo_radiancerttov_utils_mod::ufo_rttov_init_emissivity
subroutine ufo_rttov_init_emissivity(self, conf)
Definition: ufo_radiancerttov_utils_mod.F90:1314
ufo_radiancerttov_utils_mod::message
character(len=max_string), public message
Definition: ufo_radiancerttov_utils_mod.F90:45
ufo_constants_mod::g_to_kg
real(kind_real), parameter, public g_to_kg
Definition: ufo_constants_mod.F90:64
ufo_vars_mod::var_sfc_soilt
character(len=maxvarlen), parameter, public var_sfc_soilt
Definition: ufo_variables_mod.F90:65
ufo_constants_mod::zerodegc
real(kind_real), parameter, public zerodegc
Definition: ufo_constants_mod.F90:68
ufo_vars_mod::var_sfc_wspeed
character(len=maxvarlen), parameter, public var_sfc_wspeed
Definition: ufo_variables_mod.F90:61
ufo_radiancerttov_utils_mod::alloc_rttov
subroutine alloc_rttov(self, errorstatus, conf, nprofiles, nchannels, nlevels, init, asw)
Definition: ufo_radiancerttov_utils_mod.F90:1171
ufo_vars_mod::var_sfc_emiss
character(len=maxvarlen), parameter, public var_sfc_emiss
Definition: ufo_variables_mod.F90:74
ufo_radiancerttov_utils_mod::load_geom_data_rttov
subroutine, public load_geom_data_rttov(obss, profiles, prof_start1, obs_info)
Definition: ufo_radiancerttov_utils_mod.F90:779
ufo_radiancerttov_utils_mod::ufo_rttov_skip_profiles
subroutine, public ufo_rttov_skip_profiles(nprofiles, nchannels, channels, obss, Skip_Profiles)
Definition: ufo_radiancerttov_utils_mod.F90:1272
ufo_vars_mod::var_water_type_rttov
character(len=maxvarlen), parameter, public var_water_type_rttov
Definition: ufo_variables_mod.F90:91
ufo_radiancerttov_utils_mod::rttov_absorber_id
integer, dimension(ngases_max+1), parameter rttov_absorber_id
Definition: ufo_radiancerttov_utils_mod.F90:76
ufo_vars_mod::var_lvl_transmit
character(len=maxvarlen), parameter, public var_lvl_transmit
Definition: ufo_variables_mod.F90:80
ufo_vars_mod::var_sfc_wfrac
character(len=maxvarlen), parameter, public var_sfc_wfrac
Definition: ufo_variables_mod.F90:51
ufo_radiancerttov_utils_mod::rttov_errorstatus
integer, public rttov_errorstatus
Definition: ufo_radiancerttov_utils_mod.F90:48
ufo_radiancerttov_utils_mod::rttov_conf_setup
subroutine, public rttov_conf_setup(conf, f_confOpts, f_confOper)
Definition: ufo_radiancerttov_utils_mod.F90:161
ufo_vars_mod::ufo_vars_getindex
integer function, public ufo_vars_getindex(vars, varname)
Definition: ufo_variables_mod.F90:204
ufo_radiancerttov_utils_mod::xstr_diags
character(len=maxvarlen), dimension(:), allocatable xstr_diags
Definition: ufo_radiancerttov_utils_mod.F90:50
ufo_vars_mod::var_mixr
character(len=maxvarlen), parameter, public var_mixr
Definition: ufo_variables_mod.F90:21
ufo_vars_mod::var_sfc_p2m
character(len=maxvarlen), parameter, public var_sfc_p2m
Definition: ufo_variables_mod.F90:47
ufo_radiancerttov_utils_mod::null_str
character(len=maxvarlen), parameter null_str
Definition: ufo_radiancerttov_utils_mod.F90:88
ufo_geovals_mod
Definition: ufo_geovals_mod.F90:7
ufo_radiancerttov_utils_mod
Fortran module to provide code shared between nonlinear and tlm/adm radiance calculations.
Definition: ufo_radiancerttov_utils_mod.F90:8
ufo_radiancerttov_utils_mod::ufo_rttov_io
Definition: ufo_radiancerttov_utils_mod.F90:104
ufo_radiancerttov_utils_mod::get_var_name
subroutine get_var_name(n, varname)
Definition: ufo_radiancerttov_utils_mod.F90:1158
ufo_radiancerttov_utils_mod::ch_diags
integer, dimension(:), allocatable ch_diags
Definition: ufo_radiancerttov_utils_mod.F90:51
ufo_vars_mod::var_prsi
character(len=maxvarlen), parameter, public var_prsi
Definition: ufo_variables_mod.F90:26
ufo_radiancerttov_utils_mod::maxvarin
integer, parameter, public maxvarin
Definition: ufo_radiancerttov_utils_mod.F90:40
ufo_vars_mod::var_sfc_itmp
character(len=maxvarlen), parameter, public var_sfc_itmp
Definition: ufo_variables_mod.F90:57
ufo_radiancerttov_utils_mod::nchan_inst
integer, public nchan_inst
Definition: ufo_radiancerttov_utils_mod.F90:97
ufo_constants_mod
Definition: ufo_constants_mod.F90:2
ufo_vars_mod::var_sfc_wtmp
character(len=maxvarlen), parameter, public var_sfc_wtmp
Definition: ufo_variables_mod.F90:55
ufo_vars_mod::var_u
character(len=maxvarlen), parameter, public var_u
Definition: ufo_variables_mod.F90:23
ufo_vars_mod::var_sfc_q2m
character(len=maxvarlen), parameter, public var_sfc_q2m
Definition: ufo_variables_mod.F90:48
ufo_vars_mod::var_sfc_soiltyp
character(len=maxvarlen), parameter, public var_sfc_soiltyp
Definition: ufo_variables_mod.F90:68
ufo_radiancerttov_utils_mod::setup_rttov
subroutine setup_rttov(self, f_confOpts, asw)
Definition: ufo_radiancerttov_utils_mod.F90:1105
ufo_vars_mod::var_sfc_sfrac
character(len=maxvarlen), parameter, public var_sfc_sfrac
Definition: ufo_variables_mod.F90:54
ufo_radiancerttov_utils_mod::iprof
integer iprof
Definition: ufo_radiancerttov_utils_mod.F90:102
ufo_vars_mod::var_sfc_lfrac
character(len=maxvarlen), parameter, public var_sfc_lfrac
Definition: ufo_variables_mod.F90:52
ufo_constants_mod::min_q
real(kind_real), parameter, public min_q
Definition: ufo_constants_mod.F90:67
ufo_basis_mod::ufo_basis
Definition: ufo_basis_mod.F90:12
ufo_aodcrtm_mod::varin_default
character(len=maxvarlen), dimension(5), parameter varin_default
Definition: ufo_aodcrtm_mod.F90:35
ufo_vars_mod::var_sfc_landtyp
character(len=maxvarlen), parameter, public var_sfc_landtyp
Definition: ufo_variables_mod.F90:66
ufo_vars_mod
Definition: ufo_variables_mod.F90:8
ufo_constants_mod::pa_to_hpa
real(kind_real), parameter, public pa_to_hpa
Definition: ufo_constants_mod.F90:66
ufo_radiancerttov_utils_mod::varin_temp
character(len=maxvarlen), dimension(maxvarin), public varin_temp
Definition: ufo_radiancerttov_utils_mod.F90:44
ufo_vars_mod::var_tb
character(len=maxvarlen), parameter, public var_tb
Definition: ufo_variables_mod.F90:78
ufo_radiancerttov_utils_mod::varin_default_satrad
character(len=maxvarlen), dimension(11), target varin_default_satrad
Definition: ufo_radiancerttov_utils_mod.F90:64
ufo_radiancerttov_utils_mod::ops_satrad_qsplit
subroutine ops_satrad_qsplit(output_type, p, t, qtotal, q, ql, qi, UseQtSplitRain)
Definition: ufo_radiancerttov_utils_mod.F90:1809
ufo_radiancerttov_utils_mod::varin_default_crtm
character(len=maxvarlen), dimension(21), target varin_default_crtm
Definition: ufo_radiancerttov_utils_mod.F90:55
ufo_vars_mod::var_sfc_lai
character(len=maxvarlen), parameter, public var_sfc_lai
Definition: ufo_variables_mod.F90:63
ufo_geovals_mod::ufo_geovals_get_var
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Definition: ufo_geovals_mod.F90:128
ufo_radiancerttov_utils_mod::nvars_in
integer, public nvars_in
Definition: ufo_radiancerttov_utils_mod.F90:47
ufo_vars_mod::var_q
character(len=maxvarlen), parameter, public var_q
Definition: ufo_variables_mod.F90:22
ufo_vars_mod::var_tb_clr
character(len=maxvarlen), parameter, public var_tb_clr
Definition: ufo_variables_mod.F90:79
ufo_vars_mod::var_clw
character(len=maxvarlen), parameter, public var_clw
Definition: ufo_variables_mod.F90:34
ufo_constants_mod::deg2rad
real(kind_real), parameter, public deg2rad
Definition: ufo_constants_mod.F90:50
ufo_vars_mod::var_sfc_t2m
character(len=maxvarlen), parameter, public var_sfc_t2m
Definition: ufo_variables_mod.F90:49
ufo_vars_mod::var_ts
character(len=maxvarlen), parameter, public var_ts
Definition: ufo_variables_mod.F90:19
ufo_vars_mod::var_pmaxlev_weightfunc
character(len=maxvarlen), parameter, public var_pmaxlev_weightfunc
Definition: ufo_variables_mod.F90:82
ufo_vars_mod::var_sfc_stmp
character(len=maxvarlen), parameter, public var_sfc_stmp
Definition: ufo_variables_mod.F90:58
ufo_vars_mod::var_co2
character(len=maxvarlen), parameter, public var_co2
Definition: ufo_variables_mod.F90:33
ufo_radiancerttov_utils_mod::nlocs_total
integer, public nlocs_total
Definition: ufo_radiancerttov_utils_mod.F90:99
ufo_radiancerttov_utils_mod::set_defaults_rttov
subroutine set_defaults_rttov(self, default_opts_set)
Definition: ufo_radiancerttov_utils_mod.F90:1374
ufo_radiancerttov_utils_mod::populate_hofxdiags
subroutine, public populate_hofxdiags(RTProf, chanprof, hofxdiags)
Definition: ufo_radiancerttov_utils_mod.F90:1547
ufo_vars_mod::var_oz
character(len=maxvarlen), parameter, public var_oz
Definition: ufo_variables_mod.F90:32
ufo_vars_mod::var_sfc_sdepth
character(len=maxvarlen), parameter, public var_sfc_sdepth
Definition: ufo_variables_mod.F90:59
ufo_vars_mod::var_sfc_tskin
character(len=maxvarlen), parameter, public var_sfc_tskin
Definition: ufo_variables_mod.F90:50
ufo_geovals_mod::ufo_geovals
type to hold interpolated fields required by the obs operators
Definition: ufo_geovals_mod.F90:47
ufo_radiancerttov_utils_mod::ystr_diags
character(len=maxvarlen), dimension(:), allocatable ystr_diags
Definition: ufo_radiancerttov_utils_mod.F90:50
ufo_vars_mod::var_opt_depth
character(len=maxvarlen), parameter, public var_opt_depth
Definition: ufo_variables_mod.F90:76
ufo_radiancerttov_utils_mod::rttov_conf_delete
subroutine, public rttov_conf_delete(conf)
Definition: ufo_radiancerttov_utils_mod.F90:305
ufo_constants_mod::zero
real(kind_real), parameter, public zero
Definition: ufo_constants_mod.F90:24
ufo_geovals_mod::ufo_geoval
type to hold interpolated field for one variable, one observation
Definition: ufo_geovals_mod.F90:40
conf
Definition: conf.py:1
ufo_vars_mod::var_prs
character(len=maxvarlen), parameter, public var_prs
Definition: ufo_variables_mod.F90:25
ufo_vars_mod::var_sfc_ltmp
character(len=maxvarlen), parameter, public var_sfc_ltmp
Definition: ufo_variables_mod.F90:56
ufo_radiancerttov_utils_mod::parse_hofxdiags
subroutine, public parse_hofxdiags(hofxdiags, jacobian_needed)
Definition: ufo_radiancerttov_utils_mod.F90:1726
ufo_radiancerttov_utils_mod::set_options_rttov
subroutine set_options_rttov(self, f_confOpts)
Definition: ufo_radiancerttov_utils_mod.F90:852