/[lmdze]/trunk/phylmd/Interface_surf/climb_hq_down.f
ViewVC logotype

Annotation of /trunk/phylmd/Interface_surf/climb_hq_down.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 295 - (hide annotations)
Thu Jul 26 13:23:28 2018 UTC (5 years, 10 months ago) by guez
File size: 3614 byte(s)
Create procedure climb_hq_down from part of procedure clqh (following LMDZ).

1 guez 295 module climb_hq_down_m
2    
3     implicit none
4    
5     contains
6    
7     subroutine climb_hq_down(pkf, cq, dq, ch, dh, paprs, pplay, t, coef, dtime, &
8     delp, q)
9    
10     USE conf_phys_m, ONLY: iflag_pbl
11     USE dimphy, ONLY: klev
12     USE suphec_m, ONLY: rcpd, rd, rg, rkappa
13    
14     REAL, intent(out), dimension(:, :):: pkf, cq, dq, ch, dh ! (knon, klev)
15    
16     REAL, intent(in):: paprs(:, :) ! (knon, klev + 1)
17     ! pression a inter-couche (Pa)
18    
19     REAL, intent(in):: pplay(:, :) ! (knon, klev)
20     ! pression au milieu de couche (Pa)
21    
22     REAL, intent(in):: t(:, :) ! (knon, klev) temperature (K)
23    
24     REAL, intent(in):: coef(:, 2:) ! (knon, 2:klev)
25     ! Le coefficient d'echange (m**2 / s) multiplie par le cisaillement
26     ! du vent (dV / dz)
27    
28     REAL, intent(in):: dtime ! intervalle du temps (s)
29    
30     REAL, intent(in):: delp(:, :) ! (knon, klev)
31     ! epaisseur de couche en pression (Pa)
32    
33     REAL, intent(in):: q(:, :) ! (knon, klev) humidite specifique (kg / kg)
34    
35     ! Local:
36    
37     INTEGER k
38     REAL h(size(paprs, 1), klev) ! (knon, klev) enthalpie potentielle
39     REAL zx_coef(size(paprs, 1), 2:klev) ! (knon, 2:klev)
40    
41     REAL gamt(size(paprs, 1), 2:klev) ! (knon, 2:klev)
42     ! contre-gradient pour la chaleur sensible, en K m-1
43    
44     REAL gamah(size(paprs, 1), 2:klev) ! (knon, 2:klev)
45     REAL buf1(size(paprs, 1)), buf2(size(paprs, 1))
46    
47     !----------------------------------------------------------------
48    
49     forall (k = 1:klev) pkf(:, k) = (paprs(:, 1) / pplay(:, k))**RKAPPA
50     ! (La pression de r\'ef\'erence est celle au sol.)
51    
52     h = RCPD * t * pkf
53    
54     ! Convertir les coefficients en variables convenables au calcul:
55     forall (k = 2:klev) zx_coef(:, k) = coef(:, k) &
56     / (pplay(:, k - 1) - pplay(:, k)) &
57     * (paprs(:, k) * 2 / (t(:, k) + t(:, k - 1)) / RD)**2 * dtime * RG**2
58    
59     ! Preparer les flux lies aux contre-gardients
60    
61     if (iflag_pbl == 1) then
62     gamt(:, 2) = - 2.5e-3
63     gamt(:, 3:)= - 1e-3
64     forall (k = 2:klev) gamah(:, k) = gamt(:, k) * (RD * (t(:, k - 1) &
65     + t(:, k)) / 2. / RG / paprs(:, k) * (pplay(:, k - 1) &
66     - pplay(:, k))) * RCPD * (paprs(:, 1) / paprs(:, k))**RKAPPA
67     else
68     gamah = 0.
69     endif
70    
71     buf1 = zx_coef(:, klev) + delp(:, klev)
72     cq(:, klev) = q(:, klev) * delp(:, klev) / buf1
73     dq(:, klev) = zx_coef(:, klev) / buf1
74    
75     buf2 = delp(:, klev) / pkf(:, klev) + zx_coef(:, klev)
76     ch(:, klev) = (h(:, klev) / pkf(:, klev) * delp(:, klev) &
77     - zx_coef(:, klev) * gamah(:, klev)) / buf2
78     dh(:, klev) = zx_coef(:, klev) / buf2
79    
80     DO k = klev - 1, 2, - 1
81     buf1 = delp(:, k) + zx_coef(:, k) &
82     + zx_coef(:, k + 1) * (1. - dq(:, k + 1))
83     cq(:, k) = (q(:, k) * delp(:, k) &
84     + zx_coef(:, k + 1) * cq(:, k + 1)) / buf1
85     dq(:, k) = zx_coef(:, k) / buf1
86    
87     buf2 = delp(:, k) / pkf(:, k) + zx_coef(:, k) &
88     + zx_coef(:, k + 1) * (1. - dh(:, k + 1))
89     ch(:, k) = (h(:, k) / pkf(:, k) * delp(:, k) &
90     + zx_coef(:, k + 1) * ch(:, k + 1) &
91     + zx_coef(:, k + 1) * gamah(:, k + 1) &
92     - zx_coef(:, k) * gamah(:, k)) / buf2
93     dh(:, k) = zx_coef(:, k) / buf2
94     ENDDO
95    
96     buf1 = delp(:, 1) + zx_coef(:, 2) * (1. - dq(:, 2))
97     cq(:, 1) = (q(:, 1) * delp(:, 1) + zx_coef(:, 2) * cq(:, 2)) / buf1
98     dq(:, 1) = - 1. * RG / buf1
99    
100     buf2 = delp(:, 1) / pkf(:, 1) + zx_coef(:, 2) * (1. - dh(:, 2))
101     ch(:, 1) = (h(:, 1) / pkf(:, 1) * delp(:, 1) &
102     + zx_coef(:, 2) * (gamah(:, 2) + ch(:, 2))) / buf2
103     dh(:, 1) = - 1. * RG / buf2
104    
105     end subroutine climb_hq_down
106    
107     end module climb_hq_down_m

  ViewVC Help
Powered by ViewVC 1.1.21