source: LMDZ6/trunk/libf/phylmd/lmdz_thermcell_dtke.f90

Last change on this file was 5268, checked in by abarral, 9 days ago

.f90 <-> .F90 depending on cpp key use

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 3.0 KB
RevLine 
[4590]1MODULE lmdz_thermcell_dtke
2CONTAINS
3
[1312]4      subroutine thermcell_dtke(ngrid,nlay,nsrf,ptimestep,fm0,entr0,  &
5     &           rg,pplev,tke)
[2311]6      USE print_control_mod, ONLY: prt_level
[1312]7      implicit none
8
9!=======================================================================
10!
11!   Calcul du transport verticale dans la couche limite en presence
12!   de "thermiques" explicitement representes
13!   calcul du dq/dt une fois qu'on connait les ascendances
14!
15!=======================================================================
16
[4824]17      integer, intent(in) :: ngrid,nlay,nsrf
18      real, intent(in) :: ptimestep
19      real, dimension(ngrid,nlay), intent(in) :: entr0
20      real, dimension(ngrid,nlay+1), intent(in) :: fm0,pplev
21      real, intent(in) :: rg
22      real, intent(inout) :: tke(ngrid,nlay+1,nsrf)
[1312]23
[4824]24      real, dimension(ngrid,nlay) :: masse0,detr0, masse,entr,detr
25      real, dimension(ngrid,nlay+1) :: fm,wqd,q,qa
[1312]26      integer lev_out                           ! niveau pour les print
27
28
[4824]29      real :: zzm
[1312]30
31      integer ig,k
32      integer isrf
33
34
35      lev_out=0
36
[4824]37!print*,'thermcell_dtke'
[1312]38
39      if (prt_level.ge.1) print*,'Q2 THERMCEL_DQ 0'
40
41!   calcul du detrainement
42      do k=1,nlay
43         detr0(:,k)=fm0(:,k)-fm0(:,k+1)+entr0(:,k)
44         masse0(:,k)=(pplev(:,k)-pplev(:,k+1))/RG
45      enddo
46
47
48! Decalage vertical des entrainements et detrainements.
49      masse(:,1)=0.5*masse0(:,1)
50      entr(:,1)=0.5*entr0(:,1)
51      detr(:,1)=0.5*detr0(:,1)
52      fm(:,1)=0.
53      do k=1,nlay-1
54         masse(:,k+1)=0.5*(masse0(:,k)+masse0(:,k+1))
55         entr(:,k+1)=0.5*(entr0(:,k)+entr0(:,k+1))
56         detr(:,k+1)=0.5*(detr0(:,k)+detr0(:,k+1))
57         fm(:,k+1)=fm(:,k)+entr(:,k)-detr(:,k)
58      enddo
59      fm(:,nlay+1)=0.
60
[4746]61
62
63do isrf=1,nsrf
64
65   q(:,:)=tke(:,:,isrf)
[1312]66!   calcul de la valeur dans les ascendances
67      do ig=1,ngrid
68         qa(ig,1)=q(ig,1)
69      enddo
70
71
72    if (1==1) then
73      do k=2,nlay
74         do ig=1,ngrid
75            if ((fm(ig,k+1)+detr(ig,k))*ptimestep.gt.  &
76     &         1.e-5*masse(ig,k)) then
77         qa(ig,k)=(fm(ig,k)*qa(ig,k-1)+entr(ig,k)*q(ig,k))  &
78     &         /(fm(ig,k+1)+detr(ig,k))
79            else
80               qa(ig,k)=q(ig,k)
81            endif
82            if (qa(ig,k).lt.0.) then
83!               print*,'qa<0!!!'
84            endif
85            if (q(ig,k).lt.0.) then
86!               print*,'q<0!!!'
87            endif
88         enddo
89      enddo
90
91! Calcul du flux subsident
92
93      do k=2,nlay
94         do ig=1,ngrid
95            wqd(ig,k)=fm(ig,k)*q(ig,k)
96            if (wqd(ig,k).lt.0.) then
97!               print*,'wqd<0!!!'
98            endif
99         enddo
100      enddo
101      do ig=1,ngrid
102         wqd(ig,1)=0.
103         wqd(ig,nlay+1)=0.
104      enddo
105     
106
107! Calcul des tendances
108      do k=1,nlay
109         do ig=1,ngrid
110            q(ig,k)=q(ig,k)+(detr(ig,k)*qa(ig,k)-entr(ig,k)*q(ig,k)  &
111     &               -wqd(ig,k)+wqd(ig,k+1))  &
112     &               *ptimestep/masse(ig,k)
113         enddo
114      enddo
115
116 endif
117
118   tke(:,:,isrf)=q(:,:)
119
120enddo
121
122      return
123      end
[4590]124END MODULE lmdz_thermcell_dtke
Note: See TracBrowser for help on using the repository browser.