Setting input
Here we provide an explanation of the various input parameters that can be edited for a typical STREAmS run.
All the parameters are reported in a single file, named singleideal.ini, organized in different sections and using
the classical .ini format. Some ready-to-use input examples are provided in the folder examples
.
File structure of singleideal.ini
[grid]
This section allows the user to specify parameters to define the computational domain and its discretization. When the grid is directly generated in STREAmS, nodes in the x and z directions are always equally spaced.
grid_dim
=> integer (optional, default 1), used to specify Cartesian (grid_dim = 1) or 2D x-y curvilinear (grid_dim = 2) mode.domain_size_x
=> real, domain size in x directiondomain_size_y
=> real, domain size in y directiondomain_size_z
=> real, domain size in z directionnxmax
=> integer, total number of grid nodes in x direction (must be multiple of x_split)nymax
=> integer, total number of grid nodes in y directionnzmax
=> integer, total number of grid nodes in z direction (must be multiple of z_split)ng
=> integer, total number of ghost nodes (must be >= max(ep_order
,visc_order
)/2)metrics_order
=> integer, order of accuracy for metrics computation (suggested value: 2*ng
)grid_type
=> integer, specifiy mesh type1 => Read grid from files
2 => Uniform grid
3 => Channel type grid
4 => Boundary layer type grid
5 => Shock wave boundary layer interaction type grid
6 => Airfoil C-grid type
dyptarget
=> real, active forgrid_type = 3,4,5
, y spacing at the wall in wall-unitsjbgrid
=> integer, active forgrid_type = 3,4,5
, used to control transition from the inner and outer layer resolution (suggested values: 16,24,32)nywr
=> integer, active forgrid_type = 4,5
, number of nodes in the well-resolved region between the wall (y=0) and y =lywr
lywr
=> real, active forgrid_type = 4,5
, size of the well-resolved regionysmoosteps
=> integer, number of iterations to smooth wall-normal grid at boundary layer edge or channel centrelinesave_metrics_2d
=> integer (optional, default 1), used to save c2 metricsgpu_bind
=> integer (optional, default 0). Setgpu_bind = 1
only on those architectures for which mpirun associates each process to a specific GPU (eg APU).grid2d_par
=> integer (optional, default 0). Read mode of c2 grid0 => global mode (each task reads the global grid)
1 => local mode (each task reads serially the local grid)
2 => local mode (each task reads with MPI/IO the local grid)
3 => local mode (each task reads the local grid from an associated file that has been preliminarly decomposed)
Notes:
A non-uniform distribution of nodes is needed in the wall-normal direction (y) for the simulation of wall-bounded flows. In this case, STREAmS adopts the “natural” grid stretching proposed in [12], hereinafter denoted as Pirozzoli mapping function. Depending on the particular flow case under investigation, the following strategies are followed:
channel flow => Pirozzoli mapping function applied to distribute nodes between the two walls (located at \(y = \pm 1\)). To avoid odd behavior in the node distribution, the user shold use a total number of nodes in the wall-normal direction (
nymax
) close to \(N_y = \frac{8}{3} \, Re_{\tau}^{3/4}\)boundary layer and shock boundary layer interaction => Pirozzoli mapping function applied to distribute nodes in the well-resolved region, between y = 0 (lower wall) and
y = lywr
, discretized withnywr
nodes. To avoid odd behavior in the node distribution, the user should use a number of nodes in the well-resolved region (nywr
) close to \({N_y}_{wr} = \frac{4}{3} \, \left (Re_{\tau} \, {L_y}_{wr}\right)^{3/4}\). The remaining part of the computational domain (from y =lywr
up to y =domain_size_y
) is discretized usingnymax - nywr
grid nodes, distributed according to a geometric progression.
When grid_type = 1
, the grid is generated in pre-processing, and so the user must provide the code with 3 ASCII files named x.dat, y.dat, z.dat of size nxmax
,
nymax
, nzmax
respectively, each containing one column with the node coordinates in the corresponding Cartesian direction.
[bc]
This section allows the user to specify flags to impose boundary conditions.
xmin
=> integer, flag for static boundary condition on side xminxmax
=> integer, flag for static boundary condition on side xmaxymin
=> integer, flag for static boundary condition on side yminymax
=> integer, flag for static boundary condition on side ymaxzmin
=> integer, flag for static boundary condition on side zminzmax
=> integer, flag for static boundary condition on side zmaxxmin_nr
=> integer, flag for bc treatment based on characteristic decomposition on side xminxmax_nr
=> integer, flag for bc treatment based on characteristic decomposition on side xmaxymin_nr
=> integer, flag for bc treatment based on characteristic decomposition on side yminymax_nr
=> integer, flag for bc treatment based on characteristic decomposition on side ymaxzmin_nr
=> integer, flag for bc treatment based on characteristic decomposition on side zminzmax_nr
=> integer, flag for bc treatment based on characteristic decomposition on side zmaxx_recyc
=> real, specify recycling location in x directioncorrect_bound_ord
=> integer, flag to reduce the order of accuracy at boundaries and avoid the use of ghost nodes (1 -> active)
Notes:
The following flags can be specified for the static boundary condition:
0 ==> periodicity (available for all sides)
1 ==> free-stream (available for xmin)
2 ==> extrapolation (available for all sides)
4 ==> extrapolation with pressure imposed (available for xmax, ymax)
5 ==> symmetry (available for ymin)
6 ==> viscous wall (available for ymin, ymax)
7 ==> oblique shock entering the domain (available for ymax)
9 ==> compressible Blasius for laminar inflow (available for xmin).
10 ==> recycling-rescaling for turbulent inflow (available for xmin).
Concerning the characteristic treatment, the following options are available for all sides:
0 ==> characteristic treatment is not applied
1 ==> pure non-reflecting
3 ==> relaxation, see [11]
6 ==> reflective wall
[mpi]
This section allows the user to specify parameters for the MPI decomposition.
x_split
=> integer, MPI partitions in x directiony_split
=> integer, MPI partitions in y direction (should be always 1)z_split
=> integer, MPI partitions in z direction
Notes:
The product x_split * y_split * z_split
must be equal to the total number of MPI tasks.
[controls]
This section allows the user to specify parameters to control time integration and restart options.
restart_type
=> integer, 0-> new run from scratch, 1-> continue previous run, 2-> continue previous run and statistics collectioncfl
=> real, Courant number for time step selectioniter_dt_recompute
=> integer, active whencfl > 0
, frequency (number of RK steps, i.e. time iterations) for re-computation of time stepnum_iter
=> integer, number of RK steps
Notes:
When cfl > 0
, the Courant number is constant, and the time step is computed every iter_dt_recompute
cycles.
Alternatively, the time step of the simulation can be directly specified using a negative value of cfl
.
In the latter case, the time step is equal to abs(cfl
).
[numerics]
This section allows the user to specify parameters to control numerical methods for the spatial discretization
ep_order
=> integer, order of accuracy for the discretization of convective termsnkeep
=> integer, KEEP(n) or KEP scheme selectionconservative_viscous
=> integer, flag to activate computation of viscous terms in conservative form (1-> yes)visc_order
=> integer, order of accuracy for the discretization of viscous termsweno_scheme
=> integer, parameter for WENO (order of accuracy is 2 *weno_scheme
- 1)weno_version
=> integer, select WENO version (0 -> Jang & Shu, 1 -> WENO Z)flux_splitting
=> integer, select flux splitting type (0 -> Roe, 1 -> Rusanov)sensor_threshold
=> real, threshold value for shock sensorrk_type
=> integer, specify Runge-Kutta type (for future use, only the third order scheme is currently implemented, set equal to 1)rand_type
=> integer, specify random behavior, < 0 not reproducible, >0 reproducible, 0 random not activeortho
=> integer (optional, default 1), when > 0 allows to activate non-ortho terms in the curvilinear modecount_weno_control
=> integer (optional, default 0), ctive control of weno usage if >0
Notes:
The discretization of the convective terms of the Navier-Stokes equations in smooth flow regions can be selected with nkeep
.
When nkeep<0
the KEP scheme described in [13] is used. When nkeep>=0
the discretization
is based on the family of kinetic energy and entropy preserving schemes KEEP(n), described in [15].
The available orders of accuracy for ep_order
and visc_order
are 2,4,6,8.
The baseline shock capturing flux is based on the classic WENO schemes by [7] (weno_version = 0
).
Admissible values for weno_scheme
are 1,2,3,4, corresponding to first-, third-, fifth- and seventh-order schemes, respectively.
For the fifth-order scheme (weno_scheme = 3
), it is possible to select the optimized WENO-Z scheme proposed by [1] (weno_version=1
).
The conservative discretization of the viscous terms is available only in the Cartesian mode (grid_dim = 1
).
[flow]
This section allows the user to specify parameters to control the flow configuration under investigation.
flow_init
=> integer, specify flow initialization0 => Channel flow (Cartesian mode)
1 => Boundary layer (also used for SBLI, Cartesian and curvilinear mode)
4 => Curvilinear channel flow (curvilinear mode)
5 => Airfoil flow case (curvilinear mode)
Mach
=> real, Mach numberReynolds
=> real, Reynolds numbertheta_wall
=> real, wall-temperature conditionxshock_imp
=> real, oblique shock impingment locationshock_angle
=> real, oblique shock angleT_ref
=> real, reference dimensional temperatureaoa
=> real (optional, default 0), angle of attack (degrees) of the freestream flow
Notes:
For the channel flow configuration, the input Reynolds number is the (estimated) friction Reynolds number \(Re_\tau\).
Simulations are performed at constant flow rate. The bulk temperature \(T_b\) can be specified providing a value for
theta_wall
, defined as \(\Theta = \frac{T_w - T_b}{T_r - T_b}\), where \(T_r\) is the recovery
temperature based on \(T_b\) and the Mach number.
When \(\Theta \ge -1\), a bulk cooling term is added to the total energy equation to keep the bulk temperature constant.
For this type of simulation, Mach
is the Mach number based on the bulk velocity and on the speed of sound evaluated at bulk temperature,
see [10].
Alternatively, when \(\Theta < -1\), the bulk temperature is not constant and is free to evolve.
In this case, the setup is similar to that reported in [9],
with Mach
based on the bulk velocity and on the speed of sound evaluated at wall temperature.
For boundary layer and SBLI configurations, the input Reynolds number is the friction Reynolds number at the inflow station \(Re_{\tau}\),
whereas Mach
is the free-stream Mach number. Moreover, in this case \(\Theta = \frac{T_w - T_{\infty}}{T_r - T_{\infty}}\),
where \(T_r\) is the recovery temperature based on the free-stream temperature \(T_{\infty}\) and Mach number.
[fluid]
This section allows the user to specify fluid properties.
gam
=> real, specific heat ratio (for calorically-perfect gases)Prandtl
=> real, (molecular) Prandtl numbervisc_model
=> integer, viscosity-temperature relation0 => Inviscid flow
1 => Power law
2 => Sutherland law
vt_exp
=> real, exponent for power law viscosity-temperature relations_suth
=> real, Sutherland constant (dimensional temperature value, 110.4 K for air)calorically_perfect
=> integer, flag to specify perfect gas properties0 => thermally perfect gas, variable specific heats
1 => calorically perfect gas, constant specific heats
indx_cp_l
, integer, exponent low limit of \(C_p (T)\) polynomialindx_cp_r
, integer, exponent high limit of \(C_p (T)\) polynomialcp_coeff
, real, list of coefficients for \(C_p (T)\) polynomial (from lowest to highest power), see [8]. If \(T^*\) is the dimensional temperature, the dimensional specific heat at constant pressure is expressed as
[output]
This section allows the user to specify parameters to control the solver output.
print_control
=> integer, frequency for updating the output file progress.outdtsave_restart
=> real, time interval for saving restart filesdtsave
=> real, time interval for saving plot3d and vtk filesdtstat
=> real, time interval for computing statistics (2D and 3D)dtslice
=> real, time interval for writing slicesenable_stat_3d
=> integer, enable computation of 3D statistics (0 -> No, 1 -> Yes)enable_plot3d
=> integer, enable writing of plot3d files (0 -> No, 1 -> Yes)enable_vtk
=> integer, enable writing of vtr files (0 -> No, 1 -> Yes)igslice
=> integer, list of (global) indices for yz-slices, between 1 andnxmax
jgslice
=> integer, list of (global) indices for xz-slices, between 1 andnymax
kgslice
=> integer, list of (global) indices for xy-slices, between 1 andnzmax
list_aux_slice
=> integer, list of variables to print in slices (from 1 to 6)debug memory
=> integer, flag to activate CPU and GPU memory checks (0 -> No, 1 -> Yes)io_type_r
=> integer, specify IO_type for reading (0 => no IO, 1 => serial, 2 => parallel, optional default 2)io_type_w
=> integer, specify IO_type for writing (0 => no IO, 1 => serial, 2 => parallel, optional default 2)enable_slice_plot3d
=> integer, enable writing of slices in plot3d format (0 => no 1 => yes optional, default 0)dtslice_p3d
=> time interval for writing slices in P3D format (optional if not present dtslice_p3d=dtslice)igslice_p3d
=> integer, list of (global) indices for yz-slices in p3d formatjgslice_p3d
=> integer, list of (global) indices for xz-slices in p3d formatkgslice_p3d
=> integer, list of (global) indices for xy-slices in p3d formatdtprobe
=> real, time interval for writing probe signals (optional, if not present dtprobe=dtslice))save_stat_z
=> integer, enable saving of z-averaged statistics (optional, default 0)enable_stat_time
=> integer, enable statstics averaged over time (optional, default 1)enable_tspec
=> integer, enable run-time computation of spectra (0 => no spectra, 1 => restart spectra, 2=> continue spectra, optional default 0)dt_tspec
=> real, time interval for runtime spectra (mandatory if enable_tspec = 1)ndft_tspec
=> integer, total length of a single segment used for the spectra computation using the Welch methodjg_tspec
; integer, j-index (global) of the plane transformed for spectraigstart_tspec
=> integer, starting i-index (global) for spectra computationigend_tspec
=> integer, ending i-index (global) for spectra computation
[curvi]
This section allows the user to specify parameters to control a series of code features exclusivelyu associated to the curvilinear mode
Reynolds_friction
=> real, input friction Reynolds number for curved channel mesh generationl0_ducros
=> real, length scale for Ducros (optional, default l0)angle
=> real, angle of curved channel (radiant)R_curv
=> real, radius of curved channelenable_limiter
=> integer, enable limiter of density/temperature (optional, default 0)rho_lim
real, => minimum accepted value for density (optional default 0)tem_lim
real, => minimum accepted value for temperature (optional default 0)enable_forces_runtime
=> integer, enable airfoil profile force at runtime (optional, default 0)theta_threshold
=> real, bad cell theta level for diagnostics (optional default 100000)xtr1
=> real, minimum x of airfoil trippingxtr2
=> real, maximum x of airfoil trippingxtw1
=> real, minimum x of travelling waves for airfoilxtw2
=> real, maximum x of travelling waves for airfoila_tr
=> real, tripping amplitudea_tw
=> real, travelling wave amplitudedel0_tr
=> real, tripping scale (optional, default 0.001)v_bs
=> real, velocity of blowing/suction for airfoil (optional, default 0.)kx_tw
=> real, wave number of travelling wavesom_tw
=> real, angular frequency of travelling wavesthic
=> real, spatial window thickness of travelling waves/blowing-suction for airfoilteshk
=> real, weno threshold of trailing edgejweno
=> integer (optional, default 0). Force weno application for j> jweno in curvilinear computations.ifilter
=> integer (optional, default 0). Frequency of filtering activation for c2 computations. Used to manage complex transient behaviors with explosive transition to turbulence. When active setng >= 4
.jfilter
=> integer (optional, default 0). When filter is active it is applied from jfilter to nymax.enable_sponge
=> integer (optional, default 0). Sponge along y-max: 0 => no sponge, 1 => relax to w_inf, 2 => relax to initial boundary values.j_sponge
=> integer (optional, default nint(0.9375*nymax)). Minimum index along y where sponge is enabled.k_sponge
=> real (optional, default 0.1). Sponge relaxation intensity.
Notes:
In certain flow scenarios, the abrupt transition from a laminar to a turbulent flow configuration can be challenging to handle from a numerical perspective. As a result, the computation may diverge, leading to a blow-up. To robustly manage such situations, a limiter can be activated to enforce minimum thresholds on the thermodynamic variables (enable_limiter = 1).
[insitu]
enable_insitu
=> integer, 0 => in situ disabled, 1 => in situ enabledinsitu_platform
= string, only “catalyst-v2” supportedvtkpipeline
=> string, name of Python in situ pipeline filedt_insitu
=> time interval between consecutive in situ pipeline callsaux_list
=> comma-separated integer list, indexes of auxiliary variables exported to the in situ pipeline (1 => rho, 2 => u, 3 => v, 4 => w, 5 => H, 6 => T, 7 => viscosity, 8 => ducros, 9 => abs(omega), 10 => divergence)add_list
=> comma-separated integer list, indexes of additional variables exported to the in situ pipeline (1 => div, 2 => abs(omega), 3 => ducros, 4 => swirling strength, 5 => Schlieren)perc_ny_cut
=> percentage of domain along wall normal which is not included in the in situ pipeline (if negative, all domain is included)freeze_intervals
=> list of index pairs defining time iterations (time iteration start and number of freezed iterations) when the fields are freezed, “start1,interval1, start2,interval2, …”
Notes:
In situ analysis requires the code compiled activating CATALYST=”v2” in makefile.inc and linking against Catalyst2 as explained in the compilation section of the guide. Then, a typical pipeline is created by first running a few iterations using “gridwriter.py” predefined pipeline, then opening the produced files using ParaView, creating and exporting the actual pipeline. Finally, the created pipeline is specified in vtkpipeline field and the real simulation is run. In examples/curvcha_moser, the generic pipeline gridwriter.py, a simple pipeline curvchasimple.py and an input file singleideal.ini.insitu are provided to prepare a configuration.
Mesh file
The computational mesh is automatically generated by STREAmS for Cartesian flow cases
and for the curved channel flow case.
For curved boundary layers and C-shape meshes a mesh file grid2d.xyz
in plot3D format has to be provided.
The 2D mesh must have dimensions \(nxmax\times nymax\) and must be provided in the following plot3d format:
open(10,file='grid2d.xyz',form='unformatted',access='stream')
write(10) nxmax,nymax,1
do j=1,nymax
do i=1,nxmax
write(10) xc2g(i,j)
enddo
enddo
do j=1,nymax
do i=1,nxmax
write(10) yc2g(i,j)
enddo
enddo
do j=1,nymax
do i=1,nxmax
write(10) 0.d0
enddo
enddo
close(10)
For C-shape meshes the grid can be generated using Construct2D [14].