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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 76 - (hide annotations)
Fri Nov 15 18:45:49 2013 UTC (10 years, 6 months ago) by guez
Original Path: trunk/phylmd/CV3_routines/cv3_uncompress.f90
File size: 2338 byte(s)
Moved everything out of libf.
1 guez 47 SUBROUTINE cv3_uncompress(nloc, len, ncum, nd, ntra, idcum, iflag, precip, &
2     VPrecip, sig, w0, ft, fq, fu, fv, ftra, inb, Ma, upwd, dnwd, dnwd0, &
3     qcondc, wd, cape, da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, &
4     ft1, fq1, fu1, fv1, ftra1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, &
5     wd1, cape1, da1, phi1, mp1)
6    
7 guez 69 use cv3_param_m
8 guez 47
9     implicit none
10    
11     ! inputs:
12 guez 69 integer, intent(in):: len, ncum, nd, ntra, nloc
13 guez 47 integer idcum(nloc)
14     integer iflag(nloc)
15     integer inb(nloc)
16     real precip(nloc)
17     real VPrecip(nloc, nd+1)
18     real sig(nloc, nd), w0(nloc, nd)
19     real ft(nloc, nd), fq(nloc, nd), fu(nloc, nd), fv(nloc, nd)
20     real ftra(nloc, nd, ntra)
21     real Ma(nloc, nd)
22     real upwd(nloc, nd), dnwd(nloc, nd), dnwd0(nloc, nd)
23     real qcondc(nloc, nd)
24     real wd(nloc), cape(nloc)
25     real da(nloc, nd), phi(nloc, nd, nd), mp(nloc, nd)
26    
27     ! outputs:
28     integer iflag1(len)
29     integer inb1(len)
30     real precip1(len)
31     real VPrecip1(len, nd+1)
32     real sig1(len, nd), w01(len, nd)
33     real ft1(len, nd), fq1(len, nd), fu1(len, nd), fv1(len, nd)
34     real ftra1(len, nd, ntra)
35     real Ma1(len, nd)
36     real upwd1(len, nd), dnwd1(len, nd), dnwd01(len, nd)
37     real qcondc1(nloc, nd)
38     real wd1(nloc), cape1(nloc)
39     real da1(nloc, nd), phi1(nloc, nd, nd), mp1(nloc, nd)
40    
41     ! local variables:
42     integer i, k, j
43    
44     !-------------------------------------------------------------------
45    
46     do i=1, ncum
47     precip1(idcum(i))=precip(i)
48     iflag1(idcum(i))=iflag(i)
49     wd1(idcum(i))=wd(i)
50     inb1(idcum(i))=inb(i)
51     cape1(idcum(i))=cape(i)
52     end do
53    
54     do k=1, nl
55     do i=1, ncum
56     VPrecip1(idcum(i), k)=VPrecip(i, k)
57     sig1(idcum(i), k)=sig(i, k)
58     w01(idcum(i), k)=w0(i, k)
59     ft1(idcum(i), k)=ft(i, k)
60     fq1(idcum(i), k)=fq(i, k)
61     fu1(idcum(i), k)=fu(i, k)
62     fv1(idcum(i), k)=fv(i, k)
63     Ma1(idcum(i), k)=Ma(i, k)
64     upwd1(idcum(i), k)=upwd(i, k)
65     dnwd1(idcum(i), k)=dnwd(i, k)
66     dnwd01(idcum(i), k)=dnwd0(i, k)
67     qcondc1(idcum(i), k)=qcondc(i, k)
68     da1(idcum(i), k)=da(i, k)
69     mp1(idcum(i), k)=mp(i, k)
70     end do
71     end do
72    
73     do i=1, ncum
74     sig1(idcum(i), nd)=sig(i, nd)
75     end do
76    
77     do j=1, nd
78     do k=1, nd
79     do i=1, ncum
80     phi1(idcum(i), k, j)=phi(i, k, j)
81     end do
82     end do
83     end do
84    
85     end SUBROUTINE cv3_uncompress

  ViewVC Help
Powered by ViewVC 1.1.21