cuIBM
A GPU-based Immersed Boundary Method code
Public Member Functions | Private Member Functions | Private Attributes | List of all members
TairaColoniusSolver< memoryType > Class Template Reference

Solves the Navier-Stokes equations using the Immersed boundary method from Taira and Colonius (2007). More...

#include <TairaColoniusSolver.h>

Inheritance diagram for TairaColoniusSolver< memoryType >:
Inheritance graph
[legend]
Collaboration diagram for TairaColoniusSolver< memoryType >:
Collaboration graph
[legend]

Public Member Functions

 TairaColoniusSolver (parameterDB *pDB=NULL, domain *dInfo=NULL)
 Constructor. Copies the database and information about the computational grid. More...
 
virtual void initialise ()
 Initializes the solvers. More...
 
virtual void writeData ()
 Calculates and writes forces acting on each immersed body at current time. More...
 
virtual std::string name ()
 Returns the name of the solver. More...
 
- Public Member Functions inherited from NSWithBody< memoryType >
virtual void writeCommon ()
 Writes flow variables and position of body points into files. More...
 
virtual void shutDown ()
 Closes iteration file and force file. More...
 
- Public Member Functions inherited from NavierStokesSolver< memoryType >
 NavierStokesSolver (parameterDB *pDB=NULL, domain *dInfo=NULL)
 Constructor. Copies the database and information about the computational grid. More...
 
void stepTime ()
 Calculates the variables at the next time step. More...
 
bool finished ()
 Evaluates the condition required to stop the simulation. More...
 

Private Member Functions

virtual void generateQT ()
 
void updateQT ()
 
virtual void generateBC2 ()
 
virtual void updateSolverState ()
 Updates the location of the bodies and re-generates appropriate matrices. More...
 
virtual void calculateForce ()
 Calculates forces acting on each immersed body. More...
 
template<>
void generateBC2 ()
 Generates the right hand-side of the Poisson system. More...
 
template<>
void updateQT ()
 Updates the interpolation matrix using the current locations of body points. More...
 
template<>
void generateQT ()
 Generates the transposed gradient matrix and interpolation matrix. More...
 

Private Attributes

cusp::coo_matrix< int, real, memoryType > E
 Interpolation matrix from the Eulerian grid to the Lagrangian points. More...
 
cusp::coo_matrix< int, real, memoryType > ET
 Regularization matrix form the Lagrangian points to the Eulerian grid. More...
 

Additional Inherited Members

- Protected Member Functions inherited from NSWithBody< memoryType >
void initialiseBodies ()
 Stores the parameters of the simulation and initializes the location and motion of each immersed bodies. More...
 
void updateBodies ()
 Updates location and motion of each immersed body at current time. More...
 
template<>
void calculateForce ()
 Calculates forces acting on an immersed body. More...
 
- Protected Member Functions inherited from NavierStokesSolver< memoryType >
void initialiseCommon ()
 Initializes parameters common to all Navier-Stokes solvers. More...
 
void initialiseArrays (int numQ, int numLambda)
 Initializes all arrays required to solve the Navier-Stokes equations. More...
 
virtual void initialiseFluxes ()
 
virtual void initialiseFluxes (real *q)
 Initializes velocity flux vectors. More...
 
void initialiseBoundaryArrays ()
 Initializes boundary velocity arrays with values stored in the database. More...
 
void assembleMatrices ()
 Assembles matrices of the intermediate flux solver and the Poisson solver. More...
 
void generateM ()
 
virtual void generateL ()
 
virtual void generateA (real alpha)
 
void generateBN ()
 Generates approximate inverse of the matrix resulting from implicit velocity terms. More...
 
void updateQ (real gamma)
 Doing nothing. More...
 
virtual void generateC ()
 
void generateRN ()
 Generates explicit terms of the momentum equation. More...
 
void calculateExplicitQTerms ()
 
virtual void calculateExplicitLambdaTerms ()
 Doing nothing. Used in methods that use the explicit pressure term in the intermediate velocity solve step, such as FadlunEtAlSolver and DFModifiedSolver. More...
 
virtual void generateBC1 ()
 
virtual void assembleRHS1 ()
 Assembles the right hand-side of the system for the intermediate flux. More...
 
void assembleRHS2 ()
 Assembles the right hand-side of the Poisson system. More...
 
virtual void solveIntermediateVelocity ()
 Solves for the intermediate flux velocity. More...
 
virtual void solvePoisson ()
 Solves the Poisson system. More...
 
virtual void projectionStep ()
 Projects the flux onto the divergence-free field. More...
 
void updateBoundaryConditions ()
 Doing nothing. More...
 
template<>
void generateA (real alpha)
 Assembles the matrix from implicit terms of momentum equation. More...
 
template<>
void generateBC1 ()
 Generates inhomogeneous boundary conditions from the discrete Laplacian operator. More...
 
template<>
void generateBC2 ()
 Generates inhomogeneous boundary conditions from the discrete continuity equation. More...
 
template<>
void generateL ()
 Generates the discrete Laplacian matrix. More...
 
template<>
void generateM ()
 Generates the mass matrix and its inverse. More...
 
template<>
void generateQT ()
 Generates the discrete divergence matrix. More...
 
template<>
void calculateExplicitQTerms ()
 Generates explicit terms that arise from the velocity fluxes. More...
 
template<>
void initialiseFluxes ()
 Initializes velocity flux vectors (on the device). More...
 
template<>
void generateC ()
 Generates the matrix of the Poisson solver. More...
 
- Protected Attributes inherited from NSWithBody< memoryType >
bodies< memoryType > B
 bodies in the flow More...
 
real forceX
 force acting on each body in the x-direction More...
 
real forceY
 force acting on each body in the y-direction More...
 
std::ofstream forceFile
 file to write the forces More...
 
- Protected Attributes inherited from NavierStokesSolver< memoryType >
parameterDBparamDB
 database that contains all the simulation parameters More...
 
domaindomInfo
 computational grid information More...
 
integrationScheme intgSchm
 instance of the class integrationScheme More...
 
real QCoeff
 
cusp::coo_matrix< int, real, memoryType > M
 diagonal mass matrix More...
 
cusp::coo_matrix< int, real, memoryType > Minv
 inverse of the mass matrix More...
 
cusp::coo_matrix< int, real, memoryType > L
 discrete Laplacian matrix More...
 
cusp::coo_matrix< int, real, memoryType > A
 combination of mass and Laplacian matrices More...
 
cusp::coo_matrix< int, real, memoryType > Q
 gradient matrix (and regularization matrix if immersed body) More...
 
cusp::coo_matrix< int, real, memoryType > QT
 transposed gradient matrix (and interpolation matrix if immersed body) More...
 
cusp::coo_matrix< int, real, memoryType > BN
 N-th order Taylor series expansion of the inverse of A. More...
 
cusp::coo_matrix< int, real, memoryType > C
 matrix of the Poisson system More...
 
cusp::array1d< real, memoryType > q
 velocity flux vector More...
 
cusp::array1d< real, memoryType > qStar
 intermediate velocity flux vector More...
 
cusp::array1d< real, memoryType > qOld
 velocity flux vector at the previous time-step More...
 
cusp::array1d< real, memoryType > lambda
 pressure vector (and body forces if immersed body) More...
 
cusp::array1d< real, memoryType > rn
 explicit terms of the momentum equation More...
 
cusp::array1d< real, memoryType > bc1
 inhomogeneous boundary conditions from the discrete Laplacian operator More...
 
cusp::array1d< real, memoryType > bc2
 inhomogeneous boundary conditions from the discrete continuity equation More...
 
cusp::array1d< real, memoryType > rhs1
 right hand-side of the system for the intermediate velocity flux vector More...
 
cusp::array1d< real, memoryType > rhs2
 right hand-side of the Poisson system More...
 
cusp::array1d< real, memoryType > H
 hold convective terms from previous time step More...
 
cusp::array1d< real, memoryType > temp1
 temporary 1D Cusp array More...
 
cusp::array1d< real, memoryType > temp2
 temporary 1D Cusp array More...
 
cusp::array1d< real, memoryType > bc [4]
 array that contains the boundary conditions of the rectangular domain More...
 
preconditioner< cusp::coo_matrix< int, real, memoryType > > * PC1
 preconditioner for the intermediate flux solver More...
 
preconditioner< cusp::coo_matrix< int, real, memoryType > > * PC2
 preconditioner for the Poisson solver More...
 
size_t timeStep
 iteration number More...
 
size_t subStep
 number of sub-iterations More...
 
size_t iterationCount1
 number of iteration to solve the intermediate fluxes More...
 
size_t iterationCount2
 number of iteration to solve the Poisson equation More...
 
Logger logger
 instance of the class Logger to track time of different tasks More...
 
std::ofstream iterationsFile
 file that contains the number of iterations More...
 

Detailed Description

template<typename memoryType>
class TairaColoniusSolver< memoryType >

Solves the Navier-Stokes equations using the Immersed boundary method from Taira and Colonius (2007).

The immersed boundary method: a projection approach
Taira K. and Colonius T.
Journal of Computational Physics
Volume 225 Number 2 (2007).

Definition at line 24 of file TairaColoniusSolver.h.

Constructor & Destructor Documentation

template<typename memoryType >
TairaColoniusSolver< memoryType >::TairaColoniusSolver ( parameterDB pDB = NULL,
domain dInfo = NULL 
)

Constructor. Copies the database and information about the computational grid.

Definition at line 16 of file TairaColoniusSolver.cu.

Member Function Documentation

template<typename memoryType >
void TairaColoniusSolver< memoryType >::calculateForce ( )
privatevirtual

Calculates forces acting on each immersed body.

For each immersed body, it spreads the Lagrangian forces on the Eulerian grid and sum the entire field to get the x- and y- components of the force.

Reimplemented from NSWithBody< memoryType >.

Definition at line 19 of file calculateForce.inl.

References NSWithBody< memoryType >::forceX, and NSWithBody< memoryType >::forceY.

template<>
void TairaColoniusSolver< device_memory >::generateBC2 ( )
privatevirtual

Generates the right hand-side of the Poisson system.

It contains the inhomogeneous boundary conditions from the discrete divergence operator, as well as the no-slip boundary condition at the body surface.

Reimplemented from NavierStokesSolver< memoryType >.

Definition at line 19 of file generateBC2.inl.

References XMINUS, XPLUS, YMINUS, and YPLUS.

template<typename memoryType>
virtual void TairaColoniusSolver< memoryType >::generateBC2 ( )
privatevirtual

Reimplemented from NavierStokesSolver< memoryType >.

template<typename memoryType>
virtual void TairaColoniusSolver< memoryType >::generateQT ( )
privatevirtual

Reimplemented from NavierStokesSolver< memoryType >.

template<>
void TairaColoniusSolver< device_memory >::generateQT ( )
privatevirtual

Generates the transposed gradient matrix and interpolation matrix.

QT is an (np + 2*nb) x nuv matrix

Reimplemented from NavierStokesSolver< memoryType >.

Definition at line 70 of file generateQT.inl.

References kernels::generateQT(), and kernels::updateQT().

template<typename memoryType >
void TairaColoniusSolver< memoryType >::initialise ( )
virtual

Initializes the solvers.

Initializes bodies, arrays and matrices of the intermediate flux solver and the Poisson solver.

Reimplemented from NavierStokesSolver< memoryType >.

Definition at line 31 of file TairaColoniusSolver.cu.

References NavierStokesSolver< memoryType >::assembleMatrices(), NavierStokesSolver< memoryType >::initialiseArrays(), NSWithBody< memoryType >::initialiseBodies(), and NavierStokesSolver< memoryType >::initialiseCommon().

template<typename memoryType>
virtual std::string TairaColoniusSolver< memoryType >::name ( )
inlinevirtual

Returns the name of the solver.

Reimplemented from NavierStokesSolver< memoryType >.

Definition at line 59 of file TairaColoniusSolver.h.

template<>
void TairaColoniusSolver< device_memory >::updateQT ( )
private

Updates the interpolation matrix using the current locations of body points.

Typically called after the body has moved.

Definition at line 20 of file generateQT.inl.

References BLOCKSIZE.

template<typename memoryType>
void TairaColoniusSolver< memoryType >::updateQT ( )
private
template<typename memoryType >
void TairaColoniusSolver< memoryType >::updateSolverState ( )
privatevirtual

Updates the location of the bodies and re-generates appropriate matrices.

Updates only done if the bodies are moving. It re-generates the interpolation matrix, therefore the Poisson matrix too.

Reimplemented from NavierStokesSolver< memoryType >.

Definition at line 98 of file TairaColoniusSolver.cu.

References NavierStokesSolver< memoryType >::generateC(), NSWithBody< memoryType >::updateBodies(), and kernels::updateQT().

template<typename memoryType >
void TairaColoniusSolver< memoryType >::writeData ( )
virtual

Calculates and writes forces acting on each immersed body at current time.

Reimplemented from NavierStokesSolver< memoryType >.

Definition at line 61 of file TairaColoniusSolver.cu.

References NSWithBody< memoryType >::forceX, NSWithBody< memoryType >::forceY, and NSWithBody< memoryType >::writeCommon().

Member Data Documentation

template<typename memoryType>
cusp::coo_matrix<int, real, memoryType> TairaColoniusSolver< memoryType >::E
private

Interpolation matrix from the Eulerian grid to the Lagrangian points.

Definition at line 28 of file TairaColoniusSolver.h.

template<typename memoryType>
cusp::coo_matrix<int, real, memoryType> TairaColoniusSolver< memoryType >::ET
private

Regularization matrix form the Lagrangian points to the Eulerian grid.

Definition at line 28 of file TairaColoniusSolver.h.


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