1 | MODULE lmdz_thermcell_dtke |
---|
2 | CONTAINS |
---|
3 | |
---|
4 | subroutine thermcell_dtke(ngrid,nlay,nsrf,ptimestep,fm0,entr0, & |
---|
5 | & rg,pplev,tke) |
---|
6 | USE print_control_mod, ONLY: prt_level |
---|
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 | |
---|
17 | integer ngrid,nlay,nsrf |
---|
18 | |
---|
19 | real ptimestep |
---|
20 | real masse0(ngrid,nlay),fm0(ngrid,nlay+1),pplev(ngrid,nlay+1) |
---|
21 | real entr0(ngrid,nlay),rg |
---|
22 | real tke(ngrid,nlay,nsrf) |
---|
23 | real detr0(ngrid,nlay) |
---|
24 | |
---|
25 | |
---|
26 | real masse(ngrid,nlay),fm(ngrid,nlay+1) |
---|
27 | real entr(ngrid,nlay) |
---|
28 | real q(ngrid,nlay) |
---|
29 | integer lev_out ! niveau pour les print |
---|
30 | |
---|
31 | real qa(ngrid,nlay),detr(ngrid,nlay),wqd(ngrid,nlay+1) |
---|
32 | |
---|
33 | real zzm |
---|
34 | |
---|
35 | integer ig,k |
---|
36 | integer isrf |
---|
37 | |
---|
38 | |
---|
39 | lev_out=0 |
---|
40 | |
---|
41 | |
---|
42 | if (prt_level.ge.1) print*,'Q2 THERMCEL_DQ 0' |
---|
43 | |
---|
44 | ! calcul du detrainement |
---|
45 | do k=1,nlay |
---|
46 | detr0(:,k)=fm0(:,k)-fm0(:,k+1)+entr0(:,k) |
---|
47 | masse0(:,k)=(pplev(:,k)-pplev(:,k+1))/RG |
---|
48 | enddo |
---|
49 | |
---|
50 | |
---|
51 | ! Decalage vertical des entrainements et detrainements. |
---|
52 | masse(:,1)=0.5*masse0(:,1) |
---|
53 | entr(:,1)=0.5*entr0(:,1) |
---|
54 | detr(:,1)=0.5*detr0(:,1) |
---|
55 | fm(:,1)=0. |
---|
56 | do k=1,nlay-1 |
---|
57 | masse(:,k+1)=0.5*(masse0(:,k)+masse0(:,k+1)) |
---|
58 | entr(:,k+1)=0.5*(entr0(:,k)+entr0(:,k+1)) |
---|
59 | detr(:,k+1)=0.5*(detr0(:,k)+detr0(:,k+1)) |
---|
60 | fm(:,k+1)=fm(:,k)+entr(:,k)-detr(:,k) |
---|
61 | enddo |
---|
62 | fm(:,nlay+1)=0. |
---|
63 | |
---|
64 | |
---|
65 | |
---|
66 | do isrf=1,nsrf |
---|
67 | |
---|
68 | q(:,:)=tke(:,:,isrf) |
---|
69 | ! calcul de la valeur dans les ascendances |
---|
70 | do ig=1,ngrid |
---|
71 | qa(ig,1)=q(ig,1) |
---|
72 | enddo |
---|
73 | |
---|
74 | |
---|
75 | if (1==1) then |
---|
76 | do k=2,nlay |
---|
77 | do ig=1,ngrid |
---|
78 | if ((fm(ig,k+1)+detr(ig,k))*ptimestep.gt. & |
---|
79 | & 1.e-5*masse(ig,k)) then |
---|
80 | qa(ig,k)=(fm(ig,k)*qa(ig,k-1)+entr(ig,k)*q(ig,k)) & |
---|
81 | & /(fm(ig,k+1)+detr(ig,k)) |
---|
82 | else |
---|
83 | qa(ig,k)=q(ig,k) |
---|
84 | endif |
---|
85 | if (qa(ig,k).lt.0.) then |
---|
86 | ! print*,'qa<0!!!' |
---|
87 | endif |
---|
88 | if (q(ig,k).lt.0.) then |
---|
89 | ! print*,'q<0!!!' |
---|
90 | endif |
---|
91 | enddo |
---|
92 | enddo |
---|
93 | |
---|
94 | ! Calcul du flux subsident |
---|
95 | |
---|
96 | do k=2,nlay |
---|
97 | do ig=1,ngrid |
---|
98 | wqd(ig,k)=fm(ig,k)*q(ig,k) |
---|
99 | if (wqd(ig,k).lt.0.) then |
---|
100 | ! print*,'wqd<0!!!' |
---|
101 | endif |
---|
102 | enddo |
---|
103 | enddo |
---|
104 | do ig=1,ngrid |
---|
105 | wqd(ig,1)=0. |
---|
106 | wqd(ig,nlay+1)=0. |
---|
107 | enddo |
---|
108 | |
---|
109 | |
---|
110 | ! Calcul des tendances |
---|
111 | do k=1,nlay |
---|
112 | do ig=1,ngrid |
---|
113 | q(ig,k)=q(ig,k)+(detr(ig,k)*qa(ig,k)-entr(ig,k)*q(ig,k) & |
---|
114 | & -wqd(ig,k)+wqd(ig,k+1)) & |
---|
115 | & *ptimestep/masse(ig,k) |
---|
116 | enddo |
---|
117 | enddo |
---|
118 | |
---|
119 | endif |
---|
120 | |
---|
121 | tke(:,:,isrf)=q(:,:) |
---|
122 | |
---|
123 | enddo |
---|
124 | |
---|
125 | return |
---|
126 | end |
---|
127 | END MODULE lmdz_thermcell_dtke |
---|