#include #include #include #include #include #include #include #include #include #include #include #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& x, std::vector& 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; }; void plummer(double M, double b, const double3& pos, double3& acc) { double r2 = (pos[0]*pos[0] + pos[1]*pos[1] + pos[2]*pos[2] + b*b); double r = sqrt(r2); double r3_inv = 1/(r*r2); acc = -G*M*pos*r3_inv; } void miyamoto_nagai(double M, double a, double b, double phi, double theta, const double3& pos, double3& acc) { // 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::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 acc = 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, 9, 10, 11}).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<10; i++) interp[i] = Interp(t_data, data[i+1]); } void func(const std::array &y, std::array &f, const double t) { f[0] = y[3]; // vx -> x' f[1] = y[4]; // vy -> y' f[2] = y[5]; // vz -> z' double m_disk = interp[0](t); double phi = interp[1](t); double theta = interp[2](t); double a_disk = interp[3](t); double b_disk = interp[4](t); double m_halo = interp[5](t); double b_halo = interp[6](t); double3 dalembert; // for (int i=0; i<3; i++) dalembert[i] = interp[7+i](t); dalembert = {0,0,0}; double3 acc_disk, acc_halo; const double3 pos(y.data()); miyamoto_nagai(m_disk, a_disk, b_disk, phi, theta, pos, acc_disk); plummer(m_halo, b_halo, pos, acc_halo); for (int i=0; i<3; i++) f[3+i] = acc_disk[i] + acc_halo[i] + dalembert[i]; } private: Interp interp[10]; } galaxy("subhalo_411321_parameters_processed.dat"); extern "C" int integrate(const double y0[], const double t_max, const double stride_size, double y[]) { using namespace boost::numeric::odeint; using Coordinates = std::array; auto stepper = bulirsch_stoer(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 = 1./4096.; for (int i=0; i w = {8, 0.12, 0.04, -0.08, 200, -0.001}; std::cout << "Bye" << std::endl; return 0; }