Skip to content

Commit d2a33ee

Browse files
committed
Inline Harmonic Analysis
Update 4: time_ref and const_name are defined in HA_init, instead of being copied from MOM_tidal_forcing. This commit prepares for the separation of MOM_harmonic_analysis and MOM_tidal_forcing.
1 parent b3aa2bc commit d2a33ee

File tree

2 files changed

+44
-8
lines changed

2 files changed

+44
-8
lines changed

src/diagnostics/MOM_harmonic_analysis.F90

Lines changed: 43 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -1,7 +1,8 @@
11
!> Inline harmonic analysis (conventional)
22
module MOM_harmonic_analysis
33

4-
use MOM_time_manager, only : time_type, real_to_time, time_type_to_real, get_date, increment_date, &
4+
use MOM_time_manager, only : time_type, real_to_time, time_type_to_real, &
5+
set_date, get_date, increment_date, &
56
operator(+), operator(-), operator(<), operator(>), operator(>=)
67
use MOM_grid, only : ocean_grid_type
78
use MOM_unit_scaling, only : unit_scale_type
@@ -52,7 +53,7 @@ module MOM_harmonic_analysis
5253
tide_un !< Phase modulation of tides by nodal cycle [rad].
5354
integer :: nc !< The number of tidal constituents in use
5455
integer :: length !< Number of fields of which harmonic analysis is to be performed
55-
character(len=16) :: const_name(MAX_CONSTITUENTS) !< The name of each constituent
56+
character(len=2), allocatable, dimension(:) :: const_name !< The name of each constituent
5657
character(len=255) :: path !< Path to directory where output will be written
5758
type(unit_scale_type) :: US !< A dimensional unit scaling type
5859
type(HA_node), pointer :: list => NULL() !< A linked list for storing the HA info of different fields
@@ -62,20 +63,26 @@ module MOM_harmonic_analysis
6263

6364
!> This subroutine sets static variables used by this module and initializes CS%list.
6465
!! THIS MUST BE CALLED AT THE END OF tidal_forcing_init.
65-
subroutine HA_init(Time, US, param_file, time_ref, nc, freq, phase0, const_name, tide_fn, tide_un, CS)
66+
subroutine HA_init(Time, US, param_file, nc, freq, phase0, tide_fn, tide_un, CS)
6667
type(time_type), intent(in) :: Time !< The current model time
67-
type(time_type), intent(in) :: time_ref !< Reference time (t = 0) used to calculate tidal forcing
6868
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
6969
type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
7070
real, intent(in) :: freq(MAX_CONSTITUENTS) !< The frequency of a tidal constituent [T-1 ~> s-1]
7171
real, intent(in) :: phase0(MAX_CONSTITUENTS) !< The phase of a tidal constituent at time 0 [rad]
7272
real, intent(in) :: tide_fn(MAX_CONSTITUENTS) !< Amplitude modulation of tides by nodal cycle [nondim].
7373
real, intent(in) :: tide_un(MAX_CONSTITUENTS) !< Phase modulation of tides by nodal cycle [rad].
7474
integer, intent(in) :: nc !< The number of tidal constituents in use
75-
character(len=16), intent(in) :: const_name(MAX_CONSTITUENTS) !< The name of each constituent
7675
type(harmonic_analysis_CS), intent(out) :: CS !< Control structure of the MOM_harmonic_analysis module
7776

7877
! Local variables
78+
logical :: tides !< True if tidal forcing module is enabled
79+
logical :: use_eq_phase !< If true, tidal forcing is phase-shifted to match
80+
!! equilibrium tide. Set to false if providing tidal phases
81+
!! that have already been shifted by the
82+
!! astronomical/equilibrium argument.
83+
integer, dimension(3) :: tide_ref_date !< Reference date (t = 0) for tidal forcing (year, month, day).
84+
character(len=2) :: const_name(MAX_CONSTITUENTS) !< The name of each constituent
85+
7986
type(HA_type) :: ha1 !< A temporary, null field used for initializing CS%list
8087
real :: HA_start_time !< Start time of harmonic analysis [T ~> s]
8188
real :: HA_end_time !< End time of harmonic analysis [T ~> s]
@@ -84,6 +91,37 @@ subroutine HA_init(Time, US, param_file, time_ref, nc, freq, phase0, const_name,
8491
character(len=255) :: mesg
8592
integer :: year, month, day, hour, minute, second
8693

94+
call get_param(param_file, mdl, "TIDES", tides, &
95+
"If true, apply tidal momentum forcing.", default=.false., do_not_log=.true.)
96+
97+
call get_param(param_file, mdl, "TIDE_REF_DATE", tide_ref_date, &
98+
"Reference date to use for tidal calculations and equilibrium phase.", &
99+
old_name="OBC_TIDE_REF_DATE", defaults=(/0, 0, 0/), do_not_log=tides)
100+
101+
call get_param(param_file, mdl, "TIDE_USE_EQ_PHASE", use_eq_phase, &
102+
"If true, add the equilibrium phase argument to the specified tidal phases.", &
103+
old_name="OBC_TIDE_ADD_EQ_PHASE", default=.false., do_not_log=tides)
104+
105+
call get_param(param_file, mdl, "HA_CONSTITUENTS", const_name, &
106+
"Names of tidal constituents to be harmonically analyzed. "//&
107+
"They don't have to be the same as those used in MOM_tidal_forcing.", &
108+
fail_if_missing=.true.)
109+
110+
if (sum(tide_ref_date) == 0) then ! tide_ref_date defaults to 0.
111+
CS%time_ref = set_date(1, 1, 1, 0, 0, 0)
112+
else
113+
if (.not. use_eq_phase) then
114+
! Using a reference date but not using phase relative to equilibrium.
115+
! This makes sense as long as either phases are overridden, or
116+
! correctly simulating tidal phases is not desired.
117+
call MOM_mesg('Tidal phases will *not* be corrected with equilibrium arguments.')
118+
endif
119+
CS%time_ref = set_date(tide_ref_date(1), tide_ref_date(2), tide_ref_date(3), 0, 0, 0)
120+
endif
121+
122+
allocate(CS%const_name(nc))
123+
read(const_name, *) CS%const_name
124+
87125
! Determine CS%time_start and CS%time_end
88126
call get_param(param_file, mdl, "HA_START_TIME", HA_start_time, &
89127
"Start time of harmonic analysis, in units of days after "//&
@@ -137,13 +175,11 @@ subroutine HA_init(Time, US, param_file, time_ref, nc, freq, phase0, const_name,
137175
"Path to output files for runtime harmonic analysis.", default="./")
138176

139177
! Populate some parameters of the control structure
140-
CS%time_ref = time_ref
141178
CS%freq = freq
142179
CS%phase0 = phase0
143180
CS%tide_fn = tide_fn
144181
CS%tide_un = tide_un
145182
CS%nc = nc
146-
CS%const_name = const_name
147183
CS%length = 0
148184
CS%US = US
149185

src/parameterizations/lateral/MOM_tidal_forcing.F90

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -566,7 +566,7 @@ subroutine tidal_forcing_init(Time, G, US, param_file, CS, HA_CS)
566566
enddo
567567

568568
if (present(HA_CS)) then
569-
call HA_init(Time, US, param_file, CS%time_ref, CS%nc, CS%freq, CS%phase0, CS%const_name, &
569+
call HA_init(Time, US, param_file, CS%nc, CS%freq, CS%phase0, &
570570
CS%tide_fn, CS%tide_un, HA_CS)
571571
endif
572572

0 commit comments

Comments
 (0)