-
Notifications
You must be signed in to change notification settings - Fork 0
/
compute_flow_and_simulation_metrics.py
397 lines (329 loc) · 14.1 KB
/
compute_flow_and_simulation_metrics.py
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
#This script is developed by KVSlab at Simula Research Laboratory
#All rights and acknowlegement belong to KVSlab
#Please refer to https://github.com/KVSlab/VaMPy for original version
from os import path
from time import time
from glob import glob
import numpy as np
from dolfin import *
from postprocessing_common import read_command_line, epsilon
def compute_flow_and_simulation_metrics(folder, nu, dt, velocity_degree):
"""
Computes several flow field characteristics
for velocity field stored at 'folder' location
for flow_metrics given viscosity and time step
Args:
velocity_degree (int): Finite element degree of velocity
folder (str): Path to simulation results
nu (float): Viscosity
dt (float): Time step
"""
# File paths
file_path_u = path.join(folder, "u.h5")
if len(glob(folder+"/mesh.h5"))==1:
mesh_path = path.join(folder, "mesh.h5")
# Get mesh information
mesh = Mesh()
with HDF5File(MPI.comm_world, mesh_path, "r") as mesh_file:
mesh_file.read(mesh, "mesh", False)
elif len(glob(folder+"/mesh.xml.gz"))==1:
mesh_path = path.join(folder,"mesh.xml.gz")
mesh=Mesh(mesh_path)
else:
print ("No mesh found in .h5 or .xml.gz format")
print ("Exiting...")
exit(1)
f = HDF5File(MPI.comm_world, file_path_u, "r")
# Get names of data to extract
start = 0
if MPI.rank(MPI.comm_world) == 0:
print("The post processing starts from", start)
dataset_names = get_dataset_names(f, start=start)
# Function space
print ("Creating Function Spaces")
DG = FunctionSpace(mesh, 'DG', 0)
V = VectorFunctionSpace(mesh, "CG", velocity_degree)
Vv = FunctionSpace(mesh, "CG", velocity_degree)
# Mesh info
h = CellDiameter(mesh)
characteristic_edge_length = project(h, DG)
# Functions for storing values
v = TestFunction(DG)
u = Function(V)
u_mean = Function(V)
u_prime = Function(V)
# Plus-values
l_plus_avg = Function(DG)
l_plus = Function(DG)
t_plus_avg = Function(DG)
t_plus = Function(DG)
# Kolmogorov scales
length_scale = Function(DG)
length_scale_avg = Function(DG)
time_scale = Function(DG)
time_scale_avg = Function(DG)
velocity_scale = Function(DG)
velocity_scale_avg = Function(DG)
# Inner grad(u), grad(u)
turbulent_dissipation = Function(DG)
turbulent_dissipation_avg = Function(DG)
strain = Function(DG)
strain_avg = Function(DG)
dissipation = Function(DG)
dissipation_avg = Function(DG)
# Energy
kinetic_energy = Function(Vv)
kinetic_energy_avg = Function(Vv)
turbulent_kinetic_energy = Function(Vv)
turbulent_kinetic_energy_avg = Function(Vv)
# Velocity
u0 = Function(Vv)
u1 = Function(Vv)
u2 = Function(Vv)
u0_prime = Function(Vv)
u1_prime = Function(Vv)
u2_prime = Function(Vv)
# CFL
CFL = Function(DG)
CFL_avg = Function(DG)
# Create XDMF files for saving metrics
fullname = file_path_u.replace("u.h5", "%s.xdmf")
fullname = fullname.replace("Solutions", "flow_metrics")
var_name = ["u_mean", "l_plus", "t_plus", "CFL", "strain", "length_scale", "time_scale", "velocity_scale", "u_mag",
"characteristic_edge_length", "dissipation", "kinetic_energy", "turbulent_kinetic_energy",
"turbulent_dissipation", "u_prime"]
metrics = {}
for vn in var_name:
if MPI.rank(MPI.comm_world) == 0:
print(fullname % vn)
metrics[vn] = XDMFFile(MPI.comm_world, fullname % vn)
metrics[vn].parameters["rewrite_function_mesh"] = False
metrics[vn].parameters["flush_output"] = True
# Get u mean
u_mean_file_path = file_path_u.replace("u.h5", "u_mean.h5")
tmp_file = HDF5File(MPI.comm_world, u_mean_file_path, "r")
tmp_file.read(u, "u_mean/vector_0")
tmp_file.close()
assign(u_mean, u)
counter = 0
for data in dataset_names:
counter += 1
if MPI.rank(MPI.comm_world) == 0:
print(data)
# Time step and velocity
f.read(u, data)
# Compute CFL
t0 = Timer("CFL")
u_mag = project(sqrt(inner(u, u)), DG)
CFL.vector().set_local(u_mag.vector().get_local() / characteristic_edge_length.vector().get_local() * dt)
CFL.vector().apply("insert")
CFL_avg.vector().axpy(1, CFL.vector())
t0.stop()
# Compute rate-of-strain
t0 = Timer("rate of strain")
rate_of_strain(strain, u, v, mesh, h)
strain_avg.vector().axpy(1, strain.vector())
t0.stop()
# Compute l+
t0 = Timer("l plus")
u_star = np.sqrt(strain.vector().get_local() * nu)
l_plus.vector().set_local(u_star * characteristic_edge_length.vector().get_local() / nu)
l_plus.vector().apply("insert")
l_plus_avg.vector().axpy(1, l_plus.vector())
t0.stop()
# Compute t+
t0 = Timer("t plus")
t_plus.vector().set_local(u_star ** 2 * dt / nu)
t_plus.vector().apply("insert")
t_plus_avg.vector().axpy(1, t_plus.vector())
t0.stop()
# Compute Kolmogorov
t0 = Timer("dissipation")
rate_of_dissipation(dissipation, u, v, mesh, h, nu)
dissipation_avg.vector().axpy(1, dissipation.vector())
t0.stop()
# Compute u_prime
t0 = Timer("u prime")
u_prime.vector().set_local(u.vector().get_local() - u_mean.vector().get_local())
u_prime.vector().apply("insert")
t0.stop()
# Compute Turbulent dissipation
t0 = Timer("turbulent dissipation")
rate_of_dissipation(turbulent_dissipation, u_prime, v, mesh, h, nu)
turbulent_dissipation_avg.vector().axpy(1, turbulent_dissipation.vector())
eps = turbulent_dissipation.vector().get_local()
t0.stop()
# Compute length scale
t0 = Timer("length scale")
length_scale.vector().set_local((nu ** 3 / eps) ** (1. / 4))
length_scale.vector().apply("insert")
length_scale_avg.vector().axpy(1, length_scale.vector())
t0.stop()
# Compute time scale
t0 = Timer("time scale")
time_scale.vector().set_local((nu / eps) ** 0.5)
time_scale.vector().apply("insert")
time_scale_avg.vector().axpy(1, time_scale.vector())
t0.stop()
# Compute velocity scale
t0 = Timer("velocity scale")
velocity_scale.vector().set_local((eps * nu) ** (1. / 4))
velocity_scale.vector().apply("insert")
velocity_scale_avg.vector().axpy(1, velocity_scale.vector())
t0.stop()
# Compute both kinetic energy and turbulent kinetic energy
t0 = Timer("kinetic energy")
assign(u0, u.sub(0))
assign(u1, u.sub(1))
if mesh.geometry().dim() == 3:
assign(u2, u.sub(2))
kinetic_energy.vector().set_local(
0.5 * (u0.vector().get_local() ** 2 + u1.vector().get_local() ** 2 + u2.vector().get_local() ** 2))
kinetic_energy.vector().apply("insert")
kinetic_energy_avg.vector().axpy(1, kinetic_energy.vector())
t0.stop()
t0 = Timer("turbulent kinetic energy")
assign(u0_prime, u_prime.sub(0))
assign(u1_prime, u_prime.sub(1))
if mesh.geometry().dim() == 3:
assign(u2_prime, u_prime.sub(2))
turbulent_kinetic_energy.vector().set_local(
0.5 * (u0_prime.vector().get_local() ** 2
+ u1_prime.vector().get_local() ** 2
+ u2_prime.vector().get_local() ** 2))
turbulent_kinetic_energy.vector().apply("insert")
turbulent_kinetic_energy_avg.vector().axpy(1, turbulent_kinetic_energy.vector())
t0.stop()
if counter % 10 == 0:
list_timings(TimingClear.clear, [TimingType.wall])
# Get avg
N = len(dataset_names)
l_plus_avg.vector()[:] = l_plus_avg.vector()[:] / N
t_plus_avg.vector()[:] = t_plus_avg.vector()[:] / N
length_scale_avg.vector()[:] = length_scale_avg.vector()[:] / N
time_scale_avg.vector()[:] = time_scale_avg.vector()[:] / N
velocity_scale_avg.vector()[:] = velocity_scale_avg.vector()[:] / N
CFL_avg.vector()[:] = CFL_avg.vector()[:] / N
dissipation_avg.vector()[:] = dissipation_avg.vector()[:] / N
kinetic_energy_avg.vector()[:] = kinetic_energy_avg.vector()[:] / N
turbulent_kinetic_energy_avg.vector()[:] = turbulent_kinetic_energy_avg.vector()[:] / N
turbulent_dissipation_avg.vector()[:] = turbulent_dissipation_avg.vector()[:] / N
#Write in VTU Format
File(path.join(folder,"CFL.pvd"))<<CFL_avg
File(path.join(folder,"l_plus.pvd"))<<l_plus_avg
File(path.join(folder,"t_plus.pvd"))<<t_plus_avg
File(path.join(folder,"length_scale.pvd"))<<length_scale_avg
File(path.join(folder,"time_scale.pvd"))<<time_scale_avg
File(path.join(folder,"velocity_scale.pvd"))<<velocity_scale_avg
File(path.join(folder,"dissipation.pvd"))<<dissipation_avg
File(path.join(folder,"kinetic_energy.pvd"))<<kinetic_energy_avg
File(path.join(folder,"turbulent_kinetic_energy.pvd"))<<turbulent_kinetic_energy_avg
File(path.join(folder,"turbulent_dissipation.pvd"))<<turbulent_dissipation_avg
File(path.join(folder,"characteristic_edge_length.pvd"))<<characteristic_edge_length
File(path.join(folder,"strain.pvd"))<<strain_avg
File(path.join(folder,"u_mean.pvd"))<<u_mean
# Store average data
metrics["CFL"].write_checkpoint(CFL_avg, "CFL")
metrics["l_plus"].write_checkpoint(l_plus_avg, "l_plus")
metrics["t_plus"].write_checkpoint(t_plus_avg, "t_plus")
metrics["length_scale"].write_checkpoint(length_scale_avg, "length_scale")
metrics["time_scale"].write_checkpoint(time_scale_avg, "time_scale")
metrics["velocity_scale"].write_checkpoint(velocity_scale_avg, "velocity_scale")
metrics["dissipation"].write_checkpoint(dissipation_avg, "dissipation")
metrics["kinetic_energy"].write_checkpoint(kinetic_energy_avg, "kinetic_energy")
metrics["turbulent_kinetic_energy"].write_checkpoint(turbulent_kinetic_energy_avg, "turbulent_kinetic_energy")
metrics["turbulent_dissipation"].write_checkpoint(turbulent_dissipation_avg, "turbulent_dissipation")
metrics["characteristic_edge_length"].write_checkpoint(characteristic_edge_length, "characteristic_edge_length")
metrics["strain"].write_checkpoint(strain_avg, "strain_avg")
metrics["u_mean"].write_checkpoint(u_mean, "u_mean")
# Print info
flow_metrics = [("dx", characteristic_edge_length), ("l+", l_plus_avg), ("t+", t_plus_avg),
("Length scale", length_scale_avg), ("Time scale", time_scale_avg),
("Velocity scale", velocity_scale), ("CFL", CFL_avg), ("Strain", strain_avg),
("Dissipation", dissipation), ("Turbulent dissipation", turbulent_dissipation),
("Turbulent kinetic energy", turbulent_kinetic_energy), ("Kinetic energy", kinetic_energy)]
for metric_name, metric_value in flow_metrics:
sum_ = MPI.sum(MPI.comm_world, np.sum(metric_value.vector().get_local()))
num = MPI.sum(MPI.comm_world, metric_value.vector().get_local().shape[0])
mean = sum_ / num
max_ = MPI.max(MPI.comm_world, metric_value.vector().get_local().max())
min_ = MPI.min(MPI.comm_world, metric_value.vector().get_local().min())
if MPI.rank(MPI.comm_world) == 0:
print(metric_name, "mean:", mean)
print(metric_name, "max:", max_)
print(metric_name, "min:", min_)
def get_dataset_names(data_file, num_files=3000000, step=1, start=1, print_info=True,
vector_filename="/velocity/vector_%d"):
"""
Read velocity fields datasets and extract names of files
Args:
data_file (HDF5File): File object of velocity
num_files (int): Number of files
step (int): Step between each data dump
start (int): Step to start on
print_info (bool): Prints info about data if true
vector_filename (str): Name of velocity files
Returns:
names (list): List of data file names
"""
check = True
# Find start file
t0 = time()
while check:
if data_file.has_dataset(vector_filename % start):
check = False
start -= step
start += step
# Get names
names = []
for i in range(num_files):
step = 1
index = start + i * step
if data_file.has_dataset(vector_filename % index):
names.append(vector_filename % index)
t1 = time()
# Print info
if MPI.rank(MPI.comm_world) == 0 and print_info:
print()
print("=" * 6 + " Timesteps to average over " + "=" * 6)
print("Length on data set names:", len(names))
print("Start index:", start)
print("Wanted num files:", num_files)
print("Step between files:", step)
print("Time used:", t1 - t0, "s")
print()
return names
def rate_of_strain(strain, u, v, mesh, h):
"""
Computes rate of strain
Args:
strain (Function): Function to save rate of strain to
u (Function): Function for velocity field
v (Function): Test function for velocity
mesh: Mesh to compute strain rate on
h (float): Cell diameter of mesh
"""
eps = epsilon(u)
f = sqrt(inner(eps, eps))
x = assemble(inner(f, v) / h * dx(mesh))
strain.vector().set_local(x.get_local())
strain.vector().apply("insert")
def rate_of_dissipation(dissipation, u, v, mesh, h, nu):
"""
Computes rate of dissipation
Args:
dissipation (Function): Function to save rate of dissipation to
u (Function): Function for velocity field
v (Function): Test function for velocity
mesh: Mesh to compute dissipation rate on
h (float): Cell diameter of mesh
nu (float): Viscosity
"""
eps = epsilon(u)
f = 2 * nu * inner(eps, eps)
x = assemble(inner(f, v) / h * dx(mesh))
dissipation.vector().set_local(x.get_local())
dissipation.vector().apply("insert")
if __name__ == '__main__':
folder, nu, _, dt, velocity_degree, _ = read_command_line()
compute_flow_and_simulation_metrics(folder, nu, dt, velocity_degree)