2020-09-11 22:56:00 -04:00
|
|
|
#include <math.h>
|
|
|
|
|
#include <stdio.h>
|
|
|
|
|
#include <cuda_profiler_api.h>
|
|
|
|
|
|
|
|
|
|
#define NBLOCKS 1024
|
|
|
|
|
#define NTHREADS 256
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_Solid_Dirichlet_D3Q7(double *dist, double *BoundaryValue, int *BounceBackDist_list, int *BounceBackSolid_list, int count)
|
|
|
|
|
{
|
|
|
|
|
|
|
|
|
|
int idx;
|
|
|
|
|
int iq,ib;
|
|
|
|
|
double value_b,value_q;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
iq = BounceBackDist_list[idx];
|
|
|
|
|
ib = BounceBackSolid_list[idx];
|
|
|
|
|
value_b = BoundaryValue[ib];//get boundary value from a solid site
|
|
|
|
|
value_q = dist[iq];
|
|
|
|
|
dist[iq] = -1.0*value_q + value_b*0.25;//NOTE 0.25 is the speed of sound for D3Q7 lattice
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_Solid_Neumann_D3Q7(double *dist, double *BoundaryValue, int *BounceBackDist_list, int *BounceBackSolid_list, int count)
|
|
|
|
|
{
|
|
|
|
|
|
|
|
|
|
int idx;
|
|
|
|
|
int iq,ib;
|
|
|
|
|
double value_b,value_q;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
iq = BounceBackDist_list[idx];
|
|
|
|
|
ib = BounceBackSolid_list[idx];
|
|
|
|
|
value_b = BoundaryValue[ib];//get boundary value from a solid site
|
|
|
|
|
value_q = dist[iq];
|
|
|
|
|
dist[iq] = value_q + value_b;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_D3Q7_AAeven_Poisson_Potential_BC_z(int *list, double *dist, double Vin, int count, int Np)
|
|
|
|
|
{
|
|
|
|
|
int idx,n;
|
|
|
|
|
double f0,f1,f2,f3,f4,f5,f6;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
f0 = dist[n];
|
|
|
|
|
f1 = dist[2*Np+n];
|
|
|
|
|
f2 = dist[1*Np+n];
|
|
|
|
|
f3 = dist[4*Np+n];
|
|
|
|
|
f4 = dist[3*Np+n];
|
|
|
|
|
f6 = dist[5*Np+n];
|
|
|
|
|
//...................................................
|
|
|
|
|
f5 = Vin - (f0+f1+f2+f3+f4+f6);
|
|
|
|
|
dist[6*Np+n] = f5;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_D3Q7_AAeven_Poisson_Potential_BC_Z(int *list, double *dist, double Vout, int count, int Np)
|
|
|
|
|
{
|
|
|
|
|
int idx,n;
|
|
|
|
|
double f0,f1,f2,f3,f4,f5,f6;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
f0 = dist[n];
|
|
|
|
|
f1 = dist[2*Np+n];
|
|
|
|
|
f2 = dist[1*Np+n];
|
|
|
|
|
f3 = dist[4*Np+n];
|
|
|
|
|
f4 = dist[3*Np+n];
|
|
|
|
|
f5 = dist[6*Np+n];
|
|
|
|
|
//...................................................
|
|
|
|
|
f6 = Vout - (f0+f1+f2+f3+f4+f5);
|
|
|
|
|
dist[5*Np+n] = f6;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_z(int *d_neighborList, int *list, double *dist, double Vin, int count, int Np)
|
|
|
|
|
{
|
|
|
|
|
int idx, n;
|
|
|
|
|
int nread,nr5;
|
|
|
|
|
double f0,f1,f2,f3,f4,f5,f6;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
f0 = dist[n];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n];
|
|
|
|
|
f1 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+2*Np];
|
|
|
|
|
f3 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+Np];
|
|
|
|
|
f2 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+3*Np];
|
|
|
|
|
f4 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+5*Np];
|
|
|
|
|
f6 = dist[nread];
|
|
|
|
|
|
|
|
|
|
// Unknown distributions
|
|
|
|
|
nr5 = d_neighborList[n+4*Np];
|
|
|
|
|
f5 = Vin - (f0+f1+f2+f3+f4+f6);
|
|
|
|
|
dist[nr5] = f5;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
__global__ void dvc_ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_Z(int *d_neighborList, int *list, double *dist, double Vout, int count, int Np)
|
2020-09-11 22:56:00 -04:00
|
|
|
{
|
|
|
|
|
int idx, n;
|
|
|
|
|
int nread,nr6;
|
|
|
|
|
double f0,f1,f2,f3,f4,f5,f6;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
f0 = dist[n];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n];
|
|
|
|
|
f1 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+2*Np];
|
|
|
|
|
f3 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+4*Np];
|
|
|
|
|
f5 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+Np];
|
|
|
|
|
f2 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+3*Np];
|
|
|
|
|
f4 = dist[nread];
|
|
|
|
|
|
|
|
|
|
// unknown distributions
|
|
|
|
|
nr6 = d_neighborList[n+5*Np];
|
|
|
|
|
f6 = Vout - (f0+f1+f2+f3+f4+f5);
|
|
|
|
|
dist[nr6] = f6;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_Poisson_D3Q7_BC_z(int *list, int *Map, double *Psi, double Vin, int count)
|
|
|
|
|
{
|
|
|
|
|
int idx,n,nm;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
nm = Map[n];
|
|
|
|
|
Psi[nm] = Vin;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_Poisson_D3Q7_BC_Z(int *list, int *Map, double *Psi, double Vout, int count)
|
|
|
|
|
{
|
|
|
|
|
int idx,n,nm;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
nm = Map[n];
|
|
|
|
|
Psi[nm] = Vout;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_D3Q7_AAeven_Ion_Concentration_BC_z(int *list, double *dist, double Cin, int count, int Np)
|
|
|
|
|
{
|
|
|
|
|
int idx,n;
|
|
|
|
|
double f0,f1,f2,f3,f4,f5,f6;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
f0 = dist[n];
|
|
|
|
|
f1 = dist[2*Np+n];
|
|
|
|
|
f2 = dist[1*Np+n];
|
|
|
|
|
f3 = dist[4*Np+n];
|
|
|
|
|
f4 = dist[3*Np+n];
|
|
|
|
|
f6 = dist[5*Np+n];
|
|
|
|
|
//...................................................
|
|
|
|
|
f5 = Cin - (f0+f1+f2+f3+f4+f6);
|
|
|
|
|
dist[6*Np+n] = f5;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_D3Q7_AAeven_Ion_Concentration_BC_Z(int *list, double *dist, double Cout, int count, int Np)
|
|
|
|
|
{
|
|
|
|
|
int idx,n;
|
|
|
|
|
double f0,f1,f2,f3,f4,f5,f6;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
2020-09-20 11:00:36 -04:00
|
|
|
if (idx < count){
|
2020-09-11 22:56:00 -04:00
|
|
|
n = list[idx];
|
|
|
|
|
f0 = dist[n];
|
|
|
|
|
f1 = dist[2*Np+n];
|
|
|
|
|
f2 = dist[1*Np+n];
|
|
|
|
|
f3 = dist[4*Np+n];
|
|
|
|
|
f4 = dist[3*Np+n];
|
|
|
|
|
f5 = dist[6*Np+n];
|
|
|
|
|
//...................................................
|
|
|
|
|
f6 = Cout - (f0+f1+f2+f3+f4+f5);
|
|
|
|
|
dist[5*Np+n] = f6;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_D3Q7_AAodd_Ion_Concentration_BC_z(int *d_neighborList, int *list, double *dist, double Cin, int count, int Np)
|
|
|
|
|
{
|
|
|
|
|
int idx, n;
|
|
|
|
|
int nread,nr5;
|
|
|
|
|
double f0,f1,f2,f3,f4,f5,f6;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
f0 = dist[n];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n];
|
|
|
|
|
f1 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+2*Np];
|
|
|
|
|
f3 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+Np];
|
|
|
|
|
f2 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+3*Np];
|
|
|
|
|
f4 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+5*Np];
|
|
|
|
|
f6 = dist[nread];
|
|
|
|
|
|
|
|
|
|
// Unknown distributions
|
|
|
|
|
nr5 = d_neighborList[n+4*Np];
|
|
|
|
|
f5 = Cin - (f0+f1+f2+f3+f4+f6);
|
|
|
|
|
dist[nr5] = f5;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
__global__ void dvc_ScaLBL_D3Q7_AAodd_Ion_Concentration_BC_Z(int *d_neighborList, int *list, double *dist, double Cout, int count, int Np)
|
2020-09-11 22:56:00 -04:00
|
|
|
{
|
|
|
|
|
int idx, n;
|
|
|
|
|
int nread,nr6;
|
|
|
|
|
double f0,f1,f2,f3,f4,f5,f6;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx < count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
f0 = dist[n];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n];
|
|
|
|
|
f1 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+2*Np];
|
|
|
|
|
f3 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+4*Np];
|
|
|
|
|
f5 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+Np];
|
|
|
|
|
f2 = dist[nread];
|
|
|
|
|
|
|
|
|
|
nread = d_neighborList[n+3*Np];
|
|
|
|
|
f4 = dist[nread];
|
|
|
|
|
|
|
|
|
|
// unknown distributions
|
|
|
|
|
nr6 = d_neighborList[n+5*Np];
|
|
|
|
|
f6 = Cout - (f0+f1+f2+f3+f4+f5);
|
|
|
|
|
dist[nr6] = f6;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//*************************************************************************
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_Solid_Dirichlet_D3Q7(double *dist, double *BoundaryValue, int *BounceBackDist_list, int *BounceBackSolid_list, int count){
|
2020-09-11 22:56:00 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2020-09-20 11:00:36 -04:00
|
|
|
dvc_ScaLBL_Solid_Dirichlet_D3Q7<<<GRID,512>>>(dist, BoundaryValue, BounceBackDist_list, BounceBackSolid_list, count);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_Solid_Dirichlet_D3Q7 (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_Solid_Neumann_D3Q7(double *dist, double *BoundaryValue, int *BounceBackDist_list, int *BounceBackSolid_list, int count){
|
|
|
|
|
int GRID = count / 512 + 1;
|
|
|
|
|
dvc_ScaLBL_Solid_Neumann_D3Q7<<<GRID,512>>>(dist, BoundaryValue, BounceBackDist_list, BounceBackSolid_list, count);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_Solid_Neumann_D3Q7 (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_D3Q7_AAeven_Poisson_Potential_BC_z(int *list, double *dist, double Vin, int count, int Np){
|
|
|
|
|
int GRID = count / 512 + 1;
|
|
|
|
|
dvc_ScaLBL_D3Q7_AAeven_Poisson_Potential_BC_z<<<GRID,512>>>(list, dist, Vin, count, Np);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_D3Q7_AAeven_Poisson_Potential_BC_z (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_D3Q7_AAeven_Poisson_Potential_BC_Z(int *list, double *dist, double Vout, int count, int Np){
|
|
|
|
|
int GRID = count / 512 + 1;
|
|
|
|
|
dvc_ScaLBL_D3Q7_AAeven_Poisson_Potential_BC_Z<<<GRID,512>>>(list, dist, Vout, count, Np);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_D3Q7_AAeven_Poisson_Potential_BC_Z (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_z(int *d_neighborList, int *list, double *dist, double Vin, int count, int Np){
|
2020-09-11 22:56:00 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2020-09-20 11:00:36 -04:00
|
|
|
dvc_ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_z<<<GRID,512>>>(d_neighborList, list, dist, Vin, count, Np);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_z (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_Z(int *d_neighborList, int *list, double *dist, double Vout, int count, int Np){
|
2020-09-11 22:56:00 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2020-09-20 11:00:36 -04:00
|
|
|
dvc_ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_Z<<<GRID,512>>>(d_neighborList, list, dist, Vout, count, Np);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_Z (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_Poisson_D3Q7_BC_z(int *list, int *Map, double *Psi, double Vin, int count){
|
2020-09-11 22:56:00 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2020-09-20 11:00:36 -04:00
|
|
|
dvc_ScaLBL_Poisson_D3Q7_BC_z<<<GRID,512>>>(list, Map, Psi, Vin, count);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_Poisson_D3Q7_BC_z (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_Poisson_D3Q7_BC_Z(int *list, int *Map, double *Psi, double Vout, int count){
|
2020-09-11 22:56:00 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2020-09-20 11:00:36 -04:00
|
|
|
dvc_ScaLBL_Poisson_D3Q7_BC_Z<<<GRID,512>>>(list, Map, Psi, Vout, count);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_Poisson_D3Q7_BC_Z (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_D3Q7_AAeven_Ion_Concentration_BC_z(int *list, double *dist, double Cin, int count, int Np){
|
2020-09-11 22:56:00 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2020-09-20 11:00:36 -04:00
|
|
|
dvc_ScaLBL_D3Q7_AAeven_Ion_Concentration_BC_z<<<GRID,512>>>(list, dist, Cin, count, Np);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_D3Q7_AAeven_Ion_Concentration_BC_z (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_D3Q7_AAeven_Ion_Concentration_BC_Z(int *list, double *dist, double Cout, int count, int Np){
|
2020-09-11 22:56:00 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2020-09-20 11:00:36 -04:00
|
|
|
dvc_ScaLBL_D3Q7_AAeven_Ion_Concentration_BC_Z<<<GRID,512>>>(list, dist, Cout, count, Np);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_D3Q7_AAeven_Ion_Concentration_BC_Z (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_D3Q7_AAodd_Ion_Concentration_BC_z(int *d_neighborList, int *list, double *dist, double Cin, int count, int Np){
|
2020-09-11 22:56:00 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2020-09-20 11:00:36 -04:00
|
|
|
dvc_ScaLBL_D3Q7_AAodd_Ion_Concentration_BC_z<<<GRID,512>>>(d_neighborList, list, dist, Cin, count, Np);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_D3Q7_AAodd_Ion_Concentration_BC_z (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-09-20 11:00:36 -04:00
|
|
|
extern "C" void ScaLBL_D3Q7_AAodd_Ion_Concentration_BC_Z(int *d_neighborList, int *list, double *dist, double Cout, int count, int Np){
|
2020-09-11 22:56:00 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2020-09-20 11:00:36 -04:00
|
|
|
dvc_ScaLBL_D3Q7_AAodd_Ion_Concentration_BC_Z<<<GRID,512>>>(d_neighborList, list, dist, Cout, count, Np);
|
2020-09-11 22:56:00 -04:00
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
2020-09-20 11:00:36 -04:00
|
|
|
printf("CUDA error in ScaLBL_D3Q7_AAodd_Ion_Concentration_BC_Z (kernel): %s \n",cudaGetErrorString(err));
|
2020-09-11 22:56:00 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|