FLINT 0.9.9
Fortran Library for numerical INTegration of differential equations
|
FLINT Base Module. More...
Data Types | |
interface | defunc |
Interface of the user-supplied Differential Equation function. It must return a vector of double precision values with dimension same as that of the input vector Y. More... | |
type | diffeqsys |
User must extend DiffEqSys and provide the DiffEq and Events function. More... | |
interface | eventfunc |
Interface for a user-supplied procedure for computing events values during the integration. It is also called after an event is detected in order for the user to specify the event actions to be taken by the event handler. Note that if multiple events are triggered at the same location then only the event with the lowest index will be reported. More... | |
type | flint_class |
Abstract base class that must be inherited by all the numerical integrator classes. More... | |
interface | info |
Query method to get the current state information of FLINT. Can be called after each Integration to get the integration stats. More... | |
interface | init |
Interface for the Initializaition method. It must be called before Integration. More... | |
interface | integrate |
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... | |
interface | interpolate |
Method for delayed interpolation. It must be called after at least one invocation of the Init and Integrate method with the InterpOn option enabled. It can be called by the user any number of times with different points at which the solution is sought. FLINT computes the interpolation coefficients during integration and stores them internally at each natural step of the integrator. During interpolation, these coefficients are used to compute the solutions quickly. This is much less computationally intensive than using the Integrate method. More... | |
Enumerations | |
enum | { flint_success = 0 , flint_event_term = 1 , flint_stiff_problem = 2 , flint_error = 3 , flint_error_params = 4 , flint_error_maxsteps = 5 , flint_error_memalloc = 6 , flint_error_memdealloc = 7 , flint_error_dimension = 8 , flint_error_stepsz_toosmall = 9 , flint_error_intpstates = 10 , flint_error_intp_array = 11 , flint_error_intp_off = 12 , flint_error_init_reqd = 13 , flint_error_eventparams = 14 , flint_error_eventroot = 15 , flint_error_conststepsz = 16 } |
Supported Error Codes. More... | |
enum | { flint_eventaction_continue = 0 , flint_eventaction_terminate = 1 , flint_eventaction_changesoly = 2 , flint_eventaction_restartint = 4 , flint_eventaction_mask = 128 } |
Supported FLINT Event Actions. More... | |
enum | { flint_eventoption_rootfinding = 0 , flint_eventoption_stepbegin = 1 , flint_eventoption_stepend = 2 } |
Supported FLINT Event Options. More... | |
Variables | |
integer, parameter | int_maxsteps = 1024 |
Default number of maximum steps that intgerator is allowed to take. | |
character(len= *), dimension(0:16), parameter | flint_errorstrings = [character(len=64) :: 'Successfully Completed', 'Termination Event Occurred', 'Problem appears to be stiff', 'Unknown error', 'Incorrect user-provided parameters', 'Maximum integration steps reached', 'Memory allocation error', 'Memory deallocation error', 'Incorrect dimensions', 'Step size reduced to very small value', 'Incorrect interpolation states', 'Incorrect interpolation array provided', 'Interpolation is not enabled', 'Init is required', 'Incorrect event-related parameters ', 'Event root finding failed', 'Incorrect step size for the non-adaptive option' ] |
Error message strings for the supported error codes. | |
FLINT Base Module.
It includes the abstract intefaces for every solver in FLINT in addition to the intefaces for the differential equation and event functions provided by the user. FLINT error codes are defined here as well.
anonymous enum |
Supported Error Codes.
Definition at line 43 of file FLINT_base.f90.
anonymous enum |
Supported FLINT Event Actions.
Definition at line 86 of file FLINT_base.f90.
anonymous enum |
Supported FLINT Event Options.
Definition at line 107 of file FLINT_base.f90.
character(len=*), dimension(0:16), parameter flint_base::flint_errorstrings = [character(len=64) :: 'Successfully Completed', 'Termination Event Occurred', 'Problem appears to be stiff', 'Unknown error', 'Incorrect user-provided parameters', 'Maximum integration steps reached', 'Memory allocation error', 'Memory deallocation error', 'Incorrect dimensions', 'Step size reduced to very small value', 'Incorrect interpolation states', 'Incorrect interpolation array provided', 'Interpolation is not enabled', 'Init is required', 'Incorrect event-related parameters ', 'Event root finding failed', 'Incorrect step size for the non-adaptive option' ] |
Error message strings for the supported error codes.
Definition at line 64 of file FLINT_base.f90.
integer, parameter flint_base::int_maxsteps = 1024 |
Default number of maximum steps that intgerator is allowed to take.
Definition at line 40 of file FLINT_base.f90.