1   /* Contributed in the public domain.
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.models.earth;
18  
19  import org.hipparchus.Field;
20  import org.hipparchus.RealFieldElement;
21  import org.hipparchus.analysis.RealFieldUnivariateFunction;
22  import org.hipparchus.analysis.UnivariateFunction;
23  import org.hipparchus.analysis.solvers.AllowedSolution;
24  import org.hipparchus.analysis.solvers.BracketingNthOrderBrentSolver;
25  import org.hipparchus.analysis.solvers.FieldBracketingNthOrderBrentSolver;
26  import org.hipparchus.analysis.solvers.UnivariateSolver;
27  import org.hipparchus.exception.MathRuntimeException;
28  import org.hipparchus.geometry.euclidean.threed.FieldLine;
29  import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
30  import org.hipparchus.geometry.euclidean.threed.Line;
31  import org.hipparchus.geometry.euclidean.threed.Vector3D;
32  import org.hipparchus.util.FastMath;
33  import org.orekit.bodies.FieldGeodeticPoint;
34  import org.orekit.bodies.GeodeticPoint;
35  import org.orekit.errors.OrekitException;
36  import org.orekit.forces.gravity.HolmesFeatherstoneAttractionModel;
37  import org.orekit.forces.gravity.potential.NormalizedSphericalHarmonicsProvider;
38  import org.orekit.forces.gravity.potential.TideSystem;
39  import org.orekit.frames.FieldTransform;
40  import org.orekit.frames.Frame;
41  import org.orekit.frames.Transform;
42  import org.orekit.time.AbsoluteDate;
43  import org.orekit.time.FieldAbsoluteDate;
44  import org.orekit.utils.TimeStampedPVCoordinates;
45  
46  /**
47   * A geoid is a level surface of the gravity potential of a body. The gravity
48   * potential, W, is split so W = U + T, where U is the normal potential (defined
49   * by the ellipsoid) and T is the anomalous potential.[3](eq. 2-137)
50   *
51   * <p> The {@link #getIntersectionPoint(Line, Vector3D, Frame, AbsoluteDate)}
52   * method is tailored specifically for Earth's geoid. All of the other methods
53   * in this class are general and will work for an arbitrary body.
54   *
55   * <p> There are several components that are needed to define a geoid[1]:
56   *
57   * <ul> <li>Geopotential field. These are the coefficients of the spherical
58   * harmonics: S<sub>n,m</sub> and C<sub>n,m</sub></li>
59   *
60   * <li>Reference Ellipsoid. The ellipsoid is used to define the undulation of
61   * the geoid (distance between ellipsoid and geoid) and U<sub>0</sub> the value
62   * of the normal gravity potential at the surface of the ellipsoid.</li>
63   *
64   * <li>W<sub>0</sub>, the potential at the geoid. The value of the potential on
65   * the level surface. This is taken to be U<sub>0</sub>, the normal gravity
66   * potential at the surface of the {@link ReferenceEllipsoid}.</li>
67   *
68   * <li>Permanent Tide System. This implementation assumes that the geopotential
69   * field and the reference ellipsoid use the same permanent tide system. If the
70   * assumption is false it will produce errors of about 0.5 m. Conversion between
71   * tide systems is a possible improvement.[1,2]</li>
72   *
73   * <li>Topographic Masses. That is mass outside of the geoid, e.g. mountains.
74   * This implementation ignores topographic masses, which causes up to 3m error
75   * in the Himalayas, and ~ 1.5m error in the Rockies. This could be improved
76   * through the use of DTED and calculating height anomalies or using the
77   * correction coefficients.[1]</li> </ul>
78   *
79   * <p> This implementation also assumes that the normal to the reference
80   * ellipsoid is the same as the normal to the geoid. This assumption enables the
81   * equation: (height above geoid) = (height above ellipsoid) - (undulation),
82   * which is used in {@link #transform(GeodeticPoint)} and {@link
83   * #transform(Vector3D, Frame, AbsoluteDate)}.
84   *
85   * <p> In testing, the error in the undulations calculated by this class were
86   * off by less than 3 meters, which matches the assumptions outlined above.
87   *
88   * <p> References:
89   *
90   * <ol> <li>Dru A. Smith. There is no such thing as "The" EGM96 geoid: Subtle
91   * points on the use of a global geopotential model. IGeS Bulletin No. 8:17-28,
92   * 1998. <a href= "http://www.ngs.noaa.gov/PUBS_LIB/EGM96_GEOID_PAPER/egm96_geoid_paper.html"
93   * >http ://www.ngs.noaa.gov/PUBS_LIB/EGM96_GEOID_PAPER/egm96_geoid_paper.html</a></li>
94   *
95   * <li> Martin Losch, Verena Seufer. How to Compute Geoid Undulations (Geoid
96   * Height Relative to a Given Reference Ellipsoid) from Spherical Harmonic
97   * Coefficients for Satellite Altimetry Applications. , 2003. <a
98   * href="http://mitgcm.org/~mlosch/geoidcookbook.pdf">mitgcm.org/~mlosch/geoidcookbook.pdf</a>
99   * </li>
100  *
101  * <li>Weikko A. Heiskanen, Helmut Moritz. Physical Geodesy. W. H. Freeman and
102  * Company, 1967. (especially sections 2.13 and equation 2-144 Bruns
103  * Formula)</li>
104  *
105  * <li>S. A. Holmes, W. E. Featherstone. A unified approach to the Clenshaw
106  * summation and the recursive computation of very high degree and order
107  * normalised associated Legendre functions. Journal of Geodesy, 76(5):279,
108  * 2002.</li>
109  *
110  * <li>DMA TR 8350.2. 1984.</li>
111  *
112  * <li>Department of Defense World Geodetic System 1984. 2000. NIMA TR 8350.2
113  * Third Edition, Amendment 1.</li> </ol>
114  *
115  * @author Evan Ward
116  */
117 public class Geoid implements EarthShape {
118 
119     /**
120      * uid is date of last modification.
121      */
122     private static final long serialVersionUID = 20150312L;
123 
124     /**
125      * A number larger than the largest undulation. Wikipedia says the geoid
126      * height is in [-106, 85]. I chose 100 to be safe.
127      */
128     private static final double MAX_UNDULATION = 100;
129     /**
130      * A number smaller than the smallest undulation. Wikipedia says the geoid
131      * height is in [-106, 85]. I chose -150 to be safe.
132      */
133     private static final double MIN_UNDULATION = -150;
134     /**
135      * the maximum number of evaluations for the line search in {@link
136      * #getIntersectionPoint(Line, Vector3D, Frame, AbsoluteDate)}.
137      */
138     private static final int MAX_EVALUATIONS = 100;
139 
140     /**
141      * the default date to use when evaluating the {@link #harmonics}. Used when
142      * no other dates are available. Should be removed in a future release.
143      */
144     private final AbsoluteDate defaultDate;
145     /**
146      * the reference ellipsoid.
147      */
148     private final ReferenceEllipsoid referenceEllipsoid;
149     /**
150      * the geo-potential combined with an algorithm for evaluating the spherical
151      * harmonics. The Holmes and Featherstone method is very robust.
152      */
153     private final transient HolmesFeatherstoneAttractionModel harmonics;
154 
155     /**
156      * Creates a geoid from the given geopotential, reference ellipsoid and the
157      * assumptions in the comment for {@link Geoid}.
158      *
159      * @param geopotential       the gravity potential. Only the anomalous
160      *                           potential will be used. It is assumed that the
161      *                           {@code geopotential} and the {@code
162      *                           referenceEllipsoid} are defined in the same
163      *                           frame. Usually a {@link org.orekit.forces.gravity.potential.GravityFieldFactory#getConstantNormalizedProvider(int,
164      *                           int) constant geopotential} is used to define a
165      *                           time-invariant Geoid.
166      * @param referenceEllipsoid the normal gravity potential.
167      * @throws NullPointerException if {@code geopotential == null ||
168      *                              referenceEllipsoid == null}
169      */
170     public Geoid(final NormalizedSphericalHarmonicsProvider geopotential,
171                  final ReferenceEllipsoid referenceEllipsoid) {
172         // parameter check
173         if (geopotential == null || referenceEllipsoid == null) {
174             throw new NullPointerException();
175         }
176 
177         // subtract the ellipsoid from the geopotential
178         final SubtractEllipsoid potential = new SubtractEllipsoid(geopotential,
179                 referenceEllipsoid);
180 
181         // set instance parameters
182         this.referenceEllipsoid = referenceEllipsoid;
183         this.harmonics = new HolmesFeatherstoneAttractionModel(
184                 referenceEllipsoid.getBodyFrame(), potential);
185         this.defaultDate = geopotential.getReferenceDate();
186     }
187 
188     @Override
189     public Frame getBodyFrame() {
190         // same as for reference ellipsoid.
191         return this.getEllipsoid().getBodyFrame();
192     }
193 
194     /**
195      * Gets the Undulation of the Geoid, N at the given position. N is the
196      * distance between the {@link #getEllipsoid() reference ellipsoid} and the
197      * geoid. The latitude and longitude parameters are both defined with
198      * respect to the reference ellipsoid. For EGM96 and the WGS84 ellipsoid the
199      * undulation is between -107m and +86m.
200      *
201      * <p> NOTE: Restrictions are not put on the range of the arguments {@code
202      * geodeticLatitude} and {@code longitude}.
203      *
204      * @param geodeticLatitude geodetic latitude (angle between the local normal
205      *                         and the equatorial plane on the reference
206      *                         ellipsoid), in radians.
207      * @param longitude        on the reference ellipsoid, in radians.
208      * @param date             of evaluation. Used for time varying geopotential
209      *                         fields.
210      * @return the undulation in m, positive means the geoid is higher than the
211      * ellipsoid.
212      * @see Geoid
213      * @see <a href="http://en.wikipedia.org/wiki/Geoid">Geoid on Wikipedia</a>
214      */
215     public double getUndulation(final double geodeticLatitude,
216                                 final double longitude,
217                                 final AbsoluteDate date) {
218         /*
219          * equations references are to the algorithm printed in the geoid
220          * cookbook[2]. See comment for Geoid.
221          */
222         // reference ellipsoid
223         final ReferenceEllipsoid ellipsoid = this.getEllipsoid();
224 
225         // position in geodetic coordinates
226         final GeodeticPoint.html#GeodeticPoint">GeodeticPoint gp = new GeodeticPoint(geodeticLatitude, longitude, 0);
227         // position in Cartesian coordinates, is converted to geocentric lat and
228         // lon in the Holmes and Featherstone class
229         final Vector3D position = ellipsoid.transform(gp);
230 
231         // get normal gravity from ellipsoid, eq 15
232         final double normalGravity = ellipsoid
233                 .getNormalGravity(geodeticLatitude);
234 
235         // calculate disturbing potential, T, eq 30.
236         final double mu = this.harmonics.getMu();
237         final double T  = this.harmonics.nonCentralPart(date, position, mu);
238         // calculate undulation, eq 30
239         return T / normalGravity;
240     }
241 
242     @Override
243     public ReferenceEllipsoid getEllipsoid() {
244         return this.referenceEllipsoid;
245     }
246 
247     /**
248      * This class implements equations 20-24 in the geoid cook book.(Losch and
249      * Seufer) It modifies C<sub>2n,0</sub> where n = 1,2,...,5.
250      *
251      * @see "DMA TR 8350.2. 1984."
252      */
253     private static final class SubtractEllipsoid implements
254             NormalizedSphericalHarmonicsProvider {
255         /**
256          * provider of the fully normalized coefficients, includes the reference
257          * ellipsoid.
258          */
259         private final NormalizedSphericalHarmonicsProvider provider;
260         /**
261          * the reference ellipsoid to subtract from {@link #provider}.
262          */
263         private final ReferenceEllipsoid ellipsoid;
264 
265         /**
266          * @param provider  potential used for GM<sub>g</sub> and a<sub>g</sub>,
267          *                  and of course the coefficients Cnm, and Snm.
268          * @param ellipsoid Used to calculate the fully normalized
269          *                  J<sub>2n</sub>
270          */
271         private SubtractEllipsoid(
272                 final NormalizedSphericalHarmonicsProvider provider,
273                 final ReferenceEllipsoid ellipsoid) {
274             super();
275             this.provider = provider;
276             this.ellipsoid = ellipsoid;
277         }
278 
279         @Override
280         public int getMaxDegree() {
281             return this.provider.getMaxDegree();
282         }
283 
284         @Override
285         public int getMaxOrder() {
286             return this.provider.getMaxOrder();
287         }
288 
289         @Override
290         public double getMu() {
291             return this.provider.getMu();
292         }
293 
294         @Override
295         public double getAe() {
296             return this.provider.getAe();
297         }
298 
299         @Override
300         public AbsoluteDate getReferenceDate() {
301             return this.provider.getReferenceDate();
302         }
303 
304         @Override
305         public double getOffset(final AbsoluteDate date) {
306             return this.provider.getOffset(date);
307         }
308 
309         @Override
310         public NormalizedSphericalHarmonics onDate(final AbsoluteDate date) {
311             return new NormalizedSphericalHarmonics() {
312 
313                 /** the original harmonics */
314                 private final NormalizedSphericalHarmonics delegate = provider.onDate(date);
315 
316                 @Override
317                 public double getNormalizedCnm(final int n, final int m) {
318                     return getCorrectedCnm(n, m, this.delegate.getNormalizedCnm(n, m));
319                 }
320 
321                 @Override
322                 public double getNormalizedSnm(final int n, final int m) {
323                     return this.delegate.getNormalizedSnm(n, m);
324                 }
325 
326                 @Override
327                 public AbsoluteDate getDate() {
328                     return date;
329                 }
330             };
331         }
332 
333         /**
334          * Get the corrected Cnm for different GM or a values.
335          *
336          * @param n              degree
337          * @param m              order
338          * @param uncorrectedCnm uncorrected Cnm coefficient
339          * @return the corrected Cnm coefficient.
340          */
341         private double getCorrectedCnm(final int n,
342                                        final int m,
343                                        final double uncorrectedCnm) {
344             double Cnm = uncorrectedCnm;
345             // n = 2,4,6,8, or 10 and m = 0
346             if (m == 0 && n <= 10 && n % 2 == 0 && n > 0) {
347                 // correction factor for different GM and a, 1 if no difference
348                 final double gmRatio = this.ellipsoid.getGM() / this.getMu();
349                 final double aRatio = this.ellipsoid.getEquatorialRadius() /
350                         this.getAe();
351                 /*
352                  * eq 20 in the geoid cook book[2], with eq 3-61 in chapter 3 of
353                  * DMA TR 8350.2
354                  */
355                 // halfN = 1,2,3,4,5 for n = 2,4,6,8,10, respectively
356                 final int halfN = n / 2;
357                 Cnm = Cnm - gmRatio * FastMath.pow(aRatio, halfN) *
358                         this.ellipsoid.getC2n0(halfN);
359             }
360             // return is a modified Cnm
361             return Cnm;
362         }
363 
364         @Override
365         public TideSystem getTideSystem() {
366             return this.provider.getTideSystem();
367         }
368 
369     }
370 
371     /**
372      * {@inheritDoc}
373      *
374      * <p> The intersection point is computed using a line search along the
375      * specified line. This is accurate when the geoid is slowly varying.
376      */
377     @Override
378     public GeodeticPoint getIntersectionPoint(final Line lineInFrame,
379                                               final Vector3D closeInFrame,
380                                               final Frame frame,
381                                               final AbsoluteDate date) {
382         /*
383          * It is assumed that the geoid is slowly varying over it's entire
384          * surface. Therefore there will one local intersection.
385          */
386         // transform to body frame
387         final Frame bodyFrame = this.getBodyFrame();
388         final Transform frameToBody = frame.getTransformTo(bodyFrame, date);
389         final Vector3D close = frameToBody.transformPosition(closeInFrame);
390         final Line lineInBodyFrame = frameToBody.transformLine(lineInFrame);
391 
392         // set the line's direction so the solved for value is always positive
393         final Line line;
394         if (lineInBodyFrame.getAbscissa(close) < 0) {
395             line = lineInBodyFrame.revert();
396         } else {
397             line = lineInBodyFrame;
398         }
399 
400         final ReferenceEllipsoid ellipsoid = this.getEllipsoid();
401         // calculate end points
402         // distance from line to center of earth, squared
403         final double d2 = line.pointAt(0.0).getNormSq();
404         // the minimum abscissa, squared
405         final double n = ellipsoid.getPolarRadius() + MIN_UNDULATION;
406         final double minAbscissa2 = n * n - d2;
407         // smaller end point of the interval = 0.0 or intersection with
408         // min_undulation sphere
409         final double lowPoint = FastMath.sqrt(FastMath.max(minAbscissa2, 0.0));
410         // the maximum abscissa, squared
411         final double x = ellipsoid.getEquatorialRadius() + MAX_UNDULATION;
412         final double maxAbscissa2 = x * x - d2;
413         // larger end point of the interval
414         final double highPoint = FastMath.sqrt(maxAbscissa2);
415 
416         // line search function
417         final UnivariateFunction heightFunction = new UnivariateFunction() {
418             @Override
419             public double value(final double x) {
420                 try {
421                     final GeodeticPoint geodetic =
422                             transform(line.pointAt(x), bodyFrame, date);
423                     return geodetic.getAltitude();
424                 } catch (OrekitException e) {
425                     // due to frame transform -> re-throw
426                     throw new RuntimeException(e);
427                 }
428             }
429         };
430 
431         // compute answer
432         if (maxAbscissa2 < 0) {
433             // ray does not pierce bounding sphere -> no possible intersection
434             return null;
435         }
436         // solve line search problem to find the intersection
437         final UnivariateSolver solver = new BracketingNthOrderBrentSolver();
438         try {
439             final double abscissa = solver.solve(MAX_EVALUATIONS, heightFunction, lowPoint, highPoint);
440             // return intersection point
441             return this.transform(line.pointAt(abscissa), bodyFrame, date);
442         } catch (MathRuntimeException e) {
443             // no intersection
444             return null;
445         }
446     }
447 
448     @Override
449     public Vector3D projectToGround(final Vector3D point,
450                                     final AbsoluteDate date,
451                                     final Frame frame) {
452         final GeodeticPoint gp = this.transform(point, frame, date);
453         final GeodeticPoint gpZero =
454                 new GeodeticPoint(gp.getLatitude(), gp.getLongitude(), 0);
455         final Transform bodyToFrame = this.getBodyFrame().getTransformTo(frame, date);
456         return bodyToFrame.transformPosition(this.transform(gpZero));
457     }
458 
459     /**
460      * {@inheritDoc}
461      *
462      * <p> The intersection point is computed using a line search along the
463      * specified line. This is accurate when the geoid is slowly varying.
464      */
465     @Override
466     public <T extends RealFieldElement<T>> FieldGeodeticPoint<T> getIntersectionPoint(final FieldLine<T> lineInFrame,
467                                                                                       final FieldVector3D<T> closeInFrame,
468                                                                                       final Frame frame,
469                                                                                       final FieldAbsoluteDate<T> date) {
470 
471         final Field<T> field = date.getField();
472         /*
473          * It is assumed that the geoid is slowly varying over it's entire
474          * surface. Therefore there will one local intersection.
475          */
476         // transform to body frame
477         final Frame bodyFrame = this.getBodyFrame();
478         final FieldTransform<T> frameToBody = frame.getTransformTo(bodyFrame, date);
479         final FieldVector3D<T> close = frameToBody.transformPosition(closeInFrame);
480         final FieldLine<T> lineInBodyFrame = frameToBody.transformLine(lineInFrame);
481 
482         // set the line's direction so the solved for value is always positive
483         final FieldLine<T> line;
484         if (lineInBodyFrame.getAbscissa(close).getReal() < 0) {
485             line = lineInBodyFrame.revert();
486         } else {
487             line = lineInBodyFrame;
488         }
489 
490         final ReferenceEllipsoid ellipsoid = this.getEllipsoid();
491         // calculate end points
492         // distance from line to center of earth, squared
493         final T d2 = line.pointAt(0.0).getNormSq();
494         // the minimum abscissa, squared
495         final double n = ellipsoid.getPolarRadius() + MIN_UNDULATION;
496         final T minAbscissa2 = d2.negate().add(n * n);
497         // smaller end point of the interval = 0.0 or intersection with
498         // min_undulation sphere
499         final T lowPoint = minAbscissa2.getReal() < 0 ? field.getZero() : minAbscissa2.sqrt();
500         // the maximum abscissa, squared
501         final double x = ellipsoid.getEquatorialRadius() + MAX_UNDULATION;
502         final T maxAbscissa2 = d2.negate().add(x * x);
503         // larger end point of the interval
504         final T highPoint = maxAbscissa2.sqrt();
505 
506         // line search function
507         final RealFieldUnivariateFunction<T> heightFunction = z -> {
508             try {
509                 final FieldGeodeticPoint<T> geodetic =
510                         transform(line.pointAt(z), bodyFrame, date);
511                 return geodetic.getAltitude();
512             } catch (OrekitException e) {
513                 // due to frame transform -> re-throw
514                 throw new RuntimeException(e);
515             }
516         };
517 
518         // compute answer
519         if (maxAbscissa2.getReal() < 0) {
520             // ray does not pierce bounding sphere -> no possible intersection
521             return null;
522         }
523         // solve line search problem to find the intersection
524         final FieldBracketingNthOrderBrentSolver<T> solver =
525                         new FieldBracketingNthOrderBrentSolver<>(field.getZero().add(1.0e-14),
526                                                                  field.getZero().add(1.0e-6),
527                                                                  field.getZero().add(1.0e-15),
528                                                                  5);
529         try {
530             final T abscissa = solver.solve(MAX_EVALUATIONS, heightFunction, lowPoint, highPoint,
531                                             AllowedSolution.ANY_SIDE);
532             // return intersection point
533             return this.transform(line.pointAt(abscissa), bodyFrame, date);
534         } catch (MathRuntimeException e) {
535             // no intersection
536             return null;
537         }
538     }
539 
540     @Override
541     public TimeStampedPVCoordinates projectToGround(
542             final TimeStampedPVCoordinates pv,
543             final Frame frame) {
544         throw new UnsupportedOperationException();
545     }
546 
547     /**
548      * {@inheritDoc}
549      *
550      * @param date date of the conversion. Used for computing frame
551      *             transformations and for time dependent geopotential.
552      * @return The surface relative point at the same location. Altitude is
553      * orthometric height, that is height above the {@link Geoid}. Latitude and
554      * longitude are both geodetic and defined with respect to the {@link
555      * #getEllipsoid() reference ellipsoid}.
556      * @see #transform(GeodeticPoint)
557      * @see <a href="http://en.wikipedia.org/wiki/Orthometric_height">Orthometric_height</a>
558      */
559     @Override
560     public GeodeticPoint transform(final Vector3D point, final Frame frame,
561                                    final AbsoluteDate date) {
562         // convert using reference ellipsoid, altitude referenced to ellipsoid
563         final GeodeticPoint ellipsoidal = this.getEllipsoid().transform(
564                 point, frame, date);
565         // convert altitude to orthometric using the undulation.
566         final double undulation = this.getUndulation(ellipsoidal.getLatitude(),
567                 ellipsoidal.getLongitude(), date);
568         // add undulation to the altitude
569         return new GeodeticPoint(
570                 ellipsoidal.getLatitude(),
571                 ellipsoidal.getLongitude(),
572                 ellipsoidal.getAltitude() - undulation
573         );
574     }
575 
576     /**
577      * {@inheritDoc}
578      *
579      * @param date date of the conversion. Used for computing frame
580      *             transformations and for time dependent geopotential.
581      * @return The surface relative point at the same location. Altitude is
582      * orthometric height, that is height above the {@link Geoid}. Latitude and
583      * longitude are both geodetic and defined with respect to the {@link
584      * #getEllipsoid() reference ellipsoid}.
585      * @see #transform(GeodeticPoint)
586      * @see <a href="http://en.wikipedia.org/wiki/Orthometric_height">Orthometric_height</a>
587      */
588     @Override
589     public <T extends RealFieldElement<T>> FieldGeodeticPoint<T> transform(final FieldVector3D<T> point, final Frame frame,
590                                                                            final FieldAbsoluteDate<T> date) {
591         // convert using reference ellipsoid, altitude referenced to ellipsoid
592         final FieldGeodeticPoint<T> ellipsoidal = this.getEllipsoid().transform(
593                 point, frame, date);
594         // convert altitude to orthometric using the undulation.
595         final double undulation = this.getUndulation(ellipsoidal.getLatitude().getReal(),
596                                                      ellipsoidal.getLongitude().getReal(),
597                                                      date.toAbsoluteDate());
598         // add undulation to the altitude
599         return new FieldGeodeticPoint<>(
600                 ellipsoidal.getLatitude(),
601                 ellipsoidal.getLongitude(),
602                 ellipsoidal.getAltitude().subtract(undulation)
603         );
604     }
605 
606     /**
607      * {@inheritDoc}
608      *
609      * @param point The surface relative point to transform. Altitude is
610      *              orthometric height, that is height above the {@link Geoid}.
611      *              Latitude and longitude are both geodetic and defined with
612      *              respect to the {@link #getEllipsoid() reference ellipsoid}.
613      * @return point at the same location but as a Cartesian point in the {@link
614      * #getBodyFrame() body frame}.
615      * @see #transform(Vector3D, Frame, AbsoluteDate)
616      */
617     @Override
618     public Vector3D transform(final GeodeticPoint point) {
619         try {
620             // convert orthometric height to height above ellipsoid using undulation
621             // TODO pass in date to allow user to specify
622             final double undulation = this.getUndulation(
623                     point.getLatitude(),
624                     point.getLongitude(),
625                     this.defaultDate
626             );
627             final GeodeticPointdeticPoint">GeodeticPoint ellipsoidal = new GeodeticPoint(
628                     point.getLatitude(),
629                     point.getLongitude(),
630                     point.getAltitude() + undulation
631             );
632             // transform using reference ellipsoid
633             return this.getEllipsoid().transform(ellipsoidal);
634         } catch (OrekitException e) {
635             //this method, as defined in BodyShape, is not permitted to throw
636             //an OrekitException, so wrap in an exception we can throw.
637             throw new RuntimeException(e);
638         }
639     }
640 
641     /**
642      * {@inheritDoc}
643      *
644      * @param point The surface relative point to transform. Altitude is
645      *              orthometric height, that is height above the {@link Geoid}.
646      *              Latitude and longitude are both geodetic and defined with
647      *              respect to the {@link #getEllipsoid() reference ellipsoid}.
648      * @param <T> type of the field elements
649      * @return point at the same location but as a Cartesian point in the {@link
650      * #getBodyFrame() body frame}.
651      * @see #transform(Vector3D, Frame, AbsoluteDate)
652      * @since 9.0
653      */
654     @Override
655     public <T extends RealFieldElement<T>> FieldVector3D<T> transform(final FieldGeodeticPoint<T> point) {
656         try {
657             // convert orthometric height to height above ellipsoid using undulation
658             // TODO pass in date to allow user to specify
659             final double undulation = this.getUndulation(
660                     point.getLatitude().getReal(),
661                     point.getLongitude().getReal(),
662                     this.defaultDate
663             );
664             final FieldGeodeticPoint<T> ellipsoidal = new FieldGeodeticPoint<>(
665                     point.getLatitude(),
666                     point.getLongitude(),
667                     point.getAltitude().add(undulation)
668             );
669             // transform using reference ellipsoid
670             return this.getEllipsoid().transform(ellipsoidal);
671         } catch (OrekitException e) {
672             //this method, as defined in BodyShape, is not permitted to throw
673             //an OrekitException, so wrap in an exception we can throw.
674             throw new RuntimeException(e);
675         }
676     }
677 
678 }