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

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

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

trunk/Sources/phylmd/CV3_routines/cv3_prelim.f revision 134 by guez, Wed Apr 29 15:47:56 2015 UTC trunk/Sources/phylmd/CV30_routines/cv30_prelim.f revision 198 by guez, Tue May 31 16:17:35 2016 UTC
# Line 1  Line 1 
1  module cv3_prelim_m  module cv30_prelim_m
2    
3    implicit none    implicit none
4    
5  contains  contains
6    
7    SUBROUTINE cv3_prelim(len, nd, ndp1, t, q, p, ph, lv, cpn, tv, gz, h, hm, th)    SUBROUTINE cv30_prelim(t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, h1, hm1, th1)
8    
9      USE cv3_param_m, ONLY: nl      USE cv30_param_m, ONLY: nl
10      USE cvthermo, ONLY: cl, clmcpv, cpd, cpv, eps, lv0, rrd, rrv      USE cv_thermo_m, ONLY: cl, clmcpv, cpd, cpv, eps, rrd, rrv
11        USE dimphy, ONLY: klev, klon
12        use SUPHEC_M, only: rlvtt
13    
14      ! Calculate arrays of geopotential, heat capacity and static energy      ! Calculate arrays of geopotential, heat capacity and static energy
15    
16      integer, intent(in):: len, nd, ndp1      real, intent(in):: t1(klon, klev)
17      real, intent(in):: t(len, nd)      real, intent(in):: q1(klon, klev)
18      real, intent(in):: q(len, nd)      real, intent(in):: p1(klon, klev), ph1(klon, klev + 1)
     real, intent(in):: p(len, nd), ph(len, ndp1)  
19    
20      ! outputs:      ! outputs:
21      real lv(len, nd), cpn(len, nd), tv(len, nd)      real lv1(klon, klev), cpn1(klon, klev), tv1(klon, klev)
22      real gz(len, nd), h(len, nd), hm(len, nd)      real gz1(klon, klev), h1(klon, klev), hm1(klon, klev)
23      real th(len, nd)      real th1(klon, klev) ! potential temperature
24    
25      ! Local:      ! Local:
26      integer k, i      integer k, i
27      real rdcp      real rdcp
28      real tvx, tvy      real tvx, tvy
29      real cpx(len, nd)      real cpx(klon, klev)
30    
31      !--------------------------------------------------------------      !--------------------------------------------------------------
32    
33      do k=1, nl      do k = 1, nl
34         do i=1, len         do i = 1, klon
35            lv(i, k)= lv0-clmcpv*(t(i, k)-273.15)            lv1(i, k) =  rlvtt - clmcpv * (t1(i, k) - 273.15)
36            cpn(i, k)=cpd*(1.0-q(i, k)) + cpv*q(i, k)            cpn1(i, k) = cpd * (1. - q1(i, k)) + cpv * q1(i, k)
37            cpx(i, k)=cpd*(1.0-q(i, k)) + cl*q(i, k)            cpx(i, k) = cpd * (1. - q1(i, k)) + cl * q1(i, k)
38            tv(i, k)=t(i, k)*(1.0 + q(i, k)/eps-q(i, k))            tv1(i, k) = t1(i, k) * (1. + q1(i, k)/eps - q1(i, k))
39            rdcp=(rrd*(1.-q(i, k)) + q(i, k)*rrv)/cpn(i, k)            rdcp = (rrd * (1. - q1(i, k)) + q1(i, k) * rrv)/cpn1(i, k)
40            th(i, k)=t(i, k)*(1000.0/p(i, k))**rdcp            th1(i, k) = t1(i, k) * (1000./p1(i, k))**rdcp
41         end do         end do
42      end do      end do
43    
44      ! gz = phi at the full levels (same as p).      ! gz1 = phi at the full levels (same as p1).
45    
46      do i=1, len      do i = 1, klon
47         gz(i, 1)=0.0         gz1(i, 1) = 0.
48      end do      end do
49    
50      do k=2, nl      do k = 2, nl
51         do i=1, len         do i = 1, klon
52            tvx=t(i, k)*(1. + q(i, k)/eps-q(i, k))            tvx = t1(i, k) * (1. + q1(i, k)/eps - q1(i, k))
53            tvy=t(i, k-1)*(1. + q(i, k-1)/eps-q(i, k-1))            tvy = t1(i, k - 1) * (1. + q1(i, k - 1)/eps - q1(i, k - 1))
54            gz(i, k)=gz(i, k-1) + 0.5*rrd*(tvx + tvy) &            gz1(i, k) = gz1(i, k - 1) + 0.5 * rrd * (tvx + tvy) &
55                 *(p(i, k-1)-p(i, k))/ph(i, k)                 * (p1(i, k - 1) - p1(i, k))/ph1(i, k)
56         end do         end do
57      end do      end do
58    
59      ! h = phi + cpT (dry static energy).      ! h1 = phi + cpT (dry static energy).
60      ! hm = phi + cp(T-Tbase) + Lq      ! hm1 = phi + cp(T1 - Tbase) + Lq
61    
62      do k=1, nl      do k = 1, nl
63         do i=1, len         do i = 1, klon
64            h(i, k)=gz(i, k) + cpn(i, k)*t(i, k)            h1(i, k) = gz1(i, k) + cpn1(i, k) * t1(i, k)
65            hm(i, k)=gz(i, k) + cpx(i, k)*(t(i, k)-t(i, 1)) + lv(i, k)*q(i, k)            hm1(i, k) = gz1(i, k) + cpx(i, k) * (t1(i, k) - t1(i, 1)) &
66                   + lv1(i, k) * q1(i, k)
67         end do         end do
68      end do      end do
69    
70    end SUBROUTINE cv3_prelim    end SUBROUTINE cv30_prelim
71    
72  end module cv3_prelim_m  end module cv30_prelim_m

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

  ViewVC Help
Powered by ViewVC 1.1.21