/[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 97 - (show annotations)
Fri Apr 25 14:58:31 2014 UTC (10 years ago) by guez
File size: 1791 byte(s)
Module pressure_var is now only used in gcm. Created local variables
pls and p3d in etat0, added argument p3d to regr_pr_o3.

In leapfrog, moved computation of p3d and exner function immediately
after integrd, for clarity (does not change the execution).

Removed unused arguments: ntra, tra1 and tra of cv3_compress; ntra,
tra and traent of cv3_mixing; ntra, ftra, ftra1 of cv3_uncompress;
ntra, tra, trap of cv3_unsat; ntra, tra, trap, traent, ftra of
cv3_yield; tra, tvp, pbase, bbase, dtvpdt1, dtvpdq1, dplcldt,
dplcldr, ntra of concvl; ndp1, ntra, tra1 of cv_driver

Removed argument d_tra and computation of d_tra in concvl. Removed
argument ftra1 and computation of ftra1 in cv_driver. ftra1 was just
set to 0 in cv_driver, associated to d_tra in concvl, and set again to
zero in concvl.

1 module cv3_prelim_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE cv3_prelim(len, nd, ndp1, t, q, p, ph, lv, cpn, tv, gz, h, hm, th)
8
9 USE cv3_param_m, ONLY: nl
10 USE cvthermo, ONLY: cl, clmcpv, cpd, cpv, eps, lv0, rrd, rrv
11
12 ! Calculate arrays of geopotential, heat capacity and static energy
13
14 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
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:
25 integer k, i
26 real rdcp
27 real tvx, tvy
28 real cpx(len, nd)
29
30 !--------------------------------------------------------------
31
32 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