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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 288 - (hide annotations)
Tue Jul 24 16:27:12 2018 UTC (5 years, 11 months ago) by guez
File size: 3150 byte(s)
Remove tests on richum, tvirtu and opt_ec in procedure coefkz (not
used in LMDZ either).

Change the meaning of variable ecrit_ins of module clesphys from
number of seconds (integer, weird), which was modified in physiq to
make a multiple of the time step of physics, to number of physics time
steps.

1 guez 250 module coef_diff_turb_m
2    
3     implicit none
4    
5     contains
6    
7 guez 251 subroutine coef_diff_turb(dtime, nsrf, ni, paprs, pplay, u, v, q, t, ts, &
8     cdragm, zgeop, coefm, coefh, q2)
9 guez 250
10 guez 251 ! Computes coefficients for turbulent diffusion in the atmosphere.
11    
12 guez 288 use nr_util, only: assert
13    
14 guez 250 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 guez 251 USE dimphy, ONLY: klev
20 guez 250 USE suphec_m, ONLY: rd, rg, rkappa
21     use ustarhb_m, only: ustarhb
22     use yamada4_m, only: yamada4
23    
24     REAL, INTENT(IN):: dtime ! interval du temps (secondes)
25     INTEGER, INTENT(IN):: nsrf
26     INTEGER, INTENT(IN):: ni(:) ! (knon)
27 guez 251 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 guez 250 REAL, INTENT(IN):: zgeop(:, :) ! (knon, klev)
32 guez 251 REAL, intent(out):: coefm(:, 2:) ! (knon, 2:klev) coefficient, vitesse
33 guez 250
34 guez 251 real, intent(out):: coefh(:, 2:) ! (knon, 2:klev)
35 guez 250 ! coefficient, chaleur et humidité
36    
37 guez 251 real, intent(inout):: q2(:, :) ! (knon, klev + 1)
38 guez 250
39     ! Local:
40 guez 251 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 guez 250
45     !-------------------------------------------------------------------------
46    
47 guez 251 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 guez 250
54     IF (iflag_pbl == 1) THEN
55 guez 251 CALL coefkz2(nsrf, paprs, pplay, t, coefm0, coefh0)
56     coefm = max(coefm, coefm0)
57     coefh = max(coefh, coefh0)
58 guez 250 END IF
59    
60     IF (ok_kzmin) THEN
61     ! Calcul d'une diffusion minimale pour les conditions tres stables
62 guez 251 CALL coefkzmin(paprs, pplay, u, v, t, q, cdragm, coefh0)
63 guez 288 coefm = max(coefm, coefh0)
64 guez 251 coefh = max(coefh, coefh0)
65 guez 250 END IF
66    
67     IF (iflag_pbl >= 6) THEN
68     ! Mellor et Yamada adapt\'e \`a Mars, Richard Fournier et
69     ! Fr\'ed\'eric Hourdin
70 guez 251 zlay(:, 1) = rd * t(:, 1) / (0.5 * (paprs(:, 1) &
71     + pplay(:, 1))) * (paprs(:, 1) - pplay(:, 1)) / rg
72 guez 250
73     DO k = 2, klev
74 guez 251 zlay(:, k) = zlay(:, k-1) + rd * 0.5 * (t(:, k-1) + t(:, k)) &
75     / paprs(:, k) * (pplay(:, k-1) - pplay(:, k)) / rg
76 guez 250 END DO
77    
78 guez 288 forall (k = 1:klev) teta(:, k) = t(:, k) &
79     * (paprs(:, 1) / pplay(:, k))**rkappa * (1. + 0.61 * q(:, k))
80 guez 250
81 guez 251 zlev(:, 1) = 0.
82     zlev(:, klev + 1) = 2. * zlay(:, klev) - zlay(:, klev - 1)
83 guez 288 forall (k = 2:klev) zlev(:, k) = 0.5 * (zlay(:, k) + zlay(:, k-1))
84 guez 250
85 guez 288 CALL yamada4(dtime, zlev, zlay, u, v, teta, q2, coefm, coefh, &
86     ustarhb(u(:, 1), v(:, 1), cdragm))
87 guez 250 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