cuIBM
A GPU-based Immersed Boundary Method code
generateRN.h
Go to the documentation of this file.
1 
7 #pragma once
8 
9 #include "utilities/types.h"
10 
11 
16 namespace kernels
17 {
18 
19 // compute u-component of the explicit terms of the discretized momentum equation
20 // and element of the explicit convective terms
21 __global__
22 void convectionTermU(real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu);
23 
24 // compute v-component of the explicit terms of the discretized momentum equation
25 // and element of the explicit convective terms
26 __global__
27 void convectionTermV(real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu);
28 
29 // compute v-component of the explicit terms of the discretized momentum equation
30 // and element of the explicit convective terms
31 // top or bottom boundaries
32 __global__
33 void convectionTermVBottomTop(real *rn, real *H, real *q, \
34  int nx, int ny, real *dx, real *dy, \
35  real dt, real gamma, real zeta, real alpha, real nu, \
36  real *bcBottom, real *bcTop);
37 
38 // compute u-component of the explicit terms of the discretized momentum equation
39 // and element of the explicit convective terms
40 // top or bottom boundaries
41 __global__
42 void convectionTermUBottomTop(real *rn, real *H, real *q, \
43  int nx, int ny, real *dx, real *dy, \
44  real dt, real gamma, real zeta, real alpha, real nu, \
45  real *bcBottom, real *bcTop, real *bcLeft, real *bcRight);
46 
47 // compute u-component of the explicit terms of the discretized momentum equation
48 // and element of the explicit convective terms
49 // left or right boundaries
50 __global__
51 void convectionTermULeftRight(real *rn, real *H, real *q, \
52  int nx, int ny, real *dx, real *dy, \
53  real dt, real gamma, real zeta, real alpha, real nu, \
54  real *bcLeft, real *bcRight);
55 
56 // compute v-component of the explicit terms of the discretized momentum equation
57 // and element of the explicit convective terms
58 // left or right boundaries
59 __global__
60 void convectionTermVLeftRight(real *rn, real *H, real *q, \
61  int nx, int ny, real *dx, real *dy, \
62  real dt, real gamma, real zeta, real alpha, real nu, \
63  real *bcBottom, real *bcTop, real *bcLeft, real *bcRight);
64 
65 } // End of namespace kernels
Definition of custom types required by the code.
double real
Is a float or a double depending on the machine precision.
Definition: types.h:116
Contains all the custom-written CUDA kernels.
__global__ void convectionTermVLeftRight(real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu, real *bcBottom, real *bcTop, real *bcLeft, real *bcRight)
Computes v-component of the vector resulting from the explicit terms of the discretized momentum equa...
Definition: generateRN.cu:480
__global__ void convectionTermUBottomTop(real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu, real *bcBottom, real *bcTop, real *bcLeft, real *bcRight)
Computes u-component of the vector resulting from the explicit terms of the discretized momentum equa...
Definition: generateRN.cu:291
__global__ void convectionTermVBottomTop(real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu, real *bcBottom, real *bcTop)
Computes v-component of the vector resulting from the explicit terms of the discretized momentum equa...
Definition: generateRN.cu:211
__global__ void convectionTermU(real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu)
Computes u-component of the vector resulting from the explicit terms of the discretized momentum equa...
Definition: generateRN.cu:38
__global__ void convectionTermV(real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu)
Computes v-component of the vector resulting from the explicit terms of the discretized momentum equa...
Definition: generateRN.cu:122
__global__ void convectionTermULeftRight(real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu, real *bcLeft, real *bcRight)
Computes u-component of the vector resulting from the explicit terms of the discretized momentum equa...
Definition: generateRN.cu:400