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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 328 - (show annotations)
Thu Jun 13 14:40:06 2019 UTC (5 years ago) by guez
File size: 3027 byte(s)
Change all `.f` suffixes to `.f90`. (The opposite was done in revision
82.)  Because of change of philosopy in GNUmakefile: we already had a
rewritten rule for `.f`, so it does not make the makefile longer to
replace it by a rule for `.f90`. And it spares us options of
makedepf90 and of the compiler. Also we prepare the way for a simpler
`CMakeLists.txt`.

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_m, ONLY: nsplit_thermals
16 USE dimphy, ONLY: klev, klon
17 use thermcell_m, only: thermcell
18
19 REAL, intent(in):: pplay(klon, klev)
20 REAL, intent(in):: paprs(klon, klev+1)
21 REAL, intent(in):: pphi(klon, klev)
22 REAL, intent(inout):: u_seri(klon, klev), v_seri(klon, klev)
23 REAL, intent(inout):: t_seri(klon, klev)
24 real, intent(inout):: q_seri(klon, klev)
25
26 ! Update thermiques
27 REAL d_u_ajs(klon, klev), d_v_ajs(klon, klev)
28 REAL d_t_ajs(klon, klev), d_q_ajs(klon, klev)
29 real fm_therm(klon, klev+1), entr_therm(klon, klev)
30
31 ! Local:
32
33 REAL d_t_the(klon, klev), d_q_the(klon, klev)
34 REAL d_u_the(klon, klev), d_v_the(klon, klev)
35
36 real, save:: zfm_therm(klon, klev+1), zentr_therm(klon, klev)
37 real zdt
38
39 integer i, k, isplit
40
41 !----------------------------------------------------------------
42
43 ! Modele du thermique
44 print*, 'avant isplit ', nsplit_thermals
45
46 fm_therm=0.
47 entr_therm=0.
48
49 ! tests sur les valeurs negatives de l'eau
50 do k=1, klev
51 do i=1, klon
52 if (.not.q_seri(i, k) >= 0.) then
53 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
60 zdt=dtphys/real(nsplit_thermals)
61 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 zentr_therm)
65
66 ! transformation de la derivee en tendance
67 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 fm_therm=fm_therm +zfm_therm/real(nsplit_thermals)
72 entr_therm=entr_therm +zentr_therm/real(nsplit_thermals)
73 fm_therm(:, klev+1)=0.
74
75 ! 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
81 ! 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
87 ! tests sur les valeurs negatives de l'eau
88 DO k = 1, klev
89 DO i = 1, klon
90 if (.not.q_seri(i, k) >= 0.) then
91 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
99 end subroutine calltherm
100
101 end module calltherm_m

  ViewVC Help
Powered by ViewVC 1.1.21