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

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

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.21