FLINT 0.9.9
Fortran Library for numerical INTegration of differential equations
|
Abstract base class that must be inherited by all the numerical integrator classes. More...
Public Member Functions | |
procedure(init), deferred | init (me, de, maxsteps, method, atol, rtol, interpon, interpstates, minstepsize, maxstepsize, stepszparams, eventson, eventstepsz, eventoptions, eventtol) |
Init is called once and only to specify new differential equaions, events, and a few integrator options. | |
procedure(integrate), deferred | integrate (me, x0, y0, xf, yf, stepsz, useconststepsz, stepadvance, intstepson, xint, yint, eventmask, eventstates, stifftest, params) |
Integrate the solution from the initial to final value of the independent variable. | |
procedure(interpolate), deferred | interpolate (me, xarr, yarr, saveinterpcoeffs) |
If interpolation is enabled, then it outputs the solution at the user-specified grid points. | |
procedure(info), deferred | info (me, laststatus, statusmsg, nsteps, naccept, nreject, nfcalls, interpready, h0, x0, y0, hf, xf, yf) |
Current FLINT object query function. | |
Public Attributes | |
integer(kind(flint_success)), public | status = FLINT_SUCCESS |
Status of the last operation completed. | |
class(diffeqsys), pointer, public | pdiffeqsys => null() |
Pointer to the user-supplied Differential Equation object. | |
logical | isscalartol = .TRUE. |
Tolerance related constants: rtol and atol Dimension can be either 1 or n (dimension of the initial condition vector). These are used to compute the scale factor (SC) needed for error computations using the following formula: | |
real(wp), dimension(:), allocatable | rtol |
real(wp), dimension(:), allocatable | atol |
integer | totalsteps = 0 |
Total number of integrator steps. | |
integer | acceptedsteps = 0 |
Number of accepted integrator steps. | |
integer | rejectedsteps = 0 |
Number of rejected integrator steps. | |
integer | fcalls = 0 |
Total number of Diff Eq. function calls made. | |
integer | maxsteps = INT_MAXSTEPS |
Maximum number of steps integrator is allowed to take. | |
logical | interpon = .FALSE. |
True for dense output. | |
integer, dimension(:), allocatable | interpstates |
States between 1 and n for dense output. | |
real(wp), dimension(:), allocatable | xint |
Array of the natural integrator steps. | |
real(wp), dimension(:,:,:), allocatable | bip |
Array to store interpolation coefficients at integrator steps. | |
real(wp) | minstepsize |
Minimum step size of the integrator. It is always a positive number. Default value is used if user do not specifify it. | |
real(wp) | maxstepsize |
Maximum step-size of the integrator. It is always a positive number. Default value is the difference between the initial and final value of the independent variable. | |
real(wp) | initialstepsize |
Initial step size. | |
real(wp), dimension(6) | stepszparams |
Step size method-specific tuning parameters. | |
logical | eventson = .FALSE. |
Events checking off by default. | |
real(wp), dimension(:), allocatable | eventstepsz |
This specifies the maximum time after which the event function must be called for each event detection. A 0 value in EventStepSz means that event will be checked after the integrator takes a natural step always. Otherwise, events will be checked after the interval specified by EventStepSz if that interval size is smaller than the length of the current step size. In case it is bigger, then event will be checked at the integrator's natural step boundary. Setting is for each event spearately. | |
integer, dimension(:), allocatable | eventoptions |
Event options for each event separately. Must be initialized by the Init Function. Specifies what event related info to return to the user. | |
real(wp), dimension(:), allocatable | etol |
Event tolerance value for each event in case of root-finding is used for event loation. | |
real(wp) | x0 |
Initial value of the independent variable. | |
real(wp) | xf |
Final value of the independent variable. | |
real(wp) | h0 |
Initial step-size that was accepted. | |
real(wp) | hf |
The accepted step-size that was used at the last step. | |
real(wp), dimension(:), allocatable | y0 |
Initial condition. | |
real(wp), dimension(:), allocatable | yf |
Final solution. | |
Abstract base class that must be inherited by all the numerical integrator classes.
Definition at line 125 of file FLINT_base.f90.
|
pure virtual |
Current FLINT object query function.
[in,out] | me | Object of class type FLINT |
[out] | laststatus | The status code of the last operation completed |
[out] | statusmsg | Status message corresponding to the status code in LastStatus |
[out] | nsteps | Total number of steps taken |
[out] | naccept | Number of Accepted steps |
[out] | nreject | Number of rejected steps |
[out] | nfcalls | Function calls made not including event calls |
[out] | interpready | Whether Interpolation is enabled |
[out] | h0 | The initial accepted step-size |
[out] | x0 | The initial value of the independent variable |
[out] | y0 | Initial condition |
[out] | hf | Final accepted step-size |
[out] | xf | The final value of the independent variable |
[out] | yf | Final solution |
Definition at line 209 of file FLINT_base.f90.
|
pure virtual |
Init is called once and only to specify new differential equaions, events, and a few integrator options.
[in,out] | me | Object |
de | Differential Equation system object | |
[in] | maxsteps | Must specify max no. of steps The implementation can use this number to decide how and how much to allocate internal storage for storing integrator natural steps as well as the interpolation coefficients for delayed interpolation. |
[in] | method | Type of Intergrator to use. Implemented methods are:
|
[in] | atol | Absolute tolerance, size must be either 1 or n |
[in] | rtol | Relative tolerance, size must be either 1 or n |
[in] | interpon | If true, then the interpolation coefficients along with integrator natural steps are computed and stored internally. User can use Interpolate method to compute solutions at any value of the indepedent variable without doing the actual integrating again. Note that if InterpOn is true, then IntStepsOn option of Integrate must be either omitted or set to false. |
[in] | interpstates | Specify the interpolation states for which the interpolation coefficients will be computed and stored, e.g. [1,3,4,6,n]. |
[in] | minstepsize | User can specify minimum and maximum step sizes allowed. Note that these are positive quantities. Sign will be handled bu FLINT. |
[in] | maxstepsize | User can specify minimum and maximum step sizes allowed. Note that these are positive quantities. Sign will be handled bu FLINT. |
[in] | stepszparams | Step-size computation specific parameters. User can specify these for a specific step-size computation method, otherwise defaults are used. For Hairer's DOPRI5/DOP853 codes:
|
[in] | eventson | Set true for events detection |
[in] | eventstepsz | If events detection is enabled, then user can specify the maximum interval after which events must be checked. If there are even number of sign changes between X and X+EventStepSz, that event cannot be detected. If EventStepSz is not specified, then events are always checked at the integrator's natural steps. Note that these are positive quantities and FLINT will take care of the appropriate sign based on the direction of integration. |
[in] | eventoptions | Event Options that control the event information returned by FLINT Default is FLINT_EVENTOPTION_ROOTFINDING. |
[in] | eventtol | Event tolerance values used for root-finding if enabled, size must be m |
Definition at line 200 of file FLINT_base.f90.
|
pure virtual |
Integrate the solution from the initial to final value of the independent variable.
[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 203 of file FLINT_base.f90.
|
pure virtual |
If interpolation is enabled, then it outputs the solution at the user-specified grid points.
[in,out] | me | Object of class type FLINT |
[in] | xarr | Xarr is provided by the user with values of the independent variable at which the solution is sought. Xarr must be in ascending order and the min and max values must not be smaller and larger than X0 and Xf, respectively. |
[out] | yarr | Yarr provides the solution at each point specified in Xarr. |
[in] | saveinterpcoeffs | If SaveInterpCoeffs=False, then the internal memory which stores interpolation coefficients is freed. After that, Interpolate function can not be called again without calling Integrate first. By default, SaveInterpCoeffs=False. In order to be able to call Interpolate again, user must set it to Trueat for every call to Interpolate. |
Definition at line 206 of file FLINT_base.f90.
integer flint_base::flint_class::acceptedsteps = 0 |
Number of accepted integrator steps.
Definition at line 144 of file FLINT_base.f90.
real(wp), dimension(:), allocatable flint_base::flint_class::atol |
Definition at line 141 of file FLINT_base.f90.
real(wp), dimension(:,:,:), allocatable flint_base::flint_class::bip |
Array to store interpolation coefficients at integrator steps.
Definition at line 157 of file FLINT_base.f90.
real(wp), dimension(:), allocatable flint_base::flint_class::etol |
Event tolerance value for each event in case of root-finding is used for event loation.
Definition at line 187 of file FLINT_base.f90.
integer, dimension(:), allocatable flint_base::flint_class::eventoptions |
Event options for each event separately. Must be initialized by the Init Function. Specifies what event related info to return to the user.
Definition at line 184 of file FLINT_base.f90.
logical flint_base::flint_class::eventson = .FALSE. |
Events checking off by default.
Definition at line 172 of file FLINT_base.f90.
real(wp), dimension(:), allocatable flint_base::flint_class::eventstepsz |
This specifies the maximum time after which the event function must be called for each event detection. A 0 value in EventStepSz means that event will be checked after the integrator takes a natural step always. Otherwise, events will be checked after the interval specified by EventStepSz if that interval size is smaller than the length of the current step size. In case it is bigger, then event will be checked at the integrator's natural step boundary. Setting is for each event spearately.
Definition at line 180 of file FLINT_base.f90.
integer flint_base::flint_class::fcalls = 0 |
Total number of Diff Eq. function calls made.
Definition at line 147 of file FLINT_base.f90.
real(wp) flint_base::flint_class::h0 |
Initial step-size that was accepted.
Definition at line 191 of file FLINT_base.f90.
real(wp) flint_base::flint_class::hf |
The accepted step-size that was used at the last step.
Definition at line 192 of file FLINT_base.f90.
real(wp) flint_base::flint_class::initialstepsize |
Initial step size.
Definition at line 168 of file FLINT_base.f90.
logical flint_base::flint_class::interpon = .FALSE. |
True for dense output.
Definition at line 151 of file FLINT_base.f90.
integer, dimension(:), allocatable flint_base::flint_class::interpstates |
States between 1 and n for dense output.
Definition at line 152 of file FLINT_base.f90.
logical flint_base::flint_class::isscalartol = .TRUE. |
Tolerance related constants: rtol and atol Dimension can be either 1 or n (dimension of the initial condition vector). These are used to compute the scale factor (SC) needed for error computations using the following formula:
where is the absolute value of the current solution vector.
Definition at line 139 of file FLINT_base.f90.
integer flint_base::flint_class::maxsteps = INT_MAXSTEPS |
Maximum number of steps integrator is allowed to take.
Definition at line 149 of file FLINT_base.f90.
real(wp) flint_base::flint_class::maxstepsize |
Maximum step-size of the integrator. It is always a positive number. Default value is the difference between the initial and final value of the independent variable.
Definition at line 166 of file FLINT_base.f90.
real(wp) flint_base::flint_class::minstepsize |
Minimum step size of the integrator. It is always a positive number. Default value is used if user do not specifify it.
Definition at line 161 of file FLINT_base.f90.
class(diffeqsys), pointer, public flint_base::flint_class::pdiffeqsys => null() |
Pointer to the user-supplied Differential Equation object.
Definition at line 131 of file FLINT_base.f90.
integer flint_base::flint_class::rejectedsteps = 0 |
Number of rejected integrator steps.
Definition at line 145 of file FLINT_base.f90.
real(wp), dimension(:), allocatable flint_base::flint_class::rtol |
Definition at line 140 of file FLINT_base.f90.
integer(kind(flint_success)), public flint_base::flint_class::status = FLINT_SUCCESS |
Status of the last operation completed.
Definition at line 128 of file FLINT_base.f90.
real(wp), dimension(6) flint_base::flint_class::stepszparams |
Step size method-specific tuning parameters.
Definition at line 170 of file FLINT_base.f90.
integer flint_base::flint_class::totalsteps = 0 |
Total number of integrator steps.
Definition at line 143 of file FLINT_base.f90.
real(wp) flint_base::flint_class::x0 |
Initial value of the independent variable.
Definition at line 189 of file FLINT_base.f90.
real(wp) flint_base::flint_class::xf |
Final value of the independent variable.
Definition at line 190 of file FLINT_base.f90.
real(wp), dimension(:), allocatable flint_base::flint_class::xint |
Array of the natural integrator steps.
Definition at line 154 of file FLINT_base.f90.
real(wp), dimension(:), allocatable flint_base::flint_class::y0 |
Initial condition.
Definition at line 193 of file FLINT_base.f90.
real(wp), dimension(:), allocatable flint_base::flint_class::yf |
Final solution.
Definition at line 194 of file FLINT_base.f90.