Loading [MathJax]/extensions/TeX/AMSmath.js
ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages Concepts
dp3m.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/** @file
23 * P3M algorithm for long-range magnetic dipole-dipole interaction.
24 *
25 * We use here a P3M (Particle-Particle Particle-Mesh) method based
26 * on the dipolar Ewald summation. Details of the used method can be
27 * found in @cite hockney88a and @cite deserno98a @cite deserno98b.
28 *
29 * Further reading: @cite cerda08d
30 */
31
32#pragma once
33
34#include "config/config.hpp"
35
36#ifdef DP3M
37
39
40#include "p3m/common.hpp"
41#include "p3m/data_struct.hpp"
42
43#include "Particle.hpp"
44#include "ParticleRange.hpp"
45
46#include <utils/Vector.hpp>
48#include <utils/math/sqr.hpp>
49
50#include <cmath>
51#include <numbers>
52
53/** @brief Dipolar P3M solver. */
54struct DipolarP3M : public Dipoles::Actor<DipolarP3M> {
56
57public:
59
60 virtual ~DipolarP3M() = default;
61
62 [[nodiscard]] virtual bool is_tuned() const noexcept = 0;
63 [[nodiscard]] virtual bool is_gpu() const noexcept = 0;
64 [[nodiscard]] virtual bool is_double_precision() const noexcept = 0;
65
66 virtual void on_activation() = 0;
67 /** @brief Recalculate all box-length-dependent parameters. */
75 /** @brief Recalculate all derived parameters. */
76 virtual void init() = 0;
77
84
85 /**
86 * @brief Count the number of magnetic particles and calculate
87 * the sum of the squared dipole moments.
88 */
89 virtual void count_magnetic_particles() = 0;
90
91 /** Assign the physical dipoles using the tabulated assignment function. */
92 virtual void dipole_assign(ParticleRange const &particles) = 0;
93
94 /**
95 * @brief Tune dipolar P3M parameters to desired accuracy.
96 *
97 * The parameters
98 * @ref P3MParameters::mesh "mesh",
99 * @ref P3MParameters::cao "cao",
100 * @ref P3MParameters::r_cut_iL "r_cut_iL" and
101 * @ref P3MParameters::alpha_L "alpha_L" are tuned to obtain the target
102 * @ref P3MParameters::accuracy "accuracy" in optimal time.
103 * These parameters are stored in the @ref dp3m_params object.
104 *
105 * The function utilizes the analytic expression of the error estimate
106 * for the dipolar P3M method in the paper of @cite cerda08d in
107 * order to obtain the rms error in the force for a system of N randomly
108 * distributed particles in a cubic box. For the real space error, the
109 * estimate in @cite kolafa92a is used.
110 *
111 * Parameter ranges if not given explicitly in the constructor:
112 * - @p mesh is set up such that the number of mesh points is equal to the
113 * number of magnetic dipolar particles
114 * - @p cao explores all possible values
115 * - @p alpha_L is tuned for each tuple (@p r_cut_iL, @p mesh, @p cao) and
116 * calculated assuming that the error contributions of real and reciprocal
117 * space should be equal
118 *
119 * After checking if the total error lies below the target accuracy,
120 * the time needed for one force calculation is measured. Parameters
121 * that minimize the runtime are kept.
122 *
123 * The function is based on routines of the program HE_Q.cpp for charges
124 * written by M. Deserno.
125 */
126 virtual void tune() = 0;
127
128 /** Compute the k-space part of energies. */
129 virtual double long_range_energy(ParticleRange const &particles) = 0;
130
131 /** Compute the k-space part of forces. */
132 virtual void add_long_range_forces(ParticleRange const &particles) = 0;
133
134 /** Calculate real-space contribution of p3m dipolar pair forces and torques.
135 * If NPT is compiled in, update the NpT virial.
136 */
137 inline ParticleForce pair_force(Particle const &p1, Particle const &p2,
138 Utils::Vector3d const &d, double dist2,
139 double dist) const {
140 if ((p1.dipm() == 0.) || (p2.dipm() == 0.) || dist >= dp3m_params.r_cut ||
141 dist <= 0.)
142 return {};
143
144 auto const dip1 = p1.calc_dip();
145 auto const dip2 = p2.calc_dip();
146 auto const alpsq = dp3m_params.alpha * dp3m_params.alpha;
147 auto const adist = dp3m_params.alpha * dist;
148#if USE_ERFC_APPROXIMATION
149 auto const erfc_part_ri = Utils::AS_erfc_part(adist) / dist;
150#else
151 auto const erfc_part_ri = erfc(adist) / dist;
152#endif
153
154 // Calculate scalar multiplications for vectors mi, mj, rij
155 auto const mimj = dip1 * dip2;
156
157 auto const mir = dip1 * d;
158 auto const mjr = dip2 * d;
159
160 auto const coeff = 2. * dp3m_params.alpha * std::numbers::inv_sqrtpi;
161 auto const dist2i = 1. / dist2;
162 auto const exp_adist2 = exp(-Utils::sqr(adist));
163
164 auto const B_r = (dp3m_params.accuracy > 5e-06)
165 ? (erfc_part_ri + coeff) * exp_adist2 * dist2i
166 : (erfc(adist) / dist + coeff * exp_adist2) * dist2i;
167
168 auto const common_term = alpsq * coeff * exp_adist2;
169 auto const C_r = dist2i * (3. * B_r + 2. * common_term);
170 auto const D_r = dist2i * (5. * C_r + 4. * common_term * alpsq);
171
172 // Calculate real-space forces
173 auto const force = prefactor * ((mimj * d + dip1 * mjr + dip2 * mir) * C_r -
174 mir * mjr * D_r * d);
175
176 // Calculate vector multiplications for vectors mi, mj, rij
177 auto const mixmj = vector_product(dip1, dip2);
178 auto const mixr = vector_product(dip1, d);
179
180 // Calculate real-space torques
181 auto const torque = prefactor * (-mixmj * B_r + mixr * (mjr * C_r));
182#ifdef NPT
183#if USE_ERFC_APPROXIMATION
184 auto const fac = prefactor * p1.dipm() * p2.dipm() * exp_adist2;
185#else
186 auto const fac = prefactor * p1.dipm() * p2.dipm();
187#endif
188 auto const energy = fac * (mimj * B_r - mir * mjr * C_r);
190#endif // NPT
191 return ParticleForce{force, torque};
192 }
193
194 /** Calculate real-space contribution of dipolar pair energy. */
195 inline double pair_energy(Particle const &p1, Particle const &p2,
196 Utils::Vector3d const &d, double dist2,
197 double dist) const {
198 if ((p1.dipm() == 0.) || (p2.dipm() == 0.) || dist >= dp3m_params.r_cut ||
199 dist <= 0.)
200 return {};
201
202 auto const dip1 = p1.calc_dip();
203 auto const dip2 = p2.calc_dip();
204
205 auto const alpsq = dp3m_params.alpha * dp3m_params.alpha;
206 auto const adist = dp3m_params.alpha * dist;
207
208#if USE_ERFC_APPROXIMATION
209 auto const erfc_part_ri = Utils::AS_erfc_part(adist) / dist;
210#else
211 auto const erfc_part_ri = erfc(adist) / dist;
212#endif
213
214 // Calculate scalar multiplications for vectors mi, mj, rij
215 auto const mimj = dip1 * dip2;
216 auto const mir = dip1 * d;
217 auto const mjr = dip2 * d;
218
219 auto const coeff = 2. * dp3m_params.alpha * std::numbers::inv_sqrtpi;
220 auto const dist2i = 1. / dist2;
221 auto const exp_adist2 = exp(-Utils::sqr(adist));
222
223 auto const B_r = (dp3m_params.accuracy > 5e-06)
224 ? dist2i * (erfc_part_ri + coeff) * exp_adist2
225 : dist2i * (erfc(adist) / dist + coeff * exp_adist2);
226 auto const C_r = (3. * B_r + 2. * alpsq * coeff * exp_adist2) * dist2i;
227
228 return prefactor * (mimj * B_r - mir * mjr * C_r);
229 }
230
231protected:
232 /** Calculate self-energy in k-space. */
233 virtual double calc_average_self_energy_k_space() const = 0;
234
235 /** Calculate energy correction that minimizes the error.
236 * This quantity is only calculated once and is cached.
237 */
238 virtual void calc_energy_correction() = 0;
239
240 /** Calculate the influence function for the dipolar forces and torques. */
242
243 /** Calculate the influence function for the dipolar energy. */
245
246 /** Compute the dipolar surface terms */
247 virtual double calc_surface_term(bool force_flag, bool energy_flag,
248 ParticleRange const &particles) = 0;
249
250 /** Checks for correctness of the k-space cutoff. */
251 void sanity_checks_boxl() const;
252 void sanity_checks_node_grid() const;
253 void sanity_checks_periodicity() const;
254 void sanity_checks_cell_structure() const;
255
256 virtual void scaleby_box_l() = 0;
257
258#ifdef NPT
259 /** Update the NpT virial */
260 virtual void npt_add_virial_contribution(double energy) const = 0;
261#endif
262};
263
264#endif // DP3M
Vector implementation and trait types for boost qvm interoperability.
double prefactor
Magnetostatics prefactor.
A range of particles.
This file contains the defaults for ESPResSo.
__device__ void vector_product(float const *a, float const *b, float *out)
constexpr T AS_erfc_part(T d)
Approximate by applying a formula from chapter 7.
DEVICE_QUALIFIER constexpr T sqr(T x)
Calculates the SQuaRe of x.
Definition sqr.hpp:28
Common functions for dipolar and charge P3M.
Dipolar P3M solver.
Definition dp3m.hpp:54
void sanity_checks_boxl() const
Checks for correctness of the k-space cutoff.
Definition dp3m.cpp:841
virtual void init()=0
Recalculate all derived parameters.
virtual bool is_gpu() const noexcept=0
void sanity_checks() const
Definition dp3m.hpp:78
virtual void scaleby_box_l()=0
virtual void tune()=0
Tune dipolar P3M parameters to desired accuracy.
virtual void npt_add_virial_contribution(double energy) const =0
Update the NpT virial.
virtual void calc_energy_correction()=0
Calculate energy correction that minimizes the error.
void on_node_grid_change() const
Definition dp3m.hpp:69
virtual bool is_double_precision() const noexcept=0
void on_boxl_change()
Recalculate all box-length-dependent parameters.
Definition dp3m.hpp:68
void sanity_checks_cell_structure() const
Definition dp3m.cpp:875
virtual double long_range_energy(ParticleRange const &particles)=0
Compute the k-space part of energies.
virtual ~DipolarP3M()=default
virtual void dipole_assign(ParticleRange const &particles)=0
Assign the physical dipoles using the tabulated assignment function.
P3MParameters const & dp3m_params
Definition dp3m.hpp:55
virtual void on_activation()=0
virtual void calc_influence_function_force()=0
Calculate the influence function for the dipolar forces and torques.
virtual double calc_average_self_energy_k_space() const =0
Calculate self-energy in k-space.
virtual bool is_tuned() const noexcept=0
void sanity_checks_periodicity() const
Definition dp3m.cpp:867
DipolarP3M(P3MParameters const &dp3m_params)
Definition dp3m.hpp:58
void sanity_checks_node_grid() const
Definition dp3m.cpp:890
double pair_energy(Particle const &p1, Particle const &p2, Utils::Vector3d const &d, double dist2, double dist) const
Calculate real-space contribution of dipolar pair energy.
Definition dp3m.hpp:195
void on_periodicity_change() const
Definition dp3m.hpp:70
virtual void count_magnetic_particles()=0
Count the number of magnetic particles and calculate the sum of the squared dipole moments.
virtual double calc_surface_term(bool force_flag, bool energy_flag, ParticleRange const &particles)=0
Compute the dipolar surface terms.
void on_cell_structure_change()
Definition dp3m.hpp:71
virtual void calc_influence_function_energy()=0
Calculate the influence function for the dipolar energy.
virtual void add_long_range_forces(ParticleRange const &particles)=0
Compute the k-space part of forces.
ParticleForce pair_force(Particle const &p1, Particle const &p2, Utils::Vector3d const &d, double dist2, double dist) const
Calculate real-space contribution of p3m dipolar pair forces and torques.
Definition dp3m.hpp:137
Structure to hold P3M parameters and some dependent variables.
double alpha
unscaled alpha_L for use with fast inline functions only
double accuracy
accuracy of the actual parameter set.
double r_cut
unscaled r_cut_iL for use with fast inline functions only
Force information on a particle.
Definition Particle.hpp:290
Struct holding all information for one particle.
Definition Particle.hpp:395
auto calc_dip() const
Definition Particle.hpp:495
auto const & dipm() const
Definition Particle.hpp:493