Changeset 3012


Ignore:
Timestamp:
Jul 24, 2023, 9:31:31 AM (17 months ago)
Author:
emillour
Message:

Mars PCM:
Code cleanup concerning chemistry. Turn chimiedata.h into module
chemistrydata.F90 and integrate read_phototable.F in it.
Also fix an OpenMP issue with read_phototable; different OpenMP threads
should not simultaneously open/read a file. Let the master do the job
and then broadcast the result to all cores.
While at it, also turned nltecool.F and nlte_tcool.F into modules.
EM

Location:
trunk/LMDZ.MARS
Files:
1 deleted
7 edited
1 moved

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/changelog.txt

    r3008 r3012  
    41244124microphys_h.F90, and while at it also turn nuclea.F, growthrate.F90 and
    41254125massflowrateco2.F90 into modules.
     4126
     4127== 24/07/2023 == EM
     4128Code cleanup concerning chemistry. Turn chimiedata.h into module
     4129chemistrydata.F90 and integrate read_phototable.F in it.
     4130Also fix an OpenMP issue with read_phototable; different OpenMP threads
     4131should not simultaneously open/read a file. Let the master do the job
     4132and then broadcast the result to all cores.
     4133While at it, also turned nltecool.F and nlte_tcool.F into modules.
  • trunk/LMDZ.MARS/libf/aeronomars/calchim_mod.F90

    r2613 r3012  
    3131      use conc_mod, only: mmean ! mean molecular mass of the atmosphere
    3232      use comcstfi_h, only: pi
     33      use chemistrydata_mod, only: read_phototable
    3334      use photolysis_mod, only: init_photolysis, nphot
    3435      use iono_h, only: temp_elect
  • trunk/LMDZ.MARS/libf/aeronomars/chemistrydata.F90

    r3011 r3012  
     1MODULE chemistrydata_mod
    12!--------------------------------------------
    23!     data for photochemistry
    34!--------------------------------------------
    4 
     5IMPLICIT NONE
    56!--------------------------------------------
    67!     dimensions of photolysis lookup table
    78!--------------------------------------------
    89
    9       integer, parameter :: nd    = 13  ! species
    10       integer, parameter :: nz    = 143 ! altitude
    11       integer, parameter :: nozo  = 7   ! ozone
    12       integer, parameter :: nsza  = 27  ! solar zenith angle
    13       integer, parameter :: ntemp = 4   ! temperature
    14       integer, parameter :: ntau  = 8   ! dust
     10integer, parameter :: nd    = 13  ! species
     11integer, parameter :: nz    = 143 ! altitude
     12integer, parameter :: nozo  = 7   ! ozone
     13integer, parameter :: nsza  = 27  ! solar zenith angle
     14integer, parameter :: ntemp = 4   ! temperature
     15integer, parameter :: ntau  = 8   ! dust
    1516
    1617!--------------------------------------------
    1718
    18       common/chimiedata/jphot,colairtab,table_ozo
     19! tabulated solar zenith angles
     20real,parameter :: szatab(nsza) = [ 0.,  5., 10., 15., 20., 25., &
     21                                  30., 35., 40., 45., 50., 55., &
     22                                  60., 65., 70., 75., 80., 82., &
     23                                  84., 86., 88., 90., 91., 92., &
     24                                  93., 94., 95. ]
    1925
    20 !$OMP THREADPRIVATE(/chimiedata/)
     26! tabulated opacities
     27real,parameter :: tautab(ntau)=[0., 0.2, 0.4, 0.6, 0.8, 1., 2., 4.]
    2128
    22       real jphot(ntemp,nsza,nz,nozo,ntau,nd)
    23       real colairtab(nz)
    24       real szatab(nsza)
    25       real table_ozo(nozo)
    26       real tautab(ntau)
    2729
    28       data szatab/0.,  5., 10., 15., 20., 25.,                          &
    29      &            30., 35., 40., 45., 50., 55.,                         &
    30      &            60., 65., 70., 75., 80., 82.,                         &
    31      &            84., 86., 88., 90., 91., 92.,                         &
    32      &            93., 94., 95./
     30real,save,protected :: jphot(ntemp,nsza,nz,nozo,ntau,nd)
     31!$OMP THREADPRIVATE(jphot)
     32real,save,protected :: colairtab(nz)
     33!$OMP THREADPRIVATE(colairtab)
     34real,save,protected :: table_ozo(nozo)
     35!$OMP THREADPRIVATE(table_ozo)
    3336
    34       data tautab/0., 0.2, 0.4, 0.6, 0.8, 1., 2., 4./
     37CONTAINS
     38
     39  subroutine read_phototable
     40
     41!***********************************************************************
     42!
     43!   subject:
     44!   --------
     45!
     46!   read photolysis lookup table
     47!
     48!   VERSION: 8/10/2014
     49!
     50!   Author:   Franck Lefevre
     51!
     52!   Arguments:
     53!   ----------
     54!
     55!***********************************************************************
     56
     57  use ioipsl_getin_p_mod, only: getin_p
     58  use datafile_mod, only: datadir
     59  use mod_phys_lmdz_para, only: is_master
     60  use mod_phys_lmdz_transfert_para, only: bcast
     61
     62  implicit none
     63
     64!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
     65!     local:
     66!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
     67
     68  integer :: fic, ij, iozo, isza, itemp, iz, itau, ierr
     69  real    :: xsza
     70
     71  character(len = 128) :: phototable ! photolysis table file name
     72
     73!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
     74! set photolysis table input file name
     75!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
     76
     77  phototable = "jmars.20140930" ! default
     78
     79! look for a " phototable= ..." option in def files
     80
     81  call getin_p("phototable",phototable)
     82
     83  fic = 81
     84     
     85  if (is_master) then ! only the master needs to open file and read data
     86
     87    open(fic, form = 'formatted', status = 'old',                &
     88           file =trim(datadir)//"/"//trim(phototable),iostat=ierr)
     89
     90    if (ierr /= 0) THEN
     91        write(*,*)'Error : cannot open photolysis lookup table ', trim(phototable)
     92        write(*,*)'It should be in :',trim(datadir),'/'
     93        write(*,*)'1) You can change this directory in callphys.def'
     94        write(*,*)'   with:'
     95        write(*,*)'   datadir=/path/to/the/directory'
     96        write(*,*)'2) You can change the input phototable file name in'
     97        write(*,*)'   callphys.def with:'
     98        write(*,*)'   phototable=filename'
     99        call abort_physic("read_phototable","missing "//trim(phototable)//"file",1)
     100    end if
     101
     102!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
     103! read photolys table
     104!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
     105
     106    print*, 'read photolysis lookup table ',trim(phototable)
     107
     108    do itau = 1,ntau
     109       do itemp = 1,ntemp
     110          do iozo = 1,nozo
     111             do isza = 1,nsza
     112                do iz = nz,1,-1
     113                   read(fic,*) colairtab(iz), xsza, table_ozo(iozo)
     114                   read(fic,'(7e11.4)') (jphot(itemp,isza,iz,iozo,itau,ij), ij= 1,nd)
     115                   do ij = 1,nd
     116                      if (jphot(itemp,isza,iz,iozo,itau,ij) == 1.e-30) then
     117                         jphot(itemp,isza,iz,iozo,itau,ij) = 0.
     118                      end if
     119                   end do
     120                end do
     121             end do
     122          end do
     123       end do
     124    end do
     125
     126    print*, 'lookup table...ok'
     127    close(fic)
     128
     129  endif ! of if (is_master)
     130   
     131  ! broadcast the information to all cores
     132  call bcast(colairtab)
     133  call bcast(table_ozo)
     134  call bcast(jphot)
     135
     136end subroutine read_phototable
     137
     138END MODULE chemistrydata_mod
  • trunk/LMDZ.MARS/libf/aeronomars/photolysis.F90

    r2170 r3012  
    77!==========================================================================
    88
    9       use comcstfi_h
     9      use comcstfi_h, only: g
     10      use chemistrydata_mod, only: nd, nz, nozo, nsza, ntemp, ntau
     11      use chemistrydata_mod, only: szatab, tautab, colairtab, table_ozo, jphot
    1012
    1113      implicit none
    12 
    13 #include "chimiedata.h"
    1414
    1515!==========================================================================
     
    1919      integer, intent(in) :: nlayer      ! number of atmospheric layers
    2020      integer, intent(in) :: nb_phot_max ! number of processes treated numerically as photodissociations
    21       integer :: lswitch                 ! interface level between chemistries
    22       real :: press(nlayer)              ! pressure (hPa)
    23       real :: temp(nlayer)               ! temperature (K)
    24       real :: sza                        ! solar zenith angle (deg)
    25       real :: tauref                     ! optical depth at 7 hpa
    26       real :: zmmean(nlayer)             ! mean molecular mass (g)
    27       real :: dist_sol                   ! sun distance (AU)
    28       real :: rmco2(nlayer)              ! co2 mixing ratio
    29       real :: rmo3(nlayer)               ! ozone mixing ratio
     21      integer,intent(in) :: lswitch      ! interface level between chemistries
     22      real,intent(in) :: press(nlayer)   ! pressure (hPa)
     23      real,intent(in) :: temp(nlayer)    ! temperature (K)
     24      real,intent(in) :: sza             ! solar zenith angle (deg)
     25      real,intent(in) :: tauref          ! optical depth at 7 hpa
     26      real,intent(in) :: zmmean(nlayer)  ! mean molecular mass (g)
     27      real,intent(in) :: dist_sol        ! sun distance (AU)
     28      real,intent(in) :: rmco2(nlayer)   ! co2 mixing ratio
     29      real,intent(in) :: rmo3(nlayer)    ! ozone mixing ratio
    3030
    3131!==========================================================================
     
    3333!==========================================================================
    3434
    35       real (kind = 8), dimension(nlayer,nb_phot_max) :: v_phot
     35      real (kind = 8),intent(out),dimension(nlayer,nb_phot_max) :: v_phot
    3636
    3737!==========================================================================
  • trunk/LMDZ.MARS/libf/phymars/nlte_calc.F

    r2616 r3012  
    2828      subroutine MZESC110 (ig,nl_cts_real, nzy_cts_real,ierr,varerr)
    2929c***********************************************************************
     30      use nlte_tcool_mod, only: errors
    3031      implicit none
    3132
     
    710711
    711712c***********************************************************************
    712 
     713      use nlte_tcool_mod, only: errors
    713714      implicit none
    714715
     
    12571258      subroutine MZESC121
    12581259c***********************************************************************
    1259 
     1260      use nlte_tcool_mod, only: errors
    12601261      implicit none
    12611262
  • trunk/LMDZ.MARS/libf/phymars/nlte_tcool.F

    r2398 r3012  
     1      MODULE nlte_tcool_mod
     2     
     3      IMPLICIT NONE
     4     
     5      CONTAINS
    16!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    27! Fast scheme for NLTE cooling rates at 15um by CO2 in a Martian GCM !
     
    3540      include 'nlte_paramdef.h'
    3641      include 'nlte_commons.h'
    37       include "chimiedata.h"
    3842
    3943
     
    195199      enddo
    196200c     end subroutine
    197       return
    198       end
     201
     202      end subroutine nlte_tcool
    199203
    200204
     
    232236
    233237c     functions
    234       external  hrkday_convert
    235       real      hrkday_convert
     238!      external         hrkday_convert
     239!      real     hrkday_convert
    236240
    237241c***********************************************************************
     
    472476
    473477c     end
    474       return
    475       end
     478
     479      end subroutine NLTEdlvr11_ZGRID
    476480
    477481
     
    12531257
    12541258c     final
    1255       return
     1259
    12561260c     
    1257       end
     1261      end subroutine NLTEdlvr11_CZALU
    12581262
    12591263
     
    13401344
    13411345c
    1342       return
    1343       end
     1346      end subroutine NLTEdlvr11_FB626CTS
    13441347
    13451348
     
    13721375     
    13731376c     end                                           
    1374       return                                 
    1375       end        
     1377
     1378      end function hrkday_convert
    13761379
    13771380
     
    16351638      call abort_physic("nlte_tcool",
    16361639     &     'Stopped in NLTE scheme due to severe error.',1)
    1637       end
     1640
     1641      end subroutine ERRORS
     1642
     1643      END MODULE nlte_tcool_mod
  • trunk/LMDZ.MARS/libf/phymars/nltecool.F

    r2586 r3012  
     1      MODULE nltecool_mod
     2     
     3      IMPLICIT NONE
     4     
     5      CONTAINS
    16c**************************************************************************
    27c
     
    3439      implicit none
    3540
    36 #include "nltedata.h" ! (Equivalent to the reading of the "nlte_escape.dat" file)
    37 #include "chimiedata.h"
    38 #include "callkeys.h"
     41      include "nltedata.h"
     42      include "callkeys.h"
    3943
    4044c Input and output variables
    4145c
    42       integer     ngrid                           ! no. of horiz. gridpoints
    43       integer     nlayer                          ! no. of atmospheric layers
    44       integer     nq                              ! no. of tracers
    45       real        pplay(ngrid,nlayer)             ! input pressure grid
    46       real        pt(ngrid,nlayer)                ! input temperatures
    47       real        pq(ngrid,nlayer,nq)                ! input mmrs
    48       real        dtnlte(ngrid,nlayer)            ! output temp. tendencies
     46      integer,intent(in) :: ngrid               ! no. of horiz. gridpoints
     47      integer,intent(in) :: nlayer              ! no. of atmospheric layers
     48      integer,intent(in) :: nq                  ! no. of tracers
     49      real,intent(in) :: pplay(ngrid,nlayer)    ! input pressure grid (Pa)
     50      real,intent(in) :: pt(ngrid,nlayer)       ! input temperatures (K)
     51      real,intent(in) :: pq(ngrid,nlayer,nq)    ! input mmrs (kg/kg_air)
     52      real,intent(out) :: dtnlte(ngrid,nlayer)  ! output temp. tendencies (K/s)
    4953
    5054c
     
    261265c     close(7)
    262266c
    263         return
    264         end
     267
     268      end subroutine nltecool
    265269
    266270c***********************************************************************
     
    290294         endif
    291295      enddo
    292       return
    293       end
     296
     297      end subroutine interp1
    294298
    295299c***********************************************************************
     
    324328         endif
    325329      enddo
    326       return
    327       end
     330
     331      end subroutine interp3
     332
     333      END MODULE nltecool_mod
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r3006 r3012  
    2424      use calcstormfract_mod, only: calcstormfract
    2525      use topmons_mod, only: topmons,topmons_setup
     26      use nltecool_mod, only: nltecool
     27      use nlte_tcool_mod, only: nlte_tcool
    2628      use tracer_mod, only: noms, mmol, igcm_co2, igcm_n2, igcm_co2_ice,
    2729     &                      igcm_co, igcm_o, igcm_h2o_vap, igcm_h2o_ice,
Note: See TracChangeset for help on using the changeset viewer.