cuIBM
A GPU-based Immersed Boundary Method code
generateQT.inl
Go to the documentation of this file.
1 
9 
10 #define BLOCKSIZE 256
11 
12 
19 template <>
21 {
22  logger.startTimer("updateQT");
23 
24  int *QTRows = thrust::raw_pointer_cast(&(QT.row_indices[0])),
25  *QTCols = thrust::raw_pointer_cast(&(QT.column_indices[0]));
26  real *QTVals = thrust::raw_pointer_cast(&(QT.values[0]));
27 
28  int *ERows = thrust::raw_pointer_cast(&(E.row_indices[0])),
29  *ECols = thrust::raw_pointer_cast(&(E.column_indices[0]));
30  real *EVals = thrust::raw_pointer_cast(&(E.values[0]));
31 
32  real *x = thrust::raw_pointer_cast(&(domInfo->xD[0])),
33  *y = thrust::raw_pointer_cast(&(domInfo->yD[0])),
34  *dx = thrust::raw_pointer_cast(&(domInfo->dxD[0]));
35 
36  real *xB = thrust::raw_pointer_cast(&(B.x[0])),
37  *yB = thrust::raw_pointer_cast(&(B.y[0]));
38 
39  int *I = thrust::raw_pointer_cast(&(B.I[0])),
40  *J = thrust::raw_pointer_cast(&(B.J[0]));
41 
42  int nx = domInfo->nx,
43  ny = domInfo->ny;
44 
45  dim3 dimGrid(int((B.totalPoints-0.5)/BLOCKSIZE)+1, 1);
46  dim3 dimBlock(BLOCKSIZE, 1);
47 
48  kernels::updateQT <<<dimGrid, dimBlock>>>
49  (QTRows, QTCols, QTVals,
50  ERows, ECols, EVals,
51  nx, ny, x, y, dx,
52  B.totalPoints, xB, yB, I, J);
53 
54  logger.stopTimer("updateQT");
55 
56  logger.startTimer("transposeQT");
57  cusp::transpose(QT, Q);
58  logger.stopTimer("transposeQT");
59 
60  logger.startTimer("transposeE");
61  cusp::transpose(E, ET);
62  logger.stopTimer("transposeE");
63 } // updateQT
64 
65 
69 template <>
71 {
72  logger.startTimer("generateQT");
73 
74  int nx = domInfo->nx,
75  ny = domInfo->ny;
76 
77  int numU = (nx-1)*ny;
78  int numUV = numU + nx*(ny-1);
79  int numP = numU + ny;
80 
82  cooH QTHost(numP + 2*B.totalPoints, numUV, 4*numP-2*(nx+ny) + 24*B.totalPoints);
83 
84  int *QTRows = thrust::raw_pointer_cast(&(QTHost.row_indices[0])),
85  *QTCols = thrust::raw_pointer_cast(&(QTHost.column_indices[0]));
86  real *QTVals = thrust::raw_pointer_cast(&(QTHost.values[0]));
87 
88  kernels::generateQT(QTRows, QTCols, QTVals, nx, ny);
89 
90  QT = QTHost;
91 
92  logger.stopTimer("generateQT");
93 
94  updateQT();
95 } // generateQT
__global__ void updateQT(int *QTRows, int *QTCols, real *QTVals, int QTSize, int *tags, real *coeffs)
Update the divergence operator at forcing nodes.
Definition: generateQT.cu:87
double real
Is a float or a double depending on the machine precision.
Definition: types.h:116
#define BLOCKSIZE
Definition: generateQT.inl:10
void generateQT(int *QTRows, int *QTCols, real *QTVals, int nx, int ny)
Generates the divergence matrix (on the host).
Definition: generateQT.cu:110
cusp::coo_matrix< int, real, host_memory > cooH
COO matrix stored on the host.
Definition: types.h:126
Declaration of the kernels to generate gradient matrix and interpolation matrix.
virtual void generateQT()