Newer
Older

RODIER Quentin
committed
!MNH_LIC Copyright 1994-2021 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.
!-----------------------------------------------------------------

RODIER Quentin
committed
SUBROUTINE TURB(CST,CSTURB,BUCONF,TURBN,D,TLES, &

RODIER Quentin
committed
& KMI,KRR,KRRL,KRRI,HLBCX,HLBCY,KGRADIENTS,KHALO, &

RODIER Quentin
committed
& KSPLIT,KMODEL_CL,KSV,KSV_LGBEG,KSV_LGEND,HPROGRAM, &
& KSV_LIMA_NR, KSV_LIMA_NS, KSV_LIMA_NG, KSV_LIMA_NH, &
& O2D,ONOMIXLG,OFLAT,OCOUPLES,OBLOWSNOW,OIBM,OFLYER, &

RODIER Quentin
committed
& OCOMPUTE_SRC, PRSNOW, &
& OOCEAN,ODEEPOC,ODIAG_IN_RUN, &
& HTURBLEN_CL,HCLOUD, &
& PTSTEP,TPFILE, &
& PDXX,PDYY,PDZZ,PDZX,PDZY,PZZ, &
& PDIRCOSXW,PDIRCOSYW,PDIRCOSZW,PCOSSLOPE,PSINSLOPE, &

RODIER Quentin
committed
& PRHODJ,PTHVREF,PHGRAD,PZS, &
& PSFTH,PSFRV,PSFSV,PSFU,PSFV, &
& PPABST,PUT,PVT,PWT,PTKET,PSVT,PSRCT, &
& PLENGTHM,PLENGTHH,MFMOIST, &
& PBL_DEPTH,PSBL_DEPTH, &
& PCEI,PCEI_MIN,PCEI_MAX,PCOEF_AMPL_SAT, &
& PRUS,PRVS,PRWS,PRTHLS,PRRS,PRSVS,PRTKES, &
& PSIGS, &
& PFLXZTHVMF,PWTH,PWRC,PWSV,PDP,PTP,PTDIFF,PTDISS, &
& TBUDGETS, KBUDGETS, &
& PEDR,PLEM,PRTKEMS,PTPMF, &
& PDRUS_TURB,PDRVS_TURB, &

RODIER Quentin
committed
& PDRTHLS_TURB,PDRRTS_TURB,PDRSVS_TURB,PTR,PDISS, &
& PIBM_LS, PIBM_XMUT, &
& PCURRENT_TKE_DISS, PSSTFL, PSSTFL_C, PSSRFL_C, &
& PSSUFL_C, PSSVFL_C,PSSUFL,PSSVFL )
! #################################################################
!
!
!!**** *TURB* - computes the turbulent source terms for the prognostic
!!
!! PURPOSE
!! -------
!!**** The purpose of this routine is to compute the source terms in
!! the evolution equations due to the turbulent mixing.
!! The source term is computed as the divergence of the turbulent fluxes.
!! The cartesian fluxes are obtained by a one and a half order closure, based
!! on a prognostic equation for the Turbulence Kinetic Energy( TKE ). The
!! system is closed by prescribing a turbulent mixing length. Different
!! choices are available for this length.
!
!!** METHOD
!! ------
!! The dimensionality of the turbulence parameterization can be chosen by

RODIER Quentin
committed
!! means of the parameter TURBN%CTURBDIM:
!! * TURBN%CTURBDIM='1DIM' the parameterization is 1D but can be used in
!! 3D , 2D or 1D simulations. Only the sources associated to the vertical
!! turbulent fluxes are taken into account.

RODIER Quentin
committed
!! * TURBN%CTURBDIM='3DIM' the parameterization is fully 2D or 3D depending
!! on the model dimensionality. Of course, it does not make any sense to
!! activate this option with a 1D model.
!!
!! The following steps are made:
!! 1- Preliminary computations.
!! 2- The metric coefficients are recovered from the grid knowledge.
!! 3- The mixing length is computed according to its choice:

RODIER Quentin
committed
!! * TURBN%CTURBLEN='BL89' the Bougeault and Lacarrere algorithm is used.
!! The mixing length is given by the vertical displacement from its
!! original level of an air particule having an initial internal
!! energy equal to its TKE and stopped by the buoyancy forces.
!! The discrete formulation is second order accurate.

RODIER Quentin
committed
!! * TURBN%CTURBLEN='DELT' the mixing length is given by the mesh size
!! depending on the model dimensionality, this length is limited
!! with the ground distance.

RODIER Quentin
committed
!! * TURBN%CTURBLEN='DEAR' the mixing length is given by the mesh size
!! depending on the model dimensionality, this length is limited
!! with the ground distance and also by the Deardorff mixing length
!! pertinent in the stable cases.

RODIER Quentin
committed
!! * TURBN%CTURBLEN='KEPS' the mixing length is deduced from the TKE
!! dissipation, which becomes a prognostic variable of the model (
!! 3'- The cloud mixing length is computed according to HTURBLEN_CLOUD
!! and emphasized following the CEI index
!! 4- The conservative variables are computed along with Lv/Cp.
!! 5- The turbulent Prandtl numbers are computed from the resolved fields
!! 6- The sources associated to the vertical turbulent fluxes are computed
!! with a temporal scheme allowing a degree of implicitness given by

RODIER Quentin
committed
!! TURBN%XIMPL, varying from TURBN%XIMPL=0. ( purely explicit scheme) to TURBN%XIMPL=1.
!! ( purely implicit scheme)
!! The sources associated to the horizontal fluxes are computed with a
!! purely explicit temporal scheme. These sources are only computed when

RODIER Quentin
committed
!! the turbulence parameterization is 2D or 3D( TURBN%CTURBDIM='3DIM' ).
!! 7- The sources for TKE are computed, along with the dissipation of TKE

RODIER Quentin
committed
!! if TURBN%CTURBLEN='KEPS'.
!! 8- Some turbulence-related quantities are stored in the synchronous
!! FM-file.
!! 9- The non-conservative variables are retrieved.
!!
!!
!! The saving of the fields in the synchronous FM-file is controlled by:

RODIER Quentin
committed
!! * TURBN%LTURB_FLX => saves all the turbulent fluxes and correlations
!! * TURBN%LTURB_DIAG=> saves the turbulent Prandtl and Schmidt numbers, the
!! source terms of TKE and dissipation of TKE
!!
!! EXTERNAL
!! --------
!! SUBROUTINE PRANDTL : computes the turbulent Prandtl number
!! SUBROUTINE TURB_VER : computes the sources from the vertical fluxes
!! SUBROUTINE TURB_HOR : computes the sources from the horizontal fluxes
!! SUBROUTINE TKE_EPS_SOURCES : computes the sources for TKE and its
!! dissipation
!! SUBROUTINE BUDGET : computes and stores the budgets
!!
!! IMPLICIT ARGUMENTS
!! ------------------
!!
!! MODD_PARAMETERS : JPVEXT_TURB number of marginal vertical points
!!
!! MODD_CONF : CCONF model configuration (start/restart)
!! L1D switch for 1D model version
!! L2D switch for 2D model version
!!
!! MODD_CST : contains physical constants

RODIER Quentin
committed
!! CST%XG gravity constant
!! CST%XRD Gas constant for dry air
!! CST%XRV Gas constant for vapor
!!
!! MODD_CTURB : contains turbulence scheme constants
!! XCMFS,XCED to compute the dissipation mixing length
!! XTKEMIN minimum values for the TKE
!! CST%XLINI,CST%XLINF to compute Bougeault-Lacarrere mixing
!! length
!! Module MODD_BUDGET:
!! NBUMOD
!! CBUTYPE
!! LBU_RU
!! LBU_RV
!! LBU_RW
!! LBU_RTH
!! LBU_RSV1
!! LBU_RRV
!! LBU_RRC
!! LBU_RRR
!! LBU_RRI
!! LBU_RRS
!! LBU_RRG
!! LBU_RRH
!!
!! REFERENCE
!! ---------
!! Book 2 of documentation (routine TURB)
!! Book 1 of documentation (Chapter: Turbulence)
!!
!! AUTHOR
!! ------
!! Joan Cuxart * INM and Meteo-France *
!!
!! MODIFICATIONS
!! -------------
!! Original 05/10/94
!! Modifications: Feb 14, 1995 (J.Cuxart and J.Stein)
!! Doctorization and Optimization
!! Modifications: March 21, 1995 (J.M. Carriere)
!! Introduction of cloud water
!! Modifications: June 1, 1995 (J.Cuxart )
!! take min(Kz,delta)
!! Modifications: June 1, 1995 (J.Stein J.Cuxart)
!! remove unnecessary arrays and change Prandtl
!! and Schmidt numbers localizations
!! Modifications: July 20, 1995 (J.Stein) remove MODI_ground_ocean +
!! TZDTCUR + MODD_TIME because they are not used
!! change RW in RNP for the outputs
!! Modifications: August 21, 1995 (Ph. Bougeault)
!! take min(K(z-zsol),delta)
!! Modifications: Sept 14, 1995 (Ph Bougeault, J. Cuxart)
!! second order BL89 mixing length computations + add Deardorff length
!! in the Delta case for stable cases
!! Modifications: Sept 19, 1995 (J. Stein, J. Cuxart)
!! define a DEAR case for the mixing length, add MODI_BUDGET and change
!! some BUDGET calls, add LES tools
!! Modifications: Oct 16, 1995 (J. Stein) change the budget calls
!! Modifications: Feb 28, 1996 (J. Stein) optimization +
!! remove min(K(z-zsol),delta)+
!! bug in the tangential fluxes
!! Modifications: Oct 16, 1996 (J. Stein) change the subgrid condensation
!! scheme + temporal discretization
!! Modifications: Dec 19, 1996 (J.-P. Pinty) update the budget calls
!! Jun 22, 1997 (J. Stein) use the absolute pressure and
!! change the Deardorf length at the surface
!! Modifications: Apr 27, 1997 (V. Masson) BL89 mix. length computed in
!! a separate routine
!! Oct 13, 1999 (J. Stein) switch for the tgt fluxes
!! Jun 24, 1999 (P Jabouille) Add routine UPDATE_ROTATE_WIND
!! Feb 15, 2001 (J. Stein) remove tgt fluxes
!! Mar 8, 2001 (V. Masson) forces the same behaviour near the surface
!! for all mixing lengths
!! Nov 06, 2002 (V. Masson) LES budgets
!! Nov, 2002 (V. Masson) implement modifications of
!! mixing and dissipative lengths
!! near the surface (according
!! Redelsperger et al 2001)
!! Apr, 2003 (V. Masson) bug in Blackadar length
!! bug in LES in 1DIM case
!! Feb 20, 2003 (J.-P. Pinty) Add reversible ice processes
!! May,26 2004 (P Jabouille) coef for computing dissipative heating
!! Sept 2004 (M.Tomasini) Cloud Mixing length modification
!! following the instability
!! criterium CEI calculated in modeln
!! May 2006 Remove KEPS
!! Sept.2006 (I.Sandu): Modification of the stability criterion for
!! DEAR (theta_v -> theta_l)
!! Oct 2007 (J.Pergaud) Add MF contribution for vert. turb. transport
!! Oct.2009 (C.Lac) Introduction of different PTSTEP according to the
!! advection schemes
!! October 2009 (G. Tanguy) add ILENCH=LEN(YCOMMENT) after
!! change of YCOMMENT
!! 06/2011 (J.escobar ) Bypass Bug with ifort11/12 on HLBCX,HLBC
!! 2012-02 Y. Seity, add possibility to run with reversed
!! vertical levels
!! 10/2012 (J. Colin) Correct bug in DearDoff for dry simulations
!! 10/2012 J.Escobar Bypass PGI bug , redefine some allocatable array inplace of automatic
!! 2014-11 Y. Seity, add output terms for TKE DDHs budgets
!! July 2015 (Wim de Rooy) modifications to run with RACMO

RODIER Quentin
committed
!! turbulence (TURBN%LHARAT=TRUE)
!! 04/2016 (C.Lac) correction of negativity for KHKO
! P. Wautelet 05/2016-04/2018: new data structures and calls for I/O
! Q. Rodier 01/2018: introduction of RM17
! P. Wautelet 20/05/2019: add name argument to ADDnFIELD_ll + new ADD4DFIELD_ll subroutine
! Wim de Rooy 06/2019: update statistical cloud scheme
! P. Wautelet 02/2020: use the new data structures and subroutines for budgets
! B. Vie 03/2020: LIMA negativity checks after turbulence, advection and microphysics budgets
! P. Wautelet 11/06/2020: bugfix: correct PRSVS array indices
! P. Wautelet + Benoit Vié 06/2020: improve removal of negative scalar variables + adapt the corresponding budgets
! P. Wautelet 30/06/2020: move removal of negative scalar variables to Sources_neg_correct
! R. Honnert/V. Masson 02/2021: new mixing length in the grey zone
! J.L. Redelsperger 03/2021: add Ocean LES case
! --------------------------------------------------------------------------
!
!* 0. DECLARATIONS
! ------------
!
USE PARKIND1, ONLY: JPRB
USE MODE_SHUMAN_PHY, ONLY: MZF_PHY,MXF_PHY,MYF_PHY
USE YOMHOOK , ONLY: LHOOK, DR_HOOK

RODIER Quentin
committed
!
USE MODD_BUDGET, ONLY: NBUDGET_U, NBUDGET_V, NBUDGET_W, NBUDGET_TH, NBUDGET_RV, NBUDGET_RC, &
NBUDGET_RR, NBUDGET_RI, NBUDGET_RS, NBUDGET_RG, NBUDGET_RH, NBUDGET_SV1, &

RODIER Quentin
committed
TBUDGETDATA, TBUDGETCONF_t
USE MODD_CST, ONLY: CST_t
USE MODD_CTURB, ONLY: CSTURB_t
USE MODD_DIMPHYEX, ONLY: DIMPHYEX_t
USE MODD_FIELD, ONLY: TFIELDMETADATA, TYPEREAL
USE MODD_IO, ONLY: TFILEDATA
USE MODD_LES, ONLY: TLES_t
USE MODD_PARAMETERS, ONLY: JPVEXT_TURB, XUNDEF
USE MODD_TURB_n, ONLY: TURB_t
!
USE MODE_BL89, ONLY: BL89
USE MODE_BUDGET_PHY, ONLY: BUDGET_STORE_INIT_PHY, BUDGET_STORE_END_PHY
USE MODE_EMOIST, ONLY: EMOIST
USE MODE_ETHETA, ONLY: ETHETA
USE MODE_GRADIENT_U_PHY, ONLY: GZ_U_UW_PHY
USE MODE_GRADIENT_V_PHY, ONLY: GZ_V_VW_PHY
USE MODE_GRADIENT_W_PHY, ONLY: GZ_W_M_PHY
USE MODE_GRADIENT_M_PHY, ONLY: GZ_M_W_PHY
USE MODE_IBM_MIXINGLENGTH, ONLY: IBM_MIXINGLENGTH

RODIER Quentin
committed
USE MODE_IO_FIELD_WRITE_PHY, ONLY: IO_FIELD_WRITE_PHY
USE MODE_RMC01, ONLY: RMC01
USE MODE_ROTATE_WIND, ONLY: ROTATE_WIND, UPDATE_ROTATE_WIND
USE MODE_SBL_PHY, ONLY: LMO
USE MODE_SOURCES_NEG_CORRECT, ONLY: SOURCES_NEG_CORRECT_PHY
USE MODE_TM06, ONLY: TM06
USE MODE_TKE_EPS_SOURCES, ONLY: TKE_EPS_SOURCES
USE MODE_TURB_HOR_SPLT, ONLY: TURB_HOR_SPLT
USE MODE_TURB_VER, ONLY: TURB_VER
USE MODE_UPDATE_LM, ONLY: UPDATE_LM
USE MODI_SECOND_MNH, ONLY: SECOND_MNH

RODIER Quentin
committed
!
IMPLICIT NONE
!
!
!* 0.1 declarations of arguments
!
!
!

RODIER Quentin
committed
TYPE(DIMPHYEX_t), INTENT(IN) :: D ! PHYEX variables dimensions structure
TYPE(CST_t), INTENT(IN) :: CST ! modd_cst general constant structure
TYPE(CSTURB_t), INTENT(IN) :: CSTURB ! modd_csturb turb constant structure
TYPE(TBUDGETCONF_t), INTENT(IN) :: BUCONF ! budget structure
TYPE(TURB_t), INTENT(IN) :: TURBN ! modn_turbn (turb namelist) structure
TYPE(TLES_t), INTENT(INOUT) :: TLES ! modd_les structure

RODIER Quentin
committed
INTEGER, INTENT(IN) :: KGRADIENTS ! Number of stored horizontal gradients
INTEGER, INTENT(IN) :: KMI ! model index number
INTEGER, INTENT(IN) :: KRR ! number of moist var.
INTEGER, INTENT(IN) :: KRRL ! number of liquid water var.
INTEGER, INTENT(IN) :: KRRI ! number of ice water var.

RODIER Quentin
committed
INTEGER, INTENT(IN) :: KSV, KSV_LGBEG, KSV_LGEND ! number of scalar variables
INTEGER, INTENT(IN) :: KSV_LIMA_NR,KSV_LIMA_NS,KSV_LIMA_NG,KSV_LIMA_NH
CHARACTER(LEN=4),DIMENSION(2),INTENT(IN):: HLBCX, HLBCY ! X- and Y-direc LBC
INTEGER, INTENT(IN) :: KSPLIT ! number of time-splitting
INTEGER, INTENT(IN) :: KMODEL_CL ! model number for cloud mixing length
LOGICAL, INTENT(IN) :: OCOMPUTE_SRC ! flag to define dimensions of SIGS and SRCT variables
LOGICAL, INTENT(IN) :: OOCEAN ! switch for Ocean model version

RODIER Quentin
committed
LOGICAL, INTENT(IN) :: ODEEPOC ! activates sfc forcing for ideal ocean deep conv
LOGICAL, INTENT(IN) :: OFLYER ! MesoNH flyer diagnostic

RODIER Quentin
committed
LOGICAL, INTENT(IN) :: OFLAT ! Logical for zero ororography

RODIER Quentin
committed
LOGICAL, INTENT(IN) :: OCOUPLES ! switch to activate atmos-ocean LES version
LOGICAL, INTENT(IN) :: OBLOWSNOW ! switch to activate pronostic blowing snow

RODIER Quentin
committed
LOGICAL, INTENT(IN) :: ODIAG_IN_RUN ! switch to activate online diagnostics (mesonh)
LOGICAL, INTENT(IN) :: OIBM ! switch to modity mixing length near building with IBM
CHARACTER(LEN=4), INTENT(IN) :: HTURBLEN_CL ! kind of cloud mixing length
CHARACTER (LEN=4), INTENT(IN) :: HCLOUD ! Kind of microphysical scheme

RODIER Quentin
committed
INTEGER, INTENT(IN) :: KHALO ! Size of the halo for parallel distribution

RODIER Quentin
committed
REAL, INTENT(IN) :: PRSNOW ! Ratio for diffusion coeff. scalar (blowing snow)
REAL, INTENT(IN) :: PTSTEP ! timestep

RODIER Quentin
committed
TYPE(TFILEDATA), INTENT(INOUT) :: TPFILE ! Output file
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PDXX,PDYY,PDZZ,PDZX,PDZY
! metric coefficients
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PZZ ! physical distance
! between 2 succesive grid points along the K direction
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PDIRCOSXW, PDIRCOSYW, PDIRCOSZW
! Director Cosinus along x, y and z directions at surface w-point
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PCOSSLOPE ! cosinus of the angle
! between i and the slope vector
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PSINSLOPE ! sinus of the angle
! between i and the slope vector

RODIER Quentin
committed
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PZS ! orography (for LEONARD terms)
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PRHODJ ! dry density * Grid size
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: MFMOIST ! moist mass flux dual scheme
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PTHVREF ! Virtual Potential
! Temperature of the reference state
!
REAL, DIMENSION(D%NIJT), INTENT(IN) :: PSFTH,PSFRV, &
! normal surface fluxes of theta and Rv
PSFU,PSFV

RODIER Quentin
committed
! normal surface fluxes of (u,v) parallel to the orography
REAL, DIMENSION(D%NIJT,KSV), INTENT(IN) :: PSFSV
! normal surface fluxes of Scalar var.

RODIER Quentin
committed
REAL, DIMENSION(D%NIJT,D%NKT,KGRADIENTS), INTENT(IN) :: PHGRAD ! horizontal gradients
!
! prognostic variables at t- deltat
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PPABST ! Pressure at time t
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PUT,PVT,PWT ! wind components
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PTKET ! TKE
REAL, DIMENSION(D%NIJT,D%NKT,KSV), INTENT(IN) :: PSVT ! passive scal. var.
REAL, DIMENSION(MERGE(D%NIJT,0,OCOMPUTE_SRC),&
MERGE(D%NKT,0,OCOMPUTE_SRC)), INTENT(IN) :: PSRCT ! Second-order flux
! s'rc'/2Sigma_s2 at time t-1 multiplied by Lambda_3

RODIER Quentin
committed
REAL, DIMENSION(MERGE(D%NIJT,0,TURBN%CTOM=='TM06')),INTENT(INOUT) :: PBL_DEPTH ! BL height for TOMS
REAL, DIMENSION(MERGE(D%NIJT,0,TURBN%LRMC01)),INTENT(INOUT) :: PSBL_DEPTH ! SBL depth for RMC01
!
! variables for cloud mixing length
REAL, DIMENSION(MERGE(D%NIJT,0,KMODEL_CL==KMI .AND. HTURBLEN_CL/='NONE'),&
MERGE(D%NKT,0,KMODEL_CL==KMI .AND. HTURBLEN_CL/='NONE')),INTENT(IN) :: PCEI
! Cloud Entrainment instability
! index to emphasize localy
! turbulent fluxes
REAL, INTENT(IN) :: PCEI_MIN ! minimum threshold for the instability index CEI
REAL, INTENT(IN) :: PCEI_MAX ! maximum threshold for the instability index CEI
REAL, INTENT(IN) :: PCOEF_AMPL_SAT ! saturation of the amplification coefficient
!
! thermodynamical variables which are transformed in conservative var.
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(INOUT) :: PTHLT ! conservative pot. temp.
REAL, DIMENSION(D%NIJT,D%NKT,KRR), INTENT(INOUT) :: PRT ! water var. where
! PRT(:,:,:,1) is the conservative mixing ratio
! sources of momentum, conservative potential temperature, Turb. Kin. Energy,
! TKE dissipation
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(INOUT) :: PRUS,PRVS,PRWS,PRTHLS,PRTKES
! Source terms for all water kinds, PRRS(:,:,:,1) is used for the conservative
! mixing ratio
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN),OPTIONAL :: PRTKEMS
REAL, DIMENSION(D%NIJT,D%NKT,KRR), INTENT(INOUT) :: PRRS
! Source terms for all passive scalar variables
REAL, DIMENSION(D%NIJT,D%NKT,KSV), INTENT(INOUT) :: PRSVS
! Sigma_s at time t+1 : square root of the variance of the deviation to the

RODIER Quentin
committed
! saturation
REAL, DIMENSION(MERGE(D%NIJT,0,OCOMPUTE_SRC),&
MERGE(D%NKT,0,OCOMPUTE_SRC)), INTENT(OUT) :: PSIGS
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT),OPTIONAL :: PDRUS_TURB ! evolution of rhoJ*U by turbulence only
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT),OPTIONAL :: PDRVS_TURB ! evolution of rhoJ*V by turbulence only
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT),OPTIONAL :: PDRTHLS_TURB ! evolution of rhoJ*thl by turbulence only
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT),OPTIONAL :: PDRRTS_TURB ! evolution of rhoJ*rt by turbulence only
REAL, DIMENSION(D%NIJT,D%NKT,KSV), INTENT(OUT),OPTIONAL :: PDRSVS_TURB ! evolution of rhoJ*Sv by turbulence only
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PFLXZTHVMF
! MF contribution for vert. turb. transport
! used in the buoy. prod. of TKE

RODIER Quentin
committed
REAL, DIMENSION(MERGE(D%NIJT,0,OFLYER),MERGE(D%NKT,0,OFLYER)), INTENT(OUT) :: PWTH ! heat flux
REAL, DIMENSION(MERGE(D%NIJT,0,OFLYER),MERGE(D%NKT,0,OFLYER)), INTENT(OUT) :: PWRC ! cloud water flux
REAL, DIMENSION(MERGE(D%NIJT,0,OFLYER),MERGE(D%NKT,0,OFLYER),MERGE(KSV,0,OFLYER)),INTENT(OUT) :: PWSV ! scalar flux
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT) :: PTP ! Thermal TKE production
! MassFlux + turb
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT),OPTIONAL :: PTPMF ! Thermal TKE production
! MassFlux Only
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT) :: PDP ! Dynamic TKE production
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT) :: PTDIFF ! Diffusion TKE term
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT) :: PTDISS ! Dissipation TKE term
TYPE(TBUDGETDATA), DIMENSION(KBUDGETS), INTENT(INOUT) :: TBUDGETS
INTEGER, INTENT(IN) :: KBUDGETS

RODIER Quentin
committed
CHARACTER(LEN=6), INTENT(IN) :: HPROGRAM ! CPROGRAM is the program currently running (modd_conf)
LOGICAL, INTENT(IN) :: ONOMIXLG ! to use turbulence for lagrangian variables (modd_conf)
LOGICAL, INTENT(IN) :: O2D ! Logical for 2D model version (modd_conf)
!
! length scale from vdfexcu
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(IN) :: PLENGTHM, PLENGTHH
!
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT), OPTIONAL :: PEDR ! EDR
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT), OPTIONAL :: PLEM ! Mixing length
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT), OPTIONAL :: PTR ! Transport prod. of TKE
REAL, DIMENSION(D%NIJT,D%NKT), INTENT(OUT), OPTIONAL :: PDISS ! Dissipation of TKE

RODIER Quentin
committed
REAL, DIMENSION(MERGE(D%NIJT,0,ODIAG_IN_RUN),MERGE(D%NKT,0,ODIAG_IN_RUN)), INTENT(INOUT), OPTIONAL :: PCURRENT_TKE_DISS ! if ODIAG_IN_RUN in mesonh
REAL, DIMENSION(MERGE(D%NIJT,0,OCOUPLES)), INTENT(IN),OPTIONAL :: PSSTFL ! Time evol Flux of T at sea surface (LOCEAN)
REAL, DIMENSION(MERGE(D%NIJT,0,OCOUPLES)), INTENT(IN),OPTIONAL :: PSSTFL_C ! O-A interface flux for theta(LOCEAN and LCOUPLES)
REAL, DIMENSION(MERGE(D%NIJT,0,OCOUPLES)), INTENT(IN),OPTIONAL :: PSSRFL_C ! O-A interface flux for vapor (LOCEAN and LCOUPLES)
REAL, DIMENSION(MERGE(D%NIJT,0,OCOUPLES)), INTENT(IN),OPTIONAL :: PSSUFL_C ! Time evol Flux of U at sea surface (LOCEAN)
REAL, DIMENSION(MERGE(D%NIJT,0,OCOUPLES)), INTENT(IN),OPTIONAL :: PSSVFL_C !
REAL, DIMENSION(MERGE(D%NIJT,0,OCOUPLES)), INTENT(IN),OPTIONAL :: PSSUFL
REAL, DIMENSION(MERGE(D%NIJT,0,OCOUPLES)), INTENT(IN),OPTIONAL :: PSSVFL !

RODIER Quentin
committed
REAL, DIMENSION(MERGE(D%NIJT,0,OIBM),MERGE(D%NKT,0,OIBM)), INTENT(OUT), OPTIONAL :: PIBM_XMUT ! IBM turbulent viscosity
REAL, DIMENSION(MERGE(D%NIJT,0,OIBM),MERGE(D%NKT,0,OIBM)), INTENT(IN), OPTIONAL :: PIBM_LS ! IBM Level-set function
!
!-------------------------------------------------------------------------------
!
! 0.2 declaration of local variables
!
REAL, DIMENSION(D%NIJT,D%NKT) :: &
ZCP, & ! Cp at t-1
ZEXN, & ! EXN at t-1
ZT, & ! T at t-1
ZLOCPEXNM, & ! Lv/Cp/EXNREF at t-1

RODIER Quentin
committed
ZLM,ZLMW, & ! Turbulent mixing length (+ work array)
ZLEPS, & ! Dissipative length
ZATHETA,ZAMOIST, & ! coefficients for s = f (Thetal,Rnp)
ZCOEF_DISS, & ! 1/(Cph*Exner) for dissipative heating
ZFRAC_ICE, & ! ri fraction of rc+ri
ZMWTH,ZMWR,ZMTH2,ZMR2,ZMTHR,& ! 3rd order moments
ZFWTH,ZFWR,ZFTH2,ZFR2,ZFTHR,& ! opposite of verticale derivate of 3rd order moments
ZTHLM,ZRTKEMS, & ! initial potential temp; TKE advective source
ZSHEAR, ZDUDZ, ZDVDZ, & ! horizontal-wind vertical gradient
ZLVOCPEXNM,ZLSOCPEXNM, & ! Lv/Cp/EXNREF and Ls/Cp/EXNREF at t-1
ZATHETA_ICE,ZAMOIST_ICE, & ! coefficients for s = f (Thetal,Rnp)
ZRVSAT, ZDRVSATDT, & ! local array for routine compute_function_thermo

RODIER Quentin
committed
ZWORK1,ZWORK2, & ! working array syntax
ZETHETA,ZEMOIST, & ! coef ETHETA and EMOIST (for DEAR routine)
ZDTHLDZ,ZDRTDZ, & ! dtheta_l/dz, drt_dz used for computing the stablity criterion
ZCOEF_AMPL, & ! Amplification coefficient of the mixing length
! when the instability criterium is verified (routine CLOUD_MODIF_LM)
ZLM_CLOUD, & ! Turbulent mixing length in the clouds (routine CLOUD_MODIF_LM)
ZTEMP_BUD
!
!
REAL, DIMENSION(D%NIJT,D%NKT,KRR) :: ZRM ! initial mixing ratio
REAL, DIMENSION(D%NIJT) :: ZTAU11M,ZTAU12M, &
ZTAU22M,ZTAU33M, &
! tangential surface fluxes in the axes following the orography
ZUSLOPE,ZVSLOPE, &
! wind components at the first mass level parallel
! to the orography
ZCDUEFF, &
! - Cd*||u|| where ||u|| is the module of the wind tangential to
! orography (ZUSLOPE,ZVSLOPE) at the surface.
ZUSTAR, ZLMO, &
ZRVM, ZSFRV,ZWORK2D
! friction velocity, Monin Obuhkov length, work arrays for vapor
!
! Virtual Potential Temp. used
! in the Deardorff mixing length computation
!
!with LIMA, do not change rain, snow, graupel and hail concentrations (mixing ratio is not changed)
REAL, DIMENSION(D%NIJT,D%NKT,KSV) :: ZRSVS
!

RODIER Quentin
committed
REAL :: ZEXPL ! 1-TURBN%XIMPL deg of expl.
REAL :: ZRVORD ! RV/RD
REAL :: ZEPS ! XMV / XMD
REAL :: ZD ! distance to the surface (for routine DELT)
REAL :: ZVAR ! Intermediary variable (for routine DEAR)
REAL :: ZPENTE ! Slope of the amplification straight line (for routine CLOUD_MODIF_LM)
REAL :: ZCOEF_AMPL_CEI_NUL! Ordonnate at the origin of the
! amplification straight line (for routine CLOUD_MODIF_LM)
INTEGER :: IIJB,IIJE,IKB,IKE ! index value for the
INTEGER :: IINFO_ll ! return code of parallel routine
! Beginning and the End of the physical domain for the mass points
INTEGER :: IKT,IKA,IKU ! array size in k direction
INTEGER :: IKL
INTEGER :: IKTB,IKTE ! start, end of k loops in physical domain
INTEGER :: JRR,JK,JSV ! loop counters
REAL :: ZL0 ! Max. Mixing Length in Blakadar formula
REAL :: ZALPHA ! work coefficient :
! - proportionnality constant between Dz/2 and

RODIER Quentin
committed
! ! BL89 mixing length near the surface
!
REAL :: ZTIME1, ZTIME2
LOGICAL :: GOCEAN !Intermediate variable used to work around a Cray compiler bug (CCE 13.0.0)
LOGICAL :: GTURBLEN_BL89_TURBLEN_RM17_TURBLEN_HM21_ORMC01
TYPE(TFIELDMETADATA) :: TZFIELD
!
!* 1.PRELIMINARIES
! -------------
!
!* 1.1 Set the internal domains, ZEXPL

RODIER Quentin
committed
REAL(KIND=JPRB) :: ZHOOK_HANDLE,ZHOOK_HANDLE2
IF (LHOOK) CALL DR_HOOK('TURB',0,ZHOOK_HANDLE)

RODIER Quentin
committed
IF (TURBN%LHARAT .AND. TURBN%CTURBDIM /= '1DIM') THEN
CALL ABOR1('TURBN%LHARATU only implemented for option TURBN%CTURBDIM=1DIM!')

RODIER Quentin
committed
IF (TURBN%LHARAT .AND. TLES%LLES_CALL) THEN

RODIER Quentin
committed
CALL ABOR1('TURBN%LHARATU not implemented for option LLES_CALL')
IKTE=D%NKTE
IKB=D%NKB
IKE=D%NKE
IKA=D%NKA
IKU=D%NKU
IKL=D%NKL

RODIER Quentin
committed
ZEXPL = 1.- TURBN%XIMPL

RODIER Quentin
committed
ZRVORD= CST%XRV / CST%XRD
GOCEAN = OOCEAN
GTURBLEN_BL89_TURBLEN_RM17_TURBLEN_HM21_ORMC01 = &
TURBN%CTURBLEN=='BL89' .OR. TURBN%CTURBLEN=='RM17' .OR. TURBN%CTURBLEN=='HM21' .OR. TURBN%LRMC01

RODIER Quentin
committed
!Copy data into ZTHLM and ZRM only if needed
IF (GTURBLEN_BL89_TURBLEN_RM17_TURBLEN_HM21_ORMC01) THEN
!$acc kernels present_cr(ZTHLM,ZRM)

RODIER Quentin
committed
ZTHLM(:,:) = PTHLT(:,:)
ZRM(:,:,:) = PRT(:,:,:)

RODIER Quentin
committed
END IF
!Save LIMA scalar variables sources

RODIER Quentin
committed
ZRSVS(:,:,1:KSV)=PRSVS(:,:,1:KSV)
!
!
!----------------------------------------------------------------------------
!
!* 2. COMPUTE CONSERVATIVE VARIABLES AND RELATED QUANTITIES
! -----------------------------------------------------
!
!* 2.1 Cph at t
!
!$acc kernels present_cr(ZCOEF_DISS,ZTHLM,ZRM,zcp)
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZCP(:,:)=CST%XCPD

RODIER Quentin
committed
IF (KRR > 0) ZCP(:,:) = ZCP(:,:) + CST%XCPV * PRT(:,:,1)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
! PGI20.5 BUG or reproductibility problem , with pointer this loop on JRR parallelize whitout reduction
!$acc loop seq

RODIER Quentin
committed
DO JRR = 2,1+KRRL ! loop on the liquid components
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZCP(:,:) = ZCP(:,:) + CST%XCL * PRT(:,:,JRR)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
END DO
!
DO JRR = 2+KRRL,1+KRRL+KRRI ! loop on the solid components
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZCP(:,:) = ZCP(:,:) + CST%XCI * PRT(:,:,JRR)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
END DO
!
!* 2.2 Exner function at t
!
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZEXN(:,:) = 1.
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZEXN(:,:) = (PPABST(:,:)/CST%XP00) ** (CST%XRD/CST%XCPD)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!
!* 2.3 dissipative heating coeff a t
!
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
ZCOEF_DISS(:,:) = 1.0/(ZCP(:,:) * ZEXN(:,:))

RODIER Quentin
committed
ZFRAC_ICE(:,:) = 0.0
ZATHETA(:,:) = 0.0
ZAMOIST(:,:) = 0.0
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!$acc end kernels
!
IF (KRRL >=1) THEN
!
!* 2.4 Temperature at t
!
!$acc kernels present_cr(ZT)
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZT(:,:) = PTHLT(:,:) * ZEXN(:,:)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!
!* 2.5 Lv/Cph/Exn
!

RODIER Quentin
committed
IF (TURBN%LSTATNW) THEN
!wc call new functions depending on statnew
CALL COMPUTE_FUNCTION_THERMO_NEW_STAT(CST%XALPW,CST%XBETAW,CST%XGAMW,CST%XLVTT,CST%XCL,ZT,ZEXN,ZCP, &
ZLVOCPEXNM,ZAMOIST,ZATHETA)
CALL COMPUTE_FUNCTION_THERMO_NEW_STAT(CST%XALPI,CST%XBETAI,CST%XGAMI,CST%XLSTT,CST%XCI,ZT,ZEXN,ZCP, &
ZLSOCPEXNM,ZAMOIST_ICE,ZATHETA_ICE)
CALL COMPUTE_FUNCTION_THERMO(CST%XALPW,CST%XBETAW,CST%XGAMW,CST%XLVTT,CST%XCL,ZT,ZEXN,ZCP, &
ZLVOCPEXNM,ZAMOIST,ZATHETA)
CALL COMPUTE_FUNCTION_THERMO(CST%XALPI,CST%XBETAI,CST%XGAMI,CST%XLSTT,CST%XCI,ZT,ZEXN,ZCP, &
ZLSOCPEXNM,ZAMOIST_ICE,ZATHETA_ICE)
!$acc kernels present_cr( zamoist, zatheta, zlocpexnm, zlvocpexnm, zlsocpexnm, zamoist_ice, zatheta_ice )
!$mnh_expand_where(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
WHERE(PRT(:,:,2)+PRT(:,:,4)>0.0)
ZFRAC_ICE(:,:) = PRT(:,:,4) / ( PRT(:,:,2) &
+PRT(:,:,4) )
END WHERE
!$mnh_end_expand_where(JIJ=IIJB:IIJE,JK=1:IKT)
!
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZLOCPEXNM(:,:) = (1.0-ZFRAC_ICE(:,:))*ZLVOCPEXNM(:,:) &
+ZFRAC_ICE(:,:) *ZLSOCPEXNM(:,:)
ZAMOIST(:,:) = (1.0-ZFRAC_ICE(:,:))*ZAMOIST(:,:) &
+ZFRAC_ICE(:,:) *ZAMOIST_ICE(:,:)
ZATHETA(:,:) = (1.0-ZFRAC_ICE(:,:))*ZATHETA(:,:) &
+ZFRAC_ICE(:,:) *ZATHETA_ICE(:,:)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!wc call new stat functions or not

RODIER Quentin
committed
IF (TURBN%LSTATNW) THEN
CALL COMPUTE_FUNCTION_THERMO_NEW_STAT(CST%XALPW,CST%XBETAW,CST%XGAMW,CST%XLVTT,CST%XCL,ZT,ZEXN,ZCP, &
ZLOCPEXNM,ZAMOIST,ZATHETA)
CALL COMPUTE_FUNCTION_THERMO(CST%XALPW,CST%XBETAW,CST%XGAMW,CST%XLVTT,CST%XCL,ZT,ZEXN,ZCP, &
ZLOCPEXNM,ZAMOIST,ZATHETA)
ENDIF
END IF
!
!

RODIER Quentin
committed
IF ( TPFILE%LOPENED .AND. TURBN%LTURB_DIAG ) THEN
!$acc update self(ZAMOIST,ZATHETA)
TZFIELD = TFIELDMETADATA( &
CMNHNAME = 'ATHETA', &
CSTDNAME = '', &
CLONGNAME = 'ATHETA', &
CUNITS = 'm', &
CDIR = 'XY', &
CCOMMENT = 'X_Y_Z_ATHETA', &
NGRID = 1, &
NTYPE = TYPEREAL, &
NDIMS = 3, &
LTIMEDEP = .TRUE. )
CALL IO_FIELD_WRITE_PHY(D,TPFILE,TZFIELD,ZATHETA)
TZFIELD = TFIELDMETADATA( &
CMNHNAME = 'AMOIST', &
CSTDNAME = '', &
CLONGNAME = 'AMOIST', &
CUNITS = 'm', &
CDIR = 'XY', &
CCOMMENT = 'X_Y_Z_AMOIST', &
NGRID = 1, &
NTYPE = TYPEREAL, &
NDIMS = 3, &
LTIMEDEP = .TRUE. )
CALL IO_FIELD_WRITE_PHY(D,TPFILE,TZFIELD,ZAMOIST)
END IF
!
ELSE
!$acc kernels present_cr( zlocpexnm )

RODIER Quentin
committed
ZLOCPEXNM(:,:)=0.
END IF ! loop end on KRRL >= 1
!
! computes conservative variables
!
!$acc update device(PRRS,PRTHLS)
IF ( KRRL >= 1 ) THEN
!$acc kernels present_cr( zlocpexnm )
IF ( KRRI >= 1 ) THEN
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
PRT(:,:,1) = PRT(:,:,1) + PRT(:,:,2) &
+ PRT(:,:,4)
PRRS(:,:,1) = PRRS(:,:,1) + PRRS(:,:,2) &
+ PRRS(:,:,4)

RODIER Quentin
committed
PTHLT(:,:) = PTHLT(:,:) - ZLVOCPEXNM(:,:) &
* PRT(:,:,2) &
- ZLSOCPEXNM(:,:) * PRT(:,:,4)
PRTHLS(:,:) = PRTHLS(:,:) - ZLVOCPEXNM(:,:) &
* PRRS(:,:,2) &
- ZLSOCPEXNM(:,:) * PRRS(:,:,4)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
PRT(:,:,1) = PRT(:,:,1) + PRT(:,:,2)
PRRS(:,:,1) = PRRS(:,:,1) + PRRS(:,:,2)

RODIER Quentin
committed
PTHLT(:,:) = PTHLT(:,:) - ZLOCPEXNM(:,:) &
* PRT(:,:,2)
PRTHLS(:,:) = PRTHLS(:,:) - ZLOCPEXNM(:,:) &
* PRRS(:,:,2)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
END IF
!
!* stores value of conservative variables & wind before turbulence tendency (AROME diag)
IF(PRESENT(PDRUS_TURB)) THEN
PDRUS_TURB(:,:) = PRUS(:,:)
PDRVS_TURB(:,:) = PRVS(:,:)
PDRTHLS_TURB(:,:) = PRTHLS(:,:)
PDRRTS_TURB(:,:) = PRRS(:,:,1)
PDRSVS_TURB(:,:,:) = PRSVS(:,:,:)
!----------------------------------------------------------------------------
!
!* 3. MIXING LENGTH : SELECTION AND COMPUTATION
! -----------------------------------------
!
!

RODIER Quentin
committed
IF (.NOT. TURBN%LHARAT) THEN

RODIER Quentin
committed
SELECT CASE (TURBN%CTURBLEN)
!
!* 3.1 BL89 mixing length
! ------------------
CASE ('BL89')
!$acc kernels present_cr(ZSHEAR)
!$acc end kernels
CALL BL89(D,CST,CSTURB,PZZ,PDZZ,PTHVREF,ZTHLM,KRR,ZRM,PTKET,ZSHEAR,ZLM,GOCEAN,HPROGRAM)

RODIER Quentin
committed
!* 3.2 RM17 mixing length
! ------------------
CASE ('RM17')

RODIER Quentin
committed
CALL GZ_U_UW_PHY(D,PUT,PDZZ,ZWORK1)
CALL MZF_PHY(D,ZWORK1,ZWORK2)
CALL MXF_PHY(D,ZWORK2,ZDUDZ)
!
CALL GZ_V_VW_PHY(D,PVT,PDZZ,ZWORK1)
CALL MZF_PHY(D,ZWORK1,ZWORK2)
CALL MYF_PHY(D,ZWORK2,ZDVDZ)
!
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZSHEAR(:,:) = SQRT(ZDUDZ(:,:)*ZDUDZ(:,:) &
+ ZDVDZ(:,:)*ZDVDZ(:,:))
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!$acc end kernels
CALL BL89(D,CST,CSTURB,PZZ,PDZZ,PTHVREF,ZTHLM,KRR,ZRM,PTKET,ZSHEAR,ZLM,GOCEAN,HPROGRAM)

RODIER Quentin
committed
!
!* 3.3 Grey-zone combined RM17 & Deardorff mixing lengths

RODIER Quentin
committed
! --------------------------------------------------
CASE ('HM21')

RODIER Quentin
committed
CALL GZ_U_UW_PHY(D,PUT,PDZZ,ZWORK1)
CALL MZF_PHY(D,ZWORK1,ZWORK2)
CALL MXF_PHY(D,ZWORK2,ZDUDZ)
!
CALL GZ_V_VW_PHY(D,PVT,PDZZ,ZWORK1)
CALL MZF_PHY(D,ZWORK1,ZWORK2)
CALL MYF_PHY(D,ZWORK2,ZDVDZ)
!
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZSHEAR(:,:) = SQRT(ZDUDZ(:,:)*ZDUDZ(:,:) &
+ ZDVDZ(:,:)*ZDVDZ(:,:))
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!$acc end kernels
CALL BL89(D,CST,CSTURB,PZZ,PDZZ,PTHVREF,ZTHLM,KRR,ZRM,PTKET,ZSHEAR,ZLM,GOCEAN,HPROGRAM)

RODIER Quentin
committed
CALL DELT(ZLMW,ODZ=.FALSE.)
! The minimum mixing length is chosen between Horizontal grid mesh (not taking into account the vertical grid mesh) and RM17.
! For large horizontal grid meshes, this is equal to RM17
! For LES grid meshes, this is equivalent to Deardorff : the base mixing lentgh is the horizontal grid mesh,

RODIER Quentin
committed
!and it is limited by a stability-based length (RM17), as was done in Deardorff length (but taking into account shear as well)

RODIER Quentin
committed
! For grid meshes in the grey zone, then this is the smaller of the two.
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZLM(:,:) = MIN(ZLM(:,:),TURBN%XCADAP*ZLMW(:,:))
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
!
!* 3.4 Delta mixing length
! -------------------
!
CASE ('DELT')
CALL DELT(ZLM,ODZ=.TRUE.)
!* 3.5 Deardorff mixing length
! -----------------------
!
CASE ('DEAR')
CALL DEAR(ZLM)
!
!* 3.6 Blackadar mixing length
! -----------------------
!
CASE ('BLKR')
ZL0 = 100.
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZLM(:,:) = ZL0
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
ZALPHA=0.5**(-1.5)
DO JK=IKTB,IKTE

RODIER Quentin
committed
ZLM(:,JK) = ( 0.5*(PZZ(:,JK)+PZZ(:,JK+IKL)) - &
& PZZ(:,IKA+JPVEXT_TURB*IKL) ) * PDIRCOSZW(:)
ZLM(:,JK) = ZALPHA * ZLM(:,JK) * ZL0 / ( ZL0 + ZALPHA*ZLM(:,JK) )
!$mnh_end_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
ZLM(:,IKTB-1) = ZLM(:,IKTB)
ZLM(:,IKTE+1) = ZLM(:,IKTE)
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
!
!
!
END SELECT
!
!* 3.5 Mixing length modification for cloud
! -----------------------

RODIER Quentin
committed
IF (KMODEL_CL==KMI .AND. HTURBLEN_CL/='NONE') CALL CLOUD_MODIF_LM
ENDIF ! end LHARRAT
!
!* 3.6 Dissipative length
! ------------------

RODIER Quentin
committed
IF (TURBN%LHARAT) THEN
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZLEPS(:,:)=PLENGTHM(:,:)*(3.75**2.)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!$acc kernels
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZLEPS(:,:)=ZLM(:,:)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)
!$acc end kernels
ENDIF
!
!* 3.7 Correction in the Surface Boundary Layer (Redelsperger 2001)
! ----------------------------------------
!

RODIER Quentin
committed
ZLMO(:)=XUNDEF
!$mnh_end_expand_array(JIJ=IIJB:IIJE)

RODIER Quentin
committed
IF (TURBN%LRMC01) THEN

RODIER Quentin
committed
ZUSTAR(:)=(PSFU(:)**2+PSFV(:)**2)**(0.25)
!$mnh_end_expand_array(JIJ=IIJB:IIJE)
CALL LMO(D,CST,ZUSTAR,ZTHLM(:,IKB),ZRM(:,IKB,1),PSFTH,PSFRV,ZLMO)
CALL LMO(D,CST,ZUSTAR,ZTHLM(:,IKB),ZRVM,PSFTH,ZSFRV,ZLMO)

RODIER Quentin
committed
CALL RMC01(D,CST,CSTURB,TURBN%CTURBLEN,PZZ,PDXX,PDYY,PDZZ,PDIRCOSZW,PSBL_DEPTH,ZLMO,ZLM,ZLEPS)
!RMC01 is only applied on RM17 in HM21
IF (TURBN%CTURBLEN=='HM21') THEN
!$mnh_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
ZLEPS(:,:) = MIN(ZLEPS(:,:),ZLMW(:,:)*TURBN%XCADAP)
!$mnh_end_expand_array(JIJ=IIJB:IIJE,JK=1:IKT)

RODIER Quentin
committed
!* 3.8 Mixing length in external points (used if TURBN%CTURBDIM="3DIM")
! ----------------------------------------------------------
!

RODIER Quentin
committed
IF (TURBN%CTURBDIM=="3DIM") THEN
CALL UPDATE_LM(D,HLBCX,HLBCY,ZLM,ZLEPS)
!* 3.9 Mixing length correction if immersed walls
! ------------------------------------------
!
IF (OIBM) THEN
CALL IBM_MIXINGLENGTH(D,ZLM,ZLEPS,PIBM_XMUT,PIBM_LS,PTKET)
!----------------------------------------------------------------------------
!
!* 4. GO INTO THE AXES FOLLOWING THE SURFACE
! --------------------------------------
!
!
!* 4.1 rotate the wind at time t
!
!
!

RODIER Quentin
committed
IF (HPROGRAM/='AROME ') THEN
CALL ROTATE_WIND(D,PUT,PVT,PWT, &
PDIRCOSXW, PDIRCOSYW, PDIRCOSZW, &
PCOSSLOPE,PSINSLOPE, &
PDXX,PDYY,PDZZ, &
ZUSLOPE,ZVSLOPE )
CALL UPDATE_ROTATE_WIND(D,ZUSLOPE,ZVSLOPE,HLBCX,HLBCY)
!$acc kernels present_cr(ZUSLOPE,ZVSLOPE)

RODIER Quentin
committed
ZUSLOPE(:)=PUT(:,IKA)
ZVSLOPE(:)=PVT(:,IKA)
END IF
IF (GOCEAN) THEN
!$acc kernels present_cr(ZUSLOPE,ZVSLOPE)

RODIER Quentin
committed
ZUSLOPE(:)=PUT(:,IKU-1)
ZVSLOPE(:)=PVT(:,IKU-1)
END IF
!* 4.2 compute the proportionality coefficient between wind and stress

RODIER Quentin
committed
ZCDUEFF(:) =-SQRT ( (PSFU(:)**2 + PSFV(:)**2) / &
(CST%XMNH_TINY + ZUSLOPE(:)**2 + ZVSLOPE(:)**2 ) )
!
!* 4.6 compute the surface tangential fluxes
!
!$acc kernels present_cr(ZTAU22M,ZTAU33M)
IF (GOCEAN) THEN

RODIER Quentin
committed
ZTAU11M(:)=0.

RODIER Quentin
committed
ELSE

RODIER Quentin
committed
ZTAU11M(:) =2./3.*( (1.+ (PZZ(:,IKB+IKL)-PZZ(:,IKB)) &
/(PDZZ(:,IKB+IKL)+PDZZ(:,IKB)) &
) *PTKET(:,IKB) &
-0.5 *PTKET(:,IKB+IKL) &

RODIER Quentin
committed
END IF

RODIER Quentin
committed
ZTAU12M(:) =0.0
ZTAU22M(:) =ZTAU11M(:)
ZTAU33M(:) =ZTAU11M(:)
!
!* 4.7 third order terms in temperature and water fluxes and correlations
! ------------------------------------------------------------------
!
!
ZMWTH(:,:) = 0. ! w'2th'
ZMWR(:,:) = 0. ! w'2r'
ZMTH2(:,:) = 0. ! w'th'2
ZMR2(:,:) = 0. ! w'r'2
ZMTHR(:,:) = 0. ! w'th'r'

RODIER Quentin
committed
!

RODIER Quentin
committed
IF (TURBN%CTOM=='TM06') THEN
CALL TM06(D,CST,PTHVREF,PBL_DEPTH,PZZ,PSFTH,ZMWTH,ZMTH2)