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

Pre-step function 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 gpuNavierStokes3DPreStep_kernel (int npoints_grid, double gamma, const double *__restrict__ u, double *__restrict__ fast_jac)
 
int gpuNavierStokes3DPreStep (double *u, void *s, void *m, double waqt)
 

Detailed Description

Pre-step function for 3D Navier Stokes equations.

Author
Youngdae Kim

Definition in file NavierStokes3DPreStep_GPU.cu.

Function Documentation

__global__ void gpuNavierStokes3DPreStep_kernel ( int  npoints_grid,
double  gamma,
const double *__restrict__  u,
double *__restrict__  fast_jac 
)

Kernel for gpuNavierStokes3DPreStep()

Definition at line 117 of file NavierStokes3DPreStep_GPU.cu.

123 {
124  int p = blockDim.x * blockIdx.x + threadIdx.x;
125 
126  if (p < npoints_grid) {
127  double * __restrict__ A;
130  int JacSize = _MODEL_NVARS_*_MODEL_NVARS_;
131  int dir = _XDIR_;
132 
133  A = (fast_jac + _MODEL_NDIMS_*JacSize*p + dir*JacSize);
134  /* get the eigenvalues, and left & right eigenvectors */
135  _NavierStokes3DEigenvalues_ ((u+p),npoints_grid,D,gamma,dir);
136  _NavierStokes3DLeftEigenvectors_ ((u+p),npoints_grid,L,gamma,dir);
137  _NavierStokes3DRightEigenvectors_((u+p),npoints_grid,R,gamma,dir);
138  /* remove the entropy modes (corresponding to eigenvalues u) */
139  D[0] = D[12] = D[18] = 0.0;
140  /* assemble the Jacobian */
141  MatMult5(_MODEL_NVARS_,DL,D,L );
142  MatMult5(_MODEL_NVARS_,A ,R,DL);
143 
144  dir = _YDIR_;
145  A = (fast_jac + _MODEL_NDIMS_*JacSize*p + dir*JacSize);
146  /* get the eigenvalues, and left & right eigenvectors */
147  _NavierStokes3DEigenvalues_ ((u+p),npoints_grid,D,gamma,dir)
148  _NavierStokes3DLeftEigenvectors_ ((u+p),npoints_grid,L,gamma,dir);
149  _NavierStokes3DRightEigenvectors_((u+p),npoints_grid,R,gamma,dir);
150  /* remove the entropy modes (corresponding to eigenvalues v) */
151  D[0] = D[6] = D[18] = 0.0;
152  /* assemble the Jacobian */
153  MatMult5(_MODEL_NVARS_,DL,D,L );
154  MatMult5(_MODEL_NVARS_,A ,R,DL);
155 
156  dir = _ZDIR_;
157  A = (fast_jac + _MODEL_NDIMS_*JacSize*p + dir*JacSize);
158  /* get the eigenvalues, and left & right eigenvectors */
159  _NavierStokes3DEigenvalues_ ((u+p),npoints_grid,D,gamma,dir)
160  _NavierStokes3DLeftEigenvectors_ ((u+p),npoints_grid,L,gamma,dir);
161  _NavierStokes3DRightEigenvectors_((u+p),npoints_grid,R,gamma,dir);
162  /* remove the entropy modes (corresponding to eigenvalues v) */
163  D[0] = D[6] = D[12] = 0.0;
164  /* assemble the Jacobian */
165  MatMult5(_MODEL_NVARS_,DL,D,L );
166  MatMult5(_MODEL_NVARS_,A ,R,DL);
167  }
168  return;
169 }
#define _YDIR_
Definition: euler2d.h:41
#define _NavierStokes3DLeftEigenvectors_(u, stride, L, ga, dir)
#define _MODEL_NVARS_
Definition: euler1d.h:58
#define _MODEL_NDIMS_
Definition: euler1d.h:56
#define _ZDIR_
#define _NavierStokes3DEigenvalues_(u, stride, D, gamma, dir)
#define MatMult5(N, A, X, Y)
#define _NavierStokes3DRightEigenvectors_(u, stride, R, ga, dir)
#define _XDIR_
Definition: euler1d.h:75
int gpuNavierStokes3DPreStep ( double *  u,
void *  s,
void *  m,
double  waqt 
)

Pre-step function for the 3D Navier Stokes equations: This function is called at the beginning of each time step.

  • The solution at the beginning of the step is copied into NavierStokes3D::solution for the linearized flux partitioning.
  • The linearized "fast" Jacobian (representing the acoustic modes) NavierStokes3D::fast_jac is computed as:

    \begin{equation} A_f\left({\bf u}\right) = R\left({\bf u}\right)\Lambda_f\left({\bf u}\right)L\left({\bf u}\right) \end{equation}

    where \(R\) and \(L\) are the matrices of right and left eigenvectors, and,

    \begin{equation} \Lambda_f = diag\left[0,0,0,u+c,u-c \right] \end{equation}

    with \(c\) being the speed of sound.

    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.
Parameters
uSolution vector
sSolver object of type HyPar
mMPI object of type MPIVariables
waqtCurrent simulation time

Definition at line 192 of file NavierStokes3DPreStep_GPU.cu.

198 {
199  HyPar *solver = (HyPar*) s;
200  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
201 
202  int nblocks = (solver->npoints_local_wghosts-1)/GPU_THREADS_PER_BLOCK + 1;
203 
205  gpuNavierStokes3DPreStep_kernel<<<nblocks, GPU_THREADS_PER_BLOCK>>>(
206  solver->npoints_local_wghosts, param->gamma, u, param->gpu_fast_jac
207  );
208  cudaDeviceSynchronize();
209 
210  return 0;
211 }
int npoints_local_wghosts
Definition: hypar.h:42
double * gpu_solution
#define _MODEL_NVARS_
Definition: euler1d.h:58
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.
double * gpu_fast_jac
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
void gpuArrayCopy1D(const double *, double *, int)