1 | ! $Id: infotrac.F90 3945 2021-06-18 15:24:43Z acozic $ |
---|
2 | ! |
---|
3 | MODULE infotrac |
---|
4 | |
---|
5 | ! nqtot : total number of tracers and higher order of moment, water vapor and liquid included |
---|
6 | INTEGER, SAVE :: nqtot |
---|
7 | !CR: on ajoute le nombre de traceurs de l eau |
---|
8 | INTEGER, SAVE :: nqo |
---|
9 | |
---|
10 | ! nbtr : number of tracers not including higher order of moment or water vapor or liquid |
---|
11 | ! number of tracers used in the physics |
---|
12 | INTEGER, SAVE :: nbtr |
---|
13 | |
---|
14 | ! CRisi: on retranche les isotopes des traceurs habituels |
---|
15 | ! On fait un tableaux d'indices des traceurs qui passeront dans phytrac |
---|
16 | INTEGER, SAVE :: nqtottr |
---|
17 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: itr_indice |
---|
18 | |
---|
19 | ! CRisi: nb traceurs peres= directement advectes par l'air |
---|
20 | INTEGER, SAVE :: nqperes |
---|
21 | |
---|
22 | ! ThL: nb traceurs INCA |
---|
23 | INTEGER, SAVE :: nqINCA |
---|
24 | |
---|
25 | ! ThL: nb traceurs CO2 |
---|
26 | INTEGER, SAVE :: nqCO2 |
---|
27 | |
---|
28 | ! Name variables |
---|
29 | CHARACTER(len=20), ALLOCATABLE, DIMENSION(:), SAVE :: tname ! tracer short name for restart and diagnostics |
---|
30 | CHARACTER(len=23), ALLOCATABLE, DIMENSION(:), SAVE :: ttext ! tracer long name for diagnostics |
---|
31 | |
---|
32 | ! iadv : index of trasport schema for each tracer |
---|
33 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: iadv |
---|
34 | |
---|
35 | ! niadv : vector keeping the coorspondance between all tracers(nqtot) treated in the |
---|
36 | ! dynamic part of the code and the tracers (nbtr+2) used in the physics part of the code. |
---|
37 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: niadv ! equivalent dyn / physique |
---|
38 | |
---|
39 | ! CRisi: tableaux de fils |
---|
40 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: nqfils |
---|
41 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: nqdesc ! nombres de fils + nombre de tous les petits fils sur toutes les generations |
---|
42 | INTEGER, SAVE :: nqdesc_tot |
---|
43 | INTEGER, ALLOCATABLE, DIMENSION(:,:), SAVE :: iqfils |
---|
44 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: iqpere |
---|
45 | REAL :: qperemin,masseqmin,ratiomin ! MVals et CRisi |
---|
46 | PARAMETER (qperemin=1e-16,masseqmin=1e-16,ratiomin=1e-16) ! MVals |
---|
47 | |
---|
48 | ! conv_flg(it)=0 : convection desactivated for tracer number it |
---|
49 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: conv_flg |
---|
50 | ! pbl_flg(it)=0 : boundary layer diffusion desactivaded for tracer number it |
---|
51 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: pbl_flg |
---|
52 | |
---|
53 | CHARACTER(len=4),SAVE :: type_trac |
---|
54 | CHARACTER(len=8),DIMENSION(:),ALLOCATABLE, SAVE :: solsym |
---|
55 | |
---|
56 | ! CRisi: cas particulier des isotopes |
---|
57 | LOGICAL,SAVE :: ok_isotopes,ok_iso_verif,ok_isotrac,ok_init_iso |
---|
58 | INTEGER :: niso_possibles |
---|
59 | PARAMETER ( niso_possibles=5) |
---|
60 | REAL, DIMENSION (niso_possibles),SAVE :: tnat,alpha_ideal |
---|
61 | LOGICAL, DIMENSION(niso_possibles),SAVE :: use_iso |
---|
62 | INTEGER, ALLOCATABLE, DIMENSION(:,:), SAVE :: iqiso ! donne indice iq en fn de (ixt,phase) |
---|
63 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: iso_num ! donne numero iso entre 1 et niso_possibles en fn de nqtot |
---|
64 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: iso_indnum ! donne numero iso entre 1 et niso effectif en fn de nqtot |
---|
65 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: zone_num ! donne numero de la zone de tracage en fn de nqtot |
---|
66 | INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: phase_num ! donne numero de la zone de tracage en fn de nqtot |
---|
67 | INTEGER, DIMENSION(niso_possibles), SAVE :: indnum_fn_num ! donne indice entre entre 1 et niso en fonction du numero d isotope entre 1 et niso_possibles |
---|
68 | INTEGER, ALLOCATABLE, DIMENSION(:,:), SAVE :: index_trac ! numero ixt en fn izone, indnum entre 1 et niso |
---|
69 | INTEGER,SAVE :: niso,ntraceurs_zone,ntraciso |
---|
70 | |
---|
71 | #ifdef CPP_StratAer |
---|
72 | !--CK/OB for stratospheric aerosols |
---|
73 | INTEGER, SAVE :: nbtr_bin |
---|
74 | INTEGER, SAVE :: nbtr_sulgas |
---|
75 | INTEGER, SAVE :: id_OCS_strat |
---|
76 | INTEGER, SAVE :: id_SO2_strat |
---|
77 | INTEGER, SAVE :: id_H2SO4_strat |
---|
78 | INTEGER, SAVE :: id_BIN01_strat |
---|
79 | INTEGER, SAVE :: id_TEST_strat |
---|
80 | #endif |
---|
81 | |
---|
82 | CONTAINS |
---|
83 | |
---|
84 | SUBROUTINE infotrac_init |
---|
85 | USE control_mod, ONLY: planet_type, config_inca |
---|
86 | #ifdef REPROBUS |
---|
87 | USE CHEM_REP, ONLY : Init_chem_rep_trac |
---|
88 | #endif |
---|
89 | IMPLICIT NONE |
---|
90 | !======================================================================= |
---|
91 | ! |
---|
92 | ! Auteur: P. Le Van /L. Fairhead/F.Hourdin |
---|
93 | ! ------- |
---|
94 | ! Modif special traceur F.Forget 05/94 |
---|
95 | ! Modif M-A Filiberti 02/02 lecture de traceur.def |
---|
96 | ! |
---|
97 | ! Objet: |
---|
98 | ! ------ |
---|
99 | ! GCM LMD nouvelle grille |
---|
100 | ! |
---|
101 | !======================================================================= |
---|
102 | ! ... modification de l'integration de q ( 26/04/94 ) .... |
---|
103 | !----------------------------------------------------------------------- |
---|
104 | ! Declarations |
---|
105 | |
---|
106 | INCLUDE "dimensions.h" |
---|
107 | INCLUDE "iniprint.h" |
---|
108 | |
---|
109 | ! Local variables |
---|
110 | INTEGER, ALLOCATABLE, DIMENSION(:) :: hadv ! index of horizontal trasport schema |
---|
111 | INTEGER, ALLOCATABLE, DIMENSION(:) :: vadv ! index of vertical trasport schema |
---|
112 | |
---|
113 | INTEGER, ALLOCATABLE, DIMENSION(:) :: hadv_inca ! index of horizontal trasport schema |
---|
114 | INTEGER, ALLOCATABLE, DIMENSION(:) :: vadv_inca ! index of vertical trasport schema |
---|
115 | |
---|
116 | INTEGER, ALLOCATABLE, DIMENSION(:) :: conv_flg_inca |
---|
117 | INTEGER, ALLOCATABLE, DIMENSION(:) :: pbl_flg_inca |
---|
118 | CHARACTER(len=8), ALLOCATABLE, DIMENSION(:) :: solsym_inca |
---|
119 | |
---|
120 | CHARACTER(len=30), ALLOCATABLE, DIMENSION(:) :: tnom_0 ! tracer short name |
---|
121 | CHARACTER(len=30), ALLOCATABLE, DIMENSION(:) :: tnom_transp ! transporting fluid short name: CRisi |
---|
122 | CHARACTER(len=3), DIMENSION(30) :: descrq |
---|
123 | CHARACTER(len=1), DIMENSION(3) :: txts |
---|
124 | CHARACTER(len=2), DIMENSION(9) :: txtp |
---|
125 | CHARACTER(len=23) :: str1,str2 |
---|
126 | |
---|
127 | INTEGER :: nqtrue ! number of tracers read from tracer.def, without higer order of moment |
---|
128 | INTEGER :: iq, new_iq, iiq, jq, ierr,itr |
---|
129 | INTEGER :: ifils,ipere,generation ! CRisi |
---|
130 | LOGICAL :: continu,nouveau_traceurdef |
---|
131 | INTEGER :: IOstatus ! gestion de la retrocompatibilite de traceur.def |
---|
132 | CHARACTER(len=30) :: tchaine |
---|
133 | |
---|
134 | character(len=*),parameter :: modname="infotrac_init" |
---|
135 | |
---|
136 | !----------------------------------------------------------------------- |
---|
137 | ! Initialization : |
---|
138 | ! |
---|
139 | txts=(/'x','y','z'/) |
---|
140 | txtp=(/'x ','y ','z ','xx','xy','xz','yy','yz','zz'/) |
---|
141 | |
---|
142 | descrq(14)='VLH' |
---|
143 | descrq(10)='VL1' |
---|
144 | descrq(11)='VLP' |
---|
145 | descrq(12)='FH1' |
---|
146 | descrq(13)='FH2' |
---|
147 | descrq(16)='PPM' |
---|
148 | descrq(17)='PPS' |
---|
149 | descrq(18)='PPP' |
---|
150 | descrq(20)='SLP' |
---|
151 | descrq(30)='PRA' |
---|
152 | |
---|
153 | |
---|
154 | ! Coherence test between parameter type_trac, config_inca and preprocessing keys |
---|
155 | IF (type_trac=='inca') THEN |
---|
156 | WRITE(lunout,*) 'You have chosen to couple with INCA chemistry model : type_trac=', & |
---|
157 | type_trac,' config_inca=',config_inca |
---|
158 | IF (config_inca/='aero' .AND. config_inca/='aeNP' .AND. config_inca/='chem') THEN |
---|
159 | WRITE(lunout,*) 'Incoherence between type_trac and config_inca. Model stops. Modify run.def' |
---|
160 | CALL abort_gcm('infotrac_init','Incoherence between type_trac and config_inca',1) |
---|
161 | ENDIF |
---|
162 | #ifndef INCA |
---|
163 | WRITE(lunout,*) 'To run this option you must add cpp key INCA and compile with INCA code' |
---|
164 | CALL abort_gcm('infotrac_init','You must compile with cpp key INCA',1) |
---|
165 | #endif |
---|
166 | ELSE IF (type_trac=='repr') THEN |
---|
167 | WRITE(lunout,*) 'You have chosen to couple with REPROBUS chemestry model : type_trac=', type_trac |
---|
168 | #ifndef REPROBUS |
---|
169 | WRITE(lunout,*) 'To run this option you must add cpp key REPROBUS and compile with REPRPBUS code' |
---|
170 | CALL abort_gcm('infotrac_init','You must compile with cpp key REPROBUS',1) |
---|
171 | #endif |
---|
172 | ELSE IF (type_trac == 'co2i') THEN |
---|
173 | WRITE(lunout,*) 'You have chosen to run with CO2 cycle: type_trac=', type_trac |
---|
174 | ELSE IF (type_trac == 'coag') THEN |
---|
175 | WRITE(lunout,*) 'Tracers are treated for COAGULATION tests : type_trac=', type_trac |
---|
176 | #ifndef CPP_StratAer |
---|
177 | WRITE(lunout,*) 'To run this option you must add cpp key StratAer and compile with StratAer code' |
---|
178 | CALL abort_gcm('infotrac_init','You must compile with cpp key StratAer',1) |
---|
179 | #endif |
---|
180 | ELSE IF (type_trac == 'lmdz') THEN |
---|
181 | WRITE(lunout,*) 'Tracers are treated in LMDZ only : type_trac=', type_trac |
---|
182 | ELSE IF (type_trac == 'inco') THEN ! ThL |
---|
183 | WRITE(lunout,*) 'Using jointly INCA and CO2 cycle: type_trac =', type_trac |
---|
184 | IF (config_inca/='aero' .AND. config_inca/='aeNP' .AND. config_inca/='chem') THEN |
---|
185 | WRITE(lunout,*) 'Incoherence between type_trac and config_inca. Model stops. Modify run.def' |
---|
186 | CALL abort_gcm('infotrac_init','Incoherence between type_trac and config_inca',1) |
---|
187 | ENDIF |
---|
188 | #ifndef INCA |
---|
189 | WRITE(lunout,*) 'To run this option you must add cpp key INCA and compilewith INCA code' |
---|
190 | CALL abort_gcm('infotrac_init','You must compile with cpp key INCA',1) |
---|
191 | #endif |
---|
192 | ELSE |
---|
193 | WRITE(lunout,*) 'type_trac=',type_trac,' not possible. Model stops' |
---|
194 | CALL abort_gcm('infotrac_init','bad parameter',1) |
---|
195 | ENDIF |
---|
196 | |
---|
197 | ! Test if config_inca is other then none for run without INCA |
---|
198 | IF (type_trac/='inca' .AND. type_trac/='inco' .AND. config_inca/='none') THEN |
---|
199 | WRITE(lunout,*) 'config_inca will now be changed to none as you do not couple with INCA model' |
---|
200 | config_inca='none' |
---|
201 | ENDIF |
---|
202 | |
---|
203 | !----------------------------------------------------------------------- |
---|
204 | ! |
---|
205 | ! 1) Get the true number of tracers + water vapor/liquid |
---|
206 | ! Here true tracers (nqtrue) means declared tracers (only first order) |
---|
207 | ! |
---|
208 | !----------------------------------------------------------------------- |
---|
209 | IF (type_trac == 'lmdz' .OR. type_trac == 'repr' .OR. type_trac == 'coag' .OR. type_trac == 'co2i') THEN |
---|
210 | IF (type_trac=='co2i') THEN |
---|
211 | nqCO2 = 1 |
---|
212 | ELSE |
---|
213 | nqCO2 = 0 |
---|
214 | ENDIF |
---|
215 | OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr) |
---|
216 | IF(ierr.EQ.0) THEN |
---|
217 | WRITE(lunout,*) trim(modname),': Open traceur.def : ok' |
---|
218 | READ(90,*) nqtrue |
---|
219 | write(lunout,*) 'nqtrue=',nqtrue |
---|
220 | ELSE |
---|
221 | WRITE(lunout,*) trim(modname),': Problem in opening traceur.def' |
---|
222 | WRITE(lunout,*) trim(modname),': WARNING using defaut values' |
---|
223 | IF (planet_type=='earth') THEN |
---|
224 | nqtrue=4 ! Default value for Earth |
---|
225 | ELSE |
---|
226 | nqtrue=1 ! Default value for other planets |
---|
227 | ENDIF |
---|
228 | ENDIF |
---|
229 | !jyg< |
---|
230 | !! if ( planet_type=='earth') then |
---|
231 | !! ! For Earth, water vapour & liquid tracers are not in the physics |
---|
232 | !! nbtr=nqtrue-2 |
---|
233 | !! else |
---|
234 | !! ! Other planets (for now); we have the same number of tracers |
---|
235 | !! ! in the dynamics than in the physics |
---|
236 | !! nbtr=nqtrue |
---|
237 | !! endif |
---|
238 | !>jyg |
---|
239 | ELSE ! type_trac=inca or inco |
---|
240 | IF (type_trac=='inco') THEN |
---|
241 | nqCO2 = 1 |
---|
242 | ELSE |
---|
243 | nqCO2 = 0 |
---|
244 | ENDIF |
---|
245 | !jyg< |
---|
246 | ! The traceur.def file is used to define the number "nqo" of water phases |
---|
247 | ! present in the simulation. Default : nqo = 2. |
---|
248 | OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr) |
---|
249 | IF(ierr.EQ.0) THEN |
---|
250 | WRITE(lunout,*) trim(modname),': Open traceur.def : ok' |
---|
251 | READ(90,*) nqo |
---|
252 | ELSE |
---|
253 | WRITE(lunout,*) trim(modname),': Using default value for nqo' |
---|
254 | nqo=2 |
---|
255 | ENDIF |
---|
256 | IF (nqo /= 2 .AND. nqo /= 3 ) THEN |
---|
257 | IF (nqo == 4 .AND. type_trac=='inco') THEN ! ThL |
---|
258 | WRITE(lunout,*) trim(modname),': you are coupling with INCA, and also using CO2i.' |
---|
259 | nqo = 3 ! A ameliorier... je force 3 traceurs eau... ThL |
---|
260 | WRITE(lunout,*) trim(modname),': nqo = ',nqo |
---|
261 | ELSE |
---|
262 | WRITE(lunout,*) trim(modname),': nqo=',nqo, ' is not allowed. Only 2 or 3 water phases allowed' |
---|
263 | CALL abort_gcm('infotrac_init','Bad number of water phases',1) |
---|
264 | ENDIF |
---|
265 | ENDIF |
---|
266 | ! nbtr has been read from INCA by init_const_lmdz() in gcm.F |
---|
267 | #ifdef INCA |
---|
268 | CALL Init_chem_inca_trac(nqINCA) |
---|
269 | #else |
---|
270 | nqINCA=0 |
---|
271 | #endif |
---|
272 | nbtr=nqINCA+nqCO2 |
---|
273 | nqtrue=nbtr+nqo |
---|
274 | WRITE(lunout,*) trim(modname),': nqo = ',nqo |
---|
275 | WRITE(lunout,*) trim(modname),': nbtr = ',nbtr |
---|
276 | WRITE(lunout,*) trim(modname),': nqtrue = ',nqtrue |
---|
277 | WRITE(lunout,*) trim(modname),': nqCO2 = ',nqCO2 |
---|
278 | WRITE(lunout,*) trim(modname),': nqINCA = ',nqINCA |
---|
279 | ALLOCATE(hadv_inca(nqINCA), vadv_inca(nqINCA), conv_flg_inca(nqINCA), pbl_flg_inca(nqINCA), solsym_inca(nqINCA)) |
---|
280 | ENDIF ! type_trac 'inca' ou 'inco' |
---|
281 | !>jyg |
---|
282 | |
---|
283 | IF ((planet_type=="earth").and.(nqtrue < 2)) THEN |
---|
284 | WRITE(lunout,*) trim(modname),': nqtrue=',nqtrue, ' is not allowed. 2 tracers is the minimum' |
---|
285 | CALL abort_gcm('infotrac_init','Not enough tracers',1) |
---|
286 | ENDIF |
---|
287 | |
---|
288 | !jyg< |
---|
289 | |
---|
290 | ! |
---|
291 | ! Allocate variables depending on nqtrue |
---|
292 | ! |
---|
293 | ALLOCATE(tnom_0(nqtrue), hadv(nqtrue), vadv(nqtrue),tnom_transp(nqtrue)) |
---|
294 | |
---|
295 | |
---|
296 | !----------------------------------------------------------------------- |
---|
297 | ! 2) Choix des schemas d'advection pour l'eau et les traceurs |
---|
298 | ! |
---|
299 | ! iadv = 1 schema transport type "humidite specifique LMD" |
---|
300 | ! iadv = 2 schema amont |
---|
301 | ! iadv = 14 schema Van-leer + humidite specifique |
---|
302 | ! Modif F.Codron |
---|
303 | ! iadv = 10 schema Van-leer (retenu pour l'eau vapeur et liquide) |
---|
304 | ! iadv = 11 schema Van-Leer pour hadv et version PPM (Monotone) pour vadv |
---|
305 | ! iadv = 12 schema Frederic Hourdin I |
---|
306 | ! iadv = 13 schema Frederic Hourdin II |
---|
307 | ! iadv = 16 schema PPM Monotone(Collela & Woodward 1984) |
---|
308 | ! iadv = 17 schema PPM Semi Monotone (overshoots autorises) |
---|
309 | ! iadv = 18 schema PPM Positif Defini (overshoots undershoots autorises) |
---|
310 | ! iadv = 20 schema Slopes |
---|
311 | ! iadv = 30 schema Prather |
---|
312 | ! |
---|
313 | ! Dans le tableau q(ij,l,iq) : iq = 1 pour l'eau vapeur |
---|
314 | ! iq = 2 pour l'eau liquide |
---|
315 | ! Et eventuellement iq = 3,nqtot pour les autres traceurs |
---|
316 | ! |
---|
317 | ! iadv(1): choix pour l'eau vap. et iadv(2) : choix pour l'eau liq. |
---|
318 | !------------------------------------------------------------------------ |
---|
319 | ! |
---|
320 | ! Get choice of advection schema from file tracer.def or from INCA |
---|
321 | !--------------------------------------------------------------------- |
---|
322 | IF (type_trac == 'lmdz' .OR. type_trac == 'repr' .OR. type_trac == 'coag' .OR. type_trac == 'co2i') THEN |
---|
323 | |
---|
324 | IF(ierr.EQ.0) THEN |
---|
325 | |
---|
326 | ! Continue to read tracer.def |
---|
327 | DO iq=1,nqtrue |
---|
328 | |
---|
329 | write(*,*) 'infotrac 237: iq=',iq |
---|
330 | ! CRisi: ajout du nom du fluide transporteur |
---|
331 | ! mais rester retro compatible |
---|
332 | READ(90,'(I2,X,I2,X,A)',IOSTAT=IOstatus) hadv(iq),vadv(iq),tchaine |
---|
333 | write(lunout,*) 'iq,hadv(iq),vadv(iq)=',iq,hadv(iq),vadv(iq) |
---|
334 | write(lunout,*) 'tchaine=',trim(tchaine) |
---|
335 | write(*,*) 'infotrac 238: IOstatus=',IOstatus |
---|
336 | if (IOstatus.ne.0) then |
---|
337 | CALL abort_gcm('infotrac_init','Pb dans la lecture de traceur.def',1) |
---|
338 | endif |
---|
339 | ! Y-a-t-il 1 ou 2 noms de traceurs? -> On regarde s'il y a un |
---|
340 | ! espace ou pas au milieu de la chaine. |
---|
341 | continu=.true. |
---|
342 | nouveau_traceurdef=.false. |
---|
343 | iiq=1 |
---|
344 | do while (continu) |
---|
345 | if (tchaine(iiq:iiq).eq.' ') then |
---|
346 | nouveau_traceurdef=.true. |
---|
347 | continu=.false. |
---|
348 | else if (iiq.lt.LEN_TRIM(tchaine)) then |
---|
349 | iiq=iiq+1 |
---|
350 | else |
---|
351 | continu=.false. |
---|
352 | endif |
---|
353 | enddo |
---|
354 | write(*,*) 'iiq,nouveau_traceurdef=',iiq,nouveau_traceurdef |
---|
355 | if (nouveau_traceurdef) then |
---|
356 | write(lunout,*) 'C''est la nouvelle version de traceur.def' |
---|
357 | tnom_0(iq)=tchaine(1:iiq-1) |
---|
358 | tnom_transp(iq)=tchaine(iiq+1:30) |
---|
359 | else |
---|
360 | write(lunout,*) 'C''est l''ancienne version de traceur.def' |
---|
361 | write(lunout,*) 'On suppose que les traceurs sont tous d''air' |
---|
362 | tnom_0(iq)=tchaine |
---|
363 | tnom_transp(iq) = 'air' |
---|
364 | endif |
---|
365 | write(lunout,*) 'tnom_0(iq)=<',trim(tnom_0(iq)),'>' |
---|
366 | write(lunout,*) 'tnom_transp(iq)=<',trim(tnom_transp(iq)),'>' |
---|
367 | |
---|
368 | ENDDO!DO iq=1,nqtrue |
---|
369 | |
---|
370 | CLOSE(90) |
---|
371 | |
---|
372 | ELSE ! Without tracer.def, set default values |
---|
373 | |
---|
374 | if (planet_type=="earth") then |
---|
375 | ! for Earth, default is to have 4 tracers |
---|
376 | hadv(1) = 14 |
---|
377 | vadv(1) = 14 |
---|
378 | tnom_0(1) = 'H2Ov' |
---|
379 | tnom_transp(1) = 'air' |
---|
380 | hadv(2) = 10 |
---|
381 | vadv(2) = 10 |
---|
382 | tnom_0(2) = 'H2Ol' |
---|
383 | tnom_transp(2) = 'air' |
---|
384 | hadv(3) = 10 |
---|
385 | vadv(3) = 10 |
---|
386 | tnom_0(3) = 'RN' |
---|
387 | tnom_transp(3) = 'air' |
---|
388 | hadv(4) = 10 |
---|
389 | vadv(4) = 10 |
---|
390 | tnom_0(4) = 'PB' |
---|
391 | tnom_transp(4) = 'air' |
---|
392 | else ! default for other planets |
---|
393 | hadv(1) = 10 |
---|
394 | vadv(1) = 10 |
---|
395 | tnom_0(1) = 'dummy' |
---|
396 | tnom_transp(1) = 'dummy' |
---|
397 | endif ! of if (planet_type=="earth") |
---|
398 | |
---|
399 | ENDIF |
---|
400 | |
---|
401 | WRITE(lunout,*) trim(modname),': Valeur de traceur.def :' |
---|
402 | WRITE(lunout,*) trim(modname),': nombre total de traceurs ',nqtrue |
---|
403 | DO iq=1,nqtrue |
---|
404 | WRITE(lunout,*) hadv(iq),vadv(iq),tnom_0(iq),tnom_transp(iq) |
---|
405 | END DO |
---|
406 | |
---|
407 | IF ( planet_type=='earth') THEN |
---|
408 | !CR: nombre de traceurs de l eau |
---|
409 | IF (tnom_0(3) == 'H2Oi') THEN |
---|
410 | nqo=3 |
---|
411 | ELSE |
---|
412 | nqo=2 |
---|
413 | ENDIF |
---|
414 | ! For Earth, water vapour & liquid tracers are not in the physics |
---|
415 | nbtr=nqtrue-nqo |
---|
416 | ELSE |
---|
417 | ! Other planets (for now); we have the same number of tracers |
---|
418 | ! in the dynamics than in the physics |
---|
419 | nbtr=nqtrue |
---|
420 | ENDIF |
---|
421 | |
---|
422 | #ifdef CPP_StratAer |
---|
423 | IF (type_trac == 'coag') THEN |
---|
424 | nbtr_bin=0 |
---|
425 | nbtr_sulgas=0 |
---|
426 | DO iq=1,nqtrue |
---|
427 | IF (tnom_0(iq)(1:3)=='BIN') THEN !check if tracer name contains 'BIN' |
---|
428 | nbtr_bin=nbtr_bin+1 |
---|
429 | ENDIF |
---|
430 | IF (tnom_0(iq)(1:3)=='GAS') THEN !check if tracer name contains 'GAS' |
---|
431 | nbtr_sulgas=nbtr_sulgas+1 |
---|
432 | ENDIF |
---|
433 | ENDDO |
---|
434 | print*,'nbtr_bin=',nbtr_bin |
---|
435 | print*,'nbtr_sulgas=',nbtr_sulgas |
---|
436 | DO iq=1,nqtrue |
---|
437 | IF (tnom_0(iq)=='GASOCS') THEN |
---|
438 | id_OCS_strat=iq-nqo |
---|
439 | ENDIF |
---|
440 | IF (tnom_0(iq)=='GASSO2') THEN |
---|
441 | id_SO2_strat=iq-nqo |
---|
442 | ENDIF |
---|
443 | IF (tnom_0(iq)=='GASH2SO4') THEN |
---|
444 | id_H2SO4_strat=iq-nqo |
---|
445 | ENDIF |
---|
446 | IF (tnom_0(iq)=='BIN01') THEN |
---|
447 | id_BIN01_strat=iq-nqo |
---|
448 | ENDIF |
---|
449 | IF (tnom_0(iq)=='GASTEST') THEN |
---|
450 | id_TEST_strat=iq-nqo |
---|
451 | ENDIF |
---|
452 | ENDDO |
---|
453 | print*,'id_OCS_strat =',id_OCS_strat |
---|
454 | print*,'id_SO2_strat =',id_SO2_strat |
---|
455 | print*,'id_H2SO4_strat=',id_H2SO4_strat |
---|
456 | print*,'id_BIN01_strat=',id_BIN01_strat |
---|
457 | ENDIF |
---|
458 | #endif |
---|
459 | |
---|
460 | ENDIF ! (type_trac == 'lmdz' .OR. type_trac == 'repr' .OR. type_trac = 'coag' .OR. type_trac = 'co2i') |
---|
461 | !jyg< |
---|
462 | ! |
---|
463 | |
---|
464 | ! Transfert number of tracers to Reprobus |
---|
465 | IF (type_trac == 'repr') THEN |
---|
466 | #ifdef REPROBUS |
---|
467 | CALL Init_chem_rep_trac(nbtr,nqo,tnom_0) |
---|
468 | #endif |
---|
469 | ENDIF |
---|
470 | ! |
---|
471 | ! Allocate variables depending on nbtr |
---|
472 | ! |
---|
473 | ALLOCATE(conv_flg(nbtr), pbl_flg(nbtr), solsym(nbtr)) |
---|
474 | conv_flg(:) = 1 ! convection activated for all tracers |
---|
475 | pbl_flg(:) = 1 ! boundary layer activated for all tracers |
---|
476 | |
---|
477 | IF (type_trac == 'inca' .OR. type_trac == 'inco') THEN ! config_inca='aero' ou 'chem' |
---|
478 | !>jyg |
---|
479 | ! le module de chimie fournit les noms des traceurs |
---|
480 | ! et les schemas d'advection associes. excepte pour ceux lus |
---|
481 | ! dans traceur.def |
---|
482 | IF (ierr .eq. 0) then |
---|
483 | |
---|
484 | DO iq=1,nqo+nqCO2 |
---|
485 | |
---|
486 | write(*,*) 'infotrac 237: iq=',iq |
---|
487 | ! CRisi: ajout du nom du fluide transporteur |
---|
488 | ! mais rester retro compatible |
---|
489 | READ(90,'(I2,X,I2,X,A)',IOSTAT=IOstatus) hadv(iq),vadv(iq),tchaine |
---|
490 | write(lunout,*) 'iq,hadv(iq),vadv(iq)=',iq,hadv(iq),vadv(iq) |
---|
491 | write(lunout,*) 'tchaine=',trim(tchaine) |
---|
492 | write(*,*) 'infotrac 238: IOstatus=',IOstatus |
---|
493 | if (IOstatus.ne.0) then |
---|
494 | CALL abort_gcm('infotrac_init','Pb dans la lecture de traceur.def',1) |
---|
495 | endif |
---|
496 | ! Y-a-t-il 1 ou 2 noms de traceurs? -> On regarde s'il y a un |
---|
497 | ! espace ou pas au milieu de la chaine. |
---|
498 | continu=.true. |
---|
499 | nouveau_traceurdef=.false. |
---|
500 | iiq=1 |
---|
501 | |
---|
502 | do while (continu) |
---|
503 | if (tchaine(iiq:iiq).eq.' ') then |
---|
504 | nouveau_traceurdef=.true. |
---|
505 | continu=.false. |
---|
506 | else if (iiq.lt.LEN_TRIM(tchaine)) then |
---|
507 | iiq=iiq+1 |
---|
508 | else |
---|
509 | continu=.false. |
---|
510 | endif |
---|
511 | enddo |
---|
512 | |
---|
513 | write(*,*) 'iiq,nouveau_traceurdef=',iiq,nouveau_traceurdef |
---|
514 | |
---|
515 | if (nouveau_traceurdef) then |
---|
516 | write(lunout,*) 'C''est la nouvelle version de traceur.def' |
---|
517 | tnom_0(iq)=tchaine(1:iiq-1) |
---|
518 | tnom_transp(iq)=tchaine(iiq+1:30) |
---|
519 | else |
---|
520 | write(lunout,*) 'C''est l''ancienne version de traceur.def' |
---|
521 | write(lunout,*) 'On suppose que les traceurs sont tous d''air' |
---|
522 | tnom_0(iq)=tchaine |
---|
523 | tnom_transp(iq) = 'air' |
---|
524 | endif |
---|
525 | |
---|
526 | write(lunout,*) 'tnom_0(iq)=<',trim(tnom_0(iq)),'>' |
---|
527 | write(lunout,*) 'tnom_transp(iq)=<',trim(tnom_transp(iq)),'>' |
---|
528 | |
---|
529 | ENDDO !DO iq=1,nqo |
---|
530 | CLOSE(90) |
---|
531 | |
---|
532 | ELSE !! if traceur.def doesn't exist |
---|
533 | |
---|
534 | tnom_0(1)='H2Ov' |
---|
535 | tnom_transp(1) = 'air' |
---|
536 | tnom_0(2)='H2Ol' |
---|
537 | tnom_transp(2) = 'air' |
---|
538 | hadv(1) = 10 |
---|
539 | hadv(2) = 10 |
---|
540 | vadv(1) = 10 |
---|
541 | vadv(2) = 10 |
---|
542 | |
---|
543 | ENDIF |
---|
544 | |
---|
545 | #ifdef INCA |
---|
546 | CALL init_transport( & |
---|
547 | hadv_inca, & |
---|
548 | vadv_inca, & |
---|
549 | conv_flg_inca, & |
---|
550 | pbl_flg_inca, & |
---|
551 | solsym_inca) |
---|
552 | |
---|
553 | conv_flg(1+nqCO2:nbtr) = conv_flg_inca |
---|
554 | pbl_flg(1+nqCO2:nbtr) = pbl_flg_inca |
---|
555 | solsym(1+nqCO2:nbtr) = solsym_inca |
---|
556 | |
---|
557 | IF (type_trac == 'inco') THEN |
---|
558 | conv_flg(1:nqCO2) = 1 |
---|
559 | pbl_flg(1:nqCO2) = 1 |
---|
560 | solsym(1:nqCO2) = 'CO2' |
---|
561 | ENDIF |
---|
562 | #endif |
---|
563 | |
---|
564 | !jyg< |
---|
565 | DO iq = nqo+nqCO2+1, nqtrue |
---|
566 | hadv(iq) = hadv_inca(iq-nqo-nqCO2) |
---|
567 | vadv(iq) = vadv_inca(iq-nqo-nqCO2) |
---|
568 | tnom_0(iq)=solsym_inca(iq-nqo-nqCO2) |
---|
569 | tnom_transp(iq) = 'air' |
---|
570 | END DO |
---|
571 | |
---|
572 | ENDIF ! (type_trac == 'inca' or 'inco') |
---|
573 | |
---|
574 | !----------------------------------------------------------------------- |
---|
575 | ! |
---|
576 | ! 3) Verify if advection schema 20 or 30 choosen |
---|
577 | ! Calculate total number of tracers needed: nqtot |
---|
578 | ! Allocate variables depending on total number of tracers |
---|
579 | !----------------------------------------------------------------------- |
---|
580 | new_iq=0 |
---|
581 | DO iq=1,nqtrue |
---|
582 | ! Add tracers for certain advection schema |
---|
583 | IF (hadv(iq)<20 .AND. vadv(iq)<20 ) THEN |
---|
584 | new_iq=new_iq+1 ! no tracers added |
---|
585 | ELSE IF (hadv(iq)==20 .AND. vadv(iq)==20 ) THEN |
---|
586 | new_iq=new_iq+4 ! 3 tracers added |
---|
587 | ELSE IF (hadv(iq)==30 .AND. vadv(iq)==30 ) THEN |
---|
588 | new_iq=new_iq+10 ! 9 tracers added |
---|
589 | ELSE |
---|
590 | WRITE(lunout,*) trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq) |
---|
591 | CALL abort_gcm('infotrac_init','Bad choice of advection schema - 1',1) |
---|
592 | ENDIF |
---|
593 | END DO |
---|
594 | |
---|
595 | IF (new_iq /= nqtrue) THEN |
---|
596 | ! The choice of advection schema imposes more tracers |
---|
597 | ! Assigne total number of tracers |
---|
598 | nqtot = new_iq |
---|
599 | |
---|
600 | WRITE(lunout,*) trim(modname),': The choice of advection schema for one or more tracers' |
---|
601 | WRITE(lunout,*) 'makes it necessary to add tracers' |
---|
602 | WRITE(lunout,*) trim(modname)//': ',nqtrue,' is the number of true tracers' |
---|
603 | WRITE(lunout,*) trim(modname)//': ',nqtot, ' is the total number of tracers needed' |
---|
604 | |
---|
605 | ELSE |
---|
606 | ! The true number of tracers is also the total number |
---|
607 | nqtot = nqtrue |
---|
608 | ENDIF |
---|
609 | |
---|
610 | ! |
---|
611 | ! Allocate variables with total number of tracers, nqtot |
---|
612 | ! |
---|
613 | ALLOCATE(tname(nqtot), ttext(nqtot)) |
---|
614 | ALLOCATE(iadv(nqtot), niadv(nqtot)) |
---|
615 | |
---|
616 | !----------------------------------------------------------------------- |
---|
617 | ! |
---|
618 | ! 4) Determine iadv, long and short name |
---|
619 | ! |
---|
620 | !----------------------------------------------------------------------- |
---|
621 | new_iq=0 |
---|
622 | DO iq=1,nqtrue |
---|
623 | new_iq=new_iq+1 |
---|
624 | |
---|
625 | ! Verify choice of advection schema |
---|
626 | IF (hadv(iq)==vadv(iq)) THEN |
---|
627 | iadv(new_iq)=hadv(iq) |
---|
628 | ELSE IF (hadv(iq)==10 .AND. vadv(iq)==16) THEN |
---|
629 | iadv(new_iq)=11 |
---|
630 | ELSE |
---|
631 | WRITE(lunout,*)trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq) |
---|
632 | |
---|
633 | CALL abort_gcm('infotrac_init','Bad choice of advection schema - 2',1) |
---|
634 | ENDIF |
---|
635 | |
---|
636 | str1=tnom_0(iq) |
---|
637 | tname(new_iq)= tnom_0(iq) |
---|
638 | IF (iadv(new_iq)==0) THEN |
---|
639 | ttext(new_iq)=trim(str1) |
---|
640 | ELSE |
---|
641 | ttext(new_iq)=trim(tnom_0(iq))//descrq(iadv(new_iq)) |
---|
642 | ENDIF |
---|
643 | |
---|
644 | ! schemas tenant compte des moments d'ordre superieur |
---|
645 | str2=ttext(new_iq) |
---|
646 | IF (iadv(new_iq)==20) THEN |
---|
647 | DO jq=1,3 |
---|
648 | new_iq=new_iq+1 |
---|
649 | iadv(new_iq)=-20 |
---|
650 | ttext(new_iq)=trim(str2)//txts(jq) |
---|
651 | tname(new_iq)=trim(str1)//txts(jq) |
---|
652 | END DO |
---|
653 | ELSE IF (iadv(new_iq)==30) THEN |
---|
654 | DO jq=1,9 |
---|
655 | new_iq=new_iq+1 |
---|
656 | iadv(new_iq)=-30 |
---|
657 | ttext(new_iq)=trim(str2)//txtp(jq) |
---|
658 | tname(new_iq)=trim(str1)//txtp(jq) |
---|
659 | END DO |
---|
660 | ENDIF |
---|
661 | END DO |
---|
662 | |
---|
663 | ! |
---|
664 | ! Find vector keeping the correspodence between true and total tracers |
---|
665 | ! |
---|
666 | niadv(:)=0 |
---|
667 | iiq=0 |
---|
668 | DO iq=1,nqtot |
---|
669 | IF(iadv(iq).GE.0) THEN |
---|
670 | ! True tracer |
---|
671 | iiq=iiq+1 |
---|
672 | niadv(iiq)=iq |
---|
673 | ENDIF |
---|
674 | END DO |
---|
675 | |
---|
676 | |
---|
677 | WRITE(lunout,*) trim(modname),': Information stored in infotrac :' |
---|
678 | WRITE(lunout,*) trim(modname),': iadv niadv tname ttext :' |
---|
679 | |
---|
680 | DO iq=1,nqtot |
---|
681 | WRITE(lunout,*) iadv(iq),niadv(iq), ' ',trim(tname(iq)),' ',trim(ttext(iq)) |
---|
682 | END DO |
---|
683 | |
---|
684 | ! |
---|
685 | ! Test for advection schema. |
---|
686 | ! This version of LMDZ only garantees iadv=10 and iadv=14 (14 only for water vapour) . |
---|
687 | ! |
---|
688 | DO iq=1,nqtot |
---|
689 | IF (iadv(iq)/=10 .AND. iadv(iq)/=14 .AND. iadv(iq)/=0) THEN |
---|
690 | WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv(iq),' is not tested in this version of LMDZ' |
---|
691 | CALL abort_gcm('infotrac_init','In this version only iadv=10 and iadv=14 is tested!',1) |
---|
692 | ELSE IF (iadv(iq)==14 .AND. iq/=1) THEN |
---|
693 | WRITE(lunout,*)trim(modname),'STOP : The option iadv=',iadv(iq),' is not tested in this version of LMDZ' |
---|
694 | CALL abort_gcm('infotrac_init','In this version iadv=14 is only permitted for water vapour!',1) |
---|
695 | ENDIF |
---|
696 | END DO |
---|
697 | |
---|
698 | |
---|
699 | ! CRisi: quels sont les traceurs fils et les traceurs peres. |
---|
700 | ! initialiser tous les tableaux d'indices lies aux traceurs familiaux |
---|
701 | ! + verifier que tous les peres sont ecrits en premieres positions |
---|
702 | ALLOCATE(nqfils(nqtot),nqdesc(nqtot)) |
---|
703 | ALLOCATE(iqfils(nqtot,nqtot)) |
---|
704 | ALLOCATE(iqpere(nqtot)) |
---|
705 | nqperes=0 |
---|
706 | nqfils(:)=0 |
---|
707 | nqdesc(:)=0 |
---|
708 | iqfils(:,:)=0 |
---|
709 | iqpere(:)=0 |
---|
710 | nqdesc_tot=0 |
---|
711 | DO iq=1,nqtot |
---|
712 | if (tnom_transp(iq) == 'air') then |
---|
713 | ! ceci est un traceur pere |
---|
714 | WRITE(lunout,*) 'Le traceur',iq,', appele ',trim(tnom_0(iq)),', est un pere' |
---|
715 | nqperes=nqperes+1 |
---|
716 | iqpere(iq)=0 |
---|
717 | else !if (tnom_transp(iq) == 'air') then |
---|
718 | ! ceci est un fils. Qui est son pere? |
---|
719 | WRITE(lunout,*) 'Le traceur',iq,', appele ',trim(tnom_0(iq)),', est un fils' |
---|
720 | continu=.true. |
---|
721 | ipere=1 |
---|
722 | do while (continu) |
---|
723 | if (tnom_transp(iq) == tnom_0(ipere)) then |
---|
724 | ! Son pere est ipere |
---|
725 | WRITE(lunout,*) 'Le traceur',iq,'appele ', & |
---|
726 | & trim(tnom_0(iq)),' est le fils de ',ipere,'appele ',trim(tnom_0(ipere)) |
---|
727 | if (iq.eq.ipere) then |
---|
728 | CALL abort_gcm('infotrac_init','Un fils est son propre pere',1) |
---|
729 | endif |
---|
730 | nqfils(ipere)=nqfils(ipere)+1 |
---|
731 | iqfils(nqfils(ipere),ipere)=iq |
---|
732 | iqpere(iq)=ipere |
---|
733 | continu=.false. |
---|
734 | else !if (tnom_transp(iq) == tnom_0(ipere)) then |
---|
735 | ipere=ipere+1 |
---|
736 | if (ipere.gt.nqtot) then |
---|
737 | WRITE(lunout,*) 'Le traceur',iq,'appele ', & |
---|
738 | & trim(tnom_0(iq)),', est orphelin.' |
---|
739 | CALL abort_gcm('infotrac_init','Un traceur est orphelin',1) |
---|
740 | endif !if (ipere.gt.nqtot) then |
---|
741 | endif !if (tnom_transp(iq) == tnom_0(ipere)) then |
---|
742 | enddo !do while (continu) |
---|
743 | endif !if (tnom_transp(iq) == 'air') then |
---|
744 | enddo !DO iq=1,nqtot |
---|
745 | WRITE(lunout,*) 'infotrac: nqperes=',nqperes |
---|
746 | WRITE(lunout,*) 'nqfils=',nqfils |
---|
747 | WRITE(lunout,*) 'iqpere=',iqpere |
---|
748 | WRITE(lunout,*) 'iqfils=',iqfils |
---|
749 | |
---|
750 | ! Calculer le nombre de descendants a partir de iqfils et de nbfils |
---|
751 | DO iq=1,nqtot |
---|
752 | generation=0 |
---|
753 | continu=.true. |
---|
754 | ifils=iq |
---|
755 | do while (continu) |
---|
756 | ipere=iqpere(ifils) |
---|
757 | if (ipere.gt.0) then |
---|
758 | nqdesc(ipere)=nqdesc(ipere)+1 |
---|
759 | nqdesc_tot=nqdesc_tot+1 |
---|
760 | iqfils(nqdesc(ipere),ipere)=iq |
---|
761 | ifils=ipere |
---|
762 | generation=generation+1 |
---|
763 | else !if (ipere.gt.0) then |
---|
764 | continu=.false. |
---|
765 | endif !if (ipere.gt.0) then |
---|
766 | enddo !do while (continu) |
---|
767 | WRITE(lunout,*) 'Le traceur ',iq,', appele ',trim(tnom_0(iq)),' est un traceur de generation: ',generation |
---|
768 | enddo !DO iq=1,nqtot |
---|
769 | WRITE(lunout,*) 'infotrac: nqdesc=',nqdesc |
---|
770 | WRITE(lunout,*) 'iqfils=',iqfils |
---|
771 | WRITE(lunout,*) 'nqdesc_tot=',nqdesc_tot |
---|
772 | |
---|
773 | ! Interdire autres schemas que 10 pour les traceurs fils, et autres schemas |
---|
774 | ! que 10 et 14 si des peres ont des fils |
---|
775 | do iq=1,nqtot |
---|
776 | if (iqpere(iq).gt.0) then |
---|
777 | ! ce traceur a un pere qui n'est pas l'air |
---|
778 | ! Seul le schema 10 est autorise |
---|
779 | if (iadv(iq)/=10) then |
---|
780 | WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv(iq),' is not implemented for sons' |
---|
781 | CALL abort_gcm('infotrac_init','Sons should be advected by scheme 10',1) |
---|
782 | endif |
---|
783 | ! Le traceur pere ne peut etre advecte que par schema 10 ou 14: |
---|
784 | IF (iadv(iqpere(iq))/=10 .AND. iadv(iqpere(iq))/=14) THEN |
---|
785 | WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv(iq),' is not implemented for fathers' |
---|
786 | CALL abort_gcm('infotrac_init','Fathers should be advected by scheme 10 ou 14',1) |
---|
787 | endif !IF (iadv(iqpere(iq))/=10 .AND. iadv(iqpere(iq))/=14) THEN |
---|
788 | endif !if (iqpere(iq).gt.0) the |
---|
789 | enddo !do iq=1,nqtot |
---|
790 | |
---|
791 | |
---|
792 | |
---|
793 | ! detecter quels sont les traceurs isotopiques parmi des traceurs |
---|
794 | call infotrac_isoinit(tnom_0,nqtrue) |
---|
795 | |
---|
796 | ! if (ntraciso.gt.0) then |
---|
797 | ! le 18 sep 2020: on enleve la condition ntraciso.gt.0 car nqtottr doit etre |
---|
798 | ! connu meme si il n'y a pas d'isotopes! |
---|
799 | write(lunout,*) 'infotrac 702: nbtr,ntraciso=',nbtr,ntraciso |
---|
800 | ! retrancher les traceurs isotopiques de la liste des traceurs qui passent dans |
---|
801 | ! phytrac |
---|
802 | nbtr=nbtr-nqo*ntraciso |
---|
803 | |
---|
804 | ! faire un tableau d'indice des traceurs qui passeront dans phytrac |
---|
805 | nqtottr=nqtot-nqo*(1+ntraciso) |
---|
806 | write(lunout,*) 'infotrac 704: nqtottr,nqtot,nqo=',nqtottr,nqtot,nqo |
---|
807 | ! Rq: nqtottr n'est pas forcement egal a nbtr dans le cas ou new_iq /= nqtrue |
---|
808 | ALLOCATE (itr_indice(nqtottr)) |
---|
809 | itr_indice(:)=0 |
---|
810 | itr=0 |
---|
811 | do iq=nqo+1, nqtot |
---|
812 | if (iso_num(iq).eq.0) then |
---|
813 | itr=itr+1 |
---|
814 | write(*,*) 'itr=',itr |
---|
815 | itr_indice(itr)=iq |
---|
816 | endif !if (iso_num(iq).eq.0) then |
---|
817 | enddo |
---|
818 | if (itr.ne.nqtottr) then |
---|
819 | CALL abort_gcm('infotrac_init','pb dans le calcul de nqtottr',1) |
---|
820 | endif |
---|
821 | write(lunout,*) 'itr_indice=',itr_indice |
---|
822 | ! endif !if (ntraciso.gt.0) then |
---|
823 | |
---|
824 | !----------------------------------------------------------------------- |
---|
825 | ! Finalize : |
---|
826 | ! |
---|
827 | DEALLOCATE(tnom_0, hadv, vadv,tnom_transp) |
---|
828 | |
---|
829 | WRITE(lunout,*) 'infotrac init fin' |
---|
830 | |
---|
831 | END SUBROUTINE infotrac_init |
---|
832 | |
---|
833 | SUBROUTINE infotrac_isoinit(tnom_0,nqtrue) |
---|
834 | |
---|
835 | #ifdef CPP_IOIPSL |
---|
836 | use IOIPSL |
---|
837 | #else |
---|
838 | ! if not using IOIPSL, we still need to use (a local version of) getin |
---|
839 | use ioipsl_getincom |
---|
840 | #endif |
---|
841 | implicit none |
---|
842 | |
---|
843 | ! inputs |
---|
844 | INTEGER nqtrue |
---|
845 | CHARACTER(len=30) tnom_0(nqtrue) |
---|
846 | |
---|
847 | ! locals |
---|
848 | CHARACTER(len=3), DIMENSION(niso_possibles) :: tnom_iso |
---|
849 | INTEGER, ALLOCATABLE,DIMENSION(:,:) :: nb_iso,nb_traciso |
---|
850 | INTEGER, ALLOCATABLE,DIMENSION(:) :: nb_isoind |
---|
851 | INTEGER :: ntraceurs_zone_prec,iq,phase,ixt,iiso,izone |
---|
852 | CHARACTER(len=19) :: tnom_trac |
---|
853 | INCLUDE "iniprint.h" |
---|
854 | |
---|
855 | tnom_iso=(/'eau','HDO','O18','O17','HTO'/) |
---|
856 | |
---|
857 | ALLOCATE(nb_iso(niso_possibles,nqo)) |
---|
858 | ALLOCATE(nb_isoind(nqo)) |
---|
859 | ALLOCATE(nb_traciso(niso_possibles,nqo)) |
---|
860 | ALLOCATE(iso_num(nqtot)) |
---|
861 | ALLOCATE(iso_indnum(nqtot)) |
---|
862 | ALLOCATE(zone_num(nqtot)) |
---|
863 | ALLOCATE(phase_num(nqtot)) |
---|
864 | |
---|
865 | iso_num(:)=0 |
---|
866 | iso_indnum(:)=0 |
---|
867 | zone_num(:)=0 |
---|
868 | phase_num(:)=0 |
---|
869 | indnum_fn_num(:)=0 |
---|
870 | use_iso(:)=.false. |
---|
871 | nb_iso(:,:)=0 |
---|
872 | nb_isoind(:)=0 |
---|
873 | nb_traciso(:,:)=0 |
---|
874 | niso=0 |
---|
875 | ntraceurs_zone=0 |
---|
876 | ntraceurs_zone_prec=0 |
---|
877 | ntraciso=0 |
---|
878 | |
---|
879 | do iq=nqo+1,nqtot |
---|
880 | ! write(lunout,*) 'infotrac 569: iq,tnom_0(iq)=',iq,tnom_0(iq) |
---|
881 | do phase=1,nqo |
---|
882 | do ixt= 1,niso_possibles |
---|
883 | tnom_trac=trim(tnom_0(phase))//'_' |
---|
884 | tnom_trac=trim(tnom_trac)//trim(tnom_iso(ixt)) |
---|
885 | ! write(*,*) 'phase,ixt,tnom_trac=',phase,ixt,tnom_trac |
---|
886 | IF (tnom_0(iq) == tnom_trac) then |
---|
887 | ! write(lunout,*) 'Ce traceur est un isotope' |
---|
888 | nb_iso(ixt,phase)=nb_iso(ixt,phase)+1 |
---|
889 | nb_isoind(phase)=nb_isoind(phase)+1 |
---|
890 | iso_num(iq)=ixt |
---|
891 | iso_indnum(iq)=nb_isoind(phase) |
---|
892 | indnum_fn_num(ixt)=iso_indnum(iq) |
---|
893 | phase_num(iq)=phase |
---|
894 | ! write(lunout,*) 'iso_num(iq)=',iso_num(iq) |
---|
895 | ! write(lunout,*) 'iso_indnum(iq)=',iso_indnum(iq) |
---|
896 | ! write(lunout,*) 'indnum_fn_num(ixt)=',indnum_fn_num(ixt) |
---|
897 | ! write(lunout,*) 'phase_num(iq)=',phase_num(iq) |
---|
898 | goto 20 |
---|
899 | else if (iqpere(iq).gt.0) then |
---|
900 | if (tnom_0(iqpere(iq)) == tnom_trac) then |
---|
901 | ! write(lunout,*) 'Ce traceur est le fils d''un isotope' |
---|
902 | ! c'est un traceur d'isotope |
---|
903 | nb_traciso(ixt,phase)=nb_traciso(ixt,phase)+1 |
---|
904 | iso_num(iq)=ixt |
---|
905 | iso_indnum(iq)=indnum_fn_num(ixt) |
---|
906 | zone_num(iq)=nb_traciso(ixt,phase) |
---|
907 | phase_num(iq)=phase |
---|
908 | ! write(lunout,*) 'iso_num(iq)=',iso_num(iq) |
---|
909 | ! write(lunout,*) 'phase_num(iq)=',phase_num(iq) |
---|
910 | ! write(lunout,*) 'zone_num(iq)=',zone_num(iq) |
---|
911 | goto 20 |
---|
912 | endif !if (tnom_0(iqpere(iq)) == trim(tnom_0(phase))//trim(tnom_iso(ixt))) then |
---|
913 | endif !IF (tnom_0(iq) == trim(tnom_0(phase))//trim(tnom_iso(ixt))) then |
---|
914 | enddo !do ixt= niso_possibles |
---|
915 | enddo !do phase=1,nqo |
---|
916 | 20 continue |
---|
917 | enddo !do iq=1,nqtot |
---|
918 | |
---|
919 | ! write(lunout,*) 'iso_num=',iso_num |
---|
920 | ! write(lunout,*) 'iso_indnum=',iso_indnum |
---|
921 | ! write(lunout,*) 'zone_num=',zone_num |
---|
922 | ! write(lunout,*) 'phase_num=',phase_num |
---|
923 | ! write(lunout,*) 'indnum_fn_num=',indnum_fn_num |
---|
924 | |
---|
925 | do ixt= 1,niso_possibles |
---|
926 | |
---|
927 | if (nb_iso(ixt,1).eq.1) then |
---|
928 | ! on verifie que toutes les phases ont le meme nombre de |
---|
929 | ! traceurs |
---|
930 | do phase=2,nqo |
---|
931 | if (nb_iso(ixt,phase).ne.nb_iso(ixt,1)) then |
---|
932 | ! write(lunout,*) 'ixt,phase,nb_iso=',ixt,phase,nb_iso(ixt,phase) |
---|
933 | CALL abort_gcm('infotrac_init','Phases must have same number of isotopes',1) |
---|
934 | endif |
---|
935 | enddo !do phase=2,nqo |
---|
936 | |
---|
937 | niso=niso+1 |
---|
938 | use_iso(ixt)=.true. |
---|
939 | ntraceurs_zone=nb_traciso(ixt,1) |
---|
940 | |
---|
941 | ! on verifie que toutes les phases ont le meme nombre de |
---|
942 | ! traceurs |
---|
943 | do phase=2,nqo |
---|
944 | if (nb_traciso(ixt,phase).ne.ntraceurs_zone) then |
---|
945 | write(lunout,*) 'ixt,phase,nb_traciso=',ixt,phase,nb_traciso(ixt,phase) |
---|
946 | write(lunout,*) 'ntraceurs_zone=',ntraceurs_zone |
---|
947 | CALL abort_gcm('infotrac_init','Phases must have same number of tracers',1) |
---|
948 | endif |
---|
949 | enddo !do phase=2,nqo |
---|
950 | ! on verifie que tous les isotopes ont le meme nombre de |
---|
951 | ! traceurs |
---|
952 | if (ntraceurs_zone_prec.gt.0) then |
---|
953 | if (ntraceurs_zone.eq.ntraceurs_zone_prec) then |
---|
954 | ntraceurs_zone_prec=ntraceurs_zone |
---|
955 | else !if (ntraceurs_zone.eq.ntraceurs_zone_prec) then |
---|
956 | write(*,*) 'ntraceurs_zone_prec,ntraceurs_zone=',ntraceurs_zone_prec,ntraceurs_zone |
---|
957 | CALL abort_gcm('infotrac_init', & |
---|
958 | &'Isotope tracers are not well defined in traceur.def',1) |
---|
959 | endif !if (ntraceurs_zone.eq.ntraceurs_zone_prec) then |
---|
960 | endif !if (ntraceurs_zone_prec.gt.0) then |
---|
961 | |
---|
962 | else if (nb_iso(ixt,1).ne.0) then |
---|
963 | WRITE(lunout,*) 'nqo,ixt=',nqo,ixt |
---|
964 | WRITE(lunout,*) 'nb_iso(ixt,1)=',nb_iso(ixt,1) |
---|
965 | CALL abort_gcm('infotrac_init','Isotopes are not well defined in traceur.def',1) |
---|
966 | endif !if (nb_iso(ixt,1).eq.1) then |
---|
967 | enddo ! do ixt= niso_possibles |
---|
968 | |
---|
969 | ! dimensions isotopique: |
---|
970 | ntraciso=niso*(ntraceurs_zone+1) |
---|
971 | ! WRITE(lunout,*) 'niso=',niso |
---|
972 | ! WRITE(lunout,*) 'ntraceurs_zone,ntraciso=',ntraceurs_zone,ntraciso |
---|
973 | |
---|
974 | ! flags isotopiques: |
---|
975 | if (niso.gt.0) then |
---|
976 | ok_isotopes=.true. |
---|
977 | else |
---|
978 | ok_isotopes=.false. |
---|
979 | endif |
---|
980 | ! WRITE(lunout,*) 'ok_isotopes=',ok_isotopes |
---|
981 | |
---|
982 | if (ok_isotopes) then |
---|
983 | ok_iso_verif=.false. |
---|
984 | call getin('ok_iso_verif',ok_iso_verif) |
---|
985 | ok_init_iso=.false. |
---|
986 | call getin('ok_init_iso',ok_init_iso) |
---|
987 | tnat=(/1.0,155.76e-6,2005.2e-6,0.004/100.,0.0/) |
---|
988 | alpha_ideal=(/1.0,1.01,1.006,1.003,1.0/) |
---|
989 | endif !if (ok_isotopes) then |
---|
990 | ! WRITE(lunout,*) 'ok_iso_verif=',ok_iso_verif |
---|
991 | ! WRITE(lunout,*) 'ok_init_iso=',ok_init_iso |
---|
992 | |
---|
993 | if (ntraceurs_zone.gt.0) then |
---|
994 | ok_isotrac=.true. |
---|
995 | else |
---|
996 | ok_isotrac=.false. |
---|
997 | endif |
---|
998 | ! WRITE(lunout,*) 'ok_isotrac=',ok_isotrac |
---|
999 | |
---|
1000 | ! remplissage du tableau iqiso(ntraciso,phase) |
---|
1001 | ALLOCATE(iqiso(ntraciso,nqo)) |
---|
1002 | iqiso(:,:)=0 |
---|
1003 | do iq=1,nqtot |
---|
1004 | if (iso_num(iq).gt.0) then |
---|
1005 | ixt=iso_indnum(iq)+zone_num(iq)*niso |
---|
1006 | iqiso(ixt,phase_num(iq))=iq |
---|
1007 | endif |
---|
1008 | enddo |
---|
1009 | ! WRITE(lunout,*) 'iqiso=',iqiso |
---|
1010 | |
---|
1011 | ! replissage du tableau index_trac(ntraceurs_zone,niso) |
---|
1012 | ALLOCATE(index_trac(ntraceurs_zone,niso)) |
---|
1013 | if (ok_isotrac) then |
---|
1014 | do iiso=1,niso |
---|
1015 | do izone=1,ntraceurs_zone |
---|
1016 | index_trac(izone,iiso)=iiso+izone*niso |
---|
1017 | enddo |
---|
1018 | enddo |
---|
1019 | else !if (ok_isotrac) then |
---|
1020 | index_trac(:,:)=0.0 |
---|
1021 | endif !if (ok_isotrac) then |
---|
1022 | ! write(lunout,*) 'index_trac=',index_trac |
---|
1023 | |
---|
1024 | ! Finalize : |
---|
1025 | DEALLOCATE(nb_iso) |
---|
1026 | |
---|
1027 | END SUBROUTINE infotrac_isoinit |
---|
1028 | |
---|
1029 | END MODULE infotrac |
---|