OpenLB 1.8.1
Loading...
Searching...
No Matches
dynamics.h
Go to the documentation of this file.
1/* This file is part of the OpenLB library
2 *
3 * Copyright (C) 2006, 2007 Jonas Latt, Mathias J. Krause
4 * 2021 Adrian Kummerlaender
5 * E-mail contact: info@openlb.net
6 * The most recent release of OpenLB can be downloaded at
7 * <http://www.openlb.net/>
8 *
9 * This program is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU General Public License
11 * as published by the Free Software Foundation; either version 2
12 * of the License, or (at your option) any later version.
13 *
14 * This program is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 * GNU General Public License for more details.
18 *
19 * You should have received a copy of the GNU General Public
20 * License along with this program; if not, write to the Free
21 * Software Foundation, Inc., 51 Franklin Street, Fifth Floor,
22 * Boston, MA 02110-1301, USA.
23*/
24
25#ifndef DYNAMICS_DYNAMICS_H
26#define DYNAMICS_DYNAMICS_H
27
28#include "interface.h"
29#include "modifiers.h"
30
31#include "core/util.h"
34
35#include "momenta/interface.h"
36#include "momenta/aliases.h"
37
38#include "collision.h"
39#include "equilibrium.h"
40#include "forcing.h"
41
42#include "collisionModifiers.h"
43
44// TODO: Check this later
45#include "collisionLES.h"
46
47namespace olb {
48
50template <typename T, typename DESCRIPTOR>
52 T, DESCRIPTOR,
53 // Return 1 for the 0th moment, 0 for all others
59 >,
62>;
63
65template <typename T, typename DESCRIPTOR>
67 T, DESCRIPTOR,
71>;
72
74template <typename T, typename DESCRIPTOR>
76 T, DESCRIPTOR,
77 // Return 1 for the 0th moment, 0 for all others
83 >,
86>;
87
89template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
91 T, DESCRIPTOR,
92 MOMENTA,
93 std::conditional_t< (DESCRIPTOR::d == 3 && DESCRIPTOR::q == 7)
94 || (DESCRIPTOR::d == 2 && DESCRIPTOR::q == 5),
98>;
99
100//BGK collision with third order equilibrium function
101template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
103 T, DESCRIPTOR,
104 MOMENTA,
107>;
108
109// HRR collision https://hal.science/hal-02114308
110template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
112 T, DESCRIPTOR,
113 MOMENTA,
116>;
117
118template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
120 T, DESCRIPTOR,
121 MOMENTA,
124>;
125
126template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
128 T, DESCRIPTOR,
129 MOMENTA,
133>;
134
135template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
137 T, DESCRIPTOR,
138 MOMENTA,
142>;
143
144template <typename T, typename DESCRIPTOR>
146 T, DESCRIPTOR,
152 >,
155>;
156
157template <typename T, typename DESCRIPTOR>
159 T, DESCRIPTOR,
165 >,
169>;
170
171template <typename T, typename DESCRIPTOR>
173 T, DESCRIPTOR,
179 >,
182>;
183
184template <typename T, typename DESCRIPTOR>
186 T, DESCRIPTOR,
192 >,
196>;
197
199template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
201 T, DESCRIPTOR,
202 MOMENTA,
205>;
206
208
213template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
215 T, DESCRIPTOR,
216 MOMENTA,
220>;
221
222// BGK collision step with external force (Guo) for multiple component lattices
223template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::ExternalVelocityTuple>
225 T, DESCRIPTOR,
226 MOMENTA,
230>;
231
233
238template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
240 T, DESCRIPTOR,
241 MOMENTA,
245>;
246
248
253template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
255 T, DESCRIPTOR,
256 MOMENTA,
260>;
261
263
268template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
270 T, DESCRIPTOR,
271 MOMENTA,
275>;
276
278template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
280
282
285template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
287 T, DESCRIPTOR,
288 MOMENTA,
292>;
293
295template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
297 T, DESCRIPTOR,
298 MOMENTA,
302>;
303
305
311template <typename T, typename DESCRIPTOR>
313 T, DESCRIPTOR,
318>;
319
320template <typename T, typename DESCRIPTOR>
322 T, DESCRIPTOR,
327>;
328
329template <typename T, typename DESCRIPTOR>
331 T, DESCRIPTOR,
336>;
337
338template <typename T, typename DESCRIPTOR>
340 T, DESCRIPTOR,
345>;
346
347template <typename T, typename DESCRIPTOR>
349 T, DESCRIPTOR,
354>;
355
356template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::ExternalVelocityTuple>
358 T, DESCRIPTOR,
359 MOMENTA,
363>;
364
365template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::ExternalVelocityTuple>
367 T, DESCRIPTOR,
368 MOMENTA,
372>;
373
374
376
381template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
383 T, DESCRIPTOR,
384 MOMENTA,
387>;
388
390
394template <typename T, typename DESCRIPTOR, typename DYNAMICS, typename MOMENTA>
395class CombinedRLBdynamics final : public dynamics::CustomCollision<T,DESCRIPTOR,MOMENTA> {
396private:
397 // Use same MOMENTA in combined and nested (boundary) dynamics
398 using CORRECTED_DYNAMICS = typename DYNAMICS::template exchange_momenta<MOMENTA>;
399
400public:
401 using MomentaF = typename MOMENTA::template type<DESCRIPTOR>;
402 using EquilibriumF = typename CORRECTED_DYNAMICS::EquilibriumF;
403
404 using parameters = typename CORRECTED_DYNAMICS::parameters;
405
406 template <typename NEW_T>
408 NEW_T, DESCRIPTOR,
409 typename DYNAMICS::template exchange_value_type<NEW_T>,
410 MOMENTA
411 >;
412
413 template <typename M>
415
416 std::type_index id() override {
417 return typeid(CombinedRLBdynamics);
418 }
419
421 return block.template getData<OperatorParameters<CombinedRLBdynamics>>();
422 }
423
424 template <concepts::Cell CELL, typename PARAMETERS, typename V=typename CELL::value_t>
426 V rho, u[DESCRIPTOR::d], pi[util::TensorVal<DESCRIPTOR>::n];
427 MomentaF().computeAllMomenta(cell,rho,u,pi);
428 V fEq[DESCRIPTOR::q] { };
429 EquilibriumF().compute(cell, rho, u, fEq);
430
431 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
432 cell[iPop] = fEq[iPop] + equilibrium<DESCRIPTOR>::template fromPiToFneq<V>(iPop, pi);
433 }
434
435 return CORRECTED_DYNAMICS().collide(cell, parameters);
436 };
437
438 void computeEquilibrium(ConstCell<T,DESCRIPTOR>& cell, T rho, const T u[DESCRIPTOR::d], T fEq[DESCRIPTOR::q]) const override {
439 EquilibriumF().compute(cell, rho, u, fEq);
440 };
441
442 std::string getName() const override {
443 return "CombinedRLBdynamics<" + CORRECTED_DYNAMICS().getName() + ">";
444 };
445
446};
447
449template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
451 T, DESCRIPTOR,
452 MOMENTA,
455>;
456
458
461template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
463 T, DESCRIPTOR,
464 MOMENTA,
468>;
469
471
476template <typename T, typename DESCRIPTOR>
478 T, DESCRIPTOR,
484 >,
487>;
488
490
495template <typename T, typename DESCRIPTOR>
497 T, DESCRIPTOR,
503 >,
506>;
507
508template <typename T, typename DESCRIPTOR>
510 T, DESCRIPTOR,
516 >,
519>;
520
521template <typename T, typename DESCRIPTOR>
523 T, DESCRIPTOR,
529 >,
532>;
533
534template <typename T, typename DESCRIPTOR>
536 T, DESCRIPTOR,
542 >,
545>;
546
548
553template <typename T, typename DESCRIPTOR>
555 T, DESCRIPTOR,
559>;
560
565template <typename T, typename DESCRIPTOR>
567 T, DESCRIPTOR,
571>;
572
574template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::PoissonTuple>
576 T, DESCRIPTOR,
577 MOMENTA,
580>;
581
583template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::P1Tuple>
585 T, DESCRIPTOR,
586 MOMENTA,
589>;
590
592template <typename T, typename DESCRIPTOR>
594 T, DESCRIPTOR,
595 momenta::Tuple<
596 momenta::BulkDensity,
597 momenta::ZeroMomentum,
598 momenta::ZeroStress,
599 momenta::DefineSeparately
600 >
601> {
602 using EquilibriumF = typename equilibria::None::template type<DESCRIPTOR,momenta::Tuple<
607 >>;
608
610
611 template <typename NEW_T>
613
614 template <typename M>
616
617 std::type_index id() override {
618 return typeid(ZeroDistributionDynamics);
619 };
620
622 return block.template getData<OperatorParameters<ZeroDistributionDynamics>>();
623 }
624
625 template <typename CELL, typename PARAMETERS, typename V=typename CELL::value_t>
627 for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) {
628 cell[iPop] = -descriptors::t<T,DESCRIPTOR>(iPop);
629 }
630 return {-1, -1};
631 };
632
633 void computeEquilibrium(ConstCell<T,DESCRIPTOR>& cell, T rho, const T u[DESCRIPTOR::d], T fEq[DESCRIPTOR::q]) const override {
634 for ( int iPop = 0; iPop < DESCRIPTOR::q; iPop++ ) {
635 fEq[iPop] = 0;
636 }
637 };
638
639 std::string getName() const override {
640 return "ZeroDistributionDynamics";
641 };
642};
643
644template <typename T, typename DESCRIPTOR, typename MOMENTA>
646 T, DESCRIPTOR,
647 MOMENTA,
650>;
651
653
657template <typename T, typename DESCRIPTOR>
659 T, DESCRIPTOR,
663>;
664
666
670template <typename T, typename DESCRIPTOR>
672 T, DESCRIPTOR,
676>;
677
679
683template <typename T, typename DESCRIPTOR, typename MOMENTA=momenta::BulkTuple>
684struct ForcedVANSBGKdynamics final : public dynamics::CustomCollision<T,DESCRIPTOR,MOMENTA> {
685 using MomentaF = typename MOMENTA::template type<DESCRIPTOR>;
686
688
689 template <typename NEW_T>
691
692 std::type_index id() override {
693 return typeid(ForcedVANSBGKdynamics);
694 }
695
697 return block.template getData<OperatorParameters<ForcedVANSBGKdynamics>>();
698 }
699
700 void computeU(ConstCell<T,DESCRIPTOR>& cell, T u[DESCRIPTOR::d]) const override {
701 T rho;
702 T porosity = cell.template getField<descriptors::POROSITY>();
703 MomentaF().computeVANSRhoU(cell, rho, u, &porosity);
704 auto force = cell.template getFieldPointer<descriptors::FORCE>();
705 for (int iVel=0; iVel < DESCRIPTOR::d; ++iVel) {
706 u[iVel] += force[iVel] * T{0.5};
707 }
708 }
709
710 void computeRhoU(ConstCell<T,DESCRIPTOR>& cell, T& rho, T u[DESCRIPTOR::d]) const override {
711 T porosity = cell.template getField<descriptors::POROSITY>();
712 MomentaF().computeVANSRhoU(cell, rho, u, &porosity);
713 auto force = cell.template getFieldPointer<descriptors::FORCE>();
714 for (int iVel=0; iVel<DESCRIPTOR::d; ++iVel) {
715 u[iVel] += force[iVel] * T{0.5};
716 }
717 }
718
719 template <typename CELL, typename PARAMETERS, typename V=typename CELL::value_t>
720 CellStatistic<V> collide(CELL& cell, PARAMETERS& parameters) {
721 V rho, u[DESCRIPTOR::d];
722 V porosity = cell.template getField<descriptors::POROSITY>();
723 MomentaF().computeVANSRhoU(cell, rho, u, &porosity);
724 auto force = cell.template getFieldPointer<descriptors::FORCE>();
725 for (int iVel=0; iVel < DESCRIPTOR::d; ++iVel) {
726 u[iVel] += force[iVel] * V{0.5};
727 u[iVel] /= porosity;
728 }
729 rho *= porosity;
730 const auto omega = parameters.template get<descriptors::OMEGA>();
731 V uSqr = lbm<DESCRIPTOR>::bgkCollision(cell, rho, u, omega);
732 lbm<DESCRIPTOR>::addExternalForce(cell, u, omega, rho);
733 return {rho / porosity, uSqr * porosity * porosity};
734 }
735
736};
737
738}
739
740#endif
741
742#include "legacy/dynamics.h"
Regularized BGK collision followed by any other Dynamics.
Definition dynamics.h:395
AbstractParameters< T, DESCRIPTOR > & getParameters(BlockLattice< T, DESCRIPTOR > &block) override
Parameters access for legacy post processors.
Definition dynamics.h:420
typename CORRECTED_DYNAMICS::EquilibriumF EquilibriumF
Definition dynamics.h:402
std::string getName() const override
Return human-readable name.
Definition dynamics.h:442
typename MOMENTA::template type< DESCRIPTOR > MomentaF
Definition dynamics.h:401
CellStatistic< V > collide(CELL &cell, PARAMETERS &parameters) any_platform
Definition dynamics.h:425
void computeEquilibrium(ConstCell< T, DESCRIPTOR > &cell, T rho, const T u[DESCRIPTOR::d], T fEq[DESCRIPTOR::q]) const override
Return iPop equilibrium for given first and second momenta.
Definition dynamics.h:438
typename CORRECTED_DYNAMICS::parameters parameters
Definition dynamics.h:404
std::type_index id() override
Expose unique type-identifier for RTTI.
Definition dynamics.h:416
Highest-level interface to read-only Cell data.
Definition interface.h:36
Dynamics combination rule implementing the forcing scheme by Shan and Chen.
Definition forcing.h:462
Instantiation of Momenta tuples which define the computation and definition of momenta in a specific ...
Interface for post-processing steps – header file.
File provides a generic interface for the computation and definition of momenta (density,...
constexpr T t(unsigned iPop, tag::CUM) any_platform
Definition cum.h:108
Tuple< ZeroDensity, ZeroMomentum, ZeroStress, DefineSeparately > None
Definition aliases.h:255
Tuple< FixedDensity, FixedVelocityMomentumGeneric, ZeroStress, DefineSeparately > EquilibriumBoundaryTuple
Velocity and density are stored in external fields.
Definition aliases.h:118
Tuple< FixedDensity, FixedVelocityMomentum, BulkStress, DefineUSeparately > ExternalVelocityFixedDensityTuple
The Velocity is stored in descriptors::VELOCITY (and computed e.g.
Definition aliases.h:68
Top level namespace for all of OpenLB.
#define any_platform
Define preprocessor macros for device-side functions, constant storage.
Definition platform.h:77
Dynamic access interface for FIELD-valued parameters.
Return value of any collision.
Definition interface.h:45
VANS BGK collision step with external force.
Definition dynamics.h:684
void computeU(ConstCell< T, DESCRIPTOR > &cell, T u[DESCRIPTOR::d]) const override
Compute fluid velocity.
Definition dynamics.h:700
typename MOMENTA::template type< DESCRIPTOR > MomentaF
Definition dynamics.h:685
AbstractParameters< T, DESCRIPTOR > & getParameters(BlockLattice< T, DESCRIPTOR > &block) override
Parameters access for legacy post processors.
Definition dynamics.h:696
std::type_index id() override
Expose unique type-identifier for RTTI.
Definition dynamics.h:692
void computeRhoU(ConstCell< T, DESCRIPTOR > &cell, T &rho, T u[DESCRIPTOR::d]) const override
Compute fluid velocity and particle density.
Definition dynamics.h:710
CellStatistic< V > collide(CELL &cell, PARAMETERS &parameters)
Definition dynamics.h:720
Models a density sink by enforcing a zero distribution on the cell.
Definition dynamics.h:601
AbstractParameters< T, DESCRIPTOR > & getParameters(BlockLattice< T, DESCRIPTOR > &block) override
Parameters access for legacy post processors.
Definition dynamics.h:621
CellStatistic< V > collide(CELL &cell, PARAMETERS &parameters) any_platform
Definition dynamics.h:626
typename equilibria::None::template type< DESCRIPTOR, momenta::Tuple< momenta::BulkDensity, momenta::ZeroMomentum, momenta::ZeroStress, momenta::DefineSeparately > > EquilibriumF
Definition dynamics.h:602
std::string getName() const override
Return human-readable name.
Definition dynamics.h:639
void computeEquilibrium(ConstCell< T, DESCRIPTOR > &cell, T rho, const T u[DESCRIPTOR::d], T fEq[DESCRIPTOR::q]) const override
Return iPop equilibrium for given first and second momenta.
Definition dynamics.h:633
std::type_index id() override
Expose unique type-identifier for RTTI.
Definition dynamics.h:617
Nguyen-Ladd Velocity Correction using momenta-defined velocity.
Definition collision.h:323
Override COLLISION parameter OMEGA with inverse of cell field TAU_EFF.
Override COLLISION parameter PARAMETER with cell field PARAMETER.
Dynamics constructed as a tuple of momenta, equilibrium and collision.
Definition interface.h:308
Dynamics combination rule implementing the forcing scheme by Guo et al. for 3rd order RLB.
Dynamics combination rule implementing the forcing scheme by Guo et al.
Definition forcing.h:38
Dynamics combination rule implementing the forcing scheme by Kupershtokh et al.
Definition forcing.h:346
Dynamics combination rule implementing the forcing scheme by Guo et al. with barycentric velocity.
Definition forcing.h:121
Dynamics combination rule implementing the forcing scheme by A.J. Wagner, Phys. Rev....
Definition forcing.h:389
static void addExternalForce(CELL &cell, const RHO &rho, const U &u, const OMEGA &omega, const FORCE &force) any_platform
Add a force term after BGK collision.
Definition lbm.h:545
static V bgkCollision(CELL &cell, const RHO &rho, const VELOCITY &u, const OMEGA &omega) any_platform
BGK collision step.
Definition lbm.h:372
Plain wrapper for list of types.
Definition meta.h:276
Standard computation for density in the bulk as zeroth moment of the population.
Definition elements.h:150
Standard stress computation as second moment of the population.
Definition elements.h:1624
The momenta are defined one after the other.
When momenta are changed, the equilibrium part of the population is modified while the non-equilibriu...
The density is fixed and stored in the external field RHO.
Definition elements.h:237
Momentum is zero at solid material.
Definition elements.h:1442
The stress is always zero.
Definition elements.h:1779
Compute number of elements of a symmetric d-dimensional tensor.
Definition util.h:216
Set of functions commonly used in LB computations – header file.