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

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

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

revision 195 by guez, Wed May 18 17:56:44 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(ncum, iflag1, nk1, icb1, icbs1, plcl1, tnk1, qnk1, &    SUBROUTINE cv30_compress(iflag1, icb1, icbs1, plcl1, tnk1, qnk1, gznk1, &
8         gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, cpn1, &         pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, cpn1, p1, &
9         p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, iflag, nk, icb, icbs, plcl, &         ph1, tv1, tp1, tvp1, clw1, sig1, w01, icb, icbs, plcl, tnk, qnk, gznk, &
10         tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, &         pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &
11         ph, tv, tp, tvp, clw, sig, w0)         tvp, clw, sig, w0)
12    
13      ! Compress the fields (vectorization over convective gridpoints).      ! Compress the fields (vectorization over convective gridpoints).
14    
# Line 17  contains Line 17  contains
17      use nr_util, only: assert      use nr_util, only: assert
18    
19      ! inputs:      ! inputs:
20      integer, intent(in):: ncum      integer, intent(in):: iflag1(:), icb1(:), icbs1(:) ! (klon)
     integer, intent(in):: iflag1(klon), nk1(klon), icb1(klon), icbs1(klon)  
21      real, intent(in):: plcl1(klon), tnk1(klon), qnk1(klon), gznk1(klon)      real, intent(in):: plcl1(klon), tnk1(klon), qnk1(klon), gznk1(klon)
22      real pbase1(klon), buoybase1(klon)      real pbase1(klon), buoybase1(klon)
23      real, intent(in):: t1(klon, klev)      real, intent(in):: t1(klon, klev)
24      real, intent(in):: q1(klon, klev), qs1(klon, klev)      real, intent(in):: q1(klon, klev), qs1(klon, klev)
25      real, intent(in):: u1(klon, klev), v1(klon, klev)      real, intent(in):: u1(klon, klev), v1(klon, klev)
26      real gz1(klon, klev), h1(klon, klev), lv1(klon, klev), cpn1(klon, klev)      real gz1(klon, klev), h1(klon, klev), lv1(klon, klev), cpn1(klon, klev)
27      real, intent(in):: p1(klon, klev), ph1(klon, klev+1)      real, intent(in):: p1(klon, klev), ph1(klon, klev + 1)
28      real, intent(in):: tv1(klon, klev), tp1(klon, klev)      real, intent(in):: tv1(klon, klev), tp1(klon, klev)
29      real tvp1(klon, klev), clw1(klon, klev)      real tvp1(klon, klev), clw1(klon, klev)
30      real th1(klon, klev)      real th1(klon, klev)
31      real sig1(klon, klev), w01(klon, klev)      real sig1(klon, klev), w01(klon, klev)
32    
33      ! outputs:      ! outputs:
34      integer iflag(klon), nk(klon)      integer, intent(out):: icb(:) ! (ncum) {2 <= icb <= nl - 3}
     integer, intent(out):: icb(:) ! (ncum)  
35      integer icbs(klon)      integer icbs(klon)
36      real plcl(klon), tnk(klon), qnk(klon), gznk(klon)      real, intent(out):: plcl(:) ! (ncum)
37        real tnk(:), qnk(:), gznk(:) ! (klon)
38      real pbase(klon), buoybase(klon)      real pbase(klon), buoybase(klon)
39      real t(klon, klev), q(klon, klev), qs(klon, klev)      real t(klon, klev), q(klon, klev), qs(klon, klev)
40      real u(klon, klev), v(klon, klev)      real u(klon, klev), v(klon, klev)
41      real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)      real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)
42      real p(klon, klev), ph(klon, klev+1), tv(klon, klev), tp(klon, klev)      real p(klon, klev)
43        real ph(:, :) ! (klon, klev + 1)
44        real tv(klon, klev), tp(klon, klev)
45      real tvp(klon, klev), clw(klon, klev)      real tvp(klon, klev), clw(klon, klev)
46      real th(klon, klev)      real th(klon, klev)
47      real sig(klon, klev), w0(klon, klev)      real sig(klon, klev), w0(klon, klev)
48    
49      ! Local:      ! Local:
50      integer i, k, nn      integer i, k, nn, ncum
51    
52      !---------------------------------------------------------------      !---------------------------------------------------------------
53    
54      do k=1, nl+1      ncum = size(icb)
55         nn=0  
56         do i=1, klon      do k = 1, nl + 1
57           nn = 0
58           do i = 1, klon
59            if (iflag1(i) == 0) then            if (iflag1(i) == 0) then
60               nn=nn+1               nn = nn + 1
61               sig(nn, k)=sig1(i, k)               sig(nn, k) = sig1(i, k)
62               w0(nn, k)=w01(i, k)               w0(nn, k) = w01(i, k)
63               t(nn, k)=t1(i, k)               t(nn, k) = t1(i, k)
64               q(nn, k)=q1(i, k)               q(nn, k) = q1(i, k)
65               qs(nn, k)=qs1(i, k)               qs(nn, k) = qs1(i, k)
66               u(nn, k)=u1(i, k)               u(nn, k) = u1(i, k)
67               v(nn, k)=v1(i, k)               v(nn, k) = v1(i, k)
68               gz(nn, k)=gz1(i, k)               gz(nn, k) = gz1(i, k)
69               h(nn, k)=h1(i, k)               h(nn, k) = h1(i, k)
70               lv(nn, k)=lv1(i, k)               lv(nn, k) = lv1(i, k)
71               cpn(nn, k)=cpn1(i, k)               cpn(nn, k) = cpn1(i, k)
72               p(nn, k)=p1(i, k)               p(nn, k) = p1(i, k)
73               ph(nn, k)=ph1(i, k)               ph(nn, k) = ph1(i, k)
74               tv(nn, k)=tv1(i, k)               tv(nn, k) = tv1(i, k)
75               tp(nn, k)=tp1(i, k)               tp(nn, k) = tp1(i, k)
76               tvp(nn, k)=tvp1(i, k)               tvp(nn, k) = tvp1(i, k)
77               clw(nn, k)=clw1(i, k)               clw(nn, k) = clw1(i, k)
78               th(nn, k)=th1(i, k)               th(nn, k) = th1(i, k)
79            endif            endif
80         end do         end do
81      end do      end do
82    
83      call assert(nn == ncum, "cv30_compress")      nn = 0
     nn=0  
84    
85      do i=1, klon      do i = 1, klon
86         if (iflag1(i) == 0) then         if (iflag1(i) == 0) then
87            nn=nn+1            nn = nn + 1
88            pbase(nn)=pbase1(i)            pbase(nn) = pbase1(i)
89            buoybase(nn)=buoybase1(i)            buoybase(nn) = buoybase1(i)
90            plcl(nn)=plcl1(i)            plcl(nn) = plcl1(i)
91            tnk(nn)=tnk1(i)            tnk(nn) = tnk1(i)
92            qnk(nn)=qnk1(i)            qnk(nn) = qnk1(i)
93            gznk(nn)=gznk1(i)            gznk(nn) = gznk1(i)
94            nk(nn)=nk1(i)            icb(nn) = icb1(i)
95            icb(nn)=icb1(i)            icbs(nn) = icbs1(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.195  
changed lines
  Added in v.198

  ViewVC Help
Powered by ViewVC 1.1.21