-
Notifications
You must be signed in to change notification settings - Fork 0
/
FatesPlantRespPhotosynthMod.F90
1992 lines (1607 loc) · 103 KB
/
FatesPlantRespPhotosynthMod.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 FATESPlantRespPhotosynthMod
!-------------------------------------------------------------------------------------
! !DESCRIPTION:
! Calculates the plant respiration and photosynthetic fluxes for the FATES model
! This code is similar to and was originally based off of the 'photosynthesis'
! subroutine in the CLM model.
!
! Parameter for activation and deactivation energies were taken from:
! Activation energy, from:
! Bernacchi et al (2001) Plant, Cell and Environment 24:253-259
! Bernacchi et al (2003) Plant, Cell and Environment 26:1419-1430
! except TPU from: Harley et al (1992) Plant, Cell and Environment 15:271-282
! High temperature deactivation, from:
! Leuning (2002) Plant, Cell and Environment 25:1205-1210
! The factor "c" scales the deactivation to a value of 1.0 at 25C
! Photosynthesis and stomatal conductance parameters, from:
! Bonan et al (2011) JGR, 116, doi:10.1029/2010JG001593
! ------------------------------------------------------------------------------------
! !USES:
use FatesGlobals, only : endrun => fates_endrun
use FatesGlobals, only : fates_log
use FatesConstantsMod, only : r8 => fates_r8
use FatesConstantsMod, only : itrue,ifalse
use FatesConstantsMod, only : nearzero
use FatesInterfaceMod, only : hlm_use_planthydro
use FatesInterfaceMod, only : hlm_parteh_mode
use FatesInterfaceMod, only : numpft
use FatesInterfaceMod, only : nleafage
use EDTypesMod, only : maxpft
use EDTypesMod, only : nlevleaf
use EDTypesMod, only : nclmax
use EDTypesMod, only : max_nleafage
use EDTypesMod, only : do_fates_salinity
use EDParamsMod, only : q10_mr
use PRTGenericMod, only : prt_carbon_allom_hyp
use PRTGenericMod, only : prt_cnp_flex_allom_hyp
use PRTGenericMod, only : all_carbon_elements
use PRTGenericMod, only : nitrogen_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 EDParamsMod, only : ED_val_bbopt_c3, ED_val_bbopt_c4, ED_val_base_mr_20
! CIME Globals
use shr_log_mod , only : errMsg => shr_log_errMsg
implicit none
private
public :: FatesPlantRespPhotosynthDrive ! Called by the HLM-Fates interface
character(len=*), parameter, private :: sourcefile = &
__FILE__
!-------------------------------------------------------------------------------------
! maximum stomatal resistance [s/m] (used across several procedures)
real(r8),parameter :: rsmax0 = 2.e8_r8
logical :: debug = .false.
contains
!--------------------------------------------------------------------------------------
subroutine FatesPlantRespPhotosynthDrive (nsites, sites,bc_in,bc_out,dtime)
! -----------------------------------------------------------------------------------
! !DESCRIPTION:
! Leaf photosynthesis and stomatal conductance calculation as described by
! Bonan et al (2011) JGR, 116, doi:10.1029/2010JG001593 and extended to
! a multi-layer canopy
! -----------------------------------------------------------------------------------
! !USES:
use EDPftvarcon , only : EDPftvarcon_inst
use FatesSynchronizedParamsMod , only : FatesSynchronizedParamsInst
use EDTypesMod , only : ed_patch_type
use EDTypesMod , only : ed_cohort_type
use EDTypesMod , only : ed_site_type
use EDTypesMod , only : maxpft
use EDTypesMod , only : dinc_ed
use FatesInterfaceMod , only : bc_in_type
use FatesInterfaceMod , only : bc_out_type
use EDCanopyStructureMod, only : calc_areaindex
use FatesConstantsMod, only : umolC_to_kgC
use FatesConstantsMod, only : g_per_kg
use FatesConstantsMod, only : umol_per_mmol
use FatesConstantsMod, only : rgas => rgas_J_K_kmol
use FatesConstantsMod, only : tfrz => t_water_freeze_k_1atm
use FatesParameterDerivedMod, only : param_derived
use FatesAllometryMod, only : bleaf
use FatesAllometryMod, only : storage_fraction_of_target
use FatesAllometryMod, only : set_root_fraction
use FatesAllometryMod, only : i_hydro_rootprof_context
use FatesAllometryMod, only : decay_coeff_kn
! ARGUMENTS:
! -----------------------------------------------------------------------------------
integer,intent(in) :: nsites
type(ed_site_type),intent(inout),target :: sites(nsites)
type(bc_in_type),intent(in) :: bc_in(nsites)
type(bc_out_type),intent(inout) :: bc_out(nsites)
real(r8),intent(in) :: dtime
! LOCAL VARIABLES:
! -----------------------------------------------------------------------------------
type (ed_patch_type) , pointer :: currentPatch
type (ed_cohort_type), pointer :: currentCohort
! -----------------------------------------------------------------------------------
! These three arrays hold leaf-level biophysical rates that are calculated
! in one loop and then sent to the cohorts in another loop. If hydraulics are
! on, we calculate a unique solution for each level-cohort-layer combination.
! If we are not using hydraulics, we calculate a unique solution for each
! level-pft-layer combination. Thus the following three arrays are statically
! allocated for the maximum space of the two cases (numCohortsPerPatch)
! The "_z" suffix indicates these variables are discretized at the "leaf_layer"
! scale.
! Note: For these temporary arrays, we have the leaf layer dimension first
! and the canopy layer last. This order is chosen for efficiency. The arrays
! such as leaf area that are bound to the patch structure DO NOT follow this order
! as they are used in many other parts of the code with different looping, we
! are not modifying its order now.
! -----------------------------------------------------------------------------------
! leaf maintenance (dark) respiration [umol CO2/m**2/s]
real(r8) :: lmr_z(nlevleaf,maxpft,nclmax)
! stomatal resistance [s/m]
real(r8) :: rs_z(nlevleaf,maxpft,nclmax)
! net leaf photosynthesis averaged over sun and shade leaves. [umol CO2/m**2/s]
real(r8) :: anet_av_z(nlevleaf,maxpft,nclmax)
! Mask used to determine which leaf-layer biophysical rates have been
! used already
logical :: rate_mask_z(nlevleaf,maxpft,nclmax)
real(r8) :: vcmax_z ! leaf layer maximum rate of carboxylation
! (umol co2/m**2/s)
real(r8) :: jmax_z ! leaf layer maximum electron transport rate
! (umol electrons/m**2/s)
real(r8) :: tpu_z ! leaf layer triose phosphate utilization rate
! (umol CO2/m**2/s)
real(r8) :: kp_z ! leaf layer initial slope of CO2 response
! curve (C4 plants)
real(r8) :: c13disc_z(nclmax,maxpft,nlevleaf) ! carbon 13 in newly assimilated carbon at leaf level
real(r8) :: mm_kco2 ! Michaelis-Menten constant for CO2 (Pa)
real(r8) :: mm_ko2 ! Michaelis-Menten constant for O2 (Pa)
real(r8) :: co2_cpoint ! CO2 compensation point (Pa)
real(r8) :: btran_eff ! effective transpiration wetness factor (0 to 1)
real(r8) :: bbb ! Ball-Berry minimum leaf conductance (umol H2O/m**2/s)
real(r8) :: kn ! leaf nitrogen decay coefficient
real(r8) :: cf ! s m**2/umol -> s/m (ideal gas conversion) [umol/m3]
real(r8) :: gb_mol ! leaf boundary layer conductance (molar form: [umol /m**2/s])
real(r8) :: ceair ! vapor pressure of air, constrained (Pa)
real(r8) :: nscaler ! leaf nitrogen scaling coefficient
real(r8) :: leaf_frac ! ratio of to leaf biomass to total alive biomass
real(r8) :: tcsoi ! Temperature response function for root respiration.
real(r8) :: tcwood ! Temperature response function for wood
real(r8) :: elai ! exposed LAI (patch scale)
real(r8) :: live_stem_n ! Live stem (above-ground sapwood)
! nitrogen content (kgN/plant)
real(r8) :: live_croot_n ! Live coarse root (below-ground sapwood)
! nitrogen content (kgN/plant)
real(r8) :: sapw_c ! Sapwood carbon (kgC/plant)
real(r8) :: fnrt_c ! Fine root carbon (kgC/plant)
real(r8) :: fnrt_n ! Fine root nitrogen content (kgN/plant)
real(r8) :: leaf_c ! Leaf carbon (kgC/plant)
real(r8) :: leaf_n ! leaf nitrogen content (kgN/plant)
real(r8) :: g_sb_leaves ! Mean combined (stomata+boundary layer) leaf conductance [m/s]
! over all of the patch's leaves. The "sb" refers to the combined
! "s"tomatal and "b"oundary layer.
! This quantity is relevant on leaf surfaces. It does not
! have units of /m2 leaf per say, but is implicitly on leaf surfaces
real(r8) :: r_sb_leaves ! Mean leaf resistance over all the patch's leaves [s/m]
! This is the direct reciprocal of g_sb_leaves
real(r8) :: r_stomata ! Mean stomatal resistance across all leaves in the patch [s/m]
real(r8) :: maintresp_reduction_factor ! factor by which to reduce maintenance
! respiration when storage pools are low
real(r8) :: b_leaf ! leaf biomass kgC
real(r8) :: frac ! storage pool as a fraction of target leaf biomass
real(r8) :: check_elai ! This is a check on the effective LAI that is calculated
! over each cohort x layer.
real(r8) :: cohort_eleaf_area ! This is the effective leaf area [m2] reported by each cohort
real(r8) :: lnc_top ! Leaf nitrogen content per unit area at canopy top [gN/m2]
real(r8) :: lmr25top ! canopy top leaf maint resp rate at 25C
! for this plant or pft (umol CO2/m**2/s)
real(r8) :: leaf_inc ! LAI-only portion of the vegetation increment of dinc_ed
real(r8) :: lai_canopy_above ! the LAI in the canopy layers above the layer of interest
real(r8) :: lai_layers_above ! the LAI in the leaf layers, within the current canopy,
! above the leaf layer of interest
real(r8) :: lai_current ! the LAI in the current leaf layer
real(r8) :: cumulative_lai ! the cumulative LAI, top down, to the leaf layer of interest
! Junyan's addition
real(r8) :: thl ! leaf water content of the cohort
real(r8) :: thsl ! saturated leaf water content of the cohort
! -----------------------------------------------------------------------------------
! Keeping these two definitions in case they need to be added later
!
! -----------------------------------------------------------------------------------
!real(r8) :: psncanopy_pa ! patch sunlit leaf photosynthesis (umol CO2 /m**2/ s)
!real(r8) :: lmrcanopy_pa ! patch sunlit leaf maintenance respiration rate (umol CO2/m**2/s)
integer :: cl,s,iv,j,ps,ft,ifp ! indices
integer :: nv ! number of leaf layers
integer :: NCL_p ! number of canopy layers in patch
integer :: iage ! loop counter for leaf age classes
! Parameters
! -----------------------------------------------------------------------
! Base maintenance respiration rate for plant tissues base_mr_20
! M. Ryan, 1991. Effects of climate change on plant respiration.
! Ecological Applications, 1(2), 157-167.
! Original expression is br = 0.0106 molC/(molN h)
! Conversion by molecular weights of C and N gives 2.525e-6 gC/(gN s)
!
! Base rate is at 20C. Adjust to 25C using the CN Q10 = 1.5
! (gC/gN/s)
! ------------------------------------------------------------------------
! -----------------------------------------------------------------------------------
! Photosynthesis and stomatal conductance parameters, from:
! Bonan et al (2011) JGR, 116, doi:10.1029/2010JG001593
! -----------------------------------------------------------------------------------
! Ball-Berry minimum leaf conductance, unstressed (umol H2O/m**2/s)
! For C3 and C4 plants
! -----------------------------------------------------------------------------------
real(r8), dimension(0:1) :: bbbopt
associate( &
c3psn => EDPftvarcon_inst%c3psn , &
slatop => EDPftvarcon_inst%slatop , & ! specific leaf area at top of canopy,
! projected area basis [m^2/gC]
woody => EDPftvarcon_inst%woody) ! Is vegetation woody or not?
bbbopt(0) = ED_val_bbopt_c4
bbbopt(1) = ED_val_bbopt_c3
do s = 1,nsites
! Multi-layer parameters scaled by leaf nitrogen profile.
! Loop through each canopy layer to calculate nitrogen profile using
! cumulative lai at the midpoint of the layer
ifp = 0
currentpatch => sites(s)%oldest_patch
do while (associated(currentpatch))
ifp = ifp+1
NCL_p = currentPatch%NCL_p
! Part I. Zero output boundary conditions
! ---------------------------------------------------------------------------
bc_out(s)%rssun_pa(ifp) = 0._r8
bc_out(s)%rssha_pa(ifp) = 0._r8
g_sb_leaves = 0._r8
check_elai = 0._r8
! Part II. Filter out patches
! Patch level filter flag for photosynthesis calculations
! has a short memory, flags:
! 1 = patch has not been called
! 2 = patch is currently marked for photosynthesis
! 3 = patch has been called for photosynthesis already
! ---------------------------------------------------------------------------
if(bc_in(s)%filter_photo_pa(ifp)==2)then
! Part III. Calculate the number of sublayers for each pft and layer.
! And then identify which layer/pft combinations have things in them.
! Output:
! currentPatch%ncan(:,:)
! currentPatch%canopy_mask(:,:)
call UpdateCanopyNCanNRadPresent(currentPatch)
! Part IV. Identify some environmentally derived parameters:
! These quantities are biologically irrelevant
! Michaelis-Menten constant for CO2 (Pa)
! Michaelis-Menten constant for O2 (Pa)
! CO2 compensation point (Pa)
! leaf boundary layer conductance of h20
! constrained vapor pressure
call GetCanopyGasParameters(bc_in(s)%forc_pbot, & ! in
bc_in(s)%oair_pa(ifp), & ! in
bc_in(s)%t_veg_pa(ifp), & ! in
bc_in(s)%tgcm_pa(ifp), & ! in
bc_in(s)%eair_pa(ifp), & ! in
bc_in(s)%esat_tv_pa(ifp), & ! in
bc_in(s)%rb_pa(ifp), & ! in
mm_kco2, & ! out
mm_ko2, & ! out
co2_cpoint, & ! out
cf, & ! out
gb_mol, & ! out
ceair) ! out
! Part V. Pre-process some variables that are PFT dependent
! but not environmentally dependent
! ------------------------------------------------------------------------
do ft = 1,numpft
! This is probably unnecessary and already calculated
! ALSO, THIS ROOTING PROFILE IS USED TO CALCULATE RESPIRATION
! YET IT USES THE PROFILE THAT IS CONSISTENT WITH WATER UPTAKE
! AND NOT THE PROFILE WE USE FOR DECOMPOSITION
! SEEMS LIKE THE LATTER WOULD BE MORE APPROPRIATE, RIGHT? (RGK 05-2018)
call set_root_fraction(currentPatch%rootfr_ft(ft,1:bc_in(s)%nlevsoil), ft, &
bc_in(s)%zi_sisl,lowerb=lbound(bc_in(s)%zi_sisl,1), &
icontext = i_hydro_rootprof_context)
end do !ft
! ------------------------------------------------------------------------
! Part VI: Loop over all leaf layers.
! The concept of leaf layers is a result of the radiative transfer scheme.
! A leaf layer has uniform radiation environment. Leaf layers are a group
! of vegetation surfaces (stems and leaves) which inhabit the same
! canopy-layer "CL", have the same functional type "ft" and within those
! two partitions are further partitioned into vertical layers where
! downwelling radiation attenuates in order.
! In this phase we loop over the leaf layers and calculate the
! photosynthesis and respiration of the layer (since all biophysical
! properties are homogeneous). After this step, we can loop through
! our cohort list, associate each cohort with its list of leaf-layers
! and transfer these quantities to the cohort.
! With plant hydraulics, we must realize that photosynthesis and
! respiration will be different for leaves of each cohort in the leaf
! layers, as they will have there own hydraulic limitations.
! NOTE: Only need to flush mask on the number of used pfts, not the whole
! scratch space.
! ------------------------------------------------------------------------
rate_mask_z(:,1:numpft,:) = .false.
if(currentPatch%countcohorts > 0.0)then ! Ignore empty patches
currentCohort => currentPatch%tallest
do while (associated(currentCohort)) ! Cohort loop
! Identify the canopy layer (cl), functional type (ft)
! and the leaf layer (IV) for this cohort
ft = currentCohort%pft
cl = currentCohort%canopy_layer
call bleaf(currentCohort%dbh,currentCohort%pft,currentCohort%canopy_trim,b_leaf)
call storage_fraction_of_target(b_leaf, &
currentCohort%prt%GetState(store_organ, all_carbon_elements), &
frac)
call lowstorage_maintresp_reduction(frac,currentCohort%pft, &
maintresp_reduction_factor)
! are there any leaves of this pft in this layer?
if(currentPatch%canopy_mask(cl,ft) == 1)then
! Loop over leaf-layers
do iv = 1,currentCohort%nv
! ------------------------------------------------------------
! If we are doing plant hydro-dynamics (or any run-type
! where cohorts may generate different photosynthetic rates
! of other cohorts in the same canopy-pft-layer combo),
! we re-calculate the leaf biophysical rates for the
! cohort-layer combo of interest.
! but in the vanilla case, we only re-calculate if it has
! not been done yet.
! Other cases where we need to solve for every cohort
! in every leaf layer: nutrient dynamic mode, multiple leaf
! age classes
! ------------------------------------------------------------
if ( .not.rate_mask_z(iv,ft,cl) .or. &
(hlm_use_planthydro.eq.itrue) .or. &
(nleafage > 1) .or. &
(hlm_parteh_mode .ne. prt_carbon_allom_hyp ) ) then
if (hlm_use_planthydro.eq.itrue ) then
thl = currentCohort%co_hydr%th_ag(1) !Junyan
thsl = EDPftvarcon_inst%hydr_thetas_node(ft,1) ! Junyan
bbb = max( cf/rsmax0, bbbopt(nint(c3psn(ft)))*currentCohort%co_hydr%btran(1) )
btran_eff = currentCohort%co_hydr%btran(1)
! dinc_ed is the total vegetation area index of each "leaf" layer
! we convert to the leaf only portion of the increment
! ------------------------------------------------------
leaf_inc = dinc_ed * &
currentCohort%treelai/(currentCohort%treelai+currentCohort%treesai)
! Now calculate the cumulative top-down lai of the current layer's midpoint
lai_canopy_above = sum(currentPatch%canopy_layer_tlai(1:cl-1))
lai_layers_above = leaf_inc * (iv-1)
lai_current = min(leaf_inc, currentCohort%treelai - lai_layers_above)
cumulative_lai = lai_canopy_above + lai_layers_above + 0.5*lai_current
else
bbb = max( cf/rsmax0, bbbopt(nint(c3psn(ft)))*currentPatch%btran_ft(ft) )
btran_eff = currentPatch%btran_ft(ft)
! For consistency sake, we use total LAI here, and not exposed
! if the plant is under-snow, it will be effectively dormant for
! the purposes of nscaler
cumulative_lai = sum(currentPatch%canopy_layer_tlai(1:cl-1)) + &
sum(currentPatch%tlai_profile(cl,ft,1:iv-1)) + &
0.5*currentPatch%tlai_profile(cl,ft,iv)
end if
if(do_fates_salinity)then
btran_eff = btran_eff*currentPatch%bstress_sal_ft(ft)
endif
! Bonan et al (2011) JGR, 116, doi:10.1029/2010JG001593 used
! kn = 0.11. Here, derive kn from vcmax25 as in Lloyd et al
! (2010) Biogeosciences, 7, 1833-1859
kn = decay_coeff_kn(ft,currentCohort%vcmax25top)
! Scale for leaf nitrogen profile
nscaler = exp(-kn * cumulative_lai)
! Leaf maintenance respiration to match the base rate used in CN
! but with the new temperature functions for C3 and C4 plants.
! CN respiration has units: g C / g N [leaf] / s. This needs to be
! converted from g C / g N [leaf] / s to umol CO2 / m**2 [leaf] / s
! Then scale this value at the top of the canopy for canopy depth
! Leaf nitrogen concentration at the top of the canopy (g N leaf / m**2 leaf)
select case(hlm_parteh_mode)
case (prt_carbon_allom_hyp)
lnc_top = EDPftvarcon_inst%prt_nitr_stoich_p1(ft,leaf_organ)/slatop(ft)
case (prt_cnp_flex_allom_hyp)
leaf_c = currentCohort%prt%GetState(leaf_organ, all_carbon_elements)
leaf_n = currentCohort%prt%GetState(leaf_organ, all_carbon_elements)
lnc_top = leaf_n / (slatop(ft) * leaf_c )
end select
lmr25top = 2.525e-6_r8 * (1.5_r8 ** ((25._r8 - 20._r8)/10._r8))
lmr25top = lmr25top * lnc_top / (umolC_to_kgC * g_per_kg)
! Part VII: Calculate dark respiration (leaf maintenance) for this layer
call LeafLayerMaintenanceRespiration( lmr25top, & ! in
nscaler, & ! in
ft, & ! in
bc_in(s)%t_veg_pa(ifp), & ! in
lmr_z(iv,ft,cl)) ! out
! Part VII: Calculate (1) maximum rate of carboxylation (vcmax),
! (2) maximum electron transport rate, (3) triose phosphate
! utilization rate and (4) the initial slope of CO2 response curve
! (C4 plants). Earlier we calculated their base rates as dictated
! by their plant functional type and some simple scaling rules for
! nitrogen limitation baesd on canopy position (not prognostic).
! These rates are the specific rates used in the actual photosynthesis
! calculations that take localized environmental effects (temperature)
! into consideration.
call LeafLayerBiophysicalRates(currentPatch%ed_parsun_z(cl,ft,iv), & ! in
ft, & ! in
currentCohort%vcmax25top, & ! in
currentCohort%jmax25top, & ! in
currentCohort%tpu25top, & ! in
currentCohort%kp25top, & ! in
nscaler, & ! in
bc_in(s)%t_veg_pa(ifp), & ! in
btran_eff, & ! in
vcmax_z, & ! out
jmax_z, & ! out
tpu_z, & ! out
kp_z ) ! out
! Part IX: This call calculates the actual photosynthesis for the
! leaf layer, as well as the stomatal resistance and the net assimilated carbon.
call LeafLayerPhotosynthesis(currentPatch%f_sun(cl,ft,iv), & ! in
currentPatch%ed_parsun_z(cl,ft,iv), & ! in
currentPatch%ed_parsha_z(cl,ft,iv), & ! in
currentPatch%ed_laisun_z(cl,ft,iv), & ! in
currentPatch%ed_laisha_z(cl,ft,iv), & ! in
currentPatch%canopy_area_profile(cl,ft,iv), & ! in
ft, & ! in
vcmax_z, & ! in
jmax_z, & ! in
tpu_z, & ! in
kp_z, & ! in
bc_in(s)%t_veg_pa(ifp), & ! in
bc_in(s)%esat_tv_pa(ifp), & ! in
bc_in(s)%forc_pbot, & ! in
bc_in(s)%cair_pa(ifp), & ! in
bc_in(s)%oair_pa(ifp), & ! in
btran_eff, & ! in
bbb, & ! in
cf, & ! in
gb_mol, & ! in
ceair, & ! in
mm_kco2, & ! in
mm_ko2, & ! in
co2_cpoint, & ! in
lmr_z(iv,ft,cl), & ! in
thl, & ! in Junyan
thsl, & ! in Junyan
currentPatch%psn_z(cl,ft,iv), & ! out
rs_z(iv,ft,cl), & ! out
anet_av_z(iv,ft,cl), & ! out
c13disc_z(cl,ft,iv)) ! out
rate_mask_z(iv,ft,cl) = .true.
end if
end do
! Zero cohort flux accumulators.
currentCohort%npp_tstep = 0.0_r8
currentCohort%resp_tstep = 0.0_r8
currentCohort%gpp_tstep = 0.0_r8
currentCohort%rdark = 0.0_r8
currentCohort%resp_m = 0.0_r8
currentCohort%ts_net_uptake = 0.0_r8
currentCohort%c13disc_clm = 0.0_r8
! ---------------------------------------------------------------
! Part VII: Transfer leaf flux rates (like maintenance respiration,
! carbon assimilation and conductance) that are defined by the
! leaf layer (which is area independent, ie /m2) onto each cohort
! (where the rates become per cohort, ie /individual). Most likely
! a sum over layers.
! ---------------------------------------------------------------
nv = currentCohort%nv
call ScaleLeafLayerFluxToCohort(nv, & !in
currentPatch%psn_z(cl,ft,1:nv), & !in
lmr_z(1:nv,ft,cl), & !in
rs_z(1:nv,ft,cl), & !in
currentPatch%elai_profile(cl,ft,1:nv), & !in
c13disc_z(cl, ft, 1:nv), & !in
currentCohort%c_area, & !in
currentCohort%n, & !in
bc_in(s)%rb_pa(ifp), & !in
maintresp_reduction_factor, & !in
currentCohort%g_sb_laweight, & !out
currentCohort%gpp_tstep, & !out
currentCohort%rdark, & !out
currentCohort%c13disc_clm, & !out
cohort_eleaf_area) !out
! Net Uptake does not need to be scaled, just transfer directly
currentCohort%ts_net_uptake(1:nv) = anet_av_z(1:nv,ft,cl) * umolC_to_kgC
else
! In this case, the cohort had no leaves,
! so no productivity,conductance, transpiration uptake
! or dark respiration
cohort_eleaf_area = 0.0_r8
currentCohort%gpp_tstep = 0.0_r8
currentCohort%rdark = 0.0_r8
currentCohort%g_sb_laweight = 0.0_r8
currentCohort%ts_net_uptake(:) = 0.0_r8
end if ! if(currentPatch%canopy_mask(cl,ft) == 1)then
! ------------------------------------------------------------------
! Part VIII: Calculate maintenance respiration in the sapwood and
! fine root pools.
! ------------------------------------------------------------------
! Calculate the amount of nitrogen in the above and below ground
! stem and root pools, used for maint resp
! We are using the fine-root C:N ratio as an approximation for
! the sapwood pools.
! Units are in (kgN/plant)
! ------------------------------------------------------------------
sapw_c = currentCohort%prt%GetState(sapw_organ, all_carbon_elements)
fnrt_c = currentCohort%prt%GetState(fnrt_organ, all_carbon_elements)
select case(hlm_parteh_mode)
case (prt_carbon_allom_hyp)
live_stem_n = EDPftvarcon_inst%allom_agb_frac(currentCohort%pft) * &
sapw_c * EDPftvarcon_inst%prt_nitr_stoich_p1(ft,sapw_organ)
live_croot_n = (1.0_r8-EDPftvarcon_inst%allom_agb_frac(currentCohort%pft)) * &
sapw_c * EDPftvarcon_inst%prt_nitr_stoich_p1(ft,sapw_organ)
fnrt_n = fnrt_c * EDPftvarcon_inst%prt_nitr_stoich_p1(ft,fnrt_organ)
case(prt_cnp_flex_allom_hyp)
live_stem_n = EDPftvarcon_inst%allom_agb_frac(currentCohort%pft) * &
currentCohort%prt%GetState(sapw_organ, nitrogen_element)
live_croot_n = (1.0_r8-EDPftvarcon_inst%allom_agb_frac(currentCohort%pft)) * &
currentCohort%prt%GetState(sapw_organ, nitrogen_element)
fnrt_n = currentCohort%prt%GetState(fnrt_organ, nitrogen_element)
case default
end select
!------------------------------------------------------------------------------
! Calculate Whole Plant Respiration
! (this doesn't really need to be in this iteration at all, surely?)
! Response: (RGK 12-2016): I think the positioning of these calls is
! appropriate as of now. Maintenance calculations in sapwood and roots
! vary by cohort and with changing temperature at the minimum, and there are
! no sub-pools chopping up those pools any finer that need to be dealt with.
!------------------------------------------------------------------------------
! Live stem MR (kgC/plant/s) (above ground sapwood)
! ------------------------------------------------------------------
if (woody(ft) == 1) then
tcwood = q10_mr**((bc_in(s)%t_veg_pa(ifp)-tfrz - 20.0_r8)/10.0_r8)
! kgC/s = kgN * kgC/kgN/s
currentCohort%livestem_mr = live_stem_n * ED_val_base_mr_20 * tcwood * maintresp_reduction_factor
else
currentCohort%livestem_mr = 0._r8
end if
! Fine Root MR (kgC/plant/s)
! ------------------------------------------------------------------
currentCohort%froot_mr = 0._r8
do j = 1,bc_in(s)%nlevsoil
tcsoi = q10_mr**((bc_in(s)%t_soisno_sl(j)-tfrz - 20.0_r8)/10.0_r8)
currentCohort%froot_mr = currentCohort%froot_mr + &
fnrt_n * ED_val_base_mr_20 * tcsoi * currentPatch%rootfr_ft(ft,j) * maintresp_reduction_factor
enddo
! Coarse Root MR (kgC/plant/s) (below ground sapwood)
! ------------------------------------------------------------------
if (woody(ft) == 1) then
currentCohort%livecroot_mr = 0._r8
do j = 1,bc_in(s)%nlevsoil
! Soil temperature used to adjust base rate of MR
tcsoi = q10_mr**((bc_in(s)%t_soisno_sl(j)-tfrz - 20.0_r8)/10.0_r8)
currentCohort%livecroot_mr = currentCohort%livecroot_mr + &
live_croot_n * ED_val_base_mr_20 * tcsoi * &
currentPatch%rootfr_ft(ft,j) * maintresp_reduction_factor
enddo
else
currentCohort%livecroot_mr = 0._r8
end if
! ------------------------------------------------------------------
! Part IX: Perform some unit conversions (rate to integrated) and
! calcualate some fluxes that are sums and nets of the base fluxes
! ------------------------------------------------------------------
if ( debug ) write(fates_log(),*) 'EDPhoto 904 ', currentCohort%resp_m
if ( debug ) write(fates_log(),*) 'EDPhoto 905 ', currentCohort%rdark
if ( debug ) write(fates_log(),*) 'EDPhoto 906 ', currentCohort%livestem_mr
if ( debug ) write(fates_log(),*) 'EDPhoto 907 ', currentCohort%livecroot_mr
if ( debug ) write(fates_log(),*) 'EDPhoto 908 ', currentCohort%froot_mr
! add on whole plant respiration values in kgC/indiv/s-1
currentCohort%resp_m = currentCohort%livestem_mr + &
currentCohort%livecroot_mr + &
currentCohort%froot_mr
! no drought response right now.. something like:
! resp_m = resp_m * (1.0_r8 - currentPatch%btran_ft(currentCohort%pft) * &
! EDPftvarcon_inst%resp_drought_response(ft))
currentCohort%resp_m = currentCohort%resp_m + currentCohort%rdark
! convert from kgC/indiv/s to kgC/indiv/timestep
currentCohort%resp_m = currentCohort%resp_m * dtime
currentCohort%gpp_tstep = currentCohort%gpp_tstep * dtime
currentCohort%ts_net_uptake = currentCohort%ts_net_uptake * dtime
if ( debug ) write(fates_log(),*) 'EDPhoto 911 ', currentCohort%gpp_tstep
if ( debug ) write(fates_log(),*) 'EDPhoto 912 ', currentCohort%resp_tstep
if ( debug ) write(fates_log(),*) 'EDPhoto 913 ', currentCohort%resp_m
currentCohort%resp_g = EDPftvarcon_inst%grperc(ft) * &
(max(0._r8,currentCohort%gpp_tstep - &
currentCohort%resp_m))
currentCohort%resp_tstep = currentCohort%resp_m + &
currentCohort%resp_g ! kgC/indiv/ts
currentCohort%npp_tstep = currentCohort%gpp_tstep - &
currentCohort%resp_tstep ! kgC/indiv/ts
! Accumulate the combined conductance (stomatal+leaf boundary layer)
! Note that currentCohort%g_sb_laweight is weighted by the leaf area
! of each cohort and has units of [m/s] * [m2 leaf]
g_sb_leaves = g_sb_leaves + currentCohort%g_sb_laweight
! Accumulate the total effective leaf area from all cohorts
! in this patch. Normalize by canopy area outside the loop
check_elai = check_elai + cohort_eleaf_area
currentCohort => currentCohort%shorter
enddo ! end cohort loop.
end if !count_cohorts is more than zero.
check_elai = check_elai / currentPatch%total_canopy_area
elai = calc_areaindex(currentPatch,'elai')
! Normalize canopy total conductance by the effective LAI
! The value here was integrated over each cohort x leaf layer
! and was weighted by m2 of effective leaf area for each layer
if(check_elai>tiny(check_elai)) then
! Normalize the leaf-area weighted canopy conductance
! The denominator is the total effective leaf area in the canopy,
! units of [m/s]*[m2] / [m2] = [m/s]
g_sb_leaves = g_sb_leaves / (elai*currentPatch%total_canopy_area)
if( g_sb_leaves > (1._r8/rsmax0) ) then
! Combined mean leaf resistance is the inverse of mean leaf conductance
r_sb_leaves = 1.0_r8/g_sb_leaves
if (r_sb_leaves<bc_in(s)%rb_pa(ifp)) then
write(fates_log(),*) 'Combined canopy resistance was somehow smaller than'
write(fates_log(),*) 'its boundary layer resistance component'
write(fates_log(),*) 'r_sb_leaves [s/m]: ',r_sb_leaves
write(fates_log(),*) 'bc_in(s)%rb_pa(ifp) [s/m]: ',bc_in(s)%rb_pa(ifp)
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
! Mean leaf stomatal resistance for all patch leaves
r_stomata = (r_sb_leaves - bc_in(s)%rb_pa(ifp))
else
! Here we prevent super high resistances
! and use a nominal value when conductance is low
! Junyan changed that to remove the cap when hydro is on June 25, 2019
if (hlm_use_planthydro.eq.ifalse ) then
r_stomata = rsmax0
end if
end if
! This will be multiplied by scaled by effective LAI in the host model
! when it comes time to calculate a flux rate per unit ground
bc_out(s)%rssun_pa(ifp) = r_stomata
bc_out(s)%rssha_pa(ifp) = r_stomata
! This value is used for diagnostics, the molar form of conductance
! is what is used in the field usually, so we track that form
currentPatch%c_stomata = cf / r_stomata
else
! But this will prevent it from using an unintialized value
bc_out(s)%rssun_pa(ifp) = rsmax0
bc_out(s)%rssha_pa(ifp) = rsmax0
! This value is used for diagnostics, the molar form of conductance
! is what is used in the field usually, so we track that form
currentPatch%c_stomata = cf / rsmax0
end if
! This value is used for diagnostics, the molar form of conductance
! is what is used in the field usually, so we track that form
currentPatch%c_lblayer = cf / bc_in(s)%rb_pa(ifp)
end if
currentPatch => currentPatch%younger
end do
end do !site loop
end associate
end subroutine FatesPlantRespPhotosynthDrive
! =======================================================================================
subroutine LeafLayerPhotosynthesis(f_sun_lsl, & ! in
parsun_lsl, & ! in
parsha_lsl, & ! in
laisun_lsl, & ! in
laisha_lsl, & ! in
canopy_area_lsl, & ! in
ft, & ! in
vcmax, & ! in
jmax, & ! in
tpu, & ! in
co2_rcurve_islope, & ! in
veg_tempk, & ! in
veg_esat, & ! in
can_press, & ! in
can_co2_ppress, & ! in
can_o2_ppress, & ! in
btran, & ! in
bbb, & ! in
cf, & ! in
gb_mol, & ! in
ceair, & ! in
mm_kco2, & ! in
mm_ko2, & ! in
co2_cpoint, & ! in
lmr, & ! in
thl, & ! in Junyan
thsl, & ! in Junyan
psn_out, & ! out
rstoma_out, & ! out
anet_av_out, & ! out
c13disc_z) ! out
! ------------------------------------------------------------------------------------
! This subroutine calculates photosynthesis and stomatal conductance within each leaf
! sublayer.
! A note on naming conventions: As this subroutine is called for every
! leaf-sublayer, many of the arguments are specific to that "leaf sub layer"
! (LSL), those variables are given a dimension tag "_lsl"
! Other arguments or variables may be indicative of scales broader than the LSL.
! Junyan made modifications on gstomaw to deal with over estimate on transpiration by HLM
! for FATESHydrau
! ------------------------------------------------------------------------------------
use EDPftvarcon , only : EDPftvarcon_inst
! Arguments
! ------------------------------------------------------------------------------------
real(r8), intent(in) :: f_sun_lsl !
real(r8), intent(in) :: parsun_lsl ! Absorbed PAR in sunlist leaves
real(r8), intent(in) :: parsha_lsl ! Absorved PAR in shaded leaves
real(r8), intent(in) :: laisun_lsl ! LAI in sunlit leaves
real(r8), intent(in) :: laisha_lsl ! LAI in shaded leaves
real(r8), intent(in) :: canopy_area_lsl !
integer, intent(in) :: ft ! (plant) Functional Type Index
real(r8), intent(in) :: vcmax ! maximum rate of carboxylation (umol co2/m**2/s)
real(r8), intent(in) :: jmax ! maximum electron transport rate (umol electrons/m**2/s)
real(r8), intent(in) :: tpu ! triose phosphate utilization rate (umol CO2/m**2/s)
real(r8), intent(in) :: co2_rcurve_islope ! initial slope of CO2 response curve (C4 plants)
real(r8), intent(in) :: veg_tempk ! vegetation temperature
real(r8), intent(in) :: veg_esat ! saturation vapor pressure at veg_tempk (Pa)
! Important Note on the following gas pressures. This photosynthesis scheme will iteratively
! solve for the co2 partial pressure at the leaf surface (ie in the stomata). The reference
! point for these input values are NOT within that boundary layer that separates the stomata from
! the canopy air space. The reference point for these is on the outside of that boundary
! layer. This routine, which operates at the leaf scale, makes no assumptions about what the
! scale of the refernce is, it could be lower atmosphere, it could be within the canopy
! but most likely it is the closest value one can get to the edge of the leaf's boundary
! layer. We use the convention "can_" because a reference point of within the canopy
! ia a best reasonable scenario of where we can get that information from.
real(r8), intent(in) :: can_press ! Air pressure NEAR the surface of the leaf (Pa)
real(r8), intent(in) :: can_co2_ppress ! Partial pressure of CO2 NEAR the leaf surface (Pa)
real(r8), intent(in) :: can_o2_ppress ! Partial pressure of O2 NEAR the leaf surface (Pa)
real(r8), intent(in) :: btran ! transpiration wetness factor (0 to 1)
real(r8), intent(in) :: bbb ! Ball-Berry minimum leaf conductance (umol H2O/m**2/s)
real(r8), intent(in) :: cf ! s m**2/umol -> s/m (ideal gas conversion) [umol/m3]
real(r8), intent(in) :: gb_mol ! leaf boundary layer conductance (umol /m**2/s)
real(r8), intent(in) :: ceair ! vapor pressure of air, constrained (Pa)
real(r8), intent(in) :: mm_kco2 ! Michaelis-Menten constant for CO2 (Pa)
real(r8), intent(in) :: mm_ko2 ! Michaelis-Menten constant for O2 (Pa)
real(r8), intent(in) :: co2_cpoint ! CO2 compensation point (Pa)
real(r8), intent(in) :: lmr ! Leaf Maintenance Respiration (umol CO2/m**2/s)
real(r8), intent(out) :: psn_out ! carbon assimilated in this leaf layer umolC/m2/s
real(r8), intent(out) :: rstoma_out ! stomatal resistance (1/gs_lsl) (s/m)
real(r8), intent(out) :: anet_av_out ! net leaf photosynthesis (umol CO2/m**2/s)
! averaged over sun and shade leaves.
real(r8), intent(out) :: c13disc_z ! carbon 13 in newly assimilated carbon
! Locals
! ------------------------------------------------------------------------
integer :: c3c4_path_index ! Index for which photosynthetic pathway
! is active. C4 = 0, C3 = 1
integer :: sunsha ! Index for differentiating sun and shade
real(r8) :: gstoma ! Stomatal Conductance of this leaf layer (m/s)
real(r8) :: agross ! co-limited gross leaf photosynthesis (umol CO2/m**2/s)
real(r8) :: anet ! net leaf photosynthesis (umol CO2/m**2/s)
real(r8) :: je ! electron transport rate (umol electrons/m**2/s)
real(r8) :: qabs ! PAR absorbed by PS II (umol photons/m**2/s)
real(r8) :: aquad,bquad,cquad,cquadw ! terms for quadratic equations, cquadw is the same as cquad but only used when Hydrau is on
real(r8) :: r1,r2,rw1,rw2 ! roots of quadratic equation
real(r8) :: co2_inter_c ! intercellular leaf CO2 (Pa)
real(r8) :: co2_inter_c_old ! intercellular leaf CO2 (Pa) (previous iteration)
logical :: loop_continue ! Loop control variable
integer :: niter ! iteration loop index
real(r8) :: gs_mol ! leaf stomatal conductance (umol H2O/m**2/s)
real(r8) :: gs ! leaf stomatal conductance (m/s)
real(r8) :: hs ! fractional humidity at leaf surface (dimensionless)
real(r8) :: gs_mol_err ! gs_mol for error check
real(r8) :: ac ! Rubisco-limited gross photosynthesis (umol CO2/m**2/s)
real(r8) :: aj ! RuBP-limited gross photosynthesis (umol CO2/m**2/s)
real(r8) :: ap ! product-limited (C3) or CO2-limited
! (C4) gross photosynthesis (umol CO2/m**2/s)
real(r8) :: ai ! intermediate co-limited photosynthesis (umol CO2/m**2/s)
real(r8) :: leaf_co2_ppress ! CO2 partial pressure at leaf surface (Pa)
real(r8) :: init_co2_inter_c ! First guess intercellular co2 specific to C path
real(r8), dimension(0:1) :: bbbopt ! Cuticular conductance at full water potential (umol H2O /m2/s)
real(r8) :: gsw
real(r8) :: gsw_mol
real(r8) :: thl
real(r8) :: thsl
! Parameters
! ------------------------------------------------------------------------
! Fraction of light absorbed by non-photosynthetic pigments
real(r8),parameter :: fnps = 0.15_r8
! For plants with no leaves, a miniscule amount of conductance
! can happen through the stems, at a partial rate of cuticular conductance
real(r8),parameter :: stem_cuticle_loss_frac = 0.1_r8
! empirical curvature parameter for electron transport rate
real(r8),parameter :: theta_psii = 0.7_r8
! First guess on ratio between intercellular co2 and the atmosphere
! an iterator converges on actual
real(r8),parameter :: init_a2l_co2_c3 = 0.7_r8
real(r8),parameter :: init_a2l_co2_c4 = 0.4_r8
! quantum efficiency, used only for C4 (mol CO2 / mol photons) (index 0)
real(r8),parameter,dimension(0:1) :: quant_eff = [0.05_r8,0.0_r8]
! empirical curvature parameter for ac, aj photosynthesis co-limitation.
! Changed theta_cj and theta_ip to 0.999 to effectively remove smoothing logic
! following Anthony Walker's findings from MAAT.
real(r8),parameter,dimension(0:1) :: theta_cj = [0.999_r8,0.999_r8]
! empirical curvature parameter for ap photosynthesis co-limitation
real(r8),parameter :: theta_ip = 0.999_r8
associate( bb_slope => EDPftvarcon_inst%BB_slope) ! slope of BB relationship
! photosynthetic pathway: 0. = c4, 1. = c3
c3c4_path_index = nint(EDPftvarcon_inst%c3psn(ft))
bbbopt(0) = ED_val_bbopt_c4
bbbopt(1) = ED_val_bbopt_c3
if (c3c4_path_index == 1) then
init_co2_inter_c = init_a2l_co2_c3 * can_co2_ppress
else
init_co2_inter_c = init_a2l_co2_c4 * can_co2_ppress
end if
! Part III: Photosynthesis and Conductance
! ----------------------------------------------------------------------------------
if ( parsun_lsl <= 0._r8 ) then ! night time
anet_av_out = -lmr
psn_out = 0._r8
! The cuticular conductance already factored in maximum resistance as a bound
! no need to re-bound it