subroutine rain_generic(ngrid,nlayer,nq,ptimestep,pplev,pplay,pphi,t,pdt,pq,pdq,d_t,dq_rain_generic_vap,dq_rain_generic_cld,dqsrain_generic,dqssnow_generic,reevap_precip,rneb) use ioipsl_getin_p_mod, only: getin_p use generic_cloud_common_h use watercommon_h, only: T_h2O_ice_liq,T_h2O_ice_clouds,rhowater ! T_h2O_ice_clouds,rhowater are only used for precip_scheme_generic >=2 use radii_mod, only: h2o_cloudrad ! only used for precip_scheme_generic >=2 use tracer_h use comcstfi_mod, only: g, r, cpp use generic_tracer_index_mod, only: generic_tracer_index implicit none !================================================================== ! ! Purpose ! ------- ! Calculates precipitation for generic condensable tracers, using simplified microphysics. ! ! GCS : generic condensable specie ! ! Authors ! ------- ! Adapted from rain.F90 by Noé Clément (2022) ! !================================================================== ! Arguments integer,intent(in) :: ngrid ! number of atmospheric columns integer,intent(in) :: nlayer ! number of atmospheric layers integer,intent(in) :: nq ! number of tracers real,intent(in) :: ptimestep ! time interval real,intent(in) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa) real,intent(in) :: pplay(ngrid,nlayer) ! mid-layer pressure (Pa) real,intent(in) :: pphi(ngrid,nlayer) ! mid-layer geopotential real,intent(in) :: t(ngrid,nlayer) ! input temperature (K) real,intent(in) :: pdt(ngrid,nlayer) ! input tendency on temperature (K/s) real,intent(in) :: pq(ngrid,nlayer,nq) ! tracers (kg/kg) real,intent(in) :: pdq(ngrid,nlayer,nq) ! input tendency on tracers real,intent(out) :: d_t(ngrid,nlayer) ! temperature tendency (K/s) real,intent(out) :: dq_rain_generic_vap(ngrid,nlayer,nq) ! tendency of GCS tracers precipitation (kg/kg.s-1) - vapor real,intent(out) :: dq_rain_generic_cld(ngrid,nlayer,nq) ! tendency of GCS tracers precipitation (kg/kg.s-1) - cloud real,intent(out) :: dqsrain_generic(ngrid,nq) ! rain flux at the surface (kg.m-2.s-1) real,intent(out) :: dqssnow_generic(ngrid,nq) ! snow flux at the surface (kg.m-2.s-1) real,intent(out) :: reevap_precip(ngrid,nq) ! re-evaporation flux of precipitation integrated over the atmospheric column (kg.m-2.s-1) real,intent(in) :: rneb(ngrid,nlayer) ! cloud fraction real zt(ngrid,nlayer) ! working temperature (K) real ql(ngrid,nlayer) ! liquid GCS (Kg/Kg) real q(ngrid,nlayer) ! specific humidity (Kg/Kg) real d_q(ngrid,nlayer) ! GCS vapor increment real d_ql(ngrid,nlayer) ! liquid GCS / ice increment integer igcm_generic_vap, igcm_generic_ice ! index of the vap and ice fo GCS real, save :: RCPD ! equal to cpp real, save :: metallicity ! metallicity of planet --- is not used here, but necessary to call function Psat_generic !$OMP THREADPRIVATE(metallicity) ! Subroutine options real,parameter :: seuil_neb=0.001 ! Nebulosity threshold integer,save :: precip_scheme_generic ! id number for precipitaion scheme ! for simple scheme (precip_scheme_generic=1) real,save :: rainthreshold_generic ! Precipitation threshold in simple scheme ! for sundquist scheme (precip_scheme_generic=2-3) real,save :: cloud_sat_generic ! Precipitation threshold in non simple scheme real,save :: precip_timescale_generic ! Precipitation timescale ! for Boucher scheme (precip_scheme_generic=4) real,save :: Cboucher_generic ! Precipitation constant in Boucher 95 scheme real,parameter :: Kboucher=1.19E8 real,save :: c1 !$OMP THREADPRIVATE(precip_scheme_generic,rainthreshold_generic,cloud_sat_generic,precip_timescale_generic,Cboucher_generic,c1) integer,parameter :: ninter=5 ! only used for precip_scheme_generic >=2 logical,save :: evap_prec_generic ! Does the rain evaporate ? REAL,SAVE :: evap_coeff_generic ! multiplication evaporation constant. 1. gives the model of gregory et al. !$OMP THREADPRIVATE(evap_prec_generic,evap_coeff_generic) ! for simple scheme : precip_scheme_generic=1 real lconvert ! Local variables integer i, k, n, iq REAL zqs(ngrid,nlayer), dqsat(ngrid,nlayer), dlnpsat(ngrid,nlayer), Tsat(ngrid,nlayer), zdelta, zcor REAL precip_rate(ngrid), precip_rate_tmp(ngrid) ! local precipitation rate in kg of condensed GCS per m^2 per s. REAL zqev, zqevt real zoliq(ngrid) real zdz(ngrid),zrho(ngrid),ztot(ngrid), zrhol(ngrid) real zchau(ngrid),zfroi(ngrid),zfrac(ngrid),zneb(ngrid) real reffh2oliq(ngrid,nlayer),reffh2oice(ngrid,nlayer) real t_tmp, p_tmp, psat_tmp real tnext(ngrid,nlayer) real dmass(ngrid,nlayer) ! mass of air in each grid cell real dWtot ! Indices of GCS vapour and GCS ice tracers integer, save :: i_vap_generic=0 ! GCS vapour integer, save :: i_ice_generic=0 ! GCS ice !$OMP THREADPRIVATE(i_vap_generic,i_ice_generic) LOGICAL,save :: firstcall=.true. !$OMP THREADPRIVATE(firstcall) ! to call only one time the ice/vap pair of a tracer logical call_ice_vap_generic ! Online functions real fallv, fall2v, zzz ! falling speed of ice crystals fallv (zzz) = 3.29 * ((zzz)**0.16) fall2v (zzz) =10.6 * ((zzz)**0.31) !for use with radii ! Let's loop on tracers dq_rain_generic_vap(1:ngrid,1:nlayer,1:nq) = 0.0 dq_rain_generic_cld(1:ngrid,1:nlayer,1:nq) = 0.0 dqsrain_generic(:,:) = 0.0 do iq=1,nq call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic) if(call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer m=constants_mass(iq) delta_vapH=constants_delta_vapH(iq) Tref=constants_Tref(iq) Pref=constants_Pref(iq) epsi_generic=constants_epsi_generic(iq) RLVTT_generic=constants_RLVTT_generic(iq) RCPD = cpp if (firstcall) then metallicity=0.0 ! default value --- is not used here but necessary to call function Psat_generic call getin_p("metallicity",metallicity) ! --- is not used here but necessary to call function Psat_generic i_vap_generic=igcm_generic_vap i_ice_generic=igcm_generic_ice write(*,*) "rain: i_ice_generic=",i_ice_generic write(*,*) " i_vap_generic=",i_vap_generic write(*,*) "re-evaporate precipitations?" evap_prec_generic=.true. ! default value call getin_p("evap_prec_generic",evap_prec_generic) write(*,*) " evap_prec_generic = ",evap_prec_generic if (evap_prec_generic) then write(*,*) "multiplicative constant in reevaporation" evap_coeff_generic=1. ! default value call getin_p("evap_coeff_generic",evap_coeff_generic) write(*,*) " evap_coeff_generic = ",evap_coeff_generic end if write(*,*) "Precipitation scheme to use?" precip_scheme_generic=1 ! default value call getin_p("precip_scheme_generic",precip_scheme_generic) write(*,*) " precip_scheme_generic = ",precip_scheme_generic if (precip_scheme_generic.eq.1) then write(*,*) "rainthreshold_generic in simple scheme?" rainthreshold_generic=0. ! default value call getin_p("rainthreshold_generic",rainthreshold_generic) write(*,*) " rainthreshold_generic = ",rainthreshold_generic !else ! write(*,*) "precip_scheme_generic = ", precip_scheme_generic ! write(*,*) "this precip_scheme_generic has not been implemented yet," ! write(*,*) "only precip_scheme_generic = 1 has been implemented" else if (precip_scheme_generic.eq.2.or.precip_scheme_generic.eq.3) then write(*,*) "cloud GCS saturation level in non simple scheme?" cloud_sat_generic=2.6e-4 ! default value call getin_p("cloud_sat_generic",cloud_sat_generic) write(*,*) " cloud_sat_generic = ",cloud_sat_generic write(*,*) "precipitation timescale in non simple scheme?" precip_timescale_generic=3600. ! default value call getin_p("precip_timescale_generic",precip_timescale_generic) write(*,*) " precip_timescale_generic = ",precip_timescale_generic else if (precip_scheme_generic.eq.4) then write(*,*) "multiplicative constant in Boucher 95 precip scheme" Cboucher_generic=1. ! default value call getin_p("Cboucher_generic",Cboucher_generic) write(*,*) " Cboucher_generic = ",Cboucher_generic c1=1.00*1.097/rhowater*Cboucher_generic*Kboucher endif PRINT*, 'in rain_generic.F, ninter=', ninter firstcall = .false. endif ! of if (firstcall) ! GCM -----> subroutine variables do k = 1, nlayer do i = 1, ngrid zt(i,k) = t(i,k)+pdt(i,k)*ptimestep ! a big fat bug was here q(i,k) = pq(i,k,i_vap_generic)+pdq(i,k,i_vap_generic)*ptimestep ql(i,k) = pq(i,k,i_ice_generic)+pdq(i,k,i_ice_generic)*ptimestep !q(i,k) = pq(i,k,i_vap_generic)!+pdq(i,k,i_vap_generic) !ql(i,k) = pq(i,k,i_ice_generic)!+pdq(i,k,i_ice_generic) if(q(i,k).lt.0.)then ! if this is not done, we don't conserve GCS q(i,k)=0. ! vap endif if(ql(i,k).lt.0.)then ql(i,k)=0. ! ice endif enddo enddo ! Initialise the outputs d_t(1:ngrid,1:nlayer) = 0.0 d_q(1:ngrid,1:nlayer) = 0.0 d_ql(1:ngrid,1:nlayer) = 0.0 precip_rate(1:ngrid) = 0.0 precip_rate_tmp(1:ngrid) = 0.0 dq_rain_generic_vap(1:ngrid,1:nlayer,1:nq) = 0.0 dq_rain_generic_cld(1:ngrid,1:nlayer,1:nq) = 0.0 dqsrain_generic(1:ngrid,1:nq) = 0.0 dqssnow_generic(1:ngrid,1:nq) = 0.0 ! calculate saturation mixing ratio do k = 1, nlayer do i = 1, ngrid p_tmp = pplay(i,k) call Psat_generic(zt(i,k),p_tmp,metallicity,psat_tmp,zqs(i,k)) ! calculates psat_tmp & zqs(i,k) ! metallicity --- is not used here but necessary to call function Psat_generic call Lcpdqsat_generic(zt(i,k),p_tmp,psat_tmp,zqs(i,k),dqsat(i,k),dlnpsat(i,k)) ! calculates dqsat(i,k) & dlnpsat(i,k) call Tsat_generic(p_tmp,Tsat(i,k)) ! calculates Tsat(i,k) enddo enddo ! get column / layer conversion factor do k = 1, nlayer do i = 1, ngrid dmass(i,k)=(pplev(i,k)-pplev(i,k+1))/g ! mass per m² in each layer enddo enddo ! Vertical loop (from top to bottom) ! We carry the rain with us and calculate that added by warm/cold precipitation ! processes and that subtracted by evaporation at each level. ! We go from a layer to the next layer below and make the rain fall do k = nlayer, 1, -1 if (evap_prec_generic) then ! note no rneb dependence! do i = 1, ngrid if (precip_rate(i) .GT.0.) then ! if rain from upper layers has fallen in the current layer box if(zt(i,k).gt.Tsat(i,k))then ! if temperature of the layer box is greater than Tsat ! treat the case where all liquid water should boil zqev=MIN((zt(i,k)-Tsat(i,k))*RCPD*dmass(i,k)/RLVTT_generic/ptimestep,precip_rate(i)) ! on évapore precip_rate(i)=MAX(precip_rate(i)-zqev,0.) ! we withdraw from precip_rate the evaporated part d_q(i,k)=zqev/dmass(i,k)*ptimestep ! quantité évaporée d_t(i,k) = - d_q(i,k) * RLVTT_generic/RCPD else ! zqev/zqevt are the maximum amount of vapour that we are allowed to add by evaporation of rain zqev = MAX (0.0, (zqs(i,k)-q(i,k)))*dmass(i,k)/(ptimestep*(1.d0+dqsat(i,k))) !max evaporation to reach saturation implictly accounting for temperature reduction zqevt= MAX (0.0, evap_coeff_generic*2.0e-5*(1.0-q(i,k)/zqs(i,k)) & !default was 2.e-5 *sqrt(precip_rate(i))*dmass(i,k)/pplay(i,k)*zt(i,k)*R) ! BC modif here, is it R or r/(mu/1000) ? zqev = MIN (zqev, zqevt) zqev = MAX (zqev, 0.0) ! not necessary according to previous lines ! we withdraw the evaporated rain precip_rate_tmp(i)= precip_rate(i) - zqev precip_rate_tmp(i)= max(precip_rate_tmp(i),0.0) d_q(i,k) = - (precip_rate_tmp(i)-precip_rate(i))/dmass(i,k)*ptimestep d_t(i,k) = - d_q(i,k) * RLVTT_generic/RCPD precip_rate(i) = precip_rate_tmp(i) end if #ifdef MESOSCALE d_t(i,k) = d_t(i,k)+(pphi(i,k+1)-pphi(i,k))*precip_rate(i)*ptimestep/(RCPD*dmass(i,k)) ! JL22. Accounts for gravitational energy of falling precipitations (probably not to be used in the GCM ! where the counterpart is not included in the dynamics.) ! only for mesoscale simulations !!! #endif endif ! of if (precip_rate(i) .GT.0.) enddo endif ! of if (evap_prec_generic) zoliq(1:ngrid) = 0.0 if(precip_scheme_generic.eq.1)then ! for now this is the only precipitation scheme working here and is therefore the default one (NC2023) do i = 1, ngrid lconvert=rainthreshold_generic if (ql(i,k).gt.1.e-9) then zneb(i) = MAX(rneb(i,k), seuil_neb) ! in mesoscale rneb = 0 or 1 if ((ql(i,k)/zneb(i)).gt.lconvert)then ! precipitate! d_ql(i,k) = -MAX((ql(i,k)-lconvert*zneb(i)),0.0) precip_rate(i) = precip_rate(i) - d_ql(i,k)*dmass(i,k)/ptimestep endif endif enddo elseif (precip_scheme_generic.ge.2) then ! all the precipitation schemes below have not been tested and still have water-dependant routines ! these schemes should be tested, improved ... write(*,*) "Be careful, you have chosen a precip_scheme_generic equal to or greater than 2. For now it has not been tested." do i = 1, ngrid if (rneb(i,k).GT.0.0) then zoliq(i) = ql(i,k) zrho(i) = pplay(i,k) / ( zt(i,k) * R ) zdz(i) = (pplev(i,k)-pplev(i,k+1)) / (zrho(i)*g) zfrac(i) = (zt(i,k)-T_h2O_ice_clouds) / (T_h2O_ice_liq-T_h2O_ice_clouds) zfrac(i) = MAX(zfrac(i), 0.0) zfrac(i) = MIN(zfrac(i), 1.0) zneb(i) = MAX(rneb(i,k), seuil_neb) endif enddo !recalculate liquid GCS particle radii call h2o_cloudrad(ngrid,nlayer,ql,reffh2oliq,reffh2oice) SELECT CASE(precip_scheme_generic) !precip scheme from Sundquist 78 CASE(2) do n = 1, ninter do i = 1, ngrid if (rneb(i,k).GT.0.0) then ! this is the ONLY place where zneb, precip_timescale_generic and cloud_sat_generic are used zchau(i) = (ptimestep/(FLOAT(ninter)*precip_timescale_generic)) * zoliq(i) & * (1.0-EXP(-(zoliq(i)/zneb(i)/cloud_sat_generic)**2)) * zfrac(i) zrhol(i) = zrho(i) * zoliq(i) / zneb(i) zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i) & *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly... ztot(i) = zchau(i) + zfroi(i) if (zneb(i).EQ.seuil_neb) ztot(i) = 0.0 ztot(i) = MIN(MAX(ztot(i),0.0),zoliq(i)) zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0) endif enddo enddo !precip scheme modified from Sundquist 78 (in q**3) CASE(3) do n = 1, ninter do i = 1, ngrid if (rneb(i,k).GT.0.0) then ! this is the ONLY place where zneb, precip_timescale_generic and cloud_sat_generic are used zchau(i) = (ptimestep/(FLOAT(ninter)*precip_timescale_generic*cloud_sat_generic**2)) * (zoliq(i)/zneb(i))**3 zrhol(i) = zrho(i) * zoliq(i) / zneb(i) zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i) & *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly... ztot(i) = zchau(i) + zfroi(i) if (zneb(i).EQ.seuil_neb) ztot(i) = 0.0 ztot(i) = MIN(MAX(ztot(i),0.0),zoliq(i)) zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0) endif enddo enddo !precip scheme modified from Boucher 95 CASE(4) do n = 1, ninter do i = 1, ngrid if (rneb(i,k).GT.0.0) then ! this is the ONLY place where zneb and c1 are used zchau(i) = ptimestep/FLOAT(ninter) *c1* zrho(i) & *(zoliq(i)/zneb(i))**2*reffh2oliq(i,k)*zneb(i)* zfrac(i) zrhol(i) = zrho(i) * zoliq(i) / zneb(i) zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i) & *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly... ztot(i) = zchau(i) + zfroi(i) if (zneb(i).EQ.seuil_neb) ztot(i) = 0.0 ztot(i) = MIN(MAX(ztot(i),0.0),zoliq(i)) zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0) endif enddo enddo END SELECT ! precip_scheme_generic ! Change in cloud density and surface GCS values do i = 1, ngrid if (rneb(i,k).GT.0.0) then d_ql(i,k) = (zoliq(i) - ql(i,k))!/ptimestep precip_rate(i) = precip_rate(i)+ MAX(ql(i,k)-zoliq(i),0.0)*dmass(i,k)/ptimestep endif enddo endif ! if precip_scheme_generic enddo ! of do k = nlayer, 1, -1 ! Rain or snow on the ground do i = 1, ngrid if(precip_rate(i).lt.0.0)then print*,'Droplets of negative rain are falling...' call abort endif ! lines below should be uncommented and customized if you want to make a difference between rain and snow !if (t(i,1) .LT. T_h2O_ice_liq) then ! dqssnow_generic(i,i_ice_generic) = precip_rate(i) ! dqsrain_generic(i,i_ice_generic) = 0.0 !else ! dqssnow_generic(i,i_ice_generic) = 0.0 ! dqsrain_generic(i,i_ice_generic) = precip_rate(i) ! liquid water = ice for now !endif ! For now we make no distinction between rain and snow, and consider only rain: dqsrain_generic(i,i_ice_generic) = precip_rate(i) dqssnow_generic(i,i_ice_generic) = 0.0 enddo ! now subroutine -----> GCM variables if (evap_prec_generic) then dq_rain_generic_vap(1:ngrid,1:nlayer,i_vap_generic)=d_q(1:ngrid,1:nlayer)/ptimestep d_t(1:ngrid,1:nlayer)=d_t(1:ngrid,1:nlayer)/ptimestep do i=1,ngrid reevap_precip(i,i_vap_generic)=0. do k=1,nlayer reevap_precip(i,i_vap_generic)=reevap_precip(i,i_vap_generic)+dq_rain_generic_vap(i,k,i_vap_generic)*dmass(i,k) enddo enddo else dq_rain_generic_vap(1:ngrid,1:nlayer,i_vap_generic)=0.0 d_t(1:ngrid,1:nlayer)=0.0 endif dq_rain_generic_cld(1:ngrid,1:nlayer,i_ice_generic) = d_ql(1:ngrid,1:nlayer)/ptimestep end if ! if(call_ice_vap_generic) end do ! do iq=1,nq loop on tracers end subroutine rain_generic