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

Diff of /trunk/Sources/phylmd/CV30_routines/cv30_compress.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 185 by guez, Wed Mar 16 15:04:46 2016 UTC revision 195 by guez, Wed May 18 17:56:44 2016 UTC
# Line 4  module cv30_compress_m Line 4  module cv30_compress_m
4    
5  contains  contains
6    
7    SUBROUTINE cv30_compress(len, nloc, ncum, nd, iflag1, nk1, icb1, icbs1, &    SUBROUTINE cv30_compress(ncum, iflag1, nk1, icb1, icbs1, plcl1, tnk1, qnk1, &
8         plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, &         gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, cpn1, &
9         th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, iflag, &         p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, iflag, nk, icb, icbs, plcl, &
10         nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, &         tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, &
11         gz, th, h, lv, cpn, p, ph, tv, tp, tvp, clw, sig, w0)         ph, tv, tp, tvp, clw, sig, w0)
12    
13      use cv30_param_m      ! Compress the fields (vectorization over convective gridpoints).
14    
15        use cv30_param_m, only: nl
16        USE dimphy, ONLY: klev, klon
17        use nr_util, only: assert
18    
19      ! inputs:      ! inputs:
20      integer, intent(in):: len, ncum, nd, nloc      integer, intent(in):: ncum
21      integer iflag1(len), nk1(len), icb1(len), icbs1(len)      integer, intent(in):: iflag1(klon), nk1(klon), icb1(klon), icbs1(klon)
22      real plcl1(len), tnk1(len), qnk1(len), gznk1(len)      real, intent(in):: plcl1(klon), tnk1(klon), qnk1(klon), gznk1(klon)
23      real pbase1(len), buoybase1(len)      real pbase1(klon), buoybase1(klon)
24      real, intent(in):: t1(len, nd)      real, intent(in):: t1(klon, klev)
25      real, intent(in):: q1(len, nd), qs1(len, nd), u1(len, nd), v1(len, nd)      real, intent(in):: q1(klon, klev), qs1(klon, klev)
26      real gz1(len, nd), h1(len, nd), lv1(len, nd), cpn1(len, nd)      real, intent(in):: u1(klon, klev), v1(klon, klev)
27      real, intent(in):: p1(len, nd), ph1(len, nd+1), tv1(len, nd), tp1(len, nd)      real gz1(klon, klev), h1(klon, klev), lv1(klon, klev), cpn1(klon, klev)
28      real tvp1(len, nd), clw1(len, nd)      real, intent(in):: p1(klon, klev), ph1(klon, klev+1)
29      real th1(len, nd)      real, intent(in):: tv1(klon, klev), tp1(klon, klev)
30      real sig1(len, nd), w01(len, nd)      real tvp1(klon, klev), clw1(klon, klev)
31        real th1(klon, klev)
32        real sig1(klon, klev), w01(klon, klev)
33    
34      ! outputs:      ! outputs:
35      ! en fait, on a nloc=len pour l'instant (cf cv_driver)      integer iflag(klon), nk(klon)
36      integer iflag(nloc), nk(nloc), icb(nloc), icbs(nloc)      integer, intent(out):: icb(:) ! (ncum)
37      real plcl(nloc), tnk(nloc), qnk(nloc), gznk(nloc)      integer icbs(klon)
38      real pbase(nloc), buoybase(nloc)      real plcl(klon), tnk(klon), qnk(klon), gznk(klon)
39      real t(nloc, nd), q(nloc, nd), qs(nloc, nd), u(nloc, nd), v(nloc, nd)      real pbase(klon), buoybase(klon)
40      real gz(nloc, nd), h(nloc, nd), lv(nloc, nd), cpn(nloc, nd)      real t(klon, klev), q(klon, klev), qs(klon, klev)
41      real p(nloc, nd), ph(nloc, nd+1), tv(nloc, nd), tp(nloc, nd)      real u(klon, klev), v(klon, klev)
42      real tvp(nloc, nd), clw(nloc, nd)      real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)
43      real th(nloc, nd)      real p(klon, klev), ph(klon, klev+1), tv(klon, klev), tp(klon, klev)
44      real sig(nloc, nd), w0(nloc, nd)      real tvp(klon, klev), clw(klon, klev)
45        real th(klon, klev)
46        real sig(klon, klev), w0(klon, klev)
47    
48      ! local variables:      ! Local:
49      integer i, k, nn      integer i, k, nn
50    
51        !---------------------------------------------------------------
52    
53      do  k=1, nl+1      do k=1, nl+1
54         nn=0         nn=0
55         do  i=1, len         do i=1, klon
56            if(iflag1(i).eq.0)then            if (iflag1(i) == 0) then
57               nn=nn+1               nn=nn+1
58               sig(nn, k)=sig1(i, k)               sig(nn, k)=sig1(i, k)
59               w0(nn, k)=w01(i, k)               w0(nn, k)=w01(i, k)
# Line 69  contains Line 77  contains
77         end do         end do
78      end do      end do
79    
80      if (nn.ne.ncum) then      call assert(nn == ncum, "cv30_compress")
        print*, 'strange! nn not equal to ncum: ', nn, ncum  
        stop  
     endif  
   
81      nn=0      nn=0
82      do  i=1, len  
83         if(iflag1(i).eq.0)then      do i=1, klon
84           if (iflag1(i) == 0) then
85            nn=nn+1            nn=nn+1
86            pbase(nn)=pbase1(i)            pbase(nn)=pbase1(i)
87            buoybase(nn)=buoybase1(i)            buoybase(nn)=buoybase1(i)

Legend:
Removed from v.185  
changed lines
  Added in v.195

  ViewVC Help
Powered by ViewVC 1.1.21