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
ParticleDistances.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2019-2022 The ESPResSo project
3 *
4 * This file is part of ESPResSo.
5 *
6 * ESPResSo is free software: you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation, either version 3 of the License, or
9 * (at your option) any later version.
10 *
11 * ESPResSo is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program. If not, see <http://www.gnu.org/licenses/>.
18 */
19#ifndef OBSERVABLES_PARTICLEDISTANCES_HPP
20#define OBSERVABLES_PARTICLEDISTANCES_HPP
21
22#include "BoxGeometry.hpp"
23#include "PidObservable.hpp"
24#include "system/System.hpp"
25
26#include <cassert>
27#include <cstddef>
28#include <stdexcept>
29#include <utility>
30#include <vector>
31
32namespace Observables {
33
34/** Calculate distances between particles in a polymer.
35 * For @f$ n @f$ bonded particles, return the @f$ n-1 @f$ distances separating
36 * them.
37 */
39public:
41 explicit ParticleDistances(std::vector<int> ids)
42 : PidObservable(std::move(ids)) {
43 if (this->ids().size() < 2)
44 throw std::runtime_error("At least 2 particles are required");
45 }
46
47 std::vector<double>
48 evaluate(boost::mpi::communicator const &comm,
49 ParticleReferenceRange const &local_particles,
50 const ParticleObservables::traits<Particle> &traits) const override {
51 auto const positions_sorted = detail::get_all_particle_positions(
52 comm, local_particles, ids(), traits, false);
53
54 if (comm.rank() != 0) {
55 return {};
56 }
57
58 auto const &box_geo = *System::get_system().box_geo;
59 std::vector<double> res(n_values());
60
61 for (std::size_t i = 0, end = n_values(); i < end; i++) {
62 auto const v =
63 box_geo.get_mi_vector(positions_sorted[i], positions_sorted[i + 1]);
64 res[i] = v.norm();
65 }
66 return res;
67 }
68 std::vector<std::size_t> shape() const override {
69 assert(!ids().empty());
70 return {ids().size() - 1};
71 }
72};
73
74} // namespace Observables
75
76#endif
std::size_t n_values() const
Size of the flat array returned by the observable.
Calculate distances between particles in a polymer.
std::vector< std::size_t > shape() const override
std::vector< double > evaluate(boost::mpi::communicator const &comm, ParticleReferenceRange const &local_particles, const ParticleObservables::traits< Particle > &traits) const override
ParticleDistances(std::vector< int > ids)
std::vector< int > const & ids() const
std::shared_ptr< BoxGeometry > box_geo
std::vector< std::reference_wrapper< Particle const > > ParticleReferenceRange
System & get_system()