25 template <
typename memoryType>
35 cusp::coo_matrix<int, real, memoryType>
45 cusp::array1d<real, memoryType>
175 return "Navier-Stokes";
Monitors the time spent to achieve a certain task.
cusp::array1d< real, memoryType > bc1
inhomogeneous boundary conditions from the discrete Laplacian operator
Specifies the time-integration scheme used.
parameterDB * paramDB
database that contains all the simulation parameters
virtual void generateA(real alpha)
cusp::coo_matrix< int, real, memoryType > C
matrix of the Poisson system
cusp::array1d< real, memoryType > H
hold convective terms from previous time step
bool finished()
Evaluates the condition required to stop the simulation.
Definition of custom types required by the code.
virtual void generateBC1()
size_t subStep
number of sub-iterations
cusp::coo_matrix< int, real, memoryType > Q
gradient matrix (and regularization matrix if immersed body)
NavierStokesSolver(parameterDB *pDB=NULL, domain *dInfo=NULL)
Constructor. Copies the database and information about the computational grid.
cusp::array1d< real, memoryType > lambda
pressure vector (and body forces if immersed body)
double real
Is a float or a double depending on the machine precision.
cusp::coo_matrix< int, real, memoryType > Minv
inverse of the mass matrix
cusp::coo_matrix< int, real, memoryType > BN
N-th order Taylor series expansion of the inverse of A.
cusp::array1d< real, memoryType > rhs2
right hand-side of the Poisson system
void stepTime()
Calculates the variables at the next time step.
virtual void calculateExplicitLambdaTerms()
Doing nothing. Used in methods that use the explicit pressure term in the intermediate velocity solve...
cusp::array1d< real, memoryType > rhs1
right hand-side of the system for the intermediate velocity flux vector
Stores the preconditioner for a given system.
virtual void initialiseFluxes()
void calculateExplicitQTerms()
void initialiseCommon()
Initializes parameters common to all Navier-Stokes solvers.
void generateBN()
Generates approximate inverse of the matrix resulting from implicit velocity terms.
size_t iterationCount1
number of iteration to solve the intermediate fluxes
cusp::array1d< real, memoryType > temp2
temporary 1D Cusp array
cusp::coo_matrix< int, real, memoryType > L
discrete Laplacian matrix
size_t iterationCount2
number of iteration to solve the Poisson equation
Declaration of the functions of the namespace io.
size_t timeStep
iteration number
std::map< std::string, componentParameter > parameterDB
Map from a string to a componentParameter.
cusp::array1d< real, memoryType > q
velocity flux vector
void generateRN()
Generates explicit terms of the momentum equation.
cusp::array1d< real, memoryType > qOld
velocity flux vector at the previous time-step
cusp::array1d< real, memoryType > qStar
intermediate velocity flux vector
virtual std::string name()
Returns the name of the current solver.
virtual void shutDown()
Prints timing information and closes the different files.
Stores information about the computational grid.
void initialiseArrays(int numQ, int numLambda)
Initializes all arrays required to solve the Navier-Stokes equations.
cusp::array1d< real, memoryType > rn
explicit terms of the momentum equation
virtual void solvePoisson()
Solves the Poisson system.
virtual void generateBC2()
cusp::array1d< real, memoryType > temp1
temporary 1D Cusp array
domain * domInfo
computational grid information
virtual void writeCommon()
Writes numerical solution at current time-step, as well as the number of iterations performed in each...
cusp::array1d< real, memoryType > bc[4]
array that contains the boundary conditions of the rectangular domain
Definition of the class preconditioner.
integrationScheme intgSchm
instance of the class integrationScheme
Declaration of the class property.
virtual void writeData()
Writes data into files.
void initialiseBoundaryArrays()
Initializes boundary velocity arrays with values stored in the database.
virtual void solveIntermediateVelocity()
Solves for the intermediate flux velocity.
void updateBoundaryConditions()
Doing nothing.
virtual void projectionStep()
Projects the flux onto the divergence-free field.
void assembleMatrices()
Assembles matrices of the intermediate flux solver and the Poisson solver.
std::ofstream iterationsFile
file that contains the number of iterations
preconditioner< cusp::coo_matrix< int, real, memoryType > > * PC1
preconditioner for the intermediate flux solver
Definition of the class integrationScheme.
Definition of the class domain.
cusp::coo_matrix< int, real, memoryType > M
diagonal mass matrix
virtual void assembleRHS1()
Assembles the right hand-side of the system for the intermediate flux.
Logger logger
instance of the class Logger to track time of different tasks
cusp::coo_matrix< int, real, memoryType > A
combination of mass and Laplacian matrices
cusp::array1d< real, memoryType > bc2
inhomogeneous boundary conditions from the discrete continuity equation
void updateQ(real gamma)
Doing nothing.
Solves the Navier-Stokes equations in a rectangular domain.
cusp::coo_matrix< int, real, memoryType > QT
transposed gradient matrix (and interpolation matrix if immersed body)
virtual void updateSolverState()
Doing nothing. Used in immersed boundary methods when the body moves.
void assembleRHS2()
Assembles the right hand-side of the Poisson system.
virtual void initialise()
Initializes parameters, arrays and matrices required for the simulation.
virtual void generateQT()
preconditioner< cusp::coo_matrix< int, real, memoryType > > * PC2
preconditioner for the Poisson solver