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

Last change on this file since 1403 was 1403, checked in by emillour, 10 years ago

All models: Reorganizing the physics/dynamics interface.

  • makelmdz and makelmdz_fcm scripts adapted to handle the new directory settings
  • misc: (replaces what was the "bibio" directory)
  • Should only contain extremely generic (and non physics or dynamics-specific) routines
  • Therefore moved initdynav.F90, initfluxsto.F, inithist.F, writedynav.F90, write_field.F90, writehist.F to "dyn3d_common"
  • dynlonlat_phylonlat: (new interface directory)
  • This directory contains routines relevent to physics/dynamics grid interactions, e.g. routines gr_dyn_fi or gr_fi_dyn and calfis
  • Moreover the dynlonlat_phylonlat contains directories "phy*" corresponding to each physics package "phy*" to be used. These subdirectories should only contain specific interfaces (e.g. iniphysiq) or main programs (e.g. newstart)
  • phy*/dyn1d: this subdirectory contains the 1D model using physics from phy*

EM

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