added D3Q19 reflection BVC

This commit is contained in:
James McClure
2020-04-03 09:30:55 -04:00
parent e641e2e3ed
commit e64d44e438
4 changed files with 107 additions and 5 deletions

View File

@@ -1633,6 +1633,17 @@ double ScaLBL_Communicator::D3Q19_Flux_BC_z(int *neighborList, double *fq, doubl
return din; return din;
} }
void ScaLBL_Communicator::D3Q19_Reflection_BC_z(int *neighborList, double *fq){
if (kproc == 0)
ScaLBL_D3Q19_AAeven_Reflection_BC_z(dvcSendList_z, fq, sendCount_z, N);
}
void ScaLBL_Communicator::D3Q19_Reflection_BC_Z(int *neighborList, double *fq){
if (kproc == nprocz-1)
ScaLBL_D3Q19_AAeven_Reflection_BC_Z(dvcSendList_Z, fq, sendCount_Z, N);
}
void ScaLBL_Communicator::PrintD3Q19(){ void ScaLBL_Communicator::PrintD3Q19(){
printf("Printing D3Q19 communication buffer contents \n"); printf("Printing D3Q19 communication buffer contents \n");

View File

@@ -137,6 +137,10 @@ extern "C" void ScaLBL_Color_BC_z(int *list, int *Map, double *Phi, double *Den,
extern "C" void ScaLBL_Color_BC_Z(int *list, int *Map, double *Phi, double *Den, double vA, double vB, int count, int Np); extern "C" void ScaLBL_Color_BC_Z(int *list, int *Map, double *Phi, double *Den, double vA, double vB, int count, int Np);
extern "C" void ScaLBL_D3Q19_AAeven_Reflection_BC_z(int *list, double *dist, int count, int Np);
extern "C" void ScaLBL_D3Q19_AAeven_Reflection_BC_Z(int *list, double *dist, int count, int Np);
extern "C" void ScaLBL_SetSlice_z(double *Phi, double value, int Nx, int Ny, int Nz, int Slice); extern "C" void ScaLBL_SetSlice_z(double *Phi, double value, int Nx, int Ny, int Nz, int Slice);
class ScaLBL_Communicator{ class ScaLBL_Communicator{
@@ -189,6 +193,8 @@ public:
void Color_BC_Z(int *Map, double *Phi, double *Den, double vA, double vB); void Color_BC_Z(int *Map, double *Phi, double *Den, double vA, double vB);
void D3Q19_Pressure_BC_z(int *neighborList, double *fq, double din, int time); void D3Q19_Pressure_BC_z(int *neighborList, double *fq, double din, int time);
void D3Q19_Pressure_BC_Z(int *neighborList, double *fq, double dout, int time); void D3Q19_Pressure_BC_Z(int *neighborList, double *fq, double dout, int time);
void D3Q19_Reflection_BC_z(int *neighborList, double *fq);
void D3Q19_Reflection_BC_Z(int *neighborList, double *fq);
double D3Q19_Flux_BC_z(int *neighborList, double *fq, double flux, int time); double D3Q19_Flux_BC_z(int *neighborList, double *fq, double flux, int time);
// Debugging and unit testing functions // Debugging and unit testing functions

View File

@@ -448,6 +448,42 @@ extern "C" double ScaLBL_D3Q19_Flux_BC_Z(double *disteven, double *distodd, doub
return dout; return dout;
} }
extern "C" void ScaLBL_D3Q19_AAeven_Reflection_BC_z(int *list, double *dist, int count, int Np){
for (int idx=0; idx<count; idx++){
int n = list[idx];
double f5 = 0.111111111111111111111111 - dist[6*Np+n];
double f11 = 0.05555555555555555555556 - dist[12*Np+n];
double f14 = 0.05555555555555555555556 - dist[13*Np+n];
double f15 = 0.05555555555555555555556 - dist[16*Np+n];
double f18 = 0.05555555555555555555556 - dist[17*Np+n];
dist[6*Np+n] = f5;
dist[12*Np+n] = f11;
dist[13*Np+n] = f14;
dist[16*Np+n] = f15;
dist[17*Np+n] = f18;
}
}
extern "C" void ScaLBL_D3Q19_AAeven_Reflection_BC_Z(int *list, double *dist, int count, int Np){
for (int idx=0; idx<count; idx++){
int n = list[idx];
double f6 = 0.111111111111111111111111 - dist[5*Np+n];
double f12 = 0.05555555555555555555556 - dist[11*Np+n];
double f13 = 0.05555555555555555555556 - dist[14*Np+n] ;
double f16 = 0.05555555555555555555556 - dist[15*Np+n];
double f17 = 0.05555555555555555555556 - dist[18*Np+n];
dist[5*Np+n] = f6;
dist[11*Np+n] = f12;
dist[14*Np+n] = f13;
dist[15*Np+n] = f16;
dist[18*Np+n] = f17;
}
}
extern "C" void ScaLBL_D3Q19_AAeven_Pressure_BC_z(int *list, double *dist, double din, int count, int Np) extern "C" void ScaLBL_D3Q19_AAeven_Pressure_BC_z(int *list, double *dist, double din, int count, int Np)
{ {
// distributions // distributions

View File

@@ -1758,6 +1758,43 @@ __global__ void dvc_ScaLBL_D3Q19_AAeven_Pressure_BC_Z(int *list, double *dist,
//................................................... //...................................................
} }
} }
__global__ void dvc_ScaLBL_D3Q19_Reflection_BC_z(int *d_neighborList, int *list, double *dist, int count, int Np){
int idx, n;
idx = blockIdx.x*blockDim.x + threadIdx.x;
if (idx < count){
n = list[idx];
double f5 = 0.111111111111111111111111 - dist[6*Np+n];
double f11 = 0.05555555555555555555556 - dist[12*Np+n];
double f14 = 0.05555555555555555555556 - dist[13*Np+n];
double f15 = 0.05555555555555555555556 - dist[16*Np+n];
double f18 = 0.05555555555555555555556 - dist[17*Np+n];
dist[6*Np+n] = f5;
dist[12*Np+n] = f11;
dist[13*Np+n] = f14;
dist[16*Np+n] = f15;
dist[17*Np+n] = f18;
}
}
__global__ void dvc_ScaLBL_D3Q19_Reflection_BC_Z(int *d_neighborList, int *list, double *dist, int count, int Np){
int idx, n;
idx = blockIdx.x*blockDim.x + threadIdx.x;
if (idx < count){
n = list[idx];
double f6 = 0.111111111111111111111111 - dist[5*Np+n];
double f12 = 0.05555555555555555555556 - dist[11*Np+n];
double f13 = 0.05555555555555555555556 - dist[14*Np+n] ;
double f16 = 0.05555555555555555555556 - dist[15*Np+n];
double f17 = 0.05555555555555555555556 - dist[18*Np+n];
dist[5*Np+n] = f6;
dist[11*Np+n] = f12;
dist[14*Np+n] = f13;
dist[15*Np+n] = f16;
dist[18*Np+n] = f17;
}
}
__global__ void dvc_ScaLBL_D3Q19_AAodd_Pressure_BC_z(int *d_neighborList, int *list, double *dist, double din, int count, int Np) __global__ void dvc_ScaLBL_D3Q19_AAodd_Pressure_BC_z(int *d_neighborList, int *list, double *dist, double din, int count, int Np)
{ {
@@ -2652,11 +2689,23 @@ extern "C" double deviceReduce(double *in, double* out, int N) {
return sum; return sum;
} }
// extern "C" void ScaLBL_D3Q19_Reflection_BC_z(int *list, double *dist, int count, int Np){
//extern "C" void ScaLBL_D3Q19_Pressure_BC_Z(int *list, double *dist, double dout, int count, int Np){ int GRID = count / 512 + 1;
// int GRID = count / 512 + 1; dvc_ScaLBL_D3Q19_Reflection_BC_z<<<GRID,512>>>(neighborList, list, dist, count, N);
// dvc_ScaLBL_D3Q19_Pressure_BC_Z<<<GRID,512>>>(disteven, distodd, dout, Nx, Ny, Nz, outlet); cudaError_t err = cudaGetLastError();
//} if (cudaSuccess != err){
printf("CUDA error in ScaLBL_D3Q19_Reflection_BC_z (kernel): %s \n",cudaGetErrorString(err));
}
}
extern "C" void ScaLBL_D3Q19_Reflection_BC_Z(int *list, double *dist, int count, int Np){
int GRID = count / 512 + 1;
dvc_ScaLBL_D3Q19_Reflection_BC_Z<<<GRID,512>>>(neighborList, list, dist, count, N);
cudaError_t err = cudaGetLastError();
if (cudaSuccess != err){
printf("CUDA error in ScaLBL_D3Q19_Reflection_BC_Z (kernel): %s \n",cudaGetErrorString(err));
}
}
extern "C" void ScaLBL_D3Q19_AAeven_MRT(double *dist, int start, int finish, int Np, double rlx_setA, double rlx_setB, double Fx, extern "C" void ScaLBL_D3Q19_AAeven_MRT(double *dist, int start, int finish, int Np, double rlx_setA, double rlx_setB, double Fx,
double Fy, double Fz){ double Fy, double Fz){