/[lmdze]/trunk/libf/phylmd/CV3_routines/cv3_prelim.f90
ViewVC logotype

Contents of /trunk/libf/phylmd/CV3_routines/cv3_prelim.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: 2161 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 cv3_prelim(len,nd,ndp1,t,q,p,ph &
3 ,lv,cpn,tv,gz,h,hm,th)
4 use cvparam3
5 use cvthermo
6 implicit none
7
8 !=====================================================================
9 ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
10 ! "ori": from convect4.3 (vectorized)
11 ! "convect3": to be exactly consistent with convect3
12 !=====================================================================
13
14 ! inputs:
15 integer len, nd, ndp1
16 real t(len,nd), q(len,nd), p(len,nd), ph(len,ndp1)
17
18 ! outputs:
19 real lv(len,nd), cpn(len,nd), tv(len,nd)
20 real gz(len,nd), h(len,nd), hm(len,nd)
21 real th(len,nd)
22
23 ! local variables:
24 integer k, i
25 real rdcp
26 real tvx,tvy ! convect3
27 real cpx(len,nd)
28
29
30
31 ! ori do 110 k=1,nlp
32 do 110 k=1,nl ! convect3
33 do 100 i=1,len
34 !debug lv(i,k)= lv0-clmcpv*(t(i,k)-t0)
35 lv(i,k)= lv0-clmcpv*(t(i,k)-273.15)
36 cpn(i,k)=cpd*(1.0-q(i,k))+cpv*q(i,k)
37 cpx(i,k)=cpd*(1.0-q(i,k))+cl*q(i,k)
38 ! ori tv(i,k)=t(i,k)*(1.0+q(i,k)*epsim1)
39 tv(i,k)=t(i,k)*(1.0+q(i,k)/eps-q(i,k))
40 rdcp=(rrd*(1.-q(i,k))+q(i,k)*rrv)/cpn(i,k)
41 th(i,k)=t(i,k)*(1000.0/p(i,k))**rdcp
42 100 continue
43 110 continue
44 !
45 ! gz = phi at the full levels (same as p).
46 !
47 do 120 i=1,len
48 gz(i,1)=0.0
49 120 continue
50 ! ori do 140 k=2,nlp
51 do 140 k=2,nl ! convect3
52 do 130 i=1,len
53 tvx=t(i,k)*(1.+q(i,k)/eps-q(i,k)) !convect3
54 tvy=t(i,k-1)*(1.+q(i,k-1)/eps-q(i,k-1)) !convect3
55 gz(i,k)=gz(i,k-1)+0.5*rrd*(tvx+tvy) &
56 *(p(i,k-1)-p(i,k))/ph(i,k) !convect3
57
58 ! ori gz(i,k)=gz(i,k-1)+hrd*(tv(i,k-1)+tv(i,k))
59 ! ori & *(p(i,k-1)-p(i,k))/ph(i,k)
60 130 continue
61 140 continue
62 !
63 ! h = phi + cpT (dry static energy).
64 ! hm = phi + cp(T-Tbase)+Lq
65 !
66 ! ori do 170 k=1,nlp
67 do 170 k=1,nl ! convect3
68 do 160 i=1,len
69 h(i,k)=gz(i,k)+cpn(i,k)*t(i,k)
70 hm(i,k)=gz(i,k)+cpx(i,k)*(t(i,k)-t(i,1))+lv(i,k)*q(i,k)
71 160 continue
72 170 continue
73
74 return
75 end

  ViewVC Help
Powered by ViewVC 1.1.21