/[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 297 - (hide annotations)
Thu Jul 26 16:02:11 2018 UTC (5 years, 10 months ago) by guez
File size: 3529 byte(s)
Rename module interface_surf to conf_interface_m.

Move the computation of pkf out of procedure climb_hq_down into clqh.

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

  ViewVC Help
Powered by ViewVC 1.1.21