ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
FFT.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2022-2023 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 "PoissonSolver.hpp"
23
24#include <blockforest/communication/UniformBufferedScheme.h>
25#include <domain_decomposition/BlockDataID.h>
26#include <fft/Fft.h>
27#include <field/AddToStorage.h>
28#include <field/GhostLayerField.h>
29#include <field/communication/PackInfo.h>
30#include <stencil/D3Q27.h>
31
32#include <utils/constants.hpp>
33
34#include <cmath>
35#include <cstddef>
36#include <memory>
37#include <utility>
38
39namespace walberla {
40
41template <typename FloatType> class FFT : public PoissonSolver {
42private:
43 template <typename T> FloatType FloatType_c(T t) {
44 return numeric_cast<FloatType>(t);
45 }
46
47 domain_decomposition::BlockDataID m_potential_field_id;
48
49 using PotentialField = GhostLayerField<FloatType, 1>;
50
51 std::shared_ptr<fft::FourierTransform<PotentialField>> m_ft;
52 std::shared_ptr<blockforest::StructuredBlockForest> m_blocks;
53
54 using FullCommunicator = blockforest::communication::UniformBufferedScheme<
55 typename stencil::D3Q27>;
56 std::shared_ptr<FullCommunicator> m_full_communication;
57
58public:
59 FFT(std::shared_ptr<LatticeWalberla> lattice, double permittivity)
60 : PoissonSolver(std::move(lattice), permittivity) {
61 m_blocks = get_lattice().get_blocks();
62
63 Vector3<uint_t> dim(m_blocks->getNumberOfXCells(),
64 m_blocks->getNumberOfYCells(),
65 m_blocks->getNumberOfZCells());
66 auto const greens = [dim](uint_t x, uint_t y, uint_t z) -> real_t {
67 if (x == 0u && y == 0u && z == 0u)
68 return 0.;
69 return -0.5 /
70 (std::cos(2. * Utils::pi() * real_c(x) / real_c(dim[0])) +
71 std::cos(2. * Utils::pi() * real_c(y) / real_c(dim[1])) +
72 std::cos(2. * Utils::pi() * real_c(z) / real_c(dim[2])) - 3.) /
73 real_c(dim[0] * dim[1] * dim[2]);
74 };
75
76 m_potential_field_id = field::addToStorage<PotentialField>(
77 get_lattice().get_blocks(), "potential field", 0.0, field::fzyx,
78 get_lattice().get_ghost_layers());
79
80 m_ft = std::make_shared<fft::FourierTransform<PotentialField>>(
81 m_blocks, m_potential_field_id, greens);
82
83 m_full_communication =
84 std::make_shared<FullCommunicator>(get_lattice().get_blocks());
85 m_full_communication->addPackInfo(
86 std::make_shared<field::communication::PackInfo<PotentialField>>(
87 m_potential_field_id));
88 }
89 ~FFT() override = default;
90
91 void reset_charge_field() override {
92 // the FFT-solver re-uses the potential field for the charge
93 auto const potential_id = walberla::BlockDataID(get_potential_field_id());
94
95 for (auto &block : *get_lattice().get_blocks()) {
96 auto field = block.template getData<PotentialField>(potential_id);
97 WALBERLA_FOR_ALL_CELLS_XYZ(field, field->get(x, y, z) = 0.;)
98 }
99 }
100
101 void add_charge_to_field(std::size_t id, double valency,
102 bool is_double_precision) override {
103 auto const factor = FloatType_c(valency) / FloatType_c(get_permittivity());
104 // the FFT-solver re-uses the potential field for the charge
105 const auto charge_id = walberla::BlockDataID(get_potential_field_id());
106 const auto density_id = walberla::BlockDataID(id);
107 for (auto &block : *get_lattice().get_blocks()) {
108 auto charge_field = block.template getData<PotentialField>(charge_id);
109 if (is_double_precision) {
110 auto density_field =
111 block.template getData<walberla::GhostLayerField<double, 1>>(
112 density_id);
113 WALBERLA_FOR_ALL_CELLS_XYZ(
114 charge_field, charge_field->get(x, y, z) +=
115 factor * FloatType_c(density_field->get(x, y, z));)
116 } else {
117 auto density_field =
118 block.template getData<walberla::GhostLayerField<float, 1>>(
119 density_id);
120 WALBERLA_FOR_ALL_CELLS_XYZ(
121 charge_field, charge_field->get(x, y, z) +=
122 factor * FloatType_c(density_field->get(x, y, z));)
123 }
124 }
125 }
126
127 [[nodiscard]] std::size_t get_potential_field_id() const noexcept override {
128 return static_cast<std::size_t>(m_potential_field_id);
129 }
130
131 void solve() override {
132 (*m_ft)();
133 ghost_communication();
134 }
135
136private:
137 void ghost_communication() { (*m_full_communication)(); }
138};
139
140} // namespace walberla
float u[3]
void solve() override
Definition FFT.hpp:131
void reset_charge_field() override
Definition FFT.hpp:91
FFT(std::shared_ptr< LatticeWalberla > lattice, double permittivity)
Definition FFT.hpp:59
std::size_t get_potential_field_id() const noexcept override
Definition FFT.hpp:127
~FFT() override=default
void add_charge_to_field(std::size_t id, double valency, bool is_double_precision) override
Definition FFT.hpp:101
double get_permittivity() const noexcept
auto const & get_lattice() const noexcept
static double * block(double *p, std::size_t index, std::size_t size)
Definition elc.cpp:174
DEVICE_QUALIFIER constexpr T pi()
Ratio of diameter and circumference of a circle.
Definition constants.hpp:36