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

Annotation of /trunk/libf/phylmd/CV3_routines/cv3_prelim.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 52 - (hide annotations)
Fri Sep 23 12:28:01 2011 UTC (12 years, 8 months ago) by guez
File size: 2185 byte(s)
Split "conflx.f" into single-procedure files in directory "Conflx".

Split "cv_routines.f" into single-procedure files in directory
"CV_routines". Made module "cvparam" from included file
"cvparam.h". No included file other than "netcdf.inc" left in LMDZE.

1 guez 47
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 guez 52 real, intent(in):: t(len,nd)
17     real q(len,nd), p(len,nd), ph(len,ndp1)
18 guez 47
19     ! outputs:
20     real lv(len,nd), cpn(len,nd), tv(len,nd)
21     real gz(len,nd), h(len,nd), hm(len,nd)
22     real th(len,nd)
23    
24     ! local variables:
25     integer k, i
26     real rdcp
27     real tvx,tvy ! convect3
28     real cpx(len,nd)
29    
30    
31    
32     ! ori do 110 k=1,nlp
33     do 110 k=1,nl ! convect3
34     do 100 i=1,len
35     !debug lv(i,k)= lv0-clmcpv*(t(i,k)-t0)
36     lv(i,k)= lv0-clmcpv*(t(i,k)-273.15)
37     cpn(i,k)=cpd*(1.0-q(i,k))+cpv*q(i,k)
38     cpx(i,k)=cpd*(1.0-q(i,k))+cl*q(i,k)
39     ! ori tv(i,k)=t(i,k)*(1.0+q(i,k)*epsim1)
40     tv(i,k)=t(i,k)*(1.0+q(i,k)/eps-q(i,k))
41     rdcp=(rrd*(1.-q(i,k))+q(i,k)*rrv)/cpn(i,k)
42     th(i,k)=t(i,k)*(1000.0/p(i,k))**rdcp
43     100 continue
44     110 continue
45     !
46     ! gz = phi at the full levels (same as p).
47     !
48     do 120 i=1,len
49     gz(i,1)=0.0
50     120 continue
51     ! ori do 140 k=2,nlp
52     do 140 k=2,nl ! convect3
53     do 130 i=1,len
54     tvx=t(i,k)*(1.+q(i,k)/eps-q(i,k)) !convect3
55     tvy=t(i,k-1)*(1.+q(i,k-1)/eps-q(i,k-1)) !convect3
56     gz(i,k)=gz(i,k-1)+0.5*rrd*(tvx+tvy) &
57     *(p(i,k-1)-p(i,k))/ph(i,k) !convect3
58    
59     ! ori gz(i,k)=gz(i,k-1)+hrd*(tv(i,k-1)+tv(i,k))
60     ! ori & *(p(i,k-1)-p(i,k))/ph(i,k)
61     130 continue
62     140 continue
63     !
64     ! h = phi + cpT (dry static energy).
65     ! hm = phi + cp(T-Tbase)+Lq
66     !
67     ! ori do 170 k=1,nlp
68     do 170 k=1,nl ! convect3
69     do 160 i=1,len
70     h(i,k)=gz(i,k)+cpn(i,k)*t(i,k)
71     hm(i,k)=gz(i,k)+cpx(i,k)*(t(i,k)-t(i,1))+lv(i,k)*q(i,k)
72     160 continue
73     170 continue
74    
75     return
76     end

  ViewVC Help
Powered by ViewVC 1.1.21