Skip to content

Commit e52352e

Browse files
jmcvey3akeeste
authored andcommitted
Allow user to specify universal Kolmogorov constant for TKE dissipation rate function (MHKiT-Software#406)
It was pointed out to me that a universal constant of 0.5 was hardcoded into dolfyn when calculating dissipation rate using the Lumley and Terray method. According to turbulence theory and historical experimental validation, this is the value to use for the streamwise velocity spectra. According to the same, the constant is different for velocity spectra that are orthogonal to the streamwise direction - usually taken to be somewhere between 0.65 and 0.69. This PR allows the user to set this constant via function input. Pope, Stephen B. "Turbulent flows." Measurement Science and Technology 12.11 (2001): 2020-2021. pp 231-232.
1 parent c3473a7 commit e52352e

File tree

5 files changed

+62
-25
lines changed

5 files changed

+62
-25
lines changed
291 Bytes
Binary file not shown.
0 Bytes
Binary file not shown.

mhkit/dolfyn/adp/turbulence.py

Lines changed: 23 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -835,19 +835,26 @@ def check_turbulence_cascade_slope(self, psd, freq_range=[0.2, 0.4]):
835835

836836
return m, b
837837

838-
def dissipation_rate_LT83(self, psd, U_mag, freq_range=[0.2, 0.4], noise=None):
838+
def dissipation_rate_LT83(
839+
self, psd, U_mag, freq_range=[0.2, 0.4], k_constant=0.67, noise=None
840+
):
839841
"""
840842
Calculate the TKE dissipation rate from the velocity spectra.
841843
842844
Parameters
843845
----------
844846
psd : xarray.DataArray (time,f)
845-
The power spectral density from a single depth bin (range)
847+
The power spectral density from the vertical beam and depth
848+
bin (range)
846849
U_mag : xarray.DataArray (time)
847-
The bin-averaged horizontal velocity (a.k.a. speed) from a single depth bin (range)
850+
The bin-averaged horizontal velocity (a.k.a. speed) from a single
851+
depth bin (range)
848852
f_range : iterable(2)
849853
The range over which to integrate/average the spectrum, in units
850854
of the psd frequency vector (Hz or rad/s)
855+
k_constant : float or iterable(3)
856+
Kolmogorov Constant (\\alpha in Notes section below) to use. Default
857+
\\alpha is 0.67.
851858
noise : float or array-like
852859
Instrument noise level in same units as velocity. Typically
853860
found from `adp.turbulence.doppler_noise_level`.
@@ -864,10 +871,9 @@ def dissipation_rate_LT83(self, psd, U_mag, freq_range=[0.2, 0.4], noise=None):
864871
865872
.. math:: S(k) = \\alpha \\epsilon^{2/3} k^{-5/3} + N
866873
867-
where :math:`\\alpha = 0.5` (1.5 for all three velocity
868-
components), `k` is wavenumber, `S(k)` is the turbulent
869-
kinetic energy spectrum, and `N' is the doppler noise level
870-
associated with the TKE spectrum.
874+
where :math:`\\alpha is the Kolmogorov constant (0.67 for vertical direction),
875+
`k` is wavenumber, `S(k)` is the turbulent kinetic energy spectrum, and
876+
`N' is the doppler noise level associated with the TKE spectrum.
871877
872878
With :math:`k \\rightarrow \\omega / U`, then -- to preserve variance --
873879
:math:`S(k) = U S(\\omega)`, and so this becomes:
@@ -882,15 +888,19 @@ def dissipation_rate_LT83(self, psd, U_mag, freq_range=[0.2, 0.4], noise=None):
882888
by a random wave field". JPO, 1983, vol13, pp2000-2007.
883889
"""
884890

891+
if not isinstance(psd, xr.DataArray):
892+
raise TypeError("`psd` must be an instance of `xarray.DataArray`.")
885893
if len(psd.shape) != 2:
886-
raise Exception("PSD should be 2-dimensional (time, frequency)")
894+
raise Exception("`psd` should be 2-dimensional (time, frequency)")
887895
if len(U_mag.shape) != 1:
888896
raise Exception("U_mag should be 1-dimensional (time)")
889897
if not hasattr(freq_range, "__iter__") or len(freq_range) != 2:
890898
raise ValueError("`freq_range` must be an iterable of length 2.")
899+
if np.size(k_constant) != 1:
900+
raise ValueError("`k_constant` should be a single value.")
891901
if noise is not None:
892902
if np.shape(noise)[0] != np.shape(psd)[0]:
893-
raise Exception("Noise should have same first dimension as PSD")
903+
raise Exception("Noise should have same first dimension as `psd`")
894904
else:
895905
noise = np.array(0)
896906

@@ -904,12 +914,15 @@ def dissipation_rate_LT83(self, psd, U_mag, freq_range=[0.2, 0.4], noise=None):
904914
idx = np.where((freq_range[0] < freq) & (freq < freq_range[1]))
905915
idx = idx[0]
906916

917+
# Set the correct magnitude whether the frequency is in Hz or rad/s
907918
if freq.units == "Hz":
908919
U = U_mag / (2 * np.pi)
909920
else:
910921
U = U_mag
911922

912-
a = 0.5
923+
# Use the transverse value derived from the Kolmogorov constant
924+
a = k_constant
925+
# Calculate dissipation
913926
out = (psd[:, idx] * freq[idx] ** (5 / 3) / a).mean(axis=-1) ** (
914927
3 / 2
915928
) / U.values

mhkit/dolfyn/adv/turbulence.py

Lines changed: 36 additions & 12 deletions
Original file line numberDiff line numberDiff line change
@@ -339,9 +339,16 @@ def check_turbulence_cascade_slope(self, psd, freq_range=[6.28, 12.57]):
339339

340340
return m, b
341341

342-
def dissipation_rate_LT83(self, psd, U_mag, freq_range=[6.28, 12.57], noise=None):
342+
def dissipation_rate_LT83(
343+
self,
344+
psd,
345+
U_mag,
346+
freq_range=[6.28, 12.57],
347+
k_constant=[0.5, 0.67, 0.67],
348+
noise=None,
349+
):
343350
"""
344-
Calculate the dissipation rate from the PSD
351+
Calculate the dissipation rate from the power spectral density of velocity.
345352
346353
Parameters
347354
----------
@@ -353,6 +360,12 @@ def dissipation_rate_LT83(self, psd, U_mag, freq_range=[6.28, 12.57], noise=None
353360
The range over which to integrate/average the spectrum, in units
354361
of the psd frequency vector (Hz or rad/s).
355362
Default = [6.28, 12.57] rad/s
363+
k_constant : float or iterable(3)
364+
Kolmogorov Constant (\\alpha in Notes section below) to use. If a
365+
three dimensional PSD is provided, \\alpha defaults to [0.5, 0.67, 0.67];
366+
i.e. 0.5 for the streamwise PSD and 0.67 for the transverse and vertical
367+
PSDs. If the PSD is provided for a single velocity direction, \\alpha is
368+
taken to be 0.5 unless otherwise specified.
356369
noise : float or array-like
357370
Instrument noise level in same units as velocity. Typically
358371
found from `adv.turbulence.calc_doppler_noise`.
@@ -369,10 +382,9 @@ def dissipation_rate_LT83(self, psd, U_mag, freq_range=[6.28, 12.57], noise=None
369382
370383
.. math:: S(k) = \\alpha \\epsilon^{2/3} k^{-5/3} + N
371384
372-
where :math:`\\alpha = 0.5` (1.5 for all three velocity
373-
components), `k` is wavenumber, `S(k)` is the turbulent
374-
kinetic energy spectrum, and `N' is the doppler noise level
375-
associated with the TKE spectrum.
385+
where :math:`\\alpha is the Kolmogorov constant, `k` is wavenumber,
386+
`S(k)` is the turbulent kinetic energy spectrum, and `N' is the
387+
doppler noise level associated with the TKE spectrum.
376388
377389
With :math:`k \\rightarrow \\omega / U`, then -- to preserve variance --
378390
:math:`S(k) = U S(\\omega)`, and so this becomes:
@@ -390,15 +402,19 @@ def dissipation_rate_LT83(self, psd, U_mag, freq_range=[6.28, 12.57], noise=None
390402
if not isinstance(psd, xr.DataArray):
391403
raise TypeError("`psd` must be an instance of `xarray.DataArray`.")
392404
if len(U_mag.shape) != 1:
393-
raise Exception("U_mag should be 1-dimensional (time)")
405+
raise Exception("U_mag should be 1-dimensional (time).")
394406
if len(psd["time"]) != len(U_mag["time"]):
395-
raise Exception("`U_mag` should be from ensembled-averaged dataset")
407+
raise Exception("`U_mag` should be from ensembled-averaged dataset.")
396408
if not hasattr(freq_range, "__iter__") or len(freq_range) != 2:
397409
raise ValueError("`freq_range` must be an iterable of length 2.")
398-
410+
# if the spectra are 1D, then the first dimension should be time (any length)
411+
if (psd.shape[0] != 3) and (np.size(k_constant) != 1):
412+
raise ValueError("`k_constant` should be a single value.")
413+
elif (psd.shape[0] == 3) and (np.size(k_constant) != 3):
414+
raise ValueError("`k_constant` should be an iterable of length 3.")
399415
if noise is not None:
400-
if np.shape(noise)[0] != 3:
401-
raise Exception("Noise should have same first dimension as velocity")
416+
if np.shape(noise)[0] != np.shape(psd)[0]:
417+
raise Exception("Noise should have same first dimension as `psd`.")
402418
else:
403419
noise = np.array([0, 0, 0])[:, None, None]
404420

@@ -412,12 +428,20 @@ def dissipation_rate_LT83(self, psd, U_mag, freq_range=[6.28, 12.57], noise=None
412428
idx = np.where((freq_range[0] < freq) & (freq < freq_range[1]))
413429
idx = idx[0]
414430

431+
# Set the correct magnitude whether the frequency is in Hz or rad/s
415432
if freq.units == "Hz":
416433
U = U_mag / (2 * np.pi)
417434
else:
418435
U = U_mag
419436

420-
a = 0.5
437+
# Set Kolmogorov constant
438+
a = np.array(k_constant)
439+
if psd.shape[0] == 3:
440+
a = a[:, None, None] # stack properly
441+
else:
442+
a = np.squeeze(k_constant)
443+
444+
# Calculate dissipation
421445
out = (psd.isel(freq=idx) * freq.isel(freq=idx) ** (5 / 3) / a).mean(
422446
axis=-1
423447
) ** (3 / 2) / U

mhkit/tests/dolfyn/test_analysis.py

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -143,7 +143,7 @@ def test_adcp_turbulence(self):
143143
dat.velds.rotate2("beam")
144144

145145
tdat["psd"] = bnr.power_spectral_density(
146-
dat["vel"].isel(dir=2, range=len(dat["range"]) // 2), freq_units="Hz"
146+
dat["vel_b5"].isel(range_b5=len(dat["range_b5"]) // 2), freq_units="Hz"
147147
)
148148
tdat["noise"] = bnr.doppler_noise_level(tdat["psd"], pct_fN=0.8)
149149
tdat["stress_vec4"] = bnr.reynolds_stress_4beam(
@@ -179,11 +179,11 @@ def test_adcp_turbulence(self):
179179
) = bnr.dissipation_rate_SF(dat["vel"].isel(dir=2), r_range=[1, 5])
180180

181181
slope_check = bnr.check_turbulence_cascade_slope(
182-
tdat["psd"].mean("time"), freq_range=[0.4, 4]
182+
tdat["psd"].mean("time_b5"), freq_range=[0.4, 4]
183183
)
184184
# Check noise subtraction in psd function
185185
tdat["psd_noise"] = bnr.power_spectral_density(
186-
dat["vel"].isel(dir=2, range=len(dat["range"]) // 2),
186+
dat["vel_b5"].isel(range_b5=len(dat["range_b5"]) // 2),
187187
freq_units="Hz",
188188
noise=0.01,
189189
)

0 commit comments

Comments
 (0)