2018-06-11 15:19:05 -04:00
|
|
|
/*
|
|
|
|
|
Copyright 2013--2018 James E. McClure, Virginia Polytechnic & State University
|
|
|
|
|
|
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
OPM 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.
|
|
|
|
|
OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
|
*/
|
2014-03-19 09:08:39 -04:00
|
|
|
// GPU Functions for D3Q7 Lattice Boltzmann Methods
|
2020-04-07 10:43:01 -04:00
|
|
|
#include <stdio.h>
|
2014-03-18 11:55:10 -04:00
|
|
|
|
2018-01-24 10:08:43 -05:00
|
|
|
#define NBLOCKS 560
|
2014-03-18 11:55:10 -04:00
|
|
|
#define NTHREADS 128
|
2013-08-26 15:12:25 -04:00
|
|
|
|
2016-11-24 11:26:51 -05:00
|
|
|
__global__ void dvc_ScaLBL_Scalar_Pack(int *list, int count, double *sendbuf, double *Data, int N){
|
2013-08-26 15:12:25 -04:00
|
|
|
//....................................................................................
|
|
|
|
|
// Pack distribution q into the send buffer for the listed lattice sites
|
|
|
|
|
// dist may be even or odd distributions stored by stream layout
|
|
|
|
|
//....................................................................................
|
|
|
|
|
int idx,n;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx<count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
sendbuf[idx] = Data[n];
|
|
|
|
|
}
|
|
|
|
|
}
|
2016-11-23 17:03:12 -05:00
|
|
|
__global__ void dvc_ScaLBL_Scalar_Unpack(int *list, int count, double *recvbuf, double *Data, int N){
|
2013-08-26 15:12:25 -04:00
|
|
|
//....................................................................................
|
|
|
|
|
// Pack distribution q into the send buffer for the listed lattice sites
|
|
|
|
|
// dist may be even or odd distributions stored by stream layout
|
|
|
|
|
//....................................................................................
|
|
|
|
|
int idx,n;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx<count){
|
|
|
|
|
n = list[idx];
|
|
|
|
|
Data[n] = recvbuf[idx];
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2016-11-23 17:03:12 -05:00
|
|
|
__global__ void dvc_ScaLBL_PackDenD3Q7(int *list, int count, double *sendbuf, int number, double *Data, int N){
|
2013-08-26 15:12:25 -04:00
|
|
|
//....................................................................................
|
|
|
|
|
// Pack distribution into the send buffer for the listed lattice sites
|
|
|
|
|
//....................................................................................
|
|
|
|
|
int idx,n,component;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx<count){
|
2014-03-18 11:55:10 -04:00
|
|
|
for (component=0; component<number; component++){
|
2013-08-26 15:12:25 -04:00
|
|
|
n = list[idx];
|
|
|
|
|
sendbuf[idx*number+component] = Data[number*n+component];
|
|
|
|
|
Data[number*n+component] = 0.0; // Set the data value to zero once it's in the buffer!
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
2016-11-23 17:03:12 -05:00
|
|
|
__global__ void dvc_ScaLBL_UnpackDenD3Q7(int *list, int count, double *recvbuf, int number, double *Data, int N){
|
2013-08-26 15:12:25 -04:00
|
|
|
//....................................................................................
|
|
|
|
|
// Unack distribution from the recv buffer
|
|
|
|
|
// Sum to the existing density value
|
|
|
|
|
//....................................................................................
|
|
|
|
|
int idx,n,component;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx<count){
|
2014-03-17 10:14:46 -04:00
|
|
|
for (component=0; component<number; component++){
|
2013-08-26 15:12:25 -04:00
|
|
|
n = list[idx];
|
|
|
|
|
Data[number*n+component] += recvbuf[idx*number+component];
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2018-01-24 10:08:43 -05:00
|
|
|
__global__ void dvc_ScaLBL_D3Q7_Unpack(int q, int *list, int start, int count,
|
|
|
|
|
double *recvbuf, double *dist, int N){
|
|
|
|
|
//....................................................................................
|
|
|
|
|
// Unpack distribution from the recv buffer
|
|
|
|
|
// Distribution q matche Cqx, Cqy, Cqz
|
|
|
|
|
// swap rule means that the distributions in recvbuf are OPPOSITE of q
|
|
|
|
|
// dist may be even or odd distributions stored by stream layout
|
|
|
|
|
//....................................................................................
|
|
|
|
|
int n,idx;
|
|
|
|
|
idx = blockIdx.x*blockDim.x + threadIdx.x;
|
|
|
|
|
if (idx<count){
|
|
|
|
|
// Get the value from the list -- note that n is the index is from the send (non-local) process
|
|
|
|
|
n = list[idx];
|
|
|
|
|
// unpack the distribution to the proper location
|
|
|
|
|
if (!(n<0)) { dist[q*N+n] = recvbuf[start+idx];
|
|
|
|
|
//printf("%f \n",,dist[q*N+n]);
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2020-04-07 10:38:21 -04:00
|
|
|
__global__ void dvc_ScaLBL_D3Q7_Reflection_BC_z(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.222222222222222222222222 - dist[6*Np+n];
|
|
|
|
|
dist[6*Np+n] = f5;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
__global__ void dvc_ScaLBL_D3Q7_Reflection_BC_Z(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.222222222222222222222222 - dist[5*Np+n];
|
|
|
|
|
dist[5*Np+n] = f6;
|
|
|
|
|
}
|
|
|
|
|
}
|
2016-11-23 17:03:12 -05:00
|
|
|
__global__ void dvc_ScaLBL_D3Q7_Init(char *ID, double *f_even, double *f_odd, double *Den, int Nx, int Ny, int Nz)
|
2013-08-26 15:12:25 -04:00
|
|
|
{
|
2014-03-17 10:14:46 -04:00
|
|
|
int n,N;
|
|
|
|
|
N = Nx*Ny*Nz;
|
|
|
|
|
double value;
|
2015-06-30 15:35:55 -04:00
|
|
|
char id;
|
2014-03-18 11:55:10 -04:00
|
|
|
int S = N/NBLOCKS/NTHREADS + 1;
|
|
|
|
|
for (int s=0; s<S; s++){
|
2014-03-17 10:14:46 -04:00
|
|
|
//........Get 1-D index for this thread....................
|
|
|
|
|
n = S*blockIdx.x*blockDim.x + s*blockDim.x + threadIdx.x;
|
|
|
|
|
if (n<N){
|
2015-06-30 15:38:31 -04:00
|
|
|
id = ID[n];
|
2015-06-30 15:35:55 -04:00
|
|
|
if (id > 0){
|
2014-03-18 11:55:10 -04:00
|
|
|
value = Den[n];
|
|
|
|
|
f_even[n] = 0.3333333333333333*value;
|
|
|
|
|
f_odd[n] = 0.1111111111111111*value; //double(100*n)+1.f;
|
|
|
|
|
f_even[N+n] = 0.1111111111111111*value; //double(100*n)+2.f;
|
|
|
|
|
f_odd[N+n] = 0.1111111111111111*value; //double(100*n)+3.f;
|
|
|
|
|
f_even[2*N+n] = 0.1111111111111111*value; //double(100*n)+4.f;
|
|
|
|
|
f_odd[2*N+n] = 0.1111111111111111*value; //double(100*n)+5.f;
|
|
|
|
|
f_even[3*N+n] = 0.1111111111111111*value; //double(100*n)+6.f;
|
|
|
|
|
}
|
|
|
|
|
else{
|
2015-07-04 09:40:48 -04:00
|
|
|
for(int q=0; q<3; q++){
|
2014-03-18 11:55:10 -04:00
|
|
|
f_even[q*N+n] = -1.0;
|
|
|
|
|
f_odd[q*N+n] = -1.0;
|
|
|
|
|
}
|
|
|
|
|
f_even[3*N+n] = -1.0;
|
2014-03-17 10:14:46 -04:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
}
|
2013-08-26 15:12:25 -04:00
|
|
|
}
|
2014-03-17 10:14:46 -04:00
|
|
|
|
|
|
|
|
//*************************************************************************
|
2016-11-23 17:03:12 -05:00
|
|
|
__global__ void dvc_ScaLBL_D3Q7_Swap(char *ID, double *disteven, double *distodd, int Nx, int Ny, int Nz)
|
2013-08-26 15:12:25 -04:00
|
|
|
{
|
2014-03-17 10:14:46 -04:00
|
|
|
int i,j,k,n,nn,N;
|
|
|
|
|
// distributions
|
2014-03-18 11:55:10 -04:00
|
|
|
double f1,f2,f3,f4,f5,f6;
|
2015-06-26 12:14:57 -04:00
|
|
|
char id;
|
2014-03-17 10:14:46 -04:00
|
|
|
N = Nx*Ny*Nz;
|
2014-03-18 11:55:10 -04:00
|
|
|
|
|
|
|
|
int S = N/NBLOCKS/NTHREADS + 1;
|
|
|
|
|
for (int s=0; s<S; s++){
|
|
|
|
|
//........Get 1-D index for this thread....................
|
|
|
|
|
n = S*blockIdx.x*blockDim.x + s*blockDim.x + threadIdx.x;
|
|
|
|
|
|
2015-06-26 12:14:57 -04:00
|
|
|
if (n<N ){
|
|
|
|
|
id = ID[n];
|
2015-06-30 15:16:41 -04:00
|
|
|
if (id > 0){
|
2015-07-04 16:34:39 -04:00
|
|
|
//.......Back out the 3-D indices for node n..............
|
|
|
|
|
k = n/(Nx*Ny);
|
|
|
|
|
j = (n-Nx*Ny*k)/Nx;
|
2016-07-25 11:09:05 -04:00
|
|
|
i = n-Nx*Ny*k-Nx*j;
|
2015-07-04 16:34:39 -04:00
|
|
|
//........................................................................
|
|
|
|
|
// Retrieve even distributions from the local node (swap convention)
|
|
|
|
|
// f0 = disteven[n]; // Does not particupate in streaming
|
|
|
|
|
f1 = distodd[n];
|
|
|
|
|
f3 = distodd[N+n];
|
|
|
|
|
f5 = distodd[2*N+n];
|
|
|
|
|
//........................................................................
|
2014-03-18 11:55:10 -04:00
|
|
|
|
2015-07-04 16:34:39 -04:00
|
|
|
//........................................................................
|
|
|
|
|
// Retrieve odd distributions from neighboring nodes (swap convention)
|
|
|
|
|
//........................................................................
|
|
|
|
|
nn = n+1; // neighbor index (pull convention)
|
|
|
|
|
if (!(i+1<Nx)) nn -= Nx; // periodic BC along the x-boundary
|
|
|
|
|
//if (i+1<Nx){
|
|
|
|
|
f2 = disteven[N+nn]; // pull neighbor for distribution 2
|
|
|
|
|
if (!(f2 < 0.0)){
|
|
|
|
|
distodd[n] = f2;
|
|
|
|
|
disteven[N+nn] = f1;
|
|
|
|
|
}
|
|
|
|
|
//}
|
|
|
|
|
//........................................................................
|
|
|
|
|
nn = n+Nx; // neighbor index (pull convention)
|
|
|
|
|
if (!(j+1<Ny)) nn -= Nx*Ny; // Perioidic BC along the y-boundary
|
|
|
|
|
//if (j+1<Ny){
|
|
|
|
|
f4 = disteven[2*N+nn]; // pull neighbor for distribution 4
|
|
|
|
|
if (!(f4 < 0.0)){
|
|
|
|
|
distodd[N+n] = f4;
|
|
|
|
|
disteven[2*N+nn] = f3;
|
|
|
|
|
}
|
|
|
|
|
//........................................................................
|
|
|
|
|
nn = n+Nx*Ny; // neighbor index (pull convention)
|
|
|
|
|
if (!(k+1<Nz)) nn -= Nx*Ny*Nz; // Perioidic BC along the z-boundary
|
|
|
|
|
//if (k+1<Nz){
|
|
|
|
|
f6 = disteven[3*N+nn]; // pull neighbor for distribution 6
|
|
|
|
|
if (!(f6 < 0.0)){
|
|
|
|
|
distodd[2*N+n] = f6;
|
|
|
|
|
disteven[3*N+nn] = f5;
|
|
|
|
|
}
|
2015-06-26 12:14:57 -04:00
|
|
|
}
|
2014-03-17 10:14:46 -04:00
|
|
|
}
|
|
|
|
|
}
|
2013-08-26 15:12:25 -04:00
|
|
|
}
|
2014-03-17 10:14:46 -04:00
|
|
|
|
|
|
|
|
//*************************************************************************
|
2016-11-23 17:03:12 -05:00
|
|
|
__global__ void dvc_ScaLBL_D3Q7_Density(char *ID, double *disteven, double *distodd, double *Den,
|
2014-03-18 11:55:10 -04:00
|
|
|
int Nx, int Ny, int Nz)
|
2013-08-26 15:12:25 -04:00
|
|
|
{
|
2014-03-17 10:14:46 -04:00
|
|
|
char id;
|
|
|
|
|
int n;
|
|
|
|
|
double f0,f1,f2,f3,f4,f5,f6;
|
|
|
|
|
int N = Nx*Ny*Nz;
|
2014-03-18 11:55:10 -04:00
|
|
|
|
|
|
|
|
int S = N/NBLOCKS/NTHREADS + 1;
|
|
|
|
|
for (int s=0; s<S; s++){
|
|
|
|
|
//........Get 1-D index for this thread....................
|
|
|
|
|
n = S*blockIdx.x*blockDim.x + s*blockDim.x + threadIdx.x;
|
2015-06-24 20:02:16 -04:00
|
|
|
if (n<N){
|
2015-07-04 16:34:39 -04:00
|
|
|
id = ID[n];
|
|
|
|
|
if (id > 0 ){
|
|
|
|
|
// Read the distributions
|
|
|
|
|
f0 = disteven[n];
|
|
|
|
|
f2 = disteven[N+n];
|
|
|
|
|
f4 = disteven[2*N+n];
|
|
|
|
|
f6 = disteven[3*N+n];
|
|
|
|
|
f1 = distodd[n];
|
|
|
|
|
f3 = distodd[N+n];
|
|
|
|
|
f5 = distodd[2*N+n];
|
|
|
|
|
// Compute the density
|
|
|
|
|
Den[n] = f0+f1+f2+f3+f4+f5+f6;
|
|
|
|
|
}
|
2014-03-17 10:14:46 -04:00
|
|
|
}
|
|
|
|
|
}
|
2013-08-26 15:12:25 -04:00
|
|
|
}
|
2014-03-18 11:55:10 -04:00
|
|
|
|
2020-04-07 10:38:21 -04:00
|
|
|
extern "C" void ScaLBL_D3Q7_Reflection_BC_z(int *list, double *dist, int count, int Np){
|
|
|
|
|
int GRID = count / 512 + 1;
|
|
|
|
|
dvc_ScaLBL_D3Q7_Reflection_BC_z<<<GRID,512>>>(list, dist, count, Np);
|
|
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
|
|
|
|
printf("CUDA error in ScaLBL_D3Q7_Reflection_BC_z (kernel): %s \n",cudaGetErrorString(err));
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
extern "C" void ScaLBL_D3Q7_Reflection_BC_Z(int *list, double *dist, int count, int Np){
|
|
|
|
|
int GRID = count / 512 + 1;
|
|
|
|
|
dvc_ScaLBL_D3Q7_Reflection_BC_Z<<<GRID,512>>>(list, dist, count, Np);
|
|
|
|
|
cudaError_t err = cudaGetLastError();
|
|
|
|
|
if (cudaSuccess != err){
|
|
|
|
|
printf("CUDA error in ScaLBL_D3Q7_Reflection_BC_Z (kernel): %s \n",cudaGetErrorString(err));
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2018-01-24 10:08:43 -05:00
|
|
|
extern "C" void ScaLBL_D3Q7_Unpack(int q, int *list, int start, int count, double *recvbuf, double *dist, int N){
|
|
|
|
|
int GRID = count / 512 + 1;
|
|
|
|
|
dvc_ScaLBL_D3Q7_Unpack <<<GRID,512 >>>(q, list, start, count, recvbuf, dist, N);
|
|
|
|
|
}
|
|
|
|
|
|
2016-11-25 15:51:43 -05:00
|
|
|
extern "C" void ScaLBL_Scalar_Pack(int *list, int count, double *sendbuf, double *Data, int N){
|
2014-03-18 11:55:10 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2016-11-25 15:51:43 -05:00
|
|
|
dvc_ScaLBL_Scalar_Pack <<<GRID,512 >>>(list, count, sendbuf, Data, N);
|
2014-03-18 11:55:10 -04:00
|
|
|
}
|
|
|
|
|
|
2016-11-23 17:03:12 -05:00
|
|
|
extern "C" void ScaLBL_Scalar_Unpack(int *list, int count, double *recvbuf, double *Data, int N){
|
2014-03-18 11:55:10 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2016-11-23 17:03:12 -05:00
|
|
|
dvc_ScaLBL_Scalar_Unpack <<<GRID,512 >>>(list, count, recvbuf, Data, N);
|
2014-03-18 11:55:10 -04:00
|
|
|
}
|
2016-11-23 17:03:12 -05:00
|
|
|
extern "C" void ScaLBL_PackDenD3Q7(int *list, int count, double *sendbuf, int number, double *Data, int N){
|
2014-03-18 11:55:10 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2016-11-23 17:03:12 -05:00
|
|
|
dvc_ScaLBL_PackDenD3Q7 <<<GRID,512 >>>(list, count, sendbuf, number, Data, N);
|
2014-03-18 11:55:10 -04:00
|
|
|
}
|
|
|
|
|
|
2016-11-23 17:03:12 -05:00
|
|
|
extern "C" void ScaLBL_UnpackDenD3Q7(int *list, int count, double *recvbuf, int number, double *Data, int N){
|
2014-03-18 11:55:10 -04:00
|
|
|
int GRID = count / 512 + 1;
|
2016-11-23 17:03:12 -05:00
|
|
|
dvc_ScaLBL_UnpackDenD3Q7 <<<GRID,512 >>>(list, count, recvbuf, number, Data, N);
|
2014-03-18 11:55:10 -04:00
|
|
|
}
|
|
|
|
|
|
2016-11-23 17:03:12 -05:00
|
|
|
extern "C" void ScaLBL_D3Q7_Init(char *ID, double *f_even, double *f_odd, double *Den, int Nx, int Ny, int Nz){
|
|
|
|
|
dvc_ScaLBL_D3Q7_Init <<<NBLOCKS,NTHREADS >>>(ID, f_even, f_odd, Den, Nx, Ny, Nz);
|
2014-03-18 11:55:10 -04:00
|
|
|
}
|
|
|
|
|
|
2016-11-23 17:03:12 -05:00
|
|
|
extern "C" void ScaLBL_D3Q7_Swap(char *ID, double *disteven, double *distodd, int Nx, int Ny, int Nz){
|
|
|
|
|
dvc_ScaLBL_D3Q7_Swap <<<NBLOCKS,NTHREADS >>>(ID, disteven, distodd, Nx, Ny, Nz);
|
2014-03-18 11:55:10 -04:00
|
|
|
}
|
|
|
|
|
|
2016-11-23 17:03:12 -05:00
|
|
|
extern "C" void ScaLBL_D3Q7_Density(char *ID, double *disteven, double *distodd, double *Den,
|
2014-03-18 11:55:10 -04:00
|
|
|
int Nx, int Ny, int Nz){
|
2016-11-23 17:03:12 -05:00
|
|
|
dvc_ScaLBL_D3Q7_Density <<<NBLOCKS,NTHREADS >>>(ID, disteven, distodd, Den, Nx, Ny, Nz);
|
2014-03-18 11:55:10 -04:00
|
|
|
}
|
|
|
|
|
|