Files
LBPM/models/IonModel.cpp

873 lines
34 KiB
C++
Raw Normal View History

2020-08-06 15:41:40 -04:00
/*
* Dilute Ion Transport LBM Model
2020-08-06 15:41:40 -04:00
*/
#include <algorithm>
2020-08-06 16:06:52 -04:00
#include "models/IonModel.h"
2020-08-06 15:41:40 -04:00
#include "analysis/distance.h"
#include "common/ReadMicroCT.h"
ScaLBL_IonModel::ScaLBL_IonModel(int RANK, int NP, MPI_Comm COMM):
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),
tolerance(0),number_ion_species(0),Nx(0),Ny(0),Nz(0),N(0),Np(0),nprocx(0),nprocy(0),nprocz(0),
fluidVelx_dummy(0),fluidVely_dummy(0),fluidVelz_dummy(0),
BoundaryCondition(0),BoundaryConditionSolid(0),Lx(0),Ly(0),Lz(0),comm(COMM)
2020-08-06 15:41:40 -04:00
{
}
ScaLBL_IonModel::~ScaLBL_IonModel(){
}
void ScaLBL_IonModel::ReadParams(string filename,vector<int> &num_iter){
2020-08-10 12:03:28 -04:00
// read the input database
2020-08-06 15:41:40 -04:00
db = std::make_shared<Database>( filename );
domain_db = db->getDatabase( "Domain" );
ion_db = db->getDatabase( "Ions" );
2020-08-14 14:23:22 -04:00
// Universal constant
kb = 1.38e-23;//Boltzmann constant;unit [J/K]
electron_charge = 1.6e-19;//electron charge;unit [C]
//---------------------- Default model parameters --------------------------//
2020-08-10 12:03:28 -04:00
T = 300.0;//temperature; unit [K]
2020-08-14 14:23:22 -04:00
Vt = kb*T/electron_charge;//thermal voltage; unit [Vy]
k2_inv = 4.0;//speed of sound for D3Q7 lattice
2020-08-10 12:03:28 -04:00
h = 1.0;//resolution; unit: um/lu
2020-08-06 15:41:40 -04:00
tolerance = 1.0e-8;
number_ion_species = 1;
tau.push_back(1.0);
2020-08-14 14:23:22 -04:00
IonDiffusivity.push_back(1.0e-9);//user-input diffusivity has physical unit [m^2/sec]
IonValence.push_back(1);//algebraic valence charge
IonConcentration.push_back(1.0e-3);//user-input ion concentration has physical unit [mol/m^3]
Cin.push_back(1.0e-3);//user-input inlet boundary ion concentration;unit [mol/m^3]
Cout.push_back(1.0e-3);//user-input outlet boundary ion concentration;unit [mol/m^3]
//tau.push_back(0.5+k2_inv*time_conv/(h*1.0e-6)/(h*1.0e-6)*IonDiffusivity[0]);
time_conv.push_back((tau[0]-0.5)/k2_inv*(h*h*1.0e-12)/IonDiffusivity[0]);
fluidVelx_dummy = 0.0;//for debugging, unit [m/sec]
fluidVely_dummy = 0.0;//for debugging, unit [m/sec]
fluidVelz_dummy = 0.0;//for debugging, unit [m/sec]
Ex_dummy = 0.0;//for debugging, unit [V/m]
Ey_dummy = 0.0;//for debugging, unit [V/m]
Ez_dummy = 0.0;//for debugging, unit [V/m]
2020-08-14 14:23:22 -04:00
//--------------------------------------------------------------------------//
2020-08-10 12:03:28 -04:00
2020-08-20 22:47:10 -04:00
// Read domain parameters
if (domain_db->keyExists( "voxel_length" )){//default unit: um/lu
h = domain_db->getScalar<double>( "voxel_length" );
}
2020-08-10 12:03:28 -04:00
// LB-Ion Model parameters
2020-08-14 14:23:22 -04:00
//if (ion_db->keyExists( "timestepMax" )){
// timestepMax = ion_db->getScalar<int>( "timestepMax" );
//}
2020-08-10 12:03:28 -04:00
if (ion_db->keyExists( "tolerance" )){
tolerance = ion_db->getScalar<double>( "tolerance" );
}
if (ion_db->keyExists( "temperature" )){
T = ion_db->getScalar<int>( "temperature" );
//re-calculate thermal voltage
Vt = kb*T/electron_charge;//thermal voltage; unit [Vy]
2020-08-10 12:03:28 -04:00
}
if (ion_db->keyExists( "FluidVelDummy" )){
fluidVelx_dummy = ion_db->getVector<double>( "FluidVelDummy" )[0];
fluidVely_dummy = ion_db->getVector<double>( "FluidVelDummy" )[1];
fluidVelz_dummy = ion_db->getVector<double>( "FluidVelDummy" )[2];
}
if (ion_db->keyExists( "ElectricFieldDummy" )){
Ex_dummy = ion_db->getVector<double>( "ElectricFieldDummy" )[0];
Ey_dummy = ion_db->getVector<double>( "ElectricFieldDummy" )[1];
Ez_dummy = ion_db->getVector<double>( "ElectricFieldDummy" )[2];
}
if (ion_db->keyExists( "number_ion_species" )){
number_ion_species = ion_db->getScalar<int>( "number_ion_species" );
2020-08-06 15:41:40 -04:00
}
//------ Load number of iteration from multiphysics controller ------//
if (num_iter.size()!=number_ion_species){
ERROR("Error: number_ion_species and num_iter_Ion_List (from Multiphysics) must be of the same length! \n");
}
else{
timestepMax.assign(num_iter.begin(),num_iter.end());
}
//-------------------------------------------------------------------//
if (ion_db->keyExists("tauList")){
tau.clear();
tau = ion_db->getVector<double>( "tauList" );
vector<double>Di = ion_db->getVector<double>( "IonDiffusivityList" );//temp storing ion diffusivity in physical unit
if (tau.size()!=number_ion_species || Di.size()!=number_ion_species){
ERROR("Error: number_ion_species, tauList and IonDiffusivityList must be of the same length! \n");
}
else{
time_conv.clear();
for (int i=0; i<tau.size();i++){
time_conv.push_back((tau[i]-0.5)/k2_inv*(h*h*1.0e-12)/Di[i]);
}
}
}
2020-08-10 12:03:28 -04:00
//read ion related list
2020-08-14 14:23:22 -04:00
//NOTE: ion diffusivity has INPUT unit: [m^2/sec]
// it must be converted to LB unit: [lu^2/lt]
2020-08-10 12:03:28 -04:00
if (ion_db->keyExists("IonDiffusivityList")){
IonDiffusivity.clear();
IonDiffusivity = ion_db->getVector<double>( "IonDiffusivityList" );
if (IonDiffusivity.size()!=number_ion_species){
ERROR("Error: number_ion_species and IonDiffusivityList must be the same length! \n");
}
2020-08-14 14:23:22 -04:00
else{
for (int i=0; i<IonDiffusivity.size();i++){
IonDiffusivity[i] = IonDiffusivity[i]*time_conv[i]/(h*h*1.0e-12);//LB diffusivity has unit [lu^2/lt]
2020-08-14 14:23:22 -04:00
}
}
2020-08-10 12:03:28 -04:00
}
else {
for (int i=0; i<IonDiffusivity.size();i++){
//convert ion diffusivity in physical unit to LB unit
IonDiffusivity[i] = IonDiffusivity[i]*time_conv[i]/(h*h*1.0e-12);//LB diffusivity has unit [lu^2/lt]
}
}
// read time relaxation time list
2020-08-10 12:03:28 -04:00
//read ion algebric valence list
if (ion_db->keyExists("IonValenceList")){
IonValence.clear();
IonValence = ion_db->getVector<int>( "IonValenceList" );
if (IonValence.size()!=number_ion_species){
ERROR("Error: number_ion_species and IonValenceList must be the same length! \n");
}
}
2020-08-14 14:23:22 -04:00
//read initial ion concentration list; INPUT unit [mol/m^3]
//it must be converted to LB unit [mol/lu^3]
2020-08-10 12:03:28 -04:00
if (ion_db->keyExists("IonConcentrationList")){
IonConcentration.clear();
IonConcentration = ion_db->getVector<double>( "IonConcentrationList" );
if (IonConcentration.size()!=number_ion_species){
ERROR("Error: number_ion_species and IonConcentrationList must be the same length! \n");
}
2020-08-14 14:23:22 -04:00
else{
for (int i=0; i<IonConcentration.size();i++){
2020-08-14 14:23:22 -04:00
IonConcentration[i] = IonConcentration[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
2020-08-10 12:03:28 -04:00
}
else {
for (int i=0; i<IonConcentration.size();i++){
IonConcentration[i] = IonConcentration[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
2020-08-10 12:03:28 -04:00
//Read solid boundary condition specific to Ion model
BoundaryConditionSolid = 0;
if (ion_db->keyExists( "BC_Solid" )){
BoundaryConditionSolid = ion_db->getScalar<int>( "BC_Solid" );
}
// Read boundary condition for ion transport
// BC = 0: normal periodic BC
// BC = 1: fixed inlet and outlet ion concentration
BoundaryCondition = 0;
if (ion_db->keyExists( "BC" )){
BoundaryCondition = ion_db->getScalar<int>( "BC" );
}
if (BoundaryCondition==1){
//read boundary ion concentration list; INPUT unit [mol/m^3]
//it must be converted to LB unit [mol/lu^3]
//inlet
if (ion_db->keyExists("CinList")){
Cin.clear();
Cin = ion_db->getVector<double>( "CinList" );
if (Cin.size()!=number_ion_species){
ERROR("Error: number_ion_species and CinList must be the same length! \n");
}
else{
for (int i=0; i<Cin.size();i++){
Cin[i] = Cin[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
}
else {
for (int i=0; i<Cin.size();i++){
Cin[i] = Cin[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
//outlet
if (ion_db->keyExists("CoutList")){
Cout.clear();
Cout = ion_db->getVector<double>( "CoutList" );
if (Cout.size()!=number_ion_species){
ERROR("Error: number_ion_species and CoutList must be the same length! \n");
}
else{
for (int i=0; i<Cout.size();i++){
Cout[i] = Cout[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
}
else {
for (int i=0; i<Cout.size();i++){
Cout[i] = Cout[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
}
}
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
// read the input database
db = std::make_shared<Database>( filename );
domain_db = db->getDatabase( "Domain" );
ion_db = db->getDatabase( "Ions" );
2020-08-14 14:23:22 -04:00
// Universal constant
kb = 1.38e-23;//Boltzmann constant;unit [J/K]
electron_charge = 1.6e-19;//electron charge;unit [C]
//---------------------- Default model parameters --------------------------//
T = 300.0;//temperature; unit [K]
Vt = kb*T/electron_charge;//thermal voltage; unit [Vy]
k2_inv = 4.0;//speed of sound for D3Q7 lattice
h = 1.0;//resolution; unit: um/lu
tolerance = 1.0e-8;
number_ion_species = 1;
tau.push_back(1.0);
IonDiffusivity.push_back(1.0e-9);//user-input diffusivity has physical unit [m^2/sec]
IonValence.push_back(1);//algebraic valence charge
IonConcentration.push_back(1.0e-3);//user-input ion concentration has physical unit [mol/m^3]
Cin.push_back(1.0e-3);//user-input inlet boundary ion concentration;unit [mol/m^3]
Cout.push_back(1.0e-3);//user-input outlet boundary ion concentration;unit [mol/m^3]
//tau.push_back(0.5+k2_inv*time_conv/(h*1.0e-6)/(h*1.0e-6)*IonDiffusivity[0]);
time_conv.push_back((tau[0]-0.5)/k2_inv*(h*h*1.0e-12)/IonDiffusivity[0]);
fluidVelx_dummy = 0.0;//for debugging, unit [m/sec]
fluidVely_dummy = 0.0;//for debugging, unit [m/sec]
fluidVelz_dummy = 0.0;//for debugging, unit [m/sec]
Ex_dummy = 0.0;//for debugging, unit [V/m]
Ey_dummy = 0.0;//for debugging, unit [V/m]
Ez_dummy = 0.0;//for debugging, unit [V/m]
//--------------------------------------------------------------------------//
// Read domain parameters
if (domain_db->keyExists( "voxel_length" )){//default unit: um/lu
h = domain_db->getScalar<double>( "voxel_length" );
}
// LB-Ion Model parameters
//if (ion_db->keyExists( "timestepMax" )){
// timestepMax = ion_db->getScalar<int>( "timestepMax" );
//}
if (ion_db->keyExists( "tolerance" )){
tolerance = ion_db->getScalar<double>( "tolerance" );
}
if (ion_db->keyExists( "temperature" )){
T = ion_db->getScalar<int>( "temperature" );
//re-calculate thermal voltage
Vt = kb*T/electron_charge;//thermal voltage; unit [Vy]
}
if (ion_db->keyExists( "FluidVelDummy" )){
fluidVelx_dummy = ion_db->getVector<double>( "FluidVelDummy" )[0];
fluidVely_dummy = ion_db->getVector<double>( "FluidVelDummy" )[1];
fluidVelz_dummy = ion_db->getVector<double>( "FluidVelDummy" )[2];
}
if (ion_db->keyExists( "ElectricFieldDummy" )){
Ex_dummy = ion_db->getVector<double>( "ElectricFieldDummy" )[0];
Ey_dummy = ion_db->getVector<double>( "ElectricFieldDummy" )[1];
Ez_dummy = ion_db->getVector<double>( "ElectricFieldDummy" )[2];
}
if (ion_db->keyExists( "number_ion_species" )){
number_ion_species = ion_db->getScalar<int>( "number_ion_species" );
}
if (ion_db->keyExists("tauList")){
tau.clear();
tau = ion_db->getVector<double>( "tauList" );
vector<double>Di = ion_db->getVector<double>( "IonDiffusivityList" );//temp storing ion diffusivity in physical unit
if (tau.size()!=number_ion_species || Di.size()!=number_ion_species){
ERROR("Error: number_ion_species, tauList and IonDiffusivityList must be of the same length! \n");
}
else{
time_conv.clear();
for (int i=0; i<tau.size();i++){
time_conv.push_back((tau[i]-0.5)/k2_inv*(h*h*1.0e-12)/Di[i]);
}
}
}
//read ion related list
//NOTE: ion diffusivity has INPUT unit: [m^2/sec]
// it must be converted to LB unit: [lu^2/lt]
if (ion_db->keyExists("IonDiffusivityList")){
IonDiffusivity.clear();
IonDiffusivity = ion_db->getVector<double>( "IonDiffusivityList" );
if (IonDiffusivity.size()!=number_ion_species){
ERROR("Error: number_ion_species and IonDiffusivityList must be the same length! \n");
}
else{
for (int i=0; i<IonDiffusivity.size();i++){
IonDiffusivity[i] = IonDiffusivity[i]*time_conv[i]/(h*h*1.0e-12);//LB diffusivity has unit [lu^2/lt]
}
}
}
else {
for (int i=0; i<IonDiffusivity.size();i++){
//convert ion diffusivity in physical unit to LB unit
IonDiffusivity[i] = IonDiffusivity[i]*time_conv[i]/(h*h*1.0e-12);//LB diffusivity has unit [lu^2/lt]
}
}
// read time relaxation time list
//read ion algebric valence list
if (ion_db->keyExists("IonValenceList")){
IonValence.clear();
IonValence = ion_db->getVector<int>( "IonValenceList" );
if (IonValence.size()!=number_ion_species){
ERROR("Error: number_ion_species and IonValenceList must be the same length! \n");
}
}
//read initial ion concentration list; INPUT unit [mol/m^3]
//it must be converted to LB unit [mol/lu^3]
if (ion_db->keyExists("IonConcentrationList")){
IonConcentration.clear();
IonConcentration = ion_db->getVector<double>( "IonConcentrationList" );
if (IonConcentration.size()!=number_ion_species){
ERROR("Error: number_ion_species and IonConcentrationList must be the same length! \n");
}
else{
for (int i=0; i<IonConcentration.size();i++){
IonConcentration[i] = IonConcentration[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
}
else {
for (int i=0; i<IonConcentration.size();i++){
IonConcentration[i] = IonConcentration[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
2020-08-10 12:03:28 -04:00
}
//Read solid boundary condition specific to Ion model
BoundaryConditionSolid = 0;
if (ion_db->keyExists( "BC_Solid" )){
BoundaryConditionSolid = ion_db->getScalar<int>( "BC_Solid" );
}
// Read boundary condition for ion transport
// BC = 0: normal periodic BC
// BC = 1: fixed inlet and outlet ion concentration
BoundaryCondition = 0;
if (ion_db->keyExists( "BC" )){
BoundaryCondition = ion_db->getScalar<int>( "BC" );
}
if (BoundaryCondition==1){
//read boundary ion concentration list; INPUT unit [mol/m^3]
//it must be converted to LB unit [mol/lu^3]
//inlet
if (ion_db->keyExists("CinList")){
Cin.clear();
Cin = ion_db->getVector<double>( "CinList" );
if (Cin.size()!=number_ion_species){
ERROR("Error: number_ion_species and CinList must be the same length! \n");
}
else{
for (int i=0; i<Cin.size();i++){
Cin[i] = Cin[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
}
else {
for (int i=0; i<Cin.size();i++){
Cin[i] = Cin[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
//outlet
if (ion_db->keyExists("CoutList")){
Cout.clear();
Cout = ion_db->getVector<double>( "CoutList" );
if (Cout.size()!=number_ion_species){
ERROR("Error: number_ion_species and CoutList must be the same length! \n");
}
else{
for (int i=0; i<Cout.size();i++){
Cout[i] = Cout[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
}
else {
for (int i=0; i<Cout.size();i++){
Cout[i] = Cout[i]*(h*h*h*1.0e-18);//LB ion concentration has unit [mol/lu^3]
}
}
}
2020-08-06 15:41:40 -04:00
}
2020-08-10 12:03:28 -04:00
2020-08-06 15:41:40 -04:00
void ScaLBL_IonModel::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;
Distance.resize(Nx,Ny,Nz);
for (int i=0; i<Nx*Ny*Nz; i++) Dm->id[i] = 1; // initialize this way
//Averages = std::shared_ptr<TwoPhase> ( new TwoPhase(Dm) ); // TwoPhase analysis object
MPI_Barrier(comm);
Dm->CommInit();
MPI_Barrier(comm);
rank = Dm->rank();
nprocx = Dm->nprocx();
nprocy = Dm->nprocy();
nprocz = Dm->nprocz();
}
void ScaLBL_IonModel::ReadInput(){
sprintf(LocalRankString,"%05d",Dm->rank());
sprintf(LocalRankFilename,"%s%s","ID.",LocalRankString);
sprintf(LocalRestartFile,"%s%s","Restart.",LocalRankString);
if (domain_db->keyExists( "Filename" )){
auto Filename = domain_db->getScalar<std::string>( "Filename" );
Mask->Decomp(Filename);
}
else if (domain_db->keyExists( "GridFile" )){
// Read the local domain data
auto input_id = readMicroCT( *domain_db, comm );
// 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( comm, Mask->rank_info, size0, { 1, 1, 1 }, 0, 1 );
Array<signed char> id_view;
id_view.viewRaw( size1, Mask->id );
fill.copy( input_id, id_view );
fill.fill( id_view );
}
else{
Mask->ReadIDs();
}
// 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
if (Mask->id[n] > 0) id_solid(i,j,k) = 1;
else id_solid(i,j,k) = 0;
}
}
}
// 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
Distance(i,j,k) = 2.0*double(id_solid(i,j,k))-1.0;
}
}
}
// MeanFilter(Averages->SDs);
2020-08-14 14:23:22 -04:00
if (rank==0) printf("LB Ion Solver: Initialized solid phase & converting to Signed Distance function \n");
2020-08-06 15:41:40 -04:00
CalcDist(Distance,id_solid,*Dm);
2020-08-14 14:23:22 -04:00
if (rank == 0) cout << " Domain set." << endl;
2020-08-06 15:41:40 -04:00
}
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*h*1.0e-18);
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){
2020-08-20 22:47:10 -04:00
printf("LB Ion Solver: number of 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);
}
}
}
2020-08-06 15:41:40 -04:00
void ScaLBL_IonModel::Create(){
/*
* This function creates the variables needed to run a LBM
*/
int rank=Mask->rank();
//.........................................................
// 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();
//...........................................................................
2020-08-14 14:23:22 -04:00
if (rank==0) printf ("LB Ion Solver: Create ScaLBL_Communicator \n");
2020-08-06 15:41:40 -04:00
// 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));
int Npad=(Np/16 + 2)*16;
2020-08-14 14:23:22 -04:00
if (rank==0) printf ("LB Ion Solver: Set up memory efficient layout \n");
2020-08-06 15:41:40 -04:00
Map.resize(Nx,Ny,Nz); Map.fill(-2);
auto neighborList= new int[18*Npad];
Np = ScaLBL_Comm->MemoryOptimizedLayoutAA(Map,neighborList,Mask->id,Np);
MPI_Barrier(comm);
//...........................................................................
// MAIN VARIABLES ALLOCATED HERE
//...........................................................................
// LBM variables
2020-08-14 14:23:22 -04:00
if (rank==0) printf ("LB Ion Solver: Allocating distributions \n");
2020-08-06 15:41:40 -04:00
//......................device distributions.................................
int dist_mem_size = Np*sizeof(double);
int neighborSize=18*(Np*sizeof(int));
//...........................................................................
ScaLBL_AllocateDeviceMemory((void **) &NeighborList, neighborSize);
ScaLBL_AllocateDeviceMemory((void **) &fq, number_ion_species*7*dist_mem_size);
ScaLBL_AllocateDeviceMemory((void **) &Ci, number_ion_species*sizeof(double)*Np);
ScaLBL_AllocateDeviceMemory((void **) &ChargeDensity, sizeof(double)*Np);
//...........................................................................
// Update GPU data structures
2020-08-14 14:23:22 -04:00
if (rank==0) printf ("LB Ion Solver: Setting up device map and neighbor list \n");
2020-08-06 15:41:40 -04:00
// copy the neighbor list
ScaLBL_CopyToDevice(NeighborList, neighborList, neighborSize);
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;
}
2020-08-06 15:41:40 -04:00
}
void ScaLBL_IonModel::Initialize(){
/*
* This function initializes model
*/
2020-08-14 14:23:22 -04:00
if (rank==0) printf ("LB Ion Solver: initializing D3Q7 distributions\n");
2020-08-10 12:03:28 -04:00
for (int ic=0; ic<number_ion_species; ic++){
ScaLBL_D3Q7_Ion_Init(&fq[ic*Np*7],&Ci[ic*Np],IonConcentration[ic],Np);
}
2020-08-14 14:23:22 -04:00
if (rank==0) printf ("LB Ion Solver: initializing charge density\n");
for (int ic=0; ic<number_ion_species; ic++){
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);
}
if (rank==0) printf("*****************************************************\n");
if (rank==0) printf("LB Ion Transport Solver: \n");
for (int i=0; i<number_ion_species;i++){
if (rank==0) printf(" Ion %i: LB relaxation tau = %.5g\n", i+1,tau[i]);
if (rank==0) printf(" Time conversion factor: %.5g [sec/lt]\n", time_conv[i]);
if (rank==0) printf(" Internal iteration: %i [lt]\n", timestepMax[i]);
}
if (rank==0) printf("*****************************************************\n");
switch (BoundaryConditionSolid){
case 0:
if (rank==0) printf("LB Ion Solver: solid boundary: non-flux boundary is assigned\n");
break;
case 1:
if (rank==0) printf("LB Ion Solver: solid boundary: Dirichlet-type surfacen ion concentration is assigned\n");
break;
default:
if (rank==0) printf("LB Ion Solver: solid boundary: non-flux boundary is assigned\n");
break;
}
2020-08-06 15:41:40 -04:00
}
2020-08-10 12:03:28 -04:00
void ScaLBL_IonModel::Run(double *Velocity, double *ElectricField){
//Input parameter:
//1. Velocity is from StokesModel
//2. ElectricField is from Poisson model
2020-08-10 12:03:28 -04:00
//LB-related parameter
vector<double> rlx(tau.begin(),tau.end());
for (double item : rlx){
item = 1.0/item;
}
2020-08-06 15:41:40 -04:00
//.......create and start timer............
2020-08-14 14:23:22 -04:00
//double starttime,stoptime,cputime;
//ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
//starttime = MPI_Wtime();
2020-08-10 12:03:28 -04:00
for (int ic=0; ic<number_ion_species; ic++){
timestep=0;
while (timestep < timestepMax[ic]) {
//************************************************************************/
// *************ODD TIMESTEP*************//
timestep++;
//Update ion concentration and charge density
ScaLBL_Comm->SendD3Q7AA(fq, ic); //READ FROM NORMAL
2020-08-10 12:03:28 -04:00
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_DeviceBarrier();
// Set boundary conditions
if (BoundaryCondition == 1){
ScaLBL_Comm->D3Q7_Ion_Concentration_BC_z(NeighborList, &fq[ic*Np*7], Cin[ic], timestep);
ScaLBL_Comm->D3Q7_Ion_Concentration_BC_Z(NeighborList, &fq[ic*Np*7], Cout[ic], timestep);
}
//-------------------------//
2020-08-10 12:03:28 -04:00
ScaLBL_D3Q7_AAodd_IonConcentration(NeighborList, &fq[ic*Np*7],&Ci[ic*Np], 0, ScaLBL_Comm->LastExterior(), Np);
2020-08-10 12:03:28 -04:00
//LB-Ion collison
2020-08-10 12:03:28 -04:00
ScaLBL_D3Q7_AAodd_Ion(NeighborList, &fq[ic*Np*7],&Ci[ic*Np],Velocity,ElectricField,IonDiffusivity[ic],IonValence[ic],
2020-08-14 14:23:22 -04:00
rlx[ic],Vt,ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
2020-08-10 12:03:28 -04:00
ScaLBL_D3Q7_AAodd_Ion(NeighborList, &fq[ic*Np*7],&Ci[ic*Np],Velocity,ElectricField,IonDiffusivity[ic],IonValence[ic],
2020-08-14 14:23:22 -04:00
rlx[ic],Vt,0, ScaLBL_Comm->LastExterior(), Np);
if (BoundaryConditionSolid==1){
//TODO IonSolid may also be species-dependent
ScaLBL_Comm->SolidDirichletD3Q7(&fq[ic*Np*7], IonSolid);
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
}
2020-08-10 12:03:28 -04:00
// *************EVEN TIMESTEP*************//
timestep++;
//Update ion concentration and charge density
ScaLBL_Comm->SendD3Q7AA(fq, ic); //READ FORM NORMAL
2020-08-10 12:03:28 -04:00
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_DeviceBarrier();
// Set boundary conditions
if (BoundaryCondition == 1){
ScaLBL_Comm->D3Q7_Ion_Concentration_BC_z(NeighborList, &fq[ic*Np*7], Cin[ic], timestep);
ScaLBL_Comm->D3Q7_Ion_Concentration_BC_Z(NeighborList, &fq[ic*Np*7], Cout[ic], timestep);
}
//-------------------------//
2020-08-10 12:03:28 -04:00
ScaLBL_D3Q7_AAeven_IonConcentration(&fq[ic*Np*7],&Ci[ic*Np], 0, ScaLBL_Comm->LastExterior(), Np);
2020-08-10 12:03:28 -04:00
//LB-Ion collison
2020-08-10 12:03:28 -04:00
ScaLBL_D3Q7_AAeven_Ion(&fq[ic*Np*7],&Ci[ic*Np],Velocity,ElectricField,IonDiffusivity[ic],IonValence[ic],
2020-08-14 14:23:22 -04:00
rlx[ic],Vt,ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
2020-08-10 12:03:28 -04:00
ScaLBL_D3Q7_AAeven_Ion(&fq[ic*Np*7],&Ci[ic*Np],Velocity,ElectricField,IonDiffusivity[ic],IonValence[ic],
2020-08-14 14:23:22 -04:00
rlx[ic],Vt,0, ScaLBL_Comm->LastExterior(), Np);
if (BoundaryConditionSolid==1){
//TODO IonSolid may also be species-dependent
ScaLBL_Comm->SolidDirichletD3Q7(&fq[ic*Np*7], IonSolid);
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
}
}
}
//Compute charge density for Poisson equation
for (int ic=0; ic<number_ion_species; ic++){
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);
}
2020-08-06 15:41:40 -04:00
//************************************************************************/
2020-08-14 14:23:22 -04:00
//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;
2020-08-06 15:41:40 -04:00
2020-08-14 14:23:22 -04:00
//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");
2020-08-06 15:41:40 -04:00
}
void ScaLBL_IonModel::getIonConcentration(int timestep){
DoubleArray PhaseField(Nx,Ny,Nz);
for (int ic=0; ic<number_ion_species; ic++){
ScaLBL_Comm->RegularLayout(Map,&Ci[ic*Np],PhaseField);
ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
IonConcentration_LB_to_Phys(PhaseField);
FILE *OUTFILE;
sprintf(LocalRankFilename,"Ion%02i_Time_%i.%05i.raw",ic+1,timestep,rank);
OUTFILE = fopen(LocalRankFilename,"wb");
fwrite(PhaseField.data(),8,N,OUTFILE);
fclose(OUTFILE);
}
}
void ScaLBL_IonModel::IonConcentration_LB_to_Phys(DoubleArray &Den_reg){
for (int k=0;k<Nz;k++){
for (int j=0;j<Ny;j++){
for (int i=0;i<Nx;i++){
int idx=Map(i,j,k);
if (!(idx < 0)){
Den_reg(i,j,k) = Den_reg(i,j,k)/(h*h*h*1.0e-18);
}
}
}
}
}
void ScaLBL_IonModel::DummyFluidVelocity(){
double *FluidVelocity_host;
FluidVelocity_host = new double[3*Np];
for (int k=0; k<Nz; k++){
for (int j=0; j<Ny; j++){
for (int i=0; i<Nx; i++){
int idx=Map(i,j,k);
if (!(idx < 0))
FluidVelocity_host[idx+0*Np] = fluidVelx_dummy/(h*1.0e-6)*time_conv[0];
FluidVelocity_host[idx+1*Np] = fluidVely_dummy/(h*1.0e-6)*time_conv[0];
FluidVelocity_host[idx+2*Np] = fluidVelz_dummy/(h*1.0e-6)*time_conv[0];
}
}
}
ScaLBL_AllocateDeviceMemory((void **) &FluidVelocityDummy, sizeof(double)*3*Np);
ScaLBL_CopyToDevice(FluidVelocityDummy, FluidVelocity_host, sizeof(double)*3*Np);
ScaLBL_DeviceBarrier();
delete [] FluidVelocity_host;
}
void ScaLBL_IonModel::DummyElectricField(){
double *ElectricField_host;
ElectricField_host = new double[3*Np];
for (int k=0; k<Nz; k++){
for (int j=0; j<Ny; j++){
for (int i=0; i<Nx; i++){
int idx=Map(i,j,k);
if (!(idx < 0))
ElectricField_host[idx+0*Np] = Ex_dummy*(h*1.0e-6);
ElectricField_host[idx+1*Np] = Ey_dummy*(h*1.0e-6);
ElectricField_host[idx+2*Np] = Ez_dummy*(h*1.0e-6);
}
}
}
ScaLBL_AllocateDeviceMemory((void **) &ElectricFieldDummy, sizeof(double)*3*Np);
ScaLBL_CopyToDevice(ElectricFieldDummy, ElectricField_host, sizeof(double)*3*Np);
ScaLBL_DeviceBarrier();
delete [] ElectricField_host;
}
double ScaLBL_IonModel::CalIonDenConvergence(vector<double> &ci_avg_previous){
double *Ci_host;
Ci_host = new double[Np];
vector<double> error(number_ion_species,0.0);
for (int ic=0; ic<number_ion_species; ic++){
ScaLBL_CopyToHost(Ci_host,&Ci[ic*Np],Np*sizeof(double));
double count_loc=0;
double count;
double ci_avg;
double ci_loc=0.f;
for (int idx=0; idx<ScaLBL_Comm->LastExterior(); idx++){
ci_loc +=Ci_host[idx];
count_loc+=1.0;
}
for (int idx=ScaLBL_Comm->FirstInterior(); idx<ScaLBL_Comm->LastInterior(); idx++){
ci_loc +=Ci_host[idx];
count_loc+=1.0;
}
MPI_Allreduce(&ci_loc,&ci_avg,1,MPI_DOUBLE,MPI_SUM,Mask->Comm);
MPI_Allreduce(&count_loc,&count,1,MPI_DOUBLE,MPI_SUM,Mask->Comm);
ci_avg /= count;
double ci_avg_mag=ci_avg;
if (ci_avg==0.0) ci_avg_mag=1.0;
error[ic] = fabs(ci_avg-ci_avg_previous[ic])/fabs(ci_avg_mag);
ci_avg_previous[ic] = ci_avg;
}
double error_max;
error_max = *max_element(error.begin(),error.end());
if (rank==0){
printf("IonModel: error max: %.5g\n",error_max);
}
return error_max;
}
//void ScaLBL_IonModel::getIonConcentration(){
// for (int ic=0; ic<number_ion_species; ic++){
// ScaLBL_IonConcentration_Phys(Ci, h, ic, ScaLBL_Comm->FirstInterior(), ScaLBL_Comm->LastInterior(), Np);
// ScaLBL_IonConcentration_Phys(Ci, h, ic, 0, ScaLBL_Comm->LastExterior(), Np);
// }
//
// DoubleArray PhaseField(Nx,Ny,Nz);
// for (int ic=0; ic<number_ion_species; ic++){
// ScaLBL_Comm->RegularLayout(Map,&Ci[ic*Np],PhaseField);
// ScaLBL_DeviceBarrier(); MPI_Barrier(comm);
//
// FILE *OUTFILE;
// sprintf(LocalRankFilename,"Ion%02i.%05i.raw",ic+1,rank);
// OUTFILE = fopen(LocalRankFilename,"wb");
// fwrite(PhaseField.data(),8,N,OUTFILE);
// fclose(OUTFILE);
// }
//
//}