1 | ! radiation_aerosol_optics_description.F90 - Type to store aerosol optics metadata |
---|
2 | ! |
---|
3 | ! (C) Copyright 2022- ECMWF. |
---|
4 | ! |
---|
5 | ! This software is licensed under the terms of the Apache Licence Version 2.0 |
---|
6 | ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0. |
---|
7 | ! |
---|
8 | ! In applying this licence, ECMWF does not waive the privileges and immunities |
---|
9 | ! granted to it by virtue of its status as an intergovernmental organisation |
---|
10 | ! nor does it submit to any jurisdiction. |
---|
11 | ! |
---|
12 | ! Author: Robin Hogan |
---|
13 | ! Email: r.j.hogan@ecmwf.int |
---|
14 | ! |
---|
15 | |
---|
16 | module radiation_aerosol_optics_description |
---|
17 | |
---|
18 | use parkind1, only : jprb |
---|
19 | |
---|
20 | implicit none |
---|
21 | public |
---|
22 | |
---|
23 | !--------------------------------------------------------------------- |
---|
24 | ! This type holds the metadata from an aerosol optical property |
---|
25 | ! file, enabling the user to request the index to the aerosol type |
---|
26 | ! with particular properties. Note that string information is held |
---|
27 | ! in the form of 2D arrays of single characters, so comparison to |
---|
28 | ! character strings requires the to_string helper function at the |
---|
29 | ! end of this file. |
---|
30 | type aerosol_optics_description_type |
---|
31 | |
---|
32 | ! Two-character code describing the aerosol family, dimensioned |
---|
33 | ! (2,naer), e.g. |
---|
34 | ! SS: Sea salt |
---|
35 | ! OM: Organic matter |
---|
36 | ! SU: Sulfate |
---|
37 | ! OB: Secondary organic biogenic |
---|
38 | ! OA: Secondary organic anthropogenic |
---|
39 | ! AM: Fine-mode ammonium sulfate |
---|
40 | ! NI: Nitrate |
---|
41 | ! DD: Desert dust |
---|
42 | ! BC: Black carbon |
---|
43 | character(len=1), allocatable :: code_phobic(:,:) |
---|
44 | character(len=1), allocatable :: code_philic(:,:) |
---|
45 | |
---|
46 | ! Size bin, typically 1-2 or 1-3 in order fine to coarse, or zero |
---|
47 | ! if no division by size is used, dimensioned (naer) |
---|
48 | integer, allocatable :: bin_phobic(:) |
---|
49 | integer, allocatable :: bin_philic(:) |
---|
50 | |
---|
51 | ! Character string characterizing the optical model, e.g. OPAC, |
---|
52 | ! GACP, GLOMAP, Dubovik2002 etc. |
---|
53 | character(len=1), allocatable :: optical_model_phobic(:,:) |
---|
54 | character(len=1), allocatable :: optical_model_philic(:,:) |
---|
55 | |
---|
56 | ! The user can call preferred_optical_model to specify that a |
---|
57 | ! certain optical model for a certain aerosol family is to be |
---|
58 | ! preferred when get_index is called |
---|
59 | logical, allocatable :: is_preferred_phobic(:) |
---|
60 | logical, allocatable :: is_preferred_philic(:) |
---|
61 | |
---|
62 | contains |
---|
63 | procedure :: read |
---|
64 | procedure :: preferred_optical_model |
---|
65 | procedure :: get_index |
---|
66 | |
---|
67 | end type aerosol_optics_description_type |
---|
68 | |
---|
69 | contains |
---|
70 | |
---|
71 | !--------------------------------------------------------------------- |
---|
72 | ! Read optical property file file_name into an |
---|
73 | ! aerosol_optics_description_type object |
---|
74 | subroutine read(this, file_name, iverbose) |
---|
75 | |
---|
76 | use yomhook, only : lhook, dr_hook, jphook |
---|
77 | use easy_netcdf, only : netcdf_file |
---|
78 | |
---|
79 | class(aerosol_optics_description_type), intent(inout) :: this |
---|
80 | character(len=*), intent(in) :: file_name |
---|
81 | integer, intent(in), optional :: iverbose |
---|
82 | |
---|
83 | ! The NetCDF file containing the aerosol optics data |
---|
84 | type(netcdf_file) :: file |
---|
85 | |
---|
86 | real(jphook) :: hook_handle |
---|
87 | |
---|
88 | if (lhook) call dr_hook('radiation_aerosol_optics_description:load',0,hook_handle) |
---|
89 | |
---|
90 | ! Open the aerosol scattering file and configure the way it is |
---|
91 | ! read |
---|
92 | call file%open(trim(file_name), iverbose=iverbose) |
---|
93 | |
---|
94 | ! Read metadata variables |
---|
95 | call file%get('code_hydrophilic', this%code_philic) |
---|
96 | call file%get('code_hydrophobic', this%code_phobic) |
---|
97 | call file%get('bin_hydrophilic', this%bin_philic) |
---|
98 | call file%get('bin_hydrophobic', this%bin_phobic) |
---|
99 | call file%get('optical_model_hydrophilic', this%optical_model_philic) |
---|
100 | call file%get('optical_model_hydrophobic', this%optical_model_phobic) |
---|
101 | |
---|
102 | ! Allocate logical arrays of the appropriate size and set to FALSE |
---|
103 | allocate(this%is_preferred_philic(size(this%bin_philic))) |
---|
104 | allocate(this%is_preferred_phobic(size(this%bin_phobic))) |
---|
105 | this%is_preferred_philic = .false. |
---|
106 | this%is_preferred_phobic = .false. |
---|
107 | |
---|
108 | call file%close() |
---|
109 | |
---|
110 | if (lhook) call dr_hook('radiation_aerosol_optics_description:load',1,hook_handle) |
---|
111 | |
---|
112 | end subroutine read |
---|
113 | |
---|
114 | !--------------------------------------------------------------------- |
---|
115 | ! Specify the preferred optical model for a particular aerosol |
---|
116 | ! family, e.g. "call |
---|
117 | ! aer_desc%preferred_optical_model('DD','Woodward2001')" would mean |
---|
118 | ! that subsequent calls to get_index in which the optical model is |
---|
119 | ! not specified would return the Woodward model rather than the |
---|
120 | ! first matching model in the file. The check is only done on the |
---|
121 | ! first len(optical_model_str) characters, so "Woodward" and |
---|
122 | ! "Woodward2001" would both match the Woodward2001 model. |
---|
123 | subroutine preferred_optical_model(this, code_str, optical_model_str) |
---|
124 | |
---|
125 | use yomhook, only : lhook, dr_hook, jphook |
---|
126 | |
---|
127 | class(aerosol_optics_description_type), intent(inout) :: this |
---|
128 | character(len=2), intent(in) :: code_str |
---|
129 | character(len=*), intent(in) :: optical_model_str |
---|
130 | |
---|
131 | ! Aerosol loop counter |
---|
132 | integer :: ja |
---|
133 | |
---|
134 | real(jphook) :: hook_handle |
---|
135 | |
---|
136 | if (lhook) call dr_hook('radiation_aerosol_optics_description:preferred_optical_model',0,hook_handle) |
---|
137 | |
---|
138 | ! Check for empty string |
---|
139 | if (optical_model_str == ' ') then |
---|
140 | return |
---|
141 | end if |
---|
142 | |
---|
143 | ! Loop over hydrophilic types |
---|
144 | do ja = 1,size(this%bin_philic) |
---|
145 | ! Check if we have a match |
---|
146 | if (to_string(this%code_philic(:,ja)) == code_str & |
---|
147 | & .and. to_string(this%optical_model_philic(1:len(optical_model_str),ja)) & |
---|
148 | & == optical_model_str) then |
---|
149 | this%is_preferred_philic(ja) = .true. |
---|
150 | end if |
---|
151 | end do |
---|
152 | ! Repeat for the hydrophobic types |
---|
153 | do ja = 1,size(this%bin_phobic) |
---|
154 | if (to_string(this%code_phobic(:,ja)) == code_str & |
---|
155 | & .and. to_string(this%optical_model_phobic(1:len(optical_model_str),ja)) & |
---|
156 | & == optical_model_str) then |
---|
157 | this%is_preferred_phobic(ja) = .true. |
---|
158 | end if |
---|
159 | end do |
---|
160 | |
---|
161 | if (lhook) call dr_hook('radiation_aerosol_optics_description:preferred_optical_model',1,hook_handle) |
---|
162 | |
---|
163 | end subroutine preferred_optical_model |
---|
164 | |
---|
165 | !--------------------------------------------------------------------- |
---|
166 | ! Return the index to the aerosol optical properties corresponding |
---|
167 | ! to the aerosol family in code_str (e.g. SS, DD etc), whether or |
---|
168 | ! not the requested aerosol is hydrophilic in the logical |
---|
169 | ! lhydrophilic, and optionally the size bin ibin and optical model |
---|
170 | ! in optical_model_str. The return value may be used to populate the |
---|
171 | ! radiation_config%i_aerosol_map vector, where a positive number is |
---|
172 | ! a hydrophobic index, a negative number is a hydrophilic index and |
---|
173 | ! zero indicates that the aerosol type was not found in the file. |
---|
174 | ! This is a valid entry in i_aerosol_map meaning the aerosol is |
---|
175 | ! ignored, but the calling routine to get_index might wish to throw |
---|
176 | ! a warning or error. This routine works by assigning a score based |
---|
177 | ! on the closeness of the match. |
---|
178 | function get_index(this, code_str, lhydrophilic, ibin, optical_model_str) |
---|
179 | |
---|
180 | use yomhook, only : lhook, dr_hook, jphook |
---|
181 | use easy_netcdf, only : netcdf_file |
---|
182 | use radiation_io, only : nulout |
---|
183 | |
---|
184 | class(aerosol_optics_description_type), intent(in) :: this |
---|
185 | character(len=2), intent(in) :: code_str |
---|
186 | logical, intent(in) :: lhydrophilic |
---|
187 | integer, intent(in), optional :: ibin |
---|
188 | character(len=*), intent(in), optional :: optical_model_str |
---|
189 | |
---|
190 | ! Score of the currently selected aerosol index, and the score of |
---|
191 | ! the current one under consideration |
---|
192 | integer :: score, current_score |
---|
193 | |
---|
194 | ! Loop index for aerosol type |
---|
195 | integer :: ja |
---|
196 | |
---|
197 | ! Return value |
---|
198 | integer :: get_index |
---|
199 | |
---|
200 | ! Issue a warning if there is more than one equal match |
---|
201 | logical :: is_ambiguous |
---|
202 | |
---|
203 | real(jphook) :: hook_handle |
---|
204 | |
---|
205 | if (lhook) call dr_hook('radiation_aerosol_optics_description:get_index',0,hook_handle) |
---|
206 | |
---|
207 | ! Initial values |
---|
208 | get_index = 0 |
---|
209 | score = 0 |
---|
210 | is_ambiguous = .false. |
---|
211 | |
---|
212 | if (lhydrophilic) then |
---|
213 | ! Loop over hydrophilic aerosol types |
---|
214 | do ja = 1,size(this%bin_philic) |
---|
215 | current_score = 0 |
---|
216 | if (to_string(this%code_philic(:,ja)) == code_str) then |
---|
217 | ! Aerosol code matches |
---|
218 | if (present(ibin) .and. this%bin_philic(ja) > 0) then |
---|
219 | if (ibin > 0) then |
---|
220 | if (ibin == this%bin_philic(ja)) then |
---|
221 | ! Requested bin number matches |
---|
222 | current_score = 4 |
---|
223 | else |
---|
224 | ! Requested bin number does not match |
---|
225 | current_score = -1 |
---|
226 | end if |
---|
227 | else |
---|
228 | ! Bin number is zero: no request |
---|
229 | current_score = 2 |
---|
230 | end if |
---|
231 | else |
---|
232 | ! No bin number present |
---|
233 | current_score = 2 |
---|
234 | end if |
---|
235 | if (present(optical_model_str)) then |
---|
236 | if (to_string(this%optical_model_philic(1:len(optical_model_str),ja)) & |
---|
237 | & == optical_model_str) then |
---|
238 | ! Requested optical model matches |
---|
239 | if (current_score >= 0) then |
---|
240 | current_score = current_score + 4 |
---|
241 | end if |
---|
242 | else |
---|
243 | ! Requested optical model does not match |
---|
244 | current_score = -1 |
---|
245 | end if |
---|
246 | else if (current_score >= 0) then |
---|
247 | ! No requested optical model |
---|
248 | current_score = current_score + 2 |
---|
249 | end if |
---|
250 | if (current_score > 0 .and. this%is_preferred_philic(ja)) then |
---|
251 | current_score = current_score + 1 |
---|
252 | end if |
---|
253 | if (current_score > score) then |
---|
254 | ! Better score than any existing aerosol type |
---|
255 | get_index = -ja |
---|
256 | score = current_score |
---|
257 | is_ambiguous = .false. |
---|
258 | else if (current_score > 0 .and. current_score == score) then |
---|
259 | is_ambiguous = .true. |
---|
260 | end if |
---|
261 | end if |
---|
262 | end do |
---|
263 | else |
---|
264 | ! Loop over hydrophobic aerosol types |
---|
265 | do ja = 1,size(this%bin_phobic) |
---|
266 | current_score = 0 |
---|
267 | if (to_string(this%code_phobic(:,ja)) == code_str) then |
---|
268 | ! Aerosol code matches |
---|
269 | if (present(ibin) .and. this%bin_phobic(ja) > 0) then |
---|
270 | if (ibin > 0) then |
---|
271 | if (ibin == this%bin_phobic(ja)) then |
---|
272 | ! Requested bin number matches |
---|
273 | current_score = 4 |
---|
274 | else |
---|
275 | ! Requested bin number does not match |
---|
276 | current_score = -1 |
---|
277 | end if |
---|
278 | else |
---|
279 | ! Bin number is zero: no request |
---|
280 | current_score = 2 |
---|
281 | end if |
---|
282 | else |
---|
283 | ! No bin number requested or present |
---|
284 | current_score = 2 |
---|
285 | end if |
---|
286 | if (present(optical_model_str)) then |
---|
287 | if (to_string(this%optical_model_phobic(1:len(optical_model_str),ja)) & |
---|
288 | & == optical_model_str) then |
---|
289 | ! Requested optical model matches |
---|
290 | if (current_score >= 0) then |
---|
291 | current_score = current_score + 4 |
---|
292 | end if |
---|
293 | else |
---|
294 | ! Requested optical model does not match |
---|
295 | current_score = -1 |
---|
296 | end if |
---|
297 | else if (current_score >= 0) then |
---|
298 | ! No requested optical model |
---|
299 | current_score = current_score + 2 |
---|
300 | end if |
---|
301 | if (current_score > 0 .and. this%is_preferred_phobic(ja)) then |
---|
302 | current_score = current_score + 1 |
---|
303 | end if |
---|
304 | if (current_score > score) then |
---|
305 | ! Better score than any existing aerosol type |
---|
306 | get_index = ja |
---|
307 | score = current_score |
---|
308 | is_ambiguous = .false. |
---|
309 | else if (current_score > 0 .and. current_score == score) then |
---|
310 | is_ambiguous = .true. |
---|
311 | end if |
---|
312 | end if |
---|
313 | end do |
---|
314 | end if |
---|
315 | |
---|
316 | if (is_ambiguous) then |
---|
317 | write(nulout,'(a,a2,a,l1,a)') 'Warning: get_index("', code_str, '",', lhydrophilic, & |
---|
318 | & ',...) does not unambiguously identify an aerosol optical property index' |
---|
319 | end if |
---|
320 | |
---|
321 | if (lhook) call dr_hook('radiation_aerosol_optics_description:get_index',1,hook_handle) |
---|
322 | |
---|
323 | end function get_index |
---|
324 | |
---|
325 | !--------------------------------------------------------------------- |
---|
326 | ! Utility function to convert an array of single characters to a |
---|
327 | ! character string (yes Fortran's string handling is a bit rubbish) |
---|
328 | pure function to_string(arr) result(str) |
---|
329 | character, intent(in) :: arr(:) |
---|
330 | character(len=size(arr)) :: str |
---|
331 | integer :: jc |
---|
332 | do jc = 1,size(arr) |
---|
333 | str(jc:jc) = arr(jc) |
---|
334 | end do |
---|
335 | end function to_string |
---|
336 | |
---|
337 | end module radiation_aerosol_optics_description |
---|