done: add restart and write visfiles for greyscale Color
This commit is contained in:
parent
f3ba90337a
commit
d02eff3017
@ -186,6 +186,8 @@ extern "C" void ScaLBL_D3Q19_AAodd_GreyscaleColor(int *d_neighborList, int *Map,
|
|||||||
double rhoA, double rhoB, double tauA, double tauB, double tauA_eff,double tauB_eff, double alpha, double beta,
|
double rhoA, double rhoB, double tauA, double tauB, double tauA_eff,double tauB_eff, double alpha, double beta,
|
||||||
double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np);
|
double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np);
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_PhaseField_InitFromRestart(double *Den, double *Aq, double *Bq, int start, int finish, int Np);
|
||||||
|
|
||||||
// MRT MODEL
|
// MRT MODEL
|
||||||
extern "C" void ScaLBL_D3Q19_AAeven_MRT(double *dist, int start, int finish, int Np, double rlx_setA, double rlx_setB, double Fx,
|
extern "C" void ScaLBL_D3Q19_AAeven_MRT(double *dist, int start, int finish, int Np, double rlx_setA, double rlx_setB, double Fx,
|
||||||
double Fy, double Fz);
|
double Fy, double Fz);
|
||||||
|
@ -1338,6 +1338,32 @@ extern "C" void ScaLBL_D3Q19_AAeven_GreyscaleColor(int *Map, double *dist, doubl
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_PhaseField_InitFromRestart(double *Den, double *Aq, double *Bq, int start, int finish, int Np){
|
||||||
|
int idx;
|
||||||
|
double nA,nB;
|
||||||
|
|
||||||
|
for (idx=start; idx<finish; idx++){
|
||||||
|
|
||||||
|
nA = Den[idx];
|
||||||
|
nB = Den[Np+idx];
|
||||||
|
|
||||||
|
Aq[idx]=0.3333333333333333*nA;
|
||||||
|
Aq[Np+idx]=0.1111111111111111*nA;
|
||||||
|
Aq[2*Np+idx]=0.1111111111111111*nA;
|
||||||
|
Aq[3*Np+idx]=0.1111111111111111*nA;
|
||||||
|
Aq[4*Np+idx]=0.1111111111111111*nA;
|
||||||
|
Aq[5*Np+idx]=0.1111111111111111*nA;
|
||||||
|
Aq[6*Np+idx]=0.1111111111111111*nA;
|
||||||
|
|
||||||
|
Bq[idx]=0.3333333333333333*nB;
|
||||||
|
Bq[Np+idx]=0.1111111111111111*nB;
|
||||||
|
Bq[2*Np+idx]=0.1111111111111111*nB;
|
||||||
|
Bq[3*Np+idx]=0.1111111111111111*nB;
|
||||||
|
Bq[4*Np+idx]=0.1111111111111111*nB;
|
||||||
|
Bq[5*Np+idx]=0.1111111111111111*nB;
|
||||||
|
Bq[6*Np+idx]=0.1111111111111111*nB;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
//extern "C" void ScaLBL_D3Q19_GreyscaleColor_Init(double *dist, double *Porosity, int Np){
|
//extern "C" void ScaLBL_D3Q19_GreyscaleColor_Init(double *dist, double *Porosity, int Np){
|
||||||
// int n;
|
// int n;
|
||||||
|
@ -1447,6 +1447,37 @@ __global__ void dvc_ScaLBL_D3Q19_AAeven_GreyscaleColor(int *Map, double *dist,
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
__global__ void dvc_ScaLBL_PhaseField_InitFromRestart(double *Den, double *Aq, double *Bq, int start, int finish, int Np){
|
||||||
|
int idx;
|
||||||
|
double nA,nB;
|
||||||
|
|
||||||
|
int S = Np/NBLOCKS/NTHREADS + 1;
|
||||||
|
for (int s=0; s<S; s++){
|
||||||
|
//........Get 1-D index for this thread....................
|
||||||
|
idx = S*blockIdx.x*blockDim.x + s*blockDim.x + threadIdx.x + start;
|
||||||
|
if (idx<finish) {
|
||||||
|
|
||||||
|
nA = Den[idx];
|
||||||
|
nB = Den[Np+idx];
|
||||||
|
|
||||||
|
Aq[idx]=0.3333333333333333*nA;
|
||||||
|
Aq[Np+idx]=0.1111111111111111*nA;
|
||||||
|
Aq[2*Np+idx]=0.1111111111111111*nA;
|
||||||
|
Aq[3*Np+idx]=0.1111111111111111*nA;
|
||||||
|
Aq[4*Np+idx]=0.1111111111111111*nA;
|
||||||
|
Aq[5*Np+idx]=0.1111111111111111*nA;
|
||||||
|
Aq[6*Np+idx]=0.1111111111111111*nA;
|
||||||
|
|
||||||
|
Bq[idx]=0.3333333333333333*nB;
|
||||||
|
Bq[Np+idx]=0.1111111111111111*nB;
|
||||||
|
Bq[2*Np+idx]=0.1111111111111111*nB;
|
||||||
|
Bq[3*Np+idx]=0.1111111111111111*nB;
|
||||||
|
Bq[4*Np+idx]=0.1111111111111111*nB;
|
||||||
|
Bq[5*Np+idx]=0.1111111111111111*nB;
|
||||||
|
Bq[6*Np+idx]=0.1111111111111111*nB;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
////Model-2&3
|
////Model-2&3
|
||||||
//__global__ void dvc_ScaLBL_D3Q19_AAodd_GreyscaleColor(int *neighborList, int *Map, double *dist, double *Aq, double *Bq, double *Den,
|
//__global__ void dvc_ScaLBL_D3Q19_AAodd_GreyscaleColor(int *neighborList, int *Map, double *dist, double *Aq, double *Bq, double *Den,
|
||||||
// double *Phi, double *GreySolidGrad, double *Poros,double *Perm, double *Velocity,
|
// double *Phi, double *GreySolidGrad, double *Poros,double *Perm, double *Velocity,
|
||||||
@ -2959,6 +2990,13 @@ extern "C" void ScaLBL_D3Q19_AAodd_GreyscaleColor(int *d_neighborList, int *Map,
|
|||||||
//cudaProfilerStop();
|
//cudaProfilerStop();
|
||||||
}
|
}
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_PhaseField_InitFromRestart(double *Den, double *Aq, double *Bq, int start, int finish, int Np){
|
||||||
|
dvc_ScaLBL_PhaseField_InitFromRestart<<<NBLOCKS,NTHREADS >>>(Den, Aq, Bq, start, finish, Np);
|
||||||
|
cudaError_t err = cudaGetLastError();
|
||||||
|
if (cudaSuccess != err){
|
||||||
|
printf("CUDA error in ScaLBL_PhaseField_InitFromRestart: %s \n",cudaGetErrorString(err));
|
||||||
|
}
|
||||||
|
}
|
||||||
////Model-2&3
|
////Model-2&3
|
||||||
//extern "C" void ScaLBL_D3Q19_AAeven_GreyscaleColor(int *Map, double *dist, double *Aq, double *Bq, double *Den,
|
//extern "C" void ScaLBL_D3Q19_AAeven_GreyscaleColor(int *Map, double *dist, double *Aq, double *Bq, double *Den,
|
||||||
// double *Phi,double *GreySolidGrad, double *Poros,double *Perm,double *Vel,
|
// double *Phi,double *GreySolidGrad, double *Poros,double *Perm,double *Vel,
|
||||||
|
@ -673,13 +673,17 @@ void ScaLBL_GreyscaleColorModel::Create(){
|
|||||||
}
|
}
|
||||||
|
|
||||||
void ScaLBL_GreyscaleColorModel::Initialize(){
|
void ScaLBL_GreyscaleColorModel::Initialize(){
|
||||||
|
|
||||||
if (rank==0) printf ("Initializing distributions \n");
|
|
||||||
ScaLBL_D3Q19_Init(fq, Np);
|
|
||||||
//ScaLBL_D3Q19_GreyscaleColor_Init(fq, Porosity_dvc, Np);
|
|
||||||
/*
|
/*
|
||||||
* This function initializes model
|
* This function initializes model
|
||||||
*/
|
*/
|
||||||
|
if (rank==0) printf ("Initializing distributions \n");
|
||||||
|
ScaLBL_D3Q19_Init(fq, Np);
|
||||||
|
//ScaLBL_D3Q19_GreyscaleColor_Init(fq, Porosity_dvc, Np);
|
||||||
|
|
||||||
|
if (rank==0) printf ("Initializing phase field \n");
|
||||||
|
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 (Restart == true){
|
if (Restart == true){
|
||||||
if (rank==0){
|
if (rank==0){
|
||||||
printf("Reading restart file! \n");
|
printf("Reading restart file! \n");
|
||||||
@ -738,11 +742,11 @@ void ScaLBL_GreyscaleColorModel::Initialize(){
|
|||||||
ScaLBL_DeviceBarrier();
|
ScaLBL_DeviceBarrier();
|
||||||
|
|
||||||
MPI_Barrier(comm);
|
MPI_Barrier(comm);
|
||||||
}
|
|
||||||
|
|
||||||
if (rank==0) printf ("Initializing phase field \n");
|
if (rank==0) printf ("Initializing phase field from Restart\n");
|
||||||
ScaLBL_PhaseField_Init(dvcMap, Phi, Den, Aq, Bq, 0, ScaLBL_Comm->LastExterior(), Np);
|
ScaLBL_PhaseField_InitFromRestart(Den, Aq, Bq, 0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
ScaLBL_PhaseField_Init(dvcMap, Phi, Den, Aq, Bq, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
ScaLBL_PhaseField_InitFromRestart(Den, Aq, Bq, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
||||||
|
}
|
||||||
|
|
||||||
// establish reservoirs for external bC
|
// establish reservoirs for external bC
|
||||||
if (BoundaryCondition == 1 || BoundaryCondition == 2 || BoundaryCondition == 3 || BoundaryCondition == 4 ){
|
if (BoundaryCondition == 1 || BoundaryCondition == 2 || BoundaryCondition == 3 || BoundaryCondition == 4 ){
|
||||||
@ -794,6 +798,8 @@ void ScaLBL_GreyscaleColorModel::Run(){
|
|||||||
double delta_volume_target = 0.0;
|
double delta_volume_target = 0.0;
|
||||||
|
|
||||||
//TODO -------- For temporary use - should be included in the analysis framework later -------------
|
//TODO -------- For temporary use - should be included in the analysis framework later -------------
|
||||||
|
int visualization_interval = 50000;
|
||||||
|
int restart_interval = 100000;
|
||||||
if (analysis_db->keyExists( "visualization_interval" )){
|
if (analysis_db->keyExists( "visualization_interval" )){
|
||||||
visualization_interval = analysis_db->getScalar<int>( "visualization_interval" );
|
visualization_interval = analysis_db->getScalar<int>( "visualization_interval" );
|
||||||
}
|
}
|
||||||
@ -1482,11 +1488,11 @@ void ScaLBL_GreyscaleColorModel::WriteVisFiles(){
|
|||||||
PhaseVar->type = IO::VariableType::VolumeVariable;
|
PhaseVar->type = IO::VariableType::VolumeVariable;
|
||||||
PhaseVar->dim = 1;
|
PhaseVar->dim = 1;
|
||||||
PhaseVar->data.resize(Dm->Nx-2,Dm->Ny-2,Dm->Nz-2);
|
PhaseVar->data.resize(Dm->Nx-2,Dm->Ny-2,Dm->Nz-2);
|
||||||
visData[0].vars.push_back(PaseVar);
|
visData[0].vars.push_back(PhaseVar);
|
||||||
|
|
||||||
ASSERT(visData[0].vars[0]->name=="Phase");
|
ASSERT(visData[0].vars[0]->name=="Phase");
|
||||||
Array<double>& PhaseData = visData[0].vars[0]->data;
|
Array<double>& PhaseData = visData[0].vars[0]->data;
|
||||||
ScaLBL_Comm->RegularLayout(Map,Phase,DataTemp);
|
ScaLBL_CopyToHost(DataTemp.data(), Phi, sizeof(double)*Nx*Ny*Nz);
|
||||||
fillData.copy(DataTemp,PhaseData);
|
fillData.copy(DataTemp,PhaseData);
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -1546,7 +1552,7 @@ void ScaLBL_GreyscaleColorModel::WriteVisFiles(){
|
|||||||
|
|
||||||
ASSERT(visData[0].vars[5]->name=="SignDist");
|
ASSERT(visData[0].vars[5]->name=="SignDist");
|
||||||
Array<double>& SignData = visData[0].vars[5]->data;
|
Array<double>& SignData = visData[0].vars[5]->data;
|
||||||
fillData.copy(Averages.SDs,SignData);
|
fillData.copy(Averages->SDs,SignData);
|
||||||
}
|
}
|
||||||
|
|
||||||
if (vis_db->getWithDefault<bool>( "write_silo", true )){
|
if (vis_db->getWithDefault<bool>( "write_silo", true )){
|
||||||
|
@ -90,6 +90,6 @@ private:
|
|||||||
double MorphInit(const double beta, const double morph_delta);
|
double MorphInit(const double beta, const double morph_delta);
|
||||||
double SeedPhaseField(const double seed_water_in_oil);
|
double SeedPhaseField(const double seed_water_in_oil);
|
||||||
double MorphOpenConnected(double target_volume_change);
|
double MorphOpenConnected(double target_volume_change);
|
||||||
double WriteVisFiles();
|
void WriteVisFiles();
|
||||||
};
|
};
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user