FLINT 0.9.9
Fortran Library for numerical INTegration of differential equations
Loading...
Searching...
No Matches
flint_base::flint_class Type Referenceabstract

Abstract base class that must be inherited by all the numerical integrator classes. More...

+ Inheritance diagram for flint_base::flint_class:
+ Collaboration diagram for flint_base::flint_class:

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.
 

Detailed Description

Abstract base class that must be inherited by all the numerical integrator classes.

Definition at line 125 of file FLINT_base.f90.

Member Function/Subroutine Documentation

◆ info()

procedure (info), deferred flint_base::flint_class::info ( class(flint_class), intent(inout) me,
integer(kind(flint_success)), intent(out) laststatus,
character(len=:), intent(out), optional, allocatable statusmsg,
integer, intent(out), optional nsteps,
integer, intent(out), optional naccept,
integer, intent(out), optional nreject,
integer, intent(out), optional nfcalls,
logical, intent(out), optional interpready,
real(wp), intent(out), optional h0,
real(wp), intent(out), optional x0,
real(wp), dimension(:), intent(out), optional y0,
real(wp), intent(out), optional hf,
real(wp), intent(out), optional xf,
real(wp), dimension(:), intent(out), optional yf )
pure virtual

Current FLINT object query function.

Parameters
[in,out]meObject of class type FLINT
[out]laststatusThe status code of the last operation completed
[out]statusmsgStatus message corresponding to the status code in LastStatus
[out]nstepsTotal number of steps taken
[out]nacceptNumber of Accepted steps
[out]nrejectNumber of rejected steps
[out]nfcallsFunction calls made not including event calls
[out]interpreadyWhether Interpolation is enabled
[out]h0The initial accepted step-size
[out]x0The initial value of the independent variable
[out]y0Initial condition
[out]hfFinal accepted step-size
[out]xfThe final value of the independent variable
[out]yfFinal solution

Definition at line 209 of file FLINT_base.f90.

◆ init()

procedure (init), deferred flint_base::flint_class::init ( class(flint_class), intent(inout) me,
class(diffeqsys), target de,
integer, intent(in) maxsteps,
integer, intent(in), optional method,
real(wp), dimension(:), intent(in), optional atol,
real(wp), dimension(:), intent(in), optional rtol,
logical, intent(in), optional interpon,
integer, dimension(:), intent(in), optional interpstates,
real(wp), intent(in), optional minstepsize,
real(wp), intent(in), optional maxstepsize,
real(wp), dimension(6), intent(in), optional stepszparams,
logical, intent(in), optional eventson,
real(wp), dimension(:), intent(in), optional eventstepsz,
integer(kind(flint_eventoption_rootfinding)), dimension(:), intent(in), optional eventoptions,
real(wp), dimension(:), intent(in), optional eventtol )
pure virtual

Init is called once and only to specify new differential equaions, events, and a few integrator options.

Parameters
[in,out]meObject
deDifferential Equation system object
[in]maxstepsMust 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]methodType of Intergrator to use. Implemented methods are:
  • ERK_DOP853 (default)
  • ERK_DOP54
  • ERK_VERNER98R
  • ERK_VERNER65E
[in]atolAbsolute tolerance, size must be either 1 or n
[in]rtolRelative tolerance, size must be either 1 or n
[in]interponIf 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]interpstatesSpecify the interpolation states for which the interpolation coefficients will be computed and stored, e.g. [1,3,4,6,n].
[in]minstepsizeUser can specify minimum and maximum step sizes allowed. Note that these are positive quantities. Sign will be handled bu FLINT.
[in]maxstepsizeUser can specify minimum and maximum step sizes allowed. Note that these are positive quantities. Sign will be handled bu FLINT.
[in]stepszparamsStep-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:
  • StepSzParams(1): Safety-factor
  • StepSzParams(2): Minimum safety-factor
  • StepSzParams(3): Maximum safety-factor
  • StepSzParams(4): Lund Stabilization parameter beta
  • StepSzParams(5): Lund Stabilization parameter beta multiplier
  • StepSzParams(6): Lund Stabilization PI control term initial value
[in]eventsonSet true for events detection
[in]eventstepszIf 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]eventoptionsEvent Options that control the event information returned by FLINT Default is FLINT_EVENTOPTION_ROOTFINDING.
[in]eventtolEvent tolerance values used for root-finding if enabled, size must be m

Definition at line 200 of file FLINT_base.f90.

◆ integrate()

procedure (integrate), deferred flint_base::flint_class::integrate ( class(flint_class), intent(inout) me,
real(wp), intent(in) x0,
real(wp), dimension(me%pdiffeqsys%n), intent(in) y0,
real(wp), intent(inout) xf,
real(wp), dimension(me%pdiffeqsys%n), intent(out) yf,
real(wp), intent(inout) stepsz,
logical, intent(in), optional useconststepsz,
logical, intent(in), optional stepadvance,
logical, intent(in), optional intstepson,
real(wp), dimension(:), intent(out), optional, allocatable xint,
real(wp), dimension(:,:), intent(out), optional, allocatable yint,
logical, dimension(me%pdiffeqsys%m), intent(in), optional eventmask,
real(wp), dimension(:,:), intent(out), optional, allocatable eventstates,
integer, intent(inout), optional stifftest,
real(wp), dimension(:), intent(in), optional params )
pure virtual

Integrate the solution from the initial to final value of the independent variable.

Parameters
[in,out]meObject of the class type FLINT
[in]x0Initial value of the independent variable
[in]y0Initial conditions
[in,out]xfFinal value of the independent variable
[out]yfFinal Solution
[in,out]stepszUser 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]useconststepszIf 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]stepadvanceIf true, FLINT will only advance the integration by one step and return. By default, this option is set to False.
[in]intstepsonIf true, then FLINT will return the solution at the integrator's natural step size.
[out]xintIf 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]yintIf 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]eventmaskIf 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]eventstatesSolutions 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]stifftestFLINT 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.
  • StiffTest == 1: If DiffEq are detected to be stiff then -1 is returned in StiffTest and the integration halts at the current step with appropriate status code.
  • StiffTest == 2: If DiffEq are detected to be stiff then -1 is returned in StiffTest but the integration continues.
[in]paramsParameters for sensitivity analysis. Not implemented yet.

Definition at line 203 of file FLINT_base.f90.

◆ interpolate()

procedure (interpolate), deferred flint_base::flint_class::interpolate ( class(flint_class), intent(inout) me,
real(wp), dimension(1:), intent(in) xarr,
real(wp), dimension(size(me%interpstates),size(xarr)), intent(out) yarr,
logical, intent(in), optional saveinterpcoeffs )
pure virtual

If interpolation is enabled, then it outputs the solution at the user-specified grid points.

Parameters
[in,out]meObject of class type FLINT
[in]xarrXarr 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]yarrYarr provides the solution at each point specified in Xarr.
[in]saveinterpcoeffsIf 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.

Member Data Documentation

◆ acceptedsteps

integer flint_base::flint_class::acceptedsteps = 0

Number of accepted integrator steps.

Definition at line 144 of file FLINT_base.f90.

◆ atol

real(wp), dimension(:), allocatable flint_base::flint_class::atol

Definition at line 141 of file FLINT_base.f90.

◆ bip

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.

◆ etol

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.

◆ eventoptions

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.

◆ eventson

logical flint_base::flint_class::eventson = .FALSE.

Events checking off by default.

Definition at line 172 of file FLINT_base.f90.

◆ eventstepsz

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.

◆ fcalls

integer flint_base::flint_class::fcalls = 0

Total number of Diff Eq. function calls made.

Definition at line 147 of file FLINT_base.f90.

◆ h0

real(wp) flint_base::flint_class::h0

Initial step-size that was accepted.

Definition at line 191 of file FLINT_base.f90.

◆ hf

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.

◆ initialstepsize

real(wp) flint_base::flint_class::initialstepsize

Initial step size.

Definition at line 168 of file FLINT_base.f90.

◆ interpon

logical flint_base::flint_class::interpon = .FALSE.

True for dense output.

Definition at line 151 of file FLINT_base.f90.

◆ interpstates

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.

◆ isscalartol

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:

\[ SC_i = atol_i + rtol_i |Y_{i}| \quad i=1..n,\]

where $|Y|$ is the absolute value of the current solution vector.

Definition at line 139 of file FLINT_base.f90.

◆ maxsteps

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.

◆ maxstepsize

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.

◆ minstepsize

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.

◆ pdiffeqsys

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.

◆ rejectedsteps

integer flint_base::flint_class::rejectedsteps = 0

Number of rejected integrator steps.

Definition at line 145 of file FLINT_base.f90.

◆ rtol

real(wp), dimension(:), allocatable flint_base::flint_class::rtol

Definition at line 140 of file FLINT_base.f90.

◆ status

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.

◆ stepszparams

real(wp), dimension(6) flint_base::flint_class::stepszparams

Step size method-specific tuning parameters.

Definition at line 170 of file FLINT_base.f90.

◆ totalsteps

integer flint_base::flint_class::totalsteps = 0

Total number of integrator steps.

Definition at line 143 of file FLINT_base.f90.

◆ x0

real(wp) flint_base::flint_class::x0

Initial value of the independent variable.

Definition at line 189 of file FLINT_base.f90.

◆ xf

real(wp) flint_base::flint_class::xf

Final value of the independent variable.

Definition at line 190 of file FLINT_base.f90.

◆ xint

real(wp), dimension(:), allocatable flint_base::flint_class::xint

Array of the natural integrator steps.

Definition at line 154 of file FLINT_base.f90.

◆ y0

real(wp), dimension(:), allocatable flint_base::flint_class::y0

Initial condition.

Definition at line 193 of file FLINT_base.f90.

◆ yf

real(wp), dimension(:), allocatable flint_base::flint_class::yf

Final solution.

Definition at line 194 of file FLINT_base.f90.


The documentation for this type was generated from the following file: