ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
dp3m.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 DP3M
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#include <vector>
44
45template <typename FloatType>
46struct p3m_data_struct_dipoles : public p3m_data_struct<FloatType> {
47 using p3m_data_struct<FloatType>::p3m_data_struct;
48
49 /** number of dipolar particles (only on head node). */
50 int sum_dip_part = 0;
51 /** Sum of square of magnetic dipoles (only on head node). */
52 double sum_mu2 = 0.;
53
54 /** position shift for calculation of first assignment mesh point. */
55 double pos_shift = 0.;
56
57 /** cached k-space self-energy correction */
58 double energy_correction = 0.;
59 /** k-space scalar mesh for k-space calculations. */
60 std::vector<FloatType> ks_scalar;
61
63};
64
65template <typename FloatType, Arch Architecture>
66struct DipolarP3MImpl : public DipolarP3M {
67 ~DipolarP3MImpl() override = default;
68
69 /** @brief Dipolar P3M parameters. */
71
72private:
73 /** @brief Coulomb P3M meshes and FFT algorithm. */
74 std::unique_ptr<p3m_data_struct_dipoles<FloatType>> dp3m_impl;
75 int tune_timings;
76 std::pair<std::optional<int>, std::optional<int>> tune_limits;
77 bool tune_verbose;
78 bool m_is_tuned;
79
80public:
82 std::unique_ptr<p3m_data_struct_dipoles<FloatType>> &&dp3m_handle,
83 double prefactor, int tune_timings, bool tune_verbose,
84 decltype(tune_limits) tune_limits)
85 : DipolarP3M(dp3m_handle->params), dp3m{*dp3m_handle},
86 dp3m_impl{std::move(dp3m_handle)}, tune_timings{tune_timings},
87 tune_limits{std::move(tune_limits)}, tune_verbose{tune_verbose} {
88
89 if (tune_timings <= 0) {
90 throw std::domain_error("Parameter 'timings' must be > 0");
91 }
92 if (dp3m.params.mesh != Utils::Vector3i::broadcast(dp3m.params.mesh[0])) {
93 throw std::domain_error("DipolarP3M requires a cubic mesh");
94 }
95 m_is_tuned = not dp3m.params.tuning;
96 dp3m.params.tuning = false;
98 }
99
100 void init() override {
101 if constexpr (Architecture == Arch::CPU) {
103 }
104 }
105 void tune() override;
106 void count_magnetic_particles() override;
107
108 [[nodiscard]] bool is_tuned() const noexcept override { return m_is_tuned; }
109 [[nodiscard]] bool is_gpu() const noexcept override {
110 return Architecture == Arch::GPU;
111 }
112 [[nodiscard]] bool is_double_precision() const noexcept override {
113 return std::is_same_v<FloatType, double>;
114 }
115
116 void on_activation() override {
118 tune();
119 }
120
121 double long_range_energy(ParticleRange const &particles) override {
122 return long_range_kernel(false, true, particles);
123 }
124
125 void add_long_range_forces(ParticleRange const &particles) override {
126 if constexpr (Architecture == Arch::CPU) {
127 long_range_kernel(true, false, particles);
128 }
129 }
130
131 void dipole_assign(ParticleRange const &particles) override;
132
133protected:
134 /** Compute the k-space part of forces and energies. */
135 double long_range_kernel(bool force_flag, bool energy_flag,
136 ParticleRange const &particles);
137 double calc_average_self_energy_k_space() const override;
138 void calc_energy_correction() override;
139 void calc_influence_function_force() override;
140 void calc_influence_function_energy() override;
141 double calc_surface_term(bool force_flag, bool energy_flag,
142 ParticleRange const &particles) override;
143 void init_cpu_kernels();
144 void scaleby_box_l() override;
145};
146
147template <typename FloatType, Arch Architecture,
148 template <typename> class FFTBackendImpl,
149 template <typename> class P3MFFTMeshImpl, class... Args>
150std::shared_ptr<DipolarP3M> new_dp3m_handle(P3MParameters &&p3m,
151 Args &&...args) {
152 auto obj = std::make_shared<DipolarP3MImpl<FloatType, Architecture>>(
153 std::make_unique<p3m_data_struct_dipoles<FloatType>>(std::move(p3m)),
154 std::forward<Args>(args)...);
155 obj->dp3m.template make_mesh_instance<P3MFFTMeshImpl<FloatType>>();
156 obj->dp3m.template make_fft_instance<FFTBackendImpl<FloatType>>();
157 return obj;
158}
159
160#endif // DP3M
Vector implementation and trait types for boost qvm interoperability.
void set_prefactor(double new_prefactor)
double prefactor
Magnetostatics prefactor.
A range of particles.
static DEVICE_QUALIFIER constexpr Vector< T, N > broadcast(typename Base::value_type const &value)
Create a vector that has all entries set to the same value.
Definition Vector.hpp:110
Cache for interpolation weights.
This file contains the defaults for ESPResSo.
P3M algorithm for long-range magnetic dipole-dipole interaction.
std::shared_ptr< DipolarP3M > new_dp3m_handle(P3MParameters &&p3m, Args &&...args)
Common functions for dipolar and charge P3M.
Arch
P3M kernel architecture.
static SteepestDescentParameters params
Currently active steepest descent instance.
void calc_energy_correction() override
Definition dp3m.cpp:904
double long_range_kernel(bool force_flag, bool energy_flag, ParticleRange const &particles)
Compute the k-space part of forces and energies.
Definition dp3m.cpp:249
double long_range_energy(ParticleRange const &particles) override
p3m_data_struct_dipoles< FloatType > & dp3m
Dipolar P3M parameters.
Definition dp3m.impl.hpp:70
double calc_surface_term(bool force_flag, bool energy_flag, ParticleRange const &particles) override
Definition dp3m.cpp:454
void dipole_assign(ParticleRange const &particles) override
Definition dp3m.cpp:179
bool is_double_precision() const noexcept override
void init_cpu_kernels()
Definition dp3m.cpp:130
void calc_influence_function_energy() override
Definition dp3m.cpp:529
bool is_gpu() const noexcept override
void add_long_range_forces(ParticleRange const &particles) override
void count_magnetic_particles() override
Definition dp3m.cpp:85
~DipolarP3MImpl() override=default
void calc_influence_function_force() override
Definition dp3m.cpp:522
bool is_tuned() const noexcept override
DipolarP3MImpl(std::unique_ptr< p3m_data_struct_dipoles< FloatType > > &&dp3m_handle, double prefactor, int tune_timings, bool tune_verbose, decltype(tune_limits) tune_limits)
Definition dp3m.impl.hpp:81
void init() override
void on_activation() override
void tune() override
Definition dp3m.cpp:659
void scaleby_box_l() override
Definition dp3m.cpp:891
double calc_average_self_energy_k_space() const override
Definition dp3m.cpp:116
Dipolar P3M solver.
Definition dp3m.hpp:59
void sanity_checks() const
Definition dp3m.hpp:83
Structure to hold P3M parameters and some dependent variables.
double energy_correction
cached k-space self-energy correction
Definition dp3m.impl.hpp:58
int sum_dip_part
number of dipolar particles (only on head node).
Definition dp3m.impl.hpp:50
double pos_shift
position shift for calculation of first assignment mesh point.
Definition dp3m.impl.hpp:55
p3m_interpolation_cache inter_weights
Definition dp3m.impl.hpp:62
double sum_mu2
Sum of square of magnetic dipoles (only on head node).
Definition dp3m.impl.hpp:52
std::vector< FloatType > ks_scalar
k-space scalar mesh for k-space calculations.
Definition dp3m.impl.hpp:60
Base class for the electrostatics and magnetostatics P3M algorithms.