Loading [MathJax]/extensions/TeX/AMSmath.js
ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages Concepts
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, typename Boundary_T>
39class BoundaryPackInfo : public PackInfo<GhostLayerField_T> {
40protected:
41 using PackInfo<GhostLayerField_T>::bdId_;
42 /** Flag for domain cells, i.e. all cells. */
43 FlagUID const Domain_flag{"domain"};
44 /** Flag for boundary cells. */
45 FlagUID const Boundary_flag{"boundary"};
46
47public:
48 using PackInfo<GhostLayerField_T>::PackInfo;
49 using PackInfo<GhostLayerField_T>::numberOfGhostLayersToCommunicate;
50
51 ~BoundaryPackInfo() override = default;
52
53 void setup_boundary_handle(std::shared_ptr<LatticeWalberla> lattice,
54 std::shared_ptr<Boundary_T> boundary) {
55 m_lattice = std::move(lattice);
56 m_boundary = std::move(boundary);
57 }
58
59 bool constantDataExchange() const override { return false; }
60 bool threadsafeReceiving() const override { return true; }
61
62 void communicateLocal(IBlock const *sender, IBlock *receiver,
63 stencil::Direction dir) override {
64 mpi::SendBuffer sBuffer;
65 packDataImpl(sender, dir, sBuffer);
66 mpi::RecvBuffer rBuffer(sBuffer);
67 unpackData(receiver, stencil::inverseDir[dir], rBuffer);
68 }
69
70 void unpackData(IBlock *receiver, stencil::Direction dir,
71 mpi::RecvBuffer &buffer) override {
72
73 auto *flag_field = receiver->getData<GhostLayerField_T>(bdId_);
74 WALBERLA_ASSERT_NOT_NULLPTR(flag_field);
75 WALBERLA_ASSERT_NOT_NULLPTR(m_boundary);
76 WALBERLA_ASSERT_NOT_NULLPTR(m_lattice);
77
78 auto const boundary_flag = flag_field->getFlag(Boundary_flag);
79 auto const gl = numberOfGhostLayersToCommunicate(flag_field);
80 auto const begin = [gl, dir](auto const *flag_field) {
81 return flag_field->beginGhostLayerOnly(gl, dir);
82 };
83
84#ifndef NDEBUG
85 uint_t xSize, ySize, zSize, bSize;
86 buffer >> xSize >> ySize >> zSize >> bSize;
87 uint_t buf_size{0u};
88 for (auto it = begin(flag_field); it != flag_field->end(); ++it) {
89 if (isFlagSet(it, boundary_flag)) {
90 ++buf_size;
91 }
92 }
93 WALBERLA_ASSERT_EQUAL(xSize, flag_field->xSize());
94 WALBERLA_ASSERT_EQUAL(ySize, flag_field->ySize());
95 WALBERLA_ASSERT_EQUAL(zSize, flag_field->zSize());
96 WALBERLA_ASSERT_EQUAL(bSize, buf_size);
97#endif
98
99 auto const offset = m_lattice->get_block_corner(*receiver, true);
100 typename Boundary_T::value_type value;
101 for (auto it = begin(flag_field); it != flag_field->end(); ++it) {
102 if (isFlagSet(it, boundary_flag)) {
103 auto const node = offset + Utils::Vector3i{{it.x(), it.y(), it.z()}};
104 buffer >> value;
105 m_boundary->unpack_node(node, value);
106 }
107 }
108 }
109
110protected:
111 void packDataImpl(IBlock const *sender, stencil::Direction dir,
112 mpi::SendBuffer &buffer) const override {
113
114 auto const *flag_field = sender->getData<GhostLayerField_T>(bdId_);
115 WALBERLA_ASSERT_NOT_NULLPTR(flag_field);
116 WALBERLA_ASSERT_NOT_NULLPTR(m_boundary);
117 WALBERLA_ASSERT_NOT_NULLPTR(m_lattice);
118
119 auto const boundary_flag = flag_field->getFlag(Boundary_flag);
120 auto const gl = numberOfGhostLayersToCommunicate(flag_field);
121 auto const begin = [gl, dir](auto const *flag_field) {
122 return flag_field->beginSliceBeforeGhostLayer(dir, cell_idx_c(gl));
123 };
124
125#ifndef NDEBUG
126 uint_t buf_size{0u};
127 for (auto it = begin(flag_field); it != flag_field->end(); ++it) {
128 if (isFlagSet(it, boundary_flag)) {
129 ++buf_size;
130 }
131 }
132 buffer << flag_field->xSize() << flag_field->ySize() << flag_field->zSize()
133 << buf_size;
134#endif
135
136 auto const offset = m_lattice->get_block_corner(*sender, true);
137 for (auto it = begin(flag_field); it != flag_field->end(); ++it) {
138 if (isFlagSet(it, boundary_flag)) {
139 auto const node = offset + Utils::Vector3i{{it.x(), it.y(), it.z()}};
140 buffer << m_boundary->get_node_value_at_boundary(node);
141 }
142 }
143 }
144
145private:
146 std::shared_ptr<LatticeWalberla> m_lattice;
147 std::shared_ptr<Boundary_T> m_boundary;
148};
149
150} // namespace communication
151} // namespace field
152} // 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