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

Functions to compute the hyperbolic flux for 3D Navier-Stokes equations. More...

#include <basic_gpu.h>
#include <arrayfunctions_gpu.h>
#include <mathfunctions.h>
#include <matmult_native.h>
#include <physicalmodels/navierstokes3d.h>
#include <hypar.h>

Go to the source code of this file.

Functions

__global__ void gpuNavierStokes3DFlux_kernel (int npoints_grid, int dir, double gamma, const double *__restrict__ u, double *__restrict__ f)
 
int gpuNavierStokes3DFlux (double *__restrict__ f, double *__restrict__ u, int dir, void *__restrict__ s, double t)
 

Detailed Description

Functions to compute the hyperbolic flux for 3D Navier-Stokes equations.

Author
Youngdae Kim

Definition in file NavierStokes3DFlux_GPU.cu.

Function Documentation

__global__ void gpuNavierStokes3DFlux_kernel ( int  npoints_grid,
int  dir,
double  gamma,
const double *__restrict__  u,
double *__restrict__  f 
)

Kernel for gpuNavierStokes3DFlux()

Definition at line 69 of file NavierStokes3DFlux_GPU.cu.

76 {
77  int p = blockDim.x * blockIdx.x + threadIdx.x;
78 
79  if (p < npoints_grid) {
80  double rho, vx, vy, vz, e, P;
81  _NavierStokes3DGetFlowVar_((u+p),npoints_grid,rho,vx,vy,vz,e,P,gamma);
82  _NavierStokes3DSetFlux_((f+p),npoints_grid,rho,vx,vy,vz,e,P,dir);
83  }
84  return;
85 }
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)
#define _NavierStokes3DSetFlux_(f, stride, rho, vx, vy, vz, e, P, dir)
int gpuNavierStokes3DFlux ( double *__restrict__  f,
double *__restrict__  u,
int  dir,
void *__restrict__  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.

Parameters
fArray to hold the computed flux vector (same layout as u)
uArray with the solution vector
dirSpatial dimension (x, y, or z) for which to compute the flux
sSolver object of type HyPar
tCurrent simulation time

Definition at line 96 of file NavierStokes3DFlux_GPU.cu.

103 {
104  HyPar *solver = (HyPar*) s;
105  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
106 
107  int nblocks = (solver->npoints_local_wghosts-1)/GPU_THREADS_PER_BLOCK + 1;
108 
109  gpuNavierStokes3DFlux_kernel<<<nblocks, GPU_THREADS_PER_BLOCK>>>(
110  solver->npoints_local_wghosts, dir, param->gamma, u, f
111  );
112  cudaDeviceSynchronize();
113 
114  return 0;
115 }
int npoints_local_wghosts
Definition: hypar.h:42
void * physics
Definition: hypar.h:266
#define GPU_THREADS_PER_BLOCK
Definition: basic_gpu.h:7
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.
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23