1 |
guez |
295 |
module climb_hq_down_m |
2 |
|
|
|
3 |
|
|
implicit none |
4 |
|
|
|
5 |
|
|
contains |
6 |
|
|
|
7 |
guez |
299 |
subroutine climb_hq_down(pkf, cq, dq, ch, dh, paprs, pplay, t, coef, delp, q) |
8 |
guez |
295 |
|
9 |
guez |
299 |
use comconst, only: dtphys |
10 |
guez |
295 |
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):: delp(:, :) ! (knon, klev) |
30 |
|
|
! epaisseur de couche en pression (Pa) |
31 |
|
|
|
32 |
|
|
REAL, intent(in):: q(:, :) ! (knon, klev) humidite specifique (kg / kg) |
33 |
|
|
|
34 |
|
|
! Local: |
35 |
|
|
|
36 |
|
|
INTEGER k |
37 |
|
|
REAL h(size(paprs, 1), klev) ! (knon, klev) enthalpie potentielle |
38 |
|
|
REAL zx_coef(size(paprs, 1), 2:klev) ! (knon, 2:klev) |
39 |
|
|
|
40 |
|
|
REAL gamt(size(paprs, 1), 2:klev) ! (knon, 2:klev) |
41 |
|
|
! contre-gradient pour la chaleur sensible, en K m-1 |
42 |
|
|
|
43 |
|
|
REAL gamah(size(paprs, 1), 2:klev) ! (knon, 2:klev) |
44 |
|
|
REAL buf1(size(paprs, 1)), buf2(size(paprs, 1)) |
45 |
|
|
|
46 |
|
|
!---------------------------------------------------------------- |
47 |
|
|
|
48 |
|
|
h = RCPD * t * pkf |
49 |
|
|
|
50 |
|
|
! Convertir les coefficients en variables convenables au calcul: |
51 |
|
|
forall (k = 2:klev) zx_coef(:, k) = coef(:, k) & |
52 |
|
|
/ (pplay(:, k - 1) - pplay(:, k)) & |
53 |
guez |
299 |
* (paprs(:, k) * 2 / (t(:, k) + t(:, k - 1)) / RD)**2 * dtphys * RG**2 |
54 |
guez |
295 |
|
55 |
|
|
! Preparer les flux lies aux contre-gardients |
56 |
|
|
|
57 |
|
|
if (iflag_pbl == 1) then |
58 |
|
|
gamt(:, 2) = - 2.5e-3 |
59 |
|
|
gamt(:, 3:)= - 1e-3 |
60 |
|
|
forall (k = 2:klev) gamah(:, k) = gamt(:, k) * (RD * (t(:, k - 1) & |
61 |
|
|
+ t(:, k)) / 2. / RG / paprs(:, k) * (pplay(:, k - 1) & |
62 |
|
|
- pplay(:, k))) * RCPD * (paprs(:, 1) / paprs(:, k))**RKAPPA |
63 |
|
|
else |
64 |
|
|
gamah = 0. |
65 |
|
|
endif |
66 |
|
|
|
67 |
|
|
buf1 = zx_coef(:, klev) + delp(:, klev) |
68 |
|
|
cq(:, klev) = q(:, klev) * delp(:, klev) / buf1 |
69 |
|
|
dq(:, klev) = zx_coef(:, klev) / buf1 |
70 |
|
|
|
71 |
|
|
buf2 = delp(:, klev) / pkf(:, klev) + zx_coef(:, klev) |
72 |
|
|
ch(:, klev) = (h(:, klev) / pkf(:, klev) * delp(:, klev) & |
73 |
|
|
- zx_coef(:, klev) * gamah(:, klev)) / buf2 |
74 |
|
|
dh(:, klev) = zx_coef(:, klev) / buf2 |
75 |
|
|
|
76 |
|
|
DO k = klev - 1, 2, - 1 |
77 |
|
|
buf1 = delp(:, k) + zx_coef(:, k) & |
78 |
|
|
+ zx_coef(:, k + 1) * (1. - dq(:, k + 1)) |
79 |
|
|
cq(:, k) = (q(:, k) * delp(:, k) & |
80 |
|
|
+ zx_coef(:, k + 1) * cq(:, k + 1)) / buf1 |
81 |
|
|
dq(:, k) = zx_coef(:, k) / buf1 |
82 |
|
|
|
83 |
|
|
buf2 = delp(:, k) / pkf(:, k) + zx_coef(:, k) & |
84 |
|
|
+ zx_coef(:, k + 1) * (1. - dh(:, k + 1)) |
85 |
|
|
ch(:, k) = (h(:, k) / pkf(:, k) * delp(:, k) & |
86 |
|
|
+ zx_coef(:, k + 1) * ch(:, k + 1) & |
87 |
|
|
+ zx_coef(:, k + 1) * gamah(:, k + 1) & |
88 |
|
|
- zx_coef(:, k) * gamah(:, k)) / buf2 |
89 |
|
|
dh(:, k) = zx_coef(:, k) / buf2 |
90 |
|
|
ENDDO |
91 |
|
|
|
92 |
|
|
buf1 = delp(:, 1) + zx_coef(:, 2) * (1. - dq(:, 2)) |
93 |
|
|
cq(:, 1) = (q(:, 1) * delp(:, 1) + zx_coef(:, 2) * cq(:, 2)) / buf1 |
94 |
|
|
dq(:, 1) = - 1. * RG / buf1 |
95 |
|
|
|
96 |
|
|
buf2 = delp(:, 1) / pkf(:, 1) + zx_coef(:, 2) * (1. - dh(:, 2)) |
97 |
|
|
ch(:, 1) = (h(:, 1) / pkf(:, 1) * delp(:, 1) & |
98 |
|
|
+ zx_coef(:, 2) * (gamah(:, 2) + ch(:, 2))) / buf2 |
99 |
|
|
dh(:, 1) = - 1. * RG / buf2 |
100 |
|
|
|
101 |
|
|
end subroutine climb_hq_down |
102 |
|
|
|
103 |
|
|
end module climb_hq_down_m |