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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 190 - (show annotations)
Thu Apr 14 15:15:56 2016 UTC (8 years, 1 month ago) by guez
File size: 1799 byte(s)
Created module cv_thermo_m around procedure cv_thermo. Moved variables
from module cvthermo to module cv_thermo_m, where they are defined.

In ini_histins and initphysto, using part of rlon and rlat from
phyetat0_m is pretending that we do not know about the dynamical grid,
while the way we extract zx_lon(:, 1) and zx_lat(1, :) depends on
ordering inside rlon and rlat. So we might as well simplify and
clarify by using directly rlonv and rlatu.

Removed intermediary variables in write_histins and phystokenc.

1 module cv30_prelim_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE cv30_prelim(len, nd, ndp1, t, q, p, ph, lv, cpn, tv, gz, h, hm, th)
8
9 USE cv30_param_m, ONLY: nl
10 USE cv_thermo_m, 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 cv30_prelim
69
70 end module cv30_prelim_m

  ViewVC Help
Powered by ViewVC 1.1.21