find_crossing.cpp 4.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134
  1. /*
  2. * find_crossing.cpp
  3. *
  4. * Finds the energy threshold crossing for a damped oscillator.
  5. * The algorithm uses a dense out stepper with find_if to first find an
  6. * interval containing the threshold crossing and the utilizes the dense out
  7. * functionality with a bisection to further refine the interval until some
  8. * desired precision is reached.
  9. *
  10. * Copyright 2015 Mario Mulansky
  11. *
  12. * Distributed under the Boost Software License, Version 1.0.
  13. * (See accompanying file LICENSE_1_0.txt or
  14. * copy at http://www.boost.org/LICENSE_1_0.txt)
  15. */
  16. #include <iostream>
  17. #include <utility>
  18. #include <algorithm>
  19. #include <array>
  20. #include <boost/numeric/odeint/stepper/runge_kutta_dopri5.hpp>
  21. #include <boost/numeric/odeint/stepper/generation.hpp>
  22. #include <boost/numeric/odeint/iterator/adaptive_iterator.hpp>
  23. namespace odeint = boost::numeric::odeint;
  24. typedef std::array<double, 2> state_type;
  25. const double gam = 1.0; // damping strength
  26. void damped_osc(const state_type &x, state_type &dxdt, const double /*t*/)
  27. {
  28. dxdt[0] = x[1];
  29. dxdt[1] = -x[0] - gam * x[1];
  30. }
  31. struct energy_condition {
  32. // defines the threshold crossing in terms of a boolean functor
  33. double m_min_energy;
  34. energy_condition(const double min_energy)
  35. : m_min_energy(min_energy) { }
  36. double energy(const state_type &x) {
  37. return 0.5 * x[1] * x[1] + 0.5 * x[0] * x[0];
  38. }
  39. bool operator()(const state_type &x) {
  40. // becomes true if the energy becomes smaller than the threshold
  41. return energy(x) <= m_min_energy;
  42. }
  43. };
  44. template<class System, class Condition>
  45. std::pair<double, state_type>
  46. find_condition(state_type &x0, System sys, Condition cond,
  47. const double t_start, const double t_end, const double dt,
  48. const double precision = 1E-6) {
  49. // integrates an ODE until some threshold is crossed
  50. // returns time and state at the point of the threshold crossing
  51. // if no threshold crossing is found, some time > t_end is returned
  52. auto stepper = odeint::make_dense_output(1.0e-6, 1.0e-6,
  53. odeint::runge_kutta_dopri5<state_type>());
  54. auto ode_range = odeint::make_adaptive_range(std::ref(stepper), sys, x0,
  55. t_start, t_end, dt);
  56. // find the step where the condition changes
  57. auto found_iter = std::find_if(ode_range.first, ode_range.second, cond);
  58. if(found_iter == ode_range.second)
  59. {
  60. // no threshold crossing -> return time after t_end and ic
  61. return std::make_pair(t_end + dt, x0);
  62. }
  63. // the dense out stepper now covers the interval where the condition changes
  64. // improve the solution by bisection
  65. double t0 = stepper.previous_time();
  66. double t1 = stepper.current_time();
  67. double t_m;
  68. state_type x_m;
  69. // use odeint's resizing functionality to allocate memory for x_m
  70. odeint::adjust_size_by_resizeability(x_m, x0,
  71. typename odeint::is_resizeable<state_type>::type());
  72. while(std::abs(t1 - t0) > precision) {
  73. t_m = 0.5 * (t0 + t1); // get the mid point time
  74. stepper.calc_state(t_m, x_m); // obtain the corresponding state
  75. if (cond(x_m))
  76. t1 = t_m; // condition changer lies before midpoint
  77. else
  78. t0 = t_m; // condition changer lies after midpoint
  79. }
  80. // we found the interval of size eps, take it's midpoint as final guess
  81. t_m = 0.5 * (t0 + t1);
  82. stepper.calc_state(t_m, x_m);
  83. return std::make_pair(t_m, x_m);
  84. }
  85. int main(int argc, char **argv)
  86. {
  87. state_type x0 = {{10.0, 0.0}};
  88. const double t_start = 0.0;
  89. const double t_end = 10.0;
  90. const double dt = 0.1;
  91. const double threshold = 0.1;
  92. energy_condition cond(threshold);
  93. state_type x_cond;
  94. double t_cond;
  95. std::tie(t_cond, x_cond) = find_condition(x0, damped_osc, cond,
  96. t_start, t_end, dt, 1E-6);
  97. if(t_cond > t_end)
  98. {
  99. // time after t_end -> no threshold crossing within [t_start, t_end]
  100. std::cout << "No threshold crossing found." << std::endl;
  101. } else
  102. {
  103. std::cout.precision(16);
  104. std::cout << "Time of energy threshold crossing: " << t_cond << std::endl;
  105. std::cout << "State: [" << x_cond[0] << " , " << x_cond[1] << "]" << std::endl;
  106. std::cout << "Energy: " << cond.energy(x_cond) << std::endl;
  107. }
  108. }