]>
Commit | Line | Data |
---|---|---|
7c673cae FG |
1 | /* |
2 | [auto_generated] | |
3 | boost/numeric/odeint/integrate/integrate_const.hpp | |
4 | ||
5 | [begin_description] | |
6 | Constant integration of ODEs, meaning that the state of the ODE is observed on constant time intervals. | |
7 | The routines makes full use of adaptive and dense-output methods. | |
8 | [end_description] | |
9 | ||
10 | Copyright 2011-2013 Karsten Ahnert | |
11 | Copyright 2011-2015 Mario Mulansky | |
12 | ||
13 | Distributed under the Boost Software License, Version 1.0. | |
14 | (See accompanying file LICENSE_1_0.txt or | |
15 | copy at http://www.boost.org/LICENSE_1_0.txt) | |
16 | */ | |
17 | ||
18 | ||
19 | #ifndef BOOST_NUMERIC_ODEINT_INTEGRATE_INTEGRATE_CONST_HPP_INCLUDED | |
20 | #define BOOST_NUMERIC_ODEINT_INTEGRATE_INTEGRATE_CONST_HPP_INCLUDED | |
21 | ||
22 | #include <boost/type_traits/is_same.hpp> | |
23 | ||
24 | #include <boost/numeric/odeint/stepper/stepper_categories.hpp> | |
25 | #include <boost/numeric/odeint/integrate/null_observer.hpp> | |
26 | #include <boost/numeric/odeint/integrate/check_adapter.hpp> | |
27 | #include <boost/numeric/odeint/integrate/detail/integrate_const.hpp> | |
28 | #include <boost/numeric/odeint/integrate/detail/integrate_adaptive.hpp> | |
29 | ||
30 | namespace boost { | |
31 | namespace numeric { | |
32 | namespace odeint { | |
33 | ||
34 | ||
35 | /* | |
36 | * Integrates with constant time step dt. | |
37 | */ | |
38 | template<class Stepper, class System, class State, class Time, class Observer, class StepOverflowChecker> | |
39 | size_t integrate_const( | |
40 | Stepper stepper, System system, State &start_state, | |
41 | Time start_time, Time end_time, Time dt, | |
42 | Observer observer, StepOverflowChecker checker | |
43 | ) { | |
44 | typedef typename odeint::unwrap_reference<Stepper>::type::stepper_category stepper_category; | |
45 | // we want to get as fast as possible to the end | |
46 | // no overflow checks needed | |
47 | if (boost::is_same<null_observer, Observer>::value) { | |
48 | return detail::integrate_adaptive( | |
49 | stepper, system, start_state, | |
50 | start_time, end_time, dt, | |
51 | observer, stepper_category()); | |
52 | } | |
53 | else { | |
54 | // unwrap references | |
55 | typedef typename odeint::unwrap_reference< Stepper >::type stepper_type; | |
56 | typedef typename odeint::unwrap_reference< Observer >::type observer_type; | |
57 | typedef typename odeint::unwrap_reference< StepOverflowChecker >::type checker_type; | |
58 | ||
59 | return detail::integrate_const(checked_stepper<stepper_type, checker_type>(stepper, checker), | |
60 | system, start_state, | |
61 | start_time, end_time, dt, | |
62 | checked_observer<observer_type, checker_type>(observer, checker), | |
63 | stepper_category()); | |
64 | } | |
65 | } | |
66 | ||
67 | /** | |
68 | * \brief Second version to solve the forwarding problem, | |
69 | * can be called with Boost.Range as start_state. | |
70 | */ | |
71 | template<class Stepper, class System, class State, class Time, class Observer, class StepOverflowChecker > | |
72 | size_t integrate_const( | |
73 | Stepper stepper, System system, const State &start_state, | |
74 | Time start_time, Time end_time, Time dt, | |
75 | Observer observer, StepOverflowChecker checker | |
76 | ) { | |
77 | typedef typename odeint::unwrap_reference<Stepper>::type::stepper_category stepper_category; | |
78 | // we want to get as fast as possible to the end | |
79 | ||
80 | if (boost::is_same<null_observer, Observer>::value) { | |
81 | return detail::integrate_adaptive( | |
82 | stepper, system, start_state, | |
83 | start_time, end_time, dt, | |
84 | observer, stepper_category()); | |
85 | } | |
86 | else { | |
87 | typedef typename odeint::unwrap_reference< Stepper >::type stepper_type; | |
88 | typedef typename odeint::unwrap_reference< Observer >::type observer_type; | |
89 | typedef typename odeint::unwrap_reference< StepOverflowChecker >::type checker_type; | |
90 | ||
91 | return detail::integrate_const(checked_stepper<stepper_type, checker_type>(stepper, checker), | |
92 | system, start_state, | |
93 | start_time, end_time, dt, | |
94 | checked_observer<observer_type, checker_type>(observer, checker), | |
95 | stepper_category()); | |
96 | } | |
97 | } | |
98 | ||
99 | ||
100 | /** | |
101 | * \brief integrate_const without step overflow checker | |
102 | */ | |
103 | template<class Stepper, class System, class State, class Time, class Observer> | |
104 | size_t integrate_const( | |
105 | Stepper stepper, System system, State &start_state, | |
106 | Time start_time, Time end_time, Time dt, Observer observer) | |
107 | { | |
108 | typedef typename odeint::unwrap_reference<Stepper>::type::stepper_category stepper_category; | |
109 | return detail::integrate_const(stepper, system, start_state, | |
110 | start_time, end_time, dt, observer, stepper_category()); | |
111 | } | |
112 | ||
113 | /** | |
114 | * \brief Second version to solve the forwarding problem, | |
115 | * can be called with Boost.Range as start_state. | |
116 | */ | |
117 | template<class Stepper, class System, class State, class Time, class Observer> | |
118 | size_t integrate_const( | |
119 | Stepper stepper, System system, const State &start_state, | |
120 | Time start_time, Time end_time, Time dt, Observer observer | |
121 | ) { | |
122 | typedef typename odeint::unwrap_reference<Stepper>::type::stepper_category stepper_category; | |
123 | return detail::integrate_const(stepper, system, start_state, | |
124 | start_time, end_time, dt, observer, stepper_category()); | |
125 | } | |
126 | ||
127 | ||
128 | /** | |
129 | * \brief integrate_const without observer calls | |
130 | */ | |
131 | template<class Stepper, class System, class State, class Time> | |
132 | size_t integrate_const( | |
133 | Stepper stepper, System system, State &start_state, | |
134 | Time start_time, Time end_time, Time dt | |
135 | ) { | |
136 | return integrate_const(stepper, system, start_state, start_time, end_time, dt, null_observer()); | |
137 | } | |
138 | ||
139 | /** | |
140 | * \brief Second version to solve the forwarding problem, | |
141 | * can be called with Boost.Range as start_state. | |
142 | */ | |
143 | template<class Stepper, class System, class State, class Time> | |
144 | size_t integrate_const( | |
145 | Stepper stepper, System system, const State &start_state, | |
146 | Time start_time, Time end_time, Time dt | |
147 | ) { | |
148 | return integrate_const(stepper, system, start_state, start_time, end_time, dt, null_observer()); | |
149 | } | |
150 | ||
151 | ||
152 | ||
153 | ||
154 | ||
155 | ||
156 | /********* DOXYGEN *********/ | |
157 | /** | |
158 | * \fn integrate_const( Stepper stepper , System system , State &start_state , Time start_time , | |
159 | * Time end_time , Time dt , Observer observer , StepOverflowChecker checker ) | |
160 | * \brief Integrates the ODE with constant step size. | |
161 | * | |
162 | * Integrates the ODE defined by system using the given stepper. | |
163 | * This method ensures that the observer is called at constant intervals dt. | |
164 | * If the Stepper is a normal stepper without step size control, dt is also | |
165 | * used for the numerical scheme. If a ControlledStepper is provided, the | |
166 | * algorithm might reduce the step size to meet the error bounds, but it is | |
167 | * ensured that the observer is always called at equidistant time points | |
168 | * t0 + n*dt. If a DenseOutputStepper is used, the step size also may vary | |
169 | * and the dense output is used to call the observer at equidistant time | |
170 | * points. | |
171 | * If a max_step_checker is provided as StepOverflowChecker, a | |
172 | * no_progress_error is thrown if too many steps (default: 500) are performed | |
173 | * without progress, i.e. in between observer calls. If no checker is provided, | |
174 | * no such overflow check is performed. | |
175 | * | |
176 | * \param stepper The stepper to be used for numerical integration. | |
177 | * \param system Function/Functor defining the rhs of the ODE. | |
178 | * \param start_state The initial condition x0. | |
179 | * \param start_time The initial time t0. | |
180 | * \param end_time The final integration time tend. | |
181 | * \param dt The time step between observer calls, _not_ necessarily the | |
182 | * time step of the integration. | |
183 | * \param observer [optional] Function/Functor called at equidistant time intervals. | |
184 | * \param checker [optional] Functor to check for step count overflows, if no | |
185 | * checker is provided, no exception is thrown. | |
186 | * \return The number of steps performed. | |
187 | */ | |
188 | ||
189 | } // namespace odeint | |
190 | } // namespace numeric | |
191 | } // namespace boost | |
192 | ||
193 | ||
194 | ||
195 | #endif // BOOST_NUMERIC_ODEINT_INTEGRATE_INTEGRATE_CONST_HPP_INCLUDED |