mom_internal_tides module reference

Subroutines that use the ray-tracing equations to propagate the internal tide energy density.

More…

Data Types

int_tide_cs

This control structure has parameters for the MOM_internal_tides module.

loop_bounds_type

A structure with the active energy loop bounds.

Functions/Subroutines

propagate_int_tide()

Calls subroutines in this file that are needed to refract, propagate, and dissipate energy density of the internal tide.

sum_en()

Checks for energy conservation on computational domain.

itidal_lowmode_loss()

Calculates the energy lost from the propagating internal tide due to scattering over small-scale roughness along the lines of Jayne & St.

get_lowmode_loss()

This subroutine extracts the energy lost from the propagating internal which has been summed across all angles, frequencies, and modes for a given mechanism and location.

get_lowmode_diffusivity()

Returns the values of diffusivity corresponding to various mechanisms.

refract()

Implements refraction on the internal waves at a single frequency.

ppm_angular_advect()

This subroutine calculates the 1-d flux for advection in angular space using a monotonic piecewise parabolic scheme.

propagate()

Propagates internal waves at a single frequency.

propagate_x()

Propagates the internal wave energy in the logical x-direction.

propagate_y()

Propagates the internal wave energy in the logical y-direction.

zonal_flux_en()

Evaluates the zonal mass or volume fluxes in a layer.

merid_flux_en()

Evaluates the meridional mass or volume fluxes in a layer.

reflect()

Reflection of the internal waves at a single frequency.

turning_latitude()

teleport()

Moves energy across lines of partial reflection to prevent reflection of energy that is supposed to get across.

correct_halo_rotation()

Rotates points in the halos where required to accommodate changes in grid orientation, such as at the tripolar fold.

correct_halo_rotation_2d()

Rotates points in the halos where required to accommodate changes in grid orientation, such as at the tripolar fold.

ppm_reconstruction_x()

Calculates left/right edge values for PPM reconstruction in x-direction.

ppm_reconstruction_y()

Calculates left/right edge valus for PPM reconstruction in y-direction.

ppm_limit_pos()

Limits the left/right edge values of the PPM reconstruction to give a reconstruction that is positive-definite.

minmod_limiter()

Limits the left/right edge values using the simple minmod limiter written in a way that avoids branching in favor of intrinsics.

register_int_tide_restarts()

internal_tides_init()

This subroutine initializes the internal tides module.

internal_tides_end()

This subroutine deallocates the memory associated with the internal tides control structure.

Detailed Description

<undocumented>

Type Documentation

type  mom_internal_tides/int_tide_cs

This control structure has parameters for the MOM_internal_tides module.

Type fields:
  • % id_cg1 :: integer Diag handles.

  • % id_cn :: integer, dimension(:), allocatable Diag handles.

  • % id_tot_en :: integer Diag handles.

  • % id_refl_pref :: integer Diag handles.

  • % id_refl_ang :: integer Diag handles.

  • % id_land_mask :: integer Diag handles.

  • % id_trans :: integer Diag handles.

  • % id_residual :: integer Diag handles.

  • % id_dx_cv :: integer Diag handles.

  • % id_dy_cu :: integer Diag handles.

  • % id_tot_leak_loss :: integer Diag handles.

  • % id_tot_quad_loss :: integer Diag handles.

  • % id_tot_itidal_loss :: integer Diag handles.

  • % id_tot_froude_loss :: integer Diag handles.

  • % id_tot_residual_loss :: integer Diag handles.

  • % id_tot_allprocesses_loss :: integer Diag handles.

  • % id_en_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_itidal_loss_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_leak_loss_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_quad_loss_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_froude_loss_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_residual_loss_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_allprocesses_loss_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_itide_drag :: integer, dimension(:,:), allocatable Diag handles.

  • % id_ub_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_cp_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_en_ang_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_itidal_loss_ang_mode :: integer, dimension(:,:), allocatable Diag handles.

  • % id_tke_itidal_input :: integer, dimension(:), allocatable Diag handles.

  • % id_ustruct_mode :: integer, dimension(:), allocatable Diag handles.

  • % id_wstruct_mode :: integer, dimension(:), allocatable Diag handles.

  • % id_int_w2_mode :: integer, dimension(:), allocatable Diag handles.

  • % id_int_u2_mode :: integer, dimension(:), allocatable Diag handles.

  • % id_int_n2w2_mode :: integer, dimension(:), allocatable Diag handles.

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

  • % do_int_tides :: logical If true, use the internal tide code.

  • % nfreq :: integer The number of internal tide frequency bands.

  • % nmode :: integer The number of internal tide vertical modes.

  • % nangle :: integer The number of internal tide angular orientations.

  • % energized_angle :: integer If positive, only this angular band is energized for debugging purposes.

  • % dt_itides :: real The timestep for internal tides ray-tracing [T ~> s].

  • % uniform_test_cg :: real Uniform group velocity of internal tide for testing internal tides [L T-1 ~> m s-1].

  • % corner_adv :: logical If true, use a corner advection rather than PPM.

  • % upwind_1st :: logical If true, use a first-order upwind scheme.

  • % simple_2nd :: logical If true, use a simple second order (arithmetic mean) interpolation of the edge values instead of the higher order interpolation.

  • % vol_cfl :: logical If true, use the ratio of the open face lengths to the tracer cell areas when estimating CFL numbers. Without aggress_adjust, the default is false; it is always true with aggress_adjust.

  • % use_ppmang :: logical If true, use PPM for advection of energy in angular space.

  • % update_kd :: logical If true, the scheme will modify the diffusivities seen by the dynamics.

  • % apply_refraction :: logical If false, skip refraction (for debugging)

  • % apply_propagation :: logical If False, do not propagate energy (for debugging)

  • % turn_critical_lat :: logical If True, rays change direction at critical latitude instead of being trapped.

  • % reflect_critical_lat :: logical If True, rays reflect at the critical latitude instead of turning parallel to it.

  • % debug :: logical If true, use debugging prints.

  • % init_forcing_only :: logical if True, add TKE forcing only at first step (for debugging)

  • % force_posit_en :: logical if True, remove subroundoff negative values (needs enhancement)

  • % add_tke_forcing :: logical Whether to add forcing, used by init_forcing_only.

  • % fraction_tidal_input :: real, dimension(:,:), allocatable how the energy from one tidal component is distributed over the various vertical modes, 2d in frequency and mode [nondim]

  • % refl_angle :: real, dimension(:,:), allocatable local coastline/ridge/shelf angles read from file [rad]

  • % nullangle :: real placeholder value in cells with no reflection [rad]

  • % refl_pref :: real, dimension(:,:), allocatable partial reflection coeff for each “coast cell” [nondim]

  • % refl_pref_logical :: logical, dimension(:,:), allocatable true if reflecting cell with partial reflection

  • % refl_dbl :: logical, dimension(:,:), allocatable identifies reflection cells where double reflection is possible (i.e. ridge cells)

  • % trans :: real, dimension(:,:), allocatable partial transmission coeff for each “coast cell” [nondim]

  • % residual :: real, dimension(:,:), allocatable residual of reflection and transmission coeff for each “coast cell” [nondim]

  • % cp :: real, dimension(:,:,:,:), allocatable horizontal phase speed [L T-1 ~> m s-1]

  • % tke_leak_loss :: real, dimension(:,:,:,:,:), allocatable energy lost due to misc background processes [H Z2 T-3 ~> m3 s-3 or W m-2]

  • % tke_quad_loss :: real, dimension(:,:,:,:,:), allocatable energy lost due to quadratic bottom drag [H Z2 T-3 ~> m3 s-3 or W m-2]

  • % tke_froude_loss :: real, dimension(:,:,:,:,:), allocatable energy lost due to wave breaking [H Z2 T-3 ~> m3 s-3 or W m-2]

  • % tke_itidal_loss_fixed :: real, dimension(:,:), allocatable Fixed part of the energy lost due to small-scale drag [H Z2 L-2 ~> kg m-2] here. This will be multiplied by N and the squared near-bottom velocity (and by the near-bottom density in non-Boussinesq mode) to get the energy losses in [R Z4 H-1 L-2 ~> kg m-2 or m].

  • % tke_itidal_loss :: real, dimension(:,:,:,:,:), allocatable energy lost due to small-scale wave drag [H Z2 T-3 ~> m3 s-3 or W m-2]

  • % tke_residual_loss :: real, dimension(:,:,:,:,:), allocatable internal tide energy loss due to the residual at slopes [H Z2 T-3 ~> m3 s-3 or W m-2]

  • % tke_slope_loss :: real, dimension(:,:,:,:,:), allocatable internal tide energy loss due to the residual at slopes [H Z2 T-3 ~> m3 s-3 or W m-2]

  • % tke_input_glo_dt :: real, dimension(:,:), allocatable The integrated energy input to the internal waves [H Z2 L2 T-2 ~> m5 s-2 or J].

  • % tke_leak_loss_glo_dt :: real, dimension(:,:), allocatable Integrated energy lost due to misc background processes [H Z2 L2 T-2 ~> m5 s-2 or J].

  • % tke_quad_loss_glo_dt :: real, dimension(:,:), allocatable Integrated energy lost due to quadratic bottom drag [H Z2 L2 T-2 ~> m5 s-2 or J].

  • % tke_froude_loss_glo_dt :: real, dimension(:,:), allocatable Integrated energy lost due to wave breaking [H Z2 L2 T-2 ~> m5 s-2 or J].

  • % tke_itidal_loss_glo_dt :: real, dimension(:,:), allocatable energy lost due to small-scale wave drag [H Z2 T-2 ~> m3 s-2 or J m-2]

  • % tke_residual_loss_glo_dt :: real, dimension(:,:), allocatable internal tide energy loss due to the residual at slopes [H Z2 L2 T-2 ~> m5 s-2 or J]

  • % error_mode :: real, dimension(:,:), allocatable internal tide energy budget error for each mode [H Z2 L2 T-2 ~> m5 s-2 or J]

  • % tot_leak_loss :: real, dimension(:,:), allocatable Energy loss rates due to misc background processes, summed over angle, frequency and mode [H Z2 T-3 ~> m3 s-3 or W m-2].

  • % tot_quad_loss :: real, dimension(:,:), allocatable Energy loss rates due to quadratic bottom drag, summed over angle, frequency and mode [H Z2 T-3 ~> m3 s-3 or W m-2].

  • % tot_itidal_loss :: real, dimension(:,:), allocatable Energy loss rates due to small-scale drag, summed over angle, frequency and mode [H Z2 T-3 ~> m3 s-3 or W m-2].

  • % tot_froude_loss :: real, dimension(:,:), allocatable Energy loss rates due to wave breaking, summed over angle, frequency and mode [H Z2 T-3 ~> m3 s-3 or W m-2].

  • % tot_residual_loss :: real, dimension(:,:), allocatable Energy loss rates due to residual on slopes, summed over angle, frequency and mode [H Z2 T-3 ~> m3 s-3 or W m-2].

  • % tot_allprocesses_loss :: real, dimension(:,:), allocatable Energy loss rates due to all processes, summed over angle, frequency and mode [H Z2 T-3 ~> m3 s-3 or W m-2].

  • % w_struct :: real, dimension(:,:,:,:), allocatable Vertical structure of vertical velocity (normalized) for each frequency and each mode [nondim].

  • % u_struct :: real, dimension(:,:,:,:), allocatable Vertical structure of horizontal velocity (normalized and divided by layer thicknesses) for each frequency and each mode [Z-1 ~> m-1].

  • % u_struct_max :: real, dimension(:,:,:), allocatable Maximum of u_struct, for each mode [Z-1 ~> m-1].

  • % u_struct_bot :: real, dimension(:,:,:), allocatable Bottom value of u_struct, for each mode [Z-1 ~> m-1].

  • % int_w2 :: real, dimension(:,:,:), allocatable Vertical integral of w_struct squared, for each mode [H ~> m or kg m-2].

  • % int_u2 :: real, dimension(:,:,:), allocatable Vertical integral of u_struct squared, for each mode [H Z-2 ~> m-1 or kg m-4].

  • % int_n2w2 :: real, dimension(:,:,:), allocatable Depth-integrated Brunt Vaissalla freqency times vertical profile squared, for each mode [H T-2 ~> m s-2 or kg m-2 s-2].

  • % q_itides :: real fraction of local dissipation [nondim]

  • % mixing_effic :: real mixing efficiency [nondim]

  • % en_sum :: real global sum of energy for use in debugging, in MKS units [m5 s-2 or J]

  • % en_underflow :: real A minuscule amount of energy [H Z2 T-2 ~> m3 s-2 or J m-2].

  • % en_restart_power :: integer A power factor of 2 by which to multiply the energy in restart [nondim].

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

  • % pass_en :: type(group_pass_type) Pass 5d array Energy as a group of 3d arrays.

  • % inputdir :: character(len=200) directory to look for coastline angle file

  • % itides_adv_limiter :: integer The type of limiter to use for the energy advection scheme.

  • % decay_rate_2d :: real, dimension(:,:,:,:), allocatable rate at which internal tide energy is lost to the interior ocean internal wave field as a function of longitude, latitude, frequency and vertical mode [T-1 ~> s-1].

  • % cdrag :: real The bottom drag coefficient [nondim].

  • % drag_min_depth :: real The minimum total ocean thickness that will be used in the denominator of the quadratic drag terms for internal tides when INTERNAL_TIDE_QUAD_DRAG is true [H ~> m or kg m-2].

  • % gamma_osborn :: real Mixing efficiency from Osborn 1980 [nondim].

  • % kd_min :: real The minimum diapycnal diffusivity. [L2 T-1 ~> m2 s-1].

  • % max_tke_to_kd :: real Maximum allowed value for TKE_to_kd [H Z2 T-3 ~> m3 s-3 or W m-2].

  • % min_thick_layer_kd :: real minimum layer thickness allowed to use with TKE_to_kd [H ~> m or kg m-2]

  • % apply_background_drag :: logical If true, apply a drag due to background processes as a sink.

  • % apply_bottom_drag :: logical If true, apply a quadratic bottom drag as a sink.

  • % apply_wave_drag :: logical If true, apply scattering due to small-scale roughness as a sink.

  • % apply_froude_drag :: logical If true, apply wave breaking as a sink.

  • % en_check_tol :: real An energy density tolerance for flagging points with small negative internal tide energy [H Z2 T-2 ~> m3 s-2 or J m-2].

  • % apply_residual_drag :: logical If true, apply sink from residual term of reflection/transmission.

  • % use_2d_decay_rate :: logical If true, use a spatially varying decay rate for each harmonic.

  • % en :: real, dimension(:,:,:,:,:), allocatable The internal wave energy density as a function of (i,j,angle,frequency,mode) integrated within an angular and frequency band [H Z2 T-2 ~> m3 s-2 or J m-2].

  • % en_ini_glo :: real, dimension(:,:), allocatable The internal wave energy density as a function of (frequency,mode) spatially integrated within an angular and frequency band [H Z2 L2 T-2 ~> m5 s-2 or J] only at the start of the routine (for diags)

  • % en_end_glo :: real, dimension(:,:), allocatable The internal wave energy density as a function of (frequency,mode) spatially integrated within an angular and frequency band [H Z2 L2 T-2 ~> m5 s-2 or J] only at the end of the routine (for diags)

  • % en_restart_mode1 :: real, dimension(:,:,:,:), allocatable The internal wave energy density as a function of (i,j,angle,freq) for mode 1 [H Z2 T-2 ~> m3 s-2 or J m-2].

  • % en_restart_mode2 :: real, dimension(:,:,:,:), allocatable The internal wave energy density as a function of (i,j,angle,freq) for mode 2 [H Z2 T-2 ~> m3 s-2 or J m-2].

  • % en_restart_mode3 :: real, dimension(:,:,:,:), allocatable The internal wave energy density as a function of (i,j,angle,freq) for mode 3 [H Z2 T-2 ~> m3 s-2 or J m-2].

  • % en_restart_mode4 :: real, dimension(:,:,:,:), allocatable The internal wave energy density as a function of (i,j,angle,freq) for mode 4 [H Z2 T-2 ~> m3 s-2 or J m-2].

  • % en_restart_mode5 :: real, dimension(:,:,:,:), allocatable The internal wave energy density as a function of (i,j,angle,freq) for mode 5 [H Z2 T-2 ~> m3 s-2 or J m-2].

  • % frequency :: real, dimension(:), allocatable The frequency of each band [T-1 ~> s-1].

  • % int_tide_decay_scale :: real vertical decay scale for St Laurent profile [Z ~> m]

  • % int_tide_decay_scale_slope :: real vertical decay scale for St Laurent profile on slopes [Z ~> m]

  • % wave_speed :: type(wave_speed_cs) Wave speed control structure.

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

[source]

type  mom_internal_tides/loop_bounds_type

A structure with the active energy loop bounds.

Type fields:
  • % ish :: integer, private The active loop bounds.

  • % ieh :: integer, private The active loop bounds.

  • % jsh :: integer, private The active loop bounds.

  • % jeh :: integer, private The active loop bounds.

[source]

Function/Subroutine Documentation

subroutine mom_internal_tides/propagate_int_tide(h, tv, Nb, Rho_bot, dt, G, GV, US, inttide_input_CSp, CS)

Calls subroutines in this file that are needed to refract, propagate, and dissipate energy density of the internal tide.

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

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

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

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

  • tv :: tv [in] Pointer to thermodynamic variables (needed for wave structure).

  • nb :: [inout] Near-bottom buoyancy frequency [T-1 ~> s-1]. In some cases the input values are used, but in others this is set along with the wave speeds.

  • rho_bot :: [in] Near-bottom density or the Boussinesq reference density [R ~> kg m-3].

  • dt :: dt [in] Length of time over which to advance the internal tides [T ~> s].

  • inttide_input_csp :: [in] Internal tide input control structure

  • cs :: [inout] Internal tide control structure

Call to:

correct_halo_rotation mom_int_tide_input::get_barotropic_tidal_vel mom_int_tide_input::get_input_tke mom_spatial_means::global_area_integral itidal_lowmode_loss mom_error_handler::mom_error propagate refract mom_io::stdout sum_en mom_wave_speed::wave_speeds

[source]

subroutine mom_internal_tides/sum_en(G, GV, US, CS, En, label)

Checks for energy conservation on computational domain.

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

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

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

  • cs :: [inout] Internal tide control structure

  • en :: [in] The energy density of the internal tides [H Z2 T-2 ~> m3 s-2 or J m-2].

  • label :: label [in] A label to use in error messages

Call to:

mom_spatial_means::global_area_integral

Called from:

propagate propagate_int_tide

[source]

subroutine mom_internal_tides/itidal_lowmode_loss(G, GV, US, CS, Nb, Rho_bot, Ub, En, TKE_loss_fixed, TKE_loss, dt, halo_size)

Calculates the energy lost from the propagating internal tide due to scattering over small-scale roughness along the lines of Jayne & St. Laurent (2001).

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

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

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

  • cs :: [in] Internal tide control structure

  • nb :: [in] Near-bottom stratification [T-1 ~> s-1].

  • rho_bot :: [in] Near-bottom density [R ~> kg m-3].

  • ub :: [inout] RMS (over one period) near-bottom horizontal

  • tke_loss_fixed :: [in] Fixed part of energy loss [R Z4 H-1 L-2 ~> kg m-2 or m]

  • en :: [inout] Energy density of the internal waves [H Z2 T-2 ~> m3 s-2 or J m-2].

  • tke_loss :: [out] Energy loss rate [H Z2 T-3 ~> m3 s-3 or W m-2]

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

  • halo_size :: halo_size [in] The halo size over which to do the calculations

Called from:

propagate_int_tide

[source]

subroutine mom_internal_tides/get_lowmode_loss(i, j, G, CS, mechanism, TKE_loss_sum)

This subroutine extracts the energy lost from the propagating internal which has been summed across all angles, frequencies, and modes for a given mechanism and location.

It can be called from another module to get values from this module’s (private) CS.

Parameters:
  • i :: i [in] The i-index of the value to be reported.

  • j :: j [in] The j-index of the value to be reported.

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

  • cs :: [in] Internal tide control structure

  • mechanism :: mechanism [in] The named mechanism of loss to return

  • tke_loss_sum :: [out] Total energy loss rate due to specified mechanism [H Z2 T-3 ~> m3 s-3 or W m-2].

Called from:

mom_tidal_mixing::add_int_tide_diffusivity get_lowmode_diffusivity

[source]

subroutine mom_internal_tides/get_lowmode_diffusivity(G, GV, h, tv, US, h_bot, k_bot, j, N2_lay, N2_int, TKE_to_Kd, Kd_max, CS, Kd_leak, Kd_quad, Kd_itidal, Kd_Froude, Kd_slope, Kd_lay, Kd_int, profile_leak, profile_quad, profile_itidal, profile_Froude, profile_slope)

Returns the values of diffusivity corresponding to various mechanisms.

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

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

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

  • tv :: tv [in] Structure containing pointers to any available

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

  • h_bot :: h_bot [in] Bottom boundary layer thickness [H ~> m or kg m-2]

  • k_bot :: k_bot [in] Bottom boundary layer top layer index

  • j :: j [in] The j-index to work on

  • n2_lay :: [in] The squared buoyancy frequency of the layers [T-2 ~> s-2].

  • n2_int :: [in] The squared buoyancy frequency of the interfaces [T-2 ~> s-2].

  • tke_to_kd :: [in] The conversion rate between the TKE dissipated within a layer and the diapycnal diffusivity within that layer, usually (~Rho_0 / (G_Earth * dRho_lay)) [T2 Z-1 ~> s2 m-1]

  • kd_max :: [in] The maximum increment for diapycnal diffusivity due to TKE-based processes [H Z T-1 ~> m2 s-1 or kg m-1 s-1]. Set this to a negative value to have no limit. [H Z T-1 ~> m2 s-1 or kg m-1 s-1].

  • cs :: [in] The control structure for this module

  • kd_leak :: [out] Diffusivity due to background drag [H Z T-1 ~> m2 s-1 or kg m-1 s-1].

  • kd_quad :: [out] Diffusivity due to bottom drag [H Z T-1 ~> m2 s-1 or kg m-1 s-1].

  • kd_itidal :: [out] Diffusivity due to wave drag [H Z T-1 ~> m2 s-1 or kg m-1 s-1].

  • kd_froude :: [out] Diffusivity due to high Froude breaking [H Z T-1 ~> m2 s-1 or kg m-1 s-1].

  • kd_slope :: [out] Diffusivity due to critical slopes [H Z T-1 ~> m2 s-1 or kg m-1 s-1].

  • kd_lay :: [inout] The diapycnal diffusivity in layers [H Z T-1 ~> m2 s-1 or kg m-1 s-1].

  • kd_int :: [inout] The diapycnal diffusivity at interfaces [H Z T-1 ~> m2 s-1 or kg m-1 s-1].

  • profile_leak :: profile_leak [out] Normalized profile for background drag [H-1 ~> m-1 or m2 kg-1]

  • profile_quad :: profile_quad [out] Normalized profile for bottom drag [H-1 ~> m-1 or m2 kg-1]

  • profile_itidal :: profile_itidal [out] Normalized profile for wave drag [H-1 ~> m-1 or m2 kg-1]

  • profile_froude :: [out] Normalized profile for Froude drag [H-1 ~> m-1 or m2 kg-1]

  • profile_slope :: profile_slope [out] Normalized profile for critical slopes [H-1 ~> m-1 or m2 kg-1]

Call to:

get_lowmode_loss mom_error_handler::mom_error mom_io::stdout

Called from:

mom_set_diffusivity::set_diffusivity

[source]

subroutine mom_internal_tides/refract(En, cn, freq, dt, G, US, NAngle, use_PPMang)

Implements refraction on the internal waves at a single frequency.

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

  • nangle :: [in] The number of wave orientations in the discretized wave energy spectrum.

  • en :: [inout] The internal gravity wave energy density as a

  • cn :: cn [in] Baroclinic mode speed [L T-1 ~> m s-1].

  • freq :: freq [in] Wave frequency [T-1 ~> s-1].

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

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

  • use_ppmang :: [in] If true, use PPM for advection rather than upwind.

Call to:

mom_error_handler::mom_error ppm_angular_advect

Called from:

propagate_int_tide

[source]

subroutine mom_internal_tides/ppm_angular_advect(En2d, CFL_ang, Flux_En, NAngle, dt, halo_ang)

This subroutine calculates the 1-d flux for advection in angular space using a monotonic piecewise parabolic scheme. This needs to be called from within i and j spatial loops.

Parameters:
  • nangle :: [in] The number of wave orientations in the discretized wave energy spectrum [nondim]

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

  • halo_ang :: halo_ang [in] The halo size in angular space

  • en2d :: [in] The internal gravity wave energy density as a

  • cfl_ang :: [in] The CFL number of the energy advection across angles [nondim]

  • flux_en :: [out] The time integrated internal wave energy flux across angles [H Z2 T-2 ~> m3 s-2 or J m-2].

Called from:

refract

[source]

subroutine mom_internal_tides/propagate(En, cn, freq, dt, G, GV, US, CS, NAngle, test, halo_size, residual_loss)

Propagates internal waves at a single frequency.

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

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

  • nangle :: [in] The number of wave orientations in the discretized wave energy spectrum.

  • en :: [inout] The internal gravity wave energy density as a

  • cn :: cn [in] Baroclinic mode speed [L T-1 ~> m s-1].

  • freq :: freq [in] Wave frequency [T-1 ~> s-1].

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

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

  • test :: test [in] test rotation vector

  • cs :: [inout] Internal tide control structure

  • halo_size :: halo_size [in] halo size for correct rotation

  • residual_loss :: residual_loss [inout] internal tide energy loss due

Call to:

correct_halo_rotation_2d propagate_x propagate_y mom_io::stdout sum_en

Called from:

propagate_int_tide

[source]

subroutine mom_internal_tides/propagate_x(En, speed_x, Cgx_av, dCgx, dt, G, US, Nangle, CS, LB, residual_loss, freq2)

Propagates the internal wave energy in the logical x-direction.

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

  • nangle :: [in] The number of wave orientations in the discretized wave energy spectrum.

  • en :: [inout] The energy density integrated over an angular

  • speed_x :: speed_x [in] The magnitude of the group velocity at the

  • cgx_av :: [in] The average x-projection in each angular band [nondim]

  • dcgx :: [in] The difference in x-projections between the edges of each angular band [nondim].

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

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

  • cs :: [in] Internal tide control structure

  • lb :: [in] A structure with the active energy loop bounds.

  • residual_loss :: residual_loss [inout] internal tide energy loss due

  • freq2 :: freq2 [in] The square of internal tides frequency [T-2 ~> s-2].

Call to:

ppm_reconstruction_x reflect turning_latitude zonal_flux_en

Called from:

propagate

[source]

subroutine mom_internal_tides/propagate_y(En, speed_y, Cgy_av, dCgy, dt, G, US, Nangle, CS, LB, residual_loss, freq2)

Propagates the internal wave energy in the logical y-direction.

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

  • nangle :: [in] The number of wave orientations in the discretized wave energy spectrum.

  • en :: [inout] The energy density integrated over an angular

  • speed_y :: speed_y [in] The magnitude of the group velocity at the

  • cgy_av :: [in] The average y-projection in each angular band [nondim]

  • dcgy :: [in] The difference in y-projections between the edges of each angular band [nondim]

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

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

  • cs :: [in] Internal tide control structure

  • lb :: [in] A structure with the active energy loop bounds.

  • residual_loss :: residual_loss [inout] internal tide energy loss due

  • freq2 :: freq2 [in] The square of internal tides frequency [T-2 ~> s-2].

Call to:

merid_flux_en ppm_reconstruction_y reflect turning_latitude

Called from:

propagate

[source]

subroutine mom_internal_tides/zonal_flux_en(u, h, hL, hR, uh, dt, G, US, j, ish, ieh, vol_CFL)

Evaluates the zonal mass or volume fluxes in a layer.

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

  • u :: u [in] The zonal velocity [L T-1 ~> m s-1].

  • h :: h [in] Energy density used to calculate the fluxes [H Z2 T-2 ~> m3 s-2 or J m-2].

  • hl :: [in] Left- Energy densities in the reconstruction [H Z2 T-2 ~> m3 s-2 or J m-2].

  • hr :: [in] Right- Energy densities in the reconstruction [H Z2 T-2 ~> m3 s-2 or J m-2].

  • uh :: uh [out] The zonal energy transport [H Z2 L2 T-3 ~> m5 s-3 or J s-1].

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

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

  • j :: j [in] The j-index to work on.

  • ish :: ish [in] The start i-index range to work on.

  • ieh :: ieh [in] The end i-index range to work on.

  • vol_cfl :: [in] If true, rescale the ratio of face areas to the cell areas when estimating the CFL number.

Called from:

propagate_x

[source]

subroutine mom_internal_tides/merid_flux_en(v, h, hL, hR, vh, dt, G, US, J, ish, ieh, vol_CFL)

Evaluates the meridional mass or volume fluxes in a layer.

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

  • v :: v [in] The meridional velocity [L T-1 ~> m s-1].

  • h :: h [in] Energy density used to calculate the fluxes [H Z2 T-2 ~> m3 s-2 or J m-2].

  • hl :: [in] Left- Energy densities in the reconstruction [H Z2 T-2 ~> m3 s-2 or J m-2].

  • hr :: [in] Right- Energy densities in the reconstruction [H Z2 T-2 ~> m3 s-2 or J m-2].

  • vh :: vh [out] The meridional energy transport [H Z2 L2 T-3 ~> m5 s-3 or J s-1].

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

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

  • j :: [in] The j-index to work on.

  • ish :: ish [in] The start i-index range to work on.

  • ieh :: ieh [in] The end i-index range to work on.

  • vol_cfl :: [in] If true, rescale the ratio of face areas to the cell areas when estimating the CFL number.

Called from:

propagate_y

[source]

subroutine mom_internal_tides/reflect(En, NAngle, CS, G, LB)

Reflection of the internal waves at a single frequency.

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

  • nangle :: [in] The number of wave orientations in the discretized wave energy spectrum.

  • en :: [inout] The internal gravity wave energy density as a

  • cs :: [in] Internal tide control structure

  • lb :: [in] A structure with the active energy loop bounds.

Called from:

propagate_x propagate_y

[source]

subroutine mom_internal_tides/turning_latitude(En, NAngle, freq2, CS, G, LB)
Parameters:
  • g :: [in] The ocean’s grid structure

  • nangle :: [in] The number of wave orientations in the discretized wave energy spectrum.

  • en :: [inout] The internal gravity wave energy density as a

  • cs :: [in] Internal tide control structure

  • lb :: [in] A structure with the active energy loop bounds.

  • freq2 :: freq2 [in] The square of the internal tide frequency [T-2 ~> s-2]

Called from:

propagate_x propagate_y

[source]

subroutine mom_internal_tides/teleport(En, NAngle, CS, G, LB)

Moves energy across lines of partial reflection to prevent reflection of energy that is supposed to get across.

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

  • nangle :: [in] The number of wave orientations in the discretized wave energy spectrum.

  • en :: [inout] The internal gravity wave energy density as a

  • cs :: [in] Internal tide control structure

  • lb :: [in] A structure with the active energy loop bounds.

Call to:

mom_error_handler::mom_error

[source]

subroutine mom_internal_tides/correct_halo_rotation(En, test, G, NAngle, halo)

Rotates points in the halos where required to accommodate changes in grid orientation, such as at the tripolar fold.

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

  • en :: [inout] The internal gravity wave energy density as a function of space, angular orientation, frequency, and vertical mode [H Z2 T-2 ~> m3 s-2 or J m-2].

  • test :: test [in] An x-unit vector that has been passed through

  • nangle :: [in] The number of wave orientations in the discretized wave energy spectrum.

  • halo :: halo [in] The halo size over which to do the calculations

Call to:

mom_error_handler::mom_error

Called from:

propagate_int_tide

[source]

subroutine mom_internal_tides/correct_halo_rotation_2d(En, test, G, NAngle, halo)

Rotates points in the halos where required to accommodate changes in grid orientation, such as at the tripolar fold.

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

  • en :: [inout] The internal gravity wave energy density as a function of space, angular orientation, frequency, and vertical mode [H Z2 T-2 ~> m3 s-2 or J m-2].

  • test :: test [in] An x-unit vector that has been passed through

  • nangle :: [in] The number of wave orientations in the discretized wave energy spectrum.

  • halo :: halo [in] The halo size over which to do the calculations

Call to:

mom_error_handler::mom_error

Called from:

propagate

[source]

subroutine mom_internal_tides/ppm_reconstruction_x(h_in, h_l, h_r, G, LB, simple_2nd, adv_limiter)

Calculates left/right edge values for PPM reconstruction in x-direction.

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

  • h_in :: h_in [in] Energy density in a sector (2D) [H Z2 T-2 ~> m3 s-2 or J m-2]

  • h_l :: h_l [out] Left edge value of reconstruction (2D) [H Z2 T-2 ~> m3 s-2 or J m-2]

  • h_r :: h_r [out] Right edge value of reconstruction (2D) [H Z2 T-2 ~> m3 s-2 or J m-2]

  • lb :: [in] A structure with the active loop bounds.

  • simple_2nd :: simple_2nd [in] If true, use the arithmetic mean energy densities as default edge values for a simple 2nd order scheme.

  • adv_limiter :: adv_limiter [in] The type of limiter used

Call to:

limiter_adv_minmod limiter_adv_positive minmod_limiter mom_error_handler::mom_error ppm_limit_pos

Called from:

propagate_x

[source]

subroutine mom_internal_tides/ppm_reconstruction_y(h_in, h_l, h_r, G, LB, simple_2nd, adv_limiter)

Calculates left/right edge valus for PPM reconstruction in y-direction.

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

  • h_in :: h_in [in] Energy density in a sector (2D) [H Z2 T-2 ~> m3 s-2 or J m-2]

  • h_l :: h_l [out] Left edge value of reconstruction (2D) [H Z2 T-2 ~> m3 s-2 or J m-2]

  • h_r :: h_r [out] Right edge value of reconstruction (2D) [H Z2 T-2 ~> m3 s-2 or J m-2]

  • lb :: [in] A structure with the active loop bounds.

  • simple_2nd :: simple_2nd [in] If true, use the arithmetic mean energy densities as default edge values for a simple 2nd order scheme.

  • adv_limiter :: adv_limiter [in] The type of limiter used

Call to:

limiter_adv_minmod limiter_adv_positive minmod_limiter mom_error_handler::mom_error ppm_limit_pos

Called from:

propagate_y

[source]

subroutine mom_internal_tides/ppm_limit_pos(h_in, h_L, h_R, h_min, G, iis, iie, jis, jie)

Limits the left/right edge values of the PPM reconstruction to give a reconstruction that is positive-definite. Here this is reinterpreted as giving a constant value if the mean value is less than h_min, with a minimum of h_min otherwise.

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

  • h_in :: h_in [in] Energy density in each sector (2D) [H Z2 T-2 ~> m3 s-2 or J m-2]

  • h_l :: [inout] Left edge value of reconstruction [H Z2 T-2 ~> m3 s-2 or J m-2]

  • h_r :: [inout] Right edge value of reconstruction [H Z2 T-2 ~> m3 s-2 or J m-2]

  • h_min :: h_min [in] The minimum value that can be obtained by a concave parabolic fit [H Z2 T-2 ~> m3 s-2 or J m-2]

  • iis :: iis [in] Start i-index for computations

  • iie :: iie [in] End i-index for computations

  • jis :: jis [in] Start j-index for computations

  • jie :: jie [in] End j-index for computations

Called from:

ppm_reconstruction_x ppm_reconstruction_y

[source]

subroutine mom_internal_tides/minmod_limiter(h_in, h_L, h_R, G, iis, iie, jis, jie)

Limits the left/right edge values using the simple minmod limiter written in a way that avoids branching in favor of intrinsics.

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

  • h_in :: h_in [in] Energy density in each sector (2D) [H Z2 T-2 ~> m3 s-2 or J m-2]

  • h_l :: [inout] Left edge value of reconstruction [H Z2 T-2 ~> m3 s-2 or J m-2]

  • h_r :: [inout] Right edge value of reconstruction [H Z2 T-2 ~> m3 s-2 or J m-2]

  • iis :: iis [in] Start i-index for computations

  • iie :: iie [in] End i-index for computations

  • jis :: jis [in] Start j-index for computations

  • jie :: jie [in] End j-index for computations

Called from:

ppm_reconstruction_x ppm_reconstruction_y

[source]

subroutine mom_internal_tides/register_int_tide_restarts(G, GV, US, param_file, CS, restart_CS)
Parameters:
  • g :: [in] The ocean’s grid structure

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

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

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

  • cs :: Internal tide control structure

  • restart_cs :: MOM restart control structure

Call to:

mom_error_handler::mom_error mom_io::set_axis_info

[source]

subroutine mom_internal_tides/internal_tides_init(Time, G, GV, US, param_file, diag, CS)

This subroutine initializes the internal tides 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 :: param_file [in] A structure to parse for run-time parameters.

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

  • cs :: Internal tide control structure

Call to:

mom_diag_mediator::define_axes_group mom_string_functions::extract_real limiter_adv_minmod limiter_adv_minmod_string limiter_adv_positive limiter_adv_positive_string mom_error_handler::mom_error mom_error_handler::mom_mesg mom_string_functions::uppercase mom_wave_speed::wave_speed_init

Called from:

mom_diabatic_driver::diabatic_driver_init

[source]

subroutine mom_internal_tides/internal_tides_end(CS)

This subroutine deallocates the memory associated with the internal tides control structure.

Parameters:

cs :: [inout] Internal tide control structure

[source]

[source]