1 |
module cv30_compress_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
contains |
6 |
|
7 |
SUBROUTINE cv30_compress(ncum, 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, iflag, nk, icb, icbs, plcl, & |
10 |
tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, & |
11 |
ph, 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):: ncum |
20 |
integer iflag1(klon), nk1(klon), icb1(klon), icbs1(klon) |
21 |
real plcl1(klon), tnk1(klon), qnk1(klon), gznk1(klon) |
22 |
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 |
real, intent(in):: p1(klon, klev), ph1(klon, klev+1) |
28 |
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 |
|
33 |
! outputs: |
34 |
integer iflag(klon), nk(klon), icb(klon), icbs(klon) |
35 |
real plcl(klon), tnk(klon), qnk(klon), gznk(klon) |
36 |
real pbase(klon), buoybase(klon) |
37 |
real t(klon, klev), q(klon, klev), qs(klon, klev) |
38 |
real u(klon, klev), v(klon, klev) |
39 |
real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev) |
40 |
real p(klon, klev), ph(klon, klev+1), tv(klon, klev), tp(klon, klev) |
41 |
real tvp(klon, klev), clw(klon, klev) |
42 |
real th(klon, klev) |
43 |
real sig(klon, klev), w0(klon, klev) |
44 |
|
45 |
! Local: |
46 |
integer i, k, nn |
47 |
|
48 |
!--------------------------------------------------------------- |
49 |
|
50 |
do k=1, nl+1 |
51 |
nn=0 |
52 |
do i=1, klon |
53 |
if (iflag1(i) == 0) then |
54 |
nn=nn+1 |
55 |
sig(nn, k)=sig1(i, k) |
56 |
w0(nn, k)=w01(i, k) |
57 |
t(nn, k)=t1(i, k) |
58 |
q(nn, k)=q1(i, k) |
59 |
qs(nn, k)=qs1(i, k) |
60 |
u(nn, k)=u1(i, k) |
61 |
v(nn, k)=v1(i, k) |
62 |
gz(nn, k)=gz1(i, k) |
63 |
h(nn, k)=h1(i, k) |
64 |
lv(nn, k)=lv1(i, k) |
65 |
cpn(nn, k)=cpn1(i, k) |
66 |
p(nn, k)=p1(i, k) |
67 |
ph(nn, k)=ph1(i, k) |
68 |
tv(nn, k)=tv1(i, k) |
69 |
tp(nn, k)=tp1(i, k) |
70 |
tvp(nn, k)=tvp1(i, k) |
71 |
clw(nn, k)=clw1(i, k) |
72 |
th(nn, k)=th1(i, k) |
73 |
endif |
74 |
end do |
75 |
end do |
76 |
|
77 |
if (nn /= ncum) then |
78 |
print*, 'strange! nn not equal to ncum: ', nn, ncum |
79 |
stop 1 |
80 |
endif |
81 |
|
82 |
nn=0 |
83 |
do i=1, klon |
84 |
if (iflag1(i) == 0) then |
85 |
nn=nn+1 |
86 |
pbase(nn)=pbase1(i) |
87 |
buoybase(nn)=buoybase1(i) |
88 |
plcl(nn)=plcl1(i) |
89 |
tnk(nn)=tnk1(i) |
90 |
qnk(nn)=qnk1(i) |
91 |
gznk(nn)=gznk1(i) |
92 |
nk(nn)=nk1(i) |
93 |
icb(nn)=icb1(i) |
94 |
icbs(nn)=icbs1(i) |
95 |
iflag(nn)=iflag1(i) |
96 |
endif |
97 |
end do |
98 |
|
99 |
end SUBROUTINE cv30_compress |
100 |
|
101 |
end module cv30_compress_m |