19#ifndef OBSERVABLES_BONDANGLES_HPP
20#define OBSERVABLES_BONDANGLES_HPP
25#include "system/System.hpp"
47 if (this->
ids().size() < 3)
48 throw std::runtime_error(
"At least 3 particles are required");
55 auto const positions_sorted = detail::get_all_particle_positions(
56 comm, local_particles,
ids(),
traits,
false);
58 if (comm.rank() != 0) {
64 auto v1 = box_geo.get_mi_vector(positions_sorted[1], positions_sorted[0]);
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(
77 res[i] = acos(-cosine);
83 std::vector<std::size_t>
shape()
const override {
84 assert(
ids().size() >= 2);
85 return {
ids().size() - 2};
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.
Particle-based observable.
PidObservable(std::vector< int > ids)
std::vector< int > const & ids() const
std::shared_ptr< BoxGeometry > box_geo
std::vector< std::reference_wrapper< Particle const > > ParticleReferenceRange