source: branches/publications/ORCHIDEE_CN_CAN_r5698/src_stomate/stomate_data.f90 @ 8076

Last change on this file since 8076 was 5680, checked in by sebastiaan.luyssaert, 6 years ago

DEV: tested for 2 years with and without windthrow for a single pixel. The changes largely restore the calculation of wind damage and contribute to ticket #273

  • Property svn:keywords set to HeadURL Date Author Revision
File size: 25.6 KB
Line 
1! =================================================================================================================================
2! MODULE        : stomate_data
3!
4! CONTACT      : orchidee-help _at_ ipsl.jussieu.fr
5!
6! LICENCE      : IPSL (2006)
7! This software is governed by the CeCILL licence see ORCHIDEE/ORCHIDEE_CeCILL.LIC
8!
9!>\BRIEF         "stomate_data" module defines the values about the PFT parameters. It will print
10!! the values of the parameters for STOMATE in the standard outputs.
11!!
12!!\n DESCRIPTION: None
13!!
14!! RECENT CHANGE(S): Sonke Zaehle: Reich et al, 1992 find no statistically significant differences
15!!                  between broadleaved and coniferous forests, specifically, the assumption that grasses grow
16!!                  needles is not justified. Replacing the function with the one based on Reich et al. 1997.
17!!                  Given that sla=100cm2/gDW at 9 months, sla is:
18!!                  sla=exp(5.615-0.46*ln(leaflon in months))
19!!
20!! REFERENCE(S) : None
21!!
22!! SVN          :
23!! $HeadURL$
24!! $Date$
25!! $Revision$
26!! \n
27!_ ================================================================================================================================
28
29MODULE stomate_data
30
31  ! modules used:
32
33  USE constantes
34  USE pft_parameters
35  USE defprec
36 
37
38  IMPLICIT NONE
39
40 ! INTEGER(i_std), SAVE                         :: printlev_loc   !! local printlev for this module
41!!$OMP THREADPRIVATE(printlev_loc)
42
43  INTEGER(i_std),ALLOCATABLE,SAVE,DIMENSION(:) :: hori_index     !! Move to Horizontal indices
44!$OMP THREADPRIVATE(hori_index)
45
46  INTEGER(i_std),ALLOCATABLE,SAVE,DIMENSION(:) :: horipft_index  !! Horizontal + PFT indices
47!$OMP THREADPRIVATE(horipft_index)
48
49  INTEGER(i_std),ALLOCATABLE,SAVE,DIMENSION(:) :: horican_index  !! Horizontal + canopy levels indices
50!$OMP THREADPRIVATE(horican_index)
51 
52  INTEGER(i_std),ALLOCATABLE,SAVE,DIMENSION(:) :: horicut_index  !! Horizontal + cut times indices
53!$OMP THREADPRIVATE(horicut_index)
54
55  !
56  ! Land cover change
57  !
58  INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: horip_s_index   !! Horizontal + P short indices
59!$OMP THREADPRIVATE(horip_s_index)
60  INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: horip_m_index   !! Horizontal + P medium indices
61!$OMP THREADPRIVATE(horip_m_index)
62  INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: horip_l_index   !! Horizontal + P long indice
63!$OMP THREADPRIVATE(horip_l_index)
64  INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: horip_ss_index  !! Horizontal + P short indices
65!$OMP THREADPRIVATE(horip_ss_index)
66  INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: horip_mm_index  !! Horizontal + P medium indices
67!$OMP THREADPRIVATE(horip_mm_index)
68  INTEGER(i_std), ALLOCATABLE, SAVE, DIMENSION(:) :: horip_ll_index  !! Horizontal + P long indices
69!$OMP THREADPRIVATE(horip_ll_index)
70
71  INTEGER(i_std),SAVE :: itime                 !! time step
72!$OMP THREADPRIVATE(itime)
73  INTEGER(i_std),SAVE :: hist_id_stomate       !! STOMATE history file ID
74!$OMP THREADPRIVATE(hist_id_stomate)
75  INTEGER(i_std),SAVE :: hist_id_stomate_IPCC  !! STOMATE history file ID for IPCC output
76!$OMP THREADPRIVATE(hist_id_stomate_IPCC)
77  INTEGER(i_std),SAVE :: rest_id_stomate       !! STOMATE restart file ID
78!$OMP THREADPRIVATE(rest_id_stomate)
79
80  REAL(r_std),PARAMETER :: adapted_crit = 1. - ( 1. / euler ) !! critical value for being adapted (1-1/e) (unitless)
81  REAL(r_std),PARAMETER :: regenerate_crit = 1. / euler       !! critical value for being regenerative (1/e) (unitless)
82
83
84  ! private & public routines
85
86  PUBLIC data
87
88CONTAINS
89
90!! ================================================================================================================================
91!! SUBROUTINE   : data
92!!
93!>\BRIEF         This routine defines the values of the PFT parameters. It will print the values of the parameters for STOMATE
94!!               in the standard outputs of ORCHIDEE.
95!!
96!! DESCRIPTION : This routine defines PFT parameters. It initializes the pheno_crit structure by tabulated parameters.\n
97!!               Some initializations are done for parameters. The SLA is calculated according *to* Reich et al (1992).\n
98!!               Another formulation by Reich et al(1997) could be used for the computation of the SLA.
99!!               The geographical coordinates might be used for defining some additional parameters
100!!               (e.g. frequency of anthropogenic fires, irrigation of agricultural surfaces, etc.). \n
101!!               For the moment, this possibility is not used. \n
102!!               The specifc leaf area (SLA) is calculated according Reich et al, 1992 by :
103!!               \latexonly
104!!               \input{stomate_data_SLA.tex}
105!!               \endlatexonly
106!!               The sapling (young) biomass for trees and for each compartment of biomass is calculated by :
107!!               \latexonly
108!!               \input{stomate_data_sapl_tree.tex}
109!!               \endlatexonly
110!!               The sapling biomass for grasses and for each compartment of biomass is calculated by :
111!!               \latexonly
112!!               \input{stomate_data_sapl_grass.tex}
113!!               \endlatexonly
114!!               The critical stem diameter is given by the following formula :
115!!               \latexonly
116!!               \input{stomate_data_stem_diameter.tex}
117!!               \endlatexonly
118!!
119!! RECENT CHANGE(S): Sonke Zaehle: Reich et al, 1992 find no statistically significant differences
120!!                  between broadleaved and coniferous forests, specifically, the assumption that grasses grow
121!!                  needles is not justified. Replacing the function with the one based on Reich et al. 1997.
122!!                  Given that sla=100cm2/gDW at 9 months, sla is:
123!!                  sla=exp(5.615-0.46*ln(leaflon in months))
124!!                   \latexonly
125!!                   \input{stomate_data_SLA_Reich_97.tex}
126!!                   \endlatexonly
127!!
128!! MAIN OUTPUT VARIABLE(S):
129!!
130!! REFERENCE(S) :
131!! - Reich PB, Walters MB, Ellsworth DS, (1992), Leaf life-span in relation to leaf, plant and
132!! stand characteristics among diverse ecosystems. Ecological Monographs, Vol 62, pp 365-392.
133!! - Reich PB, Walters MB, Ellsworth DS (1997) From tropics to tundra: global convergence in plant
134!!  functioning. Proc Natl Acad Sci USA, 94:13730 13734
135!!
136!! FLOWCHART    :
137!! \n
138!_ ================================================================================================================================
139
140  SUBROUTINE data (npts, lalo)
141
142
143    !! 0. Variables and parameter declaration
144
145
146    !! 0.1 Input variables
147
148    INTEGER(i_std), INTENT(in)                   :: npts    !! [DISPENSABLE] Domain size (unitless)
149    REAL(r_std),DIMENSION (npts,2), INTENT (in)  :: lalo    !! [DISPENSABLE] Geographical coordinates (latitude,longitude)
150
151    !! 0.4 Local variables
152
153    INTEGER(i_std)                               :: i,j     !! Index (unitless)
154    REAL(r_std)                                  :: alpha   !! alpha's : (unitless)
155    REAL(r_std)                                  :: dia     !! stem diameter (m)
156    REAL(r_std)                                  :: csa_sap !! Crown specific area sapling @tex $(m^2.ind^{-1})$ @endtex
157    REAL(r_std)                                  :: cn_leaf         !! C to N ratio of Leaf pool (gC per gN)
158    REAL(r_std)                                  :: cn_wood         !! C to N ratio of Woody pools (gC per gN)
159    REAL(r_std)                                  :: cn_root         !! C to N ratio of Root pool (gC per gN)
160
161!_ ================================================================================================================================
162
163    IF ( printlev>=1 ) WRITE(numout,*) 'data: PFT characteristics'
164
165    ! Initialize local printlev
166    !printlev_loc=get_printlev('data')
167
168    !- pheno_gdd_crit
169    pheno_gdd_crit(:,1) = pheno_gdd_crit_c(:)
170    pheno_gdd_crit(:,2) = pheno_gdd_crit_b(:)         
171    pheno_gdd_crit(:,3) = pheno_gdd_crit_a(:) 
172    !
173    !- senescence_temp
174    senescence_temp(:,1) = senescence_temp_c(:)
175    senescence_temp(:,2) = senescence_temp_b(:)
176    senescence_temp(:,3) = senescence_temp_a(:)
177
178    !- maint_resp_slope
179    maint_resp_slope(:,1) = maint_resp_slope_c(:)
180    maint_resp_slope(:,2) = maint_resp_slope_b(:)
181    maint_resp_slope(:,3) = maint_resp_slope_a(:)
182
183    !
184    !-LC
185    LC(:,ileaf) = LC_leaf(:) 
186    LC(:,isapabove) = LC_sapabove(:) 
187    LC(:,isapbelow) = LC_sapbelow(:) 
188    LC(:,iheartabove) = LC_heartabove(:) 
189    LC(:,iheartbelow) = LC_heartbelow(:) 
190    LC(:,iroot) = LC_root(:) 
191    LC(:,ifruit) = LC_fruit(:) 
192    LC(:,icarbres) = LC_carbres(:) 
193    LC(:,ilabile) = LC_labile(:) 
194
195
196    DO j = 2,nvm ! Loop over # PFTS
197
198       IF ( printlev_loc >= 2 ) WRITE(numout,'(a,i3,a,a)') '    > PFT#',j,': ', PFT_name(j)
199
200       !
201       ! 1 tree? (true/false)
202       !
203       IF ( printlev_loc >= 2 ) WRITE(numout,*) '       tree: (::is_tree) ', is_tree(j)
204
205       !
206       ! 2 flamability (0-1, unitless)
207       !
208
209       IF ( printlev_loc >= 2 ) WRITE(numout,*) '       litter flamability (::flam) :', flam(j)
210
211       !
212       ! 3 fire resistance (unitless)
213       !
214
215       IF ( printlev_loc >= 2 ) WRITE(numout,*) '       fire resistance (::resist) :', resist(j)
216
217       !
218       ! 4 specific leaf area per mass carbon = 2 * sla / dry mass (m^2.gC^{-1})
219       !
220
221       ! S. Zaehle: Reich et al, 1992 find no statistically significant differences between broadleaved and coniferous
222       ! forests, specifically, the assumption that grasses grow needles is not justified. Replacing the function
223       ! with the one based on Reich et al. 1997. Given that sla=100cm2/gDW at 9 months, sla is:
224       ! sla=exp(5.615-0.46*ln(leaflon in months))
225
226       ! Oct 2010 : sla values are prescribed by values given by N.Viovy
227
228       ! includes conversion from
229       !!       sla(j) = 2. * 1e-4 * EXP(5.615 - 0.46 * log(12./leaflife_tab(j)))
230       !!\latexonly
231       !!\input{stomate_data_SLA.tex}
232       !!\endlatexonly
233!       IF ( leaf_tab(j) .EQ. 2 ) THEN
234!
235!          ! needle leaved tree
236!          sla(j) = 2. * ( 10. ** ( 2.29 - 0.4 * LOG10(12./leaflife_tab(j)) ) ) *1e-4
237!
238!       ELSE
239!
240!          ! broad leaved tree or grass (Reich et al 1992)
241!          sla(j) = 2. * ( 10. ** ( 2.41 - 0.38 * LOG10(12./leaflife_tab(j)) ) ) *1e-4
242!
243!       ENDIF
244
245!!!$      IF ( leaf_tab(j) .EQ. 1 ) THEN
246!!!$
247!!!$        ! broad leaved tree
248!!!$
249!!!$        sla(j) = 2. * ( 10. ** ( 2.41 - 0.38 * LOG10(12./leaflife_tab(j)) ) ) *1e-4
250!!!$
251!!!$      ELSE
252!!!$
253!!!$        ! needle leaved or grass (Reich et al 1992)
254!!!$
255!!!$        sla(j) = 2. * ( 10. ** ( 2.29 - 0.4 * LOG10(12./leaflife_tab(j)) ) ) *1e-4
256!!!$
257!!!$      ENDIF
258!!!$
259!!!$      IF ( ( leaf_tab(j) .EQ. 2 ) .AND. ( pheno_type_tab(j) .EQ. 2 ) ) THEN
260!!!$
261!!!$        ! summergreen needle leaf
262!!!$
263!!!$        sla(j) = 1.25 * sla(j)
264!!!$
265!!!$      ENDIF
266
267       IF ( printlev_loc >= 2 ) WRITE(numout,*) '       specific leaf area (m**2/gC) (::sla):', sla(j), 12./leaflife_tab(j)
268
269       !
270       ! 5 sapling characteristics
271       !
272
273       IF ( is_tree(j) ) THEN
274
275          !> 5.1 trees
276
277          !!\latexonly
278          !!\input{stomate_data_sapl_tree.tex}
279          !!\endlatexonly
280
281          alpha = alpha_tree
282
283          bm_sapl(j,ileaf,icarbon) = &
284               &     ((bm_sapl_leaf(1)*pipe_tune1(j)*(mass_ratio_heart_sap(j) *bm_sapl_leaf(2)*sla(j)/(pi*pipe_k1(j))) & 
285               &     **bm_sapl_leaf(3))/sla(j))**bm_sapl_leaf(4)
286
287          IF ( pheno_type(j) .NE. 1 ) THEN
288             ! not evergreen
289             bm_sapl(j,icarbres,icarbon) = bm_sapl_carbres * bm_sapl(j,ileaf,icarbon)
290             bm_sapl(j,ilabile,icarbon) = bm_sapl_labile * bm_sapl(j,ileaf,icarbon)
291          ELSE
292             bm_sapl(j,icarbres,icarbon) = zero
293             bm_sapl(j,ilabile,icarbon) = zero
294          ENDIF ! (pheno_type_tab(j) .NE. 1 )
295
296          csa_sap = bm_sapl(j,ileaf,icarbon) / ( pipe_k1(j) / sla(j) )
297
298          dia = (mass_ratio_heart_sap(j) * csa_sap * dia_coeff(1) / pi ) ** dia_coeff(2)
299
300          bm_sapl(j,isapabove,icarbon) = &
301               bm_sapl_sapabove * pipe_density(j) * csa_sap * pipe_tune2(j) * dia ** pipe_tune3(j)
302          bm_sapl(j,isapbelow,icarbon) = bm_sapl(j,isapabove,icarbon)
303
304          bm_sapl(j,iheartabove,icarbon) = bm_sapl_heartabove * bm_sapl(j,isapabove,icarbon)
305          bm_sapl(j,iheartbelow,icarbon) = bm_sapl_heartbelow * bm_sapl(j,isapbelow,icarbon)
306
307       ELSE
308
309          !> 5.2 grasses
310
311          !!\latexonly
312          !!\input{stomate_data_sapl_grass.tex}
313          !!\endlatexonly
314
315          alpha = alpha_grass
316
317          IF ( natural(j) ) THEN
318             bm_sapl(j,ileaf,icarbon) = init_sapl_mass_leaf_nat / sla(j)
319          ELSE
320             bm_sapl(j,ileaf,icarbon) = init_sapl_mass_leaf_agri / sla(j)
321          ENDIF
322
323          bm_sapl(j,icarbres,icarbon) = init_sapl_mass_carbres *bm_sapl(j,ileaf,icarbon)
324          bm_sapl(j,ilabile,icarbon) = init_sapl_mass_labile *bm_sapl(j,ileaf,icarbon)
325
326          bm_sapl(j,isapabove,icarbon) = zero
327          bm_sapl(j,isapbelow,icarbon) = zero
328
329          bm_sapl(j,iheartabove,icarbon) = zero
330          bm_sapl(j,iheartbelow,icarbon) = zero
331
332       ENDIF !( is_tree(j) )
333
334       bm_sapl(j,iroot,icarbon) = init_sapl_mass_root * (1./alpha) * bm_sapl(j,ileaf,icarbon)
335
336       bm_sapl(j,ifruit,icarbon) = init_sapl_mass_fruit  * bm_sapl(j,ileaf,icarbon)
337
338
339       cn_leaf=cn_leaf_init(j)
340       cn_wood=cn_leaf/fcn_wood(j)
341       cn_root=cn_leaf/fcn_root(j)
342       
343       DO i=1,nparts
344          IF (i.EQ.1) THEN
345             bm_sapl(j,i,initrogen) = bm_sapl(j,i,icarbon) / cn_leaf
346          ELSE IF (i.LT.iroot) THEN
347             bm_sapl(j,i,initrogen) = bm_sapl(j,i,icarbon) / cn_wood
348          ELSE
349             bm_sapl(j,i,initrogen) = bm_sapl(j,i,icarbon) / cn_root
350          ENDIF
351       ENDDO
352
353       IF ( printlev_loc >= 2 ) THEN
354          WRITE(numout,*) '       sapling biomass (gC):'
355          WRITE(numout,*) '         leaves: (::bm_sapl(j,ileaf,icarbon))',bm_sapl(j,ileaf,icarbon)
356          WRITE(numout,*) '         sap above ground: (::bm_sapl(j,ispabove,icarbon)):',bm_sapl(j,isapabove,icarbon)
357          WRITE(numout,*) '         sap below ground: (::bm_sapl(j,isapbelow,icarbon))',bm_sapl(j,isapbelow,icarbon)
358          WRITE(numout,*) '         heartwood above ground: (::bm_sapl(j,iheartabove,icarbon))',bm_sapl(j,iheartabove,icarbon)
359          WRITE(numout,*) '         heartwood below ground: (::bm_sapl(j,iheartbelow,icarbon))',bm_sapl(j,iheartbelow,icarbon)
360          WRITE(numout,*) '         roots: (::bm_sapl(j,iroot,icarbon))',bm_sapl(j,iroot,icarbon)
361          WRITE(numout,*) '         fruits: (::bm_sapl(j,ifruit,icarbon))',bm_sapl(j,ifruit,icarbon)
362          WRITE(numout,*) '         carbohydrate reserve: (::bm_sapl(j,icarbres,icarbon))',bm_sapl(j,icarbres,icarbon)
363          WRITE(numout,*) '         labile reserve: (::bm_sapl(j,ilabile,icarbon))',bm_sapl(j,ilabile,icarbon)
364       ENDIF
365
366       !
367       ! 6 migration speed (m/year)
368       !
369
370       IF ( is_tree(j) ) THEN
371
372          migrate(j) = migrate_tree
373
374       ELSE
375
376          ! can be any value as grasses are, per *definition*, everywhere (big leaf).
377          migrate(j) = migrate_grass
378
379       ENDIF !( is_tree(j) )
380
381       IF ( printlev_loc >= 2 ) WRITE(numout,*) '       migration speed (m/year): (::migrate(j))', migrate(j)
382
383       !
384       ! 7 critical stem diameter: beyond this diameter, the crown area no longer
385       !     increases (m)
386       !
387
388       IF ( is_tree(j) ) THEN
389
390          !!\latexonly
391          !!\input{stomate_data_stem_diameter.tex}
392          !!\endlatexonly
393
394          maxdia(j) = ( ( pipe_tune4(j) / ((pipe_tune2(j)*pipe_tune3(j))/(maxdia_coeff(1)**pipe_tune3(j))) ) &
395               ** ( un / ( pipe_tune3(j) - un ) ) ) * maxdia_coeff(2)
396          cn_sapl(j) = cn_sapl_init !crown of individual tree, first year
397
398       ELSE
399
400          maxdia(j) = undef
401          cn_sapl(j)=1
402
403       ENDIF !( is_tree(j) )
404
405       IF ( printlev_loc >= 2 ) WRITE(numout,*) '       critical stem diameter (m): (::maxdia(j))', maxdia(j)
406
407       !
408       ! 8 Coldest tolerable temperature (K)
409       !
410
411       IF ( ABS( tmin_crit(j) - undef ) .GT. min_stomate ) THEN
412          tmin_crit(j) = tmin_crit(j) + ZeroCelsius
413       ELSE
414          tmin_crit(j) = undef
415       ENDIF
416
417       IF ( printlev_loc >= 2 ) &
418            WRITE(numout,*) '       coldest tolerable temperature (K): (::tmin_crit(j))', tmin_crit(j)
419
420       !
421       ! 9 Maximum temperature of the coldest month: need to be below this temperature
422       !      for a certain time to regrow leaves next spring *(vernalization)* (K)
423       !
424
425       IF ( ABS ( tcm_crit(j) - undef ) .GT. min_stomate ) THEN
426          tcm_crit(j) = tcm_crit(j) + ZeroCelsius
427       ELSE
428          tcm_crit(j) = undef
429       ENDIF
430
431       IF ( printlev_loc >= 2 ) &
432            WRITE(numout,*) '       vernalization temperature (K): (::tcm_crit(j))', tcm_crit(j)
433
434       !
435       ! 10 critical values for phenology
436       !
437
438       ! 10.1 model used
439
440       IF ( printlev_loc >= 2 ) &
441            WRITE(numout,*) '       phenology model used: (::pheno_model(j)) ',pheno_model(j)
442
443       ! 10.2 growing degree days. What kind of gdd is meant (i.e. threshold 0 or -5 deg C
444       !        or whatever), depends on how this is used in stomate_phenology.
445
446
447       IF ( ( printlev_loc >= 2 ) .AND. ( ALL(pheno_gdd_crit(j,:) .NE. undef) ) ) THEN
448          WRITE(numout,*) '         critical GDD is a function of long term T (C): (::gdd)'
449          WRITE(numout,*) '          ',pheno_gdd_crit(j,1), &
450               ' + T *',pheno_gdd_crit(j,2), &
451               ' + T^2 *',pheno_gdd_crit(j,3)
452       ENDIF
453
454       ! consistency check
455
456       IF ( ( ( pheno_model(j) .EQ. 'moigdd' ) .OR. &
457            ( pheno_model(j) .EQ. 'humgdd' )       ) .AND. &
458            ( ANY(pheno_gdd_crit(j,:) .EQ. undef) )                      ) THEN
459          CALL ipslerr_p(3,'stomate_data','problem with phenology parameters, critical GDD. (::pheno_model)','','')
460       ENDIF
461
462       ! 10.3 number of growing days
463
464       IF ( ( printlev_loc >= 2 ) .AND. ( ngd_crit(j) .NE. undef ) ) &
465            WRITE(numout,*) '         critical NGD: (::ngd_crit(j))', ngd_crit(j)
466
467       ! 10.4 critical temperature for ncd vs. gdd function in phenology (C)
468
469       IF ( ( printlev_loc >= 2 ) .AND. ( ncdgdd_temp(j) .NE. undef ) ) &
470            WRITE(numout,*) '         critical temperature for NCD vs. GDD (C): (::ncdgdd_temp(j))', &
471            ncdgdd_temp(j)
472
473       ! 10.5 humidity fractions (0-1, unitless)
474
475       IF ( ( printlev_loc >= 2 ) .AND. ( hum_frac(j) .NE. undef ) ) &
476            WRITE(numout,*) '         critical humidity fraction: (::hum_frac(j))', &
477            &  hum_frac(j)
478
479
480       ! 10.6 minimum time elapsed since moisture minimum (days)
481
482       IF ( ( printlev_loc >= 2 ) .AND. ( hum_min_time(j) .NE. undef ) ) &
483            WRITE(numout,*) '         time to wait after moisture min (d): (::hum_min_time(j))', &
484        &    hum_min_time(j)
485
486       !
487       ! 11 critical values for senescence
488       !
489
490       ! 11.1 type of senescence
491
492       IF ( printlev_loc >= 2 ) WRITE(numout,*) '       type of senescence: (::senescence_type(j))',&
493            senescence_type(j)
494
495       ! 11.2 critical temperature for senescence (C)
496
497       IF ( ( printlev_loc >= 2 ) .AND. ( ALL(senescence_temp(j,:) .NE. undef) ) ) THEN
498          WRITE(numout,*) '         critical temperature for senescence (C) is'
499          WRITE(numout,*) '          a function of long term T (C): (::senescence_temp)'
500          WRITE(numout,*) '          ',senescence_temp(j,1), &
501               ' + T *',senescence_temp(j,2), &
502               ' + T^2 *',senescence_temp(j,3)
503       ENDIF
504
505       ! consistency check
506
507       IF ( ( ( senescence_type(j) .EQ. 'cold' ) .OR. &
508            ( senescence_type(j) .EQ. 'mixed' )      ) .AND. &
509            ( ANY(senescence_temp(j,:) .EQ. undef ) )           ) THEN
510          CALL ipslerr_p(3,'stomate_data','Problem with senescence parameters, temperature. (::senescence_type)','','')
511       ENDIF
512
513       ! 11.3 critical relative moisture availability for senescence
514
515       IF ( ( printlev_loc >= 2 ) .AND. ( senescence_hum(j) .NE. undef ) ) THEN
516          WRITE(numout,*)  ' max. critical relative moisture availability for' 
517          WRITE(numout,*)  ' senescence: (::senescence_hum(j))',  &
518               & senescence_hum(j)
519       ENDIF
520
521       ! consistency check
522
523       IF ( ( ( senescence_type(j) .EQ. 'dry' ) .OR. &
524            ( senescence_type(j) .EQ. 'mixed' )     ) .AND. &
525            ( senescence_hum(j) .EQ. undef )                   ) THEN
526          CALL ipslerr_p(3,'stomate_data','Problem with senescence parameters, humidity.(::senescence_type)','','')
527       ENDIF
528
529       ! 14.3 relative moisture availability above which there is no moisture-related
530       !      senescence (0-1, unitless)
531
532       IF ( ( printlev_loc >= 2 ) .AND. ( nosenescence_hum(j) .NE. undef ) ) THEN
533          WRITE(numout,*) '         relative moisture availability above which there is' 
534          WRITE(numout,*) '             no moisture-related senescence: (::nosenescence_hum(j))', &
535               &  nosenescence_hum(j)
536       ENDIF
537
538       ! consistency check
539
540       IF ( ( ( senescence_type(j) .EQ. 'dry' ) .OR. &
541            ( senescence_type(j) .EQ. 'mixed' )     ) .AND. &
542            ( nosenescence_hum(j) .EQ. undef )                   ) THEN
543          CALL ipslerr_p(3,'stomate_data','Problem with senescence parameters, humidity. (::senescence_type)','','')
544       ENDIF
545
546       !
547       ! 12 sapwood -> heartwood conversion time (days)
548       !
549
550       IF ( printlev_loc >= 2 ) &
551            WRITE(numout,*) '       sapwood -> heartwood conversion time (d): (::tau_sap(j))', tau_sap(j)
552
553       !
554       ! 13 fruit lifetime (days)
555       !
556
557       IF ( printlev_loc >= 2 ) WRITE(numout,*) '       fruit lifetime (d): (::tau_fruit(j))', tau_fruit(j)
558
559       !
560       ! 14 length of leaf death (days)
561       !      For evergreen trees, this variable determines the lifetime of the leaves.
562       !      Note that it is different from the value given in leaflife_tab.
563       !
564
565       IF ( printlev_loc >= 2 ) &
566            WRITE(numout,*) '       length of leaf death (d): (::leaffall(j))', leaffall(j)
567
568       !
569       ! 15 maximum lifetime of leaves (days)
570       !
571
572       IF ( ( printlev_loc >= 2 ) .AND. ( leafagecrit(j) .NE. undef ) ) &
573            WRITE(numout,*) '       critical leaf age (d): (::leafagecrit(j))', leafagecrit(j)
574
575       !
576       ! 16 time constant for leaf age discretisation (days)
577       !
578
579       leaf_timecst(j) = leafagecrit(j) / REAL( nleafages,r_std )
580
581       IF ( printlev_loc >= 2 ) &
582            WRITE(numout,*) '       time constant for leaf age discretisation (d): (::leaf_timecst(j))', &
583            leaf_timecst(j)
584
585       !
586       ! 17 minimum lai, initial (m^2.m^{-2})
587       !
588
589       IF ( is_tree(j) ) THEN
590          lai_initmin(j) = lai_initmin_tree
591       ELSE
592          lai_initmin(j) = lai_initmin_grass
593       ENDIF !( is_tree(j) )
594
595       IF ( printlev_loc >= 2 ) &
596            WRITE(numout,*) '       initial LAI: (::lai_initmin(j))', lai_initmin(j)
597
598       !
599       ! 19 maximum LAI (m^2.m^{-2})
600       !
601
602       IF ( printlev_loc >= 2 ) &
603            WRITE(numout,*) '       critical LAI above which no leaf allocation: (::lai_max(j))', lai_max(j)
604
605       !
606       ! 20 fraction of primary leaf and root allocation put into reserve (0-1, unitless)
607       !
608
609       IF ( printlev_loc >= 2 ) &
610            WRITE(numout,*) '       reserve allocation factor: (::ecureuil(j))', ecureuil(j)
611
612
613       !
614       ! 23 natural ?
615       !
616
617       IF ( printlev_loc >= 2 ) &
618            WRITE(numout,*) '       Natural: (::natural(j))', natural(j)
619
620       !
621       ! 24 Vcmax et Vjmax (umol.m^{-2}.s^{-1})
622       !
623
624       IF ( printlev_loc >= 2 ) &
625            WRITE(numout,*) '       Maximum rate of carboxylation: (::Vcmax_25(j))', vcmax25(j)
626       !
627       ! 25 constants for photosynthesis temperatures
628       !
629
630       IF ( printlev_loc >= 2 ) THEN
631
632
633          !
634          ! 26 Properties
635          !
636
637          WRITE(numout,*) '       C4 photosynthesis: (::is_c4(j))', is_c4(j)
638          WRITE(numout,*) '       Depth constant for root profile (m): (::1./humcste(j))', 1./humcste(j)
639
640       ENDIF
641
642       !
643       ! 27 extinction coefficient of the Monsi and Saeki (1953) relationship
644       !
645       IF ( printlev_loc >= 2 ) THEN
646          WRITE(numout,*) '       extinction coefficient: (::ext_coeff(j))', ext_coeff(j)
647       ENDIF
648
649       !
650       ! 30 fraction of allocatable biomass which is lost as growth respiration (0-1, unitless)
651       !
652       IF ( printlev_loc >= 2 ) &
653            WRITE(numout,*) '       growth respiration fraction: (::frac_growthresp(j))', frac_growthresp(j)
654
655    ENDDO ! Loop over # PFTS
656
657    !
658    ! 29 time scales for phenology and other processes (in days)
659    !
660
661    tau_longterm_max = coeff_tau_longterm * one_year
662
663    IF ( printlev_loc >= 2 ) THEN
664
665       WRITE(numout,*) '   > time scale for ''monthly'' moisture availability (d): (::tau_hum_month)', &
666            tau_hum_month
667       WRITE(numout,*) '   > time scale for ''weekly'' moisture availability (d): (::tau_hum_week)', &
668           tau_hum_week
669       WRITE(numout,*) '   > time scale for ''monthly'' 2 meter temperature (d): (::tau_t2m_month)', &
670            tau_t2m_month
671       WRITE(numout,*) '   > time scale for ''weekly'' 2 meter temperature (d): (::tau_t2m_week)', &
672            tau_t2m_week
673       WRITE(numout,*) '   > time scale for ''weekly'' GPP (d): (::tau_gpp_week)', &
674            tau_gpp_week
675       WRITE(numout,*) '   > time scale for ''monthly'' soil temperature (d): (::tau_tsoil_month)', &
676            tau_tsoil_month
677       WRITE(numout,*) '   > time scale for ''monthly'' soil humidity (d): (::tau_soilhum_month)', &
678            tau_soilhum_month
679       WRITE(numout,*) '   > time scale for vigour calculations (y): (::tau_longterm_max / one_year)', &
680            tau_longterm_max / one_year
681
682    ENDIF
683
684    IF (printlev >= 4) WRITE(numout,*) 'Leaving data'
685
686  END SUBROUTINE data
687
688END MODULE stomate_data
Note: See TracBrowser for help on using the repository browser.