-
Notifications
You must be signed in to change notification settings - Fork 0
/
helmholtz_backup.m
449 lines (377 loc) · 16 KB
/
helmholtz_backup.m
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
function [eps_total, eps_ratio, time_norm,...
K, K_ratio, u_solenoidal, v_solenoidal,...
u_dilatational, v_dilatational,...
u_mean, v_mean, w_mean,...
K_dilatational, K_solenoidal, chi] = helmholtz(file_location,...
file_location_nodes,...
T_ref,...
mu_ref)
% Helmholtz-Hodge decomposition of a 3D velocity field into its
% solenoidal and dilatational parts using fast Fourier transform [1].
%
% The decomposition is performed with the spectral method, which is
% only suitable for relatively smooth fields, i.e., with little power
% on small scales. The code assumes that the grid is uniform with
% dx = dy = dz.
%
% This code is based on Ref. [2] and has been rewritten in MATLAB.
%
% Notes:
% For even NX, NY, and NZ, decomposed fields can be complex,
% with the imaginary part coming from the real part of the kmode
% at Nyquist frequency. In principle, the Nyquist frequency
% kmode should be dropped when doing the first derivatives to
% maintain symmetry. See footnote on page 4 of [2]. However,
% when the field is smooth enough, the imaginary part caused by
% the Nyquist frequency kmode should be negligible.
%
% Args:
% file_location (char): Path to the data .hdf file
% file_location_nodes (char): Path to the grid .hdf file
% T_ref (float): Temperature of reference [K]
% mu_ref (float): Dynamic viscosity of reference [kg/(m-s)] or [Pa-s]
%
% Returns:
% Tuple containing
%
% * eps_total (float): Total dissipation
% * eps_ratio (float): Dissipation ratio (dilatational / solenoidal)
% * time_norm (float): Time normalized with the eddy turnover time
% * u_solenoidal (float): Solenoidal part of the velocity field in the x-axis
% * v_solenoidal (float): Solenoidal part of the velocity field in the y-axis
% * u_dilatational (float): Dilatational part of the velocity field in the x-axis
% * v_dilatational (float): Dilatational part of the velocity field in the y-axis
% * u_mean (float): Mean of the velocity field in the x-axis
% * v_mean (float): Mean of the velocity field in the y-axis
% * w_mean (float): Mean of the velocity field in the z-axis
%
% Examples:
% [eps_total, eps_ratio, time_norm,...
% K, K_ratio, u_solenoidal, v_solenoidal,...
% u_dilatational, v_dilatational,...
% u_mean, v_mean, w_mean] = helmholtz(file_location,...
% file_location_nodes,...
% T_ref,...
% mu_ref)
%
% References:
% [1] Johnson, S. G. (2011). Notes on FFT-based differentiation.
% MIT Applied Mathematics, Tech. Rep.
% Available: http://math.mit.edu/~stevenj/fft-deriv.pdf
% [2] Xun Shi, Helmholtz-Hodge decomposition using fft (Python),
% Available: https://github.com/shixun22/helmholtz
%
%
% @author: Alberto Cuadra Lara
% PhD Candidate - Group Fluid Mechanics
% Universidad Carlos III de Madrid
%
% Last update Apr 05 2023
% Definitions
FPS = 60; % Frames per second
S = 110.4; % Sutherland constant [K]
% Get coordinates
[coordinates_x, coordinates_y, coordinates_z,...
sz_coordinates] = read_3D(file_location_nodes, 'centerCoordinates');
x = reshape(coordinates_x(:, 1, 1), 1, sz_coordinates(1));
y = reshape(coordinates_y(1, :, 1), 1, sz_coordinates(2));
z = reshape(coordinates_z(1, 1, :), 1, sz_coordinates(3));
% Get domain length
L = max(x);
% Get density
rho = read_data(file_location, 'rho');
rho_mean = mean(rho, 'all');
% Get velocity components
[u, v, w, sz] = read_3D(file_location, 'velocity');
% Multiply by sqrt(rho)
u = sqrt(rho) .* u;
v = sqrt(rho) .* v;
w = sqrt(rho) .* w;
% Get mean velocity field
u_mean = mean(u, 'all');
v_mean = mean(v, 'all');
w_mean = mean(w, 'all');
% Remove mean velocity
u = u - mean(u, 'all');
v = v - mean(v, 'all');
w = w - mean(w, 'all');
% Compute ratio density deviations / velocity perturbations (for LIA)
rho_delta = rho - rho_mean;
p = read_data(file_location, 'pressure');
sound_mean = sqrt(1.4 * mean(p, 'all') ./ rho_mean); % ??
vel = sqrt(u.^2 + v.^2 + w.^2);
%chi = mean(rho_delta ./ u, 'all') * (1.1720 / sqrt(101320/1.1720)) * (347.7206 / 1.1720);
% Get turbulent Mach number
Mt = sqrt(3) * sqrt(1/3 * (u.^2 + v.^2 + w.^2)) / sound_mean;
Mt_mean = mean(Mt, 'all');
chi_ast = mean(rho_delta .* vel, 'all') / (rho_mean * sound_mean);
%chi_ast = mean(rho_delta, 'all') / rho_mean * Mt;
chi = chi_ast / Mt_mean^2;
chi = -abs(mean(rho_delta ./ Mt, 'all') / rho_mean);
fprintf('chi (Cuadra): %.2e\n', chi);
% Get temperature
% T = read_data(file_location, 'temperature');
% T_mean = mean(T, 'all');
% T_delta = T - T_mean;
%
% vel = u;
% chi = (mean(rho_delta .* vel, 'all') / (rho_mean * sound_mean)) / (mean(vel .* vel, 'all') / sound_mean^2);
% chi = -(mean(T_delta .* vel, 'all') / (T_mean * sound_mean)) / (mean(vel .* vel, 'all') / sound_mean^2);
% Turbulent Kinetic Energy
K = compute_tke(u, v, w, rho);
% Approximation of the integral length
l = L / 6;
% Velocity rms
vel_rms = sqrt(mean(u.^2 + v.^2 + w.^2, 'all') / 3);
% Eddy turnover time
eddy_turnover = l / vel_rms;
% Time
time = h5readatt(file_location, '/', 'simTime');
% Time normalized with the Eddy turnover time
time_norm = time / eddy_turnover;
% Get N-D fast Fourier transform (fft)
U = fftn(u);
V = fftn(v);
W = fftn(w);
% Get wave numbers
kx = fftfreq(sz(1));
ky = fftfreq(sz(2));
kz = fftfreq(sz(3));
% Get grid wave numbers
[KX, KY, KZ] = ndgrid(kx, ky, kz);
% Get k^2
K2 = KX.^2 + KY.^2 + KZ.^2;
% Avoid infinity value. We do not care about the k = 0 component
K2(1, 1, 1) = 1;
% Compute velocity divergence
div = (U .* KX + V .* KY + W .* KZ);
% Compute the Helmholtz decomposition (dilatational)
H = div ./ K2;
% Get dilatational contributions (curl-free)
u_dilatational = ifftn(H .* KX);
v_dilatational = ifftn(H .* KY);
w_dilatational = ifftn(H .* KZ);
% Get solenoidal contributions (divergence-free)
u_solenoidal = u - u_dilatational;
v_solenoidal = v - v_dilatational;
w_solenoidal = w - w_dilatational;
% Plot slices of the decomposed field on the X-Y plane
% set_figure();
% for slice = 1:sz(1)
% plot_slice(u_solenoidal, v_solenoidal, u_dilatational, v_dilatational, sz, slice);
% pause(1 / FPS);
% end
% Check if the solenoidal part is divergence-free
div_solenoidal = ifftn((fftn(u_solenoidal) .* KX + ...
fftn(v_solenoidal) .* KY + ...
fftn(w_solenoidal) .* KZ) * 1i * 2 * pi);
fprintf('div_solenoidal max: %.2e\n', max(abs(div_solenoidal(:))));
% Check if the dilatational part is curl-free
curl_dilatational = ifftn((fftn(w_dilatational) .* KY - fftn(v_dilatational) .* KZ + ...
fftn(u_dilatational) .* KZ - fftn(w_dilatational) .* KX + ...
fftn(v_dilatational) .* KX - fftn(u_dilatational) .* KY) * 1i * 2 * pi);
fprintf('curl_dilatational max: %.2e\n', max(abs(curl_dilatational(:))));
% Compute velocity derivatives
[dudx, dudy, dudz,...
dvdx, dvdy, dvdz,...
dwdx, dwdy, dwdz] = gradient_periodic_set(u, v, w, x, y, z);
% Compute divergence velocity | div( vf(x, y, z) ) = dudx + dvdy + dwdz
div_vf = dudx + dvdy + dwdz;
% Get temperature field
T = read_data(file_location, 'temperature');
% Compute dynamic viscosity field using the Sutherland's law
mu = compute_mu_sutherland(T, S, T_ref, mu_ref);
% Calculate shear-stress tensor tau_ij assuming Newton's relation
% and Stokes' hypothesis
tau_11 = mu .* (2 * dudx);
tau_22 = mu .* (2 * dvdy);
tau_33 = mu .* (2 * dwdz);
tau_12 = mu .* (dudy + dvdx - 2/3 * div_vf);
tau_13 = mu .* (dudz + dwdx - 2/3 * div_vf);
tau_23 = mu .* (dvdz + dwdy - 2/3 * div_vf);
% Calculate dissipation components eps_alpha = avg(tau_ij * w_i,alpha / sqrt(rho))
eps_solenoidal = mean(dissipation(u_solenoidal, v_solenoidal, w_solenoidal), 'all');
eps_dilatational = mean(dissipation(u_dilatational, v_dilatational, w_dilatational), 'all');
% Compute dissipation ratio eps_dilatational / eps_solenoidal
eps_ratio = eps_dilatational ./ eps_solenoidal;
% Compute dissipation
eps_total = mean(dissipation(u, v, w), 'all');
fprintf('Total dissipation: %.2e\n', eps_total);
% Compute dilatational and solenoidal contributions of the TKE
K_solenoidal = compute_tke(u_solenoidal, v_solenoidal, w_solenoidal, 1);
K_dilatational = compute_tke(u_dilatational, v_dilatational, w_dilatational, 1);
% Ratio dilatational to total TKE
ratio_dilatational = K_dilatational / K;
fprintf('Ratio dilatational: %.2e\n\n', ratio_dilatational);
% Compute TKE ratio
K_ratio = K_dilatational ./ K_solenoidal;
% Vorticity
% omega_x = dwdy - dvdz;
% omega_y = dudz - dwdx;
% omega_z = dvdx - dudy;
% NESTED FUNCTIONS
function value = dissipation(u, v, w)
% Compute dissipation
%
% Args:
% u (float): 3D array with the x-component of the velocity field
% v (float): 3D array with the y-component of the velocity field
% w (float): 3D array with the z-component of the velocity field
%
% Returns:
% value (float): 3D array with the dissipation
%
% Example:
% value = dissipation(u, v, w);
[dudx, dudy, dudz,...
dvdx, dvdy, dvdz,...
dwdx, dwdy, dwdz] = gradient_periodic_set(u, v, w, x, y, z);
value = ...
(tau_11 .* dudx + tau_12 .* dudy + tau_13 .* dudz + ...
tau_12 .* dvdx + tau_22 .* dvdy + tau_23 .* dvdz + ...
tau_13 .* dwdx + tau_23 .* dwdy + tau_33 .* dwdz) ./ sqrt(rho);
end
end
% SUB-PASS FUNCTIONS
function value = read_data(file_location, property)
% Read data from a .hdf file
%
% Args:
% file_location (char): Path to the .hdf file
% property (char): Name of the property to read
%
% Returns:
% value (float): 3D array with the property field
%
% Example:
% value = read_data(file_location, 'rho');
value = h5read(file_location, ['/', property]);
% Reshape
sz = size(value);
N_min = min(sz);
% if sz(1) ~= N_min
% value = value';
% end
end
function [x, y, z, sz] = read_3D(file_location, property)
% Read velocity data
%
% Args:
% file_location (char): Path to the .hdf file
% property (char): Name of the property to read
%
% Returns:
% x (float): 3D array with the x-component of the property field
% y (float): 3D array with the y-component of the property field
% z (float): 3D array with the z-component of the property field
% sz (float): Size of the 3D array
%
% Example:
% [u, v, w, sz] = read_3D(file_location, 'velocity');
value = h5read(file_location, ['/', property]);
% Reshape
sz = size(value);
value = reshape(value, sz);
% Get Dimensions
sz = sz(2:end);
% Get components
x(:, :, :) = value(1, :, :, :);
y(:, :, :) = value(2, :, :, :);
z(:, :, :) = value(3, :, :, :);
end
function mu = compute_mu_sutherland(T, S, T_ref, mu_ref)
% Compute dynamic viscosity using the Sutherland's law
%
% Args:
% T (float): 3D array with the temperature field
% S (float): Sutherland's constant [K]
% T_ref (float): Temperature of reference [K]
% mu_ref (float): Dynamic viscosity of reference [kg/(m-s)] or [Pa-s]
%
% Returns:
% mu (float): 3D array with the dynamic viscosity
%
% Example:
% mu = compute_mu_sutherland(T, S, T_ref, mu_ref);
C1 = mu_ref * ((273.15 + 110.4) / T_ref / (273.15 / T_ref)^(3/2));
mu = C1 * T.^(3/2) ./ (T + S / T_ref);
end
function K = compute_tke(u, v, w, rho)
% Compute Turbulent Kinetic Energy (TKE)
%
% Args:
% u (float): 3D array with the x-component of the velocity field
% v (float): 3D array with the y-component of the velocity field
% w (float): 3D array with the z-component of the velocity field
% rho (float): 3D array with the density field
%
% Returns:
% K (float): 3D array with the TKE
%
% Example:
% K = compute_tke(u, v, w, rho);
K = 0.5 * mean(rho .* (u.^2 + v.^2 + w.^2), 'all');
end
function [dudx, dudy, dudz,...
dvdx, dvdy, dvdz,...
dwdx, dwdy, dwdz] = gradient_periodic_set(u, v, w, x, y, z)
% Compute first derivative of a 3D field over the three axis
% considering a periodic box
%
% Args:
% u (float): 3D array with the x-component of the velocity field
% v (float): 3D array with the y-component of the velocity field
% w (float): 3D array with the z-component of the velocity field
% x (float): 3D array with the x-coordinate of the grid
% y (float): 3D array with the y-coordinate of the grid
% z (float): 3D array with the z-coordinate of the grid
%
% Returns:
% dudx (float): 3D array with the first derivative of u over x
% dudy (float): 3D array with the first derivative of u over y
% dudz (float): 3D array with the first derivative of u over z
% dvdx (float): 3D array with the first derivative of v over x
% dvdy (float): 3D array with the first derivative of v over y
% dvdz (float): 3D array with the first derivative of v over z
% dwdx (float): 3D array with the first derivative of w over x
% dwdy (float): 3D array with the first derivative of w over y
% dwdz (float): 3D array with the first derivative of w over z
%
% Example:
% [dudx, dudy, dudz,...
% dvdx, dvdy, dvdz,...
% dwdx, dwdy, dwdz] = gradient_periodic_set(u, v, w, x, y, z);
% 1. dudx_i
[dudx, dudy, dudz] = gradient_periodic(u, x, y, z);
% 2. dvdx_i
[dvdx, dvdy, dvdz] = gradient_periodic(v, x, y, z);
% 3. dwdx_i
[dwdx, dwdy, dwdz] = gradient_periodic(w, x, y, z);
end
function [grad_x, grad_y, grad_z] = gradient_periodic(f, x, y, z)
% Computes the gradient of a three-dimensional scalar field using
% second-order central finite differences considering a periodic box
%
% Args:
% f (float): Scalar field to compute the gradient of
% x (float): 3D array with the x-coordinate of the grid
% y (float): 3D array with the y-coordinate of the grid
% z (float): 3D array with the z-coordinate of the grid
%
% Returns:
% grad_x (float): 3D array with the gradient of f over x
% grad_y (float): 3D array with the gradient of f over y
% grad_z (float): 3D array with the gradient of f over z
%
% Example:
% [grad_x, grad_y, grad_z] = gradient_periodic(f, x, y, z);
% Get grid spacing in the three directions
hx = x(2) - x(1);
hy = y(2) - y(1);
hz = z(2) - z(1);
% Compute the gradient in each direction using central finite differences
% considering periodic boundary conditions
grad_x = (circshift(f, [-1, 0, 0]) - circshift(f, [1, 0, 0])) ./ (2 * hx);
grad_y = (circshift(f, [ 0, -1, 0]) - circshift(f, [0, 1, 0])) ./ (2 * hy);
grad_z = (circshift(f, [ 0, 0, -1]) - circshift(f, [0, 0, 1])) ./ (2 * hz);
end