Ignore:
Timestamp:
Sep 23, 2013, 9:56:47 AM (11 years ago)
Author:
emillour
Message:

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/aeronomars/calchim.F90

    r1036 r1047  
    1       subroutine calchim(nq,                                                &
     1      subroutine calchim(ngrid,nlayer,nq,                                   &
    22                         ptimestep,pplay,pplev,pt,pdt,dist_sol,mu0,         &
    33                         zzlev,zzlay,zday,pq,pdq,dqchim,dqschim,dqcloud,    &
     
    1313                            igcm_noplus, igcm_n2plus, igcm_hplus,         &
    1414                            igcm_hco2plus, igcm_elec, mmol
     15      use conc_mod, only: mmean ! mean molecular mass of the atmosphere
     16
    1517      implicit none
    1618
     
    3840!
    3941!    ptimestep                  timestep (s)
    40 !    pplay(ngridmx,nlayermx)    Pressure at the middle of the layers (Pa)
    41 !    pplev(ngridmx,nlayermx+1)  Intermediate pressure levels (Pa)
    42 !    pt(ngridmx,nlayermx)       Temperature (K)
    43 !    pdt(ngridmx,nlayermx)      Temperature tendency (K)
    44 !    pu(ngridmx,nlayermx)       u component of the wind (ms-1)
    45 !    pdu(ngridmx,nlayermx)      u component tendency (K)
    46 !    pv(ngridmx,nlayermx)       v component of the wind (ms-1)
    47 !    pdv(ngridmx,nlayermx)      v component tendency (K)
     42!    pplay(ngrid,nlayer)    Pressure at the middle of the layers (Pa)
     43!    pplev(ngrid,nlayer+1)  Intermediate pressure levels (Pa)
     44!    pt(ngrid,nlayer)       Temperature (K)
     45!    pdt(ngrid,nlayer)      Temperature tendency (K)
     46!    pu(ngrid,nlayer)       u component of the wind (ms-1)
     47!    pdu(ngrid,nlayer)      u component tendency (K)
     48!    pv(ngrid,nlayer)       v component of the wind (ms-1)
     49!    pdv(ngrid,nlayer)      v component tendency (K)
    4850!    dist_sol                   distance of the sun (AU)
    49 !    mu0(ngridmx)               cos of solar zenith angle (=1 when sun at zenith)
    50 !    pq(ngridmx,nlayermx,nqmx)  Advected fields, ie chemical species here
    51 !    pdq(ngridmx,nlayermx,nqmx) Previous tendencies on pq
    52 !    tauref(ngridmx)            Optical depth at 7 hPa
    53 !    co2ice(ngridmx)            co2 ice surface layer (kg.m-2)
    54 !    surfdust(ngridmx,nlayermx) dust surface area (m2/m3)
    55 !    surfice(ngridmx,nlayermx)  ice surface area (m2/m3)
     51!    mu0(ngrid)               cos of solar zenith angle (=1 when sun at zenith)
     52!    pq(ngrid,nlayer,nq)  Advected fields, ie chemical species here
     53!    pdq(ngrid,nlayer,nq) Previous tendencies on pq
     54!    tauref(ngrid)            Optical depth at 7 hPa
     55!    co2ice(ngrid)            co2 ice surface layer (kg.m-2)
     56!    surfdust(ngrid,nlayer) dust surface area (m2/m3)
     57!    surfice(ngrid,nlayer)  ice surface area (m2/m3)
    5658!
    5759!  Output:
    5860!
    59 !    dqchim(ngridmx,nlayermx,nqmx) ! tendencies on pq due to chemistry
    60 !    dqschim(ngridmx,nqmx)         ! tendencies on qsurf
    61 !
    62 !=======================================================================
    63 
    64 #include "dimensions.h"
    65 #include "dimphys.h"
     61!    dqchim(ngrid,nlayer,nq) ! tendencies on pq due to chemistry
     62!    dqschim(ngrid,nq)         ! tendencies on qsurf
     63!
     64!=======================================================================
     65
     66!#include "dimensions.h"
     67!#include "dimphys.h"
    6668#include "chimiedata.h"
    6769!#include "tracer.h"
    6870#include "comcstfi.h"
    6971#include "callkeys.h"
    70 #include "conc.h"
     72!#include "conc.h"
    7173
    7274!     input:
    7375
     76      integer,intent(in) :: ngrid ! number of atmospheric columns
     77      integer,intent(in) :: nlayer ! number of atmospheric layers
    7478      integer,intent(in) :: nq ! number of tracers
    7579      real :: ptimestep
    76       real :: pplay(ngridmx,nlayermx)    ! pressure at the middle of the layers
    77       real :: zzlay(ngridmx,nlayermx)    ! pressure at the middle of the layers
    78       real :: pplev(ngridmx,nlayermx+1)  ! intermediate pressure levels
    79       real :: zzlev(ngridmx,nlayermx+1)  ! altitude at layer boundaries
    80       real :: pt(ngridmx,nlayermx)       ! temperature
    81       real :: pdt(ngridmx,nlayermx)      ! temperature tendency
    82       real :: pu(ngridmx,nlayermx)       ! u component of the wind (m.s-1)
    83       real :: pdu(ngridmx,nlayermx)      ! u component tendency
    84       real :: pv(ngridmx,nlayermx)       ! v component of the wind (m.s-1)
    85       real :: pdv(ngridmx,nlayermx)      ! v component tendency
     80      real :: pplay(ngrid,nlayer)    ! pressure at the middle of the layers
     81      real :: zzlay(ngrid,nlayer)    ! pressure at the middle of the layers
     82      real :: pplev(ngrid,nlayer+1)  ! intermediate pressure levels
     83      real :: zzlev(ngrid,nlayer+1)  ! altitude at layer boundaries
     84      real :: pt(ngrid,nlayer)       ! temperature
     85      real :: pdt(ngrid,nlayer)      ! temperature tendency
     86      real :: pu(ngrid,nlayer)       ! u component of the wind (m.s-1)
     87      real :: pdu(ngrid,nlayer)      ! u component tendency
     88      real :: pv(ngrid,nlayer)       ! v component of the wind (m.s-1)
     89      real :: pdv(ngrid,nlayer)      ! v component tendency
    8690      real :: dist_sol                   ! distance of the sun (AU)
    87       real :: mu0(ngridmx)               ! cos of solar zenith angle (=1 when sun at zenith)
    88       real :: pq(ngridmx,nlayermx,nq)  ! tracers mass mixing ratio
    89       real :: pdq(ngridmx,nlayermx,nq) ! previous tendencies
     91      real :: mu0(ngrid)               ! cos of solar zenith angle (=1 when sun at zenith)
     92      real :: pq(ngrid,nlayer,nq)  ! tracers mass mixing ratio
     93      real :: pdq(ngrid,nlayer,nq) ! previous tendencies
    9094      real :: zday                       ! date (time since Ls=0, in martian days)
    91       real :: tauref(ngridmx)            ! optical depth at 7 hPa
    92       real :: co2ice(ngridmx)            ! co2 ice surface layer (kg.m-2)
    93       real :: surfdust(ngridmx,nlayermx) ! dust surface area (m2/m3)
    94       real :: surfice(ngridmx,nlayermx)  !  ice surface area (m2/m3)
     95      real :: tauref(ngrid)            ! optical depth at 7 hPa
     96      real :: co2ice(ngrid)            ! co2 ice surface layer (kg.m-2)
     97      real :: surfdust(ngrid,nlayer) ! dust surface area (m2/m3)
     98      real :: surfice(ngrid,nlayer)  !  ice surface area (m2/m3)
    9599
    96100!     output:
    97101
    98       real :: dqchim(ngridmx,nlayermx,nq) ! tendencies on pq due to chemistry
    99       real :: dqschim(ngridmx,nq)         ! tendencies on qsurf
    100       real :: dqcloud(ngridmx,nlayermx,nq)! tendencies on pq due to condensation
    101       real :: dqscloud(ngridmx,nq)        ! tendencies on qsurf
     102      real :: dqchim(ngrid,nlayer,nq) ! tendencies on pq due to chemistry
     103      real :: dqschim(ngrid,nq)         ! tendencies on qsurf
     104      real :: dqcloud(ngrid,nlayer,nq)! tendencies on pq due to condensation
     105      real :: dqscloud(ngrid,nq)        ! tendencies on qsurf
    102106
    103107!     local variables:
     
    143147
    144148      real    :: latvl1, lonvl1
    145       real    :: zq(ngridmx,nlayermx,nq) ! pq+pdq*ptimestep before chemistry
     149      real    :: zq(ngrid,nlayer,nq) ! pq+pdq*ptimestep before chemistry
    146150                                           ! new mole fraction after
    147       real    :: zt(ngridmx,nlayermx)      ! temperature
    148       real    :: zu(ngridmx,nlayermx)      ! u component of the wind
    149       real    :: zv(ngridmx,nlayermx)      ! v component of the wind
     151      real    :: zt(ngrid,nlayer)      ! temperature
     152      real    :: zu(ngrid,nlayer)      ! u component of the wind
     153      real    :: zv(ngrid,nlayer)      ! v component of the wind
    150154      real    :: taucol                    ! optical depth at 7 hPa
    151155
     
    155159!     for each column of atmosphere:
    156160
    157       real :: zpress(nlayermx)       !  Pressure (mbar)
    158       real :: zdens(nlayermx)        !  Density  (cm-3)
    159       real :: ztemp(nlayermx)        !  Temperature (K)
    160       real :: zlocal(nlayermx)       !  Altitude (km)
    161       real :: zycol(nlayermx,nq)   !  Composition (mole fractions)
     161      real :: zpress(nlayer)       !  Pressure (mbar)
     162      real :: zdens(nlayer)        !  Density  (cm-3)
     163      real :: ztemp(nlayer)        !  Temperature (K)
     164      real :: zlocal(nlayer)       !  Altitude (km)
     165      real :: zycol(nlayer,nq)   !  Composition (mole fractions)
    162166      real :: szacol                 !  Solar zenith angle
    163       real :: surfice1d(nlayermx)    !  Ice surface area (cm2/cm3)
    164       real :: surfdust1d(nlayermx)   !  Dust surface area (cm2/cm3)
    165       real :: jo3(nlayermx)          !  Photodissociation rate O3->O1D (s-1)
     167      real :: surfice1d(nlayer)    !  Ice surface area (cm2/cm3)
     168      real :: surfdust1d(nlayer)   !  Dust surface area (cm2/cm3)
     169      real :: jo3(nlayer)          !  Photodissociation rate O3->O1D (s-1)
    166170
    167171!     for output:
     
    169173      logical :: output                 ! to issue calls to writediagfi and stats
    170174      parameter (output = .true.)
    171       real :: jo3_3d(ngridmx,nlayermx)  ! Photodissociation rate O3->O1D (s-1)
     175      real :: jo3_3d(ngrid,nlayer)  ! Photodissociation rate O3->O1D (s-1)
    172176
    173177!=======================================================================
     
    592596!=======================================================================
    593597
    594       do ig = 1,ngridmx
     598      do ig = 1,ngrid
    595599         
    596600         foundswitch = 0
    597          do l = 1,nlayermx
     601         do l = 1,nlayer
    598602            do i = 1,nbq
    599603               iq = niq(i) ! get tracer index
     
    626630            end if
    627631            if (.not. thermochem) then
    628                lswitch = min(50,nlayermx+1)
     632               lswitch = min(50,nlayer+1)
    629633            end if
    630634
    631          end do ! of do l=1,nlayermx
     635         end do ! of do l=1,nlayer
    632636
    633637         szacol = acos(mu0(ig))*180./pi
     
    647651!        ozone photolysis, for output
    648652
    649             do l = 1,nlayermx
     653            do l = 1,nlayer
    650654               jo3_3d(ig,l) = jo3(l)
    651655            end do
     
    653657!        condensation of h2o2
    654658
    655             call perosat(ig,ptimestep,pplev,pplay,                 &
     659            call perosat(ngrid, nlayer, nq,                        &
     660                         ig,ptimestep,pplev,pplay,                 &
    656661                         ztemp,zycol,dqcloud,dqscloud)
    657662         end if
     
    667672
    668673         if (depos) then
    669             call deposition(ig, ig_vl1, pplay, pplev, zzlay, zzlev,&
     674            call deposition(ngrid, nlayer, nq,                      &
     675                            ig, ig_vl1, pplay, pplev, zzlay, zzlev, &
    670676                            zu, zv, zt, zycol, ptimestep, co2ice)
    671677         end if
     
    679685
    680686!     tendency for the most abundant species = - sum of others
    681          do l = 1,nlayermx
     687         do l = 1,nlayer
    682688            iloc=maxloc(zycol(l,:))
    683689            iqmax=iloc(1)
     
    691697               end if
    692698            end do
    693          end do ! of do l = 1,nlayermx
     699         end do ! of do l = 1,nlayer
    694700
    695701!=======================================================================
     
    697703!=======================================================================
    698704
    699       end do ! of do ig=1,ngridmx
     705      end do ! of do ig=1,ngrid
    700706
    701707!=======================================================================
     
    707713
    708714      if (photochem .and. output) then
    709          if (ngridmx > 1) then
    710             call writediagfi(ngridmx,'jo3','j o3->o1d',    &
     715         if (ngrid > 1) then
     716            call writediagfi(ngrid,'jo3','j o3->o1d',    &
    711717                             's-1',3,jo3_3d(1,1))
    712718           if (callstats) then
    713               call wstats(ngridmx,'jo3','j o3->o1d',       &
     719              call wstats(ngrid,'jo3','j o3->o1d',       &
    714720                          's-1',3,jo3_3d(1,1))
    715721           endif
    716          end if ! of if (ngridmx.gt.1)
     722         end if ! of if (ngrid.gt.1)
    717723      end if ! of if (output)
    718724
Note: See TracChangeset for help on using the changeset viewer.