ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
Dirichlet_single_precision.h
Go to the documentation of this file.
1//======================================================================================================================
2//
3// This file is part of waLBerla. waLBerla is free software: you can
4// redistribute it and/or modify it under the terms of the GNU General Public
5// License as published by the Free Software Foundation, either version 3 of
6// the License, or (at your option) any later version.
7//
8// waLBerla is distributed in the hope that it will be useful, but WITHOUT
9// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
10// FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
11// for more details.
12//
13// You should have received a copy of the GNU General Public License along
14// with waLBerla (see COPYING.txt). If not, see <http://www.gnu.org/licenses/>.
15//
16//! \\file Dirichlet_single_precision.h
17//! \\author pystencils
18//======================================================================================================================
19
20// kernel generated with pystencils v1.3.7+13.gdfd203a, lbmpy
21// v1.3.7+10.gd3f6236, sympy v1.12.1, lbmpy_walberla/pystencils_walberla from
22// waLBerla commit 191cf58b16b96d1d2f050dcbd9e88443995b2222
23
24#pragma once
25#include "core/DataTypes.h"
26#include "core/logging/Logging.h"
27
28#include "blockforest/StructuredBlockForest.h"
29#include "core/debug/Debug.h"
30#include "domain_decomposition/BlockDataID.h"
31#include "domain_decomposition/IBlock.h"
32#include "field/FlagField.h"
33#include "field/GhostLayerField.h"
34
35#include <functional>
36#include <memory>
37#include <vector>
38
39#ifdef __GNUC__
40#define RESTRICT __restrict__
41#elif _MSC_VER
42#define RESTRICT __restrict
43#else
44#define RESTRICT
45#endif
46
47#ifdef WALBERLA_BUILD_WITH_HALF_PRECISION_SUPPORT
48using walberla::half;
49#endif
50
51namespace walberla {
52namespace pystencils {
53
55public:
56 struct IndexInfo {
57 int32_t x;
58 int32_t y;
59 int32_t z;
60 int32_t dir;
61 float value;
62 IndexInfo(int32_t x_, int32_t y_, int32_t z_, int32_t dir_)
63 : x(x_), y(y_), z(z_), dir(dir_), value() {}
64 bool operator==(const IndexInfo &o) const {
65 return x == o.x && y == o.y && z == o.z && dir == o.dir &&
66 floatIsEqual(value, o.value);
67 }
68 };
69
71 public:
72 using CpuIndexVector = std::vector<IndexInfo>;
73
74 enum Type { ALL = 0, INNER = 1, OUTER = 2, NUM_TYPES = 3 };
75
76 IndexVectors() = default;
77 bool operator==(IndexVectors const &other) const {
78 return other.cpuVectors_ == cpuVectors_;
79 }
80
81 auto &indexVector(Type t) { return cpuVectors_[t]; }
82 auto const &indexVector(Type t) const { return cpuVectors_[t]; }
84 return cpuVectors_[t].empty() ? nullptr : cpuVectors_[t].data();
85 }
86
87 void syncGPU() {}
88
89 private:
90 std::vector<CpuIndexVector> cpuVectors_{NUM_TYPES};
91 };
92
94 const std::shared_ptr<StructuredBlockForest> &blocks,
95 BlockDataID fieldID_,
96 std::function<float(const Cell &,
97 const shared_ptr<StructuredBlockForest> &, IBlock &)>
98 &dirichletCallback)
99 : elementInitaliser(dirichletCallback), fieldID(fieldID_) {
100 auto createIdxVector = [](IBlock *const, StructuredBlockStorage *const) {
101 return new IndexVectors();
102 };
103 indexVectorID = blocks->addStructuredBlockData<IndexVectors>(
104 createIdxVector, "IndexField_Dirichlet_single_precision");
105 }
106
107 void run(IBlock *block);
108
109 void operator()(IBlock *block) { run(block); }
110
111 void inner(IBlock *block);
112
113 void outer(IBlock *block);
114
115 Vector3<float> getForce(IBlock * /*block*/) {
116
117 WALBERLA_ABORT(
118 "Boundary condition was not generated including force calculation.")
119 return Vector3<float>(float_c(0.0));
120 }
121
122 std::function<void(IBlock *)> getSweep() {
123 return [this](IBlock *b) { this->run(b); };
124 }
125
126 std::function<void(IBlock *)> getInnerSweep() {
127 return [this](IBlock *b) { this->inner(b); };
128 }
129
130 std::function<void(IBlock *)> getOuterSweep() {
131 return [this](IBlock *b) { this->outer(b); };
132 }
133
134 template <typename FlagField_T>
135 void fillFromFlagField(const std::shared_ptr<StructuredBlockForest> &blocks,
136 ConstBlockDataID flagFieldID, FlagUID boundaryFlagUID,
137 FlagUID domainFlagUID) {
138 for (auto blockIt = blocks->begin(); blockIt != blocks->end(); ++blockIt)
139 fillFromFlagField<FlagField_T>(blocks, &*blockIt, flagFieldID,
140 boundaryFlagUID, domainFlagUID);
141 }
142
143 template <typename FlagField_T>
144 void fillFromFlagField(const shared_ptr<StructuredBlockForest> &blocks,
145 IBlock *block, ConstBlockDataID flagFieldID,
146 FlagUID boundaryFlagUID, FlagUID domainFlagUID) {
147 auto *indexVectors = block->getData<IndexVectors>(indexVectorID);
148 auto &indexVectorAll = indexVectors->indexVector(IndexVectors::ALL);
149 auto &indexVectorInner = indexVectors->indexVector(IndexVectors::INNER);
150 auto &indexVectorOuter = indexVectors->indexVector(IndexVectors::OUTER);
151
152 auto *flagField = block->getData<FlagField_T>(flagFieldID);
153
154 if (!(flagField->flagExists(boundaryFlagUID) &&
155 flagField->flagExists(domainFlagUID)))
156 return;
157
158 auto boundaryFlag = flagField->getFlag(boundaryFlagUID);
159 auto domainFlag = flagField->getFlag(domainFlagUID);
160
161 auto inner = flagField->xyzSize();
162 inner.expand(cell_idx_t(-1));
163
164 indexVectorAll.clear();
165 indexVectorInner.clear();
166 indexVectorOuter.clear();
167
168 auto flagWithGLayers = flagField->xyzSizeWithGhostLayer();
169 for (auto it = flagField->beginWithGhostLayerXYZ(); it != flagField->end();
170 ++it) {
171
172 if (!isFlagSet(it, boundaryFlag))
173 continue;
174 if (flagWithGLayers.contains(it.x() + cell_idx_c(0),
175 it.y() + cell_idx_c(0),
176 it.z() + cell_idx_c(0)) &&
177 isFlagSet(it.neighbor(0, 0, 0, 0), domainFlag)) {
178
179 auto element = IndexInfo(it.x(), it.y(), it.z(), 0);
180 float InitialisatonAdditionalData =
181 elementInitaliser(Cell(it.x(), it.y(), it.z()), blocks, *block);
182 element.value = InitialisatonAdditionalData;
183 indexVectorAll.emplace_back(element);
184 if (inner.contains(it.x(), it.y(), it.z()))
185 indexVectorInner.emplace_back(element);
186 else
187 indexVectorOuter.emplace_back(element);
188 }
189 }
190
191 indexVectors->syncGPU();
192 }
193
194private:
195 void run_impl(IBlock *block, IndexVectors::Type type);
196
197 BlockDataID indexVectorID;
198
199 std::function<float(const Cell &, const shared_ptr<StructuredBlockForest> &,
200 IBlock &)>
201 elementInitaliser;
202
203public:
204 BlockDataID fieldID;
205};
206
207} // namespace pystencils
208} // namespace walberla
Definition Cell.hpp:96
void fillFromFlagField(const shared_ptr< StructuredBlockForest > &blocks, IBlock *block, ConstBlockDataID flagFieldID, FlagUID boundaryFlagUID, FlagUID domainFlagUID)
Dirichlet_single_precision(const std::shared_ptr< StructuredBlockForest > &blocks, BlockDataID fieldID_, std::function< float(const Cell &, const shared_ptr< StructuredBlockForest > &, IBlock &)> &dirichletCallback)
void fillFromFlagField(const std::shared_ptr< StructuredBlockForest > &blocks, ConstBlockDataID flagFieldID, FlagUID boundaryFlagUID, FlagUID domainFlagUID)
static double * block(double *p, std::size_t index, std::size_t size)
Definition elc.cpp:176
\file PackInfoPdfDoublePrecision.cpp \author pystencils
IndexInfo(int32_t x_, int32_t y_, int32_t z_, int32_t dir_)