source: trunk/LMDZ.GENERIC/libf/dyn3d/infotrac.F90 @ 2035

Last change on this file since 2035 was 2035, checked in by jvatant, 6 years ago

Fix a bug in 1D, again some (late) follow-up
of r1971, r1978 and r1980 about the size of tracer names.
--JVO

File size: 20.1 KB
Line 
1MODULE infotrac
2
3IMPLICIT NONE
4! nqtot : total number of tracers and higher order of moment, water vapor and liquid included
5  INTEGER, SAVE :: nqtot
6! CR: add number of tracers for water (for Earth model only!!)
7  INTEGER, SAVE :: nqo
8
9! nbtr : number of tracers not including higher order of moment or water vapor or liquid
10!        number of tracers used in the physics
11  INTEGER, SAVE :: nbtr
12
13! Name variables
14  CHARACTER(len=30), ALLOCATABLE, DIMENSION(:), SAVE :: tname ! tracer short name for restart and diagnostics
15  CHARACTER(len=33), ALLOCATABLE, DIMENSION(:), SAVE :: ttext ! tracer long name for diagnostics
16
17! iadv  : index of trasport schema for each tracer
18  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE    :: iadv
19
20! niadv : vector keeping the coorspondance between all tracers(nqtot) treated in the
21!         dynamic part of the code and the tracers (nbtr+2) used in the physics part of the code.
22  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE    :: niadv ! equivalent dyn / physique
23
24! conv_flg(it)=0 : convection desactivated for tracer number it
25  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE  :: conv_flg
26! pbl_flg(it)=0  : boundary layer diffusion desactivaded for tracer number it
27  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE  :: pbl_flg
28
29  CHARACTER(len=4),SAVE :: type_trac
30  CHARACTER(len=8),DIMENSION(:),ALLOCATABLE, SAVE :: solsym
31 
32CONTAINS
33
34  SUBROUTINE infotrac_init
35    USE control_mod
36#ifdef REPROBUS
37    USE CHEM_REP, ONLY : Init_chem_rep_trac
38#endif
39    IMPLICIT NONE
40!=======================================================================
41!
42!   Auteur:  P. Le Van /L. Fairhead/F.Hourdin
43!   -------
44!   Modif special traceur F.Forget 05/94
45!   Modif M-A Filiberti 02/02 lecture de traceur.def
46!
47!   Objet:
48!   ------
49!   GCM LMD nouvelle grille
50!
51!=======================================================================
52!   ... modification de l'integration de q ( 26/04/94 ) ....
53!-----------------------------------------------------------------------
54! Declarations
55
56    INCLUDE "dimensions.h"
57    INCLUDE "iniprint.h"
58
59! Local variables
60    INTEGER, ALLOCATABLE, DIMENSION(:) :: hadv  ! index of horizontal trasport schema
61    INTEGER, ALLOCATABLE, DIMENSION(:) :: vadv  ! index of vertical trasport schema
62
63    CHARACTER(len=30), ALLOCATABLE, DIMENSION(:) :: tnom_0  ! tracer short name
64    CHARACTER(len=3), DIMENSION(30) :: descrq
65    CHARACTER(len=1), DIMENSION(3)  :: txts
66    CHARACTER(len=2), DIMENSION(9)  :: txtp
67    CHARACTER(len=30)               :: str1,str2
68 
69    INTEGER :: nqtrue  ! number of tracers read from tracer.def, without higer order of moment
70    INTEGER :: iq, new_iq, iiq, jq, ierr, ierr2, ierr3
71   
72    character(len=80) :: line ! to store a line of text
73 
74    character(len=*),parameter :: modname="infotrac_init"
75!-----------------------------------------------------------------------
76! Initialization :
77!
78    txts=(/'x','y','z'/)
79    txtp=(/'x ','y ','z ','xx','xy','xz','yy','yz','zz'/)
80
81    descrq(14)='VLH'
82    descrq(10)='VL1'
83    descrq(11)='VLP'
84    descrq(12)='FH1'
85    descrq(13)='FH2'
86    descrq(16)='PPM'
87    descrq(17)='PPS'
88    descrq(18)='PPP'
89    descrq(20)='SLP'
90    descrq(30)='PRA'
91   
92    IF (planet_type=='earth') THEN
93     ! Coherence test between parameter type_trac, config_inca and preprocessing keys
94     IF (type_trac=='inca') THEN
95       WRITE(lunout,*) 'You have choosen to couple with INCA chemestry model : type_trac=', &
96            type_trac,' config_inca=',config_inca
97       IF (config_inca/='aero' .AND. config_inca/='aeNP' .AND. config_inca/='chem') THEN
98          WRITE(lunout,*) 'Incoherence between type_trac and config_inca. Model stops. Modify run.def'
99          CALL abort_gcm('infotrac_init','Incoherence between type_trac and config_inca',1)
100       END IF
101#ifndef INCA
102       WRITE(lunout,*) 'To run this option you must add cpp key INCA and compile with INCA code'
103       CALL abort_gcm('infotrac_init','You must compile with cpp key INCA',1)
104#endif
105     ELSE IF (type_trac=='repr') THEN
106       WRITE(lunout,*) 'You have choosen to couple with REPROBUS chemestry model : type_trac=', type_trac
107#ifndef REPROBUS
108       WRITE(lunout,*) 'To run this option you must add cpp key REPROBUS and compile with REPRPBUS code'
109       CALL abort_gcm('infotrac_init','You must compile with cpp key REPROBUS',1)
110#endif
111     ELSE IF (type_trac == 'lmdz') THEN
112       WRITE(lunout,*) 'Tracers are treated in LMDZ only : type_trac=', type_trac
113     ELSE
114       WRITE(lunout,*) 'type_trac=',type_trac,' not possible. Model stops'
115       CALL abort_gcm('infotrac_init','bad parameter',1)
116     END IF
117
118     ! Test if config_inca is other then none for run without INCA
119     IF (type_trac/='inca' .AND. config_inca/='none') THEN
120       WRITE(lunout,*) 'config_inca will now be changed to none as you do not couple with INCA model'
121       config_inca='none'
122     END IF
123    ELSE
124     type_trac='plnt'  ! planets... May want to dissociate between each later.
125    ENDIF ! of IF (planet_type=='earth')
126
127!-----------------------------------------------------------------------
128!
129! 1) Get the true number of tracers + water vapor/liquid
130!    Here true tracers (nqtrue) means declared tracers (only first order)
131!
132!-----------------------------------------------------------------------
133    IF (planet_type=='earth') THEN
134     IF (type_trac == 'lmdz' .OR. type_trac == 'repr') THEN
135       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
136       IF(ierr.EQ.0) THEN
137          WRITE(lunout,*) trim(modname),': Open traceur.def : ok'
138          READ(90,*) nqtrue
139       ELSE
140          WRITE(lunout,*) trim(modname),': Problem in opening traceur.def'
141          WRITE(lunout,*) trim(modname),': WARNING using defaut values'
142          nqtrue=4 ! Defaut value
143       END IF
144       ! For Earth, water vapour & liquid tracers are not in the physics
145       nbtr=nqtrue-2
146     ELSE ! type_trac=inca
147       ! nbtr has been read from INCA by init_cont_lmdz() in gcm.F
148       nqtrue=nbtr+2
149     END IF
150
151     IF (nqtrue < 2) THEN
152       WRITE(lunout,*) trim(modname),': nqtrue=',nqtrue, ' is not allowded. 2 tracers is the minimum'
153       CALL abort_gcm('infotrac_init','Not enough tracers',1)
154     END IF
155
156! Transfert number of tracers to Reprobus
157     IF (type_trac == 'repr') THEN
158#ifdef REPROBUS
159       CALL Init_chem_rep_trac(nbtr)
160#endif
161     END IF
162
163    ELSE  ! not Earth
164       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
165       IF(ierr.EQ.0) THEN
166          WRITE(lunout,*) 'Open traceur.def : ok'
167          READ(90,*) nqtrue
168       ELSE
169          WRITE(lunout,*) 'Problem in opening traceur.def'
170          WRITE(lunout,*) 'ATTENTION using defaut values: nqtrue=1'
171          nqtrue=1 ! Defaut value
172       END IF
173       ! Other planets (for now); we have the same number of tracers
174       ! in the dynamics than in the physics
175       nbtr=nqtrue
176     
177    ENDIF  ! planet_type
178!
179! Allocate variables depending on nqtrue and nbtr
180!
181    ALLOCATE(tnom_0(nqtrue), hadv(nqtrue), vadv(nqtrue))
182    ALLOCATE(conv_flg(nbtr), pbl_flg(nbtr), solsym(nbtr))
183    conv_flg(:) = 1 ! convection activated for all tracers
184    pbl_flg(:)  = 1 ! boundary layer activated for all tracers
185
186!-----------------------------------------------------------------------
187! 2)     Choix  des schemas d'advection pour l'eau et les traceurs
188!
189!     iadv = 1    schema  transport type "humidite specifique LMD"
190!     iadv = 2    schema   amont
191!     iadv = 14   schema  Van-leer + humidite specifique
192!                            Modif F.Codron
193!     iadv = 10   schema  Van-leer (retenu pour l'eau vapeur et liquide)
194!     iadv = 11   schema  Van-Leer pour hadv et version PPM (Monotone) pour vadv
195!     iadv = 12   schema  Frederic Hourdin I
196!     iadv = 13   schema  Frederic Hourdin II
197!     iadv = 16   schema  PPM Monotone(Collela & Woodward 1984)
198!     iadv = 17   schema  PPM Semi Monotone (overshoots autorisés)
199!     iadv = 18   schema  PPM Positif Defini (overshoots undershoots autorisés)
200!     iadv = 20   schema  Slopes
201!     iadv = 30   schema  Prather
202!
203!        Dans le tableau q(ij,l,iq) : iq = 1  pour l'eau vapeur
204!                                     iq = 2  pour l'eau liquide
205!       Et eventuellement             iq = 3,nqtot pour les autres traceurs
206!
207!        iadv(1): choix pour l'eau vap. et  iadv(2) : choix pour l'eau liq.
208!------------------------------------------------------------------------
209!
210!    Get choice of advection schema from file tracer.def or from INCA
211!---------------------------------------------------------------------
212    IF (planet_type=='earth') THEN
213     IF (type_trac == 'lmdz' .OR. type_trac == 'repr') THEN
214       IF(ierr.EQ.0) THEN
215          ! Continue to read tracer.def
216          DO iq=1,nqtrue
217             READ(90,*) hadv(iq),vadv(iq),tnom_0(iq)
218          END DO
219          CLOSE(90) 
220       ELSE ! Without tracer.def, set default values (for Earth!)
221         if ((nqtrue==4).and.(planet_type=="earth")) then
222          hadv(1) = 14
223          vadv(1) = 14
224          tnom_0(1) = 'H2Ov'
225          hadv(2) = 10
226          vadv(2) = 10
227          tnom_0(2) = 'H2Ol'
228          hadv(3) = 10
229          vadv(3) = 10
230          tnom_0(3) = 'RN'
231          hadv(4) = 10
232          vadv(4) = 10
233          tnom_0(4) = 'PB'
234         else
235           ! Error message, we need a traceur.def file
236           write(lunout,*) trim(modname),&
237           ': Cannot set default tracer names!'
238           write(lunout,*) trim(modname),' Make a traceur.def file!!!'
239           CALL abort_gcm('infotrac_init','Need a traceur.def file!',1)
240         endif ! of if (nqtrue==4)
241       END IF
242       
243!CR: nombre de traceurs de l eau
244       if (tnom_0(3) == 'H2Oi') then
245          nqo=3
246       else
247          nqo=2
248       endif
249
250       WRITE(lunout,*) trim(modname),': Valeur de traceur.def :'
251       WRITE(lunout,*) trim(modname),': nombre de traceurs ',nqtrue
252       DO iq=1,nqtrue
253          WRITE(lunout,*) hadv(iq),vadv(iq),tnom_0(iq)
254       END DO
255
256     ELSE  ! type_trac=inca : config_inca='aero' ou 'chem'
257! le module de chimie fournit les noms des traceurs
258! et les schemas d'advection associes.
259     
260#ifdef INCA
261       CALL init_transport( &
262            hadv, &
263            vadv, &
264            conv_flg, &
265            pbl_flg,  &
266            tracnam)
267#endif
268       tnom_0(1)='H2Ov'
269       tnom_0(2)='H2Ol'
270
271       DO iq =3,nqtrue
272          tnom_0(iq)=solsym(iq-2)
273       END DO
274       nqo = 2
275
276     END IF ! type_trac
277
278    ELSE  ! not Earth
279       IF(ierr.EQ.0) THEN
280          ! Continue to read tracer.def
281          DO iq=1,nqtrue
282             !READ(90,*) hadv(iq),vadv(iq),tnom_0(iq)
283            ! try to be smart when reading traceur.def
284            read(90,'(80a)') line ! store the line from traceur.def
285            ! assume format is hadv,vadv,tnom_0
286            read(line,*,iostat=ierr2) hadv(iq),vadv(iq),tnom_0(iq)
287            if (ierr2.ne.0) then
288              ! maybe format is tnom0,hadv,vadv
289              read(line,*,iostat=ierr3) tnom_0(iq),hadv(iq),vadv(iq)
290              if (ierr3.ne.0) then
291                ! assume only tnom0 is provided (havd and vad default to 10)
292                read(line,*) tnom_0(iq)
293                hadv(iq)=10
294                vadv(iq)=10
295              endif
296            endif ! of if(ierr2.ne.0)
297          END DO ! of DO iq=1,nqtrue
298          CLOSE(90) 
299       ELSE ! Without tracer.def
300          hadv(1) = 10
301          vadv(1) = 10
302          tnom_0(1) = 'dummy'
303       END IF
304       
305       WRITE(lunout,*) trim(modname),': Valeur de traceur.def :'
306       WRITE(lunout,*) trim(modname),': nombre de traceurs ',nqtrue
307       DO iq=1,nqtrue
308          WRITE(lunout,*) hadv(iq),vadv(iq),tnom_0(iq)
309       END DO
310
311    ENDIF  ! planet_type
312
313!-----------------------------------------------------------------------
314!
315! 3) Verify if advection schema 20 or 30 choosen
316!    Calculate total number of tracers needed: nqtot
317!    Allocate variables depending on total number of tracers
318!-----------------------------------------------------------------------
319    new_iq=0
320    DO iq=1,nqtrue
321       ! Add tracers for certain advection schema
322       IF (hadv(iq)<20 .AND. vadv(iq)<20 ) THEN
323          new_iq=new_iq+1  ! no tracers added
324       ELSE IF (hadv(iq)==20 .AND. vadv(iq)==20 ) THEN
325          new_iq=new_iq+4  ! 3 tracers added
326       ELSE IF (hadv(iq)==30 .AND. vadv(iq)==30 ) THEN
327          new_iq=new_iq+10 ! 9 tracers added
328       ELSE
329          WRITE(lunout,*) trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq)
330          CALL abort_gcm('infotrac_init','Bad choice of advection schema - 1',1)
331       END IF
332    END DO
333   
334    IF (new_iq /= nqtrue) THEN
335       ! The choice of advection schema imposes more tracers
336       ! Assigne total number of tracers
337       nqtot = new_iq
338
339       WRITE(lunout,*) trim(modname),': The choice of advection schema for one or more tracers'
340       WRITE(lunout,*) 'makes it necessary to add tracers'
341       WRITE(lunout,*) trim(modname)//': ',nqtrue,' is the number of true tracers'
342       WRITE(lunout,*) trim(modname)//': ',nqtot, ' is the total number of tracers needed'
343
344    ELSE
345       ! The true number of tracers is also the total number
346       nqtot = nqtrue
347    END IF
348
349!
350! Allocate variables with total number of tracers, nqtot
351!
352    ALLOCATE(tname(nqtot), ttext(nqtot))
353    ALLOCATE(iadv(nqtot), niadv(nqtot))
354
355!-----------------------------------------------------------------------
356!
357! 4) Determine iadv, long and short name
358!
359!-----------------------------------------------------------------------
360    new_iq=0
361    DO iq=1,nqtrue
362       new_iq=new_iq+1
363
364       ! Verify choice of advection schema
365       IF (hadv(iq)==vadv(iq)) THEN
366          iadv(new_iq)=hadv(iq)
367       ELSE IF (hadv(iq)==10 .AND. vadv(iq)==16) THEN
368          iadv(new_iq)=11
369       ELSE
370          WRITE(lunout,*)trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq)
371
372          CALL abort_gcm('infotrac_init','Bad choice of advection schema - 2',1)
373       END IF
374     
375       str1=tnom_0(iq)
376       tname(new_iq)= tnom_0(iq)
377       IF (iadv(new_iq)==0) THEN
378          ttext(new_iq)=trim(str1)
379       ELSE
380          ttext(new_iq)=trim(tnom_0(iq))//descrq(iadv(new_iq))
381       END IF
382
383       ! schemas tenant compte des moments d'ordre superieur
384       str2=ttext(new_iq)
385       IF (iadv(new_iq)==20) THEN
386          DO jq=1,3
387             new_iq=new_iq+1
388             iadv(new_iq)=-20
389             ttext(new_iq)=trim(str2)//txts(jq)
390             tname(new_iq)=trim(str1)//txts(jq)
391          END DO
392       ELSE IF (iadv(new_iq)==30) THEN
393          DO jq=1,9
394             new_iq=new_iq+1
395             iadv(new_iq)=-30
396             ttext(new_iq)=trim(str2)//txtp(jq)
397             tname(new_iq)=trim(str1)//txtp(jq)
398          END DO
399       END IF
400    END DO
401
402!
403! Find vector keeping the correspodence between true and total tracers
404!
405    niadv(:)=0
406    iiq=0
407    DO iq=1,nqtot
408       IF(iadv(iq).GE.0) THEN
409          ! True tracer
410          iiq=iiq+1
411          niadv(iiq)=iq
412       ENDIF
413    END DO
414
415
416    WRITE(lunout,*) trim(modname),': Information stored in infotrac :'
417    WRITE(lunout,*) trim(modname),': iadv  niadv tname  ttext :'
418    DO iq=1,nqtot
419       WRITE(lunout,*) iadv(iq),niadv(iq),&
420       ' ',trim(tname(iq)),' ',trim(ttext(iq))
421    END DO
422
423!
424! Test for advection schema.
425! This version of LMDZ only garantees iadv=10 and iadv=14 (14 only for water vapour) .
426!
427    DO iq=1,nqtot
428       IF (iadv(iq)/=10 .AND. iadv(iq)/=14 .AND. iadv(iq)/=0) THEN
429          WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv(iq),' is not tested in this version of LMDZ'
430          CALL abort_gcm('infotrac_init','In this version only iadv=10 and iadv=14 is tested!',1)
431       ELSE IF (iadv(iq)==14 .AND. iq/=1) THEN
432          WRITE(lunout,*)trim(modname),'STOP : The option iadv=',iadv(iq),' is not tested in this version of LMDZ'
433          CALL abort_gcm('infotrac_init','In this version iadv=14 is only permitted for water vapour!',1)
434       END IF
435    END DO
436
437!-----------------------------------------------------------------------
438! Finalize :
439!
440    DEALLOCATE(tnom_0, hadv, vadv)
441
442
443  END SUBROUTINE infotrac_init
444
445! Ehouarn: routine iniadvtrac => from Mars/generic; does essentially the
446!          same job as infotrac_init. To clean up and merge at some point...
447      subroutine iniadvtrac(nq,numvanle)
448!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
449! routine which initializes tracer names and advection schemes
450! reads these infos from file 'traceur.def' but uses default values
451! if that file is not found.
452! Ehouarn Millour. Oct. 2008  (made this LMDZ4-like) for future compatibility
453!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
454      IMPLICIT NONE
455
456!#include "dimensions.h"
457!#include "advtrac.h"
458!#include "control.h"
459
460! routine arguments:
461      INTEGER,INTENT(out) :: nq ! number of tracers
462      INTEGER,INTENT(out) :: numvanle
463
464! local variables:
465      LOGICAL :: first
466      INTEGER :: iq
467      INTEGER :: ierr
468      CHARACTER(len=3) :: qname
469
470! Look for file traceur.def
471      OPEN(90,file='traceur.def',form='formatted',status='old', &
472              iostat=ierr)
473      IF (ierr.eq.0) THEN
474        write(*,*) "iniadvtrac: Reading file traceur.def"
475        ! read number of tracers:
476        read(90,*,iostat=ierr) nq
477        if (ierr.ne.0) then
478          write(*,*) "iniadvtrac: error reading number of tracers"
479          write(*,*) "   (first line of traceur.def) "
480          stop
481        endif
482       
483        ! allocate arrays:
484        allocate(iadv(nq))
485        allocate(tname(nq))
486       
487        ! initialize advection schemes to Van-Leer for all tracers
488        do iq=1,nq
489          iadv(iq)=3 ! Van-Leer
490        enddo
491       
492        do iq=1,nq
493        ! minimal version, just read in the tracer names, 1 per line
494          read(90,*,iostat=ierr) tname(iq)
495          if (ierr.ne.0) then
496            write(*,*) 'iniadvtrac: error reading tracer names...'
497            stop
498          endif
499        enddo !of do iq=1,nq
500        close(90) ! done reading tracer names, close file
501      ENDIF ! of IF (ierr.eq.0)
502
503!  ....  Choix  des shemas d'advection pour l'eau et les traceurs  ...
504!  ...................................................................
505!
506!     iadv = 1    shema  transport type "humidite specifique LMD" 
507!     iadv = 2    shema   amont
508!     iadv = 3    shema  Van-leer
509!     iadv = 4    schema  Van-leer + humidite specifique
510!                        Modif F.Codron
511!
512!
513      DO  iq = 1, nq-1
514       IF( iadv(iq).EQ.1 ) PRINT *,' Choix du shema humidite specifique'&
515       ,' pour le traceur no ', iq
516       IF( iadv(iq).EQ.2 ) PRINT *,' Choix du shema  amont',' pour le'  &
517       ,' traceur no ', iq
518       IF( iadv(iq).EQ.3 ) PRINT *,' Choix du shema  Van-Leer ',' pour' &
519       ,'le traceur no ', iq
520
521       IF( iadv(iq).EQ.4 )  THEN
522         PRINT *,' Le shema  Van-Leer + humidite specifique ',          &
523       ' est  uniquement pour la vapeur d eau .'
524         PRINT *,' Corriger iadv( ',iq, ')  et repasser ! '
525         CALL ABORT
526       ENDIF
527
528       IF( iadv(iq).LE.0.OR.iadv(iq).GT.4 )   THEN
529        PRINT *,' Erreur dans le choix de iadv (nqtot).Corriger et '    &
530       ,' repasser car  iadv(iq) = ', iadv(iq)
531         CALL ABORT
532       ENDIF
533      ENDDO
534
535!       IF( iadv(nq).EQ.1 ) PRINT *,' Choix du shema humidite '          &
536!       ,'specifique pour la vapeur d''eau'
537!       IF( iadv(nq).EQ.2 ) PRINT *,' Choix du shema  amont',' pour la'  &
538!       ,' vapeur d''eau '
539!       IF( iadv(nq).EQ.3 ) PRINT *,' Choix du shema  Van-Leer '         &
540!       ,' pour la vapeur d''eau'
541!       IF( iadv(nq).EQ.4 ) PRINT *,' Choix du shema  Van-Leer + '       &
542!       ,' humidite specifique pour la vapeur d''eau'
543!
544!       IF( (iadv(nq).LE.0).OR.(iadv(nq).GT.4) )   THEN
545!        PRINT *,' Erreur dans le choix de iadv (nqtot).Corriger et '    &
546!       ,' repasser car  iadv(nqtot) = ', iadv(nqtot)
547!         CALL ABORT
548!       ENDIF
549
550      first = .TRUE.
551      numvanle = nq + 1
552      DO  iq = 1, nq
553        IF(((iadv(iq).EQ.3).OR.(iadv(iq).EQ.4)).AND.first ) THEN
554          numvanle = iq
555          first    = .FALSE.
556        ENDIF
557      ENDDO
558!
559      DO  iq = 1, nq
560
561      IF( (iadv(iq).NE.3.AND.iadv(iq).NE.4).AND.iq.GT.numvanle )  THEN
562          PRINT *,' Il y a discontinuite dans le choix du shema de ',   &
563          'Van-leer pour les traceurs . Corriger et repasser . '
564           CALL ABORT
565      ENDIF
566
567      ENDDO
568!
569      end subroutine iniadvtrac
570
571
572END MODULE infotrac
Note: See TracBrowser for help on using the repository browser.