ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
Dirichlet_double_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_double_precision.h
17//! \\author pystencils
18//======================================================================================================================
19
20// kernel generated with pystencils v1.3.3, lbmpy v1.3.3,
21// lbmpy_walberla/pystencils_walberla from waLBerla commit
22// b0842e1a493ce19ef1bbb8d2cf382fc343970a7f
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 <set>
36#include <vector>
37
38#ifdef __GNUC__
39#define RESTRICT __restrict__
40#elif _MSC_VER
41#define RESTRICT __restrict
42#else
43#define RESTRICT
44#endif
45
46#ifdef WALBERLA_BUILD_WITH_HALF_PRECISION_SUPPORT
47using walberla::half;
48#endif
49
50namespace walberla {
51namespace pystencils {
52
54public:
55 struct IndexInfo {
56 int32_t x;
57 int32_t y;
58 int32_t z;
59 int32_t dir;
60 double value;
61 IndexInfo(int32_t x_, int32_t y_, int32_t z_, int32_t dir_)
62 : x(x_), y(y_), z(z_), dir(dir_), value() {}
63 bool operator==(const IndexInfo &o) const {
64 return x == o.x && y == o.y && z == o.z && dir == o.dir &&
65 floatIsEqual(value, o.value);
66 }
67 };
68
70 public:
71 using CpuIndexVector = std::vector<IndexInfo>;
72
73 enum Type { ALL = 0, INNER = 1, OUTER = 2, NUM_TYPES = 3 };
74
75 IndexVectors() = default;
76 bool operator==(IndexVectors const &other) const {
77 return other.cpuVectors_ == cpuVectors_;
78 }
79
80 CpuIndexVector &indexVector(Type t) { return cpuVectors_[t]; }
81 IndexInfo *pointerCpu(Type t) { return cpuVectors_[t].data(); }
82
83 void syncGPU() {}
84
85 private:
86 std::vector<CpuIndexVector> cpuVectors_{NUM_TYPES};
87 };
88
90 const shared_ptr<StructuredBlockForest> &blocks, BlockDataID fieldID_,
91 std::function<double(const Cell &,
92 const shared_ptr<StructuredBlockForest> &, IBlock &)>
93 &dirichletCallback)
94 : elementInitaliser(dirichletCallback), fieldID(fieldID_) {
95 auto createIdxVector = [](IBlock *const, StructuredBlockStorage *const) {
96 return new IndexVectors();
97 };
98 indexVectorID = blocks->addStructuredBlockData<IndexVectors>(
99 createIdxVector, "IndexField_Dirichlet_double_precision");
100 };
101
102 void run(IBlock *block);
103
104 void operator()(IBlock *block) { run(block); }
105
106 void inner(IBlock *block);
107
108 void outer(IBlock *block);
109
110 std::function<void(IBlock *)> getSweep() {
111 return [this](IBlock *b) { this->run(b); };
112 }
113
114 std::function<void(IBlock *)> getInnerSweep() {
115 return [this](IBlock *b) { this->inner(b); };
116 }
117
118 std::function<void(IBlock *)> getOuterSweep() {
119 return [this](IBlock *b) { this->outer(b); };
120 }
121
122 template <typename FlagField_T>
123 void fillFromFlagField(const shared_ptr<StructuredBlockForest> &blocks,
124 ConstBlockDataID flagFieldID, FlagUID boundaryFlagUID,
125 FlagUID domainFlagUID) {
126 for (auto blockIt = blocks->begin(); blockIt != blocks->end(); ++blockIt)
127 fillFromFlagField<FlagField_T>(blocks, &*blockIt, flagFieldID,
128 boundaryFlagUID, domainFlagUID);
129 }
130
131 template <typename FlagField_T>
132 void fillFromFlagField(const shared_ptr<StructuredBlockForest> &blocks,
133 IBlock *block, ConstBlockDataID flagFieldID,
134 FlagUID boundaryFlagUID, FlagUID domainFlagUID) {
135 auto *indexVectors = block->getData<IndexVectors>(indexVectorID);
136 auto &indexVectorAll = indexVectors->indexVector(IndexVectors::ALL);
137 auto &indexVectorInner = indexVectors->indexVector(IndexVectors::INNER);
138 auto &indexVectorOuter = indexVectors->indexVector(IndexVectors::OUTER);
139
140 auto *flagField = block->getData<FlagField_T>(flagFieldID);
141
142 if (!(flagField->flagExists(boundaryFlagUID) &&
143 flagField->flagExists(domainFlagUID)))
144 return;
145
146 auto boundaryFlag = flagField->getFlag(boundaryFlagUID);
147 auto domainFlag = flagField->getFlag(domainFlagUID);
148
149 auto inner = flagField->xyzSize();
150 inner.expand(cell_idx_t(-1));
151
152 indexVectorAll.clear();
153 indexVectorInner.clear();
154 indexVectorOuter.clear();
155
156 auto flagWithGLayers = flagField->xyzSizeWithGhostLayer();
157 for (auto it = flagField->beginWithGhostLayerXYZ(); it != flagField->end();
158 ++it) {
159
160 if (!isFlagSet(it, boundaryFlag))
161 continue;
162 if (flagWithGLayers.contains(it.x() + cell_idx_c(0),
163 it.y() + cell_idx_c(0),
164 it.z() + cell_idx_c(0)) &&
165 isFlagSet(it.neighbor(0, 0, 0, 0), domainFlag)) {
166
167 auto element = IndexInfo(it.x(), it.y(), it.z(), 0);
168 double InitialisatonAdditionalData =
169 elementInitaliser(Cell(it.x(), it.y(), it.z()), blocks, *block);
170 element.value = InitialisatonAdditionalData;
171 indexVectorAll.push_back(element);
172 if (inner.contains(it.x(), it.y(), it.z()))
173 indexVectorInner.push_back(element);
174 else
175 indexVectorOuter.push_back(element);
176 }
177 }
178
179 indexVectors->syncGPU();
180 }
181
182private:
183 void run_impl(IBlock *block, IndexVectors::Type type);
184
185 BlockDataID indexVectorID;
186 std::function<double(const Cell &, const shared_ptr<StructuredBlockForest> &,
187 IBlock &)>
188 elementInitaliser;
189
190public:
191 BlockDataID fieldID;
192};
193
194} // namespace pystencils
195} // namespace walberla
Definition Cell.hpp:97
void fillFromFlagField(const shared_ptr< StructuredBlockForest > &blocks, IBlock *block, ConstBlockDataID flagFieldID, FlagUID boundaryFlagUID, FlagUID domainFlagUID)
Dirichlet_double_precision(const shared_ptr< StructuredBlockForest > &blocks, BlockDataID fieldID_, std::function< double(const Cell &, const shared_ptr< StructuredBlockForest > &, IBlock &)> &dirichletCallback)
void fillFromFlagField(const 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:172
\file PackInfoPdfDoublePrecision.cpp \author pystencils
IndexInfo(int32_t x_, int32_t y_, int32_t z_, int32_t dir_)