T
- the type of the field elementspublic abstract class AdaptiveStepsizeFieldIntegrator<T extends RealFieldElement<T>> extends AbstractFieldIntegrator<T>
These algorithms perform integration with stepsize control, which means the user does not specify the integration step but rather a tolerance on error. The error threshold is computed as
threshold_i = absTol_i + relTol_i * max (abs (ym), abs (ym+1))where absTol_i is the absolute tolerance for component i of the state vector and relTol_i is the relative tolerance for the same component. The user can also use only two scalar values absTol and relTol which will be used for all components.
Note that only the main part
of the state vector is used for stepsize control. The secondary parts
of the state
vector are explicitly ignored for stepsize control.
If the estimated error for ym+1 is such that
sqrt((sum (errEst_i / threshold_i)^2 ) / n) < 1(where n is the main set dimension) then the step is accepted, otherwise the step is rejected and a new attempt is made with a new stepsize.
Modifier and Type | Field and Description |
---|---|
protected int |
mainSetDimension
Main set dimension.
|
protected double |
scalAbsoluteTolerance
Allowed absolute scalar error.
|
protected double |
scalRelativeTolerance
Allowed relative scalar error.
|
protected double[] |
vecAbsoluteTolerance
Allowed absolute vectorial error.
|
protected double[] |
vecRelativeTolerance
Allowed relative vectorial error.
|
Constructor and Description |
---|
AdaptiveStepsizeFieldIntegrator(Field<T> field,
String name,
double minStep,
double maxStep,
double[] vecAbsoluteTolerance,
double[] vecRelativeTolerance)
Build an integrator with the given stepsize bounds.
|
AdaptiveStepsizeFieldIntegrator(Field<T> field,
String name,
double minStep,
double maxStep,
double scalAbsoluteTolerance,
double scalRelativeTolerance)
Build an integrator with the given stepsize bounds.
|
Modifier and Type | Method and Description |
---|---|
protected T |
filterStep(T h,
boolean forward,
boolean acceptSmall)
Filter the integration step.
|
T |
getMaxStep()
Get the maximal step.
|
T |
getMinStep()
Get the minimal step.
|
T |
initializeStep(boolean forward,
int order,
T[] scale,
FieldODEStateAndDerivative<T> state0,
FieldEquationsMapper<T> mapper)
Initialize the integration step.
|
protected void |
resetInternalState()
Reset internal state to dummy values.
|
protected void |
sanityChecks(FieldODEState<T> initialState,
T t)
Check the integration span.
|
void |
setInitialStepSize(T initialStepSize)
Set the initial step size.
|
void |
setStepSizeControl(double minimalStep,
double maximalStep,
double[] absoluteTolerance,
double[] relativeTolerance)
Set the adaptive step size control parameters.
|
void |
setStepSizeControl(double minimalStep,
double maximalStep,
double absoluteTolerance,
double relativeTolerance)
Set the adaptive step size control parameters.
|
acceptStep, addEventHandler, addEventHandler, addStepHandler, clearEventHandlers, clearStepHandlers, computeDerivatives, getCurrentSignedStepsize, getEquations, getEvaluations, getEvaluationsCounter, getEventHandlers, getField, getMaxEvaluations, getName, getStepHandlers, getStepSize, getStepStart, initIntegration, isLastStep, resetOccurred, setIsLastStep, setMaxEvaluations, setStateInitialized, setStepSize, setStepStart
clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
integrate
protected double scalAbsoluteTolerance
protected double scalRelativeTolerance
protected double[] vecAbsoluteTolerance
protected double[] vecRelativeTolerance
protected int mainSetDimension
public AdaptiveStepsizeFieldIntegrator(Field<T> field, String name, double minStep, double maxStep, double scalAbsoluteTolerance, double scalRelativeTolerance)
field
- field to which the time and state vector elements belongname
- name of the methodminStep
- 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 errorpublic AdaptiveStepsizeFieldIntegrator(Field<T> field, String name, double minStep, double maxStep, double[] vecAbsoluteTolerance, double[] vecRelativeTolerance)
field
- field to which the time and state vector elements belongname
- name of the methodminStep
- 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 errorpublic void setStepSizeControl(double minimalStep, double maximalStep, double absoluteTolerance, double relativeTolerance)
A side effect of this method is to also reset the initial
step so it will be automatically computed by the integrator
if setInitialStepSize
is not called by the user.
minimalStep
- minimal step (must be positive even for backward
integration), the last step can be smaller than thismaximalStep
- maximal step (must be positive even for backward
integration)absoluteTolerance
- allowed absolute errorrelativeTolerance
- allowed relative errorpublic void setStepSizeControl(double minimalStep, double maximalStep, double[] absoluteTolerance, double[] relativeTolerance)
A side effect of this method is to also reset the initial
step so it will be automatically computed by the integrator
if setInitialStepSize
is not called by the user.
minimalStep
- minimal step (must be positive even for backward
integration), the last step can be smaller than thismaximalStep
- maximal step (must be positive even for backward
integration)absoluteTolerance
- allowed absolute errorrelativeTolerance
- allowed relative errorpublic void setInitialStepSize(T initialStepSize)
This method allows the user to specify an initial positive step size instead of letting the integrator guess it by itself. If this method is not called before integration is started, the initial step size will be estimated by the integrator.
initialStepSize
- initial step size to use (must be positive even
for backward integration ; providing a negative value or a value
outside of the min/max step interval will lead the integrator to
ignore the value and compute the initial step size by itself)protected void sanityChecks(FieldODEState<T> initialState, T t) throws MathIllegalArgumentException
sanityChecks
in class AbstractFieldIntegrator<T extends RealFieldElement<T>>
initialState
- initial statet
- target time for the integrationMathIllegalArgumentException
- if integration span is too smallpublic T initializeStep(boolean forward, int order, T[] scale, FieldODEStateAndDerivative<T> state0, FieldEquationsMapper<T> mapper) throws MathIllegalArgumentException, MathIllegalStateException
forward
- forward integration indicatororder
- order of the methodscale
- scaling vector for the state vector (can be shorter than state vector)state0
- state at integration start timemapper
- mapper for all the equationsMathIllegalStateException
- if the number of functions evaluations is exceededMathIllegalArgumentException
- if arrays dimensions do not match equations settingsprotected T filterStep(T h, boolean forward, boolean acceptSmall) throws MathIllegalArgumentException
h
- signed stepforward
- forward integration indicatoracceptSmall
- if true, steps smaller than the minimal value
are silently increased up to this value, if false such small
steps generate an exceptionMathIllegalArgumentException
- if the step is too small and acceptSmall is falseprotected void resetInternalState()
public T getMinStep()
public T getMaxStep()
Copyright © 2016 Hipparchus.org. All rights reserved.