-
Notifications
You must be signed in to change notification settings - Fork 0
/
dfn.f90
504 lines (447 loc) · 22.3 KB
/
dfn.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
!!!
!!! Authors: Douglas P. Drob and John Emmert, NRL 7632
!!!
!!! =================================================
!!! MSIS density and composition vertical profiles
!!! =================================================
module dfn
use constants, only : nl, nsplO1, nsplNO
type dnparm
sequence
real(8) :: lnPhiF ! (Except O, H) Natural log of mixing ratio at zetaF (70 km), before chemical and dynamical corrections are applied (ln m^-3) (global term only)
real(8) :: lndref ! Natural log of number density at reference height
real(8) :: zetaM ! "Turbopause Height": Height of midpoint of effective mass transition (km)
real(8) :: HML ! Scale height of lower portion of effective mass profile (km)
real(8) :: HMU ! Scale height of upper portion of effective mass profile (km)
real(8) :: C ! Chapman term coefficient
real(8) :: zetaC ! Chapman term reference height (km)
real(8) :: HC ! Chapman term scale height (km)
real(8) :: R ! Chemical/dynamical term coefficient
real(8) :: zetaR ! Chemical/dynamical term reference height (km)
real(8) :: HR ! Chemical/dynamical term scale height (km)
real(8) :: cf(0:nsplO1+1)! Merged spline coefficients (for chemistry-dominated region of O1, NO, and (eventually), H, N)
real(8) :: zref ! Reference height for hydrostatic integral and ideal gas terms
real(8) :: Mi(0:4) ! Effective mass at nodes of piecewise mass profile (derived from zetaM, HML, HMU)
real(8) :: zetaMi(0:4) ! Height of nodes of piecewise mass profile (derived from zetaM, HML, HMU)
real(8) :: aMi(0:4) = 0d0 ! Slopes of piecewise mass profile segments (derived from zetaM, HML, HMU)
real(8) :: WMi(0:4) = 0d0 ! 2nd indefinite integral of 1/T at mass profile nodes
real(8) :: XMi(0:4) = 0d0 ! Cumulative adjustment to M/T integral due to changing effective mass
real(8) :: Izref ! Indefinite hydrostatic integral at reference height
real(8) :: Tref ! Temperature at reference height (for ideal gas law term)
real(8) :: zmin ! Minimum height of profile (missing values below)
real(8) :: zhyd ! Hydrostatic terms needed above this height
integer(8) :: massid
end type dnparm
contains
!==========================================================
! A subroutine to evaluate all composition parameters that
! are independent of altitude, but not location, time, solar
! flux, and/or mass number. These only need to be updated
! if gf (and implicitly tpro) have changed.
!==========================================================
subroutine dfnparm(massid,gf,tpro,dpro)
use constants
use tfn, only : tnparm
use msis,only: eta,sfluxmod,geomag,utdep,TN,PR,N2,O2,O1,HE,H1,AR,N1,OA,NO
implicit none
integer,intent(in) :: massid
real(8),intent(in) :: gf(0:maxnbf-1)
type(tnparm),intent(in) :: tpro
type(dnparm),intent(out) :: dpro
integer :: izf, i, i1, iz
real(8) :: Cterm
real(8) :: bc(2)
real(8) :: hbetaL,hbetaU
real(8) :: delM, delz
real(8) :: Wi !2nd indefinite integral at a piecewise mass profile node
real(8) :: Si(-5:0,2:6) !Array of b-spline values at a mass profile node
real(8) :: Mzref !Effective mass at reference altitude
real(8), external :: dilog
! ==================================================================
! Compute the species dependent profile parameters
! ==================================================================
dpro%massid = massid
select case(massid)
! Molecular Nitrogen ----------------------
case(2)
! Reference mixing ratio
dpro%lnPhiF = log(volmix(massid))
dpro%lndref = tpro%lndtotF + dpro%lnPhiF
dpro%zref = zetaF
dpro%zmin = -1d0
dpro%zhyd = zetaF
! Effective mass
dpro%zetaM = dot_product(N2%beta(0:mbf,1),gf(0:mbf))
dpro%HML = dot_product(N2%beta(0:mbf,2),gf(0:mbf))
dpro%HMU = dot_product(N2%beta(0:mbf,3),gf(0:mbf))
! Photochemical correction
dpro%R = dot_product(N2%beta(0:mbf,7),gf(0:mbf))
dpro%zetaR = dot_product(N2%beta(0:mbf,8),gf(0:mbf))
dpro%HR = dot_product(N2%beta(0:mbf,9),gf(0:mbf))
! Molecular Oxygen ------------------------
case(3)
! Reference mixing ratio
dpro%lnPhiF = log(volmix(massid))
dpro%lndref = tpro%lndtotF + dpro%lnPhiF
dpro%zref = zetaF
dpro%zmin = -1d0
dpro%zhyd = zetaF
! Effective mass
dpro%zetaM = dot_product(O2%beta(0:mbf,1),gf(0:mbf))
dpro%HML = dot_product(O2%beta(0:mbf,2),gf(0:mbf))
dpro%HMU = dot_product(O2%beta(0:mbf,3),gf(0:mbf))
! Photochemical correction
dpro%R = dot_product(O2%beta(0:mbf,7),gf(0:mbf))
dpro%R = dpro%R + geomag(O2%beta(cmag:cmag+nmag-1,7),gf(cmag:cmag+12),gf(cmag+13:cmag+26))
dpro%zetaR = dot_product(O2%beta(0:mbf,8),gf(0:mbf))
dpro%HR = dot_product(O2%beta(0:mbf,9),gf(0:mbf))
! Atomic Oxygen --------------------------
case(4)
! Reference number density
dpro%lnPhiF = 0d0
dpro%lndref = dot_product(O1%beta(0:mbf,0),gf(0:mbf))
dpro%zref = zetarefO1
dpro%zmin = nodesO1(3)
dpro%zhyd = zetarefO1
! Effective mass
dpro%zetaM = dot_product(O1%beta(0:mbf,1),gf(0:mbf))
dpro%HML = dot_product(O1%beta(0:mbf,2),gf(0:mbf))
dpro%HMU = dot_product(O1%beta(0:mbf,3),gf(0:mbf))
! Chapman correction
dpro%C = dot_product(O1%beta(0:mbf,4),gf(0:mbf))
dpro%zetaC = dot_product(O1%beta(0:mbf,5),gf(0:mbf))
dpro%HC = dot_product(O1%beta(0:mbf,6),gf(0:mbf))
! Dynamical correction
dpro%R = dot_product(O1%beta(0:mbf,7),gf(0:mbf))
dpro%R = dpro%R + sfluxmod(7,gf,O1,0d0)
dpro%R = dpro%R + geomag(O1%beta(cmag:cmag+nmag-1,7),gf(cmag:cmag+12),gf(cmag+13:cmag+26))
dpro%R = dpro%R + utdep(O1%beta(cut:cut+nut-1,7),gf(cut:cut+8))
dpro%zetaR = dot_product(O1%beta(0:mbf,8),gf(0:mbf))
dpro%HR = dot_product(O1%beta(0:mbf,9),gf(0:mbf))
! Unconstrained splines
do izf = 0,nsplO1-1
dpro%cf(izf) = dot_product(O1%beta(0:mbf,izf+10),gf(0:mbf))
enddo
!Constrained splines calculated after case statement
! Helium ----------------------
case(5)
! Reference mixing ratio
dpro%lnPhiF = log(volmix(massid))
dpro%lndref = tpro%lndtotF + dpro%lnPhiF
dpro%zref = zetaF
dpro%zmin = -1d0
dpro%zhyd = zetaF
! Effective mass
dpro%zetaM = dot_product(HE%beta(0:mbf,1),gf(0:mbf))
dpro%HML = dot_product(HE%beta(0:mbf,2),gf(0:mbf))
dpro%HMU = dot_product(HE%beta(0:mbf,3),gf(0:mbf))
! Dynamical correction
dpro%R = dot_product(HE%beta(0:mbf,7),gf(0:mbf))
dpro%R = dpro%R + sfluxmod(7,gf,HE,1d0)
dpro%R = dpro%R + geomag(HE%beta(cmag:cmag+nmag-1,7),gf(cmag:cmag+12),gf(cmag+13:cmag+26))
dpro%R = dpro%R + utdep(HE%beta(cut:cut+nut-1,7),gf(cut:cut+8))
dpro%zetaR = dot_product(HE%beta(0:mbf,8),gf(0:mbf))
dpro%HR = dot_product(HE%beta(0:mbf,9),gf(0:mbf))
! Atomic Hydrogen ----------------------
case(6)
! Reference number density
dpro%lnPhiF = 0d0
dpro%lndref = dot_product(H1%beta(0:mbf,0),gf(0:mbf))
dpro%zref = zetaA
dpro%zmin = 75d0
dpro%zhyd = zetaF
! Effective mass
dpro%zetaM = dot_product(H1%beta(0:mbf,1),gf(0:mbf))
dpro%HML = dot_product(H1%beta(0:mbf,2),gf(0:mbf))
dpro%HMU = dot_product(H1%beta(0:mbf,3),gf(0:mbf))
! Chapman correction
dpro%C = dot_product(H1%beta(0:mbf,4),gf(0:mbf))
dpro%zetaC = dot_product(H1%beta(0:mbf,5),gf(0:mbf))
dpro%HC = dot_product(H1%beta(0:mbf,6),gf(0:mbf))
! Dynamical correction
dpro%R = dot_product(H1%beta(0:mbf,7),gf(0:mbf))
dpro%R = dpro%R + sfluxmod(7,gf,H1,0d0)
dpro%R = dpro%R + geomag(H1%beta(cmag:cmag+nmag-1,7),gf(cmag:cmag+12),gf(cmag+13:cmag+26))
dpro%R = dpro%R + utdep(H1%beta(cut:cut+nut-1,7),gf(cut:cut+8))
dpro%zetaR = dot_product(H1%beta(0:mbf,8),gf(0:mbf))
dpro%HR = dot_product(H1%beta(0:mbf,9),gf(0:mbf))
! Argon ----------------------
case(7)
! Reference mixing ratio
dpro%lnPhiF = log(volmix(massid))
dpro%lndref = tpro%lndtotF + dpro%lnPhiF
dpro%zref = zetaF
dpro%zmin = -1d0
dpro%zhyd = zetaF
! Effective mass
dpro%zetaM = dot_product(AR%beta(0:mbf,1),gf(0:mbf))
dpro%HML = dot_product(AR%beta(0:mbf,2),gf(0:mbf))
dpro%HMU = dot_product(AR%beta(0:mbf,3),gf(0:mbf))
! Dynamical correction
dpro%R = dot_product(AR%beta(0:mbf,7),gf(0:mbf))
dpro%R = dpro%R + geomag(AR%beta(cmag:cmag+nmag-1,7),gf(cmag:cmag+12),gf(cmag+13:cmag+26))
dpro%R = dpro%R + utdep(AR%beta(cut:cut+nut-1,7),gf(cut:cut+8))
dpro%zetaR = dot_product(AR%beta(0:mbf,8),gf(0:mbf))
dpro%HR = dot_product(AR%beta(0:mbf,9),gf(0:mbf))
! ! Atomic Nitrogen ----------------------
case(8)
! Reference number density
dpro%lnPhiF = 0d0
dpro%lndref = dot_product(N1%beta(0:mbf,0),gf(0:mbf))
dpro%lndref = dpro%lndref + sfluxmod(0,gf,N1,0d0)
dpro%lndref = dpro%lndref + geomag(N1%beta(cmag:cmag+nmag-1,0),gf(cmag:cmag+12),gf(cmag+13:cmag+26))
dpro%lndref = dpro%lndref + utdep(N1%beta(cut:cut+nut-1,0),gf(cut:cut+8))
dpro%zref = zetaF
dpro%zmin = 110d0
dpro%zhyd = zetaF
! Effective mass
dpro%zetaM = dot_product(N1%beta(0:mbf,1),gf(0:mbf))
dpro%HML = dot_product(N1%beta(0:mbf,2),gf(0:mbf))
dpro%HMU = dot_product(N1%beta(0:mbf,3),gf(0:mbf))
! Chapman correction
dpro%C = dot_product(N1%beta(0:mbf,4),gf(0:mbf))
dpro%zetaC = dot_product(N1%beta(0:mbf,5),gf(0:mbf))
dpro%HC = dot_product(N1%beta(0:mbf,6),gf(0:mbf))
! Dynamical correction
dpro%R = dot_product(N1%beta(0:mbf,7),gf(0:mbf))
dpro%zetaR = dot_product(N1%beta(0:mbf,8),gf(0:mbf))
dpro%HR = dot_product(N1%beta(0:mbf,9),gf(0:mbf))
! ! Anomalous Oxygen ----------------------
case(9)
dpro%lndref = dot_product(OA%beta(0:mbf,0),gf(0:mbf))
dpro%lndref = dpro%lndref + geomag(OA%beta(cmag:cmag+nmag-1,0),gf(cmag:cmag+12),gf(cmag+13:cmag+26))
dpro%zref = zetarefOA
dpro%zmin = 110d0
dpro%zhyd = 0d0
dpro%C = OA%beta(0,4)
dpro%zetaC = OA%beta(0,5)
dpro%HC = OA%beta(0,6)
return !No further parameters needed for legacy anomalous oxygen profile
! ! Nitric Oxide ----------------------
case(10)
! Reference number density
dpro%lnPhiF = 0d0
dpro%lndref = dot_product(NO%beta(0:mbf,0),gf(0:mbf))
dpro%zref = zetarefNO
dpro%zmin = nodesNO(3)
dpro%zhyd = zetarefNO
! Effective mass
dpro%zetaM = dot_product(NO%beta(0:mbf,1),gf(0:mbf))
dpro%HML = dot_product(NO%beta(0:mbf,2),gf(0:mbf))
dpro%HMU = dot_product(NO%beta(0:mbf,3),gf(0:mbf))
! Chapman correction
dpro%C = dot_product(NO%beta(0:mbf,4),gf(0:mbf))
dpro%zetaC = dot_product(NO%beta(0:mbf,5),gf(0:mbf))
dpro%HC = dot_product(NO%beta(0:mbf,6),gf(0:mbf))
! Dynamical correction
dpro%R = dot_product(NO%beta(0:mbf,7),gf(0:mbf))
dpro%zetaR = dot_product(NO%beta(0:mbf,8),gf(0:mbf))
dpro%HR = dot_product(NO%beta(0:mbf,9),gf(0:mbf))
! Unconstrained splines
do izf = 0,nsplNO-1
dpro%cf(izf) = dot_product(NO%beta(0:mbf,izf+10),gf(0:mbf))
enddo
!Constrained splines calculated after case statement
! Failsafe ----- ---------------------------
case default
stop 'Species not yet implemented'
endselect
!Compute piecewise mass profile values and integration terms
dpro%zetaMi(0) = dpro%zetaM - 2d0*dpro%HML
dpro%zetaMi(1) = dpro%zetaM - dpro%HML
dpro%zetaMi(2) = dpro%zetaM
dpro%zetaMi(3) = dpro%zetaM + dpro%HMU
dpro%zetaMi(4) = dpro%zetaM + 2d0*dpro%HMU
dpro%Mi(0) = mbar
dpro%Mi(4) = amu(massid)
dpro%Mi(2) = (dpro%Mi(0) + dpro%Mi(4)) / 2d0
delM = tanh1 * (dpro%Mi(4) - dpro%Mi(0)) / 2d0
dpro%Mi(1) = dpro%Mi(2) - delM
dpro%Mi(3) = dpro%Mi(2) + delM
do i = 0, 4
i1 = i + 1
if (i .lt. 4) dpro%aMi(i) = (dpro%Mi(i1) - dpro%Mi(i)) / (dpro%zetaMi(i1) - dpro%zetaMi(i))
delz = dpro%zetaMi(i) - zetaB
if (dpro%zetaMi(i) .lt. zetaB) then
call bspline(dpro%zetaMi(i),nodes,nd+2,6,eta,Si,iz)
dpro%WMi(i) = dot_product(tpro%gamma(iz-5:iz),Si(:,6)) + tpro%cVS*delz + tpro%cWS
else
dpro%WMi(i) = (0.5d0*delz*delz + dilog(tpro%b*exp(-tpro%sigma*delz))/tpro%sigmasq)/tpro%tex &
+ tpro%cVB*delz + tpro%cWB
endif
end do
dpro%XMi(0) = -dpro%aMi(0) * dpro%WMi(0)
do i = 1, 3
dpro%XMi(i) = dpro%XMi(i-1) - dpro%WMi(i) * (dpro%aMi(i) - dpro%aMi(i-1))
end do
dpro%XMi(4) = dpro%XMi(3) + dpro%WMi(4) * dpro%aMi(3)
!Calculate hydrostatic integral at reference height, and copy temperature
if (dpro%zref .eq. zetaF) then
Mzref = mbar
dpro%Tref = tpro%TzetaF
dpro%Izref = mbar * tpro%VzetaF
else if (dpro%zref .eq. zetaB) then
Mzref = pwmp(dpro%zref,dpro%zetaMi,dpro%Mi,dpro%aMi)
dpro%Tref = tpro%Tb0
dpro%Izref = 0.0d0
if ((zetaB .gt. dpro%zetaMi(0)) .and. (zetaB .lt. dpro%zetaMi(4))) then
i = 0
do i1 = 1, 3
if (zetaB .lt. dpro%zetaMi(i1)) then
exit
else
i = i1
endif
enddo
dpro%Izref = dpro%Izref - dpro%XMi(i)
else
dpro%Izref = dpro%Izref - dpro%XMi(4)
endif
else if (dpro%zref .eq. zetaA) then
Mzref = pwmp(dpro%zref,dpro%zetaMi,dpro%Mi,dpro%aMi)
dpro%Tref = tpro%TzetaA
dpro%Izref = Mzref * tpro%VzetaA
if ((zetaA .gt. dpro%zetaMi(0)) .and. (zetaA .lt. dpro%zetaMi(4))) then
i = 0
do i1 = 1, 3
if (zetaA .lt. dpro%zetaMi(i1)) then
exit
else
i = i1
endif
enddo
dpro%Izref = dpro%Izref - (dpro%aMi(i)*tpro%WzetaA + dpro%XMi(i))
else
dpro%Izref = dpro%Izref - dpro%XMi(4)
endif
else
stop 'Integrals at reference height not available'
endif
!C1 constraint for O1 at 85 km
if (massid .eq. 4) then
Cterm = dpro%C*exp(-(dpro%zref - dpro%zetaC)/dpro%HC)
bc(1) = dpro%lndref - Cterm - dpro%cf(7)*c1o1adj(1) !Reference density, Chapman term, and subtraction of last unconstrained spline(7)
bc(2) = -Mzref*g0divR/tpro%tzetaA & !Gradient of hydrostatic term
-tpro%dlntdzA & !Gradient of ideal gas law term
+Cterm/dpro%HC & !Gradient of Chapman term
-dpro%cf(7)*c1o1adj(2) !Subtraction of gradient of last unconstrained spline(7)
! Compute coefficients for constrained splines
dpro%cf(8:9) = matmul(bc,c1o1)
endif
!C1 constraint for NO at 122.5 km
if (massid .eq. 10) then
Cterm = dpro%C*exp(-(dpro%zref - dpro%zetaC)/dpro%HC)
bc(1) = dpro%lndref - Cterm - dpro%cf(7)*c1noadj(1) !Reference density, Chapman term, and subtraction of last unconstrained spline(7)
bc(2) = -Mzref*g0divR/tpro%tb0 & !Gradient of hydrostatic term
-tpro%tgb0/tpro%tb0 & !Gradient of ideal gas law term
+Cterm/dpro%HC & !Gradient of Chapman term
-dpro%cf(7)*c1noadj(2) !Subtraction of gradient of last unconstrained spline(7)
! Compute coefficients for constrained splines
dpro%cf(8:9) = matmul(bc,c1no)
endif
return
end subroutine dfnparm
! =======================================================================
! Evaluate density from the profile parameters
! =======================================================================
real(8) function dfnx(z,tnz,lndtotz,Vz,Wz,tpro,dpro)
use constants
use msis, only : etaO1, etaNO
use tfn,only : tnparm
implicit none
real(8),intent(in) :: z
real(8),intent(in) :: tnz,lndtotz !Temperature, total number density at input z
real(8),intent(in) :: Vz,Wz !First and second indefinite integrals of 1/T at z
type(tnparm),intent(in) :: tpro
type(dnparm),intent(in) :: dpro
integer(4) :: i,i1,iz
real(8) :: Mz
real(8) :: Sz(-5:0,2:6)
real(8) :: Ihyd !Hydrostatic definite integral
! Below minimum height of profile
if (z .lt. dpro%zmin) then
dfnx = 1.0d-32
return
endif
! Anomalous Oxygen (legacy MSISE-00 formulation)
if (dpro%massid .eq. 9) then
dfnx = dpro%lndref - (z - dpro%zref)/HOA - dpro%C*exp(-(z-dpro%zetaC)/dpro%HC)
dfnx = exp(dfnx)
return !No further calculation needed for anomalous oxygen
endif
! Below height where hydrostatic terms are needed
if (z .lt. dpro%zhyd) then
select case(dpro%massid)
case(2,3,5,7) !For N2, O2, He, and Ar, apply mixing ratios and exit
dfnx = exp(lndtotz + dpro%lnPhiF)
return
case(4) !For O, evaluate splines
call bspline(z,nodesO1,ndO1,4,etaO1,Sz,iz)
dfnx = exp(dot_product(dpro%cf(iz-3:iz),Sz(-3:0,4)))
return
case(10) !For NO, evaluate splines
call bspline(z,nodesNO,ndNO,4,etaNO,Sz,iz)
dfnx = exp(dot_product(dpro%cf(iz-3:iz),Sz(-3:0,4)))
return
endselect
endif
! Calculate hydrostatic term and apply to reference density
Mz = pwmp(z,dpro%zetaMi,dpro%Mi,dpro%aMi)
Ihyd = Mz * Vz - dpro%Izref
if ((z .gt. dpro%zetaMi(0)) .and. (z .lt. dpro%zetaMi(4))) then
i = 0
do i1 = 1, 3
if (z .lt. dpro%zetaMi(i1)) then
exit
else
i = i1
endif
enddo
Ihyd = Ihyd - (dpro%aMi(i)*Wz + dpro%XMi(i))
else if (z .ge. dpro%zetaMi(4)) then
Ihyd = Ihyd - dpro%XMi(4)
endif
dfnx = dpro%lndref - Ihyd * g0divR
! Apply Chapman and logistic corrections
select case(dpro%massid)
case(2,3,5,7) !For N2, O2, He, and Ar, logistic correction only
dfnx = dfnx + dpro%R*(1+tanh((z-dpro%zetaR)/(2d0*dpro%HR)))/2d0
case(4,6,8,10) !For O, H, N, NO Chapman and logistic corrections
dfnx = dfnx - dpro%C*exp(-(z-dpro%zetaC)/dpro%HC) &
+ dpro%R*(1+tanh((z-dpro%zetaR)/(2d0*dpro%HR)))/2d0
endselect
! Apply ideal gas law
dfnx = exp(dfnx) * dpro%Tref/tnz
return
end function dfnx
! -----------------------------------------------------
! Piecewise mass profile interpolation
! -----------------------------------------------------
real(8) function pwmp(z,zm,m,dmdz)
real(8), intent(in) :: z
real(8), intent(in) :: zm(0:4)
real(8), intent(in) :: m(0:4)
real(8), intent(in) :: dmdz(0:3)
integer :: irng !Index of piecwise interval
integer :: inode
! Most probable case
if (z .ge. zm(4)) then
pwmp = m(4)
return
endif
! Second most probable case
if (z .le. zm(0)) then
pwmp = m(0)
return
endif
! None of the above
do inode = 0,3
if (z .lt. zm(inode+1)) then
pwmp = m(inode) + dmdz(inode)*(z - zm(inode))
return
endif
enddo
! If we are here this is a problem
stop 'Error in pwmp'
end function pwmp
end module dfn