Ignore:
Timestamp:
Apr 9, 2026, 11:48:44 AM (7 days ago)
Author:
flefevre
Message:

Photolysis: 1. updated H2SO4 cross-sections (UV and vibrational overtones) 2. optimisation of the wavelength grid in the H2SO4 vibrational overtones.

The updated H2SO4 cross-section file can be downloaded from this link:

https://owncloud.latmos.ipsl.fr/index.php/s/WGrRb7omYZf9ZAK

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.VENUS/libf/phyvenus/photolysis_mod.F90

    r4111 r4176  
    55! photolysis
    66
    7   integer, save :: nphot = 33             ! number of photolysis
     7  integer, save :: nphot = 34             ! number of photolysis
    88
    99!$OMP THREADPRIVATE(nphot)
     
    1313! spectral grid
    1414
    15   integer, parameter :: nw = 194          ! number of spectral intervals (low-res)
     15  integer, parameter :: nw = 208          ! number of spectral intervals (if mopt = 3)
    1616  integer, save :: mopt                   ! high-res/low-res switch
    1717
     
    5151  real, dimension(nw), save :: xsocs                                  ! cos absorption cross-section (cm2)
    5252  real, dimension(nw), save :: xscocl2                                ! cocl2 absorption cross-section (cm2)
    53   real, dimension(nw), save :: xsh2so4                                ! h2so4 absorption cross-section (cm2)
     53  real, dimension(nw), save :: xsh2so4_hso3                           ! h2so4 absorption cross-section (cm2), hso3 + oh channel
     54  real, dimension(nw), save :: xsh2so4_so3                            ! h2so4 absorption cross-section (cm2), so3 + h2o channel
    5455  real, dimension(nw), save :: xsh2                                   ! h2 absorption cross-section (cm2)
    5556  real, dimension(nw), save :: yieldh2                                ! h2 photodissociation yield
     
    168169! read and grid h2so4 cross-sections
    169170
    170   call rdxsh2so4(nw,wl,xsh2so4)
     171  call rdxsh2so4(nw,wl,xsh2so4_hso3,xsh2so4_so3)
    171172
    172173! read and grid h2 cross-sections
     
    245246!                 365-850 nm :  5.0  nm 
    246247!
    247 !     mopt = 2    low-resolution mode (162)
     248!     mopt = 2    low-resolution mode (162 intervals)
    248249!
    249250!                    0-60 nm :  6.0 nm
     
    260261!                 415-815 nm : 50.0 nm
    261262!
    262 !     mopt = 3    venusian low-resolution mode (194)
    263 !                 (205-245 nm -> 1nm of resolution)
     263!     mopt = 3   venusian low-resolution mode (208 intervals)
     264!                205-245 nm -> 1nm of resolution,
     265!                optimised for H2SO4 overtones
    264266!
    265267!                    0-60 nm :  6.0 nm
     
    274276!                 205-245 nm :  1.0 nm
    275277!                 245-415 nm : 10.0 nm
    276 !                 415-815 nm : 50.0 nm
     278!                 415-565 nm : 50.0 nm
     279!                 565-604 nm : 39.0 nm
     280!                 604-608 nm :  1.0 nm
     281!                 608-736 nm : 32.0 nm
     282!                 736-745 nm :  1.0 nm
     283!                 745-815 nm : 35.0 nm
    277284
    278285      if (mopt == 1) then   ! high-res
     
    562569      ENDDO
    563570
    564 ! define wavelength intervals of width 5.0 nm from 205 to 245 nm:
     571! define wavelength intervals of width 1.0 nm from 205 to 245 nm:
    565572
    566573      wincr = 1
     
    582589      ENDDO
    583590
    584 ! define wavelength intervals of width 50.0 nm from 415 to 815 nm:
     591! define wavelength intervals of width 50.0 nm from 415 to 565 nm:
    585592
    586593      wincr = 50.0
    587       DO iw = 415, 815, 50
     594      DO iw = 415, 515, 50
     595        kw = kw + 1
     596        wl(kw) = real(iw)
     597        wu(kw) = wl(kw) + wincr
     598        wc(kw) = (wl(kw) + wu(kw))/2.
     599      ENDDO
     600
     601! define wavelength intervals of width 39.0 nm from 565 to 604 nm:
     602
     603      wincr = 39.0
     604      DO iw = 565, 565, 39
     605        kw = kw + 1
     606        wl(kw) = real(iw)
     607        wu(kw) = wl(kw) + wincr
     608        wc(kw) = (wl(kw) + wu(kw))/2.
     609      ENDDO
     610
     611! define wavelength intervals of width 1.0 nm from 604 to 608 nm:
     612
     613      wincr = 1.0
     614      DO iw = 604, 607, 1
     615        kw = kw + 1
     616        wl(kw) = real(iw)
     617        wu(kw) = wl(kw) + wincr
     618        wc(kw) = (wl(kw) + wu(kw))/2.
     619      ENDDO
     620
     621! define wavelength intervals of width 32.0 nm from 608 to 736 nm:
     622
     623      wincr = 32.0
     624      DO iw = 608, 704, 32
     625        kw = kw + 1
     626        wl(kw) = real(iw)
     627        wu(kw) = wl(kw) + wincr
     628        wc(kw) = (wl(kw) + wu(kw))/2.
     629      ENDDO
     630
     631! define wavelength intervals of width 1.0 nm from 736 to 745 nm:
     632
     633      wincr = 1.0
     634      DO iw = 736, 744, 1
     635        kw = kw + 1
     636        wl(kw) = real(iw)
     637        wu(kw) = wl(kw) + wincr
     638        wc(kw) = (wl(kw) + wu(kw))/2.
     639      ENDDO
     640
     641! define wavelength intervals of width 35.0 nm from 745 to 815 nm:
     642
     643      wincr = 35.0
     644      DO iw = 745, 780, 35
    588645        kw = kw + 1
    589646        wl(kw) = real(iw)
     
    31273184!==============================================================================
    31283185
    3129       subroutine rdxsh2so4(nw, wl, yg)
     3186      subroutine rdxsh2so4(nw, wl, xsh2so4_hso3, xsh2so4_so3)
    31303187
    31313188!-----------------------------------------------------------------------------*
    31323189!=  PURPOSE:                                                                 =*
    3133 !=  Read H2SO4 cross-sections                                                  =*
    3134 !=  JPL 2006 recommendation                                                  =*
     3190!=  Read H2SO4 cross-sections                                                =*
     3191!=  Frandsen et al., personal communication, 2026                            =*
    31353192!-----------------------------------------------------------------------------*
    31363193!=  PARAMETERS:                                                              =*
     
    31423199      USE mod_phys_lmdz_transfert_para, ONLY: bcast
    31433200
    3144       IMPLICIT NONE
     3201      implicit none
    31453202
    31463203!     input
     
    31513208!     output
    31523209
    3153       real, dimension(nw) :: yg   ! h2so4 cross-sections (cm2)
     3210      real, dimension(nw) :: xsh2so4_hso3   ! h2so4 -> hso3 + oh cross-sections (cm2)
     3211      real, dimension(nw) :: xsh2so4_so3    ! h2so4 -> so3 + h2o cross-sections (cm2)
    31543212
    31553213!     local
    31563214
    31573215      real, parameter :: deltax = 1.e-4
    3158       integer, parameter :: kdata = 100
     3216      integer, parameter :: kdata = 10000
    31593217      real, dimension(kdata) :: x1, y1
     3218      real, dimension(nw) :: yg
    31603219      integer :: i, n, ierr
    31613220      character*100 fil
     
    31643223      kin = 10
    31653224
    3166 !*** cross sections from JPL [2006]
    3167 
    3168       fil = 'cross_sections/h2so4_cross_sections.txt'
     3225!*** cross sections
     3226
     3227      fil = 'cross_sections/h2so4_hdso4_frandsen_2026.txt'
    31693228      print*, 'section efficace h2so4: ', fil
    31703229
    3171       if(is_master) then
    3172 
    3173          n = 22
     3230      if (is_master) then
     3231
     3232      n = 9602
    31743233      OPEN(kin,FILE=fil,STATUS='OLD')
    3175       DO i = 1,3
     3234      DO i = 1,5
    31763235         READ(kin,*)
    31773236      ENDDO   
     
    31873246
    31883247      CALL inter2(nw,wl,yg,n,x1,y1,ierr)
    3189  
     3248     
    31903249      IF (ierr .NE. 0) THEN
    31913250        WRITE(*,*) ierr, fil
     
    31933252      ENDIF
    31943253
    3195       endif !is_master
    3196 
    3197       call bcast(yg)
     3254      do i = 1, nw - 1
     3255         if (wl(i) < 200.) then
     3256            xsh2so4_hso3(i) = yg(i)
     3257            xsh2so4_so3(i)  = 0.
     3258         else
     3259            xsh2so4_hso3(i) = 0.
     3260            xsh2so4_so3(i)  = yg(i)
     3261         end if
     3262      end do
     3263
     3264      end if !is_master
     3265
     3266      call bcast(xsh2so4_hso3)
     3267      call bcast(xsh2so4_so3)
    31983268 
    31993269      end subroutine rdxsh2so4
Note: See TracChangeset for help on using the changeset viewer.