cuIBM
A GPU-based Immersed Boundary Method code
NavierStokesSolver.h
Go to the documentation of this file.
1 
7 #pragma once
8 
9 #include "utilities/types.h"
10 #include "utilities/domain.h"
12 #include "utilities/parameterDB.h"
14 #include "io/io.h"
15 
16 
25 template <typename memoryType>
27 {
28 protected:
32 
34 
35  cusp::coo_matrix<int, real, memoryType>
36  M,
37  Minv,
38  L,
39  A,
40  Q,
41  QT,
42  BN,
43  C;
44 
45  cusp::array1d<real, memoryType>
46  q,
47  qStar,
48  qOld,
49  lambda,
50  rn,
51  bc1,
52  bc2,
53  rhs1,
54  rhs2,
55  H,
56  temp1,
57  temp2,
58  bc[4];
59 
61  *PC1,
62  *PC2;
63 
64  size_t timeStep,
65  subStep,
68 
70 
71  std::ofstream iterationsFile;
72 
73  // initialize parameters common to all Navier-Stokes solvers
74  void initialiseCommon();
75 
76  // initialize all arrays required to solve the Navier-Stokes equations
77  void initialiseArrays(int numQ, int numLambda);
78 
79  // initialize velocity flux vectors
80  virtual void initialiseFluxes();
81 
82  // initialize velocity flux vectors
83  virtual void initialiseFluxes(real *q);
84 
85  // initialize boundary velocity arrays with values stored in the database
87 
88  // assemble matrices of the intermediate flux solver and the Poisson solver
89  void assembleMatrices();
90 
91  // generate the mass matrix and its inverse
92  void generateM();
93 
94  // generate the discrete Laplacian matrix
95  virtual void generateL();
96 
97  // generate the matrix resulting from the implicit flux terms
98  virtual void generateA(real alpha);
99 
100  // generate approximate inverse of the matrix resulting from implicit velocity terms
101  void generateBN();
102 
103  // generate the discrete divergence matrix
104  virtual void generateQT();
105 
106  // does nothing
107  void updateQ(real gamma);
108 
109  // generate the matrix of the Poisson solver
110  virtual void generateC();
111 
112  // generate explicit terms of the momemtum equation
113  void generateRN();
114 
115  // generate explicit flux terms
117 
118  // calculates the explicit lambda terms
119  virtual void calculateExplicitLambdaTerms();
120 
121  // generate inhomogeneous boundary conditions from the discrete Laplacian operator
122  virtual void generateBC1();
123 
124  // generate inhomogeneous boundary conditions from the discrete continuity equation
125  virtual void generateBC2();
126 
127  // assemble the right hand-side of the system for the intermediate flux
128  virtual void assembleRHS1();
129 
130  // assemble the right hand-side of the Poisson system.
131  void assembleRHS2();
132 
133  // solve for the intermediate flux velocity
134  virtual void solveIntermediateVelocity();
135 
136  // solve the Poisson system for the pressure (and the body forces if immersed body)
137  virtual void solvePoisson();
138 
139  // project the flux onto the divergence-free field
140  virtual void projectionStep();
141 
142  // doing nothing
144 
145  // doing nothing
146  virtual void updateSolverState();
147 
148 public:
149  // constructor -- copy the database and information about the computational grid
150  NavierStokesSolver(parameterDB *pDB=NULL, domain *dInfo=NULL);
151 
152  // initialize parameters, arrays and matrices required for the simulation
153  virtual void initialise();
154 
155  // calculate the variables at the next time step
156  void stepTime();
157 
158  // write numerical solution and number of iterations performed in each solver.
159  virtual void writeCommon();
160 
161  // write data into files
162  virtual void writeData();
163 
164  // evaluate the condition required to stop the simulation
165  bool finished();
166 
167  // print timing information and clse the different files
168  virtual void shutDown();
169 
173  virtual std::string name()
174  {
175  return "Navier-Stokes";
176  }
177 };
Monitors the time spent to achieve a certain task.
Definition: logger.h:33
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.
Definition: types.h:116
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...
Definition: generateRN.inl:19
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.
Definition: parameterDB.h:64
cusp::array1d< real, memoryType > q
velocity flux vector
virtual void generateL()
void generateRN()
Generates explicit terms of the momentum equation.
Definition: generateRN.inl:89
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.
Definition: domain.h:16
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 generateC()
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