1   /* Copyright 2002-2021 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.propagation.semianalytical.dsst.forces;
18  
19  import java.util.List;
20  
21  import org.hipparchus.Field;
22  import org.hipparchus.CalculusFieldElement;
23  import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
24  import org.hipparchus.geometry.euclidean.threed.Vector3D;
25  import org.hipparchus.util.FastMath;
26  import org.hipparchus.util.MathArrays;
27  import org.hipparchus.util.MathUtils;
28  import org.hipparchus.util.Precision;
29  import org.orekit.forces.radiation.IsotropicRadiationSingleCoefficient;
30  import org.orekit.forces.radiation.RadiationSensitive;
31  import org.orekit.forces.radiation.SolarRadiationPressure;
32  import org.orekit.propagation.FieldSpacecraftState;
33  import org.orekit.propagation.SpacecraftState;
34  import org.orekit.propagation.events.EventDetector;
35  import org.orekit.propagation.events.FieldEventDetector;
36  import org.orekit.propagation.semianalytical.dsst.utilities.AuxiliaryElements;
37  import org.orekit.propagation.semianalytical.dsst.utilities.FieldAuxiliaryElements;
38  import org.orekit.utils.ExtendedPVCoordinatesProvider;
39  import org.orekit.utils.ParameterDriver;
40  
41  /** Solar radiation pressure contribution to the
42   *  {@link org.orekit.propagation.semianalytical.dsst.DSSTPropagator DSSTPropagator}.
43   *  <p>
44   *  The solar radiation pressure acceleration is computed through the
45   *  acceleration model of
46   *  {@link org.orekit.forces.radiation.SolarRadiationPressure SolarRadiationPressure}.
47   *  </p>
48   *
49   *  @author Pascal Parraud
50   */
51  public class DSSTSolarRadiationPressure extends AbstractGaussianContribution {
52  
53      /** Reference distance for the solar radiation pressure (m). */
54      private static final double D_REF = 149597870000.0;
55  
56      /** Reference solar radiation pressure at D_REF (N/m²). */
57      private static final double P_REF = 4.56e-6;
58  
59      /** Threshold for the choice of the Gauss quadrature order. */
60      private static final double GAUSS_THRESHOLD = 1.0e-15;
61  
62      /** Threshold for shadow equation. */
63      private static final double S_ZERO = 1.0e-6;
64  
65      /** Prefix for the coefficient keys. */
66      private static final String PREFIX = "DSST-SRP-";
67  
68      /** Sun model. */
69      private final ExtendedPVCoordinatesProvider sun;
70  
71      /** Central Body radius. */
72      private final double                ae;
73  
74      /** Spacecraft model for radiation acceleration computation. */
75      private final RadiationSensitive spacecraft;
76  
77      /**
78       * Simple constructor with default reference values and spherical spacecraft.
79       * <p>
80       * When this constructor is used, the reference values are:
81       * </p>
82       * <ul>
83       * <li>d<sub>ref</sub> = 149597870000.0 m</li>
84       * <li>p<sub>ref</sub> = 4.56 10<sup>-6</sup> N/m²</li>
85       * </ul>
86       * <p>
87       * The spacecraft has a spherical shape.
88       * </p>
89       *
90       * @param cr satellite radiation pressure coefficient (assuming total specular reflection)
91       * @param area cross sectional area of satellite
92       * @param sun Sun model
93       * @param equatorialRadius central body equatorial radius (for shadow computation)
94       * @param mu central attraction coefficient
95       * @since 9.2
96       */
97      public DSSTSolarRadiationPressure(final double cr, final double area,
98                                        final ExtendedPVCoordinatesProvider sun,
99                                        final double equatorialRadius,
100                                       final double mu) {
101         this(D_REF, P_REF, cr, area, sun, equatorialRadius, mu);
102     }
103 
104     /**
105      * Simple constructor with default reference values, but custom spacecraft.
106      * <p>
107      * When this constructor is used, the reference values are:
108      * </p>
109      * <ul>
110      * <li>d<sub>ref</sub> = 149597870000.0 m</li>
111      * <li>p<sub>ref</sub> = 4.56 10<sup>-6</sup> N/m²</li>
112      * </ul>
113      *
114      * @param sun Sun model
115      * @param equatorialRadius central body equatorial radius (for shadow computation)
116      * @param spacecraft spacecraft model
117      * @param mu central attraction coefficient
118      * @since 9.2
119      */
120     public DSSTSolarRadiationPressure(final ExtendedPVCoordinatesProvider sun,
121                                       final double equatorialRadius,
122                                       final RadiationSensitive spacecraft,
123                                       final double mu) {
124         this(D_REF, P_REF, sun, equatorialRadius, spacecraft, mu);
125     }
126 
127     /**
128      * Constructor with customizable reference values but spherical spacecraft.
129      * <p>
130      * Note that reference solar radiation pressure <code>pRef</code> in N/m² is linked
131      * to solar flux SF in W/m² using formula pRef = SF/c where c is the speed of light
132      * (299792458 m/s). So at 1UA a 1367 W/m² solar flux is a 4.56 10<sup>-6</sup>
133      * N/m² solar radiation pressure.
134      * </p>
135      *
136      * @param dRef reference distance for the solar radiation pressure (m)
137      * @param pRef reference solar radiation pressure at dRef (N/m²)
138      * @param cr satellite radiation pressure coefficient (assuming total specular reflection)
139      * @param area cross sectional area of satellite
140      * @param sun Sun model
141      * @param equatorialRadius central body equatorial radius (for shadow computation)
142      * @param mu central attraction coefficient
143      * @since 9.2
144      */
145     public DSSTSolarRadiationPressure(final double dRef, final double pRef,
146                                       final double cr, final double area,
147                                       final ExtendedPVCoordinatesProvider sun,
148                                       final double equatorialRadius,
149                                       final double mu) {
150 
151         // cR being the DSST SRP coef and assuming a spherical spacecraft,
152         // the conversion is:
153         // cR = 1 + (1 - kA) * (1 - kR) * 4 / 9
154         // with kA arbitrary sets to 0
155         this(dRef, pRef, sun, equatorialRadius, new IsotropicRadiationSingleCoefficient(area, cr), mu);
156     }
157 
158     /**
159      * Complete constructor.
160      * <p>
161      * Note that reference solar radiation pressure <code>pRef</code> in N/m² is linked
162      * to solar flux SF in W/m² using formula pRef = SF/c where c is the speed of light
163      * (299792458 m/s). So at 1UA a 1367 W/m² solar flux is a 4.56 10<sup>-6</sup>
164      * N/m² solar radiation pressure.
165      * </p>
166      *
167      * @param dRef reference distance for the solar radiation pressure (m)
168      * @param pRef reference solar radiation pressure at dRef (N/m²)
169      * @param sun Sun model
170      * @param equatorialRadius central body equatorial radius (for shadow computation)
171      * @param spacecraft spacecraft model
172      * @param mu central attraction coefficient
173      * @since 9.2
174      */
175     public DSSTSolarRadiationPressure(final double dRef, final double pRef,
176                                       final ExtendedPVCoordinatesProvider sun,
177                                       final double equatorialRadius,
178                                       final RadiationSensitive spacecraft,
179                                       final double mu) {
180 
181         //Call to the constructor from superclass using the numerical SRP model as ForceModel
182         super(PREFIX, GAUSS_THRESHOLD,
183               new SolarRadiationPressure(dRef, pRef, sun, equatorialRadius, spacecraft), mu);
184 
185         this.sun  = sun;
186         this.ae   = equatorialRadius;
187         this.spacecraft = spacecraft;
188     }
189 
190     /** Get spacecraft shape.
191      * @return the spacecraft shape.
192      */
193     public RadiationSensitive getSpacecraft() {
194         return spacecraft;
195     }
196 
197     /** {@inheritDoc} */
198     public EventDetector[] getEventsDetectors() {
199         return null;
200     }
201 
202     /** {@inheritDoc} */
203     @Override
204     public <T extends CalculusFieldElement<T>> FieldEventDetector<T>[] getFieldEventsDetectors(final Field<T> field) {
205         return null;
206     }
207 
208     /** {@inheritDoc} */
209     protected List<ParameterDriver> getParametersDriversWithoutMu() {
210         return spacecraft.getRadiationParametersDrivers();
211     }
212 
213     /** {@inheritDoc} */
214     protected double[] getLLimits(final SpacecraftState state, final AuxiliaryElements auxiliaryElements) {
215 
216         // Default bounds without shadow [-PI, PI]
217         final double[] ll = {-FastMath.PI + MathUtils.normalizeAngle(state.getLv(), 0),
218                              FastMath.PI + MathUtils.normalizeAngle(state.getLv(), 0)};
219 
220         // Direction cosines of the Sun in the equinoctial frame
221         final Vector3D sunDir = sun.getPVCoordinates(state.getDate(), state.getFrame()).getPosition().normalize();
222         final double alpha = sunDir.dotProduct(auxiliaryElements.getVectorF());
223         final double beta  = sunDir.dotProduct(auxiliaryElements.getVectorG());
224         final double gamma = sunDir.dotProduct(auxiliaryElements.getVectorW());
225 
226         // Compute limits only if the perigee is close enough from the central body to be in the shadow
227         if (FastMath.abs(gamma * auxiliaryElements.getSma() * (1. - auxiliaryElements.getEcc())) < ae) {
228 
229             // Compute the coefficients of the quartic equation in cos(L) 3.5-(2)
230             final double bet2 = beta * beta;
231             final double h2 = auxiliaryElements.getH() * auxiliaryElements.getH();
232             final double k2 = auxiliaryElements.getK() * auxiliaryElements.getK();
233             final double m  = ae / (auxiliaryElements.getSma() * auxiliaryElements.getB());
234             final double m2 = m * m;
235             final double m4 = m2 * m2;
236             final double bb = alpha * beta + m2 * auxiliaryElements.getH() * auxiliaryElements.getK();
237             final double b2 = bb * bb;
238             final double cc = alpha * alpha - bet2 + m2 * (k2 - h2);
239             final double dd = 1. - bet2 - m2 * (1. + h2);
240             final double[] a = new double[5];
241             a[0] = 4. * b2 + cc * cc;
242             a[1] = 8. * bb * m2 * auxiliaryElements.getH() + 4. * cc * m2 * auxiliaryElements.getK();
243             a[2] = -4. * b2 + 4. * m4 * h2 - 2. * cc * dd + 4. * m4 * k2;
244             a[3] = -8. * bb * m2 * auxiliaryElements.getH() - 4. * dd * m2 * auxiliaryElements.getK();
245             a[4] = -4. * m4 * h2 + dd * dd;
246             // Compute the real roots of the quartic equation 3.5-2
247             final double[] roots = new double[4];
248             final int nbRoots = realQuarticRoots(a, roots);
249             if (nbRoots > 0) {
250                 // Check for consistency
251                 boolean entryFound = false;
252                 boolean exitFound  = false;
253                 // Eliminate spurious roots
254                 for (int i = 0; i < nbRoots; i++) {
255                     final double cosL = roots[i];
256                     final double sL = FastMath.sqrt((1. - cosL) * (1. + cosL));
257                     // Check both angles: L and -L
258                     for (int j = -1; j <= 1; j += 2) {
259                         final double sinL = j * sL;
260                         final double cPhi = alpha * cosL + beta * sinL;
261                         // Is the angle on the shadow side of the central body (eq. 3.5-3) ?
262                         if (cPhi < 0.) {
263                             final double range = 1. + auxiliaryElements.getK() * cosL + auxiliaryElements.getH() * sinL;
264                             final double S  = 1. - m2 * range * range - cPhi * cPhi;
265                             // Is the shadow equation 3.5-1 satisfied ?
266                             if (FastMath.abs(S) < S_ZERO) {
267                                 // Is this the entry or exit angle ?
268                                 final double dSdL = m2 * range * (auxiliaryElements.getK() * sinL - auxiliaryElements.getH() * cosL) + cPhi * (alpha * sinL - beta * cosL);
269                                 if (dSdL > 0.) {
270                                     // Exit from shadow: 3.5-4
271                                     exitFound = true;
272                                     ll[0] = FastMath.atan2(sinL, cosL);
273                                 } else {
274                                     // Entry into shadow: 3.5-5
275                                     entryFound = true;
276                                     ll[1] = FastMath.atan2(sinL, cosL);
277                                 }
278                             }
279                         }
280                     }
281                 }
282                 // Must be one entry and one exit or none
283                 if (!(entryFound == exitFound)) {
284                     // entry or exit found but not both ! In this case, consider there is no eclipse...
285                     ll[0] = -FastMath.PI;
286                     ll[1] = FastMath.PI;
287                 }
288                 // Quadrature between L at exit and L at entry so Lexit must be lower than Lentry
289                 if (ll[0] > ll[1]) {
290                     // Keep the angles between [-2PI, 2PI]
291                     if (ll[1] < 0.) {
292                         ll[1] += 2. * FastMath.PI;
293                     } else {
294                         ll[0] -= 2. * FastMath.PI;
295                     }
296                 }
297             }
298         }
299         return ll;
300     }
301 
302     /** {@inheritDoc} */
303     protected <T extends CalculusFieldElement<T>> T[] getLLimits(final FieldSpacecraftState<T> state,
304                                                              final FieldAuxiliaryElements<T> auxiliaryElements) {
305 
306         final Field<T> field = state.getDate().getField();
307         final T zero = field.getZero();
308         final T one  = field.getOne();
309         final T pi   = one.getPi();
310 
311         // Default bounds without shadow [-PI, PI]
312         final T[] ll = MathArrays.buildArray(field, 2);
313         ll[0] = MathUtils.normalizeAngle(state.getLv(), zero).subtract(pi);
314         ll[1] = MathUtils.normalizeAngle(state.getLv(), zero).add(pi);
315 
316         // Direction cosines of the Sun in the equinoctial frame
317         final FieldVector3D<T> sunDir = sun.getPVCoordinates(state.getDate(), state.getFrame()).getPosition().normalize();
318         final T alpha = sunDir.dotProduct(auxiliaryElements.getVectorF());
319         final T beta  = sunDir.dotProduct(auxiliaryElements.getVectorG());
320         final T gamma = sunDir.dotProduct(auxiliaryElements.getVectorW());
321 
322         // Compute limits only if the perigee is close enough from the central body to be in the shadow
323         if (FastMath.abs(gamma.multiply(auxiliaryElements.getSma()).multiply(auxiliaryElements.getEcc().negate().add(one))).getReal() < ae) {
324 
325             // Compute the coefficients of the quartic equation in cos(L) 3.5-(2)
326             final T bet2 = beta.multiply(beta);
327             final T h2 = auxiliaryElements.getH().multiply(auxiliaryElements.getH());
328             final T k2 = auxiliaryElements.getK().multiply(auxiliaryElements.getK());
329             final T m  = (auxiliaryElements.getSma().multiply(auxiliaryElements.getB())).divide(ae).reciprocal();
330             final T m2 = m.multiply(m);
331             final T m4 = m2.multiply(m2);
332             final T bb = alpha.multiply(beta).add(m2.multiply(auxiliaryElements.getH()).multiply(auxiliaryElements.getK()));
333             final T b2 = bb.multiply(bb);
334             final T cc = alpha.multiply(alpha).subtract(bet2).add(m2.multiply(k2.subtract(h2)));
335             final T dd = bet2.add(m2.multiply(h2.add(1.))).negate().add(one);
336             final T[] a = MathArrays.buildArray(field, 5);
337             a[0] = b2.multiply(4.).add(cc.multiply(cc));
338             a[1] = bb.multiply(8.).multiply(m2).multiply(auxiliaryElements.getH()).add(cc.multiply(4.).multiply(m2).multiply(auxiliaryElements.getK()));
339             a[2] = m4.multiply(h2).multiply(4.).subtract(cc.multiply(dd).multiply(2.)).add(m4.multiply(k2).multiply(4.)).subtract(b2.multiply(4.));
340             a[3] = auxiliaryElements.getH().multiply(m2).multiply(bb).multiply(8.).add(auxiliaryElements.getK().multiply(m2).multiply(dd).multiply(4.)).negate();
341             a[4] = dd.multiply(dd).subtract(m4.multiply(h2).multiply(4.));
342             // Compute the real roots of the quartic equation 3.5-2
343             final T[] roots = MathArrays.buildArray(field, 4);
344             final int nbRoots = realQuarticRoots(a, roots, field);
345             if (nbRoots > 0) {
346                 // Check for consistency
347                 boolean entryFound = false;
348                 boolean exitFound  = false;
349                 // Eliminate spurious roots
350                 for (int i = 0; i < nbRoots; i++) {
351                     final T cosL = roots[i];
352                     final T sL = FastMath.sqrt((cosL.negate().add(one)).multiply(cosL.add(one)));
353                     // Check both angles: L and -L
354                     for (int j = -1; j <= 1; j += 2) {
355                         final T sinL = sL.multiply(j);
356                         final T cPhi = cosL.multiply(alpha).add(sinL.multiply(beta));
357                         // Is the angle on the shadow side of the central body (eq. 3.5-3) ?
358                         if (cPhi.getReal() < 0.) {
359                             final T range = cosL.multiply(auxiliaryElements.getK()).add(sinL.multiply(auxiliaryElements.getH())).add(one);
360                             final T S  = (range.multiply(range).multiply(m2).add(cPhi.multiply(cPhi))).negate().add(1.);
361                             // Is the shadow equation 3.5-1 satisfied ?
362                             if (FastMath.abs(S).getReal() < S_ZERO) {
363                                 // Is this the entry or exit angle ?
364                                 final T dSdL = m2.multiply(range).multiply(auxiliaryElements.getK().multiply(sinL).subtract(auxiliaryElements.getH().multiply(cosL))).add(cPhi.multiply(alpha.multiply(sinL).subtract(beta.multiply(cosL))));
365                                 if (dSdL.getReal() > 0.) {
366                                     // Exit from shadow: 3.5-4
367                                     exitFound = true;
368                                     ll[0] = FastMath.atan2(sinL, cosL);
369                                 } else {
370                                     // Entry into shadow: 3.5-5
371                                     entryFound = true;
372                                     ll[1] = FastMath.atan2(sinL, cosL);
373                                 }
374                             }
375                         }
376                     }
377                 }
378                 // Must be one entry and one exit or none
379                 if (!(entryFound == exitFound)) {
380                     // entry or exit found but not both ! In this case, consider there is no eclipse...
381                     ll[0] = pi.negate();
382                     ll[1] = pi;
383                 }
384                 // Quadrature between L at exit and L at entry so Lexit must be lower than Lentry
385                 if (ll[0].getReal() > ll[1].getReal()) {
386                     // Keep the angles between [-2PI, 2PI]
387                     if (ll[1].getReal() < 0.) {
388                         ll[1] = ll[1].add(pi.multiply(2.0));
389                     } else {
390                         ll[0] = ll[0].subtract(pi.multiply(2.0));
391                     }
392                 }
393             }
394         }
395         return ll;
396     }
397 
398     /** Get the central body equatorial radius.
399      *  @return central body equatorial radius (m)
400      */
401     public double getEquatorialRadius() {
402         return ae;
403     }
404 
405     /**
406      * Compute the real roots of a quartic equation.
407      *
408      * <pre>
409      * a[0] * x⁴ + a[1] * x³ + a[2] * x² + a[3] * x + a[4] = 0.
410      * </pre>
411      *
412      * @param a the 5 coefficients
413      * @param y the real roots
414      * @return the number of real roots
415      */
416     private int realQuarticRoots(final double[] a, final double[] y) {
417         /* Treat the degenerate quartic as cubic */
418         if (Precision.equals(a[0], 0.)) {
419             final double[] aa = new double[a.length - 1];
420             System.arraycopy(a, 1, aa, 0, aa.length);
421             return realCubicRoots(aa, y);
422         }
423 
424         // Transform coefficients
425         final double b  = a[1] / a[0];
426         final double c  = a[2] / a[0];
427         final double d  = a[3] / a[0];
428         final double e  = a[4] / a[0];
429         final double bh = b * 0.5;
430 
431         // Solve resolvant cubic
432         final double[] z3 = new double[3];
433         final int i3 = realCubicRoots(new double[] {1.0, -c, b * d - 4. * e, e * (4. * c - b * b) - d * d}, z3);
434         if (i3 == 0) {
435             return 0;
436         }
437 
438         // Largest real root of resolvant cubic
439         final double z = z3[0];
440 
441         // Compute auxiliary quantities
442         final double zh = z * 0.5;
443         final double p  = FastMath.max(z + bh * bh - c, 0.);
444         final double q  = FastMath.max(zh * zh - e, 0.);
445         final double r  = bh * z - d;
446         final double pp = FastMath.sqrt(p);
447         final double qq = FastMath.copySign(FastMath.sqrt(q), r);
448 
449         // Solve quadratic factors of quartic equation
450         final double[] y1 = new double[2];
451         final int n1 = realQuadraticRoots(new double[] {1.0, bh - pp, zh - qq}, y1);
452         final double[] y2 = new double[2];
453         final int n2 = realQuadraticRoots(new double[] {1.0, bh + pp, zh + qq}, y2);
454 
455         if (n1 == 2) {
456             if (n2 == 2) {
457                 y[0] = y1[0];
458                 y[1] = y1[1];
459                 y[2] = y2[0];
460                 y[3] = y2[1];
461                 return 4;
462             } else {
463                 y[0] = y1[0];
464                 y[1] = y1[1];
465                 return 2;
466             }
467         } else {
468             if (n2 == 2) {
469                 y[0] = y2[0];
470                 y[1] = y2[1];
471                 return 2;
472             } else {
473                 return 0;
474             }
475         }
476     }
477 
478     /**
479      * Compute the real roots of a quartic equation.
480      *
481      * <pre>
482      * a[0] * x⁴ + a[1] * x³ + a[2] * x² + a[3] * x + a[4] = 0.
483      * </pre>
484      * @param <T> the type of the field elements
485      *
486      * @param a the 5 coefficients
487      * @param y the real roots
488      * @param field field of elements
489      * @return the number of real roots
490      */
491     private <T extends CalculusFieldElement<T>> int realQuarticRoots(final T[] a, final T[] y,
492                                                                  final Field<T> field) {
493 
494         final T zero = field.getZero();
495 
496         // Treat the degenerate quartic as cubic
497         if (Precision.equals(a[0].getReal(), 0.)) {
498             final T[] aa = MathArrays.buildArray(field, a.length - 1);
499             System.arraycopy(a, 1, aa, 0, aa.length);
500             return realCubicRoots(aa, y, field);
501         }
502 
503         // Transform coefficients
504         final T b  = a[1].divide(a[0]);
505         final T c  = a[2].divide(a[0]);
506         final T d  = a[3].divide(a[0]);
507         final T e  = a[4].divide(a[0]);
508         final T bh = b.multiply(0.5);
509 
510         // Solve resolvant cubic
511         final T[] z3 = MathArrays.buildArray(field, 3);
512         final T[] i = MathArrays.buildArray(field, 4);
513         i[0] = zero.add(1.0);
514         i[1] = c.negate();
515         i[2] = b.multiply(d).subtract(e.multiply(4.0));
516         i[3] = e.multiply(c.multiply(4.).subtract(b.multiply(b))).subtract(d.multiply(d));
517         final int i3 = realCubicRoots(i, z3, field);
518         if (i3 == 0) {
519             return 0;
520         }
521 
522         // Largest real root of resolvant cubic
523         final T z = z3[0];
524 
525         // Compute auxiliary quantities
526         final T zh = z.multiply(0.5);
527         final T p  = FastMath.max(z.add(bh.multiply(bh)).subtract(c), zero);
528         final T q  = FastMath.max(zh.multiply(zh).subtract(e), zero);
529         final T r  = bh.multiply(z).subtract(d);
530         final T pp = FastMath.sqrt(p);
531         final T qq = FastMath.copySign(FastMath.sqrt(q), r);
532 
533         // Solve quadratic factors of quartic equation
534         final T[] y1 = MathArrays.buildArray(field, 2);
535         final T[] n = MathArrays.buildArray(field, 3);
536         n[0] = zero.add(1.0);
537         n[1] = bh.subtract(pp);
538         n[2] = zh.subtract(qq);
539         final int n1 = realQuadraticRoots(n, y1);
540         final T[] y2 = MathArrays.buildArray(field, 2);
541         final T[] nn = MathArrays.buildArray(field, 3);
542         nn[0] = zero.add(1.0);
543         nn[1] = bh.add(pp);
544         nn[2] = zh.add(qq);
545         final int n2 = realQuadraticRoots(nn, y2);
546 
547         if (n1 == 2) {
548             if (n2 == 2) {
549                 y[0] = y1[0];
550                 y[1] = y1[1];
551                 y[2] = y2[0];
552                 y[3] = y2[1];
553                 return 4;
554             } else {
555                 y[0] = y1[0];
556                 y[1] = y1[1];
557                 return 2;
558             }
559         } else {
560             if (n2 == 2) {
561                 y[0] = y2[0];
562                 y[1] = y2[1];
563                 return 2;
564             } else {
565                 return 0;
566             }
567         }
568     }
569 
570     /**
571      * Compute the real roots of a cubic equation.
572      *
573      * <pre>
574      * a[0] * x³ + a[1] * x² + a[2] * x + a[3] = 0.
575      * </pre>
576      *
577      * @param a the 4 coefficients
578      * @param y the real roots sorted in descending order
579      * @return the number of real roots
580      */
581     private int realCubicRoots(final double[] a, final double[] y) {
582         if (Precision.equals(a[0], 0.)) {
583             // Treat the degenerate cubic as quadratic
584             final double[] aa = new double[a.length - 1];
585             System.arraycopy(a, 1, aa, 0, aa.length);
586             return realQuadraticRoots(aa, y);
587         }
588 
589         // Transform coefficients
590         final double b  = -a[1] / (3. * a[0]);
591         final double c  =  a[2] / a[0];
592         final double d  =  a[3] / a[0];
593         final double b2 =  b * b;
594         final double p  =  b2 - c / 3.;
595         final double q  =  b * (b2 - c * 0.5) - d * 0.5;
596 
597         // Compute discriminant
598         final double disc = p * p * p - q * q;
599 
600         if (disc < 0.) {
601             // One real root
602             final double alpha  = q + FastMath.copySign(FastMath.sqrt(-disc), q);
603             final double cbrtAl = FastMath.cbrt(alpha);
604             final double cbrtBe = p / cbrtAl;
605 
606             if (p < 0.) {
607                 y[0] = b + 2. * q / (cbrtAl * cbrtAl + cbrtBe * cbrtBe - p);
608             } else if (p > 0.) {
609                 y[0] = b + cbrtAl + cbrtBe;
610             } else {
611                 y[0] = b + cbrtAl;
612             }
613 
614             return 1;
615 
616         } else if (disc > 0.) {
617             // Three distinct real roots
618             final double phi = FastMath.atan2(FastMath.sqrt(disc), q) / 3.;
619             final double sqP = 2.0 * FastMath.sqrt(p);
620 
621             y[0] = b + sqP * FastMath.cos(phi);
622             y[1] = b - sqP * FastMath.cos(FastMath.PI / 3. + phi);
623             y[2] = b - sqP * FastMath.cos(FastMath.PI / 3. - phi);
624 
625             return 3;
626 
627         } else {
628             // One distinct and two equals real roots
629             final double cbrtQ = FastMath.cbrt(q);
630             final double root1 = b + 2. * cbrtQ;
631             final double root2 = b - cbrtQ;
632             if (q < 0.) {
633                 y[0] = root2;
634                 y[1] = root2;
635                 y[2] = root1;
636             } else {
637                 y[0] = root1;
638                 y[1] = root2;
639                 y[2] = root2;
640             }
641 
642             return 3;
643         }
644     }
645 
646     /**
647      * Compute the real roots of a cubic equation.
648      *
649      * <pre>
650      * a[0] * x³ + a[1] * x² + a[2] * x + a[3] = 0.
651      * </pre>
652      *
653      * @param <T> the type of the field elements
654      * @param a the 4 coefficients
655      * @param y the real roots sorted in descending order
656      * @param field field of elements
657      * @return the number of real roots
658      */
659     private <T extends CalculusFieldElement<T>> int realCubicRoots(final T[] a, final T[] y,
660                                                                final Field<T> field) {
661 
662         if (Precision.equals(a[0].getReal(), 0.)) {
663             // Treat the degenerate cubic as quadratic
664             final T[] aa = MathArrays.buildArray(field, a.length - 1);
665             System.arraycopy(a, 1, aa, 0, aa.length);
666             return realQuadraticRoots(aa, y);
667         }
668 
669         // Transform coefficients
670         final T b  =  a[1].divide(a[0].multiply(3.)).negate();
671         final T c  =  a[2].divide(a[0]);
672         final T d  =  a[3].divide(a[0]);
673         final T b2 =  b.multiply(b);
674         final T p  =  b2.subtract(c.divide(3.));
675         final T q  =  b.multiply(b2.subtract(c.multiply(0.5))).subtract(d.multiply(0.5));
676 
677         // Compute discriminant
678         final T disc = p.multiply(p).multiply(p).subtract(q.multiply(q));
679 
680         if (disc.getReal() < 0.) {
681             // One real root
682             final T alpha  = FastMath.copySign(FastMath.sqrt(disc.negate()), q).add(q);
683             final T cbrtAl = FastMath.cbrt(alpha);
684             final T cbrtBe = p.divide(cbrtAl);
685 
686             if (p .getReal() < 0.) {
687                 y[0] = q.divide(cbrtAl.multiply(cbrtAl).add(cbrtBe.multiply(cbrtBe)).subtract(p)).multiply(2.).add(b);
688             } else if (p.getReal() > 0.) {
689                 y[0] = b.add(cbrtAl).add(cbrtBe);
690             } else {
691                 y[0] = b.add(cbrtAl);
692             }
693 
694             return 1;
695 
696         } else if (disc.getReal() > 0.) {
697             // Three distinct real roots
698             final T phi = FastMath.atan2(FastMath.sqrt(disc), q).divide(3.);
699             final T sqP = FastMath.sqrt(p).multiply(2.);
700 
701             y[0] = b.add(sqP.multiply(FastMath.cos(phi)));
702             y[1] = b.subtract(sqP.multiply(FastMath.cos(phi.add(b.getPi().divide(3.)))));
703             y[2] = b.subtract(sqP.multiply(FastMath.cos(phi.negate().add(b.getPi().divide(3.)))));
704 
705             return 3;
706 
707         } else {
708             // One distinct and two equals real roots
709             final T cbrtQ = FastMath.cbrt(q);
710             final T root1 = b.add(cbrtQ.multiply(2.0));
711             final T root2 = b.subtract(cbrtQ);
712             if (q.getReal() < 0.) {
713                 y[0] = root2;
714                 y[1] = root2;
715                 y[2] = root1;
716             } else {
717                 y[0] = root1;
718                 y[1] = root2;
719                 y[2] = root2;
720             }
721 
722             return 3;
723         }
724     }
725 
726     /**
727      * Compute the real roots of a quadratic equation.
728      *
729      * <pre>
730      * a[0] * x² + a[1] * x + a[2] = 0.
731      * </pre>
732      *
733      * @param a the 3 coefficients
734      * @param y the real roots sorted in descending order
735      * @return the number of real roots
736      */
737     private int realQuadraticRoots(final double[] a, final double[] y) {
738         if (Precision.equals(a[0], 0.)) {
739             // Degenerate quadratic
740             if (Precision.equals(a[1], 0.)) {
741                 // Degenerate linear equation: no real roots
742                 return 0;
743             }
744             // Linear equation: one real root
745             y[0] = -a[2] / a[1];
746             return 1;
747         }
748 
749         // Transform coefficients
750         final double b = -0.5 * a[1] / a[0];
751         final double c =  a[2] / a[0];
752 
753         // Compute discriminant
754         final double d =  b * b - c;
755 
756         if (d < 0.) {
757             // No real roots
758             return 0;
759         } else if (d > 0.) {
760             // Two distinct real roots
761             final double y0 = b + FastMath.copySign(FastMath.sqrt(d), b);
762             final double y1 = c / y0;
763             y[0] = FastMath.max(y0, y1);
764             y[1] = FastMath.min(y0, y1);
765             return 2;
766         } else {
767             // Discriminant is zero: two equal real roots
768             y[0] = b;
769             y[1] = b;
770             return 2;
771         }
772     }
773 
774     /**
775      * Compute the real roots of a quadratic equation.
776      *
777      * <pre>
778      * a[0] * x² + a[1] * x + a[2] = 0.
779      * </pre>
780      *
781      * @param <T> the type of the field elements
782      * @param a the 3 coefficients
783      * @param y the real roots sorted in descending order
784      * @return the number of real roots
785      */
786     private <T extends CalculusFieldElement<T>> int realQuadraticRoots(final T[] a, final T[] y) {
787 
788         if (Precision.equals(a[0].getReal(), 0.)) {
789             // Degenerate quadratic
790             if (Precision.equals(a[1].getReal(), 0.)) {
791                 // Degenerate linear equation: no real roots
792                 return 0;
793             }
794             // Linear equation: one real root
795             y[0] = a[2].divide(a[1]).negate();
796             return 1;
797         }
798 
799         // Transform coefficients
800         final T b = a[1].divide(a[0]).multiply(0.5).negate();
801         final T c =  a[2].divide(a[0]);
802 
803         // Compute discriminant
804         final T d =  b.multiply(b).subtract(c);
805 
806         if (d.getReal() < 0.) {
807             // No real roots
808             return 0;
809         } else if (d.getReal() > 0.) {
810             // Two distinct real roots
811             final T y0 = b.add(FastMath.copySign(FastMath.sqrt(d), b));
812             final T y1 = c.divide(y0);
813             y[0] = FastMath.max(y0, y1);
814             y[1] = FastMath.min(y0, y1);
815             return 2;
816         } else {
817             // Discriminant is zero: two equal real roots
818             y[0] = b;
819             y[1] = b;
820             return 2;
821         }
822     }
823 
824 }