mom_wave_interface module reference¶
Interface for surface waves.
Data Types¶
Container for all surface wave related parameters. |
Functions/Subroutines¶
Initializes parameters related to MOM_wave_interface. |
|
A ‘lite’ init subroutine to initialize a few inputs needed if using wave information with the wind-speed dependent Stokes drift formulation of LF17. |
|
Subroutine that handles updating of surface wave/Stokes drift related properties. |
|
Constructs the Stokes Drift profile on the model grid based on desired coupling options. |
|
A subroutine to fill the Stokes drift from a NetCDF file using the data_override procedures. |
|
Interface to get Langmuir number based on options stored in wave structure. |
|
Get SL averaged Stokes drift from Li/FK 17 method. |
|
Get SL Averaged Stokes drift from a Stokes drift Profile. |
|
Get SL averaged Stokes drift from the banded Spectrum method. |
|
Compute the Stokes drift at a given depth. |
|
Explicit solver for Stokes mixing. |
|
Solver to add Coriolis-Stokes to model Still in development and not meant for general use. |
|
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. |
|
Clear pointers, deallocate memory. |
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 ammend 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_surfacestokes_x
[integer,public] :: Diagnostic handles.%
id_surfacestokes_y
[integer,public] :: Diagnostic handles.%
id_3dstokes_x
[integer,public] :: Diagnostic handles.%
id_3dstokes_y
[integer,public] :: Diagnostic handles.%
id_la_turb
[integer,public] :: Diagnostic handles.%
usewaves
[logical,public] :: Flag to enable surface gravity wave feature.%
lagrangianmixing
[logical,public] :: 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,public] :: 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,public] :: This feature is in development and not ready.%
stklevelmode
[integer,public] :: Sets if Stokes drift is defined at mid-points or layer averaged. Set to 0 if mid-point and set to 1 if average value of Stokes drift over level. If advecting with Stokes transport, 1 is the correct approach.%
wavenum_cen
[real(:),allocatable, public] :: Wavenumber bands for read/coupled [m-1].%
freq_cen
[real(:),allocatable, public] :: Frequency bands for read/coupled [s-1].%
prescribedsurfstkx
[real(:),allocatable, public] :: Surface Stokes drift if prescribed [m s-1].%
prescribedsurfstky
[real(:),allocatable, public] :: Surface Stokes drift if prescribed [m s-1].%
us_x
[real(:,:,:),allocatable, public] :: 3d zonal Stokes drift profile [m s-1]%
us_y
[real(:,:,:),allocatable, public] :: 3d meridional Stokes drift profile [m s-1]%
la_sl
[real(:,:),allocatable, public] :: SL Langmuir number (directionality factored later)%
la_turb
[real(:,:),allocatable, public] :: Aligned Turbulent Langmuir number.%
us0_x
[real(:,:),allocatable, public] :: Surface Stokes Drift (zonal, m/s)%
us0_y
[real(:,:),allocatable, public] :: Surface Stokes Drift (meridional, m/s)%
stkx0
[real(:,:,:),allocatable, public] :: Stokes Drift spectrum (zonal, m/s)%
stky0
[real(:,:,:),allocatable, public] :: Stokes Drift spectrum (meridional, m/s)%
kvs
[real(:,:,:),allocatable, public] :: Viscosity for Stokes Drift shear [Z2 T-1 ~> m2 s-1].%
time
[type(time_type),pointer, public] :: A pointer to the ocean model’s clock.%
diag
[type(diag_ctrl),pointer, public] :: A structure that is used to regulate the timing of diagnostic output.%
la_min
[real] :: An arbitrary lower-bound on the Langmuir number. Run-time parameter. 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.
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
dataoverrideisinitialized
dataovr
datasource
dhh85
input
la_frachbl
la_misalignment
lf17
mdl
mom_error_handler::mom_error
numbands
partitionmode
pi
staticwaves
surfbandfilename
surfbands
testprof
tp_stkx0
tp_stky0
tp_wvl
waveage
waveagepeakfreq
wavemethod
wavewind
- Called from
mom_main
-
subroutine
mom_wave_interface/
mom_wave_interface_init_lite
(param_file)¶ A ‘lite’ init subroutine to initialize a few inputs needed if using wave information with the wind-speed dependent Stokes drift formulation of LF17.
- Parameters
param_file :: [in] Input parameter structure
- Call to
la_frachbl
lf17
mdl
null_wavemethod
pi
wavemethod
- Called from
mom_main
ocean_model_mod::ocean_model_init
-
subroutine
mom_wave_interface/
update_surface_waves
(G, GV, US, Day, dt, CS)¶ 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
day :: [in] Current model time
dt :: [in] Timestep as a time-type
- Call to
coupler
dataovr
datasource
input
numbands
surface_bands_by_data_override
surfbands
testprof
wavemethod
- Called from
mom_main
ocean_model_mod::update_ocean_model
-
subroutine
mom_wave_interface/
update_stokes_drift
(G, GV, US, CS, h, ustar)¶ 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
h :: [in] Thickness [H ~> m or kg m-2]
ustar :: [in] Wind friction velocity [Z T-1 ~> m s-1].
- Call to
dhh85
dhh85_is_set
dhh85_mid
get_langmuir_number
numbands
partitionmode
pi
staticwaves
surfbands
testprof
tp_stkx0
tp_stky0
tp_wvl
wavemethod
- Called from
-
subroutine
mom_wave_interface/
surface_bands_by_data_override
(day_center, G, GV, US, CS)¶ A subroutine to fill the Stokes drift from a NetCDF file using the data_override procedures.
- Parameters
day_center :: [in] Center of timestep
cs :: Wave structure
g :: [inout] Grid structure
gv :: [in] Vertical grid structure
us :: [in] A dimensional unit scaling type
- Call to
dataoverrideisinitialized
mom_error_handler::mom_error
numbands
partitionmode
pi
surfbandfilename
mom_domains::to_all
- Called from
-
subroutine
mom_wave_interface/
get_langmuir_number
(LA, G, GV, US, HBL, ustar, i, j, H, U_H, V_H, Override_MA, Waves)¶ 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
us :: [in] A dimensional unit scaling type
i :: [in] Meridional index of h-point
j :: [in] Zonal index of h-point
ustar :: [in] Friction velocity [Z T-1 ~> m s-1].
hbl :: [in] (Positive) thickness of boundary layer [Z ~> m].
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.
h :: [in] Grid layer thickness [H ~> m or kg m-2]
u_h :: [in] Zonal velocity at H point [m s-1]
v_h :: [in] Meridional velocity at H point [m s-1]
waves :: Surface wave control structure.
la :: [out] Langmuir number
- Call to
dhh85
get_sl_average_band
get_sl_average_prof
get_stokessl_lifoxkemper
la_frachbl
la_misalignment
lf17
mom_error_handler::mom_error
null_wavemethod
numbands
surfbands
testprof
wavemethod
- Called from
mom_energetic_pbl::epbl_column
mom_cvmix_kpp::kpp_compute_bld
update_stokes_drift
-
subroutine
mom_wave_interface/
get_stokessl_lifoxkemper
(ustar, hbl, GV, US, 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.
Computs 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
ustokes_sl :: [out] Surface layer averaged Stokes drift [m s-1]
la :: [out] Langmuir number
- Call to
- Called from
-
subroutine
mom_wave_interface/
get_sl_average_prof
(GV, AvgDepth, H, 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].
h :: [in] Grid thickness [H ~> m or kg m-2]
profile :: [in] Profile of quantity to be averaged [arbitrary]
average :: [out] Output quantity averaged over depth AvgDepth [arbitrary] (used here for Stokes drift)
- Called from
-
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 [m s-1]
average :: [out] Output average Stokes drift over depth AvgDepth [m s-1]
- Called from
-
subroutine
mom_wave_interface/
dhh85_mid
(GV, US, 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
zpt :: [in] Depth to get Stokes drift [Z ~> m].
ustokes :: [out] Stokes drift [m s-1]
- Call to
pi
waveage
waveagepeakfreq
wavewind
- Called from
-
subroutine
mom_wave_interface/
stokesmixing
(G, GV, dt, h, 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]
u :: [inout] Velocity i-component [m s-1]
v :: [inout] Velocity j-component [m s-1]
waves :: Surface wave related control structure.
-
subroutine
mom_wave_interface/
coriolisstokes
(G, GV, DT, h, u, v, WAVES, US)¶ 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 [s] CHECK IF PASSING RIGHT TIMESTEP
h :: [in] Layer thicknesses [H ~> m or kg m-2]
u :: [inout] Velocity i-component [m s-1]
v :: [inout] Velocity j-component [m s-1]
waves :: Surface wave related control structure.
us :: [in] A dimensional unit scaling type
-
subroutine
mom_wave_interface/
ust_2_u10_coare3p5
(USTair, U10, GV, US)¶ 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 [m s-1]
u10 :: [out] 10-m neutral wind speed [m s-1]
gv :: [in] vertical grid type
us :: [in] A dimensional unit scaling type
- Called from
-
subroutine
mom_wave_interface/
waves_end
(CS)¶ Clear pointers, deallocate memory.
- Parameters
cs :: Control structure