1 |
module coef_diff_turb_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
contains |
6 |
|
7 |
subroutine coef_diff_turb(nsrf, ni, paprs, pplay, u, v, q, t, ts, cdragm, & |
8 |
zgeop, coefm, coefh, q2) |
9 |
|
10 |
! Computes coefficients for turbulent diffusion in the atmosphere. |
11 |
|
12 |
use nr_util, only: assert |
13 |
|
14 |
USE clesphys, ONLY: ok_kzmin |
15 |
use coefkz_m, only: coefkz |
16 |
use coefkzmin_m, only: coefkzmin |
17 |
use coefkz2_m, only: coefkz2 |
18 |
USE conf_phys_m, ONLY: iflag_pbl |
19 |
USE dimphy, ONLY: klev |
20 |
USE suphec_m, ONLY: rd, rg, rkappa |
21 |
use ustarhb_m, only: ustarhb |
22 |
use yamada4_m, only: yamada4 |
23 |
|
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 |
integer k |
43 |
|
44 |
!------------------------------------------------------------------------- |
45 |
|
46 |
call assert(size(ni) == [size(paprs, 1), size(pplay, 1), size(u, 1), & |
47 |
size(v, 1), size(q, 1), size(t, 1), size(ts), size(cdragm), & |
48 |
size(zgeop, 1), size(coefm, 1), size(coefh, 1), size(q2, 1)], & |
49 |
"coef_diff_turb knon") |
50 |
|
51 |
IF (iflag_pbl >= 6) THEN |
52 |
! Mellor et Yamada adapt\'e \`a Mars, Richard Fournier et |
53 |
! Fr\'ed\'eric Hourdin |
54 |
zlay(:, 1) = rd * t(:, 1) / (0.5 * (paprs(:, 1) & |
55 |
+ pplay(:, 1))) * (paprs(:, 1) - pplay(:, 1)) / rg |
56 |
|
57 |
DO k = 2, klev |
58 |
zlay(:, k) = zlay(:, k-1) + rd * 0.5 * (t(:, k-1) + t(:, k)) & |
59 |
/ paprs(:, k) * (pplay(:, k-1) - pplay(:, k)) / rg |
60 |
END DO |
61 |
|
62 |
forall (k = 1:klev) teta(:, k) = t(:, k) & |
63 |
* (paprs(:, 1) / pplay(:, k))**rkappa * (1. + 0.61 * q(:, k)) |
64 |
|
65 |
zlev(:, 1) = 0. |
66 |
forall (k = 2:klev) zlev(:, k) = 0.5 * (zlay(:, k) + zlay(:, k-1)) |
67 |
zlev(:, klev + 1) = 2. * zlay(:, klev) - zlev(:, klev) |
68 |
|
69 |
CALL yamada4(zlev, zlay, u, v, teta, q2, coefm, coefh, & |
70 |
ustarhb(u(:, 1), v(:, 1), cdragm)) |
71 |
else |
72 |
CALL coefkz(nsrf, paprs, pplay, ts, u, v, t, q, zgeop, coefm, coefh) |
73 |
|
74 |
IF (iflag_pbl == 1) THEN |
75 |
CALL coefkz2(nsrf, paprs, pplay, t, coefm0, coefh0) |
76 |
coefm = max(coefm, coefm0) |
77 |
coefh = max(coefh, coefh0) |
78 |
END IF |
79 |
|
80 |
IF (ok_kzmin) THEN |
81 |
! Calcul d'une diffusion minimale pour les conditions tres stables |
82 |
CALL coefkzmin(paprs, pplay, u, v, t, q, cdragm, coefh0) |
83 |
coefm = max(coefm, coefh0) |
84 |
coefh = max(coefh, coefh0) |
85 |
END IF |
86 |
END IF |
87 |
|
88 |
end subroutine coef_diff_turb |
89 |
|
90 |
end module coef_diff_turb_m |