details_state_types_algebras_operations.qbk 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482
  1. [/============================================================================
  2. Boost.odeint
  3. Copyright 2011-2012 Karsten Ahnert
  4. Copyright 2011-2013 Mario Mulansky
  5. Copyright 2012 Sylwester Arabas
  6. Use, modification and distribution is subject to the Boost Software License,
  7. Version 1.0. (See accompanying file LICENSE_1_0.txt or copy at
  8. http://www.boost.org/LICENSE_1_0.txt)
  9. =============================================================================/]
  10. [section State types, algebras and operations]
  11. In odeint the stepper algorithms are implemented independently of the
  12. underlying fundamental mathematical operations.
  13. This is realized by giving the user full control over the state type and the
  14. mathematical operations for this state type.
  15. Technically, this is done by introducing three concepts: StateType, Algebra,
  16. Operations.
  17. Most of the steppers in odeint expect three class types fulfilling these
  18. concepts as template parameters.
  19. Note that these concepts are not fully independent of each other but rather a
  20. valid combination must be provided in order to make the steppers work.
  21. In the following we will give some examples on reasonable
  22. state_type-algebra-operations combinations.
  23. For the most common state types, like `vector<double>` or `array<double,N>`
  24. the default values range_algebra and default_operations are perfectly fine and
  25. odeint can be used as is without worrying about algebra/operations at all.
  26. [important state_type, algebra and operations are not independent, a valid
  27. combination must be provided to make odeint work properly]
  28. Moreover, as odeint handles the memory required for intermediate temporary
  29. objects itself, it also needs knowledge about how to create state_type objects
  30. and maybe how to allocate memory (resizing).
  31. All in all, the following things have to be taken care of when odeint is used
  32. with non-standard state types:
  33. * construction/destruction
  34. * resizing (if possible/required)
  35. * algebraic operations
  36. Again, odeint already provides basic interfaces for most of the usual state
  37. types.
  38. So if you use a `std::vector`, or a `boost::array` as state type no additional
  39. work is required, they just work out of the box.
  40. [section Construction/Resizing]
  41. We distinguish between two basic state types: fixed sized and dynamically
  42. sized.
  43. For fixed size state types the default constructor `state_type()` already
  44. allocates the required memory, prominent example is `boost::array<T,N>`.
  45. Dynamically sized types have to be resized to make sure enough memory is
  46. allocated, the standard constructor does not take care of the resizing.
  47. Examples for this are the STL containers like `vector<double>`.
  48. The most easy way of getting your own state type to work with odeint is to use
  49. a fixed size state, base calculations on the range_algebra and provide the
  50. following functionality:
  51. [table
  52. [[Name] [Expression] [Type] [Semantics]]
  53. [[Construct State] [`State x()`] [`void`] [Creates an instance of `State`
  54. and allocates memory.] ]
  55. [[Begin of the sequence] [boost::begin(x)] [Iterator] [Returns an iterator
  56. pointing to the begin of the sequence]]
  57. [[End of the sequence] [boost::end(x)] [Iterator] [Returns an iterator
  58. pointing to the end of the sequence]]
  59. ]
  60. [warning If your state type does not allocate memory by default construction,
  61. you [*must define it as resizeable] and provide resize functionality (see
  62. below). Otherwise segmentation faults will occur.]
  63. So fixed sized arrays supported by __boost_range immediately work with odeint.
  64. For dynamically sized arrays one has to additionally supply the resize
  65. functionality.
  66. First, the state has to be tagged as resizeable by specializing the struct
  67. `is_resizeable` which consists of one typedef and one bool value:
  68. [table
  69. [[Name] [Expression] [Type] [Semantics]]
  70. [[Resizability] [`is_resizeable<State>::type`]
  71. [`boost::true_type` or `boost::false_type`]
  72. [Determines resizeability of the state type, returns `boost::true_type` if
  73. the state is resizeable.]]
  74. [[Resizability] [`is_resizeable<State>::value`]
  75. [`bool`]
  76. [Same as above, but with `bool` value.]]
  77. ]
  78. Defining `type` to be `true_type` and `value` as `true` tells odeint that your
  79. state is resizeable.
  80. By default, odeint now expects the support of `boost::size(x)` and a
  81. `x.resize( boost::size(y) )` member function for resizing:
  82. [table
  83. [[Name] [Expression] [Type] [Semantics]]
  84. [[Get size] [`boost::size( x )`]
  85. [`size_type`] [Returns the current size of x.]]
  86. [[Resize] [`x.resize( boost::size( y ) )`]
  87. [`void`] [Resizes x to have the same size as y.]]
  88. ]
  89. [section Using the container interface]
  90. [import ../examples/my_vector.cpp]
  91. As a first example we take the most simple case and implement our own vector
  92. `my_vector` which will provide a container interface.
  93. This makes __boost_range working out-of-box.
  94. We add a little functionality to our vector which makes it allocate some
  95. default capacity by construction.
  96. This is helpful when using resizing as then a resize can be assured to not
  97. require a new allocation.
  98. [my_vector]
  99. The only thing that has to be done other than defining is thus declaring
  100. my_vector as resizeable:
  101. [my_vector_resizeable]
  102. If we wouldn't specialize the `is_resizeable` template, the code would still
  103. compile but odeint would not adjust the size of temporary internal instances
  104. of my_vector and hence try to fill zero-sized vectors resulting in
  105. segmentation faults!
  106. The full example can be found in [github_link examples/my_vector.cpp my_vector.cpp]
  107. [endsect]
  108. [section std::list]
  109. If your state type does work with __boost_range, but handles resizing
  110. differently you are required to specialize two implementations used by odeint
  111. to check a state's size and to resize:
  112. [table
  113. [[Name] [Expression] [Type] [Semantics]]
  114. [[Check size] [`same_size_impl<State,State>::same_size(x , y)`]
  115. [`bool`] [Returns true if the size of x equals the size of y.]]
  116. [[Resize] [`resize_impl<State,State>::resize(x , y)`]
  117. [`void`] [Resizes x to have the same size as y.]]
  118. ]
  119. As an example we will use a `std::list` as state type in odeint.
  120. Because `std::list` is not supported by `boost::size` we have to replace the
  121. same_size and resize implementation to get list to work with odeint.
  122. The following code shows the required template specializations:
  123. [import ../examples/list_lattice.cpp]
  124. [list_bindings]
  125. With these definitions odeint knows how to resize `std::list`s and so they can
  126. be used as state types.
  127. A complete example can be found in [github_link examples/list_lattice.cpp list_lattice.cpp].
  128. [endsect]
  129. [endsect]
  130. [section Algebras and Operations]
  131. To provide maximum flexibility odeint is implemented in a highly modularized
  132. way. This means it is possible to change the underlying mathematical
  133. operations without touching the integration algorithms.
  134. The fundamental mathematical operations are those of a vector space, that is
  135. addition of `state_types` and multiplication of `state_type`s with a scalar
  136. (`time_type`). In odeint this is realized in two concepts: _Algebra_ and
  137. _Operations_.
  138. The standard way how this works is by the range algebra which provides
  139. functions that apply a specific operation to each of the individual elements
  140. of a container based on the __boost_range library.
  141. If your state type is not supported by __boost_range there are several
  142. possibilities to tell odeint how to do algebraic operations:
  143. * Implement `boost::begin` and `boost::end` for your state type so it works
  144. with __boost_range.
  145. * Implement vector-vector addition operator `+` and scalar-vector
  146. multiplication operator `*` and use the non-standard `vector_space_algebra`.
  147. * Implement your own algebra that implements the required functions.
  148. [section GSL Vector]
  149. In the following example we will try to use the `gsl_vector` type from __gsl (GNU
  150. Scientific Library) as state type in odeint.
  151. We will realize this by implementing a wrapper around the gsl_vector that
  152. takes care of construction/destruction.
  153. Also, __boost_range is extended such that it works with `gsl_vector`s as well
  154. which required also the implementation of a new `gsl_iterator`.
  155. [note odeint already includes all the code presented here, see [github_link
  156. boost/numeric/odeint/external/gsl/gsl_wrapper.hpp gsl_wrapper.hpp], so `gsl_vector`s
  157. can be used straight out-of-box.
  158. The following description is just for educational purpose.]
  159. The GSL is a C library, so `gsl_vector` has neither constructor, nor
  160. destructor or any `begin` or `end` function, no iterators at all.
  161. So to make it work with odeint plenty of things have to be implemented.
  162. Note that all of the work shown here is already included in odeint, so using
  163. `gsl_vector`s in odeint doesn't require any further adjustments.
  164. We present it here just as an educational example.
  165. We start with defining appropriate constructors and destructors.
  166. This is done by specializing the `state_wrapper` for `gsl_vector`.
  167. State wrappers are used by the steppers internally to create and manage
  168. temporary instances of state types:
  169. ``
  170. template<>
  171. struct state_wrapper< gsl_vector* >
  172. {
  173. typedef double value_type;
  174. typedef gsl_vector* state_type;
  175. typedef state_wrapper< gsl_vector* > state_wrapper_type;
  176. state_type m_v;
  177. state_wrapper( )
  178. {
  179. m_v = gsl_vector_alloc( 1 );
  180. }
  181. state_wrapper( const state_wrapper_type &x )
  182. {
  183. resize( m_v , x.m_v );
  184. gsl_vector_memcpy( m_v , x.m_v );
  185. }
  186. ~state_wrapper()
  187. {
  188. gsl_vector_free( m_v );
  189. }
  190. };
  191. ``
  192. This `state_wrapper` specialization tells odeint how gsl_vectors are created,
  193. copied and destroyed.
  194. Next we need resizing, this is required because gsl_vectors are dynamically
  195. sized objects:
  196. ``
  197. template<>
  198. struct is_resizeable< gsl_vector* >
  199. {
  200. typedef boost::true_type type;
  201. const static bool value = type::value;
  202. };
  203. template <>
  204. struct same_size_impl< gsl_vector* , gsl_vector* >
  205. {
  206. static bool same_size( const gsl_vector* x , const gsl_vector* y )
  207. {
  208. return x->size == y->size;
  209. }
  210. };
  211. template <>
  212. struct resize_impl< gsl_vector* , gsl_vector* >
  213. {
  214. static void resize( gsl_vector* x , const gsl_vector* y )
  215. {
  216. gsl_vector_free( x );
  217. x = gsl_vector_alloc( y->size );
  218. }
  219. };
  220. ``
  221. Up to now, we defined creation/destruction and resizing, but gsl_vectors also
  222. don't support iterators, so we first implement a gsl iterator:
  223. ``
  224. /*
  225. * defines an iterator for gsl_vector
  226. */
  227. class gsl_vector_iterator
  228. : public boost::iterator_facade< gsl_vector_iterator , double ,
  229. boost::random_access_traversal_tag >
  230. {
  231. public :
  232. gsl_vector_iterator( void ): m_p(0) , m_stride( 0 ) { }
  233. explicit gsl_vector_iterator( gsl_vector *p ) : m_p( p->data ) , m_stride( p->stride ) { }
  234. friend gsl_vector_iterator end_iterator( gsl_vector * );
  235. private :
  236. friend class boost::iterator_core_access;
  237. friend class const_gsl_vector_iterator;
  238. void increment( void ) { m_p += m_stride; }
  239. void decrement( void ) { m_p -= m_stride; }
  240. void advance( ptrdiff_t n ) { m_p += n*m_stride; }
  241. bool equal( const gsl_vector_iterator &other ) const { return this->m_p == other.m_p; }
  242. bool equal( const const_gsl_vector_iterator &other ) const;
  243. double& dereference( void ) const { return *m_p; }
  244. double *m_p;
  245. size_t m_stride;
  246. };
  247. ``
  248. A similar class exists for the `const` version of the iterator.
  249. Then we have a function returning the end iterator (similarly for `const` again):
  250. ``
  251. gsl_vector_iterator end_iterator( gsl_vector *x )
  252. {
  253. gsl_vector_iterator iter( x );
  254. iter.m_p += iter.m_stride * x->size;
  255. return iter;
  256. }
  257. ``
  258. Finally, the bindings for __boost_range are added:
  259. ``
  260. // template<>
  261. inline gsl_vector_iterator range_begin( gsl_vector *x )
  262. {
  263. return gsl_vector_iterator( x );
  264. }
  265. // template<>
  266. inline gsl_vector_iterator range_end( gsl_vector *x )
  267. {
  268. return end_iterator( x );
  269. }
  270. ``
  271. Again with similar definitions for the `const` versions.
  272. This eventually makes odeint work with gsl vectors as state types.
  273. The full code for these bindings is found in [github_link
  274. boost/numeric/odeint/external/gsl/gsl_wrapper.hpp gsl_wrapper.hpp].
  275. It might look rather complicated but keep in mind that gsl is a pre-compiled C
  276. library.
  277. [endsect]
  278. [section Vector Space Algebra]
  279. As seen above, the standard way of performing algebraic operations on
  280. container-like state types in odeint is to iterate through the elements of the
  281. container and perform the operations element-wise on the underlying value type.
  282. This is realized by means of the `range_algebra` that uses __boost_range for
  283. obtaining iterators of the state types.
  284. However, there are other ways to implement the algebraic operations on
  285. containers, one of which is defining the addition/multiplication operators for
  286. the containers directly and then using the `vector_space_algebra`.
  287. If you use this algebra, the following operators have to be defined for the
  288. state_type:
  289. [table
  290. [[Name] [Expression] [Type] [Semantics]]
  291. [[Addition] [`x + y`] [`state_type`] [Calculates the vector sum 'x+y'.]]
  292. [[Assign addition] [`x += y`] [`state_type`] [Performs x+y in place.]]
  293. [[Scalar multiplication] [`a * x `] [`state_type`] [Performs multiplication of vector x with scalar a.]]
  294. [[Assign scalar multiplication] [`x *= a`] [`state_type`] [Performs in-place multiplication of vector x with scalar a.]]
  295. ]
  296. Defining these operators makes your state type work with any basic Runge-Kutta
  297. stepper.
  298. However, if you want to use step-size control, some more functionality is
  299. required.
  300. Specifically, operations like
  301. [' max[sub i]( |err[sub i]| / (alpha * |s[sub i]|) )]
  302. have to be performed.
  303. ['err] and ['s] are state_types, alpha is a scalar.
  304. As you can see, we need element wise absolute value and division as well as an
  305. reduce operation to get the maximum value.
  306. So for controlled steppers the following things have to be implemented:
  307. [table
  308. [[Name] [Expression] [Type] [Semantics]]
  309. [[Division] [`x / y`] [`state_type`] [Calculates the element-wise division 'x/y']]
  310. [[Absolute value] [`abs( x )`] [`state_type`] [Element wise absolute value]]
  311. [[Reduce] [`vector_space_reduce_impl< state_type >::reduce( state , operation , init )`] [`value_type`]
  312. [Performs the `operation` for subsequently each element of `state` and returns the aggregate value.
  313. E.g.
  314. `init = operator( init , state[0] );`
  315. `init = operator( init , state[1] )`
  316. `...`
  317. ]]
  318. ]
  319. [endsect]
  320. [/
  321. [section Boost.Ublas]
  322. As an example for the employment of the `vector_space_algebra` we will adopt
  323. `ublas::vector` from __ublas to work as a state type in odeint.
  324. This is particularly easy because `ublas::vector` supports vector-vector
  325. addition and scalar-vector multiplication described above as well as `boost::size`.
  326. It also has a resize member function so all that has to be done in this case
  327. is to declare resizability:
  328. [import ../examples/ublas/lorenz_ublas.cpp]
  329. [ublas_resizeable]
  330. Now ublas::vector can be used as state type for simple Runge-Kutta steppers
  331. in odeint by specifying the `vector_space_algebra` as algebra in the template
  332. parameter list of the stepper.
  333. The following code shows the corresponding definitions:
  334. [ublas_main]
  335. Note again, that we haven't supported the requirements for controlled steppers,
  336. but only for simple Runge-Kutta methods.
  337. You can find the full example in [github_link
  338. examples/ublas/lorenz_ublas.cpp lorenz_ublas.cpp].
  339. [endsect]
  340. /]
  341. [section Point type]
  342. [import ../examples/lorenz_point.cpp]
  343. Here we show how to implement the required operators on a state type.
  344. As example we define a new class `point3D` representing a three-dimensional
  345. vector with components x,y,z and define addition and scalar multiplication
  346. operators for it.
  347. We use __boost_operators to reduce the amount of code to be written.
  348. The class for the point type looks as follows:
  349. [point3D]
  350. By deriving from __boost_operators classes we don't have to define outer class
  351. operators like `operator+( point3D , point3D )` because that is taken care of
  352. by the operators library.
  353. Note that for simple Runge-Kutta schemes (like `runge_kutta4`) only the `+`
  354. and `*` operators are required.
  355. If, however, a controlled stepper is used one also needs to specify the
  356. division operator `/` because calculation of the error term involves an
  357. element wise division of the state types.
  358. Additionally, controlled steppers require an `abs` function calculating the
  359. element-wise absolute value for the state type:
  360. [point3D_abs_div]
  361. Finally, we have to provide a specialization to calculate the infintity norm of a state:
  362. [point3D_norm]
  363. Again, note that the two last steps were only required if you want to use
  364. controlled steppers.
  365. For simple steppers definition of the simple `+=` and `*=` operators are
  366. sufficient.
  367. Having defined such a point type, we can easily perform the integration on a Lorenz
  368. system by explicitely configuring the `vector_space_algebra` in the stepper's
  369. template argument list:
  370. [point3D_main]
  371. The whole example can be found in [github_link
  372. examples/lorenz_point.cpp lorenz_point.cpp]
  373. [note For the most `state_types`, odeint is able to automatically determine
  374. the correct algebra and operations. But if you want to use your own `state_type`, as in this
  375. example with `point3D`, you have to manually configure the right
  376. algebra/operations, unless your `state_type` works with the default choice of
  377. `range_algebra` and `default_operations`.]
  378. [endsect]
  379. [endsect]
  380. gsl_vector, gsl_matrix, ublas::matrix, blitz::matrix, thrust
  381. [section Adapt your own operations]
  382. to be continued
  383. *thrust
  384. *gsl_complex
  385. *min, max, pow
  386. [endsect]
  387. [endsect]