source: branches/publications/ORCHIDEE-GMv3.2/ORCHIDEE/src_sticslai/senescen.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! *******************************************************************
2!    version 6.0
3!    reste pb cultures fauchées à voir avec Domi
4!--------------------------------------------------
5!    calcul de la senescence de la matiere seche
6!    avec eventuellement une partie residuelle pour les cultures fauchees
7!
8!    calcul de la sénescence du LAI pour l'option LAI brut
9! derniere modif 30/05/05
10! *******************************************************************
11!> In STICS shoot senescence only concerns leaves: dry matter and LAI.  For cut crops, it also affects residual biomass after cutting.
12!> - Stics book paragraphe 3.1.2, page 44-46
13!> While in the first versions of the model senescence was implicit (Brisson et al., 1998a), it is now explicit, with a clear distinction between natural
14!! senescence due to the natural ageing of leaves, and senescence accelerated by stresses (water, nitrogen, frost). The concept of leaf lifespan,
15!! used for example by Maas (1993), is applied to the green leaf area and biomass produced. The leaf area and part of the leaf biomass produced on a given day
16!! is therefore lost through senescence once the lifetime has elapsed (Duru et al., 1995). This part corresponds to the ratiosen parameter, taking into account
17!! the part which was remobilised during its senescence.
18!!
19!! Calculation of lifespan:
20!!
21!! The natural lifespan of leaves (durage) is defined by two values: the  lifespan of early leaves, or durvieI (expressed as a proportion of durvieF ) and the
22!! lifespan of the last leaves emitted, or durvieF (assumed genotype-dependent). Until the IAMF stage, the natural lifespan, calculated for the day when the
23!! leaves are emitted (I0) is durvieI; from IAMF to ILAX, the natural lifespan increases between durvieI and durvieF as a function of the leaf development variable ULAI.
24!! Because of water or nitrogen stress, the current lifespan may be shortened if the stress undergone is more intense than previous stresses.
25!! Specific stress indices for senescence are introduced (senfac and innsenes).  Frost (fstressgel that can be either fgeljuv or fgelveg) may also reduce or
26!! even cancel  lifespan.  In the event of over-fertilisation with nitrogen (inn >1), the foliage lifespan is increased from the IAMF stage up to a maximum
27!! given by the durviesupmax parameter.
28!! The lifespan of leaves is not expressed in degree.days (like phasic development), because this has the disadvantage of stopping any progression as soon as
29!! the temperature is lower than the base temperature (tdmin).  To remedy this problem, the senescence course (somsen) is expressed by cumulated Q10 units
30!! (with Q10=2), i.e. an exponential type function.
31!! The senescence course between I0 and I is affected by the same cardinal temperatures as phasic development and can be slown down by stresses. The lifespan parameter
32!! of the leaf (durvieF) expressed in Q10 units represents about 20% of the same  lifespan expressed in degree.days.
33!!
34!! Calculation of senescence:
35!!
36!! Material produced on day I0 disappears by senescence after a period corresponding to durvie(I0). Depending on the evolution of temperature and of lifespan
37!! as a function of phenology and stresses, senescence can vary from one day to another and affect several days of production (J=I0, I0+1, …) or, on the contrary,
38!! none (durvieE(I0)>somsen(I)).  This principle is applied to the biomass (dltamsen) and leaf area index (ltaisen). In general, the leaf biomass produced
39!! does not completely disappear (remobilisation):  the ratiosen (<1) parameter enables the definition of the senescent proportion with reference to production.
40!! It is the pfeuilverte ratio which defines the proportion of leaves in the biomass produced.
41!! The cumulative senescent foliage is laisen. If the crop is a forage crop benefiting from residual dry matter from the previous cycle (msresiduel parameter),
42!! the senescence of residual dry matter (deltamsresen) starts as from cutting.  It occurs at a rate estimated from the relative daily lifespan course and
43!! weighted by the remobilisation (ratiosen).
44! *-----------------------------------------------------------------------------------------------------------------------------------------------------------* c
45!subroutine senescen(nlev,n,nbjmax,lai,P_codeinnact,P_codeh2oact,senfac,innsenes,P_codlainet,  &
46!                    P_codeperenne,nbfeuille,P_nbfgellev,P_codgellev,tcultmin,P_tletale,P_tdebgel, &
47!                    P_tgellev90,P_tgellev10,ipl,densitelev,densiteassoc,P_codgeljuv,          &
48!                    P_tgeljuv90,P_tgeljuv10,P_codgelveg,P_tgelveg90,P_tgelveg10,masecveg,         &
49!                    nstopfeuille,somcour,resperenne,ndrp,nrec,P_QNpltminINN,numcult,      &
50!                    P_codeinitprec,ulai,P_vlaimax,durvieI,P_durvieF,inn,P_durviesupmax,         &
51!                    P_codestrphot,phoi,P_phobasesen,dltams,P_msresiduel,P_ratiosen,tdevelop,    &
52!                    somtemp,pfeuilverte,deltai,P_lai0,                                    &
53!                    dernier_n,nsencour,dltaisen,dltamsen,fstressgel,fgellev,gelee,      &
54!                    densite,laisen,nlan,R_stsenlan,nsen,R_stlaxsen,namf,nlax,R_stlevamf,      &
55!                    R_stamflax,nrecbutoir,mortplante,nst2,mortplanteN,durvie,strphot,     &
56!                    msres,dltamsres,ndebsen,somsenreste,msresjaune,mafeuiljaune,        &
57!                    msneojaune,dltamstombe,QNplante,P_dltamsminsen,P_dltamsmaxsen,          &
58!                    P_alphaphot,strphotveille)
59
60subroutine senescen(nlev,n, lai, senfac,  &
61                    nbfeuille,t2m_min_daily, &
62                    densitelev,          &
63                    masecveg,         &
64                    nstopfeuille,somcour,nrec,      &
65                    ulai, &
66                    tdevelop,    &
67                    somtemp,nsenpfeuilverte,nsendltai,nsendurvie, nsendltams, nsenndurvie,         &    ! INPUTS
68                    pdlaisen, dernier_n,nsencour,dltaisen,dltamsen,fstressgel,fgellev,gelee,      & !INOUT
69                    densite,laisen,nlan,R_stsenlan,nsen,R_stlaxsen,namf,nlax,R_stlevamf,      &
70                    R_stamflax,durvie,     &
71                    ndebsen,somsenreste, mafeuiljaune, msneojaune,  & !INOUT
72                    durvieI, deltai, dltams, &
73                    histgrowth, hist_sencour, hist_latest, doyhistst,  &
74                    nbox, boxulai, boxndays, boxlai, boxlairem,boxtdev, &
75                    boxbiom,boxbiomrem, boxdurage, boxsomsenbase)  ! OUTPUT
76
77USE constantes
78USE Stics
79USE Divers_gel
80!USE Messages
81
82  implicit none
83
84! DECLARATION
85
86! 0.1 INPUTS
87
88  integer, intent(IN)    :: nlev 
89  integer, intent(IN)    :: n 
90  !integer, intent(IN)    :: nbjmax 
91  real,    intent(IN)    :: lai   !> // OUTPUT // Leaf area index (table) // m2 leafs  m-2 soil
92  !integer, intent(IN)    :: P_codeinnact  !> // PARAMETER // code activating  nitrogen stress effect on the crop: yes (1), no (2) // code 1/2 // PARAM // 0
93  !integer, intent(IN)    :: P_codeh2oact  !> // PARAMETER // code to activate  water stress effect on the crop: yes (1), no (2) // code 1/2 // PARAM // 0
94  real,    intent(IN)    :: senfac   !> // OUTPUT // Water stress index on senescence // 0-1
95  !real,    intent(IN)    :: innsenes   !> // OUTPUT // Index of nitrogen stress active on leaf death // P_innmin to 1
96  !integer, intent(IN)    :: P_codlainet  !> // PARAMETER //option of calculation of the LAI (1 : direct LAInet; 2 : LAInet = gross LAI - senescent LAI) // code 1/2 // PARPLT // 0
97  !integer, intent(IN)    :: P_codeperenne  !> // PARAMETER // option defining the annual (1) or perenial (2) character of the plant // code 1/2 // PARPLT // 0
98  integer, intent(IN)    :: nbfeuille   !> // OUTPUT // Number of leaves on main stem // SD
99  !integer, intent(IN)    :: P_nbfgellev  !> // PARAMETER // leaf number at the end of the juvenile phase (frost sensitivity)  // nb pl-1 // PARPLT // 1
100  !integer, intent(IN)    :: P_codgellev  !> // PARAMETER // activation of plantlet frost // code 1/2 // PARPLT // 0
101  real,    intent(IN)      :: t2m_min_daily    !// minimum daily temperature // celsius degree
102  !real,    intent(IN)    :: tcultmin 
103  !real,    intent(IN)    :: P_tletale  !> // PARAMETER // lethal temperature for the plant // degree C // PARPLT // 1
104  !real,    intent(IN)    :: P_tdebgel  !> // PARAMETER // temperature of frost beginning // degree C // PARPLT // 1
105  !real,    intent(IN)    :: P_tgellev90  !> // PARAMETER // temperature corresponding to 90% of frost damage on the plantlet  // degree C // PARPLT // 1
106  !real,    intent(IN)    :: P_tgellev10  !> // PARAMETER // temperature corresponding to 10% of frost damage on the plantlet  // degree C // PARPLT // 1
107  !integer, intent(IN)    :: ipl         !> number of kinds of  crops
108  real,    intent(IN)    :: densitelev 
109  !real,    intent(IN)    :: densiteassoc 
110  !integer, intent(IN)    :: P_codgeljuv  !> // PARAMETER // activation of LAI frost at the juvenile stadge // code 1/2 // PARPLT // 0
111  !real,    intent(IN)    :: P_tgeljuv90  !> // PARAMETER // temperature corresponding to 90 % of frost damage on the LAI (juvenile stage) // degree C // PARPLT // 1
112  !real,    intent(IN)    :: P_tgeljuv10  !> // PARAMETER // temperature corresponding to 10 % of frost damage on the LAI (juvenile stage) // degree C // PARPLT // 1
113  !integer, intent(IN)    :: P_codgelveg  !> // PARAMETER // activation of LAI frost at adult stage // code 1/2 // PARPLT // 0
114  !real,    intent(IN)    :: P_tgelveg90  !> // PARAMETER // temperature corresponding to 90 % of frost damage on the LAI (adult stage) // degree C // PARPLT // 1
115  !real,    intent(IN)    :: P_tgelveg10  !> // PARAMETER // temperature corresponding to 10 % of frost damage on the LAI (adult stage) // degree C // PARPLT // 1
116  real,    intent(IN)    :: masecveg   !> // OUTPUT // Vegetative dry matter // t.ha-1
117  integer, intent(IN)    :: nstopfeuille 
118  real,    intent(IN)    :: somcour   !> // OUTPUT // Cumulated units of development between two stages // degree.days
119  !real,    intent(IN)    :: resperenne   !> // OUTPUT // C crop reserve, during the cropping season, or during the intercrop period (for perenial crops) // t ha-1
120  !integer, intent(IN)    :: ndrp 
121  integer, intent(IN)    :: nrec 
122  !real,    intent(IN)    :: P_QNpltminINN  !> // PARAMETER // minimal amount of nitrogen in the plant allowing INN computing // kg ha-1 // PARAM // 1
123  !integer, intent(IN)    :: numcult 
124  !integer, intent(IN)    :: P_codeinitprec  !> // PARAMETER // reinitializing initial status in case of chaining simulations : yes (1), no (2) // code 1/2 // PARAM // 0
125  real,    intent(IN)    :: ulai  !>   // OUTPUT // Daily relative development unit for LAI // 0-3
126  !real,    intent(IN)    :: P_vlaimax  !> // PARAMETER // ULAI at inflection point of the function DELTAI=f(ULAI) // SD // PARPLT // 1
127  !real,    intent(IN)    :: P_durvieF  !> // PARAMETER // maximal  lifespan of an adult leaf expressed in summation of P_Q10=2 (2**(T-Tbase)) // P_Q10 // PARPLT // 1
128  !real,    intent(IN)    :: inn   !> // OUTPUT // Nitrogen nutrition index (satisfaction of nitrogen needs ) // 0-1
129  !real,    intent(IN)    :: P_durviesupmax  !> // PARAMETER // proportion of additional lifespan due to an overfertilization // SD // PARPLT // 1
130  !integer, intent(IN)    :: P_codestrphot  !> // PARAMETER // activation of the photoperiodic stress on lifespan : yes (1), no (2) // code 1/2 // PARPLT // 0
131  !real,    intent(IN)    :: phoi   !> // OUTPUT // Photoperiod // hours
132  !real,    intent(IN)    :: P_phobasesen  !> // PARAMETER // photoperiod under which the photoperiodic stress is activated on the leaf lifespan // heures // PARPLT // 1
133  !real,    intent(IN)    :: dltams(nbjmax)  !>   // OUTPUT // Growth rate of the plant  // t ha-1.j-1
134  !real,    intent(IN)    :: P_msresiduel  !> // PARAMETER // Residual dry matter after a cut // t ha-1 // PARTEC // 1
135  !real,    intent(IN)    :: P_ratiosen  !> // PARAMETER // fraction of senescent biomass (by ratio at the total biomass) // between 0 and 1 // PARPLT // 1
136  real,    intent(IN)    :: tdevelop
137  real,    intent(IN)    :: somtemp   !> // OUTPUT // Sum of temperatures // degree C.j
138  !real,    intent(IN)    :: pfeuilverte(nbjmax) !>  // OUTPUT // Proportion of green leaves in total non-senescent biomass // 0-1
139  !real,    intent(IN)    :: deltai(nbjmax)  !>   // OUTPUT // Daily increase of the green leaf index // m2 leafs.m-2 soil
140  !real,    intent(IN)    :: P_lai0  !> // PARAMETER // Initial leaf area index // m2 m-2 // INIT // 1
141  !real,    intent(IN)    :: dltamstombe 
142  !real,    intent(IN)    :: QNplante   !> // OUTPUT // Amount of nitrogen taken up by the plant  // kgN.ha-1
143  !real,    intent(IN)    :: P_dltamsminsen  !> // PARAMETER // threshold value of deltams from which the photoperiodic effect on senescence is maximal // t ha-1j-1 // PARAM // 1
144  !real,    intent(IN)    :: P_dltamsmaxsen  !> // PARAMETER // threshold value of deltams from which there is no more photoperiodic effect on senescence // t ha-1j-1 // PARPLT // 1
145  !real,    intent(IN)    :: P_alphaphot  !> // PARAMETER // parameter of photoperiodic effect on leaf lifespan // P_Q10 // PARAM // 1
146  real,     intent(IN)    :: nsenpfeuilverte
147  real,     intent(IN)    :: nsendltai 
148  real,     intent(IN)    :: nsendurvie 
149  real,     intent(IN)    :: nsendltams 
150  real,     intent(IN)    :: nsenndurvie 
151  real,     intent(IN)    :: deltai
152  real,     intent(IN)    :: dltams
153
154
155! 0. 2 INOUT
156
157  real,     intent(INOUT) :: pdlaisen    !>  senescence leaf area index of previous day // m2 leafs  m-2 soil
158  integer, intent(INOUT) :: dernier_n 
159  integer, intent(INOUT) :: nsencour 
160  real,    intent(INOUT) :: dltaisen   !> // OUTPUT // Daily increase of the senescent leaf index // m2.m-2 sol.j-1
161  real,    intent(INOUT) :: dltamsen   !> // OUTPUT // Senescence rate // t ha-1 j-1
162  real,    intent(INOUT) :: fstressgel   !> // OUTPUT // Frost index on the LAI // 0-1
163  real,    intent(INOUT) :: fgellev 
164  logical, intent(INOUT) :: gelee 
165  real,    intent(INOUT) :: densite   !> // OUTPUT // Actual sowing density // plants.m-2
166  real,    intent(INOUT) :: laisen    ! // OUTPUT // Leaf area index of senescent leaves // m2 leafs  m-2 soil
167  integer, intent(INOUT) :: nlan 
168  real,    intent(INOUT) :: R_stsenlan  !> // PARAMETER // Sum of development units between the stages SEN et LAN // degree.days // PARPLT // 1
169  integer, intent(INOUT) :: nsen 
170  real,    intent(INOUT) :: R_stlaxsen  !> // PARAMETER // Sum of development units between the stages LAX and SEN // degree.days // PARPLT // 1
171  integer, intent(INOUT) :: namf 
172  integer, intent(INOUT) :: nlax 
173  real,    intent(INOUT) :: R_stlevamf  !> // PARAMETER // Sum of development units between the stages LEV and AMF // degree.days // PARPLT // 1
174  real,    intent(INOUT) :: R_stamflax  !> // PARAMETER // Sum of development units between the stages AMF and LAX // degree.days // PARPLT // 1
175  !integer, intent(INOUT) :: nrecbutoir 
176  !integer, intent(INOUT) :: mortplante 
177  !integer, intent(INOUT) :: nst2 
178  !integer, intent(INOUT) :: mortplanteN 
179  real,    intent(INOUT) :: durvie !<   // OUTPUT // Actual life span of the leaf surface //  degree C
180  !real,    intent(INOUT) :: msres 
181  !real,    intent(INOUT) :: dltamsres 
182  integer, intent(INOUT) :: ndebsen 
183  real,    intent(INOUT) :: somsenreste 
184  !real,    intent(INOUT) :: msresjaune   !> // OUTPUT // Senescent residual dry matter  // t.ha-1
185  real,    intent(INOUT) :: mafeuiljaune   !> // OUTPUT // Dry matter of yellow leaves // t.ha-1
186  real,    intent(INOUT) :: msneojaune   !> // OUTPUT // Newly-formed senescent dry matter  // t.ha-1
187  !real,    intent(INOUT) :: strphotveille 
188
189  real, dimension(300,5), intent(INOUT) :: histgrowth
190  integer, intent(INOUT)                :: hist_sencour
191  integer, intent(INOUT)                :: hist_latest
192  integer, intent(INOUT)                :: doyhistst
193
194! compartment senescence module
195!boxulai, boxndays, boxlai, boxlairem,boxtdev, boxbiom, boxbiomrem, boxdurage,
196!boxsomsenbase
197       integer(i_std), intent(IN)                  :: nbox
198       real(r_std), dimension(nbox), intent(INOUT) :: boxulai
199       integer(i_std), dimension(nbox), intent(INOUT) :: boxndays
200       real(r_std), dimension(nbox), intent(INOUT) :: boxlai
201       real(r_std), dimension(nbox), intent(INOUT) :: boxlairem
202       real(r_std), dimension(nbox), intent(INOUT) :: boxtdev
203       real(r_std), dimension(nbox), intent(INOUT) :: boxbiom
204       real(r_std), dimension(nbox), intent(INOUT) :: boxbiomrem
205       real(r_std), dimension(nbox), intent(INOUT) :: boxdurage
206       real(r_std), dimension(nbox), intent(INOUT) :: boxsomsenbase
207
208
209! 0.3 OUTPUT
210
211  !real,    intent(OUT)   :: strphot 
212  real,    intent(OUT)    :: durvieI 
213
214! 0.4 LOCAL VARIABLES
215
216  integer :: i  !> 
217  integer :: ii
218  integer :: ncompteur  !> 
219  integer :: debut  !> 
220  integer :: fin 
221  !real    :: vitsenres  !> 
222  real    :: durage  !> 
223  real    :: durviesup  !> 
224  real    :: senstress  !> 
225  real    :: somsen  !> 
226  !real    :: tgel10  !>
227  !real    :: tgel90
228  !real    :: msresTMP 
229  real   :: dtdevrat
230  real   :: tempdays, templaisen, tempamsen
231  logical :: myprint
232
233
234
235! 1. whether or not go into this subroutine
236      !:On shunte ce programme si la plante n'est pas levée ou si le lai est nul
237      ! If the emergence of crop did not occur, the value of nsencour is always 0
238      myprint = .false.
239      if (nlev == 0) then
240        ! *- 07/09/06 DR et NB et IGC (la sainte famille)
241        ! *- apres un 2 jours d'errements nadine a été notre lumiere !!
242        ! *- Qu'elle en soit louée !
243        ! *- initialisation de nsencour au debut du cycle
244
245
246        nsencour = n  !initialization of nsencour
247
248        return
249      endif
250        if (P_codlainet==3) then ! accumulated part for compartment senescence
251            if (ulai < 1 .or. ulai > 3) then
252                write(*,*) 'xuhui: ulai: ', ulai, ', likely a bug'
253            elseif (ulai==3) then
254!                write(*,*) 'xuhui: ulai==3 with deltai ',deltai
255            else
256                ii = nbox
257                do while (ulai < boxulai(ii))
258                    ii = ii - 1
259                enddo
260                if (ii < 1) then
261                    write(*,*) 'xuhui: bug: ii<1'
262                    ii = 1
263                endif
264!                write(*,*) 'xuhui: ulai ',ulai,' ii ', ii
265!                write(*,*) 'xuhui: boxulai ',boxulai
266                if (boxndays(ii) == 0) then ! the 1st leaf in this pool
267                    boxsomsenbase(ii) = somtemp
268                    boxlai(ii) = boxlai(ii) + P_lai0/real(nbox) ! equally distributed the original lai to all boxes
269!                    if (ii == 1) then  ! put P_lai0 in the 1st box
270!                        boxlai(ii) = boxlai(ii) + P_lai0
271!                    endif
272                endif
273                boxndays(ii) = boxndays(ii) + 1
274                boxlai(ii) = boxlai(ii) + deltai
275                boxlairem(ii) = boxlairem(ii) + deltai
276                boxbiom(ii) = boxbiom(ii) + dltams
277                boxbiomrem(ii) = boxbiomrem(ii) + dltams
278                boxtdev(ii) = boxtdev(ii) + tdevelop
279            endif
280        endif
281 
282
283      !: PB - 08/03/2005 - si encore de la matière jaune on continue la senescence
284     
285      ! in our model, we do not consider the processes regarding the yellow leaves, so the variable mafeuiljaune is not necessary
286
287      !if (lai <= 0 .and. mafeuiljaune <= 0.) then
288      if (lai <= 0) then
289        dltaisen = 0.
290        dltamsen = 0.
291        return
292      endif
293
294
295! 2. calculation of the lifetime for leaves
296 
297! ** calcul de la durée de vie
298! *---------------------------
299! *-  a) calcul du LAI net : on utilise version 4 avec une durée de vie fixe (stlevsenms)
300! *-  b) calcul du LAI brut: la durée de vie est fonction de l'age de la plante
301! *-                         (durée de vie potentielle x facteurs de stress)
302! *-   durée de vie potentielle des feuilles (durage) calculée dans calai
303! *-   facteurs de stress eau et azote calculés à partir de senfac (transpi.for)
304! *-   et innsenes (absorn.for)
305! *-------------------------------------------------------------------
306! *- spécifiques de la sénescence
307
308
309   ! 2. 1 Comments: calculation of senstress,
310   ! if we chose the direct calculation of LAI, the senstress is 1.
311   ! calculation of senfac is implemented in stress process
312
313
314      if (P_codeinnact == 1 .and. P_codeh2oact == 1) then
315        senstress = min(senfac,P_innsenes)
316      endif
317
318      if (P_codeinnact == 2 .and. P_codeh2oact == 1) then
319        senstress = min(senfac,1.)
320      endif
321
322      if (P_codeinnact == 1 .and. P_codeh2oact == 2) then
323        senstress = min(1.,P_innsenes)
324      endif
325
326      if (P_codeinnact == 2 .and. P_codeh2oact == 2) senstress = 1.
327
328      !if (P_codlainet == 1) senstress = 1.
329
330
331
332   ! 2. 2 ** calcul du STRESS lié au GEL  (utilisation de la fonction GEL.for)
333   ! *-------------------------------------------------------------------
334   ! *- GEL entre levée et plantule diminue la densité de population uniquement pour les annuelles
335
336      if (P_codeperenne == 1)  then
337        if (nbfeuille <= P_nbfgellev) then
338          !fgellev = GEL(P_codgellev,tcultmin,P_tletale,P_tdebgel,P_tgellev90,P_tgellev10)
339          fgellev = GEL(P_codgellev, t2m_min_daily, P_tgellev90, P_tgellev10)
340          if (fgellev < 1.) gelee = .TRUE.    ! we must consider the effects of gel
341          !if (ipl == 1) then     ! only consider one kind of plant
342            densite = min(densite, densitelev*fgellev)
343          !else
344          !  densite = min(densite, (densitelev + densiteassoc)*fgellev)
345          !endif
346        endif
347      endif
348
349      ! ** GEL de l'appareil végétatif à deux stades: juvénil (avant AMF) et adulte (après AMF)
350      if (namf <= 0) then
351        fstressgel = GEL(P_codgeljuv, t2m_min_daily, P_tgeljuv90, P_tgeljuv10)
352      else
353        fstressgel = GEL(P_codgelveg, t2m_min_daily, P_tgelveg90, P_tgelveg10)
354      endif
355      if (fstressgel < 1.) gelee = .TRUE.
356
357 
358 
359      ! ** GEL LETAL pour la plante
360      if (fstressgel <= 0.) then  ! frost damage occurs
361        dltaisen  = lai
362        dltamsen  = masecveg   ! daily dry matter
363        !laisen(1) = laisen(0) + dltaisen ! consider the former day
364        laisen = pdlaisen + dltaisen
365
366
367        ! ** affectation des stades
368        if (nstopfeuille > 0 .and. nlan == 0) then
369          nlan = n
370          R_stsenlan = somcour
371          if (nsen == 0) then
372            nsen = n
373            R_stlaxsen = somcour
374            R_stsenlan = 0.
375          endif
376        else
377          ! comments: resperenne is only for perennial crops
378          ! From the reptair process, we know that the resperenne is 0 for annual determinant crops
379         
380          if (resperenne <= 0.) then  ! resperenne is global constance
381            if (namf == 0) then
382              namf = n
383              nlax = n
384              nsen = n
385              nlan = n
386              R_stlevamf = somcour
387              R_stamflax = 0.
388              R_stlaxsen = 0.
389              R_stsenlan = 0.
390            endif
391            if (nlax == 0 .and. namf > 0) then
392              nlax = n
393              nsen = n
394              nlan = n
395              R_stamflax = somcour
396              R_stlaxsen = 0.
397              R_stsenlan = 0.
398            endif
399            if (nsen == 0 .and. nlax > 0) then
400              nsen = n
401              nlan = n
402              R_stlaxsen = somcour
403              R_stsenlan = 0.
404            endif
405          endif
406        endif
407        !if (resperenne <= 0.) then
408        !  call EnvoyerMsgHistorique(1167)
409        !  !: modif - NB - 20/09 - mort le jour d'apres pour écriture rapport
410        !  !nrecbutoir = n+1
411        !  !mortplante = 1
412        !  !: pour non plantade dans bilan.for - TODO : à mettre ici ou dans bilan ?
413        !  if (ndrp == 0) nst2 = 1
414        !endif
415        goto 30
416      endif
417
418      if (n == nlev .and. P_codeperenne == 1) then
419        nsencour = nlev
420      endif
421
422      ! *- DR et NB - 25/08/08 - s'il n'y a plus d'azote dans la plante on la fait mourir
423      ! DR 11/04/2012 pour le prairie on a des pbs car quand on coupe on a plus d'azote dans la  plante
424      ! donc je mets < au lieu de <=
425     
426      ! Comments: at this moment, the crop are not N limited, so QNplante >= p_qnpltmininn
427
428      !if ((namf > 0 .and. nrec == 0) .and. QNplante < P_QNpltminINN) then
429      !  call EnvoyerMsgHistorique(2099)
430      !! modif NB 20/09 mort le jour d'apres pour écriture rapport
431      !  nrecbutoir  = n+1
432      !  mortplanteN = 1
433      !! pour non plantade dans bilan.for
434      !  if (ndrp == 0) nst2 = 1
435      !endif
436
437! *- NB - le 22/12/2003 - pb senescence perenne vu avec FR
438! *- PB - 05/01/2005 - mise en commentaire jusqu'à nouvel ordre.
439! *- Sinon cela perturbe le calcul de la sénescence lors de l'enchainement des pérennes.
440      !if (n == nlev .and. P_codeperenne == 2 .and. numcult == 1) then
441      !  nsencour = 1
442      !endif
443
444! ** calcul de la durée de vie
445! *----------------------------
446! *-  nsenscour = jour où la biomasse sénescente a été produite
447! *-  durage = durée de vie potentielle fonction de la P_phénologie
448! *-  (varie entre durvieI et P_durvieF; durvieI avant lax)
449
450      debut = nsencour
451      if (P_codeperenne == 2 .and. P_codeinitprec == 2 .and. numcult > 1 .and. nsencour > n) then
452        fin = dernier_n
453      else
454        fin = n
455        dernier_n = n
456      endif
457
458! calculation of the durvieI
459      !durvieI = P_ratiodurvieI * P_durvieF(P_variete)
460      durvieI = P_ratiodurvieI * P_durvieF   !(P_variete)
461      !do i = debut, fin
462     
463      if (P_codlainet==1) then
464            if (ulai <= P_vlaimax) then
465              durage = durvieI   ! where is the durvieI??   p(ipl)%durvieI = p(ipl)%P_ratiodurvieI * p(ipl)%P_durvieF(itk(ipl)%P_variete), it seems a parameter
466            else
467              durage = durvieI + ((ulai - P_vlaimax) / (3. - P_vlaimax)) * (P_durvieF - durvieI)
468            endif
469   
470            !! ** application des stress qui raccourcissent la durée de vie
471            !if (i < fin) then
472            !  durvie(i) = min(durvie(i), durage * min(senstress, fstressgel))
473            !else
474              durvie = durage * min(senstress, fstressgel)
475            !endif
476   
477            ! ** augmentation de durée de vie maxi en cas de surfertilisation azotée
478            !if (namf > 0 .and. i > namf .and. P_codlainet == 2) then
479            !  if (inn > 1.) then
480            !    durviesup = P_durvieF * min(P_durviesupmax, (inn - 1.))
481    !   !               write(*,*)n,'durviesup',P_durvieF,durviesup
482            !  else
483            !    durviesup = 0.
484            !  endif
485            !else
486            !  durviesup = 0.
487            !endif
488   
489            durviesup = 0.
490       
491            durvie = durvie + durviesup
492   
493            ! domi 04/01/06 pour inaki insertion du stress photoperiodique
494            ! NB le 22/06/06 introduction dans la boucle des "i"
495   
496            ! at this moment, we do not consider the stress of photoperiod on senescence---xcw
497            !if (P_codestrphot == 1 .and. nlax /= 0) then
498            !  if (phoi < P_phobasesen) then
499            !    call stressphot(dltams(n), P_dltamsminsen, P_dltamsmaxsen, P_alphaphot, strphot, strphotveille)
500            !    durvie(i) = durvie(i) * strphot
501            !  endif
502            !endif
503   
504          !end do
505   
506   
507          !if (P_codeperenne == 2 .and. P_codeinitprec == 2 .and. numcult > 1 .and. debut > n) then
508          !  debut = 1
509          !  fin = n
510          !  goto 5
511          !endif
512   
513          ! ** la base de temps de la sénescence est calculée dans develop.for
514          ! *- tdevelop(n) est un P_Q10 sommé dans la variable somtemp
515   
516          ! ** pour les cultures fauchées : senescence de la matiere seche résiduelle
517          ! *- qui demarre dès la coupe
518          ! *- PB - 25/01/2004 - ajout de msresTMP et dltamsres pour calculer le delta msres(n-1)<->msres(n)
519          ! *- et ce, pour pallier à un problème dans le calcul de msresjaune pour les pérennes enchainées fauchées
520          !if (msres > 0.) then
521          !  vitsenres = P_msresiduel / durvieI
522          !  msresTMP = msres - (P_ratiosen * vitsenres * tdevelop(n))
523          !  msresTMP = max(0.,msresTMP)
524          !  dltamsres = max(0., msres - msresTMP)
525          !  msres = msresTMP
526          !else
527          !  dltamsres = 0.
528          !endif
529   
530   
531   
532          ! ** pour les autres cultures ou pour la matière sèche néoformée
533   
534   
535          ! ** DOMI VOIR POUR PRAIRIE UTILISATION DE somcourfauche
536          ! ** 2/ senescence de la matiere seche neoformee qui demarre au stade
537          ! *- debut de senescence de la matiere seche identifie par le parcours
538          ! *- de developpement depuis la levee : stlevsenms
539          ! *- on calcule vitsen le jour du début de sénescence : ndebsen
540   
541          !if (somtemp >= durvie(n)) then
542          if (somtemp >= durvie) then
543            if (ndebsen == 0) then ! premier jour de senescence
544              ndebsen = n
545              nsencour = nlev
546   
547              !dltamsen =  dltams(nsencour) * P_ratiosen * pfeuilverte(nsencour)
548              dltamsen = nsendltams * P_ratiosen * nsenpfeuilverte     
549   
550              if (P_codlainet == 2) then
551                dltaisen = nsendltai + P_lai0
552              endif
553   
554            ! ** enlever le reste du LAI si deltai(nsencour) = 0.
555   
556            else  ! jours suivants
557              somsen = somsenreste
558              dltamsen = 0.
559              ncompteur = 0
560              if (P_codlainet == 2) dltaisen = 0.
561   
562    ! ** calcul de la somme de température depuis le dernier jour de sénescence
563    ! *- PB - 05/01/2005 - ajout d'un test pour l'enchainement des pérennes. On calcul somsen
564    ! *- à partir du nsencour du cycle précédent jusqu'à la fin du cycle précédent et
565    ! *- on l'ajoute à la somme des températures depuis le début du cycle courant jusqu'au jour n
566              !if (P_codeperenne == 2 .and. P_codeinitprec == 2 .and. numcult > 1 .and. nsencour > n) then
567              !  do i = nsencour+1, dernier_n
568              !    somsen = somsen + tdevelop(i)
569              !  end do
570              !  do i = 1, n
571              !    somsen = somsen + tdevelop(i)
572              !  end do
573              !else
574                !do i = nsencour+1, n
575                !  somsen = somsen + tdevelop(i)
576                !end do
577                if (nsencour == nlev .and. nsencour /= 0) then
578                   somsen = somsen + tdevelop
579                endif
580              !endif
581   
582   
583              !if (somsen > durvie(nsencour+1)) then
584              if (somsen > nsenndurvie) then
58520              nsencour = nsencour + 1
586   
587   
588    ! *- PB - 05/01/2004 - ajout d'un test lors de l'enchainement des pérennes pour que nsencour
589    ! *- ne dépasse pas dernier_n. Si nsencour>dernier_n alors on a atteint la fin du cycle précédent,
590    ! *- et donc on doit repartir du début du cycle courant. Soit nsencour = 1
591               ! if (P_codeperenne == 2 .and. P_codeinitprec == 2 .and. numcult > 1 .and. nsencour > dernier_n) then
592               !   nsencour = 1
593               ! endif
594   
595                !somsen = somsen - durvie(nsencour)
596                somsen = somsen - nsendurvie
597                if (somsen < 0.) then
598                  somsen = 0.
599                  nsencour = nsencour - 1
600                  goto 40 ! TODO : a priori goto 40 pourrait être remplacé par goto 30s
601                endif
602   
603                !dltamsen = dltamsen + (dltams(nsencour) * P_ratiosen * pfeuilverte(nsencour))
604                dltamsen = dltamsen + (nsendltams * P_ratiosen * nsenpfeuilverte)
605                ncompteur = ncompteur+1
606                if (P_codlainet == 2) then
607                  !dltaisen = dltaisen + deltai(nsencour)
608                  dltaisen = dltaisen + nsendltai
609                endif
610   
611                ! ** si la durée de vie est raccourcie, on peut "faire mourir" plusieurs jours de production
612                ! *- PB - 10/03/2005 - Petite modificatoin, on compare nsencour à dernier_n et non plus à n.
613     
614       
615                !if (somsen > durvie(nsencour) .and. nsencour < dernier_n) goto 20
616                if (somsen >  nsendurvie .and. nsencour < dernier_n) goto 20
617   
618                somsenreste = somsen
619              else
620                dltamsen = 0.
621                if (P_codlainet == 2) dltaisen = 0.
622              endif
62340            continue
624            endif
625          else
626            dltamsen = 0.
627            if (P_codlainet == 2) dltaisen = 0.
628          endif
629   
630    elseif (P_codlainet==2) then
631!        dltaisen = 0.
632!        write(*,*) 'new senescence process, xuhui'
633        IF (ulai <= P_vlaimax) THEN
634            durage = durvieI
635        ELSE
636            durage = durvieI + ((ulai - P_vlaimax) / (3. - P_vlaimax)) * (P_durvieF - durvieI)
637        ENDIF
638        durvie = durage * min(senstress, fstressgel)
639        ! after calculation of durvie, record it
640        histgrowth(hist_latest,5) = durvie
641        ! update all previous durvie with minimum values
642        if (hist_latest > hist_sencour) then
643            do ii = hist_sencour, hist_latest-1
644                histgrowth(ii,5) = min(histgrowth(ii,5),durvie)
645            end do
646        endif
647            !!!! this part need to be done when we have Nitrogen stress
648    !           DO ii = hist_sencour, hist_latest
649    !               if (namf > 0 .and. i > namf .and. P_codlainet == 2) then
650    !                 if (inn > 1.) then
651    !                   durviesup = P_durvieF * min(P_durviesupmax, (inn - 1.))
652    !       !           write(*,*)n,'durviesup',P_durvieF,durviesup
653    !                 else
654    !                   durviesup = 0.
655    !                 endif
656    !               else
657    !                 durviesup = 0.
658    !               endif
659    !               durvie(i) = durvie(i) + durviesup
660    !           ENDDO
661       
662        !IF (somtemp >= durvie) THEN ! histgrowth(hist_seencour+1,5)
663        IF (somtemp >= histgrowth(hist_sencour+1,5)) THEN ! histgrowth(hist_seencour+1,5)
664        ! considering durvie(1..n-1) is smaller or equal to durvie(n),
665        ! this is a conservative choice to start senescence a period of full
666        ! senescense
667        ! a radical one should be set as >= durvie(nsencour)
668            IF (ndebsen ==0) THEN
669                ndebsen = n
670    !           dltamsen = histgrowth(hist_sencour, 2) * P_ratiosen *
671    !           pfeuilverte(nsencour)
672                !!! dltamsen needs to reconsidered since pfeuilverte is not yet
673                !in  histgrowth
674                IF (P_codlainet == 2) THEN
675                    dltaisen = histgrowth(hist_sencour,1) ! + P_lai0
676                    !!! initial lai will never be senescened because we do not
677                    !know when it is initiated
678                ENDIF
679            ELSE
680                ! somsen = somsenreste
681                somsen = 0
682                dltamsen = 0.
683                IF (P_codlainet == 2) THEN
684                    dltaisen = 0.
685                ENDIF
686                DO ii = hist_sencour+1, n
687                    somsen = somsen + histgrowth(ii,3)
688                ENDDO
689               
690                IF (somsen < histgrowth(hist_sencour+1, 5) ) THEN
691                    dltamsen = 0.
692                    IF (P_codlainet == 2) THEN
693                        dltaisen = 0.
694                    ENDIF
695                ELSE
696                    DO WHILE (somsen >= histgrowth(hist_sencour+1, 5) .and. hist_sencour+1 <= hist_latest)
697                        hist_sencour = hist_sencour + 1
698                        somsen = somsen - histgrowth(hist_sencour,3)!tdevelop(hist_sencour)
699    !                   dltamsen = dltamsen + (histgrowth(hist_sencour, 2) *
700    !                   P_ratiosen * pfeuilverte(nsencour))
701                        !!! dltamsen needs to reconsidered since pfeuilverte is
702                        !not yet in  histgrowth
703                        IF (P_codlainet == 2) THEN
704                            dltaisen = dltaisen + histgrowth(hist_sencour, 1)
705                        ENDIF
706                        somsenreste = somsen
707                    ENDDO
708                    ! somsenreste = somsen - histgrowth(hist_sencour,5)
709                    !IF ( somsenreste<0 ) THEN
710                    !    write(*,*) 'xuhui: somesenreste<0, likely a bug'
711                    !    write(*,*) 'xuhui: somsen : ',somsen, 'durvie(hist_sencour) : ', histgrowth(hist_sencour,5)
712                    !ENDIF
713                ENDIF
714            ENDIF   
715        ENDIF
716
717        IF (myprint) THEN
718            write(*,*) 'xuhui: ---------------------'
719            write(*,*) 'xuhui: day, ', n
720            write(*,*) 'xuhui: hist_sencour, ',hist_sencour
721            write(*,*) 'xuhui: hist_latest, ',hist_latest
722            write(*,*) 'xuhui: dltaisen, ',dltaisen
723            write(*,*) 'xuhui: dltai, ',histgrowth(hist_latest,1)
724            write(*,*) 'xuhui: somtemp, ',somtemp
725            write(*,*) 'xuhui: durvie, ',durvie
726            write(*,*) 'xuhui: somsen, ',somsen
727            write(*,*) 'xuhui: ---------------------'
728        ENDIF
729    elseif (P_codlainet==3) then
730        IF (ulai <= P_vlaimax) THEN
731            durage = durvieI
732        ELSE
733            durage = durvieI + ((ulai - P_vlaimax) / (3. -P_vlaimax)) * (P_durvieF - durvieI)
734        ENDIF
735        durvie = durage * min(senstress, fstressgel)
736        do ii = 1, nbox
737            if (boxndays(ii) == 0) then
738                ! do nothing because of no leaf
739            else
740                if (boxdurage(ii) == 0.) then
741                    boxdurage(ii) = durvie
742                else
743                    boxdurage(ii) = min(durvie, boxdurage(ii))
744                endif
745            endif
746        enddo
747        ! senescence condition
748        dltaisen = 0.
749        dltamsen = 0.
750        do ii = 1, nbox
751            if (boxndays(ii) > 0 .and. somtemp - boxsomsenbase(ii)  >= boxdurage(ii)) then !senescence part of this compartment
752                if (ndebsen ==0) then
753                    ndebsen = n
754                    write(*,*) 'xuhui: debsen, ', ndebsen
755                endif
756!                dlairat = boxlai(ii)/real(boxndays(ii))
757!                damsrat = boxbiom(ii)/real(boxndays(ii))
758                dtdevrat = boxtdev(ii)/real(boxndays(ii))
759               
760                tempdays = (somtemp - boxsomsenbase(ii) - boxdurage(ii)) / dtdevrat
761                if (boxlairem(ii) >= 0.) then ! we have at least sth to be senesced
762                    if (tempdays >= boxndays(ii)) then
763                        templaisen = boxlai(ii)
764                    else
765                        templaisen = tempdays/real(boxndays(ii))*boxlai(ii)
766                    endif
767                    if (templaisen + boxlairem(ii) - boxlai(ii) > 0.) then
768                        dltaisen = dltaisen + templaisen + boxlairem(ii) - boxlai(ii)
769                        boxlairem(ii) = boxlairem(ii) - dltaisen
770                    else
771                        dltaisen = dltaisen + 0.
772                    endif
773                endif
774                if (boxbiomrem(ii) >= 0.) then
775                    if (tempdays >= boxndays(ii)) then 
776                        !leaf:biomass ratio (feuilverte) also change with time,
777                        ! considering including a compartment feuilverte in
778                        ! later version
779                        tempamsen = boxbiom(ii) * P_ratiosen * nsenpfeuilverte 
780                    else
781                        tempamsen = tempdays/real(boxndays(ii)) * boxbiom(ii) * P_ratiosen * nsenpfeuilverte
782                    endif
783                    if (tempamsen + (boxbiomrem(ii) - boxbiom(ii)) * P_ratiosen * nsenpfeuilverte > 0.) then
784                        dltamsen = dltamsen + tempamsen + (boxbiomrem(ii) - boxbiom(ii)) * P_ratiosen * nsenpfeuilverte
785                    else
786                        dltamsen = dltamsen + 0.
787                    endif
788                endif
789            endif
790        enddo
791
792        IF (myprint) THEN
793            write(*,*) 'xuhui: ---------------------'
794            write(*,*) 'xuhui: day, ', n
795            write(*,*) 'xuhui: ulai, ',ulai
796            write(*,*) 'xuhui: dltai, ',deltai
797            write(*,*) 'xuhui: dltaisen, ',dltaisen
798            write(*,*) 'xuhui: somtemp, ',somtemp
799            write(*,*) 'xuhui: boxndays, ',boxndays
800            write(*,*) 'xuhui: boxdurage, ',boxdurage
801            write(*,*) 'xuhui: boxlai, ',boxlai
802            write(*,*) 'xuhui: boxtdev, ',boxtdev
803            write(*,*) 'xuhui: boxlairem, ',boxlairem
804            write(*,*) 'xuhui: boxsomsenbase, ', boxsomsenbase
805    !        write(*,*) 'xuhui: somsen, ',somsen
806            write(*,*) 'xuhui: ---------------------'
807        ENDIF
808    endif
809   
81030        continue
811      ! **      msresjaune   = matière sénescente résiduelle
812      ! *-      msneojaune   = matière sénescente néoformée (limitée au végétatif)
813      ! *-      mafeuiljaune = matière sénescente cumulée
814      ! *- Hyp: msneojaune   = mafeuiljaune
815     
816      !msresjaune = msresjaune + dltamsres
817
818      ! PB&Inaki - 08/03/2005 : On enlève les feuilles tombées de la matière jaune
819     
820      ! the purpose of dltamstombe is for recycling of biomass, but at this moment we do not consider this process---xcw
821     
822      !mafeuiljaune = mafeuiljaune + dltamsen - dltamstombe
823
824      mafeuiljaune = mafeuiljaune + dltamsen
825      msneojaune = mafeuiljaune
826
827! Storage of the laisen information, assignment of the pervious day laisen
828     
829      pdlaisen = laisen
830
831return
832end subroutine senescen
833 
834 
Note: See TracBrowser for help on using the repository browser.