/[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 189 by guez, Tue Mar 29 15:20:23 2016 UTC revision 197 by guez, Tue May 24 12:25:29 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, nk1, icb1, icbs1, plcl1, tnk1, qnk1, &
8         gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, cpn1, &         gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, cpn1, &
9         p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, iflag, nk, icb, icbs, plcl, &         p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, nk, icb, icbs, plcl, tnk, &
10         tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, &         qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, &
11         ph, tv, tp, tvp, clw, sig, w0)         tv, tp, tvp, clw, sig, w0)
12    
13      ! Compress the fields (vectorization over convective gridpoints).      ! Compress the fields (vectorization over convective gridpoints).
14    
15      use cv30_param_m, only: nl      use cv30_param_m, only: nl
16      USE dimphy, ONLY: klev, klon      USE dimphy, ONLY: klev, klon
17        use nr_util, only: assert
18    
19      ! inputs:      ! inputs:
20      integer, intent(in):: ncum      integer, intent(in):: iflag1(:), nk1(:), icb1(:), icbs1(:) ! (klon)
21      integer iflag1(klon), nk1(klon), icb1(klon), icbs1(klon)      real, intent(in):: plcl1(klon), tnk1(klon), qnk1(klon), gznk1(klon)
     real 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), icb(klon), icbs(klon)      integer nk(:) ! (klon)
35      real plcl(klon), tnk(klon), qnk(klon), gznk(klon)      integer, intent(out):: icb(:) ! (ncum) {2 <= icb <= nl - 3}
36        integer icbs(klon)
37        real, intent(out):: plcl(:) ! (ncum)
38        real tnk(:), qnk(:), gznk(:) ! (klon)
39      real pbase(klon), buoybase(klon)      real pbase(klon), buoybase(klon)
40      real t(klon, klev), q(klon, klev), qs(klon, klev)      real t(klon, klev), q(klon, klev), qs(klon, klev)
41      real u(klon, klev), v(klon, klev)      real u(klon, klev), v(klon, klev)
42      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)
43      real p(klon, klev), ph(klon, klev+1), tv(klon, klev), tp(klon, klev)      real p(klon, klev)
44        real ph(:, :) ! (klon, klev + 1)
45        real tv(klon, klev), tp(klon, klev)
46      real tvp(klon, klev), clw(klon, klev)      real tvp(klon, klev), clw(klon, klev)
47      real th(klon, klev)      real th(klon, klev)
48      real sig(klon, klev), w0(klon, klev)      real sig(klon, klev), w0(klon, klev)
49    
50      ! Local:      ! Local:
51      integer i, k, nn      integer i, k, nn, ncum
52    
53      !---------------------------------------------------------------      !---------------------------------------------------------------
54    
55      do k=1, nl+1      ncum = size(icb)
56         nn=0  
57         do i=1, klon      do k = 1, nl + 1
58           nn = 0
59           do i = 1, klon
60            if (iflag1(i) == 0) then            if (iflag1(i) == 0) then
61               nn=nn+1               nn = nn + 1
62               sig(nn, k)=sig1(i, k)               sig(nn, k) = sig1(i, k)
63               w0(nn, k)=w01(i, k)               w0(nn, k) = w01(i, k)
64               t(nn, k)=t1(i, k)               t(nn, k) = t1(i, k)
65               q(nn, k)=q1(i, k)               q(nn, k) = q1(i, k)
66               qs(nn, k)=qs1(i, k)               qs(nn, k) = qs1(i, k)
67               u(nn, k)=u1(i, k)               u(nn, k) = u1(i, k)
68               v(nn, k)=v1(i, k)               v(nn, k) = v1(i, k)
69               gz(nn, k)=gz1(i, k)               gz(nn, k) = gz1(i, k)
70               h(nn, k)=h1(i, k)               h(nn, k) = h1(i, k)
71               lv(nn, k)=lv1(i, k)               lv(nn, k) = lv1(i, k)
72               cpn(nn, k)=cpn1(i, k)               cpn(nn, k) = cpn1(i, k)
73               p(nn, k)=p1(i, k)               p(nn, k) = p1(i, k)
74               ph(nn, k)=ph1(i, k)               ph(nn, k) = ph1(i, k)
75               tv(nn, k)=tv1(i, k)               tv(nn, k) = tv1(i, k)
76               tp(nn, k)=tp1(i, k)               tp(nn, k) = tp1(i, k)
77               tvp(nn, k)=tvp1(i, k)               tvp(nn, k) = tvp1(i, k)
78               clw(nn, k)=clw1(i, k)               clw(nn, k) = clw1(i, k)
79               th(nn, k)=th1(i, k)               th(nn, k) = th1(i, k)
80            endif            endif
81         end do         end do
82      end do      end do
83    
84      if (nn /= ncum) then      nn = 0
        print*, 'strange! nn not equal to ncum: ', nn, ncum  
        stop 1  
     endif  
85    
86      nn=0      do i = 1, klon
     do i=1, klon  
87         if (iflag1(i) == 0) then         if (iflag1(i) == 0) then
88            nn=nn+1            nn = nn + 1
89            pbase(nn)=pbase1(i)            pbase(nn) = pbase1(i)
90            buoybase(nn)=buoybase1(i)            buoybase(nn) = buoybase1(i)
91            plcl(nn)=plcl1(i)            plcl(nn) = plcl1(i)
92            tnk(nn)=tnk1(i)            tnk(nn) = tnk1(i)
93            qnk(nn)=qnk1(i)            qnk(nn) = qnk1(i)
94            gznk(nn)=gznk1(i)            gznk(nn) = gznk1(i)
95            nk(nn)=nk1(i)            nk(nn) = nk1(i)
96            icb(nn)=icb1(i)            icb(nn) = icb1(i)
97            icbs(nn)=icbs1(i)            icbs(nn) = icbs1(i)
           iflag(nn)=iflag1(i)  
98         endif         endif
99      end do      end do
100    
101        do i = 1, ncum
102           call assert(2 <= icb(i) .and. icb(i) <= nl - 3 .and. ph(i, icb(i) + 1) &
103                < plcl(i) .and. (plcl(i) <= ph(i, icb(i)) .or. icb(i) == 2), &
104                "cv30_compress")
105        end do
106    
107    end SUBROUTINE cv30_compress    end SUBROUTINE cv30_compress
108    
109  end module cv30_compress_m  end module cv30_compress_m

Legend:
Removed from v.189  
changed lines
  Added in v.197

  ViewVC Help
Powered by ViewVC 1.1.21