FractionalDiffEq.AtanganaSedaCFType
solve(prob::FODESystem, dt, AtanganaSedaCF())

Atagana Seda method for Caputo-Fabrizio fractional order differential equations.

Tip

Used for the Caputo Fabrizio fractional differential operators.

https://doi.org/10.1016/c2020-0-02711-8
FractionalDiffEq.BDFType

BDF

Backward Differentiation Formula generated weights fractional linear multi-steps method.

References

@article{Garrappa2015TrapezoidalMF, title={Trapezoidal methods for fractional differential equations: Theoretical and computational aspects}, author={Roberto Garrappa}, journal={ArXiv}, year={2015}, volume={abs/1912.09878} }

FractionalDiffEq.ClosedFormType

Usage

solve(prob::MultiTermsFODEProblem, ClosedForm())

Use Closed-Form solution to obtain numerical solution at zero initial condition.

Reference:

Dingyu Xue, Northeastern University, China ISBN:9787030543981

FractionalDiffEq.ClosedFormHankelMType

Usage

solve(prob::MultiTermsFODEProblem, ::ClosedFormHankelM)

Use Closed-Form Hankel matrix algorithm to obtain numerical solution at zero initial condition.

FractionalDiffEq.DODEProblemType

Defines an distributed order fractional ordinary differential equation (DODE) problem.

Mathematical Specification of an distributed order FODE problem

To define an multi-terms FODE Problem, you simply need to given the parameters, their correspoding orders, right hand side function and the initial condition $u_0$ which define an FODE:

\[\frac{du^\alpha}{d^{\alpha}t} = f(u,p,t)\]

Distributed order differential equations.

FractionalDiffEq.DOMatrixDiscreteType

Usage

solve(prob, DOMatrixDiscrete())

Use distributed order strip matrix algorithm to solve distriubted order problem.

References

@inproceedings{Jiao2012DistributedOrderDS,
  title={Distributed-Order Dynamic Systems - Stability, Simulation, Applications and Perspectives},
  author={Zhuang Jiao and Yang Quan Chen and Igor Podlubny},
  booktitle={Springer Briefs in Electrical and Computer Engineering},
  year={2012}
}
FractionalDiffEq.DelayABMType
solve(FDDE::FDDEProblem, dt, DelayABM())

Use the Adams-Bashforth-Moulton method to solve fractional delayed differential equations.

References

@inproceedings{Bhalekar2011APS,
  title={A PREDICTOR-CORRECTOR SCHEME FOR SOLVING NONLINEAR DELAY DIFFERENTIAL EQUATIONS OF FRACTIONAL ORDER},
  author={Sachin Bhalekar and Varsha Daftardar-Gejji},
  year={2011}
}
FractionalDiffEq.DelayPECEType

Usage

solve(FDDE::FDDEProblem, dt, DelayPECE())

Using the delayed predictor-corrector method to solve the delayed fractional differential equation problem in the Caputo sense.

Capable of solving both single term FDDE and multiple FDDE, support time varying lags of course😋.

References

@article{Wang2013ANM, title={A Numerical Method for Delayed Fractional-Order Differential Equations}, author={Zhen Wang}, journal={J. Appl. Math.}, year={2013}, volume={2013}, pages={256071:1-256071:7} }

@inproceedings{Nagy2014NUMERICALSF, title={NUMERICAL SIMULATIONS FOR VARIABLE-ORDER FRACTIONAL NONLINEAR DELAY DIFFERENTIAL EQUATIONS}, author={Abdelhameed M. Nagy and Taghreed Abdul Rahman Assiri}, year={2014} }

@inproceedings{Abdelmalek2019APM, title={A Predictor-Corrector Method for Fractional Delay-Differential System with Multiple Lags}, author={Salem Abdelmalek and Redouane Douaifia}, year={2019} }

FractionalDiffEq.DelayPIEXType

DelayPIEX

Explicit product integral method for initial value problems of fractional order differential equations.

FractionalDiffEq.EulerType
Euler

The classical Euler method extended for fractional order differential equations.

FractionalDiffEq.FDDEMatrixProblemType
FDDEMatrixProblem(α, τ, A, B, f, x0, tspan)

Construct a fractional matrix differential equation with delay with general form:

\[D_{t_0}^\alpha\textbf{x}(t)=\textbf{A}(t)\textbf{x}(t)+\textbf{B}(t)\textbf{x}(t-\tau)+\textbf{f}(t)\]

FractionalDiffEq.FDDEProblemType
FDDEProblem(f, ϕ, α, τ, tspan)
  • f: The function describing fractional delay differential equations.
  • ϕ: History function

Construct a fractional delayed differential equation problem.

FractionalDiffEq.FFEODEProblemType
FFEODEProblem(f, α, u0, tspan, p)

Define fractal-fractional differential equations problems with expoenntial decay kernel.

FractionalDiffEq.FFMODEProblemType
FFMODEProblem(f, α, u0, tspan, p)

Define fractal-fractional differential equations problems with Mittag Leffler kernel.

FractionalDiffEq.FFPODEProblemType
FFPODEProblem(f, α, u0, tspan, p)

Define fractal-fractional differential equations problems with power law kernel.

FractionalDiffEq.FODEProblemType

Defines an fractional ordinary differential equation (FODE) problem.

Mathematical Specification of an FODE problem

To define an FODE Problem, you simply need to given the function $f$ and the initial condition $u_0$ which define an FODE:

\[\frac{du^\alpha}{d^{\alpha}t} = f(u,p,t)\]

There are two different ways of specifying f:

  • f(du,u,p,t): in-place. Memory-efficient when avoiding allocations. Best option for most cases unless mutation is not allowed.
  • f(u,p,t): returning du. Less memory-efficient way, particularly suitable when mutation is not allowed (e.g. with certain automatic differentiation packages such as Zygote).
  • order: the fractional order of the differential equations, commensurate and non-commensurate is both supported.

-u₀ should be an AbstractArray (or number) whose geometry matches the desired geometry of u. Note that we are not limited to numbers or vectors for u₀; one is allowed to provide u₀ as arbitrary matrices / higher dimension tensors as well.

Problem Type

Constructors

FODEProblem can be constructed by first building an ODEFunction or by simply passing the FODE right-hand side to the constructor. The constructors are:

  • FODEProblem(f::ODEFunction,u0,tspan,p=NullParameters();kwargs...)
  • FODEProblem{isinplace,specialize}(f,u0,tspan,p=NullParameters();kwargs...) : Defines the FODE with the specified functions. isinplace optionally sets whether the function is inplace or not. This is determined automatically, but not inferred. specialize optionally controls the specialization level. See the specialization levels section of the SciMLBase documentation for more details. The default is AutoSpecialize.

For more details on the in-place and specialization controls, see the ODEFunction documentation.

Parameters are optional, and if not given, then a NullParameters() singleton will be used which will throw nice errors if you try to index non-existent parameters. Any extra keyword arguments are passed on to the solvers. For example, if you set a callback in the problem, then that callback will be added in every solve call.

For specifying Jacobians and mass matrices, see the ODEFunction documentation.

Fields

  • f: The function in the ODE.
  • order: The order of the FODE.
  • u0: The initial condition.
  • tspan: The timespan for the problem.
  • p: The parameters.
  • kwargs: The keyword arguments passed onto the solves.

Example Problem

using SciMLBase
function lorenz!(du,u,p,t)
 du[1] = 10.0(u[2]-u[1])
 du[2] = u[1]*(28.0-u[3]) - u[2]
 du[3] = u[1]*u[2] - (8/3)*u[3]
end
order = [0.96;0.96;0.96]
u0 = [1.0;0.0;0.0]
tspan = (0.0,100.0)
prob = FODEProblem(lorenz!,u0,tspan)

# Test that it worked
using FractionalDiffEq
sol = solve(prob,PIEX())
using Plots; plot(sol,vars=(1,2,3))
FractionalDiffEq.FractionalDiscreteSystemType
FractionalDifferenceSystem(f, α, u0)

Define Fractional Difference System with the general constructure:

\[{^G\nabla_k^\alpha x(k+1)}=f(x(k))\]

With given initial condition $x(i)$.

FractionalDiffEq.MTPECEType
solve(prob::MultiTermsFODEProblem, dt, MTPECE())

Use product integration predictor-corrector method to solve multi-terms FODE.

FractionalDiffEq.MTPIEXType

MTPIEX

Explicit product integral method for initial value problems of fractional order differential equations.

FractionalDiffEq.MTPIRectType
solve(prob::MultiTermsFODEProblem,, MTPIRect(); dt)

Use implicit product integration rectangular type method to solve multi-terms FODE.

FractionalDiffEq.MTPITrapType
solve(prob::MultiTermsFODEProblem, MTPITrap(), dt)

Use implicit product integration trapezoidal type method to solve multi-terms FODE.

FractionalDiffEq.MatrixDiscreteType

MatrixDiscrete

Triangular strip matrices to discrete fractional ordinary differential equations to simple algebra system and solve the system.

References

@inproceedings{Podlubny2000MATRIXAT, title={MATRIX APPROACH TO DISCRETE FRACTIONAL CALCULUS}, author={Igor Podlubny}, year={2000} }

FractionalDiffEq.MatrixFormType

Usage

solve(prob::FDDEMatrixProblem, h, MatrixForm())

Reference

https://github.com/mandresik/system-of-linear-fractional-differential-delayed-equations

FractionalDiffEq.MultiTermsFODEProblemType

Defines an multiple terms linear fractional ordinary differential equation (FODE) problem.

Mathematical Specification of an multi-terms FODE problem

To define an multi-terms FODE Problem, you simply need to given the parameters, their correspoding orders, right hand side function and the initial condition $u_0$ which define an FODE:

\[\frac{du^\alpha}{d^{\alpha}t} = f(u,p,t)\]

Multiple terms fractional order differential equations.

FractionalDiffEq.NewtonGregoryType
solve(prob::FODEProblem, NewtonGregory(); abstol=1e-3, maxiters=1e3)

Newton Gregory generated weights fractional linear multiple steps method.

References

@article{Garrappa2015TrapezoidalMF, title={Trapezoidal methods for fractional differential equations: Theoretical and computational aspects}, author={Roberto Garrappa}, journal={ArXiv}, year={2015}, volume={abs/1912.09878} }

FractionalDiffEq.NewtonPolynomialType
solve(prob::FODEProblem, dt, NewtonPolynomial())

Solve FODE system using Newton Polynomials methods.

Tip

Used for the Caputo Fabrizio fractional differential operators.

References

https://doi.org/10.1016/c2020-0-02711-8

FractionalDiffEq.NonLinearAlgType

Usage

solve(prob::FODEProblem, dt, NonLinearAlg())

Nonlinear algorithm for nonlinear fractional differential equations.

References

Dingyu Xue, Northeastern University, China ISBN:9787030543981

FractionalDiffEq.PECEType
Predictor-Correct scheme.

References

@inproceedings{Garrappa2018NumericalSO, title={Numerical Solution of Fractional Differential Equations: A Survey and a Software Tutorial}, author={Roberto Garrappa}, year={2018} }

FractionalDiffEq.PIEXType

PIEX

Explicit product integral method for initial value problems of fractional order differential equations.

FractionalDiffEq.TrapezoidType
solve(prob::FODEProblem, Trapezoidal(); abstol=1e-3, maxiters=1e3)

Use Trapezoidal with generating function $f(x)=\frac{1+x}{2(1-x)^\alpha}$ generated weights fractional linear multiple steps method to solve system of FODE.

References

@article{Garrappa2015TrapezoidalMF, title={Trapezoidal methods for fractional differential equations: Theoretical and computational aspects}, author={Roberto Garrappa}, journal={ArXiv}, year={2015}, volume={abs/1912.09878} }

FractionalDiffEq.DMethod
D(N, α, dt)

Using D function to construct left hand side equations.

Info

Here N is the size of discrete matrix.

FractionalDiffEq.genfunMethod

P-th precision polynomial generate function

\[g_p(z)=\sum_{k=1}^p \frac{1}{k}(1-z)^k\]

FractionalDiffEq.mittleffMethod
mittleff(α, β, γ, z)

Compute three-parametric mittleff(α, β, γ, z).

@article{2015,
   title={Numerical Evaluation of Two and Three Parameter Mittag-Leffler Functions},
   volume={53},
   ISSN={1095-7170},
   url={http://dx.doi.org/10.1137/140971191},
   DOI={10.1137/140971191},
   number={3},
   journal={SIAM Journal on Numerical Analysis},
   publisher={Society for Industrial & Applied Mathematics (SIAM)},
   author={Garrappa, Roberto},
   year={2015},
   month={Jan},
   pages={1350–1369}
}