-
Notifications
You must be signed in to change notification settings - Fork 10
/
glowbasic.f90
153 lines (138 loc) · 4.91 KB
/
glowbasic.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
program glowbasic
! This software is part of the GLOW model. Use is governed by the Open Source
! Academic Research License Agreement contained in the file glowlicense.txt.
! For more information see the file glow.txt.
! Version 0.981, 6/2017
! Adapted from glowdriver by Stan Solomon, 2/2016
! Basic single-processor driver for the GLOW model.
! Uses MSIS/IRI for input.
! Runs GLOW for designated inputs once, or multiple times.
! MPI and netCDF libraries not required.
! For definitions of use-associated variables, see subroutine GLOW and module CGLOW.
! Other definitions:
! f107p Solar 10.7 cm flux for previous day
! ap Ap index of geomagnetic activity
! z altitude array, km
! Array dimensions:
! jmax number of altitude levels
! nbins number of energetic electron energy bins
! lmax number of wavelength intervals for solar flux
! nmaj number of major species
! nst number of states produced by photoionization/dissociation
! nei number of states produced by electron impact
! nex number of ionized/excited species
! nw number of airglow emission wavelengths
! nc number of component production terms for each emission
use cglow,only: jmax,nbins,lmax,nmaj,nei,nex,nw,nc,nst
use cglow,only: idate,ut,glat,glong,f107a,f107,f107p,ap,ef,ec
use cglow,only: iscale,jlocal,kchem,xuvfac
use cglow,only: sza,dip,efrac,ierr
use cglow,only: zz,zo,zn2,zo2,zns,znd,zno,ztn,ze,zti,zte
use cglow,only: ener,del,phitop,wave1,wave2,sflux,pespec,sespec,uflx,dflx,sion
use cglow,only: photoi,photod,phono,aglw,tei,tpi,tir,ecalc,zxden,zeta,zceta,zlbh
use cglow,only: cglow_init
use cglow,only: data_dir
implicit none
character(len=1024) :: iri90_dir
real,allocatable :: z(:) ! glow height coordinate in km (jmax)
real,allocatable :: zun(:), zvn(:) ! neutral wind components (not in use)
real,allocatable :: pedcond(:), hallcond(:) ! Pederson and Hall conductivities in S/m (mho)
real,allocatable :: outf(:,:) ! iri output (11,jmax)
real :: rz12,stl,fmono,emono
real :: d(8), t(2), sw(25), oarr(30)
integer :: l,j,jj,ijf,jmag,iday,mmdd,i,ii,n,k,ix,itail
integer :: instance,iostatus
logical :: jf(12)
data sw/25*1./
!
! Initialize standard switches:
!
iscale=1
xuvfac=3.
kchem=4
jlocal=0
itail=0
fmono=0.
emono=0.
!
! Set data directories:
!
data_dir = 'data/'
iri90_dir = 'data/iri90/'
!
! Set number of altitude levels:
!
jmax = 102
!
! Allocate local arrays:
!
allocate(z(jmax))
allocate(zun(jmax))
allocate(zvn(jmax))
allocate(pedcond(jmax))
allocate(hallcond(jmax))
allocate(outf(11,jmax))
!
! Call CGLOW_INIT (module CGLOW) to set array dimensions and allocate use-associated variables:
! (This was formerly done using common blocks, including common block /cglow/.)
!
call cglow_init
!
! Call EGRID to set up electron energy grid:
!
call egrid (ener, del, nbins)
!
! Loop to call GLOW for designated inputs until end-of-file or any character on standard input:
!
do instance=1,10000
!
! Get input values:
!
write(6,"('Enter date, UT, lat, lon, F107a, F107, F107p, Ap, Ef, Ec')")
read(5,*,iostat=iostatus) idate,ut,glat,glong,f107a,f107,f107p,ap,ef,ec
if (iostatus /= 0) stop
!
! Calculate local solar time:
!
stl = ut/3600. + glong/15.
if (stl < 0.) stl = stl + 24.
if (stl >= 24.) stl = stl - 24.
!
! Call MZGRID to use MSIS/NOEM/IRI inputs on default altitude grid:
!
call mzgrid (jmax,nex,idate,ut,glat,glong,stl,f107a,f107,f107p,ap,iri90_dir, &
z,zo,zo2,zn2,zns,znd,zno,ztn,zun,zvn,ze,zti,zte,zxden)
!
! Call MAXT to put auroral electron flux specified by namelist input into phitop array:
!
phitop(:) = 0.
if (ef>.001 .and. ec>1.) call maxt (ef,ec,ener,del,nbins,itail,fmono,emono,phitop)
!
! Fill altitude array, converting to cm:
!
zz(:) = z(:) * 1.e5 ! km to cm at all jmax levels
!
! Call GLOW to calculate ionized and excited species, airglow emission rates,
! and vertical column brightnesses:
!
call glow
!
! Call CONDUCT to calculate Pederson and Hall conductivities:
!
do j=1,jmax
call conduct (glat, glong, z(j), zo(j), zo2(j), zn2(j), &
zxden(3,j), zxden(6,j), zxden(7,j), ztn(j), zti(j), zte(j), &
pedcond(j), hallcond(j))
enddo
!
! Output section:
!
write(6,"(1x,i7,9f8.1)") idate,ut,glat,glong,f107a,f107,f107p,ap,ef,ec
write(6,"(' Z Tn O N2 NO Ne(in) Ne(out) Ionrate O+ O2+ NO+ N(2D) Pederson Hall')")
write(6,"(1x,0p,f5.1,f6.0,1p,12e10.2)") (z(j),ztn(j),zo(j),zn2(j),zno(j),ze(j), &
ecalc(j),tir(j),zxden(3,j),zxden(6,j),zxden(7,j),zxden(10,j),pedcond(j),hallcond(j),j=1,jmax)
write(6,"(' Z 3371 4278 5200 5577 6300 7320 10400 3644 7774 8446 3726 LBH 1356 1493 1304')")
write(6,"(1x,f5.1,15f8.2)")(z(j),(zeta(ii,j),ii=1,15),j=1,jmax)
enddo
stop
end program glowbasic