//************************************************************************* // Lattice Boltzmann Simulator for Single Phase Flow in Porous Media // James E. McCLure //************************************************************************* #include #include #include #include "common/ScaLBL.h" #include "common/MPI_Helpers.h" using namespace std; extern void AA1_ScaLBL_D3Q19_Init(double *dist_even, double *dist_odd, int Nx, int Ny, int Nz, int iproc, int jproc, int kproc, int nprocx, int nprocy, int nprocz) { // Set of Discrete velocities for the D3Q19 Model static int D3Q19[18][3]={{1,0,0},{-1,0,0},{0,1,0},{0,-1,0},{0,0,1},{0,0,-1}, {1,1,0},{-1,-1,0},{1,-1,0},{-1,1,0}, {1,0,1},{-1,0,-1},{1,0,-1},{-1,0,1}, {0,1,1},{0,-1,-1},{0,1,-1},{0,-1,1}}; int q,i,j,k,n,N; int Cqx,Cqy,Cqz; // Discrete velocity int x,y,z; // Global indices int xn,yn,zn; // Global indices of neighbor int X,Y,Z; // Global size X = Nx*nprocx; Y = Ny*nprocy; Z = Nz*nprocz; NULL_USE(Z); N = (Nx+2)*(Ny+2)*(Nz+2); // size of the array including halo for (k=0; k ",f_even_host[(q+1)*Np + idx],idx); printf("%.02f |",f_odd_host[(q*Np+idx)]); } } } } } printf("\n"); } extern void HostToUnobtainium(IntArray Map, double * f_even_host, double * f_odd_host, double * f_even_unobtainium, double * f_odd_unobtainium, int Nx, int Ny, int Nz,int Np,int N) { int n; for (int k=1;k> nprocx; domain >> nprocy; domain >> nprocz; domain >> Nx; domain >> Ny; domain >> Nz; domain >> nspheres; domain >> Lx; domain >> Ly; domain >> Lz; } else if (nprocs==1){ nprocx=nprocy=nprocz=1; Nx=3; Ny = 1; Nz = 1; nspheres=0; Lx=Ly=Lz=1; } else if (nprocs==2){ nprocx=2; nprocy=1; nprocz=1; Nx=Ny=Nz=dim; Nx = dim; Ny = dim; Nz = dim; nspheres=0; Lx=Ly=Lz=1; } else if (nprocs==4){ nprocx=nprocy=2; nprocz=1; Nx=Ny=Nz=dim; nspheres=0; Lx=Ly=Lz=1; } else if (nprocs==8){ nprocx=nprocy=nprocz=2; Nx=Ny=Nz=dim; nspheres=0; Lx=Ly=Lz=1; } //....................................................................... } // ************************************************************** // Broadcast simulation parameters from rank 0 to all other procs MPI_Barrier(comm); //................................................. MPI_Bcast(&Nx,1,MPI_INT,0,comm); MPI_Bcast(&Ny,1,MPI_INT,0,comm); MPI_Bcast(&Nz,1,MPI_INT,0,comm); MPI_Bcast(&nprocx,1,MPI_INT,0,comm); MPI_Bcast(&nprocy,1,MPI_INT,0,comm); MPI_Bcast(&nprocz,1,MPI_INT,0,comm); MPI_Bcast(&nspheres,1,MPI_INT,0,comm); MPI_Bcast(&Lx,1,MPI_DOUBLE,0,comm); MPI_Bcast(&Ly,1,MPI_DOUBLE,0,comm); MPI_Bcast(&Lz,1,MPI_DOUBLE,0,comm); //................................................. MPI_Barrier(comm); // ************************************************************** // ************************************************************** if (nprocs != nprocx*nprocy*nprocz){ printf("nprocx = %i \n",nprocx); printf("nprocy = %i \n",nprocy); printf("nprocz = %i \n",nprocz); INSIST(nprocs == nprocx*nprocy*nprocz,"Fatal error in processor count!"); } if (rank==0){ printf("********************************************************\n"); printf("Sub-domain size = %i x %i x %i\n",Nx,Ny,Nz); printf("Process grid = %i x %i x %i\n",nprocx,nprocy,nprocz); printf("********************************************************\n"); } MPI_Barrier(comm); kproc = rank/(nprocx*nprocy); jproc = (rank-nprocx*nprocy*kproc)/nprocx; iproc = rank-nprocx*nprocy*kproc-nprocz*jproc; if (rank == 0) { printf("i,j,k proc=%d %d %d \n",iproc,jproc,kproc); } MPI_Barrier(comm); if (rank == 1){ printf("i,j,k proc=%d %d %d \n",iproc,jproc,kproc); printf("\n\n"); } double iVol_global = 1.0/Nx/Ny/Nz/nprocx/nprocy/nprocz; int BoundaryCondition=0; Domain Dm(Nx,Ny,Nz,rank,nprocx,nprocy,nprocz,Lx,Ly,Lz,BoundaryCondition); Nx += 2; Ny += 2; Nz += 2; int N = Nx*Ny*Nz; //....................................................................... // Assign the phase ID field //....................................................................... char LocalRankString[8]; sprintf(LocalRankString,"%05d",rank); char LocalRankFilename[40]; sprintf(LocalRankFilename,"ID.%05i",rank); FILE *IDFILE = fopen(LocalRankFilename,"rb"); if (IDFILE==NULL) ERROR("Error opening file: ID.xxxxx"); fread(Dm.id,1,N,IDFILE); fclose(IDFILE); MPI_Barrier(comm); Dm.CommInit(); //....................................................................... // Compute the media porosity //....................................................................... double sum; double sum_local=0.0, porosity; int Np=0; // number of local pore nodes //....................................................................... for (k=1;k 0){ sum_local+=1.0; Np++; } } } } MPI_Barrier(comm); MPI_Allreduce(&sum_local,&sum,1,MPI_DOUBLE,MPI_SUM,comm); porosity = sum*iVol_global; if (rank==0) printf("Media porosity = %f \n",porosity); MPI_Barrier(comm); if (rank == 0) cout << "Domain set." << endl; if (rank==0) printf ("Create ScaLBL_Communicator \n"); // Create a communicator for the device ScaLBL_Communicator ScaLBL_Comm(Dm); //...........device phase ID................................................. if (rank==0) printf ("Copying phase ID to device \n"); char *ID; ScaLBL_AllocateDeviceMemory((void **) &ID, N); // Allocate device memory // Copy to the device ScaLBL_CopyToDevice(ID, Dm.id, N); //........................................................................... if (rank==0){ printf("Total domain size = %i \n",N); printf("Reduced domain size = %i \n",Np); } // LBM variables if (rank==0) printf ("Set up the neighborlist \n"); int neighborSize=18*Np*sizeof(int); int *neighborList; IntArray Map(Nx,Ny,Nz); neighborList= new int[18*Np]; ScaLBL_Comm.MemoryOptimizedLayoutAA(Map,neighborList,Dm.id,Np,1); MPI_Barrier(comm); //......................device distributions................................. int dist_mem_size = Np*sizeof(double); if (rank==0) printf ("Allocating distributions \n"); int *NeighborList; // double *f_even,*f_odd; double * dist; double * Velocity; //........................................................................... ScaLBL_AllocateDeviceMemory((void **) &dist, 19*dist_mem_size); ScaLBL_AllocateDeviceMemory((void **) &NeighborList, neighborSize); ScaLBL_AllocateDeviceMemory((void **) &Velocity, 3*sizeof(double)*Np); ScaLBL_CopyToDevice(NeighborList, neighborList, neighborSize); //........................................................................... ScaLBL_D3Q19_Init(dist, Np); /************ MAIN ITERATION LOOP (timing communications)***************************************/ if (rank==0) printf("Beginning AA timesteps...\n"); if (rank==0) printf("********************************************************\n"); if (rank==0) printf("No. of timesteps for timing: %i \n", timesteps); //.......create and start timer............ double starttime,stoptime,cputime; ScaLBL_DeviceBarrier(); MPI_Barrier(comm); starttime = MPI_Wtime(); while (timestep < timesteps) { ScaLBL_Comm.SendD3Q19AA(dist); //READ FROM NORMAL ScaLBL_D3Q19_AAodd_MRT(NeighborList, dist, ScaLBL_Comm.next, Np, Np, rlx_setA, rlx_setB, Fx, Fy, Fz); ScaLBL_Comm.RecvD3Q19AA(dist); //WRITE INTO OPPOSITE ScaLBL_D3Q19_AAodd_MRT(NeighborList, dist, 0, ScaLBL_Comm.next, Np, rlx_setA, rlx_setB, Fx, Fy, Fz); ScaLBL_DeviceBarrier(); MPI_Barrier(comm); timestep++; ScaLBL_Comm.SendD3Q19AA(dist); //READ FORM NORMAL ScaLBL_D3Q19_AAeven_MRT(dist, ScaLBL_Comm.next, Np, Np, rlx_setA, rlx_setB, Fx, Fy, Fz); ScaLBL_Comm.RecvD3Q19AA(dist); //WRITE INTO OPPOSITE ScaLBL_D3Q19_AAeven_MRT(dist, 0, ScaLBL_Comm.next, Np, rlx_setA, rlx_setB, Fx, Fy, Fz); ScaLBL_DeviceBarrier(); MPI_Barrier(comm); timestep++; //************************************************************************/ } //************************************************************************/ stoptime = MPI_Wtime(); // cout << "CPU time: " << (stoptime - starttime) << " seconds" << endl; cputime = stoptime - starttime; // cout << "Lattice update rate: "<< double(Nx*Ny*Nz*timestep)/cputime/1000000 << " MLUPS" << endl; double MLUPS = double(Np*timestep)/cputime/1000000; if (rank==0) printf("********************************************************\n"); if (rank==0) printf("CPU time = %f \n", cputime); if (rank==0) printf("Lattice update rate (per process)= %f MLUPS \n", MLUPS); MLUPS *= nprocs; if (rank==0) printf("Lattice update rate (process)= %f MLUPS \n", MLUPS); if (rank==0) printf("********************************************************\n"); // Number of memory references from the swap algorithm (per timestep) // 18 reads and 18 writes for each lattice site double MemoryRefs = Np*38; // number of memory references for the swap algorithm - GigaBytes / second if (rank==0) printf("DRAM bandwidth (per process)= %f GB/sec \n",MemoryRefs*8*timestep/1e9/cputime); // Report bandwidth in Gigabits per second // communication bandwidth includes both send and recieve if (rank==0) printf("Communication bandwidth (per process)= %f Gbit/sec \n",ScaLBL_Comm.CommunicationCount*64*timestep/1e9/cputime); if (rank==0) printf("Aggregated communication bandwidth = %f Gbit/sec \n",nprocs*ScaLBL_Comm.CommunicationCount*64*timestep/1e9/cputime); double *VEL; VEL= new double [3*Np]; int SIZE=3*Np*sizeof(double); ScaLBL_D3Q19_Momentum(dist,Velocity, Np); ScaLBL_DeviceBarrier(); MPI_Barrier(comm); ScaLBL_CopyToHost(&VEL[0],&Velocity[0],SIZE); sum_local=0.f; sum = 0.f; for (k=1;k 0){\ int idx = Map(i,j,k); sum_local+=VEL[2*Np+idx]; } } } } MPI_Allreduce(&sum_local,&sum,1,MPI_DOUBLE,MPI_SUM,comm); double PoreVel = sum*iVol_global; if (rank==0) printf("Velocity = %f \n",PoreVel); } // **************************************************** comm.barrier(); Utilities::shutdown(); // **************************************************** return check; }