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
Go to the documentation of this file.
1 
6 #include <basic_gpu.h>
7 #include <arrayfunctions_gpu.h>
8 #include <mathfunctions.h>
9 #include <matmult_native.h>
11 #include <hypar.h>
12 
13 #ifdef CUDA_VAR_ORDERDING_AOS
14 
16 __global__
18  int npoints_grid,
19  int dir,
20  double gamma,
21  const double * __restrict__ u,
22  double * __restrict__ f
23 )
24 {
25  int p = blockDim.x * blockIdx.x + threadIdx.x;
26 
27  if (p < npoints_grid) {
28  double rho, vx, vy, vz, e, P;
31  }
32  return;
33 }
34 
44 extern "C" int gpuNavierStokes3DFlux(
45  double *__restrict__ f,
46  double *__restrict__ u,
47  int dir,
48  void *__restrict__ s,
49  double t
50 )
51 {
52  HyPar *solver = (HyPar*) s;
53  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
54 
55  int nblocks = (solver->npoints_local_wghosts-1)/GPU_THREADS_PER_BLOCK + 1;
56 
57  gpuNavierStokes3DFlux_kernel<<<nblocks, GPU_THREADS_PER_BLOCK>>>(
58  solver->npoints_local_wghosts, dir, param->gamma, u, f
59  );
60  cudaDeviceSynchronize();
61 
62  return 0;
63 }
64 
65 #else
66 
68 __global__
70  int npoints_grid,
71  int dir,
72  double gamma,
73  const double * __restrict__ u,
74  double * __restrict__ f
75 )
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 }
86 
96 extern "C" int gpuNavierStokes3DFlux(
97  double *__restrict__ f,
98  double *__restrict__ u,
99  int dir,
100  void *__restrict__ s,
101  double t
102 )
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 }
116 
117 #endif
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)
int npoints_local_wghosts
Definition: hypar.h:42
3D Navier Stokes equations (compressible flows)
#define _MODEL_NVARS_
Definition: euler1d.h:58
void * physics
Definition: hypar.h:266
#define _NavierStokes3DSetFlux_(f, stride, rho, vx, vy, vz, e, P, dir)
#define GPU_THREADS_PER_BLOCK
Definition: basic_gpu.h:7
int gpuNavierStokes3DFlux(double *__restrict__ f, double *__restrict__ u, int dir, void *__restrict__ s, double t)
Contains function definitions for common mathematical functions.
Contains function definitions for common array operations on GPU.
Contains structure definition for hypar.
__global__ void gpuNavierStokes3DFlux_kernel(int npoints_grid, int dir, double gamma, const double *__restrict__ u, double *__restrict__ f)
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
Contains macros and function definitions for common matrix multiplication.
static const int _NavierStokes3D_stride_