/[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 185 - (hide annotations)
Wed Mar 16 15:04:46 2016 UTC (8 years, 2 months ago) by guez
File size: 2860 byte(s)
CV3 to CV30 (following LMDZ) (continued).
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 185 SUBROUTINE cv30_compress(len, nloc, ncum, nd, iflag1, nk1, icb1, icbs1, &
8 guez 97 plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, &
9     th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, iflag, &
10     nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, &
11     gz, th, h, lv, cpn, p, ph, tv, tp, tvp, clw, sig, w0)
12 guez 47
13 guez 185 use cv30_param_m
14 guez 47
15    
16 guez 91 ! inputs:
17 guez 97 integer, intent(in):: len, ncum, nd, nloc
18 guez 91 integer iflag1(len), nk1(len), icb1(len), icbs1(len)
19     real plcl1(len), tnk1(len), qnk1(len), gznk1(len)
20     real pbase1(len), buoybase1(len)
21     real, intent(in):: t1(len, nd)
22     real, intent(in):: q1(len, nd), qs1(len, nd), u1(len, nd), v1(len, nd)
23     real gz1(len, nd), h1(len, nd), lv1(len, nd), cpn1(len, nd)
24 guez 103 real, intent(in):: p1(len, nd), ph1(len, nd+1), tv1(len, nd), tp1(len, nd)
25 guez 91 real tvp1(len, nd), clw1(len, nd)
26     real th1(len, nd)
27     real sig1(len, nd), w01(len, nd)
28 guez 47
29 guez 91 ! outputs:
30     ! en fait, on a nloc=len pour l'instant (cf cv_driver)
31     integer iflag(nloc), nk(nloc), icb(nloc), icbs(nloc)
32     real plcl(nloc), tnk(nloc), qnk(nloc), gznk(nloc)
33     real pbase(nloc), buoybase(nloc)
34     real t(nloc, nd), q(nloc, nd), qs(nloc, nd), u(nloc, nd), v(nloc, nd)
35     real gz(nloc, nd), h(nloc, nd), lv(nloc, nd), cpn(nloc, nd)
36     real p(nloc, nd), ph(nloc, nd+1), tv(nloc, nd), tp(nloc, nd)
37     real tvp(nloc, nd), clw(nloc, nd)
38     real th(nloc, nd)
39     real sig(nloc, nd), w0(nloc, nd)
40    
41     ! local variables:
42 guez 105 integer i, k, nn
43 guez 91
44    
45     do k=1, nl+1
46 guez 47 nn=0
47 guez 91 do i=1, len
48     if(iflag1(i).eq.0)then
49     nn=nn+1
50     sig(nn, k)=sig1(i, k)
51     w0(nn, k)=w01(i, k)
52     t(nn, k)=t1(i, k)
53     q(nn, k)=q1(i, k)
54     qs(nn, k)=qs1(i, k)
55     u(nn, k)=u1(i, k)
56     v(nn, k)=v1(i, k)
57     gz(nn, k)=gz1(i, k)
58     h(nn, k)=h1(i, k)
59     lv(nn, k)=lv1(i, k)
60     cpn(nn, k)=cpn1(i, k)
61     p(nn, k)=p1(i, k)
62     ph(nn, k)=ph1(i, k)
63     tv(nn, k)=tv1(i, k)
64     tp(nn, k)=tp1(i, k)
65     tvp(nn, k)=tvp1(i, k)
66     clw(nn, k)=clw1(i, k)
67     th(nn, k)=th1(i, k)
68     endif
69     end do
70     end do
71 guez 47
72 guez 91 if (nn.ne.ncum) then
73     print*, 'strange! nn not equal to ncum: ', nn, ncum
74     stop
75     endif
76 guez 47
77 guez 91 nn=0
78     do i=1, len
79     if(iflag1(i).eq.0)then
80     nn=nn+1
81     pbase(nn)=pbase1(i)
82     buoybase(nn)=buoybase1(i)
83     plcl(nn)=plcl1(i)
84     tnk(nn)=tnk1(i)
85     qnk(nn)=qnk1(i)
86     gznk(nn)=gznk1(i)
87     nk(nn)=nk1(i)
88     icb(nn)=icb1(i)
89     icbs(nn)=icbs1(i)
90     iflag(nn)=iflag1(i)
91     endif
92     end do
93 guez 47
94 guez 185 end SUBROUTINE cv30_compress
95 guez 91
96 guez 185 end module cv30_compress_m

  ViewVC Help
Powered by ViewVC 1.1.21