AbstractGaussianContribution.java

  1. /* Copyright 2002-2013 CS Systèmes d'Information
  2.  * Licensed to CS Systèmes d'Information (CS) under one or more
  3.  * contributor license agreements.  See the NOTICE file distributed with
  4.  * this work for additional information regarding copyright ownership.
  5.  * CS licenses this file to You under the Apache License, Version 2.0
  6.  * (the "License"); you may not use this file except in compliance with
  7.  * the License.  You may obtain a copy of the License at
  8.  *
  9.  *   http://www.apache.org/licenses/LICENSE-2.0
  10.  *
  11.  * Unless required by applicable law or agreed to in writing, software
  12.  * distributed under the License is distributed on an "AS IS" BASIS,
  13.  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
  14.  * See the License for the specific language governing permissions and
  15.  * limitations under the License.
  16.  */
  17. package org.orekit.propagation.semianalytical.dsst.forces;

  18. import org.apache.commons.math3.analysis.UnivariateVectorFunction;
  19. import org.apache.commons.math3.geometry.euclidean.threed.Vector3D;
  20. import org.apache.commons.math3.util.FastMath;
  21. import org.orekit.errors.OrekitException;
  22. import org.orekit.errors.OrekitExceptionWrapper;
  23. import org.orekit.propagation.SpacecraftState;
  24. import org.orekit.propagation.semianalytical.dsst.utilities.AuxiliaryElements;

  25. /** Common handling of {@link DSSTForceModel} methods for Gaussian contributions to DSST propagation.
  26.  * <p>
  27.  * This abstract class allows to provide easily a subset of {@link DSSTForceModel} methods
  28.  * for specific Gaussian contributions (i.e. atmospheric drag and solar radiation pressure).
  29.  * </p><p>
  30.  * Gaussian contributions can be expressed as: da<sub>i</sub>/dt = &delta;a<sub>i</sub>/&delta;v . q<br>
  31.  * where:
  32.  * <ul>
  33.  * <li>a<sub>i</sub> are the six equinoctial elements</li>
  34.  * <li>v is the velocity vector</li>
  35.  * <li>q is the perturbing acceleration due to the considered force</li>
  36.  * </ul>
  37.  * The averaging process and other considerations lead to integrate this contribution
  38.  * over the true longitude L possibly taking into account some limits.
  39.  * </p><p>
  40.  * Only two methods must be implemented by derived classes:
  41.  * {@link #getAcceleration(SpacecraftState, Vector3D, Vector3D)} and
  42.  * {@link #getLLimits(SpacecraftState)}.
  43.  * </p>
  44.  * @author Pascal Parraud
  45.  */
  46. public abstract class AbstractGaussianContribution implements DSSTForceModel {

  47.     /** Available orders for Gauss quadrature. */
  48.     private static final int[] GAUSS_ORDER = {12, 16, 20, 24, 32, 40, 48};

  49.     /** Max rank in Gauss quadrature orders array. */
  50.     private static final int MAX_ORDER_RANK = GAUSS_ORDER.length - 1;

  51.     // CHECKSTYLE: stop VisibilityModifierCheck

  52.     /** Retrograde factor. */
  53.     protected double I;

  54.     /** a. */
  55.     protected double a;
  56.     /** e<sub>x</sub>. */
  57.     protected double k;
  58.     /** e<sub>y</sub>. */
  59.     protected double h;
  60.     /** h<sub>x</sub>. */
  61.     protected double q;
  62.     /** h<sub>y</sub>. */
  63.     protected double p;

  64.     /** Eccentricity. */
  65.     protected double ecc;

  66.     /** Kepler mean motion: n = sqrt(&mu; / a<sup>3</sup>). */
  67.     protected double n;

  68.     /** Equinoctial frame f vector. */
  69.     protected Vector3D f;
  70.     /** Equinoctial frame g vector. */
  71.     protected Vector3D g;
  72.     /** Equinoctial frame w vector. */
  73.     protected Vector3D w;

  74.     /** A = sqrt(&mu; * a). */
  75.     protected double A;
  76.     /** B = sqrt(1 - h<sup>2</sup> - k<sup>2</sup>). */
  77.     protected double B;
  78.     /** C = 1 + p<sup>2</sup> + q<sup>2</sup>. */
  79.     protected double C;

  80.     /** 2 / (n<sup>2</sup> * a) . */
  81.     protected double ton2a;
  82.     /** 1 / A .*/
  83.     protected double ooA;
  84.     /** 1 / (A * B) .*/
  85.     protected double ooAB;
  86.     /** C / (2 * A * B) .*/
  87.     protected double co2AB;
  88.     /** 1 / (1 + B) .*/
  89.     protected double ooBpo;
  90.     /** 1 / &mu; .*/
  91.     protected double ooMu;

  92.     // CHECKSTYLE: resume VisibilityModifierCheck

  93.     /** Gauss integrator. */
  94.     private final double threshold;

  95.     /** Gauss integrator. */
  96.     private GaussQuadrature integrator;

  97.     /** Flag for Gauss order computation. */
  98.     private boolean isDirty;

  99.     /** Build a new instance.
  100.      *
  101.      *  @param threshold tolerance for the choice of the Gauss quadrature order
  102.      */
  103.     protected AbstractGaussianContribution(final double threshold) {
  104.         this.threshold  = threshold;
  105.         this.integrator = new GaussQuadrature(GAUSS_ORDER[MAX_ORDER_RANK]);
  106.         this.isDirty    = true;
  107.     }

  108.     /** {@inheritDoc} */
  109.     public void initialize(final AuxiliaryElements aux)
  110.         throws OrekitException {
  111.         // Nothing to do for gaussian contributions at the beginning of the propagation.
  112.     }

  113.     /** {@inheritDoc} */
  114.     public void initializeStep(final AuxiliaryElements aux)
  115.         throws OrekitException {

  116.         // Equinoctial elements
  117.         a  = aux.getSma();
  118.         k  = aux.getK();
  119.         h  = aux.getH();
  120.         q  = aux.getQ();
  121.         p  = aux.getP();

  122.         // Retrograde factor
  123.         I = aux.getRetrogradeFactor();

  124.         // Eccentricity
  125.         ecc = aux.getEcc();

  126.         // Equinoctial coefficients
  127.         A = aux.getA();
  128.         B = aux.getB();
  129.         C = aux.getC();

  130.         // Equinoctial frame vectors
  131.         f = aux.getVectorF();
  132.         g = aux.getVectorG();
  133.         w = aux.getVectorW();

  134.         // Kepler mean motion
  135.         n = A / (a * a);

  136.         // 1 / A
  137.         ooA = 1. / A;
  138.         // 1 / AB
  139.         ooAB = ooA / B;
  140.         // C / 2AB
  141.         co2AB = C * ooAB / 2.;
  142.         // 1 / (1 + B)
  143.         ooBpo = 1. / (1. + B);
  144.         // 2 / (n² * a)
  145.         ton2a = 2. / (n * n * a);
  146.         // 1 / mu
  147.         ooMu  = 1. / aux.getMu();
  148.     }

  149.     /** {@inheritDoc} */
  150.     public double[] getMeanElementRate(final SpacecraftState state) throws OrekitException {

  151.         double[] meanElementRate = new double[6];
  152.         // Computes the limits for the integral
  153.         final double[] ll = getLLimits(state);
  154.         // Computes integrated mean element rates if Llow < Lhigh
  155.         if (ll[0] < ll[1]) {
  156.             meanElementRate = getMeanElementRate(state, integrator, ll[0], ll[1]);
  157.             if (isDirty) {
  158.                 boolean next = true;
  159.                 for (int i = 0; i < MAX_ORDER_RANK && next; i++) {
  160.                     final double[] meanRates = getMeanElementRate(state, new GaussQuadrature(GAUSS_ORDER[i]), ll[0], ll[1]);
  161.                     if (getRatesDiff(meanElementRate, meanRates) < threshold) {
  162.                         integrator = new GaussQuadrature(GAUSS_ORDER[i]);
  163.                         next = false;
  164.                     }
  165.                 }
  166.                 isDirty = false;
  167.             }
  168.         }
  169.         return meanElementRate;
  170.     }

  171.     /** Compute the acceleration due to the non conservative perturbing force.
  172.      *
  173.      *  @param state current state information: date, kinematics, attitude
  174.      *  @param position spacecraft position
  175.      *  @param velocity spacecraft velocity
  176.      *  @return the perturbing acceleration
  177.      *  @exception OrekitException if some specific error occurs
  178.      */
  179.     protected abstract Vector3D getAcceleration(final SpacecraftState state,
  180.                                                 final Vector3D position,
  181.                                                 final Vector3D velocity) throws OrekitException;

  182.     /** Compute the limits in L, the true longitude, for integration.
  183.      *
  184.      *  @param  state current state information: date, kinematics, attitude
  185.      *  @return the integration limits in L
  186.      *  @exception OrekitException if some specific error occurs
  187.      */
  188.     protected abstract double[] getLLimits(final SpacecraftState state) throws OrekitException;

  189.     /** Computes the mean equinoctial elements rates da<sub>i</sub> / dt.
  190.      *
  191.      *  @param state current state
  192.      *  @param gauss Gauss quadrature
  193.      *  @param low lower bound of the integral interval
  194.      *  @param high upper bound of the integral interval
  195.      *  @return the mean element rates
  196.      *  @throws OrekitException if some specific error occurs
  197.      */
  198.     private double[] getMeanElementRate(final SpacecraftState state,
  199.                                         final GaussQuadrature gauss,
  200.                                         final double low,
  201.                                         final double high) throws OrekitException {
  202.         final double[] meanElementRate = gauss.integrate(new IntegrableFunction(state), low, high);
  203.         // Constant multiplier for integral
  204.         final double coef = 1. / (2. * FastMath.PI * B);
  205.         // Corrects mean element rates
  206.         for (int i = 0; i < 6; i++) {
  207.             meanElementRate[i] *= coef;
  208.         }
  209.         return meanElementRate;
  210.     }

  211.     /** Estimates the weighted magnitude of the difference between 2 sets of equinoctial elements rates.
  212.      *
  213.      *  @param meanRef reference rates
  214.      *  @param meanCur current rates
  215.      *  @return estimated magnitude of weighted differences
  216.      */
  217.     private double getRatesDiff(final double[] meanRef, final double[] meanCur) {
  218.         double maxDiff = FastMath.abs(meanRef[0] - meanCur[0]) / a;
  219.         // Corrects mean element rates
  220.         for (int i = 1; i < meanRef.length; i++) {
  221.             final double diff = FastMath.abs(meanRef[i] - meanCur[i]);
  222.             if (maxDiff < diff) maxDiff = diff;
  223.         }
  224.         return maxDiff;
  225.     }

  226.     /** Internal class for numerical quadrature. */
  227.     private class IntegrableFunction implements UnivariateVectorFunction {

  228.         /** Current state. */
  229.         private final SpacecraftState state;

  230.         /** Build a new instance.
  231.          *  @param  state current state information: date, kinematics, attitude
  232.          */
  233.         public IntegrableFunction(final SpacecraftState state) {
  234.             this.state = state;
  235.         }

  236.         /** {@inheritDoc} */
  237.         public double[] value(final double x) {
  238.             final double cosL = FastMath.cos(x);
  239.             final double sinL = FastMath.sin(x);
  240.             final double roa  = B * B / (1. + h * sinL + k * cosL);
  241.             final double roa2 = roa * roa;
  242.             final double r    = a * roa;
  243.             final double X    = r * cosL;
  244.             final double Y    = r * sinL;
  245.             final double naob = n * a / B;
  246.             final double Xdot = -naob * (h + sinL);
  247.             final double Ydot =  naob * (k + cosL);
  248.             final Vector3D pos = new Vector3D(X, f, Y, g);
  249.             final Vector3D vel = new Vector3D(Xdot, f, Ydot, g);
  250.             // Compute acceleration
  251.             Vector3D acc = Vector3D.ZERO;
  252.             try {
  253.                 acc = getAcceleration(state, pos, vel);
  254.             } catch (OrekitException oe) {
  255.                 throw new OrekitExceptionWrapper(oe);
  256.             }
  257.             // Compute mean elements rates
  258.             final double[] val = new double[6];
  259.             // da/dt
  260.             val[0] = roa2 * getAoV(vel).dotProduct(acc);
  261.             // dex/dt
  262.             val[1] = roa2 * getKoV(X, Y, Xdot, Ydot).dotProduct(acc);
  263.             // dey/dt
  264.             val[2] = roa2 * getHoV(X, Y, Xdot, Ydot).dotProduct(acc);
  265.             // dhx/dt
  266.             val[3] = roa2 * getQoV(X).dotProduct(acc);
  267.             // dhy/dt
  268.             val[4] = roa2 * getPoV(Y).dotProduct(acc);
  269.             // d&lambda;/dt
  270.             val[5] = roa2 * getLoV(X, Y, Xdot, Ydot).dotProduct(acc);

  271.             return val;
  272.         }

  273.         /** Compute &delta;a/&delta;v.
  274.          *  @param vel satellite velocity
  275.          *  @return &delta;a/&delta;v
  276.          */
  277.         private Vector3D getAoV(final Vector3D vel) {
  278.             return new Vector3D(ton2a, vel);
  279.         }

  280.         /** Compute &delta;h/&delta;v.
  281.          *  @param X satellite position component along f, equinoctial reference frame 1st vector
  282.          *  @param Y satellite position component along g, equinoctial reference frame 2nd vector
  283.          *  @param Xdot satellite velocity component along f, equinoctial reference frame 1st vector
  284.          *  @param Ydot satellite velocity component along g, equinoctial reference frame 2nd vector
  285.          *  @return &delta;h/&delta;v
  286.          */
  287.         private Vector3D getHoV(final double X, final double Y, final double Xdot, final double Ydot) {
  288.             final double kf = (2. * Xdot * Y - X * Ydot) * ooMu;
  289.             final double kg = X * Xdot * ooMu;
  290.             final double kw = k * (I * q * Y - p * X) * ooAB;
  291.             return new Vector3D(kf, f, -kg, g, kw, w);
  292.         }

  293.         /** Compute &delta;k/&delta;v.
  294.          *  @param X satellite position component along f, equinoctial reference frame 1st vector
  295.          *  @param Y satellite position component along g, equinoctial reference frame 2nd vector
  296.          *  @param Xdot satellite velocity component along f, equinoctial reference frame 1st vector
  297.          *  @param Ydot satellite velocity component along g, equinoctial reference frame 2nd vector
  298.          *  @return &delta;k/&delta;v
  299.          */
  300.         private Vector3D getKoV(final double X, final double Y, final double Xdot, final double Ydot) {
  301.             final double kf = Y * Ydot * ooMu;
  302.             final double kg = (2. * X * Ydot - Xdot * Y) * ooMu;
  303.             final double kw = h * (I * q * Y - p * X) * ooAB;
  304.             return new Vector3D(-kf, f, kg, g, -kw, w);
  305.         }

  306.         /** Compute &delta;p/&delta;v.
  307.          *  @param Y satellite position component along g, equinoctial reference frame 2nd vector
  308.          *  @return &delta;p/&delta;v
  309.          */
  310.         private Vector3D getPoV(final double Y) {
  311.             return new Vector3D(co2AB * Y, w);
  312.         }

  313.         /** Compute &delta;q/&delta;v.
  314.          *  @param X satellite position component along f, equinoctial reference frame 1st vector
  315.          *  @return &delta;q/&delta;v
  316.          */
  317.         private Vector3D getQoV(final double X) {
  318.             return new Vector3D(I * co2AB * X, w);
  319.         }

  320.         /** Compute &delta;&lambda;/&delta;v.
  321.          *  @param X satellite position component along f, equinoctial reference frame 1st vector
  322.          *  @param Y satellite position component along g, equinoctial reference frame 2nd vector
  323.          *  @param Xdot satellite velocity component along f, equinoctial reference frame 1st vector
  324.          *  @param Ydot satellite velocity component along g, equinoctial reference frame 2nd vector
  325.          *  @return &delta;&lambda;/&delta;v
  326.          */
  327.         private Vector3D getLoV(final double X, final double Y, final double Xdot, final double Ydot) {
  328.             final Vector3D pos = new Vector3D(X, f, Y, g);
  329.             final Vector3D v2  = new Vector3D(k, getHoV(X, Y, Xdot, Ydot), -h, getKoV(X, Y, Xdot, Ydot));
  330.             return new Vector3D(-2. * ooA, pos, ooBpo, v2, (I * q * Y - p * X) * ooA, w);
  331.         }

  332.     }

  333.     /** Class used to {@link #integrate(UnivariateVectorFunction, double, double) integrate}
  334.      *  a {@link org.apache.commons.math3.analysis.UnivariateVectorFunction function}
  335.      *  of the orbital elements using the Gaussian quadrature rule to get the acceleration.
  336.      */
  337.     private static class GaussQuadrature {

  338.         // CHECKSTYLE: stop NoWhitespaceAfter

  339.         // Points and weights for the available quadrature orders

  340.         /** Points for quadrature of order 12. */
  341.         private static final double[] P_12 = {
  342.             -0.98156063424671910000,
  343.             -0.90411725637047490000,
  344.             -0.76990267419430470000,
  345.             -0.58731795428661740000,
  346.             -0.36783149899818024000,
  347.             -0.12523340851146890000,
  348.             0.12523340851146890000,
  349.             0.36783149899818024000,
  350.             0.58731795428661740000,
  351.             0.76990267419430470000,
  352.             0.90411725637047490000,
  353.             0.98156063424671910000
  354.         };

  355.         /** Weights for quadrature of order 12. */
  356.         private static final double[] W_12 = {
  357.             0.04717533638651220000,
  358.             0.10693932599531830000,
  359.             0.16007832854334633000,
  360.             0.20316742672306584000,
  361.             0.23349253653835478000,
  362.             0.24914704581340286000,
  363.             0.24914704581340286000,
  364.             0.23349253653835478000,
  365.             0.20316742672306584000,
  366.             0.16007832854334633000,
  367.             0.10693932599531830000,
  368.             0.04717533638651220000
  369.         };

  370.         /** Points for quadrature of order 16. */
  371.         private static final double[] P_16 = {
  372.             -0.98940093499164990000,
  373.             -0.94457502307323260000,
  374.             -0.86563120238783160000,
  375.             -0.75540440835500310000,
  376.             -0.61787624440264380000,
  377.             -0.45801677765722737000,
  378.             -0.28160355077925890000,
  379.             -0.09501250983763745000,
  380.             0.09501250983763745000,
  381.             0.28160355077925890000,
  382.             0.45801677765722737000,
  383.             0.61787624440264380000,
  384.             0.75540440835500310000,
  385.             0.86563120238783160000,
  386.             0.94457502307323260000,
  387.             0.98940093499164990000
  388.         };

  389.         /** Weights for quadrature of order 16. */
  390.         private static final double[] W_16 = {
  391.             0.02715245941175405800,
  392.             0.06225352393864777000,
  393.             0.09515851168249283000,
  394.             0.12462897125553388000,
  395.             0.14959598881657685000,
  396.             0.16915651939500256000,
  397.             0.18260341504492360000,
  398.             0.18945061045506847000,
  399.             0.18945061045506847000,
  400.             0.18260341504492360000,
  401.             0.16915651939500256000,
  402.             0.14959598881657685000,
  403.             0.12462897125553388000,
  404.             0.09515851168249283000,
  405.             0.06225352393864777000,
  406.             0.02715245941175405800
  407.         };

  408.         /** Points for quadrature of order 20. */
  409.         private static final double[] P_20 = {
  410.             -0.99312859918509490000,
  411.             -0.96397192727791390000,
  412.             -0.91223442825132600000,
  413.             -0.83911697182221890000,
  414.             -0.74633190646015080000,
  415.             -0.63605368072651510000,
  416.             -0.51086700195082700000,
  417.             -0.37370608871541955000,
  418.             -0.22778585114164507000,
  419.             -0.07652652113349734000,
  420.             0.07652652113349734000,
  421.             0.22778585114164507000,
  422.             0.37370608871541955000,
  423.             0.51086700195082700000,
  424.             0.63605368072651510000,
  425.             0.74633190646015080000,
  426.             0.83911697182221890000,
  427.             0.91223442825132600000,
  428.             0.96397192727791390000,
  429.             0.99312859918509490000
  430.         };

  431.         /** Weights for quadrature of order 20. */
  432.         private static final double[] W_20 = {
  433.             0.01761400713915226400,
  434.             0.04060142980038684000,
  435.             0.06267204833410904000,
  436.             0.08327674157670477000,
  437.             0.10193011981724048000,
  438.             0.11819453196151844000,
  439.             0.13168863844917678000,
  440.             0.14209610931838212000,
  441.             0.14917298647260380000,
  442.             0.15275338713072600000,
  443.             0.15275338713072600000,
  444.             0.14917298647260380000,
  445.             0.14209610931838212000,
  446.             0.13168863844917678000,
  447.             0.11819453196151844000,
  448.             0.10193011981724048000,
  449.             0.08327674157670477000,
  450.             0.06267204833410904000,
  451.             0.04060142980038684000,
  452.             0.01761400713915226400
  453.         };

  454.         /** Points for quadrature of order 24. */
  455.         private static final double[] P_24 = {
  456.             -0.99518721999702130000,
  457.             -0.97472855597130950000,
  458.             -0.93827455200273270000,
  459.             -0.88641552700440100000,
  460.             -0.82000198597390300000,
  461.             -0.74012419157855440000,
  462.             -0.64809365193697550000,
  463.             -0.54542147138883950000,
  464.             -0.43379350762604520000,
  465.             -0.31504267969616340000,
  466.             -0.19111886747361634000,
  467.             -0.06405689286260563000,
  468.             0.06405689286260563000,
  469.             0.19111886747361634000,
  470.             0.31504267969616340000,
  471.             0.43379350762604520000,
  472.             0.54542147138883950000,
  473.             0.64809365193697550000,
  474.             0.74012419157855440000,
  475.             0.82000198597390300000,
  476.             0.88641552700440100000,
  477.             0.93827455200273270000,
  478.             0.97472855597130950000,
  479.             0.99518721999702130000
  480.         };

  481.         /** Weights for quadrature of order 24. */
  482.         private static final double[] W_24 = {
  483.             0.01234122979998733500,
  484.             0.02853138862893380600,
  485.             0.04427743881741981000,
  486.             0.05929858491543691500,
  487.             0.07334648141108027000,
  488.             0.08619016153195320000,
  489.             0.09761865210411391000,
  490.             0.10744427011596558000,
  491.             0.11550566805372553000,
  492.             0.12167047292780335000,
  493.             0.12583745634682825000,
  494.             0.12793819534675221000,
  495.             0.12793819534675221000,
  496.             0.12583745634682825000,
  497.             0.12167047292780335000,
  498.             0.11550566805372553000,
  499.             0.10744427011596558000,
  500.             0.09761865210411391000,
  501.             0.08619016153195320000,
  502.             0.07334648141108027000,
  503.             0.05929858491543691500,
  504.             0.04427743881741981000,
  505.             0.02853138862893380600,
  506.             0.01234122979998733500
  507.         };

  508.         /** Points for quadrature of order 32. */
  509.         private static final double[] P_32 = {
  510.             -0.99726386184948160000,
  511.             -0.98561151154526840000,
  512.             -0.96476225558750640000,
  513.             -0.93490607593773970000,
  514.             -0.89632115576605220000,
  515.             -0.84936761373256990000,
  516.             -0.79448379596794250000,
  517.             -0.73218211874028970000,
  518.             -0.66304426693021520000,
  519.             -0.58771575724076230000,
  520.             -0.50689990893222950000,
  521.             -0.42135127613063540000,
  522.             -0.33186860228212767000,
  523.             -0.23928736225213710000,
  524.             -0.14447196158279646000,
  525.             -0.04830766568773831000,
  526.             0.04830766568773831000,
  527.             0.14447196158279646000,
  528.             0.23928736225213710000,
  529.             0.33186860228212767000,
  530.             0.42135127613063540000,
  531.             0.50689990893222950000,
  532.             0.58771575724076230000,
  533.             0.66304426693021520000,
  534.             0.73218211874028970000,
  535.             0.79448379596794250000,
  536.             0.84936761373256990000,
  537.             0.89632115576605220000,
  538.             0.93490607593773970000,
  539.             0.96476225558750640000,
  540.             0.98561151154526840000,
  541.             0.99726386184948160000
  542.         };

  543.         /** Weights for quadrature of order 32. */
  544.         private static final double[] W_32 = {
  545.             0.00701861000947013600,
  546.             0.01627439473090571200,
  547.             0.02539206530926214200,
  548.             0.03427386291302141000,
  549.             0.04283589802222658600,
  550.             0.05099805926237621600,
  551.             0.05868409347853559000,
  552.             0.06582222277636193000,
  553.             0.07234579410884862000,
  554.             0.07819389578707042000,
  555.             0.08331192422694673000,
  556.             0.08765209300440380000,
  557.             0.09117387869576390000,
  558.             0.09384439908080441000,
  559.             0.09563872007927487000,
  560.             0.09654008851472784000,
  561.             0.09654008851472784000,
  562.             0.09563872007927487000,
  563.             0.09384439908080441000,
  564.             0.09117387869576390000,
  565.             0.08765209300440380000,
  566.             0.08331192422694673000,
  567.             0.07819389578707042000,
  568.             0.07234579410884862000,
  569.             0.06582222277636193000,
  570.             0.05868409347853559000,
  571.             0.05099805926237621600,
  572.             0.04283589802222658600,
  573.             0.03427386291302141000,
  574.             0.02539206530926214200,
  575.             0.01627439473090571200,
  576.             0.00701861000947013600
  577.         };

  578.         /** Points for quadrature of order 40. */
  579.         private static final double[] P_40 = {
  580.             -0.99823770971055930000,
  581.             -0.99072623869945710000,
  582.             -0.97725994998377420000,
  583.             -0.95791681921379170000,
  584.             -0.93281280827867660000,
  585.             -0.90209880696887420000,
  586.             -0.86595950321225960000,
  587.             -0.82461223083331170000,
  588.             -0.77830565142651940000,
  589.             -0.72731825518992710000,
  590.             -0.67195668461417960000,
  591.             -0.61255388966798030000,
  592.             -0.54946712509512820000,
  593.             -0.48307580168617870000,
  594.             -0.41377920437160500000,
  595.             -0.34199409082575850000,
  596.             -0.26815218500725370000,
  597.             -0.19269758070137110000,
  598.             -0.11608407067525522000,
  599.             -0.03877241750605081600,
  600.             0.03877241750605081600,
  601.             0.11608407067525522000,
  602.             0.19269758070137110000,
  603.             0.26815218500725370000,
  604.             0.34199409082575850000,
  605.             0.41377920437160500000,
  606.             0.48307580168617870000,
  607.             0.54946712509512820000,
  608.             0.61255388966798030000,
  609.             0.67195668461417960000,
  610.             0.72731825518992710000,
  611.             0.77830565142651940000,
  612.             0.82461223083331170000,
  613.             0.86595950321225960000,
  614.             0.90209880696887420000,
  615.             0.93281280827867660000,
  616.             0.95791681921379170000,
  617.             0.97725994998377420000,
  618.             0.99072623869945710000,
  619.             0.99823770971055930000
  620.         };

  621.         /** Weights for quadrature of order 40. */
  622.         private static final double[] W_40 = {
  623.             0.00452127709853309800,
  624.             0.01049828453115270400,
  625.             0.01642105838190797300,
  626.             0.02224584919416689000,
  627.             0.02793700698002338000,
  628.             0.03346019528254786500,
  629.             0.03878216797447199000,
  630.             0.04387090818567333000,
  631.             0.04869580763507221000,
  632.             0.05322784698393679000,
  633.             0.05743976909939157000,
  634.             0.06130624249292891000,
  635.             0.06480401345660108000,
  636.             0.06791204581523394000,
  637.             0.07061164739128681000,
  638.             0.07288658239580408000,
  639.             0.07472316905796833000,
  640.             0.07611036190062619000,
  641.             0.07703981816424793000,
  642.             0.07750594797842482000,
  643.             0.07750594797842482000,
  644.             0.07703981816424793000,
  645.             0.07611036190062619000,
  646.             0.07472316905796833000,
  647.             0.07288658239580408000,
  648.             0.07061164739128681000,
  649.             0.06791204581523394000,
  650.             0.06480401345660108000,
  651.             0.06130624249292891000,
  652.             0.05743976909939157000,
  653.             0.05322784698393679000,
  654.             0.04869580763507221000,
  655.             0.04387090818567333000,
  656.             0.03878216797447199000,
  657.             0.03346019528254786500,
  658.             0.02793700698002338000,
  659.             0.02224584919416689000,
  660.             0.01642105838190797300,
  661.             0.01049828453115270400,
  662.             0.00452127709853309800
  663.         };

  664.         /** Points for quadrature of order 48. */
  665.         private static final double[] P_48 = {
  666.             -0.99877100725242610000,
  667.             -0.99353017226635080000,
  668.             -0.98412458372282700000,
  669.             -0.97059159254624720000,
  670.             -0.95298770316043080000,
  671.             -0.93138669070655440000,
  672.             -0.90587913671556960000,
  673.             -0.87657202027424800000,
  674.             -0.84358826162439350000,
  675.             -0.80706620402944250000,
  676.             -0.76715903251574020000,
  677.             -0.72403413092381470000,
  678.             -0.67787237963266400000,
  679.             -0.62886739677651370000,
  680.             -0.57722472608397270000,
  681.             -0.52316097472223300000,
  682.             -0.46690290475095840000,
  683.             -0.40868648199071680000,
  684.             -0.34875588629216070000,
  685.             -0.28736248735545555000,
  686.             -0.22476379039468908000,
  687.             -0.16122235606889174000,
  688.             -0.09700469920946270000,
  689.             -0.03238017096286937000,
  690.             0.03238017096286937000,
  691.             0.09700469920946270000,
  692.             0.16122235606889174000,
  693.             0.22476379039468908000,
  694.             0.28736248735545555000,
  695.             0.34875588629216070000,
  696.             0.40868648199071680000,
  697.             0.46690290475095840000,
  698.             0.52316097472223300000,
  699.             0.57722472608397270000,
  700.             0.62886739677651370000,
  701.             0.67787237963266400000,
  702.             0.72403413092381470000,
  703.             0.76715903251574020000,
  704.             0.80706620402944250000,
  705.             0.84358826162439350000,
  706.             0.87657202027424800000,
  707.             0.90587913671556960000,
  708.             0.93138669070655440000,
  709.             0.95298770316043080000,
  710.             0.97059159254624720000,
  711.             0.98412458372282700000,
  712.             0.99353017226635080000,
  713.             0.99877100725242610000
  714.         };

  715.         /** Weights for quadrature of order 48. */
  716.         private static final double[] W_48 = {
  717.             0.00315334605230596250,
  718.             0.00732755390127620800,
  719.             0.01147723457923446900,
  720.             0.01557931572294386600,
  721.             0.01961616045735556700,
  722.             0.02357076083932435600,
  723.             0.02742650970835688000,
  724.             0.03116722783279807000,
  725.             0.03477722256477045000,
  726.             0.03824135106583080600,
  727.             0.04154508294346483000,
  728.             0.04467456085669424000,
  729.             0.04761665849249054000,
  730.             0.05035903555385448000,
  731.             0.05289018948519365000,
  732.             0.05519950369998416500,
  733.             0.05727729210040315000,
  734.             0.05911483969839566000,
  735.             0.06070443916589384000,
  736.             0.06203942315989268000,
  737.             0.06311419228625403000,
  738.             0.06392423858464817000,
  739.             0.06446616443595010000,
  740.             0.06473769681268386000,
  741.             0.06473769681268386000,
  742.             0.06446616443595010000,
  743.             0.06392423858464817000,
  744.             0.06311419228625403000,
  745.             0.06203942315989268000,
  746.             0.06070443916589384000,
  747.             0.05911483969839566000,
  748.             0.05727729210040315000,
  749.             0.05519950369998416500,
  750.             0.05289018948519365000,
  751.             0.05035903555385448000,
  752.             0.04761665849249054000,
  753.             0.04467456085669424000,
  754.             0.04154508294346483000,
  755.             0.03824135106583080600,
  756.             0.03477722256477045000,
  757.             0.03116722783279807000,
  758.             0.02742650970835688000,
  759.             0.02357076083932435600,
  760.             0.01961616045735556700,
  761.             0.01557931572294386600,
  762.             0.01147723457923446900,
  763.             0.00732755390127620800,
  764.             0.00315334605230596250
  765.         };
  766.         // CHECKSTYLE: resume NoWhitespaceAfter

  767.         /** Node points. */
  768.         private final double[] nodePoints;

  769.         /** Node weights. */
  770.         private final double[] nodeWeights;

  771.         /** Creates a Gauss integrator of the given order.
  772.          *
  773.          *  @param numberOfPoints Order of the integration rule.
  774.          */
  775.         public GaussQuadrature(final int numberOfPoints) {

  776.             switch(numberOfPoints) {
  777.             case 12 :
  778.                 this.nodePoints  = P_12.clone();
  779.                 this.nodeWeights = W_12.clone();
  780.                 break;
  781.             case 16 :
  782.                 this.nodePoints  = P_16.clone();
  783.                 this.nodeWeights = W_16.clone();
  784.                 break;
  785.             case 20 :
  786.                 this.nodePoints  = P_20.clone();
  787.                 this.nodeWeights = W_20.clone();
  788.                 break;
  789.             case 24 :
  790.                 this.nodePoints  = P_24.clone();
  791.                 this.nodeWeights = W_24.clone();
  792.                 break;
  793.             case 32 :
  794.                 this.nodePoints  = P_32.clone();
  795.                 this.nodeWeights = W_32.clone();
  796.                 break;
  797.             case 40 :
  798.                 this.nodePoints  = P_40.clone();
  799.                 this.nodeWeights = W_40.clone();
  800.                 break;
  801.             case 48 :
  802.             default :
  803.                 this.nodePoints  = P_48.clone();
  804.                 this.nodeWeights = W_48.clone();
  805.                 break;
  806.             }

  807.         }

  808.         /** Integrates a given function on the given interval.
  809.          *
  810.          *  @param f Function to integrate.
  811.          *  @param lowerBound Lower bound of the integration interval.
  812.          *  @param upperBound Upper bound of the integration interval.
  813.          *  @return the integral of the weighted function.
  814.          */
  815.         public double[] integrate(final UnivariateVectorFunction f,
  816.                                   final double lowerBound, final double upperBound) {

  817.             final double[] adaptedPoints  = nodePoints.clone();
  818.             final double[] adaptedWeights = nodeWeights.clone();
  819.             transform(adaptedPoints, adaptedWeights, lowerBound, upperBound);
  820.             return basicIntegrate(f, adaptedPoints, adaptedWeights);
  821.         }

  822.         /** Performs a change of variable so that the integration
  823.          *  can be performed on an arbitrary interval {@code [a, b]}.
  824.          *  <p>
  825.          *  It is assumed that the natural interval is {@code [-1, 1]}.
  826.          *  </p>
  827.          *
  828.          * @param points  Points to adapt to the new interval.
  829.          * @param weights Weights to adapt to the new interval.
  830.          * @param a Lower bound of the integration interval.
  831.          * @param b Lower bound of the integration interval.
  832.          */
  833.         private void transform(final double[] points, final double[] weights,
  834.                                final double a, final double b) {
  835.             // Scaling
  836.             final double scale = (b - a) / 2;
  837.             final double shift = a + scale;
  838.             for (int i = 0; i < points.length; i++) {
  839.                 points[i]   = points[i] * scale + shift;
  840.                 weights[i] *= scale;
  841.             }
  842.         }

  843.         /** Returns an estimate of the integral of {@code f(x) * w(x)},
  844.          *  where {@code w} is a weight function that depends on the actual
  845.          *  flavor of the Gauss integration scheme.
  846.          *
  847.          * @param f Function to integrate.
  848.          * @param points  Nodes.
  849.          * @param weights Nodes weights.
  850.          * @return the integral of the weighted function.
  851.          */
  852.         private double[] basicIntegrate(final UnivariateVectorFunction f,
  853.                                         final double[] points,
  854.                                         final double[] weights) {
  855.             double x = points[0];
  856.             double w = weights[0];
  857.             double[] v = f.value(x);
  858.             final double[] y = new double[v.length];
  859.             for (int j = 0; j < v.length; j++) {
  860.                 y[j] = w * v[j];
  861.             }
  862.             final double[] t = y.clone();
  863.             final double[] c = new double[v.length];
  864.             final double[] s = t.clone();
  865.             for (int i = 1; i < points.length; i++) {
  866.                 x = points[i];
  867.                 w = weights[i];
  868.                 v = f.value(x);
  869.                 for (int j = 0; j < v.length; j++) {
  870.                     y[j] = w * v[j] - c[j];
  871.                     t[j] =  s[j] + y[j];
  872.                     c[j] = (t[j] - s[j]) - y[j];
  873.                     s[j] = t[j];
  874.                 }
  875.             }
  876.             return s;
  877.         }

  878.     }
  879. }