-
Notifications
You must be signed in to change notification settings - Fork 92
/
FatesPlantHydraulicsMod.F90
6363 lines (5049 loc) · 275 KB
/
FatesPlantHydraulicsMod.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 FatesPlantHydraulicsMod
! ==============================================================================================
! This module contains the relevant code for plant hydraulics. Currently, one hydraulics module
! is available. Other methods of estimating plant hydraulics may become available in future
! releases. For now, please cite the following reference if this module is used to generate
! published research:
!
! Christoffersen, B.O., Gloor, M., Fauset, S., Fyllas, N. M., Galbraith, D. R., Baker,
! T. R., Kruijt, B., Rowland, L., Fisher, R. A., Binks, O. J., Sevanto, S., Xu, C., Jansen,
! S., Choat, B., Mencuccini, M., McDowell, N. G., Meir, P. Linking hydraulic traits to
! tropical forest function in a size-structured and trait-driven model (TFS~v.1-Hydro).
! Geoscientific Model Development, 9(11), 2016, pp: 4227-4255,
! https://www.geosci-model-dev.net/9/4227/2016/, DOI = 10.5194/gmd-9-4227-2016.
!
! WARNING WARNING WARNING WARNING WARNING WARNING WARNING WARNING WARNING WARNING WARNING
!
! PLANT HYDRAULICS IS AN EXPERIMENTAL OPTION THAT IS STILL UNDERGOING TESTING.
!
! WARNING WARNING WARNING WARNING WARNING WARNING WARNING WARNING WARNING WARNING WARNING
!
! ==============================================================================================
use FatesGlobals, only : endrun => fates_endrun
use FatesGlobals, only : fates_log
use FatesConstantsMod, only : r8 => fates_r8
use FatesConstantsMod, only : fates_huge
use FatesConstantsMod, only : denh2o => dens_fresh_liquid_water
use FatesConstantsMod, only : grav_earth
use FatesConstantsMod, only : ifalse, itrue
use FatesConstantsMod, only : pi_const
use FatesConstantsMod, only : cm2_per_m2
use FatesConstantsMod, only : g_per_kg
use FatesConstantsMod, only : nearzero
use FatesConstantsMod, only : mpa_per_pa
use FatesConstantsMod, only : m_per_mm
use FatesConstantsMod, only : mg_per_kg
use FatesConstantsMod, only : pa_per_mpa
use FatesConstantsMod, only : rsnbl_math_prec
use FatesConstantsMod, only : m3_per_mm3
use FatesConstantsMod, only : cm3_per_m3
use FatesConstantsMod, only : kg_per_g
use FatesConstantsMod, only : fates_unset_r8
use FatesConstantsMod, only : nocomp_bareground
use EDParamsMod , only : hydr_kmax_rsurf1
use EDParamsMod , only : hydr_kmax_rsurf2
use EDParamsMod , only : hydr_psi0
use EDParamsMod , only : hydr_psicap
use EDParamsMod , only : hydr_htftype_node
use EDParamsMod , only : hydr_solver
use EDTypesMod , only : ed_site_type
use FatesPatchMod , only : fates_patch_type
use FatesCohortMod , only : fates_cohort_type
use EDTypesMod , only : AREA_INV
use EDTypesMod , only : AREA
use FatesConstantsMod , only : leaves_on
use FatesInterfaceTypesMod , only : bc_in_type
use FatesInterfaceTypesMod , only : bc_out_type
use FatesInterfaceTypesMod , only : hlm_use_planthydro
use FatesInterfaceTypesMod , only : hlm_ipedof
use FatesInterfaceTypesMod , only : numpft
use FatesInterfaceTypesMod , only : nlevsclass
use FatesAllometryMod, only : bleaf
use FatesAllometryMod, only : bsap_allom
use FatesAllometryMod, only : CrownDepth
use FatesHydraulicsMemMod, only: hydr_solver_1DTaylor
use FatesHydraulicsMemMod, only: hydr_solver_2DNewton
use FatesHydraulicsMemMod, only: hydr_solver_2DPicard
use FatesHydraulicsMemMod, only: ed_site_hydr_type
use FatesHydraulicsMemMod, only: ed_cohort_hydr_type
use FatesHydraulicsMemMod, only: n_hypool_plant
use FatesHydraulicsMemMod, only: n_hypool_leaf
use FatesHydraulicsMemMod, only: n_hypool_tot
use FatesHydraulicsMemMod, only: n_hypool_stem
use FatesHydraulicsMemMod, only: n_hypool_troot
use FatesHydraulicsMemMod, only: n_hypool_aroot
use FatesHydraulicsMemMod, only: n_plant_media
use FatesHydraulicsMemMod, only: nshell
use FatesHydraulicsMemMod, only: n_hypool_ag
use FatesHydraulicsMemMod, only: stomata_p_media
use FatesHydraulicsMemMod, only: leaf_p_media
use FatesHydraulicsMemMod, only: stem_p_media
use FatesHydraulicsMemMod, only: troot_p_media
use FatesHydraulicsMemMod, only: aroot_p_media
use FatesHydraulicsMemMod, only: rhiz_p_media
use FatesHydraulicsMemMod, only: nlevsoi_hyd_max
use FatesHydraulicsMemMod, only: rwccap, rwcft
use PRTGenericMod, only : carbon12_element
use PRTGenericMod, only : leaf_organ, fnrt_organ, sapw_organ
use PRTGenericMod, only : store_organ, repro_organ, struct_organ
use PRTGenericMod, only : num_elements
use PRTGenericMod, only : element_list
use EDPftvarcon, only : EDPftvarcon_inst
use PRTParametersMod, only : prt_params
use FatesHydroWTFMod, only : wrf_arr_type
use FatesHydroWTFMod, only : wkf_arr_type
use FatesHydroWTFMod, only : wrf_type, wrf_type_vg, wrf_type_cch, wrf_type_tfs
use FatesHydroWTFMod, only : wkf_type, wkf_type_vg, wkf_type_cch, wkf_type_tfs
use FatesHydroWTFMod, only : wrf_type_smooth_cch, wkf_type_smooth_cch
! CIME Globals
use shr_log_mod , only : errMsg => shr_log_errMsg
use shr_infnan_mod , only : isnan => shr_infnan_isnan
implicit none
! 1=leaf, 2=stem, 3=troot, 4=aroot
! Several of these may be better transferred to the parameter file in due time (RGK)
integer, public :: use_ed_planthydraulics = 1 ! 0 => use vanilla btran
! 1 => use BC hydraulics;
! 2 => use CX hydraulics
! The following options are temporarily unavailable (RGK 09-06-19)
! ----------------------------------------------------------------------------------
! logical, public :: do_dqtopdth_leaf = .false. ! should a nonzero dqtopdth_leaf
! term be applied to the plant
! hydraulics numerical solution?
! logical, public :: do_dyn_xylemrefill = .false. ! should the dynamics of xylem refilling
! (i.e., non-instantaneous) be considered
! within plant hydraulics?
! logical, public :: do_kbound_upstream = .true. ! should the hydraulic conductance at the
! boundary between nodes be taken to be a
! function of the upstream loss of
! conductivity (flc)?
! DO NOT TURN THIS ON. LEAVING THIS ONLY IF THE HLMS START HAVING
! TROUBLE RESPONDING TO SUPERSATURATION
logical :: purge_supersaturation = .false. ! If for some reason the roots force water
! into a saturated soil layer, or push it slightly
! past saturation, should we attempt to help
! fix the situation by assigning some
! of the water to a runoff term?
logical, public :: do_growthrecruiteffects = .true. ! should size- or root length-dependent
! hydraulic properties and states be
! updated every day when trees grow or
! when recruitment happens?
! If this is set to true, then the conductance over a path between nodes, is defined
! by the side of the path with higher potential only.
logical, parameter :: do_upstream_k = .true.
logical :: do_parallel_stem = .true. ! If this mode is active, we treat the conduit through
! the plant (in 1D solves) as closed from root layer
! to the stomata. The effect of this, is that
! the conductances through stem and leaf surfaces
! are reduced by the fraction of active root
! conductance, and for each soil-layer, integration
! proceeds over the entire time-step.
! These switches are for developers who which to understand if there simulations
! are ever entering regimes where water contents go negative (yes physically impossible)
! or water pressures exceed that at saturation (maybe, maybe not likely)
! These situations are possible/likely due to the nature of the constant flux boundary condition
! of transpiration, due to the loosely-coupled nature of the hydro-land-energy-photosynthesis
! system
logical, parameter :: trap_neg_wc = .false.
logical, parameter :: trap_supersat_psi = .false.
real(r8), parameter :: error_thresh = 1.e-5_r8 ! site level conservation error threshold in CLM
! (mm = kg/m2)
real(r8), parameter :: thsat_buff = 0.001_r8 ! Ensure that this amount of buffer
! is left between soil moisture and saturation [m3/m3]
! (if we are going to help purge super-saturation)
logical,parameter :: debug = .false. ! flag to report warning in hydro
character(len=*), parameter, private :: sourcefile = &
__FILE__
! These index flags specify which pressure-volumen and pressure
! conductivity relationship are available.
! For plants: Users can option between useing tfs and van_genuchten
! by specifying their choice in the parameter file,
! with the model parameter hydr_htftype_node,
! the value should be 1 for TFS or 2 for VG (as shown below).
! Campbell, could technically be used, but the parameters for
! that hypothesis are not in the parameter file, so it not currently available.
! For soil: The soil hypothesis should follow the hypothesis for water transfer
! in the Host Land Model. At this time campbell is the default for both
! ELM and ALM. However, if alternatives arise (like VG), we still need to write
! interface routines to transfer over parameters. Right now we just hard-code
! the use of campbell_type for the soil (see a few lines below).
integer, public, parameter :: van_genuchten_type = 2
integer, public, parameter :: campbell_type = 3
integer, public, parameter :: smooth1_campbell_type = 31
integer, public, parameter :: smooth2_campbell_type = 32
integer, public, parameter :: tfs_type = 1
integer, parameter :: soil_wrf_type = campbell_type
integer, parameter :: soil_wkf_type = campbell_type
! Define the global object that holds the water retention functions
! for plants of each different porous media type, and plant functional type
class(wrf_arr_type),pointer :: wrf_plant(:,:)
! Define the global object that holds the water conductance functions
! for plants of each different porous media type, and plant functional type
class(wkf_arr_type), pointer :: wkf_plant(:,:)
! Testing parameters for Van Genuchten soil WRTs
! unused unless van_genuchten_type is selected, also
! it would be much better to use the native parameters passed in
! from the HLM's soil model
real(r8), parameter :: alpha_vg = 0.001_r8
real(r8), parameter :: th_sat_vg = 0.65_r8
real(r8), parameter :: th_res_vg = 0.15_r8
real(r8), parameter :: psd_vg = 2.7_r8
real(r8), parameter :: m_vg = 0.62963_r8
real(r8), parameter :: soil_tort_vg = 0.5_r8
real(r8), parameter :: plant_tort_vg = 0.0_r8
! The maximum allowable water balance error over a plant-soil continuum
! for a given step [kgs] (2 mg)
real(r8), parameter :: max_wb_step_err = 2.e-6_r8 ! original is 1.e-7_r8, Junyan changed to 2.e-6_r8
!
! !PUBLIC MEMBER FUNCTIONS:
public :: AccumulateMortalityWaterStorage
public :: RecruitWaterStorage
public :: hydraulics_drive
public :: InitHydrSites
public :: HydrSiteColdStart
public :: BTranForHLMDiagnosticsFromCohortHydr
public :: InitHydrCohort
public :: DeallocateHydrCohort
public :: UpdateH2OVeg
public :: FuseCohortHydraulics
public :: UpdateSizeDepPlantHydProps
public :: UpdateSizeDepPlantHydStates
public :: UpdatePlantPsiFTCFromTheta
public :: InitPlantHydStates
public :: UpdateSizeDepRhizHydProps
public :: RestartHydrStates
public :: SavePreviousCompartmentVolumes
public :: SavePreviousRhizVolumes
public :: UpdatePlantHydrNodes
public :: UpdatePlantHydrLenVol
public :: UpdatePlantKmax
public :: ConstrainRecruitNumber
public :: InitHydroGlobals
! RGK 12-2021: UpdateSizeDepRhizHydStates was removed
! this code can be found in tags prior to
! sci.1.52.0_api.20.0.0
!------------------------------------------------------------------------------
! 01/18/16: Created by Brad Christoffersen
! 02/xx/17: Refactoring by Ryan Knox and Brad Christoffersen
!------------------------------------------------------------------------------
contains
!------------------------------------------------------------------------------
subroutine hydraulics_drive( nsites, sites, bc_in,bc_out,dtime )
! 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
select case (use_ed_planthydraulics)
case (1)
call FillDrainRhizShells(nsites, sites, bc_in, bc_out )
call hydraulics_BC(nsites, sites,bc_in,bc_out,dtime )
case (2)
!call Hydraulics_CX()
case DEFAULT
end select
end subroutine Hydraulics_Drive
! =====================================================================================
subroutine RestartHydrStates(sites,nsites,bc_in,bc_out)
! It is assumed that the following state variables have been read in by
! the restart machinery.
!
! co_hydr%th_ag
! co_hydr%th_troot
! co_hydr%th_aroot
! si_hydr%r_node_shell
! si_hydr%v_shell
! si_hydr%h2osoi_liqvol_shell
! si_hydr%l_aroot_layer
!
! The goal of this subroutine is to call
! the correct sequence of hydraulics initializations to repopulate
! information that relies on these key states, as well as other vegetation
! states such as carbon pools and plant geometry.
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)
! locals
! ----------------------------------------------------------------------------------
! LL pointers
type(fates_patch_type),pointer :: cpatch ! current patch
type(fates_cohort_type),pointer :: ccohort ! current cohort
type(ed_cohort_hydr_type),pointer :: ccohort_hydr
type(ed_site_hydr_type),pointer :: csite_hydr
integer :: s ! site loop counter
integer :: j ! soil layer index
integer :: j_bc ! soil layer index of boundary condition
class(wrf_type_vg), pointer :: wrf_vg
class(wkf_type_vg), pointer :: wkf_vg
class(wrf_type_cch), pointer :: wrf_cch
class(wkf_type_cch), pointer :: wkf_cch
class(wrf_type_smooth_cch), pointer :: wrf_smooth_cch
class(wkf_type_smooth_cch), pointer :: wkf_smooth_cch
real(r8) :: watsat ! Mean wsat across soil layers contributing to current root layer
real(r8) :: sucsat ! Mean sucsat across soil layers contributing to current root layer
real(r8) :: bsw ! Mean bsw across soil layers contributing to current root layer
do s = 1,nsites
csite_hydr=>sites(s)%si_hydr
cpatch => sites(s)%oldest_patch
do while(associated(cpatch))
ccohort => cpatch%shortest
do while(associated(ccohort))
ccohort_hydr => ccohort%co_hydr
! This calculates node heights
call UpdatePlantHydrNodes(ccohort,ccohort%pft,ccohort%height, &
sites(s)%si_hydr)
! This calculates volumes and lengths
call UpdatePlantHydrLenVol(ccohort,csite_hydr)
! This updates the Kmax's of the plant's compartments
call UpdatePlantKmax(ccohort_hydr,ccohort,sites(s)%si_hydr)
! Since this is a newly initialized plant, we set the previous compartment-size
! equal to the ones we just calculated.
call SavePreviousCompartmentVolumes(ccohort_hydr)
ccohort => ccohort%taller
enddo
cpatch => cpatch%younger
end do
sites(s)%si_hydr%l_aroot_layer_init(:) = fates_unset_r8
sites(s)%si_hydr%r_node_shell_init(:,:) = fates_unset_r8
sites(s)%si_hydr%v_shell_init(:,:) = fates_unset_r8
! --------------------------------------------------------------------------------
! Initialize water transfer functions
! which include both water retention functions (WRFs)
! as well as the water conductance (K) functions (WKFs)
! But, this is only for soil!
! --------------------------------------------------------------------------------
! Initialize the Water Retention Functions
! -----------------------------------------------------------------------------------
select case(soil_wrf_type)
case(van_genuchten_type)
do j=1,sites(s)%si_hydr%nlevrhiz
watsat = csite_hydr%AggBCToRhiz(bc_in(s)%watsat_sisl,j,bc_in(s)%dz_sisl)
allocate(wrf_vg)
sites(s)%si_hydr%wrf_soil(j)%p => wrf_vg
call wrf_vg%set_wrf_param([alpha_vg, psd_vg, m_vg, watsat, th_res_vg])
end do
case(campbell_type)
do j=1,sites(s)%si_hydr%nlevrhiz
watsat = csite_hydr%AggBCToRhiz(bc_in(s)%watsat_sisl,j,bc_in(s)%dz_sisl)
sucsat = csite_hydr%AggBCToRhiz(bc_in(s)%sucsat_sisl,j,bc_in(s)%dz_sisl)
bsw = csite_hydr%AggBCToRhiz(bc_in(s)%bsw_sisl,j,bc_in(s)%dz_sisl)
allocate(wrf_cch)
sites(s)%si_hydr%wrf_soil(j)%p => wrf_cch
call wrf_cch%set_wrf_param([watsat, &
(-1.0_r8)*sucsat*denh2o*grav_earth*mpa_per_pa*m_per_mm , &
bsw])
end do
case(smooth1_campbell_type)
do j=1,sites(s)%si_hydr%nlevrhiz
watsat = csite_hydr%AggBCToRhiz(bc_in(s)%watsat_sisl,j,bc_in(s)%dz_sisl)
sucsat = csite_hydr%AggBCToRhiz(bc_in(s)%sucsat_sisl,j,bc_in(s)%dz_sisl)
bsw = csite_hydr%AggBCToRhiz(bc_in(s)%bsw_sisl,j,bc_in(s)%dz_sisl)
allocate(wrf_smooth_cch)
sites(s)%si_hydr%wrf_soil(j)%p => wrf_smooth_cch
call wrf_smooth_cch%set_wrf_param([watsat, &
(-1.0_r8)*sucsat*denh2o*grav_earth*mpa_per_pa*m_per_mm , &
bsw,1._r8])
end do
case(smooth2_campbell_type)
do j=1,sites(s)%si_hydr%nlevrhiz
watsat = csite_hydr%AggBCToRhiz(bc_in(s)%watsat_sisl,j,bc_in(s)%dz_sisl)
sucsat = csite_hydr%AggBCToRhiz(bc_in(s)%sucsat_sisl,j,bc_in(s)%dz_sisl)
bsw = csite_hydr%AggBCToRhiz(bc_in(s)%bsw_sisl,j,bc_in(s)%dz_sisl)
allocate(wrf_smooth_cch)
sites(s)%si_hydr%wrf_soil(j)%p => wrf_smooth_cch
call wrf_smooth_cch%set_wrf_param([watsat, &
(-1.0_r8)*sucsat*denh2o*grav_earth*mpa_per_pa*m_per_mm , &
bsw,2._r8])
end do
case(tfs_type)
write(fates_log(),*) 'TFS water retention curves not available for soil'
call endrun(msg=errMsg(sourcefile, __LINE__))
case default
write(fates_log(),*) 'undefined water retention type for soil:',soil_wrf_type
call endrun(msg=errMsg(sourcefile, __LINE__))
end select
! -----------------------------------------------------------------------------------
! Initialize the Water Conductance (K) Functions
! -----------------------------------------------------------------------------------
select case(soil_wkf_type)
case(van_genuchten_type)
do j=1,sites(s)%si_hydr%nlevrhiz
allocate(wkf_vg)
sites(s)%si_hydr%wkf_soil(j)%p => wkf_vg
call wkf_vg%set_wkf_param([alpha_vg, psd_vg, m_vg, th_sat_vg, th_res_vg, soil_tort_vg])
end do
case(campbell_type)
do j=1,sites(s)%si_hydr%nlevrhiz
watsat = csite_hydr%AggBCToRhiz(bc_in(s)%watsat_sisl,j,bc_in(s)%dz_sisl)
sucsat = csite_hydr%AggBCToRhiz(bc_in(s)%sucsat_sisl,j,bc_in(s)%dz_sisl)
bsw = csite_hydr%AggBCToRhiz(bc_in(s)%bsw_sisl,j,bc_in(s)%dz_sisl)
allocate(wkf_cch)
sites(s)%si_hydr%wkf_soil(j)%p => wkf_cch
call wkf_cch%set_wkf_param([watsat, &
(-1.0_r8)*sucsat*denh2o*grav_earth*mpa_per_pa*m_per_mm , &
bsw])
end do
case(smooth1_campbell_type)
do j=1,sites(s)%si_hydr%nlevrhiz
watsat = csite_hydr%AggBCToRhiz(bc_in(s)%watsat_sisl,j,bc_in(s)%dz_sisl)
sucsat = csite_hydr%AggBCToRhiz(bc_in(s)%sucsat_sisl,j,bc_in(s)%dz_sisl)
bsw = csite_hydr%AggBCToRhiz(bc_in(s)%bsw_sisl,j,bc_in(s)%dz_sisl)
allocate(wkf_smooth_cch)
sites(s)%si_hydr%wkf_soil(j)%p => wkf_smooth_cch
call wkf_smooth_cch%set_wkf_param([watsat, &
(-1.0_r8)*sucsat*denh2o*grav_earth*mpa_per_pa*m_per_mm , &
bsw,1._r8])
end do
case(smooth2_campbell_type)
do j=1,sites(s)%si_hydr%nlevrhiz
watsat = csite_hydr%AggBCToRhiz(bc_in(s)%watsat_sisl,j,bc_in(s)%dz_sisl)
sucsat = csite_hydr%AggBCToRhiz(bc_in(s)%sucsat_sisl,j,bc_in(s)%dz_sisl)
bsw = csite_hydr%AggBCToRhiz(bc_in(s)%bsw_sisl,j,bc_in(s)%dz_sisl)
allocate(wkf_smooth_cch)
sites(s)%si_hydr%wkf_soil(j)%p => wkf_smooth_cch
call wkf_smooth_cch%set_wkf_param([watsat, &
(-1.0_r8)*sucsat*denh2o*grav_earth*mpa_per_pa*m_per_mm , &
bsw,2._r8])
end do
case(tfs_type)
write(fates_log(),*) 'TFS conductance not used in soil'
call endrun(msg=errMsg(sourcefile, __LINE__))
case default
write(fates_log(),*) 'undefined water conductance type for soil:',soil_wkf_type
call endrun(msg=errMsg(sourcefile, __LINE__))
end select
! Update static quantities related to the rhizosphere
call UpdateSizeDepRhizVolLenCon(sites(s), bc_in(s))
! We update the "initial" values of the rhizosphere after
! the previous call to make sure that the conductances are updated
! Now we set the prevous to the current so that the water states
! are not perturbed
call SavePreviousRhizVolumes(sites(s))
call UpdateH2OVeg(sites(s),bc_out(s))
end do
return
end subroutine RestartHydrStates
! ====================================================================================
subroutine InitPlantHydStates(site, cohort)
! REQUIRED INPUTS:
!
! ccohort_hydr%z_node_troot(:)
! ccohort_hydr%z_node_aroot
! ccohort_hydr%z_node_ag
!
! !DESCRIPTION:
!
! !USES:
! !ARGUMENTS:
type(ed_site_type), intent(inout), target :: site ! current site pointer
type(fates_cohort_type), intent(inout), target :: cohort ! current cohort pointer
!
! !LOCAL VARIABLES:
type(ed_site_hydr_type), pointer :: csite_hydr
type(ed_cohort_hydr_type), pointer :: cohort_hydr
integer :: j,k ! layer and node indices
integer :: ft ! functional type index
real(r8) :: psi_rhiz1 ! pressure in first rhizosphere shell [MPa]
real(r8) :: dz ! depth of the current layer [m]
real(r8) :: h_aroot_mean ! minimum total potential of absorbing roots
real(r8), parameter :: psi_aroot_init = -0.2_r8 ! Initialize aroots with -0.2 MPa
real(r8), parameter :: dh_dz = 0.02_r8 ! amount to decrease downstream
! compartment total potentials [MPa/meter]
! In init mode = 1, set absorbing roots to -0.2 MPa
! = 2, use soil as starting point, match total potentials
! and then reduce plant compartment total potential by 1KPa
! for transporting root node, match the lowest total potential
! in absorbing roots
integer, parameter :: init_mode = 2
class(wrf_arr_type),pointer :: wrfa,wrft
class(wkf_arr_type),pointer :: wkfa,wkft
csite_hydr => site%si_hydr
cohort_hydr => cohort%co_hydr
ft = cohort%pft
wrfa => wrf_plant(aroot_p_media,ft)
wkfa => wkf_plant(aroot_p_media,ft)
wrft => wrf_plant(troot_p_media,ft)
wkft => wkf_plant(troot_p_media,ft)
! Set abosrbing root
if(init_mode == 2) then
! h_aroot_mean = 0._r8
do j=1, csite_hydr%nlevrhiz
! Checking apperance of roots. Only proceed if there are roots in that layer
if(cohort_hydr%l_aroot_layer(j) > nearzero) then
! Match the potential of the absorbing root to the inner rhizosphere shell
cohort_hydr%psi_aroot(j) = csite_hydr%wrf_soil(j)%p%psi_from_th(csite_hydr%h2osoi_liqvol_shell(j,1))
! Calculate the mean total potential (include height) of absorbing roots
! h_aroot_mean = h_aroot_mean + cohort_hydr%psi_aroot(j) +
! mpa_per_pa*denh2o*grav_earth*(-csite_hydr%zi_rhiz(j))
cohort_hydr%th_aroot(j) = max(wrfa%p%th_from_psi(cohort_hydr%psi_aroot(j)),wrfa%p%get_thmin())
cohort_hydr%ftc_aroot(j) = wkfa%p%ftc_from_psi(cohort_hydr%psi_aroot(j))
else
cohort_hydr%psi_aroot(j) = psi_aroot_init
cohort_hydr%th_aroot(j) = 0
end if
end do
else
do j=1, csite_hydr%nlevrhiz
cohort_hydr%psi_aroot(j) = psi_aroot_init
! Calculate the mean total potential (include height) of absorbing roots
! h_aroot_mean = h_aroot_mean + cohort_hydr%psi_aroot(j) +
! mpa_per_pa*denh2o*grav_earth*(-csite_hydr%zi_rhiz(j))
cohort_hydr%th_aroot(j) = max(wrfa%p%th_from_psi(cohort_hydr%psi_aroot(j)), &
wrfa%p%get_thmin())
cohort_hydr%ftc_aroot(j) = wkfa%p%ftc_from_psi(cohort_hydr%psi_aroot(j))
end do
end if
!h_aroot_mean = h_aroot_mean/real(csite_hydr%nlevrhiz,r8)
h_aroot_mean = minval(cohort_hydr%psi_aroot(:) + mpa_per_pa*denh2o*grav_earth* &
( -csite_hydr%zi_rhiz(:)+0.5*csite_hydr%dz_rhiz(:) )) ! Get layer centers
! initialize plant water potentials with slight potential gradient (or zero) (dh/dz = C)
! the assumption is made here that initial conditions for soil water will
! be in (or at least close to) hydrostatic equilibrium as well, so that
! it doesn't matter which absorbing root layer the transporting root water
! Set the transporting root to be in equilibrium with mean potential
! of the absorbing roots, minus any gradient we add
cohort_hydr%psi_troot = h_aroot_mean - &
mpa_per_pa*denh2o*grav_earth*cohort_hydr%z_node_troot - dh_dz
cohort_hydr%th_troot = max(wrft%p%th_from_psi(cohort_hydr%psi_troot), &
wrft%p%get_thmin())
cohort_hydr%ftc_troot = wkft%p%ftc_from_psi(cohort_hydr%psi_troot)
! working our way up a tree, assigning water potentials that are in
! hydrostatic equilibrium (minus dh_dz offset) with the water potential immediately below
dz = cohort_hydr%z_node_ag(n_hypool_ag) - cohort_hydr%z_node_troot
cohort_hydr%psi_ag(n_hypool_ag) = cohort_hydr%psi_troot - &
mpa_per_pa*denh2o*grav_earth*dz - dh_dz
cohort_hydr%th_ag(n_hypool_ag) = max(wrf_plant(stem_p_media,ft)%p%get_thmin(), &
wrf_plant(stem_p_media,ft)%p%th_from_psi(cohort_hydr%psi_ag(n_hypool_ag)))
cohort_hydr%ftc_ag(n_hypool_ag) = wkf_plant(stem_p_media,ft)%p%ftc_from_psi(cohort_hydr%psi_ag(n_hypool_ag))
do k=n_hypool_ag-1, 1, -1
dz = cohort_hydr%z_node_ag(k) - cohort_hydr%z_node_ag(k+1)
cohort_hydr%psi_ag(k) = cohort_hydr%psi_ag(k+1) - &
mpa_per_pa*denh2o*grav_earth*dz - &
dh_dz
cohort_hydr%th_ag(k) = max(wrf_plant(csite_hydr%pm_node(k),ft)%p%th_from_psi(cohort_hydr%psi_ag(k)), &
wrf_plant(csite_hydr%pm_node(k),ft)%p%get_thmin())
cohort_hydr%ftc_ag(k) = wkf_plant(csite_hydr%pm_node(k),ft)%p%ftc_from_psi(cohort_hydr%psi_ag(k))
end do
!initialize cohort-level btran
cohort_hydr%btran = wkf_plant(stomata_p_media,ft)%p%ftc_from_psi(cohort_hydr%psi_ag(1))
!flc_gs_from_psi(cohort_hydr%psi_ag(1),cohort%pft)
! We do allow for positive pressures.
! But starting off with positive pressures is something we try to avoid
if ( (cohort_hydr%psi_troot>0.0_r8) .or. &
any(cohort_hydr%psi_ag(:)>0._r8) .or. &
any(cohort_hydr%psi_aroot(:)>0._r8) ) then
write(fates_log(),*) 'Initialized plant compartments with positive pressure?'
write(fates_log(),*) 'psi troot: ',cohort_hydr%psi_troot
write(fates_log(),*) 'psi ag(:): ',cohort_hydr%psi_ag(:)
write(fates_log(),*) 'psi_aroot(:): ',cohort_hydr%psi_aroot(:)
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
end subroutine InitPlantHydStates
! =====================================================================================
subroutine UpdatePlantPsiFTCFromTheta(ccohort,csite_hydr)
! This subroutine updates the potential and the fractional
! of total conductivity based on the relative water
! content
! Arguments
type(fates_cohort_type),intent(inout), target :: ccohort
type(ed_site_hydr_type),intent(in), target :: csite_hydr
! Locals
integer :: ft ! Plant functional type
integer :: k ! loop index for compartments
integer :: j ! Loop index for soil layers
type(ed_cohort_hydr_type), pointer :: ccohort_hydr
ccohort_hydr => ccohort%co_hydr
ft = ccohort%pft
! Update Psi and FTC in above-ground compartments
! -----------------------------------------------------------------------------------
do k = 1,n_hypool_leaf
ccohort_hydr%psi_ag(k) = wrf_plant(leaf_p_media,ft)%p%psi_from_th(ccohort_hydr%th_ag(k))
ccohort_hydr%ftc_ag(k) = wkf_plant(leaf_p_media,ft)%p%ftc_from_psi(ccohort_hydr%psi_ag(k))
end do
do k = n_hypool_leaf+1, n_hypool_ag
ccohort_hydr%psi_ag(k) = wrf_plant(stem_p_media,ft)%p%psi_from_th(ccohort_hydr%th_ag(k))
ccohort_hydr%ftc_ag(k) = wkf_plant(stem_p_media,ft)%p%ftc_from_psi(ccohort_hydr%psi_ag(k))
end do
! Update the Psi and FTC for the transporting root compartment
ccohort_hydr%psi_troot = wrf_plant(troot_p_media,ft)%p%psi_from_th(ccohort_hydr%th_troot)
ccohort_hydr%ftc_troot = wkf_plant(troot_p_media,ft)%p%ftc_from_psi(ccohort_hydr%psi_troot)
! Update the Psi and FTC for the absorbing roots
do j = 1, csite_hydr%nlevrhiz
ccohort_hydr%psi_aroot(j) = wrf_plant(aroot_p_media,ft)%p%psi_from_th(ccohort_hydr%th_aroot(j))
ccohort_hydr%ftc_aroot(j) = wkf_plant(aroot_p_media,ft)%p%ftc_from_psi(ccohort_hydr%psi_aroot(j))
end do
return
end subroutine UpdatePlantPsiFTCFromTheta
! =====================================================================================
subroutine UpdatePlantHydrNodes(ccohort,ft,plant_height,csite_hydr)
! --------------------------------------------------------------------------------
! This subroutine calculates the nodal heights critical to hydraulics in the plant
!
! Inputs: Plant height
! Plant functional type
! Number of soil hydraulic layers
!
! Outputs: cohort_hydr%z_node_ag(:)
! %z_lower_ag(:)
! %z_upper_ag(:)
! %z_node_troot
! %z_node_aroot(:)
! --------------------------------------------------------------------------------
! Arguments
type(fates_cohort_type), intent(inout) :: ccohort
integer,intent(in) :: ft ! plant functional type index
real(r8), intent(in) :: plant_height ! [m]
type(ed_site_hydr_type), intent(in) :: csite_hydr
! Locals
type(ed_cohort_hydr_type), pointer :: ccohort_hydr
integer :: nlevrhiz ! number of rhizosphere layers
real(r8) :: roota ! root profile parameter a zeng2001_crootfr
real(r8) :: rootb ! root profile parameter b zeng2001_crootfr
real(r8) :: crown_depth ! crown depth for the plant [m]
real(r8) :: dz_canopy ! discrete crown depth intervals [m]
real(r8) :: z_stem ! the height of the plants stem below crown [m]
real(r8) :: dz_stem ! vertical stem discretization [m]
real(r8) :: dcumul_rf ! cumulative root distribution discretization [-]
real(r8) :: cumul_rf ! cumulative root distribution where depth is determined [-]
real(r8) :: z_cumul_rf ! depth at which cumul_rf occurs [m]
integer :: k ! Loop counter for compartments
real(r8) :: z_fr ! Maximum rooting depth of the plant [m]
ccohort_hydr => ccohort%co_hydr
! Crown Nodes
! in special case where n_hypool_leaf = 1, the node height of the canopy
! water pool is 1/2 the distance from the bottom of the canopy to the top of the tree
roota = prt_params%fnrt_prof_a(ft)
rootb = prt_params%fnrt_prof_b(ft)
nlevrhiz = csite_hydr%nlevrhiz
!call CrownDepth(plant_height,ft,crown_depth)
crown_depth = min(plant_height,0.1_r8)
dz_canopy = crown_depth / real(n_hypool_leaf,r8)
do k=1,n_hypool_leaf
ccohort_hydr%z_lower_ag(k) = plant_height - dz_canopy*real(k,r8)
ccohort_hydr%z_node_ag(k) = ccohort_hydr%z_lower_ag(k) + 0.5_r8*dz_canopy
ccohort_hydr%z_upper_ag(k) = ccohort_hydr%z_lower_ag(k) + dz_canopy
enddo
! Stem Nodes
! in special case where n_hypool_stem = 1, the node height of the stem water pool is
! 1/2 the height from the ground to the bottom of the canopy
z_stem = plant_height - crown_depth
dz_stem = z_stem / real(n_hypool_stem,r8)
do k=n_hypool_leaf+1,n_hypool_ag
ccohort_hydr%z_upper_ag(k) = real(n_hypool_stem - (k - 1 - n_hypool_leaf),r8)*dz_stem
ccohort_hydr%z_node_ag(k) = ccohort_hydr%z_upper_ag(k) - 0.5_r8*dz_stem
ccohort_hydr%z_lower_ag(k) = ccohort_hydr%z_upper_ag(k) - dz_stem
enddo
call MaximumRootingDepth(ccohort%dbh,ft,csite_hydr%zi_rhiz(nlevrhiz),z_fr)
! Transporting Root Node depth [m] (negative from surface)
call bisect_rootfr(roota, rootb, z_fr, 0._r8, 1.E10_r8, &
0.001_r8, 0.001_r8, 0.5_r8, z_cumul_rf)
if(z_cumul_rf > csite_hydr%zi_rhiz(nlevrhiz) ) then
write(fates_log(),*) 'z_cumul_rf > zi_rhiz(nlevrhiz)?',z_cumul_rf,csite_hydr%zi_rhiz(nlevrhiz)
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
z_cumul_rf = min(z_cumul_rf, abs(csite_hydr%zi_rhiz(nlevrhiz)))
ccohort_hydr%z_node_troot = -z_cumul_rf
return
end subroutine UpdatePlantHydrNodes
! =====================================================================================
subroutine SavePreviousCompartmentVolumes(ccohort_hydr)
type(ed_cohort_hydr_type),intent(inout) :: ccohort_hydr
! Saving the current compartment volumes into an "initial" save-space
! allows us to see how the compartments change size when plants
! change size and effect water contents
ccohort_hydr%v_ag_init(:) = ccohort_hydr%v_ag(:)
ccohort_hydr%v_troot_init = ccohort_hydr%v_troot
ccohort_hydr%v_aroot_layer_init(:) = ccohort_hydr%v_aroot_layer(:)
return
end subroutine SavePreviousCompartmentVolumes
! =====================================================================================
subroutine UpdateSizeDepPlantHydProps(currentSite,ccohort,bc_in)
! DESCRIPTION: Updates absorbing root length (total and its vertical distribution)
! as well as the consequential change in the size of the 'representative' rhizosphere
! shell radii, volumes, and compartment volumes of plant tissues
! !USES:
use shr_sys_mod , only : shr_sys_abort
! ARGUMENTS:
type(ed_site_type) , intent(in) :: currentSite ! Site stuff
type(fates_cohort_type) , intent(inout) :: ccohort ! current cohort pointer
type(bc_in_type) , intent(in) :: bc_in ! Boundary Conditions
! Locals
integer :: nlevrhiz ! Number of total soil layers
type(ed_cohort_hydr_type), pointer :: ccohort_hydr
integer :: ft
nlevrhiz = currentSite%si_hydr%nlevrhiz
ccohort_hydr => ccohort%co_hydr
ft = ccohort%pft
! Save the current vegetation compartment volumes into
! a save space so that it can be compared with the updated quantity.
call SavePreviousCompartmentVolumes(ccohort_hydr)
! This updates all of the z_node positions
call UpdatePlantHydrNodes(ccohort,ft,ccohort%height,currentSite%si_hydr)
! This updates plant compartment volumes, lengths and
! maximum conductances. Make sure for already
! initialized vegetation, that SavePreviousCompartment
! volumes, and UpdatePlantHydrNodes is called prior to this.
call UpdatePlantHydrLenVol(ccohort,currentSite%si_hydr)
! This updates the Kmax's of the plant's compartments
call UpdatePlantKmax(ccohort_hydr,ccohort,currentsite%si_hydr)
end subroutine UpdateSizeDepPlantHydProps
! =====================================================================================
subroutine UpdatePlantHydrLenVol(ccohort,csite_hydr)
! -----------------------------------------------------------------------------------
! This subroutine calculates two attributes of a plant:
! 1) the volumes of storage compartments in the plants
! 2) the lenghts of the organs
! These are not dependent on the hydraulic state of the
! plant, it is more about the structural characteristics and how much biomass
! is present in the different tissues.
!
! Inputs, plant geometries, plant carbon pools, z_node values
!
! -----------------------------------------------------------------------------------
! Arguments
type(fates_cohort_type),intent(inout) :: ccohort
type(ed_site_hydr_type),intent(in) :: csite_hydr
type(ed_cohort_hydr_type),pointer :: ccohort_hydr ! Plant hydraulics structure
integer :: j,k
integer :: ft ! Plant functional type index
real(r8) :: roota ! root profile parameter a zeng2001_crootfr
real(r8) :: rootb ! root profile parameter b zeng2001_crootfr
real(r8) :: leaf_c ! Current amount of leaf carbon in the plant [kg]
real(r8) :: leaf_c_target ! Target leaf carbon (with some conditions) [kgC]
real(r8) :: fnrt_c ! Current amount of fine-root carbon in the plant [kg]
real(r8) :: sapw_c ! Current amount of sapwood carbon in the plant [kg]
real(r8) :: struct_c ! Current amount of structural carbon in the plant [kg]
real(r8) :: woody_bg_c ! belowground woody biomass in carbon units [kgC/indiv]
real(r8) :: z_stem ! the height of the plants stem below crown [m]
real(r8) :: sla ! specific leaf area [cm2/g]
real(r8) :: v_aroot_tot ! total compartment volume of all absorbing roots for cohort [m3]
real(r8) :: l_aroot_tot ! total length of absorbing roots for cohrot [m]
real(r8) :: denleaf ! leaf dry mass per unit fresh leaf volume [kg/m3]
real(r8) :: a_sapwood ! sapwood area [m2]
real(r8) :: a_sapwood_target ! sapwood cross-section area at reference height, at target biomass [m2]
real(r8) :: sapw_c_target ! sapwood carbon, at target [kgC]
real(r8) :: v_sapwood ! sapwood volume [m3]
real(r8) :: v_troot ! transporting root volume [m3/indiv]
real(r8) :: rootfr ! mass fraction of roots in each layer [kg/kg]
real(r8) :: crown_depth ! Depth of the plant's crown [m]
real(r8) :: norm ! total root fraction used <1
integer :: nlevrhiz ! number of rhizosphere levels
real(r8) :: dbh ! the dbh of current cohort [cm]
real(r8) :: z_fr ! rooting depth of a cohort [cm]
real(r8) :: v_leaf_donate(1:n_hypool_leaf) ! the volume that leaf will donate to xylem
! We allow the transporting root to donate a fraction of its volume to the absorbing
! roots to help mitigate numerical issues due to very small volumes. This is the
! fraction the transporting roots donate to those layers
real(r8), parameter :: t2aroot_vol_donate_frac = 0.65_r8
real(r8), parameter :: l2sap_vol_donate_frac = 0.5_r8 ! Junyan added
real(r8), parameter :: min_leaf_frac = 0.1_r8 ! Fraction of maximum leaf carbon that
! we set as our lower cap on leaf volume
real(r8), parameter :: min_trim = 0.1_r8 ! The lower cap on trimming function used
! to estimate maximum leaf carbon
ccohort_hydr => ccohort%co_hydr
ft = ccohort%pft
nlevrhiz = csite_hydr%nlevrhiz
leaf_c = ccohort%prt%GetState(leaf_organ, carbon12_element)
sapw_c = ccohort%prt%GetState(sapw_organ, carbon12_element)
fnrt_c = ccohort%prt%GetState(fnrt_organ, carbon12_element)
struct_c = ccohort%prt%GetState(struct_organ, carbon12_element)
roota = prt_params%fnrt_prof_a(ft)
rootb = prt_params%fnrt_prof_b(ft)
! Leaf Volumes
! -----------------------------------------------------------------------------------
! NOTE: SLATOP currently does not use any vertical scaling functions
! but that may not be so forever. ie sla = slatop (RGK-082017)
! m2/gC * cm2/m2 -> cm2/gC
sla = prt_params%slatop(ft) * cm2_per_m2
! empirical regression data from leaves at Caxiuana (~ 8 spp)
denleaf = -2.3231_r8*sla/prt_params%c2b(ft) + 781.899_r8
! Leaf volumes
! Note: Leaf volumes of zero is problematic for two reasons. Zero volumes create
! numerical difficulties, and they could also create problems when a leaf is trying
! to re-flush.
! Therefore, if the leaf is in an "off" status, then we do not update the leaf
! volume. This way the volume is where it was when it dropped, and this is consistent
! with the theory that leaf water potentials drive growth and re-flushing, not the
! other way around. However, it is possible that we may have recruits with an
! "off" status (due to external seed rain active during a dry or cold season). If a
! cohort is newly created, we must give it a starting volume.
! We also place a lower bound on how low the leaf volume is allowed to go, which is 10%
! of the plant's carrying capacity.
! [kgC] * [kg/kgC] / [kg/m3] -> [m3]
! Get the target, or rather, maximum leaf carrying capacity of plant
! Lets also avoid super-low targets that have very low trimming functions
call bleaf(ccohort%dbh,ccohort%pft,ccohort%crowndamage, &
max(ccohort%canopy_trim,min_trim),ccohort%efleaf_coh, leaf_c_target)
if( (ccohort%status_coh == leaves_on) .or. ccohort_hydr%is_newly_recruited ) then
ccohort_hydr%v_ag(1:n_hypool_leaf) = max(leaf_c,min_leaf_frac*leaf_c_target) * &
prt_params%c2b(ft) / denleaf/ real(n_hypool_leaf,r8)
end if
! Step sapwood volume
! -----------------------------------------------------------------------------------
! BOC...may be needed for testing/comparison w/ v_sapwood
! kg / ( g cm-3 * cm3/m3 * kg/g ) -> m3
! v_stem = c_stem_biom / (prt_params%wood_density(ft) * kg_per_g * cm3_per_m3 )
! calculate the sapwood cross-sectional area
call bsap_allom(ccohort%dbh,ccohort%pft,ccohort%crowndamage, &
ccohort%canopy_trim, ccohort%efstem_coh, a_sapwood_target,sapw_c_target)