//************************************************************************* // Lattice Boltzmann Simulator for Single Phase Flow in Porous Media // James E. McCLure //************************************************************************* #include #include #include #include "common/ScaLBL.h" #include "common/MPI_Helpers.h" #include "models/ColorModel.h" std::shared_ptr loadInputs( int nprocs ) { auto db = std::make_shared( "Domain.in" ); const int dim = 50; db->putScalar( "BC", 0 ); if ( nprocs == 1 ){ db->putVector( "nproc", { 1, 1, 1 } ); db->putVector( "n", { 3, 1, 1 } ); db->putScalar( "nspheres", 0 ); db->putVector( "L", { 1, 1, 1 } ); } else if ( nprocs == 2 ) { db->putVector( "nproc", { 2, 1, 1 } ); db->putVector( "n", { dim, dim, dim } ); db->putScalar( "nspheres", 0 ); db->putVector( "L", { 1, 1, 1 } ); } else if ( nprocs == 4 ) { db->putVector( "nproc", { 2, 2, 1 } ); db->putVector( "n", { dim, dim, dim } ); db->putScalar( "nspheres", 0 ); db->putVector( "L", { 1, 1, 1 } ); } else if (nprocs==8){ db->putVector( "nproc", { 2, 2, 2 } ); db->putVector( "n", { dim, dim, dim } ); db->putScalar( "nspheres", 0 ); db->putVector( "L", { 1, 1, 1 } ); } return db; } void InitializeSquareTube(ScaLBL_ColorModel &ColorModel){ int i,j,k,n; int rank = ColorModel.Mask->rank(); int Nx = ColorModel.Mask->Nx; int Ny = ColorModel.Mask->Ny; int Nz = ColorModel.Mask->Nz; int nprocx = ColorModel.Mask->rank_info.nx; int nprocy = ColorModel.Mask->rank_info.ny; int iproc = ColorModel.Mask->rank_info.ix; int jproc = ColorModel.Mask->rank_info.jy; int kproc = ColorModel.Mask->rank_info.kz; for (k=0;kid[n]=0; } } } printf("rank=%i, %i,%i,%i \n",rank,iproc,jproc,kproc); // Initialize a square tube for (k=1;kid[n]=0; else if (iglobal > (Nx-2)*nprocx-2) ColorModel.Mask->id[n]=0; else if (jglobal < 2) ColorModel.Mask->id[n]=0; else if (jglobal > (Ny-2)*nprocy-2) ColorModel.Mask->id[n]=0; else if (kglobal < 20) ColorModel.Mask->id[n]=1; else ColorModel.Mask->id[n]=2; } } } } //*************************************************************************************** int main(int argc, char **argv) { //***************************************** // ***** MPI STUFF **************** //***************************************** // Initialize MPI int rank,nprocs; MPI_Init(&argc,&argv); MPI_Comm comm = MPI_COMM_WORLD; MPI_Comm_rank(comm,&rank); MPI_Comm_size(comm,&nprocs); int check=0; { if (rank == 0){ printf("********************************************************\n"); printf("Running Color Model: TestColor \n"); printf("********************************************************\n"); } auto filename = argv[1]; ScaLBL_ColorModel ColorModel(rank,nprocs,comm); ColorModel.ReadParams(filename); ColorModel.SetDomain(); //ColorModel.ReadInput(); InitializeSquareTube(ColorModel); ColorModel.Create(); // creating the model will create data structure to match the pore structure and allocate variables ColorModel.Initialize(); // initializing the model will set initial conditions for variables ColorModel.Run(); ColorModel.WriteDebug(); } // **************************************************** MPI_Barrier(comm); MPI_Finalize(); // **************************************************** return check; }