[4773] | 1 | ! radiation_aerosol_optics.F90 - Computing aerosol optical properties |
---|
| 2 | ! |
---|
| 3 | ! (C) Copyright 2015- 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 | ! Modifications |
---|
| 16 | ! 2018-04-15 R. Hogan Add "direct" option |
---|
| 17 | ! 2020-11-14 R. Hogan Add setup_general_aerosol_optics for ecCKD compatibility |
---|
| 18 | ! 2022-03-27 R. Hogan Add setup_general_aerosol_optics_legacy to use RRTM aerosol files with ecCKD |
---|
| 19 | ! 2022-11-22 P. Ukkonen / R. Hogan Optimizations to enhance vectorization |
---|
| 20 | |
---|
[4853] | 21 | #include "ecrad_config.h" |
---|
| 22 | |
---|
[4773] | 23 | module radiation_aerosol_optics |
---|
| 24 | |
---|
| 25 | implicit none |
---|
| 26 | public |
---|
| 27 | |
---|
| 28 | contains |
---|
| 29 | |
---|
| 30 | ! Provides the elemental function "delta_eddington_extensive" |
---|
| 31 | #include "radiation_delta_eddington.h" |
---|
| 32 | |
---|
| 33 | !--------------------------------------------------------------------- |
---|
| 34 | ! Load aerosol scattering data; this subroutine delegates to one |
---|
| 35 | ! in radiation_aerosol_optics_data.F90 |
---|
| 36 | subroutine setup_aerosol_optics(config) |
---|
| 37 | |
---|
| 38 | use parkind1, only : jprb |
---|
| 39 | use yomhook, only : lhook, dr_hook, jphook |
---|
| 40 | use radiation_config, only : config_type |
---|
| 41 | use radiation_aerosol_optics_data, only : aerosol_optics_type |
---|
| 42 | use radiation_io, only : nulerr, radiation_abort |
---|
| 43 | |
---|
| 44 | type(config_type), intent(inout) :: config |
---|
| 45 | |
---|
| 46 | real(jphook) :: hook_handle |
---|
| 47 | |
---|
| 48 | if (lhook) call dr_hook('radiation_aerosol_optics:setup_aerosol_optics',0,hook_handle) |
---|
| 49 | |
---|
| 50 | if (config%n_aerosol_types > 0) then |
---|
| 51 | ! Load data from file and prepare to map config%n_aerosol_types |
---|
| 52 | ! aerosol types |
---|
| 53 | if (config%use_general_aerosol_optics) then |
---|
| 54 | ! Read file containing high spectral resolution optical |
---|
| 55 | ! properties and average to the spectral intervals of the |
---|
| 56 | ! current gas-optics scheme |
---|
| 57 | call setup_general_aerosol_optics(config) |
---|
| 58 | else |
---|
| 59 | ! Read file containing optical properties already in the bands |
---|
| 60 | ! of the gas-optics scheme |
---|
| 61 | call config%aerosol_optics%setup(trim(config%aerosol_optics_file_name), & |
---|
| 62 | & iverbose=config%iverbosesetup) |
---|
| 63 | end if |
---|
| 64 | |
---|
| 65 | call config%aerosol_optics%initialize_types(config%n_aerosol_types) |
---|
| 66 | |
---|
| 67 | ! Check agreement in number of bands |
---|
| 68 | if (config%n_bands_lw /= config%aerosol_optics%n_bands_lw) then |
---|
| 69 | write(nulerr,'(a,i0,a,i0,a)') '*** Error: number of longwave bands (', & |
---|
| 70 | & config%n_bands_lw, ') does not match aerosol optics look-up table (', & |
---|
| 71 | & config%aerosol_optics%n_bands_lw, ')' |
---|
| 72 | call radiation_abort() |
---|
| 73 | end if |
---|
| 74 | if (config%n_bands_sw /= config%aerosol_optics%n_bands_sw) then |
---|
| 75 | write(nulerr,'(a)') '*** Error: number of shortwave bands does not match aerosol optics look-up table' |
---|
| 76 | call radiation_abort() |
---|
| 77 | end if |
---|
| 78 | |
---|
| 79 | ! Map aerosol types to those loaded from the data file |
---|
| 80 | call config%aerosol_optics%set_types(config%i_aerosol_type_map(1:config%n_aerosol_types)) |
---|
| 81 | end if |
---|
| 82 | |
---|
| 83 | if (config%iverbosesetup >= 1) then |
---|
| 84 | call config%aerosol_optics%print_description(config%i_aerosol_type_map(1:config%n_aerosol_types)) |
---|
| 85 | end if |
---|
| 86 | |
---|
| 87 | if (lhook) call dr_hook('radiation_aerosol_optics:setup_aerosol_optics',1,hook_handle) |
---|
| 88 | |
---|
| 89 | end subroutine setup_aerosol_optics |
---|
| 90 | |
---|
| 91 | |
---|
| 92 | !--------------------------------------------------------------------- |
---|
| 93 | ! Read file containing high spectral resolution optical properties |
---|
| 94 | ! and average to the spectral intervals of the current gas-optics |
---|
| 95 | ! scheme |
---|
| 96 | subroutine setup_general_aerosol_optics(config) |
---|
| 97 | |
---|
| 98 | use parkind1, only : jprb |
---|
| 99 | use yomhook, only : lhook, dr_hook, jphook |
---|
[4853] | 100 | #ifdef EASY_NETCDF_READ_MPI |
---|
| 101 | use easy_netcdf_read_mpi, only : netcdf_file |
---|
| 102 | #else |
---|
[4773] | 103 | use easy_netcdf, only : netcdf_file |
---|
[4853] | 104 | #endif |
---|
[4773] | 105 | use radiation_config, only : config_type |
---|
| 106 | use radiation_aerosol_optics_data, only : aerosol_optics_type |
---|
| 107 | use radiation_spectral_definition, only : SolarReferenceTemperature, & |
---|
| 108 | & TerrestrialReferenceTemperature |
---|
| 109 | use radiation_io, only : nulout |
---|
| 110 | |
---|
| 111 | type(config_type), intent(inout), target :: config |
---|
| 112 | |
---|
| 113 | ! The NetCDF file containing the aerosol optics data |
---|
| 114 | type(netcdf_file) :: file |
---|
| 115 | |
---|
| 116 | ! Wavenumber points in NetCDF file |
---|
| 117 | real(jprb), allocatable :: wavenumber(:) ! cm-1 |
---|
| 118 | |
---|
| 119 | ! Hydrophilic aerosol properties |
---|
| 120 | real(jprb), allocatable :: mass_ext_philic(:,:,:) ! Mass-ext coefficient (m2 kg-1) |
---|
| 121 | real(jprb), allocatable :: ssa_philic(:,:,:) ! Single-scattering albedo |
---|
| 122 | real(jprb), allocatable :: g_philic(:,:,:) ! Asymmetry factor |
---|
| 123 | real(jprb), allocatable :: lidar_ratio_philic(:,:,:) ! Lidar ratio (sr) |
---|
| 124 | |
---|
| 125 | ! Hydrophobic aerosol properties |
---|
| 126 | real(jprb), allocatable :: mass_ext_phobic(:,:) ! Mass-ext coefficient (m2 kg-1) |
---|
| 127 | real(jprb), allocatable :: ssa_phobic(:,:) ! Single-scattering albedo |
---|
| 128 | real(jprb), allocatable :: g_phobic(:,:) ! Asymmetry factor |
---|
| 129 | real(jprb), allocatable :: lidar_ratio_phobic(:,:) ! Lidar ratio (sr) |
---|
| 130 | |
---|
| 131 | ! Mapping matrix between optical properties at the wavenumbers in |
---|
| 132 | ! the file, and spectral intervals used by the gas-optics scheme |
---|
| 133 | real(jprb), allocatable :: mapping(:,:) |
---|
| 134 | |
---|
| 135 | ! Target monochromatic wavenumber for interpolation (cm-1) |
---|
| 136 | real(jprb) :: wavenumber_target |
---|
| 137 | |
---|
| 138 | ! Number of spectral points describing aerosol properties in the |
---|
| 139 | ! shortwave and longwave |
---|
| 140 | integer :: nspecsw, nspeclw |
---|
| 141 | |
---|
| 142 | ! Number of monochromatic wavelengths required |
---|
| 143 | integer :: nmono |
---|
| 144 | |
---|
| 145 | integer :: n_type_philic, n_type_phobic, nrh, nwn |
---|
| 146 | integer :: jtype, jwl, iwn |
---|
| 147 | |
---|
| 148 | ! Weight of first point in interpolation |
---|
| 149 | real(jprb) :: weight1 |
---|
| 150 | |
---|
| 151 | real(jphook) :: hook_handle |
---|
| 152 | |
---|
| 153 | if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics',0,hook_handle) |
---|
| 154 | |
---|
| 155 | associate(ao => config%aerosol_optics) |
---|
| 156 | |
---|
| 157 | call file%open(trim(config%aerosol_optics_file_name), iverbose=config%iverbosesetup) |
---|
| 158 | |
---|
| 159 | if (.not. file%exists('wavenumber')) then |
---|
| 160 | ! Assume we have an old-style aerosol optics file with optical |
---|
| 161 | ! properties provided per pre-defined band |
---|
| 162 | call file%close() |
---|
| 163 | if (config%iverbosesetup >= 2) then |
---|
| 164 | write(nulout,'(a)') 'Legacy aerosol optics file: mapping between bands' |
---|
| 165 | end if |
---|
| 166 | call setup_general_aerosol_optics_legacy(config, trim(config%aerosol_optics_file_name)) |
---|
| 167 | if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics',1,hook_handle) |
---|
| 168 | return |
---|
| 169 | end if |
---|
| 170 | |
---|
| 171 | if (file%exists('mass_ext_hydrophilic')) then |
---|
| 172 | ao%use_hydrophilic = .true. |
---|
| 173 | else |
---|
| 174 | ao%use_hydrophilic = .false. |
---|
| 175 | end if |
---|
| 176 | |
---|
| 177 | call file%get('wavenumber', wavenumber) |
---|
| 178 | nwn = size(wavenumber) |
---|
| 179 | |
---|
| 180 | ! Read the raw scattering data |
---|
| 181 | call file%get('mass_ext_hydrophobic', mass_ext_phobic) |
---|
| 182 | call file%get('ssa_hydrophobic', ssa_phobic) |
---|
| 183 | call file%get('asymmetry_hydrophobic', g_phobic) |
---|
| 184 | call file%get('lidar_ratio_hydrophobic', lidar_ratio_phobic) |
---|
| 185 | |
---|
| 186 | call file%get_global_attribute('description_hydrophobic', & |
---|
| 187 | & ao%description_phobic_str) |
---|
| 188 | |
---|
| 189 | if (ao%use_hydrophilic) then |
---|
| 190 | call file%get('mass_ext_hydrophilic', mass_ext_philic) |
---|
| 191 | call file%get('ssa_hydrophilic', ssa_philic) |
---|
| 192 | call file%get('asymmetry_hydrophilic', g_philic) |
---|
| 193 | call file%get('lidar_ratio_hydrophilic', lidar_ratio_philic) |
---|
| 194 | |
---|
| 195 | call file%get('relative_humidity1', ao%rh_lower) |
---|
| 196 | |
---|
| 197 | call file%get_global_attribute('description_hydrophilic', & |
---|
| 198 | & ao%description_philic_str) |
---|
| 199 | end if |
---|
| 200 | |
---|
| 201 | ! Close aerosol scattering file |
---|
| 202 | call file%close() |
---|
| 203 | |
---|
| 204 | n_type_phobic = size(mass_ext_phobic, 2) |
---|
| 205 | if (ao%use_hydrophilic) then |
---|
| 206 | n_type_philic = size(mass_ext_philic, 3) |
---|
| 207 | nrh = size(ao%rh_lower) |
---|
| 208 | else |
---|
| 209 | n_type_philic = 0 |
---|
| 210 | nrh = 0 |
---|
| 211 | end if |
---|
| 212 | |
---|
| 213 | if (config%do_cloud_aerosol_per_sw_g_point) then |
---|
| 214 | nspecsw = config%gas_optics_sw%spectral_def%ng |
---|
| 215 | else |
---|
| 216 | nspecsw = config%gas_optics_sw%spectral_def%nband |
---|
| 217 | end if |
---|
| 218 | |
---|
| 219 | if (config%do_cloud_aerosol_per_lw_g_point) then |
---|
| 220 | nspeclw = config%gas_optics_lw%spectral_def%ng |
---|
| 221 | else |
---|
| 222 | nspeclw = config%gas_optics_lw%spectral_def%nband |
---|
| 223 | end if |
---|
| 224 | |
---|
| 225 | if (allocated(ao%wavelength_mono)) then |
---|
| 226 | ! Monochromatic wavelengths also required |
---|
| 227 | nmono = size(ao%wavelength_mono) |
---|
| 228 | else |
---|
| 229 | nmono = 0 |
---|
| 230 | end if |
---|
| 231 | |
---|
| 232 | call ao%allocate(n_type_phobic, n_type_philic, nrh, nspeclw, nspecsw, nmono) |
---|
| 233 | |
---|
| 234 | if (config%do_sw) then |
---|
| 235 | call config%gas_optics_sw%spectral_def%calc_mapping(wavenumber, mapping, & |
---|
| 236 | & use_bands=(.not. config%do_cloud_aerosol_per_sw_g_point)) |
---|
| 237 | |
---|
| 238 | ao%mass_ext_sw_phobic = matmul(mapping, mass_ext_phobic) |
---|
| 239 | ao%ssa_sw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic) & |
---|
| 240 | & / ao%mass_ext_sw_phobic |
---|
| 241 | ao%g_sw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic*g_phobic) & |
---|
| 242 | & / (ao%mass_ext_sw_phobic*ao%ssa_sw_phobic) |
---|
| 243 | |
---|
| 244 | if (ao%use_hydrophilic) then |
---|
| 245 | do jtype = 1,n_type_philic |
---|
| 246 | ao%mass_ext_sw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype)) |
---|
| 247 | ao%ssa_sw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) & |
---|
| 248 | & *ssa_philic(:,:,jtype)) & |
---|
| 249 | & / ao%mass_ext_sw_philic(:,:,jtype) |
---|
| 250 | ao%g_sw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) & |
---|
| 251 | & *ssa_philic(:,:,jtype)*g_philic(:,:,jtype)) & |
---|
| 252 | & / (ao%mass_ext_sw_philic(:,:,jtype)*ao%ssa_sw_philic(:,:,jtype)) |
---|
| 253 | end do |
---|
| 254 | end if |
---|
| 255 | end if |
---|
| 256 | |
---|
| 257 | if (config%do_lw) then |
---|
| 258 | call config%gas_optics_lw%spectral_def%calc_mapping(wavenumber, mapping, & |
---|
| 259 | & use_bands=(.not. config%do_cloud_aerosol_per_lw_g_point)) |
---|
| 260 | |
---|
| 261 | ao%mass_ext_lw_phobic = matmul(mapping, mass_ext_phobic) |
---|
| 262 | ao%ssa_lw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic) & |
---|
| 263 | & / ao%mass_ext_lw_phobic |
---|
| 264 | ao%g_lw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic*g_phobic) & |
---|
| 265 | & / (ao%mass_ext_lw_phobic*ao%ssa_lw_phobic) |
---|
| 266 | |
---|
| 267 | if (ao%use_hydrophilic) then |
---|
| 268 | do jtype = 1,n_type_philic |
---|
| 269 | ao%mass_ext_lw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype)) |
---|
| 270 | ao%ssa_lw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) & |
---|
| 271 | & *ssa_philic(:,:,jtype)) & |
---|
| 272 | & / ao%mass_ext_lw_philic(:,:,jtype) |
---|
| 273 | ao%g_lw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) & |
---|
| 274 | & *ssa_philic(:,:,jtype)*g_philic(:,:,jtype)) & |
---|
| 275 | & / (ao%mass_ext_lw_philic(:,:,jtype)*ao%ssa_lw_philic(:,:,jtype)) |
---|
| 276 | end do |
---|
| 277 | end if |
---|
| 278 | end if |
---|
| 279 | |
---|
| 280 | if (allocated(ao%wavelength_mono)) then |
---|
| 281 | ! Monochromatic wavelengths also required |
---|
| 282 | do jwl = 1,nmono |
---|
| 283 | ! Wavelength (m) to wavenumber (cm-1) |
---|
| 284 | wavenumber_target = 0.01_jprb / ao%wavelength_mono(jwl) |
---|
| 285 | ! Find index to first interpolation point, and its weight |
---|
| 286 | if (wavenumber_target <= wavenumber(1)) then |
---|
| 287 | weight1 = 1.0_jprb |
---|
| 288 | iwn = 1 |
---|
| 289 | else if (wavenumber_target >= wavenumber(nwn)) then |
---|
| 290 | iwn = nwn-1 |
---|
| 291 | weight1 = 0.0_jprb |
---|
| 292 | else |
---|
| 293 | iwn = 1 |
---|
| 294 | do while (wavenumber(iwn+1) < wavenumber_target .and. iwn < nwn-1) |
---|
| 295 | iwn = iwn + 1 |
---|
| 296 | end do |
---|
| 297 | weight1 = (wavenumber(iwn+1)-wavenumber_target) & |
---|
| 298 | & / (wavenumber(iwn+1)-wavenumber(iwn)) |
---|
| 299 | end if |
---|
| 300 | ! Linear interpolation |
---|
| 301 | ao%mass_ext_mono_phobic(jwl,:) = weight1 * mass_ext_phobic(iwn,:) & |
---|
| 302 | & + (1.0_jprb - weight1)* mass_ext_phobic(iwn+1,:) |
---|
| 303 | ao%ssa_mono_phobic(jwl,:) = weight1 * ssa_phobic(iwn,:) & |
---|
| 304 | & + (1.0_jprb - weight1)* ssa_phobic(iwn+1,:) |
---|
| 305 | ao%g_mono_phobic(jwl,:) = weight1 * g_phobic(iwn,:) & |
---|
| 306 | & + (1.0_jprb - weight1)* g_phobic(iwn+1,:) |
---|
| 307 | ao%lidar_ratio_mono_phobic(jwl,:) = weight1 * lidar_ratio_phobic(iwn,:) & |
---|
| 308 | & + (1.0_jprb - weight1)* lidar_ratio_phobic(iwn+1,:) |
---|
| 309 | if (ao%use_hydrophilic) then |
---|
| 310 | ao%mass_ext_mono_philic(jwl,:,:) = weight1 * mass_ext_philic(iwn,:,:) & |
---|
| 311 | & + (1.0_jprb - weight1)* mass_ext_philic(iwn+1,:,:) |
---|
| 312 | ao%ssa_mono_philic(jwl,:,:) = weight1 * ssa_philic(iwn,:,:) & |
---|
| 313 | & + (1.0_jprb - weight1)* ssa_philic(iwn+1,:,:) |
---|
| 314 | ao%g_mono_philic(jwl,:,:) = weight1 * g_philic(iwn,:,:) & |
---|
| 315 | & + (1.0_jprb - weight1)* g_philic(iwn+1,:,:) |
---|
| 316 | ao%lidar_ratio_mono_philic(jwl,:,:) = weight1 * lidar_ratio_philic(iwn,:,:) & |
---|
| 317 | & + (1.0_jprb - weight1)* lidar_ratio_philic(iwn+1,:,:) |
---|
| 318 | end if |
---|
| 319 | end do |
---|
| 320 | end if |
---|
| 321 | |
---|
| 322 | ! Deallocate memory local to this routine |
---|
| 323 | deallocate(mass_ext_phobic) |
---|
| 324 | deallocate(ssa_phobic) |
---|
| 325 | deallocate(g_phobic) |
---|
| 326 | deallocate(lidar_ratio_phobic) |
---|
| 327 | if (ao%use_hydrophilic) then |
---|
| 328 | deallocate(mass_ext_philic) |
---|
| 329 | deallocate(ssa_philic) |
---|
| 330 | deallocate(g_philic) |
---|
| 331 | deallocate(lidar_ratio_philic) |
---|
| 332 | end if |
---|
| 333 | |
---|
| 334 | end associate |
---|
| 335 | |
---|
| 336 | if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics',1,hook_handle) |
---|
| 337 | |
---|
| 338 | end subroutine setup_general_aerosol_optics |
---|
| 339 | |
---|
| 340 | |
---|
| 341 | !--------------------------------------------------------------------- |
---|
| 342 | ! Read file containing legacy-style band-wise aerosol optical |
---|
| 343 | ! properties and average to the spectral intervals of the current |
---|
| 344 | ! gas-optics scheme |
---|
| 345 | subroutine setup_general_aerosol_optics_legacy(config, file_name) |
---|
| 346 | |
---|
| 347 | use parkind1, only : jprb |
---|
| 348 | use yomhook, only : lhook, dr_hook, jphook |
---|
[4853] | 349 | #ifdef EASY_NETCDF_READ_MPI |
---|
| 350 | use easy_netcdf_read_mpi, only : netcdf_file |
---|
| 351 | #else |
---|
[4773] | 352 | use easy_netcdf, only : netcdf_file |
---|
[4853] | 353 | #endif |
---|
[4773] | 354 | use radiation_config, only : config_type |
---|
| 355 | use radiation_aerosol_optics_data, only : aerosol_optics_type |
---|
| 356 | use radiation_spectral_definition, only : SolarReferenceTemperature, & |
---|
| 357 | & TerrestrialReferenceTemperature |
---|
| 358 | |
---|
| 359 | type(config_type), intent(inout), target :: config |
---|
| 360 | |
---|
| 361 | ! The NetCDF file containing the aerosol optics data |
---|
| 362 | character(len=*), intent(in) :: file_name |
---|
| 363 | |
---|
| 364 | ! Mapping matrix between optical properties at the wavenumbers in |
---|
| 365 | ! the file, and spectral intervals used by the gas-optics scheme |
---|
| 366 | real(jprb), allocatable :: mapping(:,:), mapping_transp(:,:) |
---|
| 367 | |
---|
| 368 | ! Pointer to the aerosol optics coefficients for brevity of access |
---|
| 369 | type(aerosol_optics_type), pointer :: ao |
---|
| 370 | |
---|
| 371 | ! Local copy of aerosol optical properties in the spectral |
---|
| 372 | ! intervals of the file, which is deallocated when it goes out of |
---|
| 373 | ! scope |
---|
| 374 | type(aerosol_optics_type) :: ao_legacy |
---|
| 375 | |
---|
| 376 | integer :: jtype |
---|
| 377 | |
---|
| 378 | real(jphook) :: hook_handle |
---|
| 379 | |
---|
| 380 | if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics_legacy',0,hook_handle) |
---|
| 381 | ao => config%aerosol_optics |
---|
| 382 | |
---|
| 383 | ! Load file into a local structure |
---|
| 384 | call ao_legacy%setup(file_name, iverbose=config%iverbosesetup) |
---|
| 385 | |
---|
| 386 | ! Copy over scalars and coordinate variables |
---|
| 387 | call ao%allocate(ao_legacy%n_type_phobic, ao_legacy%n_type_philic, ao_legacy%nrh, & |
---|
| 388 | & config%n_bands_lw, config%n_bands_sw, ao_legacy%n_mono_wl) |
---|
| 389 | ao%description_phobic_str = ao_legacy%description_phobic_str |
---|
| 390 | if (ao_legacy%use_hydrophilic) then |
---|
| 391 | ao%description_philic_str = ao_legacy%description_philic_str |
---|
| 392 | ao%rh_lower = ao_legacy%rh_lower |
---|
| 393 | end if |
---|
| 394 | |
---|
| 395 | ! use_hydrophilic = ao_legacy%use_hydrophilic |
---|
| 396 | ! ao%iclass = ao_legacy%iclass |
---|
| 397 | ! ao%itype = ao_legacy%itype |
---|
| 398 | ! ao%ntype = ao_legacy%ntype |
---|
| 399 | ! ao%n_type_phobic = ao_legacy%n_type_phobic |
---|
| 400 | ! ao%n_type_philic = ao_legacy%n_type_philic |
---|
| 401 | ! ao%n_mono_wl = ao_legacy%n_mono_wl |
---|
| 402 | ! ao%use_monochromatic = ao_legacy%use_monochromatic |
---|
| 403 | |
---|
| 404 | if (config%do_sw) then |
---|
| 405 | call config%gas_optics_sw%spectral_def%calc_mapping_from_wavenumber_bands( & |
---|
| 406 | & ao_legacy%wavenumber1_sw, ao_legacy%wavenumber2_sw, mapping_transp, & |
---|
| 407 | & use_bands=(.not. config%do_cloud_aerosol_per_sw_g_point)) |
---|
| 408 | if (allocated(mapping)) then |
---|
| 409 | deallocate(mapping) |
---|
| 410 | end if |
---|
| 411 | allocate(mapping(config%n_bands_sw,ao_legacy%n_bands_sw)) |
---|
| 412 | mapping = transpose(mapping_transp) |
---|
| 413 | ao%mass_ext_sw_phobic = matmul(mapping, ao_legacy%mass_ext_sw_phobic) |
---|
| 414 | ao%ssa_sw_phobic = matmul(mapping, ao_legacy%mass_ext_sw_phobic*ao_legacy%ssa_sw_phobic) & |
---|
| 415 | & / ao%mass_ext_sw_phobic |
---|
| 416 | ao%g_sw_phobic = matmul(mapping, ao_legacy%mass_ext_sw_phobic*ao_legacy%ssa_sw_phobic & |
---|
| 417 | & *ao_legacy%g_sw_phobic) & |
---|
| 418 | & / (ao%mass_ext_sw_phobic*ao%ssa_sw_phobic) |
---|
| 419 | |
---|
| 420 | if (ao%use_hydrophilic) then |
---|
| 421 | do jtype = 1,ao%n_type_philic |
---|
| 422 | ao%mass_ext_sw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_sw_philic(:,:,jtype)) |
---|
| 423 | ao%ssa_sw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_sw_philic(:,:,jtype) & |
---|
| 424 | & *ao_legacy%ssa_sw_philic(:,:,jtype)) & |
---|
| 425 | & / ao%mass_ext_sw_philic(:,:,jtype) |
---|
| 426 | ao%g_sw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_sw_philic(:,:,jtype) & |
---|
| 427 | & *ao_legacy%ssa_sw_philic(:,:,jtype)*ao_legacy%g_sw_philic(:,:,jtype)) & |
---|
| 428 | & / (ao%mass_ext_sw_philic(:,:,jtype)*ao%ssa_sw_philic(:,:,jtype)) |
---|
| 429 | end do |
---|
| 430 | end if |
---|
| 431 | end if |
---|
| 432 | |
---|
| 433 | if (config%do_lw) then |
---|
| 434 | if (allocated(mapping_transp)) then |
---|
| 435 | deallocate(mapping_transp) |
---|
| 436 | end if |
---|
| 437 | call config%gas_optics_lw%spectral_def%calc_mapping_from_wavenumber_bands( & |
---|
| 438 | & ao_legacy%wavenumber1_lw, ao_legacy%wavenumber2_lw, mapping_transp, & |
---|
| 439 | & use_bands=(.not. config%do_cloud_aerosol_per_lw_g_point)) |
---|
| 440 | if (allocated(mapping)) then |
---|
| 441 | deallocate(mapping) |
---|
| 442 | end if |
---|
| 443 | allocate(mapping(config%n_bands_lw,ao_legacy%n_bands_lw)) |
---|
| 444 | mapping = transpose(mapping_transp) |
---|
| 445 | ao%mass_ext_lw_phobic = matmul(mapping, ao_legacy%mass_ext_lw_phobic) |
---|
| 446 | ao%ssa_lw_phobic = matmul(mapping, ao_legacy%mass_ext_lw_phobic*ao_legacy%ssa_lw_phobic) & |
---|
| 447 | & / ao%mass_ext_lw_phobic |
---|
| 448 | ao%g_lw_phobic = matmul(mapping, ao_legacy%mass_ext_lw_phobic*ao_legacy%ssa_lw_phobic & |
---|
| 449 | & *ao_legacy%g_lw_phobic) & |
---|
| 450 | & / (ao%mass_ext_lw_phobic*ao%ssa_lw_phobic) |
---|
| 451 | |
---|
| 452 | if (ao%use_hydrophilic) then |
---|
| 453 | do jtype = 1,ao%n_type_philic |
---|
| 454 | ao%mass_ext_lw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_lw_philic(:,:,jtype)) |
---|
| 455 | ao%ssa_lw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_lw_philic(:,:,jtype) & |
---|
| 456 | & *ao_legacy%ssa_lw_philic(:,:,jtype)) & |
---|
| 457 | & / ao%mass_ext_lw_philic(:,:,jtype) |
---|
| 458 | ao%g_lw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_lw_philic(:,:,jtype) & |
---|
| 459 | & *ao_legacy%ssa_lw_philic(:,:,jtype)*ao_legacy%g_lw_philic(:,:,jtype)) & |
---|
| 460 | & / (ao%mass_ext_lw_philic(:,:,jtype)*ao%ssa_lw_philic(:,:,jtype)) |
---|
| 461 | end do |
---|
| 462 | end if |
---|
| 463 | end if |
---|
| 464 | |
---|
| 465 | if (allocated(ao_legacy%wavelength_mono)) then |
---|
| 466 | ao%wavelength_mono = ao_legacy%wavelength_mono |
---|
| 467 | ao%mass_ext_mono_phobic = ao_legacy%mass_ext_mono_phobic |
---|
| 468 | ao%ssa_mono_phobic = ao_legacy%ssa_mono_phobic |
---|
| 469 | ao%g_mono_phobic = ao_legacy%g_mono_phobic |
---|
| 470 | ao%lidar_ratio_mono_phobic = ao_legacy%lidar_ratio_mono_phobic |
---|
| 471 | if (ao%use_hydrophilic) then |
---|
| 472 | ao%mass_ext_mono_philic = ao_legacy%mass_ext_mono_philic |
---|
| 473 | ao%ssa_mono_philic = ao_legacy%ssa_mono_philic |
---|
| 474 | ao%g_mono_philic = ao_legacy%g_mono_philic |
---|
| 475 | ao%lidar_ratio_mono_philic = ao_legacy%lidar_ratio_mono_philic |
---|
| 476 | end if |
---|
| 477 | end if |
---|
| 478 | |
---|
| 479 | if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics_legacy',1,hook_handle) |
---|
| 480 | |
---|
| 481 | end subroutine setup_general_aerosol_optics_legacy |
---|
| 482 | |
---|
| 483 | |
---|
| 484 | !--------------------------------------------------------------------- |
---|
| 485 | ! Compute aerosol optical properties and add to existing gas optical |
---|
| 486 | ! depth and scattering properties |
---|
| 487 | subroutine add_aerosol_optics(nlev,istartcol,iendcol, & |
---|
| 488 | & config, thermodynamics, gas, aerosol, & |
---|
| 489 | & od_lw, ssa_lw, g_lw, od_sw, ssa_sw, g_sw) |
---|
| 490 | |
---|
| 491 | use parkind1, only : jprb |
---|
| 492 | use radiation_io, only : nulout, nulerr, radiation_abort |
---|
| 493 | use yomhook, only : lhook, dr_hook, jphook |
---|
| 494 | use radiation_config, only : config_type |
---|
| 495 | use radiation_thermodynamics, only : thermodynamics_type |
---|
| 496 | use radiation_gas, only : gas_type, IH2O, IMassMixingRatio |
---|
| 497 | use radiation_aerosol, only : aerosol_type |
---|
| 498 | use radiation_constants, only : AccelDueToGravity |
---|
| 499 | use radiation_aerosol_optics_data, only : aerosol_optics_type, & |
---|
| 500 | & IAerosolClassUndefined, IAerosolClassIgnored, & |
---|
| 501 | & IAerosolClassHydrophobic, IAerosolClassHydrophilic |
---|
| 502 | |
---|
| 503 | real(jprb), parameter :: OneOverAccelDueToGravity = 1.0_jprb / AccelDueToGravity |
---|
| 504 | |
---|
| 505 | integer, intent(in) :: nlev ! number of model levels |
---|
| 506 | integer, intent(in) :: istartcol, iendcol ! range of columns to process |
---|
| 507 | type(config_type), intent(in), target :: config |
---|
| 508 | type(thermodynamics_type),intent(in) :: thermodynamics |
---|
| 509 | type(gas_type), intent(in) :: gas |
---|
| 510 | type(aerosol_type), intent(in) :: aerosol |
---|
| 511 | ! Optical depth, single scattering albedo and asymmetry factor of |
---|
| 512 | ! the atmosphere (gases on input, gases and aerosols on output) |
---|
| 513 | ! for each g point. Note that longwave ssa and asymmetry and |
---|
| 514 | ! shortwave asymmetry are all zero for gases, so are not yet |
---|
| 515 | ! defined on input and are therefore intent(out). |
---|
| 516 | real(jprb), dimension(config%n_g_lw,nlev,istartcol:iendcol), & |
---|
| 517 | & intent(inout) :: od_lw |
---|
| 518 | real(jprb), dimension(config%n_g_lw_if_scattering,nlev,istartcol:iendcol), & |
---|
| 519 | & intent(out) :: ssa_lw, g_lw |
---|
| 520 | real(jprb), dimension(config%n_g_sw,nlev,istartcol:iendcol), & |
---|
| 521 | & intent(inout) :: od_sw, ssa_sw |
---|
| 522 | real(jprb), dimension(config%n_g_sw,nlev,istartcol:iendcol), & |
---|
| 523 | & intent(out) :: g_sw |
---|
| 524 | |
---|
| 525 | ! Extinction optical depth, scattering optical depth and |
---|
| 526 | ! asymmetry-times-scattering-optical-depth for all the aerosols in |
---|
| 527 | ! a column for each spectral band of the shortwave and longwave |
---|
| 528 | ! spectrum |
---|
| 529 | real(jprb), dimension(config%n_bands_sw,nlev) & |
---|
| 530 | & :: od_sw_aerosol, scat_sw_aerosol, scat_g_sw_aerosol |
---|
| 531 | real(jprb), dimension(config%n_bands_lw,nlev) & |
---|
| 532 | & :: od_lw_aerosol |
---|
| 533 | real(jprb), dimension(config%n_bands_lw_if_scattering,nlev) & |
---|
| 534 | & :: scat_lw_aerosol, scat_g_lw_aerosol |
---|
| 535 | |
---|
| 536 | real(jprb) :: local_od_sw, local_od_lw |
---|
| 537 | |
---|
| 538 | real(jprb) :: h2o_mmr(istartcol:iendcol,nlev) |
---|
| 539 | |
---|
| 540 | real(jprb) :: rh ! Relative humidity with respect to liquid water |
---|
| 541 | |
---|
| 542 | ! Factor (kg m-2) to convert mixing ratio (kg kg-1) to mass in |
---|
| 543 | ! path (kg m-2) |
---|
| 544 | real(jprb) :: factor(nlev) |
---|
| 545 | |
---|
| 546 | ! Temporary extinction and scattering optical depths of aerosol |
---|
| 547 | ! plus gas |
---|
| 548 | real(jprb) :: local_od, local_scat |
---|
| 549 | |
---|
| 550 | ! Aerosol mixing ratio as a scalar |
---|
| 551 | real(jprb) :: mixing_ratio |
---|
| 552 | |
---|
| 553 | ! Loop indices for column, level, g point, band and aerosol type |
---|
| 554 | integer :: jcol, jlev, jg, jtype, jband |
---|
| 555 | |
---|
| 556 | ! Range of levels over which aerosols are present |
---|
| 557 | integer :: istartlev, iendlev |
---|
| 558 | |
---|
| 559 | ! Indices to spectral band and relative humidity look-up table |
---|
| 560 | integer :: iband, irh, irhs(nlev) |
---|
| 561 | |
---|
| 562 | ! Short cut for ao%itype(jtype) |
---|
| 563 | integer :: itype |
---|
| 564 | |
---|
| 565 | ! Pointer to the aerosol optics coefficients for brevity of access |
---|
| 566 | type(aerosol_optics_type), pointer :: ao |
---|
| 567 | |
---|
| 568 | real(jphook) :: hook_handle |
---|
| 569 | |
---|
| 570 | if (lhook) call dr_hook('radiation_aerosol_optics:add_aerosol_optics',0,hook_handle) |
---|
| 571 | |
---|
| 572 | if (aerosol%is_direct) then |
---|
| 573 | ! Aerosol optical properties have been provided in each band |
---|
| 574 | ! directly by the user |
---|
| 575 | call add_aerosol_optics_direct(nlev,istartcol,iendcol, & |
---|
| 576 | & config, aerosol, & |
---|
| 577 | & od_lw, ssa_lw, g_lw, od_sw, ssa_sw, g_sw) |
---|
| 578 | else |
---|
| 579 | ! Aerosol mixing ratios have been provided |
---|
| 580 | |
---|
| 581 | do jtype = 1,config%n_aerosol_types |
---|
| 582 | if (config%aerosol_optics%iclass(jtype) == IAerosolClassUndefined) then |
---|
| 583 | write(nulerr,'(a)') '*** Error: not all aerosol types are defined' |
---|
| 584 | call radiation_abort() |
---|
| 585 | end if |
---|
| 586 | end do |
---|
| 587 | |
---|
| 588 | if (config%iverbose >= 2) then |
---|
| 589 | write(nulout,'(a)') 'Computing aerosol absorption/scattering properties' |
---|
| 590 | end if |
---|
| 591 | |
---|
| 592 | ao => config%aerosol_optics |
---|
| 593 | |
---|
| 594 | istartlev = lbound(aerosol%mixing_ratio,2) |
---|
| 595 | iendlev = ubound(aerosol%mixing_ratio,2) |
---|
| 596 | |
---|
| 597 | if (ubound(aerosol%mixing_ratio,3) /= config%n_aerosol_types) then |
---|
| 598 | write(nulerr,'(a,i0,a,i0)') '*** Error: aerosol%mixing_ratio contains ', & |
---|
| 599 | & ubound(aerosol%mixing_ratio,3), ' aerosol types, expected ', & |
---|
| 600 | & config%n_aerosol_types |
---|
| 601 | call radiation_abort() |
---|
| 602 | end if |
---|
| 603 | |
---|
| 604 | ! Set variables to zero that may not have been previously |
---|
| 605 | g_sw(:,:,istartcol:iendcol) = 0.0_jprb |
---|
| 606 | if (config%do_lw_aerosol_scattering) then |
---|
| 607 | ssa_lw(:,:,istartcol:iendcol) = 0.0_jprb |
---|
| 608 | g_lw(:,:,istartcol:iendcol) = 0.0_jprb |
---|
| 609 | end if |
---|
| 610 | |
---|
| 611 | call gas%get(IH2O, IMassMixingRatio, h2o_mmr, istartcol=istartcol) |
---|
| 612 | |
---|
| 613 | ! Loop over column |
---|
| 614 | do jcol = istartcol,iendcol |
---|
| 615 | |
---|
| 616 | ! Reset temporary arrays |
---|
| 617 | od_sw_aerosol = 0.0_jprb |
---|
| 618 | scat_sw_aerosol = 0.0_jprb |
---|
| 619 | scat_g_sw_aerosol = 0.0_jprb |
---|
| 620 | od_lw_aerosol = 0.0_jprb |
---|
| 621 | scat_lw_aerosol = 0.0_jprb |
---|
| 622 | scat_g_lw_aerosol = 0.0_jprb |
---|
| 623 | |
---|
| 624 | do jlev = istartlev,iendlev |
---|
| 625 | ! Compute relative humidity with respect to liquid |
---|
| 626 | ! saturation and the index to the relative-humidity index of |
---|
| 627 | ! hydrophilic-aerosol data |
---|
| 628 | rh = h2o_mmr(jcol,jlev) / thermodynamics%h2o_sat_liq(jcol,jlev) |
---|
| 629 | irhs(jlev) = ao%calc_rh_index(rh) |
---|
| 630 | |
---|
| 631 | factor(jlev) = ( thermodynamics%pressure_hl(jcol,jlev+1) & |
---|
| 632 | & -thermodynamics%pressure_hl(jcol,jlev ) ) & |
---|
| 633 | & * OneOverAccelDueToGravity |
---|
| 634 | end do |
---|
| 635 | |
---|
| 636 | do jtype = 1,config%n_aerosol_types |
---|
| 637 | itype = ao%itype(jtype) |
---|
| 638 | |
---|
| 639 | ! Add the optical depth, scattering optical depth and |
---|
| 640 | ! scattering optical depth-weighted asymmetry factor for |
---|
| 641 | ! this aerosol type to the total for all aerosols. Note |
---|
| 642 | ! that the following expressions are array-wise, the |
---|
| 643 | ! dimension being spectral band. |
---|
| 644 | if (ao%iclass(jtype) == IAerosolClassHydrophobic) then |
---|
| 645 | do jlev = istartlev,iendlev |
---|
| 646 | mixing_ratio = aerosol%mixing_ratio(jcol,jlev,jtype) |
---|
| 647 | do jband = 1,config%n_bands_sw |
---|
| 648 | local_od_sw = factor(jlev) * mixing_ratio & |
---|
| 649 | & * ao%mass_ext_sw_phobic(jband,itype) |
---|
| 650 | od_sw_aerosol(jband,jlev) = od_sw_aerosol(jband,jlev) + local_od_sw |
---|
| 651 | scat_sw_aerosol(jband,jlev) = scat_sw_aerosol(jband,jlev) & |
---|
| 652 | & + local_od_sw * ao%ssa_sw_phobic(jband,itype) |
---|
| 653 | scat_g_sw_aerosol(jband,jlev) = scat_g_sw_aerosol(jband,jlev) & |
---|
| 654 | & + local_od_sw * ao%ssa_sw_phobic(jband,itype) & |
---|
| 655 | & * ao%g_sw_phobic(jband,itype) |
---|
| 656 | end do |
---|
| 657 | if (config%do_lw_aerosol_scattering) then |
---|
| 658 | do jband = 1,config%n_bands_lw |
---|
| 659 | local_od_lw = factor(jlev) * mixing_ratio & |
---|
| 660 | & * ao%mass_ext_lw_phobic(jband,itype) |
---|
| 661 | od_lw_aerosol(jband,jlev) = od_lw_aerosol(jband,jlev) + local_od_lw |
---|
| 662 | scat_lw_aerosol(jband,jlev) = scat_lw_aerosol(jband,jlev) & |
---|
| 663 | & + local_od_lw * ao%ssa_lw_phobic(jband,itype) |
---|
| 664 | scat_g_lw_aerosol(jband,jlev) = scat_g_lw_aerosol(jband,jlev) & |
---|
| 665 | & + local_od_lw * ao%ssa_lw_phobic(jband,itype) & |
---|
| 666 | & * ao%g_lw_phobic(jband,itype) |
---|
| 667 | end do |
---|
| 668 | else |
---|
| 669 | ! If aerosol longwave scattering is not included then we |
---|
| 670 | ! weight the optical depth by the single scattering |
---|
| 671 | ! co-albedo |
---|
| 672 | do jband = 1,config%n_bands_lw |
---|
| 673 | od_lw_aerosol(jband,jlev) = od_lw_aerosol(jband,jlev) & |
---|
| 674 | & + factor(jlev) * mixing_ratio & |
---|
| 675 | & * ao%mass_ext_lw_phobic(jband,itype) & |
---|
| 676 | & * (1.0_jprb - ao%ssa_lw_phobic(jband,itype)) |
---|
| 677 | end do |
---|
| 678 | end if |
---|
| 679 | end do |
---|
| 680 | |
---|
| 681 | else if (ao%iclass(jtype) == IAerosolClassHydrophilic) then |
---|
| 682 | ! Hydrophilic aerosols require the look-up tables to |
---|
| 683 | ! be indexed with irh |
---|
| 684 | do jlev = istartlev,iendlev |
---|
| 685 | mixing_ratio = aerosol%mixing_ratio(jcol,jlev,jtype) |
---|
| 686 | irh = irhs(jlev) |
---|
| 687 | do jband = 1,config%n_bands_sw |
---|
| 688 | local_od_sw = factor(jlev) * mixing_ratio & |
---|
| 689 | & * ao%mass_ext_sw_philic(jband,irh,itype) |
---|
| 690 | od_sw_aerosol(jband,jlev) = od_sw_aerosol(jband,jlev) + local_od_sw |
---|
| 691 | scat_sw_aerosol(jband,jlev) = scat_sw_aerosol(jband,jlev) & |
---|
| 692 | & + local_od_sw * ao%ssa_sw_philic(jband,irh,itype) |
---|
| 693 | scat_g_sw_aerosol(jband,jlev) = scat_g_sw_aerosol(jband,jlev) & |
---|
| 694 | & + local_od_sw * ao%ssa_sw_philic(jband,irh,itype) & |
---|
| 695 | & * ao%g_sw_philic(jband,irh,itype) |
---|
| 696 | end do |
---|
| 697 | if (config%do_lw_aerosol_scattering) then |
---|
| 698 | do jband = 1,config%n_bands_lw |
---|
| 699 | local_od_lw = factor(jlev) * mixing_ratio & |
---|
| 700 | & * ao%mass_ext_lw_philic(jband,irh,itype) |
---|
| 701 | od_lw_aerosol(jband,jlev) = od_lw_aerosol(jband,jlev) + local_od_lw |
---|
| 702 | scat_lw_aerosol(jband,jlev) = scat_lw_aerosol(jband,jlev) & |
---|
| 703 | & + local_od_lw * ao%ssa_lw_philic(jband,irh,itype) |
---|
| 704 | scat_g_lw_aerosol(jband,jlev) = scat_g_lw_aerosol(jband,jlev) & |
---|
| 705 | & + local_od_lw * ao%ssa_lw_philic(jband,irh,itype) & |
---|
| 706 | & * ao%g_lw_philic(jband,irh,itype) |
---|
| 707 | end do |
---|
| 708 | else |
---|
| 709 | ! If aerosol longwave scattering is not included then we |
---|
| 710 | ! weight the optical depth by the single scattering |
---|
| 711 | ! co-albedo |
---|
| 712 | do jband = 1,config%n_bands_lw |
---|
| 713 | od_lw_aerosol(jband,jlev) = od_lw_aerosol(jband,jlev) & |
---|
| 714 | & + factor(jlev) * mixing_ratio & |
---|
| 715 | & * ao%mass_ext_lw_philic(jband,irh,itype) & |
---|
| 716 | & * (1.0_jprb - ao%ssa_lw_philic(jband,irh,itype)) |
---|
| 717 | end do |
---|
| 718 | end if |
---|
| 719 | end do |
---|
| 720 | |
---|
| 721 | ! Implicitly, if ao%iclass(jtype) == IAerosolClassNone, then |
---|
| 722 | ! no aerosol scattering properties are added |
---|
| 723 | end if |
---|
| 724 | |
---|
| 725 | end do ! Loop over aerosol type |
---|
| 726 | |
---|
| 727 | if (.not. config%do_sw_delta_scaling_with_gases) then |
---|
| 728 | ! Delta-Eddington scaling on aerosol only. Note that if |
---|
| 729 | ! do_sw_delta_scaling_with_gases==.true. then the delta |
---|
| 730 | ! scaling is done to the cloud-aerosol-gas mixture inside |
---|
| 731 | ! the solver |
---|
| 732 | call delta_eddington_extensive_vec(config%n_bands_sw*nlev, od_sw_aerosol, & |
---|
| 733 | & scat_sw_aerosol, scat_g_sw_aerosol) |
---|
| 734 | end if |
---|
| 735 | |
---|
| 736 | ! Combine aerosol shortwave scattering properties with gas |
---|
| 737 | ! properties (noting that any gas scattering will have an |
---|
| 738 | ! asymmetry factor of zero) |
---|
| 739 | if (config%do_cloud_aerosol_per_sw_g_point) then |
---|
| 740 | |
---|
| 741 | ! We can assume the band and g-point indices are the same |
---|
| 742 | do jlev = 1,nlev |
---|
| 743 | do jg = 1,config%n_g_sw |
---|
| 744 | local_scat = ssa_sw(jg,jlev,jcol)*od_sw(jg,jlev,jcol) + scat_sw_aerosol(jg,jlev) |
---|
| 745 | od_sw(jg,jlev,jcol) = od_sw(jg,jlev,jcol) + od_sw_aerosol(jg,jlev) |
---|
| 746 | g_sw(jg,jlev,jcol) = scat_g_sw_aerosol(jg,jlev) / max(local_scat, 1.0e-24_jprb) |
---|
| 747 | ssa_sw(jg,jlev,jcol) = min(local_scat / max(od_sw(jg,jlev,jcol), 1.0e-24_jprb), 1.0_jprb) |
---|
| 748 | end do |
---|
| 749 | end do |
---|
| 750 | |
---|
| 751 | else |
---|
| 752 | |
---|
| 753 | do jlev = 1,nlev |
---|
| 754 | do jg = 1,config%n_g_sw |
---|
| 755 | ! Need to map between bands and g-points |
---|
| 756 | iband = config%i_band_from_reordered_g_sw(jg) |
---|
| 757 | local_od = od_sw(jg,jlev,jcol) + od_sw_aerosol(iband,jlev) |
---|
| 758 | if (local_od > 0.0_jprb .and. od_sw_aerosol(iband,jlev) > 0.0_jprb) then |
---|
| 759 | local_scat = ssa_sw(jg,jlev,jcol) * od_sw(jg,jlev,jcol) & |
---|
| 760 | & + scat_sw_aerosol(iband,jlev) |
---|
| 761 | ! Note that asymmetry_sw of gases is zero so the following |
---|
| 762 | ! simply weights the aerosol asymmetry by the scattering |
---|
| 763 | ! optical depth |
---|
| 764 | if (local_scat > 0.0_jprb) then |
---|
| 765 | g_sw(jg,jlev,jcol) = scat_g_sw_aerosol(iband,jlev) / local_scat |
---|
| 766 | end if |
---|
| 767 | ssa_sw(jg,jlev,jcol) = local_scat / local_od |
---|
| 768 | od_sw (jg,jlev,jcol) = local_od |
---|
| 769 | end if |
---|
| 770 | end do |
---|
| 771 | end do |
---|
| 772 | |
---|
| 773 | end if |
---|
| 774 | |
---|
| 775 | ! Combine aerosol longwave scattering properties with gas |
---|
| 776 | ! properties, noting that in the longwave, gases do not |
---|
| 777 | ! scatter at all |
---|
| 778 | if (config%do_lw_aerosol_scattering) then |
---|
| 779 | |
---|
| 780 | call delta_eddington_extensive_vec(config%n_bands_lw*nlev, od_lw_aerosol, & |
---|
| 781 | & scat_lw_aerosol, scat_g_lw_aerosol) |
---|
| 782 | |
---|
| 783 | do jlev = istartlev,iendlev |
---|
| 784 | do jg = 1,config%n_g_lw |
---|
| 785 | iband = config%i_band_from_reordered_g_lw(jg) |
---|
| 786 | local_od = od_lw(jg,jlev,jcol) + od_lw_aerosol(iband,jlev) |
---|
| 787 | if (local_od > 0.0_jprb .and. od_lw_aerosol(iband,jlev) > 0.0_jprb) then |
---|
| 788 | ! All scattering is due to aerosols, therefore the |
---|
| 789 | ! asymmetry factor is equal to the value for aerosols |
---|
| 790 | if (scat_lw_aerosol(iband,jlev) > 0.0_jprb) then |
---|
| 791 | g_lw(jg,jlev,jcol) = scat_g_lw_aerosol(iband,jlev) & |
---|
| 792 | & / scat_lw_aerosol(iband,jlev) |
---|
| 793 | end if |
---|
| 794 | ssa_lw(jg,jlev,jcol) = scat_lw_aerosol(iband,jlev) / local_od |
---|
| 795 | od_lw (jg,jlev,jcol) = local_od |
---|
| 796 | end if |
---|
| 797 | end do |
---|
| 798 | end do |
---|
| 799 | |
---|
| 800 | else |
---|
| 801 | |
---|
| 802 | if (config%do_cloud_aerosol_per_lw_g_point) then |
---|
| 803 | ! We can assume band and g-point indices are the same |
---|
| 804 | do jlev = istartlev,iendlev |
---|
| 805 | do jg = 1,config%n_g_lw |
---|
| 806 | od_lw(jg,jlev,jcol) = od_lw(jg,jlev,jcol) + od_lw_aerosol(jg,jlev) |
---|
| 807 | end do |
---|
| 808 | end do |
---|
| 809 | else |
---|
| 810 | do jlev = istartlev,iendlev |
---|
| 811 | do jg = 1,config%n_g_lw |
---|
| 812 | od_lw(jg,jlev,jcol) = od_lw(jg,jlev,jcol) & |
---|
| 813 | & + od_lw_aerosol(config%i_band_from_reordered_g_lw(jg),jlev) |
---|
| 814 | end do |
---|
| 815 | end do |
---|
| 816 | end if |
---|
| 817 | |
---|
| 818 | end if |
---|
| 819 | |
---|
| 820 | end do ! Loop over column |
---|
| 821 | |
---|
| 822 | end if |
---|
| 823 | |
---|
| 824 | if (lhook) call dr_hook('radiation_aerosol_optics:add_aerosol_optics',1,hook_handle) |
---|
| 825 | |
---|
| 826 | end subroutine add_aerosol_optics |
---|
| 827 | |
---|
| 828 | |
---|
| 829 | !--------------------------------------------------------------------- |
---|
| 830 | ! Add precomputed optical properties to gas optical depth and |
---|
| 831 | ! scattering properties |
---|
| 832 | subroutine add_aerosol_optics_direct(nlev,istartcol,iendcol, & |
---|
| 833 | & config, aerosol, & |
---|
| 834 | & od_lw, ssa_lw, g_lw, od_sw, ssa_sw, g_sw) |
---|
| 835 | |
---|
| 836 | use parkind1, only : jprb |
---|
| 837 | use radiation_io, only : nulerr, radiation_abort |
---|
| 838 | use yomhook, only : lhook, dr_hook, jphook |
---|
| 839 | use radiation_config, only : config_type |
---|
| 840 | use radiation_aerosol, only : aerosol_type |
---|
| 841 | |
---|
| 842 | integer, intent(in) :: nlev ! number of model levels |
---|
| 843 | integer, intent(in) :: istartcol, iendcol ! range of columns to process |
---|
| 844 | type(config_type), intent(in), target :: config |
---|
| 845 | type(aerosol_type), intent(in) :: aerosol |
---|
| 846 | ! Optical depth, single scattering albedo and asymmetry factor of |
---|
| 847 | ! the atmosphere (gases on input, gases and aerosols on output) |
---|
| 848 | ! for each g point. Note that longwave ssa and asymmetry and |
---|
| 849 | ! shortwave asymmetry are all zero for gases, so are not yet |
---|
| 850 | ! defined on input and are therefore intent(out). |
---|
| 851 | real(jprb), dimension(config%n_g_lw,nlev,istartcol:iendcol), & |
---|
| 852 | & intent(inout) :: od_lw |
---|
| 853 | real(jprb), dimension(config%n_g_lw_if_scattering,nlev,istartcol:iendcol), & |
---|
| 854 | & intent(out) :: ssa_lw, g_lw |
---|
| 855 | real(jprb), dimension(config%n_g_sw,nlev,istartcol:iendcol), & |
---|
| 856 | & intent(inout) :: od_sw, ssa_sw |
---|
| 857 | real(jprb), dimension(config%n_g_sw,nlev,istartcol:iendcol), & |
---|
| 858 | & intent(out) :: g_sw |
---|
| 859 | |
---|
| 860 | ! Temporary extinction and scattering optical depths of aerosol |
---|
| 861 | ! plus gas |
---|
| 862 | real(jprb) :: local_od, local_scat |
---|
| 863 | |
---|
| 864 | ! Extinction optical depth, scattering optical depth and |
---|
| 865 | ! asymmetry-times-scattering-optical-depth for all the aerosols at |
---|
| 866 | ! a point in space for each spectral band of the shortwave and |
---|
| 867 | ! longwave spectrum |
---|
| 868 | real(jprb), dimension(config%n_bands_sw,nlev) & |
---|
| 869 | & :: od_sw_aerosol, scat_sw_aerosol, scat_g_sw_aerosol |
---|
| 870 | real(jprb), dimension(config%n_bands_lw,nlev) :: od_lw_aerosol |
---|
| 871 | real(jprb), dimension(config%n_bands_lw_if_scattering,nlev) & |
---|
| 872 | & :: scat_lw_aerosol, scat_g_lw_aerosol |
---|
| 873 | |
---|
| 874 | ! Loop indices for column, level, g point and band |
---|
| 875 | integer :: jcol, jlev, jg, jb |
---|
| 876 | |
---|
| 877 | ! Range of levels over which aerosols are present |
---|
| 878 | integer :: istartlev, iendlev |
---|
| 879 | |
---|
| 880 | ! Indices to spectral band |
---|
| 881 | integer :: iband |
---|
| 882 | |
---|
| 883 | real(jphook) :: hook_handle |
---|
| 884 | |
---|
| 885 | if (lhook) call dr_hook('radiation_aerosol_optics:add_aerosol_optics_direct',0,hook_handle) |
---|
| 886 | |
---|
| 887 | if (config%do_sw) then |
---|
| 888 | ! Check array dimensions |
---|
| 889 | if (ubound(aerosol%od_sw,1) /= config%n_bands_sw) then |
---|
| 890 | write(nulerr,'(a,i0,a,i0)') '*** Error: aerosol%od_sw contains ', & |
---|
| 891 | & ubound(aerosol%od_sw,1), ' band, expected ', & |
---|
| 892 | & config%n_bands_sw |
---|
| 893 | call radiation_abort() |
---|
| 894 | end if |
---|
| 895 | |
---|
| 896 | istartlev = lbound(aerosol%od_sw,2) |
---|
| 897 | iendlev = ubound(aerosol%od_sw,2) |
---|
| 898 | |
---|
| 899 | ! Set variables to zero that may not have been previously |
---|
| 900 | g_sw(:,:,istartcol:iendcol) = 0.0_jprb |
---|
| 901 | |
---|
| 902 | ! Loop over position |
---|
| 903 | do jcol = istartcol,iendcol |
---|
| 904 | ! Added for DWD (2020) |
---|
| 905 | !NEC$ forced_collapse |
---|
| 906 | do jlev = istartlev,iendlev |
---|
| 907 | do jb = 1,config%n_bands_sw |
---|
| 908 | od_sw_aerosol(jb,jlev) = aerosol%od_sw(jb,jlev,jcol) |
---|
| 909 | scat_sw_aerosol(jb,jlev) = aerosol%ssa_sw(jb,jlev,jcol) * od_sw_aerosol(jb,jlev) |
---|
| 910 | scat_g_sw_aerosol(jb,jlev) = aerosol%g_sw(jb,jlev,jcol) * scat_sw_aerosol(jb,jlev) |
---|
| 911 | |
---|
| 912 | if (.not. config%do_sw_delta_scaling_with_gases) then |
---|
| 913 | ! Delta-Eddington scaling on aerosol only. Note that if |
---|
| 914 | ! do_sw_delta_scaling_with_gases==.true. then the delta |
---|
| 915 | ! scaling is done to the cloud-aerosol-gas mixture |
---|
| 916 | ! inside the solver |
---|
| 917 | call delta_eddington_extensive(od_sw_aerosol(jb,jlev), scat_sw_aerosol(jb,jlev), & |
---|
| 918 | & scat_g_sw_aerosol(jb,jlev)) |
---|
| 919 | end if |
---|
| 920 | end do |
---|
| 921 | end do |
---|
| 922 | ! Combine aerosol shortwave scattering properties with gas |
---|
| 923 | ! properties (noting that any gas scattering will have an |
---|
| 924 | ! asymmetry factor of zero) |
---|
| 925 | do jlev = istartlev,iendlev |
---|
| 926 | if (od_sw_aerosol(1,jlev) > 0.0_jprb) then |
---|
| 927 | do jg = 1,config%n_g_sw |
---|
| 928 | iband = config%i_band_from_reordered_g_sw(jg) |
---|
| 929 | local_od = od_sw(jg,jlev,jcol) + od_sw_aerosol(iband,jlev) |
---|
| 930 | local_scat = ssa_sw(jg,jlev,jcol) * od_sw(jg,jlev,jcol) & |
---|
| 931 | & + scat_sw_aerosol(iband,jlev) |
---|
| 932 | ! Note that asymmetry_sw of gases is zero so the following |
---|
| 933 | ! simply weights the aerosol asymmetry by the scattering |
---|
| 934 | ! optical depth |
---|
| 935 | g_sw(jg,jlev,jcol) = scat_g_sw_aerosol(iband,jlev) / local_scat |
---|
| 936 | local_od = od_sw(jg,jlev,jcol) + od_sw_aerosol(iband,jlev) |
---|
| 937 | ssa_sw(jg,jlev,jcol) = local_scat / local_od |
---|
| 938 | od_sw (jg,jlev,jcol) = local_od |
---|
| 939 | end do |
---|
| 940 | end if |
---|
| 941 | end do |
---|
| 942 | end do |
---|
| 943 | |
---|
| 944 | end if |
---|
| 945 | |
---|
| 946 | if (config%do_lw) then |
---|
| 947 | |
---|
| 948 | if (ubound(aerosol%od_lw,1) /= config%n_bands_lw) then |
---|
| 949 | write(nulerr,'(a,i0,a,i0)') '*** Error: aerosol%od_lw contains ', & |
---|
| 950 | & ubound(aerosol%od_lw,1), ' band, expected ', & |
---|
| 951 | & config%n_bands_lw |
---|
| 952 | call radiation_abort() |
---|
| 953 | end if |
---|
| 954 | |
---|
| 955 | istartlev = lbound(aerosol%od_lw,2) |
---|
| 956 | iendlev = ubound(aerosol%od_lw,2) |
---|
| 957 | |
---|
| 958 | if (config%do_lw_aerosol_scattering) then |
---|
| 959 | ssa_lw(:,:,istartcol:iendcol) = 0.0_jprb |
---|
| 960 | g_lw(:,:,istartcol:iendcol) = 0.0_jprb |
---|
| 961 | |
---|
| 962 | ! Loop over position |
---|
| 963 | do jcol = istartcol,iendcol |
---|
| 964 | ! Added for DWD (2020) |
---|
| 965 | !NEC$ forced_collapse |
---|
| 966 | do jlev = istartlev,iendlev |
---|
| 967 | do jb = 1,config%n_bands_lw |
---|
| 968 | od_lw_aerosol(jb,jlev) = aerosol%od_lw(jb,jlev,jcol) |
---|
| 969 | scat_lw_aerosol(jb,jlev) = aerosol%ssa_lw(jb,jlev,jcol) * od_lw_aerosol(jb,jlev) |
---|
| 970 | scat_g_lw_aerosol(jb,jlev) = aerosol%g_lw(jb,jlev,jcol) * scat_lw_aerosol(jb,jlev) |
---|
| 971 | |
---|
| 972 | call delta_eddington_extensive(od_lw_aerosol(jb,jlev), scat_lw_aerosol(jb,jlev), & |
---|
| 973 | & scat_g_lw_aerosol(jb,jlev)) |
---|
| 974 | end do |
---|
| 975 | end do |
---|
| 976 | do jlev = istartlev,iendlev |
---|
| 977 | do jg = 1,config%n_g_lw |
---|
| 978 | iband = config%i_band_from_reordered_g_lw(jg) |
---|
| 979 | if (od_lw_aerosol(iband,jlev) > 0.0_jprb) then |
---|
| 980 | ! All scattering is due to aerosols, therefore the |
---|
| 981 | ! asymmetry factor is equal to the value for aerosols |
---|
| 982 | if (scat_lw_aerosol(iband,jlev) > 0.0_jprb) then |
---|
| 983 | g_lw(jg,jlev,jcol) = scat_g_lw_aerosol(iband,jlev) & |
---|
| 984 | & / scat_lw_aerosol(iband,jlev) |
---|
| 985 | end if |
---|
| 986 | local_od = od_lw(jg,jlev,jcol) + od_lw_aerosol(iband,jlev) |
---|
| 987 | ssa_lw(jg,jlev,jcol) = scat_lw_aerosol(iband,jlev) / local_od |
---|
| 988 | od_lw (jg,jlev,jcol) = local_od |
---|
| 989 | end if |
---|
| 990 | end do |
---|
| 991 | end do |
---|
| 992 | end do |
---|
| 993 | |
---|
| 994 | else ! No longwave scattering |
---|
| 995 | |
---|
| 996 | ! Loop over position |
---|
| 997 | do jcol = istartcol,iendcol |
---|
| 998 | ! Added for DWD (2020) |
---|
| 999 | !NEC$ forced_collapse |
---|
| 1000 | do jlev = istartlev,iendlev |
---|
| 1001 | ! If aerosol longwave scattering is not included then we |
---|
| 1002 | ! weight the optical depth by the single scattering |
---|
| 1003 | ! co-albedo |
---|
| 1004 | do jb = 1, config%n_bands_lw |
---|
| 1005 | od_lw_aerosol(jb,jlev) = aerosol%od_lw(jb,jlev,jcol) & |
---|
| 1006 | & * (1.0_jprb - aerosol%ssa_lw(jb,jlev,jcol)) |
---|
| 1007 | end do |
---|
| 1008 | end do |
---|
| 1009 | do jlev = istartlev,iendlev |
---|
| 1010 | do jg = 1,config%n_g_lw |
---|
| 1011 | od_lw(jg,jlev,jcol) = od_lw(jg,jlev,jcol) & |
---|
| 1012 | & + od_lw_aerosol(config%i_band_from_reordered_g_lw(jg),jlev) |
---|
| 1013 | end do |
---|
| 1014 | end do |
---|
| 1015 | end do |
---|
| 1016 | |
---|
| 1017 | end if |
---|
| 1018 | end if |
---|
| 1019 | |
---|
| 1020 | |
---|
| 1021 | if (lhook) call dr_hook('radiation_aerosol_optics:add_aerosol_optics_direct',1,hook_handle) |
---|
| 1022 | |
---|
| 1023 | end subroutine add_aerosol_optics_direct |
---|
| 1024 | |
---|
| 1025 | |
---|
| 1026 | !--------------------------------------------------------------------- |
---|
| 1027 | ! Sometimes it is useful to specify aerosol in terms of its optical |
---|
| 1028 | ! depth at a particular wavelength. This function returns the dry |
---|
| 1029 | ! mass-extinction coefficient, i.e. the extinction cross section per |
---|
| 1030 | ! unit mass, for aerosol of type "itype" at the specified wavelength |
---|
| 1031 | ! (m). For hydrophilic types, the value at the first relative |
---|
| 1032 | ! humidity bin is taken. |
---|
| 1033 | function dry_aerosol_mass_extinction(config, itype, wavelength) |
---|
| 1034 | |
---|
| 1035 | use parkind1, only : jprb |
---|
| 1036 | use radiation_io, only : nulerr, radiation_abort |
---|
| 1037 | use radiation_config, only : config_type |
---|
| 1038 | use radiation_aerosol_optics_data, only : aerosol_optics_type, & |
---|
| 1039 | & IAerosolClassUndefined, IAerosolClassIgnored, & |
---|
| 1040 | & IAerosolClassHydrophobic, IAerosolClassHydrophilic |
---|
| 1041 | |
---|
| 1042 | type(config_type), intent(in), target :: config |
---|
| 1043 | |
---|
| 1044 | ! Aerosol type |
---|
| 1045 | integer, intent(in) :: itype |
---|
| 1046 | |
---|
| 1047 | ! Wavelength (m) |
---|
| 1048 | real(jprb), intent(in) :: wavelength |
---|
| 1049 | |
---|
| 1050 | real(jprb) :: dry_aerosol_mass_extinction |
---|
| 1051 | |
---|
| 1052 | ! Index to the monochromatic wavelength requested |
---|
| 1053 | integer :: imono |
---|
| 1054 | |
---|
| 1055 | ! Pointer to the aerosol optics coefficients for brevity of access |
---|
| 1056 | type(aerosol_optics_type), pointer :: ao |
---|
| 1057 | |
---|
| 1058 | ao => config%aerosol_optics |
---|
| 1059 | |
---|
| 1060 | imono = minloc(abs(wavelength - ao%wavelength_mono), 1) |
---|
| 1061 | |
---|
| 1062 | if (abs(wavelength - ao%wavelength_mono(imono))/wavelength > 0.02_jprb) then |
---|
| 1063 | write(nulerr,'(a,e11.4,a)') '*** Error: requested wavelength ', & |
---|
| 1064 | & wavelength, ' not within 2% of stored wavelengths' |
---|
| 1065 | call radiation_abort() |
---|
| 1066 | end if |
---|
| 1067 | |
---|
| 1068 | if (ao%iclass(itype) == IAerosolClassHydrophobic) then |
---|
| 1069 | dry_aerosol_mass_extinction = ao%mass_ext_mono_phobic(imono,ao%itype(itype)) |
---|
| 1070 | else if (ao%iclass(itype) == IAerosolClassHydrophilic) then |
---|
| 1071 | ! Take the value at the first relative-humidity bin for the |
---|
| 1072 | ! "dry" aerosol value |
---|
| 1073 | dry_aerosol_mass_extinction = ao%mass_ext_mono_philic(imono,1,ao%itype(itype)) |
---|
| 1074 | else |
---|
| 1075 | dry_aerosol_mass_extinction = 0.0_jprb |
---|
| 1076 | end if |
---|
| 1077 | |
---|
| 1078 | end function dry_aerosol_mass_extinction |
---|
| 1079 | |
---|
| 1080 | |
---|
| 1081 | !--------------------------------------------------------------------- |
---|
| 1082 | ! Compute aerosol extinction coefficient at a particular wavelength |
---|
| 1083 | ! and a single height - this is useful for visibility diagnostics |
---|
| 1084 | subroutine aerosol_extinction(ncol,istartcol,iendcol, & |
---|
| 1085 | & config, wavelength, mixing_ratio, relative_humidity, extinction) |
---|
| 1086 | |
---|
| 1087 | use parkind1, only : jprb |
---|
| 1088 | use yomhook, only : lhook, dr_hook, jphook |
---|
| 1089 | use radiation_io, only : nulerr, radiation_abort |
---|
| 1090 | use radiation_config, only : config_type |
---|
| 1091 | use radiation_aerosol_optics_data, only : aerosol_optics_type, & |
---|
| 1092 | & IAerosolClassUndefined, IAerosolClassIgnored, & |
---|
| 1093 | & IAerosolClassHydrophobic, IAerosolClassHydrophilic |
---|
| 1094 | |
---|
| 1095 | integer, intent(in) :: ncol ! number of columns |
---|
| 1096 | integer, intent(in) :: istartcol, iendcol ! range of columns to process |
---|
| 1097 | type(config_type), intent(in), target :: config |
---|
| 1098 | real(jprb), intent(in) :: wavelength ! Requested wavelength (m) |
---|
| 1099 | real(jprb), intent(in) :: mixing_ratio(ncol,config%n_aerosol_types) |
---|
| 1100 | real(jprb), intent(in) :: relative_humidity(ncol) |
---|
| 1101 | real(jprb), intent(out) :: extinction(ncol) |
---|
| 1102 | |
---|
| 1103 | ! Local aerosol extinction |
---|
| 1104 | real(jprb) :: ext |
---|
| 1105 | |
---|
| 1106 | ! Index to the monochromatic wavelength requested |
---|
| 1107 | integer :: imono |
---|
| 1108 | |
---|
| 1109 | ! Pointer to the aerosol optics coefficients for brevity of access |
---|
| 1110 | type(aerosol_optics_type), pointer :: ao |
---|
| 1111 | |
---|
| 1112 | ! Loop indices for column and aerosol type |
---|
| 1113 | integer :: jcol, jtype |
---|
| 1114 | |
---|
| 1115 | ! Relative humidity index |
---|
| 1116 | integer :: irh |
---|
| 1117 | |
---|
| 1118 | real(jphook) :: hook_handle |
---|
| 1119 | |
---|
| 1120 | if (lhook) call dr_hook('radiation_aerosol_optics:aerosol_extinction',0,hook_handle) |
---|
| 1121 | |
---|
| 1122 | do jtype = 1,config%n_aerosol_types |
---|
| 1123 | if (config%aerosol_optics%iclass(jtype) == IAerosolClassUndefined) then |
---|
| 1124 | write(nulerr,'(a)') '*** Error: not all aerosol types are defined' |
---|
| 1125 | call radiation_abort() |
---|
| 1126 | end if |
---|
| 1127 | end do |
---|
| 1128 | |
---|
| 1129 | ao => config%aerosol_optics |
---|
| 1130 | |
---|
| 1131 | imono = minloc(abs(wavelength - ao%wavelength_mono), 1) |
---|
| 1132 | |
---|
| 1133 | if (abs(wavelength - ao%wavelength_mono(imono))/wavelength > 0.02_jprb) then |
---|
| 1134 | write(nulerr,'(a,e11.4,a)') '*** Error: requested wavelength ', & |
---|
| 1135 | & wavelength, ' not within 2% of stored wavelengths' |
---|
| 1136 | call radiation_abort() |
---|
| 1137 | end if |
---|
| 1138 | |
---|
| 1139 | ! Loop over position |
---|
| 1140 | do jcol = istartcol,iendcol |
---|
| 1141 | ext = 0.0_jprb |
---|
| 1142 | ! Get relative-humidity index |
---|
| 1143 | irh = ao%calc_rh_index(relative_humidity(jcol)) |
---|
| 1144 | ! Add extinction coefficients from each aerosol type |
---|
| 1145 | do jtype = 1,config%n_aerosol_types |
---|
| 1146 | if (ao%iclass(jtype) == IAerosolClassHydrophobic) then |
---|
| 1147 | ext = ext + mixing_ratio(jcol,jtype) & |
---|
| 1148 | & * ao%mass_ext_mono_phobic(imono,ao%itype(jtype)) |
---|
| 1149 | else if (ao%iclass(jtype) == IAerosolClassHydrophilic) then |
---|
| 1150 | ext = ext + mixing_ratio(jcol,jtype) & |
---|
| 1151 | & * ao%mass_ext_mono_philic(imono,irh,ao%itype(jtype)) |
---|
| 1152 | end if |
---|
| 1153 | end do |
---|
| 1154 | |
---|
| 1155 | extinction(jcol) = ext |
---|
| 1156 | end do |
---|
| 1157 | |
---|
| 1158 | if (lhook) call dr_hook('radiation_aerosol_optics:aerosol_extinction',1,hook_handle) |
---|
| 1159 | |
---|
| 1160 | end subroutine aerosol_extinction |
---|
| 1161 | |
---|
| 1162 | end module radiation_aerosol_optics |
---|