ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
BondAngles.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_BONDANGLES_HPP
20#define OBSERVABLES_BONDANGLES_HPP
21
22#include "BoxGeometry.hpp"
23#include "PidObservable.hpp"
25#include "system/System.hpp"
26
27#include <utils/Vector.hpp>
28
29#include <algorithm>
30#include <cassert>
31#include <cmath>
32#include <cstddef>
33#include <stdexcept>
34#include <utility>
35#include <vector>
36
37namespace Observables {
38
39/** Calculate bond angles between particles in a polymer.
40 * For @f$ n @f$ bonded particles, return the @f$ n-2 @f$ angles along the
41 * chain, in radians.
42 */
43class BondAngles : public PidObservable {
44public:
46 explicit BondAngles(std::vector<int> ids) : PidObservable(std::move(ids)) {
47 if (this->ids().size() < 3)
48 throw std::runtime_error("At least 3 particles are required");
49 }
50
51 std::vector<double>
52 evaluate(boost::mpi::communicator const &comm,
53 ParticleReferenceRange const &local_particles,
54 const ParticleObservables::traits<Particle> &traits) const override {
55 auto const positions_sorted = detail::get_all_particle_positions(
56 comm, local_particles, ids(), traits, false);
57
58 if (comm.rank() != 0) {
59 return {};
60 }
61
62 auto const &box_geo = *System::get_system().box_geo;
63 std::vector<double> res(n_values());
64 auto v1 = box_geo.get_mi_vector(positions_sorted[1], positions_sorted[0]);
65 auto n1 = v1.norm();
66 for (std::size_t i = 0, end = n_values(); i < end; i++) {
67 auto v2 = box_geo.get_mi_vector(positions_sorted[i + 2],
68 positions_sorted[i + 1]);
69 auto const n2 = v2.norm();
70 auto const cosine = std::clamp(
71 (v1 * v2) / (n1 * n2), -tiny_cos_angle_value, +tiny_cos_angle_value);
72 /* to reduce computational time, after calculating an angle ijk, the
73 * vector r_ij takes the value r_jk, but to orient it correctly, it has
74 * to be multiplied -1; it's cheaper to do this operation on a double
75 * than on a vector of doubles
76 */
77 res[i] = acos(-cosine);
78 v1 = v2;
79 n1 = n2;
80 }
81 return res;
82 }
83 std::vector<std::size_t> shape() const override {
84 assert(ids().size() >= 2);
85 return {ids().size() - 2};
86 }
87};
88
89} // Namespace Observables
90
91#endif
Vector implementation and trait types for boost qvm interoperability.
Common code for functions calculating angle forces.
constexpr auto tiny_cos_angle_value
Tiny angle cutoff for cosine calculations.
Calculate bond angles between particles in a polymer.
std::vector< std::size_t > shape() const override
BondAngles(std::vector< int > ids)
std::vector< double > evaluate(boost::mpi::communicator const &comm, ParticleReferenceRange const &local_particles, const ParticleObservables::traits< Particle > &traits) const override
std::size_t n_values() const
Size of the flat array returned by the observable.
std::vector< int > const & ids() const
std::shared_ptr< BoxGeometry > box_geo
std::vector< std::reference_wrapper< Particle const > > ParticleReferenceRange
System & get_system()
STL namespace.