save the work;CPU version compiled; to be tested
This commit is contained in:
@@ -857,8 +857,8 @@ int ScaLBL_Communicator::MemoryOptimizedLayoutAA(IntArray &Map, int *neighborLis
|
|||||||
void ScaLBL_Communicator::SetupBounceBackList(IntArray &Map, signed char *id, int Np)
|
void ScaLBL_Communicator::SetupBounceBackList(IntArray &Map, signed char *id, int Np)
|
||||||
{
|
{
|
||||||
|
|
||||||
int idx,i,j,k;
|
int idx,i,j,k;
|
||||||
int neighbor;
|
int neighbor;
|
||||||
// save list of bounce-back distributions and interaction sites
|
// save list of bounce-back distributions and interaction sites
|
||||||
n_bb_d3q7 = 0; n_bb_d3q19 = 0;
|
n_bb_d3q7 = 0; n_bb_d3q19 = 0;
|
||||||
|
|
||||||
@@ -929,8 +929,9 @@ void ScaLBL_Communicator::SetupBounceBackList(IntArray &Map, signed char *id, in
|
|||||||
}
|
}
|
||||||
|
|
||||||
int *bb_dist_tmp = new int [local_count];
|
int *bb_dist_tmp = new int [local_count];
|
||||||
bb_interactions = new int [local_count];
|
int *bb_interactions_tmp = new int [local_count];
|
||||||
ScaLBL_AllocateDeviceMemory((void **) &bb_dist, sizeof(int)*local_count);
|
ScaLBL_AllocateDeviceMemory((void **) &bb_dist, sizeof(int)*local_count);
|
||||||
|
ScaLBL_AllocateDeviceMemory((void **) &bb_interactions, sizeof(int)*local_count);
|
||||||
|
|
||||||
local_count=0;
|
local_count=0;
|
||||||
for (k=1;k<Nz-1;k++){
|
for (k=1;k<Nz-1;k++){
|
||||||
@@ -943,37 +944,37 @@ void ScaLBL_Communicator::SetupBounceBackList(IntArray &Map, signed char *id, in
|
|||||||
int neighbor; // cycle through the neighbors of lattice site idx
|
int neighbor; // cycle through the neighbors of lattice site idx
|
||||||
neighbor=Map(i-1,j,k);
|
neighbor=Map(i-1,j,k);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i-1) + (j)*Nx + (k)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i-1) + (j)*Nx + (k)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 2*Np;
|
bb_dist_tmp[local_count++]=idx + 2*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i+1,j,k);
|
neighbor=Map(i+1,j,k);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i+1) + (j)*Nx + (k)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i+1) + (j)*Nx + (k)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++] = idx + 1*Np;
|
bb_dist_tmp[local_count++] = idx + 1*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i,j-1,k);
|
neighbor=Map(i,j-1,k);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i) + (j-1)*Nx + (k)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i) + (j-1)*Nx + (k)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 4*Np;
|
bb_dist_tmp[local_count++]=idx + 4*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i,j+1,k);
|
neighbor=Map(i,j+1,k);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i) + (j+1)*Nx + (k)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i) + (j+1)*Nx + (k)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 3*Np;
|
bb_dist_tmp[local_count++]=idx + 3*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i,j,k-1);
|
neighbor=Map(i,j,k-1);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i) + (j)*Nx + (k-1)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i) + (j)*Nx + (k-1)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 6*Np;
|
bb_dist_tmp[local_count++]=idx + 6*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i,j,k+1);
|
neighbor=Map(i,j,k+1);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i) + (j)*Nx + (k+1)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i) + (j)*Nx + (k+1)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 5*Np;
|
bb_dist_tmp[local_count++]=idx + 5*Np;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@@ -990,73 +991,73 @@ void ScaLBL_Communicator::SetupBounceBackList(IntArray &Map, signed char *id, in
|
|||||||
|
|
||||||
neighbor=Map(i-1,j-1,k);
|
neighbor=Map(i-1,j-1,k);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i-1) + (j-1)*Nx + (k)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i-1) + (j-1)*Nx + (k)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 8*Np;
|
bb_dist_tmp[local_count++]=idx + 8*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i+1,j+1,k);
|
neighbor=Map(i+1,j+1,k);
|
||||||
if (neighbor==-1) {
|
if (neighbor==-1) {
|
||||||
bb_interactions[local_count] = (i+1) + (j+1)*Nx + (k)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i+1) + (j+1)*Nx + (k)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 7*Np;
|
bb_dist_tmp[local_count++]=idx + 7*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i-1,j+1,k);
|
neighbor=Map(i-1,j+1,k);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i-1) + (j+1)*Nx + (k)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i-1) + (j+1)*Nx + (k)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 10*Np;
|
bb_dist_tmp[local_count++]=idx + 10*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i+1,j-1,k);
|
neighbor=Map(i+1,j-1,k);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i+1) + (j-1)*Nx + (k)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i+1) + (j-1)*Nx + (k)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 9*Np;
|
bb_dist_tmp[local_count++]=idx + 9*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i-1,j,k-1);
|
neighbor=Map(i-1,j,k-1);
|
||||||
if (neighbor==-1) {
|
if (neighbor==-1) {
|
||||||
bb_interactions[local_count] = (i-1) + (j)*Nx + (k-1)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i-1) + (j)*Nx + (k-1)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 12*Np;
|
bb_dist_tmp[local_count++]=idx + 12*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i+1,j,k+1);
|
neighbor=Map(i+1,j,k+1);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i+1) + (j)*Nx + (k+1)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i+1) + (j)*Nx + (k+1)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 11*Np;
|
bb_dist_tmp[local_count++]=idx + 11*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i-1,j,k+1);
|
neighbor=Map(i-1,j,k+1);
|
||||||
if (neighbor==-1) {
|
if (neighbor==-1) {
|
||||||
bb_interactions[local_count] = (i-1) + (j)*Nx + (k+1)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i-1) + (j)*Nx + (k+1)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 14*Np;
|
bb_dist_tmp[local_count++]=idx + 14*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i+1,j,k-1);
|
neighbor=Map(i+1,j,k-1);
|
||||||
if (neighbor==-1) {
|
if (neighbor==-1) {
|
||||||
bb_interactions[local_count] = (i+1) + (j)*Nx + (k-1)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i+1) + (j)*Nx + (k-1)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 13*Np;
|
bb_dist_tmp[local_count++]=idx + 13*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i,j-1,k-1);
|
neighbor=Map(i,j-1,k-1);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i) + (j-1)*Nx + (k-1)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i) + (j-1)*Nx + (k-1)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 16*Np;
|
bb_dist_tmp[local_count++]=idx + 16*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i,j+1,k+1);
|
neighbor=Map(i,j+1,k+1);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i) + (j+1)*Nx + (k+1)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i) + (j+1)*Nx + (k+1)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 15*Np;
|
bb_dist_tmp[local_count++]=idx + 15*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i,j-1,k+1);
|
neighbor=Map(i,j-1,k+1);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i) + (j-1)*Nx + (k+1)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i) + (j-1)*Nx + (k+1)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 18*Np;
|
bb_dist_tmp[local_count++]=idx + 18*Np;
|
||||||
}
|
}
|
||||||
|
|
||||||
neighbor=Map(i,j+1,k-1);
|
neighbor=Map(i,j+1,k-1);
|
||||||
if (neighbor==-1){
|
if (neighbor==-1){
|
||||||
bb_interactions[local_count] = (i) + (j+1)*Nx + (k-1)*Nx*Ny;
|
bb_interactions_tmp[local_count] = (i) + (j+1)*Nx + (k-1)*Nx*Ny;
|
||||||
bb_dist_tmp[local_count++]=idx + 17*Np;
|
bb_dist_tmp[local_count++]=idx + 17*Np;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@@ -1065,18 +1066,26 @@ void ScaLBL_Communicator::SetupBounceBackList(IntArray &Map, signed char *id, in
|
|||||||
}
|
}
|
||||||
n_bb_d3q19 = local_count; // this gives the d3q19 distributions not part of d3q7 model
|
n_bb_d3q19 = local_count; // this gives the d3q19 distributions not part of d3q7 model
|
||||||
ScaLBL_CopyToDevice(bb_dist, bb_dist_tmp, local_count*sizeof(int));
|
ScaLBL_CopyToDevice(bb_dist, bb_dist_tmp, local_count*sizeof(int));
|
||||||
|
ScaLBL_CopyToDevice(bb_interactions, bb_interactions_tmp, local_count*sizeof(int));
|
||||||
|
ScaLBL_DeviceBarrier();
|
||||||
|
|
||||||
|
delete [] bb_dist_tmp;
|
||||||
|
delete [] bb_interactions_tmp;
|
||||||
}
|
}
|
||||||
|
|
||||||
void ScaLBL_Communicator::SolidDirichletD3Q7(double *fq, double *assignValues){
|
void ScaLBL_Communicator::SolidDirichletD3Q7(double *fq, double *BoundaryValue){
|
||||||
// fq is a D3Q7 distribution
|
// fq is a D3Q7 distribution
|
||||||
// assignValues is a list of values to assign at bounce-back sites
|
// BoundaryValues is a list of values to assign at bounce-back solid sites
|
||||||
for (int idx=0; idx<n_bb_d3q7; idx++){
|
ScaLBL_Solid_Dirichlet_D3Q7(fq,BoundaryValue,bb_dist,bb_interactions,n_bb_d3q7);
|
||||||
double value = assignValues[idx];
|
|
||||||
int iq = bb_dist[idx];
|
|
||||||
fq[iq] += value;
|
|
||||||
}
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
void ScaLBL_Communicator::SolidNeumannD3Q7(double *fq, double *BoundaryValue){
|
||||||
|
// fq is a D3Q7 distribution
|
||||||
|
// BoundaryValues is a list of values to assign at bounce-back solid sites
|
||||||
|
ScaLBL_Solid_Neumann_D3Q7(fq,BoundaryValue,bb_dist,bb_interactions,n_bb_d3q7);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
void ScaLBL_Communicator::SendD3Q19AA(double *dist){
|
void ScaLBL_Communicator::SendD3Q19AA(double *dist){
|
||||||
|
|
||||||
// NOTE: the center distribution f0 must NOT be at the start of feven, provide offset to start of f2
|
// NOTE: the center distribution f0 must NOT be at the start of feven, provide offset to start of f2
|
||||||
|
|||||||
@@ -73,16 +73,39 @@ extern "C" void ScaLBL_D3Q19_AAodd_Greyscale_IMRT(int *neighborList, double *dis
|
|||||||
double *Poros,double *Perm, double *Velocity,double Den,double *Pressure);
|
double *Poros,double *Perm, double *Velocity,double Den,double *Pressure);
|
||||||
|
|
||||||
// ION TRANSPORT MODEL
|
// ION TRANSPORT MODEL
|
||||||
extern "C" void ScaLBL_D3Q7_AAeven_Ion(double *dist, double *Velocity, int start, int finish, int Np, double rlx, double Fx, double Fy, double Fz);
|
|
||||||
|
|
||||||
extern "C" void ScaLBL_D3Q7_AAodd_Ion(int *neighborList, double *dist, double *Velocity, int start, int finish, int Np, double rlx, double Fx, double Fy, double Fz);
|
extern "C" void ScaLBL_D3Q7_AAodd_IonConcentration(int *neighborList, double *dist, double *Den, int start, int finish, int Np);
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_D3Q7_AAeven_IonConcentration(double *dist, double *Den, int start, int finish, int Np);
|
||||||
|
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_D3Q7_AAodd_Ion(int *neighborList, double *dist, double *Den, double *Velocity, double *ElectricField,
|
||||||
|
double Di, double zi, double rlx, double Vt, int start, int finish, int Np);
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_D3Q7_AAeven_Ion(double *dist, double *Den, double *Velocity, double *ElectricField,
|
||||||
|
double Di, double zi, double rlx, double Vt, int start, int finish, int Np);
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_D3Q7_Ion_Init(double *dist, double *Den, double DenInit, int Np);
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_D3Q7_Ion_ChargeDensity(double *Den, double *ChargeDensity, int IonValence, int ion_component, int start, int finish, int Np);
|
||||||
|
|
||||||
// LBM Poisson solver
|
// LBM Poisson solver
|
||||||
extern "C" void ScaLBL_D3Q7_AAeven_Poisson(double *dist, double *ChargeDensity, int start, int finish, int Np, double rlx, double Fx, double Fy, double Fz);
|
|
||||||
|
|
||||||
extern "C" void ScaLBL_D3Q7_AAodd_Poisson(int *neighborList, double *dist, double *ChargeDensity, int start, int finish, int Np, double rlx, double Fx, double Fy, double Fz);
|
extern "C" void ScaLBL_D3Q7_AAodd_Poisson(int *neighborList, double *dist, double *Den_charge, double *Psi, double *ElectricField, double tau, double epsilon_LB,double gamma,
|
||||||
|
int start, int finish, int Np);
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_D3Q7_AAeven_Poisson(double *dist, double *Den_charge, double *Psi, double *ElectricField, double tau, double epsilon_LB,double gamma,
|
||||||
|
int start, int finish, int Np);
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_D3Q7_Poisson_Init(double *dist, int Np);
|
||||||
|
|
||||||
|
// LBM Stokes Model (adapted from MRT model)
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_D3Q19_AAeven_StokesMRT(double *dist, double *Velocity, double *ChargeDensity, double *ElectricField, double rlx_setA, double rlx_setB,
|
||||||
|
double Gx, double Gy, double Gz, int start, int finish, int Np);
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_D3Q19_AAodd_StokesMRT(int *neighborList, double *dist, double *Velocity, double *ChargeDensity, double *ElectricField, double rlx_setA, double rlx_setB,
|
||||||
|
double Gx, double Gy, double Gz, 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,
|
||||||
@@ -159,6 +182,10 @@ extern "C" void ScaLBL_SetSlice_z(double *Phi, double value, int Nx, int Ny, int
|
|||||||
|
|
||||||
extern "C" void ScaLBL_CopySlice_z(double *Phi, int Nx, int Ny, int Nz, int Source, int Destination);
|
extern "C" void ScaLBL_CopySlice_z(double *Phi, int Nx, int Ny, int Nz, int Source, int Destination);
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_Solid_Dirichlet_D3Q7(double *dist,double *BoundaryValue,int *BounceBackDist_list,int *BounceBackSolid_list,int N);
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_Solid_Neumann_D3Q7(double *dist,double *BoundaryValue,int *BounceBackDist_list,int *BounceBackSolid_list,int N);
|
||||||
|
|
||||||
class ScaLBL_Communicator{
|
class ScaLBL_Communicator{
|
||||||
public:
|
public:
|
||||||
//......................................................................................
|
//......................................................................................
|
||||||
@@ -207,7 +234,8 @@ public:
|
|||||||
void RecvGrad(double *Phi, double *Gradient);
|
void RecvGrad(double *Phi, double *Gradient);
|
||||||
void RegularLayout(IntArray map, const double *data, DoubleArray ®data);
|
void RegularLayout(IntArray map, const double *data, DoubleArray ®data);
|
||||||
void SetupBounceBackList(IntArray &Map, signed char *id, int Np);
|
void SetupBounceBackList(IntArray &Map, signed char *id, int Np);
|
||||||
void SolidDirichletD3Q7(double *fq, double *assignValues);
|
void SolidDirichletD3Q7(double *fq, double *BoundaryValue);
|
||||||
|
void SolidNeumannD3Q7(double *fq, double *BoundaryValue);
|
||||||
|
|
||||||
// Routines to set boundary conditions
|
// Routines to set boundary conditions
|
||||||
void Color_BC_z(int *Map, double *Phi, double *Den, double vA, double vB);
|
void Color_BC_z(int *Map, double *Phi, double *Den, double vA, double vB);
|
||||||
|
|||||||
32
cpu/D3Q7BC.cpp
Normal file
32
cpu/D3Q7BC.cpp
Normal file
@@ -0,0 +1,32 @@
|
|||||||
|
// CPU Functions for D3Q7 Lattice Boltzmann Methods
|
||||||
|
// Boundary Conditions
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_Solid_Dirichlet_D3Q7(double *dist,double *BoundaryValue,int *BounceBackDist_list,int *BounceBackSolid_list,int N){
|
||||||
|
|
||||||
|
int idx;
|
||||||
|
int iq,ib;
|
||||||
|
double value_b,value_q;
|
||||||
|
for (idx=0; idx<N; idx++){
|
||||||
|
iq = BounceBackDist_list[idx];
|
||||||
|
ib = BounceBackSolid_list[idx];
|
||||||
|
value_b = BoundaryValue[ib];//get boundary value from a solid site
|
||||||
|
value_q = dist[iq];
|
||||||
|
dist[iq] = -1.0*value_q + value_b*2.0/9.0;//NOTE 2/9 is the speed of sound for D3Q7 lattice
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
extern "C" void ScaLBL_Solid_Neumann_D3Q7(double *dist,double *BoundaryValue,int *BounceBackDist_list,int *BounceBackSolid_list,int N){
|
||||||
|
|
||||||
|
int idx;
|
||||||
|
int iq,ib;
|
||||||
|
double value_b,value_q;
|
||||||
|
for (idx=0; idx<N; idx++){
|
||||||
|
iq = BounceBackDist_list[idx];
|
||||||
|
ib = BounceBackSolid_list[idx];
|
||||||
|
value_b = BoundaryValue[ib];//get boundary value from a solid site
|
||||||
|
value_q = dist[iq];
|
||||||
|
dist[iq] = value_q + value_b;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
23
cpu/Ion.cpp
23
cpu/Ion.cpp
@@ -219,26 +219,19 @@ extern "C" void ScaLBL_D3Q7_Ion_Init(double *dist, double *Den, double DenInit,
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
extern "C" void ScaLBL_D3Q7_Ion_ChargeDensity(double *Den, double *ChargeDensity, vector<double>& IonValence, int number_ion_species, int start, int finish, int Np){
|
extern "C" void ScaLBL_D3Q7_Ion_ChargeDensity(double *Den, double *ChargeDensity, int IonValence, int ion_component, int start, int finish, int Np){
|
||||||
|
|
||||||
int n;
|
int n;
|
||||||
int ic=number_ion_species;
|
|
||||||
double Ci;//ion concentration of species i
|
double Ci;//ion concentration of species i
|
||||||
double CD;//charge density
|
double CD;//charge density
|
||||||
|
double CD_tmp;
|
||||||
double F = 96485.0;//Faraday's constant; unit[C/mol]; F=e*Na, where Na is the Avogadro constant
|
double F = 96485.0;//Faraday's constant; unit[C/mol]; F=e*Na, where Na is the Avogadro constant
|
||||||
for (n=start; n<finish; n++){
|
|
||||||
Ci = Den[n];
|
|
||||||
CD = F*IonValence[0]*Ci;
|
|
||||||
ChargeDensity[n] = CD;
|
|
||||||
}
|
|
||||||
|
|
||||||
ic = ic - 1;
|
for (n=start; n<finish; n++){
|
||||||
while (ic>0){
|
Ci = Den[n+ion_component*Np];
|
||||||
for (n=start; n<finish; n++){
|
CD = ChargeDensity[n];
|
||||||
Ci = Den[n+ic*Np];
|
CD_tmp = F*IonValence*Ci;
|
||||||
CD = F*IonValence[ic]*Ci;
|
ChargeDensity[n] = CD*(IonValence>0) + CD_tmp;
|
||||||
ChargeDensity[n] += CD;
|
|
||||||
}
|
|
||||||
ic--;
|
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|||||||
@@ -1,14 +1,14 @@
|
|||||||
/*
|
/*
|
||||||
* Multi-relaxation time LBM Model
|
* Dilute Ion Transport LBM Model
|
||||||
*/
|
*/
|
||||||
#include "models/IonModel.h"
|
#include "models/IonModel.h"
|
||||||
#include "analysis/distance.h"
|
#include "analysis/distance.h"
|
||||||
#include "common/ReadMicroCT.h"
|
#include "common/ReadMicroCT.h"
|
||||||
|
|
||||||
ScaLBL_IonModel::ScaLBL_IonModel(int RANK, int NP, MPI_Comm COMM):
|
ScaLBL_IonModel::ScaLBL_IonModel(int RANK, int NP, MPI_Comm COMM):
|
||||||
rank(RANK), nprocs(NP), Restart(0),timestep(0),timestepMax(0),tau(0),
|
rank(RANK),nprocs(NP),timestep(0),timestepMax(0),time_conv(0),kb(0),electron_charge(0),T(0),Vt(0),k2_inv(0),h(0),
|
||||||
Fx(0),Fy(0),Fz(0),flux(0),din(0),dout(0),mu(0),
|
tolerance(0),number_ion_species(0),Nx(0),Ny(0),Nz(0),N(0),Np(0),nprocx(0),nprocy(0),nprocz(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)
|
BoundaryCondition(0),BoundaryConditionSolid(0),Lx(0),Ly(0),Lz(0),comm(COMM)
|
||||||
{
|
{
|
||||||
|
|
||||||
}
|
}
|
||||||
@@ -44,8 +44,8 @@ void ScaLBL_IonModel::ReadParams(string filename,int num_iter,int num_iter_Stoke
|
|||||||
IonValence.push_back(1);//algebraic valence charge
|
IonValence.push_back(1);//algebraic valence charge
|
||||||
IonConcentration.push_back(1.0e-3);//user-input ion concentration has physical unit [mol/m^3]
|
IonConcentration.push_back(1.0e-3);//user-input ion concentration has physical unit [mol/m^3]
|
||||||
//deltaT.push_back(1.0);
|
//deltaT.push_back(1.0);
|
||||||
//tau.push_back(0.5+k2_inv*deltaT[0]*IonDiffusivisty[0]);
|
//tau.push_back(0.5+k2_inv*deltaT[0]*IonDiffusivity[0]);
|
||||||
tau.push_back(0.5+k2_inv*time_conv/(h*1.0e-6)/(h*1.0e-6)*IonDiffusivisty[0]);
|
tau.push_back(0.5+k2_inv*time_conv/(h*1.0e-6)/(h*1.0e-6)*IonDiffusivity[0]);
|
||||||
//--------------------------------------------------------------------------//
|
//--------------------------------------------------------------------------//
|
||||||
|
|
||||||
// LB-Ion Model parameters
|
// LB-Ion Model parameters
|
||||||
@@ -74,10 +74,10 @@ void ScaLBL_IonModel::ReadParams(string filename,int num_iter,int num_iter_Stoke
|
|||||||
}
|
}
|
||||||
else{
|
else{
|
||||||
for (int i=0; i<IonDiffusivity.size();i++){
|
for (int i=0; i<IonDiffusivity.size();i++){
|
||||||
//First, re-calculate tau
|
//First, convert ion diffusivity in physical unit to LB unit
|
||||||
tau[i] = 0.5+k2_inv*time_conv/(h*h*1.0e-12)*IonDiffusivisty[i];
|
|
||||||
//Second, convert ion diffusivity in physical unit to LB unit
|
|
||||||
IonDiffusivity[i] = IonDiffusivity[i]*time_conv/(h*h*1.0e-12);//LB diffusivity has unit [lu^2/lt]
|
IonDiffusivity[i] = IonDiffusivity[i]*time_conv/(h*h*1.0e-12);//LB diffusivity has unit [lu^2/lt]
|
||||||
|
//Second, re-calculate tau
|
||||||
|
tau[i] = 0.5+k2_inv*IonDiffusivity[i];
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@@ -103,28 +103,19 @@ void ScaLBL_IonModel::ReadParams(string filename,int num_iter,int num_iter_Stoke
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
// wrong...
|
|
||||||
//if (ion_db->keyExists( "deltaT" )){
|
|
||||||
// deltaT.clear();
|
|
||||||
// tau.clear();
|
|
||||||
// deltaT = ion_db->getVector<double>( "deltaT" );
|
|
||||||
// if (deltaT.size()!=number_ion_species){
|
|
||||||
// ERROR("Error: number_ion_species and deltaT must be the same length! \n");
|
|
||||||
// }
|
|
||||||
// else{//update relaxation parameter tau
|
|
||||||
// for (int i=0;i<deltaT.size();i++){
|
|
||||||
// tau[i]=0.5+k2_inv*deltaT[i]*IonDiffusivity[i];
|
|
||||||
// }
|
|
||||||
// }
|
|
||||||
//}
|
|
||||||
|
|
||||||
// Read domain parameters
|
// Read domain parameters
|
||||||
if (domain_db->keyExists( "voxel_length" )){//default unit: um/lu
|
if (domain_db->keyExists( "voxel_length" )){//default unit: um/lu
|
||||||
h = domain_db->getScalar<double>( "voxel_length" );
|
h = domain_db->getScalar<double>( "voxel_length" );
|
||||||
}
|
}
|
||||||
|
BoundaryCondition = 0;
|
||||||
if (domain_db->keyExists( "BC" )){
|
if (domain_db->keyExists( "BC" )){
|
||||||
BoundaryCondition = domain_db->getScalar<int>( "BC" );
|
BoundaryCondition = domain_db->getScalar<int>( "BC" );
|
||||||
}
|
}
|
||||||
|
BoundaryConditionSolid = 0;
|
||||||
|
if (domain_db->keyExists( "BC_Solid" )){
|
||||||
|
BoundaryConditionSolid = domain_db->getScalar<int>( "BC_Solid" );
|
||||||
|
}
|
||||||
|
|
||||||
if (rank==0) printf("*****************************************************\n");
|
if (rank==0) printf("*****************************************************\n");
|
||||||
if (rank==0) printf("LB Ion Transport Solver: \n");
|
if (rank==0) printf("LB Ion Transport Solver: \n");
|
||||||
@@ -134,6 +125,18 @@ void ScaLBL_IonModel::ReadParams(string filename,int num_iter,int num_iter_Stoke
|
|||||||
if (rank==0) printf(" Ion %i: LB relaxation tau = %.5g\n", i+1,tau[i]);
|
if (rank==0) printf(" Ion %i: LB relaxation tau = %.5g\n", i+1,tau[i]);
|
||||||
}
|
}
|
||||||
if (rank==0) printf("*****************************************************\n");
|
if (rank==0) printf("*****************************************************\n");
|
||||||
|
|
||||||
|
switch (BoundaryConditionSolid){
|
||||||
|
case 0:
|
||||||
|
if (rank==0) printf("LB Ion Solver: solid boundary: non-flux boundary is assigned");
|
||||||
|
break;
|
||||||
|
case 1:
|
||||||
|
if (rank==0) printf("LB Ion Solver: solid boundary: Neumann-type surfacen ion concentration is assigned");
|
||||||
|
break;
|
||||||
|
default:
|
||||||
|
if (rank==0) printf("LB Ion Solver: solid boundary: non-flux boundary is assigned");
|
||||||
|
break;
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
void ScaLBL_IonModel::SetDomain(){
|
void ScaLBL_IonModel::SetDomain(){
|
||||||
@@ -220,6 +223,64 @@ void ScaLBL_IonModel::ReadInput(){
|
|||||||
if (rank == 0) cout << " Domain set." << endl;
|
if (rank == 0) cout << " Domain set." << endl;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
void ScaLBL_IonModel::AssignSolidBoundary(double *ion_solid)
|
||||||
|
{
|
||||||
|
size_t NLABELS=0;
|
||||||
|
signed char VALUE=0;
|
||||||
|
double AFFINITY=0.f;
|
||||||
|
|
||||||
|
auto LabelList = ion_db->getVector<int>( "SolidLabels" );
|
||||||
|
auto AffinityList = ion_db->getVector<double>( "SolidValues" );
|
||||||
|
|
||||||
|
NLABELS=LabelList.size();
|
||||||
|
if (NLABELS != AffinityList.size()){
|
||||||
|
ERROR("Error: LB Ion Solver: SolidLabels and SolidValues 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=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=Mask->id[n];
|
||||||
|
AFFINITY=0.f;
|
||||||
|
// 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];
|
||||||
|
//NOTE need to convert the user input phys unit to LB unit
|
||||||
|
AFFINITY = AFFINITY*(h*h*1.0e-12);
|
||||||
|
label_count[idx] += 1.0;
|
||||||
|
idx = NLABELS;
|
||||||
|
//Mask->id[n] = 0; // set mask to zero since this is an immobile component
|
||||||
|
}
|
||||||
|
}
|
||||||
|
ion_solid[n] = AFFINITY;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
for (size_t idx=0; idx<NLABELS; idx++)
|
||||||
|
label_count_global[idx]=sumReduce( Dm->Comm, label_count[idx]);
|
||||||
|
|
||||||
|
if (rank==0){
|
||||||
|
printf("LB Ion Solver: Ion Solid 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, surface ion concentration=%.3g [mol/m^2], volume fraction=%.2g\n",VALUE,AFFINITY,volume_fraction);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
void ScaLBL_IonModel::Create(){
|
void ScaLBL_IonModel::Create(){
|
||||||
/*
|
/*
|
||||||
* This function creates the variables needed to run a LBM
|
* This function creates the variables needed to run a LBM
|
||||||
@@ -262,6 +323,23 @@ void ScaLBL_IonModel::Create(){
|
|||||||
ScaLBL_CopyToDevice(NeighborList, neighborList, neighborSize);
|
ScaLBL_CopyToDevice(NeighborList, neighborList, neighborSize);
|
||||||
MPI_Barrier(comm);
|
MPI_Barrier(comm);
|
||||||
|
|
||||||
|
//Initialize solid boundary for electrical potential
|
||||||
|
//if ion concentration at solid surface is specified
|
||||||
|
if (BoundaryConditionSolid==1){
|
||||||
|
|
||||||
|
ScaLBL_AllocateDeviceMemory((void **) &IonSolid, sizeof(double)*Nx*Ny*Nz);
|
||||||
|
ScaLBL_Comm->SetupBounceBackList(Map, Mask->id, Np);
|
||||||
|
MPI_Barrier(comm);
|
||||||
|
|
||||||
|
double *IonSolid_host;
|
||||||
|
IonSolid_host = new double[Nx*Ny*Nz];
|
||||||
|
AssignSolidBoundary(IonSolid_host);
|
||||||
|
ScaLBL_CopyToDevice(IonSolid, IonSolid_host, Nx*Ny*Nz*sizeof(double));
|
||||||
|
ScaLBL_DeviceBarrier();
|
||||||
|
delete [] IonSolid_host;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
}
|
}
|
||||||
|
|
||||||
void ScaLBL_IonModel::Initialize(){
|
void ScaLBL_IonModel::Initialize(){
|
||||||
@@ -273,8 +351,10 @@ void ScaLBL_IonModel::Initialize(){
|
|||||||
ScaLBL_D3Q7_Ion_Init(&fq[ic*Np*7],&Ci[ic*Np],IonConcentration[ic],Np);
|
ScaLBL_D3Q7_Ion_Init(&fq[ic*Np*7],&Ci[ic*Np],IonConcentration[ic],Np);
|
||||||
}
|
}
|
||||||
if (rank==0) printf ("LB Ion Solver: initializing charge density\n");
|
if (rank==0) printf ("LB Ion Solver: initializing charge density\n");
|
||||||
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence, number_ion_species, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
for (int ic=0; ic<number_ion_species; ic++){
|
||||||
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence, number_ion_species, 0, ScaLBL_Comm->LastExterior(), Np);
|
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence[ic], ic, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
||||||
|
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence[ic], ic, 0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
void ScaLBL_IonModel::Run(double *Velocity, double *ElectricField){
|
void ScaLBL_IonModel::Run(double *Velocity, double *ElectricField){
|
||||||
@@ -299,10 +379,11 @@ void ScaLBL_IonModel::Run(double *Velocity, double *ElectricField){
|
|||||||
ScaLBL_Comm->SendD3Q7AA(fq, ic); //READ FROM NORMAL
|
ScaLBL_Comm->SendD3Q7AA(fq, ic); //READ FROM NORMAL
|
||||||
ScaLBL_D3Q7_AAodd_IonConcentration(NeighborList, &fq[ic*Np*7],&Ci[ic*Np],ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
ScaLBL_D3Q7_AAodd_IonConcentration(NeighborList, &fq[ic*Np*7],&Ci[ic*Np],ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
||||||
ScaLBL_Comm->RecvD3Q7AA(fq, ic); //WRITE INTO OPPOSITE
|
ScaLBL_Comm->RecvD3Q7AA(fq, ic); //WRITE INTO OPPOSITE
|
||||||
|
ScaLBL_DeviceBarrier();
|
||||||
ScaLBL_D3Q7_AAodd_IonConcentration(NeighborList, &fq[ic*Np*7],&Ci[ic*Np], 0, ScaLBL_Comm->LastExterior(), Np);
|
ScaLBL_D3Q7_AAodd_IonConcentration(NeighborList, &fq[ic*Np*7],&Ci[ic*Np], 0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
|
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence[ic], ic, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
||||||
|
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence[ic], ic, 0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
}
|
}
|
||||||
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence, number_ion_species, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
|
||||||
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence, number_ion_species, 0, ScaLBL_Comm->LastExterior(), Np);
|
|
||||||
|
|
||||||
//LB-Ion collison
|
//LB-Ion collison
|
||||||
for (int ic=0; ic<number_ion_species; ic++){
|
for (int ic=0; ic<number_ion_species; ic++){
|
||||||
@@ -317,7 +398,13 @@ void ScaLBL_IonModel::Run(double *Velocity, double *ElectricField){
|
|||||||
ScaLBL_D3Q7_AAodd_Ion(NeighborList, &fq[ic*Np*7],&Ci[ic*Np],Velocity,ElectricField,IonDiffusivity[ic],IonValence[ic],
|
ScaLBL_D3Q7_AAodd_Ion(NeighborList, &fq[ic*Np*7],&Ci[ic*Np],Velocity,ElectricField,IonDiffusivity[ic],IonValence[ic],
|
||||||
rlx[ic],Vt,0, ScaLBL_Comm->LastExterior(), Np);
|
rlx[ic],Vt,0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
}
|
}
|
||||||
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
|
if (BoundaryConditionSolid==1){
|
||||||
|
for (int ic=0; ic<number_ion_species; ic++){
|
||||||
|
//TODO IonSolid may also be species-dependent
|
||||||
|
ScaLBL_Comm->SolidNeumannD3Q7(&fq[ic*Np*7], IonSolid);
|
||||||
|
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
// *************EVEN TIMESTEP*************//
|
// *************EVEN TIMESTEP*************//
|
||||||
timestep++;
|
timestep++;
|
||||||
@@ -327,9 +414,9 @@ void ScaLBL_IonModel::Run(double *Velocity, double *ElectricField){
|
|||||||
ScaLBL_D3Q7_AAeven_IonConcentration(&fq[ic*Np*7],&Ci[ic*Np],ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
ScaLBL_D3Q7_AAeven_IonConcentration(&fq[ic*Np*7],&Ci[ic*Np],ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
||||||
ScaLBL_Comm->RecvD3Q7AA(fq, ic); //WRITE INTO OPPOSITE
|
ScaLBL_Comm->RecvD3Q7AA(fq, ic); //WRITE INTO OPPOSITE
|
||||||
ScaLBL_D3Q7_AAeven_IonConcentration(&fq[ic*Np*7],&Ci[ic*Np], 0, ScaLBL_Comm->LastExterior(), Np);
|
ScaLBL_D3Q7_AAeven_IonConcentration(&fq[ic*Np*7],&Ci[ic*Np], 0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
|
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence[ic], ic, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
||||||
|
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence[ic], ic, 0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
}
|
}
|
||||||
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence, number_ion_species, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
|
||||||
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, IonValence, number_ion_species, 0, ScaLBL_Comm->LastExterior(), Np);
|
|
||||||
|
|
||||||
//LB-Ion collison
|
//LB-Ion collison
|
||||||
for (int ic=0; ic<number_ion_species; ic++){
|
for (int ic=0; ic<number_ion_species; ic++){
|
||||||
@@ -344,7 +431,13 @@ void ScaLBL_IonModel::Run(double *Velocity, double *ElectricField){
|
|||||||
ScaLBL_D3Q7_AAeven_Ion(&fq[ic*Np*7],&Ci[ic*Np],Velocity,ElectricField,IonDiffusivity[ic],IonValence[ic],
|
ScaLBL_D3Q7_AAeven_Ion(&fq[ic*Np*7],&Ci[ic*Np],Velocity,ElectricField,IonDiffusivity[ic],IonValence[ic],
|
||||||
rlx[ic],Vt,0, ScaLBL_Comm->LastExterior(), Np);
|
rlx[ic],Vt,0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
}
|
}
|
||||||
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
|
if (BoundaryConditionSolid==1){
|
||||||
|
for (int ic=0; ic<number_ion_species; ic++){
|
||||||
|
//TODO IonSolid may also be species-dependent
|
||||||
|
ScaLBL_Comm->SolidNeumannD3Q7(&fq[ic*Np*7], IonSolid);
|
||||||
|
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
|
||||||
|
}
|
||||||
|
}
|
||||||
//************************************************************************/
|
//************************************************************************/
|
||||||
}
|
}
|
||||||
//************************************************************************/
|
//************************************************************************/
|
||||||
|
|||||||
@@ -30,9 +30,10 @@ public:
|
|||||||
void Initialize();
|
void Initialize();
|
||||||
void Run(double *Velocity, double *ElectricField);
|
void Run(double *Velocity, double *ElectricField);
|
||||||
|
|
||||||
bool Restart,pBC;
|
//bool Restart,pBC;
|
||||||
int timestep,timestepMax;
|
int timestep,timestepMax;
|
||||||
int BoundaryCondition;
|
int BoundaryCondition;
|
||||||
|
int BoundaryConditionSolid;
|
||||||
double h;//domain resolution, unit [um/lu]
|
double h;//domain resolution, unit [um/lu]
|
||||||
double time_conv;
|
double time_conv;
|
||||||
double kb,electron_charge,T,Vt;
|
double kb,electron_charge,T,Vt;
|
||||||
@@ -64,6 +65,7 @@ public:
|
|||||||
double *fq;
|
double *fq;
|
||||||
double *Ci;
|
double *Ci;
|
||||||
double *ChargeDensity;
|
double *ChargeDensity;
|
||||||
|
double *IonSolid;
|
||||||
|
|
||||||
private:
|
private:
|
||||||
MPI_Comm comm;
|
MPI_Comm comm;
|
||||||
@@ -75,4 +77,5 @@ private:
|
|||||||
|
|
||||||
//int rank,nprocs;
|
//int rank,nprocs;
|
||||||
void LoadParams(std::shared_ptr<Database> db0);
|
void LoadParams(std::shared_ptr<Database> db0);
|
||||||
|
void AssignSolidBoundary(double *ion_solid);
|
||||||
};
|
};
|
||||||
|
|||||||
@@ -1,7 +1,7 @@
|
|||||||
#include "models/MultiPhysController.h"
|
#include "models/MultiPhysController.h"
|
||||||
|
|
||||||
ScaLBL_Multiphys_Controller::ScaLBL_Multiphys_Controller(int RANK, int NP, MPI_Comm COMM):
|
ScaLBL_Multiphys_Controller::ScaLBL_Multiphys_Controller(int RANK, int NP, MPI_Comm COMM):
|
||||||
rank(RANK),Restart(0),timestepMax(0),num_iter_Stokes(0),num_iter_Ion(0),SchmidtNum(0),comm(COMM)
|
rank(RANK),nprocs(NP),Restart(0),timestepMax(0),num_iter_Stokes(0),num_iter_Ion(0),SchmidtNum(0),comm(COMM)
|
||||||
{
|
{
|
||||||
|
|
||||||
}
|
}
|
||||||
|
|||||||
@@ -30,6 +30,8 @@ public:
|
|||||||
int num_iter_Ion;
|
int num_iter_Ion;
|
||||||
double SchmidtNum;//Schmidt number = kinematic_viscosity/mass_diffusivity
|
double SchmidtNum;//Schmidt number = kinematic_viscosity/mass_diffusivity
|
||||||
|
|
||||||
|
int rank,nprocs;
|
||||||
|
|
||||||
// input database
|
// input database
|
||||||
std::shared_ptr<Database> db;
|
std::shared_ptr<Database> db;
|
||||||
std::shared_ptr<Database> study_db;
|
std::shared_ptr<Database> study_db;
|
||||||
|
|||||||
@@ -6,9 +6,9 @@
|
|||||||
#include "common/ReadMicroCT.h"
|
#include "common/ReadMicroCT.h"
|
||||||
|
|
||||||
ScaLBL_Poisson::ScaLBL_Poisson(int RANK, int NP, MPI_Comm COMM):
|
ScaLBL_Poisson::ScaLBL_Poisson(int RANK, int NP, MPI_Comm COMM):
|
||||||
rank(RANK), nprocs(NP), Restart(0),timestep(0),timestepMax(0),tau(0),
|
rank(RANK), nprocs(NP),timestep(0),timestepMax(0),tau(0),k2_inv(0),gamma(0),tolerance(0),h(0),
|
||||||
Fx(0),Fy(0),Fz(0),flux(0),din(0),dout(0),mu(0),
|
epsilon0(0),epsilon0_LB(0),epsilonR(0),epsilon_LB(0),Nx(0),Ny(0),Nz(0),N(0),Np(0),analysis_interval(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)
|
nprocx(0),nprocy(0),nprocz(0),BoundaryCondition(0),BoundaryConditionSolid(0),Lx(0),Ly(0),Lz(0),comm(COMM)
|
||||||
{
|
{
|
||||||
|
|
||||||
}
|
}
|
||||||
@@ -54,9 +54,16 @@ void ScaLBL_Poisson::ReadParams(string filename){
|
|||||||
if (domain_db->keyExists( "voxel_length" )){//default unit: um/lu
|
if (domain_db->keyExists( "voxel_length" )){//default unit: um/lu
|
||||||
h = domain_db->getScalar<double>( "voxel_length" );
|
h = domain_db->getScalar<double>( "voxel_length" );
|
||||||
}
|
}
|
||||||
|
|
||||||
|
BoundaryCondition = 0;
|
||||||
if (domain_db->keyExists( "BC" )){
|
if (domain_db->keyExists( "BC" )){
|
||||||
BoundaryCondition = domain_db->getScalar<int>( "BC" );
|
BoundaryCondition = domain_db->getScalar<int>( "BC" );
|
||||||
}
|
}
|
||||||
|
BoundaryConditionSolid = 1;
|
||||||
|
if (domain_db->keyExists( "BC_Solid" )){
|
||||||
|
BoundaryConditionSolid = domain_db->getScalar<int>( "BC_Solid" );
|
||||||
|
}
|
||||||
|
|
||||||
//Re-calcualte model parameters if user updates input
|
//Re-calcualte model parameters if user updates input
|
||||||
epsilon0_LB = epsilon0*(h*1.0e-6);//unit:[C/(V*lu)]
|
epsilon0_LB = epsilon0*(h*1.0e-6);//unit:[C/(V*lu)]
|
||||||
epsilon_LB = epsilon0_LB*epsilonR;//electrical permittivity
|
epsilon_LB = epsilon0_LB*epsilonR;//electrical permittivity
|
||||||
@@ -66,6 +73,18 @@ void ScaLBL_Poisson::ReadParams(string filename){
|
|||||||
if (rank==0) printf("LB-Poisson Solver: steady-state MaxTimeStep = %i; steady-state tolerance = %.3g \n", timestepMax,tolerance);
|
if (rank==0) printf("LB-Poisson Solver: steady-state MaxTimeStep = %i; steady-state tolerance = %.3g \n", timestepMax,tolerance);
|
||||||
if (rank==0) printf(" LB relaxation tau = %.5g \n", tau);
|
if (rank==0) printf(" LB relaxation tau = %.5g \n", tau);
|
||||||
if (rank==0) printf("***********************************************************************************\n");
|
if (rank==0) printf("***********************************************************************************\n");
|
||||||
|
|
||||||
|
switch (BoundaryConditionSolid){
|
||||||
|
case 1:
|
||||||
|
if (rank==0) printf("LB-Poisson Solver: solid boundary: Dirichlet-type surfacen potential is assigned");
|
||||||
|
break;
|
||||||
|
case 2:
|
||||||
|
if (rank==0) printf("LB-Poisson Solver: solid boundary: Neumann-type surfacen charge density is assigned");
|
||||||
|
break;
|
||||||
|
default:
|
||||||
|
if (rank==0) printf("LB-Poisson Solver: solid boundary: Dirichlet-type surfacen potential is assigned");
|
||||||
|
break;
|
||||||
|
}
|
||||||
}
|
}
|
||||||
void ScaLBL_Poisson::SetDomain(){
|
void ScaLBL_Poisson::SetDomain(){
|
||||||
Dm = std::shared_ptr<Domain>(new Domain(domain_db,comm)); // full domain for analysis
|
Dm = std::shared_ptr<Domain>(new Domain(domain_db,comm)); // full domain for analysis
|
||||||
@@ -152,6 +171,77 @@ void ScaLBL_Poisson::ReadInput(){
|
|||||||
if (rank == 0) cout << " Domain set." << endl;
|
if (rank == 0) cout << " Domain set." << endl;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
void ScaLBL_Poisson::AssignSolidBoundary(double *poisson_solid)
|
||||||
|
{
|
||||||
|
size_t NLABELS=0;
|
||||||
|
signed char VALUE=0;
|
||||||
|
double AFFINITY=0.f;
|
||||||
|
|
||||||
|
auto LabelList = electric_db->getVector<int>( "SolidLabels" );
|
||||||
|
auto AffinityList = electric_db->getVector<double>( "SolidValues" );
|
||||||
|
|
||||||
|
NLABELS=LabelList.size();
|
||||||
|
if (NLABELS != AffinityList.size()){
|
||||||
|
ERROR("Error: LB-Poisson Solver: SolidLabels and SolidValues 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=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=Mask->id[n];
|
||||||
|
AFFINITY=0.f;
|
||||||
|
// 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];
|
||||||
|
//NOTE need to convert the user input phys unit to LB unit
|
||||||
|
if (BoundaryConditionSolid==2){
|
||||||
|
//for BCS=1, i.e. Dirichlet-type, no need for unit conversion
|
||||||
|
AFFINITY = AFFINITY*(h*h*1.0e-12);
|
||||||
|
}
|
||||||
|
label_count[idx] += 1.0;
|
||||||
|
idx = NLABELS;
|
||||||
|
//Mask->id[n] = 0; // set mask to zero since this is an immobile component
|
||||||
|
}
|
||||||
|
}
|
||||||
|
poisson_solid[n] = AFFINITY;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
for (size_t idx=0; idx<NLABELS; idx++)
|
||||||
|
label_count_global[idx]=sumReduce( Dm->Comm, label_count[idx]);
|
||||||
|
|
||||||
|
if (rank==0){
|
||||||
|
printf("LB-Poisson Solver: Poisson Solid 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);
|
||||||
|
switch (BoundaryConditionSolid){
|
||||||
|
case 1:
|
||||||
|
printf(" label=%d, surface potential=%.3g [V], volume fraction=%.2g\n",VALUE,AFFINITY,volume_fraction);
|
||||||
|
break;
|
||||||
|
case 2:
|
||||||
|
printf(" label=%d, surface charge density=%.3g [C/m^2], volume fraction=%.2g\n",VALUE,AFFINITY,volume_fraction);
|
||||||
|
break;
|
||||||
|
default:
|
||||||
|
printf(" label=%d, surface potential=%.3g [V], volume fraction=%.2g\n",VALUE,AFFINITY,volume_fraction);
|
||||||
|
break;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
void ScaLBL_Poisson::Create(){
|
void ScaLBL_Poisson::Create(){
|
||||||
/*
|
/*
|
||||||
* This function creates the variables needed to run a LBM
|
* This function creates the variables needed to run a LBM
|
||||||
@@ -187,22 +277,26 @@ void ScaLBL_Poisson::Create(){
|
|||||||
ScaLBL_AllocateDeviceMemory((void **) &fq, 7*dist_mem_size);
|
ScaLBL_AllocateDeviceMemory((void **) &fq, 7*dist_mem_size);
|
||||||
ScaLBL_AllocateDeviceMemory((void **) &Psi, sizeof(double)*Np);
|
ScaLBL_AllocateDeviceMemory((void **) &Psi, sizeof(double)*Np);
|
||||||
ScaLBL_AllocateDeviceMemory((void **) &ElectricField, 3*sizeof(double)*Np);
|
ScaLBL_AllocateDeviceMemory((void **) &ElectricField, 3*sizeof(double)*Np);
|
||||||
ScaLBL_AllocateDeviceMemory((void **) &zeta, sizeof(double)*ScaLBL_Comm->n_bb_d3q7);
|
ScaLBL_AllocateDeviceMemory((void **) &PoissonSolid, sizeof(double)*Nx*Ny*Nz);
|
||||||
//...........................................................................
|
//...........................................................................
|
||||||
// initialize the zeta function (example is zeta is constant on solid surface)
|
|
||||||
double *tmpZeta = new double[ScaLBL_Comm->n_bb_d3q7];
|
|
||||||
for int (i=0; i<ScaLBL_Comm->n_bb_d3q7; i++){
|
|
||||||
tmpZeta[i] = 1.0/k2_inv; // this has to be read from input file
|
|
||||||
}
|
|
||||||
ScaLBL_CopyToDevice(zeta, tmpZeta, sizeof(double)*ScaLBL_Comm->n_bb_d3q7);
|
|
||||||
delete [] tmpZeta;
|
|
||||||
|
|
||||||
// Update GPU data structures
|
// Update GPU data structures
|
||||||
if (rank==0) printf ("LB-Poisson Solver: Setting up device map and neighbor list \n");
|
if (rank==0) printf ("LB-Poisson Solver: Setting up device map and neighbor list \n");
|
||||||
// copy the neighbor list
|
// copy the neighbor list
|
||||||
ScaLBL_CopyToDevice(NeighborList, neighborList, neighborSize);
|
ScaLBL_CopyToDevice(NeighborList, neighborList, neighborSize);
|
||||||
|
ScaLBL_DeviceBarrier();
|
||||||
MPI_Barrier(comm);
|
MPI_Barrier(comm);
|
||||||
|
|
||||||
|
//Initialize solid boundary for electrical potential
|
||||||
|
ScaLBL_Comm->SetupBounceBackList(Map, Mask->id, Np);
|
||||||
|
MPI_Barrier(comm);
|
||||||
|
|
||||||
|
double *PoissonSolid_host;
|
||||||
|
PoissonSolid_host = new double[Nx*Ny*Nz];
|
||||||
|
AssignSolidBoundary(PoissonSolid_host);
|
||||||
|
ScaLBL_CopyToDevice(PoissonSolid, PoissonSolid_host, Nx*Ny*Nz*sizeof(double));
|
||||||
|
ScaLBL_DeviceBarrier();
|
||||||
|
delete [] PoissonSolid_host;
|
||||||
}
|
}
|
||||||
|
|
||||||
void ScaLBL_Poisson::Initialize(){
|
void ScaLBL_Poisson::Initialize(){
|
||||||
@@ -233,7 +327,7 @@ void ScaLBL_Poisson::Run(double *ChargeDensity){
|
|||||||
// Set boundary conditions
|
// Set boundary conditions
|
||||||
/* ... */
|
/* ... */
|
||||||
ScaLBL_D3Q7_AAodd_Poisson(NeighborList, fq, ChargeDensity, Psi, ElectricField, tau, epsilon_LB, gamma, 0, ScaLBL_Comm->LastExterior(), Np);
|
ScaLBL_D3Q7_AAodd_Poisson(NeighborList, fq, ChargeDensity, Psi, ElectricField, tau, epsilon_LB, gamma, 0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
ScaLBL_Comm->SolidDirichletD3Q7(fq, zeta);
|
ScaLBL_Comm->SolidDirichletD3Q7(fq, PoissonSolid);
|
||||||
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
|
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
|
||||||
|
|
||||||
// *************EVEN TIMESTEP*************//
|
// *************EVEN TIMESTEP*************//
|
||||||
@@ -244,14 +338,14 @@ void ScaLBL_Poisson::Run(double *ChargeDensity){
|
|||||||
// Set boundary conditions
|
// Set boundary conditions
|
||||||
/* ... */
|
/* ... */
|
||||||
ScaLBL_D3Q7_AAeven_Poisson(fq, ChargeDensity, Psi, ElectricField, tau, epsilon_LB, gamma, 0, ScaLBL_Comm->LastExterior(), Np);
|
ScaLBL_D3Q7_AAeven_Poisson(fq, ChargeDensity, Psi, ElectricField, tau, epsilon_LB, gamma, 0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
ScaLBL_Comm->SolidDirichletD3Q7(fq, zeta);
|
ScaLBL_Comm->SolidDirichletD3Q7(fq, PoissonSolid);
|
||||||
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
|
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
|
||||||
//************************************************************************/
|
//************************************************************************/
|
||||||
|
|
||||||
// Check convergence of steady-state solution
|
// Check convergence of steady-state solution
|
||||||
if (timestep%analysis_interval==0){
|
if (timestep%analysis_interval==0){
|
||||||
|
|
||||||
ScaLBL_Comm->RegularLayout(Map,&Psi,Psi_host);
|
ScaLBL_Comm->RegularLayout(Map,Psi,Psi_host);
|
||||||
double count_loc=0;
|
double count_loc=0;
|
||||||
double count;
|
double count;
|
||||||
double psi_avg;
|
double psi_avg;
|
||||||
|
|||||||
@@ -33,9 +33,10 @@ public:
|
|||||||
int timestep,timestepMax;
|
int timestep,timestepMax;
|
||||||
int analysis_interval;
|
int analysis_interval;
|
||||||
int BoundaryCondition;
|
int BoundaryCondition;
|
||||||
|
int BoundaryConditionSolid;
|
||||||
double tau;
|
double tau;
|
||||||
double tolerance;
|
double tolerance;
|
||||||
double k2_inv,deltaT;
|
double k2_inv,gamma;
|
||||||
double epsilon0,epsilon0_LB,epsilonR,epsilon_LB;
|
double epsilon0,epsilon0_LB,epsilonR,epsilon_LB;
|
||||||
|
|
||||||
int Nx,Ny,Nz,N,Np;
|
int Nx,Ny,Nz,N,Np;
|
||||||
@@ -58,7 +59,7 @@ public:
|
|||||||
double *fq;
|
double *fq;
|
||||||
double *Psi;
|
double *Psi;
|
||||||
double *ElectricField;
|
double *ElectricField;
|
||||||
double *zeta;
|
double *PoissonSolid;
|
||||||
|
|
||||||
private:
|
private:
|
||||||
MPI_Comm comm;
|
MPI_Comm comm;
|
||||||
@@ -70,4 +71,5 @@ private:
|
|||||||
|
|
||||||
//int rank,nprocs;
|
//int rank,nprocs;
|
||||||
void LoadParams(std::shared_ptr<Database> db0);
|
void LoadParams(std::shared_ptr<Database> db0);
|
||||||
|
void AssignSolidBoundary(double *poisson_solid);
|
||||||
};
|
};
|
||||||
|
|||||||
@@ -7,7 +7,7 @@
|
|||||||
|
|
||||||
ScaLBL_StokesModel::ScaLBL_StokesModel(int RANK, int NP, MPI_Comm COMM):
|
ScaLBL_StokesModel::ScaLBL_StokesModel(int RANK, int NP, MPI_Comm COMM):
|
||||||
rank(RANK), nprocs(NP), Restart(0),timestep(0),timestepMax(0),tau(0),
|
rank(RANK), nprocs(NP), Restart(0),timestep(0),timestepMax(0),tau(0),
|
||||||
Fx(0),Fy(0),Fz(0),flux(0),din(0),dout(0),mu(0),
|
Fx(0),Fy(0),Fz(0),flux(0),din(0),dout(0),mu(0),h(0),nu_phys(0),time_conv(0),tolerance(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)
|
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)
|
||||||
{
|
{
|
||||||
|
|
||||||
|
|||||||
@@ -21,13 +21,14 @@ public:
|
|||||||
~ScaLBL_StokesModel();
|
~ScaLBL_StokesModel();
|
||||||
|
|
||||||
// functions in they should be run
|
// functions in they should be run
|
||||||
void ReadParams(string filename);
|
void ReadParams(string filename,int num_iter);
|
||||||
void ReadParams(std::shared_ptr<Database> db0);
|
void ReadParams(std::shared_ptr<Database> db0);
|
||||||
void SetDomain();
|
void SetDomain();
|
||||||
void ReadInput();
|
void ReadInput();
|
||||||
void Create();
|
void Create();
|
||||||
void Initialize();
|
void Initialize();
|
||||||
void Run();
|
void Run();
|
||||||
|
void Run_Lite(double *ChargeDensity, double *ElectricField);
|
||||||
void VelocityField();
|
void VelocityField();
|
||||||
|
|
||||||
bool Restart,pBC;
|
bool Restart,pBC;
|
||||||
@@ -37,6 +38,9 @@ public:
|
|||||||
double Fx,Fy,Fz,flux;
|
double Fx,Fy,Fz,flux;
|
||||||
double din,dout;
|
double din,dout;
|
||||||
double tolerance;
|
double tolerance;
|
||||||
|
double nu_phys;
|
||||||
|
double time_conv;
|
||||||
|
double h;//image resolution
|
||||||
|
|
||||||
int Nx,Ny,Nz,N,Np;
|
int Nx,Ny,Nz,N,Np;
|
||||||
int rank,nprocx,nprocy,nprocz,nprocs;
|
int rank,nprocx,nprocy,nprocz,nprocs;
|
||||||
|
|||||||
@@ -4,7 +4,7 @@
|
|||||||
ADD_LBPM_EXECUTABLE( lbpm_color_simulator )
|
ADD_LBPM_EXECUTABLE( lbpm_color_simulator )
|
||||||
ADD_LBPM_EXECUTABLE( lbpm_permeability_simulator )
|
ADD_LBPM_EXECUTABLE( lbpm_permeability_simulator )
|
||||||
ADD_LBPM_EXECUTABLE( lbpm_greyscale_simulator )
|
ADD_LBPM_EXECUTABLE( lbpm_greyscale_simulator )
|
||||||
ADD_LBPM_EXECUTABLE( lbpm_electrokinetic_dfh_simulator.cpp )
|
ADD_LBPM_EXECUTABLE( lbpm_electrokinetic_dfh_simulator )
|
||||||
#ADD_LBPM_EXECUTABLE( lbpm_BGK_simulator )
|
#ADD_LBPM_EXECUTABLE( lbpm_BGK_simulator )
|
||||||
#ADD_LBPM_EXECUTABLE( lbpm_color_macro_simulator )
|
#ADD_LBPM_EXECUTABLE( lbpm_color_macro_simulator )
|
||||||
ADD_LBPM_EXECUTABLE( lbpm_dfh_simulator )
|
ADD_LBPM_EXECUTABLE( lbpm_dfh_simulator )
|
||||||
|
|||||||
Reference in New Issue
Block a user