cleanup with help from valgrind
This commit is contained in:
@@ -274,6 +274,7 @@ Membrane::~Membrane() {
|
||||
ScaLBL_FreeDeviceMemory( NeighborList );
|
||||
ScaLBL_FreeDeviceMemory( MembraneLinks );
|
||||
ScaLBL_FreeDeviceMemory( MembraneCoef );
|
||||
ScaLBL_FreeDeviceMemory( MembraneDistance );
|
||||
|
||||
ScaLBL_FreeDeviceMemory( sendbuf_x );
|
||||
ScaLBL_FreeDeviceMemory( sendbuf_X );
|
||||
@@ -730,7 +731,6 @@ int Membrane::Create(std::shared_ptr <Domain> Dm, DoubleArray &Distance, IntArra
|
||||
}
|
||||
|
||||
if (rank == 0) printf(" Create device data structures... \n");
|
||||
|
||||
/* Create device copies of data structures */
|
||||
ScaLBL_AllocateDeviceMemory((void **)&MembraneLinks, 2*mlink*sizeof(int));
|
||||
ScaLBL_AllocateDeviceMemory((void **)&MembraneCoef, 2*mlink*sizeof(double));
|
||||
@@ -995,6 +995,7 @@ int Membrane::Create(std::shared_ptr <Domain> Dm, DoubleArray &Distance, IntArra
|
||||
ScaLBL_AllocateZeroCopy((void **) &coefficient_Z, (recvCount_Z - linkCount_Z[0])*sizeof(double));
|
||||
//......................................................................................
|
||||
|
||||
delete [] neighborList;
|
||||
delete [] TempBuffer;
|
||||
return mlink;
|
||||
}
|
||||
@@ -1356,36 +1357,42 @@ void Membrane::AssignCoefficients(int *Map, double *Psi, double Threshold,
|
||||
double ThresholdMassFractionOut){
|
||||
/* Assign mass transfer coefficients to the membrane data structure */
|
||||
|
||||
|
||||
if (membraneLinkCount > 0)
|
||||
ScaLBL_D3Q7_Membrane_AssignLinkCoef(MembraneLinks, Map, MembraneDistance, Psi, MembraneCoef,
|
||||
Threshold, MassFractionIn, MassFractionOut, ThresholdMassFractionIn, ThresholdMassFractionOut,
|
||||
membraneLinkCount, Nx, Ny, Nz, Np);
|
||||
|
||||
if (linkCount_X[0] < recvCount_X)
|
||||
ScaLBL_D3Q7_Membrane_AssignLinkCoef_halo(-1,0,0,Map,MembraneDistance,Psi,Threshold,
|
||||
MassFractionIn,MassFractionOut,ThresholdMassFractionIn,ThresholdMassFractionOut,
|
||||
dvcRecvDist_X,dvcRecvLinks_X,coefficient_X,0,linkCount_X[0],recvCount_X,
|
||||
Np,Nx,Ny,Nz);
|
||||
|
||||
if (linkCount_x[0] < recvCount_x)
|
||||
ScaLBL_D3Q7_Membrane_AssignLinkCoef_halo(1,0,0,Map,MembraneDistance,Psi,Threshold,
|
||||
MassFractionIn,MassFractionOut,ThresholdMassFractionIn,ThresholdMassFractionOut,
|
||||
dvcRecvDist_x,dvcRecvLinks_x,coefficient_x,0,linkCount_x[0],recvCount_x,
|
||||
Np,Nx,Ny,Nz);
|
||||
|
||||
if (linkCount_Y[0] < recvCount_Y)
|
||||
ScaLBL_D3Q7_Membrane_AssignLinkCoef_halo(0,-1,0,Map,MembraneDistance,Psi,Threshold,
|
||||
MassFractionIn,MassFractionOut,ThresholdMassFractionIn,ThresholdMassFractionOut,
|
||||
dvcRecvDist_Y,dvcRecvLinks_Y,coefficient_Y,0,linkCount_Y[0],recvCount_Y,
|
||||
Np,Nx,Ny,Nz);
|
||||
|
||||
if (linkCount_y[0]<recvCount_y)
|
||||
ScaLBL_D3Q7_Membrane_AssignLinkCoef_halo(0,1,0,Map,MembraneDistance,Psi,Threshold,
|
||||
MassFractionIn,MassFractionOut,ThresholdMassFractionIn,ThresholdMassFractionOut,
|
||||
dvcRecvDist_y,dvcRecvLinks_y,coefficient_y,0,linkCount_y[0],recvCount_y,
|
||||
Np,Nx,Ny,Nz);
|
||||
|
||||
if (linkCount_Z[0]<recvCount_Z)
|
||||
ScaLBL_D3Q7_Membrane_AssignLinkCoef_halo(0,0,-1,Map,MembraneDistance,Psi,Threshold,
|
||||
MassFractionIn,MassFractionOut,ThresholdMassFractionIn,ThresholdMassFractionOut,
|
||||
dvcRecvDist_Z,dvcRecvLinks_Z,coefficient_Z,0,linkCount_Z[0],recvCount_Z,
|
||||
Np,Nx,Ny,Nz);
|
||||
|
||||
if (linkCount_z[0]<recvCount_z)
|
||||
ScaLBL_D3Q7_Membrane_AssignLinkCoef_halo(0,0,1,Map,MembraneDistance,Psi,Threshold,
|
||||
MassFractionIn,MassFractionOut,ThresholdMassFractionIn,ThresholdMassFractionOut,
|
||||
dvcRecvDist_z,dvcRecvLinks_z,coefficient_z,0,linkCount_z[0],recvCount_z,
|
||||
|
||||
@@ -1114,6 +1114,7 @@ void ScaLBL_Communicator::SetupBounceBackList(IntArray &Map, signed char *id, in
|
||||
}
|
||||
}
|
||||
}
|
||||
if (local_count > 0){
|
||||
|
||||
int *bb_dist_tmp = new int [local_count];
|
||||
int *bb_interactions_tmp = new int [local_count];
|
||||
@@ -1410,6 +1411,16 @@ void ScaLBL_Communicator::SetupBounceBackList(IntArray &Map, signed char *id, in
|
||||
delete [] lattice_cx_tmp;
|
||||
delete [] lattice_cy_tmp;
|
||||
delete [] lattice_cz_tmp;
|
||||
}
|
||||
else {
|
||||
bb_dist = NULL;
|
||||
bb_interactions = NULL;
|
||||
fluid_boundary = NULL;
|
||||
lattice_weight = NULL;
|
||||
lattice_cx = NULL;
|
||||
lattice_cy = NULL;
|
||||
lattice_cz = NULL;
|
||||
}
|
||||
}
|
||||
|
||||
void ScaLBL_Communicator::SolidDirichletD3Q7(double *fq, double *BoundaryValue){
|
||||
|
||||
@@ -13,7 +13,20 @@ ScaLBL_IonModel::ScaLBL_IonModel(int RANK, int NP, const Utilities::MPI &COMM)
|
||||
nprocy(0), nprocz(0), fluidVelx_dummy(0), fluidVely_dummy(0),
|
||||
fluidVelz_dummy(0), BoundaryConditionInlet(0), BoundaryConditionOutlet(0),
|
||||
BoundaryConditionSolid(0), Lx(0), Ly(0), Lz(0), comm(COMM) {}
|
||||
ScaLBL_IonModel::~ScaLBL_IonModel() {}
|
||||
|
||||
ScaLBL_IonModel::~ScaLBL_IonModel() {
|
||||
|
||||
ScaLBL_FreeDeviceMemory(NeighborList);
|
||||
ScaLBL_FreeDeviceMemory(dvcMap);
|
||||
ScaLBL_FreeDeviceMemory(fq);
|
||||
ScaLBL_FreeDeviceMemory(Ci);
|
||||
ScaLBL_FreeDeviceMemory(ChargeDensity);
|
||||
ScaLBL_FreeDeviceMemory(FluxDiffusive);
|
||||
ScaLBL_FreeDeviceMemory(FluxAdvective);
|
||||
ScaLBL_FreeDeviceMemory(FluxElectrical);
|
||||
ScaLBL_FreeDeviceMemory(IonSolid);
|
||||
ScaLBL_FreeDeviceMemory(FluidVelocityDummy);
|
||||
}
|
||||
|
||||
void ScaLBL_IonModel::ReadParams(string filename, vector<int> &num_iter) {
|
||||
|
||||
@@ -50,6 +63,8 @@ void ScaLBL_IonModel::ReadParams(string filename, vector<int> &num_iter) {
|
||||
Ez_dummy = 0.0; //for debugging, unit [V/m]
|
||||
//--------------------------------------------------------------------------//
|
||||
|
||||
BoundaryConditionSolid = 0;
|
||||
|
||||
// Read domain parameters
|
||||
if (domain_db->keyExists("voxel_length")) { //default unit: um/lu
|
||||
h = domain_db->getScalar<double>("voxel_length");
|
||||
@@ -685,7 +700,6 @@ void ScaLBL_IonModel::SetDomain() {
|
||||
}
|
||||
|
||||
void ScaLBL_IonModel::SetMembrane() {
|
||||
size_t NLABELS = 0;
|
||||
|
||||
membrane_db = db->getDatabase("Membrane");
|
||||
|
||||
@@ -693,7 +707,7 @@ void ScaLBL_IonModel::SetMembrane() {
|
||||
auto MembraneLabels = membrane_db->getVector<int>("MembraneLabels");
|
||||
|
||||
IonMembrane = std::shared_ptr<Membrane>(new Membrane(Dm, NeighborList, Np));
|
||||
|
||||
size_t NLABELS = MembraneLabels.size();
|
||||
signed char LABEL = 0;
|
||||
double *label_count;
|
||||
double *label_count_global;
|
||||
@@ -1020,7 +1034,6 @@ void ScaLBL_IonModel::Create() {
|
||||
}
|
||||
ScaLBL_CopyToDevice(dvcMap, TmpMap, sizeof(int) * Np);
|
||||
ScaLBL_Comm->Barrier();
|
||||
delete[] TmpMap;
|
||||
|
||||
ScaLBL_CopyToDevice(NeighborList, neighborList, neighborSize);
|
||||
comm.barrier();
|
||||
@@ -1042,6 +1055,13 @@ void ScaLBL_IonModel::Create() {
|
||||
ScaLBL_Comm->Barrier();
|
||||
delete[] IonSolid_host;
|
||||
}
|
||||
else {
|
||||
IonSolid = NULL;
|
||||
}
|
||||
|
||||
|
||||
delete[] TmpMap;
|
||||
delete[] neighborList;
|
||||
}
|
||||
|
||||
void ScaLBL_IonModel::Initialize() {
|
||||
@@ -1444,6 +1464,7 @@ void ScaLBL_IonModel::RunMembrane(double *Velocity, double *ElectricField, doubl
|
||||
//auto t1 = std::chrono::system_clock::now();
|
||||
|
||||
for (size_t ic = 0; ic < number_ion_species; ic++) {
|
||||
|
||||
/* set the mass transfer coefficients for the membrane */
|
||||
IonMembrane->AssignCoefficients(dvcMap, Psi, ThresholdVoltage[ic],MassFractionIn[ic],
|
||||
MassFractionOut[ic],ThresholdMassFractionIn[ic],ThresholdMassFractionOut[ic]);
|
||||
@@ -1525,8 +1546,6 @@ void ScaLBL_IonModel::RunMembrane(double *Velocity, double *ElectricField, doubl
|
||||
//Compute charge density for Poisson equation
|
||||
for (size_t ic = 0; ic < number_ion_species; ic++) {
|
||||
int Valence = IonValence[ic];
|
||||
if (rank==0) printf("compute charge density for ion %i, Valence =%i \n", ic,Valence);
|
||||
|
||||
ScaLBL_D3Q7_Ion_ChargeDensity(Ci, ChargeDensity, Valence, ic,
|
||||
ScaLBL_Comm->FirstInterior(),
|
||||
ScaLBL_Comm->LastInterior(), Np);
|
||||
|
||||
@@ -32,6 +32,14 @@ ScaLBL_Poisson::ScaLBL_Poisson(int RANK, int NP, const Utilities::MPI& COMM):
|
||||
}
|
||||
ScaLBL_Poisson::~ScaLBL_Poisson()
|
||||
{
|
||||
ScaLBL_FreeDeviceMemory(NeighborList);
|
||||
ScaLBL_FreeDeviceMemory(dvcMap);
|
||||
ScaLBL_FreeDeviceMemory(Psi);
|
||||
ScaLBL_FreeDeviceMemory(Psi_BCLabel);
|
||||
ScaLBL_FreeDeviceMemory(ElectricField);
|
||||
ScaLBL_FreeDeviceMemory(ResidualError);
|
||||
ScaLBL_FreeDeviceMemory(fq);
|
||||
|
||||
if ( TIMELOG )
|
||||
fclose( TIMELOG );
|
||||
}
|
||||
@@ -416,6 +424,7 @@ void ScaLBL_Poisson::Create(){
|
||||
//ScaLBL_Comm->Barrier();
|
||||
|
||||
//Initialize solid boundary for electric potential
|
||||
// DON'T USE WITH CELLULAR SYSTEM (NO SOLID -- NEED Membrane SOLUTION)
|
||||
ScaLBL_Comm->SetupBounceBackList(Map, Mask->id.data(), Np);
|
||||
comm.barrier();
|
||||
}
|
||||
|
||||
@@ -23,7 +23,9 @@ using namespace std;
|
||||
int main(int argc, char **argv)
|
||||
{
|
||||
// Initialize MPI and error handlers
|
||||
Utilities::startup( argc, argv );
|
||||
//Utilities::startup( argc, argv );
|
||||
Utilities::startup( argc, argv, true );
|
||||
|
||||
Utilities::MPI comm( MPI_COMM_WORLD );
|
||||
int rank = comm.getRank();
|
||||
int nprocs = comm.getSize();
|
||||
|
||||
Reference in New Issue
Block a user