HyPar
1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
|
Initialization of the physics-related variables and function pointers for the 3D Navier-Stokes system. More...
#include <float.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <boundaryconditions.h>
#include <physicalmodels/navierstokes3d.h>
#include <mpivars.h>
#include <hypar.h>
Go to the source code of this file.
Functions | |
double | NavierStokes3DComputeCFL (void *, void *, double, double) |
int | NavierStokes3DFlux (double *, double *, int, void *, double) |
int | NavierStokes3DStiffFlux (double *, double *, int, void *, double) |
int | NavierStokes3DNonStiffFlux (double *, double *, int, void *, double) |
int | NavierStokes3DRoeAverage (double *, double *, double *, void *) |
int | NavierStokes3DParabolicFunction (double *, double *, void *, void *, double) |
int | NavierStokes3DSource (double *, double *, void *, void *, double) |
int | NavierStokes3DJacobian (double *, double *, void *, int, int, int) |
int | NavierStokes3DStiffJacobian (double *, double *, void *, int, int, int) |
int | NavierStokes3DLeftEigenvectors (double *, double *, void *, int) |
int | NavierStokes3DRightEigenvectors (double *, double *, void *, int) |
int | NavierStokes3DUpwindRoe (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | NavierStokes3DUpwindRF (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | NavierStokes3DUpwindLLF (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | NavierStokes3DUpwindRusanov (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | NavierStokes3DUpwinddFRoe (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | NavierStokes3DUpwindFdFRoe (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | NavierStokes3DUpwindRusanovModified (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | NavierStokes3DUpwinddFRusanovModified (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | NavierStokes3DUpwindFdFRusanovModified (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | NavierStokes3DGravityField (void *, void *) |
int | NavierStokes3DModifiedSolution (double *, double *, int, void *, void *, double) |
int | NavierStokes3DIBAdiabatic (void *, void *, double *, double) |
int | NavierStokes3DIBIsothermal (void *, void *, double *, double) |
int | NavierStokes3DPreStep (double *, void *, void *, double) |
int | NavierStokes3DIBForces (void *, void *, double) |
int | gpuNavierStokes3DInitialize (void *, void *) |
int | gpuNavierStokes3DFlux (double *, double *, int, void *, double) |
int | gpuNavierStokes3DParabolicFunction (double *, double *, void *, void *, double) |
int | gpuNavierStokes3DSource (double *, double *, void *, void *, double) |
int | gpuNavierStokes3DModifiedSolution (double *, double *, int, void *, void *, double) |
int | gpuNavierStokes3DUpwindRusanov (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | gpuNavierStokes3DUpwindRoe (double *, double *, double *, double *, double *, double *, int, void *, double) |
int | gpuNavierStokes3DPreStep (double *, void *, void *, double) |
int | NavierStokes3DInitialize (void *s, void *m) |
Initialization of the physics-related variables and function pointers for the 3D Navier-Stokes system.
Definition in file NavierStokes3DInitialize.c.
double NavierStokes3DComputeCFL | ( | void * | s, |
void * | m, | ||
double | dt, | ||
double | t | ||
) |
Computes the maximum CFL number over the domain. Note that the CFL is computed over the local domain on this processor only.
s | Solver object of type HyPar |
m | MPI object of type MPIVariables |
dt | Time step size for which to compute the CFL |
t | Time |
Definition at line 16 of file NavierStokes3DComputeCFL.c.
int NavierStokes3DFlux | ( | double * | f, |
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
Compute the hyperbolic flux function for the 3D Navier-Stokes equations:
\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}
Note: the flux function needs to be computed at the ghost points as well.
f | Array to hold the computed flux vector (same layout as u) |
u | Array with the solution vector |
dir | Spatial dimension (x, y, or z) for which to compute the flux |
s | Solver object of type HyPar |
t | Current simulation time |
Definition at line 23 of file NavierStokes3DFlux.c.
int NavierStokes3DStiffFlux | ( | double * | f, |
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
Compute the stiff flux, given the solution vector. The stiff flux is the component of the total flux that represents the acoustic modes (see _NavierStokes3DSetStiffFlux_). Here, the linearized approximation to the stiff flux is computed as:
\begin{equation} {\bf f}\left({\bf u}\right) = A_f{\bf u} \end{equation}
where \(A_f = A_f\left({\bf u}_{ref}\right)\) is the fast Jacobian (NavierStokes3D::fast_jac) evaluated for the solution at the beginning of each time step ( \({\bf u}_{ref}\) is NavierStokes3D::solution). This is done in NavierStokes3DPreStep().
Note: the flux function needs to be computed at the ghost points as well.
f | Array to hold the computed flux vector (same layout as u) |
u | Array with the solution vector |
dir | Spatial dimension (x,y, or z) for which to compute the flux |
s | Solver object of type HyPar |
t | Current simulation time |
Definition at line 70 of file NavierStokes3DFlux.c.
int NavierStokes3DNonStiffFlux | ( | double * | f, |
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
Compute the non-stiff flux, given the solution vector. The non-stiff flux is the component of the total flux that represents the entropy modes (see _NavierStokes3DSetNonStiffFlux_). Here, the linearized approximation to the non-stiff flux is computed as:
\begin{equation} {\bf f}\left({\bf u}\right) - A_f{\bf u} \end{equation}
where \({\bf f}\left({\bf u}\right)\) is the total flux computed in NavierStokes3DFlux(), and \(A_f{\bf u}\) is the linearized stiff flux computed in NavierStokes3DStiffFlux().
Note: the flux function needs to be computed at the ghost points as well.
f | Array to hold the computed flux vector (same layout as u) |
u | Array with the solution vector |
dir | Spatial dimension (x,y, or z) for which to compute the flux |
s | Solver object of type HyPar |
t | Current simulation time |
Definition at line 111 of file NavierStokes3DFlux.c.
int NavierStokes3DRoeAverage | ( | double * | uavg, |
double * | uL, | ||
double * | uR, | ||
void * | p | ||
) |
Compute the Roe-averaged state for the 3D Navier Stokes equations. This function just calls the macro _NavierStokes3DRoeAverage_ and is not used by any functions within the 3D Navier Stokes module. However, it's necessary to define it and provide it to the the solver object (HyPar) so that it can then send it to interpolation functions for a characteristic-based reconstruction.
uavg | The computed Roe-averaged state |
uL | Left state (conserved variables) |
uR | Right state (conserved variables) |
p | Object of type NavierStokes3D with physics-related variables |
Definition at line 18 of file NavierStokes3DFunctions.c.
int NavierStokes3DParabolicFunction | ( | double * | par, |
double * | u, | ||
void * | s, | ||
void * | m, | ||
double | t | ||
) |
Compute the viscous terms in the 3D Navier Stokes equations: this function computes the following:
\begin{equation} \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] \end{equation}
where
\begin{align} \tau_{xx} &= \frac{2}{3}\left(\frac{\mu}{Re}\right)\left(2\frac{\partial u}{\partial x} - \frac{\partial v}{\partial y} - \frac{\partial w}{\partial z}\right),\\ \tau_{xy} &= \left(\frac{\mu}{Re}\right)\left(\frac{\partial u}{\partial y} + \frac{\partial v}{\partial x}\right),\\ \tau_{xz} &= \left(\frac{\mu}{Re}\right)\left(\frac{\partial u}{\partial z} + \frac{\partial w}{\partial x}\right),\\ \tau_{yx} &= \left(\frac{\mu}{Re}\right)\left(\frac{\partial u}{\partial y} + \frac{\partial v}{\partial x}\right),\\ \tau_{yy} &= \frac{2}{3}\left(\frac{\mu}{Re}\right)\left(-\frac{\partial u}{\partial x} +2\frac{\partial v}{\partial y} - \frac{\partial w}{\partial z}\right),\\ \tau_{yz} &= \left(\frac{\mu}{Re}\right)\left(\frac{\partial w}{\partial y} + \frac{\partial v}{\partial z}\right),\\ \tau_{zx} &= \left(\frac{\mu}{Re}\right)\left(\frac{\partial u}{\partial z} + \frac{\partial w}{\partial x}\right),\\ \tau_{zy} &= \left(\frac{\mu}{Re}\right)\left(\frac{\partial v}{\partial z} + \frac{\partial w}{\partial y}\right),\\ \tau_{zz} &= \frac{2}{3}\left(\frac{\mu}{Re}\right)\left(-\frac{\partial u}{\partial x} - \frac{\partial v}{\partial y} + 2\frac{\partial v}{\partial y}\right),\\ q_x &= -\frac{mu}{\left(\gamma-1\right)Re Pr}\frac{\partial T}{\partial x}, \\ q_y &= -\frac{mu}{\left(\gamma-1\right)Re Pr}\frac{\partial T}{\partial y}, \\ q_z &= -\frac{mu}{\left(\gamma-1\right)Re Pr}\frac{\partial T}{\partial z} \end{align}
and the temperature is \(T = \gamma p/\rho\). \(Re\) and \(Pr\) are the Reynolds and Prandtl numbers, respectively. Note that this function computes the entire parabolic term, and thus bypasses HyPar's parabolic function calculation interfaces. NavierStokes3DInitialize() assigns this function to HyPar::ParabolicFunction.
Reference:
par | Array to hold the computed viscous terms |
u | Solution vector array |
s | Solver object of type HyPar |
m | MPI object of type MPIVariables |
t | Current simulation time |
Definition at line 50 of file NavierStokes3DParabolicFunction.c.
int NavierStokes3DSource | ( | double * | source, |
double * | u, | ||
void * | s, | ||
void * | m, | ||
double | t | ||
) |
Computes the gravitational source term using a well-balanced formulation: the source term is rewritten as follows:
\begin{equation} \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] = \left[\begin{array}{ccccc} 0 & p_0 \varrho^{-1} & 0 & 0 & p_0 u \varrho^{-1} \end{array}\right] \cdot \frac{\partial}{\partial x}\left[\begin{array}{c} 0 \\ \varphi \\ 0 \\ 0 \\ \varphi \end{array}\right] + \left[\begin{array}{ccccc} 0 & 0 & p_0 \varrho^{-1} & 0 & p_0 v \varrho^{-1} \end{array}\right] \cdot \frac{\partial}{\partial y}\left[\begin{array}{c} 0 \\ 0 \\ \varphi \\ 0 \\ \varphi \end{array}\right] + \left[\begin{array}{ccccc} 0 & 0 & 0 & p_0 \varrho^{-1} & p_0 w \varrho^{-1} \end{array}\right] \cdot \frac{\partial}{\partial y}\left[\begin{array}{c} 0 \\ 0 \\ 0 \\ \varphi \\ \varphi \end{array}\right] \end{equation}
where \(\varphi = \varphi\left(x,y\right)\) and \(\varrho = \varrho\left(x,y\right)\) are computed in NavierStokes3DGravityField(). The derivatives are computed in an exactly identical manner as the hyperbolic flux (i.e., using a conservative finite-difference formulation) (see HyperbolicFunction()).
References:
source | Array to hold the computed source |
u | Solution vector array |
s | Solver object of type HyPar |
m | MPI object of type MPIVariables |
t | Current simulation time |
Definition at line 38 of file NavierStokes3DSource.c.
int NavierStokes3DJacobian | ( | double * | Jac, |
double * | u, | ||
void * | p, | ||
int | dir, | ||
int | nvars, | ||
int | upw | ||
) |
Function to compute the flux Jacobian of the 3D Navier-Stokes equations, given the solution at a grid point. The Jacobian is square matrix of size nvar=5, and is returned as a 1D array (double) of 25 elements in row-major format.
Jac | Jacobian matrix: 1D array of size nvar^2 = 25 |
u | solution at a grid point (array of size nvar = 5) |
p | object containing the physics-related parameters |
dir | dimension (0 -> x, 1 -> y, 2 -> z) |
nvars | number of vector components |
upw | 0 -> send back complete Jacobian, 1 -> send back Jacobian of right(+)-moving flux, -1 -> send back Jacobian of left(-)-moving flux |
Definition at line 15 of file NavierStokes3DJacobian.c.
int NavierStokes3DStiffJacobian | ( | double * | Jac, |
double * | u, | ||
void * | p, | ||
int | dir, | ||
int | nvars, | ||
int | upw | ||
) |
Function to compute the Jacobian of the fast flux (representing the acoustic waves) of the 3D Navier-Stokes equations, given the solution at a grid point. The Jacobian is square matrix of size nvar=5, and is returned as a 1D array (double) of 25 elements in row-major format.
Jac | Jacobian matrix: 1D array of size nvar^2 = 25 |
u | solution at a grid point (array of size nvar = 5) |
p | object containing the physics-related parameters |
dir | dimension (0 -> x, 1 -> y, 2 -> z) |
nvars | number of vector components |
upw | 0 -> send back complete Jacobian, 1 -> send back Jacobian of right(+)-moving flux, -1 -> send back Jacobian of left(-)-moving flux |
Definition at line 53 of file NavierStokes3DJacobian.c.
int NavierStokes3DLeftEigenvectors | ( | double * | u, |
double * | L, | ||
void * | p, | ||
int | dir | ||
) |
Compute the left eigenvections for the 3D Navier Stokes equations. This function just calls the macro _NavierStokes3DLeftEigenvectors_ and is not used by any functions within the 3D Navier Stokes module. However, it's necessary to define it and provide it to the the solver object (HyPar) so that it can then send it to interpolation functions for a characteristic-based reconstruction.
u | Conserved solution at a grid point |
L | Array of size nvar^2 = 5^2 to save the matrix of left eigenvectors in (row-major format). |
p | Object of type NavierStokes3D with physics-related variables |
dir | Spatial dimension (x, y, or z) |
Definition at line 20 of file NavierStokes3DEigen.c.
int NavierStokes3DRightEigenvectors | ( | double * | u, |
double * | R, | ||
void * | p, | ||
int | dir | ||
) |
Compute the right eigenvections for the 3D Navier Stokes equations. This function just calls the macro _NavierStokes3DRightEigenvectors_ and is not used by any functions within the 3D Navier Stokes module. However, it's necessary to define it and provide it to the the solver object (HyPar) so that it can then send it to interpolation functions for a characteristic-based reconstruction.
u | Conserved solution at a grid point |
R | Array of size nvar^2 = 5^2 to save the matrix of right eigenvectors in (row-major format). |
p | Object of type NavierStokes3D with physics-related variables |
dir | Spatial dimension (x, y, or z) |
Definition at line 39 of file NavierStokes3DEigen.c.
int NavierStokes3DUpwindRoe | ( | double * | fI, |
double * | fL, | ||
double * | fR, | ||
double * | uL, | ||
double * | uR, | ||
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
Roe's upwinding scheme.
\begin{equation} {\bf f}_{j+1/2} = \frac{1}{2}\left[ {\bf f}_{j+1/2}^L + {\bf f}_{j+1/2}^R - \left| A\left({\bf u}_{j+1/2}^L,{\bf u}_{j+1/2}^R\right) \right| \left( {\bf u}_{j+1/2}^R - {\bf u}_{j+1/2}^L \right)\right] \end{equation}
This upwinding scheme is modified for the balanced discretization of the 3D Navier Stokes equations when there is a non-zero gravitational force. See the reference below. For flows without any gravitational forces, it reduces to its original form.
fI | Computed upwind interface flux |
fL | Left-biased reconstructed interface flux |
fR | Right-biased reconstructed interface flux |
uL | Left-biased reconstructed interface solution |
uR | Right-biased reconstructed interface solution |
u | Cell-centered solution |
dir | Spatial dimension (x, y, or z) |
s | Solver object of type HyPar |
t | Current solution time |
Definition at line 40 of file NavierStokes3DUpwind.c.
int NavierStokes3DUpwindRF | ( | double * | fI, |
double * | fL, | ||
double * | fR, | ||
double * | uL, | ||
double * | uR, | ||
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
Characteristic-based Roe-fixed upwinding scheme.
\begin{align} \alpha_{j+1/2}^{k,L} &= \sum_{k=1}^3 {\bf l}_{j+1/2}^k \cdot {\bf f}_{j+1/2}^{k,L}, \\ \alpha_{j+1/2}^{k,R} &= \sum_{k=1}^3 {\bf l}_{j+1/2}^k \cdot {\bf f}_{j+1/2}^{k,R}, \\ v_{j+1/2}^{k,L} &= \sum_{k=1}^3 {\bf l}_{j+1/2}^k \cdot {\bf u}_{j+1/2}^{k,L}, \\ v_{j+1/2}^{k,R} &= \sum_{k=1}^3 {\bf l}_{j+1/2}^k \cdot {\bf u}_{j+1/2}^{k,R}, \\ \alpha_{j+1/2}^k &= \left\{ \begin{array}{cc} \alpha_{j+1/2}^{k,L} & {\rm if}\ \lambda_{j,j+1/2,j+1} > 0 \\ \alpha_{j+1/2}^{k,R} & {\rm if}\ \lambda_{j,j+1/2,j+1} < 0 \\ \frac{1}{2}\left[ \alpha_{j+1/2}^{k,L} + \alpha_{j+1/2}^{k,R} - \left(\max_{\left[j,j+1\right]} \lambda\right) \left( v_{j+1/2}^{k,R} - v_{j+1/2}^{k,L} \right) \right] & {\rm otherwise} \end{array}\right., \\ {\bf f}_{j+1/2} &= \sum_{k=1}^3 \alpha_{j+1/2}^k {\bf r}_{j+1/2}^k \end{align}
where \({\bf l}\), \({\bf r}\), and \(\lambda\) are the left-eigenvectors, right-eigenvectors and eigenvalues. The subscripts denote the grid locations.
Note that this upwinding scheme cannot be used for solving flows with non-zero gravitational forces.
fI | Computed upwind interface flux |
fL | Left-biased reconstructed interface flux |
fR | Right-biased reconstructed interface flux |
uL | Left-biased reconstructed interface solution |
uR | Right-biased reconstructed interface solution |
u | Cell-centered solution |
dir | Spatial dimension (x, y, or z) |
s | Solver object of type HyPar |
t | Current solution time |
Definition at line 140 of file NavierStokes3DUpwind.c.
int NavierStokes3DUpwindLLF | ( | double * | fI, |
double * | fL, | ||
double * | fR, | ||
double * | uL, | ||
double * | uR, | ||
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
Characteristic-based local Lax-Friedrich upwinding scheme.
\begin{align} \alpha_{j+1/2}^{k,L} &= \sum_{k=1}^3 {\bf l}_{j+1/2}^k \cdot {\bf f}_{j+1/2}^{k,L}, \\ \alpha_{j+1/2}^{k,R} &= \sum_{k=1}^3 {\bf l}_{j+1/2}^k \cdot {\bf f}_{j+1/2}^{k,R}, \\ v_{j+1/2}^{k,L} &= \sum_{k=1}^3 {\bf l}_{j+1/2}^k \cdot {\bf u}_{j+1/2}^{k,L}, \\ v_{j+1/2}^{k,R} &= \sum_{k=1}^3 {\bf l}_{j+1/2}^k \cdot {\bf u}_{j+1/2}^{k,R}, \\ \alpha_{j+1/2}^k &= \frac{1}{2}\left[ \alpha_{j+1/2}^{k,L} + \alpha_{j+1/2}^{k,R} - \left(\max_{\left[j,j+1\right]} \lambda\right) \left( v_{j+1/2}^{k,R} - v_{j+1/2}^{k,L} \right) \right], \\ {\bf f}_{j+1/2} &= \sum_{k=1}^3 \alpha_{j+1/2}^k {\bf r}_{j+1/2}^k \end{align}
where \({\bf l}\), \({\bf r}\), and \(\lambda\) are the left-eigenvectors, right-eigenvectors and eigenvalues. The subscripts denote the grid locations.
This upwinding scheme is modified for the balanced discretization of the 3D Navier Stokes equations when there is a non-zero gravitational force. See the reference below. For flows without any gravitational forces, it reduces to its original form.
fI | Computed upwind interface flux |
fL | Left-biased reconstructed interface flux |
fR | Right-biased reconstructed interface flux |
uL | Left-biased reconstructed interface solution |
uR | Right-biased reconstructed interface solution |
u | Cell-centered solution |
dir | Spatial dimension (x, y, or z) |
s | Solver object of type HyPar |
t | Current solution time |
Definition at line 244 of file NavierStokes3DUpwind.c.
int NavierStokes3DUpwindRusanov | ( | double * | fI, |
double * | fL, | ||
double * | fR, | ||
double * | uL, | ||
double * | uR, | ||
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
Rusanov's upwinding scheme.
\begin{equation} {\bf f}_{j+1/2} = \frac{1}{2}\left[ {\bf f}_{j+1/2}^L + {\bf f}_{j+1/2}^R - \max_{j,j+1} \nu_j \left( {\bf u}_{j+1/2}^R - {\bf u}_{j+1/2}^L \right)\right] \end{equation}
where \(\nu = c + \left|u\right|\).
This upwinding scheme is modified for the balanced discretization of the 3D Navier Stokes equations when there is a non-zero gravitational force. See the reference below. For flows without any gravitational forces, it reduces to its original form.
fI | Computed upwind interface flux |
fL | Left-biased reconstructed interface flux |
fR | Right-biased reconstructed interface flux |
uL | Left-biased reconstructed interface solution |
uR | Right-biased reconstructed interface solution |
u | Cell-centered solution |
dir | Spatial dimension (x,y, or z) |
s | Solver object of type HyPar |
t | Current solution time |
Definition at line 349 of file NavierStokes3DUpwind.c.
int NavierStokes3DUpwinddFRoe | ( | double * | fI, |
double * | fL, | ||
double * | fR, | ||
double * | uL, | ||
double * | uR, | ||
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
The Roe upwinding scheme (NavierStokes3DUpwindRoe) for the partitioned hyperbolic flux that comprises of the acoustic waves only (see NavierStokes3DStiffFlux, _NavierStokes3DSetStiffFlux_). Thus, only the characteristic fields / eigen-modes corresponding to \( u\pm a\) are used. Reference:
fI | Computed upwind interface flux |
fL | Left-biased reconstructed interface flux |
fR | Right-biased reconstructed interface flux |
uL | Left-biased reconstructed interface solution |
uR | Right-biased reconstructed interface solution |
u | Cell-centered solution |
dir | Spatial dimension (x, y, or z) |
s | Solver object of type HyPar |
t | Current solution time |
Definition at line 430 of file NavierStokes3DUpwind.c.
int NavierStokes3DUpwindFdFRoe | ( | double * | fI, |
double * | fL, | ||
double * | fR, | ||
double * | uL, | ||
double * | uR, | ||
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
The Roe upwinding scheme (NavierStokes3DUpwindRoe) for the partitioned hyperbolic flux that comprises of the entropy waves only (see NavierStokes3DNonStiffFlux, _NavierStokes3DSetStiffFlux_). Thus, only the characteristic fields / eigen-modes corresponding to \(u\) are used. Reference:
fI | Computed upwind interface flux |
fL | Left-biased reconstructed interface flux |
fR | Right-biased reconstructed interface flux |
uL | Left-biased reconstructed interface solution |
uR | Right-biased reconstructed interface solution |
u | Cell-centered solution |
dir | Spatial dimension (x, y, or z) |
s | Solver object of type HyPar |
t | Current solution time |
Definition at line 527 of file NavierStokes3DUpwind.c.
int NavierStokes3DUpwindRusanovModified | ( | double * | fI, |
double * | fL, | ||
double * | fR, | ||
double * | uL, | ||
double * | uR, | ||
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
Modified Rusanov's upwinding scheme: NavierStokes3DUpwindRusanov() modified as described in the following paper (for consistent characteristic-based splitting):
fI | Computed upwind interface flux |
fL | Left-biased reconstructed interface flux |
fR | Right-biased reconstructed interface flux |
uL | Left-biased reconstructed interface solution |
uR | Right-biased reconstructed interface solution |
u | Cell-centered solution |
dir | Spatial dimension (x,y, or z) |
s | Solver object of type HyPar |
t | Current solution time |
Definition at line 638 of file NavierStokes3DUpwind.c.
int NavierStokes3DUpwinddFRusanovModified | ( | double * | fI, |
double * | fL, | ||
double * | fR, | ||
double * | uL, | ||
double * | uR, | ||
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
The modified Rusanov upwinding scheme (NavierStokes3DUpwindRusanovModified()) for the partitioned hyperbolic flux that comprises of the acoustic waves only (see NavierStokes3DStiffFlux, _NavierStokes3DSetStiffFlux_). Thus, only the characteristic fields / eigen-modes corresponding to \( u\pm a\) are used. Reference:
fI | Computed upwind interface flux |
fL | Left-biased reconstructed interface flux |
fR | Right-biased reconstructed interface flux |
uL | Left-biased reconstructed interface solution |
uR | Right-biased reconstructed interface solution |
u | Cell-centered solution |
dir | Spatial dimension (x,y, or z) |
s | Solver object of type HyPar |
t | Current solution time |
Definition at line 756 of file NavierStokes3DUpwind.c.
int NavierStokes3DUpwindFdFRusanovModified | ( | double * | fI, |
double * | fL, | ||
double * | fR, | ||
double * | uL, | ||
double * | uR, | ||
double * | u, | ||
int | dir, | ||
void * | s, | ||
double | t | ||
) |
The modified Rusanov upwinding scheme (NavierStokes3DUpwindRusanovModified()) for the partitioned hyperbolic flux that comprises of the entropy waves only (see NavierStokes3DNonStiffFlux, _NavierStokes3DSetNonStiffFlux_). Thus, only the characteristic fields / eigen-modes corresponding to \(u\) are used. Reference:
fI | Computed upwind interface flux |
fL | Left-biased reconstructed interface flux |
fR | Right-biased reconstructed interface flux |
uL | Left-biased reconstructed interface solution |
uR | Right-biased reconstructed interface solution |
u | Cell-centered solution |
dir | Spatial dimension (x,y, or z) |
s | Solver object of type HyPar |
t | Current solution time |
Definition at line 862 of file NavierStokes3DUpwind.c.
int NavierStokes3DGravityField | ( | void * | s, |
void * | m | ||
) |
This function computes the pressure and density variation functions for the hydrostatic balance of the type specified by the user (NavierStokes3D:HB). The pressure and density in hydrostatic balance are given by
\begin{equation} \rho = \rho_0\varrho\left(x,y\right),\ p = p_0\varphi\left(x,y\right) \end{equation}
where \(\rho_0\) and \(p_0\) are the reference density (NavierStokes3D::rho0) and pressure (NavierStokes3D::p0). This function computes \(\varrho\) (NavierStokes3D::grav_field_f) and \(\varphi\) (NavierStokes3D::grav_field_g). For flows without gravity, \(\varrho = \varphi = 1\).
References:
s | Solver object of type HyPar |
m | MPI object of type MPIVariables |
Definition at line 33 of file NavierStokes3DGravityField.c.
int NavierStokes3DModifiedSolution | ( | double * | uC, |
double * | u, | ||
int | d, | ||
void * | s, | ||
void * | m, | ||
double | waqt | ||
) |
This function computes the modified solution for the well-balanced treatment of the gravitational source terms. The modified solution vector is given by
\begin{equation} {\bf u}^* = \left[\begin{array}{c} \rho \varrho^{-1}\left(x,y\right) \\ \rho u \varrho^{-1}\left(x,y\right) \\ \rho v \varrho^{-1}\left(x,y\right) \\ \rho w \varrho^{-1}\left(x,y\right) \\ e^* \end{array}\right] \end{equation}
where
\begin{equation} e^* = \frac{p \varphi^{-1}\left(x,y\right)}{\gamma-1} + \frac{1}{2}\rho \varrho^{-1}\left(x,y\right) \left(u^2+v^2+w^2\right) \end{equation}
\(\varrho\) and \(\varphi\) are computed in NavierStokes3DGravityField(). For flows without gravity, \({\bf u}^* = {\bf u}\).
References:
uC | Array to hold the computed modified solution |
u | Solution vector array |
d | spatial dimension (not used) |
s | Solver object of type HyPar |
m | MPI object of time MPIVariables |
waqt | Current simulation time |
Definition at line 31 of file NavierStokes3DModifiedSolution.c.
int NavierStokes3DIBAdiabatic | ( | void * | s, |
void * | m, | ||
double * | u, | ||
double | t | ||
) |
Apply no-slip adiabatic wall boundary conditions on the immersed boundary points (grid points within the immersed body that are within stencil-width distance of interior points, i.e., points in the interior of the computational domain).
s | Solver object of type HyPar |
m | Solver object of type HyPar |
u | Array with the solution vector |
t | Current simulation time |
Definition at line 19 of file NavierStokes3DImmersedBoundary.c.
int NavierStokes3DIBIsothermal | ( | void * | s, |
void * | m, | ||
double * | u, | ||
double | t | ||
) |
Apply no-slip isothermal wall boundary conditions on the immersed boundary points (grid points within the immersed body that are within stencil-width distance of interior points, i.e., points in the interior of the computational domain).
s | Solver object of type HyPar |
m | Solver object of type HyPar |
u | Array with the solution vector |
t | Current simulation time |
Definition at line 119 of file NavierStokes3DImmersedBoundary.c.
int NavierStokes3DPreStep | ( | double * | u, |
void * | s, | ||
void * | m, | ||
double | waqt | ||
) |
Pre-step function for the 3D Navier Stokes equations: This function is called at the beginning of each time step.
\begin{equation} A_f\left({\bf u}\right) = R\left({\bf u}\right)\Lambda_f\left({\bf u}\right)L\left({\bf u}\right) \end{equation}
where \(R\) and \(L\) are the matrices of right and left eigenvectors, and,\begin{equation} \Lambda_f = diag\left[0,0,0,u+c,u-c \right] \end{equation}
with \(c\) being the speed of sound.u | Solution vector |
s | Solver object of type HyPar |
m | MPI object of type MPIVariables |
waqt | Current simulation time |
Definition at line 32 of file NavierStokes3DPreStep.c.
int NavierStokes3DIBForces | ( | void * | s, |
void * | m, | ||
double | a_t | ||
) |
Calculate the aerodynamic forces on the immersed body surface and write them to file
s | Solver object of type HyPar |
m | MPI object of type MPIVariables |
a_t | Current simulation time |
Definition at line 264 of file NavierStokes3DIBForces.c.
int gpuNavierStokes3DInitialize | ( | void * | s, |
void * | m | ||
) |
Initialize GPU-related arrays.
s | Solver object of type HyPar |
m | MPI object of type MPIVariables |
Definition at line 10 of file NavierStokes3DInitialize_GPU.c.
int gpuNavierStokes3DFlux | ( | double * | , |
double * | , | ||
int | , | ||
void * | , | ||
double | |||
) |
int gpuNavierStokes3DParabolicFunction | ( | double * | , |
double * | , | ||
void * | , | ||
void * | , | ||
double | |||
) |
int gpuNavierStokes3DSource | ( | double * | , |
double * | , | ||
void * | , | ||
void * | , | ||
double | |||
) |
int gpuNavierStokes3DModifiedSolution | ( | double * | , |
double * | , | ||
int | , | ||
void * | , | ||
void * | , | ||
double | |||
) |
int gpuNavierStokes3DUpwindRusanov | ( | double * | , |
double * | , | ||
double * | , | ||
double * | , | ||
double * | , | ||
double * | , | ||
int | , | ||
void * | , | ||
double | |||
) |
int gpuNavierStokes3DUpwindRoe | ( | double * | , |
double * | , | ||
double * | , | ||
double * | , | ||
double * | , | ||
double * | , | ||
int | , | ||
void * | , | ||
double | |||
) |
int gpuNavierStokes3DPreStep | ( | double * | , |
void * | , | ||
void * | , | ||
double | |||
) |
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.