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

Last change on this file since 2334 was 2307, checked in by mvals, 6 years ago

Mars GCM:
follow-up of the commit regarding the dynamical transport of isotopes: new variables for the thresholds for zq(pere) and masseq in the transport of
the isotopic Ratio
MV

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