build pass; TODO: verify the slippingBC model
This commit is contained in:
parent
fee3d9eadc
commit
cd2bcfba85
@ -973,7 +973,7 @@ int ScaLBL_Communicator::MemoryOptimizedLayoutAA(IntArray &Map, int *neighborLis
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
void ScaLBL_Communicator::SetupBounceBackList(IntArray &Map, signed char *id, int Np, bool SlippingVelBC=false)
|
void ScaLBL_Communicator::SetupBounceBackList(IntArray &Map, signed char *id, int Np, bool SlippingVelBC)
|
||||||
{
|
{
|
||||||
|
|
||||||
int idx,i,j,k;
|
int idx,i,j,k;
|
||||||
@ -1362,10 +1362,10 @@ void ScaLBL_Communicator::SolidNeumannD3Q7(double *fq, double *BoundaryValue){
|
|||||||
}
|
}
|
||||||
|
|
||||||
void ScaLBL_Communicator::SolidSlippingVelocityBCD3Q19(double *fq, double *zeta_potential, double *ElectricField, double *SolidGrad,
|
void ScaLBL_Communicator::SolidSlippingVelocityBCD3Q19(double *fq, double *zeta_potential, double *ElectricField, double *SolidGrad,
|
||||||
double epslion_LB, double tau, double rho0, double den_scale,double h, double time_conv){
|
double epsilon_LB, double tau, double rho0, double den_scale,double h, double time_conv){
|
||||||
// fq is a D3Q19 distribution
|
// fq is a D3Q19 distribution
|
||||||
// BoundaryValues is a list of values to assign at bounce-back solid sites
|
// BoundaryValues is a list of values to assign at bounce-back solid sites
|
||||||
ScaLBL_Solid_SlippingVelocityBC_D3Q19(fq,zeta_potential,Electricfield,SolidGrad,epsilon_LB,tau,rho0,den_scale,h,time_conv,
|
ScaLBL_Solid_SlippingVelocityBC_D3Q19(fq,zeta_potential,ElectricField,SolidGrad,epsilon_LB,tau,rho0,den_scale,h,time_conv,
|
||||||
bb_dist,bb_interactions,fluid_boundary,lattice_weight,lattice_cx,lattice_cy,lattice_cz,n_bb_d3q19,N);
|
bb_dist,bb_interactions,fluid_boundary,lattice_weight,lattice_cx,lattice_cy,lattice_cz,n_bb_d3q19,N);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -49,6 +49,7 @@ __global__ void dvc_ScaLBL_Solid_SlippingVelocityBC_D3Q19(double *dist, double *
|
|||||||
double value_b,value_q;
|
double value_b,value_q;
|
||||||
double Ex,Ey,Ez;
|
double Ex,Ey,Ez;
|
||||||
double Etx,Ety,Etz;//tangential part of electric field
|
double Etx,Ety,Etz;//tangential part of electric field
|
||||||
|
double E_mag_normal;
|
||||||
double nsx,nsy,nsz;//unit normal solid gradient
|
double nsx,nsy,nsz;//unit normal solid gradient
|
||||||
double ubx,uby,ubz;//slipping velocity at fluid boundary nodes
|
double ubx,uby,ubz;//slipping velocity at fluid boundary nodes
|
||||||
float cx,cy,cz;//lattice velocity (D3Q19)
|
float cx,cy,cz;//lattice velocity (D3Q19)
|
||||||
@ -75,9 +76,9 @@ __global__ void dvc_ScaLBL_Solid_SlippingVelocityBC_D3Q19(double *dist, double *
|
|||||||
Etx = Ex - E_mag_normal*nsx;
|
Etx = Ex - E_mag_normal*nsx;
|
||||||
Ety = Ey - E_mag_normal*nsy;
|
Ety = Ey - E_mag_normal*nsy;
|
||||||
Etz = Ez - E_mag_normal*nsz;
|
Etz = Ez - E_mag_normal*nsz;
|
||||||
ubx = -eplison_LB*value_b*Etx/(nu_LB*rho0)*time_conv*time_conv/h/h/den_scale;
|
ubx = -epsilon_LB*value_b*Etx/(nu_LB*rho0)*time_conv*time_conv/h/h/den_scale;
|
||||||
uby = -eplison_LB*value_b*Ety/(nu_LB*rho0)*time_conv*time_conv/h/h/den_scale;
|
uby = -epsilon_LB*value_b*Ety/(nu_LB*rho0)*time_conv*time_conv/h/h/den_scale;
|
||||||
ubz = -eplison_LB*value_b*Etz/(nu_LB*rho0)*time_conv*time_conv/h/h/den_scale;
|
ubz = -epsilon_LB*value_b*Etz/(nu_LB*rho0)*time_conv*time_conv/h/h/den_scale;
|
||||||
|
|
||||||
//compute bounce-back distribution
|
//compute bounce-back distribution
|
||||||
LB_weight = lattice_weight[idx];
|
LB_weight = lattice_weight[idx];
|
||||||
|
@ -95,8 +95,8 @@ void ScaLBL_StokesModel::ReadParams(string filename,int num_iter){
|
|||||||
if (stokes_db->keyExists( "UseElectroosmoticVelocityBC" )){
|
if (stokes_db->keyExists( "UseElectroosmoticVelocityBC" )){
|
||||||
UseSlippingVelBC = stokes_db->getScalar<bool>( "UseElectroosmoticVelocityBC" );
|
UseSlippingVelBC = stokes_db->getScalar<bool>( "UseElectroosmoticVelocityBC" );
|
||||||
}
|
}
|
||||||
if (electric_db->keyExists( "epsilonR" )){
|
if (stokes_db->keyExists( "epsilonR" )){
|
||||||
epsilonR = electric_db->getScalar<double>( "epsilonR" );
|
epsilonR = stokes_db->getScalar<double>( "epsilonR" );
|
||||||
}
|
}
|
||||||
|
|
||||||
// Re-calculate model parameters due to parameter read
|
// Re-calculate model parameters due to parameter read
|
||||||
@ -182,8 +182,8 @@ void ScaLBL_StokesModel::ReadParams(string filename){
|
|||||||
if (stokes_db->keyExists( "UseElectroosmoticVelocityBC" )){
|
if (stokes_db->keyExists( "UseElectroosmoticVelocityBC" )){
|
||||||
UseSlippingVelBC = stokes_db->getScalar<bool>( "UseElectroosmoticVelocityBC" );
|
UseSlippingVelBC = stokes_db->getScalar<bool>( "UseElectroosmoticVelocityBC" );
|
||||||
}
|
}
|
||||||
if (electric_db->keyExists( "epsilonR" )){
|
if (stokes_db->keyExists( "epsilonR" )){
|
||||||
epsilonR = electric_db->getScalar<double>( "epsilonR" );
|
epsilonR = stokes_db->getScalar<double>( "epsilonR" );
|
||||||
}
|
}
|
||||||
|
|
||||||
// Re-calculate model parameters due to parameter read
|
// Re-calculate model parameters due to parameter read
|
||||||
@ -282,14 +282,14 @@ void ScaLBL_StokesModel::ReadInput(){
|
|||||||
if (rank == 0) cout << " Domain set." << endl;
|
if (rank == 0) cout << " Domain set." << endl;
|
||||||
}
|
}
|
||||||
|
|
||||||
void ScaLBL_IonModel::AssignZetaPotentialSolid(double *zeta_potential_solid)
|
void ScaLBL_StokesModel::AssignZetaPotentialSolid(double *zeta_potential_solid)
|
||||||
{
|
{
|
||||||
size_t NLABELS=0;
|
size_t NLABELS=0;
|
||||||
signed char VALUE=0;
|
signed char VALUE=0;
|
||||||
double AFFINITY=0.f;
|
double AFFINITY=0.f;
|
||||||
|
|
||||||
auto LabelList = ion_db->getVector<int>( "SolidLabels" );
|
auto LabelList = stokes_db->getVector<int>( "SolidLabels" );
|
||||||
auto AffinityList = ion_db->getVector<double>( "ZetaPotentialSolidList" );
|
auto AffinityList = stokes_db->getVector<double>( "ZetaPotentialSolidList" );
|
||||||
|
|
||||||
NLABELS=LabelList.size();
|
NLABELS=LabelList.size();
|
||||||
if (NLABELS != AffinityList.size()){
|
if (NLABELS != AffinityList.size()){
|
||||||
@ -335,9 +335,8 @@ void ScaLBL_IonModel::AssignZetaPotentialSolid(double *zeta_potential_solid)
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
void ScaLBL_IonModel::AssignSolidGrad(double *solid_grad)
|
void ScaLBL_StokesModel::AssignSolidGrad(double *solid_grad)
|
||||||
{
|
{
|
||||||
//TODO need to normalize the computed solid grad!!!
|
|
||||||
double *Dst;
|
double *Dst;
|
||||||
Dst = new double [3*3*3];
|
Dst = new double [3*3*3];
|
||||||
for (int kk=0; kk<3; kk++){
|
for (int kk=0; kk<3; kk++){
|
||||||
@ -424,9 +423,12 @@ void ScaLBL_IonModel::AssignSolidGrad(double *solid_grad)
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
solid_grad[idx+0*Np] = phi_x;
|
//solid_grad normalization
|
||||||
solid_grad[idx+1*Np] = phi_y;
|
double phi_mag=sqrt(phi_x*phi_x+phi_y*phi_y+phi_z*phi_z);
|
||||||
solid_grad[idx+2*Np] = phi_z;
|
if (phi_mag==0.0) phi_mag=1.0;
|
||||||
|
solid_grad[idx+0*Np] = phi_x/phi_mag;
|
||||||
|
solid_grad[idx+1*Np] = phi_y/phi_mag;
|
||||||
|
solid_grad[idx+2*Np] = phi_z/phi_mag;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -543,7 +545,7 @@ void ScaLBL_StokesModel::Run_Lite(double *ChargeDensity, double *ElectricField){
|
|||||||
|
|
||||||
if (UseSlippingVelBC==true){
|
if (UseSlippingVelBC==true){
|
||||||
ScaLBL_Comm->SolidSlippingVelocityBCD3Q19(fq, ZetaPotentialSolid, ElectricField, SolidGrad,
|
ScaLBL_Comm->SolidSlippingVelocityBCD3Q19(fq, ZetaPotentialSolid, ElectricField, SolidGrad,
|
||||||
epslion_LB, 1.0/rlx_setA, rho0, den_scale, h, time_conv);
|
epsilon_LB, 1.0/rlx_setA, rho0, den_scale, h, time_conv);
|
||||||
}
|
}
|
||||||
ScaLBL_Comm->Barrier(); comm.barrier();
|
ScaLBL_Comm->Barrier(); comm.barrier();
|
||||||
|
|
||||||
@ -570,7 +572,7 @@ void ScaLBL_StokesModel::Run_Lite(double *ChargeDensity, double *ElectricField){
|
|||||||
0, ScaLBL_Comm->LastExterior(), Np);
|
0, ScaLBL_Comm->LastExterior(), Np);
|
||||||
if (UseSlippingVelBC==true){
|
if (UseSlippingVelBC==true){
|
||||||
ScaLBL_Comm->SolidSlippingVelocityBCD3Q19(fq, ZetaPotentialSolid, ElectricField, SolidGrad,
|
ScaLBL_Comm->SolidSlippingVelocityBCD3Q19(fq, ZetaPotentialSolid, ElectricField, SolidGrad,
|
||||||
epslion_LB, 1.0/rlx_setA, rho0, den_scale, h, time_conv);
|
epsilon_LB, 1.0/rlx_setA, rho0, den_scale, h, time_conv);
|
||||||
}
|
}
|
||||||
ScaLBL_Comm->Barrier(); comm.barrier();
|
ScaLBL_Comm->Barrier(); comm.barrier();
|
||||||
//************************************************************************/
|
//************************************************************************/
|
||||||
|
@ -73,6 +73,7 @@ public:
|
|||||||
double *Velocity;
|
double *Velocity;
|
||||||
double *Pressure;
|
double *Pressure;
|
||||||
double *ZetaPotentialSolid;
|
double *ZetaPotentialSolid;
|
||||||
|
double *SolidGrad;
|
||||||
|
|
||||||
//Minkowski Morphology;
|
//Minkowski Morphology;
|
||||||
DoubleArray Velocity_x;
|
DoubleArray Velocity_x;
|
||||||
@ -91,5 +92,7 @@ private:
|
|||||||
void LoadParams(std::shared_ptr<Database> db0);
|
void LoadParams(std::shared_ptr<Database> db0);
|
||||||
void Velocity_LB_to_Phys(DoubleArray &Vel_reg);
|
void Velocity_LB_to_Phys(DoubleArray &Vel_reg);
|
||||||
vector<double> computeElectricForceAvg(double *ChargeDensity, double *ElectricField);
|
vector<double> computeElectricForceAvg(double *ChargeDensity, double *ElectricField);
|
||||||
|
void AssignSolidGrad(double *solid_grad);
|
||||||
|
void AssignZetaPotentialSolid(double *zeta_potential_solid);
|
||||||
};
|
};
|
||||||
#endif
|
#endif
|
||||||
|
Loading…
Reference in New Issue
Block a user