source: trunk/LMDZ.COMMON/libf/dyn3d_common/infotrac.F90 @ 2572

Last change on this file since 2572 was 2468, checked in by yjaziri, 4 years ago
File size: 46.2 KB
Line 
1! $Id$
2!
3MODULE infotrac
4
5! nqtot : total number of tracers and higher order of moment, water vapor and liquid included
6  INTEGER, SAVE :: nqtot
7! CR: add number of tracers for water (for Earth model only!!)
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: nb of father tracers (i.e. directly advected by air)
15  INTEGER, SAVE :: nqperes
16
17! Name variables
18  CHARACTER(len=30), ALLOCATABLE, DIMENSION(:), SAVE :: tname ! tracer short name for restart and diagnostics
19  CHARACTER(len=33), ALLOCATABLE, DIMENSION(:), SAVE :: ttext ! tracer long name for diagnostics
20
21! iadv  : index of trasport schema for each tracer
22  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE    :: iadv
23
24! niadv : vector keeping the coorspondance between all tracers(nqtot) treated in the
25!         dynamic part of the code and the tracers (nbtr+2) used in the physics part of the code.
26  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE    :: niadv ! equivalent dyn / physique
27
28! CRisi: arrays for sons
29  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE    :: nqfils
30  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE    :: nqdesc ! number of sons + all gran-sons over all generations
31  INTEGER, SAVE :: nqdesc_tot
32  INTEGER, ALLOCATABLE, DIMENSION(:,:), SAVE    :: iqfils
33  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE    :: iqpere
34  REAL :: qperemin,masseqmin ! MVals: thresholds for zq(pere) and masseq in the transport of the isotopic Ratio (vlsplt_p.F)
35  PARAMETER (qperemin=1.e-16)
36  PARAMETER (masseqmin=1.e-16)
37
38! conv_flg(it)=0 : convection desactivated for tracer number it
39  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE  :: conv_flg
40! pbl_flg(it)=0  : boundary layer diffusion desactivaded for tracer number it
41  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE  :: pbl_flg
42
43  CHARACTER(len=4),SAVE :: type_trac
44  CHARACTER(len=8),DIMENSION(:),ALLOCATABLE, SAVE :: solsym
45
46    ! CRisi: specific stuff for isotopes
47    LOGICAL,SAVE :: ok_isotopes,ok_iso_verif,ok_isotrac,ok_init_iso
48    INTEGER :: niso_possibles   
49    PARAMETER ( niso_possibles=5) ! 5 possible water isotopes
50    REAL, DIMENSION (niso_possibles),SAVE :: tnat,alpha_ideal
51    LOGICAL, DIMENSION(niso_possibles),SAVE ::  use_iso
52    INTEGER, ALLOCATABLE, DIMENSION(:,:), SAVE ::  iqiso ! donne indice iq en fn de (ixt,phase)
53    INTEGER, ALLOCATABLE, DIMENSION(:), SAVE ::  iso_num ! donne numéro iso entre 1 et niso_possibles en fn de nqtot
54    INTEGER, ALLOCATABLE, DIMENSION(:), SAVE ::  iso_indnum ! donne numéro iso entre 1 et niso effectif en fn de nqtot
55    INTEGER, ALLOCATABLE, DIMENSION(:), SAVE ::  zone_num ! donne numéro de la zone de tracage en fn de nqtot
56    INTEGER, ALLOCATABLE, DIMENSION(:), SAVE ::  phase_num ! donne numéro de la zone de tracage en fn de nqtot
57    INTEGER, DIMENSION(niso_possibles), SAVE :: indnum_fn_num ! donne indice entre entre 1 et niso en fonction du numéro d isotope entre 1 et niso_possibles
58    INTEGER, ALLOCATABLE, DIMENSION(:,:), SAVE ::  index_trac ! numéro ixt en fn izone, indnum entre 1 et niso
59    INTEGER,SAVE :: niso,ntraceurs_zone,ntraciso
60
61#ifdef CPP_StratAer
62!--CK/OB for stratospheric aerosols
63  INTEGER, SAVE :: nbtr_bin
64  INTEGER, SAVE :: nbtr_sulgas
65  INTEGER, SAVE :: id_OCS_strat
66  INTEGER, SAVE :: id_SO2_strat
67  INTEGER, SAVE :: id_H2SO4_strat
68  INTEGER, SAVE :: id_BIN01_strat
69  INTEGER, SAVE :: id_TEST_strat
70#endif
71
72CONTAINS
73
74  SUBROUTINE infotrac_init
75    USE control_mod, ONLY: planet_type, config_inca
76#ifdef REPROBUS
77    USE CHEM_REP, ONLY : Init_chem_rep_trac
78#endif
79    IMPLICIT NONE
80!=======================================================================
81!
82!   Auteur:  P. Le Van /L. Fairhead/F.Hourdin
83!   -------
84!   Modif special traceur F.Forget 05/94
85!   Modif M-A Filiberti 02/02 lecture de traceur.def
86!
87!   Objet:
88!   ------
89!   GCM LMD nouvelle grille
90!
91!=======================================================================
92!   ... modification de l'integration de q ( 26/04/94 ) ....
93!-----------------------------------------------------------------------
94! Declarations
95
96    INCLUDE "dimensions.h"
97    INCLUDE "iniprint.h"
98
99! Local variables
100    INTEGER, ALLOCATABLE, DIMENSION(:) :: hadv  ! index of horizontal trasport schema
101    INTEGER, ALLOCATABLE, DIMENSION(:) :: vadv  ! index of vertical trasport schema
102
103    INTEGER, ALLOCATABLE, DIMENSION(:) :: hadv_inca  ! index of horizontal trasport schema
104    INTEGER, ALLOCATABLE, DIMENSION(:) :: vadv_inca  ! index of vertical trasport schema
105
106    CHARACTER(len=30), ALLOCATABLE, DIMENSION(:) :: tnom_0  ! tracer short name
107    CHARACTER(len=30), ALLOCATABLE, DIMENSION(:) :: tnom_transp ! transporting fluid short name: CRisi
108    CHARACTER(len=3), DIMENSION(30) :: descrq
109    CHARACTER(len=1), DIMENSION(3)  :: txts
110    CHARACTER(len=2), DIMENSION(9)  :: txtp
111    CHARACTER(len=30)               :: str1,str2
112 
113    INTEGER :: nqtrue  ! number of tracers read from tracer.def, without higer order of moment
114    INTEGER :: iq, new_iq, iiq, jq, ierr, ierr1, ierr2, ierr3, ierr4 !MVals: add ierr1 et ierr4 for the isotope case
115    INTEGER :: ifils,ipere,generation ! CRisi
116    LOGICAL :: continu,nouveau_traceurdef
117    LOGICAL :: moderntracdef=.false. ! JVO, YJ : For planets modern traceur.def
118    INTEGER :: IOstatus ! gestion de la retrocompatibilite de traceur.def
119    CHARACTER(len=30) :: tchaine   
120   
121    character(len=500) :: line ! to store a line of text
122 
123    character(len=*),parameter :: modname="infotrac_init"
124!-----------------------------------------------------------------------
125! Initialization :
126!
127    txts=(/'x','y','z'/)
128    txtp=(/'x ','y ','z ','xx','xy','xz','yy','yz','zz'/)
129
130    descrq(14)='VLH'
131    descrq(10)='VL1'
132    descrq(11)='VLP'
133    descrq(12)='FH1'
134    descrq(13)='FH2'
135    descrq(16)='PPM'
136    descrq(17)='PPS'
137    descrq(18)='PPP'
138    descrq(20)='SLP'
139    descrq(30)='PRA'
140   
141    IF (planet_type=='earth') THEN
142     ! Coherence test between parameter type_trac, config_inca and preprocessing keys
143     IF (type_trac=='inca') THEN
144       WRITE(lunout,*) 'You have choosen to couple with INCA chemestry model : type_trac=', &
145            type_trac,' config_inca=',config_inca
146       IF (config_inca/='aero' .AND. config_inca/='aeNP' .AND. config_inca/='chem') THEN
147          WRITE(lunout,*) 'Incoherence between type_trac and config_inca. Model stops. Modify run.def'
148          CALL abort_gcm('infotrac_init','Incoherence between type_trac and config_inca',1)
149       END IF
150#ifndef INCA
151       WRITE(lunout,*) 'To run this option you must add cpp key INCA and compile with INCA code'
152       CALL abort_gcm('infotrac_init','You must compile with cpp key INCA',1)
153#endif
154     ELSE IF (type_trac=='repr') THEN
155       WRITE(lunout,*) 'You have choosen to couple with REPROBUS chemestry model : type_trac=', type_trac
156#ifndef REPROBUS
157       WRITE(lunout,*) 'To run this option you must add cpp key REPROBUS and compile with REPRPBUS code'
158       CALL abort_gcm('infotrac_init','You must compile with cpp key REPROBUS',1)
159#endif
160    ELSE IF (type_trac == 'coag') THEN
161       WRITE(lunout,*) 'Tracers are treated for COAGULATION tests : type_trac=', type_trac
162#ifndef CPP_StratAer
163       WRITE(lunout,*) 'To run this option you must add cpp key StratAer and compile with StratAer code'
164       CALL abort_gcm('infotrac_init','You must compile with cpp key StratAer',1)
165#endif
166     ELSE IF (type_trac == 'lmdz') THEN
167       WRITE(lunout,*) 'Tracers are treated in LMDZ only : type_trac=', type_trac
168     ELSE
169       WRITE(lunout,*) 'type_trac=',type_trac,' not possible. Model stops'
170       CALL abort_gcm('infotrac_init','bad parameter',1)
171     END IF
172
173     ! Test if config_inca is other then none for run without INCA
174     IF (type_trac/='inca' .AND. config_inca/='none') THEN
175       WRITE(lunout,*) 'config_inca will now be changed to none as you do not couple with INCA model'
176       config_inca='none'
177     END IF
178    ELSE
179     type_trac='plnt'  ! planets... May want to dissociate between each later.
180    ENDIF ! of IF (planet_type=='earth')
181
182!-----------------------------------------------------------------------
183!
184! 1) Get the true number of tracers + water vapor/liquid
185!    Here true tracers (nqtrue) means declared tracers (only first order)
186!
187!-----------------------------------------------------------------------
188    IF (planet_type=='earth') THEN
189     IF (type_trac == 'lmdz' .OR. type_trac == 'repr' .OR. type_trac == 'coag') THEN
190       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
191       IF(ierr.EQ.0) THEN
192          WRITE(lunout,*) trim(modname),': Open traceur.def : ok'
193          READ(90,*) nqtrue
194          WRITE(lunout,*) trim(modname),' nqtrue=',nqtrue
195       ELSE
196          WRITE(lunout,*) trim(modname),': Problem in opening traceur.def'
197          WRITE(lunout,*) trim(modname),': WARNING using defaut values'
198          nqtrue=4 ! Defaut value
199       END IF
200       ! For Earth, water vapour & liquid tracers are not in the physics
201       nbtr=nqtrue-2
202     ELSE ! type_trac=inca
203       ! The traceur.def file is used to define the number "nqo" of water phases
204       ! present in the simulation. Default : nqo = 2.
205       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
206       IF(ierr.EQ.0) THEN
207          WRITE(lunout,*) trim(modname),': Open traceur.def : ok'
208          READ(90,*) nqo
209       ELSE
210          WRITE(lunout,*) trim(modname),': Using default value for nqo'
211          nqo=2
212       ENDIF
213       IF (nqo /= 2 .AND. nqo /= 3 ) THEN
214          WRITE(lunout,*) trim(modname),': nqo=',nqo, ' is not allowded. Only 2 or 3 water phases allowed'
215          CALL abort_gcm('infotrac_init','Bad number of water phases',1)
216       END IF
217       ! nbtr has been read from INCA by init_const_lmdz() in gcm.F
218#ifdef INCA
219       CALL Init_chem_inca_trac(nbtr)
220#endif
221       nqtrue=nbtr+nqo
222
223       ALLOCATE(hadv_inca(nbtr), vadv_inca(nbtr))
224
225     END IF   ! type_trac
226
227     IF (nqtrue < 2) THEN
228       WRITE(lunout,*) trim(modname),': nqtrue=',nqtrue, ' is not allowded. 2 tracers is the minimum'
229       CALL abort_gcm('infotrac_init','Not enough tracers',1)
230     END IF
231
232!jyg<
233! Transfert number of tracers to Reprobus
234!!    IF (type_trac == 'repr') THEN
235!!#ifdef REPROBUS
236!!       CALL Init_chem_rep_trac(nbtr)
237!!#endif
238!!    END IF
239!>jyg
240
241    ELSE  ! not Earth
242       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
243       IF (ierr.EQ.0) THEN
244          WRITE(lunout,*) 'Open traceur.def : ok'
245!! - Modif. by JVO and YJ for modern planetary traceur.def ---------------
246          READ(90,'(A)') line
247          IF (trim(line).ne.'#ModernTrac-v1') THEN ! Test modern traceur.def
248             READ(line,*) nqtrue ! Try standard traceur.def
249          ELSE
250             moderntracdef = .true.
251             DO
252               READ(90,'(A)',iostat=ierr1) line
253               IF (ierr1==0) THEN
254                 IF (index(line,'#').ne.1) THEN ! Allows arbitary number of comments lines in the header
255                   READ(line,*) nbtr
256                   WRITE(lunout,*) trim(modname),': Number of tracers in the physics: nbtr=',nbtr
257                   EXIT
258                   ! JVO 20 : Here nbtr corresponds to ALL what we consider as modern 'tracers' (species)
259                   ! from physics, including non-dynamic species (radiative, compo..), and nqtrue will be set
260                   ! below, according to the flags in traceur.def, by read and rewind.
261                 ENDIF
262               ELSE ! If pb, or if reached EOF without having found nbtr
263                 CALL abort_gcm('infotrac_init','Unable to read numbers of tracers in traceur.def',1)
264               ENDIF
265             ENDDO
266          ENDIF ! if modern or standard traceur.def
267!! -----------------------------------------------------------------------
268       ELSE
269          WRITE(lunout,*) 'Problem in opening traceur.def'
270          WRITE(lunout,*) 'ATTENTION using defaut values: nqtrue=1'
271          nqtrue=1 ! Defaut value
272       END IF
273
274       nqo=0
275     
276    ENDIF  ! planet_type
277!
278! Allocate variables depending on nqtrue and nbtr
279!
280    IF (.NOT.moderntracdef) ALLOCATE(tnom_0(nqtrue), hadv(nqtrue), vadv(nqtrue),tnom_transp(nqtrue))
281!
282!jyg<
283!!    ALLOCATE(conv_flg(nbtr), pbl_flg(nbtr), solsym(nbtr))
284!!    conv_flg(:) = 1 ! convection activated for all tracers
285!!    pbl_flg(:)  = 1 ! boundary layer activated for all tracers
286!>jyg
287
288!-----------------------------------------------------------------------
289! 2)     Choix  des schemas d'advection pour l'eau et les traceurs
290!
291!     iadv = 1    schema  transport type "humidite specifique LMD"
292!     iadv = 2    schema   amont
293!     iadv = 14   schema  Van-leer + humidite specifique
294!                            Modif F.Codron
295!     iadv = 10   schema  Van-leer (retenu pour l'eau vapeur et liquide)
296!     iadv = 11   schema  Van-Leer pour hadv et version PPM (Monotone) pour vadv
297!     iadv = 12   schema  Frederic Hourdin I
298!     iadv = 13   schema  Frederic Hourdin II
299!     iadv = 16   schema  PPM Monotone(Collela & Woodward 1984)
300!     iadv = 17   schema  PPM Semi Monotone (overshoots autorisés)
301!     iadv = 18   schema  PPM Positif Defini (overshoots undershoots autorisés)
302!     iadv = 20   schema  Slopes
303!     iadv = 30   schema  Prather
304!
305!        Dans le tableau q(ij,l,iq) : iq = 1  pour l'eau vapeur
306!                                     iq = 2  pour l'eau liquide
307!       Et eventuellement             iq = 3,nqtot pour les autres traceurs
308!
309!        iadv(1): choix pour l'eau vap. et  iadv(2) : choix pour l'eau liq.
310!------------------------------------------------------------------------
311!
312!    Get choice of advection schema from file tracer.def or from INCA
313!---------------------------------------------------------------------
314    IF (planet_type=='earth') THEN
315     IF (type_trac == 'lmdz' .OR. type_trac == 'repr' .OR. type_trac == 'coag') THEN
316       IF(ierr.EQ.0) THEN
317          ! Continue to read tracer.def
318          DO iq=1,nqtrue
319
320             write(*,*) 'infotrac 237: iq=',iq
321             ! CRisi: ajout du nom du fluide transporteur
322             ! mais rester retro compatible
323             READ(90,'(I2,X,I2,X,A)',IOSTAT=IOstatus) hadv(iq),vadv(iq),tchaine
324             write(lunout,*) 'iq,hadv(iq),vadv(iq)=',iq,hadv(iq),vadv(iq)
325             write(lunout,*) 'tchaine=',trim(tchaine)
326!             write(*,*) 'infotrac 238: IOstatus=',IOstatus
327             if (IOstatus.ne.0) then
328                CALL abort_gcm('infotrac_init','Pb dans la lecture de traceur.def',1)
329             endif
330             ! Y-a-t-il 1 ou 2 noms de traceurs? -> On regarde s'il y a un
331             ! espace ou pas au milieu de la chaine.
332             continu=.true.
333             nouveau_traceurdef=.false.
334             iiq=1
335             do while (continu)
336                if (tchaine(iiq:iiq).eq.' ') then
337                  nouveau_traceurdef=.true.
338                  continu=.false.
339                else if (iiq.lt.LEN_TRIM(tchaine)) then
340                  iiq=iiq+1
341                else
342                  continu=.false.     
343                endif
344             enddo
345             write(*,*) 'iiq,nouveau_traceurdef=',iiq,nouveau_traceurdef
346             if (nouveau_traceurdef) then
347                write(lunout,*) 'C''est la nouvelle version de traceur.def'
348                tnom_0(iq)=tchaine(1:iiq-1)
349                tnom_transp(iq)=tchaine(iiq+1:15)
350             else
351                write(lunout,*) 'C''est l''ancienne version de traceur.def'
352                write(lunout,*) 'On suppose que les traceurs sont tous d''air'
353                tnom_0(iq)=tchaine
354                tnom_transp(iq) = 'air'
355             endif
356             write(lunout,*) 'tnom_0(iq)=<',trim(tnom_0(iq)),'>'
357             write(lunout,*) 'tnom_transp(iq)=<',trim(tnom_transp(iq)),'>'
358
359          END DO ! DO iq=1,nqtrue
360          CLOSE(90) 
361
362       ELSE ! Without tracer.def, set default values (for Earth!)
363         if ((nqtrue==4).and.(planet_type=="earth")) then
364          hadv(1) = 14
365          vadv(1) = 14
366          tnom_0(1) = 'H2Ov'
367          tnom_transp(1) = 'air'
368          hadv(2) = 10
369          vadv(2) = 10
370          tnom_0(2) = 'H2Ol'
371          tnom_transp(2) = 'air'
372          hadv(3) = 10
373          vadv(3) = 10
374          tnom_0(3) = 'RN'
375          tnom_transp(3) = 'air'
376          hadv(4) = 10
377          vadv(4) = 10
378          tnom_0(4) = 'PB'
379          tnom_transp(4) = 'air'
380         else
381           ! Error message, we need a traceur.def file
382           write(lunout,*) trim(modname),&
383           ': Cannot set default tracer names!'
384           write(lunout,*) trim(modname),' Make a traceur.def file!!!'
385           CALL abort_gcm('infotrac_init','Need a traceur.def file!',1)
386         endif ! of if (nqtrue==4)
387       END IF ! of IF(ierr.EQ.0)
388       
389       WRITE(lunout,*) trim(modname),': Valeur de traceur.def :'
390       WRITE(lunout,*) trim(modname),': nombre de traceurs ',nqtrue
391       DO iq=1,nqtrue
392          WRITE(lunout,*) hadv(iq),vadv(iq),tnom_0(iq),tnom_transp(iq)
393       END DO
394
395!       IF ( planet_type=='earth') THEN
396         !CR: nombre de traceurs de l eau
397         IF (tnom_0(3) == 'H2Oi') then
398            nqo=3
399         ELSE
400            nqo=2
401         ENDIF
402         ! For Earth, water vapour & liquid tracers are not in the physics
403         nbtr=nqtrue-nqo
404!       ELSE
405!         ! Other planets (for now); we have the same number of tracers
406!         ! in the dynamics than in the physics
407!         nbtr=nqtrue
408!       ENDIF
409
410#ifdef CPP_StratAer
411       IF (type_trac == 'coag') THEN
412         nbtr_bin=0
413         nbtr_sulgas=0
414         DO iq=1,nqtrue
415           IF (tnom_0(iq)(1:3)=='BIN') THEN !check if tracer name contains 'BIN'
416             nbtr_bin=nbtr_bin+1
417           ENDIF
418           IF (tnom_0(iq)(1:3)=='GAS') THEN !check if tracer name contains 'GAS'
419             nbtr_sulgas=nbtr_sulgas+1
420           ENDIF
421         ENDDO
422         print*,'nbtr_bin=',nbtr_bin
423         print*,'nbtr_sulgas=',nbtr_sulgas
424         DO iq=1,nqtrue
425           IF (tnom_0(iq)=='GASOCS') THEN
426             id_OCS_strat=iq-nqo
427           ENDIF
428           IF (tnom_0(iq)=='GASSO2') THEN
429             id_SO2_strat=iq-nqo
430           ENDIF
431           IF (tnom_0(iq)=='GASH2SO4') THEN
432             id_H2SO4_strat=iq-nqo
433           ENDIF
434           IF (tnom_0(iq)=='BIN01') THEN
435             id_BIN01_strat=iq-nqo
436           ENDIF
437           IF (tnom_0(iq)=='GASTEST') THEN
438             id_TEST_strat=iq-nqo
439           ENDIF
440         ENDDO
441         print*,'id_OCS_strat  =',id_OCS_strat
442         print*,'id_SO2_strat  =',id_SO2_strat
443         print*,'id_H2SO4_strat=',id_H2SO4_strat
444         print*,'id_BIN01_strat=',id_BIN01_strat
445       ENDIF
446#endif
447
448     ENDIF  ! (type_trac == 'lmdz' .OR. type_trac == 'repr' .OR. type_trac = 'coag')
449!jyg<
450!
451! Transfert number of tracers to Reprobus
452    IF (type_trac == 'repr') THEN
453#ifdef REPROBUS
454       CALL Init_chem_rep_trac(nbtr)
455#endif
456    END IF
457!
458! Allocate variables depending on nbtr
459!
460    ALLOCATE(conv_flg(nbtr), pbl_flg(nbtr), solsym(nbtr))
461    conv_flg(:) = 1 ! convection activated for all tracers
462    pbl_flg(:)  = 1 ! boundary layer activated for all tracers
463!
464!!    ELSE  ! type_trac=inca : config_inca='aero' ou 'chem'
465!
466    IF (type_trac == 'inca') THEN   ! config_inca='aero' ou 'chem'
467!>jyg
468! le module de chimie fournit les noms des traceurs
469! et les schemas d'advection associes. excepte pour ceux lus
470! dans traceur.def
471       IF (ierr .eq. 0) then
472          DO iq=1,nqo
473
474             write(*,*) 'infotrac 237: iq=',iq
475             ! CRisi: ajout du nom du fluide transporteur
476             ! mais rester retro compatible
477             READ(90,'(I2,X,I2,X,A)',IOSTAT=IOstatus) hadv(iq),vadv(iq),tchaine
478             write(lunout,*) 'iq,hadv(iq),vadv(iq)=',iq,hadv(iq),vadv(iq)
479             write(lunout,*) 'tchaine=',trim(tchaine)
480             write(*,*) 'infotrac 238: IOstatus=',IOstatus
481             if (IOstatus.ne.0) then
482                CALL abort_gcm('infotrac_init','Pb dans la lecture de traceur.def',1)
483             endif
484             ! Y-a-t-il 1 ou 2 noms de traceurs? -> On regarde s'il y a un
485             ! espace ou pas au milieu de la chaine.
486             continu=.true.
487             nouveau_traceurdef=.false.
488             iiq=1
489             do while (continu)
490                if (tchaine(iiq:iiq).eq.' ') then
491                  nouveau_traceurdef=.true.
492                  continu=.false.
493                else if (iiq.lt.LEN_TRIM(tchaine)) then
494                  iiq=iiq+1
495                else
496                  continu=.false.
497                endif
498             enddo
499             write(*,*) 'iiq,nouveau_traceurdef=',iiq,nouveau_traceurdef
500             if (nouveau_traceurdef) then
501                write(lunout,*) 'C''est la nouvelle version de traceur.def'
502                tnom_0(iq)=tchaine(1:iiq-1)
503                tnom_transp(iq)=tchaine(iiq+1:15)
504             else
505                write(lunout,*) 'C''est l''ancienne version de traceur.def'
506                write(lunout,*) 'On suppose que les traceurs sont tous d''air'
507                tnom_0(iq)=tchaine
508                tnom_transp(iq) = 'air'
509             endif
510             write(lunout,*) 'tnom_0(iq)=<',trim(tnom_0(iq)),'>'
511             write(lunout,*) 'tnom_transp(iq)=<',trim(tnom_transp(iq)),'>'
512
513          END DO !DO iq=1,nqtrue
514          CLOSE(90) 
515       ELSE  !! if traceur.def doesn't exist
516          tnom_0(1)='H2Ov'
517          tnom_transp(1) = 'air'
518          tnom_0(2)='H2Ol'
519          tnom_transp(2) = 'air'
520          hadv(1) = 10
521          hadv(2) = 10
522          vadv(1) = 10
523          vadv(2) = 10
524       ENDIF
525     
526#ifdef INCA
527       CALL init_transport( &
528            hadv_inca, &
529            vadv_inca, &
530            conv_flg, &
531            pbl_flg,  &
532            tracnam)
533#endif
534
535!jyg<
536       DO iq = nqo+1, nqtrue
537          tnom_0(iq)=solsym(iq-nqo)
538          tnom_transp(iq) = 'air'
539       END DO
540!!       DO iq =3,nqtrue
541!!          tnom_0(iq)=solsym(iq-2)
542!!       END DO
543!!       nqo = 2
544!>jyg
545
546     END IF ! (type_trac == 'inca')
547
548    ELSE  ! not Earth
549
550       IF (.not.moderntracdef) THEN ! Standard traceur.def or no traceur.def file
551
552         ! Other planets (for now); we have the same number of tracers
553         ! in the dynamics than in the physics
554         nbtr=nqtrue
555         ! NB: Reading a traceur.def with isotopes remains to be done...
556         IF(ierr.EQ.0) THEN
557            ! Continue to read tracer.def
558            DO iq=1,nqtrue
559               !READ(90,*) hadv(iq),vadv(iq),tnom_0(iq)
560              ! try to be smart when reading traceur.def
561              read(90,'(80a)') line ! store the line from traceur.def
562              ! if format is hadv, vadv, tnom_0, tnom_transp
563              read(line,*,iostat=ierr1) hadv(iq),vadv(iq),tnom_0(iq),tnom_transp(iq)
564              if (ierr1.ne.0) then
565              ! assume format is hadv,vadv,tnom_0
566               read(line,*,iostat=ierr2) hadv(iq),vadv(iq),tnom_0(iq)
567               if (ierr2.ne.0) then
568                ! maybe format is tnom0,hadv,vadv
569                read(line,*,iostat=ierr3) tnom_0(iq),hadv(iq),vadv(iq)
570                if (ierr3.ne.0) then
571                  !assuming values of hadv et vadv
572                  hadv(iq)=10
573                  vadv(iq)=10
574                  read(line,*, iostat=ierr4) tnom_0(iq), tnom_transp(iq)
575                  if (ierr4.ne.0) then
576                  ! assume only tnom0 is provided (hadv and vadv default to 10)
577                    read(line,*) tnom_0(iq)
578                    tnom_transp(iq)='air'
579                  endif
580                else
581                  !format is tnom0,hadv,vadv
582                  tnom_transp(iq)='air' ! no isotopes
583                endif ! of if (ierr3.ne.0)
584               else
585                 ! format is hadv,vadv,tnom_0
586                 tnom_transp(iq)='air' ! no isotopes
587               endif ! of if (ierr2.ne.0)
588              endif ! of if(ierr1.ne.0)
589            END DO ! of DO iq=1,nqtrue
590            CLOSE(90) 
591         ELSE ! Without tracer.def
592            hadv(1) = 10
593            vadv(1) = 10
594            tnom_0(1) = 'dummy'
595            tnom_transp(1)='air'
596         END IF
597       
598       ELSE ! Modern traceur.def (moderntracdef=true) - Authors : JVO, YJ - 2020
599
600         nqtrue = 0
601         DO iq=1,nbtr ! This loops on ALL the 'tracers' including the non-dyn ones.
602           READ(90,'(A)') line ! store the line from traceur.def
603           IF (index(line,'is_dyn=0').ne.0) CYCLE ! Only if explicitly not dyn., otherwise default behaviour
604           nqtrue = nqtrue + 1
605         ENDDO
606         WRITE(lunout,*) trim(modname),': Number of tracers in the dynamics: nqtrue=',nqtrue
607         ALLOCATE(tnom_0(nqtrue), hadv(nqtrue), vadv(nqtrue),tnom_transp(nqtrue))
608
609         ! Now that we know the dyn. tracers and allocated variables, let's read again the file.
610         REWIND(90)
611         DO
612           READ(90,'(A)') line ! no need to check on iostat, already done at first read
613           IF (index(line,'#').ne.1) THEN ! Skip arbitary number of commented lines in the header
614             READ(line,*) ! This is the line of nbtr (see above)
615             EXIT
616           ENDIF
617         ENDDO
618         iiq=0
619         DO iq=1,nbtr ! This loops on ALL the 'tracers' including the non-dyn ones.
620           
621           READ(90,'(A)') line ! store the line from traceur.def
622           IF (index(line,'is_dyn=0').ne.0) CYCLE ! Skip the non-dyn ones
623           iiq = iiq +  1 ! This ensure that iiq is on 1:nqtrue and not 1:nbtr
624           ! Name must be first parameter in this version, but all others can be in any order
625           read(line,*) tnom_0(iiq)
626           write(*,*)"infotrac_init: iq=",iiq,"noms(iq)=",trim(tnom_0(iiq))
627
628           if (index(line,'vadv=') /= 0) then
629             read(line(index(line,'vadv='   )+len('vadv='):),*) vadv(iiq)
630             write(*,*) ' Parameter value (traceur.def) : vadv=', vadv(iiq)
631           else
632             vadv(iiq)=10
633             write(*,*) ' Parameter value (default)     : vadv=',vadv(iiq)
634           end if
635
636           if (index(line,'hadv=') /= 0) then
637             read(line(index(line,'hadv='   )+len('hadv='):),*) hadv(iiq)
638             write(*,*) ' Parameter value (traceur.def) : hadv=', hadv(iiq)
639           else
640             hadv(iiq)=10
641             write(*,*) ' Parameter value (default)     : hadv=',hadv(iiq)
642           end if
643
644           if (index(line,'tnom_transp=' ) /= 0) then
645             read(line(index(line,'tnom_transp='   )+len('tnom_transp='):),*) tnom_transp(iiq)
646             write(*,*) ' Parameter value (traceur.def) : tnom_transp=', tnom_transp(iiq)
647           else
648             tnom_transp(iiq)='air'
649             write(*,*) ' Parameter value (default)     : tnom_transp=',tnom_transp(iiq)
650           end if
651
652         ENDDO
653         CLOSE(90)
654
655       ENDIF ! if (moderntracdef)
656
657       WRITE(lunout,*) trim(modname),': Valeur de traceur.def :'
658       WRITE(lunout,*) trim(modname),': nombre de traceurs ',nqtrue
659       DO iq=1,nqtrue
660          WRITE(lunout,*) hadv(iq),vadv(iq),tnom_0(iq),tnom_transp(iq)
661       END DO
662
663    ENDIF  ! planet_type
664
665!-----------------------------------------------------------------------
666!
667! 3) Verify if advection schema 20 or 30 choosen
668!    Calculate total number of tracers needed: nqtot
669!    Allocate variables depending on total number of tracers
670!-----------------------------------------------------------------------
671    new_iq=0
672    DO iq=1,nqtrue
673       ! Add tracers for certain advection schema
674       IF (hadv(iq)<20 .AND. vadv(iq)<20 ) THEN
675          new_iq=new_iq+1  ! no tracers added
676       ELSE IF (hadv(iq)==20 .AND. vadv(iq)==20 ) THEN
677          new_iq=new_iq+4  ! 3 tracers added
678       ELSE IF (hadv(iq)==30 .AND. vadv(iq)==30 ) THEN
679          new_iq=new_iq+10 ! 9 tracers added
680       ELSE
681          WRITE(lunout,*) trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq)
682          CALL abort_gcm('infotrac_init','Bad choice of advection schema - 1',1)
683       END IF
684    END DO
685   
686    IF (new_iq /= nqtrue) THEN
687       ! The choice of advection schema imposes more tracers
688       ! Assigne total number of tracers
689       nqtot = new_iq
690
691       WRITE(lunout,*) trim(modname),': The choice of advection schema for one or more tracers'
692       WRITE(lunout,*) 'makes it necessary to add tracers'
693       WRITE(lunout,*) trim(modname)//': ',nqtrue,' is the number of true tracers'
694       WRITE(lunout,*) trim(modname)//': ',nqtot, ' is the total number of tracers needed'
695
696    ELSE
697       ! The true number of tracers is also the total number
698       nqtot = nqtrue
699    END IF
700
701!
702! Allocate variables with total number of tracers, nqtot
703!
704    ALLOCATE(tname(nqtot), ttext(nqtot))
705    ALLOCATE(iadv(nqtot), niadv(nqtot))
706
707!-----------------------------------------------------------------------
708!
709! 4) Determine iadv, long and short name
710!
711!-----------------------------------------------------------------------
712    new_iq=0
713    DO iq=1,nqtrue
714       new_iq=new_iq+1
715
716       ! Verify choice of advection schema
717       IF (hadv(iq)==vadv(iq)) THEN
718          iadv(new_iq)=hadv(iq)
719       ELSE IF (hadv(iq)==10 .AND. vadv(iq)==16) THEN
720          iadv(new_iq)=11
721       ELSE
722          WRITE(lunout,*)trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq)
723
724          CALL abort_gcm('infotrac_init','Bad choice of advection schema - 2',1)
725       END IF
726     
727       str1=tnom_0(iq)
728       tname(new_iq)= tnom_0(iq)
729       IF (iadv(new_iq)==0) THEN
730          ttext(new_iq)=trim(str1)
731       ELSE
732          ttext(new_iq)=trim(tnom_0(iq))//descrq(iadv(new_iq))
733       END IF
734
735       ! schemas tenant compte des moments d'ordre superieur
736       str2=ttext(new_iq)
737       IF (iadv(new_iq)==20) THEN
738          DO jq=1,3
739             new_iq=new_iq+1
740             iadv(new_iq)=-20
741             ttext(new_iq)=trim(str2)//txts(jq)
742             tname(new_iq)=trim(str1)//txts(jq)
743          END DO
744       ELSE IF (iadv(new_iq)==30) THEN
745          DO jq=1,9
746             new_iq=new_iq+1
747             iadv(new_iq)=-30
748             ttext(new_iq)=trim(str2)//txtp(jq)
749             tname(new_iq)=trim(str1)//txtp(jq)
750          END DO
751       END IF
752    END DO
753
754!
755! Find vector keeping the correspodence between true and total tracers
756!
757    niadv(:)=0
758    iiq=0
759    DO iq=1,nqtot
760       IF(iadv(iq).GE.0) THEN
761          ! True tracer
762          iiq=iiq+1
763          niadv(iiq)=iq
764       ENDIF
765    END DO
766
767
768    WRITE(lunout,*) trim(modname),': Information stored in infotrac :'
769    WRITE(lunout,*) trim(modname),': iadv  niadv tname  ttext :'
770    DO iq=1,nqtot
771       WRITE(lunout,*) iadv(iq),niadv(iq),&
772       ' ',trim(tname(iq)),' ',trim(ttext(iq))
773    END DO
774
775!
776! Test for advection schema.
777! This version of LMDZ only garantees iadv=10 and iadv=14 (14 only for water vapour) .
778!
779    DO iq=1,nqtot
780       IF (iadv(iq)/=10 .AND. iadv(iq)/=14 .AND. iadv(iq)/=0) THEN
781          WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv(iq),' is not tested in this version of LMDZ'
782          CALL abort_gcm('infotrac_init','In this version only iadv=10 and iadv=14 is tested!',1)
783       ELSE IF (iadv(iq)==14 .AND. iq/=1) THEN
784          WRITE(lunout,*)trim(modname),'STOP : The option iadv=',iadv(iq),' is not tested in this version of LMDZ'
785          CALL abort_gcm('infotrac_init','In this version iadv=14 is only permitted for water vapour!',1)
786       END IF
787    END DO
788
789!-----------------------------------------------------------------------
790!
791! 5) Determine father/son relations for isotopes and carrying fluid
792!
793!-----------------------------------------------------------------------
794
795! CRisi: quels sont les traceurs fils et les traceurs pères.
796! initialiser tous les tableaux d'indices liés aux traceurs familiaux
797! + vérifier que tous les pères sont écrits en premières positions
798    ALLOCATE(nqfils(nqtot),nqdesc(nqtot))   
799    ALLOCATE(iqfils(nqtot,nqtot))   
800    ALLOCATE(iqpere(nqtot))
801    nqperes=0
802    nqfils(:)=0
803    nqdesc(:)=0
804    iqfils(:,:)=0
805    iqpere(:)=0
806    nqdesc_tot=0   
807    DO iq=1,nqtot
808      if (tnom_transp(iq) == 'air') then
809        ! ceci est un traceur père
810        WRITE(lunout,*) 'Le traceur',iq,', appele ',trim(tnom_0(iq)),', est un pere'
811        nqperes=nqperes+1
812        iqpere(iq)=0
813      else !if (tnom_transp(iq) == 'air') then
814        ! ceci est un fils. Qui est son père?
815        WRITE(lunout,*) 'Le traceur',iq,', appele ',trim(tnom_0(iq)),', est un fils'
816        continu=.true.
817        ipere=1
818        do while (continu)           
819          if (tnom_transp(iq) == tnom_0(ipere)) then
820            ! Son père est ipere
821            WRITE(lunout,*) 'Le traceur',iq,'appele ', &
822      &          trim(tnom_0(iq)),' est le fils de ',ipere,'appele ',trim(tnom_0(ipere))
823            nqfils(ipere)=nqfils(ipere)+1 
824            iqfils(nqfils(ipere),ipere)=iq
825            iqpere(iq)=ipere         
826            continu=.false.
827          else !if (tnom_transp(iq) == tnom_0(ipere)) then
828            ipere=ipere+1
829            if (ipere.gt.nqtot) then
830                WRITE(lunout,*) 'Le traceur',iq,'appele ', &
831      &          trim(tnom_0(iq)),', est orpelin.'
832                CALL abort_gcm('infotrac_init','Un traceur est orphelin',1)
833            endif !if (ipere.gt.nqtot) then
834          endif !if (tnom_transp(iq) == tnom_0(ipere)) then
835        enddo !do while (continu)
836      endif !if (tnom_transp(iq) == 'air') then
837    enddo !DO iq=1,nqtot
838    WRITE(lunout,*) 'infotrac: nqperes=',nqperes   
839    WRITE(lunout,*) 'nqfils=',nqfils
840    WRITE(lunout,*) 'iqpere=',iqpere
841    WRITE(lunout,*) 'iqfils=',iqfils
842
843! Calculer le nombre de descendants à partir de iqfils et de nbfils
844    DO iq=1,nqtot   
845      generation=0
846      continu=.true.
847      ifils=iq
848      do while (continu)
849        ipere=iqpere(ifils)
850        if (ipere.gt.0) then
851         nqdesc(ipere)=nqdesc(ipere)+1   
852         nqdesc_tot=nqdesc_tot+1     
853         iqfils(nqdesc(ipere),ipere)=iq
854         ifils=ipere
855         generation=generation+1
856        else !if (ipere.gt.0) then
857         continu=.false.
858        endif !if (ipere.gt.0) then
859      enddo !do while (continu)   
860      WRITE(lunout,*) 'Le traceur ',iq,', appele ',trim(tnom_0(iq)),' est un traceur de generation: ',generation
861    enddo !DO iq=1,nqtot
862    WRITE(lunout,*) 'infotrac: nqdesc=',nqdesc
863    WRITE(lunout,*) 'iqfils=',iqfils
864    WRITE(lunout,*) 'nqdesc_tot=',nqdesc_tot
865
866! Interdire autres schémas que 10 pour les traceurs fils, et autres schémas
867! que 10 et 14 si des pères ont des fils
868    do iq=1,nqtot
869      if (iqpere(iq).gt.0) then
870        ! ce traceur a un père qui n'est pas l'air
871        ! Seul le schéma 10 est autorisé
872        if (iadv(iq)/=10) then
873           WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv(iq),' is not implemented for sons'
874          CALL abort_gcm('infotrac_init','Sons should be advected by scheme 10',1)
875        endif
876        ! Le traceur père ne peut être advecté que par schéma 10 ou 14:
877        IF (iadv(iqpere(iq))/=10 .AND. iadv(iqpere(iq))/=14) THEN
878          WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv(iq),' is not implemented for fathers'
879          CALL abort_gcm('infotrac_init','Fathers should be advected by scheme 10 ou 14',1)
880        endif !IF (iadv(iqpere(iq))/=10 .AND. iadv(iqpere(iq))/=14) THEN
881     endif !if (iqpere(iq).gt.0) the
882    enddo !do iq=1,nqtot
883
884
885! detecter quels sont les traceurs isotopiques parmi des traceurs
886    call infotrac_isoinit(tnom_0,nqtrue)
887       
888!-----------------------------------------------------------------------
889! Finalize :
890!
891    DEALLOCATE(tnom_0, hadv, vadv,tnom_transp)
892
893
894  END SUBROUTINE infotrac_init
895
896!-----------------------------------------------------------------------
897
898  SUBROUTINE infotrac_isoinit(tnom_0,nqtrue)
899
900#ifdef CPP_IOIPSL
901  use IOIPSL
902#else
903  ! if not using IOIPSL, we still need to use (a local version of) getin
904  use ioipsl_getincom
905#endif
906  implicit none
907 
908    ! inputs
909    INTEGER nqtrue
910    CHARACTER(len=*) tnom_0(nqtrue)
911   
912    ! locals   
913    CHARACTER(len=3), DIMENSION(niso_possibles) :: tnom_iso
914    INTEGER, ALLOCATABLE,DIMENSION(:,:) :: nb_iso,nb_traciso
915    INTEGER, ALLOCATABLE,DIMENSION(:) :: nb_isoind
916    INTEGER :: ntraceurs_zone_prec,iq,phase,ixt,iiso,izone
917    CHARACTER(len=30) :: tnom_trac
918    INCLUDE "iniprint.h"
919
920    tnom_iso=(/'eau','HDO','O18','O17','HTO'/)
921
922    ALLOCATE(nb_iso(niso_possibles,nqo))
923    ALLOCATE(nb_isoind(nqo))
924    ALLOCATE(nb_traciso(niso_possibles,nqo))
925    ALLOCATE(iso_num(nqtot))
926    ALLOCATE(iso_indnum(nqtot))
927    ALLOCATE(zone_num(nqtot))
928    ALLOCATE(phase_num(nqtot))
929     
930    iso_num(:)=0
931    iso_indnum(:)=0
932    zone_num(:)=0
933    phase_num(:)=0
934    indnum_fn_num(:)=0
935    use_iso(:)=.false. 
936    nb_iso(:,:)=0 
937    nb_isoind(:)=0     
938    nb_traciso(:,:)=0
939    niso=0
940    ntraceurs_zone=0 
941    ntraceurs_zone_prec=0
942    ntraciso=0
943
944    do iq=nqo+1,nqtot
945!       write(lunout,*) 'infotrac 569: iq,tnom_0(iq)=',iq,tnom_0(iq)
946       do phase=1,nqo   
947        do ixt= 1,niso_possibles   
948         tnom_trac=trim(tnom_0(phase))//'_'
949         tnom_trac=trim(tnom_trac)//trim(tnom_iso(ixt))
950!         write(*,*) 'phase,ixt,tnom_trac=',phase,ixt,tnom_trac     
951         IF (tnom_0(iq) == tnom_trac) then
952!          write(lunout,*) 'Ce traceur est un isotope'
953          nb_iso(ixt,phase)=nb_iso(ixt,phase)+1   
954          nb_isoind(phase)=nb_isoind(phase)+1   
955          iso_num(iq)=ixt
956          iso_indnum(iq)=nb_isoind(phase)
957          indnum_fn_num(ixt)=iso_indnum(iq)
958          phase_num(iq)=phase
959!          write(lunout,*) 'iso_num(iq)=',iso_num(iq)
960!          write(lunout,*) 'iso_indnum(iq)=',iso_indnum(iq)
961!          write(lunout,*) 'indnum_fn_num(ixt)=',indnum_fn_num(ixt)
962!          write(lunout,*) 'phase_num(iq)=',phase_num(iq)
963          goto 20
964         else if (iqpere(iq).gt.0) then         
965          if (tnom_0(iqpere(iq)) == tnom_trac) then
966!           write(lunout,*) 'Ce traceur est le fils d''un isotope'
967           ! c'est un traceur d'isotope
968           nb_traciso(ixt,phase)=nb_traciso(ixt,phase)+1
969           iso_num(iq)=ixt
970           iso_indnum(iq)=indnum_fn_num(ixt)
971           zone_num(iq)=nb_traciso(ixt,phase)
972           phase_num(iq)=phase
973!           write(lunout,*) 'iso_num(iq)=',iso_num(iq)
974!           write(lunout,*) 'phase_num(iq)=',phase_num(iq)
975!           write(lunout,*) 'zone_num(iq)=',zone_num(iq)
976           goto 20
977          endif !if (tnom_0(iqpere(iq)) == trim(tnom_0(phase))//trim(tnom_iso(ixt))) then
978         endif !IF (tnom_0(iq) == trim(tnom_0(phase))//trim(tnom_iso(ixt))) then
979        enddo !do ixt= niso_possibles
980       enddo !do phase=1,nqo
981  20   continue
982      enddo !do iq=1,nqtot
983
984!      write(lunout,*) 'iso_num=',iso_num
985!      write(lunout,*) 'iso_indnum=',iso_indnum
986!      write(lunout,*) 'zone_num=',zone_num 
987!      write(lunout,*) 'phase_num=',phase_num
988!      write(lunout,*) 'indnum_fn_num=',indnum_fn_num
989
990      do ixt= 1,niso_possibles 
991       
992       if (nqo.gt.0) then ! Ehouarn: because tests below only valid if nqo>=1,
993
994        if (nb_iso(ixt,1).eq.1) then
995          ! on vérifie que toutes les phases ont le même nombre de
996          ! traceurs
997          do phase=2,nqo
998            if (nb_iso(ixt,phase).ne.nb_iso(ixt,1)) then
999!              write(lunout,*) 'ixt,phase,nb_iso=',ixt,phase,nb_iso(ixt,phase)
1000              CALL abort_gcm('infotrac_init','Phases must have same number of isotopes',1)
1001            endif
1002          enddo !do phase=2,nqo
1003
1004          niso=niso+1
1005          use_iso(ixt)=.true.
1006          ntraceurs_zone=nb_traciso(ixt,1)
1007
1008          ! on vérifie que toutes les phases ont le même nombre de
1009          ! traceurs
1010          do phase=2,nqo
1011            if (nb_traciso(ixt,phase).ne.ntraceurs_zone) then
1012              write(lunout,*) 'ixt,phase,nb_traciso=',ixt,phase,nb_traciso(ixt,phase)
1013              write(lunout,*) 'ntraceurs_zone=',ntraceurs_zone
1014              CALL abort_gcm('infotrac_init','Phases must have same number of tracers',1)
1015            endif 
1016          enddo  !do phase=2,nqo
1017          ! on vérifie que tous les isotopes ont le même nombre de
1018          ! traceurs
1019          if (ntraceurs_zone_prec.gt.0) then               
1020            if (ntraceurs_zone.eq.ntraceurs_zone_prec) then
1021              ntraceurs_zone_prec=ntraceurs_zone
1022            else !if (ntraceurs_zone.eq.ntraceurs_zone_prec) then
1023              write(*,*) 'ntraceurs_zone_prec,ntraceurs_zone=',ntraceurs_zone_prec,ntraceurs_zone   
1024              CALL abort_gcm('infotrac_init', &
1025               &'Isotope tracers are not well defined in traceur.def',1)           
1026            endif !if (ntraceurs_zone.eq.ntraceurs_zone_prec) then
1027           endif !if (ntraceurs_zone_prec.gt.0) then
1028
1029        else if (nb_iso(ixt,1).ne.0) then
1030           WRITE(lunout,*) 'nqo,ixt=',nqo,ixt
1031           WRITE(lunout,*) 'nb_iso(ixt,1)=',nb_iso(ixt,1)   
1032           CALL abort_gcm('infotrac_init','Isotopes are not well defined in traceur.def',1)     
1033        endif   !if (nb_iso(ixt,1).eq.1) then       
1034
1035     endif ! of if (nqo.gt.0)
1036
1037    enddo ! do ixt= niso_possibles
1038
1039    ! dimensions isotopique:
1040    ntraciso=niso*(ntraceurs_zone+1)
1041!    WRITE(lunout,*) 'niso=',niso
1042!    WRITE(lunout,*) 'ntraceurs_zone,ntraciso=',ntraceurs_zone,ntraciso   
1043 
1044    ! flags isotopiques:
1045    if (niso.gt.0) then
1046        ok_isotopes=.true.
1047    else
1048        ok_isotopes=.false.
1049    endif
1050!    WRITE(lunout,*) 'ok_isotopes=',ok_isotopes
1051 
1052    if (ok_isotopes) then
1053        ok_iso_verif=.false.
1054        call getin('ok_iso_verif',ok_iso_verif)
1055        ok_init_iso=.false.
1056        call getin('ok_init_iso',ok_init_iso)
1057        tnat=(/1.0,155.76e-6,2005.2e-6,0.004/100.,0.0/)
1058        alpha_ideal=(/1.0,1.01,1.006,1.003,1.0/)
1059    endif !if (ok_isotopes) then 
1060!    WRITE(lunout,*) 'ok_iso_verif=',ok_iso_verif
1061!    WRITE(lunout,*) 'ok_init_iso=',ok_init_iso
1062
1063    if (ntraceurs_zone.gt.0) then
1064        ok_isotrac=.true.
1065    else
1066        ok_isotrac=.false.
1067    endif   
1068!    WRITE(lunout,*) 'ok_isotrac=',ok_isotrac
1069
1070    ! remplissage du tableau iqiso(ntraciso,phase)
1071    ALLOCATE(iqiso(ntraciso,nqo))   
1072    iqiso(:,:)=0     
1073    do iq=1,nqtot
1074        if (iso_num(iq).gt.0) then
1075          ixt=iso_indnum(iq)+zone_num(iq)*niso
1076          iqiso(ixt,phase_num(iq))=iq
1077        endif
1078    enddo
1079!    WRITE(lunout,*) 'iqiso=',iqiso
1080
1081    ! replissage du tableau index_trac(ntraceurs_zone,niso)
1082    ALLOCATE(index_trac(ntraceurs_zone,niso)) 
1083    if (ok_isotrac) then
1084        do iiso=1,niso
1085          do izone=1,ntraceurs_zone
1086             index_trac(izone,iiso)=iiso+izone*niso
1087          enddo
1088        enddo
1089    else !if (ok_isotrac) then     
1090        index_trac(:,:)=0.0
1091    endif !if (ok_isotrac) then
1092!    write(lunout,*) 'index_trac=',index_trac   
1093
1094! Finalize :
1095    DEALLOCATE(nb_iso)
1096
1097  END SUBROUTINE infotrac_isoinit
1098
1099!-----------------------------------------------------------------------
1100
1101! Ehouarn: routine iniadvtrac => from Mars/generic; does essentially the
1102!          same job as infotrac_init. To clean up and merge at some point...
1103      subroutine iniadvtrac(nq,numvanle)
1104!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1105! routine which initializes tracer names and advection schemes
1106! reads these infos from file 'traceur.def' but uses default values
1107! if that file is not found.
1108! Ehouarn Millour. Oct. 2008  (made this LMDZ4-like) for future compatibility
1109!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1110      IMPLICIT NONE
1111
1112!#include "dimensions.h"
1113!#include "advtrac.h"
1114!#include "control.h"
1115
1116! routine arguments:
1117      INTEGER,INTENT(out) :: nq ! number of tracers
1118      INTEGER,INTENT(out) :: numvanle
1119
1120! local variables:
1121      LOGICAL :: first
1122      INTEGER :: iq
1123      INTEGER :: ierr
1124      CHARACTER(len=3) :: qname
1125
1126! Look for file traceur.def
1127      OPEN(90,file='traceur.def',form='formatted',status='old', &
1128              iostat=ierr)
1129      IF (ierr.eq.0) THEN
1130        write(*,*) "iniadvtrac: Reading file traceur.def"
1131        ! read number of tracers:
1132        read(90,*,iostat=ierr) nq
1133        if (ierr.ne.0) then
1134          write(*,*) "iniadvtrac: error reading number of tracers"
1135          write(*,*) "   (first line of traceur.def) "
1136          stop
1137        endif
1138       
1139        ! allocate arrays:
1140        allocate(iadv(nq))
1141        allocate(tname(nq))
1142       
1143        ! initialize advection schemes to Van-Leer for all tracers
1144        do iq=1,nq
1145          iadv(iq)=3 ! Van-Leer
1146        enddo
1147       
1148        do iq=1,nq
1149        ! minimal version, just read in the tracer names, 1 per line
1150          read(90,*,iostat=ierr) tname(iq)
1151          if (ierr.ne.0) then
1152            write(*,*) 'iniadvtrac: error reading tracer names...'
1153            stop
1154          endif
1155        enddo !of do iq=1,nq
1156        close(90) ! done reading tracer names, close file
1157      ENDIF ! of IF (ierr.eq.0)
1158
1159!  ....  Choix  des shemas d'advection pour l'eau et les traceurs  ...
1160!  ...................................................................
1161!
1162!     iadv = 1    shema  transport type "humidite specifique LMD" 
1163!     iadv = 2    shema   amont
1164!     iadv = 3    shema  Van-leer
1165!     iadv = 4    schema  Van-leer + humidite specifique
1166!                        Modif F.Codron
1167!
1168!
1169      DO  iq = 1, nq-1
1170       IF( iadv(iq).EQ.1 ) PRINT *,' Choix du shema humidite specifique'&
1171       ,' pour le traceur no ', iq
1172       IF( iadv(iq).EQ.2 ) PRINT *,' Choix du shema  amont',' pour le'  &
1173       ,' traceur no ', iq
1174       IF( iadv(iq).EQ.3 ) PRINT *,' Choix du shema  Van-Leer ',' pour' &
1175       ,'le traceur no ', iq
1176
1177       IF( iadv(iq).EQ.4 )  THEN
1178         PRINT *,' Le shema  Van-Leer + humidite specifique ',          &
1179       ' est  uniquement pour la vapeur d eau .'
1180         PRINT *,' Corriger iadv( ',iq, ')  et repasser ! '
1181         CALL ABORT
1182       ENDIF
1183
1184       IF( iadv(iq).LE.0.OR.iadv(iq).GT.4 )   THEN
1185        PRINT *,' Erreur dans le choix de iadv (nqtot).Corriger et '    &
1186       ,' repasser car  iadv(iq) = ', iadv(iq)
1187         CALL ABORT
1188       ENDIF
1189      ENDDO
1190
1191       IF( iadv(nq).EQ.1 ) PRINT *,' Choix du shema humidite '          &
1192       ,'specifique pour la vapeur d''eau'
1193       IF( iadv(nq).EQ.2 ) PRINT *,' Choix du shema  amont',' pour la'  &
1194       ,' vapeur d''eau '
1195       IF( iadv(nq).EQ.3 ) PRINT *,' Choix du shema  Van-Leer '         &
1196       ,' pour la vapeur d''eau'
1197       IF( iadv(nq).EQ.4 ) PRINT *,' Choix du shema  Van-Leer + '       &
1198       ,' humidite specifique pour la vapeur d''eau'
1199!
1200       IF( (iadv(nq).LE.0).OR.(iadv(nq).GT.4) )   THEN
1201        PRINT *,' Erreur dans le choix de iadv (nqtot).Corriger et '    &
1202       ,' repasser car  iadv(nqtot) = ', iadv(nqtot)
1203         CALL ABORT
1204       ENDIF
1205
1206      first = .TRUE.
1207      numvanle = nq + 1
1208      DO  iq = 1, nq
1209        IF(((iadv(iq).EQ.3).OR.(iadv(iq).EQ.4)).AND.first ) THEN
1210          numvanle = iq
1211          first    = .FALSE.
1212        ENDIF
1213      ENDDO
1214!
1215      DO  iq = 1, nq
1216
1217      IF( (iadv(iq).NE.3.AND.iadv(iq).NE.4).AND.iq.GT.numvanle )  THEN
1218          PRINT *,' Il y a discontinuite dans le choix du shema de ',   &
1219          'Van-leer pour les traceurs . Corriger et repasser . '
1220           CALL ABORT
1221      ENDIF
1222
1223      ENDDO
1224!
1225      end subroutine iniadvtrac
1226
1227
1228END MODULE infotrac
Note: See TracBrowser for help on using the repository browser.