Chris@16
|
1 /*
|
Chris@16
|
2 [auto_generated]
|
Chris@16
|
3 boost/numeric/odeint/stepper/adams_moulton.hpp
|
Chris@16
|
4
|
Chris@16
|
5 [begin_description]
|
Chris@16
|
6 Implementation of the Adams-Moulton method. This is method is not a real stepper, it is more a helper class
|
Chris@16
|
7 which computes the corrector step in the Adams-Bashforth-Moulton method.
|
Chris@16
|
8 [end_description]
|
Chris@16
|
9
|
Chris@101
|
10 Copyright 2011-2012 Karsten Ahnert
|
Chris@101
|
11 Copyright 2011-2013 Mario Mulansky
|
Chris@101
|
12 Copyright 2012 Christoph Koke
|
Chris@16
|
13
|
Chris@16
|
14 Distributed under the Boost Software License, Version 1.0.
|
Chris@16
|
15 (See accompanying file LICENSE_1_0.txt or
|
Chris@16
|
16 copy at http://www.boost.org/LICENSE_1_0.txt)
|
Chris@16
|
17 */
|
Chris@16
|
18
|
Chris@16
|
19
|
Chris@16
|
20 #ifndef BOOST_NUMERIC_ODEINT_STEPPER_ADAMS_MOULTON_HPP_INCLUDED
|
Chris@16
|
21 #define BOOST_NUMERIC_ODEINT_STEPPER_ADAMS_MOULTON_HPP_INCLUDED
|
Chris@16
|
22
|
Chris@16
|
23
|
Chris@16
|
24 #include <boost/numeric/odeint/util/bind.hpp>
|
Chris@16
|
25
|
Chris@16
|
26 #include <boost/numeric/odeint/algebra/range_algebra.hpp>
|
Chris@16
|
27 #include <boost/numeric/odeint/algebra/default_operations.hpp>
|
Chris@101
|
28 #include <boost/numeric/odeint/algebra/algebra_dispatcher.hpp>
|
Chris@101
|
29 #include <boost/numeric/odeint/algebra/operations_dispatcher.hpp>
|
Chris@16
|
30
|
Chris@16
|
31 #include <boost/numeric/odeint/util/state_wrapper.hpp>
|
Chris@16
|
32 #include <boost/numeric/odeint/util/is_resizeable.hpp>
|
Chris@16
|
33 #include <boost/numeric/odeint/util/resizer.hpp>
|
Chris@16
|
34
|
Chris@16
|
35 #include <boost/numeric/odeint/stepper/stepper_categories.hpp>
|
Chris@16
|
36 #include <boost/numeric/odeint/stepper/runge_kutta4_classic.hpp>
|
Chris@16
|
37
|
Chris@16
|
38 #include <boost/numeric/odeint/stepper/detail/adams_moulton_call_algebra.hpp>
|
Chris@16
|
39 #include <boost/numeric/odeint/stepper/detail/adams_moulton_coefficients.hpp>
|
Chris@16
|
40 #include <boost/numeric/odeint/stepper/detail/rotating_buffer.hpp>
|
Chris@16
|
41
|
Chris@16
|
42
|
Chris@16
|
43
|
Chris@16
|
44
|
Chris@16
|
45 namespace boost {
|
Chris@16
|
46 namespace numeric {
|
Chris@16
|
47 namespace odeint {
|
Chris@16
|
48
|
Chris@16
|
49
|
Chris@16
|
50 /*
|
Chris@16
|
51 * Static implicit Adams-Moulton multistep-solver without step size control and without dense output.
|
Chris@16
|
52 */
|
Chris@16
|
53 template<
|
Chris@16
|
54 size_t Steps ,
|
Chris@16
|
55 class State ,
|
Chris@16
|
56 class Value = double ,
|
Chris@16
|
57 class Deriv = State ,
|
Chris@16
|
58 class Time = Value ,
|
Chris@101
|
59 class Algebra = typename algebra_dispatcher< State >::algebra_type ,
|
Chris@101
|
60 class Operations = typename operations_dispatcher< State >::operations_type ,
|
Chris@16
|
61 class Resizer = initially_resizer
|
Chris@16
|
62 >
|
Chris@16
|
63 class adams_moulton
|
Chris@16
|
64 {
|
Chris@16
|
65 private:
|
Chris@16
|
66
|
Chris@16
|
67
|
Chris@16
|
68 public :
|
Chris@16
|
69
|
Chris@16
|
70 typedef State state_type;
|
Chris@16
|
71 typedef state_wrapper< state_type > wrapped_state_type;
|
Chris@16
|
72 typedef Value value_type;
|
Chris@16
|
73 typedef Deriv deriv_type;
|
Chris@16
|
74 typedef state_wrapper< deriv_type > wrapped_deriv_type;
|
Chris@16
|
75 typedef Time time_type;
|
Chris@16
|
76 typedef Algebra algebra_type;
|
Chris@16
|
77 typedef Operations operations_type;
|
Chris@16
|
78 typedef Resizer resizer_type;
|
Chris@16
|
79 typedef stepper_tag stepper_category;
|
Chris@16
|
80
|
Chris@16
|
81 typedef adams_moulton< Steps , State , Value , Deriv , Time , Algebra , Operations , Resizer > stepper_type;
|
Chris@16
|
82
|
Chris@16
|
83 static const size_t steps = Steps;
|
Chris@16
|
84
|
Chris@16
|
85 typedef unsigned short order_type;
|
Chris@16
|
86 static const order_type order_value = steps + 1;
|
Chris@16
|
87
|
Chris@16
|
88 typedef detail::rotating_buffer< wrapped_deriv_type , steps > step_storage_type;
|
Chris@16
|
89
|
Chris@16
|
90 adams_moulton( )
|
Chris@16
|
91 : m_coefficients() , m_dxdt() , m_resizer() ,
|
Chris@16
|
92 m_algebra_instance() , m_algebra( m_algebra_instance )
|
Chris@16
|
93 { }
|
Chris@16
|
94
|
Chris@16
|
95 adams_moulton( algebra_type &algebra )
|
Chris@16
|
96 : m_coefficients() , m_dxdt() , m_resizer() ,
|
Chris@16
|
97 m_algebra_instance() , m_algebra( algebra )
|
Chris@16
|
98 { }
|
Chris@16
|
99
|
Chris@16
|
100 adams_moulton& operator=( const adams_moulton &stepper )
|
Chris@16
|
101 {
|
Chris@16
|
102 m_dxdt = stepper.m_dxdt;
|
Chris@16
|
103 m_resizer = stepper.m_resizer;
|
Chris@16
|
104 m_algebra = stepper.m_algebra;
|
Chris@16
|
105 return *this;
|
Chris@16
|
106 }
|
Chris@16
|
107
|
Chris@16
|
108 order_type order( void ) const { return order_value; }
|
Chris@16
|
109
|
Chris@16
|
110
|
Chris@16
|
111 /*
|
Chris@16
|
112 * Version 1 : do_step( system , x , t , dt , buf );
|
Chris@16
|
113 *
|
Chris@16
|
114 * solves the forwarding problem
|
Chris@16
|
115 */
|
Chris@101
|
116 template< class System , class StateInOut , class StateIn , class ABBuf >
|
Chris@101
|
117 void do_step( System system , StateInOut &x , StateIn const & pred , time_type t , time_type dt , const ABBuf &buf )
|
Chris@16
|
118 {
|
Chris@101
|
119 do_step( system , x , pred , t , x , dt , buf );
|
Chris@16
|
120 }
|
Chris@16
|
121
|
Chris@101
|
122 template< class System , class StateInOut , class StateIn , class ABBuf >
|
Chris@101
|
123 void do_step( System system , const StateInOut &x , StateIn const & pred , time_type t , time_type dt , const ABBuf &buf )
|
Chris@16
|
124 {
|
Chris@101
|
125 do_step( system , x , pred , t , x , dt , buf );
|
Chris@16
|
126 }
|
Chris@16
|
127
|
Chris@16
|
128
|
Chris@16
|
129
|
Chris@16
|
130 /*
|
Chris@16
|
131 * Version 2 : do_step( system , in , t , out , dt , buf );
|
Chris@16
|
132 *
|
Chris@16
|
133 * solves the forwarding problem
|
Chris@16
|
134 */
|
Chris@101
|
135 template< class System , class StateIn , class PredIn , class StateOut , class ABBuf >
|
Chris@101
|
136 void do_step( System system , const StateIn &in , const PredIn &pred , time_type t , StateOut &out , time_type dt , const ABBuf &buf )
|
Chris@16
|
137 {
|
Chris@101
|
138 do_step_impl( system , in , pred , t , out , dt , buf );
|
Chris@16
|
139 }
|
Chris@16
|
140
|
Chris@101
|
141 template< class System , class StateIn , class PredIn , class StateOut , class ABBuf >
|
Chris@101
|
142 void do_step( System system , const StateIn &in , const PredIn &pred , time_type t , const StateOut &out , time_type dt , const ABBuf &buf )
|
Chris@16
|
143 {
|
Chris@101
|
144 do_step_impl( system , in , pred , t , out , dt , buf );
|
Chris@16
|
145 }
|
Chris@16
|
146
|
Chris@16
|
147
|
Chris@16
|
148
|
Chris@16
|
149 template< class StateType >
|
Chris@16
|
150 void adjust_size( const StateType &x )
|
Chris@16
|
151 {
|
Chris@16
|
152 resize_impl( x );
|
Chris@16
|
153 }
|
Chris@16
|
154
|
Chris@16
|
155 algebra_type& algebra()
|
Chris@16
|
156 { return m_algebra; }
|
Chris@16
|
157
|
Chris@16
|
158 const algebra_type& algebra() const
|
Chris@16
|
159 { return m_algebra; }
|
Chris@16
|
160
|
Chris@16
|
161
|
Chris@16
|
162 private:
|
Chris@16
|
163
|
Chris@101
|
164
|
Chris@101
|
165 template< class System , class StateIn , class PredIn , class StateOut , class ABBuf >
|
Chris@101
|
166 void do_step_impl( System system , const StateIn &in , const PredIn &pred , time_type t , StateOut &out , time_type dt , const ABBuf &buf )
|
Chris@101
|
167 {
|
Chris@101
|
168 typename odeint::unwrap_reference< System >::type &sys = system;
|
Chris@101
|
169 m_resizer.adjust_size( in , detail::bind( &stepper_type::template resize_impl<StateIn> , detail::ref( *this ) , detail::_1 ) );
|
Chris@101
|
170 sys( pred , m_dxdt.m_v , t );
|
Chris@101
|
171 detail::adams_moulton_call_algebra< steps , algebra_type , operations_type >()( m_algebra , in , out , m_dxdt.m_v , buf , m_coefficients , dt );
|
Chris@101
|
172 }
|
Chris@101
|
173
|
Chris@101
|
174
|
Chris@16
|
175 template< class StateIn >
|
Chris@16
|
176 bool resize_impl( const StateIn &x )
|
Chris@16
|
177 {
|
Chris@16
|
178 return adjust_size_by_resizeability( m_dxdt , x , typename is_resizeable<deriv_type>::type() );
|
Chris@16
|
179 }
|
Chris@16
|
180
|
Chris@16
|
181
|
Chris@16
|
182 const detail::adams_moulton_coefficients< value_type , steps > m_coefficients;
|
Chris@16
|
183 wrapped_deriv_type m_dxdt;
|
Chris@16
|
184 resizer_type m_resizer;
|
Chris@16
|
185
|
Chris@16
|
186 protected:
|
Chris@16
|
187
|
Chris@16
|
188 algebra_type m_algebra_instance;
|
Chris@16
|
189 algebra_type &m_algebra;
|
Chris@16
|
190 };
|
Chris@16
|
191
|
Chris@16
|
192
|
Chris@16
|
193
|
Chris@16
|
194
|
Chris@16
|
195 } // odeint
|
Chris@16
|
196 } // numeric
|
Chris@16
|
197 } // boost
|
Chris@16
|
198
|
Chris@16
|
199
|
Chris@16
|
200
|
Chris@16
|
201 #endif // BOOST_NUMERIC_ODEINT_STEPPER_ADAMS_MOULTON_HPP_INCLUDED
|