cuIBM
A GPU-based Immersed Boundary Method code
generateBC2.cu
Go to the documentation of this file.
1 
8 #include "generateBC2.h"
9 
10 
15 namespace kernels
16 {
17 
29 __global__
30 void fillBC2_v(real *bc2, real *yminus, real *yplus, real *dx, int nx, int ny)
31 {
32  int i = threadIdx.x + blockIdx.x*blockDim.x;
33  if(i>=nx)
34  return;
35  bc2[i] -= yminus[i+nx-1]*dx[i];
36  bc2[(ny-1)*nx + i] += yplus[i+nx-1]*dx[i];
37 } // fillBC2_v
38 
39 
51 __global__
52 void fillBC2_u(real *bc2, real *xminus, real *xplus, real *dy, int nx, int ny)
53 {
54  int j = threadIdx.x + blockIdx.x*blockDim.x;
55  if(j>=ny)
56  return;
57  bc2[j*nx] -= xminus[j]*dy[j];
58  bc2[j*nx+nx-1] += xplus[j]*dy[j];
59 } // fillBC2_u
60 
61 
73 __global__
74 void fillBC2_uvB(real *bc2, real *uB, real *vB, int totalPoints, int nx, int ny)
75 {
76  int k = threadIdx.x + blockIdx.x*blockDim.x;
77  if(k>=totalPoints)
78  return;
79  bc2[nx*ny + k] = uB[k];
80  bc2[nx*ny + k + totalPoints] = vB[k];
81 } // fillBC2_uvB
82 
83 } // End of namespace kernels
__global__ void fillBC2_v(real *bc2, real *yminus, real *yplus, real *dx, int nx, int ny)
Computes inhomogeneous terms of the discrete divergence operator from the bottom and top boundaries a...
Definition: generateBC2.cu:30
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 fillBC2_u(real *bc2, real *xminus, real *xplus, real *dy, int nx, int ny)
Computes inhomogeneous terms of the discrete divergence operator from the left and right boundaries a...
Definition: generateBC2.cu:52
__global__ void fillBC2_uvB(real *bc2, real *uB, real *vB, int totalPoints, int nx, int ny)
Computes inhomogeneous terms of the discrete divergence operator from the no-slip constraint at the b...
Definition: generateBC2.cu:74
Declaration of the kernels to generate elements of the right hand-side of the Poisson solver...