source: LMDZ5/branches/testing/libf/dyn3d_common/infotrac.F90 @ 2453

Last change on this file since 2453 was 2408, checked in by Laurent Fairhead, 9 years ago

Merged trunk changes r2298:2396 into testing branch

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