/* color lattice boltzmann model */ #include "models/ColorModel.h" #include "analysis/distance.h" #include "analysis/morphology.h" ScaLBL_ColorModel::ScaLBL_ColorModel(int RANK, int NP, MPI_Comm COMM): rank(RANK), nprocs(NP), Restart(0),timestep(0),timestepMax(0),tauA(0),tauB(0),rhoA(0),rhoB(0),alpha(0),beta(0), Fx(0),Fy(0),Fz(0),flux(0),din(0),dout(0),inletA(0),inletB(0),outletA(0),outletB(0), Nx(0),Ny(0),Nz(0),N(0),Np(0),nprocx(0),nprocy(0),nprocz(0),BoundaryCondition(0),Lx(0),Ly(0),Lz(0),comm(COMM) { } ScaLBL_ColorModel::~ScaLBL_ColorModel(){ } /*void ScaLBL_ColorModel::WriteCheckpoint(const char *FILENAME, const double *cPhi, const double *cfq, int Np) { int q,n; double value; ofstream File(FILENAME,ios::binary); for (n=0; n( filename ); domain_db = db->getDatabase( "Domain" ); color_db = db->getDatabase( "Color" ); analysis_db = db->getDatabase( "Analysis" ); // Color Model parameters timestepMax = color_db->getScalar( "timestepMax" ); tauA = color_db->getScalar( "tauA" ); tauB = color_db->getScalar( "tauB" ); rhoA = color_db->getScalar( "rhoA" ); rhoB = color_db->getScalar( "rhoB" ); Fx = color_db->getVector( "F" )[0]; Fy = color_db->getVector( "F" )[1]; Fz = color_db->getVector( "F" )[2]; alpha = color_db->getScalar( "alpha" ); beta = color_db->getScalar( "beta" ); Restart = color_db->getScalar( "Restart" ); din = color_db->getScalar( "din" ); dout = color_db->getScalar( "dout" ); flux = color_db->getScalar( "flux" ); inletA=1.f; inletB=0.f; outletA=0.f; outletB=1.f; // Read domain parameters auto L = domain_db->getVector( "L" ); auto size = domain_db->getVector( "n" ); auto nproc = domain_db->getVector( "nproc" ); BoundaryCondition = domain_db->getScalar( "BC" ); Nx = size[0]; Ny = size[1]; Nz = size[2]; Lx = L[0]; Ly = L[1]; Lz = L[2]; nprocx = nproc[0]; nprocy = nproc[1]; nprocz = nproc[2]; if (BoundaryCondition==4) flux *= rhoA; // mass flux must adjust for density (see formulation for details) } void ScaLBL_ColorModel::SetDomain(){ Dm = std::shared_ptr(new Domain(domain_db,comm)); // full domain for analysis Mask = std::shared_ptr(new Domain(domain_db,comm)); // mask domain removes immobile phases Nx+=2; Ny+=2; Nz += 2; N = Nx*Ny*Nz; id = new char [N]; for (int i=0; iid[i] = 1; // initialize this way Averages = std::shared_ptr ( new TwoPhase(Dm) ); // TwoPhase analysis object MPI_Barrier(comm); Dm->CommInit(); MPI_Barrier(comm); rank = Dm->rank(); } void ScaLBL_ColorModel::ReadInput(){ size_t readID; Mask->ReadIDs(); for (int i=0; iid[i]; // save what was read sprintf(LocalRankString,"%05d",rank); sprintf(LocalRankFilename,"%s%s","ID.",LocalRankString); sprintf(LocalRestartFile,"%s%s","Restart.",LocalRankString); // Generate the signed distance map // Initialize the domain and communication Array id_solid(Nx,Ny,Nz); int count = 0; // Solve for the position of the solid phase for (int k=0;kid[n] > 0) id_solid(i,j,k) = 1; else id_solid(i,j,k) = 0; } } } // Initialize the signed distance function for (int k=0;kSDs(i,j,k) = 2.0*double(id_solid(i,j,k))-1.0; } } } // MeanFilter(Averages->SDs); if (rank==0) printf("Initialized solid phase -- Converting to Signed Distance function \n"); CalcDist(Averages->SDs,id_solid,*Mask); if (rank == 0) cout << "Domain set." << endl; Averages->SetParams(rhoA,rhoB,tauA,tauB,Fx,Fy,Fz,alpha); } void ScaLBL_ColorModel::AssignComponentLabels(double *phase) { size_t NLABELS=0; char VALUE=0; double AFFINITY=0.f; auto LabelList = color_db->getVector( "ComponentLabels" ); auto AffinityList = color_db->getVector( "ComponentAffinity" ); NLABELS=LabelList.size(); if (NLABELS != AffinityList.size()){ ERROR("Error: ComponentLabels and ComponentAffinity must be the same length! \n"); } double label_count[NLABELS]; double label_count_global[NLABELS]; // Assign the labels for (int idx=0; idxid[n] = 0; // set mask to zero since this is an immobile component } } // fluid labels are reserved if (VALUE == 1) AFFINITY=1.0; else if (VALUE == 2) AFFINITY=-1.0; phase[n] = AFFINITY; } } } // Set Dm to match Mask for (int i=0; iid[i] = Mask->id[i]; for (int idx=0; idxComm, label_count[idx]); if (rank==0){ printf("Components labels: %lu \n",NLABELS); for (unsigned int idx=0; idxid[i] = Mask->id[i]; Mask->CommInit(); Np=Mask->PoreCount(); //........................................................................... if (rank==0) printf ("Create ScaLBL_Communicator \n"); // Create a communicator for the device (will use optimized layout) // ScaLBL_Communicator ScaLBL_Comm(Mask); // original ScaLBL_Comm = std::shared_ptr(new ScaLBL_Communicator(Mask)); ScaLBL_Comm_Regular = std::shared_ptr(new ScaLBL_Communicator(Mask)); int Npad=(Np/16 + 2)*16; if (rank==0) printf ("Set up memory efficient layout, %i | %i | %i \n", Np, Npad, N); Map.resize(Nx,Ny,Nz); Map.fill(-2); auto neighborList= new int[18*Npad]; Np = ScaLBL_Comm->MemoryOptimizedLayoutAA(Map,neighborList,Mask->id,Np); MPI_Barrier(comm); //........................................................................... // MAIN VARIABLES ALLOCATED HERE //........................................................................... // LBM variables if (rank==0) printf ("Allocating distributions \n"); //......................device distributions................................. dist_mem_size = Np*sizeof(double); neighborSize=18*(Np*sizeof(int)); //........................................................................... 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 **) &Velocity, 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"); fflush(stdout); int *TmpMap; TmpMap=new int[Np]; for (int k=1; kLastExterior(); idx++){ int n = TmpMap[idx]; if (n > Nx*Ny*Nz){ printf("Bad value! idx=%i \n"); TmpMap[idx] = Nx*Ny*Nz-1; } } for (int idx=ScaLBL_Comm->FirstInterior(); idxLastInterior(); idx++){ int n = TmpMap[idx]; if (n > Nx*Ny*Nz){ printf("Bad value! idx=%i \n"); TmpMap[idx] = Nx*Ny*Nz-1; } } ScaLBL_CopyToDevice(dvcMap, TmpMap, sizeof(int)*Np); ScaLBL_DeviceBarrier(); delete [] TmpMap; // copy the neighbor list ScaLBL_CopyToDevice(NeighborList, neighborList, neighborSize); // initialize phi based on PhaseLabel (include solid component labels) double *PhaseLabel; PhaseLabel = new double[N]; AssignComponentLabels(PhaseLabel); ScaLBL_CopyToDevice(Phi, PhaseLabel, N*sizeof(double)); } /******************************************************** * AssignComponentLabels * ********************************************************/ void ScaLBL_ColorModel::Initialize(){ if (rank==0) printf ("Initializing distributions \n"); ScaLBL_D3Q19_Init(fq, Np); /* * This function initializes model */ if (Restart == true){ if (rank==0){ printf("Reading restart file! \n"); ifstream restart("Restart.txt"); if (restart.is_open()){ restart >> timestep; printf("Restarting from timestep =%i \n",timestep); } else{ printf("WARNING:No Restart.txt file, setting timestep=0 \n"); timestep=0; } } MPI_Bcast(×tep,1,MPI_INT,0,comm); // Read in the restart file to CPU buffers int *TmpMap; TmpMap = new int[Np]; double *cPhi, *cDist, *cDen; cPhi = new double[N]; cDen = new double[2*Np]; cDist = new double[19*Np]; ScaLBL_CopyToHost(TmpMap, dvcMap, Np*sizeof(int)); ScaLBL_CopyToHost(cPhi, Phi, N*sizeof(double)); ifstream File(LocalRestartFile,ios::binary); int idx; double value,va,vb; for (int n=0; nLastExterior(); n++){ va = cDen[n]; vb = cDen[Np + n]; value = (va-vb)/(va+vb); idx = TmpMap[n]; if (!(idx < 0) && idxFirstInterior(); nLastInterior(); n++){ va = cDen[n]; vb = cDen[Np + n]; value = (va-vb)/(va+vb); idx = TmpMap[n]; if (!(idx < 0) && idxLastExterior(), Np); ScaLBL_PhaseField_Init(dvcMap, Phi, Den, Aq, Bq, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np); if (BoundaryCondition >0 ){ if (Dm->kproc()==0){ ScaLBL_SetSlice_z(Phi,1.0,Nx,Ny,Nz,0); ScaLBL_SetSlice_z(Phi,1.0,Nx,Ny,Nz,1); ScaLBL_SetSlice_z(Phi,1.0,Nx,Ny,Nz,2); } if (Dm->kproc() == nprocz-1){ ScaLBL_SetSlice_z(Phi,-1.0,Nx,Ny,Nz,Nz-1); ScaLBL_SetSlice_z(Phi,-1.0,Nx,Ny,Nz,Nz-2); ScaLBL_SetSlice_z(Phi,-1.0,Nx,Ny,Nz,Nz-3); } } } void ScaLBL_ColorModel::Run(){ int nprocs=nprocx*nprocy*nprocz; const RankInfoStruct rank_info(rank,nprocx,nprocy,nprocz); bool SET_CAPILLARY_NUMBER = false; bool MORPH_ADAPT = false; bool USE_MORPH = false; int analysis_interval = 1000; // number of timesteps in between in situ analysis int MAX_MORPH_TIMESTEPS = 50000; // maximum number of LBM timesteps to spend in morphological adaptation routine int RAMP_TIMESTEPS = 50000; // number of timesteps to run initially (to get a reasonable velocity field before other pieces kick in) int morph_interval = 1000000; int CURRENT_MORPH_TIMESTEPS=0; // counter for number of timesteps spent in morphological adaptation routine (reset each time) double morph_delta = 0.0; int morph_timesteps = 0; double capillary_number = 0.0; double tolerance = 0.01; double Ca_previous = 0.f; double delta_volume = 0.0; double delta_volume_target = 0.0; double RESIDUAL_ENDPOINT_THRESHOLD = 0.04; if (color_db->keyExists( "residual_endpoint_threshold" )){ RESIDUAL_ENDPOINT_THRESHOLD = color_db->getScalar( "residual_endpoint_threshold" ); } if (color_db->keyExists( "capillary_number" )){ capillary_number = color_db->getScalar( "capillary_number" ); SET_CAPILLARY_NUMBER=true; } if (BoundaryCondition != 0 && SET_CAPILLARY_NUMBER==true){ if (rank == 0) printf("WARINING: capillary number target only supported for BC = 0 \n"); SET_CAPILLARY_NUMBER=false; } if (analysis_db->keyExists( "morph_delta" )){ morph_delta = analysis_db->getScalar( "morph_delta" ); } if (analysis_db->keyExists( "morph_interval" )){ morph_interval = analysis_db->getScalar( "morph_interval" ); USE_MORPH = true; } if (analysis_db->keyExists( "tolerance" )){ tolerance = analysis_db->getScalar( "tolerance" ); } if (analysis_db->keyExists( "analysis_interval" )){ analysis_interval = analysis_db->getScalar( "analysis_interval" ); } if (rank==0){ printf("********************************************************\n"); printf("No. of timesteps: %i \n", timestepMax); fflush(stdout); } //.......create and start timer............ double starttime,stoptime,cputime; ScaLBL_DeviceBarrier(); MPI_Barrier(comm); starttime = MPI_Wtime(); //......................................... //************ MAIN ITERATION LOOP ***************************************/ PROFILE_START("Loop"); //std::shared_ptr analysis_db; bool Regular = false; runAnalysis analysis( analysis_db, rank_info, ScaLBL_Comm, Dm, Np, Regular, beta, Map ); //analysis.createThreads( analysis_method, 4 ); while (timestep < timestepMax ) { //if ( rank==0 ) { printf("Running timestep %i (%i MB)\n",timestep+1,(int)(Utilities::getMemoryUsage()/1048576)); } PROFILE_START("Update"); // *************ODD TIMESTEP************* timestep++; // Compute the Phase indicator field // Read for Aq, Bq happens in this routine (requires communication) ScaLBL_Comm->BiSendD3Q7AA(Aq,Bq); //READ FROM NORMAL ScaLBL_D3Q7_AAodd_PhaseField(NeighborList, dvcMap, Aq, Bq, Den, Phi, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np); ScaLBL_Comm->BiRecvD3Q7AA(Aq,Bq); //WRITE INTO OPPOSITE ScaLBL_DeviceBarrier(); ScaLBL_D3Q7_AAodd_PhaseField(NeighborList, dvcMap, Aq, Bq, Den, Phi, 0, ScaLBL_Comm->LastExterior(), Np); // Perform the collision operation ScaLBL_Comm->SendD3Q19AA(fq); //READ FROM NORMAL if (BoundaryCondition > 0){ ScaLBL_Comm->Color_BC_z(dvcMap, Phi, Den, inletA, inletB); ScaLBL_Comm->Color_BC_Z(dvcMap, Phi, Den, outletA, outletB); } // Halo exchange for phase field ScaLBL_Comm_Regular->SendHalo(Phi); ScaLBL_D3Q19_AAodd_Color(NeighborList, dvcMap, fq, Aq, Bq, Den, Phi, Velocity, rhoA, rhoB, tauA, tauB, alpha, beta, Fx, Fy, Fz, Nx, Nx*Ny, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np); ScaLBL_Comm_Regular->RecvHalo(Phi); ScaLBL_Comm->RecvD3Q19AA(fq); //WRITE INTO OPPOSITE ScaLBL_DeviceBarrier(); // Set BCs if (BoundaryCondition == 3){ ScaLBL_Comm->D3Q19_Pressure_BC_z(NeighborList, fq, din, timestep); ScaLBL_Comm->D3Q19_Pressure_BC_Z(NeighborList, fq, dout, timestep); } if (BoundaryCondition == 4){ din = ScaLBL_Comm->D3Q19_Flux_BC_z(NeighborList, fq, flux, timestep); ScaLBL_Comm->D3Q19_Pressure_BC_Z(NeighborList, fq, dout, timestep); } ScaLBL_D3Q19_AAodd_Color(NeighborList, dvcMap, fq, Aq, Bq, Den, Phi, Velocity, rhoA, rhoB, tauA, tauB, alpha, beta, Fx, Fy, Fz, Nx, Nx*Ny, 0, ScaLBL_Comm->LastExterior(), Np); ScaLBL_DeviceBarrier(); MPI_Barrier(comm); // *************EVEN TIMESTEP************* timestep++; // Compute the Phase indicator field ScaLBL_Comm->BiSendD3Q7AA(Aq,Bq); //READ FROM NORMAL ScaLBL_D3Q7_AAeven_PhaseField(dvcMap, Aq, Bq, Den, Phi, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np); ScaLBL_Comm->BiRecvD3Q7AA(Aq,Bq); //WRITE INTO OPPOSITE ScaLBL_DeviceBarrier(); ScaLBL_D3Q7_AAeven_PhaseField(dvcMap, Aq, Bq, Den, Phi, 0, ScaLBL_Comm->LastExterior(), Np); // Perform the collision operation ScaLBL_Comm->SendD3Q19AA(fq); //READ FORM NORMAL // Halo exchange for phase field if (BoundaryCondition > 0){ ScaLBL_Comm->Color_BC_z(dvcMap, Phi, Den, inletA, inletB); ScaLBL_Comm->Color_BC_Z(dvcMap, Phi, Den, outletA, outletB); } ScaLBL_Comm_Regular->SendHalo(Phi); ScaLBL_D3Q19_AAeven_Color(dvcMap, fq, Aq, Bq, Den, Phi, Velocity, rhoA, rhoB, tauA, tauB, alpha, beta, Fx, Fy, Fz, Nx, Nx*Ny, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np); ScaLBL_Comm_Regular->RecvHalo(Phi); ScaLBL_Comm->RecvD3Q19AA(fq); //WRITE INTO OPPOSITE ScaLBL_DeviceBarrier(); // Set boundary conditions if (BoundaryCondition == 3){ ScaLBL_Comm->D3Q19_Pressure_BC_z(NeighborList, fq, din, timestep); ScaLBL_Comm->D3Q19_Pressure_BC_Z(NeighborList, fq, dout, timestep); } else if (BoundaryCondition == 4){ din = ScaLBL_Comm->D3Q19_Flux_BC_z(NeighborList, fq, flux, timestep); ScaLBL_Comm->D3Q19_Pressure_BC_Z(NeighborList, fq, dout, timestep); } ScaLBL_D3Q19_AAeven_Color(dvcMap, fq, Aq, Bq, Den, Phi, Velocity, rhoA, rhoB, tauA, tauB, alpha, beta, Fx, Fy, Fz, Nx, Nx*Ny, 0, ScaLBL_Comm->LastExterior(), Np); ScaLBL_DeviceBarrier(); MPI_Barrier(comm); //************************************************************************ MPI_Barrier(comm); PROFILE_STOP("Update"); // Run the analysis analysis.run( timestep, *Averages, Phi, Pressure, Velocity, fq, Den ); // allow initial ramp-up to get closer to steady state if (timestep > RAMP_TIMESTEPS && timestep%analysis_interval == 0 && USE_MORPH){ analysis.finish(); if ( morph_timesteps > morph_interval ){ double volB = Averages->Volume_w(); double volA = Averages->Volume_n(); double vA_x = Averages->van_global(0); double vA_y = Averages->van_global(1); double vA_z = Averages->van_global(2); double vB_x = Averages->vaw_global(0); double vB_y = Averages->vaw_global(1); double vB_z = Averages->vaw_global(2); double muA = rhoA*(tauA-0.5)/3.f; double muB = rhoB*(tauB-0.5)/3.f; double flow_rate_A = sqrt(vA_x*vA_x + vA_y*vA_y + vA_z*vA_z); double flow_rate_B = sqrt(vB_x*vB_x + vB_y*vB_y + vB_z*vB_z); double current_saturation = volB/(volA+volB); double Ca = fabs(volA*muA*flow_rate_A + volB*muB*flow_rate_B)/(5.796*alpha*double((Nx-2)*(Ny-2)*(Nz-2)*nprocs)); double force_magnitude = sqrt(Fx*Fx + Fy*Fy + Fz*Fz); if (fabs((Ca - Ca_previous)/Ca) < tolerance ){ MORPH_ADAPT = true; CURRENT_MORPH_TIMESTEPS=0; delta_volume_target = (volA + volB)*morph_delta; // set target volume chnage if (rank==0){ printf("** WRITE STEADY POINT *** "); printf("Ca = %f, (previous = %f) \n",Ca,Ca_previous); volA /= double((Nx-2)*(Ny-2)*(Nz-2)*nprocs); volB /= double((Nx-2)*(Ny-2)*(Nz-2)*nprocs); FILE * kr_log_file = fopen("relperm.csv","a"); fprintf(kr_log_file,"%i %.5g %.5g %.5g %.5g %.5g %.5g ",timestep-analysis_interval+20,muA,muB,5.796*alpha,Fx,Fy,Fz); fprintf(kr_log_file,"%.5g %.5g %.5g %.5g %.5g %.5g %.5g %.5g\n",volA,volB,vA_x,vA_y,vA_z,vB_x,vB_y,vB_z); fclose(kr_log_file); // flow reversal criteria based on fractional flow if (delta_volume_target < 0.0 && volA*sqrt(vA_x*vA_x + vA_y*vA_y + vA_z*vA_z)/(volB*sqrt(vB_x*vB_x + vB_y*vB_y + vB_z*vB_z)) < RESIDUAL_ENDPOINT_THRESHOLD){ delta_volume_target *= (-1.0); } else if (delta_volume_target > 0.0 && (volB*sqrt(vB_x*vB_x + vB_y*vB_y + vB_z*vB_z)) / (volA*sqrt(vA_x*vA_x + vA_y*vA_y + vA_z*vA_z)) < RESIDUAL_ENDPOINT_THRESHOLD){ delta_volume_target *= (-1.0); } printf(" Measured capillary number %f \n ",Ca); } if (SET_CAPILLARY_NUMBER ){ Fx *= capillary_number / Ca; Fy *= capillary_number / Ca; Fz *= capillary_number / Ca; if (force_magnitude > 1e-3){ Fx *= 1e-3/force_magnitude; // impose ceiling for stability Fy *= 1e-3/force_magnitude; Fz *= 1e-3/force_magnitude; } if (force_magnitude < 1e-6){ Fx *= 1e-6/force_magnitude; // impose floor Fy *= 1e-6/force_magnitude; Fz *= 1e-6/force_magnitude; } if (rank == 0) printf(" -- adjust force by factor %f \n ",capillary_number / Ca); Averages->SetParams(rhoA,rhoB,tauA,tauB,Fx,Fy,Fz,alpha); } } else{ if (rank==0){ printf("** Continue to simulate steady *** \n "); printf("Ca = %f, (previous = %f) \n",Ca,Ca_previous); } morph_timesteps=0; } Ca_previous = Ca; } if (MORPH_ADAPT ){ CURRENT_MORPH_TIMESTEPS += analysis_interval; if (rank==0) printf("***Morphological step with target volume change %f ***\n", delta_volume_target); //double delta_volume_target = volB - (volA + volB)*TARGET_SATURATION; // change in volume to A delta_volume += MorphInit(beta,delta_volume_target-delta_volume); if ( (delta_volume - delta_volume_target)/delta_volume_target > 0.0 ){ MORPH_ADAPT = false; delta_volume = 0.0; } else if (CURRENT_MORPH_TIMESTEPS > MAX_MORPH_TIMESTEPS) { MORPH_ADAPT = false; } MPI_Barrier(comm); } morph_timesteps += analysis_interval; } } analysis.finish(); PROFILE_STOP("Loop"); PROFILE_SAVE("lbpm_color_simulator",1); //************************************************************************ ScaLBL_DeviceBarrier(); MPI_Barrier(comm); stoptime = MPI_Wtime(); if (rank==0) printf("-------------------------------------------------------------------\n"); // Compute the walltime per timestep cputime = (stoptime - starttime)/timestep; // Performance obtained from each node double MLUPS = double(Np)/cputime/1000000; if (rank==0) printf("********************************************************\n"); if (rank==0) printf("CPU time = %f \n", cputime); if (rank==0) printf("Lattice update rate (per core)= %f MLUPS \n", MLUPS); MLUPS *= nprocs; if (rank==0) printf("Lattice update rate (total)= %f MLUPS \n", MLUPS); if (rank==0) printf("********************************************************\n"); // ************************************************************************ } double ScaLBL_ColorModel::MorphInit(const double beta, const double target_delta_volume){ const RankInfoStruct rank_info(rank,nprocx,nprocy,nprocz); double vF = 0.f; double vS = 0.f; double delta_volume; DoubleArray phase(Nx,Ny,Nz); IntArray phase_label(Nx,Ny,Nz);; DoubleArray phase_distance(Nx,Ny,Nz); Array phase_id(Nx,Ny,Nz); fillHalo fillDouble(Dm->Comm,Dm->rank_info,{Nx-2,Ny-2,Nz-2},{1,1,1},0,1); // Basic algorithm to // 1. Copy phase field to CPU ScaLBL_CopyToHost(phase.data(), Phi, N*sizeof(double)); double count,count_global,volume_initial,volume_final,volume_connected; count = 0.f; for (int k=1; k 0.f && Averages->SDs(i,j,k) > 0.f) count+=1.f; } } } volume_initial = sumReduce( Dm->Comm, count); /* sprintf(LocalRankFilename,"phi_initial.%05i.raw",rank); FILE *INPUT = fopen(LocalRankFilename,"wb"); fwrite(phase.data(),8,N,INPUT); fclose(INPUT); */ // 2. Identify connected components of phase field -> phase_label BlobIDstruct new_index; ComputeGlobalBlobIDs(Nx-2,Ny-2,Nz-2,rank_info,phase,Averages->SDs,vF,vS,phase_label,comm); MPI_Barrier(comm); // only operate on component "0" count = 0.0; for (int k=0; kComm, count); // 3. Generate a distance map to the largest object -> phase_distance CalcDist(phase_distance,phase_id,*Dm); double temp,value; double factor=0.5/beta; for (int k=0; k 1.f) value=1.f; if (value < -1.f) value=-1.f; // temp -- distance based on analytical form McClure, Prins et al, Comp. Phys. Comm. temp = -factor*log((1.0+value)/(1.0-value)); /// use this approximation close to the object if (fabs(value) < 0.8 && Averages->SDs(i,j,k) > 1.f ){ phase_distance(i,j,k) = temp; } // erase the original object phase(i,j,k) = -1.0; } } } } if (volume_connected < 0.025*volume_initial){ // if connected volume is less than 2.5% just delete the whole thing if (rank==0) printf("Connected region has shrunk to less than 2.5% of total fluid volume (remove the whole thing) \n"); } else { if (rank==0) printf("MorphGrow with target volume fraction change %f \n", target_delta_volume/volume_initial); double target_delta_volume_incremental = target_delta_volume; if (fabs(target_delta_volume) > 0.01*volume_initial) target_delta_volume_incremental = 0.01*volume_initial*target_delta_volume/fabs(target_delta_volume); delta_volume = MorphGrow(Averages->SDs,phase_distance,phase_id,Averages->Dm, target_delta_volume_incremental); for (int k=0; kSDs(i,j,k) > 0.f){ if (d < 3.f){ //phase(i,j,k) = -1.0; phase(i,j,k) = (2.f*(exp(-2.f*beta*d))/(1.f+exp(-2.f*beta*d))-1.f); } } } } } fillDouble.fill(phase); } count = 0.f; for (int k=1; k 0.f && Averages->SDs(i,j,k) > 0.f) count+=1.f; } } } volume_final= sumReduce( Dm->Comm, count); delta_volume = (volume_final-volume_initial); if (rank == 0) printf("MorphInit: change fluid volume fraction by %f \n", delta_volume/volume_initial); if (rank == 0) printf(" new saturation = %f \n", volume_final/(0.238323*double((Nx-2)*(Ny-2)*(Nz-2)*nprocs))); // 6. copy back to the device //if (rank==0) printf("MorphInit: copy data back to device\n"); ScaLBL_CopyToDevice(Phi,phase.data(),N*sizeof(double)); /* sprintf(LocalRankFilename,"dist_final.%05i.raw",rank); FILE *DIST = fopen(LocalRankFilename,"wb"); fwrite(phase_distance.data(),8,N,DIST); fclose(DIST); sprintf(LocalRankFilename,"phi_final.%05i.raw",rank); FILE *PHI = fopen(LocalRankFilename,"wb"); fwrite(phase.data(),8,N,PHI); fclose(PHI); */ // 7. Re-initialize phase field and density ScaLBL_PhaseField_Init(dvcMap, Phi, Den, Aq, Bq, 0, ScaLBL_Comm->LastExterior(), Np); ScaLBL_PhaseField_Init(dvcMap, Phi, Den, Aq, Bq, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np); if (BoundaryCondition >0 ){ if (Dm->kproc()==0){ ScaLBL_SetSlice_z(Phi,1.0,Nx,Ny,Nz,0); ScaLBL_SetSlice_z(Phi,1.0,Nx,Ny,Nz,1); ScaLBL_SetSlice_z(Phi,1.0,Nx,Ny,Nz,2); } if (Dm->kproc() == nprocz-1){ ScaLBL_SetSlice_z(Phi,-1.0,Nx,Ny,Nz,Nz-1); ScaLBL_SetSlice_z(Phi,-1.0,Nx,Ny,Nz,Nz-2); ScaLBL_SetSlice_z(Phi,-1.0,Nx,Ny,Nz,Nz-3); } } return delta_volume; } void ScaLBL_ColorModel::WriteDebug(){ // Copy back final phase indicator field and convert to regular layout DoubleArray PhaseField(Nx,Ny,Nz); //ScaLBL_Comm->RegularLayout(Map,Phi,PhaseField); ScaLBL_CopyToHost(PhaseField.data(), Phi, sizeof(double)*N); FILE *OUTFILE; sprintf(LocalRankFilename,"Phase.%05i.raw",rank); OUTFILE = fopen(LocalRankFilename,"wb"); fwrite(PhaseField.data(),8,N,OUTFILE); fclose(OUTFILE); }