Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:42:57

0001 /*
0002  [auto_generated]
0003  boost/numeric/odeint/stepper/rosenbrock4.hpp
0004 
0005  [begin_description]
0006  Implementation of the Rosenbrock 4 method for solving stiff ODEs. Note, that a
0007  controller and a dense-output stepper exist for this method,
0008  [end_description]
0009 
0010  Copyright 2011-2013 Karsten Ahnert
0011  Copyright 2011-2012 Mario Mulansky
0012  Copyright 2012 Christoph Koke
0013 
0014  Distributed under the Boost Software License, Version 1.0.
0015  (See accompanying file LICENSE_1_0.txt or
0016  copy at http://www.boost.org/LICENSE_1_0.txt)
0017  */
0018 
0019 
0020 #ifndef BOOST_NUMERIC_ODEINT_STEPPER_ROSENBROCK4_HPP_INCLUDED
0021 #define BOOST_NUMERIC_ODEINT_STEPPER_ROSENBROCK4_HPP_INCLUDED
0022 
0023 
0024 #include <boost/numeric/odeint/util/bind.hpp>
0025 #include <boost/numeric/odeint/util/unwrap_reference.hpp>
0026 #include <boost/numeric/ublas/vector.hpp>
0027 #include <boost/numeric/ublas/matrix.hpp>
0028 #include <boost/numeric/ublas/lu.hpp>
0029 
0030 #include <boost/numeric/odeint/stepper/stepper_categories.hpp>
0031 
0032 #include <boost/numeric/odeint/util/ublas_wrapper.hpp>
0033 #include <boost/numeric/odeint/util/is_resizeable.hpp>
0034 #include <boost/numeric/odeint/util/resizer.hpp>
0035 
0036 #include <boost/numeric/ublas/vector.hpp>
0037 #include <boost/numeric/ublas/matrix.hpp>
0038 #include <boost/numeric/ublas/lu.hpp>
0039 
0040 
0041 namespace boost {
0042 namespace numeric {
0043 namespace odeint {
0044 
0045 
0046 /*
0047  * ToDo:
0048  *
0049  * 2. Interfacing for odeint, check if controlled_error_stepper can be used
0050  * 3. dense output
0051  */
0052 
0053 
0054 
0055 template< class Value >
0056 struct default_rosenbrock_coefficients
0057 {
0058     typedef Value value_type;
0059     typedef unsigned short order_type;
0060 
0061     default_rosenbrock_coefficients( void )
0062     : gamma ( static_cast< value_type >( 0.25 ) ) ,
0063       d1 ( static_cast< value_type >( 0.25 ) ) ,
0064       d2 ( static_cast< value_type >( -0.1043 ) ) ,
0065       d3 ( static_cast< value_type >( 0.1035 ) ) ,
0066       d4 ( static_cast< value_type >( 0.3620000000000023e-01 ) ) ,
0067       c2 ( static_cast< value_type >( 0.386 ) ) ,
0068       c3 ( static_cast< value_type >( 0.21 ) ) ,
0069       c4 ( static_cast< value_type >( 0.63 ) ) ,
0070       c21 ( static_cast< value_type >( -0.5668800000000000e+01 ) ) ,
0071       a21 ( static_cast< value_type >( 0.1544000000000000e+01 ) ) ,
0072       c31 ( static_cast< value_type >( -0.2430093356833875e+01 ) ) ,
0073       c32 ( static_cast< value_type >( -0.2063599157091915e+00 ) ) ,
0074       a31 ( static_cast< value_type >( 0.9466785280815826e+00 ) ) ,
0075       a32 ( static_cast< value_type >( 0.2557011698983284e+00 ) ) ,
0076       c41 ( static_cast< value_type >( -0.1073529058151375e+00 ) ) ,
0077       c42 ( static_cast< value_type >( -0.9594562251023355e+01 ) ) ,
0078       c43 ( static_cast< value_type >( -0.2047028614809616e+02 ) ) ,
0079       a41 ( static_cast< value_type >( 0.3314825187068521e+01 ) ) ,
0080       a42 ( static_cast< value_type >( 0.2896124015972201e+01 ) ) ,
0081       a43 ( static_cast< value_type >( 0.9986419139977817e+00 ) ) ,
0082       c51 ( static_cast< value_type >( 0.7496443313967647e+01 ) ) ,
0083       c52 ( static_cast< value_type >( -0.1024680431464352e+02 ) ) ,
0084       c53 ( static_cast< value_type >( -0.3399990352819905e+02 ) ) ,
0085       c54 ( static_cast< value_type >(  0.1170890893206160e+02 ) ) ,
0086       a51 ( static_cast< value_type >( 0.1221224509226641e+01 ) ) ,
0087       a52 ( static_cast< value_type >( 0.6019134481288629e+01 ) ) ,
0088       a53 ( static_cast< value_type >( 0.1253708332932087e+02 ) ) ,
0089       a54 ( static_cast< value_type >( -0.6878860361058950e+00 ) ) ,
0090       c61 ( static_cast< value_type >( 0.8083246795921522e+01 ) ) ,
0091       c62 ( static_cast< value_type >( -0.7981132988064893e+01 ) ) ,
0092       c63 ( static_cast< value_type >( -0.3152159432874371e+02 ) ) ,
0093       c64 ( static_cast< value_type >( 0.1631930543123136e+02 ) ) ,
0094       c65 ( static_cast< value_type >( -0.6058818238834054e+01 ) ) ,
0095       d21 ( static_cast< value_type >( 0.1012623508344586e+02 ) ) ,
0096       d22 ( static_cast< value_type >( -0.7487995877610167e+01 ) ) ,
0097       d23 ( static_cast< value_type >( -0.3480091861555747e+02 ) ) ,
0098       d24 ( static_cast< value_type >( -0.7992771707568823e+01 ) ) ,
0099       d25 ( static_cast< value_type >( 0.1025137723295662e+01 ) ) ,
0100       d31 ( static_cast< value_type >( -0.6762803392801253e+00 ) ) ,
0101       d32 ( static_cast< value_type >( 0.6087714651680015e+01 ) ) ,
0102       d33 ( static_cast< value_type >( 0.1643084320892478e+02 ) ) ,
0103       d34 ( static_cast< value_type >( 0.2476722511418386e+02 ) ) ,
0104       d35 ( static_cast< value_type >( -0.6594389125716872e+01 ) )
0105     {}
0106 
0107     const value_type gamma;
0108     const value_type d1 , d2 , d3 , d4;
0109     const value_type c2 , c3 , c4;
0110     const value_type c21 ;
0111     const value_type a21;
0112     const value_type c31 , c32;
0113     const value_type a31 , a32;
0114     const value_type c41 , c42 , c43;
0115     const value_type a41 , a42 , a43;
0116     const value_type c51 , c52 , c53 , c54;
0117     const value_type a51 , a52 , a53 , a54;
0118     const value_type c61 , c62 , c63 , c64 , c65;
0119     const value_type d21 , d22 , d23 , d24 , d25;
0120     const value_type d31 , d32 , d33 , d34 , d35;
0121 
0122     static const order_type stepper_order = 4;
0123     static const order_type error_order = 3;
0124 };
0125 
0126 
0127 
0128 template< class Value , class Coefficients = default_rosenbrock_coefficients< Value > , class Resizer = initially_resizer >
0129 class rosenbrock4
0130 {
0131 private:
0132 
0133 public:
0134 
0135     typedef Value value_type;
0136     typedef boost::numeric::ublas::vector< value_type > state_type;
0137     typedef state_type deriv_type;
0138     typedef value_type time_type;
0139     typedef boost::numeric::ublas::matrix< value_type > matrix_type;
0140     typedef boost::numeric::ublas::permutation_matrix< size_t > pmatrix_type;
0141     typedef Resizer resizer_type;
0142     typedef Coefficients rosenbrock_coefficients;
0143     typedef stepper_tag stepper_category;
0144     typedef unsigned short order_type;
0145 
0146     typedef state_wrapper< state_type > wrapped_state_type;
0147     typedef state_wrapper< deriv_type > wrapped_deriv_type;
0148     typedef state_wrapper< matrix_type > wrapped_matrix_type;
0149     typedef state_wrapper< pmatrix_type > wrapped_pmatrix_type;
0150 
0151     typedef rosenbrock4< Value , Coefficients , Resizer > stepper_type;
0152 
0153     const static order_type stepper_order = rosenbrock_coefficients::stepper_order;
0154     const static order_type error_order = rosenbrock_coefficients::error_order;
0155 
0156     rosenbrock4( void )
0157     : m_resizer() , m_x_err_resizer() ,
0158       m_jac() , m_pm() ,
0159       m_dfdt() , m_dxdt() , m_dxdtnew() ,
0160       m_g1() , m_g2() , m_g3() , m_g4() , m_g5() ,
0161       m_cont3() , m_cont4() , m_xtmp() , m_x_err() ,
0162       m_coef()
0163     { }
0164 
0165 
0166     order_type order() const { return stepper_order; } 
0167 
0168     template< class System >
0169     void do_step( System system , const state_type &x , time_type t , state_type &xout , time_type dt , state_type &xerr )
0170     {
0171         // get the system and jacobi function
0172         typedef typename odeint::unwrap_reference< System >::type system_type;
0173         typedef typename odeint::unwrap_reference< typename system_type::first_type >::type deriv_func_type;
0174         typedef typename odeint::unwrap_reference< typename system_type::second_type >::type jacobi_func_type;
0175         system_type &sys = system;
0176         deriv_func_type &deriv_func = sys.first;
0177         jacobi_func_type &jacobi_func = sys.second;
0178 
0179         const size_t n = x.size();
0180 
0181         m_resizer.adjust_size( x , detail::bind( &stepper_type::template resize_impl<state_type> , detail::ref( *this ) , detail::_1 ) );
0182 
0183         for( size_t i=0 ; i<n ; ++i )
0184             m_pm.m_v( i ) = i;
0185 
0186         deriv_func( x , m_dxdt.m_v , t );
0187         jacobi_func( x , m_jac.m_v , t , m_dfdt.m_v );
0188 
0189         m_jac.m_v *= -1.0;
0190         m_jac.m_v += 1.0 / m_coef.gamma / dt * boost::numeric::ublas::identity_matrix< value_type >( n );
0191         boost::numeric::ublas::lu_factorize( m_jac.m_v , m_pm.m_v );
0192 
0193         for( size_t i=0 ; i<n ; ++i )
0194             m_g1.m_v[i] = m_dxdt.m_v[i] + dt * m_coef.d1 * m_dfdt.m_v[i];
0195         boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , m_g1.m_v );
0196 
0197 
0198         for( size_t i=0 ; i<n ; ++i )
0199             m_xtmp.m_v[i] = x[i] + m_coef.a21 * m_g1.m_v[i];
0200         deriv_func( m_xtmp.m_v , m_dxdtnew.m_v , t + m_coef.c2 * dt );
0201         for( size_t i=0 ; i<n ; ++i )
0202             m_g2.m_v[i] = m_dxdtnew.m_v[i] + dt * m_coef.d2 * m_dfdt.m_v[i] + m_coef.c21 * m_g1.m_v[i] / dt;
0203         boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , m_g2.m_v );
0204 
0205 
0206         for( size_t i=0 ; i<n ; ++i )
0207             m_xtmp.m_v[i] = x[i] + m_coef.a31 * m_g1.m_v[i] + m_coef.a32 * m_g2.m_v[i];
0208         deriv_func( m_xtmp.m_v , m_dxdtnew.m_v , t + m_coef.c3 * dt );
0209         for( size_t i=0 ; i<n ; ++i )
0210             m_g3.m_v[i] = m_dxdtnew.m_v[i] + dt * m_coef.d3 * m_dfdt.m_v[i] + ( m_coef.c31 * m_g1.m_v[i] + m_coef.c32 * m_g2.m_v[i] ) / dt;
0211         boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , m_g3.m_v );
0212 
0213 
0214         for( size_t i=0 ; i<n ; ++i )
0215             m_xtmp.m_v[i] = x[i] + m_coef.a41 * m_g1.m_v[i] + m_coef.a42 * m_g2.m_v[i] + m_coef.a43 * m_g3.m_v[i];
0216         deriv_func( m_xtmp.m_v , m_dxdtnew.m_v , t + m_coef.c4 * dt );
0217         for( size_t i=0 ; i<n ; ++i )
0218             m_g4.m_v[i] = m_dxdtnew.m_v[i] + dt * m_coef.d4 * m_dfdt.m_v[i] + ( m_coef.c41 * m_g1.m_v[i] + m_coef.c42 * m_g2.m_v[i] + m_coef.c43 * m_g3.m_v[i] ) / dt;
0219         boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , m_g4.m_v );
0220 
0221 
0222         for( size_t i=0 ; i<n ; ++i )
0223             m_xtmp.m_v[i] = x[i] + m_coef.a51 * m_g1.m_v[i] + m_coef.a52 * m_g2.m_v[i] + m_coef.a53 * m_g3.m_v[i] + m_coef.a54 * m_g4.m_v[i];
0224         deriv_func( m_xtmp.m_v , m_dxdtnew.m_v , t + dt );
0225         for( size_t i=0 ; i<n ; ++i )
0226             m_g5.m_v[i] = m_dxdtnew.m_v[i] + ( m_coef.c51 * m_g1.m_v[i] + m_coef.c52 * m_g2.m_v[i] + m_coef.c53 * m_g3.m_v[i] + m_coef.c54 * m_g4.m_v[i] ) / dt;
0227         boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , m_g5.m_v );
0228 
0229         for( size_t i=0 ; i<n ; ++i )
0230             m_xtmp.m_v[i] += m_g5.m_v[i];
0231         deriv_func( m_xtmp.m_v , m_dxdtnew.m_v , t + dt );
0232         for( size_t i=0 ; i<n ; ++i )
0233             xerr[i] = m_dxdtnew.m_v[i] + ( m_coef.c61 * m_g1.m_v[i] + m_coef.c62 * m_g2.m_v[i] + m_coef.c63 * m_g3.m_v[i] + m_coef.c64 * m_g4.m_v[i] + m_coef.c65 * m_g5.m_v[i] ) / dt;
0234         boost::numeric::ublas::lu_substitute( m_jac.m_v , m_pm.m_v , xerr );
0235 
0236         for( size_t i=0 ; i<n ; ++i )
0237             xout[i] = m_xtmp.m_v[i] + xerr[i];
0238     }
0239 
0240     template< class System >
0241     void do_step( System system , state_type &x , time_type t , time_type dt , state_type &xerr )
0242     {
0243         do_step( system , x , t , x , dt , xerr );
0244     }
0245 
0246     /*
0247      * do_step without error output - just calls above functions with and neglects the error estimate
0248      */
0249     template< class System >
0250     void do_step( System system , const state_type &x , time_type t , state_type &xout , time_type dt )
0251     {
0252         m_x_err_resizer.adjust_size( x , detail::bind( &stepper_type::template resize_x_err<state_type> , detail::ref( *this ) , detail::_1 ) );
0253         do_step( system , x , t , xout , dt , m_x_err.m_v );
0254     }
0255 
0256     template< class System >
0257     void do_step( System system , state_type &x , time_type t , time_type dt )
0258     {
0259         m_x_err_resizer.adjust_size( x , detail::bind( &stepper_type::template resize_x_err<state_type> , detail::ref( *this ) , detail::_1 ) );
0260         do_step( system , x , t , dt , m_x_err.m_v );
0261     }
0262 
0263     void prepare_dense_output()
0264     {
0265         const size_t n = m_g1.m_v.size();
0266         for( size_t i=0 ; i<n ; ++i )
0267         {
0268             m_cont3.m_v[i] = m_coef.d21 * m_g1.m_v[i] + m_coef.d22 * m_g2.m_v[i] + m_coef.d23 * m_g3.m_v[i] + m_coef.d24 * m_g4.m_v[i] + m_coef.d25 * m_g5.m_v[i];
0269             m_cont4.m_v[i] = m_coef.d31 * m_g1.m_v[i] + m_coef.d32 * m_g2.m_v[i] + m_coef.d33 * m_g3.m_v[i] + m_coef.d34 * m_g4.m_v[i] + m_coef.d35 * m_g5.m_v[i];
0270         }
0271     }
0272 
0273 
0274     void calc_state( time_type t , state_type &x ,
0275             const state_type &x_old , time_type t_old ,
0276             const state_type &x_new , time_type t_new )
0277     {
0278         const size_t n = m_g1.m_v.size();
0279         time_type dt = t_new - t_old;
0280         time_type s = ( t - t_old ) / dt;
0281         time_type s1 = 1.0 - s;
0282         for( size_t i=0 ; i<n ; ++i )
0283             x[i] = x_old[i] * s1 + s * ( x_new[i] + s1 * ( m_cont3.m_v[i] + s * m_cont4.m_v[i] ) );
0284     }
0285 
0286 
0287 
0288     template< class StateType >
0289     void adjust_size( const StateType &x )
0290     {
0291         resize_impl( x );
0292         resize_x_err( x );
0293     }
0294 
0295 
0296 protected:
0297 
0298     template< class StateIn >
0299     bool resize_impl( const StateIn &x )
0300     {
0301         bool resized = false;
0302         resized |= adjust_size_by_resizeability( m_dxdt , x , typename is_resizeable<deriv_type>::type() );
0303         resized |= adjust_size_by_resizeability( m_dfdt , x , typename is_resizeable<deriv_type>::type() );
0304         resized |= adjust_size_by_resizeability( m_dxdtnew , x , typename is_resizeable<deriv_type>::type() );
0305         resized |= adjust_size_by_resizeability( m_xtmp , x , typename is_resizeable<state_type>::type() );
0306         resized |= adjust_size_by_resizeability( m_g1 , x , typename is_resizeable<state_type>::type() );
0307         resized |= adjust_size_by_resizeability( m_g2 , x , typename is_resizeable<state_type>::type() );
0308         resized |= adjust_size_by_resizeability( m_g3 , x , typename is_resizeable<state_type>::type() );
0309         resized |= adjust_size_by_resizeability( m_g4 , x , typename is_resizeable<state_type>::type() );
0310         resized |= adjust_size_by_resizeability( m_g5 , x , typename is_resizeable<state_type>::type() );
0311         resized |= adjust_size_by_resizeability( m_cont3 , x , typename is_resizeable<state_type>::type() );
0312         resized |= adjust_size_by_resizeability( m_cont4 , x , typename is_resizeable<state_type>::type() );
0313         resized |= adjust_size_by_resizeability( m_jac , x , typename is_resizeable<matrix_type>::type() );
0314         resized |= adjust_size_by_resizeability( m_pm , x , typename is_resizeable<pmatrix_type>::type() );
0315         return resized;
0316     }
0317 
0318     template< class StateIn >
0319     bool resize_x_err( const StateIn &x )
0320     {
0321         return adjust_size_by_resizeability( m_x_err , x , typename is_resizeable<state_type>::type() );
0322     }
0323 
0324 private:
0325 
0326 
0327     resizer_type m_resizer;
0328     resizer_type m_x_err_resizer;
0329 
0330     wrapped_matrix_type m_jac;
0331     wrapped_pmatrix_type m_pm;
0332     wrapped_deriv_type m_dfdt , m_dxdt , m_dxdtnew;
0333     wrapped_state_type m_g1 , m_g2 , m_g3 , m_g4 , m_g5;
0334     wrapped_state_type m_cont3 , m_cont4;
0335     wrapped_state_type m_xtmp;
0336     wrapped_state_type m_x_err;
0337 
0338     const rosenbrock_coefficients m_coef;
0339 };
0340 
0341 
0342 } // namespace odeint
0343 } // namespace numeric
0344 } // namespace boost
0345 
0346 #endif // BOOST_NUMERIC_ODEINT_STEPPER_ROSENBROCK4_HPP_INCLUDED