source: branches/publications/ORCHIDEE-GMv3.2/ORCHIDEE/src_sticslai/Stics_Lai_Developpement.f90 @ 6940

Last change on this file since 6940 was 6940, checked in by jinfeng.chang, 4 years ago

add missing files for ORCHIDEE-GMv3.2

File size: 41.6 KB
Line 
1!*-----------------------------------------------------------------------------------------------------------------------------------------------------------* c!
2!> - Stics book paragraphe 3.1.1, page 40-44
3!> - In STICS, leaf area growth is driven by phasic development, temperature and stresses. An empirical plant density-dependent function represents inter-plant
4!! competition. For indeterminate plants, trophic competition is taken into account through a trophic stress index, while for determinate plants a maximal
5!! expansion rate threshold is calculated to avoid unrealistic leaf expansion.
6!! The calculation of leaf growth rate (deltai in m2m-2 d-1) is broken down: a first calculation of the LAI growth rate (in m2 plant-1 degree-day-1) describes
7!! a logistic curve, related to the ILEV, IAMF and ILAX phenological stages. This value is then multiplied by the effective crop temperature, the plant density
8!! combined with a density factor, supposed to stand for inter-plant competition, that is characteristic for the variety, and the water and nitrogen stress indices.
9!!
10!! At the end of the juvenile stage (IAMF), it is equal to vlaimax when the inflection of the dynamics (point of maximal rate) also occurs.
11!! Between the stages ILEV, IAMF and ILAX, the model performs linear interpolation based on development units (upvt) which include all the environmental effects
12!! on phasic development. As the ILAX stage approaches, it is possible to introduce a gradual decline in growth rate using the udlaimax parameter
13!! - the ULAI value beyond which there is a decline in the leaf growth rate. If udlaimax=3 it has no effect and the leaf stops growing at once at ILAX.
14!!
15!! The thermal function relies on crop temperature and cardinal temperatures (tcmin and tcmax) which differ from the ones used for the phasic development.
16!! The extreme threshold tcxstop is the same as for development.
17!!
18!! The density function is active solely after a given LAI threshold occurs (laicomp parameter) if the plant density (densite in plant m-2 possibly decreased
19!! by early frost) is greater than the bdens threshold, below which plant leaf area is assumed independent of density.  Beyond this density value, leaf area
20!! per plant decreases exponentially.  The adens parameter represents the ability of a plant to withstand increasing densities.  It depends on the species
21!! and may depend on the variety.  For branching or tillering plants, adens represents the plant’s branching or tillering ability (e. g. wheat or pea).
22!! For single-stem plants, adens represents competition between plant leaves within a given stand (e.g. maize or sunflower).
23!!
24!! Water and nitrogen affect leaf growth as limiting factors, i.e. stress indices whose values vary between 0 and 1. Water (turfac) and nitrogen deficits (innlai)
25!! are assumed to interact, justifying the use of the more severe of the two stresses. Meanwhile at the whole plant level the water-logging stress index is assumed
26!! to act independently
27!
28!
29!
30!> -    Features of determinate crops
31!!   Failure to account for trophic aspects in the calculation of leaf growth may cause problems when the radiation intercepted by the crop is insufficient to
32!!   ensure leaf expansion (e.g. for crops under a tree canopy or crops growing in winter).  Consequently, from the IAMF stage, we have introduced a trophic effect
33!!  to calculate the definitive LAI growth rate in the form of a maximum threshold for leaf expansion (deltaimaxi in m2m-2d-1) using the notion of the maximum
34!!   leaf expansion allowed per unit of biomass accumulated in the plant (sbvmax in cm2 g-1) and the daily biomass accumulation (dltams in t.ha-1day-1 possibly
35!!   complemented by remobilized reserve remobilj). sbvmax is calculated using the slamax and tigefeuil parameters.
36!!
37!> -    Features of indeterminate crops
38!!   It has been possible to test the robustness of the above formalisation on a variety of crops, including crops where there is an overlap between the
39!!   vegetative phase and the reproductive phase (soybean and flax for example).  However, when trophic competition between leaves and fruits is a driving force
40!!   for the production and management of the crop (for example tomato, sugarbeet), this formalisation is unsuitable.  We therefore calculate the deltai variable
41!!   so as to take trophic monitoring into consideration in the case of crops described as ‘indeterminate’, by introducing a trophic stress index (splai).
42!!   As a consequence, the LAI can decrease markedly during the theoretical growth phase if the crop is experiencing severe stresses during the harvested
43!!   organs filling phase.
44! *-----------------------------------------------------------------------------------------------------------------------------------------------------------* c
45!TODO: modify the structure into the same manner as ORCHIDEE
46
47! XCW 23/05/2013 modified
48! subroutine laidev(sc,p,pg,itk,t,soil,c,sta)
49
50 subroutine laidev(n,                               &  ! IN
51                   in_cycle,                        &  ! INout
52                   f_crop_recycle,                  &  ! INOUT
53                   f_sen_lai,                        &  ! INout
54                   tair,                            &  ! IN 
55                   t2m_min_daily,                            &  ! IN 
56                   gdh_daily,                       &  ! IN
57                   phoi,                            &  ! IN
58                   onarretesomcourdrp,              &  ! IN                         
59                   stempdiag_cm_daily,              &  ! IN
60                   shumdiag_cm_day,                 &  ! IN
61                   nlevobs,                         &  ! IN
62                   namfobs,                         &  ! IN
63                   nfloobs,                         &  ! IN
64                   nlanobs,                         &  ! IN
65                   nlaxobs,                         &  ! IN
66                   nmatobs,                         &  ! IN
67                   nrecobs,                         &  ! IN
68                   nsenobs,                         &  ! IN
69                   ndrpobs,                         &  ! IN
70                   dltams,                          &  ! IN
71                   eop,                             &  ! IN
72                   masec,                           &  ! IN
73                   masecveg,                        &  ! IN
74                   nsendltams,                      &  ! INOUT
75                   nsendltai,                       &  ! INOUT
76                   nsenpfeuilverte,                 &  ! INOU&
77                   nsendurvie,                      &  ! INOU&
78                   nsenndurvie,                     &  ! INOU&
79                   densiteequiv,                    &  ! INOU&
80                   nplt,                            &  ! INOU&
81                   tursla,                          &  ! INOU&
82                   sla,                             &  ! INOU&
83                   pfeuilverte,                     &  ! INOU&
84                   bsenlai,                         &  ! INOU&
85                   zrac,                            &  ! INOU&
86                   nrec,                            &  ! INOU&
87                   nlan,                            &  ! INOU&
88                   tcult,                           &  ! INOU&
89                   udevair,                         &  ! INOU&
90                   udevcult,                        &  ! INOU&
91                   ndrp,                            &  ! INOU&
92                   rfvi,                            &  ! INOU&
93                   nlev,                            &  ! INOU&
94                   nger,                            &  ! INOU&
95                   etatvernal,                      &  ! INOU&
96                   caljvc,                          &  ! INOU&
97                   rfpi,                            &  ! INOU&
98                   upvt,                            &  ! INOU&
99                   utp,                             &  ! INOU&
100                   somcour,                         &  ! INOU&
101                   somcourdrp,                      &  ! INOU&
102                   somcourutp,                      &  ! INOU&
103                   tdevelop,                        &  ! INOU&
104                   somtemp,                         &  ! INOU&
105                   somcourfauche,                   &  ! INOU&
106                   stpltger,                        &  ! INOU&
107                   R_stamflax,                      &  ! INOU&
108                   R_stlaxsen,                      &  ! INOU&
109                   R_stsenlan,                      &  ! INOU&
110                   stlevflo,                        &  ! INOU&
111                   nflo,                            &  ! INOU&
112                   R_stlevdrp,                      &  ! INOU&
113                   R_stflodrp,                      &  ! INOU&
114                   R_stdrpmat,                      &  ! INOU&
115                   nmat,                            &  ! INOU&
116                   nlax,                            &  ! INOU&
117                   nrecbutoir,                      &  ! INOU&
118                   group,                           &  ! INOU&
119                   ndebdes,                         &  ! INOU&
120                   R_stdrpdes,                      &  ! INOU&
121                   densite,                         &  ! INOU&
122                   densitelev,                      &  ! INOU&
123                   coeflev,                         &  ! INOU&
124                   densiteger,                      &  ! INOU&
125                   somelong,                        &  ! INOU&
126                   somger,                          &  ! INOU&
127                   humectation,                     &  ! INOU&
128                   nbjhumec,                        &  ! INOU&
129                   somtemphumec,                    &  ! INOU&
130                   stpltlev,                        &  ! INOU&
131                   namf,                            &  ! INOU&
132                   stmatrec,                        &  ! INOU&
133                   tustress,                        &  ! INOU&
134                   lai,                             &  ! INOU&
135                   somfeuille,                      &  ! INOU&
136                   pdlai,                           &  ! INOU&
137                   nbfeuille,                       &  ! INOU&
138                   reajust,                         &  ! INOU&
139                   ulai,                            &  ! INOU&
140                   pdulai,                          &  ! INOU&
141                   efdensite,                       &  ! INOU&
142                   tempeff,                         &  ! INOU&
143                   nstopfeuille,                    &  ! INOU&
144                   deltai,                          &  ! INOU&
145                   vmax,                            &  ! INOU&
146                   nsen,                            &  ! INOU&
147                   laisen,                          &  ! INOU&
148                   dltaisenat,                      &  ! INOU&
149                   nsencour,                        &  ! INOU&
150                   dltamsen,                        &  ! INOU&
151                   dltaisen,                        &  ! INOU&
152                   fgellev,                         &  ! INOU&
153                   gelee,                           &  ! INOU&
154                   fstressgel,                      &  ! INOU&
155                   pdlaisen,                        &  ! INOU&
156                   R_stlevamf,                      &  ! INOU&
157                   dernier_n,                       &  ! INOU&
158                   durvieI,                         &  ! INOU&
159                   durvie,                          &  ! INOU&
160                   ndebsen,                         &  ! INOU&
161                   somsenreste,                     &  ! INOU&
162                   shumrel,                          &  ! INOU&
163                   humrel,                          &  ! INOU&
164                   swfac,                           &  ! INOU&
165                   turfac,                          &  ! INOU&
166                   senfac,                          &  ! INOUT
167                   mafeuiljaune,                    &  ! INOUT
168                   msneojaune,                      &
169                   gslen,                           &
170                   drylen,                           &
171                   vswc, &
172                   histgrowth, &
173                   hist_sencour, &
174                   hist_latest, &
175                   doyhistst, &
176                   nbox, boxulai, boxndays, boxlai, boxlairem,boxtdev, &
177                   boxbiom,boxbiomrem, boxdurage, boxsomsenbase)           ! INOUT
178
179
180  USE Stics 
181  USE Besoins_en_froid
182  USE netcdf
183  USE ioipsl
184
185  IMPLICIT NONE
186
187 
188  ! declaration of variables
189
190  ! 0.1 input
191 
192  integer, intent(IN)                           :: n 
193  logical, intent(inout)                           :: in_cycle
194  logical, intent(inout)                           :: f_crop_recycle
195 
196  logical, intent(inout)                           :: f_sen_lai
197  real,    intent(IN)                           :: tair !> / Mean air temperature of the day // degree C
198  real,    intent(IN)                           :: t2m_min_daily !> / Minimum air temperature of the day // degree C
199  real,    intent(IN)                           :: gdh_daily  !> // daily gdh calculated according to halfhourly temperature // transmitted from stomate.f90 gdh_daily
200  real,    intent(IN)                           :: phoi   !> // OUTPUT // Photoperiod // hours
201  logical, intent(IN)                           :: onarretesomcourdrp 
202  real,    intent(IN), dimension(3)             :: stempdiag_cm_daily !> / soil temperature at 1 cm resolution for the sowing depth and neighbour layers // Degree C
203  real,    intent(IN), dimension(3)             :: shumdiag_cm_day     !> /soil moisture at 1 cm resolution for the sowing depth and neighbour layers // unit m3 m-3 with values ranging 0-1
204  integer(i_std), intent(IN)                           :: nlevobs 
205  integer, intent(IN)                           :: namfobs 
206  integer, intent(IN)                           :: nfloobs 
207  integer, intent(IN)                           :: nlanobs 
208  integer, intent(IN)                           :: nlaxobs 
209  integer, intent(IN)                           :: nmatobs 
210  integer, intent(IN)                           :: nrecobs 
211  integer, intent(IN)                           :: nsenobs 
212  integer, intent(IN)                           :: ndrpobs 
213  real,    intent(IN)                           :: dltams     !> // growth rate of the plant // t ha-1 day-1 (difference of total biomass)
214  real,    intent(IN)                           :: eop        !> // potential evanpotranspiration
215  real,    intent(INOUT)                        :: masec      !> // above ground drymatter // t ha-1
216  real,    intent(IN)                           :: masecveg   !> // vegetative dry matter  // t ha-1
217  real,    intent(IN)                           :: vswc   !> daily humrel data
218  real,    intent(IN)                           :: humrel
219
220! 0.2 inout
221 
222  ! these variables are for laidev specifically
223  real, intent(INOUT)        :: nsendltams
224  real, intent(INOUT)        :: nsendltai
225  real, intent(INOUT)        :: nsenpfeuilverte
226  real, intent(INOUT)        :: nsendurvie
227  real, intent(INOUT)        :: nsenndurvie
228  real, intent(INOUT)        :: densiteequiv
229  integer, intent(INOUT)        :: nplt
230  real, intent(INOUT)        :: tursla
231  real, intent(INOUT)        :: sla
232  real, intent(INOUT)        :: pfeuilverte
233  real, intent(INOUT)        :: bsenlai
234 
235  ! variables are involved in DEVELOPMENT
236
237  real,    intent(INOUT)        :: zrac
238  integer, intent(INOUT)        :: nrec
239  integer, intent(INOUT)        :: nlan
240  real,    intent(INOUT)        :: tcult
241  real,    intent(INOUT)        :: udevair
242  real,    intent(INOUT)        :: udevcult
243  integer,    intent(INOUT)        :: ndrp
244  real,    intent(INOUT)        :: rfvi
245  integer,    intent(INOUT)        :: nlev
246  integer,    intent(INOUT)        :: nger
247  logical,    intent(INOUT)        :: etatvernal
248  real,    intent(INOUT)        :: caljvc
249  real,    intent(INOUT)        :: rfpi
250  real,    intent(INOUT)        :: upvt
251  real,    intent(INOUT)        :: utp
252  real,    intent(INOUT)        :: somcour
253  real,    intent(INOUT)        :: somcourdrp
254  real,    intent(INOUT)        :: somcourutp
255  real,    intent(INOUT)        :: tdevelop
256  real,    intent(INOUT)        :: somtemp
257  real,    intent(INOUT)        :: somcourfauche
258  real,    intent(INOUT)        :: stpltger
259  real,    intent(INOUT)        :: R_stamflax
260  real,    intent(INOUT)        :: R_stlaxsen
261  real,    intent(INOUT)        :: R_stsenlan
262  real,    intent(INOUT)        :: stlevflo
263  integer,    intent(INOUT)        :: nflo
264  real,    intent(INOUT)        :: R_stlevdrp
265  real,    intent(INOUT)        :: R_stflodrp
266  real,    intent(INOUT)        :: R_stdrpmat
267  integer,    intent(INOUT)        :: nmat
268  integer,    intent(INOUT)        :: nlax
269  integer,    intent(INOUT)        :: nrecbutoir
270  real,    intent(INOUT)        :: group
271  integer,    intent(INOUT)        :: ndebdes
272  real,    intent(INOUT)        :: R_stdrpdes
273  real,    intent(INOUT)        :: densite
274  real,    intent(INOUT)        :: densitelev
275  real,    intent(INOUT)        :: coeflev
276  real,    intent(INOUT)        :: densiteger
277  real,    intent(INOUT)        :: somelong
278  real,    intent(INOUT)        :: somger
279  logical,    intent(INOUT)        :: humectation
280  integer,    intent(INOUT)        :: nbjhumec
281  real,    intent(INOUT)        :: somtemphumec
282  real,    intent(INOUT)        :: stpltlev
283  integer,    intent(INOUT)        :: namf
284  real,    intent(INOUT)        :: stmatrec
285 
286  ! these variables are involved in Lai_calculation
287   
288  real,    intent(INOUT)        :: tustress
289  real,    intent(INOUT)        :: lai
290  real,    intent(INOUT)        :: somfeuille
291  real,    intent(INOUT)        :: pdlai
292  integer,    intent(INOUT)        :: nbfeuille
293  real,    intent(INOUT)        :: reajust
294  real,    intent(INOUT)        :: ulai
295  real,    intent(INOUT)        :: pdulai
296  real,    intent(INOUT)        :: efdensite
297  real,    intent(INOUT)        :: tempeff
298  integer,    intent(INOUT)        :: nstopfeuille
299  real,    intent(INOUT)        :: deltai
300  real,    intent(INOUT)        :: vmax
301  integer,    intent(INOUT)        :: nsen
302  real,    intent(INOUT)        :: laisen
303  real,    intent(INOUT)        :: dltaisenat
304
305  ! these variables are involved in the LAIsenescence
306
307  integer,    intent(INOUT)      :: nsencour
308  real,    intent(INOUT)      :: dltamsen
309  real,    intent(INOUT)      :: dltaisen
310  real,    intent(INOUT)      :: fgellev
311  logical,    intent(INOUT)      :: gelee
312  real,    intent(INOUT)      :: fstressgel
313  real,    intent(INOUT)      :: pdlaisen
314  real,    intent(INOUT)      :: R_stlevamf
315  integer,    intent(INOUT)      :: dernier_n
316  real,    intent(INOUT)      :: durvieI
317  real,    intent(INOUT)      :: durvie
318  integer,    intent(INOUT)      :: ndebsen
319  real,    intent(INOUT)      :: somsenreste
320
321  real, dimension(300,5), intent(INOUT) :: histgrowth
322  integer, intent(INOUT)                :: hist_sencour
323  integer, intent(INOUT)                :: hist_latest
324  integer, intent(INOUT)                :: doyhistst
325
326! compartment senescence module
327!boxulai, boxndays, boxlai, boxlairem,boxtdev, boxbiom, boxbiomrem, boxdurage,
328!boxsomsenbase
329       integer(i_std), intent(IN)                  :: nbox
330       real(r_std), dimension(nbox), intent(INOUT) :: boxulai
331       integer(i_std), dimension(nbox), intent(INOUT) :: boxndays
332       real(r_std), dimension(nbox), intent(INOUT) :: boxlai
333       real(r_std), dimension(nbox), intent(INOUT) :: boxlairem
334       real(r_std), dimension(nbox), intent(INOUT) :: boxtdev
335       real(r_std), dimension(nbox), intent(INOUT) :: boxbiom
336       real(r_std), dimension(nbox), intent(INOUT) :: boxbiomrem
337       real(r_std), dimension(nbox), intent(INOUT) :: boxdurage
338       real(r_std), dimension(nbox), intent(INOUT) :: boxsomsenbase
339
340
341  ! these variables are involved in STRESS calculation
342 
343  real,    intent(INOUT)      :: shumrel
344  real,    intent(INOUT)      :: swfac
345  real,    intent(INOUT)      :: turfac
346  real,    intent(INOUT)      :: senfac
347
348  ! these variables are involved in CARBON ALLOCATION PROCESSES
349 
350  real,    intent(INOUT)      :: mafeuiljaune      !  Dry matter of yellow leaves // t.ha-1
351  real,    intent(INOUT)      :: msneojaune        ! Newly-formed senescent dry matter  // t.ha-1
352  integer,    intent(INOUT)      :: gslen
353  integer,    intent(INOUT)      :: drylen
354   
355
356! 0.4 local Variables
357  !real                          :: tdev 
358  real                          :: upvtutil
359  real                          :: tempdev
360
361 !  do i = 1, P_nbplantes   !possible multi species for intercropping.
362 !     on pourrait utiliser un pointeur pour raccourcir les écritures.
363 !     Ex : p = plante(i). Et au lieu de p(i), on aurait juste p.
364 !   
365 !     comments: adjustment of LAI, but it is only effective for the option of forcing LAI.
366 !     there are two options for treating lai, 1. forced, 2. calculated.  ---xcwu       
367   
368 !     if (lge(P_codesimul,'feuille') .eqv. .TRUE.) then
369 !         call recalcullevee(n, nlev, nlevobs, nlaxobs, lai(1:n+1), &
370 !                            ! tauxcouv(1:n+1),
371 !                            & P_codelaitr &
372 !                            !,estDominante
373 !                            & )
374 !     endif
375
376
377    if (P_codeplante /= 'snu') then    ! 'snu' is bare soil.
378    !ens = AS   ! AS is a scalar, the value is 1.
379
380
381    ! 1. In senescence process, we need the dltams and dltai at date of  nsencour (specific date).
382    !    but we should know that the variable nsencour is calculated in the latter process in our laidev module.
383    !    So here , we can only use the data at the date of nsencour + 1
384     
385         if (P_culturean /= 1 .AND. nsencour /= 0 .AND. nsencour > 365) then
386   
387            if (n == nsencour - 365 + 1) then
388               nsendltams = dltams
389               nsendltai  = deltai
390               nsenpfeuilverte = pfeuilverte 
391               nsendurvie = durvie
392            endif
393     
394 
395            ! We need the durvie of the following day after nsencour
396            if (n == nsencour- 365 + 2) then
397               nsenndurvie = durvie
398            endif
399         endif
400         
401         if (P_culturean == 1 .AND. nsencour /= 0) then
402            if (n == nsencour + 1) then
403               nsendltams = dltams
404               nsendltai  = deltai
405               nsenpfeuilverte = pfeuilverte 
406               nsendurvie = durvie
407            endif
408 
409            ! We need the durvie of the following day after nsencour
410            if (n == nsencour + 2) then
411               nsenndurvie = durvie
412            endif
413         endif
414
415    ! 2. DEVELOPMENT OF PLANTS   
416    !>>>> DEVELOPMENT PROCESSES<<<<<<<<<<<<<!
417    !>>>> MODIFIED 05/17/13<<<<<<<<<<<<<<<<<!
418         
419         ! SET THE NPLT
420         ! P_iwater //julian day of the beginning of the simulation // jour julien
421
422         nplt = P_iplt0 - P_iwater + 1  ! although nplt is changed each day, but it is actually a constance if the sowing date is determined.
423
424
425         call develop2(n, in_cycle, nplt, tair,  gdh_daily, turfac,  phoi, onarretesomcourdrp,  stempdiag_cm_daily, shumdiag_cm_day,lai,           &  !> INPUTS
426                       nlevobs, namfobs, nfloobs, nlanobs, nlaxobs, nmatobs, nrecobs, nsenobs, ndrpobs, nrecbutoir,                  &  !> INPUTS
427                       masec, namf,  ndrp, nflo, nger, nlan, nlev, nrec, etatvernal, caljvc,                                         &  !> INOUT
428                       upvt, utp, somcour, somcourdrp, somcourfauche, somcourutp, somtemp, zrac,                                                      &  !> INOUT
429                       coeflev, somelong, somger, humectation, nbjhumec, somtemphumec, densite, densitelev, nlax, nmat, nsen, stlevflo, ndebdes,      &  !> INOUT
430                       R_stlevamf, R_stamflax, R_stsenlan, R_stlaxsen, R_stlevdrp, R_stflodrp, R_stdrpmat, R_stdrpdes, densiteger,                    &  !> INOUT
431                       udevair, udevcult, rfvi, rfpi, tdevelop, stpltger, upvtutil, stmatrec, group, tcult, stpltlev,                                 &
432                       f_crop_recycle, gslen, drylen)                                   !> OUT
433
434   
435       ! DR 12/09/2012 on pense que il y avait un pb sur le calcul de deltalai qui pour les CAS etait calculé à partir de la densite equivalente
436       ! cette densite equivalente ne doit etre utilisée que pour le calcul de efdensite (en cas de cas et pour la plante 2 ...)
437       ! si les 2 plantes sont arrivées à levee on recalcule la densite equivalente
438
439       ! comments: at this moment, we do not consider the intercropping system---xcwu
440       !
441       ! if(sc%P_nbplantes.gt.1 .and. i==2)then
442       !       if((p(1)%nlev>0 .and. p(2)%nlev>0) .and. (p(1)%nlev==sc%n.or.p(2)%nlev==sc%n))then
443              !TODO: commenté ou pas ça ?
444       !            if (p(1)%nplt > p(2)%nplt) then
445       !             call EnvoyerMsgHistorique('culture principale semée après culture associée')
446       !           endif
447       !    call F_densite_equiv(i,p(1)%densite,p(2)%densite,p(i)%estDominante,p(1)%P_bdens ,p(2)%P_bdens,&
448       !     p(1)%densiteequiv,p(2)%densiteequiv)
449       ! le jour de la levee la desnite de la plante 2 est celle qu'on vient de recalculer (densite equivalente)
450       !     p(2)%densitelev=p(2)%densite
451       ! endif
452       ! endif
453       ! si on a une seule plante on reaffecte la densite equivalente au cas ou il y aurait eu une modif sur la densite courante (gel ou autre )
454       
455 
456       ! We consider that the crop system is homogenous and is without intercropping system, So the densiteequiv is equal to densite
457   
458        densiteequiv = densite
459
460
461       ! on décompose les valeurs AOAS calculées dans develop en AO et AS
462       !  p(i)%rdtint(sc%AO,p(i)%nbrecolte-1) = p(i)%rdtint(sc%AOAS,p(i)%nbrecolte-1) * p(i)%surf(sc%AO)
463       !  p(i)%rdtint(sc%AS,p(i)%nbrecolte-1) = p(i)%rdtint(sc%AOAS,p(i)%nbrecolte-1) * p(i)%surf(sc%AS)
464
465       !  p(i)%nfruit(sc%AO,p(i)%P_nboite) = p(i)%nfruit(sc%AOAS,p(i)%P_nboite) * p(i)%surf(sc%AO)
466       !  p(i)%nfruit(sc%AS,p(i)%P_nboite) = p(i)%nfruit(sc%AOAS,p(i)%P_nboite) * p(i)%surf(sc%AS)
467
468       !  p(i)%pdsfruit(sc%AO,p(i)%P_nboite) = p(i)%pdsfruit(sc%AOAS,p(i)%P_nboite) * p(i)%surf(sc%AO)
469       !  p(i)%pdsfruit(sc%AS,p(i)%P_nboite) = p(i)%pdsfruit(sc%AOAS,p(i)%P_nboite) * p(i)%surf(sc%AS)
470
471       !#if DEBUG == 1
472       !        if (iand(sc%develop,4) >0) call develop_debug_read_output(1232,sc,p,itk,i)
473       !        if (iand(sc%develop,8) >0) call develop_debug_write_output(1233,sc,p,itk,i)
474       !        if (iand(sc%develop,16) >0) call develop_debug_test_output(1234,sc,p,itk,i)
475       !#endif
476
477       ! dr - 17/11/05: pb de recolte
478       ! dr et ml - 06/12/05 : c'etait idiot pour les perennes de forcer la recolte des 2 cultures en meme temps!!
479     
480
481       ! this part for the harvest of the cropping system,
482       ! at this moment, we do not consider the cropping system, so we comment them
483       !    if (itk(i)%P_coderecolteassoc == 1) then
484       !      if (p(i)%nrec /= 0 .and. ALL(p(:)%P_codeperenne == 1))then
485       !        p(:)%nrec = p(i)%nrec
486       !      endif
487       !    endif
488   
489    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
490    ! 3. Record the lai at nsen
491    !    After all of the processes relating to the calculation of daily LAI, we got the daily LAI.
492    !    Specially, we need record the lai value on the date when the sen stage begin. 
493       
494          if (nsen > 0 .and. nlan == 0 .and. f_sen_lai) then
495                         
496             bsenlai = lai     ! bsenlai is the lai at the date when the sen starts
497           
498             f_sen_lai = .FALSE. 
499          endif
500
501    ! 4. CALCULATION OF LAI   
502    !>>>> MODIFIED 05/17/13<<<<<<<<<<<<<<<<<!
503
504        if (lge(P_codesimul,'feuille') .eqv. .FALSE.) then    ! If the LAI is not forced.
505          if (P_codelaitr == 1) then  ! if the calculation of LAI is true LAI instead of soil cover
506
507           
508            !do ens = sc%AS,sc%AO
509            !  if (surf > 0.) then  ! here we do not consider whether or not there are crop!
510            !    sc%ens = ens       ! we also do not need the intercropping system
511
512            ! Calcul the  LAI by calling calai subroutine
513
514            !call calai_(P_codeinnact,P_codeh2oact,P_codehypo,P_codegermin,P_codcueille,P_codlainet,codeulaivernal,      & !IN
515            !            P_codeindetermin,P_codeinitprec,P_codeperenne,P_codetemp,turfac,innlai,                         &
516            !            P_laiplantule,P_phyllotherme,udevair,udevcult,P_dlaimaxbrut,P_udlaimax,n,numcult,nbjmax,        &
517            !            nplt,nlev,namf,nlax,ndrp,nrec,upvt,upvtutil,P_vlaimax,P_laicomp,somcour,somcourutp,             &
518            !            P_stlevamf,P_stamflax,densite,P_adens,P_bdens,P_tcxstop,tcult,P_tcmin,P_tcmax,P_pentlaimax,     &
519            !            P_dlaimin,exolai,sourcepuits,P_splaimin,P_splaimax,P_slamin,P_slamax,P_tigefeuil,P_tustressmin, &
520            !            fstressgel,dltamsen,sla,remobilj,dltams,                                                        & !IN
521            !            tustress,efdensite,nstopfeuille,deltai,splai,fpv,P_stlaxsen,tempeff,                            & !OUT
522            !            lai,somfeuille,nbfeuille,nsen,nlan,P_dlaimax,reajust,ulai,vmax,dltaisenat,laisen,               & !INOUT
523            !            dltaisen,P_stsenlan,densiteequiv)
524
525            call calai_(turfac,                                                                     & ! IN
526                        udevair,udevcult,n,                                                         & ! IN
527                        nplt,nlev,namf,nlax,ndrp,nrec,upvt,upvtutil,somcour,somcourutp,             & ! IN
528                        densite,tcult,                                                              & ! IN
529                        fstressgel,dltamsen,sla,                                                    & ! IN
530                        dltams, densiteequiv, bsenlai,                                              & ! IN
531                        pdlaisen, lai, pdlai, pdulai, somfeuille, nbfeuille, nsen, nlan,            & ! INOUT
532                        reajust,  ulai, vmax,                                                       & ! INOUT
533                        dltaisenat, laisen, dltaisen, R_stsenlan,                                   & ! INOUT
534                        tustress, efdensite, nstopfeuille, deltai, R_stlaxsen, tempeff,             &
535                        histgrowth, hist_sencour, hist_latest, doyhistst)!,&
536!                        nbox, boxulai, boxndays, boxlai, boxlairem,boxtdev, &
537!                        boxbiom,boxbiomrem, boxdurage, boxsomsenbase)               ! OUT
538
539       
540            ! Rognage & Effeuillage, seulement si lai > 0
541            ! ATTENTION! on travaille sur masec(n-1) car à ce moment de la boucle journalière masec(n) n'a pas encore été calculé...
542            !     if (p(i)%lai(ens,sc%n) > 0.0) then
543            !       if (itk(i)%P_codrognage == 2) then
544            !         call rognage(itk(i)%P_codcalrogne,p(i)%hauteur(ens),p(i)%largeur,   & ! IN
545            !                      itk(i)%P_hautrogne,itk(i)%P_margerogne,itk(i)%P_largrogne, &
546            !                      p(i)%dfol,p(i)%P_hautbase,itk(i)%P_interrang,            &
547            !                      p(i)%sla(ens),p(i)%P_tigefeuil,itk(i)%P_biorognem,sc%n,  &
548            !           p(i)%nrogne,p(i)%CNplante(ens),         & DR 20/07/2012 on a plus besoin du Cnplante il est calculé apres
549            !                      p(i)%nrogne,                                             &
550            !                      p(i)%lairogne(ens),p(i)%biorogne(ens),                   & ! INOUT & OUT
551            !                      p(i)%lai(ens,sc%n),p(i)%masec(ens,sc%n-1),               &
552            !                      p(i)%varrapforme,p(i)%P_forme,p(i)%biorognecum(ens),     &
553            !                      p(i)%lairognecum(ens))
554
555            !       ! Cumuls AO/AS
556            !         !p(i)%lairogne(ens)
557            !         !p(i)%biorogne(ens)
558            !         p(i)%lai(sc%aoas,sc%n) = p(i)%lai(sc%as,sc%n) * p(i)%surf(sc%as)              &
559            !                                + p(i)%lai(sc%ao,sc%n) * p(i)%surf(sc%ao)
560            !         p(i)%masec(sc%aoas,sc%n-1) = p(i)%masec(sc%as,sc%n-1) * p(i)%surf(sc%as)      &
561            !                                    + p(i)%masec(sc%ao,sc%n-1) * p(i)%surf(sc%ao)
562            !         !p(i)%varrapforme
563            !         !p(i)%P_forme
564            !         !p(i)%biorognecum(ens)
565            !         !p(i)%lairognecum(ens)
566
567            !         ! NB le 23/05 recyclage de la biomasse rognée
568            !         ! PB - 08/08/2009 - Pas besoin d'apporter des résidus de rognage tous les jours.
569            !         !                   On conditionne l'apport en résidus, si il y a de la masse rognée.
570
571            !         if (sc%n == p(i)%nrogne .and. p(i)%biorogne(ens) > 0.) then
572            !           soil%itrav1 = 1
573            !           soil%itrav2 = 1
574            !           sc%ires = 2
575            !           Cfeupc  =  42.
576            !           CNrogne = Cfeupc / p(i)%CNplante(ens)
577            !         
578            !           !Bruno 06/2012 - ajout quantites de C et N tombees au sol apres rognage
579            !           Crogne=p(i)%biorogne(ens)*Cfeupc*10.
580            !           Nrogne=Crogne/CNrogne
581            !           p(i)%QCrogne = p(i)%QCrogne + Crogne
582            !           p(i)%QNrogne = p(i)%QNrogne + Nrogne
583
584
585            !          !EC 0!6/08/2012 Ajout du code ires des feuilles rognées
586            !          sc%ires = 2
587            !          call ResidusApportSurfaceSol(p(i)%biorogne(ens),Cfeupc,-1.e-10,sc%ires,pg%P_CNresmin(sc%ires),  & ! IN
588            !                 pg%P_CNresmax(sc%ires),1.,1.,pg%P_Qmulchdec(sc%ires),sc%nbCouchesSol,pg%nbResidus,   &
589            !                 itk(i)%nap,sc%airg(sc%n+1),CNrogne,sc%Cnondec(1:10),sc%Nnondec(1:10),                &
590            !                 sc%Cres(1:sc%nbCouchesSol,1:pg%nbResidus),sc%Nres(1:sc%nbCouchesSol,1:pg%nbResidus), &
591            !                 sc%QCapp,sc%QNapp,sc%QCresorg,sc%QNresorg)                                                 ! INOUT
592
593            !         Wr=0.
594            !         if (CNrogne > 0.) Wr=1./CNrogne
595            !         call ResiduParamDec(pg%P_awb(sc%ires),pg%P_bwb(sc%ires),pg%P_cwb(sc%ires),pg%P_CroCo(sc%ires),  &
596            !            pg%P_akres(sc%ires),pg%P_bkres(sc%ires),pg%P_ahres(sc%ires),pg%P_bhres(sc%ires),             &
597            !            sc%Wb(sc%ires),sc%kres(sc%ires),sc%hres(sc%ires),Wr)
598            !         endif
599
600            !       endif
601
602
603            !      Comments----xcwu
604            !      The removal of leaves seems only effective for indeterminate crops, so at this moment we do not consider this processes,
605            !      because we mainly focus on the determinate crops.
606            !      So, this process is commented.
607
608            !       if (P_codeffeuil /= 1) then  ! if there is no thinning management! it seems to be only effective for indeterminate crops, such as grape etc
609            !         call effeuill(P_codcaleffeuil,P_laidebeff,deltai(ens,sc%n), &
610            !                       P_effeuil,sla(ens),n,neffeuil,  &
611            !                       P_laieffeuil, P_codetransrad,P_hautmax,P_khaut,     &
612            !                       dfol,P_largrogne,P_codhauteff,largeur,      &
613            !                       lai(ens,sc%n),masec(ens,sc%n-1),P_hautbase,        & ! INOUT
614            !                       varrapforme,bioeffcum(ens),laieffcum(ens))
615
616            !       ! Cumuls AO/AS
617            !         p(i)%lai(sc%aoas,sc%n) = p(i)%lai(sc%as,sc%n) * p(i)%surf(sc%as)              &
618            !                                  + p(i)%lai(sc%ao,sc%n) * p(i)%surf(sc%ao)
619            !         p(i)%masec(sc%aoas,sc%n-1) = p(i)%masec(sc%as,sc%n-1) * p(i)%surf(sc%as)      &
620            !                                      + p(i)%masec(sc%ao,sc%n-1) * p(i)%surf(sc%ao)
621
622            !       endif
623                endif
624              endif
625      !      end do
626      !     else
627      !      call TauxRecouvrement(sc,pg,p(i),itk(i))
628      !    endif
629      !  endif
630      !endif
631
632
633    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
634    ! 4. Senescence processes for LAI
635    !    This process is belong to GROWTH module in STICS
636    !    In this subroutine, we mainly calculate the senescenced lai and biomass
637    !    In detail, see senescence.f90 and the red book
638   
639         call senescen(nlev,n, lai, senfac,  &
640                       nbfeuille,t2m_min_daily, &
641                       densitelev,          &
642                       masecveg,         &
643                       nstopfeuille,somcour,nrec,      &
644                       ulai, &
645                       tdevelop,    &
646                       somtemp,nsenpfeuilverte,nsendltai,nsendurvie, nsendltams, nsenndurvie,          &    ! INPUTS
647                       pdlaisen, dernier_n,nsencour,dltaisen,dltamsen,fstressgel,fgellev,gelee,      & !INOUT
648                       densite,laisen,nlan,R_stsenlan,nsen,R_stlaxsen,namf,nlax,R_stlevamf,      &
649                       R_stamflax,durvie,     &
650                       ndebsen,somsenreste, mafeuiljaune, msneojaune,        &
651                       durvieI,deltai,dltams, &
652                       histgrowth, hist_sencour, hist_latest, doyhistst,  &
653                       nbox, boxulai, boxndays, boxlai, boxlairem,boxtdev, &
654                       boxbiom,boxbiomrem, boxdurage, boxsomsenbase)               ! OUT
655!)  ! OUTPUT
656   
657    ! >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
658    ! 5. Calculation of sla ---xcw
659    !    In stics, the calculation of sla is in reptair.f90, The process reptair is belong to GROWTH in STICS
660    !    the reptair process is called after the senescence
661   
662         ! calculation of tursla, the slowing effects of water stress on sla
663         ! turstress is calculated in Stics_Calcul_Lai.f90
664   
665         tursla = (tursla+tustress)/2.
666         IF (P_codesla==1) THEN
667             sla = min(tursla * P_slamax, P_slamin)
668             if ((lge(P_codesimul,'feuille') .eqv. .FALSE.).and.(sla > P_slamax  )) then
669                sla = P_slamax
670             endif
671             if (sla < P_slamin) sla = P_slamin
672         ELSEIF (P_codesla==2) THEN
673            if ((n <= nlev+1) .or. (nlev==0)) then ! before or start leafing
674               sla = P_slamax
675            else if ( ndrp == 0 ) then
676               tempdev = min(1. + (3. * somcourdrp / (P_stlevamf + P_stamflax)), 4.)
677               if (tursla>1) then
678                   tursla = 1
679               else if (tursla<0) then
680                   tursla = 0
681               endif
682               tursla = (1-tursla)/2 + 0.5 !tursla ~ [0.5 1]
683               sla = (P_slamax-P_slamin) * (exp(-1.5 * (tempdev - 1.))) * tursla + P_slamin
684            else if (ndrp >0 .and. nmat==0) then
685               sla = sla ! after ndrp, there is no imposed change of sla
686            else
687               ! do nothing, keep the same sla until harvest
688            endif
689             
690         ELSE
691             WRITE(numout,*) 'codesla not recognized: ',P_codesla
692             stop
693         ENDIF
694
695    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
696    ! 6. Calculation of pfeuilverte
697    !    In stics, the calculation of pfeuilverte is in reptair.f90
698    !    the reptair process is called after the senescence
699    !    masec == aboveground dry matter (t ha-1)
700    !    dltams == growth rate of plant (t ha-1 day-1)   
701   
702     
703         if (masec <= 0.0) then
704            pfeuilverte = 0.0
705         else
706            if (dltams > 0.) then
707               pfeuilverte = (deltai/P_slamax*1e2) / dltams
708               pfeuilverte = min(1.0,pfeuilverte)
709            else
710               pfeuilverte = 0.0
711            endif
712         endif
713
714    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
715    ! 8.  Calculation of the water stress on LAI, turfac
716          call stress(n,                  &  ! IN
717                      nrec,               & 
718                      lai,                & 
719                      eop,                &
720                      shumdiag_cm_day,    &  ! IN
721                      vswc,               &  ! IN                   
722                      humrel,             &  ! IN
723                      shumrel,             &  ! OUT
724                      swfac,              &  ! INOUT)   ! calculate the turfac stress (water stress)
725                      turfac,             &
726                      senfac)                 
727
728    ! This part is for debugiing
729
730   ! write(*,*) 'wu: DOY is',n
731   ! write(*,*) 'wu: nplt, nger, nlev, ndrp, nlax, nsen is ', nplt, nger, nlev, ndrp, nlax, nsen
732   ! write(*,*) 'wu: upvt is ', upvt
733   ! write(*,*) 'wu: somcour is ', somcour
734   
735   ! write(*,*) 'wu: udevair is ', udevair
736   ! write(*,*) 'wu: rfvi and rfpi is ', rfvi, rfpi
737   ! write(*,*) 'wu: efdensity and tempeff and nstopfeille is ', efdensite, tempeff,nstopfeuille
738   ! write(*,*) 'wu: dltai and dltaisen is ', deltai, dltaisen
739   ! write(*,*) 'wu: ulai and tustress  is ', ulai, tustress
740
741
742
743    endif
744return
745end
Note: See TracBrowser for help on using the repository browser.