Skip to content

Commit ee42fdb

Browse files
awallcraftkshedstrom
authored andcommitted
Tidal angular frequency has units [rad s-1] (#764)
* Tidal angular frequency has units [rad s-1] Tidal frequencies are always angular frequencies to simplify applying sine and cosine. These have MKS units [rad s-1] but they are all currently listed as [s-1]. Updated dOxygen comments for variables, e.g. [T-1 ~> s-1] becomes [rad T-1 ~> rad s-1]. Updated get_param units. e.g. units="s-1" becomes units="rad s-1". No answers are changed, but the logged parameter units are different. There are frequencies in MOM_internal_tides.F90 but these have not been updated because they may be specified incorrectly. They are used as if they are [T-1] but they are calculated as 2PI/period [rad T-1]. real, allocatable, dimension(:) :: frequency !< The frequency of each band [T-1 ~> s-1]. real :: period ! A tidal period read from namelist [T ~> s] ! The periods of the tidal constituents for internal tides raytracing call read_param(param_file, "TIDAL_PERIODS", periods) do fr=1,num_freq period = US%s_to_T*extract_real(periods, " ,", fr, 0.) CS%frequency(fr) = 8.0*atan(1.0)/period enddo All MOM6-examples cases have INTERNAL_TIDES=False and so can't resolve this issue. * fixed too-long line
1 parent 8834f68 commit ee42fdb

File tree

4 files changed

+13
-11
lines changed

4 files changed

+13
-11
lines changed

src/core/MOM_barotropic.F90

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -5323,7 +5323,7 @@ subroutine register_barotropic_restarts(HI, GV, US, param_file, CS, restart_CS)
53235323
character(len=40) :: mdl = "MOM_barotropic" ! This module's name.
53245324
integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
53255325
real :: am2, ak1 !< Bandwidth parameters of the M2 and K1 streaming filters [nondim]
5326-
real :: om2, ok1 !< Target frequencies of the M2 and K1 streaming filters [T-1 ~> s-1]
5326+
real :: om2, ok1 !< Target frequencies of the M2 and K1 streaming filters [rad T-1 ~> rad s-1]
53275327

53285328
isd = HI%isd ; ied = HI%ied ; jsd = HI%jsd ; jed = HI%jed
53295329
IsdB = HI%IsdB ; IedB = HI%IedB ; JsdB = HI%JsdB ; JedB = HI%JedB
@@ -5354,13 +5354,13 @@ subroutine register_barotropic_restarts(HI, GV, US, param_file, CS, restart_CS)
53545354
"Frequency of the M2 tidal constituent. "//&
53555355
"This is only used if TIDES and TIDE_M2"// &
53565356
" are true, or if OBC_TIDE_N_CONSTITUENTS > 0 and M2"// &
5357-
" is in OBC_TIDE_CONSTITUENTS.", units="s-1", default=tidal_frequency("M2"), &
5357+
" is in OBC_TIDE_CONSTITUENTS.", units="rad s-1", default=tidal_frequency("M2"), &
53585358
scale=US%T_to_s, do_not_log=.true.)
53595359
call get_param(param_file, mdl, "TIDE_K1_FREQ", ok1, &
53605360
"Frequency of the K1 tidal constituent. "//&
53615361
"This is only used if TIDES and TIDE_K1"// &
53625362
" are true, or if OBC_TIDE_N_CONSTITUENTS > 0 and K1"// &
5363-
" is in OBC_TIDE_CONSTITUENTS.", units="s-1", default=tidal_frequency("K1"), &
5363+
" is in OBC_TIDE_CONSTITUENTS.", units="rad s-1", default=tidal_frequency("K1"), &
53645364
scale=US%T_to_s, do_not_log=.true.)
53655365

53665366
ALLOC_(CS%ubtav(IsdB:IedB,jsd:jed)) ; CS%ubtav(:,:) = 0.0

src/core/MOM_open_boundary.F90

Lines changed: 3 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -320,7 +320,8 @@ module MOM_open_boundary
320320
logical :: add_tide_constituents = .false. !< If true, add tidal constituents to the boundary elevation
321321
!! and velocity. Will be set to true if n_tide_constituents > 0.
322322
character(len=2), allocatable, dimension(:) :: tide_names !< Names of tidal constituents to add to the boundary data.
323-
real, allocatable, dimension(:) :: tide_frequencies !< Angular frequencies of chosen tidal constituents [T-1 ~> s-1].
323+
real, allocatable, dimension(:) :: tide_frequencies !< Angular frequencies of chosen tidal
324+
!! constituents [rad T-1 ~> rad s-1].
324325
real, allocatable, dimension(:) :: tide_eq_phases !< Equilibrium phases of chosen tidal constituents [rad].
325326
real, allocatable, dimension(:) :: tide_fn !< Amplitude modulation of boundary tides by nodal cycle [nondim].
326327
real, allocatable, dimension(:) :: tide_un !< Phase modulation of boundary tides by nodal cycle [rad].
@@ -1234,7 +1235,7 @@ subroutine initialize_obc_tides(OBC, US, param_file)
12341235
"This is only used if TIDES and TIDE_"//trim(OBC%tide_names(c))// &
12351236
" are true, or if OBC_TIDE_N_CONSTITUENTS > 0 and "//trim(OBC%tide_names(c))//&
12361237
" is in OBC_TIDE_CONSTITUENTS.", &
1237-
units="s-1", default=tidal_frequency(trim(OBC%tide_names(c))), scale=US%T_to_s)
1238+
units="rad s-1", default=tidal_frequency(trim(OBC%tide_names(c))), scale=US%T_to_s)
12381239

12391240
! Find equilibrium phase if needed
12401241
if (OBC%add_eq_phase) then

src/parameterizations/lateral/MOM_streaming_filter.F90

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -15,7 +15,7 @@ module MOM_streaming_filter
1515
!> The control structure for storing the filter infomation of a particular field
1616
type, public :: Filter_CS ; private
1717
real :: a, & !< Parameter that determines the bandwidth [nondim]
18-
om, & !< Target frequency of the filter [T-1 ~> s-1]
18+
om, & !< Target frequency of the filter [rad T-1 ~> rad s-1]
1919
old_time = -1.0 !< The time of the previous accumulating step [T ~> s]
2020
real, allocatable, dimension(:,:) :: s1, & !< Dummy variable [A]
2121
u1 !< Filtered data [A]
@@ -29,7 +29,7 @@ module MOM_streaming_filter
2929
!> This subroutine registers each of the fields to be filtered.
3030
subroutine Filt_register(a, om, grid, HI, CS)
3131
real, intent(in) :: a !< Parameter that determines the bandwidth [nondim]
32-
real, intent(in) :: om !< Target frequency of the filter [T-1 ~> s-1]
32+
real, intent(in) :: om !< Target frequency of the filter [rad T-1 ~> rad s-1]
3333
character(len=*), intent(in) :: grid !< Horizontal grid location: h, u, or v
3434
type(hor_index_type), intent(in) :: HI !< Horizontal index type structure
3535
type(Filter_CS), intent(out) :: CS !< Control structure for the current field

src/parameterizations/lateral/MOM_tidal_forcing.F90

Lines changed: 5 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -50,7 +50,7 @@ module MOM_tidal_forcing
5050
!! and bottom geopotential anomalies [nondim].
5151
integer :: nc !< The number of tidal constituents in use.
5252
real, dimension(MAX_CONSTITUENTS) :: &
53-
freq, & !< The frequency of a tidal constituent [T-1 ~> s-1].
53+
freq, & !< The frequency of a tidal constituent [rad T-1 ~> rad s-1].
5454
phase0, & !< The phase of a tidal constituent at time 0 [rad].
5555
amp, & !< The amplitude of a tidal constituent at time 0 [Z ~> m].
5656
love_no !< The Love number of a tidal constituent at time 0 [nondim].
@@ -151,7 +151,7 @@ end function eq_phase
151151
!! Values used here are from previous versions of MOM.
152152
function tidal_frequency(constit)
153153
character (len=2), intent(in) :: constit !> Constituent to look up
154-
real :: tidal_frequency !> Angular frequency [s-1]
154+
real :: tidal_frequency !> Angular frequency [rad s-1]
155155

156156
select case (constit)
157157
case ("M2")
@@ -246,7 +246,7 @@ subroutine tidal_forcing_init(Time, G, US, param_file, CS, HA_CS)
246246
phase, & ! The phase of some tidal constituent [radians].
247247
lat_rad, lon_rad ! Latitudes and longitudes of h-points [radians].
248248
real :: deg_to_rad ! A conversion factor from degrees to radians [radian degree-1]
249-
real, dimension(MAX_CONSTITUENTS) :: freq_def ! Default frequency for each tidal constituent [s-1]
249+
real, dimension(MAX_CONSTITUENTS) :: freq_def ! Default frequency for each tidal constituent [rad s-1]
250250
real, dimension(MAX_CONSTITUENTS) :: phase0_def ! Default reference phase for each tidal constituent [rad]
251251
real, dimension(MAX_CONSTITUENTS) :: amp_def ! Default amplitude for each tidal constituent [m]
252252
real, dimension(MAX_CONSTITUENTS) :: love_def ! Default love number for each constituent [nondim]
@@ -480,7 +480,8 @@ subroutine tidal_forcing_init(Time, G, US, param_file, CS, HA_CS)
480480
"Frequency of the "//trim(CS%const_name(c))//" tidal constituent. "//&
481481
"This is only used if TIDES and TIDE_"//trim(CS%const_name(c))// &
482482
" are true, or if OBC_TIDE_N_CONSTITUENTS > 0 and "//trim(CS%const_name(c))// &
483-
" is in OBC_TIDE_CONSTITUENTS.", units="s-1", default=freq_def(c), scale=US%T_to_s)
483+
" is in OBC_TIDE_CONSTITUENTS.", units="rad s-1", default=freq_def(c), &
484+
scale=US%T_to_s)
484485
call get_param(param_file, mdl, "TIDE_"//trim(CS%const_name(c))//"_AMP", CS%amp(c), &
485486
"Amplitude of the "//trim(CS%const_name(c))//" tidal constituent. "//&
486487
"This is only used if TIDES and TIDE_"//trim(CS%const_name(c))// &

0 commit comments

Comments
 (0)