OpenLB 1.7
Loading...
Searching...
No Matches
operator.h
Go to the documentation of this file.
1/* This file is part of the OpenLB library
2 *
3 * Copyright (C) 2021 Adrian Kummerlaender
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
24#ifndef SIMD_OPERATOR_H
25#define SIMD_OPERATOR_H
26
27#include "mask.h"
28
29#include "core/meta.h"
30#include "core/operator.h"
31
33
34#include "dynamics/dynamics.h"
35
37
38namespace olb {
39
40template <typename T>
41struct CellStatistic<cpu::simd::Pack<T>> {
44
45 operator bool() const {
46 return true;
47 }
48};
49
50namespace cpu {
51
53namespace simd {
54
56template <typename T, unsigned D>
57class FieldPtr : public ScalarVector<Pack<T>,D,FieldPtr<T,D>> {
58private:
59 ColumnVector<Column<T>,D>& _data;
60 std::size_t _index;
61
62 friend typename ScalarVector<Pack<T>,D,FieldPtr>::type;
63
64 std::array<Pack<T>,D> _packs;
65
66protected:
67 const Pack<T>* getComponentPointer(unsigned iDim) const
68 {
69 return &_packs[iDim];
70 }
72 {
73 return &_packs[iDim];
74 }
75
76public:
77 FieldPtr(ColumnVector<Column<T>,D>& columns, std::size_t index):
78 _data(columns),
79 _index(index) {
80 meta::call_n_times<D>([&](unsigned iD) {
81 _packs[iD] = &_data[iD][_index];
82 });
83 }
84
86 _data(rhs._data),
87 _index(rhs._index),
88 _packs(rhs._packs) { }
89
90 template <typename U, typename IMPL>
92 {
93 for (unsigned iDim=0; iDim < D; ++iDim) {
94 this->operator[](iDim) = rhs[iDim];
95 }
96 return *this;
97 }
98
99};
100
102
108template <typename T, typename DESCRIPTOR, typename V, typename... RW_FIELDS>
109class Cell {
110private:
111 using rw_fields = meta::list<RW_FIELDS...>;
112
114 std::size_t _iCell;
115 Mask<T>& _mask;
116
118
121 std::tuple<FieldD<V,DESCRIPTOR,RW_FIELDS>...> _fields;
122
123public:
124 using value_t = V;
125 using descriptor_t = DESCRIPTOR;
126
129 _lattice(lattice),
130 _iCell(iCell),
131 _mask(mask)
132 {
133 rw_fields::for_each([&](auto field) {
134 using FIELD = typename decltype(field)::type;
135 auto& array = _lattice.template getField<FIELD>();
136 auto& pack = std::get<(rw_fields::template index<FIELD>())>(_fields);
137 //meta::call_n_times<(DESCRIPTOR::template size<FIELD>())>([&](unsigned iD) {
138 for (unsigned iD=0; iD < DESCRIPTOR::template size<FIELD>(); ++iD) {
139 pack[iD] = cpu::simd::Pack<T>(&array[iD][_iCell]);
140 }
141 });
142 }
143
146 {
147 rw_fields::for_each([&](auto field) {
148 using FIELD = typename decltype(field)::type;
149 auto& array = _lattice.template getField<FIELD>();
150 auto& pack = std::get<(rw_fields::template index<FIELD>())>(_fields);
151 //meta::call_n_times<(DESCRIPTOR::template size<FIELD>())>([&](unsigned iD) {
152 for (unsigned iD=0; iD < DESCRIPTOR::template size<FIELD>(); ++iD) {
153 cpu::simd::maskstore(&array[iD][_iCell], _mask, pack[iD]);
154 }
155 });
156 }
157
159 V& operator[](unsigned iPop) {
160 return std::get<(rw_fields::template index<descriptors::POPULATION>())>(_fields)[iPop];
161 }
162
164 template <typename FIELD>
165 auto getField() const {
166 if constexpr (rw_fields::template contains<FIELD>()) {
167 if constexpr (DESCRIPTOR::template size<FIELD>() == 1) {
168 return std::get<(rw_fields::template index<FIELD>())>(_fields)[0];
169 } else {
170 return std::get<(rw_fields::template index<FIELD>())>(_fields);
171 }
172 } else {
173 auto& fieldArray = _lattice.template getField<FIELD>();
174 if constexpr (DESCRIPTOR::template size<FIELD>() == 1) {
175 return cpu::simd::Pack<T>(&fieldArray[0][_iCell]);
176 } else {
177 return FieldD<V,DESCRIPTOR,FIELD>([&](unsigned iD) {
178 return &fieldArray[iD][_iCell];
179 });
180 }
181 }
182 __builtin_unreachable();
183 }
184
186 template <typename FIELD>
188 if constexpr (rw_fields::template contains<FIELD>()) {
189 std::get<(rw_fields::template index<FIELD>())>(_fields) = value;
190 } else {
191 auto& array = _lattice.template getField<FIELD>();
192 for (unsigned iD=0; iD < DESCRIPTOR::template size<FIELD>(); ++iD) {
193 cpu::simd::maskstore(&array[iD][_iCell], _mask, value[iD]);
194 }
195 }
196 }
197
199 template <typename FIELD>
200 std::enable_if_t<rw_fields::template contains<FIELD>(), FieldD<V,DESCRIPTOR,FIELD>&>
202 return std::get<(rw_fields::template index<FIELD>())>(_fields);
203 }
204
206 template <typename FIELD>
207 std::enable_if_t<!rw_fields::template contains<FIELD>(), FieldD<V,DESCRIPTOR,FIELD>>
209 return getField<FIELD>();
210 }
211
213 template <typename FIELD>
214 std::enable_if_t<rw_fields::template contains<FIELD>(),V&>
215 getFieldComponent(unsigned iD) {
216 return std::get<(rw_fields::template index<FIELD>())>(_fields)[iD];
217 }
218
220 template <typename FIELD>
221 std::enable_if_t<!rw_fields::template contains<FIELD>(),V>
222 getFieldComponent(unsigned iD) {
223 return &_lattice.template getField<FIELD>()[iD][_iCell];
224 }
225
226};
227
228
230template <typename T, typename DESCRIPTOR, typename DYNAMICS>
231class ConcreteDynamics final : public cpu::Dynamics<T,DESCRIPTOR,Platform::CPU_SIMD> {
232private:
234
235public:
237 _parameters{parameters} {
238 }
239
241 return DYNAMICS().apply(cell, *_parameters);
242 }
243
245 return typename DYNAMICS::MomentaF().computeRho(cell);
246 }
248 typename DYNAMICS::MomentaF().computeU(cell, u);
249 }
251 typename DYNAMICS::MomentaF().computeJ(cell, j);
252 }
253 void computeRhoU(cpu::Cell<T,DESCRIPTOR,Platform::CPU_SIMD>& cell, T& rho, T* u) override {
254 typename DYNAMICS::MomentaF().computeRhoU(cell, rho, u);
255 }
256 void computeStress(cpu::Cell<T,DESCRIPTOR,Platform::CPU_SIMD>& cell, T& rho, T* u, T* pi) override {
257 typename DYNAMICS::MomentaF().computeStress(cell, rho, u, pi);
258 }
259 void computeAllMomenta(cpu::Cell<T,DESCRIPTOR,Platform::CPU_SIMD>& cell, T& rho, T* u, T* pi) override {
260 typename DYNAMICS::MomentaF().computeAllMomenta(cell, rho, u, pi);
261 }
262
263 T getOmegaOrFallback(T fallback) override {
264 if constexpr (DYNAMICS::parameters::template contains<descriptors::OMEGA>()) {
265 return _parameters->template get<descriptors::OMEGA>();
266 } else {
267 return fallback;
268 }
269 __builtin_unreachable();
270 }
271
272 T computeEquilibrium(int iPop, T rho, T* u) override {
273 return DYNAMICS().computeEquilibrium(iPop, rho, u);
274 }
275
277 typename DYNAMICS::MomentaF().defineRho(cell, rho);
278 }
279
281 typename DYNAMICS::MomentaF().defineU(cell, u);
282 }
283
284 void defineRhoU(cpu::Cell<T,DESCRIPTOR,Platform::CPU_SIMD>& cell, T& rho, T* u) override {
285 typename DYNAMICS::MomentaF().defineRhoU(cell, rho, u);
286 }
287
288 void defineAllMomenta(cpu::Cell<T,DESCRIPTOR,Platform::CPU_SIMD>& cell, T& rho, T* u, T* pi) override {
289 typename DYNAMICS::MomentaF().defineAllMomenta(cell, rho, u, pi);
290 }
291
293 typename DYNAMICS::MomentaF().inverseShiftRhoU(cell, rho, u);
294 }
295};
296
297}
298
299}
300
301
303
309template <typename T, typename DESCRIPTOR, typename DYNAMICS>
310class ConcreteBlockCollisionO<T,DESCRIPTOR,Platform::CPU_SIMD,DYNAMICS> final
311 : public BlockCollisionO<T,DESCRIPTOR,Platform::CPU_SIMD> {
312private:
313 std::unique_ptr<DYNAMICS> _dynamics;
314 std::unique_ptr<cpu::Dynamics<T,DESCRIPTOR,Platform::CPU_SIMD>> _concreteDynamics;
315
318
320
323 typename LatticeStatistics<T>::Aggregatable& statistics,
324 std::size_t iCell)
325 {
327 if (auto cellStatistic = _dynamicsOfCells[iCell]->collide(cell)) {
328 statistics.increment(cellStatistic.rho, cellStatistic.uSqr);
329 }
330 }
331
337 typename LatticeStatistics<T>::Aggregatable& statistics,
338 std::size_t iCell)
339 {
340 if constexpr (dynamics::is_vectorizable_v<DYNAMICS>) {
341 if (cpu::simd::Mask<T> m = {mask.raw(), iCell}) {
343 auto cellStatistic = DYNAMICS().apply(cell, parameters);
344 for (unsigned i=0; i < cpu::simd::Pack<T>::size; ++i) {
345 if (mask[iCell+i]) {
346 if (cellStatistic.rho[i] != T{-1}) {
347 statistics.increment(cellStatistic.rho[i], cellStatistic.uSqr[i]);
348 }
349 } else if (subdomain[iCell+i]) {
350 applyOther(block, statistics, iCell+i);
351 }
352 }
353 } else {
354 for (std::size_t i=iCell; i < iCell+cpu::simd::Pack<T>::size; ++i) {
355 if (subdomain[i]) {
356 applyOther(block, statistics, i);
357 }
358 }
359 }
360 }
361 }
362
363public:
365 _dynamics(new DYNAMICS()),
366 _parameters(nullptr),
367 _mask(nullptr)
368 { }
369
370 std::type_index id() const override
371 {
372 return typeid(DYNAMICS);
373 }
374
375 std::size_t weight() const override
376 {
377 return _mask->weight();
378 }
379
380 void set(CellID iCell, bool state, bool overlap) override
381 {
383 if constexpr (!std::is_same_v<DYNAMICS,NoDynamics<T,DESCRIPTOR>>) {
384 if (!overlap) {
385 _mask->set(iCell, state);
386 }
387 }
388 if (state) {
389 _dynamicsOfCells[iCell] = _concreteDynamics.get();
390 }
391 }
392
394 {
395 return _dynamics.get();
396 }
397
399 {
400 _parameters = &block.template getData<OperatorParameters<DYNAMICS>>().parameters;
401 _mask = &block.template getData<DynamicsMask<DYNAMICS>>();
402 if constexpr (dynamics::has_parametrized_momenta_v<DYNAMICS>) {
403 _dynamics->setMomentaParameters(_parameters);
404 }
405
406 _concreteDynamics.reset(new cpu::simd::ConcreteDynamics<T,DESCRIPTOR,DYNAMICS>(_parameters));
407 // Fetch pointer to concretized dynamic-dispatch field
408 _dynamicsOfCells = block.template getField<cpu::DYNAMICS<T,DESCRIPTOR,Platform::CPU_SIMD>>()[0].data();
409 }
410
414 CollisionDispatchStrategy strategy) override
415 {
416 if (strategy != CollisionDispatchStrategy::Dominant) {
417 throw std::runtime_error("Platform::CPU_SIMD currently only support CollisionDispatchStrategy::Dominant");
418 }
419
420 auto& mask = *_mask;
421 typename LatticeStatistics<T>::Aggregatable statistics{};
422 #ifdef PARALLEL_MODE_OMP
423 #pragma omp declare reduction(+ : typename LatticeStatistics<T>::Aggregatable : omp_out += omp_in) initializer (omp_priv={})
424 #endif
425
426 if constexpr (dynamics::is_vectorizable_v<DYNAMICS>) {
427 // Ensure that serialized mask storage is up-to-date
428 mask.setProcessingContext(ProcessingContext::Simulation);
429 // Apply collision to cells
430 #ifdef PARALLEL_MODE_OMP
431 #pragma omp parallel for schedule(static) reduction(+ : statistics)
432 #endif
433 for (CellID iCell=0; iCell < block.getNcells(); iCell += cpu::simd::Pack<T>::size) {
434 apply(block, subdomain, mask, *_parameters, statistics, iCell);
435 }
436 } else { // Fallback for non-vectorizable collision operators
437 #ifdef PARALLEL_MODE_OMP
438 #pragma omp parallel for schedule(static) reduction(+ : statistics)
439 #endif
440 for (std::size_t iCell=0; iCell < block.getNcells(); ++iCell) {
441 if (mask[iCell]) {
443 if (auto cellStatistic = DYNAMICS().apply(cell, *_parameters)) {
444 statistics.increment(cellStatistic.rho, cellStatistic.uSqr);
445 }
446 } else if (subdomain[iCell]) {
447 applyOther(block, statistics, iCell);
448 }
449 }
450 }
451
452 block.getStatistics().incrementStats(statistics);
453 }
454
455};
456
457
459template <typename T, typename DESCRIPTOR, CONCEPT(CellOperator) OPERATOR>
460class ConcreteBlockO<T,DESCRIPTOR,Platform::CPU_SIMD,OPERATOR,OperatorScope::PerCell> final
461 : public BlockO<T,DESCRIPTOR,Platform::CPU_SIMD> {
462private:
463 std::vector<CellID> _cells;
464 bool _modified;
465
466public:
467 ConcreteBlockO() = default;
468
469 std::type_index id() const override
470 {
471 return typeid(OPERATOR);
472 }
473
474 void set(CellID iCell, bool state) override
475 {
476 if (state) {
477 _cells.emplace_back(iCell);
478 _modified = true;
479 }
480 }
481
484
486 {
487 if (_modified) {
488 std::sort(_cells.begin(), _cells.end());
489 _cells.erase(std::unique(_cells.begin(), _cells.end()), _cells.end());
490 _modified = false;
491 }
492 if (_cells.size() > 0) {
494 #ifdef PARALLEL_MODE_OMP
495 #pragma omp parallel for schedule(static) firstprivate(cell)
496 #endif
497 for (CellID iCell : _cells) {
498 cell.setCellId(iCell);
499 OPERATOR().apply(cell);
500 }
501 }
502 }
503
504};
505
506
507template <typename T, typename DESCRIPTOR, CONCEPT(CellOperator) OPERATOR>
509 : public BlockO<T,DESCRIPTOR,Platform::CPU_SIMD> {
510private:
511 std::vector<CellID> _cells;
512 bool _modified;
513
515
516public:
517 ConcreteBlockO() = default;
518
519 std::type_index id() const override
520 {
521 return typeid(OPERATOR);
522 }
523
524 void set(CellID iCell, bool state) override
525 {
526 if (state) {
527 _cells.emplace_back(iCell);
528 _modified = true;
529 }
530 }
531
533 {
534 _parameters = &block.template getData<OperatorParameters<OPERATOR>>().parameters;
535 }
536
538 {
539 if (_modified) {
540 std::sort(_cells.begin(), _cells.end());
541 _cells.erase(std::unique(_cells.begin(), _cells.end()), _cells.end());
542 _modified = false;
543 }
544 if (_cells.size() > 0) {
546 #ifdef PARALLEL_MODE_OMP
547 #pragma omp parallel for schedule(static) firstprivate(cell)
548 #endif
549 for (CellID iCell : _cells) {
550 cell.setCellId(iCell);
551 OPERATOR().apply(cell, *_parameters);
552 }
553 }
554 }
555
556};
557
558
560
563template <typename T, typename DESCRIPTOR, CONCEPT(BlockOperator) OPERATOR>
564class ConcreteBlockO<T,DESCRIPTOR,Platform::CPU_SIMD,OPERATOR,OperatorScope::PerBlock> final
565 : public BlockO<T,DESCRIPTOR,Platform::CPU_SIMD> {
566public:
567 std::type_index id() const override
568 {
569 return typeid(OPERATOR);
570 }
571
572 void set(CellID iCell, bool state) override
573 {
574 throw std::logic_error("BlockO::set not supported for OperatorScope::PerBlock");
575 }
576
578 {
579 OPERATOR().setup(block);
580 }
581
583 {
584 OPERATOR().apply(block);
585 }
586
587};
588
589}
590
591#endif
LatticeStatistics< T > & getStatistics()
Return a handle to the LatticeStatistics object.
std::size_t getNcells() const
Get number of cells.
Vector of columns.
std::size_t weight() const override
Returns number of assigned cells.
Definition operator.h:375
void setup(ConcreteBlockLattice< T, DESCRIPTOR, Platform::CPU_SIMD > &block) override
Definition operator.h:398
void apply(ConcreteBlockLattice< T, DESCRIPTOR, Platform::CPU_SIMD > &block, ConcreteBlockMask< T, Platform::CPU_SIMD > &subdomain, CollisionDispatchStrategy strategy) override
Apply collision on entire block.
Definition operator.h:412
void set(CellID iCell, bool state, bool overlap) override
Set whether iCell is covered by the present collision step.
Definition operator.h:380
Collision operation of concrete DYNAMICS on concrete block lattices of PLATFORM.
Implementation of BlockLattice on a concrete PLATFORM.
void apply(ConcreteBlockLattice< T, DESCRIPTOR, Platform::CPU_SIMD > &block) override
Definition operator.h:537
void setup(ConcreteBlockLattice< T, DESCRIPTOR, Platform::CPU_SIMD > &block) override
Definition operator.h:532
void set(CellID iCell, bool state) override
Set whether iCell is covered by the operator (optional)
Definition operator.h:524
void setup(ConcreteBlockLattice< T, DESCRIPTOR, Platform::CPU_SIMD > &block) override
Definition operator.h:577
void set(CellID iCell, bool state) override
Set whether iCell is covered by the operator (optional)
Definition operator.h:572
void apply(ConcreteBlockLattice< T, DESCRIPTOR, Platform::CPU_SIMD > &block) override
Definition operator.h:582
void setup(ConcreteBlockLattice< T, DESCRIPTOR, Platform::CPU_SIMD > &block) override
Definition operator.h:482
void set(CellID iCell, bool state) override
Set whether iCell is covered by the operator (optional)
Definition operator.h:474
void apply(ConcreteBlockLattice< T, DESCRIPTOR, Platform::CPU_SIMD > &block) override
Definition operator.h:485
Block application of concrete OPERATOR called using SCOPE on PLATFORM.
Definition operator.h:65
Plain old scalar vector.
Definition vector.h:47
Cell concept for concrete block lattices on CPU platforms.
Definition cell.h:87
void setCellId(std::size_t iCell)
Definition cell.h:108
Implementation of the Cell concept for vectorized collision operators.
Definition operator.h:109
auto getField() const
Return pack-valued copy of FIELD.
Definition operator.h:165
V & operator[](unsigned iPop)
Return reference to iPop population pack.
Definition operator.h:159
std::enable_if_t< rw_fields::template contains< FIELD >(), FieldD< V, DESCRIPTOR, FIELD > & > getFieldPointer()
Return reference to pack-valued interim storage vector of r/w field.
Definition operator.h:201
std::enable_if_t<!rw_fields::template contains< FIELD >(), FieldD< V, DESCRIPTOR, FIELD > > getFieldPointer()
Return pack-valued copy of non r/w field.
Definition operator.h:208
DESCRIPTOR descriptor_t
Definition operator.h:125
void setField(FieldD< V, DESCRIPTOR, FIELD > &&value)
Set compoents of FIELD from pack-valued vector.
Definition operator.h:187
std::enable_if_t<!rw_fields::template contains< FIELD >(), V > getFieldComponent(unsigned iD)
Return pack-valued copy of non r/w field component.
Definition operator.h:222
std::enable_if_t< rw_fields::template contains< FIELD >(), V & > getFieldComponent(unsigned iD)
Return reference to pack-valued interim storage component of r/w field.
Definition operator.h:215
Cell(ConcreteBlockLattice< T, DESCRIPTOR, Platform::CPU_SIMD > &lattice, std::size_t iCell, Mask< T > &mask)
Load r/w fields into SIMD packs.
Definition operator.h:128
~Cell()
Store modified r/w fields back into lattice taking into account the mask.
Definition operator.h:145
Plain column for SIMD CPU targets.
Definition column.h:53
Implementation of cpu::Dynamics for concrete DYNAMICS on SIMD blocks.
Definition operator.h:231
void computeJ(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell, T *j) override
Definition operator.h:250
CellStatistic< T > collide(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell) override
Definition operator.h:240
void inverseShiftRhoU(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell, T &rho, T *u) override
Definition operator.h:292
T computeRho(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell) override
Definition operator.h:244
void computeRhoU(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell, T &rho, T *u) override
Definition operator.h:253
void computeStress(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell, T &rho, T *u, T *pi) override
Definition operator.h:256
T computeEquilibrium(int iPop, T rho, T *u) override
Definition operator.h:272
void defineAllMomenta(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell, T &rho, T *u, T *pi) override
Definition operator.h:288
void defineU(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell, T *u) override
Definition operator.h:280
void computeAllMomenta(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell, T &rho, T *u, T *pi) override
Definition operator.h:259
ConcreteDynamics(ParametersOfOperatorD< T, DESCRIPTOR, DYNAMICS > *parameters)
Definition operator.h:236
T getOmegaOrFallback(T fallback) override
Definition operator.h:263
void defineRhoU(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell, T &rho, T *u) override
Definition operator.h:284
void computeU(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell, T *u) override
Definition operator.h:247
void defineRho(cpu::Cell< T, DESCRIPTOR, Platform::CPU_SIMD > &cell, T &rho) override
Definition operator.h:276
SIMD-specific pointer to a pack of rows of a D-dimensional field.
Definition operator.h:57
Pack< T > * getComponentPointer(unsigned iDim)
Definition operator.h:71
FieldPtr< T, D > & operator=(const GenericVector< U, D, IMPL > &rhs)
Definition operator.h:91
FieldPtr(ColumnVector< Column< T >, D > &columns, std::size_t index)
Definition operator.h:77
FieldPtr(FieldPtr< T, D > &&rhs)
Definition operator.h:85
const Pack< T > * getComponentPointer(unsigned iDim) const
Definition operator.h:67
Interface for post-processing steps – header file.
void maskstore(T *target, Mask< T > mask, Pack< T > value)
Top level namespace for all of OpenLB.
std::uint32_t CellID
Type for sequential block-local cell indices.
@ Simulation
Data available on host for e.g. functor evaluation.
Platform
OpenLB execution targets.
Definition platform.h:36
@ CPU_SIMD
Basic scalar CPU.
CollisionDispatchStrategy
Collision dispatch strategy.
Definition operator.h:88
@ Dominant
Apply dominant dynamics using mask and fallback to virtual dispatch for others.
typename ParametersD< T, DESCRIPTOR >::template include< typename OPERATOR::parameters > ParametersOfOperatorD
Deduce ParametersD of OPERATOR w.r.t. T and DESCRIPTOR.
@ PerBlock
Per-block application, i.e. OPERATOR::apply is passed a ConcreteBlockLattice.
@ PerCell
Per-cell application, i.e. OPERATOR::apply is passed a CELL concept implementation.
@ PerCellWithParameters
Per-cell application with parameters, i.e. OPERATOR::apply is passed a CELL concept implementation an...
Collision operation on concrete blocks of PLATFORM.
Definition operator.h:97
Base of block-wide operators such as post processors.
Definition operator.h:41
Return value of any collision.
Definition interface.h:43
Interface for per-cell dynamics.
Definition interface.h:54
Generic vector of values supporting basic arithmetic.
constexpr const T & operator[](unsigned iDim) const any_platform
Vector of scalars.
GenericVector< T, D, FieldPtr< T, D > > type
CPU specific field mirroring BlockDynamicsMap.
Definition cell.h:80
Virtual interface for dynamically-dispatched dynamics access on CPU targets.
Definition cell.h:39
Plain wrapper for list of types.
Definition meta.h:276
static constexpr void for_each(F f)
Calls f for each type of TYPES by-value (in reversed order!)
Definition meta.h:331