ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
dp3m.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-2022 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#include "p3m/fft.hpp"
43#include "p3m/interpolation.hpp"
44#include "p3m/send_mesh.hpp"
45
46#include "Particle.hpp"
47#include "ParticleRange.hpp"
48
49#include <utils/Vector.hpp>
50#include <utils/constants.hpp>
52
53#include <array>
54#include <cmath>
55#include <vector>
56
57#ifdef NPT
58/** Update the NpT virial */
59void npt_add_virial_magnetic_contribution(double energy);
60#endif
61
63 explicit dp3m_data_struct(P3MParameters &&parameters)
64 : p3m_data_struct_base{std::move(parameters)} {}
65
66 /** local mesh. */
68 /** real space mesh (local) for CA/FFT. */
70 /** real space mesh (local) for CA/FFT of the dipolar field. */
71 std::array<fft_vector<double>, 3> rs_mesh_dip;
72 /** k-space mesh (local) for k-space calculation and FFT. */
73 std::vector<double> ks_mesh;
74
75 /** number of dipolar particles (only on head node). */
76 int sum_dip_part = 0;
77 /** Sum of square of magnetic dipoles (only on head node). */
78 double sum_mu2 = 0.;
79
80 /** position shift for calculation of first assignment mesh point. */
81 double pos_shift = 0.;
82
84
85 /** send/recv mesh sizes */
87
88 /** cached k-space self-energy correction */
89 double energy_correction = 0.;
90
92};
93
94/** @brief Dipolar P3M solver. */
95struct DipolarP3M : public Dipoles::Actor<DipolarP3M> {
96 /** Dipolar P3M parameters. */
98
99 /** Magnetostatics prefactor. */
102
103 DipolarP3M(P3MParameters &&parameters, double prefactor, int tune_timings,
104 bool tune_verbose);
105
108 tune();
109 }
110 /** @brief Recalculate all box-length-dependent parameters. */
111 void on_boxl_change() { scaleby_box_l(); }
112 void on_node_grid_change() const { sanity_checks_node_grid(); }
113 void on_periodicity_change() const { sanity_checks_periodicity(); }
115 sanity_checks_cell_structure();
116 init();
117 }
118 /** @brief Recalculate all derived parameters. */
119 void init();
120
121 void sanity_checks() const {
122 sanity_checks_boxl();
123 sanity_checks_node_grid();
124 sanity_checks_periodicity();
125 sanity_checks_cell_structure();
126 }
127
128 /**
129 * @brief Count the number of magnetic particles and calculate
130 * the sum of the squared dipole moments.
131 */
133
134 /** Assign the physical dipoles using the tabulated assignment function. */
135 void dipole_assign(ParticleRange const &particles);
136
137 /**
138 * @brief Tune dipolar P3M parameters to desired accuracy.
139 *
140 * The parameters
141 * @ref P3MParameters::mesh "mesh",
142 * @ref P3MParameters::cao "cao",
143 * @ref P3MParameters::r_cut_iL "r_cut_iL" and
144 * @ref P3MParameters::alpha_L "alpha_L" are tuned to obtain the target
145 * @ref P3MParameters::accuracy "accuracy" in optimal time.
146 * These parameters are stored in the @ref dp3m object.
147 *
148 * The function utilizes the analytic expression of the error estimate
149 * for the dipolar P3M method in the paper of @cite cerda08d in
150 * order to obtain the rms error in the force for a system of N randomly
151 * distributed particles in a cubic box. For the real space error, the
152 * estimate in @cite kolafa92a is used.
153 *
154 * Parameter ranges if not given explicitly in the constructor:
155 * - @p mesh is set up such that the number of mesh points is equal to the
156 * number of magnetic dipolar particles
157 * - @p cao explores all possible values
158 * - @p alpha_L is tuned for each tuple (@p r_cut_iL, @p mesh, @p cao) and
159 * calculated assuming that the error contributions of real and reciprocal
160 * space should be equal
161 *
162 * After checking if the total error lies below the target accuracy,
163 * the time needed for one force calculation is measured. Parameters
164 * that minimize the runtime are kept.
165 *
166 * The function is based on routines of the program HE_Q.cpp for charges
167 * written by M. Deserno.
168 */
169 void tune();
170 bool is_tuned() const { return m_is_tuned; }
171
172 /** Compute the k-space part of energies. */
173 double long_range_energy(ParticleRange const &particles) {
174 return long_range_kernel(false, true, particles);
175 }
176
177 /** Compute the k-space part of forces. */
178 void add_long_range_forces(ParticleRange const &particles) {
179 long_range_kernel(true, false, particles);
180 }
181
182 /** Compute the k-space part of forces and energies. */
183 double long_range_kernel(bool force_flag, bool energy_flag,
184 ParticleRange const &particles);
185
186 /** Calculate real-space contribution of p3m dipolar pair forces and torques.
187 * If NPT is compiled in, update the NpT virial.
188 */
189 inline ParticleForce pair_force(Particle const &p1, Particle const &p2,
190 Utils::Vector3d const &d, double dist2,
191 double dist) const {
192 if ((p1.dipm() == 0.) || (p2.dipm() == 0.) || dist >= dp3m.params.r_cut ||
193 dist <= 0.)
194 return {};
195
196 auto const dip1 = p1.calc_dip();
197 auto const dip2 = p2.calc_dip();
198 auto const alpsq = dp3m.params.alpha * dp3m.params.alpha;
199 auto const adist = dp3m.params.alpha * dist;
200#if USE_ERFC_APPROXIMATION
201 auto const erfc_part_ri = Utils::AS_erfc_part(adist) / dist;
202#else
203 auto const erfc_part_ri = erfc(adist) / dist;
204#endif
205
206 // Calculate scalar multiplications for vectors mi, mj, rij
207 auto const mimj = dip1 * dip2;
208
209 auto const mir = dip1 * d;
210 auto const mjr = dip2 * d;
211
212 auto const coeff = 2. * dp3m.params.alpha * Utils::sqrt_pi_i();
213 auto const dist2i = 1. / dist2;
214 auto const exp_adist2 = exp(-Utils::sqr(adist));
215
216 auto const B_r = (dp3m.params.accuracy > 5e-06)
217 ? (erfc_part_ri + coeff) * exp_adist2 * dist2i
218 : (erfc(adist) / dist + coeff * exp_adist2) * dist2i;
219
220 auto const common_term = alpsq * coeff * exp_adist2;
221 auto const C_r = dist2i * (3. * B_r + 2. * common_term);
222 auto const D_r = dist2i * (5. * C_r + 4. * common_term * alpsq);
223
224 // Calculate real-space forces
225 auto const force = prefactor * ((mimj * d + dip1 * mjr + dip2 * mir) * C_r -
226 mir * mjr * D_r * d);
227
228 // Calculate vector multiplications for vectors mi, mj, rij
229 auto const mixmj = vector_product(dip1, dip2);
230 auto const mixr = vector_product(dip1, d);
231
232 // Calculate real-space torques
233 auto const torque = prefactor * (-mixmj * B_r + mixr * (mjr * C_r));
234#ifdef NPT
235#if USE_ERFC_APPROXIMATION
236 auto const fac = prefactor * p1.dipm() * p2.dipm() * exp_adist2;
237#else
238 auto const fac = prefactor * p1.dipm() * p2.dipm();
239#endif
240 auto const energy = fac * (mimj * B_r - mir * mjr * C_r);
242#endif // NPT
243 return ParticleForce{force, torque};
244 }
245
246 /** Calculate real-space contribution of dipolar pair energy. */
247 inline double pair_energy(Particle const &p1, Particle const &p2,
248 Utils::Vector3d const &d, double dist2,
249 double dist) const {
250 if ((p1.dipm() == 0.) || (p2.dipm() == 0.) || dist >= dp3m.params.r_cut ||
251 dist <= 0.)
252 return {};
253
254 auto const dip1 = p1.calc_dip();
255 auto const dip2 = p2.calc_dip();
256
257 auto const alpsq = dp3m.params.alpha * dp3m.params.alpha;
258 auto const adist = dp3m.params.alpha * dist;
259
260#if USE_ERFC_APPROXIMATION
261 auto const erfc_part_ri = Utils::AS_erfc_part(adist) / dist;
262#else
263 auto const erfc_part_ri = erfc(adist) / dist;
264#endif
265
266 // Calculate scalar multiplications for vectors mi, mj, rij
267 auto const mimj = dip1 * dip2;
268 auto const mir = dip1 * d;
269 auto const mjr = dip2 * d;
270
271 auto const coeff = 2. * dp3m.params.alpha * Utils::sqrt_pi_i();
272 auto const dist2i = 1. / dist2;
273 auto const exp_adist2 = exp(-Utils::sqr(adist));
274
275 auto const B_r = (dp3m.params.accuracy > 5e-06)
276 ? dist2i * (erfc_part_ri + coeff) * exp_adist2
277 : dist2i * (erfc(adist) / dist + coeff * exp_adist2);
278 auto const C_r = (3. * B_r + 2. * alpsq * coeff * exp_adist2) * dist2i;
279
280 return prefactor * (mimj * B_r - mir * mjr * C_r);
281 }
282
283private:
284 bool m_is_tuned;
285
286 /** Calculate self-energy in k-space. */
287 double calc_average_self_energy_k_space() const;
288
289 /** Calculate energy correction that minimizes the error.
290 * This quantity is only calculated once and is cached.
291 */
292 void calc_energy_correction();
293
294 /** Calculate the influence function for the dipolar forces and torques. */
295 void calc_influence_function_force();
296
297 /** Calculate the influence function for the dipolar energy. */
298 void calc_influence_function_energy();
299
300 /** Compute the dipolar surface terms */
301 double calc_surface_term(bool force_flag, bool energy_flag,
302 ParticleRange const &particles);
303
304 /** Checks for correctness of the k-space cutoff. */
305 void sanity_checks_boxl() const;
306 void sanity_checks_node_grid() const;
307 void sanity_checks_periodicity() const;
308 void sanity_checks_cell_structure() const;
309
310 void scaleby_box_l();
311};
312
313#endif // DP3M
Vector implementation and trait types for boost qvm interoperability.
__global__ float float * torque
__global__ float * force
double prefactor
Magnetostatics prefactor.
A range of particles.
Cache for interpolation weights.
Structure for send/recv meshes.
Definition send_mesh.hpp:38
Common functions for dipolar and charge P3M.
This file contains the defaults for ESPResSo.
__device__ void vector_product(float const *a, float const *b, float *out)
void npt_add_virial_magnetic_contribution(double energy)
Update the NpT virial.
Definition dp3m.cpp:965
Routines, row decomposition, data structures and communication for the 3D-FFT.
std::vector< T, fft_allocator< T > > fft_vector
Definition fft.hpp:87
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:26
DEVICE_QUALIFIER constexpr T sqrt_pi_i()
One over square root of pi.
Definition constants.hpp:43
Dipolar P3M solver.
Definition dp3m.hpp:95
void init()
Recalculate all derived parameters.
Definition dp3m.cpp:119
int tune_timings
Magnetostatics prefactor.
Definition dp3m.hpp:100
void sanity_checks() const
Definition dp3m.hpp:121
void dipole_assign(ParticleRange const &particles)
Assign the physical dipoles using the tabulated assignment function.
Definition dp3m.cpp:190
void on_node_grid_change() const
Definition dp3m.hpp:112
dp3m_data_struct dp3m
Dipolar P3M parameters.
Definition dp3m.hpp:97
void on_boxl_change()
Recalculate all box-length-dependent parameters.
Definition dp3m.hpp:111
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:257
void count_magnetic_particles()
Count the number of magnetic particles and calculate the sum of the squared dipole moments.
Definition dp3m.cpp:76
double long_range_energy(ParticleRange const &particles)
Compute the k-space part of energies.
Definition dp3m.hpp:173
void add_long_range_forces(ParticleRange const &particles)
Compute the k-space part of forces.
Definition dp3m.hpp:178
void tune()
Tune dipolar P3M parameters to desired accuracy.
Definition dp3m.cpp:719
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:247
bool is_tuned() const
Definition dp3m.hpp:170
void on_periodicity_change() const
Definition dp3m.hpp:113
bool tune_verbose
Definition dp3m.hpp:101
void on_cell_structure_change()
Definition dp3m.hpp:114
void on_activation()
Definition dp3m.hpp:106
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:189
Structure for local mesh parameters.
Definition common.hpp:183
Structure to hold P3M parameters and some dependent variables.
Definition common.hpp:67
double alpha
unscaled alpha_L for use with fast inline functions only
Definition common.hpp:96
double accuracy
accuracy of the actual parameter set.
Definition common.hpp:84
double r_cut
unscaled r_cut_iL for use with fast inline functions only
Definition common.hpp:99
Force information on a particle.
Definition Particle.hpp:290
Struct holding all information for one particle.
Definition Particle.hpp:393
auto calc_dip() const
Definition Particle.hpp:493
auto const & dipm() const
Definition Particle.hpp:491
p3m_send_mesh sm
send/recv mesh sizes
Definition dp3m.hpp:86
P3MLocalMesh local_mesh
local mesh.
Definition dp3m.hpp:67
fft_data_struct fft
Definition dp3m.hpp:91
dp3m_data_struct(P3MParameters &&parameters)
Definition dp3m.hpp:63
fft_vector< double > rs_mesh
real space mesh (local) for CA/FFT.
Definition dp3m.hpp:69
double sum_mu2
Sum of square of magnetic dipoles (only on head node).
Definition dp3m.hpp:78
std::array< fft_vector< double >, 3 > rs_mesh_dip
real space mesh (local) for CA/FFT of the dipolar field.
Definition dp3m.hpp:71
int sum_dip_part
number of dipolar particles (only on head node).
Definition dp3m.hpp:76
p3m_interpolation_cache inter_weights
Definition dp3m.hpp:83
double pos_shift
position shift for calculation of first assignment mesh point.
Definition dp3m.hpp:81
std::vector< double > ks_mesh
k-space mesh (local) for k-space calculation and FFT.
Definition dp3m.hpp:73
double energy_correction
cached k-space self-energy correction
Definition dp3m.hpp:89
Information about the three one dimensional FFTs and how the nodes have to communicate inbetween.
Definition fft.hpp:152
P3MParameters params