/[lmdze]/trunk/phylmd/Thermcell/calltherm.f90
ViewVC logotype

Annotation of /trunk/phylmd/Thermcell/calltherm.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 324 - (hide annotations)
Wed Feb 6 15:58:03 2019 UTC (5 years, 4 months ago) by guez
Original Path: trunk/phylmd/Thermcell/calltherm.f
File size: 3027 byte(s)
Rename variable zmasq of module phyetat0_m to masque, which was
already its name in "restartphy.nc". Rename variable fraclic of
procedure etat0 to landice, which was already its name in
"landiceref.nc". Style guide: we try to have the same names for
identical data objects across the program.

In procedure interfsurf_hq, in case is_sic, define tsurf instead of
tsurf_new, avoiding a copy from tsurf_new to tsurf.

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

  ViewVC Help
Powered by ViewVC 1.1.21