Ignore:
Timestamp:
Feb 21, 2023, 3:26:41 PM (20 months ago)
Author:
idelkadi
Message:

Update of the ECRAD radiative code version implemented in the LMDZ model.
Upgrade to the :
https://github.com/ecmwf/ecrad/trunk
Version svn : 749
UUID du dépôt : 44b0ca93-0ed8-356e-d663-ce57b7db7bff

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ6/branches/LMDZ_ECRad/libf/phylmd/ecrad/radiation_aerosol_optics_data.F90

    r3908 r4444  
    5757     integer, allocatable, dimension(:) :: itype
    5858
     59     ! Wavenumber (cm-1) upper and lower bounds of each spectral
     60     ! interval, which if used in the RRTMG gas optics scheme should
     61     ! match its band bounds
     62     real(jprb), allocatable, dimension(:) :: wavenumber1_sw, wavenumber2_sw
     63     real(jprb), allocatable, dimension(:) :: wavenumber1_lw, wavenumber2_lw
     64
    5965     ! Scattering properties are provided separately in the shortwave
    6066     ! and longwave for hydrophobic and hydrophilic aerosols.
     
    6470          &  ssa_sw_phobic,      & ! Single scattering albedo
    6571          &  g_sw_phobic,        & ! Asymmetry factor
     72!          &  ssa_g_sw_phobic,    & ! ssa*g
    6673          &  mass_ext_lw_phobic, & ! Mass-extinction coefficient (m2 kg-1)
     74!          &  mass_abs_lw_phobic, & ! Mass-absorption coefficient (m2 kg-1)
    6775          &  ssa_lw_phobic,      & ! Single scattering albedo
    6876          &  g_lw_phobic           ! Asymmetry factor
    6977
    70      ! Hydrophilic aerosols are dimensioned (nband, nrh, n_type_philic):
     78     ! Hydrophilic aerosols are dimensioned (nband,nrh,n_type_philic):
    7179     real(jprb), allocatable, dimension(:,:,:) :: &
    7280          &  mass_ext_sw_philic, & ! Mass-extinction coefficient (m2 kg-1)
    7381          &  ssa_sw_philic,      & ! Single scattering albedo
    7482          &  g_sw_philic,        & ! Asymmetry factor
     83 !         &  ssa_g_sw_philic,    & ! ssa*g
    7584          &  mass_ext_lw_philic, & ! Mass-extinction coefficient (m2 kg-1)
     85 !         &  mass_abs_lw_philic, & ! Mass-absorption coefficient (m2 kg-1)
    7686          &  ssa_lw_philic,      & ! Single scattering albedo
    7787          &  g_lw_philic           ! Asymmetry factor
    7888
    79      ! Scattering properties at selected wavelengths
    80      ! (n_mono_wl,n_type_phobic/philic)
     89     ! Wavelengths at which monochromatic properties are stored,
     90     ! dimensioned (n_mono_wl), units metres
     91     real(jprb), allocatable :: wavelength_mono(:)
     92
     93     ! Scattering properties at selected monochromatic wavelengths
     94     ! (n_mono_wl,n_type_phobic)
    8195     real(jprb), allocatable, dimension(:,:) :: &
    8296          &  mass_ext_mono_phobic, & ! Mass-extinction coefficient (m2 kg-1)
     
    8498          &  g_mono_phobic,        & ! Asymmetry factor
    8599          &  lidar_ratio_mono_phobic ! Lidar Ratio
     100     ! ...hydrophilic aerosols dimensioned (n_mono_wl,nrh,n_type_philic):
    86101     real(jprb), allocatable, dimension(:,:,:) :: &
    87102          &  mass_ext_mono_philic, & ! Mass-extinction coefficient (m2 kg-1)
     
    104119     ! The number of hydrophobic and hydrophilic types read from the
    105120     ! aerosol optics file
    106      integer :: n_type_phobic, n_type_philic
     121     integer :: n_type_phobic = 0
     122     integer :: n_type_philic = 0
    107123
    108124     ! Number of relative humidity bins
    109      integer :: nrh
     125     integer :: nrh = 0
    110126
    111127     ! Number of longwave and shortwave bands of the data in the file,
     
    121137   contains
    122138     procedure :: setup => setup_aerosol_optics
     139     procedure :: save  => save_aerosol_optics
     140     procedure :: allocate
     141     procedure :: initialize_types
    123142     procedure :: set_hydrophobic_type
    124143     procedure :: set_hydrophilic_type
     
    135154  !---------------------------------------------------------------------
    136155  ! Setup aerosol optics coefficients by reading them from a file
    137   subroutine setup_aerosol_optics(this, file_name, ntype, iverbose)
     156  subroutine setup_aerosol_optics(this, file_name, iverbose)
    138157
    139158    use yomhook,              only : lhook, dr_hook
     
    143162    class(aerosol_optics_type), intent(inout) :: this
    144163    character(len=*), intent(in)              :: file_name
    145     integer, intent(in)                       :: ntype
    146164    integer, intent(in), optional             :: iverbose
    147165
    148166    ! The NetCDF file containing the aerosol optics data
    149167    type(netcdf_file)  :: file
     168
     169    real(jprb), allocatable :: wavelength_tmp(:)
    150170    integer            :: iverb
     171
    151172    real(jprb)         :: hook_handle
    152173
     
    168189      this%use_hydrophilic = .false.
    169190    end if
     191
     192    ! Read the wavenumber bounds
     193    call file%get('wavenumber1_sw', this%wavenumber1_sw)
     194    call file%get('wavenumber2_sw', this%wavenumber2_sw)
     195    call file%get('wavenumber1_lw', this%wavenumber1_lw)
     196    call file%get('wavenumber2_lw', this%wavenumber2_lw)
    170197
    171198    ! Read the raw scattering data
     
    180207         &                         this%description_phobic_str)
    181208
     209    ! Precompute ssa*g for the shortwave and mass-absorption for the
     210    ! longwave - TBD FIX
     211    !allocate(this%ssa_g_sw_phobic(...
     212
    182213    if (this%use_hydrophilic) then
    183214      call file%get('mass_ext_sw_hydrophilic', this%mass_ext_sw_philic)
     
    194225    end if
    195226
    196     ! Read the raw scattering data at selected wavelengths if
    197     ! available in the input file
     227    ! Read the monochromatic scattering data at selected wavelengths
     228    ! if available in the input file
    198229    if (file%exists('mass_ext_mono_hydrophobic')) then
    199230      this%use_monochromatic = .true.
     231
     232      if (allocated(this%wavelength_mono)) then
     233        ! User has provided required monochromatic wavelengths, which
     234        ! must match those in the file (in the more recent "general"
     235        ! aerosol optics, interpolation provides optical properties at
     236        ! the requested wavelengths)
     237        call file%get('wavelength_mono', wavelength_tmp)
     238        if (size(wavelength_tmp) /= size(this%wavelength_mono)) then
     239          write(nulerr,'(a,i0,a,i0,a)') '*** Error: ', size(this%wavelength_mono), &
     240               &  ' monochromatic wavelengths requested but ', &
     241               &  size(wavelength_tmp), ' in file'
     242          call radiation_abort('Radiation configuration error')
     243        end if
     244        if (any(abs(this%wavelength_mono-wavelength_tmp) &
     245               &  / this%wavelength_mono > 0.01_jprb)) then
     246          write(nulerr,'(a,a)') '*** Error: requested monochromatic wavelengths', &
     247               &  'must all be within 1% of values in file'
     248          call radiation_abort('Radiation configuration error')
     249        end if
     250      else
     251        ! User has not provided required wavelengths, so we save the
     252        ! monochromatic wavelengths in the file
     253        call file%get('wavelength_mono', this%wavelength_mono)
     254      end if
     255
    200256      call file%get('mass_ext_mono_hydrophobic', this%mass_ext_mono_phobic)
    201257      call file%get('ssa_mono_hydrophobic',      this%ssa_mono_phobic)
     
    233289        write(nulerr,'(a,a)') '*** Error: mass extinction for hydrophilic and hydrophobic ', &
    234290             &                'aerosol have different numbers of longwave bands'
    235         call radiation_abort()
     291        call radiation_abort('Radiation configuration error')
    236292      end if
    237293      if (size(this%mass_ext_sw_philic,1) /= this%n_bands_sw) then
    238294        write(nulerr,'(a,a)') '*** Error: mass extinction for hydrophilic and hydrophobic ', &
    239295             &                'aerosol have different numbers of shortwave bands'
    240         call radiation_abort()
     296        call radiation_abort('Radiation configuration error')
    241297      end if
    242298      if (size(this%rh_lower) /= this%nrh) then
    243299        write(nulerr,'(a)') '*** Error: size(relative_humidity1) /= size(mass_ext_sw_hydrophilic,2)'
    244         call radiation_abort()
     300        call radiation_abort('Radiation configuration error')
    245301      end if
    246302
     
    250306    end if
    251307
     308    if (lhook) call dr_hook('radiation_aerosol_optics_data:setup',1,hook_handle)
     309
     310  end subroutine setup_aerosol_optics
     311
     312
     313  !---------------------------------------------------------------------
     314  ! Initialize the arrays describing the user's aerosol types
     315  subroutine initialize_types(this, ntype)
     316
     317    class(aerosol_optics_type), intent(inout) :: this
     318    integer,                    intent(in)    :: ntype
     319   
    252320    ! Allocate memory for mapping arrays
    253321    this%ntype = ntype
     
    258326    this%itype  = 0
    259327
    260     if (lhook) call dr_hook('radiation_aerosol_optics_data:setup',1,hook_handle)
    261 
    262   end subroutine setup_aerosol_optics
     328  end subroutine initialize_types
     329
     330  !---------------------------------------------------------------------
     331  ! Allocate arrays for aerosol optics data type
     332  subroutine allocate(this, n_type_phobic, n_type_philic, nrh, &
     333       &              n_bands_lw, n_bands_sw, n_mono_wl)
     334
     335    use yomhook,     only : lhook, dr_hook
     336
     337    class(aerosol_optics_type), intent(inout) :: this
     338    integer, intent(in) :: n_type_phobic, n_type_philic, nrh
     339    integer, intent(in) :: n_bands_lw, n_bands_sw, n_mono_wl
     340
     341    real(jprb) :: hook_handle
     342
     343    if (lhook) call dr_hook('radiation_aerosol_optics_data:allocate',0,hook_handle)
     344
     345    this%n_type_phobic = n_type_phobic
     346    this%n_type_philic = n_type_philic
     347    this%nrh           = nrh
     348    this%n_bands_lw    = n_bands_lw
     349    this%n_bands_sw    = n_bands_sw
     350    this%n_mono_wl     = n_mono_wl
     351
     352    if (n_type_philic > 0) then
     353      this%use_hydrophilic = .true.
     354    else
     355      this%use_hydrophilic = .false.
     356    end if
     357
     358    if (n_bands_sw > 0) then
     359      allocate(this%mass_ext_sw_phobic(n_bands_sw, n_type_phobic))
     360      allocate(this%ssa_sw_phobic(n_bands_sw, n_type_phobic))
     361      allocate(this%g_sw_phobic(n_bands_sw, n_type_phobic))
     362    end if
     363    if (n_bands_lw > 0) then
     364      allocate(this%mass_ext_lw_phobic(n_bands_lw, n_type_phobic))
     365      allocate(this%ssa_lw_phobic(n_bands_lw, n_type_phobic))
     366      allocate(this%g_lw_phobic(n_bands_lw, n_type_phobic))
     367    end if
     368    if (n_mono_wl > 0) then
     369      allocate(this%mass_ext_mono_phobic(n_mono_wl, n_type_phobic))
     370      allocate(this%ssa_mono_phobic(n_mono_wl, n_type_phobic))
     371      allocate(this%g_mono_phobic(n_mono_wl, n_type_phobic))
     372      allocate(this%lidar_ratio_mono_phobic(n_mono_wl, n_type_phobic))
     373    end if
     374
     375    if (n_type_philic > 0 .and. nrh > 0) then
     376      if (n_bands_sw > 0) then
     377        allocate(this%mass_ext_sw_philic(n_bands_sw, nrh, n_type_philic))
     378        allocate(this%ssa_sw_philic(n_bands_sw, nrh, n_type_philic))
     379        allocate(this%g_sw_philic(n_bands_sw, nrh, n_type_philic))
     380      end if
     381      if (n_bands_lw > 0) then
     382        allocate(this%mass_ext_lw_philic(n_bands_lw, nrh, n_type_philic))
     383        allocate(this%ssa_lw_philic(n_bands_lw, nrh, n_type_philic))
     384        allocate(this%g_lw_philic(n_bands_lw, nrh, n_type_philic))
     385      end if
     386      if (n_mono_wl > 0) then
     387        allocate(this%mass_ext_mono_philic(n_mono_wl, nrh, n_type_philic))
     388        allocate(this%ssa_mono_philic(n_mono_wl, nrh, n_type_philic))
     389        allocate(this%g_mono_philic(n_mono_wl, nrh, n_type_philic))
     390        allocate(this%lidar_ratio_mono_philic(n_mono_wl, nrh, n_type_philic))
     391      end if
     392    end if
     393
     394    if (lhook) call dr_hook('radiation_aerosol_optics_data:allocate',1,hook_handle)
     395
     396  end subroutine allocate
     397
     398
     399  !---------------------------------------------------------------------
     400  subroutine save_aerosol_optics(this, file_name, iverbose)
     401
     402    use yomhook,     only : lhook, dr_hook
     403    use easy_netcdf, only : netcdf_file
     404
     405    class(aerosol_optics_type), intent(inout) :: this
     406    character(len=*),           intent(in)    :: file_name
     407    integer,          optional, intent(in)    :: iverbose
     408
     409    ! Object for output NetCDF file
     410    type(netcdf_file) :: out_file
     411
     412    real(jprb) :: hook_handle
     413
     414    if (lhook) call dr_hook('radiation_aerosol_optics_data:save',0,hook_handle)
     415
     416    ! Create the file
     417    call out_file%create(trim(file_name), iverbose=iverbose)
     418
     419    ! Define dimensions
     420    call out_file%define_dimension("band_lw", this%n_bands_lw)
     421    call out_file%define_dimension("band_sw", this%n_bands_sw)
     422    call out_file%define_dimension("hydrophilic", this%n_type_philic)
     423    call out_file%define_dimension("hydrophobic", this%n_type_phobic)
     424    call out_file%define_dimension("relative_humidity", this%nrh)
     425    !if (this%use_monochromatic) then
     426    !  call out_file%define_dimension("wavelength_mono", this%n_mono_wl)
     427    !end if
     428
     429    ! Put global attributes
     430    call out_file%put_global_attributes( &
     431         &   title_str="Aerosol optical properties in the spectral intervals of the gas-optics scheme for ecRad", &
     432         &   source_str="ecRad offline radiation model")
     433    call out_file%put_global_attribute( &
     434         &  "description_hydrophobic", this%description_phobic_str)
     435    call out_file%put_global_attribute( &
     436         &  "description_hydrophilic", this%description_philic_str)
     437
     438    ! Define variables
     439    call out_file%define_variable("mass_ext_sw_hydrophobic", units_str="m2 kg-1", &
     440         &  long_name="Shortwave mass-extinction coefficient of hydrophobic aerosols", &
     441         &  dim2_name="hydrophobic", dim1_name="band_sw")
     442    call out_file%define_variable("ssa_sw_hydrophobic", units_str="1", &
     443         &  long_name="Shortwave single scattering albedo of hydrophobic aerosols", &
     444         &  dim2_name="hydrophobic", dim1_name="band_sw")
     445    call out_file%define_variable("asymmetry_sw_hydrophobic", units_str="1", &
     446         &  long_name="Shortwave asymmetry factor of hydrophobic aerosols", &
     447         &  dim2_name="hydrophobic", dim1_name="band_sw")
     448
     449    call out_file%define_variable("mass_ext_lw_hydrophobic", units_str="m2 kg-1", &
     450         &  long_name="Longwave mass-extinction coefficient of hydrophobic aerosols", &
     451         &  dim2_name="hydrophobic", dim1_name="band_lw")
     452    call out_file%define_variable("ssa_lw_hydrophobic", units_str="1", &
     453         &  long_name="Longwave single scattering albedo of hydrophobic aerosols", &
     454         &  dim2_name="hydrophobic", dim1_name="band_lw")
     455    call out_file%define_variable("asymmetry_lw_hydrophobic", units_str="1", &
     456         &  long_name="Longwave asymmetry factor of hydrophobic aerosols", &
     457         &  dim2_name="hydrophobic", dim1_name="band_lw")
     458
     459    call out_file%define_variable("mass_ext_sw_hydrophilic", units_str="m2 kg-1", &
     460         &  long_name="Shortwave mass-extinction coefficient of hydrophilic aerosols", &
     461         &  dim3_name="hydrophilic", dim2_name="relative_humidity", dim1_name="band_sw")
     462    call out_file%define_variable("ssa_sw_hydrophilic", units_str="1", &
     463         &  long_name="Shortwave single scattering albedo of hydrophilic aerosols", &
     464         &  dim3_name="hydrophilic", dim2_name="relative_humidity", dim1_name="band_sw")
     465    call out_file%define_variable("asymmetry_sw_hydrophilic", units_str="1", &
     466         &  long_name="Shortwave asymmetry factor of hydrophilic aerosols", &
     467         &  dim3_name="hydrophilic", dim2_name="relative_humidity", dim1_name="band_sw")
     468
     469    call out_file%define_variable("mass_ext_lw_hydrophilic", units_str="m2 kg-1", &
     470         &  long_name="Longwave mass-extinction coefficient of hydrophilic aerosols", &
     471         &  dim3_name="hydrophilic", dim2_name="relative_humidity", dim1_name="band_lw")
     472    call out_file%define_variable("ssa_lw_hydrophilic", units_str="1", &
     473         &  long_name="Longwave single scattering albedo of hydrophilic aerosols", &
     474         &  dim3_name="hydrophilic", dim2_name="relative_humidity", dim1_name="band_lw")
     475    call out_file%define_variable("asymmetry_lw_hydrophilic", units_str="1", &
     476         &  long_name="Longwave asymmetry factor of hydrophilic aerosols", &
     477         &  dim3_name="hydrophilic", dim2_name="relative_humidity", dim1_name="band_lw")
     478
     479    ! Write variables
     480    call out_file%put("mass_ext_sw_hydrophobic", this%mass_ext_sw_phobic)
     481    call out_file%put("ssa_sw_hydrophobic", this%ssa_sw_phobic)
     482    call out_file%put("asymmetry_sw_hydrophobic", this%g_sw_phobic)
     483    call out_file%put("mass_ext_lw_hydrophobic", this%mass_ext_lw_phobic)
     484    call out_file%put("ssa_lw_hydrophobic", this%ssa_lw_phobic)
     485    call out_file%put("asymmetry_lw_hydrophobic", this%g_lw_phobic)
     486    call out_file%put("mass_ext_sw_hydrophilic", this%mass_ext_sw_philic)
     487    call out_file%put("ssa_sw_hydrophilic", this%ssa_sw_philic)
     488    call out_file%put("asymmetry_sw_hydrophilic", this%g_sw_philic)
     489    call out_file%put("mass_ext_lw_hydrophilic", this%mass_ext_lw_philic)
     490    call out_file%put("ssa_lw_hydrophilic", this%ssa_lw_philic)
     491    call out_file%put("asymmetry_lw_hydrophilic", this%g_lw_philic)
     492
     493    call out_file%close()
     494
     495    if (lhook) call dr_hook('radiation_aerosol_optics_data:save',1,hook_handle)
     496
     497  end subroutine save_aerosol_optics
    263498
    264499
Note: See TracChangeset for help on using the changeset viewer.