/[lmdze]/trunk/phylmd/CV3_routines/cv3_prelim.f
ViewVC logotype

Contents of /trunk/phylmd/CV3_routines/cv3_prelim.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 69 - (show annotations)
Mon Feb 18 16:33:12 2013 UTC (11 years, 3 months ago) by guez
Original Path: trunk/libf/phylmd/CV3_routines/cv3_prelim.f90
File size: 2212 byte(s)
Deleted files cvparam3.f90 and nuagecom.f90. Moved variables from
module cvparam3 to module cv3_param_m. Moved variables rad_chau1 and
rad_chau2 from module nuagecom to module conf_phys_m.

Read clesphys2_nml from conf_phys instead of gcm.

Removed argument iflag_con from several procedures. Access module
variable instead.

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

  ViewVC Help
Powered by ViewVC 1.1.21