224 lines
No EOL
7.3 KiB
C++
224 lines
No EOL
7.3 KiB
C++
#include <algorithm>
|
|
#include <array>
|
|
#include <boost/numeric/odeint.hpp>
|
|
#include <Eigen/Geometry>
|
|
#include <fstream>
|
|
#include <gsl/gsl_spline.h>
|
|
#include <iostream>
|
|
#include <numeric>
|
|
#include <string>
|
|
#include <stdexcept>
|
|
#include <vector>
|
|
|
|
#include "loadtxt.h"
|
|
|
|
using double3 = Eigen::Vector3d;
|
|
|
|
// Our units are {kiloparsec, solar mass, gigayear}
|
|
constexpr double G = 4.498317481097514e-06;
|
|
|
|
class Interp {
|
|
// This is a wrapper around GSL spline interpolation. I tried to use
|
|
// boost::math::interpolators but as of version 1.72 none were suitable:
|
|
// barycentric_rational is the one suitalbe for non-uniform sampling but it
|
|
// is very slow. I also tried to resample the data uniformly using
|
|
// barycentric rational interpolation and then using cardinal_cubic_b_spline
|
|
// on the result, but was still slower than GSL.
|
|
public:
|
|
Interp(std::vector<double>& x, std::vector<double>& y)
|
|
{
|
|
acc = gsl_interp_accel_alloc();
|
|
spline = gsl_spline_alloc(gsl_interp_cspline, x.size());
|
|
gsl_spline_init(spline, x.data(), y.data(), x.size());
|
|
}
|
|
Interp() {}
|
|
inline double operator()(double x) const
|
|
{
|
|
return gsl_spline_eval(spline, x, acc);
|
|
}
|
|
|
|
private:
|
|
gsl_interp_accel *acc;
|
|
gsl_spline *spline;
|
|
};
|
|
|
|
double3 plummer(const double M, const double b, const double3& pos)
|
|
{
|
|
double r2 = pos.squaredNorm() + b*b;
|
|
double r = sqrt(r2);
|
|
double r3_inv = 1/(r*r2);
|
|
return -G*M*pos*r3_inv;
|
|
}
|
|
|
|
double3 nfw(const double rho_0, const double b, const double3& pos)
|
|
{
|
|
double r2 = pos.squaredNorm();
|
|
double r = sqrt(r2);
|
|
double tmp = -4*M_PI*G*rho_0*b*b*b*(log(1+r/b) - r/(b+r))/(r2*r);
|
|
return pos*tmp;
|
|
}
|
|
|
|
double3 miyamoto_nagai(const double M, const double a, const double b, const double phi, const double theta, const double3& pos)
|
|
{
|
|
// Construct new z-axis from the angles
|
|
double cos_theta = cos(theta);
|
|
double sin_theta = sqrt(1-cos_theta*cos_theta);
|
|
double cos_phi = cos(phi);
|
|
double sin_phi = sqrt(1-cos_phi*cos_phi);
|
|
Eigen::Vector3d old_z_axis = {cos_phi*sin_theta, sin_phi*sin_theta, cos_theta};
|
|
|
|
// Construct rotation
|
|
auto rot = Eigen::Quaternion<double>::FromTwoVectors(old_z_axis, Eigen::Vector3d::UnitZ());
|
|
|
|
// Rotate the position vector
|
|
auto new_pos = rot * pos;
|
|
|
|
// Calculate acceleration in new frame
|
|
Eigen::Vector3d acc_in_rotated_frame;
|
|
auto z_tmp = sqrt(new_pos[2]*new_pos[2] + b*b);
|
|
auto r2_tmp = new_pos[0]*new_pos[0] + new_pos[1]*new_pos[1] + (z_tmp + a)*(z_tmp + a);
|
|
auto r_tmp = sqrt(r2_tmp);
|
|
auto tmp = G*M / (r2_tmp*r_tmp);
|
|
acc_in_rotated_frame[0] = - tmp * new_pos[0];
|
|
acc_in_rotated_frame[1] = - tmp * new_pos[1];
|
|
acc_in_rotated_frame[2] = - tmp * new_pos[2] * (z_tmp + a)/z_tmp;
|
|
|
|
// Return to original frame
|
|
return rot.inverse() * acc_in_rotated_frame;
|
|
}
|
|
|
|
class Galaxy {
|
|
public:
|
|
Galaxy(std::string file_name)
|
|
{
|
|
auto data = Loadtxt(file_name, {1, 2, 3, 4, 5, 6, 7, 8}).get_cols();
|
|
auto& t_data = data[0];
|
|
std::transform(begin(t_data), end(t_data), begin(t_data), [t0=t_data[0]](const double& x){ return x-t0; });
|
|
for (int i=0; i<7; i++) interp[i] = Interp(t_data, data[i+1]);
|
|
}
|
|
void func(const std::array<double, 6> &y, std::array<double, 6> &f, const double t)
|
|
{
|
|
f[0] = y[3]; // vx -> x'
|
|
f[1] = y[4]; // vy -> y'
|
|
f[2] = y[5]; // vz -> z'
|
|
double phi = interp[0](t);
|
|
double theta = interp[1](t);
|
|
double m_disk = interp[2](t);
|
|
double a_disk = interp[3](t);
|
|
double b_disk = interp[4](t);
|
|
double rho_0_nfw = interp[5](t);
|
|
double b_nfw = interp[6](t);
|
|
|
|
double3 pos(y.data());
|
|
double3 acc_disk = miyamoto_nagai(m_disk, a_disk, b_disk, phi, theta, pos);
|
|
double3 acc_halo = nfw(rho_0_nfw, b_nfw, pos);
|
|
for (int i=0; i<3; i++) f[3+i] = acc_disk[i] + acc_halo[i];
|
|
}
|
|
|
|
|
|
private:
|
|
Interp interp[10];
|
|
} galaxy("subhalo_411321_parameters_processed.dat");
|
|
|
|
extern "C"
|
|
void integrate(const double y0[], const double t_max, const double stride_size, double y[])
|
|
{
|
|
using namespace boost::numeric::odeint;
|
|
using Coordinates = std::array<double, 6>;
|
|
auto stepper = bulirsch_stoer<Coordinates>(1E-7, 0);
|
|
auto function_wrapper = [](const Coordinates &x, Coordinates &dxdt, const double t) { return galaxy.func(x, dxdt, t); };
|
|
const int stride_count = t_max / stride_size;
|
|
Coordinates y_current;
|
|
std::copy(y0, y0+6, begin(y_current));
|
|
std::copy(y0, y0+6, y);
|
|
double t = 0;
|
|
const double h = std::min(1./4096., stride_size);
|
|
for (int i=0; i<stride_count; i++) {
|
|
integrate_adaptive(stepper, function_wrapper, y_current, t, t+stride_size, h);
|
|
// NOTE h here is just a recommended initial step size for the stepper,
|
|
// the actual step is adapted as needed. Since the integration is
|
|
// interrupted here in order for the data to be saved, the result
|
|
// somewhat depends on stride_size.
|
|
std::copy(begin(y_current), end(y_current), y+(i+1)*6);
|
|
t += stride_size;
|
|
}
|
|
}
|
|
|
|
|
|
#include <random>
|
|
#include <chrono>
|
|
#include <unistd.h>
|
|
int mainXXX()
|
|
{
|
|
std::cout << "Hi" << std::endl;
|
|
|
|
double x_max = 100;
|
|
double v_max = 200;
|
|
double t_max = 10.;
|
|
|
|
double3 target_pos = {-9.1817914423447837E+03, -2.3161972143758221E+03, 4.2321813890923086E+03}; // pc
|
|
double3 target_vel = { 1.6995805137819784E+02, 1.3476658021349482E+02, -1.3342192079702110E+02}; // km/s
|
|
|
|
// Unit adjustment
|
|
target_pos *= 0.001; // Now in kpc
|
|
target_vel *= 1.0226911647958985; // Now in kpc/Gyr
|
|
|
|
double tolrance = 2;
|
|
double max_dx = target_pos.norm()*tolrance;
|
|
double max_dv = target_vel.norm()*tolrance;
|
|
|
|
const auto current_time = std::chrono::system_clock::now().time_since_epoch().count();
|
|
const unsigned long big_prime = 840580612873;
|
|
const auto pid = getpid();
|
|
const auto seed = current_time + big_prime*pid; // TODO add thread number
|
|
|
|
std::default_random_engine generator(seed);
|
|
|
|
std::uniform_real_distribution<double> position_distribution(-x_max, x_max);
|
|
auto random_position_component = [&]() { return position_distribution(generator); };
|
|
|
|
std::uniform_real_distribution<double> velocity_distribution(-v_max, v_max);
|
|
auto random_velocity_component = [&]() { return velocity_distribution(generator); };
|
|
|
|
for (int i=0; i<100000; i++) {
|
|
if (i%1000==0) std::cout << i << "\n";
|
|
double y[12], y0[6];
|
|
for (int i=0; i<3; i++) y0[i] = random_position_component();
|
|
for (int i=3; i<6; i++) y0[i] = random_velocity_component();
|
|
|
|
integrate(y0, t_max, t_max, y);
|
|
double3 new_pos(y+6);
|
|
double3 new_vel(y+9);
|
|
|
|
|
|
if (((new_pos - target_pos).norm() < max_dx) && ((new_vel - target_vel).norm() < max_dv))
|
|
std::cout << "*";
|
|
}
|
|
|
|
|
|
|
|
|
|
// double y[12];
|
|
// double y0[] = {80,0,0,0,80,0};
|
|
// integrate(y0, 10, 10, y);
|
|
// for (int i=0; i<12; i++) std::cout << y[i] << std::endl;
|
|
|
|
|
|
|
|
// double y[12];
|
|
// double y0[] = {80,0,0,0,80,0};
|
|
// for (int i=0; i<8000; i++)
|
|
// integrate(y0, 10, 10, y);
|
|
|
|
/*double y0[] = {80,0,0,0,80,0};
|
|
double y[12];
|
|
integrate(y0, 10, 10, y);
|
|
for (int i=0; i<12; i++) std::cout << y[i] << std::endl;*/
|
|
|
|
//std::vector<double> w = {8, 0.12, 0.04, -0.08, 200, -0.001};
|
|
|
|
|
|
|
|
std::cout << "Bye" << std::endl;
|
|
return 0;
|
|
} |