#include #include #include #include #include #include #include //#include "common/pmmc.h" #include "common/Domain.h" #include "common/SpherePack.h" #include "common/MPI.h" #include "common/Communication.h" /* * Pre-Processor to generate signed distance function from sphere packing * to use as an input domain for lattice Boltzmann simulator * James E. McClure 2014 */ using namespace std; inline void PackID(const 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 maxdist) maxdist=SignDist(i,j,k); } } } } // total Global is the number of nodes in the pore-space totalGlobal = Dm.Comm.sumReduce( count ); maxdistGlobal = Dm.Comm.sumReduce( maxdist ); double volume=double(nprocx*nprocy*nprocz)*double(nx-2)*double(ny-2)*double(nz-2); double porosity=totalGlobal/volume; if (rank==0) printf("Media Porosity: %f \n",porosity); if (rank==0) printf("Maximum pore size: %f \n",maxdistGlobal); // Generate the NWP configuration //if (rank==0) printf("Initializing morphological distribution with critical radius %f \n", Rcrit); if (rank==0) printf("Performing morphological opening with target saturation %f \n", SW); // GenerateResidual(id,nx,ny,nz,Saturation); // Communication buffers char *sendID_x, *sendID_y, *sendID_z, *sendID_X, *sendID_Y, *sendID_Z; char *sendID_xy, *sendID_yz, *sendID_xz, *sendID_Xy, *sendID_Yz, *sendID_xZ; char *sendID_xY, *sendID_yZ, *sendID_Xz, *sendID_XY, *sendID_YZ, *sendID_XZ; char *recvID_x, *recvID_y, *recvID_z, *recvID_X, *recvID_Y, *recvID_Z; char *recvID_xy, *recvID_yz, *recvID_xz, *recvID_Xy, *recvID_Yz, *recvID_xZ; char *recvID_xY, *recvID_yZ, *recvID_Xz, *recvID_XY, *recvID_YZ, *recvID_XZ; // send buffers sendID_x = new char [Dm.sendCount("x")]; sendID_y = new char [Dm.sendCount("y")]; sendID_z = new char [Dm.sendCount("z")]; sendID_X = new char [Dm.sendCount("X")]; sendID_Y = new char [Dm.sendCount("Y")]; sendID_Z = new char [Dm.sendCount("Z")]; sendID_xy = new char [Dm.sendCount("xy")]; sendID_yz = new char [Dm.sendCount("yz")]; sendID_xz = new char [Dm.sendCount("xz")]; sendID_Xy = new char [Dm.sendCount("Xy")]; sendID_Yz = new char [Dm.sendCount("Yz")]; sendID_xZ = new char [Dm.sendCount("xZ")]; sendID_xY = new char [Dm.sendCount("xY")]; sendID_yZ = new char [Dm.sendCount("yZ")]; sendID_Xz = new char [Dm.sendCount("Xz")]; sendID_XY = new char [Dm.sendCount("XY")]; sendID_YZ = new char [Dm.sendCount("YZ")]; sendID_XZ = new char [Dm.sendCount("XZ")]; //...................................................................................... // recv buffers recvID_x = new char [Dm.recvCount("x")]; recvID_y = new char [Dm.recvCount("y")]; recvID_z = new char [Dm.recvCount("z")]; recvID_X = new char [Dm.recvCount("X")]; recvID_Y = new char [Dm.recvCount("Y")]; recvID_Z = new char [Dm.recvCount("Z")]; recvID_xy = new char [Dm.recvCount("xy")]; recvID_yz = new char [Dm.recvCount("yz")]; recvID_xz = new char [Dm.recvCount("xz")]; recvID_Xy = new char [Dm.recvCount("Xy")]; recvID_xZ = new char [Dm.recvCount("xZ")]; recvID_xY = new char [Dm.recvCount("xY")]; recvID_yZ = new char [Dm.recvCount("yZ")]; recvID_Yz = new char [Dm.recvCount("Yz")]; recvID_Xz = new char [Dm.recvCount("Xz")]; recvID_XY = new char [Dm.recvCount("XY")]; recvID_YZ = new char [Dm.recvCount("YZ")]; recvID_XZ = new char [Dm.recvCount("XZ")]; //...................................................................................... int sendtag,recvtag; sendtag = recvtag = 7; int ii,jj,kk; int Nx = nx; int Ny = ny; int Nz = nz; double sw_old=1.0; double sw_new=1.0; double sw_diff_old = 1.0; double sw_diff_new = 1.0; // Increase the critical radius until the target saturation is met double deltaR=0.05; // amount to change the radius in voxel units double Rcrit_old; double Rcrit_new; int imin,jmin,kmin,imax,jmax,kmax; Rcrit_new = maxdistGlobal; while (sw_new > SW) { sw_diff_old = sw_diff_new; sw_old = sw_new; Rcrit_old = Rcrit_new; Rcrit_new -= deltaR*Rcrit_old; int Window=round(Rcrit_new); if (Window == 0) Window = 1; // If Window = 0 at the begining, after the following process will have sw=1.0 // and sw Rcrit_new){ // loop over the window and update imin=max(1,i-Window); jmin=max(1,j-Window); kmin=max(1,k-Window); imax=min(Nx-1,i+Window); jmax=min(Ny-1,j+Window); kmax=min(Nz-1,k+Window); for (kk=kmin; kk( filename ); auto domain_db = db->getDatabase( "Domain" ); auto Dm = std::shared_ptr(new Domain(domain_db,comm)); // full domain for analysis Nx = Dm->Nx; Ny = Dm->Ny; Nz = Dm->Nz; Lx = Dm->Lx; Ly = Dm->Ly; Lz = Dm->Lz; iproc = Dm->iproc(); jproc = Dm->jproc(); kproc = Dm->kproc(); nprocx = Dm->nprocx(); nprocy = Dm->nprocy(); nprocz = Dm->nprocz(); nspheres = domain_db->getScalar( "nspheres"); //printf("Set domain \n"); //int BoundaryCondition=1; //Nz += 2; //Nx = Ny = Nz; // Cubic domain int N = Nx*Ny*Nz; // Define Dm.Communication sub-domain -- everywhere for (int k=0; kid[n] = 1; } } } Dm->CommInit(); //....................................................................... // Filenames used char LocalRankString[8]; char LocalRankFilename[40]; char LocalRestartFile[40]; sprintf(LocalRankString,"%05d",rank); sprintf(LocalRankFilename,"%s%s","ID.",LocalRankString); sprintf(LocalRestartFile,"%s%s","Restart.",LocalRankString); // printf("Local File Name = %s \n",LocalRankFilename); // .......... READ THE INPUT FILE ....................................... // char value; char *id; id = new char[N]; int sum = 0; double sum_local; double iVol_global = 1.0/(1.0*(Nx-2)*(Ny-2)*(Nz-2)*nprocs); double porosity; //........................................................................... DoubleArray SignDist(Nx,Ny,Nz); //....................................................................... // Read in sphere pack //if (rank==1) printf("nspheres =%i \n",nspheres); //....................................................................... double *cx,*cy,*cz,*rad; cx = new double[nspheres]; cy = new double[nspheres]; cz = new double[nspheres]; rad = new double[nspheres]; //....................................................................... if (rank == 0) printf("Reading the sphere packing \n"); if (rank == 0) ReadSpherePacking(nspheres,cx,cy,cz,rad); comm.barrier(); // Broadcast the sphere packing to all processes comm.bcast(cx,nspheres,0); comm.bcast(cy,nspheres,0); comm.bcast(cz,nspheres,0); comm.bcast(rad,nspheres,0); //........................................................................... comm.barrier(); if (rank == 0) cout << "Domain set." << endl; if (rank == 0){ // Compute the Sauter mean diameter double totVol = 0.0; double totArea = 0.0; // Compute the total volume and area of all spheres for (i=0; i 0.0){ id[n] = 2; } // compute the porosity (actual interface location used) if (SignDist(n) > 0.0){ sum++; } } } } sum_local = 1.0*sum; porosity = comm.sumReduce(sum_local); porosity = porosity*iVol_global; if (rank==0) printf("Media porosity = %f \n",porosity); // Run Morphological opening to initialize 50% saturation double SW=0.50; if (rank==0) printf("MorphOpen: Initializing with saturation %f \n",SW); MorphOpen(SignDist, id, *Dm, Nx, Ny, Nz, rank, SW); //......................................................... // don't perform computations at the eight corners id[0] = id[Nx-1] = id[(Ny-1)*Nx] = id[(Ny-1)*Nx + Nx-1] = 0; id[(Nz-1)*Nx*Ny] = id[(Nz-1)*Nx*Ny+Nx-1] = id[(Nz-1)*Nx*Ny+(Ny-1)*Nx] = id[(Nz-1)*Nx*Ny+(Ny-1)*Nx + Nx-1] = 0; //......................................................... /* //....................................................................... sprintf(LocalRankString,"%05d",rank); sprintf(LocalRankFilename,"%s%s","SignDist.",LocalRankString); FILE *DIST = fopen(LocalRankFilename,"wb"); if (DIST==NULL) ERROR("Error opening file: ID.xxxxx"); fwrite(SignDist.data(),1,N*sizeof(double),DIST); fclose(DIST); //...................................................................... */ //....................................................................... sprintf(LocalRankFilename,"%s%s","ID.",LocalRankString); FILE *IDFILE = fopen(LocalRankFilename,"wb"); if (IDFILE==NULL) ERROR("Error opening file: ID.xxxxx"); fwrite(id,1,N,IDFILE); fclose(IDFILE); //...................................................................... } Utilities::shutdown(); }