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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 205 - (hide annotations)
Tue Jun 21 15:16:03 2016 UTC (7 years, 11 months ago) by guez
File size: 2812 byte(s)
dnwd0 is just - mp. Compute it simply in concvl.

da, phi and mp were set to 0 in physiq before the call to
concvl. Clearer to set da1, phi1 and mp1 to 0 in cv_driver so they are
intent out.

qcheck was debugging, printed to standard output and was called
several times per time step of physics.

zxtsol was a duplicate of ztsol.

1 guez 185 module cv30_uncompress_m
2 guez 47
3     implicit none
4    
5 guez 97 contains
6 guez 47
7 guez 185 SUBROUTINE cv30_uncompress(idcum, iflag, precip, VPrecip, sig, w0, ft, fq, &
8 guez 205 fu, fv, inb, Ma, upwd, dnwd, qcondc, cape, da, phi, mp, iflag1, &
9 guez 189 precip1, VPrecip1, sig1, w01, ft1, fq1, fu1, fv1, inb1, Ma1, upwd1, &
10 guez 205 dnwd1, qcondc1, cape1, da1, phi1, mp1)
11 guez 47
12 guez 196 ! UNCOMPRESS THE FIELDS
13    
14 guez 185 USE cv30_param_m, ONLY: nl
15 guez 103 use dimphy, only: klon, klev
16 guez 47
17 guez 103 integer, intent(in):: idcum(:) ! (ncum)
18 guez 196 integer, intent(in):: iflag(:) ! (ncum)
19 guez 103 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 guez 201 integer, intent(in):: inb(:) ! (ncum)
24 guez 103 real, intent(in):: Ma(klon, klev)
25 guez 205 real, intent(in):: upwd(klon, klev), dnwd(klon, klev)
26 guez 103 real, intent(in):: qcondc(klon, klev)
27 guez 189 real, intent(in):: cape(klon)
28 guez 201 real, intent(in):: da(:, :) ! (klon, klev)
29     real, intent(in):: phi(:, :, :) ! (klon, klev, klev)
30     real, intent(in):: mp(:, :) ! (ncum, nl)
31 guez 47
32 guez 97 ! outputs:
33 guez 196 integer, intent(out):: iflag1(:) ! (klon)
34 guez 103 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 guez 205 real upwd1(klon, klev), dnwd1(klon, klev)
41 guez 103 real qcondc1(klon, klev)
42 guez 189 real cape1(klon)
43 guez 103 real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
44     real, intent(inout):: mp1(klon, klev)
45 guez 47
46 guez 180 ! Local:
47 guez 103 integer ncum, i, k, j
48 guez 47
49 guez 97 !-------------------------------------------------------------------
50 guez 47
51 guez 103 ncum = size(idcum)
52 guez 196 iflag1 = 42 ! for non convective points
53 guez 103
54 guez 205 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 guez 97 end do
60 guez 47
61 guez 205 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     qcondc1(idcum(i), k) = qcondc(i, k)
74     da1(idcum(i), k) = da(i, k)
75     mp1(idcum(i), k) = mp(i, k)
76 guez 97 end do
77     end do
78    
79 guez 205 do i = 1, ncum
80     sig1(idcum(i), klev) = sig(i, klev)
81 guez 97 end do
82    
83 guez 205 do j = 1, klev
84     do k = 1, klev
85     do i = 1, ncum
86     phi1(idcum(i), k, j) = phi(i, k, j)
87 guez 97 end do
88     end do
89     end do
90    
91 guez 185 end SUBROUTINE cv30_uncompress
92 guez 97
93 guez 185 end module cv30_uncompress_m

  ViewVC Help
Powered by ViewVC 1.1.21