Class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
- java.lang.Object
-
- org.orekit.orbits.FieldOrbit<T>
-
- All Implemented Interfaces:
FieldTimeInterpolable<FieldOrbit<T>,T>
,FieldTimeShiftable<FieldOrbit<T>,T>
,FieldTimeStamped<T>
,FieldPVCoordinatesProvider<T>
- Direct Known Subclasses:
FieldCartesianOrbit
,FieldCircularOrbit
,FieldEquinoctialOrbit
,FieldKeplerianOrbit
public abstract class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>> extends Object implements FieldPVCoordinatesProvider<T>, FieldTimeStamped<T>, FieldTimeShiftable<FieldOrbit<T>,T>, FieldTimeInterpolable<FieldOrbit<T>,T>
This class handles orbital parameters.For user convenience, both the Cartesian and the equinoctial elements are provided by this class, regardless of the canonical representation implemented in the derived class (which may be classical Keplerian elements for example).
The parameters are defined in a frame specified by the user. It is important to make sure this frame is consistent: it probably is inertial and centered on the central body. This information is used for example by some force models.
Instance of this class are guaranteed to be immutable.
- Since:
- 9.0
- Author:
- Luc Maisonobe, Guylaine Prat, Fabien Maussion, Véronique Pommier-Maurussane
-
-
Constructor Summary
Constructors Modifier Constructor Description protected
FieldOrbit(Frame frame, FieldAbsoluteDate<T> date, double mu)
Default constructor.protected
FieldOrbit(TimeStampedFieldPVCoordinates<T> FieldPVCoordinates, Frame frame, double mu)
Set the orbit from Cartesian parameters.
-
Method Summary
All Methods Static Methods Instance Methods Abstract Methods Concrete Methods Modifier and Type Method Description abstract void
addKeplerContribution(PositionAngle type, double gm, T[] pDot)
Add the contribution of the Keplerian motion to parameters derivativesprotected abstract T[][]
computeJacobianEccentricWrtCartesian()
Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.protected abstract T[][]
computeJacobianMeanWrtCartesian()
Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.protected abstract T[][]
computeJacobianTrueWrtCartesian()
Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.protected void
fillHalfRow(T a, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v, T[] row, int j)
Fill a Jacobian half row with a single vector.protected void
fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.protected void
fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.protected void
fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T a4, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4, T[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.protected void
fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T a4, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4, T a5, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v5, T[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.protected void
fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T a4, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4, T a5, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v5, T a6, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v6, T[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.abstract T
getA()
Get the semi-major axis.abstract T
getADot()
Get the semi-major axis derivative.FieldAbsoluteDate<T>
getDate()
Get the date of orbital parameters.abstract T
getE()
Get the eccentricity.abstract T
getEDot()
Get the eccentricity derivative.abstract T
getEquinoctialEx()
Get the first component of the equinoctial eccentricity vector.abstract T
getEquinoctialExDot()
Get the first component of the equinoctial eccentricity vector.abstract T
getEquinoctialEy()
Get the second component of the equinoctial eccentricity vector.abstract T
getEquinoctialEyDot()
Get the second component of the equinoctial eccentricity vector.Frame
getFrame()
Get the frame in which the orbital parameters are defined.abstract T
getHx()
Get the first component of the inclination vector.abstract T
getHxDot()
Get the first component of the inclination vector derivative.abstract T
getHy()
Get the second component of the inclination vector.abstract T
getHyDot()
Get the second component of the inclination vector derivative.abstract T
getI()
Get the inclination.abstract T
getIDot()
Get the inclination derivative.void
getJacobianWrtCartesian(PositionAngle type, T[][] jacobian)
Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.void
getJacobianWrtParameters(PositionAngle type, T[][] jacobian)
Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.T
getKeplerianMeanMotion()
Get the Keplerian mean motion.T
getKeplerianPeriod()
Get the Keplerian period.abstract T
getLE()
Get the eccentric longitude argument.abstract T
getLEDot()
Get the eccentric longitude argument derivative.abstract T
getLM()
Get the mean longitude argument.abstract T
getLMDot()
Get the mean longitude argument derivative.abstract T
getLv()
Get the true longitude argument.abstract T
getLvDot()
Get the true longitude argument derivative.double
getMu()
TimeStampedFieldPVCoordinates<T>
getPVCoordinates()
Get theTimeStampedPVCoordinates
in definition frame.TimeStampedFieldPVCoordinates<T>
getPVCoordinates(Frame outputFrame)
Get theTimeStampedPVCoordinates
in a specified frame.TimeStampedFieldPVCoordinates<T>
getPVCoordinates(FieldAbsoluteDate<T> otherDate, Frame otherFrame)
Get theFieldPVCoordinates
of the body in the selected frame.abstract OrbitType
getType()
Get the orbit type.abstract boolean
hasDerivatives()
Check if orbit includes derivatives.protected static <T extends org.hipparchus.RealFieldElement<T>>
booleanhasNonKeplerianAcceleration(FieldPVCoordinates<T> pva, double mu)
Check if Cartesian coordinates include non-Keplerian acceleration.protected abstract TimeStampedFieldPVCoordinates<T>
initPVCoordinates()
Compute the position/velocity coordinates from the canonical parameters.abstract FieldOrbit<T>
shiftedBy(T dt)
Get a time-shifted orbit.abstract Orbit
toOrbit()
Transforms the FieldOrbit instance into an Orbit instance.-
Methods inherited from class java.lang.Object
clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
-
Methods inherited from interface org.orekit.time.FieldTimeInterpolable
interpolate, interpolate
-
Methods inherited from interface org.orekit.time.FieldTimeShiftable
shiftedBy
-
-
-
-
Constructor Detail
-
FieldOrbit
protected FieldOrbit(Frame frame, FieldAbsoluteDate<T> date, double mu) throws IllegalArgumentException
Default constructor. Build a new instance with arbitrary default elements.- Parameters:
frame
- the frame in which the parameters are defined (must be apseudo-inertial frame
)date
- date of the orbital parametersmu
- central attraction coefficient (m^3/s^2)- Throws:
IllegalArgumentException
- if frame is not apseudo-inertial frame
-
FieldOrbit
protected FieldOrbit(TimeStampedFieldPVCoordinates<T> FieldPVCoordinates, Frame frame, double mu) throws IllegalArgumentException
Set the orbit from Cartesian parameters.The acceleration provided in
FieldPVCoordinates
is accessible usinggetPVCoordinates()
andgetPVCoordinates(Frame)
. All other methods usemu
and the position to compute the acceleration, includingshiftedBy(RealFieldElement)
andgetPVCoordinates(FieldAbsoluteDate, Frame)
.- Parameters:
FieldPVCoordinates
- the position and velocity in the inertial frameframe
- the frame in which theTimeStampedPVCoordinates
are defined (must be apseudo-inertial frame
)mu
- central attraction coefficient (m^3/s^2)- Throws:
IllegalArgumentException
- if frame is not apseudo-inertial frame
-
-
Method Detail
-
hasNonKeplerianAcceleration
protected static <T extends org.hipparchus.RealFieldElement<T>> boolean hasNonKeplerianAcceleration(FieldPVCoordinates<T> pva, double mu)
Check if Cartesian coordinates include non-Keplerian acceleration.- Type Parameters:
T
- type of the field elements- Parameters:
pva
- Cartesian coordinatesmu
- central attraction coefficient- Returns:
- true if Cartesian coordinates include non-Keplerian acceleration
-
getType
public abstract OrbitType getType()
Get the orbit type.- Returns:
- orbit type
-
getFrame
public Frame getFrame()
Get the frame in which the orbital parameters are defined.- Returns:
- frame in which the orbital parameters are defined
-
toOrbit
public abstract Orbit toOrbit()
Transforms the FieldOrbit instance into an Orbit instance.- Returns:
- Orbit instance with same properties
-
getA
public abstract T getA()
Get the semi-major axis.Note that the semi-major axis is considered negative for hyperbolic orbits.
- Returns:
- semi-major axis (m)
-
getADot
public abstract T getADot()
Get the semi-major axis derivative.Note that the semi-major axis is considered negative for hyperbolic orbits.
If the orbit was created without derivatives, the value returned is null.
- Returns:
- semi-major axis derivative (m/s)
-
getEquinoctialEx
public abstract T getEquinoctialEx()
Get the first component of the equinoctial eccentricity vector.- Returns:
- first component of the equinoctial eccentricity vector
-
getEquinoctialExDot
public abstract T getEquinoctialExDot()
Get the first component of the equinoctial eccentricity vector.If the orbit was created without derivatives, the value returned is null.
- Returns:
- first component of the equinoctial eccentricity vector
-
getEquinoctialEy
public abstract T getEquinoctialEy()
Get the second component of the equinoctial eccentricity vector.- Returns:
- second component of the equinoctial eccentricity vector
-
getEquinoctialEyDot
public abstract T getEquinoctialEyDot()
Get the second component of the equinoctial eccentricity vector.If the orbit was created without derivatives, the value returned is null.
- Returns:
- second component of the equinoctial eccentricity vector
-
getHx
public abstract T getHx()
Get the first component of the inclination vector.- Returns:
- first component of the inclination vector
-
getHxDot
public abstract T getHxDot()
Get the first component of the inclination vector derivative.If the orbit was created without derivatives, the value returned is null.
- Returns:
- first component of the inclination vector derivative
-
getHy
public abstract T getHy()
Get the second component of the inclination vector.- Returns:
- second component of the inclination vector
-
getHyDot
public abstract T getHyDot()
Get the second component of the inclination vector derivative.- Returns:
- second component of the inclination vector derivative
-
getLE
public abstract T getLE()
Get the eccentric longitude argument.- Returns:
- E + ω + Ω eccentric longitude argument (rad)
-
getLEDot
public abstract T getLEDot()
Get the eccentric longitude argument derivative.If the orbit was created without derivatives, the value returned is null.
- Returns:
- d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s)
-
getLv
public abstract T getLv()
Get the true longitude argument.- Returns:
- v + ω + Ω true longitude argument (rad)
-
getLvDot
public abstract T getLvDot()
Get the true longitude argument derivative.If the orbit was created without derivatives, the value returned is null.
- Returns:
- d(v + ω + Ω)/dt true longitude argument derivative (rad/s)
-
getLM
public abstract T getLM()
Get the mean longitude argument.- Returns:
- M + ω + Ω mean longitude argument (rad)
-
getLMDot
public abstract T getLMDot()
Get the mean longitude argument derivative.If the orbit was created without derivatives, the value returned is null.
- Returns:
- d(M + ω + Ω)/dt mean longitude argument derivative (rad/s)
-
getE
public abstract T getE()
Get the eccentricity.- Returns:
- eccentricity
-
getEDot
public abstract T getEDot()
Get the eccentricity derivative.If the orbit was created without derivatives, the value returned is null.
- Returns:
- eccentricity derivative
-
getI
public abstract T getI()
Get the inclination.If the orbit was created without derivatives, the value returned is null.
- Returns:
- inclination (rad)
-
getIDot
public abstract T getIDot()
Get the inclination derivative.- Returns:
- inclination derivative (rad/s)
-
hasDerivatives
public abstract boolean hasDerivatives()
Check if orbit includes derivatives.- Returns:
- true if orbit includes derivatives
- Since:
- 9.0
- See Also:
getADot()
,getEquinoctialExDot()
,getEquinoctialEyDot()
,getHxDot()
,getHyDot()
,getLEDot()
,getLvDot()
,getLMDot()
,getEDot()
,getIDot()
-
getMu
public double getMu()
-
getKeplerianPeriod
public T getKeplerianPeriod()
Get the Keplerian period.The Keplerian period is computed directly from semi major axis and central acceleration constant.
- Returns:
- Keplerian period in seconds, or positive infinity for hyperbolic orbits
-
getKeplerianMeanMotion
public T getKeplerianMeanMotion()
Get the Keplerian mean motion.The Keplerian mean motion is computed directly from semi major axis and central acceleration constant.
- Returns:
- Keplerian mean motion in radians per second
-
getDate
public FieldAbsoluteDate<T> getDate()
Get the date of orbital parameters.- Specified by:
getDate
in interfaceFieldTimeStamped<T extends org.hipparchus.RealFieldElement<T>>
- Returns:
- date of the orbital parameters
-
getPVCoordinates
public TimeStampedFieldPVCoordinates<T> getPVCoordinates(Frame outputFrame)
Get theTimeStampedPVCoordinates
in a specified frame.- Parameters:
outputFrame
- frame in which the position/velocity coordinates shall be computed- Returns:
- FieldPVCoordinates in the specified output frame
- See Also:
getPVCoordinates()
-
getPVCoordinates
public TimeStampedFieldPVCoordinates<T> getPVCoordinates(FieldAbsoluteDate<T> otherDate, Frame otherFrame)
Get theFieldPVCoordinates
of the body in the selected frame.- Specified by:
getPVCoordinates
in interfaceFieldPVCoordinatesProvider<T extends org.hipparchus.RealFieldElement<T>>
- Parameters:
otherDate
- current dateotherFrame
- the frame where to define the position- Returns:
- time-stamped position/velocity of the body (m and m/s)
-
getPVCoordinates
public TimeStampedFieldPVCoordinates<T> getPVCoordinates()
Get theTimeStampedPVCoordinates
in definition frame.- Returns:
- FieldPVCoordinates in the definition frame
- See Also:
getPVCoordinates(Frame)
-
initPVCoordinates
protected abstract TimeStampedFieldPVCoordinates<T> initPVCoordinates()
Compute the position/velocity coordinates from the canonical parameters.- Returns:
- computed position/velocity coordinates
-
shiftedBy
public abstract FieldOrbit<T> shiftedBy(T dt)
Get a time-shifted orbit.The orbit can be slightly shifted to close dates. This shift is based on a simple Keplerian model. It is not intended as a replacement for proper orbit and attitude propagation but should be sufficient for small time shifts or coarse accuracy.
- Specified by:
shiftedBy
in interfaceFieldTimeShiftable<FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>,T extends org.hipparchus.RealFieldElement<T>>
- Parameters:
dt
- time shift in seconds- Returns:
- a new orbit, shifted with respect to the instance (which is immutable)
-
getJacobianWrtCartesian
public void getJacobianWrtCartesian(PositionAngle type, T[][] jacobian)
Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.Element
jacobian[i][j]
is the derivative of parameter i of the orbit with respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.- Parameters:
type
- type of the position angle to usejacobian
- placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix is larger than 6x6, only the 6x6 upper left corner will be modified
-
getJacobianWrtParameters
public void getJacobianWrtParameters(PositionAngle type, T[][] jacobian)
Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.Element
jacobian[i][j]
is the derivative of Cartesian coordinate i of the orbit with respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters.- Parameters:
type
- type of the position angle to usejacobian
- placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix is larger than 6x6, only the 6x6 upper left corner will be modified
-
computeJacobianMeanWrtCartesian
protected abstract T[][] computeJacobianMeanWrtCartesian()
Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.Element
jacobian[i][j]
is the derivative of parameter i of the orbit with respect to Cartesian coordinate j. This means each row correspond to one orbital parameter whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.- Returns:
- 6x6 Jacobian matrix
- See Also:
computeJacobianEccentricWrtCartesian()
,computeJacobianTrueWrtCartesian()
-
computeJacobianEccentricWrtCartesian
protected abstract T[][] computeJacobianEccentricWrtCartesian()
Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.Element
jacobian[i][j]
is the derivative of parameter i of the orbit with respect to Cartesian coordinate j. This means each row correspond to one orbital parameter whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.- Returns:
- 6x6 Jacobian matrix
- See Also:
computeJacobianMeanWrtCartesian()
,computeJacobianTrueWrtCartesian()
-
computeJacobianTrueWrtCartesian
protected abstract T[][] computeJacobianTrueWrtCartesian()
Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.Element
jacobian[i][j]
is the derivative of parameter i of the orbit with respect to Cartesian coordinate j. This means each row correspond to one orbital parameter whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.- Returns:
- 6x6 Jacobian matrix
- See Also:
computeJacobianMeanWrtCartesian()
,computeJacobianEccentricWrtCartesian()
-
addKeplerContribution
public abstract void addKeplerContribution(PositionAngle type, double gm, T[] pDot)
Add the contribution of the Keplerian motion to parameters derivativesThis method is used by integration-based propagators to evaluate the part of Keplerian motion to evolution of the orbital state.
- Parameters:
type
- type of the position angle in the stategm
- attraction coefficient to usepDot
- array containing orbital state derivatives to update (the Keplerian part must be added to the array components, as the array may already contain some non-zero elements corresponding to non-Keplerian parts)
-
fillHalfRow
protected void fillHalfRow(T a, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v, T[] row, int j)
Fill a Jacobian half row with a single vector.- Parameters:
a
- coefficient of the vectorv
- vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
fillHalfRow
protected void fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.- Parameters:
a1
- coefficient of the first vectorv1
- first vectora2
- coefficient of the second vectorv2
- second vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
fillHalfRow
protected void fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.- Parameters:
a1
- coefficient of the first vectorv1
- first vectora2
- coefficient of the second vectorv2
- second vectora3
- coefficient of the third vectorv3
- third vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
fillHalfRow
protected void fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T a4, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4, T[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.- Parameters:
a1
- coefficient of the first vectorv1
- first vectora2
- coefficient of the second vectorv2
- second vectora3
- coefficient of the third vectorv3
- third vectora4
- coefficient of the fourth vectorv4
- fourth vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
fillHalfRow
protected void fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T a4, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4, T a5, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v5, T[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.- Parameters:
a1
- coefficient of the first vectorv1
- first vectora2
- coefficient of the second vectorv2
- second vectora3
- coefficient of the third vectorv3
- third vectora4
- coefficient of the fourth vectorv4
- fourth vectora5
- coefficient of the fifth vectorv5
- fifth vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
fillHalfRow
protected void fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T a4, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4, T a5, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v5, T a6, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v6, T[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.- Parameters:
a1
- coefficient of the first vectorv1
- first vectora2
- coefficient of the second vectorv2
- second vectora3
- coefficient of the third vectorv3
- third vectora4
- coefficient of the fourth vectorv4
- fourth vectora5
- coefficient of the fifth vectorv5
- fifth vectora6
- coefficient of the sixth vectorv6
- sixth vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
-