Neko  0.9.99
A portable framework for high-order spectral element flow simulations
Statistics guide

Statistics in the context of Neko, is the common name for fields that are averaged in time and possible also space.

The statistics module in Neko computes the temporal average of a wide range of fields.

In this page we use the following convention for a field

  • \( u \), the instantaneous field.
  • \( \langle u \rangle_t \), the temporal average of \( u \).
  • \( u = \langle u \rangle + u' \), the Reynolds decomposition of \( u \), where \( u' \) is the fluctuation of \( u \) around the mean field.

The temporal average of a field \(u\) is the approximation of the integral

$$ \langle u \rangle_t = \int_{T_0}^{T_N} u dt $$

In Neko, this is computed as

$$ \langle u \rangle_t = \sum_{i=0}^N u_i \Delta t_i $$ where \( u_0 \) is the fields value at \( T_0 \) and \( N \) is the number of time steps needed to reach \( T_N \), \( T_N = T_0 + \sum_{i=0}^N \Delta t_i \).

In the statistics in Neko, various averages of the the different velocity components, derivatives and pressure are computed. In total, 44 "raw statistics" are computed that are required to compute the Reynolds stress budgets, mean fields, and the different terms in the turbulent kinetic energy equation.

Using statistics

Statistics are enable in the the case file as a simcomp with the added argument avg_direction, set_of_stats, and start_time:

Name Description Admissible values Default value
start_time Time at which to start gathering statistics. Positive real 0
avg_direction Directions to compute spatial average. x,y,z,xy,xz,yz No spatial average
set_of_stats What set of stats to compute. basic, full full
compute_value Interval, in timesteps or simulationtime, depending on compute_control, for sampling the flow fields for statistics. Positive real or int Not set (but recommended with every 50 timesteps or so

In addition to the usual controls for the output, which then outputs the averages computes from the last time the statistics were written to file.

For example, if one wants to compute only the basic statistics and sample the fields every 4 time steps and compute and output batches every 20 time units and have an initial transient of 60 time units the following would work:

"simulation_components":
[
{
"type": "fluid_stats",
"compute_control": "tsteps",
"compute_value": 4,
"output_control": "simulationtime",
"output_value": 20,
"start_time":60.0,
"avg_direction":"xz",
"set_of_stats":"basic"
}1
]

Preferably set the initial transient to a multiple of output_value as otherwise the first output will be slightly shorter than the rest. The code related to fluid statistics are located in fluid_stats and fluid_stats_simcomp.

The argument "avg_direction" is optional and if ignored we output 3d fields. The statistics are saved in a fld file according to the following in 2D and 3D. Observe that in 2D the mean Z-velocity is stored in a last scalar field. All other fields are kept the same. This is due to a limitation of the fld file format.

For 1D statistics a CSV file is outputted. The first column is the time at which the statistics are collected, the second column the spatial coordinate, and the rest of the data is stored in the order below. In this case all statistics are kept in the same order as in 3D.

List of fields in output files

When only the basic set of stats is enabled, only stats 1-11 are computed. When 2D stats are enabled \( \langle w \rangle \) is stored in s7 for basic stats and in s40 for the full set of statistics.

Number Statistic Stored in variable (for fld files)
1 \( \langle p \rangle \) Pressure
2 \( \langle u \rangle \) X-Velocity
3 \( \langle v \rangle \) Y-Velocity
4 \( \langle w \rangle \) Z-Velocity
5 \( \langle pp \rangle \) Temperature
6 \( \langle uu \rangle \) Scalar 1 (s1)
7 \( \langle vv \rangle \) Scalar 2 (s2)
8 \( \langle ww \rangle \) s3
9 \( \langle uv \rangle \) s4
10 \( \langle uw \rangle \) s5
11 \( \langle vw \rangle \) s6
12 \( \langle uuu \rangle \) s7
13 \( \langle vvv \rangle \) s8
14 \( \langle www \rangle \) s9
15 \( \langle uuv \rangle \) s10
16 \( \langle uuw \rangle \) s11
17 \( \langle uvv \rangle \) s12
18 \( \langle uvw \rangle \) s13
19 \( \langle vvw \rangle \) s14
20 \( \langle uww \rangle \) s15
21 \( \langle vww \rangle \) s16
22 \( \langle uuuu \rangle \) s17
23 \( \langle vvvv \rangle \) s18
24 \( \langle wwww \rangle \) s19
25 \( \langle ppp \rangle \) s20
26 \( \langle pppp \rangle \) s21
27 \( \langle pu \rangle \) s22
28 \( \langle pv \rangle \) s23
29 \( \langle pw \rangle \) s24
30 \( \langle p \frac{\partial u} {\partial x} \rangle \) s25
31 \( \langle p \frac{\partial u} {\partial y}\rangle \) s26
32 \( \langle p \frac{\partial u} {\partial z}\rangle \) s27
33 \( \langle p \frac{\partial v} {\partial x}\rangle \) s28
34 \( \langle p \frac{\partial v} {\partial y}\rangle \) s29
35 \( \langle p \frac{\partial v} {\partial z}\rangle \) s30
36 \( \langle p \frac{\partial w} {\partial x}\rangle \) s31
37 \( \langle p \frac{\partial w} {\partial y}\rangle \) s32
38 \( \langle p \frac{\partial w} {\partial z}\rangle \) s33
39 \( \langle e11 \rangle \) s34
40 \( \langle e22 \rangle \) s35
41 \( \langle e33 \rangle \) s36
42 \( \langle e12 \rangle \) s37
43 \( \langle e13 \rangle \) s38
44 \( \langle e23 \rangle \) s39

where \(e11,e22...\) is computed as: $$ \begin{aligned} e11 &= \left(\frac{\partial u}{\partial x}\right)^2 + \left(\frac{\partial u}{\partial y}\right)^2 + \left(\frac{\partial u}{\partial z}\right)^2 \\ e22 &= \left(\frac{\partial v}{\partial x}\right)^2 + \left(\frac{\partial v}{\partial y}\right)^2 + \left(\frac{\partial v}{\partial z}\right)^2 \\ e33 &= \left(\frac{\partial w}{\partial x}\right)^2 + \left(\frac{\partial w}{\partial y}\right)^2 + \left(\frac{\partial w}{\partial z}\right)^2 \\ e12 &= \frac{\partial u}{\partial x} \frac{\partial v}{\partial x} + \frac{\partial u}{\partial y}\frac{\partial v}{\partial y}+ \frac{\partial u}{\partial z}\frac{\partial v}{\partial z} \\ e13 &= \frac{\partial u}{\partial x} \frac{\partial w}{\partial x} + \frac{\partial u}{\partial y}\frac{\partial w}{\partial y}+ \frac{\partial u}{\partial z}\frac{\partial w}{\partial z} \\ e23 &= \frac{\partial v}{\partial x} \frac{\partial w}{\partial x} + \frac{\partial v}{\partial y}\frac{\partial w}{\partial y}+ \frac{\partial v}{\partial z}\frac{\partial w}{\partial z} \\ \end{aligned} $$

Postprocessing

These statistics are only the "raw statistics" in the sense that in general we are not interested in \( \langle uu\rangle \), but rather say \( \langle u'u'\rangle\). For this we need to postprocess the statistics.

There is some rudimentary postprocessing to compute the spatial averages of fld files and also to combine the statistics collected from several runs (compute average in time) and also compute both the mean velocity gradient and the Reynolds stresses available in the contrib scripts. By running the contrib scripts without any arguments one gets a hint on their usage, (e.g. by running ./average_fields_in_time will give the options).

To further postprocess the statistics it is suggested to look into PyNekTools which introduces convenient functions for postprocessing, largely based on the PyMech library with the addition of an expanding set of tools for interpolation, computing derivatives and more advanced functionality. PyNekTools is entirely parallelized in MPI and can also handle large data sets for postprocessing.