HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
navierstokes3d.h File Reference

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
 

Detailed Description

3D Navier Stokes equations (compressible flows)

Author
Debojyoti Ghosh

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):-

  • Tannehill, Anderson and Pletcher, Computational Fluid Mechanics and Heat Transfer, Chapter 5, Section 5.1.7 (However, the non-dimensional velocity and the Reynolds number is based on speed of sound, instead of the freestream velocity).

Reference for the well-balanced treatment of gravitational source term:

  • Ghosh, D., Constantinescu, E.M., Well-Balanced Formulation of Gravitational Source Terms for Conservative Finite-Difference Atmospheric Flow Solvers, AIAA Paper 2015-2889, 7th AIAA Atmospheric and Space Environments Conference, June 22-26, 2015, Dallas, TX, http://dx.doi.org/10.2514/6.2015-2889
  • Ghosh, D., Constantinescu, E.M., A Well-Balanced, Conservative Finite-Difference Algorithm for Atmospheric Flows, AIAA Journal, 54 (4), 2016, pp. 1370-1385, http://dx.doi.org/10.2514/1.J054580.

Reference for the partitioning of the flux into its stiff (acoustic) and non-stiff (convective) components:

  • Ghosh, D., Constantinescu, E. M., Semi-Implicit Time Integration of Atmospheric Flows with Characteristic-Based Flux Partitioning, SIAM Journal on Scientific Computing, 38 (3), 2016, A1848-A1875, http://dx.doi.org/10.1137/15M1044369.

Definition in file navierstokes3d.h.


Data Structure Documentation

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

Macro Definition Documentation

#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 
)
Value:
{ \
double vsq; \
rho = u[0]; \
vx = (rho==0) ? 0 : u[stride] / rho; \
vy = (rho==0) ? 0 : u[2*stride] / rho; \
vz = (rho==0) ? 0 : u[3*stride] / rho; \
e = u[4*stride]; \
vsq = vx*vx + vy*vy + vz*vz; \
P = (e - 0.5*rho*vsq) * (gamma-1.0); \
}

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 
)
Value:
{ \
if (dir == _XDIR_) { \
f[0] = rho * vx; \
f[stride] = rho * vx * vx + P; \
f[2*stride] = rho * vx * vy; \
f[3*stride] = rho * vx * vz; \
f[4*stride] = (e + P) * vx; \
} else if (dir == _YDIR_) { \
f[0] = rho * vy; \
f[stride] = rho * vy * vx; \
f[2*stride] = rho * vy * vy + P; \
f[3*stride] = rho * vy * vz; \
f[4*stride] = (e + P) * vy; \
} else if (dir == _ZDIR_) { \
f[0] = rho * vz; \
f[stride] = rho * vz * vx; \
f[2*stride] = rho * vz * vy; \
f[3*stride] = rho * vz * vz + P; \
f[4*stride] = (e + P) * vz; \
} \
}
#define _YDIR_
#define _ZDIR_
#define _XDIR_

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 
)
Value:
{ \
double rho ,vx, vy, vz, e ,P ,H; \
double rhoL,vxL,vyL,vzL,eL,PL,HL,cLsq; \
double rhoR,vxR,vyR,vzR,eR,PR,HR,cRsq; \
_NavierStokes3DGetFlowVar_(uL,stride,rhoL,vxL,vyL,vzL,eL,PL,gamma); \
cLsq = gamma * PL/rhoL; \
HL = 0.5*(vxL*vxL+vyL*vyL+vzL*vzL) + cLsq / (gamma-1.0); \
_NavierStokes3DGetFlowVar_(uR,stride,rhoR,vxR,vyR,vzR,eR,PR,gamma); \
cRsq = gamma * PR/rhoR; \
HR = 0.5*(vxR*vxR+vyR*vyR+vzR*vzR) + cRsq / (gamma-1.0); \
double tL = sqrt(rhoL); \
double tR = sqrt(rhoR); \
rho = tL * tR; \
vx = (tL*vxL + tR*vxR) / (tL + tR); \
vy = (tL*vyL + tR*vyR) / (tL + tR); \
vz = (tL*vzL + tR*vzR) / (tL + tR); \
H = (tL*HL + tR*HR ) / (tL + tR); \
P = (H - 0.5* (vx*vx+vy*vy+vz*vz)) * (rho*(gamma-1.0))/gamma; \
e = P/(gamma-1.0) + 0.5*rho*(vx*vx+vy*vy+vz*vz); \
uavg[0] = rho; \
uavg[1] = rho*vx; \
uavg[2] = rho*vy; \
uavg[3] = rho*vz; \
uavg[4] = e; \
}
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, 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 
)
Value:
{ \
double gamma_inv = 1.0/gamma; \
if (dir == _XDIR_) { \
f[0*stride] = gamma_inv * rho * vx; \
f[1*stride] = gamma_inv * rho * vx * vx + P; \
f[2*stride] = gamma_inv * rho * vx * vy; \
f[3*stride] = gamma_inv * rho * vx * vz; \
f[4*stride] = (e + P) * vx - 0.5 * gamma_inv * (gamma-1.0) * rho * (vx*vx+vy*vy+vz*vz) * vx; \
} else if (dir == _YDIR_) { \
f[0*stride] = gamma_inv * rho * vy; \
f[1*stride] = gamma_inv * rho * vy * vx; \
f[2*stride] = gamma_inv * rho * vy * vy + P; \
f[3*stride] = gamma_inv * rho * vy * vz; \
f[4*stride] = (e + P) * vy - 0.5 * gamma_inv * (gamma-1.0) * rho * (vx*vx+vy*vy+vz*vz) * vy; \
} else if (dir == _ZDIR_) { \
f[0*stride] = gamma_inv * rho * vz; \
f[1*stride] = gamma_inv * rho * vz * vx; \
f[2*stride] = gamma_inv * rho * vz * vy; \
f[3*stride] = gamma_inv * rho * vz * vz + P; \
f[4*stride] = (e + P) * vz - 0.5 * gamma_inv * (gamma-1.0) * rho * (vx*vx+vy*vy+vz*vz) * vz; \
} \
}
#define _YDIR_
#define _ZDIR_
#define _XDIR_

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:

  • Ghosh, D., Constantinescu, E. M., Semi-Implicit Time Integration of Atmospheric Flows with Characteristic-Based Flux Partitioning, SIAM Journal on Scientific Computing, 38 (3), 2016, A1848-A1875, http://dx.doi.org/10.1137/15M1044369.

Definition at line 183 of file navierstokes3d.h.

#define _NavierStokes3DSetNonStiffFlux_ (   f,
  rho,
  vx,
  vy,
  vz,
  e,
  P,
  dir,
  gamma 
)
Value:
{ \
double gamma_inv = 1.0/gamma; \
if (dir == _XDIR_) { \
f[0*stride] = (gamma-1.0) * gamma_inv * rho * vx; \
f[1*stride] = (gamma-1.0) * gamma_inv * rho * vx * vx; \
f[2*stride] = (gamma-1.0) * gamma_inv * rho * vx * vy; \
f[3*stride] = (gamma-1.0) * gamma_inv * rho * vx * vz; \
f[4*stride] = 0.5 * gamma_inv * (gamma-1.0) * rho * (vx*vx+vy*vy+vz*vz) * vx; \
} else if (dir == _YDIR_) { \
f[0*stride] = (gamma-1.0) * gamma_inv * rho * vy; \
f[1*stride] = (gamma-1.0) * gamma_inv * rho * vy * vx; \
f[2*stride] = (gamma-1.0) * gamma_inv * rho * vy * vy; \
f[3*stride] = (gamma-1.0) * gamma_inv * rho * vy * vz; \
f[4*stride] = 0.5 * gamma_inv * (gamma-1.0) * rho * (vx*vx+vy*vy+vz*vz) * vy; \
} else if (dir == _ZDIR_) { \
f[0*stride] = (gamma-1.0) * gamma_inv * rho * vz; \
f[1*stride] = (gamma-1.0) * gamma_inv * rho * vz * vx; \
f[2*stride] = (gamma-1.0) * gamma_inv * rho * vz * vy; \
f[3*stride] = (gamma-1.0) * gamma_inv * rho * vz * vz; \
f[4*stride] = 0.5 * gamma_inv * (gamma-1.0) * rho * (vx*vx+vy*vy+vz*vz) * vz; \
} \
}
#define _YDIR_
#define _ZDIR_
#define _XDIR_

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:

  • Ghosh, D., Constantinescu, E. M., Semi-Implicit Time Integration of Atmospheric Flows with Characteristic-Based Flux Partitioning, SIAM Journal on Scientific Computing, 38 (3), 2016, A1848-A1875, http://dx.doi.org/10.1137/15M1044369.

Definition at line 219 of file navierstokes3d.h.

#define _NavierStokes3DEigenvalues_ (   u,
  stride,
  D,
  gamma,
  dir 
)
Value:
{ \
double rho,vx,vy,vz,e,P,c; \
_NavierStokes3DGetFlowVar_(u,stride,rho,vx,vy,vz,e,P,gamma); \
c = sqrt(gamma*P/rho); \
if (dir == _XDIR_) { \
D[0*_MODEL_NVARS_+0] = vx; D[0*_MODEL_NVARS_+1] = 0; D[0*_MODEL_NVARS_+2] = 0; D[0*_MODEL_NVARS_+3] = 0; D[0*_MODEL_NVARS_+4] = 0; \
D[1*_MODEL_NVARS_+0] = 0; D[1*_MODEL_NVARS_+1] = vx-c; D[1*_MODEL_NVARS_+2] = 0; D[1*_MODEL_NVARS_+3] = 0; D[1*_MODEL_NVARS_+4] = 0; \
D[2*_MODEL_NVARS_+0] = 0; D[2*_MODEL_NVARS_+1] = 0; D[2*_MODEL_NVARS_+2] = vx; D[2*_MODEL_NVARS_+3] = 0; D[2*_MODEL_NVARS_+4] = 0; \
D[3*_MODEL_NVARS_+0] = 0; D[3*_MODEL_NVARS_+1] = 0; D[3*_MODEL_NVARS_+2] = 0; D[3*_MODEL_NVARS_+3] = vx; D[3*_MODEL_NVARS_+4] = 0; \
D[4*_MODEL_NVARS_+0] = 0; D[4*_MODEL_NVARS_+1] = 0; D[4*_MODEL_NVARS_+2] = 0; D[4*_MODEL_NVARS_+3] = 0; D[4*_MODEL_NVARS_+4] = vx+c;\
} else if (dir == _YDIR_) { \
D[0*_MODEL_NVARS_+0] = vy; D[0*_MODEL_NVARS_+1] = 0; D[0*_MODEL_NVARS_+2] = 0; D[0*_MODEL_NVARS_+3] = 0; D[0*_MODEL_NVARS_+4] = 0; \
D[1*_MODEL_NVARS_+0] = 0; D[1*_MODEL_NVARS_+1] = vy; D[1*_MODEL_NVARS_+2] = 0; D[1*_MODEL_NVARS_+3] = 0; D[1*_MODEL_NVARS_+4] = 0; \
D[2*_MODEL_NVARS_+0] = 0; D[2*_MODEL_NVARS_+1] = 0; D[2*_MODEL_NVARS_+2] = vy-c; D[2*_MODEL_NVARS_+3] = 0; D[2*_MODEL_NVARS_+4] = 0; \
D[3*_MODEL_NVARS_+0] = 0; D[3*_MODEL_NVARS_+1] = 0; D[3*_MODEL_NVARS_+2] = 0; D[3*_MODEL_NVARS_+3] = vy; D[3*_MODEL_NVARS_+4] = 0; \
D[4*_MODEL_NVARS_+0] = 0; D[4*_MODEL_NVARS_+1] = 0; D[4*_MODEL_NVARS_+2] = 0; D[4*_MODEL_NVARS_+3] = 0; D[4*_MODEL_NVARS_+4] = vy+c;\
} else if (dir == _ZDIR_) { \
D[0*_MODEL_NVARS_+0] = vz; D[0*_MODEL_NVARS_+1] = 0; D[0*_MODEL_NVARS_+2] = 0; D[0*_MODEL_NVARS_+3] = 0; D[0*_MODEL_NVARS_+4] = 0; \
D[1*_MODEL_NVARS_+0] = 0; D[1*_MODEL_NVARS_+1] = vz; D[1*_MODEL_NVARS_+2] = 0; D[1*_MODEL_NVARS_+3] = 0; D[1*_MODEL_NVARS_+4] = 0; \
D[2*_MODEL_NVARS_+0] = 0; D[2*_MODEL_NVARS_+1] = 0; D[2*_MODEL_NVARS_+2] = vz; D[2*_MODEL_NVARS_+3] = 0; D[2*_MODEL_NVARS_+4] = 0; \
D[3*_MODEL_NVARS_+0] = 0; D[3*_MODEL_NVARS_+1] = 0; D[3*_MODEL_NVARS_+2] = 0; D[3*_MODEL_NVARS_+3] = vz-c; D[3*_MODEL_NVARS_+4] = 0; \
D[4*_MODEL_NVARS_+0] = 0; D[4*_MODEL_NVARS_+1] = 0; D[4*_MODEL_NVARS_+2] = 0; D[4*_MODEL_NVARS_+3] = 0; D[4*_MODEL_NVARS_+4] = vz+c;\
} \
}
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)
#define _YDIR_
#define _ZDIR_
#define _XDIR_
#define _MODEL_NVARS_

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.

Function Documentation

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_)

Note: "physics.inp" is optional; if absent, default values will be used.

Parameters
sSolver object of type HyPar
mMPI object of type MPIVariables

Definition at line 116 of file NavierStokes3DInitialize.c.

119 {
120  HyPar *solver = (HyPar*) s;
121  MPIVariables *mpi = (MPIVariables*) m;
122  NavierStokes3D *physics = (NavierStokes3D*) solver->physics;
123  int ferr = 0;
124 
125  static int count = 0;
126 
127  if (solver->nvars != _MODEL_NVARS_) {
128  fprintf(stderr,"Error in NavierStokes3DInitialize(): nvars has to be %d.\n",_MODEL_NVARS_);
129  return(1);
130  }
131  if (solver->ndims != _MODEL_NDIMS_) {
132  fprintf(stderr,"Error in NavierStokes3DInitialize(): ndims has to be %d.\n",_MODEL_NDIMS_);
133  return(1);
134  }
135 
136  /* default values */
137  physics->gamma = 1.4;
138  physics->Pr = 0.72;
139  physics->Re = -1;
140  physics->Minf = 1.0;
141  physics->C1 = 1.458e-6;
142  physics->C2 = 110.4;
143  physics->grav_x = 0.0;
144  physics->grav_y = 0.0;
145  physics->grav_z = 0.0;
146  physics->rho0 = 1.0;
147  physics->p0 = 1.0;
148  physics->HB = 1;
149  physics->R = 1.0;
150  physics->N_bv = 0.0;
151  physics->T_ib_wall = -DBL_MAX;
152  physics->t_ib_ramp = -1.0;
153  physics->t_ib_width= 0.05;
154  physics->ib_T_tol = 5;
155  strcpy(physics->upw_choice,"roe");
156  strcpy(physics->ib_write_surface_data,"yes");
157  strcpy(physics->ib_wall_type,"adiabatic");
158  strcpy(physics->ib_ramp_type,"linear");
159 
160  /* reading physical model specific inputs - all processes */
161  if (!mpi->rank) {
162  FILE *in;
163  if (!count) printf("Reading physical model inputs from file \"physics.inp\".\n");
164  in = fopen("physics.inp","r");
165  if (!in) printf("Warning: File \"physics.inp\" not found. Using default values.\n");
166  else {
167  char word[_MAX_STRING_SIZE_];
168  ferr = fscanf(in,"%s",word);
169  if (ferr != 1) {
170  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
171  return 1;
172  }
173  if (!strcmp(word, "begin")){
174  while (strcmp(word, "end")){
175  ferr = fscanf(in,"%s",word);
176  if (ferr != 1) {
177  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
178  return 1;
179  }
180  if (!strcmp(word, "gamma")) {
181  ferr = fscanf(in,"%lf",&physics->gamma);
182  if (ferr != 1) {
183  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
184  return 1;
185  }
186  } else if (!strcmp(word,"upwinding")) {
187  ferr = fscanf(in,"%s",physics->upw_choice);
188  if (ferr != 1) {
189  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
190  return 1;
191  }
192  } else if (!strcmp(word,"Pr")) {
193  ferr = fscanf(in,"%lf",&physics->Pr);
194  if (ferr != 1) {
195  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
196  return 1;
197  }
198  } else if (!strcmp(word,"Re")) {
199  ferr = fscanf(in,"%lf",&physics->Re);
200  if (ferr != 1) {
201  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
202  return 1;
203  }
204  } else if (!strcmp(word,"Minf")) {
205  ferr = fscanf(in,"%lf",&physics->Minf);
206  if (ferr != 1) {
207  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
208  return 1;
209  }
210  } else if (!strcmp(word,"gravity")) {
211  ferr = fscanf(in,"%lf",&physics->grav_x);
212  if (ferr != 1) {
213  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
214  return 1;
215  }
216  ferr = fscanf(in,"%lf",&physics->grav_y);
217  if (ferr != 1) {
218  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
219  return 1;
220  }
221  ferr = fscanf(in,"%lf",&physics->grav_z);
222  if (ferr != 1) {
223  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
224  return 1;
225  }
226  } else if (!strcmp(word,"rho_ref")) {
227  ferr = fscanf(in,"%lf",&physics->rho0);
228  if (ferr != 1) {
229  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
230  return 1;
231  }
232  } else if (!strcmp(word,"p_ref")) {
233  ferr = fscanf(in,"%lf",&physics->p0);
234  if (ferr != 1) {
235  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
236  return 1;
237  }
238  } else if (!strcmp(word,"HB")) {
239  ferr = fscanf(in,"%d",&physics->HB);
240  if (ferr != 1) {
241  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
242  return 1;
243  }
244  if (physics->HB==3) {
245  ferr = fscanf(in,"%lf",&physics->N_bv);
246  if (ferr != 1) {
247  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
248  return 1;
249  }
250  }
251  } else if (!strcmp(word,"R")) {
252  ferr = fscanf(in,"%lf",&physics->R);
253  if (ferr != 1) {
254  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
255  return 1;
256  }
257  } else if (!strcmp(word,"ib_surface_data")) {
258  ferr = fscanf(in,"%s",physics->ib_write_surface_data);
259  if (ferr != 1) {
260  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
261  return 1;
262  }
263  } else if (!strcmp(word,"ib_wall_type")) {
264  ferr = fscanf(in,"%s",physics->ib_wall_type);
265  if (ferr != 1) {
266  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
267  return 1;
268  }
269  if (!strcmp(physics->ib_wall_type,_IB_ISOTHERMAL_)) {
270  ferr = fscanf(in,"%lf",&physics->T_ib_wall);
271  if (ferr != 1) {
272  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
273  return 1;
274  }
275  }
276  if (!solver->flag_ib) {
277  printf("Warning: in NavierStokes3DInitialize().\n");
278  printf("Warning: no immersed body present; specification of ib_wall_type unnecessary.\n");
279  }
280  } else if (!strcmp(word,"ib_ramp_time")) {
281  ferr = fscanf(in,"%lf",&physics->t_ib_ramp);
282  if (ferr != 1) {
283  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
284  return 1;
285  }
286  if (!solver->flag_ib) {
287  printf("Warning: in NavierStokes3DInitialize().\n");
288  printf("Warning: no immersed body present; specification of ib_ramp_time unnecessary.\n");
289  }
290  } else if (!strcmp(word,"ib_ramp_width")) {
291  ferr = fscanf(in,"%lf",&physics->t_ib_width);
292  if (ferr != 1) {
293  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
294  return 1;
295  }
296  if (!solver->flag_ib) {
297  printf("Warning: in NavierStokes3DInitialize().\n");
298  printf("Warning: no immersed body present; specification of ib_ramp_width unnecessary.\n");
299  }
300  } else if (!strcmp(word,"ib_ramp_type")) {
301  ferr = fscanf(in,"%s",physics->ib_ramp_type);
302  if (ferr != 1) {
303  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
304  return 1;
305  }
306  if (!solver->flag_ib) {
307  printf("Warning: in NavierStokes3DInitialize().\n");
308  printf("Warning: no immersed body present; specification of ib_ramp_type unnecessary.\n");
309  }
310  } else if (!strcmp(word,"ib_T_tolerance")) {
311  ferr = fscanf(in,"%lf",&physics->ib_T_tol);
312  if (ferr != 1) {
313  fprintf(stderr, "Read error while reading physics.inp in NavierStokes3DInitialize().\n");
314  return 1;
315  }
316  if (!solver->flag_ib) {
317  printf("Warning: in NavierStokes3DInitialize().\n");
318  printf("Warning: no immersed body present; specification of ib_T_tolerance unnecessary.\n");
319  }
320  } else if (strcmp(word,"end")) {
321  char useless[_MAX_STRING_SIZE_];
322  ferr = fscanf(in,"%s",useless); if (ferr != 1) return(ferr);
323  printf("Warning: keyword %s in file \"physics.inp\" with value %s not ",word,useless);
324  printf("recognized or extraneous. Ignoring.\n");
325  }
326  }
327  } else {
328  fprintf(stderr,"Error: Illegal format in file \"physics.inp\".\n");
329  return(1);
330  }
331  }
332  fclose(in);
333  }
334 
339  IERR MPIBroadcast_double (&physics->gamma ,1 ,0,&mpi->world); CHECKERR(ierr);
340  IERR MPIBroadcast_double (&physics->Pr ,1 ,0,&mpi->world); CHECKERR(ierr);
341  IERR MPIBroadcast_double (&physics->Re ,1 ,0,&mpi->world); CHECKERR(ierr);
342  IERR MPIBroadcast_double (&physics->Minf ,1 ,0,&mpi->world); CHECKERR(ierr);
343  IERR MPIBroadcast_double (&physics->grav_x ,1 ,0,&mpi->world); CHECKERR(ierr);
344  IERR MPIBroadcast_double (&physics->grav_y ,1 ,0,&mpi->world); CHECKERR(ierr);
345  IERR MPIBroadcast_double (&physics->grav_z ,1 ,0,&mpi->world); CHECKERR(ierr);
346  IERR MPIBroadcast_double (&physics->rho0 ,1 ,0,&mpi->world); CHECKERR(ierr);
347  IERR MPIBroadcast_double (&physics->p0 ,1 ,0,&mpi->world); CHECKERR(ierr);
348  IERR MPIBroadcast_double (&physics->R ,1 ,0,&mpi->world); CHECKERR(ierr);
349  IERR MPIBroadcast_double (&physics->N_bv ,1 ,0,&mpi->world); CHECKERR(ierr);
350  IERR MPIBroadcast_double (&physics->T_ib_wall ,1 ,0,&mpi->world); CHECKERR(ierr);
351  IERR MPIBroadcast_double (&physics->t_ib_ramp ,1 ,0,&mpi->world); CHECKERR(ierr);
352  IERR MPIBroadcast_double (&physics->t_ib_width ,1 ,0,&mpi->world); CHECKERR(ierr);
353  IERR MPIBroadcast_double (&physics->ib_T_tol ,1 ,0,&mpi->world); CHECKERR(ierr);
354  IERR MPIBroadcast_integer (&physics->HB ,1 ,0,&mpi->world); CHECKERR(ierr);
355 
356  /* if file output is disabled in HyPar, respect that */
357  if (!strcmp(solver->op_file_format,"none")) {
358  if (!strcmp(physics->ib_write_surface_data,"yes")) {
359  if (!mpi->rank) {
360  printf("Warning from NavierStokes3DInitialize(): solver->op_file_format is set to \"none\", thus ");
361  printf("setting physics->ib_write_surface_data to \"no\" (no solution files will be written).\n");
362  }
363  }
364  strcpy(physics->ib_write_surface_data,"no");
365  }
366 
367  /* Scaling Re by M_inf */
368  physics->Re /= physics->Minf;
369 
370  /* check that a well-balanced upwinding scheme is being used for cases with gravity */
371  if ( ((physics->grav_x != 0.0) || (physics->grav_y != 0.0) || (physics->grav_z != 0.0) )
372  && (strcmp(physics->upw_choice,_RUSANOV_)) ) {
373  if (!mpi->rank) {
374  fprintf(stderr,"Error in NavierStokes3DInitialize: %s upwinding is needed for flows ",_RUSANOV_);
375  fprintf(stderr,"with gravitational forces.\n");
376  }
377  return(1);
378  }
379  /* check that solver has the correct choice of diffusion formulation, if viscous flow */
380  if (strcmp(solver->spatial_type_par,_NC_2STAGE_) && (physics->Re > 0)) {
381  if (!mpi->rank) {
382  fprintf(stderr,"Error in NavierStokes3DInitialize(): Parabolic term spatial discretization must be \"%s\"\n",_NC_2STAGE_);
383  }
384  return(1);
385  }
386 
387  /* initializing physical model-specific functions */
388 #if defined(HAVE_CUDA)
389  if (solver->use_gpu) {
394  } else {
395 #endif
396  solver->PreStep = NavierStokes3DPreStep;
398  solver->FFunction = NavierStokes3DFlux;
404 #if defined(HAVE_CUDA)
405  }
406 #endif
407 
408  if (solver->flag_ib) {
409  if (!strcmp(physics->ib_wall_type,_IB_ADIABATIC_)) {
411  } else if (!strcmp(physics->ib_wall_type,_IB_ISOTHERMAL_)) {
413  } else {
414  fprintf(stderr, "Error in NavierStokes3DInitialize()\n");
415  fprintf(stderr, " invalid value for IB wall type (%s).\n",
416  physics->ib_wall_type );
417  }
418  if (!strcmp(physics->ib_write_surface_data,"yes")) {
420  }
421  }
422 
423 #if defined(HAVE_CUDA)
424  if (solver->use_gpu) {
425 
426  if (!strcmp(solver->SplitHyperbolicFlux,"yes")) {
427  if (!mpi->rank) {
428  fprintf(stderr,"Error in NavierStokes3DInitialize(): Not available on GPU!");
429  }
430  return 1;
431  } else {
433  if (!strcmp(physics->upw_choice,_ROE_ )) {
435  } else if (!strcmp(physics->upw_choice,_RUSANOV_)) {
437  } else {
438  if (!mpi->rank) {
439  fprintf(stderr,"Error in NavierStokes3DInitialize(): %s is not a valid upwinding scheme on GPU. ",
440  physics->upw_choice);
441  fprintf(stderr,"Choices are %s and %s.\n",_ROE_,_RUSANOV_);
442  }
443  return(1);
444  }
445  }
446 
447  } else {
448 #endif
449 
450  if (!strcmp(solver->SplitHyperbolicFlux,"yes")) {
454  if (!strcmp(physics->upw_choice,_ROE_)) {
458  } else if (!strcmp(physics->upw_choice,_RUSANOV_)) {
462  } else {
463  if (!mpi->rank) {
464  fprintf(stderr,"Error in NavierStokes3DInitialize(): %s is not a valid upwinding scheme ",
465  physics->upw_choice);
466  fprintf(stderr,"for use with split hyperbolic flux form. Use %s or %s.\n",
467  _ROE_,_RUSANOV_);
468  }
469  return(1);
470  }
471  } else {
473  if (!strcmp(physics->upw_choice,_ROE_ )) solver->Upwind = NavierStokes3DUpwindRoe;
474  else if (!strcmp(physics->upw_choice,_RF_ )) solver->Upwind = NavierStokes3DUpwindRF;
475  else if (!strcmp(physics->upw_choice,_LLF_ )) solver->Upwind = NavierStokes3DUpwindLLF;
476  else if (!strcmp(physics->upw_choice,_RUSANOV_)) solver->Upwind = NavierStokes3DUpwindRusanov;
477  else {
478  if (!mpi->rank) {
479  fprintf(stderr,"Error in NavierStokes3DInitialize(): %s is not a valid upwinding scheme. ",
480  physics->upw_choice);
481  fprintf(stderr,"Choices are %s, %s, %s, and %s.\n",_ROE_,_RF_,_LLF_,_RUSANOV_);
482  }
483  return(1);
484  }
485  }
486 
487 #if defined(HAVE_CUDA)
488  }
489 #endif
490 
491  /* set the value of gamma in all the boundary objects */
492  int n;
493  DomainBoundary *boundary = (DomainBoundary*) solver->boundary;
494  for (n = 0; n < solver->nBoundaryZones; n++) boundary[n].gamma = physics->gamma;
495 
496  /* finally, hijack the main solver's dissipation function pointer
497  * to this model's own function, since it's difficult to express
498  * the dissipation terms in the general form */
499 #if defined(HAVE_CUDA)
500  if (solver->use_gpu) {
502  } else {
503 #endif
505 #if defined(HAVE_CUDA)
506  }
507 #endif
508 
509  /* allocate array to hold the gravity field */
510  int *dim = solver->dim_local;
511  int ghosts = solver->ghosts;
512  int d, size = 1; for (d=0; d<_MODEL_NDIMS_; d++) size *= (dim[d] + 2*ghosts);
513  physics->grav_field_f = (double*) calloc (size, sizeof(double));
514  physics->grav_field_g = (double*) calloc (size, sizeof(double));
515  /* allocate arrays to hold the fast Jacobian for split form of the hyperbolic flux */
516  physics->fast_jac = (double*) calloc (_MODEL_NDIMS_*size*_MODEL_NVARS_*_MODEL_NVARS_,sizeof(double));
517  physics->solution = (double*) calloc (size*_MODEL_NVARS_,sizeof(double));
518  /* initialize the gravity fields */
519  IERR NavierStokes3DGravityField(solver,mpi); CHECKERR(ierr);
520 
521 #if defined(HAVE_CUDA)
522  if (solver->use_gpu) gpuNavierStokes3DInitialize(s,m);
523 #endif
524 
525  count++;
526  return(0);
527 }
int NavierStokes3DJacobian(double *, double *, void *, int, int, int)
int NavierStokes3DUpwindRoe(double *, double *, double *, double *, double *, double *, int, void *, double)
int(* JFunction)(double *, double *, void *, int, int, int)
Definition: hypar.h:326
int(* dFFunction)(double *, double *, int, void *, double)
Definition: hypar.h:280
double(* ComputeCFL)(void *, void *, double, double)
Definition: hypar.h:269
int NavierStokes3DStiffJacobian(double *, double *, void *, int, int, int)
int MPIBroadcast_double(double *, int, int, void *)
Definition: MPIBroadcast.c:9
char ib_ramp_type[_MAX_STRING_SIZE_]
#define _MODEL_NVARS_
Definition: euler1d.h:58
int NavierStokes3DUpwindRusanov(double *, double *, double *, double *, double *, double *, int, void *, double)
char ib_wall_type[_MAX_STRING_SIZE_]
char spatial_type_par[_MAX_STRING_SIZE_]
Definition: hypar.h:96
void * physics
Definition: hypar.h:266
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int(* GetRightEigenvectors)(double *, double *, void *, int)
Definition: hypar.h:359
#define _LLF_
Definition: euler1d.h:66
int NavierStokes3DSource(double *, double *, void *, void *, double)
int flag_ib
Definition: hypar.h:441
int(* UpwinddF)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:300
int gpuNavierStokes3DParabolicFunction(double *, double *, void *, void *, double)
int * dim_local
Definition: hypar.h:37
int NavierStokes3DNonStiffFlux(double *f, double *u, int dir, void *s, double t)
int(* UpwindFdF)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:307
char op_file_format[_MAX_STRING_SIZE_]
Definition: hypar.h:186
int NavierStokes3DUpwindFdFRusanovModified(double *, double *, double *, double *, double *, double *, int, void *, double)
int(* PhysicsOutput)(void *, void *, double)
Definition: hypar.h:347
int(* FdFFunction)(double *, double *, int, void *, double)
Definition: hypar.h:286
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
int gpuNavierStokes3DInitialize(void *, void *)
#define _IB_ISOTHERMAL_
#define _MODEL_NDIMS_
Definition: euler1d.h:56
double NavierStokes3DComputeCFL(void *s, void *m, double dt, double t)
#define _IB_ADIABATIC_
int(* ParabolicFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:256
int ghosts
Definition: hypar.h:52
int gpuNavierStokes3DUpwindRoe(double *, double *, double *, double *, double *, double *, int, void *, double)
int NavierStokes3DUpwinddFRoe(double *, double *, double *, double *, double *, double *, int, void *, double)
int NavierStokes3DRoeAverage(double *uavg, double *uL, double *uR, void *p)
int NavierStokes3DGravityField(void *s, void *m)
MPI_Comm world
#define _NC_2STAGE_
Definition: hypar.h:477
#define _MAX_STRING_SIZE_
Definition: basic.h:14
int gpuNavierStokes3DSource(double *, double *, void *, void *, double)
int NavierStokes3DFlux(double *f, double *u, int dir, void *s, double t)
int(* AveragingFunction)(double *, double *, double *, void *)
Definition: hypar.h:354
int NavierStokes3DUpwindRusanovModified(double *, double *, double *, double *, double *, double *, int, void *, double)
double * grav_field_f
void * boundary
Definition: hypar.h:159
int gpuNavierStokes3DFlux(double *__restrict__ f, double *__restrict__ u, int dir, void *__restrict__ s, double t)
Structure containing the variables and function pointers defining a boundary.
int NavierStokes3DModifiedSolution(double *, double *, int, void *, void *, double)
#define _RUSANOV_
Definition: euler1d.h:70
int NavierStokes3DPreStep(double *, void *, void *, double)
char ib_write_surface_data[_MAX_STRING_SIZE_]
int gpuNavierStokes3DPreStep(double *, void *, void *, double)
char upw_choice[_MAX_STRING_SIZE_]
int(* GetLeftEigenvectors)(double *, double *, void *, int)
Definition: hypar.h:357
int NavierStokes3DRightEigenvectors(double *u, double *R, void *p, int dir)
#define _ROE_
Definition: euler1d.h:62
int nBoundaryZones
Definition: hypar.h:157
int(* SFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:317
int gpuNavierStokes3DModifiedSolution(double *, double *, int, void *, void *, double)
int NavierStokes3DUpwindLLF(double *, double *, double *, double *, double *, double *, int, void *, double)
int nvars
Definition: hypar.h:29
#define CHECKERR(ierr)
Definition: basic.h:18
int(* UFunction)(double *, double *, int, void *, void *, double)
Definition: hypar.h:321
int NavierStokes3DIBIsothermal(void *s, void *m, double *u, double t)
int use_gpu
Definition: hypar.h:449
int MPIBroadcast_character(char *, int, int, void *)
Definition: MPIBroadcast.c:37
int NavierStokes3DUpwindRF(double *, double *, double *, double *, double *, double *, int, void *, double)
int NavierStokes3DIBAdiabatic(void *s, void *m, double *u, double t)
double * grav_field_g
int NavierStokes3DLeftEigenvectors(double *u, double *L, void *p, int dir)
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.
char SplitHyperbolicFlux[_MAX_STRING_SIZE_]
Definition: hypar.h:92
int ndims
Definition: hypar.h:26
#define _RF_
Definition: euler1d.h:64
int gpuNavierStokes3DUpwindRusanov(double *, double *, double *, double *, double *, double *, int, void *, double)
int NavierStokes3DParabolicFunction(double *, double *, void *, void *, double)
#define IERR
Definition: basic.h:16
int(* IBFunction)(void *, void *, double *, double)
Definition: hypar.h:446
int NavierStokes3DUpwinddFRusanovModified(double *, double *, double *, double *, double *, double *, int, void *, double)
int(* PreStep)(double *, void *, void *, double)
Definition: hypar.h:339
int NavierStokes3DUpwindFdFRoe(double *, double *, double *, double *, double *, double *, int, void *, double)
Structure of MPI-related variables.
int NavierStokes3DStiffFlux(double *f, double *u, int dir, void *s, double t)
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
int NavierStokes3DIBForces(void *s, void *m, double a_t)
int MPIBroadcast_integer(int *, int, int, void *)
Definition: MPIBroadcast.c:23
int NavierStokes3DCleanup ( void *  s)

Function to clean up all allocations in the 3D Navier Stokes module.

Parameters
sObject of type NavierStokes3D

Definition at line 11 of file NavierStokes3DCleanup.c.

12 {
13  NavierStokes3D *param = (NavierStokes3D*) s;
14 
15  free(param->grav_field_f);
16  free(param->grav_field_g);
17  free(param->fast_jac);
18  free(param->solution);
19  return(0);
20 }
double * grav_field_f
double * grav_field_g
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.
int gpuNavierStokes3DCleanup ( void *  s)

Function to clean up all allocations in the 3D Navier Stokes module.

Parameters
sObject of type NavierStokes3D

Definition at line 13 of file NavierStokes3DCleanup_GPU.c.

14 {
15  NavierStokes3D *param = (NavierStokes3D*) s;
16 
17  gpuFree(param->gpu_Q);
18  gpuFree(param->gpu_QDerivX);
19  gpuFree(param->gpu_QDerivY);
20  gpuFree(param->gpu_QDerivZ);
21  gpuFree(param->gpu_FViscous);
22  gpuFree(param->gpu_FDeriv);
23  gpuFree(param->gpu_grav_field_f);
24  gpuFree(param->gpu_grav_field_g);
25  gpuFree(param->gpu_fast_jac);
26  gpuFree(param->gpu_solution);
27 
28  return(0);
29 }
double * gpu_solution
double * gpu_FDeriv
double * gpu_QDerivY
void gpuFree(void *)
double * gpu_grav_field_f
double * gpu_QDerivX
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.
double * gpu_fast_jac
double * gpu_QDerivZ
double * gpu_grav_field_g
double * gpu_FViscous

Variable Documentation

const int _NavierStokes3D_stride_ = 1
static

Definition at line 559 of file navierstokes3d.h.