source: LMDZ6/branches/LMDZ-INCA-Dyn/libf/dyn3d_common/infotrac.F90 @ 4251

Last change on this file since 4251 was 3850, checked in by acozic, 4 years ago

arguments order was changed in inca init_transport routine

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