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

Annotation of /trunk/Sources/phylmd/calltherm.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 52 - (hide annotations)
Fri Sep 23 12:28:01 2011 UTC (12 years, 8 months ago) by guez
Original Path: trunk/libf/phylmd/calltherm.f90
File size: 3111 byte(s)
Split "conflx.f" into single-procedure files in directory "Conflx".

Split "cv_routines.f" into single-procedure files in directory
"CV_routines". Made module "cvparam" from included file
"cvparam.h". No included file other than "netcdf.inc" left in LMDZE.

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

  ViewVC Help
Powered by ViewVC 1.1.21