source: branches/publications/ORCHIDEE_GLUC_r6545/src_sticslai/reprac_calc.f90 @ 6737

Last change on this file since 6737 was 3833, checked in by albert.jornet, 8 years ago

Fix: all getin_in calls must include a default value before being called. An exception might raise.
Fix: nstoprac variable is changed from Real to Integer to be consistent. A real is passed as an argument when an integer is expected in the subroutine signature.
Fix: masec variable changed from IN to INOUT subroutines declaration. There is an inconsistency when in the same subroutine is passed as INOUT (can be modified) to another subroutine. No variable modification was expected.
Clean: delete print message in hydrol

File size: 10.5 KB
Line 
1!>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>!
2!> Purpose of this subroutine: addressing the shoot/root partitioning <!
3!> Author: Xiuchen Wu                                                 <!
4!> Date  : 31/07/2013                                                 <!
5
6
7subroutine reprac_calc(n, nsen, nlax, nflo, nmat, nrec, nplt, in_cycle, nger, tcult, nlev, &      ! IN
8                     zrac, repracmax, repracmin, kreprac, somtemprac, urac, reprac, nstoprac, &     ! inout
9                  P_stoprac, P_codeperenne, P_codehypo, P_codegermin,P_zracplantule, P_profsem, P_codtrophrac,  &
10                  P_repracpermax, P_repracpermin, P_krepracperm, P_repracseumax, P_repracseumin, P_krepracseu,  &
11                  P_codetemprac, P_tcmin, P_tcmax, P_codedyntalle, P_tcxstop, P_stlevamf, P_stamflax)     ! parameter
12
13! USE Stics
14
15! DECLARATION PART
16!
17! 0.0  INPUT 
18
19
20integer, intent(IN)                :: n
21integer, intent(IN)                :: nsen
22integer, intent(IN)                :: nlax
23integer, intent(IN)                :: nflo
24integer, intent(IN)                :: nmat
25integer, intent(IN)                :: nrec
26integer, intent(IN)                :: nplt
27logical, intent(IN)                :: in_cycle
28integer, intent(IN)                :: nger
29real,    intent(IN)                :: tcult
30integer, intent(IN)                :: nlev
31character, intent(in)              :: P_stoprac !//parameter, stage when root stop growth
32integer, intent(in)                :: P_codeperenne   ! // parameter
33integer, intent(in)                :: P_codehypo   ! // parameter
34
35integer, intent(in)                :: P_codegermin   ! // parameter
36real,    intent(in)                :: P_zracplantule ! // parameter
37real,    intent(in)                :: P_profsem      ! // parameter
38integer, intent(in)                :: P_codtrophrac   ! // parameter
39
40real,    intent(in)                :: P_repracpermax     ! // parameter
41real,    intent(in)                :: P_repracpermin      ! // parameter
42real,    intent(in)                :: P_krepracperm      ! // parameter
43
44real,    intent(in)                :: P_repracseumax     ! // parameter
45real,    intent(in)                :: P_repracseumin      ! // parameter
46real,    intent(in)                :: P_krepracseu      ! // parameter
47
48integer, intent(in)                :: P_codetemprac   ! // parameter
49real,    intent(in)                :: P_tcmin      ! // parameter
50real,    intent(in)                :: P_tcmax      ! // parameter
51integer, intent(in)                :: P_codedyntalle   ! // parameter
52real,    intent(in)                :: P_tcxstop      ! // parameter
53
54real,    intent(in)                :: P_stlevamf      ! // parameter
55real,    intent(in)                :: P_stamflax     ! // parameter
56
57
58
59!
60! 1.0  INOUT
61
62real,    intent(INOUT)             :: zrac
63real,    intent(INOUT)             :: repracmax
64real,    intent(INOUT)             :: repracmin
65real,    intent(INOUT)             :: kreprac
66real,    intent(INOUT)             :: somtemprac
67real,    intent(INOUT)             :: urac
68real,    intent(INOUT)             :: reprac
69integer, intent(INOUT)             :: nstoprac
70
71
72
73! 2.0  OUT
74
75! 3.0  LOCAL
76
77real :: znonli
78real :: dtj
79
80
81! zrac should be initialized as 0
82! when nger ~= 0, meaning that the germination occurs, we assign the zrac as P_profsem
83
84
85     !>> FIRST
86     ! we determine the nstoprac, this variable is dependent on P_stoprac
87
88     if (P_stoprac == 'sen' .or. P_stoprac == 'SEN') nstoprac = nsen
89     if (P_stoprac == 'lax' .or. P_stoprac == 'LAX') nstoprac = nlax
90     if (P_stoprac == 'flo' .or. P_stoprac == 'FLO') nstoprac = nflo
91     if (P_stoprac == 'mat' .or. P_stoprac == 'MAT') nstoprac = nmat
92     if (P_stoprac == 'rec' .or. P_stoprac == 'REC') nstoprac = nrec
93
94
95     !>> SECOND
96     ! adjustment of the root depth accoring to the stages
97     ! but note that in this version, we do not calculate the growth of root
98
99     if (P_codeperenne == 1 .and. (n >= nplt .or. in_cycle == .true.)) then ! annual crop and in the growing cycle
100       if (P_codehypo == 2 .and. P_codegermin == 1) then
101         znonli = P_zracplantule + P_profsem
102         zrac = znonli
103         !return
104       else
105         znonli = P_profsem
106         zrac = znonli
107         !return
108       endif
109     else
110       znonli = 0.
111       zrac = znonli
112     endif
113     ! now we adjust the zrac
114
115     
116     !>> THIRD
117     ! we get the three parameters for calculating the REPRAC variable
118     ! There are three cases:
119     ! 1. codtrophrac is 1
120     ! 2. codtrophrac is 2
121     ! 3. codtrophrac is 3
122
123     if (P_codtrophrac == 1) then
124       repracmax = P_repracpermax
125       repracmin = P_repracpermin
126       kreprac   = P_krepracperm
127     endif
128     if (P_codtrophrac == 2) then
129       repracmax = P_repracseumax
130       repracmin = P_repracseumin
131       kreprac   = P_krepracseu
132     endif
133     if (P_codtrophrac == 3) then
134       repracmax = 0.
135       repracmin = 0.
136       kreprac   = 0.
137     endif
138
139 
140     !>> FOURTH
141     ! if the germination did not occur, we do not need the shoot/root partitioning calculation, so the reprac keeps the initialization value
142     if (nger == 0) return
143   
144
145     
146     !>> FIFTH
147     ! we calculate the reprac
148
149
150     !if (n >= nger .and. nger > 0 .and. zrac > 0.) then
151     if (in_cycle == .true. .and. nger > 0 .and. zrac > 0) then ! if in growing cycle, and already germinate and depth of root front is larger than 0
152       if (P_codetemprac == 1) then ! using the crop temperature
153         dtj = max(tcult - P_tcmin, 0.)
154         dtj = min(dtj, P_tcmax - P_tcmin)
155
156       ! DR et ML et SYL 15/06/09
157       ! ************************
158       ! introduction de la fin des modifications de Sylvain (nadine et FR)
159       ! dans le cadre du projet PERMED
160         if (P_codedyntalle == 1) then
161       ! ####
162       ! NB le 06/03/2008 introduction de l'effet négatif des températures élevées
163       ! sur les fonctions racinaires qui utilisent "dtj"
164           if (tcult > P_tcmax .and. P_tcxstop < 100) then
165             dtj = (P_tcmax-P_tcmin)/(-P_tcxstop+P_tcmax)*(tcult-P_tcxstop)
166             dtj = max(dtj,0.)
167           endif
168         endif
169       ! ####
170
171       !else ! using the soil temperature
172       !  dtj(n) = max(tsol(int(zrac)) - P_tgmin, 0.)
173       !  dtj(n) = min(dtj(n) ,P_tcmax - P_tgmin)
174
175       !! ####
176       !! NB le 06/03/2008 introduction de l'effet négatif des températures élevées
177       !! sur les fonctions racinaires qui utilisent "dtj"
178       !  if(P_codedyntalle == 1)then
179       !    if (tsol(int(zrac)) > P_tcmax .and. P_tcxstop < 100) then
180       !      dtj(n) = (P_tcmax-P_tcmin)/(-P_tcxstop+P_tcmax)*(tsol(int(zrac))-P_tcxstop)
181       !      dtj(n) = max(dtj(n),0.)
182       !    endif
183       !  endif
184       !! ####
185
186       endif
187     else
188       dtj = 0.
189     endif
190
191      ! addressing the cumulated effective temperature of root
192
193
194      somtemprac = somtemprac + dtj
195      if (n == nlev) somtemprac = 0.
196      if (nstoprac == 0) then
197        !: Définition d'une unité de dl racinaire urac
198        !- somme de degré.jours racine
199        urac = min(1. + (2. * somtemprac / (P_stlevamf + P_stamflax)), 3.)
200        if (nlev == 0) urac = 1.
201        !: la longueur totale de racine émise par jour (rlj)
202        !- suit une logistique de façon similaire au LAI
203        !- NB - le 28/10/01:
204        !- 07/02/08: on integre ENFIN les modifs de Samuel (vive lui)
205        !- essai de rendre nouvrac continu pour Samuel
206        !- 20/02/07:
207        !--if (int(zrac)-int(zrac-deltaz) <= int(deltaz)) then
208        !--  multlvfront = int(deltaz)
209        !--else
210        !--  multlvfront = int(deltaz)+1
211        !--endif
212        !--if (nlev == 0) multlvfront = max0(1,multlvfront)
213        !--nouvrac = P_lvfront*multlvfront
214        ! nouvrac = P_lvfront*deltaz
215
216        !: Introduction de l'indice de stress de densité racinaire idzrac / NB - le 06/06
217        !efanoxd = 1.-(1.-idzrac)*P_sensanox
218
219        !: NB - 10/03/02:
220        !- Calcul du facteur densité efdensite actif à partir de P_laicomp
221        !- Calcul de l'effet densité sur la mise en place du LAI pour les stics-plante
222        !efdensite = 1.
223        !if (urac >= 1.) then
224        !  if (lai_veille < P_laicomp) then
225        !    efdensite = 1.
226        !  else
227        !    !: domi - 02/07/2002: pb si GEL total densite = 0 et pb de log(0)
228        !    if (densite == 0) then
229        !      efdensite = 0.
230        !    else
231        !      efdensite = min(exp(P_adens * (log(densite / P_bdens))), 1.)
232        !    endif
233        !  endif
234        !else
235        !  efdensite = min(exp(P_adens * (log(densite / P_bdens))),1.)
236        !endif
237
238        !rlj = (P_draclong / (1. + exp(5.5 * (P_vlaimax - urac))) * efdensite * densite *dtj(n) * efanoxd) + (nouvrac * 1.e4)
239
240        !: De la germination à la levée, il n'y a qu'une croissance du front racinaire
241        !- Nb - le 17/02/2003: suppression de coderac
242        if (P_codtrophrac /= 3) then
243          !: Option trophique
244          ! *- par calcul d'une fonction de répartition reprac = souterrain/total
245          ! *- la longueur de racine au niveau du front est soustraite de la
246          ! *- longueur produite
247          ! ** ML le 25/05/2007 - les paramÚtres initiaux repracmin, repracmax et kreprac
248          ! *- ont changé de noms dans le fichier plante: P_repracpermin, P_repracpermax et P_krepracperm
249          ! *- si la liaison trophique est permanente (P_codtrophrac = 1) et
250          ! *- P_repracseumin, P_repracseumax et P_krepracseu si la liaison est par seuils (P_codtrophrac = 2)
251!!!          if (P_codtrophrac == 1) then
252!!!            repracmax = P_repracpermax
253!!!            repracmin = P_repracpermin
254!!!            kreprac = P_krepracperm
255!!!          endif
256!!!
257!!!          if (P_codtrophrac == 2) then
258!!!            repracmax = P_repracseumax
259!!!            repracmin = P_repracseumin
260!!!            kreprac = P_krepracseu
261!!!          endif
262
263          reprac = (repracmax-repracmin) * (exp(-kreprac * (urac - 1.))) + repracmin
264
265          !: PB - 06/01/2005: On multiplie par P_longsperac plutot que diviser (modifs Nadine du 30/12/2004)
266          !rlj1 = reprac/(1.-reprac) * P_longsperac * 1.e2 * dltams
267
268          !if (rlj1 < nouvrac * 1.e4) nouvrac = rlj1 * 1.e-4
269          !if (P_codtrophrac == 1) then
270          !  rlj = rlj1
271          !else
272          !  if (rlj >= rlj1) rlj = rlj1
273          !endif
274        endif
275
276        !if (nlev == 0) rlj = nouvrac*1.e4
277
278      !else
279      !  rlj = 0.
280      endif
281     
282     
283      ! Here we get the reprac, which is an important variable addressing the above/below ground biomass partitioning.   
284      ! we end our subroutine
285
286end subroutine reprac_calc
Note: See TracBrowser for help on using the repository browser.