Ignore:
Timestamp:
May 6, 2016, 12:30:29 PM (9 years ago)
Author:
emillour
Message:

All GCMs:
Further adaptations to keep up with changes in LMDZ5 concerning
physics/dynamics separation (up to rev r2420 of LMDZ5)

  • all physics packages:
  • added module callphysiq_mod.F90 in dynphy_lonlat/phy* which contains the routine "call_physiq" which is called by calfis* and calls the physics. This way different "physiq" routine from different physics packages may be called: The calfis* routines now exposes all available fields that might be transmitted to physiq but which is actually send (ie: expected/needed by physiq) is decided in call_physiq.
  • turned "physiq.F[90]" into module "physiq_mod.F[90]" for better control of "physiq" arguments. for phyvenus/phytitan, extracted gr_fi_ecrit from physiq.F as gr_fi_ecrit.F90 (note that it can only work in serial).
  • misc:
  • updated wxios.F90 to keep up with LMDZ5 modifications.
  • dyn3d_common:
  • infotrac.F90 keep up with LMDZ5 modifications (cosmetics)
  • dyn3d:
  • gcm.F90: cosmetic cleanup.
  • leapfrog.F90: fix computation of date as function of itau.
  • dyn3dpar:
  • gcm.F: cosmetic cleanup.
  • leapfrog_p.F90: fix computation of date as function of itau.

NB: physics are given the date corresponding to the end of the
physics step.

  • dynphy_lonlat:
  • calfis.F : added computation of relative wind vorticity.
  • calfis_p.F: added computation of relative wind vorticity (input required by Earth physics)

EM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.COMMON/libf/dynphy_lonlat/calfis.F

    r1459 r1549  
    3333      USE write_field
    3434      USE cpdet_mod, only: t2tpot,tpot2t
    35  
     35#ifdef CPP_PHYS
     36      USE callphysiq_mod, ONLY: call_physiq
     37#endif
     38
    3639! used only for zonal averages
    3740      USE moyzon_mod
     
    145148      REAL zphi(ngridmx,llm),zphis(ngridmx)
    146149
     150      REAL zrot(iip1,jjm,llm) ! AdlC May 2014
    147151      REAL zufi(ngridmx,llm), zvfi(ngridmx,llm)
     152      REAL zrfi(ngridmx,llm) ! relative wind vorticity
    148153      REAL ztfi(ngridmx,llm),zqfi(ngridmx,llm,nqtot)
    149154! ADAPTATION GCM POUR CP(T)
     
    572577
    573578
     579C  Alvaro de la Camara (May 2014)
     580C  46.1 Calcul de la vorticite et passage sur la grille physique
     581C  --------------------------------------------------------------
     582      DO l=1,llm
     583        do i=1,iim
     584          do j=1,jjm
     585            zrot(i,j,l) = (pvcov(i+1,j,l) - pvcov(i,j,l)
     586     $                   + pucov(i,j+1,l) - pucov(i,j,l))
     587     $                   / (cu(i,j)+cu(i,j+1))
     588     $                   / (cv(i+1,j)+cv(i,j)) *4
     589          enddo
     590        enddo
     591      ENDDO
     592
    574593c   46.champ v:
    575594c   -----------
     
    584603c     $         ( pdvcov(i,j-1,l)/cv(i,j-1) + pdvcov(i,j,l)/cv(i,j) )
    585604            ENDDO
     605               zrfi(ig0 + 1,l)= 0.25 *(zrot(iim,j-1,l)+zrot(iim,j,l)
     606     &                                +zrot(1,j-1,l)+zrot(1,j,l))
     607            DO i=2,iim
     608               zrfi(ig0 + i,l)= 0.25 *(zrot(i-1,j-1,l)+zrot(i-1,j,l)
     609     $                   +zrot(i,j-1,l)+zrot(i,j,l))   !  AdlC MAY 2014
     610            ENDDO
    586611         ENDDO
    587612      ENDDO
     
    613638         zvfi(1,l)  = SSUM(iim,zsin,1)/pi
    614639!         pcvgv(1,l) = SSUM(iim,zsinbis,1)/pi
    615 
     640         zrfi(1, l) = 0.
    616641      ENDDO
    617642
     
    642667         zvfi(ngridmx,l)  = SSUM(iim,zsin,1)/pi
    643668!         pcvgv(ngridmx,l) = SSUM(iim,zsinbis,1)/pi
    644 
     669         zrfi(ngridmx, l) = 0.
    645670      ENDDO
    646671c
     
    674699         lafin_split=lafin.and.isplit==nsplit_phys
    675700
    676       if (planet_type.eq."earth") then
    677          CALL physiq (ngridmx,
    678      .             llm,
    679      .             debut_split,
    680      .             lafin_split,
    681      .             jD_cur,
    682      .             jH_cur_split,
    683      .             zdt_split,
    684      .             zplev,
    685      .             zplay,
    686      .             zphi,
    687      .             zphis,
    688      .             presnivs,
    689      .             zufi,
    690      .             zvfi,
    691      .             ztfi,
    692      .             zqfi,
    693      .             flxwfi,
    694      .             zdufi,
    695      .             zdvfi,
    696      .             zdtfi,
    697      .             zdqfi,
    698      .             zdpsrf,
    699      .             pducov)
    700 
    701       else if ( planet_type=="generic" ) then
    702 
    703          CALL physiq (ngridmx,     !! ngrid
    704      .             llm,            !! nlayer
    705      .             nqtot,          !! nq
    706      .             tname,          !! tracer names from dynamical core (given in infotrac)
    707      .             debut_split,    !! firstcall
    708      .             lafin_split,    !! lastcall
    709      .             jD_cur,         !! pday. see leapfrog
    710      .             jH_cur_split,   !! ptime "fraction of day"
    711      .             zdt_split,      !! ptimestep
    712      .             zplev,          !! pplev
    713      .             zplay,          !! pplay
    714      .             zphi,           !! pphi
    715      .             zufi,           !! pu
    716      .             zvfi,           !! pv
    717      .             ztfi,           !! pt
    718      .             zqfi,           !! pq
    719      .             flxwfi,         !! pw !! or 0. anyway this is for diagnostic. not used in physiq.
    720      .             zdufi,          !! pdu
    721      .             zdvfi,          !! pdv
    722      .             zdtfi,          !! pdt
    723      .             zdqfi,          !! pdq
    724      .             zdpsrf,         !! pdpsrf
    725      .             tracerdyn)      !! tracerdyn <-- utilite ???
    726 
    727       else if ( planet_type=="mars" ) then
    728 
    729         CALL physiq (ngridmx,    ! ngrid
    730      .             llm,          ! nlayer
    731      .             nqtot,        ! nq
    732      .             debut_split,  ! firstcall
    733      .             lafin_split,  ! lastcall
    734      .             jD_cur,       ! pday
    735      .             jH_cur_split, ! ptime
    736      .             zdt_split,    ! ptimestep
    737      .             zplev,        ! pplev
    738      .             zplay,        ! pplay
    739      .             zphi,         ! pphi
    740      .             zufi,         ! pu
    741      .             zvfi,         ! pv
    742      .             ztfi,         ! pt
    743      .             zqfi,         ! pq
    744      .             flxwfi,       ! pw
    745      .             zdufi,        ! pdu
    746      .             zdvfi,        ! pdv
    747      .             zdtfi,        ! pdt
    748      .             zdqfi,        ! pdq
    749      .             zdpsrf,       ! pdpsrf
    750      .             tracerdyn)    ! tracerdyn (somewhat obsolete)
    751 
    752       else if ((planet_type=="titan").or.(planet_type=="venus")) then
    753 
    754          CALL physiq (ngridmx,
    755      .             llm,
    756      .             nqtot,
    757      .             debut_split,
    758      .             lafin_split,
    759      .             jD_cur,
    760      .             jH_cur_split,
    761      .             zdt_split,
    762      .             zplev,
    763      .             zplay,
    764      .             zpk,
    765      .             zphi,
    766      .             zphis,
    767      .             presnivs,
    768      .             zufi,
    769      .             zvfi,
    770      .             ztfi,
    771      .             zqfi,
    772      .             flxwfi,
    773      .             zdufi,
    774      .             zdvfi,
    775      .             zdtfi,
    776      .             zdqfi,
    777      .             zdpsrf)
    778 
    779       else ! unknown "planet_type"
    780 
    781         write(lunout,*) "calfis_p: error, unknown planet_type: ",
    782      &                  trim(planet_type)
    783         stop
    784 
    785       endif ! planet_type
     701        CALL call_physiq(ngridmx,llm,nqtot,tname,
     702     &                   debut_split,lafin_split,
     703     &                   jD_cur,jH_cur_split,zdt_split,
     704     &                   zplev,zplay,
     705     &                   zpk,zphi,zphis,
     706     &                   presnivs,
     707     &                   zufi,zvfi,zrfi,ztfi,zqfi,
     708     &                   flxwfi,pducov,
     709     &                   zdufi,zdvfi,zdtfi,zdqfi,zdpsrf,
     710     &                   tracerdyn)
     711
     712!      if (planet_type.eq."earth") then
     713!         CALL physiq (ngridmx,
     714!     .             llm,
     715!     .             debut_split,
     716!     .             lafin_split,
     717!     .             jD_cur,
     718!     .             jH_cur_split,
     719!     .             zdt_split,
     720!     .             zplev,
     721!     .             zplay,
     722!     .             zphi,
     723!     .             zphis,
     724!     .             presnivs,
     725!     .             zufi,
     726!     .             zvfi,
     727!     .             ztfi,
     728!     .             zqfi,
     729!     .             flxwfi,
     730!     .             zdufi,
     731!     .             zdvfi,
     732!     .             zdtfi,
     733!     .             zdqfi,
     734!     .             zdpsrf,
     735!     .             pducov)
     736!
     737!      else if ( planet_type=="generic" ) then
     738!
     739!         CALL physiq (ngridmx,     !! ngrid
     740!     .             llm,            !! nlayer
     741!     .             nqtot,          !! nq
     742!     .             tname,          !! tracer names from dynamical core (given in infotrac)
     743!     .             debut_split,    !! firstcall
     744!     .             lafin_split,    !! lastcall
     745!     .             jD_cur,         !! pday. see leapfrog
     746!     .             jH_cur_split,   !! ptime "fraction of day"
     747!     .             zdt_split,      !! ptimestep
     748!     .             zplev,          !! pplev
     749!     .             zplay,          !! pplay
     750!     .             zphi,           !! pphi
     751!     .             zufi,           !! pu
     752!     .             zvfi,           !! pv
     753!     .             ztfi,           !! pt
     754!     .             zqfi,           !! pq
     755!     .             flxwfi,         !! pw !! or 0. anyway this is for diagnostic. not used in physiq.
     756!     .             zdufi,          !! pdu
     757!     .             zdvfi,          !! pdv
     758!     .             zdtfi,          !! pdt
     759!     .             zdqfi,          !! pdq
     760!     .             zdpsrf,         !! pdpsrf
     761!     .             tracerdyn)      !! tracerdyn <-- utilite ???
     762!
     763!      else if ( planet_type=="mars" ) then
     764!
     765!        CALL physiq (ngridmx,    ! ngrid
     766!     .             llm,          ! nlayer
     767!     .             nqtot,        ! nq
     768!     .             debut_split,  ! firstcall
     769!     .             lafin_split,  ! lastcall
     770!     .             jD_cur,       ! pday
     771!     .             jH_cur_split, ! ptime
     772!     .             zdt_split,    ! ptimestep
     773!     .             zplev,        ! pplev
     774!     .             zplay,        ! pplay
     775!     .             zphi,         ! pphi
     776!     .             zufi,         ! pu
     777!     .             zvfi,         ! pv
     778!     .             ztfi,         ! pt
     779!     .             zqfi,         ! pq
     780!     .             flxwfi,       ! pw
     781!     .             zdufi,        ! pdu
     782!     .             zdvfi,        ! pdv
     783!     .             zdtfi,        ! pdt
     784!     .             zdqfi,        ! pdq
     785!     .             zdpsrf,       ! pdpsrf
     786!     .             tracerdyn)    ! tracerdyn (somewhat obsolete)
     787!
     788!      else if ((planet_type=="titan").or.(planet_type=="venus")) then
     789!
     790!         CALL physiq (ngridmx,
     791!     .             llm,
     792!     .             nqtot,
     793!     .             debut_split,
     794!     .             lafin_split,
     795!     .             jD_cur,
     796!     .             jH_cur_split,
     797!     .             zdt_split,
     798!     .             zplev,
     799!     .             zplay,
     800!     .             zpk,
     801!     .             zphi,
     802!     .             zphis,
     803!     .             presnivs,
     804!     .             zufi,
     805!     .             zvfi,
     806!     .             ztfi,
     807!     .             zqfi,
     808!     .             flxwfi,
     809!     .             zdufi,
     810!     .             zdvfi,
     811!     .             zdtfi,
     812!     .             zdqfi,
     813!     .             zdpsrf)
     814!
     815!      else ! unknown "planet_type"
     816!
     817!        write(lunout,*) "calfis_p: error, unknown planet_type: ",
     818!     &                  trim(planet_type)
     819!        stop
     820!
     821!      endif ! planet_type
    786822
    787823         zufi(:,:)=zufi(:,:)+zdufi(:,:)*zdt_split
Note: See TracChangeset for help on using the changeset viewer.