/[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 195 - (show annotations)
Wed May 18 17:56:44 2016 UTC (8 years ago) by guez
File size: 3138 byte(s)
In cv30_feed, iflag1 is 0 on entry so we can simplify the test for
iflag1 = 7.

In cv30_feed, for the computation of icb, replaced sequential search
(with a useless end of loop on k) by a call to locate.

In CV30 routines, replaced len, nloc, nd, na by klon or
klev. Philosophy: no more generality than actually necessary.

Converted as many variables as possible to named constants in
cv30_param_m and downgraded pbcrit, ptcrit, dtovsh, dpbase, dttrig,
tau, delta to local objects in procedures. spfac, betad and omtrain
are useless and removed.

Instead of filling the array sigp with the constant spfac in
cv30_undilute2, just made sigp a constant in cv30_unsat.

In cv_driver, define as allocatable variables that are only
used on the range (ncum, nl).

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 use nr_util, only: assert
18
19 ! inputs:
20 integer, intent(in):: ncum
21 integer, intent(in):: iflag1(klon), nk1(klon), icb1(klon), icbs1(klon)
22 real, intent(in):: plcl1(klon), tnk1(klon), qnk1(klon), gznk1(klon)
23 real pbase1(klon), buoybase1(klon)
24 real, intent(in):: t1(klon, klev)
25 real, intent(in):: q1(klon, klev), qs1(klon, klev)
26 real, intent(in):: u1(klon, klev), v1(klon, klev)
27 real gz1(klon, klev), h1(klon, klev), lv1(klon, klev), cpn1(klon, klev)
28 real, intent(in):: p1(klon, klev), ph1(klon, klev+1)
29 real, intent(in):: tv1(klon, klev), tp1(klon, klev)
30 real tvp1(klon, klev), clw1(klon, klev)
31 real th1(klon, klev)
32 real sig1(klon, klev), w01(klon, klev)
33
34 ! outputs:
35 integer iflag(klon), nk(klon)
36 integer, intent(out):: icb(:) ! (ncum)
37 integer icbs(klon)
38 real plcl(klon), tnk(klon), qnk(klon), gznk(klon)
39 real pbase(klon), buoybase(klon)
40 real t(klon, klev), q(klon, klev), qs(klon, klev)
41 real u(klon, klev), v(klon, klev)
42 real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)
43 real p(klon, klev), ph(klon, klev+1), tv(klon, klev), tp(klon, klev)
44 real tvp(klon, klev), clw(klon, klev)
45 real th(klon, klev)
46 real sig(klon, klev), w0(klon, klev)
47
48 ! Local:
49 integer i, k, nn
50
51 !---------------------------------------------------------------
52
53 do k=1, nl+1
54 nn=0
55 do i=1, klon
56 if (iflag1(i) == 0) then
57 nn=nn+1
58 sig(nn, k)=sig1(i, k)
59 w0(nn, k)=w01(i, k)
60 t(nn, k)=t1(i, k)
61 q(nn, k)=q1(i, k)
62 qs(nn, k)=qs1(i, k)
63 u(nn, k)=u1(i, k)
64 v(nn, k)=v1(i, k)
65 gz(nn, k)=gz1(i, k)
66 h(nn, k)=h1(i, k)
67 lv(nn, k)=lv1(i, k)
68 cpn(nn, k)=cpn1(i, k)
69 p(nn, k)=p1(i, k)
70 ph(nn, k)=ph1(i, k)
71 tv(nn, k)=tv1(i, k)
72 tp(nn, k)=tp1(i, k)
73 tvp(nn, k)=tvp1(i, k)
74 clw(nn, k)=clw1(i, k)
75 th(nn, k)=th1(i, k)
76 endif
77 end do
78 end do
79
80 call assert(nn == ncum, "cv30_compress")
81 nn=0
82
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

  ViewVC Help
Powered by ViewVC 1.1.21