ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
p3m.impl.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-2024 The ESPResSo project
3 * Copyright (C) 2002,2003,2004,2005,2006,2007,2008,2009,2010
4 * Max-Planck-Institute for Polymer Research, Theory Group
5 *
6 * This file is part of ESPResSo.
7 *
8 * ESPResSo is free software: you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License as published by
10 * the Free Software Foundation, either version 3 of the License, or
11 * (at your option) any later version.
12 *
13 * ESPResSo is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public License
19 * along with this program. If not, see <http://www.gnu.org/licenses/>.
20 */
21
22#pragma once
23
24#include "config/config.hpp"
25
26#ifdef P3M
27
29
30#include "p3m/common.hpp"
31#include "p3m/data_struct.hpp"
32#include "p3m/interpolation.hpp"
33
34#include "ParticleRange.hpp"
35
36#include <utils/Vector.hpp>
37
38#include <memory>
39#include <optional>
40#include <stdexcept>
41#include <type_traits>
42#include <utility>
43
44template <typename FloatType>
45struct p3m_data_struct_coulomb : public p3m_data_struct<FloatType> {
46 using p3m_data_struct<FloatType>::p3m_data_struct;
47
48 /** number of charged particles (only on head node). */
49 int sum_qpart = 0;
50 /** Sum of square of charges (only on head node). */
51 double sum_q2 = 0.;
52 /** square of sum of charges (only on head node). */
53 double square_sum_q = 0.;
54
56};
57
58#ifdef CUDA
59struct P3MGpuParams;
60#endif
61
62template <typename FloatType, Arch Architecture>
63struct CoulombP3MImpl : public CoulombP3M {
64 ~CoulombP3MImpl() override = default;
65
66 /** @brief Coulomb P3M parameters. */
68
69private:
70 /** @brief Coulomb P3M meshes and FFT algorithm. */
71 std::unique_ptr<p3m_data_struct_coulomb<FloatType>> p3m_impl;
72 int tune_timings;
73 std::pair<std::optional<int>, std::optional<int>> tune_limits;
74 bool tune_verbose;
75 bool check_complex_residuals;
76 bool m_is_tuned;
77
78public:
80 std::unique_ptr<p3m_data_struct_coulomb<FloatType>> &&p3m_handle,
81 double prefactor, int tune_timings, bool tune_verbose,
82 decltype(tune_limits) tune_limits, bool check_complex_residuals)
83 : CoulombP3M(p3m_handle->params), p3m{*p3m_handle},
84 p3m_impl{std::move(p3m_handle)}, tune_timings{tune_timings},
85 tune_limits{std::move(tune_limits)}, tune_verbose{tune_verbose},
86 check_complex_residuals{check_complex_residuals} {
87
88 if (tune_timings <= 0) {
89 throw std::domain_error("Parameter 'timings' must be > 0");
90 }
91 m_is_tuned = not p3m.params.tuning;
92 p3m.params.tuning = false;
94 }
95
96 void init() override {
97 if constexpr (Architecture == Arch::CPU) {
99 }
100#ifdef CUDA
101 if constexpr (Architecture == Arch::GPU) {
103 }
104#endif
105 }
106 void tune() override;
107 void count_charged_particles() override;
108 void count_charged_particles_elc(int n, double sum_q2,
109 double square_sum_q) override {
110 p3m.sum_qpart = n;
111 p3m.sum_q2 = sum_q2;
112 p3m.square_sum_q = square_sum_q;
113 }
114 void adapt_epsilon_elc() override {
115 p3m.params.epsilon = P3M_EPSILON_METALLIC;
116 }
117
118 [[nodiscard]] bool is_tuned() const noexcept override { return m_is_tuned; }
119 [[nodiscard]] bool is_gpu() const noexcept override {
120 return Architecture == Arch::GPU;
121 }
122 [[nodiscard]] bool is_double_precision() const noexcept override {
123 return std::is_same_v<FloatType, double>;
124 }
125
126 void on_activation() override {
127#ifdef CUDA
128 if constexpr (Architecture == Arch::GPU) {
129 request_gpu();
130 }
131#endif
133 tune();
134#ifdef CUDA
135 if constexpr (Architecture == Arch::GPU) {
136 if (is_tuned()) {
138 }
139 }
140#endif
141 }
142
143 double long_range_energy(ParticleRange const &particles) override {
144 return long_range_kernel(false, true, particles);
145 }
146
147 void add_long_range_forces(ParticleRange const &particles) override {
148 if constexpr (Architecture == Arch::CPU) {
149 long_range_kernel(true, false, particles);
150 }
151#ifdef CUDA
152 if constexpr (Architecture == Arch::GPU) {
153 add_long_range_forces_gpu(particles);
154 }
155#endif
156 }
157
158 Utils::Vector9d long_range_pressure(ParticleRange const &particles) override;
159
160 void charge_assign(ParticleRange const &particles) override;
161 void assign_charge(double q, Utils::Vector3d const &real_pos,
162 bool skip_cache) override;
163 void prepare_fft_mesh(bool reset_weights) override {
164 if (reset_weights) {
165 p3m.inter_weights.reset(p3m.params.cao);
166 }
167 for (int i = 0; i < p3m.local_mesh.size; i++) {
168 p3m.mesh.rs_scalar[i] = FloatType(0);
169 }
170 }
171
172protected:
173 /** Compute the k-space part of forces and energies. */
174 double long_range_kernel(bool force_flag, bool energy_flag,
175 ParticleRange const &particles);
176 void calc_influence_function_force() override;
177 void calc_influence_function_energy() override;
178 void scaleby_box_l() override;
179 void init_cpu_kernels();
180#ifdef CUDA
181 void init_gpu_kernels();
182 void add_long_range_forces_gpu(ParticleRange const &particles);
183 std::shared_ptr<P3MGpuParams> m_gpu_data = nullptr;
184 void request_gpu() const;
185#endif
186};
187
188template <typename FloatType, Arch Architecture,
189 template <typename> class FFTBackendImpl,
190 template <typename> class P3MFFTMeshImpl, class... Args>
191std::shared_ptr<CoulombP3M> new_p3m_handle(P3MParameters &&p3m,
192 Args &&...args) {
193 auto obj = std::make_shared<CoulombP3MImpl<FloatType, Architecture>>(
194 std::make_unique<p3m_data_struct_coulomb<FloatType>>(std::move(p3m)),
195 std::forward<Args>(args)...);
196 obj->p3m.template make_mesh_instance<P3MFFTMeshImpl<FloatType>>();
197 obj->p3m.template make_fft_instance<FFTBackendImpl<FloatType>>();
198 return obj;
199}
200
201#endif // P3M
Vector implementation and trait types for boost qvm interoperability.
void set_prefactor(double new_prefactor)
double prefactor
Electrostatics prefactor.
A range of particles.
Cache for interpolation weights.
This file contains the defaults for ESPResSo.
Common functions for dipolar and charge P3M.
Arch
P3M kernel architecture.
auto constexpr P3M_EPSILON_METALLIC
This value indicates metallic boundary conditions.
P3M algorithm for long-range Coulomb interaction.
std::shared_ptr< CoulombP3M > new_p3m_handle(P3MParameters &&p3m, Args &&...args)
Definition p3m.impl.hpp:191
static SteepestDescentParameters params
Currently active steepest descent instance.
void charge_assign(ParticleRange const &particles) override
Definition p3m.cpp:322
double long_range_energy(ParticleRange const &particles) override
Definition p3m.impl.hpp:143
void tune() override
Definition p3m.cpp:804
bool is_gpu() const noexcept override
Definition p3m.impl.hpp:119
bool is_tuned() const noexcept override
Definition p3m.impl.hpp:118
void assign_charge(double q, Utils::Vector3d const &real_pos, bool skip_cache) override
Definition p3m.cpp:334
void calc_influence_function_energy() override
Calculate the influence function optimized for the energy and the self energy correction.
Definition p3m.cpp:139
void add_long_range_forces(ParticleRange const &particles) override
Definition p3m.impl.hpp:147
void count_charged_particles_elc(int n, double sum_q2, double square_sum_q) override
Definition p3m.impl.hpp:108
void scaleby_box_l() override
Definition p3m.cpp:902
void calc_influence_function_force() override
Calculate the optimal influence function of .
Definition p3m.cpp:128
Utils::Vector9d long_range_pressure(ParticleRange const &particles) override
Definition p3m.cpp:395
void init_cpu_kernels()
Definition p3m.cpp:247
bool is_double_precision() const noexcept override
Definition p3m.impl.hpp:122
~CoulombP3MImpl() override=default
void on_activation() override
Definition p3m.impl.hpp:126
CoulombP3MImpl(std::unique_ptr< p3m_data_struct_coulomb< FloatType > > &&p3m_handle, double prefactor, int tune_timings, bool tune_verbose, decltype(tune_limits) tune_limits, bool check_complex_residuals)
Definition p3m.impl.hpp:79
void count_charged_particles() override
Definition p3m.cpp:99
void init() override
Definition p3m.impl.hpp:96
double long_range_kernel(bool force_flag, bool energy_flag, ParticleRange const &particles)
Compute the k-space part of forces and energies.
Definition p3m.cpp:453
std::shared_ptr< P3MGpuParams > m_gpu_data
Definition p3m.impl.hpp:183
void request_gpu() const
Definition p3m.cpp:953
void prepare_fft_mesh(bool reset_weights) override
Definition p3m.impl.hpp:163
void init_gpu_kernels()
Definition p3m.cpp:940
void add_long_range_forces_gpu(ParticleRange const &particles)
Definition p3m.cpp:915
p3m_data_struct_coulomb< FloatType > & p3m
Coulomb P3M parameters.
Definition p3m.impl.hpp:67
void adapt_epsilon_elc() override
Definition p3m.impl.hpp:114
P3M solver.
Definition p3m.hpp:55
void sanity_checks() const
Definition p3m.hpp:79
Structure to hold P3M parameters and some dependent variables.
double sum_q2
Sum of square of charges (only on head node).
Definition p3m.impl.hpp:51
double square_sum_q
square of sum of charges (only on head node).
Definition p3m.impl.hpp:53
int sum_qpart
number of charged particles (only on head node).
Definition p3m.impl.hpp:49
p3m_interpolation_cache inter_weights
Definition p3m.impl.hpp:55
Base class for the electrostatics and magnetostatics P3M algorithms.