HyPar
1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
|
3D Navier Stokes equations (compressible flows) More...
#include <basic.h>
Go to the source code of this file.
Data Structures | |
struct | NavierStokes3D |
Structure containing variables and parameters specific to the 3D Navier Stokes equations. This structure contains the physical parameters, variables, and function pointers specific to the 3D Navier-Stokes equations. More... | |
Macros | |
#define | _NAVIER_STOKES_3D_ "navierstokes3d" |
#define | _MODEL_NDIMS_ 3 |
#define | _MODEL_NVARS_ 5 |
#define | _ROE_ "roe" |
#define | _RF_ "rf-char" |
#define | _LLF_ "llf-char" |
#define | _RUSANOV_ "rusanov" |
#define | _XDIR_ 0 |
#define | _YDIR_ 1 |
#define | _ZDIR_ 2 |
#define | _IB_ADIABATIC_ "adiabatic" |
#define | _IB_ISOTHERMAL_ "isothermal" |
#define | _IB_RAMP_LINEAR_ "linear" |
#define | _IB_RAMP_SMOOTHEDSLAB_ "smoothed_slab" |
#define | _IB_RAMP_DISABLE_ "no_ib" |
#define | _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma) |
#define | _NavierStokes3DSetFlux_(f, stride, rho, vx, vy, vz, e, P, dir) |
#define | _NavierStokes3DRoeAverage_(uavg, stride, uL, uR, gamma) |
#define | _NavierStokes3DSetStiffFlux_(f, stride, rho, vx, vy, vz, e, P, dir, gamma) |
#define | _NavierStokes3DSetNonStiffFlux_(f, rho, vx, vy, vz, e, P, dir, gamma) |
#define | _NavierStokes3DEigenvalues_(u, stride, D, gamma, dir) |
#define | _NavierStokes3DLeftEigenvectors_(u, stride, L, ga, dir) |
#define | _NavierStokes3DRightEigenvectors_(u, stride, R, ga, dir) |
Functions | |
int | NavierStokes3DInitialize (void *, void *) |
int | NavierStokes3DCleanup (void *) |
int | gpuNavierStokes3DCleanup (void *) |
Variables | |
static const int | _NavierStokes3D_stride_ = 1 |
3D Navier Stokes equations (compressible flows)
3D Navier-Stokes equations for viscous and inviscid compressible flows (with gravitational terms)
\begin{equation} \frac {\partial} {\partial t} \left[\begin{array}{c} \rho \\ \rho u \\ \rho v \\ \rho w \\ e \end{array}\right] + \frac {\partial} {\partial x} \left[\begin{array}{c} \rho u \\ \rho u^2 + p \\ \rho u v \\ \rho u w \\ (e+p) u\end{array}\right] + \frac {\partial} {\partial y} \left[\begin{array}{c} \rho v \\ \rho u v \\ \rho v^2 + p \\ \rho v w \\ (e+p) v \end{array}\right] + \frac {\partial} {\partial z} \left[\begin{array}{c} \rho w \\ \rho u w \\ \rho v w \\ \rho w^2 + p \\ (e+p) w \end{array}\right] = \frac {\partial} {\partial x} \left[\begin{array}{c} 0 \\ \tau_{xx} \\ \tau_{yx} \\ \tau_{zx} \\ u \tau_{xx} + v \tau_{yx} + w \tau_{zx} - q_x \end{array}\right] + \frac {\partial} {\partial y} \left[\begin{array}{c} 0 \\ \tau_{xy} \\ \tau_{yy} \\ \tau_{zy} \\ u \tau_{xy} + v \tau_{yy} + w \tau_{zy} - q_y \end{array}\right] + \frac {\partial} {\partial z} \left[\begin{array}{c} 0 \\ \tau_{xz} \\ \tau_{yz} \\ \tau_{zz} \\ u \tau_{xz} + v \tau_{yz} + w \tau_{zz} - q_z \end{array}\right] + \left[\begin{array}{c} 0 \\ -\rho {\bf g}\cdot{\bf \hat{i}} \\ -\rho {\bf g}\cdot{\bf \hat{j}} \\ -\rho {\bf g}\cdot{\bf \hat{k}} \\ -\rho u {\bf g}\cdot{\bf \hat{i}} - \rho v {\bf g}\cdot{\bf \hat{j}} - \rho w {\bf g}\cdot{\bf \hat{k}} \end{array}\right] \end{equation}
where \({\bf g}\) is the gravitational force vector per unit mass, \({\bf \hat{i}},{\bf \hat{j}}, {\bf \hat{k}}\) are the unit vectors along the x, y, and z, the viscous terms are given by
\begin{align} \tau_{ij} &= \frac{\mu}{Re_\infty} \left[ \left( \frac{\partial u_i}{\partial x_j} + \frac{\partial u_j}{\partial x_i}\right) - \frac{2}{3}\frac{\partial u_k}{\partial x_k} \delta_{ij} \right], \\ q_i &= - \frac{\mu}{\left(\gamma-1\right)Re_\infty Pr} \frac{\partial T}{\partial x_i} \end{align}
with \(\mu\) being the viscosity coefficient (computed using Sutherland's law), and the equation of state is
\begin{equation} e = \frac {p} {\gamma-1} + \frac{1}{2} \rho \left(u^2 + v^2 + w^2\right) \end{equation}
References for the governing equations (as well as non-dimensional form):-
Reference for the well-balanced treatment of gravitational source term:
Reference for the partitioning of the flux into its stiff (acoustic) and non-stiff (convective) components:
Definition in file navierstokes3d.h.
struct NavierStokes3D |
Structure containing variables and parameters specific to the 3D Navier Stokes equations. This structure contains the physical parameters, variables, and function pointers specific to the 3D Navier-Stokes equations.
Definition at line 482 of file navierstokes3d.h.
Data Fields | ||
---|---|---|
double | gamma |
Ratio of heat capacities |
double | Re |
Reynolds number |
double | Pr |
Prandtl number |
double | Minf |
Freestream Mach number |
double | C1 |
Sutherlands law constant |
double | C2 |
Sutherlands law constant |
double | grav_x |
acceleration due to gravity in x |
double | grav_y |
acceleration due to gravity in y |
double | grav_z |
acceleration due to gravity in z |
double | rho0 |
reference density at zero altitude for flows with gravity |
double | p0 |
reference pressure at zero altitude for flows with gravity |
double | R |
universal Gas constant |
char | upw_choice[_MAX_STRING_SIZE_] |
choice of upwinding |
double * | grav_field_f |
density variation function ( \(\varrho\)) for hydrostatic equilibrium for flows with gravity |
double * | grav_field_g |
pressure variation function ( \(\varrho\)) for hydrostatic equilibrium for flows with gravity |
double * | fast_jac |
"Fast" Jacobian of the flux function (comprising the acoustic modes) |
double * | solution |
array to store the solution at the beginning of each time step |
char | ib_wall_type[_MAX_STRING_SIZE_] |
Type of immersed boundary wall: isothermal or adiabatic |
double | T_ib_wall |
Immersed body wall temperature, if isothermal |
int | HB |
< Choice of hydrostatic balance for flows with gravity (1 - isothermal equilibrium, 2 - constant potential temperature 3 - stratified atmosphere with a Brunt-Vaisala frequency) |
double | N_bv |
the Brunt-Vaisala frequency for NavierStokes3D::HB = 3 |
char | ib_write_surface_data[_MAX_STRING_SIZE_] |
Flag to indicate whether to analyze and write surface data for immersed body, if present. Applicable only if HyPar::flag_ib is 1 |
double | t_ib_ramp |
Time scale to ramp up the application of immersed boundary conditions, applicable only if HyPar::flag_ib is 1 |
double | t_ib_width |
The "gentleness" with which to ramp up the application of immersed boundary conditions, applicable only if HyPar::flag_ib is 1 |
char | ib_ramp_type[_MAX_STRING_SIZE_] |
Type of ramp up the application of immersed boundary conditions (linear, exponential, etc.), applicable only if HyPar::flag_ib is 1 |
double | ib_T_tol |
Isothermal immersed boundary temperature tolerance: if ghost point temperature differs from wall temperature by more than this factor, set it to the wall temperature - sort of a limiting |
double * | gpu_Q | |
double * | gpu_QDerivX | |
double * | gpu_QDerivY | |
double * | gpu_QDerivZ | |
double * | gpu_FViscous | |
double * | gpu_FDeriv | |
double * | gpu_grav_field_f | |
double * | gpu_grav_field_g | |
double * | gpu_fast_jac | |
double * | gpu_solution |
#define _NAVIER_STOKES_3D_ "navierstokes3d" |
3D Navier-Stokes equations
Definition at line 50 of file navierstokes3d.h.
#define _MODEL_NDIMS_ 3 |
Number of spatial dimensions
Definition at line 56 of file navierstokes3d.h.
#define _MODEL_NVARS_ 5 |
Number of vector components at each grid point
Definition at line 58 of file navierstokes3d.h.
#define _ROE_ "roe" |
Roe's upwinding scheme
Definition at line 62 of file navierstokes3d.h.
#define _RF_ "rf-char" |
Characteristic-based Roe-fixed scheme
Definition at line 64 of file navierstokes3d.h.
#define _LLF_ "llf-char" |
Characteristic-based local Lax-Friedrich scheme
Definition at line 66 of file navierstokes3d.h.
#define _RUSANOV_ "rusanov" |
Rusanov's upwinding scheme
Definition at line 68 of file navierstokes3d.h.
#define _XDIR_ 0 |
dimension corresponding to the x spatial dimension
Definition at line 72 of file navierstokes3d.h.
#define _YDIR_ 1 |
dimension corresponding to the y spatial dimension
Definition at line 74 of file navierstokes3d.h.
#define _ZDIR_ 2 |
dimension corresponding to the z spatial dimension
Definition at line 76 of file navierstokes3d.h.
#define _IB_ADIABATIC_ "adiabatic" |
adiabatic immersed body wall
Definition at line 80 of file navierstokes3d.h.
#define _IB_ISOTHERMAL_ "isothermal" |
isothermal immersed body wall
Definition at line 82 of file navierstokes3d.h.
#define _IB_RAMP_LINEAR_ "linear" |
linear ramping
Definition at line 86 of file navierstokes3d.h.
#define _IB_RAMP_SMOOTHEDSLAB_ "smoothed_slab" |
smoothed slab ramp (like tanh)
Definition at line 88 of file navierstokes3d.h.
#define _IB_RAMP_DISABLE_ "no_ib" |
disable the immersed boundaries
Definition at line 90 of file navierstokes3d.h.
#define _NavierStokes3DGetFlowVar_ | ( | u, | |
stride, | |||
rho, | |||
vx, | |||
vy, | |||
vz, | |||
e, | |||
P, | |||
gamma | |||
) |
Get the flow variables from the conserved solution vector.
\begin{equation} {\bf u} = \left[\begin{array}{c} \rho \\ \rho u \\ \rho v \\ \rho w \\ e \end{array}\right] \end{equation}
Definition at line 98 of file navierstokes3d.h.
#define _NavierStokes3DSetFlux_ | ( | f, | |
stride, | |||
rho, | |||
vx, | |||
vy, | |||
vz, | |||
e, | |||
P, | |||
dir | |||
) |
Compute the flux vector, given the flow variables
\begin{eqnarray} dir = x, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \rho u \\ \rho u^2 + p \\ \rho u v \\ \rho u w \\ (e+p)u \end{array}\right], \\ dir = y, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \rho v \\ \rho u v \\ \rho v^2 + p \\ \rho v w \\ (e+p)v \end{array}\right], \\ dir = z, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \rho w \\ \rho u w \\ \rho v w \\ \rho w^2 + p \\ (e+p)w \end{array}\right] \end{eqnarray}
Definition at line 118 of file navierstokes3d.h.
#define _NavierStokes3DRoeAverage_ | ( | uavg, | |
stride, | |||
uL, | |||
uR, | |||
gamma | |||
) |
Compute the Roe-average of two solutions.
Definition at line 144 of file navierstokes3d.h.
#define _NavierStokes3DSetStiffFlux_ | ( | f, | |
stride, | |||
rho, | |||
vx, | |||
vy, | |||
vz, | |||
e, | |||
P, | |||
dir, | |||
gamma | |||
) |
Compute the stiff flux vector (comprising the acoustic modes only), given the flow variables
\begin{eqnarray} dir = x, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \frac{1}{\gamma}\rho u \\ \frac{1}{\gamma}\rho u^2 + p \\ \frac{1}{\gamma}\rho u v \\ \frac{1}{\gamma}\rho u w \\ (e+p)u - \frac{1}{2} \frac{\gamma-1}{\gamma}\rho\left(u^2+v^2+w^2\right)u \end{array}\right], \\ dir = y, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \frac{1}{\gamma}\rho v \\ \frac{1}{\gamma}\rho u v \\ \frac{1}{\gamma}\rho v^2 + p \\ \frac{1}{\gamma}\rho v w \\ (e+p)v - \frac{1}{2} \frac{\gamma-1}{\gamma}\rho\left(u^2+v^2+w^2\right)v \end{array}\right], \\ dir = z, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \frac{1}{\gamma}\rho w \\ \frac{1}{\gamma}\rho u w \\ \frac{1}{\gamma}\rho v w \\ \frac{1}{\gamma}\rho w^2 + p \\ (e+p)w - \frac{1}{2} \frac{\gamma-1}{\gamma}\rho\left(u^2+v^2+w^2\right)w \end{array}\right], \\ \end{eqnarray}
Reference:
Definition at line 183 of file navierstokes3d.h.
#define _NavierStokes3DSetNonStiffFlux_ | ( | f, | |
rho, | |||
vx, | |||
vy, | |||
vz, | |||
e, | |||
P, | |||
dir, | |||
gamma | |||
) |
Compute the non-stiff flux vector (comprising the entropy modes only), given the flow variables
\begin{eqnarray} dir = x, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \frac{\gamma-1}{\gamma}\rho u \\ \frac{\gamma-1}{\gamma}\rho u^2 \\ \frac{\gamma-11}{\gamma}\rho u v \\ \frac{\gamma-1}{\gamma}\rho u w \\ \frac{1}{2} \frac{\gamma-1}{\gamma}\rho\left(u^2+v^2+w^2\right)u \end{array}\right], \\ dir = y, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \frac{\gamma-1}{\gamma}\rho v \\ \frac{\gamma-1}{\gamma}\rho u v \\ \frac{\gamma-1}{\gamma}\rho v^2 \\ \frac{\gamma-1}{\gamma}\rho v w \\ \frac{1}{2} \frac{\gamma-1}{\gamma}\rho\left(u^2+v^2+w^2\right)v \end{array}\right], \\ dir = z, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \frac{\gamma-1}{\gamma}\rho w \\ \frac{\gamma-1}{\gamma}\rho u w \\ \frac{\gamma-1}{\gamma}\rho v w \\ \frac{\gamma-1}{\gamma}\rho w^2 \\ \frac{1}{2} \frac{\gamma-1}{\gamma}\rho\left(u^2+v^2+w^2\right)w \end{array}\right], \\ \end{eqnarray}
Reference:
Definition at line 219 of file navierstokes3d.h.
#define _NavierStokes3DEigenvalues_ | ( | u, | |
stride, | |||
D, | |||
gamma, | |||
dir | |||
) |
Compute the eigenvalues, given a solution vector in terms of the conserved variables. The eigenvalues are returned as a matrix D whose diagonal values are the eigenvalues. Admittedly, this is inefficient. The matrix D is stored in a row-major format.
The eigenvalues are \( {\bf v}\cdot\hat{\bf n}, {\bf v}\cdot\hat{\bf n}, {\bf v}\cdot\hat{\bf n}, {\bf v}\cdot\hat{\bf n} \pm c\) where \({\bf v} = u\hat{\bf i} + v\hat{\bf j} + w\hat{\bf k}\) is the velocity vector, \(\hat{\bf n} = \hat{\bf i}, \hat{\bf j}, \hat{\bf k}\) is the unit vector along the spatial dimension, and \(c\) is the speed of sound. Note that the order of the eigenvalues (and therefore, the order of the eigenvectors in _NavierStokes3DLeftEigenvectors_ and _NavierStokes3DRightEigenvectors_ has been chosen to avoid left eigen-matrices with zero on the diagonals.
Definition at line 254 of file navierstokes3d.h.
#define _NavierStokes3DLeftEigenvectors_ | ( | u, | |
stride, | |||
L, | |||
ga, | |||
dir | |||
) |
Compute the left eigenvectors, given a solution vector in terms of the conserved variables. The eigenvectors are returned as a matrix L whose rows correspond to each eigenvector. The matrix L is stored in the row-major format.
Reference:
Definition at line 288 of file navierstokes3d.h.
#define _NavierStokes3DRightEigenvectors_ | ( | u, | |
stride, | |||
R, | |||
ga, | |||
dir | |||
) |
Compute the right eigenvectors, given a solution vector in terms of the conserved variables. The eigenvectors are returned as a matrix R whose columns correspond to each eigenvector. The matrix R is stored in the row-major format.
Reference:
Definition at line 384 of file navierstokes3d.h.
int NavierStokes3DInitialize | ( | void * | s, |
void * | m | ||
) |
Initialize the 3D Navier-Stokes (NavierStokes3D) module: Sets the default parameters, read in and set physics-related parameters, and set the physics-related function pointers in HyPar.
This file reads the file "physics.inp" that must have the following format:
begin <keyword> <value> <keyword> <value> <keyword> <value> ... <keyword> <value> end
where the list of keywords are:
Keyword name | Type | Variable | Default value |
---|---|---|---|
gamma | double | NavierStokes3D::gamma | 1.4 |
Pr | double | NavierStokes3D::Pr | 0.72 |
Re | double | NavierStokes3D::Re | -1 |
Minf | double | NavierStokes3D::Minf | 1.0 |
gravity | double,double,double | NavierStokes3D::grav_x,NavierStokes3D::grav_y,NavierStokes3D::grav_z | 0.0,0.0,0.0 |
rho_ref | double | NavierStokes3D::rho0 | 1.0 |
p_ref | double | NavierStokes3D::p0 | 1.0 |
HB | int | NavierStokes3D::HB | 1 |
R | double | NavierStokes3D::R | 1.0 |
upwinding | char[] | NavierStokes3D::upw_choice | "roe" (_ROE_) |
ib_wall_type | char[] | NavierStokes3D::ib_wall_type | "adiabatic" (_IB_ADIABATIC_) |
ib_ramp_time | double | NavierStokes3D::t_ib_ramp | -1 |
ib_ramp_width | double | NavierStokes3D::t_ib_width | 0.05 |
ib_ramp_type | char[] | NavierStokes3D::ib_ramp_type | "linear" (_IB_RAMP_LINEAR_) |
if "ib_wall_type" (NavierStokes3D::ib_wall_type) is specified as "isothermal", it should be followed by the wall temperature (#NavierStokes3D::T_ib_wall), i.e,
begin ... ib_wall_type isothermal 1.0 ... end
If "HB" (NavierStokes3D::HB) is specified as 3, it should be followed by the the Brunt-Vaisala frequency (NavierStokes3D::N_bv), i.e.
begin ... HB 3 0.01 ... end
Note: "physics.inp" is optional; if absent, default values will be used.
s | Solver object of type HyPar |
m | MPI object of type MPIVariables |
Definition at line 116 of file NavierStokes3DInitialize.c.
int NavierStokes3DCleanup | ( | void * | s | ) |
Function to clean up all allocations in the 3D Navier Stokes module.
s | Object of type NavierStokes3D |
Definition at line 11 of file NavierStokes3DCleanup.c.
int gpuNavierStokes3DCleanup | ( | void * | s | ) |
Function to clean up all allocations in the 3D Navier Stokes module.
s | Object of type NavierStokes3D |
Definition at line 13 of file NavierStokes3DCleanup_GPU.c.
|
static |
Definition at line 559 of file navierstokes3d.h.