UFO
ufo_crtm_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 iso_c_binding
12 use kinds
13 
14 use crtm_module
15 
16 use ufo_vars_mod
18 use obsspace_mod
19 
20 implicit none
21 private
22 
24 public crtm_conf
25 public crtm_conf_setup
26 public crtm_conf_delete
27 public get_var_name
28 public load_atm_data
29 public load_sfc_data
30 public load_geom_data
32 
33 PUBLIC load_aerosol_data
36 PUBLIC :: qsmith
37 PUBLIC :: upper2lower
38 
39 PUBLIC :: grav,rv_rd,&
41 
42 REAL(kind_real), PARAMETER :: &
43  &rdgas = 2.8704e+2_kind_real,&
44  &rvgas = 4.6150e+2_kind_real,&
45  &rv_rd = rvgas/rdgas,&
46  &esl = 0.621971831,&
47  &zvir = rv_rd - 1_kind_real,&
48  &tice = 273.16_kind_real,&
49  &grav = 9.81_kind_real,&
50  &aerosol_concentration_minvalue=1.e-16_kind_real,&
51  &aerosol_concentration_minvalue_layer=tiny(rdgas),&
52  &ozone_default_value=1.e-3_kind_real ! in ppmv in crtm
53 
54 integer, parameter, public :: max_string=800
55 
56 !Type for general config
58  integer :: n_sensors
59  integer :: n_absorbers
60  integer :: n_clouds
61  integer :: n_aerosols
62  integer :: n_surfaces
63  character(len=MAXVARLEN), allocatable :: absorbers(:)
64  integer, allocatable :: absorber_id(:)
65  integer, allocatable :: absorber_units(:)
66  character(len=MAXVARLEN), allocatable :: clouds(:,:)
67  integer, allocatable :: cloud_id(:)
68  character(len=MAXVARLEN), allocatable :: surfaces(:)
69 
70  character(len=255), allocatable :: sensor_id(:)
71  character(len=255) :: endian_type
72  character(len=255) :: coefficient_path
73  character(len=255) :: &
74  irwatercoeff_file, irlandcoeff_file, irsnowcoeff_file, iricecoeff_file, &
75  viswatercoeff_file, vislandcoeff_file, vissnowcoeff_file, visicecoeff_file, &
76  mwwatercoeff_file
77  integer, allocatable :: land_wsi(:)
78  real(kind_real) :: cloud_fraction = -1.0_kind_real
79  integer :: inspect
80  character(len=MAXVARLEN) :: aerosol_option
81  character(len=255) :: salinity_option
82 end type crtm_conf
83 
85 
88 
89 END INTERFACE
90 
91 INTERFACE qsmith
92 
93  MODULE PROCEDURE qsmith_atm,qsmith_profiles
94 
95 END INTERFACE qsmith
96 
97 
98  ! Add more UFO_Absorbers as needed
99  ! Note: must have same ordering as CRTM_Absorbers,
100  ! CRTM_Absorber_Id, and CRTM_Absorber_Units
101  character(len=MAXVARLEN), parameter :: &
102  ufo_absorbers(3) = &
103  [ var_mixr, var_co2, var_oz ]
104 
105  ! copy of ABSORBER_ID_NAME defined in CRTM_Atmosphere_Define
106  character(len=*), parameter :: &
107  crtm_absorbers(n_valid_absorber_ids) = &
108  absorber_id_name(1:n_valid_absorber_ids)
109  integer, parameter :: &
110  crtm_absorber_id(n_valid_absorber_ids) = &
111  [ h2o_id, co2_id, o3_id, n2o_id, &
112  co_id, ch4_id, o2_id, no_id, &
113  so2_id, no2_id, nh3_id, hno3_id, &
114  oh_id, hf_id, hcl_id, hbr_id, &
115  hi_id, clo_id, ocs_id, h2co_id, &
116  hocl_id, n2_id, hcn_id, ch3l_id, &
117  h2o2_id, c2h2_id, c2h6_id, ph3_id, &
118  cof2_id, sf6_id, h2s_id,hcooh_id ]
119  integer, parameter :: &
120  crtm_absorber_units(3) = [ &
121  mass_mixing_ratio_units & !H2O
122  , volume_mixing_ratio_units & !CO2
123  , volume_mixing_ratio_units & !O3
124  ]
125 
126  character(len=MAXVARLEN), parameter :: &
127  ufo_clouds(n_valid_cloud_categories,2) = &
128  reshape( &
131  , [n_valid_cloud_categories,2] )
132 
133  ! copy of CLOUD_CATEGORY_NAME defined in CRTM_Cloud_Define
134  character(len=*), parameter :: &
135  crtm_clouds(n_valid_cloud_categories) = &
136  cloud_category_name(1:n_valid_cloud_categories)
137  integer, parameter :: &
138  crtm_cloud_id(n_valid_cloud_categories) = &
139  [ water_cloud, &
140  ice_cloud, &
141  rain_cloud, &
142  snow_cloud, &
143  graupel_cloud, &
144  hail_cloud ]
145 
146 ! Surface Variables
147 
148  character(len=MAXVARLEN), parameter :: &
149  ufo_surfaces(4) = &
151 
152  character(len=MAXVARLEN), parameter :: &
153  crtm_surfaces(4) = &
154  [ character(len=maxvarlen):: 'Water_Temperature', 'Wind_Speed', 'Wind_Direction', 'Salinity' ]
155 
156 contains
157 
158 ! ------------------------------------------------------------------------------
159 
160 subroutine crtm_conf_setup(conf, f_confOpts, f_confOper)
161 
162 implicit none
163 type(crtm_conf), intent(inout) :: conf
164 type(fckit_configuration), intent(in) :: f_confopts
165 type(fckit_configuration), intent(in) :: f_confoper
166 
167 character(*), PARAMETER :: routine_name = 'crtm_conf_setup'
168 character(len=255) :: irwatercoeff, viswatercoeff, &
169  irvislandcoeff, irvissnowcoeff, irvisicecoeff, &
170  mwwatercoeff
171 integer :: jspec, ivar
172 character(len=max_string) :: message
173 character(len=:), allocatable :: str
174 character(len=:), allocatable :: str_array(:)
175 
176 CHARACTER(len=MAXVARLEN), ALLOCATABLE :: var_aerosols(:)
177 
178 
179  !Some config needs to come from user
180  !-----------------------------------
181 
182  !Number of sensors, each call to CRTM will be for a single sensor
183  !type (zenith/scan angle will be different)
184  conf%n_Sensors = 1
185 
186  ! absorbers, clouds and aerosols (should match what model will provide)
187 
188  ! Absorbers
189  !----------
190  conf%n_Absorbers = 0
191  if (f_confoper%has("Absorbers")) &
192  conf%n_Absorbers = conf%n_Absorbers + f_confoper%get_size("Absorbers")
193 
194  allocate( conf%Absorbers ( conf%n_Absorbers ), &
195  conf%Absorber_Id ( conf%n_Absorbers ), &
196  conf%Absorber_Units( conf%n_Absorbers ) )
197 
198  if (conf%n_Absorbers > 0) then
199  call f_confoper%get_or_die("Absorbers",str_array)
200  conf%Absorbers(1:conf%n_Absorbers) = str_array
201  end if
202 
203  ! check for duplications
204  do jspec = 2, conf%n_Absorbers
205  if ( any(conf%Absorbers(jspec-1) == conf%Absorbers(jspec:conf%n_Absorbers)) ) then
206  write(message,*) trim(routine_name),' error: ',trim(conf%Absorbers(jspec)),' is duplicated in Absorbers'
207  call abor1_ftn(message)
208  end if
209  end do
210 
211  ! convert from CRTM names to UFO CF names and define Id and Units
212  do jspec = 1, conf%n_Absorbers
213  ivar = ufo_vars_getindex(crtm_absorbers, conf%Absorbers(jspec))
214  if (ivar < 1 .or. ivar > size(ufo_absorbers)) then
215  write(message,*) trim(routine_name),' error: ',trim(conf%Absorbers(jspec)),' not supported by UFO_Absorbers'
216  call abor1_ftn(message)
217  end if
218  conf%Absorbers(jspec) = ufo_absorbers(ivar)
219  conf%Absorber_Id(jspec) = crtm_absorber_id(ivar)
220  conf%Absorber_Units(jspec) = crtm_absorber_units(ivar)
221  end do
222 
223 
224  ! Clouds
225  !-------
226  conf%n_Clouds = 0
227  if (f_confoper%has("Clouds")) &
228  conf%n_Clouds = f_confoper%get_size("Clouds")
229  allocate( conf%Clouds ( conf%n_Clouds,2), &
230  conf%Cloud_Id( conf%n_Clouds ) )
231  if (conf%n_Clouds > 0) then
232  call f_confoper%get_or_die("Clouds",str_array)
233  conf%Clouds(1:conf%n_Clouds,1) = str_array
234 
235  if (f_confoper%has("Cloud_Fraction")) then
236  call f_confoper%get_or_die("Cloud_Fraction",conf%Cloud_Fraction)
237  if ( conf%Cloud_Fraction < 0.0 .or. &
238  conf%Cloud_Fraction > 1.0 ) then
239  write(message,*) trim(routine_name),' error: must specify ' // &
240  ' 0.0 <= Cloud_Fraction <= 1.0' // &
241  ' or remove Cloud_Fraction from conf' // &
242  ' and provide as a geoval'
243  call abor1_ftn(message)
244  end if
245  else
246  message = trim(routine_name) // &
247  ': Cloud_Fraction is not provided in conf.' // &
248  ' Will request as a geoval.'
249  CALL display_message(routine_name, trim(message), warning )
250  end if
251  end if
252 
253  ! check for duplications
254  do jspec = 2, conf%n_Clouds
255  if ( any(conf%Clouds(jspec-1,1) == conf%Clouds(jspec:conf%n_Clouds,1)) ) then
256  write(message,*) trim(routine_name),' error: ',trim(conf%Clouds(jspec,1)), &
257  ' is duplicated in Clouds'
258  call abor1_ftn(message)
259  end if
260  end do
261 
262  ! convert from CRTM names to UFO CF names and define Id
263  do jspec = 1, conf%n_Clouds
264  ivar = ufo_vars_getindex(crtm_clouds, conf%Clouds(jspec,1))
265  if (ivar < 1 .or. ivar > size(ufo_clouds)) then
266  write(message,*) trim(routine_name),' error: ',trim(conf%Clouds(jspec,1)),' not supported by UFO_Clouds'
267  call abor1_ftn(message)
268  end if
269  conf%Clouds(jspec,1:2) = ufo_clouds(ivar,1:2)
270  conf%Cloud_Id(jspec) = crtm_cloud_id(ivar)
271  end do
272 
273  ! Aerosols
274  !---------
275  IF (f_confopts%has("AerosolOption")) THEN
276  call f_confopts%get_or_die("AerosolOption",str)
277  conf%aerosol_option = str
278  conf%aerosol_option = upper2lower(conf%aerosol_option)
279  CALL assign_aerosol_names(conf%aerosol_option,var_aerosols)
280  conf%n_Aerosols=SIZE(var_aerosols)
281  ELSE
282  conf%n_Aerosols = 0
283  conf%aerosol_option = ""
284  ENDIF
285 
286  ! Surface variables
287  !----------
288  conf%n_Surfaces = 0
289  if (f_confoper%has("Surfaces")) &
290  conf%n_Surfaces = conf%n_Surfaces + f_confoper%get_size("Surfaces")
291 
292  allocate( conf%Surfaces ( conf%n_Surfaces ))
293 
294  if (conf%n_Surfaces > 0) then
295  call f_confoper%get_or_die("Surfaces",str_array)
296  conf%Surfaces(1:conf%n_Surfaces) = str_array
297  end if
298 
299  ! check for duplications
300  do jspec = 2, conf%n_Surfaces
301  if ( any(conf%Surfaces(jspec-1) == conf%Surfaces(jspec:conf%n_Surfaces)) ) then
302  write(message,*) 'crtm_conf_setup error: ',trim(conf%Surfaces(jspec)),' is duplicated in Surfaces'
303  call abor1_ftn(message)
304  end if
305  end do
306 
307  ! convert from CRTM names to UFO CF names and define Id and Units
308  do jspec = 1, conf%n_Surfaces
309  ivar = ufo_vars_getindex(crtm_surfaces, conf%Surfaces(jspec))
310  if (ivar < 1 .or. ivar > size(ufo_surfaces)) then
311  write(message,*) 'crtm_conf_setup error: ',trim(conf%Surfaces(jspec)),' not supported by UFO_Surfaces'
312  call abor1_ftn(message)
313  end if
314  conf%Surfaces(jspec) = ufo_surfaces(ivar)
315 
316  end do
317 
318  ! Sea_Surface_Salinity
319  !---------
320  IF (f_confopts%get("Salinity",str)) THEN
321  conf%salinity_option = str
322  ELSE
323  conf%salinity_option = 'off'
324  END IF
325 
326  !Allocate SENSOR_ID
327  allocate(conf%SENSOR_ID(conf%n_Sensors))
328 
329  !Get sensor ID from config
330  call f_confopts%get_or_die("Sensor_ID",str)
331  conf%SENSOR_ID(conf%n_Sensors) = str
332 
333  !ENDIAN type
334  call f_confopts%get_or_die("EndianType",str)
335  conf%ENDIAN_TYPE = str
336 
337  !Path to coefficient files
338  call f_confopts%get_or_die("CoefficientPath",str)
339  conf%COEFFICIENT_PATH = str
340 
341  ! Coefficient file prefixes
342  irwatercoeff = "Nalli"
343  if (f_confopts%has("IRwaterCoeff")) then
344  call f_confopts%get_or_die("IRwaterCoeff",str)
345  irwatercoeff = str
346  end if
347  viswatercoeff = "NPOESS"
348  if (f_confopts%has("VISwaterCoeff")) then
349  call f_confopts%get_or_die("VISwaterCoeff",str)
350  viswatercoeff = str
351  end if
352  irvislandcoeff = "NPOESS"
353  if (f_confopts%has("IRVISlandCoeff")) then
354  call f_confopts%get_or_die("IRVISlandCoeff",str)
355  irvislandcoeff = str
356  end if
357  irvissnowcoeff = "NPOESS"
358  if (f_confopts%has("IRVISsnowCoeff")) then
359  call f_confopts%get_or_die("IRVISsnowCoeff",str)
360  irvissnowcoeff = str
361  end if
362  irvisicecoeff = "NPOESS"
363  if (f_confopts%has("IRVISiceCoeff")) then
364  call f_confopts%get_or_die("IRVISiceCoeff",str)
365  irvisicecoeff = str
366  end if
367  mwwatercoeff = "FASTEM6"
368  if (f_confopts%has("MWwaterCoeff")) then
369  call f_confopts%get_or_die("MWwaterCoeff",str)
370  mwwatercoeff = str
371  end if
372 
373  ! Define water, snow, ice (WSI) categories
374  select case (trim(irvislandcoeff))
375  case ("USGS")
376  allocate(conf%Land_WSI(2))
377  conf%Land_WSI(1:2) = (/16,24/)
378  case ("IGBP")
379  allocate(conf%Land_WSI(2))
380  conf%Land_WSI(1:2) = (/15,17/)
381  case default
382  allocate(conf%Land_WSI(1))
383  conf%Land_WSI(1) = -1
384  end select
385 
386  ! IR emissivity coeff files
387  conf%IRwaterCoeff_File = trim(irwatercoeff)//".IRwater.EmisCoeff.bin"
388  conf%IRlandCoeff_File = trim(irvislandcoeff)//".IRland.EmisCoeff.bin"
389  conf%IRsnowCoeff_File = trim(irvissnowcoeff)//".IRsnow.EmisCoeff.bin"
390  conf%IRiceCoeff_File = trim(irvisicecoeff)//".IRice.EmisCoeff.bin"
391 
392  !VIS emissivity coeff files
393  conf%VISwaterCoeff_File = trim(viswatercoeff)//".VISwater.EmisCoeff.bin"
394  conf%VISlandCoeff_File = trim(irvislandcoeff)//".VISland.EmisCoeff.bin"
395  conf%VISsnowCoeff_File = trim(irvissnowcoeff)//".VISsnow.EmisCoeff.bin"
396  conf%VISiceCoeff_File = trim(irvisicecoeff)//".VISice.EmisCoeff.bin"
397 
398  ! MW water emissivity coeff file
399  conf%MWwaterCoeff_File = trim(mwwatercoeff)//".MWwater.EmisCoeff.bin"
400 
401  conf%inspect = 0
402  if (f_confopts%has("InspectProfileNumber")) then
403  call f_confopts%get_or_die("InspectProfileNumber",conf%inspect)
404  endif
405 
406 end subroutine crtm_conf_setup
407 
408 ! -----------------------------------------------------------------------------
409 
410 subroutine crtm_conf_delete(conf)
411 
412 implicit none
413 type(crtm_conf), intent(inout) :: conf
414 
415  deallocate(conf%SENSOR_ID)
416  deallocate(conf%Land_WSI)
417  deallocate(conf%Absorbers)
418  deallocate(conf%Absorber_Id)
419  deallocate(conf%Absorber_Units)
420  deallocate(conf%Clouds)
421  deallocate(conf%Cloud_Id)
422 
423 end subroutine crtm_conf_delete
424 
425 ! ------------------------------------------------------------------------------
426 
427 subroutine crtm_comm_stat_check(stat, PROGRAM_NAME, message, f_comm)
428 use fckit_mpi_module, only: fckit_mpi_comm
429 
430 implicit none
431 
432 integer, intent(in) :: stat
433 character(*), intent(in) :: program_name
434 character(*), intent(in) :: message
435 type(fckit_mpi_comm), intent(in) :: f_comm
436 
437 character(max_string) :: rank_message
438 
439  if ( stat /= success ) THEN
440  write(rank_message,*) trim(message)," on rank ",f_comm%rank()
441  call display_message( program_name, rank_message, failure )
442  call abor1_ftn("Abort from "//program_name)
443  end if
444 
445 end subroutine crtm_comm_stat_check
446 
447 ! ------------------------------------------------------------------------------
448 
449 subroutine ufo_crtm_skip_profiles(n_Profiles,n_Channels,channels,obss,Skip_Profiles)
450 ! Profiles are skipped when the ObsValue of all channels is missing.
451 ! TODO: Use complete QC information
452 ! It would be more comprehensive to use EffectiveQC or EffectiveError. That
453 ! would require those ObsSpace values to be initialized before calls to
454 ! this subroutine within ufo_radiancecrtm_simobs+ufo_radiancecrtm_tlad_settraj.
455 use missing_values_mod
456 
457 implicit none
458 integer, intent(in) :: n_profiles, n_channels
459 type(c_ptr), value, intent(in) :: obss
460 integer(c_int), intent(in) :: channels(:)
461 logical, intent(inout) :: skip_profiles(:)
462 
463 integer :: jprofile, jchannel
464 character(len=MAXVARLEN) :: varname
465 real(kind_real) :: obsval(n_profiles,n_channels)
466 !real(kind_real) :: EffObsErr(n_Profiles,n_Channels)
467 !integer :: EffQC(n_Profiles,n_Channels)
468 
469 real(c_double) :: missing
470 
471  ! Set missing value
472  missing = missing_value(missing)
473 
474  obsval = missing
475 ! EffObsErr = missing
476 ! EffQC = 0
477 
478  do jchannel = 1, n_channels
479  call get_var_name(channels(jchannel),varname)
480  call obsspace_get_db(obss, "ObsValue", varname, obsval(:,jchannel))
481 ! call obsspace_get_db(obss, "EffectiveError", varname, EffObsErr(:,jchannel))
482 ! call obsspace_get_db(obss, "EffectiveQC{iter}", varname, EffQC(:,jchannel))
483  enddo
484 
485  !Loop over all n_Profiles, i.e. number of locations
486  do jprofile = 1, n_profiles
487  skip_profiles(jprofile) = all(obsval(jprofile,:) == missing)
488 ! .OR. all(EffObsErr(jprofile,:) == missing) &
489 ! .OR. all(EffQC(jprofile,:) /= 0)
490  end do
491 
492 end subroutine ufo_crtm_skip_profiles
493 
494 ! ------------------------------------------------------------------------------
495 
496 SUBROUTINE load_atm_data(n_Profiles, n_Layers, geovals, atm, conf)
497 
498 implicit none
499 
500 integer, intent(in) :: n_profiles, n_layers
501 type(ufo_geovals), intent(in) :: geovals
502 type(crtm_atmosphere_type), intent(inout) :: atm(:)
503 type(crtm_conf) :: conf
504 
505 ! Local variables
506 integer :: k1, jspec
507 type(ufo_geoval), pointer :: geoval
508 character(max_string) :: err_msg
509 
510  ! Populate the atmosphere structures for CRTM
511  ! -------------------------------------------
512 
513  call ufo_geovals_get_var(geovals, var_ts, geoval)
514  ! Check model levels is consistent in geovals & crtm
515  if (geoval%nval /= n_layers) then
516  write(err_msg,*) 'Load_Atm_Data error: layers inconsistent!'
517  call abor1_ftn(err_msg)
518  endif
519 
520  do k1 = 1, n_profiles
521  atm(k1)%Temperature(1:n_layers) = geoval%vals(:, k1)
522  end do
523 
524  call ufo_geovals_get_var(geovals, var_prs, geoval)
525  do k1 = 1, n_profiles
526  atm(k1)%Pressure(1:n_layers) = geoval%vals(:, k1) * 0.01 ! to hPa
527  end do
528 
529  call ufo_geovals_get_var(geovals, var_prsi, geoval)
530  do k1 = 1, n_profiles
531  atm(k1)%Level_Pressure(:) = geoval%vals(:, k1) * 0.01 ! to hPa
532  atm(k1)%Climatology = us_standard_atmosphere
533  end do
534 
535  do jspec = 1, conf%n_Absorbers
536  ! O3 Absorber has special treatment for Aerosols
537  if ( trim(conf%Absorbers(jspec)) == trim(var_oz) .AND. &
538  ufo_vars_getindex(geovals%variables, var_oz) < 0 .AND. &
539  trim(conf%aerosol_option) /= "" ) then
540  do k1 = 1, n_profiles
541  atm(k1)%Absorber(1:n_layers, jspec) = ozone_default_value
542  end do
543  else
544  call ufo_geovals_get_var(geovals, conf%Absorbers(jspec), geoval)
545  do k1 = 1, n_profiles
546  atm(k1)%Absorber(1:n_layers, jspec) = geoval%vals(:, k1)
547  end do
548  end if
549  do k1 = 1, n_profiles
550  atm(k1)%Absorber_Id(jspec) = conf%Absorber_Id(jspec)
551  atm(k1)%Absorber_Units(jspec) = conf%Absorber_Units(jspec)
552  end do
553  end do
554 
555  do jspec = 1, conf%n_Clouds
556  ! cloud species content
557  CALL ufo_geovals_get_var(geovals, conf%Clouds(jspec,1), geoval)
558  do k1 = 1, n_profiles
559  atm(k1)%Cloud(jspec)%Water_Content = geoval%vals(:, k1)
560  atm(k1)%Cloud(jspec)%Type = conf%Cloud_Id(jspec)
561  end do
562 
563  ! effective radius
564  CALL ufo_geovals_get_var(geovals, conf%Clouds(jspec,2), geoval)
565  do k1 = 1, n_profiles
566  atm(k1)%Cloud(jspec)%Effective_Radius = geoval%vals(:, k1)
567  end do
568  end do
569 
570  ! When n_Clouds>0, Cloud_Fraction must either be provided as geoval or in conf
571  if (conf%n_Clouds > 0) then
572  if ( ufo_vars_getindex(geovals%variables, var_cldfrac) > 0 ) then
573  CALL ufo_geovals_get_var(geovals, var_cldfrac, geoval)
574  do k1 = 1, n_profiles
575  atm(k1)%Cloud_Fraction(:) = geoval%vals(:, k1)
576  end do
577  else
578  do k1 = 1, n_profiles
579  atm(k1)%Cloud_Fraction(:) = conf%Cloud_Fraction
580  end do
581  end if
582  end if
583 
584 end subroutine load_atm_data
585 
586 ! ------------------------------------------------------------------------------
587 
588 subroutine load_sfc_data(n_Profiles, n_Channels, channels, geovals, sfc, chinfo, obss, conf)
589 
590 implicit none
591 
592 integer, intent(in) :: n_profiles, n_channels
593 type(ufo_geovals), intent(in) :: geovals
594 type(crtm_surface_type), intent(inout) :: sfc(:)
595 type(crtm_channelinfo_type), intent(in) :: chinfo(:)
596 type(c_ptr), value, intent(in) :: obss
597 integer(c_int), intent(in) :: channels(:)
598 type(crtm_conf), intent(in) :: conf
599 
600 type(ufo_geoval), pointer :: geoval
601 integer :: k1, n1
602 integer :: iland
603 
604 ! Surface type definitions for default SfcOptics definitions
605 ! for IR and VIS, this is the NPOESS reflectivities.
606 integer, parameter :: tundra_surface_type = 10 ! NPOESS Land surface type for IR/VIS Land SfcOptics
607 integer, parameter :: scrub_surface_type = 7 ! NPOESS Land surface type for IR/VIS Land SfcOptics
608 integer, parameter :: coarse_soil_type = 1 ! Soil type for MW land SfcOptics
609 integer, parameter :: groundcover_vegetation_type = 7 ! Vegetation type for MW Land SfcOptics
610 integer, parameter :: bare_soil_vegetation_type = 11 ! Vegetation type for MW Land SfcOptics
611 integer, parameter :: sea_water_type = 1 ! Water type for all SfcOptics
612 integer, parameter :: fresh_snow_type = 2 ! NPOESS Snow type for IR/VIS SfcOptics
613 integer, parameter :: fresh_ice_type = 1 ! NPOESS Ice type for IR/VIS SfcOptics
614 
615 character(len=MAXVARLEN) :: varname
616 
617 real(kind_real), allocatable :: obstb(:,:)
618 
619  allocate(obstb(n_profiles, n_channels))
620  obstb = 0.0_kind_real
621 
622  do n1 = 1, n_channels
623  call get_var_name(channels(n1), varname)
624  call obsspace_get_db(obss, "ObsValue", varname, obstb(:, n1))
625  enddo
626 
627  do k1 = 1, n_profiles
628  !Pass sensor information
629  sfc(k1)%sensordata%sensor_id = chinfo(1)%sensor_id
630  sfc(k1)%sensordata%wmo_sensor_id = chinfo(1)%wmo_sensor_id
631  sfc(k1)%sensordata%wmo_satellite_id = chinfo(1)%wmo_satellite_id
632  sfc(k1)%sensordata%sensor_channel = channels
633 
634  !Pass observation value
635  do n1 = 1, n_channels
636  sfc(k1)%sensordata%tb(n1) = obstb(k1, n1)
637  enddo
638 
639  !Water_type
640  sfc(k1)%Water_Type = sea_water_type !** NOTE: need to check how to determine fresh vs sea water types (salinity???)
641  end do
642  deallocate(obstb)
643 
644  !Wind_Speed
645  call ufo_geovals_get_var(geovals, var_sfc_wspeed, geoval)
646  do k1 = 1, n_profiles
647  sfc(k1)%Wind_Speed = geoval%vals(1, k1)
648  end do
649 
650  !Wind_Direction
651  call ufo_geovals_get_var(geovals, var_sfc_wdir, geoval)
652  do k1 = 1, n_profiles
653  sfc(k1)%Wind_Direction = geoval%vals(1, k1)
654  end do
655 
656  !Water_Coverage
657  call ufo_geovals_get_var(geovals, var_sfc_wfrac, geoval)
658  do k1 = 1, n_profiles
659  sfc(k1)%Water_Coverage = geoval%vals(1, k1)
660  end do
661 
662  !Water_Temperature
663  call ufo_geovals_get_var(geovals, var_sfc_wtmp, geoval)
664  do k1 = 1, n_profiles
665  sfc(k1)%Water_Temperature = geoval%vals(1, k1)
666  end do
667 
668  !Ice_Coverage
669  call ufo_geovals_get_var(geovals, var_sfc_ifrac, geoval)
670  do k1 = 1, n_profiles
671  sfc(k1)%Ice_Coverage = geoval%vals(1, k1)
672  end do
673 
674  !Ice_Temperature
675  call ufo_geovals_get_var(geovals, var_sfc_itmp, geoval)
676  do k1 = 1, n_profiles
677  sfc(k1)%Ice_Temperature = geoval%vals(1, k1)
678  end do
679 
680  !Snow_Coverage
681  call ufo_geovals_get_var(geovals, var_sfc_sfrac, geoval)
682  do k1 = 1, n_profiles
683  sfc(k1)%Snow_Coverage = geoval%vals(1, k1)
684  end do
685 
686  !Snow_Temperature
687  call ufo_geovals_get_var(geovals, var_sfc_stmp, geoval)
688  do k1 = 1, n_profiles
689  sfc(k1)%Snow_Temperature = geoval%vals(1, k1)
690  end do
691 
692  !Snow_Depth
693  call ufo_geovals_get_var(geovals, var_sfc_sdepth, geoval)
694  do k1 = 1, n_profiles
695  sfc(k1)%Snow_Depth = geoval%vals(1, k1)
696  end do
697 
698  !Land_Coverage
699  call ufo_geovals_get_var(geovals, var_sfc_lfrac, geoval)
700  do k1 = 1, n_profiles
701  sfc(k1)%Land_Coverage = geoval%vals(1, k1)
702  end do
703 
704  !Land_Type
705  ! + used to lookup land sfc emiss. for IR and VIS
706  ! + land sfc emiss. undefined over water/snow/ice
707  call ufo_geovals_get_var(geovals, var_sfc_landtyp, geoval)
708  do k1 = 1, n_profiles
709  iland = int(geoval%vals(1, k1))
710  if (.not.any(iland == conf%Land_WSI)) then
711  sfc(k1)%Land_Type = iland
712  end if
713  end do
714 
715  !Land_Temperature
716  call ufo_geovals_get_var(geovals, var_sfc_ltmp, geoval)
717  do k1 = 1, n_profiles
718  sfc(k1)%Land_Temperature = geoval%vals(1, k1)
719  end do
720 
721  !Lai
722  call ufo_geovals_get_var(geovals, var_sfc_lai, geoval)
723  do k1 = 1, n_profiles
724  sfc(k1)%Lai = geoval%vals(1, k1)
725  end do
726 
727  !Vegetation_Fraction
728  call ufo_geovals_get_var(geovals, var_sfc_vegfrac, geoval)
729  do k1 = 1, n_profiles
730  sfc(k1)%Vegetation_Fraction = geoval%vals(1, k1)
731  end do
732 
733  !Vegetation_Type
734  call ufo_geovals_get_var(geovals, var_sfc_vegtyp, geoval)
735  do k1 = 1, n_profiles
736  sfc(k1)%Vegetation_Type = int(geoval%vals(1, k1))
737  end do
738 
739  !Soil_Type
740  call ufo_geovals_get_var(geovals, var_sfc_soiltyp, geoval)
741  do k1 = 1, n_profiles
742  sfc(k1)%Soil_Type = int(geoval%vals(1, k1))
743  end do
744 
745  !Soil_Moisture_Content
746  call ufo_geovals_get_var(geovals, var_sfc_soilm, geoval)
747  do k1 = 1, n_profiles
748  sfc(k1)%Soil_Moisture_Content = geoval%vals(1, k1)
749  end do
750 
751  !Soil_Temperature
752  call ufo_geovals_get_var(geovals, var_sfc_soilt, geoval)
753  do k1 = 1, n_profiles
754  sfc(k1)%Soil_Temperature = geoval%vals(1, k1)
755  end do
756 
757  !Sea_Surface_Salinity
758  if (trim(conf%salinity_option) == "on") THEN
759  call ufo_geovals_get_var(geovals, var_sfc_sss, geoval)
760  do k1 = 1, n_profiles
761  sfc(k1)%Salinity = geoval%vals(1, k1)
762  end do
763  end if
764 
765 end subroutine load_sfc_data
766 
767 ! ------------------------------------------------------------------------------
768 
769 subroutine load_geom_data(obss,geo)
770 
771 implicit none
772 type(c_ptr), value, intent(in) :: obss
773 type(crtm_geometry_type), intent(inout) :: geo(:)
774 real(kind_real), allocatable :: tmpvar(:)
775 integer :: nlocs
776 
777  nlocs = obsspace_get_nlocs(obss)
778  allocate(tmpvar(nlocs))
779 
780  call obsspace_get_db(obss, "MetaData", "sensor_zenith_angle", tmpvar)
781  geo(:)%Sensor_Zenith_Angle = abs(tmpvar(:)) ! needs to be absolute value
782 
783  call obsspace_get_db(obss, "MetaData", "solar_zenith_angle", tmpvar)
784  geo(:)%Source_Zenith_Angle = tmpvar(:)
785 
786  call obsspace_get_db(obss, "MetaData", "sensor_azimuth_angle", tmpvar)
787  geo(:)%Sensor_Azimuth_Angle = tmpvar(:)
788 
789  call obsspace_get_db(obss, "MetaData", "solar_azimuth_angle", tmpvar)
790  geo(:)%Source_Azimuth_Angle = tmpvar(:)
791 
792 ! For some microwave instruments the solar and sensor azimuth angles can be
793 ! missing (given a value of 10^11). Set these to zero to get past CRTM QC.
794  where (geo(:)%Source_Azimuth_Angle < 0.0_kind_real .or. &
795  geo(:)%Source_Azimuth_Angle > 360.0_kind_real) &
796  geo(:)%Source_Azimuth_Angle = 0.0_kind_real
797  where (geo(:)%Sensor_Azimuth_Angle < 0.0_kind_real .or. &
798  geo(:)%Sensor_Azimuth_Angle > 360.0_kind_real) &
799  geo(:)%Sensor_Azimuth_Angle = 0.0_kind_real
800 
801  where (abs(geo(:)%Source_Zenith_Angle) > 180.0_kind_real) &
802  geo(:)%Source_Zenith_Angle = 100.0_kind_real
803 
804  call obsspace_get_db(obss, "MetaData", "scan_position", tmpvar)
805  geo(:)%Ifov = tmpvar(:)
806 
807  call obsspace_get_db(obss, "MetaData", "sensor_view_angle", tmpvar) !The Sensor_Scan_Angle is optional
808  geo(:)%Sensor_Scan_Angle = tmpvar(:)
809 
810  where (abs(geo(:)%Sensor_Scan_Angle) > 80.0_kind_real) &
811  geo(:)%Sensor_Scan_Angle = 0.0_kind_real
812 
813  deallocate(tmpvar)
814 
815 end subroutine load_geom_data
816 
817 ! ------------------------------------------------------------------------------
818 
819 subroutine get_var_name(n,varname)
820 
821 integer, intent(in) :: n
822 character(len=*), intent(out) :: varname
823 
824 character(len=6) :: chan
825 
826  write(chan, '(I0)') n
827  varname = 'brightness_temperature_' // trim(chan)
828 
829 end subroutine get_var_name
830 
831 ! -----------------------------------------------------------------------------
832 
833 
834 SUBROUTINE load_aerosol_data(n_profiles,n_layers,geovals,&
835  &aerosol_option,atm)
836 
837  USE crtm_aerosolcoeff, ONLY: aeroc
838 
839  INTEGER, INTENT(in) :: n_profiles,n_layers
840  TYPE(ufo_geovals), INTENT(in) :: geovals
841  TYPE(crtm_atmosphere_type), INTENT(inout) :: atm(:)
842 
843  CHARACTER(*) :: aerosol_option
844  CHARACTER(max_string) :: message
845  CHARACTER(len=MAXVARLEN) :: varname
846 
847  TYPE(ufo_geoval), POINTER :: geoval
848 
849  CHARACTER(*), PARAMETER :: routine_name = 'Load_Aerosol_Data'
850 
851  REAL(kind_real), DIMENSION(n_layers,n_profiles) :: rh
852  INTEGER :: ivar
853 
854  IF (trim(aerosol_option) == "aerosols_gocart_default") THEN
855  varname=var_rh
856  CALL ufo_geovals_get_var(geovals, varname, geoval)
857  rh(1:n_layers,1:n_profiles)=geoval%vals(1:n_layers,1:n_profiles)
858  WHERE (rh > 1_kind_real) rh=1_kind_real
860  ELSEIF (trim(aerosol_option) == "aerosols_gocart_merra_2") THEN
861  varname=var_rh
862  CALL ufo_geovals_get_var(geovals, varname, geoval)
863  rh(1:n_layers,1:n_profiles)=geoval%vals(1:n_layers,1:n_profiles)
864  WHERE (rh > 1_kind_real) rh=1_kind_real
866  ELSEIF (trim(aerosol_option) == "aerosols_other") THEN
867  CALL assign_other
868  ELSE
869  message = 'this aerosol not implemented - check later'
870  CALL display_message( aerosol_option, message, failure )
871  stop
872  ENDIF
873 
874  CONTAINS
875 
877 
878  INTEGER, PARAMETER :: ndust_bins=5, nseas_bins=4
879  REAL(kind_real), DIMENSION(ndust_bins), PARAMETER :: dust_radii=[&
880  &0.55_kind_real,1.4_kind_real,2.4_kind_real,4.5_kind_real,8.0_kind_real]
881 
882  INTEGER, DIMENSION(nseas_bins), PARAMETER :: seas_types=[&
883  seasalt_ssam_aerosol,seasalt_sscm1_aerosol,seasalt_sscm2_aerosol,seasalt_sscm3_aerosol]
884 
885  REAL(kind_real), DIMENSION(n_layers) :: ugkg_kgm2
886 
887  INTEGER :: i,k,m
888 
889  CHARACTER(len=MAXVARLEN) :: varname
890 
891  DO m=1,n_profiles
892 
893  CALL calculate_aero_layer_factor(atm(m),ugkg_kgm2)
894 
896 
897  varname=var_aerosols_gocart_default(i)
898  CALL ufo_geovals_get_var(geovals,varname, geoval)
899 
900  atm(m)%aerosol(i)%Concentration(1:n_layers)=&
901  &max(geoval%vals(:,m)*ugkg_kgm2,aerosol_concentration_minvalue_layer)
902 
903  SELECT CASE ( trim(varname))
904  CASE ('sulf')
905  atm(m)%aerosol(i)%type = sulfate_aerosol
906  DO k=1,n_layers
907  atm(m)%aerosol(i)%effective_radius(k)=&
908  &gocart_aerosol_size(atm(m)%aerosol(i)%type, &
909  &rh(k,m))
910  ENDDO
911 
912  CASE ('bc1')
913  atm(m)%aerosol(i)%type = black_carbon_aerosol
914  atm(m)%aerosol(i)%effective_radius(:)=&
915  &aeroc%Reff(1,atm(m)%aerosol(i)%type)
916  CASE ('bc2')
917  atm(m)%aerosol(i)%type = black_carbon_aerosol
918  DO k=1,n_layers
919  atm(m)%aerosol(i)%effective_radius(k)=&
920  &gocart_aerosol_size(atm(m)%aerosol(i)%type, &
921  &rh(k,m))
922  ENDDO
923 
924  CASE ('oc1')
925  atm(m)%aerosol(i)%type = organic_carbon_aerosol
926  atm(m)%aerosol(i)%effective_radius(:)=&
927  &aeroc%Reff(1,atm(m)%aerosol(i)%type)
928  CASE ('oc2')
929  atm(m)%aerosol(i)%type = organic_carbon_aerosol
930  DO k=1,n_layers
931  atm(m)%aerosol(i)%effective_radius(k)=&
932  &gocart_aerosol_size(atm(m)%aerosol(i)%type, &
933  &rh(k,m))
934  ENDDO
935 
936  CASE ('dust1')
937  atm(m)%aerosol(i)%type = dust_aerosol
938  atm(m)%aerosol(i)%effective_radius(:)=dust_radii(1)
939  CASE ('dust2')
940  atm(m)%aerosol(i)%type = dust_aerosol
941  atm(m)%aerosol(i)%effective_radius(:)=dust_radii(2)
942  CASE ('dust3')
943  atm(m)%aerosol(i)%type = dust_aerosol
944  atm(m)%aerosol(i)%effective_radius(:)=dust_radii(3)
945  CASE ('dust4')
946  atm(m)%aerosol(i)%type = dust_aerosol
947  atm(m)%aerosol(i)%effective_radius(:)=dust_radii(4)
948  CASE ('dust5')
949  atm(m)%aerosol(i)%type = dust_aerosol
950  atm(m)%aerosol(i)%effective_radius(:)=dust_radii(5)
951 
952  CASE ('seas1')
953  atm(m)%aerosol(i)%type = seas_types(1)
954  DO k=1,n_layers
955  atm(m)%aerosol(i)%effective_radius(k)=&
956  &gocart_aerosol_size(atm(m)%aerosol(i)%type, &
957  &rh(k,m))
958  ENDDO
959  CASE ('seas2')
960  atm(m)%aerosol(i)%type = seas_types(2)
961  DO k=1,n_layers
962  atm(m)%aerosol(i)%effective_radius(k)=&
963  &gocart_aerosol_size(atm(m)%aerosol(i)%type, &
964  &rh(k,m))
965  ENDDO
966  CASE ('seas3')
967  atm(m)%aerosol(i)%type = seas_types(3)
968  DO k=1,n_layers
969  atm(m)%aerosol(i)%effective_radius(k)=&
970  &gocart_aerosol_size(atm(m)%aerosol(i)%type, &
971  &rh(k,m))
972  ENDDO
973  CASE ('seas4')
974  atm(m)%aerosol(i)%type = seas_types(4)
975  DO k=1,n_layers
976  atm(m)%aerosol(i)%effective_radius(k)=&
977  &gocart_aerosol_size(atm(m)%aerosol(i)%type, &
978  &rh(k,m))
979  ENDDO
980  END SELECT
981 
982  ENDDO
983 
984  ENDDO
985 
986 
987  END SUBROUTINE assign_gocart_default
988 
990 
991  message = 'this aerosol not implemented in the CRTM - check later'
992  CALL display_message( aerosol_option, message, failure )
993  stop
994 
995  END SUBROUTINE assign_gocart_merra_2
996 
997  SUBROUTINE assign_other
998 
999  message = 'this aerosol not implemented - check later'
1000  CALL display_message( aerosol_option, message, failure )
1001  stop
1002 
1003  END SUBROUTINE assign_other
1004 
1005  END SUBROUTINE load_aerosol_data
1006 
1007  SUBROUTINE assign_aerosol_names(aerosol_option,var_aerosols)
1008 
1009  CHARACTER(*), INTENT(in) :: aerosol_option
1010  CHARACTER(len=MAXVARLEN), ALLOCATABLE, INTENT(out) :: var_aerosols(:)
1011 
1012  CHARACTER(max_string) :: err_msg
1013 
1014  IF (aerosol_option == "aerosols_gocart_default") THEN
1015  ALLOCATE(var_aerosols(SIZE(var_aerosols_gocart_default)))
1016  var_aerosols=var_aerosols_gocart_default
1017  ELSEIF (aerosol_option == "aerosols_gocart_merra_2") THEN
1018  ALLOCATE(var_aerosols(SIZE(var_aerosols_gocart_merra_2)))
1019  var_aerosols=var_aerosols_gocart_merra_2
1020  ELSEIF (aerosol_option == "var_aerosols_other") THEN
1021  ALLOCATE(var_aerosols(SIZE(var_aerosols_other)))
1022  var_aerosols=var_aerosols_other
1023  ELSE
1024  WRITE(err_msg,*) 'assign_aerosol_names: aerosol_option not implemented '//trim(aerosol_option)
1025  call abor1_ftn(err_msg)
1026  END IF
1027 
1028  END SUBROUTINE assign_aerosol_names
1029 
1031 
1032  TYPE(crtm_atmosphere_type), INTENT(in) :: atm
1033  REAL(kind_real), INTENT(out) :: ugkg_kgm2(:)
1034 
1035  INTEGER :: k
1036 
1037 !rh, ug2kg need to be from top to bottom
1038 !ug2kg && hPa2Pa
1039  DO k=1,SIZE(ugkg_kgm2)
1040 !correct for mixing ratio factor ugkg_kgm2
1041 !being calculated from dry pressure, cotton eq. (2.4)
1042 !p_dry=p_total/(1+1.61*mixing_ratio)
1043  ugkg_kgm2(k)=1.0e-9_kind_real*(atm%Level_Pressure(k)-&
1044  &atm%Level_Pressure(k-1))*100_kind_real/grav/&
1045  &(1_kind_real+rv_rd*atm%Absorber(k,1)*1e-3_kind_real)
1046  ENDDO
1047 
1049 
1050  SUBROUTINE calculate_aero_layer_factor_atm(atm,ugkg_kgm2)
1051 
1052  TYPE(crtm_atmosphere_type), INTENT(in) :: atm(:)
1053  REAL(kind_real), INTENT(out) :: ugkg_kgm2(:,:)
1054 
1055  INTEGER :: k,m
1056 
1057 !rh, ug2kg need to be from top to bottom
1058 !ug2kg && hPa2Pa
1059  DO k=1,SIZE(ugkg_kgm2,1)
1060  DO m=1,SIZE(ugkg_kgm2,2)
1061 !correct for mixing ratio factor ugkg_kgm2
1062 !being calculated from dry pressure, cotton eq. (2.4)
1063 !p_dry=p_total/(1+1.61*mixing_ratio)
1064  ugkg_kgm2(k,m)=1.0e-9_kind_real*(atm(m)%Level_Pressure(k)-&
1065  &atm(m)%Level_Pressure(k-1))*100_kind_real/grav/&
1066  &(1_kind_real+rv_rd*atm(m)%Absorber(k,1)*1.e-3_kind_real)
1067  ENDDO
1068  ENDDO
1069 
1070  END SUBROUTINE calculate_aero_layer_factor_atm
1071 
1072  FUNCTION gocart_aerosol_size( itype, rh ) & ! rh input in 0-1
1073  &result(r_eff) ! in micrometer
1074 
1075  USE crtm_aerosolcoeff, ONLY: aeroc
1076  IMPLICIT NONE
1077 
1078 !
1079 ! modified from a function provided by quanhua liu
1080 !
1081  INTEGER ,INTENT(in) :: itype
1082  REAL(kind_real) ,INTENT(in) :: rh
1083 
1084  INTEGER :: j1,j2,m
1085  REAL(kind_real) :: h1
1086  REAL(kind_real) :: r_eff
1087 
1088  j2 = 0
1089  j1 = 1
1090  IF ( rh <= aeroc%rh(1) ) THEN
1091  j1 = 1
1092  ELSE IF ( rh >= aeroc%rh(aeroc%n_rh) ) THEN
1093  j1 = aeroc%n_rh
1094  ELSE
1095  DO m = 1, aeroc%n_rh-1
1096  IF ( rh < aeroc%rh(m+1) .AND. rh > aeroc%rh(m) ) THEN
1097  j1 = m
1098  j2 = m+1
1099  h1 = (rh-aeroc%rh(m))/(aeroc%rh(m+1)-aeroc%rh(m))
1100  EXIT
1101  ENDIF
1102  ENDDO
1103  ENDIF
1104 
1105  IF ( j2 == 0 ) THEN
1106  r_eff = aeroc%reff(j1,itype )
1107  ELSE
1108  r_eff = (1_kind_real-h1)*aeroc%reff(j1,itype ) + h1*aeroc%reff(j2,itype )
1109  ENDIF
1110 
1111  END FUNCTION gocart_aerosol_size
1112 
1113 
1114  FUNCTION upper2lower(str) RESULT(string)
1115 
1116  IMPLICIT NONE
1117 
1118  CHARACTER(*), INTENT(in) :: str
1119  CHARACTER(LEN(str)) :: string
1120 
1121  INTEGER :: ic, i
1122 
1123  CHARACTER(26), PARAMETER :: upper = 'ABCDEFGHIJKLMNOPQRSTUVWXYZ'
1124  CHARACTER(26), PARAMETER :: lower = 'abcdefghijklmnopqrstuvwxyz'
1125 
1126 ! lowcase each letter if it is lowecase
1127  string = str
1128  DO i = 1, len_trim(str)
1129  ic = index(upper, str(i:i))
1130  IF (ic > 0) string(i:i) = lower(ic:ic)
1131  END DO
1132 
1133  END FUNCTION upper2lower
1134 
1135  INTEGER FUNCTION getindex(names,usrname)
1136  IMPLICIT NONE
1137  CHARACTER(len=*),INTENT(in) :: names(:)
1138  CHARACTER(len=*),INTENT(in) :: usrname
1139  INTEGER i
1140  getindex=-1
1141  DO i=1,SIZE(names)
1142  IF(trim(usrname)==trim(names(i))) THEN
1143  getindex=i
1144  EXIT
1145  ENDIF
1146  ENDDO
1147  END FUNCTION getindex
1148 
1149 !from fv3
1150 
1151  SUBROUTINE qsmith_atm(atm,rh)
1152 
1153  TYPE(crtm_atmosphere_type), INTENT(in) :: atm(:)
1154  REAL(kind_real), INTENT(out),DIMENSION(:,:):: rh
1155 
1156  REAL, ALLOCATABLE :: table(:),des(:)
1157 
1158  REAL es, qs, q
1159  REAL ap1, eps10
1160  REAL Tmin
1161  INTEGER i, k, it, n_layers, n_profiles
1162 
1163  n_layers=SIZE(rh,1)
1164  n_profiles=SIZE(rh,2)
1165 
1166  tmin = tice-160.
1167  eps10 = 10.*esl
1168 
1169  IF( .NOT. ALLOCATED(table) ) CALL qsmith_init(table,des)
1170 
1171  DO k=1,n_layers
1172  DO i=1,n_profiles
1173  ap1 = 10.*dim(atm(i)%Temperature(k), tmin) + 1.
1174  ap1 = min(2621., ap1)
1175  it = ap1
1176  es = table(it) + (ap1-it)*des(it)
1177  q=atm(i)%Absorber(k,1)*1.e-3/(1.+atm(i)%Absorber(k,1)*1.e-3)
1178  qs = esl*es*(1.+zvir*q)/(atm(i)%Pressure(k)*100.)
1179  rh(k,i) = q/qs
1180  ENDDO
1181  ENDDO
1182 
1183  END SUBROUTINE qsmith_atm
1184 
1185  SUBROUTINE qsmith_profiles(t,sphum,p,rh)
1186 
1187  REAL(kind_real), DIMENSION(:,:), INTENT(in) :: t,sphum,p
1188  REAL(kind_real), DIMENSION(:,:), INTENT(out) :: rh
1189 
1190  REAL, ALLOCATABLE :: table(:),des(:)
1191 
1192  REAL es, qs, q
1193  REAL ap1, eps10
1194  REAL Tmin
1195  INTEGER i, k, it, n_layers, n_profiles
1196 
1197  n_layers=SIZE(t,1)
1198  n_profiles=SIZE(t,2)
1199 
1200  tmin = tice-160.
1201  eps10 = 10.*esl
1202 
1203  IF ( .NOT. ALLOCATED(table) ) CALL qsmith_init(table,des)
1204 
1205  DO k=1,n_layers
1206  DO i=1,n_profiles
1207  ap1 = 10.*dim(t(k,i), tmin) + 1.
1208  ap1 = min(2621., ap1)
1209  it = ap1
1210  es = table(it) + (ap1-it)*des(it)
1211  q=sphum(k,i)
1212  qs = esl*es*(1.+zvir*q)/p(k,i)
1213  rh(k,i) = q/qs
1214  ENDDO
1215  ENDDO
1216 
1217  END SUBROUTINE qsmith_profiles
1218 
1219  SUBROUTINE qsmith_init(table,des)
1220 
1221  REAL, ALLOCATABLE, INTENT(out) :: table(:),des(:)
1222  INTEGER, PARAMETER:: length=2621
1223  INTEGER i
1224 
1225  IF( .NOT. ALLOCATED(table) ) THEN
1226 ! Generate es table (dT = 0.1 deg. C)
1227 
1228  ALLOCATE ( table(length) )
1229  ALLOCATE ( des(length) )
1230 
1231  CALL qs_table(length, table)
1232 
1233  DO i=1,length-1
1234  des(i) = table(i+1) - table(i)
1235  ENDDO
1236  des(length) = des(length-1)
1237  ENDIF
1238 
1239  END SUBROUTINE qsmith_init
1240 
1241  SUBROUTINE qs_table(n,table)
1242  INTEGER, INTENT(in):: n
1243  REAL table (n)
1244  REAL :: dt=0.1
1245  REAL esbasw, tbasw, tbasi, Tmin, tem, aa, b, c, d, e
1246  INTEGER i
1247 ! Constants
1248  esbasw = 1013246.0
1249  tbasw = 373.16
1250  tbasi = 273.16
1251  tmin = tbasi - 160.
1252 ! Compute es over water
1253 ! see smithsonian meteorological tables page 350.
1254  DO i=1,n
1255  tem = tmin+dt*real(i-1)
1256  aa = -7.90298*(tbasw/tem-1)
1257  b = 5.02808*alog10(tbasw/tem)
1258  c = -1.3816e-07*(10**((1-tem/tbasw)*11.344)-1)
1259  d = 8.1328e-03*(10**((tbasw/tem-1)*(-3.49149))-1)
1260  e = alog10(esbasw)
1261  table(i) = 0.1*10**(aa+b+c+d+e)
1262  ENDDO
1263 
1264  END SUBROUTINE qs_table
1265 
1266 END MODULE ufo_crtm_utils_mod
1267 
ufo_vars_mod::var_aerosols_gocart_default
character(len=maxvarlen), dimension(n_aerosols_gocart_default), parameter, public var_aerosols_gocart_default
Definition: ufo_variables_mod.F90:95
ufo_crtm_utils_mod::load_atm_data
subroutine, public load_atm_data(n_Profiles, n_Layers, geovals, atm, conf)
Definition: ufo_crtm_utils_mod.F90:497
ufo_avgkernel_mod::max_string
integer, parameter max_string
Definition: ufo_avgkernel_mod.F90:17
ufo_vars_mod::var_aerosols_gocart_merra_2
character(len=maxvarlen), dimension(n_aerosols_gocart_merra_2), parameter, public var_aerosols_gocart_merra_2
Definition: ufo_variables_mod.F90:115
ufo_vars_mod::var_cldfrac
character(len=maxvarlen), parameter, public var_cldfrac
Definition: ufo_variables_mod.F90:46
ufo_crtm_utils_mod::qsmith_atm
subroutine qsmith_atm(atm, rh)
Definition: ufo_crtm_utils_mod.F90:1152
ufo_vars_mod::var_sfc_vegtyp
character(len=maxvarlen), parameter, public var_sfc_vegtyp
Definition: ufo_variables_mod.F90:67
ufo_crtm_utils_mod::crtm_absorber_units
integer, dimension(3), parameter crtm_absorber_units
Definition: ufo_crtm_utils_mod.F90:119
ufo_crtm_utils_mod::crtm_surfaces
character(len=maxvarlen), dimension(4), parameter crtm_surfaces
Definition: ufo_crtm_utils_mod.F90:152
ufo_vars_mod::var_sfc_vegfrac
character(len=maxvarlen), parameter, public var_sfc_vegfrac
Definition: ufo_variables_mod.F90:60
ufo_crtm_utils_mod::aerosol_concentration_minvalue_layer
real(kind_real), parameter, public aerosol_concentration_minvalue_layer
Definition: ufo_crtm_utils_mod.F90:42
ufo_vars_mod::var_sfc_soilm
character(len=maxvarlen), parameter, public var_sfc_soilm
Definition: ufo_variables_mod.F90:64
ufo_vars_mod::var_cliefr
character(len=maxvarlen), parameter, public var_cliefr
Definition: ufo_variables_mod.F90:41
ufo_crtm_utils_mod::assign_aerosol_names
subroutine, public assign_aerosol_names(aerosol_option, var_aerosols)
Definition: ufo_crtm_utils_mod.F90:1008
ufo_vars_mod::var_cls
character(len=maxvarlen), parameter, public var_cls
Definition: ufo_variables_mod.F90:37
ufo_vars_mod::var_sfc_wdir
character(len=maxvarlen), parameter, public var_sfc_wdir
Definition: ufo_variables_mod.F90:62
ufo_crtm_utils_mod::crtm_comm_stat_check
subroutine, public crtm_comm_stat_check(stat, PROGRAM_NAME, message, f_comm)
Definition: ufo_crtm_utils_mod.F90:428
ufo_vars_mod::var_clg
character(len=maxvarlen), parameter, public var_clg
Definition: ufo_variables_mod.F90:38
ufo_crtm_utils_mod::load_geom_data
subroutine, public load_geom_data(obss, geo)
Definition: ufo_crtm_utils_mod.F90:770
ufo_gnssro_bndnbam_tlad_mod::missing
real(c_double) missing
Definition: ufo_gnssro_bndnbam_tlad_mod.F90:23
ufo_vars_mod::var_sfc_ifrac
character(len=maxvarlen), parameter, public var_sfc_ifrac
Definition: ufo_variables_mod.F90:53
ufo_vars_mod::var_rh
character(len=maxvarlen), parameter, public var_rh
Definition: ufo_variables_mod.F90:90
ufo_crtm_utils_mod::load_aerosol_data
subroutine, public load_aerosol_data(n_profiles, n_layers, geovals, aerosol_option, atm)
Definition: ufo_crtm_utils_mod.F90:836
assign_other
subroutine assign_other
Definition: ufo_crtm_utils_mod.F90:998
ufo_crtm_utils_mod::crtm_clouds
character(len= *), dimension(n_valid_cloud_categories), parameter crtm_clouds
Definition: ufo_crtm_utils_mod.F90:134
ufo_crtm_utils_mod::crtm_absorbers
character(len= *), dimension(n_valid_absorber_ids), parameter crtm_absorbers
Definition: ufo_crtm_utils_mod.F90:106
assign_gocart_default
subroutine assign_gocart_default
Definition: ufo_crtm_utils_mod.F90:877
ufo_vars_mod::var_clsefr
character(len=maxvarlen), parameter, public var_clsefr
Definition: ufo_variables_mod.F90:43
ufo_crtm_utils_mod::gocart_aerosol_size
real(kind_real) function gocart_aerosol_size(itype, rh)
Definition: ufo_crtm_utils_mod.F90:1074
ufo_vars_mod::var_clgefr
character(len=maxvarlen), parameter, public var_clgefr
Definition: ufo_variables_mod.F90:44
ufo_crtm_utils_mod::qsmith_profiles
subroutine qsmith_profiles(t, sphum, p, rh)
Definition: ufo_crtm_utils_mod.F90:1186
ufo_vars_mod::var_sfc_sss
character(len=maxvarlen), parameter, public var_sfc_sss
Definition: ufo_variables_mod.F90:75
ufo_vars_mod::var_sfc_soilt
character(len=maxvarlen), parameter, public var_sfc_soilt
Definition: ufo_variables_mod.F90:65
ufo_vars_mod::var_sfc_wspeed
character(len=maxvarlen), parameter, public var_sfc_wspeed
Definition: ufo_variables_mod.F90:61
ufo_crtm_utils_mod::ufo_absorbers
character(len=maxvarlen), dimension(3), parameter ufo_absorbers
Definition: ufo_crtm_utils_mod.F90:101
ufo_crtm_utils_mod::rv_rd
real(kind_real), parameter, public rv_rd
Definition: ufo_crtm_utils_mod.F90:42
ufo_vars_mod::var_cli
character(len=maxvarlen), parameter, public var_cli
Definition: ufo_variables_mod.F90:35
ufo_crtm_utils_mod::load_sfc_data
subroutine, public load_sfc_data(n_Profiles, n_Channels, channels, geovals, sfc, chinfo, obss, conf)
Definition: ufo_crtm_utils_mod.F90:589
ufo_vars_mod::var_sfc_wfrac
character(len=maxvarlen), parameter, public var_sfc_wfrac
Definition: ufo_variables_mod.F90:51
ufo_crtm_utils_mod::getindex
integer function getindex(names, usrname)
Definition: ufo_crtm_utils_mod.F90:1136
ufo_vars_mod::ufo_vars_getindex
integer function, public ufo_vars_getindex(vars, varname)
Definition: ufo_variables_mod.F90:204
ufo_vars_mod::var_mixr
character(len=maxvarlen), parameter, public var_mixr
Definition: ufo_variables_mod.F90:21
ufo_crtm_utils_mod
Fortran module to provide code shared between nonlinear and tlm/adm radiance calculations.
Definition: ufo_crtm_utils_mod.F90:8
ufo_crtm_utils_mod::grav
real(kind_real), parameter, public grav
Definition: ufo_crtm_utils_mod.F90:42
ufo_crtm_utils_mod::calculate_aero_layer_factor_atm_profile
subroutine calculate_aero_layer_factor_atm_profile(atm, ugkg_kgm2)
Definition: ufo_crtm_utils_mod.F90:1031
ufo_vars_mod::n_aerosols_gocart_default
integer, parameter, public n_aerosols_gocart_default
Definition: ufo_variables_mod.F90:14
ufo_geovals_mod
Definition: ufo_geovals_mod.F90:7
ufo_vars_mod::var_clhefr
character(len=maxvarlen), parameter, public var_clhefr
Definition: ufo_variables_mod.F90:45
ufo_crtm_utils_mod::tice
real(kind_real), parameter tice
Definition: ufo_crtm_utils_mod.F90:42
ufo_vars_mod::var_prsi
character(len=maxvarlen), parameter, public var_prsi
Definition: ufo_variables_mod.F90:26
ufo_vars_mod::var_sfc_itmp
character(len=maxvarlen), parameter, public var_sfc_itmp
Definition: ufo_variables_mod.F90:57
ufo_crtm_utils_mod::crtm_conf_delete
subroutine, public crtm_conf_delete(conf)
Definition: ufo_crtm_utils_mod.F90:411
ufo_vars_mod::var_clwefr
character(len=maxvarlen), parameter, public var_clwefr
Definition: ufo_variables_mod.F90:40
ufo_crtm_utils_mod::crtm_absorber_id
integer, dimension(n_valid_absorber_ids), parameter crtm_absorber_id
Definition: ufo_crtm_utils_mod.F90:109
ufo_vars_mod::var_sfc_wtmp
character(len=maxvarlen), parameter, public var_sfc_wtmp
Definition: ufo_variables_mod.F90:55
ufo_vars_mod::var_sfc_soiltyp
character(len=maxvarlen), parameter, public var_sfc_soiltyp
Definition: ufo_variables_mod.F90:68
ufo_crtm_utils_mod::ufo_surfaces
character(len=maxvarlen), dimension(4), parameter ufo_surfaces
Definition: ufo_crtm_utils_mod.F90:148
ufo_vars_mod::var_sfc_sfrac
character(len=maxvarlen), parameter, public var_sfc_sfrac
Definition: ufo_variables_mod.F90:54
ufo_crtm_utils_mod::ufo_crtm_skip_profiles
subroutine, public ufo_crtm_skip_profiles(n_Profiles, n_Channels, channels, obss, Skip_Profiles)
Definition: ufo_crtm_utils_mod.F90:450
ufo_vars_mod::var_clrefr
character(len=maxvarlen), parameter, public var_clrefr
Definition: ufo_variables_mod.F90:42
ufo_vars_mod::var_sfc_lfrac
character(len=maxvarlen), parameter, public var_sfc_lfrac
Definition: ufo_variables_mod.F90:52
ufo_crtm_utils_mod::esl
real(kind_real), parameter esl
Definition: ufo_crtm_utils_mod.F90:42
ufo_crtm_utils_mod::qs_table
subroutine qs_table(n, table)
Definition: ufo_crtm_utils_mod.F90:1242
ufo_crtm_utils_mod::qsmith_init
subroutine qsmith_init(table, des)
Definition: ufo_crtm_utils_mod.F90:1220
ufo_crtm_utils_mod::calculate_aero_layer_factor
Definition: ufo_crtm_utils_mod.F90:84
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_crtm_utils_mod::crtm_conf
Definition: ufo_crtm_utils_mod.F90:57
ufo_crtm_utils_mod::calculate_aero_layer_factor_atm
subroutine calculate_aero_layer_factor_atm(atm, ugkg_kgm2)
Definition: ufo_crtm_utils_mod.F90:1051
ufo_crtm_utils_mod::ozone_default_value
real(kind_real), parameter ozone_default_value
Definition: ufo_crtm_utils_mod.F90:42
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_vars_mod::var_clw
character(len=maxvarlen), parameter, public var_clw
Definition: ufo_variables_mod.F90:34
ufo_vars_mod::var_ts
character(len=maxvarlen), parameter, public var_ts
Definition: ufo_variables_mod.F90:19
ufo_crtm_utils_mod::qsmith
Definition: ufo_crtm_utils_mod.F90:91
ufo_crtm_utils_mod::upper2lower
character(len(str)) function, public upper2lower(str)
Definition: ufo_crtm_utils_mod.F90:1115
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_vars_mod::var_clh
character(len=maxvarlen), parameter, public var_clh
Definition: ufo_variables_mod.F90:39
ufo_crtm_utils_mod::crtm_conf_setup
subroutine, public crtm_conf_setup(conf, f_confOpts, f_confOper)
Definition: ufo_crtm_utils_mod.F90:161
ufo_vars_mod::var_clr
character(len=maxvarlen), parameter, public var_clr
Definition: ufo_variables_mod.F90:36
ufo_vars_mod::maxvarlen
integer, parameter, public maxvarlen
Definition: ufo_variables_mod.F90:17
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_geovals_mod::ufo_geovals
type to hold interpolated fields required by the obs operators
Definition: ufo_geovals_mod.F90:47
ufo_crtm_utils_mod::crtm_cloud_id
integer, dimension(n_valid_cloud_categories), parameter crtm_cloud_id
Definition: ufo_crtm_utils_mod.F90:137
ufo_crtm_utils_mod::aerosol_concentration_minvalue
real(kind_real), parameter, public aerosol_concentration_minvalue
Definition: ufo_crtm_utils_mod.F90:42
ufo_vars_mod::var_aerosols_other
character(len=maxvarlen), dimension(n_aerosols_other), parameter, public var_aerosols_other
Definition: ufo_variables_mod.F90:133
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_crtm_utils_mod::ufo_clouds
character(len=maxvarlen), dimension(n_valid_cloud_categories, 2), parameter ufo_clouds
Definition: ufo_crtm_utils_mod.F90:126
ufo_crtm_utils_mod::zvir
real(kind_real), parameter zvir
Definition: ufo_crtm_utils_mod.F90:42
ufo_vars_mod::var_prs
character(len=maxvarlen), parameter, public var_prs
Definition: ufo_variables_mod.F90:25
ufo_crtm_utils_mod::get_var_name
subroutine, public get_var_name(n, varname)
Definition: ufo_crtm_utils_mod.F90:820
ufo_vars_mod::var_sfc_ltmp
character(len=maxvarlen), parameter, public var_sfc_ltmp
Definition: ufo_variables_mod.F90:56
assign_gocart_merra_2
subroutine assign_gocart_merra_2
Definition: ufo_crtm_utils_mod.F90:990