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.4+1.ge851f4e, lbmpy v1.4+1.ge9efe34,
21// sympy v1.12.1, lbmpy_walberla/pystencils_walberla from waLBerla commit
22// 007e77e077ad9d22b5eed6f3d3118240993e553c
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#else
42#define RESTRICT
43#endif
44
45#ifdef WALBERLA_BUILD_WITH_HALF_PRECISION_SUPPORT
46using walberla::half;
47#endif
48
49namespace walberla {
50namespace pystencils {
51
53public:
54 struct IndexInfo {
55 int32_t x;
56 int32_t y;
57 int32_t z;
58 int32_t dir;
59 float value;
60 IndexInfo(int32_t x_, int32_t y_, int32_t z_, int32_t dir_)
61 : x(x_), y(y_), z(z_), dir(dir_), value() {}
62 bool operator==(const IndexInfo &o) const {
63 return x == o.x && y == o.y && z == o.z && dir == o.dir &&
64 floatIsEqual(value, o.value);
65 }
66 };
67
69 public:
70 using CpuIndexVector = std::vector<IndexInfo>;
71
72 enum Type { ALL = 0, INNER = 1, OUTER = 2, NUM_TYPES = 3 };
73
74 IndexVectors() = default;
75 bool operator==(IndexVectors const &other) const {
76 return other.cpuVectors_ == cpuVectors_;
77 }
78
79 auto &indexVector(Type t) { return cpuVectors_[t]; }
80 auto const &indexVector(Type t) const { return cpuVectors_[t]; }
82 return cpuVectors_[t].empty() ? nullptr : cpuVectors_[t].data();
83 }
84
85 void syncGPU() {}
86
87 private:
88 std::vector<CpuIndexVector> cpuVectors_{NUM_TYPES};
89 };
90
92 const std::shared_ptr<StructuredBlockForest> &blocks,
93 BlockDataID fieldID_,
94 std::function<float(const Cell &,
95 const shared_ptr<StructuredBlockForest> &, IBlock &)>
96 &dirichletCallback)
97 : elementInitaliser(dirichletCallback), fieldID(fieldID_) {
98 auto createIdxVector = [](IBlock *const, StructuredBlockStorage *const) {
99 return new IndexVectors();
100 };
101 indexVectorID = blocks->addStructuredBlockData<IndexVectors>(
102 createIdxVector, "IndexField_Dirichlet_single_precision");
103 }
104
105 void run(IBlock *block);
106
107 void operator()(IBlock *block) { run(block); }
108
109 void inner(IBlock *block);
110
111 void outer(IBlock *block);
112
113 Vector3<float> getForce(IBlock * /*block*/) {
114
115 WALBERLA_ABORT(
116 "Boundary condition was not generated including force calculation.")
117 return Vector3<float>(float_c(0.0));
118 }
119
120 std::function<void(IBlock *)> getSweep() {
121 return [this](IBlock *b) { this->run(b); };
122 }
123
124 std::function<void(IBlock *)> getInnerSweep() {
125 return [this](IBlock *b) { this->inner(b); };
126 }
127
128 std::function<void(IBlock *)> getOuterSweep() {
129 return [this](IBlock *b) { this->outer(b); };
130 }
131
132 template <typename FlagField_T>
133 void fillFromFlagField(const std::shared_ptr<StructuredBlockForest> &blocks,
134 ConstBlockDataID flagFieldID, FlagUID boundaryFlagUID,
135 FlagUID domainFlagUID) {
136 for (auto &block : *blocks)
137 fillFromFlagField<FlagField_T>(blocks, &block, flagFieldID,
138 boundaryFlagUID, domainFlagUID);
139 }
140
141 template <typename FlagField_T>
142 void fillFromFlagField(const shared_ptr<StructuredBlockForest> &blocks,
143 IBlock *block, ConstBlockDataID flagFieldID,
144 FlagUID boundaryFlagUID, FlagUID domainFlagUID) {
145 auto *indexVectors = block->getData<IndexVectors>(indexVectorID);
146 auto &indexVectorAll = indexVectors->indexVector(IndexVectors::ALL);
147 auto &indexVectorInner = indexVectors->indexVector(IndexVectors::INNER);
148 auto &indexVectorOuter = indexVectors->indexVector(IndexVectors::OUTER);
149
150 auto *flagField = block->getData<FlagField_T>(flagFieldID);
151
152 if (!(flagField->flagExists(boundaryFlagUID) &&
153 flagField->flagExists(domainFlagUID)))
154 return;
155
156 auto boundaryFlag = flagField->getFlag(boundaryFlagUID);
157 auto domainFlag = flagField->getFlag(domainFlagUID);
158
159 auto inner = flagField->xyzSize();
160 inner.expand(cell_idx_t(-1));
161
162 indexVectorAll.clear();
163 indexVectorInner.clear();
164 indexVectorOuter.clear();
165
166 auto flagWithGLayers = flagField->xyzSizeWithGhostLayer();
167 for (auto it = flagField->beginWithGhostLayerXYZ(); it != flagField->end();
168 ++it) {
169
170 if (!isFlagSet(it, boundaryFlag))
171 continue;
172 if (flagWithGLayers.contains(it.x() + cell_idx_c(0),
173 it.y() + cell_idx_c(0),
174 it.z() + cell_idx_c(0)) &&
175 isFlagSet(it.neighbor(0, 0, 0, 0), domainFlag)) {
176
177 auto element = IndexInfo(it.x(), it.y(), it.z(), 0);
178 float InitialisatonAdditionalData =
179 elementInitaliser(Cell(it.x(), it.y(), it.z()), blocks, *block);
180 element.value = InitialisatonAdditionalData;
181 indexVectorAll.emplace_back(element);
182 if (inner.contains(it.x(), it.y(), it.z()))
183 indexVectorInner.emplace_back(element);
184 else
185 indexVectorOuter.emplace_back(element);
186 }
187 }
188
189 indexVectors->syncGPU();
190 }
191
192private:
193 void run_impl(IBlock *block, IndexVectors::Type type);
194
195 BlockDataID indexVectorID;
196
197 std::function<float(const Cell &, const shared_ptr<StructuredBlockForest> &,
198 IBlock &)>
199 elementInitaliser;
200
201public:
202 BlockDataID fieldID;
203};
204
205} // namespace pystencils
206} // 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_)