#pragma once #include #include #include #include #include #include "pcg-cpp/include/pcg_random.hpp" #include "randutils/randutils.hpp" using Rng = randutils::random_generator; template using Vector = Eigen::Matrix; template using Matrix = Eigen::Matrix; template class Position : public Vector { private: void recenter() { for (unsigned i = 0; i < D; i++) { double& vi = Vector::operator()(i); signed n = floor(vi); vi -= n; if (i == D - 2) { Vector::operator()(D - 1) += γ * n; } } } public: const double& γ; Position(const double& γ) : γ(γ) {} Position(const Position& p) : Vector(p), γ(p.γ) {} Position(const Vector& v, const double& γ) : Vector(v), γ(γ) {} Position& operator=(const Position& p) { Vector::operator=(p); return *this; } Position& operator+=(const Vector& u) { Vector::operator+=(u); recenter(); return *this; } friend Position operator+(Position v, const Vector& u) { v += u; return v; } Position& operator-=(const Vector& u) { Vector::operator-=(u); recenter(); return *this; } friend Vector operator-(const Position& v, const Position& u) { std::vector> ds; ds.reserve(9); unsigned iCur = 0; unsigned iMin; double dMin = std::numeric_limits::infinity(); Position<2> v2(v.γ); v2(0) = v(D - 2); v2(1) = v(D - 1); for (signed i : {-1, 0, 1}) { for (signed j : {-1, 0, 1}) { Vector<2> u2; u2(0) = u(0) + i; u2(1) = u(1) + j + v.γ * i; ds.push_back(v2 - u2); double d = ds.back().squaredNorm(); if (d < dMin) { iMin = iCur; dMin = d; } iCur++; } } Vector w = v - (Vector)u; for (unsigned i = 0; i < D - 2; i++) { double& wi = w(i); if (wi > 0.5) { wi = 1 - wi; } } w(D - 2) = ds[iMin](0); w(D - 1) = ds[iMin](1); return (Vector)w; } }; template class Euclidean { public: Vector t; Matrix r; Euclidean() { t.setZero(); r.setIdentity(); } Euclidean(const Vector& t0, const Matrix& r0) : t(t0), r(r0) {} Euclidean inverse() const { return Euclidean(-r.transpose() * t, r.transpose()); } Vector act(const Vector& x) const { return r * x + t; } Position act(const Position& x) const { return Position(r * x, x.γ) + t; } Euclidean act(const Euclidean& R) const { return {act(R.t), r * R.r}; } }; template concept Particle = requires(T v, const Position& x, double r) { { v.x } -> std::same_as>; { v.r } -> std::same_as; { v.m } -> std::same_as; { v.U(x, r) } -> std::same_as; }; template T> class NeighborLists { private: unsigned n; const double& γ; std::vector> lists; unsigned mod(signed a, unsigned b) const { if (a >= 0) { return a % b; } else { return (a + b * ((b - a) / b)) % b; } } int iPow(int x, unsigned p) const { if (p == 0) return 1; if (p == 1) return x; int tmp = iPow(x, p / 2); if (p % 2 == 0) { return tmp * tmp; } else { return x * tmp * tmp; } } unsigned index(Position x) const { unsigned ind = 0; for (unsigned i = 0; i < D; i++) { ind += pow(n, i) * std::floor(n * x(i)); } assert(ind < lists.size()); return ind; } std::pair::iterator, bool> insert(T& p, unsigned ind) { return lists[ind].insert(&p); } size_t erase(T& p, unsigned ind) { return lists[ind].erase(&p); } void neighborsOfHelper(signed i0, unsigned d, std::set& ns, signed overEdge = 0) const { if (d == 0) { const std::set& a = lists[i0]; ns.insert(a.begin(), a.end()); } else { for (signed j : {-1, 0, 1}) { signed i1 = iPow(n, d) * (i0 / iPow(n, d)) + mod(i0 + (signed)round(overEdge * n * γ) + j * iPow(n, d-1), iPow(n, d)); signed newEdge = 0; if (d == 2) { if (i0 - i1 == n) newEdge = -1; if (i1 - i0 == n) { newEdge = 1; } } neighborsOfHelper(i1, d - 1, ns, newEdge); } } } public: NeighborLists(const double& γ, unsigned m) : n(m), γ(γ), lists(pow(n, D)) {} NeighborLists(const double& γ, double l) : NeighborLists(γ, (unsigned)floor(1 / l)) {} std::pair::iterator, bool> insert(T& p, const Position& x) { return insert(p, index(x)); } std::pair::iterator, bool> insert(T& p) { return insert(p, p.x); } size_t erase(T& p, const Position& x) { return erase(p, index(x)); } size_t erase(T& p) { return erase(p, p.x); } void move(T& p, const Position& xNew) { unsigned iOld = index(p.x); unsigned iNew = index(xNew); if (iNew != iOld) { erase(p, iOld); insert(p, iNew); } } std::set neighborsOf(const Position& x) const { std::set neighbors; neighborsOfHelper(index(x), D, neighbors); return neighbors; } }; template T> class Model { private: bool metropolisAccept(double β, double ΔE, Rng& r) const { return ΔE < 0 || exp(-β * ΔE) > r.uniform(0.0, 1.0); } bool swapAccept(const T& p1, const T& p2) const { double rat = p1.r / p2.r; return (rat < 1.2 && rat > 0.83); } public: double γ; NeighborLists nl; std::vector particles; Euclidean orientation; Model(unsigned N, double l, std::function g, Rng& r) : γ(0), nl(γ, l) { particles.reserve(N); for (unsigned i = 0; i < N; i++) { Position x(γ); for (double& xi : x) { xi = r.uniform(0.0, 1.0); } particles.push_back(T(x, g(r))); nl.insert(particles.back()); } }; void shear(double Δγ) { γ += Δγ; } void move(T& p, const Position& x) { nl.move(p, x); p.x = x; } bool move(double β, T& p, const Vector& Δx, Rng& r) { Position xNew = p.x + Δx; double ΔE = 0; for (const T* neighbor : nl.neighborsOf(p.x)) { if (neighbor != &p) { ΔE -= neighbor->U(p.x, p.r); } } for (const T* neighbor : nl.neighborsOf(xNew)) { if (neighbor != &p) { ΔE += neighbor->U(xNew, p.r); } } if (metropolisAccept(β, ΔE, r)) { move(p, xNew); return true; } else { return false; } } bool swap(double β, T& p1, T& p2, Rng& r) { if (!swapAccept(p1, p2)) return false; double ΔE = 0; for (const T* neighbor : nl.neighborsOf(p1.x)) { if (neighbor != &p1 && neighbor != &p2) { ΔE += neighbor->U(p1.x, p2.r) - neighbor->U(p1.x, p1.r); } } for (const T* neighbor : nl.neighborsOf(p2.x)) { if (neighbor != &p1 && neighbor != &p2) { ΔE += neighbor->U(p2.x, p1.r) - neighbor->U(p2.x, p2.r); } } if (metropolisAccept(β, ΔE, r)) { std::swap(p1.r, p2.r); return true; } else { return false; } } bool randomSwap(double β, Rng& r) { return swap(β, r.pick(particles), r.pick(particles), r); } bool randomMove(double β, double δ, Rng& r) { Vector Δx; for (double& Δxi : Δx) { Δxi = r.variate(0, δ); } return move(β, r.pick(particles), Δx, r); } unsigned clusterFlip(double β, const Euclidean& R, T& p0, Rng& r) { unsigned n = 0; std::queue q; q.push(&p0); while (!q.empty()) { T* p = q.front(); q.pop(); if (!p->m) { p->m = true; Position xNew = R.act(p->x); for (T* neighbor : nl.neighborsOf(xNew)) { if (!neighbor->m) { double ΔE = neighbor->U(xNew, p->r) - neighbor->U(p->x, p->r); if (1 - exp(-β * ΔE) > r.uniform(0.0, 1.0)) { q.push(neighbor); } } } move(*p, xNew); n++; } } for (T& p : particles) { p.m = false; } if (n > particles.size() / 2) { orientation = R.act(orientation); } return n; } unsigned clusterSwap(double β, T& s1, T& s2, Rng& r) { if (!swapAccept(s1, s2)) return 0; Vector<2> t1 = s1.x; Vector<2> t2 = s2.x; Vector<2> t = (t1 + t2) / 2; Matrix<2> mat; mat(0, 0) = -1; mat(1, 1) = -1; mat(0, 1) = 0; mat(1, 0) = 0; Euclidean<2> g(t - mat * t, mat); return clusterFlip(β, g, s1, r); } unsigned randomClusterSwap(double β, Rng& r) { return clusterSwap(β, r.pick(particles), r.pick(particles), r); } }; template class Sphere { public: Position x; bool m; double r; Sphere() : m(false) {}; Sphere(const Position& x, double r) : x(x), m(false), r(r) {}; double regularizedDistanceSquared(const Position& y, double ry) const { return (x - y).squaredNorm() / pow(r + ry, 2); } bool intersectsBoundary() const { for (double xi : x) { if (xi < r || 1 - xi < r) { return true; } } return false; } }; std::function gContinuousPolydisperse(double Rmax) { return [Rmax] (Rng& r) -> double { return pow(r.uniform(pow(2.219, -2), 1.0), -0.5) * Rmax / 2.219; }; }