ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
statistics.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-2026 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/** \file
22 * Statistical tools to analyze simulations.
23 *
24 * The corresponding header file is statistics.hpp.
25 */
26
28
29#include "BoxGeometry.hpp"
30#include "Particle.hpp"
32#include "communication.hpp"
34#include "system/System.hpp"
35
36#include <utils/Vector.hpp>
37#include <utils/contains.hpp>
38#include <utils/math/sqr.hpp>
40
41#include <boost/mpi/collectives/broadcast.hpp>
42#include <boost/mpi/collectives/reduce.hpp>
43
44#include <cassert>
45#include <cmath>
46#include <cstdlib>
47#include <functional>
48#include <limits>
49#include <numbers>
50#include <tuple>
51#include <utility>
52#include <vector>
53
55/** @brief Serialize @c std::tuple. */
56template <typename Archive, typename... T>
57void serialize(Archive &ar, std::tuple<T...> &pack, unsigned int const) {
58 std::apply([&](auto &...item) { ((ar & item), ...); }, pack);
59}
60} // namespace boost::serialization
61
62/** @brief Decay a tuple of only 1 type to that type. */
63template <typename... F> struct DecayTupleResult {
64 using type = std::tuple<std::invoke_result_t<F, Particle const &>...>;
65};
66
67template <typename F> struct DecayTupleResult<F> {
68 using type = std::invoke_result_t<F, Particle const &>;
69};
70
71/**
72 * @brief Gather particle traits to MPI rank 0.
73 * When only one trait is requested, return a vector of type T.
74 * When multiple traits are requested, return a vector of tuples.
75 */
76template <class... Trait>
78 std::vector<int> const &p_types,
79 Trait &&...trait) {
80 std::vector<typename DecayTupleResult<Trait...>::type> buffer{};
81
82 for (auto const &p : system.cell_structure->local_particles()) {
83 if (Utils::contains(p_types, p.type())) {
84 buffer.emplace_back(trait(p)...);
85 }
86 }
87
89 if (::comm_cart.rank() != 0) {
90 buffer.clear();
91 }
92 return buffer;
93}
94
95double mindist(System::System const &system, std::vector<int> const &set1,
96 std::vector<int> const &set2) {
97 using Utils::contains;
98
99 std::vector<int> buf_type{};
100 std::vector<Utils::Vector3d> buf_pos{};
101
102 auto const &box_geo = *system.box_geo;
103 auto const accept_all = set1.empty() or set2.empty();
104 for (auto const &p : system.cell_structure->local_particles()) {
105 if (accept_all or contains(set1, p.type()) or contains(set2, p.type())) {
106 buf_type.emplace_back(p.type());
107 buf_pos.emplace_back(box_geo.unfolded_position(p.pos(), p.image_box()));
108 }
109 }
110
113 if (::comm_cart.rank() != 0) {
114 buf_type.clear();
115 buf_pos.clear();
116 }
117
118 auto mindist_sq = std::numeric_limits<double>::infinity();
119 for (std::size_t j = 0ul; j < buf_type.size(); ++j) {
120 /* check which sets particle j belongs to (bit 0: set1, bit1: set2) */
121 auto in_set = 0;
122 if (set1.empty() || contains(set1, buf_type[j]))
123 in_set = 1;
124 if (set2.empty() || contains(set2, buf_type[j]))
125 in_set |= 2;
126 assert(in_set != 0);
127
128 for (auto i = j + 1ul; i != buf_type.size(); ++i) {
129 /* accept a pair if particle j is in set1 and particle i in set2 or vice
130 * versa. */
131 if (((in_set & 1) and (set2.empty() or contains(set2, buf_type[i]))) or
132 ((in_set & 2) and (set1.empty() or contains(set1, buf_type[i])))) {
133 mindist_sq = std::min(
134 mindist_sq, box_geo.get_mi_vector(buf_pos[j], buf_pos[i]).norm2());
135 }
136 }
137 }
138
139 return std::sqrt(mindist_sq);
140}
141
144 bool include_lbfluid) {
146 if (include_particles) {
148 *(system.cell_structure),
149 [](Utils::Vector3d &acc, Particle const &p) {
150 acc += p.mass() * p.v();
151 },
152 [](Utils::Vector3d &acc, Utils::Vector3d const &v) { acc = acc + v; });
153 }
154 if (include_lbfluid and system.lb.is_solver_set()) {
155 momentum += system.lb.get_momentum() * system.lb.get_lattice_speed();
156 }
157 return momentum;
158}
159
161 auto const &box_geo = *system.box_geo;
162 auto const &cell_structure = *system.cell_structure;
164 double local_mass = 0.;
165
166 for (auto const &p : cell_structure.local_particles()) {
167 if ((p.type() == p_type or p_type == -1) and not p.is_virtual()) {
168 local_com += box_geo.unfolded_position(p.pos(), p.image_box()) * p.mass();
169 local_mass += p.mass();
170 }
171 }
173 double mass = 1.; // placeholder value to avoid division by zero
174 boost::mpi::reduce(::comm_cart, local_com, com, std::plus<>(), 0);
175 boost::mpi::reduce(::comm_cart, local_mass, mass, std::plus<>(), 0);
176 return com / mass;
177}
178
180 auto const &box_geo = *system.box_geo;
181 auto const &cell_structure = *system.cell_structure;
183
184 for (auto const &p : cell_structure.local_particles()) {
185 if ((p.type() == p_type or p_type == -1) and not p.is_virtual()) {
186 auto const pos = box_geo.unfolded_position(p.pos(), p.image_box());
187 am += p.mass() * vector_product(pos, p.v());
188 }
189 }
190 return am;
191}
192
194 std::vector<int> const &p_types) {
195 auto const &box_geo = *system.box_geo;
196 auto const trait_pos = [&box_geo](Particle const &p) {
197 return box_geo.unfolded_position(p.pos(), p.image_box());
198 };
200
202 if (::comm_cart.rank() == 0) {
203 auto const center =
204 std::accumulate(buf_pos.begin(), buf_pos.end(), Utils::Vector3d{}) /
205 static_cast<double>(buf_pos.size());
206 // compute covariance matrix
207 for (unsigned int i = 0u; i < 3u; ++i) {
208 for (unsigned int j = 0u; j < 3u; ++j) {
209 if (i > j) {
210 mat[i * 3u + j] = mat[j * 3u + i];
211 } else {
212 mat[i * 3u + j] = std::accumulate(
213 buf_pos.begin(), buf_pos.end(), 0.,
214 [i, j, &center](double acc, Utils::Vector3d const &pos) {
215 return acc + (pos[i] - center[i]) * (pos[j] - center[j]);
216 });
217 }
218 }
219 }
220 mat /= static_cast<double>(buf_pos.size());
221 }
222 return mat;
223}
224
226 int p_type) {
227 auto const &box_geo = *system.box_geo;
228 auto const &cell_structure = *system.cell_structure;
231 boost::mpi::broadcast(::comm_cart, com, 0);
232
233 for (auto const &p : cell_structure.local_particles()) {
234 if (p.type() == p_type and not p.is_virtual()) {
235 auto const pos = box_geo.unfolded_position(p.pos(), p.image_box()) - com;
236 auto const mass = p.mass();
237 mat[0] += mass * (pos[1] * pos[1] + pos[2] * pos[2]);
238 mat[4] += mass * (pos[0] * pos[0] + pos[2] * pos[2]);
239 mat[8] += mass * (pos[0] * pos[0] + pos[1] * pos[1]);
240 mat[1] -= mass * (pos[0] * pos[1]);
241 mat[2] -= mass * (pos[0] * pos[2]);
242 mat[5] -= mass * (pos[1] * pos[2]);
243 }
244 }
245 /* use symmetry */
246 mat[3] = mat[1];
247 mat[6] = mat[2];
248 mat[7] = mat[5];
249 return mat;
250}
251
252std::vector<int> nbhood(System::System const &system,
253 Utils::Vector3d const &pos, double dist) {
254 std::vector<int> buf_pid{};
255 auto const dist_sq = dist * dist;
256 auto const &box_geo = *system.box_geo;
257
258 for (auto const &p : system.cell_structure->local_particles()) {
259 auto const r_sq = box_geo.get_mi_vector(pos, p.pos()).norm2();
260 if (r_sq < dist_sq) {
261 buf_pid.push_back(p.id());
262 }
263 }
264
266 if (::comm_cart.rank() != 0) {
267 buf_pid.clear();
268 }
269
270 return buf_pid;
271}
272
273std::vector<std::vector<double>>
275 std::vector<int> const &p1_types,
276 std::vector<int> const &p2_types, double r_min,
277 double r_max, int r_bins, bool log_flag, bool int_flag) {
278
279 auto const &box_geo = *system.box_geo;
280 auto const trait_id = [](Particle const &p) { return p.id(); };
281 auto const trait_pos = [&box_geo](Particle const &p) {
282 return box_geo.unfolded_position(p.pos(), p.image_box());
283 };
284 auto const buf1 =
286 auto const buf2 =
288 auto const r_max2 = Utils::sqr(r_max);
289 auto const r_min2 = Utils::sqr(r_min);
290 auto const start_dist2 = Utils::sqr(r_max + 1.);
291 auto const inv_bin_width =
292 (log_flag) ? static_cast<double>(r_bins) / std::log(r_max / r_min)
293 : static_cast<double>(r_bins) / (r_max - r_min);
294
295 long cnt = 0;
296 double low = 0.0;
297 std::vector<double> distribution(r_bins);
298
299 for (auto const &[pid1, pos1] : buf1) {
300 auto min_dist2 = start_dist2;
301 /* particle loop: p2_types */
302 for (auto const &[pid2, pos2] : buf2) {
303 if (pid1 != pid2) {
304 auto const act_dist2 = box_geo.get_mi_vector(pos1, pos2).norm2();
305 if (act_dist2 < min_dist2) {
307 }
308 }
309 }
310 if (min_dist2 <= r_max2) {
311 if (min_dist2 >= r_min2) {
312 auto const min_dist = std::sqrt(min_dist2);
313 /* calculate bin index */
314 auto const ind = static_cast<int>(
315 ((log_flag) ? std::log(min_dist / r_min) : (min_dist - r_min)) *
317 if (ind >= 0 and ind < r_bins) {
318 distribution[ind] += 1.0;
319 }
320 } else {
321 low += 1.0;
322 }
323 }
324 cnt++;
325 }
326
327 if (cnt != 0) {
328 // normalization
329 low /= static_cast<double>(cnt);
330 for (int i = 0; i < r_bins; i++) {
331 distribution[i] /= static_cast<double>(cnt);
332 }
333
334 // integration
335 if (int_flag) {
336 distribution[0] += low;
337 for (int i = 0; i < r_bins - 1; i++)
338 distribution[i + 1] += distribution[i];
339 }
340 }
341
342 std::vector<double> radii(r_bins);
343 if (log_flag) {
344 auto const log_fac = std::pow(r_max / r_min, 1. / r_bins);
345 radii[0] = r_min * std::sqrt(log_fac);
346 for (int i = 1; i < r_bins; ++i) {
347 radii[i] = radii[i - 1] * log_fac;
348 }
349 } else {
350 auto const bin_width = (r_max - r_min) / static_cast<double>(r_bins);
351 for (int i = 0; i < r_bins; ++i) {
352 radii[i] = r_min + bin_width / 2. + static_cast<double>(i) * bin_width;
353 }
354 }
355
356 return {radii, distribution};
357}
358
359std::vector<std::vector<double>>
360structure_factor(System::System const &system, std::vector<int> const &p_types,
361 int order) {
362 auto const &box_geo = *system.box_geo;
363 auto const trait_pos = [&box_geo](Particle const &p) {
364 return box_geo.unfolded_position(p.pos(), p.image_box());
365 };
367 auto const order_sq = Utils::sqr(static_cast<std::size_t>(order));
368 auto const twoPI_L = 2. * std::numbers::pi * system.box_geo->length_inv()[0];
369 std::vector<double> ff(2ul * order_sq + 1ul);
370 std::vector<double> wavevectors;
371 std::vector<double> intensities;
372
373 if (::comm_cart.rank() == 0) {
374 for (int i = 0; i <= order; i++) {
375 for (int j = -order; j <= order; j++) {
376 for (int k = -order; k <= order; k++) {
377 auto const n = i * i + j * j + k * k;
378 if ((static_cast<std::size_t>(n) <= order_sq) && (n >= 1)) {
379 double C_sum = 0.0, S_sum = 0.0;
380 for (auto const &pos : buf_pos) {
381 auto const qr = twoPI_L * (Utils::Vector3i{{i, j, k}} * pos);
382 C_sum += cos(qr);
383 S_sum += sin(qr);
384 }
385 ff[2 * n - 2] += C_sum * C_sum + S_sum * S_sum;
386 ff[2 * n - 1]++;
387 }
388 }
389 }
390 }
391
392 std::size_t length = 0;
393 for (std::size_t qi = 0; qi < order_sq; qi++) {
394 if (ff[2 * qi + 1] != 0) {
395 ff[2 * qi] /= static_cast<double>(buf_pos.size()) * ff[2 * qi + 1];
396 ++length;
397 }
398 }
399
400 wavevectors.resize(length);
401 intensities.resize(length);
402
403 int cnt = 0;
404 for (std::size_t i = 0; i < order_sq; i++) {
405 if (ff[2 * i + 1] != 0) {
406 wavevectors[cnt] = twoPI_L * std::sqrt(static_cast<double>(i + 1));
407 intensities[cnt] = ff[2 * i];
408 cnt++;
409 }
410 }
411 }
412
413 return {std::move(wavevectors), std::move(intensities)};
414}
Vector implementation and trait types for boost qvm interoperability.
Main system class.
cudaStream_t stream[1]
CUDA streams for parallel computing on CPU and GPU.
boost::mpi::communicator comm_cart
The communicator.
__device__ void vector_product(float const *a, float const *b, float *out)
void gather_buffer(std::vector< T, Allocator > &buffer, boost::mpi::communicator const &comm, int root=0)
Gather buffer with different size on each node.
DEVICE_QUALIFIER constexpr T sqr(T x)
Calculates the SQuaRe of x.
Definition sqr.hpp:28
bool contains(Range &&rng, T const &value)
Check whether a range contains a value.
Definition contains.hpp:36
void serialize(Archive &ar, std::tuple< T... > &pack, unsigned int const)
Serialize std::tuple.
Utils::Vector3d center_of_mass(System::System const &system, int p_type)
Calculate the center of mass of particles of a certain type.
Utils::Vector3d angular_momentum(System::System const &system, int p_type)
Calculate the angular momentum of particles of a certain type.
std::vector< int > nbhood(System::System const &system, Utils::Vector3d const &pos, double dist)
Find all particles within a given radius dist around a position pos.
Utils::Vector9d moment_of_inertia_matrix(System::System const &system, int p_type)
Calculate the moment of inertia of particles of a certain type.
Utils::Vector9d gyration_tensor(System::System const &system, std::vector< int > const &p_types)
Calculate the gyration tensor of particles of certain types.
Utils::Vector3d calc_linear_momentum(System::System const &system, bool include_particles, bool include_lbfluid)
Calculate total momentum of the system (particles & LB fluid).
std::vector< std::vector< double > > structure_factor(System::System const &system, std::vector< int > const &p_types, int order)
Calculate the spherically averaged structure factor.
double mindist(System::System const &system, std::vector< int > const &set1, std::vector< int > const &set2)
Calculate the minimal distance of two particles with types in set1 and set2, respectively.
static auto gather_traits_for_types(System::System const &system, std::vector< int > const &p_types, Trait &&...trait)
Gather particle traits to MPI rank 0.
std::vector< std::vector< double > > calc_part_distribution(System::System const &system, std::vector< int > const &p1_types, std::vector< int > const &p2_types, double r_min, double r_max, int r_bins, bool log_flag, bool int_flag)
Calculate the distribution of particles around others.
Statistical tools to analyze simulations.
std::invoke_result_t< F, Particle const & > type
Decay a tuple of only 1 type to that type.
std::tuple< std::invoke_result_t< F, Particle const & >... > type
Struct holding all information for one particle.
Definition Particle.hpp:435