!$Id: infotrac.F90 4082 2022-02-28 08:37:41Z jyg $ ! MODULE infotrac USE strings_mod, ONLY: msg, find, strIdx, strFind, strParse, dispTable, int2str, reduceExpr, & cat, fmsg, test, strTail, strHead, strStack, strReduce, bool2str, maxlen, testFile USE readTracFiles_mod, ONLY: trac_type, readTracersFiles, addPhase, phases_sep, nphases, ancestor, & isot_type, readIsotopesFile, delPhase, old_phases, getKey_init, tran0, & keys_type, initIsotopes, indexUpdate, known_phases, getKey, setGeneration, & new2oldPhase IMPLICIT NONE PRIVATE !=== FOR TRACERS: PUBLIC :: infotrac_init !--- Initialization of the tracers PUBLIC :: tracers, type_trac !--- Full tracers database, tracers type keyword PUBLIC :: nqtot, nbtr, nqo, nqCO2, nqtottr !--- Main dimensions PUBLIC :: solsym, conv_flg, pbl_flg !--- Tracers names + convection & boundary layer activation keys !=== FOR ISOTOPES: General PUBLIC :: isotopes, nbIso !--- Derived type, full isotopes families database + nb of families PUBLIC :: isoSelect, ixIso !--- Isotopes family selection tool + selected family index PUBLIC :: min_qParent, min_qMass, min_ratio !--- Min. values for various isotopic quantities !=== FOR ISOTOPES: Specific to water PUBLIC :: iH2O, tnat, alpha_ideal !--- H2O isotopes index, natural abundance, fractionning coeff. !=== FOR ISOTOPES: Depending on the selected isotopes family PUBLIC :: isotope, isoKeys !--- Selected isotopes database + associated keys (cf. getKey) PUBLIC :: isoName, isoZone, isoPhas !--- Isotopes and tagging zones names, phases PUBLIC :: niso, nzone, nphas, ntiso !--- " " numbers + isotopes & tagging tracers number PUBLIC :: iZonIso, iTraPha !--- 2D index tables to get "iq" index PUBLIC :: isoCheck !--- Run isotopes checking routines !=== FOR BOTH TRACERS AND ISOTOPES PUBLIC :: getKey !--- Get a key from "tracers" or "isotope" PUBLIC :: ntraciso, ntraceurs_zone, iqiso PUBLIC :: ok_isotopes, ok_iso_verif, ok_isotrac, ok_init_iso, use_iso PUBLIC :: index_trac, iso_indnum, indnum_fn_num, niso_possibles PUBLIC :: qperemin, masseqmin, ratiomin INTERFACE isoSelect; MODULE PROCEDURE isoSelectByIndex, isoSelectByName; END INTERFACE isoSelect !=== CONVENTIONS FOR TRACERS NUMBERS: ! |--------------------+-----------------------+-----------------+---------------+----------------------------| ! | water in different | water tagging | water isotopes | other tracers | additional tracers moments | ! | phases: H2O_[gls] | isotopes | | | for higher order schemes | ! |--------------------+-----------------------+-----------------+---------------+----------------------------| ! | | | | | | ! |<-- nqo -->|<-- nqo*niso* nzone -->|<-- nqo*niso -->|<-- nbtr -->|<-- (nmom) -->| ! | | | | ! | |<-- nqo*niso*(nzone+1) = nqo*ntiso -->|<-- nqtottr = nbtr + nmom -->| ! | = nqtot - nqo*(ntiso+1) | ! | | ! |<-- nqtrue = nbtr + nqo*(ntiso+1) -->| | ! | | ! |<-- nqtot = nqtrue + nmom -->| ! | | ! |-----------------------------------------------------------------------------------------------------------| ! NOTES FOR THIS TABLE: ! * Used "niso", "nzone" and "ntiso" are components of "isotopes(ip)" for water (isotopes(ip)%parent == 'H2O'), ! since water is so far the sole tracers family, except passive CO2, removed from the main tracers table. ! * For water, "nqo" is equal to the more general field "isotopes(ip)%nphas". ! * "niso", "nzone", "ntiso", "nphas" are defined for other isotopic tracers families, if any. ! !=== DERIVED TYPE EMBEDDING MOST OF THE TRACERS-RELATED QUANTITIES (LENGTH: nqtot) ! Each entry is accessible using "%" sign. ! |-------------+------------------------------------------------------+-------------+------------------------+ ! | entry | Meaning | Former name | Possible values | ! |-------------+------------------------------------------------------+-------------+------------------------+ ! | name | Name (short) | tname | | ! | gen0Name | Name of the 1st generation ancestor | / | | ! | parent | Name of the parent | / | | ! | longName | Long name (with adv. scheme suffix) for outputs | ttext | | ! | type | Type (so far: tracer or tag) | / | tracer,tag | ! | phase | Phases list ("g"as / "l"iquid / "s"olid) | / | [g][l][s] | ! | component | Name(s) of the merged/cumulated section(s) | / | coma-separated names | ! | iadv | Advection scheme number | iadv | 1-20,30 exc. 3-9,15,19 | ! | iGeneration | Generation (>=1) | / | | ! | isAdvected | advected tracers flag (.TRUE. if iadv >= 0) | / | nqtrue .TRUE. values | ! | isInPhysics | tracers not extracted from the main table in physics | / | nqtottr .TRUE. values | ! | iqParent | Index of the parent tracer | iqpere | 1:nqtot | ! | iqDescen | Indexes of the childs (all generations) | iqfils | 1:nqtot | ! | nqDescen | Number of the descendants (all generations) | nqdesc | 1:nqtot | ! | nqChilds | Number of childs (1st generation only) | nqfils | 1:nqtot | ! | iso_iGroup | Isotopes group index in isotopes(:) | / | 1:nbIso | ! | iso_iName | Isotope name index in isotopes(iso_iGroup)%trac(:) | iso_indnum | 1:niso | ! | iso_iZone | Isotope zone index in isotopes(iso_iGroup)%zone(:) | zone_num | 1:nzone | ! | iso_iPhas | Isotope phase index in isotopes(iso_iGroup)%phas(:) | phase_num | 1:nphas | ! | keys | key/val pairs accessible with "getKey" routine | / | | ! +-------------+------------------------------------------------------+-------------+------------------------+ ! !=== DERIVED TYPE EMBEDDING MOST OF THE ISOTOPES-RELATED QUANTITIES (LENGTH: nbIso, NUMBER OF ISOTOPES FAMILIES) ! Each entry is accessible using "%" sign. ! |-----------------+--------------------------------------------------+----------------+-----------------+ ! | entry | length | Meaning | Former name | Possible values | ! |-----------------+--------------------------------------------------+----------------+-----------------+ ! | parent | Parent tracer (isotopes family name) | | | ! | keys | niso | Isotopes keys/values pairs list + number | | | ! | trac | ntiso | Isotopes + tagging tracers list + number | | | ! | zone | nzone | Geographic tagging zones list + number | | | ! | phase | nphas | Phases list + number | | [g][l][s], 1:3 | ! | niso | Number of isotopes, excluding tagging tracers | | | ! | ntiso | Number of isotopes, including tagging tracers | ntraciso | | ! | nzone | Number of geographic tagging zones | ntraceurs_zone | | ! | nphas | Number of phases | | | ! | iTraPha | Index in "trac(1:niso)" = f(name(1:ntiso)),phas) | iqiso | 1:niso | ! | iZonIso | Index in "trac(1:ntiso)" = f(zone, name(1:niso)) | index_trac | 1:nzone | ! |-----------------+--------------------------------------------------+----------------+-----------------+ REAL, PARAMETER :: min_qParent = 1.e-30, min_qMass = 1.e-18, min_ratio = 1.e-16 ! MVals et CRisi !=== DIMENSIONS OF THE TRACERS TABLES AND OTHER SCALAR VARIABLES INTEGER, SAVE :: nqtot, & !--- Tracers nb in dynamics (incl. higher moments + H2O) nbtr, & !--- Tracers nb in physics (excl. higher moments + H2O) nqo, & !--- Number of water phases nbIso, & !--- Number of available isotopes family nqtottr, & !--- Number of tracers passed to phytrac (TO BE DELETED ?) nqCO2 !--- Number of tracers of CO2 (ThL) CHARACTER(LEN=maxlen), SAVE :: type_trac !--- Keyword for tracers type !=== DERIVED TYPES EMBEDDING MOST INFORMATIONS ABOUT TRACERS AND ISOTOPES FAMILIES TYPE(trac_type), TARGET, SAVE, ALLOCATABLE :: tracers(:) !=== TRACERS DESCRIPTORS VECTOR TYPE(isot_type), TARGET, SAVE, ALLOCATABLE :: isotopes(:) !=== ISOTOPES PARAMETERS VECTOR !=== ALIASES FOR CURRENTLY SELECTED ISOTOPES FAMILY OF VARIABLES EMBEDDED IN THE VECTOR "isotopes" TYPE(isot_type), SAVE, POINTER :: isotope !--- CURRENTLY SELECTED ISOTOPES FAMILY DESCRIPTOR INTEGER, SAVE :: ixIso, iH2O !--- Index of the selected isotopes family and H2O family LOGICAL, SAVE :: isoCheck !--- Flag to trigger the checking routines TYPE(keys_type), SAVE, POINTER :: isoKeys(:) !--- ONE SET OF KEYS FOR EACH ISOTOPE (LISTED IN isoName) CHARACTER(LEN=maxlen), SAVE, POINTER :: isoName(:), & !--- ISOTOPES NAMES FOR THE CURRENTLY SELECTED FAMILY isoZone(:), & !--- TAGGING ZONES FOR THE CURRENTLY SELECTED FAMILY isoPhas !--- USED PHASES FOR THE CURRENTLY SELECTED FAMILY INTEGER, TARGET, SAVE :: niso, nzone, & !--- NUMBER OF ISOTOPES, TAGGING ZONES AND PHASES nphas, ntiso !--- NUMBER OF PHASES AND ISOTOPES + ISOTOPIC TAGGING TRACERS INTEGER, SAVE, POINTER :: iZonIso(:,:) !--- INDEX IN "isoTrac" AS f(tagging zone, isotope) INTEGER, SAVE, POINTER :: iTraPha(:,:) !--- INDEX IN "isoTrac" AS f(isotopic tracer, phase) INTEGER, ALLOCATABLE, SAVE :: index_trac(:,:) ! numero ixt en fn izone, indnum entre 1 et niso INTEGER, ALLOCATABLE, SAVE :: iqiso(:,:) ! donne indice iq en fn de (ixt,phase) !--- Aliases for older names INTEGER, POINTER, SAVE :: ntraciso, ntraceurs_zone REAL, SAVE :: qperemin, masseqmin, ratiomin ! CRisi: cas particulier des isotopes INTEGER, PARAMETER :: niso_possibles = 5 LOGICAL, SAVE :: ok_isotopes, ok_iso_verif, ok_isotrac, ok_init_iso LOGICAL, SAVE, ALLOCATABLE :: use_iso(:) INTEGER, SAVE, ALLOCATABLE :: iso_indnum(:) !--- Gives 1<=idx<=niso_possibles as function(1<=iq <=nqtot) INTEGER, SAVE, ALLOCATABLE :: indnum_fn_num(:) !--- Gives 1<=idx<=niso as function(1<=idx<=niso_possibles) !=== VARIABLES FOR ISOTOPES INITIALIZATION AND FOR INCA REAL, SAVE, ALLOCATABLE :: tnat(:), & !--- Natural relative abundance of water isotope (niso) alpha_ideal(:) !--- Ideal fractionning coefficient (for initial state) (niso) INTEGER, SAVE, ALLOCATABLE :: conv_flg(:), & !--- Convection activation ; needed for INCA (nbtr) pbl_flg(:) !--- Boundary layer activation ; needed for INCA (nbtr) CHARACTER(LEN=8), SAVE, ALLOCATABLE :: solsym(:) !--- Names from INCA (nbtr) LOGICAL, PARAMETER :: lOldCode = .TRUE. CONTAINS SUBROUTINE infotrac_init USE control_mod, ONLY: planet_type, config_inca #ifdef REPROBUS USE CHEM_REP, ONLY: Init_chem_rep_trac #endif IMPLICIT NONE !============================================================================================================================== ! ! Auteur: P. Le Van /L. Fairhead/F.Hourdin ! ------- ! ! Modifications: ! -------------- ! 05/94: F.Forget Modif special traceur ! 02/02: M-A Filiberti Lecture de traceur.def ! 01/22: D. Cugnet Nouveaux tracer.def et tracer_*.def + encapsulation (types tr et iso) ! ! Objet: ! ------ ! GCM LMD nouvelle grille ! !============================================================================================================================== ! ... modification de l'integration de q ( 26/04/94 ) .... !------------------------------------------------------------------------------------------------------------------------------ ! Declarations: INCLUDE "dimensions.h" INCLUDE "iniprint.h" !------------------------------------------------------------------------------------------------------------------------------ ! Local variables INTEGER, ALLOCATABLE :: hadv(:), vadv(:) !--- Horizontal/vertical transport scheme number #ifdef INCA INTEGER, ALLOCATABLE :: had (:), hadv_inca(:), conv_flg_inca(:), &!--- Variables specific to INCA vad (:), vadv_inca(:), pbl_flg_inca(:) CHARACTER(LEN=8), ALLOCATABLE :: solsym_inca(:) !--- Tracers names for INCA INTEGER :: nqINCA #endif CHARACTER(LEN=2) :: suff(9) !--- Suffixes for schemes of order 3 or 4 (Prather) CHARACTER(LEN=3) :: descrq(30) !--- Advection scheme description tags CHARACTER(LEN=maxlen) :: msg1 !--- String for messages CHARACTER(LEN=maxlen), ALLOCATABLE :: str(:) !--- Temporary storage INTEGER :: fType !--- Tracers description file type ; 0: none !--- 1/2/3: "traceur.def"/"tracer.def"/"tracer_*.def" INTEGER :: nqtrue !--- Tracers nb from tracer.def (no higher order moments) INTEGER :: iad !--- Advection scheme number INTEGER :: ic, ip, np, iq, jq, it, nt, im, nm, ix, iz, nz, k !--- Indexes and temporary variables LOGICAL :: lerr, ll CHARACTER(LEN=1) :: p TYPE(trac_type), ALLOCATABLE, TARGET :: ttr(:) TYPE(trac_type), POINTER :: t1, t(:) TYPE(isot_type), POINTER :: iso CHARACTER(LEN=maxlen), ALLOCATABLE :: tnom_0(:), tnom_transp(:) !--- Tracer short name + transporting fluid name CHARACTER(LEN=maxlen) :: tchaine INTEGER :: ierr LOGICAL :: lINCA CHARACTER(LEN=*), PARAMETER :: modname="infotrac_init" !------------------------------------------------------------------------------------------------------------------------------ ! Initialization : !------------------------------------------------------------------------------------------------------------------------------ suff = ['x ','y ','z ','xx','xy','xz','yy','yz','zz'] descrq( 1: 2) = ['LMV','BAK'] descrq(10:20) = ['VL1','VLP','FH1','FH2','VLH',' ','PPM','PPS','PPP',' ','SLP'] descrq(30) = 'PRA' CALL msg('type_trac = "'//TRIM(type_trac)//'"', modname) IF(lOldCode) THEN str = [type_trac]; nt = 1 ELSE IF(strParse(type_trac, ',', str, n=nt)) CALL abort_gcm(modname,'can''t parse "type_trac = '//TRIM(type_trac)//'"',1) END IF !--------------------------------------------------------------------------------------------------------------------------- DO it = 1, nt !--- nt>1=> "type_trac": coma-separated keywords list !--------------------------------------------------------------------------------------------------------------------------- !--- MESSAGE ABOUT THE CHOSEN CONFIGURATION msg1 = 'For type_trac = "'//TRIM(str(it))//'":' SELECT CASE(type_trac) CASE('inca'); CALL msg(TRIM(msg1)//' coupling with INCA chemistry model, config_inca='//config_inca, modname) CASE('inco'); CALL msg(TRIM(msg1)//' coupling jointly with INCA and CO2 cycle', modname) CASE('repr'); CALL msg(TRIM(msg1)//' coupling with REPROBUS chemistry model', modname) CASE('co2i'); CALL msg(TRIM(msg1)//' you have chosen to run with CO2 cycle', modname) CASE('coag'); CALL msg(TRIM(msg1)//' tracers are treated for COAGULATION tests', modname) CASE('lmdz'); CALL msg(TRIM(msg1)//' tracers are treated in LMDZ only', modname) CASE DEFAULT; CALL abort_gcm(modname,'type_trac='//TRIM(str(it))//' not possible yet.',1) END SELECT !--- COHERENCE TEST BETWEEN "type_trac" AND "config_inca" IF(ANY(['inca', 'inco'] == str(it)) .AND. ALL(['aero', 'aeNP', 'chem'] /= config_inca)) & CALL abort_gcm(modname, 'Incoherence between type_trac and config_inca. Please modify "run.def"',1) !--- COHERENCE TEST BETWEEN "type_trac" AND PREPROCESSING KEYS SELECT CASE(str(it)) CASE('inca','inco') #ifndef INCA CALL abort_gcm(modname, 'You must add cpp key INCA and compile with INCA code', 1) #endif CASE('repr') #ifndef REPROBUS CALL abort_gcm(modname, 'You must add cpp key REPROBUS and compile with REPROBUS code', 1) #endif CASE('coag') #ifndef CPP_StratAer CALL abort_gcm(modname, 'You must add cpp key StratAer and compile with StratAer code', 1) #endif END SELECT !--- DISABLE "config_inca" OPTION FOR A RUN WITHOUT "INCA" IF IT DIFFERS FROM "none" IF(fmsg('Setting config_inca="none" as you do not couple with INCA model', & modname, ALL(['inca', 'inco'] /= str(it)) .AND. config_inca /= 'none')) config_inca = 'none' !--------------------------------------------------------------------------------------------------------------------------- END DO !--------------------------------------------------------------------------------------------------------------------------- nqCO2 = 0; IF(ANY(str == 'inco')) nqCO2 = 1 !============================================================================================================================== ! 1) Get the numbers of: true (first order only) tracers "nqtrue", water tracers "nqo" (vapor/liquid/solid) !============================================================================================================================== !++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ IF(lOldCode) THEN !------------------------------------------------------------------------------------------------------------------------------ !--- Determine nqtrue and (INCA only) nqo, nbtr OPEN(90, FILE='traceur.def', FORM='formatted', STATUS='old', IOSTAT=ierr) IF(ierr /= 0) CALL abort_gcm(modname, 'file "traceur.def" not found !', 1) CALL msg('File "traceur.def" successfully opened.', modname) lINCA = ANY(['inca','inco'] == type_trac) IF(lINCA) THEN #ifdef INCA READ(90,*) nqo IF(nqCO2==1 .AND. nqo==4) nqo = 3 CALL Init_chem_inca_trac(nqINCA) nbtr = nqINCA + nqCO2 nqtrue = nbtr + nqo IF(ALL([2,3] /= nqo)) CALL abort_gcm(modname, 'Only 2 or 3 water phases allowed ; found nqo='//TRIM(int2str(nqo)), 1) CALL msg('nqo = '//TRIM(int2str(nqo)), modname) CALL msg('nbtr = '//TRIM(int2str(nbtr)), modname) CALL msg('nqtrue = '//TRIM(int2str(nqtrue)), modname) CALL msg('nqCO2 = '//TRIM(int2str(nqCO2)), modname) CALL msg('nqINCA = '//TRIM(int2str(nqINCA)), modname) ALLOCATE(hadv_inca(nqINCA), conv_flg_inca(nqINCA), solsym_inca(nqINCA)) ALLOCATE(vadv_inca(nqINCA), pbl_flg_inca(nqINCA)) CALL init_transport(hadv_inca, vadv_inca, conv_flg_inca, pbl_flg_inca, solsym_inca) ! DC passive CO2 tracer is at position 1: H2O was removed ; nqCO2/=0 in "inco" case only ALLOCATE(conv_flg(nbtr),pbl_flg(nbtr),solsym(nbtr)) conv_flg = [( 1, ic=1, nqCO2),conv_flg_inca] pbl_flg = [( 1, ic=1, nqCO2), pbl_flg_inca] solsym = [('CO2 ', ic=1, nqCO2), solsym_inca] DEALLOCATE(conv_flg_inca, pbl_flg_inca) #endif ELSE READ(90,*) nqtrue END IF IF (planet_type=="earth" .AND. nqtrue < 2) & CALL abort_gcm('infotrac_init', 'Not enough tracers: nqtrue='//TRIM(int2str(nqtrue))//', 2 tracers is the minimum', 1) !--- Allocate variables depending on nqtrue ALLOCATE(hadv(nqtrue), vadv(nqtrue), tnom_0(nqtrue), tnom_transp(nqtrue), tracers(nqtrue)) !--- Continue to read tracer.def it = 0 DO iq = 1, nqtrue #ifdef INCA IF(iq > nqo+nqCO2) THEN it = it+1 hadv (iq) = hadv_inca (it) vadv (iq) = vadv_inca (it) tnom_0(iq) = solsym_inca(it) tnom_transp(iq) = 'air' CYCLE END IF #endif CALL msg('237: iq='//TRIM(int2str(iq)), modname) READ(90,'(I2,X,I2,X,A)',IOSTAT=ierr) hadv(iq),vadv(iq),tchaine WRITE(msg1,'("hadv(",i0,"), vadv(",i0,") = ",i0,", ",i0)')iq, iq, hadv(iq), vadv(iq) CALL msg(TRIM(msg1), modname) CALL msg('tchaine = "'//TRIM(tchaine)//'"', modname) CALL msg('infotrac 238: IOstatus='//TRIM(int2str(ierr)), modname) IF(ierr/=0) CALL abort_gcm('infotrac_init', 'Pb dans la lecture de traceur.def', 1) jq = INDEX(tchaine(1:LEN_TRIM(tchaine)),' ') CALL msg("Ancienne version de traceur.def: traceurs d'air uniquement", modname, iq==1 .AND. jq==0) CALL msg("Nouvelle version de traceur.def", modname, iq==1 .AND. jq/=0) IF(jq /= 0) THEN !--- Space in the string chain => new format tnom_0 (iq) = tchaine(1:jq-1) tnom_transp(iq) = tchaine(jq+1:) ELSE tnom_0 (iq) = tchaine tnom_transp(iq) = 'air' END IF CALL msg( 'tnom_0(iq)=<'//TRIM(tnom_0(iq)) //'>', modname) CALL msg('tnom_transp(iq)=<'//TRIM(tnom_transp(iq))//'>', modname) END DO #ifdef INCA DEALLOCATE(solsym_inca) #endif CLOSE(90) #ifndef INCA CALL msg('Valeur de traceur.def :', modname) CALL msg('nombre total de traceurs '//TRIM(int2str(nqtrue)), modname) DO iq = 1, nqtrue CALL msg(strStack([int2str(hadv(iq)), int2str(vadv(iq)), tnom_0(iq), tnom_transp(iq)]), modname) END DO IF(planet_type /= 'earth') nqo = 0 !--- Same number of tracers in dynamics and physics IF(planet_type == 'earth') nqo = COUNT(delPhase(tnom_0) == 'H2O') !--- for all planets except for Earth nbtr = nqtrue - nqo ALLOCATE(conv_flg(nbtr),pbl_flg(nbtr),solsym(nbtr)) conv_flg(1:nbtr) = 1 !--- Convection activated for all tracers pbl_flg(1:nbtr) = 1 !--- Boundary layer activated for all tracers #endif !--- SET FIELDS %name, %parent, %phase, %component tracers(:)%name = tnom_0 tracers(:)%parent = tnom_transp tracers(:)%phase = 'g' tracers(:)%component = type_trac DO iq = 1, nqtrue ip = strIdx([(addPhase('H2O',old_phases(ix:ix),''), ix=1, nphases)], strHead(tracers(iq)%name,'_',.TRUE.)) IF(ip == 0) CYCLE tracers(iq)%phase = known_phases(ip:ip) tracers(iq)%component = 'lmdz' END DO IF(lINCA) tracers(1+nqo:nqCO2+nqo)%component = 'co2i' CALL setGeneration(tracers) !--- SET FIELDS %iGeneration, %gen0Name ! manque "type" !++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ ELSE !++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ IF(readTracersFiles(type_trac, fType, tracers)) CALL abort_gcm(modname, 'problem with tracers file(s)',1) IF(fType == 0) CALL abort_gcm(modname, 'Missing tracers file: "traceur.def", "tracer.def" or "tracer_.def file.',1) !--------------------------------------------------------------------------------------------------------------------------- IF(fType == 1) THEN !=== FOUND AN OLD STYLE "traceur.def" !--------------------------------------------------------------------------------------------------------------------------- #ifdef INCA nqo = SIZE(tracers) IF(nqCO2==1 .AND. nqo==4) nqo = 3 !--- Force nqo to 3 (ThL) CALL Init_chem_inca_trac(nqINCA) !--- Get nqINCA from INCA nbtr = nqINCA + nqCO2 !--- Number of tracers passed to phytrac nqtrue = nbtr + nqo !--- Total number of "true" tracers IF(ALL([2,3] /= nqo)) CALL abort_gcm(modname, 'Only 2 or 3 water phases allowed ; found nqo='//TRIM(int2str(nqo)), 1) CALL msg('nqo = '//TRIM(int2str(nqo)), modname) CALL msg('nbtr = '//TRIM(int2str(nbtr)), modname) CALL msg('nqtrue = '//TRIM(int2str(nqtrue)), modname) CALL msg('nqCO2 = '//TRIM(int2str(nqCO2)), modname) CALL msg('nqINCA = '//TRIM(int2str(nqINCA)), modname) ALLOCATE(hadv(nqtrue), hadv_inca(nqINCA), conv_flg_inca(nqINCA), solsym_inca(nqINCA)) ALLOCATE(vadv(nqtrue), vadv_inca(nqINCA), pbl_flg_inca(nqINCA)) CALL init_transport(hadv_inca, vadv_inca, conv_flg_inca, pbl_flg_inca, solsym_inca) ! DC passive CO2 tracer is at position 1: H2O was removed ; nqCO2/=0 in "inco" case only conv_flg = [( 1 , k=1, nqCO2), conv_flg_inca] pbl_flg = [( 1 , k=1, nqCO2), pbl_flg_inca] solsym = [('CO2 ', k=1, nqCO2), solsym_inca] DEALLOCATE(conv_flg_inca, pbl_flg_inca, solsym_inca) ALLOCATE(ttr(nqtrue)) ttr(1:nqo+nqCO2) = tracers ttr(1 : nqo )%component = 'lmdz' ttr(1+nqo:nqCO2+nqo )%component = 'co2i' ttr(1+nqo+nqCO2:nqtrue)%component = 'inca' ttr(1+nqo :nqtrue)%name = solsym ttr(1+nqo+nqCO2:nqtrue)%parent = tran0 ttr(1+nqo+nqCO2:nqtrue)%phase = 'g' lerr = getKey('hadv', had, ky=ttr(:)%keys); hadv(:) = [had, hadv_inca] lerr = getKey('vadv', vad, ky=ttr(:)%keys); vadv(:) = [vad, vadv_inca] DEALLOCATE(had, hadv_inca, vad, vadv_inca) CALL MOVE_ALLOC(FROM=ttr, TO=tracers) CALL setGeneration(tracers) !--- SET FIELDS %iGeneration, %gen0Name #else nqo = COUNT(delPhase(tracers(:)%name) == 'H2O') !--- Number of water phases nqtrue = SIZE(tracers) !--- Total number of "true" tracers nbtr = nqtrue - nqo !--- Number of tracers passed to phytrac lerr = getKey('hadv', hadv, ky=tracers(:)%keys) lerr = getKey('vadv', vadv, ky=tracers(:)%keys) ALLOCATE(solsym(nbtr)) conv_flg(1:nbtr)=1 !--- Convection activated for all tracers pbl_flg(1:nbtr)=1 !--- Boundary layer activated for all tracers #endif !--------------------------------------------------------------------------------------------------------------------------- ELSE !=== FOUND NEW STYLE TRACERS CONFIGURATION FILE(S) !--------------------------------------------------------------------------------------------------------------------------- nqo = COUNT(delPhase(tracers(:)%name) == 'H2O') !--- Number of water phases nqtrue = SIZE(tracers) !--- Total number of "true" tracers nbtr = nqtrue - nqo !--- Number of tracers passed to phytrac lerr = getKey('hadv', hadv, ky=tracers(:)%keys) lerr = getKey('vadv', vadv, ky=tracers(:)%keys) ALLOCATE(solsym(nbtr)) conv_flg(1:nbtr)=1 !--- Convection activated for all tracers pbl_flg(1:nbtr)=1 !--- Boundary layer activated for all tracers !--------------------------------------------------------------------------------------------------------------------------- END IF !--------------------------------------------------------------------------------------------------------------------------- END IF !++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ CALL getKey_init(tracers) !--- Transfert the number of tracers to Reprobus #ifdef REPROBUS CALL Init_chem_rep_trac(nbtr, nqo, tracers(:)%name) #endif !============================================================================================================================== ! 2) Calculate nqtot, number of tracers needed (greater if advection schemes 20 or 30 have been chosen). !============================================================================================================================== DO iq = 1, nqtrue IF( hadv(iq)<20 .OR. (ANY(hadv(iq)==[20,30]) .AND. hadv(iq)==vadv(iq)) ) CYCLE WRITE(msg1,'("The choice hadv=",i0,", vadv=",i0,a)')hadv(iq),vadv(iq),' for "'//TRIM(tracers(iq)%name)//'" is not available' CALL abort_gcm(modname, TRIM(msg1), 1) END DO nqtot = COUNT( hadv< 20 .AND. vadv< 20 ) & !--- No additional tracer + 4*COUNT( hadv==20 .AND. vadv==20 ) & !--- 3 additional tracers + 10*COUNT( hadv==30 .AND. vadv==30 ) !--- 9 additional tracers !--- More tracers due to the choice of advection scheme => assign total number of tracers IF( nqtot /= nqtrue ) THEN CALL msg('The choice of advection scheme for one or more tracers makes it necessary to add tracers') CALL msg('The number of true tracers is '//TRIM(int2str(nqtrue))) CALL msg('The total number of tracers needed is '//TRIM(int2str(nqtot))) END IF CALL msg('nqo = '//TRIM(int2str(nqo)), modname) CALL msg('nbtr = '//TRIM(int2str(nbtr)), modname) CALL msg('nqtrue = '//TRIM(int2str(nqtrue)), modname) CALL msg('nqtot = '//TRIM(int2str(nqtot)), modname) !============================================================================================================================== ! 3) Determine the advection scheme choice for water and tracers "iadv" and the fields long name, isAdvected. ! iadv = 1 "LMDZ-specific humidity transport" (for H2O vapour) LMV ! iadv = 2 backward (for H2O liquid) BAK ! iadv = 14 Van-Leer + specific humidity, modified by Francis Codron VLH ! iadv = 10 Van-Leer (chosen for vapour and liquid water) VL1 ! iadv = 11 Van-Leer for hadv and PPM version (Monotonic) for vadv VLP ! iadv = 12 Frederic Hourdin I FH1 ! iadv = 13 Frederic Hourdin II FH2 ! iadv = 16 Monotonic PPM (Collela & Woodward 1984) PPM ! iadv = 17 Semi-monotonic PPM (overshoots allowed) PPS ! iadv = 18 Definite positive PPM (overshoots and undershoots allowed) PPP ! iadv = 20 Slopes SLP ! iadv = 30 Prather PRA ! ! In array q(ij,l,iq) : iq = 1/2[/3] for vapour/liquid[/ice] water ! And optionaly: iq = 3[4],nqtot for other tracers !============================================================================================================================== ALLOCATE(ttr(nqtot)) jq = nqtrue+1; tracers(:)%iadv = -1 DO iq = 1, nqtrue t1 => tracers(iq) !--- VERIFY THE CHOICE OF ADVECTION SCHEME iad = -1 IF(hadv(iq) == vadv(iq) ) iad = hadv(iq) IF(hadv(iq)==10 .AND. vadv(iq)==16) iad = 11 WRITE(msg1,'("Bad choice of advection scheme for ",a,": hadv = ",i0,", vadv = ",i0)')TRIM(t1%name), hadv(iq), vadv(iq) IF(iad == -1) CALL abort_gcm(modname, msg1, 1) !--- SET FIELDS %longName, %iadv, %isAdvected, %isInPhysics t1%longName = t1%name; IF(iad > 0) t1%longName=TRIM(t1%name)//descrq(iad) t1%iadv = iad t1%isAdvected = iad >= 0 t1%isInPhysics= .not. (delPhase(t1%gen0Name) == 'H2O' .and. t1%component=='lmdz') !=== TO BE COMPLETED WITH OTHER EXCEPTIONS: CO2i, SURSATURATED CLOUDS... ttr(iq) = t1 !--- DEFINE THE HIGHER ORDER TRACERS, IF ANY nm = 0 IF(iad == 20) nm = 3 !--- 2nd order scheme IF(iad == 30) nm = 9 !--- 3rd order scheme IF(nm == 0) CYCLE !--- No higher moments ttr(jq+1:jq+nm) = t1 ttr(jq+1:jq+nm)%name = [(TRIM(t1%name) //'-'//TRIM(suff(im)), im=1, nm) ] ttr(jq+1:jq+nm)%parent = [(TRIM(t1%parent) //'-'//TRIM(suff(im)), im=1, nm) ] ttr(jq+1:jq+nm)%longName = [(TRIM(t1%longName)//'-'//TRIM(suff(im)), im=1, nm) ] ttr(jq+1:jq+nm)%iadv = [(-iad, im=1, nm) ] ttr(jq+1:jq+nm)%isAdvected = [(.FALSE., im=1, nm) ] jq = jq + nm END DO DEALLOCATE(hadv, vadv) CALL MOVE_ALLOC(FROM=ttr, TO=tracers) !--- SET FIELDS %iqParent, %nqChilds, %iGeneration, %iqDescen, %nqDescen CALL indexUpdate(tracers) CALL msg('Information stored in infotrac :', modname) CALL msg('iadv name long_name :', modname) !=== TEST ADVECTION SCHEME DO iq=1,nqtot ; t1 => tracers(iq); iad = t1%iadv !--- ONLY TESTED VALUES FOR TRACERS FOR NOW: iadv = 14, 10 (and 0 for non-transported tracers) IF(ALL([10,14,0] /= iad)) & CALL abort_gcm(modname, 'Not tested for iadv='//TRIM(int2str(iad))//' ; 10 or 14 only are allowed !', 1) !--- ONLY TESTED VALUES FOR PARENTS HAVING CHILDS FOR NOW: iadv = 14, 10 (PARENTS: TRACERS OF GENERATION 1) IF(ALL([10,14] /= iad) .AND. t1%iGeneration == 1 .AND. ANY(tracers(:)%iGeneration > 1)) & CALL abort_gcm(modname, 'iadv='//TRIM(int2str(iad))//' not implemented for parents ; 10 or 14 only are allowed !', 1) !--- ONLY TESTED VALUES FOR CHILDS FOR NOW: iadv = 10 (CHILDS: TRACERS OF GENERATION GREATER THAN 1) IF(fmsg('WARNING ! iadv='//TRIM(int2str(iad))//' not implemented for childs. Setting iadv=10 for "'//TRIM(t1%name)//'"',& modname, iad /= 10 .AND. t1%iGeneration > 1)) t1%iadv = 10 !--- ONLY VALID SCHEME NUMBER FOR WATER VAPOUR: iadv = 14 ll = t1%name /= addPhase('H2O','g'); IF(lOldCode) ll = t1%name /= 'H2Ov' IF(fmsg('WARNING ! iadv=14 is valid for water vapour only. Setting iadv=10 for "'//TRIM(t1%name)//'".', & modname, iad == 14 .AND. ll)) t1%iadv = 10 END DO IF(lOldCode) THEN CALL infotrac_setHeredity !--- SET FIELDS %iqParent, %nqChilds, %iGeneration, %gen0Name, %iqDescen, %nqDescen CALL infotrac_isoinit !--- SET FIELDS %type, %iso_iName, %iso_iZone, %iso_iPhase CALL getKey_init(tracers, isotopes) IF(isoSelect('H2O')) RETURN !--- Select water isotopes ; finished if no water isotopes iH2O = ixIso !--- Keep track of water family index !--- Remove the isotopic tracers from the tracers list passed to phytrac nbtr = nbtr -nqo* ntiso !--- ISOTOPIC TAGGING TRACERS ARE NOT PASSED TO THE PHYSICS nqtottr = nqtot-nqo*(1+ntiso) !--- NO H2O-FAMILY TRACER IS PASSED TO THE PHYSICS CALL msg('702: nbtr, ntiso='//strStack(int2str([nbtr, ntiso])), modname) CALL msg('704: nqtottr, nqtot, nqo = '//strStack(int2str([nqtottr, nqtot, nqo])), modname) ! Rq: nqtottr n'est pas forcement egal a nbtr dans le cas ou nmom/=0 IF(COUNT(tracers%iso_iName == 0) - COUNT(delPhase(tracers(:)%name) == 'H2O' .AND. tracers(:)%component=='lmdz') /= nqtottr) & CALL abort_gcm('infotrac_init', 'pb dans le calcul de nqtottr', 1) !--- Finalize : DEALLOCATE(tnom_0, tnom_transp) ELSE CALL initIsotopes(tracers, isotopes) nbIso = SIZE(isotopes); IF(nbIso==0) RETURN !--- No isotopes: finished. !=== READ PHYSICAL PARAMETERS FROM "isotopes_params.def" FILE SPECIFIC TO WATER ISOTOPES ! DONE HERE, AND NOT ONLY IN "infotrac_phy", BECAUSE SOME PHYSICAL PARAMS ARE NEEDED FOR RESTARTS (tnat, alpha_ideal) CALL getKey_init(tracers, isotopes) IF(isoSelect('H2O')) RETURN !--- Select water isotopes ; finished if no water isotopes iH2O = ixIso !--- Keep track of water family index IF(getKey('tnat' , tnat, isoName(1:niso))) CALL abort_gcm(modname, 'can''t read "tnat"', 1) IF(getKey('alpha', alpha_ideal, isoName(1:niso))) CALL abort_gcm(modname, 'can''t read "alpha_ideal"', 1) !=== ENSURE THE MEMBERS OF AN ISOTOPES FAMILY ARE PRESENT IN THE SAME PHASES DO ix = 1, nbIso iso => isotopes(ix) !--- Check whether each isotope and tagging isotopic tracer is present in the same number of phases DO it = 1, iso%ntiso np = SUM([(COUNT(tracers(:)%name == addPhase(iso%trac(it), iso%phase(ip:ip))), ip=1, iso%nphas)]) IF(np == iso%nphas) CYCLE WRITE(msg1,'("Found ",i0," phases for ",s," instead of ",i0)')np, iso%trac(it), iso%nphas CALL abort_gcm(modname, msg1, 1) END DO DO it = 1, iso%niso nz = SUM([(COUNT(iso%trac == iso%trac(it)//'_'//iso%zone(iz)), iz=1, iso%nzone)]) IF(nz == iso%nzone) CYCLE WRITE(msg1,'("Found ",i0," tagging zones for ",s," instead of ",i0)')nz, iso%trac(it), iso%nzone CALL abort_gcm(modname, msg1, 1) END DO END DO nqtottr = COUNT(tracers%iso_iName == 0) END IF !=== DISPLAY THE RESULTING LIST t => tracers CALL msg('Information stored in infotrac :') IF(dispTable('isssssssssiiiiiiiii', & ['iq ', 'name ', 'longN. ', 'gen0N. ', 'parent ', 'type ', 'phase ', 'compon. ', 'isAdv. ', 'isPhy. '& ,'iadv ', 'iGen. ', 'iqPar. ', 'nqDes. ', 'nqChil. ', 'iso_iG. ', 'iso_iN. ', 'iso_iZ. ', 'iso_iP. '], & cat(t%name, t%longName, t%gen0Name, t%parent, t%type, t%phase, & t%component, bool2str(t%isAdvected), bool2str(t%isInPhysics)), & cat([(iq, iq=1, nqtot)], t%iadv, t%iGeneration, t%iqParent, t%nqDescen, & t%nqChilds, t%iso_iGroup, t%iso_iName, t%iso_iZone, t%iso_iPhase))) & CALL abort_gcm(modname, "problem with the tracers table content", 1) !--- Some aliases to be removed later ntraciso => isotope%ntiso ntraceurs_zone => isotope%nzone qperemin = min_qParent masseqmin = min_qMass ratiomin = min_ratio CALL msg('end', modname) END SUBROUTINE infotrac_init SUBROUTINE infotrac_setHeredity !--- Purpose: Set fields %iqParent, %nqChilds, %iGeneration, %iqDescen, %nqDescen (old method) USE strings_mod, ONLY: strIdx INTEGER :: iq, ipere, ifils INTEGER, ALLOCATABLE :: iqfils(:,:) CHARACTER(LEN=maxlen) :: msg1, modname='infotrac_init' INCLUDE "iniprint.h" !=== SET FIELDS %iqParent, %nqChilds ALLOCATE(iqfils(nqtot,nqtot)); iqfils(:,:) = 0 DO iq = 1, nqtot msg1 = 'Tracer nr. '//TRIM(int2str(iq))//', called "'//TRIM(tracers(iq)%name)//'" is ' !--- IS IT A GENERATION 0 TRACER ? IF SO, tracers(iq)%iqParent KEEPS ITS DEFAULT VALUE (0) IF(fmsg(TRIM(msg1)//' a parent', modname, tracers(iq)%parent == tran0)) CYCLE !--- TRACERS OF GENERATION > 0 ; FIND ITS PARENT INDEX ipere = strIdx(tracers(:)%name, tracers(iq)%parent) IF(ipere == 0) CALL abort_gcm('infotrac_init', TRIM(msg1)//' an orphan', 1) IF(iq == ipere) CALL abort_gcm('infotrac_init', TRIM(msg1)//' its own parent',1) CALL msg(TRIM(msg1)//' the child of '//TRIM(tracers(ipere)%name), modname) tracers(iq)%iqParent = ipere tracers(ipere)%nqChilds = tracers(ipere)%nqChilds+1 iqfils(tracers(ipere)%nqChilds,ipere) = iq END DO CALL msg('nqGen0 = '//int2str(COUNT(tracers(:)%parent == 'air')), modname) CALL msg('nqChilds = '//strStack(int2str(tracers(:)%nqChilds)), modname) CALL msg('iqParent = '//strStack(int2str(tracers(:)%iqParent)), modname) CALL msg('iqChilds = '//strStack(int2str(PACK(iqfils,MASK=.TRUE.))),modname) !=== SET FIELDS %iGeneration, %iqDescen, %nqDescen tracers(:)%iGeneration = 0 DO iq = 1, nqtot ifils = iq DO WHILE(tracers(ifils)%iqParent > 0) ipere = tracers(ifils)%iqParent tracers(ipere)%nqDescen = tracers(ipere)%nqDescen+1 tracers(iq)%iGeneration = tracers(iq)%iGeneration+1 iqfils(tracers(ipere)%nqDescen,ipere) = iq ifils = ipere END DO msg1 = 'Tracer nr. '//TRIM(int2str(iq))//', called "'//TRIM(tracers(iq)%name)//'" is ' CALL msg(TRIM(msg1)//' of generation '//TRIM(int2str(tracers(iq)%iGeneration)), modname) END DO DO iq=1,nqtot tracers(iq)%iqDescen = iqfils(1:tracers(iq)%nqDescen,iq) END DO CALL msg('nqDescen = '//TRIM(strStack(int2str(tracers(:)%nqDescen))), modname) CALL msg('nqDescen_tot = ' //TRIM(int2str(SUM(tracers(:)%nqDescen))), modname) CALL msg('iqChilds = '//strStack(int2str(PACK(iqfils, MASK=.TRUE.))), modname) END SUBROUTINE infotrac_setHeredity SUBROUTINE infotrac_isoinit #ifdef CPP_IOIPSL USE IOIPSL #else USE ioipsl_getincom #endif IMPLICIT NONE CHARACTER(LEN=3) :: tnom_iso(niso_possibles) INTEGER, ALLOCATABLE :: nb_iso(:,:), nb_traciso(:,:), nb_isoind(:) INTEGER :: ii, ip, iq, it, iz, ixt, nzone_prec TYPE(isot_type), POINTER :: i TYPE(trac_type), POINTER :: t(:) CHARACTER(LEN=maxlen) :: tnom_trac CHARACTER(LEN=maxlen), ALLOCATABLE :: str(:) LOGICAL, DIMENSION(:), ALLOCATABLE :: mask INCLUDE "iniprint.h" tnom_iso = ['eau', 'HDO', 'O18', 'O17', 'HTO'] ALLOCATE(nb_iso (niso_possibles,nqo)) ALLOCATE(nb_traciso (niso_possibles,nqo)) ALLOCATE(use_iso (niso_possibles)) ALLOCATE(indnum_fn_num(niso_possibles)) ALLOCATE(iso_indnum(nqtot)) ALLOCATE(nb_isoind(nqo)) iso_indnum (:) = 0 use_iso (:) = .FALSE. indnum_fn_num(:) = 0 nb_iso (:,:) = 0 nb_traciso (:,:) = 0 nb_isoind (:) = 0 DO iq=1, nqtot IF(delPhase(tracers(iq)%name) == 'H2O' .OR. .NOT.tracers(iq)%isAdvected) CYCLE outer:DO ip = 1, nqo DO ixt= 1,niso_possibles tnom_trac = 'H2O'//old_phases(ip:ip)//'_'//TRIM(tnom_iso(ixt)) IF (tracers(iq)%name == tnom_trac) THEN nb_iso(ixt,ip) = nb_iso(ixt,ip)+1 nb_isoind (ip) = nb_isoind (ip)+1 tracers(iq)%type = 'tracer' tracers(iq)%iso_iGroup = 1 tracers(iq)%iso_iName = ixt iso_indnum(iq) = nb_isoind(ip) indnum_fn_num(ixt) = iso_indnum(iq) tracers(iq)%iso_iPhase = ip EXIT outer ELSE IF(tracers(iq)%iqParent> 0) THEN IF(tracers(tracers(iq)%iqParent)%name == tnom_trac) THEN nb_traciso(ixt,ip) = nb_traciso(ixt,ip)+1 iso_indnum(iq) = indnum_fn_num(ixt) tracers(iq)%type = 'tag' tracers(iq)%iso_iGroup = 1 tracers(iq)%iso_iName = ixt tracers(iq)%iso_iZone = nb_traciso(ixt,ip) tracers(iq)%iso_iPhase = ip EXIT outer END IF END IF END DO END DO outer END DO niso = 0; nzone_prec = nb_traciso(1,1) DO ixt = 1, niso_possibles IF(nb_iso(ixt,1) == 0) CYCLE IF(nb_iso(ixt,1) /= 1) CALL abort_gcm('infotrac_init', 'Isotopes are not well defined in traceur.def', 1) ! on verifie que toutes les phases ont le meme nombre d'isotopes IF(ANY(nb_iso(ixt,:) /= 1)) CALL abort_gcm('infotrac_init', 'Phases must have same number of isotopes', 1) niso = niso+1 use_iso(ixt) = .TRUE. nzone = nb_traciso(ixt,1) ! on verifie que toutes les phases ont le meme nombre de traceurs d'isotopes IF(ANY(nb_traciso(ixt,2:nqo) /= nzone)) CALL abort_gcm('infotrac_init','Phases must have same number of tracers',1) ! on verifie que tous les isotopes ont le meme nombre de traceurs d'isotopes IF(nzone /= nzone_prec) CALL abort_gcm('infotrac_init','Isotope tracers are not well defined in traceur.def',1) nzone_prec = nzone END DO ! dimensions et flags isotopiques: ntiso = niso*(nzone+1) ok_isotopes = niso > 0 ok_isotrac = nzone > 0 IF(ok_isotopes) THEN ok_iso_verif = .FALSE.; CALL getin('ok_iso_verif', ok_iso_verif) ok_init_iso = .FALSE.; CALL getin('ok_init_iso', ok_init_iso) END IF tnat = [1.0, 155.76e-6, 2005.2e-6, 0.004/100., 0.0] alpha_ideal = [1.0, 1.01, 1.006, 1.003, 1.0] ! END IF ! remplissage du tableau iqiso(ntiso,phase) ALLOCATE(iqiso(ntiso,nqo)) iqiso(:,:)=0 DO iq = 1, nqtot IF(tracers(iq)%iso_iName <= 0) CYCLE ixt = iso_indnum(iq) + tracers(iq)%iso_iZone*niso iqiso(ixt, tracers(iq)%iso_iPhase) = iq END DO ! remplissage du tableau index_trac(nzone,niso) ALLOCATE(index_trac(nzone, niso)) IF(ok_isotrac) then DO ii = 1, niso; index_trac(:, ii) = ii + niso*[(iz, iz=1, nzone)]; END DO ELSE index_trac(:,:)=0.0 END IF ALLOCATE(isotopes(1)) !--- Only water nbIso = 1 t => tracers i => isotopes(1) i%parent = 'H2O' !--- Isotopes names list (embedded in the "keys" field) i%niso = niso ALLOCATE(i%keys(i%niso)) mask = t%type=='tracer' .AND. delPhase(t%gen0Name)=='H2O' .AND. t%phase == 'g' .AND. t%iGeneration==1 str = strTail(PACK(delPhase(t%name), MASK=mask), '_') CALL strReduce(str) i%keys(:)%name = str !--- Full isotopes list, with isotopes tagging tracers (if any) following the previous list i%ntiso = ntiso ALLOCATE(i%trac(i%ntiso)) mask = t%type=='tag' .AND. delPhase(t%gen0Name)=='H2O' .AND. t%phase == 'g' .AND. t%iGeneration==2 str = PACK(delPhase(t%name), MASK=mask) CALL strReduce(str) i%trac(:) = [i%keys(:)%name, str] !--- Tagging zones names list i%nzone = nzone i%zone = strTail(str, '_', .TRUE.) !--- Effective phases list i%nphas = nqo i%phase = '' 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 !--- Table: index in "qx" of an isotope, knowing its indices "it","ip" in "isotope%iName,%iPhase" i%iTraPha = RESHAPE([((strIdx(t%name, TRIM(addPhase('H2O', new2oldPhase(i%phase(ip:ip)), ''))//'_'//TRIM(i%trac(it))), & it=1,i%ntiso), ip=1,i%nphas)], [i%ntiso,i%nphas]) !--- Table: index in "isotope%tracs(:)%name" of an isotopic tagging tracer, knowing its indices "iz","ip" in "isotope%iZone,%iName" 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 ]) DO it=1,ntiso WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': iqiso (',it,',:) = '//strStack(int2str(iqiso(it,:))) WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': iTraPha(',it,',:) = '//strStack(int2str(i%iTraPha(it,:))) END DO DO iz=1,nzone WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': index_trac(',iz,',:) = '//strStack(int2str(index_trac(iz,:))) WRITE(lunout,'(a,i0,a)')TRIM('infotrac_init')//': iZonIso (',iz,',:) = '//strStack(int2str(i%iZonIso(iz,:))) END DO ! Finalize : DEALLOCATE(nb_iso) END SUBROUTINE infotrac_isoinit !============================================================================================================================== !=== THE ROUTINE isoSelect IS USED TO SWITCH FROM AN ISOTOPE FAMILY TO ANOTHER: ISOTOPES DEPENDENT PARAMETERS ARE UPDATED ! Single generic "isoSelect" routine, using the predefined index of the parent (fast version) or its name (first call). !============================================================================================================================== LOGICAL FUNCTION isoSelectByName(iName, lVerbose) RESULT(lerr) IMPLICIT NONE CHARACTER(LEN=*), INTENT(IN) :: iName LOGICAL, OPTIONAL, INTENT(IN) :: lVerbose INTEGER :: iIso LOGICAL :: lV lV = .FALSE.; IF(PRESENT(lVerbose)) lV = lVerbose iIso = strIdx(isotopes(:)%parent, iName) lerr = iIso == 0 CALL msg('no isotope family named "'//TRIM(iName)//'"', ll=lerr .AND. lV) IF(lerr) RETURN lerr = isoSelectByIndex(iIso, lV) END FUNCTION isoSelectByName !============================================================================================================================== LOGICAL FUNCTION isoSelectByIndex(iIso, lVerbose) RESULT(lerr) IMPLICIT NONE INTEGER, INTENT(IN) :: iIso LOGICAL, OPTIONAL, INTENT(IN) :: lVerbose LOGICAL :: lv lv = .FALSE.; IF(PRESENT(lVerbose)) lv = lVerbose lerr = .FALSE. IF(iIso == ixIso) RETURN !--- Nothing to do if the index is already OK lerr = iIso<=0 .OR. iIso>nbIso CALL msg('Inconsistent isotopes family index '//TRIM(int2str(iIso))//': should be > 0 and <= '//TRIM(int2str(nbIso))//'"',& ll=lerr .AND. lV) IF(lerr) RETURN ixIso = iIso !--- Update currently selected family index isotope => isotopes(ixIso) !--- Select corresponding component isoKeys => isotope%keys; niso = isotope%niso isoName => isotope%trac; ntiso = isotope%ntiso isoZone => isotope%zone; nzone = isotope%nzone isoPhas => isotope%phase; nphas = isotope%nphas iZonIso => isotope%iZonIso; isoCheck = isotope%check iTraPha => isotope%iTraPha END FUNCTION isoSelectByIndex !============================================================================================================================== END MODULE infotrac