add voltage BC for d3q19 poisson solver
This commit is contained in:
parent
95233d61a7
commit
111529ff5e
@ -2472,6 +2472,29 @@ void ScaLBL_Communicator::D3Q7_Poisson_Potential_BC_Z(int *neighborList, double
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
void ScaLBL_Communicator::D3Q19_Poisson_Potential_BC_z(int *neighborList, double *fq, double Vin, int time){
|
||||
if (kproc == 0) {
|
||||
if (time%2==0){
|
||||
ScaLBL_D3Q19_AAeven_Poisson_Potential_BC_z(dvcSendList_z, fq, Vin, sendCount_z, N);
|
||||
}
|
||||
else{
|
||||
ScaLBL_D3Q19_AAodd_Poisson_Potential_BC_z(neighborList, dvcSendList_z, fq, Vin, sendCount_z, N);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
void ScaLBL_Communicator::D3Q19_Poisson_Potential_BC_Z(int *neighborList, double *fq, double Vout, int time){
|
||||
if (kproc == nprocz-1){
|
||||
if (time%2==0){
|
||||
ScaLBL_D3Q19_AAeven_Poisson_Potential_BC_Z(dvcSendList_Z, fq, Vout, sendCount_Z, N);
|
||||
}
|
||||
else{
|
||||
ScaLBL_D3Q19_AAodd_Poisson_Potential_BC_Z(neighborList, dvcSendList_Z, fq, Vout, sendCount_Z, N);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
void ScaLBL_Communicator::Poisson_D3Q7_BC_z(int *Map, double *Psi, double Vin){
|
||||
if (kproc == 0) {
|
||||
ScaLBL_Poisson_D3Q7_BC_z(dvcSendList_z, Map, Psi, Vin, sendCount_z);
|
||||
|
@ -396,6 +396,12 @@ extern "C" void ScaLBL_D3Q19_AAeven_Poisson(int *Map, double *dist,
|
||||
|
||||
extern "C" void ScaLBL_D3Q19_Poisson_getElectricField(double *dist, double *ElectricField, double tau, int Np);
|
||||
|
||||
|
||||
extern "C" void ScaLBL_D3Q19_AAodd_Poisson_Potential_BC_Z(int *d_neighborList, int *list, double *dist, double Vin, int count, int Np);
|
||||
extern "C" void ScaLBL_D3Q19_AAodd_Poisson_Potential_BC_z(int *d_neighborList, int *list, double *dist, double Vin, int count, int Np);
|
||||
extern "C" void ScaLBL_D3Q19_AAeven_Poisson_Potential_BC_Z(int *list, double *dist, double Vout, int count, int Np);
|
||||
extern "C" void ScaLBL_D3Q19_AAeven_Poisson_Potential_BC_z(int *list, double *dist, double Vout, int count, 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,double rho0, double den_scale, double h, double time_conv, bool UseSlippingVelBC, int start, int finish, int Np);
|
||||
@ -806,6 +812,8 @@ public:
|
||||
double D3Q19_Flux_BC_z(int *neighborList, double *fq, double flux, int time);
|
||||
void D3Q7_Poisson_Potential_BC_z(int *neighborList, double *fq, double Vin, int time);
|
||||
void D3Q7_Poisson_Potential_BC_Z(int *neighborList, double *fq, double Vout, int time);
|
||||
void D3Q19_Poisson_Potential_BC_z(int *neighborList, double *fq, double Vin, int time);
|
||||
void D3Q19_Poisson_Potential_BC_Z(int *neighborList, double *fq, double Vout, int time);
|
||||
void Poisson_D3Q7_BC_z(int *Map, double *Psi, double Vin);
|
||||
void Poisson_D3Q7_BC_Z(int *Map, double *Psi, double Vout);
|
||||
void D3Q7_Ion_Concentration_BC_z(int *neighborList, double *fq, double Cin, int time);
|
||||
|
@ -181,6 +181,7 @@ extern "C" void ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_z(int *d_neighborList,
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
extern "C" void ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_Z(int *d_neighborList,
|
||||
int *list,
|
||||
double *dist,
|
||||
@ -213,6 +214,8 @@ extern "C" void ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_Z(int *d_neighborList,
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
|
||||
extern "C" void ScaLBL_Poisson_D3Q7_BC_z(int *list, int *Map, double *Psi,
|
||||
double Vin, int count) {
|
||||
int idx, n, nm;
|
||||
|
@ -631,6 +631,8 @@ extern "C" void ScaLBL_D3Q19_AAeven_Poisson_ElectricPotential(
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
|
||||
extern "C" void ScaLBL_D3Q19_AAodd_Poisson(int *neighborList, int *Map,
|
||||
double *dist, double *Den_charge,
|
||||
double *Psi, double *ElectricField,
|
||||
@ -915,6 +917,92 @@ extern "C" void ScaLBL_D3Q19_AAeven_Poisson(int *Map, double *dist,
|
||||
}
|
||||
}
|
||||
|
||||
extern "C" void ScaLBL_D3Q19_AAeven_Poisson_Potential_BC_z(int *list, double *dist, double Vin, int count, int Np) {
|
||||
double W1 = 1.0/24.0;
|
||||
double W2 = 1.0/48.0;
|
||||
int nread, nr5;
|
||||
|
||||
for (int idx = 0; idx < count; idx++) {
|
||||
int n = list[idx];
|
||||
|
||||
dist[6 * Np + n] = W1*Vin;
|
||||
dist[12 * Np + n] = W2*Vin;
|
||||
dist[13 * Np + n] = W2*Vin;
|
||||
dist[16 * Np + n] = W2*Vin;
|
||||
dist[17 * Np + n] = W2*Vin;
|
||||
}
|
||||
}
|
||||
|
||||
extern "C" void ScaLBL_D3Q19_AAeven_Poisson_Potential_BC_Z(int *list,
|
||||
double *dist,
|
||||
double Vout,
|
||||
int count, int Np) {
|
||||
|
||||
double W1 = 1.0/24.0;
|
||||
double W2 = 1.0/48.0;
|
||||
|
||||
for (int idx = 0; idx < count; idx++) {
|
||||
|
||||
int n = list[idx];
|
||||
dist[5 * Np + n] = W1*Vout;
|
||||
dist[11 * Np + n] = W2*Vout;
|
||||
dist[14 * Np + n] = W2*Vout;
|
||||
dist[15 * Np + n] = W2*Vout;
|
||||
dist[18 * Np + n] = W2*Vout;
|
||||
}
|
||||
}
|
||||
|
||||
extern "C" void ScaLBL_D3Q19_AAodd_Poisson_Potential_BC_z(int *d_neighborList,
|
||||
int *list,
|
||||
double *dist,
|
||||
double Vin, int count,
|
||||
int Np) {
|
||||
double W1 = 1.0/24.0;
|
||||
double W2 = 1.0/48.0;
|
||||
int nr5, nr11, nr14, nr15, nr18;
|
||||
|
||||
for (int idx = 0; idx < count; idx++) {
|
||||
int n = list[idx];
|
||||
|
||||
// Unknown distributions
|
||||
nr5 = d_neighborList[n + 4 * Np];
|
||||
nr11 = d_neighborList[n + 10 * Np];
|
||||
nr15 = d_neighborList[n + 14 * Np];
|
||||
nr14 = d_neighborList[n + 13 * Np];
|
||||
nr18 = d_neighborList[n + 17 * Np];
|
||||
|
||||
dist[nr5] = W1*Vin;
|
||||
dist[nr11] = W2*Vin;
|
||||
dist[nr15] = W2*Vin;
|
||||
dist[nr14] = W2*Vin;
|
||||
dist[nr18] = W2*Vin;
|
||||
}
|
||||
}
|
||||
|
||||
extern "C" void ScaLBL_D3Q19_AAodd_Poisson_Potential_BC_Z(int *d_neighborList, int *list, double *dist, double Vout, int count, int Np) {
|
||||
|
||||
double W1 = 1.0/24.0;
|
||||
double W2 = 1.0/48.0;
|
||||
int nr6, nr12, nr13, nr16, nr17;
|
||||
|
||||
for (int idx = 0; idx < count; idx++) {
|
||||
|
||||
int n = list[idx];
|
||||
// unknown distributions
|
||||
nr6 = d_neighborList[n + 5 * Np];
|
||||
nr12 = d_neighborList[n + 11 * Np];
|
||||
nr16 = d_neighborList[n + 15 * Np];
|
||||
nr17 = d_neighborList[n + 16 * Np];
|
||||
nr13 = d_neighborList[n + 12 * Np];
|
||||
|
||||
dist[nr6] = W1*Vout;
|
||||
dist[nr12] = W2*Vout;
|
||||
dist[nr16] = W2*Vout;
|
||||
dist[nr17] = W2*Vout;
|
||||
dist[nr13] = W2*Vout;
|
||||
}
|
||||
}
|
||||
|
||||
extern "C" void ScaLBL_D3Q19_Poisson_Init(int *Map, double *dist, double *Psi,
|
||||
int start, int finish, int Np) {
|
||||
int n;
|
||||
|
@ -300,7 +300,8 @@ void ScaLBL_IonModel::ReadParams(string filename, vector<int> &num_iter) {
|
||||
void ScaLBL_IonModel::ReadParams(string filename) {
|
||||
//NOTE: the maximum iteration timesteps for ions are left unspecified
|
||||
// it relies on the multiphys controller to compute the max timestep
|
||||
USE_MEMBRANE = true;
|
||||
USE_MEMBRANE = false;
|
||||
Restart = false;
|
||||
// read the input database
|
||||
db = std::make_shared<Database>(filename);
|
||||
domain_db = db->getDatabase("Domain");
|
||||
|
@ -726,13 +726,13 @@ void ScaLBL_Poisson::Run(double *ChargeDensity, bool UseSlippingVelBC, int times
|
||||
// *************ODD TIMESTEP*************//
|
||||
timestep++;
|
||||
SolveElectricPotentialAAodd(timestep_from_Study);//update electric potential
|
||||
SolvePoissonAAodd(ChargeDensity, UseSlippingVelBC);//perform collision
|
||||
SolvePoissonAAodd(ChargeDensity, UseSlippingVelBC,timestep);//perform collision
|
||||
ScaLBL_Comm->Barrier(); comm.barrier();
|
||||
|
||||
// *************EVEN TIMESTEP*************//
|
||||
timestep++;
|
||||
SolveElectricPotentialAAeven(timestep_from_Study);//update electric potential
|
||||
SolvePoissonAAeven(ChargeDensity, UseSlippingVelBC);//perform collision
|
||||
SolvePoissonAAeven(ChargeDensity, UseSlippingVelBC,timestep);//perform collision
|
||||
ScaLBL_Comm->Barrier(); comm.barrier();
|
||||
//************************************************************************/
|
||||
|
||||
@ -797,14 +797,14 @@ void ScaLBL_Poisson::Run(double *ChargeDensity, bool UseSlippingVelBC, int times
|
||||
//************************************************************************/
|
||||
// *************ODD TIMESTEP*************//
|
||||
timestep++;
|
||||
//SolveElectricPotentialAAodd(timestep_from_Study,ChargeDensity, UseSlippingVelBC);//update electric potential
|
||||
SolvePoissonAAodd(ChargeDensity, UseSlippingVelBC);//perform collision
|
||||
//SolveElectricPotentialAAodd(timestep_from_Study);//,ChargeDensity, UseSlippingVelBC);//update electric potential
|
||||
SolvePoissonAAodd(ChargeDensity, UseSlippingVelBC,timestep);//perform collision
|
||||
ScaLBL_Comm->Barrier(); comm.barrier();
|
||||
|
||||
// *************EVEN TIMESTEP*************//
|
||||
timestep++;
|
||||
//SolveElectricPotentialAAeven(timestep_from_Study,ChargeDensity, UseSlippingVelBC);//update electric potential
|
||||
SolvePoissonAAeven(ChargeDensity, UseSlippingVelBC);//perform collision
|
||||
//SolveElectricPotentialAAeven(timestep_from_Study);//,ChargeDensity, UseSlippingVelBC);//update electric potential
|
||||
SolvePoissonAAeven(ChargeDensity, UseSlippingVelBC,timestep);//perform collision
|
||||
ScaLBL_Comm->Barrier(); comm.barrier();
|
||||
//************************************************************************/
|
||||
|
||||
@ -884,13 +884,13 @@ void ScaLBL_Poisson::Run(double *ChargeDensity, DoubleArray MembraneDistance, bo
|
||||
// *************ODD TIMESTEP*************//
|
||||
timestep++;
|
||||
//SolveElectricPotentialAAodd(timestep_from_Study,ChargeDensity, UseSlippingVelBC);//update electric potential
|
||||
SolvePoissonAAodd(ChargeDensity, UseSlippingVelBC);//perform collision
|
||||
SolvePoissonAAodd(ChargeDensity, UseSlippingVelBC,timestep);//perform collision
|
||||
ScaLBL_Comm->Barrier(); comm.barrier();
|
||||
|
||||
// *************EVEN TIMESTEP*************//
|
||||
timestep++;
|
||||
//SolveElectricPotentialAAeven(timestep_from_Study,ChargeDensity, UseSlippingVelBC);//update electric potential
|
||||
SolvePoissonAAeven(ChargeDensity, UseSlippingVelBC);//perform collision
|
||||
SolvePoissonAAeven(ChargeDensity, UseSlippingVelBC,timestep);//perform collision
|
||||
ScaLBL_Comm->Barrier(); comm.barrier();
|
||||
//************************************************************************/
|
||||
|
||||
@ -913,7 +913,7 @@ void ScaLBL_Poisson::Run(double *ChargeDensity, DoubleArray MembraneDistance, bo
|
||||
}
|
||||
error=Dm->Comm.maxReduce(max_error);
|
||||
|
||||
if (error > tolerance & SET_THRESHOLD){
|
||||
if (error > tolerance && SET_THRESHOLD){
|
||||
/* don't use this with an external BC */
|
||||
// cpompute the far-field electric potential
|
||||
double inside_local = 0.0;
|
||||
@ -1146,7 +1146,7 @@ void ScaLBL_Poisson::SolveElectricPotentialAAeven(int timestep_from_Study){
|
||||
}
|
||||
}
|
||||
|
||||
void ScaLBL_Poisson::SolvePoissonAAodd(double *ChargeDensity, bool UseSlippingVelBC){
|
||||
void ScaLBL_Poisson::SolvePoissonAAodd(double *ChargeDensity, bool UseSlippingVelBC, int timestep){
|
||||
|
||||
if (lattice_scheme.compare("D3Q7")==0){
|
||||
ScaLBL_D3Q7_AAodd_Poisson(NeighborList, dvcMap, fq, ChargeDensity, Psi, ElectricField, tau, epsilon_LB, UseSlippingVelBC, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
||||
@ -1167,6 +1167,30 @@ void ScaLBL_Poisson::SolvePoissonAAodd(double *ChargeDensity, bool UseSlippingVe
|
||||
ScaLBL_D3Q19_AAodd_Poisson(NeighborList, dvcMap, fq, ChargeDensity, Psi, ElectricField, tau, epsilon_LB, UseSlippingVelBC, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
||||
//ScaLBL_Comm->RecvD3Q19AA(fq); //WRITE INTO OPPOSITE
|
||||
ScaLBL_Comm->RecvD3Q19AA(fq); //WRITE INTO OPPOSITE
|
||||
|
||||
// Set boundary conditions
|
||||
if (BoundaryConditionInlet > 0){
|
||||
switch (BoundaryConditionInlet){
|
||||
case 1:
|
||||
ScaLBL_Comm->D3Q7_Poisson_Potential_BC_z(NeighborList, fq, Vin, timestep);
|
||||
break;
|
||||
case 2:
|
||||
Vin = getBoundaryVoltagefromPeriodicBC(Vin0,freqIn,PhaseShift_In,timestep);
|
||||
ScaLBL_Comm->D3Q7_Poisson_Potential_BC_z(NeighborList, fq, Vin, timestep);
|
||||
break;
|
||||
}
|
||||
}
|
||||
if (BoundaryConditionOutlet > 0){
|
||||
switch (BoundaryConditionOutlet){
|
||||
case 1:
|
||||
ScaLBL_Comm->D3Q7_Poisson_Potential_BC_Z(NeighborList, fq, Vout, timestep);
|
||||
break;
|
||||
case 2:
|
||||
Vout = getBoundaryVoltagefromPeriodicBC(Vout0,freqOut,PhaseShift_Out,timestep);
|
||||
ScaLBL_Comm->D3Q7_Poisson_Potential_BC_Z(NeighborList, fq, Vout, timestep);
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
ScaLBL_D3Q19_AAodd_Poisson(NeighborList, dvcMap, fq, ChargeDensity, Psi, ElectricField, tau, epsilon_LB, UseSlippingVelBC, 0, ScaLBL_Comm->LastExterior(), Np);
|
||||
ScaLBL_Comm->Barrier();
|
||||
@ -1176,7 +1200,7 @@ void ScaLBL_Poisson::SolvePoissonAAodd(double *ChargeDensity, bool UseSlippingVe
|
||||
}
|
||||
}
|
||||
|
||||
void ScaLBL_Poisson::SolvePoissonAAeven(double *ChargeDensity, bool UseSlippingVelBC){
|
||||
void ScaLBL_Poisson::SolvePoissonAAeven(double *ChargeDensity, bool UseSlippingVelBC, int timestep){
|
||||
|
||||
if (lattice_scheme.compare("D3Q7")==0){
|
||||
ScaLBL_D3Q7_AAeven_Poisson(dvcMap, fq, ChargeDensity, Psi, ElectricField, tau, epsilon_LB, UseSlippingVelBC, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
||||
@ -1194,6 +1218,31 @@ void ScaLBL_Poisson::SolvePoissonAAeven(double *ChargeDensity, bool UseSlippingV
|
||||
ScaLBL_D3Q19_AAeven_Poisson(dvcMap, fq, ChargeDensity, Psi, ElectricField, ResidualError, tau, epsilon_LB, UseSlippingVelBC, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
|
||||
ScaLBL_Comm->RecvD3Q19AA(fq); //WRITE INTO OPPOSITE
|
||||
// ScaLBL_Comm->RecvD3Q19AA(fq); //WRITE INTO OPPOSITE
|
||||
|
||||
// Set boundary conditions
|
||||
if (BoundaryConditionInlet > 0){
|
||||
switch (BoundaryConditionInlet){
|
||||
case 1:
|
||||
ScaLBL_Comm->D3Q7_Poisson_Potential_BC_z(NeighborList, fq, Vin, timestep);
|
||||
break;
|
||||
case 2:
|
||||
Vin = getBoundaryVoltagefromPeriodicBC(Vin0,freqIn,PhaseShift_In,timestep);
|
||||
ScaLBL_Comm->D3Q7_Poisson_Potential_BC_z(NeighborList, fq, Vin, timestep);
|
||||
break;
|
||||
}
|
||||
}
|
||||
if (BoundaryConditionOutlet > 0){
|
||||
switch (BoundaryConditionOutlet){
|
||||
case 1:
|
||||
ScaLBL_Comm->D3Q7_Poisson_Potential_BC_Z(NeighborList, fq, Vout, timestep);
|
||||
break;
|
||||
case 2:
|
||||
Vout = getBoundaryVoltagefromPeriodicBC(Vout0,freqOut,PhaseShift_Out,timestep);
|
||||
ScaLBL_Comm->D3Q7_Poisson_Potential_BC_Z(NeighborList, fq, Vout, timestep);
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
ScaLBL_D3Q19_AAeven_Poisson(dvcMap, fq, ChargeDensity, Psi, ElectricField, ResidualError, tau, epsilon_LB, UseSlippingVelBC, 0, ScaLBL_Comm->LastExterior(), Np);
|
||||
ScaLBL_Comm->Barrier();
|
||||
|
||||
@ -1323,7 +1372,6 @@ void ScaLBL_Poisson::ElectricField_LB_to_Phys(DoubleArray &Efield_reg){
|
||||
void ScaLBL_Poisson::WriteVis( int timestep) {
|
||||
|
||||
auto vis_db = db->getDatabase("Visualization");
|
||||
char VisName[40];
|
||||
auto format = vis_db->getWithDefault<string>( "format", "hdf5" );
|
||||
|
||||
DoubleArray ElectricalPotential(Nx, Ny, Nz);
|
||||
|
@ -118,8 +118,8 @@ private:
|
||||
void SolveElectricPotentialAAeven(int timestep_from_Study);
|
||||
void SolveElectricPotentialAAodd(int timestep_from_Study);
|
||||
//void SolveElectricField();
|
||||
void SolvePoissonAAodd(double *ChargeDensity, bool UseSlippingVelBC);
|
||||
void SolvePoissonAAeven(double *ChargeDensity, bool UseSlippingVelBC);
|
||||
void SolvePoissonAAodd(double *ChargeDensity, bool UseSlippingVelBC, int timestep);
|
||||
void SolvePoissonAAeven(double *ChargeDensity, bool UseSlippingVelBC, int timestep);
|
||||
void getConvergenceLog(int timestep,double error);
|
||||
double getBoundaryVoltagefromPeriodicBC(double V0,double freq,double t0,int time_step);
|
||||
|
||||
|
@ -104,7 +104,6 @@ int main(int argc, char **argv)
|
||||
//StokesModel.Run_Lite(IonModel.ChargeDensity, PoissonSolver.ElectricField);// Solve the N-S equations to get velocity
|
||||
IonModel.Run(IonModel.FluidVelocityDummy,PoissonSolver.ElectricField); //solve for ion transport and electric potential
|
||||
|
||||
|
||||
//if (timestep%Study.analysis_interval==0){
|
||||
// Analysis.Basic(IonModel,PoissonSolver,StokesModel,timestep);
|
||||
//}
|
||||
|
Loading…
Reference in New Issue
Block a user