1 | ! ================================================================================================================================= |
---|
2 | ! MODULE : pft_parameters |
---|
3 | ! |
---|
4 | ! CONTACT : orchidee-help _at_ ipsl.jussieu.fr |
---|
5 | ! |
---|
6 | ! LICENCE : IPSL (2011) |
---|
7 | ! This software is governed by the CeCILL licence see ORCHIDEE/ORCHIDEE_CeCILL.LIC |
---|
8 | ! |
---|
9 | !>\BRIEF This module initializes all the pft parameters in function of the |
---|
10 | !! number of vegetation types and of the values chosen by the user. |
---|
11 | !! |
---|
12 | !!\n DESCRIPTION: This module allocates and initializes the pft parameters in function of the number of pfts |
---|
13 | !! and the values of the parameters. \n |
---|
14 | !! The number of PFTs is read in intersurf.f90 (subroutine intsurf_config). \n |
---|
15 | !! Then we can initialize the parameters. \n |
---|
16 | !! This module is the result of the merge of constantes_co2, constantes_veg, stomate_constants.\n |
---|
17 | !! |
---|
18 | !! RECENT CHANGE(S): None |
---|
19 | !! |
---|
20 | !! REFERENCE(S) : None |
---|
21 | !! |
---|
22 | !! SVN : |
---|
23 | !! $HeadURL: $ |
---|
24 | !! $Date$ |
---|
25 | !! $Revision$ |
---|
26 | !! \n |
---|
27 | !_ ================================================================================================================================ |
---|
28 | |
---|
29 | MODULE pft_parameters |
---|
30 | |
---|
31 | USE constantes_mtc |
---|
32 | USE constantes |
---|
33 | USE ioipsl |
---|
34 | USE parallel |
---|
35 | USE defprec |
---|
36 | |
---|
37 | IMPLICIT NONE |
---|
38 | |
---|
39 | |
---|
40 | ! |
---|
41 | ! PFT GLOBAL |
---|
42 | ! |
---|
43 | INTEGER(i_std), SAVE :: nvm = 13 !! Number of vegetation types (2-N, unitless) |
---|
44 | |
---|
45 | INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pft_to_mtc !! Table of conversion : we associate one pft to one metaclass |
---|
46 | !! (1-13, unitless) |
---|
47 | |
---|
48 | CHARACTER(LEN=34), ALLOCATABLE, SAVE, DIMENSION(:) :: PFT_name !! Description of the PFT (unitless) |
---|
49 | |
---|
50 | LOGICAL, SAVE :: l_first_define_pft = .TRUE. !! To keep first call trace of the module (true/false) |
---|
51 | |
---|
52 | |
---|
53 | ! |
---|
54 | ! VEGETATION STRUCTURE |
---|
55 | ! |
---|
56 | !- |
---|
57 | ! 1. Sechiba |
---|
58 | ! |
---|
59 | ! 1.1 Labels - Characteristics |
---|
60 | !- |
---|
61 | LOGICAL, ALLOCATABLE, SAVE, DIMENSION(:) :: is_tree !! Is the vegetation type a tree ? (true/false) |
---|
62 | |
---|
63 | LOGICAL, ALLOCATABLE, SAVE, DIMENSION(:) :: is_deciduous !! Is PFT deciduous ? (true/false) |
---|
64 | |
---|
65 | LOGICAL, ALLOCATABLE, SAVE, DIMENSION(:) :: is_evergreen !! Is PFT evegreen ? (true/false) |
---|
66 | |
---|
67 | LOGICAL, ALLOCATABLE, SAVE, DIMENSION(:) :: is_c3 !! Is PFT c3 ? (true/false) |
---|
68 | |
---|
69 | CHARACTER(len=5), ALLOCATABLE, SAVE, DIMENSION(:) :: type_of_lai !! Type of behaviour of the LAI evolution algorithm |
---|
70 | !! for each vegetation type. |
---|
71 | !! Value of type_of_lai, one for each vegetation type : |
---|
72 | !! mean or interp |
---|
73 | !- |
---|
74 | ! 1.2 Prescribed Values |
---|
75 | !- |
---|
76 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: veget_ori_fixed_test_1 !! Value for veget_ori for tests in 0-dim simulations |
---|
77 | !! (0-1, unitless) |
---|
78 | |
---|
79 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: llaimax !! laimax for maximum lai see also type of lai |
---|
80 | !! interpolation (m^2.m^{-2}) |
---|
81 | |
---|
82 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: llaimin !! laimin for minimum lai see also type of lai |
---|
83 | !! interpolation (m^2.m^{-2}) |
---|
84 | |
---|
85 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: height_presc !! prescribed height of vegetation.(m) |
---|
86 | !! Value for height_presc : one for each vegetation type |
---|
87 | |
---|
88 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: rveg_pft !! Potentiometer to set vegetation resistance (unitless) |
---|
89 | !! Nathalie on March 28th, 2006 - from Fred Hourdin, |
---|
90 | !- |
---|
91 | ! 2. Stomate |
---|
92 | ! |
---|
93 | ! 2.1 Labels - Characteristics |
---|
94 | !- |
---|
95 | LOGICAL, ALLOCATABLE, SAVE, DIMENSION(:) :: natural !! natural? (true/false) |
---|
96 | |
---|
97 | INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: leaf_tab !! leaf type (1-4, unitless) |
---|
98 | !! 1=broad leaved tree, 2=needle leaved tree, |
---|
99 | !! 3=grass 4=bare ground |
---|
100 | !- |
---|
101 | ! 2.2 Prescribed Values |
---|
102 | !- |
---|
103 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: sla !! specif leaf area (m^2.gC^{-1}) |
---|
104 | |
---|
105 | |
---|
106 | ! |
---|
107 | ! EVAPOTRANSPIRATION (sechiba) |
---|
108 | ! |
---|
109 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: rstruct_const !! Structural resistance. (s.m^{-1}) |
---|
110 | !! Value for rstruct_const : one for each vegetation type |
---|
111 | |
---|
112 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: kzero !! A vegetation dependent constant used in the calculation |
---|
113 | !! of the surface resistance. (kg.m^2.s^{-1}) |
---|
114 | !! Value for kzero one for each vegetation type |
---|
115 | |
---|
116 | |
---|
117 | ! |
---|
118 | ! WATER (sechiba) |
---|
119 | ! |
---|
120 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: wmax_veg !! Maximum field capacity for each of the vegetations (Temporary). |
---|
121 | !! Value of wmax_veg : max quantity of water : |
---|
122 | !! one for each vegetation type (kg.m^{-3}) |
---|
123 | |
---|
124 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: humcste !! Root profile description for the different vegetation types. |
---|
125 | !! These are the factor in the exponential which gets |
---|
126 | !! the root density as a function of depth (m^{-1}) |
---|
127 | |
---|
128 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: throughfall_by_pft !! Fraction of rain intercepted by the canopy (0-100, unitless) |
---|
129 | |
---|
130 | |
---|
131 | ! |
---|
132 | ! ALBEDO (sechiba) |
---|
133 | ! |
---|
134 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: snowa_ini !! Initial snow albedo value for each vegetation type |
---|
135 | !! as it will be used in condveg_snow (unitless) |
---|
136 | !! Source : Values are from the Thesis of S. Chalita (1992) |
---|
137 | |
---|
138 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: snowa_dec !! Decay rate of snow albedo value for each vegetation type |
---|
139 | !! as it will be used in condveg_snow (unitless) |
---|
140 | !! Source : Values are from the Thesis of S. Chalita (1992) |
---|
141 | |
---|
142 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: alb_leaf_vis !! leaf albedo of vegetation type, visible albedo (unitless) |
---|
143 | |
---|
144 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: alb_leaf_nir !! leaf albedo of vegetation type, near infrared albedo (unitless) |
---|
145 | |
---|
146 | |
---|
147 | ! |
---|
148 | ! SOIL - VEGETATION |
---|
149 | ! |
---|
150 | INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: pref_soil_veg !! Table which contains the correlation between the soil |
---|
151 | !! types and vegetation type. Two modes exist : |
---|
152 | !! 1) pref_soil_veg = 0 then we have an equidistribution |
---|
153 | !! of vegetation on soil types |
---|
154 | !! 2) Else for each pft the prefered soil type is given : |
---|
155 | !! 1=sand, 2=loan, 3=clay |
---|
156 | !! This variable is initialized in slowproc.(1-3, unitless) |
---|
157 | |
---|
158 | INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pref_soil_veg_sand !! Prefered soil type for the first layer of soil : |
---|
159 | !! 1=sand, 2=loan, 3=clay (1-3, unitless) |
---|
160 | |
---|
161 | INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pref_soil_veg_loan !! Prefered soil type for the second layer of soil : |
---|
162 | !! 1=sand, 2=loan, 3=clay (1-3, unitless) |
---|
163 | |
---|
164 | INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pref_soil_veg_clay !! Prefered soil type for the third layer of soil : |
---|
165 | !! 1=sand, 2=loan, 3=clay (1-3, unitless) |
---|
166 | |
---|
167 | ! |
---|
168 | ! PHOTOSYNTHESIS |
---|
169 | ! |
---|
170 | !- |
---|
171 | ! 1. CO2 |
---|
172 | !- |
---|
173 | LOGICAL, ALLOCATABLE, SAVE, DIMENSION(:) :: is_c4 !! flag for C4 vegetation types (true/false) |
---|
174 | |
---|
175 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: gsslope !! Slope of the gs/A relation (Ball & al.) (unitless) |
---|
176 | |
---|
177 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: gsoffset !! intercept of the gs/A relation (Ball & al.) |
---|
178 | |
---|
179 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: vcmax_fix !! values used for vcmax when STOMATE is not activated |
---|
180 | !! (µmol.m^{-2}.s^{-1}) |
---|
181 | |
---|
182 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: vjmax_fix !! values used for vjmax when STOMATE is not activated |
---|
183 | !! (µmol.m^{-2}.s^{-1}) |
---|
184 | |
---|
185 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: co2_tmin_fix !! values used for photosynthesis tmin when STOMATE |
---|
186 | !! is not activated (C) |
---|
187 | |
---|
188 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: co2_topt_fix !! values used for photosynthesis topt when STOMATE |
---|
189 | !! is not activated (C) |
---|
190 | |
---|
191 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: co2_tmax_fix !! values used for photosynthesis tmax when STOMATE |
---|
192 | !! is not activated (C) |
---|
193 | !- |
---|
194 | ! 2. Stomate |
---|
195 | !- |
---|
196 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: ext_coeff !! extinction coefficient of the Monsi&Saeki relationship (1953) |
---|
197 | !! (unitless) |
---|
198 | |
---|
199 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: vcmax_opt !! Maximum rate of carboxylation (µmol.m^{-2}.s^{-1}) |
---|
200 | |
---|
201 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: vjmax_opt !! Maximum rate of RUbp regeneration (µmol.m^{-2}.s^{-1}) |
---|
202 | |
---|
203 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tphoto_min_a !! minimum photosynthesis temperature, |
---|
204 | !! constant a of ax^2+bx+c (deg C),tabulated (unitless) |
---|
205 | |
---|
206 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tphoto_min_b !! minimum photosynthesis temperature, |
---|
207 | !! constant b of ax^2+bx+c (deg C),tabulated (unitless) |
---|
208 | |
---|
209 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tphoto_min_c !! minimum photosynthesis temperature, |
---|
210 | !! constant c of ax^2+bx+c (deg C),tabulated (unitless) |
---|
211 | |
---|
212 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tphoto_opt_a !! optimum photosynthesis temperature, |
---|
213 | !! constant a of ax^2+bx+c (deg C),tabulated (unitless) |
---|
214 | |
---|
215 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tphoto_opt_b !! optimum photosynthesis temperature, |
---|
216 | !! constant b of ax^2+bx+c (deg C),tabulated (unitless) |
---|
217 | |
---|
218 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tphoto_opt_c !! optimum photosynthesis temperature, |
---|
219 | !! constant c of ax^2+bx+c (deg C),tabulated (unitless) |
---|
220 | |
---|
221 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tphoto_max_a !! maximum photosynthesis temperature, |
---|
222 | !! constant a of ax^2+bx+c (deg C), tabulated (unitless) |
---|
223 | |
---|
224 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tphoto_max_b !! maximum photosynthesis temperature, |
---|
225 | !! constant b of ax^2+bx+c (deg C), tabulated (unitless) |
---|
226 | |
---|
227 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tphoto_max_c !! maximum photosynthesis temperature, |
---|
228 | !! constant c of ax^2+bx+c (deg C), tabulated (unitless) |
---|
229 | |
---|
230 | |
---|
231 | ! |
---|
232 | ! RESPIRATION (stomate) |
---|
233 | ! |
---|
234 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: maint_resp_slope !! slope of maintenance respiration coefficient |
---|
235 | !! (1/K, 1/K^2, 1/K^3), used in the code |
---|
236 | |
---|
237 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: maint_resp_slope_c !! slope of maintenance respiration coefficient (1/K), |
---|
238 | !! constant c of aT^2+bT+c , tabulated |
---|
239 | |
---|
240 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: maint_resp_slope_b !! slope of maintenance respiration coefficient (1/K), |
---|
241 | !! constant b of aT^2+bT+c , tabulated |
---|
242 | |
---|
243 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: maint_resp_slope_a !! slope of maintenance respiration coefficient (1/K), |
---|
244 | !! constant a of aT^2+bT+c , tabulated |
---|
245 | |
---|
246 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: coeff_maint_zero !! maintenance respiration coefficient at 0 deg C, |
---|
247 | !! used in the code (gC.gC^{-1}.day^{-1}) |
---|
248 | |
---|
249 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: cm_zero_leaf !! maintenance respiration coefficient at 0 deg C, |
---|
250 | !! for leaves, tabulated (gC.gC^{-1}.day^{-1}) |
---|
251 | |
---|
252 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: cm_zero_sapabove !! maintenance respiration coefficient at 0 deg C, |
---|
253 | !! for sapwood above, tabulated (gC.gC^{-1}.day^{-1}) |
---|
254 | |
---|
255 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: cm_zero_sapbelow !! maintenance respiration coefficient at 0 deg C, |
---|
256 | !! for sapwood below, tabulated (gC.gC^{-1}.day^{-1}) |
---|
257 | |
---|
258 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: cm_zero_heartabove !! maintenance respiration coefficient at 0 deg C |
---|
259 | !! for heartwood above, tabulated (gC.gC^{-1}.day^{-1}) |
---|
260 | |
---|
261 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: cm_zero_heartbelow !! maintenance respiration coefficient at 0 deg C, |
---|
262 | !! for heartwood below, tabulated (gC.gC^{-1}.day^{-1}) |
---|
263 | |
---|
264 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: cm_zero_root !! maintenance respiration coefficient at 0 deg C, |
---|
265 | !! for roots, tabulated (gC.gC^{-1}.day^{-1}) |
---|
266 | |
---|
267 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: cm_zero_fruit !! maintenance respiration coefficient at 0 deg C, |
---|
268 | !! for fruits, tabulated (gC.gC^{-1}.day^{-1}) |
---|
269 | |
---|
270 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: cm_zero_carbres !! maintenance respiration coefficient at 0 deg C, |
---|
271 | !! for carbohydrate reserve, tabulated (gC.gC^{-1}.day^{-1}) |
---|
272 | |
---|
273 | |
---|
274 | ! |
---|
275 | ! FIRE (stomate) |
---|
276 | ! |
---|
277 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: flam !! flamability : critical fraction of water holding |
---|
278 | !! capacity (0-1, unitless) |
---|
279 | |
---|
280 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: resist !! fire resistance (0-1, unitless) |
---|
281 | |
---|
282 | |
---|
283 | ! |
---|
284 | ! FLUX - LUC (Land Use Change) |
---|
285 | ! |
---|
286 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: coeff_lcchange_1 !! Coeff of biomass export for the year (unitless) |
---|
287 | |
---|
288 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: coeff_lcchange_10 !! Coeff of biomass export for the decade (unitless) |
---|
289 | |
---|
290 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: coeff_lcchange_100 !! Coeff of biomass export for the century (unitless) |
---|
291 | |
---|
292 | |
---|
293 | ! |
---|
294 | ! PHENOLOGY |
---|
295 | ! |
---|
296 | !- |
---|
297 | ! 1. Stomate |
---|
298 | !- |
---|
299 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: lai_max !! maximum LAI, PFT-specific (m^2.m^{-2}) |
---|
300 | |
---|
301 | CHARACTER(len=6), ALLOCATABLE, SAVE, DIMENSION(:) :: pheno_model !! which phenology model is used? (tabulated) (unitless) |
---|
302 | |
---|
303 | INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pheno_type !! type of phenology (0-4, unitless) |
---|
304 | !! 0=bare ground 1=evergreen, 2=summergreen, |
---|
305 | !! 3=raingreen, 4=perennial |
---|
306 | !! For the moment, the bare ground phenotype is not managed, |
---|
307 | !! so it is considered as "evergreen" |
---|
308 | !- |
---|
309 | ! 2. Leaf Onset |
---|
310 | !- |
---|
311 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: pheno_gdd_crit !! critical gdd,tabulated (C), used in the code |
---|
312 | |
---|
313 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pheno_gdd_crit_c !! critical gdd,tabulated (C), |
---|
314 | !! constant c of aT^2+bT+c (unitless) |
---|
315 | |
---|
316 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pheno_gdd_crit_b !! critical gdd,tabulated (C), |
---|
317 | !! constant b of aT^2+bT+c (unitless) |
---|
318 | |
---|
319 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pheno_gdd_crit_a !! critical gdd,tabulated (C), |
---|
320 | !! constant a of aT^2+bT+c (unitless) |
---|
321 | |
---|
322 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: ngd_crit !! critical ngd,tabulated. Threshold -5 degrees (days) |
---|
323 | |
---|
324 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: ncdgdd_temp !! critical temperature for the ncd vs. gdd function |
---|
325 | !! in phenology (C) |
---|
326 | |
---|
327 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: hum_frac !! critical humidity (relative to min/max) for phenology |
---|
328 | !! (0-1, unitless) |
---|
329 | |
---|
330 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: lowgpp_time !! minimum duration of dormance (days) |
---|
331 | |
---|
332 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: hum_min_time !! minimum time elapsed since moisture minimum (days) |
---|
333 | |
---|
334 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tau_sap !! sapwood -> heartwood conversion time (days) |
---|
335 | |
---|
336 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tau_fruit !! fruit lifetime (days) |
---|
337 | |
---|
338 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: ecureuil !! fraction of primary leaf and root allocation put |
---|
339 | !! into reserve (0-1, unitless) |
---|
340 | |
---|
341 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: alloc_min !! NEW - allocation above/below = f(age) - 30/01/04 NV/JO/PF |
---|
342 | |
---|
343 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: alloc_max !! NEW - allocation above/below = f(age) - 30/01/04 NV/JO/PF |
---|
344 | |
---|
345 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: demi_alloc !! NEW - allocation above/below = f(age) - 30/01/04 NV/JO/PF |
---|
346 | |
---|
347 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: leaflife_tab !! leaf longevity, tabulated (??units??) |
---|
348 | !- |
---|
349 | ! 3. Senescence |
---|
350 | !- |
---|
351 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: leaffall !! length of death of leaves,tabulated (days) |
---|
352 | |
---|
353 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: leafagecrit !! critical leaf age,tabulated (days) |
---|
354 | |
---|
355 | CHARACTER(len=6), ALLOCATABLE, SAVE, DIMENSION(:) :: senescence_type !! type of senescence,tabulated (unitless) |
---|
356 | !! List of avaible types of senescence : |
---|
357 | !! 'cold ', 'dry ', 'mixed ', 'none ' |
---|
358 | |
---|
359 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: senescence_hum !! critical relative moisture availability for senescence |
---|
360 | !! (0-1, unitless) |
---|
361 | |
---|
362 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: nosenescence_hum !! relative moisture availability above which there is |
---|
363 | !! no humidity-related senescence (0-1, unitless) |
---|
364 | |
---|
365 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: max_turnover_time !! maximum turnover time for grasses (days) |
---|
366 | |
---|
367 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: min_turnover_time !! minimum turnover time for grasses (days) |
---|
368 | |
---|
369 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: min_leaf_age_for_senescence !! minimum leaf age to allow senescence g (days) |
---|
370 | |
---|
371 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: senescence_temp !! critical temperature for senescence (C), |
---|
372 | !! used in the code |
---|
373 | |
---|
374 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: senescence_temp_c !! critical temperature for senescence (C), |
---|
375 | !! constant c of aT^2+bT+c , tabulated (unitless) |
---|
376 | |
---|
377 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: senescence_temp_b !! critical temperature for senescence (C), |
---|
378 | !! constant b of aT^2+bT+c , tabulated (unitless) |
---|
379 | |
---|
380 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: senescence_temp_a !! critical temperature for senescence (C), |
---|
381 | !! constant a of aT^2+bT+c , tabulated (unitless) |
---|
382 | |
---|
383 | |
---|
384 | ! |
---|
385 | ! DGVM |
---|
386 | ! |
---|
387 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: residence_time !! residence time of trees (y) |
---|
388 | |
---|
389 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tmin_crit !! critical tmin, tabulated (C) |
---|
390 | |
---|
391 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: tcm_crit !! critical tcm, tabulated (C) |
---|
392 | |
---|
393 | ! |
---|
394 | ! INTERNAL PARAMETERS USED IN STOMATE_DATA |
---|
395 | ! |
---|
396 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: lai_initmin !! Initial lai for trees/grass (m^2.m^{-2}) |
---|
397 | |
---|
398 | LOGICAL, ALLOCATABLE, SAVE, DIMENSION(:) :: tree !! is pft a tree? (used for consistency with is_tree) |
---|
399 | !! (true/false) |
---|
400 | |
---|
401 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: bm_sapl !! sapling biomass (gC.ind^{-1}) |
---|
402 | |
---|
403 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: migrate !! migration speed (m.year^{-1}) |
---|
404 | |
---|
405 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: maxdia !! maximum stem diameter from which on crown area no longer |
---|
406 | !! increases (m)m |
---|
407 | |
---|
408 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: cn_sapl !! crown of tree when sapling (m^2) |
---|
409 | |
---|
410 | REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: leaf_timecst !! time constant for leaf age discretisation (days) |
---|
411 | |
---|
412 | |
---|
413 | CONTAINS |
---|
414 | ! |
---|
415 | |
---|
416 | !! ================================================================================================================================ |
---|
417 | !! SUBROUTINE : pft_parameters_main |
---|
418 | !! |
---|
419 | !>\BRIEF This subroutine initializes all the pft parameters in function of the |
---|
420 | !! number of vegetation types chosen by the user. |
---|
421 | !! |
---|
422 | !! DESCRIPTION : This subroutine is called after the reading of the number of PFTS and the options |
---|
423 | !! activated by the user in the configuration files. (structure active_flags) \n |
---|
424 | !! The allocation is done just before reading the correspondence table between PFTs and MTCs |
---|
425 | !! defined by the user in the configuration file.\n |
---|
426 | !! With the correspondence table, the subroutine can initialize the pft parameters in function |
---|
427 | !! of the flags activated (ok_sechiba, ok_stomate, ok_co2, routing, new_hydrol...) in order to |
---|
428 | !! optimize the memory allocation. \n |
---|
429 | !! If the number of PFTs and pft_to_mtc are not found, the standard configuration will be used |
---|
430 | !! (13 PFTs, PFT = MTC). \n |
---|
431 | !! Some restrictions : the pft 1 can only be the bare soil and it is unique. \n |
---|
432 | !! |
---|
433 | !! RECENT CHANGE(S): None |
---|
434 | !! |
---|
435 | !! MAIN OUTPUT VARIABLE(S): None |
---|
436 | !! |
---|
437 | !! REFERENCE(S) : None |
---|
438 | !! |
---|
439 | !! FLOWCHART : None |
---|
440 | !! \n |
---|
441 | !_ ================================================================================================================================ |
---|
442 | |
---|
443 | SUBROUTINE pft_parameters_main(active_flags) |
---|
444 | |
---|
445 | IMPLICIT NONE |
---|
446 | |
---|
447 | !! 0. Variables and parameters declaration |
---|
448 | |
---|
449 | !! 0.1 Input variables |
---|
450 | |
---|
451 | TYPE(control_type),INTENT(in) :: active_flags !! What parts of the code are activated ? (true/false) |
---|
452 | |
---|
453 | !! 0.4 Local variables |
---|
454 | |
---|
455 | INTEGER(i_std) :: i !! Index (unitless) |
---|
456 | |
---|
457 | !_ ================================================================================================================================ |
---|
458 | |
---|
459 | ! |
---|
460 | ! PFT global |
---|
461 | ! |
---|
462 | |
---|
463 | IF(l_first_define_pft) THEN |
---|
464 | |
---|
465 | !! 1. First time step |
---|
466 | IF(long_print) THEN |
---|
467 | WRITE(numout,*) 'l_first_define_pft :we read the parameters from the def files' |
---|
468 | ENDIF |
---|
469 | |
---|
470 | !! 2. Memory allocation for the pfts-parameters |
---|
471 | CALL pft_parameters_alloc(active_flags) |
---|
472 | |
---|
473 | !! 3. Correspondance table |
---|
474 | |
---|
475 | !! 3.1 Initialisation of the correspondance table |
---|
476 | !! Initialisation of the correspondance table |
---|
477 | IF (nvm == nvmc) THEN |
---|
478 | WRITE(numout,*) 'Message to the user : we will use ORCHIDEE to its standard configuration' |
---|
479 | pft_to_mtc = (/ 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13 /) |
---|
480 | ELSE |
---|
481 | pft_to_mtc(:) = undef_int |
---|
482 | ENDIF !(nvm == nvmc) |
---|
483 | |
---|
484 | !! 3.2 Reading of the conrrespondance table in the .def file |
---|
485 | ! |
---|
486 | !Config Key = PFT_TO_MTC |
---|
487 | !Config Desc = correspondance array linking a PFT to MTC |
---|
488 | !Config if = OK_SECHIBA or OK_STOMATE |
---|
489 | !Config Def = 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13 |
---|
490 | !Config Help = |
---|
491 | !Config Units = [-] |
---|
492 | CALL getin_p('PFT_TO_MTC',pft_to_mtc) |
---|
493 | |
---|
494 | !! 3.3 If the user want to use the standard configuration, he needn't to fill the correspondance array |
---|
495 | !! If the configuration is wrong, send a error message to the user. |
---|
496 | IF(nvm /= nvmc ) THEN |
---|
497 | ! |
---|
498 | IF(pft_to_mtc(1) == undef_int) THEN |
---|
499 | STOP ' The array PFT_TO_MTC is empty : we stop' |
---|
500 | ENDIF !(pft_to_mtc(1) == undef_int) |
---|
501 | ! |
---|
502 | ENDIF !(nvm /= nvmc ) |
---|
503 | |
---|
504 | !! 3.4 Some error messages |
---|
505 | |
---|
506 | !! 3.4.1 What happened if pft_to_mtc(i) > nvmc or pft_to_mtc(i) <=0 (if the mtc doesn't exist)? |
---|
507 | DO i = 1, nvm ! Loop over # PFTs |
---|
508 | ! |
---|
509 | IF( (pft_to_mtc(i) > nvmc) .OR. (pft_to_mtc(i) <= 0) ) THEN |
---|
510 | WRITE(numout,*) 'the metaclass chosen does not exist' |
---|
511 | STOP 'we stop reading pft_to_mtc' |
---|
512 | ENDIF !( (pft_to_mtc(i) > nvmc) .OR. (pft_to_mtc(i) <= 0) ) |
---|
513 | ! |
---|
514 | ENDDO ! Loop over # PFTs |
---|
515 | |
---|
516 | |
---|
517 | !! 3.4.2 Check if pft_to_mtc(1) = 1 |
---|
518 | IF(pft_to_mtc(1) /= 1) THEN |
---|
519 | ! |
---|
520 | WRITE(numout,*) 'the first pft has to be the bare soil' |
---|
521 | STOP 'we stop reading next values of pft_to_mtc' |
---|
522 | ! |
---|
523 | ELSE |
---|
524 | ! |
---|
525 | DO i = 2,nvm ! Loop over # PFTs different from bare soil |
---|
526 | ! |
---|
527 | IF(pft_to_mtc(i) == 1) THEN |
---|
528 | WRITE(numout,*) 'only pft_to_mtc(1) has to be the bare soil' |
---|
529 | STOP 'we stop reading pft_to_mtc' |
---|
530 | ENDIF ! (pft_to_mtc(i) == 1) |
---|
531 | ! |
---|
532 | ENDDO ! Loop over # PFTs different from bare soil |
---|
533 | ! |
---|
534 | ENDIF !(pft_to_mtc(1) /= 1) |
---|
535 | |
---|
536 | |
---|
537 | !! 4.Initialisation of the pfts-parameters |
---|
538 | CALL pft_parameters_init(active_flags) |
---|
539 | |
---|
540 | !! 5. Useful data |
---|
541 | |
---|
542 | !! 5.1 Read the name of the PFTs given by the user |
---|
543 | ! |
---|
544 | !Config Key = PFT_NAME |
---|
545 | !Config Desc = Name of a PFT |
---|
546 | !Config if = OK_SECHIBA or OK_STOMATE |
---|
547 | !Config Def = bare ground, tropical broad-leaved evergreen, tropical broad-leaved raingreen, |
---|
548 | !Config temperate needleleaf evergreen, temperate broad-leaved evergreen temperate broad-leaved summergreen, |
---|
549 | !Config boreal needleleaf evergreen, boreal broad-leaved summergreen, boreal needleleaf summergreen, |
---|
550 | !Config C3 grass, C4 grass, C3 agriculture, C4 agriculture |
---|
551 | !Config Help = the user can name the new PFTs he/she introducing for new species |
---|
552 | !Config Units = [-] |
---|
553 | CALL getin_p('PFT_NAME',pft_name) |
---|
554 | |
---|
555 | !! 5.2 A useful message to the user: correspondance between the number of the pft |
---|
556 | !! and the name of the associated mtc |
---|
557 | DO i = 1,nvm ! Loop over # PFTs |
---|
558 | |
---|
559 | WRITE(numout,*) 'the PFT',i, 'called ', PFT_name(i),'corresponds to the MTC : ',MTC_name(pft_to_mtc(i)) |
---|
560 | |
---|
561 | ENDDO ! Loop over # PFTs |
---|
562 | |
---|
563 | !! 6. Initialisation of 2D arrays used in the code |
---|
564 | |
---|
565 | !- pref_soil_veg |
---|
566 | pref_soil_veg(:,:) = undef_int |
---|
567 | |
---|
568 | IF (active_flags%ok_stomate) THEN |
---|
569 | !- pheno_gdd_crit |
---|
570 | pheno_gdd_crit(:,:) = zero |
---|
571 | ! |
---|
572 | !- senescence_temp |
---|
573 | senescence_temp(:,:) = zero |
---|
574 | ! |
---|
575 | !- maint_resp_slope |
---|
576 | maint_resp_slope(:,:) = zero |
---|
577 | ! |
---|
578 | !-coeff_maint_zero |
---|
579 | coeff_maint_zero(:,:) = zero |
---|
580 | ENDIF !(active_flags%ok_stomate) |
---|
581 | |
---|
582 | !! 7. End message |
---|
583 | IF(long_print) THEN |
---|
584 | WRITE(numout,*) 'pft_parameters_done' |
---|
585 | ENDIF |
---|
586 | |
---|
587 | !! 8. Reset flag |
---|
588 | l_first_define_pft = .FALSE. |
---|
589 | |
---|
590 | ELSE |
---|
591 | |
---|
592 | RETURN |
---|
593 | |
---|
594 | ENDIF !(l_first_define_pft) |
---|
595 | |
---|
596 | END SUBROUTINE pft_parameters_main |
---|
597 | ! |
---|
598 | != |
---|
599 | ! |
---|
600 | |
---|
601 | !! ================================================================================================================================ |
---|
602 | !! SUBROUTINE : pft_parameters_init |
---|
603 | !! |
---|
604 | !>\BRIEF This subroutine initializes all the pft parameters by the default values |
---|
605 | !! of the corresponding metaclasse. |
---|
606 | !! |
---|
607 | !! DESCRIPTION : This subroutine is called after the reading of the number of PFTS and the correspondence |
---|
608 | !! table defined by the user in the configuration files. \n |
---|
609 | !! With the correspondence table, the subroutine can search the default values for the parameter |
---|
610 | !! even if the PFTs are classified in a random order (except bare soil). \n |
---|
611 | !! With the correspondence table, the subroutine can initialize the pft parameters in function |
---|
612 | !! of the flags activated (ok_sechiba, ok_stomate, ok_co2, routing, new_hydrol...).\n |
---|
613 | !! |
---|
614 | !! RECENT CHANGE(S): None |
---|
615 | !! |
---|
616 | !! MAIN OUTPUT VARIABLE(S): None |
---|
617 | !! |
---|
618 | !! REFERENCE(S) : None |
---|
619 | !! |
---|
620 | !! FLOWCHART : None |
---|
621 | !! \n |
---|
622 | !_ ================================================================================================================================ |
---|
623 | |
---|
624 | SUBROUTINE pft_parameters_init(active_flags) |
---|
625 | |
---|
626 | IMPLICIT NONE |
---|
627 | |
---|
628 | !! 0. Variables and parameters declaration |
---|
629 | |
---|
630 | !! 0.1 Input variables |
---|
631 | |
---|
632 | TYPE(control_type),INTENT(in) :: active_flags !! What parts of the code are activated ? (true/false) |
---|
633 | |
---|
634 | !! 0.4 Local variables |
---|
635 | |
---|
636 | INTEGER(i_std) :: j !! Index (unitless) |
---|
637 | |
---|
638 | !_ ================================================================================================================================ |
---|
639 | |
---|
640 | ! |
---|
641 | ! 1. Correspondance between the PFTs values and thes MTCs values |
---|
642 | ! |
---|
643 | |
---|
644 | |
---|
645 | ! 1.1 For parameters used anytime |
---|
646 | |
---|
647 | DO j= 1, nvm ! Loop over PFTs |
---|
648 | !- |
---|
649 | PFT_name(j) = MTC_name(pft_to_mtc(j)) |
---|
650 | !- |
---|
651 | ! Vegetation structure |
---|
652 | !- |
---|
653 | ! |
---|
654 | ! 1 .Sechiba |
---|
655 | ! |
---|
656 | veget_ori_fixed_test_1(j) = veget_ori_fixed_mtc(pft_to_mtc(j)) |
---|
657 | llaimax(j) = llaimax_mtc(pft_to_mtc(j)) |
---|
658 | llaimin(j) = llaimin_mtc(pft_to_mtc(j)) |
---|
659 | height_presc(j) = height_presc_mtc(pft_to_mtc(j)) |
---|
660 | type_of_lai(j) = type_of_lai_mtc(pft_to_mtc(j)) |
---|
661 | is_tree(j) = is_tree_mtc(pft_to_mtc(j)) |
---|
662 | natural(j) = natural_mtc(pft_to_mtc(j)) |
---|
663 | ! |
---|
664 | ! 2 .Stomate |
---|
665 | ! |
---|
666 | is_deciduous(j) = is_deciduous_mtc(pft_to_mtc(j)) |
---|
667 | is_evergreen(j) = is_evergreen_mtc(pft_to_mtc(j)) |
---|
668 | is_c3(j) = is_c3_mtc(pft_to_mtc(j)) |
---|
669 | !- |
---|
670 | ! Water - sechiba |
---|
671 | !- |
---|
672 | humcste(j) = humcste_mtc(pft_to_mtc(j)) |
---|
673 | !- |
---|
674 | ! Soil - vegetation |
---|
675 | !- |
---|
676 | pref_soil_veg_sand(j) = pref_soil_veg_sand_mtc(pft_to_mtc(j)) |
---|
677 | pref_soil_veg_loan(j) = pref_soil_veg_loan_mtc(pft_to_mtc(j)) |
---|
678 | pref_soil_veg_clay(j) = pref_soil_veg_clay_mtc(pft_to_mtc(j)) |
---|
679 | !- |
---|
680 | ! Photosynthesis |
---|
681 | !- |
---|
682 | is_c4(j) = is_c4_mtc(pft_to_mtc(j)) |
---|
683 | gsslope(j) = gsslope_mtc(pft_to_mtc(j)) |
---|
684 | gsoffset(j) = gsoffset_mtc(pft_to_mtc(j)) |
---|
685 | vcmax_fix(j) = vcmax_fix_mtc(pft_to_mtc(j)) |
---|
686 | vjmax_fix(j) = vjmax_fix_mtc(pft_to_mtc(j)) |
---|
687 | co2_tmin_fix(j) = co2_tmin_fix_mtc(pft_to_mtc(j)) |
---|
688 | co2_topt_fix(j) = co2_topt_fix_mtc(pft_to_mtc(j)) |
---|
689 | co2_tmax_fix(j) = co2_tmax_fix_mtc(pft_to_mtc(j)) |
---|
690 | ext_coeff(j) = ext_coeff_mtc(pft_to_mtc(j)) |
---|
691 | !- |
---|
692 | ENDDO ! Loop over PFTs |
---|
693 | |
---|
694 | ! 1.2 For sechiba parameters |
---|
695 | |
---|
696 | IF (active_flags%ok_sechiba) THEN |
---|
697 | !- |
---|
698 | DO j =1 ,nvm ! Loop over PFTs |
---|
699 | !- |
---|
700 | ! Vegetation structure - sechiba |
---|
701 | !- |
---|
702 | rveg_pft(j) = rveg_mtc(pft_to_mtc(j)) |
---|
703 | !- |
---|
704 | ! Evapotranspiration - sechiba |
---|
705 | !- |
---|
706 | rstruct_const(j) = rstruct_const_mtc(pft_to_mtc(j)) |
---|
707 | kzero(j) = kzero_mtc(pft_to_mtc(j)) |
---|
708 | !- |
---|
709 | ! Water - sechiba |
---|
710 | !- |
---|
711 | wmax_veg(j) = wmax_veg_mtc(pft_to_mtc(j)) |
---|
712 | throughfall_by_pft(j) = throughfall_by_mtc(pft_to_mtc(j)) |
---|
713 | !- |
---|
714 | ! Albedo - sechiba |
---|
715 | !- |
---|
716 | snowa_ini(j) = snowa_ini_mtc(pft_to_mtc(j)) |
---|
717 | snowa_dec(j) = snowa_dec_mtc(pft_to_mtc(j)) |
---|
718 | alb_leaf_vis(j) = alb_leaf_vis_mtc(pft_to_mtc(j)) |
---|
719 | alb_leaf_nir(j) = alb_leaf_nir_mtc(pft_to_mtc(j)) |
---|
720 | ENDDO ! Loop over PFTs |
---|
721 | !- |
---|
722 | ENDIF !(active_flags%ok_sechiba) |
---|
723 | |
---|
724 | ! 1.3 For stomate parameters |
---|
725 | |
---|
726 | IF (active_flags%ok_stomate) THEN |
---|
727 | !- |
---|
728 | DO j = 1,nvm ! Loop over PFTs |
---|
729 | !- |
---|
730 | ! Vegetation structure - stomate |
---|
731 | !- |
---|
732 | leaf_tab(j) = leaf_tab_mtc(pft_to_mtc(j)) |
---|
733 | sla(j) = sla_mtc(pft_to_mtc(j)) |
---|
734 | !- |
---|
735 | ! Photosynthesis |
---|
736 | !- |
---|
737 | vcmax_opt(j) = vcmax_opt_mtc(pft_to_mtc(j)) |
---|
738 | vjmax_opt(j) = vjmax_opt_mtc(pft_to_mtc(j)) |
---|
739 | tphoto_min_a(j) = tphoto_min_a_mtc(pft_to_mtc(j)) |
---|
740 | tphoto_min_b(j) = tphoto_min_b_mtc(pft_to_mtc(j)) |
---|
741 | tphoto_min_c(j) = tphoto_min_c_mtc(pft_to_mtc(j)) |
---|
742 | tphoto_opt_a(j) = tphoto_opt_a_mtc(pft_to_mtc(j)) |
---|
743 | tphoto_opt_b(j) = tphoto_opt_b_mtc(pft_to_mtc(j)) |
---|
744 | tphoto_opt_c(j) = tphoto_opt_c_mtc(pft_to_mtc(j)) |
---|
745 | tphoto_max_a(j) = tphoto_max_a_mtc(pft_to_mtc(j)) |
---|
746 | tphoto_max_b(j) = tphoto_max_b_mtc(pft_to_mtc(j)) |
---|
747 | tphoto_max_c(j) = tphoto_max_c_mtc(pft_to_mtc(j)) |
---|
748 | !- |
---|
749 | ! Respiration - stomate |
---|
750 | !- |
---|
751 | maint_resp_slope_c(j) = maint_resp_slope_c_mtc(pft_to_mtc(j)) |
---|
752 | maint_resp_slope_b(j) = maint_resp_slope_b_mtc(pft_to_mtc(j)) |
---|
753 | maint_resp_slope_a(j) = maint_resp_slope_a_mtc(pft_to_mtc(j)) |
---|
754 | cm_zero_leaf(j) = cm_zero_leaf_mtc(pft_to_mtc(j)) |
---|
755 | cm_zero_sapabove(j) = cm_zero_sapabove_mtc(pft_to_mtc(j)) |
---|
756 | cm_zero_sapbelow(j) = cm_zero_sapbelow_mtc(pft_to_mtc(j)) |
---|
757 | cm_zero_heartabove(j) = cm_zero_heartabove_mtc(pft_to_mtc(j)) |
---|
758 | cm_zero_heartbelow(j) = cm_zero_heartbelow_mtc(pft_to_mtc(j)) |
---|
759 | cm_zero_root(j) = cm_zero_root_mtc(pft_to_mtc(j)) |
---|
760 | cm_zero_fruit(j) = cm_zero_fruit_mtc(pft_to_mtc(j)) |
---|
761 | cm_zero_carbres(j) = cm_zero_carbres_mtc(pft_to_mtc(j)) |
---|
762 | !- |
---|
763 | ! Fire - stomate |
---|
764 | !- |
---|
765 | flam(j) = flam_mtc(pft_to_mtc(j)) |
---|
766 | resist(j) = resist_mtc(pft_to_mtc(j)) |
---|
767 | !- |
---|
768 | ! Flux - LUC |
---|
769 | !- |
---|
770 | coeff_lcchange_1(j) = coeff_lcchange_1_mtc(pft_to_mtc(j)) |
---|
771 | coeff_lcchange_10(j) = coeff_lcchange_10_mtc(pft_to_mtc(j)) |
---|
772 | coeff_lcchange_100(j) = coeff_lcchange_100_mtc(pft_to_mtc(j)) |
---|
773 | !- |
---|
774 | ! Phenology |
---|
775 | !- |
---|
776 | ! |
---|
777 | ! 1 .Stomate |
---|
778 | ! |
---|
779 | lai_max(j) = lai_max_mtc(pft_to_mtc(j)) |
---|
780 | pheno_model(j) = pheno_model_mtc(pft_to_mtc(j)) |
---|
781 | pheno_type(j) = pheno_type_mtc(pft_to_mtc(j)) |
---|
782 | ! |
---|
783 | ! 2. Leaf Onset |
---|
784 | ! |
---|
785 | pheno_gdd_crit_c(j) = pheno_gdd_crit_c_mtc(pft_to_mtc(j)) |
---|
786 | pheno_gdd_crit_b(j) = pheno_gdd_crit_b_mtc(pft_to_mtc(j)) |
---|
787 | pheno_gdd_crit_a(j) = pheno_gdd_crit_a_mtc(pft_to_mtc(j)) |
---|
788 | ngd_crit(j) = ngd_crit_mtc(pft_to_mtc(j)) |
---|
789 | ncdgdd_temp(j) = ncdgdd_temp_mtc(pft_to_mtc(j)) |
---|
790 | hum_frac(j) = hum_frac_mtc(pft_to_mtc(j)) |
---|
791 | lowgpp_time(j) = lowgpp_time_mtc(pft_to_mtc(j)) |
---|
792 | hum_min_time(j) = hum_min_time_mtc(pft_to_mtc(j)) |
---|
793 | tau_sap(j) = tau_sap_mtc(pft_to_mtc(j)) |
---|
794 | tau_fruit(j) = tau_fruit_mtc(pft_to_mtc(j)) |
---|
795 | ecureuil(j) = ecureuil_mtc(pft_to_mtc(j)) |
---|
796 | alloc_min(j) = alloc_min_mtc(pft_to_mtc(j)) |
---|
797 | alloc_max(j) = alloc_max_mtc(pft_to_mtc(j)) |
---|
798 | demi_alloc(j) = demi_alloc_mtc(pft_to_mtc(j)) |
---|
799 | leaflife_tab(j) = leaflife_mtc(pft_to_mtc(j)) |
---|
800 | ! |
---|
801 | ! 3. Senescence |
---|
802 | ! |
---|
803 | leaffall(j) = leaffall_mtc(pft_to_mtc(j)) |
---|
804 | leafagecrit(j) = leafagecrit_mtc(pft_to_mtc(j)) |
---|
805 | senescence_type(j) = senescence_type_mtc(pft_to_mtc(j)) |
---|
806 | senescence_hum(j) = senescence_hum_mtc(pft_to_mtc(j)) |
---|
807 | nosenescence_hum(j) = nosenescence_hum_mtc(pft_to_mtc(j)) |
---|
808 | max_turnover_time(j) = max_turnover_time_mtc(pft_to_mtc(j)) |
---|
809 | min_turnover_time(j) = min_turnover_time_mtc(pft_to_mtc(j)) |
---|
810 | min_leaf_age_for_senescence(j) = min_leaf_age_for_senescence_mtc(pft_to_mtc(j)) |
---|
811 | senescence_temp_c(j) = senescence_temp_c_mtc(pft_to_mtc(j)) |
---|
812 | senescence_temp_b(j) = senescence_temp_b_mtc(pft_to_mtc(j)) |
---|
813 | senescence_temp_a(j) = senescence_temp_a_mtc(pft_to_mtc(j)) |
---|
814 | !- |
---|
815 | ! DGVM |
---|
816 | residence_time(j) = residence_time_mtc(pft_to_mtc(j)) |
---|
817 | tmin_crit(j) = tmin_crit_mtc(pft_to_mtc(j)) |
---|
818 | tcm_crit(j) = tcm_crit_mtc(pft_to_mtc(j)) |
---|
819 | ENDDO ! Loop over PFTs |
---|
820 | !- |
---|
821 | ENDIF !(active_flags%ok_stomate) |
---|
822 | |
---|
823 | END SUBROUTINE pft_parameters_init |
---|
824 | ! |
---|
825 | != |
---|
826 | ! |
---|
827 | |
---|
828 | !! ================================================================================================================================ |
---|
829 | !! SUBROUTINE : pft_parameters_alloc |
---|
830 | !! |
---|
831 | !>\BRIEF This subroutine allocates memory needed for the PFT parameters |
---|
832 | !! in function of the flags activated. |
---|
833 | !! |
---|
834 | !! DESCRIPTION : None |
---|
835 | !! |
---|
836 | !! RECENT CHANGE(S): None |
---|
837 | !! |
---|
838 | !! MAIN OUTPUT VARIABLE(S): None |
---|
839 | !! |
---|
840 | !! REFERENCE(S) : None |
---|
841 | !! |
---|
842 | !! FLOWCHART : None |
---|
843 | !! \n |
---|
844 | !_ ================================================================================================================================ |
---|
845 | |
---|
846 | SUBROUTINE pft_parameters_alloc(active_flags) |
---|
847 | |
---|
848 | IMPLICIT NONE |
---|
849 | |
---|
850 | !! 0. Variables and parameters declaration |
---|
851 | |
---|
852 | !! 0.1 Input variables |
---|
853 | |
---|
854 | TYPE(control_type),INTENT(in) :: active_flags !! What parts of the code are activated ? (true/false) |
---|
855 | |
---|
856 | !! 0.4 Local variables |
---|
857 | |
---|
858 | LOGICAL :: l_error !! Diagnostic boolean for error allocation (true/false) |
---|
859 | INTEGER :: ier !! Return value for memory allocation (0-N, unitless) |
---|
860 | |
---|
861 | !_ ================================================================================================================================ |
---|
862 | |
---|
863 | ! |
---|
864 | ! 1. Parameters used anytime |
---|
865 | ! |
---|
866 | l_error = .FALSE. |
---|
867 | !- |
---|
868 | ALLOCATE(pft_to_mtc(nvm),stat=ier) |
---|
869 | l_error = l_error .OR. (ier .NE. 0) |
---|
870 | ALLOCATE(PFT_name(nvm),stat=ier) |
---|
871 | l_error = l_error .OR. (ier .NE. 0) |
---|
872 | ALLOCATE(height_presc(nvm),stat=ier) |
---|
873 | l_error = l_error .OR. (ier .NE. 0) |
---|
874 | ALLOCATE(is_tree(nvm),stat=ier) |
---|
875 | l_error = l_error .OR. (ier .NE. 0) |
---|
876 | ALLOCATE(natural(nvm),stat=ier) |
---|
877 | l_error = l_error .OR. (ier .NE. 0) |
---|
878 | ALLOCATE(is_c4(nvm),stat=ier) |
---|
879 | l_error = l_error .OR. (ier .NE. 0) |
---|
880 | ALLOCATE(gsslope(nvm),stat=ier) |
---|
881 | l_error = l_error .OR. (ier .NE. 0) |
---|
882 | ALLOCATE(gsoffset(nvm),stat=ier) |
---|
883 | l_error = l_error .OR. (ier .NE. 0) |
---|
884 | ALLOCATE(humcste(nvm),stat=ier) |
---|
885 | l_error = l_error .OR. (ier .NE. 0) |
---|
886 | ALLOCATE(ext_coeff(nvm),stat=ier) |
---|
887 | l_error = l_error .OR. (ier .NE. 0) |
---|
888 | !- |
---|
889 | ALLOCATE(veget_ori_fixed_test_1(nvm),stat=ier) |
---|
890 | l_error = l_error .OR. (ier .NE. 0) |
---|
891 | ALLOCATE(llaimax(nvm),stat=ier) |
---|
892 | l_error = l_error .OR. (ier .NE. 0) |
---|
893 | ALLOCATE(llaimin(nvm),stat=ier) |
---|
894 | l_error = l_error .OR. (ier .NE. 0) |
---|
895 | ALLOCATE(type_of_lai(nvm),stat=ier) |
---|
896 | l_error = l_error .OR. (ier .NE. 0) |
---|
897 | ALLOCATE(vcmax_fix(nvm),stat=ier) |
---|
898 | l_error = l_error .OR. (ier .NE. 0) |
---|
899 | ALLOCATE(vjmax_fix(nvm),stat=ier) |
---|
900 | l_error = l_error .OR. (ier .NE. 0) |
---|
901 | ALLOCATE(co2_tmin_fix(nvm),stat=ier) |
---|
902 | l_error = l_error .OR. (ier .NE. 0) |
---|
903 | ALLOCATE(co2_topt_fix(nvm),stat=ier) |
---|
904 | l_error = l_error .OR. (ier .NE. 0) |
---|
905 | ALLOCATE(co2_tmax_fix(nvm),stat=ier) |
---|
906 | l_error = l_error .OR. (ier .NE. 0) |
---|
907 | !- |
---|
908 | ALLOCATE(pref_soil_veg_sand(nvm),stat=ier) |
---|
909 | l_error = l_error .OR. (ier .NE. 0) |
---|
910 | ALLOCATE(pref_soil_veg_loan(nvm),stat=ier) |
---|
911 | l_error = l_error .OR. (ier .NE. 0) |
---|
912 | ALLOCATE(pref_soil_veg_clay(nvm),stat=ier) |
---|
913 | l_error = l_error .OR. (ier .NE. 0) |
---|
914 | ALLOCATE(pref_soil_veg(nvm,nstm),stat=ier) |
---|
915 | l_error = l_error .OR. (ier .NE. 0) |
---|
916 | |
---|
917 | !- |
---|
918 | ALLOCATE(is_deciduous(nvm),stat=ier) |
---|
919 | l_error = l_error .OR. (ier .NE. 0) |
---|
920 | ALLOCATE(is_evergreen(nvm),stat=ier) |
---|
921 | l_error = l_error .OR. (ier .NE. 0) |
---|
922 | ALLOCATE(is_c3(nvm),stat=ier) |
---|
923 | l_error = l_error .OR. (ier .NE. 0) |
---|
924 | |
---|
925 | IF(l_error) THEN |
---|
926 | STOP 'pft _alloc : error in memory allocation of pft parameters' |
---|
927 | ENDIF |
---|
928 | |
---|
929 | ! |
---|
930 | ! 2. Parameters used if ok_sechiba only |
---|
931 | ! |
---|
932 | IF (active_flags%ok_sechiba) THEN |
---|
933 | l_error=.FALSE. |
---|
934 | !- |
---|
935 | ALLOCATE(rstruct_const(nvm),stat=ier) |
---|
936 | l_error = l_error .OR. (ier .NE. 0) |
---|
937 | ALLOCATE(kzero(nvm),stat=ier) |
---|
938 | l_error = l_error .OR. (ier .NE. 0) |
---|
939 | ALLOCATE (rveg_pft(nvm),stat=ier) |
---|
940 | l_error = l_error .OR. (ier .NE. 0) |
---|
941 | ALLOCATE(wmax_veg(nvm),stat=ier) |
---|
942 | l_error = l_error .OR. (ier .NE. 0) |
---|
943 | ALLOCATE(throughfall_by_pft(nvm),stat=ier) |
---|
944 | l_error = l_error .OR. (ier .NE. 0) |
---|
945 | !- |
---|
946 | ALLOCATE(snowa_ini(nvm),stat=ier) |
---|
947 | l_error = l_error .OR. (ier .NE. 0) |
---|
948 | ALLOCATE(snowa_dec(nvm),stat=ier) |
---|
949 | l_error = l_error .OR. (ier .NE. 0) |
---|
950 | ALLOCATE(alb_leaf_vis(nvm),stat=ier) |
---|
951 | l_error = l_error .OR. (ier .NE. 0) |
---|
952 | ALLOCATE(alb_leaf_nir(nvm),stat=ier) |
---|
953 | l_error = l_error .OR. (ier .NE. 0) |
---|
954 | ENDIF !(active_flags%ok_sechiba) |
---|
955 | |
---|
956 | IF(l_error) THEN |
---|
957 | WRITE(numout,*) 'Error in memory allocation for sechiba pft parameters' |
---|
958 | ENDIF |
---|
959 | |
---|
960 | ! |
---|
961 | ! 3. Parameters used if ok_stomate only |
---|
962 | ! |
---|
963 | IF (active_flags%ok_stomate) THEN |
---|
964 | l_error=.FALSE. |
---|
965 | !- |
---|
966 | ALLOCATE(leaf_tab(nvm),stat=ier) |
---|
967 | l_error = l_error .OR. (ier .NE. 0) |
---|
968 | ALLOCATE(sla(nvm),stat=ier) |
---|
969 | l_error = l_error .OR. (ier .NE. 0) |
---|
970 | ALLOCATE(vcmax_opt(nvm),stat=ier) |
---|
971 | l_error = l_error .OR. (ier .NE. 0) |
---|
972 | ALLOCATE(vjmax_opt(nvm),stat=ier) |
---|
973 | l_error = l_error .OR. (ier .NE. 0) |
---|
974 | !- |
---|
975 | ALLOCATE(tphoto_min_a(nvm),stat=ier) |
---|
976 | l_error = l_error .OR. (ier .NE. 0) |
---|
977 | ALLOCATE(tphoto_min_b(nvm),stat=ier) |
---|
978 | l_error = l_error .OR. (ier .NE. 0) |
---|
979 | ALLOCATE(tphoto_min_c(nvm),stat=ier) |
---|
980 | l_error = l_error .OR. (ier .NE. 0) |
---|
981 | ALLOCATE(tphoto_opt_a(nvm),stat=ier) |
---|
982 | l_error = l_error .OR. (ier .NE. 0) |
---|
983 | ALLOCATE(tphoto_opt_b(nvm),stat=ier) |
---|
984 | l_error = l_error .OR. (ier .NE. 0) |
---|
985 | ALLOCATE(tphoto_opt_c(nvm),stat=ier) |
---|
986 | l_error = l_error .OR. (ier .NE. 0) |
---|
987 | ALLOCATE(tphoto_max_a(nvm),stat=ier) |
---|
988 | l_error = l_error .OR. (ier .NE. 0) |
---|
989 | ALLOCATE(tphoto_max_b(nvm),stat=ier) |
---|
990 | l_error = l_error .OR. (ier .NE. 0) |
---|
991 | ALLOCATE(tphoto_max_c(nvm),stat=ier) |
---|
992 | l_error = l_error .OR. (ier .NE. 0) |
---|
993 | !- |
---|
994 | ALLOCATE(pheno_gdd_crit_c(nvm),stat=ier) |
---|
995 | l_error = l_error .OR. (ier .NE. 0) |
---|
996 | ALLOCATE(pheno_gdd_crit_b(nvm),stat=ier) |
---|
997 | l_error = l_error .OR. (ier .NE. 0) |
---|
998 | ALLOCATE(pheno_gdd_crit_a(nvm),stat=ier) |
---|
999 | l_error = l_error .OR. (ier .NE. 0) |
---|
1000 | ALLOCATE(pheno_gdd_crit(nvm,3),stat=ier) |
---|
1001 | l_error = l_error .OR. (ier .NE. 0) |
---|
1002 | ALLOCATE(ngd_crit(nvm),stat=ier) |
---|
1003 | l_error = l_error .OR. (ier .NE. 0) |
---|
1004 | ALLOCATE(ncdgdd_temp(nvm),stat=ier) |
---|
1005 | l_error = l_error .OR. (ier .NE. 0) |
---|
1006 | ALLOCATE(hum_frac(nvm),stat=ier) |
---|
1007 | l_error = l_error .OR. (ier .NE. 0) |
---|
1008 | ALLOCATE(lowgpp_time(nvm),stat=ier) |
---|
1009 | l_error = l_error .OR. (ier .NE. 0) |
---|
1010 | ALLOCATE(hum_min_time(nvm),stat=ier) |
---|
1011 | l_error = l_error .OR. (ier .NE. 0) |
---|
1012 | ALLOCATE(tau_sap(nvm),stat=ier) |
---|
1013 | l_error = l_error .OR. (ier .NE. 0) |
---|
1014 | ALLOCATE(tau_fruit(nvm),stat=ier) |
---|
1015 | l_error = l_error .OR. (ier .NE. 0) |
---|
1016 | ALLOCATE(ecureuil(nvm),stat=ier) |
---|
1017 | l_error = l_error .OR. (ier .NE. 0) |
---|
1018 | ALLOCATE(alloc_min(nvm),stat=ier) |
---|
1019 | l_error = l_error .OR. (ier .NE. 0) |
---|
1020 | ALLOCATE(alloc_max(nvm),stat=ier) |
---|
1021 | l_error = l_error .OR. (ier .NE. 0) |
---|
1022 | ALLOCATE(demi_alloc(nvm),stat=ier) |
---|
1023 | l_error = l_error .OR. (ier .NE. 0) |
---|
1024 | !- |
---|
1025 | ALLOCATE(maint_resp_slope(nvm,3),stat=ier) |
---|
1026 | l_error = l_error .OR. (ier .NE. 0) |
---|
1027 | ALLOCATE(maint_resp_slope_c(nvm),stat=ier) |
---|
1028 | l_error = l_error .OR. (ier .NE. 0) |
---|
1029 | ALLOCATE(maint_resp_slope_b(nvm),stat=ier) |
---|
1030 | l_error = l_error .OR. (ier .NE. 0) |
---|
1031 | ALLOCATE(maint_resp_slope_a(nvm),stat=ier) |
---|
1032 | l_error = l_error .OR. (ier .NE. 0) |
---|
1033 | ALLOCATE(coeff_maint_zero(nvm,nparts),stat=ier) |
---|
1034 | l_error = l_error .OR. (ier .NE. 0) |
---|
1035 | ALLOCATE(cm_zero_leaf(nvm),stat=ier) |
---|
1036 | l_error = l_error .OR. (ier .NE. 0) |
---|
1037 | ALLOCATE(cm_zero_sapabove(nvm),stat=ier) |
---|
1038 | l_error = l_error .OR. (ier .NE. 0) |
---|
1039 | ALLOCATE(cm_zero_sapbelow(nvm),stat=ier) |
---|
1040 | l_error = l_error .OR. (ier .NE. 0) |
---|
1041 | ALLOCATE(cm_zero_heartabove(nvm),stat=ier) |
---|
1042 | l_error = l_error .OR. (ier .NE. 0) |
---|
1043 | ALLOCATE(cm_zero_heartbelow(nvm),stat=ier) |
---|
1044 | l_error = l_error .OR. (ier .NE. 0) |
---|
1045 | ALLOCATE(cm_zero_root(nvm),stat=ier) |
---|
1046 | l_error = l_error .OR. (ier .NE. 0) |
---|
1047 | ALLOCATE(cm_zero_fruit(nvm),stat=ier) |
---|
1048 | l_error = l_error .OR. (ier .NE. 0) |
---|
1049 | ALLOCATE(cm_zero_carbres(nvm),stat=ier) |
---|
1050 | l_error = l_error .OR. (ier .NE. 0) |
---|
1051 | !- |
---|
1052 | ALLOCATE(flam(nvm),stat=ier) |
---|
1053 | l_error = l_error .OR. (ier .NE. 0) |
---|
1054 | ALLOCATE(resist(nvm),stat=ier) |
---|
1055 | l_error = l_error .OR. (ier .NE. 0) |
---|
1056 | !- |
---|
1057 | ALLOCATE(coeff_lcchange_1(nvm),stat=ier) |
---|
1058 | l_error = l_error .OR. (ier .NE. 0) |
---|
1059 | ALLOCATE(coeff_lcchange_10(nvm),stat=ier) |
---|
1060 | l_error = l_error .OR. (ier .NE. 0) |
---|
1061 | ALLOCATE(coeff_lcchange_100(nvm),stat=ier) |
---|
1062 | l_error = l_error .OR. (ier .NE. 0) |
---|
1063 | !- |
---|
1064 | ALLOCATE(lai_max(nvm),stat=ier) |
---|
1065 | l_error = l_error .OR. (ier .NE. 0) |
---|
1066 | ALLOCATE(pheno_model(nvm),stat=ier) |
---|
1067 | l_error = l_error .OR. (ier .NE. 0) |
---|
1068 | ALLOCATE(pheno_type(nvm),stat=ier) |
---|
1069 | l_error = l_error .OR. (ier .NE. 0) |
---|
1070 | !- |
---|
1071 | ALLOCATE(leaffall(nvm),stat=ier) |
---|
1072 | l_error = l_error .OR. (ier .NE. 0) |
---|
1073 | ALLOCATE(leafagecrit(nvm),stat=ier) |
---|
1074 | l_error = l_error .OR. (ier .NE. 0) |
---|
1075 | ALLOCATE(senescence_type(nvm),stat=ier) |
---|
1076 | l_error = l_error .OR. (ier .NE. 0) |
---|
1077 | ALLOCATE(senescence_hum(nvm),stat=ier) |
---|
1078 | l_error = l_error .OR. (ier .NE. 0) |
---|
1079 | ALLOCATE(nosenescence_hum(nvm),stat=ier) |
---|
1080 | l_error = l_error .OR. (ier .NE. 0) |
---|
1081 | ALLOCATE(max_turnover_time(nvm),stat=ier) |
---|
1082 | l_error = l_error .OR. (ier .NE. 0) |
---|
1083 | ALLOCATE(min_turnover_time(nvm),stat=ier) |
---|
1084 | l_error = l_error .OR. (ier .NE. 0) |
---|
1085 | ALLOCATE(min_leaf_age_for_senescence(nvm),stat=ier) |
---|
1086 | l_error = l_error .OR. (ier .NE. 0) |
---|
1087 | ALLOCATE(senescence_temp_c(nvm),stat=ier) |
---|
1088 | l_error = l_error .OR. (ier .NE. 0) |
---|
1089 | ALLOCATE(senescence_temp_b(nvm),stat=ier) |
---|
1090 | l_error = l_error .OR. (ier .NE. 0) |
---|
1091 | ALLOCATE(senescence_temp_a(nvm),stat=ier) |
---|
1092 | l_error = l_error .OR. (ier .NE. 0) |
---|
1093 | ALLOCATE(senescence_temp(nvm,3),stat=ier) |
---|
1094 | l_error = l_error .OR. (ier .NE. 0) |
---|
1095 | !- |
---|
1096 | ALLOCATE(residence_time(nvm),stat=ier) |
---|
1097 | l_error = l_error .OR. (ier .NE. 0) |
---|
1098 | ALLOCATE(tmin_crit(nvm),stat=ier) |
---|
1099 | l_error = l_error .OR. (ier .NE. 0) |
---|
1100 | ALLOCATE(tcm_crit(nvm),stat=ier) |
---|
1101 | l_error = l_error .OR. (ier .NE. 0) |
---|
1102 | !- |
---|
1103 | ALLOCATE(lai_initmin(nvm),stat=ier) |
---|
1104 | l_error = l_error .OR. (ier .NE. 0) |
---|
1105 | ALLOCATE(tree(nvm),stat=ier) |
---|
1106 | l_error = l_error .OR. (ier .NE. 0) |
---|
1107 | ALLOCATE(bm_sapl(nvm,nparts),stat=ier) |
---|
1108 | l_error = l_error .OR. (ier .NE. 0) |
---|
1109 | ALLOCATE(migrate(nvm),stat=ier) |
---|
1110 | l_error = l_error .OR. (ier .NE. 0) |
---|
1111 | ALLOCATE(maxdia(nvm),stat=ier) |
---|
1112 | l_error = l_error .OR. (ier .NE. 0) |
---|
1113 | ALLOCATE(cn_sapl(nvm),stat=ier) |
---|
1114 | l_error = l_error .OR. (ier .NE. 0) |
---|
1115 | ALLOCATE(leaf_timecst(nvm),stat=ier) |
---|
1116 | l_error = l_error .OR. (ier .NE. 0) |
---|
1117 | ALLOCATE(leaflife_tab(nvm),stat=ier) |
---|
1118 | l_error = l_error .OR. (ier .NE. 0) |
---|
1119 | ENDIF ! (active_flags%ok_stomate) |
---|
1120 | ! |
---|
1121 | IF (l_error) THEN |
---|
1122 | STOP 'pft _alloc : error in memory allocation of stomate pft parameters' |
---|
1123 | ENDIF |
---|
1124 | |
---|
1125 | END SUBROUTINE pft_parameters_alloc |
---|
1126 | ! |
---|
1127 | != |
---|
1128 | ! |
---|
1129 | |
---|
1130 | !! ================================================================================================================================ |
---|
1131 | !! SUBROUTINE : config_pft_parameters |
---|
1132 | !! |
---|
1133 | !>\BRIEF This subroutine will read the imposed values for the global pft |
---|
1134 | !! parameters (sechiba + stomate). It is not called if IMPOSE_PARAM is set to NO. |
---|
1135 | !! |
---|
1136 | !! DESCRIPTION : None |
---|
1137 | !! |
---|
1138 | !! RECENT CHANGE(S): None |
---|
1139 | !! |
---|
1140 | !! MAIN OUTPUT VARIABLE(S): None |
---|
1141 | !! |
---|
1142 | !! REFERENCE(S) : None |
---|
1143 | !! |
---|
1144 | !! FLOWCHART : None |
---|
1145 | !! \n |
---|
1146 | !_ ================================================================================================================================ |
---|
1147 | |
---|
1148 | SUBROUTINE config_pft_parameters |
---|
1149 | |
---|
1150 | IMPLICIT NONE |
---|
1151 | |
---|
1152 | !! 0. Variables and parameters declaration |
---|
1153 | |
---|
1154 | !! 0.4 Local variable |
---|
1155 | |
---|
1156 | LOGICAL, SAVE :: first_call = .TRUE. !! To keep first call trace (true/false) |
---|
1157 | |
---|
1158 | !_ ================================================================================================================================ |
---|
1159 | |
---|
1160 | IF (first_call) THEN |
---|
1161 | |
---|
1162 | !- |
---|
1163 | ! Vegetation structure |
---|
1164 | !- |
---|
1165 | ! |
---|
1166 | !Config Key = SECHIBA_LAI |
---|
1167 | !Config Desc = laimax for maximum lai(see also type of lai interpolation) |
---|
1168 | !Config if = OK_SECHIBA or IMPOSE_VEG |
---|
1169 | !Config Def = 0., 8., 8., 4., 4.5, 4.5, 4., 4.5, 4., 2., 2., 2., 2. |
---|
1170 | !Config Help = Maximum values of lai used for interpolation of the lai map |
---|
1171 | !Config Units = [m^2/m^2] |
---|
1172 | CALL getin_p('SECHIBA_LAI',llaimax) |
---|
1173 | ! |
---|
1174 | !Config Key = LLAIMIN |
---|
1175 | !Config Desc = laimin for minimum lai(see also type of lai interpolation) |
---|
1176 | !Config if = OK_SECHIBA or IMPOSE_VEG |
---|
1177 | !Config Def = 0., 8., 0., 4., 4.5, 0., 4., 0., 0., 0., 0., 0., 0. |
---|
1178 | !Config Help = Minimum values of lai used for interpolation of the lai map |
---|
1179 | !Config Units = [m^2/m^2] |
---|
1180 | CALL getin_p('LLAIMIN',llaimin) |
---|
1181 | ! |
---|
1182 | !Config Key = SLOWPROC_HEIGHT |
---|
1183 | !Config Desc = prescribed height of vegetation |
---|
1184 | !Config if = OK_SECHIBA |
---|
1185 | !Config Def = 0., 30., 30., 20., 20., 20., 15., 15., 15., .5, .6, 1., 1. |
---|
1186 | !Config Help = |
---|
1187 | !Config Units = [m] |
---|
1188 | CALL getin_p('SLOWPROC_HEIGHT', height_presc) |
---|
1189 | ! |
---|
1190 | !Config Key = TYPE_OF_LAI |
---|
1191 | !Config Desc = Type of behaviour of the LAI evolution algorithm |
---|
1192 | !Config if = OK_SECHIBA |
---|
1193 | !Config Def = inter, inter, inter, inter, inter, inter, inter, inter, inter, inter, inter, inter, inter |
---|
1194 | !Config Help = |
---|
1195 | !Config Units = [-] |
---|
1196 | CALL getin_p('TYPE_OF_LAI',type_of_lai) |
---|
1197 | ! |
---|
1198 | !Config Key = IS_TREE |
---|
1199 | !Config Desc = Is the vegetation type a tree ? |
---|
1200 | !Config if = OK_SECHIBA |
---|
1201 | !Config Def = n, y, y, y, y, y, y, y, y, n, n, n, n |
---|
1202 | !Config Help = |
---|
1203 | !Config Units = [BOOLEAN] |
---|
1204 | CALL getin_p('IS_TREE',is_tree) |
---|
1205 | ! |
---|
1206 | !Config Key = NATURAL |
---|
1207 | !Config Desc = natural? |
---|
1208 | !Config if = OK_SECHIBA, OK_STOMATE |
---|
1209 | !Config Def = y, y, y, y, y, y, y, y, y, y, y, n, n |
---|
1210 | !Config Help = |
---|
1211 | !Config Units = [BOOLEAN] |
---|
1212 | CALL getin_p('NATURAL',natural) |
---|
1213 | ! |
---|
1214 | !Config Key = IS_DECIDUOUS |
---|
1215 | !Config Desc = is PFT deciduous ? |
---|
1216 | !Config if = OK_STOMATE |
---|
1217 | !Config Def = n, n, y, n, n, y, n, y, y, n, n, n, n |
---|
1218 | !Config Help = |
---|
1219 | !Config Units = [BOOLEAN] |
---|
1220 | CALL getin_p('IS_DECIDUOUS',is_deciduous) |
---|
1221 | ! |
---|
1222 | !Config Key = IS_EVERGREEN |
---|
1223 | !Config Desc = is PFT evergreen ? |
---|
1224 | !Config if = OK_STOMATE |
---|
1225 | !Config Def = n, y, n, y, y, n, y, n, n, n, n, n, n |
---|
1226 | !Config Help = |
---|
1227 | !Config Units = [BOOLEAN] |
---|
1228 | CALL getin_p('IS_EVERGREEN',is_evergreen) |
---|
1229 | ! |
---|
1230 | !Config Key = IS_C3 |
---|
1231 | !Config Desc = is PFT C3 ? |
---|
1232 | !Config if = OK_STOMATE |
---|
1233 | !Config Def = n, n, n, n, n, n, n, n, n, n, y, n, y, n |
---|
1234 | !Config Help = |
---|
1235 | !Config Units = [BOOLEAN] |
---|
1236 | CALL getin_p('IS_C3',is_c3) |
---|
1237 | |
---|
1238 | !- |
---|
1239 | ! Photosynthesis |
---|
1240 | !- |
---|
1241 | ! |
---|
1242 | !Config Key = IS_C4 |
---|
1243 | !Config Desc = flag for C4 vegetation types |
---|
1244 | !Config if = OK_SECHIBA or OK_STOMATE |
---|
1245 | !Config Def = n, n, n, n, n, n, n, n, n, n, n, y, n, y |
---|
1246 | !Config Help = |
---|
1247 | !Config Units = [BOOLEAN] |
---|
1248 | CALL getin_p('IS_C4',is_c4) |
---|
1249 | ! |
---|
1250 | !Config Key = GSSLOPE |
---|
1251 | !Config Desc = Slope of the gs/A relation (Ball & al.) |
---|
1252 | !Config if = OK_CO2 |
---|
1253 | !Config Def = 0., 9., 9., 9., 9., 9., 9., 9., 9., 9., 3., 9., 3. |
---|
1254 | !Config Help = |
---|
1255 | !Config Units = [-] |
---|
1256 | CALL getin_p('GSSLOPE',gsslope) |
---|
1257 | ! |
---|
1258 | !Config Key = GSOFFSET |
---|
1259 | !Config Desc = intercept of the gs/A relation (Ball & al.) |
---|
1260 | !Config if = OK_CO2 or OK_STOMATE |
---|
1261 | !Config Def = 0.0, 0.01, 0.01, 0.01, 0.01, 0.01, 0.01, 0.01, 0.01, 0.01, 0.03, 0.01, 0.03 |
---|
1262 | !Config Help = |
---|
1263 | !Config Units = [-] |
---|
1264 | CALL getin_p('GSOFFSET',gsoffset) |
---|
1265 | ! |
---|
1266 | !Config Key = VCMAX_FIX |
---|
1267 | !Config Desc = values used for vcmax when STOMATE is not activated |
---|
1268 | !Config if = OK_SECHIBA and NOT(OK_STOMATE) |
---|
1269 | !Config Def = 0., 40., 50., 30., 35., 40.,30., 40., 35., 60., 60., 70., 70. |
---|
1270 | !Config Help = |
---|
1271 | !Config Units = [micromol/m^2/s] |
---|
1272 | CALL getin_p('VCMAX_FIX',vcmax_fix) |
---|
1273 | ! |
---|
1274 | !Config Key = VJMAX_FIX |
---|
1275 | !Config Desc = values used for vjmax when STOMATE is not activated |
---|
1276 | !Config if = OK_SECHIBA and NOT(OK_STOMATE) |
---|
1277 | !Config Def = 0., 80., 100., 60., 70., 80., 60., 80., 70., 120., 120., 140., 140. |
---|
1278 | !Config Help = |
---|
1279 | !Config Units = [micromol/m^2/s] |
---|
1280 | CALL getin_p('VJMAX_FIX',vjmax_fix) |
---|
1281 | ! |
---|
1282 | !Config Key = CO2_TMIN_FIX |
---|
1283 | !Config Desc = values used for photosynthesis tmin when STOMATE is not activated |
---|
1284 | !Config if = OK_SECHIBA and NOT(OK_STOMATE) |
---|
1285 | !Config Def = 0., 2., 2., -4., -3., -2., -4., -4., -4., -5., 6., -5., 6. |
---|
1286 | !Config Help = |
---|
1287 | !Config Units = [C] |
---|
1288 | CALL getin_p('CO2_TMIN_FIX',co2_tmin_fix) |
---|
1289 | ! |
---|
1290 | !Config Key = CO2_TOPT_FIX |
---|
1291 | !Config Desc = values used for photosynthesis topt when STOMATE is not activated |
---|
1292 | !Config if = OK_SECHIBA and NOT(OK_STOMATE) |
---|
1293 | !Config Def = 0., 27.5, 27.5, 17.5, 25., 20.,17.5, 17.5, 17.5, 20., 32.5, 20., 32.5 |
---|
1294 | !Config Help = |
---|
1295 | !Config Units = [C] |
---|
1296 | CALL getin_p('CO2_TOPT_FIX',co2_topt_fix) |
---|
1297 | ! |
---|
1298 | !Config Key = CO2_TMAX_FIX |
---|
1299 | !Config Desc = values used for photosynthesis tmax when STOMATE is not activated |
---|
1300 | !Config if = OK_SECHIBA and NOT(OK_STOMATE) |
---|
1301 | !Config Def = 0., 55., 55., 38., 48., 38.,38., 38., 38., 45., 55., 45., 55. |
---|
1302 | !Config Help = |
---|
1303 | !Config Units = [C] |
---|
1304 | CALL getin_p('CO2_TMAX_FIX',co2_tmax_fix) |
---|
1305 | ! |
---|
1306 | !Config Key = EXT_COEFF |
---|
1307 | !Config Desc = extinction coefficient of the Monsi&Seaki relationship (1953) |
---|
1308 | !Config if = OK_SECHIBA or OK_STOMATE |
---|
1309 | !Config Def = .5, .5, .5, .5, .5, .5, .5, .5, .5, .5, .5, .5, .5 |
---|
1310 | !Config Help = |
---|
1311 | !Config Units = [-] |
---|
1312 | CALL getin_p('EXT_COEFF',ext_coeff) |
---|
1313 | |
---|
1314 | !- |
---|
1315 | ! Water-hydrology - sechiba |
---|
1316 | !- |
---|
1317 | ! |
---|
1318 | !Config Key = HYDROL_HUMCSTE |
---|
1319 | !Config Desc = Root profile |
---|
1320 | !Config Def = 5., .8, .8, 1., .8, .8, 1., 1., .8, 4., 4., 4., 4. |
---|
1321 | !Config if = OK_SECHIBA |
---|
1322 | !Config Help = Default values were defined for 2 meters soil depth. |
---|
1323 | !Config For 4 meters soil depth, you may use those ones : |
---|
1324 | !Config 5., .4, .4, 1., .8, .8, 1., 1., .8, 4., 1., 4., 1. |
---|
1325 | !Config Units = [m] |
---|
1326 | CALL getin_p('HYDROL_HUMCSTE',humcste) |
---|
1327 | |
---|
1328 | !- |
---|
1329 | ! Soil - vegetation |
---|
1330 | !- |
---|
1331 | ! |
---|
1332 | !Config Key = PREF_SOIL_VEG_SAND |
---|
1333 | !Config Desc = Table which contains the correlation between the soil types and vegetation type |
---|
1334 | !Config if = OK_SECHIBA or OK_STOMATE |
---|
1335 | !Config Def = 1, 3, 3, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2 |
---|
1336 | !Config Help = first layer of the soil |
---|
1337 | !Config Units = [-] |
---|
1338 | CALL getin_p('PREF_SOIL_VEG_SAND',pref_soil_veg_sand) |
---|
1339 | ! |
---|
1340 | !Config Key = PREF_SOIL_VEG_LOAN |
---|
1341 | !Config Desc = Table which contains the correlation between the soil types and vegetation type |
---|
1342 | !Config if = OK_SECHIBA or OK_STOMATE |
---|
1343 | !Config Def = 2, 2, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3 |
---|
1344 | !Config Help = second layer of the soil |
---|
1345 | !Config Units = [-] |
---|
1346 | CALL getin_p('PREF_SOIL_VEG_LOAN',pref_soil_veg_loan) |
---|
1347 | ! |
---|
1348 | !Config Key = PREF_SOIL_VEG_CLAY |
---|
1349 | !Config Desc = Table which contains the correlation between the soil types and vegetation type |
---|
1350 | !Config if = OK_SECHIBA or OK_STOMATE |
---|
1351 | !Config Def = 3, 1, 1, 1, 1, 1 ,1 ,1 ,1 ,1 ,1 ,1, 1 |
---|
1352 | !Config Help = third layer of the soil |
---|
1353 | !Config Units = [-] |
---|
1354 | CALL getin_p('PREF_SOIL_VEG_CLAY',pref_soil_veg_clay) |
---|
1355 | |
---|
1356 | first_call = .FALSE. |
---|
1357 | |
---|
1358 | ENDIF !(first_call) |
---|
1359 | |
---|
1360 | END SUBROUTINE config_pft_parameters |
---|
1361 | ! |
---|
1362 | != |
---|
1363 | ! |
---|
1364 | |
---|
1365 | !! ================================================================================================================================ |
---|
1366 | !! SUBROUTINE : config_sechiba_pft_parameters |
---|
1367 | !! |
---|
1368 | !>\BRIEF This subroutine will read the imposed values for the sechiba pft |
---|
1369 | !! parameters. It is not called if IMPOSE_PARAM is set to NO. |
---|
1370 | !! |
---|
1371 | !! DESCRIPTION : None |
---|
1372 | !! |
---|
1373 | !! RECENT CHANGE(S): None |
---|
1374 | !! |
---|
1375 | !! MAIN OUTPUT VARIABLE(S): None |
---|
1376 | !! |
---|
1377 | !! REFERENCE(S) : None |
---|
1378 | !! |
---|
1379 | !! FLOWCHART : None |
---|
1380 | !! \n |
---|
1381 | !_ ================================================================================================================================ |
---|
1382 | |
---|
1383 | SUBROUTINE config_sechiba_pft_parameters |
---|
1384 | |
---|
1385 | IMPLICIT NONE |
---|
1386 | |
---|
1387 | !! 0. Variables and parameters declaration |
---|
1388 | |
---|
1389 | !! 0.4 Local variable |
---|
1390 | |
---|
1391 | LOGICAL, SAVE :: first_call = .TRUE. !! To keep first call trace (true/false) |
---|
1392 | |
---|
1393 | !_ ================================================================================================================================ |
---|
1394 | |
---|
1395 | IF (first_call) THEN |
---|
1396 | !- |
---|
1397 | ! Evapotranspiration - sechiba |
---|
1398 | !- |
---|
1399 | ! |
---|
1400 | !Config Key = RSTRUCT_CONST |
---|
1401 | !Config Desc = Structural resistance |
---|
1402 | !Config if = OK_SECHIBA |
---|
1403 | !Config Def = 0.0, 25.0, 25.0, 25.0, 25.0, 25.0, 25.0, 25.0, 25.0, 2.5, 2.0, 2.0, 2.0 |
---|
1404 | !Config Help = |
---|
1405 | !Config Units = [s/m] |
---|
1406 | CALL getin_p('RSTRUCT_CONST',rstruct_const) |
---|
1407 | ! |
---|
1408 | !Config Key = KZERO |
---|
1409 | !Config Desc = A vegetation dependent constant used in the calculation of the surface resistance. |
---|
1410 | !Config if = OK_SECHIBA |
---|
1411 | !Config Def = 0.0, 12.E-5, 12.E-5, 12.e-5, 12.e-5, 25.e-5, 12.e-5,25.e-5, 25.e-5, 30.e-5, 30.e-5, 30.e-5, 30.e-5 |
---|
1412 | !Config Help = |
---|
1413 | !Config Units = [kg/m^2/s] |
---|
1414 | CALL getin_p('KZERO',kzero) |
---|
1415 | ! |
---|
1416 | ! Ajouts Nathalie - le 28 Mars 2006 - sur conseils Fred Hourdin |
---|
1417 | ! |
---|
1418 | !Config Key = RVEG_PFT |
---|
1419 | !Config Desc = Artificial parameter to increase or decrease canopy resistance. |
---|
1420 | !Config if = OK_SECHIBA |
---|
1421 | !Config Def = 1., 1., 1., 1., 1., 1., 1., 1., 1., 1., 1., 1., 1. |
---|
1422 | !Config Help = This parameter is set by PFT. |
---|
1423 | !Config Units = [-] |
---|
1424 | CALL getin_p('RVEG_PFT',rveg_pft) |
---|
1425 | |
---|
1426 | !- |
---|
1427 | ! Water-hydrology - sechiba |
---|
1428 | !- |
---|
1429 | ! |
---|
1430 | !Config Key = WMAX_VEG |
---|
1431 | !Config Desc = Maximum field capacity for each of the vegetations (Temporary): max quantity of water |
---|
1432 | !Config if = OK_SECHIBA |
---|
1433 | !Config Def = 150., 150., 150., 150., 150., 150., 150.,150., 150., 150., 150., 150., 150. |
---|
1434 | !Config Help = |
---|
1435 | !Config Units = [kg/m^3] |
---|
1436 | CALL getin_p('WMAX_VEG',wmax_veg) |
---|
1437 | ! |
---|
1438 | !Config Key = PERCENT_THROUGHFALL_PFT |
---|
1439 | !Config Desc = Percent by PFT of precip that is not intercepted by the canopy |
---|
1440 | !Config if = OK_SECHIBA OR OK_CWRR |
---|
1441 | !Config Def = 30. 30. 30. 30. 30. 30. 30. 30. 30. 30. 30. 30. 30. |
---|
1442 | !Config Help = During one rainfall event, PERCENT_THROUGHFALL_PFT% of the incident rainfall |
---|
1443 | !Config will get directly to the ground without being intercepted, for each PFT. |
---|
1444 | !Config Units = [%] |
---|
1445 | CALL getin_p('PERCENT_THROUGHFALL_PFT',throughfall_by_pft) |
---|
1446 | |
---|
1447 | !- |
---|
1448 | ! Albedo - sechiba |
---|
1449 | !- |
---|
1450 | ! |
---|
1451 | !Config Key = SNOWA_INI |
---|
1452 | !Config Desc = Initial snow albedo value for each vegetation type as it will be used in condveg_snow |
---|
1453 | !Config if = OK_SECHIBA |
---|
1454 | !Config Def = 0.35, 0., 0., 0.14, 0.14, 0.14, 0.14, 0.14, 0.14, 0.18, 0.18, 0.18, 0.18 |
---|
1455 | !Config Help = Values are from the Thesis of S. Chalita (1992) |
---|
1456 | !Config Units = [-] |
---|
1457 | CALL getin_p('SNOWA_INI',snowa_ini) |
---|
1458 | ! |
---|
1459 | !Config Key = SNOWA_DEC |
---|
1460 | !Config Desc = Decay rate of snow albedo value for each vegetation type as it will be used in condveg_snow |
---|
1461 | !Config if = OK_SECHIBA |
---|
1462 | !Config Def = 0.45, 0., 0., 0.06, 0.06, 0.11, 0.06, 0.11, 0.11, 0.52,0.52, 0.52, 0.52 |
---|
1463 | !Config Help = Values are from the Thesis of S. Chalita (1992) |
---|
1464 | !Config Units = [-] |
---|
1465 | CALL getin_p('SNOWA_DEC',snowa_dec) |
---|
1466 | ! |
---|
1467 | !Config Key = ALB_LEAF_VIS |
---|
1468 | !Config Desc = leaf albedo of vegetation type, visible albedo |
---|
1469 | !Config if = OK_SECHIBA |
---|
1470 | !Config Def = .00, .04, .06, .06, .06,.06, .06, .06, .06, .10, .10, .10, .10 |
---|
1471 | !Config Help = |
---|
1472 | !Config Units = [-] |
---|
1473 | CALL getin_p('ALB_LEAF_VIS',alb_leaf_vis) |
---|
1474 | ! |
---|
1475 | !Config Key = ALB_LEAF_NIR |
---|
1476 | !Config Desc = leaf albedo of vegetation type, near infrared albedo |
---|
1477 | !Config if = OK_SECHIBA |
---|
1478 | !Config Def = .00, .20, .22, .22, .22,.22, .22, .22, .22, .30, .30, .30, .30 |
---|
1479 | !Config Help = |
---|
1480 | !Config Units = [-] |
---|
1481 | CALL getin_p('ALB_LEAF_NIR',alb_leaf_nir) |
---|
1482 | |
---|
1483 | first_call = .FALSE. |
---|
1484 | |
---|
1485 | ENDIF !(first_call) |
---|
1486 | |
---|
1487 | END SUBROUTINE config_sechiba_pft_parameters |
---|
1488 | ! |
---|
1489 | != |
---|
1490 | ! |
---|
1491 | |
---|
1492 | !! ================================================================================================================================ |
---|
1493 | !! SUBROUTINE : config_stomate_pft_parameters |
---|
1494 | !! |
---|
1495 | !>\BRIEF This subroutine will read the imposed values for the stomate pft |
---|
1496 | !! parameters. It is not called if IMPOSE_PARAM is set to NO. |
---|
1497 | !! |
---|
1498 | !! DESCRIPTION : None |
---|
1499 | !! |
---|
1500 | !! RECENT CHANGE(S): None |
---|
1501 | !! |
---|
1502 | !! MAIN OUTPUT VARIABLE(S): None |
---|
1503 | !! |
---|
1504 | !! REFERENCE(S) : None |
---|
1505 | !! |
---|
1506 | !! FLOWCHART : None |
---|
1507 | !! \n |
---|
1508 | !_ ================================================================================================================================ |
---|
1509 | |
---|
1510 | SUBROUTINE config_stomate_pft_parameters |
---|
1511 | |
---|
1512 | IMPLICIT NONE |
---|
1513 | |
---|
1514 | !! 0. Variables and parameters declaration |
---|
1515 | |
---|
1516 | !! 0.4 Local variable |
---|
1517 | |
---|
1518 | LOGICAL, SAVE :: first_call = .TRUE. !! To keep first call trace (true/false) |
---|
1519 | |
---|
1520 | !_ ================================================================================================================================ |
---|
1521 | |
---|
1522 | IF (first_call) THEN |
---|
1523 | |
---|
1524 | !- |
---|
1525 | ! Vegetation structure |
---|
1526 | !- |
---|
1527 | ! |
---|
1528 | !Config Key = LEAF_TAB |
---|
1529 | !Config Desc = leaf type : 1=broad leaved tree, 2=needle leaved tree, 3=grass 4=bare ground |
---|
1530 | !Config if = OK_STOMATE |
---|
1531 | !Config Def = 4, 1, 1, 2, 1, 1, 2, 1, 2, 3, 3, 3, 3 |
---|
1532 | !Config Help = |
---|
1533 | !Config Units = [-] |
---|
1534 | CALL getin_p('LEAF_TAB',leaf_tab) |
---|
1535 | ! |
---|
1536 | !Config Key = SLA |
---|
1537 | !Config Desc = specif leaf area |
---|
1538 | !Config if = OK_STOMATE |
---|
1539 | !Config Def = 1.5E-2, 1.53E-2, 2.6E-2, 9.26E-3, 2E-2, 2.6E-2, 9.26E-3, 2.6E-2, 1.9E-2, 2.6E-2, 2.6E-2, 2.6E-2, 2.6E-2 |
---|
1540 | !Config Help = |
---|
1541 | !Config Units = [m^2/gC] |
---|
1542 | CALL getin_p('SLA',sla) |
---|
1543 | |
---|
1544 | !- |
---|
1545 | ! Photosynthesis |
---|
1546 | !- |
---|
1547 | ! |
---|
1548 | !Config Key = VCMAX_OPT |
---|
1549 | !Config Desc = Maximum rate of carboxylation |
---|
1550 | !Config if = OK_STOMATE |
---|
1551 | !Config Def = undef, 65., 65., 35., 45., 55., 35., 45., 35., 70., 70., 70., 70. |
---|
1552 | !Config Help = |
---|
1553 | !Config Units = [micromol/m^2/s] |
---|
1554 | CALL getin_p('VCMAX_OPT',vcmax_opt) |
---|
1555 | ! |
---|
1556 | !Config Key = VJMAX_OPT |
---|
1557 | !Config Desc = Maximum rate of RUbp regeneration |
---|
1558 | !Config if = OK_STOMATE |
---|
1559 | !Config Def = undef, 130., 130., 70., 80., 110., 70., 90., 70., 160., 160., 200., 200. |
---|
1560 | !Config Help = |
---|
1561 | !Config Units = [micromol/m^2/s] |
---|
1562 | CALL getin_p('VJMAX_OPT',vjmax_opt) |
---|
1563 | ! |
---|
1564 | !Config Key = TPHOTO_MIN_A |
---|
1565 | !Config Desc = minimum photosynthesis temperature, constant a of ax^2+bx+c (deg C), tabulated |
---|
1566 | !Config if = OK_STOMATE |
---|
1567 | !Config Def = undef, 0., 0., 0., 0., 0., 0., 0., 0., 0.0025, 0., 0., 0. |
---|
1568 | !Config Help = |
---|
1569 | !Config Units = [-] |
---|
1570 | CALL getin_p('TPHOTO_MIN_A',tphoto_min_a) |
---|
1571 | ! |
---|
1572 | !Config Key = TPHOTO_MIN_B |
---|
1573 | !Config Desc = minimum photosynthesis temperature, constant b of ax^2+bx+c (deg C), tabulated |
---|
1574 | !Config if = OK_STOMATE |
---|
1575 | !Config Def = undef, 0., 0., 0., 0., 0., 0., 0., 0., 0.1, 0.,0.,0. |
---|
1576 | !Config Help = |
---|
1577 | !Config Units = [-] |
---|
1578 | CALL getin_p('TPHOTO_MIN_B',tphoto_min_b) |
---|
1579 | ! |
---|
1580 | !Config Key = TPHOTO_MIN_C |
---|
1581 | !Config Desc = minimum photosynthesis temperature, constant c of ax^2+bx+c (deg C), tabulated |
---|
1582 | !Config if = OK_STOMATE |
---|
1583 | !Config Def = undef, 2., 2., -4., -3.,-2.,-4., -4.,-4.,-3.25, 13.,-5.,13. |
---|
1584 | !Config Help = |
---|
1585 | !Config Units = [-] |
---|
1586 | CALL getin_p('TPHOTO_MIN_C',tphoto_min_c) |
---|
1587 | ! |
---|
1588 | !Config Key = TPHOTO_OPT_A |
---|
1589 | !Config Desc = optimum photosynthesis temperature, constant a of ax^2+bx+c (deg C), tabulated |
---|
1590 | !Config if = OK_STOMATE |
---|
1591 | !Config Def = undef, 0., 0., 0., 0., 0., 0., 0., 0., 0.0025, 0., 0., 0. |
---|
1592 | !Config Help = |
---|
1593 | !Config Units = [-] |
---|
1594 | CALL getin_p('TPHOTO_OPT_A',tphoto_opt_a) |
---|
1595 | ! |
---|
1596 | !Config Key = TPHOTO_OPT_B |
---|
1597 | !Config Desc = optimum photosynthesis temperature, constant b of ax^2+bx+c (deg C), tabulated |
---|
1598 | !Config if = OK_STOMATE |
---|
1599 | !Config Def = undef, 0., 0., 0., 0., 0., 0., 0., 0., 0.25, 0., 0., 0. |
---|
1600 | !Config Help = |
---|
1601 | !Config Units = [-] |
---|
1602 | CALL getin_p('TPHOTO_OPT_B',tphoto_opt_b) |
---|
1603 | ! |
---|
1604 | !Config Key = TPHOTO_OPT_C |
---|
1605 | !Config Desc = optimum photosynthesis temperature, constant c of ax^2+bx+c (deg C), tabulated |
---|
1606 | !Config if = OK_STOMATE |
---|
1607 | !Config Def = undef, 37., 37., 25., 32., 26., 25., 25., 25., 27.25, 36., 30., 36. |
---|
1608 | !Config Help = |
---|
1609 | !Config Units = [-] |
---|
1610 | CALL getin_p('TPHOTO_OPT_C',tphoto_opt_c) |
---|
1611 | ! |
---|
1612 | !Config Key = TPHOTO_MAX_A |
---|
1613 | !Config Desc = maximum photosynthesis temperature, constant a of ax^2+bx+c (deg C), tabulated |
---|
1614 | !Config if = OK_STOMATE |
---|
1615 | !Config Def = undef, 0., 0., 0., 0., 0., 0., 0., 0., 0.00375, 0., 0., 0. |
---|
1616 | !Config Help = |
---|
1617 | !Config Units = [-] |
---|
1618 | CALL getin_p('TPHOTO_MAX_A',tphoto_max_a) |
---|
1619 | ! |
---|
1620 | !Config Key = TPHOTO_MAX_B |
---|
1621 | !Config Desc = maximum photosynthesis temperature, constant b of ax^2+bx+c (deg C), tabulated |
---|
1622 | !Config if = OK_STOMATE |
---|
1623 | !Config Def = undef, 0., 0., 0., 0., 0., 0., 0., 0.,0.35, 0., 0., 0. |
---|
1624 | !Config Help = |
---|
1625 | !Config Units = [-] |
---|
1626 | CALL getin_p('TPHOTO_MAX_B',tphoto_max_b) |
---|
1627 | ! |
---|
1628 | !Config Key = TPHOTO_MAX_C |
---|
1629 | !Config Desc = maximum photosynthesis temperature, constant c of ax^2+bx+c (deg C), tabulated |
---|
1630 | !Config if = OK_STOMATE |
---|
1631 | !Config Def = undef, 55., 55.,38., 48.,38.,38., 38., 38., 41.125, 55., 45., 55. |
---|
1632 | !Config Help = |
---|
1633 | !Config Units = [-] |
---|
1634 | CALL getin_p('TPHOTO_MAX_C',tphoto_max_c) |
---|
1635 | |
---|
1636 | !- |
---|
1637 | ! Respiration - stomate |
---|
1638 | !- |
---|
1639 | ! |
---|
1640 | !Config Key = MAINT_RESP_SLOPE_C |
---|
1641 | !Config Desc = slope of maintenance respiration coefficient (1/K), constant c of aT^2+bT+c , tabulated |
---|
1642 | !Config if = OK_STOMATE |
---|
1643 | !Config Def = undef, .12, .12, .16, .16, .16, .16, .16, .16, .16, .12, .16, .12 |
---|
1644 | !Config Help = |
---|
1645 | !Config Units = [-] |
---|
1646 | CALL getin_p('MAINT_RESP_SLOPE_C',maint_resp_slope_c) |
---|
1647 | ! |
---|
1648 | !Config Key = MAINT_RESP_SLOPE_B |
---|
1649 | !Config Desc = slope of maintenance respiration coefficient (1/K), constant b of aT^2+bT+c , tabulated |
---|
1650 | !Config if = OK_STOMATE |
---|
1651 | !Config Def = undef, .0, .0, .0, .0, .0, .0, .0, .0, -.00133, .0, -.00133, .0 |
---|
1652 | !Config Help = |
---|
1653 | !Config Units = [-] |
---|
1654 | CALL getin_p('MAINT_RESP_SLOPE_B',maint_resp_slope_b) |
---|
1655 | ! |
---|
1656 | !Config Key = MAINT_RESP_SLOPE_A |
---|
1657 | !Config Desc = slope of maintenance respiration coefficient (1/K), constant a of aT^2+bT+c , tabulated |
---|
1658 | !Config if = OK_STOMATE |
---|
1659 | !Config Def = undef, .0, .0, .0, .0, .0, .0, .0, .0, .0, .0, .0, .0 |
---|
1660 | !Config Help = |
---|
1661 | !Config Units = [-] |
---|
1662 | CALL getin_p('MAINT_RESP_SLOPE_A',maint_resp_slope_a) |
---|
1663 | ! |
---|
1664 | !Config Key = CM_ZERO_LEAF |
---|
1665 | !Config Desc = maintenance respiration coefficient at 0 deg C, for leaves, tabulated |
---|
1666 | !Config if = OK_STOMATE |
---|
1667 | !Config Def = undef, 2.35E-3, 2.62E-3, 1.01E-3, 2.35E-3, 2.62E-3, 1.01E-3,2.62E-3, 2.05E-3, 2.62E-3, 2.62E-3, 2.62E-3, 2.62E-3 |
---|
1668 | !Config Help = |
---|
1669 | !Config Units = [g/g/day] |
---|
1670 | CALL getin_p('CM_ZERO_LEAF',cm_zero_leaf) |
---|
1671 | ! |
---|
1672 | !Config Key = CM_ZERO_SAPABOVE |
---|
1673 | !Config Desc = maintenance respiration coefficient at 0 deg C,for sapwood above, tabulated |
---|
1674 | !Config if = OK_STOMATE |
---|
1675 | !Config Def = undef, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4 |
---|
1676 | !Config Help = |
---|
1677 | !Config Units = [g/g/day] |
---|
1678 | CALL getin_p('CM_ZERO_SAPABOVE',cm_zero_sapabove) |
---|
1679 | ! |
---|
1680 | !Config Key = CM_ZERO_SAPBELOW |
---|
1681 | !Config Desc = maintenance respiration coefficient at 0 deg C, for sapwood below, tabulated |
---|
1682 | !Config if = OK_STOMATE |
---|
1683 | !Config Def = undef, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4 |
---|
1684 | !Config Help = |
---|
1685 | !Config Units = [g/g/day] |
---|
1686 | CALL getin_p('CM_ZERO_SAPBELOW',cm_zero_sapbelow) |
---|
1687 | ! |
---|
1688 | !Config Key = CM_ZERO_HEARTABOVE |
---|
1689 | !Config Desc = maintenance respiration coefficient at 0 deg C, for heartwood above, tabulated |
---|
1690 | !Config if = OK_STOMATE |
---|
1691 | !Config Def = undef, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0. |
---|
1692 | !Config Help = |
---|
1693 | !Config Units = [g/g/day] |
---|
1694 | CALL getin_p('CM_ZERO_HEARTABOVE',cm_zero_heartabove) |
---|
1695 | ! |
---|
1696 | !Config Key = CM_ZERO_HEARTBELOW |
---|
1697 | !Config Desc = maintenance respiration coefficient at 0 deg C,for heartwood below, tabulated |
---|
1698 | !Config if = OK_STOMATE |
---|
1699 | !Config Def = undef, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0. |
---|
1700 | !Config Help = |
---|
1701 | !Config Units = [g/g/day] |
---|
1702 | CALL getin_p('CM_ZERO_HEARTBELOW',cm_zero_heartbelow) |
---|
1703 | ! |
---|
1704 | !Config Key = CM_ZERO_ROOT |
---|
1705 | !Config Desc = maintenance respiration coefficient at 0 deg C, for roots, tabulated |
---|
1706 | !Config if = OK_STOMATE |
---|
1707 | !Config Def = undef,1.67E-3, 1.67E-3, 1.67E-3, 1.67E-3, 1.67E-3, 1.67E-3,1.67E-3, 1.67E-3, 1.67E-3, 1.67E-3, 1.67E-3, 1.67E-3 |
---|
1708 | !Config Help = |
---|
1709 | !Config Units = [g/g/day] |
---|
1710 | CALL getin_p('CM_ZERO_ROOT',cm_zero_root) |
---|
1711 | ! |
---|
1712 | !Config Key = CM_ZERO_FRUIT |
---|
1713 | !Config Desc = maintenance respiration coefficient at 0 deg C, for fruits, tabulated |
---|
1714 | !Config if = OK_STOMATE |
---|
1715 | !Config Def = undef, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4,1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4 |
---|
1716 | !Config Help = |
---|
1717 | !Config Units = [g/g/day] |
---|
1718 | CALL getin_p('CM_ZERO_FRUIT',cm_zero_fruit) |
---|
1719 | ! |
---|
1720 | !Config Key = CM_ZERO_CARBRES |
---|
1721 | !Config Desc = maintenance respiration coefficient at 0 deg C, for carbohydrate reserve, tabulated |
---|
1722 | !Config if = OK_STOMATE |
---|
1723 | !Config Def = undef, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4,1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4, 1.19E-4 |
---|
1724 | !Config Help = |
---|
1725 | !Config Units = [g/g/day] |
---|
1726 | CALL getin_p('CM_ZERO_CARBRES',cm_zero_carbres) |
---|
1727 | |
---|
1728 | !- |
---|
1729 | ! Fire - stomate |
---|
1730 | !- |
---|
1731 | ! |
---|
1732 | !Config Key = FLAM |
---|
1733 | !Config Desc = flamability: critical fraction of water holding capacity |
---|
1734 | !Config if = OK_STOMATE |
---|
1735 | !Config Def = undef, .15, .25, .25, .25, .25, .25, .25, .25, .25, .25, .35, .35 |
---|
1736 | !Config Help = |
---|
1737 | !Config Units = [-] |
---|
1738 | CALL getin_p('FLAM',flam) |
---|
1739 | ! |
---|
1740 | !Config Key = RESIST |
---|
1741 | !Config Desc = fire resistance |
---|
1742 | !Config if = OK_STOMATE |
---|
1743 | !Config Def = undef, .95, .90, .12, .50, .12, .12, .12, .12, .0, .0, .0, .0 |
---|
1744 | !Config Help = |
---|
1745 | !Config Units = [-] |
---|
1746 | CALL getin_p('RESIST',resist) |
---|
1747 | |
---|
1748 | !- |
---|
1749 | ! Flux - LUC |
---|
1750 | !- |
---|
1751 | ! |
---|
1752 | !Config Key = COEFF_LCCHANGE_1 |
---|
1753 | !Config Desc = Coeff of biomass export for the year |
---|
1754 | !Config if = OK_STOMATE |
---|
1755 | !Config Def = undef, 0.597, 0.597, 0.597, 0.597, 0.597, 0.597, 0.597, 0.597, 0.597, 0.597, 0.597, 0.597 |
---|
1756 | !Config Help = |
---|
1757 | !Config Units = [-] |
---|
1758 | CALL getin_p('COEFF_LCCHANGE_1',coeff_lcchange_1) |
---|
1759 | ! |
---|
1760 | !Config Key = COEFF_LCCHANGE_10 |
---|
1761 | !Config Desc = Coeff of biomass export for the decade |
---|
1762 | !Config if = OK_STOMATE |
---|
1763 | !Config Def = undef, 0.403, 0.403, 0.299, 0.299, 0.299, 0.299, 0.299, 0.299, 0.299, 0.403, 0.299, 0.403 |
---|
1764 | !Config Help = |
---|
1765 | !Config Units = [-] |
---|
1766 | CALL getin_p('COEFF_LCCHANGE_10',coeff_lcchange_10) |
---|
1767 | ! |
---|
1768 | !Config Key = COEFF_LCCHANGE_100 |
---|
1769 | !Config Desc = Coeff of biomass export for the century |
---|
1770 | !Config if = OK_STOMATE |
---|
1771 | !Config Def = undef, 0., 0., 0.104, 0.104, 0.104, 0.104, 0.104, 0.104, 0.104, 0., 0.104, 0. |
---|
1772 | !Config Help = |
---|
1773 | !Config Units = [-] |
---|
1774 | CALL getin_p('COEFF_LCCHANGE_100',coeff_lcchange_100) |
---|
1775 | |
---|
1776 | !- |
---|
1777 | ! Phenology |
---|
1778 | !- |
---|
1779 | ! |
---|
1780 | !Config Key = LAI_MAX |
---|
1781 | !Config Desc = maximum LAI, PFT-specific |
---|
1782 | !Config if = OK_STOMATE |
---|
1783 | !Config Def = undef, 7., 7., 5., 5., 5., 4.5, 4.5, 3.0, 2.5, 2.5, 5.,5. |
---|
1784 | !Config Help = |
---|
1785 | !Config Units = [m^2/m^2] |
---|
1786 | CALL getin_p('LAI_MAX',lai_max) |
---|
1787 | ! |
---|
1788 | !Config Key = PHENO_MODEL |
---|
1789 | !Config Desc = which phenology model is used? (tabulated) |
---|
1790 | !Config if = OK_STOMATE |
---|
1791 | !Config Def = none, none, moi, none, none, ncdgdd, none, ncdgdd, ngd, moigdd, moigdd, moigdd, moigdd |
---|
1792 | !Config Help = |
---|
1793 | !Config Units = [-] |
---|
1794 | CALL getin_p('PHENO_MODEL',pheno_model) |
---|
1795 | ! |
---|
1796 | !Config Key = PHENO_TYPE |
---|
1797 | !Config Desc = type of phenology, 0=bare ground 1=evergreen, 2=summergreen, 3=raingreen, 4=perennial |
---|
1798 | !Config if = OK_STOMATE |
---|
1799 | !Config Def = 0, 1, 3, 1, 1, 2, 1, 2, 2, 4, 4, 2, 3 |
---|
1800 | !Config Help = |
---|
1801 | !Config Units = [-] |
---|
1802 | CALL getin_p('PHENO_TYPE',pheno_type) |
---|
1803 | |
---|
1804 | !- |
---|
1805 | ! Phenology : Leaf Onset |
---|
1806 | !- |
---|
1807 | ! |
---|
1808 | !Config Key = PHENO_GDD_CRIT_C |
---|
1809 | !Config Desc = critical gdd, tabulated (C), constant c of aT^2+bT+c |
---|
1810 | !Config if = OK_STOMATE |
---|
1811 | !Config Def = undef, undef, undef, undef, undef, undef, undef, undef, undef, 270., 400., 125., 400. |
---|
1812 | !Config Help = |
---|
1813 | !Config Units = [-] |
---|
1814 | CALL getin_p('PHENO_GDD_CRIT_C',pheno_gdd_crit_c) |
---|
1815 | ! |
---|
1816 | !Config Key = PHENO_GDD_CRIT_B |
---|
1817 | !Config Desc = critical gdd, tabulated (C), constant b of aT^2+bT+c |
---|
1818 | !Config if = OK_STOMATE |
---|
1819 | !Config Def = undef, undef, undef, undef, undef, undef, undef,undef, undef, 6.25, 0., 0., 0. |
---|
1820 | !Config Help = |
---|
1821 | !Config Units = [-] |
---|
1822 | CALL getin_p('PHENO_GDD_CRIT_B',pheno_gdd_crit_b) |
---|
1823 | ! |
---|
1824 | !Config Key = PHENO_GDD_CRIT_A |
---|
1825 | !Config Desc = critical gdd, tabulated (C), constant a of aT^2+bT+c |
---|
1826 | !Config if = OK_STOMATE |
---|
1827 | !Config Def = undef, undef, undef, undef, undef, undef, undef, undef, undef, 0.03125, 0., 0., 0. |
---|
1828 | !Config Help = |
---|
1829 | !Config Units = [-] |
---|
1830 | CALL getin_p('PHENO_GDD_CRIT_A',pheno_gdd_crit_a) |
---|
1831 | ! |
---|
1832 | !Config Key = NGD_CRIT |
---|
1833 | !Config Desc = critical ngd, tabulated. Threshold -5 degrees |
---|
1834 | !Config if = OK_STOMATE |
---|
1835 | !Config Def = undef, undef, undef, undef, undef, undef, undef, 0., undef, undef, undef, undef, undef |
---|
1836 | !Config Help = NGD : Number of Growing Days. |
---|
1837 | !Config Units = [days] |
---|
1838 | CALL getin_p('NGD_CRIT',ngd_crit) |
---|
1839 | ! |
---|
1840 | !Config Key = NCDGDD_TEMP |
---|
1841 | !Config Desc = critical temperature for the ncd vs. gdd function in phenology |
---|
1842 | !Config if = OK_STOMATE |
---|
1843 | !Config Def = undef, undef, undef, undef, undef, 5., undef, 0., undef, undef, undef, undef, undef |
---|
1844 | !Config Help = |
---|
1845 | !Config Units = [C] |
---|
1846 | CALL getin_p('NCDGDD_TEMP', ncdgdd_temp) |
---|
1847 | ! |
---|
1848 | !Config Key = HUM_FRAC |
---|
1849 | !Config Desc = critical humidity (relative to min/max) for phenology |
---|
1850 | !Config if = OK_STOMATE |
---|
1851 | !Config Def = undef, undef, .5, undef, undef, undef, undef, undef, undef, .5, .5, .5,.5 |
---|
1852 | !Config Help = |
---|
1853 | !Config Units = [%] |
---|
1854 | CALL getin_p('HUM_FRAC', hum_frac) |
---|
1855 | ! |
---|
1856 | !Config Key = LOWGPP_TIME |
---|
1857 | !Config Desc = minimum duration of dormance for phenology |
---|
1858 | !Config if = OK_STOMATE |
---|
1859 | !Config Def = undef, undef, 30., undef, undef, 30., undef, 30., 30., 30., 30., 30., 30. |
---|
1860 | !Config Help = |
---|
1861 | !Config Units = [days] |
---|
1862 | CALL getin_p('LOWGPP_TIME', lowgpp_time) |
---|
1863 | ! |
---|
1864 | !Config Key = HUM_MIN_TIME |
---|
1865 | !Config Desc = minimum time elapsed since moisture minimum |
---|
1866 | !Config if = OK_STOMATE |
---|
1867 | !Config Def = undef, undef, 50., undef, undef, undef, undef, undef, undef, 35., 35., 75., 75. |
---|
1868 | !Config Help = |
---|
1869 | !Config Units = [days] |
---|
1870 | CALL getin_p('HUM_MIN_TIME', hum_min_time) |
---|
1871 | ! |
---|
1872 | !Config Key = TAU_SAP |
---|
1873 | !Config Desc = sapwood -> heartwood conversion time |
---|
1874 | !Config if = OK_STOMATE |
---|
1875 | !Config Def = undef, 730., 730., 730., 730., 730., 730., 730., 730., undef, undef, undef, undef |
---|
1876 | !Config Help = |
---|
1877 | !Config Units = [days] |
---|
1878 | CALL getin_p('TAU_SAP',tau_sap) |
---|
1879 | ! |
---|
1880 | !Config Key = TAU_FRUIT |
---|
1881 | !Config Desc = fruit lifetime |
---|
1882 | !Config if = OK_STOMATE |
---|
1883 | !Config Def = undef, 90., 90., 90., 90., 90., 90., 90., 90., undef, undef, undef, undef |
---|
1884 | !Config Help = |
---|
1885 | !Config Units = [days] |
---|
1886 | CALL getin_p('TAU_FRUIT',tau_fruit) |
---|
1887 | ! |
---|
1888 | !Config Key = ECUREUIL |
---|
1889 | !Config Desc = fraction of primary leaf and root allocation put into reserve |
---|
1890 | !Config if = OK_STOMATE |
---|
1891 | !Config Def = undef, .0, 1., .0, .0, 1., .0, 1., 1., 1., 1., 1., 1. |
---|
1892 | !Config Help = |
---|
1893 | !Config Units = [-] |
---|
1894 | CALL getin_p('ECUREUIL',ecureuil) |
---|
1895 | ! |
---|
1896 | !Config Key = ALLOC_MIN |
---|
1897 | !Config Desc = minimum allocation above/below = f(age) - 30/01/04 NV/JO/PF |
---|
1898 | !Config if = OK_STOMATE |
---|
1899 | !Config Def = undef, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2, undef, undef, undef, undef |
---|
1900 | !Config Help = |
---|
1901 | !Config Units = [-] |
---|
1902 | CALL getin_p('ALLOC_MIN',alloc_min) |
---|
1903 | ! |
---|
1904 | !Config Key = ALLOC_MAX |
---|
1905 | !Config Desc = maximum allocation above/below = f(age) - 30/01/04 NV/JO/PF |
---|
1906 | !Config if = OK_STOMATE |
---|
1907 | !Config Def = undef, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, undef, undef, undef, undef |
---|
1908 | !Config Help = |
---|
1909 | !Config Units = [-] |
---|
1910 | CALL getin_p('ALLOC_MAX',alloc_max) |
---|
1911 | ! |
---|
1912 | !Config Key = DEMI_ALLOC |
---|
1913 | !Config Desc = mean allocation above/below = f(age) - 30/01/04 NV/JO/PF |
---|
1914 | !Config if = OK_STOMATE |
---|
1915 | !Config Def = undef, 5., 5., 5., 5., 5., 5., 5., 5., undef, undef, undef, undef |
---|
1916 | !Config Help = |
---|
1917 | !Config Units = [-] |
---|
1918 | CALL getin_p('DEMI_ALLOC',demi_alloc) |
---|
1919 | ! |
---|
1920 | !Config Key = LEAFLIFE_TAB |
---|
1921 | !Config Desc = leaf longevity |
---|
1922 | !Config if = OK_STOMATE |
---|
1923 | !Config Def = undef, .5, 2., .33, 1., 2., .33, 2., 2., 2., 2., 2., 2. |
---|
1924 | !Config Help = |
---|
1925 | !Config Units = [years] |
---|
1926 | CALL getin_p('LEAFLIFE_TAB',leaflife_tab) |
---|
1927 | |
---|
1928 | !- |
---|
1929 | ! Phenology : Senescence |
---|
1930 | !- |
---|
1931 | ! |
---|
1932 | !Config Key = LEAFFALL |
---|
1933 | !Config Desc = length of death of leaves, tabulated |
---|
1934 | !Config if = OK_STOMATE |
---|
1935 | !Config Def = undef, undef, 10., undef, undef, 10., undef, 10., 10., 10., 10., 10., 10. |
---|
1936 | !Config Help = |
---|
1937 | !Config Units = [days] |
---|
1938 | CALL getin_p('LEAFFALL',leaffall) |
---|
1939 | ! |
---|
1940 | !Config Key = LEAFAGECRIT |
---|
1941 | !Config Desc = critical leaf age, tabulated |
---|
1942 | !Config if = OK_STOMATE |
---|
1943 | !Config Def = undef, 730., 180., 910., 730., 180., 910., 180., 180., 120., 120., 90., 90. |
---|
1944 | !Config Help = |
---|
1945 | !Config Units = [days] |
---|
1946 | CALL getin_p('LEAFAGECRIT',leafagecrit) |
---|
1947 | ! |
---|
1948 | !Config Key = SENESCENCE_TYPE |
---|
1949 | !Config Desc = type of senescence, tabulated |
---|
1950 | !Config if = OK_STOMATE |
---|
1951 | !Config Def = none, none, dry, none, none, cold, none, cold, cold, mixed, mixed, mixed, mixed |
---|
1952 | !Config Help = |
---|
1953 | !Config Units = [-] |
---|
1954 | CALL getin_p('SENESCENCE_TYPE',senescence_type) |
---|
1955 | ! |
---|
1956 | !Config Key = SENESCENCE_HUM |
---|
1957 | !Config Desc = critical relative moisture availability for senescence |
---|
1958 | !Config if = OK_STOMATE |
---|
1959 | !Config Def = undef, undef, .3, undef, undef, undef, undef, undef, undef, .2, .2, .3, .2 |
---|
1960 | !Config Help = |
---|
1961 | !Config Units = [-] |
---|
1962 | CALL getin_p('SENESCENCE_HUM',senescence_hum) |
---|
1963 | ! |
---|
1964 | !Config Key = NOSENESCENCE_HUM |
---|
1965 | !Config Desc = relative moisture availability above which there is no humidity-related senescence |
---|
1966 | !Config if = OK_STOMATE |
---|
1967 | !Config Def = undef, undef, .8, undef, undef, undef, undef, undef, undef, .3, .3, .3, .3 |
---|
1968 | !Config Help = |
---|
1969 | !Config Units = [-] |
---|
1970 | CALL getin_p('NOSENESCENCE_HUM',nosenescence_hum) |
---|
1971 | ! |
---|
1972 | !Config Key = MAX_TURNOVER_TIME |
---|
1973 | !Config Desc = maximum turnover time for grasse |
---|
1974 | !Config if = OK_STOMATE |
---|
1975 | !Config Def = undef, undef, undef, undef, undef, undef, undef, undef, undef, 80., 80., 80., 80. |
---|
1976 | !Config Help = |
---|
1977 | !Config Units = [days] |
---|
1978 | CALL getin_p('MAX_TURNOVER_TIME',max_turnover_time) |
---|
1979 | ! |
---|
1980 | !Config Key = MIN_TURNOVER_TIME |
---|
1981 | !Config Desc = minimum turnover time for grasse |
---|
1982 | !Config if = OK_STOMATE |
---|
1983 | !Config Def = undef, undef, undef, undef, undef, undef, undef, undef, undef, 10., 10., 10., 10. |
---|
1984 | !Config Help = |
---|
1985 | !Config Units = [days] |
---|
1986 | CALL getin_p('MIN_TURNOVER_TIME',min_turnover_time) |
---|
1987 | ! |
---|
1988 | !Config Key = MIN_LEAF_AGE_FOR_SENESCENCE |
---|
1989 | !Config Desc = minimum leaf age to allow senescence g |
---|
1990 | !Config if = OK_STOMATE |
---|
1991 | !Config Def = undef, undef, 90., undef, undef, 90., undef, 60., 60., 30., 30., 30., 30. |
---|
1992 | !Config Help = |
---|
1993 | !Config Units = [days] |
---|
1994 | CALL getin_p('MIN_LEAF_AGE_FOR_SENESCENCE',min_leaf_age_for_senescence) |
---|
1995 | ! |
---|
1996 | !Config Key = SENESCENCE_TEMP_C |
---|
1997 | !Config Desc = critical temperature for senescence (C), constant c of aT^2+bT+c, tabulated |
---|
1998 | !Config if = OK_STOMATE |
---|
1999 | !Config Def = undef, undef, undef, undef, undef, 12., undef, 7., 2., -1.375, 5., 5., 10. |
---|
2000 | !Config Help = |
---|
2001 | !Config Units = [-] |
---|
2002 | CALL getin_p('SENESCENCE_TEMP_C',senescence_temp_c) |
---|
2003 | ! |
---|
2004 | !Config Key = SENESCENCE_TEMP_B |
---|
2005 | !Config Desc = critical temperature for senescence (C), constant b of aT^2+bT+c ,tabulated |
---|
2006 | !Config if = OK_STOMATE |
---|
2007 | !Config Def = undef, undef, undef, undef, undef, 0., undef, 0., 0., .1, 0., 0., 0. |
---|
2008 | !Config Help = |
---|
2009 | !Config Units = [-] |
---|
2010 | CALL getin_p('SENESCENCE_TEMP_B',senescence_temp_b) |
---|
2011 | ! |
---|
2012 | !Config Key = SENESCENCE_TEMP_A |
---|
2013 | !Config Desc = critical temperature for senescence (C), constant a of aT^2+bT+c , tabulated |
---|
2014 | !Config if = OK_STOMATE |
---|
2015 | !Config Def = undef, undef, undef, undef, undef, 0., undef, 0., 0.,.00375, 0., 0., 0. |
---|
2016 | !Config Help = |
---|
2017 | !Config Units = [-] |
---|
2018 | CALL getin_p('SENESCENCE_TEMP_A',senescence_temp_a) |
---|
2019 | |
---|
2020 | !- |
---|
2021 | ! DGVM |
---|
2022 | !- |
---|
2023 | ! |
---|
2024 | !Config Key = RESIDENCE_TIME |
---|
2025 | !Config Desc = residence time of trees |
---|
2026 | !Config if = OK_DGVM and NOT(LPJ_GAP_CONST_MORT) |
---|
2027 | !Config Def = undef, 30.0, 30.0, 40.0, 40.0, 40.0, 80.0, 80.0, 80.0, 0.0, 0.0, 0.0, 0.0 |
---|
2028 | !Config Help = |
---|
2029 | !Config Units = [years] |
---|
2030 | CALL getin_p('RESIDENCE_TIME',residence_time) |
---|
2031 | ! |
---|
2032 | !Config Key = TMIN_CRIT |
---|
2033 | !Config Desc = critical tmin, tabulated |
---|
2034 | !Config if = OK_STOMATE |
---|
2035 | !Config Def = undef, 0.0, 0.0, -30.0, -14.0, -30.0, -45.0, -45.0, undef, undef, undef, undef, undef |
---|
2036 | !Config Help = |
---|
2037 | !Config Units = [C] |
---|
2038 | CALL getin_p('TMIN_CRIT',tmin_crit) |
---|
2039 | ! |
---|
2040 | !Config Key = TCM_CRIT |
---|
2041 | !Config Desc = critical tcm, tabulated |
---|
2042 | !Config if = OK_STOMATE |
---|
2043 | !Config Def = undef, undef, undef, 5.0, 15.5, 15.5, -8.0, -8.0, -8.0, undef, undef, undef, undef |
---|
2044 | !Config Help = |
---|
2045 | !Config Units = [C] |
---|
2046 | CALL getin_p('TCM_CRIT',tcm_crit) |
---|
2047 | |
---|
2048 | first_call = .FALSE. |
---|
2049 | |
---|
2050 | ENDIF !(first_call) |
---|
2051 | |
---|
2052 | END SUBROUTINE config_stomate_pft_parameters |
---|
2053 | ! |
---|
2054 | != |
---|
2055 | ! |
---|
2056 | |
---|
2057 | !! ================================================================================================================================ |
---|
2058 | !! SUBROUTINE : pft_parameters_clear |
---|
2059 | !! |
---|
2060 | !>\BRIEF This subroutine deallocates memory at the end of the simulation. |
---|
2061 | !! |
---|
2062 | !! DESCRIPTION : None |
---|
2063 | !! |
---|
2064 | !! RECENT CHANGE(S): None |
---|
2065 | !! |
---|
2066 | !! MAIN OUTPUT VARIABLE(S): None |
---|
2067 | !! |
---|
2068 | !! REFERENCE(S) : None |
---|
2069 | !! |
---|
2070 | !! FLOWCHART : None |
---|
2071 | !! \n |
---|
2072 | !_ ================================================================================================================================ |
---|
2073 | |
---|
2074 | SUBROUTINE pft_parameters_clear |
---|
2075 | |
---|
2076 | l_first_define_pft = .TRUE. |
---|
2077 | |
---|
2078 | IF(ALLOCATED(pft_to_mtc))DEALLOCATE(pft_to_mtc) |
---|
2079 | IF(ALLOCATED(PFT_name))DEALLOCATE(PFT_name) |
---|
2080 | !- |
---|
2081 | IF(ALLOCATED(veget_ori_fixed_test_1))DEALLOCATE(veget_ori_fixed_test_1) |
---|
2082 | IF(ALLOCATED(llaimax))DEALLOCATE(llaimax) |
---|
2083 | IF(ALLOCATED(llaimin))DEALLOCATE(llaimin) |
---|
2084 | IF(ALLOCATED(height_presc))DEALLOCATE(height_presc) |
---|
2085 | IF(ALLOCATED(type_of_lai))DEALLOCATE(type_of_lai) |
---|
2086 | IF(ALLOCATED(is_tree))DEALLOCATE(is_tree) |
---|
2087 | IF(ALLOCATED(natural))DEALLOCATE(natural) |
---|
2088 | !- |
---|
2089 | IF(ALLOCATED(is_deciduous))DEALLOCATE(is_deciduous) |
---|
2090 | IF(ALLOCATED(is_evergreen))DEALLOCATE(is_evergreen) |
---|
2091 | IF(ALLOCATED(is_c3))DEALLOCATE(is_c3) |
---|
2092 | !- |
---|
2093 | IF(ALLOCATED(humcste))DEALLOCATE(humcste) |
---|
2094 | !- |
---|
2095 | IF(ALLOCATED(pref_soil_veg_sand))DEALLOCATE(pref_soil_veg_sand) |
---|
2096 | IF(ALLOCATED(pref_soil_veg_loan))DEALLOCATE(pref_soil_veg_loan) |
---|
2097 | IF(ALLOCATED(pref_soil_veg_clay))DEALLOCATE(pref_soil_veg_clay) |
---|
2098 | IF(ALLOCATED(pref_soil_veg))DEALLOCATE(pref_soil_veg) |
---|
2099 | !- |
---|
2100 | IF(ALLOCATED(is_c4))DEALLOCATE(is_c4) |
---|
2101 | IF(ALLOCATED(gsslope))DEALLOCATE(gsslope) |
---|
2102 | IF(ALLOCATED(gsoffset))DEALLOCATE(gsoffset) |
---|
2103 | IF(ALLOCATED(vcmax_fix))DEALLOCATE(vcmax_fix) |
---|
2104 | IF(ALLOCATED(vjmax_fix))DEALLOCATE(vjmax_fix) |
---|
2105 | IF(ALLOCATED(co2_tmin_fix))DEALLOCATE(co2_tmin_fix) |
---|
2106 | IF(ALLOCATED(co2_topt_fix))DEALLOCATE(co2_topt_fix) |
---|
2107 | IF(ALLOCATED(co2_tmax_fix))DEALLOCATE(co2_tmax_fix) |
---|
2108 | IF(ALLOCATED(ext_coeff))DEALLOCATE(ext_coeff) |
---|
2109 | !- |
---|
2110 | IF(ALLOCATED(rveg_pft))DEALLOCATE(rveg_pft) |
---|
2111 | !- |
---|
2112 | IF(ALLOCATED(rstruct_const))DEALLOCATE(rstruct_const) |
---|
2113 | IF(ALLOCATED(kzero))DEALLOCATE(kzero) |
---|
2114 | !- |
---|
2115 | IF(ALLOCATED(wmax_veg))DEALLOCATE(wmax_veg) |
---|
2116 | IF(ALLOCATED(throughfall_by_pft))DEALLOCATE(throughfall_by_pft) |
---|
2117 | !- |
---|
2118 | IF(ALLOCATED(snowa_ini))DEALLOCATE(snowa_ini) |
---|
2119 | IF(ALLOCATED(snowa_dec))DEALLOCATE(snowa_dec) |
---|
2120 | IF(ALLOCATED(alb_leaf_vis))DEALLOCATE(alb_leaf_vis) |
---|
2121 | IF(ALLOCATED(alb_leaf_nir))DEALLOCATE(alb_leaf_nir) |
---|
2122 | !- |
---|
2123 | IF(ALLOCATED(leaf_tab))DEALLOCATE(leaf_tab) |
---|
2124 | IF(ALLOCATED(sla))DEALLOCATE(sla) |
---|
2125 | !- |
---|
2126 | IF(ALLOCATED(vcmax_opt))DEALLOCATE(vcmax_opt) |
---|
2127 | IF(ALLOCATED(vjmax_opt))DEALLOCATE(vjmax_opt) |
---|
2128 | IF(ALLOCATED(tphoto_min_a))DEALLOCATE(tphoto_min_a) |
---|
2129 | IF(ALLOCATED(tphoto_min_b))DEALLOCATE(tphoto_min_b) |
---|
2130 | IF(ALLOCATED(tphoto_min_c))DEALLOCATE(tphoto_min_c) |
---|
2131 | IF(ALLOCATED(tphoto_opt_a))DEALLOCATE(tphoto_opt_a) |
---|
2132 | IF(ALLOCATED(tphoto_opt_b))DEALLOCATE(tphoto_opt_b) |
---|
2133 | IF(ALLOCATED(tphoto_opt_c))DEALLOCATE(tphoto_opt_c) |
---|
2134 | IF(ALLOCATED(tphoto_max_a))DEALLOCATE(tphoto_max_a) |
---|
2135 | IF(ALLOCATED(tphoto_max_b))DEALLOCATE(tphoto_max_b) |
---|
2136 | IF(ALLOCATED(tphoto_max_c))DEALLOCATE(tphoto_max_c) |
---|
2137 | !- |
---|
2138 | IF(ALLOCATED(maint_resp_slope))DEALLOCATE(maint_resp_slope) |
---|
2139 | IF(ALLOCATED(maint_resp_slope_c))DEALLOCATE(maint_resp_slope_c) |
---|
2140 | IF(ALLOCATED(maint_resp_slope_b))DEALLOCATE(maint_resp_slope_b) |
---|
2141 | IF(ALLOCATED(maint_resp_slope_a))DEALLOCATE(maint_resp_slope_a) |
---|
2142 | IF(ALLOCATED(coeff_maint_zero))DEALLOCATE(coeff_maint_zero) |
---|
2143 | IF(ALLOCATED(cm_zero_leaf))DEALLOCATE(cm_zero_leaf) |
---|
2144 | IF(ALLOCATED(cm_zero_sapabove))DEALLOCATE(cm_zero_sapabove) |
---|
2145 | IF(ALLOCATED(cm_zero_sapbelow))DEALLOCATE(cm_zero_sapbelow) |
---|
2146 | IF(ALLOCATED(cm_zero_heartabove))DEALLOCATE(cm_zero_heartabove) |
---|
2147 | IF(ALLOCATED(cm_zero_heartbelow))DEALLOCATE(cm_zero_heartbelow) |
---|
2148 | IF(ALLOCATED(cm_zero_root))DEALLOCATE(cm_zero_root) |
---|
2149 | IF(ALLOCATED(cm_zero_fruit))DEALLOCATE(cm_zero_fruit) |
---|
2150 | IF(ALLOCATED(cm_zero_carbres))DEALLOCATE(cm_zero_carbres) |
---|
2151 | !- |
---|
2152 | IF(ALLOCATED(flam))DEALLOCATE(flam) |
---|
2153 | IF(ALLOCATED(resist))DEALLOCATE(resist) |
---|
2154 | !- |
---|
2155 | IF(ALLOCATED(coeff_lcchange_1))DEALLOCATE(coeff_lcchange_1) |
---|
2156 | IF(ALLOCATED(coeff_lcchange_10))DEALLOCATE(coeff_lcchange_10) |
---|
2157 | IF(ALLOCATED(coeff_lcchange_100))DEALLOCATE(coeff_lcchange_100) |
---|
2158 | !- |
---|
2159 | IF(ALLOCATED(lai_max)) DEALLOCATE(lai_max) |
---|
2160 | IF(ALLOCATED(pheno_model))DEALLOCATE(pheno_model) |
---|
2161 | IF(ALLOCATED(pheno_type))DEALLOCATE(pheno_type) |
---|
2162 | !- |
---|
2163 | IF(ALLOCATED(pheno_gdd_crit_c))DEALLOCATE(pheno_gdd_crit_c) |
---|
2164 | IF(ALLOCATED(pheno_gdd_crit_b))DEALLOCATE(pheno_gdd_crit_b) |
---|
2165 | IF(ALLOCATED(pheno_gdd_crit_a))DEALLOCATE(pheno_gdd_crit_a) |
---|
2166 | IF(ALLOCATED(pheno_gdd_crit))DEALLOCATE(pheno_gdd_crit) |
---|
2167 | IF(ALLOCATED(ngd_crit))DEALLOCATE(ngd_crit) |
---|
2168 | IF(ALLOCATED(ncdgdd_temp))DEALLOCATE(ncdgdd_temp) |
---|
2169 | IF(ALLOCATED(hum_frac))DEALLOCATE(hum_frac) |
---|
2170 | IF(ALLOCATED(lowgpp_time))DEALLOCATE(lowgpp_time) |
---|
2171 | IF(ALLOCATED(hum_min_time))DEALLOCATE(hum_min_time) |
---|
2172 | IF(ALLOCATED(tau_sap))DEALLOCATE(tau_sap) |
---|
2173 | IF(ALLOCATED(tau_fruit))DEALLOCATE(tau_fruit) |
---|
2174 | IF(ALLOCATED(ecureuil))DEALLOCATE(ecureuil) |
---|
2175 | IF(ALLOCATED(alloc_min))DEALLOCATE(alloc_min) |
---|
2176 | IF(ALLOCATED(alloc_max))DEALLOCATE(alloc_max) |
---|
2177 | IF(ALLOCATED(demi_alloc))DEALLOCATE(demi_alloc) |
---|
2178 | IF(ALLOCATED(leaflife_tab))DEALLOCATE(leaflife_tab) |
---|
2179 | !- |
---|
2180 | IF(ALLOCATED(leaffall))DEALLOCATE(leaffall) |
---|
2181 | IF(ALLOCATED(leafagecrit))DEALLOCATE(leafagecrit) |
---|
2182 | IF(ALLOCATED(senescence_type))DEALLOCATE(senescence_type) |
---|
2183 | IF(ALLOCATED(senescence_hum))DEALLOCATE(senescence_hum) |
---|
2184 | IF(ALLOCATED(nosenescence_hum))DEALLOCATE(nosenescence_hum) |
---|
2185 | IF(ALLOCATED(max_turnover_time))DEALLOCATE(max_turnover_time) |
---|
2186 | IF(ALLOCATED(min_turnover_time))DEALLOCATE(min_turnover_time) |
---|
2187 | IF(ALLOCATED(min_leaf_age_for_senescence))DEALLOCATE(min_leaf_age_for_senescence) |
---|
2188 | !- |
---|
2189 | IF(ALLOCATED(senescence_temp_c))DEALLOCATE(senescence_temp_c) |
---|
2190 | IF(ALLOCATED(senescence_temp_b))DEALLOCATE(senescence_temp_b) |
---|
2191 | IF(ALLOCATED(senescence_temp_a))DEALLOCATE(senescence_temp_a) |
---|
2192 | IF(ALLOCATED(senescence_temp))DEALLOCATE(senescence_temp) |
---|
2193 | !- |
---|
2194 | IF(ALLOCATED(residence_time))DEALLOCATE(residence_time) |
---|
2195 | IF(ALLOCATED(tmin_crit))DEALLOCATE(tmin_crit) |
---|
2196 | IF(ALLOCATED(tcm_crit))DEALLOCATE(tcm_crit) |
---|
2197 | !- |
---|
2198 | IF(ALLOCATED(lai_initmin))DEALLOCATE(lai_initmin) |
---|
2199 | IF(ALLOCATED(tree))DEALLOCATE(tree) |
---|
2200 | IF(ALLOCATED(bm_sapl))DEALLOCATE(bm_sapl) |
---|
2201 | IF(ALLOCATED(migrate))DEALLOCATE(migrate) |
---|
2202 | IF(ALLOCATED(maxdia))DEALLOCATE(maxdia) |
---|
2203 | IF(ALLOCATED(cn_sapl))DEALLOCATE(cn_sapl) |
---|
2204 | IF(ALLOCATED(leaf_timecst))DEALLOCATE(leaf_timecst) |
---|
2205 | |
---|
2206 | END SUBROUTINE pft_parameters_clear |
---|
2207 | |
---|
2208 | END MODULE pft_parameters |
---|