Scalable Lattice Boltzmann Library (ScalBL)

Functions

void ScaLBL_AllocateDeviceMemory(void **address, size_t size)

Allocate memory.

Parameters
  • address – memory address

  • size – size in bytes

void ScaLBL_AllocateZeroCopy(void **address, size_t size)

Allocate zero copy memory buffer (i.e. shared memory)

=*

Parameters
  • address – memory address

  • size – size in bytes

void ScaLBL_Color_BC_z(int *list, int *Map, double *Phi, double *Den, double vA, double vB, int count, int Np)
void ScaLBL_Color_BC_Z(int *list, int *Map, double *Phi, double *Den, double vA, double vB, int count, int Np)
void ScaLBL_CopySlice_z(double *Phi, int Nx, int Ny, int Nz, int Source, int Destination)
void ScaLBL_CopyToDevice(void *dest, const void *source, size_t size)

Copy memory from host to device.

Device memory should be close to simulation (based on NUMA cost) Host memory may be a shared memory region (with possibly higher NUMA cost for simulation) Analysis routine should minimize NUMA for host memory (based on process placement)

Parameters
  • dest – memory location to copy to

  • source – memory region to copy from

  • size – size of the region to copy in bytes

void ScaLBL_CopyToHost(void *dest, const void *source, size_t size)

Copy memory from device to host.

Device memory should be close to simulation (based on NUMA cost) Host memory may be a shared memory region (with possibly higher NUMA cost for simulation) Analysis routine should minimize NUMA for host memory (based on process placement)

Parameters
  • dest – memory location to copy to

  • source – memory region to copy from

  • size – size of the region to copy in bytes

void ScaLBL_CopyToZeroCopy(void *dest, const void *source, size_t size)

Copy memory from host to zero copy buffer.

Device memory should be close to simulation (based on NUMA cost) Host memory may be a shared memory region (with possibly higher NUMA cost for simulation) Analysis routine should minimize NUMA for host memory (based on process placement)

Parameters
  • dest – memory location to copy to

  • source – memory region to copy from

  • size – size of the region to copy in bytes

void ScaLBL_D3Q19_AAeven_BGK(double *dist, int start, int finish, int Np, double rlx, double Fx, double Fy, double Fz)

BGK collision based on AA even access pattern for D3Q19.

Parameters
  • dist – - D3Q19 distributions

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

  • rlx – - relaxation parameter

  • Fx – - force in x direction

  • Fy – - force in y direction

  • Fz – - force in z direction

void ScaLBL_D3Q19_AAeven_Color(int *Map, double *dist, double *Aq, double *Bq, double *Den, double *Phi, double *Vel, double rhoA, double rhoB, double tauA, double tauB, double alpha, double beta, double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np)

Color model collision based on AA even access pattern for D3Q19.

Parameters
  • Map – - mapping between sparse and dense data structures

  • dist – - D3Q19 distributions

  • Aq – - D3Q7 distribution for component A

  • Bq – - D3Q7 distribution for component B

  • Den – - density field

  • Phi – - phase indicator field

  • Vel – - velocity field

  • rhoA – - density of component A

  • rhoB – - density of component B

  • tauA – - relaxation time for component A

  • tauB – - relaxation time for component B

  • alpha – - parameter to control interfacial tension

  • beta – - parameter to control interface width

  • Fx – - force in x direction

  • Fy – - force in y direction

  • Fz – - force in z direction

  • strideY – - stride in y-direction for gradient computation

  • strideZ – - stride in z-direction for gradient computation

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q19_AAeven_Compact(double *d_dist, int Np)
void ScaLBL_D3Q19_AAeven_DFH(int *neighborList, double *dist, double *Aq, double *Bq, double *Den, double *Phi, double *Gradient, double *SolidForce, double rhoA, double rhoB, double tauA, double tauB, double alpha, double beta, double Fx, double Fy, double Fz, int start, int finish, int Np)
double ScaLBL_D3Q19_AAeven_Flux_BC_z(int *list, double *dist, double flux, double area, int count, int N)
void ScaLBL_D3Q19_AAeven_FreeLeeModel(int *Map, double *dist, double *Den, double *Phi, double *mu_phi, double *Vel, double *Pressure, double *ColorGrad, double rhoA, double rhoB, double tauA, double tauB, double kappa, double beta, double W, double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_D3Q19_AAeven_FreeLeeModel_Combined(int *Map, double *dist, double *hq, double *Den, double *Phi, double *mu_phi, double *Vel, double *Pressure, double *ColorGrad, double rhoA, double rhoB, double tauA, double tauB, double tauM, double kappa, double beta, double W, double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_D3Q19_AAeven_FreeLeeModel_Combined_HigherOrder(int *Map, double *dist, double *hq, double *Den, double *Phi, double *mu_phi, double *Vel, double *Pressure, double *ColorGrad, double rhoA, double rhoB, double tauA, double tauB, double tauM, double kappa, double beta, double W, double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_D3Q19_AAeven_FreeLeeModel_SingleFluid_BGK(double *dist, double *Vel, double *Pressure, double tau, double rho0, double Fx, double Fy, double Fz, int start, int finish, int Np)
void ScaLBL_D3Q19_AAeven_Greyscale(double *dist, int start, int finish, int Np, double rlx, double rlx_eff, double Fx, double Fy, double Fz, double *Poros, double *Perm, double *Velocity, double *Pressure)
void ScaLBL_D3Q19_AAeven_Greyscale_IMRT(double *dist, int start, int finish, int Np, double rlx, double rlx_eff, double Fx, double Fy, double Fz, double *Poros, double *Perm, double *Velocity, double Den, double *Pressure)
void ScaLBL_D3Q19_AAeven_Greyscale_MRT(double *dist, int start, int finish, int Np, double rlx, double rlx_eff, double Fx, double Fy, double Fz, double *Poros, double *Perm, double *Velocity, double Den, double *Pressure)
void ScaLBL_D3Q19_AAeven_GreyscaleColor(int *Map, double *dist, double *Aq, double *Bq, double *Den, double *Phi, double *GreySolidGrad, double *Poros, double *Perm, double *Vel, double *Pressure, double rhoA, double rhoB, double tauA, double tauB, double tauA_eff, double tauB_eff, double alpha, double beta, double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_D3Q19_AAeven_GreyscaleColor_CP(int *Map, double *dist, double *Aq, double *Bq, double *Den, double *Phi, double *GreySolidW, double *GreySn, double *GreySw, double *GreyKn, double *GreyKw, double *Poros, double *Perm, double *Vel, double *MobilityRatio, double *Pressure, double rhoA, double rhoB, double tauA, double tauB, double tauA_eff, double tauB_eff, double alpha, double beta, double Fx, double Fy, double Fz, bool RecoloringOff, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_D3Q19_AAeven_MRT(double *dist, int start, int finish, int Np, double rlx_setA, double rlx_setB, double Fx, double Fy, double Fz)

MRT collision based on AA even access pattern for D3Q19.

Parameters
  • dist – - D3Q19 distributions

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

  • rlx_setA – - relaxation parameter for viscous modes

  • rlx_setB – - relaxation parameter for non-viscous modes

  • Fx – - force in x direction

  • Fy – - force in y direction

  • Fz – - force in z direction

void ScaLBL_D3Q19_AAeven_Poisson(int *Map, double *dist, double *Den_charge, double *Psi, double *ElectricField, double *Error, double tau, double epsilon_LB, bool UseSlippingVelBC, int start, int finish, int Np)
void ScaLBL_D3Q19_AAeven_Pressure_BC_z(int *list, double *dist, double din, int count, int Np)
void ScaLBL_D3Q19_AAeven_Pressure_BC_Z(int *list, double *dist, double dout, int count, int Np)
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)
void ScaLBL_D3Q19_AAodd_BGK(int *neighborList, double *dist, int start, int finish, int Np, double rlx, double Fx, double Fy, double Fz)

BGK collision based on AA odd access pattern for D3Q19.

Parameters
  • neighborList – - neighbors based on D3Q19 lattice structure

  • dist – - D3Q19 distributions

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

  • rlx – - relaxation parameter

  • Fx – - force in x direction

  • Fy – - force in y direction

  • Fz – - force in z direction

void ScaLBL_D3Q19_AAodd_Color(int *NeighborList, int *Map, double *dist, double *Aq, double *Bq, double *Den, double *Phi, double *Vel, double rhoA, double rhoB, double tauA, double tauB, double alpha, double beta, double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np)

Color model collision based on AA even access pattern for D3Q19.

Parameters
  • NeighborList – - neighbors based on D3Q19 lattice structure

  • Map – - mapping between sparse and dense data structures

  • dist – - D3Q19 distributions

  • Aq – - D3Q7 distribution for component A

  • Bq – - D3Q7 distribution for component B

  • Den – - density field

  • Phi – - phase indicator field

  • Vel – - velocity field

  • rhoA – - density of component A

  • rhoB – - density of component B

  • tauA – - relaxation time for component A

  • tauB – - relaxation time for component B

  • alpha – - parameter to control interfacial tension

  • beta – - parameter to control interface width

  • Fx – - force in x direction

  • Fy – - force in y direction

  • Fz – - force in z direction

  • strideY – - stride in y-direction for gradient computation

  • strideZ – - stride in z-direction for gradient computation

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q19_AAodd_Compact(int *d_neighborList, double *d_dist, int Np)
void ScaLBL_D3Q19_AAodd_DFH(int *neighborList, double *dist, double *Aq, double *Bq, double *Den, double *Phi, double *Gradient, double *SolidForce, double rhoA, double rhoB, double tauA, double tauB, double alpha, double beta, double Fx, double Fy, double Fz, int start, int finish, int Np)
double ScaLBL_D3Q19_AAodd_Flux_BC_z(int *neighborList, int *list, double *dist, double flux, double area, int count, int N)
void ScaLBL_D3Q19_AAodd_FreeLeeModel(int *neighborList, int *Map, double *dist, double *Den, double *Phi, double *mu_phi, double *Vel, double *Pressure, double *ColorGrad, double rhoA, double rhoB, double tauA, double tauB, double kappa, double beta, double W, double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_D3Q19_AAodd_FreeLeeModel_Combined(int *neighborList, int *Map, double *dist, double *hq, double *Den, double *Phi, double *mu_phi, double *Vel, double *Pressure, double *ColorGrad, double rhoA, double rhoB, double tauA, double tauB, double tauM, double kappa, double beta, double W, double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_D3Q19_AAodd_FreeLeeModel_Combined_HigherOrder(int *neighborList, int *Map, double *dist, double *hq, double *Den, double *Phi, double *mu_phi, double *Vel, double *Pressure, double *ColorGrad, double rhoA, double rhoB, double tauA, double tauB, double tauM, double kappa, double beta, double W, double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_D3Q19_AAodd_FreeLeeModel_SingleFluid_BGK(int *neighborList, double *dist, double *Vel, double *Pressure, double tau, double rho0, double Fx, double Fy, double Fz, int start, int finish, int Np)
void ScaLBL_D3Q19_AAodd_Greyscale(int *neighborList, double *dist, int start, int finish, int Np, double rlx, double rlx_eff, double Fx, double Fy, double Fz, double *Poros, double *Perm, double *Velocity, double *Pressure)
void ScaLBL_D3Q19_AAodd_Greyscale_IMRT(int *neighborList, double *dist, int start, int finish, int Np, double rlx, double rlx_eff, double Fx, double Fy, double Fz, double *Poros, double *Perm, double *Velocity, double Den, double *Pressure)
void ScaLBL_D3Q19_AAodd_Greyscale_MRT(int *neighborList, double *dist, int start, int finish, int Np, double rlx, double rlx_eff, double Fx, double Fy, double Fz, double *Poros, double *Perm, double *Velocity, double Den, double *Pressure)
void ScaLBL_D3Q19_AAodd_GreyscaleColor(int *d_neighborList, int *Map, double *dist, double *Aq, double *Bq, double *Den, double *Phi, double *GreySolidGrad, double *Poros, double *Perm, double *Vel, double *Pressure, double rhoA, double rhoB, double tauA, double tauB, double tauA_eff, double tauB_eff, double alpha, double beta, double Fx, double Fy, double Fz, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_D3Q19_AAodd_GreyscaleColor_CP(int *d_neighborList, int *Map, double *dist, double *Aq, double *Bq, double *Den, double *Phi, double *GreySolidW, double *GreySn, double *GreySw, double *GreyKn, double *GreyKw, double *Poros, double *Perm, double *Vel, double *MobilityRatio, double *Pressure, double rhoA, double rhoB, double tauA, double tauB, double tauA_eff, double tauB_eff, double alpha, double beta, double Fx, double Fy, double Fz, bool RecoloringOff, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_D3Q19_AAodd_MRT(int *neighborList, double *dist, int start, int finish, int Np, double rlx_setA, double rlx_setB, double Fx, double Fy, double Fz)

MRT collision based on AA even access pattern for D3Q19.

Parameters
  • neighborList – - index of neighbors based on D3Q19 lattice structure

  • dist – - D3Q19 distributions

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

  • rlx_setA – - relaxation parameter for viscous modes

  • rlx_setB – - relaxation parameter for non-viscous modes

  • Fx – - force in x direction

  • Fy – - force in y direction

  • Fz – - force in z direction

void ScaLBL_D3Q19_AAodd_Poisson(int *neighborList, int *Map, double *dist, double *Den_charge, double *Psi, double *ElectricField, double tau, double epsilon_LB, bool UseSlippingVelBC, int start, int finish, int Np)
void ScaLBL_D3Q19_AAodd_Pressure_BC_z(int *neighborList, int *list, double *dist, double din, int count, int Np)
void ScaLBL_D3Q19_AAodd_Pressure_BC_Z(int *neighborList, int *list, double *dist, double dout, int count, int Np)
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, double rho0, double den_scale, double h, double time_conv, bool UseSlippingVelBC, int start, int finish, int Np)
void ScaLBL_D3Q19_FreeLeeModel_SingleFluid_Init(double *gqbar, double Fx, double Fy, double Fz, int Np)
void ScaLBL_D3Q19_FreeLeeModel_TwoFluid_Init(double *gqbar, double *mu_phi, double *ColorGrad, double Fx, double Fy, double Fz, int Np)
void ScaLBL_D3Q19_Gradient(int *Map, double *Phi, double *ColorGrad, int start, int finish, int Np, int Nx, int Ny, int Nz)
void ScaLBL_D3Q19_Gradient_DFH(int *NeighborList, double *Phi, double *ColorGrad, int start, int finish, int Np)
void ScaLBL_D3Q19_GreyIMRT_Init(double *Dist, int Np, double Den)
void ScaLBL_D3Q19_Init(double *Dist, int Np)

Initialize D3Q19 distributions.

Parameters
  • Dist – - D3Q19 distributions

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q19_MixedGradient(int *Map, double *Phi, double *Gradient, int start, int finish, int Np, int Nx, int Ny, int Nz)
void ScaLBL_D3Q19_Momentum(double *dist, double *vel, int Np)

Compute momentum from D3Q19 distribution.

Parameters
  • dist – - D3Q19 distributions

  • vel – - memory buffer to store the momentum that is computed

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q19_Pack(int q, int *list, int start, int count, double *sendbuf, double *dist, int N)

Pack D3Q19 distributions for communication.

Parameters
  • q – - index for distribution based on D3Q19 discrete velocity structure

  • list – - list of distributions to communicate

  • start – - index to start parsing the list

  • count – - number of values to pack

  • sendbuf – - memory buffer to hold values that will be sent

  • dist – - memory buffer to hold the distributions

  • N – - size of the distributions (derived from Domain structure)

void ScaLBL_D3Q19_Poisson_getElectricField(double *dist, double *ElectricField, double tau, int Np)
void ScaLBL_D3Q19_Poisson_Init(int *Map, double *dist, double *Psi, int start, int finish, int Np)
void ScaLBL_D3Q19_Pressure(double *dist, double *press, int Np)

compute pressure from D3Q19 distribution

Parameters
  • dist – - D3Q19 distributions

  • press – - memory buffer to store the pressure field that is computed

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q19_Reflection_BC_z(int *list, double *dist, int count, int Np)
void ScaLBL_D3Q19_Reflection_BC_Z(int *list, double *dist, int count, int Np)
void ScaLBL_D3Q19_Unpack(int q, int *list, int start, int count, double *recvbuf, double *dist, int N)

Unpack D3Q19 distributions after communication.

Parameters
  • q – - index for distribution based on D3Q19 discrete velocity structure

  • list – - list of distributions to communicate

  • start – - index to start parsing the list

  • count – - number of values to unppack

  • recvbuf – - memory buffer where recieved values have been stored

  • dist – - memory buffer to hold the distributions

  • N – - size of the distributions (derived from Domain structure)

void ScaLBL_D3Q7_AAeven_Color(int *Map, double *Aq, double *Bq, double *Den, double *Phi, double *ColorGrad, double *Vel, double rhoA, double rhoB, double beta, int start, int finish, int Np)
void ScaLBL_D3Q7_AAeven_DFH(double *Aq, double *Bq, double *Den, double *Phi, int start, int finish, int Np)
void ScaLBL_D3Q7_AAeven_FreeLee_PhaseField(int *Map, double *hq, double *Den, double *Phi, double *ColorGrad, double *Vel, double rhoA, double rhoB, double tauM, double W, int start, int finish, int Np)
void ScaLBL_D3Q7_AAeven_FreeLeeModel_PhaseField(int *Map, double *hq, double *Den, double *Phi, double rhoA, double rhoB, int start, int finish, int Np)
void ScaLBL_D3Q7_AAeven_Ion(double *dist, double *Den, double *FluxDiffusive, double *FluxAdvective, double *FluxElectrical, double *Velocity, double *ElectricField, double Di, int zi, double rlx, double Vt, int start, int finish, int Np)
void ScaLBL_D3Q7_AAeven_Ion_Concentration_BC_z(int *list, double *dist, double Cin, int count, int Np)
void ScaLBL_D3Q7_AAeven_Ion_Concentration_BC_Z(int *list, double *dist, double Cout, int count, int Np)
void ScaLBL_D3Q7_AAeven_Ion_Flux_Diff_BC_z(int *list, double *dist, double Cin, double tau, double *VelocityZ, int count, int Np)
void ScaLBL_D3Q7_AAeven_Ion_Flux_Diff_BC_Z(int *list, double *dist, double Cout, double tau, double *VelocityZ, int count, int Np)
void ScaLBL_D3Q7_AAeven_Ion_Flux_DiffAdvc_BC_z(int *list, double *dist, double Cin, double tau, double *VelocityZ, int count, int Np)
void ScaLBL_D3Q7_AAeven_Ion_Flux_DiffAdvc_BC_Z(int *list, double *dist, double Cout, double tau, double *VelocityZ, int count, int Np)
void ScaLBL_D3Q7_AAeven_Ion_Flux_DiffAdvcElec_BC_z(int *list, double *dist, double Cin, double tau, double *VelocityZ, double *ElectricField, double Di, double zi, double Vt, int count, int Np)
void ScaLBL_D3Q7_AAeven_Ion_Flux_DiffAdvcElec_BC_Z(int *list, double *dist, double Cout, double tau, double *VelocityZ, double *ElectricField, double Di, double zi, double Vt, int count, int Np)
void ScaLBL_D3Q7_AAeven_Ion_v0(double *dist, double *Den, double *FluxDiffusive, double *FluxAdvective, double *FluxElectrical, double *Velocity, double *ElectricField, double Di, int zi, double rlx, double Vt, int start, int finish, int Np)
void ScaLBL_D3Q7_AAeven_IonConcentration(double *dist, double *Den, int start, int finish, int Np)
void ScaLBL_D3Q7_AAeven_PhaseField(int *Map, double *Aq, double *Bq, double *Den, double *Phi, int start, int finish, int Np)

Compute phase field based on AA even access pattern for D3Q19.

Parameters
  • Map – - mapping between sparse and dense data structures

  • Aq – - D3Q7 distribution for component A

  • Bq – - D3Q7 distribution for component B

  • Den – - density field

  • Phi – - phase indicator field

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q7_AAeven_Poisson(int *Map, double *dist, double *Den_charge, double *Psi, double *ElectricField, double tau, double epsilon_LB, bool UseSlippingVelBC, int start, int finish, int Np)

Poisson-Boltzmann collision based on AA even access pattern for D3Q7.

Parameters
  • Map – - mapping between sparse and dense representations

  • dist – - D3Q7 distributions

  • Den_charge – - charge density

  • Psi – -

  • ElectricField – - electric field

  • tau – - relaxation time

  • epsilon_LB – - dielectric constant of medium

  • UseSlippingVelBC – - flag indicating the use of Helmholtz-Smoluchowski slipping velocity equation when EDL is too small to resolve

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q7_AAeven_Poisson_ElectricPotential(int *Map, double *dist, double *Psi, int start, int finish, int Np)

Poisson-Boltzmann solver / solve electric potential based on AA odd access pattern for D3Q7.

Parameters
  • Map – - mapping between sparse and dense representations

  • dist – - D3Q7 distributions

  • Psi – -

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q7_AAeven_Poisson_Potential_BC_z(int *list, double *dist, double Vin, int count, int Np)
void ScaLBL_D3Q7_AAeven_Poisson_Potential_BC_Z(int *list, double *dist, double Vout, int count, int Np)
void ScaLBL_D3Q7_AAodd_Color(int *neighborList, int *Map, double *Aq, double *Bq, double *Den, double *Phi, double *ColorGrad, double *Vel, double rhoA, double rhoB, double beta, int start, int finish, int Np)
void ScaLBL_D3Q7_AAodd_DFH(int *NeighborList, double *Aq, double *Bq, double *Den, double *Phi, int start, int finish, int Np)
void ScaLBL_D3Q7_AAodd_FreeLee_PhaseField(int *neighborList, int *Map, double *hq, double *Den, double *Phi, double *ColorGrad, double *Vel, double rhoA, double rhoB, double tauM, double W, int start, int finish, int Np)
void ScaLBL_D3Q7_AAodd_FreeLeeModel_PhaseField(int *neighborList, int *Map, double *hq, double *Den, double *Phi, double rhoA, double rhoB, int start, int finish, int Np)
void ScaLBL_D3Q7_AAodd_Ion(int *neighborList, double *dist, double *Den, double *FluxDiffusive, double *FluxAdvective, double *FluxElectrical, double *Velocity, double *ElectricField, double Di, int zi, double rlx, double Vt, int start, int finish, int Np)
void ScaLBL_D3Q7_AAodd_Ion_Concentration_BC_z(int *d_neighborList, int *list, double *dist, double Cin, int count, int Np)
void ScaLBL_D3Q7_AAodd_Ion_Concentration_BC_Z(int *d_neighborList, int *list, double *dist, double Cout, int count, int Np)
void ScaLBL_D3Q7_AAodd_Ion_Flux_Diff_BC_z(int *d_neighborList, int *list, double *dist, double Cin, double tau, double *VelocityZ, int count, int Np)
void ScaLBL_D3Q7_AAodd_Ion_Flux_Diff_BC_Z(int *d_neighborList, int *list, double *dist, double Cout, double tau, double *VelocityZ, int count, int Np)
void ScaLBL_D3Q7_AAodd_Ion_Flux_DiffAdvc_BC_z(int *d_neighborList, int *list, double *dist, double Cin, double tau, double *VelocityZ, int count, int Np)
void ScaLBL_D3Q7_AAodd_Ion_Flux_DiffAdvc_BC_Z(int *d_neighborList, int *list, double *dist, double Cout, double tau, double *VelocityZ, int count, int Np)
void ScaLBL_D3Q7_AAodd_Ion_Flux_DiffAdvcElec_BC_z(int *d_neighborList, int *list, double *dist, double Cin, double tau, double *VelocityZ, double *ElectricField, double Di, double zi, double Vt, int count, int Np)
void ScaLBL_D3Q7_AAodd_Ion_Flux_DiffAdvcElec_BC_Z(int *d_neighborList, int *list, double *dist, double Cout, double tau, double *VelocityZ, double *ElectricField, double Di, double zi, double Vt, int count, int Np)
void ScaLBL_D3Q7_AAodd_Ion_v0(int *neighborList, double *dist, double *Den, double *FluxDiffusive, double *FluxAdvective, double *FluxElectrical, double *Velocity, double *ElectricField, double Di, int zi, double rlx, double Vt, int start, int finish, int Np)
void ScaLBL_D3Q7_AAodd_IonConcentration(int *neighborList, double *dist, double *Den, int start, int finish, int Np)
void ScaLBL_D3Q7_AAodd_PhaseField(int *NeighborList, int *Map, double *Aq, double *Bq, double *Den, double *Phi, int start, int finish, int Np)

Compute phase field based on AA odd access pattern for D3Q19.

Parameters
  • NeighborList – - neighbors based on D3Q19 lattice structure

  • Map – - mapping between sparse and dense data structures

  • Aq – - D3Q7 distribution for component A

  • Bq – - D3Q7 distribution for component B

  • Den – - density field

  • Phi – - phase indicator field

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q7_AAodd_Poisson(int *neighborList, int *Map, double *dist, double *Den_charge, double *Psi, double *ElectricField, double tau, double epsilon_LB, bool UseSlippingVelBC, int start, int finish, int Np)

Poisson-Boltzmann collision based on AA odd access pattern for D3Q19.

Parameters
  • neighborList – - neighbors based on D3Q19 lattice structure

  • Map – - mapping between sparse and dense representations

  • dist – - D3Q7 distributions

  • Den_charge – - charge density

  • Psi – -

  • ElectricField – - electric field

  • tau – - relaxation time

  • epsilon_LB – -

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q7_AAodd_Poisson_ElectricPotential(int *neighborList, int *Map, double *dist, double *Psi, int start, int finish, int Np)

Poisson-Boltzmann solver / solve electric potential based on AA odd access pattern for D3Q7.

Parameters
  • neighborList – - neighbors based on D3Q19 lattice structure

  • Map – - mapping between sparse and dense representations

  • dist – - D3Q7 distributions

  • Psi – -

  • epsilon_LB – - dielectric constant of medium

  • UseSlippingVelBC – - flag indicating the use of Helmholtz-Smoluchowski slipping velocity equation when EDL is too small to resolve

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_z(int *d_neighborList, int *list, double *dist, double Vin, int count, int Np)
void ScaLBL_D3Q7_AAodd_Poisson_Potential_BC_Z(int *d_neighborList, int *list, double *dist, double Vout, int count, int Np)
void ScaLBL_D3Q7_ComputePhaseField(int *Map, double *hq, double *Den, double *Phi, double rhoA, double rhoB, int start, int finish, int Np)
void ScaLBL_D3Q7_Ion_ChargeDensity(double *Den, double *ChargeDensity, double IonValence, int ion_component, int start, int finish, int Np)
void ScaLBL_D3Q7_Ion_Init(double *dist, double *Den, double DenInit, int Np)
void ScaLBL_D3Q7_Ion_Init_FromFile(double *dist, double *Den, int Np)
void ScaLBL_D3Q7_Membrane_AssignLinkCoef(int *membrane, int *Map, double *Distance, double *Psi, double *coef, double Threshold, double MassFractionIn, double MassFractionOut, double ThresholdMassFractionIn, double ThresholdMassFractionOut, int memLinks, int Nx, int Ny, int Nz, int Np)
void ScaLBL_D3Q7_Membrane_AssignLinkCoef_halo(const int Cqx, const int Cqy, int const Cqz, int *Map, double *Distance, double *Psi, double Threshold, double MassFractionIn, double MassFractionOut, double ThresholdMassFractionIn, double ThresholdMassFractionOut, int *d3q7_recvlist, int *d3q7_linkList, double *coef, int start, int nlinks, int count, const int N, const int Nx, const int Ny, const int Nz)
void ScaLBL_D3Q7_Membrane_IonTransport(int *membrane, double *coef, double *dist, double *Den, int memLinks, int Np)
void ScaLBL_D3Q7_Membrane_Unpack(int q, int *d3q7_recvlist, double *recvbuf, int count, double *dist, int N, double *coef)
void ScaLBL_D3Q7_Poisson_Init(int *Map, double *dist, double *Psi, int start, int finish, int Np)

Initialize Poisson-Boltzmann solver.

Parameters
  • Map – - mapping between sparse and dense representations

  • dist – - D3Q7 distributions

  • Psi – -

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_D3Q7_Reflection_BC_z(int *list, double *dist, int count, int Np)
void ScaLBL_D3Q7_Reflection_BC_Z(int *list, double *dist, int count, int Np)
void ScaLBL_D3Q7_Unpack(int q, int *list, int start, int count, double *recvbuf, double *dist, int N)

Unpack D3Q7 distributions after communication.

Parameters
  • q – - index for distribution based on D3Q19 discrete velocity structure

  • list – - list of distributions to communicate

  • start – - index to start parsing the list

  • count – - number of values to unppack

  • recvbuf – - memory buffer where recieved values have been stored

  • dist – - memory buffer to hold the distributions

  • N – - size of the distributions (derived from Domain structure)

void ScaLBL_D3Q9_MGTest(int *Map, double *Phi, double *ColorGrad, int strideY, int strideZ, int start, int finish, int Np)
void ScaLBL_DeviceBarrier()

Device barrier routine.

void ScaLBL_DFH_Init(double *Phi, double *Den, double *Aq, double *Bq, int start, int finish, int Np)
void ScaLBL_FreeDeviceMemory(void *pointer)

Free memory.

Parameters

pointer – pointer to memory to free

void ScaLBL_FreeLeeModel_PhaseField_Init(int *Map, double *Phi, double *Den, double *hq, double *ColorGrad, double rhonA, double rhoB, double tauM, double W, int start, int finish, int Np)
void ScaLBL_Gradient_Unpack(double weight, double Cqx, double Cqy, double Cqz, int *list, int start, int count, double *recvbuf, double *phi, double *grad, int N)

Unpack values and compute Shan-Chen type of gradient.

Parameters
  • weight – - weight value for gradient sum

  • Cqx – - contribution to x-part of gradient

  • Cqy – - contribution to y-part of gradient

  • Cqz – - contribution to z-part of gradient

  • list – - list of distributions to communicate

  • start – - location to start reading the list

  • count – - number of values to unppack

  • recvbuf – - memory buffer where recieved values have been stored

  • phi – - scalar field

  • grad – - gradient

  • N – - size of local sub-domain (derived from Domain structure)

void ScaLBL_PhaseField_Init(int *Map, double *Phi, double *Den, double *Aq, double *Bq, int start, int finish, int Np)

Initialize phase field for color model.

Parameters
  • Map – - mapping between sparse and dense data structures

  • Phi – - phase indicator field

  • Den – - density field

  • Aq – - D3Q7 distribution for component A

  • Bq – - D3Q7 distribution for component B

  • start – - lattice node to start loop

  • finish – - lattice node to finish loop

  • Np – - size of local sub-domain (derived from Domain structure)

void ScaLBL_PhaseField_InitFromRestart(double *Den, double *Aq, double *Bq, int start, int finish, int Np)
void ScaLBL_Poisson_D3Q7_BC_z(int *list, int *Map, double *Psi, double Vin, int count)
void ScaLBL_Poisson_D3Q7_BC_Z(int *list, int *Map, double *Psi, double Vout, int count)
void ScaLBL_Scalar_Pack(int *list, int count, double *sendbuf, double *Data, int N)

Pack halo for scalar field to be prepare for communication.

Parameters
  • list – - list of distributions to communicate

  • count – - number of values to ppack

  • sendbuf – - memory buffer to pack values into

  • Data – - scalar field

  • N – - size of local sub-domain (derived from Domain structure)

void ScaLBL_Scalar_Unpack(int *list, int count, double *recvbuf, double *Data, int N)

Pack halo for scalar field to be prepare for communication.

Parameters
  • list – - list of distributions to communicate

  • count – - number of values to unppack

  • recvbuf – - memory buffer where recieved values have been stored

  • Data – - scalar field

  • N – - size of local sub-domain (derived from Domain structure)

int ScaLBL_SetDevice(int rank)

Set compute device.

Parameters

rank – rank of MPI process

void ScaLBL_SetSlice_z(double *Phi, double value, int Nx, int Ny, int Nz, int Slice)
void ScaLBL_Solid_Dirichlet_D3Q7(double *dist, double *BoundaryValue, int *BounceBackDist_list, int *BounceBackSolid_list, int N)
void ScaLBL_Solid_DirichletAndNeumann_D3Q7(double *dist, double *BoundaryValue, int *BoundaryLabel, int *BounceBackDist_list, int *BounceBackSolid_list, int N)
void ScaLBL_Solid_Neumann_D3Q7(double *dist, double *BoundaryValue, int *BounceBackDist_list, int *BounceBackSolid_list, int N)
void ScaLBL_Solid_SlippingVelocityBC_D3Q19(double *dist, double *zeta_potential, double *ElectricField, double *SolidGrad, double epsilon_LB, double tau, double rho0, double den_scale, double h, double time_conv, int *BounceBackDist_list, int *BounceBackSolid_list, int *FluidBoundary_list, double *lattice_weight, float *lattice_cx, float *lattice_cy, float *lattice_cz, int count, int Np)
class ScaLBL_Communicator
#include <ScaLBL.h>

Generalized communication routines for lattice Boltzmann methods.

Public Functions

inline void Barrier()

Create membrane data structure.

  • cut lattice links based on distance map

Parameters
  • Distance – - signed distance to membrane

  • neighborList – - data structure that retains lattice links

  • Np – - number of lattice sites

  • width – - halo width for the model

void BiRecvD3Q7AA(double *Aq, double *Bq)
void BiSendD3Q7AA(double *Aq, double *Bq)
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)
int copyRecvList(const char *dir, int *buffer)
int copySendList(const char *dir, int *buffer)
double D3Q19_Flux_BC_z(int *neighborList, double *fq, double flux, int time)
void D3Q19_Pressure_BC_z(int *neighborList, double *fq, double din, int time)
void D3Q19_Pressure_BC_Z(int *neighborList, double *fq, double dout, int time)
void D3Q19_Reflection_BC_z(double *fq)
void D3Q19_Reflection_BC_Z(double *fq)
void D3Q7_Ion_Concentration_BC_z(int *neighborList, double *fq, double Cin, int time)
void D3Q7_Ion_Concentration_BC_Z(int *neighborList, double *fq, double Cout, int time)
void D3Q7_Ion_Flux_Diff_BC_z(int *neighborList, double *fq, double Cin, double tau, double *VelocityZ, int time)
void D3Q7_Ion_Flux_Diff_BC_Z(int *neighborList, double *fq, double Cout, double tau, double *VelocityZ, int time)
void D3Q7_Ion_Flux_DiffAdvc_BC_z(int *neighborList, double *fq, double Cin, double tau, double *VelocityZ, int time)
void D3Q7_Ion_Flux_DiffAdvc_BC_Z(int *neighborList, double *fq, double Cout, double tau, double *VelocityZ, int time)
void D3Q7_Ion_Flux_DiffAdvcElec_BC_z(int *neighborList, double *fq, double Cin, double tau, double *VelocityZ, double *ElectricField_Z, double Di, double zi, double Vt, int time)
void D3Q7_Ion_Flux_DiffAdvcElec_BC_Z(int *neighborList, double *fq, double Cout, double tau, double *VelocityZ, double *ElectricField_Z, double Di, double zi, double Vt, 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)
int FirstInterior()
double GetPerformance(int *NeighborList, double *fq, int Np)
void GreyscaleSC_BC_z(int *Map, double *DenA, double *DenB, double vA, double vB)
void GreyscaleSC_BC_Z(int *Map, double *DenA, double *DenB, double vA, double vB)
void GreyscaleSC_Pressure_BC_z(int *neighborList, double *fqA, double *fqB, double dinA, double dinB, int time)
void GreyscaleSC_Pressure_BC_Z(int *neighborList, double *fqA, double *fqB, double doutA, double doutB, int time)
int LastExterior()
int LastInterior()
int MemoryOptimizedLayoutAA(IntArray &Map, int *neighborList, signed char *id, int Np, int width)
void Poisson_D3Q7_BC_z(int *Map, double *Psi, double Vin)
void Poisson_D3Q7_BC_Z(int *Map, double *Psi, double Vout)
void PrintD3Q19()
void RecvD3Q19AA(double *dist)
void RecvD3Q7AA(double *fq, int Component)
void RecvGrad(double *Phi, double *Gradient)
void RecvHalo(double *data)
void RegularLayout(IntArray map, const double *data, DoubleArray &regdata)
ScaLBL_Communicator(std::shared_ptr<Domain> Dm)

Constructor.

Parameters

Dm – - Domain information

void SendD3Q19AA(double *dist)
void SendD3Q7AA(double *fq, int Component)
void SendHalo(double *data)
void SetupBounceBackList(IntArray &Map, signed char *id, int Np, bool SlippingVelBC = false)
void SolidDirichletAndNeumannD3Q7(double *fq, double *BoundaryValue, int *BoundaryLabel)
void SolidDirichletD3Q7(double *fq, double *BoundaryValue)
void SolidNeumannD3Q7(double *fq, double *BoundaryValue)
void 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)
void TriRecvD3Q7AA(double *Aq, double *Bq, double *Cq)
void TriSendD3Q7AA(double *Aq, double *Bq, double *Cq)
~ScaLBL_Communicator()

Destructor.

Public Members

int BoundaryCondition
unsigned long int CommunicationCount
int first_interior
int last_interior
Utilities::MPI MPI_COMM_SCALBL
int N
int n_bb_d3q19
int n_bb_d3q7
int next
int Nx
int Ny
int Nz
int rank
int rank_x
int rank_X
int rank_xy
int rank_XY
int rank_xY
int rank_Xy
int rank_Xz
int rank_xz
int rank_XZ
int rank_xZ
int rank_y
int rank_Y
int rank_yz
int rank_YZ
int rank_yZ
int rank_Yz
int rank_Z
int rank_z
double *recvbuf_x
double *recvbuf_X
double *recvbuf_xy
double *recvbuf_Xy
double *recvbuf_xY
double *recvbuf_XY
double *recvbuf_xZ
double *recvbuf_xz
double *recvbuf_Xz
double *recvbuf_XZ
double *recvbuf_y
double *recvbuf_Y
double *recvbuf_yz
double *recvbuf_Yz
double *recvbuf_yZ
double *recvbuf_YZ
double *recvbuf_z
double *recvbuf_Z
unsigned long int RecvCount
int recvCount_X
int recvCount_x
int recvCount_xy
int recvCount_Xy
int recvCount_XY
int recvCount_xY
int recvCount_xz
int recvCount_xZ
int recvCount_XZ
int recvCount_Xz
int recvCount_Y
int recvCount_y
int recvCount_yz
int recvCount_Yz
int recvCount_YZ
int recvCount_yZ
int recvCount_z
int recvCount_Z
double *sendbuf_x
double *sendbuf_X
double *sendbuf_xy
double *sendbuf_Xy
double *sendbuf_xY
double *sendbuf_XY
double *sendbuf_xZ
double *sendbuf_xz
double *sendbuf_Xz
double *sendbuf_XZ
double *sendbuf_y
double *sendbuf_Y
double *sendbuf_Yz
double *sendbuf_yz
double *sendbuf_YZ
double *sendbuf_yZ
double *sendbuf_z
double *sendbuf_Z
unsigned long int SendCount
int sendCount_X
int sendCount_x
int sendCount_xy
int sendCount_Xy
int sendCount_XY
int sendCount_xY
int sendCount_xZ
int sendCount_xz
int sendCount_Xz
int sendCount_XZ
int sendCount_y
int sendCount_Y
int sendCount_yz
int sendCount_Yz
int sendCount_YZ
int sendCount_yZ
int sendCount_z
int sendCount_Z

Private Functions

void D3Q19_MapRecv(int Cqx, int Cqy, int Cqz, const int *list, int start, int count, int *d3q19_recvlist)
void start(std::vector<std::shared_ptr<MPI_Request>> &requests)
void wait(std::vector<std::shared_ptr<MPI_Request>> &requests)

Private Members

int *bb_dist
int *bb_interactions
int *dvcRecvDist_x
int *dvcRecvDist_X
int *dvcRecvDist_xy
int *dvcRecvDist_Xy
int *dvcRecvDist_xY
int *dvcRecvDist_XY
int *dvcRecvDist_xZ
int *dvcRecvDist_xz
int *dvcRecvDist_Xz
int *dvcRecvDist_XZ
int *dvcRecvDist_y
int *dvcRecvDist_Y
int *dvcRecvDist_yz
int *dvcRecvDist_Yz
int *dvcRecvDist_yZ
int *dvcRecvDist_YZ
int *dvcRecvDist_z
int *dvcRecvDist_Z
int *dvcRecvList_X
int *dvcRecvList_x
int *dvcRecvList_xy
int *dvcRecvList_Xy
int *dvcRecvList_xY
int *dvcRecvList_XY
int *dvcRecvList_xz
int *dvcRecvList_xZ
int *dvcRecvList_XZ
int *dvcRecvList_Xz
int *dvcRecvList_y
int *dvcRecvList_Y
int *dvcRecvList_yz
int *dvcRecvList_Yz
int *dvcRecvList_yZ
int *dvcRecvList_YZ
int *dvcRecvList_z
int *dvcRecvList_Z
int *dvcSendList_X
int *dvcSendList_x
int *dvcSendList_xy
int *dvcSendList_Xy
int *dvcSendList_XY
int *dvcSendList_xY
int *dvcSendList_xz
int *dvcSendList_xZ
int *dvcSendList_XZ
int *dvcSendList_Xz
int *dvcSendList_Y
int *dvcSendList_y
int *dvcSendList_yz
int *dvcSendList_Yz
int *dvcSendList_YZ
int *dvcSendList_yZ
int *dvcSendList_z
int *dvcSendList_Z
int *fluid_boundary
int i
int iproc
int j
int jproc
int k
int kproc
float *lattice_cx
float *lattice_cy
float *lattice_cz
double *lattice_weight
bool Lock
int n
int nprocx
int nprocy
int nprocz
RankInfoStruct rank_info
int recvtag
MPI_Request req1[18]
MPI_Request req2[18]
std::vector<std::shared_ptr<MPI_Request>> req_BiD3Q19AA
std::vector<std::shared_ptr<MPI_Request>> req_D3Q19AA
std::vector<std::shared_ptr<MPI_Request>> req_TriD3Q19AA
int sendtag