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