Skip to content

Commit

Permalink
Move entrainment fact to climaparams method call
Browse files Browse the repository at this point in the history
  • Loading branch information
charleskawczynski committed Aug 8, 2021
1 parent 81c9f0f commit 5cfd332
Show file tree
Hide file tree
Showing 4 changed files with 16 additions and 16 deletions.
2 changes: 2 additions & 0 deletions integration_tests/utils/parameter_set.jl
Original file line number Diff line number Diff line change
Expand Up @@ -9,13 +9,15 @@ end

CLIMAParameters.Planet.MSLP(ps::EarthParameterSet) = ps.nt.MSLP
CLIMAParameters.Atmos.Microphysics.τ_cond_evap(ps::EarthParameterSet) = ps.nt.τ_cond_evap
CLIMAParameters.Atmos.EDMF.c_ε(ps::EarthParameterSet) = ps.nt.c_ε

#! format: off
function create_parameter_set(namelist)
TC = TurbulenceConvection
nt = (;
MSLP = 100000.0, # or grab from, e.g., namelist[""][...]
τ_cond_evap = TC.parse_namelist(namelist, "microphysics", "τ_cond_evap"; default = 10.0),
c_ε = TC.parse_namelist(namelist, "turbulence", "EDMF_PrognosticTKE", "entrainment_factor"),
)
return EarthParameterSet(nt)
end
Expand Down
4 changes: 2 additions & 2 deletions src/Turbulence_PrognosticTKE.jl
Original file line number Diff line number Diff line change
Expand Up @@ -934,6 +934,7 @@ function compute_updraft_closures(self::EDMF_PrognosticTKE, GMV::GridMeanVariabl
sa = eos_struct()
quadrature_order = 3
grid = get_grid(self)
param_set = parameter_set(self)

upd_cloud_diagnostics(self.UpdVar, reference_state(self))

Expand All @@ -942,7 +943,6 @@ function compute_updraft_closures(self::EDMF_PrognosticTKE, GMV::GridMeanVariabl
input.dz = get_grid(self).dz
input.zbl = compute_zbl_qt_grad(GMV)
input.sort_pow = self.sorting_power
input.c_ent = self.entrainment_factor
input.c_det = self.detrainment_factor
input.c_div = self.entrainment_Mdiv_factor
input.c_mu = self.entrainment_sigma
Expand Down Expand Up @@ -999,7 +999,7 @@ function compute_updraft_closures(self::EDMF_PrognosticTKE, GMV::GridMeanVariabl
input.dMdz = (input.M - Mm) * grid.dzi
input.tke = self.EnvVar.TKE.values[k]

ret = self.entr_detr_fp(input)
ret = self.entr_detr_fp(param_set, input)
self.entr_sc[i, k] = ret.entr_sc
self.detr_sc[i, k] = ret.detr_sc
self.sorting_function[i, k] = ret.sorting_function
Expand Down
23 changes: 12 additions & 11 deletions src/turbulence_functions.jl
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
# Entrainment Rates

function entr_detr_env_moisture_deficit_b_ED_MF(entr_in::entr_struct)
function entr_detr_env_moisture_deficit_b_ED_MF(param_set::APS, entr_in::entr_struct)
_ret = entr_struct()

moisture_deficit_d =
Expand Down Expand Up @@ -39,13 +39,14 @@ function entr_detr_env_moisture_deficit_b_ED_MF(entr_in::entr_struct)
fabs(entr_in.buoy_ed_flux) /
(fabs(entr_in.a_upd * entr_in.a_env * (entr_in.w_upd - entr_in.w_env) * (entr_in.b_upd - entr_in.b_env)) + 1e-8)
logistic_e *= (1.0 / (1.0 + exp(entr_in.c_ed_mf * (ed_mf_ratio - 1.0))))
_ret.entr_sc = inv_timescale / dw * (entr_in.c_ent * logistic_e + c_det * moisture_deficit_e)
_ret.detr_sc = inv_timescale / dw * (entr_in.c_ent * logistic_d + c_det * moisture_deficit_d)
c_ε = CPEDMF.c_ε(param_set)
_ret.entr_sc = inv_timescale / dw * (c_ε * logistic_e + c_det * moisture_deficit_e)
_ret.detr_sc = inv_timescale / dw * (c_ε * logistic_d + c_det * moisture_deficit_d)

return _ret
end

function entr_detr_env_moisture_deficit(entr_in::entr_in_struct)
function entr_detr_env_moisture_deficit(param_set::APS, entr_in::entr_in_struct)
_ret = entr_struct()
l = zeros(2)

Expand Down Expand Up @@ -83,13 +84,14 @@ function entr_detr_env_moisture_deficit(entr_in::entr_in_struct)
l[1] = entr_in.tke_coef * fabs(db / sqrt(entr_in.tke + 1e-8))
l[2] = fabs(db / dw)
inv_timescale = lamb_smooth_minimum(l, 0.1, 0.0005)
_ret.entr_sc = inv_timescale / dw * (entr_in.c_ent * logistic_e + c_det * moisture_deficit_e)
_ret.detr_sc = inv_timescale / dw * (entr_in.c_ent * logistic_d + c_det * moisture_deficit_d)
c_ε = CPEDMF.c_ε(param_set)
_ret.entr_sc = inv_timescale / dw * (c_ε * logistic_e + c_det * moisture_deficit_e)
_ret.detr_sc = inv_timescale / dw * (c_ε * logistic_d + c_det * moisture_deficit_d)

return _ret
end

function entr_detr_env_moisture_deficit_div(entr_in::entr_struct)
function entr_detr_env_moisture_deficit_div(param_set::APS, entr_in::entr_struct)
_ret = entr_struct()
l = zeros(2)

Expand Down Expand Up @@ -132,10 +134,9 @@ function entr_detr_env_moisture_deficit_div(entr_in::entr_struct)
l[2] = fabs(db / dw)
inv_timescale = lamb_smooth_minimum(l, 0.1, 0.0005)

_ret.entr_sc =
inv_timescale / dw * (entr_in.c_ent * logistic_e + c_det * moisture_deficit_e) + entr_MdMdz * entr_in.c_div
_ret.detr_sc =
inv_timescale / dw * (entr_in.c_ent * logistic_d + c_det * moisture_deficit_d) + detr_MdMdz * entr_in.c_div
c_ε = CPEDMF.c_ε(param_set)
_ret.entr_sc = inv_timescale / dw * (c_ε * logistic_e + c_det * moisture_deficit_e) + entr_MdMdz * entr_in.c_div
_ret.detr_sc = inv_timescale / dw * (c_ε * logistic_d + c_det * moisture_deficit_d) + detr_MdMdz * entr_in.c_div


return _ret
Expand Down
3 changes: 0 additions & 3 deletions src/types.jl
Original file line number Diff line number Diff line change
Expand Up @@ -73,7 +73,6 @@ Base.@kwdef mutable struct entr_in_struct
dz::Float64 = 0
w_upd::Float64 = 0
b_upd::Float64 = 0
c_ent::Float64 = 0
dt::Float64 = 0
b_env::Float64 = 0
a_upd::Float64 = 0
Expand Down Expand Up @@ -761,7 +760,6 @@ mutable struct EDMF_PrognosticTKE{PS, A1, A2}
extrapolate_buoyancy::Bool
surface_area::Float64
max_area::Float64
entrainment_factor::Float64
entrainment_Mdiv_factor::Float64
updraft_mixing_frac::Float64
entrainment_sigma::Float64
Expand Down Expand Up @@ -1077,7 +1075,6 @@ mutable struct EDMF_PrognosticTKE{PS, A1, A2}
extrapolate_buoyancy,
surface_area,
max_area,
entrainment_factor,
entrainment_Mdiv_factor,
updraft_mixing_frac,
entrainment_sigma,
Expand Down

0 comments on commit 5cfd332

Please sign in to comment.