// Sequential blob analysis // Reads parallel simulation data and performs connectivity analysis // and averaging on a blob-by-blob basis // James E. McClure 2014 #include #include #include "common/Communication.h" #include "analysis/analysis.h" #include "common/TwoPhase.h" //#include "Domain.h" using namespace std; void readRankData( int proc, int nx, int ny, int nz, DoubleArray& Phase, DoubleArray& SignDist ) { Phase.resize(nx,ny,nz); SignDist.resize(nx,ny,nz); char file1[40], file2[40]; sprintf(file1,"SignDist.%05d",proc); sprintf(file2,"Phase.%05d",proc); ReadBinaryFile(file1, Phase.get(), nx*ny*nz); ReadBinaryFile(file2, SignDist.get(), nx*ny*nz); } inline void WriteBlobs(TwoPhase Averages){ printf("Writing the blob list \n"); FILE *BLOBLOG; BLOBLOG=fopen("blobs.tcat","w"); fprintf(BLOBLOG,"%.5g %.5g %.5g\n",Averages.vol_w_global,Averages.paw_global,Averages.aws_global); for (int b=0; b<(int)Averages.ComponentAverages_NWP.size(1); b++){ if (Averages.ComponentAverages_NWP(0,b) > 0.0){ double Vn,pn,awn,ans,Jwn,Kwn,lwns,cwns; Vn = Averages.ComponentAverages_NWP(1,b); pn = Averages.ComponentAverages_NWP(2,b); awn = Averages.ComponentAverages_NWP(3,b); ans = Averages.ComponentAverages_NWP(4,b); Jwn = Averages.ComponentAverages_NWP(5,b); Kwn = Averages.ComponentAverages_NWP(6,b); lwns = Averages.ComponentAverages_NWP(7,b); cwns = Averages.ComponentAverages_NWP(8,b); fprintf(BLOBLOG,"%.5g ", Vn); //Vn fprintf(BLOBLOG,"%.5g ", pn); //pn fprintf(BLOBLOG,"%.5g ", awn); //awn fprintf(BLOBLOG,"%.5g ", ans); //ans fprintf(BLOBLOG,"%.5g ", Jwn); //Jwn fprintf(BLOBLOG,"%.5g ", Kwn); //Kwn fprintf(BLOBLOG,"%.5g ", lwns); //lwns fprintf(BLOBLOG,"%.5g\n",cwns); //cwns } } fclose(BLOBLOG); } inline void WriteBlobStates(TwoPhase TCAT, double D, double porosity){ int a; double iVol=1.0/TCAT.Dm.Volume; double PoreVolume; double nwp_volume,vol_n,pan,pn,pw,pawn,pwn,awn,ans,aws,Jwn,Kwn,lwns,cwns,clwns; double sw,awnD,awsD,ansD,lwnsDD,JwnD,pc; nwp_volume=vol_n=pan=awn=ans=Jwn=Kwn=lwns=clwns=pawn=0.0; sw = TCAT.sat_w; pw = TCAT.paw_global; aws = TCAT.aws; // Compute the averages over the entire non-wetting phase printf("Writing blobstates.tcat for %i components \n",TCAT.NumberComponents_NWP); FILE *BLOBSTATES; BLOBSTATES = fopen("./blobstates.tcat","w"); if (BLOBSTATES==NULL) ERROR("Cannot open blobstates.tcat for writing"); for (a=0; a0; a--){ // Subtract the features one-by-one vol_n -= TCAT.ComponentAverages_NWP(0,a); pan -= TCAT.ComponentAverages_NWP(2,a)*TCAT.ComponentAverages_NWP(0,a); awn -= TCAT.ComponentAverages_NWP(3,a); ans -= TCAT.ComponentAverages_NWP(4,a); Jwn -= TCAT.ComponentAverages_NWP(5,a)*TCAT.ComponentAverages_NWP(3,a); Kwn -= TCAT.ComponentAverages_NWP(6,a)*TCAT.ComponentAverages_NWP(3,a); lwns -= TCAT.ComponentAverages_NWP(7,a); clwns -= TCAT.ComponentAverages_NWP(8,a)*TCAT.ComponentAverages_NWP(7,a); nwp_volume -= TCAT.ComponentAverages_NWP(1,a); pawn -= TCAT.ComponentAverages_NWP(2,a)*TCAT.ComponentAverages_NWP(3,a); // Update wetting phase averages aws += TCAT.ComponentAverages_NWP(4,a); if (vol_n > 64){ // Only consider systems with "large enough" blobs -- 4^3 if (fabs(1.0 - nwp_volume/PoreVolume - sw) > 0.005 || a == 1){ sw = 1.0 - nwp_volume/PoreVolume; JwnD = Jwn*D/awn; //trJwnD = -trJwn*D/trawn; cwns = clwns / lwns; pwn = (pawn/awn-pw)*D/0.058; pn = pan/vol_n; awnD = awn*D*iVol; awsD = aws*D*iVol; ansD = ans*D*iVol; lwnsDD = lwns*D*D*iVol; pc = (pn-pw)*D/0.058; // hard-coded surface tension due to being lazy fprintf(BLOBSTATES,"%.5g %.5g %.5g ",sw,pn,pw); fprintf(BLOBSTATES,"%.5g %.5g %.5g %.5g ",awnD,awsD,ansD,lwnsDD); fprintf(BLOBSTATES,"%.5g %.5g %.5g %.5g %i\n",pc,pwn,JwnD,cwns,a); } } } fclose(BLOBSTATES); } int main(int argc, char **argv) { // Initialize MPI int rank, nprocs; MPI_Init(&argc,&argv); MPI_Comm_rank(MPI_COMM_WORLD,&rank); MPI_Comm_size(MPI_COMM_WORLD,&nprocs); { // Limit scope so variables that contain communicators will free before MPI_Finialize if ( rank==0 ) { printf("-----------------------------------------------------------\n"); printf("Labeling Blobs from Two-Phase Lattice Boltzmann Simulation \n"); printf("-----------------------------------------------------------\n"); } //....................................................................... // Reading the domain information file //....................................................................... int nprocx, nprocy, nprocz, nx, ny, nz, nspheres; double Lx, Ly, Lz; int Nx,Ny,Nz; int i,j,k,n; if (rank==0){ ifstream domain("Domain.in"); domain >> nprocx; domain >> nprocy; domain >> nprocz; domain >> nx; domain >> ny; domain >> nz; domain >> nspheres; domain >> Lx; domain >> Ly; domain >> Lz; } MPI_Barrier(MPI_COMM_WORLD); // Computational domain MPI_Bcast(&nx,1,MPI_INT,0,MPI_COMM_WORLD); MPI_Bcast(&ny,1,MPI_INT,0,MPI_COMM_WORLD); MPI_Bcast(&nz,1,MPI_INT,0,MPI_COMM_WORLD); MPI_Bcast(&nprocx,1,MPI_INT,0,MPI_COMM_WORLD); MPI_Bcast(&nprocy,1,MPI_INT,0,MPI_COMM_WORLD); MPI_Bcast(&nprocz,1,MPI_INT,0,MPI_COMM_WORLD); MPI_Bcast(&nspheres,1,MPI_INT,0,MPI_COMM_WORLD); MPI_Bcast(&Lx,1,MPI_DOUBLE,0,MPI_COMM_WORLD); MPI_Bcast(&Ly,1,MPI_DOUBLE,0,MPI_COMM_WORLD); MPI_Bcast(&Lz,1,MPI_DOUBLE,0,MPI_COMM_WORLD); //................................................. MPI_Barrier(MPI_COMM_WORLD); // 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"); // Set up the domain int BC=0; // Get the rank info Domain Dm(nx,ny,nz,rank,nprocx,nprocy,nprocz,Lx,Ly,Lz,BC); // const RankInfoStruct rank_info(rank,nprocx,nprocy,nprocz); TwoPhase Averages(Dm); int N = (nx+2)*(ny+2)*(nz+2); Nx = nx+2; Ny = ny+2; Nz = nz+2; if (rank == 0) cout << "Domain set." << endl; //....................................................................... for ( k=1;k 0.0){ Dm.id[n] = 2; } else{ Dm.id[n] = 1; } Averages.SDn(i,j,k) = -Averages.Phase(i,j,k); Averages.Phase(i,j,k) = Averages.SDn(i,j,k); Averages.Phase_tplus(i,j,k) = Averages.SDn(i,j,k); Averages.Phase_tminus(i,j,k) = Averages.SDn(i,j,k); Averages.DelPhi(i,j,k) = 0.0; Averages.Press(i,j,k) = 0.0; Averages.Vel_x(i,j,k) = 0.0; Averages.Vel_y(i,j,k) = 0.0; Averages.Vel_z(i,j,k) = 0.0; } } } double vF,vS; vF = vS = 0.0; double beta = 0.95; if (rank==0) printf("initializing the system \n"); Averages.UpdateSolid(); Averages.Initialize(); Averages.UpdateMeshValues(); Dm.CommunicateMeshHalo(Averages.Phase); Dm.CommunicateMeshHalo(Averages.SDn); // if (rank==0) printf("computing blobs \n"); // int nblobs_global = ComputeGlobalBlobIDs(Dm.Nx-2,Dm.Ny-2,Dm.Nz-2,Dm.rank_info, // Averages.Phase,Averages.SDs,vF,vS,Averages.BlobLabel); // if (Dm.rank==0) printf("Number of blobs is %i \n",nblobs_global); if (rank==0) printf("computing local averages \n"); Averages.AssignComponentLabels(); Averages.ComponentAverages(); Averages.PrintComponents(int(5)); if (rank==0) printf("reducing averages \n"); // Averages.Reduce(); } // Limit scope so variables that contain communicators will free before MPI_Finialize MPI_Barrier(MPI_COMM_WORLD); MPI_Finalize(); return 0; }