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

Last change on this file since 4071 was 4071, checked in by dcugnet, 2 years ago
  • Fix for unadvected tracers (iadv==0)
  • The key %isH2Ofamily, from the derived type "trac_type", is replaced with the more general

key %isInPhysics, which is TRUE for tracers both in "qx" and "tr_seri".

Currently, FALSE for tracers descending on H2O (isotopes and tagging tracers included). Could be set to FALSE
for interactive CO2 (type_trac=='inco') or ice supersaturated cloud content (tranfered to "rneb_seri")

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