#include #include #include #include #include #include #include #include "Domain.h" #include "Extras.h" #include "ScaLBL.h" #include "Communication.h" #include "TwoPhase.h" #include "common/MPI_Helpers.h" //#define WRITE_SURFACES /* * Simulator for two-phase flow in porous media * James E. McClure 2013-2014 */ using namespace std; //************************************************************************* // Implementation of Steady State Single-Phase LBM for permeability measurement //************************************************************************* 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> tau; // Viscosity parameter // Line 2: External force components (Fx,Fy, Fz) input >> Fx; input >> Fy; input >> Fz; // Line 3: Pressure Boundary conditions input >> Restart; input >> pBC; input >> din; input >> dout; // Line 4: time-stepping criteria input >> timestepMax; // max no. of timesteps input >> interval; // restart interval input >> tol; // error tolerance //............................................................. //....................................................................... // Reading the domain information file //....................................................................... ifstream domain("Domain.in"); domain >> nprocx; domain >> nprocy; domain >> nprocz; domain >> Nx; domain >> Ny; domain >> Nz; domain >> nspheres; domain >> Lx; domain >> Ly; domain >> Lz; //....................................................................... } // ************************************************************** // Broadcast simulation parameters from rank 0 to all other procs MPI_Barrier(MPI_COMM_WORLD); //................................................. MPI_Bcast(&tau,1,MPI_DOUBLE,0,MPI_COMM_WORLD); MPI_Bcast(&pBC,1,MPI_LOGICAL,0,MPI_COMM_WORLD); MPI_Bcast(&Restart,1,MPI_LOGICAL,0,MPI_COMM_WORLD); MPI_Bcast(&din,1,MPI_DOUBLE,0,MPI_COMM_WORLD); MPI_Bcast(&dout,1,MPI_DOUBLE,0,MPI_COMM_WORLD); MPI_Bcast(&Fx,1,MPI_DOUBLE,0,MPI_COMM_WORLD); MPI_Bcast(&Fy,1,MPI_DOUBLE,0,MPI_COMM_WORLD); MPI_Bcast(&Fz,1,MPI_DOUBLE,0,MPI_COMM_WORLD); MPI_Bcast(×tepMax,1,MPI_INT,0,MPI_COMM_WORLD); MPI_Bcast(&interval,1,MPI_INT,0,MPI_COMM_WORLD); MPI_Bcast(&tol,1,MPI_DOUBLE,0,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); RESTART_INTERVAL=interval; // ************************************************************** // ************************************************************** double rlxA = 1.f/tau; double rlxB = 8.f*(2.f-rlxA)/(8.f-rlxA); 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("tau = %f \n", tau); printf("Force(x) = %f \n", Fx); printf("Force(y) = %f \n", Fy); printf("Force(z) = %f \n", Fz); printf("Sub-domain size = %i x %i x %i\n",Nz,Nz,Nz); printf("Parallel domain size = %i x %i x %i\n",nprocx,nprocy,nprocz); printf("********************************************************\n"); } double viscosity=(tau-0.5)/3.0; // Initialized domain and averaging framework for Two-Phase Flow int BC=pBC; Domain Dm(Nx,Ny,Nz,rank,nprocx,nprocy,nprocz,Lx,Ly,Lz,BC); TwoPhase Averages(Dm); InitializeRanks( rank, nprocx, nprocy, nprocz, iproc, jproc, kproc, rank_x, rank_y, rank_z, rank_X, rank_Y, rank_Z, rank_xy, rank_XY, rank_xY, rank_Xy, rank_xz, rank_XZ, rank_xZ, rank_Xz, rank_yz, rank_YZ, rank_yZ, rank_Yz ); MPI_Barrier(MPI_COMM_WORLD); Nz += 2; Nx = Ny = Nz; // Cubic domain int N = Nx*Ny*Nz; int dist_mem_size = N*sizeof(double); //....................................................................... if (rank == 0) printf("Read input media... \n"); //....................................................................... //....................................................................... // Filenames used char LocalRankString[8]; char LocalRankFilename[40]; char LocalRestartFile[40]; char tmpstr[10]; 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); if (pBC) iVol_global = 1.0/(1.0*(Nx-2)*nprocx*(Ny-2)*nprocy*((Nz-2)*nprocz-6)); double porosity, pore_vol; //........................................................................... if (rank == 0) cout << "Reading in domain from signed distance function..." << endl; //....................................................................... sprintf(LocalRankString,"%05d",rank); // sprintf(LocalRankFilename,"%s%s","ID.",LocalRankString); // WriteLocalSolidID(LocalRankFilename, id, N); sprintf(LocalRankFilename,"%s%s","SignDist.",LocalRankString); ReadBinaryFile(LocalRankFilename, Averages.SDs.get(), N); MPI_Barrier(MPI_COMM_WORLD); if (rank == 0) cout << "Domain set." << endl; //....................................................................... // Assign the phase ID field based on the signed distance //....................................................................... for (k=0;k 0.0){ id[n] = 2; } // compute the porosity (actual interface location used) if (Averages.SDs(n) > 0.0){ sum++; } } } } // Set up kstart, kfinish so that the reservoirs are excluded from averaging int kstart,kfinish; kstart = 1; kfinish = Nz-1; if (pBC && kproc==0) kstart = 4; if (pBC && kproc==nprocz-1) kfinish = Nz-4; // Compute the pore volume sum_local = 0.0; for ( k=kstart;k 0){ sum_local += 1.0; } } } } MPI_Allreduce(&sum_local,&pore_vol,1,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD); // MPI_Allreduce(&sum_local,&porosity,1,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD); porosity = pore_vol*iVol_global; if (rank==0) printf("Media porosity = %f \n",porosity); //......................................................... // If pressure boundary conditions are applied remove solid if (pBC && kproc == 0){ for (k=0; k<3; k++){ for (j=0;j tol ){ //************************************************************************* // Fused Color Gradient and Collision //************************************************************************* MRT( ID,f_even,f_odd,rlxA,rlxB,Fx,Fy,Fz,Nx,Ny,Nz); //************************************************************************* //................................................................................... PackDist(1,dvcSendList_x,0,sendCount_x,sendbuf_x,f_even,N); PackDist(4,dvcSendList_x,sendCount_x,sendCount_x,sendbuf_x,f_even,N); PackDist(5,dvcSendList_x,2*sendCount_x,sendCount_x,sendbuf_x,f_even,N); PackDist(6,dvcSendList_x,3*sendCount_x,sendCount_x,sendbuf_x,f_even,N); PackDist(7,dvcSendList_x,4*sendCount_x,sendCount_x,sendbuf_x,f_even,N); //...Packing for X face(1,7,9,11,13)................................ PackDist(0,dvcSendList_X,0,sendCount_X,sendbuf_X,f_odd,N); PackDist(3,dvcSendList_X,sendCount_X,sendCount_X,sendbuf_X,f_odd,N); PackDist(4,dvcSendList_X,2*sendCount_X,sendCount_X,sendbuf_X,f_odd,N); PackDist(5,dvcSendList_X,3*sendCount_X,sendCount_X,sendbuf_X,f_odd,N); PackDist(6,dvcSendList_X,4*sendCount_X,sendCount_X,sendbuf_X,f_odd,N); //...Packing for y face(4,8,9,16,18)................................. PackDist(2,dvcSendList_y,0,sendCount_y,sendbuf_y,f_even,N); PackDist(4,dvcSendList_y,sendCount_y,sendCount_y,sendbuf_y,f_even,N); PackDist(4,dvcSendList_y,2*sendCount_y,sendCount_y,sendbuf_y,f_odd,N); PackDist(8,dvcSendList_y,3*sendCount_y,sendCount_y,sendbuf_y,f_even,N); PackDist(9,dvcSendList_y,4*sendCount_y,sendCount_y,sendbuf_y,f_even,N); //...Packing for Y face(3,7,10,15,17)................................. PackDist(1,dvcSendList_Y,0,sendCount_Y,sendbuf_Y,f_odd,N); PackDist(3,dvcSendList_Y,sendCount_Y,sendCount_Y,sendbuf_Y,f_odd,N); PackDist(5,dvcSendList_Y,2*sendCount_Y,sendCount_Y,sendbuf_Y,f_even,N); PackDist(7,dvcSendList_Y,3*sendCount_Y,sendCount_Y,sendbuf_Y,f_odd,N); PackDist(8,dvcSendList_Y,4*sendCount_Y,sendCount_Y,sendbuf_Y,f_odd,N); //...Packing for z face(6,12,13,16,17)................................ PackDist(3,dvcSendList_z,0,sendCount_z,sendbuf_z,f_even,N); PackDist(6,dvcSendList_z,sendCount_z,sendCount_z,sendbuf_z,f_even,N); PackDist(6,dvcSendList_z,2*sendCount_z,sendCount_z,sendbuf_z,f_odd,N); PackDist(8,dvcSendList_z,3*sendCount_z,sendCount_z,sendbuf_z,f_even,N); PackDist(8,dvcSendList_z,4*sendCount_z,sendCount_z,sendbuf_z,f_odd,N); //...Packing for Z face(5,11,14,15,18)................................ PackDist(2,dvcSendList_Z,0,sendCount_Z,sendbuf_Z,f_odd,N); PackDist(5,dvcSendList_Z,sendCount_Z,sendCount_Z,sendbuf_Z,f_odd,N); PackDist(7,dvcSendList_Z,2*sendCount_Z,sendCount_Z,sendbuf_Z,f_even,N); PackDist(7,dvcSendList_Z,3*sendCount_Z,sendCount_Z,sendbuf_Z,f_odd,N); PackDist(9,dvcSendList_Z,4*sendCount_Z,sendCount_Z,sendbuf_Z,f_even,N); //...Pack the xy edge (8)................................ PackDist(4,dvcSendList_xy,0,sendCount_xy,sendbuf_xy,f_even,N); //...Pack the Xy edge (9)................................ PackDist(4,dvcSendList_Xy,0,sendCount_Xy,sendbuf_Xy,f_odd,N); //...Pack the xY edge (10)................................ PackDist(5,dvcSendList_xY,0,sendCount_xY,sendbuf_xY,f_even,N); //...Pack the XY edge (7)................................ PackDist(3,dvcSendList_XY,0,sendCount_XY,sendbuf_XY,f_odd,N); //...Pack the xz edge (12)................................ PackDist(6,dvcSendList_xz,0,sendCount_xz,sendbuf_xz,f_even,N); //...Pack the xZ edge (14)................................ PackDist(7,dvcSendList_xZ,0,sendCount_xZ,sendbuf_xZ,f_even,N); //...Pack the Xz edge (13)................................ PackDist(6,dvcSendList_Xz,0,sendCount_Xz,sendbuf_Xz,f_odd,N); //...Pack the XZ edge (11)................................ PackDist(5,dvcSendList_XZ,0,sendCount_XZ,sendbuf_XZ,f_odd,N); //...Pack the xz edge (12)................................ //...Pack the yz edge (16)................................ PackDist(8,dvcSendList_yz,0,sendCount_yz,sendbuf_yz,f_even,N); //...Pack the yZ edge (18)................................ PackDist(9,dvcSendList_yZ,0,sendCount_yZ,sendbuf_yZ,f_even,N); //...Pack the Yz edge (17)................................ PackDist(8,dvcSendList_Yz,0,sendCount_Yz,sendbuf_Yz,f_odd,N); //...Pack the YZ edge (15)................................ PackDist(7,dvcSendList_YZ,0,sendCount_YZ,sendbuf_YZ,f_odd,N); //................................................................................... //................................................................................... // Send all the distributions MPI_Isend(sendbuf_x, 5*sendCount_x,MPI_DOUBLE,rank_x,sendtag,MPI_COMM_WORLD,&req1[0]); MPI_Irecv(recvbuf_X, 5*recvCount_X,MPI_DOUBLE,rank_X,recvtag,MPI_COMM_WORLD,&req2[0]); MPI_Isend(sendbuf_X, 5*sendCount_X,MPI_DOUBLE,rank_X,sendtag,MPI_COMM_WORLD,&req1[1]); MPI_Irecv(recvbuf_x, 5*recvCount_x,MPI_DOUBLE,rank_x,recvtag,MPI_COMM_WORLD,&req2[1]); MPI_Isend(sendbuf_y, 5*sendCount_y,MPI_DOUBLE,rank_y,sendtag,MPI_COMM_WORLD,&req1[2]); MPI_Irecv(recvbuf_Y, 5*recvCount_Y,MPI_DOUBLE,rank_Y,recvtag,MPI_COMM_WORLD,&req2[2]); MPI_Isend(sendbuf_Y, 5*sendCount_Y,MPI_DOUBLE,rank_Y,sendtag,MPI_COMM_WORLD,&req1[3]); MPI_Irecv(recvbuf_y, 5*recvCount_y,MPI_DOUBLE,rank_y,recvtag,MPI_COMM_WORLD,&req2[3]); MPI_Isend(sendbuf_z, 5*sendCount_z,MPI_DOUBLE,rank_z,sendtag,MPI_COMM_WORLD,&req1[4]); MPI_Irecv(recvbuf_Z, 5*recvCount_Z,MPI_DOUBLE,rank_Z,recvtag,MPI_COMM_WORLD,&req2[4]); MPI_Isend(sendbuf_Z, 5*sendCount_Z,MPI_DOUBLE,rank_Z,sendtag,MPI_COMM_WORLD,&req1[5]); MPI_Irecv(recvbuf_z, 5*recvCount_z,MPI_DOUBLE,rank_z,recvtag,MPI_COMM_WORLD,&req2[5]); MPI_Isend(sendbuf_xy, sendCount_xy,MPI_DOUBLE,rank_xy,sendtag,MPI_COMM_WORLD,&req1[6]); MPI_Irecv(recvbuf_XY, recvCount_XY,MPI_DOUBLE,rank_XY,recvtag,MPI_COMM_WORLD,&req2[6]); MPI_Isend(sendbuf_XY, sendCount_XY,MPI_DOUBLE,rank_XY,sendtag,MPI_COMM_WORLD,&req1[7]); MPI_Irecv(recvbuf_xy, recvCount_xy,MPI_DOUBLE,rank_xy,recvtag,MPI_COMM_WORLD,&req2[7]); MPI_Isend(sendbuf_Xy, sendCount_Xy,MPI_DOUBLE,rank_Xy,sendtag,MPI_COMM_WORLD,&req1[8]); MPI_Irecv(recvbuf_xY, recvCount_xY,MPI_DOUBLE,rank_xY,recvtag,MPI_COMM_WORLD,&req2[8]); MPI_Isend(sendbuf_xY, sendCount_xY,MPI_DOUBLE,rank_xY,sendtag,MPI_COMM_WORLD,&req1[9]); MPI_Irecv(recvbuf_Xy, recvCount_Xy,MPI_DOUBLE,rank_Xy,recvtag,MPI_COMM_WORLD,&req2[9]); MPI_Isend(sendbuf_xz, sendCount_xz,MPI_DOUBLE,rank_xz,sendtag,MPI_COMM_WORLD,&req1[10]); MPI_Irecv(recvbuf_XZ, recvCount_XZ,MPI_DOUBLE,rank_XZ,recvtag,MPI_COMM_WORLD,&req2[10]); MPI_Isend(sendbuf_XZ, sendCount_XZ,MPI_DOUBLE,rank_XZ,sendtag,MPI_COMM_WORLD,&req1[11]); MPI_Irecv(recvbuf_xz, recvCount_xz,MPI_DOUBLE,rank_xz,recvtag,MPI_COMM_WORLD,&req2[11]); MPI_Isend(sendbuf_Xz, sendCount_Xz,MPI_DOUBLE,rank_Xz,sendtag,MPI_COMM_WORLD,&req1[12]); MPI_Irecv(recvbuf_xZ, recvCount_xZ,MPI_DOUBLE,rank_xZ,recvtag,MPI_COMM_WORLD,&req2[12]); MPI_Isend(sendbuf_xZ, sendCount_xZ,MPI_DOUBLE,rank_xZ,sendtag,MPI_COMM_WORLD,&req1[13]); MPI_Irecv(recvbuf_Xz, recvCount_Xz,MPI_DOUBLE,rank_Xz,recvtag,MPI_COMM_WORLD,&req2[13]); MPI_Isend(sendbuf_yz, sendCount_yz,MPI_DOUBLE,rank_yz,sendtag,MPI_COMM_WORLD,&req1[14]); MPI_Irecv(recvbuf_YZ, recvCount_YZ,MPI_DOUBLE,rank_YZ,recvtag,MPI_COMM_WORLD,&req2[14]); MPI_Isend(sendbuf_YZ, sendCount_YZ,MPI_DOUBLE,rank_YZ,sendtag,MPI_COMM_WORLD,&req1[15]); MPI_Irecv(recvbuf_yz, recvCount_yz,MPI_DOUBLE,rank_yz,recvtag,MPI_COMM_WORLD,&req2[15]); MPI_Isend(sendbuf_Yz, sendCount_Yz,MPI_DOUBLE,rank_Yz,sendtag,MPI_COMM_WORLD,&req1[16]); MPI_Irecv(recvbuf_yZ, recvCount_yZ,MPI_DOUBLE,rank_yZ,recvtag,MPI_COMM_WORLD,&req2[16]); MPI_Isend(sendbuf_yZ, sendCount_yZ,MPI_DOUBLE,rank_yZ,sendtag,MPI_COMM_WORLD,&req1[17]); MPI_Irecv(recvbuf_Yz, recvCount_Yz,MPI_DOUBLE,rank_Yz,recvtag,MPI_COMM_WORLD,&req2[17]); //................................................................................... //************************************************************************* // Swap the distributions for momentum transport //************************************************************************* SwapD3Q19(ID, f_even, f_odd, Nx, Ny, Nz); //************************************************************************* //................................................................................... // Wait for completion of D3Q19 communication MPI_Waitall(18,req1,stat1); MPI_Waitall(18,req2,stat2); //................................................................................... // Unpack the distributions on the device //................................................................................... //...Map recieve list for the X face: q=2,8,10,12,13 ................................. UnpackDist(0,-1,0,0,dvcRecvList_X,0,recvCount_X,recvbuf_X,f_odd,Nx,Ny,Nz); UnpackDist(3,-1,-1,0,dvcRecvList_X,recvCount_X,recvCount_X,recvbuf_X,f_odd,Nx,Ny,Nz); UnpackDist(4,-1,1,0,dvcRecvList_X,2*recvCount_X,recvCount_X,recvbuf_X,f_odd,Nx,Ny,Nz); UnpackDist(5,-1,0,-1,dvcRecvList_X,3*recvCount_X,recvCount_X,recvbuf_X,f_odd,Nx,Ny,Nz); UnpackDist(6,-1,0,1,dvcRecvList_X,4*recvCount_X,recvCount_X,recvbuf_X,f_odd,Nx,Ny,Nz); //................................................................................... //...Map recieve list for the x face: q=1,7,9,11,13.................................. UnpackDist(1,1,0,0,dvcRecvList_x,0,recvCount_x,recvbuf_x,f_even,Nx,Ny,Nz); UnpackDist(4,1,1,0,dvcRecvList_x,recvCount_x,recvCount_x,recvbuf_x,f_even,Nx,Ny,Nz); UnpackDist(5,1,-1,0,dvcRecvList_x,2*recvCount_x,recvCount_x,recvbuf_x,f_even,Nx,Ny,Nz); UnpackDist(6,1,0,1,dvcRecvList_x,3*recvCount_x,recvCount_x,recvbuf_x,f_even,Nx,Ny,Nz); UnpackDist(7,1,0,-1,dvcRecvList_x,4*recvCount_x,recvCount_x,recvbuf_x,f_even,Nx,Ny,Nz); //................................................................................... //...Map recieve list for the y face: q=4,8,9,16,18 ................................... UnpackDist(1,0,-1,0,dvcRecvList_Y,0,recvCount_Y,recvbuf_Y,f_odd,Nx,Ny,Nz); UnpackDist(3,-1,-1,0,dvcRecvList_Y,recvCount_Y,recvCount_Y,recvbuf_Y,f_odd,Nx,Ny,Nz); UnpackDist(5,1,-1,0,dvcRecvList_Y,2*recvCount_Y,recvCount_Y,recvbuf_Y,f_even,Nx,Ny,Nz); UnpackDist(7,0,-1,-1,dvcRecvList_Y,3*recvCount_Y,recvCount_Y,recvbuf_Y,f_odd,Nx,Ny,Nz); UnpackDist(8,0,-1,1,dvcRecvList_Y,4*recvCount_Y,recvCount_Y,recvbuf_Y,f_odd,Nx,Ny,Nz); //................................................................................... //...Map recieve list for the Y face: q=3,7,10,15,17 .................................. UnpackDist(2,0,1,0,dvcRecvList_y,0,recvCount_y,recvbuf_y,f_even,Nx,Ny,Nz); UnpackDist(4,1,1,0,dvcRecvList_y,recvCount_y,recvCount_y,recvbuf_y,f_even,Nx,Ny,Nz); UnpackDist(4,-1,1,0,dvcRecvList_y,2*recvCount_y,recvCount_y,recvbuf_y,f_odd,Nx,Ny,Nz); UnpackDist(8,0,1,1,dvcRecvList_y,3*recvCount_y,recvCount_y,recvbuf_y,f_even,Nx,Ny,Nz); UnpackDist(9,0,1,-1,dvcRecvList_y,4*recvCount_y,recvCount_y,recvbuf_y,f_even,Nx,Ny,Nz); //................................................................................... //...Map recieve list for the z face<<<6,12,13,16,17).............................................. UnpackDist(2,0,0,-1,dvcRecvList_Z,0,recvCount_Z,recvbuf_Z,f_odd,Nx,Ny,Nz); UnpackDist(5,-1,0,-1,dvcRecvList_Z,recvCount_Z,recvCount_Z,recvbuf_Z,f_odd,Nx,Ny,Nz); UnpackDist(7,1,0,-1,dvcRecvList_Z,2*recvCount_Z,recvCount_Z,recvbuf_Z,f_even,Nx,Ny,Nz); UnpackDist(7,0,-1,-1,dvcRecvList_Z,3*recvCount_Z,recvCount_Z,recvbuf_Z,f_odd,Nx,Ny,Nz); UnpackDist(9,0,1,-1,dvcRecvList_Z,4*recvCount_Z,recvCount_Z,recvbuf_Z,f_even,Nx,Ny,Nz); //...Map recieve list for the Z face<<<5,11,14,15,18).............................................. UnpackDist(3,0,0,1,dvcRecvList_z,0,recvCount_z,recvbuf_z,f_even,Nx,Ny,Nz); UnpackDist(6,1,0,1,dvcRecvList_z,recvCount_z,recvCount_z,recvbuf_z,f_even,Nx,Ny,Nz); UnpackDist(6,-1,0,1,dvcRecvList_z,2*recvCount_z,recvCount_z,recvbuf_z,f_odd,Nx,Ny,Nz); UnpackDist(8,0,1,1,dvcRecvList_z,3*recvCount_z,recvCount_z,recvbuf_z,f_even,Nx,Ny,Nz); UnpackDist(8,0,-1,1,dvcRecvList_z,4*recvCount_z,recvCount_z,recvbuf_z,f_odd,Nx,Ny,Nz); //.................................................................................. //...Map recieve list for the xy edge <<<8)................................ UnpackDist(3,-1,-1,0,dvcRecvList_XY,0,recvCount_XY,recvbuf_XY,f_odd,Nx,Ny,Nz); //...Map recieve list for the Xy edge <<<9)................................ UnpackDist(5,1,-1,0,dvcRecvList_xY,0,recvCount_xY,recvbuf_xY,f_even,Nx,Ny,Nz); //...Map recieve list for the xY edge <<<10)................................ UnpackDist(4,-1,1,0,dvcRecvList_Xy,0,recvCount_Xy,recvbuf_Xy,f_odd,Nx,Ny,Nz); //...Map recieve list for the XY edge <<<7)................................ UnpackDist(4,1,1,0,dvcRecvList_xy,0,recvCount_xy,recvbuf_xy,f_even,Nx,Ny,Nz); //...Map recieve list for the xz edge <<<12)................................ UnpackDist(5,-1,0,-1,dvcRecvList_XZ,0,recvCount_XZ,recvbuf_XZ,f_odd,Nx,Ny,Nz); //...Map recieve list for the xZ edge <<<14)................................ UnpackDist(6,-1,0,1,dvcRecvList_Xz,0,recvCount_Xz,recvbuf_Xz,f_odd,Nx,Ny,Nz); //...Map recieve list for the Xz edge <<<13)................................ UnpackDist(7,1,0,-1,dvcRecvList_xZ,0,recvCount_xZ,recvbuf_xZ,f_even,Nx,Ny,Nz); //...Map recieve list for the XZ edge <<<11)................................ UnpackDist(6,1,0,1,dvcRecvList_xz,0,recvCount_xz,recvbuf_xz,f_even,Nx,Ny,Nz); //...Map recieve list for the yz edge <<<16)................................ UnpackDist(7,0,-1,-1,dvcRecvList_YZ,0,recvCount_YZ,recvbuf_YZ,f_odd,Nx,Ny,Nz); //...Map recieve list for the yZ edge <<<18)....................`............ UnpackDist(8,0,-1,1,dvcRecvList_Yz,0,recvCount_Yz,recvbuf_Yz,f_odd,Nx,Ny,Nz); //...Map recieve list for the Yz edge <<<17)................................ UnpackDist(9,0,1,-1,dvcRecvList_yZ,0,recvCount_yZ,recvbuf_yZ,f_even,Nx,Ny,Nz); //...Map recieve list for the YZ edge <<<15)................................ UnpackDist(8,0,1,1,dvcRecvList_yz,0,recvCount_yz,recvbuf_yz,f_even,Nx,Ny,Nz); //................................................................................... if (pBC && kproc == 0) { PressureBC_inlet(f_even,f_odd,din,Nx,Ny,Nz); } if (pBC && kproc == nprocz-1){ PressureBC_outlet(f_even,f_odd,dout,Nx,Ny,Nz,Nx*Ny*(Nz-2)); } //................................................................................... MPI_Barrier(MPI_COMM_WORLD); // Timestep completed! timestep++; if (timestep%500 == 0){ //........................................................................... // Copy the data for for the analysis timestep //........................................................................... // Copy the phase from the GPU -> CPU //........................................................................... DeviceBarrier(); ComputePressureD3Q19(ID,f_even,f_odd,Pressure,Nx,Ny,Nz); ComputeVelocityD3Q19(ID,f_even,f_odd,Velocity,Nx,Ny,Nz); CopyToHost(Averages.Press.get(),Pressure,N*sizeof(double)); CopyToHost(Averages.Vel_x.get(),&Velocity[0],N*sizeof(double)); CopyToHost(Averages.Vel_y.get(),&Velocity[N],N*sizeof(double)); CopyToHost(Averages.Vel_z.get(),&Velocity[2*N],N*sizeof(double)); // copy the velocity CopyToHost(Vel,Velocity,3*N*sizeof(double)); // Compute the macroscopic fluid velocity double vawx,vawy,vawz,Vw; double vawx_global,vawy_global,vawz_global,Vw_global; vawx=vawy=vawz=0.0; for (k=1; k0){ vawx+=Vel[n]; vawy+=Vel[N+n]; vawz+=Vel[2*N+n]; Vw+=1.0; } } } } // Reduce to get the global averages MPI_Barrier(MPI_COMM_WORLD); MPI_Allreduce(&vawx,&vawx_global,1,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD); MPI_Allreduce(&vawy,&vawy_global,1,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD); MPI_Allreduce(&vawz,&vawz_global,1,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD); MPI_Allreduce(&Vw,&Vw_global,1,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD); MPI_Barrier(MPI_COMM_WORLD); // Compute volume averaged velocity vawx = vawx_global/Vw_global; vawy = vawy_global/Vw_global; vawz = vawz_global/Vw_global; // Way more work than necessary -- this is just to get the solid interfacial area!! Averages.Initialize(); Averages.UpdateMeshValues(); Averages.ComputeLocal(); Averages.Reduce(); // Check for steady macroscopic velocity velocity = vawx*Fx/mag_force + vawy*Fy/mag_force + vawz*Fz/mag_force; err=fabs(velocity-vel_prev)/velocity; if (rank==0){ // ************* DIMENSIONLESS FORCHEIMER EQUATION ************************* // Dye, A.L., McClure, J.E., Gray, W.G. and C.T. Miller // Description of Non-Darcy Flows in Porous Medium Systems // Physical Review E 87 (3), 033012 // Fo := density*D32^3*(density*force) / (viscosity^2) // Re := density*D32*velocity / viscosity // Fo = a*Re + b*Re^2 // ************************************************************************* //viscosity = (tau-0.5)*0.333333333333333333; D32 = 6.0*(Dm.Volume-Vw_global)/Averages.As_global; printf("Sauter Mean Diameter = %f \n",D32); mag_force = sqrt(Fx*Fx+Fy*Fy+Fz*Fz); Fo = D32*D32*D32*mag_force/viscosity/viscosity; // .... 1-D flow should be aligned with force ... velocity = vawx*Fx/mag_force + vawy*Fy/mag_force + vawz*Fz/mag_force; err1D = fabs(velocity-sqrt(vawx*vawx+vawy*vawy+vawz*vawz))/velocity; //.......... Computation of the Reynolds number Re .............. Re = D32*velocity/viscosity; printf("Relative error for 1D representation: %f \n",err1D); printf("Dimensionless force: %f \n", Fo); printf("Reynolds number: %f \n", Re); printf("Permeability: %f \n", Re/Fo); } } } //************************************************************************/ DeviceBarrier(); MPI_Barrier(MPI_COMM_WORLD); 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(Nx*Ny*Nz)/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"); // **************************************************** MPI_Barrier(MPI_COMM_WORLD); MPI_Finalize(); // **************************************************** }