OpenLB 1.7
Loading...
Searching...
No Matches
operator.hh
Go to the documentation of this file.
1/* This file is part of the OpenLB library
2 *
3 * Copyright (C) 2022 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 GPU_CUDA_OPERATOR_HH
25#define GPU_CUDA_OPERATOR_HH
26
27#include "operator.h"
28
29#include "context.hh"
30#include "dynamics.hh"
31
32namespace olb {
33
35struct CollisionSubdomainMask;
36
38namespace gpu {
39
41namespace cuda {
42
44template <typename T, typename DESCRIPTOR, typename DYNAMICS>
46private:
50 bool* _mask;
51
52 CellStatistic<T> apply(DeviceContext<T,DESCRIPTOR>& lattice, CellID iCell) __device__ {
53 DataOnlyCell<T,DESCRIPTOR> cell(lattice, iCell);
54 return DYNAMICS().apply(cell, _parameters);
55 }
56
57public:
59
63 _parameters{parameters},
64 _mask{mask}
65 { }
66
68
72 bool operator()(DeviceContext<T,DESCRIPTOR>& lattice, CellID iCell) __device__ {
73 if (_mask[iCell]) {
74 apply(lattice, iCell);
75 return true;
76 }
77 return false;
78 }
79
81 bool operator()(DeviceContext<T,DESCRIPTOR>& lattice, CellID iCell, CellStatistic<T>& statistic) __device__ {
82 if (_mask[iCell]) {
83 statistic = apply(lattice, iCell);
84 return true;
85 }
86 return false;
87 }
88
89};
90
92template <typename T, typename DESCRIPTOR, typename DYNAMICS>
94private:
96
97public:
99 _parameters{parameters}
100 { }
101
102 bool operator()(DeviceContext<T,DESCRIPTOR>& lattice, CellID iCell) __device__ {
103 DataOnlyCell<T,DESCRIPTOR> cell(lattice, iCell);
104 DYNAMICS().apply(cell, _parameters);
105 return true;
106 }
107
108 bool operator()(DeviceContext<T,DESCRIPTOR>& lattice, CellID iCell, CellStatistic<T>& statistic) __device__ {
109 DataOnlyCell<T,DESCRIPTOR> cell(lattice, iCell);
110 statistic = DYNAMICS().apply(cell, _parameters);
111 return true;
112 }
113
114};
115
116
118template <typename OPERATOR>
120private:
122 bool* _mask;
123
124public:
126 _mask{mask}
127 { }
128
129 template <typename T, typename DESCRIPTOR>
130 bool operator()(DeviceBlockLattice<T,DESCRIPTOR>& lattice, CellID iCell) __device__ {
131 if (_mask[iCell]) {
132 Cell<T,DESCRIPTOR> cell(lattice, iCell);
133 OPERATOR().apply(cell);
134 return true;
135 }
136 return false;
137 }
138
139};
140
142
145template <typename OPERATOR>
147 template <typename T, typename DESCRIPTOR>
148 bool operator()(DeviceBlockLattice<T,DESCRIPTOR>& lattice, CellID iCell) __device__ {
149 Cell<T,DESCRIPTOR> cell(lattice, iCell);
150 OPERATOR().apply(cell);
151 return true;
152 }
153
154};
155
157template <typename T, typename DESCRIPTOR, typename OPERATOR>
159private:
161
162public:
164 _parameters{parameters}
165 { }
166
167 bool operator()(DeviceBlockLattice<T,DESCRIPTOR>& lattice, CellID iCell) __device__ {
168 Cell<T,DESCRIPTOR> cell(lattice, iCell);
169 OPERATOR().apply(cell, _parameters);
170 return true;
171 }
172
173};
174
176template <typename COUPLER>
178 template <typename CONTEXT>
179 bool operator()(CONTEXT& lattices,
180 CellID iCell) __device__ {
181 auto cells = lattices.exchange_values([&](auto name) -> auto {
182 return Cell{lattices.get(name), iCell};
183 });
184 COUPLER().apply(cells);
185 return true;
186 }
187};
188
190template <typename COUPLER, typename COUPLEES>
192private:
193 typename COUPLER::parameters::template decompose_into<
195 > _parameters;
196
197public:
198 template <typename PARAMETERS>
200 _parameters{parameters}
201 { }
202
203 template <typename CONTEXT>
204 bool operator()(CONTEXT& lattices,
205 CellID iCell) __device__ {
206 auto cells = lattices.exchange_values([&](auto name) -> auto {
207 return Cell{lattices.get(name), iCell};
208 });
209 COUPLER().apply(cells, _parameters);
210 return true;
211 }
212
213};
214
216
231template <typename T, typename DESCRIPTOR, typename... DYNAMICS>
235 bool* subdomain = block.template getData<CollisionSubdomainMask>().deviceData();
236 DeviceContext<T,DESCRIPTOR> lattice(block);
237 if (block.statisticsEnabled()) {
239 lattice,
240 subdomain,
242 block.template getData<OperatorParameters<DYNAMICS>>().parameters,
243 block.template getData<DynamicsMask<DYNAMICS>>().deviceData()
244 }...,
246 } else {
248 lattice,
249 subdomain,
251 block.template getData<OperatorParameters<DYNAMICS>>().parameters,
252 block.template getData<DynamicsMask<DYNAMICS>>().deviceData()
253 }...,
255 }
256 };
257}
258
259
261namespace kernel {
262
264template <typename CONTEXT, typename... OPERATORS>
265void call_operators(CONTEXT lattice, bool* subdomain, OPERATORS... ops) __global__ {
266 const CellID iCell = blockIdx.x * blockDim.x + threadIdx.x;
267 if (!(iCell < lattice.getNcells()) || !subdomain[iCell]) {
268 return;
269 }
270 (ops(lattice, iCell) || ... );
271}
272
274
277template <typename CONTEXT, typename... OPERATORS>
278void call_operators_with_statistics(CONTEXT lattice, bool* subdomain, OPERATORS... ops) __global__ {
279 const CellID iCell = blockIdx.x * blockDim.x + threadIdx.x;
280 if (!(iCell < lattice.getNcells()) || !subdomain[iCell]) {
281 return;
282 }
283 typename CONTEXT::value_t** statistic = lattice.template getField<descriptors::STATISTIC>();
284 int* statisticGenerated = lattice.template getField<descriptors::STATISTIC_GENERATED>()[0];
285 CellStatistic<typename CONTEXT::value_t> cellStatistic{-1, -1};
286 if ((ops(lattice, iCell, cellStatistic) || ... )) {
287 if (cellStatistic) {
288 statisticGenerated[iCell] = 1;
289 statistic[0][iCell] = cellStatistic.rho;
290 statistic[1][iCell] = cellStatistic.uSqr;
291 } else {
292 statisticGenerated[iCell] = 0;
293 statistic[0][iCell] = 0;
294 statistic[1][iCell] = 0;
295 }
296 }
297}
298
300template <typename CONTEXT, typename... OPERATORS>
301void call_list_operators(CONTEXT lattice,
302 const CellID* indices, std::size_t nIndices,
303 OPERATORS... ops) __global__ {
304 const std::size_t iIndex = blockIdx.x * blockDim.x + threadIdx.x;
305 if (!(iIndex < nIndices)) {
306 return;
307 }
308 (ops(lattice, indices[iIndex]) || ... );
309}
310
312
315template <typename CONTEXT, typename... OPERATORS>
317 const CellID* indices, std::size_t nIndices,
318 OPERATORS... ops) __global__ {
319 const std::size_t iIndex = blockIdx.x * blockDim.x + threadIdx.x;
320 if (!(iIndex < nIndices)) {
321 return;
322 }
323 typename CONTEXT::value_t** statistic = lattice.template getField<descriptors::STATISTIC>();
324 int* statisticGenerated = lattice.template getField<descriptors::STATISTIC_GENERATED>()[0];
325 CellStatistic<typename CONTEXT::value_t> cellStatistic{-1, -1};
326 if ((ops(lattice, indices[iIndex], cellStatistic) || ... )) {
327 if (cellStatistic) {
328 statisticGenerated[indices[iIndex]] = 1;
329 statistic[0][indices[iIndex]] = cellStatistic.rho;
330 statistic[1][indices[iIndex]] = cellStatistic.uSqr;
331 } else {
332 statisticGenerated[indices[iIndex]] = 0;
333 statistic[0][indices[iIndex]] = 0;
334 statistic[1][indices[iIndex]] = 0;
335 }
336 }
337}
338
340template <typename CONTEXTS, typename... OPERATORS>
341void call_coupling_operators(CONTEXTS lattices, bool* subdomain, OPERATORS... ops) __global__ {
342 const CellID iCell = blockIdx.x * blockDim.x + threadIdx.x;
343 const auto nCells = lattices.template get<0>().getNcells();
344 if (!(iCell < nCells) || !subdomain[iCell]) {
345 return;
346 }
347 (ops(lattices, iCell) || ... );
348}
349
351template <typename T, typename DESCRIPTOR, typename DYNAMICS, typename PARAMETERS=typename DYNAMICS::ParametersD>
352void construct_dynamics(void* target, PARAMETERS* parameters) __global__ {
353 new (target) ConcreteDynamics<T,DESCRIPTOR,DYNAMICS>(parameters);
354}
355
356}
357
359
362template <typename CONTEXT, typename... ARGS>
363void call_operators(CONTEXT& lattice, bool* subdomain, ARGS&&... args) {
364 const auto block_size = 32;
365 const auto block_count = (lattice.getNcells() + block_size - 1) / block_size;
366 kernel::call_operators<CONTEXT,ARGS...><<<block_count,block_size>>>(
367 lattice, subdomain, std::forward<decltype(args)>(args)...);
369}
370
372template <typename CONTEXT, typename... ARGS>
373void async_call_operators(cudaStream_t stream, CONTEXT& lattice, bool* subdomain, ARGS&&... args) {
374 const auto block_size = 32;
375 const auto block_count = (lattice.getNcells() + block_size - 1) / block_size;
376 kernel::call_operators<CONTEXT,ARGS...><<<block_count,block_size,0,stream>>>(
377 lattice, subdomain, std::forward<decltype(args)>(args)...);
379}
380
382
385template <typename CONTEXT, typename... ARGS>
386void call_operators_with_statistics(CONTEXT& lattice, bool* subdomain, ARGS&&... args) {
387 const auto block_size = 32;
388 const auto block_count = (lattice.getNcells() + block_size - 1) / block_size;
389 kernel::call_operators_with_statistics<CONTEXT,ARGS...><<<block_count,block_size>>>(
390 lattice, subdomain, std::forward<decltype(args)>(args)...);
392}
393
395template <typename CONTEXT, typename... ARGS>
396void async_call_operators_with_statistics(cudaStream_t stream, CONTEXT& lattice, bool* subdomain, ARGS&&... args) {
397 const auto block_size = 32;
398 const auto block_count = (lattice.getNcells() + block_size - 1) / block_size;
399 kernel::call_operators_with_statistics<CONTEXT,ARGS...><<<block_count,block_size,0,stream>>>(
400 lattice, subdomain, std::forward<decltype(args)>(args)...);
402}
403
405
408template <typename CONTEXT, typename... ARGS>
409void call_list_operators(CONTEXT& lattice,
410 const gpu::cuda::Column<CellID>& cells,
411 ARGS&&... args) {
412 const auto block_size = 32;
413 const auto block_count = (cells.size() + block_size - 1) / block_size;
414 kernel::call_list_operators<CONTEXT,ARGS...><<<block_count, block_size>>>(
415 lattice,
416 cells.deviceData(), cells.size(),
417 std::forward<decltype(args)>(args)...);
419}
420
422template <typename CONTEXT, typename... ARGS>
423void async_call_list_operators(cudaStream_t stream,
424 CONTEXT& lattice,
425 const gpu::cuda::Column<CellID>& cells,
426 ARGS&&... args) {
427 const auto block_size = 32;
428 const auto block_count = (cells.size() + block_size - 1) / block_size;
429 kernel::call_list_operators<<<block_count,block_size,0,stream>>>(
430 lattice,
431 cells.deviceData(), cells.size(),
432 std::forward<decltype(args)>(args)...);
434}
435
437template <typename CONTEXT, typename... ARGS>
439 CONTEXT& lattice,
440 const gpu::cuda::Column<CellID>& cells,
441 ARGS&&... args) {
442 const auto block_size = 32;
443 const auto block_count = (cells.size() + block_size - 1) / block_size;
444 kernel::call_list_operators_with_statistics<<<block_count,block_size,0,stream>>>(
445 lattice,
446 cells.deviceData(), cells.size(),
447 std::forward<decltype(args)>(args)...);
449}
450
452template <typename CONTEXT, typename... ARGS>
453void call_coupling_operators(CONTEXT& lattices, bool* subdomain, ARGS&&... args) {
454 const auto nCells = lattices.template get<0>().getNcells();
455 const auto block_size = 32;
456 const auto block_count = (nCells + block_size - 1) / block_size;
457 kernel::call_coupling_operators<CONTEXT,ARGS...><<<block_count,block_size>>>(
458 lattices, subdomain, std::forward<decltype(args)>(args)...);
460}
461
462}
463
464}
465
466template <typename T, typename DESCRIPTOR, typename DYNAMICS>
468 _dynamics(new DYNAMICS()),
469 _parameters(nullptr),
470 _mask(nullptr),
471 _cells(0),
472 _modified(true),
473 _stream(cudaStreamDefault)
474{ }
475
476template <typename T, typename DESCRIPTOR, typename DYNAMICS>
480{
481 using namespace gpu::cuda;
482 DeviceContext<T,DESCRIPTOR> lattice(block);
483 if (block.statisticsEnabled()) {
484 call_operators_with_statistics(
485 lattice,
486 subdomain.deviceData(),
487 MaskedCollision<T,DESCRIPTOR,DYNAMICS>{_parameters->parameters, _mask->deviceData()},
488 DynamicDispatchCollision{});
489 } else {
490 call_operators(
491 lattice,
492 subdomain.deviceData(),
493 MaskedCollision<T,DESCRIPTOR,DYNAMICS>{_parameters->parameters, _mask->deviceData()},
494 DynamicDispatchCollision{});
495 }
496}
497
498template <typename T, typename DESCRIPTOR, typename DYNAMICS>
499void ConcreteBlockCollisionO<T,DESCRIPTOR,Platform::GPU_CUDA,DYNAMICS>::applyIndividual(
500 ConcreteBlockLattice<T,DESCRIPTOR,Platform::GPU_CUDA>& block,
501 ConcreteBlockMask<T,Platform::GPU_CUDA>& subdomain)
502{
503 using namespace gpu::cuda;
504 DeviceContext<T,DESCRIPTOR> lattice(block);
505 // Primitive heuristic for preferring mask-based to list-based dispatch
506 if (_mask->weight() > 0.5*subdomain.weight()) {
507 if (block.statisticsEnabled()) {
509 _stream.get(),
510 lattice,
511 subdomain.deviceData(),
512 MaskedCollision<T,DESCRIPTOR,DYNAMICS>{_parameters->parameters, _mask->deviceData()});
513 } else {
515 _stream.get(),
516 lattice,
517 subdomain.deviceData(),
518 MaskedCollision<T,DESCRIPTOR,DYNAMICS>{_parameters->parameters, _mask->deviceData()});
519 }
520
521 // Use list of cell indices
522 } else {
523 // Update cell list from mask
524 if (_modified) {
525 _cells.clear();
526 for (CellID iCell=0; iCell < block.getNcells(); ++iCell) {
527 if (_mask->operator[](iCell)) {
528 _cells.push_back(iCell);
529 }
530 }
531 _cells.setProcessingContext(ProcessingContext::Simulation);
532 _modified = false;
533 }
534
535 if (block.statisticsEnabled()) {
537 _stream.get(),
538 lattice,
539 _cells,
540 ListedCollision<T,DESCRIPTOR,DYNAMICS>{_parameters->parameters});
541 } else {
543 _stream.get(),
544 lattice,
545 _cells,
546 ListedCollision<T,DESCRIPTOR,DYNAMICS>{_parameters->parameters});
547 }
548 }
549}
550
551template <typename T, typename DESCRIPTOR, typename DYNAMICS>
554{
555 // Fetch pointers to DYNAMICS-specific parameter and mask data
556 _parameters = &block.template getData<OperatorParameters<DYNAMICS>>();
557 _mask = &block.template getData<DynamicsMask<DYNAMICS>>();
558
559 {
560 // Construct on-device dynamics proxy for dynamic dispatch
561 _deviceDynamics = gpu::cuda::device::malloc<gpu::cuda::ConcreteDynamics<T,DESCRIPTOR,DYNAMICS>>(1);
562 gpu::cuda::kernel::construct_dynamics<T,DESCRIPTOR,DYNAMICS><<<1,1>>>(
563 _deviceDynamics.get(),
564 _parameters->deviceData());
566
567 // Fetch pointer to on-device dynamic-dispatch field
568 _dynamicsOfCells = block.template getField<gpu::cuda::DYNAMICS<T,DESCRIPTOR>>()[0].data();
569 }
570}
571
572template <typename T, typename DESCRIPTOR, typename OPERATOR>
574 _cells(0),
575 _modified{false},
576 _stream{cudaStreamDefault}
577{ }
578
579template <typename T, typename DESCRIPTOR, typename OPERATOR>
582{
583 if (_cells.size() > 0) {
584 if (_modified) {
585 _cells.deduplicate();
586 _cells.setProcessingContext(ProcessingContext::Simulation);
587 _modified = false;
588 }
591 lattice,
592 _cells,
594 }
595}
596
597template <typename T, typename DESCRIPTOR, typename OPERATOR>
599 _cells(0),
600 _modified{false},
601 _stream{cudaStreamDefault}
602{ }
603
604template <typename T, typename DESCRIPTOR, typename OPERATOR>
607{
608 if (_cells.size() > 0) {
609 if (_modified) {
610 _cells.deduplicate();
611 _cells.setProcessingContext(ProcessingContext::Simulation);
612 _modified = false;
613 }
614 using namespace gpu::cuda;
615 DeviceBlockLattice<T,DESCRIPTOR> lattice(block);
616 async_call_list_operators(_stream.get(),
617 lattice,
618 _cells,
619 ListedPostProcessorWithParameters<T,DESCRIPTOR,OPERATOR>{_parameters->parameters});
620 }
621}
622
623
624template <typename COUPLER, typename COUPLEES>
626 CellID iCell, bool state)
627{
628 if (!_mask) {
629 _mask = std::make_unique<ConcreteBlockMask<typename COUPLEES::values_t::template get<0>::value_t,
631 _lattices.template get<0>()->template getData<CollisionSubdomainMask>()
632 );
633 }
634 _mask->set(iCell, state);
635}
636
637template <typename COUPLER, typename COUPLEES>
639{
640 auto deviceLattice = _lattices.exchange_values([&](auto name) -> auto {
641 return gpu::cuda::DeviceBlockLattice{*_lattices.get(name)};
642 });
643 if (_mask) {
644 _mask->setProcessingContext(ProcessingContext::Simulation);
646 deviceLattice, _mask->deviceData(),
648 } else {
649 auto& mask = _lattices.template get<0>()->template getData<CollisionSubdomainMask>();
651 deviceLattice, mask.deviceData(),
653 }
654}
655
656template <typename COUPLER, typename COUPLEES>
658 CellID iCell, bool state)
659{
660 if (!_mask) {
661 _mask = std::make_unique<ConcreteBlockMask<typename COUPLEES::values_t::template get<0>::value_t,
663 _lattices.template get<0>()->template getData<CollisionSubdomainMask>()
664 );
665 }
666 _mask->set(iCell, state);
667}
668
669template <typename COUPLER, typename COUPLEES>
671{
672 auto deviceLattice = _lattices.exchange_values([&](auto name) -> auto {
673 return gpu::cuda::DeviceBlockLattice{*_lattices.get(name)};
674 });
675 if (_mask) {
676 _mask->setProcessingContext(ProcessingContext::Simulation);
678 deviceLattice, _mask->deviceData(),
680 } else {
681 auto& mask = _lattices.template get<0>()->template getData<CollisionSubdomainMask>();
683 deviceLattice, mask.deviceData(),
685 }
686}
687
688}
689
690#endif
bool statisticsEnabled() const
Collision operation of concrete DYNAMICS on concrete block lattices of PLATFORM.
Coupling of COUPLEES using concrete OPERATOR with SCOPE on PLATFORM lattices.
Definition operator.h:142
Implementation of BlockLattice on a concrete PLATFORM.
Block application of concrete OPERATOR called using SCOPE on PLATFORM.
Definition operator.h:65
Device-side implementation of the Cell concept for post processors.
Definition context.hh:248
Plain column for CUDA GPU targets.
Definition column.h:49
const T * deviceData() const
Definition column.hh:146
std::size_t size() const
Definition column.hh:128
Implementation of gpu::cuda::Dynamics for concrete DYNAMICS.
Definition dynamics.hh:81
Device-side implementation of the data-only Cell concept for collision steps.
Definition context.hh:140
Device-side view of a block lattice.
Definition context.hh:207
Cell< T, DESCRIPTOR > get(CellID iCell) __device__
Definition context.hh:233
Structure for passing pointers to on-device data into CUDA kernels.
Definition context.hh:55
List-based application of DYNAMICS::apply for use in kernel::call_list_operators.
Definition operator.hh:93
bool operator()(DeviceContext< T, DESCRIPTOR > &lattice, CellID iCell) __device__
Definition operator.hh:102
bool operator()(DeviceContext< T, DESCRIPTOR > &lattice, CellID iCell, CellStatistic< T > &statistic) __device__
Definition operator.hh:108
ListedCollision(ParametersOfOperatorD< T, DESCRIPTOR, DYNAMICS > &parameters) __host__
Definition operator.hh:98
List-based application of OPERATOR::apply with parameters.
Definition operator.hh:158
bool operator()(DeviceBlockLattice< T, DESCRIPTOR > &lattice, CellID iCell) __device__
Definition operator.hh:167
ListedPostProcessorWithParameters(ParametersOfOperatorD< T, DESCRIPTOR, OPERATOR > &parameters) __host__
Definition operator.hh:163
Masked application of DYNAMICS::apply for use in kernel::call_operators.
Definition operator.hh:45
bool operator()(DeviceContext< T, DESCRIPTOR > &lattice, CellID iCell) __device__
Chainable call operator for use in kernel::call_operators.
Definition operator.hh:72
MaskedCollision(ParametersOfOperatorD< T, DESCRIPTOR, DYNAMICS > &parameters, bool *mask) any_platform
Constructor (commonly called on the host side)
Definition operator.hh:62
bool operator()(DeviceContext< T, DESCRIPTOR > &lattice, CellID iCell, CellStatistic< T > &statistic) __device__
Chainable call operator with statistics storage.
Definition operator.hh:81
Masked application of OPERATOR::apply.
Definition operator.hh:119
MaskedPostProcessor(bool *mask) any_platform
Definition operator.hh:125
bool operator()(DeviceBlockLattice< T, DESCRIPTOR > &lattice, CellID iCell) __device__
Definition operator.hh:130
Unrestricted application of COUPLING::apply with parameters.
Definition operator.hh:191
bool operator()(CONTEXT &lattices, CellID iCell) __device__
Definition operator.hh:204
UnmaskedCouplingWithParameters(PARAMETERS &parameters) any_platform
Definition operator.hh:199
void check()
Check errors.
Definition device.hh:48
void call_coupling_operators(CONTEXTS lattices, bool *subdomain, OPERATORS... ops) __global__
CUDA kernel for applying UnmaskedCoupling(WithParameters)
Definition operator.hh:341
void call_list_operators(CONTEXT lattice, const CellID *indices, std::size_t nIndices, OPERATORS... ops) __global__
CUDA kernel for applying generic OPERATORS with OperatorScope::PerCell or ListedCollision.
Definition operator.hh:301
void construct_dynamics(void *target, PARAMETERS *parameters) __global__
CUDA kernel for constructing on-device ConcreteDynamics.
Definition operator.hh:352
void call_operators(CONTEXT lattice, bool *subdomain, OPERATORS... ops) __global__
CUDA kernel for applying purely local collision steps.
Definition operator.hh:265
void call_list_operators_with_statistics(CONTEXT lattice, const CellID *indices, std::size_t nIndices, OPERATORS... ops) __global__
CUDA kernel for applying ListedCollision.
Definition operator.hh:316
void call_operators_with_statistics(CONTEXT lattice, bool *subdomain, OPERATORS... ops) __global__
CUDA kernel for applying purely local collision steps while tracking statistics.
Definition operator.hh:278
void async_call_list_operators(cudaStream_t stream, CONTEXT &lattice, const gpu::cuda::Column< CellID > &cells, ARGS &&... args)
Apply operators to listed cell indices (async version)
Definition operator.hh:423
void async_call_operators(cudaStream_t stream, CONTEXT &lattice, bool *subdomain, ARGS &&... args)
Apply masked collision operators to lattice (async)
Definition operator.hh:373
void call_operators(CONTEXT &lattice, bool *subdomain, ARGS &&... args)
Apply masked collision operators to lattice.
Definition operator.hh:363
void async_call_list_operators_with_statistics(cudaStream_t stream, CONTEXT &lattice, const gpu::cuda::Column< CellID > &cells, ARGS &&... args)
Apply ListedCollision with statistics (async version)
Definition operator.hh:438
void async_call_operators_with_statistics(cudaStream_t stream, CONTEXT &lattice, bool *subdomain, ARGS &&... args)
Apply masked collision operators to lattice while tracking statistics (async)
Definition operator.hh:396
void call_list_operators(CONTEXT &lattice, const gpu::cuda::Column< CellID > &cells, ARGS &&... args)
Apply operators to listed cell indices.
Definition operator.hh:409
std::function< void(ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > &) getFusedCollisionO)()
Helper for constructing fused collision operators.
Definition operator.hh:233
void call_operators_with_statistics(CONTEXT &lattice, bool *subdomain, ARGS &&... args)
Apply masked collision operators to lattice while tracking statistics.
Definition operator.hh:386
void call_coupling_operators(CONTEXT &lattices, bool *subdomain, ARGS &&... args)
Apply coupling on subdomain.
Definition operator.hh:453
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.
@ GPU_CUDA
Vector CPU (AVX2 / AVX-512 collision)
typename ParametersD< T, DESCRIPTOR >::template include< typename OPERATOR::parameters > ParametersOfOperatorD
Deduce ParametersD of OPERATOR w.r.t. T and DESCRIPTOR.
#define any_platform
Define preprocessor macros for device-side functions, constant storage.
Definition platform.h:78
Base of block-wide coupling operators executed by SuperLatticeCoupling.
Definition operator.h:115
Return value of any collision.
Definition interface.h:43
Describe mask of DYNAMICS in Data.
Definition data.h:62
On-device field mirroring BlockDynamicsMap.
Definition dynamics.h:39
Last node in a MaskedDynamics chain in kernel::call_operators.
Definition dynamics.hh:148
List-based application of OPERATOR::apply.
Definition operator.hh:146
bool operator()(DeviceBlockLattice< T, DESCRIPTOR > &lattice, CellID iCell) __device__
Definition operator.hh:148
Unrestricted application of COUPLING::apply.
Definition operator.hh:177
bool operator()(CONTEXT &lattices, CellID iCell) __device__
Definition operator.hh:179