#ifdef ISO ! $Id: $ MODULE isotopes_mod USE strings_mod, ONLY: msg, real2str, int2str, bool2str, maxlen, strIdx, strStack USE infotrac_phy, ONLY: isoName, isoSelect, niso, ntiso, nbIso, isoFamilies USE iso_params_mod IMPLICIT NONE INTERFACE get_in; MODULE PROCEDURE getinp_s, getinp_i, getinp_r, getinp_l; END INTERFACE get_in SAVE !--- Contains all isotopic variables + their initialization !--- Isotopes-specific routines are in isotopes_routines_mod to avoid circular dependencies with isotopes_verif_mod. !--- Negligible lower thresholds: no need to check for absurd values under these lower limits REAL, PARAMETER :: & ridicule = 1e-12, & ! For mixing ratios ridicule_rain = 1e-8, & ! For rain fluxes (rain, zrfl...) in kg/s <-> 1e-3 mm/day ridicule_evap = ridicule_rain*1e-2, & ! For evaporations in kg/s <-> 1e-3 mm/day ridicule_qsol = ridicule_rain, & ! For qsol in kg <-> 1e-8 kg ridicule_snow = ridicule_qsol ! For snow in kg <-> 1e-8 kg REAL, PARAMETER :: expb_max = 30.0 !--- Fractionation coefficients for H217O REAL, PARAMETER :: fac_coeff_eq17_liq = 0.529, & fac_coeff_eq17_ice = 0.529 !--- H218O reference REAL, PARAMETER :: fac_enrichoce18 = 0.0005, alpha_liq_sol_O18 = 1.00291, & talph1_O18 = 1137., talps1_O18 = 11.839, tkcin0_O18 = 0.006, & talph2_O18 = -0.4156, talps2_O18 = -0.028244, tkcin1_O18 = 0.000285, & talph3_O18 = -2.0667E-3, tdifrel_O18 = 1./0.9723, tkcin2_O18 = 0.00082 !--- Parameters that do not depend on the nature of water isotopes: REAL, PARAMETER :: pxtmelt = 273.15 !--- temperature at which ice formation starts REAL, PARAMETER :: pxtice = 273.15 - 10.0 !--- temperature at which all condensate is ice: REAL, PARAMETER :: pxtmin = 273.15 - 120.0 !--- computation done only under -120°C REAL, PARAMETER :: pxtmax = 273.15 + 60.0 !--- computation done only over +60°C REAL, PARAMETER :: tdifexp = 0.58 !--- a constant for alpha_eff for equilibrium below cloud base: REAL, PARAMETER :: tv0cin = 7.0 !--- wind threshold (m/s) for smooth/rough regime (Merlivat and Jouzel 1979) REAL, PARAMETER :: musi = 1.0 !--- facteurs lambda et mu dans Si=musi-lambda*T REAL, PARAMETER :: Kd = 2.5e-9 ! m2/s !--- diffusion in soil REAL, PARAMETER :: rh_cste_surf_cond = 0.6 !--- cste_surf_cond case: rhs and/or Ts set to constants REAL, PARAMETER :: T_cste_surf_cond = 288.0 !--- Isotopes indices (in [1,niso] ; non-existing => 0 index) INTEGER, SAVE :: iso_eau, iso_HDO, iso_O18, iso_O17, iso_HTO !$OMP THREADPRIVATE(iso_eau, iso_HDO, iso_O18, iso_O17, iso_HTO) INTEGER, SAVE :: ntracisoOR !$OMP THREADPRIVATE(ntracisoOR) !--- Variables not depending on isotopes REAL, SAVE :: thumxt1 !$OMP THREADPRIVATE(thumxt1) INTEGER, SAVE :: ntot !$OMP THREADPRIVATE(ntot) REAL, SAVE :: h_land_ice !$OMP THREADPRIVATE(h_land_ice) REAL, SAVE :: P_veg !$OMP THREADPRIVATE(P_veg) REAL, SAVE :: lambda_sursat !$OMP THREADPRIVATE(lambda_sursat) LOGICAL, SAVE :: bidouille_anti_divergence ! T: regularly, xteau <- q to avoid slow drifts !$OMP THREADPRIVATE(bidouille_anti_divergence) LOGICAL, SAVE :: essai_convergence ! F: as in LMDZ without isotopes (bad for isotopes) !$OMP THREADPRIVATE(essai_convergence) INTEGER, SAVE :: initialisation_iso ! 0: file ; 1: R=0 ; 2: R=distill. Rayleigh ; 3: R=Rsmow !$OMP THREADPRIVATE(initialisation_iso) INTEGER, SAVE :: modif_SST ! 0: default ; 1: modified SST ; 2, 3: SST profiles !$OMP THREADPRIVATE(modif_SST) REAL, SAVE :: deltaTtest ! Uniform modification of the SST !$OMP THREADPRIVATE(deltaTtest) INTEGER, SAVE :: modif_sic ! Holes in the Sea Ice !$OMP THREADPRIVATE(modif_sic) REAL, SAVE :: deltasic ! Minimal holes fraction !$OMP THREADPRIVATE(deltasic) REAL, SAVE :: deltaTtestpoles !$OMP THREADPRIVATE(deltaTtestpoles) REAL, SAVE :: sstlatcrit, dsstlatcrit !$OMP THREADPRIVATE(sstlatcrit, dsstlatcrit) INTEGER, SAVE :: albedo_prescrit ! 0: default ; 1: constant albedo !$OMP THREADPRIVATE(albedo_prescrit) REAL, SAVE :: lon_min_albedo, lon_max_albedo, lat_min_albedo, lat_max_albedo !$OMP THREADPRIVATE(lon_min_albedo, lon_max_albedo, lat_min_albedo, lat_max_albedo) REAL, SAVE :: deltaP_BL,tdifexp_sol !$OMP THREADPRIVATE(deltaP_BL,tdifexp_sol) INTEGER, SAVE :: ruissellement_pluie, alphak_stewart !$OMP THREADPRIVATE(ruissellement_pluie, alphak_stewart) INTEGER, SAVE :: calendrier_guide !$OMP THREADPRIVATE(calendrier_guide) INTEGER, SAVE :: cste_surf_cond !$OMP THREADPRIVATE(cste_surf_cond) REAL, SAVE :: mixlen !$OMP THREADPRIVATE(mixlen) INTEGER, SAVE :: evap_cont_cste !$OMP THREADPRIVATE(evap_cont_cste) REAL, SAVE :: deltaO18_evap_cont, d_evap_cont !$OMP THREADPRIVATE(deltaO18_evap_cont, d_evap_cont) INTEGER, SAVE :: nudge_qsol, region_nudge_qsol !$OMP THREADPRIVATE(nudge_qsol, region_nudge_qsol) INTEGER, SAVE :: nlevmaxO17 !$OMP THREADPRIVATE(nlevmaxO17) INTEGER, SAVE :: no_pce !$OMP THREADPRIVATE(no_pce) REAL, SAVE :: A_satlim !$OMP THREADPRIVATE(A_satlim) INTEGER, SAVE :: ok_restrict_A_satlim, modif_ratqs !$OMP THREADPRIVATE(ok_restrict_A_satlim, modif_ratqs) REAL, SAVE :: Pcrit_ratqs, ratqsbasnew !$OMP THREADPRIVATE(Pcrit_ratqs, ratqsbasnew) REAL, SAVE :: fac_modif_evaoce !$OMP THREADPRIVATE(fac_modif_evaoce) REAL, SAVE :: deltaO18_oce !$OMP THREADPRIVATE(deltaO18_oce) INTEGER, SAVE :: ok_bidouille_wake !$OMP THREADPRIVATE(ok_bidouille_wake) LOGICAL, SAVE :: cond_temp_env !$OMP THREADPRIVATE(cond_temp_env) !--- Vectors of length "niso" REAL, ALLOCATABLE, DIMENSION(:), SAVE :: & alpha, tnat, toce, tcorr, tdifrel !$OMP THREADPRIVATE(alpha, tnat, toce, tcorr, tdifrel) REAL, ALLOCATABLE, DIMENSION(:), SAVE :: & talph1, talph2, talph3, talps1, talps2 !$OMP THREADPRIVATE(talph1, talph2, talph3, talps1, talps2) REAL, ALLOCATABLE, DIMENSION(:), SAVE :: & tkcin0, tkcin1, tkcin2 !$OMP THREADPRIVATE(tkcin0, tkcin1, tkcin2) REAL, ALLOCATABLE, DIMENSION(:), SAVE :: & alpha_liq_sol, Rdefault, Rmethox !$OMP THREADPRIVATE(alpha_liq_sol, Rdefault, Rmethox) !--- Specific to HTO: LOGICAL, SAVE :: ok_prod_nucl_tritium !--- TRUE => HTO production by nuclear tests !$OMP THREADPRIVATE(ok_prod_nucl_tritium) INTEGER, PARAMETER :: nessai = 486 INTEGER, DIMENSION(nessai), SAVE :: & day_nucl, month_nucl, year_nucl !$OMP THREADPRIVATE(day_nucl, month_nucl, year_nucl) REAL, DIMENSION(nessai), SAVE :: & lat_nucl, lon_nucl, zmin_nucl, zmax_nucl, HTO_nucl !$OMP THREADPRIVATE(lat_nucl, lon_nucl, zmin_nucl, zmax_nucl, HTO_nucl) CONTAINS SUBROUTINE iso_init() IMPLICIT NONE !=== Local variables: INTEGER :: ixt, ii, is LOGICAL :: ltnat1 CHARACTER(LEN=maxlen) :: modname, sxt !--- For H2[17]O REAL :: fac_kcin, pente_MWL !--- Sensitivity tests LOGICAL, PARAMETER :: ok_nocinsfc = .FALSE. ! if T: no kinetic effect in sfc evap LOGICAL, PARAMETER :: ok_nocinsat = .FALSE. ! if T: no sursaturation effect for ice LOGICAL, PARAMETER :: Rdefault_smow = .FALSE. ! if T: Rdefault=smow; if F: nul !--- For [3]H INTEGER :: iessai modname = 'iso_init' CALL msg('219: entree', modname) !-------------------------------------------------------------- ! General: !-------------------------------------------------------------- !--- Check number of isotopes CALL msg('64: niso = '//TRIM(int2str(niso)), modname) DO ii = 1, nbIso CALL msg('Can''t select isotopes class "'//TRIM(isoFamilies(ii))//'"', modname, isoSelect(ii, lVerbose=.TRUE.)) !============================================================================================================================== IF(isoFamilies(ii) == 'H2O') THEN !============================================================================================================================== !--- Init de ntracisoOR: on ecrasera en cas de traceurs de tagging isotopiques ! (nzone>0) si complications avec ORCHIDEE ntracisoOR = ntiso !--- Type of water isotopes: iso_eau = strIdx(isoName, 'H216O'); CALL msg('iso_eau='//int2str(iso_eau), modname) iso_HDO = strIdx(isoName, 'HDO'); CALL msg('iso_HDO='//int2str(iso_HDO), modname) iso_O18 = strIdx(isoName, 'H218O'); CALL msg('iso_O18='//int2str(iso_O18), modname) iso_O17 = strIdx(isoName, 'H217O'); CALL msg('iso_O17='//int2str(iso_O17), modname) iso_HTO = strIdx(isoName, 'HTO'); CALL msg('iso_HTO='//int2str(iso_HTO), modname) !--- Initialisation: reading the isotopic parameters. CALL get_in('lambda', lambda_sursat, 0.004); IF(ok_nocinsat) lambda_sursat = 0. CALL get_in('thumxt1', thumxt1, 0.75*1.2) CALL get_in('ntot', ntot, 20, .FALSE.) CALL get_in('h_land_ice', h_land_ice, 20., .FALSE.) CALL get_in('P_veg', P_veg, 1.0, .FALSE.) CALL get_in('bidouille_anti_divergence', bidouille_anti_divergence, .FALSE.) CALL get_in('essai_convergence', essai_convergence, .FALSE.) CALL get_in('initialisation_iso', initialisation_iso, 0) ! IF(nzone>0 .AND. initialisation_iso==0) & ! CALL get_in('initialisation_isotrac',initialisation_isotrac) CALL get_in('modif_sst', modif_sst, 0) CALL get_in('deltaTtest', deltaTtest, 0.0) !--- For modif_sst>=1 CALL get_in('deltaTtestpoles',deltaTtestpoles, 0.0) !--- For modif_sst>=2 CALL get_in( 'sstlatcrit', sstlatcrit, 30.0) !--- For modif_sst>=3 CALL get_in('dsstlatcrit', dsstlatcrit, 0.0) !--- For modif_sst>=3 #ifdef ISOVERIF CALL msg('iso_init 270: sstlatcrit='//real2str( sstlatcrit), modname, sstlatcrit < 0.0) !--- For modif_sst>=2 CALL msg('iso_init 279: dsstlatcrit='//real2str(dsstlatcrit), modname, sstlatcrit < 0.0) !--- For modif_sst>=3 IF(modif_sst >= 2 .AND. sstlatcrit < 0.0) STOP #endif CALL get_in('modif_sic', modif_sic, 0) IF(modif_sic >= 1) & CALL get_in('deltasic', deltasic, 0.1) CALL get_in('albedo_prescrit', albedo_prescrit, 0) IF(albedo_prescrit == 1) THEN CALL get_in('lon_min_albedo', lon_min_albedo, -200.) CALL get_in('lon_max_albedo', lon_max_albedo, 200.) CALL get_in('lat_min_albedo', lat_min_albedo, -100.) CALL get_in('lat_max_albedo', lat_max_albedo, 100.) END IF CALL get_in('deltaO18_oce', deltaO18_oce, 0.0) CALL get_in('deltaP_BL', deltaP_BL, 10.0) CALL get_in('ruissellement_pluie', ruissellement_pluie, 0) CALL get_in('alphak_stewart', alphak_stewart, 1) CALL get_in('tdifexp_sol', tdifexp_sol, 0.67) CALL get_in('calendrier_guide', calendrier_guide, 0) CALL get_in('cste_surf_cond', cste_surf_cond, 0) CALL get_in('mixlen', mixlen, 35.0) CALL get_in('evap_cont_cste', evap_cont_cste, 0) CALL get_in('deltaO18_evap_cont', deltaO18_evap_cont,0.0) CALL get_in('d_evap_cont', d_evap_cont, 0.0) CALL get_in('nudge_qsol', nudge_qsol, 0) CALL get_in('region_nudge_qsol', region_nudge_qsol, 1) nlevmaxO17 = 50 CALL msg('nlevmaxO17='//TRIM(int2str(nlevmaxO17))) CALL get_in('no_pce', no_pce, 0) CALL get_in('A_satlim', A_satlim, 1.0) CALL get_in('ok_restrict_A_satlim', ok_restrict_A_satlim, 0) #ifdef ISOVERIF CALL msg(' 315: A_satlim='//real2str(A_satlim), modname, A_satlim > 1.0) IF(A_satlim > 1.0) STOP #endif ! CALL get_in('slope_limiterxy', slope_limiterxy, 2.0) ! CALL get_in('slope_limiterz', slope_limiterz, 2.0) CALL get_in('modif_ratqs', modif_ratqs, 0) CALL get_in('Pcrit_ratqs', Pcrit_ratqs, 500.0) CALL get_in('ratqsbasnew', ratqsbasnew, 0.05) CALL get_in('fac_modif_evaoce', fac_modif_evaoce, 1.0) CALL get_in('ok_bidouille_wake', ok_bidouille_wake, 0) ! si oui, la temperature de cond est celle de l'environnement, pour eviter ! bugs quand temperature dans ascendances convs est mal calculee CALL get_in('cond_temp_env', cond_temp_env, .FALSE.) IF(ANY(isoName == 'HTO')) & CALL get_in('ok_prod_nucl_tritium', ok_prod_nucl_tritium, .FALSE., .FALSE.) CALL get_in('tnateq1', ltnat1, .TRUE.) CALL msg('iso_O18, iso_HDO, iso_eau = '//TRIM(strStack(int2str([iso_O18, iso_HDO, iso_eau]))), modname) !-------------------------------------------------------------- ! Parameters that depend on the nature of water isotopes: !-------------------------------------------------------------- ALLOCATE(tnat (niso), talph1(niso), talps1(niso), tkcin0(niso), tdifrel (niso), alpha (niso)) ALLOCATE(toce (niso), talph2(niso), talps2(niso), tkcin1(niso), Rdefault(niso), alpha_liq_sol(niso)) ALLOCATE(tcorr(niso), talph3(niso), tkcin2(niso), Rmethox (niso)) !=== H216O is = iso_eau IF(is /= 0) THEN tdifrel (is) = 1.0 alpha (is) = alpha_ideal_H216O tnat (is) = tnat_H216O; IF(ltnat1) tnat(is) = 1.0 toce (is) = tnat(is) tcorr (is) = 1.0 talph1 (is) = 0.0; talps1(is) = 0.0; tkcin0(is) = 0.0 talph2 (is) = 0.0; talps2(is) = 0.0; tkcin1(is) = 0.0 talph3 (is) = 0.0; tkcin2(is) = 0.0 Rdefault(is) = tnat(is)*1.0 Rmethox (is) = 1.0 alpha_liq_sol(is) = 1.0 END IF !=== H217O is = iso_O17 IF(is /= 0) THEN; pente_MWL = 0.528 tdifrel (is) = 1./0.98555 ! used in 1D and in LdG model ; tdifrel=1./0.985452: from Amaelle alpha (is) = alpha_ideal_H217O tnat (is) = tnat_H217O; IF(ltnat1) tnat(is) = 1.0 toce (is) = tnat(is)*(1.0+deltaO18_oce/1000.0)**pente_MWL tcorr (is) = 1.0+fac_enrichoce18*pente_MWL fac_kcin = (tdifrel(is)-1.0)/(tdifrel_O18-1.0) ! fac_kcin=0.5145: from Amaelle talph1 (is) = talph1_O18; talps1(is) = talps1_O18; tkcin0(is) = tkcin0_O18*fac_kcin talph2 (is) = talph2_O18; talps2(is) = talps2_O18; tkcin1(is) = tkcin1_O18*fac_kcin talph3 (is) = talph3_O18; tkcin2(is) = tkcin2_O18*fac_kcin Rdefault(is) = tnat(is)*(1.0-3.15/1000.) Rmethox (is) = tnat(is)*(1.0+230./1000.) alpha_liq_sol(is) = alpha_liq_sol_O18**fac_coeff_eq17_liq END IF !=== H218O is = iso_O18 IF(is /= 0) THEN tdifrel (is) = tdifrel_O18 alpha (is) = alpha_ideal_H218O tnat (is) = tnat_H218O; IF(ltnat1) tnat(is) = 1.0 toce (is) = tnat(is)*(1.0+deltaO18_oce/1000.0) tcorr (is) = 1.0+fac_enrichoce18 talph1 (is) = talph1_O18; talps1(is) = talps1_O18; tkcin0(is) = tkcin0_O18 talph2 (is) = talph2_O18; talps2(is) = talps2_O18; tkcin1(is) = tkcin1_O18 talph3 (is) = talph3_O18; tkcin2(is) = tkcin2_O18 Rdefault(is) = tnat(is)*(1.0-6.00/1000.) Rmethox (is) = tnat(is)*(1.0+130./1000.) ! Zahn & al. 2006 alpha_liq_sol(is) = alpha_liq_sol_O18 END IF !=== HDO is = iso_HDO IF(is /= 0) THEN; pente_MWL = 8.0 tdifrel (is) = 1./0.9755 ! fac_kcin=0.88 alpha (is) = alpha_ideal_HDO tnat (is) = tnat_HDO; IF(ltnat1) tnat(is) = 1.0 toce (is) = tnat(is)*(1.0+deltaO18_oce/1000.0*pente_MWL) tcorr (is) = 1.0+fac_enrichoce18*pente_MWL fac_kcin = (tdifrel(is)-1.0)/(tdifrel_O18-1.0) talph1 (is) = 24844.; talps1(is) = 16288.; tkcin0(is) = tkcin0_O18*fac_kcin talph2 (is) = -76.248; talps2(is) = -0.0934; tkcin1(is) = tkcin1_O18*fac_kcin talph3 (is) = 52.612E-3; tkcin2(is) = tkcin2_O18*fac_kcin Rdefault(is) = tnat(is)*(1.0+(10.0-6.0*pente_MWL)/1000.) Rmethox (is) = tnat(is)*(1.0-25.0/1000.) alpha_liq_sol(is) = 1.0212 ! Lehmann & Siegenthaler, 1991, Jo. of Glaciology, vol 37, p 23 ! alpha_liq_sol=1.0192: Weston, Ralph, 1955 END IF !=== HTO is = iso_HTO IF(is /= 0) THEN tdifrel (is) = 1./0.968 alpha (is) = alpha_ideal_HTO tnat (is) = tnat_HTO; IF(ltnat1) tnat(is) = 1.0 toce (is) = 4.0E-19 ! ratio T/H = 0.2 TU Dreisigacker & Roether 1978 tcorr (is) = 1.0 talph1 (is) = 46480.; talps1(is) = 46480.; tkcin0(is) = 0.01056 talph2 (is) = -103.87; talps2(is) = -103.87; tkcin1(is) = 0.0005016 talph3 (is) = 0.0; tkcin2(is) = 0.0014432 Rdefault(is) = 0.0 Rmethox (is) = 0.0 alpha_liq_sol(is) = 1.0 END IF IF(.NOT. Rdefault_smow) THEN Rdefault(:) = 0.0; IF(iso_eau > 0) Rdefault(iso_eau) = 1.0 END IF WRITE(*,*) 'Rdefault = ',Rdefault WRITE(*,*) 'toce = ', toce !--- Sensitivity test: no kinetic effect in sfc evaporation IF(ok_nocinsfc) THEN tkcin0(1:niso) = 0.0 tkcin1(1:niso) = 0.0 tkcin2(1:niso) = 0.0 END IF CALL msg('285: verif initialisation:', modname) DO ixt=1,niso sxt=int2str(ixt) CALL msg(' * isoName('//TRIM(sxt)//') = <'//TRIM(isoName(ixt))//'>', modname) CALL msg( ' tnat('//TRIM(sxt)//') = '//TRIM(real2str(tnat(ixt))), modname) ! CALL msg(' alpha_liq_sol('//TRIM(sxt)//') = '//TRIM(real2str(alpha_liq_sol(ixt))), modname) ! CALL msg( ' tkcin0('//TRIM(sxt)//') = '//TRIM(real2str(tkcin0(ixt))), modname) ! CALL msg( ' tdifrel('//TRIM(sxt)//') = '//TRIM(real2str(tdifrel(ixt))), modname) END DO CALL msg('69: lambda = '//TRIM(real2str(lambda_sursat)), modname) CALL msg('69: thumxt1 = '//TRIM(real2str(thumxt1)), modname) CALL msg('69: h_land_ice = '//TRIM(real2str(h_land_ice)), modname) CALL msg('69: P_veg = '//TRIM(real2str(P_veg)), modname) !============================================================================================================================== ELSE !============================================================================================================================== CALL abort_physic('"isotopes_mod" is not set up yet for isotopes family "'//TRIM(isoFamilies(ii))//'"', modname, 1) !============================================================================================================================== END IF !============================================================================================================================== END DO END SUBROUTINE iso_init SUBROUTINE getinp_s(nam, val, def, lDisp) USE ioipsl_getincom, ONLY: getin USE mod_phys_lmdz_mpi_data, ONLY : is_mpi_root USE mod_phys_lmdz_omp_data, ONLY : is_omp_root USE mod_phys_lmdz_transfert_para, ONLY : bcast CHARACTER(LEN=*), INTENT(IN) :: nam CHARACTER(LEN=*), INTENT(INOUT) :: val CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: def LOGICAL, OPTIONAL, INTENT(IN) :: lDisp LOGICAL :: lD !$OMP BARRIER IF(is_mpi_root.AND.is_omp_root) THEN IF(PRESENT(def)) val=def; CALL getin(nam,val) lD=.TRUE.; IF(PRESENT(lDisp)) lD=lDisp IF(lD) CALL msg(TRIM(nam)//' = '//TRIM(val)) END IF CALL bcast(val) END SUBROUTINE getinp_s SUBROUTINE getinp_i(nam, val, def, lDisp) USE ioipsl_getincom, ONLY: getin USE mod_phys_lmdz_mpi_data, ONLY : is_mpi_root USE mod_phys_lmdz_omp_data, ONLY : is_omp_root USE mod_phys_lmdz_transfert_para, ONLY : bcast CHARACTER(LEN=*), INTENT(IN) :: nam INTEGER, INTENT(INOUT) :: val INTEGER, OPTIONAL, INTENT(IN) :: def LOGICAL, OPTIONAL, INTENT(IN) :: lDisp LOGICAL :: lD !$OMP BARRIER IF(is_mpi_root.AND.is_omp_root) THEN IF(PRESENT(def)) val=def; CALL getin(nam,val) lD=.TRUE.; IF(PRESENT(lDisp)) lD=lDisp IF(lD) CALL msg(TRIM(nam)//' = '//TRIM(int2str(val))) END IF CALL bcast(val) END SUBROUTINE getinp_i SUBROUTINE getinp_r(nam, val, def, lDisp) USE ioipsl_getincom, ONLY: getin USE mod_phys_lmdz_mpi_data, ONLY : is_mpi_root USE mod_phys_lmdz_omp_data, ONLY : is_omp_root USE mod_phys_lmdz_transfert_para, ONLY : bcast CHARACTER(LEN=*), INTENT(IN) :: nam REAL, INTENT(INOUT) :: val REAL, OPTIONAL, INTENT(IN) :: def LOGICAL, OPTIONAL, INTENT(IN) :: lDisp LOGICAL :: lD !$OMP BARRIER IF(is_mpi_root.AND.is_omp_root) THEN IF(PRESENT(def)) val=def; CALL getin(nam,val) lD=.TRUE.; IF(PRESENT(lDisp)) lD=lDisp IF(lD) CALL msg(TRIM(nam)//' = '//TRIM(real2str(val))) END IF CALL bcast(val) END SUBROUTINE getinp_r SUBROUTINE getinp_l(nam, val, def, lDisp) USE ioipsl_getincom, ONLY: getin USE mod_phys_lmdz_mpi_data, ONLY : is_mpi_root USE mod_phys_lmdz_omp_data, ONLY : is_omp_root USE mod_phys_lmdz_transfert_para, ONLY : bcast CHARACTER(LEN=*), INTENT(IN) :: nam LOGICAL, INTENT(INOUT) :: val LOGICAL, OPTIONAL, INTENT(IN) :: def LOGICAL, OPTIONAL, INTENT(IN) :: lDisp LOGICAL :: lD !$OMP BARRIER IF(is_mpi_root.AND.is_omp_root) THEN IF(PRESENT(def)) val=def; CALL getin(nam,val) lD=.TRUE.; IF(PRESENT(lDisp)) lD=lDisp IF(lD) CALL msg(TRIM(nam)//' = '//TRIM(bool2str(val))) END IF CALL bcast(val) END SUBROUTINE getinp_l END MODULE isotopes_mod #endif