cuIBM
A GPU-based Immersed Boundary Method code
generateA.cu
Go to the documentation of this file.
1 
8 #include "generateA.h"
9 
10 
15 namespace kernels
16 {
17 
37 __global__
38 void generateA(int *ARows, int *ACols, real *AVals,
39  real *MVals,
40  int *LRows, int *LCols, real *LVals,
41  int ASize, real alpha)
42 {
43  for (int I=threadIdx.x + blockIdx.x*blockDim.x; I<ASize; I += blockDim.x*gridDim.x)
44  {
45  ARows[I] = LRows[I];
46  ACols[I] = LCols[I];
47  AVals[I] = -alpha*LVals[I] + (LRows[I]==LCols[I])*MVals[LRows[I]];
48  }
49 } // generateA
50 
51 
76 __global__
77 void generateADirectForcing(int *ARows, int *ACols, real *AVals,
78  real *MVals,
79  int *LRows, int *LCols, real *LVals,
80  int ASize, real alpha, int *tags)
81 {
82  for(int I=threadIdx.x + blockIdx.x*blockDim.x; I<ASize; I += blockDim.x*gridDim.x)
83  {
84  ARows[I] = LRows[I];
85  ACols[I] = LCols[I];
86  AVals[I] = (tags[LRows[I]] == -1)*(-alpha*LVals[I]) // if the current location is untagged, add -alpha*L
87  + (tags[LRows[I]] != -1)*(-LVals[I]) // if the current location is tagged, add -L
88  + (LRows[I]==LCols[I])*MVals[LRows[I]]; // if it is a diagonal, add M
89  }
90 } // generateADirectForcing
91 
92 } // End of namespace kernels
double real
Is a float or a double depending on the machine precision.
Definition: types.h:116
__global__ void generateADirectForcing(int *ARows, int *ACols, real *AVals, real *MVals, int *LRows, int *LCols, real *LVals, int ASize, real alpha, int *tags)
Generates a block of the matrix resulting from implicit terms in the momentum equation for the direct...
Definition: generateA.cu:77
Contains all the custom-written CUDA kernels.
__global__ void generateA(int *ARows, int *ACols, real *AVals, real *MVals, int *LRows, int *LCols, real *LVals, int ASize, real alpha)
Generates a block of the matrix resulting from implicit terms in the momentum equation.
Definition: generateA.cu:38
Declaration of the kernels to generate the matrix resulting from the implicit terms in the momentum e...