MOM6
mom_tidal_forcing Module Reference

Detailed Description

Tidal contributions to geopotential.

Data Types

type  astro_longitudes
 Simple type to store astronomical longitudes used to calculate tidal phases. More...
 
type  tidal_forcing_cs
 The control structure for the MOM_tidal_forcing module. More...
 

Functions/Subroutines

subroutine, public astro_longitudes_init (time_ref, longitudes)
 Finds astronomical longitudes s, h, p, and N, the mean longitude of the moon, sun, lunar perigee, and ascending node, respectively, at the specified reference time time_ref. These formulas were obtained from Kowalik and Luick, "Modern Theory and Practice of Tide Analysis and Tidal Power", 2019 (their Equation I.71), which are based on Schureman, 1958. For simplicity, the time associated with time_ref should be at midnight. These formulas also only make sense if the calendar is gregorian.
 
real function, public eq_phase (constit, longitudes)
 Calculates the equilibrium phase argument for the given tidal constituent constit and the astronomical longitudes and the reference time. These formulas follow Table I.4 of Kowalik and Luick, "Modern Theory and Practice of Tide Analysis and Tidal Power", 2019.
 
real function, public tidal_frequency (constit)
 Looks up angular frequencies for the main tidal constituents. Values used here are from previous versions of MOM.
 
subroutine, public nodal_fu (constit, N, fn, un)
 Find amplitude (f) and phase (u) modulation of tidal constituents by the 18.6 year nodal cycle. Values here follow Table I.6 in Kowalik and Luick, "Modern Theory and Practice of Tide Analysis and Tidal Power", 2019.
 
subroutine, public tidal_forcing_init (Time, G, param_file, CS)
 This subroutine allocates space for the static variables used by this module. The metrics may be effectively 0, 1, or 2-D arrays, while fields like the background viscosities are 2-D arrays. ALLOC is a macro defined in MOM_memory.h for allocate or nothing with static memory. More...
 
subroutine find_in_files (filenames, varname, array, G)
 This subroutine finds a named variable in a list of files and reads its values into a domain-decomposed 2-d array. More...
 
subroutine, public tidal_forcing_sensitivity (G, CS, deta_tidal_deta)
 This subroutine calculates returns the partial derivative of the local geopotential height with the input sea surface height due to self-attraction and loading. More...
 
subroutine, public calc_tidal_forcing (Time, eta, eta_tidal, G, CS, deta_tidal_deta, m_to_Z)
 This subroutine calculates the geopotential anomalies that drive the tides, including self-attraction and loading. Optionally, it also returns the partial derivative of the local geopotential height with the input sea surface height. For now, eta and eta_tidal are both geopotential heights in depth units, but probably the input for eta should really be replaced with the column mass anomalies. More...
 
subroutine, public tidal_forcing_end (CS)
 This subroutine deallocates memory associated with the tidal forcing module. More...
 

Variables

integer, parameter max_constituents = 10
 The maximum number of tidal constituents that could be used.
 
integer id_clock_tides
 CPU clock for tides.
 

Function/Subroutine Documentation

◆ calc_tidal_forcing()

subroutine, public mom_tidal_forcing::calc_tidal_forcing ( type(time_type), intent(in)  Time,
real, dimension(szi_(g),szj_(g)), intent(in)  eta,
real, dimension(szi_(g),szj_(g)), intent(out)  eta_tidal,
type(ocean_grid_type), intent(in)  G,
type(tidal_forcing_cs), pointer  CS,
real, intent(out), optional  deta_tidal_deta,
real, intent(in), optional  m_to_Z 
)

This subroutine calculates the geopotential anomalies that drive the tides, including self-attraction and loading. Optionally, it also returns the partial derivative of the local geopotential height with the input sea surface height. For now, eta and eta_tidal are both geopotential heights in depth units, but probably the input for eta should really be replaced with the column mass anomalies.

Parameters
[in]gThe ocean's grid structure.
[in]timeThe time for the caluculation.
[in]etaThe sea surface height anomaly from a time-mean geoid [Z ~> m].
[out]eta_tidalThe tidal forcing geopotential height anomalies [Z ~> m].
csThe control structure returned by a previous call to tidal_forcing_init.
[out]deta_tidal_detaThe partial derivative of eta_tidal with the local value of eta [nondim].
[in]m_to_zA scaling factor from m to the units of eta.

Definition at line 583 of file MOM_tidal_forcing.F90.

583  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
584  type(time_type), intent(in) :: time !< The time for the caluculation.
585  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: eta !< The sea surface height anomaly from
586  !! a time-mean geoid [Z ~> m].
587  real, dimension(SZI_(G),SZJ_(G)), intent(out) :: eta_tidal !< The tidal forcing geopotential height
588  !! anomalies [Z ~> m].
589  type(tidal_forcing_cs), pointer :: cs !< The control structure returned by a
590  !! previous call to tidal_forcing_init.
591  real, optional, intent(out) :: deta_tidal_deta !< The partial derivative of
592  !! eta_tidal with the local value of
593  !! eta [nondim].
594  real, optional, intent(in) :: m_to_z !< A scaling factor from m to the units of eta.
595 
596  ! Local variables
597  real :: now ! The relative time in seconds.
598  real :: amp_cosomegat, amp_sinomegat
599  real :: cosomegat, sinomegat
600  real :: m_z ! A scaling factor from m to depth units.
601  real :: eta_prop ! The nondimenional constant of proportionality beteen eta and eta_tidal.
602  integer :: i, j, c, m, is, ie, js, je, isq, ieq, jsq, jeq
603  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
604  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
605 
606  if (.not.associated(cs)) return
607 
608  call cpu_clock_begin(id_clock_tides)
609 
610  if (cs%nc == 0) then
611  do j=jsq,jeq+1 ; do i=isq,ieq+1 ; eta_tidal(i,j) = 0.0 ; enddo ; enddo
612  return
613  endif
614 
615  now = time_type_to_real(time - cs%time_ref)
616 
617  if (cs%USE_SAL_SCALAR .and. cs%USE_PREV_TIDES) then
618  eta_prop = 2.0*cs%SAL_SCALAR
619  elseif (cs%USE_SAL_SCALAR .or. cs%USE_PREV_TIDES) then
620  eta_prop = cs%SAL_SCALAR
621  else
622  eta_prop = 0.0
623  endif
624 
625  if (present(deta_tidal_deta)) then
626  deta_tidal_deta = eta_prop
627  do j=jsq,jeq+1 ; do i=isq,ieq+1 ; eta_tidal(i,j) = 0.0 ; enddo ; enddo
628  else
629  do j=jsq,jeq+1 ; do i=isq,ieq+1
630  eta_tidal(i,j) = eta_prop*eta(i,j)
631  enddo ; enddo
632  endif
633 
634  m_z = 1.0 ; if (present(m_to_z)) m_z = m_to_z
635 
636  do c=1,cs%nc
637  m = cs%struct(c)
638  amp_cosomegat = m_z*cs%amp(c)*cs%love_no(c) * cos(cs%freq(c)*now + cs%phase0(c))
639  amp_sinomegat = m_z*cs%amp(c)*cs%love_no(c) * sin(cs%freq(c)*now + cs%phase0(c))
640  do j=jsq,jeq+1 ; do i=isq,ieq+1
641  eta_tidal(i,j) = eta_tidal(i,j) + (amp_cosomegat*cs%cos_struct(i,j,m) + &
642  amp_sinomegat*cs%sin_struct(i,j,m))
643  enddo ; enddo
644  enddo
645 
646  if (cs%tidal_sal_from_file) then ; do c=1,cs%nc
647  cosomegat = cos(cs%freq(c)*now)
648  sinomegat = sin(cs%freq(c)*now)
649  do j=jsq,jeq+1 ; do i=isq,ieq+1
650  eta_tidal(i,j) = eta_tidal(i,j) + m_z*cs%ampsal(i,j,c) * &
651  (cosomegat*cs%cosphasesal(i,j,c) + sinomegat*cs%sinphasesal(i,j,c))
652  enddo ; enddo
653  enddo ; endif
654 
655  if (cs%USE_PREV_TIDES) then ; do c=1,cs%nc
656  cosomegat = cos(cs%freq(c)*now)
657  sinomegat = sin(cs%freq(c)*now)
658  do j=jsq,jeq+1 ; do i=isq,ieq+1
659  eta_tidal(i,j) = eta_tidal(i,j) - m_z*cs%SAL_SCALAR*cs%amp_prev(i,j,c) * &
660  (cosomegat*cs%cosphase_prev(i,j,c) + sinomegat*cs%sinphase_prev(i,j,c))
661  enddo ; enddo
662  enddo ; endif
663 
664  call cpu_clock_end(id_clock_tides)
665 

◆ find_in_files()

subroutine mom_tidal_forcing::find_in_files ( character(len=*), dimension(:), intent(in)  filenames,
character(len=*), intent(in)  varname,
real, dimension(szi_(g),szj_(g)), intent(out)  array,
type(ocean_grid_type), intent(in)  G 
)
private

This subroutine finds a named variable in a list of files and reads its values into a domain-decomposed 2-d array.

Parameters
[in]filenamesThe names of the files to search for the named variable
[in]varnameThe name of the variable to read
[in]gThe ocean's grid structure
[out]arrayThe array to fill with the data

Definition at line 530 of file MOM_tidal_forcing.F90.

530  character(len=*), dimension(:), intent(in) :: filenames !< The names of the files to search for the named variable
531  character(len=*), intent(in) :: varname !< The name of the variable to read
532  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
533  real, dimension(SZI_(G),SZJ_(G)), intent(out) :: array !< The array to fill with the data
534  ! Local variables
535  integer :: nf
536 
537  do nf=1,size(filenames)
538  if (len_trim(filenames(nf)) == 0) cycle
539  if (field_exists(filenames(nf), varname, g%Domain%mpp_domain)) then
540  call mom_read_data(filenames(nf), varname, array, g%Domain)
541  return
542  endif
543  enddo
544 
545  do nf=size(filenames),1,-1
546  if (file_exists(filenames(nf), g%Domain)) then
547  call mom_error(fatal, "MOM_tidal_forcing.F90: Unable to find "// &
548  trim(varname)//" in any of the tidal input files, last tried "// &
549  trim(filenames(nf)))
550  endif
551  enddo
552 
553  call mom_error(fatal, "MOM_tidal_forcing.F90: Unable to find any of the "// &
554  "tidal input files, including "//trim(filenames(1)))
555 

◆ tidal_forcing_end()

subroutine, public mom_tidal_forcing::tidal_forcing_end ( type(tidal_forcing_cs), pointer  CS)

This subroutine deallocates memory associated with the tidal forcing module.

Parameters
csThe control structure returned by a previous call to tidal_forcing_init; it is deallocated here.

Definition at line 670 of file MOM_tidal_forcing.F90.

670  type(tidal_forcing_cs), pointer :: cs !< The control structure returned by a previous call
671  !! to tidal_forcing_init; it is deallocated here.
672 
673  if (associated(cs%sin_struct)) deallocate(cs%sin_struct)
674  if (associated(cs%cos_struct)) deallocate(cs%cos_struct)
675 
676  if (associated(cs%cosphasesal)) deallocate(cs%cosphasesal)
677  if (associated(cs%sinphasesal)) deallocate(cs%sinphasesal)
678  if (associated(cs%ampsal)) deallocate(cs%ampsal)
679 
680  if (associated(cs%cosphase_prev)) deallocate(cs%cosphase_prev)
681  if (associated(cs%sinphase_prev)) deallocate(cs%sinphase_prev)
682  if (associated(cs%amp_prev)) deallocate(cs%amp_prev)
683 
684  if (associated(cs)) deallocate(cs)
685 

◆ tidal_forcing_init()

subroutine, public mom_tidal_forcing::tidal_forcing_init ( type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(inout)  G,
type(param_file_type), intent(in)  param_file,
type(tidal_forcing_cs), pointer  CS 
)

This subroutine allocates space for the static variables used by this module. The metrics may be effectively 0, 1, or 2-D arrays, while fields like the background viscosities are 2-D arrays. ALLOC is a macro defined in MOM_memory.h for allocate or nothing with static memory.

Parameters
[in]timeThe current model time.
[in,out]gThe ocean's grid structure.
[in]param_fileA structure to parse for run-time parameters.
csA pointer that is set to point to the control structure for this module.

Definition at line 230 of file MOM_tidal_forcing.F90.

230  type(time_type), intent(in) :: time !< The current model time.
231  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure.
232  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters.
233  type(tidal_forcing_cs), pointer :: cs !< A pointer that is set to point to the control
234  !! structure for this module.
235  ! Local variables
236  real, dimension(SZI_(G), SZJ_(G)) :: &
237  phase, & ! The phase of some tidal constituent.
238  lat_rad, lon_rad ! Latitudes and longitudes of h-points in radians.
239  real :: deg_to_rad
240  real, dimension(MAX_CONSTITUENTS) :: freq_def, phase0_def, amp_def, love_def
241  integer, dimension(3) :: tide_ref_date !< Reference date (t = 0) for tidal forcing.
242  logical :: use_const ! True if a constituent is being used.
243  logical :: use_m2, use_s2, use_n2, use_k2, use_k1, use_o1, use_p1, use_q1
244  logical :: use_mf, use_mm
245  logical :: tides ! True if a tidal forcing is to be used.
246  logical :: fail_if_missing = .true.
247 ! This include declares and sets the variable "version".
248 #include "version_variable.h"
249  character(len=40) :: mdl = "MOM_tidal_forcing" ! This module's name.
250  character(len=128) :: mesg
251  character(len=200) :: tidal_input_files(4*max_constituents)
252  integer :: i, j, c, is, ie, js, je, isd, ied, jsd, jed, nc
253  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
254  isd = g%isd ; ied = g%ied ; jsd = g%jsd; jed = g%jed
255 
256  if (associated(cs)) then
257  call mom_error(warning, "tidal_forcing_init called with an associated "// &
258  "control structure.")
259  return
260  endif
261 
262  ! Read all relevant parameters and write them to the model log.
263  call log_version(param_file, mdl, version, "")
264  call get_param(param_file, mdl, "TIDES", tides, &
265  "If true, apply tidal momentum forcing.", default=.false.)
266 
267  if (.not.tides) return
268 
269  allocate(cs)
270 
271  ! Set up the spatial structure functions for the diurnal, semidiurnal, and
272  ! low-frequency tidal components.
273  allocate(cs%sin_struct(isd:ied,jsd:jed,3)) ; cs%sin_struct(:,:,:) = 0.0
274  allocate(cs%cos_struct(isd:ied,jsd:jed,3)) ; cs%cos_struct(:,:,:) = 0.0
275  deg_to_rad = 4.0*atan(1.0)/180.0
276  do j=js-1,je+1 ; do i=is-1,ie+1
277  lat_rad(i,j) = g%geoLatT(i,j)*deg_to_rad
278  lon_rad(i,j) = g%geoLonT(i,j)*deg_to_rad
279  enddo ; enddo
280  do j=js-1,je+1 ; do i=is-1,ie+1
281  cs%sin_struct(i,j,1) = -sin(2.0*lat_rad(i,j)) * sin(lon_rad(i,j))
282  cs%cos_struct(i,j,1) = sin(2.0*lat_rad(i,j)) * cos(lon_rad(i,j))
283  cs%sin_struct(i,j,2) = -cos(lat_rad(i,j))**2 * sin(2.0*lon_rad(i,j))
284  cs%cos_struct(i,j,2) = cos(lat_rad(i,j))**2 * cos(2.0*lon_rad(i,j))
285  cs%sin_struct(i,j,3) = 0.0
286  cs%cos_struct(i,j,3) = (0.5-1.5*sin(lat_rad(i,j))**2)
287  enddo ; enddo
288 
289  call get_param(param_file, mdl, "TIDE_M2", use_m2, &
290  "If true, apply tidal momentum forcing at the M2 "//&
291  "frequency. This is only used if TIDES is true.", &
292  default=.false.)
293  call get_param(param_file, mdl, "TIDE_S2", use_s2, &
294  "If true, apply tidal momentum forcing at the S2 "//&
295  "frequency. This is only used if TIDES is true.", &
296  default=.false.)
297  call get_param(param_file, mdl, "TIDE_N2", use_n2, &
298  "If true, apply tidal momentum forcing at the N2 "//&
299  "frequency. This is only used if TIDES is true.", &
300  default=.false.)
301  call get_param(param_file, mdl, "TIDE_K2", use_k2, &
302  "If true, apply tidal momentum forcing at the K2 "//&
303  "frequency. This is only used if TIDES is true.", &
304  default=.false.)
305  call get_param(param_file, mdl, "TIDE_K1", use_k1, &
306  "If true, apply tidal momentum forcing at the K1 "//&
307  "frequency. This is only used if TIDES is true.", &
308  default=.false.)
309  call get_param(param_file, mdl, "TIDE_O1", use_o1, &
310  "If true, apply tidal momentum forcing at the O1 "//&
311  "frequency. This is only used if TIDES is true.", &
312  default=.false.)
313  call get_param(param_file, mdl, "TIDE_P1", use_p1, &
314  "If true, apply tidal momentum forcing at the P1 "//&
315  "frequency. This is only used if TIDES is true.", &
316  default=.false.)
317  call get_param(param_file, mdl, "TIDE_Q1", use_q1, &
318  "If true, apply tidal momentum forcing at the Q1 "//&
319  "frequency. This is only used if TIDES is true.", &
320  default=.false.)
321  call get_param(param_file, mdl, "TIDE_MF", use_mf, &
322  "If true, apply tidal momentum forcing at the MF "//&
323  "frequency. This is only used if TIDES is true.", &
324  default=.false.)
325  call get_param(param_file, mdl, "TIDE_MM", use_mm, &
326  "If true, apply tidal momentum forcing at the MM "//&
327  "frequency. This is only used if TIDES is true.", &
328  default=.false.)
329 
330  ! Determine how many tidal components are to be used.
331  nc = 0
332  if (use_m2) nc=nc+1 ; if (use_s2) nc=nc+1
333  if (use_n2) nc=nc+1 ; if (use_k2) nc=nc+1
334  if (use_k1) nc=nc+1 ; if (use_o1) nc=nc+1
335  if (use_p1) nc=nc+1 ; if (use_q1) nc=nc+1
336  if (use_mf) nc=nc+1 ; if (use_mm) nc=nc+1
337  cs%nc = nc
338 
339  if (nc == 0) then
340  call mom_error(fatal, "tidal_forcing_init: "// &
341  "TIDES are defined, but no tidal constituents are used.")
342  return
343  endif
344 
345  call get_param(param_file, mdl, "TIDAL_SAL_FROM_FILE", cs%tidal_sal_from_file, &
346  "If true, read the tidal self-attraction and loading "//&
347  "from input files, specified by TIDAL_INPUT_FILE. "//&
348  "This is only used if TIDES is true.", default=.false.)
349  call get_param(param_file, mdl, "USE_PREVIOUS_TIDES", cs%use_prev_tides, &
350  "If true, use the SAL from the previous iteration of the "//&
351  "tides to facilitate convergent iteration. "//&
352  "This is only used if TIDES is true.", default=.false.)
353  call get_param(param_file, mdl, "TIDE_USE_SAL_SCALAR", cs%use_sal_scalar, &
354  "If true and TIDES is true, use the scalar approximation "//&
355  "when calculating self-attraction and loading.", &
356  default=.not.cs%tidal_sal_from_file)
357  ! If it is being used, sal_scalar MUST be specified in param_file.
358  if (cs%use_sal_scalar .or. cs%use_prev_tides) &
359  call get_param(param_file, mdl, "TIDE_SAL_SCALAR_VALUE", cs%sal_scalar, &
360  "The constant of proportionality between sea surface "//&
361  "height (really it should be bottom pressure) anomalies "//&
362  "and bottom geopotential anomalies. This is only used if "//&
363  "TIDES and TIDE_USE_SAL_SCALAR are true.", units="m m-1", &
364  fail_if_missing=.true.)
365 
366  if (nc > max_constituents) then
367  write(mesg,'("Increase MAX_CONSTITUENTS in MOM_tidal_forcing.F90 to at least",I3, &
368  &"to accommodate all the registered tidal constituents.")') nc
369  call mom_error(fatal, "MOM_tidal_forcing"//mesg)
370  endif
371 
372  do c=1,4*max_constituents ; tidal_input_files(c) = "" ; enddo
373 
374  if (cs%tidal_sal_from_file .or. cs%use_prev_tides) then
375  call get_param(param_file, mdl, "TIDAL_INPUT_FILE", tidal_input_files, &
376  "A list of input files for tidal information.", &
377  default = "", fail_if_missing=.true.)
378  endif
379 
380  call get_param(param_file, mdl, "TIDE_REF_DATE", tide_ref_date, &
381  "Year,month,day to use as reference date for tidal forcing. "//&
382  "If not specified, defaults to 0.", &
383  default=0)
384 
385  call get_param(param_file, mdl, "TIDE_USE_EQ_PHASE", cs%use_eq_phase, &
386  "Correct phases by calculating equilibrium phase arguments for TIDE_REF_DATE. ", &
387  default=.false., fail_if_missing=.false.)
388 
389  if (sum(tide_ref_date) == 0) then ! tide_ref_date defaults to 0.
390  cs%time_ref = set_date(1, 1, 1)
391  else
392  if(.not. cs%use_eq_phase) then
393  ! Using a reference date but not using phase relative to equilibrium.
394  ! This makes sense as long as either phases are overridden, or
395  ! correctly simulating tidal phases is not desired.
396  call mom_mesg('Tidal phases will *not* be corrected with equilibrium arguments.')
397  endif
398  cs%time_ref = set_date(tide_ref_date(1), tide_ref_date(2), tide_ref_date(3))
399  endif
400  ! Set the parameters for all components that are in use.
401  ! Initialize reference time for tides and
402  ! find relevant lunar and solar longitudes at the reference time.
403  if (cs%use_eq_phase) call astro_longitudes_init(cs%time_ref, cs%tidal_longitudes)
404  c=0
405  if (use_m2) then
406  c=c+1 ; cs%const_name(c) = "M2" ; cs%struct(c) = 2
407  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.242334
408  endif
409 
410  if (use_s2) then
411  c=c+1 ; cs%const_name(c) = "S2" ; cs%struct(c) = 2
412  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.112743
413  endif
414 
415  if (use_n2) then
416  c=c+1 ; cs%const_name(c) = "N2" ; cs%struct(c) = 2
417  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.046397
418  endif
419 
420  if (use_k2) then
421  c=c+1 ; cs%const_name(c) = "K2" ; cs%struct(c) = 2
422  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.030684
423  endif
424 
425  if (use_k1) then
426  c=c+1 ; cs%const_name(c) = "K1" ; cs%struct(c) = 1
427  cs%love_no(c) = 0.736 ; cs%amp(c) = 0.141565
428  endif
429 
430  if (use_o1) then
431  c=c+1 ; cs%const_name(c) = "O1" ; cs%struct(c) = 1
432  cs%love_no(c) = 0.695 ; cs%amp(c) = 0.100661
433  endif
434 
435  if (use_p1) then
436  c=c+1 ; cs%const_name(c) = "P1" ; cs%struct(c) = 1
437  cs%love_no(c) = 0.706 ; cs%amp(c) = 0.046848
438  endif
439 
440  if (use_q1) then
441  c=c+1 ; cs%const_name(c) = "Q1" ; cs%struct(c) = 1
442  cs%love_no(c) = 0.695 ; cs%amp(c) = 0.019273
443  endif
444 
445  if (use_mf) then
446  c=c+1 ; cs%const_name(c) = "MF" ; cs%struct(c) = 3
447  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.042041
448  endif
449 
450  if (use_mm) then
451  c=c+1 ; cs%const_name(c) = "MM" ; cs%struct(c) = 3
452  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.022191
453  endif
454 
455  ! Set defaults for all included constituents
456  ! and things that can be set by functions
457  do c=1,nc
458  cs%freq(c) = tidal_frequency(cs%const_name(c))
459  freq_def(c) = cs%freq(c)
460  love_def(c) = cs%love_no(c)
461  amp_def(c) = cs%amp(c)
462  cs%phase0(c) = 0.0
463  if (cs%use_eq_phase) then
464  phase0_def(c) = eq_phase(cs%const_name(c), cs%tidal_longitudes)
465  else
466  phase0_def(c) = 0.0
467  endif
468  enddo
469 
470  ! Parse the input file to potentially override the default values for the
471  ! frequency, amplitude and initial phase of each constituent, and log the
472  ! values that are actually used.
473  do c=1,nc
474  call get_param(param_file, mdl, "TIDE_"//trim(cs%const_name(c))//"_FREQ", cs%freq(c), &
475  "Frequency of the "//trim(cs%const_name(c))//" tidal constituent. "//&
476  "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
477  " are true, or if OBC_TIDE_N_CONSTITUENTS > 0 and "//trim(cs%const_name(c))// &
478  " is in OBC_TIDE_CONSTITUENTS.", units="s-1", default=freq_def(c))
479  call get_param(param_file, mdl, "TIDE_"//trim(cs%const_name(c))//"_AMP", cs%amp(c), &
480  "Amplitude of the "//trim(cs%const_name(c))//" tidal constituent. "//&
481  "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
482  " are true.", units="m", default=amp_def(c))
483  call get_param(param_file, mdl, "TIDE_"//trim(cs%const_name(c))//"_PHASE_T0", cs%phase0(c), &
484  "Phase of the "//trim(cs%const_name(c))//" tidal constituent at time 0. "//&
485  "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
486  " are true.", units="radians", default=phase0_def(c))
487  enddo
488 
489  if (cs%tidal_sal_from_file) then
490  allocate(cs%cosphasesal(isd:ied,jsd:jed,nc))
491  allocate(cs%sinphasesal(isd:ied,jsd:jed,nc))
492  allocate(cs%ampsal(isd:ied,jsd:jed,nc))
493  do c=1,nc
494  ! Read variables with names like PHASE_SAL_M2 and AMP_SAL_M2.
495  call find_in_files(tidal_input_files,"PHASE_SAL_"//trim(cs%const_name(c)),phase,g)
496  call find_in_files(tidal_input_files,"AMP_SAL_"//trim(cs%const_name(c)),cs%ampsal(:,:,c),g)
497  call pass_var(phase, g%domain,complete=.false.)
498  call pass_var(cs%ampsal(:,:,c),g%domain,complete=.true.)
499  do j=js-1,je+1 ; do i=is-1,ie+1
500  cs%cosphasesal(i,j,c) = cos(phase(i,j)*deg_to_rad)
501  cs%sinphasesal(i,j,c) = sin(phase(i,j)*deg_to_rad)
502  enddo ; enddo
503  enddo
504  endif
505 
506  if (cs%USE_PREV_TIDES) then
507  allocate(cs%cosphase_prev(isd:ied,jsd:jed,nc))
508  allocate(cs%sinphase_prev(isd:ied,jsd:jed,nc))
509  allocate(cs%amp_prev(isd:ied,jsd:jed,nc))
510  do c=1,nc
511  ! Read variables with names like PHASE_PREV_M2 and AMP_PREV_M2.
512  call find_in_files(tidal_input_files,"PHASE_PREV_"//trim(cs%const_name(c)),phase,g)
513  call find_in_files(tidal_input_files,"AMP_PREV_"//trim(cs%const_name(c)),cs%amp_prev(:,:,c),g)
514  call pass_var(phase, g%domain,complete=.false.)
515  call pass_var(cs%amp_prev(:,:,c),g%domain,complete=.true.)
516  do j=js-1,je+1 ; do i=is-1,ie+1
517  cs%cosphase_prev(i,j,c) = cos(phase(i,j)*deg_to_rad)
518  cs%sinphase_prev(i,j,c) = sin(phase(i,j)*deg_to_rad)
519  enddo ; enddo
520  enddo
521  endif
522 
523  id_clock_tides = cpu_clock_id('(Ocean tides)', grain=clock_module)
524 

◆ tidal_forcing_sensitivity()

subroutine, public mom_tidal_forcing::tidal_forcing_sensitivity ( type(ocean_grid_type), intent(in)  G,
type(tidal_forcing_cs), pointer  CS,
real, intent(out)  deta_tidal_deta 
)

This subroutine calculates returns the partial derivative of the local geopotential height with the input sea surface height due to self-attraction and loading.

Parameters
[in]gThe ocean's grid structure.
csThe control structure returned by a previous call to tidal_forcing_init.
[out]deta_tidal_detaThe partial derivative of eta_tidal with the local value of eta [nondim].

Definition at line 562 of file MOM_tidal_forcing.F90.

562  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
563  type(tidal_forcing_cs), pointer :: cs !< The control structure returned by a previous call to tidal_forcing_init.
564  real, intent(out) :: deta_tidal_deta !< The partial derivative of eta_tidal with
565  !! the local value of eta [nondim].
566 
567  if (cs%USE_SAL_SCALAR .and. cs%USE_PREV_TIDES) then
568  deta_tidal_deta = 2.0*cs%SAL_SCALAR
569  elseif (cs%USE_SAL_SCALAR .or. cs%USE_PREV_TIDES) then
570  deta_tidal_deta = cs%SAL_SCALAR
571  else
572  deta_tidal_deta = 0.0
573  endif