namespace mom_wave_speed¶
Overview¶
Routines for calculating baroclinic wave speeds. More…
namespace mom_wave_speed { // global functions subroutine, public wave_speed( h h, tv tv, G G, GV GV, US US, cg1 cg1, CS CS, full_halos full_halos, use_ebt_mode use_ebt_mode, mono_N2_column_fraction mono_N2_column_fraction, mono_N2_depth mono_N2_depth, modal_structure modal_structure ); subroutine, public wave_speeds( h h, tv tv, G G, GV GV, US US, nmodes nmodes, cn cn, CS CS, full_halos full_halos ); subroutine, public wave_speed_init( CS CS, use_ebt_mode use_ebt_mode, mono_N2_column_fraction mono_N2_column_fraction, mono_N2_depth mono_N2_depth ); subroutine, public wave_speed_set_param( CS CS, use_ebt_mode use_ebt_mode, mono_N2_column_fraction mono_N2_column_fraction, mono_N2_depth mono_N2_depth ); } // namespace mom_wave_speed
Detailed Documentation¶
Routines for calculating baroclinic wave speeds.
Subroutine wave_speed() solves for the first baroclinic mode wave speed. (It could solve for all the wave speeds, but the iterative approach taken here means that this is not particularly efficient.)
If e(k)
is the perturbation interface height, this means solving for the smallest eigenvalue (lam
= 1/c^2) of the system
-Igu(k)*e(k-1) + (Igu(k)+Igl(k)-lam)*e(k) - Igl(k)*e(k+1) = 0.0
with rigid lid boundary conditions e(1) = e(nz+1) = 0.0 giving
- (Igu(2)+Igl(2)-lam)*e(2) - Igl(2)*e(3) = 0.0
-Igu(nz)*e(nz-1) + (Igu(nz)+Igl(nz)-lam)*e(nz) = 0.0
Here Igl(k) = 1.0/(gprime(k)*h(k)) ; Igu(k) = 1.0/(gprime(k)*h(k-1))
Alternately, these same eigenvalues can be found from the second smallest eigenvalue of the Montgomery potential (M(k)) calculation:
-Igl(k)*M(k-1) + (Igl(k)+Igu(k+1)-lam)*M(k) - Igu(k+1)*M(k+1) = 0.0
with rigid lid and flat bottom boundary conditions
- (Igu(2)-lam)*M(1) - Igu(2)*M(2) = 0.0
-Igl(nz)*M(nz-1) + (Igl(nz)-lam)*M(nz) = 0.0
Note that the barotropic mode has been eliminated from the rigid lid interface height equations, hence the matrix is one row smaller. Without the rigid lid, the top boundary condition is simpler to implement with the M equations.
Global Functions¶
subroutine, public wave_speed( h h, tv tv, G G, GV GV, US US, cg1 cg1, CS CS, full_halos full_halos, use_ebt_mode use_ebt_mode, mono_N2_column_fraction mono_N2_column_fraction, mono_N2_depth mono_N2_depth, modal_structure modal_structure )
Calculates the wave speed of the first baroclinic mode.
Parameters:
g |
Ocean grid structure |
gv |
Vertical grid structure |
us |
A dimensional unit scaling type |
h |
Layer thickness [H ~> m or kg m-2] |
tv |
Thermodynamic variables |
cg1 |
First mode internal wave speed [m s-1] |
cs |
Control structure for MOM_wave_speed |
full_halos |
If true, do the calculation over the entire computational domain. |
use_ebt_mode |
If true, use the equivalent barotropic mode instead of the first baroclinic mode. |
mono_n2_column_fraction |
The lower fraction of water column over which N2 is limited as monotonic for the purposes of calculating vertical modal structure. |
mono_n2_depth |
A depth below which N2 is limited as monotonic for the purposes of calculating vertical modal structure [m]. |
modal_structure |
Normalized model structure [nondim] |
subroutine, public wave_speeds( h h, tv tv, G G, GV GV, US US, nmodes nmodes, cn cn, CS CS, full_halos full_halos )
Calculates the wave speeds for the first few barolinic modes.
Parameters:
g |
Ocean grid structure |
gv |
Vertical grid structure |
us |
A dimensional unit scaling type |
h |
Layer thickness [H ~> m or kg m-2] |
tv |
Thermodynamic variables |
nmodes |
Number of modes |
cn |
Waves speeds [L T-1 ~> m s-1] |
cs |
Control structure for MOM_wave_speed |
full_halos |
If true, do the calculation over the entire computational domain. |
subroutine, public wave_speed_init( CS CS, use_ebt_mode use_ebt_mode, mono_N2_column_fraction mono_N2_column_fraction, mono_N2_depth mono_N2_depth )
Initialize control structure for MOM_wave_speed.
Parameters:
cs |
Control structure for MOM_wave_speed |
use_ebt_mode |
If true, use the equivalent barotropic mode instead of the first baroclinic mode. |
mono_n2_column_fraction |
The lower fraction of water column over which N2 is limited as monotonic for the purposes of calculating the vertical modal structure. |
mono_n2_depth |
The depth below which N2 is limited as monotonic for the purposes of calculating the vertical modal structure. |
subroutine, public wave_speed_set_param( CS CS, use_ebt_mode use_ebt_mode, mono_N2_column_fraction mono_N2_column_fraction, mono_N2_depth mono_N2_depth )
Sets internal parameters for MOM_wave_speed.
Parameters:
cs |
Control structure for MOM_wave_speed |
use_ebt_mode |
If true, use the equivalent barotropic mode instead of the first baroclinic mode. |
mono_n2_column_fraction |
The lower fraction of water column over which N2 is limited as monotonic for the purposes of calculating the vertical modal structure. |
mono_n2_depth |
The depth below which N2 is limited as monotonic for the purposes of calculating the vertical modal structure. |