cuIBM
A GPU-based Immersed Boundary Method code
generateBC1.cu
Go to the documentation of this file.
1 
8 #include "generateBC1.h"
9 
10 
15 namespace kernels
16 {
17 
34 __global__
35 void bc1DirichletU(real *bc1, int N, int nx, int offset, int stride, real *dx, real C, real *bc)
36 {
37  int idx = threadIdx.x + blockIdx.x*blockDim.x;
38 
39  if ( idx >= N ) return;
40 
41  int I = (offset + idx*stride),
42  i = I % (nx-1);
43  bc1[I] += bc[idx] * C * 0.5*(dx[i] + dx[i+1]);
44 } // bc1DirichletU
45 
46 
65 __global__
66 void bc1DirichletV(real *bc1, int N, int nx, int numU, int offset, int stride, real *dy, real C, real *bc, int numUbc)
67 {
68  int idx = threadIdx.x + blockIdx.x*blockDim.x;
69 
70  if ( idx >= N ) return;
71 
72  int I = (offset + idx*stride),
73  j = I / nx;
74  bc1[numU + I] += bc[idx+numUbc] * C * 0.5*(dy[j] + dy[j+1]);
75 } // bc1DirichletV
76 
77 
97 __global__
98 void bc1ConvectiveU(real *bc1, int N, int nx, int offset, int stride, real *dx, real *dy, real C, real *bc, real *q, real beta)
99 {
100  int idx = threadIdx.x + blockIdx.x*blockDim.x;
101 
102  if ( idx >= N ) return;
103 
104  int I = (offset + idx*stride),
105  i = I % (nx-1),
106  j = I / (nx-1);
107 
108  bc[idx] = (1.0-beta)*bc[idx] + beta*q[I]/dy[j];
109 
110  bc1[I] += bc[idx] * C * 0.5*(dx[i] + dx[i+1]);
111 } // bc1ConvectiveU
112 
113 
135 __global__
136 void bc1ConvectiveV(real *bc1, int N, int nx, int numU, int offset, int stride, real *dx, real *dy, real C, real *bc, int numUbc, real *q, real beta)
137 {
138  int idx = threadIdx.x + blockIdx.x*blockDim.x;
139 
140  if ( idx >= N ) return;
141 
142  int I = (offset + idx*stride),
143  i = I % nx,
144  j = I / nx;
145 
146  bc[idx+numUbc] = (1.0-beta)*bc[idx+numUbc] + beta*q[numU + I]/dx[i];
147 
148  bc1[numU + I] += bc[idx+numUbc] * C * 0.5*(dy[j] + dy[j+1]);
149 } // bc1ConvectiveV
150 
151 
170 __global__
171 void bc1SpecialU(real *bc1, int N, int nx, int offset, int stride, real *dx, real C, real *bc, real time)
172 {
173  int idx = threadIdx.x + blockIdx.x*blockDim.x;
174 
175  if ( idx >= N ) return;
176 
177  int I = (offset + idx*stride),
178  i = I % (nx-1);
179 
180  const real T = 10.0;
181 
182  bc[idx] = sin(M_PI*time/T);
183 
184  bc1[I] += bc[idx] * C * 0.5*(dx[i] + dx[i+1]);
185 } // bc1SpecialU
186 
187 } // End of namespace kernels
Declaration of the kernels to generate right hand-side terms of the intermediate velocity flux solver...
__global__ void bc1DirichletV(real *bc1, int N, int nx, int numU, int offset, int stride, real *dy, real C, real *bc, int numUbc)
Computes inhomogeneous term from the discrete Laplacian operator for the v-velocity at a given bounda...
Definition: generateBC1.cu:66
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 bc1DirichletU(real *bc1, int N, int nx, int offset, int stride, real *dx, real C, real *bc)
Computes inhomogeneous term from the discrete Laplacian operator for the u-velocity at a given bounda...
Definition: generateBC1.cu:35
__global__ void bc1SpecialU(real *bc1, int N, int nx, int offset, int stride, real *dx, real C, real *bc, real time)
Computes inhomogeneous term from the discrete Laplacian operator for the u-velocity at a given bounda...
Definition: generateBC1.cu:171
__global__ void bc1ConvectiveV(real *bc1, int N, int nx, int numU, int offset, int stride, real *dx, real *dy, real C, real *bc, int numUbc, real *q, real beta)
Computes inhomogeneous term from the discrete Laplacian operator for the v-velocity at a given bounda...
Definition: generateBC1.cu:136
__global__ void bc1ConvectiveU(real *bc1, int N, int nx, int offset, int stride, real *dx, real *dy, real C, real *bc, real *q, real beta)
Computes inhomogeneous term from the discrete Laplacian operator for the u-velocity at a given bounda...
Definition: generateBC1.cu:98