1 | ! |
---|
2 | ! $Header: /home/cvsroot/LMDZ4/libf/phylmd/phytrac.F,v 1.16 2006/03/24 15:06:23 lmdzadmin Exp $ |
---|
3 | ! |
---|
4 | c |
---|
5 | c |
---|
6 | SUBROUTINE phytrac_chimie ( |
---|
7 | I debutphy, |
---|
8 | I gmtime, |
---|
9 | I nqmax, |
---|
10 | I n_lon, |
---|
11 | I lat, |
---|
12 | I lon, |
---|
13 | I n_lev, |
---|
14 | I pdtphys, |
---|
15 | I temp, |
---|
16 | I pplev, |
---|
17 | O trac) |
---|
18 | c====================================================================== |
---|
19 | c Auteur(s) FH |
---|
20 | c Objet: Moniteur general des tendances traceurs |
---|
21 | c |
---|
22 | cAA Remarques en vrac: |
---|
23 | cAA-------------------- |
---|
24 | cAA 1/ le call phytrac se fait avec nqmax |
---|
25 | c====================================================================== |
---|
26 | USE chemparam_mod |
---|
27 | use conc, only: mmean |
---|
28 | IMPLICIT none |
---|
29 | |
---|
30 | #include "clesphys.h" |
---|
31 | c#include "temps.h" |
---|
32 | c#include "paramet.h" |
---|
33 | c#include "comcstfi.h" !me permet de recuperer mugaz et d'autres constantes comme rad,pi etc |
---|
34 | #include "YOMCST.h" |
---|
35 | c====================================================================== |
---|
36 | c Arguments: |
---|
37 | c |
---|
38 | c |
---|
39 | c EN ENTREE: |
---|
40 | c ========== |
---|
41 | c |
---|
42 | c divers: |
---|
43 | c ------- |
---|
44 | c |
---|
45 | |
---|
46 | REAL sza_local |
---|
47 | REAL gmtime |
---|
48 | c INTEGER, SAVE :: cpt_cloudIO !un compteur pour fichier sortie cloud_parameter en 1D |
---|
49 | INTEGER iq |
---|
50 | INTEGER i |
---|
51 | INTEGER ilon, ilev |
---|
52 | integer n_lon ! nombre de points horizontaux |
---|
53 | INTEGER n_lev ! nombre de couches verticales |
---|
54 | INTEGER nqmax ! nombre de traceurs auxquels on applique la physique |
---|
55 | |
---|
56 | real pdtphys ! pas d'integration pour la physique (seconde) |
---|
57 | real lat(n_lon), lat_local(n_lon) |
---|
58 | real lon(n_lon), lon_local(n_lon) |
---|
59 | real temp(n_lon,n_lev) ! temp |
---|
60 | real trac(n_lon,n_lev,nqmax) ! traceur |
---|
61 | real trac_sav(n_lon,n_lev,nqmax) |
---|
62 | real trac_sum(n_lon,n_lev) |
---|
63 | real pplev(n_lon,n_lev) ! pression pour le mileu de chaque couche (en Pa) |
---|
64 | real lon_sun |
---|
65 | |
---|
66 | logical debutphy ! le flag de l'initialisation de la physique |
---|
67 | |
---|
68 | C Auxilary variables: |
---|
69 | |
---|
70 | REAL, DIMENSION(n_lon,n_lev) :: mrtwv,mrtsa, |
---|
71 | + mrwv,mrsa |
---|
72 | |
---|
73 | C ps_sa: satur pressure pure SA |
---|
74 | C satps_sa: satur pres over mixture in dyne/cm2=Pa/10 |
---|
75 | C---------------------------------------------------------------------------- |
---|
76 | |
---|
77 | C Variables liees a l'ecriture de la bande histoire : phytrac.nc |
---|
78 | |
---|
79 | logical ok_sync |
---|
80 | parameter (ok_sync = .true.) |
---|
81 | |
---|
82 | c modname = 'phytrac' |
---|
83 | |
---|
84 | c PRINT*,'DEBUT subroutine PHYTRAC' |
---|
85 | |
---|
86 | c---------------------------------- |
---|
87 | c debutphy: Initiation des traceurs |
---|
88 | c---------------------------------- |
---|
89 | |
---|
90 | if (debutphy) then |
---|
91 | |
---|
92 | if (n_lon .EQ. 1) then |
---|
93 | PRINT*,'n_lon 1D: ',n_lon |
---|
94 | end if |
---|
95 | |
---|
96 | c if ((n_lon .GT. 1) .AND. ok_chem) then |
---|
97 | c !!! DONC 3D !!! |
---|
98 | c CALL chemparam_ini() |
---|
99 | c endif |
---|
100 | |
---|
101 | c if ((n_lon .GT. 1) .AND. ok_cloud) then |
---|
102 | c !!! DONC 3D !!! |
---|
103 | c CALL cloud_ini(n_lon,n_lev) |
---|
104 | c endif |
---|
105 | |
---|
106 | IF (reinit_trac) THEN |
---|
107 | PRINT*,'REINIT MIXING RATIO TRACEURS' |
---|
108 | |
---|
109 | c ============================================================= |
---|
110 | c Passage de Rm à Rv |
---|
111 | c ============================================================= |
---|
112 | c Necessaire si on reprend les start.nc qui sont en MMR |
---|
113 | |
---|
114 | DO iq=1,nqmax |
---|
115 | trac(:,:,iq)=trac(:,:,iq)*mmean(:,:)/M_tr(iq) |
---|
116 | END DO |
---|
117 | c ============================================================= |
---|
118 | |
---|
119 | c============================================================= |
---|
120 | c Initialisation des profils traceurs en Rv |
---|
121 | c============================================================= |
---|
122 | c initialisation sert a mettre les valeurs voulues par utilisateur pour |
---|
123 | c chaque traceur |
---|
124 | c exemple: trac(ilon,ilev,q)=xx |
---|
125 | |
---|
126 | c trac_sav sert a sauver les valeurs initiales du start.nc |
---|
127 | trac_sav=trac |
---|
128 | |
---|
129 | c On initialise les traceurs a zero obligatoire pour la chimie |
---|
130 | trac(:,:,:)=1.0E-30 |
---|
131 | |
---|
132 | c !!! AVEC NUAGE !!! |
---|
133 | trac(:,1:22,i_ocs)=3.E-6 |
---|
134 | trac(:,:,i_hcl)=0.4E-6 |
---|
135 | trac(:,1:22,i_so2)=10.E-6 |
---|
136 | trac(:,1:22,i_h2o)=30.0E-6 |
---|
137 | |
---|
138 | c remettre tous les traceurs du start => trac(:,:,:)=trac_sav(:,:,:) |
---|
139 | |
---|
140 | c N2 n est pas encore une espece chimique du modele chimique |
---|
141 | c traceur passif pour la chimie-transport |
---|
142 | trac(:,:,i_n2)=0.35d-1 |
---|
143 | |
---|
144 | !!!! GG: Initialization CO2 = 1 - qtot |
---|
145 | !! It assures that vmr_tot = 1 |
---|
146 | c On a donc le CO2 qui est le restant d atmosphere Venus |
---|
147 | trac_sum(:,:)=0.0 |
---|
148 | DO iq=2,nqmax |
---|
149 | trac_sum(:,:)= trac_sum(:,:) + trac(:,:,iq) |
---|
150 | END DO |
---|
151 | |
---|
152 | trac(:,:,i_co2)= 1-trac_sum(:,:) |
---|
153 | |
---|
154 | c============================================================= |
---|
155 | |
---|
156 | c ============================================================= |
---|
157 | c Passage de Rv à Rm |
---|
158 | c ============================================================= |
---|
159 | DO iq=1,nqmax |
---|
160 | trac(:,:,iq)=trac(:,:,iq)*M_tr(iq)/mmean(:,:) |
---|
161 | END DO |
---|
162 | c ============================================================= |
---|
163 | |
---|
164 | ENDIF !FIN REINIT TRAC |
---|
165 | |
---|
166 | |
---|
167 | c------------- |
---|
168 | c fin debutphy |
---|
169 | c------------- |
---|
170 | |
---|
171 | ENDIF ! fin debutphy |
---|
172 | |
---|
173 | c ============================================================= |
---|
174 | c Passage de Rm à Rv |
---|
175 | c ============================================================= |
---|
176 | DO iq=1,nqmax |
---|
177 | trac(:,:,iq)=MAX(trac(:,:,iq)*mmean(:,:)/M_tr(iq),1.E-30) |
---|
178 | END DO |
---|
179 | c ============================================================= |
---|
180 | |
---|
181 | |
---|
182 | c ============================================================= |
---|
183 | c Appel Microphysique (sans nucleation) |
---|
184 | c Volume Mixing Ratio |
---|
185 | c ============================================================= |
---|
186 | |
---|
187 | IF (ok_cloud) THEN |
---|
188 | |
---|
189 | c PRINT*,'DEBUT CLOUD' |
---|
190 | c On remet tout le RM liq dans la partie gaz |
---|
191 | c !!! On reforme un nuage à chaque fois !!! |
---|
192 | |
---|
193 | DO ilev=1, n_lev |
---|
194 | DO ilon=1, n_lon |
---|
195 | mrtwv(ilon,ilev)=trac(ilon,ilev,i_h2o) + |
---|
196 | & trac(ilon,ilev,i_h2oliq) |
---|
197 | mrtsa(ilon,ilev)=trac(ilon,ilev,i_h2so4) + |
---|
198 | & trac(ilon,ilev,i_h2so4liq) |
---|
199 | mrwv(ilon,ilev)=mrtwv(ilon,ilev) |
---|
200 | mrsa(ilon,ilev)=mrtsa(ilon,ilev) |
---|
201 | ENDDO |
---|
202 | ENDDO |
---|
203 | |
---|
204 | CALL new_cloud_venus(n_lev, n_lon, |
---|
205 | e temp,pplev, |
---|
206 | e mrtwv,mrtsa, |
---|
207 | e mrwv,mrsa) |
---|
208 | |
---|
209 | c ========================================= |
---|
210 | c Actualisation des mixing ratio liq et gaz |
---|
211 | c ========================================= |
---|
212 | c Si la routine new_cloud_venus n'a pas actualisé mrwv et mrsa |
---|
213 | c on a alors bien mr=mrt pour sa et wv, donc les parties liq sont=0 hors du nuage |
---|
214 | c ou si on ne condense pas |
---|
215 | |
---|
216 | c PRINT*,'DEBUT ACTUALISATION OUTPUT CLOUD' |
---|
217 | c si tout se passe bien, mrtwv et mrtsa ne changent pas |
---|
218 | DO ilev=1, n_lev |
---|
219 | DO ilon=1, n_lon |
---|
220 | trac(ilon,ilev,i_h2o) = mrwv(ilon,ilev) |
---|
221 | trac(ilon,ilev,i_h2oliq) = mrtwv(ilon,ilev) - |
---|
222 | & trac(ilon,ilev,i_h2o) |
---|
223 | |
---|
224 | trac(ilon,ilev,i_h2so4) = mrsa(ilon,ilev) |
---|
225 | trac(ilon,ilev,i_h2so4liq) = mrtsa(ilon,ilev) - |
---|
226 | & trac(ilon,ilev,i_h2so4) |
---|
227 | ENDDO |
---|
228 | ENDDO |
---|
229 | |
---|
230 | c ============================================================= |
---|
231 | c PRINT*,'FIN CLOUD' |
---|
232 | ENDIF |
---|
233 | |
---|
234 | c============================================================= |
---|
235 | c CHIMIE: Boucle sur les lon, lat (n_lon) |
---|
236 | c============================================================= |
---|
237 | |
---|
238 | c AS: |
---|
239 | c Ici, la longitude au midi local se deplace vers l'Ouest |
---|
240 | c c'est le sens terrestre |
---|
241 | c pour Vénus on prend juste l'opposé de la longitude et on a la rotation |
---|
242 | c de Vénus et donc le midi local qui se déplace vers l'Est |
---|
243 | |
---|
244 | lon_sun = (0.5 - gmtime) * 2.0 * RPI |
---|
245 | lon_local = lon * RPI/180.0E+0 |
---|
246 | lat_local = lat * RPI/180.0E+0 |
---|
247 | |
---|
248 | DO ilon=1, n_lon |
---|
249 | |
---|
250 | c calcul sza_local pour obtenir des sza_local > 90, utile pour la chimie |
---|
251 | sza_local = acos(cos(lat_local(ilon))*cos(lon_local(ilon))* |
---|
252 | & cos(lon_sun) + cos(lat_local(ilon))*sin(lon_local(ilon)) |
---|
253 | & *sin(lon_sun))* 180.0E+0/RPI |
---|
254 | |
---|
255 | c PRINT*,'sza_local :', sza_local |
---|
256 | |
---|
257 | IF (ok_chem) THEN |
---|
258 | c PRINT*,'DEBUT CHEMISTRY' |
---|
259 | c ============================================================= |
---|
260 | c Appel Photochimie |
---|
261 | c ============================================================= |
---|
262 | c Pression en hPa => pplev/100. |
---|
263 | |
---|
264 | CALL new_photochemistry_venus(n_lev, n_lon, pdtphys, |
---|
265 | e pplev(ilon,:)/100., |
---|
266 | e temp(ilon,:), |
---|
267 | e trac(ilon,:,:), |
---|
268 | e mmean(ilon,:), |
---|
269 | e sza_local, nqmax) |
---|
270 | c ============================================================= |
---|
271 | c PRINT*,'FIN CHEMISTRY' |
---|
272 | |
---|
273 | END IF |
---|
274 | |
---|
275 | END DO |
---|
276 | c ============================================================= |
---|
277 | c Passage de Rv à Rm |
---|
278 | c ============================================================= |
---|
279 | DO iq=1,nqmax |
---|
280 | c trac(:,:,iq)=trac(:,:,iq)*M_tr(iq)/RMD |
---|
281 | trac(:,:,iq)=trac(:,:,iq)*M_tr(iq)/mmean(:,:) |
---|
282 | |
---|
283 | END DO |
---|
284 | c ============================================================= |
---|
285 | C PRINT*,'FIN PHYTRAC' |
---|
286 | RETURN |
---|
287 | END |
---|