/[lmdze]/trunk/phylmd/Interface_surf/coef_diff_turb.f
ViewVC logotype

Contents of /trunk/phylmd/Interface_surf/coef_diff_turb.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 298 - (show annotations)
Thu Jul 26 16:45:51 2018 UTC (5 years, 10 months ago) by guez
File size: 3155 byte(s)
Use directly dtphys from module comconst when possible instead of
having it trickle down through procedure arguments.

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 comconst, only: dtphys
16 use coefkz_m, only: coefkz
17 use coefkzmin_m, only: coefkzmin
18 use coefkz2_m, only: coefkz2
19 USE conf_phys_m, ONLY: iflag_pbl
20 USE dimphy, ONLY: klev
21 USE suphec_m, ONLY: rd, rg, rkappa
22 use ustarhb_m, only: ustarhb
23 use yamada4_m, only: yamada4
24
25 INTEGER, INTENT(IN):: nsrf
26 INTEGER, INTENT(IN):: ni(:) ! (knon)
27 REAL, INTENT(IN):: paprs(:, :) ! (knon, klev + 1)
28 REAL, INTENT(IN):: pplay(:, :) ! (knon, klev)
29 REAL, INTENT(IN), dimension(:, :):: u, v, q, t ! (knon, klev)
30 REAL, INTENT(IN):: ts(:), cdragm(:) ! (knon)
31 REAL, INTENT(IN):: zgeop(:, :) ! (knon, klev)
32 REAL, intent(out):: coefm(:, 2:) ! (knon, 2:klev) coefficient, vitesse
33
34 real, intent(out):: coefh(:, 2:) ! (knon, 2:klev)
35 ! coefficient, chaleur et humidité
36
37 real, intent(inout):: q2(:, :) ! (knon, klev + 1)
38
39 ! Local:
40 REAL coefm0(size(ni), 2:klev), coefh0(size(ni), 2:klev) ! (knon, 2:klev)
41 REAL zlay(size(ni), klev), teta(size(ni), klev) ! (knon, klev)
42 real zlev(size(ni), klev + 1) ! (knon, klev + 1)
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 IF (iflag_pbl >= 6) THEN
53 ! Mellor et Yamada adapt\'e \`a Mars, Richard Fournier et
54 ! Fr\'ed\'eric Hourdin
55 zlay(:, 1) = rd * t(:, 1) / (0.5 * (paprs(:, 1) &
56 + pplay(:, 1))) * (paprs(:, 1) - pplay(:, 1)) / rg
57
58 DO k = 2, klev
59 zlay(:, k) = zlay(:, k-1) + rd * 0.5 * (t(:, k-1) + t(:, k)) &
60 / paprs(:, k) * (pplay(:, k-1) - pplay(:, k)) / rg
61 END DO
62
63 forall (k = 1:klev) teta(:, k) = t(:, k) &
64 * (paprs(:, 1) / pplay(:, k))**rkappa * (1. + 0.61 * q(:, k))
65
66 zlev(:, 1) = 0.
67 zlev(:, klev + 1) = 2. * zlay(:, klev) - zlay(:, klev - 1)
68 forall (k = 2:klev) zlev(:, k) = 0.5 * (zlay(:, k) + zlay(:, k-1))
69
70 CALL yamada4(dtphys, zlev, zlay, u, v, teta, q2, coefm, coefh, &
71 ustarhb(u(:, 1), v(:, 1), cdragm))
72 else
73 CALL coefkz(nsrf, paprs, pplay, ts, u, v, t, q, zgeop, coefm, coefh)
74
75 IF (iflag_pbl == 1) THEN
76 CALL coefkz2(nsrf, paprs, pplay, t, coefm0, coefh0)
77 coefm = max(coefm, coefm0)
78 coefh = max(coefh, coefh0)
79 END IF
80
81 IF (ok_kzmin) THEN
82 ! Calcul d'une diffusion minimale pour les conditions tres stables
83 CALL coefkzmin(paprs, pplay, u, v, t, q, cdragm, coefh0)
84 coefm = max(coefm, coefh0)
85 coefh = max(coefh, coefh0)
86 END IF
87 END IF
88
89 end subroutine coef_diff_turb
90
91 end module coef_diff_turb_m

  ViewVC Help
Powered by ViewVC 1.1.21