Changeset 2060


Ignore:
Timestamp:
Jan 9, 2019, 3:25:30 PM (6 years ago)
Author:
aboissinot
Message:

Add the thermal plume model (cf. Rio et al. 2010) extended to gas giant.
Specific parameters are set in thermcell_mod.

Location:
trunk/LMDZ.GENERIC
Files:
12 added
5 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/README

    r2056 r2060  
    14191419== 07/01/2018 == AB
    14201420- Planck step function is replaced by a piecewise linear function in gfluxi.F and sfluxi.F
    1421 in the computation of B1, B0, PLTOP, PLANCKSUM and BSURF.
     1421  in the computation of B1, B0, PLTOP, PLANCKSUM and BSURF.
     1422
     1423== 09/01/2018 == AB
     1424- Add the thermal plume model (cf. Rio et al. 2010) extended to gas giant. Specific parameters are set in thermcell_mod.
     1425
  • trunk/LMDZ.GENERIC/libf/phystd/callkeys_mod.F90

    r1797 r2060  
    44      logical,save :: callrad,corrk,calldifv,UseTurbDiff
    55!$OMP THREADPRIVATE(callrad,corrk,calldifv,UseTurbDiff)
    6       logical,save :: calladj,co2cond,callsoil
    7 !$OMP THREADPRIVATE(calladj,co2cond,callsoil)
     6      logical,save :: calladj,calltherm,co2cond,callsoil
     7!$OMP THREADPRIVATE(calladj,calltherm,co2cond,callsoil)
    88      logical,save :: season,diurnal,tlocked,rings_shadow,lwrite
    99!$OMP THREADPRIVATE(season,diurnal,tlocked,rings_shadow,lwrite)
  • trunk/LMDZ.GENERIC/libf/phystd/inifis_mod.F90

    r1832 r2060  
    263263     call getin_p("calladj",calladj)
    264264     write(*,*) " calladj = ",calladj
     265
     266     write(*,*) "call thermal plume model ?"
     267     calltherm=.false. ! default value
     268     call getin_p("calltherm",calltherm)
     269     write(*,*) " calltherm = ",calltherm
    265270
    266271     write(*,*) "call CO2 condensation ?"
  • trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90

    r2058 r2060  
    1515 
    1616      use radinc_h, only : L_NSPECTI,L_NSPECTV,naerkind
    17       use watercommon_h, only : RLVTT, Psat_water,epsi,su_watercycle
     17      use watercommon_h, only : RLVTT, Psat_water,epsi,su_watercycle, RV, T_h2o_ice_liq
     18      use thermcell_mod, only: nbsrf, init_thermcell_mod
    1819      use gases_h, only: gnom, gfrac
    1920      use radcommon_h, only: sigma, glat, grav, BWNV
     
    99100!      III. Vertical diffusion (turbulent mixing) :
    100101!
    101 !      IV. Dry Convective adjustment :
     102!      IV. Convection :
     103!         IV.a Thermal plume model
     104!         IV.b Dry convective adjusment
    102105!
    103106!      V. Condensation and sublimation of gases (currently just CO2).
     
    248251      real zzlev(ngrid,nlayer+1)     ! Altitude at the atmospheric layer boundaries.
    249252
     253! VARIABLES for the thermal plume model
     254     
     255      INTEGER lmin(ngrid)                 ! Plume base level
     256      INTEGER lmix(ngrid)                 ! Plume maximal velocity level
     257      INTEGER lalim(ngrid)                ! Plume alimentation top level
     258      INTEGER lmax(ngrid)                 ! Maximal level reached by the plume
     259     
     260! AB : Integers used only in alp for diagnoses
     261      INTEGER lalim_conv(ngrid)
     262      INTEGER zlcl(ngrid)
     263     
     264      real f0(ngrid)                      ! Mass flux norm
     265      real fm0(ngrid, nlayer+1)           ! Mass flux
     266      real entr0(ngrid, nlayer)           ! Entrainment
     267      real detr0(ngrid, nlayer)           ! Detrainment
     268      real zmax0(ngrid)                   ! Plume height
     269      real fraca(ngrid, nlayer+1)         ! Fraction of the surface that plumes occupies
     270      real zw2(ngrid, nlayer+1)           ! Squared vertical speed or its square root
     271      real zqsatth(ngrid, nlayer)         ! Water vapor mixing ratio at saturation
     272      real zqta(ngrid, nlayer)            ! Total water mass mixing ratio in the plume
     273      real zqla(ngrid, nlayer)            ! Liquid water mass mixing ratio in the plume
     274      real ztv(ngrid, nlayer)             ! Virtual potential temperature in the environment considering large scale condensation
     275      real ztva(ngrid, nlayer)            ! Virtual potential temperature in the plume
     276      real zthl(ngrid, nlayer)            ! Potential temperature in the environment without considering condensation
     277      real ztla(ngrid, nlayer)            ! Potential temperature in the plume
     278      real ratqscth(ngrid, nlayer)        !
     279      real ratqsdiff(ngrid, nlayer)       !
     280     
     281! AB : Reals only used in alp for diagnoses
     282      real Ale_bl(ngrid), Alp_bl(ngrid)
     283      real therm_tke_max0(ngrid), env_tke_max0(ngrid)
     284      real n2(ngrid), s2(ngrid)
     285      real ale_bl_stat(ngrid)
     286      real therm_tke_max(ngrid, nlayer), env_tke_max(ngrid, nlayer)
     287      real alp_bl_det(ngrid), alp_bl_fluct_m(ngrid), alp_bl_fluct_tke(ngrid), alp_bl_conv(ngrid), alp_bl_stat(ngrid)
     288      real wght_th(ngrid, nlayer)
     289      real pbl_tke(ngrid,nlayer+1,nbsrf)
     290      real pctsrf(ngrid, nbsrf)
     291! AB : omega already defined, do we have to fusion them ?
     292      real omega_therm(ngrid, nlayer)
     293      real airephy(ngrid)
     294      real w0(ngrid)
     295      real w_conv(ngrid)
     296      real fraca0(ngrid)
     297     
     298! AB : variables used for outputs
     299      REAL pmax(ngrid)                       ! Maximal pressure reached by the plume
     300      REAL pmin(ngrid)                       ! Minimal pressure reached by the plume
     301
    250302! TENDENCIES due to various processes :
    251303
     
    260312      real dtlscale(ngrid,nlayer)                             ! Largescale routine.
    261313      real zdtc(ngrid,nlayer)                                 ! Condense_co2 routine.
    262       real zdtdif(ngrid,nlayer)                                      ! Turbdiff/vdifc routines.
     314      real zdtdif(ngrid,nlayer)                               ! Turbdiff/vdifc routines.
     315      real zdttherm(ngrid,nlayer)                             ! Calltherm routine.
    263316      real zdtmr(ngrid,nlayer)                                ! Mass_redistribution routine.
    264317      real zdtrain(ngrid,nlayer)                              ! Rain routine.
     
    282335      real zdqdif(ngrid,nlayer,nq)    ! Turbdiff/vdifc routines.
    283336      real zdqevap(ngrid,nlayer)      ! Turbdiff routine.
     337      real zdotherm(ngrid,nlayer)     ! Calltherm routine.
    284338      real zdqsed(ngrid,nlayer,nq)    ! Callsedim routine.
    285339      real zdqmr(ngrid,nlayer,nq)     ! Mass_redistribution routine.
     
    295349                       
    296350      ! For Winds : (m/s/s)
    297       real zdvadj(ngrid,nlayer),zduadj(ngrid,nlayer) ! Convadj routine.
    298       real zdumr(ngrid,nlayer),zdvmr(ngrid,nlayer)   ! Mass_redistribution routine.
    299       real zdvdif(ngrid,nlayer),zdudif(ngrid,nlayer) ! Turbdiff/vdifc routines.
    300       real zdhdif(ngrid,nlayer)                      ! Turbdiff/vdifc routines.
    301       real zdhadj(ngrid,nlayer)                      ! Convadj routine.
     351      real zdvadj(ngrid,nlayer), zduadj(ngrid,nlayer)       ! Convadj routine.
     352      real zdutherm(ngrid,nlayer), zdvtherm(ngrid,nlayer)   ! Calltherm routine.
     353      real zdumr(ngrid,nlayer), zdvmr(ngrid,nlayer)         ! Mass_redistribution routine.
     354      real zdvdif(ngrid,nlayer), zdudif(ngrid,nlayer)       ! Turbdiff/vdifc routines.
     355      real zdhdif(ngrid,nlayer)                             ! Turbdiff/vdifc routines.
     356      real zdhadj(ngrid,nlayer)                             ! Convadj routine.
    302357     
    303358      ! For Pressure and Mass :
     
    644699         endif
    645700         
     701!        Set some parameters for the thermal plume model
     702!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
     703         if (calltherm) then
     704            call init_thermcell_mod(g, rcp, r, pi, T_h2o_ice_liq, RV)
     705         endif
     706         
    646707         call su_watercycle ! even if we don't have a water cycle, we might
    647708                            ! need epsi for the wvp definitions in callcorrk.F
     709                            ! or RETV, RLvCp for the thermal plume model
    648710#ifndef MESOSCALE
    649711         if (ngrid.ne.1) then ! Note : no need to create a restart file in 1d.
     
    11141176
    11151177
    1116 !----------------------------------
    1117 !   IV. Dry convective adjustment :
    1118 !----------------------------------
     1178!-------------------
     1179!   IV. Convection :
     1180!-------------------
     1181     
     1182! ~~~~~~~~~~~~~~~~~~~~~~~~~~
     1183! IV.a Thermal plume model :
     1184! ~~~~~~~~~~~~~~~~~~~~~~~~~~
     1185     
     1186      IF (calltherm) THEN
     1187! AB:  WARNING: if a plume stops, the parametrization never look above if somewhere the atmosphere is still unstable!
     1188!               As is, there cannot be more than one plume by grid point by time step.
     1189         CALL thermcell_main(icount, ngrid, nlayer, ptimestep,                   &
     1190                             pplay, pplev, pphi, firstcall,                      &
     1191                             pu, pv, pt, pq(:,:,igcm_h2o_vap),                   &
     1192                             zdutherm, zdvtherm, zdttherm, zdotherm,             &
     1193                             f0, fm0, entr0, detr0,                              &
     1194                             zqta, zqla, ztv, ztva, ztla, zthl, zqsatth,         &
     1195                             zmax0, zw2, fraca,                                  &
     1196                             lmin, lmix, lalim, lmax,                            &
     1197                             zpopsk, ratqscth, ratqsdiff,                        &
     1198! AB : next variables are only used for diagnoses
     1199                             Ale_bl,Alp_bl,lalim_conv,wght_th,                   &
     1200                             pbl_tke,pctsrf,omega_therm,airephy,                 &
     1201                             zlcl,fraca0,w0,w_conv,                              &
     1202                             therm_tke_max0,env_tke_max0,                        &
     1203                             n2,s2,ale_bl_stat,                                  &
     1204                             therm_tke_max,env_tke_max,                          &
     1205                             alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke,         &
     1206                             alp_bl_conv,alp_bl_stat)
     1207         
     1208         DO ig=1,ngrid
     1209            pmin(ig) = pplev(ig,lmin(ig))
     1210            pmax(ig) = pplev(ig,lmax(ig))
     1211         ENDDO
     1212         
     1213         pdu(1:ngrid,1:nlayer) = pdu(1:ngrid,1:nlayer) + zdutherm(1:ngrid,1:nlayer)
     1214         pdv(1:ngrid,1:nlayer) = pdv(1:ngrid,1:nlayer) + zdvtherm(1:ngrid,1:nlayer)
     1215         pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer) + zdttherm(1:ngrid,1:nlayer)
     1216         
     1217         IF(tracer) THEN
     1218            pdq(1:ngrid,1:nlayer,igcm_h2o_vap) = pdq(1:ngrid,1:nlayer,igcm_h2o_vap) + zdotherm(1:ngrid,1:nlayer)
     1219         ENDIF
     1220         
     1221      ENDIF ! end of 'calltherm'
     1222     
     1223! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
     1224! IV.b Dry convective adjustment :
     1225! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    11191226
    11201227      if(calladj) then
     
    11691276         
    11701277      endif ! end of 'calladj'
    1171 
     1278     
    11721279!-----------------------------------------------
    11731280!   V. Carbon dioxide condensation-sublimation :
     
    20232130      endif
    20242131     
    2025 
     2132      ! Thermal plume model
     2133      if (calltherm) then
     2134         call writediagfi(ngrid,'entr','Entrainment','kg m$^{-2}$ s$^{-1}$', 3, entr0)
     2135         call writediagfi(ngrid,'detr','Detrainment','kg m$^{-2}$ s$^{-1}$', 3, detr0)
     2136         call writediagfi(ngrid,'fm','Mass flux','kg m$^{-2}$ s$^{-1}$', 3, fm0)
     2137         call writediagfi(ngrid,'w_plm','Squared vertical velocity','m s$^{-1}$', 3, zw2)
     2138         call writediagfi(ngrid,'fraca','Updraft fraction','', 3, fraca)
     2139!         call writediagfi(ngrid,'pmin', 'pmin', 'Pa', 2, pmin)
     2140!         call writediagfi(ngrid,'pmax', 'pmax', 'Pa', 2, pmax)
     2141      endif
     2142     
    20262143      ! Total energy balance diagnostics
    20272144      if(callrad.and.(.not.newtonian))then
     
    22082325      CALL send_xios_field("ps",ps)
    22092326      CALL send_xios_field("area",cell_area)
    2210 
     2327     
    22112328      CALL send_xios_field("temperature",zt)
    22122329      CALL send_xios_field("u",zu)
    22132330      CALL send_xios_field("v",zv)
    22142331      CALL send_xios_field("omega",omega)
    2215 
     2332     
    22162333      CALL send_xios_field("ISR",fluxtop_dn)
    22172334      CALL send_xios_field("OLR",fluxtop_lw)
    2218 
     2335     
    22192336      if (lastcall.and.is_omp_master) then
    22202337        write(*,*) "physiq: call xios_context_finalize"
     
    22222339      endif
    22232340#endif
    2224 
     2341     
    22252342      icount=icount+1
    2226 
    2227     end subroutine physiq
    2228    
    2229     end module physiq_mod
     2343     
     2344      end subroutine physiq
     2345     
     2346   end module physiq_mod
  • trunk/LMDZ.GENERIC/libf/phystd/watercommon_h.F90

    r1993 r2060  
    1919
    2020      real, save :: epsi, RCPD, RCPV, RV, RVTMP2
     21      real, save :: RETV
     22      real, save :: RLvCp
    2123!$OMP THREADPRIVATE(epsi,RCPD,RCPV,RV,RVTMP2)
    2224     
     
    2628      subroutine su_watercycle
    2729
    28       use comcstfi_mod, only: cpp, mugaz
     30      use comcstfi_mod, only: r, cpp, mugaz
    2931      implicit none
    3032
     
    5254
    5355      RVTMP2 = RCPV/RCPD-1. ! not currently used...
    54 
     56     
     57! AB : initializations added for the thermal plume model
     58      RETV = RV / r - 1.
     59      RLvCp = RLVTT / RCPD
     60     
    5561      end subroutine su_watercycle
    5662     
Note: See TracChangeset for help on using the changeset viewer.