/*
This OpenCL kernel implements the second order HLL flux
Copyright (C) 2016 SINTEF ICT
This program is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with this program. If not, see .
*/
#include "common.cu"
/**
* Computes the flux along the x axis for all faces
*/
__device__
void computeFluxF(float Q[3][block_height+4][block_width+4],
float Qx[3][block_height+2][block_width+2],
float F[3][block_height+1][block_width+1],
const float g_, const float dx_, const float dt_) {
//Index of thread within block
const int tx = get_local_id(0);
const int ty = get_local_id(1);
for (int j=ty; j evolve x first, then y
if (step_ == 0) {
//Compute fluxes along the x axis and evolve
minmodSlopeX(Q, Qx, theta_);
__syncthreads();
computeFluxF(Q, Qx, F, g_, dx_, dt_);
__syncthreads();
evolveF2(Q, F, nx_, ny_, dx_, dt_);
__syncthreads();
//Set boundary conditions
noFlowBoundary2(Q, nx_, ny_);
__syncthreads();
//Compute fluxes along the y axis and evolve
minmodSlopeY(Q, Qx, theta_);
__syncthreads();
computeFluxG(Q, Qx, F, g_, dy_, dt_);
__syncthreads();
evolveG2(Q, F, nx_, ny_, dy_, dt_);
__syncthreads();
}
//Step 1 => evolve y first, then x
else {
//Compute fluxes along the y axis and evolve
minmodSlopeY(Q, Qx, theta_);
__syncthreads();
computeFluxG(Q, Qx, F, g_, dy_, dt_);
__syncthreads();
evolveG2(Q, F, nx_, ny_, dy_, dt_);
__syncthreads();
//Set boundary conditions
noFlowBoundary2(Q, nx_, ny_);
__syncthreads();
//Compute fluxes along the x axis and evolve
minmodSlopeX(Q, Qx, theta_);
__syncthreads();
computeFluxF(Q, Qx, F, g_, dx_, dt_);
__syncthreads();
evolveF2(Q, F, nx_, ny_, dx_, dt_);
__syncthreads();
}
// Write to main memory for all internal cells
writeBlock2(h1_ptr_, h1_pitch_,
hu1_ptr_, hu1_pitch_,
hv1_ptr_, hv1_pitch_,
Q, nx_, ny_);
}