1   /* Copyright 2013-2017 CS Systèmes d'Information
2    * Licensed to CS Systèmes d'Information (CS) under one or more
3    * contributor license agreements.  See the NOTICE file distributed with
4    * this work for additional information regarding copyright ownership.
5    * CS licenses this file to You under the Apache License, Version 2.0
6    * (the "License"); you may not use this file except in compliance with
7    * the License.  You may obtain a copy of the License at
8    *
9    *   http://www.apache.org/licenses/LICENSE-2.0
10   *
11   * Unless required by applicable law or agreed to in writing, software
12   * distributed under the License is distributed on an "AS IS" BASIS,
13   * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14   * See the License for the specific language governing permissions and
15   * limitations under the License.
16   */
17  package org.orekit.rugged.linesensor;
18  
19  import java.util.ArrayList;
20  import java.util.List;
21  import java.util.stream.Stream;
22  
23  import org.hipparchus.analysis.UnivariateFunction;
24  import org.hipparchus.analysis.solvers.BracketingNthOrderBrentSolver;
25  import org.hipparchus.analysis.solvers.UnivariateSolver;
26  import org.hipparchus.exception.MathIllegalArgumentException;
27  import org.hipparchus.exception.MathIllegalStateException;
28  import org.hipparchus.geometry.euclidean.threed.Vector3D;
29  import org.hipparchus.linear.Array2DRowRealMatrix;
30  import org.hipparchus.linear.ArrayRealVector;
31  import org.hipparchus.linear.DecompositionSolver;
32  import org.hipparchus.linear.MatrixUtils;
33  import org.hipparchus.linear.QRDecomposition;
34  import org.hipparchus.linear.RealMatrix;
35  import org.hipparchus.linear.RealVector;
36  import org.hipparchus.linear.SingularValueDecomposition;
37  import org.hipparchus.util.FastMath;
38  import org.hipparchus.util.Precision;
39  import org.orekit.frames.Transform;
40  import org.orekit.rugged.errors.RuggedException;
41  import org.orekit.rugged.errors.RuggedExceptionWrapper;
42  import org.orekit.rugged.utils.SpacecraftToObservedBody;
43  import org.orekit.time.AbsoluteDate;
44  import org.orekit.utils.Constants;
45  import org.orekit.utils.PVCoordinates;
46  
47  /** Class dedicated to find when ground point crosses mean sensor plane.
48   * <p>
49   * This class is used in the first stage of inverse location.
50   * </p>
51   * @author Luc Maisonobe
52   */
53  public class SensorMeanPlaneCrossing {
54  
55      /** Number of cached results for guessing the line number. */
56      private static final int CACHED_RESULTS = 6;
57  
58      /** Converter between spacecraft and body. */
59      private final SpacecraftToObservedBody scToBody;
60  
61      /** Middle line. */
62      private final double midLine;
63  
64      /** Body to inertial transform for middle line. */
65      private final Transform midBodyToInert;
66  
67      /** Spacecraft to inertial transform for middle line. */
68      private final Transform midScToInert;
69  
70      /** Minimum line number in the search interval. */
71      private final int minLine;
72  
73      /** Maximum line number in the search interval. */
74      private final int maxLine;
75  
76      /** Flag for light time correction. */
77      private final boolean lightTimeCorrection;
78  
79      /** Flag for aberration of light correction. */
80      private final boolean aberrationOfLightCorrection;
81  
82      /** Line sensor. */
83      private final LineSensor sensor;
84  
85      /** Sensor mean plane normal. */
86      private final Vector3D meanPlaneNormal;
87  
88      /** Maximum number of evaluations. */
89      private final int maxEval;
90  
91      /** Accuracy to use for finding crossing line number. */
92      private final double accuracy;
93  
94      /** Cached results. */
95      private final List<CrossingResult> cachedResults;
96  
97      /** Simple constructor.
98       * @param sensor sensor to consider
99       * @param scToBody converter between spacecraft and body
100      * @param minLine minimum line number
101      * @param maxLine maximum line number
102      * @param lightTimeCorrection flag for light time correction
103      * @param aberrationOfLightCorrection flag for aberration of light correction.
104      * @param maxEval maximum number of evaluations
105      * @param accuracy accuracy to use for finding crossing line number
106      * @exception RuggedException if some frame conversion fails
107      */
108     public SensorMeanPlaneCrossing(final LineSensor sensor,
109                                    final SpacecraftToObservedBody scToBody,
110                                    final int minLine, final int maxLine,
111                                    final boolean lightTimeCorrection,
112                                    final boolean aberrationOfLightCorrection,
113                                    final int maxEval, final double accuracy)
114         throws RuggedException {
115         this(sensor, scToBody, minLine, maxLine, lightTimeCorrection, aberrationOfLightCorrection,
116              maxEval, accuracy, computeMeanPlaneNormal(sensor, minLine, maxLine),
117              Stream.<CrossingResult>empty());
118     }
119 
120     /** Simple constructor.
121      * @param sensor sensor to consider
122      * @param scToBody converter between spacecraft and body
123      * @param minLine minimum line number
124      * @param maxLine maximum line number
125      * @param lightTimeCorrection flag for light time correction
126      * @param aberrationOfLightCorrection flag for aberration of light correction.
127      * @param maxEval maximum number of evaluations
128      * @param accuracy accuracy to use for finding crossing line number
129      * @param meanPlaneNormal mean plane normal
130      * @param cachedResults cached results
131      * @exception RuggedException if some frame conversion fails
132      */
133     public SensorMeanPlaneCrossing(final LineSensor sensor,
134                                    final SpacecraftToObservedBody scToBody,
135                                    final int minLine, final int maxLine,
136                                    final boolean lightTimeCorrection,
137                                    final boolean aberrationOfLightCorrection,
138                                    final int maxEval, final double accuracy,
139                                    final Vector3D meanPlaneNormal,
140                                    final Stream<CrossingResult> cachedResults)
141         throws RuggedException {
142 
143         this.sensor                      = sensor;
144         this.minLine                     = minLine;
145         this.maxLine                     = maxLine;
146         this.lightTimeCorrection         = lightTimeCorrection;
147         this.aberrationOfLightCorrection = aberrationOfLightCorrection;
148         this.maxEval                     = maxEval;
149         this.accuracy                    = accuracy;
150         this.scToBody                    = scToBody;
151 
152         this.midLine                     = 0.5 * (minLine + maxLine);
153         final AbsoluteDate midDate       = sensor.getDate(midLine);
154         this.midBodyToInert              = scToBody.getBodyToInertial(midDate);
155         this.midScToInert                = scToBody.getScToInertial(midDate);
156 
157         this.meanPlaneNormal             = meanPlaneNormal;
158 
159         this.cachedResults               = new ArrayList<CrossingResult>(CACHED_RESULTS);
160         cachedResults.forEach(crossingResult -> {
161             if (crossingResult != null && this.cachedResults.size() < CACHED_RESULTS) {
162                 this.cachedResults.add(crossingResult);
163             }
164         });
165 
166     }
167 
168     /** Compute the plane containing origin that best fits viewing directions point cloud.
169      * @param sensor line sensor
170      * @param minLine minimum line number
171      * @param maxLine maximum line number
172      * <p>
173      * The normal is oriented such that traversing pixels in increasing indices
174      * order corresponds to trigonometric order (i.e. counterclockwise).
175      * </p>
176      * @return normal of the mean plane
177      * @exception RuggedException if mid date cannot be handled
178      */
179     private static Vector3D computeMeanPlaneNormal(final LineSensor sensor, final int minLine, final int maxLine)
180         throws RuggedException {
181 
182         final AbsoluteDate midDate = sensor.getDate(0.5 * (minLine + maxLine));
183 
184         // build a centered data matrix
185         // (for each viewing direction, we add both the direction and its
186         //  opposite, thus ensuring the plane will contain origin)
187         final RealMatrix matrix = MatrixUtils.createRealMatrix(3, 2 * sensor.getNbPixels());
188         for (int i = 0; i < sensor.getNbPixels(); ++i) {
189             final Vector3D l = sensor.getLOS(midDate, i);
190             matrix.setEntry(0, 2 * i,      l.getX());
191             matrix.setEntry(1, 2 * i,      l.getY());
192             matrix.setEntry(2, 2 * i,      l.getZ());
193             matrix.setEntry(0, 2 * i + 1, -l.getX());
194             matrix.setEntry(1, 2 * i + 1, -l.getY());
195             matrix.setEntry(2, 2 * i + 1, -l.getZ());
196         }
197 
198         // compute Singular Value Decomposition
199         final SingularValueDecomposition svd = new SingularValueDecomposition(matrix);
200 
201         // extract the left singular vector corresponding to least singular value
202         // (i.e. last vector since Hipparchus returns the values
203         //  in non-increasing order)
204         final Vector3D singularVector = new Vector3D(svd.getU().getColumn(2)).normalize();
205 
206         // check rotation order
207         final Vector3D first = sensor.getLOS(midDate, 0);
208         final Vector3D last  = sensor.getLOS(midDate, sensor.getNbPixels() - 1);
209         if (Vector3D.dotProduct(singularVector, Vector3D.crossProduct(first, last)) >= 0) {
210             return singularVector;
211         } else {
212             return singularVector.negate();
213         }
214 
215     }
216 
217     /** Get the underlying sensor.
218      * @return underlying sensor
219      */
220     public LineSensor getSensor() {
221         return sensor;
222     }
223 
224     /** Get converter between spacecraft and body.
225      * @return converter between spacecraft and body
226      */
227     public SpacecraftToObservedBody getScToBody() {
228         return scToBody;
229     }
230 
231     /** Get the minimum line number in the search interval.
232      * @return minimum line number in the search interval
233      */
234     public int getMinLine() {
235         return minLine;
236     }
237 
238     /** Get the maximum line number in the search interval.
239      * @return maximum line number in the search interval
240      */
241     public int getMaxLine() {
242         return maxLine;
243     }
244 
245     /** Get the maximum number of evaluations.
246      * @return maximum number of evaluations
247      */
248     public int getMaxEval() {
249         return maxEval;
250     }
251 
252     /** Get the accuracy to use for finding crossing line number.
253      * @return accuracy to use for finding crossing line number
254      */
255     public double getAccuracy() {
256         return accuracy;
257     }
258 
259     /** Get the mean plane normal.
260      * <p>
261      * The normal is oriented such traversing pixels in increasing indices
262      * order corresponds is consistent with trigonometric order (i.e.
263      * counterclockwise).
264      * </p>
265      * @return mean plane normal
266      */
267     public Vector3D getMeanPlaneNormal() {
268         return meanPlaneNormal;
269     }
270 
271     /** Get cached previous results.
272      * @return cached previous results
273      */
274     public Stream<CrossingResult> getCachedResults() {
275         return cachedResults.stream();
276     }
277 
278     /** Container for mean plane crossing result. */
279     public static class CrossingResult {
280 
281         /** Crossing date. */
282         private final AbsoluteDate crossingDate;
283 
284         /** Crossing line. */
285         private final double crossingLine;
286 
287         /** Target ground point. */
288         private final Vector3D target;
289 
290         /** Target direction in spacecraft frame. */
291         private final Vector3D targetDirection;
292 
293         /** Derivative of the target direction in spacecraft frame with respect to line number. */
294         private final Vector3D targetDirectionDerivative;
295 
296         /** Simple constructor.
297          * @param crossingDate crossing date
298          * @param crossingLine crossing line
299          * @param target target ground point
300          * @param targetDirection target direction in spacecraft frame
301          * @param targetDirectionDerivative derivative of the target direction in spacecraft frame with respect to line number
302          */
303         public CrossingResult(final AbsoluteDate crossingDate, final double crossingLine,
304                               final Vector3D target,
305                               final Vector3D targetDirection,
306                               final Vector3D targetDirectionDerivative) {
307             this.crossingDate              = crossingDate;
308             this.crossingLine              = crossingLine;
309             this.target                    = target;
310             this.targetDirection           = targetDirection;
311             this.targetDirectionDerivative = targetDirectionDerivative;
312         }
313 
314         /** Get the crossing date.
315          * @return crossing date
316          */
317         public AbsoluteDate getDate() {
318             return crossingDate;
319         }
320 
321         /** Get the crossing line.
322          * @return crossing line
323          */
324         public double getLine() {
325             return crossingLine;
326         }
327 
328         /** Get the target ground point.
329          * @return target ground point
330          */
331         public Vector3D getTarget() {
332             return target;
333         }
334 
335         /** Get the normalized target direction in spacecraft frame at crossing.
336          * @return normalized target direction in spacecraft frame at crossing
337          * with respect to line number
338          */
339         public Vector3D getTargetDirection() {
340             return targetDirection;
341         }
342 
343         /** Get the derivative of the normalized target direction in spacecraft frame at crossing.
344          * @return derivative of the normalized target direction in spacecraft frame at crossing
345          * with respect to line number
346          * @since 2.0
347          */
348         public Vector3D getTargetDirectionDerivative() {
349             return targetDirectionDerivative;
350         }
351 
352     }
353 
354     /** Find mean plane crossing.
355      * @param target target ground point
356      * @return line number and target direction at mean plane crossing,
357      * or null if search interval does not bracket a solution
358      * @exception RuggedException if geometry cannot be computed for some line or
359      * if the maximum number of evaluations is exceeded
360      */
361     public CrossingResult find(final Vector3D target)
362         throws RuggedException {
363 
364         double crossingLine     = midLine;
365         Transform bodyToInert   = midBodyToInert;
366         Transform scToInert     = midScToInert;
367 
368         // count the number of available results
369         if (cachedResults.size() >= 4) {
370             // we already have computed at lest 4 values, we attempt to build a linear
371             // model to guess a better start line
372             final double guessedCrossingLine = guessStartLine(target);
373             if (guessedCrossingLine >= minLine && guessedCrossingLine <= maxLine) {
374                 crossingLine = guessedCrossingLine;
375                 final AbsoluteDate date = sensor.getDate(crossingLine);
376                 bodyToInert = scToBody.getBodyToInertial(date);
377                 scToInert   = scToBody.getScToInertial(date);
378             }
379         }
380 
381         final PVCoordinates targetPV = new PVCoordinates(target, Vector3D.ZERO);
382 
383         // we don't use an Hipparchus solver here because we are more
384         // interested in reducing the number of evaluations than being accurate,
385         // as we know the solution is improved in the second stage of inverse location.
386         // We expect two or three evaluations only. Each new evaluation shows up quickly in
387         // the performances as it involves frames conversions
388         final double[]  crossingLineHistory = new double[maxEval];
389         final double[]  betaHistory         = new double[maxEval];
390         final double[]  betaDerHistory      = new double[maxEval];
391         boolean         atMin               = false;
392         boolean         atMax               = false;
393         for (int i = 0; i < maxEval; ++i) {
394 
395             crossingLineHistory[i] = crossingLine;
396             final Vector3D[] targetDirection =
397                     evaluateLine(crossingLine, targetPV, bodyToInert, scToInert);
398             betaHistory[i] = FastMath.acos(Vector3D.dotProduct(targetDirection[0], meanPlaneNormal));
399             final double s = Vector3D.dotProduct(targetDirection[1], meanPlaneNormal);
400             betaDerHistory[i] = -s / FastMath.sqrt(1 - s * s);
401 
402             final double deltaL;
403             if (i == 0) {
404                 // simple Newton iteration
405                 deltaL        = (0.5 * FastMath.PI - betaHistory[i]) / betaDerHistory[i];
406                 crossingLine += deltaL;
407             } else {
408                 // inverse cubic iteration
409                 final double a0    = betaHistory[i - 1] - 0.5 * FastMath.PI;
410                 final double l0    = crossingLineHistory[i - 1];
411                 final double d0    = betaDerHistory[i - 1];
412                 final double a1    = betaHistory[i]     - 0.5 * FastMath.PI;
413                 final double l1    = crossingLineHistory[i];
414                 final double d1    = betaDerHistory[i];
415                 final double a1Ma0 = a1 - a0;
416                 crossingLine = ((l0 * (a1 - 3 * a0) - a0 * a1Ma0 / d0) * a1 * a1 +
417                                 (l1 * (3 * a1 - a0) - a1 * a1Ma0 / d1) * a0 * a0
418                                ) / (a1Ma0 * a1Ma0 * a1Ma0);
419                 deltaL = crossingLine - l1;
420             }
421             if (FastMath.abs(deltaL) <= accuracy) {
422                 // return immediately, without doing any additional evaluation!
423                 final CrossingResult crossingResult =
424                     new CrossingResult(sensor.getDate(crossingLine), crossingLine, target,
425                                        targetDirection[0], targetDirection[1]);
426                 boolean isNew = true;
427                 for (final CrossingResult existing : cachedResults) {
428                     isNew = isNew && FastMath.abs(crossingLine - existing.crossingLine) > accuracy;
429                 }
430                 if (isNew) {
431                     // this result is different from the existing ones,
432                     // it brings new sampling data to the cache
433                     if (cachedResults.size() >= CACHED_RESULTS) {
434                         cachedResults.remove(cachedResults.size() - 1);
435                     }
436                     cachedResults.add(0, crossingResult);
437                 }
438                 return crossingResult;
439             }
440             for (int j = 0; j < i; ++j) {
441                 if (FastMath.abs(crossingLine - crossingLineHistory[j]) <= 1.0) {
442                     // rare case: we are stuck in a loop!
443                     // switch to a more robust (but slower) algorithm in this case
444                     final CrossingResult slowResult = slowFind(targetPV, crossingLine);
445                     if (slowResult == null) {
446                         return null;
447                     }
448                     if (cachedResults.size() >= CACHED_RESULTS) {
449                         cachedResults.remove(cachedResults.size() - 1);
450                     }
451                     cachedResults.add(0, slowResult);
452                     return cachedResults.get(0);
453                 }
454             }
455 
456             if (crossingLine < minLine) {
457                 if (atMin) {
458                     // we were already trying at minLine and we need to go below that
459                     // give up as the solution is out of search interval
460                     return null;
461                 }
462                 atMin        = true;
463                 crossingLine = minLine;
464             } else if (crossingLine > maxLine) {
465                 if (atMax) {
466                     // we were already trying at maxLine and we need to go above that
467                     // give up as the solution is out of search interval
468                     return null;
469                 }
470                 atMax        = true;
471                 crossingLine = maxLine;
472             } else {
473                 // the next evaluation will be a regular point
474                 atMin = false;
475                 atMax = false;
476             }
477 
478             final AbsoluteDate date = sensor.getDate(crossingLine);
479             bodyToInert = scToBody.getBodyToInertial(date);
480             scToInert   = scToBody.getScToInertial(date);
481         }
482 
483         return null;
484 
485     }
486 
487     /** Guess a start line using the last four results.
488      * @param target target ground point
489      * @return guessed start line
490      */
491     private double guessStartLine(final Vector3D target) {
492         try {
493 
494             // assume a linear model of the form: l = ax + by + cz + d
495             final int        n = cachedResults.size();
496             final RealMatrix m = new Array2DRowRealMatrix(n, 4);
497             final RealVector v = new ArrayRealVector(n);
498             for (int i = 0; i < n; ++i) {
499                 final CrossingResult crossingResult = cachedResults.get(i);
500                 m.setEntry(i, 0, crossingResult.getTarget().getX());
501                 m.setEntry(i, 1, crossingResult.getTarget().getY());
502                 m.setEntry(i, 2, crossingResult.getTarget().getZ());
503                 m.setEntry(i, 3, 1.0);
504                 v.setEntry(i, crossingResult.getLine());
505             }
506 
507             final DecompositionSolver solver = new QRDecomposition(m, Precision.SAFE_MIN).getSolver();
508             final RealVector coefficients = solver.solve(v);
509 
510             // apply the linear model
511             return target.getX() * coefficients.getEntry(0) +
512                    target.getY() * coefficients.getEntry(1) +
513                    target.getZ() * coefficients.getEntry(2) +
514                    coefficients.getEntry(3);
515 
516         } catch (MathIllegalStateException sme) {
517             // the points are not independent
518             return Double.NaN;
519         }
520     }
521 
522     /** Find mean plane crossing using a slow but robust method.
523      * @param targetPV target ground point
524      * @param initialGuess initial guess for the crossing line
525      * @return line number and target direction at mean plane crossing,
526      * or null if search interval does not bracket a solution
527      * @exception RuggedException if geometry cannot be computed for some line or
528      * if the maximum number of evaluations is exceeded
529      */
530     private CrossingResult slowFind(final PVCoordinates targetPV, final double initialGuess)
531         throws RuggedException {
532         try {
533 
534             // safety check
535             final double startValue;
536             if (initialGuess < minLine || initialGuess > maxLine) {
537                 startValue = midLine;
538             } else {
539                 startValue = initialGuess;
540             }
541 
542             final UnivariateSolver solver = new BracketingNthOrderBrentSolver(accuracy, 5);
543             double crossingLine = solver.solve(maxEval, new UnivariateFunction() {
544                 /** {@inheritDoc} */
545                 @Override
546                 public double value(final double x) throws RuggedExceptionWrapper {
547                     try {
548                         final AbsoluteDate date = sensor.getDate(x);
549                         final Vector3D[] targetDirection =
550                                 evaluateLine(x, targetPV, scToBody.getBodyToInertial(date), scToBody.getScToInertial(date));
551                         return 0.5 * FastMath.PI - FastMath.acos(Vector3D.dotProduct(targetDirection[0], meanPlaneNormal));
552                     } catch (RuggedException re) {
553                         throw new RuggedExceptionWrapper(re);
554                     }
555                 }
556             }, minLine, maxLine, startValue);
557 
558             final AbsoluteDate date = sensor.getDate(crossingLine);
559             final Vector3D[] targetDirection =
560                     evaluateLine(crossingLine, targetPV, scToBody.getBodyToInertial(date), scToBody.getScToInertial(date));
561             return new CrossingResult(sensor.getDate(crossingLine), crossingLine, targetPV.getPosition(),
562                                       targetDirection[0], targetDirection[1]);
563 
564         } catch (MathIllegalArgumentException nbe) {
565             return null;
566         } catch (RuggedExceptionWrapper rew) {
567             throw rew.getException();
568         }
569     }
570 
571     /** Evaluate geometry for a given line number.
572      * @param lineNumber current line number
573      * @param targetPV target ground point
574      * @param bodyToInert transform from observed body to inertial frame, for current line
575      * @param scToInert transform from inertial frame to spacecraft frame, for current line
576      * @return target direction in spacecraft frame, with its first derivative
577      * with respect to line number
578      */
579     private Vector3D[] evaluateLine(final double lineNumber, final PVCoordinates targetPV,
580                                     final Transform bodyToInert, final Transform scToInert) {
581 
582         // compute the transform between spacecraft and observed body
583         final PVCoordinates refInert =
584                 scToInert.transformPVCoordinates(new PVCoordinates(sensor.getPosition(), Vector3D.ZERO));
585 
586         final PVCoordinates targetInert;
587         if (lightTimeCorrection) {
588             // apply light time correction
589             final Vector3D iT     = bodyToInert.transformPosition(targetPV.getPosition());
590             final double   deltaT = refInert.getPosition().distance(iT) / Constants.SPEED_OF_LIGHT;
591             targetInert           = bodyToInert.shiftedBy(-deltaT).transformPVCoordinates(targetPV);
592         } else {
593             // don't apply light time correction
594             targetInert = bodyToInert.transformPVCoordinates(targetPV);
595         }
596 
597         final PVCoordinates lInert    = new PVCoordinates(refInert, targetInert);
598         final Transform     inertToSc = scToInert.getInverse();
599         final Vector3D obsLInert;
600         final Vector3D obsLInertDot;
601         if (aberrationOfLightCorrection) {
602 
603             // apply aberration of light correction
604             // as the spacecraft velocity is small with respect to speed of light,
605             // we use classical velocity addition and not relativistic velocity addition
606             // we have: c * lInert + vsat = k * obsLInert
607             final PVCoordinates spacecraftPV = scToInert.transformPVCoordinates(PVCoordinates.ZERO);
608             final Vector3D  l        = lInert.getPosition().normalize();
609             final Vector3D  lDot     = normalizedDot(lInert.getPosition(), lInert.getVelocity());
610             final Vector3D  kObs     = new Vector3D(Constants.SPEED_OF_LIGHT, l, +1.0, spacecraftPV.getVelocity());
611             obsLInert = kObs.normalize();
612 
613             // the following derivative is computed under the assumption the spacecraft velocity
614             // is constant in inertial frame ... It is obviously not true, but as this velocity
615             // is very small with respect to speed of light, the error is expected to remain small
616             obsLInertDot = normalizedDot(kObs, new Vector3D(Constants.SPEED_OF_LIGHT, lDot));
617 
618         } else {
619 
620             // don't apply aberration of light correction
621             obsLInert    = lInert.getPosition().normalize();
622             obsLInertDot = normalizedDot(lInert.getPosition(), lInert.getVelocity());
623 
624         }
625         final Vector3D dir    = inertToSc.transformVector(obsLInert);
626         final Vector3D dirDot = new Vector3D(+1.0, inertToSc.transformVector(obsLInertDot),
627                                              -1.0, Vector3D.crossProduct(inertToSc.getRotationRate(), dir));
628 
629         // combine vector value and derivative
630         final double rate = sensor.getRate(lineNumber);
631         return new Vector3D[] {
632             dir, dirDot.scalarMultiply(1.0 / rate)
633         };
634 
635     }
636 
637     /** Compute the derivative of normalized vector.
638      * @param u base vector
639      * @param uDot derivative of the base vector
640      * @return vDot, where v = u / ||u||
641      */
642     private Vector3D normalizedDot(final Vector3D u, final Vector3D uDot) {
643         final double n = u.getNorm();
644         return new Vector3D(1.0 / n, uDot, -Vector3D.dotProduct(u, uDot) / (n * n * n), u);
645     }
646 
647 }