1 /*! @file NavierStokes2DFlux_GPU.cu
3 @brief Functions to compute the hyperbolic flux for 2D Navier-Stokes equations
8 #include <arrayfunctions_gpu.h>
9 #include <mathfunctions.h>
10 #include <matmult_native.h>
11 #include <physicalmodels/navierstokes2d.h>
14 #ifdef CUDA_VAR_ORDERDING_AOS
16 /*! Kernel for gpuNavierStokes2DFlux() */
18 void NavierStokes2DFlux_kernel(
26 int p = threadIdx.x + (blockDim.x * blockIdx.x);
27 if (p < ngrid_points) {
28 double rho, vx, vy, e, P;
29 _NavierStokes2DGetFlowVar_((u+_MODEL_NVARS_*p),rho,vx,vy,e,P,gamma);
30 _NavierStokes2DSetFlux_((f+_MODEL_NVARS_*p),rho,vx,vy,e,P,dir);
36 Compute the hyperbolic flux function for the 2D Navier-Stokes equations:
38 dir = x, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \rho u \\ \rho u^2 + p \\ \rho u v \\ (e+p)u \end{array}\right], \\
39 dir = y, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \rho v \\ \rho u v \\ \rho v^2 + p \\ (e+p)v \end{array}\right]
41 Note: the flux function needs to be computed at the ghost points as well.
43 extern "C" int gpuNavierStokes2DFlux(
44 double *f, /*!< Array to hold the computed flux vector (same layout as u) */
45 double *u, /*!< Array with the solution vector */
46 int dir, /*!< Spatial dimension (x or y) for which to compute the flux */
47 void *s, /*!< Solver object of type #HyPar */
48 double t /*!< Current simulation time */
51 HyPar *solver = (HyPar*) s;
52 NavierStokes2D *param = (NavierStokes2D*) solver->physics;
53 int *dim = solver->dim_local;
54 int ghosts = solver->ghosts;
56 int ngrid_points = 1; for (int i=0; i<_MODEL_NDIMS_; i++) ngrid_points *= (dim[i]+2*ghosts);
57 int nblocks = (ngrid_points - 1) / GPU_THREADS_PER_BLOCK + 1;
59 clock_t cpu_start, cpu_end;
60 double cpu_time = 0.0;
63 NavierStokes2DFlux_kernel<<<nblocks, GPU_THREADS_PER_BLOCK>>>(
64 ngrid_points,dir,param->gamma,u,f
66 cudaDeviceSynchronize();
68 cpu_time += (double)(cpu_end - cpu_start) / CLOCKS_PER_SEC;