cuIBM
A GPU-based Immersed Boundary Method code
generateM.cu
Go to the documentation of this file.
1 
7 #include "generateM.h"
8 
9 
14 namespace kernels
15 {
16 
32 __global__
33 void fillM_u(int *MRows, int *MCols, real *MVals, int *MinvRows, int *MinvCols, real *MinvVals, int nx, int ny, real *dx, real *dy, real dt)
34 {
35  int I = threadIdx.x + blockIdx.x*blockDim.x;
36 
37  if(I >= (nx-1)*ny) return;
38 
39  int i = I % (nx-1);
40  int j = I / (nx-1);
41  real value = 0.5*(dx[i]+dx[i+1])/dy[j]/dt;
42 
43  MRows[I] = I;
44  MCols[I] = I;
45  MVals[I] = value;
46 
47  MinvRows[I] = I;
48  MinvCols[I] = I;
49  MinvVals[I] = 1.0/value;
50 } // fillM_u
51 
52 
68 __global__
69 void fillM_v(int *MRows, int *MCols, real *MVals, int *MinvRows, int *MinvCols, real *MinvVals, int nx, int ny, real *dx, real *dy, real dt)
70 {
71  int I = threadIdx.x + blockIdx.x*blockDim.x;
72 
73  if(I >= nx*(ny-1)) return;
74 
75  int numU = (nx-1)*ny;
76  int i = I % nx;
77  int j = I / nx;
78  real value = 0.5*(dy[j]+dy[j+1])/dx[i]/dt;
79 
80  MRows[I+numU] = I+numU;
81  MCols[I+numU] = I+numU;
82  MVals[I+numU] = value;
83 
84  MinvRows[I+numU] = I+numU;
85  MinvCols[I+numU] = I+numU;
86  MinvVals[I+numU] = 1.0/value;
87 } // fillM_v
88 
89 } // End of namespace kernels
__global__ void fillM_v(int *MRows, int *MCols, real *MVals, int *MinvRows, int *MinvCols, real *MinvVals, int nx, int ny, real *dx, real *dy, real dt)
Computes an element of the mass matrix and its inverse for a y-velocity node.
Definition: generateM.cu:69
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 fillM_u(int *MRows, int *MCols, real *MVals, int *MinvRows, int *MinvCols, real *MinvVals, int nx, int ny, real *dx, real *dy, real dt)
Computes an element of the mass matrix and its inverse for a x-velocity node.
Definition: generateM.cu:33
Declaration of the kernels to generate the mass matrix and its inverse.