cuIBM
A GPU-based Immersed Boundary Method code
DirectForcingSolver.h
Go to the documentation of this file.
1 
7 #pragma once
8 
9 #include "NSWithBody.h"
10 
11 
26 template <typename memoryType>
27 class DirectForcingSolver : public NSWithBody<memoryType>
28 {
29 protected:
30  // for the 1D interpolations
31  cusp::array1d<int, host_memory>
32  tags,
33  tags2;
34  cusp::array1d<int, device_memory>
36  tags2D;
38  coeffs2;
40  coeffs2D;
41 
42  vecH uv;
44  cusp::array1d<real, memoryType> pressure;
45 
46  // tag forcing points
47  void tagPoints();
48  void tagPoints(real *bx, real *by, real *uB, real *vB);
49 
50  // assemble the LHS matrix of the velocity system
51  virtual void generateA(real alpha);
52 
53  // assemble the Laplacian operator
54  virtual void generateL();
55 
56  // assemble the divergence operator
57  void generateQT(int *QTRows, int *QTCols, real *QTVals){} // check this!
58  virtual void generateQT();
59 
60  // update the gradient operator
61  void updateQ();
62 
63  // compute the Poisson matrix
64  virtual void generateC();
65 
66  // assemble the RHS of the velocity system
67  virtual void assembleRHS1();
68  // update the RHS of the velocity system
69  void updateRHS1();
70 
71  // update the solver
72  virtual void updateSolverState();
73 
74  // projection velocity onto divergence-free space
75  virtual void projectionStep();
76 
77  // write info about mass flux
78  void writeMassFluxInfo();
79 
80 public:
81  // constructor
82  DirectForcingSolver(parameterDB *pDB=NULL, domain *dInfo=NULL);
83 
84  // initialize the direct-forcing solver
85  virtual void initialise();
86 
87  // write data
88  virtual void writeData();
89 
93  virtual std::string name()
94  {
95  return "Direct Forcing";
96  }
97 
98 }; // DirectForcingSolver
cusp::array1d< int, host_memory > tags
vector used to differentiate forcing points from regular ones (host)
virtual void writeData()
Writes the velocity, pressure, force and mass flux data at every save point.
Declaration of the class NSWithBody.
double real
Is a float or a double depending on the machine precision.
Definition: types.h:116
cusp::array1d< int, device_memory > tagsD
vector used to differentiate forcing points from regular ones (device)
virtual void initialise()
Initialize the vectors used in the simulation.
virtual void generateQT()
Compute the divergence operator (transpose of the gradient one).
Definition: generateQT.inl:46
vecD uvD
velocity field (device)
vecH coeffs2
other coefficients of interpolation; quadratic interpolation (host)
void writeMassFluxInfo()
Prints the min, max and sum of the divergences of the velocity field in every cell of the domain...
A fully discrete formulation of the direct forcing method.
vecD coeffsD
coefficients of interpolation (device)
std::map< std::string, componentParameter > parameterDB
Map from a string to a componentParameter.
Definition: parameterDB.h:64
cusp::array1d< real, host_memory > vecH
Cusp 1D array stored in the host.
Definition: types.h:154
Stores information about the computational grid.
Definition: domain.h:16
DirectForcingSolver(parameterDB *pDB=NULL, domain *dInfo=NULL)
Constructor. Initializes the simulation parameters and the domain info.
cusp::array1d< real, memoryType > pressure
pressure field
virtual void generateA(real alpha)
vecD coeffs2D
other coefficients of interpolation; quadratic interpolation (device)
virtual void projectionStep()
Projects the pressure gradient on to the intermediate velocity field to obtain the divergence-free ve...
virtual void generateL()
cusp::array1d< int, host_memory > tags2
vector used to differentiate forcing points from regular ones (host)
void generateQT(int *QTRows, int *QTCols, real *QTVals)
cusp::array1d< int, device_memory > tags2D
vector used to differentiate forcing points from regular ones (device)
virtual void assembleRHS1()
Assembles the matrix rhs1 for DirectForcingSolver.
virtual std::string name()
Returns the name of the solver as a string.
vecH coeffs
coefficients of interpolation (host)
vecH uv
velocity field (host)
Generic Navier-Stokes solver in the presence of immersed boundaries.
Definition: NSWithBody.h:21
virtual void updateSolverState()
Updates the matrices every time the body is moved.
virtual void generateC()
Generates the right-hand side matrix in the Poisson step.
cusp::array1d< real, device_memory > vecD
Cusp 1D array stored in the device.
Definition: types.h:161