Saras
Finite   Difference   Solver   for   Fluid   Dynamics   Simulations
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
scalar Class Reference

The base class scalar to solve the incompressible Navier-Stokes equations with energy equation. More...

#include "lib/scalar.h"

Inheritance diagram for scalar:
hydro scalar_d2 scalar_d3

Public Member Functions

 scalar (const grid &mesh, const parser &solParam, parallel &mpiParam)
 Constructor of the base scalar class. More...
 
- Public Member Functions inherited from hydro
 hydro (const grid &mesh, const parser &solParam, parallel &mpiParam)
 Constructor of the base hydro class. More...
 
virtual void solvePDE ()
 The core publicly accessible function of the hydro class to solve the Navier-Stokes equations. More...
 
virtual real testPeriodic ()
 Function to test whether periodic BC is being implemented properly. More...
 

Public Attributes

sfield T
 The scalar field that stores the temperature field.
 
forcetForcing
 Instance of force class to handle temperature field forcing.
 
real nu
 
real kappa
 
- Public Attributes inherited from hydro
vfield V
 The vector field that stores the velocity field.
 
sfield P
 The scalar field that stores the pressure field.
 
forcevForcing
 Instance of force class to handle velocity field forcing.
 

Protected Member Functions

void initTBC ()
 Function to initialize the boundary conditions for temperature. More...
 
void imposeTBCs ()
 Function to impose the boundary conditions for temperature. More...
 
void initVForcing ()
 Function to initialize the forcing terms for velocity. More...
 
void initTForcing ()
 Function to initialize the forcing terms for temperature. More...
 
virtual void solveT ()
 Function to solve the implicit equation for scalar field. More...
 
- Protected Member Functions inherited from hydro
void checkPeriodic ()
 Function to enable/disable periodic data transfer as per the problem. More...
 
void setCoefficients ()
 Function to set the coefficients used for solving the implicit equations of U, V and W. More...
 
void initVBC ()
 Function to initialize the boundary conditions for velocity. More...
 
void imposeUBCs ()
 Function to impose the boundary conditions for the X-component of velocity. More...
 
void imposeVBCs ()
 Function to impose the boundary conditions for the Y-component of velocity. More...
 
void imposeWBCs ()
 Function to impose the boundary conditions for the Z-component of velocity. More...
 
void initVForcing ()
 Function to initialize the forcing terms for velocity. More...
 
virtual void solveVx ()
 Function to solve the implicit equation for x-velocity. More...
 
virtual void solveVy ()
 Function to solve the implicit equation for y-velocity. More...
 
virtual void solveVz ()
 Function to solve the implicit equation for z-velocity. More...
 
virtual void timeAdvance ()
 The subroutine to solve the NS equations using the implicit Crank-Nicholson method. More...
 

Protected Attributes

plainsf tmpRHS
 
plainsf guessedScalar
 
plainsf scalarLaplacian
 
boundarytLft
 
boundarytRgt
 
boundarytFrn
 
boundarytBak
 
boundarytTop
 
boundarytBot
 
- Protected Attributes inherited from hydro
int timeStepCount
 Integer value for the number of time-steps elapsed - it is incremented by 1 in each time-step.
 
int maxIterations
 Maximum number of iterations for the iterative solvers hydro::solveVx, hydro::solveVy and hydro::solveVz.
 
real time
 
real dt
 
real hx
 
real hy
 
real hz
 
real hx2
 
real hz2
 
real hz2hx2
 
real hx2hy2
 
real hy2hz2
 
real hx2hy2hz2
 
const gridmesh
 
const parserinputParams
 
const real inverseRe
 
probesdataProbe
 Instance of the probe class to collect data from probes in the domain.
 
boundaryuLft
 Instances of the boundary class to impose boundary conditions on all the 6 walls for the 3 components of the velocity field.
 
boundaryuRgt
 
boundaryuFrn
 
boundaryuBak
 
boundaryuTop
 
boundaryuBot
 
boundaryvLft
 
boundaryvRgt
 
boundaryvFrn
 
boundaryvBak
 
boundaryvTop
 
boundaryvBot
 
boundarywLft
 
boundarywRgt
 
boundarywFrn
 
boundarywBak
 
boundarywTop
 
boundarywBot
 
parallelmpiData
 Instance of the parallel class that holds the MPI-related data like rank, xRank, etc.
 
plainsf Pp
 Plain scalar field into which the pressure correction is calculated and written by the Poisson solver.
 
plainsf mgRHS
 Plain scalar field into which the RHS for pressure Poisson equation is written and passed to the Poisson solver.
 
plainvf nseRHS
 Plain vector field into which the RHS of the Navier-Stokes equation is written and stored.
 
plainvf velocityLaplacian
 Plain vector field into which the RHS of the implicit equation for velocities are calculated during iterative solving.
 
plainvf pressureGradient
 Plain vector field which stores the pressure gradient term.
 
plainvf guessedVelocity
 Plain vector field which serves as a temporary array during iterative solution procedure for velocity terms.
 

Detailed Description

The base class scalar to solve the incompressible Navier-Stokes equations with energy equation.

The class initializes and stores the velocity vector field and the pressure scalar field along with a few auxilliary fields to solve the PDE. In addition to the fields in hydro class, the scalar equation is also solved here. It solves the NSE using the solvePDE function from within which the implicit Crank-Nicholson method is used to solve the PDE.

Constructor & Destructor Documentation

◆ scalar()

scalar::scalar ( const grid mesh,
const parser solParam,
parallel mpiParam 
)

Constructor of the base scalar class.

     The short base constructor of the scalar class merely assigns the const references to the grid and parser
     class instances being used in the solver.
     Also, the maximum allowable number of iterations for the Jacobi iterative solver being used to solve for the
     velocities implicitly is set as \f$ N_{max} = N_x \times N_y \times N_z \f$, where \f$N_x\f$, \f$N_y\f$ and \f$N_z\f$
     are the number of grid points in the collocated grid at the local sub-domains along x, y and z directions
     respectively.
Parameters
meshis a const reference to the global data contained in the grid class
solParamis a const reference to the user-set parameters contained in the parser class
mpiParamis a reference to the object of parallel class containing the necessary rank data

Member Function Documentation

◆ imposeTBCs()

void scalar::imposeTBCs ( )
protected

Function to impose the boundary conditions for temperature.

     The function first calls the syncData() function of the temperature field to update the sub-domain pads.
     Then the boundary conditions are applied at the full domain boundaries by calling the imposeBC()
     of each boundary class object assigned to each wall.

◆ initTBC()

void scalar::initTBC ( )
protected

Function to initialize the boundary conditions for temperature.

     The temperature boundary conditions for all the 6 walls (4 in case of 2D simulations) are initialized here.
     Out of the different boundary conditions available in the boundary class,
     the appropriate BCs are chosen according to the type of problem being solved.

◆ initTForcing()

void scalar::initTForcing ( )
protected

Function to initialize the forcing terms for temperature.

     The forcing terms for the temperature field are initialized here.
     Out of the different forcings available in the force class,
     the appropriate forcing is chosen according to the parameters set by the user.

◆ initVForcing()

void scalar::initVForcing ( )
protected

Function to initialize the forcing terms for velocity.

     The forcing terms for the velocity field are initialized here.
     Out of the different forcings available in the force class,
     the appropriate forcing is chosen according to the parameters set by the user.

◆ solveT()

void scalar::solveT ( )
protectedvirtual

Function to solve the implicit equation for scalar field.

     The implicit equation for \f$ \theta' \f$ of the implicit Crank-Nicholson method is solved using the Jacobi
     iterative method here.

     The loop exits when the global maximum of the error in computed solution obtained using the \ref plainsf#fxMax "fxMax" function
     of scalar fields in sfield.h falls below the specified tolerance.
     If the solution doesn't converge even after an internally assigned maximum number for iterations, the solver
     aborts with an error message.

     Note that this function uses the blitz index place holders firstIndex, secondIndex and thirdIndex.
     They are declared as i, j, and k respectively.
     Hence the variables i, j and k are not scalars in this function.

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