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