1 /*! @file NavierStokes2DModifiedSolution_GPU.cu
3 @brief Compute the modified solution for the well-balanced treatment of gravitational source terms.
6 #include <arrayfunctions.h>
7 #include <physicalmodels/navierstokes2d.h>
11 #ifdef CUDA_VAR_ORDERDING_AOS
13 /*! Kernel for gpuNavierStokes2DModifiedSolution() */
15 void NavierStokes2DModifiedSolution_kernel(
19 const double *grav_field_f,
20 const double *grav_field_g,
25 int p = threadIdx.x + (blockDim.x * blockIdx.x);
28 double rho, uvel, vvel, E, P;
29 _NavierStokes2DGetFlowVar_((u+_MODEL_NVARS_*p),rho,uvel,vvel,E,P,gamma);
30 uC[_MODEL_NVARS_*p+0] = u[_MODEL_NVARS_*p+0] * grav_field_f[p];
31 uC[_MODEL_NVARS_*p+1] = u[_MODEL_NVARS_*p+1] * grav_field_f[p];
32 uC[_MODEL_NVARS_*p+2] = u[_MODEL_NVARS_*p+2] * grav_field_f[p];
33 uC[_MODEL_NVARS_*p+3] = (P*inv_gamma_m1)*(1.0/grav_field_g[p]) + (0.5*rho*(uvel*uvel+vvel*vvel))*grav_field_f[p];
40 This function computes the modified solution for the well-balanced treatment of the
41 gravitational source terms. The modified solution vector is given by
43 {\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) \\ e^* \end{array}\right]
47 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\right)
49 \f$\varrho\f$ and \f$\varphi\f$ are computed in #NavierStokes2DGravityField(). For flows without gravity, \f${\bf u}^* = {\bf u}\f$.
52 + Ghosh, D., Constantinescu, E.M., Well-Balanced Formulation of Gravitational Source
53 Terms for Conservative Finite-Difference Atmospheric Flow Solvers, AIAA Paper 2015-2889,
54 7th AIAA Atmospheric and Space Environments Conference, June 22-26, 2015, Dallas, TX,
55 http://dx.doi.org/10.2514/6.2015-2889
56 + Ghosh, D., Constantinescu, E.M., A Well-Balanced, Conservative Finite-Difference Algorithm
57 for Atmospheric Flows, AIAA Journal, 54 (4), 2016, pp. 1370-1385, http://dx.doi.org/10.2514/1.J054580.
59 extern "C" int gpuNavierStokes2DModifiedSolution(
60 double *uC, /*!< Array to hold the computed modified solution */
61 double *u, /*!< Solution vector array */
62 int d, /*!< spatial dimension (not used) */
63 void *s, /*!< Solver object of type #HyPar */
64 void *m, /*!< MPI object of time #MPIVariables */
65 double waqt /*!< Current simulation time */
68 HyPar *solver = (HyPar*) s;
69 NavierStokes2D *param = (NavierStokes2D*) solver->physics;
71 int ghosts = solver->ghosts;
72 int *dim = solver->dim_local;
73 int ndims = solver->ndims;
75 int size = 1; for (int i=0; i <ndims; i++) size *= (dim[i] + 2*ghosts);
76 int nblocks = (size - 1) / GPU_THREADS_PER_BLOCK + 1;
77 double inv_gamma_m1 = 1.0 / (param->gamma-1.0);
79 double cpu_time = 0.0;
80 clock_t cpu_start, cpu_end;
83 NavierStokes2DModifiedSolution_kernel<<<nblocks, GPU_THREADS_PER_BLOCK>>>(
84 size, param->gamma, inv_gamma_m1,
85 param->gpu_grav_field_f, param->gpu_grav_field_g,
87 cudaDeviceSynchronize();
89 cpu_time += (double)(cpu_end - cpu_start) / CLOCKS_PER_SEC;