ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
PoissonSolver.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
24
25#include <cstddef>
26#include <memory>
27#include <optional>
28#include <utility>
29
30namespace walberla {
31
33private:
34 std::shared_ptr<LatticeWalberla> m_lattice;
35 double m_permittivity;
36
37public:
38 PoissonSolver(std::shared_ptr<LatticeWalberla> lattice, double permittivity)
39 : m_lattice(std::move(lattice)), m_permittivity(permittivity) {}
40 ~PoissonSolver() override = default;
41
42 virtual void reset_charge_field() = 0;
43
44 virtual void add_charge_to_field(std::size_t id, double valency) = 0;
45
46 [[nodiscard]] virtual std::size_t get_potential_field_id() const noexcept = 0;
47
49 m_permittivity = permittivity;
50 }
51
53 return m_permittivity;
54 }
55
57 return *m_lattice;
58 }
59
60 virtual void solve() = 0;
61
62 [[nodiscard]] virtual std::optional<double>
64 bool consider_ghosts = false) = 0;
65
66 [[nodiscard]] virtual std::vector<double>
68 Utils::Vector3i const &upper_corner) const = 0;
69
70 void register_vtk_field_writers(walberla::vtk::VTKOutput &,
72 int) override {}
73
74 virtual void setup_fft(bool use_gpu_aware) = 0;
75
76 [[nodiscard]] virtual bool is_gpu() const noexcept = 0;
78
80 void integrate_vtk_writers() override {}
81 void register_vtk_field_filters(walberla::vtk::VTKOutput &) override {}
82};
83
84} // namespace walberla
Abstract representation of a lattice-based model.
std::unordered_map< std::string, double > units_map
Class that runs and controls the BlockForest in waLBerla.
void register_vtk_field_writers(walberla::vtk::VTKOutput &, LatticeModel::units_map const &, int) override
~PoissonSolver() override=default
virtual std::size_t get_potential_field_id() const noexcept=0
virtual void set_permittivity(double permittivity) noexcept
virtual void add_charge_to_field(std::size_t id, double valency)=0
PoissonSolver(std::shared_ptr< LatticeWalberla > lattice, double permittivity)
virtual std::optional< double > get_node_potential(Utils::Vector3i const &node, bool consider_ghosts=false)=0
virtual void setup_fft(bool use_gpu_aware)=0
virtual void solve()=0
void register_vtk_field_filters(walberla::vtk::VTKOutput &) override
virtual std::vector< double > get_slice_potential(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const =0
virtual bool is_double_precision() const noexcept=0
void integrate_vtk_writers() override
LatticeWalberla const & get_lattice() const noexcept override
virtual void reset_charge_field()=0
virtual double get_permittivity() const noexcept
virtual bool is_gpu() const noexcept=0
cudaStream_t stream[1]
CUDA streams for parallel computing on CPU and GPU.
STL namespace.
\file PackInfoPdfDoublePrecision.cpp \author pystencils