Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members

Abstract base class for all kinds of algorithms for integrating differential equations (ODEs or DAEs). More...

#include <integrator.hpp>

Inheritance diagram for Integrator:
Inheritance graph
[legend]

Public Member Functions

virtual BooleanType canHandleImplicitSwitches () const
 
virtual Integratorclone () const =0
 
virtual returnValue deleteAllSeeds ()
 
virtual returnValue evaluateSensitivities ()=0
 
virtual returnValue freezeAll ()=0
 
virtual returnValue freezeMesh ()=0
 
returnValue getBackwardSensitivities (DVector &Dx_x0, DVector &Dx_p, DVector &Dx_u, DVector &Dx_w, int order) const
 
virtual double getDifferentialEquationSampleTime () const
 
returnValue getForwardSensitivities (DVector &Dx, int order) const
 
returnValue getForwardSensitivities (VariablesGrid &Dx, int order) const
 
returnValue getI (VariablesGrid &I) const
 
virtual int getNumberOfRejectedSteps () const =0
 
virtual int getNumberOfSteps () const =0
 
virtual double getStepSize () const =0
 
returnValue getX (DVector &xEnd) const
 
returnValue getX (VariablesGrid &X) const
 
returnValue getXA (DVector &xaEnd) const
 
returnValue getXA (VariablesGrid &XA) const
 
virtual returnValue init (const DifferentialEquation &rhs)=0
 
returnValue init (const DifferentialEquation &rhs, const Transition &trs)
 
returnValue integrate (double t0, double tend, double *x0, double *xa=0, double *p=0, double *u=0, double *w=0)
 
returnValue integrate (const Grid &t, double *x0, double *xa=0, double *p=0, double *u=0, double *w=0)
 
returnValue integrate (double t0, double tend, const DVector &x0, const DVector &xa=emptyVector, const DVector &p=emptyVector, const DVector &u=emptyVector, const DVector &w=emptyVector)
 
returnValue integrate (const Grid &t, const DVector &x0, const DVector &xa=emptyVector, const DVector &p=emptyVector, const DVector &u=emptyVector, const DVector &w=emptyVector)
 
returnValue integrateSensitivities ()
 
 Integrator ()
 
 Integrator (const Integrator &arg)
 
virtual BooleanType isDifferentialEquationAffine () const
 
virtual BooleanType isDifferentialEquationDefined () const
 
virtual returnValue printRunTimeProfile () const
 
returnValue setBackwardSeed (const int &order, const DVector &seed)
 
virtual returnValue setDxInitialization (double *dx0)=0
 
returnValue setForwardSeed (const int &order, const DVector &xSeed, const DVector &pSeed=emptyVector, const DVector &uSeed=emptyVector, const DVector &wSeed=emptyVector)
 
returnValue setTransition (const Transition &trs)
 
virtual returnValue unfreeze ()=0
 
virtual ~Integrator ()
 
- Public Member Functions inherited from AlgorithmicBase
int addLogRecord (LogRecord &_record)
 
returnValue addOption (OptionsName name, int value)
 
returnValue addOption (OptionsName name, double value)
 
returnValue addOption (uint idx, OptionsName name, int value)
 
returnValue addOption (uint idx, OptionsName name, double value)
 
returnValue addOptionsList ()
 
 AlgorithmicBase ()
 
 AlgorithmicBase (UserInteraction *_userInteraction)
 
 AlgorithmicBase (const AlgorithmicBase &rhs)
 
returnValue get (OptionsName name, int &value) const
 
returnValue get (OptionsName name, double &value) const
 
returnValue get (OptionsName name, std::string &value) const
 
returnValue get (uint idx, OptionsName name, int &value) const
 
returnValue get (uint idx, OptionsName name, double &value) const
 
returnValue getAll (LogName _name, MatrixVariablesGrid &values) const
 
returnValue getFirst (LogName _name, DMatrix &firstValue) const
 
returnValue getFirst (LogName _name, VariablesGrid &firstValue) const
 
returnValue getLast (LogName _name, DMatrix &lastValue) const
 
returnValue getLast (LogName _name, VariablesGrid &lastValue) const
 
Options getOptions (uint idx) const
 
BooleanType haveOptionsChanged () const
 
BooleanType haveOptionsChanged (uint idx) const
 
AlgorithmicBaseoperator= (const AlgorithmicBase &rhs)
 
returnValue plot (PlotFrequency _frequency=PLOT_IN_ANY_CASE)
 
returnValue printLogRecord (std::ostream &_stream, int idx, LogPrintMode _mode=PRINT_ITEM_BY_ITEM) const
 
returnValue replot (PlotFrequency _frequency=PLOT_IN_ANY_CASE)
 
returnValue set (OptionsName name, int value)
 
returnValue set (OptionsName name, double value)
 
returnValue set (OptionsName name, const std::string &value)
 
returnValue set (uint idx, OptionsName name, int value)
 
returnValue set (uint idx, OptionsName name, double value)
 
returnValue setAll (LogName _name, const MatrixVariablesGrid &values)
 
returnValue setLast (LogName _name, int lastValue, double time=-INFTY)
 
returnValue setLast (LogName _name, double lastValue, double time=-INFTY)
 
returnValue setLast (LogName _name, const DVector &lastValue, double time=-INFTY)
 
returnValue setLast (LogName _name, const DMatrix &lastValue, double time=-INFTY)
 
returnValue setLast (LogName _name, const VariablesGrid &lastValue, double time=-INFTY)
 
returnValue setOptions (const Options &arg)
 
returnValue setOptions (uint idx, const Options &arg)
 
virtual ~AlgorithmicBase ()
 

Protected Member Functions

virtual returnValue diffTransitionBackward (DVector &DX, DVector &DP, DVector &DU, DVector &DW, int &order)
 
virtual returnValue diffTransitionForward (DVector &DX, const DVector &DP, const DVector &DU, const DVector &DW, const int &order)
 
virtual returnValue evaluate (const DVector &x0, const DVector &xa, const DVector &p, const DVector &u, const DVector &w, const Grid &t_)=0
 
virtual returnValue evaluateTransition (const double time, DVector &xd, const DVector &xa, const DVector &p, const DVector &u, const DVector &w)
 
virtual int getDim () const =0
 
virtual int getDimX () const
 
virtual returnValue getProtectedBackwardSensitivities (DVector &Dx_x0, DVector &Dx_p, DVector &Dx_u, DVector &Dx_w, int order) const =0
 
virtual returnValue getProtectedForwardSensitivities (DMatrix *Dx, int order) const =0
 
virtual returnValue getProtectedX (DVector *xEnd) const =0
 
void initializeOptions ()
 
virtual returnValue setProtectedBackwardSeed (const DVector &seed, const int &order)=0
 
virtual returnValue setProtectedForwardSeed (const DVector &xSeed, const DVector &pSeed, const DVector &uSeed, const DVector &wSeed, const int &order)=0
 
virtual returnValue setupLogging ()
 
virtual returnValue setupOptions ()
 

Protected Attributes

int * alg_index
 
int * control_index
 
int count
 
int count2
 
int count3
 
int * ddiff_index
 
VariablesGrid ddxStore
 
int * diff_index
 
DVector diff_scale
 
int * disturbance_index
 
DVector dP
 
DVector dPb
 
DVector dU
 
DVector dUb
 
DVector dW
 
DVector dWb
 
DVector dX
 
DVector dXb
 
VariablesGrid dxStore
 
RealClock functionEvaluation
 
double * h
 
double hini
 
double hmax
 
double hmin
 
int * int_control_index
 
int * int_parameter_index
 
VariablesGrid iStore
 
int las
 
short int m
 
short int ma
 
int maxNumberOfSteps
 
short int md
 
short int mdx
 
short int mn
 
short int mp
 
short int mpi
 
short int mu
 
short int mui
 
short int mw
 
int nBDirs
 
int nBDirs2
 
int nFcnEvaluations
 
int nFDirs
 
int nFDirs2
 
int * parameter_index
 
int PrintLevel
 
DifferentialEquationrhs
 
StateOfAggregation soa
 
int time_index
 
Grid timeInterval
 
double TOL
 
RealClock totalTime
 
Transitiontransition
 
double tune
 
DVector xE
 
VariablesGrid xStore
 
- Protected Attributes inherited from AlgorithmicBase
int outputLoggingIdx
 
BooleanType useModuleStandalone
 
UserInteractionuserInteraction
 

Friends

class ShootingMethod
 
class SimulationByIntegration
 

Detailed Description

Abstract base class for all kinds of algorithms for integrating differential equations (ODEs or DAEs).

The class Integrator serves as an abstract base class for all kinds of algorithms for integrating differential equations (ODEs or DAEs).

Author
Boris Houska, Hans Joachim Ferreau

Definition at line 61 of file integrator.hpp.

Constructor & Destructor Documentation

BEGIN_NAMESPACE_ACADO Integrator::Integrator ( )

Default constructor.

Definition at line 56 of file integrator.cpp.

Integrator::Integrator ( const Integrator arg)

Copy constructor.

Definition at line 131 of file integrator.cpp.

Integrator::~Integrator ( )
virtual

Default destructor.

Definition at line 139 of file integrator.cpp.

Member Function Documentation

BooleanType Integrator::canHandleImplicitSwitches ( ) const
virtual

Returns if integrator is able to handle implicit switches.

Returns
BT_TRUE: if integrator can handle implicit switches.
BT_FALSE: otherwise

Definition at line 459 of file integrator.cpp.

virtual Integrator* Integrator::clone ( ) const
pure virtual
returnValue Integrator::deleteAllSeeds ( )
virtual

Deletes all seeds that have been set with the methods above.
This function will also give the corresponding memory free.

Returns
SUCCESSFUL_RETURN
RET_NO_SEED_ALLOCATED

Definition at line 488 of file integrator.cpp.

returnValue Integrator::diffTransitionBackward ( DVector DX,
DVector DP,
DVector DU,
DVector DW,
int &  order 
)
protectedvirtual

Definition at line 542 of file integrator.cpp.

returnValue Integrator::diffTransitionForward ( DVector DX,
const DVector DP,
const DVector DU,
const DVector DW,
const int &  order 
)
protectedvirtual

Definition at line 521 of file integrator.cpp.

virtual returnValue Integrator::evaluate ( const DVector x0,
const DVector xa,
const DVector p,
const DVector u,
const DVector w,
const Grid t_ 
)
protectedpure virtual

Starts integration: cf. integrate(...) for
more details.

Parameters
x0the initial state
xathe initial algebraic state
pthe parameters
uthe controls
wthe disturbance
t_the time interval

Implemented in IntegratorLYAPUNOV, IntegratorRK, and IntegratorBDF.

virtual returnValue Integrator::evaluateSensitivities ( )
pure virtual
returnValue Integrator::evaluateTransition ( const double  time,
DVector xd,
const DVector xa,
const DVector p,
const DVector u,
const DVector w 
)
protectedvirtual

Evaluates the transtion (protected, only for internal use).

Parameters
timethe time
xdthe state
xathe algebraic state
pthe parameters
uthe controls
wthe disturbance

Definition at line 505 of file integrator.cpp.

virtual returnValue Integrator::freezeAll ( )
pure virtual

Freezes the mesh as well as all intermediate values. This function
is necessary for the case that automatic differentiation in backward
mode should is used. (Note: This function might for large right hand
sides lead to memory problems as all intemediate values will be
stored!)

Returns
SUCCESSFUL_RETURN
RET_ALREADY_FROZEN

Implemented in IntegratorRK, IntegratorLYAPUNOV, and IntegratorBDF.

virtual returnValue Integrator::freezeMesh ( )
pure virtual

Freezes the mesh: Storage of the step sizes. If the integrator is
freezed, the mesh will be stored when calling the function integrate
for the first time. If the function integrate is called more than
once the same mesh will be reused (i.e. the step size control will
be turned off). Note that the mesh should be frozen if any kind of
sensitivity generation is used.

Returns
SUCCESSFUL_RETURN
RET_ALREADY_FROZEN

Implemented in IntegratorRK, IntegratorLYAPUNOV, and IntegratorBDF.

returnValue Integrator::getBackwardSensitivities ( DVector Dx_x0,
DVector Dx_p,
DVector Dx_u,
DVector Dx_w,
int  order 
) const

Returns the result for the backward sensitivities at the time tend.

Parameters
Dx_x0backward sensitivities w.r.t. the initial states
Dx_pbackward sensitivities w.r.t. the parameters
Dx_ubackward sensitivities w.r.t. the controls
Dx_wbackward sensitivities w.r.t. the disturbance

Returns
SUCCESSFUL_RETURN
RET_INPUT_OUT_OF_RANGE

Definition at line 424 of file integrator.cpp.

double Integrator::getDifferentialEquationSampleTime ( ) const
virtual

Returns the

Definition at line 478 of file integrator.cpp.

virtual int Integrator::getDim ( ) const
protectedpure virtual

Returns the dimension of the Differential Equation

Implemented in IntegratorBDF, IntegratorLYAPUNOV, and IntegratorRK.

int Integrator::getDimX ( ) const
protectedvirtual

Returns the number of Dynamic Equations in the Differential Equation

Reimplemented in IntegratorBDF.

Definition at line 618 of file integrator.cpp.

returnValue Integrator::getForwardSensitivities ( DVector Dx,
int  order 
) const

Returns the result for the forward sensitivities at the time tend.

Parameters
Dxthe result for the forward sensitivities.
orderthe order.

Returns
SUCCESSFUL_RETURN
RET_INPUT_OUT_OF_RANGE

Definition at line 406 of file integrator.cpp.

returnValue Integrator::getForwardSensitivities ( VariablesGrid Dx,
int  order 
) const

Returns the result for the forward sensitivities on the grid.

Parameters
Dxthe result for the forward sensitivities.
orderthe order.

Returns
SUCCESSFUL_RETURN
RET_INPUT_OUT_OF_RANGE

Definition at line 414 of file integrator.cpp.

returnValue Integrator::getI ( VariablesGrid I) const
inline

Returns the requested output on the specified grid.
The intermediate states constructed based on linear
interpolation.

Parameters
Ithe intermediates states on the grid.

Returns
SUCCESSFUL_RETURN
RET_TRAJECTORY_NOT_FROZEN
virtual int Integrator::getNumberOfRejectedSteps ( ) const
pure virtual

Returns the number of rejected Steps.

Returns
The requested number of rejected steps.

Implemented in IntegratorBDF, IntegratorRK, and IntegratorLYAPUNOV.

virtual int Integrator::getNumberOfSteps ( ) const
pure virtual

Returns the number of accepted Steps.

Returns
The requested number of accepted steps.

Implemented in IntegratorBDF, IntegratorRK, and IntegratorLYAPUNOV.

virtual returnValue Integrator::getProtectedBackwardSensitivities ( DVector Dx_x0,
DVector Dx_p,
DVector Dx_u,
DVector Dx_w,
int  order 
) const
protectedpure virtual

Returns the result for the backward sensitivities at the time tend.

Parameters
Dx_x0backward sensitivities w.r.t. the initial states
Dx_pbackward sensitivities w.r.t. the parameters
Dx_ubackward sensitivities w.r.t. the controls
Dx_wbackward sensitivities w.r.t. the disturbance
orderthe order of the derivative

Returns
SUCCESSFUL_RETURN
RET_INPUT_OUT_OF_RANGE

Implemented in IntegratorLYAPUNOV, IntegratorRK, and IntegratorBDF.

virtual returnValue Integrator::getProtectedForwardSensitivities ( DMatrix Dx,
int  order 
) const
protectedpure virtual

Returns the result for the forward sensitivities at the time tend.

Returns
SUCCESSFUL_RETURN
RET_INPUT_OUT_OF_RANGE
Parameters
Dxthe result for the forward sensitivi- ties
orderthe order

Implemented in IntegratorLYAPUNOV, IntegratorRK, and IntegratorBDF.

virtual returnValue Integrator::getProtectedX ( DVector xEnd) const
protectedpure virtual

Returns the result for the state at the time tend.

Returns
SUCCESSFUL_RETURN
Parameters
xEndthe result for the states at the time tend.

Implemented in IntegratorLYAPUNOV, IntegratorRK, and IntegratorBDF.

virtual double Integrator::getStepSize ( ) const
pure virtual

Returns the current step size

Implemented in IntegratorBDF, IntegratorRK, and IntegratorLYAPUNOV.

returnValue Integrator::getX ( DVector xEnd) const
inline

Returns the result for the differential states at the
time tend.

Parameters
xEndthe result for the states at the time tend.

Returns
SUCCESSFUL_RETURN
returnValue Integrator::getX ( VariablesGrid X) const
inline

Returns the requested output on the specified grid. Note
that the output X will be evaluated based on polynomial
interpolation depending on the order of the integrator.

Parameters
Xthe differential states on the grid.

Returns
SUCCESSFUL_RETURN
RET_TRAJECTORY_NOT_FROZEN
returnValue Integrator::getXA ( DVector xaEnd) const
inline

Returns the result for the algebraic states at the
time tend.

Parameters
xaEndthe result for the algebraic states at the time tend.

Returns
SUCCESSFUL_RETURN
returnValue Integrator::getXA ( VariablesGrid XA) const
inline

Returns the requested output on the specified grid. Note
that the output X will be evaluated based on polynomial
interpolation depending on the order of the integrator.

Parameters
XAthe algebraic states on the grid.

Returns
SUCCESSFUL_RETURN
RET_TRAJECTORY_NOT_FROZEN
virtual returnValue Integrator::init ( const DifferentialEquation rhs)
pure virtual

The initialization routine which takes the right-hand side of
the differential equation to be integrated.

Parameters
rhsthe right-hand side of the ODE/DAE.

Returns
SUCCESSFUL_RETURN if all dimension checks succeed.
otherwise: integrator dependent error message.

Implemented in IntegratorLYAPUNOV, IntegratorRK, IntegratorBDF, IntegratorDiscretizedODE, and IntegratorRK12.

returnValue Integrator::init ( const DifferentialEquation rhs,
const Transition trs 
)

The initialization routine which takes the right-hand side of
the differential equation to be integrated. In addition a
transition function can be set which is evaluated at the end
of the integration interval.

Parameters
rhsthe right-hand side of the ODE/DAE.
trsthe transition to be evaluated at the end.

Returns
SUCCESSFUL_RETURN if all dimension checks succeed.
otherwise: integrator dependent error message.

Definition at line 185 of file integrator.cpp.

void Integrator::initializeOptions ( )
protected

Initializes the options. This routine should be called by the
integrators derived from this class in order to initialize
the following members based on the user options:
maxNumberOfSteps hmin hmax tune TOL las (cf. SETTINGS for more details)

Definition at line 587 of file integrator.cpp.

returnValue Integrator::integrate ( double  t0,
double  tend,
double *  x0,
double *  xa = 0,
double *  p = 0,
double *  u = 0,
double *  w = 0 
)

Starts the integration of the right hand side at time t0.
If neither the maximum number of allowed iteration is
exceeded nor any other error occurs the functions stops the
integration at time tend.

Returns
SUCCESFUL_RETURN or
a error message that depends on the specific
integration routine. (cf. the corresponding header
file that implements the integration routine)
Parameters
t0the start time
tendthe end time
x0the initial state
xathe initial algebraic state
pthe parameters
uthe controls
wthe disturbance

Definition at line 207 of file integrator.cpp.

returnValue Integrator::integrate ( const Grid t,
double *  x0,
double *  xa = 0,
double *  p = 0,
double *  u = 0,
double *  w = 0 
)

Starts the integration of the right hand side at time t0.
If neither the maximum number of allowed iteration is
exceeded nor any other error occurs the functions stops the
integration at time tend.
In addition, results at intermediate grid points can be
stored. Note that these grid points are for storage only and
have nothing to do the integrator steps.

Returns
SUCCESFUL_RETURN or
a error message that depends on the specific
integration routine. (cf. the corresponding header
file that implements the integration routine)
Parameters
tthe grid [t0,tend]
x0the initial state
xathe initial algebraic state
pthe parameters
uthe controls
wthe disturbance

Definition at line 220 of file integrator.cpp.

returnValue Integrator::integrate ( double  t0,
double  tend,
const DVector x0,
const DVector xa = emptyVector,
const DVector p = emptyVector,
const DVector u = emptyVector,
const DVector w = emptyVector 
)

Starts the integration of the right hand side at time t0.
If neither the maximum number of allowed iteration is
exceeded nor any other error occurs the functions stops the
integration at time tend.

Returns
SUCCESFUL_RETURN or
a error message that depends on the specific
integration routine. (cf. the corresponding header
file that implements the integration routine)
Parameters
t0the start time
tendthe end time
x0the initial state
xathe initial algebraic state
pthe parameters
uthe controls
wthe disturbance

Definition at line 241 of file integrator.cpp.

returnValue Integrator::integrate ( const Grid t,
const DVector x0,
const DVector xa = emptyVector,
const DVector p = emptyVector,
const DVector u = emptyVector,
const DVector w = emptyVector 
)

Starts the integration of the right hand side at time t0.
If neither the maximum number of allowed iteration is
exceeded nor any other error occurs the functions stops the
integration at time tend.
In addition, results at intermediate grid points can be
stored. Note that these grid points are for storage only and
have nothing to do the integrator steps.

Returns
SUCCESFUL_RETURN or
a error message that depends on the specific
integration routine. (cf. the corresponding header
file that implements the integration routine)
Parameters
tthe grid [t0,tend]
x0the initial state
xathe initial algebraic state
pthe parameters
uthe controls
wthe disturbance

Definition at line 255 of file integrator.cpp.

returnValue Integrator::integrateSensitivities ( )

< Integrates forward and/or backward depending on the specified seeds.

Returns
SUCCESSFUL_RETURN
RET_NO_SEED_ALLOCATED

Definition at line 357 of file integrator.cpp.

BooleanType Integrator::isDifferentialEquationAffine ( ) const
virtual

Returns if the differential equation of the integrator is affine.

Returns
BT_TRUE: if differential equation is affine.
BT_FALSE: otherwise

Definition at line 471 of file integrator.cpp.

BooleanType Integrator::isDifferentialEquationDefined ( ) const
virtual

Returns if the differential equation of the integrator is defined.

Returns
BT_TRUE: if differential equation is defined.
BT_FALSE: otherwise

Definition at line 465 of file integrator.cpp.

returnValue Integrator::printRunTimeProfile ( ) const
virtual

Prints the run-time profile. This routine
can be used after an integration run in
order to assess the performance.
Integrates forward and/or backward depending on the specified seeds.

Returns
SUCCESSFUL_RETURN
RET_NO_SEED_ALLOCATED

Definition at line 624 of file integrator.cpp.

returnValue Integrator::setBackwardSeed ( const int &  order,
const DVector seed 
)

Define a backward seed

Returns
SUCCESFUL_RETURN
RET_INPUT_OUT_OF_RANGE
Parameters
orderthe order of the seed.
seedthe backward seed

Definition at line 338 of file integrator.cpp.

virtual returnValue Integrator::setDxInitialization ( double *  dx0)
pure virtual

Sets an initial guess for the differential state derivatives
(consistency condition)

Returns
SUCCESSFUL_RETURN
Parameters
dx0initial guess for the differential state derivatives

Implemented in IntegratorRK, IntegratorLYAPUNOV, and IntegratorBDF.

returnValue Integrator::setForwardSeed ( const int &  order,
const DVector xSeed,
const DVector pSeed = emptyVector,
const DVector uSeed = emptyVector,
const DVector wSeed = emptyVector 
)

Define a forward seed matrix.

Returns
SUCCESFUL RETURN
RET_INPUT_OUT_OF_RANGE
Parameters
orderthe order of the seed.
xSeedthe seed w.r.t states
pSeedthe seed w.r.t parameters
uSeedthe seed w.r.t controls
wSeedthe seed w.r.t disturbances

Definition at line 308 of file integrator.cpp.

virtual returnValue Integrator::setProtectedBackwardSeed ( const DVector seed,
const int &  order 
)
protectedpure virtual

Define a backward seed

Returns
SUCCESFUL_RETURN
RET_INPUT_OUT_OF_RANGE
Parameters
seedthe seed matrix
orderthe order of the seed.

Implemented in IntegratorLYAPUNOV, IntegratorRK, and IntegratorBDF.

virtual returnValue Integrator::setProtectedForwardSeed ( const DVector xSeed,
const DVector pSeed,
const DVector uSeed,
const DVector wSeed,
const int &  order 
)
protectedpure virtual

Define a forward seed.

Returns
SUCCESFUL RETURN
RET_INPUT_OUT_OF_RANGE
Parameters
xSeedthe seed w.r.t the initial states
pSeedthe seed w.r.t the parameters
uSeedthe seed w.r.t the controls
wSeedthe seed w.r.t the disturbances
orderthe order of the seed.

Implemented in IntegratorLYAPUNOV, IntegratorRK, and IntegratorBDF.

returnValue Integrator::setTransition ( const Transition trs)

Definition at line 198 of file integrator.cpp.

returnValue Integrator::setupLogging ( )
protectedvirtual

Definition at line 599 of file integrator.cpp.

returnValue Integrator::setupOptions ( )
protectedvirtual

Definition at line 566 of file integrator.cpp.

virtual returnValue Integrator::unfreeze ( )
pure virtual

Unfreezes the mesh: Gives the memory free that has previously
been allocated by "freeze". If you use the function
integrate after unfreezing the usual step size control will be
switched on.

Returns
SUCCESSFUL_RETURN
RET_MESH_ALREADY_UNFREEZED

Implemented in IntegratorRK, IntegratorLYAPUNOV, and IntegratorBDF.

Friends And Related Function Documentation

friend class ShootingMethod
friend

Definition at line 65 of file integrator.hpp.

friend class SimulationByIntegration
friend

Definition at line 64 of file integrator.hpp.

Member Data Documentation

int* Integrator::alg_index
protected

the index list for the algebraic states

Definition at line 655 of file integrator.hpp.

int* Integrator::control_index
protected

the index list for the controls

Definition at line 656 of file integrator.hpp.

int Integrator::count
protected

a counter of the (accepted) steps

Definition at line 667 of file integrator.hpp.

int Integrator::count2
protected

number of steps after integration

Definition at line 668 of file integrator.hpp.

int Integrator::count3
protected

number of steps after integration

Definition at line 669 of file integrator.hpp.

int* Integrator::ddiff_index
protected

the index list for the differential state derivatives

Definition at line 654 of file integrator.hpp.

VariablesGrid Integrator::ddxStore
protected

Definition at line 713 of file integrator.hpp.

int* Integrator::diff_index
protected

the index list for the differential states

Definition at line 653 of file integrator.hpp.

DVector Integrator::diff_scale
protected

the scale of the differential states

Definition at line 670 of file integrator.hpp.

int* Integrator::disturbance_index
protected

the index list for the disturbances

Definition at line 660 of file integrator.hpp.

DVector Integrator::dP
protected

Definition at line 702 of file integrator.hpp.

DVector Integrator::dPb
protected

Definition at line 707 of file integrator.hpp.

DVector Integrator::dU
protected

Definition at line 703 of file integrator.hpp.

DVector Integrator::dUb
protected

Definition at line 708 of file integrator.hpp.

DVector Integrator::dW
protected

Definition at line 704 of file integrator.hpp.

DVector Integrator::dWb
protected

Definition at line 709 of file integrator.hpp.

DVector Integrator::dX
protected

Definition at line 701 of file integrator.hpp.

DVector Integrator::dXb
protected

Definition at line 706 of file integrator.hpp.

VariablesGrid Integrator::dxStore
protected

Definition at line 712 of file integrator.hpp.

RealClock Integrator::functionEvaluation
protected

Definition at line 695 of file integrator.hpp.

double* Integrator::h
protected

the initial step size = h[0]

Definition at line 640 of file integrator.hpp.

double Integrator::hini
protected

storage of the initial step size

Definition at line 641 of file integrator.hpp.

double Integrator::hmax
protected

the maximum step size

Definition at line 643 of file integrator.hpp.

double Integrator::hmin
protected

the minimum step size

Definition at line 642 of file integrator.hpp.

int* Integrator::int_control_index
protected

the index list for the integer controls

Definition at line 658 of file integrator.hpp.

int* Integrator::int_parameter_index
protected

the index list for the integer parameters

Definition at line 659 of file integrator.hpp.

VariablesGrid Integrator::iStore
protected

Definition at line 714 of file integrator.hpp.

int Integrator::las
protected

Definition at line 646 of file integrator.hpp.

short int Integrator::m
protected

the dimension of the right-hand side

Definition at line 620 of file integrator.hpp.

short int Integrator::ma
protected

the number of algebraic states

Definition at line 621 of file integrator.hpp.

int Integrator::maxNumberOfSteps
protected

the maximum number of integrator steps.

Definition at line 666 of file integrator.hpp.

short int Integrator::md
protected

number of differential states

Definition at line 629 of file integrator.hpp.

short int Integrator::mdx
protected

the dimension of differential state derivatives

Definition at line 622 of file integrator.hpp.

short int Integrator::mn
protected

the number of intermediate states in the rhs

Definition at line 623 of file integrator.hpp.

short int Integrator::mp
protected

the number of parameters

Definition at line 626 of file integrator.hpp.

short int Integrator::mpi
protected

the number of integer parameters

Definition at line 627 of file integrator.hpp.

short int Integrator::mu
protected

the number of controls

Definition at line 624 of file integrator.hpp.

short int Integrator::mui
protected

the number of integer controls

Definition at line 625 of file integrator.hpp.

short int Integrator::mw
protected

the number of disturbances

Definition at line 628 of file integrator.hpp.

int Integrator::nBDirs
protected

The number of backward directions

Definition at line 681 of file integrator.hpp.

int Integrator::nBDirs2
protected

The number of second order backward directions

Definition at line 683 of file integrator.hpp.

int Integrator::nFcnEvaluations
protected

Definition at line 696 of file integrator.hpp.

int Integrator::nFDirs
protected

The number of forward directions

Definition at line 680 of file integrator.hpp.

int Integrator::nFDirs2
protected

The number of second order forward directions

Definition at line 682 of file integrator.hpp.

int* Integrator::parameter_index
protected

the index list for the parameters

Definition at line 657 of file integrator.hpp.

int Integrator::PrintLevel
protected

The PrintLevel (default: LOW)

Definition at line 675 of file integrator.hpp.

DifferentialEquation* Integrator::rhs
protected

the right-hand side to be integrated

Definition at line 619 of file integrator.hpp.

StateOfAggregation Integrator::soa
protected

the state of aggregation

Definition at line 688 of file integrator.hpp.

int Integrator::time_index
protected

the time index

Definition at line 661 of file integrator.hpp.

Grid Integrator::timeInterval
protected

the type of linear algebra solver to be used the time interval

Definition at line 648 of file integrator.hpp.

double Integrator::TOL
protected

the integration tolerance

Definition at line 645 of file integrator.hpp.

RealClock Integrator::totalTime
protected

Definition at line 694 of file integrator.hpp.

Transition* Integrator::transition
protected

the transition to be evaluated at switches

Definition at line 634 of file integrator.hpp.

double Integrator::tune
protected

tuning parameter for the step size control.

Definition at line 644 of file integrator.hpp.

DVector Integrator::xE
protected

Definition at line 699 of file integrator.hpp.

VariablesGrid Integrator::xStore
protected

Definition at line 711 of file integrator.hpp.


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


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Mon Jun 10 2019 12:35:25