Loading [MathJax]/extensions/tex2jax.js
MAIA bb96820c
Multiphysics at AIA
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Structured Finite Volume (FVSTRCTRD)

Overview

m-AIA is divided into two distinct parts, one uses a hierarchical Cartesian cut-cell approach with automatic grid generation techniques. The other part is the structured block environment for solving the Navier-Stokes equation (NSE) on curvilinear boundary fitted grid. The Structured grid needs to be provided beforehand as the solver is not able to create them. This page gives introduction on how to configure your properties files for the structured finite volume (FVSTRCTRD) solver environment.

At present, the FVSTRCTRD solver can be used for both RANS and LES, also in coupled simulation using RANS for some blocks and LES for others. Simulations can be conducted in 2D and 3D for static and moving meshes.

General settings

Some common properties, described in General application control, are also necessary for the FVSTRUCTURED solver, including:

### These are the global properties
nDim = 3
flowSolver = true
timeSteps = 100
restartFile = false
outputDir = 'out/'
solutionOutput = './out/'
gridInputFileName = 'grid.hdf5'
scratchSize = 30
maxNoCells = 10000
### The properties below specify which solvers and couplers are used
multiSolver = 1
noSolvers = 1
executionRecipe = "RECIPE_BASE"
solvertype.default = "MAIA_UNIFIED"
solvertype.0 = "MAIA_STRUCTURED"
### The following are not required for FVSTRUCTURED but the general framework still reads them
gridGenerator = false
maxRfnmntLvl = 8

After this, you can choose LES or RANS

DNS/LES

DNS/Implicit LES is the default case, you do not need to set any properties to turn it on.

To deactivate the viscous flux computation you can use:

govEqs = "EULER"

RANS

To activate RANS, it is necessary to make the following settings

fullRANS = true

Choose the desired RANS method:

ransMethod = "RANS_SA_DV"

Initial condition

The initial condition is specified via a number

initialCondition = 1233

where the specific IC is found in the following table.

initialCondition Description Applicable for RANS
0 Parallel inflow field. Yes
2 Shear flow No
11 Point source in the middle No
43/333 Parallel inflow field with pressure peak in the middle of the domain. No
314 Stagnating flow field \( u = v = w = 0 \), \(\rho = \rho_\infty\), \(p = p_\infty\) No
315 Poiseuille flow No
101 Taylor-Green Vortex \(\rho = \rho_\infty\) No
1234 Laminar channel flow No
1233 Turbulent channel flow with perturbations No
1236 Pipe flow with perturbations No
79092 Approximate mean turbulent boundary layer profile Yes
79091 Turbulent plate with log-law Yes
111/112/113/44 Three-dimensional non-reflective boundary conditions No
777 Falkner Skan Cooke initial condition (incompressible) No
999 Blasius laminar boundary layer No

Boundary condition

At the boundaries of the domain suitable boundary conditions need to be applied. These are specified in the grid file by assigning numbers to each window for each block:

BC

Description

1000

No-slip wall

  • \( u_{wall} = v_{wall} = w_{wall} = 0\)
  • \(\left.\frac{\partial \rho}{\partial n}\right|_{wall} = 0 \)
  • \(\left.\frac{\partial p}{\partial n}\right|_{wall} = 0 \)

When moving grid is switched on this will become 1004.

1003

Isothermal no-slip wall

  • \( u_{wall} = v_{wall} = w_{wall} = 0\)
  • \( T_{wall} = \alpha T_\infty \)
  • \(\left.\frac{\partial p}{\partial n}\right|_{wall} = 0 \)

When moving grid is switched on this will become 1006. For this BC, isothermalWallTemperature in the properties files sets the value of \( \alpha \) which controls the wall temperature as a ratio of the freestream temperature.

1004

Moving no-slip wall

  • \( u_{wall} = u_{mesh}(t)\)
  • \( v_{wall} = v_{mesh}(t)\)
  • \( w_{wall} = w_{mesh}(t)\)
  • \( \left.\frac{\partial \rho}{\partial n}\right|_{wall} = 0 \)
  • \( \left.\frac{\partial p}{\partial n}\right|_{wall} = 0 \).

That is, the velocity of the boundary grid points, according to the prescribed boundary motion, are applied to the fluid.

1006

Moving isothermal no-slip wall

  • \( u_{wall} = u_{mesh}(t)\)
  • \( v_{wall} = v_{mesh}(t)\)
  • \( w_{wall} = w_{mesh}(t)\)
  • \( T_{wall} = \alpha * T_\infty \)
  • \( \left.\frac{\partial p}{\partial n}\right|_{wall} = 0 \).

That is, the velocity of the boundary grid points, according to the prescribed boundary motion, are applied to the fluid. For this BC, isothermalWallTemperature in the properties files sets the value of \( \alpha \) which controls the wall temperature as a ratio of the freestream temperature.

2002

Supersonic inflow

  • \( u = u_{\infty}, v = v_{\infty}, w = w_{\infty}\)
  • \( \rho = \rho_{\infty}, p = p_{\infty}.\)

2003

Simples subsonic in/outflow

Simple extrapolation and prescription

For inflows, i.e., \(M < 0\), the following conditions are used

  • \( u = u_\infty\)
  • \( v = v_\infty\)
  • \( w = w_\infty\)
  • \(\frac{\partial \rho}{\partial n} = 0 \)
  • \(\frac{\partial p}{\partial n} = 0 \)

On outflow boundaries, i.e., where \(M > 0\), the conditions

  • \(\frac{\partial u}{\partial n} = 0 \)
  • \(\frac{\partial v}{\partial n} = 0 \)
  • \(\frac{\partial w}{\partial n} = 0 \)
  • \(\frac{\partial \rho}{\partial n} = 0 \)
  • \( p = p_\infty \)

2004

Subsonic in/outflow

Depending on the flow direction this acts as an in- or outflow, simplified characteristics approach, Whitfield 1984

For inflows, i.e., \(M < 0\), the following conditions are used

\begin{align} p|_{BC} &= \frac{1}{2} \left[p|_i + p_\infty + \frac{\rho_{old} a_{old}}{\left| \nabla \vec{\xi}\right|}\left(\vec{\xi}\cdot\left(\vec{u}|_i - \vec{u}_\infty\right)\right)\right]\\ \rho|_{BC} &= \rho_\infty + \frac{1}{a_{old}^2}\left(p|_{BC} - p_\infty\right)\\ \vec{u}|_{BC} &= \vec{u}_\infty + \frac{\vec{\xi}}{\left|\nabla \vec{\xi}\right|\rho_{old}a_{old}}\left(p|_{BC} - p_\infty\right)\,, \end{align}

On outflow boundaries, i.e., where \(M > 0\), the conditions

\begin{align} p|_{BC} &= p_\infty\\ \rho|_{BC} &= \rho|_i + \frac{1}{a_{old}^2}\left(p|_{BC} - p|_i \right)\\ \vec{u}|_{BC} &= \vec{u}|_i - \frac{\vec{\xi}}{\left|\nabla \vec{\xi}\right|\rho_{old}a_{old}}\left(p|_{BC} - p|_{i}\right) \end{align}

2005

Supersonic outflow

  • \(\frac{\partial u}{\partial n} = 0 \)
  • \(\frac{\partial v}{\partial n} = 0 \)
  • \(\frac{\partial w}{\partial n} = 0 \)
  • \(\frac{\partial p}{\partial n} = 0 \)
  • \(\frac{\partial \rho}{\partial n} = 0 \)

2500,2501

Recycling-rescaling

A combined rescycling-rescaling boundary condition for compressible wall-bounded flows by El-Askary et al.. Takes a slice from the window connected to BC 2501 (needs to be specified in the grid file, typically at some distance downstream of the inflow window), rescales it to the desired boundary layer thickness, and prescribes it to the inflow window with BC 2500 to achieve a turbulent boundary layer flow. The property rescalingBLT in the properties file is used to rescale the boundary layer at the inflow which means the smaller this value is the thinner the boundary layer is. The heart of this method is a means of estimating the velocity at the inlet plane based on the solution downstream. The velocity field is extracted from a plane near the domain exit, which is uesed to rescale, and then it is reintroduced as a boundary condition at the inlet.

2600,2601

External profile

Prescribe a given profile to the window from an external HDF5 file.

4401,4403,4405

Periodic boundary

Periodic boundary condition for two opposing windows. Use the same odd number for two opposing windows, for multiple periodic directions use different numbers, e.g. 4401 for the first two windows (first periodic direction), 4403 for other two opposing windows (second periodic direction), etc.

6000

Multiblock connection

Enables Multiblock communication between user defined blocks in the grid file

7909

Synthetic turbulence generation (STG)

A modified version of the reformulated synthetic turbulence generation(RSTG) method by Roidl et al. is used. Coherent structures are generated by superimposing the velocity fluctuations of virtual eddy cores, which populate a virtual volume, onto the inflow plane.

The synthetic turbulence generation inflow boundary condition is used by setting the inflow window boundary condition number to 7909 and setting the property

useSTG = true

At the initial startup, which is indicated by setting

stgInitialStartup = true

The averaged primitive variables \( (\overline{u}, \overline{v}, \overline{w}, \overline{p}, \overline{\rho})\) and the turbulent viscosity \( \overline{\nu_t}\) need to be in the first layers of cells near the inflow. This can be achieved, e.g., by interpolating the result from a precursor RANS simulation to the flow field of the current LES case. The boundary layer thickness \( \delta_{99} \) of the average velocity profile at the inflow needs to be set via

deltaIn = 7.08

and the number of eddies is set by

stgMaxNoEddies = 100

where the number of eddies should be \( n_{eddies} \approx L_z/\delta_{99}*100 \). The property

stgCreateNewEddies = false

can normally be left to false, only when the generation of a completely new distribution is desired this property can be switched on once.

The properties

BLT1 = 1.0
BLT2 = 1.0
BLT3 = 1.0

control the size of the virtual box by specifying a factor in each space dimension by which the given \(\delta_{99}\) is multiplied to compute the size in the \(x\)- and \(y\)-direction, in the spanwise direction \(z\) the factor is multiplied with the spanwise size of the domain. For instance, using the factors above the virtual box extends one \(\delta_{99}\)-thickness in the \(x\)-drection, centered around the inflow plane, one \(\delta_{99}\)-thickness in the \(y\)-direction, starting at the bottom wall, and covers exactly the spanwise domain extent.

The property

stgExple = 0.5

controls the exponent in the length scale computation, typically a value of 0.5 delivers reasonable results.

Output

In MAIA there are different options to write out your solutions: full solution output, box output, interpolated output, and force field output.

Solution/Restart output

The normal solution output saves the solutions of the whole flow field, these files can also be used as restart files. That is, in the FVSTRCTRD solver, any solution file can be used as a restart file.

property Definition
outputDir Directory to which the files will be written
solutionInterval Interval of output in time steps

To restart the computation use the following properties:

property Definition
restartFile Set true for restart
restartVariablesFileName Name of the restart file

Boxes output

This output can be used to write out 3D subvolumes of the domain. The offsets and sizes of the subvolume are given in computational indices \((i,j,k)\). An example is found in the 3D_prescribingBC_LES testcase. The properties are:

property Definition
boxOutputInterval Number of timesteps between writing out the files
boxBlock Number of boxes to be written out
boxOffsetK cell-offset in k-direction
boxOffsetJ cell-offset in j-direction
boxOffsetI cell-offset in i-direction
boxSizeK size in cells in k-direction
boxSizeJ size in cells in j-direction
boxSizeI size in cells in i-direction
boxWriteCoordinates set true to write out cell-center coordinates
Note
If there is no box output wanted, set the boxOutputInterval to 0.

Forces output

The forces on all surfaces with wall boundary conditions, i.e., any boundary condition \( 1000 \geq 1000 < 2000\), can be written instantaneously to HDF5 files (full 2D resolution) with name auxData????.hdf, where ???? is the time step, or integrated to scalar coefficients to a text file with name forces.?.dat, where ? is the number of wall for which the coefficients were computed. That is, if you have three surfaces with a wall boundary condition three text files forces.0.dat, forces.1.dat, and forces.2.dat are written.

To trigger the 2D output to HDF5 files set the following property:

computeForce = true

This will trigger the force calculation, inclusive skin friction and pressure forces. For moving walls the power consumption can also be computed by setting

computePower = true

The force coefficients, which are computed by integrated over the forces over the wall surface, are written as attributes to the HDF5 files. Averaging the coeffient along \(i\), \(j\), or \(k\) can be triggered by

computeCpLineAverage = false
cpAveragingDir = 0

where cpAveragingDir is set to \(0\rightarrow i\), \(1\rightarrow j\), or \(2\rightarrow k\).

The parameter below will define the output interval of the auxdata HDF5 files.

forceOutputInterval = 10000

To trigger the writeout of the integrated force coefficients to a text file you can set

forceAsciiOutputInterval = 50
forceAsciiComputeInterval = 1

where forceAsciiComputeInterval controls how often the coefficients are computed, e.g., 1 for every time step, and forceAsciiOutputInterval controls how often the data is written out to the text files. For the values above the forces are computed in every time step, but saved internally, and are written out only every 50 time steps to the text files, this procedure will save computational time.

The region where the forces are computed can be limited by

auxDataCoordinateLimits = true
auxDataLimits = [155.0, 205.0, -30.0, 130.0]

The first two number of auxDataLimits control the lower and upper interval boundary in the \(x\)-direction and the latter two control the lower and upper interval boundary in the \(z\)-direction.

Interpolated points output

This output interpolates the cell-centered variables to arbitrary points which are located equidistantly on given lines. Currently, one or two directions can be specified, i.e., 1D and 2D distribution of points can be realized. An example of this can be found in the 3D_cube_periodic_LES testcase. The parameters are:

property Definition
intpPointsOutputInterval Number of timesteps between writing out the files
intpPointsStartX x-Coordinate of the starting point of the line
intpPointsStartY y-Coordinate of the starting point of the line
intpPointsStartZ z-Coordinate of the starting point of the line
intpPointsDeltaX delta in x-Direction between two points on the line
intpPointsDeltaY delta in x-Direction between two points on the line
intpPointsDeltaZ delta in x-Direction between two points on the line
intpPointsDeltaX2D delta in x-Direction between two points on the line for the second direction
intpPointsDeltaY2D delta in x-Direction between two points on the line for the second direction
intpPointsDeltaZ2D delta in x-Direction between two points on the line for the second direction
intpPointsNoPoints delta in x-Direction between two points on the line
intpPointsNoPoints2D delta in x-Direction between two points on the line for the second direction

Moving grid methods

To activate moving grids set

movingGrid = true

The following list of moving grid methods uses hard-coded parameters, which cannot be set in the property file:

gridMovingMethod Description
1 Oscillation in \(i\)-direction, generating a Stokes layer
2 Channel with oscillating indentation
3 Piston moving in \(i\)-direction
4 Inner grid movement. The grid points in the center of the domain are actuated in the \(x\)- and \(y\)-direction with a gradual weakening of the movement towards the boundaries. The mesh at the boundaries is static.

In the following table a number of more complex moving grid functions is described which can be controlled via the property file.

Traveling wave

Spanwise traveling transversal surface wave with parameters defined in inner units (gridMovingMethod = 9) or outer units (gridMovingMethod = 10)

properties Description
gridMovingMethod 9/10
waveLengthPlus/waveLength The wavelength of traveling wave with inner scale/outer scale
waveAmplitudePlus/waveAmplitude The amplitude of traveling wave with inner scale/outer scale
waveTimePlus/waveTime The wavetime of traveling wave with inner scale/outer scale
waveBeginTransition The start of wave spatial transition along i direction, from flat surface to wave surface
waveEndTransition The end of wave spatial transition along i direction, from flat surface to wave surface
waveOutBeginTransition The start of wave spatial transition along i direction, from wave surface to flat surface
waveOutEndTransition The end of wave spatial transition along i direction, from wave surface to flat surface
waveAngle The angle of wave
waveYBeginTransition The start of wave spatial transition along j direction, from wave to flat
waveYEndTransition The end of wave spatial transition along j direction, from wave to flat
waveTemporalTransition The time setting for wave temporal trasition

Oscillating cylinder

This controls the up- and down oscillation of a cylinder, i.e., a circular mesh around a round cylinder is required.

properties Description
gridMovingMethod 14
oscAmplitude The amplitude of oscillation
oscSr ???
oscFreqFactor The frequency factor used to control the oscillating frequency

Sponge

To reduce the numerial reflections from boundaries a damping sponge layer zone can be introduced to absorb outgoing waves.

In the FVSTRCTRD solver the sponge is realized by an additional local damping term, the strength and type of the damping as well as the region where the sponge is applied can be controlled via properties.

Activate the sponge by

useSponge = true

The location where sponges are applied is controlled by

readSpongeFromBC = false
spongeBndryCndIds = [7909]
spongeWindowIds = [1, 2, 4]

If readSpongeFromBC is true the sponge will be added inwards from all windows which have a BC in the list spongeBndryCndsIds. Otherwise, if readSpongeFromBC is false, the list spongeWindowIds will be used to apply the sponge inwards from these specific window ids (as given in the grid).

The type and strength of the sponge are controlled with the following properties

properties Description
spongeLayerType = 6

Type of sponge layer:

  • 1: sponge to \(\rho_\infty\) and \((\rho E)_\infty\)
  • 2: sponge to \(\overline{\rho}\) and \(\overline{\rho E}\) from the fields spongeRho and spongeRhoE from the restart file
  • 3: sponge to \(\rho_\infty\) and \(p_\infty\)
  • 4: sponge to \(\overline{\rho}\) from the field spongeRho from the restart file and \( p_\infty\)
  • 5: sponge to \(p\) from the Falkner-Skan Cooke pressure function
  • 6: sponge to \(p_\infty\)

spongeLayerThickness = [50.0, 30.0, 20.0] Thickness of the sponge layer normal to the window
sigmaSponge = [1.0, 1.0, 1.0] Sponge strength
betaSponge = [2.0, 2.0, 2.0] Controls the decay exponent of the sponge strength inwards from the window. Linear: spongeBeta = 1. Quadratic: spongeBeta = 2.
computeSpongeFactor = false Disables the sponge strength computation. If set to false, the sponge values will be read from the restart file which may be much faster due to the slow sponge computation. That is, the sponge needs to be computed once (first solver run, set computeSpongeFactor = true) and is then saved to the restart file. Following simulation runs (set computeSpongeFactor = false) the sponge factor is read from the restart file.

Postprocessing

The FVSTRCTRD solver uses a standalone postprocessing part, independent from the other Cartesian-grid solvers.

To activate postprocessing ste

postprocessing = true

Averaging of the flow variables can be done while the simulation is running (in solve) without the need to write out snapshots. Otherwise, a set of snapshots written to disk can be read in successively by the solver and averaged. Furthermore, other postprocessing options are available.

The general postprocessing operations which are desired can be controlled by specifiying

postprocessingOps = ["PP_AVERAGE_IN", "PP_TAUW_PRE"]

In the following table are the postprocessing operations which are currently available. Not all operations can be used together with all other operations.

postprocessingOps Description
PP_AVERAGE_PRE Averaging of existing solution files before the simulation start
PP_AVERAGE_POST Averaging of existing solution files after the simulation ends
PP_AVERAGE_IN Averaging of the flow variables during the simulation (in solve)
PP_TAUW_PRE Compute the average skin friction from the primitive variables
PP_LOAD_AVERAGED_SOLUTION_PRE Loads the averaged variables from the Mean*.hdf5 file into the solver and into the primitive variables
PP_SUBTRACT_MEAN Subtracts the mean from the currently loaded solution file, mean file needs to be available
PP_SUBTRACT_PERIODIC_FLUCTUATIONS Special option for traveling wave actuation. The periodic fluctuations (from a triple decomposition) are subtracted from the currently loaded solution file
PP_WRITE_GRADIENTS Computes the spatial gradients of the primitive variables and writes them to an HDF5 file
PP_DECOMPOSE_CF Perform an FIK decomposition
PP_COMPUTE_PRODUCTION_TERMS_PRE Compute production terms
PP_COMPUTE_DISSIPATION_TERMS_PRE Compute dissipation terms from a set of solution files and an existing mean file

General averaging is controlled by the properties

Property Description
pp_averageStartTimestep Timestep at which averaging is started
pp_averageStopTimestep Timestep at which averaging is ended
pp_averageInterval Interval at which probes of the flow variables are added onto the average
pp_fileName Input file name from which a mean file is read if the postprocessing operation requires it

Interpolation

The flow data for grids which do not match the current grid can be interpolated to the current grid with a built-in method. The interpolation procedure is performed at the initial startup. The interpolation is second-order accurate if a cell center on the current grid is surrounded by 8 cell centers of the donor grid. Then a trilinear interpolation, modified for curvilinear grids, is performed. If no surrounding donor points can be found the fallback procedure resorts to a simple nearest neighbor interpolation. Therefore, for a high quality of the interpolation the goal is to have a donor grid which is at least as big as the current grid, such that all points on the current grid are covered.

It is also possible to interpolate from 2D grids to 3D grids.

Here are properties to be set:

Property Description
restartInterpolation = true Turn on this function
donorGridName = "donorGrid.hdf5" The donor grid file
donorVars = "donorVars.hdf5" The donor solution file
donorScale = [1.0,1.0,1.0] Scaling factor of the donor grid
donorTranslation = [0.0,0.0,0.0] Translation of the donor grid

Zonal computations

High-resolution LES of a turbulent flow at high Reynolds numbers consumes high computational power and time due to the necessary fine grid resolution. However, a turbulent flow problem in non-equilibrium state occurs only in a specific area. Applying LES throughout the domain is therefore not an efficient way to solve the turbulent flow. Therefore, it would be elegant to perform a RANS simulation for noninteresting areas in the flow field, while using the LES simultaneously in the important regions where a higher grid resolution is required. The FVSTRCTRD solver offers fully coupled zonal RANS-LES simulations using multiple RANS and LES zones. The grid should contain at least two blocks which should overlap by a certain amount. In the property file the type of closure model (RANS or LES) is then specified for each block. The transition from RANS to LES in the near-wall region is typically done by using the STG boundary condition.

Here are relevant properties:

Property Description
zonal = true To turn on zonal computations
fullRANS = false If you use the zonal approach, fullRANS needs to be false
noRansZones = 1 Number of blocks with RANS
ransZone = [0] Block ids (in the grid file) which should be RANS
ransMethod = "RANS_SA_DV" RANS model to be used
zonalExchangeInterval = 50 Timesteps Interval for exchanging data between RANS and LES zones
zonalAvgFactor = 128 Averaging factor for the moving average in the LES zones

Test cases

For STRCTRDFV solver, there are complete **workshops...** to start with. Furthermore, some ready-to-use testcases are found in the STRUCTURED_FV subfolder of the testcase repository:

cd ~/scratch
svn co http://svn/maia/testcases/STRUCTURED_FV/

References

  • W. El-Aksary, W. Schröder, M. Meinke. LES of compressible wall-bounded flows. AIAA 2003-3554 (2003. 10.2514/6.2003-3554
  • B. Roidl, M. Meinke, W. Schröder. A reformulated synthetic turbulence generation method for a zonal RANS-LES method and its application to zero-pressure gradient boundary layers. Int. J. of Heat Fluid Flow, 44 (2013) 28-40. https://doi.org/10.1016/j.ijheatfluidflow.2013.03.017