FLINT 0.9.9
Fortran Library for numerical INTegration of differential equations
Loading...
Searching...
No Matches
flint_base::integrate Interface Reference

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...

+ Collaboration diagram for flint_base::integrate:

Public Member Functions

subroutine integrate (me, x0, y0, xf, yf, stepsz, useconststepsz, stepadvance, intstepson, xint, yint, eventmask, eventstates, stifftest, params)
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ integrate()

subroutine flint_base::integrate::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 )
virtual
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 400 of file FLINT_base.f90.


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