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