/[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 198 - (hide annotations)
Tue May 31 16:17:35 2016 UTC (8 years ago) by guez
File size: 2051 byte(s)
Removed variables nk1 and nk in cv_driver and below. These arrays were
just equal to the constant minorig. (This is also the case in LMDZ.)

In cv_thermo, removed some variables which were copies of variables of
suphec_m. Changed some variables to named constants.

1 guez 185 module cv30_prelim_m
2 guez 47
3 guez 97 implicit none
4 guez 47
5 guez 97 contains
6 guez 47
7 guez 198 SUBROUTINE cv30_prelim(t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, h1, hm1, th1)
8 guez 47
9 guez 185 USE cv30_param_m, ONLY: nl
10 guez 198 USE cv_thermo_m, ONLY: cl, clmcpv, cpd, cpv, eps, rrd, rrv
11     USE dimphy, ONLY: klev, klon
12     use SUPHEC_M, only: rlvtt
13 guez 47
14 guez 97 ! Calculate arrays of geopotential, heat capacity and static energy
15 guez 47
16 guez 198 real, intent(in):: t1(klon, klev)
17     real, intent(in):: q1(klon, klev)
18     real, intent(in):: p1(klon, klev), ph1(klon, klev + 1)
19 guez 47
20 guez 97 ! outputs:
21 guez 198 real lv1(klon, klev), cpn1(klon, klev), tv1(klon, klev)
22     real gz1(klon, klev), h1(klon, klev), hm1(klon, klev)
23     real th1(klon, klev) ! potential temperature
24 guez 47
25 guez 97 ! Local:
26     integer k, i
27     real rdcp
28     real tvx, tvy
29 guez 198 real cpx(klon, klev)
30 guez 47
31 guez 97 !--------------------------------------------------------------
32 guez 47
33 guez 198 do k = 1, nl
34     do i = 1, klon
35     lv1(i, k) = rlvtt - clmcpv * (t1(i, k) - 273.15)
36     cpn1(i, k) = cpd * (1. - q1(i, k)) + cpv * q1(i, k)
37     cpx(i, k) = cpd * (1. - q1(i, k)) + cl * q1(i, k)
38     tv1(i, k) = t1(i, k) * (1. + q1(i, k)/eps - q1(i, k))
39     rdcp = (rrd * (1. - q1(i, k)) + q1(i, k) * rrv)/cpn1(i, k)
40     th1(i, k) = t1(i, k) * (1000./p1(i, k))**rdcp
41 guez 97 end do
42     end do
43    
44 guez 198 ! gz1 = phi at the full levels (same as p1).
45 guez 97
46 guez 198 do i = 1, klon
47     gz1(i, 1) = 0.
48 guez 97 end do
49    
50 guez 198 do k = 2, nl
51     do i = 1, klon
52     tvx = t1(i, k) * (1. + q1(i, k)/eps - q1(i, k))
53     tvy = t1(i, k - 1) * (1. + q1(i, k - 1)/eps - q1(i, k - 1))
54     gz1(i, k) = gz1(i, k - 1) + 0.5 * rrd * (tvx + tvy) &
55     * (p1(i, k - 1) - p1(i, k))/ph1(i, k)
56 guez 97 end do
57     end do
58    
59 guez 198 ! h1 = phi + cpT (dry static energy).
60     ! hm1 = phi + cp(T1 - Tbase) + Lq
61 guez 97
62 guez 198 do k = 1, nl
63     do i = 1, klon
64     h1(i, k) = gz1(i, k) + cpn1(i, k) * t1(i, k)
65     hm1(i, k) = gz1(i, k) + cpx(i, k) * (t1(i, k) - t1(i, 1)) &
66     + lv1(i, k) * q1(i, k)
67 guez 97 end do
68     end do
69    
70 guez 185 end SUBROUTINE cv30_prelim
71 guez 97
72 guez 185 end module cv30_prelim_m

  ViewVC Help
Powered by ViewVC 1.1.21