1 |
module coef_diff_turb_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
contains |
6 |
|
7 |
subroutine coef_diff_turb(dtime, nsrf, ni, ypaprs, ypplay, yu, yv, yq, yt, & |
8 |
yts, ycdragm, zgeop, ycoefm, ycoefh, yq2) |
9 |
|
10 |
USE clesphys, ONLY: ok_kzmin |
11 |
use coefkz_m, only: coefkz |
12 |
use coefkzmin_m, only: coefkzmin |
13 |
use coefkz2_m, only: coefkz2 |
14 |
USE conf_phys_m, ONLY: iflag_pbl |
15 |
USE dimphy, ONLY: klev, klon |
16 |
USE suphec_m, ONLY: rd, rg, rkappa |
17 |
use ustarhb_m, only: ustarhb |
18 |
use yamada4_m, only: yamada4 |
19 |
|
20 |
REAL, INTENT(IN):: dtime ! interval du temps (secondes) |
21 |
INTEGER, INTENT(IN):: nsrf |
22 |
INTEGER, INTENT(IN):: ni(:) ! (knon) |
23 |
REAL, INTENT(IN):: ypaprs(:, :) ! (klon, klev + 1) |
24 |
REAL, INTENT(IN):: ypplay(:, :) ! (klon, klev) |
25 |
REAL, INTENT(IN), dimension(:, :):: yu, yv, yq, yt ! (klon, klev) |
26 |
REAL, INTENT(IN):: yts(:), ycdragm(:) ! (knon) |
27 |
REAL, INTENT(IN):: zgeop(:, :) ! (knon, klev) |
28 |
REAL, intent(out):: ycoefm(:, 2:) ! (knon, 2:klev) coefficient, vitesse |
29 |
|
30 |
real, intent(out):: ycoefh(:, 2:) ! (knon, 2:klev) |
31 |
! coefficient, chaleur et humidité |
32 |
|
33 |
real, intent(inout):: yq2(:, :) ! (klon, klev + 1) |
34 |
|
35 |
! Local: |
36 |
REAL ycoefm0(klon, 2:klev), ycoefh0(klon, 2:klev) |
37 |
REAL yzlay(klon, klev), zlev(klon, klev + 1), yteta(klon, klev) |
38 |
REAL ustar(klon) |
39 |
integer k, knon |
40 |
|
41 |
!------------------------------------------------------------------------- |
42 |
|
43 |
knon = size(ni) |
44 |
CALL coefkz(nsrf, ypaprs(:knon, :), ypplay(:knon, :), yts(:knon), & |
45 |
yu(:knon, :), yv(:knon, :), yt(:knon, :), yq(:knon, :), zgeop, & |
46 |
ycoefm, ycoefh) |
47 |
|
48 |
IF (iflag_pbl == 1) THEN |
49 |
CALL coefkz2(nsrf, knon, ypaprs, ypplay, yt, ycoefm0(:knon, :), & |
50 |
ycoefh0(:knon, :)) |
51 |
ycoefm = max(ycoefm, ycoefm0(:knon, :)) |
52 |
ycoefh = max(ycoefh, ycoefh0(:knon, :)) |
53 |
END IF |
54 |
|
55 |
IF (ok_kzmin) THEN |
56 |
! Calcul d'une diffusion minimale pour les conditions tres stables |
57 |
CALL coefkzmin(knon, ypaprs, ypplay, yu, yv, yt, yq, ycdragm(:knon), & |
58 |
ycoefh0(:knon, :)) |
59 |
ycoefm0(:knon, :) = ycoefh0(:knon, :) |
60 |
ycoefm = max(ycoefm, ycoefm0(:knon, :)) |
61 |
ycoefh = max(ycoefh, ycoefh0(:knon, :)) |
62 |
END IF |
63 |
|
64 |
IF (iflag_pbl >= 6) THEN |
65 |
! Mellor et Yamada adapt\'e \`a Mars, Richard Fournier et |
66 |
! Fr\'ed\'eric Hourdin |
67 |
yzlay(:knon, 1) = rd * yt(:knon, 1) / (0.5 * (ypaprs(:knon, 1) & |
68 |
+ ypplay(:knon, 1))) & |
69 |
* (ypaprs(:knon, 1) - ypplay(:knon, 1)) / rg |
70 |
|
71 |
DO k = 2, klev |
72 |
yzlay(:knon, k) = yzlay(:knon, k-1) & |
73 |
+ rd * 0.5 * (yt(1:knon, k-1) + yt(1:knon, k)) & |
74 |
/ ypaprs(1:knon, k) & |
75 |
* (ypplay(1:knon, k-1) - ypplay(1:knon, k)) / rg |
76 |
END DO |
77 |
|
78 |
DO k = 1, klev |
79 |
yteta(1:knon, k) = yt(1:knon, k) * (ypaprs(1:knon, 1) & |
80 |
/ ypplay(1:knon, k))**rkappa * (1. + 0.61 * yq(1:knon, k)) |
81 |
END DO |
82 |
|
83 |
zlev(:knon, 1) = 0. |
84 |
zlev(:knon, klev + 1) = 2. * yzlay(:knon, klev) & |
85 |
- yzlay(:knon, klev - 1) |
86 |
|
87 |
DO k = 2, klev |
88 |
zlev(:knon, k) = 0.5 * (yzlay(:knon, k) + yzlay(:knon, k-1)) |
89 |
END DO |
90 |
|
91 |
ustar(:knon) = ustarhb(yu(:knon, 1), yv(:knon, 1), ycdragm(:knon)) |
92 |
CALL yamada4(dtime, rg, zlev(:knon, :), yzlay(:knon, :), & |
93 |
yu(:knon, :), yv(:knon, :), yteta(:knon, :), yq2(:knon, :), & |
94 |
ycoefm, ycoefh, ustar(:knon)) |
95 |
END IF |
96 |
|
97 |
end subroutine coef_diff_turb |
98 |
|
99 |
end module coef_diff_turb_m |