1 |
module cv30_prelim_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
contains |
6 |
|
7 |
SUBROUTINE cv30_prelim(t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, h1, hm1, th1) |
8 |
|
9 |
USE cv30_param_m, ONLY: nl |
10 |
USE cv_thermo_m, ONLY: clmcpv, eps |
11 |
USE dimphy, ONLY: klev, klon |
12 |
use SUPHEC_M, only: rcw, rlvtt, rcpd, rcpv, rd, rv |
13 |
|
14 |
! Calculate arrays of geopotential, heat capacity and static energy |
15 |
|
16 |
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 |
|
21 |
! outputs: |
22 |
|
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 |
real gz1(klon, klev), h1(klon, klev), hm1(klon, klev) |
31 |
real, intent(out):: th1(:, :) ! (klon, nl) potential temperature, in K |
32 |
|
33 |
! Local: |
34 |
integer k, i |
35 |
real kappa |
36 |
real tvx, tvy |
37 |
real cpx(klon, klev) |
38 |
|
39 |
!-------------------------------------------------------------- |
40 |
|
41 |
do k = 1, nl |
42 |
do i = 1, klon |
43 |
lv1(i, k) = rlvtt - clmcpv * (t1(i, k) - 273.15) |
44 |
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 |
end do |
50 |
end do |
51 |
|
52 |
! gz1 = phi at the full levels (same as p1). |
53 |
|
54 |
do i = 1, klon |
55 |
gz1(i, 1) = 0. |
56 |
end do |
57 |
|
58 |
do k = 2, nl |
59 |
do i = 1, klon |
60 |
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 |
end do |
65 |
end do |
66 |
|
67 |
! h1 = phi + cpT (dry static energy). |
68 |
! hm1 = phi + cp(T1 - Tbase) + Lq |
69 |
|
70 |
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 |
end do |
76 |
end do |
77 |
|
78 |
end SUBROUTINE cv30_prelim |
79 |
|
80 |
end module cv30_prelim_m |