MOM6
mom_int_tide_input Module Reference

Detailed Description

Calculates energy input to the internal tides.

Data Types

type  int_tide_input_cs
 This control structure holds parameters that regulate internal tide energy inputs. More...
 
type  int_tide_input_type
 This type is used to exchange fields related to the internal tides. More...
 

Functions/Subroutines

subroutine, public set_int_tide_input (u, v, h, tv, fluxes, itide, dt, G, GV, US, CS)
 Sets the model-state dependent internal tide energy sources. More...
 
subroutine find_n2_bottom (h, tv, T_f, S_f, h2, fluxes, G, GV, US, N2_bot)
 Estimates the near-bottom buoyancy frequency (N^2). More...
 
subroutine, public int_tide_input_init (Time, G, GV, US, param_file, diag, CS, itide)
 Initializes the data related to the internal tide input module. More...
 
subroutine, public int_tide_input_end (CS)
 Deallocates any memory related to the internal tide input module. More...
 

Function/Subroutine Documentation

◆ find_n2_bottom()

subroutine mom_int_tide_input::find_n2_bottom ( real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
type(thermo_var_ptrs), intent(in)  tv,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  T_f,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  S_f,
real, dimension(szi_(g),szj_(g)), intent(in)  h2,
type(forcing), intent(in)  fluxes,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  N2_bot 
)
private

Estimates the near-bottom buoyancy frequency (N^2).

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]usA dimensional unit scaling type
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]tvA structure containing pointers to the thermodynamic fields
[in]t_fTemperature after vertical filtering to smooth out the values in thin layers [degC].
[in]s_fSalinity after vertical filtering to smooth out the values in thin layers [ppt].
[in]h2Bottom topographic roughness [Z2 ~> m2].
[in]fluxesA structure of thermodynamic surface fluxes
[out]n2_botThe squared buoyancy freqency at the ocean bottom [T-2 ~> s-2].

Definition at line 151 of file MOM_internal_tide_input.F90.

151  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
152  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
153  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
154  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
155  type(thermo_var_ptrs), intent(in) :: tv !< A structure containing pointers to the
156  !! thermodynamic fields
157  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: T_f !< Temperature after vertical filtering to
158  !! smooth out the values in thin layers [degC].
159  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: S_f !< Salinity after vertical filtering to
160  !! smooth out the values in thin layers [ppt].
161  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: h2 !< Bottom topographic roughness [Z2 ~> m2].
162  type(forcing), intent(in) :: fluxes !< A structure of thermodynamic surface fluxes
163  type(int_tide_input_CS), pointer :: CS !< This module's control structure.
164  real, dimension(SZI_(G),SZJ_(G)), intent(out) :: N2_bot !< The squared buoyancy freqency at the
165  !! ocean bottom [T-2 ~> s-2].
166  ! Local variables
167  real, dimension(SZI_(G),SZK_(G)+1) :: &
168  dRho_int ! The unfiltered density differences across interfaces [R ~> kg m-3].
169  real, dimension(SZI_(G)) :: &
170  pres, & ! The pressure at each interface [R L2 T-2 ~> Pa].
171  Temp_int, & ! The temperature at each interface [degC].
172  Salin_int, & ! The salinity at each interface [ppt].
173  drho_bot, & ! The density difference at the bottom of a layer [R ~> kg m-3]
174  h_amp, & ! The amplitude of topographic roughness [Z ~> m].
175  hb, & ! The depth below a layer [Z ~> m].
176  z_from_bot, & ! The height of a layer center above the bottom [Z ~> m].
177  dRho_dT, & ! The partial derivative of density with temperature [R degC-1 ~> kg m-3 degC-1]
178  dRho_dS ! The partial derivative of density with salinity [R ppt-1 ~> kg m-3 ppt-1].
179 
180  real :: dz_int ! The thickness associated with an interface [Z ~> m].
181  real :: G_Rho0 ! The gravitation acceleration divided by the Boussinesq
182  ! density [Z T-2 R-1 ~> m4 s-2 kg-1].
183  logical :: do_i(SZI_(G)), do_any
184  integer, dimension(2) :: EOSdom ! The i-computational domain for the equation of state
185  integer :: i, j, k, is, ie, js, je, nz
186 
187  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
188  g_rho0 = (us%L_to_Z**2*gv%g_Earth) / gv%Rho0
189  eosdom(:) = eos_domain(g%HI)
190 
191  ! Find the (limited) density jump across each interface.
192  do i=is,ie
193  drho_int(i,1) = 0.0 ; drho_int(i,nz+1) = 0.0
194  enddo
195 !$OMP parallel do default(none) shared(is,ie,js,je,nz,tv,fluxes,G,GV,US,h,T_f,S_f, &
196 !$OMP h2,N2_bot,G_Rho0,EOSdom) &
197 !$OMP private(pres,Temp_Int,Salin_Int,dRho_dT,dRho_dS, &
198 !$OMP hb,dRho_bot,z_from_bot,do_i,h_amp, &
199 !$OMP do_any,dz_int) &
200 !$OMP firstprivate(dRho_int)
201  do j=js,je
202  if (associated(tv%eqn_of_state)) then
203  if (associated(fluxes%p_surf)) then
204  do i=is,ie ; pres(i) = fluxes%p_surf(i,j) ; enddo
205  else
206  do i=is,ie ; pres(i) = 0.0 ; enddo
207  endif
208  do k=2,nz
209  do i=is,ie
210  pres(i) = pres(i) + (gv%g_Earth*gv%H_to_RZ)*h(i,j,k-1)
211  temp_int(i) = 0.5 * (t_f(i,j,k) + t_f(i,j,k-1))
212  salin_int(i) = 0.5 * (s_f(i,j,k) + s_f(i,j,k-1))
213  enddo
214  call calculate_density_derivs(temp_int, salin_int, pres, drho_dt(:), drho_ds(:), &
215  tv%eqn_of_state, eosdom)
216  do i=is,ie
217  drho_int(i,k) = max(drho_dt(i)*(t_f(i,j,k) - t_f(i,j,k-1)) + &
218  drho_ds(i)*(s_f(i,j,k) - s_f(i,j,k-1)), 0.0)
219  enddo
220  enddo
221  else
222  do k=2,nz ; do i=is,ie
223  drho_int(i,k) = (gv%Rlay(k) - gv%Rlay(k-1))
224  enddo ; enddo
225  endif
226 
227  ! Find the bottom boundary layer stratification.
228  do i=is,ie
229  hb(i) = 0.0 ; drho_bot(i) = 0.0
230  z_from_bot(i) = 0.5*gv%H_to_Z*h(i,j,nz)
231  do_i(i) = (g%mask2dT(i,j) > 0.5)
232  h_amp(i) = sqrt(h2(i,j))
233  enddo
234 
235  do k=nz,2,-1
236  do_any = .false.
237  do i=is,ie ; if (do_i(i)) then
238  dz_int = 0.5*gv%H_to_Z*(h(i,j,k) + h(i,j,k-1))
239  z_from_bot(i) = z_from_bot(i) + dz_int ! middle of the layer above
240 
241  hb(i) = hb(i) + dz_int
242  drho_bot(i) = drho_bot(i) + drho_int(i,k)
243 
244  if (z_from_bot(i) > h_amp(i)) then
245  if (k>2) then
246  ! Always include at least one full layer.
247  hb(i) = hb(i) + 0.5*gv%H_to_Z*(h(i,j,k-1) + h(i,j,k-2))
248  drho_bot(i) = drho_bot(i) + drho_int(i,k-1)
249  endif
250  do_i(i) = .false.
251  else
252  do_any = .true.
253  endif
254  endif ; enddo
255  if (.not.do_any) exit
256  enddo
257 
258  do i=is,ie
259  if (hb(i) > 0.0) then
260  n2_bot(i,j) = (g_rho0 * drho_bot(i)) / hb(i)
261  else ; n2_bot(i,j) = 0.0 ; endif
262  enddo
263  enddo
264 

◆ int_tide_input_end()

subroutine, public mom_int_tide_input::int_tide_input_end ( type(int_tide_input_cs), pointer  CS)

Deallocates any memory related to the internal tide input module.

Parameters
csThis module's control structure, which is deallocated here.

Definition at line 426 of file MOM_internal_tide_input.F90.

426  type(int_tide_input_CS), pointer :: CS !< This module's control structure, which is deallocated here.
427 
428  if (associated(cs)) deallocate(cs)
429 

◆ int_tide_input_init()

subroutine, public mom_int_tide_input::int_tide_input_init ( type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(inout), target  diag,
type(int_tide_input_cs), pointer  CS,
type(int_tide_input_type), pointer  itide 
)

Initializes the data related to the internal tide input module.

Parameters
[in]timeThe current model time
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]usA dimensional unit scaling type
[in]param_fileA structure to parse for run-time parameters
[in,out]diagstructure used to regulate diagnostic output.
csThis module's control structure, which is initialized here.
itideA structure containing fields related to the internal tide sources.

Definition at line 269 of file MOM_internal_tide_input.F90.

269  type(time_type), intent(in) :: Time !< The current model time
270  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
271  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
272  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
273  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
274  type(diag_ctrl), target, intent(inout) :: diag !< structure used to regulate diagnostic output.
275  type(int_tide_input_CS), pointer :: CS !< This module's control structure, which is initialized here.
276  type(int_tide_input_type), pointer :: itide !< A structure containing fields related
277  !! to the internal tide sources.
278  ! Local variables
279  type(vardesc) :: vd
280  logical :: read_tideamp
281  ! This include declares and sets the variable "version".
282 # include "version_variable.h"
283  character(len=40) :: mdl = "MOM_int_tide_input" ! This module's name.
284  character(len=20) :: tmpstr
285  character(len=200) :: filename, tideamp_file, h2_file
286 
287  real :: mask_itidal ! A multiplicative land mask, 0 or 1 [nondim]
288  real :: max_frac_rough ! The fraction relating the maximum topographic roughness
289  ! to the mean depth [nondim]
290  real :: utide ! constant tidal amplitude [L T-1 ~> m s-1] to be used if
291  ! tidal amplitude file is not present.
292  real :: kappa_h2_factor ! factor for the product of wavenumber * rms sgs height [nondim].
293  real :: kappa_itides ! topographic wavenumber and non-dimensional scaling [L-1 ~> m-1]
294  real :: min_zbot_itides ! Minimum ocean depth for internal tide conversion [Z ~> m].
295  integer :: tlen_days !< Time interval from start for adding wave source
296  !! for testing internal tides (BDM)
297  integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
298 
299  if (associated(cs)) then
300  call mom_error(warning, "int_tide_input_init called with an associated "// &
301  "control structure.")
302  return
303  endif
304  if (associated(itide)) then
305  call mom_error(warning, "int_tide_input_init called with an associated "// &
306  "internal tide input type.")
307  return
308  endif
309  allocate(cs)
310  allocate(itide)
311 
312  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
313  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
314 
315  cs%diag => diag
316 
317  ! Read all relevant parameters and write them to the model log.
318  call log_version(param_file, mdl, version, "")
319 
320  call get_param(param_file, mdl, "INPUTDIR", cs%inputdir, default=".")
321  cs%inputdir = slasher(cs%inputdir)
322 
323  call get_param(param_file, mdl, "DEBUG", cs%debug, default=.false., do_not_log=.true.)
324 
325  call get_param(param_file, mdl, "MIN_ZBOT_ITIDES", min_zbot_itides, &
326  "Turn off internal tidal dissipation when the total "//&
327  "ocean depth is less than this value.", units="m", default=0.0, scale=us%m_to_Z)
328  call get_param(param_file, mdl, "KD_SMOOTH", cs%kappa_fill, &
329  "A diapycnal diffusivity that is used to interpolate "//&
330  "more sensible values of T & S into thin layers.", &
331  units="m2 s-1", default=1.0e-6, scale=us%m2_s_to_Z2_T)
332 
333  call get_param(param_file, mdl, "UTIDE", utide, &
334  "The constant tidal amplitude used with INT_TIDE_DISSIPATION.", &
335  units="m s-1", default=0.0, scale=us%m_s_to_L_T)
336 
337  allocate(itide%Nb(isd:ied,jsd:jed)) ; itide%Nb(:,:) = 0.0
338  allocate(itide%h2(isd:ied,jsd:jed)) ; itide%h2(:,:) = 0.0
339  allocate(itide%TKE_itidal_input(isd:ied,jsd:jed)) ; itide%TKE_itidal_input(:,:) = 0.0
340  allocate(itide%tideamp(isd:ied,jsd:jed)) ; itide%tideamp(:,:) = utide
341  allocate(cs%TKE_itidal_coef(isd:ied,jsd:jed)) ; cs%TKE_itidal_coef(:,:) = 0.0
342 
343  call get_param(param_file, mdl, "KAPPA_ITIDES", kappa_itides, &
344  "A topographic wavenumber used with INT_TIDE_DISSIPATION. "//&
345  "The default is 2pi/10 km, as in St.Laurent et al. 2002.", &
346  units="m-1", default=8.e-4*atan(1.0), scale=us%L_to_m)
347 
348  call get_param(param_file, mdl, "KAPPA_H2_FACTOR", kappa_h2_factor, &
349  "A scaling factor for the roughness amplitude with n"//&
350  "INT_TIDE_DISSIPATION.", units="nondim", default=1.0)
351  call get_param(param_file, mdl, "TKE_ITIDE_MAX", cs%TKE_itide_max, &
352  "The maximum internal tide energy source available to mix "//&
353  "above the bottom boundary layer with INT_TIDE_DISSIPATION.", &
354  units="W m-2", default=1.0e3, scale=us%W_m2_to_RZ3_T3)
355 
356  call get_param(param_file, mdl, "READ_TIDEAMP", read_tideamp, &
357  "If true, read a file (given by TIDEAMP_FILE) containing "//&
358  "the tidal amplitude with INT_TIDE_DISSIPATION.", default=.false.)
359  if (read_tideamp) then
360  call get_param(param_file, mdl, "TIDEAMP_FILE", tideamp_file, &
361  "The path to the file containing the spatially varying "//&
362  "tidal amplitudes with INT_TIDE_DISSIPATION.", default="tideamp.nc")
363  filename = trim(cs%inputdir) // trim(tideamp_file)
364  call log_param(param_file, mdl, "INPUTDIR/TIDEAMP_FILE", filename)
365  call mom_read_data(filename, 'tideamp', itide%tideamp, g%domain, timelevel=1, scale=us%m_s_to_L_T)
366  endif
367 
368  call get_param(param_file, mdl, "H2_FILE", h2_file, &
369  "The path to the file containing the sub-grid-scale "//&
370  "topographic roughness amplitude with INT_TIDE_DISSIPATION.", &
371  fail_if_missing=.true.)
372  filename = trim(cs%inputdir) // trim(h2_file)
373  call log_param(param_file, mdl, "INPUTDIR/H2_FILE", filename)
374  call mom_read_data(filename, 'h2', itide%h2, g%domain, timelevel=1, scale=us%m_to_Z**2)
375 
376  call get_param(param_file, mdl, "FRACTIONAL_ROUGHNESS_MAX", max_frac_rough, &
377  "The maximum topographic roughness amplitude as a fraction of the mean depth, "//&
378  "or a negative value for no limitations on roughness.", &
379  units="nondim", default=0.1)
380 
381  ! The following parameters are used in testing the internal tide code.
382  call get_param(param_file, mdl, "INTERNAL_TIDE_SOURCE_TEST", cs%int_tide_source_test, &
383  "If true, apply an arbitrary generation site for internal tide testing", &
384  default=.false.)
385  if (cs%int_tide_source_test)then
386  call get_param(param_file, mdl, "INTERNAL_TIDE_SOURCE_X", cs%int_tide_source_x, &
387  "X Location of generation site for internal tide", default=1.)
388  call get_param(param_file, mdl, "INTERNAL_TIDE_SOURCE_Y", cs%int_tide_source_y, &
389  "Y Location of generation site for internal tide", default=1.)
390  call get_param(param_file, mdl, "INTERNAL_TIDE_SOURCE_TLEN_DAYS", tlen_days, &
391  "Time interval from start of experiment for adding wave source", &
392  units="days", default=0)
393  cs%time_max_source = time + set_time(0, days=tlen_days)
394  endif
395 
396  do j=js,je ; do i=is,ie
397  mask_itidal = 1.0
398  if (g%bathyT(i,j) < min_zbot_itides) mask_itidal = 0.0
399 
400  itide%tideamp(i,j) = itide%tideamp(i,j) * mask_itidal * g%mask2dT(i,j)
401 
402  ! Restrict rms topo to a fraction (often 10 percent) of the column depth.
403  if (max_frac_rough >= 0.0) &
404  itide%h2(i,j) = min((max_frac_rough*g%bathyT(i,j))**2, itide%h2(i,j))
405 
406  ! Compute the fixed part of internal tidal forcing; units are [R Z3 T-2 ~> J m-2] here.
407  cs%TKE_itidal_coef(i,j) = 0.5*us%L_to_Z*kappa_h2_factor*gv%Rho0*&
408  kappa_itides * itide%h2(i,j) * itide%tideamp(i,j)**2
409  enddo ; enddo
410 
411 
412  cs%id_TKE_itidal = register_diag_field('ocean_model','TKE_itidal_itide',diag%axesT1,time, &
413  'Internal Tide Driven Turbulent Kinetic Energy', &
414  'W m-2', conversion=us%RZ3_T3_to_W_m2)
415 
416  cs%id_Nb = register_diag_field('ocean_model','Nb_itide',diag%axesT1,time, &
417  'Bottom Buoyancy Frequency', 's-1', conversion=us%s_to_T)
418 
419  cs%id_N2_bot = register_diag_field('ocean_model','N2_b_itide',diag%axesT1,time, &
420  'Bottom Buoyancy frequency squared', 's-2', conversion=us%s_to_T**2)
421 

◆ set_int_tide_input()

subroutine, public mom_int_tide_input::set_int_tide_input ( real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(in)  u,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(in)  v,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(thermo_var_ptrs), intent(in)  tv,
type(forcing), intent(in)  fluxes,
type(int_tide_input_type), intent(inout)  itide,
real, intent(in)  dt,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(int_tide_input_cs), pointer  CS 
)

Sets the model-state dependent internal tide energy sources.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]usA dimensional unit scaling type
[in]uThe zonal velocity [L T-1 ~> m s-1]
[in]vThe meridional velocity [L T-1 ~> m s-1]
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]tvA structure containing pointers to the thermodynamic fields
[in]fluxesA structure of thermodynamic surface fluxes
[in,out]itideA structure containing fields related to the internal tide sources.
[in]dtThe time increment [T ~> s].
csThis module's control structure.

Definition at line 75 of file MOM_internal_tide_input.F90.

75  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
76  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
77  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
78  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1]
79  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1]
80  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
81  type(thermo_var_ptrs), intent(in) :: tv !< A structure containing pointers to the
82  !! thermodynamic fields
83  type(forcing), intent(in) :: fluxes !< A structure of thermodynamic surface fluxes
84  type(int_tide_input_type), intent(inout) :: itide !< A structure containing fields related
85  !! to the internal tide sources.
86  real, intent(in) :: dt !< The time increment [T ~> s].
87  type(int_tide_input_CS), pointer :: CS !< This module's control structure.
88  ! Local variables
89  real, dimension(SZI_(G),SZJ_(G)) :: &
90  N2_bot ! The bottom squared buoyancy frequency [T-2 ~> s-2].
91 
92  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
93  T_f, S_f ! The temperature and salinity in [degC] and [ppt] with the values in
94  ! the massless layers filled vertically by diffusion.
95  logical :: use_EOS ! If true, density is calculated from T & S using an
96  ! equation of state.
97  logical :: avg_enabled ! for testing internal tides (BDM)
98  type(time_type) :: time_end !< For use in testing internal tides (BDM)
99 
100  integer :: i, j, k, is, ie, js, je, nz, isd, ied, jsd, jed
101 
102  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
103  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
104 
105  if (.not.associated(cs)) call mom_error(fatal,"set_diffusivity: "//&
106  "Module must be initialized before it is used.")
107 
108  use_eos = associated(tv%eqn_of_state)
109 
110  ! Smooth the properties through massless layers.
111  if (use_eos) then
112  call vert_fill_ts(h, tv%T, tv%S, cs%kappa_fill*dt, t_f, s_f, g, gv, larger_h_denom=.true.)
113  endif
114 
115  call find_n2_bottom(h, tv, t_f, s_f, itide%h2, fluxes, g, gv, us, n2_bot)
116 
117  !$OMP parallel do default(shared)
118  do j=js,je ; do i=is,ie
119  itide%Nb(i,j) = g%mask2dT(i,j) * sqrt(n2_bot(i,j))
120  itide%TKE_itidal_input(i,j) = min(cs%TKE_itidal_coef(i,j)*itide%Nb(i,j), cs%TKE_itide_max)
121  enddo ; enddo
122 
123  if (cs%int_tide_source_test) then
124  itide%TKE_itidal_input(:,:) = 0.0
125  avg_enabled = query_averaging_enabled(cs%diag, time_end=time_end)
126  if (time_end <= cs%time_max_source) then
127  do j=js,je ; do i=is,ie
128  ! Input an arbitrary energy point source.id_
129  if (((g%geoLonCu(i-1,j)-cs%int_tide_source_x) * (g%geoLonBu(i,j)-cs%int_tide_source_x) <= 0.0) .and. &
130  ((g%geoLatCv(i,j-1)-cs%int_tide_source_y) * (g%geoLatCv(i,j)-cs%int_tide_source_y) <= 0.0)) then
131  itide%TKE_itidal_input(i,j) = 1.0*us%kg_m3_to_R*us%m_to_Z**3*us%T_to_s**3
132  endif
133  enddo ; enddo
134  endif
135  endif
136 
137  if (cs%debug) then
138  call hchksum(n2_bot,"N2_bot",g%HI,haloshift=0, scale=us%s_to_T**2)
139  call hchksum(itide%TKE_itidal_input,"TKE_itidal_input",g%HI,haloshift=0, &
140  scale=us%RZ3_T3_to_W_m2)
141  endif
142 
143  if (cs%id_TKE_itidal > 0) call post_data(cs%id_TKE_itidal, itide%TKE_itidal_input, cs%diag)
144  if (cs%id_Nb > 0) call post_data(cs%id_Nb, itide%Nb, cs%diag)
145  if (cs%id_N2_bot > 0 ) call post_data(cs%id_N2_bot, n2_bot, cs%diag)
146