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

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

Clean: simplify some stics related variables
New: introduce an interface for restput and restget for integer type variables

File size: 24.7 KB
Line 
1! INITIALIZATION OF VARIABLES REGARDING THE CROP CYCLE
2! 27/06/2013---xcw
3
4subroutine Stics_init(&
5               kjpindex        ,&
6               !nvm             ,&
7               f_crop_init     ,&   
8               f_crop_recycle      ,&   
9               in_cycle                ,&   
10               f_sen_lai                ,&   
11               onarretesomcourdrp      ,&   
12!               nlevobs                 ,&
13!               namfobs                 ,&
14!               nfloobs                 ,&
15!               nlanobs                 ,&
16!               nlaxobs                 ,&
17!               nmatobs                 ,&
18!               nrecobs                 ,&
19!               nsenobs                 ,&
20!               ndrpobs                 ,&
21               nsendltams              ,&
22               nsendltai               ,&
23               nsenpfeuilverte         ,&
24               nsendurvie              ,&
25               nsenndurvie             ,&
26               densiteequiv            ,&
27               nplt                    ,&
28               tursla                  ,&
29               ssla                     ,&
30               pfeuilverte             ,&
31               bsenlai                 ,&
32               zrac                    ,&
33               nrec                    ,& 
34               nlan                    ,&
35               tcult                   ,&
36               udevair                 ,&
37               udevcult                ,&
38               ndrp                    ,&
39               rfvi                    ,&
40               nlev                    ,&
41               nger                    ,&
42               etatvernal              ,&
43               caljvc                  ,&
44               rfpi                    ,&
45               upvt                    ,&
46               utp                     ,&
47               somcour                 ,&
48               somcourdrp              ,&
49               somcourutp              ,&
50               tdevelop                ,&
51               somtemp                 ,&
52               somcourfauche           ,&
53               stpltger                ,&
54               R_stamflax              ,&
55               R_stlaxsen              ,&
56               R_stsenlan              ,&
57               stlevflo                ,&
58               nflo                    ,&
59               R_stlevdrp              ,&
60               R_stflodrp              ,&
61               R_stdrpmat              ,&
62               nmat                    ,&
63               nlax                    ,&
64               nrecbutoir              ,&
65               group                   ,&
66               ndebdes                 ,&
67               R_stdrpdes              ,&
68               densite                 ,&
69               densitelev              ,&
70               coeflev                 ,&
71               densiteger              ,&
72               somelong                 ,&
73               somger                  ,&
74               humectation             ,&
75               nbjhumec                ,&
76               somtemphumec            ,&
77               stpltlev                ,&
78               namf                    ,&
79               stmatrec                ,&
80               tustress                ,&
81               slai                     ,&
82               somfeuille              ,&
83               pdlai                   ,&
84               nbfeuille               ,&
85               reajust                 ,&
86               ulai                    ,&
87               pdulai                  ,&
88               efdensite               ,&
89               tempeff                 ,&
90               nstopfeuille            ,&
91               deltai                  ,&
92               svmax                    ,&
93               nsen                    ,&
94               laisen                  ,&
95               pdlaisen                ,&
96               dltaisenat              ,&
97               nsencour                ,&
98               dltamsen                ,&
99               dltaisen                ,&
100               fgellev                 ,&
101               gelee                   ,&
102               fstressgel              ,&
103               R_stlevamf              ,&
104               dernier_n               ,&
105               durvieI                 ,&
106               durvie                  ,&
107               ndebsen                 ,&
108               somsenreste             ,&
109               shumrel                  ,&
110               swfac                   ,&
111               turfac                  ,&
112               senfac                  ,& 
113               mafeuiljaune            ,&
114               msneojaune              ,&             
115               v_dltams                ,&
116               fgelflo                 ,&
117               pdircarb                ,&
118               ircarb                  ,&
119               nbgrains                ,&
120               pgrain                  ,&
121               vitmoy                  ,&
122               nbgraingel              ,&
123               pgraingel               ,&
124               dltags                  ,&
125               ftempremp               ,&
126               magrain                 ,&
127               pdmagrain               ,&
128               nbj0remp                ,&
129               pdsfruittot             ,&
130               repracmax               ,&
131               repracmin               ,&
132               kreprac                 ,&
133               somtemprac              ,&
134               urac                    ,&
135               reprac                  ,&
136               nstoprac                ,&
137               c_reserve               ,&
138               c_leafb                 ,&
139               !biomass                 ,&
140               deltgrain               ,&
141               gslen                   ,&
142               drylen, &
143               histgrowthset, hist_sencourset, hist_latestset, doyhiststset, &
144               nboxmax,box_ulai, box_ndays, box_lai, box_lairem, box_tdev, box_biom, box_biomrem, box_durage, box_somsenbase)               
145
146
147  !USE stomate_io
148  USE pft_parameters
149  USE constantes
150
151
152  implicit none
153
154  ! DECLARATION
155  integer,intent(in)                             ::         kjpindex
156  !integer,intent(in)                             ::         nvm
157  integer,intent(in)                             ::         nboxmax
158  logical, intent(in)::         f_crop_init       
159  logical,dimension(kjpindex, nvm), intent(inout)::         f_crop_recycle         
160  logical,dimension(kjpindex, nvm), intent(out)::         in_cycle         
161  logical,dimension(kjpindex, nvm), intent(out)::         f_sen_lai         
162  logical,dimension(kjpindex, nvm), intent(out)::         onarretesomcourdrp         
163!  integer,dimension(kjpindex, nvm), intent(out)::         nlevobs                 
164!  integer,dimension(kjpindex, nvm), intent(out)::         namfobs                 
165!  integer,dimension(kjpindex, nvm), intent(out)::         nfloobs                 
166!  integer,dimension(kjpindex, nvm), intent(out)::         nlanobs                 
167!  integer,dimension(kjpindex, nvm), intent(out)::         nlaxobs                 
168!  integer,dimension(kjpindex, nvm), intent(out)::         nmatobs                 
169!  integer,dimension(kjpindex, nvm), intent(out)::         nrecobs                 
170!  integer,dimension(kjpindex, nvm), intent(out)::         nsenobs                 
171!  integer,dimension(kjpindex, nvm), intent(out)::         ndrpobs                 
172  !
173  real,dimension(kjpindex, nvm), intent(out)::         nsendltams             
174  real,dimension(kjpindex, nvm), intent(out)::         nsendltai               
175  real,dimension(kjpindex, nvm), intent(out)::         nsenpfeuilverte         
176  real,dimension(kjpindex, nvm), intent(out)::         nsendurvie             
177  real,dimension(kjpindex, nvm), intent(out)::         nsenndurvie             
178  real,dimension(kjpindex, nvm), intent(out)::         densiteequiv           
179  integer,dimension(kjpindex, nvm), intent(out)::         nplt                   
180  real,dimension(kjpindex, nvm), intent(out)::         tursla                 
181  real,dimension(kjpindex, nvm), intent(out)::         ssla                     
182  real,dimension(kjpindex, nvm), intent(out)::         pfeuilverte             
183  real,dimension(kjpindex, nvm), intent(out)::         bsenlai                 
184 
185  ! STICS::LAIdev::DEVELOPMENT
186  real,dimension(kjpindex, nvm), intent(out)::         zrac                   
187  integer,dimension(kjpindex, nvm), intent(out)::         nrec                     
188  integer,dimension(kjpindex, nvm), intent(out)::         nlan                   
189  real,dimension(kjpindex, nvm), intent(out)::         tcult                   
190  real,dimension(kjpindex, nvm), intent(out)::         udevair                 
191  real,dimension(kjpindex, nvm), intent(out)::         udevcult               
192  integer,dimension(kjpindex, nvm), intent(out)::         ndrp                   
193  real,dimension(kjpindex, nvm), intent(out)::         rfvi                   
194  integer,dimension(kjpindex, nvm), intent(out)::         nlev                   
195  integer,dimension(kjpindex, nvm), intent(out)::         nger                   
196  logical,dimension(kjpindex, nvm), intent(out)::         etatvernal             
197  real,dimension(kjpindex, nvm), intent(out)::         caljvc                 
198  real,dimension(kjpindex, nvm), intent(out)::         rfpi                   
199  real,dimension(kjpindex, nvm), intent(out)::         upvt                   
200  real,dimension(kjpindex, nvm), intent(out)::         utp                     
201  real,dimension(kjpindex, nvm), intent(out)::         somcour                 
202  real,dimension(kjpindex, nvm), intent(out)::         somcourdrp             
203  real,dimension(kjpindex, nvm), intent(out)::         somcourutp             
204  real,dimension(kjpindex, nvm), intent(out)::         tdevelop               
205  real,dimension(kjpindex, nvm), intent(out)::         somtemp                 
206  real,dimension(kjpindex, nvm), intent(out)::         somcourfauche           
207  real,dimension(kjpindex, nvm), intent(out)::         stpltger               
208  real,dimension(kjpindex, nvm), intent(out)::         R_stamflax             
209  real,dimension(kjpindex, nvm), intent(out)::         R_stlaxsen             
210  real,dimension(kjpindex, nvm), intent(out)::         R_stsenlan             
211  real,dimension(kjpindex, nvm), intent(out)::         stlevflo               
212  integer,dimension(kjpindex, nvm), intent(out)::         nflo                   
213  real,dimension(kjpindex, nvm), intent(out)::         R_stlevdrp             
214  real,dimension(kjpindex, nvm), intent(out)::         R_stflodrp             
215  real,dimension(kjpindex, nvm), intent(out)::         R_stdrpmat             
216  integer,dimension(kjpindex, nvm), intent(out)::         nmat                   
217  integer,dimension(kjpindex, nvm), intent(out)::         nlax                   
218  integer,dimension(kjpindex, nvm), intent(out)::         nrecbutoir             
219  real,dimension(kjpindex, nvm), intent(out)::         group                   
220  integer,dimension(kjpindex, nvm), intent(out)::         ndebdes                 
221  real,dimension(kjpindex, nvm), intent(out)::         R_stdrpdes             
222  real,dimension(kjpindex, nvm), intent(out)::         densite                 
223  real,dimension(kjpindex, nvm), intent(out)::         densitelev
224  real,dimension(kjpindex, nvm), intent(out)::         coeflev
225             
226  real,dimension(kjpindex, nvm), intent(out)::         densiteger             
227  real,dimension(kjpindex, nvm), intent(out)::         somelong                 
228  real,dimension(kjpindex, nvm), intent(out)::         somger                 
229  logical,dimension(kjpindex, nvm), intent(out)::         humectation             
230  integer,dimension(kjpindex, nvm), intent(out)::         nbjhumec               
231  real,dimension(kjpindex, nvm), intent(out)::         somtemphumec           
232  real,dimension(kjpindex, nvm), intent(out)::         stpltlev               
233  integer,dimension(kjpindex, nvm), intent(out)::         namf                   
234  real,dimension(kjpindex, nvm), intent(out)::         stmatrec               
235  ! STICS::LAIdev:: LAI calculation
236  real,dimension(kjpindex, nvm), intent(out)::         tustress               
237  real,dimension(kjpindex, nvm), intent(out)::         slai                     
238  real,dimension(kjpindex, nvm), intent(out)::         somfeuille             
239  real,dimension(kjpindex, nvm), intent(out)::         pdlai                   
240  integer,dimension(kjpindex, nvm), intent(out)::         nbfeuille               
241  real,dimension(kjpindex, nvm), intent(out)::         reajust                 
242  real,dimension(kjpindex, nvm), intent(out)::         ulai                   
243  real,dimension(kjpindex, nvm), intent(out)::         pdulai                 
244  real,dimension(kjpindex, nvm), intent(out)::         efdensite               
245  real,dimension(kjpindex, nvm), intent(out)::         tempeff                 
246  integer,dimension(kjpindex, nvm), intent(out)::         nstopfeuille           
247  real,dimension(kjpindex, nvm), intent(out)::         deltai                 
248  real,dimension(kjpindex, nvm), intent(out)::         svmax                   
249  integer,dimension(kjpindex, nvm), intent(out)::         nsen                   
250  real,dimension(kjpindex, nvm), intent(out)::         laisen                 
251  real,dimension(kjpindex, nvm), intent(out)::         pdlaisen               
252  real,dimension(kjpindex, nvm), intent(out)::         dltaisenat             
253  ! STICS:: LAIdev:: LAI senescence
254  integer,dimension(kjpindex, nvm), intent(out)::         nsencour               
255  real,dimension(kjpindex, nvm), intent(out)::         dltamsen               
256  real,dimension(kjpindex, nvm), intent(out)::         dltaisen               
257  real,dimension(kjpindex, nvm), intent(out)::         fgellev                 
258  logical,dimension(kjpindex, nvm), intent(out)::         gelee                   
259  real,dimension(kjpindex, nvm), intent(out)::         fstressgel             
260  real,dimension(kjpindex, nvm), intent(out)::         R_stlevamf             
261  integer,dimension(kjpindex, nvm), intent(out)::         dernier_n               
262  real,dimension(kjpindex, nvm), intent(out)::         durvieI                 
263  real,dimension(kjpindex, nvm), intent(out)::         durvie                 
264  integer,dimension(kjpindex, nvm), intent(out)::         ndebsen                 
265  real,dimension(kjpindex, nvm), intent(out)::         somsenreste
266  ! STICS:: LAIdev:: STRESS             
267  real,dimension(kjpindex, nvm), intent(out)::         shumrel                 
268  real,dimension(kjpindex, nvm), intent(out)::         swfac                   
269  real,dimension(kjpindex, nvm), intent(out)::         turfac                 
270  real,dimension(kjpindex, nvm), intent(out)::         senfac               
271 
272  real,dimension(kjpindex, nvm), intent(out)::         mafeuiljaune               
273  real,dimension(kjpindex, nvm), intent(out)::         msneojaune               
274  ! STICS:: CARBON ALLOCATION
275
276  real,dimension(kjpindex, nvm, 60), intent(out)::         v_dltams               
277  real,dimension(kjpindex, nvm), intent(out)::         fgelflo               
278  real,dimension(kjpindex, nvm), intent(out)::         pdircarb               
279  real,dimension(kjpindex, nvm), intent(out)::         ircarb               
280  real,dimension(kjpindex, nvm), intent(out)::         nbgrains               
281  real,dimension(kjpindex, nvm), intent(out)::         pgrain             
282  real,dimension(kjpindex, nvm), intent(out)::         vitmoy               
283  real,dimension(kjpindex, nvm), intent(out)::         nbgraingel               
284  real,dimension(kjpindex, nvm), intent(out)::         pgraingel               
285  real,dimension(kjpindex, nvm), intent(out)::         dltags               
286  real,dimension(kjpindex, nvm), intent(out)::         ftempremp               
287  real,dimension(kjpindex, nvm), intent(out)::         magrain               
288  real,dimension(kjpindex, nvm), intent(out)::         pdmagrain               
289  integer,dimension(kjpindex, nvm), intent(out)::         nbj0remp                 
290  real,dimension(kjpindex, nvm), intent(out)::         pdsfruittot               
291  real,dimension(kjpindex, nvm), intent(out)::         repracmax               
292  real,dimension(kjpindex, nvm), intent(out)::         repracmin               
293  real,dimension(kjpindex, nvm), intent(out)::         kreprac               
294  real,dimension(kjpindex, nvm), intent(out)::         somtemprac               
295  real,dimension(kjpindex, nvm), intent(out)::         urac             
296  real,dimension(kjpindex, nvm), intent(out)::         reprac               
297
298  integer,dimension(kjpindex, nvm), intent(out)::         nstoprac                 
299  real,dimension(kjpindex, nvm), intent(out)::         c_reserve               
300  real,dimension(kjpindex, nvm), intent(out)::         c_leafb           
301  !real,dimension(kjpindex, nvm, nparts ), intent(out)::         biomass           
302  real,dimension(kjpindex, nvm), intent(out)::         deltgrain         
303  integer,dimension(kjpindex, nvm), intent(out)::         gslen                 
304  integer,dimension(kjpindex, nvm), intent(out)::         drylen                 
305  real,dimension(kjpindex, nvm, 300, 5), intent(out) :: histgrowthset
306  integer, dimension(kjpindex,nvm), intent(out) :: hist_sencourset
307  integer, dimension(kjpindex,nvm), intent(out) :: hist_latestset
308  integer, dimension(kjpindex,nvm), intent(out) :: doyhiststset
309
310  integer, dimension(kjpindex,nvm,nboxmax), intent(out) :: box_ndays
311  real, dimension(kjpindex,nvm,nboxmax), intent(out)    :: box_lai
312  real, dimension(kjpindex,nvm,nboxmax), intent(out)    :: box_lairem
313  real, dimension(kjpindex,nvm,nboxmax), intent(out)    :: box_tdev
314  real, dimension(kjpindex,nvm,nboxmax), intent(out)    :: box_biom
315  real, dimension(kjpindex,nvm,nboxmax), intent(out)    :: box_biomrem
316  real, dimension(kjpindex,nvm,nboxmax), intent(out)    :: box_durage
317  real, dimension(kjpindex,nvm,nboxmax), intent(out)    :: box_somsenbase
318  real, dimension(nvm,nboxmax), intent(out)             :: box_ulai
319  ! local variables
320  integer    :: j, ip, k, mid
321  real       :: uinflex
322
323 
324  !
325  DO j= 1, nvm
326     
327     DO ip = 1, kjpindex 
328       IF (f_crop_init .OR. f_crop_recycle(ip, j)) THEN 
329         
330
331          ! LAIdev :: FORCED RUN VARIABLES
332          in_cycle(ip, j) = .FALSE.
333          f_sen_lai(ip, j) = .TRUE.
334          onarretesomcourdrp(ip, j) = .TRUE.
335!          nlevobs(ip, j) = 999
336!          namfobs(ip, j) = 999   
337!          nfloobs(ip, j) = 999
338!          nlanobs(ip, j) = 999
339!          nlaxobs(ip, j) = 999
340!          nmatobs(ip, j) = 999
341!          nrecobs(ip, j) = 999
342!          nsenobs(ip, j) = 999
343!          ndrpobs(ip, j) = 999
344       
345          ! LAIdev ::  SPECIFIC
346          nsendltams(ip, j) = 0.
347          nsendltai(ip, j) = 0.
348          nsenpfeuilverte(ip, j) = 0.
349          nsendurvie(ip, j) = 0.
350          nsenndurvie(ip, j) = 0.
351          densiteequiv(ip, j) = 0.
352          nplt(ip, j) = 0
353          tursla(ip, j) = 1.
354          ssla(ip, j) = 0.
355          pfeuilverte(ip, j) = 0.
356          bsenlai(ip, j) = 0.
357
358          ! LAIdev ::  DEVELOPMENT
359          zrac(ip, j) = 0.
360          nrec(ip, j) = 0
361          nlan(ip, j) = 0
362          tcult(ip, j) = 0.
363          udevair(ip, j) = 0.
364          udevcult(ip, j) = 0.
365          ndrp(ip, j) = 0
366          rfvi(ip, j) = 0.
367          nlev(ip, j) = 0
368          nger(ip, j) = 0
369          etatvernal(ip, j) = .FALSE.
370          caljvc(ip, j) = 0.
371          rfpi(ip, j) = 0.
372          upvt(ip, j) = 0.
373          utp(ip, j) = 0.
374          somcour(ip, j) = 0.
375          somcourdrp(ip, j) = 0.
376          somcourutp(ip, j) = 0.
377          tdevelop(ip, j) = 0.
378          somtemp(ip, j) = 0.
379          somcourfauche(ip, j) = 0.
380          stpltger(ip, j) = SP_stpltger(j)
381          R_stamflax(ip, j) = SP_stamflax(j) 
382          R_stlaxsen(ip, j) = SP_stlaxsen(j)
383          R_stsenlan(ip, j) = SP_stsenlan(j)
384          stlevflo(ip, j) = SP_stlevdrp(j) - SP_stflodrp(j)
385          nflo(ip, j) = 0
386          R_stlevdrp(ip, j) = SP_stlevdrp(j)
387          R_stflodrp(ip, j) = SP_stflodrp(j)
388          R_stdrpmat(ip, j) = SP_stdrpmat(j)
389          nmat(ip, j) = 0
390          nlax(ip, j) = 0
391          nrecbutoir(ip, j) = 999
392          group(ip, j) = 0.
393          ndebdes(ip, j) = 0
394          R_stdrpdes(ip, j) = SP_stdrpdes(j)
395          densite(ip, j) = 0.
396          densitelev(ip, j) = 0.
397          coeflev(ip, j) = 1.
398          densiteger(ip, j) = 0.
399          somelong(ip, j) = 0.
400          somger(ip, j) = 0.
401          humectation(ip, j) = .FALSE.
402          nbjhumec(ip, j) = 0
403          somtemphumec(ip, j) = 0.
404          stpltlev(ip, j) = 0.
405          namf(ip, j) = 0
406          stmatrec(ip, j) = 0.
407         
408          ! LAIdev ::  LAI calculation
409          tustress(ip, j) = 1.
410          slai(ip, j) = 0.
411          somfeuille(ip, j) = 0.
412          pdlai(ip, j) = 0.
413          nbfeuille(ip, j) = 0
414          reajust(ip, j) = 0.
415          ulai(ip, j) = 0.
416          pdulai(ip, j) = 0.
417          efdensite(ip, j) = 0.
418          tempeff(ip, j) = 0.
419          nstopfeuille(ip, j) = 0
420          deltai(ip, j) = 0.
421          svmax(ip, j) = 0.
422          nsen(ip, j) = 0
423          laisen(ip, j) = 0.
424          pdlaisen(ip, j) = 0.
425          dltaisenat(ip, j) = 0.
426         
427          ! LAIdev :: LAI senescence
428          nsencour(ip, j) = 0
429          dltamsen(ip, j) = 0.
430          dltaisen(ip, j) = 0.
431          fgellev(ip, j) = 0.
432          gelee(ip, j) = .FALSE.
433          fstressgel(ip, j) = 0.
434          R_stlevamf(ip, j) = SP_stlevamf(j)
435          dernier_n(ip, j) = 0
436          durvieI(ip, j) = 0.
437          durvie(ip, j) = 0.
438          ndebsen(ip, j) = 0
439          somsenreste(ip, j) = 0.
440         
441          IF (any(ok_LAIdev(:)) .AND. any(SP_codlainet(:) == 3)) THEN
442!              write(*,*) 'xuhui: box module to be initialize ', j
443              box_ndays(ip,j,:) = 0
444              box_lai(ip,j,:) = 0.
445              box_lairem(ip,j,:) = 0.
446              box_tdev(ip,j,:) = 0.
447              box_biom(ip,j,:) = 0.
448              box_biomrem(ip,j,:) = 0.
449              box_durage(ip,j,:) = 0.
450              box_somsenbase(ip,j,:) = 0.
451!              write(*,*) 'xuhui: box module initialized ', j
452          ENDIF
453          ! LAIdev :: STRESS
454           
455          shumrel(ip, j) = 0.
456          swfac(ip, j) = 1.
457          turfac(ip, j) = 1.
458          senfac(ip, j) = 1.
459
460          mafeuiljaune(ip, j) = 0.
461          msneojaune(ip, j) = 0.
462          ! STICS: CARBON ALLOCATION
463
464          v_dltams(ip, j, :) = 0.
465          fgelflo(ip, j) = 1.
466          pdircarb(ip, j) = 0.
467          ircarb(ip, j) = 0.
468          nbgrains(ip, j) = 0.
469          pgrain(ip, j) = 0.
470          vitmoy(ip, j) = 0.
471          nbgraingel(ip, j) = 0.
472          pgraingel(ip, j) = 0.
473          dltags(ip, j) = 0.
474          ftempremp(ip, j) = 0.
475          magrain(ip, j) = 0.
476          pdmagrain(ip, j) = 0.
477          nbj0remp(ip, j) = 0
478          pdsfruittot(ip, j) = 0.
479          repracmax(ip, j) = 0.
480          repracmin(ip, j) = 0.
481          kreprac(ip, j) = 0.
482          somtemprac(ip, j) = 0.
483          urac(ip, j) = 0.
484          reprac(ip, j) = 0.
485          nstoprac(ip, j) = 0
486          c_reserve(ip, j) = 0. 
487          c_leafb(ip, j)= 0. 
488          deltgrain(ip, j) = 0.
489          gslen(ip, j) = 0 
490          drylen(ip, j) = 0 
491          ! FINISH THE INITIALIZATION AND RESET THE VALUE
492          IF (any(ok_LAIdev(:)) .AND. any(SP_codlainet(:) == 2)) THEN
493              histgrowthset(ip,j,:,:) = 0.
494              hist_sencourset(ip,j) = 0
495              hist_latestset(ip,j) = 0
496              doyhiststset(ip,j) = 0
497          ENDIF
498          f_crop_recycle(ip, j) = .FALSE.
499
500        ENDIF
501     ENDDO ! ip
502  ENDDO ! j
503   
504!    write(*,*) 'xuhui: initialized boxulai'
505    IF (any(ok_LAIdev(:)) .AND. any(SP_codlainet(:)==3)) THEN
506        box_ulai(:,:) = 0.
507        DO j=1,nvm
508    !        IF (f_crop_init) THEN
509!            IF (f_crop_init .OR. any(f_crop_recycle(:, j))) THEN
510                if (SP_nbox(j) < 2) then
511!                    write(*,*) 'xuhui: box_ulai not allocated for ', j, 'with nbox ', SP_nbox(j)
512                else   
513                    uinflex = SP_vlaimax(j)
514                    mid = floor(real(SP_nbox(j))/2)
515                    if (mid<1) then
516                        mid=1
517                    elseif (mid .eq. SP_nbox(j)) then
518                        mid=SP_nbox(j)-1
519                    endif
520                    do k=1,mid
521                        box_ulai(j,k) = 1+(k-1)*((uinflex-1)/real(mid))
522                    enddo
523                    do k=mid+1,SP_nbox(j)
524                        box_ulai(j,k) = uinflex + (k-mid-1)*(3-uinflex)/real(SP_nbox(j)-mid)
525                    enddo
526!                    write(*,*) 'xuhui: box_ulai for pft ', j, ': ',box_ulai(j,:)
527                endif
528               
529 !           ENDIF
530        ENDDO
531    ENDIF
532
533end subroutine Stics_init
534
Note: See TracBrowser for help on using the repository browser.