Newer
Older
!MNH_LIC Copyright 2004-2019 CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC This is part of the Meso-NH software governed by the CeCILL-C licence
!MNH_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
!MNH_LIC for details. version 1.
!-----------------------------------------------------------------
! ######spl
MODULE MODE_COMPUTE_UPDRAFT
! ###########################
!
IMPLICIT NONE
CONTAINS
SUBROUTINE COMPUTE_UPDRAFT(D,CST,NEBN,PARAMMF,TURBN,CSTURB, &

RIETTE Sébastien
committed
KSV, &

RODIER Quentin
committed
OENTR_DETR, &
ONOMIXLG,KSV_LGBEG,KSV_LGEND, &
PZZ,PDZZ, &
PSFTH,PSFRV, &
PPABSM,PRHODREF,PUM,PVM, PTKEM, &
PTHM,PRVM,PTHLM,PRTM, &
PSVM,PTHL_UP,PRT_UP, &
PRV_UP,PRC_UP,PRI_UP,PTHV_UP, &
PW_UP,PU_UP, PV_UP, PSV_UP, &
PFRAC_UP,PFRAC_ICE_UP,PRSAT_UP, &
PEMF,PDETR,PENTR, &
PBUO_INTEG,KKLCL,KKETL,KKCTL, &
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
! #################################################################
!!
!!**** *COMPUTE_UPDRAFT* - calculates caracteristics of the updraft
!!
!!
!! PURPOSE
!! -------
!!**** The purpose of this routine is to build the updraft model
!!
!
!!** METHOD
!! ------
!!
!! EXTERNAL
!! --------
!!
!! IMPLICIT ARGUMENTS
!! ------------------
!!
!! !! REFERENCE
!! ---------
!! Book 1 of Meso-NH documentation (chapter Turbulence)
!! Soares et al. 2004 QJ
!!
!! AUTHOR
!! ------
!! J.Pergaud
!! V.Masson : Optimization 07/2010
!! S. Riette : 07/2010 : modification for reproducibility
!! S. Riette may 2011: ice added, interface modified
!! S. Riette Jan 2012: support for both order of vertical levels
!! V.Masson, C.Lac : 02/2011 : SV_UP initialized by a non-zero value
!! S. Riette Apr 2013: improvement of continuity at the condensation level
!! R.Honnert Oct 2016 : Add ZSURF and Update with AROME
!! Q.Rodier 01/2019 : support RM17 mixing length
!! R.Honnert 01/2019 : add LGZ (reduction of the mass-flux surface closure with the resolution)
!! S. Riette 06/2022: compute_entr_detr is inlined
!! --------------------------------------------------------------------------
!
!* 0. DECLARATIONS
! ------------
!
USE MODD_DIMPHYEX, ONLY: DIMPHYEX_t
USE MODD_CST, ONLY: CST_t
USE MODD_PARAM_MFSHALL_n, ONLY: PARAM_MFSHALL_t
USE MODD_TURB_n, ONLY: TURB_t

RIETTE Sébastien
committed
USE MODD_CTURB, ONLY: CSTURB_t
USE MODI_SHUMAN_MF, ONLY: MZM_MF, MZF_MF, GZ_M_W_MF
USE MODE_COMPUTE_BL89_ML, ONLY: COMPUTE_BL89_ML
USE MODE_MSG, ONLY: PRINT_MSG, NVERB_FATAL
USE YOMHOOK , ONLY : LHOOK, DR_HOOK, JPHOOK
IMPLICIT NONE
!* 1.1 Declaration of Arguments
!
!
!
TYPE(DIMPHYEX_t), INTENT(IN) :: D
TYPE(CST_t), INTENT(IN) :: CST
TYPE(PARAM_MFSHALL_t), INTENT(IN) :: PARAMMF

RODIER Quentin
committed
TYPE(TURB_t), INTENT(IN) :: TURBN

RIETTE Sébastien
committed
TYPE(CSTURB_t), INTENT(IN) :: CSTURB
INTEGER, INTENT(IN) :: KSV
LOGICAL, INTENT(IN) :: OENTR_DETR! flag to recompute entrainment, detrainment and mass flux
LOGICAL, INTENT(IN) :: ONOMIXLG ! False if mixing of lagrangian tracer
INTEGER, INTENT(IN) :: KSV_LGBEG ! first index of lag. tracer
INTEGER, INTENT(IN) :: KSV_LGEND ! last index of lag. tracer
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PZZ ! Height at the flux point
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PDZZ ! Metrics coefficient
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PSFTH,PSFRV
! normal surface fluxes of theta,rv,(u,v) parallel to the orography
!
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PPABSM ! Pressure at t-dt
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PRHODREF ! dry density of the
! reference state
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PUM ! u mean wind
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PVM ! v mean wind
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PTKEM ! TKE at t-dt
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PTHM ! liquid pot. temp. at t-dt
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PRVM ! vapor mixing ratio at t-dt
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PTHLM,PRTM ! cons. var. at t-dt
REAL, DIMENSION(D%NIJT,D%NKT,KSV), INTENT(IN) :: PSVM ! scalar var. at t-dt
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT) :: PTHL_UP,PRT_UP ! updraft properties
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT) :: PU_UP, PV_UP ! updraft wind components
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(INOUT):: PRV_UP,PRC_UP, & ! updraft rv, rc
PRI_UP,PTHV_UP,& ! updraft ri, THv
PW_UP,PFRAC_UP,& ! updraft w, fraction
PFRAC_ICE_UP,& ! liquid/solid fraction in updraft
PRSAT_UP ! Rsat
REAL, DIMENSION(D%NIJT,D%NKT,KSV), INTENT(OUT) :: PSV_UP ! updraft scalar var.
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(INOUT):: PEMF,PDETR,PENTR ! Mass_flux,
! detrainment,entrainment
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(INOUT) :: PBUO_INTEG ! Integrated Buoyancy
INTEGER, DIMENSION(D%NIJT), INTENT(INOUT) :: KKLCL,KKETL,KKCTL! LCL, ETL, CTL
REAL, DIMENSION(D%NIJT), INTENT(OUT) :: PDEPTH ! Deepness of cloud
! 1.2 Declaration of local variables
!
!
! Mean environment variables at t-dt at flux point
REAL, DIMENSION(D%NIJT,D%NKT) :: &
ZTHM_F,ZRVM_F ! Theta,rv of
! updraft environnement
REAL, DIMENSION(D%NIJT,D%NKT) :: &
ZRTM_F, ZTHLM_F, ZTKEM_F,& ! rt, thetal,TKE,pressure,
ZUM_F,ZVM_F,ZRHO_F, & ! density,momentum
ZPRES_F,ZTHVM_F,ZTHVM, & ! interpolated at the flux point
ZG_O_THVREF, & ! g*ThetaV ref
ZW_UP2, & ! w**2 of the updraft
ZBUO_INTEG_DRY, ZBUO_INTEG_CLD,&! Integrated Buoyancy
ZENTR_CLD,ZDETR_CLD ! wet entrainment and detrainment
REAL, DIMENSION(D%NIJT,D%NKT,KSV) :: &
ZSVM_F ! scalar variables
REAL, DIMENSION(D%NIJT,D%NKT) :: &
ZTH_UP, & ! updraft THETA
ZRC_MIX, ZRI_MIX ! guess of Rc and Ri for KF mixture
REAL, DIMENSION(D%NIJT,D%NKT) :: ZCOEF ! diminution coefficient for too high clouds
REAL, DIMENSION(D%NIJT) :: ZWTHVSURF ! Surface w'thetav'
REAL :: ZRDORV ! RD/RV
REAL :: ZRVORD ! RV/RD
REAL, DIMENSION(D%NIJT) :: ZMIX1,ZMIX2,ZMIX3_CLD,ZMIX2_CLD
REAL, DIMENSION(D%NIJT) :: ZLUP ! Upward Mixing length from the ground
INTEGER :: JK,JIJ,JSV ! loop counters
LOGICAL, DIMENSION(D%NIJT) :: GTEST,GTESTLCL,GTESTETL
! Test if the ascent continue, if LCL or ETL is reached
LOGICAL :: GLMIX
! To choose upward or downward mixing length
LOGICAL, DIMENSION(D%NIJT) :: GWORK1
LOGICAL, DIMENSION(D%NIJT,D%NKT) :: GWORK2
INTEGER :: ITEST
REAL, DIMENSION(D%NIJT) :: ZRC_UP, ZRI_UP, ZRV_UP,&
ZRSATW, ZRSATI,&
ZPART_DRY
REAL :: ZDEPTH_MAX1, ZDEPTH_MAX2 ! control auto-extinction process
REAL :: ZTMAX,ZRMAX ! control value
REAL, DIMENSION(D%NIJT) :: ZSURF
REAL, DIMENSION(D%NIJT,D%NKT) :: ZSHEAR,ZDUDZ,ZDVDZ ! vertical wind shear

RODIER Quentin
committed
REAL, DIMENSION(D%NIJT,D%NKT) :: ZWK, KDEPTH
REAL, DIMENSION(D%NIJT,16) :: ZBUF
REAL(KIND=JPHOOK) :: ZHOOK_HANDLE
!
! 1.3 Declaration of additional local variables for compute_entr_detr
!
! Variables for cloudy part
REAL, DIMENSION(D%NIJT) :: ZKIC, ZKIC_F2 ! fraction of env. mass in the muxtures
REAL, DIMENSION(D%NIJT) :: ZEPSI,ZDELTA ! factor entrainment detrainment
REAL :: ZEPSI_CLOUD ! factor entrainment detrainment
REAL :: ZCOEFFMF_CLOUD ! factor for compputing entr. detr.
REAL, DIMENSION(D%NIJT) :: ZMIXTHL,ZMIXRT ! Thetal and rt in the mixtures
REAL, DIMENSION(D%NIJT) :: ZTHMIX ! Theta and Thetav of mixtures
REAL, DIMENSION(D%NIJT) :: ZRVMIX,ZRCMIX,ZRIMIX ! mixing ratios in mixtures
REAL, DIMENSION(D%NIJT) :: ZTHVMIX, ZTHVMIX_F2 ! Theta and Thetav of mixtures
REAL, DIMENSION(D%NIJT) :: ZTHV_UP_F2 ! thv_up at flux point kk+kkl
REAL, DIMENSION(D%NIJT) :: ZRSATW_ED, ZRSATI_ED ! working arrays (mixing ratio at saturation)
REAL, DIMENSION(D%NIJT) :: ZTHV ! theta V of environment at the bottom of cloudy part
REAL :: ZKIC_INIT !Initial value of ZKIC
REAL :: ZCOTHVU ! Variation of Thvup between bottom and top of cloudy part
! Variables for dry part
REAL :: ZFOESW, ZFOESI ! saturating vapor pressure
REAL :: ZDRSATODP ! d.Rsat/dP
REAL :: ZT ! Temperature
REAL :: ZWK0D ! Work array
! Variables for dry and cloudy parts
REAL, DIMENSION(D%NIJT) :: ZCOEFF_MINUS_HALF,& ! Variation of Thv between mass points kk-kkl and kk
ZCOEFF_PLUS_HALF ! Variation of Thv between mass points kk and kk+kkl
REAL, DIMENSION(D%NIJT) :: ZPRE ! pressure at the bottom of the cloudy part
REAL, DIMENSION(D%NIJT) :: ZG_O_THVREF_ED
REAL, DIMENSION(D%NIJT) :: ZFRAC_ICE ! fraction of ice
REAL, DIMENSION(D%NIJT) :: ZDZ_STOP,& ! Exact Height of the LCL above flux level KK
ZTHV_MINUS_HALF,& ! Thv at flux point(kk)
ZTHV_PLUS_HALF ! Thv at flux point(kk+kkl)
REAL :: ZDZ ! Delta Z used in computations
INTEGER :: JKLIM
INTEGER :: IIJB,IIJE ! physical horizontal domain indices
INTEGER :: IKT,IKB,IKE,IKL
IF (LHOOK) CALL DR_HOOK('COMPUTE_UPDRAFT',0,ZHOOK_HANDLE)
!
IIJE=D%NIJE
IIJB=D%NIJB
IKT=D%NKT
IKB=D%NKB
IKE=D%NKE
IKL=D%NKL
!
! Thresholds for the perturbation of
! theta_l and r_t at the first level of the updraft
ZTMAX=2.0
ZRMAX=1.E-3
!------------------------------------------------------------------------
! INITIALISATION
! Initialisation of the constants
ZRDORV = CST%XRD / CST%XRV !=0.622
ZRVORD = (CST%XRV / CST%XRD)
ZDEPTH_MAX1=3000. ! clouds with depth inferior to this value are keeped untouched
ZDEPTH_MAX2=4000. ! clouds with depth superior to this value are suppressed
! Local variables, internal domain
IF (OENTR_DETR) THEN
! Initialisation of intersesting Level :LCL,ETL,CTL
KKLCL(:)=IKE
KKETL(:)=IKE
KKCTL(:)=IKE
!
! Initialisation
!* udraft governing variables
PEMF(:,:)=0.
PDETR(:,:)=0.
PENTR(:,:)=0.
! Initialisation
!* updraft core variables
PRV_UP(:,:)=0.
PRC_UP(:,:)=0.
PRI_UP(:,:)=0.
PW_UP(:,:)=0.
ZTH_UP(:,:)=0.
PFRAC_UP(:,:)=0.
PTHV_UP(:,:)=0.
PBUO_INTEG=0.
PFRAC_ICE_UP(:,:)=0.
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
PRSAT_UP(:,:)=PRVM(:,:) ! should be initialised correctly but is (normaly) not used
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!cloud/dry air mixture cloud content
ZRC_MIX = 0.
ZRI_MIX = 0.
END IF
! Initialisation of environment variables at t-dt
! variables at flux level
CALL MZM_MF(D, PTHLM(:,:), ZTHLM_F(:,:))
CALL MZM_MF(D, PRTM(:,:), ZRTM_F (:,:))
CALL MZM_MF(D, PUM(:,:), ZUM_F (:,:))
CALL MZM_MF(D, PVM(:,:), ZVM_F (:,:))
CALL MZM_MF(D, PTKEM(:,:), ZTKEM_F(:,:))

RIETTE Sébastien
committed
DO JSV=1,KSV
IF (ONOMIXLG .AND. JSV >= KSV_LGBEG .AND. JSV<= KSV_LGEND) CYCLE
CALL MZM_MF(D, PSVM(:,:,JSV), ZSVM_F(:,:,JSV))
END DO
!
! Initialisation of updraft characteristics
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
PTHL_UP(:,:)=ZTHLM_F(:,:)
PRT_UP(:,:)=ZRTM_F(:,:)
PU_UP(:,:)=ZUM_F(:,:)
PV_UP(:,:)=ZVM_F(:,:)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT,JSV=1:KSV)

RODIER Quentin
committed
PSV_UP(:,:,:)=ZSVM_F(:,:,:)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT,JSV=1:KSV)
! Computation or initialisation of updraft characteristics at the KKB level
! thetal_up,rt_up,thetaV_up, w2,Buoyancy term and mass flux (PEMF)
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
PTHL_UP(:,IKB)= ZTHLM_F(:,IKB)+ &
& MAX(0.,MIN(ZTMAX,(PSFTH(:)/SQRT(ZTKEM_F(:,IKB)))* PARAMMF%XALP_PERT))
PRT_UP(:,IKB) = ZRTM_F(:,IKB)+ &
& MAX(0.,MIN(ZRMAX,(PSFRV(:)/SQRT(ZTKEM_F(:,IKB)))* PARAMMF%XALP_PERT))
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
IF (OENTR_DETR) THEN
CALL MZM_MF(D, PTHM (:,:), ZTHM_F (:,:))
CALL MZM_MF(D, PPABSM(:,:), ZPRES_F(:,:))
CALL MZM_MF(D, PRHODREF(:,:), ZRHO_F (:,:))
CALL MZM_MF(D, PRVM(:,:), ZRVM_F (:,:))
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
! thetav at mass and flux levels

RODIER Quentin
committed
ZTHVM_F(:,:)=ZTHM_F(:,:)* &
&((1.+ZRVORD*ZRVM_F(:,:))/(1.+ZRTM_F(:,:)))
ZTHVM(:,:)=PTHM(:,:)* &
&((1.+ZRVORD*PRVM(:,:))/(1.+PRTM(:,:)))

RODIER Quentin
committed
PTHV_UP(:,:)=ZTHVM_F(:,:)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
ZW_UP2(:,:)=0.
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZW_UP2(:,IKB) = MAX(0.0001,(2./3.)*ZTKEM_F(:,IKB))
! Computation of non conservative variable for the KKB level of the updraft
! (all or nothing ajustement)
PRC_UP(:,IKB)=0.
PRI_UP(:,IKB)=0.
!$mnh_end_expand_array(JIJ=IIJB:IIJE)

RIETTE Sébastien
committed
CALL TH_R_FROM_THL_RT(CST, NEBN, D%NIJT, NEBN%CFRAC_ICE_SHALLOW_MF,PFRAC_ICE_UP(:,IKB),ZPRES_F(:,IKB), &
PTHL_UP(:,IKB),PRT_UP(:,IKB),ZTH_UP(:,IKB), &
PRV_UP(:,IKB),PRC_UP(:,IKB),PRI_UP(:,IKB),ZRSATW(:),ZRSATI(:), OOCEAN=.FALSE., &
PBUF=ZBUF(:,:), KB=D%NIJB, KE=D%NIJE)
!$mnh_expand_array(JIJ=IIJB:IIJE)
! compute updraft thevav and buoyancy term at KKB level

RODIER Quentin
committed
PTHV_UP(:,IKB) = ZTH_UP(:,IKB)*&
& ((1+ZRVORD*PRV_UP(:,IKB))/(1+PRT_UP(:,IKB)))
! compute mean rsat in updraft

RODIER Quentin
committed
PRSAT_UP(:,IKB) = ZRSATW(:)*(1-PFRAC_ICE_UP(:,IKB)) + &
& ZRSATI(:)*PFRAC_ICE_UP(:,IKB)
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
! Closure assumption for mass flux at KKB level
!
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZG_O_THVREF(:,:)=CST%XG/ZTHVM_F(:,:)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
! compute L_up
GLMIX=.TRUE.
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZTKEM_F(:,IKB)=0.
!$mnh_end_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
IF(TURBN%CTURBLEN=='RM17') THEN
CALL GZ_M_W_MF(D, PUM, PDZZ, ZWK)
CALL MZF_MF(D, ZWK, ZDUDZ)
CALL GZ_M_W_MF(D, PVM, PDZZ, ZWK)
CALL MZF_MF(D, ZWK, ZDVDZ)
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZSHEAR(:,:) = SQRT(ZDUDZ(:,:)**2 + ZDVDZ(:,:)**2)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
ELSE
ZSHEAR = 0. !no shear in bl89 mixing length
END IF
CALL COMPUTE_BL89_ML(D, CST, CSTURB, PDZZ,ZTKEM_F(:,IKB),&
&ZG_O_THVREF(:,IKB),ZTHVM,IKB,GLMIX,.TRUE.,ZSHEAR,ZLUP)

RIETTE Sébastien
committed
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZLUP(:)=MAX(ZLUP(:),1.E-10)
! Compute Buoyancy flux at the ground

RODIER Quentin
committed
ZWTHVSURF(:) = (ZTHVM_F(:,IKB)/ZTHM_F(:,IKB))*PSFTH(:)+ &
(0.61*ZTHM_F(:,IKB))*PSFRV(:)

RIETTE Sébastien
committed
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
! Mass flux at KKB level (updraft triggered if PSFTH>0.)
IF (PARAMMF%LGZ) THEN
IF(PDX==0. .OR. PDY==0.) THEN
CALL PRINT_MSG(NVERB_FATAL, 'GEN', 'COMPUTE_UPDRAFT', 'PDX or PDY is NULL with option LGZ!')

RIETTE Sébastien
committed
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZSURF(:)=TANH(PARAMMF%XGZ*SQRT(PDX*PDY)/ZLUP(:))

RIETTE Sébastien
committed
!$mnh_end_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZSURF(:)=1.

RIETTE Sébastien
committed
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE (ZWTHVSURF(:)>0.)
PEMF(:,IKB) = PARAMMF%XCMF * ZSURF(:) * ZRHO_F(:,IKB) * &
((ZG_O_THVREF(:,IKB))*ZWTHVSURF(:)*ZLUP(:))**(1./3.)
PFRAC_UP(:,IKB)=MIN(PEMF(:,IKB)/(SQRT(ZW_UP2(:,IKB))*ZRHO_F(:,IKB)), &
&PARAMMF%XFRAC_UP_MAX)

RODIER Quentin
committed
ZW_UP2(:,IKB)=(PEMF(:,IKB)/(PFRAC_UP(:,IKB)*ZRHO_F(:,IKB)))**2
GTEST(:)=.TRUE.
ELSEWHERE

RODIER Quentin
committed
PEMF(:,IKB) =0.
GTEST(:)=.FALSE.
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
GTEST(:)=PEMF(:,IKB+IKL)>0.
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
END IF
!--------------------------------------------------------------------------
! 3. Vertical ascending loop
! -----------------------
!
! If GTEST = T the updraft starts from the KKB level and stops when GTEST becomes F
!
!
GTESTLCL(:)=.FALSE.
GTESTETL(:)=.FALSE.
! Loop on vertical level
DO JK=IKB,IKE-IKL,IKL
! IF the updraft top is reached for all column, stop the loop on levels
ITEST=0
DO JIJ=IIJB,IIJE
IF(GTEST(JIJ)) ITEST = ITEST + 1
END DO
IF (ITEST==0) CYCLE
! Computation of entrainment and detrainment with KF90
! parameterization in clouds and LR01 in subcloud layer
! to find the LCL (check if JK is LCL or not)
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE ((PRC_UP(:,JK)+PRI_UP(:,JK)>0.).AND.(.NOT.(GTESTLCL(:))))
KKLCL(:) = JK
GTESTLCL(:)=.TRUE.
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
! COMPUTE PENTR and PDETR at mass level JK
IF (OENTR_DETR) THEN
IF(JK/=IKB) THEN
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZRC_MIX(:,JK) = ZRC_MIX(:,JK-IKL) ! guess of Rc of mixture
ZRI_MIX(:,JK) = ZRI_MIX(:,JK-IKL) ! guess of Ri of mixture
!$mnh_end_expand_array(JIJ=IIJB:IIJE)

RIETTE Sébastien
committed
CALL COMPUTE_ENTR_DETR(D, CST, NEBN, PARAMMF, JK,IKB,IKE,IKL,GTEST,GTESTLCL,PFRAC_ICE_UP(:,JK),&
PRHODREF(:,JK),ZPRES_F(:,JK),ZPRES_F(:,JK+IKL),&
PZZ(:,:),PDZZ(:,:),ZTHVM(:,:), &
PTHLM(:,:),PRTM(:,:),ZW_UP2(:,:),ZTH_UP(:,JK), &
PTHL_UP(:,JK),PRT_UP(:,JK),ZLUP(:), &
PRC_UP(:,JK),PRI_UP(:,JK),PTHV_UP(:,JK),&
PRSAT_UP(:,JK),ZRC_MIX(:,JK),ZRI_MIX(:,JK), &
PENTR(:,JK),PDETR(:,JK),ZENTR_CLD(:,JK),ZDETR_CLD(:,JK),&
ZBUO_INTEG_DRY(:,JK), ZBUO_INTEG_CLD(:,JK), &
ZPART_DRY(:) )
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
PBUO_INTEG(:,JK)=ZBUO_INTEG_DRY(:,JK)+ZBUO_INTEG_CLD(:,JK)
IF (JK==IKB) THEN

RODIER Quentin
committed
PDETR(:,JK)=0.
ZDETR_CLD(:,JK)=0.
ENDIF
! Computation of updraft characteristics at level JK+KKL

RODIER Quentin
committed
WHERE(GTEST(:))
ZMIX1(:)=0.5*(PZZ(:,JK+IKL)-PZZ(:,JK))*&
&(PENTR(:,JK)-PDETR(:,JK))
PEMF(:,JK+IKL)=PEMF(:,JK)*EXP(2*ZMIX1(:))
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
ELSE !OENTR_DETR
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
GTEST(:) = (PEMF(:,JK+IKL)>0.)
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
END IF !OENTR_DETR
! stop the updraft if MF becomes negative
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE (GTEST(:).AND.(PEMF(:,JK+IKL)<=0.))
PEMF(:,JK+IKL)=0.
KKCTL(:) = JK+IKL
GTEST(:)=.FALSE.
PFRAC_ICE_UP(:,JK+IKL)=PFRAC_ICE_UP(:,JK)
PRSAT_UP(:,JK+IKL)=PRSAT_UP(:,JK)
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
! If the updraft did not stop, compute cons updraft characteritics at jk+KKL

RODIER Quentin
committed
DO JIJ=IIJB,IIJE
IF(GTEST(JIJ)) THEN
ZMIX2(JIJ) = (PZZ(JIJ,JK+IKL)-PZZ(JIJ,JK))*PENTR(JIJ,JK) !&
ZMIX3_CLD(JIJ) = (PZZ(JIJ,JK+IKL)-PZZ(JIJ,JK))*(1.-ZPART_DRY(JIJ))*ZDETR_CLD(JIJ,JK) !&
ZMIX2_CLD(JIJ) = (PZZ(JIJ,JK+IKL)-PZZ(JIJ,JK))*(1.-ZPART_DRY(JIJ))*ZENTR_CLD(JIJ,JK)
PTHL_UP(JIJ,JK+IKL)=(PTHL_UP(JIJ,JK)*(1.-0.5*ZMIX2(JIJ)) + PTHLM(JIJ,JK)*ZMIX2(JIJ)) &
/(1.+0.5*ZMIX2(JIJ))
PRT_UP(JIJ,JK+IKL) =(PRT_UP (JIJ,JK)*(1.-0.5*ZMIX2(JIJ)) + PRTM(JIJ,JK)*ZMIX2(JIJ)) &
/(1.+0.5*ZMIX2(JIJ))

RODIER Quentin
committed
IF(PARAMMF%LMIXUV) THEN
IF(JK/=IKB) THEN
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE(GTEST(:))
PU_UP(:,JK+IKL) = (PU_UP(:,JK)*(1-0.5*ZMIX2(:)) + &
&PUM(:,JK)*ZMIX2(:)+ &
0.5*PARAMMF%XPRES_UV*(PZZ(:,JK+IKL)-PZZ(:,JK))*&
((PUM(:,JK+IKL)-PUM(:,JK))/PDZZ(:,JK+IKL)+&
(PUM(:,JK)-PUM(:,JK-IKL))/PDZZ(:,JK)) ) &
/(1+0.5*ZMIX2(:))
PV_UP(:,JK+IKL) = (PV_UP(:,JK)*(1-0.5*ZMIX2(:)) + &
&PVM(:,JK)*ZMIX2(:)+ &
0.5*PARAMMF%XPRES_UV*(PZZ(:,JK+IKL)-PZZ(:,JK))*&
((PVM(:,JK+IKL)-PVM(:,JK))/PDZZ(:,JK+IKL)+&
(PVM(:,JK)-PVM(:,JK-IKL))/PDZZ(:,JK)) ) &
/(1+0.5*ZMIX2(:))
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE(GTEST(:))
PU_UP(:,JK+IKL) = (PU_UP(:,JK)*(1-0.5*ZMIX2(:)) + &
&PUM(:,JK)*ZMIX2(:)+ &
0.5*PARAMMF%XPRES_UV*(PZZ(:,JK+IKL)-PZZ(:,JK))*&
((PUM(:,JK+IKL)-PUM(:,JK))/PDZZ(:,JK+IKL)) ) &
/(1+0.5*ZMIX2(:))
PV_UP(:,JK+IKL) = (PV_UP(:,JK)*(1-0.5*ZMIX2(:)) + &
&PVM(:,JK)*ZMIX2(:)+ &
0.5*PARAMMF%XPRES_UV*(PZZ(:,JK+IKL)-PZZ(:,JK))*&
((PVM(:,JK+IKL)-PVM(:,JK))/PDZZ(:,JK+IKL)) ) &
/(1+0.5*ZMIX2(:))
!$mnh_end_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ENDIF !PARAMMF%LMIXUV

RIETTE Sébastien
committed
DO JSV=1,KSV
IF (ONOMIXLG .AND. JSV >= KSV_LGBEG .AND. JSV<= KSV_LGEND) CYCLE
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE(GTEST(:))
PSV_UP(:,JK+IKL,JSV) = (PSV_UP(:,JK,JSV)*(1-0.5*ZMIX2(:)) + &
PSVM(:,JK,JSV)*ZMIX2(:)) /(1+0.5*ZMIX2(:))
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
END DO
IF (OENTR_DETR) THEN
! Compute non cons. var. at level JK+KKL
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZRC_UP(:)=PRC_UP(:,JK) ! guess = level just below
ZRI_UP(:)=PRI_UP(:,JK) ! guess = level just below
!$mnh_end_expand_array(JIJ=IIJB:IIJE)

RIETTE Sébastien
committed
CALL TH_R_FROM_THL_RT(CST, NEBN, D%NIJT, NEBN%CFRAC_ICE_SHALLOW_MF,PFRAC_ICE_UP(:,JK+IKL),ZPRES_F(:,JK+IKL), &
PTHL_UP(:,JK+IKL),PRT_UP(:,JK+IKL),ZTH_UP(:,JK+IKL), &
ZRV_UP(:),ZRC_UP(:),ZRI_UP(:),ZRSATW(:),ZRSATI(:), OOCEAN=.FALSE., &
PBUF=ZBUF(:,:), KB=D%NIJB, KE=D%NIJE)
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE(GTEST(:))
PRC_UP(:,JK+IKL)=ZRC_UP(:)
PRV_UP(:,JK+IKL)=ZRV_UP(:)
PRI_UP(:,JK+IKL)=ZRI_UP(:)
PRSAT_UP(:,JK+IKL) = ZRSATW(:)*(1-PFRAC_ICE_UP(:,JK+IKL)) + &
& ZRSATI(:)*PFRAC_ICE_UP(:,JK+IKL)
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
! Compute the updraft theta_v, buoyancy and w**2 for level JK+KKL
DO JIJ=IIJB,IIJE
IF(GTEST(JIJ)) THEN
PTHV_UP(JIJ,JK+IKL) = ZTH_UP(JIJ,JK+IKL)* &
& ((1+ZRVORD*PRV_UP(JIJ,JK+IKL))/(1+PRT_UP(JIJ,JK+IKL)))
IF (ZBUO_INTEG_DRY(JIJ,JK)>0.) THEN
ZW_UP2(JIJ,JK+IKL) = ZW_UP2(JIJ,JK) + 2.*(PARAMMF%XABUO-PARAMMF%XBENTR*PARAMMF%XENTR_DRY)* &
&ZBUO_INTEG_DRY(JIJ,JK)
ELSE
ZW_UP2(JIJ,JK+IKL) = ZW_UP2(JIJ,JK) + 2.*PARAMMF%XABUO* ZBUO_INTEG_DRY(JIJ,JK)
END IF
ZW_UP2(JIJ,JK+IKL) = ZW_UP2(JIJ,JK+IKL)*(1.-(PARAMMF%XBDETR*ZMIX3_CLD(JIJ)+ &
&PARAMMF%XBENTR*ZMIX2_CLD(JIJ)))&
/(1.+(PARAMMF%XBDETR*ZMIX3_CLD(JIJ)+PARAMMF%XBENTR*ZMIX2_CLD(JIJ))) &
+2.*(PARAMMF%XABUO)*ZBUO_INTEG_CLD(JIJ,JK)/ &
&(1.+(PARAMMF%XBDETR*ZMIX3_CLD(JIJ)+PARAMMF%XBENTR*ZMIX2_CLD(JIJ)))
END IF
END DO
! Test if the updraft has reach the ETL
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE (GTEST(:).AND.(PBUO_INTEG(:,JK)<=0.))
KKETL(:) = JK+IKL
GTESTETL(:)=.TRUE.

RODIER Quentin
committed
GTESTETL(:)=.FALSE.
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
! Test is we have reached the top of the updraft
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE (GTEST(:).AND.((ZW_UP2(:,JK+IKL)<=0.).OR.(PEMF(:,JK+IKL)<=0.)))
ZW_UP2(:,JK+IKL)=0.
PEMF(:,JK+IKL)=0.
GTEST(:)=.FALSE.
PTHL_UP(:,JK+IKL)=ZTHLM_F(:,JK+IKL)
PRT_UP(:,JK+IKL)=ZRTM_F(:,JK+IKL)
PRC_UP(:,JK+IKL)=0.
PRI_UP(:,JK+IKL)=0.
PRV_UP(:,JK+IKL)=0.
PTHV_UP(:,JK+IKL)=ZTHVM_F(:,JK+IKL)
PFRAC_UP(:,JK+IKL)=0.
KKCTL(:)=JK+IKL
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
! compute frac_up at JK+KKL
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE (GTEST(:))
PFRAC_UP(:,JK+IKL)=PEMF(:,JK+IKL)/&
&(SQRT(ZW_UP2(:,JK+IKL))*ZRHO_F(:,JK+IKL))
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
! Updraft fraction must be smaller than XFRAC_UP_MAX
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE (GTEST(:))
PFRAC_UP(:,JK+IKL)=MIN(PARAMMF%XFRAC_UP_MAX,PFRAC_UP(:,JK+IKL))
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
! When cloudy and non-buoyant, updraft fraction must decrease
!$mnh_expand_where(JIJ=IIJB:IIJE)

RODIER Quentin
committed
WHERE ((GTEST(:).AND.GTESTETL(:)).AND.GTESTLCL(:))
PFRAC_UP(:,JK+IKL)=MIN(PFRAC_UP(:,JK+IKL),PFRAC_UP(:,JK))
!$mnh_end_expand_where(JIJ=IIJB:IIJE)
! Mass flux is updated with the new updraft fraction
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
IF (OENTR_DETR) PEMF(:,JK+IKL)=PFRAC_UP(:,JK+IKL)*SQRT(ZW_UP2(:,JK+IKL))* &
&ZRHO_F(:,JK+IKL)
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
END IF !OENTR_DETR
ENDDO
IF(OENTR_DETR) THEN
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
PW_UP(:,:)=SQRT(ZW_UP2(:,:))
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
PEMF(:,IKB) =0.
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
! Limits the shallow convection scheme when cloud heigth is higher than 3000m.
! To do this, mass flux is multiplied by a coefficient decreasing linearly
! from 1 (for clouds of ZDEPTH_MAX1 m of depth) to 0 (for clouds of ZDEPTH_MAX2 m of depth).
! This way, all MF fluxes are diminished by this amount.
! Diagnosed cloud fraction is also multiplied by the same coefficient.
!

RODIER Quentin
committed
DO JIJ=IIJB,IIJE
PDEPTH(JIJ) = MAX(0., PZZ(JIJ,KKCTL(JIJ)) - PZZ(JIJ,KKLCL(JIJ)) )

RODIER Quentin
committed
IF(PARAMMF%LVERLIMUP) THEN
DO JK=1,IKT
DO JIJ=IIJB,IIJE
KDEPTH(JIJ,JK) = MIN(MAX(0., PZZ(JIJ,JK) - PZZ(JIJ,KKLCL(JIJ)) ), PDEPTH(JIJ))
END DO
END DO
END IF
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
GWORK1(:)= (GTESTLCL(:) .AND. (PDEPTH(:) > ZDEPTH_MAX1) )
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
DO JK=1,IKT
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
GWORK2(:,JK) = GWORK1(:)

RODIER Quentin
committed
IF(PARAMMF%LVERLIMUP) THEN

RODIER Quentin
committed
ZCOEF(:,JK) = (1.-(KDEPTH(:,JK)-ZDEPTH_MAX1)/(ZDEPTH_MAX2-ZDEPTH_MAX1))

RODIER Quentin
committed
ELSE

RODIER Quentin
committed
ZCOEF(:,JK) = (1.-(PDEPTH(:)-ZDEPTH_MAX1)/(ZDEPTH_MAX2-ZDEPTH_MAX1))

RODIER Quentin
committed
END IF

RODIER Quentin
committed
ZCOEF(:,JK)=MIN(MAX(ZCOEF(:,JK),0.),1.)
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
!$mnh_expand_where(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
WHERE (GWORK2(:,:))
PEMF(:,:) = PEMF(:,:) * ZCOEF(:,:)
PFRAC_UP(:,:) = PFRAC_UP(:,:) * ZCOEF(:,:)
!$mnh_end_expand_where(JIJ=IIJB:IIJE,JK=1:IKT)
ENDIF
IF (LHOOK) CALL DR_HOOK('COMPUTE_UPDRAFT',1,ZHOOK_HANDLE)
SUBROUTINE COMPUTE_ENTR_DETR(D, CST, NEBN, PARAMMF,&

RIETTE Sébastien
committed
KK,KKB,KKE,KKL,OTEST,OTESTLCL,&

RIETTE Sébastien
committed
PFRAC_ICE,PRHODREF,&

RIETTE Sébastien
committed
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
PPRE_MINUS_HALF,&
PPRE_PLUS_HALF,PZZ,PDZZ,&
PTHVM,PTHLM,PRTM,PW_UP2,PTH_UP,&
PTHL_UP,PRT_UP,PLUP,&
PRC_UP,PRI_UP,PTHV_UP,&
PRSAT_UP,PRC_MIX,PRI_MIX, &
PENTR,PDETR,PENTR_CLD,PDETR_CLD,&
PBUO_INTEG_DRY,PBUO_INTEG_CLD,&
PPART_DRY)
! #############################################################
!!
!!***COMPUTE_ENTR_DETR* - calculates caracteristics of the updraft or downdraft
!! using model of the EDMF scheme
!!
!! PURPOSE
!! -------
!!**** The purpose of this routine is to compute entrainement and
!! detrainement at one level of the updraft
!
!!** METHOD
!! ------
!!
!! EXTERNAL
!! --------
!!
!! IMPLICIT ARGUMENTS
!! ------------------
!!
!! REFERENCE
!! ---------
!! Book 1 of Meso-NH documentation (chapter Convection)
!!
!!
!! AUTHOR
!! ------
!! J.Pergaud : 2009
!!
!! MODIFICATIONS
!! -------------
!! Y.Seity (06/2010) Bug correction
!! V.Masson (09/2010) Optimization
!! S. Riette april 2011 : ice added, protection against zero divide by Yves Bouteloup
!! protection against too big ZPART_DRY, interface modified
!! S. Riette Jan 2012: support for both order of vertical levels
!! S. Riette & J. Escobar (11/2013) : remove div by 0 on real*4 case
!! P.Marguinaud Jun 2012: fix uninitialized variable
!! P.Marguinaud Nov 2012: fix gfortran bug
!! S. Riette Apr 2013: bugs correction, rewriting (for optimisation) and
!! improvement of continuity at the condensation level
!! S. Riette Nov 2013: protection against zero divide for min value of dry PDETR
!! R.Honnert Oct 2016 : Update with AROME
! P. Wautelet 08/02/2019: bugfix: compute ZEPSI_CLOUD only once and only when it is needed
!! R. El Khatib 29-Apr-2019 portability fix : compiler may get confused by embricked WHERE statements
!! eventually breaking tests with NaN initializations at compile time.
!! Replace by IF conditions and traditional DO loops can only improve the performance.
! P. Wautelet 10/02/2021: bugfix: initialized PPART_DRY everywhere
!! --------------------------------------------------------------------------
!
!* 0. DECLARATIONS
! ------------
!
USE MODD_DIMPHYEX, ONLY: DIMPHYEX_t
USE MODD_CST, ONLY: CST_t

RIETTE Sébastien
committed
USE MODD_PARAM_MFSHALL_n, ONLY: PARAM_MFSHALL_t
!
IMPLICIT NONE
!
!
!* 1.1 Declaration of Arguments
!
!
TYPE(DIMPHYEX_t), INTENT(IN) :: D
TYPE(CST_t), INTENT(IN) :: CST

RIETTE Sébastien
committed
TYPE(PARAM_MFSHALL_t), INTENT(IN) :: PARAMMF
!
INTEGER, INTENT(IN) :: KK
INTEGER, INTENT(IN) :: KKB ! near ground physical index
INTEGER, INTENT(IN) :: KKE ! uppest atmosphere physical index
INTEGER, INTENT(IN) :: KKL ! +1 if grid goes from ground to atmosphere top, -1 otherwise
LOGICAL,DIMENSION(D%NIJT), INTENT(IN) :: OTEST ! test to see if updraft is running
LOGICAL,DIMENSION(D%NIJT), INTENT(IN) :: OTESTLCL !test of condensation
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PFRAC_ICE ! fraction of ice

RIETTE Sébastien
committed
!
! prognostic variables at t- deltat
!
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PRHODREF !rhodref
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PPRE_MINUS_HALF ! Pressure at flux level KK
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PPRE_PLUS_HALF ! Pressure at flux level KK+KKL
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PZZ ! Height at the flux point
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PDZZ ! metrics coefficient
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PTHVM ! ThetaV environment

RIETTE Sébastien
committed
!
! thermodynamical variables which are transformed in conservative var.
!
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PTHLM ! Thetal
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PRTM ! total mixing ratio
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PW_UP2 ! Vertical velocity^2
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PTH_UP,PTHL_UP,PRT_UP ! updraft properties
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PLUP ! LUP compute from the ground
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PRC_UP,PRI_UP ! Updraft cloud content
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PTHV_UP ! Thetav of updraft
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PRSAT_UP ! Mixing ratio at saturation in updraft
REAL, DIMENSION(D%NIJT), INTENT(INOUT) :: PRC_MIX, PRI_MIX ! Mixture cloud content
REAL, DIMENSION(D%NIJT), INTENT(OUT) :: PENTR ! Mass flux entrainment of the updraft
REAL, DIMENSION(D%NIJT), INTENT(OUT) :: PDETR ! Mass flux detrainment of the updraft
REAL, DIMENSION(D%NIJT), INTENT(OUT) :: PENTR_CLD ! Mass flux entrainment of the updraft in cloudy part
REAL, DIMENSION(D%NIJT), INTENT(OUT) :: PDETR_CLD ! Mass flux detrainment of the updraft in cloudy part
REAL, DIMENSION(D%NIJT), INTENT(OUT) :: PBUO_INTEG_DRY, PBUO_INTEG_CLD! Integral Buoyancy
REAL, DIMENSION(D%NIJT), INTENT(OUT) :: PPART_DRY ! ratio of dry part at the transition level

RIETTE Sébastien
committed
!
!
! 1.2 Declaration of local variables
!
! Local array declaration must be put in the compute_updraft subroutine
! For simplicity all local variables (including scalars) are moved in the compute_updraft subroutine

RIETTE Sébastien
committed
!
!----------------------------------------------------------------------------------
! 1.3 Initialisation
! ------------------
ZCOEFFMF_CLOUD=PARAMMF%XENTR_MF * CST%XG / PARAMMF%XCRAD_MF
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZG_O_THVREF_ED(:)=CST%XG/PTHVM(:,KK)

RODIER Quentin
committed
ZFRAC_ICE(:)=PFRAC_ICE(:) ! to not modify fraction of ice

RODIER Quentin
committed
ZPRE(:)=PPRE_MINUS_HALF(:)
!$mnh_end_expand_array(JIJ=IIJB:IIJE)

RIETTE Sébastien
committed
! 1.4 Estimation of PPART_DRY

RODIER Quentin
committed
DO JIJ=IIJB,IIJE
IF(OTEST(JIJ) .AND. OTESTLCL(JIJ)) THEN
!No dry part when condensation level is reached

RODIER Quentin
committed
PPART_DRY(JIJ)=0.
ZDZ_STOP(JIJ)=0.
ZPRE(JIJ)=PPRE_MINUS_HALF(JIJ)
ELSE IF (OTEST(JIJ) .AND. .NOT. OTESTLCL(JIJ)) THEN
!Temperature at flux level KK

RODIER Quentin
committed
ZT=PTH_UP(JIJ)*(PPRE_MINUS_HALF(JIJ)/CST%XP00) ** (CST%XRD/CST%XCPD)
!Saturating vapor pressure at flux level KK

RODIER Quentin
committed
ZFOESW = MIN(EXP( CST%XALPW - CST%XBETAW/ZT - CST%XGAMW*LOG(ZT) ), 0.99*PPRE_MINUS_HALF(JIJ))
ZFOESI = MIN(EXP( CST%XALPI - CST%XBETAI/ZT - CST%XGAMI*LOG(ZT) ), 0.99*PPRE_MINUS_HALF(JIJ))
!Computation of d.Rsat / dP (partial derivations with respect to P and T
!and use of T=Theta*(P/P0)**(R/Cp) to transform dT into dP with theta_up
!constant at the vertical)

RODIER Quentin
committed
ZDRSATODP=(CST%XBETAW/ZT-CST%XGAMW)*(1-ZFRAC_ICE(JIJ))+(CST%XBETAI/ZT-CST%XGAMI)*ZFRAC_ICE(JIJ)
ZDRSATODP=((CST%XRD/CST%XCPD)*ZDRSATODP-1.)*PRSAT_UP(JIJ)/ &
&(PPRE_MINUS_HALF(JIJ)-(ZFOESW*(1-ZFRAC_ICE(JIJ)) + ZFOESI*ZFRAC_ICE(JIJ)))
!Use of d.Rsat / dP and pressure at flux level KK to find pressure (ZPRE)
!where Rsat is equal to PRT_UP

RODIER Quentin
committed
ZPRE(JIJ)=PPRE_MINUS_HALF(JIJ)+(PRT_UP(JIJ)-PRSAT_UP(JIJ))/ZDRSATODP
!Fraction of dry part (computed with pressure and used with heights, no
!impact found when using log function here and for pressure on flux levels
!computation)

RODIER Quentin
committed
PPART_DRY(JIJ)=MAX(0., MIN(1., (PPRE_MINUS_HALF(JIJ)-ZPRE(JIJ))/(PPRE_MINUS_HALF(JIJ)-PPRE_PLUS_HALF(JIJ))))
!Height above flux level KK of the cloudy part

RODIER Quentin
committed
ZDZ_STOP(JIJ) = (PZZ(JIJ,KK+KKL)-PZZ(JIJ,KK))*PPART_DRY(JIJ)

RODIER Quentin
committed
PPART_DRY(JIJ)=0. ! value does not matter, here
END IF
END DO

RIETTE Sébastien
committed
! 1.5 Gradient and flux values of thetav
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZCOEFF_MINUS_HALF(:)=((PTHVM(:,KK)-PTHVM(:,KK-KKL))/PDZZ(:,KK))
ZTHV_MINUS_HALF(:) = PTHVM(:,KK) - &
& ZCOEFF_MINUS_HALF(:)*0.5*(PZZ(:,KK+KKL)-PZZ(:,KK))

RODIER Quentin
committed
ZCOEFF_MINUS_HALF(:)=0.
ZTHV_MINUS_HALF(:) = PTHVM(:,KK)

RODIER Quentin
committed
ZCOEFF_PLUS_HALF(:) = ((PTHVM(:,KK+KKL)-PTHVM(:,KK))/PDZZ(:,KK+KKL))
ZTHV_PLUS_HALF(:) = PTHVM(:,KK) + &
& ZCOEFF_PLUS_HALF(:)*0.5*(PZZ(:,KK+KKL)-PZZ(:,KK))
!$mnh_end_expand_array(JIJ=IIJB:IIJE)

RIETTE Sébastien
committed
! 2 Dry part computation:
! Integral buoyancy and computation of PENTR and PDETR for dry part
! --------------------------------------------------------------------

RODIER Quentin
committed
DO JIJ=IIJB,IIJE
IF (OTEST(JIJ) .AND. PPART_DRY(JIJ)>0.) THEN

RIETTE Sébastien
committed
!Buoyancy computation in two parts to use change of gradient of theta v of environment
!Between flux level KK and min(mass level, bottom of cloudy part)

RODIER Quentin
committed
ZDZ=MIN(ZDZ_STOP(JIJ),(PZZ(JIJ,KK+KKL)-PZZ(JIJ,KK))*0.5)
PBUO_INTEG_DRY(JIJ) = ZG_O_THVREF_ED(JIJ)*ZDZ*&
(0.5 * ( - ZCOEFF_MINUS_HALF(JIJ))*ZDZ &
- ZTHV_MINUS_HALF(JIJ) + PTHV_UP(JIJ) )

RIETTE Sébastien
committed
!Between mass flux KK and bottom of cloudy part (if above mass flux)

RODIER Quentin
committed
ZDZ=MAX(0., ZDZ_STOP(JIJ)-(PZZ(JIJ,KK+KKL)-PZZ(JIJ,KK))*0.5)
PBUO_INTEG_DRY(JIJ) = PBUO_INTEG_DRY(JIJ) + ZG_O_THVREF_ED(JIJ)*ZDZ*&
(0.5 * ( - ZCOEFF_PLUS_HALF(JIJ))*ZDZ &
- PTHVM(JIJ,KK) + PTHV_UP(JIJ) )

RIETTE Sébastien
committed
!Entr//Detr. computation

RODIER Quentin
committed
IF (PBUO_INTEG_DRY(JIJ)>=0.) THEN
PENTR(JIJ) = 0.5/(PARAMMF%XABUO-PARAMMF%XBENTR*PARAMMF%XENTR_DRY)*&
LOG(1.+ (2.*(PARAMMF%XABUO-PARAMMF%XBENTR*PARAMMF%XENTR_DRY)/PW_UP2(JIJ,KK))* &
PBUO_INTEG_DRY(JIJ))
PDETR(JIJ) = 0.

RIETTE Sébastien
committed
ELSE

RODIER Quentin
committed
PENTR(JIJ) = 0.
PDETR(JIJ) = 0.5/(PARAMMF%XABUO)*&
LOG(1.+ (2.*(PARAMMF%XABUO)/PW_UP2(JIJ,KK))* &
(-PBUO_INTEG_DRY(JIJ)))

RIETTE Sébastien
committed
ENDIF

RODIER Quentin
committed
PENTR(JIJ) = PARAMMF%XENTR_DRY*PENTR(JIJ)/(PZZ(JIJ,KK+KKL)-PZZ(JIJ,KK))
PDETR(JIJ) = PARAMMF%XDETR_DRY*PDETR(JIJ)/(PZZ(JIJ,KK+KKL)-PZZ(JIJ,KK))

RIETTE Sébastien
committed
!Minimum value of detrainment

RODIER Quentin
committed
ZWK0D=PLUP(JIJ)-0.5*(PZZ(JIJ,KK)+PZZ(JIJ,KK+KKL))
ZWK0D=SIGN(MAX(1., ABS(ZWK0D)), ZWK0D) ! ZWK0D must not be zero

RODIER Quentin
committed
PDETR(JIJ) = MAX(PPART_DRY(JIJ)*PARAMMF%XDETR_LUP/ZWK0D, PDETR(JIJ))

RIETTE Sébastien
committed
ELSE
!No dry part, condensation reached (OTESTLCL)

RODIER Quentin
committed
PBUO_INTEG_DRY(JIJ) = 0.
PENTR(JIJ)=0.
PDETR(JIJ)=0.

RIETTE Sébastien
committed
ENDIF
ENDDO
! 3 Wet part computation
! -----------------------
! 3.1 Integral buoyancy for cloudy part
! Compute theta_v of updraft at flux level KK+KKL
!MIX variables are used to avoid declaring new variables
!but we are dealing with updraft and not mixture
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZRCMIX(:)=PRC_UP(:)
ZRIMIX(:)=PRI_UP(:)
!$mnh_end_expand_array(JIJ=IIJB:IIJE)

RIETTE Sébastien
committed
CALL TH_R_FROM_THL_RT(CST,NEBN,D%NIJT,NEBN%CFRAC_ICE_SHALLOW_MF,ZFRAC_ICE,&
PPRE_PLUS_HALF,PTHL_UP,PRT_UP,&
ZTHMIX,ZRVMIX,ZRCMIX,ZRIMIX,&
ZRSATW_ED, ZRSATI_ED,OOCEAN=.FALSE.,&
PBUF=ZBUF, KB=D%NIJB, KE=D%NIJE)
!$mnh_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZTHV_UP_F2(:) = ZTHMIX(:)*(1.+ZRVORD*ZRVMIX(:))/(1.+PRT_UP(:))
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
! Integral buoyancy for cloudy part

RODIER Quentin
committed
DO JIJ=IIJB,IIJE
IF(OTEST(JIJ) .AND. PPART_DRY(JIJ)<1.) THEN
!Gradient of Theta V updraft over the cloudy part, assuming that thetaV updraft don't change
!between flux level KK and bottom of cloudy part

RODIER Quentin
committed
ZCOTHVU=(ZTHV_UP_F2(JIJ)-PTHV_UP(JIJ))/((PZZ(JIJ,KK+KKL)-PZZ(JIJ,KK))*(1-PPART_DRY(JIJ)))
!Computation in two parts to use change of gradient of theta v of environment
!Between bottom of cloudy part (if under mass level) and mass level KK

RODIER Quentin
committed
ZDZ=MAX(0., 0.5*(PZZ(JIJ,KK+KKL)-PZZ(JIJ,KK))-ZDZ_STOP(JIJ))
PBUO_INTEG_CLD(JIJ) = ZG_O_THVREF_ED(JIJ)*ZDZ*&
(0.5*( ZCOTHVU - ZCOEFF_MINUS_HALF(JIJ))*ZDZ &
- (PTHVM(JIJ,KK)-ZDZ*ZCOEFF_MINUS_HALF(JIJ)) + PTHV_UP(JIJ) )
!Between max(mass level, bottom of cloudy part) and flux level KK+KKL

RODIER Quentin
committed
ZDZ=(PZZ(JIJ,KK+KKL)-PZZ(JIJ,KK))-MAX(ZDZ_STOP(JIJ),0.5*(PZZ(JIJ,KK+KKL)-PZZ(JIJ,KK)))
PBUO_INTEG_CLD(JIJ) = PBUO_INTEG_CLD(JIJ)+ZG_O_THVREF_ED(JIJ)*ZDZ*&
(0.5*( ZCOTHVU - ZCOEFF_PLUS_HALF(JIJ))*ZDZ&
- (PTHVM(JIJ,KK)+(0.5*((PZZ(JIJ,KK+KKL)-PZZ(JIJ,KK)))-ZDZ)*ZCOEFF_PLUS_HALF(JIJ)) +&
PTHV_UP(JIJ) )

RIETTE Sébastien
committed
ELSE
!No cloudy part

RODIER Quentin
committed
PBUO_INTEG_CLD(JIJ)=0.
END IF
END DO