ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
EKReactionImplIndexed.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2022-2023 The ESPResSo project
3 *
4 * This file is part of ESPResSo.
5 *
6 * ESPResSo is free software: you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation, either version 3 of the License, or
9 * (at your option) any later version.
10 *
11 * ESPResSo is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program. If not, see <http://www.gnu.org/licenses/>.
18 */
19
20#pragma once
21
23
29
30#include <blockforest/StructuredBlockForest.h>
31#include <domain_decomposition/BlockDataID.h>
32#include <domain_decomposition/IBlock.h>
33#include <field/AddToStorage.h>
34#include <field/FlagField.h>
35#include <field/FlagUID.h>
36
37#include <utils/Vector.hpp>
38
39#include <cassert>
40#include <cstddef>
41#include <memory>
42#include <optional>
43#include <tuple>
44#include <vector>
45
46namespace walberla {
47template <lbmpy::Arch Architecture = lbmpy::Arch::CPU>
49private:
50 BlockDataID m_flagfield_id;
51 BlockDataID m_indexvector_id;
52 bool m_pending_changes;
53
54public:
55 /** Flag for domain cells, i.e. all cells. */
56 FlagUID const Domain_flag{"domain"};
57 /** Flag for boundary cells. */
58 FlagUID const Boundary_flag{"boundary"};
59
60 using FlagField = field::FlagField<uint8_t>;
61#if defined(__CUDACC__)
62 using IndexVectors =
63 std::conditional<Architecture == lbmpy::Arch::CPU,
64 detail::ReactionKernelIndexedSelector::KernelTrait<>::
65 ReactionKernelIndexed::IndexVectors,
66 detail::ReactionKernelIndexedSelector::KernelTraitGPU<>::
67 ReactionKernelIndexedGPU::IndexVectors>::type;
68 using IndexInfo =
69 std::conditional<Architecture == lbmpy::Arch::CPU,
70 detail::ReactionKernelIndexedSelector::KernelTrait<>::
71 ReactionKernelIndexed::IndexInfo,
72 detail::ReactionKernelIndexedSelector::KernelTraitGPU<>::
73 ReactionKernelIndexedGPU::IndexInfo>::type;
74#else
75 using IndexVectors = detail::ReactionKernelIndexedSelector::KernelTrait<>::
76 ReactionKernelIndexed::IndexVectors;
77 using IndexInfo = detail::ReactionKernelIndexedSelector::KernelTrait<>::
78 ReactionKernelIndexed::IndexInfo;
79#endif
80
81private:
82 auto get_flag_field_and_flag(IBlock *block, BlockDataID const &flagfield_id) {
83 auto const flag_field =
84 block->template uncheckedFastGetData<FlagField>(flagfield_id);
85 auto const boundary_flag = flag_field->getFlag(Boundary_flag);
86 return std::make_tuple(flag_field, boundary_flag);
87 }
88
89public:
90 EKReactionImplIndexed(std::shared_ptr<LatticeWalberla> const &lattice,
91 reactants_type const &reactants, double coefficient)
92 : EKReactionBaseIndexed(lattice, reactants, coefficient),
93 m_pending_changes(false) {
94 m_flagfield_id = field::addFlagFieldToStorage<FlagField>(
95 get_lattice()->get_blocks(), "flag field reaction",
96 get_lattice()->get_ghost_layers());
97
98 auto createIdxVector = [](IBlock *const, StructuredBlockStorage *const) {
99 return new IndexVectors();
100 };
101 m_indexvector_id = get_lattice()
102 ->get_blocks()
103 ->template addStructuredBlockData<IndexVectors>(
104 createIdxVector, "IndexField");
105
106 for (auto &block : *get_lattice()->get_blocks()) {
107 auto flag_field = block.template getData<FlagField>(m_flagfield_id);
108 // register flags
109 flag_field->registerFlag(Domain_flag);
110 flag_field->registerFlag(Boundary_flag);
111 // mark all cells as domain cells and fluid cells
112 auto domain_flag = flag_field->getFlag(Domain_flag);
113 auto boundary_flag = flag_field->getFlag(Boundary_flag);
114 for (auto it = flag_field->begin(); it != flag_field->end(); ++it) {
115 flag_field->addFlag(it.x(), it.y(), it.z(), domain_flag);
116 flag_field->removeFlag(it.x(), it.y(), it.z(), boundary_flag);
117 }
118 }
119 }
120 ~EKReactionImplIndexed() override = default;
121
125
126 void perform_reaction() override {
128 std::function<void(IBlock *)> kernel;
129 if (Architecture == lbmpy::Arch::CPU) {
130 kernel = detail::ReactionKernelIndexedSelector::get_kernel(
131 get_reactants(), get_coefficient(), m_indexvector_id);
132 } else {
133#if defined(__CUDACC__)
134 kernel = detail::ReactionKernelIndexedSelector::get_kernel_gpu(
135 get_reactants(), get_coefficient(), m_indexvector_id);
136#endif
137 }
138 for (auto &block : *get_lattice()->get_blocks()) {
139 kernel(&block);
140 }
141 }
142
144 bool is_boundary) override {
145 if (auto bc = get_block_and_cell(*get_lattice(), node, true)) {
146 auto const [flag_field, boundary_flag] =
147 get_flag_field_and_flag(bc->block, m_flagfield_id);
148 if (is_boundary) {
149 flag_field->addFlag(bc->cell, boundary_flag);
150 } else {
151 flag_field->removeFlag(bc->cell, boundary_flag);
152 }
153 m_pending_changes = true;
154 }
155 }
156
157 [[nodiscard]] std::optional<bool>
158 get_node_is_boundary(Utils::Vector3i const &node) override {
159 if (auto bc = get_block_and_cell(*get_lattice(), node, true)) {
160 auto const [flag_field, boundary_flag] =
161 get_flag_field_and_flag(bc->block, m_flagfield_id);
162 return {flag_field->isFlagSet(bc->cell, boundary_flag)};
163 }
164 return std::nullopt;
165 }
166
168 if (m_pending_changes) {
169 for (auto &block : *get_lattice()->get_blocks()) {
170 fillFromFlagField(block);
171 }
172 m_pending_changes = false;
173 }
174 }
175
176private:
177 void fillFromFlagField(IBlock &block) {
178 auto *indexVectors =
179 block.uncheckedFastGetData<IndexVectors>(m_indexvector_id);
180 auto &indexVectorAll = indexVectors->indexVector(IndexVectors::ALL);
181 auto &indexVectorInner = indexVectors->indexVector(IndexVectors::INNER);
182 auto &indexVectorOuter = indexVectors->indexVector(IndexVectors::OUTER);
183
184 auto *flagField = block.getData<FlagField>(m_flagfield_id);
185
186 assert(flagField->flagExists(Boundary_flag) and
187 flagField->flagExists(Domain_flag));
188
189 auto boundaryFlag = flagField->getFlag(Boundary_flag);
190 auto domainFlag = flagField->getFlag(Domain_flag);
191
192 auto inner = flagField->xyzSize();
193 inner.expand(cell_idx_t(-1));
194
195 indexVectorAll.clear();
196 indexVectorInner.clear();
197 indexVectorOuter.clear();
198
199 auto flagWithGLayers = flagField->xyzSizeWithGhostLayer();
200 for (auto it = flagField->beginWithGhostLayerXYZ(); it != flagField->end();
201 ++it) {
202 if (!isFlagSet(it, boundaryFlag))
203 continue;
204
205 if (flagWithGLayers.contains(it.x(), it.y(), it.z()) &&
206 isFlagSet(it.neighbor(0, 0, 0, 0), domainFlag)) {
207 auto element = IndexInfo(it.x(), it.y(), it.z());
208 indexVectorAll.push_back(element);
209 if (inner.contains(it.x(), it.y(), it.z())) {
210 indexVectorInner.push_back(element);
211 } else {
212 indexVectorOuter.push_back(element);
213 }
214 }
215 }
216
217 indexVectors->syncGPU();
218 }
219};
220
221} // namespace walberla
Vector implementation and trait types for boost qvm interoperability.
std::vector< std::shared_ptr< EKReactant > > reactants_type
auto get_lattice() const noexcept
double get_coefficient() const noexcept
auto const & get_reactants() const noexcept
std::optional< bool > get_node_is_boundary(Utils::Vector3i const &node) override
FlagUID const Boundary_flag
Flag for boundary cells.
FlagUID const Domain_flag
Flag for domain cells, i.e.
field::FlagField< uint8_t > FlagField
detail::ReactionKernelIndexedSelector::KernelTrait<>::ReactionKernelIndexed::IndexInfo IndexInfo
~EKReactionImplIndexed() override=default
detail::ReactionKernelIndexedSelector::KernelTrait<>::ReactionKernelIndexed::IndexVectors IndexVectors
EKReactionImplIndexed(std::shared_ptr< LatticeWalberla > const &lattice, reactants_type const &reactants, double coefficient)
void set_node_is_boundary(Utils::Vector3i const &node, bool is_boundary) override
static double * block(double *p, std::size_t index, std::size_t size)
Definition elc.cpp:176
\file PackInfoPdfDoublePrecision.cpp \author pystencils
std::optional< BlockAndCell > get_block_and_cell(::LatticeWalberla const &lattice, signed_integral_vector auto const &node, bool consider_ghost_layers)