/[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 198 by guez, Tue May 31 16:17:35 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(iflag1, icb1, icbs1, plcl1, tnk1, qnk1, gznk1, &
8         plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, &         pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, cpn1, p1, &
9         th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, iflag, &         ph1, tv1, tp1, tvp1, clw1, sig1, w01, icb, icbs, plcl, tnk, qnk, gznk, &
10         nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, &         pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &
11         gz, th, h, lv, cpn, p, ph, tv, tp, tvp, clw, sig, w0)         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):: iflag1(:), icb1(:), icbs1(:) ! (klon)
21      integer iflag1(len), nk1(len), icb1(len), icbs1(len)      real, intent(in):: plcl1(klon), tnk1(klon), qnk1(klon), gznk1(klon)
22      real plcl1(len), tnk1(len), qnk1(len), gznk1(len)      real pbase1(klon), buoybase1(klon)
23      real pbase1(len), buoybase1(len)      real, intent(in):: t1(klon, klev)
24      real, intent(in):: t1(len, nd)      real, intent(in):: q1(klon, klev), qs1(klon, klev)
25      real, intent(in):: q1(len, nd), qs1(len, nd), u1(len, nd), v1(len, nd)      real, intent(in):: u1(klon, klev), v1(klon, klev)
26      real gz1(len, nd), h1(len, nd), lv1(len, nd), cpn1(len, nd)      real gz1(klon, klev), h1(klon, klev), lv1(klon, klev), cpn1(klon, klev)
27      real, intent(in):: p1(len, nd), ph1(len, nd+1), tv1(len, nd), tp1(len, nd)      real, intent(in):: p1(klon, klev), ph1(klon, klev + 1)
28      real tvp1(len, nd), clw1(len, nd)      real, intent(in):: tv1(klon, klev), tp1(klon, klev)
29      real th1(len, nd)      real tvp1(klon, klev), clw1(klon, klev)
30      real sig1(len, nd), w01(len, nd)      real th1(klon, klev)
31        real sig1(klon, klev), w01(klon, klev)
32    
33      ! outputs:      ! outputs:
34      ! en fait, on a nloc=len pour l'instant (cf cv_driver)      integer, intent(out):: icb(:) ! (ncum) {2 <= icb <= nl - 3}
35      integer iflag(nloc), nk(nloc), icb(nloc), icbs(nloc)      integer icbs(klon)
36      real plcl(nloc), tnk(nloc), qnk(nloc), gznk(nloc)      real, intent(out):: plcl(:) ! (ncum)
37      real pbase(nloc), buoybase(nloc)      real tnk(:), qnk(:), gznk(:) ! (klon)
38      real t(nloc, nd), q(nloc, nd), qs(nloc, nd), u(nloc, nd), v(nloc, nd)      real pbase(klon), buoybase(klon)
39      real gz(nloc, nd), h(nloc, nd), lv(nloc, nd), cpn(nloc, nd)      real t(klon, klev), q(klon, klev), qs(klon, klev)
40      real p(nloc, nd), ph(nloc, nd+1), tv(nloc, nd), tp(nloc, nd)      real u(klon, klev), v(klon, klev)
41      real tvp(nloc, nd), clw(nloc, nd)      real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)
42      real th(nloc, nd)      real p(klon, klev)
43      real sig(nloc, nd), w0(nloc, nd)      real ph(:, :) ! (klon, klev + 1)
44        real tv(klon, klev), tp(klon, klev)
45      ! local variables:      real tvp(klon, klev), clw(klon, klev)
46      integer i, k, nn      real th(klon, klev)
47        real sig(klon, klev), w0(klon, klev)
48    
49      do  k=1, nl+1      ! Local:
50         nn=0      integer i, k, nn, ncum
51         do  i=1, len  
52            if(iflag1(i).eq.0)then      !---------------------------------------------------------------
53               nn=nn+1  
54               sig(nn, k)=sig1(i, k)      ncum = size(icb)
55               w0(nn, k)=w01(i, k)  
56               t(nn, k)=t1(i, k)      do k = 1, nl + 1
57               q(nn, k)=q1(i, k)         nn = 0
58               qs(nn, k)=qs1(i, k)         do i = 1, klon
59               u(nn, k)=u1(i, k)            if (iflag1(i) == 0) then
60               v(nn, k)=v1(i, k)               nn = nn + 1
61               gz(nn, k)=gz1(i, k)               sig(nn, k) = sig1(i, k)
62               h(nn, k)=h1(i, k)               w0(nn, k) = w01(i, k)
63               lv(nn, k)=lv1(i, k)               t(nn, k) = t1(i, k)
64               cpn(nn, k)=cpn1(i, k)               q(nn, k) = q1(i, k)
65               p(nn, k)=p1(i, k)               qs(nn, k) = qs1(i, k)
66               ph(nn, k)=ph1(i, k)               u(nn, k) = u1(i, k)
67               tv(nn, k)=tv1(i, k)               v(nn, k) = v1(i, k)
68               tp(nn, k)=tp1(i, k)               gz(nn, k) = gz1(i, k)
69               tvp(nn, k)=tvp1(i, k)               h(nn, k) = h1(i, k)
70               clw(nn, k)=clw1(i, k)               lv(nn, k) = lv1(i, k)
71               th(nn, k)=th1(i, k)               cpn(nn, k) = cpn1(i, k)
72                 p(nn, k) = p1(i, k)
73                 ph(nn, k) = ph1(i, k)
74                 tv(nn, k) = tv1(i, k)
75                 tp(nn, k) = tp1(i, k)
76                 tvp(nn, k) = tvp1(i, k)
77                 clw(nn, k) = clw1(i, k)
78                 th(nn, k) = th1(i, k)
79            endif            endif
80         end do         end do
81      end do      end do
82    
83      if (nn.ne.ncum) then      nn = 0
84         print*, 'strange! nn not equal to ncum: ', nn, ncum  
85         stop      do i = 1, klon
86      endif         if (iflag1(i) == 0) then
87              nn = nn + 1
88      nn=0            pbase(nn) = pbase1(i)
89      do  i=1, len            buoybase(nn) = buoybase1(i)
90         if(iflag1(i).eq.0)then            plcl(nn) = plcl1(i)
91            nn=nn+1            tnk(nn) = tnk1(i)
92            pbase(nn)=pbase1(i)            qnk(nn) = qnk1(i)
93            buoybase(nn)=buoybase1(i)            gznk(nn) = gznk1(i)
94            plcl(nn)=plcl1(i)            icb(nn) = icb1(i)
95            tnk(nn)=tnk1(i)            icbs(nn) = icbs1(i)
           qnk(nn)=qnk1(i)  
           gznk(nn)=gznk1(i)  
           nk(nn)=nk1(i)  
           icb(nn)=icb1(i)  
           icbs(nn)=icbs1(i)  
           iflag(nn)=iflag1(i)  
96         endif         endif
97      end do      end do
98    
99        do i = 1, ncum
100           call assert(2 <= icb(i) .and. icb(i) <= nl - 3 .and. ph(i, icb(i) + 1) &
101                < plcl(i) .and. (plcl(i) <= ph(i, icb(i)) .or. icb(i) == 2), &
102                "cv30_compress")
103        end do
104    
105    end SUBROUTINE cv30_compress    end SUBROUTINE cv30_compress
106    
107  end module cv30_compress_m  end module cv30_compress_m

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

  ViewVC Help
Powered by ViewVC 1.1.21