GlobalMappingFunctionModel.java

  1. /* Copyright 2002-2024 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.models.earth.troposphere;

  18. import org.hipparchus.CalculusFieldElement;
  19. import org.hipparchus.Field;
  20. import org.hipparchus.util.FastMath;
  21. import org.hipparchus.util.FieldSinCos;
  22. import org.hipparchus.util.MathArrays;
  23. import org.hipparchus.util.SinCos;
  24. import org.orekit.annotation.DefaultDataContext;
  25. import org.orekit.bodies.FieldGeodeticPoint;
  26. import org.orekit.bodies.GeodeticPoint;
  27. import org.orekit.data.DataContext;
  28. import org.orekit.models.earth.weather.FieldPressureTemperatureHumidity;
  29. import org.orekit.models.earth.weather.PressureTemperatureHumidity;
  30. import org.orekit.time.AbsoluteDate;
  31. import org.orekit.time.DateTimeComponents;
  32. import org.orekit.time.FieldAbsoluteDate;
  33. import org.orekit.time.TimeScale;
  34. import org.orekit.utils.FieldLegendrePolynomials;
  35. import org.orekit.utils.FieldTrackingCoordinates;
  36. import org.orekit.utils.LegendrePolynomials;
  37. import org.orekit.utils.TrackingCoordinates;

  38. /** The Global Mapping Function  model for radio techniques.
  39.  *  This model is an empirical mapping function. It only needs the
  40.  *  values of the station latitude, longitude, height and the
  41.  *  date for the computations.
  42.  *  <p>
  43.  *  The Global Mapping Function is based on spherical harmonics up
  44.  *  to degree and order of 9. It was developed to be consistent
  45.  *  with the {@link ViennaOneModel Vienna1} mapping function model.
  46.  *  </p>
  47.  *
  48.  *  @see "Boehm, J., A.E. Niell, P. Tregoning, H. Schuh (2006),
  49.  *       Global Mapping Functions (GMF): A new empirical mapping function based
  50.  *       on numerical weather model data, Geoph. Res. Letters, Vol. 33, L07304,
  51.  *       doi:10.1029/2005GL025545."
  52.  *
  53.  *  @see "Petit, G. and Luzum, B. (eds.), IERS Conventions (2010),
  54.  *       IERS Technical Note No. 36, BKG (2010)"
  55.  *
  56.  *  @author Bryan Cazabonne
  57.  *
  58.  */
  59. @SuppressWarnings("deprecation")
  60. public class GlobalMappingFunctionModel implements MappingFunction, TroposphereMappingFunction {

  61.     /** Multiplication factor for mapping function coefficients. */
  62.     private static final double FACTOR = 1.0e-5;

  63.     /** UTC time scale. */
  64.     private final TimeScale utc;

  65.     /** Build a new instance.
  66.      *
  67.      * <p>This constructor uses the {@link DataContext#getDefault() default data context}.
  68.      *
  69.      * @see #GlobalMappingFunctionModel(TimeScale)
  70.      */
  71.     @DefaultDataContext
  72.     public GlobalMappingFunctionModel() {
  73.         this(DataContext.getDefault().getTimeScales().getUTC());
  74.     }

  75.     /** Build a new instance.
  76.      * @param utc UTC time scale.
  77.      * @since 10.1
  78.      */
  79.     public GlobalMappingFunctionModel(final TimeScale utc) {
  80.         this.utc = utc;
  81.     }

  82.     /** {@inheritDoc} */
  83.     @Override
  84.     @Deprecated
  85.     public double[] mappingFactors(final double elevation, final GeodeticPoint point,
  86.                                    final AbsoluteDate date) {
  87.         return mappingFactors(new TrackingCoordinates(0.0, elevation, 0.0), point,
  88.                               TroposphericModelUtils.STANDARD_ATMOSPHERE,
  89.                               date);
  90.     }

  91.     /** {@inheritDoc} */
  92.     @Override
  93.     public double[] mappingFactors(final TrackingCoordinates trackingCoordinates,
  94.                                    final GeodeticPoint point,
  95.                                    final PressureTemperatureHumidity weather,
  96.                                    final AbsoluteDate date) {
  97.         // Day of year computation
  98.         final DateTimeComponents dtc = date.getComponents(utc);
  99.         final int dofyear = dtc.getDate().getDayOfYear();

  100.         // bh and ch constants (Boehm, J et al, 2006) | HYDROSTATIC PART
  101.         final double bh  = 0.0029;
  102.         final double c0h = 0.062;
  103.         final double c10h;
  104.         final double c11h;
  105.         final double psi;

  106.         // Latitude and longitude of the station
  107.         final double latitude  = point.getLatitude();
  108.         final double longitude = point.getLongitude();

  109.         if (FastMath.sin(latitude) > 0) {
  110.             // northern hemisphere case
  111.             c10h = 0.001;
  112.             c11h = 0.005;
  113.             psi  = 0.0;
  114.         } else {
  115.             // southern hemisphere case
  116.             c10h = 0.002;
  117.             c11h = 0.007;
  118.             psi  = FastMath.PI;
  119.         }

  120.         double t0 = 28;
  121.         if (latitude < 0) {
  122.             // southern hemisphere: t0 = 28 + an integer half of year
  123.             t0 += 183;
  124.         }
  125.         final double coef = ((dofyear + 1 - t0) / 365.25) * 2 * FastMath.PI + psi;
  126.         final double ch = c0h + ((FastMath.cos(coef) + 1) * (c11h / 2.0) + c10h) * (1.0 - FastMath.cos(latitude));

  127.         // bw and cw constants (Boehm, J et al, 2006) | WET PART
  128.         final double bw = 0.00146;
  129.         final double cw = 0.04391;

  130.         // Compute coefficients ah and aw with spherical harmonics Eq. 3 (Ref 1)

  131.         // Compute Legendre Polynomials Pnm(sin(phi))
  132.         final int degree = 9;
  133.         final int order  = 9;
  134.         final LegendrePolynomials p = new LegendrePolynomials(degree, order, FastMath.sin(latitude));

  135.         double a0Hydro   = 0.;
  136.         double amplHydro = 0.;
  137.         double a0Wet   = 0.;
  138.         double amplWet = 0.;
  139.         final ABCoefficients abCoef = new ABCoefficients();
  140.         int j = 0;
  141.         for (int n = 0; n <= 9; n++) {
  142.             for (int m = 0; m <= n; m++) {
  143.                 // Sine and cosine of m * longitude
  144.                 final SinCos sc = FastMath.sinCos(m * longitude);
  145.                 // Compute coefficients
  146.                 a0Hydro   = a0Hydro + (abCoef.getAHMean(j) * p.getPnm(n, m) * sc.cos() +
  147.                                        abCoef.getBHMean(j) * p.getPnm(n, m) * sc.sin()) * FACTOR;

  148.                 a0Wet     = a0Wet + (abCoef.getAWMean(j) * p.getPnm(n, m) * sc.cos() +
  149.                                      abCoef.getBWMean(j) * p.getPnm(n, m) * sc.sin()) * FACTOR;

  150.                 amplHydro = amplHydro + (abCoef.getAHAmplitude(j) * p.getPnm(n, m) * sc.cos() +
  151.                                          abCoef.getBHAmplitude(j) * p.getPnm(n, m) * sc.sin()) * FACTOR;

  152.                 amplWet   = amplWet + (abCoef.getAWAmplitude(j) * p.getPnm(n, m) * sc.cos() +
  153.                                        abCoef.getBWAmplitude(j) * p.getPnm(n, m) * sc.sin()) * FACTOR;

  154.                 j = j + 1;
  155.             }
  156.         }

  157.         // Eq. 2 (Ref 1)
  158.         final double ah = a0Hydro + amplHydro * FastMath.cos(coef - psi);
  159.         final double aw = a0Wet + amplWet * FastMath.cos(coef - psi);

  160.         final double[] function = new double[2];
  161.         function[0] = TroposphericModelUtils.mappingFunction(ah, bh, ch,
  162.                                                              trackingCoordinates.getElevation());
  163.         function[1] = TroposphericModelUtils.mappingFunction(aw, bw, cw,
  164.                                                              trackingCoordinates.getElevation());

  165.         // Apply height correction
  166.         final double correction = TroposphericModelUtils.computeHeightCorrection(trackingCoordinates.getElevation(),
  167.                                                                                  point.getAltitude());
  168.         function[0] = function[0] + correction;

  169.         return function;
  170.     }

  171.     /** {@inheritDoc} */
  172.     @Override
  173.     @Deprecated
  174.     public <T extends CalculusFieldElement<T>> T[] mappingFactors(final T elevation, final FieldGeodeticPoint<T> point,
  175.                                                                   final FieldAbsoluteDate<T> date) {
  176.         return mappingFactors(new FieldTrackingCoordinates<>(date.getField().getZero(), elevation, date.getField().getZero()),
  177.                               point,
  178.                               new FieldPressureTemperatureHumidity<>(date.getField(),
  179.                                                                      TroposphericModelUtils.STANDARD_ATMOSPHERE),
  180.                               date);
  181.     }

  182.     /** {@inheritDoc} */
  183.     @Override
  184.     public <T extends CalculusFieldElement<T>> T[] mappingFactors(final FieldTrackingCoordinates<T> trackingCoordinates,
  185.                                                                   final FieldGeodeticPoint<T> point,
  186.                                                                   final FieldPressureTemperatureHumidity<T> weather,
  187.                                                                   final FieldAbsoluteDate<T> date) {
  188.         // Day of year computation
  189.         final DateTimeComponents dtc = date.getComponents(utc);
  190.         final int dofyear = dtc.getDate().getDayOfYear();

  191.         final Field<T> field = date.getField();
  192.         final T zero = field.getZero();

  193.         // bh and ch constants (Boehm, J et al, 2006) | HYDROSTATIC PART
  194.         final T bh  = zero.newInstance(0.0029);
  195.         final T c0h = zero.newInstance(0.062);
  196.         final T c10h;
  197.         final T c11h;
  198.         final T psi;

  199.         // Latitude and longitude of the station
  200.         final T latitude  = point.getLatitude();
  201.         final T longitude = point.getLongitude();

  202.         // sin(latitude) > 0 -> northern hemisphere
  203.         if (FastMath.sin(latitude.getReal()) > 0) {
  204.             c10h = zero.newInstance(0.001);
  205.             c11h = zero.newInstance(0.005);
  206.             psi  = zero;
  207.         } else {
  208.             c10h = zero.newInstance(0.002);
  209.             c11h = zero.newInstance(0.007);
  210.             psi  = zero.getPi();
  211.         }

  212.         double t0 = 28;
  213.         if (latitude.getReal() < 0) {
  214.             // southern hemisphere: t0 = 28 + an integer half of year
  215.             t0 += 183;
  216.         }
  217.         final T coef = psi.add(zero.getPi().multiply(2.0).multiply((dofyear + 1 - t0) / 365.25));
  218.         final T ch = c11h.divide(2.0).multiply(FastMath.cos(coef).add(1.0)).add(c10h).multiply(FastMath.cos(latitude).negate().add(1.0)).add(c0h);

  219.         // bw and cw constants (Boehm, J et al, 2006) | WET PART
  220.         final T bw = zero.newInstance(0.00146);
  221.         final T cw = zero.newInstance(0.04391);

  222.         // Compute coefficients ah and aw with spherical harmonics Eq. 3 (Ref 1)

  223.         // Compute Legendre Polynomials Pnm(sin(phi))
  224.         final int degree = 9;
  225.         final int order  = 9;
  226.         final FieldLegendrePolynomials<T> p = new FieldLegendrePolynomials<>(degree, order, FastMath.sin(latitude));

  227.         T a0Hydro   = zero;
  228.         T amplHydro = zero;
  229.         T a0Wet     = zero;
  230.         T amplWet   = zero;
  231.         final ABCoefficients abCoef = new ABCoefficients();
  232.         int j = 0;
  233.         for (int n = 0; n <= 9; n++) {
  234.             for (int m = 0; m <= n; m++) {
  235.                 // Sine and cosine of m * longitude
  236.                 final FieldSinCos<T> sc = FastMath.sinCos(longitude.multiply(m));

  237.                 // Compute coefficients
  238.                 a0Hydro   = a0Hydro.add(p.getPnm(n, m).multiply(abCoef.getAHMean(j)).multiply(sc.cos()).
  239.                                     add(p.getPnm(n, m).multiply(abCoef.getBHMean(j)).multiply(sc.sin())).
  240.                                     multiply(FACTOR));

  241.                 a0Wet     = a0Wet.add(p.getPnm(n, m).multiply(abCoef.getAWMean(j)).multiply(sc.cos()).
  242.                                   add(p.getPnm(n, m).multiply(abCoef.getBWMean(j)).multiply(sc.sin())).
  243.                                   multiply(FACTOR));

  244.                 amplHydro = amplHydro.add(p.getPnm(n, m).multiply(abCoef.getAHAmplitude(j)).multiply(sc.cos()).
  245.                                       add(p.getPnm(n, m).multiply(abCoef.getBHAmplitude(j)).multiply(sc.sin())).
  246.                                       multiply(FACTOR));

  247.                 amplWet   = amplWet.add(p.getPnm(n, m).multiply(abCoef.getAWAmplitude(j)).multiply(sc.cos()).
  248.                                     add(p.getPnm(n, m).multiply(abCoef.getBWAmplitude(j)).multiply(sc.sin())).
  249.                                     multiply(FACTOR));

  250.                 j = j + 1;
  251.             }
  252.         }

  253.         // Eq. 2 (Ref 1)
  254.         final T ah = a0Hydro.add(amplHydro.multiply(FastMath.cos(coef.subtract(psi))));
  255.         final T aw = a0Wet.add(amplWet.multiply(FastMath.cos(coef.subtract(psi))));

  256.         final T[] function = MathArrays.buildArray(field, 2);
  257.         function[0] = TroposphericModelUtils.mappingFunction(ah, bh, ch,
  258.                                                              trackingCoordinates.getElevation());
  259.         function[1] = TroposphericModelUtils.mappingFunction(aw, bw, cw,
  260.                                                              trackingCoordinates.getElevation());

  261.         // Apply height correction
  262.         final T correction = TroposphericModelUtils.computeHeightCorrection(trackingCoordinates.getElevation(),
  263.                                                                             point.getAltitude(),
  264.                                                                             field);
  265.         function[0] = function[0].add(correction);

  266.         return function;
  267.     }

  268.     private static class ABCoefficients {

  269.         /** Mean hydrostatic coefficients a.*/
  270.         private static final double[] AH_MEAN = {
  271.             1.2517e02,
  272.             8.503e-01,
  273.             6.936e-02,
  274.             -6.760e+00,
  275.             1.771e-01,
  276.             1.130e-02,
  277.             5.963e-01,
  278.             1.808e-02,
  279.             2.801e-03,
  280.             -1.414e-03,
  281.             -1.212e+00,
  282.             9.300e-02,
  283.             3.683e-03,
  284.             1.095e-03,
  285.             4.671e-05,
  286.             3.959e-01,
  287.             -3.867e-02,
  288.             5.413e-03,
  289.             -5.289e-04,
  290.             3.229e-04,
  291.             2.067e-05,
  292.             3.000e-01,
  293.             2.031e-02,
  294.             5.900e-03,
  295.             4.573e-04,
  296.             -7.619e-05,
  297.             2.327e-06,
  298.             3.845e-06,
  299.             1.182e-01,
  300.             1.158e-02,
  301.             5.445e-03,
  302.             6.219e-05,
  303.             4.204e-06,
  304.             -2.093e-06,
  305.             1.540e-07,
  306.             -4.280e-08,
  307.             -4.751e-01,
  308.             -3.490e-02,
  309.             1.758e-03,
  310.             4.019e-04,
  311.             -2.799e-06,
  312.             -1.287e-06,
  313.             5.468e-07,
  314.             7.580e-08,
  315.             -6.300e-09,
  316.             -1.160e-01,
  317.             8.301e-03,
  318.             8.771e-04,
  319.             9.955e-05,
  320.             -1.718e-06,
  321.             -2.012e-06,
  322.             1.170e-08,
  323.             1.790e-08,
  324.             -1.300e-09,
  325.             1.000e-10
  326.         };

  327.         /** Mean hydrostatic coefficients b.*/
  328.         private static final double[] BH_MEAN = {
  329.             0.000e+00,
  330.             0.000e+00,
  331.             3.249e-02,
  332.             0.000e+00,
  333.             3.324e-02,
  334.             1.850e-02,
  335.             0.000e+00,
  336.             -1.115e-01,
  337.             2.519e-02,
  338.             4.923e-03,
  339.             0.000e+00,
  340.             2.737e-02,
  341.             1.595e-02,
  342.             -7.332e-04,
  343.             1.933e-04,
  344.             0.000e+00,
  345.             -4.796e-02,
  346.             6.381e-03,
  347.             -1.599e-04,
  348.             -3.685e-04,
  349.             1.815e-05,
  350.             0.000e+00,
  351.             7.033e-02,
  352.             2.426e-03,
  353.             -1.111e-03,
  354.             -1.357e-04,
  355.             -7.828e-06,
  356.             2.547e-06,
  357.             0.000e+00,
  358.             5.779e-03,
  359.             3.133e-03,
  360.             -5.312e-04,
  361.             -2.028e-05,
  362.             2.323e-07,
  363.             -9.100e-08,
  364.             -1.650e-08,
  365.             0.000e+00,
  366.             3.688e-02,
  367.             -8.638e-04,
  368.             -8.514e-05,
  369.             -2.828e-05,
  370.             5.403e-07,
  371.             4.390e-07,
  372.             1.350e-08,
  373.             1.800e-09,
  374.             0.000e+00,
  375.             -2.736e-02,
  376.             -2.977e-04,
  377.             8.113e-05,
  378.             2.329e-07,
  379.             8.451e-07,
  380.             4.490e-08,
  381.             -8.100e-09,
  382.             -1.500e-09,
  383.             2.000e-10
  384.         };

  385.         /** Amplitude for hydrostatic coefficients a.*/
  386.         private static final double[] AH_AMPL = {
  387.             -2.738e-01,
  388.             -2.837e+00,
  389.             1.298e-02,
  390.             -3.588e-01,
  391.             2.413e-02,
  392.             3.427e-02,
  393.             -7.624e-01,
  394.             7.272e-02,
  395.             2.160e-02,
  396.             -3.385e-03,
  397.             4.424e-01,
  398.             3.722e-02,
  399.             2.195e-02,
  400.             -1.503e-03,
  401.             2.426e-04,
  402.             3.013e-01,
  403.             5.762e-02,
  404.             1.019e-02,
  405.             -4.476e-04,
  406.             6.790e-05,
  407.             3.227e-05,
  408.             3.123e-01,
  409.             -3.535e-02,
  410.             4.840e-03,
  411.             3.025e-06,
  412.             -4.363e-05,
  413.             2.854e-07,
  414.             -1.286e-06,
  415.             -6.725e-01,
  416.             -3.730e-02,
  417.             8.964e-04,
  418.             1.399e-04,
  419.             -3.990e-06,
  420.             7.431e-06,
  421.             -2.796e-07,
  422.             -1.601e-07,
  423.             4.068e-02,
  424.             -1.352e-02,
  425.             7.282e-04,
  426.             9.594e-05,
  427.             2.070e-06,
  428.             -9.620e-08,
  429.             -2.742e-07,
  430.             -6.370e-08,
  431.             -6.300e-09,
  432.             8.625e-02,
  433.             -5.971e-03,
  434.             4.705e-04,
  435.             2.335e-05,
  436.             4.226e-06,
  437.             2.475e-07,
  438.             -8.850e-08,
  439.             -3.600e-08,
  440.             -2.900e-09,
  441.             0.000e+00
  442.         };

  443.         /** Amplitude for hydrostatic coefficients b.*/
  444.         private static final double[] BH_AMPL = {
  445.             0.000e+00,
  446.             0.000e+00,
  447.             -1.136e-01,
  448.             0.000e+00,
  449.             -1.868e-01,
  450.             -1.399e-02,
  451.             0.000e+00,
  452.             -1.043e-01,
  453.             1.175e-02,
  454.             -2.240e-03,
  455.             0.000e+00,
  456.             -3.222e-02,
  457.             1.333e-02,
  458.             -2.647e-03,
  459.             -2.316e-05,
  460.             0.000e+00,
  461.             5.339e-02,
  462.             1.107e-02,
  463.             -3.116e-03,
  464.             -1.079e-04,
  465.             -1.299e-05,
  466.             0.000e+00,
  467.             4.861e-03,
  468.             8.891e-03,
  469.             -6.448e-04,
  470.             -1.279e-05,
  471.             6.358e-06,
  472.             -1.417e-07,
  473.             0.000e+00,
  474.             3.041e-02,
  475.             1.150e-03,
  476.             -8.743e-04,
  477.             -2.781e-05,
  478.             6.367e-07,
  479.             -1.140e-08,
  480.             -4.200e-08,
  481.             0.000e+00,
  482.             -2.982e-02,
  483.             -3.000e-03,
  484.             1.394e-05,
  485.             -3.290e-05,
  486.             -1.705e-07,
  487.             7.440e-08,
  488.             2.720e-08,
  489.             -6.600e-09,
  490.             0.000e+00,
  491.             1.236e-02,
  492.             -9.981e-04,
  493.             -3.792e-05,
  494.             -1.355e-05,
  495.             1.162e-06,
  496.             -1.789e-07,
  497.             1.470e-08,
  498.             -2.400e-09,
  499.             -4.000e-10
  500.         };

  501.         /** Mean wet coefficients a.*/
  502.         private static final double[] AW_MEAN = {
  503.             5.640e+01,
  504.             1.555e+00,
  505.             -1.011e+00,
  506.             -3.975e+00,
  507.             3.171e-02,
  508.             1.065e-01,
  509.             6.175e-01,
  510.             1.376e-01,
  511.             4.229e-02,
  512.             3.028e-03,
  513.             1.688e+00,
  514.             -1.692e-01,
  515.             5.478e-02,
  516.             2.473e-02,
  517.             6.059e-04,
  518.             2.278e+00,
  519.             6.614e-03,
  520.             -3.505e-04,
  521.             -6.697e-03,
  522.             8.402e-04,
  523.             7.033e-04,
  524.             -3.236e+00,
  525.             2.184e-01,
  526.             -4.611e-02,
  527.             -1.613e-02,
  528.             -1.604e-03,
  529.             5.420e-05,
  530.             7.922e-05,
  531.             -2.711e-01,
  532.             -4.406e-01,
  533.             -3.376e-02,
  534.             -2.801e-03,
  535.             -4.090e-04,
  536.             -2.056e-05,
  537.             6.894e-06,
  538.             2.317e-06,
  539.             1.941e+00,
  540.             -2.562e-01,
  541.             1.598e-02,
  542.             5.449e-03,
  543.             3.544e-04,
  544.             1.148e-05,
  545.             7.503e-06,
  546.             -5.667e-07,
  547.             -3.660e-08,
  548.             8.683e-01,
  549.             -5.931e-02,
  550.             -1.864e-03,
  551.             -1.277e-04,
  552.             2.029e-04,
  553.             1.269e-05,
  554.             1.629e-06,
  555.             9.660e-08,
  556.             -1.015e-07,
  557.             -5.000e-10
  558.         };

  559.         /** Mean wet coefficients b.*/
  560.         private static final double[] BW_MEAN = {
  561.             0.000e+00,
  562.             0.000e+00,
  563.             2.592e-01,
  564.             0.000e+00,
  565.             2.974e-02,
  566.             -5.471e-01,
  567.             0.000e+00,
  568.             -5.926e-01,
  569.             -1.030e-01,
  570.             -1.567e-02,
  571.             0.000e+00,
  572.             1.710e-01,
  573.             9.025e-02,
  574.             2.689e-02,
  575.             2.243e-03,
  576.             0.000e+00,
  577.             3.439e-01,
  578.             2.402e-02,
  579.             5.410e-03,
  580.             1.601e-03,
  581.             9.669e-05,
  582.             0.000e+00,
  583.             9.502e-02,
  584.             -3.063e-02,
  585.             -1.055e-03,
  586.             -1.067e-04,
  587.             -1.130e-04,
  588.             2.124e-05,
  589.             0.000e+00,
  590.             -3.129e-01,
  591.             8.463e-03,
  592.             2.253e-04,
  593.             7.413e-05,
  594.             -9.376e-05,
  595.             -1.606e-06,
  596.             2.060e-06,
  597.             0.000e+00,
  598.             2.739e-01,
  599.             1.167e-03,
  600.             -2.246e-05,
  601.             -1.287e-04,
  602.             -2.438e-05,
  603.             -7.561e-07,
  604.             1.158e-06,
  605.             4.950e-08,
  606.             0.000e+00,
  607.             -1.344e-01,
  608.             5.342e-03,
  609.             3.775e-04,
  610.             -6.756e-05,
  611.             -1.686e-06,
  612.             -1.184e-06,
  613.             2.768e-07,
  614.             2.730e-08,
  615.             5.700e-09
  616.         };

  617.         /** Amplitude for wet coefficients a.*/
  618.         private static final double[] AW_AMPL = {
  619.             1.023e-01,
  620.             -2.695e+00,
  621.             3.417e-01,
  622.             -1.405e-01,
  623.             3.175e-01,
  624.             2.116e-01,
  625.             3.536e+00,
  626.             -1.505e-01,
  627.             -1.660e-02,
  628.             2.967e-02,
  629.             3.819e-01,
  630.             -1.695e-01,
  631.             -7.444e-02,
  632.             7.409e-03,
  633.             -6.262e-03,
  634.             -1.836e+00,
  635.             -1.759e-02,
  636.             -6.256e-02,
  637.             -2.371e-03,
  638.             7.947e-04,
  639.             1.501e-04,
  640.             -8.603e-01,
  641.             -1.360e-01,
  642.             -3.629e-02,
  643.             -3.706e-03,
  644.             -2.976e-04,
  645.             1.857e-05,
  646.             3.021e-05,
  647.             2.248e+00,
  648.             -1.178e-01,
  649.             1.255e-02,
  650.             1.134e-03,
  651.             -2.161e-04,
  652.             -5.817e-06,
  653.             8.836e-07,
  654.             -1.769e-07,
  655.             7.313e-01,
  656.             -1.188e-01,
  657.             1.145e-02,
  658.             1.011e-03,
  659.             1.083e-04,
  660.             2.570e-06,
  661.             -2.140e-06,
  662.             -5.710e-08,
  663.             2.000e-08,
  664.             -1.632e+00,
  665.             -6.948e-03,
  666.             -3.893e-03,
  667.             8.592e-04,
  668.             7.577e-05,
  669.             4.539e-06,
  670.             -3.852e-07,
  671.             -2.213e-07,
  672.             -1.370e-08,
  673.             5.800e-09
  674.         };

  675.         /** Amplitude for wet coefficients b.*/
  676.         private static final double[] BW_AMPL = {
  677.             0.000e+00,
  678.             0.000e+00,
  679.             -8.865e-02,
  680.             0.000e+00,
  681.             -4.309e-01,
  682.             6.340e-02,
  683.             0.000e+00,
  684.             1.162e-01,
  685.             6.176e-02,
  686.             -4.234e-03,
  687.             0.000e+00,
  688.             2.530e-01,
  689.             4.017e-02,
  690.             -6.204e-03,
  691.             4.977e-03,
  692.             0.000e+00,
  693.             -1.737e-01,
  694.             -5.638e-03,
  695.             1.488e-04,
  696.             4.857e-04,
  697.             -1.809e-04,
  698.             0.000e+00,
  699.             -1.514e-01,
  700.             -1.685e-02,
  701.             5.333e-03,
  702.             -7.611e-05,
  703.             2.394e-05,
  704.             8.195e-06,
  705.             0.000e+00,
  706.             9.326e-02,
  707.             -1.275e-02,
  708.             -3.071e-04,
  709.             5.374e-05,
  710.             -3.391e-05,
  711.             -7.436e-06,
  712.             6.747e-07,
  713.             0.000e+00,
  714.             -8.637e-02,
  715.             -3.807e-03,
  716.             -6.833e-04,
  717.             -3.861e-05,
  718.             -2.268e-05,
  719.             1.454e-06,
  720.             3.860e-07,
  721.             -1.068e-07,
  722.             0.000e+00,
  723.             -2.658e-02,
  724.             -1.947e-03,
  725.             7.131e-04,
  726.             -3.506e-05,
  727.             1.885e-07,
  728.             5.792e-07,
  729.             3.990e-08,
  730.             2.000e-08,
  731.             -5.700e-09
  732.         };

  733.         /** Build a new instance. */
  734.         ABCoefficients() {

  735.         }

  736.         /** Get the value of the mean hydrostatique coefficient ah for the given index.
  737.          * @param index index
  738.          * @return the mean hydrostatique coefficient ah for the given index
  739.          */
  740.         public double getAHMean(final int index) {
  741.             return AH_MEAN[index];
  742.         }

  743.         /** Get the value of the mean hydrostatique coefficient bh for the given index.
  744.          * @param index index
  745.          * @return the mean hydrostatique coefficient bh for the given index
  746.          */
  747.         public double getBHMean(final int index) {
  748.             return BH_MEAN[index];
  749.         }

  750.         /** Get the value of the mean wet coefficient aw for the given index.
  751.          * @param index index
  752.          * @return the mean wet coefficient aw for the given index
  753.          */
  754.         public double getAWMean(final int index) {
  755.             return AW_MEAN[index];
  756.         }

  757.         /** Get the value of the mean wet coefficient bw for the given index.
  758.          * @param index index
  759.          * @return the mean wet coefficient bw for the given index
  760.          */
  761.         public double getBWMean(final int index) {
  762.             return BW_MEAN[index];
  763.         }

  764.         /** Get the value of the amplitude of the hydrostatique coefficient ah for the given index.
  765.          * @param index index
  766.          * @return the amplitude of the hydrostatique coefficient ah for the given index
  767.          */
  768.         public double getAHAmplitude(final int index) {
  769.             return AH_AMPL[index];
  770.         }

  771.         /** Get the value of the amplitude of the hydrostatique coefficient bh for the given index.
  772.          * @param index index
  773.          * @return the amplitude of the hydrostatique coefficient bh for the given index
  774.          */
  775.         public double getBHAmplitude(final int index) {
  776.             return BH_AMPL[index];
  777.         }

  778.         /** Get the value of the amplitude of the wet coefficient aw for the given index.
  779.          * @param index index
  780.          * @return the amplitude of the wet coefficient aw for the given index
  781.          */
  782.         public double getAWAmplitude(final int index) {
  783.             return AW_AMPL[index];
  784.         }

  785.         /** Get the value of the amplitude of the wet coefficient bw for the given index.
  786.          * @param index index
  787.          * @return the amplitude of the wet coefficient bw for the given index
  788.          */
  789.         public double getBWAmplitude(final int index) {
  790.             return BW_AMPL[index];
  791.         }
  792.     }
  793. }