-
Notifications
You must be signed in to change notification settings - Fork 92
/
Copy pathEDCohortDynamicsMod.F90
1857 lines (1524 loc) · 88 KB
/
EDCohortDynamicsMod.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 EDCohortDynamicsMod
!
! DESCRIPTION:
! Cohort stuctures in FATES
!
! USES:
use FatesGlobals , only : endrun => fates_endrun
use FatesGlobals , only : fates_log
use FatesInterfaceTypesMod , only : hlm_freq_day
use FatesInterfaceTypesMod , only : bc_in_type
use FatesInterfaceTypesMod , only : hlm_use_planthydro
use FatesInterfaceTypesMod , only : hlm_use_sp
use FatesInterfaceTypesMod , only : hlm_use_cohort_age_tracking
use FatesInterfaceTypesMod , only : hlm_use_tree_damage
use FatesInterfaceTypesMod , only : hlm_is_restart
use FatesConstantsMod , only : r8 => fates_r8
use FatesConstantsMod , only : fates_unset_int
use FatesConstantsMod , only : itrue,ifalse
use FatesConstantsMod , only : fates_unset_r8
use FatesConstantsMod , only : nearzero
use FatesConstantsMod , only : calloc_abs_error
use FatesInterfaceTypesMod , only : hlm_days_per_year
use FatesInterfaceTypesMod , only : nleafage
use SFParamsMod , only : SF_val_CWD_frac
use EDPftvarcon , only : EDPftvarcon_inst
use EDPftvarcon , only : GetDecompyFrac
use PRTParametersMod , only : prt_params
use FatesParameterDerivedMod, only : param_derived
use EDTypesMod , only : ed_site_type
use FatesPatchMod, only : fates_patch_type
use FatesCohortMod , only : fates_cohort_type
use EDParamsMod , only : nclmax
use PRTGenericMod , only : element_list
use PRTGenericMod , only : StorageNutrientTarget
use FatesLitterMod , only : ncwd
use FatesLitterMod , only : ndcmpy
use FatesLitterMod , only : litter_type
use FatesLitterMod , only : adjust_SF_CWD_frac
use EDParamsMod , only : max_cohort_per_patch
use EDTypesMod , only : AREA
use EDTypesMod , only : min_npm2, min_nppatch
use EDTypesMod , only : min_n_safemath
use EDParamsMod , only : nlevleaf
use PRTGenericMod , only : max_nleafage
use FatesConstantsMod , only : ican_upper
use EDTypesMod , only : site_fluxdiags_type
use EDTypesMod , only : elem_diag_type
use PRTGenericMod , only : num_elements
use FatesConstantsMod , only : leaves_on
use FatesConstantsMod , only : leaves_off
use FatesConstantsMod , only : leaves_shedding
use FatesConstantsMod , only : ihard_stress_decid
use FatesConstantsMod , only : isemi_stress_decid
use EDParamsMod , only : ED_val_cohort_age_fusion_tol
use FatesInterfaceTypesMod , only : hlm_use_planthydro
use FatesInterfaceTypesMod , only : hlm_parteh_mode
use FatesPlantHydraulicsMod, only : FuseCohortHydraulics
use FatesPlantHydraulicsMod, only : UpdateSizeDepPlantHydProps
use FatesPlantHydraulicsMod, only : InitPlantHydStates
use FatesPlantHydraulicsMod, only : InitHydrCohort
use FatesPlantHydraulicsMod, only : DeallocateHydrCohort
use FatesPlantHydraulicsMod, only : AccumulateMortalityWaterStorage
use FatesPlantHydraulicsMod, only : UpdatePlantHydrNodes
use FatesPlantHydraulicsMod, only : UpdatePlantHydrLenVol
use FatesPlantHydraulicsMod, only : UpdatePlantKmax
use FatesPlantHydraulicsMod, only : SavePreviousCompartmentVolumes
use FatesPlantHydraulicsMod, only : ConstrainRecruitNumber
use FatesSizeAgeTypeIndicesMod, only : sizetype_class_index
use FatesSizeAgeTypeIndicesMod, only : coagetype_class_index
use FatesAllometryMod , only : bleaf
use FatesAllometryMod , only : bfineroot
use FatesAllometryMod , only : bsap_allom
use FatesAllometryMod , only : bagw_allom
use FatesAllometryMod , only : bbgw_allom
use FatesAllometryMod , only : bdead_allom
use FatesAllometryMod , only : h_allom
use FatesAllometryMod , only : carea_allom
use FatesAllometryMod , only : bstore_allom
use FatesAllometryMod , only : ForceDBH
use FatesAllometryMod , only : set_root_fraction
use PRTGenericMod, only : prt_carbon_allom_hyp
use PRTGenericMod, only : prt_cnp_flex_allom_hyp
use PRTGenericMod, only : prt_vartypes
use PRTGenericMod, only : carbon12_element
use PRTGenericMod, only : nitrogen_element
use PRTGenericMod, only : phosphorus_element
use PRTGenericMod, only : leaf_organ
use PRTGenericMod, only : fnrt_organ
use PRTGenericMod, only : sapw_organ
use PRTGenericMod, only : store_organ
use PRTGenericMod, only : repro_organ
use PRTGenericMod, only : struct_organ
use PRTGenericMod, only : SetState
use PRTAllometricCarbonMod, only : callom_prt_vartypes
use PRTAllometricCarbonMod, only : ac_bc_inout_id_netdc
use PRTAllometricCarbonMod, only : ac_bc_in_id_pft
use PRTAllometricCarbonMod, only : ac_bc_in_id_ctrim
use PRTAllometricCarbonMod, only : ac_bc_inout_id_dbh
use PRTAllometricCarbonMod, only : ac_bc_in_id_lstat
use PRTAllometricCarbonMod, only : ac_bc_in_id_cdamage
use PRTAllometricCarbonMod, only : ac_bc_in_id_efleaf
use PRTAllometricCarbonMod, only : ac_bc_in_id_effnrt
use PRTAllometricCarbonMod, only : ac_bc_in_id_efstem
use PRTAllometricCNPMod, only : cnp_allom_prt_vartypes
use PRTAllometricCNPMod, only : acnp_bc_in_id_pft, acnp_bc_in_id_ctrim
use PRTAllometricCNPMod, only : acnp_bc_in_id_lstat, acnp_bc_inout_id_dbh
use PRTAllometricCNPMod, only : acnp_bc_in_id_efleaf
use PRTAllometricCNPMod, only : acnp_bc_in_id_effnrt
use PRTAllometricCNPMod, only : acnp_bc_in_id_efstem
use PRTAllometricCNPMod, only : acnp_bc_inout_id_l2fr
use PRTAllometricCNPMod, only : acnp_bc_inout_id_cx_int
use PRTAllometricCNPMod, only : acnp_bc_inout_id_cx0
use PRTAllometricCNPMod, only : acnp_bc_inout_id_emadcxdt
use PRTAllometricCNPMod, only : acnp_bc_in_id_nc_repro
use PRTAllometricCNPMod, only : acnp_bc_in_id_pc_repro
use PRTAllometricCNPMod, only : acnp_bc_inout_id_resp_excess, acnp_bc_in_id_netdc
use PRTAllometricCNPMod, only : acnp_bc_inout_id_netdn, acnp_bc_inout_id_netdp
use PRTAllometricCNPMod, only : acnp_bc_out_id_cefflux, acnp_bc_out_id_nefflux
use PRTAllometricCNPMod, only : acnp_bc_out_id_pefflux, acnp_bc_out_id_limiter
use PRTAllometricCNPMod, only : acnp_bc_in_id_cdamage
use DamageMainMod, only : undamaged_class
use FatesConstantsMod, only : n_term_mort_types
use FatesConstantsMod, only : i_term_mort_type_cstarv
use FatesConstantsMod, only : i_term_mort_type_canlev
use FatesConstantsMod, only : i_term_mort_type_numdens
use shr_infnan_mod, only : nan => shr_infnan_nan, assignment(=)
use shr_log_mod, only : errMsg => shr_log_errMsg
!
implicit none
private
!
public :: create_cohort
public :: terminate_cohorts
public :: terminate_cohort
public :: fuse_cohorts
public :: insert_cohort
public :: sort_cohorts
public :: count_cohorts
public :: InitPRTObject
public :: SendCohortToLitter
public :: EvaluateAndCorrectDBH
public :: DamageRecovery
logical, parameter :: debug = .false. ! local debug flag
character(len=*), parameter, private :: sourcefile = &
__FILE__
integer, parameter, private :: conserve_crownarea_and_number_not_dbh = 1
integer, parameter, private :: conserve_dbh_and_number_not_crownarea = 2
integer, parameter, private :: cohort_fusion_conservation_method = conserve_crownarea_and_number_not_dbh
! 10/30/09: Created by Rosie Fisher
!-------------------------------------------------------------------------------------!
contains
!-------------------------------------------------------------------------------------!
subroutine create_cohort(currentSite, patchptr, pft, nn, height, coage, dbh, &
prt, elongf_leaf, elongf_fnrt, elongf_stem, status, recruitstatus, ctrim, &
carea, clayer, crowndamage, spread, bc_in)
!
! DESCRIPTION:
! create new cohort
! There are 4 places this is called
! 1) Initializing new cohorts at the beginning of a cold-start simulation
! 2) Initializing new recruits during dynamics
! 3) Initializing new cohorts at the beginning of a inventory read
! 4) Initializing new cohorts during restart
!
! It is assumed that in the first 3, this is called with a reasonable amount of starter information.
!
! ARGUMENTS:
type(ed_site_type), intent(inout), target :: currentSite ! site object
type(fates_patch_type), intent(inout), pointer :: patchptr ! pointer to patch object
integer, intent(in) :: pft ! cohort Plant Functional Type
integer, intent(in) :: crowndamage ! cohort damage class
integer, intent(in) :: clayer ! canopy status of cohort [1=canopy; 2=understorey]
integer, intent(in) :: status ! growth status of plant [1=leaves off; 2=leaves on]
integer, intent(in) :: recruitstatus ! recruit status of plant [1 = recruitment , 0 = other]
real(r8), intent(in) :: nn ! number of individuals in cohort [/m2]
real(r8), intent(in) :: height ! cohort height [m]
real(r8), intent(in) :: coage ! cohort age [m]
real(r8), intent(in) :: dbh ! cohort diameter at breast height [cm]
real(r8), intent(in) :: elongf_leaf ! leaf elongation factor [fraction] - 0: fully abscissed; 1: fully flushed
real(r8), intent(in) :: elongf_fnrt ! fine-root "elongation factor" [fraction]
real(r8), intent(in) :: elongf_stem ! stem "elongation factor" [fraction]
class(prt_vartypes), intent(inout), pointer :: prt ! allocated PARTEH object
real(r8), intent(in) :: ctrim ! fraction of the maximum leaf biomass we are targeting
real(r8), intent(in) :: spread ! how spread crowns are in horizontal space
real(r8), intent(in) :: carea ! area of cohort - ONLY USED IN SP MODE [m2]
type(bc_in_type), intent(in) :: bc_in ! external boundary conditions
! LOCAL VARIABLES:
type(fates_cohort_type), pointer :: newCohort ! pointer to New Cohort structure.
type(fates_cohort_type), pointer :: storesmallcohort
type(fates_cohort_type), pointer :: storebigcohort
real(r8) :: rmean_temp ! running mean temperature
integer :: tnull, snull ! are the tallest and shortest cohorts allocate
integer :: nlevrhiz ! number of rhizosphere layers
!----------------------------------------------------------------------
! create new cohort
allocate(newCohort)
call newCohort%Create(prt, pft, nn, height, coage, dbh, status, ctrim, carea, &
clayer, crowndamage, spread, patchptr%canopy_layer_tlai, elongf_leaf, elongf_fnrt, &
elongf_stem)
! Put cohort at the right place in the linked list
storebigcohort => patchptr%tallest
storesmallcohort => patchptr%shortest
if (associated(patchptr%tallest)) then
tnull = 0
else
tnull = 1
patchptr%tallest => newCohort
endif
if (associated(patchptr%shortest)) then
snull = 0
else
snull = 1
patchptr%shortest => newCohort
endif
! Allocate running mean functions
! (Keeping as an example)
!! allocate(newCohort%tveg_lpa)
!! call newCohort%tveg_lpa%InitRMean(ema_lpa,init_value=patchptr%tveg_lpa%GetMean())
if (hlm_use_planthydro .eq. itrue) then
nlevrhiz = currentSite%si_hydr%nlevrhiz
! This allocates array spaces
call InitHydrCohort(currentSite, newCohort)
! zero out the water balance error
newCohort%co_hydr%errh2o = 0._r8
! This calculates node heights
call UpdatePlantHydrNodes(newCohort, newCohort%pft, &
newCohort%height,currentSite%si_hydr)
! This calculates volumes and lengths
call UpdatePlantHydrLenVol(newCohort,currentSite%si_hydr)
! This updates the Kmax's of the plant's compartments
call UpdatePlantKmax(newCohort%co_hydr,newCohort,currentSite%si_hydr)
! Since this is a newly initialized plant, we set the previous compartment-size
! equal to the ones we just calculated.
call SavePreviousCompartmentVolumes(newCohort%co_hydr)
! This comes up with starter suctions and then water contents
! based on the soil values
call InitPlantHydStates(currentSite,newCohort)
if(recruitstatus==1)then
newCohort%co_hydr%is_newly_recruited = .true.
! If plant hydraulics is active, we must constrain the
! number density of the new recruits based on the moisture
! available to be subsumed in the new plant tissues.
! So we go through the process of pre-initializing the hydraulic
! states in the temporary cohort, to calculate this new number density
rmean_temp = patchptr%tveg24%GetMean()
call ConstrainRecruitNumber(currentSite, newCohort, patchptr, &
bc_in, rmean_temp)
endif
endif
call insert_cohort(patchptr, newCohort, patchptr%tallest, patchptr%shortest, tnull, snull, &
storebigcohort, storesmallcohort)
patchptr%tallest => storebigcohort
patchptr%shortest => storesmallcohort
end subroutine create_cohort
! ------------------------------------------------------------------------------------!
subroutine InitPRTObject(prt)
! -----------------------------------------------------------------------------------
!
! This routine allocates the PARTEH object that is associated with each cohort.
! The argument that is passed in is a pointer that is then associated with this
! newly allocated object.
! The object that is allocated is the specific extended class for the hypothesis
! of choice.
! Following this, the object and its internal mappings are initialized.
! This routine does NOT set any of the initial conditions, or boundary conditions
! such as the organ/element masses. Those are handled after this call.
!
! -----------------------------------------------------------------------------------
! Argument
class(prt_vartypes), pointer :: prt
! Potential Extended types
class(callom_prt_vartypes), pointer :: c_allom_prt
class(cnp_allom_prt_vartypes), pointer :: cnp_allom_prt
select case(hlm_parteh_mode)
case (prt_carbon_allom_hyp)
allocate(c_allom_prt)
prt => c_allom_prt
case (prt_cnp_flex_allom_hyp)
allocate(cnp_allom_prt)
prt => cnp_allom_prt
case DEFAULT
write(fates_log(),*) 'You specified an unknown PRT module'
write(fates_log(),*) 'Aborting'
call endrun(msg=errMsg(sourcefile, __LINE__))
end select
! This is the call to allocate the data structures in the PRT object
! This call will be extended to each specific class.
call prt%InitPRTVartype()
return
end subroutine InitPRTObject
!-------------------------------------------------------------------------------------!
subroutine terminate_cohorts( currentSite, currentPatch, level , call_index, bc_in)
!
! !DESCRIPTION:
! terminates all cohorts when they get too small
!
! !USES:
!
! !ARGUMENTS
type (ed_site_type) , intent(inout) :: currentSite
type (fates_patch_type), intent(inout) :: currentPatch
integer , intent(in) :: level
integer :: call_index
type(bc_in_type), intent(in) :: bc_in
! Important point regarding termination levels. Termination is typically
! called after fusion. We do this so that we can re-capture the biomass that would
! otherwise be lost from termination. The biomass of a fused plant remains in the
! live pool. However, some plant number densities can be so low that they
! can cause numerical instabilities. Thus, we call terminate_cohorts at level=1
! before fusion to get rid of these cohorts that are so incredibly sparse, and then
! terminate the remainder at level 2 for various other reasons.
!
! !LOCAL VARIABLES:
type (fates_cohort_type) , pointer :: currentCohort
type (fates_cohort_type) , pointer :: shorterCohort
type (fates_cohort_type) , pointer :: tallerCohort
real(r8) :: leaf_c ! leaf carbon [kg]
real(r8) :: store_c ! storage carbon [kg]
real(r8) :: sapw_c ! sapwood carbon [kg]
real(r8) :: fnrt_c ! fineroot carbon [kg]
real(r8) :: repro_c ! reproductive carbon [kg]
real(r8) :: struct_c ! structural carbon [kg]
integer :: terminate ! do we terminate (itrue) or not (ifalse)
integer :: istat ! return status code
character(len=255) :: smsg
integer :: termination_type
!----------------------------------------------------------------------
currentCohort => currentPatch%shortest
do while (associated(currentCohort))
terminate = ifalse
termination_type = 0
tallerCohort => currentCohort%taller
leaf_c = currentCohort%prt%GetState(leaf_organ, carbon12_element)
store_c = currentCohort%prt%GetState(store_organ, carbon12_element)
sapw_c = currentCohort%prt%GetState(sapw_organ, carbon12_element)
fnrt_c = currentCohort%prt%GetState(fnrt_organ, carbon12_element)
struct_c = currentCohort%prt%GetState(struct_organ, carbon12_element)
repro_c = currentCohort%prt%GetState(repro_organ, carbon12_element)
! Check if number density is so low is breaks math (level 1)
if (currentcohort%n < min_n_safemath .and. level == 1) then
terminate = itrue
termination_type = i_term_mort_type_numdens
if ( debug ) then
write(fates_log(),*) 'terminating cohorts 0',currentCohort%n/currentPatch%area,currentCohort%dbh,currentCohort%pft,call_index
endif
endif
! The rest of these are only allowed if we are not dealing with a recruit (level 2)
if (.not.currentCohort%isnew .and. level == 2) then
! Not enough n or dbh
if (currentCohort%n/currentPatch%area <= min_npm2 .or. & !
currentCohort%n <= min_nppatch .or. &
(currentCohort%dbh < 0.00001_r8 .and. store_c < 0._r8) ) then
terminate = itrue
termination_type = i_term_mort_type_numdens
if ( debug ) then
write(fates_log(),*) 'terminating cohorts 1',currentCohort%n/currentPatch%area,currentCohort%dbh,currentCohort%pft,call_index
endif
endif
! Outside the maximum canopy layer
if (currentCohort%canopy_layer > nclmax ) then
terminate = itrue
termination_type = i_term_mort_type_canlev
if ( debug ) then
write(fates_log(),*) 'terminating cohorts 2', currentCohort%canopy_layer,currentCohort%pft,call_index
endif
endif
! live biomass pools are terminally depleted
if ( ( sapw_c+leaf_c+fnrt_c ) < 1e-10_r8 .or. &
store_c < 1e-10_r8) then
terminate = itrue
termination_type = i_term_mort_type_cstarv
if ( debug ) then
write(fates_log(),*) 'terminating cohorts 3', &
sapw_c,leaf_c,fnrt_c,store_c,currentCohort%pft,call_index
endif
endif
! Total cohort biomass is negative
if ( ( struct_c+sapw_c+leaf_c+fnrt_c+store_c ) < 0._r8) then
terminate = itrue
termination_type = i_term_mort_type_cstarv
if ( debug ) then
write(fates_log(),*) 'terminating cohorts 4', &
struct_c,sapw_c,leaf_c,fnrt_c,store_c,currentCohort%pft,call_index
endif
endif
endif ! if (.not.currentCohort%isnew .and. level == 2) then
if (terminate == itrue) then
call terminate_cohort(currentSite, currentPatch, currentCohort, bc_in, termination_type)
deallocate(currentCohort, stat=istat, errmsg=smsg)
if (istat/=0) then
write(fates_log(),*) 'dealloc001: fail on terminate_cohorts:deallocate(currentCohort):'//trim(smsg)
call endrun(msg=errMsg(sourcefile, __LINE__))
endif
endif
currentCohort => tallerCohort
enddo
end subroutine terminate_cohorts
!-------------------------------------------------------------------------------------!
subroutine terminate_cohort(currentSite, currentPatch, currentCohort, bc_in, termination_type)
!
! !DESCRIPTION:
! Terminates an individual cohort and updates the site-level
! updates the carbon flux and nuber of individuals appropriately
!
! !USES:
!
! !ARGUMENTS
type (ed_site_type) , intent(inout), target :: currentSite
type (fates_patch_type) , intent(inout), target :: currentPatch
type (fates_cohort_type), intent(inout), target :: currentCohort
type(bc_in_type), intent(in) :: bc_in
integer, intent(in) :: termination_type
! !LOCAL VARIABLES:
type (fates_cohort_type) , pointer :: shorterCohort
type (fates_cohort_type) , pointer :: tallerCohort
real(r8) :: leaf_c ! leaf carbon [kg]
real(r8) :: store_c ! storage carbon [kg]
real(r8) :: sapw_c ! sapwood carbon [kg]
real(r8) :: fnrt_c ! fineroot carbon [kg]
real(r8) :: repro_c ! reproductive carbon [kg]
real(r8) :: struct_c ! structural carbon [kg]
integer :: terminate ! do we terminate (itrue) or not (ifalse)
integer :: c ! counter for litter size class.
integer :: levcan ! canopy level
!----------------------------------------------------------------------
! check termination_type; it should not be 0
if (termination_type == 0) then
write(fates_log(),*) 'termination_type=0'
call endrun(msg=errMsg(sourcefile, __LINE__))
endif
leaf_c = currentCohort%prt%GetState(leaf_organ, carbon12_element)
store_c = currentCohort%prt%GetState(store_organ, carbon12_element)
sapw_c = currentCohort%prt%GetState(sapw_organ, carbon12_element)
fnrt_c = currentCohort%prt%GetState(fnrt_organ, carbon12_element)
struct_c = currentCohort%prt%GetState(struct_organ, carbon12_element)
repro_c = currentCohort%prt%GetState(repro_organ, carbon12_element)
! preserve a record of the to-be-terminated cohort for mortality accounting
levcan = currentCohort%canopy_layer
if( hlm_use_planthydro == itrue ) &
call AccumulateMortalityWaterStorage(currentSite,currentCohort,currentCohort%n)
! Update the site-level carbon flux and individuals count for the appropriate canopy layer
if(levcan==ican_upper) then
currentSite%term_nindivs_canopy(termination_type,currentCohort%size_class,currentCohort%pft) = &
currentSite%term_nindivs_canopy(termination_type,currentCohort%size_class,currentCohort%pft) + currentCohort%n
currentSite%term_carbonflux_canopy(termination_type,currentCohort%pft) = currentSite%term_carbonflux_canopy(termination_type,currentCohort%pft) + &
currentCohort%n * (struct_c+sapw_c+leaf_c+fnrt_c+store_c+repro_c)
else
currentSite%term_nindivs_ustory(termination_type,currentCohort%size_class,currentCohort%pft) = &
currentSite%term_nindivs_ustory(termination_type,currentCohort%size_class,currentCohort%pft) + currentCohort%n
currentSite%term_carbonflux_ustory(termination_type,currentCohort%pft) = currentSite%term_carbonflux_ustory(termination_type,currentCohort%pft) + &
currentCohort%n * (struct_c+sapw_c+leaf_c+fnrt_c+store_c+repro_c)
end if
currentSite%term_abg_flux(currentCohort%size_class, currentCohort%pft) = &
currentSite%term_abg_flux(currentCohort%size_class, currentCohort%pft) + &
currentCohort%n * ( (struct_c+sapw_c+store_c) * prt_params%allom_agb_frac(currentCohort%pft) + &
leaf_c )
! put the litter from the terminated cohorts
! straight into the fragmenting pools
if (currentCohort%n.gt.0.0_r8) then
call SendCohortToLitter(currentSite,currentPatch, &
currentCohort,currentCohort%n,bc_in)
end if
! Set pointers and deallocate the current cohort from the list
shorterCohort => currentCohort%shorter
tallerCohort => currentCohort%taller
if (.not. associated(tallerCohort)) then
currentPatch%tallest => shorterCohort
if(associated(shorterCohort)) shorterCohort%taller => null()
else
tallerCohort%shorter => shorterCohort
endif
if (.not. associated(shorterCohort)) then
currentPatch%shortest => tallerCohort
if(associated(tallerCohort)) tallerCohort%shorter => null()
else
shorterCohort%taller => tallerCohort
endif
call currentCohort%FreeMemory()
end subroutine terminate_cohort
! =====================================================================================
subroutine SendCohortToLitter(csite,cpatch,ccohort,nplant,bc_in)
! -----------------------------------------------------------------------------------
! This routine transfers the existing mass in all pools and all elements
! on a vegetation cohort, into the litter pool.
!
! Important: (1) This IS NOT turnover, this is not a partial transfer.
! (2) This is from a select number of plants in the cohort. ie this is
! not a "whole-sale" sending of all plants to litter.
! (3) This does not affect the PER PLANT mass pools, so
! do not update any PARTEH structures.
! (4) The change in plant number density (due to death or termination)
! IS NOT handled here.
! (5) This routine is NOT used for disturbance, mostly
! because this routine assumes a cohort lands in its patch
! Whereas the disturbance scheme does NOT assume that.
! -----------------------------------------------------------------------------------
! Arguments
type (ed_site_type) , target :: csite
type (fates_patch_type) , target :: cpatch
type (fates_cohort_type) , target :: ccohort
real(r8) :: nplant ! Number (absolute)
! of plants to transfer
type(bc_in_type), intent(in) :: bc_in
type(litter_type), pointer :: litt ! Litter object for each element
type(elem_diag_type),pointer :: elflux_diags
real(r8) :: leaf_m ! leaf mass [kg]
real(r8) :: store_m ! storage mass [kg]
real(r8) :: sapw_m ! sapwood mass [kg]
real(r8) :: fnrt_m ! fineroot mass [kg]
real(r8) :: repro_m ! reproductive mass [kg]
real(r8) :: struct_m ! structural mass [kg]
real(r8) :: plant_dens! plant density [/m2]
real(r8) :: dcmpy_frac! fraction of mass going to each decomposability partition
integer :: el ! loop index for elements
integer :: c ! loop index for CWD
integer :: pft ! pft index of the cohort
integer :: crowndamage ! the crown damage class of the cohort
integer :: sl ! loop index for soil layers
integer :: dcmpy ! loop index for decomposability
real(r8) :: SF_val_CWD_frac_adj(4) !Updated wood partitioning to CWD based on dbh
!----------------------------------------------------------------------
pft = ccohort%pft
plant_dens = nplant/cpatch%area
call set_root_fraction(csite%rootfrac_scr, pft, csite%zi_soil, &
bc_in%max_rooting_depth_index_col)
do el=1,num_elements
store_m = ccohort%prt%GetState(store_organ, element_list(el))
fnrt_m = ccohort%prt%GetState(fnrt_organ, element_list(el))
repro_m = ccohort%prt%GetState(repro_organ, element_list(el))
if (prt_params%woody(ccohort%pft) == itrue) then
leaf_m = ccohort%prt%GetState(leaf_organ, element_list(el))
sapw_m = ccohort%prt%GetState(sapw_organ, element_list(el))
struct_m = ccohort%prt%GetState(struct_organ, element_list(el))
else
leaf_m = ccohort%prt%GetState(leaf_organ, element_list(el)) + &
ccohort%prt%GetState(sapw_organ, element_list(el)) + &
ccohort%prt%GetState(struct_organ, element_list(el))
sapw_m = 0._r8
struct_m = 0._r8
endif
litt => cpatch%litter(el)
elflux_diags => csite%flux_diags%elem(el)
!adjust how wood is partitioned between the cwd classes based on cohort dbh
call adjust_SF_CWD_frac(ccohort%dbh,ncwd,SF_val_CWD_frac,SF_val_CWD_frac_adj)
do c=1,ncwd
! above ground CWD
litt%ag_cwd(c) = litt%ag_cwd(c) + plant_dens * &
(struct_m+sapw_m) * SF_val_CWD_frac_adj(c) * &
prt_params%allom_agb_frac(pft)
! below ground CWD
do sl=1,csite%nlevsoil
litt%bg_cwd(c,sl) = litt%bg_cwd(c,sl) + plant_dens * &
(struct_m+sapw_m) * SF_val_CWD_frac_adj(c) * &
(1.0_r8 - prt_params%allom_agb_frac(pft)) * &
csite%rootfrac_scr(sl)
enddo
! above ground
elflux_diags%cwd_ag_input(c) = elflux_diags%cwd_ag_input(c) + &
(struct_m+sapw_m) * SF_val_CWD_frac_adj(c) * &
prt_params%allom_agb_frac(pft) * nplant
! below ground
elflux_diags%cwd_bg_input(c) = elflux_diags%cwd_bg_input(c) + &
(struct_m + sapw_m) * SF_val_CWD_frac_adj(c) * &
(1.0_r8 - prt_params%allom_agb_frac(pft)) * nplant
enddo
do dcmpy=1,ndcmpy
dcmpy_frac = GetDecompyFrac(pft,leaf_organ,dcmpy)
litt%leaf_fines(dcmpy) = litt%leaf_fines(dcmpy) + &
plant_dens * (leaf_m+repro_m) * dcmpy_frac
dcmpy_frac = GetDecompyFrac(pft,fnrt_organ,dcmpy)
do sl=1,csite%nlevsoil
litt%root_fines(dcmpy,sl) = litt%root_fines(dcmpy,sl) + &
plant_dens * (fnrt_m+store_m) * csite%rootfrac_scr(sl) * dcmpy_frac
end do
end do
elflux_diags%surf_fine_litter_input(pft) = &
elflux_diags%surf_fine_litter_input(pft) + &
(leaf_m+repro_m) * nplant
elflux_diags%root_litter_input(pft) = &
elflux_diags%root_litter_input(pft) + &
(fnrt_m+store_m) * nplant
end do
return
end subroutine SendCohortToLitter
!--------------------------------------------------------------------------------------
subroutine fuse_cohorts(currentSite, currentPatch, bc_in)
!
! !DESCRIPTION:
! Join similar cohorts to reduce total number
!
! !USES:
use EDParamsMod , only : ED_val_cohort_size_fusion_tol
use EDParamsMod , only : ED_val_cohort_age_fusion_tol
use FatesInterfaceTypesMod , only : hlm_use_cohort_age_tracking
use FatesConstantsMod , only : itrue
use FatesConstantsMod, only : days_per_year
!
! !ARGUMENTS
type (ed_site_type), intent(inout) :: currentSite
type (fates_patch_type), intent(inout), pointer :: currentPatch
type (bc_in_type), intent(in) :: bc_in
!
! !LOCAL VARIABLES:
type (fates_cohort_type) , pointer :: currentCohort
type (fates_cohort_type) , pointer :: nextc
type (fates_cohort_type) , pointer :: nextnextc
type (fates_cohort_type) , pointer :: shorterCohort
type (fates_cohort_type) , pointer :: tallerCohort
integer :: i
integer :: fusion_took_place
integer :: iterate ! do we need to keep fusing to get below maxcohorts?
integer :: nocohorts
real(r8) :: newn
real(r8) :: diff
real(r8) :: coage_diff
real(r8) :: leaf_c_next ! Leaf carbon * plant density of current (for weighting)
real(r8) :: leaf_c_curr ! Leaf carbon * plant density of next (for weighting)
real(r8) :: leaf_c_target
real(r8) :: dynamic_size_fusion_tolerance
real(r8) :: dynamic_age_fusion_tolerance
real(r8) :: dbh
real(r8) :: leaf_c ! leaf carbon [kg]
real(r8) :: target_storec ! Target storage C
integer :: largersc, smallersc, sc_i ! indices for tracking the growth flux caused by fusion
real(r8) :: larger_n, smaller_n
integer :: oldercacls, youngercacls, cacls_i ! indices for tracking the age flux caused by fusion
real(r8) :: older_n, younger_n
real(r8) :: crown_reduction
logical, parameter :: fuse_debug = .false. ! This debug is over-verbose
! and gets its own flag
integer :: istat ! return status code
character(len=255) :: smsg
!----------------------------------------------------------------------
!set initial fusion tolerance (in cm)
dynamic_size_fusion_tolerance = ED_val_cohort_size_fusion_tol
! set the cohort age fusion tolerance (in fraction of years)
dynamic_age_fusion_tolerance = ED_val_cohort_age_fusion_tol
!This needs to be a function of the canopy layer, because otherwise, at canopy closure
!the number of cohorts doubles and very dissimilar cohorts are fused together
!because c_area and biomass are non-linear with dbh, this causes several mass inconsistancies
!in theory, all of this routine therefore causes minor losses of C and area, but these are below
!detection limit normally.
iterate = 1
fusion_took_place = 0
!---------------------------------------------------------------------!
! Keep doing this until nocohorts <= maxcohorts !
!---------------------------------------------------------------------!
if (associated(currentPatch%shortest)) then
do while(iterate == 1)
currentCohort => currentPatch%tallest
! The following logic continues the loop while the current cohort is not the shortest cohort
! if they point to the same target (ie equivalence), then the loop ends.
! This loop is different than the simple "continue while associated" loop in that
! it omits the last cohort (because it has already been compared by that point)
do while ( .not.associated(currentCohort,currentPatch%shortest) )
nextc => currentPatch%tallest
do while (associated(nextc))
nextnextc => nextc%shorter
diff = abs((currentCohort%dbh - nextc%dbh)/(0.5_r8*(currentCohort%dbh + nextc%dbh)))
!Criteria used to divide up the height continuum into different cohorts.
if (diff < dynamic_size_fusion_tolerance) then
! Only fuse if the cohorts are within x years of each other
! if they are the same age we make diff 0- to avoid errors divding by zero
!NB if cohort age tracking is off then the age of both should be 0
! and hence the age fusion criterion is met
if (abs(currentCohort%coage - nextc%coage)<nearzero ) then
coage_diff = 0.0_r8
else
coage_diff = abs((currentCohort%coage - nextc%coage)/ &
(0.5_r8*(currentCohort%coage + nextc%coage)))
end if
if (coage_diff <= dynamic_age_fusion_tolerance ) then
! Don't fuse a cohort with itself!
if (.not.associated(currentCohort,nextc) ) then
if (currentCohort%pft == nextc%pft) then
! check cohorts have same damage class before fusing
if (currentCohort%crowndamage == nextc%crowndamage) then
! check cohorts in same c. layer. before fusing
if (currentCohort%canopy_layer == nextc%canopy_layer) then
! Note: because newly recruited cohorts that have not experienced
! a day yet will have un-known flux quantities or change rates
! we don't want them fusing with non-new cohorts. We allow them
! to fuse with other new cohorts to keep the total number of cohorts
! down.
if( currentCohort%isnew.eqv.nextc%isnew ) then
newn = currentCohort%n + nextc%n
fusion_took_place = 1
if ( fuse_debug .and. currentCohort%isnew ) then
write(fates_log(),*) 'Fusing Two Cohorts'
write(fates_log(),*) 'newn: ',newn
write(fates_log(),*) 'Cohort I, Cohort II'
write(fates_log(),*) 'n:',currentCohort%n,nextc%n
write(fates_log(),*) 'isnew:',currentCohort%isnew,nextc%isnew
write(fates_log(),*) 'height:',currentCohort%height,nextc%height
write(fates_log(),*) 'coage:',currentCohort%coage,nextc%coage
write(fates_log(),*) 'dbh:',currentCohort%dbh,nextc%dbh
write(fates_log(),*) 'pft:',currentCohort%pft,nextc%pft
write(fates_log(),*) 'crowndamage:',currentCohort%crowndamage,nextc%crowndamage
write(fates_log(),*) 'canopy_trim:',currentCohort%canopy_trim,nextc%canopy_trim
write(fates_log(),*) 'canopy_layer_yesterday:', &
currentCohort%canopy_layer_yesterday,nextc%canopy_layer_yesterday
do i=1, nlevleaf
write(fates_log(),*) 'leaf level: ',i,'year_net_uptake', &
currentCohort%year_net_uptake(i),nextc%year_net_uptake(i)
end do
end if
! new cohort age is weighted mean of two cohorts
currentCohort%coage = &
(currentCohort%coage * (currentCohort%n/(currentCohort%n + nextc%n))) + &
(nextc%coage * (nextc%n/(currentCohort%n + nextc%n)))
! update the cohort age again
if (hlm_use_cohort_age_tracking .eq.itrue) then
call coagetype_class_index(currentCohort%coage, currentCohort%pft, &
currentCohort%coage_class, currentCohort%coage_by_pft_class)
end if
! Fuse all mass pools
call currentCohort%prt%WeightedFusePRTVartypes(nextc%prt, &
currentCohort%n/newn )
! Leaf biophysical rates (use leaf mass weighting)
! -----------------------------------------------------------------
call currentCohort%UpdateCohortBioPhysRates()
currentCohort%l2fr = (currentCohort%n*currentCohort%l2fr &
+ nextc%n*nextc%l2fr)/newn
currentCohort%canopy_trim = (currentCohort%n*currentCohort%canopy_trim &
+ nextc%n*nextc%canopy_trim)/newn
! c13disc_acc calculation; weighted mean by GPP
if ((currentCohort%n * currentCohort%gpp_acc + nextc%n * nextc%gpp_acc) .eq. 0.0_r8) then
currentCohort%c13disc_acc = 0.0_r8
else
currentCohort%c13disc_acc = (currentCohort%n * currentCohort%gpp_acc * currentCohort%c13disc_acc + &
nextc%n * nextc%gpp_acc * nextc%c13disc_acc)/ &
(currentCohort%n * currentCohort%gpp_acc + nextc%n * nextc%gpp_acc)
endif
select case(cohort_fusion_conservation_method)
!
! -----------------------------------------------------------------
! Because cohort fusion is an unavoidable but non-physical process,
! and because of the various nonlinear allometric relationships,
! it isn't possible to simultaneously conserve all of the allometric
! relationships during cohort fusion. We will always conserve carbon,
! but there are choices to made about what else to conserve or not.
! In particular, there is a choice to be made of conservation amongst
! the number density, stem diameter, and crown area. Below,
! some different conservation relationships can be chosen during fusion.
! -----------------------------------------------------------------
!
case(conserve_crownarea_and_number_not_dbh)
!
! -----------------------------------------------------------------
! conserve total crown area during the fusion step, and then calculate
! dbh of the fused cohort as that which conserves both crown area and
! the dbh to crown area allometry. dbh will be updated in the next
! growth step in the (likely) event that dbh to structural iomass
! allometry is exceeded. if using a capped crown area allometry and
! above the cap, then calculate as the weighted average of fusing
! cohorts' dbh
! -----------------------------------------------------------------
!
call carea_allom(currentCohort%dbh,currentCohort%n, &
currentSite%spread,currentCohort%pft,&
currentCohort%crowndamage, &
currentCohort%c_area,inverse=.false.)
call carea_allom(nextc%dbh,nextc%n, &
currentSite%spread,nextc%pft,&
nextc%crowndamage, &
nextc%c_area,inverse=.false.)
currentCohort%c_area = currentCohort%c_area + nextc%c_area
!
dbh = currentCohort%dbh
call carea_allom(dbh,newn,currentSite%spread,currentCohort%pft,&
currentCohort%crowndamage,currentCohort%c_area,inverse=.true.)
!
if (abs(dbh-fates_unset_r8)<nearzero) then
currentCohort%dbh = (currentCohort%n*currentCohort%dbh &
+ nextc%n*nextc%dbh)/newn
if( prt_params%woody(currentCohort%pft) == itrue ) then
call ForceDBH( currentCohort%pft, currentCohort%crowndamage, &
currentCohort%canopy_trim, &
currentCohort%efleaf_coh, currentCohort%efstem_coh, &
currentCohort%dbh, currentCohort%height, &
bdead = currentCohort%prt%GetState(struct_organ,carbon12_element))
end if
!
call carea_allom(currentCohort%dbh,newn,currentSite%spread,currentCohort%pft,&
currentCohort%crowndamage, currentCohort%c_area,inverse=.false.)
else
currentCohort%dbh = dbh
endif
!
call h_allom(currentCohort%dbh,currentCohort%pft,currentCohort%height)
!
case(conserve_dbh_and_number_not_crownarea)
!
! -----------------------------------------------------------------
! Here we conserve the mean stem diameter of the trees in the cohorts
! rather than the crown area of the cohort
! -----------------------------------------------------------------
!
currentCohort%dbh = (currentCohort%n*currentCohort%dbh &
+ nextc%n*nextc%dbh)/newn
!
call h_allom(currentCohort%dbh,currentCohort%pft,currentCohort%height)
!
! -----------------------------------------------------------------
! If fusion pushed structural biomass to be larger than
! the allometric target value derived by diameter, we
! then increase diameter and height until the allometric
! target matches actual bdead. (if it is the other way around
! we then just let the carbon pools grow to fill out allometry)
! -----------------------------------------------------------------
!
if( prt_params%woody(currentCohort%pft) == itrue ) then
call ForceDBH( currentCohort%pft, currentCohort%crowndamage, &
currentCohort%canopy_trim, &
currentCohort%efleaf_coh, currentCohort%efstem_coh, &
currentCohort%dbh, currentCohort%height, &
bdead = currentCohort%prt%GetState(struct_organ,carbon12_element))
end if