ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
PoissonSolverNone.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2022-2026 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
25
26#include "greens_function.hpp"
27
30#if defined(__CUDACC__)
33#endif
34
35#include <utils/Vector.hpp>
36
37#include <blockforest/communication/UniformBufferedScheme.h>
38#include <domain_decomposition/BlockDataID.h>
39#include <field/AddToStorage.h>
40#include <field/GhostLayerField.h>
41#include <field/vtk/VTKWriter.h>
42#include <stencil/D3Q27.h>
43#include <waLBerlaDefinitions.h>
44#if defined(__CUDACC__)
45#include <gpu/AddGPUFieldToStorage.h>
46#include <gpu/FieldAccessor.h>
47#include <gpu/FieldIndexing.h>
48#include <gpu/GPUField.h>
49#include <gpu/HostFieldAllocator.h>
50#include <gpu/Kernel.h>
51#include <gpu/communication/UniformGPUScheme.h>
52#endif
53
54#if defined(__clang__)
55#pragma clang diagnostic push
56#pragma clang diagnostic ignored "-Wfloat-conversion"
57#pragma clang diagnostic ignored "-Wimplicit-float-conversion"
58#elif defined(__GNUC__) or defined(__GNUG__)
59#pragma GCC diagnostic push
60#pragma GCC diagnostic ignored "-Wfloat-conversion"
61#endif
62
63#if defined(__clang__)
64#pragma clang diagnostic pop
65#elif defined(__GNUC__) or defined(__GNUG__)
66#pragma GCC diagnostic pop
67#endif
68
69#include <algorithm>
70#include <array>
71#include <complex>
72#include <cstddef>
73#include <functional>
74#include <memory>
75#include <optional>
76#include <ranges>
77#include <string>
78#include <type_traits>
79#include <utility>
80#include <vector>
81
82namespace walberla {
83template <typename FloatType, lbmpy::Arch Architecture>
85private:
86 template <typename T> FloatType FloatType_c(T t) {
88 }
89
90protected:
91 template <typename FT, lbmpy::Arch AT = lbmpy::Arch::CPU> struct FieldTrait {
92 using PotentialField = field::GhostLayerField<FT, 1u>;
93 };
94
95#if defined(__CUDACC__)
96 template <typename FT> struct FieldTrait<FT, lbmpy::Arch::GPU> {
97 using PotentialField = gpu::GPUField<FT>;
98 };
99#endif
100
101public:
103
104private:
105 BlockDataID m_potential_field_id;
106
107public:
108 ~PoissonSolverNone() override = default;
109 explicit PoissonSolverNone(std::shared_ptr<LatticeWalberla> lattice)
110 : PoissonSolver(std::move(lattice), 0.0) {
111 auto blocks = get_lattice().get_blocks();
112#if defined(__CUDACC__)
113 m_potential_field_id = gpu::addGPUFieldToStorage<PotentialField>(
114 blocks, "potential field", 1u, field::fzyx,
115 get_lattice().get_ghost_layers());
116 for (auto &block : *blocks) {
117 auto field = block.template getData<PotentialField>(m_potential_field_id);
118 ek::accessor::Scalar::initialize(field, FloatType{0});
119 }
120#else // __CUDACC__
121 m_potential_field_id = field::addToStorage<PotentialField>(
122 blocks, "potential field", FloatType{0}, field::fzyx,
124#endif // __CUDACC__
125 }
126
127 void setup_fft(bool) override {}
128
129 [[nodiscard]] bool is_gpu() const noexcept override {
131 }
132
134 return std::is_same_v<FloatType, double>;
135 }
136
137 std::size_t get_potential_field_id() const noexcept override {
138 return static_cast<std::size_t>(m_potential_field_id);
139 }
140
141 [[nodiscard]] std::optional<double>
143 bool consider_ghosts = false) override {
145
146 if (not bc or get_potential_field_id() == 0u)
147 return std::nullopt;
148
149 auto const potential_field =
150 bc->block->template getData<PotentialField>(m_potential_field_id);
151 return {double_c(
153 }
154
156 double potential) override {
157 auto bc = get_block_and_cell(get_lattice(), node, false);
158 if (!bc) {
159 return false;
160 }
161 auto potential_field =
162 bc->block->template getData<PotentialField>(m_potential_field_id);
164 bc->cell);
165 return true;
166 }
167
168 [[nodiscard]] std::vector<double>
170 Utils::Vector3i const &upper_corner) const override {
171 std::vector<double> out;
172#ifndef NDEBUG
174#endif
175 auto const &lattice = get_lattice();
176 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
177 out = std::vector<double>(ci->numCells());
178 for (auto &block : *lattice.get_blocks()) {
179 auto const block_offset = lattice.get_block_corner(block, true);
180 if (auto const bci = get_block_interval(
182 auto const potential_field =
183 block.template getData<PotentialField>(m_potential_field_id);
185 assert(values.size() == bci->numCells());
186#ifndef NDEBUG
187 values_size += bci->numCells();
188#endif
189 auto kernel = [&values, &out](unsigned const block_index,
190 unsigned const local_index,
191 Utils::Vector3i const &) {
193 };
194
196 }
197 }
198 assert(values_size == ci->numCells());
199 }
200 return out;
201 }
202
205 std::vector<double> const &potential) override {
206 auto const &lattice = get_lattice();
207 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
208 assert(potential.size() == ci->numCells());
209 for (auto &block : *lattice.get_blocks()) {
210 auto const block_offset = lattice.get_block_corner(block, true);
211 if (auto const bci = get_block_interval(
213 auto potential_field =
214 block.template getData<PotentialField>(m_potential_field_id);
215 std::vector<FloatType> values(bci->numCells());
216
217 auto kernel = [&values, &potential](unsigned const block_index,
218 unsigned const local_index,
219 Utils::Vector3i const &) {
222 };
223
226 }
227 }
228 }
229 }
230 void ghost_communication() override {}
231
232 void solve() override { integrate_vtk_writers(); }
233
234 void add_charge_to_field(std::size_t id, double valency) override {}
235 void reset_charge_field() override {}
236
237protected:
238 void integrate_vtk_writers() override {
239 for (auto const &vtk_handle : m_vtk_auto | std::views::values) {
240 if (vtk_handle->enabled) {
241 vtk::writeFiles(vtk_handle->ptr)();
242 vtk_handle->execution_count++;
243 }
244 }
245 }
246
247protected:
248 template <typename VecType, uint_t F_SIZE_ARG, typename OutputType>
249 class VTKWriter : public vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG> {
250 public:
251 VTKWriter(ConstBlockDataID const &block_id, std::string const &id,
252 FloatType unit_conversion)
253 : vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG>(id),
255
256 protected:
257 void configure() override { WALBERLA_ASSERT_NOT_NULLPTR(this->block_); }
258
259 std::size_t get_first_index(cell_idx_t const x, cell_idx_t const y,
260 cell_idx_t const z) {
261 return (static_cast<std::size_t>(x) * m_dims[2] * m_dims[1] +
262 static_cast<std::size_t>(y) * m_dims[2] +
263 static_cast<std::size_t>(z)) *
265 }
266
267 FloatType m_conversion;
270
271 public:
273
274 void set_dims(Vector3<uint_t> dims) { m_dims = dims; }
275 };
276
277 template <typename OutputType = float>
279 : public VTKWriter<std::vector<FloatType>, 1u, OutputType> {
280 public:
282 using Base::Base;
283 using Base::evaluate;
284
285 protected:
287 cell_idx_t const z, cell_idx_t const) override {
288 WALBERLA_ASSERT(!this->m_content.empty());
289 auto const potential = this->m_content[this->get_first_index(x, y, z)];
291 }
292 };
293
294public:
295 void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj,
297 int flag_observables) override {
298 if (flag_observables & static_cast<int>(EKPoissonOutputVTK::potential)) {
299 auto const unit_conversion = FloatType_c(units.at("potential"));
300 auto const blocks = get_lattice().get_blocks();
303 m_potential_field_id, "potential", unit_conversion);
304 auto before_function = [this, blocks, potential_writer]() {
305 for (auto &block : *blocks) {
306 auto *potential_field =
307 block.template getData<PotentialField>(m_potential_field_id);
308 auto const bci = potential_field->xyzSize();
309 potential_writer->set_content(
312 uint_c(bci.xSize()), uint_c(bci.ySize()), uint_c(bci.zSize())));
313 }
314 };
315 vtk_obj.addBeforeFunction(std::move(before_function));
316 vtk_obj.addCellDataWriter(potential_writer);
317 }
318 }
319};
320
321} // namespace walberla
Vector implementation and trait types for boost qvm interoperability.
std::map< std::string, std::shared_ptr< VTKHandle > > m_vtk_auto
VTK writers that are executed automatically.
std::unordered_map< std::string, double > units_map
OutputType evaluate(cell_idx_t const x, cell_idx_t const y, cell_idx_t const z, cell_idx_t const) override
VTKWriter(ConstBlockDataID const &block_id, std::string const &id, FloatType unit_conversion)
std::size_t get_first_index(cell_idx_t const x, cell_idx_t const y, cell_idx_t const z)
bool is_double_precision() const noexcept override
~PoissonSolverNone() override=default
bool is_gpu() const noexcept override
std::size_t get_potential_field_id() const noexcept override
void set_slice_potential(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< double > const &potential) override
std::vector< double > get_slice_potential(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
PoissonSolverNone(std::shared_ptr< LatticeWalberla > lattice)
void add_charge_to_field(std::size_t id, double valency) override
FieldTrait< FloatType, Architecture >::PotentialField PotentialField
bool set_node_potential(Utils::Vector3i const &node, double potential) override
std::optional< double > get_node_potential(Utils::Vector3i const &node, bool consider_ghosts=false) override
void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj, LatticeModel::units_map const &units, int flag_observables) override
LatticeWalberla const & get_lattice() const noexcept override
cudaStream_t stream[1]
CUDA streams for parallel computing on CPU and GPU.
static double * block(double *p, std::size_t index, std::size_t size)
Definition elc.cpp:175
STL namespace.
void initialize(GhostLayerField< double, 1u > *scalar_field, double const &value)
void set(GhostLayerField< double, 1u > *scalar_field, double const &value, Cell const &cell)
auto get(GhostLayerField< double, 1u > const *scalar_field, Cell const &cell)
\file PackInfoPdfDoublePrecision.cpp \author pystencils
void copy_block_buffer(CellInterval const &bci, CellInterval const &ci, Utils::Vector3i const &block_offset, Utils::Vector3i const &lower_corner, auto &&kernel)
Synchronize data between a sliced block and a container.
std::optional< BlockAndCell > get_block_and_cell(::LatticeWalberla const &lattice, signed_integral_vector auto const &node, bool consider_ghost_layers)
std::optional< walberla::cell::CellInterval > get_block_interval(::LatticeWalberla const &lattice, Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, Utils::Vector3i const &block_offset, IBlock const &block)
std::optional< walberla::cell::CellInterval > get_interval(::LatticeWalberla const &lattice, Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner)
field::GhostLayerField< FT, 1u > PotentialField