1 | !$Id: infotrac.F90 4082 2022-02-28 08:37:41Z fhourdin $ |
---|
2 | ! |
---|
3 | MODULE 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 | |
---|
160 | CONTAINS |
---|
161 | |
---|
162 | SUBROUTINE 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, k !--- 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 | ALLOCATE(conv_flg(nbtr),pbl_flg(nbtr),solsym(nbtr)) |
---|
311 | conv_flg = [( 1, ic=1, nqCO2),conv_flg_inca] |
---|
312 | pbl_flg = [( 1, ic=1, nqCO2), pbl_flg_inca] |
---|
313 | solsym = [('CO2 ', ic=1, nqCO2), solsym_inca] |
---|
314 | DEALLOCATE(conv_flg_inca, pbl_flg_inca) |
---|
315 | #endif |
---|
316 | ELSE |
---|
317 | READ(90,*) nqtrue |
---|
318 | END IF |
---|
319 | |
---|
320 | IF (planet_type=="earth" .AND. nqtrue < 2) & |
---|
321 | CALL abort_gcm('infotrac_init', 'Not enough tracers: nqtrue='//TRIM(int2str(nqtrue))//', 2 tracers is the minimum', 1) |
---|
322 | |
---|
323 | !--- Allocate variables depending on nqtrue |
---|
324 | ALLOCATE(hadv(nqtrue), vadv(nqtrue), tnom_0(nqtrue), tnom_transp(nqtrue), tracers(nqtrue)) |
---|
325 | |
---|
326 | !--- Continue to read tracer.def |
---|
327 | it = 0 |
---|
328 | DO iq = 1, nqtrue |
---|
329 | #ifdef INCA |
---|
330 | IF(iq > nqo+nqCO2) THEN |
---|
331 | it = it+1 |
---|
332 | hadv (iq) = hadv_inca (it) |
---|
333 | vadv (iq) = vadv_inca (it) |
---|
334 | tnom_0(iq) = solsym_inca(it) |
---|
335 | tnom_transp(iq) = 'air' |
---|
336 | CYCLE |
---|
337 | END IF |
---|
338 | #endif |
---|
339 | CALL msg('237: iq='//TRIM(int2str(iq)), modname) |
---|
340 | READ(90,'(I2,X,I2,X,A)',IOSTAT=ierr) hadv(iq),vadv(iq),tchaine |
---|
341 | WRITE(msg1,'("hadv(",i0,"), vadv(",i0,") = ",i0,", ",i0)')iq, iq, hadv(iq), vadv(iq) |
---|
342 | CALL msg(TRIM(msg1), modname) |
---|
343 | CALL msg('tchaine = "'//TRIM(tchaine)//'"', modname) |
---|
344 | CALL msg('infotrac 238: IOstatus='//TRIM(int2str(ierr)), modname) |
---|
345 | IF(ierr/=0) CALL abort_gcm('infotrac_init', 'Pb dans la lecture de traceur.def', 1) |
---|
346 | jq = INDEX(tchaine(1:LEN_TRIM(tchaine)),' ') |
---|
347 | CALL msg("Ancienne version de traceur.def: traceurs d'air uniquement", modname, iq==1 .AND. jq==0) |
---|
348 | CALL msg("Nouvelle version de traceur.def", modname, iq==1 .AND. jq/=0) |
---|
349 | IF(jq /= 0) THEN !--- Space in the string chain => new format |
---|
350 | tnom_0 (iq) = tchaine(1:jq-1) |
---|
351 | tnom_transp(iq) = tchaine(jq+1:) |
---|
352 | ELSE |
---|
353 | tnom_0 (iq) = tchaine |
---|
354 | tnom_transp(iq) = 'air' |
---|
355 | END IF |
---|
356 | CALL msg( 'tnom_0(iq)=<'//TRIM(tnom_0(iq)) //'>', modname) |
---|
357 | CALL msg('tnom_transp(iq)=<'//TRIM(tnom_transp(iq))//'>', modname) |
---|
358 | END DO |
---|
359 | #ifdef INCA |
---|
360 | DEALLOCATE(solsym_inca) |
---|
361 | #endif |
---|
362 | |
---|
363 | CLOSE(90) |
---|
364 | |
---|
365 | #ifndef INCA |
---|
366 | CALL msg('Valeur de traceur.def :', modname) |
---|
367 | CALL msg('nombre total de traceurs '//TRIM(int2str(nqtrue)), modname) |
---|
368 | DO iq = 1, nqtrue |
---|
369 | CALL msg(strStack([int2str(hadv(iq)), int2str(vadv(iq)), tnom_0(iq), tnom_transp(iq)]), modname) |
---|
370 | END DO |
---|
371 | IF(planet_type /= 'earth') nqo = 0 !--- Same number of tracers in dynamics and physics |
---|
372 | IF(planet_type == 'earth') nqo = COUNT(delPhase(tnom_0) == 'H2O') !--- for all planets except for Earth |
---|
373 | nbtr = nqtrue - nqo |
---|
374 | ALLOCATE(conv_flg(nbtr),pbl_flg(nbtr),solsym(nbtr)) |
---|
375 | conv_flg(1:nbtr) = 1 !--- Convection activated for all tracers |
---|
376 | pbl_flg(1:nbtr) = 1 !--- Boundary layer activated for all tracers |
---|
377 | #endif |
---|
378 | |
---|
379 | !--- SET FIELDS %name, %parent, %phase, %component |
---|
380 | tracers(:)%name = tnom_0 |
---|
381 | tracers(:)%parent = tnom_transp |
---|
382 | tracers(:)%phase = 'g' |
---|
383 | tracers(:)%component = type_trac |
---|
384 | |
---|
385 | |
---|
386 | DO iq = 1, nqtrue |
---|
387 | ip = strIdx([(addPhase('H2O',old_phases(ix:ix),''), ix=1, nphases)], strHead(tracers(iq)%name,'_',.TRUE.)) |
---|
388 | IF(ip == 0) CYCLE |
---|
389 | tracers(iq)%phase = known_phases(ip:ip) |
---|
390 | tracers(iq)%component = 'lmdz' |
---|
391 | END DO |
---|
392 | IF(lINCA) tracers(1+nqo:nqCO2+nqo)%component = 'co2i' |
---|
393 | CALL setGeneration(tracers) !--- SET FIELDS %iGeneration, %gen0Name |
---|
394 | ! manque "type" |
---|
395 | |
---|
396 | !++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
397 | ELSE |
---|
398 | !++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
399 | IF(readTracersFiles(type_trac, fType, tracers)) CALL abort_gcm(modname, 'problem with tracers file(s)',1) |
---|
400 | IF(fType == 0) CALL abort_gcm(modname, 'Missing tracers file: "traceur.def", "tracer.def" or "tracer_<keyword>.def file.',1) |
---|
401 | !--------------------------------------------------------------------------------------------------------------------------- |
---|
402 | IF(fType == 1) THEN !=== FOUND AN OLD STYLE "traceur.def" |
---|
403 | !--------------------------------------------------------------------------------------------------------------------------- |
---|
404 | #ifdef INCA |
---|
405 | nqo = SIZE(tracers) |
---|
406 | IF(nqCO2==1 .AND. nqo==4) nqo = 3 !--- Force nqo to 3 (ThL) |
---|
407 | CALL Init_chem_inca_trac(nqINCA) !--- Get nqINCA from INCA |
---|
408 | nbtr = nqINCA + nqCO2 !--- Number of tracers passed to phytrac |
---|
409 | nqtrue = nbtr + nqo !--- Total number of "true" tracers |
---|
410 | IF(ALL([2,3] /= nqo)) CALL abort_gcm(modname, 'Only 2 or 3 water phases allowed ; found nqo='//TRIM(int2str(nqo)), 1) |
---|
411 | CALL msg('nqo = '//TRIM(int2str(nqo)), modname) |
---|
412 | CALL msg('nbtr = '//TRIM(int2str(nbtr)), modname) |
---|
413 | CALL msg('nqtrue = '//TRIM(int2str(nqtrue)), modname) |
---|
414 | CALL msg('nqCO2 = '//TRIM(int2str(nqCO2)), modname) |
---|
415 | CALL msg('nqINCA = '//TRIM(int2str(nqINCA)), modname) |
---|
416 | ALLOCATE(hadv(nqtrue), hadv_inca(nqINCA), conv_flg_inca(nqINCA), solsym_inca(nqINCA)) |
---|
417 | ALLOCATE(vadv(nqtrue), vadv_inca(nqINCA), pbl_flg_inca(nqINCA)) |
---|
418 | CALL init_transport(hadv_inca, vadv_inca, conv_flg_inca, pbl_flg_inca, solsym_inca) |
---|
419 | ! DC passive CO2 tracer is at position 1: H2O was removed ; nqCO2/=0 in "inco" case only |
---|
420 | |
---|
421 | conv_flg = [( 1 , k=1, nqCO2), conv_flg_inca] |
---|
422 | pbl_flg = [( 1 , k=1, nqCO2), pbl_flg_inca] |
---|
423 | solsym = [('CO2 ', k=1, nqCO2), solsym_inca] |
---|
424 | DEALLOCATE(conv_flg_inca, pbl_flg_inca, solsym_inca) |
---|
425 | ALLOCATE(ttr(nqtrue)) |
---|
426 | ttr(1:nqo+nqCO2) = tracers |
---|
427 | ttr(1 : nqo )%component = 'lmdz' |
---|
428 | ttr(1+nqo:nqCO2+nqo )%component = 'co2i' |
---|
429 | ttr(1+nqo+nqCO2:nqtrue)%component = 'inca' |
---|
430 | ttr(1+nqo :nqtrue)%name = solsym |
---|
431 | ttr(1+nqo+nqCO2:nqtrue)%parent = tran0 |
---|
432 | ttr(1+nqo+nqCO2:nqtrue)%phase = 'g' |
---|
433 | lerr = getKey('hadv', had, ky=ttr(:)%keys); hadv(:) = [had, hadv_inca] |
---|
434 | lerr = getKey('vadv', vad, ky=ttr(:)%keys); vadv(:) = [vad, vadv_inca] |
---|
435 | DEALLOCATE(had, hadv_inca, vad, vadv_inca) |
---|
436 | CALL MOVE_ALLOC(FROM=ttr, TO=tracers) |
---|
437 | CALL setGeneration(tracers) !--- SET FIELDS %iGeneration, %gen0Name |
---|
438 | #else |
---|
439 | nqo = COUNT(delPhase(tracers(:)%name) == 'H2O') !--- Number of water phases |
---|
440 | nqtrue = SIZE(tracers) !--- Total number of "true" tracers |
---|
441 | nbtr = nqtrue - nqo !--- Number of tracers passed to phytrac |
---|
442 | lerr = getKey('hadv', hadv, ky=tracers(:)%keys) |
---|
443 | lerr = getKey('vadv', vadv, ky=tracers(:)%keys) |
---|
444 | ALLOCATE(solsym(nbtr)) |
---|
445 | conv_flg(1:nbtr)=1 !--- Convection activated for all tracers |
---|
446 | pbl_flg(1:nbtr)=1 !--- Boundary layer activated for all tracers |
---|
447 | #endif |
---|
448 | !--------------------------------------------------------------------------------------------------------------------------- |
---|
449 | ELSE !=== FOUND NEW STYLE TRACERS CONFIGURATION FILE(S) |
---|
450 | !--------------------------------------------------------------------------------------------------------------------------- |
---|
451 | nqo = COUNT(delPhase(tracers(:)%name) == 'H2O') !--- Number of water phases |
---|
452 | nqtrue = SIZE(tracers) !--- Total number of "true" tracers |
---|
453 | nbtr = nqtrue - nqo !--- Number of tracers passed to phytrac |
---|
454 | lerr = getKey('hadv', hadv, ky=tracers(:)%keys) |
---|
455 | lerr = getKey('vadv', vadv, ky=tracers(:)%keys) |
---|
456 | ALLOCATE(solsym(nbtr)) |
---|
457 | conv_flg(1:nbtr)=1 !--- Convection activated for all tracers |
---|
458 | pbl_flg(1:nbtr)=1 !--- Boundary layer activated for all tracers |
---|
459 | !--------------------------------------------------------------------------------------------------------------------------- |
---|
460 | END IF |
---|
461 | !--------------------------------------------------------------------------------------------------------------------------- |
---|
462 | END IF |
---|
463 | !++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
464 | |
---|
465 | CALL getKey_init(tracers) |
---|
466 | |
---|
467 | !--- Transfert the number of tracers to Reprobus |
---|
468 | #ifdef REPROBUS |
---|
469 | CALL Init_chem_rep_trac(nbtr, nqo, tracers(:)%name) |
---|
470 | #endif |
---|
471 | |
---|
472 | !============================================================================================================================== |
---|
473 | ! 2) Calculate nqtot, number of tracers needed (greater if advection schemes 20 or 30 have been chosen). |
---|
474 | !============================================================================================================================== |
---|
475 | DO iq = 1, nqtrue |
---|
476 | IF( hadv(iq)<20 .OR. (ANY(hadv(iq)==[20,30]) .AND. hadv(iq)==vadv(iq)) ) CYCLE |
---|
477 | WRITE(msg1,'("The choice hadv=",i0,", vadv=",i0,a)')hadv(iq),vadv(iq),' for "'//TRIM(tracers(iq)%name)//'" is not available' |
---|
478 | CALL abort_gcm(modname, TRIM(msg1), 1) |
---|
479 | END DO |
---|
480 | nqtot = COUNT( hadv< 20 .AND. vadv< 20 ) & !--- No additional tracer |
---|
481 | + 4*COUNT( hadv==20 .AND. vadv==20 ) & !--- 3 additional tracers |
---|
482 | + 10*COUNT( hadv==30 .AND. vadv==30 ) !--- 9 additional tracers |
---|
483 | |
---|
484 | !--- More tracers due to the choice of advection scheme => assign total number of tracers |
---|
485 | IF( nqtot /= nqtrue ) THEN |
---|
486 | CALL msg('The choice of advection scheme for one or more tracers makes it necessary to add tracers') |
---|
487 | CALL msg('The number of true tracers is '//TRIM(int2str(nqtrue))) |
---|
488 | CALL msg('The total number of tracers needed is '//TRIM(int2str(nqtot))) |
---|
489 | END IF |
---|
490 | CALL msg('nqo = '//TRIM(int2str(nqo)), modname) |
---|
491 | CALL msg('nbtr = '//TRIM(int2str(nbtr)), modname) |
---|
492 | CALL msg('nqtrue = '//TRIM(int2str(nqtrue)), modname) |
---|
493 | CALL msg('nqtot = '//TRIM(int2str(nqtot)), modname) |
---|
494 | |
---|
495 | !============================================================================================================================== |
---|
496 | ! 3) Determine the advection scheme choice for water and tracers "iadv" and the fields long name, isAdvected. |
---|
497 | ! iadv = 1 "LMDZ-specific humidity transport" (for H2O vapour) LMV |
---|
498 | ! iadv = 2 backward (for H2O liquid) BAK |
---|
499 | ! iadv = 14 Van-Leer + specific humidity, modified by Francis Codron VLH |
---|
500 | ! iadv = 10 Van-Leer (chosen for vapour and liquid water) VL1 |
---|
501 | ! iadv = 11 Van-Leer for hadv and PPM version (Monotonic) for vadv VLP |
---|
502 | ! iadv = 12 Frederic Hourdin I FH1 |
---|
503 | ! iadv = 13 Frederic Hourdin II FH2 |
---|
504 | ! iadv = 16 Monotonic PPM (Collela & Woodward 1984) PPM |
---|
505 | ! iadv = 17 Semi-monotonic PPM (overshoots allowed) PPS |
---|
506 | ! iadv = 18 Definite positive PPM (overshoots and undershoots allowed) PPP |
---|
507 | ! iadv = 20 Slopes SLP |
---|
508 | ! iadv = 30 Prather PRA |
---|
509 | ! |
---|
510 | ! In array q(ij,l,iq) : iq = 1/2[/3] for vapour/liquid[/ice] water |
---|
511 | ! And optionaly: iq = 3[4],nqtot for other tracers |
---|
512 | !============================================================================================================================== |
---|
513 | ALLOCATE(ttr(nqtot)) |
---|
514 | jq = nqtrue+1; tracers(:)%iadv = -1 |
---|
515 | DO iq = 1, nqtrue |
---|
516 | t1 => tracers(iq) |
---|
517 | |
---|
518 | !--- VERIFY THE CHOICE OF ADVECTION SCHEME |
---|
519 | iad = -1 |
---|
520 | IF(hadv(iq) == vadv(iq) ) iad = hadv(iq) |
---|
521 | IF(hadv(iq)==10 .AND. vadv(iq)==16) iad = 11 |
---|
522 | WRITE(msg1,'("Bad choice of advection scheme for ",a,": hadv = ",i0,", vadv = ",i0)')TRIM(t1%name), hadv(iq), vadv(iq) |
---|
523 | IF(iad == -1) CALL abort_gcm(modname, msg1, 1) |
---|
524 | |
---|
525 | !--- SET FIELDS %longName, %iadv, %isAdvected, %isInPhysics |
---|
526 | t1%longName = t1%name; IF(iad > 0) t1%longName=TRIM(t1%name)//descrq(iad) |
---|
527 | t1%iadv = iad |
---|
528 | t1%isAdvected = iad >= 0 |
---|
529 | t1%isInPhysics= .not. (delPhase(t1%gen0Name) == 'H2O' .and. t1%component=='lmdz') !=== TO BE COMPLETED WITH OTHER EXCEPTIONS: CO2i, SURSATURATED CLOUDS... |
---|
530 | ttr(iq) = t1 |
---|
531 | |
---|
532 | !--- DEFINE THE HIGHER ORDER TRACERS, IF ANY |
---|
533 | nm = 0 |
---|
534 | IF(iad == 20) nm = 3 !--- 2nd order scheme |
---|
535 | IF(iad == 30) nm = 9 !--- 3rd order scheme |
---|
536 | IF(nm == 0) CYCLE !--- No higher moments |
---|
537 | ttr(jq+1:jq+nm) = t1 |
---|
538 | ttr(jq+1:jq+nm)%name = [(TRIM(t1%name) //'-'//TRIM(suff(im)), im=1, nm) ] |
---|
539 | ttr(jq+1:jq+nm)%parent = [(TRIM(t1%parent) //'-'//TRIM(suff(im)), im=1, nm) ] |
---|
540 | ttr(jq+1:jq+nm)%longName = [(TRIM(t1%longName)//'-'//TRIM(suff(im)), im=1, nm) ] |
---|
541 | ttr(jq+1:jq+nm)%iadv = [(-iad, im=1, nm) ] |
---|
542 | ttr(jq+1:jq+nm)%isAdvected = [(.FALSE., im=1, nm) ] |
---|
543 | jq = jq + nm |
---|
544 | END DO |
---|
545 | DEALLOCATE(hadv, vadv) |
---|
546 | CALL MOVE_ALLOC(FROM=ttr, TO=tracers) |
---|
547 | |
---|
548 | !--- SET FIELDS %iqParent, %nqChilds, %iGeneration, %iqDescen, %nqDescen |
---|
549 | CALL indexUpdate(tracers) |
---|
550 | |
---|
551 | CALL msg('Information stored in infotrac :', modname) |
---|
552 | CALL msg('iadv name long_name :', modname) |
---|
553 | |
---|
554 | !=== TEST ADVECTION SCHEME |
---|
555 | DO iq=1,nqtot ; t1 => tracers(iq); iad = t1%iadv |
---|
556 | |
---|
557 | !--- ONLY TESTED VALUES FOR TRACERS FOR NOW: iadv = 14, 10 (and 0 for non-transported tracers) |
---|
558 | IF(ALL([10,14,0] /= iad)) & |
---|
559 | CALL abort_gcm(modname, 'Not tested for iadv='//TRIM(int2str(iad))//' ; 10 or 14 only are allowed !', 1) |
---|
560 | |
---|
561 | !--- ONLY TESTED VALUES FOR PARENTS HAVING CHILDS FOR NOW: iadv = 14, 10 (PARENTS: TRACERS OF GENERATION 1) |
---|
562 | IF(ALL([10,14] /= iad) .AND. t1%iGeneration == 1 .AND. ANY(tracers(:)%iGeneration > 1)) & |
---|
563 | CALL abort_gcm(modname, 'iadv='//TRIM(int2str(iad))//' not implemented for parents ; 10 or 14 only are allowed !', 1) |
---|
564 | |
---|
565 | !--- ONLY TESTED VALUES FOR CHILDS FOR NOW: iadv = 10 (CHILDS: TRACERS OF GENERATION GREATER THAN 1) |
---|
566 | IF(fmsg('WARNING ! iadv='//TRIM(int2str(iad))//' not implemented for childs. Setting iadv=10 for "'//TRIM(t1%name)//'"',& |
---|
567 | modname, iad /= 10 .AND. t1%iGeneration > 1)) t1%iadv = 10 |
---|
568 | |
---|
569 | !--- ONLY VALID SCHEME NUMBER FOR WATER VAPOUR: iadv = 14 |
---|
570 | ll = t1%name /= addPhase('H2O','g'); IF(lOldCode) ll = t1%name /= 'H2Ov' |
---|
571 | IF(fmsg('WARNING ! iadv=14 is valid for water vapour only. Setting iadv=10 for "'//TRIM(t1%name)//'".', & |
---|
572 | modname, iad == 14 .AND. ll)) t1%iadv = 10 |
---|
573 | END DO |
---|
574 | |
---|
575 | IF(lOldCode) THEN |
---|
576 | |
---|
577 | CALL infotrac_setHeredity !--- SET FIELDS %iqParent, %nqChilds, %iGeneration, %gen0Name, %iqDescen, %nqDescen |
---|
578 | CALL infotrac_isoinit !--- SET FIELDS %type, %iso_iName, %iso_iZone, %iso_iPhase |
---|
579 | CALL getKey_init(tracers, isotopes) |
---|
580 | IF(isoSelect('H2O')) RETURN !--- Select water isotopes ; finished if no water isotopes |
---|
581 | iH2O = ixIso !--- Keep track of water family index |
---|
582 | |
---|
583 | !--- Remove the isotopic tracers from the tracers list passed to phytrac |
---|
584 | nbtr = nbtr -nqo* ntiso !--- ISOTOPIC TAGGING TRACERS ARE NOT PASSED TO THE PHYSICS |
---|
585 | nqtottr = nqtot-nqo*(1+ntiso) !--- NO H2O-FAMILY TRACER IS PASSED TO THE PHYSICS |
---|
586 | CALL msg('702: nbtr, ntiso='//strStack(int2str([nbtr, ntiso])), modname) |
---|
587 | CALL msg('704: nqtottr, nqtot, nqo = '//strStack(int2str([nqtottr, nqtot, nqo])), modname) |
---|
588 | ! Rq: nqtottr n'est pas forcement egal a nbtr dans le cas ou nmom/=0 |
---|
589 | IF(COUNT(tracers%iso_iName == 0) - COUNT(delPhase(tracers(:)%name) == 'H2O' .AND. tracers(:)%component=='lmdz') /= nqtottr) & |
---|
590 | CALL abort_gcm('infotrac_init', 'pb dans le calcul de nqtottr', 1) |
---|
591 | |
---|
592 | !--- Finalize : |
---|
593 | DEALLOCATE(tnom_0, tnom_transp) |
---|
594 | |
---|
595 | ELSE |
---|
596 | |
---|
597 | CALL initIsotopes(tracers, isotopes) |
---|
598 | nbIso = SIZE(isotopes); IF(nbIso==0) RETURN !--- No isotopes: finished. |
---|
599 | |
---|
600 | !=== READ PHYSICAL PARAMETERS FROM "isotopes_params.def" FILE SPECIFIC TO WATER ISOTOPES |
---|
601 | ! DONE HERE, AND NOT ONLY IN "infotrac_phy", BECAUSE SOME PHYSICAL PARAMS ARE NEEDED FOR RESTARTS (tnat, alpha_ideal) |
---|
602 | CALL getKey_init(tracers, isotopes) |
---|
603 | IF(isoSelect('H2O')) RETURN !--- Select water isotopes ; finished if no water isotopes |
---|
604 | iH2O = ixIso !--- Keep track of water family index |
---|
605 | IF(getKey('tnat' , tnat, isoName(1:niso))) CALL abort_gcm(modname, 'can''t read "tnat"', 1) |
---|
606 | IF(getKey('alpha', alpha_ideal, isoName(1:niso))) CALL abort_gcm(modname, 'can''t read "alpha_ideal"', 1) |
---|
607 | |
---|
608 | !=== ENSURE THE MEMBERS OF AN ISOTOPES FAMILY ARE PRESENT IN THE SAME PHASES |
---|
609 | DO ix = 1, nbIso |
---|
610 | iso => isotopes(ix) |
---|
611 | !--- Check whether each isotope and tagging isotopic tracer is present in the same number of phases |
---|
612 | DO it = 1, iso%ntiso |
---|
613 | np = SUM([(COUNT(tracers(:)%name == addPhase(iso%trac(it), iso%phase(ip:ip))), ip=1, iso%nphas)]) |
---|
614 | IF(np == iso%nphas) CYCLE |
---|
615 | WRITE(msg1,'("Found ",i0," phases for ",s," instead of ",i0)')np, iso%trac(it), iso%nphas |
---|
616 | CALL abort_gcm(modname, msg1, 1) |
---|
617 | END DO |
---|
618 | DO it = 1, iso%niso |
---|
619 | nz = SUM([(COUNT(iso%trac == iso%trac(it)//'_'//iso%zone(iz)), iz=1, iso%nzone)]) |
---|
620 | IF(nz == iso%nzone) CYCLE |
---|
621 | WRITE(msg1,'("Found ",i0," tagging zones for ",s," instead of ",i0)')nz, iso%trac(it), iso%nzone |
---|
622 | CALL abort_gcm(modname, msg1, 1) |
---|
623 | END DO |
---|
624 | END DO |
---|
625 | nqtottr = COUNT(tracers%iso_iName == 0) |
---|
626 | |
---|
627 | END IF |
---|
628 | |
---|
629 | !=== DISPLAY THE RESULTING LIST |
---|
630 | t => tracers |
---|
631 | CALL msg('Information stored in infotrac :') |
---|
632 | IF(dispTable('isssssssssiiiiiiiii', & |
---|
633 | ['iq ', 'name ', 'longN. ', 'gen0N. ', 'parent ', 'type ', 'phase ', 'compon. ', 'isAdv. ', 'isPhy. '& |
---|
634 | ,'iadv ', 'iGen. ', 'iqPar. ', 'nqDes. ', 'nqChil. ', 'iso_iG. ', 'iso_iN. ', 'iso_iZ. ', 'iso_iP. '], & |
---|
635 | cat(t%name, t%longName, t%gen0Name, t%parent, t%type, t%phase, & |
---|
636 | t%component, bool2str(t%isAdvected), bool2str(t%isInPhysics)), & |
---|
637 | cat([(iq, iq=1, nqtot)], t%iadv, t%iGeneration, t%iqParent, t%nqDescen, & |
---|
638 | t%nqChilds, t%iso_iGroup, t%iso_iName, t%iso_iZone, t%iso_iPhase))) & |
---|
639 | CALL abort_gcm(modname, "problem with the tracers table content", 1) |
---|
640 | |
---|
641 | !--- Some aliases to be removed later |
---|
642 | ntraciso => isotope%ntiso |
---|
643 | ntraceurs_zone => isotope%nzone |
---|
644 | qperemin = min_qParent |
---|
645 | masseqmin = min_qMass |
---|
646 | ratiomin = min_ratio |
---|
647 | CALL msg('end', modname) |
---|
648 | |
---|
649 | END SUBROUTINE infotrac_init |
---|
650 | |
---|
651 | |
---|
652 | |
---|
653 | SUBROUTINE infotrac_setHeredity |
---|
654 | !--- Purpose: Set fields %iqParent, %nqChilds, %iGeneration, %iqDescen, %nqDescen (old method) |
---|
655 | USE strings_mod, ONLY: strIdx |
---|
656 | INTEGER :: iq, ipere, ifils |
---|
657 | INTEGER, ALLOCATABLE :: iqfils(:,:) |
---|
658 | CHARACTER(LEN=maxlen) :: msg1, modname='infotrac_init' |
---|
659 | INCLUDE "iniprint.h" |
---|
660 | |
---|
661 | !=== SET FIELDS %iqParent, %nqChilds |
---|
662 | ALLOCATE(iqfils(nqtot,nqtot)); iqfils(:,:) = 0 |
---|
663 | |
---|
664 | DO iq = 1, nqtot |
---|
665 | msg1 = 'Tracer nr. '//TRIM(int2str(iq))//', called "'//TRIM(tracers(iq)%name)//'" is ' |
---|
666 | |
---|
667 | !--- IS IT A GENERATION 0 TRACER ? IF SO, tracers(iq)%iqParent KEEPS ITS DEFAULT VALUE (0) |
---|
668 | IF(fmsg(TRIM(msg1)//' a parent', modname, tracers(iq)%parent == tran0)) CYCLE |
---|
669 | |
---|
670 | !--- TRACERS OF GENERATION > 0 ; FIND ITS PARENT INDEX |
---|
671 | ipere = strIdx(tracers(:)%name, tracers(iq)%parent) |
---|
672 | IF(ipere == 0) CALL abort_gcm('infotrac_init', TRIM(msg1)//' an orphan', 1) |
---|
673 | IF(iq == ipere) CALL abort_gcm('infotrac_init', TRIM(msg1)//' its own parent',1) |
---|
674 | |
---|
675 | CALL msg(TRIM(msg1)//' the child of '//TRIM(tracers(ipere)%name), modname) |
---|
676 | tracers(iq)%iqParent = ipere |
---|
677 | tracers(ipere)%nqChilds = tracers(ipere)%nqChilds+1 |
---|
678 | iqfils(tracers(ipere)%nqChilds,ipere) = iq |
---|
679 | END DO |
---|
680 | |
---|
681 | CALL msg('nqGen0 = '//int2str(COUNT(tracers(:)%parent == 'air')), modname) |
---|
682 | CALL msg('nqChilds = '//strStack(int2str(tracers(:)%nqChilds)), modname) |
---|
683 | CALL msg('iqParent = '//strStack(int2str(tracers(:)%iqParent)), modname) |
---|
684 | CALL msg('iqChilds = '//strStack(int2str(PACK(iqfils,MASK=.TRUE.))),modname) |
---|
685 | |
---|
686 | |
---|
687 | !=== SET FIELDS %iGeneration, %iqDescen, %nqDescen |
---|
688 | tracers(:)%iGeneration = 0 |
---|
689 | DO iq = 1, nqtot |
---|
690 | ifils = iq |
---|
691 | DO WHILE(tracers(ifils)%iqParent > 0) |
---|
692 | ipere = tracers(ifils)%iqParent |
---|
693 | tracers(ipere)%nqDescen = tracers(ipere)%nqDescen+1 |
---|
694 | tracers(iq)%iGeneration = tracers(iq)%iGeneration+1 |
---|
695 | iqfils(tracers(ipere)%nqDescen,ipere) = iq |
---|
696 | ifils = ipere |
---|
697 | END DO |
---|
698 | msg1 = 'Tracer nr. '//TRIM(int2str(iq))//', called "'//TRIM(tracers(iq)%name)//'" is ' |
---|
699 | CALL msg(TRIM(msg1)//' of generation '//TRIM(int2str(tracers(iq)%iGeneration)), modname) |
---|
700 | END DO |
---|
701 | DO iq=1,nqtot |
---|
702 | tracers(iq)%iqDescen = iqfils(1:tracers(iq)%nqDescen,iq) |
---|
703 | END DO |
---|
704 | |
---|
705 | CALL msg('nqDescen = '//TRIM(strStack(int2str(tracers(:)%nqDescen))), modname) |
---|
706 | CALL msg('nqDescen_tot = ' //TRIM(int2str(SUM(tracers(:)%nqDescen))), modname) |
---|
707 | CALL msg('iqChilds = '//strStack(int2str(PACK(iqfils, MASK=.TRUE.))), modname) |
---|
708 | |
---|
709 | |
---|
710 | END SUBROUTINE infotrac_setHeredity |
---|
711 | |
---|
712 | |
---|
713 | |
---|
714 | SUBROUTINE infotrac_isoinit |
---|
715 | |
---|
716 | #ifdef CPP_IOIPSL |
---|
717 | USE IOIPSL |
---|
718 | #else |
---|
719 | USE ioipsl_getincom |
---|
720 | #endif |
---|
721 | IMPLICIT NONE |
---|
722 | CHARACTER(LEN=3) :: tnom_iso(niso_possibles) |
---|
723 | INTEGER, ALLOCATABLE :: nb_iso(:,:), nb_traciso(:,:), nb_isoind(:) |
---|
724 | INTEGER :: ii, ip, iq, it, iz, ixt, nzone_prec |
---|
725 | TYPE(isot_type), POINTER :: i |
---|
726 | TYPE(trac_type), POINTER :: t(:) |
---|
727 | CHARACTER(LEN=maxlen) :: tnom_trac |
---|
728 | CHARACTER(LEN=maxlen), ALLOCATABLE :: str(:) |
---|
729 | LOGICAL, DIMENSION(:), ALLOCATABLE :: mask |
---|
730 | INCLUDE "iniprint.h" |
---|
731 | |
---|
732 | tnom_iso = ['eau', 'HDO', 'O18', 'O17', 'HTO'] |
---|
733 | ALLOCATE(nb_iso (niso_possibles,nqo)) |
---|
734 | ALLOCATE(nb_traciso (niso_possibles,nqo)) |
---|
735 | ALLOCATE(use_iso (niso_possibles)) |
---|
736 | ALLOCATE(indnum_fn_num(niso_possibles)) |
---|
737 | ALLOCATE(iso_indnum(nqtot)) |
---|
738 | ALLOCATE(nb_isoind(nqo)) |
---|
739 | |
---|
740 | iso_indnum (:) = 0 |
---|
741 | use_iso (:) = .FALSE. |
---|
742 | indnum_fn_num(:) = 0 |
---|
743 | nb_iso (:,:) = 0 |
---|
744 | nb_traciso (:,:) = 0 |
---|
745 | nb_isoind (:) = 0 |
---|
746 | |
---|
747 | DO iq=1, nqtot |
---|
748 | IF(delPhase(tracers(iq)%name) == 'H2O' .OR. .NOT.tracers(iq)%isAdvected) CYCLE |
---|
749 | outer:DO ip = 1, nqo |
---|
750 | DO ixt= 1,niso_possibles |
---|
751 | tnom_trac = 'H2O'//old_phases(ip:ip)//'_'//TRIM(tnom_iso(ixt)) |
---|
752 | IF (tracers(iq)%name == tnom_trac) THEN |
---|
753 | nb_iso(ixt,ip) = nb_iso(ixt,ip)+1 |
---|
754 | nb_isoind (ip) = nb_isoind (ip)+1 |
---|
755 | tracers(iq)%type = 'tracer' |
---|
756 | tracers(iq)%iso_iGroup = 1 |
---|
757 | tracers(iq)%iso_iName = ixt |
---|
758 | iso_indnum(iq) = nb_isoind(ip) |
---|
759 | indnum_fn_num(ixt) = iso_indnum(iq) |
---|
760 | tracers(iq)%iso_iPhase = ip |
---|
761 | EXIT outer |
---|
762 | ELSE IF(tracers(iq)%iqParent> 0) THEN |
---|
763 | IF(tracers(tracers(iq)%iqParent)%name == tnom_trac) THEN |
---|
764 | nb_traciso(ixt,ip) = nb_traciso(ixt,ip)+1 |
---|
765 | iso_indnum(iq) = indnum_fn_num(ixt) |
---|
766 | tracers(iq)%type = 'tag' |
---|
767 | tracers(iq)%iso_iGroup = 1 |
---|
768 | tracers(iq)%iso_iName = ixt |
---|
769 | tracers(iq)%iso_iZone = nb_traciso(ixt,ip) |
---|
770 | tracers(iq)%iso_iPhase = ip |
---|
771 | EXIT outer |
---|
772 | END IF |
---|
773 | END IF |
---|
774 | END DO |
---|
775 | END DO outer |
---|
776 | END DO |
---|
777 | |
---|
778 | niso = 0; nzone_prec = nb_traciso(1,1) |
---|
779 | DO ixt = 1, niso_possibles |
---|
780 | IF(nb_iso(ixt,1) == 0) CYCLE |
---|
781 | IF(nb_iso(ixt,1) /= 1) CALL abort_gcm('infotrac_init', 'Isotopes are not well defined in traceur.def', 1) |
---|
782 | |
---|
783 | ! on verifie que toutes les phases ont le meme nombre d'isotopes |
---|
784 | IF(ANY(nb_iso(ixt,:) /= 1)) CALL abort_gcm('infotrac_init', 'Phases must have same number of isotopes', 1) |
---|
785 | |
---|
786 | niso = niso+1 |
---|
787 | use_iso(ixt) = .TRUE. |
---|
788 | nzone = nb_traciso(ixt,1) |
---|
789 | |
---|
790 | ! on verifie que toutes les phases ont le meme nombre de traceurs d'isotopes |
---|
791 | IF(ANY(nb_traciso(ixt,2:nqo) /= nzone)) CALL abort_gcm('infotrac_init','Phases must have same number of tracers',1) |
---|
792 | |
---|
793 | ! on verifie que tous les isotopes ont le meme nombre de traceurs d'isotopes |
---|
794 | IF(nzone /= nzone_prec) CALL abort_gcm('infotrac_init','Isotope tracers are not well defined in traceur.def',1) |
---|
795 | nzone_prec = nzone |
---|
796 | END DO |
---|
797 | |
---|
798 | ! dimensions et flags isotopiques: |
---|
799 | ntiso = niso*(nzone+1) |
---|
800 | ok_isotopes = niso > 0 |
---|
801 | ok_isotrac = nzone > 0 |
---|
802 | |
---|
803 | IF(ok_isotopes) THEN |
---|
804 | ok_iso_verif = .FALSE.; CALL getin('ok_iso_verif', ok_iso_verif) |
---|
805 | ok_init_iso = .FALSE.; CALL getin('ok_init_iso', ok_init_iso) |
---|
806 | END IF |
---|
807 | tnat = [1.0, 155.76e-6, 2005.2e-6, 0.004/100., 0.0] |
---|
808 | alpha_ideal = [1.0, 1.01, 1.006, 1.003, 1.0] |
---|
809 | ! END IF |
---|
810 | |
---|
811 | ! remplissage du tableau iqiso(ntiso,phase) |
---|
812 | ALLOCATE(iqiso(ntiso,nqo)) |
---|
813 | iqiso(:,:)=0 |
---|
814 | DO iq = 1, nqtot |
---|
815 | IF(tracers(iq)%iso_iName <= 0) CYCLE |
---|
816 | ixt = iso_indnum(iq) + tracers(iq)%iso_iZone*niso |
---|
817 | iqiso(ixt, tracers(iq)%iso_iPhase) = iq |
---|
818 | END DO |
---|
819 | |
---|
820 | ! remplissage du tableau index_trac(nzone,niso) |
---|
821 | ALLOCATE(index_trac(nzone, niso)) |
---|
822 | IF(ok_isotrac) then |
---|
823 | DO ii = 1, niso; index_trac(:, ii) = ii + niso*[(iz, iz=1, nzone)]; END DO |
---|
824 | ELSE |
---|
825 | index_trac(:,:)=0.0 |
---|
826 | END IF |
---|
827 | |
---|
828 | ALLOCATE(isotopes(1)) !--- Only water |
---|
829 | nbIso = 1 |
---|
830 | t => tracers |
---|
831 | i => isotopes(1) |
---|
832 | i%parent = 'H2O' |
---|
833 | |
---|
834 | !--- Isotopes names list (embedded in the "keys" field) |
---|
835 | i%niso = niso |
---|
836 | ALLOCATE(i%keys(i%niso)) |
---|
837 | mask = t%type=='tracer' .AND. delPhase(t%gen0Name)=='H2O' .AND. t%phase == 'g' .AND. t%iGeneration==1 |
---|
838 | str = strTail(PACK(delPhase(t%name), MASK=mask), '_') |
---|
839 | CALL strReduce(str) |
---|
840 | i%keys(:)%name = str |
---|
841 | |
---|
842 | !--- Full isotopes list, with isotopes tagging tracers (if any) following the previous list |
---|
843 | i%ntiso = ntiso |
---|
844 | ALLOCATE(i%trac(i%ntiso)) |
---|
845 | mask = t%type=='tag' .AND. delPhase(t%gen0Name)=='H2O' .AND. t%phase == 'g' .AND. t%iGeneration==2 |
---|
846 | str = PACK(delPhase(t%name), MASK=mask) |
---|
847 | CALL strReduce(str) |
---|
848 | i%trac(:) = [i%keys(:)%name, str] |
---|
849 | |
---|
850 | !--- Tagging zones names list |
---|
851 | i%nzone = nzone |
---|
852 | i%zone = strTail(str, '_', .TRUE.) |
---|
853 | |
---|
854 | !--- Effective phases list |
---|
855 | i%nphas = nqo |
---|
856 | i%phase = '' |
---|
857 | 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 |
---|
858 | |
---|
859 | !--- Table: index in "qx" of an isotope, knowing its indices "it","ip" in "isotope%iName,%iPhase" |
---|
860 | i%iTraPha = RESHAPE([((strIdx(t%name, TRIM(addPhase('H2O', new2oldPhase(i%phase(ip:ip)), ''))//'_'//TRIM(i%trac(it))), & |
---|
861 | it=1,i%ntiso), ip=1,i%nphas)], [i%ntiso,i%nphas]) |
---|
862 | |
---|
863 | !--- Table: index in "isotope%tracs(:)%name" of an isotopic tagging tracer, knowing its indices "iz","ip" in "isotope%iZone,%iName" |
---|
864 | 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 ]) |
---|
865 | DO it=1,ntiso |
---|
866 | WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': iqiso (',it,',:) = '//strStack(int2str(iqiso(it,:))) |
---|
867 | WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': iTraPha(',it,',:) = '//strStack(int2str(i%iTraPha(it,:))) |
---|
868 | END DO |
---|
869 | DO iz=1,nzone |
---|
870 | WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': index_trac(',iz,',:) = '//strStack(int2str(index_trac(iz,:))) |
---|
871 | WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': iZonIso (',iz,',:) = '//strStack(int2str(i%iZonIso(iz,:))) |
---|
872 | END DO |
---|
873 | |
---|
874 | ! Finalize : |
---|
875 | DEALLOCATE(nb_iso) |
---|
876 | |
---|
877 | END SUBROUTINE infotrac_isoinit |
---|
878 | |
---|
879 | |
---|
880 | !============================================================================================================================== |
---|
881 | !=== THE ROUTINE isoSelect IS USED TO SWITCH FROM AN ISOTOPE FAMILY TO ANOTHER: ISOTOPES DEPENDENT PARAMETERS ARE UPDATED |
---|
882 | ! Single generic "isoSelect" routine, using the predefined index of the parent (fast version) or its name (first call). |
---|
883 | !============================================================================================================================== |
---|
884 | LOGICAL FUNCTION isoSelectByName(iName, lVerbose) RESULT(lerr) |
---|
885 | IMPLICIT NONE |
---|
886 | CHARACTER(LEN=*), INTENT(IN) :: iName |
---|
887 | LOGICAL, OPTIONAL, INTENT(IN) :: lVerbose |
---|
888 | INTEGER :: iIso |
---|
889 | LOGICAL :: lV |
---|
890 | lV = .FALSE.; IF(PRESENT(lVerbose)) lV = lVerbose |
---|
891 | iIso = strIdx(isotopes(:)%parent, iName) |
---|
892 | lerr = iIso == 0 |
---|
893 | CALL msg('no isotope family named "'//TRIM(iName)//'"', ll=lerr .AND. lV) |
---|
894 | IF(lerr) RETURN |
---|
895 | lerr = isoSelectByIndex(iIso, lV) |
---|
896 | END FUNCTION isoSelectByName |
---|
897 | !============================================================================================================================== |
---|
898 | LOGICAL FUNCTION isoSelectByIndex(iIso, lVerbose) RESULT(lerr) |
---|
899 | IMPLICIT NONE |
---|
900 | INTEGER, INTENT(IN) :: iIso |
---|
901 | LOGICAL, OPTIONAL, INTENT(IN) :: lVerbose |
---|
902 | LOGICAL :: lv |
---|
903 | lv = .FALSE.; IF(PRESENT(lVerbose)) lv = lVerbose |
---|
904 | lerr = .FALSE. |
---|
905 | IF(iIso == ixIso) RETURN !--- Nothing to do if the index is already OK |
---|
906 | lerr = iIso<=0 .OR. iIso>nbIso |
---|
907 | CALL msg('Inconsistent isotopes family index '//TRIM(int2str(iIso))//': should be > 0 and <= '//TRIM(int2str(nbIso))//'"',& |
---|
908 | ll=lerr .AND. lV) |
---|
909 | IF(lerr) RETURN |
---|
910 | ixIso = iIso !--- Update currently selected family index |
---|
911 | isotope => isotopes(ixIso) !--- Select corresponding component |
---|
912 | isoKeys => isotope%keys; niso = isotope%niso |
---|
913 | isoName => isotope%trac; ntiso = isotope%ntiso |
---|
914 | isoZone => isotope%zone; nzone = isotope%nzone |
---|
915 | isoPhas => isotope%phase; nphas = isotope%nphas |
---|
916 | iZonIso => isotope%iZonIso; isoCheck = isotope%check |
---|
917 | iTraPha => isotope%iTraPha |
---|
918 | END FUNCTION isoSelectByIndex |
---|
919 | !============================================================================================================================== |
---|
920 | |
---|
921 | END MODULE infotrac |
---|