namespace mom_neutral_diffusion_aux¶
Overview¶
A column-wise toolbox for implementing neutral diffusion. More…
namespace mom_neutral_diffusion_aux { // global functions subroutine, public set_ndiff_aux_params( CS CS, deg deg, max_iter max_iter, drho_tol drho_tol, xtol xtol, ref_pres ref_pres, force_brent force_brent, EOS EOS, debug debug ); subroutine, public mark_unstable_cells( nk nk, dRdT dRdT, dRdS dRdS, T T, S S, stable_cell stable_cell, ns ns ); subroutine, public increment_interface( nk nk, kl kl, ki ki, stable stable, reached_bottom reached_bottom, searching_this_column searching_this_column, searching_other_column searching_other_column ); real function, public calc_drho( T1 T1, S1 S1, dRdT1 dRdT1, dRdS1 dRdS1, T2 T2, S2 S2, dRdT2 dRdT2, dRdS2 dRdS2 ); subroutine, public drho_at_pos( CS CS, T_ref T_ref, S_ref S_ref, alpha_ref alpha_ref, beta_ref beta_ref, P_top P_top, P_bot P_bot, ppoly_T ppoly_T, ppoly_S ppoly_S, x0 x0, delta_rho delta_rho, P_out P_out, T_out T_out, S_out S_out, alpha_avg_out alpha_avg_out, beta_avg_out beta_avg_out, delta_T_out delta_T_out, delta_S_out delta_S_out ); subroutine, public search_other_column( dRhoTop dRhoTop, dRhoBot dRhoBot, Ptop Ptop, Pbot Pbot, lastP lastP, lastK lastK, kl kl, kl_0 kl_0, ki ki, top_connected top_connected, bot_connected bot_connected, out_P out_P, out_K out_K, search_layer search_layer ); real function, public interpolate_for_nondim_position(dRhoNeg dRhoNeg, Pneg Pneg, dRhoPos dRhoPos, Ppos Ppos); real function, public refine_nondim_position( CS CS, T_ref T_ref, S_ref S_ref, alpha_ref alpha_ref, beta_ref beta_ref, P_top P_top, P_bot P_bot, ppoly_T ppoly_T, ppoly_S ppoly_S, drho_top drho_top, drho_bot drho_bot, min_bound min_bound ); subroutine, public check_neutral_positions( CS CS, Ptop_l Ptop_l, Pbot_l Pbot_l, Ptop_r Ptop_r, Pbot_r Pbot_r, PoL PoL, PoR PoR, Tl_coeffs Tl_coeffs, Tr_coeffs Tr_coeffs, Sl_coeffs Sl_coeffs, Sr_coeffs Sr_coeffs ); subroutine, public kahan_sum(sum sum, summand summand, c c); } // namespace mom_neutral_diffusion_aux
Detailed Documentation¶
A column-wise toolbox for implementing neutral diffusion.
Global Functions¶
subroutine, public set_ndiff_aux_params( CS CS, deg deg, max_iter max_iter, drho_tol drho_tol, xtol xtol, ref_pres ref_pres, force_brent force_brent, EOS EOS, debug debug )
Initialize the parameters used to iteratively find the neutral direction.
Parameters:
cs |
Control structure for refine_pos |
deg |
Degree of polynommial used in reconstruction |
max_iter |
Maximum number of iterations |
drho_tol |
Tolerance for function convergence |
xtol |
Tolerance for change in position |
ref_pres |
Reference pressure to use |
force_brent |
Force Brent method for linear, TEOS-10, and WRIGHT |
debug |
If true, print output use to help debug neutral diffusion |
eos |
Equation of state |
subroutine, public mark_unstable_cells( nk nk, dRdT dRdT, dRdS dRdS, T T, S S, stable_cell stable_cell, ns ns )
Given the reconsturcitons of dRdT, dRdS, T, S mark the cells which are stably stratified parts of the water column For an layer to be unstable the top interface must be denser than the bottom or the bottom interface of the layer.
Parameters:
nk |
Number of levels in a column |
drdt |
drho/dT [kg m-3 degC-1] at interfaces |
drds |
drho/dS [kg m-3 ppt-1] at interfaces |
t |
Temperature [degC] at interfaces |
s |
Salinity [ppt] at interfaces |
stable_cell |
True if this cell is unstably stratified |
ns |
Number of neutral surfaces in unmasked part of the column |
subroutine, public increment_interface( nk nk, kl kl, ki ki, stable stable, reached_bottom reached_bottom, searching_this_column searching_this_column, searching_other_column searching_other_column )
Increments the interface which was just connected and also set flags if the bottom is reached.
Parameters:
nk |
Number of vertical levels |
kl |
Current layer (potentially updated) |
ki |
Current interface |
stable |
True if the cell is stably stratified |
reached_bottom |
Updated when kl == nk and ki == 2 |
searching_this_column |
Updated when kl == nk and ki == 2 |
searching_other_column |
Updated when kl == nk and ki == 2 |
real function, public calc_drho( T1 T1, S1 S1, dRdT1 dRdT1, dRdS1 dRdS1, T2 T2, S2 S2, dRdT2 dRdT2, dRdS2 dRdS2 )
Calculates difference in density at two points (rho1-rho2) with known density derivatives, T, and S.
Parameters:
t1 |
Temperature at point 1 |
s1 |
Salinity at point 1 |
drdt1 |
dRhodT at point 1 |
drds1 |
dRhodS at point 1 |
t2 |
Temperature at point 2 |
s2 |
Salinity at point 2 |
drdt2 |
dRhodT at point 2 |
drds2 |
dRhodS at point |
subroutine, public drho_at_pos( CS CS, T_ref T_ref, S_ref S_ref, alpha_ref alpha_ref, beta_ref beta_ref, P_top P_top, P_bot P_bot, ppoly_T ppoly_T, ppoly_S ppoly_S, x0 x0, delta_rho delta_rho, P_out P_out, T_out T_out, S_out S_out, alpha_avg_out alpha_avg_out, beta_avg_out beta_avg_out, delta_T_out delta_T_out, delta_S_out delta_S_out )
Calculate the difference in neutral density between a reference T, S, alpha, and beta at a point on the polynomial reconstructions of T, S.
Parameters:
cs |
Control structure with parameters for this module |
t_ref |
Temperature at reference surface |
s_ref |
Salinity at reference surface |
alpha_ref |
dRho/dT at reference surface |
beta_ref |
dRho/dS at reference surface |
p_top |
Pressure (Pa) at top interface of layer to be searched |
p_bot |
Pressure (Pa) at bottom interface |
ppoly_t |
Coefficients of T reconstruction |
ppoly_s |
Coefficients of S reconstruciton |
x0 |
Nondimensional position to evaluate |
delta_rho |
The density difference from a reference value |
p_out |
Pressure at point x0 |
t_out |
Temperature at point x0 |
s_out |
Salinity at point x0 |
alpha_avg_out |
Average of alpha between reference and x0 |
beta_avg_out |
Average of beta between reference and x0 |
delta_t_out |
Difference in temperature between reference and x0 |
delta_s_out |
Difference in salinity between reference and x0 |
subroutine, public search_other_column( dRhoTop dRhoTop, dRhoBot dRhoBot, Ptop Ptop, Pbot Pbot, lastP lastP, lastK lastK, kl kl, kl_0 kl_0, ki ki, top_connected top_connected, bot_connected bot_connected, out_P out_P, out_K out_K, search_layer search_layer )
Searches the “other” (searched) column for the position of the neutral surface.
Parameters:
drhotop |
Density difference across top interface |
drhobot |
Density difference across top interface |
ptop |
Pressure at top interface |
pbot |
Pressure at bottom interface |
lastp |
Last position connected in the searched column |
lastk |
Last layer connected in the searched column |
kl |
Layer in the searched column |
kl_0 |
Layer in the searched column |
ki |
Interface of the searched column |
top_connected |
True if the top interface was pointed to |
bot_connected |
True if the top interface was pointed to |
out_p |
Position within searched column |
out_k |
Layer within searched column |
search_layer |
Neutral surface within cell |
real function, public interpolate_for_nondim_position( dRhoNeg dRhoNeg, Pneg Pneg, dRhoPos dRhoPos, Ppos Ppos )
Returns the non-dimensional position between Pneg and Ppos where the interpolated density difference equals zero. The result is always bounded to be between 0 and 1.
Parameters:
drhoneg |
Negative density difference |
pneg |
Position of negative density difference |
drhopos |
Positive density difference |
ppos |
Position of positive density difference |
real function, public refine_nondim_position( CS CS, T_ref T_ref, S_ref S_ref, alpha_ref alpha_ref, beta_ref beta_ref, P_top P_top, P_bot P_bot, ppoly_T ppoly_T, ppoly_S ppoly_S, drho_top drho_top, drho_bot drho_bot, min_bound min_bound )
Use root-finding methods to find where dRho = 0, based on the equation of state and the polynomial reconstructions of temperature, salinity. Initial guess is based on the zero crossing of based on linear profiles of dRho, T, and S, between the top and bottom interface. If second derivatives of the EOS are available, it starts with a Newton’s method. However, Newton’s method is not guaranteed to be bracketed, a check is performed to see if it it diverges outside the interval. In that case (or in the case that second derivatives are not available), Brent’s method is used following the implementation found at https://people.sc.fsu.edu/~jburkardt/f_src/brent/brent.f90.
Parameters:
cs |
Control structure with parameters for this module |
t_ref |
Temperature of the neutral surface at the searched from interface |
s_ref |
Salinity of the neutral surface at the searched from interface |
alpha_ref |
dRho/dT of the neutral surface at the searched from interface |
beta_ref |
dRho/dS of the neutral surface at the searched from interface |
p_top |
Pressure at the top interface in the layer to be searched |
p_bot |
Pressure at the bottom interface in the layer to be searched |
ppoly_t |
Coefficients of the order N polynomial reconstruction of T within the layer to be searched. |
ppoly_s |
Coefficients of the order N polynomial reconstruction of T within the layer to be searched. |
drho_top |
Delta rho at top interface (or previous position in cell |
drho_bot |
Delta rho at bottom interface |
min_bound |
Lower bound of position, also serves as the initial guess |
subroutine, public check_neutral_positions( CS CS, Ptop_l Ptop_l, Pbot_l Pbot_l, Ptop_r Ptop_r, Pbot_r Pbot_r, PoL PoL, PoR PoR, Tl_coeffs Tl_coeffs, Tr_coeffs Tr_coeffs, Sl_coeffs Sl_coeffs, Sr_coeffs Sr_coeffs )
Parameters:
cs |
Control structure with parameters for this module |
ptop_l |
Pressure at top interface of left cell |
pbot_l |
Pressure at bottom interface of left cell |
ptop_r |
Pressure at top interface of right cell |
pbot_r |
Pressure at bottom interface of right cell |
pol |
Nondim position in left cell |
por |
Nondim position in right cell |
tl_coeffs |
T polynomial coefficients of left cell |
tr_coeffs |
T polynomial coefficients of right cell |
sl_coeffs |
S polynomial coefficients of left cell |
sr_coeffs |
S polynomial coefficients of right cell |
subroutine, public kahan_sum(sum sum, summand summand, c c)
Do a compensated sum to account for roundoff level.
Parameters:
sum |
Running sum |
summand |
Term to be added |
c |
Keep track of roundoff |