OpenLB 1.7
Loading...
Searching...
No Matches
porousForcedBGKDynamics.h
Go to the documentation of this file.
1/* This file is part of the OpenLB library
2 *
3 * Copyright (C) 2012 Asher Zarth, Lukas Baron, Mathias J. Krause, Jonas Latt, Jan E. Marquardt
4 * E-mail contact: info@openlb.net
5 * The most recent release of OpenLB can be downloaded at
6 * <http://www.openlb.net/>
7 *
8 * This program is free software; you can redistribute it and/or
9 * modify it under the terms of the GNU General Public License
10 * as published by the Free Software Foundation; either version 2
11 * of the License, or (at your option) any later version.
12 *
13 * This program is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public
19 * License along with this program; if not, write to the Free
20 * Software Foundation, Inc., 51 Franklin Street, Fifth Floor,
21 * Boston, MA 02110-1301, USA.
22*/
23
27#ifndef POROUS_FORCED_BGK_DYNAMICS_H
28#define POROUS_FORCED_BGK_DYNAMICS_H
29
30#include "dynamics/dynamics.h"
31#include "core/cell.h"
32
33namespace olb {
34
36template<typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
38 T,DESCRIPTOR,
43>;
44
46
47/* Implementations of the BGK collision for moving porous media (HLBM approach).
48 * As this scheme requires additionla data stored in an external field,
49 * it is meant to be used along with a PorousParticle descriptor.
50 * \param omega Lattice relaxation frequency
51 * \param momenta A standard object for the momenta computation
52 */
54template<typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
56 T, DESCRIPTOR,
61>;
62
64template<typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
66 T, DESCRIPTOR,
71>;
72
74template<typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
76 T, DESCRIPTOR,
81>;
82
84template<typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
86 T, DESCRIPTOR,
91>;
92
93namespace forcing {
94
95template <bool isStatic = false>
97 static std::string getName() {
98 return "PorousParticleKupershtokhForcing";
99 }
100
101 template <typename DESCRIPTOR, typename MOMENTA>
102 using combined_momenta = typename momenta::Forced<momenta::PorousParticle<MOMENTA>>::template type<DESCRIPTOR>;
103
104 template <typename DESCRIPTOR, typename MOMENTA, typename EQUILIBRIUM>
105 using combined_equilibrium = typename EQUILIBRIUM::template type<DESCRIPTOR,momenta::PorousParticle<MOMENTA>>;
106
107 template <typename DESCRIPTOR, typename MOMENTA, typename EQUILIBRIUM, typename COLLISION>
109 using MomentaF = typename MOMENTA::template type<DESCRIPTOR>;
111 using CollisionO = typename COLLISION::template type<DESCRIPTOR,MOMENTA,EQUILIBRIUM>;
112
113 constexpr static bool is_vectorizable = false;
114
115 template <typename CELL, typename VELOCITY, typename V=typename CELL::value_t>
116 void calculate(CELL& cell, VELOCITY& u) {
117 if constexpr (isStatic) {
118 for (int i=0; i<DESCRIPTOR::d; ++i) {
119 u[i] -= (V{1} - cell.template getField<descriptors::POROSITY>()) * u[i];
120 }
121 } else {
122 for (int i=0; i<DESCRIPTOR::d; ++i) {
123 u[i] += (V{1} - cell.template getField<descriptors::POROSITY>())
124 * ( cell.template getFieldComponent<descriptors::VELOCITY_NUMERATOR>(i)
125 / cell.template getField<descriptors::VELOCITY_DENOMINATOR>() - u[i]);
126 }
127 }
128 }
129
130 template <typename CELL, typename PARAMETERS, typename V=typename CELL::value_t>
131 CellStatistic<V> apply(CELL& cell, PARAMETERS& parameters) any_platform {
132 V u[DESCRIPTOR::d];
133 MomentaF().computeU(cell, u);
134 const auto velDenominator = cell.template getField<descriptors::VELOCITY_DENOMINATOR>();
135 const auto statistic = CollisionO().apply(cell, parameters);
136 const auto force = cell.template getField<descriptors::FORCE>();
137 V uPlusDeltaU[DESCRIPTOR::d];
138 for (int iVel=0; iVel < DESCRIPTOR::d; ++iVel) {
139 uPlusDeltaU[iVel] = u[iVel] + force[iVel];
140 }
141 if (velDenominator > std::numeric_limits<V>::epsilon()) {
142 calculate(cell, uPlusDeltaU);
143 if constexpr (!isStatic) {
144 // reset external field for next timestep
145 cell.template setField<descriptors::POROSITY>(1.);
146 cell.template setField<descriptors::VELOCITY_DENOMINATOR>(0.);
147 cell.template setField<descriptors::VELOCITY_NUMERATOR>({0.,0.,0.});
148 }
149 }
150 for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
151 cell[iPop] += EquilibriumF().compute(iPop, statistic.rho, uPlusDeltaU)
152 - EquilibriumF().compute(iPop, statistic.rho, u);
153 }
154
155 // Reset contact helper if utilized
156 if constexpr (DESCRIPTOR::template provides<descriptors::CONTACT_DETECTION>()) {
157 cell.template setField<descriptors::CONTACT_DETECTION>(0);
158 }
159
160 return statistic;
161 };
162 };
163
164 template <typename DESCRIPTOR, typename MOMENTA, typename EQUILIBRIUM, typename COLLISION>
165 using combined_parameters = typename COLLISION::parameters;
166
167};
168
169}
170
172template<typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
174 T, DESCRIPTOR,
175 MOMENTA,
179>;
180
182template<typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
184 T, DESCRIPTOR,
185 MOMENTA,
189>;
190
191} // olb
192
193#endif
Definition of a LB cell – header file.
Dynamics combination rule implementing the forcing scheme by Shan and Chen.
Definition forcing.h:188
Top level namespace for all of OpenLB.
#define any_platform
Define preprocessor macros for device-side functions, constant storage.
Definition platform.h:78
Return value of any collision.
Definition interface.h:43
Dynamics constructed as a tuple of momenta, equilibrium and collision.
Definition interface.h:182
Dynamics combination rule implementing the forcing scheme by Guo et al.
Definition forcing.h:38
typename COLLISION::template type< DESCRIPTOR, MOMENTA, EQUILIBRIUM > CollisionO
combined_equilibrium< DESCRIPTOR, MOMENTA, EQUILIBRIUM > EquilibriumF
CellStatistic< V > apply(CELL &cell, PARAMETERS &parameters) any_platform
typename EQUILIBRIUM::template type< DESCRIPTOR, momenta::PorousParticle< MOMENTA > > combined_equilibrium
typename COLLISION::parameters combined_parameters
typename momenta::Forced< momenta::PorousParticle< MOMENTA > >::template type< DESCRIPTOR > combined_momenta