MaMiCo 1.2
Loading...
Searching...
No Matches
coupling::solvers::LBCouetteSolver Class Reference

implements a three-dimensional Lattice-Boltzmann Couette flow solver. More...

#include <LBCouetteSolver.h>

Inheritance diagram for coupling::solvers::LBCouetteSolver:
Collaboration diagram for coupling::solvers::LBCouetteSolver:

Public Member Functions

 LBCouetteSolver (const double channelheight, tarch::la::Vector< 3, double > wallVelocity, const double kinVisc, const double dx, const double dt, const int plotEveryTimestep, const std::string filestem, const tarch::la::Vector< 3, unsigned int > processes, const unsigned int numThreads=1, const Scenario *scen=nullptr)
 a simple constructor
 
virtual ~LBCouetteSolver ()
 a simple destructor
 
void advance (double dt) override
 advances one time step dt in time and triggers vtk plot if required
 
void setMDBoundaryValues (coupling::datastructures::FlexibleCellContainer< 3 > &md2macroBuffer) override
 applies the values received from the MD-solver within the conntinuum solver
 
tarch::la::Vector< 3, double > getVelocity (tarch::la::Vector< 3, double > pos) const override
 returns velocity at a certain position
 
double getDensity (tarch::la::Vector< 3, double > pos) const override
 returns density at a certain position
 
virtual void setWallVelocity (const tarch::la::Vector< 3, double > wallVelocity) override
 changes the velocity at the moving wall (z=0)
 
std::unique_ptr< StategetState () override
 
void setState (const std::unique_ptr< State > &input, int cycle) override
 
std::unique_ptr< Stateoperator() (const std::unique_ptr< State > &input, int cycle) override
 
Mode getMode () const override
 
std::unique_ptr< PintableMacroSolver > getSupervisor (int num_cycles, double visc_multiplier) const override
 
void print (std::ostream &os) const override
 
double get_avg_vel (const std::unique_ptr< State > &state) const override
 
- Public Member Functions inherited from coupling::solvers::NumericalSolver
 NumericalSolver (const double channelheight, const double dx, const double dt, const double kinVisc, const int plotEveryTimestep, const std::string filestem, const tarch::la::Vector< 3, unsigned int > processes, const Scenario *scen=nullptr)
 a simple constructor
 
virtual ~NumericalSolver ()
 a simple destructor
 
void setMDBoundary (tarch::la::Vector< 3, double > mdDomainOffset, tarch::la::Vector< 3, double > mdDomainSize, unsigned int overlapStrip)
 flags the domain boundary cells.
 
tarch::la::Vector< 3, unsigned int > getNumberProcesses () const
 returns the number of process, regards parallel runs
 
tarch::la::Vector< 3, unsigned int > getAvgNumberLBCells () const
 returns the average number of cells on each process
 
- Public Member Functions inherited from coupling::solvers::AbstractCouetteSolver< 3 >
virtual ~AbstractCouetteSolver ()
 a dummy destructor
 
virtual tarch::la::Vector< dim, double > getVelocity (tarch::la::Vector< dim, double > pos) const=0
 returns the current velocity at the given position
 
virtual void setWallVelocity (const tarch::la::Vector< dim, double > wallVelocity)=0
 changes the velocity at the moving for, refers to Couette scenario
 

Private Member Functions

void computeDensityAndVelocityEverywhere ()
 
void collidestream ()
 collide-stream algorithm for the Lattice-Boltzmann method
 
void stream (int index)
 the stream part of the LB algorithm (from pdf1 to pdf2)
 
void collide (int index, int x, int y, int z)
 
void boundary (double *const pdf, int index, int x, int y, int z, int q, const Flag &flag, int nbIndex)
 takes care of the correct boundary treatment for the LB method
 
void communicatePart (double *pdf, double *sendBuffer, double *recvBuffer, NbFlag nbFlagTo, NbFlag nbFlagFrom, tarch::la::Vector< 3, int > startSend, tarch::la::Vector< 3, int > endSend, tarch::la::Vector< 3, int > startRecv, tarch::la::Vector< 3, int > endRecv)
 
void communicate ()
 comunicates the boundary field data between the different processes
 

Static Private Member Functions

static void computeDensityAndVelocity (double *const vel, double &density, const double *const pdf)
 refers to the LB method; computes density and velocity on pdf
 

Private Attributes

Mode _mode
 
double _dt_pint
 
const double _omega
 relaxation frequency
 
tarch::la::Vector< 3, double > _wallVelocity
 velocity of moving wall of Couette flow
 
int _pdfsize {0}
 
double * _pdf1 {NULL}
 partical distribution function field
 
double * _pdf2 {NULL}
 partial distribution function field (stores the old time step)
 
const int _C [19][3]
 lattice velocities
 
const double _W [19]
 lattice weights
 

Additional Inherited Members

- Public Types inherited from coupling::interface::PintableMacroSolver
enum class  Mode { supervising , coupling }
 
using State = PintableMacroSolverState
 
- Static Public Member Functions inherited from coupling::solvers::NumericalSolver
static int getAvgDomainSize (double channelheight, double dx, tarch::la::Vector< 3, unsigned int > processes, int d)
 
- Protected Types inherited from coupling::solvers::NumericalSolver
enum  Flag {
  FLUID = 0 , NO_SLIP = 1 , MOVING_WALL = 2 , PERIODIC = 3 ,
  MD_BOUNDARY = 4 , PARALLEL_BOUNDARY = 5
}
 for every cell exists a flag entry, upon this is defined how the cell is handled More...
 
enum  NbFlag {
  LEFT = 0 , RIGHT = 1 , BACK = 2 , FRONT = 3 ,
  BOTTOM = 4 , TOP = 5
}
 The flags are used on parallel boundaries to define in which direction the boundary goes. More...
 
- Protected Member Functions inherited from coupling::solvers::NumericalSolver
int get (int i) const
 returns i and performs checks in debug mode
 
int get (int x, int y, int z) const
 returns linearized index and performs checks in debug mode
 
int getParBuf (int x, int y, int lengthx, int lengthy) const
 returns index in 2D parallel buffer with buffer dimensions lengthx+2,lengthy+2. Performs checks in debug mode
 
void plot () const
 
void plot (std::string filename) const
 create vtk plot if required
 
bool skipRank () const
 returns true, if this rank is not of relevance for the LB simulation
 
- Protected Attributes inherited from coupling::solvers::NumericalSolver
const double _channelheight
 the height and width of the channel in z and y direction
 
const double _dx
 mesh size, dx=dy=dz
 
const double _dt
 time step
 
const double _kinVisc
 kinematic viscosity of the fluid
 
tarch::la::Vector< 3, unsigned int > _processes
 domain decomposition on MPI rank basis; total number is given by multipling all entries
 
const int _plotEveryTimestep
 number of time steps between vtk plots
 
const std::string _filestem
 file stem for vtk plot
 
const int _domainSizeX {getDomainSize(_channelheight, _dx, _processes, 0)}
 domain size in x-direction
 
const int _domainSizeY {getDomainSize(_channelheight, _dx, _processes, 1)}
 domain size in y-direction
 
const int _domainSizeZ {getDomainSize(_channelheight, _dx, _processes, 2)}
 domain size in z-direction
 
const int _avgDomainSizeX {getAvgDomainSize(_channelheight, _dx, _processes, 0)}
 avg. domain size in MPI-parallel simulation in x-direction
 
const int _avgDomainSizeY {getAvgDomainSize(_channelheight, _dx, _processes, 1)}
 avg. domain size in MPI-parallel simulation in y-direction
 
const int _avgDomainSizeZ {getAvgDomainSize(_channelheight, _dx, _processes, 2)}
 avg. domain size in MPI-parallel simulation in z-direction
 
const tarch::la::Vector< 3, unsigned int > _coords {getProcessCoordinates()}
 coordinates of this process (=1,1,1, unless parallel run of the solver )
 
int _counter {0}
 time step counter
 
double * _vel {NULL}
 velocity field
 
double * _density {NULL}
 density field
 
Flag_flag {NULL}
 flag field
 
double * _sendBufferX {NULL}
 buffer to send data from left/right to right/left neighbour
 
double * _recvBufferX {NULL}
 buffer to receive data from from left/right neighbour
 
double * _sendBufferY {NULL}
 buffer to send data from front/back to front/back neighbour
 
double * _recvBufferY {NULL}
 buffer to receive data from from front/back neighbour
 
double * _sendBufferZ {NULL}
 buffer to send data from top/buttom to top/buttom neighbour
 
double * _recvBufferZ {NULL}
 buffer to receive data from from top/buttom neighbour
 
const int _xO {_domainSizeX + 2}
 offset for y-direction (lexicographic grid ordering)
 
const int _yO {(_domainSizeX + 2) * (_domainSizeY + 2)}
 offset for z-direction
 
tarch::la::Vector< 6, int > _parallelNeighbours {(-1)}
 neighbour ranks
 
tarch::la::Vector< 3, int > _offset {(-1)}
 offset of the md domain
 
tarch::la::Vector< 3, int > _globalNumberCouplingCells {(-1)}
 the total number of coupling cells of the coupled simulation
 
const Scenario_scen
 

Detailed Description

implements a three-dimensional Lattice-Boltzmann Couette flow solver.

In our scenario, the lower wall is accelerated and the upper wall stands still. The lower wall is located at zero height.

Author
Philipp Neumann

Constructor & Destructor Documentation

◆ LBCouetteSolver()

coupling::solvers::LBCouetteSolver::LBCouetteSolver ( const double channelheight,
tarch::la::Vector< 3, double > wallVelocity,
const double kinVisc,
const double dx,
const double dt,
const int plotEveryTimestep,
const std::string filestem,
const tarch::la::Vector< 3, unsigned int > processes,
const unsigned int numThreads = 1,
const Scenario * scen = nullptr )
inline

a simple constructor

Parameters
channelheightthe width and height of the channel in y and z direction
wallVelocityvelocity at the moving wall, refers to Couette scenario
dxthe spacial step size, and equidistant grid is applied
dtthe time step
kinViscthe kinematic viscosity of the fluid
plotEveryTimestepthe time step interval for plotting data; 4 means, every 4th time step is plotted
filestemthe name of the plotted file
processesdefines on how many processes the solver will run; 1,1,1 - sequential run - 1,2,2 = 1*2*2 = 4 processes
numThreadsnumber of OpenMP threads

Member Function Documentation

◆ advance()

void coupling::solvers::LBCouetteSolver::advance ( double dt)
inlineoverridevirtual

advances one time step dt in time and triggers vtk plot if required

Implements coupling::solvers::AbstractCouetteSolver< 3 >.

◆ boundary()

void coupling::solvers::LBCouetteSolver::boundary ( double *const pdf,
int index,
int x,
int y,
int z,
int q,
const Flag & flag,
int nbIndex )
inlineprivate

takes care of the correct boundary treatment for the LB method

Parameters
pdfparticle distribution function
indexstart index for current cell in pdf-array
xthe position in x direction of the cell
ythe position in y direction of the cell
zthe position in z direction of the cell
qdistribution function number
flagboundary flag of neighbouring cell
nbIndexindex of neighbouring cell

◆ collide()

void coupling::solvers::LBCouetteSolver::collide ( int index,
int x,
int y,
int z )
inlineprivate

@brieff the collide step within pdf2

◆ collidestream()

void coupling::solvers::LBCouetteSolver::collidestream ( )
inlineprivate

collide-stream algorithm for the Lattice-Boltzmann method

------------------------------------------------------------------------------------------------------ End of Pint methods —— End of Pint methods —— End of Pint methods —— End of Pint methods —————————————————————————————————— calls stream() and collide() and swaps the fields

◆ communicatePart()

void coupling::solvers::LBCouetteSolver::communicatePart ( double * pdf,
double * sendBuffer,
double * recvBuffer,
NbFlag nbFlagTo,
NbFlag nbFlagFrom,
tarch::la::Vector< 3, int > startSend,
tarch::la::Vector< 3, int > endSend,
tarch::la::Vector< 3, int > startRecv,
tarch::la::Vector< 3, int > endRecv )
inlineprivate

takes care of communication across one face in one direction.

Parameters
pdfpartial distribution function
sendBuffersend buffer
recvBufferreceive buffer
nbFlagTodirection into which message is sent
nbFlagFromdirection from which message is received
startSend3d coordinates that define the start of the data to be sent to neighbouring process
endSend3d coordinates that define the end of the data to to be sent to neighbouring process
startRecv3d coordinates that define the start of the data to be received from neighbouring process
endRecv3d coordinates that define the end of the data to be received from neighbouring process

◆ computeDensityAndVelocity()

static void coupling::solvers::LBCouetteSolver::computeDensityAndVelocity ( double *const vel,
double & density,
const double *const pdf )
inlinestaticprivate

refers to the LB method; computes density and velocity on pdf

Parameters
velvelocity
densitydensity
pdfpartial distribution function

◆ get_avg_vel()

double coupling::solvers::LBCouetteSolver::get_avg_vel ( const std::unique_ptr< State > & state) const
inlineoverridevirtual

Useful to test if a solver behaves as expected, similar to a "det(state)": should return a single numeric value that characterises the state object, by computing the mean flow velocity

Reimplemented from coupling::interface::PintableMacroSolver.

◆ getDensity()

double coupling::solvers::LBCouetteSolver::getDensity ( tarch::la::Vector< 3, double > pos) const
inlineoverridevirtual

returns density at a certain position

Parameters
posposition for which the density will be returned
Returns
the density vector for the position

Implements coupling::solvers::NumericalSolver.

◆ getMode()

Mode coupling::solvers::LBCouetteSolver::getMode ( ) const
inlineoverridevirtual

returns either Mode.supervising or Mode.coupling

Implements coupling::interface::PintableMacroSolver.

◆ getState()

std::unique_ptr< State > coupling::solvers::LBCouetteSolver::getState ( )
inlineoverridevirtual

------------------------------------------------------------------------------------------------------ Pint methods —— Pint methods —— Pint methods —— Pint methods —— Pint methods ——————————————————————————————————

Implements coupling::interface::PintableMacroSolver.

◆ getSupervisor()

std::unique_ptr< PintableMacroSolver > coupling::solvers::LBCouetteSolver::getSupervisor ( int num_cycles,
double visc_multiplier ) const
inlineoverridevirtual

This will create a new instance of this LBCouetteSolver In the supervisor, setMDBoundary() has not been called, independently of state of couette solver in coupling mode The supervisor can run with a modified viscosity The supervisor's operator() advances many coupling cycles at once, interval is stored in _dt_pint

Implements coupling::interface::PintableMacroSolver.

◆ getVelocity()

tarch::la::Vector< 3, double > coupling::solvers::LBCouetteSolver::getVelocity ( tarch::la::Vector< 3, double > pos) const
inlineoverride

returns velocity at a certain position

Parameters
posposition for which the velocity will be returned
Returns
the velocity vector for the position

◆ operator()()

std::unique_ptr< State > coupling::solvers::LBCouetteSolver::operator() ( const std::unique_ptr< State > & ,
int cycle )
inlineoverridevirtual

This advances the macroSolver by one time step, starting from the given state

Implements coupling::interface::PintableMacroSolver.

◆ print()

void coupling::solvers::LBCouetteSolver::print ( std::ostream & os) const
inlineoverridevirtual

this is expected to print information about the mode and type of this solver object to os

Implements coupling::interface::PintableMacroSolver.

◆ setMDBoundaryValues()

void coupling::solvers::LBCouetteSolver::setMDBoundaryValues ( coupling::datastructures::FlexibleCellContainer< 3 > & md2macroBuffer)
inlineoverridevirtual

applies the values received from the MD-solver within the conntinuum solver

Parameters
md2macroBufferholds the data from the md solver coupling cells

Implements coupling::solvers::NumericalSolver.

◆ setState()

void coupling::solvers::LBCouetteSolver::setState ( const std::unique_ptr< State > & ,
int cycle )
inlineoverridevirtual

This method fully reset / initialize the CFD solver to the given state

Implements coupling::interface::PintableMacroSolver.

◆ setWallVelocity()

virtual void coupling::solvers::LBCouetteSolver::setWallVelocity ( const tarch::la::Vector< 3, double > wallVelocity)
inlineoverridevirtual

changes the velocity at the moving wall (z=0)

Parameters
wallVelocitythe velocity will be set at the moving wall

Implements coupling::solvers::NumericalSolver.

Member Data Documentation

◆ _C

const int coupling::solvers::LBCouetteSolver::_C[19][3]
private
Initial value:
{{0, -1, -1}, {-1, 0, -1}, {0, 0, -1}, {1, 0, -1}, {0, 1, -1}, {-1, -1, 0}, {0, -1, 0}, {1, -1, 0}, {-1, 0, 0}, {0, 0, 0},
{1, 0, 0}, {-1, 1, 0}, {0, 1, 0}, {1, 1, 0}, {0, -1, 1}, {-1, 0, 1}, {0, 0, 1}, {1, 0, 1}, {0, 1, 1}}

lattice velocities

◆ _W

const double coupling::solvers::LBCouetteSolver::_W[19]
private
Initial value:
{1.0 / 36.0, 1.0 / 36.0, 1.0 / 18.0, 1.0 / 36.0, 1.0 / 36.0, 1.0 / 36.0, 1.0 / 18.0, 1.0 / 36.0, 1.0 / 18.0, 1.0 / 3.0,
1.0 / 18.0, 1.0 / 36.0, 1.0 / 18.0, 1.0 / 36.0, 1.0 / 36.0, 1.0 / 36.0, 1.0 / 18.0, 1.0 / 36.0, 1.0 / 36.0}

lattice weights


The documentation for this class was generated from the following file: