/[lmdze]/trunk/libf/phylmd/Thermcell/dqthermcell2.f90
ViewVC logotype

Contents of /trunk/libf/phylmd/Thermcell/dqthermcell2.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 47 - (show annotations)
Fri Jul 1 15:00:48 2011 UTC (12 years, 10 months ago) by guez
File size: 2061 byte(s)
Split "thermcell.f" and "cv3_routines.f".
Removed copies of files that are now in "L_util".
Moved "mva9" and "diagetpq" to their own files.
Unified variable names across procedures.

1
2 subroutine dqthermcell2(ngrid,nlay,ptimestep,fm,entr,masse,frac &
3 ,q,dq,qa)
4 use dimens_m
5 use dimphy
6 implicit none
7
8 !=======================================================================
9 !
10 ! Calcul du transport verticale dans la couche limite en presence
11 ! de "thermiques" explicitement representes
12 ! calcul du dq/dt une fois qu'on connait les ascendances
13 !
14 !=======================================================================
15
16
17 integer ngrid,nlay
18
19 real ptimestep
20 real masse(ngrid,nlay),fm(ngrid,nlay+1)
21 real entr(ngrid,nlay),frac(ngrid,nlay)
22 real q(ngrid,nlay)
23 real dq(ngrid,nlay)
24
25 real qa(klon,klev),detr(klon,klev),wqd(klon,klev+1)
26 real qe(klon,klev),zf,zf2
27
28 integer ig,k
29
30 ! calcul du detrainement
31
32 do k=1,nlay
33 do ig=1,ngrid
34 detr(ig,k)=fm(ig,k)-fm(ig,k+1)+entr(ig,k)
35 enddo
36 enddo
37
38 ! calcul de la valeur dans les ascendances
39 do ig=1,ngrid
40 qa(ig,1)=q(ig,1)
41 qe(ig,1)=q(ig,1)
42 enddo
43
44 do k=2,nlay
45 do ig=1,ngrid
46 if ((fm(ig,k+1)+detr(ig,k))*ptimestep.gt. &
47 1.e-5*masse(ig,k)) then
48 zf=0.5*(frac(ig,k)+frac(ig,k+1))
49 zf2=1./(1.-zf)
50 qa(ig,k)=(fm(ig,k)*qa(ig,k-1)+zf2*entr(ig,k)*q(ig,k)) &
51 /(fm(ig,k+1)+detr(ig,k)+entr(ig,k)*zf*zf2)
52 qe(ig,k)=(q(ig,k)-zf*qa(ig,k))*zf2
53 else
54 qa(ig,k)=q(ig,k)
55 qe(ig,k)=q(ig,k)
56 endif
57 enddo
58 enddo
59
60 do k=2,nlay
61 do ig=1,ngrid
62 ! wqd(ig,k)=fm(ig,k)*0.5*(q(ig,k-1)+q(ig,k))
63 wqd(ig,k)=fm(ig,k)*qe(ig,k)
64 enddo
65 enddo
66 do ig=1,ngrid
67 wqd(ig,1)=0.
68 wqd(ig,nlay+1)=0.
69 enddo
70
71 do k=1,nlay
72 do ig=1,ngrid
73 dq(ig,k)=(detr(ig,k)*qa(ig,k)-entr(ig,k)*qe(ig,k) &
74 -wqd(ig,k)+wqd(ig,k+1)) &
75 /masse(ig,k)
76 enddo
77 enddo
78
79 return
80 end

  ViewVC Help
Powered by ViewVC 1.1.21