source: branches/ORCHIDEE_2_2/ORCHIDEE/src_parameters/constantes_var.f90 @ 7547

Last change on this file since 7547 was 7547, checked in by agnes.ducharne, 2 years ago

Create new flag IMPOSE_SLOPE to control the possibility to impose a uniform reinf_slope for surface runoff reinfiltration. orchidee.default is manually updated.

  • Property svn:keywords set to Date Revision
File size: 56.9 KB
Line 
1! =================================================================================================================================
2! MODULE       : constantes_var
3!
4! CONTACT      : orchidee-help _at_ listes.ipsl.fr
5!
6! LICENCE      : IPSL (2006)
7! This software is governed by the CeCILL licence see ORCHIDEE/ORCHIDEE_CeCILL.LIC
8!
9!>\BRIEF        constantes_var module contains most constantes like pi, Earth radius, etc...
10!!              and all externalized parameters except pft-dependent constants.
11!!
12!!\n DESCRIPTION: This module contains most constantes and the externalized parameters of ORCHIDEE which
13!!                are not pft-dependent.\n
14!!                In this module, you can set the flag diag_qsat in order to detect the pixel where the
15!!                temperature is out of range (look qsatcalc and dev_qsatcalc in qsat_moisture.f90).\n
16!!                The Earth radius is approximated by the Equatorial radius.The Earth's equatorial radius a,
17!!                or semi-major axis, is the distance from its center to the equator and equals 6,378.1370 km.
18!!                The equatorial radius is often used to compare Earth with other planets.\n
19!!                The meridional mean is well approximated by the semicubic mean of the two axe yielding
20!!                6367.4491 km or less accurately by the quadratic mean of the two axes about 6,367.454 km
21!!                or even just the mean of the two axes about 6,367.445 km.\n
22!!                This module is already USE in module constantes. Therefor no need to USE it seperatly except
23!!                if the subroutines in module constantes are not needed.\n
24!!               
25!! RECENT CHANGE(S):
26!!
27!! REFERENCE(S) :
28!! - Louis, Jean-Francois (1979), A parametric model of vertical eddy fluxes in the atmosphere.
29!! Boundary Layer Meteorology, 187-202.\n
30!!
31!! SVN          :
32!! $HeadURL: $
33!! $Date$
34!! $Revision$
35!! \n
36!_ ================================================================================================================================
37
38MODULE constantes_var
39
40  USE defprec
41
42  IMPLICIT NONE
43!-
44
45                         !-----------------------!
46                         !  ORCHIDEE CONSTANTS   !
47                         !-----------------------!
48
49  !
50  ! FLAGS
51  !
52  LOGICAL :: river_routing      !! activate river routing
53!$OMP THREADPRIVATE(river_routing)
54  LOGICAL, SAVE :: ok_nudge_mc  !! Activate nudging of soil moisture
55!$OMP THREADPRIVATE(ok_nudge_mc)
56  LOGICAL, SAVE :: ok_nudge_snow!! Activate nudging of snow variables
57!$OMP THREADPRIVATE(ok_nudge_snow)
58  LOGICAL, SAVE :: nudge_interpol_with_xios  !! Activate reading and interpolation with XIOS for nudging fields
59!$OMP THREADPRIVATE(nudge_interpol_with_xios)
60  LOGICAL :: do_floodplains     !! activate flood plains
61!$OMP THREADPRIVATE(do_floodplains)
62  LOGICAL :: do_irrigation      !! activate computation of irrigation flux
63!$OMP THREADPRIVATE(do_irrigation)
64  LOGICAL :: ok_sechiba         !! activate physic of the model
65!$OMP THREADPRIVATE(ok_sechiba)
66  LOGICAL :: ok_stomate         !! activate carbon cycle
67!$OMP THREADPRIVATE(ok_stomate)
68  LOGICAL :: ok_dgvm            !! activate dynamic vegetation
69!$OMP THREADPRIVATE(ok_dgvm)
70  LOGICAL :: do_wood_harvest    !! activate wood harvest
71!$OMP THREADPRIVATE(do_wood_harvest)
72  LOGICAL :: ok_pheno           !! activate the calculation of lai using stomate rather than a prescription
73!$OMP THREADPRIVATE(ok_pheno)
74  LOGICAL :: ok_bvoc            !! activate biogenic volatile organic coumpounds
75!$OMP THREADPRIVATE(ok_bvoc)
76  LOGICAL :: ok_leafage         !! activate leafage
77!$OMP THREADPRIVATE(ok_leafage)
78  LOGICAL :: ok_radcanopy       !! use canopy radiative transfer model
79!$OMP THREADPRIVATE(ok_radcanopy)
80  LOGICAL :: ok_multilayer      !! use canopy radiative transfer model with multi-layers
81!$OMP THREADPRIVATE(ok_multilayer)
82  LOGICAL :: ok_pulse_NOx       !! calculate NOx emissions with pulse
83!$OMP THREADPRIVATE(ok_pulse_NOx)
84  LOGICAL :: ok_bbgfertil_NOx   !! calculate NOx emissions with bbg fertilizing effect
85!$OMP THREADPRIVATE(ok_bbgfertil_NOx)
86  LOGICAL :: ok_cropsfertil_NOx !! calculate NOx emissions with fertilizers use
87!$OMP THREADPRIVATE(ok_cropsfertil_NOx)
88
89  LOGICAL :: ok_co2bvoc_poss    !! CO2 inhibition on isoprene activated following Possell et al. (2005) model
90!$OMP THREADPRIVATE(ok_co2bvoc_poss)
91  LOGICAL :: ok_co2bvoc_wilk    !! CO2 inhibition on isoprene activated following Wilkinson et al. (2006) model
92!$OMP THREADPRIVATE(ok_co2bvoc_wilk)
93 
94  LOGICAL, SAVE :: OFF_LINE_MODE = .FALSE.  !! ORCHIDEE detects if it is coupled with a GCM or
95                                            !! just use with one driver in OFF-LINE. (true/false)
96!$OMP THREADPRIVATE(OFF_LINE_MODE) 
97  LOGICAL, SAVE :: impose_param = .TRUE.    !! Flag impos_param : read all the parameters in the run.def file
98!$OMP THREADPRIVATE(impose_param)
99  CHARACTER(LEN=80), SAVE     :: restname_in       = 'NONE'                 !! Input Restart files name for Sechiba component 
100!$OMP THREADPRIVATE(restname_in)
101  CHARACTER(LEN=80), SAVE     :: restname_out      = 'sechiba_rest_out.nc'  !! Output Restart files name for Sechiba component
102!$OMP THREADPRIVATE(restname_out)
103  CHARACTER(LEN=80), SAVE     :: stom_restname_in  = 'NONE'                 !! Input Restart files name for Stomate component
104!$OMP THREADPRIVATE(stom_restname_in)
105  CHARACTER(LEN=80), SAVE     :: stom_restname_out = 'stomate_rest_out.nc'  !! Output Restart files name for Stomate component
106!$OMP THREADPRIVATE(stom_restname_out)
107  INTEGER, SAVE :: printlev=2       !! Standard level for text output [0, 1, 2, 3]
108!$OMP THREADPRIVATE(printlev)
109
110  !
111  ! TIME
112  !
113  INTEGER(i_std), PARAMETER  :: spring_days_max = 40  !! Maximum number of days during which we watch for possible spring frost damage
114  !
115  ! SPECIAL VALUES
116  !
117  INTEGER(i_std), PARAMETER :: undef_int = 999999999     !! undef integer for integer arrays (unitless)
118  !-
119  REAL(r_std), SAVE :: val_exp = 999999.                 !! Specific value if no restart value  (unitless)
120!$OMP THREADPRIVATE(val_exp)
121  REAL(r_std), PARAMETER :: undef = -9999.               !! Special value for stomate (unitless)
122 
123  REAL(r_std), PARAMETER :: min_sechiba = 1.E-8_r_std    !! Epsilon to detect a near zero floating point (unitless)
124  REAL(r_std), PARAMETER :: undef_sechiba = 1.E+20_r_std !! The undef value used in SECHIBA (unitless)
125 
126  REAL(r_std), PARAMETER :: min_stomate = 1.E-8_r_std    !! Epsilon to detect a near zero floating point (unitless)
127  REAL(r_std), PARAMETER :: large_value = 1.E33_r_std    !! some large value (for stomate) (unitless)
128
129  REAL(r_std), SAVE :: alpha_nudge_mc                    !! Nudging constant for soil moisture
130!$OMP THREADPRIVATE(alpha_nudge_mc)
131  REAL(r_std), SAVE :: alpha_nudge_snow                  !! Nudging constant for snow variables
132!$OMP THREADPRIVATE(alpha_nudge_snow)
133
134  !
135  !  DIMENSIONING AND INDICES PARAMETERS 
136  !
137  INTEGER(i_std), PARAMETER :: ibare_sechiba = 1 !! Index for bare soil in Sechiba (unitless)
138  INTEGER(i_std), PARAMETER :: ivis = 1          !! index for albedo in visible range (unitless)
139  INTEGER(i_std), PARAMETER :: inir = 2          !! index for albeod i near-infrared range (unitless)
140  INTEGER(i_std), PARAMETER :: nnobio = 1        !! Number of other surface types: land ice (lakes,cities, ...) (unitless)
141  INTEGER(i_std), PARAMETER :: iice = 1          !! Index for land ice (see nnobio) (unitless)
142  !-
143  !! Soil
144  INTEGER(i_std), PARAMETER :: classnb = 9       !! Levels of soil colour classification (unitless)
145  !-
146  INTEGER(i_std), PARAMETER :: nleafages = 4     !! leaf age discretisation ( 1 = no discretisation )(unitless)
147  !-
148  !! litter fractions: indices (unitless)
149  INTEGER(i_std), PARAMETER :: ileaf = 1         !! Index for leaf compartment (unitless)
150  INTEGER(i_std), PARAMETER :: isapabove = 2     !! Index for sapwood above compartment (unitless)
151  INTEGER(i_std), PARAMETER :: isapbelow = 3     !! Index for sapwood below compartment (unitless)
152  INTEGER(i_std), PARAMETER :: iheartabove = 4   !! Index for heartwood above compartment (unitless)
153  INTEGER(i_std), PARAMETER :: iheartbelow = 5   !! Index for heartwood below compartment (unitless)
154  INTEGER(i_std), PARAMETER :: iroot = 6         !! Index for roots compartment (unitless)
155  INTEGER(i_std), PARAMETER :: ifruit = 7        !! Index for fruits compartment (unitless)
156  INTEGER(i_std), PARAMETER :: icarbres = 8      !! Index for reserve compartment (unitless)
157  INTEGER(i_std), PARAMETER :: nparts = 8        !! Number of biomass compartments (unitless)
158  !-
159  !! indices for assimilation parameters
160  INTEGER(i_std), PARAMETER :: ivcmax = 1        !! Index for vcmax (assimilation parameters) (unitless)
161  INTEGER(i_std), PARAMETER :: npco2 = 1         !! Number of assimilation parameters (unitless)
162  !-
163  !! trees and litter: indices for the parts of heart-
164  !! and sapwood above and below the ground
165  INTEGER(i_std), PARAMETER :: iabove = 1       !! Index for above part (unitless)
166  INTEGER(i_std), PARAMETER :: ibelow = 2       !! Index for below part (unitless)
167  INTEGER(i_std), PARAMETER :: nlevs = 2        !! Number of levels for trees and litter (unitless)
168  !-
169  !! litter: indices for metabolic and structural part
170  INTEGER(i_std), PARAMETER :: imetabolic = 1   !! Index for metabolic litter (unitless)
171  INTEGER(i_std), PARAMETER :: istructural = 2  !! Index for structural litter (unitless)
172  INTEGER(i_std), PARAMETER :: nlitt = 2        !! Number of levels for litter compartments (unitless)
173  !-
174  !! carbon pools: indices
175  INTEGER(i_std), PARAMETER :: iactive = 1      !! Index for active carbon pool (unitless)
176  INTEGER(i_std), PARAMETER :: islow = 2        !! Index for slow carbon pool (unitless)
177  INTEGER(i_std), PARAMETER :: ipassive = 3     !! Index for passive carbon pool (unitless)
178  INTEGER(i_std), PARAMETER :: ncarb = 3        !! Number of soil carbon pools (unitless)
179  !-
180  !! For isotopes and nitrogen
181  INTEGER(i_std), PARAMETER :: nelements = 1    !! Number of isotopes considered
182  INTEGER(i_std), PARAMETER :: icarbon = 1      !! Index for carbon
183  !
184  !! Indices used for analytical spin-up
185  INTEGER(i_std), PARAMETER :: nbpools = 7              !! Total number of carbon pools (unitless)
186  INTEGER(i_std), PARAMETER :: istructural_above = 1    !! Index for structural litter above (unitless)
187  INTEGER(i_std), PARAMETER :: istructural_below = 2    !! Index for structural litter below (unitless)
188  INTEGER(i_std), PARAMETER :: imetabolic_above = 3     !! Index for metabolic litter above (unitless)
189  INTEGER(i_std), PARAMETER :: imetabolic_below = 4     !! Index for metabolic litter below (unitless)
190  INTEGER(i_std), PARAMETER :: iactive_pool = 5         !! Index for active carbon pool (unitless)
191  INTEGER(i_std), PARAMETER :: islow_pool   = 6         !! Index for slow carbon pool (unitless)
192  INTEGER(i_std), PARAMETER :: ipassive_pool = 7        !! Index for passive carbon pool (unitless)
193  !
194  !! Indicies used for output variables on Landuse tiles defined according to LUMIP project
195  !! Note that ORCHIDEE do not represent pasture and urban land. Therefor the variables will have
196  !! val_exp as missing value for these tiles.
197  INTEGER(i_std), PARAMETER :: nlut=4                   !! Total number of landuse tiles according to LUMIP
198  INTEGER(i_std), PARAMETER :: id_psl=1                 !! Index for primary and secondary land
199  INTEGER(i_std), PARAMETER :: id_crp=2                 !! Index for crop land
200  INTEGER(i_std), PARAMETER :: id_pst=3                 !! Index for pasture land
201  INTEGER(i_std), PARAMETER :: id_urb=4                 !! Index for urban land
202
203
204  !
205  ! NUMERICAL AND PHYSICS CONSTANTS
206  !
207  !
208
209  !-
210  ! 1. Mathematical and numerical constants
211  !-
212  REAL(r_std), PARAMETER :: pi = 3.141592653589793238   !! pi souce : http://mathworld.wolfram.com/Pi.html (unitless)
213  REAL(r_std), PARAMETER :: euler = 2.71828182845904523 !! e source : http://mathworld.wolfram.com/e.html (unitless)
214  REAL(r_std), PARAMETER :: zero = 0._r_std             !! Numerical constant set to 0 (unitless)
215  REAL(r_std), PARAMETER :: undemi = 0.5_r_std          !! Numerical constant set to 1/2 (unitless)
216  REAL(r_std), PARAMETER :: un = 1._r_std               !! Numerical constant set to 1 (unitless)
217  REAL(r_std), PARAMETER :: moins_un = -1._r_std        !! Numerical constant set to -1 (unitless)
218  REAL(r_std), PARAMETER :: deux = 2._r_std             !! Numerical constant set to 2 (unitless)
219  REAL(r_std), PARAMETER :: trois = 3._r_std            !! Numerical constant set to 3 (unitless)
220  REAL(r_std), PARAMETER :: quatre = 4._r_std           !! Numerical constant set to 4 (unitless)
221  REAL(r_std), PARAMETER :: cinq = 5._r_std             !![DISPENSABLE] Numerical constant set to 5 (unitless)
222  REAL(r_std), PARAMETER :: six = 6._r_std              !![DISPENSABLE] Numerical constant set to 6 (unitless)
223  REAL(r_std), PARAMETER :: huit = 8._r_std             !! Numerical constant set to 8 (unitless)
224  REAL(r_std), PARAMETER :: mille = 1000._r_std         !! Numerical constant set to 1000 (unitless)
225
226  !-
227  ! 2 . Physics
228  !-
229  REAL(r_std), PARAMETER :: R_Earth = 6378000.              !! radius of the Earth : Earth radius ~= Equatorial radius (m)
230  REAL(r_std), PARAMETER :: mincos  = 0.0001                !! Minimum cosine value used for interpolation (unitless)
231  REAL(r_std), PARAMETER :: pb_std = 1013.                  !! standard pressure (hPa)
232  REAL(r_std), PARAMETER :: ZeroCelsius = 273.15            !! 0 degre Celsius in degre Kelvin (K)
233  REAL(r_std), PARAMETER :: tp_00 = 273.15                  !! 0 degre Celsius in degre Kelvin (K)
234  REAL(r_std), PARAMETER :: chalsu0 = 2.8345E06             !! Latent heat of sublimation (J.kg^{-1})
235  REAL(r_std), PARAMETER :: chalev0 = 2.5008E06             !! Latent heat of evaporation (J.kg^{-1})
236  REAL(r_std), PARAMETER :: chalfu0 = chalsu0-chalev0       !! Latent heat of fusion (J.kg^{-1})
237  REAL(r_std), PARAMETER :: c_stefan = 5.6697E-8            !! Stefan-Boltzman constant (W.m^{-2}.K^{-4})
238  REAL(r_std), PARAMETER :: cp_air = 1004.675               !! Specific heat of dry air (J.kg^{-1}.K^{-1})
239  REAL(r_std), PARAMETER :: cte_molr = 287.05               !! Specific constant of dry air (kg.mol^{-1})
240  REAL(r_std), PARAMETER :: kappa = cte_molr/cp_air         !! Kappa : ratio between specific constant and specific heat
241                                                            !! of dry air (unitless)
242  REAL(r_std), PARAMETER :: msmlr_air = 28.964E-03          !! Molecular weight of dry air (kg.mol^{-1})
243  REAL(r_std), PARAMETER :: msmlr_h2o = 18.02E-03           !! Molecular weight of water vapor (kg.mol^{-1})
244  REAL(r_std), PARAMETER :: cp_h2o = &                      !! Specific heat of water vapor (J.kg^{-1}.K^{-1})
245       & cp_air*(quatre*msmlr_air)/( 3.5_r_std*msmlr_h2o) 
246  REAL(r_std), PARAMETER :: cte_molr_h2o = cte_molr/quatre  !! Specific constant of water vapor (J.kg^{-1}.K^{-1})
247  REAL(r_std), PARAMETER :: retv = msmlr_air/msmlr_h2o-un   !! Ratio between molecular weight of dry air and water
248                                                            !! vapor minus 1(unitless) 
249  REAL(r_std), PARAMETER :: rvtmp2 = cp_h2o/cp_air-un       !! Ratio between specific heat of water vapor and dry air
250                                                            !! minus 1 (unitless)
251  REAL(r_std), PARAMETER :: cepdu2 = (0.1_r_std)**2         !! Squared wind shear (m^2.s^{-2})
252  REAL(r_std), PARAMETER :: ct_karman = 0.41_r_std          !! Van Karmann Constant (unitless)
253  REAL(r_std), PARAMETER :: cte_grav = 9.80665_r_std        !! Acceleration of the gravity (m.s^{-2})
254  REAL(r_std), PARAMETER :: pa_par_hpa = 100._r_std         !! Transform pascal into hectopascal (unitless)
255  REAL(r_std), PARAMETER :: RR = 8.314                      !! Ideal gas constant (J.mol^{-1}.K^{-1})
256  REAL(r_std), PARAMETER :: Sct = 1370.                     !! Solar constant (W.m^{-2})
257
258
259  INTEGER(i_std), SAVE :: testpft = 6
260!$OMP THREADPRIVATE(testpft)
261  !-
262  ! 3. Climatic constants
263  !-
264  !! Constantes of the Louis scheme
265  REAL(r_std), SAVE :: cb = 5._r_std              !! Constant of the Louis scheme (unitless);
266                                                  !! reference to Louis (1979)
267!$OMP THREADPRIVATE(cb)
268  REAL(r_std), SAVE :: cc = 5._r_std              !! Constant of the Louis scheme (unitless);
269                                                  !! reference to Louis (1979)
270!$OMP THREADPRIVATE(cc)
271  REAL(r_std), SAVE :: cd = 5._r_std              !! Constant of the Louis scheme (unitless);
272                                                  !! reference to Louis (1979)
273!$OMP THREADPRIVATE(cd)
274  REAL(r_std), SAVE :: rayt_cste = 125.           !! Constant in the computation of surface resistance (W.m^{-2})
275!$OMP THREADPRIVATE(rayt_cste)
276  REAL(r_std), SAVE :: defc_plus = 23.E-3         !! Constant in the computation of surface resistance (K.W^{-1})
277!$OMP THREADPRIVATE(defc_plus)
278  REAL(r_std), SAVE :: defc_mult = 1.5            !! Constant in the computation of surface resistance (K.W^{-1})
279!$OMP THREADPRIVATE(defc_mult)
280
281  !-
282  ! 4. Soil thermodynamics constants
283  !-
284  ! Look at constantes_soil.f90
285
286
287  !
288  ! OPTIONAL PARTS OF THE MODEL
289  !
290  LOGICAL,PARAMETER :: diag_qsat = .TRUE.         !! One of the most frequent problems is a temperature out of range
291                                                  !! we provide here a way to catch that in the calling procedure.
292                                                  !! (from Jan Polcher)(true/false)
293  LOGICAL, SAVE     :: almaoutput =.FALSE.        !! Selects the type of output for the model.(true/false)
294                                                  !! Value is read from run.def in intersurf_history
295!$OMP THREADPRIVATE(almaoutput)
296
297  !
298  ! DIVERSE
299  !
300  CHARACTER(LEN=100), SAVE :: stomate_forcing_name='NONE'  !! NV080800 Name of STOMATE forcing file (unitless)
301                                                           ! Compatibility with Nicolas Viovy driver.
302!$OMP THREADPRIVATE(stomate_forcing_name)
303  CHARACTER(LEN=100), SAVE :: stomate_Cforcing_name='NONE' !! NV080800 Name of soil forcing file (unitless)
304                                                           ! Compatibility with Nicolas Viovy driver.
305!$OMP THREADPRIVATE(stomate_Cforcing_name)
306  INTEGER(i_std), SAVE :: forcing_id                 !! Index of the forcing file (unitless)
307!$OMP THREADPRIVATE(forcing_id)
308  LOGICAL, SAVE :: allow_forcing_write=.TRUE.        !! Allow writing of stomate_forcing file.
309                                                     !! This variable will be set to false for teststomate.
310!$OMP THREADPRIVATE(allow_forcing_write)
311
312
313
314                         !------------------------!
315                         !  SECHIBA PARAMETERS    !
316                         !------------------------!
317 
318
319  !
320  ! GLOBAL PARAMETERS   
321  !
322  REAL(r_std), SAVE :: min_wind = 0.1      !! The minimum wind (m.s^{-1})
323!$OMP THREADPRIVATE(min_wind)
324  REAL(r_std), PARAMETER :: min_qc = 1.e-4 !! The minimum value for qc (qc=drag*wind) used in coupled(enerbil) and forced mode (enerbil and diffuco)
325  REAL(r_std), SAVE :: snowcri = 1.5       !! Sets the amount above which only sublimation occures (kg.m^{-2})
326!$OMP THREADPRIVATE(snowcri)
327
328
329  !
330  ! FLAGS ACTIVATING SUB-MODELS
331  !
332  LOGICAL, SAVE :: treat_expansion = .FALSE.   !! Do we treat PFT expansion across a grid point after introduction? (true/false)
333!$OMP THREADPRIVATE(treat_expansion)
334  LOGICAL, SAVE :: ok_herbivores = .FALSE.     !! flag to activate herbivores (true/false)
335!$OMP THREADPRIVATE(ok_herbivores)
336  LOGICAL, SAVE :: harvest_agri = .TRUE.       !! flag to harvest aboveground biomass from agricultural PFTs)(true/false)
337!$OMP THREADPRIVATE(harvest_agri)
338  LOGICAL, SAVE :: lpj_gap_const_mort          !! constant moratlity (true/false). Default value depend on OK_DGVM.
339!$OMP THREADPRIVATE(lpj_gap_const_mort)
340  LOGICAL, SAVE :: disable_fire = .TRUE.       !! flag that disable fire (true/false)
341!$OMP THREADPRIVATE(disable_fire)
342  LOGICAL, SAVE :: spinup_analytic = .FALSE.   !! Flag to activate analytical resolution for spinup (true/false)
343!$OMP THREADPRIVATE(spinup_analytic)
344
345  !
346  ! CONFIGURATION VEGETATION
347  !
348  LOGICAL, SAVE :: agriculture = .TRUE.    !! allow agricultural PFTs (true/false)
349!$OMP THREADPRIVATE(agriculture)
350  LOGICAL, SAVE :: impveg = .FALSE.        !! Impose vegetation ? (true/false)
351!$OMP THREADPRIVATE(impveg)
352  LOGICAL, SAVE :: impsoilt = .FALSE.      !! Impose soil ? (true/false)
353!$OMP THREADPRIVATE(impsoilt)
354  LOGICAL, SAVE :: impslope = .FALSE.      !! Impose reinf_slope ? (true/false)
355!$OMP THREADPRIVATE(impslope)
356  LOGICAL, SAVE :: do_now_stomate_lcchange = .FALSE.  !! Time to call lcchange in stomate_lpj
357!$OMP THREADPRIVATE(do_now_stomate_lcchange)
358  LOGICAL, SAVE :: do_now_stomate_woodharvest = .FALSE.  !! Time to call woodharvest in stomate_lpj
359!$OMP THREADPRIVATE(do_now_stomate_woodharvest)
360  LOGICAL, SAVE :: done_stomate_lcchange = .FALSE.    !! If true, call lcchange in stomate_lpj has just been done.
361!$OMP THREADPRIVATE(done_stomate_lcchange)
362  LOGICAL, SAVE :: read_lai = .FALSE.      !! Flag to read a map of LAI if STOMATE is not activated (true/false)
363!$OMP THREADPRIVATE(read_lai)
364  LOGICAL, SAVE :: veget_reinit = .TRUE.   !! To change LAND USE file in a run. (true/false)
365!$OMP THREADPRIVATE(veget_reinit)
366  LOGICAL, SAVE :: vegetmap_reset = .FALSE.!! Reset the vegetation map and reset carbon related variables
367!$OMP THREADPRIVATE(vegetmap_reset)
368  INTEGER(i_std) , SAVE :: veget_update    !! Update frequency in years for landuse (nb of years)
369!$OMP THREADPRIVATE(veget_update)
370  !
371  ! PARAMETERS USED BY BOTH HYDROLOGY MODELS
372  !
373  REAL(r_std), SAVE :: max_snow_age = 50._r_std !! Maximum period of snow aging (days)
374!$OMP THREADPRIVATE(max_snow_age)
375  REAL(r_std), SAVE :: snow_trans = 0.2_r_std   !! Transformation time constant for snow (m), reduced from the value 0.3 (04/07/2016)
376!$OMP THREADPRIVATE(snow_trans)
377  REAL(r_std), SAVE :: sneige                   !! Lower limit of snow amount (kg.m^{-2})
378!$OMP THREADPRIVATE(sneige)
379  REAL(r_std), SAVE :: maxmass_snow = 3000.     !! The maximum mass of snow (kg.m^{-2})
380!$OMP THREADPRIVATE(maxmass_snow)
381
382  !! Heat capacity
383  REAL(r_std), PARAMETER :: rho_water = 1000.           !! Density of water (kg/m3)
384  REAL(r_std), PARAMETER :: rho_ice = 920.              !! Density of ice (kg/m3)
385  REAL(r_std), PARAMETER :: rho_soil = 2700.            !! Density of soil particles (kg/m3), value from Peters-Lidard et al. 1998
386
387  !! Thermal conductivities
388  REAL(r_std), PARAMETER :: cond_water = 0.6            !! Thermal conductivity of liquid water (W/m/K)
389  REAL(r_std), PARAMETER :: cond_ice = 2.2              !! Thermal conductivity of ice (W/m/K)
390  REAL(r_std), PARAMETER :: cond_solid = 2.32           !! Thermal conductivity of mineral soil particles (W/m/K)
391
392  !! Time constant of long-term soil humidity (s)
393  REAL(r_std), PARAMETER :: lhf = 0.3336*1.E6           !! Latent heat of fusion (J/kg)
394
395  INTEGER(i_std), PARAMETER :: nsnow=3                  !! Number of levels in the snow for explicit snow scheme   
396  REAL(r_std), PARAMETER    :: XMD    = 28.9644E-3 
397  REAL(r_std), PARAMETER    :: XBOLTZ      = 1.380658E-23 
398  REAL(r_std), PARAMETER    :: XAVOGADRO   = 6.0221367E+23 
399  REAL(r_std), PARAMETER    :: XRD    = XAVOGADRO * XBOLTZ / XMD 
400  REAL(r_std), PARAMETER    :: XCPD   = 7.* XRD /2. 
401  REAL(r_std), PARAMETER    :: phigeoth = 0.057 ! 0. DKtest
402  REAL(r_std), PARAMETER    :: thick_min_snow = .01 
403
404  !! The maximum snow density and water holding characterisicts
405  REAL(r_std), SAVE         :: xrhosmax = 750.  ! (kg m-3)
406!$OMP THREADPRIVATE(xrhosmax)
407  REAL(r_std), SAVE         :: xwsnowholdmax1   = 0.03  ! (-)
408!$OMP THREADPRIVATE(xwsnowholdmax1)
409  REAL(r_std), SAVE         :: xwsnowholdmax2   = 0.10  ! (-)
410!$OMP THREADPRIVATE(xwsnowholdmax2)
411  REAL(r_std), SAVE         :: xsnowrhohold     = 200.0 ! (kg/m3)
412!$OMP THREADPRIVATE(xsnowrhohold)
413  REAL(r_std), SAVE         :: xrhosmin = 50. 
414!$OMP THREADPRIVATE(xrhosmin)
415  REAL(r_std), PARAMETER    :: xci = 2.106e+3 
416  REAL(r_std), PARAMETER    :: xrv = 6.0221367e+23 * 1.380658e-23 /18.0153e-3 
417
418  !! ISBA-ES Critical snow depth at which snow grid thicknesses constant
419  REAL(r_std), PARAMETER    :: xsnowcritd = 0.03  ! (m)
420
421  !! The threshold of snow depth used for preventing numerical problem in thermal calculations
422  REAL(r_std), PARAMETER    :: snowcritd_thermal = 0.01  ! (m) 
423 
424  !! ISBA-ES CROCUS (Pahaut 1976): snowfall density coefficients:
425  REAL(r_std), PARAMETER       :: snowfall_a_sn = 109.0  !! (kg/m3)
426  REAL(r_std), PARAMETER       :: snowfall_b_sn =   6.0  !! (kg/m3/K)
427  REAL(r_std), PARAMETER       :: snowfall_c_sn =  26.0  !! [kg/(m7/2 s1/2)]
428
429  REAL(r_std), PARAMETER       :: dgrain_new_max=  2.0e-4!! (m) : Maximum grain size of new snowfall
430 
431  !! Used in explicitsnow to prevent numerical problems as snow becomes vanishingly thin.
432  REAL(r_std), PARAMETER                :: psnowdzmin = .0001   ! m
433  REAL(r_std), PARAMETER                :: xsnowdmin = .000001  ! m
434
435  REAL(r_std), PARAMETER                :: ph2o = 1000.         !! Water density [kg/m3]
436 
437  ! ISBA-ES Thermal conductivity coefficients from Anderson (1976):
438  ! see Boone, Meteo-France/CNRM Note de Centre No. 70 (2002)
439  REAL(r_std), SAVE                     :: ZSNOWTHRMCOND1 = 0.02    ! [W/m/K]
440!$OMP THREADPRIVATE(ZSNOWTHRMCOND1)
441  REAL(r_std), SAVE                     :: ZSNOWTHRMCOND2 = 2.5E-6  ! [W m5/(kg2 K)]
442!$OMP THREADPRIVATE(ZSNOWTHRMCOND2)
443 
444  ! ISBA-ES Thermal conductivity: Implicit vapor diffn effects
445  ! (sig only for new snow OR high altitudes)
446  ! from Sun et al. (1999): based on data from Jordan (1991)
447  ! see Boone, Meteo-France/CNRM Note de Centre No. 70 (2002)
448  !
449  REAL(r_std), SAVE                       :: ZSNOWTHRMCOND_AVAP  = -0.06023 ! (W/m/K)
450!$OMP THREADPRIVATE(ZSNOWTHRMCOND_AVAP)
451  REAL(r_std), SAVE                       :: ZSNOWTHRMCOND_BVAP  = -2.5425  ! (W/m)
452!$OMP THREADPRIVATE(ZSNOWTHRMCOND_BVAP)
453  REAL(r_std), SAVE                       :: ZSNOWTHRMCOND_CVAP  = -289.99  ! (K)
454!$OMP THREADPRIVATE(ZSNOWTHRMCOND_CVAP)
455 
456  REAL(r_std),SAVE :: xansmax = 0.85      !! Maxmimum snow albedo
457!$OMP THREADPRIVATE(xansmax)
458  REAL(r_std),SAVE :: xansmin = 0.50      !! Miniumum snow albedo
459!$OMP THREADPRIVATE(xansmin)
460  REAL(r_std),SAVE :: xans_todry = 0.008  !! Albedo decay rate for dry snow
461!$OMP THREADPRIVATE(xans_todry)
462  REAL(r_std),SAVE :: xans_t = 0.240      !! Albedo decay rate for wet snow
463!$OMP THREADPRIVATE(xans_t)
464
465  ! ISBA-ES Thermal conductivity coefficients from Anderson (1976):
466  ! see Boone, Meteo-France/CNRM Note de Centre No. 70 (2002)
467  REAL(r_std), PARAMETER                  :: XP00 = 1.E5
468
469  ! ISBA-ES Thermal conductivity: Implicit vapor diffn effects
470  ! (sig only for new snow OR high altitudes)
471  ! from Sun et al. (1999): based on data from Jordan (1991)
472  ! see Boone, Meteo-France/CNRM Note de Centre No. 70 (2002)
473  !
474  REAL(r_std), SAVE          :: ZSNOWCMPCT_RHOD  = 150.0        !! (kg/m3)
475!$OMP THREADPRIVATE(ZSNOWCMPCT_RHOD)
476  REAL(r_std), SAVE          :: ZSNOWCMPCT_ACM   = 2.8e-6       !! (1/s
477!$OMP THREADPRIVATE(ZSNOWCMPCT_ACM)
478  REAL(r_std), SAVE          :: ZSNOWCMPCT_BCM   = 0.04         !! (1/K)
479!$OMP THREADPRIVATE(ZSNOWCMPCT_BCM)
480  REAL(r_std), SAVE          :: ZSNOWCMPCT_CCM   = 460.         !! (m3/kg)
481!$OMP THREADPRIVATE(ZSNOWCMPCT_CCM)
482  REAL(r_std), SAVE          :: ZSNOWCMPCT_V0    = 3.7e7        !! (Pa/s)
483!$OMP THREADPRIVATE(ZSNOWCMPCT_V0)
484  REAL(r_std), SAVE          :: ZSNOWCMPCT_VT    = 0.081        !! (1/K)
485!$OMP THREADPRIVATE(ZSNOWCMPCT_VT)
486  REAL(r_std), SAVE          :: ZSNOWCMPCT_VR    = 0.018        !! (m3/kg)
487!$OMP THREADPRIVATE(ZSNOWCMPCT_VR)
488
489  !
490  ! BVOC : Biogenic activity  for each age class
491  !
492  REAL(r_std), SAVE, DIMENSION(nleafages) :: iso_activity = (/0.5, 1.5, 1.5, 0.5/)     !! Biogenic activity for each
493                                                                                       !! age class : isoprene (unitless)
494!$OMP THREADPRIVATE(iso_activity)
495  REAL(r_std), SAVE, DIMENSION(nleafages) :: methanol_activity = (/1., 1., 0.5, 0.5/)  !! Biogenic activity for each
496                                                                                       !! age class : methanol (unnitless)
497!$OMP THREADPRIVATE(methanol_activity)
498
499  !
500  ! condveg.f90
501  !
502
503  ! 1. Scalar
504
505  ! 1.1 Flags used inside the module
506
507  LOGICAL, SAVE :: alb_bare_model = .FALSE. !! Switch for choosing values of bare soil
508                                            !! albedo (see header of subroutine)
509                                            !! (true/false)
510!$OMP THREADPRIVATE(alb_bare_model)
511  LOGICAL, SAVE :: alb_bg_modis = .TRUE.    !! Switch for choosing values of bare soil
512                                            !! albedo read from file
513                                            !! (true/false)
514!$OMP THREADPRIVATE(alb_bg_modis)
515  LOGICAL, SAVE :: impaze = .FALSE.         !! Switch for choosing surface parameters
516                                            !! (see header of subroutine). 
517                                            !! (true/false)
518!$OMP THREADPRIVATE(impaze)
519  LOGICAL, SAVE :: rough_dyn = .TRUE.       !! Chooses between two methods to calculate the
520                                            !! the roughness height : static or dynamic (varying with LAI)
521                                            !! (true/false)
522!$OMP THREADPRIVATE(rough_dyn)
523
524  LOGICAL, SAVE :: new_watstress = .FALSE.
525!$OMP THREADPRIVATE(new_watstress)
526
527  REAL(r_std), SAVE :: alpha_watstress = 1.
528!$OMP THREADPRIVATE(alpha_watstress)
529
530  ! 1.2 Others
531
532
533  REAL(r_std), SAVE :: height_displacement = 0.66        !! Factor to calculate the zero-plane displacement
534                                                         !! height from vegetation height (m)
535!$OMP THREADPRIVATE(height_displacement)
536  REAL(r_std), SAVE :: z0_bare = 0.01                    !! bare soil roughness length (m)
537!$OMP THREADPRIVATE(z0_bare)
538  REAL(r_std), SAVE :: z0_ice = 0.001                    !! ice roughness length (m)
539!$OMP THREADPRIVATE(z0_ice)
540  REAL(r_std), SAVE :: tcst_snowa = 10.0                 !! Time constant of the albedo decay of snow (days), increased from the value 5.0 (04/07/2016)
541!$OMP THREADPRIVATE(tcst_snowa)
542  REAL(r_std), SAVE :: snowcri_alb = 10.                 !! Critical value for computation of snow albedo (cm)
543!$OMP THREADPRIVATE(snowcri_alb)
544  REAL(r_std), SAVE :: fixed_snow_albedo = undef_sechiba !! To choose a fixed snow albedo value (unitless)
545!$OMP THREADPRIVATE(fixed_snow_albedo)
546  REAL(r_std), SAVE :: z0_scal = 0.15                    !! Surface roughness height imposed (m)
547!$OMP THREADPRIVATE(z0_scal)
548  REAL(r_std), SAVE :: roughheight_scal = zero           !! Effective roughness Height depending on zero-plane
549                                                         !! displacement height (m) (imposed)
550!$OMP THREADPRIVATE(roughheight_scal)
551  REAL(r_std), SAVE :: emis_scal = 1.0                   !! Surface emissivity imposed (unitless)
552!$OMP THREADPRIVATE(emis_scal)
553
554  REAL(r_std), SAVE :: c1 = 0.32                         !! Constant used in the formulation of the ratio of
555!$OMP THREADPRIVATE(c1)                                  !! friction velocity to the wind speed at the canopy top
556                                                         !! see Ershadi et al. (2015) for more info
557  REAL(r_std), SAVE :: c2 = 0.264                        !! Constant used in the formulation of the ratio of
558!$OMP THREADPRIVATE(c2)                                  !! friction velocity to the wind speed at the canopy top
559                                                         !! see Ershadi et al. (2015) for more info
560  REAL(r_std), SAVE :: c3 = 15.1                         !! Constant used in the formulation of the ratio of
561!$OMP THREADPRIVATE(c3)                                  !! friction velocity to the wind speed at the canopy top
562                                                         !! see Ershadi et al. (2015) for more info
563  REAL(r_std), SAVE :: Cdrag_foliage = 0.2               !! Drag coefficient of the foliage
564!$OMP THREADPRIVATE(Cdrag_foliage)                       !! See Ershadi et al. (2015) and Su et. al (2001) for more info
565  REAL(r_std), SAVE :: Ct = 0.01                         !! Heat transfer coefficient of the leaf
566!$OMP THREADPRIVATE(Ct)                                  !! See Ershadi et al. (2015) and Su et. al (2001) for more info
567  REAL(r_std), SAVE :: Prandtl = 0.71                    !! Prandtl number used in the calculation of Ct_star
568!$OMP THREADPRIVATE(Prandtl)                             !! See Su et. al (2001) for more info
569
570
571
572  ! 2. Arrays
573
574  REAL(r_std), SAVE, DIMENSION(2) :: alb_deadleaf = (/ .12, .35/)    !! albedo of dead leaves, VIS+NIR (unitless)
575!$OMP THREADPRIVATE(alb_deadleaf)
576  REAL(r_std), SAVE, DIMENSION(2) :: alb_ice = (/ .60, .20/)         !! albedo of ice, VIS+NIR (unitless)
577!$OMP THREADPRIVATE(alb_ice)
578  REAL(r_std), SAVE, DIMENSION(2) :: albedo_scal = (/ 0.25, 0.25 /)  !! Albedo values for visible and near-infrared
579                                                                     !! used imposed (unitless)
580!$OMP THREADPRIVATE(albedo_scal)
581  REAL(r_std) , SAVE, DIMENSION(classnb) :: vis_dry = (/0.24,&
582       &0.22, 0.20, 0.18, 0.16, 0.14, 0.12, 0.10, 0.27/)  !! Soil albedo values to soil colour classification:
583                                                          !! dry soil albedo values in visible range
584!$OMP THREADPRIVATE(vis_dry)
585  REAL(r_std), SAVE, DIMENSION(classnb) :: nir_dry = (/0.48,&
586       &0.44, 0.40, 0.36, 0.32, 0.28, 0.24, 0.20, 0.55/)  !! Soil albedo values to soil colour classification:
587                                                          !! dry soil albedo values in near-infrared range
588!$OMP THREADPRIVATE(nir_dry)
589  REAL(r_std), SAVE, DIMENSION(classnb) :: vis_wet = (/0.12,&
590       &0.11, 0.10, 0.09, 0.08, 0.07, 0.06, 0.05, 0.15/)  !! Soil albedo values to soil colour classification:
591                                                          !! wet soil albedo values in visible range
592!$OMP THREADPRIVATE(vis_wet)
593  REAL(r_std), SAVE, DIMENSION(classnb) :: nir_wet = (/0.24,&
594       &0.22, 0.20, 0.18, 0.16, 0.14, 0.12, 0.10, 0.31/)  !! Soil albedo values to soil colour classification:
595                                                          !! wet soil albedo values in near-infrared range
596!$OMP THREADPRIVATE(nir_wet)
597  REAL(r_std), SAVE, DIMENSION(classnb) :: albsoil_vis = (/ &
598       &0.18, 0.16, 0.16, 0.15, 0.12, 0.105, 0.09, 0.075, 0.25/)   !! Soil albedo values to soil colour classification:
599                                                                   !! Averaged of wet and dry soil albedo values
600                                                                   !! in visible and near-infrared range
601!$OMP THREADPRIVATE(albsoil_vis)
602  REAL(r_std), SAVE, DIMENSION(classnb) :: albsoil_nir = (/ &
603       &0.36, 0.34, 0.34, 0.33, 0.30, 0.25, 0.20, 0.15, 0.45/)  !! Soil albedo values to soil colour classification:
604                                                                !! Averaged of wet and dry soil albedo values
605                                                                !! in visible and near-infrared range
606!$OMP THREADPRIVATE(albsoil_nir)
607
608  !
609  ! diffuco.f90
610  !
611
612  ! 0. Constants
613
614  REAL(r_std), PARAMETER :: Tetens_1 = 0.622         !! Ratio between molecular weight of water vapor and molecular weight 
615                                                     !! of dry air (unitless)
616  REAL(r_std), PARAMETER :: Tetens_2 = 0.378         !!
617  REAL(r_std), PARAMETER :: ratio_H2O_to_CO2 = 1.6   !! Ratio of water vapor diffusivity to the CO2 diffusivity (unitless)
618  REAL(r_std), PARAMETER :: mol_to_m_1 = 0.0244      !!
619  REAL(r_std), PARAMETER :: RG_to_PAR = 0.5          !!
620  REAL(r_std), PARAMETER :: W_to_mol = 4.6          !! W_to_mmol * RG_to_PAR = 2.3
621
622  ! 1. Scalar
623
624  INTEGER(i_std), SAVE :: nlai = 20             !! Number of LAI levels (unitless)
625!$OMP THREADPRIVATE(nlai)
626  LOGICAL, SAVE :: ldq_cdrag_from_gcm = .FALSE. !! Set to .TRUE. if you want q_cdrag coming from GCM
627!$OMP THREADPRIVATE(ldq_cdrag_from_gcm)
628  REAL(r_std), SAVE :: laimax = 12.             !! Maximal LAI used for splitting LAI into N layers (m^2.m^{-2})
629!$OMP THREADPRIVATE(laimax)
630  LOGICAL, SAVE :: downregulation_co2 = .TRUE.             !! Set to .TRUE. if you want CO2 downregulation version used for CMIP6 6.1.0-6.1.10
631!$OMP THREADPRIVATE(downregulation_co2)
632  LOGICAL, SAVE :: downregulation_co2_new = .FALSE.        !! Set to .TRUE. if you want CO2 downregulation version revised for CMIP6 6.1.11
633!$OMP THREADPRIVATE(downregulation_co2_new)
634  REAL(r_std), SAVE :: downregulation_co2_baselevel = 380. !! CO2 base level (ppm)
635!$OMP THREADPRIVATE(downregulation_co2_baselevel)
636  REAL(r_std), SAVE :: downregulation_co2_minimum = 280.   !! CO2 value above which downregulation is taken into account
637!$OMP THREADPRIVATE(downregulation_co2_minimum)
638
639  REAL(r_std), SAVE :: gb_ref = 1./25.                     !! Leaf bulk boundary layer resistance (s m-1)
640!$OMP THREADPRIVATE(gb_ref)
641
642  ! 3. Coefficients of equations
643
644  REAL(r_std), SAVE :: lai_level_depth = 0.15  !!
645!$OMP THREADPRIVATE(lai_level_depth)
646!
647  REAL(r_std), SAVE, DIMENSION(6) :: dew_veg_poly_coeff = &            !! coefficients of the 5 degree polynomomial used
648  & (/ 0.887773, 0.205673, 0.110112, 0.014843,  0.000824,  0.000017 /) !! in the equation of coeff_dew_veg
649!$OMP THREADPRIVATE(dew_veg_poly_coeff)
650!
651  REAL(r_std), SAVE               :: Oi=210000.    !! Intercellular oxygen partial pressure (ubar)
652!$OMP THREADPRIVATE(Oi)
653  !
654  ! slowproc.f90
655  !
656
657  ! 1. Scalar
658
659  INTEGER(i_std), SAVE :: veget_year_orig = 0        !!  first year for landuse (number)
660!$OMP THREADPRIVATE(veget_year_orig)
661 
662  REAL(r_std), SAVE :: min_vegfrac = 0.001           !! Minimal fraction of mesh a vegetation type can occupy (0-1, unitless)
663!$OMP THREADPRIVATE(min_vegfrac)
664  REAL(r_std), SAVE :: frac_nobio_fixed_test_1 = 0.0 !! Value for frac_nobio for tests in 0-dim simulations (0-1, unitless)
665!$OMP THREADPRIVATE(frac_nobio_fixed_test_1)
666 
667  REAL(r_std), SAVE :: stempdiag_bid = 280.          !! only needed for an initial LAI if there is no restart file
668!$OMP THREADPRIVATE(stempdiag_bid)
669
670
671                           !-----------------------------!
672                           !  STOMATE AND LPJ PARAMETERS !
673                           !-----------------------------!
674
675
676  !
677  ! lpj_constraints.f90
678  !
679 
680  ! 1. Scalar
681
682  REAL(r_std), SAVE  :: too_long = 5.      !! longest sustainable time without
683                                           !! regeneration (vernalization) (years)
684!$OMP THREADPRIVATE(too_long)
685
686
687  !
688  ! lpj_establish.f90
689  !
690
691  ! 1. Scalar
692
693  REAL(r_std), SAVE :: estab_max_tree = 0.12   !! Maximum tree establishment rate (ind/m2/dt_stomate)
694!$OMP THREADPRIVATE(estab_max_tree)
695  REAL(r_std), SAVE :: estab_max_grass = 0.12  !! Maximum grass establishment rate (ind/m2/dt_stomate)
696!$OMP THREADPRIVATE(estab_max_grass)
697 
698  ! 3. Coefficients of equations
699
700  REAL(r_std), SAVE :: establish_scal_fact = 5.  !!
701!$OMP THREADPRIVATE(establish_scal_fact)
702  REAL(r_std), SAVE :: max_tree_coverage = 0.98  !! (0-1, unitless)
703!$OMP THREADPRIVATE(max_tree_coverage)
704  REAL(r_std), SAVE :: ind_0_estab = 0.2         !! = ind_0 * 10.
705!$OMP THREADPRIVATE(ind_0_estab)
706
707
708  !
709  ! lpj_fire.f90
710  !
711
712  ! 1. Scalar
713
714  REAL(r_std), SAVE :: tau_fire = 30.           !! Time scale for memory of the fire index (days).
715!$OMP THREADPRIVATE(tau_fire)
716  REAL(r_std), SAVE :: litter_crit = 200.       !! Critical litter quantity for fire
717                                                !! below which iginitions extinguish
718                                                !! @tex $(gC m^{-2})$ @endtex
719!$OMP THREADPRIVATE(litter_crit)
720  REAL(r_std), SAVE :: fire_resist_struct = 0.5 !!
721!$OMP THREADPRIVATE(fire_resist_struct)
722  ! 2. Arrays
723
724  REAL(r_std), SAVE, DIMENSION(nparts) :: co2frac = &    !! The fraction of the different biomass
725       & (/ .95, .95, 0., 0.3, 0., 0., .95, .95 /)       !! compartments emitted to the atmosphere
726!$OMP THREADPRIVATE(co2frac)                                                         !! when burned (unitless, 0-1) 
727
728  ! 3. Coefficients of equations
729
730  REAL(r_std), SAVE, DIMENSION(3) :: bcfrac_coeff = (/ .3,  1.3,  88.2 /)         !! (unitless)
731!$OMP THREADPRIVATE(bcfrac_coeff)
732  REAL(r_std), SAVE, DIMENSION(4) :: firefrac_coeff = (/ 0.45, 0.8, 0.6, 0.13 /)  !! (unitless)
733!$OMP THREADPRIVATE(firefrac_coeff)
734
735  !
736  ! lpj_gap.f90
737  !
738
739  ! 1. Scalar
740
741  REAL(r_std), SAVE :: ref_greff = 0.035         !! Asymptotic maximum mortality rate
742                                                 !! @tex $(year^{-1})$ @endtex
743!$OMP THREADPRIVATE(ref_greff)
744
745  !               
746  ! lpj_light.f90
747  !             
748
749  ! 1. Scalar
750 
751  LOGICAL, SAVE :: annual_increase = .TRUE. !! for diagnosis of fpc increase, compare today's fpc to last year's maximum (T) or
752                                            !! to fpc of last time step (F)? (true/false)
753!$OMP THREADPRIVATE(annual_increase)
754  REAL(r_std), SAVE :: min_cover = 0.05     !! For trees, minimum fraction of crown area occupied
755                                            !! (due to its branches etc.) (0-1, unitless)
756                                            !! This means that only a small fraction of its crown area
757                                            !! can be invaded by other trees.
758!$OMP THREADPRIVATE(min_cover)
759  !
760  ! lpj_pftinout.f90
761  !
762
763  ! 1. Scalar
764
765  REAL(r_std), SAVE :: min_avail = 0.01         !! minimum availability
766!$OMP THREADPRIVATE(min_avail)
767  REAL(r_std), SAVE :: ind_0 = 0.02             !! initial density of individuals
768!$OMP THREADPRIVATE(ind_0)
769  ! 3. Coefficients of equations
770 
771  REAL(r_std), SAVE :: RIP_time_min = 1.25      !! test whether the PFT has been eliminated lately (years)
772!$OMP THREADPRIVATE(RIP_time_min)
773  REAL(r_std), SAVE :: npp_longterm_init = 10.  !! Initialisation value for npp_longterm (gC.m^{-2}.year^{-1})
774!$OMP THREADPRIVATE(npp_longterm_init)
775  REAL(r_std), SAVE :: everywhere_init = 0.05   !!
776!$OMP THREADPRIVATE(everywhere_init)
777
778
779  !
780  ! stomate_alloc.f90
781  !
782
783  ! 0. Constants
784
785  REAL(r_std), PARAMETER :: max_possible_lai = 10. !! (m^2.m^{-2})
786  REAL(r_std), PARAMETER :: Nlim_Q10 = 10.         !!
787
788  ! 1. Scalar
789
790  LOGICAL, SAVE :: ok_minres = .TRUE.              !! [DISPENSABLE] Do we try to reach a minimum reservoir even if
791                                                   !! we are severely stressed? (true/false)
792!$OMP THREADPRIVATE(ok_minres)
793  REAL(r_std), SAVE :: reserve_time_tree = 30.     !! Maximum number of days during which
794                                                   !! carbohydrate reserve may be used for
795                                                   !! trees (days)
796!$OMP THREADPRIVATE(reserve_time_tree)
797  REAL(r_std), SAVE :: reserve_time_grass = 20.    !! Maximum number of days during which
798                                                   !! carbohydrate reserve may be used for
799                                                   !! grasses (days)
800!$OMP THREADPRIVATE(reserve_time_grass)
801
802  REAL(r_std), SAVE :: f_fruit = 0.1               !! Default fruit allocation (0-1, unitless)
803!$OMP THREADPRIVATE(f_fruit)
804  REAL(r_std), SAVE :: alloc_sap_above_grass = 1.0 !! fraction of sapwood allocation above ground
805                                                   !! for grass (0-1, unitless)
806!$OMP THREADPRIVATE(alloc_sap_above_grass)
807  REAL(r_std), SAVE :: min_LtoLSR = 0.2            !! Prescribed lower bounds for leaf
808                                                   !! allocation (0-1, unitless)
809!$OMP THREADPRIVATE(min_LtoLSR)
810  REAL(r_std), SAVE :: max_LtoLSR = 0.5            !! Prescribed upper bounds for leaf
811                                                   !! allocation (0-1, unitless)
812!$OMP THREADPRIVATE(max_LtoLSR)
813  REAL(r_std), SAVE :: z_nitrogen = 0.2            !! Curvature of the root profile (m)
814!$OMP THREADPRIVATE(z_nitrogen)
815
816  ! 3. Coefficients of equations
817
818  REAL(r_std), SAVE :: Nlim_tref = 25.             !! (C)
819!$OMP THREADPRIVATE(Nlim_tref)
820
821
822  !
823  ! stomate_data.f90
824  !
825
826  ! 1. Scalar
827
828  ! 1.1 Parameters for the pipe model
829
830  REAL(r_std), SAVE :: pipe_tune1 = 100.0        !! crown area = pipe_tune1. stem diameter**(1.6) (Reinicke's theory) (unitless)
831!$OMP THREADPRIVATE(pipe_tune1)
832  REAL(r_std), SAVE :: pipe_tune2 = 40.0         !! height=pipe_tune2 * diameter**pipe_tune3 (unitless)
833!$OMP THREADPRIVATE(pipe_tune2)
834  REAL(r_std), SAVE :: pipe_tune3 = 0.5          !! height=pipe_tune2 * diameter**pipe_tune3 (unitless)
835!$OMP THREADPRIVATE(pipe_tune3)
836  REAL(r_std), SAVE :: pipe_tune4 = 0.3          !! needed for stem diameter (unitless)
837!$OMP THREADPRIVATE(pipe_tune4)
838  REAL(r_std), SAVE :: pipe_density = 2.e5       !! Density
839!$OMP THREADPRIVATE(pipe_density)
840  REAL(r_std), SAVE :: pipe_k1 = 8.e3            !! one more SAVE
841!$OMP THREADPRIVATE(pipe_k1)
842  REAL(r_std), SAVE :: pipe_tune_exp_coeff = 1.6 !! pipe tune exponential coeff (unitless)
843!$OMP THREADPRIVATE(pipe_tune_exp_coeff)
844
845  ! 1.2 climatic parameters
846
847  REAL(r_std), SAVE :: precip_crit = 100.        !! minimum precip, in (mm/year)
848!$OMP THREADPRIVATE(precip_crit)
849  REAL(r_std), SAVE :: gdd_crit_estab = 150.     !! minimum gdd for establishment of saplings
850!$OMP THREADPRIVATE(gdd_crit_estab)
851  REAL(r_std), SAVE :: fpc_crit = 0.95           !! critical fpc, needed for light competition and establishment (0-1, unitless)
852!$OMP THREADPRIVATE(fpc_crit)
853
854  ! 1.3 sapling characteristics
855
856  REAL(r_std), SAVE :: alpha_grass = 0.5         !! alpha coefficient for grasses (unitless)
857!$OMP THREADPRIVATE(alpha_grass)
858  REAL(r_std), SAVE :: alpha_tree = 1.           !! alpha coefficient for trees (unitless)
859!$OMP THREADPRIVATE(alpha_tree)
860  REAL(r_std), SAVE :: mass_ratio_heart_sap = 3. !! mass ratio (heartwood+sapwood)/sapwood (unitless)
861!$OMP THREADPRIVATE(mass_ratio_heart_sap)
862
863  ! 1.4  time scales for phenology and other processes (in days)
864
865  REAL(r_std), SAVE :: tau_hum_month = 20.        !! (days)       
866!$OMP THREADPRIVATE(tau_hum_month)
867  REAL(r_std), SAVE :: tau_hum_week = 7.          !! (days) 
868!$OMP THREADPRIVATE(tau_hum_week)
869  REAL(r_std), SAVE :: tau_t2m_month = 20.        !! (days)     
870!$OMP THREADPRIVATE(tau_t2m_month)
871  REAL(r_std), SAVE :: tau_t2m_week = 7.          !! (days) 
872!$OMP THREADPRIVATE(tau_t2m_week)
873  REAL(r_std), SAVE :: tau_tsoil_month = 20.      !! (days)     
874!$OMP THREADPRIVATE(tau_tsoil_month)
875  REAL(r_std), SAVE :: tau_soilhum_month = 20.    !! (days)     
876!$OMP THREADPRIVATE(tau_soilhum_month)
877  REAL(r_std), SAVE :: tau_gpp_week = 7.          !! (days) 
878!$OMP THREADPRIVATE(tau_gpp_week)
879  REAL(r_std), SAVE :: tau_gdd = 40.              !! (days) 
880!$OMP THREADPRIVATE(tau_gdd)
881  REAL(r_std), SAVE :: tau_ngd = 50.              !! (days) 
882!$OMP THREADPRIVATE(tau_ngd)
883  REAL(r_std), SAVE :: coeff_tau_longterm = 3.    !! (unitless)
884!$OMP THREADPRIVATE(coeff_tau_longterm)
885  REAL(r_std), SAVE :: tau_longterm_max           !! (days) 
886!$OMP THREADPRIVATE(tau_longterm_max)
887
888  ! 3. Coefficients of equations
889
890  REAL(r_std), SAVE :: bm_sapl_carbres = 5.             !!
891!$OMP THREADPRIVATE(bm_sapl_carbres)
892  REAL(r_std), SAVE :: bm_sapl_sapabove = 0.5           !!
893!$OMP THREADPRIVATE(bm_sapl_sapabove)
894  REAL(r_std), SAVE :: bm_sapl_heartabove = 2.          !!
895!$OMP THREADPRIVATE(bm_sapl_heartabove)
896  REAL(r_std), SAVE :: bm_sapl_heartbelow = 2.          !!
897!$OMP THREADPRIVATE(bm_sapl_heartbelow)
898  REAL(r_std), SAVE :: init_sapl_mass_leaf_nat = 0.1    !!
899!$OMP THREADPRIVATE(init_sapl_mass_leaf_nat)
900  REAL(r_std), SAVE :: init_sapl_mass_leaf_agri = 1.    !!
901!$OMP THREADPRIVATE(init_sapl_mass_leaf_agri)
902  REAL(r_std), SAVE :: init_sapl_mass_carbres = 5.      !!
903!$OMP THREADPRIVATE(init_sapl_mass_carbres)
904  REAL(r_std), SAVE :: init_sapl_mass_root = 0.1        !!
905!$OMP THREADPRIVATE(init_sapl_mass_root)
906  REAL(r_std), SAVE :: init_sapl_mass_fruit = 0.3       !! 
907!$OMP THREADPRIVATE(init_sapl_mass_fruit)
908  REAL(r_std), SAVE :: cn_sapl_init = 0.5               !!
909!$OMP THREADPRIVATE(cn_sapl_init)
910  REAL(r_std), SAVE :: migrate_tree = 10.*1.E3          !!
911!$OMP THREADPRIVATE(migrate_tree)
912  REAL(r_std), SAVE :: migrate_grass = 10.*1.E3         !!
913!$OMP THREADPRIVATE(migrate_grass)
914  REAL(r_std), SAVE :: lai_initmin_tree = 0.3           !!
915!$OMP THREADPRIVATE(lai_initmin_tree)
916  REAL(r_std), SAVE :: lai_initmin_grass = 0.1          !!
917!$OMP THREADPRIVATE(lai_initmin_grass)
918  REAL(r_std), SAVE, DIMENSION(2) :: dia_coeff = (/ 4., 0.5 /)            !!
919!$OMP THREADPRIVATE(dia_coeff)
920  REAL(r_std), SAVE, DIMENSION(2) :: maxdia_coeff =(/ 100., 0.01/)        !!
921!$OMP THREADPRIVATE(maxdia_coeff)
922  REAL(r_std), SAVE, DIMENSION(4) :: bm_sapl_leaf = (/ 4., 4., 0.8, 5./)  !!
923!$OMP THREADPRIVATE(bm_sapl_leaf)
924
925
926
927  !
928  ! stomate_litter.f90
929  !
930
931  ! 0. Constants
932
933  REAL(r_std), PARAMETER :: Q10 = 10.               !!
934
935  ! 1. Scalar
936
937  REAL(r_std), SAVE :: z_decomp = 0.2               !!  Maximum depth for soil decomposer's activity (m)
938!$OMP THREADPRIVATE(z_decomp)
939
940  ! 2. Arrays
941
942  REAL(r_std), SAVE :: frac_soil_struct_aa = 0.55   !! corresponding to frac_soil(istructural,iactive,iabove)
943!$OMP THREADPRIVATE(frac_soil_struct_aa)
944  REAL(r_std), SAVE :: frac_soil_struct_ab = 0.45   !! corresponding to frac_soil(istructural,iactive,ibelow)
945!$OMP THREADPRIVATE(frac_soil_struct_ab)
946  REAL(r_std), SAVE :: frac_soil_struct_sa = 0.7    !! corresponding to frac_soil(istructural,islow,iabove)
947!$OMP THREADPRIVATE(frac_soil_struct_sa)
948  REAL(r_std), SAVE :: frac_soil_struct_sb = 0.7    !! corresponding to frac_soil(istructural,islow,ibelow)
949!$OMP THREADPRIVATE(frac_soil_struct_sb)
950  REAL(r_std), SAVE :: frac_soil_metab_aa = 0.45    !! corresponding to frac_soil(imetabolic,iactive,iabove)
951!$OMP THREADPRIVATE(frac_soil_metab_aa)
952  REAL(r_std), SAVE :: frac_soil_metab_ab = 0.45    !! corresponding to frac_soil(imetabolic,iactive,ibelow)
953!$OMP THREADPRIVATE(frac_soil_metab_ab)
954  REAL(r_std), SAVE, DIMENSION(nparts) :: CN = &    !! C/N ratio of each plant pool (0-100, unitless)
955       & (/ 40., 40., 40., 40., 40., 40., 40., 40. /) 
956!$OMP THREADPRIVATE(CN)
957  REAL(r_std), SAVE, DIMENSION(nparts) :: LC = &    !! Lignin/C ratio of different plant parts (0,22-0,35, unitless)
958       & (/ 0.22, 0.35, 0.35, 0.35, 0.35, 0.22, 0.22, 0.22 /)
959!$OMP THREADPRIVATE(LC)
960
961  ! 3. Coefficients of equations
962
963  REAL(r_std), SAVE :: metabolic_ref_frac = 0.85    !! used by litter and soilcarbon (0-1, unitless)
964!$OMP THREADPRIVATE(metabolic_ref_frac)
965  REAL(r_std), SAVE :: metabolic_LN_ratio = 0.018   !! (0-1, unitless)   
966!$OMP THREADPRIVATE(metabolic_LN_ratio)
967  REAL(r_std), SAVE :: tau_metabolic = 0.066        !!
968!$OMP THREADPRIVATE(tau_metabolic)
969  REAL(r_std), SAVE :: tau_struct = 0.245           !!
970!$OMP THREADPRIVATE(tau_struct)
971  REAL(r_std), SAVE :: soil_Q10 = 0.69              !!= ln 2
972!$OMP THREADPRIVATE(soil_Q10)
973  REAL(r_std), SAVE :: tsoil_ref = 30.              !!
974!$OMP THREADPRIVATE(tsoil_ref)
975  REAL(r_std), SAVE :: litter_struct_coef = 3.      !!
976!$OMP THREADPRIVATE(litter_struct_coef)
977  REAL(r_std), SAVE, DIMENSION(3) :: moist_coeff = (/ 1.1,  2.4,  0.29 /) !!
978!$OMP THREADPRIVATE(moist_coeff)
979  REAL(r_std), SAVE :: moistcont_min = 0.25  !! minimum soil wetness to limit the heterotrophic respiration
980!$OMP THREADPRIVATE(moistcont_min)
981
982
983  !
984  ! stomate_lpj.f90
985  !
986
987  ! 1. Scalar
988
989  REAL(r_std), SAVE :: frac_turnover_daily = 0.55  !! (0-1, unitless)
990!$OMP THREADPRIVATE(frac_turnover_daily)
991
992
993  !
994  ! stomate_npp.f90
995  !
996
997  ! 1. Scalar
998
999  REAL(r_std), SAVE :: tax_max = 0.8 !! Maximum fraction of allocatable biomass used
1000                                     !! for maintenance respiration (0-1, unitless)
1001!$OMP THREADPRIVATE(tax_max)
1002
1003
1004  !
1005  ! stomate_phenology.f90
1006  !
1007
1008  ! 1. Scalar
1009
1010  REAL(r_std), SAVE :: min_growthinit_time = 300.  !! minimum time since last beginning of a growing season (days)
1011!$OMP THREADPRIVATE(min_growthinit_time)
1012  REAL(r_std), SAVE :: moiavail_always_tree = 1.0  !! moisture monthly availability above which moisture tendency doesn't matter
1013                                                   !!  - for trees (0-1, unitless)
1014!$OMP THREADPRIVATE(moiavail_always_tree)
1015  REAL(r_std), SAVE :: moiavail_always_grass = 0.6 !! moisture monthly availability above which moisture tendency doesn't matter
1016                                                   !! - for grass (0-1, unitless)
1017!$OMP THREADPRIVATE(moiavail_always_grass)
1018  REAL(r_std), SAVE :: t_always                    !! monthly temp. above which temp. tendency doesn't matter
1019!$OMP THREADPRIVATE(t_always)
1020  REAL(r_std), SAVE :: t_always_add = 10.          !! monthly temp. above which temp. tendency doesn't matter (C)
1021!$OMP THREADPRIVATE(t_always_add)
1022
1023  ! 3. Coefficients of equations
1024 
1025  REAL(r_std), SAVE :: gddncd_ref = 603.           !!
1026!$OMP THREADPRIVATE(gddncd_ref)
1027  REAL(r_std), SAVE :: gddncd_curve = 0.0091       !!
1028!$OMP THREADPRIVATE(gddncd_curve)
1029  REAL(r_std), SAVE :: gddncd_offset = 64.         !!
1030!$OMP THREADPRIVATE(gddncd_offset)
1031
1032
1033  !
1034  ! stomate_prescribe.f90
1035  !
1036
1037  ! 3. Coefficients of equations
1038
1039  REAL(r_std), SAVE :: bm_sapl_rescale = 40.       !!
1040!$OMP THREADPRIVATE(bm_sapl_rescale)
1041
1042
1043  !
1044  ! stomate_resp.f90
1045  !
1046
1047  ! 3. Coefficients of equations
1048
1049  REAL(r_std), SAVE :: maint_resp_min_vmax = 0.3   !!
1050!$OMP THREADPRIVATE(maint_resp_min_vmax)
1051  REAL(r_std), SAVE :: maint_resp_coeff = 1.4      !!
1052!$OMP THREADPRIVATE(maint_resp_coeff)
1053
1054
1055  !
1056  ! stomate_soilcarbon.f90
1057  !
1058
1059  ! 2. Arrays
1060
1061  ! 2.1 frac_carb_coefficients
1062
1063  REAL(r_std), SAVE :: frac_carb_ap = 0.004  !! from active pool: depends on clay content  (0-1, unitless)
1064                                             !! corresponding to frac_carb(:,iactive,ipassive)
1065!$OMP THREADPRIVATE(frac_carb_ap)
1066  REAL(r_std), SAVE :: frac_carb_sa = 0.42   !! from slow pool (0-1, unitless)
1067                                             !! corresponding to frac_carb(:,islow,iactive)
1068!$OMP THREADPRIVATE(frac_carb_sa)
1069  REAL(r_std), SAVE :: frac_carb_sp = 0.03   !! from slow pool (0-1, unitless)
1070                                             !! corresponding to frac_carb(:,islow,ipassive)
1071!$OMP THREADPRIVATE(frac_carb_sp)
1072  REAL(r_std), SAVE :: frac_carb_pa = 0.45   !! from passive pool (0-1, unitless)
1073                                             !! corresponding to frac_carb(:,ipassive,iactive)
1074!$OMP THREADPRIVATE(frac_carb_pa)
1075  REAL(r_std), SAVE :: frac_carb_ps = 0.0    !! from passive pool (0-1, unitless)
1076                                             !! corresponding to frac_carb(:,ipassive,islow)
1077!$OMP THREADPRIVATE(frac_carb_ps)
1078
1079  ! 3. Coefficients of equations
1080
1081  REAL(r_std), SAVE :: active_to_pass_clay_frac = 0.68  !! (0-1, unitless)
1082!$OMP THREADPRIVATE(active_to_pass_clay_frac)
1083  !! residence times in carbon pools (days)
1084  REAL(r_std), SAVE :: carbon_tau_iactive = 0.149   !! residence times in active pool (days)
1085!$OMP THREADPRIVATE(carbon_tau_iactive)
1086  REAL(r_std), SAVE :: carbon_tau_islow = 7.0       !! residence times in slow pool (days)
1087!$OMP THREADPRIVATE(carbon_tau_islow)
1088  REAL(r_std), SAVE :: carbon_tau_ipassive = 300.   !! residence times in passive pool (days)
1089!$OMP THREADPRIVATE(carbon_tau_ipassive)
1090  REAL(r_std), SAVE, DIMENSION(3) :: flux_tot_coeff = (/ 1.2, 1.4, .75/)
1091!$OMP THREADPRIVATE(flux_tot_coeff)
1092
1093  !
1094  ! stomate_turnover.f90
1095  !
1096
1097  ! 3. Coefficients of equations
1098
1099  REAL(r_std), SAVE :: new_turnover_time_ref = 20. !!(days)
1100!$OMP THREADPRIVATE(new_turnover_time_ref)
1101  REAL(r_std), SAVE :: leaf_age_crit_tref = 20.    !! (C)
1102!$OMP THREADPRIVATE(leaf_age_crit_tref)
1103  REAL(r_std), SAVE, DIMENSION(3) :: leaf_age_crit_coeff = (/ 1.5, 0.75, 10./) !! (unitless)
1104!$OMP THREADPRIVATE(leaf_age_crit_coeff)
1105
1106
1107  !
1108  ! stomate_vmax.f90
1109  !
1110 
1111  ! 1. Scalar
1112
1113  REAL(r_std), SAVE :: vmax_offset = 0.3        !! minimum leaf efficiency (unitless)
1114!$OMP THREADPRIVATE(vmax_offset)
1115  REAL(r_std), SAVE :: leafage_firstmax = 0.03  !! relative leaf age at which efficiency
1116                                                !! reaches 1 (unitless)
1117!$OMP THREADPRIVATE(leafage_firstmax)
1118  REAL(r_std), SAVE :: leafage_lastmax = 0.5    !! relative leaf age at which efficiency
1119                                                !! falls below 1 (unitless)
1120!$OMP THREADPRIVATE(leafage_lastmax)
1121  REAL(r_std), SAVE :: leafage_old = 1.         !! relative leaf age at which efficiency
1122                                                !! reaches its minimum (vmax_offset)
1123                                                !! (unitless)
1124!$OMP THREADPRIVATE(leafage_old)
1125  !
1126  ! stomate_season.f90
1127  !
1128
1129  ! 1. Scalar
1130
1131  REAL(r_std), SAVE :: gppfrac_dormance = 0.2  !! report maximal GPP/GGP_max for dormance (0-1, unitless)
1132!$OMP THREADPRIVATE(gppfrac_dormance)
1133  REAL(r_std), SAVE :: tau_climatology = 20.   !! tau for "climatologic variables (years)
1134!$OMP THREADPRIVATE(tau_climatology)
1135  REAL(r_std), SAVE :: hvc1 = 0.019            !! parameters for herbivore activity (unitless)
1136!$OMP THREADPRIVATE(hvc1)
1137  REAL(r_std), SAVE :: hvc2 = 1.38             !! parameters for herbivore activity (unitless)
1138!$OMP THREADPRIVATE(hvc2)
1139  REAL(r_std), SAVE :: leaf_frac_hvc = 0.33    !! leaf fraction (0-1, unitless)
1140!$OMP THREADPRIVATE(leaf_frac_hvc)
1141  REAL(r_std), SAVE :: tlong_ref_max = 303.1   !! maximum reference long term temperature (K)
1142!$OMP THREADPRIVATE(tlong_ref_max)
1143  REAL(r_std), SAVE :: tlong_ref_min = 253.1   !! minimum reference long term temperature (K)
1144!$OMP THREADPRIVATE(tlong_ref_min)
1145
1146  ! 3. Coefficients of equations
1147
1148  REAL(r_std), SAVE :: ncd_max_year = 3.
1149!$OMP THREADPRIVATE(ncd_max_year)
1150  REAL(r_std), SAVE :: gdd_threshold = 5.
1151!$OMP THREADPRIVATE(gdd_threshold)
1152  REAL(r_std), SAVE :: green_age_ever = 2.
1153!$OMP THREADPRIVATE(green_age_ever)
1154  REAL(r_std), SAVE :: green_age_dec = 0.5
1155!$OMP THREADPRIVATE(green_age_dec)
1156
1157END MODULE constantes_var
Note: See TracBrowser for help on using the repository browser.