mom_energetic_pbl module reference

Energetically consistent planetary boundary layer parameterization.

More…

Data Types

energetic_pbl_cs

This control structure holds parameters for the MOM_energetic_PBL() module.

epbl_column_diags

A type for conveniently passing around ePBL diagnostics for a column.

Functions/Subroutines

energetic_pbl()

This subroutine determines the diffusivities from the integrated energetics mixed layer model.

epbl_column()

This subroutine determines the diffusivities from the integrated energetics mixed layer model for a single column of water.

find_pe_chg()

This subroutine calculates the change in potential energy and or derivatives for several changes in an interface’s diapycnal diffusivity times a timestep.

find_pe_chg_orig()

This subroutine calculates the change in potential energy and or derivatives for several changes in an interface’s diapycnal diffusivity times a timestep using the original form used in the first version of ePBL.

find_mstar()

This subroutine finds the Mstar value for ePBL.

mstar_langmuir()

This subroutine modifies the Mstar value if the Langmuir number is present.

energetic_pbl_get_mld()

Copies the ePBL active mixed layer depth into MLD, in units of [Z ~> m] unless other units are specified.

energetic_pbl_init()

This subroutine initializes the energetic_PBL module.

energetic_pbl_end()

Clean up and deallocate memory associated with the energetic_PBL module.

Detailed Description

Energetically consistent planetary boundary layer parameterization.

Type Documentation

type mom_energetic_pbl/energetic_pbl_cs

This control structure holds parameters for the MOM_energetic_PBL() module. module.

Type fields:
  • % id_ml_depth [integer] :: Diagnostic IDs.

  • % id_hml_depth [integer] :: Diagnostic IDs.

  • % id_tke_wind [integer] :: Diagnostic IDs.

  • % id_tke_mixing [integer] :: Diagnostic IDs.

  • % id_tke_mke [integer] :: Diagnostic IDs.

  • % id_tke_conv [integer] :: Diagnostic IDs.

  • % id_tke_forcing [integer] :: Diagnostic IDs.

  • % id_tke_mech_decay [integer] :: Diagnostic IDs.

  • % id_tke_conv_decay [integer] :: Diagnostic IDs.

  • % id_mixing_length [integer] :: Diagnostic IDs.

  • % id_velocity_scale [integer] :: Diagnostic IDs.

  • % id_mstar_mix [integer] :: Diagnostic IDs.

  • % id_la_mod [integer] :: Diagnostic IDs.

  • % id_la [integer] :: Diagnostic IDs.

  • % id_mstar_lt [integer] :: Diagnostic IDs.

  • % initialized [logical] :: True if this control structure has been initialized.

  • % vonkar [real] :: The von Karman coefficient as used in the ePBL module [nondim].

  • % omega [real] :: The Earth’s rotation rate [T-1 ~> s-1].

  • % omega_frac [real] :: When setting the decay scale for turbulence, use this fraction of the absolute rotation rate blended with the local value of f, as sqrt((1-omega_frac)*f^2 + omega_frac*4*omega^2) [nondim].

  • % nstar [real] :: The fraction of the TKE input to the mixed layer available to drive entrainment [nondim]. This quantity is the vertically integrated buoyancy production minus the vertically integrated dissipation of TKE produced by buoyancy.

  • % use_mld_iteration [logical] :: If true, use the proximity to the bottom of the actively turbulent surface boundary layer to constrain the mixing lengths.

  • % mld_iteration_guess [logical] :: False to default to guessing half the ocean depth for the first iteration.

  • % mld_bisection [logical] :: 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 bound have been evaluated and the returned value from the previous guess or bisection before this.

  • % max_mld_its [integer] :: The maximum number of iterations that can be used to find a self-consistent mixed layer depth with Use_MLD_iteration.

  • % mixlenexponent [real] :: Exponent in the mixing length shape-function [nondim]. 1 is law-of-the-wall at top and bottom, 2 is more KPP like.

  • % mke_to_tke_effic [real] :: The efficiency with which mean kinetic energy released by mechanically forced entrainment of the mixed layer is converted to TKE [nondim].

  • % ustar_min [real] :: A minimum value of ustar to avoid numerical problems [Z T-1 ~> m s-1]. If the value is small enough, this should not affect the solution.

  • % ekman_scale_coef [real] :: A nondimensional scaling factor controlling the inhibition of the diffusive length scale by rotation [nondim]. Making this larger decreases the diffusivity in the planetary boundary layer.

  • % translay_scale [real] :: A scale for the mixing length in the transition layer at the edge of the boundary layer as a fraction of the boundary layer thickness [nondim]. The default is 0, but a value of 0.1 might be better justified by observations.

  • % mld_tol [real] :: A tolerance for determining the boundary layer thickness when Use_MLD_iteration is true [H ~> m or kg m-2].

  • % min_mix_len [real] :: The minimum mixing length scale that will be used by ePBL [Z ~> m]. The default (0) does not set a minimum.

  • % wt_scheme [integer] :: An enumerated value indicating the method for finding the turbulent velocity scale. There are currently two options: wT_mwT_from_cRoot_TKE is the original (TKE_remaining)^1/3 wT_from_RH18 is the version described by Reichl and Hallberg, 2018.

  • % wstar_ustar_coef [real] :: A ratio relating the efficiency with which convectively released energy is converted to a turbulent velocity, relative to mechanically forced turbulent kinetic energy [nondim]. Making this larger increases the diffusivity.

  • % vstar_surf_fac [real] :: If (wT_scheme == wT_from_RH18) this is the proportionality coefficient between ustar and the surface mechanical contribution to vstar [nondim].

  • % vstar_scale_fac [real] :: An overall nondimensional scaling factor for vstar times a unit conversion factor [Z s T-1 m-1 ~> nondim]. Making this larger increases the diffusivity.

  • % mstar_scheme [integer] :: An encoded integer to determine which formula is used to set mstar.

  • % mstar_flatcap [logical] :: Set false to use asymptotic mstar cap.

  • % mstar_cap [real] :: Since MSTAR is restoring undissipated energy to mixing, there must be a cap on how large it can be [nondim]. This is definitely a function of latitude (Ekman limit), but will be taken as constant for now.

  • % tke_decay [real] :: The ratio of the natural Ekman depth to the TKE decay scale [nondim].

  • % fixed_mstar [real] :: Mstar is the ratio of the friction velocity cubed to the TKE available to drive entrainment, nondimensional. This quantity is the vertically integrated shear production minus the vertically integrated dissipation of TKE produced by shear. This value is used if the option for using a fixed mstar is used.

  • % c_ek [real] :: MSTAR Coefficient in rotation limit for mstar_scheme=OM4 [nondim].

  • % mstar_coef [real] :: MSTAR coefficient in rotation/stabilizing balance for mstar_scheme=OM4 [nondim].

  • % rh18_mstar_cn1 [real] :: MSTAR_N coefficient 1 (outer-most coefficient for fit) [nondim]. Value of 0.275 in RH18. Increasing this coefficient increases mechanical mixing for all values of Hf/ust, but is most effective at low values (weakly developed OSBLs).

  • % rh18_mstar_cn2 [real] :: MSTAR_N coefficient 2 (coefficient outside of exponential decay) [nondim]. Value of 8.0 in RH18. Increasing this coefficient increases MSTAR for all values of HF/ust, with a consistent affect across a wide range of Hf/ust.

  • % rh18_mstar_cn3 [real] :: MSTAR_N coefficient 3 (exponential decay coefficient) [nondim]. Value of -5.0 in RH18. Increasing this increases how quickly the value of MSTAR decreases as Hf/ust increases.

  • % rh18_mstar_cs1 [real] :: MSTAR_S coefficient for RH18 in stabilizing limit [nondim]. Value of 0.2 in RH18.

  • % rh18_mstar_cs2 [real] :: MSTAR_S exponent for RH18 in stabilizing limit [nondim]. Value of 0.4 in RH18.

  • % mstar_convect_coef [real] :: Factor to reduce mstar when statically unstable [nondim].

  • % use_lt [logical] :: Flag for using LT in Energy calculation.

  • % lt_enhance_form [integer] :: Integer for Enhancement functional form (various options)

  • % lt_enhance_coef [real] :: Coefficient in fit for Langmuir Enhancement [nondim].

  • % lt_enhance_exp [real] :: Exponent in fit for Langmuir Enhancement [nondim].

  • % lac_mldoek [real] :: Coefficient for Langmuir number modification based on the ratio of the mixed layer depth over the Ekman depth [nondim].

  • % lac_mldoob_stab [real] :: Coefficient for Langmuir number modification based on the ratio of the mixed layer depth over the Obukhov depth with stabilizing forcing [nondim].

  • % lac_ekoob_stab [real] :: Coefficient for Langmuir number modification based on the ratio of the Ekman depth over the Obukhov depth with stabilizing forcing [nondim].

  • % lac_mldoob_un [real] :: Coefficient for Langmuir number modification based on the ratio of the mixed layer depth over the Obukhov depth with destabilizing forcing [nondim].

  • % lac_ekoob_un [real] :: Coefficient for Langmuir number modification based on the ratio of the Ekman depth over the Obukhov depth with destabilizing forcing [nondim].

  • % max_enhance_m [real] :: The maximum allowed LT enhancement to the mixing [nondim].

  • % time [type(time_type),pointer] :: A pointer to the ocean model’s clock.

  • % tke_diagnostics [logical] :: If true, diagnostics of the TKE budget are being calculated.

  • % answer_date [integer] :: The vintage of the order of arithmetic and expressions in the ePBL calculations. Values below 20190101 recover the answers from the end of 2018, while higher values use updated and more robust forms of the same expressions.

  • % orig_pe_calc [logical] :: If true, the ePBL code uses the original form of the potential energy change code. Otherwise, it uses a newer version that can work with successive increments to the diffusivity in upward or downward passes.

  • % diag [type( diag_ctrl ),pointer] :: A structure that is used to regulate the timing of diagnostic output.

  • % ml_depth [real(:,:),allocatable] :: The mixed layer depth determined by active mixing in ePBL [Z ~> m].

  • % diag_tke_wind [real(:,:),allocatable] :: The wind source of TKE [R Z3 T-3 ~> W m-2].

  • % diag_tke_mke [real(:,:),allocatable] :: The resolved KE source of TKE [R Z3 T-3 ~> W m-2].

  • % diag_tke_conv [real(:,:),allocatable] :: The convective source of TKE [R Z3 T-3 ~> W m-2].

  • % diag_tke_forcing [real(:,:),allocatable] :: The TKE sink required to mix surface penetrating shortwave heating.

  • % diag_tke_mech_decay [real(:,:),allocatable] :: The decay of mechanical TKE [R Z3 T-3 ~> W m-2].

  • % diag_tke_conv_decay [real(:,:),allocatable] :: The decay of convective TKE [R Z3 T-3 ~> W m-2].

  • % diag_tke_mixing [real(:,:),allocatable] :: The work done by TKE to deepen the mixed layer [R Z3 T-3 ~> W m-2].

  • % mstar_mix [real(:,:),allocatable] :: Mstar used in EPBL [nondim].

  • % mstar_lt [real(:,:),allocatable] :: Mstar due to Langmuir turbulence [nondim].

  • % la [real(:,:),allocatable] :: Langmuir number [nondim].

  • % la_mod [real(:,:),allocatable] :: Modified Langmuir number [nondim].

  • % sum_its [type( efp_type )(2)] :: The total number of iterations and columns worked on.

  • % velocity_scale [real(:,:,:),allocatable] :: The velocity scale used in getting Kd [Z T-1 ~> m s-1].

  • % mixing_length [real(:,:,:),allocatable] :: The length scale used in getting Kd [Z ~> m].

type mom_energetic_pbl/epbl_column_diags

A type for conveniently passing around ePBL diagnostics for a column.

Type fields:
  • % dtke_conv [real] :: Local column copies of energy change diagnostics, all in [R Z3 T-3 ~> W m-2].

  • % dtke_forcing [real] :: Local column copies of energy change diagnostics, all in [R Z3 T-3 ~> W m-2].

  • % dtke_wind [real] :: Local column copies of energy change diagnostics, all in [R Z3 T-3 ~> W m-2].

  • % dtke_mixing [real] :: Local column copies of energy change diagnostics, all in [R Z3 T-3 ~> W m-2].

  • % dtke_mke [real] :: Local column copies of energy change diagnostics, all in [R Z3 T-3 ~> W m-2].

  • % dtke_mech_decay [real] :: Local column copies of energy change diagnostics, all in [R Z3 T-3 ~> W m-2].

  • % dtke_conv_decay [real] :: Local column copies of energy change diagnostics, all in [R Z3 T-3 ~> W m-2].

  • % la [real] :: The value of the Langmuir number [nondim].

  • % lamod [real] :: The modified Langmuir number by convection [nondim].

  • % mstar [real] :: The value of mstar used in ePBL [nondim].

  • % mstar_lt [real] :: The portion of mstar due to Langmuir turbulence [nondim].

Function/Subroutine Documentation

subroutine mom_energetic_pbl/energetic_pbl(h_3d, u_3d, v_3d, tv, fluxes, dt, Kd_int, G, GV, US, CS, stoch_CS, dSV_dT, dSV_dS, TKE_forced, buoy_flux, Waves)

This subroutine determines the diffusivities from the integrated energetics mixed layer model. It assumes that heating, cooling and freshwater fluxes have already been applied. All calculations are done implicitly, and there is no stability limit on the time step.

Parameters:
  • g :: [inout] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

  • us :: [in] A dimensional unit scaling type

  • h_3d :: [inout] Layer thicknesses [H ~> m or kg m-2].

  • u_3d :: [in] Zonal velocities interpolated to h points

  • v_3d :: [in] Zonal velocities interpolated to h points

  • dsv_dt :: [in] The partial derivative of in-situ specific

  • dsv_ds :: [in] The partial derivative of in-situ specific

  • tke_forced :: [in] The forcing requirements to homogenize the

  • tv :: [inout] A structure containing pointers to any available thermodynamic fields. Absent fields have NULL ptrs.

  • fluxes :: [inout] A structure containing pointers to any possible forcing fields. Unused fields have NULL ptrs.

  • dt :: [in] Time increment [T ~> s].

  • kd_int :: [out] The diagnosed diffusivities at interfaces

  • cs :: [inout] Energetic PBL control structure

  • buoy_flux :: [in] The surface buoyancy flux [Z2 T-3 ~> m2 s-3].

  • waves :: Waves control structure for Langmuir turbulence

  • stoch_cs :: The control structure returned by a previous

Call to:

epbl_column mom_error_handler::mom_error

subroutine mom_energetic_pbl/epbl_column(h, u, v, T0, S0, dSV_dT, dSV_dS, TKE_forcing, B_flux, absf, u_star, u_star_mean, dt, MLD_io, Kd, mixvel, mixlen, GV, US, CS, eCD, Waves, G, i, j, TKE_gen_stoch, TKE_diss_stoch)

This subroutine determines the diffusivities from the integrated energetics mixed layer model for a single column of water.

Parameters:
  • gv :: [in] The ocean’s vertical grid structure.

  • us :: [in] A dimensional unit scaling type

  • h :: [in] Layer thicknesses [H ~> m or kg m-2].

  • u :: [in] Zonal velocities interpolated to h points [L T-1 ~> m s-1].

  • v :: [in] Zonal velocities interpolated to h points [L T-1 ~> m s-1].

  • t0 :: [in] The initial layer temperatures [C ~> degC].

  • s0 :: [in] The initial layer salinities [S ~> ppt].

  • dsv_dt :: [in] The partial derivative of in-situ specific volume with potential temperature [R-1 C-1 ~> m3 kg-1 degC-1].

  • dsv_ds :: [in] The partial derivative of in-situ specific volume with salinity [R-1 S-1 ~> m3 kg-1 ppt-1].

  • tke_forcing :: [in] The forcing requirements to homogenize the forcing that has been applied to each layer [R Z3 T-2 ~> J m-2].

  • b_flux :: [in] The surface buoyancy flux [Z2 T-3 ~> m2 s-3]

  • absf :: [in] The absolute value of the Coriolis parameter [T-1 ~> s-1].

  • u_star :: [in] The surface friction velocity [Z T-1 ~> m s-1].

  • u_star_mean :: [in] The surface friction velocity without any contribution from unresolved gustiness [Z T-1 ~> m s-1].

  • mld_io :: [inout] A first guess at the mixed layer depth on input, and the calculated mixed layer depth on output [Z ~> m].

  • dt :: [in] Time increment [T ~> s].

  • kd :: [out] The diagnosed diffusivities at interfaces

  • mixvel :: [out] The mixing velocity scale used in Kd

  • mixlen :: [out] The mixing length scale used in Kd [Z ~> m].

  • cs :: [inout] Energetic PBL control structure

  • ecd :: [inout] A container for passing around diagnostics.

  • waves :: Waves control structure for Langmuir turbulence

  • g :: [inout] The ocean’s grid structure.

  • tke_gen_stoch :: [in] random factor used to perturb TKE generation [nondim]

  • tke_diss_stoch :: [in] random factor used to perturb TKE dissipation [nondim]

  • i :: [in] The i-index to work on (used for Waves)

  • j :: [in] The i-index to work on (used for Waves)

Call to:

find_mstar find_pe_chg find_pe_chg_orig mom_wave_interface::get_langmuir_number mom_coms::real_to_efp report_avg_its use_fixed_mstar wt_from_croot_tke wt_from_rh18

Called from:

energetic_pbl

subroutine mom_energetic_pbl/find_pe_chg(Kddt_h0, dKddt_h, hp_a, hp_b, Th_a, Sh_a, Th_b, Sh_b, dT_to_dPE_a, dS_to_dPE_a, dT_to_dPE_b, dS_to_dPE_b, pres_Z, dT_to_dColHt_a, dS_to_dColHt_a, dT_to_dColHt_b, dS_to_dColHt_b, PE_chg, dPEc_dKd, dPE_max, dPEc_dKd_0, PE_ColHt_cor)

This subroutine calculates the change in potential energy and or derivatives for several changes in an interface’s diapycnal diffusivity times a timestep.

Parameters:
  • kddt_h0 :: [in] The previously used diffusivity at an interface times the time step and divided by the average of the thicknesses around the interface [H ~> m or kg m-2].

  • dkddt_h :: [in] The trial change in the diffusivity at an interface times the time step and divided by the average of the thicknesses around the interface [H ~> m or kg m-2].

  • hp_a :: [in] The effective pivot thickness of the layer above the interface, given by h_k plus a term that is a fraction (determined from the tridiagonal solver) of Kddt_h for the interface above [H ~> m or kg m-2].

  • hp_b :: [in] The effective pivot thickness of the layer below the interface, given by h_k plus a term that is a fraction (determined from the tridiagonal solver) of Kddt_h for the interface above [H ~> m or kg m-2].

  • th_a :: [in] An effective temperature times a thickness in the layer above, including implicit mixing effects with other yet higher layers [C H ~> degC m or degC kg m-2].

  • sh_a :: [in] An effective salinity times a thickness in the layer above, including implicit mixing effects with other yet higher layers [S H ~> ppt m or ppt kg m-2].

  • th_b :: [in] An effective temperature times a thickness in the layer below, including implicit mixing effects with other yet lower layers [C H ~> degC m or degC kg m-2].

  • sh_b :: [in] An effective salinity times a thickness in the layer below, including implicit mixing effects with other yet lower layers [S H ~> ppt m or ppt kg m-2].

  • dt_to_dpe_a :: [in] A factor (pres_lay*mass_lay*dSpec_vol/dT) relating a layer’s temperature change to the change in column potential energy, including all implicit diffusive changes in the temperatures of all the layers above [R Z3 T-2 C-1 ~> J m-2 degC-1].

  • ds_to_dpe_a :: [in] A factor (pres_lay*mass_lay*dSpec_vol/dS) relating a layer’s salinity change to the change in column potential energy, including all implicit diffusive changes in the salinities of all the layers above [R Z3 T-2 S-1 ~> J m-2 ppt-1].

  • dt_to_dpe_b :: [in] A factor (pres_lay*mass_lay*dSpec_vol/dT) relating a layer’s temperature change to the change in column potential energy, including all implicit diffusive changes in the temperatures of all the layers below [R Z3 T-2 C-1 ~> J m-2 degC-1].

  • ds_to_dpe_b :: [in] A factor (pres_lay*mass_lay*dSpec_vol/dS) relating a layer’s salinity change to the change in column potential energy, including all implicit diffusive changes in the salinities of all the layers below [R Z3 T-2 S-1 ~> J m-2 ppt-1].

  • pres_z :: [in] The rescaled hydrostatic interface pressure, which relates the changes in column thickness to the energy that is radiated as gravity waves and unavailable to drive mixing [R Z2 T-2 ~> J m-3].

  • dt_to_dcolht_a :: [in] A factor (mass_lay*dSColHtc_vol/dT) relating a layer’s temperature change to the change in column height, including all implicit diffusive changes in the temperatures of all the layers above [Z C-1 ~> m degC-1].

  • ds_to_dcolht_a :: [in] A factor (mass_lay*dSColHtc_vol/dS) relating a layer’s salinity change to the change in column height, including all implicit diffusive changes in the salinities of all the layers above [Z S-1 ~> m ppt-1].

  • dt_to_dcolht_b :: [in] A factor (mass_lay*dSColHtc_vol/dT) relating a layer’s temperature change to the change in column height, including all implicit diffusive changes in the temperatures of all the layers below [Z C-1 ~> m degC-1].

  • ds_to_dcolht_b :: [in] A factor (mass_lay*dSColHtc_vol/dS) relating a layer’s salinity change to the change in column height, including all implicit diffusive changes in the salinities of all the layers below [Z S-1 ~> m ppt-1].

  • pe_chg :: [out] The change in column potential energy from applying Kddt_h at the present interface [R Z3 T-2 ~> J m-2].

  • dpec_dkd :: [out] The partial derivative of PE_chg with Kddt_h [R Z3 T-2 H-1 ~> J m-3 or J kg-1].

  • dpe_max :: [out] The maximum change in column potential energy that could be realized by applying a huge value of Kddt_h at the present interface [R Z3 T-2 ~> J m-2].

  • dpec_dkd_0 :: [out] The partial derivative of PE_chg with Kddt_h in the limit where Kddt_h = 0 [R Z3 T-2 H-1 ~> J m-3 or J kg-1].

  • pe_colht_cor :: [out] The correction to PE_chg that is made due to a net change in the column height [R Z3 T-2 ~> J m-2].

Called from:

epbl_column

subroutine mom_energetic_pbl/find_pe_chg_orig(Kddt_h, h_k, b_den_1, dTe_term, dSe_term, dT_km1_t2, dS_km1_t2, dT_to_dPE_k, dS_to_dPE_k, dT_to_dPEa, dS_to_dPEa, pres_Z, dT_to_dColHt_k, dS_to_dColHt_k, dT_to_dColHta, dS_to_dColHta, PE_chg, dPEc_dKd, dPE_max, dPEc_dKd_0)

This subroutine calculates the change in potential energy and or derivatives for several changes in an interface’s diapycnal diffusivity times a timestep using the original form used in the first version of ePBL.

Parameters:
  • kddt_h :: [in] The diffusivity at an interface times the time step and divided by the average of the thicknesses around the interface [H ~> m or kg m-2].

  • h_k :: [in] The thickness of the layer below the interface [H ~> m or kg m-2].

  • b_den_1 :: [in] The first term in the denominator of the pivot for the tridiagonal solver, given by h_k plus a term that is a fraction (determined from the tridiagonal solver) of Kddt_h for the interface above [H ~> m or kg m-2].

  • dte_term :: [in] A diffusivity-independent term related to the temperature change in the layer below the interface [C H ~> degC m or degC kg m-2].

  • dse_term :: [in] A diffusivity-independent term related to the salinity change in the layer below the interface [S H ~> ppt m or ppt kg m-2].

  • dt_km1_t2 :: [in] A diffusivity-independent term related to the temperature change in the layer above the interface [C ~> degC].

  • ds_km1_t2 :: [in] A diffusivity-independent term related to the salinity change in the layer above the interface [S ~> ppt].

  • pres_z :: [in] The rescaled hydrostatic interface pressure, which relates the changes in column thickness to the energy that is radiated as gravity waves and unavailable to drive mixing [R Z2 T-2 ~> J m-3].

  • dt_to_dpe_k :: [in] A factor (pres_lay*mass_lay*dSpec_vol/dT) relating a layer’s temperature change to the change in column potential energy, including all implicit diffusive changes in the temperatures of all the layers below [R Z3 T-2 C-1 ~> J m-2 degC-1].

  • ds_to_dpe_k :: [in] A factor (pres_lay*mass_lay*dSpec_vol/dS) relating a layer’s salinity change to the change in column potential energy, including all implicit diffusive changes in the in the salinities of all the layers below [R Z3 T-2 S-1 ~> J m-2 ppt-1].

  • dt_to_dpea :: [in] A factor (pres_lay*mass_lay*dSpec_vol/dT) relating a layer’s temperature change to the change in column potential energy, including all implicit diffusive changes in the temperatures of all the layers above [R Z3 T-2 C-1 ~> J m-2 degC-1].

  • ds_to_dpea :: [in] A factor (pres_lay*mass_lay*dSpec_vol/dS) relating a layer’s salinity change to the change in column potential energy, including all implicit diffusive changes in the salinities of all the layers above [R Z3 T-2 S-1 ~> J m-2 ppt-1].

  • dt_to_dcolht_k :: [in] A factor (mass_lay*dSColHtc_vol/dT) relating a layer’s temperature change to the change in column height, including all implicit diffusive changes in the temperatures of all the layers below [Z C-1 ~> m degC-1].

  • ds_to_dcolht_k :: [in] A factor (mass_lay*dSColHtc_vol/dS) relating a layer’s salinity change to the change in column height, including all implicit diffusive changes in the salinities of all the layers below [Z S-1 ~> m ppt-1].

  • dt_to_dcolhta :: [in] A factor (mass_lay*dSColHtc_vol/dT) relating a layer’s temperature change to the change in column height, including all implicit diffusive changes in the temperatures of all the layers above [Z C-1 ~> m degC-1].

  • ds_to_dcolhta :: [in] A factor (mass_lay*dSColHtc_vol/dS) relating a layer’s salinity change to the change in column height, including all implicit diffusive changes in the salinities of all the layers above [Z S-1 ~> m ppt-1].

  • pe_chg :: [out] The change in column potential energy from applying Kddt_h at the present interface [R Z3 T-2 ~> J m-2].

  • dpec_dkd :: [out] The partial derivative of PE_chg with Kddt_h [R Z3 T-2 H-1 ~> J m-3 or J kg-1].

  • dpe_max :: [out] The maximum change in column potential energy that could be realized by applying a huge value of Kddt_h at the present interface [R Z3 T-2 ~> J m-2].

  • dpec_dkd_0 :: [out] The partial derivative of PE_chg with Kddt_h in the limit where Kddt_h = 0 [R Z3 T-2 H-1 ~> J m-3 or J kg-1].

Called from:

epbl_column

subroutine mom_energetic_pbl/find_mstar(CS, US, Buoyancy_Flux, UStar, UStar_Mean, BLD, Abs_Coriolis, MStar, Langmuir_Number, MStar_LT, Convect_Langmuir_Number)

This subroutine finds the Mstar value for ePBL.

Parameters:
  • cs :: [in] Energetic PBL control structure

  • us :: [in] A dimensional unit scaling type

  • ustar :: [in] ustar w/ gustiness [Z T-1 ~> m s-1]

  • ustar_mean :: [in] ustar w/o gustiness [Z T-1 ~> m s-1]

  • abs_coriolis :: [in] absolute value of the Coriolis parameter [T-1 ~> s-1]

  • buoyancy_flux :: [in] Buoyancy flux [Z2 T-3 ~> m2 s-3]

  • bld :: [in] boundary layer depth [Z ~> m]

  • mstar :: [out] Output mstar (Mixing/ustar**3) [nondim]

  • langmuir_number :: [in] Langmuir number [nondim]

  • mstar_lt :: [out] Mstar increase due to Langmuir turbulence [nondim]

  • convect_langmuir_number :: [out] Langmuir number including buoyancy flux [nondim]

Call to:

mstar_from_ekman mstar_from_rh18 mstar_langmuir use_fixed_mstar

Called from:

epbl_column

subroutine mom_energetic_pbl/mstar_langmuir(CS, US, Abs_Coriolis, Buoyancy_Flux, UStar, BLD, Langmuir_Number, Mstar, MStar_LT, Convect_Langmuir_Number)

This subroutine modifies the Mstar value if the Langmuir number is present.

Parameters:
  • cs :: [in] Energetic PBL control structure

  • us :: [in] A dimensional unit scaling type

  • abs_coriolis :: [in] Absolute value of the Coriolis parameter [T-1 ~> s-1]

  • buoyancy_flux :: [in] Buoyancy flux [Z2 T-3 ~> m2 s-3]

  • ustar :: [in] Surface friction velocity with? gustiness [Z T-1 ~> m s-1]

  • bld :: [in] boundary layer depth [Z ~> m]

  • mstar :: [inout] Input/output mstar (Mixing/ustar**3) [nondim]

  • langmuir_number :: [in] Langmuir number [nondim]

  • mstar_lt :: [out] Mstar increase due to Langmuir turbulence [nondim]

  • convect_langmuir_number :: [out] Langmuir number including buoyancy flux [nondim]

Call to:

langmuir_add langmuir_rescale no_langmuir

Called from:

find_mstar

subroutine mom_energetic_pbl/energetic_pbl_get_mld(CS, MLD, G, US, m_to_MLD_units)

Copies the ePBL active mixed layer depth into MLD, in units of [Z ~> m] unless other units are specified.

Parameters:
  • cs :: [in] Energetic PBL control structure

  • g :: [in] Grid structure

  • us :: [in] A dimensional unit scaling type

  • mld :: [out] Depth of ePBL active mixing layer [Z ~> m] or other units

  • m_to_mld_units :: [in] A conversion factor from meters to the desired units for MLD, sometimes [m Z-1 ~> 1]

Called from:

mom_diabatic_driver::diabatic_ale mom_diabatic_driver::diabatic_ale_legacy mom_hor_bnd_diffusion::hor_bnd_diffusion mom_neutral_diffusion::neutral_diffusion_calc_coeffs mom_dynamics_split_rk2::step_mom_dyn_split_rk2

subroutine mom_energetic_pbl/energetic_pbl_init(Time, G, GV, US, param_file, diag, CS)

This subroutine initializes the energetic_PBL module.

Parameters:
  • time :: [in] The current model time

  • g :: [in] The ocean’s grid structure

  • gv :: [in] The ocean’s vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • param_file :: [in] A structure to parse for run-time parameters

  • diag :: [inout] A structure that is used to regulate diagnostic output

  • cs :: [inout] Energetic PBL control structure

Call to:

additive_string constant_string langmuir_add langmuir_rescale mom_error_handler::mom_error mom_error_handler::mom_mesg mstar_from_ekman mstar_from_rh18 no_langmuir none_string om4_string mom_coms::real_to_efp report_avg_its rescaled_string rh18_string root_tke_string mom_string_functions::uppercase use_fixed_mstar wt_from_croot_tke wt_from_rh18

subroutine mom_energetic_pbl/energetic_pbl_end(CS)

Clean up and deallocate memory associated with the energetic_PBL module.

Parameters:

cs :: [inout] Energetic_PBL control structure

Call to:

mom_coms::efp_to_real mom_error_handler::mom_mesg report_avg_its