black_hole.cpp 4.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151
  1. /*
  2. [auto_generated]
  3. libs/numeric/odeint/examples/black_hole.cpp
  4. [begin_description]
  5. This example shows how the __float128 from gcc libquadmath can be used with odeint.
  6. [end_description]
  7. Copyright 2012 Karsten Ahnert
  8. Copyright 2012 Lee Hodgkinson
  9. Copyright 2012 Mario Mulansky
  10. Distributed under the Boost Software License, Version 1.0.
  11. (See accompanying file LICENSE_1_0.txt or
  12. copy at http://www.boost.org/LICENSE_1_0.txt)
  13. */
  14. #include <cstdlib>
  15. #include <cmath>
  16. #include <iostream>
  17. #include <iterator>
  18. #include <utility>
  19. #include <algorithm>
  20. #include <cassert>
  21. #include <vector>
  22. #include <complex>
  23. extern "C" {
  24. #include <quadmath.h>
  25. }
  26. const __float128 zero =strtoflt128 ("0.0", NULL);
  27. namespace std {
  28. inline __float128 abs( __float128 x )
  29. {
  30. return fabsq( x );
  31. }
  32. inline __float128 sqrt( __float128 x )
  33. {
  34. return sqrtq( x );
  35. }
  36. inline __float128 pow( __float128 x , __float128 y )
  37. {
  38. return powq( x , y );
  39. }
  40. inline __float128 abs( std::complex< __float128 > x )
  41. {
  42. return sqrtq( x.real() * x.real() + x.imag() * x.imag() );
  43. }
  44. inline std::complex< __float128 > pow( std::complex< __float128> x , __float128 y )
  45. {
  46. __float128 r = pow( abs(x) , y );
  47. __float128 phi = atanq( x.imag() / x.real() );
  48. return std::complex< __float128 >( r * cosq( y * phi ) , r * sinq( y * phi ) );
  49. }
  50. }
  51. inline std::ostream& operator<< (std::ostream& os, const __float128& f) {
  52. char* y = new char[1000];
  53. quadmath_snprintf(y, 1000, "%.30Qg", f) ;
  54. os.precision(30);
  55. os<<y;
  56. delete[] y;
  57. return os;
  58. }
  59. #include <boost/array.hpp>
  60. #include <boost/range/algorithm.hpp>
  61. #include <boost/range/adaptor/filtered.hpp>
  62. #include <boost/range/numeric.hpp>
  63. #include <boost/numeric/odeint.hpp>
  64. using namespace boost::numeric::odeint;
  65. using namespace std;
  66. typedef __float128 my_float;
  67. typedef std::vector< std::complex < my_float > > state_type;
  68. struct radMod
  69. {
  70. my_float m_om;
  71. my_float m_l;
  72. radMod( my_float om , my_float l )
  73. : m_om( om ) , m_l( l ) { }
  74. void operator()( const state_type &x , state_type &dxdt , my_float r ) const
  75. {
  76. dxdt[0] = x[1];
  77. dxdt[1] = -(2*(r-1)/(r*(r-2)))*x[1]-((m_om*m_om*r*r/((r-2)*(r-2)))-(m_l*(m_l+1)/(r*(r-2))))*x[0];
  78. }
  79. };
  80. int main( int argc , char **argv )
  81. {
  82. state_type x(2);
  83. my_float re0 = strtoflt128( "-0.00008944230755601224204687038354994353820468" , NULL );
  84. my_float im0 = strtoflt128( "0.00004472229441850588228136889483397204368247" , NULL );
  85. my_float re1 = strtoflt128( "-4.464175354293244250869336196695966076150E-6 " , NULL );
  86. my_float im1 = strtoflt128( "-8.950483248390306670770345406051469584488E-6" , NULL );
  87. x[0] = complex< my_float >( re0 ,im0 );
  88. x[1] = complex< my_float >( re1 ,im1 );
  89. const my_float dt =strtoflt128 ("-0.001", NULL);
  90. const my_float start =strtoflt128 ("10000.0", NULL);
  91. const my_float end =strtoflt128 ("9990.0", NULL);
  92. const my_float omega =strtoflt128 ("2.0", NULL);
  93. const my_float ell =strtoflt128 ("1.0", NULL);
  94. my_float abs_err = strtoflt128( "1.0E-15" , NULL ) , rel_err = strtoflt128( "1.0E-10" , NULL );
  95. my_float a_x = strtoflt128( "1.0" , NULL ) , a_dxdt = strtoflt128( "1.0" , NULL );
  96. typedef runge_kutta_dopri5< state_type, my_float > dopri5_type;
  97. typedef controlled_runge_kutta< dopri5_type > controlled_dopri5_type;
  98. typedef dense_output_runge_kutta< controlled_dopri5_type > dense_output_dopri5_type;
  99. dense_output_dopri5_type dopri5( controlled_dopri5_type( default_error_checker< my_float >( abs_err , rel_err , a_x , a_dxdt ) ) );
  100. std::for_each( make_adaptive_time_iterator_begin(dopri5 , radMod(omega , ell) , x , start , end , dt) ,
  101. make_adaptive_time_iterator_end(dopri5 , radMod(omega , ell) , x ) ,
  102. []( const std::pair< state_type&, my_float > &x ) {
  103. std::cout << x.second << ", " << x.first[0].real() << "\n"; }
  104. );
  105. return 0;
  106. }