/[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 197 - (hide annotations)
Tue May 24 12:25:29 2016 UTC (8 years ago) by guez
File size: 3405 byte(s)
Clarified the computation of sigt in cv30_unsat. Replacing pr2 by 1 -
pr1 changes the results.

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 196 SUBROUTINE cv30_compress(iflag1, nk1, icb1, icbs1, plcl1, tnk1, qnk1, &
8 guez 189 gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, cpn1, &
9 guez 196 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 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 196 integer, intent(in):: iflag1(:), nk1(:), 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 196 integer nk(:) ! (klon)
35 guez 197 integer, intent(out):: icb(:) ! (ncum) {2 <= icb <= nl - 3}
36 guez 195 integer icbs(klon)
37 guez 196 real, intent(out):: plcl(:) ! (ncum)
38     real tnk(:), qnk(:), gznk(:) ! (klon)
39 guez 189 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 guez 196 real p(klon, klev)
44     real ph(:, :) ! (klon, klev + 1)
45     real tv(klon, klev), tp(klon, klev)
46 guez 189 real tvp(klon, klev), clw(klon, klev)
47     real th(klon, klev)
48     real sig(klon, klev), w0(klon, klev)
49 guez 91
50 guez 189 ! Local:
51 guez 196 integer i, k, nn, ncum
52 guez 91
53 guez 189 !---------------------------------------------------------------
54 guez 91
55 guez 196 ncum = size(icb)
56    
57     do k = 1, nl + 1
58     nn = 0
59     do i = 1, klon
60 guez 189 if (iflag1(i) == 0) then
61 guez 196 nn = nn + 1
62     sig(nn, k) = sig1(i, k)
63     w0(nn, k) = w01(i, k)
64     t(nn, k) = t1(i, k)
65     q(nn, k) = q1(i, k)
66     qs(nn, k) = qs1(i, k)
67     u(nn, k) = u1(i, k)
68     v(nn, k) = v1(i, k)
69     gz(nn, k) = gz1(i, k)
70     h(nn, k) = h1(i, k)
71     lv(nn, k) = lv1(i, k)
72     cpn(nn, k) = cpn1(i, k)
73     p(nn, k) = p1(i, k)
74     ph(nn, k) = ph1(i, k)
75     tv(nn, k) = tv1(i, k)
76     tp(nn, k) = tp1(i, k)
77     tvp(nn, k) = tvp1(i, k)
78     clw(nn, k) = clw1(i, k)
79     th(nn, k) = th1(i, k)
80 guez 91 endif
81     end do
82     end do
83 guez 47
84 guez 196 nn = 0
85 guez 47
86 guez 196 do i = 1, klon
87 guez 189 if (iflag1(i) == 0) then
88 guez 196 nn = nn + 1
89     pbase(nn) = pbase1(i)
90     buoybase(nn) = buoybase1(i)
91     plcl(nn) = plcl1(i)
92     tnk(nn) = tnk1(i)
93     qnk(nn) = qnk1(i)
94     gznk(nn) = gznk1(i)
95     nk(nn) = nk1(i)
96     icb(nn) = icb1(i)
97     icbs(nn) = icbs1(i)
98 guez 91 endif
99     end do
100 guez 47
101 guez 197 do i = 1, ncum
102     call assert(2 <= icb(i) .and. icb(i) <= nl - 3 .and. ph(i, icb(i) + 1) &
103     < plcl(i) .and. (plcl(i) <= ph(i, icb(i)) .or. icb(i) == 2), &
104     "cv30_compress")
105     end do
106    
107 guez 185 end SUBROUTINE cv30_compress
108 guez 91
109 guez 185 end module cv30_compress_m

  ViewVC Help
Powered by ViewVC 1.1.21