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

(*)Internal thickness variable unit correction #139

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
34 changes: 17 additions & 17 deletions src/parameterizations/vertical/MOM_bulk_mixed_layer.F90
Original file line number Diff line number Diff line change
Expand Up @@ -230,7 +230,7 @@ subroutine bulkmixedlayer(h_3d, u_3d, v_3d, tv, fluxes, dt, ea, eb, G, GV, US, C
integer, dimension(SZI_(G),SZK_(GV)) :: &
ksort ! The sorted k-index that each original layer goes to.
real, dimension(SZI_(G),SZJ_(G)) :: &
h_miss ! The summed absolute mismatch [Z ~> m].
h_miss ! The summed absolute mismatch [H ~> m or kg m-2].
real, dimension(SZI_(G)) :: &
TKE, & ! The turbulent kinetic energy available for mixing over a
! time step [Z L2 T-2 ~> m3 s-2].
Expand Down Expand Up @@ -299,12 +299,12 @@ subroutine bulkmixedlayer(h_3d, u_3d, v_3d, tv, fluxes, dt, ea, eb, G, GV, US, C
! adjustment [Z L2 T-2 ~> m3 s-2].
real, dimension(SZI_(G),SZJ_(G)) :: &
Hsfc_max, & ! The thickness of the surface region (mixed and buffer layers)
! after entrainment but before any buffer layer detrainment [Z ~> m].
! after entrainment but before any buffer layer detrainment [H ~> m or kg m-2].
Hsfc_used, & ! The thickness of the surface region after buffer layer
! detrainment [Z ~> m].
! detrainment [H ~> m or kg m-2].
Hsfc_min, & ! The minimum thickness of the surface region based on the
! new mixed layer depth and the previous thickness of the
! neighboring water columns [Z ~> m].
! neighboring water columns [H ~> m or kg m-2].
h_sum, & ! The total thickness of the water column [H ~> m or kg m-2].
hmbl_prev ! The previous thickness of the mixed and buffer layers [H ~> m or kg m-2].
real, dimension(SZI_(G)) :: &
Expand Down Expand Up @@ -538,7 +538,7 @@ subroutine bulkmixedlayer(h_3d, u_3d, v_3d, tv, fluxes, dt, ea, eb, G, GV, US, C
h(i,0) = htot(i)
endif ; enddo
if (write_diags .and. allocated(CS%ML_depth)) then ; do i=is,ie
CS%ML_depth(i,j) = h(i,0) * GV%H_to_m ! Rescale the diagnostic.
CS%ML_depth(i,j) = h(i,0) ! Store the diagnostic.
enddo ; endif
if (associated(Hml)) then ; do i=is,ie
Hml(i,j) = G%mask2dT(i,j) * (h(i,0) * GV%H_to_Z) ! Rescale the diagnostic for output.
Expand Down Expand Up @@ -573,14 +573,14 @@ subroutine bulkmixedlayer(h_3d, u_3d, v_3d, tv, fluxes, dt, ea, eb, G, GV, US, C
hmbl_prev(i,j-1) - dHD*min(h_sum(i,j),h_sum(i,j-1)), &
hmbl_prev(i,j+1) - dHD*min(h_sum(i,j),h_sum(i,j+1))) )

Hsfc_min(i,j) = GV%H_to_Z * max(h(i,0), min(Hsfc(i), H_nbr))
Hsfc_min(i,j) = max(h(i,0), min(Hsfc(i), H_nbr))

if (CS%limit_det) max_BL_det(i) = max(0.0, Hsfc(i)-H_nbr)
enddo
endif

if (CS%id_Hsfc_max > 0) then ; do i=is,ie
Hsfc_max(i,j) = GV%H_to_Z * Hsfc(i)
Hsfc_max(i,j) = Hsfc(i)
enddo ; endif
endif

Expand All @@ -601,9 +601,9 @@ subroutine bulkmixedlayer(h_3d, u_3d, v_3d, tv, fluxes, dt, ea, eb, G, GV, US, C
endif

if (CS%id_Hsfc_used > 0) then
do i=is,ie ; Hsfc_used(i,j) = GV%H_to_Z * h(i,0) ; enddo
do i=is,ie ; Hsfc_used(i,j) = h(i,0) ; enddo
do k=CS%nkml+1,nkmb ; do i=is,ie
Hsfc_used(i,j) = Hsfc_used(i,j) + GV%H_to_Z * h(i,k)
Hsfc_used(i,j) = Hsfc_used(i,j) + h(i,k)
enddo ; enddo
endif

Expand Down Expand Up @@ -686,15 +686,15 @@ subroutine bulkmixedlayer(h_3d, u_3d, v_3d, tv, fluxes, dt, ea, eb, G, GV, US, C

if (CS%id_h_mismatch > 0) then
do i=is,ie
h_miss(i,j) = GV%H_to_Z * abs(h_3d(i,j,1) - (h_orig(i,1) + &
h_miss(i,j) = abs(h_3d(i,j,1) - (h_orig(i,1) + &
(eaml(i,1) + (ebml(i,1) - eaml(i,1+1)))))
enddo
do k=2,nz-1 ; do i=is,ie
h_miss(i,j) = h_miss(i,j) + GV%H_to_Z * abs(h_3d(i,j,k) - (h_orig(i,k) + &
h_miss(i,j) = h_miss(i,j) + abs(h_3d(i,j,k) - (h_orig(i,k) + &
((eaml(i,k) - ebml(i,k-1)) + (ebml(i,k) - eaml(i,k+1)))))
enddo ; enddo
do i=is,ie
h_miss(i,j) = h_miss(i,j) + GV%H_to_Z * abs(h_3d(i,j,nz) - (h_orig(i,nz) + &
h_miss(i,j) = h_miss(i,j) + abs(h_3d(i,j,nz) - (h_orig(i,nz) + &
((eaml(i,nz) - ebml(i,nz-1)) + ebml(i,nz))))
enddo
endif
Expand Down Expand Up @@ -3501,7 +3501,7 @@ subroutine bulkmixedlayer_init(Time, G, GV, US, param_file, diag, CS)
"during mixedlayer convection.", default=.false.)

CS%id_ML_depth = register_diag_field('ocean_model', 'h_ML', diag%axesT1, &
Time, 'Surface mixed layer depth', 'm')
Time, 'Surface mixed layer depth', 'm', conversion=GV%H_to_m)
CS%id_TKE_wind = register_diag_field('ocean_model', 'TKE_wind', diag%axesT1, &
Time, 'Wind-stirring source of mixed layer TKE', &
'm3 s-3', conversion=US%Z_to_m*(US%L_to_m**2)*(US%s_to_T**3))
Expand Down Expand Up @@ -3533,13 +3533,13 @@ subroutine bulkmixedlayer_init(Time, G, GV, US, param_file, diag, CS)
Time, 'Spurious source of potential energy from mixed layer only detrainment', &
'W m-2', conversion=US%RZ3_T3_to_W_m2*US%L_to_Z**2)
CS%id_h_mismatch = register_diag_field('ocean_model', 'h_miss_ML', diag%axesT1, &
Time, 'Summed absolute mismatch in entrainment terms', 'm', conversion=US%Z_to_m)
Time, 'Summed absolute mismatch in entrainment terms', 'm', conversion=GV%H_to_m)
CS%id_Hsfc_used = register_diag_field('ocean_model', 'Hs_used', diag%axesT1, &
Time, 'Surface region thickness that is used', 'm', conversion=US%Z_to_m)
Time, 'Surface region thickness that is used', 'm', conversion=GV%H_to_m)
CS%id_Hsfc_max = register_diag_field('ocean_model', 'Hs_max', diag%axesT1, &
Time, 'Maximum surface region thickness', 'm', conversion=US%Z_to_m)
Time, 'Maximum surface region thickness', 'm', conversion=GV%H_to_m)
CS%id_Hsfc_min = register_diag_field('ocean_model', 'Hs_min', diag%axesT1, &
Time, 'Minimum surface region thickness', 'm', conversion=US%Z_to_m)
Time, 'Minimum surface region thickness', 'm', conversion=GV%H_to_m)
!CS%lim_det_dH_sfc = 0.5 ; CS%lim_det_dH_bathy = 0.2 ! Technically these should not get used if limit_det is false?
if (CS%limit_det .or. (CS%id_Hsfc_min > 0)) then
call get_param(param_file, mdl, "LIMIT_BUFFER_DET_DH_SFC", CS%lim_det_dH_sfc, &
Expand Down
31 changes: 16 additions & 15 deletions src/parameterizations/vertical/MOM_diapyc_energy_req.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1046,7 +1046,8 @@ subroutine find_PE_chg(Kddt_h0, dKddt_h, hp_a, hp_b, Th_a, Sh_a, Th_b, Sh_b, &
real :: ColHt_core ! The diffusivity-independent core term in the expressions
! for the column height changes [R L2 T-2 ~> J m-3].
real :: ColHt_chg ! The change in the column height [Z ~> m].
real :: y1 ! A local temporary term, in [H-3] or [H-4] in various contexts.
real :: y1_3 ! A local temporary term in [H-3 ~> m-3 or m6 kg-3].
real :: y1_4 ! A local temporary term in [H-4 ~> m-4 or m8 kg-4].

! The expression for the change in potential energy used here is derived
! from the expression for the final estimates of the changes in temperature
Expand All @@ -1068,37 +1069,37 @@ subroutine find_PE_chg(Kddt_h0, dKddt_h, hp_a, hp_b, Th_a, Sh_a, Th_b, Sh_b, &
if (present(PE_chg)) then
! Find the change in column potential energy due to the change in the
! diffusivity at this interface by dKddt_h.
y1 = dKddt_h / (bdt1 * (bdt1 + dKddt_h * hps))
PE_chg = PEc_core * y1
ColHt_chg = ColHt_core * y1
y1_3 = dKddt_h / (bdt1 * (bdt1 + dKddt_h * hps))
PE_chg = PEc_core * y1_3
ColHt_chg = ColHt_core * y1_3
if (ColHt_chg < 0.0) PE_chg = PE_chg - pres_Z * ColHt_chg
if (present(ColHt_cor)) ColHt_cor = -pres_Z * min(ColHt_chg, 0.0)
elseif (present(ColHt_cor)) then
y1 = dKddt_h / (bdt1 * (bdt1 + dKddt_h * hps))
ColHt_cor = -pres_Z * min(ColHt_core * y1, 0.0)
y1_3 = dKddt_h / (bdt1 * (bdt1 + dKddt_h * hps))
ColHt_cor = -pres_Z * min(ColHt_core * y1_3, 0.0)
endif

if (present(dPEc_dKd)) then
! Find the derivative of the potential energy change with dKddt_h.
y1 = 1.0 / (bdt1 + dKddt_h * hps)**2
dPEc_dKd = PEc_core * y1
ColHt_chg = ColHt_core * y1
y1_4 = 1.0 / (bdt1 + dKddt_h * hps)**2
dPEc_dKd = PEc_core * y1_4
ColHt_chg = ColHt_core * y1_4
if (ColHt_chg < 0.0) dPEc_dKd = dPEc_dKd - pres_Z * ColHt_chg
endif

if (present(dPE_max)) then
! This expression is the limit of PE_chg for infinite dKddt_h.
y1 = 1.0 / (bdt1 * hps)
dPE_max = PEc_core * y1
ColHt_chg = ColHt_core * y1
y1_3 = 1.0 / (bdt1 * hps)
dPE_max = PEc_core * y1_3
ColHt_chg = ColHt_core * y1_3
if (ColHt_chg < 0.0) dPE_max = dPE_max - pres_Z * ColHt_chg
endif

if (present(dPEc_dKd_0)) then
! This expression is the limit of dPEc_dKd for dKddt_h = 0.
y1 = 1.0 / bdt1**2
dPEc_dKd_0 = PEc_core * y1
ColHt_chg = ColHt_core * y1
y1_4 = 1.0 / bdt1**2
dPEc_dKd_0 = PEc_core * y1_4
ColHt_chg = ColHt_core * y1_4
if (ColHt_chg < 0.0) dPEc_dKd_0 = dPEc_dKd_0 - pres_Z * ColHt_chg
endif

Expand Down
50 changes: 26 additions & 24 deletions src/parameterizations/vertical/MOM_energetic_PBL.F90
Original file line number Diff line number Diff line change
Expand Up @@ -76,7 +76,7 @@ module MOM_energetic_PBL
!! boundary layer thickness. The default is 0, but a
!! value of 0.1 might be better justified by observations.
real :: MLD_tol !< A tolerance for determining the boundary layer thickness when
!! Use_MLD_iteration is true [Z ~> m].
!! Use_MLD_iteration is true [H ~> m or kg m-2].
real :: min_mix_len !< The minimum mixing length scale that will be used by ePBL [Z ~> m].
!! The default (0) does not set a minimum.

Expand Down Expand Up @@ -634,9 +634,9 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
real :: dt_h ! The timestep divided by the averages of the thicknesses around
! a layer, times a thickness conversion factor [H T Z-2 ~> s m-1 or kg s m-4].
real :: h_bot ! The distance from the bottom [H ~> m or kg m-2].
real :: h_rsum ! The running sum of h from the top [Z ~> m].
real :: h_rsum ! The running sum of h from the top [H ~> m or kg m-2].
real :: I_hs ! The inverse of h_sum [H-1 ~> m-1 or m2 kg-1].
real :: I_MLD ! The inverse of the current value of MLD [Z-1 ~> m-1].
real :: I_MLD ! The inverse of the current value of MLD [H-1 ~> m-1 or m2 kg-1].
real :: h_tt ! The distance from the surface or up to the next interface
! that did not exhibit turbulent mixing from this scheme plus
! a surface mixing roughness length given by h_tt_min [H ~> m or kg m-2].
Expand All @@ -648,7 +648,7 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
real :: vstar ! An in-situ turbulent velocity [Z T-1 ~> m s-1].
real :: mstar_total ! The value of mstar used in ePBL [nondim]
real :: mstar_LT ! An addition to mstar due to Langmuir turbulence [nondim] (output for diagnostic)
real :: MLD_output ! The mixed layer depth output from this routine [Z ~> m].
real :: MLD_output ! The mixed layer depth output from this routine [H ~> m or kg m-2].
real :: LA ! The value of the Langmuir number [nondim]
real :: LAmod ! The modified Langmuir number by convection [nondim]
real :: hbs_here ! The local minimum of hb_hs and MixLen_shape, times a
Expand Down Expand Up @@ -706,8 +706,9 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
!----------------------------------------------------------------------
!/BGR added Aug24,2016 for adding iteration to get boundary layer depth
! - needed to compute new mixing length.
real :: MLD_guess, MLD_found ! Mixing Layer depth guessed/found for iteration [Z ~> m].
real :: min_MLD ! Iteration bounds [Z ~> m], which are adjusted at each step
real :: MLD_guess, MLD_found ! Mixing Layer depth guessed/found for iteration [H ~> m or kg m-2].
real :: MLD_guess_Z ! A guessed mixed layer depth, converted to height units [Z ~> m]
real :: min_MLD ! Iteration bounds [H ~> m or kg m-2], which are adjusted at each step
real :: max_MLD ! - These are initialized based on surface/bottom
! 1. The iteration guesses a value (possibly from prev step or neighbor).
! 2. The iteration checks if value is converged, too shallow, or too deep.
Expand All @@ -720,8 +721,8 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
! manner giving a usable guess. When it does fail, it is due to convection
! within the boundary layer. Likely, a new method e.g. surface_disconnect,
! can improve this.
real :: dMLD_min ! The change in diagnosed mixed layer depth when the guess is min_MLD [Z ~> m]
real :: dMLD_max ! The change in diagnosed mixed layer depth when the guess is max_MLD [Z ~> m]
real :: dMLD_min ! The change in diagnosed mixed layer depth when the guess is min_MLD [H ~> m or kg m-2]
real :: dMLD_max ! The change in diagnosed mixed layer depth when the guess is max_MLD [H ~> m or kg m-2]
logical :: OBL_converged ! Flag for convergence of MLD
integer :: OBL_it ! Iteration counter

Expand Down Expand Up @@ -754,7 +755,7 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
I_dtrho = 0.0 ; if (dt*GV%Rho0 > 0.0) I_dtrho = (US%Z_to_m**3*US%s_to_T**3) / (dt*GV%Rho0)
vstar_unit_scale = US%m_to_Z * US%T_to_s

MLD_guess = MLD_io
MLD_guess = MLD_io*GV%Z_to_H

! Determine the initial mech_TKE and conv_PErel, including the energy required
! to mix surface heating through the topmost cell, the energy released by mixing
Expand Down Expand Up @@ -787,15 +788,15 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
hb_hs(K) = h_bot * I_hs
enddo

MLD_output = h(1)*GV%H_to_Z
MLD_output = h(1)

!/The following lines are for the iteration over MLD
! max_MLD will initialized as ocean bottom depth
max_MLD = 0.0 ; do k=1,nz ; max_MLD = max_MLD + h(k)*GV%H_to_Z ; enddo
max_MLD = 0.0 ; do k=1,nz ; max_MLD = max_MLD + h(k) ; enddo
! min_MLD will be initialized to 0.
min_MLD = 0.0
! Set values of the wrong signs to indicate that these changes are not based on valid estimates
dMLD_min = -1.0*US%m_to_Z ; dMLD_max = 1.0*US%m_to_Z
dMLD_min = -1.0*GV%m_to_H ; dMLD_max = 1.0*GV%m_to_H

! If no first guess is provided for MLD, try the middle of the water column
if (MLD_guess <= min_MLD) MLD_guess = 0.5 * (min_MLD + max_MLD)
Expand All @@ -811,18 +812,19 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
if (debug) then ; mech_TKE_k(:) = 0.0 ; conv_PErel_k(:) = 0.0 ; endif

! Reset ML_depth
MLD_output = h(1)*GV%H_to_Z
MLD_output = h(1)
sfc_connected = .true.

!/ Here we get MStar, which is the ratio of convective TKE driven mixing to UStar**3
MLD_guess_z = GV%H_to_Z*MLD_guess ! Convert MLD from thickness to height coordinates for these calls
if (CS%Use_LT) then
call get_Langmuir_Number(LA, G, GV, US, abs(MLD_guess), u_star_mean, i, j, h, Waves, &
call get_Langmuir_Number(LA, G, GV, US, abs(MLD_guess_z), u_star_mean, i, j, h, Waves, &
U_H=u, V_H=v)
call find_mstar(CS, US, B_flux, u_star, u_star_Mean, MLD_Guess, absf, &
call find_mstar(CS, US, B_flux, u_star, u_star_Mean, MLD_guess_z, absf, &
MStar_total, Langmuir_Number=La, Convect_Langmuir_Number=LAmod,&
mstar_LT=mstar_LT)
else
call find_mstar(CS, US, B_flux, u_star, u_star_mean, MLD_guess, absf, mstar_total)
call find_mstar(CS, US, B_flux, u_star, u_star_mean, MLD_guess_z, absf, mstar_total)
endif

!/ Apply MStar to get mech_TKE
Expand Down Expand Up @@ -879,7 +881,7 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
h_rsum = 0.0
MixLen_shape(1) = 1.0
do K=2,nz+1
h_rsum = h_rsum + h(k-1)*GV%H_to_Z
h_rsum = h_rsum + h(k-1)
if (CS%MixLenExponent==2.0) then
MixLen_shape(K) = CS%transLay_scale + (1.0 - CS%transLay_scale) * &
(max(0.0, (MLD_guess - h_rsum)*I_MLD) )**2 ! CS%MixLenExponent
Expand Down Expand Up @@ -1076,7 +1078,7 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
if (CS%wT_scheme==wT_from_cRoot_TKE) then
vstar = CS%vstar_scale_fac * vstar_unit_scale * (I_dtrho*TKE_here)**C1_3
elseif (CS%wT_scheme==wT_from_RH18) then
Surface_Scale = max(0.05, 1.0 - htot * GV%H_to_Z / MLD_guess)
Surface_Scale = max(0.05, 1.0 - htot / MLD_guess)
vstar = CS%vstar_scale_fac * Surface_Scale * (CS%vstar_surf_fac*u_star + &
vstar_unit_scale * (CS%wstar_ustar_coef*conv_PErel*I_dtrho)**C1_3)
endif
Expand Down Expand Up @@ -1125,7 +1127,7 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
if (CS%wT_scheme==wT_from_cRoot_TKE) then
vstar = CS%vstar_scale_fac * vstar_unit_scale * (I_dtrho*TKE_here)**C1_3
elseif (CS%wT_scheme==wT_from_RH18) then
Surface_Scale = max(0.05, 1. - htot * GV%H_to_Z / MLD_guess)
Surface_Scale = max(0.05, 1. - htot / MLD_guess)
vstar = CS%vstar_scale_fac * Surface_Scale * (CS%vstar_surf_fac*u_star + &
vstar_unit_scale * (CS%wstar_ustar_coef*conv_PErel*I_dtrho)**C1_3)
endif
Expand Down Expand Up @@ -1178,7 +1180,7 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
eCD%dTKE_MKE = eCD%dTKE_MKE + MKE_src * I_dtdiag
endif
if (sfc_connected) then
MLD_output = MLD_output + GV%H_to_Z * h(k)
MLD_output = MLD_output + h(k)
endif

Kddt_h(K) = Kd(K) * dt_h
Expand All @@ -1202,7 +1204,7 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
mech_TKE = TKE_reduc*(mech_TKE + MKE_src)
conv_PErel = TKE_reduc*conv_PErel
if (sfc_connected) then
MLD_output = MLD_output + GV%H_to_Z * h(k)
MLD_output = MLD_output + h(k)
endif

elseif (tot_TKE == 0.0) then
Expand Down Expand Up @@ -1303,7 +1305,7 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
endif

if (sfc_connected) MLD_output = MLD_output + &
(PE_chg / (PE_chg_g0)) * GV%H_to_Z * h(k)
(PE_chg / (PE_chg_g0)) * h(k)

tot_TKE = 0.0 ; mech_TKE = 0.0 ; conv_PErel = 0.0
sfc_disconnect = .true.
Expand Down Expand Up @@ -1422,7 +1424,7 @@ subroutine ePBL_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, abs
eCD%LA = 0.0 ; eCD%LAmod = 0.0 ; eCD%mstar = mstar_total ; eCD%mstar_LT = 0.0
endif

MLD_io = MLD_output
MLD_io = GV%H_to_Z*MLD_output

end subroutine ePBL_column

Expand Down Expand Up @@ -2125,7 +2127,7 @@ subroutine energetic_PBL_init(Time, G, GV, US, param_file, diag, CS)
call get_param(param_file, mdl, "EPBL_MLD_TOLERANCE", CS%MLD_tol, &
"The tolerance for the iteratively determined mixed "//&
"layer depth. This is only used with USE_MLD_ITERATION.", &
units="meter", default=1.0, scale=US%m_to_Z, do_not_log=.not.CS%Use_MLD_iteration)
units="meter", default=1.0, scale=GV%m_to_H, do_not_log=.not.CS%Use_MLD_iteration)
call get_param(param_file, mdl, "EPBL_MLD_BISECTION", CS%MLD_bisection, &
"If true, use bisection with the iterative determination of the self-consistent "//&
"mixed layer depth. Otherwise use the false position after a maximum and minimum "//&
Expand Down
Loading