FLINT 0.9.9
Fortran Library for numerical INTegration of differential equations
|
Interface for the main Integrate method. It must be called after initialization and can be called multiple times with different IC and options without calling the init routine first every time. More...
Public Member Functions | |
subroutine | integrate (me, x0, y0, xf, yf, stepsz, useconststepsz, stepadvance, intstepson, xint, yint, eventmask, eventstates, stifftest, params) |
Interface for the main Integrate method. It must be called after initialization and can be called multiple times with different IC and options without calling the init routine first every time.
Definition at line 400 of file FLINT_base.f90.
|
virtual |
[in,out] | me | Object of the class type FLINT |
[in] | x0 | Initial value of the independent variable |
[in] | y0 | Initial conditions |
[in,out] | xf | Final value of the independent variable |
[out] | yf | Final Solution |
[in,out] | stepsz | User can specify the initial step-size to use. If it is 0, then FLINT will compute the initial step-size internally. After the completion of the integration, the last accepted step size will be returned back. Note that user must provide the correct sign of the desired step size. |
[in] | useconststepsz | If true, then FLINT will not use adaptive step size algorithm, rather it will use the initial step size value given in StepSz as the integrator's natural step. StepSz must take into account the sign of the step otherwise FLINT will return an error. By default, this option is set to False. |
[in] | stepadvance | If true, FLINT will only advance the integration by one step and return. By default, this option is set to False. |
[in] | intstepson | If true, then FLINT will return the solution at the integrator's natural step size. |
[out] | xint | If IntStepsOn=true, user must provide this parameter. FLINT will allocate the memory and return the integrator's natural steps taken in Xint array. FLINT will not explicitly deallocate this memory. |
[out] | yint | If IntStepsOn=true, user must provide this parameter. FLINT will allocate the memory and return the solution at the integrator's natural steps in Yint. The first column is the initial condition Y0 and the final column is Yf. FLINT will not explicitly deallocate this memory. In case of an event action changing the solution at the location of the event, the new solution will be returned in Yint. |
[in] | eventmask | If EventMask(i)=false, then the event "i" will Not be detected. Default is true. If all the events are masked, then the event detection will be turned off. |
[out] | eventstates | Solutions at the event locations for each of the unmasked event. If events are enabled, then user must specify this parameter. FLINT will internally allocate the storage and it will not explicitly deallocate this memory. EventStates(1:(n+2),i) = [Xevent, Yevent, EventID], where EventID is the index of the triggered event (between 1 and m) to which Xevent and Yevent corresponds. In case of an event action changing the solution at the location of the event, EventStates will still return the orginal state. |
[in,out] | stifftest | FLINT will use Hairer's DOP853 alogorithm to test for stiffness of the differential equations if StffTest == 1 or StffTest == 2. Test is disabled by default or if Stifftest == 0.
|
[in] | params | Parameters for sensitivity analysis. Not implemented yet. |
Definition at line 400 of file FLINT_base.f90.