/* * Pre-processor to generate signed distance function from segmented data * segmented data should be stored in a raw binary file as 1-byte integer (type char) * will output distance functions for phases */ #include #include #include #include #include #include #include "common/Array.h" #include "common/Domain.h" //************************************************************************* // Morpohological drainage pre-processor // Generate states on primary drainage using morphological approach // Signed distance function is used to determine fluid configuration //************************************************************************* 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 1){ filename=argv[1]; } else ERROR("No input database provided\n"); // read the input database auto db = std::make_shared( filename ); auto domain_db = db->getDatabase( "Domain" ); // Read domain parameters auto L = domain_db->getVector( "L" ); auto size = domain_db->getVector( "n" ); auto nproc = domain_db->getVector( "nproc" ); double SW = domain_db->getScalar( "Sw" ); nx = size[0]; ny = size[1]; nz = size[2]; nprocx = nproc[0]; nprocy = nproc[1]; nprocz = nproc[2]; // Check that the number of processors >= the number of ranks if ( rank==0 ) { printf("Number of MPI ranks required: %i \n", nprocx*nprocy*nprocz); printf("Number of MPI ranks used: %i \n", nprocs); printf("Full domain size: %i x %i x %i \n",nx*nprocx,ny*nprocy,nz*nprocz); } if ( nprocs < nprocx*nprocy*nprocz ){ ERROR("Insufficient number of processors"); } char LocalRankFilename[40]; int BoundaryCondition=1; Domain Dm(nx,ny,nz,rank,nprocx,nprocy,nprocz,Lx,Ly,Lz,BoundaryCondition); nx+=2; ny+=2; nz+=2; int N = nx*ny*nz; char *id; id = new char[N]; // Define communication sub-domain -- everywhere for (int k=0; k fillData(Dm.Comm, Dm.rank_info,{xdim,ydim,zdim},{1,1,1},0,1); DoubleArray SignDist(nx,ny,nz); // Read the signed distance from file sprintf(LocalRankFilename,"SignDist.%05i",rank); FILE *DIST = fopen(LocalRankFilename,"rb"); size_t ReadSignDist; ReadSignDist=fread(SignDist.data(),8,N,DIST); if (ReadSignDist != size_t(N)) printf("lbpm_morphdrain_pp: Error reading signed distance function (rank=%i)\n",rank); fclose(DIST); fillData.fill(SignDist); sprintf(LocalRankFilename,"ID.%05i",rank); size_t readID; FILE *IDFILE = fopen(LocalRankFilename,"rb"); if (IDFILE==NULL) ERROR("Error opening file: ID.xxxxx"); readID=fread(id,1,N,IDFILE); if (readID != size_t(N)) printf("lbpm_segmented_pp: Error reading ID (rank=%i) \n",rank); fclose(IDFILE); Dm.CommInit(); int iproc = Dm.iproc(); int jproc = Dm.jproc(); int kproc = Dm.kproc(); // Generate the NWP configuration //if (rank==0) printf("Performing morphological drainage with critical radius %f \n", Rcrit); if (rank==0) printf("Performing morphological drainage 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 x,y,z; int ii,jj,kk; int Nx = nx; int Ny = ny; int Nz = nz; double GlobalNumber = 1.f; double count,countGlobal,totalGlobal; count = 0.f; for (int k=0; k radius){ radius=SignDist(i,j,0); } } } } MPI_Allreduce(&radius,&Rcrit_new,1,MPI_DOUBLE,MPI_MAX,comm); if (rank==0) printf("Starting morhpological drainage with critical radius = %f \n",Rcrit_new); int imin,jmin,kmin,imax,jmax,kmax; // Decrease the critical radius until the target saturation is met double deltaR=0.01; // amount to change the radius in voxel units double Rcrit_old; double sw_old=1.0; // initial WP saturation set to one double sw_new=1.0; // initial WP saturation set to one double sw_diff_old = 1.0; double sw_diff_new = 1.0; while (sw_new > SW && Rcrit_new > 0.99){ Rcrit_old = Rcrit_new; Rcrit_new -= deltaR*Rcrit_old;// decrease critical radius sw_old = sw_new; sw_diff_old = sw_diff_new; int Window=round(Rcrit_new); GlobalNumber = 1.0; while (GlobalNumber > 0){ //if (rank==0) printf("GlobalNumber=%f \n",GlobalNumber); double LocalNumber=GlobalNumber=0.f; for(k=0; k 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