This repository has been archived by the owner on Sep 14, 2018. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 12
/
Copy pathEDPatchDynamicsMod.F90
1940 lines (1564 loc) · 91.8 KB
/
EDPatchDynamicsMod.F90
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
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
module EDPatchDynamicsMod
! ============================================================================
! Controls formation, creation, fusing and termination of patch level processes.
! ============================================================================
use FatesGlobals , only : fates_log
use FatesInterfaceMod , only : hlm_freq_day
use EDPftvarcon , only : EDPftvarcon_inst
use EDCohortDynamicsMod , only : fuse_cohorts, sort_cohorts, insert_cohort
use EDtypesMod , only : ncwd, n_dbh_bins, area, patchfusion_dbhbin_loweredges
use EDtypesMod , only : force_patchfuse_min_biomass
use EDTypesMod , only : maxPatchesPerSite
use EDTypesMod , only : ed_site_type, ed_patch_type, ed_cohort_type
use EDTypesMod , only : min_patch_area
use EDTypesMod , only : nclmax
use EDTypesMod , only : maxpft
use EDTypesMod , only : dtype_ifall
use EDTypesMod , only : dtype_ilog
use EDTypesMod , only : dtype_ifire
use FatesInterfaceMod , only : hlm_use_planthydro
use FatesInterfaceMod , only : hlm_numSWb
use FatesInterfaceMod , only : bc_in_type
use FatesInterfaceMod , only : hlm_days_per_year
use FatesInterfaceMod , only : numpft
use FatesGlobals , only : endrun => fates_endrun
use FatesConstantsMod , only : r8 => fates_r8
use FatesConstantsMod , only : itrue
use FatesPlantHydraulicsMod, only : InitHydrCohort
use FatesPlantHydraulicsMod, only : AccumulateMortalityWaterStorage
use FatesPlantHydraulicsMod, only : DeallocateHydrCohort
use EDLoggingMortalityMod, only : logging_litter_fluxes
use EDLoggingMortalityMod, only : logging_time
use EDParamsMod , only : fates_mortality_disturbance_fraction
use FatesAllometryMod , only : carea_allom
use FatesConstantsMod , only : g_per_kg
use FatesConstantsMod , only : ha_per_m2
use FatesConstantsMod , only : days_per_sec
use FatesConstantsMod , only : years_per_day
use FatesConstantsMod , only : nearzero
! CIME globals
use shr_infnan_mod , only : nan => shr_infnan_nan, assignment(=)
use shr_log_mod , only : errMsg => shr_log_errMsg
!
implicit none
private
!
public :: create_patch
public :: spawn_patches
public :: zero_patch
public :: fuse_patches
public :: terminate_patches
public :: patch_pft_size_profile
public :: disturbance_rates
public :: check_patch_area
public :: set_patchno
private:: fuse_2_patches
character(len=*), parameter, private :: sourcefile = &
__FILE__
logical, parameter :: debug = .false.
! 10/30/09: Created by Rosie Fisher
! ============================================================================
contains
! ============================================================================
subroutine disturbance_rates( site_in, bc_in)
!
! !DESCRIPTION:
! Calculates the fire and mortality related disturbance rates for each patch,
! and then determines which is the larger at the patch scale (for now, there an only
! be one disturbance type for each timestep.
! all disturbance rates here are per daily timestep.
! 2016-2017
! Modify to add logging disturbance
! !USES:
use EDMortalityFunctionsMod , only : mortality_rates
! loging flux
use EDLoggingMortalityMod , only : LoggingMortality_frac
! !ARGUMENTS:
type(ed_site_type) , intent(inout), target :: site_in
type(bc_in_type) , intent(in) :: bc_in
!
! !LOCAL VARIABLES:
type (ed_patch_type) , pointer :: currentPatch
type (ed_cohort_type), pointer :: currentCohort
real(r8) :: cmort
real(r8) :: bmort
real(r8) :: hmort
real(r8) :: frmort
real(r8) :: lmort_direct
real(r8) :: lmort_collateral
real(r8) :: lmort_infra
integer :: threshold_sizeclass
!----------------------------------------------------------------------------------------------
! Calculate Mortality Rates (these were previously calculated during growth derivatives)
! And the same rates in understory plants have already been applied to %dndt
!----------------------------------------------------------------------------------------------
currentPatch => site_in%oldest_patch
do while (associated(currentPatch))
currentCohort => currentPatch%shortest
do while(associated(currentCohort))
! Mortality for trees in the understorey.
currentCohort%patchptr => currentPatch
call mortality_rates(currentCohort,bc_in,cmort,hmort,bmort,frmort)
currentCohort%dmort = cmort+hmort+bmort+frmort
call carea_allom(currentCohort%dbh,currentCohort%n,site_in%spread,currentCohort%pft, &
currentCohort%c_area)
! Initialize diagnostic mortality rates
currentCohort%cmort = cmort
currentCohort%bmort = bmort
currentCohort%hmort = hmort
currentCohort%frmort = frmort
currentCohort%fmort = 0.0_r8 ! Fire mortality is initialized as zero, but may be changed
call LoggingMortality_frac(currentCohort%pft, currentCohort%dbh, &
lmort_direct,lmort_collateral,lmort_infra )
currentCohort%lmort_direct = lmort_direct
currentCohort%lmort_collateral = lmort_collateral
currentCohort%lmort_infra = lmort_infra
currentCohort => currentCohort%taller
end do
currentPatch => currentPatch%younger
end do
! ---------------------------------------------------------------------------------------------
! Calculate Disturbance Rates based on the mortality rates just calculated
! ---------------------------------------------------------------------------------------------
currentPatch => site_in%oldest_patch
do while (associated(currentPatch))
currentPatch%disturbance_rates(dtype_ifall) = 0.0_r8
currentPatch%disturbance_rates(dtype_ilog) = 0.0_r8
currentPatch%disturbance_rates(dtype_ifire) = 0.0_r8
currentCohort => currentPatch%shortest
do while(associated(currentCohort))
if(currentCohort%canopy_layer == 1)then
! Treefall Disturbance Rate
currentPatch%disturbance_rates(dtype_ifall) = currentPatch%disturbance_rates(dtype_ifall) + &
fates_mortality_disturbance_fraction * &
min(1.0_r8,currentCohort%dmort)*hlm_freq_day*currentCohort%c_area/currentPatch%area
! Logging Disturbance Rate
currentPatch%disturbance_rates(dtype_ilog) = currentPatch%disturbance_rates(dtype_ilog) + &
min(1.0_r8, currentCohort%lmort_direct + &
currentCohort%lmort_collateral + &
currentCohort%lmort_infra ) * &
currentCohort%c_area/currentPatch%area
endif
currentCohort => currentCohort%taller
enddo !currentCohort
! Fire Disturbance Rate
! Fudge - fires can't burn the whole patch, as this causes /0 errors.
! This is accumulating the daily fires over the whole 30 day patch generation phase.
currentPatch%disturbance_rates(dtype_ifire) = &
min(0.99_r8,currentPatch%disturbance_rates(dtype_ifire) + currentPatch%frac_burnt)
if (currentPatch%disturbance_rates(dtype_ifire) > 0.98_r8)then
write(fates_log(),*) 'very high fire areas', &
currentPatch%disturbance_rates(dtype_ifire),currentPatch%frac_burnt
endif
! ------------------------------------------------------------------------------------------
! Determine which disturbance is dominant, and force mortality diagnostics in the upper
! canopy to be zero for the non-dominant mode. Note: upper-canopy tree-fall mortality is
! not always disturbance generating, so when tree-fall mort is non-dominant, make sure
! to still diagnose and track the non-disturbance rate
! ------------------------------------------------------------------------------------------
if (currentPatch%disturbance_rates(dtype_ilog) > currentPatch%disturbance_rates(dtype_ifall) .and. &
currentPatch%disturbance_rates(dtype_ilog) > currentPatch%disturbance_rates(dtype_ifire) ) then
currentPatch%disturbance_rate = currentPatch%disturbance_rates(dtype_ilog)
! Update diagnostics
currentCohort => currentPatch%shortest
do while(associated(currentCohort))
if(currentCohort%canopy_layer == 1)then
currentCohort%fmort = 0.0_r8
currentCohort%cmort = currentCohort%cmort*(1.0_r8 - fates_mortality_disturbance_fraction)
currentCohort%hmort = currentCohort%hmort*(1.0_r8 - fates_mortality_disturbance_fraction)
currentCohort%bmort = currentCohort%bmort*(1.0_r8 - fates_mortality_disturbance_fraction)
currentCohort%dmort = currentCohort%dmort*(1.0_r8 - fates_mortality_disturbance_fraction)
currentCohort%frmort = currentCohort%frmort*(1.0_r8 - fates_mortality_disturbance_fraction)
end if
currentCohort => currentCohort%taller
enddo !currentCohort
elseif (currentPatch%disturbance_rates(dtype_ifire) > currentPatch%disturbance_rates(dtype_ifall) .and. &
currentPatch%disturbance_rates(dtype_ifire) > currentPatch%disturbance_rates(dtype_ilog) ) then ! DISTURBANCE IS FIRE
currentPatch%disturbance_rate = currentPatch%disturbance_rates(dtype_ifire)
! Update diagnostics, zero non-fire mortality rates
currentCohort => currentPatch%shortest
do while(associated(currentCohort))
if(currentCohort%canopy_layer == 1)then
currentCohort%cmort = currentCohort%cmort*(1.0_r8 - fates_mortality_disturbance_fraction)
currentCohort%hmort = currentCohort%hmort*(1.0_r8 - fates_mortality_disturbance_fraction)
currentCohort%bmort = currentCohort%bmort*(1.0_r8 - fates_mortality_disturbance_fraction)
currentCohort%dmort = currentCohort%dmort*(1.0_r8 - fates_mortality_disturbance_fraction)
currentCohort%frmort = currentCohort%frmort*(1.0_r8 - fates_mortality_disturbance_fraction)
currentCohort%lmort_direct = 0.0_r8
currentCohort%lmort_collateral = 0.0_r8
currentCohort%lmort_infra = 0.0_r8
end if
! This may be counter-intuitive, but the diagnostic fire-mortality rate
! will stay zero in the patch that undergoes fire, this is because
! the actual cohorts who experience the fire are only those in the
! newly created patch so currentCohort%fmort = 0.0_r8
! Don't worry, the cohorts in the newly created patch will reflect burn
currentCohort => currentCohort%taller
enddo !currentCohort
else ! If fire and loggin are not greater than treefall, just set disturbance rate to tree-fall
! which is most likely a 0.0
currentPatch%disturbance_rate = currentPatch%disturbance_rates(dtype_ifall)
! Update diagnostics, zero non-treefall mortality rates
currentCohort => currentPatch%shortest
do while(associated(currentCohort))
if(currentCohort%canopy_layer == 1)then
currentCohort%lmort_direct = 0.0_r8
currentCohort%lmort_collateral = 0.0_r8
currentCohort%lmort_infra = 0.0_r8
currentCohort%fmort = 0.0_r8
end if
currentCohort => currentCohort%taller
enddo !currentCohort
endif
currentPatch => currentPatch%younger
enddo !patch loop
end subroutine disturbance_rates
! ============================================================================
subroutine spawn_patches( currentSite, bc_in)
!
! !DESCRIPTION:
! In this subroutine, the following happens
! 1) the total area disturbed is calculated
! 2) a new patch is created
! 3) properties are averaged
! 4) litter fluxes from fire and mortality are added
! 5) For mortality, plants in existing patch canopy are killed.
! 6) For mortality, Plants in new and existing understorey are killed
! 7) For fire, burned plants are killed, and unburned plants are added to new patch.
! 8) New cohorts are added to new patch and sorted.
! 9) New patch is added into linked list
! 10) Area checked, and patchno recalculated.
!
! !USES:
use EDParamsMod , only : ED_val_understorey_death, logging_coll_under_frac
use EDCohortDynamicsMod , only : zero_cohort, copy_cohort, terminate_cohorts
!
! !ARGUMENTS:
type (ed_site_type), intent(inout), target :: currentSite
type (bc_in_type), intent(in) :: bc_in
!
! !LOCAL VARIABLES:
type (ed_patch_type) , pointer :: new_patch
type (ed_patch_type) , pointer :: currentPatch
type (ed_cohort_type), pointer :: currentCohort
type (ed_cohort_type), pointer :: nc
type (ed_cohort_type), pointer :: storesmallcohort
type (ed_cohort_type), pointer :: storebigcohort
real(r8) :: site_areadis ! total area disturbed in m2 per site per day
real(r8) :: patch_site_areadis ! total area disturbed in m2 per patch per day
real(r8) :: age ! notional age of this patch in years
integer :: tnull ! is there a tallest cohort?
integer :: snull ! is there a shortest cohort?
real(r8) :: root_litter_local(maxpft) ! initial value of root litter. KgC/m2
real(r8) :: leaf_litter_local(maxpft) ! initial value of leaf litter. KgC/m2
real(r8) :: cwd_ag_local(ncwd) ! initial value of above ground coarse woody debris. KgC/m2
real(r8) :: cwd_bg_local(ncwd) ! initial value of below ground coarse woody debris. KgC/m2
!---------------------------------------------------------------------
storesmallcohort => null() ! storage of the smallest cohort for insertion routine
storebigcohort => null() ! storage of the largest cohort for insertion routine
! calculate area of disturbed land, in this timestep, by summing contributions from each existing patch.
currentPatch => currentSite%youngest_patch
! zero site-level fire fluxes
currentSite%cwd_ag_burned = 0.0_r8
currentSite%leaf_litter_burned = 0.0_r8
currentSite%total_burn_flux_to_atm = 0.0_r8
site_areadis = 0.0_r8
do while(associated(currentPatch))
!FIX(RF,032414) Does using the max(fire,mort) actually make sense here?
if(currentPatch%disturbance_rate>1.0_r8) then
write(fates_log(),*) 'patch disturbance rate > 1 ?',currentPatch%disturbance_rate
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
site_areadis = site_areadis + currentPatch%area * currentPatch%disturbance_rate
currentPatch => currentPatch%older
enddo ! end loop over patches. sum area disturbed for all patches.
if (site_areadis > 0.0_r8) then
cwd_ag_local = 0.0_r8
cwd_bg_local = 0.0_r8
leaf_litter_local = 0.0_r8
root_litter_local = 0.0_r8
age = 0.0_r8
allocate(new_patch)
call create_patch(currentSite, new_patch, age, site_areadis, &
cwd_ag_local, cwd_bg_local, leaf_litter_local, &
root_litter_local, bc_in%nlevsoil)
new_patch%tallest => null()
new_patch%shortest => null()
currentPatch => currentSite%oldest_patch
! loop round all the patches that contribute surviving indivduals and litter pools to the new patch.
do while(associated(currentPatch))
! This is the amount of patch area that is disturbed, and donated by the donor
patch_site_areadis = currentPatch%area * currentPatch%disturbance_rate
call average_patch_properties(currentPatch, new_patch, patch_site_areadis)
if (currentPatch%disturbance_rates(dtype_ilog) > currentPatch%disturbance_rates(dtype_ifall) .and. &
currentPatch%disturbance_rates(dtype_ilog) > currentPatch%disturbance_rates(dtype_ifire) ) then
call logging_litter_fluxes(currentSite, currentPatch, new_patch, patch_site_areadis)
elseif (currentPatch%disturbance_rates(dtype_ifire) > currentPatch%disturbance_rates(dtype_ifall) .and. &
currentPatch%disturbance_rates(dtype_ifire) > currentPatch%disturbance_rates(dtype_ilog) ) then
call fire_litter_fluxes(currentSite, currentPatch, new_patch, patch_site_areadis)
else
call mortality_litter_fluxes(currentSite, currentPatch, new_patch, patch_site_areadis)
endif
!INSERT SURVIVORS FROM DISTURBANCE INTO NEW PATCH
currentCohort => currentPatch%shortest
do while(associated(currentCohort))
allocate(nc)
if(hlm_use_planthydro.eq.itrue) call InitHydrCohort(CurrentSite,nc)
call zero_cohort(nc)
! nc is the new cohort that goes in the disturbed patch (new_patch)... currentCohort
! is the curent cohort that stays in the donor patch (currentPatch)
call copy_cohort(currentCohort, nc)
!this is the case as the new patch probably doesn't have a closed canopy, and
! even if it does, that will be sorted out in canopy_structure.
nc%canopy_layer = 1
nc%canopy_layer_yesterday = 1._r8
! treefall mortality is the dominant disturbance
if(currentPatch%disturbance_rates(dtype_ifall) > currentPatch%disturbance_rates(dtype_ifire) .and. &
currentPatch%disturbance_rates(dtype_ifall) > currentPatch%disturbance_rates(dtype_ilog))then
if(currentCohort%canopy_layer == 1)then
! In the donor patch we are left with fewer trees because the area has decreased
! the plant density for large trees does not actually decrease in the donor patch
! because this is the part of the original patch where no trees have actually fallen
! The diagnostic cmort,bmort,hmort, and frmort rates have already been saved
currentCohort%n = currentCohort%n * (1.0_r8 - fates_mortality_disturbance_fraction * &
min(1.0_r8,currentCohort%dmort * hlm_freq_day))
nc%n = 0.0_r8 ! kill all of the trees who caused the disturbance.
nc%cmort = nan ! The mortality diagnostics are set to nan because the cohort should dissappear
nc%hmort = nan
nc%bmort = nan
nc%fmort = nan
nc%frmort = nan
nc%lmort_direct = nan
nc%lmort_collateral = nan
nc%lmort_infra = nan
else
! small trees
if(EDPftvarcon_inst%woody(currentCohort%pft) == 1)then
! Survivorship of undestory woody plants. Two step process.
! Step 1: Reduce current number of plants to reflect the change in area.
! The number density per square are doesn't change, but since the patch is smaller
! and cohort counts are absolute, reduce this number.
nc%n = currentCohort%n * patch_site_areadis/currentPatch%area
! because the mortality rate due to impact for the cohorts which had been in the understory and are now in the newly-
! disturbed patch is very high, passing the imort directly to history results in large numerical errors, on account
! of the sharply reduced number densities. so instead pass this info via a site-level diagnostic variable before reducing
! the number density.
currentSite%imort_rate(currentCohort%size_class, currentCohort%pft) = &
currentSite%imort_rate(currentCohort%size_class, currentCohort%pft) + &
nc%n * ED_val_understorey_death / hlm_freq_day
currentSite%imort_carbonflux = currentSite%imort_carbonflux + &
(nc%n * ED_val_understorey_death / hlm_freq_day ) * &
currentCohort%b_total() * g_per_kg * days_per_sec * years_per_day * ha_per_m2
! Step 2: Apply survivor ship function based on the understory death fraction
! remaining of understory plants of those that are knocked over by the overstorey trees dying...
nc%n = nc%n * (1.0_r8 - ED_val_understorey_death)
! since the donor patch split and sent a fraction of its members
! to the new patch and a fraction to be preserved in itself,
! when reporting diagnostic rates, we must carry over the mortality rates from
! the donor that were applied before the patch split. Remember this is only
! for diagnostics. But think of it this way, the rates are weighted by
! number density in EDCLMLink, and the number density of this new patch is donated
! so with the number density must come the effective mortality rates.
nc%fmort = 0.0_r8 ! Should had also been zero in the donor
nc%cmort = currentCohort%cmort
nc%hmort = currentCohort%hmort
nc%bmort = currentCohort%bmort
nc%frmort = currentCohort%frmort
nc%dmort = currentCohort%dmort
nc%lmort_direct = currentCohort%lmort_direct
nc%lmort_collateral = currentCohort%lmort_collateral
nc%lmort_infra = currentCohort%lmort_infra
! understory trees that might potentially be knocked over in the disturbance.
! The existing (donor) patch should not have any impact mortality, it should
! only lose cohorts due to the decrease in area. This is not mortality.
! Besides, the current and newly created patch sum to unity
currentCohort%n = currentCohort%n * (1._r8 - patch_site_areadis/currentPatch%area)
else
! grass is not killed by mortality disturbance events. Just move it into the new patch area.
! Just split the grass into the existing and new patch structures
nc%n = currentCohort%n * patch_site_areadis/currentPatch%area
! Those remaining in the existing
currentCohort%n = currentCohort%n * (1._r8 - patch_site_areadis/currentPatch%area)
nc%fmort = 0.0_r8
nc%cmort = currentCohort%cmort
nc%hmort = currentCohort%hmort
nc%bmort = currentCohort%bmort
nc%frmort = currentCohort%frmort
nc%dmort = currentCohort%dmort
nc%lmort_direct = currentCohort%lmort_direct
nc%lmort_collateral = currentCohort%lmort_collateral
nc%lmort_infra = currentCohort%lmort_infra
endif
endif
! Fire is the dominant disturbance
elseif (currentPatch%disturbance_rates(dtype_ifire) > currentPatch%disturbance_rates(dtype_ifall) .and. &
currentPatch%disturbance_rates(dtype_ifire) > currentPatch%disturbance_rates(dtype_ilog)) then !fire
! Number of members in the new patch, before we impose fire survivorship
nc%n = currentCohort%n * patch_site_areadis/currentPatch%area
! loss of individuals from source patch due to area shrinking
currentCohort%n = currentCohort%n * (1._r8 - patch_site_areadis/currentPatch%area)
! loss of individual from fire in new patch.
nc%n = nc%n * (1.0_r8 - currentCohort%fire_mort)
nc%fmort = currentCohort%fire_mort/hlm_freq_day
nc%cmort = currentCohort%cmort
nc%hmort = currentCohort%hmort
nc%bmort = currentCohort%bmort
nc%frmort = currentCohort%frmort
nc%dmort = currentCohort%dmort
nc%lmort_direct = currentCohort%lmort_direct
nc%lmort_collateral = currentCohort%lmort_collateral
nc%lmort_infra = currentCohort%lmort_infra
! Logging is the dominant disturbance
elseif (currentPatch%disturbance_rates(dtype_ilog) > currentPatch%disturbance_rates(dtype_ifall) .and. &
currentPatch%disturbance_rates(dtype_ilog) > currentPatch%disturbance_rates(dtype_ifire)) then ! Logging
! If this cohort is in the upper canopy. It generated
if(currentCohort%canopy_layer == 1)then
! Trees generating this disturbance are not there by definition
nc%n = 0.0_r8
! Reduce counts in the existing/donor patch according to the logging rate
currentCohort%n = currentCohort%n * (1.0_r8 - min(1.0_r8,(currentCohort%lmort_direct + &
currentCohort%lmort_collateral + &
currentCohort%lmort_infra)))
! The mortality diagnostics are set to nan because the cohort should dissappear
nc%cmort = nan
nc%hmort = nan
nc%bmort = nan
nc%frmort = nan
nc%fmort = nan
nc%lmort_direct = nan
nc%lmort_collateral = nan
nc%lmort_infra = nan
else
! WHat to do with cohorts in the understory of a logging generated
! disturbance patch?
if(EDPftvarcon_inst%woody(currentCohort%pft) == 1)then
! Survivorship of undestory woody plants. Two step process.
! Step 1: Reduce current number of plants to reflect the change in area.
! The number density per square are doesn't change, but since the patch is smaller
! and cohort counts are absolute, reduce this number.
nc%n = currentCohort%n * patch_site_areadis/currentPatch%area
! because the mortality rate due to impact for the cohorts which had been in the understory and are now in the newly-
! disturbed patch is very high, passing the imort directly to history results in large numerical errors, on account
! of the sharply reduced number densities. so instead pass this info via a site-level diagnostic variable before reducing
! the number density.
currentSite%imort_rate(currentCohort%size_class, currentCohort%pft) = &
currentSite%imort_rate(currentCohort%size_class, currentCohort%pft) + &
nc%n * logging_coll_under_frac / hlm_freq_day
currentSite%imort_carbonflux = currentSite%imort_carbonflux + &
(nc%n * logging_coll_under_frac/ hlm_freq_day ) * &
currentCohort%b_total() * g_per_kg * days_per_sec * years_per_day * ha_per_m2
! Step 2: Apply survivor ship function based on the understory death fraction
! remaining of understory plants of those that are knocked over by the overstorey trees dying...
! LOGGING SURVIVORSHIP OF UNDERSTORY PLANTS IS SET AS A NEW PARAMETER in the fatesparameter files
nc%n = nc%n * (1.0_r8 - logging_coll_under_frac)
! Step 3: Reduce the number count of cohorts in the original/donor/non-disturbed patch
! to reflect the area change
currentCohort%n = currentCohort%n * (1._r8 - patch_site_areadis/currentPatch%area)
nc%fmort = 0.0_r8
nc%cmort = currentCohort%cmort
nc%hmort = currentCohort%hmort
nc%bmort = currentCohort%bmort
nc%frmort = currentCohort%frmort
nc%dmort = currentCohort%dmort
nc%lmort_direct = currentCohort%lmort_direct
nc%lmort_collateral = currentCohort%lmort_collateral
nc%lmort_infra = currentCohort%lmort_infra
else
! grass is not killed by mortality disturbance events. Just move it into the new patch area.
! Just split the grass into the existing and new patch structures
nc%n = currentCohort%n * patch_site_areadis/currentPatch%area
! Those remaining in the existing
currentCohort%n = currentCohort%n * (1._r8 - patch_site_areadis/currentPatch%area)
! No grass impact mortality imposed on the newly created patch
nc%fmort = 0.0_r8
nc%cmort = currentCohort%cmort
nc%hmort = currentCohort%hmort
nc%bmort = currentCohort%bmort
nc%frmort = currentCohort%frmort
nc%dmort = currentCohort%dmort
nc%lmort_direct = currentCohort%lmort_direct
nc%lmort_collateral = currentCohort%lmort_collateral
nc%lmort_infra = currentCohort%lmort_infra
endif ! is/is-not woody
endif ! Select canopy layer
end if ! Select disturbance mode
if (nc%n > 0.0_r8) then
storebigcohort => new_patch%tallest
storesmallcohort => new_patch%shortest
if(associated(new_patch%tallest))then
tnull = 0
else
tnull = 1
new_patch%tallest => nc
nc%taller => null()
endif
if(associated(new_patch%shortest))then
snull = 0
else
snull = 1
new_patch%shortest => nc
nc%shorter => null()
endif
nc%patchptr => new_patch
call insert_cohort(nc, new_patch%tallest, new_patch%shortest, tnull, snull, storebigcohort, storesmallcohort)
new_patch%tallest => storebigcohort
new_patch%shortest => storesmallcohort
else
if(hlm_use_planthydro.eq.itrue) call DeallocateHydrCohort(nc)
deallocate(nc) !get rid of the new memory.
endif
currentCohort => currentCohort%taller
enddo ! currentCohort
call sort_cohorts(currentPatch)
!zero disturbance accumulators
currentPatch%disturbance_rate = 0._r8
currentPatch%disturbance_rates = 0._r8
!update area of donor patch
currentPatch%area = currentPatch%area - patch_site_areadis
! sort out the cohorts, since some of them may be so small as to need removing.
! the first call to terminate cohorts removes sparse number densities,
! the second call removes for all other reasons (sparse culling must happen
! before fusion)
call terminate_cohorts(currentSite, currentPatch, 1)
call fuse_cohorts(currentSite,currentPatch, bc_in)
call terminate_cohorts(currentSite, currentPatch, 2)
call sort_cohorts(currentPatch)
currentPatch => currentPatch%younger
enddo ! currentPatch patch loop.
!*************************/
!** INSERT NEW PATCH INTO LINKED LIST
!**********`***************/
currentPatch => currentSite%youngest_patch
new_patch%older => currentPatch
new_patch%younger => NULL()
currentPatch%younger => new_patch
currentSite%youngest_patch => new_patch
! sort out the cohorts, since some of them may be so small as to need removing.
! the first call to terminate cohorts removes sparse number densities,
! the second call removes for all other reasons (sparse culling must happen
! before fusion)
call terminate_cohorts(currentSite, new_patch, 1)
call fuse_cohorts(currentSite,new_patch, bc_in)
call terminate_cohorts(currentSite, new_patch, 2)
call sort_cohorts(new_patch)
endif !end new_patch area
call check_patch_area(currentSite)
call set_patchno(currentSite)
end subroutine spawn_patches
! ============================================================================
subroutine check_patch_area( currentSite )
!
! !DESCRIPTION:
! Check to see that total area is not exceeded.
!
! !USES:
!
! !ARGUMENTS:
type(ed_site_type), intent(in), target :: currentSite
!
! !LOCAL VARIABLES:
real(r8) :: areatot
type(ed_patch_type), pointer :: currentPatch
type(ed_patch_type), pointer :: largestPatch
real(r8) :: largest_area
real(r8), parameter :: area_error_fail = 1.0e-6_r8
!---------------------------------------------------------------------
areatot = 0._r8
largest_area = 0._r8
largestPatch => null()
currentPatch => currentSite%oldest_patch
do while(associated(currentPatch))
areatot = areatot + currentPatch%area
if(currentPatch%area>largest_area) then
largestPatch => currentPatch
largest_area = currentPatch%area
end if
currentPatch => currentPatch%younger
end do
if ( abs( areatot - area ) > nearzero ) then
if ( abs(areatot-area) > area_error_fail ) then
write(fates_log(),*) 'Patch areas do not sum to 10000 within tolerance'
write(fates_log(),*) 'Total area: ',areatot,'absolute error: ',areatot-area
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
if(debug) then
write(fates_log(),*) 'Total patch area precision being fixed, adjusting'
write(fates_log(),*) 'largest patch. This may have slight impacts on carbon balance.'
end if
largestPatch%area = largestPatch%area + (area-areatot)
endif
return
end subroutine check_patch_area
! ============================================================================
subroutine set_patchno( currentSite )
!
! !DESCRIPTION:
! Give patches an order number from the oldest to youngest.
!
! !USES:
!
! !ARGUMENTS:
type(ed_site_type),intent(in), target :: currentSite
!
! !LOCAL VARIABLES:
type(ed_patch_type), pointer :: currentPatch
integer patchno
!---------------------------------------------------------------------
patchno = 1
currentPatch => currentSite%oldest_patch
do while(associated(currentPatch))
currentPatch%patchno = patchno
patchno = patchno + 1
currentPatch => currentPatch%younger
enddo
end subroutine set_patchno
! ============================================================================
subroutine average_patch_properties( currentPatch, newPatch, patch_site_areadis )
!
! !DESCRIPTION:
! Average together the state properties of all of the donor patches that
! make up the new patch.
!
! !USES:
!
! !ARGUMENTS:
type(ed_patch_type) , intent(in), target :: currentPatch
type(ed_patch_type) , intent(inout) :: newPatch
real(r8) , intent(out) :: patch_site_areadis ! amount of land disturbed in this patch. m2
!
! !LOCAL VARIABLES:
integer :: c,p ! counters for PFT and litter size class.
!---------------------------------------------------------------------
patch_site_areadis = currentPatch%area * currentPatch%disturbance_rate ! how much land is disturbed in this donor patch?
do c = 1,ncwd !move litter pool en mass into the new patch.
newPatch%cwd_ag(c) = newPatch%cwd_ag(c) + currentPatch%cwd_ag(c) * patch_site_areadis/newPatch%area
newPatch%cwd_bg(c) = newPatch%cwd_bg(c) + currentPatch%cwd_bg(c) * patch_site_areadis/newPatch%area
enddo
do p = 1,numpft !move litter pool en mass into the new patch
newPatch%root_litter(p) = newPatch%root_litter(p) + &
currentPatch%root_litter(p) * patch_site_areadis/newPatch%area
newPatch%leaf_litter(p) = newPatch%leaf_litter(p) + &
currentPatch%leaf_litter(p) * patch_site_areadis/newPatch%area
! The fragmentation/decomposition flux from donor patches has already occured in existing patches. However
! some of their area has been carved out for this new patches which is receiving donations.
! Lets maintain conservation on that pre-existing mass flux in these newly disturbed patches
newPatch%root_litter_out(p) = newPatch%root_litter_out(p) + &
currentPatch%root_litter_out(p) * patch_site_areadis/newPatch%area
newPatch%leaf_litter_out(p) = newPatch%leaf_litter_out(p) + &
currentPatch%leaf_litter_out(p) * patch_site_areadis/newPatch%area
enddo
end subroutine average_patch_properties
! ============================================================================
subroutine fire_litter_fluxes(currentSite, cp_target, new_patch_target, patch_site_areadis)
!
! !DESCRIPTION:
! CWD pool burned by a fire.
! Carbon going from burned trees into CWD pool
! Burn parts of trees that don't die in fire
! Burn live grasses and kill them.
!
! !USES:
use SFParamsMod, only : SF_VAL_CWD_FRAC
use EDtypesMod , only : dl_sf
!
! !ARGUMENTS:
type(ed_site_type) , intent(inout), target :: currentSite
type(ed_patch_type) , intent(inout), target :: cp_target
type(ed_patch_type) , intent(inout), target :: new_patch_target
real(r8) , intent(inout) :: patch_site_areadis
!
! !LOCAL VARIABLES:
type(ed_patch_type) , pointer :: currentPatch
type(ed_patch_type) , pointer :: new_patch
type(ed_cohort_type), pointer :: currentCohort
real(r8) :: bcroot ! amount of below ground coarse root per cohort kgC. (goes into CWD_BG)
real(r8) :: bstem ! amount of above ground stem biomass per cohort kgC.(goes into CWG_AG)
real(r8) :: dead_tree_density ! no trees killed by fire per m2
reaL(r8) :: burned_litter ! amount of each litter pool burned by fire. kgC/m2/day
real(r8) :: burned_leaves ! amount of tissue consumed by fire for grass. KgC/individual/day
integer :: c, p
!---------------------------------------------------------------------
!check that total area is not exceeded.
currentPatch => cp_target
new_patch => new_patch_target
if ( currentPatch%fire == 1 ) then !only do this if there was a fire in this actual patch.
patch_site_areadis = currentPatch%area * currentPatch%disturbance_rate ! how much land is disturbed in this donor patch?
!************************************/
!PART 1) Burn the fractions of existing litter in the new patch that were consumed by the fire.
!************************************/
do c = 1,ncwd
burned_litter = new_patch%cwd_ag(c) * patch_site_areadis/new_patch%area * &
currentPatch%burnt_frac_litter(c+1) !kG/m2/day
new_patch%cwd_ag(c) = new_patch%cwd_ag(c) - burned_litter
currentSite%flux_out = currentSite%flux_out + burned_litter * new_patch%area !kG/site/day
currentSite%total_burn_flux_to_atm = currentSite%total_burn_flux_to_atm + &
burned_litter * new_patch%area !kG/site/day
enddo
do p = 1,numpft
burned_litter = new_patch%leaf_litter(p) * patch_site_areadis/new_patch%area * &
currentPatch%burnt_frac_litter(dl_sf)
new_patch%leaf_litter(p) = new_patch%leaf_litter(p) - burned_litter
currentSite%flux_out = currentSite%flux_out + burned_litter * new_patch%area !kG/site/day
currentSite%total_burn_flux_to_atm = currentSite%total_burn_flux_to_atm + &
burned_litter * new_patch%area !kG/site/day
enddo
!************************************/
!PART 2) Put unburned parts of plants that died in the fire into the litter pool of new and old patches
! This happens BEFORE the plant numbers have been updated. So we are working with the
! pre-fire population of plants, which is the right way round.
!************************************/
currentCohort => currentPatch%shortest
do while(associated(currentCohort))
p = currentCohort%pft
if(EDPftvarcon_inst%woody(p) == 1)then !DEAD (FROM FIRE) TREES
!************************************/
! Number of trees that died because of the fire, per m2 of ground.
! Divide their litter into the four litter streams, and spread evenly across ground surface.
!************************************/
! stem biomass per tree
bstem = (currentCohort%bsw + currentCohort%bdead) * EDPftvarcon_inst%allom_agb_frac(p)
! coarse root biomass per tree
bcroot = (currentCohort%bsw + currentCohort%bdead) * (1.0_r8 - EDPftvarcon_inst%allom_agb_frac(p) )
! density of dead trees per m2.
dead_tree_density = (currentCohort%fire_mort * currentCohort%n*patch_site_areadis/currentPatch%area) / AREA
if( hlm_use_planthydro == itrue ) then
call AccumulateMortalityWaterStorage(currentSite,currentCohort,dead_tree_density*AREA)
end if
! Unburned parts of dead tree pool.
! Unburned leaves and roots
new_patch%leaf_litter(p) = new_patch%leaf_litter(p) + dead_tree_density * (currentCohort%bl) &
* (1.0_r8-currentCohort%cfa)
new_patch%root_litter(p) = new_patch%root_litter(p) + dead_tree_density * (currentCohort%br+currentCohort%bstore)
currentPatch%leaf_litter(p) = currentPatch%leaf_litter(p) + dead_tree_density * &
(currentCohort%bl) * (1.0_r8-currentCohort%cfa)
currentPatch%root_litter(p) = currentPatch%root_litter(p) + dead_tree_density * &
(currentCohort%br+currentCohort%bstore)
! track as diagnostic fluxes
currentSite%leaf_litter_diagnostic_input_carbonflux(p) = currentSite%leaf_litter_diagnostic_input_carbonflux(p) + &
(currentCohort%bl) * (1.0_r8-currentCohort%cfa) * currentCohort%fire_mort * currentCohort%n * &
hlm_days_per_year / AREA
currentSite%root_litter_diagnostic_input_carbonflux(p) = currentSite%root_litter_diagnostic_input_carbonflux(p) + &
(currentCohort%br+currentCohort%bstore) * (1.0_r8-currentCohort%cfa) * currentCohort%fire_mort * &
currentCohort%n * hlm_days_per_year / AREA
! below ground coarse woody debris from burned trees
do c = 1,ncwd
new_patch%cwd_bg(c) = new_patch%cwd_bg(c) + dead_tree_density * SF_val_CWD_frac(c) * bcroot
currentPatch%cwd_bg(c) = currentPatch%cwd_bg(c) + dead_tree_density * SF_val_CWD_frac(c) * bcroot
! track as diagnostic fluxes
currentSite%CWD_BG_diagnostic_input_carbonflux(c) = currentSite%CWD_BG_diagnostic_input_carbonflux(c) + &
SF_val_CWD_frac(c) * bcroot * currentCohort%fire_mort * currentCohort%n * &
hlm_days_per_year / AREA
enddo
! above ground coarse woody debris from unburned twigs and small branches
do c = 1,2
new_patch%cwd_ag(c) = new_patch%cwd_ag(c) + dead_tree_density * SF_val_CWD_frac(c) * bstem &
* (1.0_r8-currentCohort%cfa)
currentPatch%cwd_ag(c) = currentPatch%cwd_ag(c) + dead_tree_density * SF_val_CWD_frac(c) * &
bstem * (1.0_r8-currentCohort%cfa)
! track as diagnostic fluxes
currentSite%CWD_AG_diagnostic_input_carbonflux(c) = currentSite%CWD_AG_diagnostic_input_carbonflux(c) + &
SF_val_CWD_frac(c) * bstem * (1.0_r8-currentCohort%cfa) * currentCohort%fire_mort * currentCohort%n * &
hlm_days_per_year / AREA
enddo
! above ground coarse woody debris from large branches and stems: these do not burn in crown fires.
do c = 3,4
new_patch%cwd_ag(c) = new_patch%cwd_ag(c) + dead_tree_density * SF_val_CWD_frac(c) * bstem
currentPatch%cwd_ag(c) = currentPatch%cwd_ag(c) + dead_tree_density * SF_val_CWD_frac(c) * bstem
! track as diagnostic fluxes
currentSite%CWD_AG_diagnostic_input_carbonflux(c) = currentSite%CWD_AG_diagnostic_input_carbonflux(c) + &
SF_val_CWD_frac(c) * bstem * currentCohort%fire_mort * currentCohort%n * &
hlm_days_per_year / AREA
enddo
! Burned parts of dead tree pool.
! Burned twigs and small branches.
do c = 1,2
currentSite%cwd_ag_burned(c) = currentSite%cwd_ag_burned(c) + dead_tree_density * &
SF_val_CWD_frac(c) * bstem * currentCohort%cfa
currentSite%flux_out = currentSite%flux_out + dead_tree_density * &
AREA * SF_val_CWD_frac(c) * bstem * currentCohort%cfa
currentSite%total_burn_flux_to_atm = currentSite%total_burn_flux_to_atm + dead_tree_density * &
AREA * SF_val_CWD_frac(c) * bstem * currentCohort%cfa
enddo
!burned leaves.
do p = 1,numpft
currentSite%leaf_litter_burned(p) = currentSite%leaf_litter_burned(p) + &
dead_tree_density * currentCohort%bl * currentCohort%cfa
currentSite%flux_out = currentSite%flux_out + &
dead_tree_density * AREA * currentCohort%bl * currentCohort%cfa
currentSite%total_burn_flux_to_atm = currentSite%total_burn_flux_to_atm + &
dead_tree_density * AREA * currentCohort%bl * currentCohort%cfa
enddo
endif
currentCohort => currentCohort%taller
enddo ! currentCohort
!************************************/
! PART 3) Burn parts of trees that did *not* die in the fire.
! PART 4) Burn parts of grass that are consumed by the fire.
! grasses are not killed directly by fire. They die by losing all of their leaves and starving.
!************************************/
currentCohort => new_patch%shortest
do while(associated(currentCohort))
call carea_allom(currentCohort%dbh,currentCohort%n,currentSite%spread,currentCohort%pft,currentCohort%c_area)
if(EDPftvarcon_inst%woody(currentCohort%pft) == 1)then
burned_leaves = min(currentCohort%bl, (currentCohort%bl+currentCohort%bsw) * currentCohort%cfa)
else
burned_leaves = min(currentCohort%bl, (currentCohort%bl+currentCohort%bsw) * currentPatch%burnt_frac_litter(6))