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