NRLMSISE00.java

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

  18. import org.hipparchus.CalculusFieldElement;
  19. import org.hipparchus.Field;
  20. import org.hipparchus.exception.LocalizedCoreFormats;
  21. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  22. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  23. import org.hipparchus.util.FastMath;
  24. import org.hipparchus.util.FieldSinCos;
  25. import org.hipparchus.util.MathArrays;
  26. import org.hipparchus.util.SinCos;
  27. import org.orekit.annotation.DefaultDataContext;
  28. import org.orekit.bodies.BodyShape;
  29. import org.orekit.bodies.FieldGeodeticPoint;
  30. import org.orekit.bodies.GeodeticPoint;
  31. import org.orekit.data.DataContext;
  32. import org.orekit.errors.OrekitException;
  33. import org.orekit.errors.OrekitMessages;
  34. import org.orekit.frames.Frame;
  35. import org.orekit.time.AbsoluteDate;
  36. import org.orekit.time.DateTimeComponents;
  37. import org.orekit.time.FieldAbsoluteDate;
  38. import org.orekit.time.TimeComponents;
  39. import org.orekit.time.TimeScale;
  40. import org.orekit.utils.IERSConventions;
  41. import org.orekit.utils.ExtendedPositionProvider;

  42. import java.util.Arrays;


  43. /** This class implements the mathematical representation of the 2001
  44.  *  Naval Research Laboratory Mass Spectrometer and Incoherent Scatter
  45.  *  Radar Exosphere (NRLMSISE-00) of the MSIS® class model.
  46.  *  <p>
  47.  *  NRLMSISE-00 calculates the neutral atmosphere empirical model from the surface
  48.  *  to lower exosphere (0 to 1000 km) and provides:
  49.  *  <ul>
  50.  *  <li>Exospheric Temperature above Input Position (K)</li>
  51.  *  <li>Local Temperature at Input Position (K)</li>
  52.  *  <li>Total Mass-Density at Input Position (kg/m³)</li>
  53.  *  <li>Partial Densities at Input Position (1/m³) for:
  54.  *  <ul>
  55.  *      <li>He,</li>
  56.  *      <li>H,</li>
  57.  *      <li>N,</li>
  58.  *      <li>O,</li>
  59.  *      <li>Ar,</li>
  60.  *      <li>N2,</li>
  61.  *      <li>O2,</li>
  62.  *      <li>anomalous oxygen.</li>
  63.  *  </ul>
  64.  *  </li>
  65.  *  </ul>
  66.  *  <p>
  67.  *  The model needs geographical and time information to compute general values,
  68.  *  but also needs space weather data:
  69.  *  <ul>
  70.  *  <li>mean and daily solar flux,</li>
  71.  *  <li>geomagnetic indices.</li>
  72.  *  </ul>
  73.  *  <p>
  74.  *  Switches can be used to turn on and off particular variations:<br>
  75.  *  0 is off, 1 is on, and 2 is main effects off but cross terms on.<br>
  76.  *  The standard value is 1 for all the 23 available switches.<br>
  77.  *  Function of each switch according to its number:
  78.  *  <ul>
  79.  *  <li>#1 - F10.7 effect on mean</li>
  80.  *  <li>#2 - Independent of time</li>
  81.  *  <li>#3 - Symmetrical annual</li>
  82.  *  <li>#4 - Symmetrical semiannual</li>
  83.  *  <li>#5 - Asymmetrical annual</li>
  84.  *  <li>#6 - Asymmetrical semiannual</li>
  85.  *  <li>#7 - Diurnal</li>
  86.  *  <li>#8 - Semidiurnal</li>
  87.  *  <li>#9 - Daily Ap [**]</li>
  88.  *  <li>#10 - All UT, longitudinal effects</li>
  89.  *  <li>#11 - Longitudinal</li>
  90.  *  <li>#12 - UT and mixed UT, longitudinal</li>
  91.  *  <li>#13 - Mixed AP, UT, longitudinal</li>
  92.  *  <li>#14 - Terdiurnal</li>
  93.  *  <li>#15 - Departures from diffusive equilibrium</li>
  94.  *  <li>#16 - All exospheric temperature variations</li>
  95.  *  <li>#17 - All variations from 120 km temperature (TLB)</li>
  96.  *  <li>#18 - All lower thermosphere (TN1) temperature variations</li>
  97.  *  <li>#19 - All 120 km gradient (S) variations</li>
  98.  *  <li>#20 - All upper stratosphere (TN2) temperature variations</li>
  99.  *  <li>#21 - All variations from 120 km values (ZLB)</li>
  100.  *  <li>#22 - All lower mesosphere temperature (TN3) variations</li>
  101.  *  <li>#23 - Turbopause scale height variations</li>
  102.  *  </ul>
  103.  *  [**] Switch #9 is a bit specific:
  104.  *  <ul>
  105.  *  <li>set to  1, the daily Ap only is used (first element of ap array),</li>
  106.  *  <li>set to -1, the entire array of ap is used, including 3 hr ap indices.</li>
  107.  *  </ul>
  108.  *  <p>
  109.  *  The NRLMSISE-00 model was developed by Mike Picone, Alan Hedin, and Doug Drob.<br>
  110.  *  They also wrote a NRLMSISE-00 distribution package in FORTRAN available at:<br>
  111.  *  ftp://hanna.ccmc.gsfc.nasa.gov/pub/modelweb/atmospheric/msis/nrlmsise00/<br>
  112.  *  <br>
  113.  *  Dominik Brodowski implemented a C version of the NRLMSISE-00 model available at:<br>
  114.  *  http://www.brodo.de/space/nrlmsise/index.html
  115.  *  <p>
  116.  *  Instances of this class are immutable.
  117.  *  </p>
  118.  *
  119.  *  @author Mike Picone &amp; al (Naval Research Laboratory), 2001: FORTRAN routine
  120.  *  @author Dominik Brodowski, 2004: C routine
  121.  *  @author Pascal Parraud, 2016: Java translation
  122.  *  @since 8.1
  123.  */
  124. public class NRLMSISE00 extends AbstractSunInfluencedAtmosphere {
  125.     // Constants

  126.     /** Identifier for helium density. */
  127.     private static final int HELIUM = 0;

  128.     /** Identifier for atomic oxygen density. */
  129.     private static final int ATOMIC_OXYGEN = 1;

  130.     /** Identifier for molecular nitrogen density. */
  131.     private static final int MOLECULAR_NITROGEN = 2;

  132.     /** Identifier for molecular oxygen density. */
  133.     private static final int MOLECULAR_OXYGEN = 3;

  134.     /** Identifier for argon density. */
  135.     private static final int ARGON = 4;

  136.     /** Identifier for atomic nitrogen density. */
  137.     private static final int TOTAL_MASS = 5;

  138.     /** Identifier for hydrogen density. */
  139.     private static final int HYDROGEN = 6;

  140.     /** Identifier for atomic nitrogen density. */
  141.     private static final int ATOMIC_NITROGEN = 7;

  142.     /** Identifier for anomalous oxygen density. */
  143.     private static final int ANOMALOUS_OXYGEN = 8;

  144.     /** Identifier for exospheric temperature. */
  145.     private static final int EXOSPHERIC = 0;

  146.     /** Identifier for temperature at altitude. */
  147.     private static final int ALTITUDE = 1;

  148.     // CONVERSION CONSTANTS

  149.     /** Conversion from degree to radian. */
  150.     private static final double DEG_TO_RAD = 1.74533e-2;

  151.     /** Conversion from day to radian. */
  152.     private static final double DAY_TO_RAD = 1.72142e-2;

  153.     /** Conversion from hour to radian. */
  154.     private static final double HOUR_TO_RAD = 0.2618;

  155.     /** Conversion from second to radian. */
  156.     private static final double SEC_TO_RAD = 7.2722e-5;

  157.     // EARTH GEOPHYSICAL CONSTANTS

  158.     /** Reference latitude (°). */
  159.     private static final double LAT_REF = 45.;

  160.     /** Reference gravity on Earth surface at reference latitude (cm/s2). */
  161.     private static final double G_REF = 980.616;

  162.     // CHEMICAL CONSTANTS

  163.     /** Unified atomic mass unit (kg). */
  164.     private static final double AMU = 1.66e-27;

  165.     /** Gas constant (inverse of). */
  166.     private static final double R_GAS = 831.4;

  167.     /** Hydrogen atomic mass. */
  168.     private static final double H_MASS = 1.;

  169.     /** Helium atomic mass. */
  170.     private static final double HE_MASS = 4.;

  171.     /** Nitrogen atomic mass. */
  172.     private static final double N_MASS = 14.;

  173.     /** N2 molecular mass. */
  174.     private static final double N2_MASS = 2. * N_MASS;

  175.     /** Oxygen atomic mass. */
  176.     private static final double O_MASS = 16.;

  177.     /** O2 molecular mass. */
  178.     private static final double O2_MASS = 2. * O_MASS;

  179.     /** Argon atomic mass. */
  180.     private static final double AR_MASS = 40.;

  181.     // NRL MSISE 2000 SPECIFIC CONSTANTS

  182.     /** Reference average flux. */
  183.     private static final double FLUX_REF = 150.;

  184.     /** Array of altitudes #1. */
  185.     private static final double[] ZN1 = {123.435, 110.0, 100.0, 90.0, 72.5};

  186.     /** Array of altitudes #2. */
  187.     private static final double[] ZN2 = {72.5, 55.0, 45.0, 32.5};

  188.     /** Array of altitudes #3. */
  189.     private static final double[] ZN3 = {32.5, 20.0, 15.0, 10.0, 0.0};

  190.     /** Mix altitude (km). */
  191.     private static final double ZMIX = 62.5;

  192.     /** NRLMSISE-00 data: temperature pt[150]. */
  193.     private static final double[] PT = {
  194.         9.86573e-01, 1.62228e-02, 1.55270e-02, -1.04323e-01, -3.75801e-03,
  195.         -1.18538e-03, -1.24043e-01, 4.56820e-03, 8.76018e-03, -1.36235e-01,
  196.         -3.52427e-02, 8.84181e-03, -5.92127e-03, -8.61650e+00, 0.00000e+00,
  197.         1.28492e-02, 0.00000e+00, 1.30096e+02, 1.04567e-02, 1.65686e-03,
  198.         -5.53887e-06, 2.97810e-03, 0.00000e+00, 5.13122e-03, 8.66784e-02,
  199.         1.58727e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, -7.27026e-06,
  200.         0.00000e+00, 6.74494e+00, 4.93933e-03, 2.21656e-03, 2.50802e-03,
  201.         0.00000e+00, 0.00000e+00, -2.08841e-02, -1.79873e+00, 1.45103e-03,
  202.         2.81769e-04, -1.44703e-03, -5.16394e-05, 8.47001e-02, 1.70147e-01,
  203.         5.72562e-03, 5.07493e-05, 4.36148e-03, 1.17863e-04, 4.74364e-03,
  204.         6.61278e-03, 4.34292e-05, 1.44373e-03, 2.41470e-05, 2.84426e-03,
  205.         8.56560e-04, 2.04028e-03, 0.00000e+00, -3.15994e+03, -2.46423e-03,
  206.         1.13843e-03, 4.20512e-04, 0.00000e+00, -9.77214e+01, 6.77794e-03,
  207.         5.27499e-03, 1.14936e-03, 0.00000e+00, -6.61311e-03, -1.84255e-02,
  208.         -1.96259e-02, 2.98618e+04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  209.         6.44574e+02, 8.84668e-04, 5.05066e-04, 0.00000e+00, 4.02881e+03,
  210.         -1.89503e-03, 0.00000e+00, 0.00000e+00, 8.21407e-04, 2.06780e-03,
  211.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  212.         -1.20410e-02, -3.63963e-03, 9.92070e-05, -1.15284e-04, -6.33059e-05,
  213.         -6.05545e-01, 8.34218e-03, -9.13036e+01, 3.71042e-04, 0.00000e+00,
  214.         4.19000e-04, 2.70928e-03, 3.31507e-03, -4.44508e-03, -4.96334e-03,
  215.         -1.60449e-03, 3.95119e-03, 2.48924e-03, 5.09815e-04, 4.05302e-03,
  216.         2.24076e-03, 0.00000e+00, 6.84256e-03, 4.66354e-04, 0.00000e+00,
  217.         -3.68328e-04, 0.00000e+00, 0.00000e+00, -1.46870e+02, 0.00000e+00,
  218.         0.00000e+00, 1.09501e-03, 4.65156e-04, 5.62583e-04, 3.21596e+00,
  219.         6.43168e-04, 3.14860e-03, 3.40738e-03, 1.78481e-03, 9.62532e-04,
  220.         5.58171e-04, 3.43731e+00, -2.33195e-01, 5.10289e-04, 0.00000e+00,
  221.         0.00000e+00, -9.25347e+04, 0.00000e+00, -1.99639e-03, 0.00000e+00,
  222.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  223.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  224.     };

  225.     /** NRLMSISE-00 data: density pd[9][150]. */
  226.     private static final double[][] PD = {
  227.         // HE DENSITY
  228.         {
  229.             1.09979e+00, -4.88060e-02, -1.97501e-01, -9.10280e-02, -6.96558e-03,
  230.             2.42136e-02, 3.91333e-01, -7.20068e-03, -3.22718e-02, 1.41508e+00,
  231.             1.68194e-01, 1.85282e-02, 1.09384e-01, -7.24282e+00, 0.00000e+00,
  232.             2.96377e-01, -4.97210e-02, 1.04114e+02, -8.61108e-02, -7.29177e-04,
  233.             1.48998e-06, 1.08629e-03, 0.00000e+00, 0.00000e+00, 8.31090e-02,
  234.             1.12818e-01, -5.75005e-02, -1.29919e-02, -1.78849e-02, -2.86343e-06,
  235.             0.00000e+00, -1.51187e+02, -6.65902e-03, 0.00000e+00, -2.02069e-03,
  236.             0.00000e+00, 0.00000e+00, 4.32264e-02, -2.80444e+01, -3.26789e-03,
  237.             2.47461e-03, 0.00000e+00, 0.00000e+00, 9.82100e-02, 1.22714e-01,
  238.             -3.96450e-02, 0.00000e+00, -2.76489e-03, 0.00000e+00, 1.87723e-03,
  239.             -8.09813e-03, 4.34428e-05, -7.70932e-03, 0.00000e+00, -2.28894e-03,
  240.             -5.69070e-03, -5.22193e-03, 6.00692e-03, -7.80434e+03, -3.48336e-03,
  241.             -6.38362e-03, -1.82190e-03, 0.00000e+00, -7.58976e+01, -2.17875e-02,
  242.             -1.72524e-02, -9.06287e-03, 0.00000e+00, 2.44725e-02, 8.66040e-02,
  243.             1.05712e-01, 3.02543e+04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  244.             -6.01364e+03, -5.64668e-03, -2.54157e-03, 0.00000e+00, 3.15611e+02,
  245.             -5.69158e-03, 0.00000e+00, 0.00000e+00, -4.47216e-03, -4.49523e-03,
  246.             4.64428e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  247.             4.51236e-02, 2.46520e-02, 6.17794e-03, 0.00000e+00, 0.00000e+00,
  248.             -3.62944e-01, -4.80022e-02, -7.57230e+01, -1.99656e-03, 0.00000e+00,
  249.             -5.18780e-03, -1.73990e-02, -9.03485e-03, 7.48465e-03, 1.53267e-02,
  250.             1.06296e-02, 1.18655e-02, 2.55569e-03, 1.69020e-03, 3.51936e-02,
  251.             -1.81242e-02, 0.00000e+00, -1.00529e-01, -5.10574e-03, 0.00000e+00,
  252.             2.10228e-03, 0.00000e+00, 0.00000e+00, -1.73255e+02, 5.07833e-01,
  253.             -2.41408e-01, 8.75414e-03, 2.77527e-03, -8.90353e-05, -5.25148e+00,
  254.             -5.83899e-03, -2.09122e-02, -9.63530e-03, 9.77164e-03, 4.07051e-03,
  255.             2.53555e-04, -5.52875e+00, -3.55993e-01, -2.49231e-03, 0.00000e+00,
  256.             0.00000e+00, 2.86026e+01, 0.00000e+00, 3.42722e-04, 0.00000e+00,
  257.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  258.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  259.         },
  260.         // O DENSITY
  261.         {
  262.             1.02315e+00, -1.59710e-01, -1.06630e-01, -1.77074e-02, -4.42726e-03,
  263.             3.44803e-02, 4.45613e-02, -3.33751e-02, -5.73598e-02, 3.50360e-01,
  264.             6.33053e-02, 2.16221e-02, 5.42577e-02, -5.74193e+00, 0.00000e+00,
  265.             1.90891e-01, -1.39194e-02, 1.01102e+02, 8.16363e-02, 1.33717e-04,
  266.             6.54403e-06, 3.10295e-03, 0.00000e+00, 0.00000e+00, 5.38205e-02,
  267.             1.23910e-01, -1.39831e-02, 0.00000e+00, 0.00000e+00, -3.95915e-06,
  268.             0.00000e+00, -7.14651e-01, -5.01027e-03, 0.00000e+00, -3.24756e-03,
  269.             0.00000e+00, 0.00000e+00, 4.42173e-02, -1.31598e+01, -3.15626e-03,
  270.             1.24574e-03, -1.47626e-03, -1.55461e-03, 6.40682e-02, 1.34898e-01,
  271.             -2.42415e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00, 6.13666e-04,
  272.             -5.40373e-03, 2.61635e-05, -3.33012e-03, 0.00000e+00, -3.08101e-03,
  273.             -2.42679e-03, -3.36086e-03, 0.00000e+00, -1.18979e+03, -5.04738e-02,
  274.             -2.61547e-03, -1.03132e-03, 1.91583e-04, -8.38132e+01, -1.40517e-02,
  275.             -1.14167e-02, -4.08012e-03, 1.73522e-04, -1.39644e-02, -6.64128e-02,
  276.             -6.85152e-02, -1.34414e+04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  277.             6.07916e+02, -4.12220e-03, -2.20996e-03, 0.00000e+00, 1.70277e+03,
  278.             -4.63015e-03, 0.00000e+00, 0.00000e+00, -2.25360e-03, -2.96204e-03,
  279.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  280.             3.92786e-02, 1.31186e-02, -1.78086e-03, 0.00000e+00, 0.00000e+00,
  281.             -3.90083e-01, -2.84741e-02, -7.78400e+01, -1.02601e-03, 0.00000e+00,
  282.             -7.26485e-04, -5.42181e-03, -5.59305e-03, 1.22825e-02, 1.23868e-02,
  283.             6.68835e-03, -1.03303e-02, -9.51903e-03, 2.70021e-04, -2.57084e-02,
  284.             -1.32430e-02, 0.00000e+00, -3.81000e-02, -3.16810e-03, 0.00000e+00,
  285.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  286.             0.00000e+00, -9.05762e-04, -2.14590e-03, -1.17824e-03, 3.66732e+00,
  287.             -3.79729e-04, -6.13966e-03, -5.09082e-03, -1.96332e-03, -3.08280e-03,
  288.             -9.75222e-04, 4.03315e+00, -2.52710e-01, 0.00000e+00, 0.00000e+00,
  289.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  290.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  291.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  292.         },
  293.         // N2 DENSITY
  294.         {
  295.             1.16112e+00, 0.00000e+00, 0.00000e+00, 3.33725e-02, 0.00000e+00,
  296.             3.48637e-02, -5.44368e-03, 0.00000e+00, -6.73940e-02, 1.74754e-01,
  297.             0.00000e+00, 0.00000e+00, 0.00000e+00, 1.74712e+02, 0.00000e+00,
  298.             1.26733e-01, 0.00000e+00, 1.03154e+02, 5.52075e-02, 0.00000e+00,
  299.             0.00000e+00, 8.13525e-04, 0.00000e+00, 0.00000e+00, 8.66784e-02,
  300.             1.58727e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  301.             0.00000e+00, -2.50482e+01, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  302.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -2.48894e-03,
  303.             6.16053e-04, -5.79716e-04, 2.95482e-03, 8.47001e-02, 1.70147e-01,
  304.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  305.             0.00000e+00, 2.47425e-05, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  306.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  307.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  308.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  309.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  310.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  311.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  312.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  313.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  314.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  315.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  316.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  317.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  318.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  319.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  320.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  321.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  322.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  323.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  324.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  325.         },
  326.         // TOTAL MASS
  327.         {
  328.             9.44846e-01, 0.00000e+00, 0.00000e+00, -3.08617e-02, 0.00000e+00,
  329.             -2.44019e-02, 6.48607e-03, 0.00000e+00, 3.08181e-02, 4.59392e-02,
  330.             0.00000e+00, 0.00000e+00, 0.00000e+00, 1.74712e+02, 0.00000e+00,
  331.             2.13260e-02, 0.00000e+00, -3.56958e+02, 0.00000e+00, 1.82278e-04,
  332.             0.00000e+00, 3.07472e-04, 0.00000e+00, 0.00000e+00, 8.66784e-02,
  333.             1.58727e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  334.             0.00000e+00, 0.00000e+00, 3.83054e-03, 0.00000e+00, 0.00000e+00,
  335.             -1.93065e-03, -1.45090e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  336.             0.00000e+00, -1.23493e-03, 1.36736e-03, 8.47001e-02, 1.70147e-01,
  337.             3.71469e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  338.             5.10250e-03, 2.47425e-05, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  339.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  340.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  341.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  342.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  343.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  344.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  345.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  346.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  347.             0.00000e+00, 3.68756e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  348.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  349.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  350.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  351.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  352.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  353.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  354.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  355.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  356.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  357.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  358.         },
  359.         // O2 DENSITY
  360.         {
  361.             1.35580e+00, 1.44816e-01, 0.00000e+00, 6.07767e-02, 0.00000e+00,
  362.             2.94777e-02, 7.46900e-02, 0.00000e+00, -9.23822e-02, 8.57342e-02,
  363.             0.00000e+00, 0.00000e+00, 0.00000e+00, 2.38636e+01, 0.00000e+00,
  364.             7.71653e-02, 0.00000e+00, 8.18751e+01, 1.87736e-02, 0.00000e+00,
  365.             0.00000e+00, 1.49667e-02, 0.00000e+00, 0.00000e+00, 8.66784e-02,
  366.             1.58727e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  367.             0.00000e+00, -3.67874e+02, 5.48158e-03, 0.00000e+00, 0.00000e+00,
  368.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  369.             0.00000e+00, 0.00000e+00, 0.00000e+00, 8.47001e-02, 1.70147e-01,
  370.             1.22631e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  371.             8.17187e-03, 3.71617e-05, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  372.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  373.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -2.10826e-03,
  374.             -3.13640e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  375.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  376.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  377.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  378.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  379.             -7.35742e-02, -5.00266e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  380.             0.00000e+00, 1.94965e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  381.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  382.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  383.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  384.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  385.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  386.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  387.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  388.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  389.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  390.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  391.         },
  392.         // AR DENSITY
  393.         {
  394.             1.04761e+00, 2.00165e-01, 2.37697e-01, 3.68552e-02, 0.00000e+00,
  395.             3.57202e-02, -2.14075e-01, 0.00000e+00, -1.08018e-01, -3.73981e-01,
  396.             0.00000e+00, 3.10022e-02, -1.16305e-03, -2.07596e+01, 0.00000e+00,
  397.             8.64502e-02, 0.00000e+00, 9.74908e+01, 5.16707e-02, 0.00000e+00,
  398.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 8.66784e-02,
  399.             1.58727e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  400.             0.00000e+00, 3.46193e+02, 1.34297e-02, 0.00000e+00, 0.00000e+00,
  401.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -3.48509e-03,
  402.             -1.54689e-04, 0.00000e+00, 0.00000e+00, 8.47001e-02, 1.70147e-01,
  403.             1.47753e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  404.             1.89320e-02, 3.68181e-05, 1.32570e-02, 0.00000e+00, 0.00000e+00,
  405.             3.59719e-03, 7.44328e-03, -1.00023e-03, -6.50528e+03, 0.00000e+00,
  406.             1.03485e-02, -1.00983e-03, -4.06916e-03, -6.60864e+01, -1.71533e-02,
  407.             1.10605e-02, 1.20300e-02, -5.20034e-03, 0.00000e+00, 0.00000e+00,
  408.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  409.             -2.62769e+03, 7.13755e-03, 4.17999e-03, 0.00000e+00, 1.25910e+04,
  410.             0.00000e+00, 0.00000e+00, 0.00000e+00, -2.23595e-03, 4.60217e-03,
  411.             5.71794e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  412.             -3.18353e-02, -2.35526e-02, -1.36189e-02, 0.00000e+00, 0.00000e+00,
  413.             0.00000e+00, 2.03522e-02, -6.67837e+01, -1.09724e-03, 0.00000e+00,
  414.             -1.38821e-02, 1.60468e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  415.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 1.51574e-02,
  416.             -5.44470e-04, 0.00000e+00, 7.28224e-02, 6.59413e-02, 0.00000e+00,
  417.             -5.15692e-03, 0.00000e+00, 0.00000e+00, -3.70367e+03, 0.00000e+00,
  418.             0.00000e+00, 1.36131e-02, 5.38153e-03, 0.00000e+00, 4.76285e+00,
  419.             -1.75677e-02, 2.26301e-02, 0.00000e+00, 1.76631e-02, 4.77162e-03,
  420.             0.00000e+00, 5.39354e+00, 0.00000e+00, -7.51710e-03, 0.00000e+00,
  421.             0.00000e+00, -8.82736e+01, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  422.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  423.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  424.         },
  425.         // H DENSITY
  426.         {
  427.             1.26376e+00, -2.14304e-01, -1.49984e-01, 2.30404e-01, 2.98237e-02,
  428.             2.68673e-02, 2.96228e-01, 2.21900e-02, -2.07655e-02, 4.52506e-01,
  429.             1.20105e-01, 3.24420e-02, 4.24816e-02, -9.14313e+00, 0.00000e+00,
  430.             2.47178e-02, -2.88229e-02, 8.12805e+01, 5.10380e-02, -5.80611e-03,
  431.             2.51236e-05, -1.24083e-02, 0.00000e+00, 0.00000e+00, 8.66784e-02,
  432.             1.58727e-01, -3.48190e-02, 0.00000e+00, 0.00000e+00, 2.89885e-05,
  433.             0.00000e+00, 1.53595e+02, -1.68604e-02, 0.00000e+00, 1.01015e-02,
  434.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.84552e-04,
  435.             -1.22181e-03, 0.00000e+00, 0.00000e+00, 8.47001e-02, 1.70147e-01,
  436.             -1.04927e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00, -5.91313e-03,
  437.             -2.30501e-02, 3.14758e-05, 0.00000e+00, 0.00000e+00, 1.26956e-02,
  438.             8.35489e-03, 3.10513e-04, 0.00000e+00, 3.42119e+03, -2.45017e-03,
  439.             -4.27154e-04, 5.45152e-04, 1.89896e-03, 2.89121e+01, -6.49973e-03,
  440.             -1.93855e-02, -1.48492e-02, 0.00000e+00, -5.10576e-02, 7.87306e-02,
  441.             9.51981e-02, -1.49422e+04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  442.             2.65503e+02, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  443.             0.00000e+00, 0.00000e+00, 0.00000e+00, 6.37110e-03, 3.24789e-04,
  444.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  445.             6.14274e-02, 1.00376e-02, -8.41083e-04, 0.00000e+00, 0.00000e+00,
  446.             0.00000e+00, -1.27099e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  447.             -3.94077e-03, -1.28601e-02, -7.97616e-03, 0.00000e+00, 0.00000e+00,
  448.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  449.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  450.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  451.             0.00000e+00, -6.71465e-03, -1.69799e-03, 1.93772e-03, 3.81140e+00,
  452.             -7.79290e-03, -1.82589e-02, -1.25860e-02, -1.04311e-02, -3.02465e-03,
  453.             2.43063e-03, 3.63237e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  454.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  455.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  456.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  457.         },
  458.         // N DENSITY
  459.         {
  460.             7.09557e+01, -3.26740e-01, 0.00000e+00, -5.16829e-01, -1.71664e-03,
  461.             9.09310e-02, -6.71500e-01, -1.47771e-01, -9.27471e-02, -2.30862e-01,
  462.             -1.56410e-01, 1.34455e-02, -1.19717e-01, 2.52151e+00, 0.00000e+00,
  463.             -2.41582e-01, 5.92939e-02, 4.39756e+00, 9.15280e-02, 4.41292e-03,
  464.             0.00000e+00, 8.66807e-03, 0.00000e+00, 0.00000e+00, 8.66784e-02,
  465.             1.58727e-01, 9.74701e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  466.             0.00000e+00, 6.70217e+01, -1.31660e-03, 0.00000e+00, -1.65317e-02,
  467.             0.00000e+00, 0.00000e+00, 8.50247e-02, 2.77428e+01, 4.98658e-03,
  468.             6.15115e-03, 9.50156e-03, -2.12723e-02, 8.47001e-02, 1.70147e-01,
  469.             -2.38645e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00, 1.37380e-03,
  470.             -8.41918e-03, 2.80145e-05, 7.12383e-03, 0.00000e+00, -1.66209e-02,
  471.             1.03533e-04, -1.68898e-02, 0.00000e+00, 3.64526e+03, 0.00000e+00,
  472.             6.54077e-03, 3.69130e-04, 9.94419e-04, 8.42803e+01, -1.16124e-02,
  473.             -7.74414e-03, -1.68844e-03, 1.42809e-03, -1.92955e-03, 1.17225e-01,
  474.             -2.41512e-02, 1.50521e+04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  475.             1.60261e+03, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  476.             0.00000e+00, 0.00000e+00, 0.00000e+00, -3.54403e-04, -1.87270e-02,
  477.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  478.             2.76439e-02, 6.43207e-03, -3.54300e-02, 0.00000e+00, 0.00000e+00,
  479.             0.00000e+00, -2.80221e-02, 8.11228e+01, -6.75255e-04, 0.00000e+00,
  480.             -1.05162e-02, -3.48292e-03, -6.97321e-03, 0.00000e+00, 0.00000e+00,
  481.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  482.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  483.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  484.             0.00000e+00, -1.45546e-03, -1.31970e-02, -3.57751e-03, -1.09021e+00,
  485.             -1.50181e-02, -7.12841e-03, -6.64590e-03, -3.52610e-03, -1.87773e-02,
  486.             -2.22432e-03, -3.93895e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  487.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  488.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  489.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  490.         },
  491.         // HOT O DENSITY
  492.         {
  493.             6.04050e-02, 1.57034e+00, 2.99387e-02, 0.00000e+00, 0.00000e+00,
  494.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -1.51018e+00,
  495.             0.00000e+00, 0.00000e+00, 0.00000e+00, -8.61650e+00, 1.26454e-02,
  496.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  497.             0.00000e+00, 5.50878e-03, 0.00000e+00, 0.00000e+00, 8.66784e-02,
  498.             1.58727e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  499.             0.00000e+00, 0.00000e+00, 6.23881e-02, 0.00000e+00, 0.00000e+00,
  500.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  501.             0.00000e+00, 0.00000e+00, 0.00000e+00, 8.47001e-02, 1.70147e-01,
  502.             -9.45934e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  503.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  504.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  505.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  506.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  507.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  508.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  509.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  510.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  511.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  512.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  513.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  514.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  515.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  516.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  517.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  518.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  519.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  520.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  521.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  522.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  523.         }
  524.     };

  525.     /** NRLMSISE-00 data: ps[150]. */
  526.     private static final double[] PS = {
  527.         9.56827e-01, 6.20637e-02, 3.18433e-02, 0.00000e+00, 0.00000e+00,
  528.         3.94900e-02, 0.00000e+00, 0.00000e+00, -9.24882e-03, -7.94023e-03,
  529.         0.00000e+00, 0.00000e+00, 0.00000e+00, 1.74712e+02, 0.00000e+00,
  530.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  531.         0.00000e+00, 2.74677e-03, 0.00000e+00, 1.54951e-02, 8.66784e-02,
  532.         1.58727e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  533.         0.00000e+00, 0.00000e+00, 0.00000e+00, -6.99007e-04, 0.00000e+00,
  534.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  535.         0.00000e+00, 1.24362e-02, -5.28756e-03, 8.47001e-02, 1.70147e-01,
  536.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  537.         0.00000e+00, 2.47425e-05, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  538.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  539.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  540.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  541.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  542.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  543.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  544.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  545.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  546.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  547.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  548.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  549.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  550.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  551.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  552.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  553.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  554.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  555.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  556.         0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  557.     };

  558.     /** NRLMSISE-00 data: TURBO pdl[2][25]. */
  559.     private static final double[][] PDL = {
  560.         {
  561.             1.09930e+00, 3.90631e+00, 3.07165e+00, 9.86161e-01, 1.63536e+01,
  562.             4.63830e+00, 1.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  563.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  564.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  565.             0.00000e+00, 0.00000e+00, 1.28840e+00, 3.10302e-02, 1.18339e-01
  566.         },
  567.         {
  568.             1.00000e+00, 7.00000e-01, 1.15020e+00, 3.44689e+00, 1.28840e+00,
  569.             1.00000e+00, 1.08738e+00, 1.22947e+00, 1.10016e+00, 7.34129e-01,
  570.             1.15241e+00, 2.22784e+00, 7.95046e-01, 4.01612e+00, 4.47749e+00,
  571.             1.23435e+02, -7.60535e-02, 1.68986e-06, 7.44294e-01, 1.03604e+00,
  572.             1.72783e+02, 1.15020e+00, 3.44689e+00, -7.46230e-01, 9.49154e-01
  573.         }
  574.     };

  575.     /** NRLMSISE-00 data: LOWER BOUNDARY ptm[10]. */
  576.     private static final double[] PTM = {
  577.         1.04130e+03, 3.86000e+02, 1.95000e+02, 1.66728e+01, 2.13000e+02,
  578.         1.20000e+02, 2.40000e+02, 1.87000e+02, -2.00000e+00, 0.00000e+00
  579.     };

  580.     /** NRLMSISE-00 data: pdm[8][10]. */
  581.     private static final double[][] PDM = {
  582.         {
  583.             2.45600e+07, 6.71072e-06, 1.00000e+02, 0.00000e+00, 1.10000e+02,
  584.             1.00000e+01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  585.         },
  586.         {
  587.             8.59400E+10, 1.00000e+00, 1.05000e+02, -8.00000e+00, 1.10000e+02,
  588.             1.00000e+01, 9.00000e+01, 2.00000e+00, 0.00000e+00, 0.00000e+00
  589.         },
  590.         {
  591.             2.81000E+11, 0.00000e+00, 1.05000e+02, 2.80000e+01, 2.89500e+01,
  592.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  593.         },
  594.         {
  595.             3.30000E+10, 2.68270e-01, 1.05000e+02, 1.00000e+00, 1.10000e+02,
  596.             1.00000e+01, 1.10000e+02, -1.00000e+01, 0.00000e+00, 0.00000e+00
  597.         },
  598.         {
  599.             1.33000e+09, 1.19615e-02, 1.05000e+02, 0.00000e+00, 1.10000e+02,
  600.             1.00000e+01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00
  601.         },
  602.         {
  603.             1.76100e+05, 1.00000e+00, 9.50000e+01, -8.00000e+00, 1.10000e+02,
  604.             1.00000e+01, 9.00000e+01, 2.00000e+00, 0.00000e+00, 0.00000e+00,
  605.         },
  606.         {
  607.             1.00000e+07, 1.00000e+00, 1.05000e+02, -8.00000e+00, 1.10000e+02,
  608.             1.00000e+01, 9.00000e+01, 2.00000e+00, 0.00000e+00, 0.00000e+00
  609.         },
  610.         {
  611.             1.00000e+06, 1.00000e+00, 1.05000e+02, -8.00000e+00, 5.50000e+02,
  612.             7.60000e+01, 9.00000e+01, 2.00000e+00, 0.00000e+00, 4.00000e+03
  613.         }
  614.     };

  615.     /** NRLMSISE-00 data: ptl[4][100]. */
  616.     private static final double[][] PTL = {
  617.         // TN1(2)
  618.         {
  619.             1.00858e+00, 4.56011e-02, -2.22972e-02, -5.44388e-02, 5.23136e-04,
  620.             -1.88849e-02, 5.23707e-02, -9.43646e-03, 6.31707e-03, -7.80460e-02,
  621.             -4.88430e-02, 0.00000e+00, 0.00000e+00, -7.60250e+00, 0.00000e+00,
  622.             -1.44635e-02, -1.76843e-02, -1.21517e+02, 2.85647e-02, 0.00000e+00,
  623.             0.00000e+00, 6.31792e-04, 0.00000e+00, 5.77197e-03, 8.66784e-02,
  624.             1.58727e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  625.             0.00000e+00, -8.90272e+03, 3.30611e-03, 3.02172e-03, 0.00000e+00,
  626.             -2.13673e-03, -3.20910e-04, 0.00000e+00, 0.00000e+00, 2.76034e-03,
  627.             2.82487e-03, -2.97592e-04, -4.21534e-03, 8.47001e-02, 1.70147e-01,
  628.             8.96456e-03, 0.00000e+00, -1.08596e-02, 0.00000e+00, 0.00000e+00,
  629.             5.57917e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  630.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  631.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  632.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  633.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  634.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  635.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  636.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  637.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  638.             0.00000e+00, 9.65405e-03, 0.00000e+00, 0.00000e+00, 2.00000e+00
  639.         },
  640.         // TN1(3)
  641.         {
  642.             9.39664e-01, 8.56514e-02, -6.79989e-03, 2.65929e-02, -4.74283e-03,
  643.             1.21855e-02, -2.14905e-02, 6.49651e-03, -2.05477e-02, -4.24952e-02,
  644.             0.00000e+00, 0.00000e+00, 0.00000e+00, 1.19148e+01, 0.00000e+00,
  645.             1.18777e-02, -7.28230e-02, -8.15965e+01, 1.73887e-02, 0.00000e+00,
  646.             0.00000e+00, 0.00000e+00, -1.44691e-02, 2.80259e-04, 8.66784e-02,
  647.             1.58727e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  648.             0.00000e+00, 2.16584e+02, 3.18713e-03, 7.37479e-03, 0.00000e+00,
  649.             -2.55018e-03, -3.92806e-03, 0.00000e+00, 0.00000e+00, -2.89757e-03,
  650.             -1.33549e-03, 1.02661e-03, 3.53775e-04, 8.47001e-02, 1.70147e-01,
  651.             -9.17497e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  652.             3.56082e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  653.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  654.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  655.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  656.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  657.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  658.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  659.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  660.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  661.             0.00000e+00, -1.00902e-02, 0.00000e+00, 0.00000e+00, 2.00000e+00
  662.         },
  663.         // TN1(4)
  664.         {
  665.             9.85982e-01, -4.55435e-02, 1.21106e-02, 2.04127e-02, -2.40836e-03,
  666.             1.11383e-02, -4.51926e-02, 1.35074e-02, -6.54139e-03, 1.15275e-01,
  667.             1.28247e-01, 0.00000e+00, 0.00000e+00, -5.30705e+00, 0.00000e+00,
  668.             -3.79332e-02, -6.24741e-02, 7.71062e-01, 2.96315e-02, 0.00000e+00,
  669.             0.00000e+00, 0.00000e+00, 6.81051e-03, -4.34767e-03, 8.66784e-02,
  670.             1.58727e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  671.             0.00000e+00, 1.07003e+01, -2.76907e-03, 4.32474e-04, 0.00000e+00,
  672.             1.31497e-03, -6.47517e-04, 0.00000e+00, -2.20621e+01, -1.10804e-03,
  673.             -8.09338e-04, 4.18184e-04, 4.29650e-03, 8.47001e-02, 1.70147e-01,
  674.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  675.             -4.04337e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  676.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  677.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -9.52550e-04,
  678.             8.56253e-04, 4.33114e-04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  679.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 1.21223e-03,
  680.             2.38694e-04, 9.15245e-04, 1.28385e-03, 8.67668e-04, -5.61425e-06,
  681.             1.04445e+00, 3.41112e+01, 0.00000e+00, -8.40704e-01, -2.39639e+02,
  682.             7.06668e-01, -2.05873e+01, -3.63696e-01, 2.39245e+01, 0.00000e+00,
  683.             -1.06657e-03, -7.67292e-04, 1.54534e-04, 0.00000e+00, 0.00000e+00,
  684.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  685.         },
  686.         // TN1(5) TN2(1)
  687.         {
  688.             1.00320e+00, 3.83501e-02, -2.38983e-03, 2.83950e-03, 4.20956e-03,
  689.             5.86619e-04, 2.19054e-02, -1.00946e-02, -3.50259e-03, 4.17392e-02,
  690.             -8.44404e-03, 0.00000e+00, 0.00000e+00, 4.96949e+00, 0.00000e+00,
  691.             -7.06478e-03, -1.46494e-02, 3.13258e+01, -1.86493e-03, 0.00000e+00,
  692.             -1.67499e-02, 0.00000e+00, 0.00000e+00, 5.12686e-04, 8.66784e-02,
  693.             1.58727e-01, -4.64167e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  694.             4.37353e-03, -1.99069e+02, 0.00000e+00, -5.34884e-03, 0.00000e+00,
  695.             1.62458e-03, 2.93016e-03, 2.67926e-03, 5.90449e+02, 0.00000e+00,
  696.             0.00000e+00, -1.17266e-03, -3.58890e-04, 8.47001e-02, 1.70147e-01,
  697.             0.00000e+00, 0.00000e+00, 1.38673e-02, 0.00000e+00, 0.00000e+00,
  698.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  699.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  700.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 1.60571e-03,
  701.             6.28078e-04, 5.05469e-05, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  702.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -1.57829e-03,
  703.             -4.00855e-04, 5.04077e-05, -1.39001e-03, -2.33406e-03, -4.81197e-04,
  704.             1.46758e+00, 6.20332e+00, 0.00000e+00, 3.66476e-01, -6.19760e+01,
  705.             3.09198e-01, -1.98999e+01, 0.00000e+00, -3.29933e+02, 0.00000e+00,
  706.             -1.10080e-03, -9.39310e-05, 1.39638e-04, 0.00000e+00, 0.00000e+00,
  707.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  708.         }
  709.     };

  710.     /** NRLMSISE-00 data: pma[10][100]. */
  711.     private static final double[][] PMA = {
  712.         // TN2(2)
  713.         {
  714.             9.81637e-01, -1.41317e-03, 3.87323e-02, 0.00000e+00, 0.00000e+00,
  715.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -3.58707e-02,
  716.             -8.63658e-03, 0.00000e+00, 0.00000e+00, -2.02226e+00, 0.00000e+00,
  717.             -8.69424e-03, -1.91397e-02, 8.76779e+01, 4.52188e-03, 0.00000e+00,
  718.             2.23760e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  719.             0.00000e+00, -7.07572e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  720.             -4.11210e-03, 3.50060e+01, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  721.             0.00000e+00, 0.00000e+00, -8.36657e-03, 1.61347e+01, 0.00000e+00,
  722.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  723.             0.00000e+00, 0.00000e+00, -1.45130e-02, 0.00000e+00, 0.00000e+00,
  724.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  725.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  726.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 1.24152e-03,
  727.             6.43365e-04, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  728.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 1.33255e-03,
  729.             2.42657e-03, 1.60666e-03, -1.85728e-03, -1.46874e-03, -4.79163e-06,
  730.             1.22464e+00, 3.53510e+01, 0.00000e+00, 4.49223e-01, -4.77466e+01,
  731.             4.70681e-01, 8.41861e+00, -2.88198e-01, 1.67854e+02, 0.00000e+00,
  732.             7.11493e-04, 6.05601e-04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  733.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  734.         },
  735.         // TN2(3)
  736.         {
  737.             1.00422e+00, -7.11212e-03, 5.24480e-03, 0.00000e+00, 0.00000e+00,
  738.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -5.28914e-02,
  739.             -2.41301e-02, 0.00000e+00, 0.00000e+00, -2.12219e+01, -1.03830e-02,
  740.             -3.28077e-03, 1.65727e-02, 1.68564e+00, -6.68154e-03, 0.00000e+00,
  741.             1.45155e-02, 0.00000e+00, 8.42365e-03, 0.00000e+00, 0.00000e+00,
  742.             0.00000e+00, -4.34645e-03, 0.00000e+00, 0.00000e+00, 2.16780e-02,
  743.             0.00000e+00, -1.38459e+02, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  744.             0.00000e+00, 0.00000e+00, 7.04573e-03, -4.73204e+01, 0.00000e+00,
  745.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  746.             0.00000e+00, 0.00000e+00, 1.08767e-02, 0.00000e+00, 0.00000e+00,
  747.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  748.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -8.08279e-03,
  749.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 5.21769e-04,
  750.             -2.27387e-04, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  751.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 3.26769e-03,
  752.             3.16901e-03, 4.60316e-04, -1.01431e-04, 1.02131e-03, 9.96601e-04,
  753.             1.25707e+00, 2.50114e+01, 0.00000e+00, 4.24472e-01, -2.77655e+01,
  754.             3.44625e-01, 2.75412e+01, 0.00000e+00, 7.94251e+02, 0.00000e+00,
  755.             2.45835e-03, 1.38871e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  756.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  757.         },
  758.         // TN2(4) TN3(1)
  759.         {
  760.             1.01890e+00, -2.46603e-02, 1.00078e-02, 0.00000e+00, 0.00000e+00,
  761.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -6.70977e-02,
  762.             -4.02286e-02, 0.00000e+00, 0.00000e+00, -2.29466e+01, -7.47019e-03,
  763.             2.26580e-03, 2.63931e-02, 3.72625e+01, -6.39041e-03, 0.00000e+00,
  764.             9.58383e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  765.             0.00000e+00, -1.85291e-03, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  766.             0.00000e+00, 1.39717e+02, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  767.             0.00000e+00, 0.00000e+00, 9.19771e-03, -3.69121e+02, 0.00000e+00,
  768.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  769.             0.00000e+00, 0.00000e+00, -1.57067e-02, 0.00000e+00, 0.00000e+00,
  770.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  771.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -7.07265e-03,
  772.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -2.92953e-03,
  773.             -2.77739e-03, -4.40092e-04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  774.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.47280e-03,
  775.             2.95035e-04, -1.81246e-03, 2.81945e-03, 4.27296e-03, 9.78863e-04,
  776.             1.40545e+00, -6.19173e+00, 0.00000e+00, 0.00000e+00, -7.93632e+01,
  777.             4.44643e-01, -4.03085e+02, 0.00000e+00, 1.15603e+01, 0.00000e+00,
  778.             2.25068e-03, 8.48557e-04, -2.98493e-04, 0.00000e+00, 0.00000e+00,
  779.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  780.         },
  781.         // TN3(2)
  782.         {
  783.             9.75801e-01, 3.80680e-02, -3.05198e-02, 0.00000e+00, 0.00000e+00,
  784.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 3.85575e-02,
  785.             5.04057e-02, 0.00000e+00, 0.00000e+00, -1.76046e+02, 1.44594e-02,
  786.             -1.48297e-03, -3.68560e-03, 3.02185e+01, -3.23338e-03, 0.00000e+00,
  787.             1.53569e-02, 0.00000e+00, -1.15558e-02, 0.00000e+00, 0.00000e+00,
  788.             0.00000e+00, 4.89620e-03, 0.00000e+00, 0.00000e+00, -1.00616e-02,
  789.             -8.21324e-03, -1.57757e+02, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  790.             0.00000e+00, 0.00000e+00, 6.63564e-03, 4.58410e+01, 0.00000e+00,
  791.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  792.             0.00000e+00, 0.00000e+00, -2.51280e-02, 0.00000e+00, 0.00000e+00,
  793.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  794.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 9.91215e-03,
  795.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -8.73148e-04,
  796.             -1.29648e-03, -7.32026e-05, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  797.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -4.68110e-03,
  798.             -4.66003e-03, -1.31567e-03, -7.39390e-04, 6.32499e-04, -4.65588e-04,
  799.             -1.29785e+00, -1.57139e+02, 0.00000e+00, 2.58350e-01, -3.69453e+01,
  800.             4.10672e-01, 9.78196e+00, -1.52064e-01, -3.85084e+03, 0.00000e+00,
  801.             -8.52706e-04, -1.40945e-03, -7.26786e-04, 0.00000e+00, 0.00000e+00,
  802.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  803.         },
  804.         // TN3(3)
  805.         {
  806.             9.60722e-01, 7.03757e-02, -3.00266e-02, 0.00000e+00, 0.00000e+00,
  807.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.22671e-02,
  808.             4.10423e-02, 0.00000e+00, 0.00000e+00, -1.63070e+02, 1.06073e-02,
  809.             5.40747e-04, 7.79481e-03, 1.44908e+02, 1.51484e-04, 0.00000e+00,
  810.             1.97547e-02, 0.00000e+00, -1.41844e-02, 0.00000e+00, 0.00000e+00,
  811.             0.00000e+00, 5.77884e-03, 0.00000e+00, 0.00000e+00, 9.74319e-03,
  812.             0.00000e+00, -2.88015e+03, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  813.             0.00000e+00, 0.00000e+00, -4.44902e-03, -2.92760e+01, 0.00000e+00,
  814.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  815.             0.00000e+00, 0.00000e+00, 2.34419e-02, 0.00000e+00, 0.00000e+00,
  816.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  817.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 5.36685e-03,
  818.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -4.65325e-04,
  819.             -5.50628e-04, 3.31465e-04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  820.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -2.06179e-03,
  821.             -3.08575e-03, -7.93589e-04, -1.08629e-04, 5.95511e-04, -9.05050e-04,
  822.             1.18997e+00, 4.15924e+01, 0.00000e+00, -4.72064e-01, -9.47150e+02,
  823.             3.98723e-01, 1.98304e+01, 0.00000e+00, 3.73219e+03, 0.00000e+00,
  824.             -1.50040e-03, -1.14933e-03, -1.56769e-04, 0.00000e+00, 0.00000e+00,
  825.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  826.         },
  827.         // TN3(4)
  828.         {
  829.             1.03123e+00, -7.05124e-02, 8.71615e-03, 0.00000e+00, 0.00000e+00,
  830.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -3.82621e-02,
  831.             -9.80975e-03, 0.00000e+00, 0.00000e+00, 2.89286e+01, 9.57341e-03,
  832.             0.00000e+00, 0.00000e+00, 8.66153e+01, 7.91938e-04, 0.00000e+00,
  833.             0.00000e+00, 0.00000e+00, 4.68917e-03, 0.00000e+00, 0.00000e+00,
  834.             0.00000e+00, 7.86638e-03, 0.00000e+00, 0.00000e+00, 9.90827e-03,
  835.             0.00000e+00, 6.55573e+01, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  836.             0.00000e+00, 0.00000e+00, 0.00000e+00, -4.00200e+01, 0.00000e+00,
  837.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  838.             0.00000e+00, 0.00000e+00, 7.07457e-03, 0.00000e+00, 0.00000e+00,
  839.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  840.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 5.72268e-03,
  841.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -2.04970e-04,
  842.             1.21560e-03, -8.05579e-06, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  843.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -2.49941e-03,
  844.             -4.57256e-04, -1.59311e-04, 2.96481e-04, -1.77318e-03, -6.37918e-04,
  845.             1.02395e+00, 1.28172e+01, 0.00000e+00, 1.49903e-01, -2.63818e+01,
  846.             0.00000e+00, 4.70628e+01, -2.22139e-01, 4.82292e-02, 0.00000e+00,
  847.             -8.67075e-04, -5.86479e-04, 5.32462e-04, 0.00000e+00, 0.00000e+00,
  848.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  849.         },
  850.         // TN3(5) SURFACE TEMP TSL
  851.         {
  852.             1.00828e+00, -9.10404e-02, -2.26549e-02, 0.00000e+00, 0.00000e+00,
  853.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -2.32420e-02,
  854.             -9.08925e-03, 0.00000e+00, 0.00000e+00, 3.36105e+01, 0.00000e+00,
  855.             0.00000e+00, 0.00000e+00, -1.24957e+01, -5.87939e-03, 0.00000e+00,
  856.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  857.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  858.             0.00000e+00, 2.79765e+01, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  859.             0.00000e+00, 0.00000e+00, 0.00000e+00, 2.01237e+03, 0.00000e+00,
  860.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  861.             0.00000e+00, 0.00000e+00, -1.75553e-02, 0.00000e+00, 0.00000e+00,
  862.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  863.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  864.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 3.29699e-03,
  865.             1.26659e-03, 2.68402e-04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  866.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 1.17894e-03,
  867.             1.48746e-03, 1.06478e-04, 1.34743e-04, -2.20939e-03, -6.23523e-04,
  868.             6.36539e-01, 1.13621e+01, 0.00000e+00, -3.93777e-01, 2.38687e+03,
  869.             0.00000e+00, 6.61865e+02, -1.21434e-01, 9.27608e+00, 0.00000e+00,
  870.             1.68478e-04, 1.24892e-03, 1.71345e-03, 0.00000e+00, 0.00000e+00,
  871.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  872.         },
  873.         // TGN3(2) SURFACE GRAD TSLG
  874.         {
  875.             1.57293e+00, -6.78400e-01, 6.47500e-01, 0.00000e+00, 0.00000e+00,
  876.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -7.62974e-02,
  877.             -3.60423e-01, 0.00000e+00, 0.00000e+00, 1.28358e+02, 0.00000e+00,
  878.             0.00000e+00, 0.00000e+00, 4.68038e+01, 0.00000e+00, 0.00000e+00,
  879.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  880.             0.00000e+00, -1.67898e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  881.             0.00000e+00, 2.90994e+04, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  882.             0.00000e+00, 0.00000e+00, 0.00000e+00, 3.15706e+01, 0.00000e+00,
  883.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  884.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  885.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  886.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  887.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  888.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  889.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  890.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  891.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  892.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  893.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  894.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  895.         },
  896.         // TGN2(1) TGN1(2)
  897.         {
  898.             8.60028e-01, 3.77052e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  899.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -1.17570e+00,
  900.             0.00000e+00, 0.00000e+00, 0.00000e+00, 7.77757e-03, 0.00000e+00,
  901.             0.00000e+00, 0.00000e+00, 1.01024e+02, 0.00000e+00, 0.00000e+00,
  902.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  903.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  904.             0.00000e+00, 6.54251e+02, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  905.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  906.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  907.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  908.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  909.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  910.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  911.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  912.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, -1.56959e-02,
  913.             1.91001e-02, 3.15971e-02, 1.00982e-02, -6.71565e-03, 2.57693e-03,
  914.             1.38692e+00, 2.82132e-01, 0.00000e+00, 0.00000e+00, 3.81511e+02,
  915.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  916.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  917.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  918.         },
  919.         // TGN3(1) TGN2(2)
  920.         {
  921.             1.06029e+00, -5.25231e-02, 3.73034e-01, 0.00000e+00, 0.00000e+00,
  922.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 3.31072e-02,
  923.             -3.88409e-01, 0.00000e+00, 0.00000e+00, -1.65295e+02, -2.13801e-01,
  924.             -4.38916e-02, -3.22716e-01, -8.82393e+01, 1.18458e-01, 0.00000e+00,
  925.             -4.35863e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  926.             0.00000e+00, -1.19782e-01, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  927.             0.00000e+00, 2.62229e+01, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  928.             0.00000e+00, 0.00000e+00, 0.00000e+00, -5.37443e+01, 0.00000e+00,
  929.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  930.             0.00000e+00, 0.00000e+00, -4.55788e-01, 0.00000e+00, 0.00000e+00,
  931.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  932.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  933.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 3.84009e-02,
  934.             3.96733e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  935.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 5.05494e-02,
  936.             7.39617e-02, 1.92200e-02, -8.46151e-03, -1.34244e-02, 1.96338e-02,
  937.             1.50421e+00, 1.88368e+01, 0.00000e+00, 0.00000e+00, -5.13114e+01,
  938.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  939.             5.11923e-02, 3.61225e-02, 0.00000e+00, 0.00000e+00, 0.00000e+00,
  940.             0.00000e+00, 0.00000e+00, 0.00000e+00, 0.00000e+00, 2.00000e+00
  941.         }
  942.     };

  943.     /**  NRLMSISE-00 data: MIDDLE ATMOSPHERE AVERAGES pavgm[10]. */
  944.     private static final double[] PAVGM = {
  945.         2.61000e+02, 2.64000e+02, 2.29000e+02, 2.17000e+02, 2.17000e+02,
  946.         2.23000e+02, 2.86760e+02, -2.93940e+00, 2.50000e+00, 0.00000e+00
  947.     };

  948.     /** NRLMSISE-00 minimum temperature, used in many cases in density computation. */
  949.     private static final double MIN_TEMP = 50.;

  950.     // Fields

  951.     /** External data container. */
  952.     private final NRLMSISE00InputParameters inputParams;

  953.     /** Earth body shape. */
  954.     private final BodyShape earth;

  955.     /** Switches for main effects. */
  956.     private final int[] sw;

  957.     /** Switches for cross effects. */
  958.     private final int[] swc;

  959.     /** UT time scale. */
  960.     private final TimeScale ut;

  961.     /** Constructor.
  962.      * <p>
  963.      * The model is constructed with all switches set to 1.
  964.      * </p>
  965.      * <p>
  966.      * Parameters are mandatory only for the
  967.      * {@link #getDensity(AbsoluteDate, Vector3D, Frame) getDensity()} and
  968.      * {@link #getVelocity(AbsoluteDate, Vector3D, Frame) getVelocity()} methods.
  969.      * </p>
  970.      *
  971.      * <p>This constructor uses the {@link DataContext#getDefault() default data context}.
  972.      *
  973.      * @param parameters the solar and magnetic activity data
  974.      * @param sun the Sun position
  975.      * @param earth the Earth body shape
  976.      * @see #NRLMSISE00(NRLMSISE00InputParameters, ExtendedPositionProvider, BodyShape,
  977.      * TimeScale)
  978.      */
  979.     @DefaultDataContext
  980.     public NRLMSISE00(final NRLMSISE00InputParameters parameters,
  981.                       final ExtendedPositionProvider sun,
  982.                       final BodyShape earth) {
  983.         this(parameters, sun, earth,
  984.                 DataContext.getDefault().getTimeScales()
  985.                         .getUT1(IERSConventions.IERS_2010, true));
  986.     }

  987.     /** Constructor.
  988.      * <p>
  989.      * The model is constructed with all switches set to 1.
  990.      * </p>
  991.      * <p>
  992.      * Parameters are mandatory only for the
  993.      * {@link #getDensity(AbsoluteDate, Vector3D, Frame) getDensity()} and
  994.      * {@link #getVelocity(AbsoluteDate, Vector3D, Frame) getVelocity()} methods.
  995.      * </p>
  996.      * @param parameters the solar and magnetic activity data
  997.      * @param sun the Sun position
  998.      * @param earth the Earth body shape
  999.      * @param ut UT time scale. The original documentation for NRLMSISE00 does not
  1000.      *           distinguish between UTC and UT1. In Orekit 10.0 {@code
  1001.      *           TimeScalesFactory.getUT1(IERSConventions.IERS_2010, true)} was used.
  1002.      * @since 10.1
  1003.      */
  1004.     public NRLMSISE00(final NRLMSISE00InputParameters parameters,
  1005.                       final ExtendedPositionProvider sun,
  1006.                       final BodyShape earth,
  1007.                       final TimeScale ut) {
  1008.         this(parameters, sun, earth, allOnes(), allOnes(), ut);
  1009.     }

  1010.     /** Constructor.
  1011.      * <p>
  1012.      * The model is constructed with all switches set to 1.
  1013.      * </p>
  1014.      * <p>
  1015.      * Parameters are mandatory only for the
  1016.      * {@link #getDensity(AbsoluteDate, Vector3D, Frame) getDensity()} and
  1017.      * {@link #getVelocity(AbsoluteDate, Vector3D, Frame) getVelocity()} methods.
  1018.      * </p>
  1019.      * @param parameters the solar and magnetic activity data
  1020.      * @param sun the Sun position
  1021.      * @param earth the Earth body shape
  1022.      * @param sw switches for main effects
  1023.      * @param swc switches for cross effects
  1024.      * @param ut UT time scale.
  1025.      */
  1026.     private NRLMSISE00(final NRLMSISE00InputParameters parameters,
  1027.                        final ExtendedPositionProvider sun,
  1028.                        final BodyShape earth,
  1029.                        final int[] sw,
  1030.                        final int[] swc,
  1031.                        final TimeScale ut) {
  1032.         super(sun);
  1033.         this.inputParams = parameters;
  1034.         this.earth       = earth;
  1035.         this.sw          = sw;
  1036.         this.swc         = swc;
  1037.         this.ut = ut;
  1038.     }

  1039.     /** Change a switch.
  1040.      * <p>
  1041.      * This method creates a new instance, the current instance is
  1042.      * not changed at all!
  1043.      * </p>
  1044.      * @param number switch number between 1 and 23
  1045.      * @param value switch value
  1046.      * @return a <em>new</em> instance, with switch changed
  1047.      */
  1048.     public NRLMSISE00 withSwitch(final int number, final int value) {
  1049.         if (number < 1 || number > 23) {
  1050.             throw new OrekitException(LocalizedCoreFormats.OUT_OF_RANGE_SIMPLE, number, 1, 23);
  1051.         }

  1052.         final int[] newSw       = sw.clone();
  1053.         final int[] newSwc      = swc.clone();
  1054.         if (number != 9) {
  1055.             newSw[number]  = (value == 1) ? 1 : 0;
  1056.             newSwc[number] = (value > 0) ? 1 : 0;
  1057.         } else {
  1058.             if (value == -1 || value == 1) {
  1059.                 newSw[number] = value;
  1060.             } else {
  1061.                 newSw[number] = 0;
  1062.             }
  1063.             newSwc[number] = newSw[number];
  1064.         }

  1065.         return new NRLMSISE00(inputParams, getSun(), earth, newSwc, newSwc, ut);

  1066.     }

  1067.     /** Create an array of switches set to 1.
  1068.      * @return array of switches
  1069.      */
  1070.     private static int[] allOnes() {
  1071.         final int[] array = new int[24];
  1072.         Arrays.fill(array, 1);
  1073.         return array;
  1074.     }

  1075.     /** {@inheritDoc} */
  1076.     @Override
  1077.     public Frame getFrame() {
  1078.         return earth.getBodyFrame();
  1079.     }

  1080.     /** {@inheritDoc} */
  1081.     @Override
  1082.     public double getDensity(final AbsoluteDate date,
  1083.                              final Vector3D position,
  1084.                              final Frame frame) {

  1085.         // check if data are available :
  1086.         if (!date.isBetweenOrEqualTo(inputParams.getMinDate(), inputParams.getMaxDate())) {
  1087.             throw new OrekitException(OrekitMessages.NO_SOLAR_ACTIVITY_AT_DATE,
  1088.                                       date, inputParams.getMinDate(), inputParams.getMaxDate());
  1089.         }

  1090.         // compute day number in current year and the seconds within the day
  1091.         final DateTimeComponents dtc = date.getComponents(ut);
  1092.         final int    doy = dtc.getDate().getDayOfYear();
  1093.         final double sec = dtc.getTime().getSecondsInLocalDay();

  1094.         // compute geodetic position (km and °)
  1095.         final GeodeticPoint inBody = earth.transform(position, frame, date);
  1096.         final double alt = inBody.getAltitude() / 1000.;
  1097.         final double lon = FastMath.toDegrees(inBody.getLongitude());
  1098.         final double lat = FastMath.toDegrees(inBody.getLatitude());

  1099.         // compute local solar time
  1100.         final double lst = localSolarTime(date, position, frame);

  1101.         // get solar activity data and compute
  1102.         final Output out = new Output(doy, sec, lat, lon, lst, inputParams.getAverageFlux(date),
  1103.                                       inputParams.getDailyFlux(date), inputParams.getAp(date));
  1104.         out.gtd7d(alt);

  1105.         // return the local density
  1106.         return out.getDensity(TOTAL_MASS);

  1107.     }

  1108.     /** {@inheritDoc} */
  1109.     @Override
  1110.     public <T extends CalculusFieldElement<T>> T getDensity(final FieldAbsoluteDate<T> date,
  1111.                                                         final FieldVector3D<T> position,
  1112.                                                         final Frame frame) {
  1113.         // check if data are available :
  1114.         final AbsoluteDate dateD = date.toAbsoluteDate();
  1115.         if (!dateD.isBetweenOrEqualTo(inputParams.getMinDate(), inputParams.getMaxDate())) {
  1116.             throw new OrekitException(OrekitMessages.NO_SOLAR_ACTIVITY_AT_DATE,
  1117.                                       dateD, inputParams.getMinDate(), inputParams.getMaxDate());
  1118.         }

  1119.         // compute day number in current year and the seconds within the day
  1120.         final DateTimeComponents dtc = dateD.getComponents(ut);
  1121.         final int    doy = dtc.getDate().getDayOfYear();
  1122.         final T sec = date.durationFrom(new AbsoluteDate(dtc.getDate(), TimeComponents.H00, ut));

  1123.         // compute geodetic position (km and °)
  1124.         final FieldGeodeticPoint<T> inBody = earth.transform(position, frame, date);
  1125.         final T alt = inBody.getAltitude().divide(1000.);
  1126.         final T lon = FastMath.toDegrees(inBody.getLongitude());
  1127.         final T lat = FastMath.toDegrees(inBody.getLatitude());

  1128.         // compute local solar time
  1129.         final T lst = localSolarTime(date, position, frame);

  1130.         // get solar activity data and compute
  1131.         final FieldOutput<T> out = new FieldOutput<>(doy, sec, lat, lon, lst,
  1132.                                                      inputParams.getAverageFlux(dateD),
  1133.                                                      inputParams.getDailyFlux(dateD), inputParams.getAp(dateD));
  1134.         out.gtd7d(alt);

  1135.         // return the local density
  1136.         return out.getDensity(TOTAL_MASS);

  1137.     }

  1138.     /** Get local solar time.
  1139.      * @param date current date
  1140.      * @param position current position in frame
  1141.      * @param frame the frame in which is defined the position
  1142.      * @return the local solar time (hour in [0, 24[)
  1143.      */
  1144.     private double localSolarTime(final AbsoluteDate date,
  1145.                                   final Vector3D position,
  1146.                                   final Frame frame) {
  1147.         final Vector3D sunPos = getSunPosition(date, frame);
  1148.         final double lst = FastMath.PI + FastMath.atan2(
  1149.                 sunPos.getX() * position.getY() - sunPos.getY() * position.getX(),
  1150.                 sunPos.getX() * position.getX() + sunPos.getY() * position.getY());
  1151.         return lst * 12. / FastMath.PI;
  1152.     }

  1153.     /** Get local solar time.
  1154.      * @param date current date
  1155.      * @param position current position in frame
  1156.      * @param frame the frame in which is defined the position
  1157.      * @param <T> type of the filed elements
  1158.      * @return the local solar time (hour in [0, 24[)
  1159.      */
  1160.     private <T extends CalculusFieldElement<T>> T localSolarTime(final FieldAbsoluteDate<T> date,
  1161.                                                              final FieldVector3D<T> position,
  1162.                                                              final Frame frame) {
  1163.         final FieldVector3D<T> sunPos = getSunPosition(date, frame);
  1164.         final T y  = position.getY().multiply(sunPos.getX()).subtract(position.getX().multiply(sunPos.getY()));
  1165.         final T x  = position.getX().multiply(sunPos.getX()).add(position.getY().multiply(sunPos.getY()));
  1166.         final T hl = y.atan2(x).add(y.getPi());

  1167.         return hl.divide(y.getPi()).multiply(12.);

  1168.     }

  1169.     /**
  1170.      * This class is a placeholder for the computed densities and temperatures.
  1171.      * <p>
  1172.      * Densities are provided as an array d such as:
  1173.      * <ul>
  1174.      * <li>d[0] = He number density (1/m³)</li>
  1175.      * <li>d[1] = O number density (1/m³)</li>
  1176.      * <li>d[2] = N2 number density (1/m³)</li>
  1177.      * <li>d[3] = O2 number density (1/m³)</li>
  1178.      * <li>d[4] = Ar number density (1/m³)</li>
  1179.      * <li>d[5] = total mass density (kg/m³) (*)</li>
  1180.      * <li>d[6] = H number density (1/m³)</li>
  1181.      * <li>d[7] = N number density (1/m³)</li>
  1182.      * <li>d[8] = anomalous oxygen number density (1/m³)
  1183.      * </ul>
  1184.      * Total mass density, d[5], is NOT the same for methods gtd7 and gtd7d:
  1185.      * <ul>
  1186.      * <li>For gtd7: d[5] is the sum of the mass densities of the species
  1187.      * He, O, N2, O2, Ar, H and N but does NOT include anomalous oxygen.</li>
  1188.      * <li>For gtd7d: d[5] is the "effective total mass density for drag" and is the sum
  1189.      * of the mass densities of all species in this model, INCLUDING anomalous oxygen.</li>
  1190.      * </ul>
  1191.      * O, H, and N are set to zero below 72.5 km.
  1192.      * </p>
  1193.      * <p>
  1194.      * Temperatures are provided as an array t such as:
  1195.      * <ul>
  1196.      * <li>t[0] = exospheric temperature (K)</li>
  1197.      * <li>t[1] = temperature at altitude (K)</li>
  1198.      * </ul>
  1199.      * t[0] is set to global average for altitudes below 120 km.<br>
  1200.      * The 120 km gradient is left at global average value for altitudes below 72 km.
  1201.      * </p>
  1202.      */
  1203.     private class Output {

  1204.         /** Day of year (from 1 to 365 or 366). */
  1205.         private final int doy;

  1206.         /** Seconds in day (UT scale). */
  1207.         private final double sec;

  1208.         /** Geodetic latitude (°). */
  1209.         private final double lat;

  1210.         /** Geodetic longitude (°). */
  1211.         private final double lon;

  1212.         /** Local apparent solar time (hours). */
  1213.         private final double hl;

  1214.         /** 81 day average of F10.7 flux (centered on day). */
  1215.         private final double f107a;

  1216.         /** Daily F10.7 flux for previous day. */
  1217.         private final double f107;

  1218.         /** Array containing:
  1219.         *  <ul>
  1220.         *  <li>0: daily Ap</li>
  1221.         *  <li>1: 3 hr ap index for current time</li>
  1222.         *  <li>2: 3 hr ap index for 3 hrs before current time</li>
  1223.         *  <li>3: 3 hr ap index for 6 hrs before current time</li>
  1224.         *  <li>4: 3 hr ap index for FOR 9 hrs before current time</li>
  1225.         *  <li>5: average of eight 3 hr ap indices from 12 to 33 hrs prior to current time</li>
  1226.         *  <li>6: average of eight 3 hr ap indices from 36 to 57 hrs prior to current time</li>
  1227.         *  </ul>. */
  1228.         private final double[] ap;

  1229.         /** Gravity at latitude (cm/s2). */
  1230.         private final double glat;

  1231.         /** Effective Earth radius at latitude (km). */
  1232.         private final double rlat;

  1233.         /** N2 mixed density at alt. */
  1234.         private double dm28;

  1235.         /** Legendre polynomials. */
  1236.         private final double[][] plg;

  1237.         /** Cosinus of local solar time. */
  1238.         private final double ctloc;
  1239.         /** Sinus of local solar time. */
  1240.         private final double stloc;
  1241.         /** Square of ctloc. */
  1242.         private final double c2tloc;
  1243.         /** Square of stloc. */
  1244.         private final double s2tloc;
  1245.         /** Cube of ctloc. */
  1246.         private final double c3tloc;
  1247.         /** Cube of stloc. */
  1248.         private final double s3tloc;

  1249.         /** Magnetic activity based on daily ap. */
  1250.         private double apdf;

  1251.         /** Magnetic activity based on daily ap. */
  1252.         private double apt;

  1253.         /** Temperature at nodes for ZN1 scale. */
  1254.         private final double[] meso_tn1;

  1255.         /** Temperature at nodes for ZN2 scale. */
  1256.         private final double[] meso_tn2;

  1257.         /** Temperature at nodes for ZN3 scale. */
  1258.         private final double[] meso_tn3;

  1259.         /** Temperature gradients at end nodes for ZN1 scale. */
  1260.         private final double[] meso_tgn1;

  1261.         /** Temperature gradients at end nodes for ZN2 scale. */
  1262.         private final double[] meso_tgn2;

  1263.         /** Temperature gradients at end nodes for ZN3 scale. */
  1264.         private final double[] meso_tgn3;

  1265.         /** Densities. */
  1266.         private final double[] densities;

  1267.         /** Temperatures. */
  1268.         private final double[] temperatures;

  1269.         /** Simple constructor.
  1270.          *  @param doy day of year (from 1 to 365 or 366)
  1271.          *  @param sec seconds in day (UT scale)
  1272.          *  @param lat geodetic latitude (°)
  1273.          *  @param lon geodetic longitude (°)
  1274.          *  @param hl local apparent solar time (hours)
  1275.          *  @param f107a 81 day average of F10.7 flux (centered on day)
  1276.          *  @param f107 daily F10.7 flux for previous day
  1277.          *  @param ap array containing:
  1278.          *  <ul>
  1279.          *  <li>0: daily Ap</li>
  1280.          *  <li>1: 3 hr ap index for current time</li>
  1281.          *  <li>2: 3 hr ap index for 3 hrs before current time</li>
  1282.          *  <li>3: 3 hr ap index for 6 hrs before current time</li>
  1283.          *  <li>4: 3 hr ap index for FOR 9 hrs before current time</li>
  1284.          *  <li>5: average of eight 3 hr ap indices from 12 to 33 hrs prior to current time</li>
  1285.          *  <li>6: average of eight 3 hr ap indices from 36 to 57 hrs prior to current time</li>
  1286.          *  </ul>
  1287.          */
  1288.         Output(final int doy, final double sec,
  1289.                final double lat, final double lon, final double hl,
  1290.                final double f107a, final double f107, final double[] ap) {

  1291.             this.doy   = doy;
  1292.             this.sec   = sec;
  1293.             this.lat   = lat;
  1294.             this.lon   = lon;
  1295.             this.hl    = hl;
  1296.             this.f107a = f107a;
  1297.             this.f107  = f107;
  1298.             this.ap    = ap.clone();

  1299.             this.plg       = new double[4][8];

  1300.             this.meso_tn1  = new double[ZN1.length];
  1301.             this.meso_tn2  = new double[ZN2.length];
  1302.             this.meso_tn3  = new double[ZN3.length];
  1303.             this.meso_tgn1 = new double[2];
  1304.             this.meso_tgn2 = new double[2];
  1305.             this.meso_tgn3 = new double[2];

  1306.             densities       = new double[9];
  1307.             temperatures    = new double[2];

  1308.             // Calculates latitude variable gravity and effective radius
  1309.             final double xlat = (sw[2] == 0) ? LAT_REF : lat;
  1310.             final double c2   = FastMath.cos(2 * DEG_TO_RAD * xlat);
  1311.             glat = G_REF * (1. - .0026373 * c2);
  1312.             rlat = 2. * glat / (3.085462e-6 + 2.27e-9 * c2) * 1.e-5;

  1313.             // Convert latitude into radians
  1314.             final double latr = DEG_TO_RAD * lat;

  1315.             // Calculate legendre polynomials
  1316.             final SinCos scLatr = FastMath.sinCos(latr);
  1317.             final double c      = scLatr.sin();
  1318.             final double s      = scLatr.cos();

  1319.             plg[0][1] = c;
  1320.             plg[0][2] = ( 3.0 * c * plg[0][1] - 1.0) / 2.0;
  1321.             plg[0][3] = ( 5.0 * c * plg[0][2] - 2.0 * plg[0][1]) / 3.0;
  1322.             plg[0][4] = ( 7.0 * c * plg[0][3] - 3.0 * plg[0][2]) / 4.0;
  1323.             plg[0][5] = ( 9.0 * c * plg[0][4] - 4.0 * plg[0][3]) / 5.0;
  1324.             plg[0][6] = (11.0 * c * plg[0][5] - 5.0 * plg[0][4]) / 6.0;

  1325.             plg[1][1] = s;
  1326.             plg[1][2] =   3.0 * c * plg[1][1];
  1327.             plg[1][3] = ( 5.0 * c * plg[1][2] - 3.0 * plg[1][1]) / 2.0;
  1328.             plg[1][4] = ( 7.0 * c * plg[1][3] - 4.0 * plg[1][2]) / 3.0;
  1329.             plg[1][5] = ( 9.0 * c * plg[1][4] - 5.0 * plg[1][3]) / 4.0;
  1330.             plg[1][6] = (11.0 * c * plg[1][5] - 6.0 * plg[1][4]) / 5.0;

  1331.             plg[2][2] = 3.0 * s * plg[1][1];
  1332.             plg[2][3] =   5.0 * c * plg[2][2];
  1333.             plg[2][4] = ( 7.0 * c * plg[2][3] - 5.0 * plg[2][2]) / 2.0;
  1334.             plg[2][5] = ( 9.0 * c * plg[2][4] - 6.0 * plg[2][3]) / 3.0;
  1335.             plg[2][6] = (11.0 * c * plg[2][5] - 7.0 * plg[2][4]) / 4.0;
  1336.             plg[2][7] = (13.0 * c * plg[2][6] - 8.0 * plg[2][5]) / 5.0;

  1337.             plg[3][3] = 5.0 * s * plg[2][2];
  1338.             plg[3][4] =   7.0 * c * plg[3][3];
  1339.             plg[3][5] = ( 9.0 * c * plg[3][4] - 7.0 * plg[3][3]) / 2.0;
  1340.             plg[3][6] = (11.0 * c * plg[3][5] - 8.0 * plg[3][4]) / 3.0;

  1341.             // Calculate additional data
  1342.             if (!(sw[7] == 0 && sw[8] == 0 && sw[14] == 0)) {
  1343.                 final double tloc = HOUR_TO_RAD * hl;
  1344.                 final SinCos sc  = FastMath.sinCos(tloc);
  1345.                 final SinCos sc2 = SinCos.sum(sc, sc);
  1346.                 final SinCos sc3 = SinCos.sum(sc, sc2);
  1347.                 stloc  = sc.sin();
  1348.                 ctloc  = sc.cos();
  1349.                 s2tloc = sc2.sin();
  1350.                 c2tloc = sc2.cos();
  1351.                 s3tloc = sc3.sin();
  1352.                 c3tloc = sc3.cos();
  1353.             } else {
  1354.                 stloc  = 0;
  1355.                 ctloc  = 0;
  1356.                 s2tloc = 0;
  1357.                 c2tloc = 0;
  1358.                 s3tloc = 0;
  1359.                 c3tloc = 0;
  1360.             }

  1361.         }

  1362.         /** Calculate temperatures and densities not including anomalous oxygen.
  1363.          *  <p>
  1364.          *  This method is the thermospheric portion of NRLMSISE-00 for alt > 72.5 km.
  1365.          *  </p>
  1366.          *  <p>NOTES ON INPUT VARIABLES:<br>
  1367.          *  Seconds, Local Time, and Longitude are used independently in the
  1368.          *  model and are not of equal importance for every situation.<br>
  1369.          *  For the most physically realistic calculation these three
  1370.          *  variables should be consistent (lst=sec/3600 + lon/15).<br>
  1371.          *  The Equation of Time departures from the above formula
  1372.          *  for apparent local time can be included if available but
  1373.          *  are of minor importance.<br><br>
  1374.          *
  1375.          *  f107 and f107A values used to generate the model correspond
  1376.          *  to the 10.7 cm radio flux at the actual distance of the Earth
  1377.          *  from the Sun rather than the radio flux at 1 AU. The following
  1378.          *  site provides both classes of values:<br>
  1379.          *  ftp://ftp.ngdc.noaa.gov/STP/SOLAR_DATA/SOLAR_RADIO/FLUX/<br><br>
  1380.          *
  1381.          *  f107, f107A, and ap effects are neither large nor well established below 80 km
  1382.          *  and these parameters should be set to 150., 150., and 4. respectively.
  1383.          *  </p>
  1384.          *  @param alt altitude (km)
  1385.          */
  1386.         void gts7(final double alt) {

  1387.             // Thermal diffusion coefficients for species
  1388.             final double[] alpha = {-0.38, 0.0, 0.0, 0.0, 0.17, 0.0, -0.38, 0.0, 0.0};
  1389.             // Altitude limits for net density computation for species
  1390.             final double[] altl  = {200.0, 300.0, 160.0, 250.0, 240.0, 450.0, 320.0, 450.0};
  1391.             // N2 mixed density
  1392.             final double xmm = PDM[2][4];

  1393.             /**** Exospheric temperature ****/
  1394.             double tinf = PTM[0] * PT[0];
  1395.             // Tinf variations not important below ZA or ZN[0]
  1396.             if (alt > ZN1[0]) {
  1397.                 tinf *= 1.0 + sw[16] * globe7(PT);
  1398.             }
  1399.             setTemperature(EXOSPHERIC, tinf);

  1400.             // Gradient variations not important below ZN[4]
  1401.             double g0 = PTM[3] * PS[0];
  1402.             if (alt > ZN1[4]) {
  1403.                 g0 *= 1.0 + sw[19] * globe7(PS);
  1404.             }

  1405.             // Temperature at lower boundary
  1406.             double tlb = PTM[1] * PD[3][0];
  1407.             tlb *= 1.0 + sw[17] * globe7(PD[3]);

  1408.             // Slope
  1409.             final double s = g0 / (tinf - tlb);

  1410.             // Lower thermosphere temp variations not significant for density above 300 km
  1411.             meso_tn1[1]  = PTM[6] * PTL[0][0];
  1412.             meso_tn1[2]  = PTM[2] * PTL[1][0];
  1413.             meso_tn1[3]  = PTM[7] * PTL[2][0];
  1414.             meso_tn1[4]  = PTM[4] * PTL[3][0];
  1415.             meso_tgn1[1] = PTM[8] * PMA[8][0];
  1416.             if (alt < 300.0) {
  1417.                 final double r = PTM[4] * PTL[3][0];
  1418.                 meso_tn1[1]  /= 1.0 - sw[18] * glob7s(PTL[0]);
  1419.                 meso_tn1[2]  /= 1.0 - sw[18] * glob7s(PTL[1]);
  1420.                 meso_tn1[3]  /= 1.0 - sw[18] * glob7s(PTL[2]);
  1421.                 meso_tn1[4]  /= 1.0 - sw[18] * sw[20] * glob7s(PTL[3]);
  1422.                 meso_tgn1[1] *= 1.0 + sw[18] * sw[20] * glob7s(PMA[8]);
  1423.                 meso_tgn1[1] *= meso_tn1[4] * meso_tn1[4] / (r * r);
  1424.             }

  1425.             /**** Temperature at altitude ****/
  1426.             setTemperature(ALTITUDE, densu(alt, 1.0, tinf, tlb, 0.0, 0.0, PTM[5], s));

  1427.             /**** N2 density ****/
  1428.             /*   Density variation factor at Zlb */
  1429.             final double g28 = sw[21] * globe7(PD[2]);
  1430.             /* Diffusive density at Zlb */
  1431.             final double db28 = PDM[2][0] * FastMath.exp(g28) * PD[2][0];
  1432.             /* Diffusive density at Alt */
  1433.             double diffusiveDensity = densu(alt, db28, tinf, tlb, N2_MASS, alpha[2], PTM[5], s);
  1434.             setDensity(MOLECULAR_NITROGEN, diffusiveDensity);
  1435.             // Variation of turbopause height
  1436.             final double zhf = PDL[1][24] * (1.0 + sw[5] * PDL[0][24] *
  1437.                                        FastMath.sin(DEG_TO_RAD * lat) *
  1438.                                        FastMath.cos(DAY_TO_RAD * (doy - PT[13])));
  1439.             /* Turbopause */
  1440.             final double zh28  = PDM[2][2] * zhf;
  1441.             final double zhm28 = PDM[2][3] * PDL[1][5];
  1442.             /* Mixed density at Zlb */
  1443.             final double b28 = densu(zh28, db28, tinf, tlb, N2_MASS - xmm, alpha[2] - 1.0, PTM[5], s);
  1444.             if (sw[15] != 0 && alt <= altl[2]) {
  1445.                 /*  Mixed density at Alt */
  1446.                 dm28 = densu(alt, b28, tinf, tlb, xmm, alpha[2], PTM[5], s);
  1447.                 /*  Net density at Alt */
  1448.                 setDensity(MOLECULAR_NITROGEN, dnet(diffusiveDensity, dm28, zhm28, xmm, N2_MASS));
  1449.             }

  1450.             /**** He density ****/
  1451.             /*   Density variation factor at Zlb */
  1452.             final double g4 = sw[21] * globe7(PD[0]);
  1453.             /*  Diffusive density at Zlb */
  1454.             final double db04 = PDM[0][0] * FastMath.exp(g4) * PD[0][0];
  1455.             /*  Diffusive density at Alt */
  1456.             diffusiveDensity = densu(alt, db04, tinf, tlb, HE_MASS, alpha[0], PTM[5], s);
  1457.             setDensity(HELIUM, diffusiveDensity);
  1458.             if (sw[15] != 0 && alt < altl[0]) {
  1459.                 /*  Turbopause */
  1460.                 final double zh04 = PDM[0][2];
  1461.                 /*  Mixed density at Zlb */
  1462.                 final double b04 = densu(zh04, db04, tinf, tlb, HE_MASS - xmm, alpha[0] - 1., PTM[5], s);
  1463.                 /*  Mixed density at Alt */
  1464.                 final double dm04 = densu(alt, b04, tinf, tlb, xmm, 0., PTM[5], s);
  1465.                 final double zhm04 = zhm28;
  1466.                 /*  Net density at Alt */
  1467.                 diffusiveDensity = dnet(diffusiveDensity, dm04, zhm04, xmm, HE_MASS);
  1468.                 /*  Correction to specified mixing ratio at ground */
  1469.                 final double rl = FastMath.log(b28 * PDM[0][1] / b04);
  1470.                 final double zc04 = PDM[0][4] * PDL[1][0];
  1471.                 final double hc04 = PDM[0][5] * PDL[1][1];
  1472.                 /*  Net density corrected at Alt */
  1473.                 setDensity(HELIUM, diffusiveDensity * ccor(alt, rl, hc04, zc04));
  1474.             }

  1475.             /**** O density ****/
  1476.             /* Density variation factor at Zlb */
  1477.             final double g16 = sw[21] * globe7(PD[1]);
  1478.             /* Diffusive density at Zlb */
  1479.             final double db16 = PDM[1][0] * FastMath.exp(g16) * PD[1][0];
  1480.             /* Diffusive density at Alt */
  1481.             diffusiveDensity = densu(alt, db16, tinf, tlb, O_MASS, alpha[1], PTM[5], s);
  1482.             setDensity(ATOMIC_OXYGEN, diffusiveDensity);
  1483.             if (sw[15] != 0 && alt < altl[1]) {
  1484.                 /* Turbopause */
  1485.                 final double zh16 = PDM[1][2];
  1486.                 /* Mixed density at Zlb */
  1487.                 final double b16 = densu(zh16, db16, tinf, tlb, O_MASS - xmm, alpha[1] - 1.0, PTM[5], s);
  1488.                 /* Mixed density at Alt */
  1489.                 final double dm16 = densu(alt, b16, tinf, tlb, xmm, 0., PTM[5], s);
  1490.                 final double zhm16 = zhm28;
  1491.                 /* Net density at Alt */
  1492.                 diffusiveDensity = dnet(diffusiveDensity, dm16, zhm16, xmm, O_MASS);
  1493.                 final double rl = PDM[1][1] * PDL[1][16] * (1.0 + sw[1] * PDL[0][23] * (f107a - FLUX_REF));
  1494.                 final double hc16 = PDM[1][5] * PDL[1][3];
  1495.                 final double zc16 = PDM[1][4] * PDL[1][2];
  1496.                 final double hc216 = PDM[1][5] * PDL[1][4];
  1497.                 diffusiveDensity *= ccor2(alt, rl, hc16, zc16, hc216);
  1498.                 /* Chemistry correction */
  1499.                 final double hcc16 = PDM[1][7] * PDL[1][13];
  1500.                 final double zcc16 = PDM[1][6] * PDL[1][12];
  1501.                 final double rc16  = PDM[1][3] * PDL[1][14];
  1502.                 /* Net density corrected at Alt */
  1503.                 setDensity(ATOMIC_OXYGEN, diffusiveDensity * ccor(alt, rc16, hcc16, zcc16));
  1504.             }

  1505.             /**** O2 density ****/
  1506.             /* Density variation factor at Zlb */
  1507.             final double g32 = sw[21] * globe7(PD[4]);
  1508.             /* Diffusive density at Zlb */
  1509.             final double db32 = PDM[3][0] * FastMath.exp(g32) * PD[4][0];
  1510.             /* Diffusive density at Alt */
  1511.             diffusiveDensity = densu(alt, db32, tinf, tlb, O2_MASS, alpha[3], PTM[5], s);
  1512.             setDensity(MOLECULAR_OXYGEN, diffusiveDensity);
  1513.             if (sw[15] != 0) {
  1514.                 if (alt <= altl[3]) {
  1515.                     /* Turbopause */
  1516.                     final double zh32 = PDM[3][2];
  1517.                     /* Mixed density at Zlb */
  1518.                     final double b32 = densu(zh32, db32, tinf, tlb, O2_MASS - xmm, alpha[3] - 1., PTM[5], s);
  1519.                     /* Mixed density at Alt */
  1520.                     final double dm32 = densu(alt, b32, tinf, tlb, xmm, 0., PTM[5], s);
  1521.                     final double zhm32 = zhm28;
  1522.                     /* Net density at Alt */
  1523.                     diffusiveDensity = dnet(diffusiveDensity, dm32, zhm32, xmm, O2_MASS);
  1524.                     /* Correction to specified mixing ratio at ground */
  1525.                     final double rl = FastMath.log(b28 * PDM[3][1] / b32);
  1526.                     final double hc32 = PDM[3][5] * PDL[1][7];
  1527.                     final double zc32 = PDM[3][4] * PDL[1][6];
  1528.                     diffusiveDensity *= ccor(alt, rl, hc32, zc32);
  1529.                 }
  1530.                 /* Correction for general departure from diffusive equilibrium above Zlb */
  1531.                 final double hcc32  = PDM[3][7] * PDL[1][22];
  1532.                 final double hcc232 = PDM[3][7] * PDL[0][22];
  1533.                 final double zcc32  = PDM[3][6] * PDL[1][21];
  1534.                 final double rc32   = PDM[3][3] * PDL[1][23] * (1. + sw[1] * PDL[0][23] * (f107a - FLUX_REF));
  1535.                 /* Net density corrected at Alt */
  1536.                 setDensity(MOLECULAR_OXYGEN, diffusiveDensity * ccor2(alt, rc32, hcc32, zcc32, hcc232));
  1537.             }

  1538.             /**** Ar density ****/
  1539.             /* Density variation factor at Zlb */
  1540.             final double g40 = sw[21] * globe7(PD[5]);
  1541.             /* Diffusive density at Zlb */
  1542.             final double db40 = PDM[4][0] * FastMath.exp(g40) * PD[5][0];
  1543.             /* Diffusive density at Alt */
  1544.             diffusiveDensity = densu(alt, db40, tinf, tlb, AR_MASS, alpha[4], PTM[5], s);
  1545.             setDensity(ARGON, diffusiveDensity);
  1546.             if (sw[15] != 0 && alt <= altl[4]) {
  1547.                 /* Turbopause */
  1548.                 final double zh40 = PDM[4][2];
  1549.                 /* Mixed density at Zlb */
  1550.                 final double b40 = densu(zh40, db40, tinf, tlb, AR_MASS - xmm, alpha[4] - 1., PTM[5], s);
  1551.                 /* Mixed density at Alt */
  1552.                 final double dm40 = densu(alt, b40, tinf, tlb, xmm, 0., PTM[5], s);
  1553.                 final double zhm40 = zhm28;
  1554.                 /* Net density at Alt */
  1555.                 diffusiveDensity = dnet(diffusiveDensity, dm40, zhm40, xmm, AR_MASS);
  1556.                 /* Correction to specified mixing ratio at ground */
  1557.                 final double rl = FastMath.log(b28 * PDM[4][1] / b40);
  1558.                 final double hc40 = PDM[4][5] * PDL[1][9];
  1559.                 final double zc40 = PDM[4][4] * PDL[1][8];
  1560.                 /* Net density corrected at Alt */
  1561.                 setDensity(ARGON, diffusiveDensity * ccor(alt, rl, hc40, zc40));
  1562.             }

  1563.             /**** H density ****/
  1564.             /* Density variation factor at Zlb */
  1565.             final double g1 = sw[21] * globe7(PD[6]);
  1566.             /* Diffusive density at Zlb */
  1567.             final double db01 = PDM[5][0] * FastMath.exp(g1) * PD[6][0];
  1568.             /* Diffusive density at Alt */
  1569.             diffusiveDensity = densu(alt, db01, tinf, tlb, H_MASS, alpha[6], PTM[5], s);
  1570.             setDensity(HYDROGEN, diffusiveDensity);
  1571.             if (sw[15] != 0 && alt <= altl[6]) {
  1572.                 /* Turbopause */
  1573.                 final double zh01 = PDM[5][2];
  1574.                 /* Mixed density at Zlb */
  1575.                 final double b01 = densu(zh01, db01, tinf, tlb, H_MASS - xmm, alpha[6] - 1., PTM[5], s);
  1576.                 /* Mixed density at Alt */
  1577.                 final double dm01 = densu(alt, b01, tinf, tlb, xmm, 0., PTM[5], s);
  1578.                 final double zhm01 = zhm28;
  1579.                 /* Net density at Alt */
  1580.                 diffusiveDensity = dnet(diffusiveDensity, dm01, zhm01, xmm, H_MASS);
  1581.                 /* Correction to specified mixing ratio at ground */
  1582.                 final double rl = FastMath.log(b28 * PDM[5][1] * FastMath.sqrt(PDL[1][17] * PDL[1][17]) / b01);
  1583.                 final double hc01 = PDM[5][5] * PDL[1][11];
  1584.                 final double zc01 = PDM[5][4] * PDL[1][10];
  1585.                 diffusiveDensity *= ccor(alt, rl, hc01, zc01);
  1586.                 /* Chemistry correction */
  1587.                 final double hcc01 = PDM[5][7] * PDL[1][19];
  1588.                 final double zcc01 = PDM[5][6] * PDL[1][18];
  1589.                 final double rc01 = PDM[5][3] * PDL[1][20];
  1590.                 /* Net density corrected at Alt */
  1591.                 setDensity(HYDROGEN, diffusiveDensity * ccor(alt, rc01, hcc01, zcc01));
  1592.             }

  1593.             /**** N density ****/
  1594.             /* Density variation factor at Zlb */
  1595.             final double g14 = sw[21] * globe7(PD[7]);
  1596.             /* Diffusive density at Zlb */
  1597.             final double db14 = PDM[6][0] * FastMath.exp(g14) * PD[7][0];
  1598.             /* Diffusive density at Alt */
  1599.             diffusiveDensity = densu(alt, db14, tinf, tlb, N_MASS, alpha[7], PTM[5], s);
  1600.             setDensity(ATOMIC_NITROGEN, diffusiveDensity);
  1601.             if (sw[15] != 0 && alt <= altl[7]) {
  1602.                 /* Turbopause */
  1603.                 final double zh14 = PDM[6][2];
  1604.                 /* Mixed density at Zlb */
  1605.                 final double b14 = densu(zh14, db14, tinf, tlb, N_MASS - xmm, alpha[7] - 1., PTM[5], s);
  1606.                 /* Mixed density at Alt */
  1607.                 final double dm14 = densu(alt, b14, tinf, tlb, xmm, 0., PTM[5], s);
  1608.                 final double zhm14 = zhm28;
  1609.                 /* Net density at Alt */
  1610.                 diffusiveDensity = dnet(diffusiveDensity, dm14, zhm14, xmm, N_MASS);
  1611.                 /* Correction to specified mixing ratio at ground */
  1612.                 final double rl = FastMath.log(b28 * PDM[6][1] * PDL[0][2] / b14);
  1613.                 final double hc14 = PDM[6][5] * PDL[0][1];
  1614.                 final double zc14 = PDM[6][4] * PDL[0][0];
  1615.                 diffusiveDensity *= ccor(alt, rl, hc14, zc14);
  1616.                 /* Chemistry correction */
  1617.                 final double hcc14 = PDM[6][7] * PDL[0][4];
  1618.                 final double zcc14 = PDM[6][6] * PDL[0][3];
  1619.                 final double rc14 = PDM[6][3] * PDL[0][5];
  1620.                 /* Net density corrected at Alt */
  1621.                 setDensity(ATOMIC_NITROGEN, diffusiveDensity * ccor(alt, rc14, hcc14, zcc14));
  1622.             }

  1623.             /**** Anomalous O density ****/
  1624.             final double g16h  = sw[21] * globe7(PD[8]);
  1625.             final double db16h = PDM[7][0] * FastMath.exp(g16h) * PD[8][0];
  1626.             final double tho   = PDM[7][9] * PDL[0][6];
  1627.             diffusiveDensity = densu(alt, db16h, tho, tho, O_MASS, alpha[8], PTM[5], s);
  1628.             final double zsht = PDM[7][5];
  1629.             final double zmho = PDM[7][4];
  1630.             final double zsho = scalh(zmho, O_MASS, tho);
  1631.             diffusiveDensity *= FastMath.exp(-zsht / zsho * (FastMath.exp((zmho - alt ) / zsht) - 1.));
  1632.             setDensity(ANOMALOUS_OXYGEN, diffusiveDensity);

  1633.             // Convert densities from cm-3 to m-3
  1634.             for (int i = 0; i < 9; i++) {
  1635.                 setDensity(i, getDensity(i) * 1.0e+06);
  1636.             }

  1637.             /**** Total mass density ****/
  1638.             final double tmd = AMU * (HE_MASS * getDensity(HELIUM) +
  1639.                                       O_MASS  * getDensity(ATOMIC_OXYGEN) +
  1640.                                       N2_MASS * getDensity(MOLECULAR_NITROGEN) +
  1641.                                       O2_MASS * getDensity(MOLECULAR_OXYGEN) +
  1642.                                       AR_MASS * getDensity(ARGON) +
  1643.                                       H_MASS  * getDensity(HYDROGEN) +
  1644.                                       N_MASS  * getDensity(ATOMIC_NITROGEN));
  1645.             setDensity(TOTAL_MASS, tmd);

  1646.         }

  1647.         /** Calculate temperatures and densities not including anomalous oxygen.
  1648.          *  <p>NOTES ON INPUT VARIABLES:<br>
  1649.          *  Seconds, Local Time, and Longitude are used independently in the
  1650.          *  model and are not of equal importance for every situation.<br>
  1651.          *  For the most physically realistic calculation these three
  1652.          *  variables should be consistent (lst=sec/3600 + lon/15).<br>
  1653.          *  The Equation of Time departures from the above formula
  1654.          *  for apparent local time can be included if available but
  1655.          *  are of minor importance.<br><br>
  1656.          *
  1657.          *  f107 and f107A values used to generate the model correspond
  1658.          *  to the 10.7 cm radio flux at the actual distance of the Earth
  1659.          *  from the Sun rather than the radio flux at 1 AU. The following
  1660.          *  site provides both classes of values:<br>
  1661.          *  ftp://ftp.ngdc.noaa.gov/STP/SOLAR_DATA/SOLAR_RADIO/FLUX/<br><br>
  1662.          *
  1663.          *  f107, f107A, and ap effects are neither large nor well established below 80 km
  1664.          *  and these parameters should be set to 150., 150., and 4. respectively.
  1665.          *  </p>
  1666.          *  @param alt altitude (km)
  1667.          */
  1668.         void gtd7(final double alt) {

  1669.             // Calculates for thermosphere/mesosphere (above ZN2[0])
  1670.             final double altt = (alt > ZN2[0]) ? alt : ZN2[0];
  1671.             gts7(altt);
  1672.             if (alt >= ZN2[0]) {
  1673.                 return;
  1674.             }

  1675.             // Calculates for lower mesosphere/upper stratosphere (between ZN2[0] and ZN3[0]):
  1676.             // Temperature at nodes and gradients at end nodes
  1677.             // Inverse temperature a linear function of spherical harmonics
  1678.             final double r = PMA[2][0] * PAVGM[2];
  1679.             meso_tgn2[0] = meso_tgn1[1];
  1680.             meso_tn2[0]  = meso_tn1[4];
  1681.             meso_tn2[1]  = PMA[0][0] * PAVGM[0] / (1.0 - sw[20] * glob7s(PMA[0]));
  1682.             meso_tn2[2]  = PMA[1][0] * PAVGM[1] / (1.0 - sw[20] * glob7s(PMA[1]));
  1683.             meso_tn2[3]  = PMA[2][0] * PAVGM[2] / (1.0 - sw[20] * sw[22] * glob7s(PMA[2]));
  1684.             meso_tgn2[1] = PMA[9][0] * PAVGM[8] * (1.0 + sw[20] * sw[22] * glob7s(PMA[9])) *
  1685.                            meso_tn2[3] * meso_tn2[3] / (r * r);
  1686.             meso_tn3[0]  = meso_tn2[3];

  1687.             // Calculates for lower stratosphere and troposphere (below ZN3[0])
  1688.             // Temperature at nodes and gradients at end nodes
  1689.             // Inverse temperature a linear function of spherical harmonics
  1690.             if (alt <= ZN3[0]) {
  1691.                 final double q = PMA[6][0] * PAVGM[6];
  1692.                 meso_tgn3[0] = meso_tgn2[1];
  1693.                 meso_tn3[1]  = PMA[3][0] * PAVGM[3] / (1.0 - sw[22] * glob7s(PMA[3]));
  1694.                 meso_tn3[2]  = PMA[4][0] * PAVGM[4] / (1.0 - sw[22] * glob7s(PMA[4]));
  1695.                 meso_tn3[3]  = PMA[5][0] * PAVGM[5] / (1.0 - sw[22] * glob7s(PMA[5]));
  1696.                 meso_tn3[4]  = PMA[6][0] * PAVGM[6] / (1.0 - sw[22] * glob7s(PMA[6]));
  1697.                 meso_tgn3[1] = PMA[7][0] * PAVGM[7] * (1.0 + sw[22] * glob7s(PMA[7])) *
  1698.                                meso_tn3[4] * meso_tn3[4] / (q * q);

  1699.             }

  1700.             // Linear transition to full mixing below ZN2[0]
  1701.             final double dmc = (alt > ZMIX) ? 1.0 - (ZN2[0] - alt) / (ZN2[0] - ZMIX) : 0.;
  1702.             final double dz28 = getDensity(MOLECULAR_NITROGEN);

  1703.             // N2 density
  1704.             final double dm28m = dm28 * 1.0e+06;
  1705.             double dmr = dz28 / dm28m - 1.0;
  1706.             double dst = densm(alt, dm28m, PDM[2][4]) * (1.0 + dmr * dmc);
  1707.             setDensity(MOLECULAR_NITROGEN, dst);

  1708.             // HE density
  1709.             dmr = getDensity(HELIUM) / (dz28 * PDM[0][1]) - 1.0;
  1710.             dst = getDensity(MOLECULAR_NITROGEN) * PDM[0][1] * (1.0 + dmr * dmc);
  1711.             setDensity(HELIUM, dst);

  1712.             // O density
  1713.             setDensity(ATOMIC_OXYGEN, 0.);
  1714.             setDensity(ANOMALOUS_OXYGEN, 0.);

  1715.             // O2 density
  1716.             dmr = getDensity(MOLECULAR_OXYGEN) / (dz28 * PDM[3][1]) - 1.0;
  1717.             dst = getDensity(MOLECULAR_NITROGEN) * PDM[3][1] * (1.0 + dmr * dmc);
  1718.             setDensity(MOLECULAR_OXYGEN, dst);

  1719.             // AR density
  1720.             dmr = getDensity(ARGON) / (dz28 * PDM[4][1]) - 1.0;
  1721.             dst = getDensity(MOLECULAR_NITROGEN) * PDM[4][1] * (1.0 + dmr * dmc);
  1722.             setDensity(ARGON, dst);

  1723.             // H density
  1724.             setDensity(HYDROGEN, 0.);

  1725.             // N density
  1726.             setDensity(ATOMIC_NITROGEN, 0.);

  1727.             // Total mass density
  1728.             final double tmd = AMU * (HE_MASS * getDensity(HELIUM) +
  1729.                                       O_MASS  * getDensity(ATOMIC_OXYGEN) +
  1730.                                       N2_MASS * getDensity(MOLECULAR_NITROGEN) +
  1731.                                       O2_MASS * getDensity(MOLECULAR_OXYGEN) +
  1732.                                       AR_MASS * getDensity(ARGON) +
  1733.                                       H_MASS  * getDensity(HYDROGEN) +
  1734.                                       N_MASS  * getDensity(ATOMIC_NITROGEN));
  1735.             setDensity(TOTAL_MASS, tmd);

  1736.             // Temperature at altitude
  1737.             setTemperature(ALTITUDE, densm(alt, 1.0, 0));

  1738.         }

  1739.         /** Calculate temperatures and densities including anomalous oxygen.
  1740.          *  <p></p>
  1741.          *  <p>NOTES ON INPUT VARIABLES:<br>
  1742.          *  Seconds, Local Time, and Longitude are used independently in the
  1743.          *  model and are not of equal importance for every situation.<br>
  1744.          *  For the most physically realistic calculation these three
  1745.          *  variables should be consistent (lst=sec/3600 + lon/15).<br>
  1746.          *  The Equation of Time departures from the above formula
  1747.          *  for apparent local time can be included if available but
  1748.          *  are of minor importance.<br>
  1749.          *  <br>
  1750.          *  f107 and f107A values used to generate the model correspond
  1751.          *  to the 10.7 cm radio flux at the actual distance of the Earth
  1752.          *  from the Sun rather than the radio flux at 1 AU. The following
  1753.          *  site provides both classes of values:<br>
  1754.          *  ftp://ftp.ngdc.noaa.gov/STP/SOLAR_DATA/SOLAR_RADIO/FLUX/<br>
  1755.          *  <br>
  1756.          *  f107, f107A, and ap effects are neither large nor well established below 80 km
  1757.          *  and these parameters should be set to 150., 150., and 4. respectively.
  1758.          *  </p>
  1759.          *  @param alt altitude (km)
  1760.          */
  1761.         void gtd7d(final double alt) {

  1762.             // Compute densities and temperatures
  1763.             gtd7(alt);

  1764.             // Update the total mass density with anomalous oxygen contribution
  1765.             final double dTot = getDensity(TOTAL_MASS) + AMU * O_MASS * getDensity(ANOMALOUS_OXYGEN);
  1766.             setDensity(TOTAL_MASS, dTot);

  1767.         }

  1768.         /** Set one density.
  1769.          * @param index one of the nine elements :
  1770.          * <ul>
  1771.          * <li>{@link #HELIUM}</li>
  1772.          * <li>{@link #ATOMIC_OXYGEN}</li>
  1773.          * <li>{@link #MOLECULAR_NITROGEN}</li>
  1774.          * <li>{@link #MOLECULAR_OXYGEN}</li>
  1775.          * <li>{@link #ARGON}</li>
  1776.          * <li>{@link #TOTAL_MASS}</li>
  1777.          * <li>{@link #HYDROGEN}</li>
  1778.          * <li>{@link #ATOMIC_NITROGEN}</li>
  1779.          * <li>{@link #ATOMIC_NITROGEN}</li>
  1780.          * </ul>
  1781.          * @param d the value of density to set
  1782.          */
  1783.         void setDensity(final int index, final double d) {
  1784.             densities[index] = d;
  1785.         }

  1786.         /** Set one temperature.
  1787.          * @param index one of the two elements :
  1788.          * <ul>
  1789.          * <li>{@link #EXOSPHERIC}</li>
  1790.          * <li>{@link #ALTITUDE}</li>
  1791.          * </ul>
  1792.          * @param t the value of temperature to set
  1793.          */
  1794.         void setTemperature(final int index, final double t) {
  1795.             temperatures[index] = t;
  1796.         }

  1797.         /** Get one of the stored densities.
  1798.          * @param index one of the nine elements :
  1799.          * <ul>
  1800.          * <li>{@link #HELIUM}</li>
  1801.          * <li>{@link #ATOMIC_OXYGEN}</li>
  1802.          * <li>{@link #MOLECULAR_NITROGEN}</li>
  1803.          * <li>{@link #MOLECULAR_OXYGEN}</li>
  1804.          * <li>{@link #ARGON}</li>
  1805.          * <li>{@link #TOTAL_MASS}</li>
  1806.          * <li>{@link #HYDROGEN}</li>
  1807.          * <li>{@link #ATOMIC_NITROGEN}</li>
  1808.          * <li>{@link #ATOMIC_NITROGEN}</li>
  1809.          * </ul>
  1810.          * @return the requested density
  1811.          */
  1812.         public double getDensity(final int index) {
  1813.             return densities[index];
  1814.         }

  1815.         /** Calculate G(L) function with upper thermosphere parameters.
  1816.          *  @param p array of parameters
  1817.          *  @return G(L) value
  1818.          */
  1819.         private double globe7(final double[] p) {

  1820.             final double[] t = new double[14];
  1821.             final double cd32 = FastMath.cos(DAY_TO_RAD * (doy - p[31]));
  1822.             final double cd18 = FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[17]));
  1823.             final double cd14 = FastMath.cos(DAY_TO_RAD * (doy - p[13]));
  1824.             final double cd39 = FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[38]));

  1825.             // F10.7 effect
  1826.             final double df  = f107  - f107a;
  1827.             final double dfa = f107a - FLUX_REF;
  1828.             t[0] = p[19] * df * (1.0 + p[59] * dfa) + p[20] * df * df + p[21] * dfa + p[29] * dfa * dfa;

  1829.             final double f1 = 1.0 + (p[47] * dfa + p[19] * df + p[20] * df * df) * swc[1];
  1830.             final double f2 = 1.0 + (p[49] * dfa + p[19] * df + p[20] * df * df) * swc[1];

  1831.             // Time independent
  1832.             t[1] = (p[1]  * plg[0][2] + p[2] * plg[0][4] + p[22] * plg[0][6]) +
  1833.                    (p[14] * plg[0][2]) * dfa * swc[1] + p[26] * plg[0][1];

  1834.             // Symmetrical annual
  1835.             t[2] = p[18] * cd32;

  1836.             // Symmetrical semiannual
  1837.             t[3] = (p[15] + p[16] * plg[0][2]) * cd18;

  1838.             // Asymmetrical annual
  1839.             t[4] = f1 * (p[9] * plg[0][1] + p[10] * plg[0][3]) * cd14;

  1840.             // Asymmetrical semiannual
  1841.             t[5] = p[37] * plg[0][1] * cd39;

  1842.             // Diurnal
  1843.             if (sw[7] != 0) {
  1844.                 final double t71 = (p[11] * plg[1][2]) * cd14 * swc[5];
  1845.                 final double t72 = (p[12] * plg[1][2]) * cd14 * swc[5];
  1846.                 t[6] = f2 * ((p[3] * plg[1][1] + p[4] * plg[1][3] + p[27] * plg[1][5] + t71) * ctloc +
  1847.                              (p[6] * plg[1][1] + p[7] * plg[1][3] + p[28] * plg[1][5] + t72) * stloc);
  1848.             }

  1849.             // Semidiurnal
  1850.             if (sw[8] != 0) {
  1851.                 final double t81 = (p[23] * plg[2][3] + p[35] * plg[2][5]) * cd14 * swc[5];
  1852.                 final double t82 = (p[33] * plg[2][3] + p[36] * plg[2][5]) * cd14 * swc[5];
  1853.                 t[7] = f2 * ((p[5] * plg[2][2] + p[41] * plg[2][4] + t81) * c2tloc +
  1854.                              (p[8] * plg[2][2] + p[42] * plg[2][4] + t82) * s2tloc);
  1855.             }

  1856.             // Terdiurnal
  1857.             if (sw[14] != 0) {
  1858.                 t[13] = f2 * ((p[39] * plg[3][3] + (p[93] * plg[3][4] + p[46] * plg[3][6]) * cd14 * swc[5]) * s3tloc +
  1859.                               (p[40] * plg[3][3] + (p[94] * plg[3][4] + p[48] * plg[3][6]) * cd14 * swc[5]) * c3tloc);
  1860.             }

  1861.             // magnetic activity based on daily ap
  1862.             if (sw[9] == -1) {
  1863.                 if (p[51] != 0) {
  1864.                     final double exp1 = FastMath.exp(-10800.0 * FastMath.abs(p[51]) /
  1865.                                                      (1.0 + p[138] * (LAT_REF - FastMath.abs(lat))));
  1866.                     final double p24 = FastMath.max(p[24], 1.0e-4);
  1867.                     apt = sg0(FastMath.min(exp1, 0.99999), p24, p[25]);
  1868.                     t[8] = apt * (p[50] + p[96] * plg[0][2] + p[54] * plg[0][4] +
  1869.                                   (p[125] * plg[0][1] + p[126] * plg[0][3] + p[127] * plg[0][5]) * cd14 * swc[5] +
  1870.                                   (p[128] * plg[1][1] + p[129] * plg[1][3] + p[130] * plg[1][5]) * swc[7] *
  1871.                                   FastMath.cos(HOUR_TO_RAD * (hl - p[131])));
  1872.                 }
  1873.             } else {
  1874.                 final double apd = ap[0] - 4.0;
  1875.                 final double p44 = (p[43] < 0.) ? 1.0E-5 : p[43];
  1876.                 final double p45 = p[44];
  1877.                 apdf = apd + (p45 - 1.0) * (apd + (FastMath.exp(-p44 * apd) - 1.0) / p44);
  1878.                 if (sw[9] != 0) {
  1879.                     t[8] = apdf * (p[32] + p[45] * plg[0][2] + p[34] * plg[0][4] +
  1880.                                    (p[100] * plg[0][1] + p[101] * plg[0][3] + p[102] * plg[0][5]) * cd14 * swc[5] +
  1881.                                    (p[121] * plg[1][1] + p[122] * plg[1][3] + p[123] * plg[1][5]) * swc[7] *
  1882.                                    FastMath.cos(HOUR_TO_RAD * (hl - p[124])));
  1883.                 }
  1884.             }

  1885.             if (sw[10] != 0) {
  1886.                 final double lonr   = DEG_TO_RAD * lon;
  1887.                 final SinCos scLonr = FastMath.sinCos(lonr);
  1888.                 // Longitudinal
  1889.                 if (sw[11] != 0) {
  1890.                     t[10] = (1.0 + p[80] * dfa * swc[1]) *
  1891.                             ((p[64]  * plg[1][2] + p[65]  * plg[1][4] + p[66]  * plg[1][6] +
  1892.                               p[103] * plg[1][1] + p[104] * plg[1][3] + p[105] * plg[1][5] +
  1893.                              (p[109] * plg[1][1] + p[110] * plg[1][3] + p[111] * plg[1][5]) * swc[5] * cd14) *
  1894.                              scLonr.cos() +
  1895.                              (p[90]  * plg[1][2] + p[91]  * plg[1][4] + p[92]  * plg[1][6] +
  1896.                               p[106] * plg[1][1] + p[107] * plg[1][3] + p[108] * plg[1][5] +
  1897.                              (p[112] * plg[1][1] + p[113] * plg[1][3] + p[114] * plg[1][5]) * swc[5] * cd14) *
  1898.                              scLonr.sin());
  1899.                 }

  1900.                 // ut and mixed ut, longitude
  1901.                 if (sw[12] != 0) {
  1902.                     t[11] = (1.0 + p[95]  * plg[0][1]) * (1.0 + p[81] * dfa * swc[1]) *
  1903.                             (1.0 + p[119] * plg[0][1] * swc[5] * cd14) *
  1904.                             (p[68] * plg[0][1] + p[69] * plg[0][3] + p[70] * plg[0][5]) *
  1905.                             FastMath.cos(SEC_TO_RAD * (sec - p[71]));
  1906.                     t[11] += swc[11] * (1.0 + p[137] * dfa * swc[1]) *
  1907.                             (p[76] * plg[2][3] + p[77] * plg[2][5] + p[78] * plg[2][7]) *
  1908.                             FastMath.cos(SEC_TO_RAD * (sec - p[79]) + 2.0 * lonr);
  1909.                 }

  1910.                 /* ut, longitude magnetic activity */
  1911.                 if (sw[13] != 0) {
  1912.                     if (sw[9] == -1) {
  1913.                         if (p[51] != 0.) {
  1914.                             t[12] = apt * swc[11] * (1. + p[132] * plg[0][1]) *
  1915.                                     (p[52] * plg[1][2] + p[98] * plg[1][4] + p[67] * plg[1][6]) *
  1916.                                     FastMath.cos(DEG_TO_RAD * (lon - p[97])) +
  1917.                                     apt * swc[11] * swc[5] * cd14 *
  1918.                                     (p[133] * plg[1][1] + p[134] * plg[1][3] + p[135] * plg[1][5]) *
  1919.                                     FastMath.cos(DEG_TO_RAD * (lon - p[136])) +
  1920.                                     apt * swc[12] *
  1921.                                     (p[55] * plg[0][1] + p[56] * plg[0][3] + p[57] * plg[0][5]) *
  1922.                                     FastMath.cos(SEC_TO_RAD * (sec - p[58]));
  1923.                         }
  1924.                     } else {
  1925.                         t[12] = apdf * swc[11] * (1.0 + p[120] * plg[0][1]) *
  1926.                                 ((p[60] * plg[1][2] + p[61] * plg[1][4] + p[62] * plg[1][6]) *
  1927.                                 FastMath.cos(DEG_TO_RAD * (lon - p[63]))) +
  1928.                                 apdf * swc[11] * swc[5] * cd14 *
  1929.                                 (p[115] * plg[1][1] + p[116] * plg[1][3] + p[117] * plg[1][5]) *
  1930.                                 FastMath.cos(DEG_TO_RAD * (lon - p[118])) +
  1931.                                 apdf * swc[12] *
  1932.                                 (p[83] * plg[0][1] + p[84] * plg[0][3] + p[85] * plg[0][5]) *
  1933.                                 FastMath.cos(SEC_TO_RAD * (sec - p[75]));
  1934.                     }
  1935.                 }
  1936.             }

  1937.             // Sum all effects (params not used: 82, 89, 99, 139-149)
  1938.             double tinf = p[30];
  1939.             for (int i = 0; i < 14; i++) {
  1940.                 tinf += FastMath.abs(sw[i + 1]) * t[i];
  1941.             }

  1942.             // Return G(L)
  1943.             return tinf;

  1944.         }

  1945.         /** Calculate G(L) function with lower atmosphere parameters.
  1946.          *  @param p array of parameters
  1947.          *  @return G(L) value
  1948.          */
  1949.         private double glob7s(final double[] p) {

  1950.             final double[] t = new double[14];
  1951.             final double cd32 = FastMath.cos(DAY_TO_RAD * (doy - p[31]));
  1952.             final double cd18 = FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[17]));
  1953.             final double cd14 = FastMath.cos(DAY_TO_RAD * (doy - p[13]));
  1954.             final double cd39 = FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[38]));

  1955.             // F10.7 effect
  1956.             t[0] = p[21] * (f107a - FLUX_REF);

  1957.             // Time independent
  1958.             t[1] = p[1]  * plg[0][2] + p[2]  * plg[0][4] + p[22] * plg[0][6] +
  1959.                    p[26] * plg[0][1] + p[14] * plg[0][3] + p[59] * plg[0][5];

  1960.             // Symmetrical annual
  1961.             t[2] = (p[18] + p[47] * plg[0][2] + p[29] * plg[0][4]) * cd32;

  1962.             // Symmetrical semiannual
  1963.             t[3] = (p[15] + p[16] * plg[0][2] + p[30] * plg[0][4]) * cd18;

  1964.             // Asymmetrical annual
  1965.             t[4] = (p[9] * plg[0][1] + p[10] * plg[0][3] + p[20] * plg[0][5]) * cd14;

  1966.             // Asymmetrical semiannual
  1967.             t[5] = (p[37] * plg[0][1]) * cd39;

  1968.             // Diurnal
  1969.             if (sw[7] != 0) {
  1970.                 final double t71 = p[11] * plg[1][2] * cd14 * swc[5];
  1971.                 final double t72 = p[12] * plg[1][2] * cd14 * swc[5];
  1972.                 t[6] = (p[3] * plg[1][1] + p[4] * plg[1][3] + t71) * ctloc +
  1973.                        (p[6] * plg[1][1] + p[7] * plg[1][3] + t72) * stloc;
  1974.             }

  1975.             // Semidiurnal
  1976.             if (sw[8] != 0) {
  1977.                 final double t81 = (p[23] * plg[2][3] + p[35] * plg[2][5]) * cd14 * swc[5];
  1978.                 final double t82 = (p[33] * plg[2][3] + p[36] * plg[2][5]) * cd14 * swc[5];
  1979.                 t[7] = (p[5] * plg[2][2] + p[41] * plg[2][4] + t81) * c2tloc +
  1980.                        (p[8] * plg[2][2] + p[42] * plg[2][4] + t82) * s2tloc;
  1981.             }

  1982.             // Terdiurnal
  1983.             if (sw[14] != 0) {
  1984.                 t[13] = p[39] * plg[3][3] * s3tloc + p[40] * plg[3][3] * c3tloc;
  1985.             }

  1986.             // Magnetic activity
  1987.             if (sw[9] == 1) {
  1988.                 t[8] = apdf * (p[32] + p[45] * plg[0][2] * swc[2]);
  1989.             } else if (sw[9] == -1) {
  1990.                 t[8] = apt  * (p[50] + p[96] * plg[0][2] * swc[2]);
  1991.             }

  1992.             // Longitudinal
  1993.             if (!(sw[10] == 0 || sw[11] == 0)) {
  1994.                 final double lonr   = DEG_TO_RAD * lon;
  1995.                 final SinCos scLonr = FastMath.sinCos(lonr);
  1996.                 t[10] = (1.0 + plg[0][1] * (p[80] * swc[5] * FastMath.cos(DAY_TO_RAD * (doy - p[81])) +
  1997.                                             p[85] * swc[6] * FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[86]))) +
  1998.                                p[83] * swc[3] * FastMath.cos(DAY_TO_RAD * (doy - p[84])) +
  1999.                                p[87] * swc[4] * FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[88]))) *
  2000.                         ((p[64] * plg[1][2] + p[65] * plg[1][4] + p[66] * plg[1][6] +
  2001.                           p[74] * plg[1][1] + p[75] * plg[1][3] + p[76] * plg[1][5]) * scLonr.cos() +
  2002.                          (p[90] * plg[1][2] + p[91] * plg[1][4] + p[92] * plg[1][6] +
  2003.                           p[77] * plg[1][1] + p[78] * plg[1][3] + p[79] * plg[1][5]) * scLonr.sin());
  2004.             }

  2005.             // Sum all effects
  2006.             double gl = 0;
  2007.             for (int i = 0; i < 14; i++) {
  2008.                 gl += FastMath.abs(sw[i + 1]) * t[i];
  2009.             }

  2010.             // Return G(L)
  2011.             return gl;
  2012.         }

  2013.         /** Implements sg0 function (Eq. A24a).
  2014.          * @param ex ex
  2015.          * @param p24 abs(p[24])
  2016.          * @param p25 p[25]
  2017.          * @return sg0
  2018.          */
  2019.         private double sg0(final double ex, final double p24, final double p25) {
  2020.             final double g01 = g0(ap[1], p24, p25);
  2021.             final double g02 = g0(ap[2], p24, p25);
  2022.             final double g03 = g0(ap[3], p24, p25);
  2023.             final double g04 = g0(ap[4], p24, p25);
  2024.             final double g05 = g0(ap[5], p24, p25);
  2025.             final double g06 = g0(ap[6], p24, p25);
  2026.             final double ex2 = ex * ex;
  2027.             final double ex3 = ex * ex2;
  2028.             final double ex4 = ex2 * ex2;
  2029.             final double ex8 = ex4 * ex4;
  2030.             final double ex12 = ex4 * ex8;
  2031.             final double g234 = g02 * ex + g03 * ex2 + g04 * ex3;
  2032.             final double g56  = g05 * ex4 + g06 * ex12;
  2033.             final double ex19 = ex3 * ex4 * ex12;
  2034.             final double omex = 1.0 - ex;
  2035.             final double sumex = 1.0 + (1.0 - ex19) / omex * FastMath.sqrt(ex);
  2036.             return (g01 + (g234 + g56 * (1.0 - ex8) / omex)) / sumex;
  2037.         }

  2038.         /** Implements go function (Eq. A24d).
  2039.          * @param apI 3 hrs ap
  2040.          * @param p24 abs(p[24])
  2041.          * @param p25 p[25]
  2042.          * @return go
  2043.          */
  2044.         private double g0(final double apI, final double p24, final double p25) {
  2045.             final double am4 = apI - 4.0;
  2046.             return am4 + (p25 - 1.0) * (am4 + (FastMath.exp(-p24 * am4) - 1.0) / p24);
  2047.         }

  2048.         /** Calculates chemistry/dissociation correction for MSIS models.
  2049.          * @param alt altitude
  2050.          * @param r target ratio
  2051.          * @param h1 transition scale length
  2052.          * @param zh altitude of 1/2 R
  2053.          * @return correction
  2054.          */
  2055.         private double ccor(final double alt, final double r, final double h1, final double zh) {
  2056.             final double e = (alt - zh) / h1;
  2057.             if (e > 70.) {
  2058.                 return 1.;
  2059.             } else if (e < -70.) {
  2060.                 return FastMath.exp(r);
  2061.             } else {
  2062.                 return FastMath.exp(r / (1.0 + FastMath.exp(e)));
  2063.             }
  2064.         }


  2065.         /** Calculates O & O2 chemistry/dissociation correction for MSIS models.
  2066.          * @param alt altitude
  2067.          * @param r target ratio
  2068.          * @param h1 transition scale length
  2069.          * @param zh altitude of 1/2 R
  2070.          * @param h2 transition scale length
  2071.          * @return correction
  2072.          */
  2073.         private double ccor2(final double alt, final double r,
  2074.                              final double h1, final double zh, final double h2) {
  2075.             final double e1 = (alt - zh) / h1;
  2076.             final double e2 = (alt - zh) / h2;
  2077.             if (e1 > 70. || e2 > 70.) {
  2078.                 return 1.;
  2079.             } else if (e1 < -70. && e2 < -70.) {
  2080.                 return FastMath.exp(r);
  2081.             } else {
  2082.                 final double ex1 = FastMath.exp(e1);
  2083.                 final double ex2 = FastMath.exp(e2);
  2084.                 return FastMath.exp(r / (1.0 + 0.5 * (ex1 + ex2)));
  2085.             }
  2086.         }

  2087.         /** Calculates scale height.
  2088.          * @param alt altitude
  2089.          * @param xm species molecular weight
  2090.          * @param temp temperature
  2091.          * @return scale height (km)
  2092.          */
  2093.         private double scalh(final double alt, final double xm, final double temp) {
  2094.             // Gravity at altitude
  2095.             final double denom = 1.0 + alt / rlat;
  2096.             final double galt = glat / (denom * denom);
  2097.             return R_GAS * temp / (galt * xm);
  2098.         }

  2099.         /** Calculates turbopause correction for MSIS models.
  2100.          * @param dd diffusive density
  2101.          * @param dm full mixed density
  2102.          * @param zhm transition scale length
  2103.          * @param xmm full mixed molecular weight
  2104.          * @param xm species molecular weight
  2105.          * @return combined density
  2106.          */
  2107.         private double dnet(final double dd, final double dm,
  2108.                             final double zhm, final double xmm, final double xm) {
  2109.             if (!(dm > 0 && dd > 0)) {
  2110.                 double ddd = dd;
  2111.                 if (dd == 0 && dm == 0) {
  2112.                     ddd = 1;
  2113.                 }
  2114.                 if (dm == 0) {
  2115.                     return ddd;
  2116.                 }
  2117.                 if (dd == 0) {
  2118.                     return dm;
  2119.                 }
  2120.             }

  2121.             final double a  = zhm / (xmm - xm);
  2122.             final double ylog = a * FastMath.log(dm / dd);
  2123.             if (ylog < -10.) {
  2124.                 return dd;
  2125.             } else if (ylog > 10.) {
  2126.                 return dm;
  2127.             } else {
  2128.                 return dd * FastMath.pow(1.0 + FastMath.exp(ylog), 1.0 / a);
  2129.             }
  2130.         }

  2131.         /** Integrate cubic spline function from xa[0] to x.
  2132.          * <p>ADAPTED FROM NUMERICAL RECIPES</p>
  2133.          * @param xa array of abscissas in ascending order
  2134.          * @param ya array of ordinates in ascending order by xa
  2135.          * @param y2a array of second derivatives in ascending order by xa
  2136.          * @param x abscissa end point
  2137.          * @return integral value
  2138.          */
  2139.         private double splini(final double[] xa, final double[] ya, final double[] y2a, final double x) {
  2140.             final int n = xa.length;
  2141.             double yi = 0;
  2142.             int klo = 0;
  2143.             int khi = 1;
  2144.             while (x > xa[klo] && khi < n) {
  2145.                 double xx = x;
  2146.                 if (khi < n - 1) {
  2147.                     xx = (x < xa[khi]) ? x : xa[khi];
  2148.                 }
  2149.                 final double h = xa[khi] - xa[klo];
  2150.                 final double a = (xa[khi] - xx) / h;
  2151.                 final double b = (xx - xa[klo]) / h;
  2152.                 final double a2 = a * a;
  2153.                 final double b2 = b * b;
  2154.                 yi += ((1.0 - a2) * ya[klo] / 2.0 + b2 * ya[khi] / 2.0 +
  2155.                        ((-(1.0 + a2 * a2) / 4.0 + a2 / 2.0) * y2a[klo] +
  2156.                           (b2 * b2 / 4.0 - b2 / 2.0) * y2a[khi]) * h * h / 6.0) * h;
  2157.                 klo++;
  2158.                 khi++;
  2159.             }
  2160.             return yi;
  2161.         }

  2162.         /** Calculate cubic spline interpolated value.
  2163.          * <p>ADAPTED FROM NUMERICAL RECIPES</p>
  2164.          * @param xa array of abscissas in ascending order
  2165.          * @param ya array of ordinates in ascending order by xa
  2166.          * @param y2a array of second derivatives in ascending order by xa
  2167.          * @param x abscissa for interpolation
  2168.          * @return interpolated value
  2169.          */
  2170.         private double splint(final double[] xa, final double[] ya, final double[] y2a, final double x) {
  2171.             final int n = xa.length;
  2172.             int klo = 0;
  2173.             int khi = n - 1;
  2174.             while (khi - klo > 1) {
  2175.                 final int k = (khi + klo) >>> 1;
  2176.                 if (xa[k] > x) {
  2177.                     khi = k;
  2178.                 } else {
  2179.                     klo = k;
  2180.                 }
  2181.             }
  2182.             final double h = xa[khi] - xa[klo];
  2183.             final double a = (xa[khi] - x) / h;
  2184.             final double b = (x - xa[klo]) / h;
  2185.             return a * ya[klo] + b * ya[khi] +
  2186.                     ((a * a * a - a) * y2a[klo] + (b * b * b - b) * y2a[khi]) * h * h / 6.0;
  2187.         }

  2188.         /** Calculate 2nd derivatives of cubic spline interpolation function.
  2189.          * <p>ADAPTED FROM NUMERICAL RECIPES</p>
  2190.          * @param x array of abscissas in ascending order
  2191.          * @param y array of ordinates in ascending order by x
  2192.          * @param yp1 derivative at x[0] (2nd derivatives null if > 1E30)
  2193.          * @param ypn derivative at x[n-1] (2nd derivatives null if > 1E30)
  2194.          * @return array of second derivatives
  2195.          */
  2196.         private double[] spline(final double[] x, final double[] y, final double yp1, final double ypn) {
  2197.             final int n = x.length;
  2198.             final double[] y2 = new double[n];
  2199.             final double[] u  = new double[n];

  2200.             if (yp1 < 1e+30) {
  2201.                 y2[0] = -0.5;
  2202.                 u[0]  = (3.0 / (x[1] - x[0])) * ((y[1] - y[0]) / (x[1] - x[0]) - yp1);
  2203.             }
  2204.             for (int i = 1; i < n - 1; i++) {
  2205.                 final double sig = (x[i] - x[i - 1]) / (x[i + 1] - x[i - 1]);
  2206.                 final double p = sig * y2[i - 1] + 2.0;
  2207.                 y2[i] = (sig - 1.0) / p;
  2208.                 u[i] = (6.0 * ((y[i + 1] - y[i]) / (x[i + 1] - x[i]) - (y[i] - y[i - 1]) / (x[i] - x[i - 1])) /
  2209.                         (x[i + 1] - x[i - 1]) - sig * u[i - 1]) / p;
  2210.             }

  2211.             double qn = 0;
  2212.             double un = 0;
  2213.             if (ypn < 1e+30) {
  2214.                 qn = 0.5;
  2215.                 un = (3.0 / (x[n - 1] - x[n - 2])) * (ypn - (y[n - 1] - y[n - 2]) / (x[n - 1] - x[n - 2]));
  2216.             }

  2217.             y2[n - 1] = (un - qn * u[n - 2]) / (qn * y2[n - 2] + 1.0);
  2218.             for (int k = n - 2; k >= 0; k--) {
  2219.                 y2[k] = y2[k] * y2[k + 1] + u[k];
  2220.             }

  2221.             return y2;
  2222.         }

  2223.         /** Calculate Temperature and Density Profiles for lower atmosphere.
  2224.          * @param alt altitude
  2225.          * @param d0 density
  2226.          * @param xm mixed density
  2227.          * @return temperature or density profile
  2228.          */
  2229.         private double densm(final double alt, final double d0, final double xm) {

  2230.             double densm = d0;

  2231.             // stratosphere/mesosphere temperature
  2232.             int mn = ZN2.length;
  2233.             double z = (alt > ZN2[mn - 1]) ? alt : ZN2[mn - 1];

  2234.             double z1 = ZN2[0];
  2235.             double z2 = ZN2[mn - 1];
  2236.             double t1 = meso_tn2[0];
  2237.             double t2 = meso_tn2[mn - 1];
  2238.             double zg  = zeta(z, z1);
  2239.             double zgdif = zeta(z2, z1);

  2240.             /* set up spline nodes */
  2241.             double[] xs = new double[mn];
  2242.             double[] ys = new double[mn];
  2243.             for (int k = 0; k < mn; k++) {
  2244.                 xs[k] = zeta(ZN2[k], z1) / zgdif;
  2245.                 ys[k] = 1.0 / meso_tn2[k];
  2246.             }
  2247.             final double qSM = (rlat + z2) / (rlat + z1);
  2248.             double yd1 = -meso_tgn2[0] / (t1 * t1) * zgdif;
  2249.             double yd2 = -meso_tgn2[1] / (t2 * t2) * zgdif * qSM * qSM;

  2250.             /* calculate spline coefficients */
  2251.             double[] y2out = spline(xs, ys, yd1, yd2);
  2252.             double x = zg / zgdif;
  2253.             double y = splint(xs, ys, y2out, x);

  2254.             /* temperature at altitude */
  2255.             double tz = 1.0 / y;

  2256.             if (xm != 0.0) {
  2257.                 /* calculate stratosphere / mesospehere density */
  2258.                 final double glb  = galt(z1);
  2259.                 final double gamm = xm * glb * zgdif / R_GAS;

  2260.                 /* Integrate temperature profile */
  2261.                 final double yi = splini(xs, ys, y2out, x);
  2262.                 final double expl = FastMath.min(MIN_TEMP, gamm * yi);

  2263.                 /* Density at altitude */
  2264.                 densm *= (t1 / tz) * FastMath.exp(-expl);
  2265.             }

  2266.             if (alt > ZN3[0]) {
  2267.                 return (xm == 0.0) ? tz : densm;
  2268.             }

  2269.             // troposhere/stratosphere temperature
  2270.             z = alt;
  2271.             mn = ZN3.length;
  2272.             z1 = ZN3[0];
  2273.             z2 = ZN3[mn - 1];
  2274.             t1 = meso_tn3[0];
  2275.             t2 = meso_tn3[mn - 1];
  2276.             zg = zeta(z, z1);
  2277.             zgdif = zeta(z2, z1);

  2278.             /* set up spline nodes */
  2279.             xs = new double[mn];
  2280.             ys = new double[mn];
  2281.             for (int k = 0; k < mn; k++) {
  2282.                 xs[k] = zeta(ZN3[k], z1) / zgdif;
  2283.                 ys[k] = 1.0 / meso_tn3[k];
  2284.             }
  2285.             final double qTS = (rlat + z2) / (rlat + z1);
  2286.             yd1 = -meso_tgn3[0] / (t1 * t1) * zgdif;
  2287.             yd2 = -meso_tgn3[1] / (t2 * t2) * zgdif * qTS * qTS;

  2288.             /* calculate spline coefficients */
  2289.             y2out = spline(xs, ys, yd1, yd2);
  2290.             x = zg / zgdif;
  2291.             y = splint(xs, ys, y2out, x);

  2292.             /* temperature at altitude */
  2293.             tz = 1.0 / y;

  2294.             if (xm != 0.0) {
  2295.                 /* calculate tropospheric / stratosphere density */
  2296.                 final double glb = galt(z1);
  2297.                 final double gamm = xm * glb * zgdif / R_GAS;

  2298.                 /* Integrate temperature profile */
  2299.                 final double yi = splini(xs, ys, y2out, x);
  2300.                 final double expl = FastMath.min(MIN_TEMP, gamm * yi);

  2301.                 /* Density at altitude */
  2302.                 densm *= (t1 / tz) * FastMath.exp(-expl);
  2303.             }

  2304.             return (xm == 0.0) ? tz : densm;
  2305.         }

  2306.         /** Calculate temperature and density profiles according to new lower thermo polynomial.
  2307.          * @param alt altitude
  2308.          * @param dlb density at lower boundary
  2309.          * @param tinf exospheric temperature
  2310.          * @param tlb temperature at lower boundary
  2311.          * @param xm species molecular weight
  2312.          * @param alpha thermal diffusion coefficient
  2313.          * @param zlb altitude of the lower boundary
  2314.          * @param s2 slope
  2315.          * @return temperature or density profile
  2316.          */
  2317.         private double densu(final double alt, final double dlb, final double tinf,
  2318.                              final double tlb, final double xm, final double alpha,
  2319.                              final double zlb, final double s2) {
  2320.             /* joining altitudes of Bates and spline */
  2321.             double z = (alt > ZN1[0]) ? alt : ZN1[0];

  2322.             /* geopotential altitude difference from ZLB */
  2323.             final double zg2 = zeta(z, zlb);

  2324.             /* Bates temperature */
  2325.             final double tt = tinf - (tinf - tlb) * FastMath.exp(-s2 * zg2);
  2326.             final double ta = tt;
  2327.             double tz = tt;

  2328.             final int mn = ZN1.length;
  2329.             final double[] xs = new double[mn];
  2330.             final double[] ys = new double[mn];
  2331.             double x = 0.;
  2332.             double[] y2out =  new double[mn];
  2333.             double zgdif = 0.;
  2334.             if (alt < ZN1[0]) {
  2335.                 /* calculate temperature below ZA
  2336.                  * temperature gradient at ZA from Bates profile */
  2337.                 final double p = (rlat + zlb) / (rlat + ZN1[0]);
  2338.                 final double dta = (tinf - ta) * s2 * p * p;
  2339.                 meso_tgn1[0] = dta;
  2340.                 meso_tn1[0] = ta;
  2341.                 z = (alt > ZN1[mn - 1]) ? alt : ZN1[mn - 1];

  2342.                 final double t1 = meso_tn1[0];
  2343.                 final double t2 = meso_tn1[mn - 1];
  2344.                 /* geopotental difference from z1 */
  2345.                 final double zg = zeta(z, ZN1[0]);
  2346.                 zgdif = zeta(ZN1[mn - 1], ZN1[0]);
  2347.                 /* set up spline nodes */
  2348.                 for (int k = 0; k < mn; k++) {
  2349.                     xs[k] = zeta(ZN1[k], ZN1[0]) / zgdif;
  2350.                     ys[k] = 1.0 / meso_tn1[k];
  2351.                 }
  2352.                 /* end node derivatives */
  2353.                 final double q   = (rlat + ZN1[mn - 1]) / (rlat + ZN1[0]);
  2354.                 final double yd1 = -meso_tgn1[0] / (t1 * t1) * zgdif;
  2355.                 final double yd2 = -meso_tgn1[1] / (t2 * t2) * zgdif * q * q;
  2356.                 /* calculate spline coefficients */
  2357.                 y2out = spline(xs, ys, yd1, yd2);
  2358.                 x = zg / zgdif;
  2359.                 final double y = splint(xs, ys, y2out, x);
  2360.                 /* temperature at altitude */
  2361.                 tz = 1.0 / y;
  2362.             }

  2363.             if (xm == 0) {
  2364.                 return tz;
  2365.             }

  2366.             /* calculate density above za */
  2367.             double glb   = galt(zlb);
  2368.             double gamma = xm * glb / (R_GAS * s2 * tinf);
  2369.             double expl  = (tt <= 0) ? MIN_TEMP : FastMath.min(MIN_TEMP, FastMath.exp(-s2 * gamma * zg2));
  2370.             double densu = dlb * expl * FastMath.pow(tlb / tt, 1.0 + alpha + gamma);

  2371.             // Correction for issue 1365 - protection against "densu" being infinite
  2372.             if (!Double.isFinite(densu)) {
  2373.                 if (expl < MIN_TEMP) {
  2374.                     densu = dlb * FastMath.exp(FastMath.log(tlb / tt) * (1.0 + alpha + gamma) - s2 * gamma * zg2);
  2375.                 } else {
  2376.                     throw new OrekitException( OrekitMessages.INFINITE_NRLMSISE00_DENSITY);
  2377.                 }
  2378.             }

  2379.             /* calculate density below za */
  2380.             if (alt < ZN1[0]) {
  2381.                 glb   = galt(ZN1[0]);
  2382.                 gamma = xm * glb * zgdif / R_GAS;
  2383.                 /* integrate spline temperatures */
  2384.                 expl  = (tz <= 0) ? MIN_TEMP : FastMath.min(MIN_TEMP, gamma * splini(xs, ys, y2out, x));
  2385.                 /* correct density at altitude */
  2386.                 densu *= FastMath.pow(meso_tn1[0] / tz, 1.0 + alpha) * FastMath.exp(-expl);
  2387.             }

  2388.             /* Return density at altitude */
  2389.             return densu;
  2390.         }

  2391.         /** Calculate gravity at altitude.
  2392.          * @param alt altitude (km)
  2393.          * @return gravity at altitude (cm/s2)
  2394.          */
  2395.         private double galt(final double alt) {
  2396.             final double r = 1.0 + alt / rlat;
  2397.             return glat / (r * r);
  2398.         }

  2399.         /** Calculate zeta function.
  2400.          * @param zz zz value
  2401.          * @param zl zl value
  2402.          * @return value of zeta function
  2403.          */
  2404.         private double zeta(final double zz, final double zl) {
  2405.             return (zz - zl) * (rlat + zl) / (rlat + zz);
  2406.         }

  2407.     }

  2408.     /**
  2409.      * This class is a placeholder for the computed densities and temperatures.
  2410.      * <p>
  2411.      * Densities are provided as an array d such as:
  2412.      * <ul>
  2413.      * <li>d[0] = He number density (1/m³)</li>
  2414.      * <li>d[1] = O number density (1/m³)</li>
  2415.      * <li>d[2] = N2 number density (1/m³)</li>
  2416.      * <li>d[3] = O2 number density (1/m³)</li>
  2417.      * <li>d[4] = Ar number density (1/m³)</li>
  2418.      * <li>d[5] = total mass density (kg/m³) (*)</li>
  2419.      * <li>d[6] = H number density (1/m³)</li>
  2420.      * <li>d[7] = N number density (1/m³)</li>
  2421.      * <li>d[8] = anomalous oxygen number density (1/m³)
  2422.      * </ul>
  2423.      * Total mass density, d[5], is NOT the same for methods gtd7 and gtd7d:
  2424.      * <ul>
  2425.      * <li>For gtd7: d[5] is the sum of the mass densities of the species
  2426.      * He, O, N2, O2, Ar, H and N but does NOT include anomalous oxygen.</li>
  2427.      * <li>For gtd7d: d[5] is the "effective total mass density for drag" and is the sum
  2428.      * of the mass densities of all species in this model, INCLUDING anomalous oxygen.</li>
  2429.      * </ul>
  2430.      * O, H, and N are set to zero below 72.5 km.
  2431.      * <p>
  2432.      * Temperatures are provided as an array t such as:
  2433.      * <ul>
  2434.      * <li>t[0] = exospheric temperature (K)</li>
  2435.      * <li>t[1] = temperature at altitude (K)</li>
  2436.      * </ul>
  2437.      * <p>
  2438.      * t[0] is set to global average for altitudes below 120 km.<br>
  2439.      * The 120 km gradient is left at global average value for altitudes below 72 km.
  2440.      * </p>
  2441.      * @param <T> type of the field elements
  2442.      * @since 9.0
  2443.      */
  2444.     public class FieldOutput<T extends CalculusFieldElement<T>> {

  2445.         /** Type of the field elements. */
  2446.         private final Field<T> field;

  2447.         /** Zero for the field. */
  2448.         private final T zero;

  2449.         /** Day of year (from 1 to 365 or 366). */
  2450.         private final int doy;

  2451.         /** Seconds in day (UT scale). */
  2452.         private final T sec;

  2453.         /** Geodetic latitude (°). */
  2454.         private final T lat;

  2455.         /** Geodetic longitude (°). */
  2456.         private final T lon;

  2457.         /** Local apparent solar time (hours). */
  2458.         private final T hl;

  2459.         /** 81 day average of F10.7 flux (centered on day). */
  2460.         private final double f107a;

  2461.         /** Daily F10.7 flux for previous day. */
  2462.         private final double f107;

  2463.         /** Array containing:
  2464.         *  <ul>
  2465.         *  <li>0: daily Ap</li>
  2466.         *  <li>1: 3 hr ap index for current time</li>
  2467.         *  <li>2: 3 hr ap index for 3 hrs before current time</li>
  2468.         *  <li>3: 3 hr ap index for 6 hrs before current time</li>
  2469.         *  <li>4: 3 hr ap index for FOR 9 hrs before current time</li>
  2470.         *  <li>5: average of eight 3 hr ap indices from 12 to 33 hrs prior to current time</li>
  2471.         *  <li>6: average of eight 3 hr ap indices from 36 to 57 hrs prior to current time</li>
  2472.         *  </ul>. */
  2473.         private final double[] ap;

  2474.         /** Gravity at latitude (cm/s2). */
  2475.         private final T glat;

  2476.         /** Effective Earth radius at latitude (km). */
  2477.         private final T rlat;

  2478.         /** N2 mixed density at alt. */
  2479.         private T dm28;

  2480.         /** Legendre polynomials. */
  2481.         private final T[][] plg;

  2482.         /** Cosinus of local solar time. */
  2483.         private final T ctloc;
  2484.         /** Sinus of local solar time. */
  2485.         private final T stloc;
  2486.         /** Square of ctloc. */
  2487.         private final T c2tloc;
  2488.         /** Square of stloc. */
  2489.         private final T s2tloc;
  2490.         /** Cube of ctloc. */
  2491.         private final T c3tloc;
  2492.         /** Cube of stloc. */
  2493.         private final T s3tloc;

  2494.         /** Magnetic activity based on daily ap. */
  2495.         private double apdf;

  2496.         /** Magnetic activity based on daily ap. */
  2497.         private T apt;

  2498.         /** Temperature at nodes for ZN1 scale. */
  2499.         private final T[] meso_tn1;

  2500.         /** Temperature at nodes for ZN2 scale. */
  2501.         private final T[] meso_tn2;

  2502.         /** Temperature at nodes for ZN3 scale. */
  2503.         private final T[] meso_tn3;

  2504.         /** Temperature gradients at end nodes for ZN1 scale. */
  2505.         private final T[] meso_tgn1;

  2506.         /** Temperature gradients at end nodes for ZN2 scale. */
  2507.         private final T[] meso_tgn2;

  2508.         /** Temperature gradients at end nodes for ZN3 scale. */
  2509.         private final T[] meso_tgn3;

  2510.         /** Densities. */
  2511.         private final T[] densities;

  2512.         /** Temperatures. */
  2513.         private final T[] temperatures;

  2514.         /** Simple constructor.
  2515.          *  @param doy day of year (from 1 to 365 or 366)
  2516.          *  @param sec seconds in day (UT scale)
  2517.          *  @param lat geodetic latitude (°)
  2518.          *  @param lon geodetic longitude (°)
  2519.          *  @param hl local apparent solar time (hours)
  2520.          *  @param f107a 81 day average of F10.7 flux (centered on day)
  2521.          *  @param f107 daily F10.7 flux for previous day
  2522.          *  @param ap array containing:
  2523.          *  <ul>
  2524.          *  <li>0: daily Ap</li>
  2525.          *  <li>1: 3 hr ap index for current time</li>
  2526.          *  <li>2: 3 hr ap index for 3 hrs before current time</li>
  2527.          *  <li>3: 3 hr ap index for 6 hrs before current time</li>
  2528.          *  <li>4: 3 hr ap index for FOR 9 hrs before current time</li>
  2529.          *  <li>5: average of eight 3 hr ap indices from 12 to 33 hrs prior to current time</li>
  2530.          *  <li>6: average of eight 3 hr ap indices from 36 to 57 hrs prior to current time</li>
  2531.          *  </ul>
  2532.          */
  2533.         FieldOutput(final int doy, final T sec,
  2534.                     final T lat, final T lon, final T hl,
  2535.                     final double f107a, final double f107, final double[] ap) {

  2536.             this.field = sec.getField();
  2537.             this.zero = field.getZero();

  2538.             this.doy   = doy;
  2539.             this.sec   = sec;
  2540.             this.lat   = lat;
  2541.             this.lon   = lon;
  2542.             this.hl    = hl;
  2543.             this.f107a = f107a;
  2544.             this.f107  = f107;
  2545.             this.ap    = ap.clone();

  2546.             this.plg       = MathArrays.buildArray(field, 4, 8);

  2547.             this.meso_tn1  = MathArrays.buildArray(field, ZN1.length);
  2548.             this.meso_tn2  = MathArrays.buildArray(field, ZN2.length);
  2549.             this.meso_tn3  = MathArrays.buildArray(field, ZN3.length);
  2550.             this.meso_tgn1 = MathArrays.buildArray(field, 2);
  2551.             this.meso_tgn2 = MathArrays.buildArray(field, 2);
  2552.             this.meso_tgn3 = MathArrays.buildArray(field, 2);

  2553.             densities       = MathArrays.buildArray(field, 9);
  2554.             temperatures    = MathArrays.buildArray(field, 2);

  2555.             // Calculates latitude variable gravity and effective radius
  2556.             final T xlat = (sw[2] == 0) ? zero.newInstance(LAT_REF) : lat;
  2557.             final T c2   = xlat.multiply(2 * DEG_TO_RAD).cos();
  2558.             glat = c2.multiply(-0.0026373).add(1).multiply(G_REF);
  2559.             rlat = glat.multiply(2).divide(c2.multiply(2.27e-9).add(3.085462e-6)).multiply(1.e-5);

  2560.             // Convert latitude into radians
  2561.             final T latr = lat.multiply(DEG_TO_RAD);

  2562.             // Calculate legendre polynomials
  2563.             final FieldSinCos<T> scLatr = FastMath.sinCos(latr);
  2564.             final T c = scLatr.sin();
  2565.             final T s = scLatr.cos();

  2566.             plg[0][1] = c;
  2567.             plg[0][2] = c.multiply( 3.0).multiply(plg[0][1]).subtract(1.0).divide(2.0);
  2568.             plg[0][3] = c.multiply( 5.0).multiply(plg[0][2]).subtract(plg[0][1].multiply(2.0)).divide(3.0);
  2569.             plg[0][4] = c.multiply( 7.0).multiply(plg[0][3]).subtract(plg[0][2].multiply(3.0)).divide(4.0);
  2570.             plg[0][5] = c.multiply( 9.0).multiply(plg[0][4]).subtract(plg[0][3].multiply(4.0)).divide(5.0);
  2571.             plg[0][6] = c.multiply(11.0).multiply(plg[0][5]).subtract(plg[0][4].multiply(5.0)).divide(6.0);

  2572.             plg[1][1] = s;
  2573.             plg[1][2] = c.multiply( 3.0).multiply(plg[1][1]);
  2574.             plg[1][3] = c.multiply( 5.0).multiply(plg[1][2]).subtract(plg[1][1].multiply(3.0)).divide(2.0);
  2575.             plg[1][4] = c.multiply( 7.0).multiply(plg[1][3]).subtract(plg[1][2].multiply(4.0)).divide(3.0);
  2576.             plg[1][5] = c.multiply( 9.0).multiply(plg[1][4]).subtract(plg[1][3].multiply(5.0)).divide(4.0);
  2577.             plg[1][6] = c.multiply(11.0).multiply(plg[1][5]).subtract(plg[1][4].multiply(6.0)).divide(5.0);

  2578.             plg[2][2] = s.multiply( 3.0).multiply(plg[1][1]);
  2579.             plg[2][3] = c.multiply( 5.0).multiply(plg[2][2]);
  2580.             plg[2][4] = c.multiply( 7.0).multiply(plg[2][3]).subtract(plg[2][2].multiply(5.0)).divide(2.0);
  2581.             plg[2][5] = c.multiply( 9.0).multiply(plg[2][4]).subtract(plg[2][3].multiply(6.0)).divide(3.0);
  2582.             plg[2][6] = c.multiply(11.0).multiply(plg[2][5]).subtract(plg[2][4].multiply(7.0)).divide(4.0);
  2583.             plg[2][7] = c.multiply(13.0).multiply(plg[2][6]).subtract(plg[2][5].multiply(8.0)).divide(5.0);

  2584.             plg[3][3] = s.multiply( 5.0).multiply(plg[2][2]);
  2585.             plg[3][4] = c.multiply( 7.0).multiply(plg[3][3]);
  2586.             plg[3][5] = c.multiply( 9.0).multiply(plg[3][4]).subtract(plg[3][3].multiply(7.0)).divide(2.0);
  2587.             plg[3][6] = c.multiply(11.0).multiply(plg[3][5]).subtract(plg[3][4].multiply(8.0)).divide(3.0);

  2588.             // Calculate additional data
  2589.             if (!(sw[7] == 0 && sw[8] == 0 && sw[14] == 0)) {
  2590.                 final T tloc = hl.multiply(HOUR_TO_RAD);
  2591.                 final FieldSinCos<T> sc  = FastMath.sinCos(tloc);
  2592.                 final FieldSinCos<T> sc2 = FieldSinCos.sum(sc, sc);
  2593.                 final FieldSinCos<T> sc3 = FieldSinCos.sum(sc, sc2);
  2594.                 stloc  = sc.sin();
  2595.                 ctloc  = sc.cos();
  2596.                 s2tloc = sc2.sin();
  2597.                 c2tloc = sc2.cos();
  2598.                 s3tloc = sc3.sin();
  2599.                 c3tloc = sc3.cos();
  2600.             } else {
  2601.                 stloc  = zero;
  2602.                 ctloc  = zero;
  2603.                 s2tloc = zero;
  2604.                 c2tloc = zero;
  2605.                 s3tloc = zero;
  2606.                 c3tloc = zero;
  2607.             }

  2608.         }

  2609.         /** Calculate temperatures and densities not including anomalous oxygen.
  2610.          *  <p>
  2611.          *  This method is the thermospheric portion of NRLMSISE-00 for alt > 72.5 km.
  2612.          *  </p>
  2613.          *  <p>NOTES ON INPUT VARIABLES:<br>
  2614.          *  Seconds, Local Time, and Longitude are used independently in the
  2615.          *  model and are not of equal importance for every situation.<br>
  2616.          *  For the most physically realistic calculation these three
  2617.          *  variables should be consistent (lst=sec/3600 + lon/15).<br>
  2618.          *  The Equation of Time departures from the above formula
  2619.          *  for apparent local time can be included if available but
  2620.          *  are of minor importance.<br><br>
  2621.          *
  2622.          *  f107 and f107A values used to generate the model correspond
  2623.          *  to the 10.7 cm radio flux at the actual distance of the Earth
  2624.          *  from the Sun rather than the radio flux at 1 AU. The following
  2625.          *  site provides both classes of values:<br>
  2626.          *  ftp://ftp.ngdc.noaa.gov/STP/SOLAR_DATA/SOLAR_RADIO/FLUX/<br><br>
  2627.          *
  2628.          *  f107, f107A, and ap effects are neither large nor well established below 80 km
  2629.          *  and these parameters should be set to 150., 150., and 4. respectively.
  2630.          *  </p>
  2631.          *  @param alt altitude (km)
  2632.          */
  2633.         void gts7(final T alt) {

  2634.             // Thermal diffusion coefficients for species
  2635.             final double[] alpha = {-0.38, 0.0, 0.0, 0.0, 0.17, 0.0, -0.38, 0.0, 0.0};
  2636.             // Altitude limits for net density computation for species
  2637.             final double[] altl  = {200.0, 300.0, 160.0, 250.0, 240.0, 450.0, 320.0, 450.0};
  2638.             // N2 mixed density
  2639.             final double xmm = PDM[2][4];

  2640.             /**** Exospheric temperature ****/
  2641.             T tinf = zero.newInstance(PTM[0] * PT[0]);
  2642.             // Tinf variations not important below ZA or ZN[0]
  2643.             if (alt.getReal() > ZN1[0]) {
  2644.                 tinf = tinf.multiply(globe7(PT).multiply(sw[16]).add(1));
  2645.             }
  2646.             setTemperature(EXOSPHERIC, tinf);

  2647.             // Gradient variations not important below ZN[4]
  2648.             T g0 = zero.newInstance(PTM[3] * PS[0]);
  2649.             if (alt.getReal() > ZN1[4]) {
  2650.                 g0 = g0.multiply(globe7(PS).multiply(sw[19]).add(1));
  2651.             }

  2652.             // Temperature at lower boundary
  2653.             T tlb = zero.newInstance(PTM[1] * PD[3][0]);
  2654.             tlb = tlb.multiply(globe7(PD[3]).multiply(sw[17]).add(1));

  2655.             // Slope
  2656.             final T s = g0.divide(tinf.subtract(tlb));

  2657.             // Lower thermosphere temp variations not significant for density above 300 km
  2658.             meso_tn1[1]  = zero.newInstance(PTM[6] * PTL[0][0]);
  2659.             meso_tn1[2]  = zero.newInstance(PTM[2] * PTL[1][0]);
  2660.             meso_tn1[3]  = zero.newInstance(PTM[7] * PTL[2][0]);
  2661.             meso_tn1[4]  = zero.newInstance(PTM[4] * PTL[3][0]);
  2662.             meso_tgn1[1] = zero.newInstance(PTM[8] * PMA[8][0]);
  2663.             if (alt.getReal() < 300.0) {
  2664.                 final double r = PTM[4] * PTL[3][0];
  2665.                 meso_tn1[1]  =  meso_tn1[1].divide(glob7s(PTL[0]).multiply(sw[18]         ).negate().add(1));
  2666.                 meso_tn1[2]  =  meso_tn1[2].divide(glob7s(PTL[1]).multiply(sw[18]         ).negate().add(1));
  2667.                 meso_tn1[3]  =  meso_tn1[3].divide(glob7s(PTL[2]).multiply(sw[18]         ).negate().add(1));
  2668.                 meso_tn1[4]  =  meso_tn1[4].divide(glob7s(PTL[3]).multiply(sw[18] * sw[20]).negate().add(1));
  2669.                 meso_tgn1[1] =  meso_tgn1[1].multiply(glob7s(PMA[8]).multiply(sw[18] * sw[20]).add(1));
  2670.                 meso_tgn1[1] =  meso_tgn1[1].multiply(meso_tn1[4].multiply(meso_tn1[4]).divide(r * r));
  2671.             }

  2672.             /**** Temperature at altitude ****/
  2673.             setTemperature(ALTITUDE, densu(alt, zero.newInstance(1.0), tinf, tlb, 0, 0, PTM[5], s));

  2674.             /**** N2 density ****/
  2675.             /*   Density variation factor at Zlb */
  2676.             final T g28 = globe7(PD[2]).multiply(sw[21]);
  2677.             /* Diffusive density at Zlb */
  2678.             final T db28 = g28.exp().multiply(PDM[2][0] * PD[2][0]);
  2679.             /* Diffusive density at Alt */
  2680.             T diffusiveDensity = densu(alt, db28, tinf, tlb, N2_MASS, alpha[2], PTM[5], s);
  2681.             setDensity(MOLECULAR_NITROGEN, diffusiveDensity);
  2682.             // Variation of turbopause height
  2683.             final T zhf = lat.multiply(DEG_TO_RAD).sin().
  2684.                             multiply(sw[5] * PDL[0][24] * FastMath.cos(DAY_TO_RAD * (doy - PT[13]))).
  2685.                             add(1).
  2686.                             multiply(PDL[1][24]);
  2687.             /* Turbopause */
  2688.             final T zh28  = zhf.multiply(PDM[2][2]);
  2689.             final double zhm28 = PDM[2][3] * PDL[1][5];
  2690.             /* Mixed density at Zlb */
  2691.             final T b28 = densu(zh28, db28, tinf, tlb, N2_MASS - xmm, alpha[2] - 1.0, PTM[5], s);
  2692.             if (sw[15] != 0 && alt.getReal() <= altl[2]) {
  2693.                 /*  Mixed density at Alt */
  2694.                 dm28 = densu(alt, b28, tinf, tlb, xmm, alpha[2], PTM[5], s);
  2695.                 /*  Net density at Alt */
  2696.                 setDensity(MOLECULAR_NITROGEN, dnet(diffusiveDensity, dm28, zhm28, xmm, N2_MASS));
  2697.             } else {
  2698.                 dm28 = zero;
  2699.             }

  2700.             /**** He density ****/
  2701.             /*   Density variation factor at Zlb */
  2702.             final T g4 = globe7(PD[0]).multiply(sw[21]);
  2703.             /*  Diffusive density at Zlb */
  2704.             final T db04 = g4.exp().multiply(PDM[0][0] * PD[0][0]);
  2705.             /*  Diffusive density at Alt */
  2706.             diffusiveDensity = densu(alt, db04, tinf, tlb, HE_MASS, alpha[0], PTM[5], s);
  2707.             setDensity(HELIUM, diffusiveDensity);
  2708.             if (sw[15] != 0 && alt.getReal() < altl[0]) {
  2709.                 /*  Turbopause */
  2710.                 final double zh04 = PDM[0][2];
  2711.                 /*  Mixed density at Zlb */
  2712.                 final T b04 = densu(zero.newInstance(zh04), db04, tinf, tlb, HE_MASS - xmm, alpha[0] - 1., PTM[5], s);
  2713.                 /*  Mixed density at Alt */
  2714.                 final T dm04 = densu(alt, b04, tinf, tlb, xmm, 0., PTM[5], s);
  2715.                 final double zhm04 = zhm28;
  2716.                 /*  Net density at Alt */
  2717.                 diffusiveDensity = dnet(diffusiveDensity, dm04, zhm04, xmm, HE_MASS);
  2718.                 /*  Correction to specified mixing ratio at ground */
  2719.                 final T rl = b28.multiply(PDM[0][1]).divide(b04).log();
  2720.                 final double zc04 = PDM[0][4] * PDL[1][0];
  2721.                 final double hc04 = PDM[0][5] * PDL[1][1];
  2722.                 /*  Net density corrected at Alt */
  2723.                 setDensity(HELIUM, diffusiveDensity.multiply(ccor(alt, rl, hc04, zc04)));
  2724.             }

  2725.             /**** O density ****/
  2726.             /* Density variation factor at Zlb */
  2727.             final T g16 = globe7(PD[1]).multiply(sw[21]);
  2728.             /* Diffusive density at Zlb */
  2729.             final T db16 = g16.exp().multiply(PDM[1][0] * PD[1][0]);
  2730.             /* Diffusive density at Alt */
  2731.             diffusiveDensity = densu(alt, db16, tinf, tlb, O_MASS, alpha[1], PTM[5], s);
  2732.             setDensity(ATOMIC_OXYGEN, diffusiveDensity);
  2733.             if (sw[15] != 0 && alt.getReal() < altl[1]) {
  2734.                 /* Turbopause */
  2735.                 final double zh16 = PDM[1][2];
  2736.                 /* Mixed density at Zlb */
  2737.                 final T b16 = densu(zero.newInstance(zh16), db16, tinf, tlb, O_MASS - xmm, alpha[1] - 1.0, PTM[5], s);
  2738.                 /* Mixed density at Alt */
  2739.                 final T dm16 = densu(alt, b16, tinf, tlb, xmm, 0., PTM[5], s);
  2740.                 final double zhm16 = zhm28;
  2741.                 /* Net density at Alt */
  2742.                 diffusiveDensity = dnet(diffusiveDensity, dm16, zhm16, xmm, O_MASS);
  2743.                 final double rl = PDM[1][1] * PDL[1][16] * (1.0 + sw[1] * PDL[0][23] * (f107a - FLUX_REF));
  2744.                 final double hc16 = PDM[1][5] * PDL[1][3];
  2745.                 final double zc16 = PDM[1][4] * PDL[1][2];
  2746.                 final double hc216 = PDM[1][5] * PDL[1][4];
  2747.                 diffusiveDensity = diffusiveDensity.multiply(ccor2(alt, rl, hc16, zc16, hc216));
  2748.                 /* Chemistry correction */
  2749.                 final double hcc16 = PDM[1][7] * PDL[1][13];
  2750.                 final double zcc16 = PDM[1][6] * PDL[1][12];
  2751.                 final double rc16  = PDM[1][3] * PDL[1][14];
  2752.                 /* Net density corrected at Alt */
  2753.                 setDensity(ATOMIC_OXYGEN, diffusiveDensity.multiply(ccor(alt, zero.newInstance(rc16), hcc16, zcc16)));
  2754.             }

  2755.             /**** O2 density ****/
  2756.             /* Density variation factor at Zlb */
  2757.             final T g32 = globe7(PD[4]).multiply(sw[21]);
  2758.             /* Diffusive density at Zlb */
  2759.             final T db32 = g32.exp().multiply(PDM[3][0] * PD[4][0]);
  2760.             /* Diffusive density at Alt */
  2761.             diffusiveDensity = densu(alt, db32, tinf, tlb, O2_MASS, alpha[3], PTM[5], s);
  2762.             setDensity(MOLECULAR_OXYGEN, diffusiveDensity);
  2763.             if (sw[15] != 0) {
  2764.                 if (alt.getReal() <= altl[3]) {
  2765.                     /* Turbopause */
  2766.                     final double zh32 = PDM[3][2];
  2767.                     /* Mixed density at Zlb */
  2768.                     final T b32 = densu(zero.newInstance(zh32), db32, tinf, tlb, O2_MASS - xmm, alpha[3] - 1., PTM[5], s);
  2769.                     /* Mixed density at Alt */
  2770.                     final T dm32 = densu(alt, b32, tinf, tlb, xmm, 0., PTM[5], s);
  2771.                     final double zhm32 = zhm28;
  2772.                     /* Net density at Alt */
  2773.                     diffusiveDensity = dnet(diffusiveDensity, dm32, zhm32, xmm, O2_MASS);
  2774.                     /* Correction to specified mixing ratio at ground */
  2775.                     final T rl = b28.multiply(PDM[3][1]).divide(b32).log();
  2776.                     final double hc32 = PDM[3][5] * PDL[1][7];
  2777.                     final double zc32 = PDM[3][4] * PDL[1][6];
  2778.                     diffusiveDensity = diffusiveDensity.multiply(ccor(alt, rl, hc32, zc32));
  2779.                 }
  2780.                 /* Correction for general departure from diffusive equilibrium above Zlb */
  2781.                 final double hcc32  = PDM[3][7] * PDL[1][22];
  2782.                 final double hcc232 = PDM[3][7] * PDL[0][22];
  2783.                 final double zcc32  = PDM[3][6] * PDL[1][21];
  2784.                 final double rc32   = PDM[3][3] * PDL[1][23] * (1. + sw[1] * PDL[0][23] * (f107a - FLUX_REF));
  2785.                 /* Net density corrected at Alt */
  2786.                 setDensity(MOLECULAR_OXYGEN, diffusiveDensity.multiply(ccor2(alt, rc32, hcc32, zcc32, hcc232)));
  2787.             }

  2788.             /**** Ar density ****/
  2789.             /* Density variation factor at Zlb */
  2790.             final T g40 = globe7(PD[5]).multiply(sw[21]);
  2791.             /* Diffusive density at Zlb */
  2792.             final T db40 = g40.exp().multiply(PDM[4][0] * PD[5][0]);
  2793.             /* Diffusive density at Alt */
  2794.             diffusiveDensity = densu(alt, db40, tinf, tlb, AR_MASS, alpha[4], PTM[5], s);
  2795.             setDensity(ARGON, diffusiveDensity);
  2796.             if (sw[15] != 0 && alt.getReal() <= altl[4]) {
  2797.                 /* Turbopause */
  2798.                 final double zh40 = PDM[4][2];
  2799.                 /* Mixed density at Zlb */
  2800.                 final T b40 = densu(zero.newInstance(zh40), db40, tinf, tlb, AR_MASS - xmm, alpha[4] - 1., PTM[5], s);
  2801.                 /* Mixed density at Alt */
  2802.                 final T dm40 = densu(alt, b40, tinf, tlb, xmm, 0., PTM[5], s);
  2803.                 final double zhm40 = zhm28;
  2804.                 /* Net density at Alt */
  2805.                 diffusiveDensity = dnet(diffusiveDensity, dm40, zhm40, xmm, AR_MASS);
  2806.                 /* Correction to specified mixing ratio at ground */
  2807.                 final T rl = b28.multiply(PDM[4][1]).divide(b40).log();
  2808.                 final double hc40 = PDM[4][5] * PDL[1][9];
  2809.                 final double zc40 = PDM[4][4] * PDL[1][8];
  2810.                 /* Net density corrected at Alt */
  2811.                 setDensity(ARGON, diffusiveDensity.multiply(ccor(alt, rl, hc40, zc40)));
  2812.             }

  2813.             /**** H density ****/
  2814.             /* Density variation factor at Zlb */
  2815.             final T g1 = globe7(PD[6]).multiply(sw[21]);
  2816.             /* Diffusive density at Zlb */
  2817.             final T db01 = g1.exp().multiply(PDM[5][0] * PD[6][0]);
  2818.             /* Diffusive density at Alt */
  2819.             diffusiveDensity = densu(alt, db01, tinf, tlb, H_MASS, alpha[6], PTM[5], s);
  2820.             setDensity(HYDROGEN, diffusiveDensity);
  2821.             if (sw[15] != 0 && alt.getReal() <= altl[6]) {
  2822.                 /* Turbopause */
  2823.                 final double zh01 = PDM[5][2];
  2824.                 /* Mixed density at Zlb */
  2825.                 final T b01 = densu(zero.newInstance(zh01), db01, tinf, tlb, H_MASS - xmm, alpha[6] - 1., PTM[5], s);
  2826.                 /* Mixed density at Alt */
  2827.                 final T dm01 = densu(alt, b01, tinf, tlb, xmm, 0., PTM[5], s);
  2828.                 final double zhm01 = zhm28;
  2829.                 /* Net density at Alt */
  2830.                 diffusiveDensity = dnet(diffusiveDensity, dm01, zhm01, xmm, H_MASS);
  2831.                 /* Correction to specified mixing ratio at ground */
  2832.                 final T rl = b28.multiply(PDM[5][1] * FastMath.sqrt(PDL[1][17] * PDL[1][17])).divide(b01).log();
  2833.                 final double hc01 = PDM[5][5] * PDL[1][11];
  2834.                 final double zc01 = PDM[5][4] * PDL[1][10];
  2835.                 diffusiveDensity = diffusiveDensity.multiply(ccor(alt, rl, hc01, zc01));
  2836.                 /* Chemistry correction */
  2837.                 final double hcc01 = PDM[5][7] * PDL[1][19];
  2838.                 final double zcc01 = PDM[5][6] * PDL[1][18];
  2839.                 final double rc01 = PDM[5][3] * PDL[1][20];
  2840.                 /* Net density corrected at Alt */
  2841.                 setDensity(HYDROGEN, diffusiveDensity.multiply(ccor(alt, zero.newInstance(rc01), hcc01, zcc01)));
  2842.             }

  2843.             /**** N density ****/
  2844.             /* Density variation factor at Zlb */
  2845.             final T g14 = globe7(PD[7]).multiply(sw[21]);
  2846.             /* Diffusive density at Zlb */
  2847.             final T db14 = g14.exp().multiply(PDM[6][0] * PD[7][0]);
  2848.             /* Diffusive density at Alt */
  2849.             diffusiveDensity = densu(alt, db14, tinf, tlb, N_MASS, alpha[7], PTM[5], s);
  2850.             setDensity(ATOMIC_NITROGEN, diffusiveDensity);
  2851.             if (sw[15] != 0 && alt.getReal() <= altl[7]) {
  2852.                 /* Turbopause */
  2853.                 final double zh14 = PDM[6][2];
  2854.                 /* Mixed density at Zlb */
  2855.                 final T b14 = densu(zero.newInstance(zh14), db14, tinf, tlb, N_MASS - xmm, alpha[7] - 1., PTM[5], s);
  2856.                 /* Mixed density at Alt */
  2857.                 final T dm14 = densu(alt, b14, tinf, tlb, xmm, 0., PTM[5], s);
  2858.                 final double zhm14 = zhm28;
  2859.                 /* Net density at Alt */
  2860.                 diffusiveDensity = dnet(diffusiveDensity, dm14, zhm14, xmm, N_MASS);
  2861.                 /* Correction to specified mixing ratio at ground */
  2862.                 final T rl = b28.multiply(PDM[6][1] * PDL[0][2]).divide(b14).log();
  2863.                 final double hc14 = PDM[6][5] * PDL[0][1];
  2864.                 final double zc14 = PDM[6][4] * PDL[0][0];
  2865.                 diffusiveDensity = diffusiveDensity.multiply(ccor(alt, rl, hc14, zc14));
  2866.                 /* Chemistry correction */
  2867.                 final double hcc14 = PDM[6][7] * PDL[0][4];
  2868.                 final double zcc14 = PDM[6][6] * PDL[0][3];
  2869.                 final double rc14 = PDM[6][3] * PDL[0][5];
  2870.                 /* Net density corrected at Alt */
  2871.                 setDensity(ATOMIC_NITROGEN, diffusiveDensity.multiply(ccor(alt, zero.newInstance(rc14), hcc14, zcc14)));
  2872.             }

  2873.             /**** Anomalous O density ****/
  2874.             final T g16h = globe7(PD[8]).multiply(sw[21]);
  2875.             final T db16h = g16h.exp().multiply(PDM[7][0] * PD[8][0]);
  2876.             final double tho   = PDM[7][9] * PDL[0][6];
  2877.             diffusiveDensity = densu(alt, db16h, zero.newInstance(tho), zero.newInstance(tho), O_MASS, alpha[8], PTM[5], s);
  2878.             final double zsht = PDM[7][5];
  2879.             final double zmho = PDM[7][4];
  2880.             final T zsho = scalh(zmho, O_MASS, tho);
  2881.             diffusiveDensity = diffusiveDensity.multiply(alt.negate().add(zmho).divide(zsht).exp().subtract(1).multiply(-zsht).divide(zsho).exp());
  2882.             setDensity(ANOMALOUS_OXYGEN, diffusiveDensity);

  2883.             // Convert densities from cm-3 to m-3
  2884.             for (int i = 0; i < 9; i++) {
  2885.                 setDensity(i, getDensity(i).multiply(1.0e+06));
  2886.             }

  2887.             /**** Total mass density ****/
  2888.             final T tmd =     getDensity(HELIUM)            .multiply(HE_MASS).
  2889.                           add(getDensity(ATOMIC_OXYGEN)     .multiply( O_MASS)).
  2890.                           add(getDensity(MOLECULAR_NITROGEN).multiply(N2_MASS)).
  2891.                           add(getDensity(MOLECULAR_OXYGEN)  .multiply(O2_MASS)).
  2892.                           add(getDensity(ARGON)             .multiply(AR_MASS)).
  2893.                           add(getDensity(HYDROGEN)          .multiply( H_MASS)).
  2894.                           add(getDensity(ATOMIC_NITROGEN)   .multiply( N_MASS)).
  2895.                           multiply(AMU);
  2896.             setDensity(TOTAL_MASS, tmd);

  2897.         }

  2898.         /** Calculate temperatures and densities not including anomalous oxygen.
  2899.          *  <p>NOTES ON INPUT VARIABLES:<br>
  2900.          *  Seconds, Local Time, and Longitude are used independently in the
  2901.          *  model and are not of equal importance for every situation.<br>
  2902.          *  For the most physically realistic calculation these three
  2903.          *  variables should be consistent (lst=sec/3600 + lon/15).<br>
  2904.          *  The Equation of Time departures from the above formula
  2905.          *  for apparent local time can be included if available but
  2906.          *  are of minor importance.<br><br>
  2907.          *
  2908.          *  f107 and f107A values used to generate the model correspond
  2909.          *  to the 10.7 cm radio flux at the actual distance of the Earth
  2910.          *  from the Sun rather than the radio flux at 1 AU. The following
  2911.          *  site provides both classes of values:<br>
  2912.          *  ftp://ftp.ngdc.noaa.gov/STP/SOLAR_DATA/SOLAR_RADIO/FLUX/<br><br>
  2913.          *
  2914.          *  f107, f107A, and ap effects are neither large nor well established below 80 km
  2915.          *  and these parameters should be set to 150., 150., and 4. respectively.
  2916.          *  </p>
  2917.          *  @param alt altitude (km)
  2918.          */
  2919.         void gtd7(final T alt) {

  2920.             // Calculates for thermosphere/mesosphere (above ZN2[0])
  2921.             final T altt = (alt.getReal() > ZN2[0]) ? alt : zero.newInstance(ZN2[0]);
  2922.             gts7(altt);
  2923.             if (alt.getReal() >= ZN2[0]) {
  2924.                 return;
  2925.             }

  2926.             // Calculates for lower mesosphere/upper stratosphere (between ZN2[0] and ZN3[0]):
  2927.             // Temperature at nodes and gradients at end nodes
  2928.             // Inverse temperature a linear function of spherical harmonics
  2929.             final double r = PMA[2][0] * PAVGM[2];
  2930.             meso_tgn2[0] = meso_tgn1[1];
  2931.             meso_tn2[0]  = meso_tn1[4];
  2932.             meso_tn2[1]  = glob7s(PMA[0]).multiply(sw[20]         ).negate().add(1).reciprocal().multiply(PMA[0][0] * PAVGM[0]);
  2933.             meso_tn2[2]  = glob7s(PMA[1]).multiply(sw[20]         ).negate().add(1).reciprocal().multiply(PMA[1][0] * PAVGM[1]);
  2934.             meso_tn2[3]  = glob7s(PMA[2]).multiply(sw[20] * sw[22]).negate().add(1).reciprocal().multiply(PMA[2][0] * PAVGM[2]);
  2935.             meso_tgn2[1] = glob7s(PMA[9]).multiply(sw[20] * sw[22]).add(1).multiply(PMA[9][0] * PAVGM[8]).
  2936.                            multiply(meso_tn2[3]).multiply(meso_tn2[3]).divide(r * r);
  2937.             meso_tn3[0]  = meso_tn2[3];

  2938.             // Calculates for lower stratosphere and troposphere (below ZN3[0])
  2939.             // Temperature at nodes and gradients at end nodes
  2940.             // Inverse temperature a linear function of spherical harmonics
  2941.             if (alt.getReal() <= ZN3[0]) {
  2942.                 final double q = PMA[6][0] * PAVGM[6];
  2943.                 meso_tgn3[0] = meso_tgn2[1];
  2944.                 meso_tn3[1]  = glob7s(PMA[3]).multiply(sw[22]).negate().add(1).reciprocal().multiply(PMA[3][0] * PAVGM[3]);
  2945.                 meso_tn3[2]  = glob7s(PMA[4]).multiply(sw[22]).negate().add(1).reciprocal().multiply(PMA[4][0] * PAVGM[4]);
  2946.                 meso_tn3[3]  = glob7s(PMA[5]).multiply(sw[22]).negate().add(1).reciprocal().multiply(PMA[5][0] * PAVGM[5]);
  2947.                 meso_tn3[4]  = glob7s(PMA[6]).multiply(sw[22]).negate().add(1).reciprocal().multiply(PMA[6][0] * PAVGM[6]);
  2948.                 meso_tgn3[1] = glob7s(PMA[7]).multiply(sw[22])         .add(1).multiply(PMA[7][0] * PAVGM[7]).
  2949.                                multiply(meso_tn3[4]).multiply(meso_tn3[4]).divide(q * q);

  2950.             }

  2951.             // Linear transition to full mixing below ZN2[0]
  2952.             final T dmc = (alt.getReal() > ZMIX) ?
  2953.                            alt.subtract(ZN2[0]).divide(ZN2[0] - ZMIX).add(1) :
  2954.                            zero;
  2955.             final T dz28 = getDensity(MOLECULAR_NITROGEN);

  2956.             // N2 density
  2957.             final T dm28m = dm28.multiply(1.0e+06);
  2958.             T dmr = dz28.divide(dm28m).subtract(1);
  2959.             T dst = densm(alt, dm28m, PDM[2][4]).multiply(dmr.multiply(dmc).add(1));
  2960.             setDensity(MOLECULAR_NITROGEN, dst);

  2961.             // HE density
  2962.             dmr = getDensity(HELIUM).divide(dz28.multiply(PDM[0][1])).subtract(1);
  2963.             dst = getDensity(MOLECULAR_NITROGEN).multiply(PDM[0][1]).multiply(dmr.multiply(dmc).add(1));
  2964.             setDensity(HELIUM, dst);

  2965.             // O density
  2966.             setDensity(ATOMIC_OXYGEN, zero);
  2967.             setDensity(ANOMALOUS_OXYGEN, zero);

  2968.             // O2 density
  2969.             dmr = getDensity(MOLECULAR_OXYGEN).divide(dz28.multiply(PDM[3][1])).subtract(1);
  2970.             dst = getDensity(MOLECULAR_NITROGEN).multiply(PDM[3][1]).multiply(dmr.multiply(dmc).add(1));
  2971.             setDensity(MOLECULAR_OXYGEN, dst);

  2972.             // AR density
  2973.             dmr = getDensity(ARGON).divide(dz28.multiply(PDM[4][1])).subtract(1);
  2974.             dst = getDensity(MOLECULAR_NITROGEN).multiply(PDM[4][1]).multiply(dmr.multiply(dmc).add(1));
  2975.             setDensity(ARGON, dst);

  2976.             // H density
  2977.             setDensity(HYDROGEN, zero);

  2978.             // N density
  2979.             setDensity(ATOMIC_NITROGEN, zero);

  2980.             // Total mass density
  2981.             final T tmd =       getDensity(HELIUM)            .multiply(HE_MASS).
  2982.                             add(getDensity(ATOMIC_OXYGEN)     .multiply( O_MASS)).
  2983.                             add(getDensity(MOLECULAR_NITROGEN).multiply(N2_MASS)).
  2984.                             add(getDensity(MOLECULAR_OXYGEN)  .multiply(O2_MASS)).
  2985.                             add(getDensity(ARGON)             .multiply(AR_MASS)).
  2986.                             add(getDensity(HYDROGEN)          .multiply( H_MASS)).
  2987.                             add(getDensity(ATOMIC_NITROGEN)   .multiply( N_MASS)).
  2988.                             multiply(AMU);
  2989.             setDensity(TOTAL_MASS, tmd);

  2990.             // Temperature at altitude
  2991.             setTemperature(ALTITUDE, densm(alt, field.getOne(), 0));

  2992.         }

  2993.         /** Calculate temperatures and densities including anomalous oxygen.
  2994.          *  <p></p>
  2995.          *  <p>NOTES ON INPUT VARIABLES:<br>
  2996.          *  Seconds, Local Time, and Longitude are used independently in the
  2997.          *  model and are not of equal importance for every situation.<br>
  2998.          *  For the most physically realistic calculation these three
  2999.          *  variables should be consistent (lst=sec/3600 + lon/15).<br>
  3000.          *  The Equation of Time departures from the above formula
  3001.          *  for apparent local time can be included if available but
  3002.          *  are of minor importance.<br>
  3003.          *  <br>
  3004.          *  f107 and f107A values used to generate the model correspond
  3005.          *  to the 10.7 cm radio flux at the actual distance of the Earth
  3006.          *  from the Sun rather than the radio flux at 1 AU. The following
  3007.          *  site provides both classes of values:<br>
  3008.          *  ftp://ftp.ngdc.noaa.gov/STP/SOLAR_DATA/SOLAR_RADIO/FLUX/<br>
  3009.          *  <br>
  3010.          *  f107, f107A, and ap effects are neither large nor well established below 80 km
  3011.          *  and these parameters should be set to 150., 150., and 4. respectively.
  3012.          *  </p>
  3013.          *  @param alt altitude (km)
  3014.          */
  3015.         void gtd7d(final T alt) {

  3016.             // Compute densities and temperatures
  3017.             gtd7(alt);

  3018.             // Update the total mass density with anomalous oxygen contribution
  3019.             final T dTot = getDensity(TOTAL_MASS).add(getDensity(ANOMALOUS_OXYGEN).multiply( AMU * O_MASS));
  3020.             setDensity(TOTAL_MASS, dTot);

  3021.         }

  3022.         /** Set one density.
  3023.          * @param index one of the nine elements :
  3024.          * <ul>
  3025.          * <li>{@link #HELIUM}</li>
  3026.          * <li>{@link #ATOMIC_OXYGEN}</li>
  3027.          * <li>{@link #MOLECULAR_NITROGEN}</li>
  3028.          * <li>{@link #MOLECULAR_OXYGEN}</li>
  3029.          * <li>{@link #ARGON}</li>
  3030.          * <li>{@link #TOTAL_MASS}</li>
  3031.          * <li>{@link #HYDROGEN}</li>
  3032.          * <li>{@link #ATOMIC_NITROGEN}</li>
  3033.          * <li>{@link #ATOMIC_NITROGEN}</li>
  3034.          * </ul>
  3035.          * @param d the value of density to set
  3036.          */
  3037.         void setDensity(final int index, final T d) {
  3038.             densities[index] = d;
  3039.         }

  3040.         /** Set one temperature.
  3041.          * @param index one of the two elements :
  3042.          * <ul>
  3043.          * <li>{@link #EXOSPHERIC}</li>
  3044.          * <li>{@link #ALTITUDE}</li>
  3045.          * </ul>
  3046.          * @param t the value of temperature to set
  3047.          */
  3048.         void setTemperature(final int index, final T t) {
  3049.             temperatures[index] = t;
  3050.         }

  3051.         /** Get one of the stored densities.
  3052.          * @param index one of the nine elements :
  3053.          * <ul>
  3054.          * <li>{@link #HELIUM}</li>
  3055.          * <li>{@link #ATOMIC_OXYGEN}</li>
  3056.          * <li>{@link #MOLECULAR_NITROGEN}</li>
  3057.          * <li>{@link #MOLECULAR_OXYGEN}</li>
  3058.          * <li>{@link #ARGON}</li>
  3059.          * <li>{@link #TOTAL_MASS}</li>
  3060.          * <li>{@link #HYDROGEN}</li>
  3061.          * <li>{@link #ATOMIC_NITROGEN}</li>
  3062.          * <li>{@link #ATOMIC_NITROGEN}</li>
  3063.          * </ul>
  3064.          * @return the requested density
  3065.          */
  3066.         public T getDensity(final int index) {
  3067.             return densities[index];
  3068.         }

  3069.         /** Calculate G(L) function with upper thermosphere parameters.
  3070.          *  @param p array of parameters
  3071.          *  @return G(L) value
  3072.          */
  3073.         private T globe7(final double[] p) {

  3074.             final T[] t = MathArrays.buildArray(field, 14);
  3075.             final double cd32 = FastMath.cos(DAY_TO_RAD * (doy - p[31]));
  3076.             final double cd18 = FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[17]));
  3077.             final double cd14 = FastMath.cos(DAY_TO_RAD * (doy - p[13]));
  3078.             final double cd39 = FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[38]));

  3079.             // F10.7 effect
  3080.             final double df  = f107  - f107a;
  3081.             final double dfa = f107a - FLUX_REF;
  3082.             t[0] = zero.newInstance(p[19] * df * (1.0 + p[59] * dfa) +
  3083.                                     p[20] * df * df +
  3084.                                     p[21] * dfa +
  3085.                                     p[29] * dfa * dfa);

  3086.             final double f1 = 1.0 + (p[47] * dfa + p[19] * df + p[20] * df * df) * swc[1];
  3087.             final double f2 = 1.0 + (p[49] * dfa + p[19] * df + p[20] * df * df) * swc[1];

  3088.             // Time independent
  3089.             t[1] =     plg[0][2].multiply(p[ 1]).
  3090.                    add(plg[0][4].multiply(p[ 2])).
  3091.                    add(plg[0][6].multiply(p[22])).
  3092.                    add(plg[0][2].multiply(p[14] * dfa * swc[1])).
  3093.                    add(plg[0][1].multiply(p[26]));

  3094.             // Symmetrical annual
  3095.             t[2] = zero.newInstance(p[18] * cd32);

  3096.             // Symmetrical semiannual
  3097.             t[3] = plg[0][2].multiply(p[16]).add(p[15]).multiply(cd18);

  3098.             // Asymmetrical annual
  3099.             t[4] = plg[0][1].multiply(p[9]).add(plg[0][3].multiply(p[10])).multiply(f1 * cd14);

  3100.             // Asymmetrical semiannual
  3101.             t[5] = plg[0][1].multiply(p[37] * cd39);

  3102.             // Diurnal
  3103.             if (sw[7] != 0) {
  3104.                 final T t71 = plg[1][2].multiply(p[11] * cd14 * swc[5]);
  3105.                 final T t72 = plg[1][2].multiply(p[12] * cd14 * swc[5]);
  3106.                 t[6] =      plg[1][1].multiply(p[3]).add(plg[1][3].multiply(p[4])).add(plg[1][5].multiply(p[27])).add(t71).multiply(ctloc).
  3107.                         add(plg[1][1].multiply(p[6]).add(plg[1][3].multiply(p[7])).add(plg[1][5].multiply(p[28])).add(t72).multiply(stloc)).
  3108.                         multiply(f2);
  3109.             }

  3110.             // Semidiurnal
  3111.             if (sw[8] != 0) {
  3112.                 final T t81 = plg[2][3].multiply(p[23]).add(plg[2][5].multiply(p[35])).multiply(cd14 * swc[5]);
  3113.                 final T t82 = plg[2][3].multiply(p[33]).add(plg[2][5].multiply(p[36])).multiply(cd14 * swc[5]);
  3114.                 t[7] =     plg[2][2].multiply(p[5]).add(plg[2][4].multiply(p[41])).add(t81).multiply(c2tloc).
  3115.                        add(plg[2][2].multiply(p[8]).add(plg[2][4].multiply(p[42])).add(t82).multiply(s2tloc)).
  3116.                        multiply(f2);
  3117.             }

  3118.             // Terdiurnal
  3119.             if (sw[14] != 0) {
  3120.                 t[13] =     plg[3][3].multiply(p[39]).add(plg[3][4].multiply(p[93]).add(plg[3][6].multiply(p[46])).multiply(cd14 * swc[5])).multiply(s3tloc).
  3121.                         add(plg[3][3].multiply(p[40]).add(plg[3][4].multiply(p[94]).add(plg[3][6].multiply(p[48])).multiply(cd14 * swc[5])).multiply(c3tloc)).
  3122.                         multiply(f2);
  3123.             }

  3124.             // magnetic activity based on daily ap
  3125.             if (sw[9] == -1) {
  3126.                 if (p[51] != 0) {
  3127.                     final T exp1 = lat.abs().negate().add(LAT_REF).multiply(p[138]).add(1).
  3128.                                     reciprocal().multiply(-10800.0 * FastMath.abs(p[51])).
  3129.                                     exp();
  3130.                     final double p24 = FastMath.max(p[24], 1.0e-4);
  3131.                     apt = sg0(min(0.99999, exp1), p24, p[25]);
  3132.                     t[8] =      plg[0][2].multiply(p[96]).add(plg[0][4].multiply(p[54])).add(p[50]).
  3133.                            add((plg[0][1].multiply(p[125]).add(plg[0][3].multiply(p[126])).add(plg[0][5].multiply(p[127]))).multiply(cd14 * swc[5])).
  3134.                            add((plg[1][1].multiply(p[128]).add(plg[1][3].multiply(p[129])).add(plg[1][5].multiply(p[130]))).multiply(swc[7]).multiply(hl.subtract(p[131]).multiply(HOUR_TO_RAD).cos())).
  3135.                            multiply(apt);
  3136.                 }
  3137.             } else {
  3138.                 final double apd = ap[0] - 4.0;
  3139.                 final double p44 = (p[43] < 0.) ? 1.0E-5 : p[43];
  3140.                 final double p45 = p[44];
  3141.                 apdf = apd + (p45 - 1.0) * (apd + (FastMath.exp(-p44 * apd) - 1.0) / p44);
  3142.                 if (sw[9] != 0) {
  3143.                     t[8] =      plg[0][2].multiply(p[45]).add(plg[0][4].multiply(p[34])).add(p[32]).
  3144.                            add((plg[0][1].multiply(p[100]).add(plg[0][3].multiply(p[101])).add(plg[0][5].multiply(p[102]))).multiply(cd14 * swc[5])).
  3145.                            add((plg[1][1].multiply(p[121]).add(plg[1][3].multiply(p[122])).add(plg[1][5].multiply(p[123]))).multiply(swc[7]).multiply(hl.subtract(p[124]).multiply(HOUR_TO_RAD).cos())).
  3146.                            multiply(apdf);
  3147.                 }
  3148.             }

  3149.             if (sw[10] != 0) {
  3150.                 final T lonr = lon.multiply(DEG_TO_RAD);
  3151.                 final FieldSinCos<T> scLonr = FastMath.sinCos(lonr);
  3152.                 // Longitudinal
  3153.                 if (sw[11] != 0) {
  3154.                     t[10] =         plg[1][2].multiply(p[ 64]) .add(plg[1][4].multiply(p[ 65])).add(plg[1][6].multiply(p[ 66])).
  3155.                                 add(plg[1][1].multiply(p[103])).add(plg[1][3].multiply(p[104])).add(plg[1][5].multiply(p[105])).
  3156.                                 add((plg[1][1].multiply(p[109])).add(plg[1][3].multiply(p[110])).add(plg[1][5].multiply(p[111])).multiply(swc[5] * cd14)).
  3157.                                 multiply(scLonr.cos()).
  3158.                             add(    plg[1][2].multiply(p[ 90]) .add(plg[1][4].multiply(p[ 91])).add(plg[1][6].multiply(p[ 92])).
  3159.                                 add(plg[1][1].multiply(p[106])).add(plg[1][3].multiply(p[107])).add(plg[1][5].multiply(p[108])).
  3160.                                 add((plg[1][1].multiply(p[112])).add(plg[1][3].multiply(p[113])).add(plg[1][5].multiply(p[114])).multiply(swc[5] * cd14)).
  3161.                                 multiply(scLonr.sin())).
  3162.                             multiply(1.0 + p[80] * dfa * swc[1]);
  3163.                 }

  3164.                 // ut and mixed ut, longitude
  3165.                 if (sw[12] != 0) {
  3166.                     t[11] =          plg[0][1].multiply(p[95]).add(1).multiply(1.0 + p[81] * dfa * swc[1]).
  3167.                             multiply(plg[0][1].multiply(p[119] * swc[5] * cd14).add(1)).
  3168.                             multiply(plg[0][1].multiply(p[68]).add(plg[0][3].multiply(p[69])).add(plg[0][5].multiply(p[70]))).
  3169.                             multiply(sec.subtract(p[71]).multiply(SEC_TO_RAD).cos());
  3170.                     t[11] = t[11].
  3171.                             add(plg[2][3].multiply(p[76]).add(plg[2][5].multiply(p[77])).add(plg[2][7].multiply(p[78])).
  3172.                                 multiply(swc[11] * (1.0 + p[137] * dfa * swc[1])).
  3173.                                 multiply(sec.subtract(p[79]).multiply(SEC_TO_RAD).add(lonr.multiply(2)).cos()));
  3174.                 }

  3175.                 /* ut, longitude magnetic activity */
  3176.                 if (sw[13] != 0) {
  3177.                     if (sw[9] == -1) {
  3178.                         if (p[51] != 0.) {
  3179.                             t[12] = apt.multiply(swc[11]).multiply(plg[0][1].multiply(p[132]).add(1)).
  3180.                                     multiply(plg[1][2].multiply(p[52]).add(plg[1][4].multiply(p[98])).add(plg[1][6].multiply(p[67]))).
  3181.                                     multiply(lon.subtract(p[97]).multiply(DEG_TO_RAD).cos()).
  3182.                                     add(apt.multiply(swc[11] * swc[5] * cd14).
  3183.                                         multiply(plg[1][1].multiply(p[133]).add(plg[1][3].multiply(p[134])).add(plg[1][5].multiply(p[135]))).
  3184.                                         multiply(lon.subtract(p[136]).multiply(DEG_TO_RAD).cos())).
  3185.                                     add(apt.multiply(swc[12]).
  3186.                                         multiply(plg[0][1].multiply(p[55]).add(plg[0][3].multiply(p[56])).add(plg[0][5].multiply(p[57]))).
  3187.                                         multiply(sec.subtract(p[58]).multiply(SEC_TO_RAD).cos()));
  3188.                         }
  3189.                     } else {
  3190.                         t[12] = plg[0][1].multiply(p[120]).add(1).multiply(apdf * swc[11]).
  3191.                                 multiply(plg[1][2].multiply(p[60]).add(plg[1][4].multiply(p[61])).add(plg[1][6].multiply(p[62]))).
  3192.                                 multiply(lon.subtract(p[63]).multiply(DEG_TO_RAD).cos()).
  3193.                                 add(plg[1][1].multiply(p[115]).add(plg[1][3].multiply(p[116])).add(plg[1][5].multiply(p[117])).
  3194.                                     multiply(apdf * swc[11] * swc[5] * cd14).
  3195.                                     multiply(lon.subtract(p[118]).multiply(DEG_TO_RAD).cos())).
  3196.                                 add(plg[0][1].multiply(p[83]).add(plg[0][3].multiply(p[84])).add(plg[0][5].multiply(p[85])).
  3197.                                     multiply(apdf * swc[12]).
  3198.                                     multiply(sec.subtract(p[75]).multiply(SEC_TO_RAD).cos()));
  3199.                     }
  3200.                 }
  3201.             }

  3202.             // Sum all effects (params not used: 82, 89, 99, 139-149)
  3203.             T tinf = zero.newInstance(p[30]);
  3204.             for (int i = 0; i < 14; i++) {
  3205.                 tinf = tinf.add(t[i].multiply(FastMath.abs(sw[i + 1])));
  3206.             }

  3207.             // Return G(L)
  3208.             return tinf;

  3209.         }

  3210.         /** Calculate G(L) function with lower atmosphere parameters.
  3211.          *  @param p array of parameters
  3212.          *  @return G(L) value
  3213.          */
  3214.         private T glob7s(final double[] p) {

  3215.             final T[] t = MathArrays.buildArray(field, 14);
  3216.             final double cd32 = FastMath.cos(DAY_TO_RAD * (doy - p[31]));
  3217.             final double cd18 = FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[17]));
  3218.             final double cd14 = FastMath.cos(DAY_TO_RAD * (doy - p[13]));
  3219.             final double cd39 = FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[38]));

  3220.             // F10.7 effect
  3221.             t[0] = zero.newInstance(p[21] * (f107a - FLUX_REF));

  3222.             // Time independent
  3223.             t[1] =     plg[0][2].multiply(p[1]).
  3224.                    add(plg[0][4].multiply(p[2])).
  3225.                    add(plg[0][6].multiply(p[22])).
  3226.                    add(plg[0][1].multiply(p[26])).
  3227.                    add(plg[0][3].multiply(p[14])).
  3228.                    add(plg[0][5].multiply(p[59]));

  3229.             // Symmetrical annual
  3230.             t[2] = plg[0][2].multiply(p[47]).add(plg[0][4].multiply(p[29])).add(p[18]).multiply(cd32);

  3231.             // Symmetrical semiannual
  3232.             t[3] = plg[0][2].multiply(p[16]).add(plg[0][4].multiply(p[30])).add(p[15]).multiply(cd18);

  3233.             // Asymmetrical annual
  3234.             t[4] = plg[0][1].multiply(p[9]).add(plg[0][3].multiply(p[10])).add(plg[0][5].multiply(p[20])).multiply(cd14);

  3235.             // Asymmetrical semiannual
  3236.             t[5] = plg[0][1].multiply(p[37]).multiply(cd39);

  3237.             // Diurnal
  3238.             if (sw[7] != 0) {
  3239.                 final T t71 = plg[1][2].multiply(p[11]).multiply(cd14 * swc[5]);
  3240.                 final T t72 = plg[1][2].multiply(p[12]).multiply(cd14 * swc[5]);
  3241.                 t[6] =     plg[1][1].multiply(p[3]).add(plg[1][3].multiply(p[4])).add(t71).multiply(ctloc).
  3242.                        add(plg[1][1].multiply(p[6]).add(plg[1][3].multiply(p[7])).add(t72).multiply(stloc));
  3243.             }

  3244.             // Semidiurnal
  3245.             if (sw[8] != 0) {
  3246.                 final T t81 = plg[2][3].multiply(p[23]).add(plg[2][5].multiply(p[35])).multiply(cd14 * swc[5]);
  3247.                 final T t82 = plg[2][3].multiply(p[33]).add(plg[2][5].multiply(p[36])).multiply(cd14 * swc[5]);
  3248.                 t[7] =     plg[2][2].multiply(p[5]).add(plg[2][4].multiply(p[41])).add(t81).multiply(c2tloc).
  3249.                        add(plg[2][2].multiply(p[8]).add(plg[2][4].multiply(p[42])).add(t82).multiply(s2tloc));
  3250.             }

  3251.             // Terdiurnal
  3252.             if (sw[14] != 0) {
  3253.                 t[13] = plg[3][3].multiply(p[39]).multiply(s3tloc).add(plg[3][3].multiply(p[40]).multiply(c3tloc));
  3254.             }

  3255.             // Magnetic activity
  3256.             if (sw[9] == 1) {
  3257.                 t[8] = plg[0][2].multiply(p[45] * swc[2]).add(p[32]).multiply(apdf);
  3258.             } else if (sw[9] == -1) {
  3259.                 t[8] = plg[0][2].multiply(p[96] * swc[2]).add(p[50]).multiply(apt);
  3260.             }

  3261.             // Longitudinal
  3262.             if (!(sw[10] == 0 || sw[11] == 0)) {
  3263.                 final T lonr = lon.multiply(DEG_TO_RAD);
  3264.                 final FieldSinCos<T> scLonr = FastMath.sinCos(lonr);
  3265.                 t[10] = plg[0][1].multiply(p[80] * swc[5] * FastMath.cos(DAY_TO_RAD * (doy - p[81])) +
  3266.                                            p[85] * swc[6] * FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[86]))).
  3267.                        add(1.0 +
  3268.                            p[83] * swc[3] * FastMath.cos(DAY_TO_RAD * (doy - p[84])) +
  3269.                            p[87] * swc[4] * FastMath.cos(2.0 * DAY_TO_RAD * (doy - p[88]))).
  3270.                        multiply(    plg[1][2].multiply(p[64]).
  3271.                                 add(plg[1][4].multiply(p[65])).
  3272.                                 add(plg[1][6].multiply(p[66])).
  3273.                                 add(plg[1][1].multiply(p[74])).
  3274.                                 add(plg[1][3].multiply(p[75])).
  3275.                                 add(plg[1][5].multiply(p[76])).multiply(scLonr.cos()).
  3276.                           add(      plg[1][2].multiply(p[90]).
  3277.                                 add(plg[1][4].multiply(p[91])).
  3278.                                 add(plg[1][6].multiply(p[92])).
  3279.                                 add(plg[1][1].multiply(p[77])).
  3280.                                 add(plg[1][3].multiply(p[78])).
  3281.                                 add(plg[1][5].multiply(p[79])).multiply(scLonr.sin())));
  3282.             }

  3283.             // Sum all effects
  3284.             T gl = zero;
  3285.             for (int i = 0; i < 14; i++) {
  3286.                 gl = gl.add(t[i].multiply(FastMath.abs(sw[i + 1])));
  3287.             }

  3288.             // Return G(L)
  3289.             return gl;
  3290.         }

  3291.         /** Implements sg0 function (Eq. A24a).
  3292.          * @param ex ex
  3293.          * @param p24 abs(p[24])
  3294.          * @param p25 p[25]
  3295.          * @return sg0
  3296.          */
  3297.         private T sg0(final T ex, final double p24, final double p25) {
  3298.             final double g01 = g0(ap[1], p24, p25);
  3299.             final double g02 = g0(ap[2], p24, p25);
  3300.             final double g03 = g0(ap[3], p24, p25);
  3301.             final double g04 = g0(ap[4], p24, p25);
  3302.             final double g05 = g0(ap[5], p24, p25);
  3303.             final double g06 = g0(ap[6], p24, p25);
  3304.             final T ex2      = ex.square();
  3305.             final T ex3      = ex.multiply(ex2);
  3306.             final T ex4      = ex2.square();
  3307.             final T ex8      = ex4.square();
  3308.             final T ex12     = ex4.multiply(ex8);
  3309.             final T g234     = ex.multiply(g02).add(ex2.multiply(g03)).add(ex3.multiply(g04));
  3310.             final T g56      = ex4.multiply(g05).add(ex12.multiply(g06));
  3311.             final T ex19     = ex3.multiply(ex4).multiply(ex12);
  3312.             final T omex     = ex.negate().add(1);
  3313.             final T sumex    = ex19.negate().add(1).divide(omex).multiply(ex.sqrt()).add(1);
  3314.             return ex8.negate().add(1).multiply(g56).divide(omex).add(g234).add(g01).divide(sumex);
  3315.         }

  3316.         /** Implements go function (Eq. A24d).
  3317.          * @param apI 3 hrs ap
  3318.          * @param p24 abs(p[24])
  3319.          * @param p25 p[25]
  3320.          * @return go
  3321.          */
  3322.         private double g0(final double apI, final double p24, final double p25) {
  3323.             final double am4 = apI - 4.0;
  3324.             return am4 + (p25 - 1.0) * (am4 + (FastMath.exp(-p24 * am4) - 1.0) / p24);
  3325.         }

  3326.         /** Calculates chemistry/dissociation correction for MSIS models.
  3327.          * @param alt altitude
  3328.          * @param r target ratio
  3329.          * @param h1 transition scale length
  3330.          * @param zh altitude of 1/2 R
  3331.          * @return correction
  3332.          */
  3333.         private T ccor(final T alt, final T r, final double h1, final double zh) {
  3334.             final T e = alt.subtract(zh).divide(h1);
  3335.             if (e.getReal() > 70.) {
  3336.                 return field.getOne();
  3337.             } else if (e.getReal() < -70.) {
  3338.                 return r.exp();
  3339.             } else {
  3340.                 return r.divide(e.exp().add(1)).exp();
  3341.             }
  3342.         }


  3343.         /** Calculates O & O2 chemistry/dissociation correction for MSIS models.
  3344.          * @param alt altitude
  3345.          * @param r target ratio
  3346.          * @param h1 transition scale length
  3347.          * @param zh altitude of 1/2 R
  3348.          * @param h2 transition scale length
  3349.          * @return correction
  3350.          */
  3351.         private T ccor2(final T alt, final double r, final double h1, final double zh, final double h2) {
  3352.             final T e1 = alt.subtract(zh).divide(h1);
  3353.             final T e2 = alt.subtract(zh).divide(h2);
  3354.             if (e1.getReal() > 70. || e2.getReal() > 70.) {
  3355.                 return field.getOne();
  3356.             } else if (e1.getReal() < -70. && e2.getReal() < -70.) {
  3357.                 return zero.newInstance(FastMath.exp(r));
  3358.             } else {
  3359.                 final T ex1 = e1.exp();
  3360.                 final T ex2 = e2.exp();
  3361.                 return ex1.add(ex2).multiply(0.5).add(1).reciprocal().multiply(r).exp();
  3362.             }
  3363.         }

  3364.         /** Calculates scale height.
  3365.          * @param alt altitude
  3366.          * @param xm species molecular weight
  3367.          * @param temp temperature
  3368.          * @return scale height (km)
  3369.          */
  3370.         private T scalh(final double alt, final double xm, final double temp) {
  3371.             // Gravity at altitude
  3372.             final T denom = rlat.reciprocal().multiply(alt).add(1);
  3373.             final T galt = glat.divide(denom.square());
  3374.             return galt.reciprocal().multiply(R_GAS * temp / xm);
  3375.         }

  3376.         /** Calculates turbopause correction for MSIS models.
  3377.          * @param dd diffusive density
  3378.          * @param dm full mixed density
  3379.          * @param zhm transition scale length
  3380.          * @param xmm full mixed molecular weight
  3381.          * @param xm species molecular weight
  3382.          * @return combined density
  3383.          */
  3384.         private T dnet(final T dd, final T dm, final double zhm, final double xmm, final double xm) {
  3385.             if (!(dm.getReal() > 0 && dd.getReal() > 0)) {
  3386.                 T ddd = dd;
  3387.                 if (dd.getReal() == 0 && dm.getReal() == 0) {
  3388.                     ddd = field.getOne();
  3389.                 }
  3390.                 if (dm.getReal() == 0) {
  3391.                     return ddd;
  3392.                 }
  3393.                 if (dd.getReal() == 0) {
  3394.                     return dm;
  3395.                 }
  3396.             }

  3397.             final double a  = zhm / (xmm - xm);
  3398.             final T ylog = dm.divide(dd).log().multiply(a);
  3399.             if (ylog.getReal() < -10.) {
  3400.                 return dd;
  3401.             } else if (ylog.getReal() > 10.) {
  3402.                 return dm;
  3403.             } else {
  3404.                 return ylog.exp().add(1).pow(1.0 / a).multiply(dd);
  3405.             }
  3406.         }

  3407.         /** Integrate cubic spline function from xa[0] to x.
  3408.          * <p>ADAPTED FROM NUMERICAL RECIPES</p>
  3409.          * @param xa array of abscissas in ascending order
  3410.          * @param ya array of ordinates in ascending order by xa
  3411.          * @param y2a array of second derivatives in ascending order by xa
  3412.          * @param x abscissa end point
  3413.          * @return integral value
  3414.          */
  3415.         private T splini(final T[] xa, final T[] ya, final T[] y2a, final T x) {
  3416.             final int n = xa.length;
  3417.             T yi = zero;
  3418.             int klo = 0;
  3419.             int khi = 1;
  3420.             while (x.getReal() > xa[klo].getReal() && khi < n) {
  3421.                 T xx = x;
  3422.                 if (khi < n - 1) {
  3423.                     xx = (x.getReal() < xa[khi].getReal()) ? x : xa[khi];
  3424.                 }
  3425.                 final T h = xa[khi].subtract(xa[klo]);
  3426.                 final T a = xa[khi].subtract(xx).divide(h);
  3427.                 final T b = xx.subtract(xa[klo]).divide(h);
  3428.                 final T a2 = a.square();
  3429.                 final T b2 = b.square();

  3430.                 final T z =
  3431.                            a2.divide(2).subtract(a2.square().add(1).divide(4)).multiply(y2a[klo]).
  3432.                            add(b2.multiply(b2).divide(4).subtract(b2.divide(2)).multiply(y2a[khi]));
  3433.                 yi = yi.add(    a2.negate().add(1).multiply(ya[klo]).divide(2).
  3434.                             add(b2.multiply(ya[khi]).divide(2)).
  3435.                             add(z.multiply(h).multiply(h).divide(6)).
  3436.                             multiply(h));
  3437.                 klo++;
  3438.                 khi++;
  3439.             }
  3440.             return yi;
  3441.         }

  3442.         /** Calculate cubic spline interpolated value.
  3443.          * <p>ADAPTED FROM NUMERICAL RECIPES</p>
  3444.          * @param xa array of abscissas in ascending order
  3445.          * @param ya array of ordinates in ascending order by xa
  3446.          * @param y2a array of second derivatives in ascending order by xa
  3447.          * @param x abscissa for interpolation
  3448.          * @return interpolated value
  3449.          */
  3450.         private T splint(final T[] xa, final T[] ya, final T[] y2a, final T x) {
  3451.             final int n = xa.length;
  3452.             int klo = 0;
  3453.             int khi = n - 1;
  3454.             while (khi - klo > 1) {
  3455.                 final int k = (khi + klo) >>> 1;
  3456.                 if (xa[k].getReal() > x.getReal()) {
  3457.                     khi = k;
  3458.                 } else {
  3459.                     klo = k;
  3460.                 }
  3461.             }
  3462.             final T h = xa[khi].subtract(xa[klo]);
  3463.             final T a = xa[khi].subtract(x).divide(h);
  3464.             final T b = x.subtract(xa[klo]).divide(h);
  3465.             return a.multiply(ya[klo]).add(b.multiply(ya[khi])).
  3466.                    add((    a.square().multiply(a).subtract(a).multiply(y2a[klo]).
  3467.                         add(b.multiply(b).multiply(b).subtract(b).multiply(y2a[khi]))
  3468.                        ).multiply(h).multiply(h).divide(6));
  3469.         }

  3470.         /** Calculate 2nd derivatives of cubic spline interpolation function.
  3471.          * <p>ADAPTED FROM NUMERICAL RECIPES</p>
  3472.          * @param x array of abscissas in ascending order
  3473.          * @param y array of ordinates in ascending order by x
  3474.          * @param yp1 derivative at x[0] (2nd derivatives null if > 1E30)
  3475.          * @param ypn derivative at x[n-1] (2nd derivatives null if > 1E30)
  3476.          * @return array of second derivatives
  3477.          */
  3478.         private T[] spline(final T[] x, final T[] y, final T yp1, final T ypn) {
  3479.             final int n = x.length;
  3480.             final T[] y2 = MathArrays.buildArray(field, n);
  3481.             final T[] u  = MathArrays.buildArray(field, n);

  3482.             if (yp1.getReal() < 1e+30) {
  3483.                 y2[0] = zero.newInstance(-0.5);
  3484.                 final T dx = x[1].subtract(x[0]);
  3485.                 final T dy = y[1].subtract(y[0]);
  3486.                 u[0]  = dx.reciprocal().multiply(3.0).multiply(dy.divide(dx).subtract(yp1));
  3487.             }
  3488.             for (int i = 1; i < n - 1; i++) {
  3489.                 final T dx0m = x[i].subtract(x[i - 1]);
  3490.                 final T dy0m = y[i].subtract(y[i - 1]);
  3491.                 final T dxpm = x[i + 1].subtract(x[i - 1]);
  3492.                 final T dxp0 = x[i + 1].subtract(x[i]);
  3493.                 final T dyp0 = y[i + 1].subtract(y[i]);
  3494.                 final T sig = dx0m.divide(dxpm);
  3495.                 final T p = sig.multiply(y2[i - 1]).add(2.0);
  3496.                 y2[i] = sig.subtract(1.0).divide(p);
  3497.                 u[i] = dyp0.divide(dxp0).subtract(dy0m.divide(dx0m)).multiply(6).divide(dxpm).subtract(sig.multiply(u[i - 1])).divide(p);
  3498.             }

  3499.             double qn = 0;
  3500.             T un = zero;
  3501.             if (ypn.getReal() < 1e+30) {
  3502.                 final T dx12 = x[n - 1].subtract(x[n - 2]);
  3503.                 final T dy12 = y[n - 1].subtract(y[n - 2]);
  3504.                 qn = 0.5;
  3505.                 un = dx12.reciprocal().multiply(3.0).multiply(ypn.subtract(dy12.divide(dx12)));
  3506.             }

  3507.             y2[n - 1] = un.subtract(u[n - 2].multiply(qn)).divide(y2[n - 2].multiply(qn).add(1.0));
  3508.             for (int k = n - 2; k >= 0; k--) {
  3509.                 y2[k] = y2[k].multiply(y2[k + 1]).add(u[k]);
  3510.             }

  3511.             return y2;

  3512.         }

  3513.         /** Calculate Temperature and Density Profiles for lower atmosphere.
  3514.          * @param alt altitude
  3515.          * @param d0 density
  3516.          * @param xm mixed density
  3517.          * @return temperature or density profile
  3518.          */
  3519.         private T densm(final T alt, final T d0, final double xm) {

  3520.             T densm = d0;

  3521.             // stratosphere/mesosphere temperature
  3522.             int mn = ZN2.length;
  3523.             T z = (alt.getReal() > ZN2[mn - 1]) ? alt : zero.newInstance(ZN2[mn - 1]);

  3524.             double z1 = ZN2[0];
  3525.             double z2 = ZN2[mn - 1];
  3526.             T t1 = meso_tn2[0];
  3527.             T t2 = meso_tn2[mn - 1];
  3528.             T zg  = zeta(z, z1);
  3529.             T zgdif = zeta(zero.newInstance(z2), z1);

  3530.             /* set up spline nodes */
  3531.             T[] xs = MathArrays.buildArray(field, mn);
  3532.             T[] ys = MathArrays.buildArray(field, mn);
  3533.             for (int k = 0; k < mn; k++) {
  3534.                 xs[k] = zeta(zero.newInstance(ZN2[k]), z1).divide(zgdif);
  3535.                 ys[k] = meso_tn2[k].reciprocal();
  3536.             }
  3537.             final T qSM = rlat.add(z2).divide(rlat.add(z1));
  3538.             T yd1 = meso_tgn2[0].negate().divide(t1.square()).multiply(zgdif);
  3539.             T yd2 = meso_tgn2[1].negate().divide(t2.square()).multiply(zgdif).multiply(qSM.square());

  3540.             /* calculate spline coefficients */
  3541.             T[] y2out = spline(xs, ys, yd1, yd2);
  3542.             T x = zg.divide(zgdif);
  3543.             T y = splint(xs, ys, y2out, x);

  3544.             /* temperature at altitude */
  3545.             T tz = y.reciprocal();

  3546.             if (xm != 0.0) {
  3547.                 /* calculate stratosphere / mesospehere density */
  3548.                 final T glb  = galt(zero.newInstance(z1));
  3549.                 final T gamm = glb.multiply(zgdif).multiply(xm / R_GAS);

  3550.                 /* Integrate temperature profile */
  3551.                 final T yi = splini(xs, ys, y2out, x);
  3552.                 final T expl = min(MIN_TEMP, gamm.multiply(yi));

  3553.                 /* Density at altitude */
  3554.                 densm = densm.multiply(t1.divide(tz).multiply(expl.negate().exp()));
  3555.             }

  3556.             if (alt.getReal() > ZN3[0]) {
  3557.                 return (xm == 0.0) ? tz : densm;
  3558.             }

  3559.             // troposhere/stratosphere temperature
  3560.             z = alt;
  3561.             mn = ZN3.length;
  3562.             z1 = ZN3[0];
  3563.             z2 = ZN3[mn - 1];
  3564.             t1 = meso_tn3[0];
  3565.             t2 = meso_tn3[mn - 1];
  3566.             zg = zeta(z, z1);
  3567.             zgdif = zeta(zero.newInstance(z2), z1);

  3568.             /* set up spline nodes */
  3569.             xs = MathArrays.buildArray(field, mn);
  3570.             ys = MathArrays.buildArray(field, mn);
  3571.             for (int k = 0; k < mn; k++) {
  3572.                 xs[k] = zeta(zero.newInstance(ZN3[k]), z1).divide(zgdif);
  3573.                 ys[k] = meso_tn3[k].reciprocal();
  3574.             }
  3575.             final T qTS = rlat.add(z2) .divide(rlat.add(z1));
  3576.             yd1 = meso_tgn3[0].negate().divide(t1.multiply(t1)).multiply(zgdif);
  3577.             yd2 = meso_tgn3[1].negate().divide(t2.multiply(t2)).multiply(zgdif).multiply(qTS).multiply(qTS);

  3578.             /* calculate spline coefficients */
  3579.             y2out = spline(xs, ys, yd1, yd2);
  3580.             x = zg.divide(zgdif);
  3581.             y = splint(xs, ys, y2out, x);

  3582.             /* temperature at altitude */
  3583.             tz = y.reciprocal();

  3584.             if (xm != 0.0) {
  3585.                 /* calculate tropospheric / stratosphere density */
  3586.                 final T glb = galt(zero.newInstance(z1));
  3587.                 final T gamm = glb.multiply(zgdif).multiply(xm / R_GAS);

  3588.                 /* Integrate temperature profile */
  3589.                 final T yi = splini(xs, ys, y2out, x);
  3590.                 final T expl = min(MIN_TEMP, gamm.multiply(yi));

  3591.                 /* Density at altitude */
  3592.                 densm = densm.multiply(t1.divide(tz).multiply(expl.negate().exp()));
  3593.             }

  3594.             return (xm == 0.0) ? tz : densm;
  3595.         }

  3596.         /** Calculate temperature and density profiles according to new lower thermo polynomial.
  3597.          * @param alt altitude
  3598.          * @param dlb density at lower boundary
  3599.          * @param tinf exospheric temperature
  3600.          * @param tlb temperature at lower boundary
  3601.          * @param xm species molecular weight
  3602.          * @param alpha thermal diffusion coefficient
  3603.          * @param zlb altitude of the lower boundary
  3604.          * @param s2 slope
  3605.          * @return temperature or density profile
  3606.          */
  3607.         private T densu(final T alt, final T dlb, final T tinf,
  3608.                         final T tlb, final double xm,  final double alpha,
  3609.                         final double zlb, final T s2) {
  3610.             /* joining altitudes of Bates and spline */
  3611.             T z = (alt.getReal() > ZN1[0]) ? alt : zero.newInstance(ZN1[0]);

  3612.             /* geopotential altitude difference from ZLB */
  3613.             final T zg2 = zeta(z, zlb);

  3614.             /* Bates temperature */
  3615.             final T tt = tinf.subtract(tinf.subtract(tlb).multiply(s2.negate().multiply(zg2).exp()));
  3616.             final T ta = tt;
  3617.             T tz = tt;

  3618.             final int mn = ZN1.length;
  3619.             final T[] xs = MathArrays.buildArray(field, mn);
  3620.             final T[] ys = MathArrays.buildArray(field, mn);
  3621.             T x = zero;
  3622.             T[] y2out =  MathArrays.buildArray(field, mn);
  3623.             T zgdif = zero;
  3624.             if (alt.getReal() < ZN1[0]) {
  3625.                 /* calculate temperature below ZA
  3626.                  * temperature gradient at ZA from Bates profile */
  3627.                 final T p = rlat.add(zlb).divide(rlat.add(ZN1[0]));
  3628.                 final T dta = tinf.subtract(ta).multiply(s2).multiply(p.square());
  3629.                 meso_tgn1[0] = dta;
  3630.                 meso_tn1[0] = ta;
  3631.                 final T tzn1mn1 = zero.newInstance(ZN1[mn - 1]);
  3632.                 z = (alt.getReal() > ZN1[mn - 1]) ? alt : tzn1mn1;

  3633.                 final T t1 = meso_tn1[0];
  3634.                 final T t2 = meso_tn1[mn - 1];
  3635.                 /* geopotental difference from z1 */
  3636.                 final T zg = zeta(z, ZN1[0]);
  3637.                 zgdif = zeta(tzn1mn1, ZN1[0]);
  3638.                 /* set up spline nodes */
  3639.                 for (int k = 0; k < mn; k++) {
  3640.                     xs[k] = zeta(zero.newInstance(ZN1[k]), ZN1[0]).divide(zgdif);
  3641.                     ys[k] =  meso_tn1[k].reciprocal();
  3642.                 }
  3643.                 /* end node derivatives */
  3644.                 final T q   = rlat.add(ZN1[mn - 1]).divide(rlat.add(ZN1[0]));
  3645.                 final T yd1 = meso_tgn1[0].negate().divide(t1.square()).multiply(zgdif);
  3646.                 final T yd2 = meso_tgn1[1].negate().divide(t2.square()).multiply(zgdif).multiply(q.square());
  3647.                 /* calculate spline coefficients */
  3648.                 y2out = spline(xs, ys, yd1, yd2);
  3649.                 x = zg.divide(zgdif);
  3650.                 final T y = splint(xs, ys, y2out, x);
  3651.                 /* temperature at altitude */
  3652.                 tz = y.reciprocal();
  3653.             }

  3654.             if (xm == 0) {
  3655.                 return tz;
  3656.             }

  3657.             /* calculate density above za */
  3658.             T glb   = galt(zero.newInstance(zlb));
  3659.             T gamma = glb.divide(s2.multiply(tinf)).multiply(xm / R_GAS);
  3660.             T expl = tt.getReal() <= 0 ?
  3661.                      zero.newInstance(MIN_TEMP) :
  3662.                      min(MIN_TEMP, s2.negate().multiply(gamma).multiply(zg2).exp());
  3663.             T densu = dlb.multiply(expl).multiply(tlb.divide(tt).pow(gamma.add(alpha + 1)));

  3664.             // Correction for issue 1365 - protection against "densu" being infinite
  3665.             if (!Double.isFinite(densu.getReal())) {
  3666.                 if (expl.getReal() < MIN_TEMP) {
  3667.                     densu = dlb.multiply(FastMath.exp((FastMath.log(tlb.divide(tt)).multiply(gamma.add(alpha + 1))).
  3668.                                                       subtract(s2.multiply(gamma).multiply(zg2))));
  3669.                 } else {
  3670.                     throw new OrekitException(OrekitMessages.INFINITE_NRLMSISE00_DENSITY);
  3671.                 }
  3672.             }

  3673.             /* calculate density below za */
  3674.             if (alt.getReal() < ZN1[0]) {
  3675.                 glb   = galt(zero.newInstance(ZN1[0]));
  3676.                 gamma = glb.multiply(zgdif).multiply(xm / R_GAS);
  3677.                 /* integrate spline temperatures */
  3678.                 expl = tz.getReal() <= 0 ?
  3679.                        zero.newInstance(MIN_TEMP) :
  3680.                        min(MIN_TEMP, gamma.multiply(splini(xs, ys, y2out, x)));
  3681.                 /* correct density at altitude */
  3682.                 densu = densu.multiply(meso_tn1[0].divide(tz).pow(alpha + 1).multiply(expl.negate().exp()));
  3683.             }

  3684.             /* Return density at altitude */
  3685.             return densu;
  3686.         }

  3687.         /** Compute min of two values, one double and one field element.
  3688.          * @param d double value
  3689.          * @param f field element
  3690.          * @return min value
  3691.          */
  3692.         private T min(final double d, final T f) {
  3693.             return (f.getReal() > d) ? zero.newInstance(d) : f;
  3694.         }

  3695.         /** Calculate gravity at altitude.
  3696.          * @param alt altitude (km)
  3697.          * @return gravity at altitude (cm/s2)
  3698.          */
  3699.         private T galt(final T alt) {
  3700.             final T r = alt.divide(rlat).add(1);
  3701.             return glat.divide(r.square());
  3702.         }

  3703.         /** Calculate zeta function.
  3704.          * @param zz zz value
  3705.          * @param zl zl value
  3706.          * @return value of zeta function
  3707.          */
  3708.         private T zeta(final T zz, final double zl) {
  3709.             return zz.subtract(zl).multiply(rlat.add(zl)).divide(rlat.add(zz));
  3710.         }

  3711.     }

  3712. }