/[lmdze]/trunk/Sources/phylmd/CV30_routines/cv30_prelim.f
ViewVC logotype

Annotation of /trunk/Sources/phylmd/CV30_routines/cv30_prelim.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 134 - (hide annotations)
Wed Apr 29 15:47:56 2015 UTC (9 years, 1 month ago) by guez
Original Path: trunk/Sources/phylmd/CV3_routines/cv3_prelim.f
File size: 1791 byte(s)
Sources inside, compilation outside.
1 guez 97 module cv3_prelim_m
2 guez 47
3 guez 97 implicit none
4 guez 47
5 guez 97 contains
6 guez 47
7 guez 97 SUBROUTINE cv3_prelim(len, nd, ndp1, t, q, p, ph, lv, cpn, tv, gz, h, hm, th)
8 guez 47
9 guez 97 USE cv3_param_m, ONLY: nl
10     USE cvthermo, ONLY: cl, clmcpv, cpd, cpv, eps, lv0, rrd, rrv
11 guez 47
12 guez 97 ! Calculate arrays of geopotential, heat capacity and static energy
13 guez 47
14 guez 97 integer, intent(in):: len, nd, ndp1
15     real, intent(in):: t(len, nd)
16     real, intent(in):: q(len, nd)
17     real, intent(in):: p(len, nd), ph(len, ndp1)
18 guez 47
19 guez 97 ! 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 guez 47
24 guez 97 ! Local:
25     integer k, i
26     real rdcp
27     real tvx, tvy
28     real cpx(len, nd)
29 guez 47
30 guez 97 !--------------------------------------------------------------
31 guez 47
32 guez 97 do k=1, nl
33     do i=1, len
34     lv(i, k)= lv0-clmcpv*(t(i, k)-273.15)
35     cpn(i, k)=cpd*(1.0-q(i, k)) + cpv*q(i, k)
36     cpx(i, k)=cpd*(1.0-q(i, k)) + cl*q(i, k)
37     tv(i, k)=t(i, k)*(1.0 + q(i, k)/eps-q(i, k))
38     rdcp=(rrd*(1.-q(i, k)) + q(i, k)*rrv)/cpn(i, k)
39     th(i, k)=t(i, k)*(1000.0/p(i, k))**rdcp
40     end do
41     end do
42    
43     ! gz = phi at the full levels (same as p).
44    
45     do i=1, len
46     gz(i, 1)=0.0
47     end do
48    
49     do k=2, nl
50     do i=1, len
51     tvx=t(i, k)*(1. + q(i, k)/eps-q(i, k))
52     tvy=t(i, k-1)*(1. + q(i, k-1)/eps-q(i, k-1))
53     gz(i, k)=gz(i, k-1) + 0.5*rrd*(tvx + tvy) &
54     *(p(i, k-1)-p(i, k))/ph(i, k)
55     end do
56     end do
57    
58     ! h = phi + cpT (dry static energy).
59     ! hm = phi + cp(T-Tbase) + Lq
60    
61     do k=1, nl
62     do i=1, len
63     h(i, k)=gz(i, k) + cpn(i, k)*t(i, k)
64     hm(i, k)=gz(i, k) + cpx(i, k)*(t(i, k)-t(i, 1)) + lv(i, k)*q(i, k)
65     end do
66     end do
67    
68     end SUBROUTINE cv3_prelim
69    
70     end module cv3_prelim_m

  ViewVC Help
Powered by ViewVC 1.1.21