cuIBM
A GPU-based Immersed Boundary Method code
integrationScheme.h
Go to the documentation of this file.
1 
7 #pragma once
8 
9 #include "types.h"
10 
11 
17 {
18 public:
19  int subSteps;
20 
22  zeta,
25 
26  // Convection term = gamma*H(n) + zeta*H(n-1)
27  // DIffusion term = alphaImplicit*D(n+1) + alphaExplicit*D(n)
28 
29 
36  void initialise(timeScheme convScheme, timeScheme diffScheme)
37  {
38  // set the number of substeps required for the specified scheme
39  switch(convScheme)
40  {
41  case EULER_EXPLICIT:
42  subSteps = 1;
43  break;
44  case ADAMS_BASHFORTH_2:
45  subSteps = 1;
46  break;
47  case RUNGE_KUTTA_3:
48  subSteps = 3;
49  break;
50  default:
51  break;
52  }
53 
54  // resize the arrays with the number of sub-iterations
55  gamma.resize(subSteps);
56  zeta.resize(subSteps);
57  alphaExplicit.resize(subSteps);
58  alphaImplicit.resize(subSteps);
59 
60  // set the coefficients of the convection terms
61  switch(convScheme)
62  {
63  case EULER_EXPLICIT:
64  gamma[0] = 1;
65  zeta[0] = 0;
66  break;
67  case ADAMS_BASHFORTH_2:
68  gamma[0] = 1.5;
69  zeta[0] = -0.5;
70  break;
71  case RUNGE_KUTTA_3:
72  gamma[0] = 8.0/15;
73  zeta[0] = 0;
74  gamma[1] = 5.0/12;
75  zeta[1] = -17.0/60;
76  gamma[2] = 3.0/4;
77  zeta[2] = -5.0/12;
78  break;
79  default:
80  break;
81  }
82 
83  // set the coefficients of the diffusion terms
84  real aI, aE;
85  switch(diffScheme)
86  {
87  case EULER_EXPLICIT:
88  aI = 0.0;
89  aE = 1.0;
90  break;
91  case EULER_IMPLICIT:
92  aI = 1.0;
93  aE = 0.0;
94  break;
95  case CRANK_NICOLSON:
96  aI = 0.5;
97  aE = 0.5;
98  default:
99  break;
100  }
101  for(int i=0; i<subSteps; i++)
102  {
103  alphaExplicit[i] = aE*(gamma[i]+zeta[i]);
104  alphaImplicit[i] = aI*(gamma[i]+zeta[i]);
105  }
106  } // initialise
107 
108 }; // integrationScheme
Specifies the time-integration scheme used.
Definition of custom types required by the code.
double real
Is a float or a double depending on the machine precision.
Definition: types.h:116
vecH zeta
coefficient of the convection term in the previous time step
vecH alphaImplicit
coefficient of the implicit diffusion term
implicit Euler method (first order)
Definition: types.h:63
timeScheme
Specifies the numerical scheme used for time-integration.
Definition: types.h:60
cusp::array1d< real, host_memory > vecH
Cusp 1D array stored in the host.
Definition: types.h:154
vecH alphaExplicit
coefficient of the explicit diffusion term
second-order Adams-Bashforth scheme
Definition: types.h:64
third-order low storage Runge-Kutta method
Definition: types.h:65
Crank-Nicolson scheme (second order)
Definition: types.h:66
int subSteps
number of substeps inside each time step
void initialise(timeScheme convScheme, timeScheme diffScheme)
Initializes the coefficients of the time-integration scheme.
vecH gamma
coefficient of the convection term in the current time step
explicit Euler method (first order)
Definition: types.h:62