CAPD::DynSys Library
6.0.0
|
PoicareMap class rigorously computes Poincare Map. More...
#include <capd/poincare/PoincareMap.h>
Public Types | |
typedef SolverT | Solver |
ODE solver. More... | |
typedef Solver::VectorFieldType | VectorFieldType |
vector field type More... | |
typedef Solver::MatrixType | MatrixType |
typedef Solver::VectorType | VectorType |
typedef Solver::ScalarType | ScalarType |
typedef ScalarType::BoundType | BoundType |
a type of the end for intervals More... | |
typedef Solver::SolutionCurve | CurveType |
typedef VectorType::size_type | size_type |
integral type used to index containers (vectors, matrices, etc) More... | |
typedef dynset::AffineCoordinateChange< MatrixType > | AffineCoordinateChange |
typedef SectionT | SectionType |
type of function that represents Poincare section. More... | |
typedef SectionType::SectionDerivativesEnclosureType | SectionDerivativesEnclosureType |
typedef capd::dynset::C0Set< MatrixType > | C0Set |
type of abstract base class for all C0 sets More... | |
typedef capd::dynset::C1Set< MatrixType > | C1Set |
type of abstract base class for all C1 sets More... | |
typedef capd::dynset::C2Set< MatrixType > | C2Set |
type of abstract base class for all C2 sets More... | |
typedef capd::dynset::CnSet< MatrixType, 0 > | CnSet |
type of abstract base class for all C2 sets More... | |
typedef C2Set::HessianType | HessianType |
type of data structure for storing of Hessians of multivariate mappings. More... | |
typedef CnSet::JetType | JetType |
typedef BasicPoincareMap< Solver, SectionT >::CrossingDirection | CrossingDirection |
section crossing direction; possible values are PlusMinus, Both(default) and MinusPlus. More... | |
typedef TypeTraits< ScalarType >::Real | RealType |
floating point type: for floating points is equal to ScalarType, for intervals is a type of their bounds. More... | |
typedef capd::diffAlgebra::SolutionCurve< CurveType > | SolutionCurve |
Public Member Functions | |
PoincareMap (Solver &solver, SectionType §ion, CrossingDirection direction=Both, const BoundType &errorTolerance=0.1) | |
Constructs PoincareMap for given dynamical system and section. More... | |
template<typename T > | |
VectorType | operator() (T &theSet, const VectorType &x0, const MatrixType &A, ScalarType &out_returnTime, int n=1) |
Computes value of n-th iterate of the Poincare Map. The result is returned in given affine coordinate system. More... | |
template<typename T > | |
VectorType | operator() (T &theSet, ScalarType &out_returnTime, int n=1) |
Computes value of n-th iterate of the Poincare Map (used for C^0 sets) and returns bound of return time. More... | |
template<typename T > | |
VectorType | operator() (T &theSet, int n=1) |
Computes value of n-th iterate of the Poincare Map (used for C^0 sets) More... | |
template<typename T > | |
VectorType | operator() (T &theSet, MatrixType &dF, ScalarType &out_returnTime, int n=1) |
Computes value of n-th iterate of Poincare Map and derivatives of the flow (used for C^1 sets) More... | |
template<typename T > | |
VectorType | operator() (T &theSet, MatrixType &dF, int n=1) |
Computes value of n-th iterate of Poincare Map and derivatives of the flow (used for C^1 sets) More... | |
template<typename T > | |
VectorType | operator() (T &theSet, MatrixType &dF, HessianType &hessian, int n=1) |
Computes value of n-th iterate of Poincare Map, derivatives and hessian of the flow (used for C^2 sets) More... | |
template<typename T > | |
VectorType | operator() (T &theSet, MatrixType &dF, HessianType &hessian, ScalarType &out_returnTime, int n=1) |
Computes value of n-th iterate of Poincare Map, derivatives and hessian of the flow (used for C^2 sets) More... | |
template<typename T > | |
VectorType | operator() (T &theSet, typename T::JetType &result, int n=1) |
Computes n-th iterate of Poincare Map and derivatives of the flow to given order (used for C^n sets) More... | |
template<typename T > | |
VectorType | operator() (T &theSet, typename T::JetType &result, ScalarType &out_returnTime, int n=1) |
Computes n-th iterate of Poincare Map and derivatives of the flow to given order (used for C^n sets) More... | |
const SectionDerivativesEnclosureType & | getSectionDerivativesEnclosure () const |
VectorType | operator() (const VectorType &v) |
Computes value of Poincare Map. More... | |
VectorType | operator() (VectorType v, ScalarType &in_out_time) |
Computes value of Poincare Map. More... | |
VectorType | operator() (const VectorType &v, VectorType &afterSection) |
Computes value of Poincare Map. More... | |
VectorType | operator() (VectorType v, VectorType &afterSection, ScalarType &in_out_time) |
Computes value of Poincare Map. More... | |
VectorType | operator() (const VectorType &v, MatrixType &dF) |
Computes value of Poincare Map and derivativeof the flow. More... | |
VectorType | operator() (VectorType v, MatrixType &dF, ScalarType &in_out_time) |
Computes value of Poincare Map and derivative of the flow. More... | |
VectorType | operator() (const VectorType &v, MatrixType &dF, HessianType &h) |
Computes value of Poincare Map, derivative and hessian of the flow. More... | |
VectorType | operator() (VectorType v, MatrixType &dF, HessianType &h, ScalarType &in_out_time) |
Computes value of Poincare Map, derivative and hessian of the flow. More... | |
VectorType | operator() (JetType &x) |
Computes Poincare Map and derivatives of the flow to given order evaluated at return time. More... | |
VectorType | operator() (JetType &x, ScalarType &in_out_time) |
Computes Poincare Map and derivatives of the flow to given order evaluated at return time. More... | |
MatrixType | computeDP (const VectorType &Px, const MatrixType &derivativeOfFlow, VectorType &dT, ScalarType returnTime=TypeTraits< ScalarType >::zero()) |
Simultaneous computation of gradient of return time and derivative of Poincare Map dP. More... | |
MatrixType | computeDP (const VectorType &Px, const MatrixType &derivativeOfFlow, ScalarType returnTime=TypeTraits< ScalarType >::zero()) |
Computes derivative of Poincare Map dP. More... | |
void | computeDP (const VectorType &Px, const MatrixType &derivativeOfFlow, const HessianType &hessianOfFlow, MatrixType &DP, HessianType &D2P, VectorType &dT, MatrixType &d2T, ScalarType returnTime=TypeTraits< ScalarType >::zero()) |
Simultaneous computation of first and second Taylor coefficients of return time and Poincare map. More... | |
void | computeDP (const VectorType &Px, const MatrixType &derivativeOfFlow, const HessianType &hessianOfFlow, MatrixType &DP, HessianType &D2P, ScalarType returnTime=TypeTraits< ScalarType >::zero()) |
Simultaneous computation of first and second Taylor coefficients of return time and Poincare map. More... | |
JetT | computeDP (const JetT &DPhi) |
Recomputes Taylor expansion of the flow into Taylor expansion of Poincare map. More... | |
const Solver & | getSolver () const |
Returns read-only reference to solver used to integrate the system. More... | |
Solver & | getSolver () |
Returns reference to solver used to integrate the system. More... | |
const Solver & | getDynamicalSystem () const |
Returns read-only reference to solver used to integrate the system. More... | |
Solver & | getDynamicalSystem () |
Returns reference to solver used to integrate the system. More... | |
VectorFieldType & | getVectorField () |
Returns reference to solver used to integrate the system. More... | |
const VectorFieldType & | getVectorField () const |
Returns read-only reference to solver used to integrate the system. More... | |
const SectionType & | getSection () const |
Returns reference to Poincare section object. More... | |
size_type | getOrder () const |
Returns order of numerical method used. More... | |
ScalarType | getStep () const |
Returns time step. More... | |
void | setOrder (size_type newOrder) |
Sets order of Taylor method. More... | |
void | setSection (const SectionType &newSection) |
Sets new section function. More... | |
void | setStep (const ScalarType &newStep) |
Sets time step for integration of DS. More... | |
void | setFactor (double newFactor) |
void | setCrossingDirection (CrossingDirection cs) |
void | turnOnStepControl () |
Disables automatic step control. More... | |
void | turnOffStepControl () |
Enables automatic step control. Step control strategy is built in the solver. More... | |
void | onOffStepControl (bool sc) |
Disables or enables automatic step control. More... | |
void | setMaxReturnTime (double maxReturnTime) |
Sets maximal return time to Poincare section. If the trajectory does not reach the section during time [0,maxReturnTime] then an exception is thrown. This prevents looping of the procedure computing Poincare map in the case when the trajectory is captured by an attractor (like sink or attracting periodic orbit) and never intersect Poincare section. More... | |
void | setBlowUpMaxNorm (double blowUpMaxNorm) |
Sets threshold value of norm that is considered as a blow up of the solution. A trajectory may escape to infinity in finite time without crossing Poincare section. If the norm of a point on trajectory is bigger than this threshold value then an exception is thrown. More... | |
Public Attributes | |
BoundType | timeStepCorrectionFactor |
constant timeStepCorrectionFactor is used as a correction factor to multiply the time step in procedure reachSection when coming close to section and actual time step is to large. The value 0.9 is just from numerical simulations as a good working. It can be changed to a value from the interval (0,1) More... | |
int | maxCorrectionAttempts |
maximal number of correction attempts More... | |
ScalarType | minCrossingTimeStep |
minimal time step during section crossing More... | |
Protected Member Functions | |
template<class T > | |
VectorType | computePoincareMap (T &originalSet, int n) |
template<typename T > | |
void | checkTransversability (const T &theSet, const VectorType &bounds) |
Function checks if we crossed section and then returned in one step. In this case it throws an exception of PoincareException<T> type. More... | |
template<typename T > | |
ScalarType | getSign (const T &theSet, VectorType &position, bool updatePosition, const VectorType &bounds) |
Function returns section sign. It also checks If we crossed section and returned in one step. In this case it throws an exception. More... | |
template<typename T > | |
ScalarType | getSign (const T &theSet, VectorType &position, bool updatePosition) |
Function returns section sign. It also checks for possible crossing the section and returning in one step. In this case it throws an exception. More... | |
template<typename T > | |
ScalarType | getSign (const T &theSet) |
Function returns section sign. It also checks If we crossed section and returned in one step. In this case it throws an exception. More... | |
template<typename T > | |
void | integrateUntilSectionCrossing (T &theSet, T &nextSet, int n=1) |
The function propagates the set of initial conditions by the flow until its -th intersection with the Poincare section. More... | |
template<typename T > | |
void | getCloserToSection (T &theSet, T &nextSet) |
On input: theSet is 'one step' before the section On output theSet is before the section and closer than sizeFactor*diam(theSet) in perpendicular direction to the section. More... | |
template<typename T > | |
VectorType | crossSection (T &theSet, T &nextSet) |
Crosses the section and returns the value of Poincare Map. It updates also derivatives. We expect theSet to be just before the section. After this function theSet is on the section or just after the section. This function is common for all types of C^0, C^1, C^2 and C^n sets. More... | |
template<typename T > | |
bool | crossSectionInOneStep (T &theSet, T &nextSet, ScalarType &oneStepReturnTime, VectorType &bound) |
If it is possible to cross the section in one step, the function does this and returns bounf Computes return time by means of the Newton method. If not possible, returns false. More... | |
T | integrateUntilSectionCrossing (T &x, int n=1, ScalarType *lastStep=0) |
Given an initial condition the function computes its trajetory until -th intersection with the Poincare section. The point just before -th intersection is returned. More... | |
ScalarType | findRelativeCrossingTime (ScalarType timeJustBeforeSection, const CurveType &curve) |
This function computes time needed to cross the section from a point which is just before the section. Given solves for satisfying by a Newton like scheme. This is just scalar equation with one unknown. Then, computed return time is used to integrate the point and compute Poincare map. More... | |
Protected Attributes | |
MatrixType | sectionCoordinateSystem |
SectionDerivativesEnclosureType | sectionDerivativesEnclosure |
ScalarType | m_signBeforeSection |
ScalarType | m_lastTimeStep |
Solver & | m_solver |
dynamical system (e.g. Taylor, C2Taylor, CnTaylor, FadTaylor, etc) More... | |
const SectionType * | m_section |
section function More... | |
CrossingDirection | m_crossingDirection |
requested direction of section crossing More... | |
RealType | m_sectionFactor |
sectionFactor is used in the procedures reachSection and crossSection. IN NONRIGOROUS : we want to be closer than m_sectionFactor after section crossing IN RIGOROUS VERSION: Before and after section crossing we want the set to be closer to section than sectionFactor * sizeOfSet More... | |
double | m_blowUpMaxNorm |
double | m_maxReturnTime |
PoicareMap class rigorously computes Poincare Map.
For given Dynamical System and section it rigorously computes first return map to section (Poincare Map) and its derivatives.
It uses BasicPoincareMap and adds error estimates to the results.
Let
In the following we denote by
Parameters:
typedef dynset::AffineCoordinateChange<MatrixType> capd::poincare::PoincareMap< SolverT, SectionT >::AffineCoordinateChange |
typedef ScalarType::BoundType capd::poincare::PoincareMap< SolverT, SectionT >::BoundType |
a type of the end for intervals
typedef capd::dynset::C0Set<MatrixType> capd::poincare::PoincareMap< SolverT, SectionT >::C0Set |
type of abstract base class for all C0 sets
typedef capd::dynset::C1Set<MatrixType> capd::poincare::PoincareMap< SolverT, SectionT >::C1Set |
type of abstract base class for all C1 sets
typedef capd::dynset::C2Set<MatrixType> capd::poincare::PoincareMap< SolverT, SectionT >::C2Set |
type of abstract base class for all C2 sets
typedef capd::dynset::CnSet<MatrixType,0> capd::poincare::PoincareMap< SolverT, SectionT >::CnSet |
type of abstract base class for all C2 sets
typedef BasicPoincareMap<Solver,SectionT>::CrossingDirection capd::poincare::PoincareMap< SolverT, SectionT >::CrossingDirection |
section crossing direction; possible values are PlusMinus, Both(default) and MinusPlus.
typedef Solver::SolutionCurve capd::poincare::PoincareMap< SolverT, SectionT >::CurveType |
typedef C2Set::HessianType capd::poincare::PoincareMap< SolverT, SectionT >::HessianType |
type of data structure for storing of Hessians of multivariate mappings.
typedef CnSet::JetType capd::poincare::PoincareMap< SolverT, SectionT >::JetType |
typedef Solver::MatrixType capd::poincare::PoincareMap< SolverT, SectionT >::MatrixType |
|
inherited |
floating point type: for floating points is equal to ScalarType, for intervals is a type of their bounds.
typedef Solver::ScalarType capd::poincare::PoincareMap< SolverT, SectionT >::ScalarType |
typedef SectionType::SectionDerivativesEnclosureType capd::poincare::PoincareMap< SolverT, SectionT >::SectionDerivativesEnclosureType |
typedef SectionT capd::poincare::PoincareMap< SolverT, SectionT >::SectionType |
type of function that represents Poincare section.
typedef VectorType::size_type capd::poincare::PoincareMap< SolverT, SectionT >::size_type |
integral type used to index containers (vectors, matrices, etc)
|
inherited |
typedef SolverT capd::poincare::PoincareMap< SolverT, SectionT >::Solver |
ODE solver.
typedef Solver::VectorFieldType capd::poincare::PoincareMap< SolverT, SectionT >::VectorFieldType |
vector field type
typedef Solver::VectorType capd::poincare::PoincareMap< SolverT, SectionT >::VectorType |
|
inherited |
Recomputes Taylor expansion of the flow into Taylor expansion of Poincare map.
Computes all partial derivatives of a Poincare map to an arbitrarily order.
[in] | Px | partial derivatives of the flow to a given order |
|
inlineinherited |
Simultaneous computation of first and second Taylor coefficients of return time and Poincare map.
[in] | Px | - value of Poincare map |
[in] | derivativeOfFlow | - solution to first variational equation computed at return time |
[in] | hessianOfFlow | - solution to first variational equation computed at return time |
[out] | DP | - computed derivative of Poincare map |
[out] | D2P | - computed second order Taylor coefficients of Poincare map |
[in] | returnTime | - return time to the section |
|
inherited |
Simultaneous computation of first and second Taylor coefficients of return time and Poincare map.
[in] | Px | - value of Poincare map |
[in] | derivativeOfFlow | - solution to first variational equation computed at return time |
[in] | hessianOfFlow | - solution to first variational equation computed at return time |
[out] | DP | - computed derivative of Poincare map |
[out] | D2P | - computed second order Taylor coefficients of Poincare map |
[out] | dT | - computed gradient of return time |
[out] | d2T | - computed second order Taylor coefficients of return time |
[in] | returnTime | - return time to the section |
|
inlineinherited |
Computes derivative of Poincare Map dP.
[in] | Px | - value of Poincare map |
[in] | derivativeOfFlow | - solution to first variational equation computed at return time |
|
inlineinherited |
Simultaneous computation of gradient of return time and derivative of Poincare Map dP.
[in] | Px | - value of Poincare map |
[in] | derivativeOfFlow | - solution to first order variational equation computed at return time |
[out] | dT | - computed gradient of return time |
[in] | returnTime | - computed return time |
|
protectedinherited |
This function computes time needed to cross the section from a point which is just before the section. Given solves for satisfying by a Newton like scheme. This is just scalar equation with one unknown. Then, computed return time is used to integrate the point and compute Poincare map.
Crosses the section and returns the value of Poincare Map.
|
inlineinherited |
Returns reference to solver used to integrate the system.
|
inlineinherited |
Returns read-only reference to solver used to integrate the system.
|
inlineinherited |
Returns order of numerical method used.
|
inlineinherited |
Returns reference to Poincare section object.
|
inline |
|
inlineinherited |
Returns reference to solver used to integrate the system.
|
inlineinherited |
Returns read-only reference to solver used to integrate the system.
|
inlineinherited |
Returns time step.
|
inlineinherited |
Returns reference to solver used to integrate the system.
|
inlineinherited |
Returns read-only reference to solver used to integrate the system.
|
protectedinherited |
Given an initial condition the function computes its trajetory until -th intersection with the Poincare section. The point just before -th intersection is returned.
////////////////////////////////////////////////////////////////////////////////
This function moves theSet with the given flow and stops just before the section (for n-th iterate of Poincare Map).
Parameter theSet contains on return the set just before the section.
Common function for different set types. //////////////////////////////////////////////////////////////////////////////
@TODO : implement for n>1
theSet | ||
[in,out] | theSet | on input: initial set position, on return: the set just before the section. |
n | number of iterates | |
lastStep | if given |
|
inlineinherited |
Disables or enables automatic step control.
|
inherited |
Computes value of Poincare Map.
[in] | v | - argument of Poincare map |
|
inherited |
Computes value of Poincare Map and derivativeof the flow.
[in] | v | - vector at which we want compute Poincare map |
[out] | dF | - derivative of the flow evaluated at (v,returnTime). |
|
inherited |
Computes value of Poincare Map, derivative and hessian of the flow.
[in] | v | - vector at which we want compute Poincare map |
[out] | dF | - derivative of the flow evaluated at (v,returnTime). |
[out] | h | - second order Taylor coefficients of flow evaluated at (v,returnTime). |
|
inherited |
Computes value of Poincare Map.
[in] | v | - argument of Poincare map |
[out] | afterSection | - on output it contains a point just after first intersection with Poincare section that might be used for further integration. |
|
inherited |
Computes Poincare Map and derivatives of the flow to given order evaluated at return time.
Nonrigorous Poincare map.
This routine is valid for an autonomous flows, only.
A point just after the section on the nonrigorous trajectory is returned The result contains also approximate values of the partial derivatives of the flow
|
inherited |
Computes Poincare Map and derivatives of the flow to given order evaluated at return time.
Nonrigorous Poincare map.
This routine is valid for both autonomous and time-dependent vector fields.
A point just after the section on the nonrigorous trajectory is returned The result contains also approximate values of the partial derivatives of the flow
|
inherited |
Computes value of Poincare Map, derivative and hessian of the flow.
[in] | v | - vector at which we want compute Poincare map |
[out] | dF | - derivative of the flow evaluated at (v,returnTime). |
[out] | h | - second order Taylor coefficients of flow evaluated at (v,returnTime). |
[in,out] | in_out_time | - on input it should specify initial time for integration for nonautonomous flows. For autonomous systems it should be initialized by zero. On output it contains approximate time at which intersection with Poincare section occurred. |
|
inherited |
Computes value of Poincare Map and derivative of the flow.
[in] | v | - vector at which we want compute Poincare map |
[out] | dF | - derivative of the flow evaluated at (v,returnTime). |
[in,out] | in_out_time | - on input it should specify initial time for integration for nonautonomous flows. For autonomous systems it should be initialized by zero. On output it contains approximate time at which intersection with Poincare section occurred. |
|
inherited |
Computes value of Poincare Map.
[in] | v | - argument of Poincare map |
[in,out] | in_out_time | - on input it should specify initial time for integration for nonautonomous flows. For autonomous systems it should be initialized by zero. On output it contains approximate time at which intersection with Poincare section occurred. |
|
inherited |
Computes value of Poincare Map.
[in] | v | - argument of Poincare map |
[out] | afterSection | - on output it contains a point just after first intersection with Poincare section that might be used for further integration. |
[in,out] | in_out_time | - on input it should specify initial time for integration for nonautonomous flows. For autonomous systems it should be initialized by zero. On output it contains approximate time at which intersection with Poincare section occurred. |
|
inlineinherited |
Sets threshold value of norm that is considered as a blow up of the solution. A trajectory may escape to infinity in finite time without crossing Poincare section. If the norm of a point on trajectory is bigger than this threshold value then an exception is thrown.
Default value is blowUpMaxNorm = 1e+10.
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
Sets maximal return time to Poincare section. If the trajectory does not reach the section during time [0,maxReturnTime] then an exception is thrown. This prevents looping of the procedure computing Poincare map in the case when the trajectory is captured by an attractor (like sink or attracting periodic orbit) and never intersect Poincare section.
Default value is maxReturnTime = 1000.
|
inlineinherited |
Sets order of Taylor method.
|
inlineinherited |
Sets new section function.
|
inlineinherited |
Sets time step for integration of DS.
|
inlineinherited |
Enables automatic step control. Step control strategy is built in the solver.
|
inlineinherited |
Disables automatic step control.
|
protectedinherited |
|
protectedinherited |
requested direction of section crossing
|
protected |
|
protectedinherited |
|
protectedinherited |
section function
|
protectedinherited |
sectionFactor is used in the procedures reachSection and crossSection. IN NONRIGOROUS : we want to be closer than m_sectionFactor after section crossing IN RIGOROUS VERSION: Before and after section crossing we want the set to be closer to section than sectionFactor * sizeOfSet
|
protected |
|
protectedinherited |
dynamical system (e.g. Taylor, C2Taylor, CnTaylor, FadTaylor, etc)
int capd::poincare::PoincareMap< SolverT, SectionT >::maxCorrectionAttempts |
maximal number of correction attempts
ScalarType capd::poincare::PoincareMap< SolverT, SectionT >::minCrossingTimeStep |
minimal time step during section crossing
|
protected |
|
protected |
BoundType capd::poincare::PoincareMap< SolverT, SectionT >::timeStepCorrectionFactor |
constant timeStepCorrectionFactor is used as a correction factor to multiply the time step in procedure reachSection when coming close to section and actual time step is to large. The value 0.9 is just from numerical simulations as a good working. It can be changed to a value from the interval (0,1)