Changeset 2945 for LMDZ5


Ignore:
Timestamp:
Jul 12, 2017, 4:20:24 PM (7 years ago)
Author:
jbmadeleine
Message:
  • Added a new output called rneblsvol which is the cloud fraction by volume

computed in the thermals (see cloudth_vert in cloudth_mod.F90)

  • Added an option called iflag_rain_incloud_vol that computes the conversion

of cloud water to rain using the cloud fraction by volume instead of the cloud
fraction by area, which is larger and otherwise erroneously reduces the in-cloud
water content; iflag_rain_incloud_vol can only be used for iflag_cloudth_vert>=3

Location:
LMDZ5/trunk
Files:
14 edited

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/DefLists/field_def_lmdz.xml

    r2928 r2945  
    515515        <field id="rnebcon"    long_name="Convective Cloud Fraction"    unit="-" />
    516516        <field id="rnebls"    long_name="LS Cloud fraction"    unit="-" />
     517        <field id="rneblsvol" long_name="LS Cloud fraction by volume"    unit="-" />
    517518        <field id="rhum"    long_name="Relative humidity"    unit="-" />
    518519        <field id="ozone"    long_name="Ozone mole fraction"    unit="-" />
  • LMDZ5/trunk/DefLists/file_def_histLES_lmdz.xml

    r2928 r2945  
    474474                <field field_ref="rnebcon" level="10" />
    475475                <field field_ref="rnebls" level="10" />
     476                <field field_ref="rneblsvol" level="10" />
    476477                <field field_ref="rhum" level="10" />
    477478                <field field_ref="ozone" level="10" />
  • LMDZ5/trunk/DefLists/file_def_histday_lmdz.xml

    r2928 r2945  
    481481                <field field_ref="rnebcon" level="5" />
    482482                <field field_ref="rnebls" level="5" />
     483                <field field_ref="rneblsvol" level="5" />
    483484                <field field_ref="rhum" level="5" />
    484485                <field field_ref="ozone" level="10" />
  • LMDZ5/trunk/DefLists/file_def_histhf_lmdz.xml

    r2928 r2945  
    504504                <field field_ref="rnebcon" level="10" />
    505505                <field field_ref="rnebls" level="10" />
     506                <field field_ref="rneblsvol" level="10" />
    506507                <field field_ref="rhum" level="10" />
    507508                <field field_ref="ozone" level="10" />
  • LMDZ5/trunk/DefLists/file_def_histins_lmdz.xml

    r2928 r2945  
    474474                <field field_ref="rnebcon" level="10" />
    475475                <field field_ref="rnebls" level="10" />
     476                <field field_ref="rneblsvol" level="10" />
    476477                <field field_ref="rhum" level="10" />
    477478                <field field_ref="ozone" level="10" />
  • LMDZ5/trunk/DefLists/file_def_histmth_lmdz.xml

    r2928 r2945  
    514514                <field field_ref="rnebcon" level="2" />
    515515                <field field_ref="rnebls" level="2" />
     516                <field field_ref="rneblsvol" level="2" />
    516517                <field field_ref="rhum" level="2" />
    517518                <field field_ref="ozone" level="2" />
  • LMDZ5/trunk/DefLists/file_def_histstn_lmdz.xml

    r2928 r2945  
    474474                <field field_ref="rnebcon" level="10" />
    475475                <field field_ref="rnebls" level="10" />
     476                <field field_ref="rneblsvol" level="10" />
    476477                <field field_ref="rhum" level="10" />
    477478                <field field_ref="ozone" level="10" />
  • LMDZ5/trunk/libf/phylmd/cloudth_mod.F90

    r2911 r2945  
    587587       SUBROUTINE cloudth_v3(ngrid,klev,ind2,  &
    588588     &           ztv,po,zqta,fraca, &
    589      &           qcloud,ctot,zpspsk,paprs,ztla,zthl, &
     589     &           qcloud,ctot,ctot_vol,zpspsk,paprs,ztla,zthl, &
    590590     &           ratqs,zqs,t)
    591591
     
    624624      REAL zqsatenv(ngrid,klev)
    625625     
    626      
    627       REAL sigma1(ngrid,klev)
     626      REAL sigma1(ngrid,klev)                                                         
    628627      REAL sigma2(ngrid,klev)
    629628      REAL qlth(ngrid,klev)
    630629      REAL qlenv(ngrid,klev)
    631630      REAL qltot(ngrid,klev)
    632       REAL cth(ngrid,klev) 
     631      REAL cth(ngrid,klev)
    633632      REAL cenv(ngrid,klev)   
    634633      REAL ctot(ngrid,klev)
    635       REAL rneb(ngrid,klev)
     634      REAL cth_vol(ngrid,klev)
     635      REAL cenv_vol(ngrid,klev)
     636      REAL ctot_vol(ngrid,klev)
     637      REAL rneb(ngrid,klev)     
    636638      REAL t(ngrid,klev)
    637639      REAL qsatmmussig1,qsatmmussig2,sqrt2pi,sqrt2,sqrtpi,pi
     
    652654      CALL cloudth_vert_v3(ngrid,klev,ind2,  &
    653655     &           ztv,po,zqta,fraca, &
    654      &           qcloud,ctot,zpspsk,paprs,ztla,zthl, &
     656     &           qcloud,ctot,ctot_vol,zpspsk,paprs,ztla,zthl, &
    655657     &           ratqs,zqs,t)
    656658      RETURN
     
    672674      cenv(:,:)=0.
    673675      ctot(:,:)=0.
     676      cth_vol(:,:)=0.
     677      cenv_vol(:,:)=0.
     678      ctot_vol(:,:)=0.
    674679      qsatmmussig1=0.
    675680      qsatmmussig2=0.
     
    747752      cenv(ind1,ind2)=0.5*(1.+1.*erf(xenv))     !4.18 p 111, l.7 p115 & 4.20 p 119 thesis Arnaud Jam
    748753      ctot(ind1,ind2)=fraca(ind1,ind2)*cth(ind1,ind2)+(1.-1.*fraca(ind1,ind2))*cenv(ind1,ind2)
     754      ctot_vol(ind1,ind2)=ctot(ind1,ind2)
    749755
    750756      qlth(ind1,ind2)=sigma2s*((exp(-1.*xth**2)/sqrt2pi)+xth*sqrt2*cth(ind1,ind2))
     
    780786      xenv=senv/(sqrt2*sigma1s)
    781787      ctot(ind1,ind2)=0.5*(1.+1.*erf(xenv))
     788      ctot_vol(ind1,ind2)=ctot(ind1,ind2)
    782789      qltot(ind1,ind2)=sigma1s*((exp(-1.*xenv**2)/sqrt2pi)+xenv*sqrt2*cenv(ind1,ind2))
    783790
     
    801808     SUBROUTINE cloudth_vert_v3(ngrid,klev,ind2,  &
    802809     &           ztv,po,zqta,fraca, &
    803      &           qcloud,ctot,zpspsk,paprs,ztla,zthl, &
     810     &           qcloud,ctot,ctot_vol,zpspsk,paprs,ztla,zthl, &
    804811     &           ratqs,zqs,t)
    805812
     
    838845      REAL zqsatenv(ngrid,klev)
    839846     
    840      
     847   
     848
     849
    841850      REAL sigma1(ngrid,klev)                                                         
    842851      REAL sigma2(ngrid,klev)
     
    844853      REAL qlenv(ngrid,klev)
    845854      REAL qltot(ngrid,klev)
    846       REAL cth(ngrid,klev) 
     855      REAL cth(ngrid,klev)
    847856      REAL cenv(ngrid,klev)   
    848857      REAL ctot(ngrid,klev)
     858      REAL cth_vol(ngrid,klev)
     859      REAL cenv_vol(ngrid,klev)
     860      REAL ctot_vol(ngrid,klev)
    849861      REAL rneb(ngrid,klev)
    850862      REAL t(ngrid,klev)                                                                 
     
    854866      REAL sth,senv,sigma1s,sigma2s,xth,xenv,exp_xenv1,exp_xenv2,exp_xth1,exp_xth2
    855867      REAL xth1,xth2,xenv1,xenv2,deltasth, deltasenv
    856       REAL IntJ,IntI1,IntI2,IntI3,coeffqlenv,coeffqlth
     868      REAL IntJ,IntI1,IntI2,IntI3,IntJ_CF,IntI1_CF,IntI3_CF,coeffqlenv,coeffqlth
    857869      REAL Tbef,zdelta,qsatbef,zcor
    858870      REAL qlbef 
     
    883895      cenv(:,:)=0.
    884896      ctot(:,:)=0.
     897      cth_vol(:,:)=0.
     898      cenv_vol(:,:)=0.
     899      ctot_vol(:,:)=0.
    885900      qsatmmussig1=0.
    886901      qsatmmussig2=0.
     
    10321047      exp_xth2 = exp(-1.*xth2**2)
    10331048     
     1049      !CF_surfacique
    10341050      cth(ind1,ind2)=0.5*(1.-1.*erf(xth1))
    10351051      cenv(ind1,ind2)=0.5*(1.-1.*erf(xenv1))
    1036       ctot(ind1,ind2)=fraca(ind1,ind2)*cth(ind1,ind2)+(1.-1.*fraca(ind1,ind2))*cenv(ind1,ind2)
    1037 
    1038      
     1052      ctot(ind1,ind2)=fraca(ind1,ind2)*cth(ind1,ind2)+(1.-1.*fraca(ind1,ind2))*cenv(ind1,ind2)
     1053
     1054
     1055      !CF_volumique & eau condense
    10391056      !environnement
    10401057      IntJ=0.5*senv*(1-erf(xenv2))+(sigma1s/sqrt2pi)*exp_xenv2
     1058      IntJ_CF=0.5*(1.-1.*erf(xenv2))
    10411059      if (deltasenv .lt. 1.e-10) then
    10421060      qlenv(ind1,ind2)=IntJ
     1061      cenv_vol(ind1,ind2)=IntJ_CF
    10431062      else
    10441063      IntI1=(((senv+deltasenv)**2+(sigma1s)**2)/(8*deltasenv))*(erf(xenv2)-erf(xenv1))
    10451064      IntI2=(sigma1s**2/(4*deltasenv*sqrtpi))*(xenv1*exp_xenv1-xenv2*exp_xenv2)
    10461065      IntI3=((sqrt2*sigma1s*(senv+deltasenv))/(4*sqrtpi*deltasenv))*(exp_xenv1-exp_xenv2)
     1066      IntI1_CF=((senv+deltasenv)*(erf(xenv2)-erf(xenv1)))/(4*deltasenv)
     1067      IntI3_CF=(sqrt2*sigma1s*(exp_xenv1-exp_xenv2))/(4*sqrtpi*deltasenv)
    10471068      qlenv(ind1,ind2)=IntJ+IntI1+IntI2+IntI3
     1069      cenv_vol(ind1,ind2)=IntJ_CF+IntI1_CF+IntI3_CF
    10481070      endif
    1049 
    10501071
    10511072      !thermique
    10521073      IntJ=0.5*sth*(1-erf(xth2))+(sigma2s/sqrt2pi)*exp_xth2
    1053       if (deltasth .lt. 1.e-10) then ! Seuil a choisir !!!
     1074      IntJ_CF=0.5*(1.-1.*erf(xth2))
     1075      if (deltasth .lt. 1.e-10) then
    10541076      qlth(ind1,ind2)=IntJ
     1077      cth_vol(ind1,ind2)=IntJ_CF
    10551078      else
    10561079      IntI1=(((sth+deltasth)**2+(sigma2s)**2)/(8*deltasth))*(erf(xth2)-erf(xth1))
    10571080      IntI2=(sigma2s**2/(4*deltasth*sqrtpi))*(xth1*exp_xth1-xth2*exp_xth2)
    10581081      IntI3=((sqrt2*sigma2s*(sth+deltasth))/(4*sqrtpi*deltasth))*(exp_xth1-exp_xth2)
     1082      IntI1_CF=((sth+deltasth)*(erf(xth2)-erf(xth1)))/(4*deltasth)
     1083      IntI3_CF=(sqrt2*sigma2s*(exp_xth1-exp_xth2))/(4*sqrtpi*deltasth)
    10591084      qlth(ind1,ind2)=IntJ+IntI1+IntI2+IntI3
     1085      cth_vol(ind1,ind2)=IntJ_CF+IntI1_CF+IntI3_CF
    10601086      endif
    10611087
     1088
    10621089      qltot(ind1,ind2)=fraca(ind1,ind2)*qlth(ind1,ind2)+(1.-1.*fraca(ind1,ind2))*qlenv(ind1,ind2)
    1063 
     1090      ctot_vol(ind1,ind2)=fraca(ind1,ind2)*cth_vol(ind1,ind2)+(1.-1.*fraca(ind1,ind2))*cenv_vol(ind1,ind2)
    10641091
    10651092      ENDIF ! of if (iflag_cloudth_vert==1 or 3 or 4)
     1093
    10661094
    10671095      if (cenv(ind1,ind2).lt.1.e-10.or.cth(ind1,ind2).lt.1.e-10) then
    10681096      ctot(ind1,ind2)=0.
     1097      ctot_vol(ind1,ind2)=0.
    10691098      qcloud(ind1)=zqsatenv(ind1,ind2)
    10701099
     
    11001129      xenv=senv/(sqrt2*sigma1s)
    11011130      ctot(ind1,ind2)=0.5*(1.+1.*erf(xenv))
     1131      ctot_vol(ind1,ind2)=ctot(ind1,ind2)
    11021132      qltot(ind1,ind2)=sigma1s*((exp(-1.*xenv**2)/sqrt2pi)+xenv*sqrt2*cenv(ind1,ind2))
    11031133     
  • LMDZ5/trunk/libf/phylmd/conf_phys_m.F90

    r2915 r2945  
    169169    INTEGER,SAVE :: iflag_t_glace_omp
    170170    INTEGER,SAVE :: iflag_cloudth_vert_omp
     171    INTEGER,SAVE :: iflag_rain_incloud_vol_omp
    171172    REAL,SAVE :: rad_froid_omp, rad_chau1_omp, rad_chau2_omp
    172173    REAL,SAVE :: t_glace_min_omp, t_glace_max_omp
     
    12481249
    12491250    !
     1251    !Config Key  = iflag_rain_incloud_vol
     1252    !Config Desc = 
     1253    !Config Def  = 0
     1254    !Config Help =
     1255    !
     1256    iflag_rain_incloud_vol_omp = 0
     1257    CALL getin('iflag_rain_incloud_vol',iflag_rain_incloud_vol_omp)
     1258
     1259    !
    12501260    !Config Key  = iflag_ice_thermo
    12511261    !Config Desc = 
     
    21632173    iflag_t_glace = iflag_t_glace_omp
    21642174    iflag_cloudth_vert=iflag_cloudth_vert_omp
     2175    iflag_rain_incloud_vol=iflag_rain_incloud_vol_omp
    21652176    iflag_ice_thermo = iflag_ice_thermo_omp
    21662177    rei_min = rei_min_omp
     
    25112522    write(lunout,*)' iflag_t_glace = ',iflag_t_glace
    25122523    write(lunout,*)' iflag_cloudth_vert = ',iflag_cloudth_vert
     2524    write(lunout,*)' iflag_rain_incloud_vol = ',iflag_rain_incloud_vol
    25132525    write(lunout,*)' iflag_ice_thermo = ',iflag_ice_thermo
    25142526    write(lunout,*)' rei_min = ',rei_min
  • LMDZ5/trunk/libf/phylmd/fisrtilp.F90

    r2923 r2945  
    1818  USE ioipsl_getin_p_mod, ONLY : getin_p
    1919  USE phys_local_var_mod, ONLY: ql_seri,qs_seri
     20  USE phys_local_var_mod, ONLY: rneblsvol
    2021  ! flag to include modifications to ensure energy conservation (if flag >0)
    2122  USE add_phys_tend_mod, only : fl_cor_ebil
     
    132133  INTEGER ncoreczq 
    133134  REAL ctot(klon,klev)
     135  REAL ctot_vol(klon,klev)
    134136  REAL zqs(klon), zdqs(klon), zdelta, zcor, zcvm5 
    135137  REAL zdqsdT_raw(klon)
     
    217219!  ice_thermo = iflag_ice_thermo .GE. 1
    218220
    219 itap=itap+1
    220 znebprecip(:)=0.
    221 
    222    ice_thermo = (iflag_ice_thermo .EQ. 1).OR.(iflag_ice_thermo .GE. 3)
     221  itap=itap+1
     222  znebprecip(:)=0.
     223
     224  ice_thermo = (iflag_ice_thermo .EQ. 1).OR.(iflag_ice_thermo .GE. 3)
    223225  zdelq=0.0
     226  ctot_vol(1:klon,1:klev)=0.0
     227  rneblsvol(1:klon,1:klev)=0.0
    224228
    225229  if (prt_level>9)write(lunout,*)'NUAGES4 A. JAM'
     
    727731               call cloudth_v3(klon,klev,k,ztv, &
    728732                   zq,zqta,fraca, &
    729                    qcloud,ctot,zpspsk,paprs,ztla,zthl, &
     733                   qcloud,ctot,ctot_vol,zpspsk,paprs,ztla,zthl, &
    730734                   ratqs,zqs,t)
    731735              endif
    732736              do i=1,klon
    733737                 rneb(i,k)=ctot(i,k)
     738                 rneblsvol(i,k)=ctot_vol(i,k)
    734739                 zqn(i)=qcloud(i)
    735740              enddo
     
    11481153                         *fallvs(zrhol(i)) * zfice(i)
    11491154                 endif
    1150                  zchau    = zct   *dtime/REAL(ninter) * zoliq(i) &
    1151                       *(1.0-EXP(-(zoliq(i)/zneb(i)/zcl   )**2)) *(1.-zfice(i))
     1155
     1156                 ! si l'heterogeneite verticale est active, on utilise
     1157                 ! la fraction volumique "vraie" plutot que la fraction
     1158                 ! surfacique modifiee, qui est plus grande et reduit
     1159                 ! sinon l'eau in-cloud de facon artificielle
     1160                 if ((iflag_cloudth_vert>=3).AND.(iflag_rain_incloud_vol==1)) then
     1161                    zchau    = zct   *dtime/REAL(ninter) * zoliq(i) &
     1162                         *(1.0-EXP(-(zoliq(i)/ctot_vol(i,k)/zcl   )**2)) *(1.-zfice(i))
     1163                 else
     1164                    zchau    = zct   *dtime/REAL(ninter) * zoliq(i) &
     1165                         *(1.0-EXP(-(zoliq(i)/zneb(i)/zcl   )**2)) *(1.-zfice(i))
     1166                 endif
    11521167!AJ<
    11531168                 IF (.NOT. ice_thermo) THEN
  • LMDZ5/trunk/libf/phylmd/nuage.h

    r2606 r2945  
    1010
    1111      INTEGER iflag_t_glace, iflag_cloudth_vert, iflag_cld_cv
     12      INTEGER iflag_rain_incloud_vol
    1213
    1314      common /nuagecom/ rad_froid,rad_chau1, rad_chau2,t_glace_max,     &
     
    1516     &                  tau_cld_cv,coefw_cld_cv,                        &
    1617     &                  tmax_fonte_cv,                                  &
    17      &                  iflag_t_glace,iflag_cloudth_vert,iflag_cld_cv
     18     &                  iflag_t_glace,iflag_cloudth_vert,iflag_cld_cv,  &
     19     &                  iflag_rain_incloud_vol
    1820!$OMP THREADPRIVATE(/nuagecom/)
  • LMDZ5/trunk/libf/phylmd/phys_local_var_mod.F90

    r2898 r2945  
    417417      REAL,ALLOCATABLE,SAVE,DIMENSION(:,:) :: beta_prec
    418418!$OMP THREADPRIVATE(beta_prec)
    419       REAL, ALLOCATABLE, SAVE, DIMENSION(:,:) :: rneb,rnebjn
    420 !$OMP THREADPRIVATE(rneb,rnebjn)
     419      REAL, ALLOCATABLE, SAVE, DIMENSION(:,:) :: rneb,rnebjn,rneblsvol
     420!$OMP THREADPRIVATE(rneb,rnebjn,rneblsvol)
    421421
    422422! variables de sorties MM
     
    732732      ALLOCATE(wdtrainA(klon,klev),wdtrainM(klon,klev))
    733733      ALLOCATE(beta_prec(klon,klev))
    734       ALLOCATE(rneb(klon,klev),rnebjn(klon,klev))
     734      ALLOCATE(rneb(klon,klev),rnebjn(klon,klev),rneblsvol(klon,klev))
    735735
    736736
  • LMDZ5/trunk/libf/phylmd/phys_output_ctrlout_mod.F90

    r2876 r2945  
    13281328  TYPE(ctrl_out), SAVE :: o_rnebls = ctrl_out((/ 2, 5, 10, 10, 10, 10, 11, 11, 11, 11/), &
    13291329    'rnebls', 'LS Cloud fraction', '-', (/ ('', i=1, 10) /))
     1330  TYPE(ctrl_out), SAVE :: o_rneblsvol = ctrl_out((/ 2, 5, 10, 10, 10, 10, 11, 11, 11, 11/), &
     1331    'rneblsvol', 'LS Cloud fraction by volume', '-', (/ ('', i=1, 10) /))
    13301332  TYPE(ctrl_out), SAVE :: o_rhum = ctrl_out((/ 2, 5, 10, 10, 10, 10, 11, 11, 11, 11/), &
    13311333    'rhum', 'Relative humidity', '-', (/ ('', i=1, 10) /))
  • LMDZ5/trunk/libf/phylmd/phys_output_write_mod.F90

    r2854 r2945  
    124124         o_vitu, o_vitv, o_vitw, o_pres, o_paprs, &
    125125         o_zfull, o_zhalf, o_rneb, o_rnebjn, o_rnebcon, &
    126          o_rnebls, o_rhum, o_ozone, o_ozone_light, &
     126         o_rnebls, o_rneblsvol, o_rhum, o_ozone, o_ozone_light, &
    127127         o_duphy, o_dtphy, o_dqphy, o_dqphy2d, o_dqlphy, o_dqlphy2d, &
    128128         o_dqsphy, o_dqsphy2d, o_albe_srf, o_z0m_srf, o_z0h_srf, &
     
    269269         ql_seri, qs_seri, tr_seri, &
    270270         zphi, u_seri, v_seri, omega, cldfra, &
    271          rneb, rnebjn, zx_rh, d_t_dyn,  &
     271         rneb, rnebjn, rneblsvol, zx_rh, d_t_dyn,  &
    272272         d_q_dyn,  d_ql_dyn, d_qs_dyn, &
    273273         d_q_dyn2d,  d_ql_dyn2d, d_qs_dyn2d, &
     
    13441344       CALL histwrite_phy(o_rnebcon, rnebcon)
    13451345       CALL histwrite_phy(o_rnebls, rneb)
     1346       CALL histwrite_phy(o_rneblsvol, rneblsvol)
    13461347       IF (vars_defined)  THEN
    13471348          DO k=1, klev
Note: See TracChangeset for help on using the changeset viewer.