added flexibility
This commit is contained in:
64
src/motions/base.cpp
Normal file
64
src/motions/base.cpp
Normal file
@ -0,0 +1,64 @@
|
||||
|
||||
#include "base.h"
|
||||
#include "coordinates.h"
|
||||
#include "bimodalangle.h"
|
||||
#include "isosmallangle.h"
|
||||
#include "random.h"
|
||||
#include "tetrahedral.h"
|
||||
|
||||
#include <stdexcept>
|
||||
|
||||
Motion::Motion(std::string name, const double delta, const double eta, std::mt19937_64& rng) : m_name(std::move(name)), m_delta(delta), m_eta(eta), m_rng(rng) {
|
||||
m_uni_dist = std::uniform_real_distribution(0., 1.);
|
||||
}
|
||||
|
||||
Motion::Motion(std::string name, std::mt19937_64& rng) : m_name(std::move(name)), m_rng(rng) {
|
||||
m_uni_dist = std::uniform_real_distribution(0., 1.);
|
||||
}
|
||||
|
||||
double Motion::omega_q(const double cos_theta, const double phi) const {
|
||||
const double cos_theta_square = cos_theta * cos_theta;
|
||||
const double sin_theta_square = 1. - cos_theta_square;
|
||||
|
||||
return M_PI * m_delta * (3. * cos_theta_square - 1. - m_eta * sin_theta_square * std::cos(2.*phi));
|
||||
}
|
||||
|
||||
double Motion::omega_q(const SphericalPos& pos) const {
|
||||
auto [cos_theta, phi] = pos;
|
||||
|
||||
return omega_q(cos_theta, phi);
|
||||
}
|
||||
|
||||
SphericalPos Motion::draw_position() {
|
||||
const double cos_theta = 1 - 2 * m_uni_dist(m_rng);
|
||||
const double phi = 2.0 * M_PI * m_uni_dist(m_rng);
|
||||
|
||||
return {cos_theta, phi};
|
||||
}
|
||||
|
||||
Motion* Motion::createFromInput(const std::string& input, std::mt19937_64& rng) {
|
||||
if (input == "TetrahedralJump")
|
||||
return new TetrahedralJump(rng);
|
||||
|
||||
if (input == "IsotropicAngle")
|
||||
return new SmallAngle(rng);
|
||||
|
||||
if (input == "RandomJump")
|
||||
return new RandomJump(rng);
|
||||
|
||||
if (input == "BimodalAngle")
|
||||
return new BimodalAngle(rng);
|
||||
|
||||
throw std::invalid_argument("Invalid input " + input);
|
||||
}
|
||||
|
||||
void Motion::setParameters(const std::unordered_map<std::string, double> ¶meters) {
|
||||
m_delta = parameters.at("delta");
|
||||
m_eta = parameters.at("eta");
|
||||
}
|
||||
|
||||
|
||||
std::ostream& operator<<(std::ostream& os, const Motion& m) {
|
||||
os << m.getName();
|
||||
return os;
|
||||
}
|
43
src/motions/base.h
Normal file
43
src/motions/base.h
Normal file
@ -0,0 +1,43 @@
|
||||
#ifndef RWSIM_MOTIONBASE_H
|
||||
#define RWSIM_MOTIONBASE_H
|
||||
|
||||
#include "coordinates.h"
|
||||
|
||||
#include <random>
|
||||
#include <unordered_map>
|
||||
|
||||
class Motion {
|
||||
public:
|
||||
virtual ~Motion() = default;
|
||||
|
||||
Motion(std::string, double, double, std::mt19937_64&);
|
||||
explicit Motion(std::string, std::mt19937_64&);
|
||||
|
||||
SphericalPos draw_position();
|
||||
[[nodiscard]] double omega_q(double, double) const;
|
||||
[[nodiscard]] double omega_q(const SphericalPos&) const;
|
||||
|
||||
virtual void initialize() = 0;
|
||||
virtual double jump() = 0;
|
||||
|
||||
virtual void setParameters(const std::unordered_map<std::string, double>&);
|
||||
|
||||
[[nodiscard]] double getDelta() const { return m_delta; }
|
||||
void setDelta(const double delta) { m_delta = delta; }
|
||||
[[nodiscard]] double getEta() const { return m_eta; }
|
||||
void setEta(const double eta) { m_eta = eta; }
|
||||
[[nodiscard]] std::string getName() const { return m_name; }
|
||||
|
||||
static Motion* createFromInput(const std::string& input, std::mt19937_64& rng);
|
||||
|
||||
protected:
|
||||
std::string m_name{"BaseMotion"};
|
||||
double m_delta{1.};
|
||||
double m_eta{0.};
|
||||
std::mt19937_64& m_rng;
|
||||
std::uniform_real_distribution<> m_uni_dist;
|
||||
};
|
||||
|
||||
std::ostream& operator<<(std::ostream& os, const Motion& m);
|
||||
|
||||
#endif //RWSIM_MOTIONBASE_H
|
30
src/motions/bimodalangle.cpp
Normal file
30
src/motions/bimodalangle.cpp
Normal file
@ -0,0 +1,30 @@
|
||||
|
||||
#include "bimodalangle.h"
|
||||
#include "base.h"
|
||||
|
||||
BimodalAngle::BimodalAngle(const double delta, const double eta, const double angle1, const double angle2, const double prob, std::mt19937_64 &rng) :
|
||||
Motion(std::string("BimodalAngle"), delta, eta, rng),
|
||||
m_angle1(angle1 * M_PI / 180.0),
|
||||
m_angle2(angle2 * M_PI / 180.0),
|
||||
m_prob(prob) {};
|
||||
BimodalAngle::BimodalAngle(std::mt19937_64 &rng) : Motion(std::string("BimodalAngle"), rng) {}
|
||||
|
||||
void BimodalAngle::initialize() {
|
||||
m_prev_pos = draw_position();
|
||||
};
|
||||
|
||||
double BimodalAngle::jump() {
|
||||
const double angle = m_uni_dist(m_rng) < m_prob ? m_angle1 : m_angle2;
|
||||
const double gamma{2 * M_PI * m_uni_dist(m_rng)};
|
||||
m_prev_pos = rotate(m_prev_pos, angle, gamma);
|
||||
|
||||
return omega_q(m_prev_pos);
|
||||
}
|
||||
|
||||
void BimodalAngle::setParameters(const std::unordered_map<std::string, double> ¶meter) {
|
||||
Motion::setParameters(parameter);
|
||||
|
||||
m_angle1 = parameter.at("angle1") * M_PI / 180.;
|
||||
m_angle2 = parameter.at("angle2") * M_PI / 180.;
|
||||
m_prob = parameter.at("probability1");
|
||||
}
|
23
src/motions/bimodalangle.h
Normal file
23
src/motions/bimodalangle.h
Normal file
@ -0,0 +1,23 @@
|
||||
|
||||
#ifndef BIMODALANGLE_H
|
||||
#define BIMODALANGLE_H
|
||||
|
||||
#include "base.h"
|
||||
|
||||
class BimodalAngle : public Motion {
|
||||
public:
|
||||
BimodalAngle(double, double, double, double, double, std::mt19937_64& );
|
||||
explicit BimodalAngle(std::mt19937_64&);
|
||||
|
||||
void initialize() override;
|
||||
double jump() override;
|
||||
void setParameters(const std::unordered_map<std::string, double> &) override;
|
||||
|
||||
protected:
|
||||
double m_angle1{0};
|
||||
double m_angle2{0};
|
||||
double m_prob{0};
|
||||
SphericalPos m_prev_pos{0., 0.};
|
||||
};
|
||||
|
||||
#endif //BIMODALANGLE_H
|
39
src/motions/coordinates.cpp
Normal file
39
src/motions/coordinates.cpp
Normal file
@ -0,0 +1,39 @@
|
||||
#include "coordinates.h"
|
||||
|
||||
#include <cmath>
|
||||
#include <iostream>
|
||||
|
||||
SphericalPos rotate(const SphericalPos& old_pos, const double alpha, const double beta) {
|
||||
const double sin_alpha{std::sin(alpha)};
|
||||
const double cos_alpha{std::cos(alpha)};
|
||||
|
||||
const double sin_beta{std::sin(beta)};
|
||||
const double cos_beta{std::cos(beta)};
|
||||
|
||||
const double cos_theta{old_pos.cos_theta};
|
||||
|
||||
if (std::abs(cos_theta) == 1) {
|
||||
return xyz_to_spherical(CartesianPos{cos_alpha * cos_beta, cos_alpha * sin_beta, cos_alpha * cos_theta});
|
||||
}
|
||||
|
||||
const double norm{std::sqrt(1 - cos_theta * cos_theta) + 1e-15};
|
||||
|
||||
auto [x, y , z] = spherical_to_xyz(old_pos);
|
||||
|
||||
const auto new_pos = CartesianPos{
|
||||
cos_alpha * x + sin_alpha * (-x * z * sin_beta - y * cos_beta) / norm,
|
||||
cos_alpha * y + sin_alpha * (-y * z * sin_beta + x * cos_beta) / norm,
|
||||
cos_alpha * z + sin_alpha * norm * sin_beta
|
||||
};
|
||||
|
||||
return xyz_to_spherical(new_pos);
|
||||
}
|
||||
|
||||
CartesianPos spherical_to_xyz(const SphericalPos& pos) {
|
||||
const double sin_theta = std::sin(std::acos(pos.cos_theta));
|
||||
return {sin_theta * std::cos(pos.phi), sin_theta * std::sin(pos.phi), pos.cos_theta};
|
||||
}
|
||||
|
||||
SphericalPos xyz_to_spherical(const CartesianPos& pos) {
|
||||
return {pos.z, std::atan2(pos.y, pos.x)};
|
||||
}
|
20
src/motions/coordinates.h
Normal file
20
src/motions/coordinates.h
Normal file
@ -0,0 +1,20 @@
|
||||
|
||||
#ifndef COORDINATES_H
|
||||
#define COORDINATES_H
|
||||
|
||||
struct SphericalPos {
|
||||
double cos_theta;
|
||||
double phi;
|
||||
};
|
||||
|
||||
struct CartesianPos {
|
||||
double x;
|
||||
double y;
|
||||
double z;
|
||||
};
|
||||
|
||||
SphericalPos rotate(const SphericalPos&, double, double);
|
||||
CartesianPos spherical_to_xyz(const SphericalPos&);
|
||||
SphericalPos xyz_to_spherical(const CartesianPos&);
|
||||
|
||||
#endif //COORDINATES_H
|
26
src/motions/isosmallangle.cpp
Normal file
26
src/motions/isosmallangle.cpp
Normal file
@ -0,0 +1,26 @@
|
||||
#include "isosmallangle.h"
|
||||
#include "coordinates.h"
|
||||
|
||||
#include <iostream>
|
||||
|
||||
|
||||
|
||||
SmallAngle::SmallAngle(const double delta, const double eta, const double chi, std::mt19937_64 &rng) :
|
||||
Motion(std::string("IsotropicAngle"), delta, eta, rng), m_chi(chi * M_PI / 180.0) {};
|
||||
SmallAngle::SmallAngle(std::mt19937_64 &rng) : Motion(std::string("IsotropicAngle"), rng) {}
|
||||
|
||||
void SmallAngle::initialize() {
|
||||
m_prev_pos = draw_position();
|
||||
};
|
||||
|
||||
double SmallAngle::jump() {
|
||||
const double gamma{2 * M_PI * m_uni_dist(m_rng)};
|
||||
m_prev_pos = rotate(m_prev_pos, m_chi, gamma);
|
||||
|
||||
return omega_q(m_prev_pos);
|
||||
}
|
||||
|
||||
void SmallAngle::setParameters(const std::unordered_map<std::string, double> ¶meters) {
|
||||
m_chi = parameters.at("angle") * M_PI / 180.0;
|
||||
Motion::setParameters(parameters);
|
||||
}
|
22
src/motions/isosmallangle.h
Normal file
22
src/motions/isosmallangle.h
Normal file
@ -0,0 +1,22 @@
|
||||
|
||||
#ifndef RWSIM_MOTIONISOSMALLANGLE_H
|
||||
#define RWSIM_MOTIONISOSMALLANGLE_H
|
||||
|
||||
#include "base.h"
|
||||
#include "coordinates.h"
|
||||
|
||||
class SmallAngle final : public Motion {
|
||||
public:
|
||||
SmallAngle(double, double, double, std::mt19937_64& );
|
||||
explicit SmallAngle(std::mt19937_64&);
|
||||
|
||||
void initialize() override;
|
||||
double jump() override;
|
||||
void setParameters(const std::unordered_map<std::string, double> &) override;
|
||||
|
||||
private:
|
||||
double m_chi{0};
|
||||
SphericalPos m_prev_pos{0., 0.};
|
||||
};
|
||||
|
||||
#endif //RWSIM_MOTIONISOSMALLANGLE_H
|
13
src/motions/random.cpp
Normal file
13
src/motions/random.cpp
Normal file
@ -0,0 +1,13 @@
|
||||
|
||||
#include "random.h"
|
||||
|
||||
|
||||
RandomJump::RandomJump(const double delta, const double eta, std::mt19937_64 &rng) : Motion(std::string("RandomJump"), delta, eta, rng) {}
|
||||
|
||||
RandomJump::RandomJump(std::mt19937_64 &rng) : Motion(std::string("RandomJump"), rng) {}
|
||||
|
||||
void RandomJump::initialize() {}
|
||||
|
||||
double RandomJump::jump() {
|
||||
return omega_q(draw_position());
|
||||
}
|
17
src/motions/random.h
Normal file
17
src/motions/random.h
Normal file
@ -0,0 +1,17 @@
|
||||
#ifndef RWSIM_MOTIONRANDOMJUMP_H
|
||||
#define RWSIM_MOTIONRANDOMJUMP_H
|
||||
|
||||
|
||||
#include "base.h"
|
||||
#include <random>
|
||||
|
||||
class RandomJump final : public Motion {
|
||||
public:
|
||||
RandomJump(double, double, std::mt19937_64&);
|
||||
explicit RandomJump(std::mt19937_64&);
|
||||
|
||||
void initialize() override;
|
||||
double jump() override;
|
||||
};
|
||||
|
||||
#endif //RWSIM_MOTIONRANDOMJUMP_H
|
29
src/motions/tetrahedral.cpp
Normal file
29
src/motions/tetrahedral.cpp
Normal file
@ -0,0 +1,29 @@
|
||||
#include "tetrahedral.h"
|
||||
|
||||
#include <random>
|
||||
|
||||
#include "tetrahedral.h"
|
||||
|
||||
TetrahedralJump::TetrahedralJump(const double delta, const double eta, std::mt19937_64& rng) :
|
||||
Motion(std::string{"FourSiteTetrahedral"}, delta, eta, rng) {}
|
||||
|
||||
TetrahedralJump::TetrahedralJump(std::mt19937_64& rng) : Motion(std::string{"FourSiteTetrahedral"}, rng) {}
|
||||
|
||||
void TetrahedralJump::initialize() {
|
||||
const auto pos = draw_position();
|
||||
m_corners[0] = omega_q(pos);
|
||||
|
||||
const double alpha = 2. * M_PI * m_uni_dist(m_rng);
|
||||
|
||||
for (int i = 1; i<4; i++) {
|
||||
auto corner_pos = rotate(pos, m_beta, alpha + (i-1) * 2*M_PI/3.);
|
||||
m_corners[i] = omega_q(corner_pos);
|
||||
}
|
||||
}
|
||||
|
||||
double TetrahedralJump::jump() {
|
||||
m_corner_idx += m_chooser(m_rng);
|
||||
m_corner_idx %= 4;
|
||||
|
||||
return m_corners[m_corner_idx];
|
||||
}
|
28
src/motions/tetrahedral.h
Normal file
28
src/motions/tetrahedral.h
Normal file
@ -0,0 +1,28 @@
|
||||
#ifndef RWSIM_MOTIONTETRAHEDRAL_H
|
||||
#define RWSIM_MOTIONTETRAHEDRAL_H
|
||||
|
||||
#include "base.h"
|
||||
#include <random>
|
||||
#include <cmath>
|
||||
#include <array>
|
||||
|
||||
class TetrahedralJump final : public Motion {
|
||||
public:
|
||||
TetrahedralJump(double, double, std::mt19937_64&);
|
||||
explicit TetrahedralJump(std::mt19937_64&);
|
||||
|
||||
void initialize() override;
|
||||
double jump() override;
|
||||
|
||||
private:
|
||||
const double m_beta{std::acos(-1/3.)};
|
||||
|
||||
std::array<double, 4> m_corners{};
|
||||
int m_corner_idx{0};
|
||||
|
||||
std::uniform_int_distribution<> m_chooser{1, 3};
|
||||
|
||||
};
|
||||
|
||||
|
||||
#endif //RWSIM_MOTIONTETRAHEDRAL_H
|
Reference in New Issue
Block a user