CAPD::DynSys Library  6.0.0
capd::poincare::PoincareMap< SolverT, SectionT > Class Template Reference

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< MatrixTypeAffineCoordinateChange
 
typedef SectionT SectionType
 type of function $ s: R^n \to R $ that represents Poincare section. More...
 
typedef SectionType::SectionDerivativesEnclosureType SectionDerivativesEnclosureType
 
typedef capd::dynset::C0Set< MatrixTypeC0Set
 type of abstract base class for all C0 sets More...
 
typedef capd::dynset::C1Set< MatrixTypeC1Set
 type of abstract base class for all C1 sets More...
 
typedef capd::dynset::C2Set< MatrixTypeC2Set
 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< CurveTypeSolutionCurve
 

Public Member Functions

 PoincareMap (Solver &solver, SectionType &section, 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 SectionDerivativesEnclosureTypegetSectionDerivativesEnclosure () 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 SolvergetSolver () const
 Returns read-only reference to solver used to integrate the system. More...
 
SolvergetSolver ()
 Returns reference to solver used to integrate the system. More...
 
const SolvergetDynamicalSystem () const
 Returns read-only reference to solver used to integrate the system. More...
 
SolvergetDynamicalSystem ()
 Returns reference to solver used to integrate the system. More...
 
VectorFieldTypegetVectorField ()
 Returns reference to solver used to integrate the system. More...
 
const VectorFieldTypegetVectorField () const
 Returns read-only reference to solver used to integrate the system. More...
 
const SectionTypegetSection () 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 $ L_1$ 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 $ L_1 $ 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 $ theSet $ by the flow until its $n$-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...
 
integrateUntilSectionCrossing (T &x, int n=1, ScalarType *lastStep=0)
 Given an initial condition $ x $ the function computes its trajetory until $n$-th intersection with the Poincare section. The point just before $n$-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 $x$ solves for $t$ satisfying $\alpha(\varphi(t,x))=0$ by a Newton like scheme. This is just scalar equation with one unknown. Then, computed return time is used to integrate the point $x$ and compute Poincare map. More...
 

Protected Attributes

MatrixType sectionCoordinateSystem
 
SectionDerivativesEnclosureType sectionDerivativesEnclosure
 
ScalarType m_signBeforeSection
 
ScalarType m_lastTimeStep
 
Solverm_solver
 dynamical system (e.g. Taylor, C2Taylor, CnTaylor, FadTaylor, etc) More...
 
const SectionTypem_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
 

Detailed Description

template<typename SolverT, typename SectionT = AbstractSection<typename SolverT::MatrixType>>
class capd::poincare::PoincareMap< SolverT, SectionT >

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

  • $ \varphi(t, x): R \times R^n -> R^n $ be a dynamical system generated by given vector field.
  • $ s: R^n \to R $ be a section function,
  • $ S = \{x \in R^n : s(x) = 0\} $
  • $ P: S \to S $ be a Poincare Map
  • for given point $ x \in S $ let T(x) be first return time (in given crossing direction) i.e. $ P(x) = \varphi(T(x), x) \in S $

In the following we denote by

  • dP the derivative of Poincare Map P : $ dP(x) = \frac{\partial P(x)}{\partial x} $
  • dT the derivative of T(x) : $ dT(x) = \frac{\partial T(x)}{\partial x} $
  • dF the derivative of the flow : $ dF(x) = \frac{\partial \varphi(T(x), x)}{\partial x} $ Then $ dP = dF + \frac{\partial \varphi}{\partial t} dT $

Parameters:

  • DS dynamical system
    • DS::MapType vector field type
  • SectionT scalar function describing section :

Member Typedef Documentation

◆ AffineCoordinateChange

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef dynset::AffineCoordinateChange<MatrixType> capd::poincare::PoincareMap< SolverT, SectionT >::AffineCoordinateChange

◆ BoundType

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef ScalarType::BoundType capd::poincare::PoincareMap< SolverT, SectionT >::BoundType

a type of the end for intervals

◆ C0Set

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef capd::dynset::C0Set<MatrixType> capd::poincare::PoincareMap< SolverT, SectionT >::C0Set

type of abstract base class for all C0 sets

◆ C1Set

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef capd::dynset::C1Set<MatrixType> capd::poincare::PoincareMap< SolverT, SectionT >::C1Set

type of abstract base class for all C1 sets

◆ C2Set

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef capd::dynset::C2Set<MatrixType> capd::poincare::PoincareMap< SolverT, SectionT >::C2Set

type of abstract base class for all C2 sets

◆ CnSet

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef capd::dynset::CnSet<MatrixType,0> capd::poincare::PoincareMap< SolverT, SectionT >::CnSet

type of abstract base class for all C2 sets

◆ CrossingDirection

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef BasicPoincareMap<Solver,SectionT>::CrossingDirection capd::poincare::PoincareMap< SolverT, SectionT >::CrossingDirection

section crossing direction; possible values are PlusMinus, Both(default) and MinusPlus.

◆ CurveType

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef Solver::SolutionCurve capd::poincare::PoincareMap< SolverT, SectionT >::CurveType

◆ HessianType

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef C2Set::HessianType capd::poincare::PoincareMap< SolverT, SectionT >::HessianType

type of data structure for storing of Hessians of multivariate mappings.

◆ JetType

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef CnSet::JetType capd::poincare::PoincareMap< SolverT, SectionT >::JetType

◆ MatrixType

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef Solver::MatrixType capd::poincare::PoincareMap< SolverT, SectionT >::MatrixType

◆ RealType

typedef TypeTraits<ScalarType>::Real capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::RealType
inherited

floating point type: for floating points is equal to ScalarType, for intervals is a type of their bounds.

◆ ScalarType

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef Solver::ScalarType capd::poincare::PoincareMap< SolverT, SectionT >::ScalarType

◆ SectionDerivativesEnclosureType

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef SectionType::SectionDerivativesEnclosureType capd::poincare::PoincareMap< SolverT, SectionT >::SectionDerivativesEnclosureType

◆ SectionType

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef SectionT capd::poincare::PoincareMap< SolverT, SectionT >::SectionType

type of function $ s: R^n \to R $ that represents Poincare section.

◆ size_type

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef VectorType::size_type capd::poincare::PoincareMap< SolverT, SectionT >::size_type

integral type used to index containers (vectors, matrices, etc)

◆ SolutionCurve

typedef capd::diffAlgebra::SolutionCurve<CurveType> capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::SolutionCurve
inherited

◆ Solver

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef SolverT capd::poincare::PoincareMap< SolverT, SectionT >::Solver

ODE solver.

◆ VectorFieldType

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef Solver::VectorFieldType capd::poincare::PoincareMap< SolverT, SectionT >::VectorFieldType

vector field type

◆ VectorType

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
typedef Solver::VectorType capd::poincare::PoincareMap< SolverT, SectionT >::VectorType

Member Function Documentation

◆ computeDP() [1/5]

JetT capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::computeDP ( const JetT &  PhiSeries)
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.

Parameters
[in]Pxpartial derivatives of the flow to a given order
Remarks
Computations are valid only for linear sections.

◆ computeDP() [2/5]

void capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::computeDP ( const VectorType Px,
const MatrixType derivativeOfFlow,
const HessianType hessianOfFlow,
MatrixType DP,
HessianType D2P,
ScalarType  returnTime = TypeTraits<ScalarType>::zero() 
)
inlineinherited

Simultaneous computation of first and second Taylor coefficients of return time and Poincare map.

Parameters
[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
Note
all input and output parameters are Taylor coefficients, not derivatives!
returnTime must be specified for nonautonomous flows only. Otherwise, default value 0 is valid.

◆ computeDP() [3/5]

void capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::computeDP ( const VectorType Px,
const MatrixType derivativeOfFlow,
const HessianType hessianOfFlow,
MatrixType DP,
HessianType D2P,
VectorType dT,
MatrixType d2T,
ScalarType  returnTime = TypeTraits<ScalarType>::zero() 
)
inherited

Simultaneous computation of first and second Taylor coefficients of return time and Poincare map.

Parameters
[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
Note
all input and output parameters are Taylor coefficients, not derivatives!
returnTime must be specified for nonautonomous flows only. Otherwise, default value 0 is valid.

◆ computeDP() [4/5]

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::MatrixType capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::computeDP ( const VectorType Px,
const MatrixType derivativeOfFlow,
ScalarType  returnTime = TypeTraits<ScalarType>::zero() 
)
inlineinherited

Computes derivative of Poincare Map dP.

Parameters
[in]Px- value of Poincare map
[in]derivativeOfFlow- solution to first variational equation computed at return time
Note
returnTime must be specified for nonautonomous flows only. Otherwise, default value 0 is valid.

◆ computeDP() [5/5]

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::MatrixType capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::computeDP ( const VectorType Px,
const MatrixType derivativeOfFlow,
VectorType dT,
ScalarType  returnTime = TypeTraits<ScalarType>::zero() 
)
inlineinherited

Simultaneous computation of gradient of return time and derivative of Poincare Map dP.

Parameters
[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
Note
returnTime must be specified for nonautonomous flows only. Otherwise, default value 0 is valid.

◆ findRelativeCrossingTime()

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::ScalarType capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::findRelativeCrossingTime ( ScalarType  timeJustBeforeSection,
const CurveType curve 
)
protectedinherited

This function computes time needed to cross the section from a point which is just before the section. Given $x$ solves for $t$ satisfying $\alpha(\varphi(t,x))=0$ by a Newton like scheme. This is just scalar equation with one unknown. Then, computed return time is used to integrate the point $x$ and compute Poincare map.

Crosses the section and returns the value of Poincare Map.

◆ getDynamicalSystem() [1/2]

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::Solver & capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::getDynamicalSystem
inlineinherited

Returns reference to solver used to integrate the system.

◆ getDynamicalSystem() [2/2]

const BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::Solver & capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::getDynamicalSystem
inlineinherited

Returns read-only reference to solver used to integrate the system.

◆ getOrder()

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::size_type capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::getOrder
inlineinherited

Returns order of numerical method used.

◆ getSection()

const BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::SectionType & capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::getSection
inlineinherited

Returns reference to Poincare section object.

◆ getSectionDerivativesEnclosure()

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
const SectionDerivativesEnclosureType& capd::poincare::PoincareMap< SolverT, SectionT >::getSectionDerivativesEnclosure ( ) const
inline

◆ getSolver() [1/2]

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::Solver & capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::getSolver
inlineinherited

Returns reference to solver used to integrate the system.

◆ getSolver() [2/2]

const BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::Solver & capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::getSolver
inlineinherited

Returns read-only reference to solver used to integrate the system.

◆ getStep()

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::ScalarType capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::getStep
inlineinherited

Returns time step.

◆ getVectorField() [1/2]

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::VectorFieldType & capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::getVectorField
inlineinherited

Returns reference to solver used to integrate the system.

◆ getVectorField() [2/2]

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::VectorFieldType const & capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::getVectorField
inlineinherited

Returns read-only reference to solver used to integrate the system.

◆ integrateUntilSectionCrossing()

T capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::integrateUntilSectionCrossing ( T &  theSet,
int  n = 1,
ScalarType lastStep = 0 
)
protectedinherited

Given an initial condition $ x $ the function computes its trajetory until $n$-th intersection with the Poincare section. The point just before $n$-th intersection is returned.

////////////////////////////////////////////////////////////////////////////////

Note
: An exception is thrown if either the maximal accepted return time is exceeded or the norm of a point on the trajectory exceeds maximal accpeted value. This may occur when the trajectory escapes to infinity.
: An initial condition is a template parameter. It can be just a vector or a structure holding also intitial conditions for (higher order) variational equations.

This function moves theSet with the given flow and stops just before the section (for n-th iterate of Poincare Map).

Returns
set on the section or just after the section (suitable for a computation of next Poincare Map iteration).

Parameter theSet contains on return the set just before the section.

Common function for different set types. //////////////////////////////////////////////////////////////////////////////

@TODO : implement for n>1

Parameters
theSet
[in,out]theSeton input: initial set position, on return: the set just before the section.
nnumber of iterates
lastStepif given

◆ onOffStepControl()

void capd::poincare::BasicPoincareMap< SolverT, FunT >::onOffStepControl ( bool  sc)
inlineinherited

Disables or enables automatic step control.

◆ operator()() [1/10]

BasicPoincareMap< SolverT, FunT >::VectorType capd::poincare::BasicPoincareMap< SolverT, FunT >::operator() ( const VectorType v)
inherited

Computes value of Poincare Map.

Parameters
[in]v- argument of Poincare map
Returns
P(v) - first intersection of the trajectory of v with Poincare section.

◆ operator()() [2/10]

BasicPoincareMap< SolverT, FunT >::VectorType capd::poincare::BasicPoincareMap< SolverT, FunT >::operator() ( const VectorType v,
MatrixType dF 
)
inherited

Computes value of Poincare Map and derivativeof the flow.

Parameters
[in]v- vector at which we want compute Poincare map
[out]dF- derivative of the flow evaluated at (v,returnTime).
Note
This routine is valid for an autonomous flows, only.

◆ operator()() [3/10]

BasicPoincareMap< SolverT, FunT >::VectorType capd::poincare::BasicPoincareMap< SolverT, FunT >::operator() ( const VectorType v,
MatrixType dF,
HessianType h 
)
inherited

Computes value of Poincare Map, derivative and hessian of the flow.

Parameters
[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).
Note
This routine is valid for an autonomous flows, only.

◆ operator()() [4/10]

VectorType capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::operator() ( const VectorType v,
VectorType afterSection 
)
inherited

Computes value of Poincare Map.

Parameters
[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.
Returns
P(v) - first intersection of the trajectory of v with Poincare section.

◆ operator()() [5/10]

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::VectorType capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::operator() ( JetType x)
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

◆ operator()() [6/10]

BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::VectorType capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::operator() ( JetType x,
ScalarType in_out_time 
)
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

◆ operator()() [7/10]

BasicPoincareMap< SolverT, FunT >::VectorType capd::poincare::BasicPoincareMap< SolverT, FunT >::operator() ( VectorType  v,
MatrixType dF,
HessianType h,
ScalarType in_out_time 
)
inherited

Computes value of Poincare Map, derivative and hessian of the flow.

Parameters
[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.

◆ operator()() [8/10]

BasicPoincareMap< SolverT, FunT >::VectorType capd::poincare::BasicPoincareMap< SolverT, FunT >::operator() ( VectorType  v,
MatrixType dF,
ScalarType in_out_time 
)
inherited

Computes value of Poincare Map and derivative of the flow.

Parameters
[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.

◆ operator()() [9/10]

BasicPoincareMap< SolverT, FunT >::VectorType capd::poincare::BasicPoincareMap< SolverT, FunT >::operator() ( VectorType  v,
ScalarType in_out_time 
)
inherited

Computes value of Poincare Map.

Parameters
[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.
Returns
P(v) - first intersection of the trajectory of v with Poincare section.

◆ operator()() [10/10]

BasicPoincareMap< SolverT, FunT >::VectorType capd::poincare::BasicPoincareMap< SolverT, FunT >::operator() ( VectorType  v,
VectorType afterSection,
ScalarType in_out_time 
)
inherited

Computes value of Poincare Map.

Parameters
[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.
Returns
P(v) - first intersection of the trajectory of v with Poincare section.

◆ setBlowUpMaxNorm()

void capd::poincare::BasicPoincareMap< SolverT, FunT >::setBlowUpMaxNorm ( double  blowUpMaxNorm)
inlineinherited

Sets threshold value of $ L_1$ 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 $ L_1 $ norm of a point on trajectory is bigger than this threshold value then an exception is thrown.

Default value is blowUpMaxNorm = 1e+10.

◆ setCrossingDirection()

void capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::setCrossingDirection ( CrossingDirection  cs)
inlineinherited

◆ setFactor()

void capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::setFactor ( double  newFactor)
inlineinherited

◆ setMaxReturnTime()

void capd::poincare::BasicPoincareMap< SolverT, FunT >::setMaxReturnTime ( double  maxReturnTime)
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.

◆ setOrder()

void capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::setOrder ( size_type  newOrder)
inlineinherited

Sets order of Taylor method.

◆ setSection()

void capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::setSection ( const SectionType newSection)
inlineinherited

Sets new section function.

◆ setStep()

void capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::setStep ( const ScalarType newStep)
inlineinherited

Sets time step for integration of DS.

◆ turnOffStepControl()

void capd::poincare::BasicPoincareMap< SolverT, FunT >::turnOffStepControl
inlineinherited

Enables automatic step control. Step control strategy is built in the solver.

◆ turnOnStepControl()

void capd::poincare::BasicPoincareMap< SolverT, FunT >::turnOnStepControl
inlineinherited

Disables automatic step control.

Member Data Documentation

◆ m_blowUpMaxNorm

double capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::m_blowUpMaxNorm
protectedinherited

◆ m_crossingDirection

CrossingDirection capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::m_crossingDirection
protectedinherited

requested direction of section crossing

◆ m_lastTimeStep

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
ScalarType capd::poincare::PoincareMap< SolverT, SectionT >::m_lastTimeStep
protected

◆ m_maxReturnTime

double capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::m_maxReturnTime
protectedinherited

◆ m_section

const SectionType* capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::m_section
protectedinherited

section function

◆ m_sectionFactor

RealType capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::m_sectionFactor
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

◆ m_signBeforeSection

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
ScalarType capd::poincare::PoincareMap< SolverT, SectionT >::m_signBeforeSection
protected

◆ m_solver

Solver& capd::poincare::BasicPoincareMap< SolverT, AbstractSection< typename SolverT::MatrixType > >::m_solver
protectedinherited

dynamical system (e.g. Taylor, C2Taylor, CnTaylor, FadTaylor, etc)

◆ maxCorrectionAttempts

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
int capd::poincare::PoincareMap< SolverT, SectionT >::maxCorrectionAttempts

maximal number of correction attempts

◆ minCrossingTimeStep

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
ScalarType capd::poincare::PoincareMap< SolverT, SectionT >::minCrossingTimeStep

minimal time step during section crossing

◆ sectionCoordinateSystem

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
MatrixType capd::poincare::PoincareMap< SolverT, SectionT >::sectionCoordinateSystem
protected

◆ sectionDerivativesEnclosure

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
SectionDerivativesEnclosureType capd::poincare::PoincareMap< SolverT, SectionT >::sectionDerivativesEnclosure
protected

◆ timeStepCorrectionFactor

template<typename SolverT , typename SectionT = AbstractSection<typename SolverT::MatrixType>>
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)


The documentation for this class was generated from the following files: