implicit_euler_mtl4.hpp 4.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161
  1. /*
  2. [begin_description]
  3. Modification of the implicit Euler method, works with the MTL4 matrix library only.
  4. [end_description]
  5. Copyright 2012-2013 Andreas Angelopoulos
  6. Copyright 2012-2013 Karsten Ahnert
  7. Copyright 2012-2013 Mario Mulansky
  8. Distributed under the Boost Software License, Version 1.0.
  9. (See accompanying file LICENSE_1_0.txt or
  10. copy at http://www.boost.org/LICENSE_1_0.txt)
  11. */
  12. #ifndef BOOST_NUMERIC_ODEINT_EXTERNAL_IMPLICIT_EULER_MTL4_HPP_INCLUDED
  13. #define BOOST_NUMERIC_ODEINT_EXTERNAL_IMPLICIT_EULER_MTL4_HPP_INCLUDED
  14. #include <utility>
  15. #include <boost/numeric/odeint/util/bind.hpp>
  16. #include <boost/numeric/odeint/util/unwrap_reference.hpp>
  17. #include <boost/numeric/odeint/stepper/stepper_categories.hpp>
  18. #include <boost/numeric/odeint/external/mtl4/mtl4_resize.hpp>
  19. #include <boost/numeric/mtl/mtl.hpp>
  20. #include <boost/numeric/itl/itl.hpp>
  21. namespace boost {
  22. namespace numeric {
  23. namespace odeint {
  24. template< class ValueType , class Resizer = initially_resizer >
  25. class implicit_euler_mtl4
  26. {
  27. public:
  28. typedef ValueType value_type;
  29. typedef value_type time_type;
  30. typedef mtl::dense_vector<value_type> state_type;
  31. typedef state_wrapper< state_type > wrapped_state_type;
  32. typedef state_type deriv_type;
  33. typedef state_wrapper< deriv_type > wrapped_deriv_type;
  34. typedef mtl::compressed2D< value_type > matrix_type;
  35. typedef state_wrapper< matrix_type > wrapped_matrix_type;
  36. typedef Resizer resizer_type;
  37. typedef stepper_tag stepper_category;
  38. typedef implicit_euler_mtl4< ValueType , Resizer > stepper_type;
  39. implicit_euler_mtl4( const value_type epsilon = 1E-6 )
  40. : m_epsilon( epsilon ) , m_resizer() ,
  41. m_dxdt() , m_x() ,
  42. m_identity() , m_jacobi()
  43. { }
  44. template< class System >
  45. void do_step( System system , state_type &x , time_type t , time_type dt )
  46. {
  47. typedef typename odeint::unwrap_reference< System >::type system_type;
  48. typedef typename odeint::unwrap_reference< typename system_type::first_type >::type deriv_func_type;
  49. typedef typename odeint::unwrap_reference< typename system_type::second_type >::type jacobi_func_type;
  50. system_type &sys = system;
  51. deriv_func_type &deriv_func = sys.first;
  52. jacobi_func_type &jacobi_func = sys.second;
  53. m_resizer.adjust_size(x, [this](auto&& arg) { return this->resize_impl<StateIn>(std::forward<decltype(arg)>(arg)); });
  54. m_identity.m_v = 1;
  55. t += dt;
  56. m_x.m_v = x;
  57. deriv_func( x , m_dxdt.m_v , t );
  58. jacobi_func( x , m_jacobi.m_v , t );
  59. m_dxdt.m_v *= -dt;
  60. m_jacobi.m_v *= dt;
  61. m_jacobi.m_v -= m_identity.m_v ;
  62. // using ilu_0 preconditioning -incomplete LU factorisation
  63. // itl::pc::diagonal<matrix_type,double> L(m_jacobi.m_v);
  64. itl::pc::ilu_0<matrix_type> L( m_jacobi.m_v );
  65. solve( m_jacobi.m_v , m_x.m_v , m_dxdt.m_v , L );
  66. x+= m_x.m_v;
  67. }
  68. template< class StateType >
  69. void adjust_size( const StateType &x )
  70. {
  71. resize_impl( x );
  72. }
  73. private:
  74. /*
  75. Applying approximate iterative linear solvers
  76. default solver is Biconjugate gradient stabilized method
  77. itl::bicgstab(A, x, b, L, iter);
  78. */
  79. template < class LinearOperator, class HilbertSpaceX, class HilbertSpaceB, class Preconditioner>
  80. void solve(const LinearOperator& A, HilbertSpaceX& x, const HilbertSpaceB& b,
  81. const Preconditioner& L, int max_iteractions =500)
  82. {
  83. // Termination criterion: r < 1e-6 * b or N iterations
  84. itl::basic_iteration< double > iter( b , max_iteractions , 1e-6 );
  85. itl::bicgstab( A , x , b , L , iter );
  86. }
  87. template< class StateIn >
  88. bool resize_impl( const StateIn &x )
  89. {
  90. bool resized = false;
  91. resized |= adjust_size_by_resizeability( m_dxdt , x , typename is_resizeable<deriv_type>::type() );
  92. resized |= adjust_size_by_resizeability( m_x , x , typename is_resizeable<state_type>::type() );
  93. resized |= adjust_size_by_resizeability( m_identity , x , typename is_resizeable<matrix_type>::type() );
  94. resized |= adjust_size_by_resizeability( m_jacobi , x , typename is_resizeable<matrix_type>::type() );
  95. return resized;
  96. }
  97. private:
  98. value_type m_epsilon;
  99. resizer_type m_resizer;
  100. wrapped_deriv_type m_dxdt;
  101. wrapped_state_type m_x;
  102. wrapped_matrix_type m_identity;
  103. wrapped_matrix_type m_jacobi;
  104. };
  105. } // odeint
  106. } // numeric
  107. } // boost
  108. #endif // BOOST_NUMERIC_ODEINT_EXTERNAL_IMPLICIT_EULER_MTL4_HPP_INCLUDED