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

Last change on this file since 1575 was 1575, checked in by emillour, 9 years ago

All GCMs:
Further adaptations to keep up with changes in LMDZ5 concerning
physics/dynamics separation (up to rev r2575 of LMDZ5)

  • dyn3d_common:
  • infotrac.F90 : propagate initialisations for INCA (Earth GCM)
  • misc:
  • wxios.F90: updates to use the XIOS2 library
  • dynphy_lonlat:
  • grid_atob_m.F90: fix for some zoomed grid interpolation cases

EM

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