ODEBND
|
00001 // Copyright (C) 2012-2014 Benoit Chachuat, Imperial College London. 00002 // All Rights Reserved. 00003 // This code is published under the Eclipse Public License. 00004 00005 #ifndef MC__ODESLV_GSL_HPP 00006 #define MC__ODESLV_GSL_HPP 00007 00008 #undef MC__ODESLV_GSL_SAMPLE_DEBUG 00009 #undef MC__ODESLV_GSL_DEBUG_INTERP 00010 #undef MC__ODESLV_GSL_DEBUG_ADJOINT 00011 #define MC__ODESLV_GSL_CHECK 00012 00013 #include <stdexcept> 00014 #include <cassert> 00015 #include <cmath> 00016 #include <fstream> 00017 #include <vector> 00018 #include <sys/time.h> 00019 00020 #include "base_gsl.hpp" 00021 00022 namespace mc 00023 { 00029 template <typename T> 00030 class ODESLV_GSL: public BASE_GSL 00031 { 00032 private: 00034 gsl_odeiv2_system _sys_traj; 00035 00037 gsl_odeiv2_driver *_driver_traj; 00038 00040 gsl_odeiv2_system _sys_adj; 00041 00043 std::vector<gsl_odeiv2_driver*> _driver_adj; 00044 00046 std::vector<gsl_interp_accel*> _accel_interp; 00047 00049 std::vector<gsl_interp*> _driver_interp; 00050 00052 unsigned int _nVAR; 00053 00055 FFVar* _pVAR; 00056 00058 double *_dVAR; 00059 00061 std::list<const FFOp*> _opRHS; 00062 00064 std::list<const FFOp*> _opJAC; 00065 00067 std::list<const FFOp*> _opIC; 00068 00070 std::list<const FFOp*> _opADJ; 00071 00073 std::list<const FFOp*> _opTC; 00074 00076 const FFVar* _pRHS; 00077 00079 double* _dRHS; 00080 00082 const FFVar* _pJAC; 00083 00085 double* _dJAC; 00086 00088 const FFVar* _pIC; 00089 00091 const FFVar* _pADJ; 00092 00094 double* _dADJ; 00095 00097 const FFVar* _pTC; 00098 00100 double* _dTC; 00101 00102 public: 00106 00107 ODESLV_GSL() 00108 : BASE_GSL(), _driver_traj(0), _nVAR(0), _pVAR(0), _dVAR(0), 00109 _pRHS(0), _dRHS(0), _pJAC(0), _dJAC(0), _pIC(0), _pADJ(0), 00110 _dADJ(0), _pTC(0), _dTC(0), _mesh(0), _meshndx(0), _ifct(0) 00111 {} 00112 00114 virtual ~ODESLV_GSL() 00115 { 00116 if( _driver_traj ) gsl_odeiv2_driver_free( _driver_traj ); 00117 for( unsigned i=0; i<_driver_adj.size(); i++ ) 00118 if( _driver_adj[i] ) gsl_odeiv2_driver_free( _driver_adj[i] ); 00119 for( unsigned i=0; i<_driver_interp.size(); i++ ) 00120 if( _driver_interp[i] ) gsl_interp_free( _driver_interp[i] ); 00121 for( unsigned i=0; i<_accel_interp.size(); i++ ) 00122 if( _accel_interp[i] ) gsl_interp_accel_free( _accel_interp[i] ); 00123 delete[] _dVAR; 00124 delete[] _dRHS; 00125 delete[] _dJAC; 00126 delete[] _dADJ; 00127 delete[] _dTC; 00128 delete[] _pVAR; 00129 /* DO NOT FREE _pRHS, _pIC */ 00130 delete[] _pJAC; 00131 delete[] _pADJ; 00132 delete[] _pTC; 00133 delete[] _mesh; 00134 delete[] _meshndx; 00135 } 00136 00138 struct Options: public BASE_GSL::Options 00139 { 00141 Options(): 00142 BASE_GSL::Options(), INTERPMETH(CSPLINE), MESHPREALLOC(0), DISPLAY(0), RESRECORD(false) 00143 {} 00145 template <typename U> Options& operator= 00146 ( U&options ){ 00147 BASE_GSL::Options::operator=(options); 00148 INTERPMETH = options.INTERPMETH; 00149 MESHPREALLOC = options.MESHPREALLOC; 00150 DISPLAY = options.DISPLAY; 00151 RESRECORD = options.RESRECORD; 00152 return *this; 00153 } 00155 enum INTERPOLATION_METHOD{ 00156 LINEAR=0, 00157 CSPLINE, 00158 AKIMA 00159 } options; 00161 INTERPOLATION_METHOD INTERPMETH; 00163 bool MESHPREALLOC; 00165 int DISPLAY; 00167 bool RESRECORD; 00168 } options; 00169 00171 class Exceptions 00172 { 00173 public: 00175 enum TYPE{ 00176 UNDEF=-33 00177 }; 00179 Exceptions( TYPE ierr ) : _ierr( ierr ){} 00181 int ierr(){ return _ierr; } 00183 std::string what(){ 00184 switch( _ierr ){ 00185 case UNDEF: default: 00186 return "ODESLV_GSL::Exceptions Error due to calling a not yet implemented feature"; 00187 } 00188 } 00189 private: 00190 TYPE _ierr; 00191 }; 00192 00194 struct Results 00195 { 00197 Results 00198 ( const double tk, const unsigned nxk, const T*Ixk ): 00199 t( tk ), nx( nxk ) 00200 { X = new T[nx]; 00201 for( unsigned ix=0; ix<nx; ix++ ) X[ix] = Ixk[ix]; } 00202 Results 00203 ( const double tk, const unsigned nxk, const T*Ixk, 00204 const unsigned nyk, const T*Iyk ): 00205 t( tk ), nx( nxk+nyk ) 00206 { X = new T[nx]; unsigned i(0); 00207 for( unsigned ix=0; ix<nxk; ix++, i++ ) X[i] = Ixk[ix]; 00208 for( unsigned iy=0; iy<nyk; iy++, i++ ) X[i] = Iyk[iy]; } 00209 Results 00210 ( const Results&res ): 00211 t( res.t ), nx( res.nx ) 00212 { X = new T[nx]; 00213 for( unsigned int ix=0; ix<nx; ix++ ) X[ix] = res.X[ix]; } 00215 ~Results() 00216 { delete[] X; } 00218 double t; 00220 unsigned int nx; 00222 T* X; 00223 }; 00224 00226 Stats stats_traj; 00227 00229 Stats stats_adj; 00230 00232 STATUS states 00233 ( const unsigned int ns, const double*tk, const double*p, double**xk, 00234 double*f, std::ostream&os=std::cout ); 00235 00237 STATUS states_ASA 00238 ( const unsigned int ns, const double*tk, const double*p, double**xk, 00239 double*f, double**lk, double*df, std::ostream&os=std::cout ); 00240 00242 STATUS bounds 00243 ( const unsigned int ns, const double*tk, const T*Ip, T**Ixk, 00244 T*If, const unsigned int nsamp, std::ostream&os=std::cout ); 00245 00247 STATUS bounds_ASA 00248 ( const unsigned int ns, const double*tk, const T*Ip, T**Ixk, 00249 T*If, T**Ilk, T*Idf, const unsigned int nsamp, std::ostream&os=std::cout ); 00250 00252 void record 00253 ( std::ofstream&bndrec, const unsigned int iprec=5 ) const; 00256 00257 static ODESLV_GSL<T> *pODESLV_GSL; 00258 00259 private: 00261 std::vector< Results > _results; 00262 00264 std::vector<double>* _mesh; 00265 00267 unsigned* _meshndx; 00268 00270 unsigned int _ifct; 00271 00273 std::vector<double> _h_adj; 00274 00276 void _GSL_term 00277 ( Stats &stats ); 00278 00280 void _GSL_init 00281 ( gsl_odeiv2_system &sys, gsl_odeiv2_driver *&driver ); 00282 00284 void _GSL_init 00285 ( gsl_odeiv2_system &sys, std::vector<gsl_odeiv2_driver*>&driver ); 00286 00288 void _GSL_init 00289 ( std::vector<gsl_interp*>&driver, const unsigned size ); 00290 00292 void _GSL_init 00293 ( std::vector<gsl_interp_accel*>&accel ); 00294 00296 void _GSL_init 00297 ( const double* p, const bool store, const unsigned ns ); 00298 00300 void _GSL_init 00301 ( const double* p ); 00302 00304 bool _set_RHS 00305 ( const unsigned int iRHS ); 00306 00308 bool _set_ADJRHS 00309 ( const unsigned int iRHS ); 00310 00312 void _mesh_init 00313 ( const unsigned ns ); 00314 00316 void _mesh_add 00317 ( const double&t, const double*x, const bool index ); 00318 00320 bool _mesh_interp 00321 ( const unsigned ns ); 00322 00324 bool _mesh_eval 00325 ( const double t, double*x ); 00326 00328 STATUS _states_traj 00329 ( const double tf, double*xf, const bool store ); 00330 00332 STATUS _states 00333 ( const unsigned int ns, const double*tk, const double*p, double**xk, 00334 double*f, const bool store, std::ostream&os ); 00335 00337 STATUS _adjoints_traj 00338 ( const double tf, double&h ); 00339 00341 STATUS _adjoints 00342 ( const unsigned int ns, const double*tk, const double*p, 00343 const double*const*xk, double**lk, double*df, std::ostream&os ); 00344 00346 static int MC_GSLRHS__ 00347 ( double t, const double* x, double* xdot, void* user_data ); 00348 00350 int _ODESLV_GSL_RHS 00351 ( double t, const double* x, double* xdot, void* user_data ); 00352 00354 static int MC_GSLJAC__ 00355 ( double t, const double* x, double* jac, double* xdot, void* user_data ); 00356 00358 int _ODESLV_GSL_JAC 00359 ( double t, const double* x, double* jac, double* xdot, void* user_data ); 00360 00362 static int MC_GSLADJRHS__ 00363 ( double t, const double* l, double* ldot, void* user_data ); 00364 00366 int _ODESLV_GSL_ADJRHS 00367 ( double t, const double* l, double* ldot, void* user_data ); 00368 00370 static int MC_GSLADJJAC__ 00371 ( double t, const double* l, double* jac, double* xdot, void* user_data ); 00372 00374 int _ODESLV_GSL_ADJJAC 00375 ( double t, const double* l, double* jac, double* xdot, void* user_data ); 00376 00378 STATUS _states 00379 ( const unsigned int ns, const double*tk, const T*Ip, T**Ixk, T*If, 00380 const unsigned int nsamp, unsigned int* vsamp, const unsigned int ip, 00381 double*p, double**xk, double*f, std::ostream&os ); 00382 00384 STATUS _states_ASA 00385 ( const unsigned int ns, const double*tk, const T*Ip, T**Ixk, T*If, 00386 T**Ilk, T*Idf, const unsigned int nsamp, unsigned int* vsamp, 00387 const unsigned int ip, double*p, double**xk, double*f, double**lk, 00388 double*df, std::ostream&os ); 00389 00391 template<typename U> static void _print_interm 00392 ( const double t, const unsigned nx, const U*x, const std::string&var, 00393 std::ostream&os=std::cout ); 00394 00396 template<typename U> static void _print_interm 00397 ( const std::string&header, const unsigned nx, const U*x, const std::string&var, 00398 std::ostream&os=std::cout ); 00399 00401 ODESLV_GSL(const ODESLV_GSL&); 00402 ODESLV_GSL& operator=(const ODESLV_GSL&); 00403 }; 00404 00405 template <typename T> 00406 ODESLV_GSL<T>* ODESLV_GSL<T>::pODESLV_GSL = 0; 00407 00408 template <typename T> inline void 00409 ODESLV_GSL<T>::_GSL_init 00410 ( const double*p, const bool store, const unsigned ns ) 00411 { 00412 // Define ODE system in GSL format 00413 _sys_traj.function = MC_GSLRHS__; 00414 _sys_traj.jacobian = MC_GSLJAC__; 00415 _sys_traj.dimension = _nx; 00416 _sys_traj.params = 0; 00417 00418 // Set GSL drivers for ODE integration 00419 _GSL_init( _sys_traj, _driver_traj ); 00420 00421 // Initialize mesh 00422 if( store ) _mesh_init( ns ); 00423 00424 // Size and set DAG evaluation arrays 00425 if( _nVAR < _nx+_np+1 ){ 00426 _nVAR = _nx+_np+1; 00427 delete[] _pVAR; _pVAR = new FFVar[_nVAR]; 00428 delete[] _dVAR; _dVAR = new double[_nVAR]; 00429 } 00430 for( unsigned ix=0; ix<_nx; ix++ ) _pVAR[ix] = _pX[ix]; 00431 for( unsigned ip=0; ip<_np; ip++ ) _pVAR[_nx+ip] = _pP[ip]; 00432 _pVAR[_nx+_np] = (_pT? *_pT: 0. ); 00433 for( unsigned ip=0; ip<_np; ip++ ) _dVAR[_nx+ip] = p[ip]; 00434 00435 // Initialize statistics 00436 _init_stats( stats_traj ); 00437 00438 return; 00439 } 00440 template <typename T> inline void 00441 ODESLV_GSL<T>::_GSL_init 00442 ( const double*p ) 00443 { 00444 // Define adjoint ODE system in GSL format 00445 _sys_adj.function = MC_GSLADJRHS__; 00446 _sys_adj.jacobian = MC_GSLADJJAC__; 00447 _sys_adj.dimension = _nx+_np; 00448 _sys_adj.params = 0; 00449 00450 // Set GSL drivers for adjoint ODE integration 00451 _GSL_init( _sys_adj, _driver_adj ); 00452 00453 // Set GSL lookup accelerator for interpolation 00454 _GSL_init( _accel_interp ); 00455 00456 // Size and set DAG evaluation arrays 00457 if( _nVAR < 2*_nx+_np+1 ){ 00458 _nVAR = 2*_nx+_np+1; 00459 delete[] _pVAR; _pVAR = new FFVar[_nVAR]; 00460 delete[] _dVAR; _dVAR = new double[_nVAR]; 00461 } 00462 00463 set_adjoint(); 00464 for( unsigned ix=0; ix<_nx; ix++ ) _pVAR[ix] = _pY[ix]; 00465 for( unsigned ix=0; ix<_nx; ix++ ) _pVAR[_nx+ix] = _pX[ix]; 00466 for( unsigned ip=0; ip<_np; ip++ ) _pVAR[2*_nx+ip] = _pP[ip]; 00467 _pVAR[2*_nx+_np] = (_pT? *_pT: 0. ); 00468 00469 for( unsigned ip=0; ip<_np; ip++ ) _dVAR[2*_nx+ip] = p[ip]; 00470 00471 // Initialize statistics 00472 _init_stats( stats_adj ); 00473 00474 return; 00475 } 00476 00477 template <typename T> inline void 00478 ODESLV_GSL<T>::_GSL_init 00479 ( gsl_odeiv2_system &sys, gsl_odeiv2_driver *&driver ) 00480 { 00481 // Set GSL numerical integration driver 00482 if( driver ) gsl_odeiv2_driver_free( driver ); 00483 switch( options.INTMETH ){ 00484 case Options::RK8PD: 00485 driver = gsl_odeiv2_driver_alloc_y_new( &sys, gsl_odeiv2_step_rk8pd, 00486 options.H0, options.ATOL, options.RTOL ); 00487 break; 00488 case Options::MSADAMS: 00489 driver = gsl_odeiv2_driver_alloc_y_new( &sys, gsl_odeiv2_step_msadams, 00490 options.H0, options.ATOL, options.RTOL ); 00491 break; 00492 case Options::MSBDF: 00493 driver = gsl_odeiv2_driver_alloc_y_new( &sys, gsl_odeiv2_step_msbdf, 00494 options.H0, options.ATOL, options.RTOL ); 00495 break; 00496 case Options::RKF45: default: 00497 driver = gsl_odeiv2_driver_alloc_y_new( &sys, gsl_odeiv2_step_rkf45, 00498 options.H0, options.ATOL, options.RTOL ); 00499 break; 00500 } 00501 gsl_odeiv2_driver_set_hmin( driver, options.HMIN ); 00502 gsl_odeiv2_driver_set_nmax( driver, options.NMAX ); 00503 00504 // Set state storage 00505 delete [] _vec_state; 00506 _vec_state = new double[ sys.dimension ]; 00507 00508 return; 00509 } 00510 00511 template <typename T> inline void 00512 ODESLV_GSL<T>::_GSL_init 00513 ( gsl_odeiv2_system &sys, std::vector<gsl_odeiv2_driver*>&driver ) 00514 { 00515 // Reset GSL numerical integration drivers 00516 for( unsigned i=0; i<driver.size(); i++ ) 00517 gsl_odeiv2_driver_free( driver[i] ); 00518 driver.clear(); 00519 00520 for( unsigned i=0; i<_nf; i++ ){ 00521 switch( options.INTMETH ){ 00522 case Options::RK8PD: 00523 driver.push_back( gsl_odeiv2_driver_alloc_y_new( &sys, 00524 gsl_odeiv2_step_rk8pd, options.H0, options.ATOL, options.RTOL ) ); 00525 break; 00526 case Options::MSADAMS: 00527 driver.push_back( gsl_odeiv2_driver_alloc_y_new( &sys, 00528 gsl_odeiv2_step_msadams, options.H0, options.ATOL, options.RTOL ) ); 00529 break; 00530 case Options::MSBDF: 00531 driver.push_back( gsl_odeiv2_driver_alloc_y_new( &sys, 00532 gsl_odeiv2_step_msbdf, options.H0, options.ATOL, options.RTOL ) ); 00533 break; 00534 case Options::RKF45: default: 00535 driver.push_back( gsl_odeiv2_driver_alloc_y_new( &sys, 00536 gsl_odeiv2_step_rkf45, options.H0, options.ATOL, options.RTOL ) ); 00537 break; 00538 } 00539 gsl_odeiv2_driver_set_hmin( driver[i], options.HMIN ); 00540 gsl_odeiv2_driver_set_nmax( driver[i], options.NMAX ); 00541 } 00542 00543 // Set intermediate storage 00544 delete [] _vec_state; _vec_state = new double[ sys.dimension ]; 00545 delete [] _dTC; _dTC = new double[ sys.dimension*_nf ]; 00546 00547 return; 00548 } 00549 00550 template <typename T> inline void 00551 ODESLV_GSL<T>::_GSL_init 00552 ( std::vector<gsl_interp_accel*>&accel ) 00553 { 00554 // Set GSL interpolation lookup accelerator - expand vector if necessary 00555 for( unsigned i=accel.size(); i<_nx; i++ ) 00556 accel.push_back( gsl_interp_accel_alloc() ); 00557 } 00558 00559 template <typename T> inline void 00560 ODESLV_GSL<T>::_GSL_init 00561 ( std::vector<gsl_interp*>&driver, const unsigned size ) 00562 { 00563 // Reset GSL interpolation drivers 00564 for( unsigned i=0; i<driver.size(); i++ ) 00565 gsl_interp_free( driver[i] ); 00566 driver.clear(); 00567 00568 for( unsigned i=0; i<_nx; i++ ) 00569 switch( options.INTERPMETH ){ 00570 case Options::AKIMA: 00571 if( size >= 5 ){ 00572 driver.push_back( gsl_interp_alloc( gsl_interp_akima, size ) ); 00573 break; 00574 } 00575 case Options::CSPLINE: 00576 if( size >= 3 ){ 00577 driver.push_back( gsl_interp_alloc( gsl_interp_cspline, size ) ); 00578 break; 00579 } 00580 case Options::LINEAR: default: 00581 driver.push_back( gsl_interp_alloc( gsl_interp_linear, size ) ); 00582 break; 00583 } 00584 } 00585 00586 template <typename T> inline void 00587 ODESLV_GSL<T>::_mesh_init 00588 ( const unsigned ns ) 00589 { 00590 // Allocate mesh for holding intermediate states 00591 if( !_mesh ) _mesh = new std::vector<double>[_sys_traj.dimension+1]; 00592 for( unsigned ix=0; options.MESHPREALLOC && ix<=_sys_traj.dimension; ix++ ) 00593 { _mesh[ix].clear(); _mesh[ix].reserve( options.MESHPREALLOC ); } 00594 00595 // Allocate array for holding index at stage times 00596 delete[] _meshndx; _meshndx = new unsigned[ns]; 00597 } 00598 00599 template <typename T> inline bool 00600 ODESLV_GSL<T>::_mesh_interp 00601 ( const unsigned ns ) 00602 { 00603 const unsigned offset = _meshndx[_istg-1]; 00604 const unsigned size = ( _istg==ns? _mesh[0].size(): _meshndx[_istg] ) - offset; 00605 00606 // Initialize interpolation driver 00607 _GSL_init( _driver_interp, size ); 00608 00609 // Interpolate every state in current stage 00610 for( unsigned i=0; i<_nx; i++ ){ 00611 #ifdef MC__ODESLV_GSL_DEBUG_INTERP 00612 std::cout << "Stage #" << _istg << " State #" << i << std::endl; 00613 std::cout << _mesh[0].data() << " " << _mesh[i+1].data() << std::endl; 00614 for( unsigned k=0; k<size; k++ ) 00615 std::cout << (_mesh[0].data()+offset)[k] << " " << (_mesh[i+1].data()+offset)[k] << std::endl; 00616 { int dum; std::cin >> dum; } 00617 #endif 00618 int flag = gsl_interp_init( _driver_interp[i], _mesh[0].data()+offset, 00619 _mesh[i+1].data()+offset, size ); 00620 #ifdef MC__ODESLV_GSL_DEBUG_INTERP 00621 std::cout << "flag: " << flag << std::endl; 00622 #endif 00623 if( flag != GSL_SUCCESS ) return false; 00624 } 00625 // !!!gsl_interp_init returns an int flag, but it is unclear what it is!!! 00626 // Function: int gsl_interp_init (gsl_interp * interp, const double xa[], const double ya[], size_t size) 00627 00628 // Reset lookup accelerator 00629 #ifdef MC__ODESLV_GSL_CHECK 00630 if( _accel_interp.size() < _nx ) return false; 00631 #endif 00632 for( unsigned i=0; i<_nx; i++ ) 00633 if( gsl_interp_accel_reset( _accel_interp[i] ) != GSL_SUCCESS ) return false; 00634 00635 return true; 00636 } 00637 00638 template <typename T> inline void 00639 ODESLV_GSL<T>::_mesh_add 00640 ( const double&t, const double*x, const bool index ) 00641 { 00642 // Index current location 00643 if( index ) _meshndx[_istg] = _mesh[0].size(); 00644 00645 // Store time and state into mesh 00646 _mesh[0].push_back( t ); 00647 for( unsigned ix=0; ix<_sys_traj.dimension; ix++ ) 00648 _mesh[ix+1].push_back( x[ix] ); 00649 } 00650 00651 template <typename T> inline bool 00652 ODESLV_GSL<T>::_mesh_eval 00653 ( const double t, double*x ) 00654 { 00655 const unsigned offset = _meshndx[_istg-1]; 00656 #ifdef MC__ODESLV_GSL_CHECK 00657 if( _driver_interp.size() != _nx ) return false; 00658 #endif 00659 for( unsigned i=0; i<_nx; i++ ) 00660 if( gsl_interp_eval_e( _driver_interp[i], _mesh[0].data()+offset, 00661 _mesh[i+1].data()+offset, t, _accel_interp[i], x+i ) == GSL_EDOM ) 00662 return false; 00663 return true; 00664 } 00665 00666 template <typename T> inline void 00667 ODESLV_GSL<T>::_GSL_term 00668 ( Stats& stats ) 00669 { 00670 // Get final CPU time 00671 _final_stats( stats ); 00672 } 00673 00674 template <typename T> inline int 00675 ODESLV_GSL<T>::_ODESLV_GSL_RHS 00676 ( double t, const double* x, double* xdot, void* user_data ) 00677 { 00678 stats_traj.numRHS++; // increment RHS counter 00679 if( !_pRHS ) return GSL_EBADFUNC; // **error** RHS not defined 00680 for( unsigned i=0; i<_nx; i++ ) _dVAR[i] = x[i]; // set current state values 00681 _dVAR[_nx+_np] = t; // set current time 00682 _pDAG->eval( _opRHS, _dRHS, _nx, _pRHS, xdot, _nVAR, _pVAR, _dVAR ); 00683 00684 return GSL_SUCCESS; 00685 } 00686 00687 template <typename T> inline int 00688 ODESLV_GSL<T>::MC_GSLRHS__ 00689 ( double t, const double* x, double* xdot, void* user_data ) 00690 { 00691 ODESLV_GSL<T> *pODESLV_GSL = ODESLV_GSL<T>::pODESLV_GSL; 00692 int flag = pODESLV_GSL->_ODESLV_GSL_RHS( t, x, xdot, user_data ); 00693 ODESLV_GSL<T>::pODESLV_GSL = pODESLV_GSL; 00694 return flag; 00695 } 00696 00697 template <typename T> inline int 00698 ODESLV_GSL<T>::_ODESLV_GSL_JAC 00699 ( double t, const double* x, double* jac, double* xdot, void* user_data ) 00700 { 00701 stats_traj.numJAC++; 00702 if( !_pRHS || !_pJAC ) return GSL_EBADFUNC; // **error** RHS or JAC not defined 00703 for( unsigned i=0; i<_nx; i++ ) _dVAR[i] = x[i]; // set current state values 00704 _dVAR[_nx+_np] = t; // set current time 00705 _pDAG->eval( _opRHS, _dRHS, _nx, _pRHS, xdot, _nVAR, _pVAR, _dVAR ); 00706 _pDAG->eval( _opJAC, _dJAC, _nx*_nx, _pJAC, jac, _nVAR, _pVAR, _dVAR ); 00707 00708 return GSL_SUCCESS; 00709 } 00710 00711 template <typename T> inline int 00712 ODESLV_GSL<T>::MC_GSLJAC__ 00713 ( double t, const double* x, double* jac, double* xdot, void* user_data ) 00714 { 00715 ODESLV_GSL<T> *pODESLV_GSL = ODESLV_GSL<T>::pODESLV_GSL; 00716 int flag = pODESLV_GSL->_ODESLV_GSL_JAC( t, x, jac, xdot, user_data ); 00717 ODESLV_GSL<T>::pODESLV_GSL = pODESLV_GSL; 00718 return flag; 00719 } 00720 00721 template <typename T> inline bool 00722 ODESLV_GSL<T>::_set_RHS 00723 ( const unsigned int iRHS ) 00724 { 00725 if( _vRHS.size() <= iRHS ) return false; 00726 00727 _pRHS = _vRHS.at( iRHS ); 00728 _opRHS = _pDAG->subgraph( _nx, _pRHS ); 00729 delete[] _dRHS; _dRHS = new double[ _opRHS.size() ]; 00730 00731 switch( options.INTMETH ){ 00732 case Options::MSADAMS: 00733 case Options::MSBDF: 00734 delete[] _pJAC; _pJAC = _pDAG->FAD( _nx, _pRHS, _nx, _pX ); 00735 _opJAC = _pDAG->subgraph( _nx*_nx, _pJAC ); 00736 delete[] _dJAC; _dJAC = new double[ _opJAC.size() ]; 00737 break; 00738 case Options::RK8PD: 00739 case Options::RKF45: 00740 default: 00741 delete[] _pJAC; _pJAC = 0; 00742 _opJAC.clear(); 00743 delete[] _dJAC; _dJAC = 0; 00744 break; 00745 } 00746 00747 return true; 00748 } 00749 00750 template <typename T> inline int 00751 ODESLV_GSL<T>::_ODESLV_GSL_ADJRHS 00752 ( double t, const double* y, double* ydot, void* user_data ) 00753 { 00754 stats_adj.numRHS++; // increment RHS counter 00755 if( !_pADJ ) return GSL_EBADFUNC; // **error** RHS not defined 00756 if( !_mesh_eval( -t, _dVAR+_nx ) ) return GSL_EBADFUNC; // set interpolated state values 00757 for( unsigned i=0; i<_nx; i++ ) _dVAR[i] = y[i]; // set current adjoint/quadrature values 00758 _dVAR[2*_nx+_np] = -t; // set current time 00759 _pDAG->eval( _opADJ, _dADJ, _nx+_np, _pADJ, ydot, _nVAR, _pVAR, _dVAR ); 00760 00761 return GSL_SUCCESS; 00762 } 00763 00764 template <typename T> inline int 00765 ODESLV_GSL<T>::MC_GSLADJRHS__ 00766 ( double t, const double* y, double* ydot, void* user_data ) 00767 { 00768 ODESLV_GSL<T> *pODESLV_GSL = ODESLV_GSL<T>::pODESLV_GSL; 00769 int flag = pODESLV_GSL->_ODESLV_GSL_ADJRHS( t, y, ydot, user_data ); 00770 ODESLV_GSL<T>::pODESLV_GSL = pODESLV_GSL; 00771 return flag; 00772 } 00773 00774 template <typename T> inline int 00775 ODESLV_GSL<T>::_ODESLV_GSL_ADJJAC 00776 ( double t, const double* y, double* jac, double* ydot, void* user_data ) 00777 { 00778 stats_adj.numJAC++; 00779 if( !_pADJ || !_pJAC ) return GSL_EBADFUNC; // **error** RHS or JAC not defined 00780 if( !_mesh_eval( -t, _dVAR+_nx ) ) return GSL_EBADFUNC; // set interpolated state values 00781 for( unsigned i=0; i<_nx; i++ ) _dVAR[i] = y[i]; // set current adjoint/quadrature values 00782 _dVAR[2*_nx+_np] = -t; // set current time 00783 _pDAG->eval( _opADJ, _dADJ, _nx+_np, _pADJ, ydot, _nVAR, _pVAR, _dVAR ); 00784 _pDAG->eval( _opJAC, _dJAC, (_nx+_np)*(_nx+_np), _pJAC, jac, _nVAR, _pVAR, _dVAR ); 00785 00786 return GSL_SUCCESS; 00787 } 00788 00789 template <typename T> inline int 00790 ODESLV_GSL<T>::MC_GSLADJJAC__ 00791 ( double t, const double* y, double* jac, double* ydot, void* user_data ) 00792 { 00793 ODESLV_GSL<T> *pODESLV_GSL = ODESLV_GSL<T>::pODESLV_GSL; 00794 int flag = pODESLV_GSL->_ODESLV_GSL_ADJJAC( t, y, jac, ydot, user_data ); 00795 ODESLV_GSL<T>::pODESLV_GSL = pODESLV_GSL; 00796 return flag; 00797 } 00798 00799 template <typename T> inline bool 00800 ODESLV_GSL<T>::_set_ADJRHS 00801 ( const unsigned int iRHS ) 00802 { 00803 if( _vRHS.size() <= iRHS ) return false; 00804 00805 _pRHS = _vRHS.at( iRHS ); 00806 FFVar pHAM( 0. ); 00807 for( unsigned ix=0; ix<_nx; ix++ ) pHAM += _pVAR[ix] * _pRHS[ix]; 00808 delete[] _pADJ; _pADJ = _pDAG->BAD( 1, &pHAM, _nx+_np, _pVAR+_nx ); 00809 _opADJ = _pDAG->subgraph( _nx+_np, _pADJ ); 00810 00811 switch( options.INTMETH ){ 00812 case Options::MSADAMS: 00813 case Options::MSBDF: 00814 delete[] _pJAC; _pJAC = _pDAG->FAD( _nx+_np, _pADJ, _nx+_np, _pVAR+_nx ); 00815 _opJAC = _pDAG->subgraph( (_nx+_np)*(_nx+_np), _pJAC ); 00816 delete[] _dJAC; _dJAC = new double[ _opJAC.size() ]; 00817 break; 00818 case Options::RK8PD: 00819 case Options::RKF45: 00820 default: 00821 delete[] _pJAC; _pJAC = 0; 00822 _opJAC.clear(); 00823 delete[] _dJAC; _dJAC = 0; 00824 break; 00825 } 00826 00827 delete[] _dADJ; _dADJ = new double[ _opADJ.size() ]; 00828 00829 return true; 00830 } 00831 00832 template <typename T> inline typename ODESLV_GSL<T>::STATUS 00833 ODESLV_GSL<T>::_states_traj 00834 ( const double tf, double*xf, const bool store ) 00835 { 00836 // Store initial point 00837 if( store ) _mesh_add( _t, xf, true ); 00838 00839 // integrate till end of time stage 00840 while( _t < tf ){ 00841 if( gsl_odeiv2_evolve_apply( _driver_traj->e, _driver_traj->c, _driver_traj->s, 00842 &_sys_traj, &_t, tf, &_h, xf ) != GSL_SUCCESS 00843 || _h < options.HMIN 00844 || (options.NMAX && stats_traj.numSteps > options.NMAX) ) return FAILURE; 00845 stats_traj.numSteps++; 00846 if( options.HMAX > 0 && _h > options.HMAX ) _h = options.HMAX; 00847 00848 // Store intermediate point 00849 if( store ) _mesh_add( _t, xf, false ); 00850 } 00851 return NORMAL; 00852 } 00853 00854 template <typename T> inline typename ODESLV_GSL<T>::STATUS 00855 ODESLV_GSL<T>::_states 00856 ( const unsigned int ns, const double*tk, const double*p, double**xk, 00857 double*f, const bool store, std::ostream&os ) 00858 { 00859 // Initialize trajectory integration with GSL 00860 _GSL_init( p, store, ns ); 00861 00862 // States initial conditions 00863 if( !_vIC.size() ) 00864 { _GSL_term( stats_traj ); return FATAL; } 00865 _pIC = _vIC.at(0); 00866 _opIC = _pDAG->subgraph( _nx, _pIC ); 00867 _dVAR[_nx+_np] = tk[0]; // set initial time 00868 _pDAG->eval( _opIC, _nx, _pIC, xk[0], _np+1, _pVAR+_nx, _dVAR+_nx ); 00869 if( options.DISPLAY >= 1 ) _print_interm( tk[0], _nx, xk[0], "x", os ); 00870 00871 // Integrate ODEs through each stage using GSL 00872 _t = tk[0]; 00873 _h = options.H0; 00874 pODESLV_GSL = this; 00875 00876 for( _istg=0; _istg<ns; _istg++ ){ 00877 00878 // Reinitialize states at intermediate times 00879 if( _istg && _vIC.size() > 1 ){ 00880 _pIC = _vIC.at(_istg); 00881 _opIC = _pDAG->subgraph( _nx, _pIC ); 00882 for( unsigned i=0; i<_nx; i++ ) _dVAR[i] = xk[_istg][i]; // set current state values 00883 _dVAR[_nx+_np] = tk[_istg+1]; // set current time 00884 _pDAG->eval( _opIC, _nx, _pIC, xk[_istg+1], _nVAR, _pVAR, _dVAR ); 00885 // Reset ODE solver - needed to account for discontinuity 00886 gsl_odeiv2_driver_reset( _driver_traj ); 00887 } 00888 else 00889 for( unsigned int ix=0; ix<_nx; ix++ ) xk[_istg+1][ix] = xk[_istg][ix]; 00890 00891 // Update list of operations in RHS and JAC 00892 const unsigned int iRHS = ( _vRHS.size()<2? 0: _istg ); 00893 if( (!_istg || iRHS) && !_set_RHS( iRHS ) ) 00894 { _GSL_term( stats_traj ); return FATAL; } 00895 00896 // Integrate until end of time stage 00897 if( _states_traj( tk[_istg+1], xk[_istg+1], store ) == FAILURE ) 00898 { _GSL_term( stats_traj ); return FAILURE; } 00899 if( options.DISPLAY >= 1 ) _print_interm( tk[_istg+1], _nx, xk[_istg+1], "x", os ); 00900 } 00901 00902 // Evaluate functions 00903 if( f && _vFCT.size() == 1 ){ 00904 const FFVar* pFCT = _vFCT.at( 0 ); 00905 for( unsigned i=0; i<_nx; i++ ) _dVAR[i] = xk[ns][i]; // set final state values 00906 _dVAR[_nx+_np] = tk[ns]; // set final time 00907 _pDAG->eval( _nf, pFCT, f, _nx+_np+1, _pVAR, _dVAR, false ); 00908 if( options.DISPLAY >= 1 ) _print_interm( "functions", _nf, f, "f", os ); 00909 } 00910 else if( f && _vFCT.size() >= ns ){ 00911 for( unsigned int istg=0; istg<ns; istg++ ){ 00912 const FFVar* pFCT = _vFCT.at( istg ); 00913 for( unsigned i=0; i<_nx; i++ ) _dVAR[i] = xk[istg+1][i]; // set intermediate state values 00914 _dVAR[_nx+_np] = tk[istg+1]; // set stage time 00915 _pDAG->eval( _nf, pFCT, f, _nx+_np+1, _pVAR, _dVAR, istg?true:false ); 00916 } 00917 if( options.DISPLAY >= 1 ) _print_interm( "function", _nf, f, "f", os ); 00918 } 00919 else if( f ){ _GSL_term( stats_traj ); return FATAL; } 00920 00921 _GSL_term( stats_traj ); 00922 if( options.DISPLAY >= 1 ) _print_stats( stats_traj, os ); 00923 return NORMAL; 00924 } 00925 00926 template <typename T> inline typename ODESLV_GSL<T>::STATUS 00927 ODESLV_GSL<T>::_adjoints_traj 00928 ( const double tf, double&h ) 00929 { 00930 // integrate till end of time stage 00931 #ifdef MC__ODESLV_GSL_DEBUG_ADJOINT 00932 std::cout << "adjoint #" << _ifct << " " << "stage #" << _istg << " times:\n"; 00933 #endif 00934 while( _t < tf ){ 00935 if( gsl_odeiv2_evolve_apply( _driver_adj[_ifct]->e, _driver_adj[_ifct]->c, 00936 _driver_adj[_ifct]->s, &_sys_adj, &_t, tf, &h, _vec_state ) != GSL_SUCCESS 00937 || h < options.HMIN 00938 || (options.NMAX && stats_adj.numSteps > options.NMAX) ) return FAILURE; 00939 stats_adj.numSteps++; 00940 if( options.HMAX > 0 && h > options.HMAX ) h = options.HMAX; 00941 #ifdef MC__ODESLV_GSL_DEBUG_ADJOINT 00942 std::cout << -_t << std::endl; 00943 #endif 00944 } 00945 return NORMAL; 00946 } 00947 00948 template <typename T> inline typename ODESLV_GSL<T>::STATUS 00949 ODESLV_GSL<T>::_adjoints 00950 ( const unsigned int ns, const double*tk, const double*p, const double*const*xk, 00951 double**lk, double*q, std::ostream&os ) 00952 { 00953 // Initialize trajectory integration with GSL 00954 _GSL_init( p ); 00955 00956 // Adjoint terminal conditions 00957 if( _vFCT.size() != 1 && _vFCT.size() < ns ) 00958 { _GSL_term( stats_adj ); return FATAL; } 00959 const FFVar* pFCT = _vFCT.size()>1? _vFCT.at(ns-1): _vFCT.at(0); 00960 delete[] _pTC; _pTC = _pDAG->FAD( _nf, pFCT, _nx+_np, _pVAR+_nx ); 00961 _opTC = _pDAG->subgraph( (_nx+_np)*_nf, _pTC ); 00962 00963 for( unsigned i=0; i<_nx; i++ ) _dVAR[_nx+i] = xk[ns][i]; // set current state values 00964 _dVAR[2*_nx+_np] = tk[ns]; // set current time 00965 _pDAG->eval( _opTC, (_nx+_np)*_nf, _pTC, _dTC, _nx+_np+1, _pVAR+_nx, _dVAR+_nx ); 00966 00967 for( unsigned k=0, ik=0; k<_nf; k++ ){ 00968 for( unsigned i=0; i<_nx; i++, ik++ ) lk[ns][k*_nx+i] = _dTC[ik]; 00969 for( unsigned i=0; i<_np; i++, ik++ ) q[k*_np+i] = _dTC[ik]; 00970 } 00971 if( options.DISPLAY >= 1 ){ 00972 _print_interm( tk[ns], _nf*_nx, lk[ns], "l", os ); 00973 //_print_interm( tk[ns], _nf*_np, q, "q", os ); 00974 } 00975 00976 // Integrate adjoint ODEs backward through each stage using GSL 00977 _h_adj.assign( _nf, options.H0 ); 00978 pODESLV_GSL = this; 00979 00980 for( _istg=ns; _istg>0; _istg-- ){ 00981 00982 // Adjoint discontinuity at stage times 00983 if( _istg<ns ){ 00984 00985 if( _vFCT.size() > 1 ){ 00986 pFCT = _vFCT.at(_istg-1); 00987 delete[] _pTC; _pTC = _pDAG->FAD( _nf, pFCT, _nx+_np, _pVAR+_nx ); 00988 _opTC = _pDAG->subgraph( (_nx+_np)*_nf, _pTC ); 00989 00990 _dVAR[2*_nx+_np] = tk[_istg]; // set current time 00991 for( unsigned i=0; i<_nx; i++ ) _dVAR[_nx+i] = xk[_istg][i]; // set current state values 00992 _pDAG->eval( _opTC, (_nx+_np)*_nf, _pTC, _dTC, _nx+_np+1, _pVAR+_nx, _dVAR+_nx ); 00993 // Reset adjoint ODE solver 00994 // for( _ifct=0; _ifct<_nf; _ifct++ ) gsl_odeiv2_driver_reset( _driver_adj[_ifct] ); 00995 } 00996 else 00997 for( unsigned i=0; i<(_nx+_np)*_nf; i++ ) _dTC[i] = 0.; 00998 00999 // State discontinuity contribution 01000 if( _vIC.size() > 1 ){ 01001 _pIC = _vIC.at(_istg); 01002 FFVar pHAM( 0. ); 01003 for( unsigned ix=0; ix<_nx; ix++ ) pHAM += _pVAR[ix] * _pIC[ix]; 01004 delete[] _pTC; _pTC = _pDAG->BAD( 1, &pHAM, _nx+_np, _pVAR+_nx ); 01005 _opTC = _pDAG->subgraph( _nx+_np, _pTC ); 01006 01007 _dVAR[2*_nx+_np] = tk[_istg]; // set current time 01008 for( _ifct=0; _ifct<_nf; _ifct++ ){ 01009 for( unsigned i=0; i<_nx; i++ ){ 01010 _dVAR[i] = lk[_istg][_ifct*_nx+i] + _dTC[_ifct*(_nx+_np)+i]; // l(ti+)+dphi/dx(ti) 01011 _dTC[_ifct*(_nx+_np)+i] = 0.; 01012 _dVAR[_nx+i] = xk[_istg][i]; // set current state values 01013 } 01014 for( unsigned i=0; i<_np; i++ ) 01015 _dTC[_ifct*(_nx+_np)+_nx+i] += q[_ifct*_np+i]; 01016 _pDAG->eval( _opTC, _nx+_np, _pTC, _dTC+_ifct*(_nx+_np), _nVAR, _pVAR, _dVAR, true ); // add result to _dTC 01017 // Reset adjoint ODE solver 01018 // gsl_odeiv2_driver_reset( _driver_adj[_ifct] ); 01019 } 01020 } 01021 else 01022 for( _ifct=0; _ifct<_nf; _ifct++ ){ 01023 for( unsigned i=0; i<_nx; i++ ) _dTC[_ifct*(_nx+_np)+i] += lk[_istg][_ifct*_nx+i]; 01024 for( unsigned i=0; i<_np; i++ ) _dTC[_ifct*(_nx+_np)+_nx+i] += q[_ifct*_np+i]; 01025 } 01026 } 01027 01028 // Update list of operations in adjoint RHS and JAC 01029 const unsigned int iRHS = ( _vRHS.size()<2? 0: _istg-1 ); 01030 if( (_istg==ns || iRHS) && !_set_ADJRHS( iRHS ) ) 01031 { _GSL_term( stats_adj ); return FATAL; } 01032 01033 // Interpolate state mesh 01034 if( !_mesh_interp( ns ) ) 01035 { _GSL_term( stats_adj ); return FAILURE; } 01036 01037 // Integrate backward along time stage for each function 01038 for( _ifct=0; _ifct<_nf; _ifct++ ){ 01039 _t = -tk[_istg]; // need to reinitialize stepsize too??? 01040 for( unsigned i=0; i<_nx+_np; i++ ) _vec_state[i] = _dTC[_ifct*(_nx+_np)+i]; 01041 if( _adjoints_traj( -tk[_istg-1], _h_adj[_ifct] ) == FAILURE ) 01042 { _GSL_term( stats_adj ); return FAILURE; } 01043 for( unsigned i=0; i<_nx; i++ ) lk[_istg-1][_ifct*_nx+i] = _vec_state[i]; 01044 for( unsigned i=0; i<_np; i++ ) q[_ifct*_np+i] = _vec_state[_nx+i]; 01045 } 01046 01047 // Initial state contribution 01048 if( _istg==1 ){ 01049 _pIC = _vIC.at(0); 01050 FFVar pHAM( 0. ); 01051 for( unsigned ix=0; ix<_nx; ix++ ) pHAM += _pVAR[ix] * _pIC[ix]; 01052 delete[] _pTC; _pTC = _pDAG->BAD( 1, &pHAM, _np, _pVAR+2*_nx ); 01053 _opTC = _pDAG->subgraph( _np, _pTC ); 01054 01055 _dVAR[2*_nx+_np] = tk[0]; // set initial time 01056 for( _ifct=0; _ifct<_nf; _ifct++ ){ 01057 for( unsigned i=0; i<_nx; i++ ) _dVAR[i] = lk[0][_ifct*_nx+i]; // set initial adjoint 01058 _pDAG->eval( _opTC, _np, _pTC, q+_ifct*_np, _nVAR, _pVAR, _dVAR, true ); // add result to q 01059 } 01060 } 01061 01062 if( options.DISPLAY >= 1 ){ 01063 _print_interm( tk[_istg-1], _nf*_nx, lk[_istg-1], "l", os ); 01064 //_print_interm( tk[_istg-1], _nf*_np, q, "q", os ); 01065 } 01066 } 01067 01068 if( options.DISPLAY >= 1 ) 01069 _print_interm( "gradient", _nf*_np, q, "df", os ); 01070 01071 _GSL_term( stats_adj ); 01072 if( options.DISPLAY >= 1 ) _print_stats( stats_adj, os ); 01073 return NORMAL; 01074 } 01075 01076 01077 template <typename T> inline typename ODESLV_GSL<T>::STATUS 01078 ODESLV_GSL<T>::_states 01079 ( const unsigned int ns, const double*tk, const T*Ip, T**Ixk, T*If, 01080 const unsigned int nsamp, unsigned int* vsamp, const unsigned int ip, 01081 double*p, double**xk, double*f, std::ostream&os ) 01082 { 01083 STATUS flag = NORMAL; 01084 01085 // Update bounds for all sampling points 01086 for( unsigned int isamp=0; isamp<nsamp; isamp++ ){ 01087 vsamp[ip] = isamp; 01088 01089 // Continue recursive call 01090 if( ip+1 < _np ){ 01091 flag = _states( ns, tk, Ip, Ixk, If, nsamp, vsamp, ip+1, p, xk, f, os ); 01092 if( flag != NORMAL ) return flag; 01093 continue; 01094 } 01095 01096 // Update bounds for current point 01097 #ifdef MC__ODESLV_GSL_SAMPLE_DEBUG 01098 std::cout << "Sample: "; 01099 #endif 01100 for( unsigned int ip=0; ip<_np; ip++ ){ 01101 p[ip] = Op<T>::l( Ip[ip] ) + vsamp[ip]/(nsamp-1.) * Op<T>::diam( Ip[ip] ); 01102 #ifdef MC__ODESLV_GSL_SAMPLE_DEBUG 01103 std::cout << p[ip] << " "; 01104 #endif 01105 } 01106 #ifdef MC__ODESLV_GSL_SAMPLE_DEBUG 01107 std::cout << std::endl; 01108 #endif 01109 flag = states( ns, tk, p, xk, f, os ); 01110 if( flag != NORMAL ) return flag; 01111 for( unsigned int is=0; is<=ns; is++ ) 01112 for( unsigned int ix=0; ix<_nx; ix++ ) 01113 Ixk[is][ix] = Op<T>::hull( xk[is][ix], Ixk[is][ix] ); 01114 for( unsigned int ifn=0; If && ifn<_nf; ifn++ ) 01115 If[ifn] = Op<T>::hull( f[ifn], If[ifn] ); 01116 } 01117 01118 return flag; 01119 } 01120 01121 template <typename T> inline typename ODESLV_GSL<T>::STATUS 01122 ODESLV_GSL<T>::_states_ASA 01123 ( const unsigned int ns, const double*tk, const T*Ip, T**Ixk, T*If, 01124 T**Ilk, T*Idf, const unsigned int nsamp, unsigned int* vsamp, 01125 const unsigned int ip, double*p, double**xk, double*f, double**lk, 01126 double*df, std::ostream&os ) 01127 { 01128 STATUS flag = NORMAL; 01129 01130 // Update bounds for all sampling points 01131 for( unsigned int isamp=0; isamp<nsamp; isamp++ ){ 01132 vsamp[ip] = isamp; 01133 01134 // Continue recursive call 01135 if( ip+1 < _np ){ 01136 flag = _states_ASA( ns, tk, Ip, Ixk, If, Ilk, Idf, nsamp, vsamp, ip+1, p, xk, f, lk, df, os ); 01137 if( flag != NORMAL ) return flag; 01138 continue; 01139 } 01140 01141 // Update bounds for current point 01142 #ifdef MC__ODESLV_GSL_SAMPLE_DEBUG 01143 std::cout << "Sample: "; 01144 #endif 01145 for( unsigned int ip=0; ip<_np; ip++ ){ 01146 p[ip] = Op<T>::l( Ip[ip] ) + vsamp[ip]/(nsamp-1.) * Op<T>::diam( Ip[ip] ); 01147 #ifdef MC__ODESLV_GSL_SAMPLE_DEBUG 01148 std::cout << p[ip] << " "; 01149 #endif 01150 } 01151 #ifdef MC__ODESLV_GSL_SAMPLE_DEBUG 01152 std::cout << std::endl; 01153 #endif 01154 flag = states_ASA( ns, tk, p, xk, f, lk, df, os ); 01155 if( flag != NORMAL ) return flag; 01156 for( unsigned int is=0; is<=ns; is++ ){ 01157 for( unsigned int ix=0; ix<_nx; ix++ ) 01158 Ixk[is][ix] = Op<T>::hull( xk[is][ix], Ixk[is][ix] ); 01159 for( unsigned int ix=0; ix<_nf*_nx; ix++ ) 01160 Ilk[is][ix] = Op<T>::hull( lk[is][ix], Ilk[is][ix] ); 01161 } 01162 for( unsigned int ifn=0; If && ifn<_nf; ifn++ ){ 01163 If[ifn] = Op<T>::hull( f[ifn], If[ifn] ); 01164 for( unsigned int ip=0; ip<_np; ip++ ) 01165 Idf[ifn*_np+ip] = Op<T>::hull( df[ifn*_np+ip], Idf[ifn*_np+ip] ); 01166 } 01167 } 01168 01169 return flag; 01170 } 01171 01185 template <typename T> inline typename ODESLV_GSL<T>::STATUS 01186 ODESLV_GSL<T>::states 01187 ( const unsigned int ns, const double*tk, const double*p, double**xk, 01188 double*f, std::ostream&os ) 01189 { 01190 return _states( ns, tk, p, xk, f, false, os ); 01191 } 01192 01208 template <typename T> inline typename ODESLV_GSL<T>::STATUS 01209 ODESLV_GSL<T>::states_ASA 01210 ( const unsigned int ns, const double*tk, const double*p, double**xk, 01211 double*f, double**lk, double*df, std::ostream&os ) 01212 { 01213 STATUS flag = NORMAL; 01214 01215 flag = _states( ns, tk, p, xk, f, true, os ); 01216 if( flag != NORMAL ) return flag; 01217 01218 flag = _adjoints( ns, tk, p, xk, lk, df, os ); 01219 return flag; 01220 } 01221 01238 template <typename T> inline typename ODESLV_GSL<T>::STATUS 01239 ODESLV_GSL<T>::bounds 01240 ( const unsigned int ns, const double*tk, const T*Ip, T**Ixk, 01241 T*If, const unsigned int nsamp, std::ostream&os ) 01242 { 01243 int DISPLAY_SAVE = options.DISPLAY; 01244 options.DISPLAY = 0; 01245 STATUS flag = NORMAL; 01246 01247 // Initialization of sampled bounds at parameter lower bound 01248 double *p = new double[_np]; 01249 for( unsigned int ip=0; ip<_np; ip++ ) 01250 p[ip] = Op<T>::l(Ip[ip]); 01251 double **xk = new double*[ns+1]; 01252 for( unsigned int is=0; is<=ns; is++ ) 01253 xk[is] = new double[_nx]; 01254 double *f = If? new double[_nf]: 0; 01255 flag = states( ns, tk, p, xk, f, os ); 01256 if( flag != NORMAL || nsamp <= 1 ){ 01257 delete[] p; delete[] f; 01258 for( unsigned int is=0; is<=ns; is++ ) delete[] xk[is]; delete[] xk; 01259 return flag; 01260 } 01261 for( unsigned int is=0; is<=ns; is++ ) 01262 for( unsigned int ix=0; ix<_nx; ix++ ) 01263 Ixk[is][ix] = xk[is][ix]; 01264 for( unsigned int ifn=0; If && ifn<_nf; ifn++ ) 01265 If[ifn] = f[ifn]; 01266 01267 // Start sampling process 01268 unsigned int* vsamp = new unsigned int[_np]; 01269 flag = _states( ns, tk, Ip, Ixk, If, nsamp, vsamp, 0, p, xk, f, os ); 01270 01271 // Display results 01272 options.DISPLAY = DISPLAY_SAVE; 01273 if( options.DISPLAY >= 1 ){ 01274 for( unsigned int is=0; is<=ns; is++ ) 01275 _print_interm( tk[is], _nx, Ixk[is], "x", os ); 01276 if( If ) _print_interm( "function", _nf, If, "f", os ); 01277 } 01278 01279 // Record intermediate results 01280 _results.clear(); 01281 if( options.RESRECORD ) 01282 for( unsigned int is=0; is<=ns; is++ ) 01283 _results.push_back( Results( tk[is], _nx, Ixk[is] ) ); 01284 01285 // Clean-up 01286 delete[] p; delete[] f; 01287 for( unsigned int is=0; is<=ns; is++ ) delete[] xk[is]; delete[] xk; 01288 delete[] vsamp; 01289 01290 return flag; 01291 } 01292 01311 template <typename T> inline typename ODESLV_GSL<T>::STATUS 01312 ODESLV_GSL<T>::bounds_ASA 01313 ( const unsigned int ns, const double*tk, const T*Ip, T**Ixk, 01314 T*If, T**Ilk, T*Idf, const unsigned int nsamp, std::ostream&os ) 01315 { 01316 int DISPLAY_SAVE = options.DISPLAY; 01317 options.DISPLAY = 0; 01318 STATUS flag = NORMAL; 01319 01320 // Initialization of sampled bounds at parameter lower bound 01321 double *p = new double[_np]; 01322 for( unsigned int ip=0; ip<_np; ip++ ) 01323 p[ip] = Op<T>::l(Ip[ip]); 01324 double **xk = new double*[ns+1]; 01325 double **lk = new double*[ns+1]; 01326 for( unsigned int is=0; is<=ns; is++ ){ 01327 xk[is] = new double[_nx]; 01328 lk[is] = new double[_nf*_nx]; 01329 } 01330 double *f = If? new double[_nf]: 0; 01331 double *df = new double[_nf*_nx]; 01332 flag = states_ASA( ns, tk, p, xk, f, lk, df, os ); 01333 if( flag != NORMAL || nsamp <= 1 ){ 01334 delete[] p; 01335 for( unsigned int is=0; is<=ns; is++ ) delete[] xk[is]; delete[] xk; 01336 for( unsigned int is=0; is<=ns; is++ ) delete[] lk[is]; delete[] lk; 01337 delete[] f; delete[] df; 01338 return flag; 01339 } 01340 for( unsigned int is=0; is<=ns; is++ ){ 01341 for( unsigned int ix=0; ix<_nx; ix++ ) 01342 Ixk[is][ix] = xk[is][ix]; 01343 for( unsigned int ix=0; ix<_nf*_nx; ix++ ) 01344 Ilk[is][ix] = lk[is][ix]; 01345 } 01346 for( unsigned int ifn=0; If && ifn<_nf; ifn++ ){ 01347 If[ifn] = f[ifn]; 01348 for( unsigned int ip=0; ip<_np; ip++ ) 01349 Idf[ifn*_np+ip] = df[ifn*_np+ip]; 01350 } 01351 01352 // Start sampling process 01353 unsigned int* vsamp = new unsigned int[_np]; 01354 flag = _states_ASA( ns, tk, Ip, Ixk, If, Ilk, Idf, nsamp, vsamp, 0, p, xk, f, lk, df, os ); 01355 01356 // Display results 01357 options.DISPLAY = DISPLAY_SAVE; 01358 if( options.DISPLAY >= 1 ){ 01359 for( unsigned int is=0; is<=ns; is++ ) 01360 _print_interm( tk[is], _nx, Ixk[is], "x", os ); 01361 if( If ) _print_interm( "function", _nf, If, "f", os ); 01362 for( unsigned int is=0; is<=ns; is++ ) 01363 _print_interm( tk[ns-is], _nx*_nf, Ilk[ns-is], "l", os ); 01364 _print_interm( "gradient", _nf*_np, Idf, "df", os ); 01365 } 01366 01367 // Record intermediate results 01368 _results.clear(); 01369 if( options.RESRECORD ) 01370 for( unsigned int is=0; is<=ns; is++ ) 01371 _results.push_back( Results( tk[is], _nx, Ixk[is], _nf*_nx, Ilk[is] ) ); 01372 01373 // Clean-up 01374 delete[] p; delete[] f; delete[] df; 01375 for( unsigned int is=0; is<=ns; is++ ) delete[] xk[is]; delete[] xk; 01376 for( unsigned int is=0; is<=ns; is++ ) delete[] lk[is]; delete[] lk; 01377 delete[] vsamp; 01378 01379 return flag; 01380 } 01381 01382 template <typename T> inline void 01383 ODESLV_GSL<T>::record 01384 ( std::ofstream&bndrec, const unsigned int iprec ) const 01385 { 01386 if( !bndrec ) return; 01387 01388 // Specify format 01389 bndrec << std::right << std::scientific << std::setprecision(iprec); 01390 01391 // Record computed interval bounds at stage times 01392 typename std::vector< Results >::const_iterator it = _results.begin(); 01393 for( ; it != _results.end(); ++it ){ 01394 bndrec << std::setw(iprec+9) << (*it).t; 01395 for( unsigned int ix=0; ix<(*it).nx; ix++ ) 01396 bndrec << std::setw(iprec+9) << mc::Op<T>::l( (*it).X[ix] ) 01397 << std::setw(iprec+9) << mc::Op<T>::u( (*it).X[ix] ); 01398 bndrec << std::endl; 01399 } 01400 } 01401 01402 template <typename T> template<typename U> inline void 01403 ODESLV_GSL<T>::_print_interm 01404 ( const double t, const unsigned nx, const U*x, const std::string&var, 01405 std::ostream&os ) 01406 { 01407 os << " @t = " << std::scientific << std::setprecision(6) 01408 << std::left << t << " :" << std::endl; 01409 for( unsigned int ix=0; ix<nx; ix++ ) 01410 os << " " << var.c_str() << "[" << ix << "] = " << x[ix] << std::endl; 01411 return; 01412 } 01413 01414 template <typename T> template<typename U> inline void 01415 ODESLV_GSL<T>::_print_interm 01416 ( const std::string&header, const unsigned nx, const U*x, const std::string&var, 01417 std::ostream&os ) 01418 { 01419 os << header.c_str() << " :" << std::endl 01420 << std::scientific << std::setprecision(6) << std::left; 01421 for( unsigned int ix=0; ix<nx; ix++ ) 01422 os << " " << var.c_str() << "[" << ix << "] = " << x[ix] << std::endl; 01423 return; 01424 } 01425 01426 } // end namescape mc 01427 01428 #endif 01429