Files
LBPM/models/FreeLeeModel.cpp

799 lines
26 KiB
C++
Raw Normal View History

2020-09-25 16:18:54 -04:00
/*
color lattice boltzmann model
*/
2020-09-29 11:05:09 -04:00
#include "models/FreeLeeModel.h"
2020-09-25 16:18:54 -04:00
#include "analysis/distance.h"
#include "analysis/morphology.h"
#include "common/Communication.h"
#include "common/ReadMicroCT.h"
#include <stdlib.h>
#include <time.h>
2021-01-13 22:06:08 -05:00
ScaLBL_FreeLeeModel::ScaLBL_FreeLeeModel(int RANK, int NP, const Utilities::MPI& COMM):
2020-09-25 16:18:54 -04:00
rank(RANK), nprocs(NP), Restart(0),timestep(0),timestepMax(0),tauA(0),tauB(0),rhoA(0),rhoB(0),W(0),gamma(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)
{
2020-09-29 13:37:02 -04:00
}
2020-09-25 16:18:54 -04:00
ScaLBL_FreeLeeModel::~ScaLBL_FreeLeeModel(){
}
void ScaLBL_FreeLeeModel::ReadParams(string filename){
// read the input database
db = std::make_shared<Database>( filename );
domain_db = db->getDatabase( "Domain" );
freelee_db = db->getDatabase( "FreeLee" );
analysis_db = db->getDatabase( "Analysis" );
vis_db = db->getDatabase( "Visualization" );
// set defaults
timestepMax = 100000;
tauA = tauB = 1.0;
rhoA = rhoB = 1.0;
Fx = Fy = Fz = 0.0;
2021-01-18 21:30:27 -05:00
gamma=1e-3;//surface tension
W=5.0;//interfacial thickness
2020-09-25 16:18:54 -04:00
Restart=false;
din=dout=1.0;
flux=0.0;
// Color Model parameters
if (freelee_db->keyExists( "timestepMax" )){
timestepMax = freelee_db->getScalar<int>( "timestepMax" );
}
if (freelee_db->keyExists( "tauA" )){
tauA = freelee_db->getScalar<double>( "tauA" );
}
if (freelee_db->keyExists( "tauB" )){
tauB = freelee_db->getScalar<double>( "tauB" );
}
if (freelee_db->keyExists( "rhoA" )){
rhoA = freelee_db->getScalar<double>( "rhoA" );
}
if (freelee_db->keyExists( "rhoB" )){
rhoB = freelee_db->getScalar<double>( "rhoB" );
}
if (freelee_db->keyExists( "F" )){
Fx = freelee_db->getVector<double>( "F" )[0];
Fy = freelee_db->getVector<double>( "F" )[1];
Fz = freelee_db->getVector<double>( "F" )[2];
}
if (freelee_db->keyExists( "gamma" )){
gamma = freelee_db->getScalar<double>( "gamma" );
}
if (freelee_db->keyExists( "W" )){
W = freelee_db->getScalar<double>( "W" );
}
if (freelee_db->keyExists( "Restart" )){
Restart = freelee_db->getScalar<bool>( "Restart" );
}
if (freelee_db->keyExists( "din" )){
din = freelee_db->getScalar<double>( "din" );
}
if (freelee_db->keyExists( "dout" )){
dout = freelee_db->getScalar<double>( "dout" );
}
if (freelee_db->keyExists( "flux" )){
flux = freelee_db->getScalar<double>( "flux" );
}
inletA=1.f;
inletB=0.f;
outletA=0.f;
outletB=1.f;
//if (BoundaryCondition==4) flux *= rhoA; // mass flux must adjust for density (see formulation for details)
BoundaryCondition = 0;
if (domain_db->keyExists( "BC" )){
BoundaryCondition = domain_db->getScalar<int>( "BC" );
}
}
void ScaLBL_FreeLeeModel::SetDomain(){
Dm = std::shared_ptr<Domain>(new Domain(domain_db,comm)); // full domain for analysis
Mask = std::shared_ptr<Domain>(new Domain(domain_db,comm)); // mask domain removes immobile phases
// domain parameters
Nx = Dm->Nx;
Ny = Dm->Ny;
Nz = Dm->Nz;
Lx = Dm->Lx;
Ly = Dm->Ly;
Lz = Dm->Lz;
N = Nx*Ny*Nz;
Nxh = Nx+2;
Nyh = Ny+2;
Nzh = Nz+2;
Nh = Nxh*Nyh*Nzh;
id = new signed char [N];
for (int i=0; i<Nx*Ny*Nz; i++) Dm->id[i] = 1; // initialize this way
2020-09-29 13:37:02 -04:00
2021-01-15 16:33:57 -05:00
comm.barrier();
2020-09-25 16:18:54 -04:00
Dm->CommInit();
2021-01-15 16:33:57 -05:00
comm.barrier();
2020-09-25 16:18:54 -04:00
// Read domain parameters
rank = Dm->rank();
nprocx = Dm->nprocx();
nprocy = Dm->nprocy();
nprocz = Dm->nprocz();
}
void ScaLBL_FreeLeeModel::ReadInput(){
sprintf(LocalRankString,"%05d",rank);
sprintf(LocalRankFilename,"%s%s","ID.",LocalRankString);
sprintf(LocalRestartFile,"%s%s","Restart.",LocalRankString);
if (freelee_db->keyExists( "image_sequence" )){
auto ImageList = freelee_db->getVector<std::string>( "image_sequence");
int IMAGE_INDEX = freelee_db->getWithDefault<int>( "image_index", 0 );
std::string first_image = ImageList[IMAGE_INDEX];
Mask->Decomp(first_image);
IMAGE_INDEX++;
}
else if (domain_db->keyExists( "GridFile" )){
// Read the local domain data
auto input_id = readMicroCT( *domain_db, MPI_COMM_WORLD );
// Fill the halo (assuming GCW of 1)
array<int,3> size0 = { (int) input_id.size(0), (int) input_id.size(1), (int) input_id.size(2) };
ArraySize size1 = { (size_t) Mask->Nx, (size_t) Mask->Ny, (size_t) Mask->Nz };
ASSERT( (int) size1[0] == size0[0]+2 && (int) size1[1] == size0[1]+2 && (int) size1[2] == size0[2]+2 );
fillHalo<signed char> fill( MPI_COMM_WORLD, Mask->rank_info, size0, { 1, 1, 1 }, 0, 1 );
Array<signed char> id_view;
2021-01-15 16:33:57 -05:00
id_view.viewRaw( size1, Mask->id.data() );
2020-09-25 16:18:54 -04:00
fill.copy( input_id, id_view );
fill.fill( id_view );
}
else if (domain_db->keyExists( "Filename" )){
auto Filename = domain_db->getScalar<std::string>( "Filename" );
Mask->Decomp(Filename);
}
else{
Mask->ReadIDs();
}
for (int i=0; i<Nx*Ny*Nz; i++) id[i] = Mask->id[i]; // save what was read
// Generate the signed distance map
// Initialize the domain and communication
Array<char> id_solid(Nx,Ny,Nz);
// Solve for the position of the solid phase
for (int k=0;k<Nz;k++){
for (int j=0;j<Ny;j++){
for (int i=0;i<Nx;i++){
int n = k*Nx*Ny+j*Nx+i;
// Initialize the solid phase
signed char label = Mask->id[n];
if (label > 0) id_solid(i,j,k) = 1;
else id_solid(i,j,k) = 0;
}
}
}
2020-09-29 13:37:02 -04:00
SignDist.resize(Nx,Ny,Nz);
2020-09-25 16:18:54 -04:00
// Initialize the signed distance function
for (int k=0;k<Nz;k++){
for (int j=0;j<Ny;j++){
for (int i=0;i<Nx;i++){
// Initialize distance to +/- 1
2020-09-29 13:37:02 -04:00
SignDist(i,j,k) = 2.0*double(id_solid(i,j,k))-1.0;
2020-09-25 16:18:54 -04:00
}
}
}
if (rank==0) printf("Initialized solid phase -- Converting to Signed Distance function \n");
2020-09-29 13:37:02 -04:00
CalcDist(SignDist,id_solid,*Mask);
2020-09-25 16:18:54 -04:00
if (rank == 0) cout << "Domain set." << endl;
2020-09-29 13:40:41 -04:00
2020-09-25 16:18:54 -04:00
}
void ScaLBL_FreeLeeModel::Create(){
/*
* This function creates the variables needed to run a LBM
*/
//.........................................................
// Initialize communication structures in averaging domain
for (int i=0; i<Nx*Ny*Nz; i++) Dm->id[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<ScaLBL_Communicator>(new ScaLBL_Communicator(Mask));
ScaLBL_Comm_Regular = std::shared_ptr<ScaLBL_Communicator>(new ScaLBL_Communicator(Mask));
2020-09-29 16:35:09 -04:00
ScaLBL_Comm_WideHalo = std::shared_ptr<ScaLBLWideHalo_Communicator>(new ScaLBLWideHalo_Communicator(Mask,2));
2020-09-25 16:18:54 -04:00
// create the layout for the LBM
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];
2021-01-13 22:06:08 -05:00
Np = ScaLBL_Comm->MemoryOptimizedLayoutAA(Map,neighborList,Mask->id.data(),Np,2);
2021-01-15 16:33:57 -05:00
comm.barrier();
2020-09-25 16:18:54 -04:00
//...........................................................................
// 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);
2021-01-18 21:30:27 -05:00
ScaLBL_AllocateDeviceMemory((void **) &gqbar, 19*dist_mem_size);
2020-09-25 16:18:54 -04:00
ScaLBL_AllocateDeviceMemory((void **) &hq, 7*dist_mem_size);
ScaLBL_AllocateDeviceMemory((void **) &mu_phi, dist_mem_size);
ScaLBL_AllocateDeviceMemory((void **) &Den, dist_mem_size);
ScaLBL_AllocateDeviceMemory((void **) &Phi, sizeof(double)*Nh);
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; k<Nz-1; k++){
for (int j=1; j<Ny-1; j++){
for (int i=1; i<Nx-1; i++){
int idx=Map(i,j,k);
if (!(idx < 0))
2021-01-18 21:30:27 -05:00
TmpMap[idx] = ScaLBL_Comm_WideHalo->Map(i,j,k);
2020-09-25 16:18:54 -04:00
}
}
}
// check that TmpMap is valid
for (int idx=0; idx<ScaLBL_Comm->LastExterior(); idx++){
auto n = TmpMap[idx];
2021-01-18 23:37:08 -05:00
if (n > Nxh*Nyh*Nzh){
2020-09-25 16:18:54 -04:00
printf("Bad value! idx=%i \n", n);
2021-01-18 23:37:08 -05:00
TmpMap[idx] = Nxh*Nyh*Nzh-1;
2020-09-25 16:18:54 -04:00
}
}
for (int idx=ScaLBL_Comm->FirstInterior(); idx<ScaLBL_Comm->LastInterior(); idx++){
auto n = TmpMap[idx];
2021-01-18 23:37:08 -05:00
if ( n > Nxh*Nyh*Nzh ){
2020-09-25 16:18:54 -04:00
printf("Bad value! idx=%i \n",n);
2021-01-18 23:37:08 -05:00
TmpMap[idx] = Nxh*Nyh*Nzh-1;
2020-09-25 16:18:54 -04:00
}
}
ScaLBL_CopyToDevice(dvcMap, TmpMap, sizeof(int)*Np);
ScaLBL_DeviceBarrier();
delete [] TmpMap;
// copy the neighbor list
ScaLBL_CopyToDevice(NeighborList, neighborList, neighborSize);
}
2021-01-18 21:30:27 -05:00
void ScaLBL_FreeLeeModel::AssignComponentLabels()
{
double *phase;
phase = new double[Nh];
2020-09-25 16:18:54 -04:00
2021-01-18 21:30:27 -05:00
size_t NLABELS=0;
signed char VALUE=0;
double AFFINITY=0.f;
auto LabelList = greyscaleColor_db->getVector<int>( "ComponentLabels" );
auto AffinityList = greyscaleColor_db->getVector<double>( "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 (size_t idx=0; idx<NLABELS; idx++) label_count[idx]=0;
for (int k=1;k<Nzh-1;k++){
for (int j=1;j<Nyh-1;j++){
for (int i=1;i<Nxh-1;i++){
int n = (k-1)*Nx*Ny+(j-1)*Nx+(i-1);
int nh = k*Nxh*Nyh+j*Nxh+i;
VALUE=id[n];
// Assign the affinity from the paired list
for (unsigned int idx=0; idx < NLABELS; idx++){
//printf("idx=%i, value=%i, %i, \n",idx, VALUE,LabelList[idx]);
if (VALUE == LabelList[idx]){
AFFINITY=AffinityList[idx];
label_count[idx] += 1.0;
idx = NLABELS;
//Mask->id[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; i<Nx*Ny*Nz; i++) Dm->id[i] = Mask->id[i];
2020-09-25 16:18:54 -04:00
2021-01-18 21:30:27 -05:00
for (size_t idx=0; idx<NLABELS; idx++)
label_count_global[idx] = sumReduce( Dm->Comm, label_count[idx]);
if (rank==0){
printf("Number of component labels: %lu \n",NLABELS);
for (unsigned int idx=0; idx<NLABELS; idx++){
VALUE=LabelList[idx];
AFFINITY=AffinityList[idx];
double volume_fraction = double(label_count_global[idx])/double((Nx-2)*(Ny-2)*(Nz-2)*nprocs);
printf(" label=%d, affinity=%f, volume fraction==%f\n",VALUE,AFFINITY,volume_fraction);
}
}
//compute color gradient and laplacian of phase field
//copy all data to device
ScaLBL_CopyToDevice(Phi, phase, N*sizeof(double));
ScaLBL_DeviceBarrier();
MPI_Barrier(ScaLBL_Comm->MPI_COMM_SCALBL);
delete [] phase;
}
void ScaLBL_FreeLeeModel::AssignChemPotential_ColorGrad()
{
double *SolidPotential_host = new double [Nx*Ny*Nz];
double *GreySolidGrad_host = new double [3*Np];
size_t NLABELS=0;
signed char VALUE=0;
double AFFINITY=0.f;
auto LabelList = greyscaleColor_db->getVector<int>( "GreySolidLabels" );
auto AffinityList = greyscaleColor_db->getVector<double>( "GreySolidAffinity" );
NLABELS=LabelList.size();
if (NLABELS != AffinityList.size()){
ERROR("Error: GreySolidLabels and GreySolidAffinity must be the same length! \n");
}
for (int k=0;k<Nz;k++){
for (int j=0;j<Ny;j++){
for (int i=0;i<Nx;i++){
int n = k*Nx*Ny+j*Nx+i;
VALUE=id[n];
AFFINITY=0.f;//all nodes except the specified grey nodes have grey-solid affinity = 0.0
// Assign the affinity from the paired list
for (unsigned int idx=0; idx < NLABELS; idx++){
//printf("idx=%i, value=%i, %i, \n",idx, VALUE,LabelList[idx]);
if (VALUE == LabelList[idx]){
AFFINITY=AffinityList[idx];
idx = NLABELS;
//Mask->id[n] = 0; // set mask to zero since this is an immobile component
}
}
SolidPotential_host[n] = AFFINITY;
}
}
}
// Calculate grey-solid color-gradient
double *Dst;
Dst = new double [3*3*3];
for (int kk=0; kk<3; kk++){
for (int jj=0; jj<3; jj++){
for (int ii=0; ii<3; ii++){
int index = kk*9+jj*3+ii;
Dst[index] = sqrt(double(ii-1)*double(ii-1) + double(jj-1)*double(jj-1)+ double(kk-1)*double(kk-1));
}
}
}
double w_face = 1.f;
double w_edge = 0.5;
double w_corner = 0.f;
//local
Dst[13] = 0.f;
//faces
Dst[4] = w_face;
Dst[10] = w_face;
Dst[12] = w_face;
Dst[14] = w_face;
Dst[16] = w_face;
Dst[22] = w_face;
// corners
Dst[0] = w_corner;
Dst[2] = w_corner;
Dst[6] = w_corner;
Dst[8] = w_corner;
Dst[18] = w_corner;
Dst[20] = w_corner;
Dst[24] = w_corner;
Dst[26] = w_corner;
// edges
Dst[1] = w_edge;
Dst[3] = w_edge;
Dst[5] = w_edge;
Dst[7] = w_edge;
Dst[9] = w_edge;
Dst[11] = w_edge;
Dst[15] = w_edge;
Dst[17] = w_edge;
Dst[19] = w_edge;
Dst[21] = w_edge;
Dst[23] = w_edge;
Dst[25] = w_edge;
for (int k=1; k<Nz-1; k++){
for (int j=1; j<Ny-1; j++){
for (int i=1; i<Nx-1; i++){
int idx=Map(i,j,k);
if (!(idx < 0)){
double phi_x = 0.f;
double phi_y = 0.f;
double phi_z = 0.f;
for (int kk=0; kk<3; kk++){
for (int jj=0; jj<3; jj++){
for (int ii=0; ii<3; ii++){
int index = kk*9+jj*3+ii;
double weight= Dst[index];
int idi=i+ii-1;
int idj=j+jj-1;
int idk=k+kk-1;
if (idi < 0) idi=0;
if (idj < 0) idj=0;
if (idk < 0) idk=0;
if (!(idi < Nx)) idi=Nx-1;
if (!(idj < Ny)) idj=Ny-1;
if (!(idk < Nz)) idk=Nz-1;
int nn = idk*Nx*Ny + idj*Nx + idi;
double vec_x = double(ii-1);
double vec_y = double(jj-1);
double vec_z = double(kk-1);
double GWNS=SolidPotential_host[nn];
phi_x += GWNS*weight*vec_x;
phi_y += GWNS*weight*vec_y;
phi_z += GWNS*weight*vec_z;
}
}
}
if (Averages->SDs(i,j,k)<2.0){
GreySolidGrad_host[idx+0*Np] = phi_x;
GreySolidGrad_host[idx+1*Np] = phi_y;
GreySolidGrad_host[idx+2*Np] = phi_z;
}
else{
GreySolidGrad_host[idx+0*Np] = 0.0;
GreySolidGrad_host[idx+1*Np] = 0.0;
GreySolidGrad_host[idx+2*Np] = 0.0;
}
}
}
}
}
if (rank==0){
printf("Number of Grey-solid labels: %lu \n",NLABELS);
for (unsigned int idx=0; idx<NLABELS; idx++){
VALUE=LabelList[idx];
AFFINITY=AffinityList[idx];
printf(" grey-solid label=%d, grey-solid affinity=%f\n",VALUE,AFFINITY);
}
}
ScaLBL_CopyToDevice(GreySolidGrad, GreySolidGrad_host, 3*Np*sizeof(double));
ScaLBL_DeviceBarrier();
delete [] SolidPotential_host;
delete [] GreySolidGrad_host;
delete [] Dst;
}
void ScaLBL_FreeLeeModel::Initialize(){
2020-09-25 16:18:54 -04:00
/*
* This function initializes model
*/
2021-01-18 21:30:27 -05:00
if (rank==0) printf ("Initializing phase field, chemical potential and color gradient\n");
AssignComponentLabels_ChemPotential_ColorGrad();//initialize phase field Phi
//if (rank==0) printf ("Initializing chemical potential and color gradient \n");
//AssignChemPotential_ColorGrad();
if (rank==0) printf ("Initializing distributions for momentum transport\n");
ScaLBL_D3Q19_FreeLeeModel_Init(gqbar, mu_phi, ColorGrad, Fx, Fy, Fz, Np);
if (rank==0) printf ("Initializing density field and distributions for phase-field transport\n");
ScaLBL_FreeLeeModel_PhaseField_Init(dvcMap, Phi, Den, hq, 0, ScaLBL_Comm->LastExterior(), Np);
ScaLBL_FreeLeeModel_PhaseField_Init(dvcMap, Phi, Den, hq, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
2020-09-25 16:18:54 -04:00
if (Restart == true){
2021-01-18 21:30:27 -05:00
//TODO need to revise this function
2020-09-25 16:18:54 -04:00
if (rank==0){
printf("Reading restart file! \n");
}
// 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));
2021-01-18 21:30:27 -05:00
//ScaLBL_CopyToHost(cPhi, Phi, N*sizeof(double));
2020-09-25 16:18:54 -04:00
ifstream File(LocalRestartFile,ios::binary);
int idx;
double value,va,vb;
for (int n=0; n<Np; n++){
File.read((char*) &va, sizeof(va));
File.read((char*) &vb, sizeof(vb));
cDen[n] = va;
cDen[Np+n] = vb;
}
for (int n=0; n<Np; n++){
// Read the distributions
for (int q=0; q<19; q++){
File.read((char*) &value, sizeof(value));
cDist[q*Np+n] = value;
}
}
File.close();
for (int n=0; n<ScaLBL_Comm->LastExterior(); n++){
va = cDen[n];
vb = cDen[Np + n];
value = (va-vb)/(va+vb);
idx = TmpMap[n];
if (!(idx < 0) && idx<N)
cPhi[idx] = value;
}
for (int n=ScaLBL_Comm->FirstInterior(); n<ScaLBL_Comm->LastInterior(); n++){
va = cDen[n];
vb = cDen[Np + n];
value = (va-vb)/(va+vb);
idx = TmpMap[n];
if (!(idx < 0) && idx<N)
cPhi[idx] = value;
}
// Copy the restart data to the GPU
ScaLBL_CopyToDevice(Den,cDen,2*Np*sizeof(double));
ScaLBL_CopyToDevice(fq,cDist,19*Np*sizeof(double));
ScaLBL_CopyToDevice(Phi,cPhi,N*sizeof(double));
2021-01-18 22:51:16 -05:00
ScaLBL_Comm->Barrier();
2021-01-15 16:33:57 -05:00
comm.barrier();
2020-09-25 16:18:54 -04:00
2021-01-18 21:30:27 -05:00
if (rank==0) printf ("Initializing phase and density fields on device from Restart\n");
ScaLBL_FreeLeeModel_PhaseField_InitFromRestart(Den, hq, 0, ScaLBL_Comm->LastExterior(), Np);
ScaLBL_FreeLeeModel_PhaseField_InitFromRestart(Den, hq, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
}
2020-09-25 16:18:54 -04:00
// establish reservoirs for external bC
if (BoundaryCondition == 1 || BoundaryCondition == 2 || BoundaryCondition == 3 || BoundaryCondition == 4 ){
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);
}
}
2020-09-29 13:51:23 -04:00
//ScaLBL_CopyToHost(Averages->Phi.data(),Phi,N*sizeof(double));
2020-09-25 16:18:54 -04:00
}
void ScaLBL_FreeLeeModel::Run(){
int nprocs=nprocx*nprocy*nprocz;
const RankInfoStruct rank_info(rank,nprocx,nprocy,nprocz);
if (rank==0){
printf("********************************************************\n");
printf("No. of timesteps: %i \n", timestepMax);
fflush(stdout);
}
//.......create and start timer............
double starttime,stoptime,cputime;
ScaLBL_DeviceBarrier();
2021-01-15 16:33:57 -05:00
comm.barrier();
2020-09-25 16:18:54 -04:00
starttime = MPI_Wtime();
//.........................................
//************ MAIN ITERATION LOOP ***************************************/
PROFILE_START("Loop");
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++;
2021-01-18 21:30:27 -05:00
//-------------------------------------------------------------------------------------------------------------------
// Compute the Phase indicator field
2020-09-25 16:18:54 -04:00
// Read for hq, Bq happens in this routine (requires communication)
2021-01-18 21:30:27 -05:00
//ScaLBL_Comm->SendD3Q7AA(hq); //READ FROM NORMAL
ScaLBL_Comm->SendD3Q7AA(hq); //READ FROM NORMAL
ScaLBL_D3Q7_AAodd_PhaseField(NeighborList, dvcMap, hq, Den, Phi, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
ScaLBL_Comm->RecvD3Q7AA(hq); //WRITE INTO OPPOSITE
2020-09-25 16:18:54 -04:00
ScaLBL_DeviceBarrier();
2021-01-18 21:30:27 -05:00
ScaLBL_D3Q7_AAodd_PhaseField(NeighborList, dvcMap, hq, Den, Phi, 0, ScaLBL_Comm->LastExterior(), Np);
2020-09-25 16:18:54 -04:00
// Perform the collision operation
2021-01-18 21:30:27 -05:00
ScaLBL_Comm->SendD3Q19AA(gqbar); //READ FROM NORMAL
2020-09-25 16:18:54 -04:00
if (BoundaryCondition > 0 && BoundaryCondition < 5){
2021-01-18 21:30:27 -05:00
//TODO to be revised
2020-09-25 16:18:54 -04:00
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
2021-01-18 21:30:27 -05:00
ScaLBL_Comm_WideHalo->Send(Phi);
2020-09-25 16:18:54 -04:00
2021-01-18 21:30:27 -05:00
ScaLBL_D3Q19_AAodd_FreeLeeModel(NeighborList, dvcMap, fq, hq, Bq, Den, Phi, Velocity, rhoA, rhoB, tauA, tauB,
2020-09-25 16:18:54 -04:00
alpha, beta, Fx, Fy, Fz, Nx, Nx*Ny, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
2021-01-18 21:30:27 -05:00
ScaLBL_Comm_WideHalo->Recv(Phi);
ScaLBL_Comm->RecvD3Q19AA(gqbar); //WRITE INTO OPPOSITE
2020-09-25 16:18:54 -04:00
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);
}
else if (BoundaryCondition == 5){
ScaLBL_Comm->D3Q19_Reflection_BC_z(fq);
ScaLBL_Comm->D3Q19_Reflection_BC_Z(fq);
}
2021-01-18 21:30:27 -05:00
ScaLBL_D3Q19_AAodd_FreeLeeModel(NeighborList, dvcMap, fq, hq, Bq, Den, Phi, Velocity, rhoA, rhoB, tauA, tauB,
2020-09-25 16:18:54 -04:00
alpha, beta, Fx, Fy, Fz, Nx, Nx*Ny, 0, ScaLBL_Comm->LastExterior(), Np);
ScaLBL_DeviceBarrier();
MPI_Barrier(ScaLBL_Comm->MPI_COMM_SCALBL);
// *************EVEN TIMESTEP*************
timestep++;
// Compute the Phase indicator field
2021-01-18 21:30:27 -05:00
ScaLBL_Comm->SendD3Q7AA(hq); //READ FROM NORMAL
ScaLBL_D3Q7_AAeven_PhaseField(dvcMap, hq, Den, Phi, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
ScaLBL_Comm->RecvD3Q7AA(hq); //WRITE INTO OPPOSITE
2020-09-25 16:18:54 -04:00
ScaLBL_DeviceBarrier();
2021-01-18 21:30:27 -05:00
ScaLBL_D3Q7_AAeven_PhaseField(dvcMap, hq, Den, Phi, 0, ScaLBL_Comm->LastExterior(), Np);
2020-09-25 16:18:54 -04:00
// Perform the collision operation
2021-01-18 21:30:27 -05:00
ScaLBL_Comm->SendD3Q19AA(gqbar); //READ FORM NORMAL
2020-09-25 16:18:54 -04:00
// Halo exchange for phase field
if (BoundaryCondition > 0 && BoundaryCondition < 5){
ScaLBL_Comm->Color_BC_z(dvcMap, Phi, Den, inletA, inletB);
ScaLBL_Comm->Color_BC_Z(dvcMap, Phi, Den, outletA, outletB);
}
2021-01-18 21:30:27 -05:00
ScaLBL_Comm_WideHalo->Send(Phi);
2020-09-25 16:18:54 -04:00
ScaLBL_D3Q19_AAeven_Color(dvcMap, fq, hq, Bq, Den, Phi, Velocity, rhoA, rhoB, tauA, tauB,
alpha, beta, Fx, Fy, Fz, Nx, Nx*Ny, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
2021-01-18 21:30:27 -05:00
ScaLBL_Comm_WideHalo->Recv(Phi);
ScaLBL_Comm->RecvD3Q19AA(gqbar); //WRITE INTO OPPOSITE
2020-09-25 16:18:54 -04:00
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);
}
else if (BoundaryCondition == 5){
ScaLBL_Comm->D3Q19_Reflection_BC_z(fq);
ScaLBL_Comm->D3Q19_Reflection_BC_Z(fq);
}
ScaLBL_D3Q19_AAeven_Color(dvcMap, fq, hq, Bq, Den, Phi, Velocity, rhoA, rhoB, tauA, tauB,
alpha, beta, Fx, Fy, Fz, Nx, Nx*Ny, 0, ScaLBL_Comm->LastExterior(), Np);
2021-01-15 16:33:57 -05:00
ScaLBL_Comm->Barrier();
2020-09-25 16:18:54 -04:00
//************************************************************************
PROFILE_STOP("Update");
}
PROFILE_STOP("Loop");
PROFILE_SAVE("lbpm_color_simulator",1);
//************************************************************************
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");
// ************************************************************************
}
void ScaLBL_FreeLeeModel::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);
ScaLBL_Comm->RegularLayout(Map,&Den[0],PhaseField);
FILE *AFILE;
sprintf(LocalRankFilename,"A.%05i.raw",rank);
AFILE = fopen(LocalRankFilename,"wb");
fwrite(PhaseField.data(),8,N,AFILE);
fclose(AFILE);
ScaLBL_Comm->RegularLayout(Map,&Den[Np],PhaseField);
FILE *BFILE;
sprintf(LocalRankFilename,"B.%05i.raw",rank);
BFILE = fopen(LocalRankFilename,"wb");
fwrite(PhaseField.data(),8,N,BFILE);
fclose(BFILE);
ScaLBL_Comm->RegularLayout(Map,Pressure,PhaseField);
FILE *PFILE;
sprintf(LocalRankFilename,"Pressure.%05i.raw",rank);
PFILE = fopen(LocalRankFilename,"wb");
fwrite(PhaseField.data(),8,N,PFILE);
fclose(PFILE);
ScaLBL_Comm->RegularLayout(Map,&Velocity[0],PhaseField);
FILE *VELX_FILE;
sprintf(LocalRankFilename,"Velocity_X.%05i.raw",rank);
VELX_FILE = fopen(LocalRankFilename,"wb");
fwrite(PhaseField.data(),8,N,VELX_FILE);
fclose(VELX_FILE);
ScaLBL_Comm->RegularLayout(Map,&Velocity[Np],PhaseField);
FILE *VELY_FILE;
sprintf(LocalRankFilename,"Velocity_Y.%05i.raw",rank);
VELY_FILE = fopen(LocalRankFilename,"wb");
fwrite(PhaseField.data(),8,N,VELY_FILE);
fclose(VELY_FILE);
ScaLBL_Comm->RegularLayout(Map,&Velocity[2*Np],PhaseField);
FILE *VELZ_FILE;
sprintf(LocalRankFilename,"Velocity_Z.%05i.raw",rank);
VELZ_FILE = fopen(LocalRankFilename,"wb");
fwrite(PhaseField.data(),8,N,VELZ_FILE);
fclose(VELZ_FILE);
/* ScaLBL_Comm->RegularLayout(Map,&ColorGrad[0],PhaseField);
FILE *CGX_FILE;
sprintf(LocalRankFilename,"Gradient_X.%05i.raw",rank);
CGX_FILE = fopen(LocalRankFilename,"wb");
fwrite(PhaseField.data(),8,N,CGX_FILE);
fclose(CGX_FILE);
ScaLBL_Comm->RegularLayout(Map,&ColorGrad[Np],PhaseField);
FILE *CGY_FILE;
sprintf(LocalRankFilename,"Gradient_Y.%05i.raw",rank);
CGY_FILE = fopen(LocalRankFilename,"wb");
fwrite(PhaseField.data(),8,N,CGY_FILE);
fclose(CGY_FILE);
ScaLBL_Comm->RegularLayout(Map,&ColorGrad[2*Np],PhaseField);
FILE *CGZ_FILE;
sprintf(LocalRankFilename,"Gradient_Z.%05i.raw",rank);
CGZ_FILE = fopen(LocalRankFilename,"wb");
fwrite(PhaseField.data(),8,N,CGZ_FILE);
fclose(CGZ_FILE);
*/
}