MOM6
mom_tidal_forcing Module Reference

Data Types

type  tidal_forcing_cs
 

Functions/Subroutines

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 (tidal_input_files, varname, array, G)
 
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)
 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 m, but probably the input for eta should really be replaced with the column mass anomalies. More...
 
subroutine, public tidal_forcing_end (CS)
 

Variables

integer, parameter max_constituents = 10
 
integer id_clock_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 
)

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 m, 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 in m.
[out]eta_tidalThe tidal forcing geopotential anomalies, in 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.

Definition at line 467 of file MOM_tidal_forcing.F90.

References id_clock_tides.

Referenced by mom_pressureforce_afv::pressureforce_afv_bouss(), mom_pressureforce_afv::pressureforce_afv_nonbouss(), mom_pressureforce_mont::pressureforce_mont_bouss(), and mom_pressureforce_mont::pressureforce_mont_nonbouss().

467  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
468  type(time_type), intent(in) :: time !< The time for the caluculation.
469  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: eta !< The sea surface height anomaly from
470  !! a time-mean geoid in m.
471  real, dimension(SZI_(G),SZJ_(G)), intent(out) :: eta_tidal !< The tidal forcing geopotential
472  !! anomalies, in m.
473  type(tidal_forcing_cs), pointer :: cs !< The control structure returned by a
474  !! previous call to tidal_forcing_init.
475  real, optional, intent(out) :: deta_tidal_deta !< The partial derivative of
476  !! eta_tidal with the local value of
477  !! eta, nondim.
478 
479 ! This subroutine calculates the geopotential anomalies that drive the tides,
480 ! including self-attraction and loading. Optionally, it also returns the
481 ! partial derivative of the local geopotential height with the input sea surface
482 ! height. For now, eta and eta_tidal are both geopotential heights in m, but
483 ! probably the input for eta should really be replaced with the column mass
484 ! anomalies.
485 !
486 ! Arguments: Time - The time for the caluculation.
487 ! (in) eta - The sea surface height anomaly from a time-mean geoid in m.
488 ! (out) eta_tidal - The tidal forcing geopotential anomalies, in m.
489 ! (in) G - The ocean's grid structure.
490 ! (in) CS - The control structure returned by a previous call to
491 ! tidal_forcing_init.
492 ! (out) deta_tidal_deta - the partial derivative of eta_tidal with the
493 ! local value of eta, nondim.
494 
495  real :: eta_astro(szi_(g),szj_(g))
496  real :: eta_sal(szi_(g),szj_(g))
497  real :: now ! The relative time in seconds.
498  real :: amp_cosomegat, amp_sinomegat
499  real :: cosomegat, sinomegat
500  real :: eta_prop
501  integer :: i, j, c, m, is, ie, js, je, isq, ieq, jsq, jeq
502  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
503  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
504 
505  if (.not.associated(cs)) return
506 
507  call cpu_clock_begin(id_clock_tides)
508 
509  if (cs%nc == 0) then
510  do j=jsq,jeq+1 ; do i=isq,ieq+1 ; eta_tidal(i,j) = 0.0 ; enddo ; enddo
511  return
512  endif
513 
514  now = time_type_to_real(time)
515 
516  if (cs%USE_SAL_SCALAR .and. cs%USE_PREV_TIDES) then
517  eta_prop = 2.0*cs%SAL_SCALAR
518  elseif (cs%USE_SAL_SCALAR .or. cs%USE_PREV_TIDES) then
519  eta_prop = cs%SAL_SCALAR
520  else
521  eta_prop = 0.0
522  endif
523 
524  if (present(deta_tidal_deta)) then
525  deta_tidal_deta = eta_prop
526  do j=jsq,jeq+1 ; do i=isq,ieq+1 ; eta_tidal(i,j) = 0.0 ; enddo ; enddo
527  else
528  do j=jsq,jeq+1 ; do i=isq,ieq+1
529  eta_tidal(i,j) = eta_prop*eta(i,j)
530  enddo ; enddo
531  endif
532 
533  do c=1,cs%nc
534  m = cs%struct(c)
535  amp_cosomegat = cs%amp(c)*cs%love_no(c)*cos(cs%freq(c)*now + cs%phase0(c))
536  amp_sinomegat = cs%amp(c)*cs%love_no(c)*sin(cs%freq(c)*now + cs%phase0(c))
537  do j=jsq,jeq+1 ; do i=isq,ieq+1
538  eta_tidal(i,j) = eta_tidal(i,j) + (amp_cosomegat*cs%cos_struct(i,j,m) + &
539  amp_sinomegat*cs%sin_struct(i,j,m))
540  enddo ; enddo
541  enddo
542 
543  if (cs%tidal_sal_from_file) then ; do c=1,cs%nc
544  cosomegat = cos(cs%freq(c)*now)
545  sinomegat = sin(cs%freq(c)*now)
546  do j=jsq,jeq+1 ; do i=isq,ieq+1
547  eta_tidal(i,j) = eta_tidal(i,j) + cs%ampsal(i,j,c) * &
548  (cosomegat*cs%cosphasesal(i,j,c) + sinomegat*cs%sinphasesal(i,j,c))
549  enddo ; enddo
550  enddo ; endif
551 
552  if (cs%USE_PREV_TIDES) then ; do c=1,cs%nc
553  cosomegat = cos(cs%freq(c)*now)
554  sinomegat = sin(cs%freq(c)*now)
555  do j=jsq,jeq+1 ; do i=isq,ieq+1
556  eta_tidal(i,j) = eta_tidal(i,j) - cs%SAL_SCALAR*cs%amp_prev(i,j,c) * &
557  (cosomegat*cs%cosphase_prev(i,j,c)+sinomegat*cs%sinphase_prev(i,j,c))
558  enddo ; enddo
559  enddo ; endif
560 
561  call cpu_clock_end(id_clock_tides)
562 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the caller graph for this function:

◆ find_in_files()

subroutine mom_tidal_forcing::find_in_files ( character(len=*), dimension(:), intent(in)  tidal_input_files,
character(len=*), intent(in)  varname,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  array,
type(ocean_grid_type), intent(in)  G 
)
private
Parameters
[in]gThe ocean's grid structure

Definition at line 405 of file MOM_tidal_forcing.F90.

References mom_error_handler::mom_error().

Referenced by tidal_forcing_init().

405  character(len=*), intent(in) :: tidal_input_files(:)
406  character(len=*), intent(in) :: varname
407  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
408  real, dimension(SZI_(G),SZJ_(G)), intent(out) :: array
409 
410  integer :: nf
411 
412  do nf=1,size(tidal_input_files)
413  if (len_trim(tidal_input_files(nf)) == 0) cycle
414  if (field_exists(tidal_input_files(nf), varname, g%Domain%mpp_domain)) then
415  call read_data(tidal_input_files(nf), varname, array, g%Domain%mpp_domain)
416  return
417  endif
418  enddo
419 
420  do nf=size(tidal_input_files),1,-1
421  if (file_exists(tidal_input_files(nf), g%Domain)) then
422  call mom_error(fatal, "MOM_tidal_forcing.F90: Unable to find "// &
423  trim(varname)//" in any of the tidal input files, last tried "// &
424  trim(tidal_input_files(nf)))
425  endif
426  enddo
427 
428  call mom_error(fatal, "MOM_tidal_forcing.F90: Unable to find any of the "// &
429  "tidal input files, including "//trim(tidal_input_files(1)))
430 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ tidal_forcing_end()

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

Definition at line 566 of file MOM_tidal_forcing.F90.

566  type(tidal_forcing_cs), pointer :: cs
567 
568  if (associated(cs%sin_struct)) deallocate(cs%sin_struct)
569  if (associated(cs%cos_struct)) deallocate(cs%cos_struct)
570 
571  if (associated(cs%cosphasesal)) deallocate(cs%cosphasesal)
572  if (associated(cs%sinphasesal)) deallocate(cs%sinphasesal)
573  if (associated(cs%ampsal)) deallocate(cs%ampsal)
574 
575  if (associated(cs%cosphase_prev)) deallocate(cs%cosphase_prev)
576  if (associated(cs%sinphase_prev)) deallocate(cs%sinphase_prev)
577  if (associated(cs%amp_prev)) deallocate(cs%amp_prev)
578 
579  if (associated(cs)) deallocate(cs)
580 

◆ 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 113 of file MOM_tidal_forcing.F90.

References find_in_files(), id_clock_tides, max_constituents, and mom_error_handler::mom_error().

113  type(time_type), intent(in) :: time !< The current model time.
114  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure.
115  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
116  !! parameters.
117  type(tidal_forcing_cs), pointer :: cs !< A pointer that is set to point to the control
118  !! structure for this module.
119 
120 ! This subroutine allocates space for the static variables used
121 ! by this module. The metrics may be effectively 0, 1, or 2-D arrays,
122 ! while fields like the background viscosities are 2-D arrays.
123 ! ALLOC is a macro defined in MOM_memory.h for allocate or nothing with
124 ! static memory.
125 !
126 ! Arguments: Time - The current model time.
127 ! (in) G - The ocean's grid structure.
128 ! (in) param_file - A structure indicating the open file to parse for
129 ! model parameter values.
130 ! (in/out) CS - A pointer that is set to point to the control structure
131 ! for this module.
132  real, dimension(SZI_(G), SZJ_(G)) :: &
133  phase, & ! The phase of some tidal constituent.
134  lat_rad, lon_rad ! Latitudes and longitudes of h-points in radians.
135  real :: deg_to_rad
136  real, dimension(MAX_CONSTITUENTS) :: freq_def, phase0_def, amp_def, love_def
137  logical :: use_const ! True if a constituent is being used.
138  logical :: use_m2, use_s2, use_n2, use_k2, use_k1, use_o1, use_p1, use_q1
139  logical :: use_mf, use_mm
140  logical :: tides ! True if a tidal forcing is to be used.
141  logical :: fail_if_missing = .true.
142 ! This include declares and sets the variable "version".
143 #include "version_variable.h"
144  character(len=40) :: mdl = "MOM_tidal_forcing" ! This module's name.
145  character(len=128) :: mesg
146  character(len=200) :: tidal_input_files(4*max_constituents)
147  integer :: i, j, c, is, ie, js, je, isd, ied, jsd, jed, nc
148  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
149  isd = g%isd ; ied = g%ied ; jsd = g%jsd; jed = g%jed
150 
151  if (associated(cs)) then
152  call mom_error(warning, "tidal_forcing_init called with an associated "// &
153  "control structure.")
154  return
155  endif
156 
157  ! Read all relevant parameters and write them to the model log.
158  call log_version(param_file, mdl, version, "")
159  call get_param(param_file, mdl, "TIDES", tides, &
160  "If true, apply tidal momentum forcing.", default=.false.)
161 
162  if (.not.tides) return
163 
164  allocate(cs)
165 
166  ! Set up the spatial structure functions for the diurnal, semidiurnal, and
167  ! low-frequency tidal components.
168  allocate(cs%sin_struct(isd:ied,jsd:jed,3)) ; cs%sin_struct(:,:,:) = 0.0
169  allocate(cs%cos_struct(isd:ied,jsd:jed,3)) ; cs%cos_struct(:,:,:) = 0.0
170  deg_to_rad = 4.0*atan(1.0)/180.0
171  do j=js-1,je+1 ; do i=is-1,ie+1
172  lat_rad(i,j) = g%geoLatT(i,j)*deg_to_rad
173  lon_rad(i,j) = g%geoLonT(i,j)*deg_to_rad
174  enddo ; enddo
175  do j=js-1,je+1 ; do i=is-1,ie+1
176  cs%sin_struct(i,j,1) = -sin(2.0*lat_rad(i,j)) * sin(lon_rad(i,j))
177  cs%cos_struct(i,j,1) = sin(2.0*lat_rad(i,j)) * cos(lon_rad(i,j))
178  cs%sin_struct(i,j,2) = -cos(lat_rad(i,j))**2 * sin(2.0*lon_rad(i,j))
179  cs%cos_struct(i,j,2) = cos(lat_rad(i,j))**2 * cos(2.0*lon_rad(i,j))
180  cs%sin_struct(i,j,3) = 0.0
181  cs%cos_struct(i,j,3) = (0.5-1.5*sin(lat_rad(i,j))**2)
182  enddo ; enddo
183 
184  call get_param(param_file, mdl, "TIDE_M2", use_m2, &
185  "If true, apply tidal momentum forcing at the M2 \n"//&
186  "frequency. This is only used if TIDES is true.", &
187  default=.false.)
188  call get_param(param_file, mdl, "TIDE_S2", use_s2, &
189  "If true, apply tidal momentum forcing at the S2 \n"//&
190  "frequency. This is only used if TIDES is true.", &
191  default=.false.)
192  call get_param(param_file, mdl, "TIDE_N2", use_n2, &
193  "If true, apply tidal momentum forcing at the N2 \n"//&
194  "frequency. This is only used if TIDES is true.", &
195  default=.false.)
196  call get_param(param_file, mdl, "TIDE_K2", use_k2, &
197  "If true, apply tidal momentum forcing at the K2 \n"//&
198  "frequency. This is only used if TIDES is true.", &
199  default=.false.)
200  call get_param(param_file, mdl, "TIDE_K1", use_k1, &
201  "If true, apply tidal momentum forcing at the K1 \n"//&
202  "frequency. This is only used if TIDES is true.", &
203  default=.false.)
204  call get_param(param_file, mdl, "TIDE_O1", use_o1, &
205  "If true, apply tidal momentum forcing at the O1 \n"//&
206  "frequency. This is only used if TIDES is true.", &
207  default=.false.)
208  call get_param(param_file, mdl, "TIDE_P1", use_p1, &
209  "If true, apply tidal momentum forcing at the P1 \n"//&
210  "frequency. This is only used if TIDES is true.", &
211  default=.false.)
212  call get_param(param_file, mdl, "TIDE_Q1", use_q1, &
213  "If true, apply tidal momentum forcing at the Q1 \n"//&
214  "frequency. This is only used if TIDES is true.", &
215  default=.false.)
216  call get_param(param_file, mdl, "TIDE_MF", use_mf, &
217  "If true, apply tidal momentum forcing at the MF \n"//&
218  "frequency. This is only used if TIDES is true.", &
219  default=.false.)
220  call get_param(param_file, mdl, "TIDE_MM", use_mm, &
221  "If true, apply tidal momentum forcing at the MM \n"//&
222  "frequency. This is only used if TIDES is true.", &
223  default=.false.)
224 
225  ! Determine how many tidal components are to be used.
226  nc = 0
227  if (use_m2) nc=nc+1 ; if (use_s2) nc=nc+1
228  if (use_n2) nc=nc+1 ; if (use_k2) nc=nc+1
229  if (use_k1) nc=nc+1 ; if (use_o1) nc=nc+1
230  if (use_p1) nc=nc+1 ; if (use_q1) nc=nc+1
231  if (use_mf) nc=nc+1 ; if (use_mm) nc=nc+1
232  cs%nc = nc
233 
234  if (nc == 0) then
235  call mom_error(fatal, "tidal_forcing_init: "// &
236  "TIDES are defined, but no tidal constituents are used.")
237  return
238  endif
239 
240  call get_param(param_file, mdl, "TIDAL_SAL_FROM_FILE", cs%tidal_sal_from_file, &
241  "If true, read the tidal self-attraction and loading \n"//&
242  "from input files, specified by TIDAL_INPUT_FILE. \n"//&
243  "This is only used if TIDES is true.", default=.false.)
244  call get_param(param_file, mdl, "USE_PREVIOUS_TIDES", cs%use_prev_tides, &
245  "If true, use the SAL from the previous iteration of the \n"//&
246  "tides to facilitate convergent iteration. \n"//&
247  "This is only used if TIDES is true.", default=.false.)
248  call get_param(param_file, mdl, "TIDE_USE_SAL_SCALAR", cs%use_sal_scalar, &
249  "If true and TIDES is true, use the scalar approximation \n"//&
250  "when calculating self-attraction and loading.", &
251  default=.not.cs%tidal_sal_from_file)
252  ! If it is being used, sal_scalar MUST be specified in param_file.
253  if (cs%use_sal_scalar .or. cs%use_prev_tides) &
254  call get_param(param_file, mdl, "TIDE_SAL_SCALAR_VALUE", cs%sal_scalar, &
255  "The constant of proportionality between sea surface \n"//&
256  "height (really it should be bottom pressure) anomalies \n"//&
257  "and bottom geopotential anomalies. This is only used if \n"//&
258  "TIDES and TIDE_USE_SAL_SCALAR are true.", units="m m-1", &
259  fail_if_missing=.true.)
260 
261  if (nc > max_constituents) then
262  write(mesg,'("Increase MAX_CONSTITUENTS in MOM_tidal_forcing.F90 to at least",I3, &
263  &"to accomodate all the registered tidal constituents.")') nc
264  call mom_error(fatal, "MOM_tidal_forcing"//mesg)
265  endif
266 
267  do c=1,4*max_constituents ; tidal_input_files(c) = "" ; enddo
268 
269  if (cs%tidal_sal_from_file .or. cs%use_prev_tides) then
270  call get_param(param_file, mdl, "TIDAL_INPUT_FILE", tidal_input_files, &
271  "A list of input files for tidal information.", &
272  default = "", fail_if_missing=.true.)
273  endif
274 
275  ! Set the parameters for all components that are in use.
276  c=0
277  if (use_m2) then
278  c=c+1 ; cs%const_name(c) = "M2" ; cs%freq(c) = 1.4051890e-4 ; cs%struct(c) = 2
279  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.242334 ; cs%phase0(c) = 0.0
280  freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
281  amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
282  endif
283 
284  if (use_s2) then
285  c=c+1 ; cs%const_name(c) = "S2" ; cs%freq(c) = 1.4544410e-4 ; cs%struct(c) = 2
286  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.112743 ; cs%phase0(c) = 0.0
287  freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
288  amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
289  endif
290 
291  if (use_n2) then
292  c=c+1 ; cs%const_name(c) = "N2" ; cs%freq(c) = 1.3787970e-4 ; cs%struct(c) = 2
293  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.046397 ; cs%phase0(c) = 0.0
294  freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
295  amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
296  endif
297 
298  if (use_k2) then
299  c=c+1 ; cs%const_name(c) = "K2" ; cs%freq(c) = 1.4584234e-4 ; cs%struct(c) = 2
300  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.030684 ; cs%phase0(c) = 0.0
301  freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
302  amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
303  endif
304 
305  if (use_k1) then
306  c=c+1 ; cs%const_name(c) = "K1" ; cs%freq(c) = 0.7292117e-4 ; cs%struct(c) = 1
307  cs%love_no(c) = 0.736 ; cs%amp(c) = 0.141565 ; cs%phase0(c) = 0.0
308  freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
309  amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
310  endif
311 
312  if (use_o1) then
313  c=c+1 ; cs%const_name(c) = "O1" ; cs%freq(c) = 0.6759774e-4 ; cs%struct(c) = 1
314  cs%love_no(c) = 0.695 ; cs%amp(c) = 0.100661 ; cs%phase0(c) = 0.0
315  freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
316  amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
317  endif
318 
319  if (use_p1) then
320  c=c+1 ; cs%const_name(c) = "P1" ; cs%freq(c) = 0.7252295e-4 ; cs%struct(c) = 1
321  cs%love_no(c) = 0.706 ; cs%amp(c) = 0.046848 ; cs%phase0(c) = 0.0
322  freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
323  amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
324  endif
325 
326  if (use_q1) then
327  c=c+1 ; cs%const_name(c) = "Q1" ; cs%freq(c) = 0.6495854e-4 ; cs%struct(c) = 1
328  cs%love_no(c) = 0.695 ; cs%amp(c) = 0.019273 ; cs%phase0(c) = 0.0
329  freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
330  amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
331  endif
332 
333  if (use_mf) then
334  c=c+1 ; cs%const_name(c) = "MF" ; cs%freq(c) = 0.053234e-4 ; cs%struct(c) = 3
335  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.042041 ; cs%phase0(c) = 0.0
336  freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
337  amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
338  endif
339 
340  if (use_mm) then
341  c=c+1 ; cs%const_name(c) = "MM" ; cs%freq(c) = 0.026392e-4 ; cs%struct(c) = 3
342  cs%love_no(c) = 0.693 ; cs%amp(c) = 0.022191 ; cs%phase0(c) = 0.0
343  freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
344  amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
345  endif
346 
347  ! Parse the input file to potentially override the default values for the
348  ! frequency, amplitude and initial phase of each constituent, and log the
349  ! values that are actually used.
350  do c=1,nc
351  call get_param(param_file, mdl, "TIDE_"//trim(cs%const_name(c))//"_FREQ", cs%freq(c), &
352  "Frequency of the "//trim(cs%const_name(c))//" tidal constituent. \n"//&
353  "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
354  " are true.", units="s-1", default=freq_def(c))
355  call get_param(param_file, mdl, "TIDE_"//trim(cs%const_name(c))//"_AMP", cs%amp(c), &
356  "Amplitude of the "//trim(cs%const_name(c))//" tidal constituent. \n"//&
357  "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
358  " are true.", units="m", default=amp_def(c))
359  call get_param(param_file, mdl, "TIDE_"//trim(cs%const_name(c))//"_PHASE_T0", cs%phase0(c), &
360  "Phase of the "//trim(cs%const_name(c))//" tidal constituent at time 0. \n"//&
361  "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
362  " are true.", units="radians", default=phase0_def(c))
363  enddo
364 
365  if (cs%tidal_sal_from_file) then
366  allocate(cs%cosphasesal(isd:ied,jsd:jed,nc))
367  allocate(cs%sinphasesal(isd:ied,jsd:jed,nc))
368  allocate(cs%ampsal(isd:ied,jsd:jed,nc))
369  do c=1,nc
370  ! Read variables with names like PHASE_SAL_M2 and AMP_SAL_M2.
371  call find_in_files(tidal_input_files,"PHASE_SAL_"//trim(cs%const_name(c)),phase,g)
372  call find_in_files(tidal_input_files,"AMP_SAL_"//trim(cs%const_name(c)),cs%ampsal(:,:,c),g)
373  call pass_var(phase, g%domain,complete=.false.)
374  call pass_var(cs%ampsal(:,:,c),g%domain,complete=.true.)
375  do j=js-1,je+1 ; do i=is-1,ie+1
376  cs%cosphasesal(i,j,c) = cos(phase(i,j)*deg_to_rad)
377  cs%sinphasesal(i,j,c) = sin(phase(i,j)*deg_to_rad)
378  enddo ; enddo
379  enddo
380  endif
381 
382  if (cs%USE_PREV_TIDES) then
383  allocate(cs%cosphase_prev(isd:ied,jsd:jed,nc))
384  allocate(cs%sinphase_prev(isd:ied,jsd:jed,nc))
385  allocate(cs%amp_prev(isd:ied,jsd:jed,nc))
386  do c=1,nc
387  ! Read variables with names like PHASE_PREV_M2 and AMP_PREV_M2.
388  call find_in_files(tidal_input_files,"PHASE_PREV_"//trim(cs%const_name(c)),phase,g)
389  call find_in_files(tidal_input_files,"AMP_PREV_"//trim(cs%const_name(c)),cs%amp_prev(:,:,c),g)
390  call pass_var(phase, g%domain,complete=.false.)
391  call pass_var(cs%amp_prev(:,:,c),g%domain,complete=.true.)
392  do j=js-1,je+1 ; do i=is-1,ie+1
393  cs%cosphase_prev(i,j,c) = cos(phase(i,j)*deg_to_rad)
394  cs%sinphase_prev(i,j,c) = sin(phase(i,j)*deg_to_rad)
395  enddo ; enddo
396  enddo
397  endif
398 
399  id_clock_tides = cpu_clock_id('(Ocean tides)', grain=clock_module)
400 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:

◆ 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 437 of file MOM_tidal_forcing.F90.

Referenced by mom_barotropic::btstep(), mom_legacy_barotropic::legacy_btstep(), and mom_barotropic::set_dtbt().

437  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
438  type(tidal_forcing_cs), pointer :: cs !< The control structure returned by a previous call to
439  !! tidal_forcing_init.
440  real, intent(out) :: deta_tidal_deta !< The partial derivative of eta_tidal with
441  !! the local value of eta, nondim.
442 ! This subroutine calculates returns the partial derivative of the local
443 ! geopotential height with the input sea surface height due to self-attraction
444 ! and loading.
445 ! (in) G - The ocean's grid structure.
446 ! (in) CS - The control structure returned by a previous call to
447 ! tidal_forcing_init.
448 ! (out) deta_tidal_deta - the partial derivative of eta_tidal with the
449 ! local value of eta, nondim.
450 
451  if (cs%USE_SAL_SCALAR .and. cs%USE_PREV_TIDES) then
452  deta_tidal_deta = 2.0*cs%SAL_SCALAR
453  elseif (cs%USE_SAL_SCALAR .or. cs%USE_PREV_TIDES) then
454  deta_tidal_deta = cs%SAL_SCALAR
455  else
456  deta_tidal_deta = 0.0
457  endif
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the caller graph for this function:

Variable Documentation

◆ id_clock_tides

integer mom_tidal_forcing::id_clock_tides

Definition at line 103 of file MOM_tidal_forcing.F90.

Referenced by calc_tidal_forcing(), and tidal_forcing_init().

103 integer :: id_clock_tides

◆ max_constituents

integer, parameter mom_tidal_forcing::max_constituents = 10
private

Definition at line 69 of file MOM_tidal_forcing.F90.

Referenced by tidal_forcing_init().

69 integer, parameter :: max_constituents = 10 ! The maximum number of tidal