-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathModInputs.f90
379 lines (281 loc) · 11.1 KB
/
ModInputs.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
module ModInputs
! ------------------------------------------------------------------------------
! Toggle Paramaters to test
! real, parameter :: PLONG = 1.5E-1 ! Pa (also 1.5 ubar): 74 km (aphelion)
! real, parameter :: PLONG = 4.0E-1 ! Pa (also 4.0 ubar): 60 km (aphelion) (on)
! real, parameter :: RPTAU = 6.1 mbar (standard) (on)
! AltMinIono = 60.0 km ! revised for NOx chemsitry (170310: S. W. Bougher)
! New Logicals for (ON) Secondary Ionization and (ON) EUVData FISM fluxes
! ------------------------------------------------------------------------------
use ModConstants
use GITM_planet, only : nSpecies, Rotation_Period
use ModIoUnit, only : UnitTmp_
use ModKind, only: Real8_
implicit none
integer, parameter :: iCharLen_ = 100
integer :: iOutputUnit_ = UnitTmp_
integer :: iInputUnit_ = UnitTmp_
integer :: iRestartUnit_ = UnitTmp_
integer :: iLogFileUnit_ = 92
logical :: IsOpenLogFile = .false.
real :: DtLogFile = 60.0
integer, parameter :: nInputMaxLines = 10000
integer :: nInputLines
character (len=iCharLen_) :: cInputText(nInputMaxLines)
character (len=iCharLen_) :: cInputFile = "UAM.in"
character (len=iCharLen_) :: cAMIEFileSouth
character (len=iCharLen_) :: cAMIEFileNorth
character (len=iCharLen_) :: PotentialModel
character (len=iCharLen_) :: AuroralModel
logical :: UseIMF = .true.
logical :: UseHpi = .true.
logical :: UseNewellAurora = .false.
logical :: UseNewellAveraged = .true.
logical :: UseNewellMono = .false.
logical :: UseNewellWave = .false.
logical :: DoNewellRemoveSpikes = .true.
logical :: DoNewellAverage = .true.
character (len=iCharLen_) :: TypeLimiter = "minmod"
integer, dimension(7) :: iStartTime
real :: CPUTimeMax = 1.0e32
logical :: UseVariableInputs = .false.
logical :: IsFramework = .false.
logical :: DoRestart = .false.
integer :: iAltTest = -1
integer :: iDebugLevel = 0
logical :: UseBarriers = .false.
integer :: nSteps = 10
real :: CFL = 0.25
integer :: nOutputTypes = 0
integer, parameter :: nMaxOutputTypes = 50
character (len=iCharLen_), dimension(nMaxOutputTypes) :: OutputType
real :: DtPlot(nMaxOutputTypes)
real :: DtPlotSave(nMaxOutputTypes)
real :: PlotTimeChangeDt(nMaxOutputTypes) = 1.0
real(Real8_) :: PlotTimeChangeStart = 0.0, PlotTimeChangeEnd = 0.0
logical :: DoAppendFiles = .false.
real :: DtRestart = 60.0*60.0
real :: DtReport = 1.0*60.0
real :: DtAurora = 60.0*1.0
real :: DtPotential = 60.0*1.0
real :: DtGlow = 60.0
real :: f107 = 150.0
real :: f107a = 150.0
integer :: iModelSolar = 0
real :: AltMin = 95.0 * 1000.0
real :: AltMax = 500.0 * 1000.0
real :: BetaLimiter = 2.0
real :: ConcentrationLatitude = 45.0
real :: StretchingPercentage = 0.0
real :: StretchingFactor = 1.0
logical :: UseTopography = .false.
! real :: AltMinIono=80.0 ! in km
real :: AltMinIono=60.0 ! in km
real :: TempMax = 1000.0
real :: TempMin = 200.0
real :: TempWidth = 25.0*1e3
real :: TempHeight = 150.0*1e3
real :: LogRho0
integer :: nBlocksLon = 2
integer :: nBlocksLat = 2
logical :: Is1D = .false.
logical :: IsFullSphere = .false.
logical :: UseGlow = .false.
real :: LonStart = 0.0
real :: LonEnd = 0.0
real :: LatStart = -pi/4.0
real :: LatEnd = pi/4.0
logical :: UseStretchedAltitude = .true.
logical :: UseApex = .true.
logical :: UseMSIS = .true.
logical :: UseIRI = .true.
logical :: UseMSISTides = .true.
logical :: UseMSISOnly = .false.
logical :: UseGSWMTides = .false.
logical :: UseWACCMTides = .false.
logical :: UseStatisticalModelsOnly = .false.
real :: DtStatisticalModels = 3600.0
logical :: UseGswmComp(4) = .true.
real :: MagneticPoleRotation = 0.0
real :: MagneticPoleTilt = 0.0
real :: xDipoleCenter = 0.0
real :: yDipoleCenter = 0.0
real :: zDipoleCenter = 0.0
!\
! These are things for the ion precipitation for the April 2002 storm:
!/
logical :: UseIonPrecipitation = .false.
character (len=iCharLen_) :: IonIonizationFilename
character (len=iCharLen_) :: IonHeatingRateFilename
logical :: UseEddyInSolver = .false.
logical :: UseNeutralFrictionInSolver = .false.
real :: MaximumVerticalVelocity = 1000.0
logical :: UseDamping = .false.
!! EddyVelocity Terms
logical :: UseBoquehoAndBlelly = .false.
logical :: UseEddyCorrection = .false.
!\
! These are logicals to turn off and on forcing terms:
!/
logical :: UsePressureGradient = .true.
logical :: UseGravity = .true.
logical :: UseIonDrag = .true.
logical :: UseViscosity = .true.
logical :: UseCoriolis = .true.
logical :: UseGravityWave = .false.
logical :: UseHorAdvection = .true.
logical :: UseVerAdvection = .true.
logical :: UseNeutralFriction = .true.
logical :: UseCrustalField = .false.
logical :: UseMHDField = .false.
logical :: crustalFieldOnly = .false.
logical :: UseIonPressureGradient = .true.
logical :: UseIonGravity = .true.
logical :: UseNeutralDrag = .true.
logical :: UseExB = .true.
logical :: UseImplicitChemistry = .false.
logical :: IsAsymmetric = .false.
Real :: BetaPointImpl = 1.0
logical :: UseEmpiricalIonization = .false.
integer :: DefaultFieldType = 4
logical :: UseDynamo = .false.
real :: DynamoHighLatBoundary = 65.0
integer :: nItersMax = 500
real :: MaxResidual = 1.0
logical :: UseSolarHeating = .true.
logical :: UseJouleHeating = .true.
logical :: UseAuroralHeating = .true.
logical :: UseNOCooling = .true.
logical :: UseOCooling = .true.
logical :: UseConduction = .true.
logical :: UseDiffusion = .false.
logical :: UseVerAdvectionT = .true.
!! New Turbulent (Eddy) Conduction Routines
logical :: UseTurbulentCond = .false.
logical :: UseUpdatedTurbulentCond = .false.
real :: EddyScaling = 1.0
character (len=iCharLen_), dimension(3) :: EddyDiffusionTypes(3)
character (len=iCharLen_) :: EddyDiffusionMethod = 'minmax'
!kMax is used as kConstant for the constant method
real :: kEddyMin = 500, kEddyMax = 1500
!! WAVE DRAG FORCINGS
logical :: UseStressHeating = .false.
real :: KappaTemp0 = 5.6e-4
real :: Kappa1DCorrectionFactor = 45.0
real :: Kappa1DCorrectionPower = 1.75
logical :: UseKappa1DCorrection = .false.
logical :: UseIonChemistry = .true.
logical :: UseNeutralChemistry = .true.
logical :: UseIonAdvection = .true.
logical :: DoCheckStopFile = .true.
real :: MaxVParallel = 100.0
!\
! Methods for completing chemistry
!/
integer, parameter :: cSteadyStateChemType_ = 1
integer, parameter :: cImplicitChemType_ = 2
integer, parameter :: cSubCycleChemType_ = 3
integer, parameter :: nChemTypes_ = 3
character (len=100), dimension(nChemTypes_) :: sChemType
character (len=100) :: sInputIonChemType
character (len=100) :: sInputNeutralChemType
integer :: iInputIonChemType, iInputNeutralChemType,nMHDFiles
real :: LogNS0(nSpecies)
logical :: UseEUVData =.false.
character (len=iCharLen_) :: cEUVFile
character (len=iCharLen_) :: cMHDFile,cMHDFilelist
character (len=100), dimension(1000) :: MHDFiles
real :: dtBMHD = 60*15 !15 minutes by default
! These are Mars Specific, but ignored by other codes:
! Some are modified in Planet.f90 (set_planet_defaults)
real :: DtLTERadiation = 100.0*Rotation_Period
!!!!!!!!!!!!!!! DUST
Logical :: UseDustDistribution = .False.
! Setting the depth to which dust is mixed based on a reference dust opacity
!This can now be read in as a horizontal distribution, or as a constant value
!using UAM.in
character (len=iCharLen_) :: cDustFile
character (len=iCharLen_) :: cConrathFile
real :: dtDust
character (len=iCharLen_) :: DustFileType
! Global mean dust opacity and CR (standard))
! real :: CONRNU_temp = 0.03 ! Standard value ~25km half-height
! real :: TAUTOT_temp = 0.3 !do not set to 0.0 or less
! Global mean dust opacity and CR (Viking)
real :: CONRNU_temp = 0.003 ! ~50 km half-height
real :: TAUTOT_temp = 0.5 !do not set to 0.0 or less
!!!!!!!!!!!!!!!!!!!!!!
! RPTAU : Reference Pressure optical depth; 6.1 mbar for now
real, parameter :: RPTAU = 6.1
! RPTAU : Reference Pressure optical depth; 7.0 mbar for Brecht tests
! real, parameter :: RPTAU = 7.0
! Top of the shortwave calculation for lower atmosphere radiation code(mbars)
real, parameter :: PRAD = 1.0E-6 !mb (~115 km)
! Top of the longwave calculation for lower atmosphere radiation code(PASCALS)
! and where radcool begins
! real, parameter :: PLONG = 1.5E-1 ! Pa (also 1.5 ubar): 74 km (aphelion)
real, parameter :: PLONG = 4.0E-1 ! Pa (also 4.0 ubar): 60 km (aphelion)
!!!!!!!!!!!!!!!!!!!!!!
! Mars Secondary Ionization Logicals
! logical :: UseSecondaryIonization = .false.
! character (len=iCharLen_) :: SecondaryIonFile
logical :: UseWValue = .false.
real :: WValue = 28.0
! Mars EUVdata Logicals
logical :: UseFluxAtPlanet = .true.
!!!!!!!!!!!!!!!!!!!!!!
! Ls Variables
real, dimension(7) :: Ls_a
real, dimension(7) :: Ls_tau
real, dimension(7) :: Ls_phi
DATA Ls_a / 0.007, 0.006, 0.004, 0.004, 0.002, 0.002, 0.002 /
DATA Ls_tau / 2.2353, 2.7543, 1.1177, 15.7866, 2.1354, &
2.4694, 32.8493 /
DATA Ls_phi / 49.409, 168.173, 191.837, 21.736, 15.704, &
95.528, 49.095 /
contains
! -------------------------------------------------------------
! set_strings
! Here we initialize all of the strings which need to be
! initialized in the code.
! -------------------------------------------------------------
subroutine set_strings
sChemType(cSteadyStateChemType_) = "steady"
sChemType(cImplicitChemType_) = "implicit"
sChemType(cSubCycleChemType_) = "subcycle"
end subroutine set_strings
! -------------------------------------------------------------
! set_defaults
! -------------------------------------------------------------
subroutine set_defaults
use ModTime
call set_strings
sInputIonChemType = sChemType(cSubCycleChemType_)
sInputNeutralChemType = sChemType(cSubCycleChemType_)
iInputIonChemType = cSubCycleChemType_
iInputNeutralChemType = cSubCycleChemType_
PotentialModel = "Weimer05"
AuroralModel = "ihp"
dTAurora = 120.0
tSimulation = 0.0
iStartTime(1) = 1999
iStartTime(2) = 3
iStartTime(3) = 21
iStartTime(4) = 0
iStartTime(5) = 0
iStartTime(6) = 0
iStartTime(7) = 0
iTimeArray = iStartTime
call time_int_to_real(iStartTime, CurrentTime)
nOutputTypes = 0
DtPlot = -1.0
DtRestart = -1.0
! Initiate Neutral Variables From MSIS90
f107a = 150.0
f107 = 150.0
iModelSolar = 0
UseApex = .false.
DoRestart = .false.
call set_planet_defaults
end subroutine set_defaults
end module ModInputs