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