ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
PoissonSolverFFT.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2025-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/communication/PackInfo.h>
42#include <field/vtk/VTKWriter.h>
43#include <stencil/D3Q27.h>
44#include <waLBerlaDefinitions.h>
45#if defined(__CUDACC__)
46#include <gpu/AddGPUFieldToStorage.h>
47#include <gpu/FieldAccessor.h>
48#include <gpu/FieldIndexing.h>
49#include <gpu/GPUField.h>
50#include <gpu/Kernel.h>
51#include <gpu/communication/MemcpyPackInfo.h>
52#include <gpu/communication/UniformGPUScheme.h>
53#endif
54
55#if defined(__clang__)
56#pragma clang diagnostic push
57#pragma clang diagnostic ignored "-Wfloat-conversion"
58#pragma clang diagnostic ignored "-Wimplicit-float-conversion"
59#elif defined(__GNUC__) or defined(__GNUG__)
60#pragma GCC diagnostic push
61#pragma GCC diagnostic ignored "-Wfloat-conversion"
62#endif
63
64#include <heffte.h>
65#include <heffte_backends.h>
66#include <heffte_geometry.h>
67
68#if defined(__clang__)
69#pragma clang diagnostic pop
70#elif defined(__GNUC__) or defined(__GNUG__)
71#pragma GCC diagnostic pop
72#endif
73
74#include <algorithm>
75#include <array>
76#include <complex>
77#include <cstddef>
78#include <functional>
79#include <memory>
80#include <optional>
81#include <string>
82#include <type_traits>
83#include <utility>
84#include <vector>
85
86namespace walberla {
87
88template <typename T, std::size_t N>
89auto to_array(Utils::Vector<T, N> const &vec) {
90 std::array<T, N> res{};
91 std::ranges::copy(vec, res.begin());
92 return res;
93}
94
95inline int pos_to_linear_index(int x, int y, int z, auto const &dim) {
96 return (z * dim[1] + y) * dim[0] + x;
97}
98
99template <typename FloatType, lbmpy::Arch Architecture>
101private:
102 template <typename T> FloatType FloatType_c(T t) {
103 return numeric_cast<FloatType>(t);
104 }
105
106protected:
107 template <typename FT, lbmpy::Arch AT = lbmpy::Arch::CPU> struct FieldTrait {
108 using ComplexType = std::complex<FloatType>;
109 using PotentialField = field::GhostLayerField<FT, 1u>;
110 template <class Field>
111 using PackInfo = field::communication::PackInfo<Field>;
112 template <class Stencil>
114 blockforest::communication::UniformBufferedScheme<Stencil>;
115 };
116
117#if defined(__CUDACC__)
118 template <typename FT> struct FieldTrait<FT, lbmpy::Arch::GPU> {
119 using ComplexType = std::conditional_t<std::is_same_v<FloatType, float>,
120 cufftComplex, cufftDoubleComplex>;
121 using PotentialField = gpu::GPUField<FT>;
122 template <class Field>
123 using PackInfo = gpu::communication::MemcpyPackInfo<Field>;
124 template <class Stencil>
125 using FullCommunicator = gpu::communication::UniformGPUScheme<Stencil>;
126 };
127#endif
128
129public:
133 FieldTrait<FloatType,
134 Architecture>::template FullCommunicator<stencil::D3Q27>;
135 template <class Field>
136 using PackInfo =
138
139protected:
140 template <typename ComplexType> struct heffte_container {
141#if defined(__CUDACC__)
142 using backend = heffte::backend::cufft;
143#else // __CUDACC__
144 using backend = heffte::backend::fftw;
145#endif // __CUDACC__
146 std::unique_ptr<heffte::box3d<>> box_in;
147 std::unique_ptr<heffte::box3d<>> box_out;
148 std::unique_ptr<heffte::fft3d<backend>> fft;
149 std::unique_ptr<heffte::fft3d<backend>::buffer_container<ComplexType>>
151
152 heffte_container(heffte::plan_options const &options,
153 auto const &grid_range) {
154 auto const offset_vec = Utils::Vector3i({1, 1, 1});
155 auto const order = Utils::Vector3i({0, 1, 2});
156 box_in = std::make_unique<heffte::box3d<>>(
157 to_array(Utils::Vector3i(grid_range.first)),
158 to_array(grid_range.second - offset_vec), to_array(order));
159 box_out = std::make_unique<heffte::box3d<>>(
160 to_array(Utils::Vector3i(grid_range.first)),
161 to_array(grid_range.second - offset_vec), to_array(order));
162 fft = std::make_unique<heffte::fft3d<backend>>(*box_in, *box_out,
163 MPI_COMM_WORLD, options);
164 buffer = std::make_unique<
165 heffte::fft3d<backend>::buffer_container<ComplexType>>(
166 fft->size_workspace());
167 }
168 };
169
170private:
171 BlockDataID m_potential_field_with_ghosts_id;
172#if defined(__CUDACC__)
173 BlockDataID m_potential_field_id;
174 BlockDataID m_greens_function_field_id;
175 BlockDataID m_potential_fourier_id;
176#endif // __CUDACC__
177
178 std::unique_ptr<heffte_container<ComplexType>> heffte;
179 std::shared_ptr<FullCommunicator> m_full_communication;
180
181#if defined(__CUDACC__)
182 using GreenFunctionField = gpu::GPUField<FloatType>;
183 using PotentialFourier = gpu::GPUField<ComplexType>;
184
185 walberla::gpu::Kernel<void (*)(walberla::gpu::FieldAccessor<ComplexType>,
186 walberla::gpu::FieldAccessor<FloatType>)>
187 kernel_greens;
188 walberla::gpu::Kernel<void (*)(walberla::gpu::FieldAccessor<FloatType>,
189 walberla::gpu::FieldAccessor<FloatType>)>
190 kernel_move_fields;
191#else // __CUDACC__
192 std::vector<FloatType> m_greens;
193 std::vector<FloatType> m_potential;
194 std::vector<ComplexType> m_potential_fourier;
195#endif // __CUDACC__
196
197public:
198 ~PoissonSolverFFT() override = default;
199 PoissonSolverFFT(std::shared_ptr<LatticeWalberla> lattice,
200 double permittivity)
201 : PoissonSolver(std::move(lattice), permittivity)
202#if defined(__CUDACC__)
203 ,
204 kernel_greens(gpu::make_kernel(
205 multiply_by_greens_function<FloatType, ComplexType>)),
206 kernel_move_fields(gpu::make_kernel(move_field<FloatType>))
207#endif
208 {
209 auto blocks = get_lattice().get_blocks();
210#if defined(__CUDACC__)
211 m_potential_field_id = gpu::addGPUFieldToStorage<PotentialField>(
212 blocks, "potential field", 1u, field::fzyx, 0u, false);
213 m_potential_field_with_ghosts_id =
214 gpu::addGPUFieldToStorage<PotentialField>(
215 blocks, "potential field with ghosts", 1u, field::fzyx,
216 get_lattice().get_ghost_layers());
217 m_greens_function_field_id = gpu::addGPUFieldToStorage<GreenFunctionField>(
218 blocks, "greens function", 1u, field::fzyx, 0u, false);
219 m_potential_fourier_id = gpu::addGPUFieldToStorage<PotentialFourier>(
220 blocks, "fourier field", 1u, field::fzyx, 0u, false);
221
222 auto const grid_range = get_lattice().get_local_grid_range();
223 auto const global_dim = get_lattice().get_grid_dimensions();
224 for (auto &block : *blocks) {
225 auto green_field = block.template getData<GreenFunctionField>(
226 m_greens_function_field_id);
227 auto kernel = gpu::make_kernel(create_greens_function<FloatType>);
228 kernel.addFieldIndexingParam(
229 gpu::FieldIndexing<FloatType>::xyz(*green_field));
230 kernel.addParam(grid_range.first[0]);
231 kernel.addParam(grid_range.first[1]);
232 kernel.addParam(grid_range.first[2]);
233 kernel.addParam(grid_range.second[0]);
234 kernel.addParam(grid_range.second[1]);
235 kernel.addParam(grid_range.second[2]);
236 kernel.addParam(global_dim[0]);
237 kernel.addParam(global_dim[1]);
238 kernel.addParam(global_dim[2]);
239 kernel();
240
241 auto potential =
242 block.template getData<PotentialField>(m_potential_field_id);
243 auto potential_ghosts = block.template getData<PotentialField>(
244 m_potential_field_with_ghosts_id);
245 auto green = block.template getData<GreenFunctionField>(
246 m_greens_function_field_id);
247 auto fourier =
248 block.template getData<PotentialFourier>(m_potential_fourier_id);
249
250 kernel_greens =
251 gpu::make_kernel(multiply_by_greens_function<FloatType, ComplexType>);
252 kernel_greens.addFieldIndexingParam(
253 gpu::FieldIndexing<ComplexType>::allInner(*fourier));
254 kernel_greens.addFieldIndexingParam(
255 gpu::FieldIndexing<FloatType>::allInner(*green));
256
257 kernel_move_fields = gpu::make_kernel(move_field<FloatType>);
258 kernel_move_fields.addFieldIndexingParam(
259 gpu::FieldIndexing<FloatType>::xyz(*potential_ghosts));
260 kernel_move_fields.addFieldIndexingParam(
261 gpu::FieldIndexing<FloatType>::xyz(*potential));
262 }
263
264#else // __CUDACC__
265 m_potential_field_with_ghosts_id = field::addToStorage<PotentialField>(
266 blocks, "potential field with ghosts", 0., field::fzyx,
267 get_lattice().get_ghost_layers());
268#endif // __CUDACC__
269
270 m_full_communication = std::make_shared<FullCommunicator>(blocks);
271 m_full_communication->addPackInfo(
272 std::make_shared<PackInfo<PotentialField>>(
273 m_potential_field_with_ghosts_id));
274 }
275
276 [[nodiscard]] bool is_gpu() const noexcept override {
277 return Architecture == lbmpy::Arch::GPU;
278 }
279
280 [[nodiscard]] bool is_double_precision() const noexcept override {
281 return std::is_same_v<FloatType, double>;
282 }
283
284 std::size_t get_potential_field_id() const noexcept override {
285 return static_cast<std::size_t>(m_potential_field_with_ghosts_id);
286 }
287
288 void setup_fft([[maybe_unused]] bool use_gpu_aware) override {
289 auto const grid_range = get_lattice().get_local_grid_range();
290 heffte::plan_options options = heffte::default_options<
292#if defined(__CUDACC__)
293 if constexpr (Architecture == lbmpy::Arch::GPU) {
294 options.use_reorder = false;
295 options.algorithm = heffte::reshape_algorithm::p2p_plined;
296 options.use_pencils = true;
297 options.use_gpu_aware = use_gpu_aware;
298 }
299#endif
300 heffte =
301 std::make_unique<heffte_container<ComplexType>>(options, grid_range);
302#if not defined(__CUDACC__)
303 if constexpr (Architecture == lbmpy::Arch::CPU) {
304 m_potential = std::vector<FloatType>(heffte->fft->size_inbox());
305 m_greens = std::vector<FloatType>(heffte->fft->size_outbox());
306 m_potential_fourier =
307 std::vector<ComplexType>(heffte->fft->size_outbox());
308 auto const dim = grid_range.second - grid_range.first;
309 auto const global_dim = get_lattice().get_grid_dimensions();
310 for (int x = 0; x < dim[0]; x++) {
311 for (int y = 0; y < dim[1]; y++) {
312 for (int z = 0; z < dim[2]; z++) {
313 m_greens[pos_to_linear_index(x, y, z, dim)] =
314 greens_function<FloatType>(x + grid_range.first[0],
315 y + grid_range.first[1],
316 z + grid_range.first[2], global_dim);
317 }
318 }
319 }
320 }
321#endif
323#if not defined(__CUDACC__)
324 if constexpr (Architecture == lbmpy::Arch::CPU) {
325 // make FFT plan
326 heffte->fft->forward(m_potential.data(), m_potential_fourier.data(),
327 heffte->buffer->data());
328 }
329#endif
330#if defined(__CUDACC__)
331 if constexpr (Architecture == lbmpy::Arch::GPU) {
332 for (auto &block : *get_lattice().get_blocks()) {
333 auto potential =
334 block.template getData<PotentialField>(m_potential_field_id);
335 auto fourier =
336 block.template getData<PotentialFourier>(m_potential_fourier_id);
337 FloatType *_data_potential = potential->dataAt(0, 0, 0, 0);
338 ComplexType *_data_fourier = fourier->dataAt(0, 0, 0, 0);
339 // make FFT plan
340 heffte->fft->forward(_data_potential, _data_fourier,
341 heffte->buffer->data());
342 }
343 }
344#endif
345 }
346
347 [[nodiscard]] std::optional<double>
349 bool consider_ghosts = false) override {
350 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
351
352 if (not bc or get_potential_field_id() == 0u)
353 return std::nullopt;
354
355 auto const potential_field = bc->block->template getData<PotentialField>(
356 m_potential_field_with_ghosts_id);
357 return {double_c(
358 walberla::ek::accessor::Scalar::get(potential_field, bc->cell))};
359 }
360
361 [[nodiscard]] std::vector<double>
363 Utils::Vector3i const &upper_corner) const override {
364 std::vector<double> out;
365#ifndef NDEBUG
366 uint_t values_size{0u};
367#endif
368 auto const &lattice = get_lattice();
369 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
370 out = std::vector<double>(ci->numCells());
371 for (auto &block : *lattice.get_blocks()) {
372 auto const block_offset = lattice.get_block_corner(block, true);
373 if (auto const bci = get_block_interval(
374 lattice, lower_corner, upper_corner, block_offset, block)) {
375 auto const potential_field = block.template getData<PotentialField>(
376 m_potential_field_with_ghosts_id);
377 auto const values = ek::accessor::Scalar::get(potential_field, *bci);
378 assert(values.size() == bci->numCells());
379#ifndef NDEBUG
380 values_size += bci->numCells();
381#endif
382 auto kernel = [&values, &out](unsigned const block_index,
383 unsigned const local_index,
384 Utils::Vector3i const &) {
385 out[local_index] = double_c(values[block_index]);
386 };
387
388 copy_block_buffer(*bci, *ci, block_offset, lower_corner, kernel);
389 }
390 }
391 assert(values_size == ci->numCells());
392 }
393 return out;
394 }
395
396 void solve() override {
397#if not defined(__CUDACC__)
398 if constexpr (Architecture == lbmpy::Arch::CPU) {
399 auto grid_range = get_lattice().get_local_grid_range();
400 auto dim = grid_range.second - grid_range.first;
401 heffte->fft->forward(m_potential.data(), m_potential_fourier.data(),
402 heffte->buffer->data());
403 std::ranges::transform(m_potential_fourier, m_greens,
404 m_potential_fourier.begin(), std::multiplies<>{});
405 heffte->fft->backward(m_potential_fourier.data(), m_potential.data(),
406 heffte->buffer->data());
407
408 for (auto &block : *get_lattice().get_blocks()) {
409 auto potential_with_ghosts = block.template getData<PotentialField>(
410 m_potential_field_with_ghosts_id);
411 for (int x = 0; x < dim[0]; x++) {
412 for (int y = 0; y < dim[1]; y++) {
413 for (int z = 0; z < dim[2]; z++) {
414 potential_with_ghosts->get(x, y, z) =
415 m_potential[pos_to_linear_index(x, y, z, dim)];
416 }
417 }
418 }
419 }
420 }
421#endif
422#if defined(__CUDACC__)
423 if constexpr (Architecture == lbmpy::Arch::GPU) {
424 for (auto &block : *get_lattice().get_blocks()) {
425 auto potential =
426 block.template getData<PotentialField>(m_potential_field_id);
427 auto fourier =
428 block.template getData<PotentialFourier>(m_potential_fourier_id);
429 FloatType *_data_potential = potential->dataAt(0, 0, 0, 0);
430 ComplexType *_data_fourier = fourier->dataAt(0, 0, 0, 0);
431 heffte->fft->forward(_data_potential, _data_fourier,
432 heffte->buffer->data());
433 kernel_greens();
434 heffte->fft->backward(_data_fourier, _data_potential,
435 heffte->buffer->data());
436 kernel_move_fields();
437 }
438 }
439#endif
442 }
443
444 void add_charge_to_field(std::size_t id, double valency) override {
445 auto const factor = FloatType_c(valency / get_permittivity());
446 auto const density_id = BlockDataID(id);
447#if not defined(__CUDACC__)
448 if constexpr (Architecture == lbmpy::Arch::CPU) {
449 auto grid_range = get_lattice().get_local_grid_range();
450 auto dim = grid_range.second - grid_range.first;
451 for (auto &block : *get_lattice().get_blocks()) {
452 auto density_field = block.template getData<PotentialField>(density_id);
453 for (int x = 0; x < dim[0]; x++) {
454 for (int y = 0; y < dim[1]; y++) {
455 for (int z = 0; z < dim[2]; z++) {
456 m_potential[pos_to_linear_index(x, y, z, dim)] +=
457 factor * density_field->get(x, y, z);
458 }
459 }
460 }
461 }
462 }
463#endif
464#if defined(__CUDACC__)
465 if constexpr (Architecture == lbmpy::Arch::GPU) {
466 for (auto &block : *get_lattice().get_blocks()) {
467 auto field =
468 block.template getData<PotentialField>(m_potential_field_id);
469 auto density_field =
470 block.template getData<gpu::GPUField<FloatType>>(density_id);
471 add_fields(field, density_field, FloatType_c(factor));
472 }
473 }
474#endif
475 }
476
477 void reset_charge_field() override {
478#if not defined(__CUDACC__)
479 if constexpr (Architecture == lbmpy::Arch::CPU) {
480 auto const grid_range = get_lattice().get_local_grid_range();
481 auto const dim = grid_range.second - grid_range.first;
482 for (int x = 0; x < dim[0]; x++) {
483 for (int i = 0; i < m_potential_fourier.size(); i++) {
484 m_potential_fourier[i] *= m_greens[i];
485 }
486 for (int y = 0; y < dim[1]; y++) {
487 for (int z = 0; z < dim[2]; z++) {
488 m_potential[pos_to_linear_index(x, y, z, dim)] = FloatType(0.0);
489 }
490 }
491 }
492 }
493#endif
494#if defined(__CUDACC__)
495 if constexpr (Architecture == lbmpy::Arch::GPU) {
496 // the FFT-solver re-uses the potential field for the charge
497 for (auto &block : *get_lattice().get_blocks()) {
498 auto field =
499 block.template getData<PotentialField>(m_potential_field_id);
500 ek::accessor::Scalar::initialize(field, FloatType_c(0.));
501 }
502 }
503#endif
504 }
505
506protected:
507 void ghost_communication() { (*m_full_communication)(); }
508
509 void integrate_vtk_writers() override {
510 for (auto const &it : m_vtk_auto) {
511 auto &vtk_handle = it.second;
512 if (vtk_handle->enabled) {
513 vtk::writeFiles(vtk_handle->ptr)();
514 vtk_handle->execution_count++;
515 }
516 }
517 }
518
519protected:
520 template <typename VecType, uint_t F_SIZE_ARG, typename OutputType>
521 class VTKWriter : public vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG> {
522 public:
523 VTKWriter(ConstBlockDataID const &block_id, std::string const &id,
524 FloatType unit_conversion)
525 : vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG>(id),
526 m_conversion(unit_conversion), m_content{} {}
527
528 protected:
529 void configure() override { WALBERLA_ASSERT_NOT_NULLPTR(this->block_); }
530
531 std::size_t get_first_index(cell_idx_t const x, cell_idx_t const y,
532 cell_idx_t const z) {
533 return (static_cast<std::size_t>(x) * m_dims[2] * m_dims[1] +
534 static_cast<std::size_t>(y) * m_dims[2] +
535 static_cast<std::size_t>(z)) *
536 F_SIZE_ARG;
537 }
538
539 FloatType m_conversion;
540 VecType m_content;
541 Vector3<uint_t> m_dims;
542
543 public:
544 void set_content(VecType content) { m_content = content; }
545
546 void set_dims(Vector3<uint_t> dims) { m_dims = dims; }
547 };
548
549 template <typename OutputType = float>
551 : public VTKWriter<std::vector<FloatType>, 1u, OutputType> {
552 public:
553 using Base = VTKWriter<std::vector<FloatType>, 1u, OutputType>;
554 using Base::Base;
555 using Base::evaluate;
556
557 protected:
558 OutputType evaluate(cell_idx_t const x, cell_idx_t const y,
559 cell_idx_t const z, cell_idx_t const) override {
560 WALBERLA_ASSERT(!this->m_content.empty());
561 auto const potential = this->m_content[this->get_first_index(x, y, z)];
562 return numeric_cast<OutputType>(this->m_conversion * potential);
563 }
564 };
565
566public:
567 void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj,
568 LatticeModel::units_map const &units,
569 int flag_observables) override {
570 if (flag_observables & static_cast<int>(EKPoissonOutputVTK::potential)) {
571 auto const unit_conversion = FloatType_c(units.at("potential"));
572 auto const blocks = get_lattice().get_blocks();
573 WALBERLA_ASSERT_NOT_NULLPTR(blocks);
574 auto potential_writer = make_shared<PotentialVTKWriter<float>>(
575 m_potential_field_with_ghosts_id, "potential", unit_conversion);
576 auto before_function = [this, blocks, potential_writer]() {
577 for (auto &block : *blocks) {
578 auto *potential_field = block.template getData<PotentialField>(
579 m_potential_field_with_ghosts_id);
580 auto const bci = potential_field->xyzSize();
581 potential_writer->set_content(
582 walberla::ek::accessor::Scalar::get(potential_field, bci));
583 potential_writer->set_dims(Vector3<uint_t>(
584 uint_c(bci.xSize()), uint_c(bci.ySize()), uint_c(bci.zSize())));
585 }
586 };
587 vtk_obj.addBeforeFunction(std::move(before_function));
588 vtk_obj.addCellDataWriter(potential_writer);
589 }
590 }
591
592private:
593#if defined(__CUDACC__)
594 void add_fields(PotentialField *field_out,
595 gpu::GPUField<FloatType> *field_add, FloatType factor) {
596 auto kernel = gpu::make_kernel(add_fields_with_factor<FloatType>);
597 kernel.addFieldIndexingParam(
598 gpu::FieldIndexing<FloatType>::xyz(*field_out));
599 kernel.addFieldIndexingParam(
600 gpu::FieldIndexing<FloatType>::xyz(*field_add));
601 kernel.addParam(factor);
602 kernel();
603 }
604#endif
605};
606
607} // 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
std::pair< Utils::Vector3i, Utils::Vector3i > get_local_grid_range(bool with_halo=false) const
DEVICE_QUALIFIER constexpr iterator begin() noexcept
Definition Array.hpp:140
OutputType evaluate(cell_idx_t const x, cell_idx_t const y, cell_idx_t const z, cell_idx_t const) override
std::size_t get_first_index(cell_idx_t const x, cell_idx_t const y, cell_idx_t const z)
void set_dims(Vector3< uint_t > dims)
VTKWriter(ConstBlockDataID const &block_id, std::string const &id, FloatType unit_conversion)
FieldTrait< FloatType, Architecture >::template FullCommunicator< stencil::D3Q27 > FullCommunicator
bool is_gpu() const noexcept override
bool is_double_precision() const noexcept override
std::vector< double > get_slice_potential(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
PoissonSolverFFT(std::shared_ptr< LatticeWalberla > lattice, double permittivity)
void setup_fft(bool use_gpu_aware) override
FieldTrait< FloatType, Architecture >::ComplexType ComplexType
std::optional< double > get_node_potential(Utils::Vector3i const &node, bool consider_ghosts=false) override
void add_charge_to_field(std::size_t id, double valency) override
FieldTrait< FloatType, Architecture >::PotentialField PotentialField
void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj, LatticeModel::units_map const &units, int flag_observables) override
FieldTrait< FloatType, Architecture >::template PackInfo< Field > PackInfo
~PoissonSolverFFT() override=default
std::size_t get_potential_field_id() const noexcept override
LatticeWalberla const & get_lattice() const noexcept override
virtual double get_permittivity() const noexcept
static double * block(double *p, std::size_t index, std::size_t size)
Definition elc.cpp:175
VectorXi< 3 > Vector3i
Definition Vector.hpp:193
Definition fft.cpp:76
STL namespace.
void initialize(GhostLayerField< double, 1u > *scalar_field, double const &value)
auto get(GhostLayerField< double, 1u > const *scalar_field, Cell const &cell)
\file PackInfoPdfDoublePrecision.cpp \author pystencils
int pos_to_linear_index(int x, int y, int z, auto const &dim)
auto to_array(Utils::Vector< T, N > const &vec)
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)
blockforest::communication::UniformBufferedScheme< Stencil > FullCommunicator
field::GhostLayerField< FT, 1u > PotentialField
field::communication::PackInfo< Field > PackInfo
std::unique_ptr< heffte::fft3d< backend > > fft
std::unique_ptr< heffte::box3d<> > box_in
std::unique_ptr< heffte::box3d<> > box_out
std::unique_ptr< heffte::fft3d< backend >::buffer_container< ComplexType > > buffer
heffte_container(heffte::plan_options const &options, auto const &grid_range)