source: LMDZ6/trunk/libf/dyn3d_common/infotrac.F90 @ 4050

Last change on this file since 4050 was 4050, checked in by dcugnet, 2 years ago

Second commit for new tracers.

  • include most of the keys in the tracers descriptor vector "tracers(:)".
  • fix in phylmdiso/cv3_routines: fq_* variables were used where their fxt_* counterparts were expected.
  • multiple IF(nqdesc(iq)>0) and IF(nqfils(iq)>0) tests suppressed, because they are not needed: "do ... enddo" loops with 0 upper bound are not executed.
  • remove French accents from comments (encoding problem) in phylmdiso/cv3_routines and phylmdiso/cv30_routines.
  • modifications in "isotopes_verif_mod", where the call to function "iso_verif_tag17_q_deltad_chn" in "iso_verif_tag17_q_deltad_chn" was not detected at linking stage, although defined in the same module (?).
  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:keywords set to Id
File size: 35.8 KB
Line 
1! $Id: infotrac.F90 4050 2021-12-23 17:54:17Z dcugnet $
2!
3MODULE infotrac
4
5  USE readTracFiles_mod, ONLY: trac_type, isot_type, maxlen, ancestor, delPhase
6
7! nqtot : total number of tracers and higher order of moment, water vapor and liquid included
8  INTEGER, SAVE :: nqtot
9!CR: on ajoute le nombre de traceurs de l eau
10  INTEGER, SAVE :: nqo
11
12! nbtr : number of tracers not including higher order of moment or water vapor or liquid
13!        number of tracers used in the physics
14  INTEGER, SAVE :: nbtr
15
16! CRisi: on retranche les isotopes des traceurs habituels
17! On fait un tableaux d'indices des traceurs qui passeront dans phytrac
18  INTEGER, SAVE :: nqtottr
19  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: itr_indice
20
21! CRisi: nb traceurs peres= directement advectes par l'air
22  INTEGER, SAVE :: nqperes
23
24! ThL: nb traceurs INCA
25  INTEGER, SAVE :: nqINCA
26
27! ThL: nb traceurs CO2
28  INTEGER, SAVE :: nqCO2
29
30! DC: derived types containing informations about tracers and isotopes
31  TYPE(trac_type), TARGET,  SAVE, ALLOCATABLE ::  tracers(:)    !=== TRACERS DESCRIPTORS VECTOR
32  TYPE(isot_type), TARGET,  SAVE, ALLOCATABLE :: isotopes(:)    !=== ISOTOPES PARAMETERS VECTOR
33
34! niadv : vector keeping the coorspondance between all tracers(nqtot) treated in the
35!         dynamic part of the code and the tracers (nbtr+2) used in the physics part of the code.
36  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE    :: niadv ! equivalent dyn / physique
37
38  REAL :: qperemin,masseqmin,ratiomin ! MVals et CRisi
39  PARAMETER (qperemin=1e-30,masseqmin=1e-18,ratiomin=1e-16) ! MVals
40
41! conv_flg(it)=0 : convection desactivated for tracer number it
42  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE  :: conv_flg
43! pbl_flg(it)=0  : boundary layer diffusion desactivaded for tracer number it
44  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE  :: pbl_flg
45
46  CHARACTER(len=4),SAVE :: type_trac
47  CHARACTER(len=8),DIMENSION(:),ALLOCATABLE, SAVE :: solsym
48
49! CRisi: cas particulier des isotopes
50  LOGICAL,SAVE :: ok_isotopes,ok_iso_verif,ok_isotrac,ok_init_iso
51  INTEGER :: niso_possibles   
52  PARAMETER ( niso_possibles=5)
53  REAL, DIMENSION (niso_possibles),SAVE :: tnat,alpha_ideal
54  LOGICAL, DIMENSION(niso_possibles),SAVE ::  use_iso
55  INTEGER, ALLOCATABLE, DIMENSION(:,:), SAVE ::  iqiso ! donne indice iq en fn de (ixt,phase)
56  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE ::  iso_indnum ! donne numero iso entre 1 et niso effectif en fn de nqtot
57  INTEGER, DIMENSION(niso_possibles), SAVE :: indnum_fn_num ! donne indice entre entre 1 et niso en fonction du numero d isotope entre 1 et niso_possibles
58  INTEGER, ALLOCATABLE, DIMENSION(:,:), SAVE ::  index_trac ! numero 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    INTEGER, ALLOCATABLE, DIMENSION(:) :: conv_flg_inca
107    INTEGER, ALLOCATABLE, DIMENSION(:) :: pbl_flg_inca
108    CHARACTER(len=8), ALLOCATABLE, DIMENSION(:) :: solsym_inca
109
110    CHARACTER(len=maxlen), ALLOCATABLE, DIMENSION(:) :: tnom_0  ! tracer short name
111    CHARACTER(len=maxlen), ALLOCATABLE, DIMENSION(:) :: tnom_transp ! transporting fluid short name: CRisi
112    CHARACTER(len=3), DIMENSION(30) :: descrq
113    CHARACTER(len=1), DIMENSION(3)  :: txts
114    CHARACTER(len=2), DIMENSION(9)  :: txtp
115    CHARACTER(len=maxlen)           :: str1,str2
116 
117    INTEGER :: nqtrue  ! number of tracers read from tracer.def, without higer order of moment
118    INTEGER :: iq, new_iq, iiq, jq, ierr,itr, iadv
119    INTEGER :: ifils,ipere ! CRisi
120    LOGICAL :: continu,nouveau_traceurdef
121    INTEGER :: IOstatus ! gestion de la retrocompatibilite de traceur.def
122    CHARACTER(len=maxlen) :: tchaine   
123    INTEGER, ALLOCATABLE  :: iqfils(:,:)
124
125    character(len=*),parameter :: modname="infotrac_init"
126
127!-----------------------------------------------------------------------
128! Initialization :
129!
130    txts=(/'x','y','z'/)
131    txtp=(/'x ','y ','z ','xx','xy','xz','yy','yz','zz'/)
132
133    descrq(14)='VLH'
134    descrq(10)='VL1'
135    descrq(11)='VLP'
136    descrq(12)='FH1'
137    descrq(13)='FH2'
138    descrq(16)='PPM'
139    descrq(17)='PPS'
140    descrq(18)='PPP'
141    descrq(20)='SLP'
142    descrq(30)='PRA'
143   
144
145    ! Coherence test between parameter type_trac, config_inca and preprocessing keys
146    IF (type_trac=='inca') THEN
147       WRITE(lunout,*) 'You have chosen to couple with INCA chemistry model : type_trac=', &
148            type_trac,' config_inca=',config_inca
149       IF (config_inca/='aero' .AND. config_inca/='aeNP' .AND. config_inca/='chem') THEN
150          WRITE(lunout,*) 'Incoherence between type_trac and config_inca. Model stops. Modify run.def'
151          CALL abort_gcm('infotrac_init','Incoherence between type_trac and config_inca',1)
152       ENDIF
153#ifndef INCA
154       WRITE(lunout,*) 'To run this option you must add cpp key INCA and compile with INCA code'
155       CALL abort_gcm('infotrac_init','You must compile with cpp key INCA',1)
156#endif
157    ELSE IF (type_trac=='repr') THEN
158       WRITE(lunout,*) 'You have chosen to couple with REPROBUS chemestry model : type_trac=', type_trac
159#ifndef REPROBUS
160       WRITE(lunout,*) 'To run this option you must add cpp key REPROBUS and compile with REPRPBUS code'
161       CALL abort_gcm('infotrac_init','You must compile with cpp key REPROBUS',1)
162#endif
163    ELSE IF (type_trac == 'co2i') THEN
164       WRITE(lunout,*) 'You have chosen to run with CO2 cycle: type_trac=', type_trac
165    ELSE IF (type_trac == 'coag') THEN
166       WRITE(lunout,*) 'Tracers are treated for COAGULATION tests : type_trac=', type_trac
167#ifndef CPP_StratAer
168       WRITE(lunout,*) 'To run this option you must add cpp key StratAer and compile with StratAer code'
169       CALL abort_gcm('infotrac_init','You must compile with cpp key StratAer',1)
170#endif
171    ELSE IF (type_trac == 'lmdz') THEN
172       WRITE(lunout,*) 'Tracers are treated in LMDZ only : type_trac=', type_trac
173    ELSE IF (type_trac == 'inco') THEN ! ThL
174       WRITE(lunout,*) 'Using jointly INCA and CO2 cycle: type_trac =', type_trac
175       IF (config_inca/='aero' .AND. config_inca/='aeNP' .AND. config_inca/='chem') THEN
176          WRITE(lunout,*) 'Incoherence between type_trac and config_inca. Model stops. Modify run.def'
177          CALL abort_gcm('infotrac_init','Incoherence between type_trac and config_inca',1)
178       ENDIF
179#ifndef INCA
180       WRITE(lunout,*) 'To run this option you must add cpp key INCA and compilewith INCA code'
181       CALL abort_gcm('infotrac_init','You must compile with cpp key INCA',1)
182#endif   
183    ELSE
184       WRITE(lunout,*) 'type_trac=',type_trac,' not possible. Model stops'
185       CALL abort_gcm('infotrac_init','bad parameter',1)
186    ENDIF
187
188    ! Test if config_inca is other then none for run without INCA
189    IF (type_trac/='inca' .AND. type_trac/='inco' .AND. config_inca/='none') THEN
190       WRITE(lunout,*) 'config_inca will now be changed to none as you do not couple with INCA model'
191       config_inca='none'
192    ENDIF
193
194!-----------------------------------------------------------------------
195!
196! 1) Get the true number of tracers + water vapor/liquid
197!    Here true tracers (nqtrue) means declared tracers (only first order)
198!
199!-----------------------------------------------------------------------
200    IF (type_trac == 'lmdz' .OR. type_trac == 'repr' .OR. type_trac == 'coag' .OR. type_trac == 'co2i') THEN
201       IF (type_trac=='co2i') THEN
202          nqCO2 = 1
203       ELSE
204          nqCO2 = 0
205       ENDIF
206       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
207       IF(ierr.EQ.0) THEN
208          WRITE(lunout,*) trim(modname),': Open traceur.def : ok'
209          READ(90,*) nqtrue
210          write(lunout,*) 'nqtrue=',nqtrue
211       ELSE
212          WRITE(lunout,*) trim(modname),': Failed opening traceur.def'
213          CALL abort_gcm(modname,"file traceur.def not found!",1)
214       ENDIF
215!jyg<
216!!       if ( planet_type=='earth') then
217!!         ! For Earth, water vapour & liquid tracers are not in the physics
218!!         nbtr=nqtrue-2
219!!       else
220!!         ! Other planets (for now); we have the same number of tracers
221!!         ! in the dynamics than in the physics
222!!         nbtr=nqtrue
223!!       endif
224!>jyg
225    ELSE ! type_trac=inca or inco
226       IF (type_trac=='inco') THEN
227          nqCO2 = 1
228       ELSE
229          nqCO2 = 0
230       ENDIF
231!jyg<
232       ! The traceur.def file is used to define the number "nqo" of water phases
233       ! present in the simulation. Default : nqo = 2.
234       OPEN(90,file='traceur.def',form='formatted',status='old', iostat=ierr)
235       IF(ierr.EQ.0) THEN
236          WRITE(lunout,*) trim(modname),': Open traceur.def : ok'
237          READ(90,*) nqo
238       ELSE
239          WRITE(lunout,*) trim(modname),': Failed opening traceur.def'
240          CALL abort_gcm(modname,"file traceur.def not found!",1)
241       ENDIF
242       IF (nqo /= 2 .AND. nqo /= 3 ) THEN
243          IF (nqo == 4 .AND. type_trac=='inco') THEN ! ThL
244             WRITE(lunout,*) trim(modname),': you are coupling with INCA, and also using CO2i.'
245             nqo = 3    ! A ameliorier... je force 3 traceurs eau...  ThL
246             WRITE(lunout,*) trim(modname),': nqo = ',nqo
247          ELSE
248          WRITE(lunout,*) trim(modname),': nqo=',nqo, ' is not allowed. Only 2 or 3 water phases allowed'
249          CALL abort_gcm('infotrac_init','Bad number of water phases',1)
250          ENDIF
251       ENDIF
252       ! nbtr has been read from INCA by init_const_lmdz() in gcm.F
253#ifdef INCA
254       CALL Init_chem_inca_trac(nqINCA)
255#else
256       nqINCA=0
257#endif
258       nbtr=nqINCA+nqCO2
259       nqtrue=nbtr+nqo
260       WRITE(lunout,*) trim(modname),': nqo = ',nqo
261       WRITE(lunout,*) trim(modname),': nbtr = ',nbtr
262       WRITE(lunout,*) trim(modname),': nqtrue = ',nqtrue
263       WRITE(lunout,*) trim(modname),': nqCO2 = ',nqCO2
264       WRITE(lunout,*) trim(modname),': nqINCA = ',nqINCA
265       ALLOCATE(hadv_inca(nqINCA), vadv_inca(nqINCA), conv_flg_inca(nqINCA), pbl_flg_inca(nqINCA), solsym_inca(nqINCA))
266    ENDIF   ! type_trac 'inca' ou 'inco'
267!>jyg
268
269    IF ((planet_type=="earth").and.(nqtrue < 2)) THEN
270       WRITE(lunout,*) trim(modname),': nqtrue=',nqtrue, ' is not allowed. 2 tracers is the minimum'
271       CALL abort_gcm('infotrac_init','Not enough tracers',1)
272    ENDIF
273   
274!jyg<
275       
276!
277! Allocate variables depending on nqtrue
278!
279    ALLOCATE(tnom_0(nqtrue), hadv(nqtrue), vadv(nqtrue),tnom_transp(nqtrue))
280
281
282!-----------------------------------------------------------------------
283! 2)     Choix  des schemas d'advection pour l'eau et les traceurs
284!
285!     iadv = 1    schema  transport type "humidite specifique LMD"
286!     iadv = 2    schema   amont
287!     iadv = 14   schema  Van-leer + humidite specifique
288!                            Modif F.Codron
289!     iadv = 10   schema  Van-leer (retenu pour l'eau vapeur et liquide)
290!     iadv = 11   schema  Van-Leer pour hadv et version PPM (Monotone) pour vadv
291!     iadv = 12   schema  Frederic Hourdin I
292!     iadv = 13   schema  Frederic Hourdin II
293!     iadv = 16   schema  PPM Monotone(Collela & Woodward 1984)
294!     iadv = 17   schema  PPM Semi Monotone (overshoots autorises)
295!     iadv = 18   schema  PPM Positif Defini (overshoots undershoots autorises)
296!     iadv = 20   schema  Slopes
297!     iadv = 30   schema  Prather
298!
299!        Dans le tableau q(ij,l,iq) : iq = 1  pour l'eau vapeur
300!                                     iq = 2  pour l'eau liquide
301!       Et eventuellement             iq = 3,nqtot pour les autres traceurs
302!
303!        iadv(1): choix pour l'eau vap. et  iadv(2) : choix pour l'eau liq.
304!------------------------------------------------------------------------
305!
306!    Get choice of advection schema from file tracer.def or from INCA
307!---------------------------------------------------------------------
308    IF (type_trac == 'lmdz' .OR. type_trac == 'repr' .OR. type_trac == 'coag' .OR. type_trac == 'co2i') THEN
309
310          ! Continue to read tracer.def
311          DO iq=1,nqtrue
312
313             write(*,*) 'infotrac 237: iq=',iq
314             ! CRisi: ajout du nom du fluide transporteur
315             ! mais rester retro compatible
316             READ(90,'(I2,X,I2,X,A)',IOSTAT=IOstatus) hadv(iq),vadv(iq),tchaine
317             write(lunout,*) 'iq,hadv(iq),vadv(iq)=',iq,hadv(iq),vadv(iq)
318             write(lunout,*) 'tchaine=',trim(tchaine)
319             write(*,*) 'infotrac 238: IOstatus=',IOstatus
320             if (IOstatus.ne.0) then
321                CALL abort_gcm('infotrac_init','Pb dans la lecture de traceur.def',1)
322             endif
323             ! Y-a-t-il 1 ou 2 noms de traceurs? -> On regarde s'il y a un
324             ! espace ou pas au milieu de la chaine.
325             continu=.true.
326             nouveau_traceurdef=.false.
327             iiq=1
328             do while (continu)
329                if (tchaine(iiq:iiq).eq.' ') then
330                  nouveau_traceurdef=.true.
331                  continu=.false.
332                else if (iiq.lt.LEN_TRIM(tchaine)) then
333                  iiq=iiq+1
334                else
335                  continu=.false.
336                endif
337             enddo
338             write(*,*) 'iiq,nouveau_traceurdef=',iiq,nouveau_traceurdef
339             if (nouveau_traceurdef) then
340                write(lunout,*) 'C''est la nouvelle version de traceur.def'
341                tnom_0(iq)=TRIM(tchaine(1:iiq-1))
342                tnom_transp(iq)=TRIM(tchaine(iiq+1:))
343             else
344                write(lunout,*) 'C''est l''ancienne version de traceur.def'
345                write(lunout,*) 'On suppose que les traceurs sont tous d''air'
346                tnom_0(iq)=tchaine
347                tnom_transp(iq) = 'air'
348             endif
349             write(lunout,*) 'tnom_0(iq)=<',trim(tnom_0(iq)),'>'
350             write(lunout,*) 'tnom_transp(iq)=<',trim(tnom_transp(iq)),'>'
351
352          ENDDO!DO iq=1,nqtrue
353
354          CLOSE(90) 
355
356       WRITE(lunout,*) trim(modname),': Valeur de traceur.def :'
357       WRITE(lunout,*) trim(modname),': nombre total de traceurs ',nqtrue
358       DO iq=1,nqtrue
359          WRITE(lunout,*) hadv(iq),vadv(iq),' ',trim(tnom_0(iq)),' ',trim(tnom_transp(iq))
360       END DO
361
362       IF ( planet_type=='earth') THEN
363         !CR: nombre de traceurs de l eau
364         IF (tnom_0(3) == 'H2Oi') THEN
365            nqo=3
366         ELSE
367            nqo=2
368         ENDIF
369         ! For Earth, water vapour & liquid tracers are not in the physics
370         nbtr=nqtrue-nqo
371       ELSE
372         ! Other planets (for now); we have the same number of tracers
373         ! in the dynamics than in the physics
374         nbtr=nqtrue
375       ENDIF
376
377#ifdef CPP_StratAer
378       IF (type_trac == 'coag') THEN
379         nbtr_bin=0
380         nbtr_sulgas=0
381         DO iq=1,nqtrue
382           IF (tnom_0(iq)(1:3)=='BIN') THEN !check if tracer name contains 'BIN'
383             nbtr_bin=nbtr_bin+1
384           ENDIF
385           IF (tnom_0(iq)(1:3)=='GAS') THEN !check if tracer name contains 'GAS'
386             nbtr_sulgas=nbtr_sulgas+1
387           ENDIF
388         ENDDO
389         print*,'nbtr_bin=',nbtr_bin
390         print*,'nbtr_sulgas=',nbtr_sulgas
391         DO iq=1,nqtrue
392           IF (tnom_0(iq)=='GASOCS') THEN
393             id_OCS_strat=iq-nqo
394           ENDIF
395           IF (tnom_0(iq)=='GASSO2') THEN
396             id_SO2_strat=iq-nqo
397           ENDIF
398           IF (tnom_0(iq)=='GASH2SO4') THEN
399             id_H2SO4_strat=iq-nqo
400           ENDIF
401           IF (tnom_0(iq)=='BIN01') THEN
402             id_BIN01_strat=iq-nqo
403           ENDIF
404           IF (tnom_0(iq)=='GASTEST') THEN
405             id_TEST_strat=iq-nqo
406           ENDIF
407         ENDDO
408         print*,'id_OCS_strat  =',id_OCS_strat
409         print*,'id_SO2_strat  =',id_SO2_strat
410         print*,'id_H2SO4_strat=',id_H2SO4_strat
411         print*,'id_BIN01_strat=',id_BIN01_strat
412       ENDIF
413#endif
414
415    ENDIF  ! (type_trac == 'lmdz' .OR. type_trac == 'repr' .OR. type_trac = 'coag' .OR. type_trac = 'co2i')
416!jyg<
417!
418
419! Transfert number of tracers to Reprobus
420    IF (type_trac == 'repr') THEN
421#ifdef REPROBUS
422       CALL Init_chem_rep_trac(nbtr,nqo,tnom_0)
423#endif
424    ENDIF
425!
426! Allocate variables depending on nbtr
427!
428    ALLOCATE(conv_flg(nbtr), pbl_flg(nbtr), solsym(nbtr))
429    conv_flg(:) = 1 ! convection activated for all tracers
430    pbl_flg(:)  = 1 ! boundary layer activated for all tracers
431
432    IF (type_trac == 'inca' .OR. type_trac == 'inco') THEN   ! config_inca='aero' ou 'chem'
433!>jyg
434! le module de chimie fournit les noms des traceurs
435! et les schemas d'advection associes. excepte pour ceux lus
436! dans traceur.def
437
438          DO iq=1,nqo+nqCO2
439
440             write(*,*) 'infotrac 237: iq=',iq
441             ! CRisi: ajout du nom du fluide transporteur
442             ! mais rester retro compatible
443             READ(90,'(I2,X,I2,X,A)',IOSTAT=IOstatus) hadv(iq),vadv(iq),tchaine
444             write(lunout,*) 'iq,hadv(iq),vadv(iq)=',iq,hadv(iq),vadv(iq)
445             write(lunout,*) 'tchaine=',trim(tchaine)
446             write(*,*) 'infotrac 238: IOstatus=',IOstatus
447             if (IOstatus.ne.0) then
448                CALL abort_gcm('infotrac_init','Pb dans la lecture de traceur.def',1)
449             endif
450             ! Y-a-t-il 1 ou 2 noms de traceurs? -> On regarde s'il y a un
451             ! espace ou pas au milieu de la chaine.
452             continu=.true.
453             nouveau_traceurdef=.false.
454             iiq=1
455
456             do while (continu)
457                if (tchaine(iiq:iiq).eq.' ') then
458                  nouveau_traceurdef=.true.
459                  continu=.false.
460                else if (iiq.lt.LEN_TRIM(tchaine)) then
461                  iiq=iiq+1
462                else
463                  continu=.false.
464                endif
465             enddo
466
467             write(*,*) 'iiq,nouveau_traceurdef=',iiq,nouveau_traceurdef
468
469             if (nouveau_traceurdef) then
470                write(lunout,*) 'C''est la nouvelle version de traceur.def'
471                tnom_0(iq)=tchaine(1:iiq-1)
472                tnom_transp(iq)=tchaine(iiq+1:)
473             else
474                write(lunout,*) 'C''est l''ancienne version de traceur.def'
475                write(lunout,*) 'On suppose que les traceurs sont tous d''air'
476                tnom_0(iq)=tchaine
477                tnom_transp(iq) = 'air'
478             endif
479
480             write(lunout,*) 'tnom_0(iq)=<',trim(tnom_0(iq)),'>'
481             write(lunout,*) 'tnom_transp(iq)=<',trim(tnom_transp(iq)),'>'
482
483          ENDDO  !DO iq=1,nqo
484          CLOSE(90) 
485
486 
487#ifdef INCA
488       CALL init_transport( &
489            hadv_inca, &
490            vadv_inca, &
491            conv_flg_inca, &
492            pbl_flg_inca,  &
493            solsym_inca)
494
495       conv_flg(1+nqCO2:nbtr) = conv_flg_inca
496       pbl_flg(1+nqCO2:nbtr) = pbl_flg_inca
497       solsym(1+nqCO2:nbtr) = solsym_inca
498
499       IF (type_trac == 'inco') THEN
500          conv_flg(1:nqCO2) = 1
501          pbl_flg(1:nqCO2) = 1
502          solsym(1:nqCO2) = 'CO2'
503       ENDIF
504#endif
505
506!jyg<
507       DO iq = nqo+nqCO2+1, nqtrue
508          hadv(iq) = hadv_inca(iq-nqo-nqCO2)
509          vadv(iq) = vadv_inca(iq-nqo-nqCO2)
510          tnom_0(iq)=solsym_inca(iq-nqo-nqCO2)
511          tnom_transp(iq) = 'air'
512       END DO
513
514    ENDIF ! (type_trac == 'inca' or 'inco')
515
516!-----------------------------------------------------------------------
517!
518! 3) Verify if advection schema 20 or 30 choosen
519!    Calculate total number of tracers needed: nqtot
520!    Allocate variables depending on total number of tracers
521!-----------------------------------------------------------------------
522    new_iq=0
523    DO iq=1,nqtrue
524       ! Add tracers for certain advection schema
525       IF (hadv(iq)<20 .AND. vadv(iq)<20 ) THEN
526          new_iq=new_iq+1  ! no tracers added
527       ELSE IF (hadv(iq)==20 .AND. vadv(iq)==20 ) THEN
528          new_iq=new_iq+4  ! 3 tracers added
529       ELSE IF (hadv(iq)==30 .AND. vadv(iq)==30 ) THEN
530          new_iq=new_iq+10 ! 9 tracers added
531       ELSE
532          WRITE(lunout,*) trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq)
533          CALL abort_gcm('infotrac_init','Bad choice of advection schema - 1',1)
534       ENDIF
535    END DO
536   
537    IF (new_iq /= nqtrue) THEN
538       ! The choice of advection schema imposes more tracers
539       ! Assigne total number of tracers
540       nqtot = new_iq
541
542       WRITE(lunout,*) trim(modname),': The choice of advection schema for one or more tracers'
543       WRITE(lunout,*) 'makes it necessary to add tracers'
544       WRITE(lunout,*) trim(modname)//': ',nqtrue,' is the number of true tracers'
545       WRITE(lunout,*) trim(modname)//': ',nqtot, ' is the total number of tracers needed'
546
547    ELSE
548       ! The true number of tracers is also the total number
549       nqtot = nqtrue
550    ENDIF
551
552!
553! Allocate variables with total number of tracers, nqtot
554!
555    ALLOCATE(tracers(nqtot))
556    ALLOCATE(niadv(nqtot))
557
558!-----------------------------------------------------------------------
559!
560! 4) Determine iadv, long and short name
561!
562!-----------------------------------------------------------------------
563    new_iq=0
564    DO iq=1,nqtrue
565       new_iq=new_iq+1
566
567       ! Verify choice of advection schema
568       IF (hadv(iq)==vadv(iq)) THEN
569          tracers(new_iq)%iadv=hadv(iq)
570       ELSE IF (hadv(iq)==10 .AND. vadv(iq)==16) THEN
571          tracers(new_iq)%iadv=11
572       ELSE
573          WRITE(lunout,*)trim(modname),': This choice of advection schema is not available',iq,hadv(iq),vadv(iq)
574
575          CALL abort_gcm('infotrac_init','Bad choice of advection schema - 2',1)
576       ENDIF
577     
578       str1=tnom_0(iq)
579       tracers(new_iq)%name = TRIM(tnom_0(iq))
580       tracers(new_iq)%parent = TRIM(tnom_transp(iq))
581       IF (tracers(new_iq)%iadv==0) THEN
582          tracers(new_iq)%longName=trim(str1)
583       ELSE
584          tracers(new_iq)%longName=trim(tnom_0(iq))//descrq(tracers(new_iq)%iadv)
585       ENDIF
586
587       ! schemas tenant compte des moments d'ordre superieur
588       str2=TRIM(tracers(new_iq)%longName)
589       IF (tracers(new_iq)%iadv==20) THEN
590          DO jq=1,3
591             new_iq=new_iq+1
592             tracers(new_iq)%iadv=-20
593             tracers(new_iq)%longName=trim(str2)//txts(jq)
594             tracers(new_iq)%name=trim(str1)//txts(jq)
595          END DO
596       ELSE IF (tracers(new_iq)%iadv==30) THEN
597          DO jq=1,9
598             new_iq=new_iq+1
599             tracers(new_iq)%iadv=-30
600             tracers(new_iq)%longName=trim(str2)//txtp(jq)
601             tracers(new_iq)%name=trim(str1)//txtp(jq)
602          END DO
603       ENDIF
604    END DO
605
606!
607! Find vector keeping the correspodence between true and total tracers
608!
609    niadv(:)=0
610    iiq=0
611    DO iq=1,nqtot
612       IF(tracers(iq)%iadv.GE.0) THEN
613          ! True tracer
614          iiq=iiq+1
615          niadv(iiq)=iq
616       ENDIF
617    END DO
618
619
620    WRITE(lunout,*) trim(modname),': Information stored in infotrac :'
621    WRITE(lunout,*) trim(modname),': iadv  niadv  name  long_name :'
622
623    DO iq=1,nqtot
624       WRITE(lunout,*) tracers(iq)%iadv,niadv(iq), ' ',trim(tracers(iq)%name),' ',trim(tracers(iq)%longName)
625    END DO
626
627!
628! Test for advection schema.
629! This version of LMDZ only garantees iadv=10 and iadv=14 (14 only for water vapour) .
630!
631    DO iq=1,nqtot
632       iadv=tracers(iq)%iadv
633       IF (ALL([10, 14, 0]/=iadv)) THEN
634          WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv,' is not tested in this version of LMDZ'
635          CALL abort_gcm('infotrac_init','In this version only iadv=10 and iadv=14 is tested!',1)
636       ELSE IF (iadv==14 .AND. iq/=1) THEN
637          WRITE(lunout,*)trim(modname),'STOP : The option iadv=',iadv,' is not tested in this version of LMDZ'
638          CALL abort_gcm('infotrac_init','In this version iadv=14 is only permitted for water vapour!',1)
639       ENDIF
640    END DO
641
642
643! CRisi: quels sont les traceurs fils et les traceurs peres.
644! initialiser tous les tableaux d'indices lies aux traceurs familiaux
645! + verifier que tous les peres sont ecrits en premieres positions
646    ALLOCATE(iqfils(nqtot,nqtot))   
647    nqperes=0
648    iqfils(:,:)=0
649    tracers(:)%iqParent=0
650    DO iq=1,nqtot
651      if (tnom_transp(iq) == 'air') then
652        ! ceci est un traceur pere
653        WRITE(lunout,*) 'Le traceur',iq,', appele ',trim(tnom_0(iq)),', est un pere'
654        nqperes=nqperes+1
655        tracers(iq)%iqParent=0
656      else !if (tnom_transp(iq) == 'air') then
657        ! ceci est un fils. Qui est son pere?
658        WRITE(lunout,*) 'Le traceur',iq,', appele ',trim(tnom_0(iq)),', est un fils'
659        continu=.true.
660        ipere=1
661        do while (continu)           
662          if (tnom_transp(iq) == tnom_0(ipere)) then
663            ! Son pere est ipere
664            WRITE(lunout,*) 'Le traceur',iq,'appele ', &
665      &          trim(tnom_0(iq)),' est le fils de ',ipere,'appele ',trim(tnom_0(ipere))
666            if (iq.eq.ipere) then
667                CALL abort_gcm('infotrac_init','Un fils est son propre pere',1)
668            endif
669            tracers(ipere)%nqChilds = tracers(ipere)%nqChilds+1 
670            iqfils(tracers(ipere)%nqChilds,ipere)=iq
671            tracers(iq)%iqParent=ipere         
672            continu=.false.
673          else !if (tnom_transp(iq) == tnom_0(ipere)) then
674            ipere=ipere+1
675            if (ipere.gt.nqtot) then
676                WRITE(lunout,*) 'Le traceur',iq,'appele ', &
677      &          trim(tnom_0(iq)),', est orphelin.'
678                CALL abort_gcm('infotrac_init','Un traceur est orphelin',1)
679            endif !if (ipere.gt.nqtot) then
680          endif !if (tnom_transp(iq) == tnom_0(ipere)) then
681        enddo !do while (continu)
682      endif !if (tnom_transp(iq) == 'air') then
683    enddo !DO iq=1,nqtot
684    WRITE(lunout,*) 'infotrac: nqperes=',nqperes   
685    WRITE(lunout,*) 'nqChilds=',tracers(:)%nqChilds
686    WRITE(lunout,*) 'iqParent=',tracers(:)%iqParent
687    WRITE(lunout,*) 'iqfils=',iqfils
688
689! Calculer le nombre de descendants a partir de iqfils et de nbfils
690    DO iq=1,nqtot   
691      tracers(iq)%iGeneration=0
692      continu=.true.
693      ifils=iq
694      do while (continu)
695        ipere=tracers(ifils)%iqParent
696        if (ipere.gt.0) then
697         tracers(ipere)%nqDescen = tracers(ipere)%nqDescen+1   
698         iqfils(tracers(ipere)%nqDescen,ipere)=iq
699         ifils=ipere
700         tracers(iq)%iGeneration=tracers(iq)%iGeneration+1
701        else !if (ipere.gt.0) then
702         continu=.false.
703        endif !if (ipere.gt.0) then
704      enddo !do while (continu)   
705      WRITE(lunout,*) 'Le traceur ',iq,', appele ',trim(tnom_0(iq)), &
706        ' est un traceur de generation: ',tracers(iq)%iGeneration
707    enddo !DO iq=1,nqtot
708    DO iq=1,nqtot
709      ALLOCATE(tracers(iq)%iqDescen(tracers(iq)%nqDescen))
710      tracers(iq)%iqDescen(:) = iqfils(1:tracers(iq)%nqDescen,iq)
711    END DO
712
713    WRITE(lunout,*) 'infotrac: nqDescen=',tracers(iq)%nqDescen
714    WRITE(lunout,*) 'iqfils=',iqfils
715    WRITE(lunout,*) 'nqDescen_tot=',SUM(tracers(:)%nqDescen)
716
717! Interdire autres schemas que 10 pour les traceurs fils, et autres schemas
718! que 10 et 14 si des peres ont des fils
719    do iq=1,nqtot
720      if (tracers(iq)%iqParent > 0) then
721        ! ce traceur a un pere qui n'est pas l'air
722        ! Seul le schema 10 est autorise
723        iadv=tracers(iq)%iadv
724        if (iadv/=10) then
725           WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv,' is not implemented for sons'
726          CALL abort_gcm('infotrac_init','Sons should be advected by scheme 10',1)
727        endif
728        ! Le traceur pere ne peut etre advecte que par schema 10 ou 14:
729        IF (ALL([10,14]/=tracers(tracers(iq)%iqParent)%iadv)) THEN
730          WRITE(lunout,*)trim(modname),' STOP : The option iadv=',iadv,' is not implemented for fathers'
731          CALL abort_gcm('infotrac_init','Fathers should be advected by scheme 10 ou 14',1)
732        endif
733     endif
734    enddo !do iq=1,nqtot
735    tracers(:)%gen0Name = ancestor(tracers)      !--- Name of the first generation ancestor
736
737
738! detecter quels sont les traceurs isotopiques parmi des traceurs
739    call infotrac_isoinit(tnom_0,nqtrue)
740
741!    if (ntraciso.gt.0) then
742! le 18 sep 2020: on enleve la condition ntraciso.gt.0 car nqtottr doit etre
743! connu meme si il n'y a pas d'isotopes!
744        write(lunout,*) 'infotrac 702: nbtr,ntraciso=',nbtr,ntraciso
745! retrancher les traceurs isotopiques de la liste des traceurs qui passent dans
746! phytrac
747        nbtr=nbtr-nqo*ntraciso
748
749! faire un tableau d'indice des traceurs qui passeront dans phytrac
750        nqtottr=nqtot-nqo*(1+ntraciso)
751        write(lunout,*) 'infotrac 704: nqtottr,nqtot,nqo=',nqtottr,nqtot,nqo
752        ! Rq: nqtottr n'est pas forcement egal a nbtr dans le cas ou new_iq /= nqtrue
753        ALLOCATE (itr_indice(nqtottr)) 
754        itr_indice(:)=0 
755        itr=0
756        do iq=nqo+1, nqtot
757          if (tracers(iq)%iso_iName.eq.0) then
758            itr=itr+1
759            write(*,*) 'itr=',itr
760            itr_indice(itr)=iq
761          endif !if (tracers(iq)%iso_iName.eq.0) then
762        enddo
763        if (itr.ne.nqtottr) then
764            CALL abort_gcm('infotrac_init','pb dans le calcul de nqtottr',1)
765        endif
766        write(lunout,*) 'itr_indice=',itr_indice
767!    endif !if (ntraciso.gt.0) then
768
769!-----------------------------------------------------------------------
770! Finalize :
771!
772    DEALLOCATE(tnom_0, hadv, vadv,tnom_transp)
773
774    WRITE(lunout,*) 'infotrac init fin'
775
776  END SUBROUTINE infotrac_init
777
778  SUBROUTINE infotrac_isoinit(tnom_0,nqtrue)
779
780#ifdef CPP_IOIPSL
781  use IOIPSL
782#else
783  ! if not using IOIPSL, we still need to use (a local version of) getin
784  use ioipsl_getincom
785#endif
786  implicit none
787 
788    ! inputs
789    INTEGER,INTENT(IN) :: nqtrue
790    CHARACTER(len=*),INTENT(IN) :: tnom_0(nqtrue)
791   
792    ! locals   
793    CHARACTER(len=3), DIMENSION(niso_possibles) :: tnom_iso
794    INTEGER, ALLOCATABLE,DIMENSION(:,:) :: nb_iso,nb_traciso
795    INTEGER, ALLOCATABLE,DIMENSION(:) :: nb_isoind
796    INTEGER :: ntraceurs_zone_prec,iq,phase,ixt,iiso,izone
797    CHARACTER(len=maxlen) :: tnom_trac
798    INCLUDE "iniprint.h"
799
800    tnom_iso=(/'eau','HDO','O18','O17','HTO'/)
801
802    ALLOCATE(nb_iso(niso_possibles,nqo))
803    ALLOCATE(nb_isoind(nqo))
804    ALLOCATE(nb_traciso(niso_possibles,nqo))
805    ALLOCATE(iso_indnum(nqtot))
806     
807    iso_indnum(:)=0
808    indnum_fn_num(:)=0
809    use_iso(:)=.false. 
810    nb_iso(:,:)=0 
811    nb_isoind(:)=0     
812    nb_traciso(:,:)=0
813    niso=0
814    ntraceurs_zone=0 
815    ntraceurs_zone_prec=0
816    ntraciso=0
817
818    do iq=nqo+1,nqtot
819!       write(lunout,*) 'infotrac 569: iq,tnom_0(iq)=',iq,tnom_0(iq)
820       do phase=1,nqo   
821        do ixt= 1,niso_possibles   
822         tnom_trac=trim(tnom_0(phase))//'_'
823         tnom_trac=trim(tnom_trac)//trim(tnom_iso(ixt))
824!         write(*,*) 'phase,ixt,tnom_trac=',phase,ixt,tnom_trac     
825         IF (tnom_0(iq) == tnom_trac) then
826!          write(lunout,*) 'Ce traceur est un isotope'
827          nb_iso(ixt,phase)=nb_iso(ixt,phase)+1   
828          nb_isoind(phase)=nb_isoind(phase)+1   
829          tracers(iq)%iso_iName=ixt
830          iso_indnum(iq)=nb_isoind(phase)
831          indnum_fn_num(ixt)=iso_indnum(iq)
832          tracers(iq)%iso_iPhase=phase
833          goto 20
834         else if ( tracers(iq)%iqParent> 0) then         
835          if (tnom_0(tracers(iq)%iqParent) == tnom_trac) then
836!           write(lunout,*) 'Ce traceur est le fils d''un isotope'
837           ! c'est un traceur d'isotope
838           nb_traciso(ixt,phase)=nb_traciso(ixt,phase)+1
839           tracers(iq)%iso_iName=ixt
840           iso_indnum(iq)=indnum_fn_num(ixt)
841           tracers(iq)%iso_iZone=nb_traciso(ixt,phase)
842           tracers(iq)%iso_iPhase=phase
843           goto 20
844          endif !if (tnom_0(tracers(iq)%iqParent) == trim(tnom_0(phase))//trim(tnom_iso(ixt))) then
845         endif !IF (tnom_0(iq) == trim(tnom_0(phase))//trim(tnom_iso(ixt))) then
846        enddo !do ixt= niso_possibles
847       enddo !do phase=1,nqo
848  20   continue
849      enddo !do iq=1,nqtot
850
851      do ixt= 1,niso_possibles 
852
853        if (nb_iso(ixt,1).eq.1) then
854          ! on verifie que toutes les phases ont le meme nombre de
855          ! traceurs
856          do phase=2,nqo
857            if (nb_iso(ixt,phase).ne.nb_iso(ixt,1)) then
858!              write(lunout,*) 'ixt,phase,nb_iso=',ixt,phase,nb_iso(ixt,phase)
859              CALL abort_gcm('infotrac_init','Phases must have same number of isotopes',1)
860            endif
861          enddo !do phase=2,nqo
862
863          niso=niso+1
864          use_iso(ixt)=.true.
865          ntraceurs_zone=nb_traciso(ixt,1)
866
867          ! on verifie que toutes les phases ont le meme nombre de
868          ! traceurs
869          do phase=2,nqo
870            if (nb_traciso(ixt,phase).ne.ntraceurs_zone) then
871              write(lunout,*) 'ixt,phase,nb_traciso=',ixt,phase,nb_traciso(ixt,phase)
872              write(lunout,*) 'ntraceurs_zone=',ntraceurs_zone
873              CALL abort_gcm('infotrac_init','Phases must have same number of tracers',1)
874            endif 
875          enddo  !do phase=2,nqo
876          ! on verifie que tous les isotopes ont le meme nombre de
877          ! traceurs
878          if (ntraceurs_zone_prec.gt.0) then               
879            if (ntraceurs_zone.eq.ntraceurs_zone_prec) then
880              ntraceurs_zone_prec=ntraceurs_zone
881            else !if (ntraceurs_zone.eq.ntraceurs_zone_prec) then
882              write(*,*) 'ntraceurs_zone_prec,ntraceurs_zone=',ntraceurs_zone_prec,ntraceurs_zone   
883              CALL abort_gcm('infotrac_init', &
884               &'Isotope tracers are not well defined in traceur.def',1)           
885            endif !if (ntraceurs_zone.eq.ntraceurs_zone_prec) then
886           endif !if (ntraceurs_zone_prec.gt.0) then
887
888        else if (nb_iso(ixt,1).ne.0) then
889           WRITE(lunout,*) 'nqo,ixt=',nqo,ixt
890           WRITE(lunout,*) 'nb_iso(ixt,1)=',nb_iso(ixt,1)   
891           CALL abort_gcm('infotrac_init','Isotopes are not well defined in traceur.def',1)     
892        endif   !if (nb_iso(ixt,1).eq.1) then       
893    enddo ! do ixt= niso_possibles
894
895    ! dimensions isotopique:
896    ntraciso=niso*(ntraceurs_zone+1)
897!    WRITE(lunout,*) 'niso=',niso
898!    WRITE(lunout,*) 'ntraceurs_zone,ntraciso=',ntraceurs_zone,ntraciso   
899 
900    ! flags isotopiques:
901    ok_isotopes = niso > 0
902!    WRITE(lunout,*) 'ok_isotopes=',ok_isotopes
903 
904    if (ok_isotopes) then
905        ok_iso_verif=.false.
906        call getin('ok_iso_verif',ok_iso_verif)
907        ok_init_iso=.false.
908        call getin('ok_init_iso',ok_init_iso)
909        tnat=(/1.0,155.76e-6,2005.2e-6,0.004/100.,0.0/)
910        alpha_ideal=(/1.0,1.01,1.006,1.003,1.0/)
911    endif !if (ok_isotopes) then 
912!    WRITE(lunout,*) 'ok_iso_verif=',ok_iso_verif
913!    WRITE(lunout,*) 'ok_init_iso=',ok_init_iso
914
915    if (ntraceurs_zone.gt.0) then
916        ok_isotrac=.true.
917    else
918        ok_isotrac=.false.
919    endif   
920!    WRITE(lunout,*) 'ok_isotrac=',ok_isotrac
921
922    ! remplissage du tableau iqiso(ntraciso,phase)
923    ALLOCATE(iqiso(ntraciso,nqo))   
924    iqiso(:,:)=0     
925    do iq=1,nqtot
926        if (tracers(iq)%iso_iName > 0) then
927          ixt=iso_indnum(iq)+tracers(iq)%iso_iZone*niso
928          iqiso(ixt,tracers(iq)%iso_iPhase)=iq
929        endif
930    enddo
931!    WRITE(lunout,*) 'iqiso=',iqiso
932
933    ! replissage du tableau index_trac(ntraceurs_zone,niso)
934    ALLOCATE(index_trac(ntraceurs_zone,niso)) 
935    if (ok_isotrac) then
936        do iiso=1,niso
937          do izone=1,ntraceurs_zone
938             index_trac(izone,iiso)=iiso+izone*niso
939          enddo
940        enddo
941    else !if (ok_isotrac) then     
942        index_trac(:,:)=0.0
943    endif !if (ok_isotrac) then
944!    write(lunout,*) 'index_trac=',index_trac   
945
946! Finalize :
947    DEALLOCATE(nb_iso)
948
949  END SUBROUTINE infotrac_isoinit
950
951END MODULE infotrac
Note: See TracBrowser for help on using the repository browser.