UFO
ufo_aodluts_mod.F90
Go to the documentation of this file.
1 ! (c) copyright 2017-2018 ucar
2 !
3 ! this software is licensed under the terms of the apache licence version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/license-2.0.
5 
6 !> fortran module to handle aod observations
7 
9 
10  USE fckit_configuration_module, ONLY: fckit_configuration
11  USE iso_c_binding
12  USE kinds
13  USE missing_values_mod
14 
16  USE ufo_vars_mod
20  USE crtm_module
21  USE crtm_spccoeff, ONLY: sc
22  USE obsspace_mod
23 
24  USE cf_mieobs_mod, ONLY: get_cf_aod
25 
26  IMPLICIT NONE
27  PRIVATE
28 
29 !> fortran derived type for aod trajectory
30  TYPE, PUBLIC :: ufo_aodluts
31  PRIVATE
32  CHARACTER(len=maxvarlen), PUBLIC, ALLOCATABLE :: varin(:) ! variablesrequested from the model
33  INTEGER, ALLOCATABLE :: channels(:)
34  REAL(kind_real), ALLOCATABLE :: wavelengths(:)
35  INTEGER :: n_aerosols
36  TYPE(luts_conf) :: conf
37  CONTAINS
38  PROCEDURE :: setup => ufo_aodluts_setup
39  PROCEDURE :: delete => ufo_aodluts_delete
40  PROCEDURE :: simobs => ufo_aodluts_simobs
41  END TYPE ufo_aodluts
42  CHARACTER(len=maxvarlen), DIMENSION(4), PARAMETER :: varin_default = &
44 
45  CHARACTER(maxvarlen), PARAMETER :: varname_tmplate="aerosol_optical_depth"
46 
47 CONTAINS
48 
49 ! ------------------------------------------------------------------------------
50 
51  SUBROUTINE ufo_aodluts_setup(self, f_confoper, channels)
52 
53  IMPLICIT NONE
54  CLASS(ufo_aodluts), INTENT(inout) :: self
55  TYPE(fckit_configuration), INTENT(in) :: f_confoper
56  INTEGER(c_int), INTENT(in) :: channels(:) !list of channels to use
57 
58  INTEGER :: nvars_in, rc
59  CHARACTER(len=max_string) :: err_msg
60  TYPE(fckit_configuration) :: f_confopts
61 
62  CHARACTER(len=maxvarlen), ALLOCATABLE :: var_aerosols(:)
63 
64  CALL f_confoper%get_or_die("obs options",f_confopts)
65 
66  CALL luts_conf_setup(self%conf, f_confopts, f_confoper)
67 
68  CALL assign_aerosol_names(self%conf%aerosol_option,var_aerosols)
69 
70  self%n_aerosols=SIZE(var_aerosols)
71  nvars_in = SIZE(varin_default)+self%n_aerosols
72 
73  ALLOCATE(self%varin(nvars_in))
74  self%varin(1:SIZE(varin_default)) = varin_default
75  self%varin(SIZE(varin_default)+1:) = var_aerosols
76 
77  ALLOCATE(self%channels(SIZE(channels)))
78  ALLOCATE(self%wavelengths(SIZE(channels)))
79 
80  self%channels(:) = channels(:)
81 
82  DEALLOCATE(var_aerosols)
83 
84  END SUBROUTINE ufo_aodluts_setup
85 
86 ! ------------------------------------------------------------------------------
87 
88  SUBROUTINE ufo_aodluts_delete(self)
89 
90  IMPLICIT NONE
91  CLASS(ufo_aodluts), INTENT(inout) :: self
92 
93  CALL luts_conf_delete(self%conf)
94 
95  IF (ALLOCATED(self%varin)) DEALLOCATE(self%varin)
96  IF (ALLOCATED(self%channels)) DEALLOCATE(self%channels)
97  IF (ALLOCATED(self%wavelengths)) DEALLOCATE(self%wavelengths)
98 
99  END SUBROUTINE ufo_aodluts_delete
100 
101 ! ------------------------------------------------------------------------------
102 
103  SUBROUTINE ufo_aodluts_simobs(self, geovals, obss, nvars, nlocs, hofx)
104 
105  IMPLICIT NONE
106  CLASS(ufo_aodluts), INTENT(inout) :: self
107  TYPE(ufo_geovals), INTENT(in) :: geovals
108  INTEGER, INTENT(in) :: nvars, nlocs
109  REAL(c_double), INTENT(inout) :: hofx(nvars, nlocs)
110  TYPE(c_ptr), VALUE, INTENT(in) :: obss
111 
112 ! local variables
113  CHARACTER(*), PARAMETER :: program_name = 'ufo_aodluts_mod.f90'
114  CHARACTER(255) :: message, version
115  INTEGER :: err_stat, alloc_stat
116  INTEGER :: l, m, n, i
117  TYPE(ufo_geoval), POINTER :: temp
118  REAL(c_double) :: missing
119 
120  INTEGER :: n_profiles
121  INTEGER :: n_layers
122  INTEGER :: n_channels
123  INTEGER :: n_aerosols
124 
125 ! define the "non-demoninational" arguments
126  TYPE(crtm_channelinfo_type) :: chinfo(self%conf%n_sensors)
127 
128  REAL(kind_real), ALLOCATABLE :: wavelengths_all(:)
129  REAL(kind_real), ALLOCATABLE :: aero_layers(:,:,:),rh(:,:)
130 
131  CHARACTER(len=maxvarlen), ALLOCATABLE :: var_aerosols(:)
132 
133  INTEGER :: rc
134 
135  CALL assign_aerosol_names(self%conf%aerosol_option,var_aerosols)
136 
137  n_profiles = geovals%nlocs
138  CALL ufo_geovals_get_var(geovals, var_aerosols(1), temp)
139  n_layers = temp%nval
140  NULLIFY(temp)
141 
142  n_aerosols=self%n_aerosols
143 
144  ALLOCATE(aero_layers(n_aerosols,n_layers,n_profiles),&
145  &rh(n_layers,n_profiles))
146 
147  err_stat = crtm_init( self%conf%sensor_id, &
148  chinfo, &
149  file_path=trim(self%conf%coefficient_path), &
150  quiet=.true.)
151 
152  IF ( err_stat /= success ) THEN
153  message = 'error initializing crtm'
154  CALL display_message( program_name, message, failure )
155  stop
156  END IF
157 
158  sensor_loop:DO n = 1, self%conf%n_sensors
159 
160  n_channels = crtm_channelinfo_n_channels(chinfo(n))
161 
162  IF (ALLOCATED(wavelengths_all)) DEALLOCATE(wavelengths_all)
163 
164  ALLOCATE(wavelengths_all(n_channels), stat = alloc_stat)
165 
166  IF ( alloc_stat /= 0 ) THEN
167  message = 'error allocating wavelengths_all'
168  CALL display_message( program_name, message, failure )
169  stop
170  END IF
171 
172  wavelengths_all=1.e7/sc(chinfo(n)%sensor_index)%wavenumber(:)
173 
174  self%wavelengths=wavelengths_all(self%channels)
175 
176  CALL calculate_aero_layers(self%conf%aerosol_option,&
177  &n_aerosols, n_profiles, n_layers,&
178  &geovals, aero_layers=aero_layers, rh=rh)
179 
180  CALL get_cf_aod(n_layers, n_profiles, nvars, n_aerosols, &
181  &self%conf%rcfile, &
182  &self%wavelengths, var_aerosols, aero_layers, rh, &
183  &aod_tot = hofx, rc = rc)
184 
185 
186  DEALLOCATE(aero_layers,rh,wavelengths_all)
187 
188  IF (rc /= 0) THEN
189  message = 'error on exit from get_cf_aod'
190  CALL display_message( program_name, message, failure )
191  stop
192  END IF
193 
194  END DO sensor_loop
195 
196  err_stat = crtm_destroy( chinfo )
197  IF ( err_stat /= success ) THEN
198  message = 'error destroying crtm (settraj)'
199  CALL display_message( program_name, message, failure )
200  stop
201  END IF
202 
203  END SUBROUTINE ufo_aodluts_simobs
204 
205 ! ------------------------------------------------------------------------------
206 
207 END MODULE ufo_aodluts_mod
fortran module to handle aod observations
character(maxvarlen), parameter varname_tmplate
subroutine ufo_aodluts_simobs(self, geovals, obss, nvars, nlocs, hofx)
subroutine ufo_aodluts_delete(self)
character(len=maxvarlen), dimension(4), parameter varin_default
subroutine ufo_aodluts_setup(self, f_confoper, channels)
Fortran module to provide code shared between nonlinear and tlm/adm radiance calculations.
integer, parameter, public max_string
subroutine, public assign_aerosol_names(aerosol_option, var_aerosols)
subroutine, public ufo_geovals_get_var(self, varname, geoval)
fortran module to provide code shared between nonlinear and tlm/adm radiance calculations
subroutine, public luts_conf_delete(conf)
subroutine, public calculate_aero_layers(aerosol_option, n_aerosols, n_profiles, n_layers, geovals, aero_layers, rh, layer_factors)
subroutine, public luts_conf_setup(conf, f_confopts, f_confoper)
character(len=maxvarlen), parameter, public var_prsi
character(len=maxvarlen), parameter, public var_prs
character(len=maxvarlen), parameter, public var_q
character(len=maxvarlen), parameter, public var_ts
fortran derived type for aod trajectory
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators