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.

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_corner_spread()

This subroutine does first-order corner advection.

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.

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.

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.

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(:),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(:,:),allocatable] :: Diag handles.

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

  • % id_int_n2w2_mode [integer(:),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.

  • % 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.

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

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

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

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

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

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

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

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

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

  • % tke_leak_loss [real(:,:,:,:,:),allocatable] :: energy lost due to misc background processes [R Z3 T-3 ~> W m-2]

  • % tke_quad_loss [real(:,:,:,:,:),allocatable] :: energy lost due to quadratic bottom drag [R Z3 T-3 ~> W m-2]

  • % tke_froude_loss [real(:,:,:,:,:),allocatable] :: energy lost due to wave breaking [R Z3 T-3 ~> W m-2]

  • % tke_itidal_loss_fixed [real(:,:),allocatable] :: Fixed part of the energy lost due to small-scale drag [R Z3 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(:,:,:,:,:),allocatable] :: energy lost due to small-scale wave drag [R Z3 T-3 ~> W m-2]

  • % tke_residual_loss [real(:,:,:,:,:),allocatable] :: internal tide energy loss due to the residual at slopes [R Z3 T-3 ~> W m-2]

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

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

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

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

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

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

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

  • % u_struct [real(:,:,:,:),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(:,:,:),allocatable] :: Maximum of u_struct, for each mode [Z-1 ~> m-1].

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

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

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

  • % int_n2w2 [real(:,:,:),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]

  • % en_sum [real] :: global sum of energy for use in debugging, in MKS units [J]

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

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

  • % decay_rate [real] :: A constant rate at which internal tide energy is lost to the interior ocean internal wave field [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].

  • % 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 an imbalance in the internal tide energy budget when apply_Froude_drag is True [R Z3 T-2 ~> J m-2].

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

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

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

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

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

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

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

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

  • % 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.

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.

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 :: [in] Layer thicknesses [H ~> m or kg m-2]

  • 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 :: [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 itidal_lowmode_loss mom_error_handler::mom_error propagate refract sum_en mom_wave_speed::wave_speeds

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

Checks for energy conservation on computational domain.

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

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

  • cs :: [inout] Internal tide control structure

  • en :: [in] The energy density of the internal tides [R Z3 T-2 ~> J m-2].

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

Call to:

mom_spatial_means::global_area_integral

Called from:

propagate_int_tide

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 [R Z3 T-2 ~> J m-2].

  • tke_loss :: [out] Energy loss rate [R Z3 T-3 ~> W m-2]

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

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

Called from:

propagate_int_tide

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 :: [in] The i-index of the value to be reported.

  • 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 :: [in] The named mechanism of loss to return

  • tke_loss_sum :: [out] Total energy loss rate due to specified mechanism [R Z3 T-3 ~> W m-2].

Called from:

mom_tidal_mixing::add_int_tide_diffusivity

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 :: [in] Baroclinic mode speed [L T-1 ~> m s-1].

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

  • 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

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 :: [in] Time increment [T ~> s].

  • 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 [R Z3 T-2 ~> J m-2].

Called from:

refract

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

Propagates internal waves at a single frequency.

Parameters:
  • g :: [inout] 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 :: [in] Baroclinic mode speed [L T-1 ~> m s-1].

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

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

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

  • cs :: [in] Internal tide control structure

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

Call to:

propagate_corner_spread propagate_x propagate_y

Called from:

propagate_int_tide

subroutine mom_internal_tides/propagate_corner_spread(En, energized_wedge, NAngle, speed, dt, G, CS, LB)

This subroutine does first-order corner advection. It was written with the hopes of smoothing out the garden sprinkler effect, but is too numerically diffusive to be of much use as of yet. It is not yet compatible with reflection schemes (BDM).

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

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

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

  • energized_wedge :: [in] Index of current ray direction.

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

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

  • cs :: [in] Internal tide control structure

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

Call to:

mom_error_handler::mom_error

Called from:

propagate

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

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 :: [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 :: [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 :: [inout] internal tide energy loss due

Call to:

ppm_reconstruction_x reflect zonal_flux_en

Called from:

propagate

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

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 :: [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 :: [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 :: [inout] internal tide energy loss due

Call to:

merid_flux_en ppm_reconstruction_y reflect

Called from:

propagate

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 :: [in] The zonal velocity [L T-1 ~> m s-1].

  • h :: [in] Energy density used to calculate the fluxes [R Z3 T-2 ~> J m-2].

  • hl :: [in] Left- Energy densities in the reconstruction [R Z3 T-2 ~> J m-2].

  • hr :: [in] Right- Energy densities in the reconstruction [R Z3 T-2 ~> J m-2].

  • uh :: [inout] The zonal energy transport [R Z3 L2 T-3 ~> J s-1].

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

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

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

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

  • 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

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 :: [in] The meridional velocity [L T-1 ~> m s-1].

  • h :: [in] Energy density used to calculate the fluxes [R Z3 T-2 ~> J m-2].

  • hl :: [in] Left- Energy densities in the reconstruction [R Z3 T-2 ~> J m-2].

  • hr :: [in] Right- Energy densities in the reconstruction [R Z3 T-2 ~> J m-2].

  • vh :: [inout] The meridional energy transport [R Z3 L2 T-3 ~> J s-1].

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

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

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

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

  • 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

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

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

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 [R Z3 T-2 ~> J m-2].

  • 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 :: [in] The halo size over which to do the calculations

Call to:

mom_error_handler::mom_error

Called from:

propagate_int_tide

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

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

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

  • h_in :: [in] Energy density in a sector (2D) [R Z3 T-2 ~> J m-2]

  • h_l :: [out] Left edge value of reconstruction (2D) [R Z3 T-2 ~> J m-2]

  • h_r :: [out] Right edge value of reconstruction (2D) [R Z3 T-2 ~> J m-2]

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

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

Call to:

mom_error_handler::mom_error ppm_limit_pos

Called from:

propagate_x

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

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

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

  • h_in :: [in] Energy density in a sector (2D) [R Z3 T-2 ~> J m-2]

  • h_l :: [out] Left edge value of reconstruction (2D) [R Z3 T-2 ~> J m-2]

  • h_r :: [out] Right edge value of reconstruction (2D) [R Z3 T-2 ~> J m-2]

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

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

Call to:

mom_error_handler::mom_error ppm_limit_pos

Called from:

propagate_y

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 :: [in] Energy density in each sector (2D) [R Z3 T-2 ~> J m-2]

  • h_l :: [inout] Left edge value of reconstruction [R Z3 T-2 ~> J m-2]

  • h_r :: [inout] Right edge value of reconstruction [R Z3 T-2 ~> J m-2]

  • h_min :: [in] The minimum value that can be obtained by a concave parabolic fit [R Z3 T-2 ~> J m-2]

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

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

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

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

Called from:

ppm_reconstruction_x ppm_reconstruction_y

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

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

  • 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

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 :: [inout] 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 :: [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 mom_error_handler::mom_error mom_error_handler::mom_mesg mom_wave_speed::wave_speed_init

Called from:

mom_diabatic_driver::diabatic_driver_init

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