summaryrefslogtreecommitdiff
path: root/src/dynamics/porousBGKdynamics.h
blob: f0f33470ca77c297575e9230d8bce056e6ea3ab0 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
/*  This file is part of the OpenLB library
 *
 *  Copyright (C) 2012 Lukas Baron, Mathias J. Krause, Jonas Latt
 *  E-mail contact: info@openlb.net
 *  The most recent release of OpenLB can be downloaded at
 *  <http://www.openlb.net/>
 *
 *  This program is free software; you can redistribute it and/or
 *  modify it under the terms of the GNU General Public License
 *  as published by the Free Software Foundation; either version 2
 *  of the License, or (at your option) any later version.
 *
 *  This program is distributed in the hope that it will be useful,
 *  but WITHOUT ANY WARRANTY; without even the implied warranty of
 *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *  GNU General Public License for more details.
 *
 *  You should have received a copy of the GNU General Public
 *  License along with this program; if not, write to the Free
 *  Software Foundation, Inc., 51 Franklin Street, Fifth Floor,
 *  Boston, MA  02110-1301, USA.
*/

/** \file
 * BGK Dynamics for porous media -- header file.
 */
#ifndef POROUS_BGK_DYNAMICS_H
#define POROUS_BGK_DYNAMICS_H

#include "dynamics/dynamics.h"
#include "core/cell.h"

namespace olb {

/// Implementation of the BGK collision step for a porosity model
template<typename T, typename DESCRIPTOR>
class PorousBGKdynamics : public BGKdynamics<T,DESCRIPTOR> {
public:
  /// Constructor
  PorousBGKdynamics(T omega_, Momenta<T,DESCRIPTOR>& momenta_);

  /// Collision step
  virtual void collide(Cell<T,DESCRIPTOR>& cell,
                       LatticeStatistics<T>& statistics_);

  /// get relaxation parameter
  T    getOmega() const;
  /// set relaxation parameter
  void setOmega(T omega_);


private:
  T omega;      ///< relaxation parameter
};

/// Implementation of the BGK collision step for a porosity model enabling
/// drag computation
template<typename T, typename DESCRIPTOR>
class ExtendedPorousBGKdynamics : public BGKdynamics<T,DESCRIPTOR> {
public:
  /// Constructor
  ExtendedPorousBGKdynamics(T omega_, Momenta<T,DESCRIPTOR>& momenta_);
  /// extended Collision step, computes local drag in a given direction
  virtual void collide(Cell<T,DESCRIPTOR>& cell,
                       LatticeStatistics<T>& statistics_);
  /// get relaxation parameter
  T    getOmega() const;
  /// set relaxation parameter
  void setOmega(T omega_);


private:
  T omega;      ///< relaxation parameter
};

/// Implementation of the BGK collision step for subgridscale particles
template<typename T, typename DESCRIPTOR>
class SubgridParticleBGKdynamics : public BGKdynamics<T,DESCRIPTOR> {
public:
  /// Constructor
  SubgridParticleBGKdynamics(T omega_, Momenta<T,DESCRIPTOR>& momenta_);
  /// extended Collision step, computes local drag in a given direction
  virtual void collide(Cell<T,DESCRIPTOR>& cell,
                       LatticeStatistics<T>& statistics_);
  /// get relaxation parameter
  T    getOmega() const;
  /// set relaxation parameter
  void setOmega(T omega_);


private:
  T omega;      ///< relaxation parameter
  T _fieldTmp[4];
};

/* Implementation of the BGK collision for moving porous media (HLBM approach).
 * As this scheme requires additionla data stored in an external field, 
 * it is meant to be used along with a PorousParticle descriptor.
 * \param omega Lattice relaxation frequency
 * \param momenta A standard object for the momenta computation
 */
template<typename T, typename DESCRIPTOR, bool isStatic=false>
class PorousParticleBGKdynamics : public BGKdynamics<T,DESCRIPTOR> {
public:
  /// Constructor
  PorousParticleBGKdynamics(T omega_, Momenta<T,DESCRIPTOR>& momenta_);
  /// extended Collision step, computes local drag in a given direction
  virtual void collide(Cell<T,DESCRIPTOR>& cell,
                       LatticeStatistics<T>& statistics_);
  /// get relaxation parameter
  T    getOmega() const;
  /// set relaxation parameter
  void setOmega(T omega_);


private:
  T omega;      ///< relaxation parameter
  /// This structure is used to emulate a "static if" to switch between static and
  /// dynamic case. It can be replaced by a "constexpr if" when switching to C++17
  /// standard.
  template<bool isStaticStruct, bool dummy = true> struct effectiveVelocity {
      static void calculate(T* pExternal, T* pVelocity);
  };
};

/// Implementation of the HBGK collision step for a porosity model enabling
/// drag computation for many particles
/// including the Krause turbulence modell

template<typename T, typename DESCRIPTOR>
class KrauseHBGKdynamics : public BGKdynamics<T,DESCRIPTOR> {
public:
  /// Constructor
  KrauseHBGKdynamics(T omega_, Momenta<T,DESCRIPTOR>& momenta_, T smagoConst_, T dx_ = 1, T dt_ = 1);
  /// extended Collision step, computes local drag in a given direction
  virtual void collide(Cell<T,DESCRIPTOR>& cell,
                       LatticeStatistics<T>& statistics_);
  /// get relaxation parameter
  T    getOmega() const;
  /// set relaxation parameter
  void setOmega(T omega_);


private:
  /// Computes a constant prefactor in order to speed up the computation
  T computePreFactor(T omega_, T smagoConst_);
  /// Computes the local smagorinsky relaxation parameter
  void computeOmega(T omega0_, Cell<T,DESCRIPTOR>& cell, T preFactor_, T rho_,
                    T u[DESCRIPTOR::d],
                    T newOmega[DESCRIPTOR::q] );

private:

  T omega;      ///< relaxation parameter
  /// effective collision time based upon Smagorisnky approach
  T tau_eff;
  /// Smagorinsky constant
  T smagoConst;
  /// Precomputed constant which speeeds up the computation
  T preFactor;
  T dx;
  T dt;

  T _fieldTmp[4];

};

/// Implementation of the BGK collision step for a porosity model enabling
/// drag computation
template<typename T, typename DESCRIPTOR>
class ParticlePorousBGKdynamics : public BGKdynamics<T,DESCRIPTOR> {
public:
  /// Constructor
  ParticlePorousBGKdynamics(T omega_, Momenta<T,DESCRIPTOR>& momenta_);
  /// extended Collision step, computes local drag in a given direction
  virtual void collide(Cell<T,DESCRIPTOR>& cell,
                       LatticeStatistics<T>& statistics_);
  /// get relaxation parameter
  T    getOmega() const;
  /// set relaxation parameter
  void setOmega(T omega_);


private:
  T omega;      ///< relaxation parameter
};

/// Implementation of the BGK collision step for a small particles enabling
/// two way coupling
template<typename T, typename DESCRIPTOR>
class SmallParticleBGKdynamics : public BGKdynamics<T,DESCRIPTOR> {
public:
  /// Constructor
  SmallParticleBGKdynamics(T omega_, Momenta<T,DESCRIPTOR>& momenta_);
  /// extended Collision step, computes local drag in a given direction
  virtual void collide(Cell<T,DESCRIPTOR>& cell,
                       LatticeStatistics<T>& statistics_);
  /// get relaxation parameter
  T    getOmega() const;
  /// set relaxation parameter
  void setOmega(T omega_);


private:
  T omega;      ///< relaxation parameter
};


enum Mode {M2, M3};

/// Implementation of the Partially Saturated Method (PSM),
/// see Krüger, Timm, et al. The Lattice Boltzmann Method. Springer, 2017. (p.447-451)
template<typename T, typename DESCRIPTOR>
class PSMBGKdynamics : public BGKdynamics<T,DESCRIPTOR> {
public:
  /// Constructor
  PSMBGKdynamics(T omega_, Momenta<T,DESCRIPTOR>& momenta_, int mode_=0);
  ///  Compute fluid velocity on the cell.
  void computeU (
    Cell<T,DESCRIPTOR> const& cell,
    T u[DESCRIPTOR::d] ) const override;
  /// Compute fluid velocity and particle density on the cell.
  void computeRhoU (
    Cell<T,DESCRIPTOR> const& cell,
    T& rho, T u[DESCRIPTOR::d]) const override;
  /// Collision step
  void collide(Cell<T,DESCRIPTOR>& cell,
                       LatticeStatistics<T>& statistics_) override;
  /// get relaxation parameter
  T    getOmega() const;
  /// set relaxation parameter
  void setOmega(T omega_);


private:
  T omega;      ///< relaxation parameter
  T paramA;      /// speed up parameter
  Mode mode;
};

/// Implementation of the Partially Saturated Method (PSM),
/// see Krüger, Timm, et al. The Lattice Boltzmann Method. Springer, 2017. (p.447-451)
template<typename T, typename DESCRIPTOR>
class ForcedPSMBGKdynamics : public ForcedBGKdynamics<T,DESCRIPTOR> {
public:
  /// Constructor
  ForcedPSMBGKdynamics(T omega_, Momenta<T,DESCRIPTOR>& momenta_, int mode_=0);
  ///  Compute fluid velocity on the cell.
  void computeU (
    Cell<T,DESCRIPTOR> const& cell,
    T u[DESCRIPTOR::d] ) const override;
  /// Compute fluid velocity and particle density on the cell.
  void computeRhoU (
    Cell<T,DESCRIPTOR> const& cell,
    T& rho, T u[DESCRIPTOR::d]) const override;
  /// Collision step
  void collide(Cell<T,DESCRIPTOR>& cell,
                       LatticeStatistics<T>& statistics_) override;
  /// get relaxation parameter
  T    getOmega() const;
  /// set relaxation parameter
  void setOmega(T omega_);


private:
  T omega;      ///< relaxation parameter
  T paramA;      /// speed up parameter
  Mode mode;
};

} // olb

#endif