/[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 201 - (hide annotations)
Mon Jun 6 17:42:15 2016 UTC (7 years, 11 months ago) by guez
File size: 3885 byte(s)
Removed intermediary objects of cv_thermo_m, access suphec_m
directly. Procedure cv_thermo disappeared, all objects are named
constants.

In cv_driver and below, limited extents of arrays to what is needed.

lv, cpn and th in cv30_compress were set at level nl + 1 but lv1, cpn1
and th1 are not defined at this level. This did not lead to an error
because values at nl + 1 were not used.

Removed test on ok_sync in phystokenc because it is not read at run
time. Printing min and max of output NetCDF variables is heavy and
archaic.

Used histwrite_phy in phytrac.

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 201 SUBROUTINE cv30_compress(idcum, iflag1, icb1, icbs1, plcl1, tnk1, qnk1, &
8     gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, &
9     cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, 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 201 integer, intent(in):: idcum(:) ! (ncum)
21 guez 198 integer, intent(in):: iflag1(:), icb1(:), icbs1(:) ! (klon)
22 guez 195 real, intent(in):: plcl1(klon), tnk1(klon), qnk1(klon), gznk1(klon)
23 guez 189 real pbase1(klon), buoybase1(klon)
24 guez 201 real, intent(in):: t1(klon, klev) ! temperature (K)
25 guez 189 real, intent(in):: q1(klon, klev), qs1(klon, klev)
26     real, intent(in):: u1(klon, klev), v1(klon, klev)
27 guez 201 real gz1(klon, klev), h1(klon, klev)
28    
29     real, intent(in):: lv1(:, :) ! (klon, nl)
30     ! specific latent heat of vaporization of water, in J kg-1
31    
32     real, intent(in):: cpn1(:, :) ! (klon, nl)
33     ! specific heat capacity at constant pressure of humid air, in J K-1 kg-1
34    
35 guez 196 real, intent(in):: p1(klon, klev), ph1(klon, klev + 1)
36 guez 189 real, intent(in):: tv1(klon, klev), tp1(klon, klev)
37     real tvp1(klon, klev), clw1(klon, klev)
38 guez 201 real, intent(in):: th1(:, :) ! (klon, nl) potential temperature, in K
39 guez 189 real sig1(klon, klev), w01(klon, klev)
40 guez 47
41 guez 91 ! outputs:
42 guez 197 integer, intent(out):: icb(:) ! (ncum) {2 <= icb <= nl - 3}
43 guez 195 integer icbs(klon)
44 guez 196 real, intent(out):: plcl(:) ! (ncum)
45     real tnk(:), qnk(:), gznk(:) ! (klon)
46 guez 189 real pbase(klon), buoybase(klon)
47 guez 201 real t(klon, klev) ! temperature (K)
48     real q(klon, klev), qs(klon, klev)
49 guez 189 real u(klon, klev), v(klon, klev)
50 guez 201 real gz(klon, klev), h(klon, klev)
51    
52     real, intent(out):: lv(:, :) ! (ncum, nl)
53     ! specific latent heat of vaporization of water, in J kg-1
54    
55     real cpn(:, :) ! (ncum, nl)
56     ! specific heat capacity at constant pressure of humid air, in J K-1 kg-1
57    
58 guez 196 real p(klon, klev)
59     real ph(:, :) ! (klon, klev + 1)
60     real tv(klon, klev), tp(klon, klev)
61 guez 189 real tvp(klon, klev), clw(klon, klev)
62 guez 201 real, intent(out):: th(:, :) ! (ncum, nl) potential temperature, in K
63 guez 189 real sig(klon, klev), w0(klon, klev)
64 guez 91
65 guez 189 ! Local:
66 guez 196 integer i, k, nn, ncum
67 guez 91
68 guez 189 !---------------------------------------------------------------
69 guez 91
70 guez 196 ncum = size(icb)
71    
72     do k = 1, nl + 1
73     nn = 0
74     do i = 1, klon
75 guez 189 if (iflag1(i) == 0) then
76 guez 196 nn = nn + 1
77     sig(nn, k) = sig1(i, k)
78     w0(nn, k) = w01(i, k)
79     t(nn, k) = t1(i, k)
80     q(nn, k) = q1(i, k)
81     qs(nn, k) = qs1(i, k)
82     u(nn, k) = u1(i, k)
83     v(nn, k) = v1(i, k)
84     gz(nn, k) = gz1(i, k)
85     h(nn, k) = h1(i, k)
86     p(nn, k) = p1(i, k)
87     ph(nn, k) = ph1(i, k)
88     tv(nn, k) = tv1(i, k)
89     tp(nn, k) = tp1(i, k)
90     tvp(nn, k) = tvp1(i, k)
91     clw(nn, k) = clw1(i, k)
92 guez 91 endif
93     end do
94     end do
95 guez 47
96 guez 201 th = th1(idcum, :)
97     lv = lv1(idcum, :)
98     cpn = cpn1(idcum, :)
99    
100 guez 196 nn = 0
101 guez 47
102 guez 196 do i = 1, klon
103 guez 189 if (iflag1(i) == 0) then
104 guez 196 nn = nn + 1
105     pbase(nn) = pbase1(i)
106     buoybase(nn) = buoybase1(i)
107     plcl(nn) = plcl1(i)
108     tnk(nn) = tnk1(i)
109     qnk(nn) = qnk1(i)
110     gznk(nn) = gznk1(i)
111     icb(nn) = icb1(i)
112     icbs(nn) = icbs1(i)
113 guez 91 endif
114     end do
115 guez 47
116 guez 197 do i = 1, ncum
117     call assert(2 <= icb(i) .and. icb(i) <= nl - 3 .and. ph(i, icb(i) + 1) &
118     < plcl(i) .and. (plcl(i) <= ph(i, icb(i)) .or. icb(i) == 2), &
119     "cv30_compress")
120     end do
121    
122 guez 185 end SUBROUTINE cv30_compress
123 guez 91
124 guez 185 end module cv30_compress_m

  ViewVC Help
Powered by ViewVC 1.1.21