//************************************************************************* // Lattice Boltzmann Simulator for Single Phase Flow in Porous Media // James E. McCLure //************************************************************************* #include #include #include #include #include "common/MPI.h" #include using namespace std; //************************************************************************* // MRT implementation of the LBM using CUDA //************************************************************************* //************************************************************************* //************************************************************************* extern "C" void dvc_AllocateDeviceMemory(void** address, size_t size); //************************************************************************* extern "C" void dvc_CopyToDevice(void* dest, void* source, size_t size); //************************************************************************* extern "C" void dvc_Barrier(); //************************************************************************* extern "C" void dvc_InitD3Q19(int nblocks, int nthreads, int S, char *ID, double *f_even, double *f_odd, int Nx, int Ny, int Nz); //************************************************************************* extern "C" void dvc_SwapD3Q19( int nblocks, int nthreads, int S, char *ID, double *f_even, double *f_odd, int Nx, int Ny, int Nz); //************************************************************************* extern "C" void dvc_MRT( int nblocks, int nthreads, int S, char *ID, double *f_even, double *f_odd, double rlxA, double rlxB, double Fx, double Fy, double Fz, int Nx, int Ny, int Nz); //************************************************************************* extern "C" void dvc_PackDist(int grid, int threads, int q, int *SendList, int start, int sendCount, double *sendbuf, double *Dist, int N); //************************************************************************* extern "C" void dvc_UnpackDist(int grid, int threads, int q, int Cqx, int Cqy, int Cqz, int *RecvList, int start, int recvCount, double *recvbuf, double *Dist, int Nx, int Ny, int Nz); //************************************************************************* inline void PackID(int *list, int count, char *sendbuf, char *ID){ // Fill in the phase ID values from neighboring processors // This packs up the values that need to be sent from one processor to another int idx,n; for (idx=0; idx> FILENAME; // name of the input file input >> Nz; // number of nodes (x,y,z) input >> nBlocks; input >> nthreads; input >> tau; // relaxation time input >> Fx; // External force components (x,y,z) input >> Fy; input >> Fz; input >> timestepMax; // max no. of timesteps input >> interval; // error interval input >> tol; // error tolerance ifstream domain("Domain.in"); domain >> nprocx; domain >> nprocy; domain >> nprocz; } // ************************************************************** // Broadcast simulation parameters from rank 0 to all other procs comm.barrier(); //................................................. comm.bcast(&Nz,1,0); comm.bcast(&nBlocks,1,0); comm.bcast(&nthreads,1,0); comm.bcast(&tau,1,0); comm.bcast(&Fx,1,0); comm.bcast(&Fy,1,0); comm.bcast(&Fz,1,0); comm.bcast(×tepMax,1,0); comm.bcast(&interval,1,0); comm.bcast(&tol,1,0); comm.bcast(&nprocx,1,0); comm.bcast(&nprocy,1,0); comm.bcast(&nprocz,1,0); //................................................. comm.barrier(); // ************************************************************** double rlx_setA = 1.f/tau; double rlx_setB = 8.f*(2.f-rlx_setA)/(8.f-rlx_setA); if (nprocs != nprocx*nprocy*nprocz){ printf("Fatal error in processor number! \n"); printf("nprocx = %i \n",nprocx); printf("nprocy = %i \n",nprocy); printf("nprocz = %i \n",nprocz); } if (rank==0){ printf("tau = %f \n", tau); printf("Set A = %f \n", rlx_setA); printf("Set B = %f \n", rlx_setB); printf("Force(x) = %f \n", Fx); printf("Force(y) = %f \n", Fy); printf("Force(z) = %f \n", Fz); printf("Sub-domain size = %i x %i x %i\n",Nz,Nz,Nz); } comm.barrier(); kproc = rank/(nprocx*nprocy); jproc = (rank-nprocx*nprocy*kproc)/nprocx; iproc = rank-nprocx*nprocy*kproc-nprocz*jproc; //.......................................... // set up the neighbor ranks //.......................................... i=iproc; j=jproc; k =kproc; i+=1; j+=0; k+=0; if (i<0) i+=nprocx; if (j<0) j+=nprocy; if (k<0) k+=nprocz; if (!(i 0) sum++; } } } PM.close(); // printf("File porosity = %f\n", double(sum)/N); //........................................................................... comm.barrier(); if (rank == 0) cout << "Domain set." << endl; //........................................................................... // Write the communcation structure into a file for debugging char LocalCommFile[40]; sprintf(LocalCommFile,"%s%s","Comm.",LocalRankString); FILE *CommFile; CommFile = fopen(LocalCommFile,"w"); fprintf(CommFile,"rank=%d, ",rank); fprintf(CommFile,"i=%d,j=%d,k=%d :",iproc,jproc,kproc); fprintf(CommFile,"x=%d, ",rank_x); fprintf(CommFile,"X=%d, ",rank_X); fprintf(CommFile,"y=%d, ",rank_y); fprintf(CommFile,"Y=%d, ",rank_Y); fprintf(CommFile,"z=%d, ",rank_z); fprintf(CommFile,"Z=%d, ",rank_Z); fprintf(CommFile,"xy=%d, ",rank_xy); fprintf(CommFile,"XY=%d, ",rank_XY); fprintf(CommFile,"xY=%d, ",rank_xY); fprintf(CommFile,"Xy=%d, ",rank_Xy); fprintf(CommFile,"xz=%d, ",rank_xz); fprintf(CommFile,"XZ=%d, ",rank_XZ); fprintf(CommFile,"xZ=%d, ",rank_xZ); fprintf(CommFile,"Xz=%d, ",rank_Xz); fprintf(CommFile,"yz=%d, ",rank_yz); fprintf(CommFile,"YZ=%d, ",rank_YZ); fprintf(CommFile,"yZ=%d, ",rank_yZ); fprintf(CommFile,"Yz=%d, ",rank_Yz); fprintf(CommFile,"\n"); fclose(CommFile); //........................................................................... // Set up MPI communication structures if (rank==0) printf ("Setting up communication control structures \n"); //...................................................................................... // Get the actual D3Q19 communication counts (based on location of solid phase) // Discrete velocity set symmetry implies the sendcount = recvcount int sendCount_x, sendCount_y, sendCount_z, sendCount_X, sendCount_Y, sendCount_Z; int sendCount_xy, sendCount_yz, sendCount_xz, sendCount_Xy, sendCount_Yz, sendCount_xZ; int sendCount_xY, sendCount_yZ, sendCount_Xz, sendCount_XY, sendCount_YZ, sendCount_XZ; sendCount_x = sendCount_y = sendCount_z = sendCount_X = sendCount_Y = sendCount_Z = 0; sendCount_xy = sendCount_yz = sendCount_xz = sendCount_Xy = sendCount_Yz = sendCount_xZ = 0; sendCount_xY = sendCount_yZ = sendCount_Xz = sendCount_XY = sendCount_YZ = sendCount_XZ = 0; //...................................................................................... for (k=0; k