Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Dynamic Patch Arrays - larger nclmax and nlevleaf #1198

Merged
merged 24 commits into from
Aug 14, 2024
Merged
Show file tree
Hide file tree
Changes from 3 commits
Commits
Show all changes
24 commits
Select commit Hold shift + click to select a range
e1fb819
removing dependency on small nclmax
rgknox May 10, 2024
5fb8e97
final pass at the first version of dynamic patch arrays
rgknox May 10, 2024
aa370ba
Added ncl_p protections
rgknox May 10, 2024
d9f25a0
Reverted some changes where we update ncl_p, but it doesn't need to b…
rgknox May 10, 2024
6dc21d1
Revert some removed b4b code, will re-remove later under its own chan…
rgknox May 10, 2024
556edd4
nclmax nlevleaf: reverting for b4b testing
rgknox May 16, 2024
b4ab283
added memory cleanup to new patch allocated arrays
rgknox May 22, 2024
40e67cf
Moved patch dynamic reallocations, zeroing and nanning into FatesPatc…
rgknox May 23, 2024
e4576eb
Merge branch 'main' into remove_nclmax
rgknox May 23, 2024
256dd62
Updated nclmax and nlevleaf
rgknox May 23, 2024
ea19532
Merge branch 'main' into remove_nclmax-merged
rgknox Jul 9, 2024
a2b3ead
changed ncan to nleaf
rgknox Jul 9, 2024
bcec046
fix on recruit l2fr
rgknox Jul 10, 2024
85c6821
Merge tag 'sci.1.77.1_api.36.0.0' into remove_nclmax
glemieux Jul 19, 2024
fc20a61
Merge branch 'remove_nclmax' of github.com:rgknox/fates into remove_n…
rgknox Jul 23, 2024
d71c9f0
reverting nclmax to 3 for tests
rgknox Jul 24, 2024
cccfac6
More zeroing of dynamic arrays
rgknox Jul 25, 2024
248595c
resetting nclmax to 2, since there is a bug that is preventing nclmax…
rgknox Jul 26, 2024
0897ba1
attempts to get b4b, unsuccesfull so far
rgknox Jul 30, 2024
d914891
zeroing psn_z
rgknox Jul 30, 2024
ed62e82
remove canopy closure clause
rgknox Jul 30, 2024
9b6b336
fixes to summed lai
rgknox Jul 31, 2024
8da0669
small additions to patch arrays
rgknox Jul 31, 2024
ae32bd4
Merge branch 'main' into remove_nclmax
rgknox Jul 31, 2024
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
371 changes: 184 additions & 187 deletions biogeochem/EDCanopyStructureMod.F90

Large diffs are not rendered by default.

1 change: 0 additions & 1 deletion biogeochem/EDPatchDynamicsMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -29,7 +29,6 @@ module EDPatchDynamicsMod
use EDTypesMod , only : site_fluxdiags_type
use EDTypesMod , only : min_patch_area
use EDTypesMod , only : min_patch_area_forced
use EDParamsMod , only : nclmax
use EDParamsMod , only : regeneration_model
use FatesInterfaceTypesMod, only : numpft
use FatesConstantsMod , only : dtype_ifall
Expand Down
2 changes: 1 addition & 1 deletion biogeochem/EDPhysiologyMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -3301,7 +3301,7 @@ subroutine UpdateRecruitL2FR(csite)
end do

! Find the daily mean for each PFT weighted by number and add it to the running mean
do cl = 1,nclmax
do cl = 1,cpatch%ncl_p
do ft = 1,numpft
if(rec_n(ft,cl)>nearzero)then
rec_l2fr0(ft,cl) = rec_l2fr0(ft,cl) / rec_n(ft,cl)
Expand Down
5 changes: 2 additions & 3 deletions biogeochem/FatesAllometryMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -98,7 +98,6 @@ module FatesAllometryMod
use FatesGlobals , only : endrun => fates_endrun
use FatesGlobals , only : FatesWarn,N2S,A2S,I2S
use EDParamsMod , only : nlevleaf,dinc_vai,dlower_vai
use EDParamsMod , only : nclmax
use DamageMainMod , only : GetCrownReduction

implicit none
Expand Down Expand Up @@ -670,7 +669,7 @@ real(r8) function tree_lai( leaf_c, pft, c_area, nplant, cl, canopy_lai, vcmax25
real(r8), intent(in) :: c_area ! areal extent of canopy (m2)
real(r8), intent(in) :: nplant ! number of individuals in cohort per ha
integer, intent(in) :: cl ! canopy layer index
real(r8), intent(in) :: canopy_lai(nclmax) ! total leaf area index of
real(r8), intent(in) :: canopy_lai(:) ! total leaf area index of
! each canopy layer
real(r8), intent(in) :: vcmax25top ! maximum carboxylation rate at canopy
! top, ref 25C
Expand Down Expand Up @@ -803,7 +802,7 @@ real(r8) function tree_sai(pft, dbh, crowndamage, canopy_trim, elongf_stem, c_ar
real(r8), intent(in) :: c_area ! crown area (m2)
real(r8), intent(in) :: nplant ! number of plants
integer, intent(in) :: cl ! canopy layer index
real(r8), intent(in) :: canopy_lai(nclmax) ! total leaf area index of
real(r8), intent(in) :: canopy_lai(:) ! total leaf area index of
! each canopy layer
real(r8), intent(in) :: treelai ! tree LAI for checking purposes only
real(r8), intent(in) :: vcmax25top ! maximum carboxylation rate at top of crown
Expand Down
3 changes: 1 addition & 2 deletions biogeochem/FatesCohortMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -6,7 +6,6 @@ module FatesCohortMod
use FatesConstantsMod, only : nearzero
use FatesConstantsMod, only : ican_upper, ican_ustory
use EDParamsMod, only : nlevleaf
use EDParamsMod, only : nclmax
use FatesGlobals, only : endrun => fates_endrun
use FatesGlobals, only : fates_log
use PRTGenericMod, only : max_nleafage
Expand Down Expand Up @@ -560,7 +559,7 @@ subroutine Create(this, prt, pft, nn, height, coage, dbh, status, &
real(r8), intent(in) :: ctrim ! fraction of the maximum leaf biomass
real(r8), intent(in) :: spread ! how spread crowns are in horizontal space
real(r8), intent(in) :: carea ! area of cohort, for SP mode [m2]
real(r8), intent(in) :: can_tlai(nclmax) ! patch-level total LAI of each leaf layer
real(r8), intent(in) :: can_tlai(:) ! patch-level total LAI of each canopy layer
real(r8), intent(in) :: elongf_leaf ! leaf elongation factor [fraction]
real(r8), intent(in) :: elongf_fnrt ! fine-root "elongation factor" [fraction]
real(r8), intent(in) :: elongf_stem ! stem "elongation factor" [fraction]
Expand Down
54 changes: 30 additions & 24 deletions biogeochem/FatesPatchMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -99,23 +99,29 @@ module FatesPatchMod
real(r8) :: total_canopy_area ! area that is covered by vegetation [m2]
real(r8) :: total_tree_area ! area that is covered by woody vegetation [m2]
real(r8) :: zstar ! height of smallest canopy tree, only meaningful in "strict PPA" mode [m]
real(r8) :: elai_profile(nclmax,maxpft,nlevleaf) ! exposed leaf area in each canopy layer, pft, and leaf layer [m2 leaf/m2 contributing crown area]
real(r8) :: esai_profile(nclmax,maxpft,nlevleaf) ! exposed stem area in each canopy layer, pft, and leaf layer [m2 leaf/m2 contributing crown area]
real(r8) :: tlai_profile(nclmax,maxpft,nlevleaf)
real(r8) :: tsai_profile(nclmax,maxpft,nlevleaf)
real(r8) :: canopy_area_profile(nclmax,maxpft,nlevleaf) ! fraction of crown area per canopy area in each layer
! they will sum to 1.0 in the fully closed canopy layers
! but only in leaf-layers that contain contributions
! from all cohorts that donate to canopy_area

! exposed leaf area in each canopy layer, pft, and leaf layer [m2 leaf/m2 contributing crown area]
real(r8), allocatable :: elai_profile(:,:,:) ! nclmax,maxpft,nlevleaf)
! exposed stem area in each canopy layer, pft, and leaf layer [m2 leaf/m2 contributing crown area]
real(r8), allocatable :: esai_profile(:,:,:) ! nclmax,maxpft,nlevleaf)
! total leaf area (includes that which is under snow-pack)
real(r8), allocatable :: tlai_profile(:,:,:) ! nclmax,maxpft,nlevleaf)
! total stem area (includes that which is under snow-pack)
real(r8), allocatable :: tsai_profile(:,:,:) ! nclmax,maxpft,nlevleaf)

real(r8), allocatable :: canopy_area_profile(:,:,:) ! nclmax,maxpft,nlevleaf) ! fraction of crown area per canopy area in each layer
! they will sum to 1.0 in the fully closed canopy layers
! but only in leaf-layers that contain contributions
! from all cohorts that donate to canopy_area

integer :: canopy_mask(nclmax,maxpft) ! is there any of this pft in this canopy layer?
integer :: nrad(nclmax,maxpft) ! number of exposed leaf layers for each canopy layer and pft
integer :: ncan(nclmax,maxpft) ! number of total leaf layers for each canopy layer and pft
real(r8) :: c_stomata ! mean stomatal conductance of all leaves in the patch [umol/m2/s]
real(r8) :: c_lblayer ! mean boundary layer conductance of all leaves in the patch [umol/m2/s]

real(r8) :: psn_z(nclmax,maxpft,nlevleaf)
real(r8) :: nrmlzd_parprof_pft_dir_z(num_rad_stream_types,nclmax,maxpft,nlevleaf)
real(r8) :: nrmlzd_parprof_pft_dif_z(num_rad_stream_types,nclmax,maxpft,nlevleaf)
real(r8),allocatable :: nrmlzd_parprof_pft_dir_z(:,:,:,:) !num_rad_stream_types,nclmax,maxpft,nlevleaf)
real(r8),allocatable :: nrmlzd_parprof_pft_dif_z(:,:,:,:) !num_rad_stream_types,nclmax,maxpft,nlevleaf)

!---------------------------------------------------------------------------

Expand All @@ -129,20 +135,20 @@ module FatesPatchMod


! organized by canopy layer, pft, and leaf layer
real(r8) :: fabd_sun_z(nclmax,maxpft,nlevleaf) ! sun fraction of direct light absorbed [0-1]
real(r8) :: fabd_sha_z(nclmax,maxpft,nlevleaf) ! shade fraction of direct light absorbed [0-1]
real(r8) :: fabi_sun_z(nclmax,maxpft,nlevleaf) ! sun fraction of indirect light absorbed [0-1]
real(r8) :: fabi_sha_z(nclmax,maxpft,nlevleaf) ! shade fraction of indirect light absorbed [0-1]
real(r8) :: ed_parsun_z(nclmax,maxpft,nlevleaf) ! PAR absorbed in the sun [W/m2]
real(r8) :: ed_parsha_z(nclmax,maxpft,nlevleaf) ! PAR absorbed in the shade [W/m2]
real(r8) :: f_sun(nclmax,maxpft,nlevleaf) ! fraction of leaves in the sun [0-1]
real(r8) :: ed_laisun_z(nclmax,maxpft,nlevleaf)
real(r8) :: ed_laisha_z(nclmax,maxpft,nlevleaf)
real(r8),allocatable :: fabd_sun_z(:,:,:) !nclmax,maxpft,nlevleaf) ! sun fraction of direct light absorbed [0-1]
real(r8),allocatable :: fabd_sha_z(:,:,:) !nclmax,maxpft,nlevleaf) ! shade fraction of direct light absorbed [0-1]
real(r8),allocatable :: fabi_sun_z(:,:,:) !nclmax,maxpft,nlevleaf) ! sun fraction of indirect light absorbed [0-1]
real(r8),allocatable :: fabi_sha_z(:,:,:) !nclmax,maxpft,nlevleaf) ! shade fraction of indirect light absorbed [0-1]
real(r8),allocatable :: ed_parsun_z(:,:,:) !nclmax,maxpft,nlevleaf) ! PAR absorbed in the sun [W/m2]
real(r8),allocatable :: ed_parsha_z(:,:,:) !nclmax,maxpft,nlevleaf) ! PAR absorbed in the shade [W/m2]
real(r8),allocatable :: f_sun(:,:,:) !nclmax,maxpft,nlevleaf) ! fraction of leaves in the sun [0-1]
real(r8),allocatable :: ed_laisun_z(:,:,:) !nclmax,maxpft,nlevleaf)
real(r8),allocatable :: ed_laisha_z(:,:,:) !nclmax,maxpft,nlevleaf)


! radiation profiles for comparison against observations
real(r8) :: parprof_pft_dir_z(nclmax,maxpft,nlevleaf) ! direct-beam PAR profile through canopy, by canopy, PFT, leaf level [W/m2]
real(r8) :: parprof_pft_dif_z(nclmax,maxpft,nlevleaf) ! diffuse PAR profile through canopy, by canopy, PFT, leaf level [W/m2]
real(r8),allocatable :: parprof_pft_dir_z(:,:,:) !nclmax,maxpft,nlevleaf) ! direct-beam PAR profile through canopy, by canopy, PFT, leaf level [W/m2]
real(r8),allocatable :: parprof_pft_dif_z(:,:,:) !nclmax,maxpft,nlevleaf) ! diffuse PAR profile through canopy, by canopy, PFT, leaf level [W/m2]

real(r8), allocatable :: tr_soil_dir(:) ! fraction of incoming direct radiation transmitted to the soil as direct, by numSWB [0-1]
real(r8), allocatable :: tr_soil_dif(:) ! fraction of incoming diffuse radiation that is transmitted to the soil as diffuse [0-1]
Expand Down Expand Up @@ -321,7 +327,7 @@ subroutine NanValues(this)
this%c_stomata = nan
this%c_lblayer = nan

this%psn_z(:,:,:) = nan
!this%psn_z(:,:,:) = nan
this%nrmlzd_parprof_pft_dir_z(:,:,:,:) = nan
this%nrmlzd_parprof_pft_dif_z(:,:,:,:) = nan

Expand Down Expand Up @@ -412,7 +418,7 @@ subroutine ZeroValues(this)
this%elai_profile(:,:,:) = 0.0_r8
this%c_stomata = 0.0_r8
this%c_lblayer = 0.0_r8
this%psn_z(:,:,:) = 0.0_r8
! this%psn_z(:,:,:) = 0.0_r8
this%nrmlzd_parprof_pft_dir_z(:,:,:,:) = 0.0_r8
this%nrmlzd_parprof_pft_dif_z(:,:,:,:) = 0.0_r8

Expand Down
13 changes: 8 additions & 5 deletions biogeophys/FatesPlantRespPhotosynthMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -186,13 +186,16 @@ subroutine FatesPlantRespPhotosynthDrive (nsites, sites,bc_in,bc_out,dtime)

! 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)

! Photosynthesis [umol /m2 /s]
real(r8) :: psn_z(nlevleaf,maxpft,nclmax)

! Mask used to determine which leaf-layer biophysical rates have been
! used already
logical :: rate_mask_z(nlevleaf,maxpft,nclmax)
Expand Down Expand Up @@ -723,7 +726,7 @@ subroutine FatesPlantRespPhotosynthDrive (nsites, sites,bc_in,bc_out,dtime)
lmr_z(iv,ft,cl), & ! in
leaf_psi, & ! in
bc_in(s)%rb_pa(ifp), & ! in
currentPatch%psn_z(cl,ft,iv), & ! out
psn_z(iv,ft,cl), & ! out
rs_z(iv,ft,cl), & ! out
anet_av_z(iv,ft,cl), & ! out
c13disc_z(cl,ft,iv)) ! out
Expand Down Expand Up @@ -755,7 +758,7 @@ subroutine FatesPlantRespPhotosynthDrive (nsites, sites,bc_in,bc_out,dtime)
if(radiation_model.eq.norman_solver) then

call ScaleLeafLayerFluxToCohort(nv, & !in
currentPatch%psn_z(cl,ft,1:nv), & !in
psn_z(1:nv,ft,cl), & !in
lmr_z(1:nv,ft,cl), & !in
rs_z(1:nv,ft,cl), & !in
currentPatch%elai_profile(cl,ft,1:nv), & !in
Expand All @@ -773,7 +776,7 @@ subroutine FatesPlantRespPhotosynthDrive (nsites, sites,bc_in,bc_out,dtime)
else

call ScaleLeafLayerFluxToCohort(nv, & !in
currentPatch%psn_z(cl,ft,1:nv), & !in
psn_z(1:nv,ft,cl), & !in
lmr_z(1:nv,ft,cl), & !in
rs_z(1:nv,ft,cl), & !in
cohort_layer_elai(1:nv), & !in
Expand Down Expand Up @@ -1991,7 +1994,7 @@ subroutine UpdateCanopyNCanNRadPresent(currentPatch)
currentPatch%nrad = currentPatch%ncan

! Now loop through and identify which layer and pft combo has scattering elements
do cl = 1,nclmax
do cl = 1,currentPatch%ncl_p
do ft = 1,numpft
currentPatch%canopy_mask(cl,ft) = 0
do iv = 1, currentPatch%nrad(cl,ft);
Expand Down
7 changes: 5 additions & 2 deletions main/EDParamsMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -104,10 +104,13 @@ module EDParamsMod

real(r8), parameter, public :: soil_tfrz_thresh = -2.0_r8 ! Soil temperature threshold below which hydraulic failure mortality is off (non-hydro only) in degrees C

integer, parameter, public :: nclmax = 2 ! Maximum number of canopy layers
integer, parameter, public :: nclmax = 5 ! Maximum number of canopy layers (used only for scratch arrays)
! We would make this even higher, but making this
! a little lower keeps the size down on some output arrays
! For large arrays at patch level we use dynamic allocation

! parameters that govern the VAI (LAI+SAI) bins used in radiative transfer code
integer, parameter, public :: nlevleaf = 30 ! number of leaf+stem layers in each canopy layer
integer, parameter, public :: nlevleaf = 50 ! number of leaf+stem layers in each canopy layer

real(r8), public :: dinc_vai(nlevleaf) = fates_unset_r8 ! VAI bin widths array
real(r8), public :: dlower_vai(nlevleaf) = fates_unset_r8 ! lower edges of VAI bins
Expand Down
2 changes: 1 addition & 1 deletion main/FatesRestartInterfaceMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -3688,7 +3688,7 @@ subroutine update_3dpatch_radiation(this, nsites, sites, bc_out)
currentpatch => sites(s)%oldest_patch
do while (associated(currentpatch))
ifp = ifp+1

currentPatch%f_sun (:,:,:) = 0._r8
currentPatch%fabd_sun_z (:,:,:) = 0._r8
currentPatch%fabd_sha_z (:,:,:) = 0._r8
Expand Down
2 changes: 1 addition & 1 deletion radiation/FatesNormanRadMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -189,7 +189,7 @@ subroutine PatchNormanRadiation (currentPatch, &
tau_layer(:,:,:,:)=0.0_r8
f_abs(:,:,:,:)=0.0_r8
f_abs_leaf(:,:,:,:)=0._r8
do L = 1,nclmax
do L = 1,currentPatch%ncl_p
do ft = 1,numpft
currentPatch%canopy_mask(L,ft) = 0
do iv = 1, currentPatch%nrad(L,ft)
Expand Down
2 changes: 1 addition & 1 deletion radiation/FatesRadiationDriveMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -471,7 +471,7 @@ subroutine FatesSunShadeFracs(nsites, sites,bc_in,bc_out)
end do do_icol

do ft = 1,numpft
do_iv: do iv = 1, nlevleaf
do_iv: do iv = 1, cpatch%ncan(cl,ft)! nlevleaf
if(area_vlpfcl(iv,ft,cl)<nearzero) exit do_iv
cpatch%parprof_pft_dir_z(cl,ft,iv) = &
cpatch%parprof_pft_dir_z(cl,ft,iv) / area_vlpfcl(iv,ft,cl)
Expand Down
3 changes: 0 additions & 3 deletions radiation/FatesTwoStreamUtilsMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -87,9 +87,6 @@ subroutine FatesConstructRadElements(site,fcansno_pa,coszen_pa)
!real(r8), parameter :: init_max_vai_diff_per_elem = 0.2_r8
!type(fates_cohort_type), pointer :: elem_co_ptrs(ncl*max_el_per_layer,100)




max_elements = -1
ifp=0
patch => site%oldest_patch
Expand Down