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

Last change on this file since 5838 was 5838, checked in by rkazeroni, 2 months ago

For GPU porting of thermcell_dtke routine:

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