Skip to content

Commit b239629

Browse files
authored
Deprecate TIDE_SAL_SCALAR_VALUE (#819)
Using the recently added `old_name` option in `get_params`.
1 parent d6f3fa0 commit b239629

File tree

2 files changed

+19
-32
lines changed

2 files changed

+19
-32
lines changed

src/parameterizations/lateral/MOM_self_attr_load.F90

Lines changed: 9 additions & 16 deletions
Original file line numberDiff line numberDiff line change
@@ -190,10 +190,9 @@ subroutine SAL_init(G, GV, US, param_file, CS)
190190
real :: rhoE ! The average density of Earth [R ~> kg m-3].
191191
character(len=200) :: filename, ebot_ref_file, inputdir ! Strings for file/path
192192
character(len=200) :: ebot_ref_varname ! Variable name in file
193-
logical :: calculate_sal=.false.
194-
logical :: tides=.false., use_tidal_sal_file=.false., bq_sal_tides_bug=.false.
195-
integer :: tides_answer_date=99991231 ! Recover old answers with tides
196-
real :: sal_scalar_value, tide_sal_scalar_value ! Scaling SAL factors [nondim]
193+
logical :: calculate_sal, tides, use_tidal_sal_file
194+
integer :: tides_answer_date ! Recover old answers with tides
195+
real :: sal_scalar_value ! Scaling SAL factors [nondim]
197196
integer :: isd, ied, jsd, jed
198197

199198
isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
@@ -245,18 +244,12 @@ subroutine SAL_init(G, GV, US, param_file, CS)
245244
if (CS%use_sal_scalar .and. CS%use_bpa) &
246245
call MOM_error(WARNING, trim(mdl) // ", SAL_init: Using bottom pressure anomaly for scalar "//&
247246
"approximation SAL is unsubstantiated.")
248-
call get_param(param_file, '', "TIDE_SAL_SCALAR_VALUE", tide_sal_scalar_value, &
249-
units="m m-1", default=0.0, do_not_log=.True.)
250-
if (tide_sal_scalar_value/=0.0) &
251-
call MOM_error(WARNING, "TIDE_SAL_SCALAR_VALUE is a deprecated parameter. "//&
252-
"Use SAL_SCALAR_VALUE instead." )
253-
call get_param(param_file, mdl, "SAL_SCALAR_VALUE", sal_scalar_value, &
254-
"The constant of proportionality between sea surface "//&
255-
"height (really it should be bottom pressure) anomalies "//&
256-
"and bottom geopotential anomalies. This is only used if "//&
257-
"USE_SAL_SCALAR is true or USE_PREVIOUS_TIDES is true.", &
258-
default=tide_sal_scalar_value, units="m m-1", &
259-
do_not_log=(.not. CS%use_sal_scalar) .and. (.not. CS%use_tidal_sal_prev))
247+
call get_param(param_file, mdl, "SAL_SCALAR_VALUE", sal_scalar_value, "The constant of "//&
248+
"proportionality between self-attraction and loading (SAL) geopotential "//&
249+
"anomaly and barotropic geopotential anomaly. This is only used if "//&
250+
"SAL_SCALAR_APPROX is true or USE_PREVIOUS_TIDES is true.", default=0.0, &
251+
units="m m-1", do_not_log=.not.(CS%use_sal_scalar .or. CS%use_tidal_sal_prev), &
252+
old_name='TIDE_SAL_SCALAR_VALUE')
260253
call get_param(param_file, mdl, "SAL_HARMONICS", CS%use_sal_sht, &
261254
"If true, use the online spherical harmonics method to calculate "//&
262255
"self-attraction and loading.", default=.false.)

src/parameterizations/lateral/MOM_tidal_forcing.F90

Lines changed: 10 additions & 16 deletions
Original file line numberDiff line numberDiff line change
@@ -45,9 +45,9 @@ module MOM_tidal_forcing
4545
!! equilibrium tide. Set to false if providing tidal phases
4646
!! that have already been shifted by the
4747
!! astronomical/equilibrium argument.
48-
real :: sal_scalar !< The constant of proportionality between sea surface
49-
!! height (really it should be bottom pressure) anomalies
50-
!! and bottom geopotential anomalies [nondim].
48+
real :: sal_scalar = 0.0 !< The constant of proportionality between self-attraction and
49+
!! loading (SAL) geopotential anomaly and total geopotential geopotential
50+
!! anomalies. This is only used if USE_PREVIOUS_TIDES is true. [nondim].
5151
integer :: nc !< The number of tidal constituents in use.
5252
real, dimension(MAX_CONSTITUENTS) :: &
5353
freq, & !< The frequency of a tidal constituent [rad T-1 ~> rad s-1].
@@ -267,7 +267,6 @@ subroutine tidal_forcing_init(Time, G, US, param_file, CS, HA_CS)
267267
character(len=40) :: mdl = "MOM_tidal_forcing" ! This module's name.
268268
character(len=128) :: mesg
269269
character(len=200) :: tidal_input_files(4*MAX_CONSTITUENTS)
270-
real :: tide_sal_scalar_value ! The constant of proportionality with the scalar approximation to SAL [nondim]
271270
integer :: i, j, c, is, ie, js, je, isd, ied, jsd, jed, nc
272271

273272
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec
@@ -362,18 +361,13 @@ subroutine tidal_forcing_init(Time, G, US, param_file, CS, HA_CS)
362361
"If true, use the SAL from the previous iteration of the "//&
363362
"tides to facilitate convergent iteration. "//&
364363
"This is only used if TIDES is true.", default=.false.)
365-
call get_param(param_file, '', "TIDE_SAL_SCALAR_VALUE", tide_sal_scalar_value, &
366-
units="m m-1", default=0.0, do_not_log=.True.)
367-
if (tide_sal_scalar_value/=0.0) &
368-
call MOM_error(WARNING, "TIDE_SAL_SCALAR_VALUE is a deprecated parameter. "//&
369-
"Use SAL_SCALAR_VALUE instead." )
370-
call get_param(param_file, mdl, "SAL_SCALAR_VALUE", CS%sal_scalar, &
371-
"The constant of proportionality between sea surface "//&
372-
"height (really it should be bottom pressure) anomalies "//&
373-
"and bottom geopotential anomalies. This is only used if "//&
374-
"USE_SAL_SCALAR is true or USE_PREVIOUS_TIDES is true.", &
375-
default=tide_sal_scalar_value, units="m m-1", &
376-
do_not_log=(.not. CS%use_tidal_sal_prev))
364+
if (CS%use_tidal_sal_prev) &
365+
call get_param(param_file, mdl, "SAL_SCALAR_VALUE", CS%sal_scalar, "The constant of "//&
366+
"proportionality between self-attraction and loading (SAL) geopotential "//&
367+
"anomaly and barotropic geopotential anomalies. This is only used if "//&
368+
"SAL_SCALAR_APPROX is true or USE_PREVIOUS_TIDES is true.", default=0.0, &
369+
units="m m-1", do_not_log=(.not.CS%use_tidal_sal_prev), &
370+
old_name='TIDE_SAL_SCALAR_VALUE')
377371

378372
if (nc > MAX_CONSTITUENTS) then
379373
write(mesg,'("Increase MAX_CONSTITUENTS in MOM_tidal_forcing.F90 to at least",I3, &

0 commit comments

Comments
 (0)