Skip to content

Commit 84bc0ad

Browse files
committed
Streaming Filter
Fixed dimensional inconsistency of the om2 and ok1 parameters.
1 parent 4c7fbfe commit 84bc0ad

File tree

2 files changed

+15
-5
lines changed

2 files changed

+15
-5
lines changed

src/core/MOM_barotropic.F90

Lines changed: 13 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -4455,7 +4455,7 @@ subroutine barotropic_init(u, v, h, eta, Time, G, GV, US, param_file, diag, CS,
44554455
real :: wave_drag_scale ! A scaling factor for the barotropic linear wave drag
44564456
! piston velocities [nondim].
44574457
real :: am2, ak1 !< Bandwidth parameters of the M2 and K1 streaming filters [nondim]
4458-
real :: om2, ok1 !< Target frequencies of the M2 and K1 streaming filters [s-1]
4458+
real :: om2, ok1 !< Target frequencies of the M2 and K1 streaming filters [T-1 ~> s-1]
44594459
character(len=200) :: inputdir ! The directory in which to find input files.
44604460
character(len=200) :: wave_drag_file ! The file from which to read the wave
44614461
! drag piston velocity.
@@ -4717,8 +4717,18 @@ subroutine barotropic_init(u, v, h, eta, Time, G, GV, US, param_file, diag, CS,
47174717
"Bandwidth parameter of the streaming filter targeting the K1 frequency. "//&
47184718
"Must be positive. To turn off filtering, set FILTER_ALPHA_K1 <= 0.0.", &
47194719
default=0.0, units="nondim")
4720-
om2 = 1.4051890e-4
4721-
ok1 = 0.7292117e-4
4720+
call get_param(param_file, mdl, "TIDE_M2_FREQ", om2, &
4721+
"Frequency of the M2 tidal constituent. "//&
4722+
"This is only used if TIDES and TIDE_M2"// &
4723+
" are true, or if OBC_TIDE_N_CONSTITUENTS > 0 and M2"// &
4724+
" is in OBC_TIDE_CONSTITUENTS.", units="s-1", default=1.4051890e-4, &
4725+
scale=US%T_to_s, do_not_log=.true.)
4726+
call get_param(param_file, mdl, "TIDE_K1_FREQ", ok1, &
4727+
"Frequency of the K1 tidal constituent. "//&
4728+
"This is only used if TIDES and TIDE_K1"// &
4729+
" are true, or if OBC_TIDE_N_CONSTITUENTS > 0 and K1"// &
4730+
" is in OBC_TIDE_CONSTITUENTS.", units="s-1", default=0.7292117e-4, &
4731+
scale=US%T_to_s, do_not_log=.true.)
47224732

47234733
call get_param(param_file, mdl, "CLIP_BT_VELOCITY", CS%clip_velocity, &
47244734
"If true, limit any velocity components that exceed "//&

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
type, private :: Filt_type
1616
character(len=16) :: key = 'none' !< Name of the current field
1717
real :: a, & !< Parameter that determines the bandwidth [nondim]
18-
om, & !< Target frequency of the filter [s-1]
18+
om, & !< Target frequency of the filter [T-1 ~> s-1]
1919
old_time = -1.0 !< The time of the previous accumulating step [T ~> s]
2020
integer :: is, ie, js, je !< Lower and upper bounds of input data
2121
real, allocatable :: s1(:,:), & !< Dummy variable [A]
@@ -54,7 +54,7 @@ end subroutine Filt_init
5454
subroutine Filt_register(key, a, om, CS)
5555
character(len=*), intent(in) :: key !< Name of the current field
5656
real, intent(in) :: a !< Parameter that determines the bandwidth [nondim]
57-
real, intent(in) :: om !< Target frequency of the filter [s-1]
57+
real, intent(in) :: om !< Target frequency of the filter [T-1 ~> s-1]
5858
type(streaming_filter_CS), intent(inout) :: CS !< Control structure of the MOM_streaming_filter module
5959

6060
! Local variables

0 commit comments

Comments
 (0)