//************************************************************************* // 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; inline void AssignComponentLabels(char *id, double *phase, int Nx, int Ny, int Nz, int rank, MPI_Comm comm) { int NLABELS=0; char VALUE=0; double AFFINITY=0.f; vector Label; vector Affinity; // Read the labels if (rank==0){ printf("Component labels:\n"); ifstream iFILE("ComponentLabels.csv");\ if (iFILE.good()){ while (!iFILE.eof()){ iFILE>>VALUE; iFILE>>AFFINITY; Label.push_back(VALUE); Affinity.push_back(AFFINITY); NLABELS++; printf("%i %f\n",VALUE,AFFINITY); } } else{ printf("Using default labels: Solid (0 --> -1.0), NWP (1 --> 1.0), WP (2 --> -1.0)\n"); // Set default values VALUE=0; AFFINITY=-1.0; Label.push_back(VALUE); Affinity.push_back(AFFINITY); NLABELS++; printf("%i %f\n",VALUE,AFFINITY); VALUE=1; AFFINITY=1.0; Label.push_back(VALUE); Affinity.push_back(AFFINITY); NLABELS++; printf("%i %f\n",VALUE,AFFINITY); VALUE=2; AFFINITY=-1.0; Label.push_back(VALUE); Affinity.push_back(AFFINITY); NLABELS++; printf("%i %f\n",VALUE,AFFINITY); } } // Broadcast the list MPI_Bcast(&NLABELS,1,MPI_INT,0,comm); // Copy into contiguous buffers char *LabelList; double * AffinityList; LabelList=new char[NLABELS]; AffinityList=new double[NLABELS]; MPI_Bcast(&LabelList,NLABELS,MPI_CHAR,0,comm); MPI_Bcast(&AffinityList,NLABELS,MPI_DOUBLE,0,comm); // Assign the labels for (int k=0;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("********************************************************\n"); } MPI_Barrier(comm); 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); for (k=0;k 0){ sum_local+=1.0; Np++; } } } } MPI_Allreduce(&sum_local,&sum,1,MPI_DOUBLE,MPI_SUM,comm); porosity = sum*iVol_global; if (rank==0) printf("Media porosity = %f \n",porosity); if (rank==0) printf ("Create ScaLBL_Communicator \n"); MPI_Barrier(comm); // Create a communicator for the device (will use optimized layout) ScaLBL_Communicator ScaLBL_Comm(Dm); //Create a second communicator based on the regular data layout ScaLBL_Communicator ScaLBL_Comm_Regular(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); MPI_Barrier(comm); //......................device distributions................................. int dist_mem_size = Np*sizeof(double); if (rank==0) printf ("Allocating distributions \n"); int *NeighborList; int *dvcMap; // double *f_even,*f_odd; double *fq, *Aq, *Bq; double *Den, *Phi; double *ColorGrad; double *Vel; double *Pressure; //........................................................................... ScaLBL_AllocateDeviceMemory((void **) &NeighborList, neighborSize); ScaLBL_AllocateDeviceMemory((void **) &dvcMap, sizeof(int)*Np); ScaLBL_AllocateDeviceMemory((void **) &fq, 19*dist_mem_size); ScaLBL_AllocateDeviceMemory((void **) &Aq, 7*dist_mem_size); ScaLBL_AllocateDeviceMemory((void **) &Bq, 7*dist_mem_size); ScaLBL_AllocateDeviceMemory((void **) &Den, 2*dist_mem_size); ScaLBL_AllocateDeviceMemory((void **) &Phi, sizeof(double)*Nx*Ny*Nz); ScaLBL_AllocateDeviceMemory((void **) &Pressure, sizeof(double)*Np); ScaLBL_AllocateDeviceMemory((void **) &Vel, 3*sizeof(double)*Np); ScaLBL_AllocateDeviceMemory((void **) &ColorGrad, 3*sizeof(double)*Np); //........................................................................... // Update GPU data structures if (rank==0) printf ("Setting up device map and neighbor list \n"); int *TmpMap; TmpMap=new int[Np]; 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); /* double *PHASE; PHASE= new double [Nx*Ny*Nz]; SIZE=Nx*Ny*Nz*sizeof(double); ScaLBL_CopyToHost(&PHASE[0],&Phi[0],SIZE); FILE *OUTFILE; sprintf(LocalRankFilename,"Phase.%05i.raw",rank); OUTFILE = fopen(LocalRankFilename,"wb"); fwrite(PHASE,8,N,OUTFILE); fclose(OUTFILE); */ } // **************************************************** MPI_Barrier(comm); MPI_Finalize(); // **************************************************** return check; }