lmdz_thermcell_dq.f90 Source File


This file depends on

sourcefile~~lmdz_thermcell_dq.f90~~EfferentGraph sourcefile~lmdz_thermcell_dq.f90 lmdz_thermcell_dq.f90 sourcefile~print_control_mod.f90 print_control_mod.f90 sourcefile~lmdz_thermcell_dq.f90->sourcefile~print_control_mod.f90 sourcefile~lmdz_thermcell_ini.f90 lmdz_thermcell_ini.f90 sourcefile~lmdz_thermcell_dq.f90->sourcefile~lmdz_thermcell_ini.f90 sourcefile~strings_mod.f90 strings_mod.f90 sourcefile~lmdz_thermcell_ini.f90->sourcefile~strings_mod.f90 sourcefile~ioipsl_getin_p_mod.f90 ioipsl_getin_p_mod.f90 sourcefile~lmdz_thermcell_ini.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~strings_mod.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~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~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~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~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~lmdz_mpi.f90 lmdz_mpi.F90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_mpi.f90 sourcefile~lmdz_cppkeys_wrapper.f90 lmdz_cppkeys_wrapper.F90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_cppkeys_wrapper.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~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 sourcefile~dimphy.f90 dimphy.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~dimphy.f90

Files dependent on this one

sourcefile~~lmdz_thermcell_dq.f90~~AfferentGraph sourcefile~lmdz_thermcell_dq.f90 lmdz_thermcell_dq.f90 sourcefile~lmdz_thermcell_old.f90 lmdz_thermcell_old.f90 sourcefile~lmdz_thermcell_old.f90->sourcefile~lmdz_thermcell_dq.f90 sourcefile~lmdz_thermcell_old.f90~2 lmdz_thermcell_old.f90 sourcefile~lmdz_thermcell_old.f90~2->sourcefile~lmdz_thermcell_dq.f90 sourcefile~lmdz_ratqs_multi.f90~2 lmdz_ratqs_multi.f90 sourcefile~lmdz_ratqs_multi.f90~2->sourcefile~lmdz_thermcell_dq.f90 sourcefile~lmdz_ratqs_multi.f90 lmdz_ratqs_multi.f90 sourcefile~lmdz_ratqs_multi.f90->sourcefile~lmdz_thermcell_dq.f90 sourcefile~phytracr_spl_mod.f90 phytracr_spl_mod.F90 sourcefile~phytracr_spl_mod.f90->sourcefile~lmdz_thermcell_dq.f90 sourcefile~phytracr_spl_mod.f90~2 phytracr_spl_mod.F90 sourcefile~phytracr_spl_mod.f90~2->sourcefile~lmdz_thermcell_dq.f90 sourcefile~phytrac_mod.f90 phytrac_mod.f90 sourcefile~phytrac_mod.f90->sourcefile~lmdz_thermcell_dq.f90 sourcefile~lmdz_thermcell_main.f90 lmdz_thermcell_main.F90 sourcefile~lmdz_thermcell_main.f90->sourcefile~lmdz_thermcell_dq.f90 sourcefile~lmdz_thermcell_main.f90~2 lmdz_thermcell_main.F90 sourcefile~lmdz_thermcell_main.f90~2->sourcefile~lmdz_thermcell_dq.f90 sourcefile~phytrac_mod.f90~2 phytrac_mod.f90 sourcefile~phytrac_mod.f90~2->sourcefile~lmdz_thermcell_dq.f90 sourcefile~phys_output_write_spl_mod.f90~2 phys_output_write_spl_mod.F90 sourcefile~phys_output_write_spl_mod.f90~2->sourcefile~phytracr_spl_mod.f90 sourcefile~phys_output_write_mod.f90 phys_output_write_mod.F90 sourcefile~phys_output_write_mod.f90->sourcefile~phytrac_mod.f90 sourcefile~lmdz_thermcell_alp.f90 lmdz_thermcell_alp.f90 sourcefile~lmdz_thermcell_alp.f90->sourcefile~lmdz_thermcell_main.f90 sourcefile~lmdz_thermcell_alp.f90~2 lmdz_thermcell_alp.f90 sourcefile~lmdz_thermcell_alp.f90~2->sourcefile~lmdz_thermcell_main.f90 sourcefile~phys_output_write_mod.f90~2 phys_output_write_mod.F90 sourcefile~phys_output_write_mod.f90~2->sourcefile~phytrac_mod.f90 sourcefile~physiq_mod.f90 physiq_mod.F90 sourcefile~physiq_mod.f90->sourcefile~phytracr_spl_mod.f90 sourcefile~physiq_mod.f90->sourcefile~phytrac_mod.f90 sourcefile~physiq_mod.f90->sourcefile~phys_output_write_mod.f90 sourcefile~lmdz_ratqs_main.f90 lmdz_ratqs_main.f90 sourcefile~physiq_mod.f90->sourcefile~lmdz_ratqs_main.f90 sourcefile~phys_output_write_spl_mod.f90 phys_output_write_spl_mod.F90 sourcefile~physiq_mod.f90->sourcefile~phys_output_write_spl_mod.f90 sourcefile~calltherm_mod.f90 calltherm_mod.F90 sourcefile~physiq_mod.f90->sourcefile~calltherm_mod.f90 sourcefile~diag_slp.f90 diag_slp.f90 sourcefile~physiq_mod.f90->sourcefile~diag_slp.f90 sourcefile~phys_output_mod.f90 phys_output_mod.F90 sourcefile~physiq_mod.f90->sourcefile~phys_output_mod.f90 sourcefile~lmdz_ratqs_main.f90~2 lmdz_ratqs_main.f90 sourcefile~lmdz_ratqs_main.f90~2->sourcefile~lmdz_ratqs_multi.f90 sourcefile~lmdz_ratqs_main.f90->sourcefile~lmdz_ratqs_multi.f90 sourcefile~calltherm_mod.f90~2 calltherm_mod.F90 sourcefile~calltherm_mod.f90~2->sourcefile~lmdz_thermcell_old.f90 sourcefile~calltherm_mod.f90~2->sourcefile~lmdz_thermcell_main.f90 sourcefile~calltherm_mod.f90~2->sourcefile~lmdz_thermcell_alp.f90 sourcefile~phys_output_write_spl_mod.f90->sourcefile~phytracr_spl_mod.f90 sourcefile~calltherm_mod.f90->sourcefile~lmdz_thermcell_old.f90 sourcefile~calltherm_mod.f90->sourcefile~lmdz_thermcell_main.f90 sourcefile~calltherm_mod.f90->sourcefile~lmdz_thermcell_alp.f90 sourcefile~physiq_mod.f90~2 physiq_mod.F90 sourcefile~physiq_mod.f90~2->sourcefile~phytracr_spl_mod.f90 sourcefile~physiq_mod.f90~2->sourcefile~phytrac_mod.f90 sourcefile~physiq_mod.f90~2->sourcefile~phys_output_write_mod.f90 sourcefile~physiq_mod.f90~2->sourcefile~lmdz_ratqs_main.f90 sourcefile~physiq_mod.f90~2->sourcefile~phys_output_write_spl_mod.f90 sourcefile~physiq_mod.f90~2->sourcefile~calltherm_mod.f90 sourcefile~physiq_mod.f90~2->sourcefile~diag_slp.f90 sourcefile~physiq_mod.f90~2->sourcefile~phys_output_mod.f90 sourcefile~diag_slp.f90->sourcefile~phys_output_write_mod.f90 sourcefile~phys_output_mod.f90->sourcefile~phys_output_write_mod.f90 sourcefile~old_lmdz1d.f90 old_lmdz1d.f90 sourcefile~old_lmdz1d.f90->sourcefile~physiq_mod.f90 sourcefile~diag_slp.f90~2 diag_slp.f90 sourcefile~diag_slp.f90~2->sourcefile~phys_output_write_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~scm.f90 scm.f90 sourcefile~scm.f90->sourcefile~physiq_mod.f90 sourcefile~phys_output_mod.f90~2 phys_output_mod.F90 sourcefile~phys_output_mod.f90~2->sourcefile~phys_output_write_mod.f90 sourcefile~recmwf_aero.f90 recmwf_aero.F90 sourcefile~recmwf_aero.f90->sourcefile~phys_output_mod.f90 sourcefile~recmwf_aero.f90~2 recmwf_aero.F90 sourcefile~recmwf_aero.f90~2->sourcefile~phys_output_mod.f90 sourcefile~sw_aeroar4.f90~2 sw_aeroAR4.f90 sourcefile~sw_aeroar4.f90~2->sourcefile~phys_output_mod.f90 sourcefile~calfis.f90 calfis.f90 sourcefile~calfis.f90->sourcefile~callphysiq_mod.f90 sourcefile~sw_aeroar4.f90 sw_aeroAR4.f90 sourcefile~sw_aeroar4.f90->sourcefile~phys_output_mod.f90

Contents

Source Code


Source Code

MODULE lmdz_thermcell_dq
CONTAINS

      subroutine thermcell_dq(ngrid,nlay,impl,ptimestep,fm,entr,  &
     &           masse,q,dq,qa,lev_out)
      USE print_control_mod, ONLY: prt_level

      implicit none

!=======================================================================
!
!   Calcul du transport verticale dans la couche limite en presence
!   de "thermiques" explicitement representes
!   calcul du dq/dt une fois qu'on connait les ascendances
!
! Modif 2013/01/04 (FH hourdin@lmd.jussieu.fr)
!  Introduction of an implicit computation of vertical advection in
!  the environment of thermal plumes in thermcell_dq
!  impl =     0 : explicit, 1 : implicit, -1 : old version
!
!=======================================================================

! arguments
      integer, intent(in) :: ngrid,nlay,impl
      real, intent(in) :: ptimestep
      real, intent(in), dimension(ngrid,nlay) :: masse
      real, intent(inout), dimension(ngrid,nlay) :: entr,q
      real, intent(in), dimension(ngrid,nlay+1) :: fm
      real, intent(out), dimension(ngrid,nlay) :: dq,qa
      integer, intent(in) :: lev_out                           ! niveau pour les print

! Local
      real, dimension(ngrid,nlay) :: detr,qold
      real, dimension(ngrid,nlay+1) :: wqd,fqa
      real zzm
      integer ig,k
      real cfl

      integer niter,iter
      CHARACTER (LEN=20), PARAMETER :: modname='thermcell_dq'
      CHARACTER (LEN=80) :: abort_message


! Old explicite scheme
if (impl<=-1) then

         call thermcell_dq_o(ngrid,nlay,impl,ptimestep,fm,entr,  &
     &           masse,q,dq,qa,lev_out)

else
  

! Calcul du critere CFL pour l'advection dans la subsidence
      cfl = 0.
      do k=1,nlay
         do ig=1,ngrid
            zzm=masse(ig,k)/ptimestep
            cfl=max(cfl,fm(ig,k)/zzm)
            if (entr(ig,k).gt.zzm) then
               print*,'entr*dt>m,1',k,entr(ig,k)*ptimestep,masse(ig,k)
               abort_message = 'entr dt > m, 1st'
               CALL abort_physic (modname,abort_message,1)
            endif
         enddo
      enddo

      qold=q


      if (prt_level.ge.1) print*,'Q2 THERMCEL_DQ 0'

!   calcul du detrainement
      do k=1,nlay
         do ig=1,ngrid
            detr(ig,k)=fm(ig,k)-fm(ig,k+1)+entr(ig,k)
!           print*,'Q2 DQ ',detr(ig,k),fm(ig,k),entr(ig,k)
!test
            if (detr(ig,k).lt.0.) then
               entr(ig,k)=entr(ig,k)-detr(ig,k)
               detr(ig,k)=0.
!               print*,'detr2<0!!!','ig=',ig,'k=',k,'f=',fm(ig,k),
!     s         'f+1=',fm(ig,k+1),'e=',entr(ig,k),'d=',detr(ig,k)
            endif
            if (fm(ig,k+1).lt.0.) then
!               print*,'fm2<0!!!'
            endif
            if (entr(ig,k).lt.0.) then
!               print*,'entr2<0!!!'
            endif
         enddo
      enddo

! Computation of tracer concentrations in the ascending plume
      do ig=1,ngrid
         qa(ig,1)=q(ig,1)
      enddo

      do k=2,nlay
         do ig=1,ngrid
            if ((fm(ig,k+1)+detr(ig,k))*ptimestep.gt.  &
     &         1.e-5*masse(ig,k)) then
         qa(ig,k)=(fm(ig,k)*qa(ig,k-1)+entr(ig,k)*q(ig,k))  &
     &         /(fm(ig,k+1)+detr(ig,k))
            else
               qa(ig,k)=q(ig,k)
            endif
            if (qa(ig,k).lt.0.) then
!               print*,'qa<0!!!'
            endif
            if (q(ig,k).lt.0.) then
!               print*,'q<0!!!'
            endif
         enddo
      enddo

! Plume vertical flux
      do k=2,nlay-1
         fqa(:,k)=fm(:,k)*qa(:,k-1)
      enddo
      fqa(:,1)=0. ; fqa(:,nlay)=0.


! Trace species evolution
   if (impl==0) then
      do k=1,nlay-1
         q(:,k)=q(:,k)+(fqa(:,k)-fqa(:,k+1)-fm(:,k)*q(:,k)+fm(:,k+1)*q(:,k+1)) &
     &               *ptimestep/masse(:,k)
      enddo
   else
      do k=nlay-1,1,-1
! FH debut de modif : le calcul ci dessous modifiait numériquement
! la concentration quand le flux de masse etait nul car on divisait
! puis multipliait par masse/ptimestep.
!        q(:,k)=(masse(:,k)*q(:,k)/ptimestep+fqa(:,k)-fqa(:,k+1)+fm(:,k+1)*q(:,k+1)) &
!    &               /(fm(:,k)+masse(:,k)/ptimestep)
         q(:,k)=(q(:,k)+ptimestep/masse(:,k)*(fqa(:,k)-fqa(:,k+1)+fm(:,k+1)*q(:,k+1))) &
      &               /(1.+fm(:,k)*ptimestep/masse(:,k))
! FH fin de modif.
      enddo
   endif

! Tendencies
      do k=1,nlay
         do ig=1,ngrid
            dq(ig,k)=(q(ig,k)-qold(ig,k))/ptimestep
            q(ig,k)=qold(ig,k)
         enddo
      enddo

endif ! impl=-1
RETURN
end subroutine thermcell_dq


!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Obsolete version kept for convergence with Cmip5 NPv3.1 simulations
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

      subroutine thermcell_dq_o(ngrid,nlay,impl,ptimestep,fm,entr,  &
     &           masse,q,dq,qa,lev_out)
      USE print_control_mod, ONLY: prt_level
      USE lmdz_thermcell_ini, ONLY : thermals_subsid_advect_scheme,thermals_subsid_advect_more_than_one
      implicit none

!=======================================================================
!
!   Calcul du transport verticale dans la couche limite en presence
!   de "thermiques" explicitement representes
!   calcul du dq/dt une fois qu'on connait les ascendances
!
!=======================================================================

      integer ngrid,nlay,impl

      real ptimestep
      real masse(ngrid,nlay),fm(ngrid,nlay+1)
      real entr(ngrid,nlay)
      real q(ngrid,nlay)
      real dq(ngrid,nlay)
      integer lev_out                           ! niveau pour les print

      real qa(ngrid,nlay),detr(ngrid,nlay),wqd(ngrid,nlay+1)

      real zzm

      integer ig,k
      real cfl

      real qold(ngrid,nlay)
      real ztimestep
      integer niter,iter
      CHARACTER (LEN=20), PARAMETER :: modname='thermcell_dq'
      CHARACTER (LEN=80) :: abort_message



! Calcul du critere CFL pour l'advection dans la subsidence
      cfl = 0.
      do k=1,nlay
         do ig=1,ngrid
            zzm=masse(ig,k)/ptimestep
            cfl=max(cfl,fm(ig,k)/zzm)
            if (entr(ig,k).gt.zzm) then
               print*,'entr*dt>m,2',k,entr(ig,k)*ptimestep,masse(ig,k)
               abort_message = 'entr dt > m, 2nd'
               CALL abort_physic (modname,abort_message,1)
            endif
         enddo
      enddo


!     niter=int(cfl)+1 ! pour tourner avec un CFL différent en splitant
      niter=1

      ztimestep=ptimestep/niter
      qold=q


do iter=1,niter
      if (prt_level.ge.1) print*,'Q2 THERMCEL_DQ 0'

!   calcul du detrainement
      do k=1,nlay
         do ig=1,ngrid
            detr(ig,k)=fm(ig,k)-fm(ig,k+1)+entr(ig,k)
!           print*,'Q2 DQ ',detr(ig,k),fm(ig,k),entr(ig,k)
!test
            if (detr(ig,k).lt.0.) then
               entr(ig,k)=entr(ig,k)-detr(ig,k)
               detr(ig,k)=0.
!               print*,'detr2<0!!!','ig=',ig,'k=',k,'f=',fm(ig,k),
!     s         'f+1=',fm(ig,k+1),'e=',entr(ig,k),'d=',detr(ig,k)
            endif
            if (fm(ig,k+1).lt.0.) then
!               print*,'fm2<0!!!'
            endif
            if (entr(ig,k).lt.0.) then
!               print*,'entr2<0!!!'
            endif
         enddo
      enddo

!   calcul de la valeur dans les ascendances
      do ig=1,ngrid
         qa(ig,1)=q(ig,1)
      enddo

      do k=2,nlay
         do ig=1,ngrid
            if ((fm(ig,k+1)+detr(ig,k))*ztimestep.gt.  &
     &         1.e-5*masse(ig,k)) then
         qa(ig,k)=(fm(ig,k)*qa(ig,k-1)+entr(ig,k)*q(ig,k))  &
     &         /(fm(ig,k+1)+detr(ig,k))
            else
               qa(ig,k)=q(ig,k)
            endif
            if (qa(ig,k).lt.0.) then
!               print*,'qa<0!!!'
            endif
            if (q(ig,k).lt.0.) then
!               print*,'q<0!!!'
            endif
         enddo
      enddo

! Calcul du flux subsident

      if ( thermals_subsid_advect_scheme == 'center' ) then
         do k=2,nlay
            do ig=1,ngrid
                wqd(ig,k)=fm(ig,k)*0.5*(q(ig,k-1)+q(ig,k))
            enddo
         enddo
      else ! upstream scheme (default and recomanded)
         do k=2,nlay
            do ig=1,ngrid
               zzm=masse(ig,k)/ztimestep
               if ( fm(ig,k)<=zzm .or. thermals_subsid_advect_more_than_one == 0 ) then
                   wqd(ig,k)=fm(ig,k)*q(ig,k)
               else
                  wqd(ig,k)=zzm*q(ig,k)+(fm(ig,k)-zzm)*q(ig,k+1)
               endif
            enddo
         enddo
      endif
      do ig=1,ngrid
         wqd(ig,1)=0.
         wqd(ig,nlay+1)=0.
      enddo
     

! Calcul des tendances
      do k=1,nlay
         do ig=1,ngrid
            q(ig,k)=q(ig,k)+(detr(ig,k)*qa(ig,k)-entr(ig,k)*q(ig,k)  &
     &               -wqd(ig,k)+wqd(ig,k+1))  &
     &               *ztimestep/masse(ig,k)
!            if (dq(ig,k).lt.0.) then
!               print*,'dq<0!!!'
!            endif
         enddo
      enddo


enddo


! Calcul des tendances
      do k=1,nlay
         do ig=1,ngrid
            dq(ig,k)=(q(ig,k)-qold(ig,k))/ptimestep
            q(ig,k)=qold(ig,k)
         enddo
      enddo

      return
    end subroutine thermcell_dq_o
END MODULE lmdz_thermcell_dq