Class AdamsBashforthFieldIntegrator<T extends CalculusFieldElement<T>>
- java.lang.Object
-
- org.hipparchus.ode.AbstractFieldIntegrator<T>
-
- org.hipparchus.ode.nonstiff.AdaptiveStepsizeFieldIntegrator<T>
-
- org.hipparchus.ode.MultistepFieldIntegrator<T>
-
- org.hipparchus.ode.nonstiff.AdamsFieldIntegrator<T>
-
- org.hipparchus.ode.nonstiff.AdamsBashforthFieldIntegrator<T>
-
- Type Parameters:
T
- the type of the field elements
- All Implemented Interfaces:
FieldODEIntegrator<T>
public class AdamsBashforthFieldIntegrator<T extends CalculusFieldElement<T>> extends AdamsFieldIntegrator<T>
This class implements explicit Adams-Bashforth integrators for Ordinary Differential Equations.Adams-Bashforth methods (in fact due to Adams alone) are explicit multistep ODE solvers. This implementation is a variation of the classical one: it uses adaptive stepsize to implement error control, whereas classical implementations are fixed step size. The value of state vector at step n+1 is a simple combination of the value at step n and of the derivatives at steps n, n-1, n-2 ... Depending on the number k of previous steps one wants to use for computing the next value, different formulas are available:
- k = 1: yn+1 = yn + h y'n
- k = 2: yn+1 = yn + h (3y'n-y'n-1)/2
- k = 3: yn+1 = yn + h (23y'n-16y'n-1+5y'n-2)/12
- k = 4: yn+1 = yn + h (55y'n-59y'n-1+37y'n-2-9y'n-3)/24
- ...
A k-steps Adams-Bashforth method is of order k.
There must be sufficient time for the
starter integrator
to take several steps between the the last reset event, and the end of integration, otherwise an exception may be thrown during integration. The user can adjust the end date of integration, or the step size of the starter integrator to ensure a sufficient number of steps can be completed before the end of integration.Implementation details
We define scaled derivatives si(n) at step n as: \[ \left\{\begin{align} s_1(n) &= h y'_n \text{ for first derivative}\\ s_2(n) &= \frac{h^2}{2} y_n'' \text{ for second derivative}\\ s_3(n) &= \frac{h^3}{6} y_n''' \text{ for third derivative}\\ &\cdots\\ s_k(n) &= \frac{h^k}{k!} y_n^{(k)} \text{ for } k^\mathrm{th} \text{ derivative} \end{align}\right. \]
The definitions above use the classical representation with several previous first derivatives. Lets define \[ q_n = [ s_1(n-1) s_1(n-2) \ldots s_1(n-(k-1)) ]^T \] (we omit the k index in the notation for clarity). With these definitions, Adams-Bashforth methods can be written: \[ \left\{\begin{align} k = 1: & y_{n+1} = y_n + s_1(n) \\ k = 2: & y_{n+1} = y_n + \frac{3}{2} s_1(n) + [ \frac{-1}{2} ] q_n \\ k = 3: & y_{n+1} = y_n + \frac{23}{12} s_1(n) + [ \frac{-16}{12} \frac{5}{12} ] q_n \\ k = 4: & y_{n+1} = y_n + \frac{55}{24} s_1(n) + [ \frac{-59}{24} \frac{37}{24} \frac{-9}{24} ] q_n \\ & \cdots \end{align}\right. \]
Instead of using the classical representation with first derivatives only (yn, s1(n) and qn), our implementation uses the Nordsieck vector with higher degrees scaled derivatives all taken at the same step (yn, s1(n) and rn) where rn is defined as: \[ r_n = [ s_2(n), s_3(n) \ldots s_k(n) ]^T \] (here again we omit the k index in the notation for clarity)
Taylor series formulas show that for any index offset i, s1(n-i) can be computed from s1(n), s2(n) ... sk(n), the formula being exact for degree k polynomials. \[ s_1(n-i) = s_1(n) + \sum_{j\gt 0} (j+1) (-i)^j s_{j+1}(n) \] The previous formula can be used with several values for i to compute the transform between classical representation and Nordsieck vector. The transform between rn and qn resulting from the Taylor series formulas above is: \[ q_n = s_1(n) u + P r_n \] where u is the [ 1 1 ... 1 ]T vector and P is the (k-1)×(k-1) matrix built with the \((j+1) (-i)^j\) terms with i being the row number starting from 1 and j being the column number starting from 1: \[ P=\begin{bmatrix} -2 & 3 & -4 & 5 & \ldots \\ -4 & 12 & -32 & 80 & \ldots \\ -6 & 27 & -108 & 405 & \ldots \\ -8 & 48 & -256 & 1280 & \ldots \\ & & \ldots\\ \end{bmatrix} \]
Using the Nordsieck vector has several advantages:
- it greatly simplifies step interpolation as the interpolator mainly applies Taylor series formulas,
- it simplifies step changes that occur when discrete events that truncate the step are triggered,
- it allows to extend the methods in order to support adaptive stepsize.
The Nordsieck vector at step n+1 is computed from the Nordsieck vector at step n as follows:
- yn+1 = yn + s1(n) + uT rn
- s1(n+1) = h f(tn+1, yn+1)
- rn+1 = (s1(n) - s1(n+1)) P-1 u + P-1 A P rn
where A is a rows shifting matrix (the lower left part is an identity matrix):
[ 0 0 ... 0 0 | 0 ] [ ---------------+---] [ 1 0 ... 0 0 | 0 ] A = [ 0 1 ... 0 0 | 0 ] [ ... | 0 ] [ 0 0 ... 1 0 | 0 ] [ 0 0 ... 0 1 | 0 ]
The P-1u vector and the P-1 A P matrix do not depend on the state, they only depend on k and therefore are precomputed once for all.
-
-
Field Summary
-
Fields inherited from class org.hipparchus.ode.MultistepFieldIntegrator
nordsieck, scaled
-
-
Constructor Summary
Constructors Constructor Description AdamsBashforthFieldIntegrator(Field<T> field, int nSteps, double minStep, double maxStep, double[] vecAbsoluteTolerance, double[] vecRelativeTolerance)
Build an Adams-Bashforth integrator with the given order and step control parameters.AdamsBashforthFieldIntegrator(Field<T> field, int nSteps, double minStep, double maxStep, double scalAbsoluteTolerance, double scalRelativeTolerance)
Build an Adams-Bashforth integrator with the given order and step control parameters.
-
Method Summary
All Methods Instance Methods Concrete Methods Modifier and Type Method Description protected double
errorEstimation(T[] previousState, T predictedTime, T[] predictedState, T[] predictedScaled, FieldMatrix<T> predictedNordsieck)
Estimate error.protected org.hipparchus.ode.nonstiff.AdamsFieldStateInterpolator<T>
finalizeStep(T stepSize, T[] predictedY, T[] predictedScaled, Array2DRowFieldMatrix<T> predictedNordsieck, boolean isForward, FieldODEStateAndDerivative<T> globalPreviousState, FieldODEStateAndDerivative<T> globalCurrentState, FieldEquationsMapper<T> equationsMapper)
Finalize the step.-
Methods inherited from class org.hipparchus.ode.nonstiff.AdamsFieldIntegrator
initializeHighOrderDerivatives, integrate, updateHighOrderDerivativesPhase1, updateHighOrderDerivativesPhase2
-
Methods inherited from class org.hipparchus.ode.MultistepFieldIntegrator
computeStepGrowShrinkFactor, getMaxGrowth, getMinReduction, getNSteps, getSafety, getStarterIntegrator, rescale, setMaxGrowth, setMinReduction, setSafety, setStarterIntegrator, start
-
Methods inherited from class org.hipparchus.ode.nonstiff.AdaptiveStepsizeFieldIntegrator
getMaxStep, getMinStep, getStepSizeHelper, initializeStep, resetInternalState, sanityChecks, setInitialStepSize, setStepSizeControl, setStepSizeControl
-
Methods inherited from class org.hipparchus.ode.AbstractFieldIntegrator
acceptStep, addEventDetector, addStepEndHandler, addStepHandler, clearEventDetectors, clearStepEndHandlers, clearStepHandlers, computeDerivatives, getCurrentSignedStepsize, getEquations, getEvaluations, getEvaluationsCounter, getEventDetectors, getField, getMaxEvaluations, getName, getStepEndHandlers, getStepHandlers, getStepSize, getStepStart, initIntegration, isLastStep, resetOccurred, setIsLastStep, setMaxEvaluations, setStateInitialized, setStepSize, setStepStart
-
-
-
-
Constructor Detail
-
AdamsBashforthFieldIntegrator
public AdamsBashforthFieldIntegrator(Field<T> field, int nSteps, double minStep, double maxStep, double scalAbsoluteTolerance, double scalRelativeTolerance) throws MathIllegalArgumentException
Build an Adams-Bashforth integrator with the given order and step control parameters.- Parameters:
field
- field to which the time and state vector elements belongnSteps
- number of steps of the method excluding the one being computedminStep
- minimal step (sign is irrelevant, regardless of integration direction, forward or backward), the last step can be smaller than thismaxStep
- maximal step (sign is irrelevant, regardless of integration direction, forward or backward), the last step can be smaller than thisscalAbsoluteTolerance
- allowed absolute errorscalRelativeTolerance
- allowed relative error- Throws:
MathIllegalArgumentException
- if order is 1 or less
-
AdamsBashforthFieldIntegrator
public AdamsBashforthFieldIntegrator(Field<T> field, int nSteps, double minStep, double maxStep, double[] vecAbsoluteTolerance, double[] vecRelativeTolerance) throws IllegalArgumentException
Build an Adams-Bashforth integrator with the given order and step control parameters.- Parameters:
field
- field to which the time and state vector elements belongnSteps
- number of steps of the method excluding the one being computedminStep
- minimal step (sign is irrelevant, regardless of integration direction, forward or backward), the last step can be smaller than thismaxStep
- maximal step (sign is irrelevant, regardless of integration direction, forward or backward), the last step can be smaller than thisvecAbsoluteTolerance
- allowed absolute errorvecRelativeTolerance
- allowed relative error- Throws:
IllegalArgumentException
- if order is 1 or less
-
-
Method Detail
-
errorEstimation
protected double errorEstimation(T[] previousState, T predictedTime, T[] predictedState, T[] predictedScaled, FieldMatrix<T> predictedNordsieck)
Estimate error.- Specified by:
errorEstimation
in classAdamsFieldIntegrator<T extends CalculusFieldElement<T>>
- Parameters:
previousState
- state vector at step startpredictedTime
- time at step endpredictedState
- predicted state vector at step endpredictedScaled
- predicted value of the scaled derivatives at step endpredictedNordsieck
- predicted value of the Nordsieck vector at step end- Returns:
- estimated normalized local discretization error
-
finalizeStep
protected org.hipparchus.ode.nonstiff.AdamsFieldStateInterpolator<T> finalizeStep(T stepSize, T[] predictedY, T[] predictedScaled, Array2DRowFieldMatrix<T> predictedNordsieck, boolean isForward, FieldODEStateAndDerivative<T> globalPreviousState, FieldODEStateAndDerivative<T> globalCurrentState, FieldEquationsMapper<T> equationsMapper)
Finalize the step.- Specified by:
finalizeStep
in classAdamsFieldIntegrator<T extends CalculusFieldElement<T>>
- Parameters:
stepSize
- step size used in the scaled and Nordsieck arrayspredictedY
- predicted state at end of steppredictedScaled
- predicted first scaled derivativepredictedNordsieck
- predicted Nordsieck vectorisForward
- integration direction indicatorglobalPreviousState
- start of the global stepglobalCurrentState
- end of the global stepequationsMapper
- mapper for ODE equations primary and secondary components- Returns:
- step interpolator
-
-