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

Last change on this file since 4127 was 4127, checked in by acozic, 2 years ago

modify lmdz code in order to transfer to inca model all infomations necessary for the coupling with dynamico

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