/[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 201 - (hide annotations)
Mon Jun 6 17:42:15 2016 UTC (7 years, 11 months ago) by guez
File size: 2461 byte(s)
Removed intermediary objects of cv_thermo_m, access suphec_m
directly. Procedure cv_thermo disappeared, all objects are named
constants.

In cv_driver and below, limited extents of arrays to what is needed.

lv, cpn and th in cv30_compress were set at level nl + 1 but lv1, cpn1
and th1 are not defined at this level. This did not lead to an error
because values at nl + 1 were not used.

Removed test on ok_sync in phystokenc because it is not read at run
time. Printing min and max of output NetCDF variables is heavy and
archaic.

Used histwrite_phy in phytrac.

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 201 USE cv_thermo_m, ONLY: clmcpv, eps
11 guez 198 USE dimphy, ONLY: klev, klon
12 guez 201 use SUPHEC_M, only: rcw, rlvtt, rcpd, rcpv, rd, rv
13 guez 47
14 guez 97 ! Calculate arrays of geopotential, heat capacity and static energy
15 guez 47
16 guez 201 real, intent(in):: t1(:, :) ! (klon, klev) temperature, in K
17     real, intent(in):: q1(:, :) ! (klon, klev) specific humidity
18     real, intent(in):: p1(:, :) ! (klon, klev) full level pressure, in hPa
19     real, intent(in):: ph1(:, :) ! (klon, klev + 1) half level pressure, in hPa
20 guez 47
21 guez 97 ! outputs:
22 guez 201
23     real, intent(out):: lv1(:, :) ! (klon, nl)
24     ! specific latent heat of vaporization of water, in J kg-1
25    
26     real, intent(out):: cpn1(:, :) ! (klon, nl)
27     ! specific heat capacity at constant pressure of humid air, in J K-1 kg-1
28    
29     real tv1(:, :) ! (klon, klev)
30 guez 198 real gz1(klon, klev), h1(klon, klev), hm1(klon, klev)
31 guez 201 real, intent(out):: th1(:, :) ! (klon, nl) potential temperature, in K
32 guez 47
33 guez 97 ! Local:
34     integer k, i
35 guez 201 real kappa
36 guez 97 real tvx, tvy
37 guez 198 real cpx(klon, klev)
38 guez 47
39 guez 97 !--------------------------------------------------------------
40 guez 47
41 guez 198 do k = 1, nl
42     do i = 1, klon
43     lv1(i, k) = rlvtt - clmcpv * (t1(i, k) - 273.15)
44 guez 201 cpn1(i, k) = rcpd * (1. - q1(i, k)) + rcpv * q1(i, k)
45     cpx(i, k) = rcpd * (1. - q1(i, k)) + rcw * q1(i, k)
46     tv1(i, k) = t1(i, k) * (1. + q1(i, k) / eps - q1(i, k))
47     kappa = (rd * (1. - q1(i, k)) + q1(i, k) * rv) / cpn1(i, k)
48     th1(i, k) = t1(i, k) * (1000. / p1(i, k))**kappa
49 guez 97 end do
50     end do
51    
52 guez 198 ! gz1 = phi at the full levels (same as p1).
53 guez 97
54 guez 198 do i = 1, klon
55     gz1(i, 1) = 0.
56 guez 97 end do
57    
58 guez 198 do k = 2, nl
59     do i = 1, klon
60 guez 201 tvx = t1(i, k) * (1. + q1(i, k) / eps - q1(i, k))
61     tvy = t1(i, k - 1) * (1. + q1(i, k - 1) / eps - q1(i, k - 1))
62     gz1(i, k) = gz1(i, k - 1) + 0.5 * rd * (tvx + tvy) &
63     * (p1(i, k - 1) - p1(i, k)) / ph1(i, k)
64 guez 97 end do
65     end do
66    
67 guez 198 ! h1 = phi + cpT (dry static energy).
68     ! hm1 = phi + cp(T1 - Tbase) + Lq
69 guez 97
70 guez 198 do k = 1, nl
71     do i = 1, klon
72     h1(i, k) = gz1(i, k) + cpn1(i, k) * t1(i, k)
73     hm1(i, k) = gz1(i, k) + cpx(i, k) * (t1(i, k) - t1(i, 1)) &
74     + lv1(i, k) * q1(i, k)
75 guez 97 end do
76     end do
77    
78 guez 185 end SUBROUTINE cv30_prelim
79 guez 97
80 guez 185 end module cv30_prelim_m

  ViewVC Help
Powered by ViewVC 1.1.21