-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathharvkillop.f
533 lines (459 loc) · 20 KB
/
harvkillop.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
subroutine harvkillop
!! ~ ~ ~ PURPOSE ~ ~ ~
!! this subroutine performs the harvest and kill operation
!! ~ ~ ~ INCOMING VARIABLES ~ ~ ~
!! name |units |definition
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! auto_eff(:) |none |fertilizer application efficiency calculated
!! |as the amount of N applied divided by the
!! |amount of N removed at harvest
!! bio_hv(:,:,:)|kg/ha |harvested biomass (dry weight)
!! bio_ms(:) |kg/ha |land cover/crop biomass (dry weight)
!! bio_yrms(:) |metric tons/ha |annual biomass (dry weight) in the HRU
!! cnop |none |SCS runoff curve number for moisture
!! |condition II
!! cnyld(:) |kg N/kg yield |fraction of nitrogen in yield
!! cpyld(:) |kg P/kg yield |fraction of phosphorus in yield
!! curyr |none |current year in simulation
!! hi_targ(:,:,:)|(kg/ha)/(kg/ha)|harvest index target of cover defined at
!! |planting
!! hru_dafr(:) |km2/km2 |fraction of watershed in HRU
!! hrupest(:) |none |pesticide use flag:
!! | 0: no pesticides used in HRU
!! | 1: pesticides used in HRU
!! hvsti(:) |(kg/ha)/(kg/ha)|harvest index: crop yield/aboveground
!! |biomass
!! hvstiadj(:) |(kg/ha)/(kg/ha)|optimal harvest index for current time during
!! |growing season
!! icr(:) |none |sequence number of crop grown within the
!! |current year
!! idplt(:) |none |land cover code from crop.dat
!! ihru |none |HRU number
!! ncrops(:,:,:)|
!! npmx |none |number of different pesticides used in
!! |the simulation
!! nro(:) |none |sequence number for year in rotation
!! nyskip |none |number of years to not summarize/print output
!! plantp(:) |kg P/ha |amount of phosphorus in plant biomass
!! plt_et(:) |mm H2O |actual ET simulated during life of plant
!! plt_pet(:) |mm H2O |potential ET simulated during life of plant
!! plt_pst(:,:)|kg/ha |pesticide on plant foliage
!! pltfr_n(:) |none |fraction of plant biomass that is nitrogen
!! rwt(:) |none |fraction of total plant biomass that is
!! |in roots
!! plantn(:) |kg N/ha |amount of nitrogen in plant biomass
!! sol_fon(:,:)|kg N/ha |amount of nitrogen stored in the fresh
!! |organic (residue) pool
!! sol_fop(:,:)|kg P/ha |amount of phosphorus stored in the fresh
!! sol_pst(:,:,1)|kg/ha |pesticide in first layer of soil
!! |organic (residue) pool0
!! sol_rsd(:,:)|kg/ha |amount of organic matter in the soil
!! |classified as residue
!! wshd_yldn |kg N/ha |amount of nitrogen removed from soil in
!! |watershed in the yield
!! wshd_yldp |kg P/ha |amount of phosphorus removed from soil in
!! |watershed in the yield
!! wsyf(:) |(kg/ha)/(kg/ha)|Value of harvest index between 0 and HVSTI
!! |which represents the lowest value expected
!! |due to water stress
!! yldanu(:) |metric tons/ha |annual yield (dry weight) in the HRU
!! yldkg(:,:,:)|kg/ha |yield (dry weight) by crop type in the HRU
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! ~ ~ ~ OUTGOING VARIABLES ~ ~ ~
!! name |units |definition
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! bio_hv(:,:,:)|kg/ha |harvested biomass (dry weight)
!! bio_ms(:) |kg/ha |land cover/crop biomass (dry weight)
!! bio_yrms(:) |metric tons/ha|annual biomass (dry weight) in the HRU
!! hvstiadj(:) |(kg/ha)/(kg/ha)|optimal harvest index for current time during
!! |growing season
!! idorm(:) |none |dormancy status code:
!! |0 land cover growing (not dormant)
!! |1 land cover dormant
!! igro(:) |none |land cover status code:
!! |0 no land cover currently growing
!! |1 land cover growing
!! laiday(:) |m**2/m**2 |leaf area index
!! ncrops(:,:,:)|
!! phuacc(:) |none |fraction of plant heat units accumulated
!! plantn(:) |kg N/ha |amount of nitrogen in plant biomass
!! plantp(:) |kg P/ha |amount of phosphorus in plant biomass
!! plt_pst(:,:)|kg/ha |pesticide on plant foliage
!! sol_fon(:,:)|kg N/ha |amount of nitrogen stored in the fresh
!! |organic (residue) pool
!! sol_fop(:,:)|kg P/ha |amount of phosphorus stored in the fresh
!! |organic (residue) pool
!! sol_pst(:,:,1)|kg/ha |pesticide in first layer of soil
!! sol_rsd(:,:)|kg/ha |amount of organic matter in the soil
!! |classified as residue
!! strsw(:) |none |fraction of potential plant growth achieved
!! |on the day where the reduction is caused by
!! |water stress
!! tnyld(:) |kg N/kg yield |modifier for autofertilization target
!! |nitrogen content for plant
!! wshd_yldn |kg N/ha |amount of nitrogen removed from soil in
!! |watershed in the yield
!! wshd_yldp |kg P/ha |amount of phosphorus removed from soil in
!! |watershed in the yield
!! yldanu(:) |metric tons/ha|annual yield (dry weight) in the HRU
!! yldkg(:,:,:)|kg/ha |yield (dry weight) by crop type in the HRU
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! ~ ~ ~ LOCAL DEFINITIONS ~ ~ ~
!! name |units |definition
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! hiad1 |
!! j |none |HRU number
!! k |none |counter
!! resnew |
!! wur |
!! yield |kg |yield (dry weight)
!! yieldn |
!! yieldp |
!! yldpst |kg pst/ha |pesticide removed in yield
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! ~ ~ ~ SUBROUTINES/FUNCTIONS CALLED ~ ~ ~
!! Intrinsic: Exp, Max, Min
!! SWAT: curno
!! ~ ~ ~ ~ ~ ~ END SPECIFICATIONS ~ ~ ~ ~ ~ ~
use parm
integer :: j, k
!! change per JGA 8/31/2011 gsm PUT YIELD IN modparm.f
!! real :: wur, hiad1, yield, yieldn, yieldp, yldpst
real :: wur, hiad1, yldpst
real :: resnew, rtresnew
!!By Zhang
!!=============
real :: BLG1, BLG2, BLG3, CLG, sf
real :: sol_min_n, resnew_n, resnew_ne
real :: LMF, LSF, LSLF, LSNF,LMNF
orgc_f = 0.
BLG1 = 0.
BLG2 = 0.
BLG3 = 0.
CLG = 0.
sf = 0.
sol_min_n = 0.
resnew = 0.
resnew_n = 0.
resnew_ne = 0.
LMF = 0.
LSF = 0.
LSLF = 0.
LSNF = 0.
LMNF = 0.
!!By Zhang
!!==================
j = 0
j = ihru
!! calculate modifier for autofertilization target nitrogen content
tnyld(j) = 0.
tnyld(j) = (1. - rwt(j)) * bio_ms(j) * pltfr_n(j) * auto_eff(j)
! if (icr(j) > 1) then
! tnyld(nro(j),icr(j)-1,j) = tnyld(nro(j),icr(j),j)
! else
! tnyld(nro(j),icr(j)+1,j) = tnyld(nro(j),icr(j),j)
! end if
hiad1 = 0.
if (hi_targ(j) > 0.) then
hiad1 = hi_targ(j)
else
if (plt_pet(j) < 10.) then
wur = 100.
else
wur = 0.
wur = 100. * plt_et(j) / plt_pet(j)
endif
hiad1 = (hvstiadj(j) - wsyf(idplt(j))) *
& (wur / (wur + Exp(6.13 - .0883 * wur))) +
& wsyf(idplt(j))
if (hiad1 > hvsti(idplt(j))) then
hiad1 = hvsti(idplt(j))
end if
end if
!! check if yield is from above or below ground
yield = 0.
resnew = 0.
rtresnew = 0.
!! stover fraction during harvkillop
xx = frac_harvk
if (xx .lt. 1.e-6) then
xx = hi_ovr
endif
!! stover fraction during harvkillop
if (hi_ovr > 1.e-6) then
yield = bio_ms(j) * hi_ovr
resnew = bio_ms(j) - yield
else
if (idc(idplt(j)) == 7) then
yield = bio_ms(j) * (1. - bio_leaf(idplt(j)))
resnew = bio_ms(j) - yield
else
if (hvsti(idplt(j)) > 1.001) then
yield = bio_ms(j) * (1. - 1. / (1. + hiad1))
resnew = bio_ms(j) / (1. + hiad1)
resnew = resnew * (1. - xx)
else
yield = (1. - rwt(j)) * bio_ms(j) * hiad1
resnew = (1. - rwt(j)) * (1. - hiad1) * bio_ms(j)
!! remove stover during harvkillop
resnew = resnew * (1. - xx)
rtresnew = rwt(j) * bio_ms(j)
endif
endif
end if
if (yield < 0.) yield = 0.
if (resnew < 0.) resnew = 0.
if (rtresnew < 0.) rtresnew = 0. ! Armen 19 May 2008
! I would avoid this check, it is
! safer to know if variable is negative
!!add by zhang
!!=================
!!use idplt(:,:,:) to calculate the crop type, then
!! decide which type of crop yield should be used.
if (cswat == 2) then
grainc_d(j) = grainc_d(j) + yield * 0.42
stoverc_d(j) = stoverc_d(j)+(bio_ms(j)-yield-rtresnew)*0.42*xx
rsdc_d(j) = rsdc_d(j) + resnew * 0.42
rsdc_d(j) = rsdc_d(j) + rtresnew * 0.42
end if
!!add by zhang
!!=================
!! calculate nutrients removed with yield
yieldn = 0.
yieldp = 0.
yieldn = yield * cnyld(idplt(j))
yieldp = yield * cpyld(idplt(j))
yieldn = Min(yieldn, 0.80 * plantn(j))
yieldp = Min(yieldp, 0.80 * plantp(j))
!! Armen 19 May 2008 / 21 January 2008
!! fraction of roots in each layer
call rootfr
!! fraction of N, P in residue (ff1) or roots (ff2)
ff1 = (1 - hiad1) / (1 - hiad1 + rwt(j))
ff2 = 1 - ff1
!! update residue, N, P on soil surface
sol_rsd(1,j) = resnew + sol_rsd(1,j)
sol_fon(1,j) = sol_fon(1,j) + FF1 * (plantn(j) - yieldn)
sol_fop(1,j) = sol_fop(1,j) + FF1 * (plantp(j) - yieldp)
sol_rsd(1,j) = Max(sol_rsd(1,j),0.)
sol_fon(1,j) = Max(sol_fon(1,j),0.)
sol_fop(1,j) = Max(sol_fop(1,j),0.)
!!insert new biomss by zhang
!!=================================
if (cswat == 2) then
!!all the lignin from STD is assigned to LSL,
!!add STDL calculation
!!
!sol_LSL(k,ihru) = sol_STDL(k,ihru)
!CLG=BLG(3,JJK)*HUI(JJK)/(HUI(JJK)+EXP(BLG(1,JJK)-BLG(2,JJK)*&HUI(JJK))
! 52 BLG1 = LIGNIN FRACTION IN PLANT AT .5 MATURITY
! 53 BLG2 = LIGNIN FRACTION IN PLANT AT MATURITY
!CROPCOM.dat BLG1 = 0.01 BLG2 = 0.10
!SUBROUTINE ASCRV(X1,X2,X3,X4)
!EPIC0810
!THIS SUBPROGRAM COMPUTES S CURVE PARMS GIVEN 2 (X,Y) POINTS.
!USE PARM
!XX=LOG(X3/X1-X3)
!X2=(XX-LOG(X4/X2-X4))/(X4-X3)
!X1=XX+X3*X2
!RETURN
!END
!HUI(JJK)=HU(JJK)/XPHU
BLG1 = 0.01/0.10
BLG2 = 0.99
BLG3 = 0.10
XX = log(0.5/BLG1-0.5)
BLG2 = (XX -log(1./BLG2-1.))/(1.-0.5)
BLG1 = XX + 0.5*BLG2
CLG=BLG3*phuacc(j)/(phuacc(j)+EXP(BLG1-BLG2*phuacc(j)))
!if (k == 1) then
sf = 0.05
!else
!sf = 0.1
!end if
!kg/ha
sol_min_n = 0.
sol_min_n = (sol_no3(1,j)+sol_nh3(1,j))
resnew = resnew
resnew_n = ff1 * (plantn(j) - yieldn)
resnew_ne = resnew_n + sf * sol_min_n
!Not sure 1000 should be here or not!
!RLN = 1000*(resnew * CLG/(resnew_n+1.E-5))
!RLN is the ratio of lignin to nitrogen in the newly added residue
RLN = (resnew * CLG/(resnew_n+1.E-5))
RLR = MIN(.8, resnew * CLG/(resnew+1.E-5))
LMF = 0.85 - 0.018 * RLN
if (LMF <0.01) then
LMF = 0.01
else
if (LMF >0.7) then
LMF = 0.7
end if
end if
!if ((resnew * CLG/(resnew_n+1.E-5)) < 47.22) then
! LMF = 0.85 - 0.018 * (resnew * CLG/(resnew_n+1.E-5))
!else
! LMF = 0.
!end if
LSF = 1 - LMF
sol_LM(1,j) = sol_LM(1,j) + LMF * resnew
sol_LS(1,j) = sol_LS(1,j) + LSF * resnew
!here a simplified assumption of 0.5 LSL
!LSLF = 0.0
!LSLF = CLG
sol_LSL(1,j) = sol_LSL(1,j) + RLR*resnew
sol_LSC(1,j) = sol_LSC(1,j) + 0.42*LSF * resnew
sol_LSLC(1,j) = sol_LSLC(1,j) + RLR*0.42*resnew
sol_LSLNC(1,j) = sol_LSC(1,j) - sol_LSLC(1,j)
!X3 = MIN(X6,0.42*LSF * resnew/150)
if (resnew_n >= (0.42 * LSF * resnew /150)) then
sol_LSN(1,j) = sol_LSN(1,j) + 0.42 * LSF * resnew / 150
sol_LMN(1,j) = sol_LMN(1,j) + resnew_n -
& (0.42 * LSF * resnew / 150) + 1.E-25
else
sol_LSN(1,j) = sol_LSN(1,j) + resnew_n
sol_LMN(1,j) = sol_LMN(1,j) + 1.E-25
end if
!LSNF = sol_LSN(1,j)/(sol_LS(1,j)+1.E-5)
sol_LMC(1,j) = sol_LMC(1,j) + 0.42 * LMF * resnew
!LMNF = sol_LMN(1,j)/(sol_LM(1,j) + 1.E-5)
!update no3 and nh3 in soil
sol_no3(1,j) = sol_no3(1,j) * (1-sf)
sol_nh3(1,j) = sol_nh3(1,j) * (1-sf)
end if
!!insert new biomss by zhang
!!===============================
!! allocate dead roots, N, P to soil layers
do l=1, sol_nly(j)
sol_rsd(l,j) = sol_rsd(l,j) + rtfr(l) *rtresnew
sol_fon(l,j) = sol_fon(l,j) + rtfr(l) *ff2 * (plantn(j) - yieldn)
sol_fop(l,j) = sol_fop(l,j) + rtfr(l) *ff2 * (plantp(j) - yieldp)
!!insert new biomss by zhang
!!==============================
if (cswat == 2) then
!!all the lignin from STD is assigned to LSL,
!!add STDL calculation
!!
!sol_LSL(k,ihru) = sol_STDL(k,ihru)
!CLG=BLG(3,JJK)*HUI(JJK)/(HUI(JJK)+EXP(BLG(1,JJK)-BLG(2,JJK)*&HUI(JJK))
! 52 BLG1 = LIGNIN FRACTION IN PLANT AT .5 MATURITY
! 53 BLG2 = LIGNIN FRACTION IN PLANT AT MATURITY
!CROPCOM.dat BLG1 = 0.01 BLG2 = 0.10
!SUBROUTINE ASCRV(X1,X2,X3,X4)
!EPIC0810
!THIS SUBPROGRAM COMPUTES S CURVE PARMS GIVEN 2 (X,Y) POINTS.
!USE PARM
!XX=LOG(X3/X1-X3)
!X2=(XX-LOG(X4/X2-X4))/(X4-X3)
!X1=XX+X3*X2
!RETURN
!END
!HUI(JJK)=HU(JJK)/XPHU
BLG1 = 0.01/0.10
BLG2 = 0.99
BLG3 = 0.10
XX = log(0.5/BLG1-0.5)
BLG2 = (XX -log(1./BLG2-1.))/(1.-0.5)
BLG1 = XX + 0.5*BLG2
CLG=BLG3*phuacc(j)/(phuacc(j)+EXP(BLG1-BLG2*phuacc(j)))
if (l == 1) then
sf = 0.05
else
sf = 0.1
end if
!kg/ha
sol_min_n = 0.
sol_min_n = (sol_no3(l,j)+sol_nh3(l,j))
resnew = rtfr(l) *rtresnew
resnew_n = rtfr(l) *ff2 * (plantn(j) - yieldn)
resnew_ne = resnew_n + sf * sol_min_n
!Not sure 1000 should be here or not!
!RLN = 1000*(resnew * CLG/(resnew_n+1.E-5))
RLN = (resnew * CLG/(resnew_n+1.E-5))
RLR = MIN(.8, resnew * CLG/1000/(resnew/1000+1.E-5))
LMF = 0.85 - 0.018 * RLN
if (LMF <0.01) then
LMF = 0.01
else
if (LMF >0.7) then
LMF = 0.7
end if
end if
!if ((resnew * CLG/(resnew_n+1.E-5)) < 47.22) then
! LMF = 0.85 - 0.018 * (resnew * CLG/(resnew_n+1.E-5))
!else
! LMF = 0.
!end if
LSF = 1 - LMF
sol_LM(l,j) = sol_LM(l,j) + LMF * resnew
sol_LS(l,j) = sol_LS(l,j) + LSF * resnew
!here a simplified assumption of 0.5 LSL
LSLF = 0.0
LSLF = CLG
sol_LSL(l,j) = sol_LSL(l,j) + RLR* LSF * resnew
sol_LSC(l,j) = sol_LSC(l,j) + 0.42*LSF * resnew
sol_LSLC(l,j) = sol_LSLC(l,j) + RLR*0.42*LSF * resnew
sol_LSLNC(l,j) = sol_LSC(l,j) - sol_LSLC(l,j)
!X3 = MIN(X6,0.42*LSF * resnew/150)
if (resnew_ne >= (0.42 * LSF * resnew /150)) then
sol_LSN(l,j) = sol_LSN(l,j) + 0.42 * LSF * resnew / 150
sol_LMN(l,j) = sol_LMN(l,j) + resnew_ne -
& (0.42 * LSF * resnew / 150) + 1.E-25
else
sol_LSN(l,j) = sol_LSN(l,j) + resnew_ne
sol_LMN(l,j) = sol_LMN(l,j) + 1.E-25
end if
!LSNF = sol_LSN(l,j)/(sol_LS(l,j)+1.E-5)
sol_LMC(l,j) = sol_LMC(l,j) + 0.42 * LMF * resnew
!LMNF = sol_LMN(l,j)/(sol_LM(l,j) + 1.E-5)
!update no3 and nh3 in soil
sol_no3(l,j) = sol_no3(l,j) * (1-sf)
sol_nh3(l,j) = sol_nh3(l,j) * (1-sf)
end if
!!insert new biomss by zhang
!!===============================
end do
!! adjust foliar pesticide for plant removal
if (hrupest(j) == 1) then
do k = 1, npmx
!! calculate amount of pesticide removed with yield
yldpst = 0.
if (hvsti(idplt(j)) > 1.001) then
yldpst = plt_pst(k,j)
plt_pst(k,j) = 0.
else
yldpst = hiad1 * plt_pst(k,j)
plt_pst(k,j) = plt_pst(k,j) - yldpst
if (plt_pst(k,j) < 0.) plt_pst(k,j) = 0.
endif
!! add pesticide in residue to soil surface
sol_pst(k,j,1) = sol_pst(k,j,1) + plt_pst(k,j)
plt_pst(k,j) = 0.
end do
end if
!! summary calculations
if (curyr > nyskip) then
wshd_yldn = wshd_yldn + yieldn * hru_dafr(j)
wshd_yldp = wshd_yldp + yieldp * hru_dafr(j)
yldkg(icr(j),j) = yldkg(icr(j),j) + yield
bio_hv(icr(j),j) = bio_ms(j) + bio_hv(icr(j),j)
yldanu(j) = yldanu(j) + yield / 1000.
bio_yrms(j) = bio_yrms(j) + bio_ms(j) / 1000.
ncrops(icr(j),j) = ncrops(icr(j),j) + 1
endif
!! update curve number
if (cnop > 0.)
* call curno(cnop,j)
!! reset variables
igro(j) = 0
idorm(j) = 0
bio_ms(j) = 0.
rwt(j) = 0.
plantn(j) = 0.
plantp(j) = 0.
strsw(j) = 1.
laiday(j) = 0.
hvstiadj(j) = 0.
rtfr = 0. ! resetting root fraction per layer array to 0
return
end