Changeset 4518 for LMDZ6/trunk/libf/phylmd/ecrad
- Timestamp:
- Apr 24, 2023, 5:10:27 PM (19 months ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ6/trunk/libf/phylmd/ecrad/radiation_aerosol_optics.F90
r4517 r4518 53 53 ! properties and average to the spectral intervals of the 54 54 ! current gas-optics scheme 55 ! call setup_general_aerosol_optics(config) 56 call setup_general_aerosol_optics_lmdz(config,trim(config%aerosol_optics_file_name)) 55 call setup_general_aerosol_optics(config) 57 56 else 58 57 ! Read file containing optical properties already in the bands … … 334 333 335 334 end subroutine setup_general_aerosol_optics 336 337 !---------------------------------------------------------------------338 ! Read LMDZ file containing high spectral resolution optical properties339 ! and average to the spectral intervals of the current gas-optics340 ! scheme341 subroutine setup_general_aerosol_optics_lmdz(config,file_name)342 343 use parkind1, only : jprb344 use yomhook, only : lhook, dr_hook345 ! use easy_netcdf, only : netcdf_file346 use radiation_config, only : config_type347 use radiation_aerosol_optics_data, only : aerosol_optics_type348 use radiation_spectral_definition, only : SolarReferenceTemperature, &349 & TerrestrialReferenceTemperature350 use radiation_io, only : nulout351 use netcdf95, only: nf95_open, nf95_inq_grp_full_ncid, nf95_close, &352 nf95_inq_dimid, nf95_inq_varid, nf95_inquire_dimension, &353 nf95_get_var, nf95_gw_var354 use netcdf, only: nf90_nowrite355 356 357 type(config_type), intent(inout), target :: config358 359 ! ! The NetCDF file containing the aerosol optics data360 ! type(netcdf_file) :: file361 362 character(len=*), intent(in):: file_name363 ! NetCDF file containing the aerosol optics data364 365 ! Wavenumber points in NetCDF file366 real(jprb), allocatable :: wavenumber(:) ! cm-1367 368 ! Hydrophilic aerosol properties369 real(jprb), allocatable :: mass_ext_philic(:,:,:) ! Mass-ext coefficient (m2 kg-1)370 real(jprb), allocatable :: ssa_philic(:,:,:) ! Single-scattering albedo371 real(jprb), allocatable :: g_philic(:,:,:) ! Asymmetry factor372 real(jprb), allocatable :: lidar_ratio_philic(:,:,:) ! Lidar ratio (sr)373 374 ! Hydrophobic aerosol properties375 real(jprb), allocatable :: mass_ext_phobic(:,:) ! Mass-ext coefficient (m2 kg-1)376 real(jprb), allocatable :: ssa_phobic(:,:) ! Single-scattering albedo377 real(jprb), allocatable :: g_phobic(:,:) ! Asymmetry factor378 real(jprb), allocatable :: lidar_ratio_phobic(:,:) ! Lidar ratio (sr)379 380 ! Mapping matrix between optical properties at the wavenumbers in381 ! the file, and spectral intervals used by the gas-optics scheme382 real(jprb), allocatable :: mapping(:,:)383 384 ! Pointer to the aerosol optics coefficients for brevity of access385 type(aerosol_optics_type), pointer :: ao386 387 ! Target monochromatic wavenumber for interpolation (cm-1)388 real(jprb) :: wavenumber_target389 390 ! Number of spectral points describing aerosol properties in the391 ! shortwave and longwave392 integer :: nspecsw, nspeclw393 394 ! Number of monochromatic wavelengths required395 integer :: nmono396 397 integer :: n_type_philic, n_type_phobic, nrh, nwn398 integer :: jtype, jwl, iwn399 400 ! Weight of first point in interpolation401 real(jprb) :: weight1402 403 real(jprb) :: hook_handle404 405 ! Local:406 integer ncid, grpid, dimid, varid407 408 if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics',0,hook_handle)409 410 ao => config%aerosol_optics411 412 ao%use_hydrophilic = .true.413 ao%use_monochromatic = .true.414 print*,'file_name= ',file_name415 call nf95_open(file_name, nf90_nowrite, ncid)416 call nf95_inq_grp_full_ncid(ncid, "Hydrophilic", grpid)417 call nf95_inq_dimid(grpid, "hur", dimid)418 call nf95_inquire_dimension(grpid, dimid, nclen = ao%nrh)419 ! allocate(ao%rh_lower(ao%nrh))420 call nf95_inq_varid(grpid, "hur_bounds", varid)421 call nf95_get_var(grpid, varid, ao%rh_lower, count_nc = [1, ao%nrh])422 423 ! Hydrophilic/LW_bands:424 call nf95_inq_grp_full_ncid(ncid, "Hydrophilic/LW_bands", grpid)425 call nf95_inq_varid(grpid, "asymmetry", varid)426 call nf95_gw_var(grpid, varid, ao%g_lw_philic)427 call nf95_inq_varid(grpid, "single_scat_alb", varid)428 call nf95_gw_var(grpid, varid, ao%ssa_lw_philic)429 call nf95_inq_varid(grpid, "mass_ext", varid)430 call nf95_gw_var(grpid, varid, ao%mass_ext_lw_philic)431 432 ! Hydrophilic/SW_bands:433 call nf95_inq_grp_full_ncid(ncid, "Hydrophilic/SW_bands", grpid)434 call nf95_inq_varid(grpid, "asymmetry", varid)435 call nf95_gw_var(grpid, varid, ao%g_sw_philic)436 ao%g_sw_philic = cshift(ao%g_sw_philic, 1)437 call nf95_inq_varid(grpid, "single_scat_alb", varid)438 call nf95_gw_var(grpid, varid, ao%ssa_sw_philic)439 ao%ssa_sw_philic = cshift(ao%ssa_sw_philic, 1)440 call nf95_inq_varid(grpid, "mass_ext", varid)441 call nf95_gw_var(grpid, varid, ao%mass_ext_sw_philic)442 ao%mass_ext_sw_philic = cshift(ao%mass_ext_sw_philic, 1)443 444 ! Hydrophilic/Monochromatic:445 call nf95_inq_grp_full_ncid(ncid, "Hydrophilic/Monochromatic", grpid)446 call nf95_inq_varid(grpid, "mass_ext", varid)447 call nf95_gw_var(grpid, varid, ao%mass_ext_mono_philic)448 449 ! Hydrophobic/LW_bands:450 call nf95_inq_grp_full_ncid(ncid, "Hydrophobic/LW_bands", grpid)451 call nf95_inq_varid(grpid, "asymmetry", varid)452 call nf95_gw_var(grpid, varid, ao%g_lw_phobic)453 call nf95_inq_varid(grpid, "single_scat_alb", varid)454 call nf95_gw_var(grpid, varid, ao%ssa_lw_phobic)455 call nf95_inq_varid(grpid, "mass_ext", varid)456 call nf95_gw_var(grpid, varid, ao%mass_ext_lw_phobic)457 458 ! Hydrophobic/SW_bands:459 call nf95_inq_grp_full_ncid(ncid, "Hydrophobic/SW_bands", grpid)460 call nf95_inq_varid(grpid, "asymmetry", varid)461 call nf95_gw_var(grpid, varid, ao%g_sw_phobic)462 ao%g_sw_phobic = cshift(ao%g_sw_phobic, 1)463 call nf95_inq_varid(grpid, "single_scat_alb", varid)464 call nf95_gw_var(grpid, varid, ao%ssa_sw_phobic)465 ao%ssa_sw_phobic = cshift(ao%ssa_sw_phobic, 1)466 call nf95_inq_varid(grpid, "mass_ext", varid)467 call nf95_gw_var(grpid, varid, ao%mass_ext_sw_phobic)468 ao%mass_ext_sw_phobic = cshift(ao%mass_ext_sw_phobic, 1)469 ! AI ATTENTION470 call nf95_inq_varid(grpid, "wavenumber", varid)471 call nf95_gw_var(grpid, varid, wavenumber)472 473 ! Hydrophobic/Monochromatic:474 call nf95_inq_grp_full_ncid(ncid, "Hydrophobic/Monochromatic", grpid)475 call nf95_inq_varid(grpid, "mass_ext", varid)476 call nf95_gw_var(grpid, varid, ao%mass_ext_mono_phobic)477 478 ! call file%get('wavenumber', wavenumber)479 ! nwn = size(wavenumber)480 481 ! call file%get_global_attribute('description_hydrophobic', &482 ! & ao%description_phobic_str)483 484 485 ! call file%get('relative_humidity1', ao%rh_lower)486 487 ! call file%get_global_attribute('description_hydrophilic', &488 ! & ao%description_philic_str)489 490 ! Close aerosol scattering file491 ! call file%close()492 493 call nf95_close(ncid)494 495 ! Get array sizes496 ! ao%n_bands_lw = size(ao%mass_ext_lw_phobic, 1)497 ! ao%n_bands_sw = size(ao%mass_ext_sw_phobic, 1)498 ! ao%n_mono_wl = size(ao%mass_ext_mono_phobic, 1)499 ! ao%n_type_phobic = size(ao%mass_ext_lw_phobic, 2)500 ! ao%n_type_philic = size(ao%mass_ext_lw_philic, 3)501 502 ! Allocate memory for mapping arrays503 ! ao%ntype = ao%n_type_phobic + ao%n_type_philic504 ! allocate(ao%iclass(ao%ntype))505 ! allocate(ao%itype(ao%ntype))506 507 ! ao%iclass = IAerosolClassUndefined508 ! ao%itype = 0509 510 n_type_phobic = size(mass_ext_phobic, 2)511 if (ao%use_hydrophilic) then512 n_type_philic = size(mass_ext_philic, 3)513 nrh = size(ao%rh_lower)514 else515 n_type_philic = 0516 nrh = 0517 end if518 519 if (config%do_cloud_aerosol_per_sw_g_point) then520 nspecsw = config%gas_optics_sw%spectral_def%ng521 else522 nspecsw = config%gas_optics_sw%spectral_def%nband523 end if524 525 if (config%do_cloud_aerosol_per_lw_g_point) then526 nspeclw = config%gas_optics_lw%spectral_def%ng527 else528 nspeclw = config%gas_optics_lw%spectral_def%nband529 end if530 531 if (allocated(ao%wavelength_mono)) then532 ! Monochromatic wavelengths also required533 nmono = size(ao%wavelength_mono)534 else535 nmono = 0536 end if537 538 call ao%allocate(n_type_phobic, n_type_philic, nrh, nspeclw, nspecsw, nmono)539 540 if (config%do_sw) then541 call config%gas_optics_sw%spectral_def%calc_mapping(SolarReferenceTemperature, &542 & wavenumber, mapping, use_bands=(.not. config%do_cloud_aerosol_per_sw_g_point))543 544 ao%mass_ext_sw_phobic = matmul(mapping, mass_ext_phobic)545 ao%ssa_sw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic) &546 & / ao%mass_ext_sw_phobic547 ao%g_sw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic*g_phobic) &548 & / (ao%mass_ext_sw_phobic*ao%ssa_sw_phobic)549 550 if (ao%use_hydrophilic) then551 do jtype = 1,n_type_philic552 ao%mass_ext_sw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype))553 ao%ssa_sw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) &554 & *ssa_philic(:,:,jtype)) &555 & / ao%mass_ext_sw_philic(:,:,jtype)556 ao%g_sw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) &557 & *ssa_philic(:,:,jtype)*g_philic(:,:,jtype)) &558 & / (ao%mass_ext_sw_philic(:,:,jtype)*ao%ssa_sw_philic(:,:,jtype))559 end do560 end if561 end if562 if (config%do_lw) then563 call config%gas_optics_lw%spectral_def%calc_mapping(TerrestrialReferenceTemperature, &564 & wavenumber, mapping, use_bands=(.not. config%do_cloud_aerosol_per_lw_g_point))565 566 ao%mass_ext_lw_phobic = matmul(mapping, mass_ext_phobic)567 ao%ssa_lw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic) &568 & / ao%mass_ext_lw_phobic569 ao%g_lw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic*g_phobic) &570 & / (ao%mass_ext_lw_phobic*ao%ssa_lw_phobic)571 572 if (ao%use_hydrophilic) then573 do jtype = 1,n_type_philic574 ao%mass_ext_lw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype))575 ao%ssa_lw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) &576 & *ssa_philic(:,:,jtype)) &577 & / ao%mass_ext_lw_philic(:,:,jtype)578 ao%g_lw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) &579 & *ssa_philic(:,:,jtype)*g_philic(:,:,jtype)) &580 & / (ao%mass_ext_lw_philic(:,:,jtype)*ao%ssa_lw_philic(:,:,jtype))581 end do582 end if583 end if584 585 if (allocated(ao%wavelength_mono)) then586 ! Monochromatic wavelengths also required587 do jwl = 1,nmono588 ! Wavelength (m) to wavenumber (cm-1)589 wavenumber_target = 0.01_jprb / ao%wavelength_mono(jwl)590 ! Find index to first interpolation point, and its weight591 if (wavenumber_target <= wavenumber(1)) then592 weight1 = 1.0_jprb593 iwn = 1594 else if (wavenumber_target >= wavenumber(nwn)) then595 iwn = nwn-1596 weight1 = 0.0_jprb597 else598 iwn = 1599 do while (wavenumber(iwn+1) < wavenumber_target .and. iwn < nwn-1)600 iwn = iwn + 1601 end do602 weight1 = (wavenumber(iwn+1)-wavenumber_target) &603 & / (wavenumber(iwn+1)-wavenumber(iwn))604 end if605 ! Linear interpolation606 ao%mass_ext_mono_phobic(jwl,:) = weight1 * mass_ext_phobic(iwn,:) &607 & + (1.0_jprb - weight1)* mass_ext_phobic(iwn+1,:)608 ao%ssa_mono_phobic(jwl,:) = weight1 * ssa_phobic(iwn,:) &609 & + (1.0_jprb - weight1)* ssa_phobic(iwn+1,:)610 ao%g_mono_phobic(jwl,:) = weight1 * g_phobic(iwn,:) &611 & + (1.0_jprb - weight1)* g_phobic(iwn+1,:)612 ao%lidar_ratio_mono_phobic(jwl,:) = weight1 * lidar_ratio_phobic(iwn,:) &613 & + (1.0_jprb - weight1)* lidar_ratio_phobic(iwn+1,:)614 if (ao%use_hydrophilic) then615 ao%mass_ext_mono_philic(jwl,:,:) = weight1 * mass_ext_philic(iwn,:,:) &616 & + (1.0_jprb - weight1)* mass_ext_philic(iwn+1,:,:)617 ao%ssa_mono_philic(jwl,:,:) = weight1 * ssa_philic(iwn,:,:) &618 & + (1.0_jprb - weight1)* ssa_philic(iwn+1,:,:)619 ao%g_mono_philic(jwl,:,:) = weight1 * g_philic(iwn,:,:) &620 & + (1.0_jprb - weight1)* g_philic(iwn+1,:,:)621 ao%lidar_ratio_mono_philic(jwl,:,:) = weight1 * lidar_ratio_philic(iwn,:,:) &622 & + (1.0_jprb - weight1)* lidar_ratio_philic(iwn+1,:,:)623 end if624 end do625 end if626 627 ! Deallocate memory local to this routine628 deallocate(mass_ext_phobic)629 deallocate(ssa_phobic)630 deallocate(g_phobic)631 deallocate(lidar_ratio_phobic)632 if (ao%use_hydrophilic) then633 deallocate(mass_ext_philic)634 deallocate(ssa_philic)635 deallocate(g_philic)636 deallocate(lidar_ratio_philic)637 end if638 639 if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics',1,hook_handle)640 641 end subroutine setup_general_aerosol_optics_lmdz642 643 335 644 336 !---------------------------------------------------------------------
Note: See TracChangeset
for help on using the changeset viewer.