/[lmdze]/trunk/phylmd/CV3_routines/cv3_uncompress.f
ViewVC logotype

Contents of /trunk/phylmd/CV3_routines/cv3_uncompress.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 103 - (show annotations)
Fri Aug 29 13:00:05 2014 UTC (9 years, 8 months ago) by guez
File size: 2754 byte(s)
Renamed module cvparam to cv_param. Deleted procedure
cv_param. Changed variables of module cv_param into parameters.

In procedures cv_driver, cv_uncompress and cv3_uncompress, removed
some arguments giving dimensions and used module variables klon and
klev instead.

In procedures gradiv2, laplacien_gam and laplacien, changed
declarations of local variables because klevel is not always klev.

Removed code for nudging surface pressure.

Removed arguments pim and pjm of tau2alpha. Added assignment of false
to variable first.

Replaced real argument del of procedures foeew and FOEDE by logical
argument.

1 module cv3_uncompress_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE cv3_uncompress(idcum, iflag, precip, VPrecip, sig, w0, ft, fq, &
8 fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, wd, cape, da, phi, mp, &
9 iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, fu1, fv1, inb1, Ma1, &
10 upwd1, dnwd1, dnwd01, qcondc1, wd1, cape1, da1, phi1, mp1)
11
12 USE cv3_param_m, ONLY: nl
13 use dimphy, only: klon, klev
14
15 ! inputs:
16 integer, intent(in):: idcum(:) ! (ncum)
17 integer, intent(in):: iflag(klon)
18 real, intent(in):: precip(klon)
19 real, intent(in):: VPrecip(klon, klev+1)
20 real, intent(in):: sig(klon, klev), w0(klon, klev)
21 real, intent(in), dimension(klon, klev):: ft, fq, fu, fv
22 integer, intent(in):: inb(klon)
23 real, intent(in):: Ma(klon, klev)
24 real, intent(in):: upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
25 real, intent(in):: qcondc(klon, klev)
26 real, intent(in):: wd(klon), cape(klon)
27 real, intent(in):: da(klon, klev), phi(klon, klev, klev), mp(klon, klev)
28
29 ! outputs:
30 integer iflag1(klon)
31 real precip1(klon)
32 real VPrecip1(klon, klev+1)
33 real sig1(klon, klev), w01(klon, klev)
34 real ft1(klon, klev), fq1(klon, klev), fu1(klon, klev), fv1(klon, klev)
35 integer, intent(inout):: inb1(klon)
36 real Ma1(klon, klev)
37 real upwd1(klon, klev), dnwd1(klon, klev), dnwd01(klon, klev)
38 real qcondc1(klon, klev)
39 real wd1(klon), cape1(klon)
40 real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
41 real, intent(inout):: mp1(klon, klev)
42
43 ! local variables:
44 integer ncum, i, k, j
45
46 !-------------------------------------------------------------------
47
48 ncum = size(idcum)
49
50 do i=1, ncum
51 precip1(idcum(i))=precip(i)
52 iflag1(idcum(i))=iflag(i)
53 wd1(idcum(i))=wd(i)
54 inb1(idcum(i))=inb(i)
55 cape1(idcum(i))=cape(i)
56 end do
57
58 do k=1, nl
59 do i=1, ncum
60 VPrecip1(idcum(i), k)=VPrecip(i, k)
61 sig1(idcum(i), k)=sig(i, k)
62 w01(idcum(i), k)=w0(i, k)
63 ft1(idcum(i), k)=ft(i, k)
64 fq1(idcum(i), k)=fq(i, k)
65 fu1(idcum(i), k)=fu(i, k)
66 fv1(idcum(i), k)=fv(i, k)
67 Ma1(idcum(i), k)=Ma(i, k)
68 upwd1(idcum(i), k)=upwd(i, k)
69 dnwd1(idcum(i), k)=dnwd(i, k)
70 dnwd01(idcum(i), k)=dnwd0(i, k)
71 qcondc1(idcum(i), k)=qcondc(i, k)
72 da1(idcum(i), k)=da(i, k)
73 mp1(idcum(i), k)=mp(i, k)
74 end do
75 end do
76
77 do i=1, ncum
78 sig1(idcum(i), klev)=sig(i, klev)
79 end do
80
81 do j=1, klev
82 do k=1, klev
83 do i=1, ncum
84 phi1(idcum(i), k, j)=phi(i, k, j)
85 end do
86 end do
87 end do
88
89 end SUBROUTINE cv3_uncompress
90
91 end module cv3_uncompress_m

  ViewVC Help
Powered by ViewVC 1.1.21