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

implements a simple one-dimensional finite-diffference solver for the Couette flow. More...

#include <FDCouetteSolver.h>

Inheritance diagram for coupling::solvers::FiniteDifferenceSolver:
Collaboration diagram for coupling::solvers::FiniteDifferenceSolver:

Public Member Functions

 FiniteDifferenceSolver (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)
 a simple constructor
 
 ~FiniteDifferenceSolver ()
 a simple destructor
 
void advance (double dt) override
 advances one time step dt in time and triggers vtk plot if required
 
tarch::la::Vector< 3, double > getVelocity (tarch::la::Vector< 3, double > pos) const override
 gets the velocity at a given position
 
virtual void setWallVelocity (const tarch::la::Vector< 3, double > wallVelocity) override
 changes the velocity at the moving wall (z=0)
 
double getDensity (tarch::la::Vector< 3, double > pos) const override
 returns density at a certain position
 
void setMDBoundaryValues (coupling::datastructures::FlexibleCellContainer< 3 > &md2macroBuffer) override
 applies the values received from the MD-solver within the conntinuum solver
 
- 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 setBeyondWall ()
 updates the velocities at boundaries
 
void update ()
 update the velocity in the cells which are flagged as FLUID
 

Private Attributes

const double _omega
 factor for the finite difference stencil = dt*kinVisc/(dx*dx)
 
tarch::la::Vector< 3, double > _wallVelocity
 velocity of moving wall of Couette flow
 
double * _velold {nullptr}
 the velocity field from the last time step
 

Additional Inherited Members

- 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 simple one-dimensional finite-diffference solver for the Couette flow.

In our scenario, the lower wall is accelerated and the upper wall stands still. The lower wall is located at zero height. The grid is just a simple cubic, equidistant mesh.

Author
Helene Wittenberg

Constructor & Destructor Documentation

◆ FiniteDifferenceSolver()

coupling::solvers::FiniteDifferenceSolver::FiniteDifferenceSolver ( 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 )
inline

a simple constructor

Parameters
channelheightthe width and height of the channel in y and z direction
kinViscthe kinematic viscosity of the fluid
wallVelocityvelocity at the moving wall, refers to Couette scenario
dxthe spacial step size, and equidistant grid is applied
dtthe time step
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::FiniteDifferenceSolver::advance ( double dt)
inlineoverridevirtual

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

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

◆ getDensity()

double coupling::solvers::FiniteDifferenceSolver::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
a density vector

Implements coupling::solvers::NumericalSolver.

◆ getVelocity()

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

gets the velocity at a given position

Parameters
posa position within the continuum domain
Returns
the velocity vector

◆ setBeyondWall()

void coupling::solvers::FiniteDifferenceSolver::setBeyondWall ( )
inlineprivate

updates the velocities at boundaries

there is no grid point on the wall, since the wall is located directly between the cells on the boundary. Therefore so set the boundary condition, the velocity boundary condition is applied by setting the velocity in the outer cell

◆ setMDBoundaryValues()

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

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

Parameters
recvBufferholds the data from the md solver
recvIndicethe indices to connect the data from the buffer with coupling cells

Implements coupling::solvers::NumericalSolver.

◆ setWallVelocity()

virtual void coupling::solvers::FiniteDifferenceSolver::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.


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