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

Contents of /trunk/libf/phylmd/calltherm.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 47 - (show annotations)
Fri Jul 1 15:00:48 2011 UTC (12 years, 10 months ago) by guez
File size: 3445 byte(s)
Split "thermcell.f" and "cv3_routines.f".
Removed copies of files that are now in "L_util".
Moved "mva9" and "diagetpq" to their own files.
Unified variable names across procedures.

1 !
2 ! $Header: /home/cvsroot/LMDZ4/libf/phylmd/calltherm.F,v 1.2 2004/12/10 11:27:46 lmdzadmin Exp $
3 !
4 subroutine calltherm(dtime
5 s ,pplay,paprs,pphi
6 s ,u_seri,v_seri,t_seri,q_seri
7 s ,d_u_ajs,d_v_ajs,d_t_ajs,d_q_ajs
8 s ,fm_therm,entr_therm)
9
10 use dimens_m
11 use dimphy
12 use ctherm
13 implicit none
14
15 REAL, intent(in):: dtime
16
17 REAL u_seri(klon,klev),v_seri(klon,klev)
18 REAL t_seri(klon,klev),q_seri(klon,klev)
19 REAL, intent(in):: paprs(klon,klev+1)
20 REAL, intent(in):: pplay(klon,klev)
21 REAL, intent(in):: pphi(klon,klev)
22
23 CFH Update Thermiques
24 REAL d_t_ajs(klon,klev), d_q_ajs(klon,klev)
25 REAL d_u_ajs(klon,klev),d_v_ajs(klon,klev)
26 real fm_therm(klon,klev+1),entr_therm(klon,klev)
27
28
29 c variables locales
30 REAL d_t_the(klon,klev), d_q_the(klon,klev)
31 REAL d_u_the(klon,klev),d_v_the(klon,klev)
32 c
33 real zfm_therm(klon,klev+1),zentr_therm(klon,klev),zdt
34 save zentr_therm,zfm_therm
35
36 integer i,k, isplit
37
38 *********************************************************
39
40 c Modele du thermique
41 c ===================
42 c print*,'thermiques: WARNING on passe t au lieu de t_seri'
43 print*,'avant isplit ',nsplit_thermals
44
45
46 fm_therm(:,:)=0.
47 entr_therm(:,:)=0.
48
49 c 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).ge.0.) then
53 print*,'WARN eau<0 avant therm i=',i,' k=',k
54 s ,' dq,q',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
61 zdt=dtime/float(nsplit_thermals)
62 do isplit=1,nsplit_thermals
63
64 CALL thermcell(klon,klev,zdt
65 s ,pplay,paprs,pphi
66 s ,u_seri,v_seri,t_seri,q_seri
67 s ,d_u_the,d_v_the,d_t_the,d_q_the
68 s ,zfm_therm,zentr_therm
69 s ,r_aspect_thermals,l_mix_thermals,w2di_thermals
70 s ,tho_thermals)
71
72 c transformation de la derivee en tendance
73 d_t_the(:,:)=d_t_the(:,:)*dtime/float(nsplit_thermals)
74 d_u_the(:,:)=d_u_the(:,:)*dtime/float(nsplit_thermals)
75 d_v_the(:,:)=d_v_the(:,:)*dtime/float(nsplit_thermals)
76 d_q_the(:,:)=d_q_the(:,:)*dtime/float(nsplit_thermals)
77 fm_therm(:,:)=fm_therm(:,:)
78 s +zfm_therm(:,:)/float(nsplit_thermals)
79 entr_therm(:,:)=entr_therm(:,:)
80 s +zentr_therm(:,:)/float(nsplit_thermals)
81 fm_therm(:,klev+1)=0.
82
83
84
85 c accumulation de la tendance
86 d_t_ajs(:,:)=d_t_ajs(:,:)+d_t_the(:,:)
87 d_u_ajs(:,:)=d_u_ajs(:,:)+d_u_the(:,:)
88 d_v_ajs(:,:)=d_v_ajs(:,:)+d_v_the(:,:)
89 d_q_ajs(:,:)=d_q_ajs(:,:)+d_q_the(:,:)
90
91 c incrementation des variables meteo
92 t_seri(:,:) = t_seri(:,:) + d_t_the(:,:)
93 u_seri(:,:) = u_seri(:,:) + d_u_the(:,:)
94 v_seri(:,:) = v_seri(:,:) + d_v_the(:,:)
95 q_seri(:,:) = q_seri(:,:) + d_q_the(:,:)
96
97 c tests sur les valeurs negatives de l'eau
98 DO k = 1, klev
99 DO i = 1, klon
100 if (.not.q_seri(i,k).ge.0.) then
101 print*,'WARN eau<0 apres therm i=',i,' k=',k
102 s ,' dq,q',d_q_the(i,k),q_seri(i,k)
103 q_seri(i,k)=1.e-15
104 endif
105 ENDDO
106 ENDDO
107
108 enddo ! isplit
109
110 return
111
112 end

  ViewVC Help
Powered by ViewVC 1.1.21