Refactoring

This commit is contained in:
André R. Brodtkorb
2018-11-01 15:17:13 +01:00
parent 064027fc0b
commit 2b899d1c80
10 changed files with 234 additions and 194 deletions

View File

@@ -46,22 +46,22 @@ __device__ __inline__ float minmodSlope(float left, float center, float right, f
/**
* Reconstructs a minmod slope for a whole block along the abscissa
*/
__device__ void minmodSlopeX(float Q[3][BLOCK_HEIGHT+4][BLOCK_WIDTH+4],
float Qx[3][BLOCK_HEIGHT+2][BLOCK_WIDTH+2],
template<int block_width, int block_height, int ghost_cells, int vars>
__device__ void minmodSlopeX(float Q[vars][block_height+2*ghost_cells][block_width+2*ghost_cells],
float Qx[vars][block_height+2*(ghost_cells-1)][block_width+2*(ghost_cells-1)],
const float theta_) {
//Index of thread within block
const int tx = threadIdx.x;
const int ty = threadIdx.y;
const int j = ty;
const int l = j + ghost_cells; //Skip ghost cells
//Reconstruct slopes along x axis
{
const int j = ty;
const int l = j + 2; //Skip ghost cells
for (int i=tx; i<BLOCK_WIDTH+2; i+=BLOCK_WIDTH) {
const int k = i + 1;
for (int p=0; p<3; ++p) {
Qx[p][j][i] = minmodSlope(Q[p][l][k-1], Q[p][l][k], Q[p][l][k+1], theta_);
}
for (int i=tx; i<block_width+2*(ghost_cells-1); i+=block_width) {
const int k = i + 1;
for (int p=0; p<vars; ++p) {
Qx[p][j][i] = minmodSlope(Q[p][l][k-1], Q[p][l][k], Q[p][l][k+1], theta_);
}
}
}
@@ -70,21 +70,22 @@ __device__ void minmodSlopeX(float Q[3][BLOCK_HEIGHT+4][BLOCK_WIDTH+4],
/**
* Reconstructs a minmod slope for a whole block along the ordinate
*/
__device__ void minmodSlopeY(float Q[3][BLOCK_HEIGHT+4][BLOCK_WIDTH+4],
float Qy[3][BLOCK_HEIGHT+2][BLOCK_WIDTH+2],
template<int block_width, int block_height, int ghost_cells, int vars>
__device__ void minmodSlopeY(float Q[vars][block_height+2*ghost_cells][block_width+2*ghost_cells],
float Qy[vars][block_height+2*(ghost_cells-1)][block_width+2*(ghost_cells-1)],
const float theta_) {
//Index of thread within block
const int tx = threadIdx.x;
const int ty = threadIdx.y;
for (int j=ty; j<BLOCK_HEIGHT+2; j+=BLOCK_HEIGHT) {
const int i = tx;
const int k = i + ghost_cells; //Skip ghost cells
//Reconstruct slopes along y axis
for (int j=ty; j<block_height+2*(ghost_cells-1); j+=block_height) {
const int l = j + 1;
{
const int i = tx;
const int k = i + 2; //Skip ghost cells
for (int p=0; p<3; ++p) {
Qy[p][j][i] = minmodSlope(Q[p][l-1][k], Q[p][l][k], Q[p][l+1][k], theta_);
}
for (int p=0; p<vars; ++p) {
Qy[p][j][i] = minmodSlope(Q[p][l-1][k], Q[p][l][k], Q[p][l+1][k], theta_);
}
}
}