source: tags/ORCHIDEE_2_0/ORCHIDEE/src_parameters/constantes_var.f90 @ 5391

Last change on this file since 5391 was 5391, checked in by josefine.ghattas, 6 years ago

Integrated changeset [5389] done on the trunk in the tag. This is needed for ESM configurations with transported carbon. P Cadule

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