]> git.mxchange.org Git - flightgear.git/blob - src/FDM/JSBSim/models/atmosphere/FGMSIS.cpp
remove a cvs conflict
[flightgear.git] / src / FDM / JSBSim / models / atmosphere / FGMSIS.cpp
1 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
2
3  Module:       FGMSIS.cpp
4  Author:       David Culp
5                (incorporated into C++ JSBSim class heirarchy, see model authors below)
6  Date started: 12/14/03
7  Purpose:      Models the MSIS-00 atmosphere
8
9  ------------- Copyright (C) 2003  David P. Culp (davidculp2@comcast.net) ------
10
11  This program is free software; you can redistribute it and/or modify it under
12  the terms of the GNU Lesser General Public License as published by the Free Software
13  Foundation; either version 2 of the License, or (at your option) any later
14  version.
15
16  This program is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18  FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public License for more
19  details.
20
21  You should have received a copy of the GNU Lesser General Public License along with
22  this program; if not, write to the Free Software Foundation, Inc., 59 Temple
23  Place - Suite 330, Boston, MA  02111-1307, USA.
24
25  Further information about the GNU Lesser General Public License can also be found on
26  the world wide web at http://www.gnu.org.
27
28 FUNCTIONAL DESCRIPTION
29 --------------------------------------------------------------------------------
30 Models the MSIS-00 atmosphere. Provides temperature and density to FGAtmosphere,
31 given day-of-year, time-of-day, altitude, latitude, longitude and local time.
32
33 HISTORY
34 --------------------------------------------------------------------------------
35 12/14/03   DPC   Created
36 01/11/04   DPC   Derived from FGAtmosphere
37
38  --------------------------------------------------------------------
39  ---------  N R L M S I S E - 0 0    M O D E L    2 0 0 1  ----------
40  --------------------------------------------------------------------
41
42  This file is part of the NRLMSISE-00  C source code package - release
43  20020503
44
45  The NRLMSISE-00 model was developed by Mike Picone, Alan Hedin, and
46  Doug Drob. They also wrote a NRLMSISE-00 distribution package in
47  FORTRAN which is available at
48  http://uap-www.nrl.navy.mil/models_web/msis/msis_home.htm
49
50  Dominik Brodowski implemented and maintains this C version. You can
51  reach him at devel@brodo.de. See the file "DOCUMENTATION" for details,
52  and check http://www.brodo.de/english/pub/nrlmsise/index.html for
53  updated releases of this package.
54 */
55
56 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
57 INCLUDES
58 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
59
60 #include "FGMSIS.h"
61 #include "models/FGAuxiliary.h"
62 #include <cmath>          /* maths functions */
63 #include <iostream>        // for cout, endl
64
65 using namespace std;
66
67 namespace JSBSim {
68
69 static const char *IdSrc = "$Id: FGMSIS.cpp,v 1.18 2011/06/21 13:54:40 jberndt Exp $";
70 static const char *IdHdr = ID_MSIS;
71
72 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
73 EXTERNAL GLOBAL DATA
74 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
75
76   /* POWER7 */
77   extern double pt[150];
78   extern double pd[9][150];
79   extern double ps[150];
80   extern double pdl[2][25];
81   extern double ptl[4][100];
82   extern double pma[10][100];
83   extern double sam[100];
84
85   /* LOWER7 */
86   extern double ptm[10];
87   extern double pdm[8][10];
88   extern double pavgm[10];
89
90 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
91 CLASS IMPLEMENTATION
92 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
93
94
95 MSIS::MSIS(FGFDMExec* fdmex) : FGAtmosphere(fdmex)
96 {
97   Name = "MSIS";
98
99   for (int i=0; i<9; i++) output.d[i] = 0.0;
100   for (int i=0; i<2; i++) output.t[i] = 0.0;
101
102   dm04 = dm16 = dm28 = dm32 = dm40 = dm01 = dm14 = dfa = 0.0;
103
104   for (int i=0; i<5; i++) meso_tn1[i] = 0.0;
105   for (int i=0; i<4; i++) meso_tn2[i] = 0.0;
106   for (int i=0; i<5; i++) meso_tn3[i] = 0.0;
107   for (int i=0; i<2; i++) meso_tgn1[i] = 0.0;
108   for (int i=0; i<2; i++) meso_tgn2[i] = 0.0;
109   for (int i=0; i<2; i++) meso_tgn3[i] = 0.0;
110
111   Debug(0);
112 }
113
114 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
115
116 MSIS::~MSIS()
117 {
118   Debug(1);
119 }
120
121 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
122
123 bool MSIS::InitModel(void)
124 {
125   unsigned int i;
126
127   flags.switches[0] = 0;
128   flags.sw[0] = 0;
129   flags.swc[0] = 0;
130   for (i=1;i<24;i++) {
131     flags.switches[i] = 1;
132     flags.sw[i] = 1;
133     flags.swc[i] = 1;
134   }
135
136   for (i=0;i<7;i++) aph.a[i] = 100.0;
137
138   // set some common magnetic flux values
139   input.f107A = 150.0;
140   input.f107 = 150.0;
141   input.ap = 4.0;
142
143 //  UseInternal();
144
145 //  SLtemperature = intTemperature = 518.0;
146 //  SLpressure    = intPressure = 2116.7;
147 //  SLdensity     = intDensity = 0.002378;
148 //  SLsoundspeed  = sqrt(2403.0832 * SLtemperature);
149 //  rSLtemperature = 1.0/intTemperature;
150 //  rSLpressure    = 1.0/intPressure;
151 //  rSLdensity     = 1.0/intDensity;
152 //  rSLsoundspeed  = 1.0/SLsoundspeed;
153
154   return true;
155 }
156
157 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
158
159 bool MSIS::Run(bool Holding)
160 {
161   if (FGModel::Run(Holding)) return true;
162   if (Holding) return false;
163
164   RunPreFunctions();
165
166   double h = FDMExec->GetPropagate()->GetAltitudeASL();
167
168   //do temp, pressure, and density first
169 //  if (!useExternal) {
170     // get sea-level values
171     Calculate(FDMExec->GetAuxiliary()->GetDayOfYear(),
172               FDMExec->GetAuxiliary()->GetSecondsInDay(),
173               0.0,
174               FDMExec->GetPropagate()->GetLocation().GetLatitudeDeg(),
175               FDMExec->GetPropagate()->GetLocation().GetLongitudeDeg());
176     SLtemperature = output.t[1] * 1.8;
177     SLdensity     = output.d[5] * 1.940321;
178     SLpressure    = 1716.488 * SLdensity * SLtemperature;
179     SLsoundspeed  = sqrt(2403.0832 * SLtemperature);
180     rSLtemperature = 1.0/SLtemperature;
181     rSLpressure    = 1.0/SLpressure;
182     rSLdensity     = 1.0/SLdensity;
183     rSLsoundspeed  = 1.0/SLsoundspeed;
184
185     // get at-altitude values
186     Calculate(FDMExec->GetAuxiliary()->GetDayOfYear(),
187               FDMExec->GetAuxiliary()->GetSecondsInDay(),
188               h,
189               FDMExec->GetPropagate()->GetLocation().GetLatitudeDeg(),
190               FDMExec->GetPropagate()->GetLocation().GetLongitudeDeg());
191 //    intTemperature = output.t[1] * 1.8;
192 //    intDensity     = output.d[5] * 1.940321;
193 //    intPressure    = 1716.488 * intDensity * intTemperature;
194     //cout << "T=" << intTemperature << " D=" << intDensity << " P=";
195     //cout << intPressure << " a=" << soundspeed << endl;
196 //  }
197
198   RunPostFunctions();
199
200   Debug(2);
201
202   return false;
203 }
204
205 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
206
207 void MSIS::Calculate(int day, double sec, double alt, double lat, double lon)
208 {
209   input.year = 2000;
210   input.doy = day;
211   input.sec = sec;
212   input.alt = alt / 3281;  //feet to kilometers
213   input.g_lat = lat;
214   input.g_long = lon;
215
216   input.lst = (sec/3600) + (lon/15);
217   if (input.lst > 24.0) input.lst -= 24.0;
218   if (input.lst < 0.0) input.lst = 24 - input.lst;
219
220   gtd7d(&input, &flags, &output);
221 }
222
223 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
224
225
226 void MSIS::UseExternal(void){
227   // do nothing, external control not allowed
228 }
229
230
231 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
232
233
234 void MSIS::tselec(struct nrlmsise_flags *flags)
235 {
236   int i;
237   for (i=0;i<24;i++) {
238     if (i!=9) {
239       if (flags->switches[i]==1)
240         flags->sw[i]=1;
241       else
242         flags->sw[i]=0;
243       if (flags->switches[i]>0)
244         flags->swc[i]=1;
245       else
246         flags->swc[i]=0;
247     } else {
248       flags->sw[i]=flags->switches[i];
249       flags->swc[i]=flags->switches[i];
250     }
251   }
252 }
253
254
255 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
256
257 void MSIS::glatf(double lat, double *gv, double *reff)
258 {
259   double dgtr = 1.74533E-2;
260   double c2;
261   c2 = cos(2.0*dgtr*lat);
262   *gv = 980.616 * (1.0 - 0.0026373 * c2);
263   *reff = 2.0 * (*gv) / (3.085462E-6 + 2.27E-9 * c2) * 1.0E-5;
264 }
265
266 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
267
268 double MSIS::ccor(double alt, double r, double h1, double zh)
269 {
270 /*        CHEMISTRY/DISSOCIATION CORRECTION FOR MSIS MODELS
271  *         ALT - altitude
272  *         R - target ratio
273  *         H1 - transition scale length
274  *         ZH - altitude of 1/2 R
275  */
276   double e;
277   double ex;
278   e = (alt - zh) / h1;
279   if (e>70)
280     return exp(0.0);
281   if (e<-70)
282     return exp(r);
283   ex = exp(e);
284   e = r / (1.0 + ex);
285   return exp(e);
286 }
287
288 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
289
290 double MSIS::ccor2(double alt, double r, double h1, double zh, double h2)
291 {
292 /*        CHEMISTRY/DISSOCIATION CORRECTION FOR MSIS MODELS
293  *         ALT - altitude
294  *         R - target ratio
295  *         H1 - transition scale length
296  *         ZH - altitude of 1/2 R
297  *         H2 - transition scale length #2 ?
298  */
299   double e1, e2;
300   double ex1, ex2;
301   double ccor2v;
302   e1 = (alt - zh) / h1;
303   e2 = (alt - zh) / h2;
304   if ((e1 > 70) || (e2 > 70))
305     return exp(0.0);
306   if ((e1 < -70) && (e2 < -70))
307     return exp(r);
308   ex1 = exp(e1);
309   ex2 = exp(e2);
310   ccor2v = r / (1.0 + 0.5 * (ex1 + ex2));
311   return exp(ccor2v);
312 }
313
314 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
315
316 double MSIS::scalh(double alt, double xm, double temp)
317 {
318   double g;
319   double rgas=831.4;
320   g = gsurf / (pow((1.0 + alt/re),2.0));
321   g = rgas * temp / (g * xm);
322   return g;
323 }
324
325 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
326
327 double MSIS::dnet (double dd, double dm, double zhm, double xmm, double xm)
328 {
329 /*       TURBOPAUSE CORRECTION FOR MSIS MODELS
330  *        Root mean density
331  *         DD - diffusive density
332  *         DM - full mixed density
333  *         ZHM - transition scale length
334  *         XMM - full mixed molecular weight
335  *         XM  - species molecular weight
336  *         DNET - combined density
337  */
338   double a;
339   double ylog;
340   a  = zhm / (xmm-xm);
341   if (!((dm>0) && (dd>0))) {
342     cerr << "dnet log error " << dm << ' ' << dd << ' ' << xm << ' ' << endl;
343     if ((dd==0) && (dm==0))
344       dd=1;
345     if (dm==0)
346       return dd;
347     if (dd==0)
348       return dm;
349   }
350   ylog = a * log(dm/dd);
351   if (ylog<-10)
352     return dd;
353   if (ylog>10)
354     return dm;
355   a = dd*pow((1.0 + exp(ylog)),(1.0/a));
356   return a;
357 }
358
359 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
360
361 void MSIS::splini (double *xa, double *ya, double *y2a, int n, double x, double *y)
362 {
363 /*      INTEGRATE CUBIC SPLINE FUNCTION FROM XA(1) TO X
364  *       XA,YA: ARRAYS OF TABULATED FUNCTION IN ASCENDING ORDER BY X
365  *       Y2A: ARRAY OF SECOND DERIVATIVES
366  *       N: SIZE OF ARRAYS XA,YA,Y2A
367  *       X: ABSCISSA ENDPOINT FOR INTEGRATION
368  *       Y: OUTPUT VALUE
369  */
370   double yi=0;
371   int klo=0;
372   int khi=1;
373   double xx=0.0, h=0.0, a=0.0, b=0.0, a2=0.0, b2=0.0;
374   while ((x>xa[klo]) && (khi<n)) {
375     xx=x;
376     if (khi<(n-1)) {
377       if (x<xa[khi])
378         xx=x;
379       else
380         xx=xa[khi];
381     }
382     h = xa[khi] - xa[klo];
383     a = (xa[khi] - xx)/h;
384     b = (xx - xa[klo])/h;
385     a2 = a*a;
386     b2 = b*b;
387     yi += ((1.0 - a2) * ya[klo] / 2.0 + b2 * ya[khi] / 2.0 + ((-(1.0+a2*a2)/4.0 + a2/2.0) * y2a[klo] + (b2*b2/4.0 - b2/2.0) * y2a[khi]) * h * h / 6.0) * h;
388     klo++;
389     khi++;
390   }
391   *y = yi;
392 }
393
394 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
395
396 void MSIS::splint (double *xa, double *ya, double *y2a, int n, double x, double *y)
397 {
398 /*      CALCULATE CUBIC SPLINE INTERP VALUE
399  *       ADAPTED FROM NUMERICAL RECIPES BY PRESS ET AL.
400  *       XA,YA: ARRAYS OF TABULATED FUNCTION IN ASCENDING ORDER BY X
401  *       Y2A: ARRAY OF SECOND DERIVATIVES
402  *       N: SIZE OF ARRAYS XA,YA,Y2A
403  *       X: ABSCISSA FOR INTERPOLATION
404  *       Y: OUTPUT VALUE
405  */
406   int klo=0;
407   int khi=n-1;
408   int k;
409   double h;
410   double a, b, yi;
411   while ((khi-klo)>1) {
412     k=(khi+klo)/2;
413     if (xa[k]>x)
414       khi=k;
415     else
416       klo=k;
417   }
418   h = xa[khi] - xa[klo];
419   if (h==0.0)
420     cerr << "bad XA input to splint" << endl;
421   a = (xa[khi] - x)/h;
422   b = (x - xa[klo])/h;
423   yi = a * ya[klo] + b * ya[khi] + ((a*a*a - a) * y2a[klo] + (b*b*b - b) * y2a[khi]) * h * h/6.0;
424   *y = yi;
425 }
426
427 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
428
429 void MSIS::spline (double *x, double *y, int n, double yp1, double ypn, double *y2)
430 {
431 /*       CALCULATE 2ND DERIVATIVES OF CUBIC SPLINE INTERP FUNCTION
432  *       ADAPTED FROM NUMERICAL RECIPES BY PRESS ET AL
433  *       X,Y: ARRAYS OF TABULATED FUNCTION IN ASCENDING ORDER BY X
434  *       N: SIZE OF ARRAYS X,Y
435  *       YP1,YPN: SPECIFIED DERIVATIVES AT X[0] AND X[N-1]; VALUES
436  *                >= 1E30 SIGNAL SIGNAL SECOND DERIVATIVE ZERO
437  *       Y2: OUTPUT ARRAY OF SECOND DERIVATIVES
438  */
439   double *u;
440   double sig, p, qn, un;
441   int i, k;
442   u=new double[n];
443   if (u==NULL) {
444     cerr << "Out Of Memory in spline - ERROR" << endl;
445     return;
446   }
447   if (yp1>0.99E30) {
448     y2[0]=0;
449     u[0]=0;
450   } else {
451     y2[0]=-0.5;
452     u[0]=(3.0/(x[1]-x[0]))*((y[1]-y[0])/(x[1]-x[0])-yp1);
453   }
454   for (i=1;i<(n-1);i++) {
455     sig = (x[i]-x[i-1])/(x[i+1] - x[i-1]);
456     p = sig * y2[i-1] + 2.0;
457     y2[i] = (sig - 1.0) / p;
458     u[i] = (6.0 * ((y[i+1] - y[i])/(x[i+1] - x[i]) -(y[i] - y[i-1]) / (x[i] - x[i-1]))/(x[i+1] - x[i-1]) - sig * u[i-1])/p;
459   }
460   if (ypn>0.99E30) {
461     qn = 0;
462     un = 0;
463   } else {
464     qn = 0.5;
465     un = (3.0 / (x[n-1] - x[n-2])) * (ypn - (y[n-1] - y[n-2])/(x[n-1] - x[n-2]));
466   }
467   y2[n-1] = (un - qn * u[n-2]) / (qn * y2[n-2] + 1.0);
468   for (k=n-2;k>=0;k--)
469     y2[k] = y2[k] * y2[k+1] + u[k];
470
471   delete u;
472 }
473
474 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
475
476 double MSIS::zeta(double zz, double zl)
477 {
478   return ((zz-zl)*(re+zl)/(re+zz));
479 }
480
481 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
482
483 double MSIS::densm(double alt, double d0, double xm, double *tz, int mn3,
484                      double *zn3, double *tn3, double *tgn3, int mn2, double *zn2,
485                      double *tn2, double *tgn2)
486 {
487 /*      Calculate Temperature and Density Profiles for lower atmos.  */
488   double xs[10] = {0,0,0,0,0,0,0,0,0,0};
489   double ys[10] = {0,0,0,0,0,0,0,0,0,0};
490   double y2out[10] = {0,0,0,0,0,0,0,0,0,0};
491   double rgas = 831.4;
492   double z=0, z1=0, z2=0, t1=0, t2=0, zg=0, zgdif=0;
493   double yd1=0, yd2=0;
494   double x=0, y=0, yi=0;
495   double expl=0, gamm=0, glb=0;
496   double densm_tmp=0;
497   int mn=0;
498   int k=0;
499   densm_tmp=d0;
500   if (alt>zn2[0]) {
501     if (xm==0.0)
502       return *tz;
503     else
504       return d0;
505   }
506
507   /* STRATOSPHERE/MESOSPHERE TEMPERATURE */
508   if (alt>zn2[mn2-1])
509     z=alt;
510   else
511     z=zn2[mn2-1];
512   mn=mn2;
513   z1=zn2[0];
514   z2=zn2[mn-1];
515   t1=tn2[0];
516   t2=tn2[mn-1];
517   zg = zeta(z, z1);
518   zgdif = zeta(z2, z1);
519
520   /* set up spline nodes */
521   for (k=0;k<mn;k++) {
522     xs[k]=zeta(zn2[k],z1)/zgdif;
523     ys[k]=1.0 / tn2[k];
524   }
525   yd1=-tgn2[0] / (t1*t1) * zgdif;
526   yd2=-tgn2[1] / (t2*t2) * zgdif * (pow(((re+z2)/(re+z1)),2.0));
527
528   /* calculate spline coefficients */
529   spline (xs, ys, mn, yd1, yd2, y2out);
530   x = zg/zgdif;
531   splint (xs, ys, y2out, mn, x, &y);
532
533   /* temperature at altitude */
534   *tz = 1.0 / y;
535   if (xm!=0.0) {
536     /* calaculate stratosphere / mesospehere density */
537     glb = gsurf / (pow((1.0 + z1/re),2.0));
538     gamm = xm * glb * zgdif / rgas;
539
540     /* Integrate temperature profile */
541     splini(xs, ys, y2out, mn, x, &yi);
542     expl=gamm*yi;
543     if (expl>50.0)
544       expl=50.0;
545
546     /* Density at altitude */
547     densm_tmp = densm_tmp * (t1 / *tz) * exp(-expl);
548   }
549
550   if (alt>zn3[0]) {
551     if (xm==0.0)
552       return *tz;
553     else
554       return densm_tmp;
555   }
556
557   /* troposhere / stratosphere temperature */
558   z = alt;
559   mn = mn3;
560   z1=zn3[0];
561   z2=zn3[mn-1];
562   t1=tn3[0];
563   t2=tn3[mn-1];
564   zg=zeta(z,z1);
565   zgdif=zeta(z2,z1);
566
567   /* set up spline nodes */
568   for (k=0;k<mn;k++) {
569     xs[k] = zeta(zn3[k],z1) / zgdif;
570     ys[k] = 1.0 / tn3[k];
571   }
572   yd1=-tgn3[0] / (t1*t1) * zgdif;
573   yd2=-tgn3[1] / (t2*t2) * zgdif * (pow(((re+z2)/(re+z1)),2.0));
574
575   /* calculate spline coefficients */
576   spline (xs, ys, mn, yd1, yd2, y2out);
577   x = zg/zgdif;
578   splint (xs, ys, y2out, mn, x, &y);
579
580   /* temperature at altitude */
581   *tz = 1.0 / y;
582   if (xm!=0.0) {
583     /* calaculate tropospheric / stratosphere density */
584     glb = gsurf / (pow((1.0 + z1/re),2.0));
585     gamm = xm * glb * zgdif / rgas;
586
587     /* Integrate temperature profile */
588     splini(xs, ys, y2out, mn, x, &yi);
589     expl=gamm*yi;
590     if (expl>50.0)
591       expl=50.0;
592
593     /* Density at altitude */
594     densm_tmp = densm_tmp * (t1 / *tz) * exp(-expl);
595   }
596   if (xm==0.0)
597     return *tz;
598   else
599     return densm_tmp;
600 }
601
602 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
603
604 double MSIS::densu(double alt, double dlb, double tinf, double tlb, double xm,
605                      double alpha, double *tz, double zlb, double s2, int mn1,
606                      double *zn1, double *tn1, double *tgn1)
607 {
608 /*      Calculate Temperature and Density Profiles for MSIS models
609  *      New lower thermo polynomial
610  */
611   double yd2=0.0, yd1=0.0, x=0.0, y=0.0;
612   double rgas=831.4;
613   double densu_temp=1.0;
614   double za=0.0, z=0.0, zg2=0.0, tt=0.0, ta=0.0;
615   double dta=0.0, z1=0.0, z2=0.0, t1=0.0, t2=0.0, zg=0.0, zgdif=0.0;
616   int mn=0;
617   int k=0;
618   double glb=0.0;
619   double expl=0.0;
620   double yi=0.0;
621   double densa=0.0;
622   double gamma=0.0, gamm=0.0;
623   double xs[5]={0.0,0.0,0.0,0.0,0.0}, ys[5]={0.0,0.0,0.0,0.0,0.0}, y2out[5]={0.0,0.0,0.0,0.0,0.0};
624   /* joining altitudes of Bates and spline */
625   za=zn1[0];
626   if (alt>za)
627     z=alt;
628   else
629     z=za;
630
631   /* geopotential altitude difference from ZLB */
632   zg2 = zeta(z, zlb);
633
634   /* Bates temperature */
635   tt = tinf - (tinf - tlb) * exp(-s2*zg2);
636   ta = tt;
637   *tz = tt;
638   densu_temp = *tz;
639
640   if (alt<za) {
641     /* calculate temperature below ZA
642      * temperature gradient at ZA from Bates profile */
643     dta = (tinf - ta) * s2 * pow(((re+zlb)/(re+za)),2.0);
644     tgn1[0]=dta;
645     tn1[0]=ta;
646     if (alt>zn1[mn1-1])
647       z=alt;
648     else
649       z=zn1[mn1-1];
650     mn=mn1;
651     z1=zn1[0];
652     z2=zn1[mn-1];
653     t1=tn1[0];
654     t2=tn1[mn-1];
655     /* geopotental difference from z1 */
656     zg = zeta (z, z1);
657     zgdif = zeta(z2, z1);
658     /* set up spline nodes */
659     for (k=0;k<mn;k++) {
660       xs[k] = zeta(zn1[k], z1) / zgdif;
661       ys[k] = 1.0 / tn1[k];
662     }
663     /* end node derivatives */
664     yd1 = -tgn1[0] / (t1*t1) * zgdif;
665     yd2 = -tgn1[1] / (t2*t2) * zgdif * pow(((re+z2)/(re+z1)),2.0);
666     /* calculate spline coefficients */
667     spline (xs, ys, mn, yd1, yd2, y2out);
668     x = zg / zgdif;
669     splint (xs, ys, y2out, mn, x, &y);
670     /* temperature at altitude */
671     *tz = 1.0 / y;
672     densu_temp = *tz;
673   }
674   if (xm==0)
675     return densu_temp;
676
677   /* calculate density above za */
678   glb = gsurf / pow((1.0 + zlb/re),2.0);
679   gamma = xm * glb / (s2 * rgas * tinf);
680   expl = exp(-s2 * gamma * zg2);
681   if (expl>50.0)
682       expl=50.0;
683   if (tt<=0)
684     expl=50.0;
685
686   /* density at altitude */
687   densa = dlb * pow((tlb/tt),((1.0+alpha+gamma))) * expl;
688   densu_temp=densa;
689   if (alt>=za)
690     return densu_temp;
691
692   /* calculate density below za */
693   glb = gsurf / pow((1.0 + z1/re),2.0);
694   gamm = xm * glb * zgdif / rgas;
695
696   /* integrate spline temperatures */
697   splini (xs, ys, y2out, mn, x, &yi);
698   expl = gamm * yi;
699   if (expl>50.0)
700     expl=50.0;
701   if (*tz<=0)
702     expl=50.0;
703
704   /* density at altitude */
705   densu_temp = densu_temp * pow ((t1 / *tz),(1.0 + alpha)) * exp(-expl);
706   return densu_temp;
707 }
708
709 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
710
711 /*    3hr Magnetic activity functions */
712 /*    Eq. A24d */
713 double MSIS::g0(double a, double *p)
714 {
715   return (a - 4.0 + (p[25] - 1.0) * (a - 4.0 + (exp(-sqrt(p[24]*p[24]) *
716                 (a - 4.0)) - 1.0) / sqrt(p[24]*p[24])));
717 }
718
719 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
720
721 /*    Eq. A24c */
722 double MSIS::sumex(double ex)
723 {
724   return (1.0 + (1.0 - pow(ex,19.0)) / (1.0 - ex) * pow(ex,0.5));
725 }
726
727 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
728
729 /*    Eq. A24a */
730 double MSIS::sg0(double ex, double *p, double *ap)
731 {
732   return (g0(ap[1],p) + (g0(ap[2],p)*ex + g0(ap[3],p)*ex*ex +
733                 g0(ap[4],p)*pow(ex,3.0)  + (g0(ap[5],p)*pow(ex,4.0) +
734                 g0(ap[6],p)*pow(ex,12.0))*(1.0-pow(ex,8.0))/(1.0-ex)))/sumex(ex);
735 }
736
737 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
738
739 double MSIS::globe7(double *p, struct nrlmsise_input *input,
740                       struct nrlmsise_flags *flags)
741 {
742 /*       CALCULATE G(L) FUNCTION
743  *       Upper Thermosphere Parameters */
744   double t[15];
745   int i,j;
746   int sw9=1;
747   double apd;
748   double xlong;
749   double tloc;
750   double c, s, c2, c4, s2;
751   double sr = 7.2722E-5;
752   double dgtr = 1.74533E-2;
753   double dr = 1.72142E-2;
754   double hr = 0.2618;
755   double cd32, cd18, cd14, cd39;
756   double p32, p18, p14, p39;
757   double df;
758   double f1, f2;
759   double tinf;
760   struct ap_array *ap;
761
762   tloc=input->lst;
763   for (j=0;j<14;j++)
764     t[j]=0;
765   if (flags->sw[9]>0)
766     sw9=1;
767   else if (flags->sw[9]<0)
768     sw9=-1;
769   xlong = input->g_long;
770
771   /* calculate legendre polynomials */
772   c = sin(input->g_lat * dgtr);
773   s = cos(input->g_lat * dgtr);
774   c2 = c*c;
775   c4 = c2*c2;
776   s2 = s*s;
777
778   plg[0][1] = c;
779   plg[0][2] = 0.5*(3.0*c2 -1.0);
780   plg[0][3] = 0.5*(5.0*c*c2-3.0*c);
781   plg[0][4] = (35.0*c4 - 30.0*c2 + 3.0)/8.0;
782   plg[0][5] = (63.0*c2*c2*c - 70.0*c2*c + 15.0*c)/8.0;
783   plg[0][6] = (11.0*c*plg[0][5] - 5.0*plg[0][4])/6.0;
784 /*      plg[0][7] = (13.0*c*plg[0][6] - 6.0*plg[0][5])/7.0; */
785   plg[1][1] = s;
786   plg[1][2] = 3.0*c*s;
787   plg[1][3] = 1.5*(5.0*c2-1.0)*s;
788   plg[1][4] = 2.5*(7.0*c2*c-3.0*c)*s;
789   plg[1][5] = 1.875*(21.0*c4 - 14.0*c2 +1.0)*s;
790   plg[1][6] = (11.0*c*plg[1][5]-6.0*plg[1][4])/5.0;
791 /*      plg[1][7] = (13.0*c*plg[1][6]-7.0*plg[1][5])/6.0; */
792 /*      plg[1][8] = (15.0*c*plg[1][7]-8.0*plg[1][6])/7.0; */
793   plg[2][2] = 3.0*s2;
794   plg[2][3] = 15.0*s2*c;
795   plg[2][4] = 7.5*(7.0*c2 -1.0)*s2;
796   plg[2][5] = 3.0*c*plg[2][4]-2.0*plg[2][3];
797   plg[2][6] =(11.0*c*plg[2][5]-7.0*plg[2][4])/4.0;
798   plg[2][7] =(13.0*c*plg[2][6]-8.0*plg[2][5])/5.0;
799   plg[3][3] = 15.0*s2*s;
800   plg[3][4] = 105.0*s2*s*c;
801   plg[3][5] =(9.0*c*plg[3][4]-7.*plg[3][3])/2.0;
802   plg[3][6] =(11.0*c*plg[3][5]-8.*plg[3][4])/3.0;
803
804   if (!(((flags->sw[7]==0)&&(flags->sw[8]==0))&&(flags->sw[14]==0))) {
805     stloc = sin(hr*tloc);
806     ctloc = cos(hr*tloc);
807     s2tloc = sin(2.0*hr*tloc);
808     c2tloc = cos(2.0*hr*tloc);
809     s3tloc = sin(3.0*hr*tloc);
810     c3tloc = cos(3.0*hr*tloc);
811   }
812
813   cd32 = cos(dr*(input->doy-p[31]));
814   cd18 = cos(2.0*dr*(input->doy-p[17]));
815   cd14 = cos(dr*(input->doy-p[13]));
816   cd39 = cos(2.0*dr*(input->doy-p[38]));
817   p32=p[31];
818   p18=p[17];
819   p14=p[13];
820   p39=p[38];
821
822   /* F10.7 EFFECT */
823   df = input->f107 - input->f107A;
824   dfa = input->f107A - 150.0;
825   t[0] =  p[19]*df*(1.0+p[59]*dfa) + p[20]*df*df + p[21]*dfa + p[29]*pow(dfa,2.0);
826   f1 = 1.0 + (p[47]*dfa +p[19]*df+p[20]*df*df)*flags->swc[1];
827   f2 = 1.0 + (p[49]*dfa+p[19]*df+p[20]*df*df)*flags->swc[1];
828
829   /*  TIME INDEPENDENT */
830   t[1] = (p[1]*plg[0][2]+ p[2]*plg[0][4]+p[22]*plg[0][6]) +
831         (p[14]*plg[0][2])*dfa*flags->swc[1] +p[26]*plg[0][1];
832
833   /*  SYMMETRICAL ANNUAL */
834   t[2] = p[18]*cd32;
835
836   /*  SYMMETRICAL SEMIANNUAL */
837   t[3] = (p[15]+p[16]*plg[0][2])*cd18;
838
839   /*  ASYMMETRICAL ANNUAL */
840   t[4] =  f1*(p[9]*plg[0][1]+p[10]*plg[0][3])*cd14;
841
842   /*  ASYMMETRICAL SEMIANNUAL */
843   t[5] =    p[37]*plg[0][1]*cd39;
844
845         /* DIURNAL */
846   if (flags->sw[7]) {
847     double t71, t72;
848     t71 = (p[11]*plg[1][2])*cd14*flags->swc[5];
849     t72 = (p[12]*plg[1][2])*cd14*flags->swc[5];
850     t[6] = f2*((p[3]*plg[1][1] + p[4]*plg[1][3] + p[27]*plg[1][5] + t71) * \
851          ctloc + (p[6]*plg[1][1] + p[7]*plg[1][3] + p[28]*plg[1][5] \
852             + t72)*stloc);
853 }
854
855   /* SEMIDIURNAL */
856   if (flags->sw[8]) {
857     double t81, t82;
858     t81 = (p[23]*plg[2][3]+p[35]*plg[2][5])*cd14*flags->swc[5];
859     t82 = (p[33]*plg[2][3]+p[36]*plg[2][5])*cd14*flags->swc[5];
860     t[7] = f2*((p[5]*plg[2][2]+ p[41]*plg[2][4] + t81)*c2tloc +(p[8]*plg[2][2] + p[42]*plg[2][4] + t82)*s2tloc);
861   }
862
863   /* TERDIURNAL */
864   if (flags->sw[14]) {
865     t[13] = f2 * ((p[39]*plg[3][3]+(p[93]*plg[3][4]+p[46]*plg[3][6])*cd14*flags->swc[5])* s3tloc +(p[40]*plg[3][3]+(p[94]*plg[3][4]+p[48]*plg[3][6])*cd14*flags->swc[5])* c3tloc);
866 }
867
868   /* magnetic activity based on daily ap */
869   if (flags->sw[9]==-1) {
870     ap = input->ap_a;
871     if (p[51]!=0) {
872       double exp1;
873       exp1 = exp(-10800.0*sqrt(p[51]*p[51])/(1.0+p[138]*(45.0-sqrt(input->g_lat*input->g_lat))));
874       if (exp1>0.99999)
875         exp1=0.99999;
876       if (p[24]<1.0E-4)
877         p[24]=1.0E-4;
878       apt[0]=sg0(exp1,p,ap->a);
879       /* apt[1]=sg2(exp1,p,ap->a);
880          apt[2]=sg0(exp2,p,ap->a);
881          apt[3]=sg2(exp2,p,ap->a);
882       */
883       if (flags->sw[9]) {
884         t[8] = apt[0]*(p[50]+p[96]*plg[0][2]+p[54]*plg[0][4]+ \
885      (p[125]*plg[0][1]+p[126]*plg[0][3]+p[127]*plg[0][5])*cd14*flags->swc[5]+ \
886      (p[128]*plg[1][1]+p[129]*plg[1][3]+p[130]*plg[1][5])*flags->swc[7]* \
887                  cos(hr*(tloc-p[131])));
888       }
889     }
890   } else {
891     double p44, p45;
892     apd=input->ap-4.0;
893     p44=p[43];
894     p45=p[44];
895     if (p44<0)
896       p44 = 1.0E-5;
897     apdf = apd + (p45-1.0)*(apd + (exp(-p44 * apd) - 1.0)/p44);
898     if (flags->sw[9]) {
899       t[8]=apdf*(p[32]+p[45]*plg[0][2]+p[34]*plg[0][4]+ \
900      (p[100]*plg[0][1]+p[101]*plg[0][3]+p[102]*plg[0][5])*cd14*flags->swc[5]+
901      (p[121]*plg[1][1]+p[122]*plg[1][3]+p[123]*plg[1][5])*flags->swc[7]*
902             cos(hr*(tloc-p[124])));
903     }
904   }
905
906   if ((flags->sw[10])&&(input->g_long>-1000.0)) {
907
908     /* longitudinal */
909     if (flags->sw[11]) {
910       t[10] = (1.0 + p[80]*dfa*flags->swc[1])* \
911      ((p[64]*plg[1][2]+p[65]*plg[1][4]+p[66]*plg[1][6]\
912       +p[103]*plg[1][1]+p[104]*plg[1][3]+p[105]*plg[1][5]\
913       +flags->swc[5]*(p[109]*plg[1][1]+p[110]*plg[1][3]+p[111]*plg[1][5])*cd14)* \
914           cos(dgtr*input->g_long) \
915       +(p[90]*plg[1][2]+p[91]*plg[1][4]+p[92]*plg[1][6]\
916       +p[106]*plg[1][1]+p[107]*plg[1][3]+p[108]*plg[1][5]\
917       +flags->swc[5]*(p[112]*plg[1][1]+p[113]*plg[1][3]+p[114]*plg[1][5])*cd14)* \
918       sin(dgtr*input->g_long));
919     }
920
921     /* ut and mixed ut, longitude */
922     if (flags->sw[12]){
923       t[11]=(1.0+p[95]*plg[0][1])*(1.0+p[81]*dfa*flags->swc[1])*\
924         (1.0+p[119]*plg[0][1]*flags->swc[5]*cd14)*\
925         ((p[68]*plg[0][1]+p[69]*plg[0][3]+p[70]*plg[0][5])*\
926         cos(sr*(input->sec-p[71])));
927       t[11]+=flags->swc[11]*\
928         (p[76]*plg[2][3]+p[77]*plg[2][5]+p[78]*plg[2][7])*\
929         cos(sr*(input->sec-p[79])+2.0*dgtr*input->g_long)*(1.0+p[137]*dfa*flags->swc[1]);
930     }
931
932     /* ut, longitude magnetic activity */
933     if (flags->sw[13]) {
934       if (flags->sw[9]==-1) {
935         if (p[51]) {
936           t[12]=apt[0]*flags->swc[11]*(1.+p[132]*plg[0][1])*\
937             ((p[52]*plg[1][2]+p[98]*plg[1][4]+p[67]*plg[1][6])*\
938              cos(dgtr*(input->g_long-p[97])))\
939             +apt[0]*flags->swc[11]*flags->swc[5]*\
940             (p[133]*plg[1][1]+p[134]*plg[1][3]+p[135]*plg[1][5])*\
941             cd14*cos(dgtr*(input->g_long-p[136])) \
942             +apt[0]*flags->swc[12]* \
943             (p[55]*plg[0][1]+p[56]*plg[0][3]+p[57]*plg[0][5])*\
944             cos(sr*(input->sec-p[58]));
945         }
946       } else {
947         t[12] = apdf*flags->swc[11]*(1.0+p[120]*plg[0][1])*\
948           ((p[60]*plg[1][2]+p[61]*plg[1][4]+p[62]*plg[1][6])*\
949           cos(dgtr*(input->g_long-p[63])))\
950           +apdf*flags->swc[11]*flags->swc[5]* \
951           (p[115]*plg[1][1]+p[116]*plg[1][3]+p[117]*plg[1][5])* \
952           cd14*cos(dgtr*(input->g_long-p[118])) \
953           + apdf*flags->swc[12]* \
954           (p[83]*plg[0][1]+p[84]*plg[0][3]+p[85]*plg[0][5])* \
955           cos(sr*(input->sec-p[75]));
956       }
957     }
958   }
959
960   /* parms not used: 82, 89, 99, 139-149 */
961   tinf = p[30];
962   for (i=0;i<14;i++)
963     tinf = tinf + fabs(flags->sw[i+1])*t[i];
964   return tinf;
965 }
966
967 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
968
969 double MSIS::glob7s(double *p, struct nrlmsise_input *input,
970                       struct nrlmsise_flags *flags)
971 {
972 /*    VERSION OF GLOBE FOR LOWER ATMOSPHERE 10/26/99
973  */
974   double pset=2.0;
975   double t[14] = {0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,};
976   double tt=0.0;
977   double cd32=0.0, cd18=0.0, cd14=0.0, cd39=0.0;
978   double p32=0.0, p18=0.0, p14=0.0, p39=0.0;
979   int i=0,j=0;
980   double dr=1.72142E-2;
981   double dgtr=1.74533E-2;
982   /* confirm parameter set */
983   if (p[99]==0)
984     p[99]=pset;
985   if (p[99]!=pset) {
986     cerr << "Wrong parameter set for glob7s" << endl;
987     return -1;
988   }
989   for (j=0;j<14;j++)
990     t[j]=0.0;
991   cd32 = cos(dr*(input->doy-p[31]));
992   cd18 = cos(2.0*dr*(input->doy-p[17]));
993   cd14 = cos(dr*(input->doy-p[13]));
994   cd39 = cos(2.0*dr*(input->doy-p[38]));
995   p32=p[31];
996   p18=p[17];
997   p14=p[13];
998   p39=p[38];
999
1000   /* F10.7 */
1001   t[0] = p[21]*dfa;
1002
1003   /* time independent */
1004   t[1]=p[1]*plg[0][2] + p[2]*plg[0][4] + p[22]*plg[0][6] + p[26]*plg[0][1] + p[14]*plg[0][3] + p[59]*plg[0][5];
1005
1006         /* SYMMETRICAL ANNUAL */
1007   t[2]=(p[18]+p[47]*plg[0][2]+p[29]*plg[0][4])*cd32;
1008
1009         /* SYMMETRICAL SEMIANNUAL */
1010   t[3]=(p[15]+p[16]*plg[0][2]+p[30]*plg[0][4])*cd18;
1011
1012         /* ASYMMETRICAL ANNUAL */
1013   t[4]=(p[9]*plg[0][1]+p[10]*plg[0][3]+p[20]*plg[0][5])*cd14;
1014
1015   /* ASYMMETRICAL SEMIANNUAL */
1016   t[5]=(p[37]*plg[0][1])*cd39;
1017
1018         /* DIURNAL */
1019   if (flags->sw[7]) {
1020     double t71, t72;
1021     t71 = p[11]*plg[1][2]*cd14*flags->swc[5];
1022     t72 = p[12]*plg[1][2]*cd14*flags->swc[5];
1023     t[6] = ((p[3]*plg[1][1] + p[4]*plg[1][3] + t71) * ctloc + (p[6]*plg[1][1] + p[7]*plg[1][3] + t72) * stloc) ;
1024   }
1025
1026   /* SEMIDIURNAL */
1027   if (flags->sw[8]) {
1028     double t81, t82;
1029     t81 = (p[23]*plg[2][3]+p[35]*plg[2][5])*cd14*flags->swc[5];
1030     t82 = (p[33]*plg[2][3]+p[36]*plg[2][5])*cd14*flags->swc[5];
1031     t[7] = ((p[5]*plg[2][2] + p[41]*plg[2][4] + t81) * c2tloc + (p[8]*plg[2][2] + p[42]*plg[2][4] + t82) * s2tloc);
1032   }
1033
1034   /* TERDIURNAL */
1035   if (flags->sw[14]) {
1036     t[13] = p[39] * plg[3][3] * s3tloc + p[40] * plg[3][3] * c3tloc;
1037   }
1038
1039   /* MAGNETIC ACTIVITY */
1040   if (flags->sw[9]) {
1041     if (flags->sw[9]==1)
1042       t[8] = apdf * (p[32] + p[45] * plg[0][2] * flags->swc[2]);
1043     if (flags->sw[9]==-1)
1044       t[8]=(p[50]*apt[0] + p[96]*plg[0][2] * apt[0]*flags->swc[2]);
1045   }
1046
1047   /* LONGITUDINAL */
1048   if (!((flags->sw[10]==0) || (flags->sw[11]==0) || (input->g_long<=-1000.0))) {
1049     t[10] = (1.0 + plg[0][1]*(p[80]*flags->swc[5]*cos(dr*(input->doy-p[81]))\
1050             +p[85]*flags->swc[6]*cos(2.0*dr*(input->doy-p[86])))\
1051       +p[83]*flags->swc[3]*cos(dr*(input->doy-p[84]))\
1052       +p[87]*flags->swc[4]*cos(2.0*dr*(input->doy-p[88])))\
1053       *((p[64]*plg[1][2]+p[65]*plg[1][4]+p[66]*plg[1][6]\
1054       +p[74]*plg[1][1]+p[75]*plg[1][3]+p[76]*plg[1][5]\
1055       )*cos(dgtr*input->g_long)\
1056       +(p[90]*plg[1][2]+p[91]*plg[1][4]+p[92]*plg[1][6]\
1057       +p[77]*plg[1][1]+p[78]*plg[1][3]+p[79]*plg[1][5]\
1058       )*sin(dgtr*input->g_long));
1059   }
1060   tt=0;
1061   for (i=0;i<14;i++)
1062     tt+=fabs(flags->sw[i+1])*t[i];
1063   return tt;
1064 }
1065
1066 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1067
1068 void MSIS::gtd7(struct nrlmsise_input *input, struct nrlmsise_flags *flags,
1069                   struct nrlmsise_output *output)
1070 {
1071   double xlat=0.0;
1072   double xmm=0.0;
1073   int mn3 = 5;
1074   double zn3[5]={32.5,20.0,15.0,10.0,0.0};
1075   int mn2 = 4;
1076   double zn2[4]={72.5,55.0,45.0,32.5};
1077   double altt=0.0;
1078   double zmix=62.5;
1079   double tmp=0.0;
1080   double dm28m=0.0;
1081   double tz=0.0;
1082   double dmc=0.0;
1083   double dmr=0.0;
1084   double dz28=0.0;
1085   struct nrlmsise_output soutput;
1086   int i;
1087
1088   for (int i=0; i<9; i++) soutput.d[i] = 0.0;
1089   for (int i=0; i<2; i++) soutput.t[i] = 0.0;
1090
1091   tselec(flags);
1092
1093   /* Latitude variation of gravity (none for sw[2]=0) */
1094   xlat=input->g_lat;
1095   if (flags->sw[2]==0)
1096     xlat=45.0;
1097   glatf(xlat, &gsurf, &re);
1098
1099   xmm = pdm[2][4];
1100
1101   /* THERMOSPHERE / MESOSPHERE (above zn2[0]) */
1102   if (input->alt>zn2[0])
1103     altt=input->alt;
1104   else
1105     altt=zn2[0];
1106
1107   tmp=input->alt;
1108   input->alt=altt;
1109   gts7(input, flags, &soutput);
1110   altt=input->alt;
1111   input->alt=tmp;
1112   if (flags->sw[0])   /* metric adjustment */
1113     dm28m=dm28*1.0E6;
1114   else
1115     dm28m=dm28;
1116   output->t[0]=soutput.t[0];
1117   output->t[1]=soutput.t[1];
1118   if (input->alt>=zn2[0]) {
1119     for (i=0;i<9;i++)
1120       output->d[i]=soutput.d[i];
1121     return;
1122   }
1123
1124 /*       LOWER MESOSPHERE/UPPER STRATOSPHERE (between zn3[0] and zn2[0])
1125  *         Temperature at nodes and gradients at end nodes
1126  *         Inverse temperature a linear function of spherical harmonics
1127  */
1128   meso_tgn2[0]=meso_tgn1[1];
1129   meso_tn2[0]=meso_tn1[4];
1130         meso_tn2[1]=pma[0][0]*pavgm[0]/(1.0-flags->sw[20]*glob7s(pma[0], input, flags));
1131         meso_tn2[2]=pma[1][0]*pavgm[1]/(1.0-flags->sw[20]*glob7s(pma[1], input, flags));
1132         meso_tn2[3]=pma[2][0]*pavgm[2]/(1.0-flags->sw[20]*flags->sw[22]*glob7s(pma[2], input, flags));
1133   meso_tgn2[1]=pavgm[8]*pma[9][0]*(1.0+flags->sw[20]*flags->sw[22]*glob7s(pma[9], input, flags))*meso_tn2[3]*meso_tn2[3]/(pow((pma[2][0]*pavgm[2]),2.0));
1134   meso_tn3[0]=meso_tn2[3];
1135
1136   if (input->alt<zn3[0]) {
1137 /*       LOWER STRATOSPHERE AND TROPOSPHERE (below zn3[0])
1138  *         Temperature at nodes and gradients at end nodes
1139  *         Inverse temperature a linear function of spherical harmonics
1140  */
1141     meso_tgn3[0]=meso_tgn2[1];
1142     meso_tn3[1]=pma[3][0]*pavgm[3]/(1.0-flags->sw[22]*glob7s(pma[3], input, flags));
1143     meso_tn3[2]=pma[4][0]*pavgm[4]/(1.0-flags->sw[22]*glob7s(pma[4], input, flags));
1144     meso_tn3[3]=pma[5][0]*pavgm[5]/(1.0-flags->sw[22]*glob7s(pma[5], input, flags));
1145     meso_tn3[4]=pma[6][0]*pavgm[6]/(1.0-flags->sw[22]*glob7s(pma[6], input, flags));
1146     meso_tgn3[1]=pma[7][0]*pavgm[7]*(1.0+flags->sw[22]*glob7s(pma[7], input, flags)) *meso_tn3[4]*meso_tn3[4]/(pow((pma[6][0]*pavgm[6]),2.0));
1147   }
1148
1149         /* LINEAR TRANSITION TO FULL MIXING BELOW zn2[0] */
1150
1151   dmc=0;
1152   if (input->alt>zmix)
1153     dmc = 1.0 - (zn2[0]-input->alt)/(zn2[0] - zmix);
1154   dz28=soutput.d[2];
1155
1156   /**** N2 density ****/
1157   dmr=soutput.d[2] / dm28m - 1.0;
1158   output->d[2]=densm(input->alt,dm28m,xmm, &tz, mn3, zn3, meso_tn3, meso_tgn3, mn2, zn2, meso_tn2, meso_tgn2);
1159   output->d[2]=output->d[2] * (1.0 + dmr*dmc);
1160
1161   /**** HE density ****/
1162   dmr = soutput.d[0] / (dz28 * pdm[0][1]) - 1.0;
1163   output->d[0] = output->d[2] * pdm[0][1] * (1.0 + dmr*dmc);
1164
1165   /**** O density ****/
1166   output->d[1] = 0;
1167   output->d[8] = 0;
1168
1169   /**** O2 density ****/
1170   dmr = soutput.d[3] / (dz28 * pdm[3][1]) - 1.0;
1171   output->d[3] = output->d[2] * pdm[3][1] * (1.0 + dmr*dmc);
1172
1173   /**** AR density ***/
1174   dmr = soutput.d[4] / (dz28 * pdm[4][1]) - 1.0;
1175   output->d[4] = output->d[2] * pdm[4][1] * (1.0 + dmr*dmc);
1176
1177   /**** Hydrogen density ****/
1178   output->d[6] = 0;
1179
1180   /**** Atomic nitrogen density ****/
1181   output->d[7] = 0;
1182
1183   /**** Total mass density */
1184   output->d[5] = 1.66E-24 * (4.0 * output->d[0] + 16.0 * output->d[1] +
1185                      28.0 * output->d[2] + 32.0 * output->d[3] + 40.0 * output->d[4]
1186                      + output->d[6] + 14.0 * output->d[7]);
1187
1188   if (flags->sw[0])
1189     output->d[5]=output->d[5]/1000;
1190
1191   /**** temperature at altitude ****/
1192   dd = densm(input->alt, 1.0, 0, &tz, mn3, zn3, meso_tn3, meso_tgn3,
1193                    mn2, zn2, meso_tn2, meso_tgn2);
1194   output->t[1]=tz;
1195
1196 }
1197
1198 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1199
1200 void MSIS::gtd7d(struct nrlmsise_input *input, struct nrlmsise_flags *flags,
1201                    struct nrlmsise_output *output)
1202 {
1203   gtd7(input, flags, output);
1204   output->d[5] = 1.66E-24 * (4.0 * output->d[0] + 16.0 * output->d[1] +
1205                    28.0 * output->d[2] + 32.0 * output->d[3] + 40.0 * output->d[4]
1206                    + output->d[6] + 14.0 * output->d[7] + 16.0 * output->d[8]);
1207   if (flags->sw[0])
1208     output->d[5]=output->d[5]/1000;
1209 }
1210
1211 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1212
1213 void MSIS::ghp7(struct nrlmsise_input *input, struct nrlmsise_flags *flags,
1214                   struct nrlmsise_output *output, double press)
1215 {
1216   double bm = 1.3806E-19;
1217   double rgas = 831.4;
1218   double test = 0.00043;
1219   double ltest = 12;
1220   double pl, p;
1221   double zi = 0.0;
1222   double z;
1223   double cl, cl2;
1224   double ca, cd;
1225   double xn, xm, diff;
1226   double g, sh;
1227   int l;
1228   pl = log10(press);
1229   if (pl >= -5.0) {
1230     if (pl>2.5)
1231       zi = 18.06 * (3.00 - pl);
1232     else if ((pl>0.075) && (pl<=2.5))
1233       zi = 14.98 * (3.08 - pl);
1234     else if ((pl>-1) && (pl<=0.075))
1235       zi = 17.80 * (2.72 - pl);
1236     else if ((pl>-2) && (pl<=-1))
1237       zi = 14.28 * (3.64 - pl);
1238     else if ((pl>-4) && (pl<=-2))
1239       zi = 12.72 * (4.32 -pl);
1240     else if (pl<=-4)
1241       zi = 25.3 * (0.11 - pl);
1242     cl = input->g_lat/90.0;
1243     cl2 = cl*cl;
1244     if (input->doy<182)
1245       cd = (1.0 - (double) input->doy) / 91.25;
1246     else
1247       cd = ((double) input->doy) / 91.25 - 3.0;
1248     ca = 0;
1249     if ((pl > -1.11) && (pl<=-0.23))
1250       ca = 1.0;
1251     if (pl > -0.23)
1252       ca = (2.79 - pl) / (2.79 + 0.23);
1253     if ((pl <= -1.11) && (pl>-3))
1254       ca = (-2.93 - pl)/(-2.93 + 1.11);
1255     z = zi - 4.87 * cl * cd * ca - 1.64 * cl2 * ca + 0.31 * ca * cl;
1256   } else
1257     z = 22.0 * pow((pl + 4.0),2.0) + 110.0;
1258
1259   /* iteration  loop */
1260   l = 0;
1261   do {
1262     l++;
1263     input->alt = z;
1264     gtd7(input, flags, output);
1265     z = input->alt;
1266     xn = output->d[0] + output->d[1] + output->d[2] + output->d[3] + output->d[4] + output->d[6] + output->d[7];
1267     p = bm * xn * output->t[1];
1268     if (flags->sw[0])
1269       p = p*1.0E-6;
1270     diff = pl - log10(p);
1271     if (sqrt(diff*diff)<test)
1272       return;
1273     if (l==ltest) {
1274       cerr << "ERROR: ghp7 not converging for press " << press << ", diff " << diff << endl;
1275       return;
1276     }
1277     xm = output->d[5] / xn / 1.66E-24;
1278     if (flags->sw[0])
1279       xm = xm * 1.0E3;
1280     g = gsurf / (pow((1.0 + z/re),2.0));
1281     sh = rgas * output->t[1] / (xm * g);
1282
1283     /* new altitude estimate using scale height */
1284     if (l <  6)
1285       z = z - sh * diff * 2.302;
1286     else
1287       z = z - sh * diff;
1288   } while (1==1);
1289 }
1290
1291 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1292
1293 void MSIS::gts7(struct nrlmsise_input *input, struct nrlmsise_flags *flags,
1294                   struct nrlmsise_output *output)
1295 {
1296 /*     Thermospheric portion of NRLMSISE-00
1297  *     See GTD7 for more extensive comments
1298  *     alt > 72.5 km!
1299  */
1300   double za=0.0;
1301   int i, j;
1302   double ddum=0.0, z=0.0;
1303   double zn1[5] = {120.0, 110.0, 100.0, 90.0, 72.5};
1304   double tinf=0.0;
1305   int mn1 = 5;
1306   double g0=0.0;
1307   double tlb=0.0;
1308   double s=0.0, z0=0.0, t0=0.0, tr12=0.0;
1309   double db01=0.0, db04=0.0, db14=0.0, db16=0.0, db28=0.0, db32=0.0, db40=0.0, db48=0.0;
1310   double zh28=0.0, zh04=0.0, zh16=0.0, zh32=0.0, zh40=0.0, zh01=0.0, zh14=0.0;
1311   double zhm28=0.0, zhm04=0.0, zhm16=0.0, zhm32=0.0, zhm40=0.0, zhm01=0.0, zhm14=0.0;
1312   double xmd=0.0;
1313   double b28=0.0, b04=0.0, b16=0.0, b32=0.0, b40=0.0, b01=0.0, b14=0.0;
1314   double tz=0.0;
1315   double g28=0.0, g4=0.0, g16=0.0, g32=0.0, g40=0.0, g1=0.0, g14=0.0;
1316   double zhf=0.0, xmm=0.0;
1317   double zc04=0.0, zc16=0.0, zc32=0.0, zc40=0.0, zc01=0.0, zc14=0.0;
1318   double hc04=0.0, hc16=0.0, hc32=0.0, hc40=0.0, hc01=0.0, hc14=0.0;
1319   double hcc16=0.0, hcc32=0.0, hcc01=0.0, hcc14=0.0;
1320   double zcc16=0.0, zcc32=0.0, zcc01=0.0, zcc14=0.0;
1321   double rc16=0.0, rc32=0.0, rc01=0.0, rc14=0.0;
1322   double rl=0.0;
1323   double g16h=0.0, db16h=0.0, tho=0.0, zsht=0.0, zmho=0.0, zsho=0.0;
1324   double dgtr=1.74533E-2;
1325   double dr=1.72142E-2;
1326   double alpha[9]={-0.38, 0.0, 0.0, 0.0, 0.17, 0.0, -0.38, 0.0, 0.0};
1327   double altl[8]={200.0, 300.0, 160.0, 250.0, 240.0, 450.0, 320.0, 450.0};
1328   double dd=0.0;
1329   double hc216=0.0, hcc232=0.0;
1330   za = pdl[1][15];
1331   zn1[0] = za;
1332   for (j=0;j<9;j++)
1333     output->d[j]=0;
1334
1335   /* TINF VARIATIONS NOT IMPORTANT BELOW ZA OR ZN1(1) */
1336   if (input->alt>zn1[0])
1337     tinf = ptm[0]*pt[0] * \
1338       (1.0+flags->sw[16]*globe7(pt,input,flags));
1339   else
1340     tinf = ptm[0]*pt[0];
1341   output->t[0]=tinf;
1342
1343   /*  GRADIENT VARIATIONS NOT IMPORTANT BELOW ZN1(5) */
1344   if (input->alt>zn1[4])
1345     g0 = ptm[3]*ps[0] * \
1346       (1.0+flags->sw[19]*globe7(ps,input,flags));
1347   else
1348     g0 = ptm[3]*ps[0];
1349   tlb = ptm[1] * (1.0 + flags->sw[17]*globe7(pd[3],input,flags))*pd[3][0];
1350   s = g0 / (tinf - tlb);
1351
1352 /*      Lower thermosphere temp variations not significant for
1353  *       density above 300 km */
1354   if (input->alt<300.0) {
1355     meso_tn1[1]=ptm[6]*ptl[0][0]/(1.0-flags->sw[18]*glob7s(ptl[0], input, flags));
1356     meso_tn1[2]=ptm[2]*ptl[1][0]/(1.0-flags->sw[18]*glob7s(ptl[1], input, flags));
1357     meso_tn1[3]=ptm[7]*ptl[2][0]/(1.0-flags->sw[18]*glob7s(ptl[2], input, flags));
1358     meso_tn1[4]=ptm[4]*ptl[3][0]/(1.0-flags->sw[18]*flags->sw[20]*glob7s(ptl[3], input, flags));
1359     meso_tgn1[1]=ptm[8]*pma[8][0]*(1.0+flags->sw[18]*flags->sw[20]*glob7s(pma[8], input, flags))*meso_tn1[4]*meso_tn1[4]/(pow((ptm[4]*ptl[3][0]),2.0));
1360   } else {
1361     meso_tn1[1]=ptm[6]*ptl[0][0];
1362     meso_tn1[2]=ptm[2]*ptl[1][0];
1363     meso_tn1[3]=ptm[7]*ptl[2][0];
1364     meso_tn1[4]=ptm[4]*ptl[3][0];
1365     meso_tgn1[1]=ptm[8]*pma[8][0]*meso_tn1[4]*meso_tn1[4]/(pow((ptm[4]*ptl[3][0]),2.0));
1366   }
1367
1368   z0 = zn1[3];
1369   t0 = meso_tn1[3];
1370   tr12 = 1.0;
1371
1372   /* N2 variation factor at Zlb */
1373   g28=flags->sw[21]*globe7(pd[2], input, flags);
1374
1375   /* VARIATION OF TURBOPAUSE HEIGHT */
1376   zhf=pdl[1][24]*(1.0+flags->sw[5]*pdl[0][24]*sin(dgtr*input->g_lat)*cos(dr*(input->doy-pt[13])));
1377   output->t[0]=tinf;
1378   xmm = pdm[2][4];
1379   z = input->alt;
1380
1381
1382         /**** N2 DENSITY ****/
1383
1384   /* Diffusive density at Zlb */
1385   db28 = pdm[2][0]*exp(g28)*pd[2][0];
1386   /* Diffusive density at Alt */
1387   output->d[2]=densu(z,db28,tinf,tlb,28.0,alpha[2],&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1388   dd=output->d[2];
1389   /* Turbopause */
1390   zh28=pdm[2][2]*zhf;
1391   zhm28=pdm[2][3]*pdl[1][5];
1392   xmd=28.0-xmm;
1393   /* Mixed density at Zlb */
1394   b28=densu(zh28,db28,tinf,tlb,xmd,(alpha[2]-1.0),&tz,ptm[5],s,mn1, zn1,meso_tn1,meso_tgn1);
1395   if ((flags->sw[15])&&(z<=altl[2])) {
1396     /*  Mixed density at Alt */
1397     dm28=densu(z,b28,tinf,tlb,xmm,alpha[2],&tz,ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1398     /*  Net density at Alt */
1399     output->d[2]=dnet(output->d[2],dm28,zhm28,xmm,28.0);
1400   }
1401
1402
1403         /**** HE DENSITY ****/
1404
1405   /*   Density variation factor at Zlb */
1406   g4 = flags->sw[21]*globe7(pd[0], input, flags);
1407   /*  Diffusive density at Zlb */
1408   db04 = pdm[0][0]*exp(g4)*pd[0][0];
1409         /*  Diffusive density at Alt */
1410   output->d[0]=densu(z,db04,tinf,tlb, 4.,alpha[0],&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1411   dd=output->d[0];
1412   if ((flags->sw[15]) && (z<altl[0])) {
1413     /*  Turbopause */
1414     zh04=pdm[0][2];
1415     /*  Mixed density at Zlb */
1416     b04=densu(zh04,db04,tinf,tlb,4.-xmm,alpha[0]-1.,&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1417     /*  Mixed density at Alt */
1418     dm04=densu(z,b04,tinf,tlb,xmm,0.,&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1419     zhm04=zhm28;
1420     /*  Net density at Alt */
1421     output->d[0]=dnet(output->d[0],dm04,zhm04,xmm,4.);
1422     /*  Correction to specified mixing ratio at ground */
1423     rl=log(b28*pdm[0][1]/b04);
1424     zc04=pdm[0][4]*pdl[1][0];
1425     hc04=pdm[0][5]*pdl[1][1];
1426     /*  Net density corrected at Alt */
1427     output->d[0]=output->d[0]*ccor(z,rl,hc04,zc04);
1428   }
1429
1430
1431         /**** O DENSITY ****/
1432
1433   /*  Density variation factor at Zlb */
1434   g16= flags->sw[21]*globe7(pd[1],input,flags);
1435   /*  Diffusive density at Zlb */
1436   db16 =  pdm[1][0]*exp(g16)*pd[1][0];
1437         /*   Diffusive density at Alt */
1438   output->d[1]=densu(z,db16,tinf,tlb, 16.,alpha[1],&output->t[1],ptm[5],s,mn1, zn1,meso_tn1,meso_tgn1);
1439   dd=output->d[1];
1440   if ((flags->sw[15]) && (z<=altl[1])) {
1441     /*   Turbopause */
1442     zh16=pdm[1][2];
1443     /*  Mixed density at Zlb */
1444     b16=densu(zh16,db16,tinf,tlb,16.0-xmm,(alpha[1]-1.0), &output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1445     /*  Mixed density at Alt */
1446     dm16=densu(z,b16,tinf,tlb,xmm,0.,&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1447     zhm16=zhm28;
1448     /*  Net density at Alt */
1449     output->d[1]=dnet(output->d[1],dm16,zhm16,xmm,16.);
1450     rl=pdm[1][1]*pdl[1][16]*(1.0+flags->sw[1]*pdl[0][23]*(input->f107A-150.0));
1451     hc16=pdm[1][5]*pdl[1][3];
1452     zc16=pdm[1][4]*pdl[1][2];
1453     hc216=pdm[1][5]*pdl[1][4];
1454     output->d[1]=output->d[1]*ccor2(z,rl,hc16,zc16,hc216);
1455     /*   Chemistry correction */
1456     hcc16=pdm[1][7]*pdl[1][13];
1457     zcc16=pdm[1][6]*pdl[1][12];
1458     rc16=pdm[1][3]*pdl[1][14];
1459     /*  Net density corrected at Alt */
1460     output->d[1]=output->d[1]*ccor(z,rc16,hcc16,zcc16);
1461   }
1462
1463
1464         /**** O2 DENSITY ****/
1465
1466         /*   Density variation factor at Zlb */
1467   g32= flags->sw[21]*globe7(pd[4], input, flags);
1468         /*  Diffusive density at Zlb */
1469   db32 = pdm[3][0]*exp(g32)*pd[4][0];
1470         /*   Diffusive density at Alt */
1471   output->d[3]=densu(z,db32,tinf,tlb, 32.,alpha[3],&output->t[1],ptm[5],s,mn1, zn1,meso_tn1,meso_tgn1);
1472   dd=output->d[3];
1473   if (flags->sw[15]) {
1474     if (z<=altl[3]) {
1475       /*   Turbopause */
1476       zh32=pdm[3][2];
1477       /*  Mixed density at Zlb */
1478       b32=densu(zh32,db32,tinf,tlb,32.-xmm,alpha[3]-1., &output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1479       /*  Mixed density at Alt */
1480       dm32=densu(z,b32,tinf,tlb,xmm,0.,&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1481       zhm32=zhm28;
1482       /*  Net density at Alt */
1483       output->d[3]=dnet(output->d[3],dm32,zhm32,xmm,32.);
1484       /*   Correction to specified mixing ratio at ground */
1485       rl=log(b28*pdm[3][1]/b32);
1486       hc32=pdm[3][5]*pdl[1][7];
1487       zc32=pdm[3][4]*pdl[1][6];
1488       output->d[3]=output->d[3]*ccor(z,rl,hc32,zc32);
1489     }
1490     /*  Correction for general departure from diffusive equilibrium above Zlb */
1491     hcc32=pdm[3][7]*pdl[1][22];
1492     hcc232=pdm[3][7]*pdl[0][22];
1493     zcc32=pdm[3][6]*pdl[1][21];
1494     rc32=pdm[3][3]*pdl[1][23]*(1.+flags->sw[1]*pdl[0][23]*(input->f107A-150.));
1495     /*  Net density corrected at Alt */
1496     output->d[3]=output->d[3]*ccor2(z,rc32,hcc32,zcc32,hcc232);
1497   }
1498
1499
1500         /**** AR DENSITY ****/
1501
1502         /*   Density variation factor at Zlb */
1503   g40= flags->sw[21]*globe7(pd[5],input,flags);
1504         /*  Diffusive density at Zlb */
1505   db40 = pdm[4][0]*exp(g40)*pd[5][0];
1506   /*   Diffusive density at Alt */
1507   output->d[4]=densu(z,db40,tinf,tlb, 40.,alpha[4],&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1508   dd=output->d[4];
1509   if ((flags->sw[15]) && (z<=altl[4])) {
1510     /*   Turbopause */
1511     zh40=pdm[4][2];
1512     /*  Mixed density at Zlb */
1513     b40=densu(zh40,db40,tinf,tlb,40.-xmm,alpha[4]-1.,&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1514     /*  Mixed density at Alt */
1515     dm40=densu(z,b40,tinf,tlb,xmm,0.,&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1516     zhm40=zhm28;
1517     /*  Net density at Alt */
1518     output->d[4]=dnet(output->d[4],dm40,zhm40,xmm,40.);
1519     /*   Correction to specified mixing ratio at ground */
1520     rl=log(b28*pdm[4][1]/b40);
1521     hc40=pdm[4][5]*pdl[1][9];
1522     zc40=pdm[4][4]*pdl[1][8];
1523     /*  Net density corrected at Alt */
1524     output->d[4]=output->d[4]*ccor(z,rl,hc40,zc40);
1525     }
1526
1527
1528         /**** HYDROGEN DENSITY ****/
1529
1530         /*   Density variation factor at Zlb */
1531   g1 = flags->sw[21]*globe7(pd[6], input, flags);
1532         /*  Diffusive density at Zlb */
1533   db01 = pdm[5][0]*exp(g1)*pd[6][0];
1534         /*   Diffusive density at Alt */
1535   output->d[6]=densu(z,db01,tinf,tlb,1.,alpha[6],&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1536   dd=output->d[6];
1537   if ((flags->sw[15]) && (z<=altl[6])) {
1538     /*   Turbopause */
1539     zh01=pdm[5][2];
1540     /*  Mixed density at Zlb */
1541     b01=densu(zh01,db01,tinf,tlb,1.-xmm,alpha[6]-1., &output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1542     /*  Mixed density at Alt */
1543     dm01=densu(z,b01,tinf,tlb,xmm,0.,&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1544     zhm01=zhm28;
1545     /*  Net density at Alt */
1546     output->d[6]=dnet(output->d[6],dm01,zhm01,xmm,1.);
1547     /*   Correction to specified mixing ratio at ground */
1548     rl=log(b28*pdm[5][1]*sqrt(pdl[1][17]*pdl[1][17])/b01);
1549     hc01=pdm[5][5]*pdl[1][11];
1550     zc01=pdm[5][4]*pdl[1][10];
1551     output->d[6]=output->d[6]*ccor(z,rl,hc01,zc01);
1552     /*   Chemistry correction */
1553     hcc01=pdm[5][7]*pdl[1][19];
1554     zcc01=pdm[5][6]*pdl[1][18];
1555     rc01=pdm[5][3]*pdl[1][20];
1556     /*  Net density corrected at Alt */
1557     output->d[6]=output->d[6]*ccor(z,rc01,hcc01,zcc01);
1558 }
1559
1560
1561         /**** ATOMIC NITROGEN DENSITY ****/
1562
1563   /*   Density variation factor at Zlb */
1564   g14 = flags->sw[21]*globe7(pd[7],input,flags);
1565         /*  Diffusive density at Zlb */
1566   db14 = pdm[6][0]*exp(g14)*pd[7][0];
1567         /*   Diffusive density at Alt */
1568   output->d[7]=densu(z,db14,tinf,tlb,14.,alpha[7],&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1569   dd=output->d[7];
1570   if ((flags->sw[15]) && (z<=altl[7])) {
1571     /*   Turbopause */
1572     zh14=pdm[6][2];
1573     /*  Mixed density at Zlb */
1574     b14=densu(zh14,db14,tinf,tlb,14.-xmm,alpha[7]-1., &output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1575     /*  Mixed density at Alt */
1576     dm14=densu(z,b14,tinf,tlb,xmm,0.,&output->t[1],ptm[5],s,mn1,zn1,meso_tn1,meso_tgn1);
1577     zhm14=zhm28;
1578     /*  Net density at Alt */
1579     output->d[7]=dnet(output->d[7],dm14,zhm14,xmm,14.);
1580     /*   Correction to specified mixing ratio at ground */
1581     rl=log(b28*pdm[6][1]*sqrt(pdl[0][2]*pdl[0][2])/b14);
1582     hc14=pdm[6][5]*pdl[0][1];
1583     zc14=pdm[6][4]*pdl[0][0];
1584     output->d[7]=output->d[7]*ccor(z,rl,hc14,zc14);
1585     /*   Chemistry correction */
1586     hcc14=pdm[6][7]*pdl[0][4];
1587     zcc14=pdm[6][6]*pdl[0][3];
1588     rc14=pdm[6][3]*pdl[0][5];
1589     /*  Net density corrected at Alt */
1590     output->d[7]=output->d[7]*ccor(z,rc14,hcc14,zcc14);
1591   }
1592
1593
1594         /**** Anomalous OXYGEN DENSITY ****/
1595
1596   g16h = flags->sw[21]*globe7(pd[8],input,flags);
1597   db16h = pdm[7][0]*exp(g16h)*pd[8][0];
1598   tho = pdm[7][9]*pdl[0][6];
1599   dd=densu(z,db16h,tho,tho,16.,alpha[8],&output->t[1],ptm[5],s,mn1, zn1,meso_tn1,meso_tgn1);
1600   zsht=pdm[7][5];
1601   zmho=pdm[7][4];
1602   zsho=scalh(zmho,16.0,tho);
1603   output->d[8]=dd*exp(-zsht/zsho*(exp(-(z-zmho)/zsht)-1.));
1604
1605
1606   /* total mass density */
1607   output->d[5] = 1.66E-24*(4.0*output->d[0]+16.0*output->d[1]+28.0*output->d[2]+32.0*output->d[3]+40.0*output->d[4]+ output->d[6]+14.0*output->d[7]);
1608   db48=1.66E-24*(4.0*db04+16.0*db16+28.0*db28+32.0*db32+40.0*db40+db01+14.0*db14);
1609
1610
1611
1612   /* temperature */
1613   z = sqrt(input->alt*input->alt);
1614   ddum = densu(z,1.0, tinf, tlb, 0.0, 0.0, &output->t[1], ptm[5], s, mn1, zn1, meso_tn1, meso_tgn1);
1615   if (flags->sw[0]) {
1616     for(i=0;i<9;i++)
1617       output->d[i]=output->d[i]*1.0E6;
1618     output->d[5]=output->d[5]/1000;
1619   }
1620 }
1621
1622
1623 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1624 //    The bitmasked value choices are as follows:
1625 //    unset: In this case (the default) JSBSim would only print
1626 //       out the normally expected messages, essentially echoing
1627 //       the config files as they are read. If the environment
1628 //       variable is not set, debug_lvl is set to 1 internally
1629 //    0: This requests JSBSim not to output any messages
1630 //       whatsoever.
1631 //    1: This value explicity requests the normal JSBSim
1632 //       startup messages
1633 //    2: This value asks for a message to be printed out when
1634 //       a class is instantiated
1635 //    4: When this value is set, a message is displayed when a
1636 //       FGModel object executes its Run() method
1637 //    8: When this value is set, various runtime state variables
1638 //       are printed out periodically
1639 //    16: When set various parameters are sanity checked and
1640 //       a message is printed out when they go out of bounds
1641
1642 void MSIS::Debug(int from)
1643 {
1644   if (debug_lvl <= 0) return;
1645
1646   if (debug_lvl & 1) { // Standard console startup message output
1647     if (from == 0) { // Constructor
1648     }
1649   }
1650   if (debug_lvl & 2 ) { // Instantiation/Destruction notification
1651     if (from == 0) cout << "Instantiated: MSIS" << endl;
1652     if (from == 1) cout << "Destroyed:    MSIS" << endl;
1653   }
1654   if (debug_lvl & 4 ) { // Run() method entry print for FGModel-derived objects
1655   }
1656   if (debug_lvl & 8 ) { // Runtime state variables
1657   }
1658   if (debug_lvl & 16) { // Sanity checking
1659   }
1660   if (debug_lvl & 32) { // Turbulence
1661   }
1662   if (debug_lvl & 64) {
1663     if (from == 0) { // Constructor
1664       cout << IdSrc << endl;
1665       cout << IdHdr << endl;
1666     }
1667   }
1668 }
1669
1670
1671
1672 } // namespace JSBSim
1673