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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 196 - (show annotations)
Mon May 23 13:50:39 2016 UTC (8 years ago) by guez
File size: 3221 byte(s)
Removed argument icbmax of cv30_feed, not used in cv_driver (not used
in LMDZ either).

Clearer to have iflag1 = 0 in cv30_feed than in cv_driver. Clearer to
have iflag1 = 42 in cv30_uncompress than in cv_driver.

Removed argument iflag of cv30_compress. Since there is iflag1 = 42 in
cv30_uncompress, iflag1 that comes out of cv_driver is disconnected
from iflag1 computed by cv30_feed and cv30_trigger.

Removed some references to convect3 and convect4 in comments. This
program derives from the convect3 version, we do not need to know
about other versions.

Bug fix in cv30_undilute1: icbs1 was not made >= 2.

1 module cv30_compress_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE cv30_compress(iflag1, nk1, icb1, icbs1, plcl1, tnk1, qnk1, &
8 gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, cpn1, &
9 p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, nk, icb, icbs, plcl, tnk, &
10 qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, &
11 tv, tp, tvp, clw, sig, w0)
12
13 ! Compress the fields (vectorization over convective gridpoints).
14
15 use cv30_param_m, only: nl
16 USE dimphy, ONLY: klev, klon
17
18 ! inputs:
19 integer, intent(in):: iflag1(:), nk1(:), icb1(:), icbs1(:) ! (klon)
20 real, intent(in):: plcl1(klon), tnk1(klon), qnk1(klon), gznk1(klon)
21 real pbase1(klon), buoybase1(klon)
22 real, intent(in):: t1(klon, klev)
23 real, intent(in):: q1(klon, klev), qs1(klon, klev)
24 real, intent(in):: u1(klon, klev), v1(klon, klev)
25 real gz1(klon, klev), h1(klon, klev), lv1(klon, klev), cpn1(klon, klev)
26 real, intent(in):: p1(klon, klev), ph1(klon, klev + 1)
27 real, intent(in):: tv1(klon, klev), tp1(klon, klev)
28 real tvp1(klon, klev), clw1(klon, klev)
29 real th1(klon, klev)
30 real sig1(klon, klev), w01(klon, klev)
31
32 ! outputs:
33 integer nk(:) ! (klon)
34
35 integer, intent(out):: icb(:) ! (ncum)
36 ! {2 <= icb <= nl - 3}
37 ! {ph(i, icb(i) + 1) < plcl(i) <= ph(i, icb(i))}
38
39 integer icbs(klon)
40 real, intent(out):: plcl(:) ! (ncum)
41 real tnk(:), qnk(:), gznk(:) ! (klon)
42 real pbase(klon), buoybase(klon)
43 real t(klon, klev), q(klon, klev), qs(klon, klev)
44 real u(klon, klev), v(klon, klev)
45 real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)
46 real p(klon, klev)
47 real ph(:, :) ! (klon, klev + 1)
48 real tv(klon, klev), tp(klon, klev)
49 real tvp(klon, klev), clw(klon, klev)
50 real th(klon, klev)
51 real sig(klon, klev), w0(klon, klev)
52
53 ! Local:
54 integer i, k, nn, ncum
55
56 !---------------------------------------------------------------
57
58 ncum = size(icb)
59
60 do k = 1, nl + 1
61 nn = 0
62 do i = 1, klon
63 if (iflag1(i) == 0) then
64 nn = nn + 1
65 sig(nn, k) = sig1(i, k)
66 w0(nn, k) = w01(i, k)
67 t(nn, k) = t1(i, k)
68 q(nn, k) = q1(i, k)
69 qs(nn, k) = qs1(i, k)
70 u(nn, k) = u1(i, k)
71 v(nn, k) = v1(i, k)
72 gz(nn, k) = gz1(i, k)
73 h(nn, k) = h1(i, k)
74 lv(nn, k) = lv1(i, k)
75 cpn(nn, k) = cpn1(i, k)
76 p(nn, k) = p1(i, k)
77 ph(nn, k) = ph1(i, k)
78 tv(nn, k) = tv1(i, k)
79 tp(nn, k) = tp1(i, k)
80 tvp(nn, k) = tvp1(i, k)
81 clw(nn, k) = clw1(i, k)
82 th(nn, k) = th1(i, k)
83 endif
84 end do
85 end do
86
87 nn = 0
88
89 do i = 1, klon
90 if (iflag1(i) == 0) then
91 nn = nn + 1
92 pbase(nn) = pbase1(i)
93 buoybase(nn) = buoybase1(i)
94 plcl(nn) = plcl1(i)
95 tnk(nn) = tnk1(i)
96 qnk(nn) = qnk1(i)
97 gznk(nn) = gznk1(i)
98 nk(nn) = nk1(i)
99 icb(nn) = icb1(i)
100 icbs(nn) = icbs1(i)
101 endif
102 end do
103
104 end SUBROUTINE cv30_compress
105
106 end module cv30_compress_m

  ViewVC Help
Powered by ViewVC 1.1.21