File kernel/src/simulationTools/Hem5OSI.hpp#
Go to the source code of this file
Hem5OSI solver (from E.
Hairer software lists)
-
class Hem5OSI : public OneStepIntegrator
- #include <Hem5OSI.hpp>
Hem5OSI solver (odepack)
Many parameters are required as input/output for LSODAR. See the documentation of this function in externals/odepack/opkdmain.f to have a full description of these parameters.
Most of them are read-only parameters (ie can not be set by user).
Except:jt: Jacobian type indicator (1 means a user-supplied full Jacobian, 2 means an internally generated full Jacobian).
Default = 2.
itol, rtol and atol
ITOL = an indicator for the type of error control.
RTOL = a relative error tolerance parameter, either a scalar or array of length NEQ.
ATOL = an absolute error tolerance parameter, either a scalar or an array of length NEQ. Input only.
Public Types
Public Functions
-
Hem5OSI()
Default and only constructor.
-
~Hem5OSI() noexcept = default
destructor
-
inline const std::vector<int> intData() const
get vector of int parameters for lsodar
- Returns:
a vector<int>
-
inline int intData(unsigned int i) const
get _intData[i]
- Parameters:
i – index
- Returns:
an int
-
inline void setIntData(unsigned int i, int newValue)
set _intData[i]
- Parameters:
i – index
newValue –
-
inline const std::vector<double> &getRtol() const
- Returns:
relative tolerance parameter for Hem5
-
inline const std::vector<double> &getAtol() const
- Returns:
absolute tolerance parameter for Hem5
-
inline int getMaxNstep() const
get the maximum number of steps for one call
- Returns:
an interger
-
inline const std::vector<double> &getRwork() const
- Returns:
real work vector parameter for lsodar
-
inline const std::vector<int> &getIwork() const
- Returns:
iwork
-
void setTol(int itol, std::vector<double> &&rtol, std::vector<double> &&atol)
set itol, rtol and atol (tolerance parameters for Hem5)
- Parameters:
itol – int (itol value)
rtol – double * (rtol)
atol – double * (atol)
-
void setTol(int itol, double rtol, double atol)
set itol, rtol and atol (scalar tolerance parameters for Hem5)
- Parameters:
itol – int (itol value)
rtol – double (rtol)
atol – double (atol)
-
void setMaxNstep(int nstepmax)
set the maximul number of steps for one call of Hem5OSI
- Parameters:
nstepmax – an int
-
void setMaxStepSize(double maxstepsize)
set the minimum and maximum step sizes
- Parameters:
maxstepsize – double (maximul step size)
-
void updateIntData()
update _intData
-
void updateData()
update doubleData and iwork memory size, when changes occur in _intData.
-
void fillqWork(int *sizex, double *x)
fill qWork with a double
- Parameters:
sizex – int*, size of x array
x – double* x:array of double
-
void fillvWork(int *sizex, double *x)
fill vWork with a double
- Parameters:
sizex – int*, size of x array
x – double* x:array of double
-
void computeRhs(double)
compute rhs(t) for all dynamical systems in the set
-
void computeJacobianRhs(double)
compute jacobian of the rhs at time t for all dynamical systems in the set
-
unsigned int numberOfConstraints()#
-
void f(int *sizeOfX, double *time, double *x, double *xdot)#
-
void g(int *nEq, double *time, double *x, int *ng, double *gOut)#
-
void jacobianfx(int*, double*, double*, int*, int*, double*, int*)#
-
virtual void initialize() override
initialization of the integrator
-
virtual void initializeWorkVectorsForDS(double t, SP::DynamicalSystem ds) override
initialization of the work vectors and matrices (properties) related to one dynamical system on the graph and needed by the osi
- Parameters:
t – time of initialization
ds – the dynamical system
-
virtual void initializeWorkVectorsForInteraction(Interaction &inter, InteractionProperties &interProp, DynamicalSystemsGraph &DSG) override
initialization of the work vectors and matrices (properties) related to one interaction on the graph and needed by the osi
- Parameters:
inter – the interaction
interProp – the properties on the graph
DSG – the dynamical systems graph
-
inline virtual unsigned int numberOfIndexSets() const override
get the number of index sets required for the simulation
- Returns:
unsigned int
-
virtual void integrate(double &tinit, double &tend, double &tout, int &idid) override
integrate the system, between tinit and tend (->iout=true), with possible stop at tout (->iout=false)
- Parameters:
tinit – initial time
tend – end time
tout – real end time
idid – in-out parameter, input: 1 for first call, else 2. Output: 2 if no root was found, else 3.
-
virtual void updateState(const unsigned int level) override
update the state of the DynamicalSystems attached to this Integrator
- Parameters:
level – level of interest for the dynamics
-
inline virtual void prepareNewtonIteration(double time) override#
-
virtual void computeFreeOutput(InteractionsGraph::VDescriptor &vertex_inter, OneStepNSProblem *osnsp) override
integrates the Interaction linked to this integrator, without taking non-smooth effects into account
- Parameters:
vertex_inter – of the interaction graph
osnsp – pointer to OneStepNSProblem
-
inline virtual SiconosVector &osnsp_rhs(InteractionsGraph::VDescriptor &vertex_inter, InteractionsGraph &indexSet) override
return the workVector corresponding to the right hand side of the OneStepNonsmooth problem
-
virtual void display() override
print the data to the screen
-
ACCEPT_STD_VISITORS()#
Public Members
-
std::shared_ptr<Hem5OSI_impl> _impl = {nullptr}#
Private Members
-
std::vector<int> _intData = {0, 0, 0, 0, 0, 0, 0, 0, 0}#
vector of integer data for the integrator _intData[0] NQ size of the position vector q _intData[1] NV size of the velocity vector v NQ >= NQ _intData[2] NU size of the external dynamic vector u _intData[3] NL size of the Lagrange multiplier vector lambda _intData[4] ITOL indicates whether RTOL and ATOL are scalar (ITOL=0), or array of dimension NQ + NV + NU (ITOL=1) _intData[5] IOUT selects the dense output formula _intData[6] LWK length of real array rwork _intData[7] LIWK length of integer array iwork See hem5.f
-
int _idid = {0}#
-
std::vector<double> rtol = {}#
relative tolerance
-
std::vector<double> atol = {}#
absolute tolerance
-
std::vector<double> rwork = {}#
real work array
-
std::vector<int> iwork = {}#
int work array
-
double _timeStep = {INITIAL_GUESS_TS}#
Private Static Attributes
-
static constexpr auto HEM5_ATOL_DEFAULT = 100 * siconos::internal::MACHINE_PREC#
-
static constexpr auto HEM5_RTOL_DEFAULT = 10 * siconos::internal::MACHINE_PREC#
-
static constexpr double INITIAL_GUESS_TS = 1.e-3#
Friends
- friend struct _NSLEffectOnFreeOutput
-
class Hem5OSI_impl#
Public Functions
-
void fprob(int *IFCN, int *NQ, int *NV, int *NU, int *NL, int *LDG, int *LDF, int *LDA, int *NBLK, int *NMRC, int *NPGP, int *NPFL, int *INDGR, int *INDGC, int *INDFLR, int *INDFLC, double *time, double *q, double *v, double *u, double *xl, double *G, double *GQ, double *F, double *GQQ, double *GT, double *FL, double *QDOT, double *UDOT, double *AM)#
-
void solout(int *MODE, int *NSTEP, int *NQ, int *NV, int *NU, int *NL, int *LDG, int *LDF, int *LDA, int *LRDO, int *LIDO, siconos::fortran::hairer::fprobpointer FPROB, double *q, double *v, double *u, double *DOWK, int *IDOWK)#
-
void fprob(int *IFCN, int *NQ, int *NV, int *NU, int *NL, int *LDG, int *LDF, int *LDA, int *NBLK, int *NMRC, int *NPGP, int *NPFL, int *INDGR, int *INDGC, int *INDFLR, int *INDFLC, double *time, double *q, double *v, double *u, double *xl, double *G, double *GQ, double *F, double *GQQ, double *GT, double *FL, double *QDOT, double *UDOT, double *AM)#