Skip to content
This repository has been archived by the owner on Sep 14, 2018. It is now read-only.

disaggregate max mort rate params and some updates to ncvarsort.py #7

Merged
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
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
9 changes: 4 additions & 5 deletions biogeochem/EDMortalityFunctionsMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -12,7 +12,6 @@ module EDMortalityFunctionsMod
use EDTypesMod , only : ed_patch_type
use FatesConstantsMod , only : itrue,ifalse
use FatesAllometryMod , only : bleaf
use EDParamsMod , only : ED_val_stress_mort
use FatesInterfaceMod , only : bc_in_type
use FatesInterfaceMod , only : hlm_use_ed_prescribed_phys
use FatesInterfaceMod , only : hlm_freq_day
Expand Down Expand Up @@ -62,7 +61,6 @@ subroutine mortality_rates( cohort_in,bc_in,cmort,hmort,bmort,frmort )
real(r8) :: hf_sm_threshold ! hydraulic failure soil moisture threshold
real(r8) :: temp_dep ! Temp. function (freezing mortality)
real(r8) :: temp_in_C ! Daily averaged temperature in Celcius
real(r8),parameter :: frost_mort_scaler = 3.0_r8 ! Scaling factor for freezing mortality
real(r8),parameter :: frost_mort_buffer = 5.0_r8 ! 5deg buffer for freezing mortality

logical, parameter :: test_zero_mortality = .false. ! Developer test which
Expand All @@ -79,7 +77,7 @@ subroutine mortality_rates( cohort_in,bc_in,cmort,hmort,bmort,frmort )
hf_sm_threshold = EDPftvarcon_inst%hf_sm_threshold(cohort_in%pft)

if(cohort_in%patchptr%btran_ft(cohort_in%pft) <= hf_sm_threshold)then
hmort = ED_val_stress_mort
hmort = EDPftvarcon_inst%mort_scalar_hydrfailure(cohort_in%pft)
else
hmort = 0.0_r8
endif
Expand All @@ -89,7 +87,8 @@ subroutine mortality_rates( cohort_in,bc_in,cmort,hmort,bmort,frmort )
call bleaf(cohort_in%dbh,cohort_in%pft,cohort_in%canopy_trim,b_leaf)
if( b_leaf > 0._r8 .and. cohort_in%bstore <= b_leaf )then
frac = cohort_in%bstore/ b_leaf
cmort = max(0.0_r8,ED_val_stress_mort*(1.0_r8 - frac))
cmort = max(0.0_r8,EDPftvarcon_inst%mort_scalar_cstarvation(cohort_in%pft) * &
(1.0_r8 - frac))
else
cmort = 0.0_r8
endif
Expand All @@ -109,7 +108,7 @@ subroutine mortality_rates( cohort_in,bc_in,cmort,hmort,bmort,frmort )
temp_in_C = bc_in%t_veg24_si - tfrz
temp_dep = max(0.0,min(1.0,1.0 - (temp_in_C - &
EDPftvarcon_inst%freezetol(cohort_in%pft))/frost_mort_buffer) )
frmort = frost_mort_scaler * temp_dep
frmort = EDPftvarcon_inst%mort_scalar_coldstress(cohort_in%pft) * temp_dep


!mortality_rates = bmort + hmort + cmort
Expand Down
10 changes: 0 additions & 10 deletions main/EDParamsMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -22,7 +22,6 @@ module EDParamsMod

real(r8),protected :: fates_mortality_disturbance_fraction ! the fraction of canopy mortality that results in disturbance
real(r8),protected :: ED_val_comp_excln
real(r8),protected :: ED_val_stress_mort
real(r8),protected :: ED_val_init_litter
real(r8),protected :: ED_val_nignitions
real(r8),protected :: ED_val_understorey_death
Expand Down Expand Up @@ -50,7 +49,6 @@ module EDParamsMod

character(len=param_string_length),parameter :: ED_name_mort_disturb_frac = "fates_mort_disturb_frac"
character(len=param_string_length),parameter :: ED_name_comp_excln = "fates_comp_excln"
character(len=param_string_length),parameter :: ED_name_stress_mort = "fates_stress_mort"
character(len=param_string_length),parameter :: ED_name_init_litter = "fates_init_litter"
character(len=param_string_length),parameter :: ED_name_nignitions = "fates_fire_nignitions"
character(len=param_string_length),parameter :: ED_name_understorey_death = "fates_mort_understorey_death"
Expand Down Expand Up @@ -128,7 +126,6 @@ subroutine FatesParamsInit()

fates_mortality_disturbance_fraction = nan
ED_val_comp_excln = nan
ED_val_stress_mort = nan
ED_val_init_litter = nan
ED_val_nignitions = nan
ED_val_understorey_death = nan
Expand Down Expand Up @@ -187,9 +184,6 @@ subroutine FatesRegisterParams(fates_params)
call fates_params%RegisterParameter(name=ED_name_comp_excln, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names)

call fates_params%RegisterParameter(name=ED_name_stress_mort, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names)

call fates_params%RegisterParameter(name=ED_name_init_litter, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names)

Expand Down Expand Up @@ -302,9 +296,6 @@ subroutine FatesReceiveParams(fates_params)
call fates_params%RetreiveParameter(name=ED_name_comp_excln, &
data=ED_val_comp_excln)

call fates_params%RetreiveParameter(name=ED_name_stress_mort, &
data=ED_val_stress_mort)

call fates_params%RetreiveParameter(name=ED_name_init_litter, &
data=ED_val_init_litter)

Expand Down Expand Up @@ -416,7 +407,6 @@ subroutine FatesReportParams(is_master)
write(fates_log(),*) '----------- FATES Scalar Parameters -----------------'
write(fates_log(),fmt0) 'fates_mortality_disturbance_fraction = ',fates_mortality_disturbance_fraction
write(fates_log(),fmt0) 'ED_val_comp_excln = ',ED_val_comp_excln
write(fates_log(),fmt0) 'ED_val_stress_mort = ',ED_val_stress_mort
write(fates_log(),fmt0) 'ED_val_init_litter = ',ED_val_init_litter
write(fates_log(),fmt0) 'ED_val_nignitions = ',ED_val_nignitions
write(fates_log(),fmt0) 'ED_val_understorey_death = ',ED_val_understorey_death
Expand Down
46 changes: 38 additions & 8 deletions main/EDPftvarcon.F90
Original file line number Diff line number Diff line change
Expand Up @@ -76,6 +76,9 @@ module EDPftvarcon
real(r8), allocatable :: maintresp_reduction_intercept(:) ! intercept of MR reduction as f(carbon storage),
! 0=no throttling, 1=max throttling
real(r8), allocatable :: bmort(:)
real(r8), allocatable :: mort_scalar_coldstress(:)
real(r8), allocatable :: mort_scalar_cstarvation(:)
real(r8), allocatable :: mort_scalar_hydrfailure(:)
real(r8), allocatable :: hf_sm_threshold(:)
real(r8), allocatable :: vcmaxha(:)
real(r8), allocatable :: jmaxha(:)
Expand Down Expand Up @@ -281,7 +284,7 @@ subroutine Register_PFT(this, fates_params)
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

name = 'fates_seed_hgt_min'
name = 'fates_recruit_hgt_min'
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

Expand All @@ -305,7 +308,7 @@ subroutine Register_PFT(this, fates_params)
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

name = 'fates_seed_initd'
name = 'fates_recruit_initd'
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

Expand Down Expand Up @@ -409,11 +412,11 @@ subroutine Register_PFT(this, fates_params)
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

name = 'fates_leafcn'
name = 'fates_leaf_c2n_ratio'
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

name = 'fates_root_frootcn'
name = 'fates_froot_c2n_ratio'
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

Expand Down Expand Up @@ -589,6 +592,18 @@ subroutine Register_PFT(this, fates_params)
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

name = 'fates_mort_scalar_coldstress'
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

name = 'fates_mort_scalar_cstarvation'
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

name = 'fates_mort_scalar_hydrfailure'
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)

name = 'fates_mort_hf_sm_threshold'
call fates_params%RegisterParameter(name=name, dimension_shape=dimension_shape_1d, &
dimension_names=dim_names, lower_bounds=dim_lower_bound)
Expand Down Expand Up @@ -697,7 +712,7 @@ subroutine Receive_PFT(this, fates_params)
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%wood_density)

name = 'fates_seed_hgt_min'
name = 'fates_recruit_hgt_min'
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%hgt_min)

Expand All @@ -721,7 +736,7 @@ subroutine Receive_PFT(this, fates_params)
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%crown_kill)

name = 'fates_seed_initd'
name = 'fates_recruit_initd'
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%initd)

Expand Down Expand Up @@ -821,11 +836,11 @@ subroutine Receive_PFT(this, fates_params)
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%vcmax25top)

name = 'fates_leafcn'
name = 'fates_leaf_c2n_ratio'
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%leafcn)

name = 'fates_root_frootcn'
name = 'fates_froot_c2n_ratio'
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%frootcn)

Expand Down Expand Up @@ -1005,6 +1020,18 @@ subroutine Receive_PFT(this, fates_params)
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%bmort)

name = 'fates_mort_scalar_coldstress'
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%mort_scalar_coldstress)

name = 'fates_mort_scalar_cstarvation'
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%mort_scalar_cstarvation)

name = 'fates_mort_scalar_hydrfailure'
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%mort_scalar_hydrfailure)

name = 'fates_mort_hf_sm_threshold'
call fates_params%RetreiveParameterAllocate(name=name, &
data=this%hf_sm_threshold)
Expand Down Expand Up @@ -1479,6 +1506,9 @@ subroutine FatesReportPFTParams(is_master)
write(fates_log(),fmt0) 'grperc = ',EDPftvarcon_inst%grperc
write(fates_log(),fmt0) 'c2b = ',EDPftvarcon_inst%c2b
write(fates_log(),fmt0) 'bmort = ',EDPftvarcon_inst%bmort
write(fates_log(),fmt0) 'mort_scalar_coldstress = ',EDPftvarcon_inst%mort_scalar_coldstress
write(fates_log(),fmt0) 'mort_scalar_cstarvation = ',EDPftvarcon_inst%mort_scalar_cstarvation
write(fates_log(),fmt0) 'mort_scalar_hydrfailure = ',EDPftvarcon_inst%mort_scalar_hydrfailure
write(fates_log(),fmt0) 'hf_sm_threshold = ',EDPftvarcon_inst%hf_sm_threshold
write(fates_log(),fmt0) 'vcmaxha = ',EDPftvarcon_inst%vcmaxha
write(fates_log(),fmt0) 'jmaxha = ',EDPftvarcon_inst%jmaxha
Expand Down
Loading