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

Last change on this file since 4064 was 4064, checked in by dcugnet, 2 years ago
  • minor fixes (unused variables suppressed, comas after a WRITE() statement, etc.)
  • parser routines taken from version 7 of https://svn.lmd.jussieu.fr/tracers-parser
  • few changes in infotrac, and few fixes of (at least) the sequential version:
    • uadv and vadv were deallocated twice (fix was lost by mistake just before last commit)
    • in [( dum(im), im=1, nm)] implicit loops, ifort evaluates "dum(im)" even if nm==0,

resulting in a crash, "im" being unitialized.

  • 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: 50.0 KB
Line 
1!$Id: infotrac.F90 4064 2022-01-25 18:14:39Z dcugnet $
2!
3MODULE infotrac
4
5   USE       strings_mod, ONLY: msg, find, strIdx,  strFind, strParse, dispTable, int2str,  reduceExpr,  &
6                          cat, fmsg, test, strTail, strHead, strStack, strReduce, bool2str, maxlen, testFile
7   USE readTracFiles_mod, ONLY: trac_type, readTracersFiles, addPhase,  phases_sep,  nphases, ancestor,  &
8                                isot_type, readIsotopesFile, delPhase,   old_phases, getKey_init, tran0, &
9                                keys_type, initIsotopes,  indexUpdate, known_phases, getKey, setGeneration
10
11   IMPLICIT NONE
12
13   PRIVATE
14
15   !=== FOR TRACERS:
16   PUBLIC :: infotrac_init                                 !--- Initialization of the tracers
17   PUBLIC :: tracers, type_trac                            !--- Full tracers database, tracers type keyword
18   PUBLIC :: nqtot,   nbtr,   nqo,   nqCO2,   nqtottr      !--- Main dimensions
19   PUBLIC :: solsym, conv_flg, pbl_flg                     !--- Tracers names + convection & boundary layer activation keys
20
21   !=== FOR ISOTOPES: General
22   PUBLIC :: isotopes,  nbIso                              !--- Derived type, full isotopes families database + nb of families
23   PUBLIC :: isoSelect, ixIso                              !--- Isotopes family selection tool + selected family index
24   PUBLIC :: min_qParent, min_qMass, min_ratio             !--- Min. values for various isotopic quantities
25   !=== FOR ISOTOPES: Specific to water
26   PUBLIC :: iH2O, tnat, alpha_ideal                       !--- H2O isotopes index, natural abundance, fractionning coeff.
27   !=== FOR ISOTOPES: Depending on the selected isotopes family
28   PUBLIC :: isotope, isoKeys                              !--- Selected isotopes database + associated keys (cf. getKey)
29   PUBLIC :: isoName, isoZone, isoPhas                     !--- Isotopes and tagging zones names, phases
30   PUBLIC :: niso,    nzone,   nphas,   ntiso              !---  " " numbers + isotopes & tagging tracers number
31   PUBLIC :: iZonIso, iTraPha                              !--- 2D index tables to get "iq" index
32   PUBLIC :: isoCheck                                      !--- Run isotopes checking routines
33   !=== FOR BOTH TRACERS AND ISOTOPES
34   PUBLIC :: getKey                                        !--- Get a key from "tracers" or "isotope"
35
36   PUBLIC :: ntraciso, ntraceurs_zone, iqiso
37   PUBLIC :: ok_isotopes, ok_iso_verif, ok_isotrac, ok_init_iso, use_iso
38   PUBLIC :: index_trac, iso_indnum, indnum_fn_num, niso_possibles
39   PUBLIC :: qperemin, masseqmin, ratiomin
40
41   INTERFACE isoSelect; MODULE PROCEDURE isoSelectByIndex, isoSelectByName; END INTERFACE isoSelect
42
43!=== CONVENTIONS FOR TRACERS NUMBERS:
44!  |--------------------+-----------------------+-----------------+---------------+----------------------------|
45!  | water in different |    water tagging      |  water isotopes | other tracers | additional tracers moments |
46!  | phases: H2O_[gls]  |      isotopes         |                 |               |  for higher order schemes  |
47!  |--------------------+-----------------------+-----------------+---------------+----------------------------|
48!  |                    |                       |                 |               |                            |
49!  |<--     nqo      -->|<-- nqo*niso* nzone -->|<-- nqo*niso  -->|<--  nbtr   -->|<--        (nmom)        -->|         
50!  |                    |                                         |                                            |
51!  |                    |<-- nqo*niso*(nzone+1)  =   nqo*ntiso -->|<--    nqtottr = nbtr + nmom             -->|
52!  |                                                                              = nqtot - nqo*(ntiso+1)      |
53!  |                                                                                                           |
54!  |<--                        nqtrue  =  nbtr + nqo*(ntiso+1)                 -->|                            |
55!  |                                                                                                           |
56!  |<--                        nqtot   =  nqtrue + nmom                                                     -->|
57!  |                                                                                                           |
58!  |-----------------------------------------------------------------------------------------------------------|
59!  NOTES FOR THIS TABLE:
60!  * Used "niso", "nzone" and "ntiso" are components of "isotopes(ip)" for water (isotopes(ip)%parent == 'H2O'),
61!    since water is so far the sole tracers family, except passive CO2, removed from the main tracers table.
62!  * For water, "nqo" is equal to the more general field "isotopes(ip)%nphas".
63!  * "niso", "nzone", "ntiso", "nphas" are defined for other isotopic tracers families, if any.
64!
65!=== DERIVED TYPE EMBEDDING MOST OF THE TRACERS-RELATED QUANTITIES (LENGTH: nqtot)
66!    Each entry is accessible using "%" sign.
67!  |-------------+------------------------------------------------------+-------------+------------------------+
68!  |  entry      | Meaning                                              | Former name | Possible values        |
69!  |-------------+------------------------------------------------------+-------------+------------------------+
70!  | name        | Name (short)                                         | tname       |                        |
71!  | gen0Name    | Name of the 1st generation ancestor                  | /           |                        |
72!  | parent      | Name of the parent                                   | /           |                        |
73!  | longName    | Long name (with adv. scheme suffix) for outputs      | ttext       |                        |
74!  | type        | Type (so far: tracer or tag)                         | /           | tracer,tag             |
75!  | phase       | Phases list ("g"as / "l"iquid / "s"olid)             | /           | [g][l][s]              |
76!  | component   | Name(s) of the merged/cumulated section(s)           | /           | coma-separated names   |
77!  | iadv        | Advection scheme number                              | iadv        | 1-20,30 exc. 3-9,15,19 |
78!  | iGeneration | Generation (>=1)                                     | /           |                        |
79!  | isAdvected  | advected tracers flag (.TRUE. if iadv > 0)           | /           | nqtrue  .TRUE. values  |
80!  | isH2Ofamily | H2O tracer/isotope/tagging tracer                    | /           | nqtottr .TRUE. values  |
81!  | iqParent    | Index of the parent tracer                           | iqpere      | 1:nqtot                |
82!  | iqDescen    | Indexes of the childs       (all generations)        | iqfils      | 1:nqtot                |
83!  | nqDescen    | Number of the descendants   (all generations)        | nqdesc      | 1:nqtot                |
84!  | nqChilds    | Number of childs            (1st generation only)    | nqfils      | 1:nqtot                |
85!  | iso_iGroup  | Isotopes group index in isotopes(:)                  | /           | 1:nbIso                |
86!  | iso_iName   | Isotope  name  index in isotopes(iso_iGroup)%trac(:) | iso_indnum  | 1:niso                 |
87!  | iso_iZone   | Isotope  zone  index in isotopes(iso_iGroup)%zone(:) | zone_num    | 1:nzone                |
88!  | iso_iPhas   | Isotope  phase index in isotopes(iso_iGroup)%phas(:) | phase_num   | 1:nphas                |
89!  | keys        | key/val pairs accessible with "getKey" routine       | /           |                        |
90!  +-------------+------------------------------------------------------+-------------+------------------------+
91!
92!=== DERIVED TYPE EMBEDDING MOST OF THE ISOTOPES-RELATED QUANTITIES (LENGTH: nbIso, NUMBER OF ISOTOPES FAMILIES)
93!    Each entry is accessible using "%" sign.
94!  |-----------------+--------------------------------------------------+----------------+-----------------+
95!  |  entry | length | Meaning                                          | Former name    | Possible values |
96!  |-----------------+--------------------------------------------------+----------------+-----------------+
97!  | parent          | Parent tracer (isotopes family name)             |                |                 |
98!  | keys   | niso   | Isotopes keys/values pairs list + number         |                |                 |
99!  | trac   | ntiso  | Isotopes + tagging tracers list + number         |                |                 |
100!  | zone   | nzone  | Geographic tagging zones   list + number         |                |                 |
101!  | phase  | nphas  | Phases                     list + number         |                | [g][l][s], 1:3  |
102!  | niso            | Number of isotopes, excluding tagging tracers    |                |                 |
103!  | ntiso           | Number of isotopes, including tagging tracers    | ntraciso       |                 |
104!  | nzone           | Number of geographic tagging zones               | ntraceurs_zone |                 |
105!  | nphas           | Number of phases                                 |                |                 |
106!  | iTraPha         | Index in "trac(1:niso)" = f(name(1:ntiso)),phas) | iqiso          | 1:niso          |
107!  | iZonIso         | Index in "trac(1:ntiso)" = f(zone, name(1:niso)) | index_trac     | 1:nzone         |
108!  |-----------------+--------------------------------------------------+----------------+-----------------+
109
110   REAL, PARAMETER :: min_qParent = 1.e-30, min_qMass = 1.e-18, min_ratio = 1.e-16 ! MVals et CRisi
111
112   !=== DIMENSIONS OF THE TRACERS TABLES AND OTHER SCALAR VARIABLES
113   INTEGER,                 SAVE :: nqtot,  &                   !--- Tracers nb in dynamics (incl. higher moments + H2O)
114                                    nbtr,   &                   !--- Tracers nb in physics  (excl. higher moments + H2O)
115                                    nqo,    &                   !--- Number of water phases
116                                    nbIso,  &                   !--- Number of available isotopes family
117                                    nqtottr, &                  !--- Number of tracers passed to phytrac (TO BE DELETED ?)
118                                    nqCO2                       !--- Number of tracers of CO2  (ThL)
119   CHARACTER(LEN=maxlen),   SAVE :: type_trac                   !--- Keyword for tracers type
120
121   !=== DERIVED TYPES EMBEDDING MOST INFORMATIONS ABOUT TRACERS AND ISOTOPES FAMILIES
122   TYPE(trac_type), TARGET, SAVE, ALLOCATABLE ::  tracers(:)    !=== TRACERS DESCRIPTORS VECTOR
123   TYPE(isot_type), TARGET, SAVE, ALLOCATABLE :: isotopes(:)    !=== ISOTOPES PARAMETERS VECTOR
124
125   !=== ALIASES FOR CURRENTLY SELECTED ISOTOPES FAMILY OF VARIABLES EMBEDDED IN THE VECTOR "isotopes"
126   TYPE(isot_type),         SAVE, POINTER :: isotope            !--- CURRENTLY SELECTED ISOTOPES FAMILY DESCRIPTOR
127   INTEGER,                 SAVE          :: ixIso, iH2O        !--- Index of the selected isotopes family and H2O family
128   LOGICAL,                 SAVE          :: isoCheck           !--- Flag to trigger the checking routines
129   TYPE(keys_type),         SAVE, POINTER :: isoKeys(:)         !--- ONE SET OF KEYS FOR EACH ISOTOPE (LISTED IN isoName)
130   CHARACTER(LEN=maxlen),   SAVE, POINTER :: isoName(:),   &    !--- ISOTOPES NAMES FOR THE CURRENTLY SELECTED FAMILY
131                                             isoZone(:),   &    !--- TAGGING ZONES  FOR THE CURRENTLY SELECTED FAMILY
132                                             isoPhas            !--- USED PHASES    FOR THE CURRENTLY SELECTED FAMILY
133   INTEGER, TARGET,         SAVE          ::  niso, nzone, &    !--- NUMBER OF ISOTOPES, TAGGING ZONES AND PHASES
134                                             nphas, ntiso       !--- NUMBER OF PHASES AND ISOTOPES + ISOTOPIC TAGGING TRACERS
135   INTEGER,                 SAVE, POINTER :: iZonIso(:,:)       !--- INDEX IN "isoTrac" AS f(tagging zone, isotope)
136   INTEGER,                 SAVE, POINTER :: iTraPha(:,:)       !--- INDEX IN "isoTrac" AS f(isotopic tracer, phase)
137   INTEGER, ALLOCATABLE,    SAVE ::  index_trac(:,:) ! numero ixt en fn izone, indnum entre 1 et niso
138   INTEGER, ALLOCATABLE,    SAVE ::  iqiso(:,:)      ! donne indice iq en fn de (ixt,phase)
139
140   !--- Aliases for older names
141   INTEGER, POINTER, SAVE :: ntraciso, ntraceurs_zone
142   REAL,             SAVE :: qperemin, masseqmin, ratiomin
143
144! CRisi: cas particulier des isotopes
145   INTEGER, PARAMETER :: niso_possibles = 5
146   LOGICAL, SAVE      :: ok_isotopes, ok_iso_verif, ok_isotrac, ok_init_iso
147   LOGICAL, SAVE, ALLOCATABLE ::       use_iso(:)
148   INTEGER, SAVE, ALLOCATABLE ::    iso_indnum(:)     !--- Gives 1<=idx<=niso_possibles as function(1<=iq <=nqtot)
149   INTEGER, SAVE, ALLOCATABLE :: indnum_fn_num(:)     !--- Gives 1<=idx<=niso           as function(1<=idx<=niso_possibles)
150
151   !=== VARIABLES FOR ISOTOPES INITIALIZATION AND FOR INCA
152   REAL,               SAVE, ALLOCATABLE ::     tnat(:),  &     !--- Natural relative abundance of water isotope        (niso)
153                                         alpha_ideal(:)         !--- Ideal fractionning coefficient (for initial state) (niso)
154   INTEGER,            SAVE, ALLOCATABLE :: conv_flg(:),  &     !--- Convection     activation ; needed for INCA        (nbtr)
155                                             pbl_flg(:)         !--- Boundary layer activation ; needed for INCA        (nbtr)
156   CHARACTER(LEN=8),   SAVE, ALLOCATABLE ::   solsym(:)         !--- Names from INCA                                    (nbtr)
157   LOGICAL, PARAMETER :: lOldCode = .TRUE.
158
159CONTAINS
160
161SUBROUTINE infotrac_init
162   USE control_mod, ONLY: planet_type, config_inca
163#ifdef REPROBUS
164   USE CHEM_REP,    ONLY: Init_chem_rep_trac
165#endif
166   IMPLICIT NONE
167!==============================================================================================================================
168!
169!   Auteur:  P. Le Van /L. Fairhead/F.Hourdin
170!   -------
171!
172!   Modifications:
173!   --------------
174!   05/94: F.Forget      Modif special traceur
175!   02/02: M-A Filiberti Lecture de traceur.def
176!   01/22: D. Cugnet     Nouveaux tracer.def et tracer_*.def + encapsulation (types tr et iso)
177!
178!   Objet:
179!   ------
180!   GCM LMD nouvelle grille
181!
182!==============================================================================================================================
183!   ... modification de l'integration de q ( 26/04/94 ) ....
184!------------------------------------------------------------------------------------------------------------------------------
185! Declarations:
186   INCLUDE "dimensions.h"
187   INCLUDE "iniprint.h"
188
189!------------------------------------------------------------------------------------------------------------------------------
190! Local variables
191   INTEGER, ALLOCATABLE :: hadv(:), vadv(:)                          !--- Horizontal/vertical transport scheme number
192#ifdef INCA
193   INTEGER, ALLOCATABLE :: had (:), hadv_inca(:), conv_flg_inca(:), &!--- Variables specific to INCA
194                           vad (:), vadv_inca(:),  pbl_flg_inca(:)
195   CHARACTER(LEN=8), ALLOCATABLE :: solsym_inca(:)                   !--- Tracers names for INCA
196   INTEGER :: nqINCA
197#endif
198   CHARACTER(LEN=2)      ::   suff(9)                                !--- Suffixes for schemes of order 3 or 4 (Prather)
199   CHARACTER(LEN=3)      :: descrq(30)                               !--- Advection scheme description tags
200   CHARACTER(LEN=maxlen) :: msg1                                     !--- String for messages
201   CHARACTER(LEN=maxlen), ALLOCATABLE :: str(:)                      !--- Temporary storage
202   INTEGER :: fType                                                  !--- Tracers description file type ; 0: none
203                                                                     !--- 1/2/3: "traceur.def"/"tracer.def"/"tracer_*.def"
204   INTEGER :: nqtrue                                                 !--- Tracers nb from tracer.def (no higher order moments)
205   INTEGER :: iad                                                    !--- Advection scheme number
206   INTEGER :: ic, ip, np, iq, jq, it, nt, im, nm, ix, iz, nz         !--- Indexes and temporary variables
207   LOGICAL :: lerr, ll
208   CHARACTER(LEN=1) :: p
209   TYPE(trac_type), ALLOCATABLE, TARGET :: ttr(:)
210   TYPE(trac_type), POINTER             :: t1, t(:)
211   TYPE(isot_type), POINTER             :: iso
212
213   CHARACTER(LEN=maxlen), ALLOCATABLE :: tnom_0(:), tnom_transp(:)        !--- Tracer short name + transporting fluid name
214   CHARACTER(LEN=maxlen)              :: tchaine
215   INTEGER :: ierr
216   LOGICAL :: lINCA
217
218   CHARACTER(LEN=*), PARAMETER :: modname="infotrac_init"
219!------------------------------------------------------------------------------------------------------------------------------
220! Initialization :
221!------------------------------------------------------------------------------------------------------------------------------
222   suff          = ['x ','y ','z ','xx','xy','xz','yy','yz','zz']
223   descrq( 1: 2) = ['LMV','BAK']
224   descrq(10:20) = ['VL1','VLP','FH1','FH2','VLH','   ','PPM','PPS','PPP','   ','SLP']
225   descrq(30)    =  'PRA'
226   
227   CALL msg('type_trac='//TRIM(type_trac), modname)
228   IF(lOldCode) THEN
229      str = [type_trac]; nt = 1
230   ELSE
231      IF(strParse(type_trac, ',', str, n=nt)) CALL abort_gcm(modname,'can''t parse "type_trac = '//TRIM(type_trac)//'"',1)
232   END IF
233
234   !---------------------------------------------------------------------------------------------------------------------------
235   DO it = 1, nt                                                          !--- nt>1=> "type_trac": coma-separated keywords list
236   !---------------------------------------------------------------------------------------------------------------------------
237      !--- MESSAGE ABOUT THE CHOSEN CONFIGURATION
238      msg1 = 'For type_trac = "'//TRIM(str(it))//'":'
239      SELECT CASE(type_trac)
240         CASE('inca'); CALL msg(TRIM(msg1)//' coupling with INCA chemistry model, config_inca='//config_inca, modname)
241         CASE('inco'); CALL msg(TRIM(msg1)//' coupling jointly with INCA and CO2 cycle',  modname)
242         CASE('repr'); CALL msg(TRIM(msg1)//' coupling with REPROBUS chemistry model',    modname)
243         CASE('co2i'); CALL msg(TRIM(msg1)//' you have chosen to run with CO2 cycle',     modname)
244         CASE('coag'); CALL msg(TRIM(msg1)//' tracers are treated for COAGULATION tests', modname)
245         CASE('lmdz'); CALL msg(TRIM(msg1)//' tracers are treated in LMDZ only',          modname)
246         CASE DEFAULT; CALL abort_gcm(modname,'type_trac='//TRIM(str(it))//' not possible yet.',1)
247      END SELECT
248
249      !--- COHERENCE TEST BETWEEN "type_trac" AND "config_inca"
250      IF(ANY(['inca', 'inco'] == str(it)) .AND. ALL(['aero', 'aeNP', 'chem'] /= config_inca)) &
251         CALL abort_gcm(modname, 'Incoherence between type_trac and config_inca. Please modify "run.def"',1)
252
253      !--- COHERENCE TEST BETWEEN "type_trac" AND PREPROCESSING KEYS
254      SELECT CASE(str(it))
255         CASE('inca','inco')
256#ifndef INCA
257            CALL abort_gcm(modname, 'You must add cpp key INCA and compile with INCA code', 1)
258#endif
259         CASE('repr')
260#ifndef REPROBUS
261            CALL abort_gcm(modname, 'You must add cpp key REPROBUS and compile with REPROBUS code', 1)
262#endif
263         CASE('coag')
264#ifndef CPP_StratAer
265            CALL abort_gcm(modname, 'You must add cpp key StratAer and compile with StratAer code', 1)
266#endif
267      END SELECT
268
269      !--- DISABLE "config_inca" OPTION FOR A RUN WITHOUT "INCA" IF IT DIFFERS FROM "none"
270      IF(fmsg('Setting config_inca="none" as you do not couple with INCA model', &
271              modname, ALL(['inca', 'inco'] /= str(it)) .AND. config_inca /= 'none')) config_inca = 'none'
272
273   !---------------------------------------------------------------------------------------------------------------------------
274   END DO
275   !---------------------------------------------------------------------------------------------------------------------------
276
277   nqCO2 = 0; IF(ANY(str == 'inco')) nqCO2 = 1
278
279!==============================================================================================================================
280! 1) Get the numbers of: true (first order only) tracers "nqtrue", water tracers "nqo" (vapor/liquid/solid)
281!==============================================================================================================================
282
283!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
284   IF(lOldCode) THEN
285!------------------------------------------------------------------------------------------------------------------------------
286   !--- Determine nqtrue and (INCA only) nqo, nbtr
287   OPEN(90, FILE='traceur.def', FORM='formatted', STATUS='old', IOSTAT=ierr)
288   IF(ierr /= 0) CALL abort_gcm(modname, 'file "traceur.def" not found !', 1)
289   CALL msg('File "traceur.def" successfully opened.', modname)
290   lINCA = ANY(['inca','inco'] == type_trac)
291
292   IF(lINCA) THEN
293#ifdef INCA
294      READ(90,*) nqo
295      IF(nqCO2==1 .AND. nqo==4) nqo = 3
296      CALL Init_chem_inca_trac(nqINCA)
297      nbtr = nqINCA + nqCO2
298      nqtrue = nbtr + nqo
299      IF(ALL([2,3] /= nqo) CALL abort_gcm(modname, 'Only 2 or 3 water phases allowed ; found nqo='//TRIM(int2str(nqo)), 1)
300      CALL msg('nqo    = '//TRIM(int2str(nqo)),    modname)
301      CALL msg('nbtr   = '//TRIM(int2str(nbtr)),   modname)
302      CALL msg('nqtrue = '//TRIM(int2str(nqtrue)), modname)
303      CALL msg('nqCO2  = '//TRIM(int2str(nqCO2)),  modname)
304      CALL msg('nqINCA = '//TRIM(int2str(nqINCA)), modname)
305      ALLOCATE(hadv_inca(nqINCA), conv_flg_inca(nqINCA), solsym_inca(nqINCA))
306      ALLOCATE(vadv_inca(nqINCA),  pbl_flg_inca(nqINCA))
307      CALL init_transport(hadv_inca, vadv_inca, conv_flg_inca, pbl_flg_inca, solsym_inca)
308      ! DC passive CO2 tracer is at position 1: H2O was removed ; nqCO2/=0 in "inco" case only
309      conv_flg = [(  1,        ic=1, nqCO2),conv_flg_inca]
310       pbl_flg = [(  1,        ic=1, nqCO2), pbl_flg_inca]
311       solsym  = [('CO2     ', ic=1, nqCO2)   solsym_inca]
312      DEALLOCATE(conv_flg_inca, pbl_flg_inca, solsym_inca)
313#endif
314   ELSE
315      READ(90,*) nqtrue
316   END IF
317
318   IF (planet_type=="earth" .AND. nqtrue < 2) &
319      CALL abort_gcm('infotrac_init', 'Not enough tracers: nqtrue='//TRIM(int2str(nqtrue))//', 2 tracers is the minimum', 1)
320
321   !--- Allocate variables depending on nqtrue
322   ALLOCATE(hadv(nqtrue), vadv(nqtrue), tnom_0(nqtrue), tnom_transp(nqtrue), tracers(nqtrue))
323
324   !--- Continue to read tracer.def
325   it = 0
326   DO iq = 1, nqtrue
327#ifdef INCA
328      IF(iq > nqo+nqCO2) THEN
329         it = it+1
330         hadv  (iq) = hadv_inca  (it)
331         vadv  (iq) = vadv_inca  (it)
332         tnom_0(iq) = solsym_inca(it)
333         tnom_transp(iq) = 'air'
334         CYCLE
335      END IF
336#endif
337      CALL msg('237: iq='//TRIM(int2str(iq)), modname)
338      READ(90,'(I2,X,I2,X,A)',IOSTAT=ierr) hadv(iq),vadv(iq),tchaine
339      WRITE(msg1,'("iq, hadv(iq), vadv(iq) = ",i0,2(", ",i0))')iq, hadv(iq), vadv(iq)
340      CALL msg(msg1, modname)
341      CALL msg('tchaine = "'//TRIM(tchaine)//'"', modname)
342      CALL msg('infotrac 238: IOstatus='//TRIM(int2str(ierr)), modname)
343      IF(ierr/=0) CALL abort_gcm('infotrac_init', 'Pb dans la lecture de traceur.def', 1)
344      jq = INDEX(tchaine(1:LEN_TRIM(tchaine)),' ')
345      CALL msg("Ancienne version de traceur.def: traceurs d'air uniquement", modname, iq==1 .AND. jq==0)
346      CALL msg("Nouvelle version de traceur.def",                            modname, iq==1 .AND. jq/=0)
347      IF(jq /= 0) THEN                                               !--- Space in the string chain => new format
348         tnom_0     (iq) = tchaine(1:jq-1)
349         tnom_transp(iq) = tchaine(jq+1:)
350      ELSE
351         tnom_0     (iq) = tchaine
352         tnom_transp(iq) = 'air'
353      END IF
354      CALL msg(     'tnom_0(iq)=<'//TRIM(tnom_0(iq))     //'>', modname)
355      CALL msg('tnom_transp(iq)=<'//TRIM(tnom_transp(iq))//'>', modname)
356   END DO
357
358   CLOSE(90)
359
360#ifndef INCA
361   conv_flg = [(1, ic=1, nbtr)]                                      !--- Convection activated for all tracers
362    pbl_flg = [(1, ic=1, nbtr)]                                      !--- Boundary layer activated for all tracers
363   ALLOCATE(solsym(nbtr))
364   CALL msg('Valeur de traceur.def :', modname)
365   CALL msg('nombre total de traceurs '//TRIM(int2str(nqtrue)), modname)
366   DO iq = 1, nqtrue
367      CALL msg(strStack([int2str(hadv(iq)), int2str(vadv(iq)), tnom_0(iq), tnom_transp(iq)]), modname)
368   END DO
369   IF(planet_type /= 'earth') nqo = 0                                !--- Same number of tracers in dynamics and physics
370   IF(planet_type == 'earth') nqo = COUNT(delPhase(tnom_0) == 'H2O') !--- for all planets except for Earth
371   nbtr = nqtrue - nqo               
372#endif
373
374   !--- SET FIELDS %name, %parent, %phase, %component
375   tracers(:)%name      = tnom_0
376   tracers(:)%parent    = tnom_transp
377   tracers(:)%phase     = 'g'
378   tracers(:)%component = type_trac
379   DO ip = 1, nphases
380      p = old_phases(ip:ip)
381      iq = strIdx(tracers(:)%name, 'H2O'//p)
382      IF(iq /= 0) CYCLE
383      tracers(iq)%phase = p
384      IF(lINCA) tracers(iq)%component = 'lmdz'
385   END DO
386   IF(lINCA) tracers(1+nqo:nqCO2+nqo)%component = 'co2i'
387   CALL setGeneration(tracers)                                       !--- SET FIELDS %iGeneration, %gen0Name
388
389! manque "type"
390
391!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
392ELSE
393!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
394   IF(readTracersFiles(type_trac, fType, tracers)) CALL abort_gcm(modname, 'problem with tracers file(s)',1)
395   IF(fType == 0) CALL abort_gcm(modname, 'Missing tracers file: "traceur.def", "tracer.def" or "tracer_<keyword>.def file.',1)
396   !---------------------------------------------------------------------------------------------------------------------------
397   IF(fType == 1) THEN                                               !=== FOUND AN OLD STYLE "traceur.def"
398   !---------------------------------------------------------------------------------------------------------------------------
399#ifdef INCA
400      nqo = SIZE(tracers)
401      IF(nqCO2==1 .AND. nqo==4) nqo = 3                              !--- Force nqo to 3 (ThL)
402      CALL Init_chem_inca_trac(nqINCA)                               !--- Get nqINCA from INCA
403      nbtr = nqINCA + nqCO2                                          !--- Number of tracers passed to phytrac
404      nqtrue = nbtr + nqo                                            !--- Total number of "true" tracers
405      IF(ALL([2,3] /= nqo) CALL abort_gcm(modname, 'Only 2 or 3 water phases allowed ; found nqo='//TRIM(int2str(nqo)), 1)
406      CALL msg('nqo    = '//TRIM(int2str(nqo)),    modname)
407      CALL msg('nbtr   = '//TRIM(int2str(nbtr)),   modname)
408      CALL msg('nqtrue = '//TRIM(int2str(nqtrue)), modname)
409      CALL msg('nqCO2  = '//TRIM(int2str(nqCO2)),  modname)
410      CALL msg('nqINCA = '//TRIM(int2str(nqINCA)), modname)
411      ALLOCATE(hadv(nqtrue), hadv_inca(nqINCA), conv_flg_inca(nqINCA), solsym_inca(nqINCA))
412      ALLOCATE(vadv(nqtrue), vadv_inca(nqINCA),  pbl_flg_inca(nqINCA))
413      CALL init_transport(hadv_inca, vadv_inca, conv_flg_inca, pbl_flg_inca, solsym_inca)
414      ! DC passive CO2 tracer is at position 1: H2O was removed ; nqCO2/=0 in "inco" case only
415      conv_flg = [(  1,        k=1, nqCO2), conv_flg_inca]
416       pbl_flg = [(  1,        k=1, nqCO2),  pbl_flg_inca]
417       solsym  = [('CO2     ', k=1, nqCO2)    solsym_inca]
418      DEALLOCATE(conv_flg_inca, pbl_flg_inca, solsym_inca)
419      ALLOCATE(ttr(nqtrue))
420      ttr(1:nqo+nqCO2)                    = tracers
421      ttr(1    :      nqo   )%component   = 'lmdz'
422      ttr(1+nqo:nqCO2+nqo   )%component   = 'co2i'
423      ttr(1+nqo+nqCO2:nqtrue)%component   = 'inca'
424      ttr(1+nqo      :nqtrue)%name        = solsym
425      ttr(1+nqo+nqCO2:nqtrue)%parent      = tran0
426      ttr(1+nqo+nqCO2:nqtrue)%phase       = 'g'
427      lerr = getKey('hadv', had, ky=ttr(:)%keys); hadv(:) = [had, hadv_inca]
428      lerr = getKey('vadv', vad, ky=ttr(:)%keys); vadv(:) = [vad, vadv_inca]
429      DEALLOCATE(had, hadv_inca, vad, vadv_inca)
430      CALL MOVE_ALLOC(FROM=ttr, TO=tracers)
431      CALL setGeneration(tracers)                                    !--- SET FIELDS %iGeneration, %gen0Name
432#else
433      nqo    = COUNT(delPhase(tracers(:)%name) == 'H2O')             !--- Number of water phases
434      nqtrue = SIZE(tracers)                                         !--- Total number of "true" tracers
435      nbtr   = nqtrue - nqo                                          !--- Number of tracers passed to phytrac
436      lerr = getKey('hadv', hadv, ky=tracers(:)%keys)
437      lerr = getKey('vadv', vadv, ky=tracers(:)%keys)
438#endif
439   !---------------------------------------------------------------------------------------------------------------------------
440   ELSE                                                              !=== FOUND NEW STYLE TRACERS CONFIGURATION FILE(S)
441   !---------------------------------------------------------------------------------------------------------------------------
442      nqo    = COUNT(delPhase(tracers(:)%name) == 'H2O')             !--- Number of water phases
443      nqtrue = SIZE(tracers)                                         !--- Total number of "true" tracers
444      nbtr   = nqtrue - nqo                                          !--- Number of tracers passed to phytrac
445      lerr = getKey('hadv', hadv, ky=tracers(:)%keys)
446      lerr = getKey('vadv', vadv, ky=tracers(:)%keys)
447      ALLOCATE(solsym(nbtr))
448      conv_flg = [(1, it=1, nbtr)]
449       pbl_flg = [(1, it=1, nbtr)]
450   !---------------------------------------------------------------------------------------------------------------------------
451   END IF
452   !---------------------------------------------------------------------------------------------------------------------------
453END IF
454!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
455
456   CALL getKey_init(tracers)
457
458   !--- Transfert the number of tracers to Reprobus
459#ifdef REPROBUS
460   CALL Init_chem_rep_trac(nbtr, nqo, tracers(:)%name)
461#endif
462
463!==============================================================================================================================
464! 2) Calculate nqtot, number of tracers needed (greater if advection schemes 20 or 30 have been chosen).
465!==============================================================================================================================
466   DO iq = 1, nqtrue
467      IF( hadv(iq)<20 .OR. (ANY(hadv(iq)==[20,30]) .AND. hadv(iq)==vadv(iq)) ) CYCLE
468      WRITE(msg1,'("The choice hadv=",i0,", vadv=",i0,a)')hadv(iq),vadv(iq),' for "'//TRIM(tracers(iq)%name)//'" is not available'
469      CALL abort_gcm(modname, TRIM(msg1), 1)
470   END DO
471   nqtot =    COUNT( hadv< 20 .AND. vadv< 20 ) &                     !--- No additional tracer
472         +  4*COUNT( hadv==20 .AND. vadv==20 ) &                     !--- 3  additional tracers
473         + 10*COUNT( hadv==30 .AND. vadv==30 )                       !--- 9  additional tracers
474
475   !--- More tracers due to the choice of advection scheme => assign total number of tracers
476   IF( nqtot /= nqtrue ) THEN
477      CALL msg('The choice of advection scheme for one or more tracers makes it necessary to add tracers')
478      CALL msg('The number of true tracers is '//TRIM(int2str(nqtrue)))
479      CALL msg('The total number of tracers needed is '//TRIM(int2str(nqtot)))
480   END IF
481   CALL msg('nqo    = '//TRIM(int2str(nqo)),    modname)
482   CALL msg('nbtr   = '//TRIM(int2str(nbtr)),   modname)
483   CALL msg('nqtrue = '//TRIM(int2str(nqtrue)), modname)
484   CALL msg('nqtot  = '//TRIM(int2str(nqtot)),  modname)
485
486!==============================================================================================================================
487! 3) Determine the advection scheme choice for water and tracers "iadv" and the fields long name, isAdvected.
488!     iadv = 1    "LMDZ-specific humidity transport" (for H2O vapour)          LMV
489!     iadv = 2    backward                           (for H2O liquid)          BAK
490!     iadv = 14   Van-Leer + specific humidity, modified by Francis Codron     VLH
491!     iadv = 10   Van-Leer (chosen for vapour and liquid water)                VL1
492!     iadv = 11   Van-Leer for hadv and PPM version (Monotonic) for vadv       VLP
493!     iadv = 12   Frederic Hourdin I                                           FH1
494!     iadv = 13   Frederic Hourdin II                                          FH2
495!     iadv = 16   Monotonic         PPM (Collela & Woodward 1984)              PPM
496!     iadv = 17   Semi-monotonic    PPM (overshoots allowed)                   PPS
497!     iadv = 18   Definite positive PPM (overshoots and undershoots allowed)   PPP
498!     iadv = 20   Slopes                                                       SLP
499!     iadv = 30   Prather                                                      PRA
500!
501!        In array q(ij,l,iq) : iq = 1/2[/3]    for vapour/liquid[/ice] water
502!        And optionaly:        iq = 3[4],nqtot for other tracers
503!==============================================================================================================================
504   ALLOCATE(ttr(nqtot))
505   jq = nqtrue+1; tracers(:)%iadv = -1
506   DO iq = 1, nqtrue
507      t1 => tracers(iq)
508
509      !--- VERIFY THE CHOICE OF ADVECTION SCHEME
510      iad = -1
511      IF(hadv(iq)     ==    vadv(iq)    ) iad = hadv(iq)
512      IF(hadv(iq)==10 .AND. vadv(iq)==16) iad = 11
513      WRITE(msg1,'("Bad choice of advection scheme for ",a,": hadv = ",i0,", vadv = ",i0)')TRIM(t1%name), hadv(iq), vadv(iq)
514      IF(iad == -1) CALL abort_gcm(modname, msg1, 1)
515
516      !--- SET FIELDS %longName, %iadv, %isAdvected, %isH2Ofamily
517      t1%longName   = t1%name; IF(iad > 0) t1%longName=TRIM(t1%name)//descrq(iad)
518      t1%iadv       = iad
519      t1%isAdvected = iad > 0
520      t1%isH2Ofamily= delPhase(t1%gen0Name) == 'H2O'
521      ttr(iq)       = t1
522
523      !--- DEFINE THE HIGHER ORDER TRACERS, IF ANY
524      nm = 0
525      IF(iad == 20) nm = 3                                             !--- 2nd order scheme
526      IF(iad == 30) nm = 9                                             !--- 3rd order scheme
527      IF(nm == 0) CYCLE                                                !--- No higher moments
528      ttr(jq+1:jq+nm)             = t1
529      ttr(jq+1:jq+nm)%name        = [(TRIM(t1%name)    //'-'//TRIM(suff(im)), im=1, nm) ]
530      ttr(jq+1:jq+nm)%parent      = [(TRIM(t1%parent)  //'-'//TRIM(suff(im)), im=1, nm) ]
531      ttr(jq+1:jq+nm)%longName    = [(TRIM(t1%longName)//'-'//TRIM(suff(im)), im=1, nm) ]
532      ttr(jq+1:jq+nm)%iadv        = [(-iad,    im=1, nm) ]
533      ttr(jq+1:jq+nm)%isAdvected  = [(.FALSE., im=1, nm) ]
534      jq = jq + nm
535   END DO
536   DEALLOCATE(hadv, vadv)
537   CALL MOVE_ALLOC(FROM=ttr, TO=tracers)
538
539   !--- SET FIELDS %iqParent, %nqChilds, %iGeneration, %iqDescen, %nqDescen
540   CALL indexUpdate(tracers)
541
542   CALL msg('Information stored in infotrac :', modname)
543   CALL msg('iadv  name  long_name :', modname)
544
545   !=== TEST ADVECTION SCHEME
546   DO iq=1,nqtot ; t1 => tracers(iq); iad = t1%iadv
547
548      !--- ONLY TESTED VALUES FOR TRACERS FOR NOW:               iadv = 14, 10 (and 0 for non-transported tracers)
549      IF(ALL([10,14,0] /= iad)) &
550         CALL abort_gcm(modname, 'Not tested for iadv='//TRIM(int2str(iad))//' ; 10 or 14 only are allowed !', 1)
551
552      !--- ONLY TESTED VALUES FOR PARENTS HAVING CHILDS FOR NOW: iadv = 14, 10 (PARENTS: TRACERS OF GENERATION 1)
553      IF(ALL([10,14] /= iad) .AND. t1%iGeneration == 1 .AND. ANY(tracers(:)%iGeneration > 1)) &
554         CALL abort_gcm(modname, 'iadv='//TRIM(int2str(iad))//' not implemented for parents ; 10 or 14 only are allowed !', 1)
555
556      !--- ONLY TESTED VALUES FOR CHILDS FOR NOW:                iadv = 10     (CHILDS:  TRACERS OF GENERATION GREATER THAN 1)
557      IF(fmsg('WARNING ! iadv='//TRIM(int2str(iad))//' not implemented for childs. Setting iadv=10 for "'//TRIM(t1%name)//'"',&
558         modname, iad /= 10 .AND. t1%iGeneration > 1)) t1%iadv = 10
559
560      !--- ONLY VALID SCHEME NUMBER FOR WATER VAPOUR:            iadv = 14
561      ll = t1%name /= addPhase('H2O','g'); IF(lOldCode) ll = t1%name /= 'H2Ov'
562      IF(fmsg('WARNING ! iadv=14 is valid for water vapour only. Setting iadv=10 for "'//TRIM(t1%name)//'".', &
563         modname, iad == 14 .AND. ll))                 t1%iadv = 10
564   END DO
565
566IF(lOldCode) THEN
567
568   CALL infotrac_setHeredity                !--- SET FIELDS %iqParent, %nqChilds, %iGeneration, %gen0Name, %iqDescen, %nqDescen
569   CALL infotrac_isoinit                    !--- SET FIELDS %type, %iso_iName, %iso_iZone, %iso_iPhase
570   CALL getKey_init(tracers, isotopes)
571   IF(isoSelect('H2O')) RETURN                                    !--- Select water isotopes ; finished if no water isotopes
572   iH2O = ixIso                                                   !--- Keep track of water family index
573
574   !--- Remove the isotopic tracers from the tracers list passed to phytrac
575   nbtr    = nbtr -nqo*   ntiso             !--- ISOTOPIC TAGGING TRACERS ARE NOT PASSED TO THE PHYSICS
576   nqtottr = nqtot-nqo*(1+ntiso)            !--- NO H2O-FAMILY    TRACER  IS      PASSED TO THE PHYSICS
577   CALL msg('702: nbtr, ntiso='//strStack(int2str([nbtr, ntiso])), modname)
578   CALL msg('704: nqtottr, nqtot, nqo = '//strStack(int2str([nqtottr, nqtot, nqo])), modname)
579   ! Rq: nqtottr n'est pas forcement egal a nbtr dans le cas ou nmom/=0
580   IF(COUNT(tracers%iso_iName == 0 .AND. delPhase(tracers(:)%name)/='H2O') /= nqtottr) &
581      CALL abort_gcm('infotrac_init', 'pb dans le calcul de nqtottr', 1)
582
583   !--- Finalize :
584   DEALLOCATE(tnom_0, tnom_transp)
585
586ELSE
587
588   CALL initIsotopes(tracers, isotopes)
589   nbIso = SIZE(isotopes); IF(nbIso==0) RETURN                    !--- No isotopes: finished.
590
591   !=== READ PHYSICAL PARAMETERS FROM "isotopes_params.def" FILE SPECIFIC TO WATER ISOTOPES
592   !    DONE HERE, AND NOT ONLY IN "infotrac_phy", BECAUSE SOME PHYSICAL PARAMS ARE NEEDED FOR RESTARTS (tnat, alpha_ideal)
593   CALL getKey_init(tracers, isotopes)
594   IF(isoSelect('H2O')) RETURN                                    !--- Select water isotopes ; finished if no water isotopes
595   iH2O = ixIso                                                   !--- Keep track of water family index
596   IF(getKey('tnat' , tnat,        isoName(1:niso))) CALL abort_gcm(modname, 'can''t read "tnat"', 1)
597   IF(getKey('alpha', alpha_ideal, isoName(1:niso))) CALL abort_gcm(modname, 'can''t read "alpha_ideal"', 1)
598
599   !=== ENSURE THE MEMBERS OF AN ISOTOPES FAMILY ARE PRESENT IN THE SAME PHASES
600   DO ix = 1, nbIso
601      iso => isotopes(ix)
602      !--- Check whether each isotope and tagging isotopic tracer is present in the same number of phases
603      DO it = 1, iso%ntiso
604         np = SUM([(COUNT(tracers(:)%name == addPhase(iso%trac(it), iso%phase(ip:ip))), ip=1, iso%nphas)])
605         IF(np == iso%nphas) CYCLE
606         WRITE(msg1,'("Found ",i0," phases for ",s," instead of ",i0)')np, iso%trac(it), iso%nphas
607         CALL abort_gcm(modname, msg1, 1)
608      END DO
609      DO it = 1, iso%niso
610         nz = SUM([(COUNT(iso%trac == iso%trac(it)//'_'//iso%zone(iz)), iz=1, iso%nzone)])
611         IF(nz == iso%nzone) CYCLE
612         WRITE(msg1,'("Found ",i0," tagging zones for ",s," instead of ",i0)')nz, iso%trac(it), iso%nzone
613         CALL abort_gcm(modname, msg1, 1)
614      END DO
615   END DO
616   nqtottr = COUNT(tracers%iso_iName == 0)
617
618END IF
619
620   !=== DISPLAY THE RESULTING LIST
621   t => tracers
622   CALL msg('Information stored in infotrac :')
623   IF(dispTable('isssssssssiiiiiiiii', &
624      ['iq      ', 'name    ', 'longN.  ', 'gen0N.  ', 'parent  ', 'type    ', 'phase   ', 'compon. ', 'isAdv.  ', 'isH2O.  '&
625      ,'iadv    ', 'iGen.   ', 'iqPar.  ', 'nqDes.  ', 'nqChil. ', 'iso_iG. ', 'iso_iN. ', 'iso_iZ. ', 'iso_iP. '],          &
626      cat(t%name,  t%longName,  t%gen0Name,  t%parent,  t%type,  t%phase, &
627          t%component, bool2str(t%isAdvected), bool2str(t%isH2Ofamily)),  &
628      cat([(iq, iq=1, nqtot)],  t%iadv,  t%iGeneration, t%iqParent, t%nqDescen, &
629         t%nqChilds, t%iso_iGroup, t%iso_iName, t%iso_iZone, t%iso_iPhase))) &
630      CALL abort_gcm(modname, "problem with the tracers table content", 1)
631
632   !--- Some aliases to be removed later
633   ntraciso       => isotope%ntiso
634   ntraceurs_zone => isotope%nzone
635   qperemin       =  min_qParent
636   masseqmin      =  min_qMass
637   ratiomin       =  min_ratio
638   CALL msg('end', modname)
639
640END SUBROUTINE infotrac_init
641
642
643
644SUBROUTINE infotrac_setHeredity
645   !--- Purpose: Set fields %iqParent, %nqChilds, %iGeneration, %iqDescen, %nqDescen (old method)
646   USE strings_mod, ONLY: strIdx
647   INTEGER               :: iq, ipere, ifils
648   INTEGER, ALLOCATABLE  :: iqfils(:,:)
649   CHARACTER(LEN=maxlen) :: msg1, modname='infotrac_init'
650   INCLUDE "iniprint.h"
651
652   !=== SET FIELDS %iqParent, %nqChilds
653   ALLOCATE(iqfils(nqtot,nqtot)); iqfils(:,:) = 0
654
655   DO iq = 1, nqtot
656      msg1 = 'Tracer nr. '//TRIM(int2str(iq))//', called "'//TRIM(tracers(iq)%name)//'" is '
657
658      !--- IS IT A GENERATION 0 TRACER ? IF SO, tracers(iq)%iqParent KEEPS ITS DEFAULT VALUE (0)
659      IF(fmsg(TRIM(msg1)//' a parent', modname, tracers(iq)%parent == tran0)) CYCLE
660
661      !--- TRACERS OF GENERATION > 0 ; FIND ITS PARENT INDEX
662      ipere = strIdx(tracers(:)%name, tracers(iq)%parent)
663      IF(ipere == 0)  CALL abort_gcm('infotrac_init', TRIM(msg1)//' an orphan', 1)
664      IF(iq == ipere) CALL abort_gcm('infotrac_init', TRIM(msg1)//' its own parent',1)
665
666      CALL msg(TRIM(msg1)//' the child of '//TRIM(tracers(ipere)%name), modname)
667      tracers(iq)%iqParent    = ipere
668      tracers(ipere)%nqChilds = tracers(ipere)%nqChilds+1
669      iqfils(tracers(ipere)%nqChilds,ipere) = iq
670   END DO
671   CALL msg('nqGen0   = '//int2str(COUNT(tracers(:)%parent == 'air')), modname)
672   CALL msg('nqChilds = '//strStack(int2str(tracers(:)%nqChilds)),     modname)
673   CALL msg('iqParent = '//strStack(int2str(tracers(:)%iqParent)),     modname)
674   CALL msg('iqChilds = '//strStack(int2str(PACK(iqfils,MASK=.TRUE.))),modname)
675
676   !=== SET FIELDS %iGeneration, %iqDescen, %nqDescen
677   tracers(:)%iGeneration = 0
678   DO iq = 1, nqtot
679      ifils = iq
680      DO WHILE(tracers(ifils)%iqParent > 0)
681         ipere = tracers(ifils)%iqParent
682         tracers(ipere)%nqDescen = tracers(ipere)%nqDescen+1   
683         tracers(iq)%iGeneration = tracers(iq)%iGeneration+1
684         iqfils(tracers(ipere)%nqDescen,ipere) = iq
685         ifils = ipere
686      END DO
687      msg1 = 'Tracer nr. '//TRIM(int2str(iq))//', called "'//TRIM(tracers(iq)%name)//'" is '
688      CALL msg(TRIM(msg1)//' of generation '//TRIM(int2str(tracers(iq)%iGeneration)), modname)
689   END DO
690   DO iq=1,nqtot
691      tracers(iq)%iqDescen = iqfils(1:tracers(iq)%nqDescen,iq)
692   END DO
693
694   CALL msg('nqDescen = '//TRIM(strStack(int2str(tracers(:)%nqDescen))), modname)
695   CALL msg('nqDescen_tot = ' //TRIM(int2str(SUM(tracers(:)%nqDescen))), modname)
696   CALL msg('iqChilds = '//strStack(int2str(PACK(iqfils, MASK=.TRUE.))), modname)
697
698END SUBROUTINE infotrac_setHeredity
699
700
701
702SUBROUTINE infotrac_isoinit
703
704#ifdef CPP_IOIPSL
705   USE IOIPSL
706#else
707   USE ioipsl_getincom
708#endif
709   IMPLICIT NONE
710   CHARACTER(LEN=3)      :: tnom_iso(niso_possibles)
711   INTEGER, ALLOCATABLE  :: nb_iso(:,:), nb_traciso(:,:)
712   INTEGER               :: ii, ip, iq, it, iz, ixt, nb_isoind, nzone_prec
713   TYPE(isot_type), POINTER :: i
714   TYPE(trac_type), POINTER :: t(:)
715   CHARACTER(LEN=maxlen)    :: tnom_trac
716   CHARACTER(LEN=maxlen), ALLOCATABLE :: str(:)
717   INCLUDE "iniprint.h"
718
719   tnom_iso = ['eau', 'HDO', 'O18', 'O17', 'HTO']
720   ALLOCATE(nb_iso       (niso_possibles,nqo))
721   ALLOCATE(nb_traciso   (niso_possibles,nqo))
722   ALLOCATE(use_iso      (niso_possibles))
723   ALLOCATE(indnum_fn_num(niso_possibles))
724   ALLOCATE(iso_indnum(nqtot))
725
726   iso_indnum   (:) = 0
727   use_iso      (:) = .FALSE.
728   indnum_fn_num(:) = 0
729   nb_iso     (:,:) = 0 
730   nb_traciso (:,:) = 0
731
732   DO iq=1, nqtot
733      IF(delPhase(tracers(iq)%name) == 'H2O' .OR. .NOT.tracers(iq)%isAdvected) CYCLE
734outer:DO ip = 1, nqo
735         nb_isoind = 0
736         DO ixt= 1,niso_possibles
737            tnom_trac = 'H2O'//old_phases(ip:ip)//'_'//TRIM(tnom_iso(ixt))
738            IF (tracers(iq)%name == tnom_trac) THEN
739               nb_iso(ixt,ip)         = nb_iso(ixt,ip)+1
740               nb_isoind              = nb_isoind+1
741               tracers(iq)%type       = 'tracer'
742               tracers(iq)%iso_iGroup = 1
743               tracers(iq)%iso_iName  = ixt
744               iso_indnum(iq)         = nb_isoind
745               indnum_fn_num(ixt)     = iso_indnum(iq)
746               tracers(iq)%iso_iPhase = ip
747               EXIT outer
748            ELSE IF(tracers(iq)%iqParent> 0) THEN
749               IF(tracers(tracers(iq)%iqParent)%name == tnom_trac) THEN
750                  nb_traciso(ixt,ip)  = nb_traciso(ixt,ip)+1
751                  iso_indnum(iq)         = indnum_fn_num(ixt)
752                  tracers(iq)%type    = 'tag'
753                  tracers(iq)%iso_iGroup = 1
754                  tracers(iq)%iso_iName  = ixt
755                  tracers(iq)%iso_iZone  = nb_traciso(ixt,ip)
756                  tracers(iq)%iso_iPhase = ip
757                  EXIT outer
758               END IF
759            END IF
760         END DO
761      END DO outer
762   END DO
763
764   niso = 0; nzone_prec = nb_traciso(1,1)
765   DO ixt = 1, niso_possibles
766      IF(nb_iso(ixt,1) == 0) CYCLE
767      IF(nb_iso(ixt,1) /= 1) CALL abort_gcm('infotrac_init', 'Isotopes are not well defined in traceur.def', 1)
768
769      ! on verifie que toutes les phases ont le meme nombre d'isotopes
770      IF(ANY(nb_iso(ixt,:) /= 1)) CALL abort_gcm('infotrac_init', 'Phases must have same number of isotopes', 1)
771
772      niso = niso+1
773      use_iso(ixt) = .TRUE.
774      nzone = nb_traciso(ixt,1)
775
776      ! on verifie que toutes les phases ont le meme nombre de traceurs d'isotopes
777      IF(ANY(nb_traciso(ixt,2:nqo) /= nzone)) CALL abort_gcm('infotrac_init','Phases must have same number of tracers',1)
778
779      ! on verifie que tous les isotopes ont le meme nombre de traceurs d'isotopes
780      IF(nzone /= nzone_prec) CALL abort_gcm('infotrac_init','Isotope tracers are not well defined in traceur.def',1)
781      nzone_prec = nzone
782   END DO
783
784   ! dimensions et flags isotopiques:
785   ntiso = niso*(nzone+1)
786   ok_isotopes = niso  > 0
787   ok_isotrac  = nzone > 0
788 
789   IF(ok_isotopes) THEN
790      ok_iso_verif = .FALSE.; CALL getin('ok_iso_verif', ok_iso_verif)
791      ok_init_iso  = .FALSE.; CALL getin('ok_init_iso',  ok_init_iso)
792   END IF
793      tnat        = [1.0, 155.76e-6, 2005.2e-6, 0.004/100., 0.0]
794      alpha_ideal = [1.0, 1.01,      1.006,     1.003,      1.0]
795!   END IF
796
797   ! remplissage du tableau iqiso(ntiso,phase)
798   ALLOCATE(iqiso(ntiso,nqo))   
799   iqiso(:,:)=0     
800   DO iq = 1, nqtot
801      IF(tracers(iq)%iso_iName <= 0) CYCLE
802      ixt = iso_indnum(iq) + tracers(iq)%iso_iZone*niso
803      iqiso(ixt, tracers(iq)%iso_iPhase) = iq
804   END DO
805
806   ! remplissage du tableau index_trac(nzone,niso)
807   ALLOCATE(index_trac(nzone, niso)) 
808   IF(ok_isotrac) then
809      DO ii = 1, niso; index_trac(:, ii) = ii + niso*[(iz, iz=1, nzone)]; END DO
810   ELSE
811      index_trac(:,:)=0.0
812   END IF
813
814   ALLOCATE(isotopes(1))                                             !--- Only water
815   nbIso = 1
816   i => isotopes(1)
817   t => tracers
818   str = PACK(delPhase(t%name), MASK = t%type=='tracer' .AND. delPhase(t%parent) == 'H2O' .AND. t%phase == 'g')
819   i%parent = 'H2O'
820   i%niso  = SIZE(str)
821   i%nzone = nzone
822   i%nphas = nqo
823   FORALL(it = 1:i%niso) i%keys(it)%name = str(it)
824   i%zone = PACK(strTail(t%name,'_',.TRUE.), MASK = t%type=='tag' .AND. delPhase(t%gen0Name)=='H2O' .AND. t%iGeneration==3)
825   CALL strReduce(i%zone)
826   i%phase = strStack([(known_phases(ip:ip), ip=1, nphases)], MASK=[(strIdx(t%name,addPhase('H2O',known_phases(ip:ip)))/=0)])
827   i%iTraPha = RESHAPE([((strIdx(t(:)%name,addPhase(i%trac(it),i%phase(ip:ip))) ,it=1,i%ntiso),ip=1,i%nphas)],[i%ntiso,i%nphas])
828   i%iZonIso = RESHAPE([((strIdx(i%trac,TRIM(i%trac(it))//'_'//TRIM(i%zone(iz))),iz=1,i%nzone),it=1,i%niso )],[i%nzone,i%niso ])
829   DO it=1,ntiso
830      WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': iqiso  (',it,',:) = '//strStack(int2str(iqiso(it,:)))
831      WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': iTraPha(',it,',:) = '//strStack(int2str(i%iTraPha(it,:)))
832   END DO
833   DO iz=1,nzone
834      WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': index_trac(',iz,',:) = '//strStack(int2str(index_trac(iz,:)))
835      WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': iZonIso   (',iz,',:) = '//strStack(int2str(i%iZonIso(iz,:)))
836   END DO
837
838   ! Finalize :
839   DEALLOCATE(nb_iso)
840
841END SUBROUTINE infotrac_isoinit
842
843
844!==============================================================================================================================
845!=== THE ROUTINE isoSelect IS USED TO SWITCH FROM AN ISOTOPE FAMILY TO ANOTHER: ISOTOPES DEPENDENT PARAMETERS ARE UPDATED
846!     Single generic "isoSelect" routine, using the predefined index of the parent (fast version) or its name (first call).
847!==============================================================================================================================
848LOGICAL FUNCTION isoSelectByName(iName, lVerbose) RESULT(lerr)
849   IMPLICIT NONE
850   CHARACTER(LEN=*),  INTENT(IN)  :: iName
851   LOGICAL, OPTIONAL, INTENT(IN) :: lVerbose
852   INTEGER :: iIso
853   LOGICAL :: lV
854   lV = .FALSE.; IF(PRESENT(lVerbose)) lV = lVerbose
855   iIso = strIdx(isotopes(:)%parent, iName)
856   lerr = iIso == 0
857   CALL msg('no isotope family named "'//TRIM(iName)//'"', ll=lerr .AND. lV)
858   IF(lerr) RETURN
859   lerr = isoSelectByIndex(iIso, lV)
860END FUNCTION isoSelectByName
861!==============================================================================================================================
862LOGICAL FUNCTION isoSelectByIndex(iIso, lVerbose) RESULT(lerr)
863   IMPLICIT NONE
864   INTEGER,           INTENT(IN) :: iIso
865   LOGICAL, OPTIONAL, INTENT(IN) :: lVerbose
866   LOGICAL :: lv
867   lv = .FALSE.; IF(PRESENT(lVerbose)) lv = lVerbose
868   lerr = .FALSE.
869   IF(iIso == ixIso) RETURN                                      !--- Nothing to do if the index is already OK
870   lerr = iIso<=0 .OR. iIso>nbIso
871   CALL msg('Inconsistent isotopes family index '//TRIM(int2str(iIso))//': should be > 0 and <= '//TRIM(int2str(nbIso))//'"',&
872            ll=lerr .AND. lV)
873   IF(lerr) RETURN
874   ixIso = iIso                                                  !--- Update currently selected family index
875   isotope => isotopes(ixIso)                                    !--- Select corresponding component
876   isoKeys => isotope%keys;    niso     = isotope%niso
877   isoName => isotope%trac;    ntiso    = isotope%ntiso
878   isoZone => isotope%zone;    nzone    = isotope%nzone
879   isoPhas => isotope%phase;   nphas    = isotope%nphas
880   iZonIso => isotope%iZonIso; isoCheck = isotope%check
881   iTraPha => isotope%iTraPha
882END FUNCTION isoSelectByIndex
883!==============================================================================================================================
884
885END MODULE infotrac
Note: See TracBrowser for help on using the repository browser.