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

Contents of /trunk/Sources/phylmd/CV30_routines/cv30_uncompress.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 201 - (show annotations)
Mon Jun 6 17:42:15 2016 UTC (7 years, 11 months ago) by guez
File size: 2856 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 module cv30_uncompress_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE cv30_uncompress(idcum, iflag, precip, VPrecip, sig, w0, ft, fq, &
8 fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, cape, da, phi, mp, iflag1, &
9 precip1, VPrecip1, sig1, w01, ft1, fq1, fu1, fv1, inb1, Ma1, upwd1, &
10 dnwd1, dnwd01, qcondc1, cape1, da1, phi1, mp1)
11
12 ! UNCOMPRESS THE FIELDS
13
14 USE cv30_param_m, ONLY: nl
15 use dimphy, only: klon, klev
16
17 integer, intent(in):: idcum(:) ! (ncum)
18 integer, intent(in):: iflag(:) ! (ncum)
19 real, intent(in):: precip(klon)
20 real, intent(in):: VPrecip(klon, klev+1)
21 real, intent(in):: sig(klon, klev), w0(klon, klev)
22 real, intent(in), dimension(klon, klev):: ft, fq, fu, fv
23 integer, intent(in):: inb(:) ! (ncum)
24 real, intent(in):: Ma(klon, klev)
25 real, intent(in):: upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
26 real, intent(in):: qcondc(klon, klev)
27 real, intent(in):: cape(klon)
28 real, intent(in):: da(:, :) ! (klon, klev)
29 real, intent(in):: phi(:, :, :) ! (klon, klev, klev)
30 real, intent(in):: mp(:, :) ! (ncum, nl)
31
32 ! outputs:
33 integer, intent(out):: iflag1(:) ! (klon)
34 real precip1(klon)
35 real VPrecip1(klon, klev+1)
36 real sig1(klon, klev), w01(klon, klev)
37 real ft1(klon, klev), fq1(klon, klev), fu1(klon, klev), fv1(klon, klev)
38 integer, intent(inout):: inb1(klon)
39 real Ma1(klon, klev)
40 real upwd1(klon, klev), dnwd1(klon, klev), dnwd01(klon, klev)
41 real qcondc1(klon, klev)
42 real cape1(klon)
43 real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
44 real, intent(inout):: mp1(klon, klev)
45
46 ! Local:
47 integer ncum, i, k, j
48
49 !-------------------------------------------------------------------
50
51 ncum = size(idcum)
52 iflag1 = 42 ! for non convective points
53
54 do i=1, ncum
55 precip1(idcum(i))=precip(i)
56 iflag1(idcum(i))=iflag(i)
57 inb1(idcum(i))=inb(i)
58 cape1(idcum(i))=cape(i)
59 end do
60
61 do k=1, nl
62 do i=1, ncum
63 VPrecip1(idcum(i), k)=VPrecip(i, k)
64 sig1(idcum(i), k)=sig(i, k)
65 w01(idcum(i), k)=w0(i, k)
66 ft1(idcum(i), k)=ft(i, k)
67 fq1(idcum(i), k)=fq(i, k)
68 fu1(idcum(i), k)=fu(i, k)
69 fv1(idcum(i), k)=fv(i, k)
70 Ma1(idcum(i), k)=Ma(i, k)
71 upwd1(idcum(i), k)=upwd(i, k)
72 dnwd1(idcum(i), k)=dnwd(i, k)
73 dnwd01(idcum(i), k)=dnwd0(i, k)
74 qcondc1(idcum(i), k)=qcondc(i, k)
75 da1(idcum(i), k)=da(i, k)
76 mp1(idcum(i), k)=mp(i, k)
77 end do
78 end do
79
80 do i=1, ncum
81 sig1(idcum(i), klev)=sig(i, klev)
82 end do
83
84 do j=1, klev
85 do k=1, klev
86 do i=1, ncum
87 phi1(idcum(i), k, j)=phi(i, k, j)
88 end do
89 end do
90 end do
91
92 end SUBROUTINE cv30_uncompress
93
94 end module cv30_uncompress_m

  ViewVC Help
Powered by ViewVC 1.1.21