ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
BoundaryPackInfo.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2024 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
22#include <core/debug/Debug.h>
23#include <core/mpi/RecvBuffer.h>
24#include <core/mpi/SendBuffer.h>
25#include <domain_decomposition/IBlock.h>
26#include <field/FlagUID.h>
27#include <field/communication/PackInfo.h>
28#include <stencil/Directions.h>
29
30#include <memory>
31#include <tuple>
32#include <utility>
33
34namespace walberla {
35namespace field {
36namespace communication {
37
38template <typename GhostLayerField_T>
39class BoundaryFlagPackInfo : public PackInfo<GhostLayerField_T> {
40
41public:
42 using PackInfo<GhostLayerField_T>::PackInfo;
43 using PackInfo<GhostLayerField_T>::numberOfGhostLayersToCommunicate;
44
45 ~BoundaryFlagPackInfo() override = default;
46
47 bool constantDataExchange() const override { return false; }
48 bool threadsafeReceiving() const override { return false; }
49};
50
51template <typename GhostLayerField_T, typename Boundary_T>
52class BoundaryPackInfo : public PackInfo<GhostLayerField_T> {
53protected:
54 using PackInfo<GhostLayerField_T>::bdId_;
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
60public:
61 using PackInfo<GhostLayerField_T>::PackInfo;
62 using PackInfo<GhostLayerField_T>::numberOfGhostLayersToCommunicate;
63
64 ~BoundaryPackInfo() override = default;
65
66 void setup_boundary_handle(std::shared_ptr<LatticeWalberla> lattice,
67 std::shared_ptr<Boundary_T> boundary) {
68 m_lattice = std::move(lattice);
69 m_boundary = std::move(boundary);
70 }
71
72 bool constantDataExchange() const override { return false; }
73 bool threadsafeReceiving() const override { return false; }
74
75 void communicateLocal(IBlock const *sender, IBlock *receiver,
76 stencil::Direction dir) override {
77#ifdef _OPENMP
78#pragma omp critical
79#endif
80 {
81 mpi::SendBuffer sBuffer;
82 packDataImpl(sender, dir, sBuffer);
83 mpi::RecvBuffer rBuffer(sBuffer);
84 unpackData(receiver, stencil::inverseDir[dir], rBuffer);
85 }
86 }
87
88 void unpackData(IBlock *receiver, stencil::Direction dir,
89 mpi::RecvBuffer &buffer) override {
90
91 auto *flag_field = receiver->getData<GhostLayerField_T>(bdId_);
92 WALBERLA_ASSERT_NOT_NULLPTR(flag_field);
93 WALBERLA_ASSERT_NOT_NULLPTR(m_boundary);
94 WALBERLA_ASSERT_NOT_NULLPTR(m_lattice);
95
96 auto const boundary_flag = flag_field->getFlag(Boundary_flag);
97 auto const gl = numberOfGhostLayersToCommunicate(flag_field);
98 auto const begin = [gl, dir](auto const *field) {
99 return field->beginGhostLayerOnly(gl, dir);
100 };
101
102#ifndef NDEBUG
103 uint_t xSize, ySize, zSize, bSize;
104 buffer >> xSize >> ySize >> zSize >> bSize;
105 uint_t buf_size{0u};
106 for (auto it = begin(flag_field); it != flag_field->end(); ++it) {
107 if (isFlagSet(it, boundary_flag)) {
108 ++buf_size;
109 }
110 }
111 WALBERLA_ASSERT_EQUAL(xSize, flag_field->xSize());
112 WALBERLA_ASSERT_EQUAL(ySize, flag_field->ySize());
113 WALBERLA_ASSERT_EQUAL(zSize, flag_field->zSize());
114 WALBERLA_ASSERT_EQUAL(bSize, buf_size);
115#endif
116
117 auto const offset = m_lattice->get_block_corner(*receiver, true);
118 typename Boundary_T::value_type value;
119 for (auto it = begin(flag_field); it != flag_field->end(); ++it) {
120 if (isFlagSet(it, boundary_flag)) {
121 auto const node = offset + Utils::Vector3i{{it.x(), it.y(), it.z()}};
122 buffer >> value;
123 m_boundary->unpack_node(node, value);
124 }
125 }
126 }
127
128protected:
129 void packDataImpl(IBlock const *sender, stencil::Direction dir,
130 mpi::SendBuffer &buffer) const override {
131
132 auto const *flag_field = sender->getData<GhostLayerField_T>(bdId_);
133 WALBERLA_ASSERT_NOT_NULLPTR(flag_field);
134 WALBERLA_ASSERT_NOT_NULLPTR(m_boundary);
135 WALBERLA_ASSERT_NOT_NULLPTR(m_lattice);
136
137 auto const boundary_flag = flag_field->getFlag(Boundary_flag);
138 auto const gl = numberOfGhostLayersToCommunicate(flag_field);
139 auto const begin = [gl, dir](auto const *field) {
140 return field->beginSliceBeforeGhostLayer(dir, cell_idx_c(gl));
141 };
142
143#ifndef NDEBUG
144 uint_t buf_size{0u};
145 for (auto it = begin(flag_field); it != flag_field->end(); ++it) {
146 if (isFlagSet(it, boundary_flag)) {
147 ++buf_size;
148 }
149 }
150 buffer << flag_field->xSize() << flag_field->ySize() << flag_field->zSize()
151 << buf_size;
152#endif
153
154 auto const offset = m_lattice->get_block_corner(*sender, true);
155 for (auto it = begin(flag_field); it != flag_field->end(); ++it) {
156 if (isFlagSet(it, boundary_flag)) {
157 auto const node = offset + Utils::Vector3i{{it.x(), it.y(), it.z()}};
158 buffer << m_boundary->get_node_value_at_boundary(node);
159 }
160 }
161 }
162
163private:
164 std::shared_ptr<LatticeWalberla> m_lattice;
165 std::shared_ptr<Boundary_T> m_boundary;
166};
167
168} // namespace communication
169} // namespace field
170} // namespace walberla
void setup_boundary_handle(std::shared_ptr< LatticeWalberla > lattice, std::shared_ptr< Boundary_T > boundary)
void unpackData(IBlock *receiver, stencil::Direction dir, mpi::RecvBuffer &buffer) override
void packDataImpl(IBlock const *sender, stencil::Direction dir, mpi::SendBuffer &buffer) const override
void communicateLocal(IBlock const *sender, IBlock *receiver, stencil::Direction dir) override
FlagUID const Domain_flag
Flag for domain cells, i.e.
FlagUID const Boundary_flag
Flag for boundary cells.
\file PackInfoPdfDoublePrecision.cpp \author pystencils