Skip to content

Commit 211cab5

Browse files
committed
Nodal modulation
This commit fixes a few (potential) inconsistencies between open boundary tidal forcing and astronomical tidal forcing. 1. There was an inconsistency in the code that the nodal modulation can be calculated in OBC tidal forcing but not in the astronomical tidal forcing. This commit fixes this bug so that nodal modulation is applied to both forcings. 2. In the previous version of MOM_open_boundary.F90, a different set of tidal parameters can be set for open boundary tidal forcing, leading to potential inconsistency with astronomical tidal forcing. This commit obsoletes those for open boundary tidal forcing. 3. Another important bug fix is that the equilibrium phase is added to the SAL term, which was missing in the original code.
1 parent 91eee52 commit 211cab5

File tree

4 files changed

+93
-32
lines changed

4 files changed

+93
-32
lines changed

src/core/MOM_open_boundary.F90

Lines changed: 12 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -1127,26 +1127,30 @@ subroutine initialize_obc_tides(OBC, US, param_file)
11271127
type(astro_longitudes) :: nodal_longitudes !< Solar and lunar longitudes for tidal forcing
11281128
type(time_type) :: nodal_time !< Model time to calculate nodal modulation for.
11291129
integer :: c !< Index to tidal constituent.
1130+
logical :: tides !< True if astronomical tides are also used.
11301131

11311132
call get_param(param_file, mdl, "OBC_TIDE_CONSTITUENTS", tide_constituent_str, &
11321133
"Names of tidal constituents being added to the open boundaries.", &
11331134
fail_if_missing=.true.)
11341135

1135-
call get_param(param_file, mdl, "OBC_TIDE_ADD_EQ_PHASE", OBC%add_eq_phase, &
1136+
call get_param(param_file, mdl, "TIDES", tides, &
1137+
"If true, apply tidal momentum forcing.", default=.false., do_not_log=.true.)
1138+
1139+
call get_param(param_file, mdl, "TIDE_USE_EQ_PHASE", OBC%add_eq_phase, &
11361140
"If true, add the equilibrium phase argument to the specified tidal phases.", &
1137-
default=.false., fail_if_missing=.false.)
1141+
default=.false., do_not_log=tides)
11381142

1139-
call get_param(param_file, mdl, "OBC_TIDE_ADD_NODAL", OBC%add_nodal_terms, &
1143+
call get_param(param_file, mdl, "TIDE_ADD_NODAL", OBC%add_nodal_terms, &
11401144
"If true, include 18.6 year nodal modulation in the boundary tidal forcing.", &
1141-
default=.false.)
1145+
default=.false., do_not_log=tides)
11421146

1143-
call get_param(param_file, mdl, "OBC_TIDE_REF_DATE", tide_ref_date, &
1147+
call get_param(param_file, mdl, "TIDE_REF_DATE", tide_ref_date, &
11441148
"Reference date to use for tidal calculations and equilibrium phase.", &
1145-
fail_if_missing=.true.)
1149+
default=0, do_not_log=tides)
11461150

1147-
call get_param(param_file, mdl, "OBC_TIDE_NODAL_REF_DATE", nodal_ref_date, &
1151+
call get_param(param_file, mdl, "TIDE_NODAL_REF_DATE", nodal_ref_date, &
11481152
"Fixed reference date to use for nodal modulation of boundary tides.", &
1149-
fail_if_missing=.false., default=0)
1153+
default=0, do_not_log=tides)
11501154

11511155
if (.not. OBC%add_eq_phase) then
11521156
! If equilibrium phase argument is not added, the input phases

src/diagnostics/MOM_harmonic_analysis.F90

Lines changed: 17 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -46,7 +46,9 @@ module MOM_harmonic_analysis
4646
time_ref !< Reference time (t = 0) used to calculate tidal forcing
4747
real, dimension(MAX_CONSTITUENTS) :: &
4848
freq, & !< The frequency of a tidal constituent [T-1 ~> s-1]
49-
phase0 !< The phase of a tidal constituent at time 0 [rad]
49+
phase0, & !< The phase of a tidal constituent at time 0 [rad]
50+
tide_fn, & !< Amplitude modulation of tides by nodal cycle [nondim].
51+
tide_un !< Phase modulation of tides by nodal cycle [rad].
5052
real, allocatable :: FtF(:,:) !< Accumulator of (F' * F) for all fields [nondim]
5153
integer :: nc !< The number of tidal constituents in use
5254
integer :: length !< Number of fields of which harmonic analysis is to be performed
@@ -60,13 +62,16 @@ module MOM_harmonic_analysis
6062

6163
!> This subroutine sets static variables used by this module and initializes CS%list.
6264
!! THIS MUST BE CALLED AT THE END OF tidal_forcing_init.
63-
subroutine HA_init(Time, US, param_file, time_ref, nc, freq, phase0, const_name, CS)
65+
subroutine HA_init(Time, US, param_file, time_ref, nc, freq, phase0, const_name, tide_fn, tide_un, CS)
6466
type(time_type), intent(in) :: Time !< The current model time
6567
type(time_type), intent(in) :: time_ref !< Reference time (t = 0) used to calculate tidal forcing
6668
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
6769
type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
68-
real, dimension(MAX_CONSTITUENTS), intent(in) :: freq !< The frequency of a tidal constituent [T-1 ~> s-1]
69-
real, dimension(MAX_CONSTITUENTS), intent(in) :: phase0 !< The phase of a tidal constituent at time 0 [rad]
70+
real, dimension(MAX_CONSTITUENTS), intent(in) :: &
71+
freq, & !< The frequency of a tidal constituent [T-1 ~> s-1]
72+
phase0, & !< The phase of a tidal constituent at time 0 [rad]
73+
tide_fn, & !< Amplitude modulation of tides by nodal cycle [nondim].
74+
tide_un !< Phase modulation of tides by nodal cycle [rad].
7075
integer, intent(in) :: nc !< The number of tidal constituents in use
7176
character(len=16), intent(in) :: const_name(MAX_CONSTITUENTS) !< The name of each constituent
7277
type(harmonic_analysis_CS), intent(out) :: CS !< Control structure of the MOM_harmonic_analysis module
@@ -135,6 +140,8 @@ subroutine HA_init(Time, US, param_file, time_ref, nc, freq, phase0, const_name,
135140
CS%time_ref = time_ref
136141
CS%freq = freq
137142
CS%phase0 = phase0
143+
CS%tide_fn = tide_fn
144+
CS%tide_un = tide_un
138145
CS%nc = nc
139146
CS%const_name = const_name
140147
CS%length = 0
@@ -196,17 +203,17 @@ subroutine HA_accum_FtF(Time, CS)
196203
do c=1,nc
197204
icos = 2*c
198205
isin = 2*c+1
199-
cosomegat = cos(CS%freq(c) * now + CS%phase0(c))
200-
sinomegat = sin(CS%freq(c) * now + CS%phase0(c))
206+
cosomegat = CS%tide_fn(c) * cos(CS%freq(c) * now + CS%phase0(c) + CS%tide_un(c))
207+
sinomegat = CS%tide_fn(c) * sin(CS%freq(c) * now + CS%phase0(c) + CS%tide_un(c))
201208
CS%FtF(icos,1) = CS%FtF(icos,1) + cosomegat
202209
CS%FtF(isin,1) = CS%FtF(isin,1) + sinomegat
203210
CS%FtF(1,icos) = CS%FtF(icos,1)
204211
CS%FtF(1,isin) = CS%FtF(isin,1)
205212
do cc=c,nc
206213
iccos = 2*cc
207214
issin = 2*cc+1
208-
ccosomegat = cos(CS%freq(cc) * now + CS%phase0(cc))
209-
ssinomegat = sin(CS%freq(cc) * now + CS%phase0(cc))
215+
ccosomegat = CS%tide_fn(cc) * cos(CS%freq(cc) * now + CS%phase0(cc) + CS%tide_un(cc))
216+
ssinomegat = CS%tide_fn(cc) * sin(CS%freq(cc) * now + CS%phase0(cc) + CS%tide_un(cc))
210217
CS%FtF(icos,iccos) = CS%FtF(icos,iccos) + cosomegat * ccosomegat
211218
CS%FtF(icos,issin) = CS%FtF(icos,issin) + cosomegat * ssinomegat
212219
CS%FtF(isin,iccos) = CS%FtF(isin,iccos) + sinomegat * ccosomegat
@@ -280,8 +287,8 @@ subroutine HA_accum_FtSSH(key, data, Time, G, CS)
280287
do c=1,nc
281288
icos = 2*c
282289
isin = 2*c+1
283-
cosomegat = cos(CS%freq(c) * now + CS%phase0(c))
284-
sinomegat = sin(CS%freq(c) * now + CS%phase0(c))
290+
cosomegat = CS%tide_fn(c) * cos(CS%freq(c) * now + CS%phase0(c) + CS%tide_un(c))
291+
sinomegat = CS%tide_fn(c) * sin(CS%freq(c) * now + CS%phase0(c) + CS%tide_un(c))
285292
do j=js,je ; do i=is,ie
286293
ha1%FtSSH(i,j,1) = ha1%FtSSH(i,j,1) + (data(i,j) - ha1%ref(i,j))
287294
ha1%FtSSH(i,j,icos) = ha1%FtSSH(i,j,icos) + (data(i,j) - ha1%ref(i,j)) * cosomegat

src/diagnostics/MOM_obsolete_params.F90

Lines changed: 5 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -155,6 +155,11 @@ subroutine find_obsolete_params(param_file)
155155
call obsolete_logical(param_file, "USE_GRID_SPACE_DIAGNOSTIC_AXES", &
156156
hint="Instead use USE_INDEX_DIAGNOSTIC_AXIS.")
157157

158+
call obsolete_logical(param_file, "OBC_TIDE_ADD_EQ_PHASE", hint="Instead use TIDE_USE_EQ_PHASE.")
159+
call obsolete_logical(param_file, "OBC_TIDE_ADD_NODAL", hint="Instead use TIDE_ADD_NODAL.")
160+
call obsolete_int(param_file, "OBC_TIDE_REF_DATE", hint="Instead use TIDE_REF_DATE.")
161+
call obsolete_int(param_file, "OBC_TIDE_NODAL_REF_DATE", hint="Instead use TIDE_NODAL_REF_DATE.")
162+
158163
! Write the file version number to the model log.
159164
call log_version(param_file, mdl, version)
160165

src/parameterizations/lateral/MOM_tidal_forcing.F90

Lines changed: 59 additions & 14 deletions
Original file line numberDiff line numberDiff line change
@@ -70,7 +70,9 @@ module MOM_tidal_forcing
7070
ampsal(:,:,:), & !< The amplitude of the SAL [Z ~> m].
7171
cosphase_prev(:,:,:), & !< The cosine of the phase of the amphidromes in the previous tidal solutions [nondim].
7272
sinphase_prev(:,:,:), & !< The sine of the phase of the amphidromes in the previous tidal solutions [nondim].
73-
amp_prev(:,:,:) !< The amplitude of the previous tidal solution [Z ~> m].
73+
amp_prev(:,:,:), & !< The amplitude of the previous tidal solution [Z ~> m].
74+
tide_fn(:), & !< Amplitude modulation of tides by nodal cycle [nondim].
75+
tide_un(:) !< Phase modulation of tides by nodal cycle [rad].
7476
end type tidal_forcing_CS
7577

7678
integer :: id_clock_tides !< CPU clock for tides
@@ -251,9 +253,14 @@ subroutine tidal_forcing_init(Time, G, US, param_file, CS, HA_CS)
251253
real, dimension(MAX_CONSTITUENTS) :: amp_def ! Default amplitude for each tidal constituent [m]
252254
real, dimension(MAX_CONSTITUENTS) :: love_def ! Default love number for each constituent [nondim]
253255
integer, dimension(3) :: tide_ref_date !< Reference date (t = 0) for tidal forcing.
256+
integer, dimension(3) :: nodal_ref_date !< Reference date for calculating nodal modulation for tidal forcing.
254257
logical :: use_M2, use_S2, use_N2, use_K2, use_K1, use_O1, use_P1, use_Q1
255258
logical :: use_MF, use_MM
256259
logical :: tides ! True if a tidal forcing is to be used.
260+
logical :: add_nodal_terms = .false. !< If true, insert terms for the 18.6 year modulation when
261+
!! calculating tidal forcing.
262+
type(time_type) :: nodal_time !< Model time to calculate nodal modulation for.
263+
type(astro_longitudes) :: nodal_longitudes !< Solar and lunar longitudes for tidal forcing
257264
logical :: HA_ssh, HA_ubt, HA_vbt
258265
! This include declares and sets the variable "version".
259266
# include "version_variable.h"
@@ -527,8 +534,46 @@ subroutine tidal_forcing_init(Time, G, US, param_file, CS, HA_CS)
527534
enddo
528535
endif
529536

537+
call get_param(param_file, mdl, "TIDE_ADD_NODAL", add_nodal_terms, &
538+
"If true, include 18.6 year nodal modulation in the astronomical tidal forcing.", &
539+
default=.false.)
540+
call get_param(param_file, mdl, "TIDE_NODAL_REF_DATE", nodal_ref_date, &
541+
"Fixed reference date to use for nodal modulation of astronomical tidal forcing.", &
542+
fail_if_missing=.false., default=0)
543+
544+
! If the nodal correction is based on a different time, initialize that.
545+
! Otherwise, it can use N from the time reference.
546+
if (add_nodal_terms) then
547+
if (sum(nodal_ref_date) /= 0) then
548+
! A reference date was provided for the nodal correction
549+
nodal_time = set_date(nodal_ref_date(1), nodal_ref_date(2), nodal_ref_date(3))
550+
call astro_longitudes_init(nodal_time, nodal_longitudes)
551+
elseif (CS%use_eq_phase) then
552+
! Astronomical longitudes were already calculated for use in equilibrium phases,
553+
! so use nodal longitude from that.
554+
nodal_longitudes = CS%tidal_longitudes
555+
else
556+
! Tidal reference time is a required parameter, so calculate the longitudes from that.
557+
call astro_longitudes_init(CS%time_ref, nodal_longitudes)
558+
endif
559+
endif
560+
561+
allocate(CS%tide_fn(nc))
562+
allocate(CS%tide_un(nc))
563+
564+
do c=1,nc
565+
! Find nodal corrections if needed
566+
if (add_nodal_terms) then
567+
call nodal_fu(trim(CS%const_name(c)), nodal_longitudes%N, CS%tide_fn(c), CS%tide_un(c))
568+
else
569+
CS%tide_fn(c) = 1.0
570+
CS%tide_un(c) = 0.0
571+
endif
572+
enddo
573+
530574
if (present(HA_CS)) then
531-
call HA_init(Time, US, param_file, CS%time_ref, CS%nc, CS%freq, CS%phase0, CS%const_name, HA_CS)
575+
call HA_init(Time, US, param_file, CS%time_ref, CS%nc, CS%freq, CS%phase0, CS%const_name, &
576+
CS%tide_fn, CS%tide_un, HA_CS)
532577
call get_param(param_file, mdl, "HA_SSH", HA_ssh, &
533578
"If true, perform harmonic analysis of sea serface height.", default=.false.)
534579
if (HA_ssh) call HA_register('ssh', 'h', HA_CS)
@@ -613,26 +658,26 @@ subroutine calc_tidal_forcing(Time, e_tide_eq, e_tide_sal, G, US, CS)
613658

614659
do c=1,CS%nc
615660
m = CS%struct(c)
616-
amp_cosomegat = CS%amp(c)*CS%love_no(c) * cos(CS%freq(c)*now + CS%phase0(c))
617-
amp_sinomegat = CS%amp(c)*CS%love_no(c) * sin(CS%freq(c)*now + CS%phase0(c))
661+
amp_cosomegat = CS%amp(c)*CS%love_no(c)*CS%tide_fn(c) * cos(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
662+
amp_sinomegat = CS%amp(c)*CS%love_no(c)*CS%tide_fn(c) * sin(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
618663
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
619664
e_tide_eq(i,j) = e_tide_eq(i,j) + (amp_cosomegat*CS%cos_struct(i,j,m) + &
620665
amp_sinomegat*CS%sin_struct(i,j,m))
621666
enddo ; enddo
622667
enddo
623668

624669
if (CS%use_tidal_sal_file) then ; do c=1,CS%nc
625-
cosomegat = cos(CS%freq(c)*now)
626-
sinomegat = sin(CS%freq(c)*now)
670+
cosomegat = CS%tide_fn(c) * cos(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
671+
sinomegat = CS%tide_fn(c) * sin(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
627672
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
628673
e_tide_sal(i,j) = e_tide_sal(i,j) + CS%ampsal(i,j,c) * &
629674
(cosomegat*CS%cosphasesal(i,j,c) + sinomegat*CS%sinphasesal(i,j,c))
630675
enddo ; enddo
631676
enddo ; endif
632677

633678
if (CS%use_tidal_sal_prev) then ; do c=1,CS%nc
634-
cosomegat = cos(CS%freq(c)*now)
635-
sinomegat = sin(CS%freq(c)*now)
679+
cosomegat = CS%tide_fn(c) * cos(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
680+
sinomegat = CS%tide_fn(c) * sin(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
636681
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
637682
e_tide_sal(i,j) = e_tide_sal(i,j) - CS%sal_scalar * CS%amp_prev(i,j,c) * &
638683
(cosomegat*CS%cosphase_prev(i,j,c) + sinomegat*CS%sinphase_prev(i,j,c))
@@ -691,8 +736,8 @@ subroutine calc_tidal_forcing_legacy(Time, e_sal, e_sal_tide, e_tide_eq, e_tide_
691736

692737
do c=1,CS%nc
693738
m = CS%struct(c)
694-
amp_cosomegat = CS%amp(c)*CS%love_no(c) * cos(CS%freq(c)*now + CS%phase0(c))
695-
amp_sinomegat = CS%amp(c)*CS%love_no(c) * sin(CS%freq(c)*now + CS%phase0(c))
739+
amp_cosomegat = CS%amp(c)*CS%love_no(c)*CS%tide_fn(c) * cos(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
740+
amp_sinomegat = CS%amp(c)*CS%love_no(c)*CS%tide_fn(c) * sin(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
696741
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
697742
amp_cossin = (amp_cosomegat*CS%cos_struct(i,j,m) + amp_sinomegat*CS%sin_struct(i,j,m))
698743
e_sal_tide(i,j) = e_sal_tide(i,j) + amp_cossin
@@ -701,8 +746,8 @@ subroutine calc_tidal_forcing_legacy(Time, e_sal, e_sal_tide, e_tide_eq, e_tide_
701746
enddo
702747

703748
if (CS%use_tidal_sal_file) then ; do c=1,CS%nc
704-
cosomegat = cos(CS%freq(c)*now)
705-
sinomegat = sin(CS%freq(c)*now)
749+
cosomegat = CS%tide_fn(c) * cos(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
750+
sinomegat = CS%tide_fn(c) * sin(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
706751
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
707752
amp_cossin = CS%ampsal(i,j,c) &
708753
* (cosomegat*CS%cosphasesal(i,j,c) + sinomegat*CS%sinphasesal(i,j,c))
@@ -712,8 +757,8 @@ subroutine calc_tidal_forcing_legacy(Time, e_sal, e_sal_tide, e_tide_eq, e_tide_
712757
enddo ; endif
713758

714759
if (CS%use_tidal_sal_prev) then ; do c=1,CS%nc
715-
cosomegat = cos(CS%freq(c)*now)
716-
sinomegat = sin(CS%freq(c)*now)
760+
cosomegat = CS%tide_fn(c) * cos(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
761+
sinomegat = CS%tide_fn(c) * sin(CS%freq(c)*now + CS%phase0(c) + CS%tide_un(c))
717762
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
718763
amp_cossin = -CS%sal_scalar * CS%amp_prev(i,j,c) &
719764
* (cosomegat*CS%cosphase_prev(i,j,c) + sinomegat*CS%sinphase_prev(i,j,c))

0 commit comments

Comments
 (0)