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.api;
18  
19  
20  import static org.junit.jupiter.api.Assertions.assertEquals;
21  import static org.junit.jupiter.api.Assertions.assertNull;
22  import static org.junit.jupiter.api.Assertions.assertTrue;
23  
24  import java.io.File;
25  import java.io.IOException;
26  import java.io.RandomAccessFile;
27  import java.lang.reflect.InvocationTargetException;
28  import java.lang.reflect.Method;
29  import java.net.URISyntaxException;
30  import java.nio.MappedByteBuffer;
31  import java.nio.channels.FileChannel;
32  import java.util.ArrayList;
33  import java.util.Collections;
34  import java.util.List;
35  import java.util.Locale;
36  
37  import org.hipparchus.analysis.differentiation.DSFactory;
38  import org.hipparchus.analysis.differentiation.FiniteDifferencesDifferentiator;
39  import org.hipparchus.analysis.differentiation.Gradient;
40  import org.hipparchus.analysis.differentiation.UnivariateDifferentiableFunction;
41  import org.hipparchus.geometry.euclidean.threed.Rotation;
42  import org.hipparchus.geometry.euclidean.threed.RotationConvention;
43  import org.hipparchus.geometry.euclidean.threed.Vector3D;
44  import org.hipparchus.stat.descriptive.StreamingStatistics;
45  import org.hipparchus.stat.descriptive.rank.Percentile;
46  import org.hipparchus.util.FastMath;
47  import org.junit.jupiter.api.Assertions;
48  import org.junit.jupiter.api.BeforeEach;
49  import org.junit.jupiter.api.Disabled;
50  import org.junit.jupiter.api.Test;
51  import org.junit.jupiter.api.io.TempDir;
52  import org.orekit.bodies.BodyShape;
53  import org.orekit.bodies.GeodeticPoint;
54  import org.orekit.data.DataContext;
55  import org.orekit.data.DirectoryCrawler;
56  import org.orekit.errors.OrekitException;
57  import org.orekit.forces.gravity.potential.NormalizedSphericalHarmonicsProvider;
58  import org.orekit.frames.Frame;
59  import org.orekit.frames.FramesFactory;
60  import org.orekit.orbits.Orbit;
61  import org.orekit.propagation.EphemerisGenerator;
62  import org.orekit.propagation.Propagator;
63  import org.orekit.rugged.TestUtils;
64  import org.orekit.rugged.adjustment.GroundOptimizationProblemBuilder;
65  import org.orekit.rugged.adjustment.measurements.Observables;
66  import org.orekit.rugged.adjustment.util.InitInterRefiningTest;
67  import org.orekit.rugged.errors.RuggedException;
68  import org.orekit.rugged.errors.RuggedMessages;
69  import org.orekit.rugged.intersection.IgnoreDEMAlgorithm;
70  import org.orekit.rugged.intersection.IntersectionAlgorithm;
71  import org.orekit.rugged.linesensor.LineDatation;
72  import org.orekit.rugged.linesensor.LineSensor;
73  import org.orekit.rugged.linesensor.LinearLineDatation;
74  import org.orekit.rugged.linesensor.SensorPixel;
75  import org.orekit.rugged.los.FixedRotation;
76  import org.orekit.rugged.los.LOSBuilder;
77  import org.orekit.rugged.los.TimeDependentLOS;
78  import org.orekit.rugged.raster.RandomLandscapeUpdater;
79  import org.orekit.rugged.raster.TileUpdater;
80  import org.orekit.rugged.raster.VolcanicConeElevationUpdater;
81  import org.orekit.rugged.refraction.AtmosphericRefraction;
82  import org.orekit.rugged.utils.DerivativeGenerator;
83  import org.orekit.rugged.utils.NormalizedGeodeticPoint;
84  import org.orekit.time.AbsoluteDate;
85  import org.orekit.time.TimeScale;
86  import org.orekit.time.TimeScalesFactory;
87  import org.orekit.utils.AngularDerivativesFilter;
88  import org.orekit.utils.CartesianDerivativesFilter;
89  import org.orekit.utils.Constants;
90  import org.orekit.utils.IERSConventions;
91  import org.orekit.utils.ParameterDriver;
92  import org.orekit.utils.TimeStampedAngularCoordinates;
93  import org.orekit.utils.TimeStampedPVCoordinates;
94  
95  
96  public class RuggedTest {
97  
98      @TempDir
99      public File tempFolder;
100 
101     // the following test is disabled by default
102     // it is only used to check timings, and also creates a large (66M) temporary file
103     @Disabled
104     @Test
105     public void testMayonVolcanoTiming()
106         throws URISyntaxException {
107 
108         long t0 = System.currentTimeMillis();
109         int dimension = 2000;
110 
111         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
112         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
113         BodyShape  earth                                  = TestUtils.createEarth();
114         NormalizedSphericalHarmonicsProvider gravityField = TestUtils.createGravityField();
115         Orbit      orbit                                  = TestUtils.createOrbit(gravityField.getMu());
116 
117         // Mayon Volcano location according to Wikipedia: 13°15′24″N 123°41′6″E
118         GeodeticPoint summit =
119                 new GeodeticPoint(FastMath.toRadians(13.25667), FastMath.toRadians(123.685), 2463.0);
120         VolcanicConeElevationUpdater updater =
121                 new VolcanicConeElevationUpdater(summit,
122                                                  FastMath.toRadians(30.0), 16.0,
123                                                  FastMath.toRadians(1.0), 1201);
124         final AbsoluteDate crossing = new AbsoluteDate("2012-01-06T02:27:16.139", TimeScalesFactory.getUTC());
125 
126         // one line sensor
127         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
128         // los: swath in the (YZ) plane, centered at +Z, ±10° aperture, 960 pixels
129         Vector3D position = new Vector3D(1.5, 0, -0.2);
130         TimeDependentLOS los = TestUtils.createLOSPerfectLine(Vector3D.PLUS_K, Vector3D.PLUS_I,
131                                                               FastMath.toRadians(10.0), dimension).
132                                                               build();
133 
134         // linear datation model: at reference time we get line 1000, and the rate is one line every 1.5ms
135         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
136         int firstLine = 0;
137         int lastLine  = dimension;
138         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
139 
140         Propagator propagator = TestUtils.createPropagator(earth, gravityField, orbit);
141         propagator.propagate(lineDatation.getDate(firstLine).shiftedBy(-1.0));
142        final EphemerisGenerator generator = propagator.getEphemerisGenerator();
143         propagator.propagate(lineDatation.getDate(lastLine).shiftedBy(+1.0));
144         Propagator ephemeris = generator.getGeneratedEphemeris();
145 
146         Rugged rugged = new RuggedBuilder().
147                 setDigitalElevationModel(updater, 8).
148                 setAlgorithm(AlgorithmId.DUVENHAGE).
149                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
150                 setTimeSpan(lineDatation.getDate(firstLine), lineDatation.getDate(lastLine), 0.001, 5.0).
151                 setTrajectory(1.0 / lineDatation.getRate(0.0),
152                               8, CartesianDerivativesFilter.USE_PV,AngularDerivativesFilter.USE_R,
153                               ephemeris).
154                 addLineSensor(lineSensor).
155                 build();
156 
157         try {
158 
159             int              size   = (lastLine - firstLine) * los.getNbPixels() * 3 * Integer.SIZE / 8;
160             RandomAccessFile out    = new RandomAccessFile(File.createTempFile("junit", null, tempFolder), "rw");
161             MappedByteBuffer buffer = out.getChannel().map(FileChannel.MapMode.READ_WRITE, 0, size);
162 
163             long t1 = System.currentTimeMillis();
164             int pixels = 0;
165             for (double line = firstLine; line < lastLine; line++) {
166                 GeodeticPoint[] gp = rugged.directLocation("line", line);
167                 for (int i = 0; i < gp.length; ++i) {
168                     final int latCode = (int) FastMath.rint(FastMath.scalb(gp[i].getLatitude(),  29));
169                     final int lonCode = (int) FastMath.rint(FastMath.scalb(gp[i].getLongitude(), 29));
170                     final int altCode = (int) FastMath.rint(FastMath.scalb(gp[i].getAltitude(),  17));
171                     buffer.putInt(latCode);
172                     buffer.putInt(lonCode);
173                     buffer.putInt(altCode);
174                 }
175                 pixels += los.getNbPixels();
176                 if  (line % 100 == 0) {
177                     System.out.format(Locale.US, "%5.0f%n", line);
178                 }
179             }
180             long t2 = System.currentTimeMillis();
181             out.close();
182             int sizeM = size / (1024 * 1024);
183             System.out.format(Locale.US,
184                               "%n%n%5dx%5d:%n" +
185                               "  Orekit initialization and DEM creation   : %5.1fs%n" +
186                               "  direct location and %3dM grid writing: %5.1fs (%.1f px/s)%n",
187                               lastLine - firstLine, los.getNbPixels(),
188                               1.0e-3 * (t1 - t0), sizeM, 1.0e-3 * (t2 - t1), pixels / (1.0e-3 * (t2 - t1)));
189         } catch (IOException ioe) {
190             Assertions.fail(ioe.getLocalizedMessage());
191         }
192 
193     }
194 
195     @Test
196     public void testLightTimeCorrection()
197         throws URISyntaxException {
198 
199         int dimension = 400;
200 
201         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
202         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
203         final BodyShape  earth = TestUtils.createEarth();
204         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
205 
206         AbsoluteDate crossing = new AbsoluteDate("2012-01-07T11:21:15.000", TimeScalesFactory.getUTC());
207 
208         // one line sensor
209         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
210         // los: swath in the (YZ) plane, centered at +Z, ±10° aperture, 960 pixels
211         Vector3D position = new Vector3D(1.5, 0, -0.2);
212         TimeDependentLOS los = TestUtils.createLOSPerfectLine(Vector3D.PLUS_K, Vector3D.PLUS_I,
213                                                               FastMath.toRadians(10.0), dimension).build();
214 
215         // linear datation model: at reference time we get line 200, and the rate is one line every 1.5ms
216         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
217         int firstLine = 0;
218         int lastLine  = dimension;
219         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
220         AbsoluteDate minDate = lineSensor.getDate(firstLine);
221         AbsoluteDate maxDate = lineSensor.getDate(lastLine);
222 
223         RuggedBuilder builder = new RuggedBuilder().
224                 setAlgorithm(AlgorithmId.IGNORE_DEM_USE_ELLIPSOID).
225                 setEllipsoid(EllipsoidId.IERS2003, BodyRotatingFrameId.ITRF).
226                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
227                 setTrajectory(InertialFrameId.EME2000,
228                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
229                               8, CartesianDerivativesFilter.USE_PV,
230                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
231                               2, AngularDerivativesFilter.USE_R).
232                 addLineSensor(lineSensor);
233 
234         Rugged rugged = builder.build();
235         Assertions.assertEquals(1, rugged.getLineSensors().size());
236         Assertions.assertTrue(lineSensor == rugged.getLineSensor("line"));
237         try {
238             rugged.getLineSensor("dummy");
239             Assertions.fail("an exception should have been thrown");
240         } catch (RuggedException re) {
241             Assertions.assertEquals(RuggedMessages.UNKNOWN_SENSOR, re.getSpecifier());
242             Assertions.assertEquals("dummy", re.getParts()[0]);
243         }
244         Assertions.assertEquals(7176419.526,
245                             rugged.getScToInertial(lineSensor.getDate(dimension / 2)).getTranslation().getNorm(),
246                             1.0e-3);
247         Assertions.assertEquals(0.0,
248                             rugged.getBodyToInertial(lineSensor.getDate(dimension / 2)).getTranslation().getNorm(),
249                             1.0e-3);
250         Assertions.assertEquals(0.0,
251                             rugged.getInertialToBody(lineSensor.getDate(dimension / 2)).getTranslation().getNorm(),
252                             1.0e-3);
253 
254         builder.setLightTimeCorrection(true);
255         builder.setAberrationOfLightCorrection(false);
256         rugged = builder.build();
257         GeodeticPoint[] gpWithLightTimeCorrection = rugged.directLocation("line", 200);
258 
259         builder.setLightTimeCorrection(false);
260         builder.setAberrationOfLightCorrection(false);
261         rugged = builder.build();
262         GeodeticPoint[] gpWithoutLightTimeCorrection = rugged.directLocation("line", 200);
263 
264         for (int i = 0; i < gpWithLightTimeCorrection.length; ++i) {
265             Vector3D pWith    = earth.transform(gpWithLightTimeCorrection[i]);
266             Vector3D pWithout = earth.transform(gpWithoutLightTimeCorrection[i]);
267             Assertions.assertTrue(Vector3D.distance(pWith, pWithout) > 1.23);
268             Assertions.assertTrue(Vector3D.distance(pWith, pWithout) < 1.27);
269         }
270 
271     }
272 
273     @Test
274     public void testAberrationOfLightCorrection()
275         throws URISyntaxException {
276 
277         int dimension = 400;
278 
279         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
280         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
281         final BodyShape  earth = TestUtils.createEarth();
282         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
283 
284         AbsoluteDate crossing = new AbsoluteDate("2012-01-07T11:46:35.000", TimeScalesFactory.getUTC());
285 
286         // one line sensor
287         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
288         // los: swath in the (YZ) plane, centered at +Z, ±10° aperture, 960 pixels
289         Vector3D position = new Vector3D(1.5, 0, -0.2);
290         TimeDependentLOS los = TestUtils.createLOSPerfectLine(Vector3D.PLUS_K, Vector3D.PLUS_I,
291                                                               FastMath.toRadians(10.0), dimension).build();
292 
293         // linear datation model: at reference time we get line 200, and the rate is one line every 1.5ms
294         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
295         int firstLine = 0;
296         int lastLine  = dimension;
297         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
298         AbsoluteDate minDate = lineSensor.getDate(firstLine);
299         AbsoluteDate maxDate = lineSensor.getDate(lastLine);
300 
301         RuggedBuilder builder = new RuggedBuilder().
302                 
303                 setAlgorithm(AlgorithmId.IGNORE_DEM_USE_ELLIPSOID).
304                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
305                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
306                 setTrajectory(InertialFrameId.EME2000,
307                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
308                               8, CartesianDerivativesFilter.USE_PV,
309                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
310                               2, AngularDerivativesFilter.USE_R).
311                 addLineSensor(lineSensor).
312                 setLightTimeCorrection(false).
313                 setAberrationOfLightCorrection(true);
314         Rugged rugged = builder.build();
315         GeodeticPoint[] gpWithAberrationOfLightCorrection = rugged.directLocation("line", 200);
316 
317         builder.setLightTimeCorrection(false);
318         builder.setAberrationOfLightCorrection(false);
319         rugged = builder.build();
320         GeodeticPoint[] gpWithoutAberrationOfLightCorrection = rugged.directLocation("line", 200);
321 
322         for (int i = 0; i < gpWithAberrationOfLightCorrection.length; ++i) {
323             Vector3D pWith    = earth.transform(gpWithAberrationOfLightCorrection[i]);
324             Vector3D pWithout = earth.transform(gpWithoutAberrationOfLightCorrection[i]);
325             Assertions.assertTrue(Vector3D.distance(pWith, pWithout) > 20.0);
326             Assertions.assertTrue(Vector3D.distance(pWith, pWithout) < 20.5);
327         }
328     }
329     
330     @Test
331     public void testFlatBodyCorrection()
332         throws URISyntaxException {
333 
334         int dimension = 200;
335 
336         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
337         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
338         final BodyShape  earth = TestUtils.createEarth();
339         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
340 
341         AbsoluteDate crossing = new AbsoluteDate("2012-01-01T12:30:00.000", TimeScalesFactory.getUTC());
342 
343         // one line sensor
344         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
345         // los: swath in the (YZ) plane, looking at 50° roll, ±1° aperture
346         Vector3D position = new Vector3D(1.5, 0, -0.2);
347         TimeDependentLOS los = TestUtils.createLOSPerfectLine(new Rotation(Vector3D.PLUS_I,
348                                                                            FastMath.toRadians(50.0),
349                                                                            RotationConvention.VECTOR_OPERATOR).applyTo(Vector3D.PLUS_K),
350                                                               Vector3D.PLUS_I, FastMath.toRadians(1.0), dimension).build();
351 
352         // linear datation model: at reference time we get line 100, and the rate is one line every 1.5ms
353         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
354         int firstLine = 0;
355         int lastLine  = dimension;
356         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
357         AbsoluteDate minDate = lineSensor.getDate(firstLine);
358         AbsoluteDate maxDate = lineSensor.getDate(lastLine);
359 
360         TileUpdater updater =
361                 new RandomLandscapeUpdater(0.0, 9000.0, 0.5, 0xf0a401650191f9f6L,
362                                            FastMath.toRadians(1.0), 257);
363 
364         RuggedBuilder builder = new RuggedBuilder().
365                 setDigitalElevationModel(updater, 8).
366                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
367                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
368                 setTrajectory(InertialFrameId.EME2000,
369                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
370                               8, CartesianDerivativesFilter.USE_PV,
371                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
372                               2, AngularDerivativesFilter.USE_R).
373                 addLineSensor(lineSensor);
374         GeodeticPoint[] gpWithFlatBodyCorrection =
375                 builder.setAlgorithm(AlgorithmId.DUVENHAGE).build().directLocation("line", 100);
376 
377         GeodeticPoint[] gpWithoutFlatBodyCorrection =
378                 builder.setAlgorithm(AlgorithmId.DUVENHAGE_FLAT_BODY).build().directLocation("line", 100);
379 
380         StreamingStatistics stats = new StreamingStatistics();
381         for (int i = 0; i < gpWithFlatBodyCorrection.length; ++i) {
382             Vector3D pWith    = earth.transform(gpWithFlatBodyCorrection[i]);
383             Vector3D pWithout = earth.transform(gpWithoutFlatBodyCorrection[i]);
384             stats.addValue(Vector3D.distance(pWith, pWithout));
385         }
386         Assertions.assertEquals( 0.004, stats.getMin(),  1.0e-3);
387         Assertions.assertEquals(28.344, stats.getMax(),  1.0e-3);
388         Assertions.assertEquals( 4.801, stats.getMean(), 1.0e-3);
389 
390     }
391 
392     @Test
393     public void testLocationSinglePoint()
394         throws URISyntaxException {
395 
396         int dimension = 200;
397 
398         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
399         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
400         final BodyShape  earth = TestUtils.createEarth();
401         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
402 
403         AbsoluteDate crossing = new AbsoluteDate("2012-01-01T12:30:00.000", TimeScalesFactory.getUTC());
404 
405         // one line sensor
406         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
407         // los: swath in the (YZ) plane, looking at 50° roll, ±1° aperture
408         Vector3D position = new Vector3D(1.5, 0, -0.2);
409         TimeDependentLOS los = TestUtils.createLOSPerfectLine(new Rotation(Vector3D.PLUS_I,
410                                                                            FastMath.toRadians(50.0),
411                                                                            RotationConvention.VECTOR_OPERATOR).applyTo(Vector3D.PLUS_K),
412                                                               Vector3D.PLUS_I, FastMath.toRadians(1.0), dimension).build();
413 
414         // linear datation model: at reference time we get line 100, and the rate is one line every 1.5ms
415         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
416         int firstLine = 0;
417         int lastLine  = dimension;
418         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
419         AbsoluteDate minDate = lineSensor.getDate(firstLine);
420         AbsoluteDate maxDate = lineSensor.getDate(lastLine);
421 
422         TileUpdater updater =
423                 new RandomLandscapeUpdater(0.0, 9000.0, 0.5, 0xf0a401650191f9f6L,
424                                            FastMath.toRadians(1.0), 257);
425 
426         Rugged rugged = new RuggedBuilder().
427                 setDigitalElevationModel(updater, 8).
428                 setAlgorithm(AlgorithmId.DUVENHAGE).
429                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
430                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
431                 setTrajectory(InertialFrameId.EME2000,
432                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
433                               8, CartesianDerivativesFilter.USE_PV,
434                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
435                               2, AngularDerivativesFilter.USE_R).
436                 addLineSensor(lineSensor).build();
437         GeodeticPoint[] gpLine = rugged.directLocation("line", 100);
438 
439         for (int i = 0; i < gpLine.length; ++i) {
440             GeodeticPoint gpPixel =
441                     rugged.directLocation(lineSensor.getDate(100), lineSensor.getPosition(),
442                                               lineSensor.getLOS(lineSensor.getDate(100), i));
443             Assertions.assertEquals(gpLine[i].getLatitude(),  gpPixel.getLatitude(),  1.0e-10);
444             Assertions.assertEquals(gpLine[i].getLongitude(), gpPixel.getLongitude(), 1.0e-10);
445             Assertions.assertEquals(gpLine[i].getAltitude(),  gpPixel.getAltitude(),  1.0e-10);
446         }
447 
448     }
449 
450     @Test
451     public void testLocationsinglePointNoCorrections()
452         throws URISyntaxException {
453 
454         int dimension = 200;
455 
456         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
457         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
458         final BodyShape  earth = TestUtils.createEarth();
459         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
460 
461         AbsoluteDate crossing = new AbsoluteDate("2012-01-01T12:30:00.000", TimeScalesFactory.getUTC());
462 
463         // one line sensor
464         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
465         // los: swath in the (YZ) plane, looking at 50° roll, ±1° aperture
466         Vector3D position = new Vector3D(1.5, 0, -0.2);
467         TimeDependentLOS los = TestUtils.createLOSPerfectLine(new Rotation(Vector3D.PLUS_I,
468                                                                            FastMath.toRadians(50.0),
469                                                                            RotationConvention.VECTOR_OPERATOR).applyTo(Vector3D.PLUS_K),
470                                                               Vector3D.PLUS_I, FastMath.toRadians(1.0), dimension).build();
471 
472         // linear datation model: at reference time we get line 100, and the rate is one line every 1.5ms
473         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
474         int firstLine = 0;
475         int lastLine  = dimension;
476         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
477         AbsoluteDate minDate = lineSensor.getDate(firstLine);
478         AbsoluteDate maxDate = lineSensor.getDate(lastLine);
479 
480         TileUpdater updater =
481                 new RandomLandscapeUpdater(0.0, 9000.0, 0.5, 0xf0a401650191f9f6L,
482                                            FastMath.toRadians(1.0), 257);
483 
484         Rugged rugged = new RuggedBuilder().
485                 setDigitalElevationModel(updater, 8).
486                 setAlgorithm(AlgorithmId.DUVENHAGE).
487                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
488                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
489                 setTrajectory(InertialFrameId.EME2000,
490                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
491                               8, CartesianDerivativesFilter.USE_PV,
492                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
493                               2, AngularDerivativesFilter.USE_R).
494                 setAberrationOfLightCorrection(false).
495                 setLightTimeCorrection(false).
496                 addLineSensor(lineSensor).
497                 build();
498         GeodeticPoint[] gpLine = rugged.directLocation("line", 100);
499 
500         for (int i = 0; i < gpLine.length; ++i) {
501             GeodeticPoint gpPixel =
502                     rugged.directLocation(lineSensor.getDate(100), lineSensor.getPosition(),
503                                               lineSensor.getLOS(lineSensor.getDate(100), i));
504             Assertions.assertEquals(gpLine[i].getLatitude(),  gpPixel.getLatitude(),  1.0e-10);
505             Assertions.assertEquals(gpLine[i].getLongitude(), gpPixel.getLongitude(), 1.0e-10);
506             Assertions.assertEquals(gpLine[i].getAltitude(),  gpPixel.getAltitude(),  1.0e-10);
507         }
508 
509     }
510 
511     @Test
512     public void testBasicScan()
513         throws URISyntaxException {
514 
515         int dimension = 200;
516 
517         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
518         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
519         final BodyShape  earth = TestUtils.createEarth();
520         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
521 
522         AbsoluteDate crossing = new AbsoluteDate("2012-01-01T12:30:00.000", TimeScalesFactory.getUTC());
523 
524         // one line sensor
525         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
526         // los: swath in the (YZ) plane, looking at 50° roll, ±1° aperture
527         Vector3D position = new Vector3D(1.5, 0, -0.2);
528         TimeDependentLOS los = TestUtils.createLOSPerfectLine(new Rotation(Vector3D.PLUS_I,
529                                                                            FastMath.toRadians(50.0),
530                                                                            RotationConvention.VECTOR_OPERATOR).applyTo(Vector3D.PLUS_K),
531                                                               Vector3D.PLUS_I, FastMath.toRadians(1.0), dimension).build();
532 
533         // linear datation model: at reference time we get line 100, and the rate is one line every 1.5ms
534         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
535         int firstLine = 0;
536         int lastLine  = dimension;
537         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
538         AbsoluteDate minDate = lineSensor.getDate(firstLine);
539         AbsoluteDate maxDate = lineSensor.getDate(lastLine);
540 
541         TileUpdater updater =
542                 new RandomLandscapeUpdater(0.0, 9000.0, 0.5, 0xe12ef744f224cf43L,
543                                            FastMath.toRadians(1.0), 257);
544 
545         RuggedBuilder builder = new RuggedBuilder().
546                 setDigitalElevationModel(updater, 8).
547                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
548                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
549                 setTrajectory(InertialFrameId.EME2000,
550                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
551                               8, CartesianDerivativesFilter.USE_PV,
552                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
553                               2, AngularDerivativesFilter.USE_R).
554                 addLineSensor(lineSensor);
555         GeodeticPoint[] gpDuvenhage =
556                 builder.setAlgorithm(AlgorithmId.DUVENHAGE).build().directLocation("line", 100);
557 
558         GeodeticPoint[] gpBasicScan =
559                 builder.setAlgorithm(AlgorithmId.BASIC_SLOW_EXHAUSTIVE_SCAN_FOR_TESTS_ONLY).build().directLocation("line", 100);
560 
561         double[] data = new double[gpDuvenhage.length]; 
562         for (int i = 0; i < gpDuvenhage.length; ++i) {
563             Vector3D pDuvenhage = earth.transform(gpDuvenhage[i]);
564             Vector3D pBasicScan = earth.transform(gpBasicScan[i]);
565             data[i] = Vector3D.distance(pDuvenhage, pBasicScan);
566         }
567         Assertions.assertEquals(0.0, new Percentile(99).evaluate(data), 5.1e-4);
568 
569     }
570 
571     // the following test is disabled by default
572     // it is only used to check timings, and also creates a large (38M) temporary file
573     @Disabled
574     @Test
575     public void testInverseLocationTiming()
576         throws URISyntaxException {
577 
578         long t0       = System.currentTimeMillis();
579         int dimension = 2000;
580         int nbSensors = 3;
581 
582         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
583         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
584         final BodyShape  earth = TestUtils.createEarth();
585         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
586 
587         AbsoluteDate crossing = new AbsoluteDate("2012-01-01T12:30:00.000", TimeScalesFactory.getUTC());
588 
589         List<LineSensor> sensors = new ArrayList<>();
590         for (int i = 0; i < nbSensors; ++i) {
591             // one line sensor
592             // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
593             // los: swath in the (YZ) plane, looking roughly at 50° roll (sensor-dependent), 5.2" per pixel
594             Vector3D position = new Vector3D(1.5, 0, -0.2);
595             TimeDependentLOS los = TestUtils.createLOSPerfectLine(new Rotation(Vector3D.PLUS_I,
596                                                                                FastMath.toRadians(50.0 - 0.001 * i),
597                                                                                RotationConvention.VECTOR_OPERATOR).applyTo(Vector3D.PLUS_K),
598                                                                   Vector3D.PLUS_I, FastMath.toRadians((dimension/2.) * 5.2 / 3600.0), dimension).build();
599 
600             // linear datation model: at reference time we get roughly middle line, and the rate is one line every 1.5ms
601             LineDatation lineDatation = new LinearLineDatation(crossing, i + dimension / 2, 1.0 / 1.5e-3);
602             sensors.add(new LineSensor("line-" + i, lineDatation, position, los));
603         }
604 
605         int firstLine = 0;
606         int lastLine  = dimension;
607         AbsoluteDate minDate = sensors.get(0).getDate(firstLine).shiftedBy(-1.0);
608         AbsoluteDate maxDate = sensors.get(sensors.size() - 1).getDate(lastLine).shiftedBy(+1.0);
609 
610         TileUpdater updater =
611                 new RandomLandscapeUpdater(0.0, 9000.0, 0.5, 0xf0a401650191f9f6L,
612                                            FastMath.toRadians(1.0), 257);
613 
614         RuggedBuilder builder = new RuggedBuilder().
615                 setDigitalElevationModel(updater, 8).
616                 setAlgorithm(AlgorithmId.DUVENHAGE).
617                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
618                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
619                 setTrajectory(InertialFrameId.EME2000,
620                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
621                               8, CartesianDerivativesFilter.USE_PV,
622                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
623                               2, AngularDerivativesFilter.USE_R);
624         builder.setLightTimeCorrection(true);
625         builder.setAberrationOfLightCorrection(true);
626         for (LineSensor lineSensor : sensors) {
627             builder.addLineSensor(lineSensor);
628         }
629         Rugged rugged = builder.build();
630 
631         double lat0  = FastMath.toRadians(-22.9);
632         double lon0  = FastMath.toRadians(142.4);
633         double delta = FastMath.toRadians(0.5);
634 
635         try {
636             long             size   = nbSensors * dimension * dimension * 2L * Integer.SIZE / 8L;
637             RandomAccessFile out    = new RandomAccessFile(File.createTempFile("junit", null, tempFolder), "rw");
638             MappedByteBuffer buffer = out.getChannel().map(FileChannel.MapMode.READ_WRITE, 0, size);
639 
640             long t1 = System.currentTimeMillis();
641             int goodPixels = 0;
642             int badPixels  = 0;
643             for (final LineSensor lineSensor : sensors) {
644                 for (int i = 0; i < dimension; ++i) {
645                     double latitude  = lat0 + (i * delta) / dimension;
646                     if (i %100 == 0) {
647                         System.out.println(lineSensor.getName() + " " + i);
648                     }
649                     for (int j = 0; j < dimension; ++j) {
650                         double longitude = lon0 + (j * delta) / dimension;
651                         SensorPixel sp = rugged.inverseLocation(lineSensor.getName(), latitude, longitude,
652                                                                 firstLine, lastLine);
653                         if (sp == null) {
654                             ++badPixels;
655                             buffer.putInt(-1);
656                             buffer.putInt(-1);
657                         } else {
658                             ++goodPixels;
659                             final int lineCode  = (int) FastMath.rint(FastMath.scalb(sp.getLineNumber(),  16));
660                             final int pixelCode = (int) FastMath.rint(FastMath.scalb(sp.getPixelNumber(), 16));
661                             buffer.putInt(lineCode);
662                             buffer.putInt(pixelCode);
663                         }
664                     }
665                 }
666             }
667 
668             long t2 = System.currentTimeMillis();
669             out.close();
670             int sizeM = (int) (size / (1024L * 1024L));
671             System.out.format(Locale.US,
672                               "%n%n%5dx%5d, %d sensors:%n" +
673                               "  Orekit initialization and DEM creation   : %5.1fs%n" +
674                               "  inverse location and %3dM grid writing: %5.1fs (%.1f px/s, %.1f%% covered)%n",
675                               dimension, dimension, nbSensors,
676                               1.0e-3 * (t1 - t0), sizeM, 1.0e-3 * (t2 - t1),
677                               (badPixels + goodPixels) / (1.0e-3 * (t2 - t1)),
678                               (100.0 * goodPixels) / (goodPixels + badPixels));
679         } catch (IOException ioe) {
680             Assertions.fail(ioe.getLocalizedMessage());
681         }
682     }
683 
684     @Test
685     public void testInverseLocation()
686         throws URISyntaxException {
687         checkInverseLocation(2000, false, false, 4.0e-7, 5.0e-6);
688         checkInverseLocation(2000, false, true,  1.0e-5, 2.0e-7);
689         checkInverseLocation(2000, true,  false, 4.0e-7, 4.0e-7);
690         checkInverseLocation(2000, true,  true,  2.0e-5, 3.0e-7);
691     }
692 
693     @Test
694     public void testDateLocation()
695         throws URISyntaxException {
696         checkDateLocation(2000, false, false, 7.0e-7);
697         checkDateLocation(2000, false, true,  2.0e-5);
698         checkDateLocation(2000, true,  false, 8.0e-7);
699         checkDateLocation(2000, true,  true,  3.0e-6);
700     }
701     
702     @Test
703     public void testLineDatation()
704         throws URISyntaxException {
705         checkLineDatation(2000, 7.0e-7);
706         checkLineDatation(10000, 8.0e-7);
707     }
708 
709 
710     @Test
711     public void testInverseLocNearLineEnd() throws URISyntaxException {
712 
713         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
714         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
715         Vector3D offset = Vector3D.ZERO;
716         TimeScale gps = TimeScalesFactory.getGPS();
717         Frame eme2000 = FramesFactory.getEME2000();
718         Frame itrf = FramesFactory.getITRF(IERSConventions.IERS_2010, true);
719         ArrayList<TimeStampedAngularCoordinates> satelliteQList = new ArrayList<>();
720         ArrayList<TimeStampedPVCoordinates> satellitePVList = new ArrayList<>();
721 
722         TestUtils.addSatelliteQ(gps, satelliteQList, "2009-12-11T16:58:42.592937", -0.340236d, 0.333952d, -0.844012d, -0.245684d);
723         TestUtils.addSatelliteQ(gps, satelliteQList, "2009-12-11T16:59:06.592937", -0.354773d, 0.329336d, -0.837871d, -0.252281d);
724         TestUtils.addSatelliteQ(gps, satelliteQList, "2009-12-11T16:59:30.592937", -0.369237d, 0.324612d, -0.831445d, -0.258824d);
725         TestUtils.addSatelliteQ(gps, satelliteQList, "2009-12-11T16:59:54.592937", -0.3836d, 0.319792d, -0.824743d, -0.265299d);
726         TestUtils.addSatelliteQ(gps, satelliteQList, "2009-12-11T17:00:18.592937", -0.397834d, 0.314883d, -0.817777d, -0.271695d);
727         TestUtils.addSatelliteQ(gps, satelliteQList, "2009-12-11T17:00:42.592937", -0.411912d, 0.309895d, -0.810561d, -0.278001d);
728         TestUtils.addSatelliteQ(gps, satelliteQList, "2009-12-11T17:01:06.592937", -0.42581d, 0.304838d, -0.803111d, -0.284206d);
729         TestUtils.addSatelliteQ(gps, satelliteQList, "2009-12-11T17:01:30.592937", -0.439505d, 0.299722d, -0.795442d, -0.290301d);
730         TestUtils.addSatelliteQ(gps, satelliteQList, "2009-12-11T17:01:54.592937", -0.452976d, 0.294556d, -0.787571d, -0.296279d);
731         TestUtils.addSatelliteQ(gps, satelliteQList, "2009-12-11T17:02:18.592937", -0.466207d, 0.28935d, -0.779516d, -0.302131d);
732 
733 
734         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T16:58:42.592937", -726361.466d, -5411878.485d, 4637549.599d, -2068.995d, -4500.601d, -5576.736d);
735         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T16:59:04.192937", -779538.267d, -5506500.533d, 4515934.894d, -2058.308d, -4369.521d, -5683.906d);
736         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T16:59:25.792937", -832615.368d, -5598184.195d, 4392036.13d, -2046.169d, -4236.279d, -5788.201d);
737         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T16:59:47.392937", -885556.748d, -5686883.696d, 4265915.971d, -2032.579d, -4100.944d, -5889.568d);
738         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T17:00:08.992937", -938326.32d, -5772554.875d, 4137638.207d, -2017.537d, -3963.59d, -5987.957d);
739         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T17:00:30.592937", -990887.942d, -5855155.21d, 4007267.717d, -2001.046d, -3824.291d, -6083.317d);
740         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T17:00:52.192937", -1043205.448d, -5934643.836d, 3874870.441d, -1983.107d, -3683.122d, -6175.6d);
741         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T17:01:13.792937", -1095242.669d, -6010981.571d, 3740513.34d, -1963.723d, -3540.157d, -6264.76d);
742         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T17:01:35.392937", -1146963.457d, -6084130.93d, 3604264.372d, -1942.899d, -3395.473d, -6350.751d);
743         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T17:01:56.992937", -1198331.706d, -6154056.146d, 3466192.446d, -1920.64d, -3249.148d, -6433.531d);
744         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2009-12-11T17:02:18.592937", -1249311.381d, -6220723.191d, 3326367.397d, -1896.952d, -3101.26d, -6513.056d);
745 
746         TileUpdater updater = new RandomLandscapeUpdater(0.0, 9000.0, 0.5, 0xf0a401650191f9f6L, FastMath.toRadians(1.0), 257);
747         RuggedBuilder builder = new RuggedBuilder().
748                 setDigitalElevationModel(updater, 8).
749                 setAlgorithm(AlgorithmId.DUVENHAGE).
750                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
751                 setTimeSpan(satellitePVList.get(0).getDate(), satellitePVList.get(satellitePVList.size() - 1).getDate(), 0.001, 5.0).
752                 setTrajectory(InertialFrameId.EME2000,
753                               satellitePVList, 8, CartesianDerivativesFilter.USE_PV,
754                               satelliteQList, 8, AngularDerivativesFilter.USE_R);
755 
756         List<Vector3D> lineOfSight = new ArrayList<>();
757         lineOfSight.add(new Vector3D(-0.011204, 0.181530, 1d).normalize());
758         lineOfSight.add(new Vector3D(-0.011204, 0.181518, 1d).normalize());
759         lineOfSight.add(new Vector3D(-0.011204, 0.181505, 1d).normalize());
760         lineOfSight.add(new Vector3D(-0.011204, 0.181492, 1d).normalize());
761         lineOfSight.add(new Vector3D(-0.011204, 0.181480, 1d).normalize());
762         lineOfSight.add(new Vector3D(-0.011204, 0.181467, 1d).normalize());
763         lineOfSight.add(new Vector3D(-0.011204, 0.181455, 1d).normalize());
764         lineOfSight.add(new Vector3D(-0.011204, 0.181442, 1d).normalize());
765         lineOfSight.add(new Vector3D(-0.011204, 0.181430, 1d).normalize());
766         lineOfSight.add(new Vector3D(-0.011204, 0.181417, 1d).normalize());
767         lineOfSight.add(new Vector3D(-0.011204, 0.181405, 1d).normalize());
768         lineOfSight.add(new Vector3D(-0.011204, 0.181392, 1d).normalize());
769 
770         lineOfSight.add(new Vector3D(-0.011204, 0.149762, 1d).normalize());
771         lineOfSight.add(new Vector3D(-0.011204, 0.149749, 1d).normalize());
772         lineOfSight.add(new Vector3D(-0.011204, 0.149737, 1d).normalize());
773         lineOfSight.add(new Vector3D(-0.011204, 0.149724, 1d).normalize());
774         lineOfSight.add(new Vector3D(-0.011204, 0.149712, 1d).normalize());
775         lineOfSight.add(new Vector3D(-0.011204, 0.149699, 1d).normalize());
776         lineOfSight.add(new Vector3D(-0.011204, 0.149686, 1d).normalize());
777         lineOfSight.add(new Vector3D(-0.011204, 0.149674, 1d).normalize());
778         lineOfSight.add(new Vector3D(-0.011204, 0.149661, 1d).normalize());
779         lineOfSight.add(new Vector3D(-0.011204, 0.149649, 1d).normalize());
780         lineOfSight.add(new Vector3D(-0.011204, 0.149636, 1d).normalize());
781         lineOfSight.add(new Vector3D(-0.011204, 0.149624, 1d).normalize());
782         lineOfSight.add(new Vector3D(-0.011204, 0.149611, 1d).normalize());
783         lineOfSight.add(new Vector3D(-0.011204, 0.149599, 1d).normalize());
784         lineOfSight.add(new Vector3D(-0.011204, 0.149586, 1d).normalize());
785 
786         AbsoluteDate absDate = new AbsoluteDate("2009-12-11T16:58:51.593", gps);
787         LinearLineDatation lineDatation = new LinearLineDatation(absDate, 1d, 638.5696040868454);
788         LineSensor lineSensor = new LineSensor("perfect-line", lineDatation, offset,
789                                                new LOSBuilder(lineOfSight).build());
790         builder.addLineSensor(lineSensor);
791 
792         Rugged rugged = builder.build();
793         GeodeticPoint point1 = new GeodeticPoint(0.7053784581520293, -1.7354535645320581, 691.856741468848);
794         SensorPixel sensorPixel1 = rugged.inverseLocation(lineSensor.getName(), point1, 1, 131328);
795         Assertions.assertEquals(   2.01474, sensorPixel1.getLineNumber(), 1.0e-5);
796         Assertions.assertEquals(   2.09271, sensorPixel1.getPixelNumber(), 1.0e-5);
797         GeodeticPoint point2 = new GeodeticPoint(0.704463899881073, -1.7303503789334154, 648.9200602492216);
798         SensorPixel sensorPixel2 = rugged.inverseLocation(lineSensor.getName(), point2, 1, 131328);
799         Assertions.assertEquals(   2.02185, sensorPixel2.getLineNumber(), 1.0e-5);
800         Assertions.assertEquals(  27.53008, sensorPixel2.getPixelNumber(), 1.0e-5);
801         GeodeticPoint point3 = new GeodeticPoint(0.7009593480939814, -1.7314283804521957, 588.3075485689468);
802         SensorPixel sensorPixel3 = rugged.inverseLocation(lineSensor.getName(), point3, 1, 131328);
803         Assertions.assertEquals(2305.26101, sensorPixel3.getLineNumber(),  1.0e-5);
804         Assertions.assertEquals(  27.18381, sensorPixel3.getPixelNumber(), 1.0e-5);
805         GeodeticPoint point4 = new GeodeticPoint(0.7018731669637096, -1.73651769725183, 611.2759403696498);
806         SensorPixel sensorPixel4 = rugged.inverseLocation(lineSensor.getName(), point4, 1, 131328);
807         Assertions.assertEquals(2305.25436, sensorPixel4.getLineNumber(), 1.0e-5);
808         Assertions.assertEquals(   1.54447, sensorPixel4.getPixelNumber(), 1.0e-5);
809 
810     }
811 
812     @Test
813     public void testInverseLoc() throws URISyntaxException {
814 
815         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
816         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
817         Vector3D offset = Vector3D.ZERO;
818         TimeScale gps = TimeScalesFactory.getGPS();
819         Frame eme2000 = FramesFactory.getEME2000();
820         Frame itrf = FramesFactory.getITRF(IERSConventions.IERS_2010, true);
821         ArrayList<TimeStampedAngularCoordinates> satelliteQList = new ArrayList<>();
822         ArrayList<TimeStampedPVCoordinates> satellitePVList = new ArrayList<>();
823 
824         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:16:27", -0.327993d, -0.715194d, -0.56313d, 0.252592d);
825         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:16:29", -0.328628d, -0.71494d, -0.562769d, 0.25329d);
826         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:16:31", -0.329263d, -0.714685d, -0.562407d, 0.253988d);
827         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:16:33", -0.329898d, -0.714429d, -0.562044d, 0.254685d);
828         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:16:35", -0.330532d, -0.714173d, -0.561681d, 0.255383d);
829         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:16:37", -0.331166d, -0.713915d, -0.561318d, 0.256079d);
830         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:16:39", -0.3318d, -0.713657d, -0.560954d, 0.256776d);
831         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:16:41", -0.332434d, -0.713397d, -0.560589d, 0.257472d);
832         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:16:43", -0.333067d, -0.713137d, -0.560224d, 0.258168d);
833         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:16:45", -0.333699d, -0.712876d, -0.559859d, 0.258864d);
834         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:18:17", -0.36244d, -0.699935d, -0.542511d, 0.290533d);
835         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:20:27", -0.401688d, -0.678574d, -0.516285d, 0.334116d);
836         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:20:29", -0.402278d, -0.678218d, -0.515866d, 0.334776d);
837         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:20:31", -0.402868d, -0.677861d, -0.515447d, 0.335435d);
838         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:20:33", -0.403457d, -0.677503d, -0.515028d, 0.336093d);
839         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:20:35", -0.404046d, -0.677144d, -0.514608d, 0.336752d);
840         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:20:37", -0.404634d, -0.676785d, -0.514187d, 0.337409d);
841         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:20:39", -0.405222d, -0.676424d, -0.513767d, 0.338067d);
842         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:20:41", -0.40581d, -0.676063d, -0.513345d, 0.338724d);
843         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:20:43", -0.406397d, -0.675701d, -0.512924d, 0.339381d);
844         TestUtils.addSatelliteQ(gps, satelliteQList, "2013-07-07T17:20:45", -0.406983d, -0.675338d, -0.512502d, 0.340038d);
845 
846         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:16:27.857531", -379110.393d, -5386317.278d, 4708158.61d, -1802.078d, -4690.847d,  -5512.223d);
847         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:16:36.857531", -398874.476d, -5428039.968d, 4658344.906d, -1801.326d, -4636.91d,  -5557.915d);
848         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:16:45.857531", -418657.992d, -5469262.453d, 4608122.145d, -1800.345d, -4582.57d,  -5603.119d);
849         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:16:54.857531", -438458.554d, -5509981.109d, 4557494.737d, -1799.136d, -4527.831d, -5647.831d);
850         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:17:03.857531", -458273.771d, -5550192.355d, 4506467.128d, -1797.697d, -4472.698d, -5692.046d);
851         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:17:12.857531", -478101.244d, -5589892.661d, 4455043.798d, -1796.029d, -4417.176d, -5735.762d);
852         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:17:21.857531", -497938.57d, -5629078.543d, 4403229.263d, -1794.131d, -4361.271d,  -5778.975d);
853         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:17:30.857531", -517783.34d, -5667746.565d, 4351028.073d, -1792.003d, -4304.987d,  -5821.679d);
854         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:17:39.857531", -537633.139d, -5705893.34d, 4298444.812d, -1789.644d, -4248.329d,  -5863.873d);
855         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:17:48.857531", -557485.549d, -5743515.53d, 4245484.097d, -1787.055d, -4191.304d,  -5905.552d);
856         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:17:57.857531", -577338.146d, -5780609.846d, 4192150.579d, -1784.234d, -4133.916d, -5946.712d);
857         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:18:06.857531", -597188.502d, -5817173.047d, 4138448.941d, -1781.183d, -4076.171d, -5987.35d);
858         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:18:15.857531", -617034.185d, -5853201.943d, 4084383.899d, -1777.899d, -4018.073d, -6027.462d);
859         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:18:24.857531", -636872.759d, -5888693.393d, 4029960.2d, -1774.385d, -3959.629d,   -6067.045d);
860         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:18:33.857531", -656701.786d, -5923644.307d, 3975182.623d, -1770.638d, -3900.844d, -6106.095d);
861         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:18:42.857531", -676518.822d, -5958051.645d, 3920055.979d, -1766.659d, -3841.723d, -6144.609d);
862         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:18:51.857531", -696321.424d, -5991912.417d, 3864585.108d, -1762.449d, -3782.271d, -6182.583d);
863         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:19:00.857531", -716107.143d, -6025223.686d, 3808774.881d, -1758.006d, -3722.495d, -6220.015d);
864         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:19:09.857531", -735873.528d, -6057982.563d, 3752630.2d, -1753.332d, -3662.399d,   -6256.9d);
865         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:19:18.857531", -755618.129d, -6090186.214d, 3696155.993d, -1748.425d, -3601.99d,  -6293.236d);
866         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:19:27.857531", -775338.49d, -6121831.854d, 3639357.221d, -1743.286d, -3541.272d,  -6329.019d);
867         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:19:36.857531", -795032.157d, -6152916.751d, 3582238.87d, -1737.915d, -3480.252d,  -6364.246d);
868         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:19:45.857531", -814696.672d, -6183438.226d, 3524805.957d, -1732.313d, -3418.935d, -6398.915d);
869         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:19:54.857531", -834329.579d, -6213393.652d, 3467063.525d, -1726.478d, -3357.327d, -6433.022d);
870         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:20:03.857531", -853928.418d, -6242780.453d, 3409016.644d, -1720.412d, -3295.433d, -6466.563d);
871         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:20:12.857531", -873490.732d, -6271596.108d, 3350670.411d, -1714.114d, -3233.259d, -6499.537d);
872         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:20:21.857531", -893014.061d, -6299838.148d, 3292029.951d, -1707.585d, -3170.811d, -6531.941d);
873         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:20:30.857531", -912495.948d, -6327504.159d, 3233100.411d, -1700.825d, -3108.095d, -6563.77d);
874         TestUtils.addSatellitePV(gps, eme2000, itrf, satellitePVList, "2013-07-07T17:20:39.857531", -931933.933d, -6354591.778d, 3173886.968d, -1693.833d, -3045.116d, -6595.024d);
875 
876         List<Vector3D> lineOfSight = new ArrayList<>();
877         lineOfSight.add(new Vector3D(0.0046536264d, -0.1851800945d, 1d));
878         lineOfSight.add(new Vector3D(0.0000001251d, -0.0002815246d, 1d));
879         lineOfSight.add(new Vector3D(0.0046694108d, 0.1853863933d, 1d));
880 
881         AbsoluteDate absDate = new AbsoluteDate("2013-07-07T17:16:36.857", gps);
882         LinearLineDatation lineDatation = new LinearLineDatation(absDate, 0.03125d, 19.95565693384045);
883         LineSensor lineSensor = new LineSensor("QUICK_LOOK", lineDatation, offset,
884                                                new LOSBuilder(lineOfSight).build());
885         
886         // in order not to have a problem when calling the pixelIsInside method (AtmosphericRefraction must be not null)
887         AtmosphericRefraction atmos = new AtmosphericRefraction() {
888             @Override
889             public NormalizedGeodeticPoint applyCorrection(Vector3D satPos, Vector3D satLos,
890                     NormalizedGeodeticPoint rawIntersection, IntersectionAlgorithm algorithm) {
891                 return rawIntersection;
892             }
893         };
894         atmos.deactivateComputation();
895         
896         Rugged rugged = new RuggedBuilder().
897                 setAlgorithm(AlgorithmId.IGNORE_DEM_USE_ELLIPSOID).
898                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
899                 setTimeSpan(satellitePVList.get(0).getDate(),
900                             satellitePVList.get(satellitePVList.size() - 1).getDate(), 0.1, 10.0).
901                 setTrajectory(InertialFrameId.EME2000,
902                               satellitePVList, 6, CartesianDerivativesFilter.USE_P,
903                               satelliteQList, 8, AngularDerivativesFilter.USE_R).
904                 addLineSensor(lineSensor).
905                 setRefractionCorrection(atmos).
906                 build();
907 
908         GeodeticPoint[] temp = rugged.directLocation("QUICK_LOOK", -250);
909         GeodeticPoint first = temp[0];
910         double minLon = first.getLongitude();
911         double minLat = first.getLatitude();
912         temp = rugged.directLocation("QUICK_LOOK", 350);
913         GeodeticPoint last = temp[temp.length - 1];
914         double maxLon = last.getLongitude();
915         double maxLat = last.getLatitude();
916 
917         GeodeticPoint gp = new GeodeticPoint((minLat + maxLat) / 2d, (minLon + maxLon) / 2d, 0d);
918         SensorPixel sensorPixel = rugged.inverseLocation("QUICK_LOOK", gp, -250, 350);
919 
920         Assertions.assertNotNull(sensorPixel);
921 
922         Assertions.assertFalse(inside(rugged, null, lineSensor));
923         Assertions.assertFalse(inside(rugged, new SensorPixel(-100, -100), lineSensor));
924         Assertions.assertFalse(inside(rugged, new SensorPixel(-100, +100), lineSensor));
925         Assertions.assertFalse(inside(rugged, new SensorPixel(+100, -100), lineSensor));
926         Assertions.assertFalse(inside(rugged, new SensorPixel(+100, +100), lineSensor));
927         Assertions.assertTrue(inside(rugged, new SensorPixel(0.2, 0.3), lineSensor));
928 
929     }
930 
931     private boolean inside(final Rugged rugged, final SensorPixel sensorPixel, LineSensor lineSensor) {
932         try {
933             final Method inside =
934                             Rugged.class.getDeclaredMethod("pixelIsInside",
935                                                            SensorPixel.class,
936                                                            LineSensor.class);
937             inside.setAccessible(true);
938             return ((Boolean) inside.invoke(rugged, sensorPixel, lineSensor));
939         } catch (NoSuchMethodException | IllegalAccessException |
940                  IllegalArgumentException | InvocationTargetException e) {
941             Assertions.fail(e.getLocalizedMessage());
942             return false;
943         }
944     }
945 
946     @Test
947     public void testInverseLocCurvedLine()
948         throws URISyntaxException {
949 
950         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
951         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
952         final BodyShape  earth = TestUtils.createEarth();
953         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
954 
955         AbsoluteDate crossing = new AbsoluteDate("2012-01-01T12:30:00.000", TimeScalesFactory.getUTC());
956         int dimension = 200;
957 
958         // one line sensor
959         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
960         // los: swath in the (YZ) plane, looking at nadir, 5.2" per pixel, 3" sagitta
961         Vector3D position = new Vector3D(1.5, 0, -0.2);
962         TimeDependentLOS los = TestUtils.createLOSCurvedLine(Vector3D.PLUS_K, Vector3D.PLUS_I,
963                                                              FastMath.toRadians((dimension/2.) * 5.2 / 3600.0),
964                                                              FastMath.toRadians(3.0 / 3600.0),
965                                                              dimension);
966 
967         // linear datation model: at reference time we get the middle line, and the rate is one line every 1.5ms
968         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
969         int firstLine = 0;
970         int lastLine  = dimension;
971         LineSensor lineSensor = new LineSensor("curved", lineDatation, position, los);
972         AbsoluteDate minDate = lineSensor.getDate(firstLine).shiftedBy(-1.0);
973         AbsoluteDate maxDate = lineSensor.getDate(lastLine).shiftedBy(+1.0);
974 
975         TileUpdater updater =
976                 new RandomLandscapeUpdater(0.0, 9000.0, 0.3, 0xf0a401650191f9f6L,
977                                            FastMath.toRadians(1.0), 257);
978 
979         Rugged rugged = new RuggedBuilder().
980                 setDigitalElevationModel(updater, 8).
981                 setAlgorithm(AlgorithmId.DUVENHAGE).
982                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
983                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
984                 setTrajectory(InertialFrameId.EME2000,
985                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
986                               8, CartesianDerivativesFilter.USE_PV,
987                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
988                               2, AngularDerivativesFilter.USE_R).
989                 addLineSensor(lineSensor).
990                 build();
991 
992         int lineNumber = 97;
993         GeodeticPoint[] gp = rugged.directLocation("curved", lineNumber);
994         for (int i = 0; i < gp.length; ++i) {
995             SensorPixel pixel = rugged.inverseLocation("curved", gp[i], firstLine, lastLine);
996             Assertions.assertEquals(lineNumber, pixel.getLineNumber(),  5.0e-4);
997             Assertions.assertEquals(i,          pixel.getPixelNumber(), 3.1e-7);
998         }
999 
1000     }
1001 
1002     private void checkInverseLocation(int dimension, boolean lightTimeCorrection, boolean aberrationOfLightCorrection,
1003                                       double maxLineError, double maxPixelError)
1004         throws URISyntaxException {
1005 
1006         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
1007         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
1008         final BodyShape  earth = TestUtils.createEarth();
1009         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
1010 
1011         AbsoluteDate crossing = new AbsoluteDate("2012-01-01T12:30:00.000", TimeScalesFactory.getUTC());
1012 
1013         // one line sensor
1014         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
1015         // los: swath in the (YZ) plane, looking at 50° roll, 5.2" per pixel
1016         Vector3D position = new Vector3D(1.5, 0, -0.2);
1017         TimeDependentLOS los = TestUtils.createLOSPerfectLine(new Rotation(Vector3D.PLUS_I,
1018                                                                            FastMath.toRadians(5.0),
1019                                                                            RotationConvention.VECTOR_OPERATOR).applyTo(Vector3D.PLUS_K),
1020                                                               Vector3D.PLUS_I,
1021                                                               FastMath.toRadians((dimension/2.) * 5.2 / 3600.0), dimension).build();
1022         // In fact the pixel size = 5.2" as we construct the LOS with the full line (dimension) instead of dimension/2
1023 
1024         // linear datation model: at reference time we get the middle line, and the rate is one line every 1.5ms
1025         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
1026         int firstLine = 0;
1027         int lastLine  = dimension;
1028         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
1029         AbsoluteDate minDate = lineSensor.getDate(firstLine).shiftedBy(-1.0);
1030         AbsoluteDate maxDate = lineSensor.getDate(lastLine).shiftedBy(+1.0);
1031 
1032         TileUpdater updater =
1033                 new RandomLandscapeUpdater(0.0, 9000.0, 0.3, 0xf0a401650191f9f6L,
1034                                            FastMath.toRadians(1.0), 257);
1035 
1036         Rugged rugged = new RuggedBuilder().
1037                 setDigitalElevationModel(updater, 8).
1038                 setAlgorithm(AlgorithmId.DUVENHAGE).
1039                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
1040                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
1041                 setTrajectory(InertialFrameId.EME2000,
1042                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
1043                               8, CartesianDerivativesFilter.USE_PV,
1044                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
1045                               2, AngularDerivativesFilter.USE_R).
1046                 setLightTimeCorrection(lightTimeCorrection).
1047                 setAberrationOfLightCorrection(aberrationOfLightCorrection).
1048                 addLineSensor(lineSensor).
1049                 build();
1050 
1051         double referenceLine = 0.87654 * dimension;
1052         GeodeticPoint[] gp = rugged.directLocation("line", referenceLine);
1053 
1054         for (double p = 0; p < gp.length - 1; p += 1.0) {
1055             int    i = (int) FastMath.floor(p);
1056             double d = p - i;
1057             SensorPixel sp = rugged.inverseLocation("line",
1058                                                         (1 - d) * gp[i].getLatitude()  + d * gp[i + 1].getLatitude(),
1059                                                         (1 - d) * gp[i].getLongitude() + d * gp[i + 1].getLongitude(),
1060                                                         0, dimension);
1061 
1062             Assertions.assertEquals(referenceLine, sp.getLineNumber(),  maxLineError);
1063             Assertions.assertEquals(p,             sp.getPixelNumber(), maxPixelError);
1064         }
1065 
1066         // point out of line (20 pixels before first pixel)
1067         Assertions.assertNull(rugged.inverseLocation("line",
1068                                                     21 * gp[0].getLatitude()  - 20 * gp[1].getLatitude(),
1069                                                     21 * gp[0].getLongitude() - 20 * gp[1].getLongitude(),
1070                                                     0, dimension));
1071 
1072         // point out of line (20 pixels after last pixel)
1073         Assertions.assertNull(rugged.inverseLocation("line",
1074                                                     -20 * gp[gp.length - 2].getLatitude()  + 21 * gp[gp.length - 1].getLatitude(),
1075                                                     -20 * gp[gp.length - 2].getLongitude() + 21 * gp[gp.length - 1].getLongitude(),
1076                                                     0, dimension));
1077 
1078         // point out of line (20 lines before first line)
1079         GeodeticPoint[] gp0 = rugged.directLocation("line", 0);
1080         GeodeticPoint[] gp1 = rugged.directLocation("line", 1);
1081         Assertions.assertNull(rugged.inverseLocation("line",
1082                                                     21 * gp0[dimension / 2].getLatitude()  - 20 * gp1[dimension / 2].getLatitude(),
1083                                                     21 * gp0[dimension / 2].getLongitude() - 20 * gp1[dimension / 2].getLongitude(),
1084                                                     0, dimension));
1085 
1086         // point out of line (20 lines after last line)
1087         GeodeticPoint[] gp2 = rugged.directLocation("line", dimension - 2);
1088         GeodeticPoint[] gp3 = rugged.directLocation("line", dimension - 1);
1089         Assertions.assertNull(rugged.inverseLocation("line",
1090                                                     -20 * gp2[dimension / 2].getLatitude()  + 21 * gp3[dimension / 2].getLatitude(),
1091                                                     -20 * gp2[dimension / 2].getLongitude() + 21 * gp3[dimension / 2].getLongitude(),
1092                                                     0, dimension));
1093 
1094     }
1095 
1096     @Test
1097     public void testInverseLocationDerivativesWithoutCorrections()
1098         {
1099         doTestInverseLocationDerivatives(2000, false, false,
1100                                          8.0e-9, 5.0e-10, 5.0e-12, 9.0e-8);
1101     }
1102 
1103     @Test
1104     public void testInverseLocationDerivativesWithLightTimeCorrection()
1105         {
1106         doTestInverseLocationDerivatives(2000, true, false,
1107                                          3.0e-9, 9.0e-9, 2.8e-12, 7e-6);
1108     }
1109 
1110     @Test
1111     public void testInverseLocationDerivativesWithAberrationOfLightCorrection()
1112         {
1113         doTestInverseLocationDerivatives(2000, false, true,
1114                                          1e-9, 3.0e-10, 3.4e-12, 8.0e-8);
1115     }
1116 
1117     @Test
1118     public void testInverseLocationDerivativesWithAllCorrections()
1119         {
1120         doTestInverseLocationDerivatives(2000, true, true,
1121                                          1e-8, 5.0e-10, 2.0e-12, 7.0e-8);
1122     }
1123 
1124     /**
1125      * @param dimension
1126      * @param lightTimeCorrection
1127      * @param aberrationOfLightCorrection
1128      * @param lineTolerance
1129      * @param pixelTolerance
1130      * @param lineDerivativeRelativeTolerance
1131      * @param pixelDerivativeRelativeTolerance
1132      */
1133     private void doTestInverseLocationDerivatives(int dimension,
1134                                                   boolean lightTimeCorrection,
1135                                                   boolean aberrationOfLightCorrection,
1136                                                   double lineTolerance,
1137                                                   double pixelTolerance,
1138                                                   double lineDerivativeRelativeTolerance,
1139                                                   double pixelDerivativeRelativeTolerance)
1140         {
1141         try {
1142 
1143             String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
1144             DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
1145             final BodyShape  earth = TestUtils.createEarth();
1146             final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
1147 
1148             AbsoluteDate crossing = new AbsoluteDate("2012-01-01T12:30:00.000", TimeScalesFactory.getUTC());
1149 
1150             // one line sensor
1151             // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
1152             // los: swath in the (YZ) plane, looking at 50° roll, 5.2" per pixel
1153             Vector3D position = new Vector3D(1.5, 0, -0.2);
1154             LOSBuilder losBuilder =
1155              TestUtils.createLOSPerfectLine(new Rotation(Vector3D.PLUS_I,
1156                                                          FastMath.toRadians(50.0),
1157                                                          RotationConvention.VECTOR_OPERATOR).applyTo(Vector3D.PLUS_K),
1158                                             Vector3D.PLUS_I,
1159                                             FastMath.toRadians((dimension/2.) * 5.2 / 3600.0), dimension);
1160             losBuilder.addTransform(new FixedRotation("roll",  Vector3D.MINUS_I, 0.0));
1161             losBuilder.addTransform(new FixedRotation("pitch", Vector3D.MINUS_J, 0.0));
1162             TimeDependentLOS los = losBuilder.build();
1163             // In fact the pixel size = 5.2" as we construct the LOS with the full line (dimension) instead of dimension/2
1164 
1165             // linear datation model: at reference time we get the middle line, and the rate is one line every 1.5ms
1166             LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
1167             int firstLine = 0;
1168             int lastLine  = dimension;
1169             final LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
1170             AbsoluteDate minDate = lineSensor.getDate(firstLine).shiftedBy(-1.0);
1171             AbsoluteDate maxDate = lineSensor.getDate(lastLine).shiftedBy(+1.0);
1172 
1173             TileUpdater updater =
1174                             new RandomLandscapeUpdater(0.0, 9000.0, 0.3, 0xf0a401650191f9f6L,
1175                                                        FastMath.toRadians(1.0), 257);
1176 
1177             Rugged rugged = new RuggedBuilder().
1178                             setDigitalElevationModel(updater, 8).
1179                             setAlgorithm(AlgorithmId.DUVENHAGE).
1180                             setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
1181                             setTimeSpan(minDate, maxDate, 0.001, 5.0).
1182                             setTrajectory(InertialFrameId.EME2000,
1183                                           TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
1184                                           8, CartesianDerivativesFilter.USE_PV,
1185                                           TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
1186                                           2, AngularDerivativesFilter.USE_R).
1187                             setLightTimeCorrection(lightTimeCorrection).
1188                             setAberrationOfLightCorrection(aberrationOfLightCorrection).
1189                             addLineSensor(lineSensor).
1190                             build();
1191 
1192             // we want to adjust sensor roll and pitch angles
1193             ParameterDriver rollDriver =
1194                             lineSensor.getParametersDrivers().
1195                             filter(driver -> driver.getName().equals("roll")).findFirst().get();
1196             rollDriver.setSelected(true);
1197             ParameterDriver pitchDriver =
1198                             lineSensor.getParametersDrivers().
1199                             filter(driver -> driver.getName().equals("pitch")).findFirst().get();
1200             pitchDriver.setSelected(true);
1201 
1202             // prepare generator
1203             final Observables measurements = new Observables(1);
1204             GroundOptimizationProblemBuilder optimizationPbBuilder = new GroundOptimizationProblemBuilder(Collections.singletonList(lineSensor),
1205                                                                                                                measurements, rugged);
1206             
1207             java.lang.reflect.Method getGenerator = GroundOptimizationProblemBuilder.class.getSuperclass().getDeclaredMethod("getGenerator");
1208             getGenerator.setAccessible(true);
1209 
1210             @SuppressWarnings("unchecked")
1211             DerivativeGenerator<Gradient> generator = (DerivativeGenerator<Gradient>) getGenerator.invoke(optimizationPbBuilder);
1212 
1213             double referenceLine = 0.87654 * dimension;
1214             GeodeticPoint[] gp = rugged.directLocation("line", referenceLine);
1215 
1216             Method inverseLoc = Rugged.class.getDeclaredMethod("inverseLocationDerivatives",
1217                                                                String.class, GeodeticPoint.class,
1218                                                                Integer.TYPE, Integer.TYPE,
1219                                                                DerivativeGenerator.class);
1220             inverseLoc.setAccessible(true);
1221             int referencePixel = (3 * dimension) / 4;
1222             Gradient[] result = 
1223                             (Gradient[]) inverseLoc.invoke(rugged,
1224                                                            "line", gp[referencePixel], 0, dimension,
1225                                                            generator);
1226             Assertions.assertEquals(referenceLine,  result[0].getValue(), lineTolerance);
1227             Assertions.assertEquals(referencePixel, result[1].getValue(), pixelTolerance);
1228             Assertions.assertEquals(2, result[0].getFreeParameters());
1229             Assertions.assertEquals(1, result[0].getOrder());
1230 
1231             // check the partial derivatives
1232             DSFactory factory = new DSFactory(1, 1);
1233             double h = 1.0e-6;
1234             FiniteDifferencesDifferentiator differentiator = new FiniteDifferencesDifferentiator(8, h);
1235 
1236             UnivariateDifferentiableFunction lineVSroll =
1237                             differentiator.differentiate((double roll) -> {
1238                                 rollDriver.setValue(roll);
1239                                 pitchDriver.setValue(0);
1240                                 return rugged.inverseLocation("line", gp[referencePixel], 0, dimension).getLineNumber();
1241                             });
1242             double dLdR = lineVSroll.value(factory.variable(0, 0.0)).getPartialDerivative(1);
1243             Assertions.assertEquals(dLdR, result[0].getPartialDerivative(1, 0), dLdR * lineDerivativeRelativeTolerance);
1244 
1245             UnivariateDifferentiableFunction lineVSpitch =
1246                             differentiator.differentiate((double pitch) -> {
1247                                 rollDriver.setValue(0);
1248                                 pitchDriver.setValue(pitch);
1249                                 return rugged.inverseLocation("line", gp[referencePixel], 0, dimension).getLineNumber();
1250                             });
1251             double dLdP = lineVSpitch.value(factory.variable(0, 0.0)).getPartialDerivative(1);
1252             Assertions.assertEquals(dLdP, result[0].getPartialDerivative(0, 1), dLdP * lineDerivativeRelativeTolerance);
1253 
1254             UnivariateDifferentiableFunction pixelVSroll =
1255                             differentiator.differentiate((double roll) -> {
1256                                 rollDriver.setValue(roll);
1257                                 pitchDriver.setValue(0);
1258                                 return rugged.inverseLocation("line", gp[referencePixel], 0, dimension).getPixelNumber();
1259                             });
1260             double dXdR = pixelVSroll.value(factory.variable(0, 0.0)).getPartialDerivative(1);
1261             Assertions.assertEquals(dXdR, result[1].getPartialDerivative(1, 0), dXdR * pixelDerivativeRelativeTolerance);
1262 
1263             UnivariateDifferentiableFunction pixelVSpitch =
1264                             differentiator.differentiate((double pitch) -> {
1265                                 rollDriver.setValue(0);
1266                                 pitchDriver.setValue(pitch);
1267                                 return rugged.inverseLocation("line", gp[referencePixel], 0, dimension).getPixelNumber();
1268                             });
1269             double dXdP = pixelVSpitch.value(factory.variable(0, 0.0)).getPartialDerivative(1);
1270             Assertions.assertEquals(dXdP, result[1].getPartialDerivative(0, 1), dXdP * pixelDerivativeRelativeTolerance);
1271 
1272         } catch (InvocationTargetException | NoSuchMethodException |
1273                 SecurityException | IllegalAccessException |
1274                 IllegalArgumentException | URISyntaxException |
1275                 OrekitException | RuggedException e) {
1276             Assertions.fail(e.getLocalizedMessage());
1277         }
1278     }
1279 
1280 
1281     private void checkDateLocation(int dimension, boolean lightTimeCorrection, boolean aberrationOfLightCorrection,
1282                                        double maxDateError)
1283         throws URISyntaxException {
1284 
1285         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
1286         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
1287         final BodyShape  earth = TestUtils.createEarth();
1288         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
1289 
1290         AbsoluteDate crossing = new AbsoluteDate("2012-01-01T12:30:00.000", TimeScalesFactory.getUTC());
1291 
1292         // one line sensor
1293         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
1294         // los: swath in the (YZ) plane, looking at 50° roll, 5.2" per pixel
1295         Vector3D position = new Vector3D(1.5, 0, -0.2);
1296         TimeDependentLOS los = TestUtils.createLOSPerfectLine(new Rotation(Vector3D.PLUS_I,
1297                                                                            FastMath.toRadians(50.0),
1298                                                                            RotationConvention.VECTOR_OPERATOR).applyTo(Vector3D.PLUS_K),
1299                                                               Vector3D.PLUS_I,
1300                                                               FastMath.toRadians((dimension/2.) * 5.2 / 3600.0), dimension).build();
1301         // In fact the pixel size = 5.2" as we construct the LOS with the full line (dimension) instead of dimension/2
1302 
1303         // linear datation model: at reference time we get line 100, and the rate is one line every 1.5ms
1304         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
1305         int firstLine = 0;
1306         int lastLine  = dimension;
1307         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
1308         AbsoluteDate minDate = lineSensor.getDate(firstLine).shiftedBy(-1.0);
1309         AbsoluteDate maxDate = lineSensor.getDate(lastLine).shiftedBy(+1.0);
1310 
1311         TileUpdater updater =
1312                 new RandomLandscapeUpdater(0.0, 9000.0, 0.5, 0xf0a401650191f9f6L,
1313                                            FastMath.toRadians(1.0), 257);
1314 
1315         Rugged rugged = new RuggedBuilder().
1316                 setDigitalElevationModel(updater, 8).
1317                 setAlgorithm(AlgorithmId.DUVENHAGE).
1318                 setEllipsoid(EllipsoidId.WGS84, BodyRotatingFrameId.ITRF).
1319                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
1320                 setTrajectory(InertialFrameId.EME2000,
1321                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
1322                               8, CartesianDerivativesFilter.USE_PV,
1323                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
1324                               2, AngularDerivativesFilter.USE_R).
1325                 setLightTimeCorrection(lightTimeCorrection).
1326                 setAberrationOfLightCorrection(aberrationOfLightCorrection).
1327                 addLineSensor(lineSensor).
1328                 build();
1329 
1330         double referenceLine = 0.87654 * dimension;
1331         GeodeticPoint[] gp = rugged.directLocation("line", referenceLine);
1332 
1333         for (double p = 0; p < gp.length - 1; p += 1.0) {
1334             int    i = (int) FastMath.floor(p);
1335             double d = p - i;
1336             AbsoluteDate date = rugged.dateLocation("line",
1337                                                         (1 - d) * gp[i].getLatitude()  + d * gp[i + 1].getLatitude(),
1338                                                         (1 - d) * gp[i].getLongitude() + d * gp[i + 1].getLongitude(),
1339                                                         0, dimension);
1340             Assertions.assertEquals(0.0, date.durationFrom(lineSensor.getDate(referenceLine)),  maxDateError);
1341         }
1342 
1343         // point out of line (20 lines before first line)
1344         GeodeticPoint[] gp0 = rugged.directLocation("line", 0);
1345         GeodeticPoint[] gp1 = rugged.directLocation("line", 1);
1346         Assertions.assertNull(rugged.dateLocation("line",
1347                                                     21 * gp0[dimension / 2].getLatitude()  - 20 * gp1[dimension / 2].getLatitude(),
1348                                                     21 * gp0[dimension / 2].getLongitude() - 20 * gp1[dimension / 2].getLongitude(),
1349                                                     0, dimension));
1350 
1351         // point out of line (20 lines after last line)
1352         GeodeticPoint[] gp2 = rugged.directLocation("line", dimension - 2);
1353         GeodeticPoint[] gp3 = rugged.directLocation("line", dimension - 1);
1354         Assertions.assertNull(rugged.dateLocation("line",
1355                                                     -20 * gp2[dimension / 2].getLatitude()  + 21 * gp3[dimension / 2].getLatitude(),
1356                                                     -20 * gp2[dimension / 2].getLongitude() + 21 * gp3[dimension / 2].getLongitude(),
1357                                                     0, dimension));
1358         
1359     }
1360 
1361     private void checkLineDatation(int dimension, double maxLineError)
1362                     throws URISyntaxException {
1363 
1364         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
1365         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
1366 
1367         AbsoluteDate crossing = new AbsoluteDate("2012-01-01T12:30:00.000", TimeScalesFactory.getUTC());
1368 
1369         // one line sensor
1370         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
1371         // los: swath in the (YZ) plane, looking at 50° roll, 2.6" per pixel
1372         Vector3D position = new Vector3D(1.5, 0, -0.2);
1373         TimeDependentLOS los = TestUtils.createLOSPerfectLine(new Rotation(Vector3D.PLUS_I,
1374                 FastMath.toRadians(50.0),
1375                 RotationConvention.VECTOR_OPERATOR).applyTo(Vector3D.PLUS_K),
1376                 Vector3D.PLUS_I,
1377                 FastMath.toRadians(dimension * 2.6 / 3600.0), dimension).build();
1378         // In fact the pixel size = 5.2" as we construct the LOS with the full line (dimension) instead of dimension/2
1379 
1380         // linear datation model: at reference time we get the middle line, and the rate is one line every 1.5ms
1381         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
1382         int firstLine = 0;
1383         int lastLine  = dimension;
1384         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
1385         AbsoluteDate minDate = lineSensor.getDate(firstLine).shiftedBy(-1.0);
1386         AbsoluteDate maxDate = lineSensor.getDate(lastLine).shiftedBy(+1.0);
1387 
1388         // Recompute the lines from the date with the appropriate shift of date
1389         double recomputedFirstLine = lineSensor.getLine(minDate.shiftedBy(+1.0));
1390         double recomputedLastLine = lineSensor.getLine(maxDate.shiftedBy(-1.0));
1391 
1392         Assertions.assertEquals(firstLine, recomputedFirstLine, maxLineError);
1393         Assertions.assertEquals(lastLine, recomputedLastLine, maxLineError);
1394     }
1395 
1396     @Test
1397     public void testDistanceBetweenLOS() {
1398         
1399         InitInterRefiningTest refiningTest = new InitInterRefiningTest();
1400         refiningTest.initRefiningTest();
1401 
1402         final SensorPixel realPixelA = new SensorPixel(2005.015883575199, 18004.968656395424);
1403         final SensorPixel realPixelB = new SensorPixel(4798.487736488162, 13952.2195710654);
1404 
1405         double[] distancesBetweenLOS = refiningTest.computeDistancesBetweenLOS(realPixelA, realPixelB);
1406         
1407         double expectedDistanceBetweenLOS = 3.88800245;
1408         double expectedDistanceToTheGround = 6368020.559109;
1409 
1410         Assertions.assertEquals(expectedDistanceBetweenLOS, distancesBetweenLOS[0], 2.e-6);
1411         Assertions.assertEquals(expectedDistanceToTheGround, distancesBetweenLOS[1], 3.e-5);
1412      }
1413 
1414     @Test
1415     public void testDistanceBetweenLOSDerivatives() throws NoSuchMethodException, SecurityException, IllegalAccessException, IllegalArgumentException, InvocationTargetException {
1416         
1417         InitInterRefiningTest refiningTest = new InitInterRefiningTest();
1418         refiningTest.initRefiningTest();
1419 
1420         final SensorPixel realPixelA = new SensorPixel(2005.015883575199, 18004.968656395424);
1421         final SensorPixel realPixelB = new SensorPixel(4798.487736488162, 13952.2195710654);
1422 
1423         // Expected distances between LOS and to the ground
1424         double expectedDistanceBetweenLOS = 3.88800245;
1425         double expectedDistanceToTheGround = 6368020.559109;
1426 
1427         // Expected derivatives for
1428         // minimum distance between LOS
1429         double[] expectedDminDerivatives = {-153874.01319097, -678866.03112033, 191294.06938169, 668600.16715270} ;
1430         // minimum distance to the ground
1431         double[] expectedDcentralBodyDerivatives = {7007767.46926062, -1577060.82402054, -6839286.39593802, 1956452.66636262};
1432 
1433         Gradient[] distancesBetweenLOSGradient = refiningTest.computeDistancesBetweenLOSGradient(realPixelA, realPixelB, expectedDistanceBetweenLOS, expectedDistanceToTheGround);
1434 
1435         // Minimum distance between LOS
1436         Gradient dMin = distancesBetweenLOSGradient[0];
1437         // Minimum distance to the ground
1438         Gradient dCentralBody = distancesBetweenLOSGradient[1];
1439 
1440         Assertions.assertEquals(expectedDistanceBetweenLOS, dMin.getValue(), 2.e-6);
1441         Assertions.assertEquals(expectedDistanceToTheGround, dCentralBody.getValue() , 4.e-5);
1442 
1443 
1444         for (int i = 0; i < dMin.getFreeParameters(); i++) {
1445             Assertions.assertEquals(expectedDminDerivatives[i], dMin.getPartialDerivative(i), 3e-5);
1446         }
1447 
1448         for (int i = 0; i < dCentralBody.getFreeParameters(); i++) {
1449             Assertions.assertEquals(expectedDcentralBodyDerivatives[i], dCentralBody.getPartialDerivative(i), 3.e-4);
1450         }
1451     }
1452 
1453     @Test
1454     public void testForCoverage() throws URISyntaxException {
1455 
1456         int dimension = 400;
1457 
1458         String path = getClass().getClassLoader().getResource("orekit-data").toURI().getPath();
1459         DataContext.getDefault().getDataProvidersManager().addProvider(new DirectoryCrawler(new File(path)));
1460         final BodyShape  earth = TestUtils.createEarth();
1461         final Orbit      orbit = TestUtils.createOrbit(Constants.EIGEN5C_EARTH_MU);
1462 
1463         AbsoluteDate crossing = new AbsoluteDate("2012-01-07T11:21:15.000", TimeScalesFactory.getUTC());
1464 
1465         // one line sensor
1466         // position: 1.5m in front (+X) and 20 cm above (-Z) of the S/C center of mass
1467         // los: swath in the (YZ) plane, centered at +Z, ±10° aperture, 960 pixels
1468         Vector3D position = new Vector3D(1.5, 0, -0.2);
1469         TimeDependentLOS los = TestUtils.createLOSPerfectLine(Vector3D.PLUS_K, Vector3D.PLUS_I,
1470                                                               FastMath.toRadians(10.0), dimension).build();
1471 
1472         // linear datation model: at reference time we get line 200, and the rate is one line every 1.5ms
1473         LineDatation lineDatation = new LinearLineDatation(crossing, dimension / 2, 1.0 / 1.5e-3);
1474         int firstLine = 0;
1475         int lastLine  = dimension;
1476         LineSensor lineSensor = new LineSensor("line", lineDatation, position, los);
1477         AbsoluteDate minDate = lineSensor.getDate(firstLine);
1478         AbsoluteDate maxDate = lineSensor.getDate(lastLine);
1479 
1480         RuggedBuilder builder = new RuggedBuilder().
1481                 setAlgorithm(AlgorithmId.IGNORE_DEM_USE_ELLIPSOID).
1482                 setEllipsoid(EllipsoidId.IERS2003, BodyRotatingFrameId.ITRF).
1483                 setTimeSpan(minDate, maxDate, 0.001, 5.0).
1484                 setTrajectory(InertialFrameId.EME2000,
1485                               TestUtils.orbitToPV(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
1486                               8, CartesianDerivativesFilter.USE_PV,
1487                               TestUtils.orbitToQ(orbit, earth, minDate.shiftedBy(-1.0), maxDate.shiftedBy(+1.0), 0.25),
1488                               2, AngularDerivativesFilter.USE_R).
1489                 addLineSensor(lineSensor);
1490 
1491         Rugged rugged = builder.build();
1492         
1493         
1494         // Check builder 
1495         assertTrue(builder.getName().equalsIgnoreCase("Rugged"));
1496         
1497         // Check a date in the range of minDate - maxDate
1498         AbsoluteDate midddleDate = lineSensor.getDate((firstLine+lastLine)/2.);
1499 
1500         assertTrue(rugged.isInRange(midddleDate));
1501         
1502         // Get the algorithm
1503         assertTrue(rugged.getAlgorithm().getClass().isInstance(new IgnoreDEMAlgorithm()));
1504         
1505         // Get the algorithm Id
1506         assertEquals(AlgorithmId.IGNORE_DEM_USE_ELLIPSOID, rugged.getAlgorithmId());
1507         
1508         // Change the min and max line in inverse location to update the SensorMeanPlaneCrossing when the planeCrossing is not null
1509         int minLine = firstLine;
1510         int maxLine = lastLine;
1511         double line = (firstLine + lastLine)/2.;
1512         double pixel = dimension/2.;
1513         AbsoluteDate date = lineSensor.getDate(line);
1514         Vector3D pixelLos = lineSensor.getLOS(date, pixel);
1515         GeodeticPoint gp = rugged.directLocation(date, position, pixelLos);
1516         
1517         SensorPixel sp = rugged.inverseLocation("line", gp, minLine, maxLine);
1518         int minLineNew = minLine + 10;
1519         int maxLineNew = maxLine - 10;
1520         SensorPixel spChangeLines = rugged.inverseLocation("line", gp, minLineNew, maxLineNew);
1521         
1522         assertEquals(sp.getPixelNumber(), spChangeLines.getPixelNumber(), 1.e-9);
1523         assertEquals(sp.getLineNumber(), spChangeLines.getLineNumber(), 1.e-9);
1524         
1525         // For computeInverseLocOnGridWithoutAtmosphere special cases
1526         try {
1527             java.lang.reflect.Method computeWithoutAtmosphere = 
1528                     rugged.getClass().getDeclaredMethod("computeInverseLocOnGridWithoutAtmosphere",
1529                                                         GeodeticPoint[][].class,
1530                                                         Integer.TYPE, Integer.TYPE,
1531                                                         LineSensor.class, Integer.TYPE, Integer.TYPE);
1532             computeWithoutAtmosphere.setAccessible(true);
1533             final int nbPixelGrid = 2; 
1534             final int nbLineGrid = 2;
1535             GeodeticPoint[][] groundGridWithAtmosphere = new GeodeticPoint[nbPixelGrid][nbLineGrid];
1536             for (int i = 0; i < nbPixelGrid; i++) {
1537                 for (int j = 0; j < nbLineGrid; j++) {
1538                     groundGridWithAtmosphere[i][j] = null;
1539                 }
1540             }
1541              
1542             SensorPixel[][] spNull = (SensorPixel[][]) computeWithoutAtmosphere.invoke(rugged, groundGridWithAtmosphere, nbPixelGrid, nbLineGrid, lineSensor, minLine, maxLine);
1543             for (int i = 0; i < nbPixelGrid; i++) {
1544                 for (int j = 0; j < nbLineGrid; j++) {
1545                     assertNull(spNull[i][j]);
1546                 }
1547             }
1548         } catch (NoSuchMethodException | SecurityException | 
1549                 IllegalAccessException | IllegalArgumentException | InvocationTargetException e) {
1550             Assertions.fail(e.getLocalizedMessage());
1551         }
1552     }
1553     
1554     
1555     @BeforeEach
1556     public void setUp() throws URISyntaxException {
1557         TestUtils.clearFactories();
1558     }
1559 }
1560