mom_wave_interface module reference

Interface for surface waves.

More…

Data Types

wave_parameters_cs

Container for all surface wave related parameters.

Functions/Subroutines

mom_wave_interface_init()

Initializes parameters related to MOM_wave_interface.

set_lf17_wave_params()

Set the parameters that are used to determine the averaged Stokes drift and Langmuir numbers.

query_wave_properties()

This interface provides the caller with information from the waves control structure.

update_surface_waves()

Subroutine that handles updating of surface wave/Stokes drift related properties.

update_stokes_drift()

Constructs the Stokes Drift profile on the model grid based on desired coupling options.

one_minus_exp_x()

Return the value of (1 - exp(-x))/x, using an accurate expression for small values of x.

surface_bands_by_data_override()

A subroutine to fill the Stokes drift from a NetCDF file using the data_override procedures.

get_langmuir_number()

Interface to get Langmuir number based on options stored in wave structure.

get_wave_method()

function to return the wave method string set in the param file

get_stokessl_lifoxkemper()

Get SL averaged Stokes drift from Li/FK 17 method.

get_sl_average_prof()

Get SL Averaged Stokes drift from a Stokes drift Profile.

get_sl_average_band()

Get SL averaged Stokes drift from the banded Spectrum method.

dhh85_mid()

Compute the Stokes drift at a given depth.

stokesmixing()

Explicit solver for Stokes mixing.

coriolisstokes()

Solver to add Coriolis-Stokes to model Still in development and not meant for general use.

stokes_pgf()

Computes tendency due to Stokes pressure gradient force anomaly including analytical integration of Stokes shear using multiple-exponential decay Stokes drift profile and vertical integration of the resulting pressure anomaly to the total pressure gradient force.

ust_2_u10_coare3p5()

Computes wind speed from ustar_air based on COARE 3.5 Cd relationship Probably doesn’t belong in this module, but it is used here to estimate wind speed for wind-wave relationships.

waves_end()

Clear pointers, deallocate memory.

waves_register_restarts()

Register wave restart fields.

Detailed Description

This module should be moved as wave coupling progresses and likely will should mirror the iceberg or sea-ice model set-up.

This module is meant to contain the routines to read in and interpret surface wave data for MOM6. In its original form, the capabilities include setting the Stokes drift in the model (from a variety of sources including prescribed, empirical, and input files). In short order, the plan is to also amend the subroutine to accept Stokes drift information from an external coupler. Eventually, it will be necessary to break this file apart so that general wave information may be stored in the control structure and the Stokes drift effect can be isolated from processes such as sea-state dependent momentum fluxes, gas fluxes, and other wave related air-sea interaction and boundary layer phenomenon.

The Stokes drift are stored on the C-grid with the conventional protocol to interpolate to the h-grid to compute Langmuir number, the primary quantity needed for Langmuir turbulence parameterizations in both the ePBL and KPP approach. This module also computes full 3d Stokes drift profiles, which will be useful if second-order type boundary layer parameterizations are implemented (perhaps via GOTM, work in progress).

Type Documentation

type mom_wave_interface/wave_parameters_cs

Container for all surface wave related parameters.

Type fields:
  • % id_pfu_stokes [integer,public] :: Diagnostic handles.

  • % id_pfv_stokes [integer,public] :: Diagnostic handles.

  • % id_3dstokes_x_from_ddt [integer,public] :: Diagnostic handles.

  • % id_3dstokes_y_from_ddt [integer,public] :: Diagnostic handles.

  • % id_p_deltastokes_l [integer] :: Diagnostic handles.

  • % id_p_deltastokes_i [integer] :: Diagnostic handles.

  • % id_surfacestokes_x [integer] :: Diagnostic handles.

  • % id_surfacestokes_y [integer] :: Diagnostic handles.

  • % id_3dstokes_x [integer] :: Diagnostic handles.

  • % id_3dstokes_y [integer] :: Diagnostic handles.

  • % id_ddt_3dstokes_x [integer] :: Diagnostic handles.

  • % id_ddt_3dstokes_y [integer] :: Diagnostic handles.

  • % id_la_turb [integer] :: Diagnostic handles.

  • % usewaves [logical,public] :: Flag to enable surface gravity wave feature.

  • % stokes_vf [logical,public] :: True if Stokes vortex force is used.

  • % passive_stokes_vf [logical,public] :: Computes Stokes VF, but doesn’t affect dynamics.

  • % stokes_pgf [logical,public] :: True if Stokes shear pressure Gradient force is used.

  • % passive_stokes_pgf [logical,public] :: Keeps Stokes_PGF on, but doesn’t affect dynamics.

  • % stokes_ddt [logical,public] :: Developmental: True if Stokes d/dt is used.

  • % passive_stokes_ddt [logical,public] :: Keeps Stokes_DDT on, but doesn’t affect dynamics.

  • % us_x [real(:,:,:),allocatable, public] :: 3d zonal Stokes drift profile [L T-1 ~> m s-1]

  • % us_y [real(:,:,:),allocatable, public] :: 3d meridional Stokes drift profile [L T-1 ~> m s-1]

  • % ddt_us_x [real(:,:,:),allocatable, public] :: 3d time tendency of zonal Stokes drift profile [L T-2 ~> m s-2]

  • % ddt_us_y [real(:,:,:),allocatable, public] :: 3d time tendency of meridional Stokes drift profile [L T-2 ~> m s-2]

  • % us_x_from_ddt [real(:,:,:),allocatable, public] :: Check of 3d zonal Stokes drift profile [L T-1 ~> m s-1].

  • % us_y_from_ddt [real(:,:,:),allocatable, public] :: Check of 3d meridional Stokes drift profile [L T-1 ~> m s-1].

  • % us_x_prev [real(:,:,:),allocatable, public] :: 3d zonal Stokes drift profile, previous dynamics call [L T-1 ~> m s-1]

  • % us_y_prev [real(:,:,:),allocatable, public] :: 3d meridional Stokes drift profile, previous dynamics call [L T-1 ~> m s-1]

  • % kvs [real(:,:,:),allocatable, public] :: Viscosity for Stokes Drift shear [H Z T-1 ~> m2 s-1 or Pa s].

  • % wavenum_cen [real(:),allocatable, public] :: Wavenumber bands for read/coupled [Z-1 ~> m-1].

  • % ustk_hb [real(:,:,:),allocatable, public] :: Surface Stokes Drift spectrum (zonal) [L T-1 ~> m s-1].

  • % vstk_hb [real(:,:,:),allocatable, public] :: Surface Stokes Drift spectrum (meridional) [L T-1 ~> m s-1].

  • % omega_w2x [real(:,:),allocatable, public] :: wind direction ccw from model x- axis [nondim radians]

  • % numbands [integer,public] :: Number of wavenumber/frequency partitions Must match the number of bands provided via either coupling or file.

  • % wavemethod [integer] :: Options for including wave information Valid (tested) choices are: 0 - Test Profile 1 - Surface Stokes Drift Bands 2 - DHH85 3 - LF17 -99 - No waves computed, but empirical Langmuir number used.

  • % lagrangianmixing [logical] :: This feature is in development and not ready True if Stokes drift is present and mixing should be applied to Lagrangian current (mean current + Stokes drift). See Reichl et al., 2016 KPP-LT approach.

  • % stokesmixing [logical] :: This feature is in development and not ready. True if vertical mixing of momentum should be applied directly to Stokes current (with separate mixing parameter for Eulerian mixing contribution). See Harcourt 2013, 2015 Second-Moment approach.

  • % coriolisstokes [logical] :: This feature is in development and not ready.

  • % stokes_min_thick_avg [real] :: A layer thickness below which the cell-center Stokes drift is used instead of the cell average [Z ~> m]. This is only used if WAVE_INTERFACE_ANSWER_DATE < 20230101.

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

  • % partitionmode [integer] :: Method for partition mode (meant to check input) 0 - wavenumbers 1 - frequencies.

  • % datasource [integer] :: Integer that specifies where the model Looks for data Valid choices are: 1 - FMS DataOverride Routine 2 - Reserved For Coupler 3 - User input (fixed values, useful for 1d testing)

  • % surfbandfilename [character (len=40)] :: Filename if using DataOverride.

  • % land_speed [real] :: A large Stokes velocity that can be used to indicate land values in a data override file [L T-1 ~> m s-1]. Stokes drift components larger than this are set to zero in data override calls for the Stokes drift.

  • % dataover_initialized [logical] :: Flag for DataOverride Initialization.

  • % la_frachbl [real] :: Fraction of OSBL for averaging Langmuir number [nondim].

  • % la_hbl_min [real] :: Minimum boundary layer depth for averaging Langmuir number [Z ~> m].

  • % la_misalignment [logical] :: Flag to use misalignment in Langmuir number.

  • % g_earth [real] :: The gravitational acceleration, equivalent to GVg_Earth but with different dimensional rescaling appropriate for deep-water gravity waves [Z T-2 ~> m s-2].

  • % i_g_earth [real] :: The inverse of the gravitational acceleration, with dimensional rescaling appropriate for deep-water gravity waves [T2 Z-1 ~> s2 m-1].

  • % freq_cen [real(:),allocatable] :: Central frequency for wave bands, including a factor of 2*pi [T-1 ~> s-1].

  • % prescribedsurfstkx [real(:),allocatable] :: Surface Stokes drift if prescribed [L T-1 ~> m s-1].

  • % prescribedsurfstky [real(:),allocatable] :: Surface Stokes drift if prescribed [L T-1 ~> m s-1].

  • % la_turb [real(:,:),allocatable] :: Aligned Turbulent Langmuir number [nondim].

  • % us0_x [real(:,:),allocatable] :: Surface Stokes Drift (zonal) [L T-1 ~> m s-1].

  • % us0_y [real(:,:),allocatable] :: Surface Stokes Drift (meridional) [L T-1 ~> m s-1].

  • % stkx0 [real(:,:,:),allocatable] :: Stokes Drift spectrum (zonal) [L T-1 ~> m s-1].

  • % stky0 [real(:,:,:),allocatable] :: Stokes Drift spectrum (meridional) [L T-1 ~> m s-1].

  • % la_min [real] :: An arbitrary lower-bound on the Langmuir number [nondim]. Langmuir number is sqrt(u_star/u_stokes). When both are small but u_star is orders of magnitude smaller, the Langmuir number could have unintended consequences. Since both are small it can be safely capped to avoid such consequences.

  • % la_stk_backgnd [real] :: A small background Stokes velocity used in the denominator of some expressions for the Langmuir number [L T-1 ~> m s-1].

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

  • % rho_air [real] :: A typical density of air at sea level, as used in wave calculations [R ~> kg m-3].

  • % nu_air [real] :: The viscosity of air, as used in wave calculations [Z2 T-1 ~> m2 s-1].

  • % rho_ocn [real] :: A typical surface density of seawater, as used in wave calculations in comparison with the density of air [R ~> kg m-3]. The default is RHO_0.

  • % swh_from_u10sq [real] :: A factor for converting the square of the 10 m wind speed to the significant wave height [Z T2 L-2 ~> s2 m-1].

  • % charnock_min [real] :: The minimum value of the Charnock coefficient, which relates the square of the air friction velocity divided by the gravitational acceleration to the wave roughness length [nondim].

  • % charnock_slope_u10 [real] :: The partial derivative of the Charnock coefficient with the 10 m wind speed [T L-1 ~> s m-1]. Note that in eq. 13 of the Edson et al. 2013 describing the COARE 3.5 bulk flux algorithm, this slope is given as 0.017. However, 0.0017 reproduces the curve in their figure 6, so that is the default value used in MOM6.

  • % charnock_intercept [real] :: The intercept of the fit for the Charnock coefficient in the limit of no wind [nondim]. Note that this can be negative because CHARNOCK_MIN will keep the final value for the Charnock coefficient from being from being negative.

  • % tp_stkx0 [real] :: Test profile x-stokes drift amplitude [L T-1 ~> m s-1].

  • % tp_stky0 [real] :: Test profile y-stokes drift amplitude [L T-1 ~> m s-1].

  • % tp_wvl [real] :: Test profile wavelength [Z ~> m].

  • % waveagepeakfreq [logical] :: Flag to use wave age to determine the peak frequency with DHH85.

  • % staticwaves [logical] :: Flag to disable updating DHH85 Stokes drift.

  • % dhh85_is_set [logical] :: The if the wave properties have been set when WaveMethod = DHH85.

  • % waveage [real] :: The fixed wave age used with the DHH85 spectrum [nondim].

  • % wavewind [real] :: Wind speed for the DHH85 spectrum [L T-1 ~> m s-1].

  • % omega_min [real] :: Minimum wave frequency with the DHH85 spectrum [T-1 ~> s-1].

  • % omega_max [real] :: Maximum wave frequency with the DHH85 spectrum [T-1 ~> s-1].

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

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

Function/Subroutine Documentation

subroutine mom_wave_interface/mom_wave_interface_init(time, G, GV, US, param_file, CS, diag)

Initializes parameters related to MOM_wave_interface.

Parameters:
  • time :: [in] Model time

  • g :: [inout] Grid structure

  • gv :: [in] Vertical grid structure

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

  • param_file :: [in] Input parameter structure

  • cs :: Wave parameter control structure

  • diag :: [inout] Diagnostic Pointer

Call to:

coupler dataovr dhh85 dhh85_string efactor efactor_string input lf17 lf17_string mom_error_handler::mom_error null_string null_wavemethod set_lf17_wave_params surfbands surfbands_string testprof testprof_string

Called from:

mom6

subroutine mom_wave_interface/set_lf17_wave_params(param_file, mdl, GV, US, CS)

Set the parameters that are used to determine the averaged Stokes drift and Langmuir numbers.

Parameters:
  • param_file :: [in] Input parameter structure

  • mdl :: [in] A module name to use in the get_param calls

  • gv :: [in] Vertical grid structure

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

  • cs :: Wave parameter control structure

Called from:

mom_wave_interface_init

subroutine mom_wave_interface/query_wave_properties(CS, NumBands, WaveNumbers, US)

This interface provides the caller with information from the waves control structure.

Parameters:
  • cs :: Wave parameter Control structure

  • numbands :: [out] If present, this returns the number of < wavenumber partitions in the wave discretization

  • wavenumbers :: [out] If present this returns the characteristic wavenumbers of the wave discretization [m-1] or [Z-1 ~> m-1]

  • us :: [in] A dimensional unit scaling type that is used to undo the dimensional scaling of the output variables, if present

Call to:

mom_error_handler::mom_error

subroutine mom_wave_interface/update_surface_waves(G, GV, US, Time_present, dt, CS, forces)

Subroutine that handles updating of surface wave/Stokes drift related properties.

Parameters:
  • cs :: Wave parameter Control structure

  • g :: [inout] Grid structure

  • gv :: [in] Vertical grid structure

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

  • time_present :: [in] Model Time

  • dt :: [in] Time increment as a time-type

  • forces :: [in] MOM_forcing_type

Call to:

coupler dataovr input mom_error_handler::mom_error surface_bands_by_data_override surfbands testprof

Called from:

mom6 ocean_model_mod::update_ocean_model

subroutine mom_wave_interface/update_stokes_drift(G, GV, US, CS, dz, ustar, dt, dynamics_step)

Constructs the Stokes Drift profile on the model grid based on desired coupling options.

Parameters:
  • cs :: Wave parameter Control structure

  • g :: [inout] Grid structure

  • gv :: [in] Vertical grid structure

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

  • dz :: [in] Thickness in height units [Z ~> m]

  • ustar :: [in] Wind friction velocity [Z T-1 ~> m s-1].

  • dt :: [in] Time-step for computing Stokes-tendency [T ~> s]

  • dynamics_step :: [in] True if this call is on a dynamics step

Call to:

dhh85 dhh85_mid efactor get_langmuir_number one_minus_exp_x surfbands testprof

Called from:

mom::step_mom

function mom_wave_interface/one_minus_exp_x(x) [real]

Return the value of (1 - exp(-x))/x, using an accurate expression for small values of x.

Parameters:

x :: [in] The argument of the function ((1 - exp(-x))/x) [nondim]

Called from:

get_stokessl_lifoxkemper update_stokes_drift

subroutine mom_wave_interface/surface_bands_by_data_override(Time, G, GV, US, CS)

A subroutine to fill the Stokes drift from a NetCDF file using the data_override procedures.

Parameters:
  • time :: [in] Time to get Stokes drift bands

  • cs :: Wave structure

  • g :: [inout] Grid structure

  • gv :: [in] Vertical grid structure

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

Call to:

mom_error_handler::mom_error

Called from:

update_surface_waves

subroutine mom_wave_interface/get_langmuir_number(LA, G, GV, US, HBL, ustar, i, j, dz, Waves, U_H, V_H, Override_MA)

Interface to get Langmuir number based on options stored in wave structure.

Note this can be called with an unallocated Waves pointer, which is okay if we want the wind-speed only dependent Langmuir number. Therefore, we need to be careful about what we try to access here.

Parameters:
  • g :: [in] Ocean grid structure

  • gv :: [in] Ocean vertical grid structure

  • la :: [out] Langmuir number [nondim]

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

  • hbl :: [in] (Positive) thickness of boundary layer [Z ~> m]

  • ustar :: [in] Friction velocity [Z T-1 ~> m s-1]

  • i :: [in] Meridional index of h-point

  • j :: [in] Zonal index of h-point

  • dz :: [in] Grid layer thickness [Z ~> m]

  • waves :: Surface wave control structure.

  • u_h :: [in] Zonal velocity at H point [L T-1 ~> m s-1] or [m s-1]

  • v_h :: [in] Meridional velocity at H point [L T-1 ~> m s-1] or [m s-1]

  • override_ma :: [in] Override to use misalignment in LA calculation. This can be used if diagnostic LA outputs are desired that are different than those used by the dynamical model.

Call to:

dhh85 get_sl_average_band get_sl_average_prof get_stokessl_lifoxkemper lf17 mom_error_handler::mom_error null_wavemethod surfbands testprof

Called from:

mom_energetic_pbl::epbl_column mom_cvmix_kpp::kpp_compute_bld update_stokes_drift

function mom_wave_interface/get_wave_method(CS) [character(:)]

function to return the wave method string set in the param file

Parameters:

cs :: Control structure

Call to:

dhh85 dhh85_string efactor efactor_string lf17 lf17_string null_string null_wavemethod surfbands surfbands_string testprof testprof_string

subroutine mom_wave_interface/get_stokessl_lifoxkemper(ustar, hbl, GV, US, CS, UStokes_SL, LA)

Get SL averaged Stokes drift from Li/FK 17 method.

Original description: * This function returns the enhancement factor, given the 10-meter wind [m s-1], friction velocity [m s-1] and the boundary layer depth [m].

Update (Jan/25): * Converted from function to subroutine, now returns Langmuir number.

  • Compute 10m wind internally, so only ustar and hbl need passed to subroutine.

Qing Li, 160606 * BGR port from CVMix to MOM6 Jan/25/2017

  • BGR change output to LA from Efactor

  • BGR remove u10 input

  • BGR note: fixed parameter values should be changed to “get_params”

Parameters:
  • ustar :: [in] water-side surface friction velocity [Z T-1 ~> m s-1].

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

  • gv :: [in] Ocean vertical grid structure

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

  • cs :: Wave parameter Control structure

  • ustokes_sl :: [out] Surface layer averaged Stokes drift [L T-1 ~> m s-1]

  • la :: [out] Langmuir number [nondim]

Call to:

one_minus_exp_x ust_2_u10_coare3p5

Called from:

get_langmuir_number

subroutine mom_wave_interface/get_sl_average_prof(GV, AvgDepth, dz, Profile, Average)

Get SL Averaged Stokes drift from a Stokes drift Profile.

Parameters:
  • gv :: [in] Ocean vertical grid structure

  • avgdepth :: [in] Depth to average over (negative) [Z ~> m]

  • dz :: [in] Grid thickness [Z ~> m]

  • profile :: [in] Profile of quantity to be averaged in arbitrary units [A]

  • average :: [out] Output quantity averaged over depth AvgDepth [A] (used here for Stokes drift)

Called from:

get_langmuir_number

subroutine mom_wave_interface/get_sl_average_band(GV, AvgDepth, NB, WaveNumbers, SurfStokes, Average)

Get SL averaged Stokes drift from the banded Spectrum method.

Parameters:
  • gv :: [in] Ocean vertical grid

  • avgdepth :: [in] Depth to average over [Z ~> m].

  • nb :: [in] Number of bands used

  • wavenumbers :: [in] Wavenumber corresponding to each band [Z-1 ~> m-1]

  • surfstokes :: [in] Surface Stokes drift for each band [L T-1 ~> m s-1]

  • average :: [out] Output average Stokes drift over depth AvgDepth [L T-1 ~> m s-1]

Called from:

get_langmuir_number

subroutine mom_wave_interface/dhh85_mid(GV, US, CS, zpt, UStokes)

Compute the Stokes drift at a given depth.

Taken from Qing Li (Brown) use for comparing MOM6 simulation to his LES computed at z mid point (I think) and not depth averaged. Should be fine to integrate in frequency from 0.1 to sqrt(-0.2*grav*2pi/dz

Parameters:
  • gv :: [in] Ocean vertical grid

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

  • cs :: Wave parameter Control structure

  • zpt :: [in] Depth to get Stokes drift [Z ~> m].

  • ustokes :: [out] Stokes drift [L T-1 ~> m s-1]

Called from:

update_stokes_drift

subroutine mom_wave_interface/stokesmixing(G, GV, dt, h, dz, u, v, Waves)

Explicit solver for Stokes mixing. Still in development do not use.

Parameters:
  • g :: [in] Ocean grid

  • gv :: [in] Ocean vertical grid

  • dt :: [in] Time step of MOM6 [T ~> s] for explicit solver

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

  • dz :: [in] Vertical distance between interfaces around a layer [Z ~> m]

  • u :: [inout] Velocity i-component [L T-1 ~> m s-1]

  • v :: [inout] Velocity j-component [L T-1 ~> m s-1]

  • waves :: Surface wave related control structure.

subroutine mom_wave_interface/coriolisstokes(G, GV, dt, h, u, v, Waves)

Solver to add Coriolis-Stokes to model Still in development and not meant for general use. Can be activated (with code intervention) for LES comparison CHECK THAT RIGHT TIMESTEP IS PASSED IF YOU USE THIS**.

Not accessed in the standard code.

Parameters:
  • g :: [in] Ocean grid

  • gv :: [in] Ocean vertical grid

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

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

  • u :: [inout] Velocity i-component [L T-1 ~> m s-1]

  • v :: [inout] Velocity j-component [L T-1 ~> m s-1]

  • waves :: Surface wave related control structure.

subroutine mom_wave_interface/stokes_pgf(G, GV, US, dz, u, v, PFu_Stokes, PFv_Stokes, CS)

Computes tendency due to Stokes pressure gradient force anomaly including analytical integration of Stokes shear using multiple-exponential decay Stokes drift profile and vertical integration of the resulting pressure anomaly to the total pressure gradient force.

Parameters:
  • g :: [in] Ocean grid

  • gv :: [in] Ocean vertical grid

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

  • dz :: [in] Layer thicknesses in height units [Z ~> m]

  • u :: [in] Lagrangian Velocity i-component [L T-1 ~> m s-1]

  • v :: [in] Lagrangian Velocity j-component [L T-1 ~> m s-1]

  • pfu_stokes :: [out] PGF Stokes-shear i-component [L T-2 ~> m s-2]

  • pfv_stokes :: [out] PGF Stokes-shear j-component [L T-2 ~> m s-2]

  • cs :: Surface wave related control structure.

Called from:

mom_dynamics_split_rk2::step_mom_dyn_split_rk2 mom_dynamics_split_rk2b::step_mom_dyn_split_rk2b

subroutine mom_wave_interface/ust_2_u10_coare3p5(USTair, U10, GV, US, CS)

Computes wind speed from ustar_air based on COARE 3.5 Cd relationship Probably doesn’t belong in this module, but it is used here to estimate wind speed for wind-wave relationships. Should be a fine way to estimate the neutral wind-speed as written here.

Parameters:
  • ustair :: [in] Wind friction velocity [Z T-1 ~> m s-1]

  • u10 :: [out] 10-m neutral wind speed [L T-1 ~> m s-1]

  • gv :: [in] vertical grid type

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

  • cs :: Wave parameter Control structure

Call to:

mom_error_handler::mom_error

Called from:

get_stokessl_lifoxkemper

subroutine mom_wave_interface/waves_end(CS)

Clear pointers, deallocate memory.

Parameters:

cs :: Control structure

subroutine mom_wave_interface/waves_register_restarts(CS, HI, GV, US, param_file, restart_CSp)

Register wave restart fields. To be called before MOM_wave_interface_init.

Parameters:
  • cs :: Wave parameter Control structure

  • hi :: [inout] Grid structure

  • gv :: [in] Vertical grid structure

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

  • param_file :: [in] Input parameter structure

  • restart_csp :: Restart structure, data intent(inout)

Call to:

mom_error_handler::mom_error mom_io::var_desc