ODEBND
|
Classes | |
struct | mc::ODEBND_GSL< T, PMT, PVT >::Results |
Integration results at a given time instant. More... | |
struct | mc::ODEBND_GSL< T, PMT, PVT >::Options |
Integrator options. More... | |
class | mc::ODEBND_GSL< T, PMT, PVT >::Exceptions |
Structure for setting up storing the solver exceptions. More... | |
Functions | |
mc::ODEBND_GSL< T, PMT, PVT >::ODEBND_GSL () | |
Default constructor. | |
virtual | mc::ODEBND_GSL< T, PMT, PVT >::~ODEBND_GSL () |
Virtual destructor. | |
STATUS | mc::ODEBND_GSL< T, PMT, PVT >::bounds (const unsigned int ns, const double *tk, const T *Ip, T **Ixk, std::ostream &os=std::cout) |
Computes interval enclosure of reachable set of parametric ODEs. | |
STATUS | mc::ODEBND_GSL< T, PMT, PVT >::bounds (const unsigned int ns, const double *tk, const PVT *PMp, PVT **PMxk, E *Exk=0, std::ostream &os=std::cout) |
Computes polynomial model enclosure of reachable set of parametric ODEs. | |
STATUS | mc::ODEBND_GSL< T, PMT, PVT >::hausdorff (const unsigned int ns, const double *tk, const T *Ip, double **Hxk, const unsigned int nsamp, std::ostream &os=std::cout) |
Computes Hausdorff distance between interval enclosure and actual reachable set of parametric ODEs, using parameter sampling. | |
STATUS | mc::ODEBND_GSL< T, PMT, PVT >::hausdorff (const unsigned int ns, const double *tk, const PVT *PMp, double **Hxk, const unsigned int nsamp, std::ostream &os=std::cout) |
Computes Hausdorff distance between polynomial model remainder enclosure and actual remainder function range, using parameter sampling. | |
void | mc::ODEBND_GSL< T, PMT, PVT >::record (std::ofstream &bndrec, const unsigned int iprec=5) const |
Record results in file bndrec, with accuracy of iprec digits. | |
double | mc::ODEBND_GSL< T, PMT, PVT >::final_time () const |
Return value of final time reached. | |
Variables | |
Stats | mc::ODEBND_GSL< T, PMT, PVT >::stats_traj |
Statistics for state bounds integration. | |
mc::BASE_GSL::BASE_GSL () | |
Default class constructor. | |
virtual | mc::BASE_GSL::~BASE_GSL () |
Class destructor. |
ODEBND_GSL< T, PMT, PVT >::STATUS mc::ODEBND_GSL< T, PMT, PVT >::bounds | ( | const unsigned int | ns, |
const double * | tk, | ||
const T * | Ip, | ||
T ** | Ixk, | ||
std::ostream & | os = std::cout |
||
) | [inline] |
This function computes an interval enclosure of the reachable set of the parametric ODEs defined in IVP using equally spaced samples:
The return value is the status.
ODEBND_GSL< T, PMT, PVT >::STATUS mc::ODEBND_GSL< T, PMT, PVT >::bounds | ( | const unsigned int | ns, |
const double * | tk, | ||
const PVT * | PMp, | ||
PVT ** | PMxk, | ||
E * | Exk = 0 , |
||
std::ostream & | os = std::cout |
||
) | [inline] |
This function computes an enclosure of the reachable set of the parametric ODEs using propagation of polynomial models with convex remainders (intervals, ellipsoids):
The return value is the status.
ODEBND_GSL< T, PMT, PVT >::STATUS mc::ODEBND_GSL< T, PMT, PVT >::hausdorff | ( | const unsigned int | ns, |
const double * | tk, | ||
const T * | Ip, | ||
double ** | Hxk, | ||
const unsigned int | nsamp, | ||
std::ostream & | os = std::cout |
||
) | [inline] |
This function computes the Hausdorff distance between the interval enclosure and the exact reachable set projected onto each variable:
The return value is the status.
ODEBND_GSL< T, PMT, PVT >::STATUS mc::ODEBND_GSL< T, PMT, PVT >::hausdorff | ( | const unsigned int | ns, |
const double * | tk, | ||
const PVT * | PMp, | ||
double ** | Hxk, | ||
const unsigned int | nsamp, | ||
std::ostream & | os = std::cout |
||
) | [inline] |
This function computes the Hausdorff distance between the polynomial model remainder and the actual (sampled) range of the remainder function in projection onto each variable and for each stage time remainder and the actual range of the remainder function: