MOM6
MOM_tidal_mixing.F90
1 !> Interface to vertical tidal mixing schemes including CVMix tidal mixing.
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_diag_mediator, only : diag_ctrl, time_type, register_diag_field
7 use mom_diag_mediator, only : safe_alloc_ptr, post_data
8 use mom_debugging, only : hchksum
9 use mom_eos, only : calculate_density
10 use mom_error_handler, only : mom_error, is_root_pe, fatal, warning, note
11 use mom_file_parser, only : openparameterblock, closeparameterblock
13 use mom_grid, only : ocean_grid_type
14 use mom_io, only : slasher, mom_read_data, field_size
15 use mom_remapping, only : remapping_cs, initialize_remapping, remapping_core_h
16 use mom_string_functions, only : uppercase, lowercase
20 use cvmix_tidal, only : cvmix_init_tidal, cvmix_compute_simmons_invariant
21 use cvmix_tidal, only : cvmix_coeffs_tidal, cvmix_tidal_params_type
22 use cvmix_tidal, only : cvmix_compute_schmittner_invariant, cvmix_compute_schmittnercoeff
23 use cvmix_tidal, only : cvmix_coeffs_tidal_schmittner
24 use cvmix_kinds_and_types, only : cvmix_global_params_type
25 use cvmix_put_get, only : cvmix_put
26 
27 implicit none ; private
28 
29 #include <MOM_memory.h>
30 
31 public tidal_mixing_init
32 public setup_tidal_diagnostics
33 public calculate_tidal_mixing
34 public post_tidal_diagnostics
35 public tidal_mixing_h_amp
36 public tidal_mixing_end
37 
38 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
39 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
40 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
41 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
42 
43 !> Containers for tidal mixing diagnostics
44 type, public :: tidal_mixing_diags ; private
45  real, pointer, dimension(:,:,:) :: &
46  kd_itidal => null(),& !< internal tide diffusivity at interfaces [Z2 T-1 ~> m2 s-1].
47  fl_itidal => null(),& !< vertical flux of tidal turbulent dissipation [Z3 T-3 ~> m3 s-3]
48  kd_niku => null(),& !< lee-wave diffusivity at interfaces [Z2 T-1 ~> m2 s-1].
49  kd_niku_work => null(),& !< layer integrated work by lee-wave driven mixing [R Z3 T-3 ~> W m-2]
50  kd_itidal_work => null(),& !< layer integrated work by int tide driven mixing [R Z3 T-3 ~> W m-2]
51  kd_lowmode_work => null(),& !< layer integrated work by low mode driven mixing [R Z3 T-3 ~> W m-2]
52  n2_int => null(),& !< Bouyancy frequency squared at interfaces [T-2 ~> s-2]
53  vert_dep_3d => null(),& !< The 3-d mixing energy deposition [W m-3]
54  schmittner_coeff_3d => null() !< The coefficient in the Schmittner et al mixing scheme, in UNITS?
55  real, pointer, dimension(:,:,:) :: tidal_qe_md => null() !< Input tidal energy dissipated locally,
56  !! interpolated to model vertical coordinate [W m-3?]
57  real, pointer, dimension(:,:,:) :: kd_lowmode => null() !< internal tide diffusivity at interfaces
58  !! due to propagating low modes [Z2 T-1 ~> m2 s-1].
59  real, pointer, dimension(:,:,:) :: fl_lowmode => null() !< vertical flux of tidal turbulent
60  !! dissipation due to propagating low modes [Z3 T-3 ~> m3 s-3]
61  real, pointer, dimension(:,:) :: &
62  tke_itidal_used => null(),& !< internal tide TKE input at ocean bottom [R Z3 T-3 ~> W m-2]
63  n2_bot => null(),& !< bottom squared buoyancy frequency [T-2 ~> s-2]
64  n2_meanz => null(),& !< vertically averaged buoyancy frequency [T-2 ~> s-2]
65  polzin_decay_scale_scaled => null(),& !< vertical scale of decay for tidal dissipation [Z ~> m]
66  polzin_decay_scale => null(),& !< vertical decay scale for tidal diss with Polzin [Z ~> m]
67  simmons_coeff_2d => null() !< The Simmons et al mixing coefficient
68 
69 end type
70 
71 !> Control structure with parameters for the tidal mixing module.
72 type, public :: tidal_mixing_cs ; private
73  logical :: debug = .true. !< If true, do more extensive debugging checks. This is hard-coded.
74 
75  ! Parameters
76  logical :: int_tide_dissipation = .false. !< Internal tide conversion (from barotropic)
77  !! with the schemes of St Laurent et al (2002) & Simmons et al (2004)
78 
79  integer :: int_tide_profile !< A coded integer indicating the vertical profile
80  !! for dissipation of the internal waves. Schemes that are
81  !! currently encoded are St Laurent et al (2002) and Polzin (2009).
82  logical :: lee_wave_dissipation = .false. !< Enable lee-wave driven mixing, following
83  !! Nikurashin (2010), with a vertical energy
84  !! deposition profile specified by Lee_wave_profile to be
85  !! St Laurent et al (2002) or Simmons et al (2004) scheme
86 
87  integer :: lee_wave_profile !< A coded integer indicating the vertical profile
88  !! for dissipation of the lee waves. Schemes that are
89  !! currently encoded are St Laurent et al (2002) and
90  !! Polzin (2009).
91  real :: int_tide_decay_scale !< decay scale for internal wave TKE [Z ~> m].
92 
93  real :: mu_itides !< efficiency for conversion of dissipation
94  !! to potential energy [nondim]
95 
96  real :: gamma_itides !< fraction of local dissipation [nondim]
97 
98  real :: gamma_lee !< fraction of local dissipation for lee waves
99  !! (Nikurashin's energy input) [nondim]
100  real :: decay_scale_factor_lee !< Scaling factor for the decay scale of lee
101  !! wave energy dissipation [nondim]
102 
103  real :: min_zbot_itides !< minimum depth for internal tide conversion [Z ~> m].
104  logical :: lowmode_itidal_dissipation = .false. !< If true, consider mixing due to breaking low
105  !! modes that have been remotely generated using an internal tidal
106  !! dissipation scheme to specify the vertical profile of the energy
107  !! input to drive diapycnal mixing, along the lines of St. Laurent
108  !! et al. (2002) and Simmons et al. (2004).
109 
110  real :: nu_polzin !< The non-dimensional constant used in Polzin form of
111  !! the vertical scale of decay of tidal dissipation [nondim]
112 
113  real :: nbotref_polzin !< Reference value for the buoyancy frequency at the
114  !! ocean bottom used in Polzin formulation of the
115  !! vertical scale of decay of tidal dissipation [T-1 ~> s-1]
116  real :: polzin_decay_scale_factor !< Scaling factor for the decay length scale
117  !! of the tidal dissipation profile in Polzin [nondim]
118  real :: polzin_decay_scale_max_factor !< The decay length scale of tidal dissipation
119  !! profile in Polzin formulation should not exceed
120  !! Polzin_decay_scale_max_factor * depth of the ocean [nondim].
121  real :: polzin_min_decay_scale !< minimum decay scale of the tidal dissipation
122  !! profile in Polzin formulation [Z ~> m].
123 
124  real :: tke_itide_max !< maximum internal tide conversion [R Z3 T-3 ~> W m-2]
125  !! available to mix above the BBL
126 
127  real :: utide !< constant tidal amplitude [Z T-1 ~> m s-1] if READ_TIDEAMP is false.
128  real :: kappa_itides !< topographic wavenumber and non-dimensional scaling [Z-1 ~> m-1].
129  real :: kappa_h2_factor !< factor for the product of wavenumber * rms sgs height
130  character(len=200) :: inputdir !< The directory in which to find input files
131 
132  logical :: use_cvmix_tidal = .false. !< true if CVMix is to be used for determining
133  !! diffusivity due to tidal mixing
134 
135  real :: min_thickness !< Minimum thickness allowed [m]
136 
137  ! CVMix-specific parameters
138  integer :: cvmix_tidal_scheme = -1 !< 1 for Simmons, 2 for Schmittner
139  type(cvmix_tidal_params_type) :: cvmix_tidal_params !< A CVMix-specific type with parameters for tidal mixing
140  type(cvmix_global_params_type) :: cvmix_glb_params !< CVMix-specific for Prandtl number only
141  real :: tidal_max_coef !< CVMix-specific maximum allowable tidal diffusivity. [m^2/s]
142  real :: tidal_diss_lim_tc !< CVMix-specific dissipation limit depth for
143  !! tidal-energy-constituent data [Z ~> m].
144  type(remapping_cs) :: remap_cs !< The control structure for remapping
145  logical :: remap_answers_2018 = .true. !< If true, use the order of arithmetic and expressions that
146  !! recover the remapping answers from 2018. If false, use more
147  !! robust forms of the same remapping expressions.
148 
149  ! Data containers
150  real, pointer, dimension(:,:) :: tke_niku => null() !< Lee wave driven Turbulent Kinetic Energy input
151  !! [R Z3 T-3 ~> W m-2]
152  real, pointer, dimension(:,:) :: tke_itidal => null() !< The internal Turbulent Kinetic Energy input divided
153  !! by the bottom stratfication [R Z3 T-2 ~> J m-2].
154  real, pointer, dimension(:,:) :: nb => null() !< The near bottom buoyancy frequency [T-1 ~> s-1].
155  real, pointer, dimension(:,:) :: mask_itidal => null() !< A mask of where internal tide energy is input
156  real, pointer, dimension(:,:) :: h2 => null() !< Squared bottom depth variance [Z2 ~> m2].
157  real, pointer, dimension(:,:) :: tideamp => null() !< RMS tidal amplitude [Z T-1 ~> m s-1]
158  real, allocatable, dimension(:) :: h_src !< tidal constituent input layer thickness [m]
159  real, allocatable, dimension(:,:) :: tidal_qe_2d !< Tidal energy input times the local dissipation
160  !! fraction, q*E(x,y), with the CVMix implementation
161  !! of Jayne et al tidal mixing [W m-2].
162  !! TODO: make this E(x,y) only
163  real, allocatable, dimension(:,:,:) :: tidal_qe_3d_in !< q*E(x,y,z) with the Schmittner parameterization [W m-3?]
164 
165  logical :: answers_2018 !< If true, use the order of arithmetic and expressions that recover the
166  !! answers from the end of 2018. Otherwise, use updated and more robust
167  !! forms of the same expressions.
168 
169  ! Diagnostics
170  type(diag_ctrl), pointer :: diag => null() !< structure to regulate diagnostic output timing
171  type(tidal_mixing_diags), pointer :: dd => null() !< A pointer to a structure of diagnostic arrays
172 
173  !>@{ Diagnostic identifiers
174  integer :: id_tke_itidal = -1
175  integer :: id_tke_leewave = -1
176  integer :: id_kd_itidal = -1
177  integer :: id_kd_niku = -1
178  integer :: id_kd_lowmode = -1
179  integer :: id_kd_itidal_work = -1
180  integer :: id_kd_niku_work = -1
181  integer :: id_kd_lowmode_work = -1
182  integer :: id_nb = -1
183  integer :: id_n2_bot = -1
184  integer :: id_n2_meanz = -1
185  integer :: id_fl_itidal = -1
186  integer :: id_fl_lowmode = -1
187  integer :: id_polzin_decay_scale = -1
188  integer :: id_polzin_decay_scale_scaled = -1
189  integer :: id_n2_int = -1
190  integer :: id_simmons_coeff = -1
191  integer :: id_schmittner_coeff = -1
192  integer :: id_tidal_qe_md = -1
193  integer :: id_vert_dep = -1
194  !>@}
195 
196 end type tidal_mixing_cs
197 
198 !>@{ Coded parmameters for specifying mixing schemes
199 character*(20), parameter :: stlaurent_profile_string = "STLAURENT_02"
200 character*(20), parameter :: polzin_profile_string = "POLZIN_09"
201 integer, parameter :: stlaurent_02 = 1
202 integer, parameter :: polzin_09 = 2
203 character*(20), parameter :: simmons_scheme_string = "SIMMONS"
204 character*(20), parameter :: schmittner_scheme_string = "SCHMITTNER"
205 integer, parameter :: simmons = 1
206 integer, parameter :: schmittner = 2
207 !>@}
208 
209 contains
210 
211 !> Initializes internal tidal dissipation scheme for diapycnal mixing
212 logical function tidal_mixing_init(Time, G, GV, US, param_file, diag, CS)
213  type(time_type), intent(in) :: Time !< The current time.
214  type(ocean_grid_type), intent(in) :: G !< Grid structure.
215  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure.
216  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
217  type(param_file_type), intent(in) :: param_file !< Run-time parameter file handle
218  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure.
219  type(tidal_mixing_cs), pointer :: CS !< This module's control structure.
220 
221  ! Local variables
222  logical :: read_tideamp
223  logical :: default_2018_answers
224  character(len=20) :: tmpstr, int_tide_profile_str
225  character(len=20) :: CVMix_tidal_scheme_str, tidal_energy_type
226  character(len=200) :: filename, h2_file, Niku_TKE_input_file
227  character(len=200) :: tidal_energy_file, tideamp_file
228  real :: utide, hamp, prandtl_tidal, max_frac_rough
229  real :: Niku_scale ! local variable for scaling the Nikurashin TKE flux data
230  integer :: i, j, is, ie, js, je
231  integer :: isd, ied, jsd, jed
232  ! This include declares and sets the variable "version".
233 # include "version_variable.h"
234  character(len=40) :: mdl = "MOM_tidal_mixing" !< This module's name.
235 
236  if (associated(cs)) then
237  call mom_error(warning, "tidal_mixing_init called when control structure "// &
238  "is already associated.")
239  return
240  endif
241  allocate(cs)
242  allocate(cs%dd)
243 
244  cs%debug = cs%debug.and.is_root_pe()
245 
246  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
247  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
248 
249  cs%diag => diag
250 
251  ! Read parameters
252  call get_param(param_file, mdl, "USE_CVMix_TIDAL", cs%use_CVMix_tidal, &
253  default=.false., do_not_log=.true.)
254  call get_param(param_file, mdl, "INT_TIDE_DISSIPATION", cs%int_tide_dissipation, &
255  default=cs%use_CVMix_tidal, do_not_log=.true.)
256  call log_version(param_file, mdl, version, &
257  "Vertical Tidal Mixing Parameterization", &
258  all_default=.not.(cs%use_CVMix_tidal .or. cs%int_tide_dissipation))
259  call get_param(param_file, mdl, "USE_CVMix_TIDAL", cs%use_CVMix_tidal, &
260  "If true, turns on tidal mixing via CVMix", &
261  default=.false.)
262 
263  call get_param(param_file, mdl, "INPUTDIR", cs%inputdir, default=".",do_not_log=.true.)
264  cs%inputdir = slasher(cs%inputdir)
265  call get_param(param_file, mdl, "INT_TIDE_DISSIPATION", cs%int_tide_dissipation, &
266  "If true, use an internal tidal dissipation scheme to "//&
267  "drive diapycnal mixing, along the lines of St. Laurent "//&
268  "et al. (2002) and Simmons et al. (2004).", default=cs%use_CVMix_tidal)
269 
270  ! return if tidal mixing is inactive
271  tidal_mixing_init = cs%int_tide_dissipation
272  if (.not. tidal_mixing_init) return
273 
274  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
275  "This sets the default value for the various _2018_ANSWERS parameters.", &
276  default=.false.)
277  call get_param(param_file, mdl, "TIDAL_MIXING_2018_ANSWERS", cs%answers_2018, &
278  "If true, use the order of arithmetic and expressions that recover the "//&
279  "answers from the end of 2018. Otherwise, use updated and more robust "//&
280  "forms of the same expressions.", default=default_2018_answers)
281  call get_param(param_file, mdl, "REMAPPING_2018_ANSWERS", cs%remap_answers_2018, &
282  "If true, use the order of arithmetic and expressions that recover the "//&
283  "answers from the end of 2018. Otherwise, use updated and more robust "//&
284  "forms of the same expressions.", default=default_2018_answers)
285 
286  if (cs%int_tide_dissipation) then
287 
288  ! Read in CVMix tidal scheme if CVMix tidal mixing is on
289  if (cs%use_CVMix_tidal) then
290  call get_param(param_file, mdl, "CVMIX_TIDAL_SCHEME", cvmix_tidal_scheme_str, &
291  "CVMIX_TIDAL_SCHEME selects the CVMix tidal mixing "//&
292  "scheme with INT_TIDE_DISSIPATION. Valid values are:\n"//&
293  "\t SIMMONS - Use the Simmons et al (2004) tidal \n"//&
294  "\t mixing scheme.\n"//&
295  "\t SCHMITTNER - Use the Schmittner et al (2014) tidal \n"//&
296  "\t mixing scheme.", &
297  default=simmons_scheme_string)
298  cvmix_tidal_scheme_str = uppercase(cvmix_tidal_scheme_str)
299 
300  select case (cvmix_tidal_scheme_str)
301  case (simmons_scheme_string) ; cs%CVMix_tidal_scheme = simmons
302  case (schmittner_scheme_string) ; cs%CVMix_tidal_scheme = schmittner
303  case default
304  call mom_error(fatal, "tidal_mixing_init: Unrecognized setting "// &
305  "#define CVMIX_TIDAL_SCHEME "//trim(cvmix_tidal_scheme_str)//" found in input file.")
306  end select
307  endif ! CS%use_CVMix_tidal
308 
309  ! Read in vertical profile of tidal energy dissipation
310  if ( cs%CVMix_tidal_scheme.eq.schmittner .or. .not. cs%use_CVMix_tidal) then
311  call get_param(param_file, mdl, "INT_TIDE_PROFILE", int_tide_profile_str, &
312  "INT_TIDE_PROFILE selects the vertical profile of energy "//&
313  "dissipation with INT_TIDE_DISSIPATION. Valid values are:\n"//&
314  "\t STLAURENT_02 - Use the St. Laurent et al exponential \n"//&
315  "\t decay profile.\n"//&
316  "\t POLZIN_09 - Use the Polzin WKB-stretched algebraic \n"//&
317  "\t decay profile.", &
318  default=stlaurent_profile_string)
319  int_tide_profile_str = uppercase(int_tide_profile_str)
320 
321  select case (int_tide_profile_str)
322  case (stlaurent_profile_string) ; cs%int_tide_profile = stlaurent_02
323  case (polzin_profile_string) ; cs%int_tide_profile = polzin_09
324  case default
325  call mom_error(fatal, "tidal_mixing_init: Unrecognized setting "// &
326  "#define INT_TIDE_PROFILE "//trim(int_tide_profile_str)//" found in input file.")
327  end select
328  endif
329 
330  elseif (cs%use_CVMix_tidal) then
331  call mom_error(fatal, "tidal_mixing_init: Cannot set INT_TIDE_DISSIPATION to False "// &
332  "when USE_CVMix_TIDAL is set to True.")
333  endif
334 
335  call get_param(param_file, mdl, "LEE_WAVE_DISSIPATION", cs%Lee_wave_dissipation, &
336  "If true, use an lee wave driven dissipation scheme to "//&
337  "drive diapycnal mixing, along the lines of Nikurashin "//&
338  "(2010) and using the St. Laurent et al. (2002) "//&
339  "and Simmons et al. (2004) vertical profile", default=.false.)
340  if (cs%lee_wave_dissipation) then
341  if (cs%use_CVMix_tidal) then
342  call mom_error(fatal, "tidal_mixing_init: Lee wave driven dissipation scheme cannot "// &
343  "be used when CVMix tidal mixing scheme is active.")
344  endif
345  call get_param(param_file, mdl, "LEE_WAVE_PROFILE", tmpstr, &
346  "LEE_WAVE_PROFILE selects the vertical profile of energy "//&
347  "dissipation with LEE_WAVE_DISSIPATION. Valid values are:\n"//&
348  "\t STLAURENT_02 - Use the St. Laurent et al exponential \n"//&
349  "\t decay profile.\n"//&
350  "\t POLZIN_09 - Use the Polzin WKB-stretched algebraic \n"//&
351  "\t decay profile.", &
352  default=stlaurent_profile_string)
353  tmpstr = uppercase(tmpstr)
354  select case (tmpstr)
355  case (stlaurent_profile_string) ; cs%lee_wave_profile = stlaurent_02
356  case (polzin_profile_string) ; cs%lee_wave_profile = polzin_09
357  case default
358  call mom_error(fatal, "tidal_mixing_init: Unrecognized setting "// &
359  "#define LEE_WAVE_PROFILE "//trim(tmpstr)//" found in input file.")
360  end select
361  endif
362 
363  call get_param(param_file, mdl, "INT_TIDE_LOWMODE_DISSIPATION", cs%Lowmode_itidal_dissipation, &
364  "If true, consider mixing due to breaking low modes that "//&
365  "have been remotely generated; as with itidal drag on the "//&
366  "barotropic tide, use an internal tidal dissipation scheme to "//&
367  "drive diapycnal mixing, along the lines of St. Laurent "//&
368  "et al. (2002) and Simmons et al. (2004).", default=.false.)
369 
370  if ((cs%Int_tide_dissipation .and. (cs%int_tide_profile == polzin_09)) .or. &
371  (cs%lee_wave_dissipation .and. (cs%lee_wave_profile == polzin_09))) then
372  if (cs%use_CVMix_tidal) then
373  call mom_error(fatal, "tidal_mixing_init: Polzin scheme cannot "// &
374  "be used when CVMix tidal mixing scheme is active.")
375  endif
376  call get_param(param_file, mdl, "NU_POLZIN", cs%Nu_Polzin, &
377  "When the Polzin decay profile is used, this is a "//&
378  "non-dimensional constant in the expression for the "//&
379  "vertical scale of decay for the tidal energy dissipation.", &
380  units="nondim", default=0.0697)
381  call get_param(param_file, mdl, "NBOTREF_POLZIN", cs%Nbotref_Polzin, &
382  "When the Polzin decay profile is used, this is the "//&
383  "reference value of the buoyancy frequency at the ocean "//&
384  "bottom in the Polzin formulation for the vertical "//&
385  "scale of decay for the tidal energy dissipation.", &
386  units="s-1", default=9.61e-4, scale=us%T_to_s)
387  call get_param(param_file, mdl, "POLZIN_DECAY_SCALE_FACTOR", &
388  cs%Polzin_decay_scale_factor, &
389  "When the Polzin decay profile is used, this is a "//&
390  "scale factor for the vertical scale of decay of the tidal "//&
391  "energy dissipation.", default=1.0, units="nondim")
392  call get_param(param_file, mdl, "POLZIN_SCALE_MAX_FACTOR", &
393  cs%Polzin_decay_scale_max_factor, &
394  "When the Polzin decay profile is used, this is a factor "//&
395  "to limit the vertical scale of decay of the tidal "//&
396  "energy dissipation to POLZIN_DECAY_SCALE_MAX_FACTOR "//&
397  "times the depth of the ocean.", units="nondim", default=1.0)
398  call get_param(param_file, mdl, "POLZIN_MIN_DECAY_SCALE", cs%Polzin_min_decay_scale, &
399  "When the Polzin decay profile is used, this is the "//&
400  "minimum vertical decay scale for the vertical profile\n"//&
401  "of internal tide dissipation with the Polzin (2009) formulation", &
402  units="m", default=0.0, scale=us%m_to_Z)
403  endif
404 
405  if (cs%Int_tide_dissipation .or. cs%Lee_wave_dissipation) then
406  call get_param(param_file, mdl, "INT_TIDE_DECAY_SCALE", cs%Int_tide_decay_scale, &
407  "The decay scale away from the bottom for tidal TKE with "//&
408  "the new coding when INT_TIDE_DISSIPATION is used.", &
409  !units="m", default=0.0)
410  units="m", default=500.0, scale=us%m_to_Z) ! TODO: confirm this new default
411  call get_param(param_file, mdl, "MU_ITIDES", cs%Mu_itides, &
412  "A dimensionless turbulent mixing efficiency used with "//&
413  "INT_TIDE_DISSIPATION, often 0.2.", units="nondim", default=0.2)
414  call get_param(param_file, mdl, "GAMMA_ITIDES", cs%Gamma_itides, &
415  "The fraction of the internal tidal energy that is "//&
416  "dissipated locally with INT_TIDE_DISSIPATION. "//&
417  "THIS NAME COULD BE BETTER.", &
418  units="nondim", default=0.3333)
419  call get_param(param_file, mdl, "MIN_ZBOT_ITIDES", cs%min_zbot_itides, &
420  "Turn off internal tidal dissipation when the total "//&
421  "ocean depth is less than this value.", units="m", default=0.0, scale=us%m_to_Z)
422  endif
423 
424  if ( (cs%Int_tide_dissipation .or. cs%Lee_wave_dissipation) .and. &
425  .not. cs%use_CVMix_tidal) then
426 
427  call safe_alloc_ptr(cs%Nb,isd,ied,jsd,jed)
428  call safe_alloc_ptr(cs%h2,isd,ied,jsd,jed)
429  call safe_alloc_ptr(cs%TKE_itidal,isd,ied,jsd,jed)
430  call safe_alloc_ptr(cs%mask_itidal,isd,ied,jsd,jed) ; cs%mask_itidal(:,:) = 1.0
431 
432  call get_param(param_file, mdl, "KAPPA_ITIDES", cs%kappa_itides, &
433  "A topographic wavenumber used with INT_TIDE_DISSIPATION. "//&
434  "The default is 2pi/10 km, as in St.Laurent et al. 2002.", &
435  units="m-1", default=8.e-4*atan(1.0), scale=us%Z_to_m)
436 
437  call get_param(param_file, mdl, "UTIDE", cs%utide, &
438  "The constant tidal amplitude used with INT_TIDE_DISSIPATION.", &
439  units="m s-1", default=0.0, scale=us%m_to_Z*us%T_to_s)
440  call safe_alloc_ptr(cs%tideamp,is,ie,js,je) ; cs%tideamp(:,:) = cs%utide
441 
442  call get_param(param_file, mdl, "KAPPA_H2_FACTOR", cs%kappa_h2_factor, &
443  "A scaling factor for the roughness amplitude with "//&
444  "INT_TIDE_DISSIPATION.", units="nondim", default=1.0)
445  call get_param(param_file, mdl, "TKE_ITIDE_MAX", cs%TKE_itide_max, &
446  "The maximum internal tide energy source available to mix "//&
447  "above the bottom boundary layer with INT_TIDE_DISSIPATION.", &
448  units="W m-2", default=1.0e3, scale=us%W_m2_to_RZ3_T3)
449 
450  call get_param(param_file, mdl, "READ_TIDEAMP", read_tideamp, &
451  "If true, read a file (given by TIDEAMP_FILE) containing "//&
452  "the tidal amplitude with INT_TIDE_DISSIPATION.", default=.false.)
453  if (read_tideamp) then
454  if (cs%use_CVMix_tidal) then
455  call mom_error(fatal, "tidal_mixing_init: Tidal amplitude files are "// &
456  "not compatible with CVMix tidal mixing. ")
457  endif
458  call get_param(param_file, mdl, "TIDEAMP_FILE", tideamp_file, &
459  "The path to the file containing the spatially varying "//&
460  "tidal amplitudes with INT_TIDE_DISSIPATION.", default="tideamp.nc")
461  filename = trim(cs%inputdir) // trim(tideamp_file)
462  call log_param(param_file, mdl, "INPUTDIR/TIDEAMP_FILE", filename)
463  call mom_read_data(filename, 'tideamp', cs%tideamp, g%domain, timelevel=1, scale=us%m_to_Z*us%T_to_s)
464  endif
465 
466  call get_param(param_file, mdl, "H2_FILE", h2_file, &
467  "The path to the file containing the sub-grid-scale "//&
468  "topographic roughness amplitude with INT_TIDE_DISSIPATION.", &
469  fail_if_missing=(.not.cs%use_CVMix_tidal))
470  filename = trim(cs%inputdir) // trim(h2_file)
471  call log_param(param_file, mdl, "INPUTDIR/H2_FILE", filename)
472  call mom_read_data(filename, 'h2', cs%h2, g%domain, timelevel=1, scale=us%m_to_Z**2)
473 
474  call get_param(param_file, mdl, "FRACTIONAL_ROUGHNESS_MAX", max_frac_rough, &
475  "The maximum topographic roughness amplitude as a fraction of the mean depth, "//&
476  "or a negative value for no limitations on roughness.", &
477  units="nondim", default=0.1)
478 
479  do j=js,je ; do i=is,ie
480  if (g%bathyT(i,j) < cs%min_zbot_itides) cs%mask_itidal(i,j) = 0.0
481  cs%tideamp(i,j) = cs%tideamp(i,j) * cs%mask_itidal(i,j) * g%mask2dT(i,j)
482 
483  ! Restrict rms topo to a fraction (often 10 percent) of the column depth.
484  if (cs%answers_2018 .and. (max_frac_rough >= 0.0)) then
485  hamp = min(max_frac_rough*g%bathyT(i,j), sqrt(cs%h2(i,j)))
486  cs%h2(i,j) = hamp*hamp
487  else
488  if (max_frac_rough >= 0.0) &
489  cs%h2(i,j) = min((max_frac_rough*g%bathyT(i,j))**2, cs%h2(i,j))
490  endif
491 
492  utide = cs%tideamp(i,j)
493  ! Compute the fixed part of internal tidal forcing.
494  ! The units here are [R Z3 T-2 ~> J m-2 = kg s-2] here.
495  cs%TKE_itidal(i,j) = 0.5 * cs%kappa_h2_factor * gv%Rho0 * &
496  cs%kappa_itides * cs%h2(i,j) * utide*utide
497  enddo ; enddo
498 
499  endif
500 
501  if (cs%Lee_wave_dissipation) then
502 
503  call get_param(param_file, mdl, "NIKURASHIN_TKE_INPUT_FILE",niku_tke_input_file, &
504  "The path to the file containing the TKE input from lee "//&
505  "wave driven mixing. Used with LEE_WAVE_DISSIPATION.", &
506  fail_if_missing=.true.)
507  call get_param(param_file, mdl, "NIKURASHIN_SCALE",niku_scale, &
508  "A non-dimensional factor by which to scale the lee-wave "//&
509  "driven TKE input. Used with LEE_WAVE_DISSIPATION.", &
510  units="nondim", default=1.0)
511 
512  filename = trim(cs%inputdir) // trim(niku_tke_input_file)
513  call log_param(param_file, mdl, "INPUTDIR/NIKURASHIN_TKE_INPUT_FILE", &
514  filename)
515  call safe_alloc_ptr(cs%TKE_Niku,is,ie,js,je) ; cs%TKE_Niku(:,:) = 0.0
516  call mom_read_data(filename, 'TKE_input', cs%TKE_Niku, g%domain, timelevel=1, & ! ??? timelevel -aja
517  scale=niku_scale*us%W_m2_to_RZ3_T3)
518 
519  call get_param(param_file, mdl, "GAMMA_NIKURASHIN",cs%Gamma_lee, &
520  "The fraction of the lee wave energy that is dissipated "//&
521  "locally with LEE_WAVE_DISSIPATION.", units="nondim", &
522  default=0.3333)
523  call get_param(param_file, mdl, "DECAY_SCALE_FACTOR_LEE",cs%Decay_scale_factor_lee, &
524  "Scaling for the vertical decay scaleof the local "//&
525  "dissipation of lee waves dissipation.", units="nondim", &
526  default=1.0)
527  else
528  cs%Decay_scale_factor_lee = -9.e99 ! This should never be used if CS%Lee_wave_dissipation = False
529  endif
530 
531  ! Configure CVMix
532  if (cs%use_CVMix_tidal) then
533 
534  ! Read in CVMix params
535  !call openParameterBlock(param_file,'CVMix_TIDAL')
536  call get_param(param_file, mdl, "TIDAL_MAX_COEF", cs%tidal_max_coef, &
537  "largest acceptable value for tidal diffusivity", &
538  units="m^2/s", default=50e-4) ! the default is 50e-4 in CVMix, 100e-4 in POP.
539  call get_param(param_file, mdl, "TIDAL_DISS_LIM_TC", cs%tidal_diss_lim_tc, &
540  "Min allowable depth for dissipation for tidal-energy-constituent data. "//&
541  "No dissipation contribution is applied above TIDAL_DISS_LIM_TC.", &
542  units="m", default=0.0, scale=us%m_to_Z)
543  call get_param(param_file, mdl, "TIDAL_ENERGY_FILE",tidal_energy_file, &
544  "The path to the file containing tidal energy "//&
545  "dissipation. Used with CVMix tidal mixing schemes.", &
546  fail_if_missing=.true.)
547  call get_param(param_file, mdl, 'MIN_THICKNESS', cs%min_thickness, default=0.001, &
548  do_not_log=.true.)
549  call get_param(param_file, mdl, "PRANDTL_TIDAL", prandtl_tidal, &
550  "Prandtl number used by CVMix tidal mixing schemes "//&
551  "to convert vertical diffusivities into viscosities.", &
552  units="nondim", default=1.0, &
553  do_not_log=.true.)
554  call cvmix_put(cs%CVMix_glb_params,'Prandtl',prandtl_tidal)
555 
556  tidal_energy_file = trim(cs%inputdir) // trim(tidal_energy_file)
557  call get_param(param_file, mdl, "TIDAL_ENERGY_TYPE",tidal_energy_type, &
558  "The type of input tidal energy flux dataset. Valid values are"//&
559  "\t Jayne\n"//&
560  "\t ER03 \n",&
561  fail_if_missing=.true.)
562  ! Check whether tidal energy input format and CVMix tidal mixing scheme are consistent
563  if ( .not. ( &
564  (uppercase(tidal_energy_type(1:4)).eq.'JAYN' .and. cs%CVMix_tidal_scheme.eq.simmons).or. &
565  (uppercase(tidal_energy_type(1:4)).eq.'ER03' .and. cs%CVMix_tidal_scheme.eq.schmittner) ) )then
566  call mom_error(fatal, "tidal_mixing_init: Tidal energy file type ("//&
567  trim(tidal_energy_type)//") is incompatible with CVMix tidal "//&
568  " mixing scheme: "//trim(cvmix_tidal_scheme_str) )
569  endif
570  cvmix_tidal_scheme_str = lowercase(cvmix_tidal_scheme_str)
571 
572  ! Set up CVMix
573  call cvmix_init_tidal(cvmix_tidal_params_user = cs%CVMix_tidal_params, &
574  mix_scheme = cvmix_tidal_scheme_str, &
575  efficiency = cs%Mu_itides, &
576  vertical_decay_scale = cs%int_tide_decay_scale*us%Z_to_m, &
577  max_coefficient = cs%tidal_max_coef, &
578  local_mixing_frac = cs%Gamma_itides, &
579  depth_cutoff = cs%min_zbot_itides*us%Z_to_m)
580 
581  call read_tidal_energy(g, us, tidal_energy_type, tidal_energy_file, cs)
582 
583  !call closeParameterBlock(param_file)
584 
585  endif ! CVMix on
586 
587  ! Register Diagnostics fields
588 
589  if (cs%Int_tide_dissipation .or. cs%Lee_wave_dissipation .or. &
590  cs%Lowmode_itidal_dissipation) then
591 
592  cs%id_Kd_itidal = register_diag_field('ocean_model','Kd_itides',diag%axesTi,time, &
593  'Internal Tide Driven Diffusivity', 'm2 s-1', conversion=us%Z2_T_to_m2_s)
594 
595  if (cs%use_CVMix_tidal) then
596  cs%id_N2_int = register_diag_field('ocean_model','N2_int',diag%axesTi,time, &
597  'Bouyancy frequency squared, at interfaces', 's-2', conversion=us%s_to_T**2)
598  !> TODO: add units
599  cs%id_Simmons_coeff = register_diag_field('ocean_model','Simmons_coeff',diag%axesT1,time, &
600  'time-invariant portion of the tidal mixing coefficient using the Simmons', '')
601  cs%id_Schmittner_coeff = register_diag_field('ocean_model','Schmittner_coeff',diag%axesTL,time, &
602  'time-invariant portion of the tidal mixing coefficient using the Schmittner', '')
603  cs%id_tidal_qe_md = register_diag_field('ocean_model','tidal_qe_md',diag%axesTL,time, &
604  'input tidal energy dissipated locally interpolated to model vertical coordinates', '')
605  cs%id_vert_dep = register_diag_field('ocean_model','vert_dep',diag%axesTi,time, &
606  'vertical deposition function needed for Simmons et al tidal mixing', '')
607 
608  else
609  cs%id_TKE_itidal = register_diag_field('ocean_model','TKE_itidal',diag%axesT1,time, &
610  'Internal Tide Driven Turbulent Kinetic Energy', &
611  'W m-2', conversion=us%RZ3_T3_to_W_m2)
612  cs%id_Nb = register_diag_field('ocean_model','Nb',diag%axesT1,time, &
613  'Bottom Buoyancy Frequency', 's-1', conversion=us%s_to_T)
614 
615  cs%id_Kd_lowmode = register_diag_field('ocean_model','Kd_lowmode',diag%axesTi,time, &
616  'Internal Tide Driven Diffusivity (from propagating low modes)', &
617  'm2 s-1', conversion=us%Z2_T_to_m2_s)
618 
619  cs%id_Fl_itidal = register_diag_field('ocean_model','Fl_itides',diag%axesTi,time, &
620  'Vertical flux of tidal turbulent dissipation', &
621  'm3 s-3', conversion=(us%Z_to_m**3*us%s_to_T**3))
622 
623  cs%id_Fl_lowmode = register_diag_field('ocean_model','Fl_lowmode',diag%axesTi,time, &
624  'Vertical flux of tidal turbulent dissipation (from propagating low modes)', &
625  'm3 s-3', conversion=(us%Z_to_m**3*us%s_to_T**3))
626 
627  cs%id_Polzin_decay_scale = register_diag_field('ocean_model','Polzin_decay_scale',diag%axesT1,time, &
628  'Vertical decay scale for the tidal turbulent dissipation with Polzin scheme', &
629  'm', conversion=us%Z_to_m)
630 
631  cs%id_Polzin_decay_scale_scaled = register_diag_field('ocean_model', &
632  'Polzin_decay_scale_scaled', diag%axesT1, time, &
633  'Vertical decay scale for the tidal turbulent dissipation with Polzin scheme, '// &
634  'scaled by N2_bot/N2_meanz', 'm', conversion=us%Z_to_m)
635 
636  cs%id_N2_bot = register_diag_field('ocean_model','N2_b',diag%axesT1,time, &
637  'Bottom Buoyancy frequency squared', 's-2', conversion=us%s_to_T**2)
638 
639  cs%id_N2_meanz = register_diag_field('ocean_model','N2_meanz', diag%axesT1, time, &
640  'Buoyancy frequency squared averaged over the water column', 's-2', conversion=us%s_to_T**2)
641 
642  cs%id_Kd_Itidal_Work = register_diag_field('ocean_model','Kd_Itidal_Work',diag%axesTL,time, &
643  'Work done by Internal Tide Diapycnal Mixing', &
644  'W m-2', conversion=us%RZ3_T3_to_W_m2)
645 
646  cs%id_Kd_Niku_Work = register_diag_field('ocean_model','Kd_Nikurashin_Work',diag%axesTL,time, &
647  'Work done by Nikurashin Lee Wave Drag Scheme', &
648  'W m-2', conversion=us%RZ3_T3_to_W_m2)
649 
650  cs%id_Kd_Lowmode_Work = register_diag_field('ocean_model','Kd_Lowmode_Work',diag%axesTL,time, &
651  'Work done by Internal Tide Diapycnal Mixing (low modes)', &
652  'W m-2', conversion=us%RZ3_T3_to_W_m2)
653 
654  if (cs%Lee_wave_dissipation) then
655  cs%id_TKE_leewave = register_diag_field('ocean_model','TKE_leewave',diag%axesT1,time, &
656  'Lee wave Driven Turbulent Kinetic Energy', &
657  'W m-2', conversion=us%RZ3_T3_to_W_m2)
658  cs%id_Kd_Niku = register_diag_field('ocean_model','Kd_Nikurashin',diag%axesTi,time, &
659  'Lee Wave Driven Diffusivity', 'm2 s-1', conversion=us%Z2_T_to_m2_s)
660  endif
661  endif ! S%use_CVMix_tidal
662  endif
663 
664 end function tidal_mixing_init
665 
666 
667 !> Depending on whether or not CVMix is active, calls the associated subroutine to compute internal
668 !! tidal dissipation and to add the effect of internal-tide-driven mixing to the layer or interface
669 !! diffusivities.
670 subroutine calculate_tidal_mixing(h, N2_bot, j, TKE_to_Kd, max_TKE, G, GV, US, CS, &
671  N2_lay, N2_int, Kd_lay, Kd_int, Kd_max, Kv)
672  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
673  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
674  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
675  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
676  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
677  real, dimension(SZI_(G)), intent(in) :: N2_bot !< The near-bottom squared buoyancy
678  !! frequency [T-2 ~> s-2].
679  real, dimension(SZI_(G),SZK_(G)), intent(in) :: N2_lay !< The squared buoyancy frequency of the
680  !! layers [T-2 ~> s-2].
681  real, dimension(SZI_(G),SZK_(G)+1), intent(in) :: N2_int !< The squared buoyancy frequency at the
682  !! interfaces [T-2 ~> s-2].
683  integer, intent(in) :: j !< The j-index to work on
684  real, dimension(SZI_(G),SZK_(G)), intent(in) :: TKE_to_Kd !< The conversion rate between the TKE
685  !! dissipated within a layer and the
686  !! diapycnal diffusivity within that layer,
687  !! usually (~Rho_0 / (G_Earth * dRho_lay))
688  !! [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
689  real, dimension(SZI_(G),SZK_(G)), intent(in) :: max_TKE !< The energy required to for a layer to entrain
690  !! to its maximum realizable thickness [Z3 T-3 ~> m3 s-3]
691  type(tidal_mixing_cs), pointer :: CS !< The control structure for this module
692  real, dimension(SZI_(G),SZK_(G)), &
693  optional, intent(inout) :: Kd_lay !< The diapycnal diffusivity in layers [Z2 T-1 ~> m2 s-1].
694  real, dimension(SZI_(G),SZK_(G)+1), &
695  optional, intent(inout) :: Kd_int !< The diapycnal diffusivity at interfaces,
696  !! [Z2 T-1 ~> m2 s-1].
697  real, intent(in) :: Kd_max !< The maximum increment for diapycnal
698  !! diffusivity due to TKE-based processes,
699  !! [Z2 T-1 ~> m2 s-1].
700  !! Set this to a negative value to have no limit.
701  real, dimension(:,:,:), pointer :: Kv !< The "slow" vertical viscosity at each interface
702  !! (not layer!) [Z2 T-1 ~> m2 s-1].
703 
704  if (cs%Int_tide_dissipation .or. cs%Lee_wave_dissipation .or. cs%Lowmode_itidal_dissipation) then
705  if (cs%use_CVMix_tidal) then
706  call calculate_cvmix_tidal(h, j, g, gv, us, cs, n2_int, kd_lay, kd_int, kv)
707  else
708  call add_int_tide_diffusivity(h, n2_bot, j, tke_to_kd, max_tke, &
709  g, gv, us, cs, n2_lay, kd_lay, kd_int, kd_max)
710  endif
711  endif
712 end subroutine calculate_tidal_mixing
713 
714 
715 !> Calls the CVMix routines to compute tidal dissipation and to add the effect of internal-tide-driven
716 !! mixing to the interface diffusivities.
717 subroutine calculate_cvmix_tidal(h, j, G, GV, US, CS, N2_int, Kd_lay, Kd_int, Kv)
718  integer, intent(in) :: j !< The j-index to work on
719  type(ocean_grid_type), intent(in) :: G !< Grid structure.
720  type(verticalgrid_type), intent(in) :: GV !< ocean vertical grid structure
721  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
722  type(tidal_mixing_cs), pointer :: CS !< This module's control structure.
723  real, dimension(SZI_(G),SZK_(G)+1), intent(in) :: N2_int !< The squared buoyancy
724  !! frequency at the interfaces [T-2 ~> s-2].
725  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
726  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
727  real, dimension(SZI_(G),SZK_(G)), &
728  optional, intent(inout) :: Kd_lay!< The diapycnal diffusivity in the layers [Z2 T-1 ~> m2 s-1].
729  real, dimension(SZI_(G),SZK_(G)+1), &
730  optional, intent(inout) :: Kd_int!< The diapycnal diffusivity at interfaces [Z2 T-1 ~> m2 s-1].
731  real, dimension(:,:,:), pointer :: Kv !< The "slow" vertical viscosity at each interface
732  !! (not layer!) [Z2 T-1 ~> m2 s-1].
733  ! Local variables
734  real, dimension(SZK_(G)+1) :: Kd_tidal ! tidal diffusivity [m2 s-1]
735  real, dimension(SZK_(G)+1) :: Kv_tidal ! tidal viscosity [m2 s-1]
736  real, dimension(SZK_(G)+1) :: vert_dep ! vertical deposition
737  real, dimension(SZK_(G)+1) :: iFaceHeight ! Height of interfaces [m]
738  real, dimension(SZK_(G)+1) :: SchmittnerSocn
739  real, dimension(SZK_(G)) :: cellHeight ! Height of cell centers [m]
740  real, dimension(SZK_(G)) :: tidal_qe_md ! Tidal dissipation energy interpolated from 3d input
741  ! to model coordinates
742  real, dimension(SZK_(G)+1) :: N2_int_i ! De-scaled interface buoyancy frequency [s-2]
743  real, dimension(SZK_(G)) :: Schmittner_coeff
744  real, dimension(SZK_(G)) :: h_m ! Cell thickness [m]
745  real, allocatable, dimension(:,:) :: exp_hab_zetar
746 
747  integer :: i, k, is, ie
748  real :: dh, hcorr, Simmons_coeff
749  real, parameter :: rho_fw = 1000.0 ! fresh water density [kg/m^3]
750  ! TODO: when coupled, get this from CESM (SHR_CONST_RHOFW)
751  type(tidal_mixing_diags), pointer :: dd => null()
752 
753  is = g%isc ; ie = g%iec
754  dd => cs%dd
755 
756  select case (cs%CVMix_tidal_scheme)
757  case (simmons)
758  do i=is,ie
759 
760  if (g%mask2dT(i,j)<1) cycle
761 
762  ifaceheight = 0.0 ! BBL is all relative to the surface
763  hcorr = 0.0
764  do k=1,g%ke
765  ! cell center and cell bottom in meters (negative values in the ocean)
766  dh = h(i,j,k) * gv%H_to_m ! Nominal thickness to use for increment, rescaled to m for use by CVMix.
767  dh = dh + hcorr ! Take away the accumulated error (could temporarily make dh<0)
768  hcorr = min( dh - cs%min_thickness, 0. ) ! If inflating then hcorr<0
769  dh = max( dh, cs%min_thickness ) ! Limit increment dh>=min_thickness
770  cellheight(k) = ifaceheight(k) - 0.5 * dh
771  ifaceheight(k+1) = ifaceheight(k) - dh
772  enddo
773 
774  call cvmix_compute_simmons_invariant( nlev = g%ke, &
775  energy_flux = cs%tidal_qe_2d(i,j), &
776  rho = rho_fw, &
777  simmonscoeff = simmons_coeff, &
778  vertdep = vert_dep, &
779  zw = ifaceheight, &
780  zt = cellheight, &
781  cvmix_tidal_params_user = cs%CVMix_tidal_params)
782 
783  ! Since we pass tidal_qe_2d=(CS%Gamma_itides)*tidal_energy_flux_2d, and not tidal_energy_flux_2d in
784  ! above subroutine call, we divide Simmons_coeff by CS%Gamma_itides as a corrective step:
785  ! TODO: (CS%Gamma_itides)*tidal_energy_flux_2d is unnecessary, directly use tidal_energy_flux_2d
786  simmons_coeff = simmons_coeff / cs%Gamma_itides
787 
788 
789  ! XXX: Temporary de-scaling of N2_int(i,:) into a temporary variable
790  do k=1,g%ke+1
791  n2_int_i(k) = us%s_to_T**2 * n2_int(i,k)
792  enddo
793 
794  call cvmix_coeffs_tidal( mdiff_out = kv_tidal, &
795  tdiff_out = kd_tidal, &
796  nsqr = n2_int_i, &
797  oceandepth = -ifaceheight(g%ke+1),&
798  simmonscoeff = simmons_coeff, &
799  vert_dep = vert_dep, &
800  nlev = g%ke, &
801  max_nlev = g%ke, &
802  cvmix_params = cs%CVMix_glb_params, &
803  cvmix_tidal_params_user = cs%CVMix_tidal_params)
804 
805  ! Update diffusivity
806  if (present(kd_lay)) then
807  do k=1,g%ke
808  kd_lay(i,k) = kd_lay(i,k) + 0.5 * us%m2_s_to_Z2_T * (kd_tidal(k) + kd_tidal(k+1))
809  enddo
810  endif
811  if (present(kd_int)) then
812  do k=1,g%ke+1
813  kd_int(i,k) = kd_int(i,k) + (us%m2_s_to_Z2_T * kd_tidal(k))
814  enddo
815  endif
816  ! Update viscosity with the proper unit conversion.
817  if (associated(kv)) then
818  do k=1,g%ke+1
819  kv(i,j,k) = kv(i,j,k) + us%m2_s_to_Z2_T * kv_tidal(k) ! Rescale from m2 s-1 to Z2 T-1.
820  enddo
821  endif
822 
823  ! diagnostics
824  if (associated(dd%Kd_itidal)) then
825  dd%Kd_itidal(i,j,:) = us%m2_s_to_Z2_T*kd_tidal(:)
826  endif
827  if (associated(dd%N2_int)) then
828  dd%N2_int(i,j,:) = n2_int(i,:)
829  endif
830  if (associated(dd%Simmons_coeff_2d)) then
831  dd%Simmons_coeff_2d(i,j) = simmons_coeff
832  endif
833  if (associated(dd%vert_dep_3d)) then
834  dd%vert_dep_3d(i,j,:) = vert_dep(:)
835  endif
836 
837  enddo ! i=is,ie
838 
839  case (schmittner)
840 
841  ! TODO: correct exp_hab_zetar shapes in CVMix_compute_Schmittner_invariant
842  ! and CVMix_compute_SchmittnerCoeff low subroutines
843 
844  allocate(exp_hab_zetar(g%ke+1,g%ke+1))
845 
846  do i=is,ie
847 
848  if (g%mask2dT(i,j)<1) cycle
849 
850  ifaceheight = 0.0 ! BBL is all relative to the surface
851  hcorr = 0.0
852  do k=1,g%ke
853  h_m(k) = h(i,j,k)*gv%H_to_m ! Rescale thicknesses to m for use by CVmix.
854  ! cell center and cell bottom in meters (negative values in the ocean)
855  dh = h_m(k) + hcorr ! Nominal thickness less the accumulated error (could temporarily make dh<0)
856  hcorr = min( dh - cs%min_thickness, 0. ) ! If inflating then hcorr<0
857  dh = max( dh, cs%min_thickness ) ! Limit increment dh>=min_thickness
858  cellheight(k) = ifaceheight(k) - 0.5 * dh
859  ifaceheight(k+1) = ifaceheight(k) - dh
860  enddo
861 
862  schmittnersocn = 0.0 ! TODO: compute this
863 
864  ! form the time-invariant part of Schmittner coefficient term
865  call cvmix_compute_schmittner_invariant(nlev = g%ke, &
866  vertdep = vert_dep, &
867  efficiency = cs%Mu_itides, &
868  rho = rho_fw, &
869  exp_hab_zetar = exp_hab_zetar, &
870  zw = ifaceheight, &
871  cvmix_tidal_params_user = cs%CVMix_tidal_params)
872  !TODO: in above call, there is no need to pass efficiency, since it gets
873  ! passed via CVMix_init_tidal and stored in CVMix_tidal_params. Change
874  ! CVMix API to prevent this redundancy.
875 
876  ! remap from input z coordinate to model coordinate:
877  tidal_qe_md = 0.0
878  call remapping_core_h(cs%remap_cs, size(cs%h_src), cs%h_src, cs%tidal_qe_3d_in(i,j,:), &
879  g%ke, h_m, tidal_qe_md)
880 
881  ! form the Schmittner coefficient that is based on 3D q*E, which is formed from
882  ! summing q_i*TidalConstituent_i over the number of constituents.
883  call cvmix_compute_schmittnercoeff( nlev = g%ke, &
884  energy_flux = tidal_qe_md(:), &
885  rho = rho_fw, &
886  schmittnercoeff = schmittner_coeff, &
887  exp_hab_zetar = exp_hab_zetar, &
888  cvmix_tidal_params_user = cs%CVMix_tidal_params)
889 
890  ! XXX: Temporary de-scaling of N2_int(i,:) into a temporary variable
891  do k=1,g%ke+1
892  n2_int_i(k) = us%s_to_T**2 * n2_int(i,k)
893  enddo
894 
895  call cvmix_coeffs_tidal_schmittner( mdiff_out = kv_tidal, &
896  tdiff_out = kd_tidal, &
897  nsqr = n2_int_i, &
898  oceandepth = -ifaceheight(g%ke+1), &
899  vert_dep = vert_dep, &
900  nlev = g%ke, &
901  max_nlev = g%ke, &
902  schmittnercoeff = schmittner_coeff, &
903  schmittnersouthernocean = schmittnersocn, &
904  cvmix_params = cs%CVMix_glb_params, &
905  cvmix_tidal_params_user = cs%CVMix_tidal_params)
906 
907  ! Update diffusivity
908  if (present(kd_lay)) then
909  do k=1,g%ke
910  kd_lay(i,k) = kd_lay(i,k) + 0.5 * us%m2_s_to_Z2_T * (kd_tidal(k) + kd_tidal(k+1))
911  enddo
912  endif
913  if (present(kd_int)) then
914  do k=1,g%ke+1
915  kd_int(i,k) = kd_int(i,k) + (us%m2_s_to_Z2_T * kd_tidal(k))
916  enddo
917  endif
918 
919  ! Update viscosity
920  if (associated(kv)) then
921  do k=1,g%ke+1
922  kv(i,j,k) = kv(i,j,k) + us%m2_s_to_Z2_T * kv_tidal(k) ! Rescale from m2 s-1 to Z2 T-1.
923  enddo
924  endif
925 
926  ! diagnostics
927  if (associated(dd%Kd_itidal)) then
928  dd%Kd_itidal(i,j,:) = us%m2_s_to_Z2_T*kd_tidal(:)
929  endif
930  if (associated(dd%N2_int)) then
931  dd%N2_int(i,j,:) = n2_int(i,:)
932  endif
933  if (associated(dd%Schmittner_coeff_3d)) then
934  dd%Schmittner_coeff_3d(i,j,:) = schmittner_coeff(:)
935  endif
936  if (associated(dd%tidal_qe_md)) then
937  dd%tidal_qe_md(i,j,:) = tidal_qe_md(:)
938  endif
939  if (associated(dd%vert_dep_3d)) then
940  dd%vert_dep_3d(i,j,:) = vert_dep(:)
941  endif
942  enddo ! i=is,ie
943 
944  deallocate(exp_hab_zetar)
945 
946  case default
947  call mom_error(fatal, "tidal_mixing_init: Unrecognized setting "// &
948  "#define CVMIX_TIDAL_SCHEME found in input file.")
949  end select
950 
951 end subroutine calculate_cvmix_tidal
952 
953 
954 !> This subroutine adds the effect of internal-tide-driven mixing to the layer diffusivities.
955 !! The mechanisms considered are (1) local dissipation of internal waves generated by the
956 !! barotropic flow ("itidal"), (2) local dissipation of internal waves generated by the propagating
957 !! low modes (rays) of the internal tide ("lowmode"), and (3) local dissipation of internal lee waves.
958 !! Will eventually need to add diffusivity due to other wave-breaking processes (e.g. Bottom friction,
959 !! Froude-number-depending breaking, PSI, etc.).
960 subroutine add_int_tide_diffusivity(h, N2_bot, j, TKE_to_Kd, max_TKE, G, GV, US, CS, &
961  N2_lay, Kd_lay, Kd_int, Kd_max)
962  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
963  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
964  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
965  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
966  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
967  real, dimension(SZI_(G)), intent(in) :: N2_bot !< The near-bottom squared buoyancy frequency
968  !! frequency [T-2 ~> s-2].
969  real, dimension(SZI_(G),SZK_(G)), intent(in) :: N2_lay !< The squared buoyancy frequency of the
970  !! layers [T-2 ~> s-2].
971  integer, intent(in) :: j !< The j-index to work on
972  real, dimension(SZI_(G),SZK_(G)), intent(in) :: TKE_to_Kd !< The conversion rate between the TKE
973  !! dissipated within a layer and the
974  !! diapycnal diffusivity within that layer,
975  !! usually (~Rho_0 / (G_Earth * dRho_lay))
976  !! [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
977  real, dimension(SZI_(G),SZK_(G)), intent(in) :: max_TKE !< The energy required to for a layer to entrain
978  !! to its maximum realizable thickness [Z3 T-3 ~> m3 s-3]
979  type(tidal_mixing_cs), pointer :: CS !< The control structure for this module
980  real, dimension(SZI_(G),SZK_(G)), &
981  optional, intent(inout) :: Kd_lay !< The diapycnal diffusivity in layers [Z2 T-1 ~> m2 s-1].
982  real, dimension(SZI_(G),SZK_(G)+1), &
983  optional, intent(inout) :: Kd_int !< The diapycnal diffusivity at interfaces
984  !! [Z2 T-1 ~> m2 s-1].
985  real, intent(in) :: Kd_max !< The maximum increment for diapycnal
986  !! diffusivity due to TKE-based processes
987  !! [Z2 T-1 ~> m2 s-1].
988  !! Set this to a negative value to have no limit.
989 
990  ! local
991 
992  real, dimension(SZI_(G)) :: &
993  htot, & ! total thickness above or below a layer, or the
994  ! integrated thickness in the BBL [Z ~> m].
995  htot_wkb, & ! WKB scaled distance from top to bottom [Z ~> m].
996  tke_itidal_bot, & ! internal tide TKE at ocean bottom [Z3 T-3 ~> m3 s-3]
997  tke_niku_bot, & ! lee-wave TKE at ocean bottom [Z3 T-3 ~> m3 s-3]
998  tke_lowmode_bot, & ! internal tide TKE at ocean bottom lost from all remote low modes [Z3 T-3 ~> m3 s-3] (BDM)
999  inv_int, & ! inverse of TKE decay for int tide over the depth of the ocean [nondim]
1000  inv_int_lee, & ! inverse of TKE decay for lee waves over the depth of the ocean [nondim]
1001  inv_int_low, & ! inverse of TKE decay for low modes over the depth of the ocean [nondim] (BDM)
1002  z0_polzin, & ! TKE decay scale in Polzin formulation [Z ~> m].
1003  z0_polzin_scaled, & ! TKE decay scale in Polzin formulation [Z ~> m].
1004  ! multiplied by N2_bot/N2_meanz to be coherent with the WKB scaled z
1005  ! z*=int(N2/N2_bot) * N2_bot/N2_meanz = int(N2/N2_meanz)
1006  ! z0_Polzin_scaled = z0_Polzin * N2_bot/N2_meanz
1007  n2_meanz, & ! vertically averaged squared buoyancy frequency [T-2] for WKB scaling
1008  tke_itidal_rem, & ! remaining internal tide TKE (from barotropic source) [Z3 T-3 ~> m3 s-3]
1009  tke_niku_rem, & ! remaining lee-wave TKE [Z3 T-3 ~> m3 s-3]
1010  tke_lowmode_rem, & ! remaining internal tide TKE (from propagating low mode source) [Z3 T-3 ~> m3 s-3] (BDM)
1011  tke_frac_top, & ! fraction of bottom TKE that should appear at top of a layer [nondim]
1012  tke_frac_top_lee, & ! fraction of bottom TKE that should appear at top of a layer [nondim]
1013  tke_frac_top_lowmode, &
1014  ! fraction of bottom TKE that should appear at top of a layer [nondim] (BDM)
1015  z_from_bot, & ! distance from bottom [Z ~> m].
1016  z_from_bot_wkb ! WKB scaled distance from bottom [Z ~> m].
1017 
1018  real :: I_rho0 ! Inverse of the Boussinesq reference density, i.e. 1 / RHO0 [R-1 ~> m3 kg-1]
1019  real :: Kd_add ! diffusivity to add in a layer [Z2 T-1 ~> m2 s-1].
1020  real :: TKE_itide_lay ! internal tide TKE imparted to a layer (from barotropic) [Z3 T-3 ~> m3 s-3]
1021  real :: TKE_Niku_lay ! lee-wave TKE imparted to a layer [Z3 T-3 ~> m3 s-3]
1022  real :: TKE_lowmode_lay ! internal tide TKE imparted to a layer (from low mode) [Z3 T-3 ~> m3 s-3] (BDM)
1023  real :: frac_used ! fraction of TKE that can be used in a layer [nondim]
1024  real :: Izeta ! inverse of TKE decay scale [Z-1 ~> m-1].
1025  real :: Izeta_lee ! inverse of TKE decay scale for lee waves [Z-1 ~> m-1].
1026  real :: z0Ps_num ! The numerator of the unlimited z0_Polzin_scaled [Z T-3 ~> m s-3].
1027  real :: z0Ps_denom ! The denominator of the unlimited z0_Polzin_scaled [T-3 ~> s-3].
1028  real :: z0_psl ! temporary variable [Z ~> m].
1029  real :: TKE_lowmode_tot ! TKE from all low modes [R Z3 T-3 ~> W m-2] (BDM)
1030 
1031  logical :: use_Polzin, use_Simmons
1032  character(len=160) :: mesg ! The text of an error message
1033  integer :: i, k, is, ie, nz
1034  integer :: a, fr, m
1035  type(tidal_mixing_diags), pointer :: dd => null()
1036 
1037  is = g%isc ; ie = g%iec ; nz = g%ke
1038  dd => cs%dd
1039 
1040  if (.not.(cs%Int_tide_dissipation .or. cs%Lee_wave_dissipation)) return
1041 
1042  do i=is,ie ; htot(i) = 0.0 ; inv_int(i) = 0.0 ; inv_int_lee(i) = 0.0 ; inv_int_low(i) = 0.0 ;enddo
1043  do k=1,nz ; do i=is,ie
1044  htot(i) = htot(i) + gv%H_to_Z*h(i,j,k)
1045  enddo ; enddo
1046 
1047  i_rho0 = 1.0 / (gv%Rho0)
1048 
1049  use_polzin = ((cs%Int_tide_dissipation .and. (cs%int_tide_profile == polzin_09)) .or. &
1050  (cs%lee_wave_dissipation .and. (cs%lee_wave_profile == polzin_09)) .or. &
1051  (cs%Lowmode_itidal_dissipation .and. (cs%int_tide_profile == polzin_09)))
1052  use_simmons = ((cs%Int_tide_dissipation .and. (cs%int_tide_profile == stlaurent_02)) .or. &
1053  (cs%lee_wave_dissipation .and. (cs%lee_wave_profile == stlaurent_02)) .or. &
1054  (cs%Lowmode_itidal_dissipation .and. (cs%int_tide_profile == stlaurent_02)))
1055 
1056  ! Calculate parameters for vertical structure of dissipation
1057  ! Simmons:
1058  if ( use_simmons ) then
1059  izeta = 1.0 / max(cs%Int_tide_decay_scale, gv%H_subroundoff*gv%H_to_Z)
1060  izeta_lee = 1.0 / max(cs%Int_tide_decay_scale*cs%Decay_scale_factor_lee, &
1061  gv%H_subroundoff*gv%H_to_Z)
1062  do i=is,ie
1063  cs%Nb(i,j) = sqrt(n2_bot(i))
1064  if (associated(dd%N2_bot)) dd%N2_bot(i,j) = n2_bot(i)
1065  if ( cs%Int_tide_dissipation ) then
1066  if (izeta*htot(i) > 1.0e-14) then ! L'Hospital's version of Adcroft's reciprocal rule.
1067  inv_int(i) = 1.0 / (1.0 - exp(-izeta*htot(i)))
1068  endif
1069  endif
1070  if ( cs%Lee_wave_dissipation ) then
1071  if (izeta_lee*htot(i) > 1.0e-14) then ! L'Hospital's version of Adcroft's reciprocal rule.
1072  inv_int_lee(i) = 1.0 / (1.0 - exp(-izeta_lee*htot(i)))
1073  endif
1074  endif
1075  if ( cs%Lowmode_itidal_dissipation) then
1076  if (izeta*htot(i) > 1.0e-14) then ! L'Hospital's version of Adcroft's reciprocal rule.
1077  inv_int_low(i) = 1.0 / (1.0 - exp(-izeta*htot(i)))
1078  endif
1079  endif
1080  z_from_bot(i) = gv%H_to_Z*h(i,j,nz)
1081  enddo
1082  endif ! Simmons
1083 
1084  ! Polzin:
1085  if ( use_polzin ) then
1086  ! WKB scaling of the vertical coordinate
1087  do i=is,ie ; n2_meanz(i) = 0.0 ; enddo
1088  do k=1,nz ; do i=is,ie
1089  n2_meanz(i) = n2_meanz(i) + n2_lay(i,k) * gv%H_to_Z * h(i,j,k)
1090  enddo ; enddo
1091  do i=is,ie
1092  n2_meanz(i) = n2_meanz(i) / (htot(i) + gv%H_subroundoff*gv%H_to_Z)
1093  if (associated(dd%N2_meanz)) dd%N2_meanz(i,j) = n2_meanz(i)
1094  enddo
1095 
1096  ! WKB scaled z*(z=H) z* at the surface using the modified Polzin WKB scaling
1097  do i=is,ie ; htot_wkb(i) = htot(i) ; enddo
1098 ! do i=is,ie ; htot_WKB(i) = 0.0 ; enddo
1099 ! do k=1,nz ; do i=is,ie
1100 ! htot_WKB(i) = htot_WKB(i) + GV%H_to_Z*h(i,j,k) * N2_lay(i,k) / N2_meanz(i)
1101 ! enddo ; enddo
1102  ! htot_WKB(i) = htot(i) ! Nearly equivalent and simpler
1103 
1104  do i=is,ie
1105  cs%Nb(i,j) = sqrt(n2_bot(i))
1106  if (cs%answers_2018) then
1107  if ((cs%tideamp(i,j) > 0.0) .and. &
1108  (cs%kappa_itides**2 * cs%h2(i,j) * cs%Nb(i,j)**3 > 1.0e-14*us%T_to_s**3) ) then
1109  z0_polzin(i) = cs%Polzin_decay_scale_factor * cs%Nu_Polzin * &
1110  cs%Nbotref_Polzin**2 * cs%tideamp(i,j) / &
1111  ( cs%kappa_itides**2 * cs%h2(i,j) * cs%Nb(i,j)**3 )
1112  if (z0_polzin(i) < cs%Polzin_min_decay_scale) &
1113  z0_polzin(i) = cs%Polzin_min_decay_scale
1114  if (n2_meanz(i) > 1.0e-14*us%T_to_s**2 ) then
1115  z0_polzin_scaled(i) = z0_polzin(i)*cs%Nb(i,j)**2 / n2_meanz(i)
1116  else
1117  z0_polzin_scaled(i) = cs%Polzin_decay_scale_max_factor * htot(i)
1118  endif
1119  if (z0_polzin_scaled(i) > (cs%Polzin_decay_scale_max_factor * htot(i)) ) &
1120  z0_polzin_scaled(i) = cs%Polzin_decay_scale_max_factor * htot(i)
1121  else
1122  z0_polzin(i) = cs%Polzin_decay_scale_max_factor * htot(i)
1123  z0_polzin_scaled(i) = cs%Polzin_decay_scale_max_factor * htot(i)
1124  endif
1125  else
1126  z0ps_num = (cs%Polzin_decay_scale_factor * cs%Nu_Polzin * cs%Nbotref_Polzin**2) * cs%tideamp(i,j)
1127  z0ps_denom = ( cs%kappa_itides**2 * cs%h2(i,j) * cs%Nb(i,j) * n2_meanz(i) )
1128  if ((cs%tideamp(i,j) > 0.0) .and. &
1129  (z0ps_num < z0ps_denom * cs%Polzin_decay_scale_max_factor * htot(i))) then
1130  z0_polzin_scaled(i) = z0ps_num / z0ps_denom
1131 
1132  if (abs(n2_meanz(i) * z0_polzin_scaled(i)) < &
1133  cs%Nb(i,j)**2 * (cs%Polzin_decay_scale_max_factor * htot(i))) then
1134  z0_polzin(i) = z0_polzin_scaled(i) * (n2_meanz(i) / cs%Nb(i,j)**2)
1135  else
1136  z0_polzin(i) = cs%Polzin_decay_scale_max_factor * htot(i)
1137  endif
1138  else
1139  z0_polzin(i) = cs%Polzin_decay_scale_max_factor * htot(i)
1140  z0_polzin_scaled(i) = cs%Polzin_decay_scale_max_factor * htot(i)
1141  endif
1142  endif
1143 
1144  if (associated(dd%Polzin_decay_scale)) &
1145  dd%Polzin_decay_scale(i,j) = z0_polzin(i)
1146  if (associated(dd%Polzin_decay_scale_scaled)) &
1147  dd%Polzin_decay_scale_scaled(i,j) = z0_polzin_scaled(i)
1148  if (associated(dd%N2_bot)) dd%N2_bot(i,j) = cs%Nb(i,j)*cs%Nb(i,j)
1149 
1150  if (cs%answers_2018) then
1151  ! These expressions use dimensional constants to avoid NaN values.
1152  if ( cs%Int_tide_dissipation .and. (cs%int_tide_profile == polzin_09) ) then
1153  if (htot_wkb(i) > 1.0e-14*us%m_to_Z) &
1154  inv_int(i) = ( z0_polzin_scaled(i) / htot_wkb(i) ) + 1.0
1155  endif
1156  if ( cs%lee_wave_dissipation .and. (cs%lee_wave_profile == polzin_09) ) then
1157  if (htot_wkb(i) > 1.0e-14*us%m_to_Z) &
1158  inv_int_lee(i) = ( z0_polzin_scaled(i)*cs%Decay_scale_factor_lee / htot_wkb(i) ) + 1.0
1159  endif
1160  if ( cs%Lowmode_itidal_dissipation .and. (cs%int_tide_profile == polzin_09) ) then
1161  if (htot_wkb(i) > 1.0e-14*us%m_to_Z) &
1162  inv_int_low(i) = ( z0_polzin_scaled(i) / htot_wkb(i) ) + 1.0
1163  endif
1164  else
1165  ! These expressions give values of Inv_int < 10^14 using a variant of Adcroft's reciprocal rule.
1166  inv_int(i) = 0.0 ; inv_int_lee(i) = 0.0 ; inv_int_low(i) = 0.0
1167  if ( cs%Int_tide_dissipation .and. (cs%int_tide_profile == polzin_09) ) then
1168  if (z0_polzin_scaled(i) < 1.0e14 * htot_wkb(i)) &
1169  inv_int(i) = ( z0_polzin_scaled(i) / htot_wkb(i) ) + 1.0
1170  endif
1171  if ( cs%lee_wave_dissipation .and. (cs%lee_wave_profile == polzin_09) ) then
1172  if (z0_polzin_scaled(i) < 1.0e14 * htot_wkb(i)) &
1173  inv_int_lee(i) = ( z0_polzin_scaled(i)*cs%Decay_scale_factor_lee / htot_wkb(i) ) + 1.0
1174  endif
1175  if ( cs%Lowmode_itidal_dissipation .and. (cs%int_tide_profile == polzin_09) ) then
1176  if (z0_polzin_scaled(i) < 1.0e14 * htot_wkb(i)) &
1177  inv_int_low(i) = ( z0_polzin_scaled(i) / htot_wkb(i) ) + 1.0
1178  endif
1179  endif
1180 
1181  z_from_bot(i) = gv%H_to_Z*h(i,j,nz)
1182  ! Use the new formulation for WKB scaling. N2 is referenced to its vertical mean.
1183  if (cs%answers_2018) then
1184  if (n2_meanz(i) > 1.0e-14*us%T_to_s**2 ) then
1185  z_from_bot_wkb(i) = gv%H_to_Z*h(i,j,nz) * n2_lay(i,nz) / n2_meanz(i)
1186  else ; z_from_bot_wkb(i) = 0 ; endif
1187  else
1188  if (gv%H_to_Z*h(i,j,nz) * n2_lay(i,nz) < n2_meanz(i) * (1.0e14 * htot_wkb(i))) then
1189  z_from_bot_wkb(i) = gv%H_to_Z*h(i,j,nz) * n2_lay(i,nz) / n2_meanz(i)
1190  else ; z_from_bot_wkb(i) = 0 ; endif
1191  endif
1192  enddo
1193  endif ! Polzin
1194 
1195  ! Calculate/get dissipation values at bottom
1196  ! Both Polzin and Simmons:
1197  do i=is,ie
1198  ! Dissipation of locally trapped internal tide (non-propagating high modes)
1199  tke_itidal_bot(i) = min(cs%TKE_itidal(i,j)*cs%Nb(i,j), cs%TKE_itide_max)
1200  if (associated(dd%TKE_itidal_used)) &
1201  dd%TKE_itidal_used(i,j) = tke_itidal_bot(i)
1202  tke_itidal_bot(i) = (i_rho0 * cs%Mu_itides * cs%Gamma_itides) * tke_itidal_bot(i)
1203  ! Dissipation of locally trapped lee waves
1204  tke_niku_bot(i) = 0.0
1205  if (cs%Lee_wave_dissipation) then
1206  tke_niku_bot(i) = (i_rho0 * cs%Mu_itides * cs%Gamma_lee) * cs%TKE_Niku(i,j)
1207  endif
1208  ! Dissipation of propagating internal tide (baroclinic low modes; rays) (BDM)
1209  tke_lowmode_tot = 0.0
1210  tke_lowmode_bot(i) = 0.0
1211  if (cs%Lowmode_itidal_dissipation) then
1212  ! get loss rate due to wave drag on low modes (already multiplied by q)
1213 
1214  ! TODO: uncomment the following call and fix it
1215  !call get_lowmode_loss(i,j,G,CS%int_tide_CSp,"WaveDrag",TKE_lowmode_tot)
1216  write (mesg,*) "========", __file__, __line__
1217  call mom_error(fatal,trim(mesg)//": this block not supported yet. (aa)")
1218 
1219  tke_lowmode_bot(i) = cs%Mu_itides * i_rho0 * tke_lowmode_tot
1220  endif
1221  ! Vertical energy flux at bottom
1222  tke_itidal_rem(i) = inv_int(i) * tke_itidal_bot(i)
1223  tke_niku_rem(i) = inv_int_lee(i) * tke_niku_bot(i)
1224  tke_lowmode_rem(i) = inv_int_low(i) * tke_lowmode_bot(i)
1225 
1226  if (associated(dd%Fl_itidal)) dd%Fl_itidal(i,j,nz) = tke_itidal_rem(i) !why is this here? BDM
1227  enddo
1228 
1229  ! Estimate the work that would be done by mixing in each layer.
1230  ! Simmons:
1231  if ( use_simmons ) then
1232  do k=nz-1,2,-1 ; do i=is,ie
1233  if (max_tke(i,k) <= 0.0) cycle
1234  z_from_bot(i) = z_from_bot(i) + gv%H_to_Z*h(i,j,k)
1235 
1236  ! Fraction of bottom flux predicted to reach top of this layer
1237  tke_frac_top(i) = inv_int(i) * exp(-izeta * z_from_bot(i))
1238  tke_frac_top_lee(i) = inv_int_lee(i) * exp(-izeta_lee * z_from_bot(i))
1239  tke_frac_top_lowmode(i) = inv_int_low(i) * exp(-izeta * z_from_bot(i))
1240 
1241  ! Actual influx at bottom of layer minus predicted outflux at top of layer to give
1242  ! predicted power expended
1243  tke_itide_lay = tke_itidal_rem(i) - tke_itidal_bot(i) * tke_frac_top(i)
1244  tke_niku_lay = tke_niku_rem(i) - tke_niku_bot(i) * tke_frac_top_lee(i)
1245  tke_lowmode_lay = tke_lowmode_rem(i) - tke_lowmode_bot(i)* tke_frac_top_lowmode(i)
1246 
1247  ! Actual power expended may be less than predicted if stratification is weak; adjust
1248  if (tke_itide_lay + tke_niku_lay + tke_lowmode_lay > (max_tke(i,k))) then
1249  frac_used = (max_tke(i,k)) / (tke_itide_lay + tke_niku_lay + tke_lowmode_lay)
1250  tke_itide_lay = frac_used * tke_itide_lay
1251  tke_niku_lay = frac_used * tke_niku_lay
1252  tke_lowmode_lay = frac_used * tke_lowmode_lay
1253  endif
1254 
1255  ! Calculate vertical flux available to bottom of layer above
1256  tke_itidal_rem(i) = tke_itidal_rem(i) - tke_itide_lay
1257  tke_niku_rem(i) = tke_niku_rem(i) - tke_niku_lay
1258  tke_lowmode_rem(i) = tke_lowmode_rem(i) - tke_lowmode_lay
1259 
1260  ! Convert power to diffusivity
1261  kd_add = tke_to_kd(i,k) * (tke_itide_lay + tke_niku_lay + tke_lowmode_lay)
1262 
1263  if (kd_max >= 0.0) kd_add = min(kd_add, kd_max)
1264  if (present(kd_lay)) then
1265  kd_lay(i,k) = kd_lay(i,k) + kd_add
1266  endif
1267 
1268  if (present(kd_int)) then
1269  kd_int(i,k) = kd_int(i,k) + 0.5 * kd_add
1270  kd_int(i,k+1) = kd_int(i,k+1) + 0.5 * kd_add
1271  endif
1272 
1273  ! diagnostics
1274  if (associated(dd%Kd_itidal)) then
1275  ! If at layers, dd%Kd_itidal is just TKE_to_Kd(i,k) * TKE_itide_lay
1276  ! The following sets the interface diagnostics.
1277  kd_add = tke_to_kd(i,k) * tke_itide_lay
1278  if (kd_max >= 0.0) kd_add = min(kd_add, kd_max)
1279  if (k>1) dd%Kd_itidal(i,j,k) = dd%Kd_itidal(i,j,k) + 0.5*kd_add
1280  if (k<nz) dd%Kd_itidal(i,j,k+1) = dd%Kd_itidal(i,j,k+1) + 0.5*kd_add
1281  endif
1282  if (associated(dd%Kd_Itidal_work)) &
1283  dd%Kd_itidal_work(i,j,k) = gv%Rho0 * tke_itide_lay
1284  if (associated(dd%Fl_itidal)) dd%Fl_itidal(i,j,k) = tke_itidal_rem(i)
1285 
1286  if (associated(dd%Kd_Niku)) then
1287  ! If at layers, dd%Kd_Niku(i,j,K) is just TKE_to_Kd(i,k) * TKE_Niku_lay
1288  ! The following sets the interface diagnostics.
1289  kd_add = tke_to_kd(i,k) * tke_niku_lay
1290  if (kd_max >= 0.0) kd_add = min(kd_add, kd_max)
1291  if (k>1) dd%Kd_Niku(i,j,k) = dd%Kd_Niku(i,j,k) + 0.5*kd_add
1292  if (k<nz) dd%Kd_Niku(i,j,k+1) = dd%Kd_Niku(i,j,k+1) + 0.5*kd_add
1293  endif
1294 ! if (associated(dd%Kd_Niku)) dd%Kd_Niku(i,j,K) = TKE_to_Kd(i,k) * TKE_Niku_lay
1295  if (associated(dd%Kd_Niku_work)) &
1296  dd%Kd_Niku_work(i,j,k) = gv%Rho0 * tke_niku_lay
1297 
1298  if (associated(dd%Kd_lowmode)) then
1299  ! If at layers, dd%Kd_lowmode is just TKE_to_Kd(i,k) * TKE_lowmode_lay
1300  ! The following sets the interface diagnostics.
1301  kd_add = tke_to_kd(i,k) * tke_lowmode_lay
1302  if (kd_max >= 0.0) kd_add = min(kd_add, kd_max)
1303  if (k>1) dd%Kd_lowmode(i,j,k) = dd%Kd_lowmode(i,j,k) + 0.5*kd_add
1304  if (k<nz) dd%Kd_lowmode(i,j,k+1) = dd%Kd_lowmode(i,j,k+1) + 0.5*kd_add
1305  endif
1306  if (associated(dd%Kd_lowmode_work)) &
1307  dd%Kd_lowmode_work(i,j,k) = gv%Rho0 * tke_lowmode_lay
1308  if (associated(dd%Fl_lowmode)) dd%Fl_lowmode(i,j,k) = tke_lowmode_rem(i)
1309 
1310  enddo ; enddo
1311  endif ! Simmons
1312 
1313  ! Polzin:
1314  if ( use_polzin ) then
1315  do k=nz-1,2,-1 ; do i=is,ie
1316  if (max_tke(i,k) <= 0.0) cycle
1317  z_from_bot(i) = z_from_bot(i) + gv%H_to_Z*h(i,j,k)
1318  if (cs%answers_2018) then
1319  if (n2_meanz(i) > 1.0e-14*us%T_to_s**2 ) then
1320  z_from_bot_wkb(i) = z_from_bot_wkb(i) &
1321  + gv%H_to_Z * h(i,j,k) * n2_lay(i,k) / n2_meanz(i)
1322  else ; z_from_bot_wkb(i) = 0 ; endif
1323  else
1324  if (gv%H_to_Z*h(i,j,k) * n2_lay(i,k) < (1.0e14 * htot_wkb(i)) * n2_meanz(i)) then
1325  z_from_bot_wkb(i) = z_from_bot_wkb(i) + &
1326  gv%H_to_Z*h(i,j,k) * n2_lay(i,k) / n2_meanz(i)
1327  endif
1328  endif
1329 
1330  ! Fraction of bottom flux predicted to reach top of this layer
1331  tke_frac_top(i) = ( inv_int(i) * z0_polzin_scaled(i) ) / &
1332  ( z0_polzin_scaled(i) + z_from_bot_wkb(i) )
1333  z0_psl = z0_polzin_scaled(i)*cs%Decay_scale_factor_lee
1334  tke_frac_top_lee(i) = (inv_int_lee(i) * z0_psl) / (z0_psl + z_from_bot_wkb(i))
1335  tke_frac_top_lowmode(i) = ( inv_int_low(i) * z0_polzin_scaled(i) ) / &
1336  ( z0_polzin_scaled(i) + z_from_bot_wkb(i) )
1337 
1338  ! Actual influx at bottom of layer minus predicted outflux at top of layer to give
1339  ! predicted power expended
1340  tke_itide_lay = tke_itidal_rem(i) - tke_itidal_bot(i) *tke_frac_top(i)
1341  tke_niku_lay = tke_niku_rem(i) - tke_niku_bot(i) * tke_frac_top_lee(i)
1342  tke_lowmode_lay = tke_lowmode_rem(i) - tke_lowmode_bot(i)*tke_frac_top_lowmode(i)
1343 
1344  ! Actual power expended may be less than predicted if stratification is weak; adjust
1345  if (tke_itide_lay + tke_niku_lay + tke_lowmode_lay > (max_tke(i,k))) then
1346  frac_used = (max_tke(i,k)) / (tke_itide_lay + tke_niku_lay + tke_lowmode_lay)
1347  tke_itide_lay = frac_used * tke_itide_lay
1348  tke_niku_lay = frac_used * tke_niku_lay
1349  tke_lowmode_lay = frac_used * tke_lowmode_lay
1350  endif
1351 
1352  ! Calculate vertical flux available to bottom of layer above
1353  tke_itidal_rem(i) = tke_itidal_rem(i) - tke_itide_lay
1354  tke_niku_rem(i) = tke_niku_rem(i) - tke_niku_lay
1355  tke_lowmode_rem(i) = tke_lowmode_rem(i) - tke_lowmode_lay
1356 
1357  ! Convert power to diffusivity
1358  kd_add = tke_to_kd(i,k) * (tke_itide_lay + tke_niku_lay + tke_lowmode_lay)
1359 
1360  if (kd_max >= 0.0) kd_add = min(kd_add, kd_max)
1361  if (present(kd_lay)) then
1362  kd_lay(i,k) = kd_lay(i,k) + kd_add
1363  endif
1364 
1365  if (present(kd_int)) then
1366  kd_int(i,k) = kd_int(i,k) + 0.5 * kd_add
1367  kd_int(i,k+1) = kd_int(i,k+1) + 0.5 * kd_add
1368  endif
1369 
1370  ! diagnostics
1371  if (associated(dd%Kd_itidal)) then
1372  ! If at layers, this is just dd%Kd_itidal(i,j,K) = TKE_to_Kd(i,k) * TKE_itide_lay
1373  ! The following sets the interface diagnostics.
1374  kd_add = tke_to_kd(i,k) * tke_itide_lay
1375  if (kd_max >= 0.0) kd_add = min(kd_add, kd_max)
1376  if (k>1) dd%Kd_itidal(i,j,k) = dd%Kd_itidal(i,j,k) + 0.5*kd_add
1377  if (k<nz) dd%Kd_itidal(i,j,k+1) = dd%Kd_itidal(i,j,k+1) + 0.5*kd_add
1378  endif
1379  if (associated(dd%Kd_Itidal_work)) &
1380  dd%Kd_itidal_work(i,j,k) = gv%Rho0 * tke_itide_lay
1381  if (associated(dd%Fl_itidal)) dd%Fl_itidal(i,j,k) = tke_itidal_rem(i)
1382 
1383  if (associated(dd%Kd_Niku)) then
1384  ! If at layers, this is just dd%Kd_Niku(i,j,K) = TKE_to_Kd(i,k) * TKE_Niku_lay
1385  ! The following sets the interface diagnostics.
1386  kd_add = tke_to_kd(i,k) * tke_niku_lay
1387  if (kd_max >= 0.0) kd_add = min(kd_add, kd_max)
1388  if (k>1) dd%Kd_Niku(i,j,k) = dd%Kd_Niku(i,j,k) + 0.5*kd_add
1389  if (k<nz) dd%Kd_Niku(i,j,k+1) = dd%Kd_Niku(i,j,k+1) + 0.5*kd_add
1390  endif
1391  ! if (associated(dd%Kd_Niku)) dd%Kd_Niku(i,j,K) = TKE_to_Kd(i,k) * TKE_Niku_lay
1392  if (associated(dd%Kd_Niku_work)) dd%Kd_Niku_work(i,j,k) = gv%Rho0 * tke_niku_lay
1393 
1394  if (associated(dd%Kd_lowmode)) then
1395  ! If at layers, dd%Kd_lowmode is just TKE_to_Kd(i,k) * TKE_lowmode_lay
1396  ! The following sets the interface diagnostics.
1397  kd_add = tke_to_kd(i,k) * tke_lowmode_lay
1398  if (kd_max >= 0.0) kd_add = min(kd_add, kd_max)
1399  if (k>1) dd%Kd_lowmode(i,j,k) = dd%Kd_lowmode(i,j,k) + 0.5*kd_add
1400  if (k<nz) dd%Kd_lowmode(i,j,k+1) = dd%Kd_lowmode(i,j,k+1) + 0.5*kd_add
1401  endif
1402  if (associated(dd%Kd_lowmode_work)) &
1403  dd%Kd_lowmode_work(i,j,k) = gv%Rho0 * tke_lowmode_lay
1404  if (associated(dd%Fl_lowmode)) dd%Fl_lowmode(i,j,k) = tke_lowmode_rem(i)
1405 
1406  enddo ; enddo
1407  endif ! Polzin
1408 
1409 end subroutine add_int_tide_diffusivity
1410 
1411 !> Sets up diagnostics arrays for tidal mixing.
1412 subroutine setup_tidal_diagnostics(G,CS)
1413  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1414  type(tidal_mixing_cs), pointer :: CS !< The control structure for this module
1415 
1416  ! local
1417  integer :: isd, ied, jsd, jed, nz
1418  type(tidal_mixing_diags), pointer :: dd => null()
1419 
1420  isd = g%isd; ied = g%ied; jsd = g%jsd; jed = g%jed; nz = g%ke
1421  dd => cs%dd
1422 
1423  if ((cs%id_Kd_itidal > 0) .or. (cs%id_Kd_Itidal_work > 0)) then
1424  allocate(dd%Kd_itidal(isd:ied,jsd:jed,nz+1)) ; dd%Kd_itidal(:,:,:) = 0.0
1425  endif
1426  if ((cs%id_Kd_lowmode > 0) .or. (cs%id_Kd_lowmode_work > 0)) then
1427  allocate(dd%Kd_lowmode(isd:ied,jsd:jed,nz+1)) ; dd%Kd_lowmode(:,:,:) = 0.0
1428  endif
1429  if ( (cs%id_Fl_itidal > 0) ) then
1430  allocate(dd%Fl_itidal(isd:ied,jsd:jed,nz+1)) ; dd%Fl_itidal(:,:,:) = 0.0
1431  endif
1432  if ( (cs%id_Fl_lowmode > 0) ) then
1433  allocate(dd%Fl_lowmode(isd:ied,jsd:jed,nz+1)) ; dd%Fl_lowmode(:,:,:) = 0.0
1434  endif
1435  if ( (cs%id_Polzin_decay_scale > 0) ) then
1436  allocate(dd%Polzin_decay_scale(isd:ied,jsd:jed))
1437  dd%Polzin_decay_scale(:,:) = 0.0
1438  endif
1439  if ( (cs%id_N2_bot > 0) ) then
1440  allocate(dd%N2_bot(isd:ied,jsd:jed)) ; dd%N2_bot(:,:) = 0.0
1441  endif
1442  if ( (cs%id_N2_meanz > 0) ) then
1443  allocate(dd%N2_meanz(isd:ied,jsd:jed)) ; dd%N2_meanz(:,:) = 0.0
1444  endif
1445  if ( (cs%id_Polzin_decay_scale_scaled > 0) ) then
1446  allocate(dd%Polzin_decay_scale_scaled(isd:ied,jsd:jed))
1447  dd%Polzin_decay_scale_scaled(:,:) = 0.0
1448  endif
1449  if ((cs%id_Kd_Niku > 0) .or. (cs%id_Kd_Niku_work > 0)) then
1450  allocate(dd%Kd_Niku(isd:ied,jsd:jed,nz+1)) ; dd%Kd_Niku(:,:,:) = 0.0
1451  endif
1452  if (cs%id_Kd_Niku_work > 0) then
1453  allocate(dd%Kd_Niku_work(isd:ied,jsd:jed,nz)) ; dd%Kd_Niku_work(:,:,:) = 0.0
1454  endif
1455  if (cs%id_Kd_Itidal_work > 0) then
1456  allocate(dd%Kd_Itidal_work(isd:ied,jsd:jed,nz))
1457  dd%Kd_Itidal_work(:,:,:) = 0.0
1458  endif
1459  if (cs%id_Kd_Lowmode_Work > 0) then
1460  allocate(dd%Kd_Lowmode_Work(isd:ied,jsd:jed,nz))
1461  dd%Kd_Lowmode_Work(:,:,:) = 0.0
1462  endif
1463  if (cs%id_TKE_itidal > 0) then
1464  allocate(dd%TKE_Itidal_used(isd:ied,jsd:jed)) ; dd%TKE_Itidal_used(:,:) = 0.
1465  endif
1466  ! additional diags for CVMix
1467  if (cs%id_N2_int > 0) then
1468  allocate(dd%N2_int(isd:ied,jsd:jed,nz+1)) ; dd%N2_int(:,:,:) = 0.0
1469  endif
1470  if (cs%id_Simmons_coeff > 0) then
1471  if (cs%CVMix_tidal_scheme .ne. simmons) then
1472  call mom_error(fatal, "setup_tidal_diagnostics: Simmons_coeff diagnostics is available "//&
1473  "only when CVMix_tidal_scheme is Simmons")
1474  endif
1475  allocate(dd%Simmons_coeff_2d(isd:ied,jsd:jed)) ; dd%Simmons_coeff_2d(:,:) = 0.0
1476  endif
1477  if (cs%id_vert_dep > 0) then
1478  allocate(dd%vert_dep_3d(isd:ied,jsd:jed,nz+1)) ; dd%vert_dep_3d(:,:,:) = 0.0
1479  endif
1480  if (cs%id_Schmittner_coeff > 0) then
1481  if (cs%CVMix_tidal_scheme .ne. schmittner) then
1482  call mom_error(fatal, "setup_tidal_diagnostics: Schmittner_coeff diagnostics is available "//&
1483  "only when CVMix_tidal_scheme is Schmittner.")
1484  endif
1485  allocate(dd%Schmittner_coeff_3d(isd:ied,jsd:jed,nz)) ; dd%Schmittner_coeff_3d(:,:,:) = 0.0
1486  endif
1487  if (cs%id_tidal_qe_md > 0) then
1488  if (cs%CVMix_tidal_scheme .ne. schmittner) then
1489  call mom_error(fatal, "setup_tidal_diagnostics: tidal_qe_md diagnostics is available "//&
1490  "only when CVMix_tidal_scheme is Schmittner.")
1491  endif
1492  allocate(dd%tidal_qe_md(isd:ied,jsd:jed,nz)) ; dd%tidal_qe_md(:,:,:) = 0.0
1493  endif
1494 end subroutine setup_tidal_diagnostics
1495 
1496 !> This subroutine offers up diagnostics of the tidal mixing.
1497 subroutine post_tidal_diagnostics(G, GV, h ,CS)
1498  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1499  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
1500  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1501  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
1502  type(tidal_mixing_cs), pointer :: CS !< The control structure for this module
1503 
1504  ! local
1505  type(tidal_mixing_diags), pointer :: dd => null()
1506 
1507  dd => cs%dd
1508 
1509  if (cs%Int_tide_dissipation .or. cs%Lee_wave_dissipation .or. cs%Lowmode_itidal_dissipation) then
1510  if (cs%id_TKE_itidal > 0) call post_data(cs%id_TKE_itidal, dd%TKE_itidal_used, cs%diag)
1511  if (cs%id_TKE_leewave > 0) call post_data(cs%id_TKE_leewave, cs%TKE_Niku, cs%diag)
1512  if (cs%id_Nb > 0) call post_data(cs%id_Nb, cs%Nb, cs%diag)
1513  if (cs%id_N2_bot > 0) call post_data(cs%id_N2_bot, dd%N2_bot, cs%diag)
1514  if (cs%id_N2_meanz > 0) call post_data(cs%id_N2_meanz,dd%N2_meanz,cs%diag)
1515 
1516  if (cs%id_Fl_itidal > 0) call post_data(cs%id_Fl_itidal, dd%Fl_itidal, cs%diag)
1517  if (cs%id_Kd_itidal > 0) call post_data(cs%id_Kd_itidal, dd%Kd_itidal, cs%diag)
1518  if (cs%id_Kd_Niku > 0) call post_data(cs%id_Kd_Niku, dd%Kd_Niku, cs%diag)
1519  if (cs%id_Kd_lowmode> 0) call post_data(cs%id_Kd_lowmode, dd%Kd_lowmode, cs%diag)
1520  if (cs%id_Fl_lowmode> 0) call post_data(cs%id_Fl_lowmode, dd%Fl_lowmode, cs%diag)
1521 
1522  if (cs%id_N2_int> 0) call post_data(cs%id_N2_int, dd%N2_int, cs%diag)
1523  if (cs%id_vert_dep> 0) call post_data(cs%id_vert_dep, dd%vert_dep_3d, cs%diag)
1524  if (cs%id_Simmons_coeff> 0) call post_data(cs%id_Simmons_coeff, dd%Simmons_coeff_2d, cs%diag)
1525  if (cs%id_Schmittner_coeff> 0) call post_data(cs%id_Schmittner_coeff, dd%Schmittner_coeff_3d, cs%diag)
1526  if (cs%id_tidal_qe_md> 0) call post_data(cs%id_tidal_qe_md, dd%tidal_qe_md, cs%diag)
1527 
1528  if (cs%id_Kd_Itidal_Work > 0) &
1529  call post_data(cs%id_Kd_Itidal_Work, dd%Kd_Itidal_Work, cs%diag)
1530  if (cs%id_Kd_Niku_Work > 0) call post_data(cs%id_Kd_Niku_Work, dd%Kd_Niku_Work, cs%diag)
1531  if (cs%id_Kd_Lowmode_Work > 0) &
1532  call post_data(cs%id_Kd_Lowmode_Work, dd%Kd_Lowmode_Work, cs%diag)
1533 
1534  if (cs%id_Polzin_decay_scale > 0 ) &
1535  call post_data(cs%id_Polzin_decay_scale, dd%Polzin_decay_scale, cs%diag)
1536  if (cs%id_Polzin_decay_scale_scaled > 0 ) &
1537  call post_data(cs%id_Polzin_decay_scale_scaled, dd%Polzin_decay_scale_scaled, cs%diag)
1538  endif
1539 
1540  if (associated(dd%Kd_itidal)) deallocate(dd%Kd_itidal)
1541  if (associated(dd%Kd_lowmode)) deallocate(dd%Kd_lowmode)
1542  if (associated(dd%Fl_itidal)) deallocate(dd%Fl_itidal)
1543  if (associated(dd%Fl_lowmode)) deallocate(dd%Fl_lowmode)
1544  if (associated(dd%Polzin_decay_scale)) deallocate(dd%Polzin_decay_scale)
1545  if (associated(dd%Polzin_decay_scale_scaled)) deallocate(dd%Polzin_decay_scale_scaled)
1546  if (associated(dd%N2_bot)) deallocate(dd%N2_bot)
1547  if (associated(dd%N2_meanz)) deallocate(dd%N2_meanz)
1548  if (associated(dd%Kd_Niku)) deallocate(dd%Kd_Niku)
1549  if (associated(dd%Kd_Niku_work)) deallocate(dd%Kd_Niku_work)
1550  if (associated(dd%Kd_Itidal_Work)) deallocate(dd%Kd_Itidal_Work)
1551  if (associated(dd%Kd_Lowmode_Work)) deallocate(dd%Kd_Lowmode_Work)
1552  if (associated(dd%TKE_itidal_used)) deallocate(dd%TKE_itidal_used)
1553  if (associated(dd%N2_int)) deallocate(dd%N2_int)
1554  if (associated(dd%vert_dep_3d)) deallocate(dd%vert_dep_3d)
1555  if (associated(dd%Simmons_coeff_2d)) deallocate(dd%Simmons_coeff_2d)
1556  if (associated(dd%Schmittner_coeff_3d)) deallocate(dd%Schmittner_coeff_3d)
1557  if (associated(dd%tidal_qe_md)) deallocate(dd%tidal_qe_md)
1558 
1559 end subroutine post_tidal_diagnostics
1560 
1561 !> This subroutine returns a zonal slice of the topographic roughness amplitudes
1562 subroutine tidal_mixing_h_amp(h_amp, G, j, CS)
1563  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1564  real, dimension(SZI_(G)), intent(out) :: h_amp !< The topographic roughness amplitude [Z ~> m]
1565  integer, intent(in) :: j !< j-index of the row to work on
1566  type(tidal_mixing_cs), pointer :: CS !< The control structure for this module
1567 
1568  integer :: i
1569 
1570  h_amp(:) = 0.0
1571  if ( cs%Int_tide_dissipation .and. .not. cs%use_CVMix_tidal ) then
1572  do i=g%isc,g%iec
1573  h_amp(i) = sqrt(cs%h2(i,j))
1574  enddo
1575  endif
1576 
1577 end subroutine tidal_mixing_h_amp
1578 
1579 ! TODO: move this subroutine to MOM_internal_tide_input module (?)
1580 !> This subroutine read tidal energy inputs from a file.
1581 subroutine read_tidal_energy(G, US, tidal_energy_type, tidal_energy_file, CS)
1582  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1583  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1584  character(len=20), intent(in) :: tidal_energy_type !< The type of tidal energy inputs to read
1585  character(len=200), intent(in) :: tidal_energy_file !< The file from which to read tidalinputs
1586  type(tidal_mixing_cs), pointer :: CS !< The control structure for this module
1587  ! local
1588  integer :: i, j, isd, ied, jsd, jed, nz
1589  real, allocatable, dimension(:,:) :: tidal_energy_flux_2d ! input tidal energy flux at T-grid points [W m-2]
1590 
1591  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = g%ke
1592 
1593  select case (uppercase(tidal_energy_type(1:4)))
1594  case ('JAYN') ! Jayne 2009
1595  if (.not. allocated(cs%tidal_qe_2d)) allocate(cs%tidal_qe_2d(isd:ied,jsd:jed))
1596  allocate(tidal_energy_flux_2d(isd:ied,jsd:jed))
1597  call mom_read_data(tidal_energy_file,'wave_dissipation',tidal_energy_flux_2d, g%domain)
1598  do j=g%jsc,g%jec ; do i=g%isc,g%iec
1599  cs%tidal_qe_2d(i,j) = cs%Gamma_itides * tidal_energy_flux_2d(i,j)
1600  enddo ; enddo
1601  deallocate(tidal_energy_flux_2d)
1602  case ('ER03') ! Egbert & Ray 2003
1603  call read_tidal_constituents(g, us, tidal_energy_file, cs)
1604  case default
1605  call mom_error(fatal, "read_tidal_energy: Unknown tidal energy file type.")
1606  end select
1607 
1608 end subroutine read_tidal_energy
1609 
1610 !> This subroutine reads tidal input energy from a file by constituent.
1611 subroutine read_tidal_constituents(G, US, tidal_energy_file, CS)
1612  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1613  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1614  character(len=200), intent(in) :: tidal_energy_file !< The file from which to read tidal energy inputs
1615  type(tidal_mixing_cs), pointer :: CS !< The control structure for this module
1616 
1617  ! local variables
1618  real, parameter :: C1_3 = 1.0/3.0
1619  real, dimension(SZI_(G),SZJ_(G)) :: &
1620  tidal_qk1, & ! qk1 coefficient used in Schmittner & Egbert
1621  tidal_qo1 ! qo1 coefficient used in Schmittner & Egbert
1622  real, allocatable, dimension(:) :: &
1623  z_t, & ! depth from surface to midpoint of input layer [Z]
1624  z_w ! depth from surface to top of input layer [Z]
1625  real, allocatable, dimension(:,:,:) :: &
1626  tc_m2, & ! input lunar semidiurnal tidal energy flux [W/m^2]
1627  tc_s2, & ! input solar semidiurnal tidal energy flux [W/m^2]
1628  tc_k1, & ! input lunar diurnal tidal energy flux [W/m^2]
1629  tc_o1 ! input lunar diurnal tidal energy flux [W/m^2]
1630  integer, dimension(4) :: nz_in
1631  integer :: k, is, ie, js, je, isd, ied, jsd, jed, i, j
1632 
1633  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1634  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1635 
1636  ! get number of input levels:
1637  call field_size(tidal_energy_file, 'z_t', nz_in)
1638 
1639  ! allocate local variables
1640  allocate(z_t(nz_in(1)), z_w(nz_in(1)) )
1641  allocate(tc_m2(isd:ied,jsd:jed,nz_in(1)), &
1642  tc_s2(isd:ied,jsd:jed,nz_in(1)), &
1643  tc_k1(isd:ied,jsd:jed,nz_in(1)), &
1644  tc_o1(isd:ied,jsd:jed,nz_in(1)) )
1645 
1646  ! allocate CS variables associated with 3d tidal energy dissipation
1647  if (.not. allocated(cs%tidal_qe_3d_in)) allocate(cs%tidal_qe_3d_in(isd:ied,jsd:jed,nz_in(1)))
1648  if (.not. allocated(cs%h_src)) allocate(cs%h_src(nz_in(1)))
1649 
1650  ! read in tidal constituents
1651  call mom_read_data(tidal_energy_file, 'M2', tc_m2, g%domain)
1652  call mom_read_data(tidal_energy_file, 'S2', tc_s2, g%domain)
1653  call mom_read_data(tidal_energy_file, 'K1', tc_k1, g%domain)
1654  call mom_read_data(tidal_energy_file, 'O1', tc_o1, g%domain)
1655  ! Note the hard-coded assumption that z_t and z_w in the file are in centimeters.
1656  call mom_read_data(tidal_energy_file, 'z_t', z_t, scale=100.0*us%m_to_Z)
1657  call mom_read_data(tidal_energy_file, 'z_w', z_w, scale=100.0*us%m_to_Z)
1658 
1659  do j=js,je ; do i=is,ie
1660  if (abs(g%geoLatT(i,j)) < 30.0) then
1661  tidal_qk1(i,j) = c1_3
1662  tidal_qo1(i,j) = c1_3
1663  else
1664  tidal_qk1(i,j) = 1.0
1665  tidal_qo1(i,j) = 1.0
1666  endif
1667  enddo ; enddo
1668 
1669  cs%tidal_qe_3d_in(:,:,:) = 0.0
1670  do k=1,nz_in(1)
1671  ! Store the input cell thickness in m for use with CVmix.
1672  cs%h_src(k) = us%Z_to_m*(z_t(k)-z_w(k))*2.0
1673  ! form tidal_qe_3d_in from weighted tidal constituents
1674  do j=js,je ; do i=is,ie
1675  if ((z_t(k) <= g%bathyT(i,j)) .and. (z_w(k) > cs%tidal_diss_lim_tc)) &
1676  cs%tidal_qe_3d_in(i,j,k) = c1_3*tc_m2(i,j,k) + c1_3*tc_s2(i,j,k) + &
1677  tidal_qk1(i,j)*tc_k1(i,j,k) + tidal_qo1(i,j)*tc_o1(i,j,k)
1678  enddo ; enddo
1679  enddo
1680 
1681  !open(unit=1905,file="out_1905.txt",access="APPEND")
1682  !do j=G%jsd,G%jed
1683  ! do i=isd,ied
1684  ! if ( i+G%idg_offset .eq. 90 .and. j+G%jdg_offset .eq. 126) then
1685  ! write(1905,*) "-------------------------------------------"
1686  ! do k=50,nz_in(1)
1687  ! write(1905,*) i,j,k
1688  ! write(1905,*) CS%tidal_qe_3d_in(i,j,k), tc_m2(i,j,k)
1689  ! write(1905,*) z_t(k), G%bathyT(i,j), z_w(k),CS%tidal_diss_lim_tc
1690  ! end do
1691  ! endif
1692  ! enddo
1693  !enddo
1694  !close(1905)
1695 
1696  ! test if qE is positive
1697  if (any(cs%tidal_qe_3d_in<0.0)) then
1698  call mom_error(fatal, "read_tidal_constituents: Negative tidal_qe_3d_in terms.")
1699  endif
1700 
1701  !! collapse 3D q*E to 2D q*E
1702  !CS%tidal_qe_2d(:,:) = 0.0
1703  !do k=1,nz_in(1) ; do j=js,je ; do i=is,ie
1704  ! if (z_t(k) <= G%bathyT(i,j)) &
1705  ! CS%tidal_qe_2d(i,j) = CS%tidal_qe_2d(i,j) + CS%tidal_qe_3d_in(i,j,k)
1706  !enddo ; enddo ; enddo
1707 
1708  ! initialize input remapping:
1709  call initialize_remapping(cs%remap_cs, remapping_scheme="PLM", &
1710  boundary_extrapolation=.false., check_remapping=cs%debug, &
1711  answers_2018=cs%remap_answers_2018)
1712 
1713  deallocate(tc_m2)
1714  deallocate(tc_s2)
1715  deallocate(tc_k1)
1716  deallocate(tc_o1)
1717  deallocate(z_t)
1718  deallocate(z_w)
1719 
1720 end subroutine read_tidal_constituents
1721 
1722 !> Clear pointers and deallocate memory
1723 subroutine tidal_mixing_end(CS)
1724  type(tidal_mixing_cs), pointer :: CS !< This module's control structure, which
1725  !! will be deallocated in this routine.
1726 
1727  if (.not.associated(cs)) return
1728 
1729  !TODO deallocate all the dynamically allocated members here ...
1730  if (allocated(cs%tidal_qe_2d)) deallocate(cs%tidal_qe_2d)
1731  if (allocated(cs%tidal_qe_3d_in)) deallocate(cs%tidal_qe_3d_in)
1732  if (allocated(cs%h_src)) deallocate(cs%h_src)
1733  deallocate(cs%dd)
1734  deallocate(cs)
1735 
1736 end subroutine tidal_mixing_end
1737 
1738 end module mom_tidal_mixing
Control structure with parameters for the tidal mixing module.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
Containers for tidal mixing diagnostics.
Provides column-wise vertical remapping functions.
This module contains I/O framework code.
Definition: MOM_io.F90:2
The MOM6 facility to parse input files for runtime parameters.
An overloaded interface to log the values of various types of parameters.
Container for remapping parameters.
Make a diagnostic available for averaging or output.
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Routines for error handling and I/O management.
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Interface to vertical tidal mixing schemes including CVMix tidal mixing.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Handy functions for manipulating strings.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Read a data field from a file.
Definition: MOM_io.F90:74
An overloaded interface to read and log the values of various types of parameters.
Provides checksumming functions for debugging.
A structure for creating arrays of pointers to 3D arrays.