/[lmdze]/trunk/phylmd/calltherm.f
ViewVC logotype

Contents of /trunk/phylmd/calltherm.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 322 - (show annotations)
Thu Jan 24 16:17:39 2019 UTC (5 years, 3 months ago) by guez
File size: 3178 byte(s)
Remove unused constants in module cv_thermo.
1 module calltherm_m
2
3 implicit none
4
5 contains
6
7 subroutine calltherm(pplay, paprs, pphi, u_seri, v_seri, t_seri, q_seri, &
8 d_u_ajs, d_v_ajs, d_t_ajs, d_q_ajs, fm_therm, entr_therm)
9
10 ! From LMDZ4/libf/phylmd/calltherm.F, version 1.2 2004/12/10 11:27:46
11
12 ! Thermiques.
13
14 use comconst, only: dtphys
15 USE ctherm, ONLY: l_mix_thermals, nsplit_thermals, r_aspect_thermals, &
16 tho_thermals, w2di_thermals
17 USE dimphy, ONLY: klev, klon
18 use thermcell_m, only: thermcell
19
20 REAL, intent(in):: pplay(klon, klev)
21 REAL, intent(in):: paprs(klon, klev+1)
22 REAL, intent(in):: pphi(klon, klev)
23 REAL, intent(inout):: u_seri(klon, klev), v_seri(klon, klev)
24 REAL, intent(inout):: t_seri(klon, klev)
25 real, intent(inout):: q_seri(klon, klev)
26
27 ! Update thermiques
28 REAL d_u_ajs(klon, klev), d_v_ajs(klon, klev)
29 REAL d_t_ajs(klon, klev), d_q_ajs(klon, klev)
30 real fm_therm(klon, klev+1), entr_therm(klon, klev)
31
32 ! Local:
33
34 REAL d_t_the(klon, klev), d_q_the(klon, klev)
35 REAL d_u_the(klon, klev), d_v_the(klon, klev)
36
37 real, save:: zfm_therm(klon, klev+1), zentr_therm(klon, klev)
38 real zdt
39
40 integer i, k, isplit
41
42 !----------------------------------------------------------------
43
44 ! Modele du thermique
45 print*, 'avant isplit ', nsplit_thermals
46
47 fm_therm=0.
48 entr_therm=0.
49
50 ! tests sur les valeurs negatives de l'eau
51 do k=1, klev
52 do i=1, klon
53 if (.not.q_seri(i, k) >= 0.) then
54 print*, 'WARN eau<0 avant therm i=', i, ' k=', k, ' dq, q', &
55 d_q_the(i, k), q_seri(i, k)
56 q_seri(i, k)=1.e-15
57 endif
58 enddo
59 enddo
60
61 zdt=dtphys/real(nsplit_thermals)
62 do isplit = 1, nsplit_thermals
63 CALL thermcell(klon, klev, zdt, pplay, paprs, pphi, u_seri, v_seri, &
64 t_seri, q_seri, d_u_the, d_v_the, d_t_the, d_q_the, zfm_therm, &
65 zentr_therm, r_aspect_thermals, l_mix_thermals, w2di_thermals, &
66 tho_thermals)
67
68 ! transformation de la derivee en tendance
69 d_t_the=d_t_the*dtphys/real(nsplit_thermals)
70 d_u_the=d_u_the*dtphys/real(nsplit_thermals)
71 d_v_the=d_v_the*dtphys/real(nsplit_thermals)
72 d_q_the=d_q_the*dtphys/real(nsplit_thermals)
73 fm_therm=fm_therm +zfm_therm/real(nsplit_thermals)
74 entr_therm=entr_therm +zentr_therm/real(nsplit_thermals)
75 fm_therm(:, klev+1)=0.
76
77 ! accumulation de la tendance
78 d_t_ajs=d_t_ajs+d_t_the
79 d_u_ajs=d_u_ajs+d_u_the
80 d_v_ajs=d_v_ajs+d_v_the
81 d_q_ajs=d_q_ajs+d_q_the
82
83 ! incrementation des variables meteo
84 t_seri = t_seri + d_t_the
85 u_seri = u_seri + d_u_the
86 v_seri = v_seri + d_v_the
87 q_seri = q_seri + d_q_the
88
89 ! tests sur les valeurs negatives de l'eau
90 DO k = 1, klev
91 DO i = 1, klon
92 if (.not.q_seri(i, k) >= 0.) then
93 print*, 'WARN eau<0 apres therm i=', i, ' k=', k, ' dq, q', &
94 d_q_the(i, k), q_seri(i, k)
95 q_seri(i, k)=1.e-15
96 endif
97 ENDDO
98 ENDDO
99 enddo
100
101 end subroutine calltherm
102
103 end module calltherm_m

  ViewVC Help
Powered by ViewVC 1.1.21