MOM6
mom_wave_interface Module Reference

Detailed Description

Interface for surface waves.

Author
Brandon Reichl, 2018.

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).

Data Types

type  wave_parameters_cs
 Container for all surface wave related parameters. More...
 

Functions/Subroutines

subroutine, public mom_wave_interface_init (time, G, GV, US, param_file, CS, diag)
 Initializes parameters related to MOM_wave_interface. More...
 
subroutine, public 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. More...
 
subroutine, public update_surface_waves (G, GV, US, Day, dt, CS)
 Subroutine that handles updating of surface wave/Stokes drift related properties. More...
 
subroutine, public update_stokes_drift (G, GV, US, CS, h, ustar)
 Constructs the Stokes Drift profile on the model grid based on desired coupling options. More...
 
subroutine 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. More...
 
subroutine, public 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. More...
 
subroutine get_stokessl_lifoxkemper (ustar, hbl, GV, US, UStokes_SL, LA)
 Get SL averaged Stokes drift from Li/FK 17 method. More...
 
subroutine get_sl_average_prof (GV, AvgDepth, H, Profile, Average)
 Get SL Averaged Stokes drift from a Stokes drift Profile. More...
 
subroutine get_sl_average_band (GV, AvgDepth, NB, WaveNumbers, SurfStokes, Average)
 Get SL averaged Stokes drift from the banded Spectrum method. More...
 
subroutine dhh85_mid (GV, US, zpt, UStokes)
 Compute the Stokes drift at a given depth. More...
 
subroutine, public stokesmixing (G, GV, dt, h, u, v, Waves)
 Explicit solver for Stokes mixing. Still in development do not use. More...
 
subroutine, public 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**. More...
 
subroutine 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. More...
 
subroutine, public waves_end (CS)
 Clear pointers, deallocate memory. More...
 

Variables

integer wavemethod =-99
 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. More...
 
integer, public numbands =0
 Number of wavenumber/frequency partitions to receive This needs to match the number of bands provided via either coupling or file. More...
 
integer, public partitionmode
 Method for partition mode (meant to check input) 0 - wavenumbers 1 - frequencies. More...
 
integer datasource
 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) More...
 
character(len=40) surfbandfilename
 Filename if using DataOverride. More...
 
logical dataoverrideisinitialized
 Flag for DataOverride Initialization. More...
 
real la_frachbl
 Fraction of OSBL for averaging Langmuir number. More...
 
logical la_misalignment = .false.
 Flag to use misalignment in Langmuir number. More...
 
character(len=40) mdl = "MOM_wave_interface"
 This module's name.
 
integer, parameter testprof = 0
 Undocumented parameters. More...
 
integer, parameter surfbands = 1
 Undocumented parameters. More...
 
integer, parameter dhh85 = 2
 Undocumented parameters. More...
 
integer, parameter lf17 = 3
 Undocumented parameters. More...
 
integer, parameter null_wavemethod =-99
 Undocumented parameters. More...
 
integer, parameter dataovr = 1
 Undocumented parameters. More...
 
integer, parameter coupler = 2
 Undocumented parameters. More...
 
integer, parameter input = 3
 Undocumented parameters. More...
 
real tp_stkx0
 Undocumented parameters. More...
 
real tp_stky0
 Undocumented parameters. More...
 
real tp_wvl
 Undocumented parameters. More...
 
logical waveagepeakfreq
 Undocumented parameters. More...
 
logical staticwaves
 Undocumented parameters. More...
 
logical dhh85_is_set
 Undocumented parameters. More...
 
real waveage
 Undocumented parameters. More...
 
real wavewind
 Undocumented parameters. More...
 
real pi
 Undocumented parameters. More...
 

Function/Subroutine Documentation

◆ coriolisstokes()

subroutine, public mom_wave_interface::coriolisstokes ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, intent(in)  DT,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
real, dimension(szib_(g),szj_(g),szk_(g)), intent(inout)  u,
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(inout)  v,
type(wave_parameters_cs), pointer  WAVES,
type(unit_scale_type), intent(in)  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
[in]gOcean grid
[in]gvOcean vertical grid
[in]dtTime step of MOM6 [s] CHECK IF PASSING RIGHT TIMESTEP
[in]hLayer thicknesses [H ~> m or kg m-2]
[in,out]uVelocity i-component [m s-1]
[in,out]vVelocity j-component [m s-1]
wavesSurface wave related control structure.
[in]usA dimensional unit scaling type

Definition at line 1288 of file MOM_wave_interface.F90.

1288  type(ocean_grid_type), &
1289  intent(in) :: g !< Ocean grid
1290  type(verticalgrid_type), &
1291  intent(in) :: gv !< Ocean vertical grid
1292  real, intent(in) :: dt !< Time step of MOM6 [s] CHECK IF PASSING RIGHT TIMESTEP
1293  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1294  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
1295  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1296  intent(inout) :: u !< Velocity i-component [m s-1]
1297  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1298  intent(inout) :: v !< Velocity j-component [m s-1]
1299  type(wave_parameters_cs), &
1300  pointer :: waves !< Surface wave related control structure.
1301  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1302  ! Local variables
1303  real :: dvel ! A rescaled velocity change [m s-1 T-1 ~> m s-2]
1304  integer :: i,j,k
1305 
1306  do k = 1, g%ke
1307  do j = g%jsc, g%jec
1308  do i = g%iscB, g%iecB
1309  dvel = 0.25*(waves%us_y(i,j+1,k)+waves%us_y(i-1,j+1,k))*g%CoriolisBu(i,j+1) + &
1310  0.25*(waves%us_y(i,j,k)+waves%us_y(i-1,j,k))*g%CoriolisBu(i,j)
1311  u(i,j,k) = u(i,j,k) + dvel*us%s_to_T*dt
1312  enddo
1313  enddo
1314  enddo
1315 
1316  do k = 1, g%ke
1317  do j = g%jscB, g%jecB
1318  do i = g%isc, g%iec
1319  dvel = 0.25*(waves%us_x(i+1,j,k)+waves%us_x(i+1,j-1,k))*g%CoriolisBu(i+1,j) + &
1320  0.25*(waves%us_x(i,j,k)+waves%us_x(i,j-1,k))*g%CoriolisBu(i,j)
1321  v(i,j,k) = v(i,j,k) - dvel*us%s_to_T*dt
1322  enddo
1323  enddo
1324  enddo

◆ dhh85_mid()

subroutine mom_wave_interface::dhh85_mid ( type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, intent(in)  zpt,
real, intent(out)  UStokes 
)
private

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
[in]gvOcean vertical grid
[in]usA dimensional unit scaling type
[in]zptDepth to get Stokes drift [Z ~> m].
[out]ustokesStokes drift [m s-1]

Definition at line 1163 of file MOM_wave_interface.F90.

1163  type(verticalgrid_type), intent(in) :: gv !< Ocean vertical grid
1164  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1165  real, intent(in) :: zpt !< Depth to get Stokes drift [Z ~> m].
1166  real, intent(out) :: ustokes !< Stokes drift [m s-1]
1167  !
1168  real :: ann, bnn, snn, cnn, dnn
1169  real :: omega_peak, omega, u10, wa, domega
1170  real :: omega_min, omega_max, wavespec, stokes
1171  real :: g_earth ! Gravitational acceleration [m s-2]
1172  integer :: nomega, oi
1173 
1174  wa = waveage
1175  u10 = wavewind
1176  g_earth = us%L_T_to_m_s**2*us%m_to_Z * gv%g_Earth
1177 
1178  !/
1179  omega_min = 0.1 ! Hz
1180  ! Cut off at 30cm for now...
1181  omega_max = 10. ! ~sqrt(0.2*g_Earth*2*pi/0.3)
1182  nomega = 1000
1183  domega = (omega_max-omega_min)/real(nomega)
1184 
1185  !
1186  if (waveagepeakfreq) then
1187  omega_peak = g_earth / (wa * u10)
1188  else
1189  omega_peak = 2. * pi * 0.13 * g_earth / u10
1190  endif
1191  !/
1192  ann = 0.006 * waveage**(-0.55)
1193  bnn = 1.0
1194  snn = 0.08 * (1.0 + 4.0 * waveage**3)
1195  cnn = 1.7
1196  if (wa < 1.) then
1197  cnn = cnn - 6.0*log10(wa)
1198  endif
1199  !/
1200  ustokes = 0.0
1201  omega = omega_min + 0.5*domega
1202  do oi = 1,nomega-1
1203  dnn = exp( -0.5 * (omega-omega_peak)**2 / (snn**2 * omega_peak**2) )
1204  ! wavespec units = m2s
1205  wavespec = (ann * g_earth**2 / (omega_peak*omega**4 ) ) * &
1206  exp(-bnn*(omega_peak/omega)**4)*cnn**dnn
1207  ! Stokes units m (multiply by frequency range for units of m/s)
1208  stokes = 2.0 * wavespec * omega**3 * &
1209  exp( 2.0 * omega**2 * us%Z_to_m*zpt / g_earth) / g_earth
1210  ustokes = ustokes + stokes*domega
1211  omega = omega + domega
1212  enddo
1213 
1214  return

◆ get_langmuir_number()

subroutine, public mom_wave_interface::get_langmuir_number ( real, intent(out)  LA,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, intent(in)  HBL,
real, intent(in)  ustar,
integer, intent(in)  i,
integer, intent(in)  j,
real, dimension( gv %ke), intent(in), optional  H,
real, dimension( gv %ke), intent(in), optional  U_H,
real, dimension( gv %ke), intent(in), optional  V_H,
logical, intent(in), optional  Override_MA,
type(wave_parameters_cs), pointer  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
[in]gOcean grid structure
[in]gvOcean vertical grid structure
[in]usA dimensional unit scaling type
[in]iMeridional index of h-point
[in]jZonal index of h-point
[in]ustarFriction velocity [Z T-1 ~> m s-1].
[in]hbl(Positive) thickness of boundary layer [Z ~> m].
[in]override_maOverride 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.
[in]hGrid layer thickness [H ~> m or kg m-2]
[in]u_hZonal velocity at H point [m s-1]
[in]v_hMeridional velocity at H point [m s-1]
wavesSurface wave control structure.
[out]laLangmuir number

Definition at line 880 of file MOM_wave_interface.F90.

880  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
881  type(verticalgrid_type), intent(in) :: gv !< Ocean vertical grid structure
882  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
883  integer, intent(in) :: i !< Meridional index of h-point
884  integer, intent(in) :: j !< Zonal index of h-point
885  real, intent(in) :: ustar !< Friction velocity [Z T-1 ~> m s-1].
886  real, intent(in) :: hbl !< (Positive) thickness of boundary layer [Z ~> m].
887  logical, optional, intent(in) :: override_ma !< Override to use misalignment in LA
888  !! calculation. This can be used if diagnostic
889  !! LA outputs are desired that are different than
890  !! those used by the dynamical model.
891  real, dimension(SZK_(GV)), optional, &
892  intent(in) :: h !< Grid layer thickness [H ~> m or kg m-2]
893  real, dimension(SZK_(GV)), optional, &
894  intent(in) :: u_h !< Zonal velocity at H point [m s-1]
895  real, dimension(SZK_(GV)), optional, &
896  intent(in) :: v_h !< Meridional velocity at H point [m s-1]
897  type(wave_parameters_cs), &
898  pointer :: waves !< Surface wave control structure.
899 
900  real, intent(out) :: la !< Langmuir number
901 
902 !Local Variables
903  real :: top, bottom, midpoint
904  real :: dpt_lasl, sheardirection, wavedirection
905  real :: la_stkx, la_stky, la_stk ! Stokes velocities in [m s-1]
906  logical :: continueloop, use_ma
907  real, dimension(SZK_(G)) :: us_h, vs_h
908  real, dimension(NumBands) :: stkband_x, stkband_y
909  integer :: kk, bb
910 
911  ! Compute averaging depth for Stokes drift (negative)
912  dpt_lasl = min(-0.1*us%m_to_Z, -la_frachbl*hbl)
913 
914  use_ma = la_misalignment
915  if (present(override_ma)) use_ma = override_ma
916 
917  ! If requesting to use misalignment in the Langmuir number compute the Shear Direction
918  if (use_ma) then
919  if (.not.(present(h).and.present(u_h).and.present(v_h))) then
920  call mom_error(fatal,'Get_LA_waves requested to consider misalignment.')
921  endif
922  continueloop = .true.
923  bottom = 0.0
924  do kk = 1,g%ke
925  top = bottom
926  midpoint = bottom + gv%H_to_Z*0.5*h(kk)
927  bottom = bottom + gv%H_to_Z*h(kk)
928  if (midpoint > dpt_lasl .and. kk > 1 .and. continueloop) then
929  sheardirection = atan2(v_h(1)-v_h(kk),u_h(1)-u_h(kk))
930  continueloop = .false.
931  endif
932  enddo
933  endif
934 
935  if (wavemethod==testprof) then
936  do kk = 1,g%ke
937  us_h(kk) = 0.5*(waves%US_X(i,j,kk)+waves%US_X(i-1,j,kk))
938  vs_h(kk) = 0.5*(waves%US_Y(i,j,kk)+waves%US_Y(i,j-1,kk))
939  enddo
940  call get_sl_average_prof( gv, dpt_lasl, h, us_h, la_stkx)
941  call get_sl_average_prof( gv, dpt_lasl, h, vs_h, la_stky)
942  la_stk = sqrt(la_stkx*la_stkx+la_stky*la_stky)
943  elseif (wavemethod==surfbands) then
944  do bb = 1,numbands
945  stkband_x(bb) = 0.5*(waves%STKx0(i,j,bb)+waves%STKx0(i-1,j,bb))
946  stkband_y(bb) = 0.5*(waves%STKy0(i,j,bb)+waves%STKy0(i,j-1,bb))
947  enddo
948  call get_sl_average_band(gv, dpt_lasl, numbands, waves%WaveNum_Cen, stkband_x, la_stkx )
949  call get_sl_average_band(gv, dpt_lasl, numbands, waves%WaveNum_Cen, stkband_y, la_stky )
950  la_stk = sqrt(la_stkx**2 + la_stky**2)
951  elseif (wavemethod==dhh85) then
952  ! Temporarily integrating profile rather than spectrum for simplicity
953  do kk = 1,gv%ke
954  us_h(kk) = 0.5*(waves%US_X(i,j,kk)+waves%US_X(i-1,j,kk))
955  vs_h(kk) = 0.5*(waves%US_Y(i,j,kk)+waves%US_Y(i,j-1,kk))
956  enddo
957  call get_sl_average_prof( gv, dpt_lasl, h, us_h, la_stkx)
958  call get_sl_average_prof( gv, dpt_lasl, h, vs_h, la_stky)
959  la_stk = sqrt(la_stkx**2 + la_stky**2)
960  elseif (wavemethod==lf17) then
961  call get_stokessl_lifoxkemper(ustar, hbl*la_frachbl, gv, us, la_stk, la)
962  elseif (wavemethod==null_wavemethod) then
963  call mom_error(fatal, "Get_Langmuir_number called without defining a WaveMethod. "//&
964  "Suggest to make sure USE_LT is set/overridden to False or "//&
965  "choose a wave method (or set USE_LA_LI2016 to use statistical "//&
966  "waves.")
967  endif
968 
969  if (.not.(wavemethod==lf17)) then
970  ! This is an arbitrary lower bound on Langmuir number.
971  ! We shouldn't expect values lower than this, but
972  ! there is also no good reason to cap it here other then
973  ! to prevent large enhancements in unconstrained parts of
974  ! the curve fit parameterizations.
975  ! Note the dimensional constant background Stokes velocity of 10^-10 m s-1.
976  la = max(waves%La_min, sqrt(us%Z_to_m*us%s_to_T*ustar / (la_stk+1.e-10)))
977  endif
978 
979  if (use_ma) then
980  wavedirection = atan2(la_stky, la_stkx)
981  la = la / sqrt(max(1.e-8, cos( wavedirection - sheardirection)))
982  endif
983 
984  return

◆ get_sl_average_band()

subroutine mom_wave_interface::get_sl_average_band ( type(verticalgrid_type), intent(in)  GV,
real, intent(in)  AvgDepth,
integer, intent(in)  NB,
real, dimension(nb), intent(in)  WaveNumbers,
real, dimension(nb), intent(in)  SurfStokes,
real, intent(out)  Average 
)
private

Get SL averaged Stokes drift from the banded Spectrum method.

Parameters
[in]gvOcean vertical grid
[in]avgdepthDepth to average over [Z ~> m].
[in]nbNumber of bands used
[in]wavenumbersWavenumber corresponding to each band [Z-1 ~> m-1]
[in]surfstokesSurface Stokes drift for each band [m s-1]
[out]averageOutput average Stokes drift over depth AvgDepth [m s-1]

Definition at line 1130 of file MOM_wave_interface.F90.

1130  type(verticalgrid_type), &
1131  intent(in) :: gv !< Ocean vertical grid
1132  real, intent(in) :: avgdepth !< Depth to average over [Z ~> m].
1133  integer, intent(in) :: nb !< Number of bands used
1134  real, dimension(NB), &
1135  intent(in) :: wavenumbers !< Wavenumber corresponding to each band [Z-1 ~> m-1]
1136  real, dimension(NB), &
1137  intent(in) :: surfstokes !< Surface Stokes drift for each band [m s-1]
1138  real, intent(out) :: average !< Output average Stokes drift over depth AvgDepth [m s-1]
1139 
1140  ! Local variables
1141  integer :: bb
1142 
1143  ! Loop over bands
1144  average = 0.0
1145  do bb = 1, nb
1146  ! Factor includes analytical integration of e(2kz)
1147  ! - divided by (-H_LA) to get average from integral.
1148  average = average + surfstokes(bb) * &
1149  (1.-exp(-abs(avgdepth * 2.0 * wavenumbers(bb)))) / &
1150  abs(avgdepth * 2.0 * wavenumbers(bb))
1151  enddo
1152 
1153  return

◆ get_sl_average_prof()

subroutine mom_wave_interface::get_sl_average_prof ( type(verticalgrid_type), intent(in)  GV,
real, intent(in)  AvgDepth,
real, dimension(szk_(gv)), intent(in)  H,
real, dimension(szk_(gv)), intent(in)  Profile,
real, intent(out)  Average 
)
private

Get SL Averaged Stokes drift from a Stokes drift Profile.

Parameters
[in]gvOcean vertical grid structure
[in]avgdepthDepth to average over (negative) [Z ~> m].
[in]hGrid thickness [H ~> m or kg m-2]
[in]profileProfile of quantity to be averaged [arbitrary]
[out]averageOutput quantity averaged over depth AvgDepth [arbitrary] (used here for Stokes drift)

Definition at line 1083 of file MOM_wave_interface.F90.

1083  type(verticalgrid_type), &
1084  intent(in) :: gv !< Ocean vertical grid structure
1085  real, intent(in) :: avgdepth !< Depth to average over (negative) [Z ~> m].
1086  real, dimension(SZK_(GV)), &
1087  intent(in) :: h !< Grid thickness [H ~> m or kg m-2]
1088  real, dimension(SZK_(GV)), &
1089  intent(in) :: profile !< Profile of quantity to be averaged [arbitrary]
1090  !! (used here for Stokes drift)
1091  real, intent(out) :: average !< Output quantity averaged over depth AvgDepth [arbitrary]
1092  !! (used here for Stokes drift)
1093  !Local variables
1094  real :: top, midpoint, bottom ! Depths, negative downward [Z ~> m].
1095  real :: sum
1096  integer :: kk
1097 
1098  ! Initializing sum
1099  sum = 0.0
1100 
1101  ! Integrate
1102  bottom = 0.0
1103  do kk = 1, gv%ke
1104  top = bottom
1105  midpoint = bottom - gv%H_to_Z * 0.5*h(kk)
1106  bottom = bottom - gv%H_to_Z * h(kk)
1107  if (avgdepth < bottom) then ! The whole cell is within H_LA
1108  sum = sum + profile(kk) * (gv%H_to_Z * h(kk))
1109  elseif (avgdepth < top) then ! A partial cell is within H_LA
1110  sum = sum + profile(kk) * (top-avgdepth)
1111  exit
1112  else
1113  exit
1114  endif
1115  enddo
1116 
1117  ! Divide by AvgDepth or the depth in the column, whichever is smaller.
1118  if (abs(avgdepth) <= abs(bottom)) then
1119  average = sum / abs(avgdepth)
1120  elseif (abs(bottom) > 0.0) then
1121  average = sum / abs(bottom)
1122  else
1123  average = 0.0
1124  endif
1125 

◆ get_stokessl_lifoxkemper()

subroutine mom_wave_interface::get_stokessl_lifoxkemper ( real, intent(in)  ustar,
real, intent(in)  hbl,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, intent(out)  UStokes_SL,
real, intent(out)  LA 
)
private

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
[in]ustarwater-side surface friction velocity [Z T-1 ~> m s-1].
[in]hblboundary layer depth [Z ~> m].
[in]gvOcean vertical grid structure
[in]usA dimensional unit scaling type
[out]ustokes_slSurface layer averaged Stokes drift [m s-1]
[out]laLangmuir number

Definition at line 1004 of file MOM_wave_interface.F90.

1004  real, intent(in) :: ustar !< water-side surface friction velocity [Z T-1 ~> m s-1].
1005  real, intent(in) :: hbl !< boundary layer depth [Z ~> m].
1006  type(verticalgrid_type), intent(in) :: gv !< Ocean vertical grid structure
1007  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1008  real, intent(out) :: ustokes_sl !< Surface layer averaged Stokes drift [m s-1]
1009  real, intent(out) :: la !< Langmuir number
1010  ! Local variables
1011  ! parameters
1012  real, parameter :: &
1013  ! ratio of U19.5 to U10 (Holthuijsen, 2007)
1014  u19p5_to_u10 = 1.075, &
1015  ! ratio of mean frequency to peak frequency for
1016  ! Pierson-Moskowitz spectrum (Webb, 2011)
1017  fm_into_fp = 1.296, &
1018  ! ratio of surface Stokes drift to U10
1019  us_to_u10 = 0.0162, &
1020  ! loss ratio of Stokes transport
1021  r_loss = 0.667
1022  real :: ustokes, hm0, fm, fp, vstokes, kphil, kstar
1023  real :: z0, z0i, r1, r2, r3, r4, tmp, lasl_sqr_i
1024  real :: u10
1025 
1026  if (ustar > 0.0) then
1027  ! Computing u10 based on u_star and COARE 3.5 relationships
1028  call ust_2_u10_coare3p5(us%Z_to_m*us%s_to_T*ustar*sqrt(us%R_to_kg_m3*gv%Rho0/1.225), u10, gv, us)
1029  ! surface Stokes drift
1030  ustokes = us_to_u10*u10
1031  !
1032  ! significant wave height from Pierson-Moskowitz
1033  ! spectrum (Bouws, 1998)
1034  hm0 = 0.0246 *u10**2
1035  !
1036  ! peak frequency (PM, Bouws, 1998)
1037  tmp = 2.0 * pi * u19p5_to_u10 * u10
1038  fp = 0.877 * us%L_T_to_m_s**2*us%m_to_Z * gv%g_Earth / tmp
1039  !
1040  ! mean frequency
1041  fm = fm_into_fp * fp
1042  !
1043  ! total Stokes transport (a factor r_loss is applied to account
1044  ! for the effect of directional spreading, multidirectional waves
1045  ! and the use of PM peak frequency and PM significant wave height
1046  ! on estimating the Stokes transport)
1047  vstokes = 0.125 * pi * r_loss * fm * hm0**2
1048  !
1049  ! the general peak wavenumber for Phillips' spectrum
1050  ! (Breivik et al., 2016) with correction of directional spreading
1051  kphil = 0.176 * ustokes / vstokes
1052  !
1053  ! surface layer averaged Stokes dirft with Stokes drift profile
1054  ! estimated from Phillips' spectrum (Breivik et al., 2016)
1055  ! the directional spreading effect from Webb and Fox-Kemper, 2015
1056  ! is also included
1057  kstar = kphil * 2.56
1058  ! surface layer
1059  z0 = abs(us%Z_to_m*hbl)
1060  z0i = 1.0 / z0
1061  ! term 1 to 4
1062  r1 = ( 0.151 / kphil * z0i -0.84 ) * &
1063  ( 1.0 - exp(-2.0 * kphil * z0) )
1064  r2 = -( 0.84 + 0.0591 / kphil * z0i ) * &
1065  sqrt( 2.0 * pi * kphil * z0 ) * &
1066  erfc( sqrt( 2.0 * kphil * z0 ) )
1067  r3 = ( 0.0632 / kstar * z0i + 0.125 ) * &
1068  (1.0 - exp(-2.0 * kstar * z0) )
1069  r4 = ( 0.125 + 0.0946 / kstar * z0i ) * &
1070  sqrt( 2.0 * pi *kstar * z0) * &
1071  erfc( sqrt( 2.0 * kstar * z0 ) )
1072  ustokes_sl = ustokes * (0.715 + r1 + r2 + r3 + r4)
1073  la = sqrt(us%Z_to_m*us%s_to_T*ustar / ustokes_sl)
1074  else
1075  ustokes_sl = 0.0
1076  la=1.e8
1077  endif
1078 

◆ mom_wave_interface_init()

subroutine, public mom_wave_interface::mom_wave_interface_init ( type(time_type), intent(in), target  time,
type(ocean_grid_type), intent(inout)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
type(wave_parameters_cs), pointer  CS,
type(diag_ctrl), intent(inout), target  diag 
)

Initializes parameters related to MOM_wave_interface.

Parameters
[in]timeModel time
[in,out]gGrid structure
[in]gvVertical grid structure
[in]usA dimensional unit scaling type
[in]param_fileInput parameter structure
csWave parameter control structure
[in,out]diagDiagnostic Pointer

Definition at line 193 of file MOM_wave_interface.F90.

193  type(time_type), target, intent(in) :: time !< Model time
194  type(ocean_grid_type), intent(inout) :: g !< Grid structure
195  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
196  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
197  type(param_file_type), intent(in) :: param_file !< Input parameter structure
198  type(wave_parameters_cs), pointer :: cs !< Wave parameter control structure
199  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostic Pointer
200  ! Local variables
201  ! I/O
202  character*(13) :: tmpstring1,tmpstring2
203  character*(5), parameter :: null_string = "EMPTY"
204  character*(12), parameter :: testprof_string = "TEST_PROFILE"
205  character*(13), parameter :: surfbands_string = "SURFACE_BANDS"
206  character*(5), parameter :: dhh85_string = "DHH85"
207  character*(4), parameter :: lf17_string = "LF17"
208  character*(12), parameter :: dataovr_string = "DATAOVERRIDE"
209  character*(7), parameter :: coupler_string = "COUPLER"
210  character*(5), parameter :: input_string = "INPUT"
211 
212  ! Dummy Check
213  if (associated(cs)) then
214  call mom_error(fatal, "wave_interface_init called with an associated"//&
215  "control structure.")
216  return
217  endif
218 
219  pi=4.0*atan(1.0)
220 
221  ! Allocate CS and set pointers
222  allocate(cs)
223 
224  cs%diag => diag
225  cs%Time => time
226 
227  ! Add any initializations needed here
228  dataoverrideisinitialized = .false.
229 
230  ! The only way to get here is with UseWaves enabled.
231  cs%UseWaves=.true.
232 
233  call log_version(param_file, mdl, version)
234 
235  ! Wave modified physics
236  ! Presently these are all in research mode
237  call get_param(param_file, mdl, "LAGRANGIAN_MIXING", cs%LagrangianMixing, &
238  "Flag to use Lagrangian Mixing of momentum", units="", &
239  default=.false.)
240  if (cs%LagrangianMixing) then
241  ! Force Code Intervention
242  call mom_error(fatal,"Should you be enabling Lagrangian Mixing? Code not ready.")
243  endif
244  call get_param(param_file, mdl, "STOKES_MIXING", cs%StokesMixing, &
245  "Flag to use Stokes Mixing of momentum", units="", &
246  default=.false.)
247  if (cs%StokesMixing) then
248  ! Force Code Intervention
249  call mom_error(fatal,"Should you be enabling Stokes Mixing? Code not ready.")
250  endif
251  call get_param(param_file, mdl, "CORIOLIS_STOKES", cs%CoriolisStokes, &
252  "Flag to use Coriolis Stokes acceleration", units="", &
253  default=.false.)
254  if (cs%CoriolisStokes) then
255  ! Force Code Intervention
256  call mom_error(fatal,"Should you be enabling Coriolis-Stokes? Code not ready.")
257  endif
258 
259  ! Get Wave Method and write to integer WaveMethod
260  call get_param(param_file,mdl,"WAVE_METHOD",tmpstring1, &
261  "Choice of wave method, valid options include: \n"// &
262  " TEST_PROFILE - Prescribed from surface Stokes drift \n"// &
263  " and a decay wavelength.\n"// &
264  " SURFACE_BANDS - Computed from multiple surface values \n"// &
265  " and decay wavelengths.\n"// &
266  " DHH85 - Uses Donelan et al. 1985 empirical \n"// &
267  " wave spectrum with prescribed values. \n"// &
268  " LF17 - Infers Stokes drift profile from wind \n"// &
269  " speed following Li and Fox-Kemper 2017.\n", &
270  units='', default=null_string)
271  select case (trim(tmpstring1))
272  case (null_string)! No Waves
273  call mom_error(fatal, "wave_interface_init called with no specified "//&
274  "WAVE_METHOD.")
275  case (testprof_string)! Test Profile
276  wavemethod = testprof
277  call get_param(param_file,mdl,"TP_STKX_SURF",tp_stkx0,&
278  'Surface Stokes (x) for test profile',&
279  units='m/s',default=0.1)
280  call get_param(param_file,mdl,"TP_STKY_SURF",tp_stky0,&
281  'Surface Stokes (y) for test profile',&
282  units='m/s',default=0.0)
283  call get_param(param_file,mdl,"TP_WVL",tp_wvl,&
284  units='m', default=50.0, scale=us%m_to_Z)
285  case (surfbands_string)! Surface Stokes Drift Bands
286  wavemethod = surfbands
287  call get_param(param_file, mdl, "SURFBAND_SOURCE",tmpstring2, &
288  "Choice of SURFACE_BANDS data mode, valid options include: \n"// &
289  " DATAOVERRIDE - Read from NetCDF using FMS DataOverride. \n"// &
290  " COUPLER - Look for variables from coupler pass \n"// &
291  " INPUT - Testing with fixed values.", &
292  units='', default=null_string)
293  select case (trim(tmpstring2))
294  case (null_string)! Default
295  call mom_error(fatal, "wave_interface_init called with SURFACE_BANDS"//&
296  " but no SURFBAND_SOURCE.")
297  case (dataovr_string)! Using Data Override
298  datasource = dataovr
299  call get_param(param_file, mdl, "SURFBAND_FILENAME", surfbandfilename, &
300  "Filename of surface Stokes drift input band data.", default="StkSpec.nc")
301  case (coupler_string)! Reserved for coupling
302  datasource = coupler
303  case (input_string)! A method to input the Stokes band (globally uniform)
304  datasource = input
305  call get_param(param_file,mdl,"SURFBAND_NB",numbands, &
306  "Prescribe number of wavenumber bands for Stokes drift. "// &
307  "Make sure this is consistnet w/ WAVENUMBERS, STOKES_X, and "// &
308  "STOKES_Y, there are no safety checks in the code.", &
309  units='', default=1)
310  allocate( cs%WaveNum_Cen(1:numbands) )
311  cs%WaveNum_Cen(:) = 0.0
312  allocate( cs%PrescribedSurfStkX(1:numbands))
313  cs%PrescribedSurfStkX(:) = 0.0
314  allocate( cs%PrescribedSurfStkY(1:numbands))
315  cs%PrescribedSurfStkY(:) = 0.0
316  allocate( cs%STKx0(g%isdB:g%iedB,g%jsd:g%jed,1:numbands))
317  cs%STKx0(:,:,:) = 0.0
318  allocate( cs%STKy0(g%isd:g%ied,g%jsdB:g%jedB,1:numbands))
319  cs%STKy0(:,:,:) = 0.0
320  partitionmode=0
321  call get_param(param_file,mdl,"SURFBAND_WAVENUMBERS",cs%WaveNum_Cen, &
322  "Central wavenumbers for surface Stokes drift bands.",units='rad/m', &
323  default=0.12566)
324  call get_param(param_file,mdl,"SURFBAND_STOKES_X",cs%PrescribedSurfStkX, &
325  "X-direction surface Stokes drift for bands.",units='m/s', &
326  default=0.15)
327  call get_param(param_file,mdl,"SURFBAND_STOKES_Y",cs%PrescribedSurfStkY, &
328  "Y-direction surface Stokes drift for bands.",units='m/s', &
329  default=0.0)
330  case default! No method provided
331  call mom_error(fatal,'Check WAVE_METHOD.')
332  end select
333 
334  case (dhh85_string)!Donelan et al., 1985 spectrum
335  wavemethod = dhh85
336  call mom_error(warning,"DHH85 only ever set-up for uniform cases w/"//&
337  " Stokes drift in x-direction.")
338  call get_param(param_file,mdl,"DHH85_AGE_FP",waveagepeakfreq, &
339  "Choose true to use waveage in peak frequency.", &
340  units='', default=.false.)
341  call get_param(param_file,mdl,"DHH85_AGE",waveage, &
342  "Wave Age for DHH85 spectrum.", &
343  units='', default=1.2)
344  call get_param(param_file,mdl,"DHH85_WIND",wavewind, &
345  "Wind speed for DHH85 spectrum.", &
346  units='', default=10.0)
347  call get_param(param_file,mdl,"STATIC_DHH85",staticwaves, &
348  "Flag to disable updating DHH85 Stokes drift.", &
349  default=.false.)
350  case (lf17_string)!Li and Fox-Kemper 17 wind-sea Langmuir number
351  wavemethod = lf17
352  case default
353  call mom_error(fatal,'Check WAVE_METHOD.')
354  end select
355 
356  ! Langmuir number Options
357  call get_param(param_file, mdl, "LA_DEPTH_RATIO", la_frachbl, &
358  "The depth (normalized by BLD) to average Stokes drift over in "//&
359  "Langmuir number calculation, where La = sqrt(ust/Stokes).", &
360  units="nondim",default=0.04)
361  call get_param(param_file, mdl, "LA_MISALIGNMENT", la_misalignment, &
362  "Flag (logical) if using misalignment bt shear and waves in LA",&
363  default=.false.)
364  call get_param(param_file, mdl, "MIN_LANGMUIR", cs%La_min, &
365  "A minimum value for all Langmuir numbers that is not physical, "//&
366  "but is likely only encountered when the wind is very small and "//&
367  "therefore its effects should be mostly benign.",units="nondim",&
368  default=0.05)
369 
370  ! Allocate and initialize
371  ! a. Stokes driftProfiles
372  allocate(cs%Us_x(g%isdB:g%IedB,g%jsd:g%jed,g%ke))
373  cs%Us_x(:,:,:) = 0.0
374  allocate(cs%Us_y(g%isd:g%Ied,g%jsdB:g%jedB,g%ke))
375  cs%Us_y(:,:,:) = 0.0
376  ! b. Surface Values
377  allocate(cs%US0_x(g%isdB:g%iedB,g%jsd:g%jed))
378  cs%US0_x(:,:) = 0.0
379  allocate(cs%US0_y(g%isd:g%ied,g%jsdB:g%jedB))
380  cs%US0_y(:,:) = 0.0
381  ! c. Langmuir number
382  allocate(cs%La_SL(g%isc:g%iec,g%jsc:g%jec))
383  allocate(cs%La_turb(g%isc:g%iec,g%jsc:g%jec))
384  cs%La_SL(:,:) = 0.0
385  cs%La_turb (:,:) = 0.0
386  ! d. Viscosity for Stokes drift
387  if (cs%StokesMixing) then
388  allocate(cs%KvS(g%isd:g%Ied,g%jsd:g%jed,g%ke))
389  cs%KvS(:,:,:) = 0.0
390  endif
391 
392  ! Initialize Wave related outputs
393  cs%id_surfacestokes_y = register_diag_field('ocean_model','surface_stokes_y', &
394  cs%diag%axesCu1,time,'Surface Stokes drift (y)','m s-1')
395  cs%id_surfacestokes_x = register_diag_field('ocean_model','surface_stokes_x', &
396  cs%diag%axesCv1,time,'Surface Stokes drift (x)','m s-1')
397  cs%id_3dstokes_y = register_diag_field('ocean_model','3d_stokes_y', &
398  cs%diag%axesCvL,time,'3d Stokes drift (y)','m s-1')
399  cs%id_3dstokes_x = register_diag_field('ocean_model','3d_stokes_x', &
400  cs%diag%axesCuL,time,'3d Stokes drift (y)','m s-1')
401  cs%id_La_turb = register_diag_field('ocean_model','La_turbulent',&
402  cs%diag%axesT1,time,'Surface (turbulent) Langmuir number','nondim')
403 
404  return

◆ mom_wave_interface_init_lite()

subroutine, public mom_wave_interface::mom_wave_interface_init_lite ( type(param_file_type), intent(in)  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
[in]param_fileInput parameter structure

Definition at line 410 of file MOM_wave_interface.F90.

410  type(param_file_type), intent(in) :: param_file !< Input parameter structure
411  character*(5), parameter :: null_string = "EMPTY"
412  character*(4), parameter :: lf17_string = "LF17"
413  character*(13) :: tmpstring1
414  logical :: statisticalwaves
415 
416  ! Langmuir number Options
417  call get_param(param_file, mdl, "LA_DEPTH_RATIO", la_frachbl, &
418  "The depth (normalized by BLD) to average Stokes drift over in "//&
419  "Langmuir number calculation, where La = sqrt(ust/Stokes).", &
420  units="nondim",default=0.04)
421 
422  ! Check if using LA_LI2016
423  call get_param(param_file,mdl,"USE_LA_LI2016",statisticalwaves, &
424  do_not_log=.true.,default=.false.)
425  if (statisticalwaves) then
426  wavemethod = lf17
427  pi=4.0*atan(1.0)
428  else
429  wavemethod = null_wavemethod
430  end if
431 
432  return

◆ stokesmixing()

subroutine, public mom_wave_interface::stokesmixing ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, intent(in)  dt,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
real, dimension(szib_(g),szj_(g),szk_(g)), intent(inout)  u,
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(inout)  v,
type(wave_parameters_cs), pointer  Waves 
)

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

Parameters
[in]gOcean grid
[in]gvOcean vertical grid
[in]dtTime step of MOM6 [T ~> s] for explicit solver
[in]hLayer thicknesses [H ~> m or kg m-2]
[in,out]uVelocity i-component [m s-1]
[in,out]vVelocity j-component [m s-1]
wavesSurface wave related control structure.

Definition at line 1220 of file MOM_wave_interface.F90.

1220  type(ocean_grid_type), &
1221  intent(in) :: g !< Ocean grid
1222  type(verticalgrid_type), &
1223  intent(in) :: gv !< Ocean vertical grid
1224  real, intent(in) :: dt !< Time step of MOM6 [T ~> s] for explicit solver
1225  real, dimension(SZI_(G),SZJ_(G),SZK_(G)),&
1226  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
1227  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1228  intent(inout) :: u !< Velocity i-component [m s-1]
1229  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1230  intent(inout) :: v !< Velocity j-component [m s-1]
1231  type(wave_parameters_cs), &
1232  pointer :: waves !< Surface wave related control structure.
1233  ! Local variables
1234  real :: dtauup, dtaudn ! Vertical momentum fluxes [Z T-1 m s-1]
1235  real :: h_lay ! The layer thickness at a velocity point [Z ~> m].
1236  integer :: i,j,k
1237 
1238 ! This is a template to think about down-Stokes mixing.
1239 ! This is not ready for use...
1240 
1241  do k = 1, g%ke
1242  do j = g%jsc, g%jec
1243  do i = g%iscB, g%iecB
1244  h_lay = gv%H_to_Z*0.5*(h(i,j,k)+h(i+1,j,k))
1245  dtauup = 0.0
1246  if (k > 1) &
1247  dtauup = 0.5*(waves%Kvs(i,j,k)+waves%Kvs(i+1,j,k)) * &
1248  (waves%us_x(i,j,k-1)-waves%us_x(i,j,k)) / &
1249  (0.5*(h_lay + gv%H_to_Z*0.5*(h(i,j,k-1)+h(i+1,j,k-1)) ))
1250  dtaudn = 0.0
1251  if (k < g%ke-1) &
1252  dtaudn = 0.5*(waves%Kvs(i,j,k+1)+waves%Kvs(i+1,j,k+1)) * &
1253  (waves%us_x(i,j,k)-waves%us_x(i,j,k+1)) / &
1254  (0.5*(h_lay + gv%H_to_Z*0.5*(h(i,j,k+1)+h(i+1,j,k+1)) ))
1255  u(i,j,k) = u(i,j,k) + dt * (dtauup-dtaudn) / h_lay
1256  enddo
1257  enddo
1258  enddo
1259 
1260  do k = 1, g%ke
1261  do j = g%jscB, g%jecB
1262  do i = g%isc, g%iec
1263  h_lay = gv%H_to_Z*0.5*(h(i,j,k)+h(i,j+1,k))
1264  dtauup = 0.
1265  if (k > 1) &
1266  dtauup = 0.5*(waves%Kvs(i,j,k)+waves%Kvs(i,j+1,k)) * &
1267  (waves%us_y(i,j,k-1)-waves%us_y(i,j,k)) / &
1268  (0.5*(h_lay + gv%H_to_Z*0.5*(h(i,j,k-1)+h(i,j+1,k-1)) ))
1269  dtaudn = 0.0
1270  if (k < g%ke-1) &
1271  dtaudn =0.5*(waves%Kvs(i,j,k+1)+waves%Kvs(i,j+1,k+1)) * &
1272  (waves%us_y(i,j,k)-waves%us_y(i,j,k+1)) / &
1273  (0.5*(h_lay + gv%H_to_Z*0.5*(h(i,j,k+1)+h(i,j+1,k+1)) ))
1274  v(i,j,k) = v(i,j,k) + dt * (dtauup-dtaudn) / h_lay
1275  enddo
1276  enddo
1277  enddo
1278 

◆ surface_bands_by_data_override()

subroutine mom_wave_interface::surface_bands_by_data_override ( type(time_type), intent(in)  day_center,
type(ocean_grid_type), intent(inout)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(wave_parameters_cs), pointer  CS 
)
private

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

Parameters
[in]day_centerCenter of timestep
csWave structure
[in,out]gGrid structure
[in]gvVertical grid structure
[in]usA dimensional unit scaling type

Definition at line 710 of file MOM_wave_interface.F90.

710  use netcdf
711  type(time_type), intent(in) :: day_center !< Center of timestep
712  type(wave_parameters_cs), pointer :: cs !< Wave structure
713  type(ocean_grid_type), intent(inout) :: g !< Grid structure
714  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
715  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
716  ! Local variables
717  real :: temp_x(szi_(g),szj_(g)) ! Pseudo-zonal Stokes drift of band at h-points [m s-1]
718  real :: temp_y(szi_(g),szj_(g)) ! Psuedo-meridional Stokes drift of band at h-points [m s-1]
719  real :: top, midpoint
720  integer :: b
721  integer :: i, j
722  integer, dimension(4) :: start, counter, dims, dim_id
723  character(len=12) :: dim_name(4)
724  character(20) :: varname, varread1, varread2
725  integer :: rcode_fr, rcode_wn, ncid, varid_fr, varid_wn, id, ndims
726 
727  if (.not.dataoverrideisinitialized) then
728  call data_override_init(ocean_domain_in=g%Domain%mpp_domain)
729  dataoverrideisinitialized = .true.
730 
731  ! Read in number of wavenumber bands in file to set number to be read in
732  ! Hardcoded filename/variables
733  varread1 = 'wavenumber' !Old method gives wavenumber
734  varread2 = 'frequency' !New method gives frequency
735  rcode_wn = nf90_open(trim(surfbandfilename), nf90_nowrite, ncid)
736  if (rcode_wn /= 0) then
737  call mom_error(fatal,"error opening file "//trim(surfbandfilename)//&
738  " in MOM_wave_interface.")
739  endif
740 
741  ! Check if rcode_wn or rcode_fr is 0 (checks if input has wavenumber or frequency)
742  rcode_wn = nf90_inq_varid(ncid, varread1, varid_wn)
743  rcode_fr = nf90_inq_varid(ncid, varread2, varid_fr)
744 
745  if (rcode_wn /= 0 .and. rcode_fr /= 0) then
746  call mom_error(fatal,"error finding variable "//trim(varread1)//&
747  " or "//trim(varread2)//" in file "//trim(surfbandfilename)//" in MOM_wave_interface.")
748 
749  elseif (rcode_wn == 0) then
750  ! wavenumbers found:
751  partitionmode = 0
752  rcode_wn = nf90_inquire_variable(ncid, varid_wn, ndims=ndims, &
753  dimids=dims)
754  if (rcode_wn /= 0) then
755  call mom_error(fatal, &
756  'error inquiring dimensions MOM_wave_interface.')
757  endif
758  rcode_wn = nf90_inquire_dimension(ncid, dims(1), dim_name(1), len=id)
759  if (rcode_wn /= 0) then
760  call mom_error(fatal,"error reading dimension 1 data for "// &
761  trim(varread1)//" in file "// trim(surfbandfilename)// &
762  " in MOM_wave_interface.")
763  endif
764  rcode_wn = nf90_inq_varid(ncid, dim_name(1), dim_id(1))
765  if (rcode_wn /= 0) then
766  call mom_error(fatal,"error finding variable "//trim(dim_name(1))//&
767  " in file "//trim(surfbandfilename)//" in MOM_wave_interace.")
768  endif
769  ! Allocating size of wavenumber bins
770  allocate( cs%WaveNum_Cen(1:id) )
771  cs%WaveNum_Cen(:) = 0.0
772  elseif (rcode_fr == 0) then
773  ! frequencies found:
774  partitionmode = 1
775  rcode_fr = nf90_inquire_variable(ncid, varid_fr, ndims=ndims, &
776  dimids=dims)
777  if (rcode_fr /= 0) then
778  call mom_error(fatal,&
779  'error inquiring dimensions MOM_wave_interface.')
780  endif
781  rcode_fr = nf90_inquire_dimension(ncid, dims(1), dim_name(1), len=id)
782  if (rcode_fr /= 0) then
783  call mom_error(fatal,"error reading dimension 1 data for "// &
784  trim(varread2)//" in file "// trim(surfbandfilename)// &
785  " in MOM_wave_interface.")
786  endif
787  rcode_fr = nf90_inq_varid(ncid, dim_name(1), dim_id(1))
788  if (rcode_fr /= 0) then
789  call mom_error(fatal,"error finding variable "//trim(dim_name(1))//&
790  " in file "//trim(surfbandfilename)//" in MOM_wave_interace.")
791  endif
792  ! Allocating size of frequency bins
793  allocate( cs%Freq_Cen(1:id) )
794  cs%Freq_Cen(:) = 0.0
795  ! Allocating size of wavenumber bins
796  allocate( cs%WaveNum_Cen(1:id) )
797  cs%WaveNum_Cen(:) = 0.0
798  allocate( cs%STKx0(g%isdB:g%iedB,g%jsd:g%jed,1:id))
799  cs%STKx0(:,:,:) = 0.0
800  allocate( cs%STKy0(g%isd:g%ied,g%jsdB:g%jedB,1:id))
801  cs%STKy0(:,:,:) = 0.0
802  endif
803 
804  ! Reading wavenumber bins/Frequencies
805  start(:) = 1 ! Set all start to 1
806  counter(:) = 1 ! Set all counter to 1
807  counter(1) = id ! Set counter(1) to id (number of frequency bins)
808  if (partitionmode==0) then
809  rcode_wn = nf90_get_var(ncid, dim_id(1), cs%WaveNum_Cen, start, counter)
810  if (rcode_wn /= 0) then
811  call mom_error(fatal,&
812  "error reading dimension 1 values for var_name "// &
813  trim(varread1)//",dim_name "//trim(dim_name(1))// &
814  " in file "// trim(surfbandfilename)//" in MOM_wave_interface")
815  endif
816  numbands = id
817  do b = 1,numbands ; cs%WaveNum_Cen(b) = us%Z_to_m*cs%WaveNum_Cen(b) ; enddo
818  elseif (partitionmode==1) then
819  rcode_fr = nf90_get_var(ncid, dim_id(1), cs%Freq_Cen, start, counter)
820  if (rcode_fr /= 0) then
821  call mom_error(fatal,&
822  "error reading dimension 1 values for var_name "// &
823  trim(varread2)//",dim_name "//trim(dim_name(1))// &
824  " in file "// trim(surfbandfilename)//" in MOM_wave_interface")
825  endif
826  numbands = id
827  do b = 1,numbands
828  cs%WaveNum_Cen(b) = (2.*pi*cs%Freq_Cen(b)*us%T_to_s)**2 / (us%L_to_Z**2*gv%g_Earth)
829  enddo
830  endif
831 
832  endif
833 
834  do b = 1,numbands
835  temp_x(:,:) = 0.0
836  temp_y(:,:) = 0.0
837  varname = ' '
838  write(varname,"(A3,I0)")'Usx',b
839  call data_override('OCN',trim(varname), temp_x, day_center)
840  varname = ' '
841  write(varname,'(A3,I0)')'Usy',b
842  call data_override('OCN',trim(varname), temp_y, day_center)
843  ! Disperse into halo on h-grid
844  call pass_vector(temp_x, temp_y, g%Domain, to_all, agrid)
845  !Filter land values
846  do j = g%jsd,g%jed
847  do i = g%Isd,g%Ied
848  if (abs(temp_x(i,j)) > 10. .or. abs(temp_y(i,j)) > 10.) then
849  ! Assume land-mask and zero out
850  temp_x(i,j) = 0.0
851  temp_y(i,j) = 0.0
852  endif
853  enddo
854  enddo
855 
856  ! Interpolate to u/v grids
857  do j = g%jsc,g%jec
858  do i = g%IscB,g%IecB
859  cs%STKx0(i,j,b) = 0.5 * (temp_x(i,j) + temp_x(i+1,j))
860  enddo
861  enddo
862  do j = g%JscB,g%JecB
863  do i = g%isc,g%iec
864  cs%STKy0(i,j,b) = 0.5 * (temp_y(i,j) + temp_y(i,j+1))
865  enddo
866  enddo
867  ! Disperse into halo on u/v grids
868  call pass_vector(cs%STKx0(:,:,b),cs%STKy0(:,:,b), g%Domain, to_all)
869  enddo !Closes b-loop
870 

◆ update_stokes_drift()

subroutine, public mom_wave_interface::update_stokes_drift ( type(ocean_grid_type), intent(inout)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(wave_parameters_cs), pointer  CS,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in)  ustar 
)

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

Parameters
csWave parameter Control structure
[in,out]gGrid structure
[in]gvVertical grid structure
[in]usA dimensional unit scaling type
[in]hThickness [H ~> m or kg m-2]
[in]ustarWind friction velocity [Z T-1 ~> m s-1].

Definition at line 479 of file MOM_wave_interface.F90.

479  type(wave_parameters_cs), pointer :: cs !< Wave parameter Control structure
480  type(ocean_grid_type), intent(inout) :: g !< Grid structure
481  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
482  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
483  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
484  intent(in) :: h !< Thickness [H ~> m or kg m-2]
485  real, dimension(SZI_(G),SZJ_(G)), &
486  intent(in) :: ustar !< Wind friction velocity [Z T-1 ~> m s-1].
487  ! Local Variables
488  real :: top, midpoint, bottom, one_cm
489  real :: decayscale
490  real :: cmn_fac, wn, ustokes
491  real :: la
492  integer :: ii, jj, kk, b, iim1, jjm1
493 
494  one_cm = 0.01*us%m_to_Z
495 
496  ! 1. If Test Profile Option is chosen
497  ! Computing mid-point value from surface value and decay wavelength
498  if (wavemethod==testprof) then
499  decayscale = 4.*pi / tp_wvl !4pi
500  do ii = g%isdB,g%iedB
501  do jj = g%jsd,g%jed
502  iim1 = max(1,ii-1)
503  bottom = 0.0
504  midpoint = 0.0
505  do kk = 1,g%ke
506  top = bottom
507  midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(iim1,jj,kk))
508  bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(iim1,jj,kk))
509  cs%Us_x(ii,jj,kk) = tp_stkx0*exp(midpoint*decayscale)
510  enddo
511  enddo
512  enddo
513  do ii = g%isd,g%ied
514  do jj = g%jsdB,g%jedB
515  jjm1 = max(1,jj-1)
516  bottom = 0.0
517  midpoint = 0.0
518  do kk = 1,g%ke
519  top = bottom
520  midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(ii,jjm1,kk))
521  bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(ii,jjm1,kk))
522  cs%Us_y(ii,jj,kk) = tp_stky0*exp(midpoint*decayscale)
523  enddo
524  enddo
525  enddo
526  ! 2. If Surface Bands is chosen
527  ! In wavenumber mode compute integral for layer averaged Stokes drift.
528  ! In frequency mode compuate value at midpoint.
529  elseif (wavemethod==surfbands) then
530  cs%Us_x(:,:,:) = 0.0
531  cs%Us_y(:,:,:) = 0.0
532  cs%Us0_x(:,:) = 0.0
533  cs%Us0_y(:,:) = 0.0
534  ! Computing X direction Stokes drift
535  do ii = g%isdB,g%iedB
536  do jj = g%jsd,g%jed
537  ! 1. First compute the surface Stokes drift
538  ! by integrating over the partitionas.
539  do b = 1,numbands
540  if (partitionmode==0) then
541  ! In wavenumber we are averaging over (small) level
542  cmn_fac = (1.0-exp(-one_cm*2.*cs%WaveNum_Cen(b))) / &
543  (one_cm*2.*cs%WaveNum_Cen(b))
544  elseif (partitionmode==1) then
545  ! In frequency we are not averaging over level and taking top
546  cmn_fac = 1.0
547  endif
548  cs%US0_x(ii,jj) = cs%US0_x(ii,jj) + cs%STKx0(ii,jj,b)*cmn_fac
549  enddo
550  ! 2. Second compute the level averaged Stokes drift
551  bottom = 0.0
552  do kk = 1,g%ke
553  top = bottom
554  iim1 = max(ii-1,1)
555  midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(iim1,jj,kk))
556  bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(iim1,jj,kk))
557  do b = 1,numbands
558  if (partitionmode==0) then
559  ! In wavenumber we are averaging over level
560  cmn_fac = (exp(top*2.*cs%WaveNum_Cen(b))-exp(bottom*2.*cs%WaveNum_Cen(b)))&
561  / ((top-bottom)*(2.*cs%WaveNum_Cen(b)))
562  elseif (partitionmode==1) then
563  if (cs%StkLevelMode==0) then
564  ! Take the value at the midpoint
565  cmn_fac = exp(midpoint*2.*(2.*pi*cs%Freq_Cen(b)*us%T_to_s)**2/(us%L_to_Z**2*gv%g_Earth))
566  elseif (cs%StkLevelMode==1) then
567  ! Use a numerical integration and then
568  ! divide by layer thickness
569  wn = (2.*pi*cs%Freq_Cen(b)*us%T_to_s)**2 / (us%L_to_Z**2*gv%g_Earth) !bgr bug-fix missing g
570  cmn_fac = (exp(2.*wn*top)-exp(2.*wn*bottom)) / (2.*wn*(top-bottom))
571  endif
572  endif
573  cs%US_x(ii,jj,kk) = cs%US_x(ii,jj,kk) + cs%STKx0(ii,jj,b)*cmn_fac
574  enddo
575  enddo
576  enddo
577  enddo
578  ! Computing Y direction Stokes drift
579  do ii = g%isd,g%ied
580  do jj = g%jsdB,g%jedB
581  ! Compute the surface values.
582  do b = 1,numbands
583  if (partitionmode==0) then
584  ! In wavenumber we are averaging over (small) level
585  cmn_fac = (1.0-exp(-one_cm*2.*cs%WaveNum_Cen(b))) / &
586  (one_cm*2.*cs%WaveNum_Cen(b))
587  elseif (partitionmode==1) then
588  ! In frequency we are not averaging over level and taking top
589  cmn_fac = 1.0
590  endif
591  cs%US0_y(ii,jj) = cs%US0_y(ii,jj) + cs%STKy0(ii,jj,b)*cmn_fac
592  enddo
593  ! Compute the level averages.
594  bottom = 0.0
595  do kk = 1,g%ke
596  top = bottom
597  jjm1 = max(jj-1,1)
598  midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(ii,jjm1,kk))
599  bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(ii,jjm1,kk))
600  do b = 1,numbands
601  if (partitionmode==0) then
602  ! In wavenumber we are averaging over level
603  cmn_fac = (exp(top*2.*cs%WaveNum_Cen(b)) - &
604  exp(bottom*2.*cs%WaveNum_Cen(b))) / &
605  ((top-bottom)*(2.*cs%WaveNum_Cen(b)))
606  elseif (partitionmode==1) then
607  if (cs%StkLevelMode==0) then
608  ! Take the value at the midpoint
609  cmn_fac = exp(midpoint*2.*(2.*pi*cs%Freq_Cen(b)*us%T_to_s)**2/(us%L_to_Z**2*gv%g_Earth))
610  elseif (cs%StkLevelMode==1) then
611  ! Use a numerical integration and then
612  ! divide by layer thickness
613  wn = (2.*pi*cs%Freq_Cen(b)*us%T_to_s)**2 / (us%L_to_Z**2*gv%g_Earth)
614  cmn_fac = (exp(2.*wn*top)-exp(2.*wn*bottom)) / (2.*wn*(top-bottom))
615  endif
616  endif
617  cs%US_y(ii,jj,kk) = cs%US_y(ii,jj,kk) + cs%STKy0(ii,jj,b)*cmn_fac
618  enddo
619  enddo
620  enddo
621  enddo
622  elseif (wavemethod==dhh85) then
623  if (.not.(staticwaves .and. dhh85_is_set)) then
624  do ii = g%isdB,g%iedB
625  do jj = g%jsd,g%jed
626  bottom = 0.0
627  do kk = 1,g%ke
628  top = bottom
629  iim1 = max(ii-1,1)
630  midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(iim1,jj,kk))
631  bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(iim1,jj,kk))
632  !bgr note that this is using a u-point ii on h-point ustar
633  ! this code has only been previous used for uniform
634  ! grid cases. This needs fixed if DHH85 is used for non
635  ! uniform cases.
636  call dhh85_mid(gv, us, midpoint, ustokes)
637  ! Putting into x-direction (no option for direction
638  cs%US_x(ii,jj,kk) = ustokes
639  enddo
640  enddo
641  enddo
642  do ii = g%isd,g%ied
643  do jj = g%jsdB,g%jedB
644  bottom = 0.0
645  do kk=1, g%ke
646  top = bottom
647  jjm1 = max(jj-1,1)
648  midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(ii,jjm1,kk))
649  bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(ii,jjm1,kk))
650  !bgr note that this is using a v-point jj on h-point ustar
651  ! this code has only been previous used for uniform
652  ! grid cases. This needs fixed if DHH85 is used for non
653  ! uniform cases.
654  ! call DHH85_mid(GV, US, Midpoint, UStokes)
655  ! Putting into x-direction, so setting y direction to 0
656  cs%US_y(ii,jj,kk) = 0.0
657  ! For rotational symmetry there should be the option for this to become = UStokes
658  ! bgr - see note above, but this is true
659  ! if this is used for anything
660  ! other than simple LES comparison
661  enddo
662  enddo
663  enddo
664  dhh85_is_set = .true.
665  endif
666  else! Keep this else, fallback to 0 Stokes drift
667  do kk= 1,g%ke
668  do ii = g%isdB,g%iedB
669  do jj = g%jsd,g%jed
670  cs%Us_x(ii,jj,kk) = 0.
671  enddo
672  enddo
673  do ii = g%isd,g%ied
674  do jj = g%jsdB,g%jedB
675  cs%Us_y(ii,jj,kk) = 0.
676  enddo
677  enddo
678  enddo
679  endif
680 
681  ! Turbulent Langmuir number is computed here and available to use anywhere.
682  ! SL Langmuir number requires mixing layer depth, and therefore is computed
683  ! in the routine it is needed by (e.g. KPP or ePBL).
684  do ii = g%isc,g%iec
685  do jj = g%jsc, g%jec
686  top = h(ii,jj,1)*gv%H_to_Z
687  call get_langmuir_number( la, g, gv, us, top, ustar(ii,jj), ii, jj, &
688  h(ii,jj,:),override_ma=.false.,waves=cs)
689  cs%La_turb(ii,jj) = la
690  enddo
691  enddo
692 
693  ! Output any desired quantities
694  if (cs%id_surfacestokes_y>0) &
695  call post_data(cs%id_surfacestokes_y, cs%us0_y, cs%diag)
696  if (cs%id_surfacestokes_x>0) &
697  call post_data(cs%id_surfacestokes_x, cs%us0_x, cs%diag)
698  if (cs%id_3dstokes_y>0) &
699  call post_data(cs%id_3dstokes_y, cs%us_y, cs%diag)
700  if (cs%id_3dstokes_x>0) &
701  call post_data(cs%id_3dstokes_x, cs%us_x, cs%diag)
702  if (cs%id_La_turb>0) &
703  call post_data(cs%id_La_turb, cs%La_turb, cs%diag)
704 

◆ update_surface_waves()

subroutine, public mom_wave_interface::update_surface_waves ( type(ocean_grid_type), intent(inout)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(time_type), intent(in)  Day,
type(time_type), intent(in)  dt,
type(wave_parameters_cs), pointer  CS 
)

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

Parameters
csWave parameter Control structure
[in,out]gGrid structure
[in]gvVertical grid structure
[in]usA dimensional unit scaling type
[in]dayCurrent model time
[in]dtTimestep as a time-type

Definition at line 437 of file MOM_wave_interface.F90.

437  type(wave_parameters_cs), pointer :: cs !< Wave parameter Control structure
438  type(ocean_grid_type), intent(inout) :: g !< Grid structure
439  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
440  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
441  type(time_type), intent(in) :: day !< Current model time
442  type(time_type), intent(in) :: dt !< Timestep as a time-type
443  ! Local variables
444  integer :: ii, jj, kk, b
445  type(time_type) :: day_center
446 
447  ! Computing central time of time step
448  day_center = day + dt/2
449 
450  if (wavemethod == testprof) then
451  ! Do nothing
452  elseif (wavemethod==surfbands) then
453  if (datasource==dataovr) then
454  call surface_bands_by_data_override(day_center, g, gv, us, cs)
455  elseif (datasource==coupler) then
456  ! Reserve for coupler hooks
457  elseif (datasource==input) then
458  do b=1,numbands
459  do ii=g%isdB,g%iedB
460  do jj=g%jsd,g%jed
461  cs%STKx0(ii,jj,b) = cs%PrescribedSurfStkX(b)
462  enddo
463  enddo
464  do ii=g%isd,g%ied
465  do jj=g%jsdB, g%jedB
466  cs%STKY0(ii,jj,b) = cs%PrescribedSurfStkY(b)
467  enddo
468  enddo
469  enddo
470  endif
471  endif
472 
473  return

◆ ust_2_u10_coare3p5()

subroutine mom_wave_interface::ust_2_u10_coare3p5 ( real, intent(in)  USTair,
real, intent(out)  U10,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US 
)
private

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
[in]ustairWind friction velocity [m s-1]
[out]u1010-m neutral wind speed [m s-1]
[in]gvvertical grid type
[in]usA dimensional unit scaling type

Definition at line 1332 of file MOM_wave_interface.F90.

1332  real, intent(in) :: ustair !< Wind friction velocity [m s-1]
1333  real, intent(out) :: u10 !< 10-m neutral wind speed [m s-1]
1334  type(verticalgrid_type), intent(in) :: gv !< vertical grid type
1335  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1336 
1337  ! Local variables
1338  real, parameter :: vonkar = 0.4 ! Should access a get_param von karman
1339  real, parameter :: nu=1e-6 ! Should access a get_param air-viscosity
1340  real :: z0sm, z0, z0rough, u10a, alpha, cd
1341  integer :: ct
1342 
1343  ! Uses empirical formula for z0 to convert ustar_air to u10 based on the
1344  ! COARE 3.5 paper (Edson et al., 2013)
1345  ! alpha=m*U10+b
1346  ! Note in Edson et al. 2013, eq. 13 m is given as 0.017. However,
1347  ! m=0.0017 reproduces the curve in their figure 6.
1348 
1349  z0sm = 0.11 * nu * us%m_to_Z / ustair !Compute z0smooth from ustar guess
1350  u10 = ustair/sqrt(0.001) !Guess for u10
1351  u10a = 1000
1352 
1353  ct=0
1354  do while (abs(u10a/u10-1.) > 0.001)
1355  ct=ct+1
1356  u10a = u10
1357  alpha = min(0.028, 0.0017 * u10 - 0.005)
1358  z0rough = alpha * (us%m_s_to_L_T*ustair)**2 / gv%g_Earth ! Compute z0rough from ustar guess
1359  z0 = z0sm + z0rough
1360  cd = ( vonkar / log(10.*us%m_to_Z / z0) )**2 ! Compute CD from derived roughness
1361  u10 = ustair/sqrt(cd) ! Compute new u10 from derived CD, while loop
1362  ! ends and checks for convergence...CT counter
1363  ! makes sure loop doesn't run away if function
1364  ! doesn't converge. This code was produced offline
1365  ! and converged rapidly (e.g. 2 cycles)
1366  ! for ustar=0.0001:0.0001:10.
1367  if (ct>20) then
1368  u10 = ustair/sqrt(0.0015) ! I don't expect to get here, but just
1369  ! in case it will output a reasonable value.
1370  exit
1371  endif
1372  enddo
1373  return

◆ waves_end()

subroutine, public mom_wave_interface::waves_end ( type(wave_parameters_cs), pointer  CS)

Clear pointers, deallocate memory.

Parameters
csControl structure

Definition at line 1378 of file MOM_wave_interface.F90.

1378  type(wave_parameters_cs), pointer :: cs !< Control structure
1379 
1380  if (allocated(cs%WaveNum_Cen)) deallocate( cs%WaveNum_Cen )
1381  if (allocated(cs%Freq_Cen)) deallocate( cs%Freq_Cen )
1382  if (allocated(cs%Us_x)) deallocate( cs%Us_x )
1383  if (allocated(cs%Us_y)) deallocate( cs%Us_y )
1384  if (allocated(cs%La_SL)) deallocate( cs%La_SL )
1385  if (allocated(cs%La_turb)) deallocate( cs%La_turb )
1386  if (allocated(cs%STKx0)) deallocate( cs%STKx0 )
1387  if (allocated(cs%STKy0)) deallocate( cs%STKy0 )
1388  if (allocated(cs%KvS)) deallocate( cs%KvS )
1389  if (allocated(cs%Us0_y)) deallocate( cs%Us0_y )
1390  if (allocated(cs%Us0_x)) deallocate( cs%Us0_x )
1391 
1392  deallocate( cs )
1393 
1394  return

Variable Documentation

◆ coupler

integer, parameter mom_wave_interface::coupler = 2
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 177 of file MOM_wave_interface.F90.

◆ dataoverrideisinitialized

logical mom_wave_interface::dataoverrideisinitialized
private

Flag for DataOverride Initialization.

Todo:
Module variable! Move into a control structure.

Definition at line 159 of file MOM_wave_interface.F90.

159 logical :: dataoverrideisinitialized !< Flag for DataOverride Initialization

◆ dataovr

integer, parameter mom_wave_interface::dataovr = 1
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 177 of file MOM_wave_interface.F90.

◆ datasource

integer mom_wave_interface::datasource
private

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)

Todo:
Module variable! Move into a control structure.

Definition at line 149 of file MOM_wave_interface.F90.

149 integer :: datasource !< Integer that specifies where the Model Looks for Data

◆ dhh85

integer, parameter mom_wave_interface::dhh85 = 2
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 177 of file MOM_wave_interface.F90.

◆ dhh85_is_set

logical mom_wave_interface::dhh85_is_set
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 184 of file MOM_wave_interface.F90.

◆ input

integer, parameter mom_wave_interface::input = 3
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 177 of file MOM_wave_interface.F90.

◆ la_frachbl

real mom_wave_interface::la_frachbl
private

Fraction of OSBL for averaging Langmuir number.

Todo:
Module variable! Move into a control structure.

Definition at line 163 of file MOM_wave_interface.F90.

163 real :: la_frachbl !< Fraction of OSBL for averaging Langmuir number

◆ la_misalignment

logical mom_wave_interface::la_misalignment = .false.
private

Flag to use misalignment in Langmuir number.

Todo:
Module variable! Move into a control structure.

Definition at line 165 of file MOM_wave_interface.F90.

165 logical :: la_misalignment = .false. !< Flag to use misalignment in Langmuir number

◆ lf17

integer, parameter mom_wave_interface::lf17 = 3
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 177 of file MOM_wave_interface.F90.

◆ null_wavemethod

integer, parameter mom_wave_interface::null_wavemethod =-99
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 177 of file MOM_wave_interface.F90.

◆ numbands

integer, public mom_wave_interface::numbands =0

Number of wavenumber/frequency partitions to receive This needs to match the number of bands provided via either coupling or file.

Todo:
Module variable! Move into a control structure.

Definition at line 141 of file MOM_wave_interface.F90.

141 integer, public :: numbands =0 !< Number of wavenumber/frequency partitions to receive

◆ partitionmode

integer, public mom_wave_interface::partitionmode

Method for partition mode (meant to check input) 0 - wavenumbers 1 - frequencies.

Todo:
Module variable! Move into a control structure.

Definition at line 145 of file MOM_wave_interface.F90.

145 integer, public :: partitionmode !< Method for partition mode (meant to check input)

◆ pi

real mom_wave_interface::pi
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 186 of file MOM_wave_interface.F90.

186 real :: pi

◆ staticwaves

logical mom_wave_interface::staticwaves
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 184 of file MOM_wave_interface.F90.

184 logical :: staticwaves, dhh85_is_set

◆ surfbandfilename

character(len=40) mom_wave_interface::surfbandfilename
private

Filename if using DataOverride.

Todo:
Module variable! Move into a control structure.

Definition at line 157 of file MOM_wave_interface.F90.

157 character(len=40) :: SurfBandFileName !< Filename if using DataOverride

◆ surfbands

integer, parameter mom_wave_interface::surfbands = 1
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 177 of file MOM_wave_interface.F90.

◆ testprof

integer, parameter mom_wave_interface::testprof = 0
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 177 of file MOM_wave_interface.F90.

177 integer, parameter :: testprof = 0, surfbands = 1, &
178  dhh85 = 2, lf17 = 3, null_wavemethod=-99, &
179  dataovr = 1, coupler = 2, input = 3

◆ tp_stkx0

real mom_wave_interface::tp_stkx0
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 182 of file MOM_wave_interface.F90.

182 Real :: tp_stkx0, tp_stky0, tp_wvl

◆ tp_stky0

real mom_wave_interface::tp_stky0
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 182 of file MOM_wave_interface.F90.

◆ tp_wvl

real mom_wave_interface::tp_wvl
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 182 of file MOM_wave_interface.F90.

◆ waveage

real mom_wave_interface::waveage
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 185 of file MOM_wave_interface.F90.

185 real :: waveage, wavewind

◆ waveagepeakfreq

logical mom_wave_interface::waveagepeakfreq
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 183 of file MOM_wave_interface.F90.

183 logical :: waveagepeakfreq ! Flag to use W

◆ wavemethod

integer mom_wave_interface::wavemethod =-99

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.

Todo:
Module variable! Move into a control structure.

Definition at line 131 of file MOM_wave_interface.F90.

131 integer :: wavemethod=-99 !< Options for including wave information

◆ wavewind

real mom_wave_interface::wavewind
private

Undocumented parameters.

Todo:
These module variables need to be documented as static/private variables or moved into a control structure.

Definition at line 185 of file MOM_wave_interface.F90.