source: LMDZ6/trunk/libf/phylmdiso/isotopes_mod.F90 @ 4158

Last change on this file since 4158 was 4158, checked in by dcugnet, 2 years ago
  • include new type "tracer.def" tracers description files in DefLists/?
  • Remove the old code (lOldCode=.TRUE.)
  • Remove the hard-coded part in isotopes_mod that was ensuring the convergence => this version will not converge with previous isotopic physics, but the

difference, only due to round off errors on "tdifrel", is acceptable.

  • Property svn:executable set to *
File size: 18.0 KB
Line 
1#ifdef ISO
2! $Id: $
3
4MODULE isotopes_mod
5   USE strings_mod,  ONLY: msg, real2str, int2str, bool2str, maxlen, strIdx, strStack
6   USE infotrac_phy, ONLY: isoName
7   IMPLICIT NONE
8   INTERFACE get_in; MODULE PROCEDURE getinp_s, getinp_i, getinp_r, getinp_l;  END INTERFACE get_in
9   SAVE
10
11  !--- Contains all isotopic variables + their initialization
12  !--- Isotopes-specific routines are in isotopes_routines_mod to avoid circular dependencies with isotopes_verif_mod.
13
14   !--- Isotopes indices (in [1,niso] ; non-existing => 0 index)
15   INTEGER, SAVE :: iso_eau, iso_HDO, iso_O18, iso_O17, iso_HTO
16!$OMP THREADPRIVATE(iso_eau, iso_HDO, iso_O18, iso_O17, iso_HTO)
17
18   INTEGER, SAVE :: ntracisoOR
19!$OMP THREADPRIVATE(ntracisoOR)
20
21   !--- Variables not depending on isotopes
22   REAL,    SAVE :: pxtmelt, pxtice, pxtmin, pxtmax
23!$OMP THREADPRIVATE(pxtmelt, pxtice, pxtmin, pxtmax)
24   REAL,    SAVE :: tdifexp, tv0cin, thumxt1
25!$OMP THREADPRIVATE(tdifexp, tv0cin, thumxt1)
26   INTEGER, SAVE :: ntot
27!$OMP THREADPRIVATE(ntot)
28   REAL,    SAVE :: h_land_ice
29!$OMP THREADPRIVATE(h_land_ice)
30   REAL,    SAVE :: P_veg
31!$OMP THREADPRIVATE(P_veg)
32   REAL,    SAVE :: musi, lambda_sursat
33!$OMP THREADPRIVATE(musi, lambda_sursat)
34   REAL,    SAVE :: Kd
35!$OMP THREADPRIVATE(Kd)
36   REAL,    SAVE :: rh_cste_surf_cond, T_cste_surf_cond
37!$OMP THREADPRIVATE(rh_cste_surf_cond, T_cste_surf_cond)
38   LOGICAL, SAVE :: bidouille_anti_divergence    ! T: regularly, xteau <- q to avoid slow drifts
39!$OMP THREADPRIVATE(bidouille_anti_divergence)
40   LOGICAL, SAVE :: essai_convergence            ! F: as in LMDZ without isotopes (bad for isotopes)
41!$OMP THREADPRIVATE(essai_convergence)
42   INTEGER, SAVE :: initialisation_iso           ! 0: file ; 1: R=0 ; 2: R=distill. Rayleigh ; 3: R=Rsmow
43!$OMP THREADPRIVATE(initialisation_iso)
44   INTEGER, SAVE :: modif_SST                    ! 0: default ; 1: modified SST ; 2, 3: SST profiles
45!$OMP THREADPRIVATE(modif_SST)
46   REAL,    SAVE :: deltaTtest                   ! Uniform modification of the SST
47!$OMP THREADPRIVATE(deltaTtest)
48   INTEGER, SAVE :: modif_sic                    ! Holes in the Sea Ice
49!$OMP THREADPRIVATE(modif_sic)
50   REAL,    SAVE :: deltasic                     ! Minimal holes fraction
51!$OMP THREADPRIVATE(deltasic)
52   REAL,    SAVE :: deltaTtestpoles
53!$OMP THREADPRIVATE(deltaTtestpoles)
54   REAL,    SAVE :: sstlatcrit, dsstlatcrit
55!$OMP THREADPRIVATE(sstlatcrit, dsstlatcrit)
56   REAL,    SAVE :: deltaO18_oce
57!$OMP THREADPRIVATE(deltaO18_oce)
58   INTEGER, SAVE :: albedo_prescrit              ! 0: default ; 1: constant albedo
59!$OMP THREADPRIVATE(albedo_prescrit)
60   REAL,    SAVE :: lon_min_albedo, lon_max_albedo, lat_min_albedo, lat_max_albedo
61!$OMP THREADPRIVATE(lon_min_albedo, lon_max_albedo, lat_min_albedo, lat_max_albedo)
62   REAL,    SAVE :: deltaP_BL,tdifexp_sol
63!$OMP THREADPRIVATE(deltaP_BL,tdifexp_sol)
64   INTEGER, SAVE :: ruissellement_pluie, alphak_stewart
65!$OMP THREADPRIVATE(ruissellement_pluie, alphak_stewart)
66   INTEGER, SAVE :: calendrier_guide
67!$OMP THREADPRIVATE(calendrier_guide)
68   INTEGER, SAVE :: cste_surf_cond
69!$OMP THREADPRIVATE(cste_surf_cond)
70   REAL,    SAVE :: mixlen
71!$OMP THREADPRIVATE(mixlen)
72   INTEGER, SAVE :: evap_cont_cste
73!$OMP THREADPRIVATE(evap_cont_cste)
74   REAL,    SAVE :: deltaO18_evap_cont, d_evap_cont
75!$OMP THREADPRIVATE(deltaO18_evap_cont, d_evap_cont)
76   INTEGER, SAVE :: nudge_qsol, region_nudge_qsol
77!$OMP THREADPRIVATE(nudge_qsol, region_nudge_qsol)
78   INTEGER, SAVE :: nlevmaxO17
79!$OMP THREADPRIVATE(nlevmaxO17)
80   INTEGER, SAVE :: no_pce
81!$OMP THREADPRIVATE(no_pce)
82   REAL,    SAVE :: A_satlim
83!$OMP THREADPRIVATE(A_satlim)
84   INTEGER, SAVE :: ok_restrict_A_satlim, modif_ratqs
85!$OMP THREADPRIVATE(ok_restrict_A_satlim, modif_ratqs)
86   REAL,    SAVE :: Pcrit_ratqs, ratqsbasnew
87!$OMP THREADPRIVATE(Pcrit_ratqs, ratqsbasnew)
88   REAL,    SAVE :: fac_modif_evaoce
89!$OMP THREADPRIVATE(fac_modif_evaoce)
90   INTEGER, SAVE :: ok_bidouille_wake
91!$OMP THREADPRIVATE(ok_bidouille_wake)
92   LOGICAL, SAVE :: cond_temp_env
93!$OMP THREADPRIVATE(cond_temp_env)
94
95   !--- Vectors of length "niso"
96   REAL, ALLOCATABLE, DIMENSION(:), SAVE :: &
97                    tnat, toce, tcorr, tdifrel
98!$OMP THREADPRIVATE(tnat, toce, tcorr, tdifrel)
99   REAL, ALLOCATABLE, DIMENSION(:), SAVE :: &
100                    talph1, talph2, talph3, talps1, talps2
101!$OMP THREADPRIVATE(talph1, talph2, talph3, talps1, talps2)
102   REAL, ALLOCATABLE, DIMENSION(:), SAVE :: &
103                    tkcin0, tkcin1, tkcin2
104!$OMP THREADPRIVATE(tkcin0, tkcin1, tkcin2)
105   REAL, ALLOCATABLE, DIMENSION(:), SAVE :: &
106                    alpha_liq_sol, Rdefault, Rmethox
107!$OMP THREADPRIVATE(alpha_liq_sol, Rdefault, Rmethox)
108   REAL, SAVE ::    fac_coeff_eq17_liq, fac_coeff_eq17_ice
109!$OMP THREADPRIVATE(fac_coeff_eq17_liq, fac_coeff_eq17_ice)
110
111   !--- Negligible lower thresholds: no need to check for absurd values under these lower limits
112   REAL, PARAMETER :: &
113      ridicule      = 1e-12,              & ! For mixing ratios
114      ridicule_rain = 1e-8,               & ! For rain fluxes (rain, zrfl...) in kg/s <-> 1e-3 mm/day
115      ridicule_evap = ridicule_rain*1e-2, & ! For evaporations                in kg/s <-> 1e-3 mm/day
116      ridicule_qsol = ridicule_rain,      & ! For qsol                        in kg <-> 1e-8 kg
117      ridicule_snow = ridicule_qsol         ! For snow                        in kg <-> 1e-8 kg
118   REAL, PARAMETER :: expb_max = 30.0
119
120   !--- Specific to HTO:
121   LOGICAL, SAVE :: ok_prod_nucl_tritium    !--- TRUE => HTO production by nuclear tests
122!$OMP THREADPRIVATE(ok_prod_nucl_tritium)
123   INTEGER, PARAMETER :: nessai = 486
124   INTEGER, DIMENSION(nessai), SAVE :: &
125                    day_nucl, month_nucl, year_nucl
126!$OMP THREADPRIVATE(day_nucl, month_nucl, year_nucl)
127   REAL,    DIMENSION(nessai), SAVE :: &
128                    lat_nucl, lon_nucl, zmin_nucl, zmax_nucl, HTO_nucl
129!$OMP THREADPRIVATE(lat_nucl, lon_nucl, zmin_nucl, zmax_nucl, HTO_nucl)
130 
131 
132CONTAINS
133
134SUBROUTINE iso_init()
135   USE ioipsl_getin_p_mod, ONLY: getin_p
136   USE infotrac_phy,       ONLY: ntiso, niso, getKey
137    USE strings_mod,       ONLY: maxlen
138   IMPLICIT NONE
139
140   !=== Local variables:
141   INTEGER :: ixt
142
143   !--- H2[18]O reference
144   REAL :: fac_enrichoce18, alpha_liq_sol_O18, &
145           talph1_O18, talph2_O18, talph3_O18, talps1_O18, talps2_O18, &
146           tkcin0_O18, tkcin1_O18, tkcin2_O18, tdifrel_O18 
147
148   !--- For H2[17]O
149   REAL    :: fac_kcin, pente_MWL
150     
151   !--- Sensitivity tests
152   LOGICAL, PARAMETER ::   ok_nocinsfc = .FALSE. ! if T: no kinetic effect in sfc evap
153   LOGICAL, PARAMETER ::   ok_nocinsat = .FALSE. ! if T: no sursaturation effect for ice
154   LOGICAL, PARAMETER :: Rdefault_smow = .FALSE. ! if T: Rdefault=smow; if F: nul
155
156   !--- For [3]H
157   INTEGER :: iessai
158
159   CHARACTER(LEN=maxlen) :: modname, sxt
160   REAL, ALLOCATABLE :: tmp(:)
161
162   modname = 'iso_init'
163   CALL msg('219: entree', modname)
164
165   !--------------------------------------------------------------
166   ! General:
167   !--------------------------------------------------------------
168
169   !--- Check number of isotopes
170   CALL msg('64: niso = '//TRIM(int2str(niso)), modname)
171
172   !--- Init de ntracisoOR: on ecrasera en cas de traceurs de tagging isotopiques
173   !                     (nzone>0) si complications avec ORCHIDEE
174   ntracisoOR = ntiso 
175
176   !--- Type of water isotopes:
177   iso_eau = strIdx(isoName, 'H2[16]O'); CALL msg('iso_eau='//int2str(iso_eau), modname)
178   iso_HDO = strIdx(isoName, 'H[2]HO');  CALL msg('iso_HDO='//int2str(iso_HDO), modname)
179   iso_O18 = strIdx(isoName, 'H2[18]O'); CALL msg('iso_O18='//int2str(iso_O18), modname)
180   iso_O17 = strIdx(isoName, 'H2[17]O'); CALL msg('iso_O17='//int2str(iso_O17), modname)
181   iso_HTO = strIdx(isoName, 'H[3]HO');  CALL msg('iso_HTO='//int2str(iso_HTO), modname)
182
183   ! initialisation
184   ! lecture des parametres isotopiques:
185   ! pour que ca marche en openMP, il faut utiliser getin_p. Car le getin ne peut
186   ! etre appele que par un thread a la fois, et ca pose tout un tas de problemes,
187   ! d'ou tout un tas de magouilles comme dans conf_phys_m. A terme, tout le monde
188   ! lira par getin_p.
189   CALL get_in('lambda',     lambda_sursat, 0.004)
190   CALL get_in('thumxt1',    thumxt1,       0.75*1.2)
191   CALL get_in('ntot',       ntot,          20,  .FALSE.)
192   CALL get_in('h_land_ice', h_land_ice,    20., .FALSE.)
193   CALL get_in('P_veg',      P_veg,         1.0, .FALSE.)
194   CALL get_in('bidouille_anti_divergence', bidouille_anti_divergence, .FALSE.)
195   CALL get_in('essai_convergence',         essai_convergence,         .FALSE.)
196   CALL get_in('initialisation_iso',        initialisation_iso,        0)
197
198!  IF(nzone>0 .AND. initialisation_iso==0) &
199!      CALL get_in('initialisation_isotrac',initialisation_isotrac)
200   CALL get_in('modif_sst',      modif_sst,         0)
201   CALL get_in('deltaTtest',     deltaTtest,      0.0)     !--- For modif_sst>=1
202   CALL get_in('deltaTtestpoles',deltaTtestpoles, 0.0)     !--- For modif_sst>=2
203   CALL get_in( 'sstlatcrit',    sstlatcrit,     30.0)     !--- For modif_sst>=3
204   CALL get_in('dsstlatcrit',   dsstlatcrit,      0.0)     !--- For modif_sst>=3
205#ifdef ISOVERIF
206   CALL msg('iso_init 270:  sstlatcrit='//real2str( sstlatcrit), modname, sstlatcrit < 0.0) !--- For modif_sst>=2
207   CALL msg('iso_init 279: dsstlatcrit='//real2str(dsstlatcrit), modname, sstlatcrit < 0.0) !--- For modif_sst>=3
208   IF(modif_sst >= 2 .AND. sstlatcrit < 0.0) STOP
209#endif             
210
211   CALL get_in('modif_sic', modif_sic,  0)
212   IF(modif_sic >= 1) &
213   CALL get_in('deltasic',  deltasic, 0.1)
214
215   CALL get_in('albedo_prescrit', albedo_prescrit, 0)
216   IF(albedo_prescrit == 1) THEN
217      CALL get_in('lon_min_albedo', lon_min_albedo, -200.)
218      CALL get_in('lon_max_albedo', lon_max_albedo,  200.)
219      CALL get_in('lat_min_albedo', lat_min_albedo, -100.)
220      CALL get_in('lat_max_albedo', lat_max_albedo,  100.)
221   END IF
222   deltaO18_oce=0.0
223   CALL get_in('deltaP_BL',           deltaP_BL,     10.0)
224   CALL get_in('ruissellement_pluie', ruissellement_pluie, 0)
225   CALL get_in('alphak_stewart',      alphak_stewart,      1)
226   CALL get_in('tdifexp_sol',         tdifexp_sol,      0.67)
227   CALL get_in('calendrier_guide',    calendrier_guide,    0)
228   CALL get_in('cste_surf_cond',      cste_surf_cond,      0)
229   CALL get_in('mixlen',              mixlen,           35.0)
230   CALL get_in('evap_cont_cste',      evap_cont_cste,      0)
231   CALL get_in('deltaO18_evap_cont',  deltaO18_evap_cont,0.0)
232   CALL get_in('d_evap_cont',         d_evap_cont,       0.0)
233   CALL get_in('nudge_qsol',          nudge_qsol,          0)
234   CALL get_in('region_nudge_qsol',   region_nudge_qsol,   1)
235   nlevmaxO17 = 50
236   CALL msg('nlevmaxO17='//TRIM(int2str(nlevmaxO17)))
237   CALL get_in('no_pce',   no_pce,     0)
238   CALL get_in('A_satlim', A_satlim, 1.0)
239   CALL get_in('ok_restrict_A_satlim', ok_restrict_A_satlim, 0)
240#ifdef ISOVERIF
241   CALL msg(' 315: A_satlim='//real2str(A_satlim), modname, A_satlim > 1.0)
242   IF(A_satlim > 1.0) STOP
243#endif
244!  CALL get_in('slope_limiterxy',   slope_limiterxy,  2.0)
245!  CALL get_in('slope_limiterz',    slope_limiterz,   2.0)
246   CALL get_in('modif_ratqs',       modif_ratqs,        0)
247   CALL get_in('Pcrit_ratqs',       Pcrit_ratqs,    500.0)
248   CALL get_in('ratqsbasnew',       ratqsbasnew,     0.05)
249   CALL get_in('fac_modif_evaoce',  fac_modif_evaoce, 1.0)
250   CALL get_in('ok_bidouille_wake', ok_bidouille_wake,  0)
251   ! si oui, la temperature de cond est celle de l'environnement, pour eviter
252   ! bugs quand temperature dans ascendances convs est mal calculee
253   CALL get_in('cond_temp_env',        cond_temp_env,        .FALSE.)
254   IF(ANY(isoName == 'H[3]HO')) &
255   CALL get_in('ok_prod_nucl_tritium', ok_prod_nucl_tritium, .FALSE., .FALSE.)
256
257   !--------------------------------------------------------------
258   ! Parameters that do not depend on the nature of water isotopes:
259   !--------------------------------------------------------------
260   ! -- temperature at which ice condensate starts to form (valeur ECHAM?):
261   pxtmelt = 273.15
262
263   ! -- temperature at which all condensate is ice:
264   pxtice  = 273.15-10.0
265
266   !- -- test PHASE
267!   pxtmelt = 273.15 - 10.0
268!   pxtice  = 273.15 - 30.0
269
270   ! -- minimum temperature to calculate fractionation coeff
271   pxtmin = 273.15 - 120.0   ! On ne calcule qu'au dessus de -120°C
272   pxtmax = 273.15 +  60.0   ! On ne calcule qu'au dessus de +60°C
273   !    Remarque: les coeffs ont ete mesures seulement jusq'à -40!
274
275   ! -- a constant for alpha_eff for equilibrium below cloud base:
276   tdifexp = 0.58
277   tv0cin  = 7.0
278
279   ! facteurs lambda et mu dans Si=musi-lambda*T
280   musi=1.0
281   if (ok_nocinsat) lambda_sursat = 0.0          ! no sursaturation effect
282
283   ! diffusion dans le sol
284   Kd=2.5e-9 ! m2/s   
285
286   ! cas où cste_surf_cond: on met rhs ou/et Ts cste pour voir
287   rh_cste_surf_cond = 0.6
288    T_cste_surf_cond = 288.0
289   
290   CALL msg('iso_O18, iso_HDO, iso_eau = '//TRIM(strStack(int2str([iso_O18, iso_HDO, iso_eau]))), modname)
291
292   !--------------------------------------------------------------
293   ! Parameters that depend on the nature of water isotopes:
294   !--------------------------------------------------------------
295   IF(getKey('tnat',    tnat,    isoName)) CALL abort_physic(modname, 'can''t get tnat',    1)
296   IF(getKey('toce',    toce,    isoName)) CALL abort_physic(modname, 'can''t get toce',    1)
297   IF(getKey('tcorr',   tcorr,   isoName)) CALL abort_physic(modname, 'can''t get tcorr',   1)
298   IF(getKey('talph1',  talph1,  isoName)) CALL abort_physic(modname, 'can''t get talph1',  1)
299   IF(getKey('talph2',  talph2,  isoName)) CALL abort_physic(modname, 'can''t get talph2',  1)
300   IF(getKey('talph3',  talph3,  isoName)) CALL abort_physic(modname, 'can''t get talph3',  1)
301   IF(getKey('talps1',  talps1,  isoName)) CALL abort_physic(modname, 'can''t get talps1',  1)
302   IF(getKey('talps2',  talps2,  isoName)) CALL abort_physic(modname, 'can''t get talps2',  1)
303   IF(getKey('tkcin0',  tkcin0,  isoName)) CALL abort_physic(modname, 'can''t get tkcin0',  1)
304   IF(getKey('tkcin1',  tkcin1,  isoName)) CALL abort_physic(modname, 'can''t get tkcin1',  1)
305   IF(getKey('tkcin2',  tkcin2,  isoName)) CALL abort_physic(modname, 'can''t get tkcin2',  1)
306   IF(getKey('tdifrel', tdifrel, isoName)) CALL abort_physic(modname, 'can''t get tdifrel', 1)
307   IF(getKey('alpha_liq_sol', alpha_liq_sol, isoName)) CALL abort_physic(modname, 'can''t get alpha_liq_sol',  1)
308   IF(getKey('Rdefault',Rdefault,isoName)) CALL abort_physic(modname, 'can''t get Rdefault',1)
309   IF(getKey('Rmethox', Rmethox, isoName)) CALL abort_physic(modname, 'can''t get Rmethox', 1)
310   IF(.NOT.Rdefault_smow) Rdefault(:) = 0.0
311
312   !--- Sensitivity test: no kinetic effect in sfc evaporation
313   IF(ok_nocinsfc) THEN
314      tkcin0(1:niso) = 0.0
315      tkcin1(1:niso) = 0.0
316      tkcin2(1:niso) = 0.0
317   END IF
318
319   CALL msg('285: verif initialisation:', modname)
320   DO ixt=1,niso
321      sxt=int2str(ixt)
322      CALL msg(' * isoName('//TRIM(sxt)//') = <'//TRIM(isoName(ixt))//'>',  modname)
323      CALL msg(  '    tnat('//TRIM(sxt)//') = '//TRIM(real2str(tnat(ixt))), modname)
324!     CALL msg('    alpha_liq_sol('//TRIM(sxt)//') = '//TRIM(real2str(alpha_liq_sol(ixt))), modname)
325!     CALL msg(        '   tkcin0('//TRIM(sxt)//') = '//TRIM(real2str(tkcin0(ixt))),        modname)
326!     CALL msg(       '   tdifrel('//TRIM(sxt)//') = '//TRIM(real2str(tdifrel(ixt))),       modname)
327   END DO
328   CALL msg('69:     lambda = '//TRIM(real2str(lambda_sursat)), modname)
329   CALL msg('69:    thumxt1 = '//TRIM(real2str(thumxt1)),       modname)
330   CALL msg('69: h_land_ice = '//TRIM(real2str(h_land_ice)),    modname)
331   CALL msg('69:      P_veg = '//TRIM(real2str(P_veg)),         modname)
332
333END SUBROUTINE iso_init
334
335
336SUBROUTINE getinp_s(nam, val, def, lDisp)
337   USE ioipsl_getincom, ONLY: getin
338   USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
339   USE mod_phys_lmdz_omp_data, ONLY :  is_omp_root
340   USE mod_phys_lmdz_transfert_para, ONLY : bcast
341   CHARACTER(LEN=*),  INTENT(IN)    :: nam
342   CHARACTER(LEN=*),  INTENT(INOUT) :: val
343   CHARACTER(LEN=*),  INTENT(IN)    :: def
344   LOGICAL, OPTIONAL, INTENT(IN)    :: lDisp
345   LOGICAL :: lD
346!$OMP BARRIER
347   IF(.NOT.(is_mpi_root.AND.is_omp_root)) RETURN
348   val=def; CALL getin(nam,val); CALL bcast(val)
349   lD=.TRUE.; IF(PRESENT(lDisp)) lD=lDisp
350   IF(lD) CALL msg(TRIM(nam)//' = '//TRIM(val))
351END SUBROUTINE getinp_s
352
353SUBROUTINE getinp_i(nam, val, def, lDisp)
354   USE ioipsl_getincom, ONLY: getin
355   USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
356   USE mod_phys_lmdz_omp_data, ONLY :  is_omp_root
357   USE mod_phys_lmdz_transfert_para, ONLY : bcast
358   CHARACTER(LEN=*),  INTENT(IN)    :: nam
359   INTEGER,           INTENT(INOUT) :: val
360   INTEGER,           INTENT(IN)    :: def
361   LOGICAL, OPTIONAL, INTENT(IN)    :: lDisp
362   LOGICAL :: lD
363!$OMP BARRIER
364   IF(.NOT.(is_mpi_root.AND.is_omp_root)) RETURN
365   val=def; CALL getin(nam,val); CALL bcast(val)
366   lD=.TRUE.; IF(PRESENT(lDisp)) lD=lDisp
367   IF(lD) CALL msg(TRIM(nam)//' = '//TRIM(int2str(val)))
368END SUBROUTINE getinp_i
369
370SUBROUTINE getinp_r(nam, val, def, lDisp)
371   USE ioipsl_getincom, ONLY: getin
372   USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
373   USE mod_phys_lmdz_omp_data, ONLY :  is_omp_root
374   USE mod_phys_lmdz_transfert_para, ONLY : bcast
375   CHARACTER(LEN=*),  INTENT(IN)    :: nam
376   REAL,              INTENT(INOUT) :: val
377   REAL,              INTENT(IN)    :: def
378   LOGICAL, OPTIONAL, INTENT(IN)    :: lDisp
379   LOGICAL :: lD
380!$OMP BARRIER
381   IF(.NOT.(is_mpi_root.AND.is_omp_root)) RETURN
382   val=def; CALL getin(nam,val); CALL bcast(val)
383   lD=.TRUE.; IF(PRESENT(lDisp)) lD=lDisp
384   IF(lD) CALL msg(TRIM(nam)//' = '//TRIM(real2str(val)))
385END SUBROUTINE getinp_r
386
387SUBROUTINE getinp_l(nam, val, def, lDisp)
388   USE ioipsl_getincom, ONLY: getin
389   USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
390   USE mod_phys_lmdz_omp_data, ONLY :  is_omp_root
391   USE mod_phys_lmdz_transfert_para, ONLY : bcast
392   CHARACTER(LEN=*),  INTENT(IN)    :: nam
393   LOGICAL,           INTENT(INOUT) :: val
394   LOGICAL,           INTENT(IN)    :: def
395   LOGICAL, OPTIONAL, INTENT(IN)    :: lDisp
396   LOGICAL :: lD
397!$OMP BARRIER
398   IF(.NOT.(is_mpi_root.AND.is_omp_root)) RETURN
399   val=def; CALL getin(nam,val); CALL bcast(val)
400   lD=.TRUE.; IF(PRESENT(lDisp)) lD=lDisp
401   IF(lD) CALL msg(TRIM(nam)//' = '//TRIM(bool2str(val)))
402END SUBROUTINE getinp_l
403
404END MODULE isotopes_mod
405#endif
406
407
Note: See TracBrowser for help on using the repository browser.