lmdz_thermcell_alp.f90 Source File


This file depends on

sourcefile~~lmdz_thermcell_alp.f90~~EfferentGraph sourcefile~lmdz_thermcell_alp.f90 lmdz_thermcell_alp.f90 sourcefile~indice_sol_mod.f90 indice_sol_mod.f90 sourcefile~lmdz_thermcell_alp.f90->sourcefile~indice_sol_mod.f90 sourcefile~yomcst_mod_h.f90 yomcst_mod_h.f90 sourcefile~lmdz_thermcell_alp.f90->sourcefile~yomcst_mod_h.f90 sourcefile~lmdz_thermcell_main.f90 lmdz_thermcell_main.F90 sourcefile~lmdz_thermcell_alp.f90->sourcefile~lmdz_thermcell_main.f90 sourcefile~yoethf_mod_h.f90 yoethf_mod_h.f90 sourcefile~lmdz_thermcell_alp.f90->sourcefile~yoethf_mod_h.f90 sourcefile~alpale_mod.f90 alpale_mod.f90 sourcefile~lmdz_thermcell_alp.f90->sourcefile~alpale_mod.f90 sourcefile~lmdz_thermcell_down.f90 lmdz_thermcell_down.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_down.f90 sourcefile~lmdz_thermcell_height.f90 lmdz_thermcell_height.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_height.f90 sourcefile~lmdz_thermcell_plume.f90 lmdz_thermcell_plume.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_plume.f90 sourcefile~lmdz_thermcell_flux2.f90 lmdz_thermcell_flux2.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_flux2.f90 sourcefile~lmdz_thermcell_ini.f90 lmdz_thermcell_ini.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_ini.f90 sourcefile~lmdz_thermcell_dq.f90 lmdz_thermcell_dq.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_dq.f90 sourcefile~lmdz_thermcell_dv2.f90 lmdz_thermcell_dv2.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_dv2.f90 sourcefile~lmdz_thermcell_qsat.f90 lmdz_thermcell_qsat.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_qsat.f90 sourcefile~lmdz_thermcell_dry.f90 lmdz_thermcell_dry.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_dry.f90 sourcefile~lmdz_thermcell_plume_6a.f90 lmdz_thermcell_plume_6A.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_plume_6a.f90 sourcefile~lmdz_thermcell_env.f90 lmdz_thermcell_env.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_env.f90 sourcefile~lmdz_thermcell_closure.f90 lmdz_thermcell_closure.f90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_closure.f90 sourcefile~alpale_mod.f90->sourcefile~yomcst_mod_h.f90 sourcefile~alpale_mod.f90->sourcefile~yoethf_mod_h.f90 sourcefile~dimphy.f90 dimphy.f90 sourcefile~alpale_mod.f90->sourcefile~dimphy.f90 sourcefile~print_control_mod.f90 print_control_mod.f90 sourcefile~alpale_mod.f90->sourcefile~print_control_mod.f90 sourcefile~ioipsl_getin_p_mod.f90 ioipsl_getin_p_mod.f90 sourcefile~alpale_mod.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~phys_local_var_mod.f90 phys_local_var_mod.F90 sourcefile~alpale_mod.f90->sourcefile~phys_local_var_mod.f90 sourcefile~lmdz_thermcell_down.f90->sourcefile~lmdz_thermcell_ini.f90 sourcefile~lmdz_thermcell_plume.f90->sourcefile~lmdz_thermcell_ini.f90 sourcefile~lmdz_thermcell_plume.f90->sourcefile~lmdz_thermcell_qsat.f90 sourcefile~lmdz_thermcell_alim.f90 lmdz_thermcell_alim.f90 sourcefile~lmdz_thermcell_plume.f90->sourcefile~lmdz_thermcell_alim.f90 sourcefile~lmdz_thermcell_flux2.f90->sourcefile~lmdz_thermcell_ini.f90 sourcefile~lmdz_thermcell_ini.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~strings_mod.f90 strings_mod.f90 sourcefile~lmdz_thermcell_ini.f90->sourcefile~strings_mod.f90 sourcefile~lmdz_thermcell_dq.f90->sourcefile~lmdz_thermcell_ini.f90 sourcefile~lmdz_thermcell_dq.f90->sourcefile~print_control_mod.f90 sourcefile~lmdz_thermcell_dv2.f90->sourcefile~print_control_mod.f90 sourcefile~lmdz_thermcell_qsat.f90->sourcefile~yomcst_mod_h.f90 sourcefile~lmdz_thermcell_qsat.f90->sourcefile~yoethf_mod_h.f90 sourcefile~lmdz_thermcell_dry.f90->sourcefile~lmdz_thermcell_ini.f90 sourcefile~lmdz_thermcell_plume_6a.f90->sourcefile~lmdz_thermcell_ini.f90 sourcefile~lmdz_thermcell_plume_6a.f90->sourcefile~lmdz_thermcell_qsat.f90 sourcefile~lmdz_thermcell_plume_6a.f90->sourcefile~lmdz_thermcell_alim.f90 sourcefile~mod_phys_lmdz_para.f90 mod_phys_lmdz_para.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~strings_mod.f90 sourcefile~mod_phys_lmdz_transfert_para.f90 mod_phys_lmdz_transfert_para.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~phys_local_var_mod.f90->sourcefile~indice_sol_mod.f90 sourcefile~phys_local_var_mod.f90->sourcefile~dimphy.f90 sourcefile~phys_output_var_mod.f90 phys_output_var_mod.f90 sourcefile~phys_local_var_mod.f90->sourcefile~phys_output_var_mod.f90 sourcefile~phys_state_var_mod.f90 phys_state_var_mod.F90 sourcefile~phys_local_var_mod.f90->sourcefile~phys_state_var_mod.f90 sourcefile~lmdz_cppkeys_wrapper.f90 lmdz_cppkeys_wrapper.F90 sourcefile~phys_local_var_mod.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~aero_mod.f90 aero_mod.f90 sourcefile~phys_local_var_mod.f90->sourcefile~aero_mod.f90 sourcefile~infotrac_phy.f90 infotrac_phy.F90 sourcefile~phys_local_var_mod.f90->sourcefile~infotrac_phy.f90 sourcefile~lmdz_thermcell_env.f90->sourcefile~lmdz_thermcell_ini.f90 sourcefile~lmdz_thermcell_env.f90->sourcefile~lmdz_thermcell_qsat.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~print_control_mod.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~mod_phys_lmdz_mpi_data.f90 mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_grid_phy_lmdz.f90 mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_omp_data.f90 mod_phys_lmdz_omp_data.F90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_omp_data.f90 sourcefile~phys_output_var_mod.f90->sourcefile~dimphy.f90 sourcefile~phys_output_var_mod.f90->sourcefile~strings_mod.f90 sourcefile~clesphys_mod_h.f90 clesphys_mod_h.f90 sourcefile~phys_output_var_mod.f90->sourcefile~clesphys_mod_h.f90 sourcefile~config_ocean_skin_m.f90 config_ocean_skin_m.F90 sourcefile~phys_output_var_mod.f90->sourcefile~config_ocean_skin_m.f90 sourcefile~phys_state_var_mod.f90->sourcefile~indice_sol_mod.f90 sourcefile~phys_state_var_mod.f90->sourcefile~dimphy.f90 sourcefile~phys_state_var_mod.f90->sourcefile~aero_mod.f90 sourcefile~phys_state_var_mod.f90->sourcefile~infotrac_phy.f90 sourcefile~surface_data.f90 surface_data.f90 sourcefile~phys_state_var_mod.f90->sourcefile~surface_data.f90 sourcefile~phys_state_var_mod.f90->sourcefile~clesphys_mod_h.f90 sourcefile~phys_state_var_mod.f90->sourcefile~config_ocean_skin_m.f90 sourcefile~dimsoil_mod_h.f90 dimsoil_mod_h.f90 sourcefile~phys_state_var_mod.f90->sourcefile~dimsoil_mod_h.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90 mod_phys_lmdz_omp_transfert.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_omp_transfert.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90 mod_phys_lmdz_mpi_transfert.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_mpi_transfert.f90 sourcefile~infotrac_phy.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~infotrac_phy.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~infotrac_phy.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~infotrac_phy.f90->sourcefile~strings_mod.f90 sourcefile~iniprint_mod_h.f90 iniprint_mod_h.f90 sourcefile~infotrac_phy.f90->sourcefile~iniprint_mod_h.f90 sourcefile~readtracfiles_mod.f90 readTracFiles_mod.f90 sourcefile~infotrac_phy.f90->sourcefile~readtracfiles_mod.f90 sourcefile~lmdz_reprobus_wrappers.f90 lmdz_reprobus_wrappers.F90 sourcefile~infotrac_phy.f90->sourcefile~lmdz_reprobus_wrappers.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90->sourcefile~mod_phys_lmdz_omp_data.f90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~print_control_mod.f90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~lmdz_mpi.f90 lmdz_mpi.F90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_mpi.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~lmdz_mpi.f90 sourcefile~readtracfiles_mod.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~readtracfiles_mod.f90->sourcefile~strings_mod.f90 sourcefile~lmdz_reprobus_wrappers.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~dimphy.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~print_control_mod.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~mod_phys_lmdz_mpi_data.f90

Files dependent on this one

sourcefile~~lmdz_thermcell_alp.f90~~AfferentGraph sourcefile~lmdz_thermcell_alp.f90 lmdz_thermcell_alp.f90 sourcefile~calltherm_mod.f90 calltherm_mod.F90 sourcefile~calltherm_mod.f90->sourcefile~lmdz_thermcell_alp.f90 sourcefile~calltherm_mod.f90~2 calltherm_mod.F90 sourcefile~calltherm_mod.f90~2->sourcefile~lmdz_thermcell_alp.f90 sourcefile~physiq_mod.f90 physiq_mod.F90 sourcefile~physiq_mod.f90->sourcefile~calltherm_mod.f90 sourcefile~physiq_mod.f90~2 physiq_mod.F90 sourcefile~physiq_mod.f90~2->sourcefile~calltherm_mod.f90 sourcefile~old_lmdz1d.f90 old_lmdz1d.f90 sourcefile~old_lmdz1d.f90->sourcefile~physiq_mod.f90 sourcefile~scm.f90 scm.f90 sourcefile~scm.f90->sourcefile~physiq_mod.f90 sourcefile~callphysiq_mod.f90 callphysiq_mod.f90 sourcefile~callphysiq_mod.f90->sourcefile~physiq_mod.f90 sourcefile~callphysiq_mod.f90~2 callphysiq_mod.f90 sourcefile~callphysiq_mod.f90~2->sourcefile~physiq_mod.f90 sourcefile~calfis.f90 calfis.f90 sourcefile~calfis.f90->sourcefile~callphysiq_mod.f90

Contents


Source Code

MODULE lmdz_thermcell_alp
! $Id: thermcell_main.F90 2351 2015-08-25 15:14:59Z emillour $
!
CONTAINS

      SUBROUTINE thermcell_alp(ngrid,nlay,ptimestep  &                         ! in
     &                  ,pplay,pplev  &                                        ! in
     &                  ,fm0,entr0,lmax  &                                     ! in
     &                  ,pbl_tke,pctsrf,omega,airephy &                        ! in
     &                  ,zw2,fraca &                                           ! in
     &                  ,pcon,rhobarz,wth3,wmax_sec,lalim,fm,alim_star,zmax &  ! in
!
     &                  ,zcong,ale_bl,alp_bl,lalim_conv,wght_th &                    ! out
     &                  ,zlcl,fraca0,w0,w_conv,therm_tke_max0,env_tke_max0 &   ! out
     &                  ,n2,s2,strig,ale_bl_stat &                                   ! out
     &                  ,therm_tke_max,env_tke_max &                           ! out
     &                  ,alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke &          ! out
     &                  ,alp_bl_conv,alp_bl_stat &                             ! out
     &)

USE yoethf_mod_h
      USE alpale_mod
            USE yomcst_mod_h
      USE indice_sol_mod
      USE lmdz_thermcell_main, ONLY : thermcell_tke_transport
      IMPLICIT NONE

!=======================================================================
!
!   Auteurs: Catherine Rio
!   Modifications :
!   Nicolas Rochetin et Jean-Yves Grandpeix
!         pour la fermeture stochastique. 2012
!   Fr�d�ric Hourdin :
!         netoyage informatique. 2022
!
!=======================================================================
!-----------------------------------------------------------------------
!   declarations:
!   -------------

      INCLUDE "FCTTRE.h"

!   arguments:
!   ----------

!------Entrees
      integer, intent(in) :: ngrid,nlay
      real, intent(in) :: ptimestep
      real, intent(in) :: pplay(ngrid,nlay),pplev(ngrid,nlay+1)
      integer, intent(in), dimension(ngrid) ::lmax,lalim
      real, intent(in), dimension(ngrid) :: zmax, zcong
      real, intent(in), dimension(ngrid,nlay+1) :: zw2
      real, intent(in), dimension(ngrid,nlay+1) :: fraca
      real, intent(in), dimension(ngrid,nlay) :: wth3
      real, intent(in), dimension(ngrid,nlay) :: rhobarz
      real, intent(in), dimension(ngrid) :: wmax_sec
      real, intent(in), dimension(ngrid,nlay) :: entr0
      real, intent(in), dimension(ngrid,nlay+1) :: fm0,fm
      real, intent(in), dimension(ngrid) :: pcon
      real, intent(in), dimension(ngrid,nlay) :: alim_star
      real, intent(in), dimension(ngrid,nlay+1,nbsrf) :: pbl_tke
      real, intent(in), dimension(ngrid,nbsrf) :: pctsrf
      real, intent(in), dimension(ngrid,nlay) :: omega
      real, intent(in), dimension(ngrid) :: airephy
!------Sorties
      real, intent(out), dimension(ngrid) :: ale_bl,alp_bl
      real, intent(out), dimension(ngrid,nlay) :: wght_th
      integer, intent(out), dimension(ngrid) :: lalim_conv
      real, intent(out), dimension(ngrid) :: zlcl,fraca0,w0,w_conv
      real, intent(out), dimension(ngrid) :: therm_tke_max0,env_tke_max0,n2,s2,ale_bl_stat,strig
      real, intent(out), dimension(ngrid,nlay) :: therm_tke_max,env_tke_max
      real, intent(out), dimension(ngrid) :: alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke
      real, intent(out), dimension(ngrid) :: alp_bl_conv,alp_bl_stat

!=============================================================================================
!------Local
!=============================================================================================

      REAL susqr2pi, reuler
      INTEGER ig,k,l
      integer nsrf
      real rhobarz0(ngrid)                    ! Densit� au LCL
      logical ok_lcl(ngrid)                   ! Existence du LCL des thermiques
      integer klcl(ngrid)                     ! Niveau du LCL
      real interp(ngrid)                      ! Coef d'interpolation pour le LCL
!--Triggering
      real, parameter :: su_cst=4e4              ! Surface unite: celle d'un updraft �l�mentaire
      real, parameter :: hcoef=1             ! Coefficient directeur pour le calcul de s2
      real, parameter :: hmincoef=0.3        ! Coefficient directeur pour l'ordonn�e � l'origine pour le calcul de s2 
      real, parameter :: eps1=0.3            ! Fraction de surface occup�e par la population 1 : eps1=n1*s1/(fraca0*Sd)
      real, dimension(ngrid) :: hmin         ! Ordonn�e � l'origine pour le calcul de s2
      real, dimension(ngrid) :: zmax_moy     ! Hauteur moyenne des thermiques : zmax_moy = zlcl + 0.33 (zmax-zlcl)
      real, parameter :: zmax_moy_coef=0.33
      real, dimension(ngrid) :: depth        ! Epaisseur moyenne du cumulus
      real, dimension(ngrid) :: zcong_moy 
      real, dimension(ngrid) ::  w_max                 ! Vitesse max statistique 
      real, dimension(ngrid) ::  s_max(ngrid)
!--Closure
      real, dimension(ngrid,nlay) :: pbl_tke_max       ! Profil de TKE moyenne 
      real, dimension(ngrid) :: pbl_tke_max0           ! TKE moyenne au LCL
      real, dimension(ngrid,nlay) :: w_ls              ! Vitesse verticale grande �chelle (m/s)
      real, parameter :: coef_m=1.            ! On consid�re un rendement pour alp_bl_fluct_m
      real, parameter :: coef_tke=1.          ! On consid�re un rendement pour alp_bl_fluct_tke
      real :: zdp
      real, dimension(ngrid) :: alp_int,dp_int
      real, dimension(ngrid) :: fm_tot

!------------------------------------------------------------
!  Initialize output arrays related to stochastic triggering
!------------------------------------------------------------
  DO ig = 1,ngrid
     zlcl(ig) = 0.
     fraca0(ig) = 0.
     w0(ig) = 0.
     w_conv(ig) = 0.
     therm_tke_max0(ig) = 0.
     env_tke_max0(ig) = 0.
     n2(ig) = 0.
     s2(ig) = 0.
     ale_bl_stat(ig) = 0.
     strig(ig) = 0.
     alp_bl_det(ig) = 0.
     alp_bl_fluct_m(ig) = 0.
     alp_bl_fluct_tke(ig) = 0.
     alp_bl_conv(ig) = 0.
     alp_bl_stat(ig) = 0.
  ENDDO
  DO l = 1,nlay
    DO ig = 1,ngrid
     therm_tke_max(ig,l) = 0.
     env_tke_max(ig,l) = 0.
    ENDDO
  ENDDO

!------------Test sur le LCL des thermiques
    do ig=1,ngrid
      ok_lcl(ig)=.false.
      if ( (pcon(ig) .gt. pplay(ig,nlay-1)) .and. (pcon(ig) .lt. pplay(ig,1)) ) ok_lcl(ig)=.true.
    enddo

!------------Localisation des niveaux entourant le LCL et du coef d'interpolation
    do l=1,nlay-1
      do ig=1,ngrid
        if (ok_lcl(ig)) then 
!ATTENTION,zw2 calcule en pplev
!          if ((pplay(ig,l) .ge. pcon(ig)) .and. (pplay(ig,l+1) .le. pcon(ig))) then
!          klcl(ig)=l
!          interp(ig)=(pcon(ig)-pplay(ig,klcl(ig)))/(pplay(ig,klcl(ig)+1)-pplay(ig,klcl(ig)))
!          endif
          if ((pplev(ig,l) .ge. pcon(ig)) .and. (pplev(ig,l+1) .le. pcon(ig))) then
          klcl(ig)=l
          interp(ig)=(pcon(ig)-pplev(ig,klcl(ig)))/(pplev(ig,klcl(ig)+1)-pplev(ig,klcl(ig)))
          endif
        endif
      enddo
    enddo

    do ig =1,ngrid
!CR:REHABILITATION ZMAX CONTINU
     if (ok_lcl(ig)) then 
      rhobarz0(ig)=rhobarz(ig,klcl(ig))+(rhobarz(ig,klcl(ig)+1) &
 &               -rhobarz(ig,klcl(ig)))*interp(ig)
      zlcl(ig)=(pplev(ig,1)-pcon(ig))/(rhobarz0(ig)*RG)
      zlcl(ig)=min(zlcl(ig),zmax(ig))   ! Si zlcl > zmax alors on pose zlcl = zmax
     else
      rhobarz0(ig)=0.
      zlcl(ig)=zmax(ig)
     endif
    enddo
!!jyg fin

!------------Calcul des propri�t�s du thermique au LCL 
  IF ( (iflag_trig_bl.ge.1) .or. (iflag_clos_bl.ge.1) ) THEN 

  !-----Initialisation de la TKE moyenne 
   do l=1,nlay
    do ig=1,ngrid
     pbl_tke_max(ig,l)=0.
    enddo
   enddo

!-----Calcul de la TKE moyenne 
   do nsrf=1,nbsrf
    do l=1,nlay
     do ig=1,ngrid
     pbl_tke_max(ig,l)=pctsrf(ig,nsrf)*pbl_tke(ig,l,nsrf)+pbl_tke_max(ig,l)
     enddo
    enddo
   enddo

!-----Initialisations des TKE dans et hors des thermiques 
   do l=1,nlay
    do ig=1,ngrid
    therm_tke_max(ig,l)=pbl_tke_max(ig,l)
    env_tke_max(ig,l)=pbl_tke_max(ig,l)
    enddo
   enddo

!-----Calcul de la TKE transport�e par les thermiques : therm_tke_max
   call thermcell_tke_transport(ngrid,nlay,ptimestep,fm0,entr0,  &  ! in
  &           rg,pplev,therm_tke_max)                               ! out
!   print *,' thermcell_tke_transport -> '   !!jyg

!-----Calcul des profils verticaux de TKE hors thermiques : env_tke_max, et de la vitesse verticale grande �chelle : W_ls
   do l=1,nlay
    do ig=1,ngrid
     pbl_tke_max(ig,l)=fraca(ig,l)*therm_tke_max(ig,l)+(1.-fraca(ig,l))*env_tke_max(ig,l)         !  Recalcul de TKE moyenne apr�s transport de TKE_TH
     env_tke_max(ig,l)=(pbl_tke_max(ig,l)-fraca(ig,l)*therm_tke_max(ig,l))/(1.-fraca(ig,l))       !  Recalcul de TKE dans  l'environnement apr�s transport de TKE_TH
     w_ls(ig,l)=-1.*omega(ig,l)/(RG*rhobarz(ig,l))                                                !  Vitesse verticale de grande �chelle
    enddo
   enddo
!    print *,' apres w_ls = '   !!jyg

  do ig=1,ngrid
   if (ok_lcl(ig)) then
     fraca0(ig)=fraca(ig,klcl(ig))+(fraca(ig,klcl(ig)+1) &
 &             -fraca(ig,klcl(ig)))*interp(ig)
     w0(ig)=zw2(ig,klcl(ig))+(zw2(ig,klcl(ig)+1) &
 &         -zw2(ig,klcl(ig)))*interp(ig)
     w_conv(ig)=w_ls(ig,klcl(ig))+(w_ls(ig,klcl(ig)+1) &
 &             -w_ls(ig,klcl(ig)))*interp(ig)
     therm_tke_max0(ig)=therm_tke_max(ig,klcl(ig)) &
 &                     +(therm_tke_max(ig,klcl(ig)+1)-therm_tke_max(ig,klcl(ig)))*interp(ig)
     env_tke_max0(ig)=env_tke_max(ig,klcl(ig))+(env_tke_max(ig,klcl(ig)+1) &
 &                   -env_tke_max(ig,klcl(ig)))*interp(ig)
     pbl_tke_max0(ig)=pbl_tke_max(ig,klcl(ig))+(pbl_tke_max(ig,klcl(ig)+1) &
 &                   -pbl_tke_max(ig,klcl(ig)))*interp(ig)
     if (therm_tke_max0(ig).ge.20.) therm_tke_max0(ig)=20.
     if (env_tke_max0(ig).ge.20.) env_tke_max0(ig)=20.
     if (pbl_tke_max0(ig).ge.20.) pbl_tke_max0(ig)=20.
   else 
     fraca0(ig)=0.
     w0(ig)=0.
!!jyg le 27/04/2012
!!     zlcl(ig)=0.
!!
   endif
  enddo

  ENDIF ! IF ( (iflag_trig_bl.ge.1) .or. (iflag_clos_bl.ge.1) )
!  print *,'ENDIF  ( (iflag_trig_bl.ge.1) .or. (iflag_clos_bl.ge.1) ) '    !!jyg

!------------Triggering------------------
  IF (iflag_trig_bl.ge.1) THEN 

!-----Initialisations
   depth(:)=0.
   n2(:)=0.
   s2(:)=100. ! some low value, arbitrary
   s_max(:)=0.



!-----Epaisseur du nuage (depth) et d�termination de la queue du spectre de panaches (n2,s2) et du panache le plus gros (s_max)
   do ig=1,ngrid
     zmax_moy(ig)=zlcl(ig)+zmax_moy_coef*(zmax(ig)-zlcl(ig))
     depth(ig)=zmax_moy(ig)-zlcl(ig)
     hmin(ig)=hmincoef*zlcl(ig)
     if (depth(ig).ge.10.) then 
       s2(ig)=(hcoef*depth(ig)+hmin(ig))**2
       n2(ig)=(1.-eps1)*fraca0(ig)*airephy(ig)/s2(ig)
!!
!!jyg le 27/04/2012
!!       s_max(ig)=s2(ig)*log(n2(ig))
!!       if (n2(ig) .lt. 1) s_max(ig)=0.
       s_max(ig)=s2(ig)*log(max(n2(ig),1.))
!!fin jyg
     else
       n2(ig)=0.
       s_max(ig)=0.
     endif
   enddo
!   print *,'avant Calcul de Wmax '    !!jyg

!CR: calcul de strig
   if (iflag_strig.eq.0) then
      strig(:)=s_trig
   else if (iflag_strig.eq.1) then
      do ig=1,ngrid
!         zcong_moy(ig)=zlcl(ig)+zmax_moy_coef*(zcong(ig)-zlcl(ig))
!         strig(ig)=(hcoef*(zcong_moy(ig)-zlcl(ig))+hmin(ig))**2  
         strig(ig)=(zcong(ig)-zlcl(ig))**2  
      enddo
   else if (iflag_strig.eq.2) then
      do ig=1,ngrid 
         if (h_trig.gt.zlcl(ig)) then
         strig(ig)=(h_trig-zlcl(ig))**2
         else
         strig(ig)=s_trig
         endif
      enddo   
   endif

   susqr2pi=su_cst*sqrt(2.*Rpi)
   reuler=exp(1.)
   do ig=1,ngrid
     if ( (depth(ig).ge.10.) .and. (s_max(ig).gt.susqr2pi*reuler) ) then
      w_max(ig)=w0(ig)*(1.+sqrt(2.*log(s_max(ig)/susqr2pi)-log(2.*log(s_max(ig)/susqr2pi))))
      ale_bl_stat(ig)=0.5*w_max(ig)**2
     else
      w_max(ig)=0.
      ale_bl_stat(ig)=0.
     endif
   enddo

  ENDIF ! iflag_trig_bl
!  print *,'ENDIF  iflag_trig_bl'    !!jyg

!------------Closure------------------

  IF (iflag_clos_bl.ge.2) THEN 

!-----Calcul de ALP_BL_STAT
  do ig=1,ngrid
  alp_bl_det(ig)=0.5*coef_m*rhobarz0(ig)*(w0(ig)**3)*fraca0(ig)*(1.-2.*fraca0(ig))/((1.-fraca0(ig))**2)
  alp_bl_fluct_m(ig)=1.5*rhobarz0(ig)*fraca0(ig)*(w_conv(ig)+coef_m*w0(ig))* &
 &                   (w0(ig)**2)
  alp_bl_fluct_tke(ig)=3.*coef_m*rhobarz0(ig)*w0(ig)*fraca0(ig)*(therm_tke_max0(ig)-env_tke_max0(ig)) &
 &                    +3.*rhobarz0(ig)*w_conv(ig)*pbl_tke_max0(ig)
    if (iflag_clos_bl.ge.2) then 
    alp_bl_conv(ig)=1.5*coef_m*rhobarz0(ig)*fraca0(ig)*(fraca0(ig)/(1.-fraca0(ig)))*w_conv(ig)* &
 &                   (w0(ig)**2)
    else
    alp_bl_conv(ig)=0.
    endif
  alp_bl_stat(ig)=alp_bl_det(ig)+alp_bl_fluct_m(ig)+alp_bl_fluct_tke(ig)+alp_bl_conv(ig)
  enddo

!-----S�curit� ALP infinie
  do ig=1,ngrid
   if (fraca0(ig).gt.0.98) alp_bl_stat(ig)=2.
  enddo

  ENDIF ! (iflag_clos_bl.ge.2)

!!! fin nrlmd le 10/04/2012

!      print*,'avant calcul ale et alp' 
!calcul de ALE et ALP pour la convection
      alp_bl(:)=0.
      ale_bl(:)=0.
!          print*,'ALE,ALP ,l,zw2(ig,l),ale_bl(ig),alp_bl(ig)'
      do l=1,nlay
      do ig=1,ngrid
           alp_bl(ig)=max(alp_bl(ig),0.5*rhobarz(ig,l)*wth3(ig,l) )
           ale_bl(ig)=max(ale_bl(ig),0.5*zw2(ig,l)**2)
!          print*,'ALE,ALP',l,zw2(ig,l),ale_bl(ig),alp_bl(ig)
      enddo
      enddo

! ale sec (max de wmax/2 sous la zone d'inhibition) dans
! le cas iflag_trig_bl=3
      IF (iflag_trig_bl==3) ale_bl(:)=0.5*wmax_sec(:)**2

!test:calcul de la ponderation des couches pour KE
!initialisations

      fm_tot(:)=0.
      wght_th(:,:)=1.
      lalim_conv(:)=lalim(:)

      do k=1,nlay
         do ig=1,ngrid
            if (k<=lalim_conv(ig)) fm_tot(ig)=fm_tot(ig)+fm(ig,k)
         enddo
      enddo

! assez bizarre car, si on est dans la couche d'alim et que alim_star et
! plus petit que 1.e-10, on prend wght_th=1.
      do k=1,nlay
         do ig=1,ngrid
            if (k<=lalim_conv(ig).and.alim_star(ig,k)>1.e-10) then
               wght_th(ig,k)=alim_star(ig,k)
            endif
         enddo
      enddo

!      print*,'apres wght_th'
!test pour prolonger la convection
      do ig=1,ngrid
!v1d  if ((alim_star(ig,1).lt.1.e-10).and.(therm)) then
      if ((alim_star(ig,1).lt.1.e-10)) then
      lalim_conv(ig)=1
      wght_th(ig,1)=1.
!      print*,'lalim_conv ok',lalim_conv(ig),wght_th(ig,1)
      endif
      enddo

!------------------------------------------------------------------------
! Modif CR/FH 20110310 : alp integree sur la verticale.
! Integrale verticale de ALP.
! wth3 etant aux niveaux inter-couches, on utilise d play comme masse des
! couches
!------------------------------------------------------------------------

      alp_int(:)=0.
      dp_int(:)=0.
      do l=2,nlay
        do ig=1,ngrid
           if(l.LE.lmax(ig)) THEN
           zdp=pplay(ig,l-1)-pplay(ig,l)
           alp_int(ig)=alp_int(ig)+0.5*rhobarz(ig,l)*wth3(ig,l)*zdp
           dp_int(ig)=dp_int(ig)+zdp
           endif
        enddo
      enddo

      if (iflag_coupl>=3 .and. iflag_coupl<=5) then
      do ig=1,ngrid
!valeur integree de alp_bl * 0.5:
        if (dp_int(ig)>0.) then
        alp_bl(ig)=alp_int(ig)/dp_int(ig)
        endif
      enddo!
      endif


! Facteur multiplicatif sur alp_bl
      alp_bl(:)=alp_bl_k*alp_bl(:)

!------------------------------------------------------------------------



      return
    END SUBROUTINE thermcell_alp
END MODULE lmdz_thermcell_alp