-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathcalc_rates.Jupiter.f90
428 lines (325 loc) · 14.7 KB
/
calc_rates.Jupiter.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
subroutine calc_rates(iBlock)
! ----------------------------------------------------------------------------
! Modified (01/29/07) : SWB : Add cp, kt, km formulation from vtgcm2d codes.
! Variables from Mars GITM used. CGS to MKS conversion.
! Modified (02/01/07) : SWB : Add eddy thermal conductivity (reduced by Prandtl)
! Modified (02/06/07) : SWB : Dimensioned array math redone correctly!
! Modified (02/07/07) : SWB : RGAS specified (cgs); access arrays explicitly with
! all indices inside loops
! Modified (02/20/07) : SWB : Conversion from cgs to mks for Kt and Km corrected.
! Comment out write statements (no longer needed)
! ----------------------------------------------------------------------------
use ModGITM
use ModRates
use ModConstants
use GITM_planet
use ModInputs
use ModEUV, only : SunOrbitEccentricity
implicit none
integer, intent(in) :: iBlock
integer :: iAlt, iIon, iSpecies, iError, iiAlt, iLat,iLon
real, dimension(nLons, nLats, nAlts) :: &
Tn, Ti, TWork1, TWork2, TWork3, NO2
real, dimension(-1:nLons+2, -1:nLats+2, -1:nAlts+2) :: &
Ne, mnd, Te, tmp
real :: ScaleHeight(nLons, nLats)
real :: e2
! ------------------------------------------------------------------------------
! cpktkm.F add-ons
integer :: is
real :: rrco2,cpco2,crn,prco
real, dimension (8) :: cmrf, com
real, dimension(-1:nLons+2, -1:nLats+2, -1:nAlts+2) :: &
po, pco, pco2, pn2, cpmix, ktmix, kmmix, &
ttot,cokm,co2kt,cokt,tt,co2km
! ------------------------------------------------------------------------------
logical :: trouble
call report("calc_rates",2)
call start_timing("calc_rates")
trouble = .false.
! ------------------------------------------------------------------------------
cmrf = (/135.8,185.7,230.4,271.0,308.3,343.3,373.7,406.1/)
com = (/0.3966,0.7692,1.0776,1.340,1.574,1.787,1.986,2.172/)
! -------------------------------------------------------------------------------
! write(*,*) '1st, Some Preliminary Diagnostics for Calc_Rates ===----------+'
do iLon = -1,nLons+2
do iLat = -1,nLats+2
do iAlt = -1,nAlts+2
if( .not. (Temperature(iLon,iLat,iAlt,iBlock) < 0.0) .and. &
.not. (Temperature(iLon,iLat,iAlt,iBlock) > 0.0) ) then
write(*,*)'Temperature(',iLon,iLat,iAlt,iBlock,') = ',&
Temperature(iLon,iLat,iAlt,iBlock)
trouble = .true.
elseif( Temperature(iLon,iLat,iAlt,iBlock) <= -1.0e+300) then
write(*,*)'Temperature(',iLon,iLat,iAlt,iBlock,') = ',&
Temperature(iLon,iLat,iAlt,iBlock)
trouble = .true.
elseif(Temperature(iLon,iLat,iAlt,iBlock) >= 1.0e+300) then
write(*,*)'Temperature(',iLon,iLat,iAlt,iBlock,') = ',&
Temperature(iLon,iLat,iAlt,iBlock)
trouble = .true.
elseif(Temperature(iLon,iLat,iAlt,iBlock) <= 0.0) then
write(*,*)'Temperature(',iLon,iLat,iAlt,iBlock,') = ',&
Temperature(iLon,iLat,iAlt,iBlock)
trouble = .true.
endif
do iSpecies = 1,nSpeciesTotal
if( .not. (NDensityS(iLon,iLat,iAlt,iSpecies,iBlock) < 0.0) .and. &
.not. (NDensityS(iLon,iLat,iAlt,iSpecies,iBlock) > 0.0) ) then
write(*,*)'NDensityS(',iLon,iLat,iAlt,iSpecies,iBlock,') = ',&
NDensityS(iLon,iLat,iAlt,iSpecies,iBlock)
trouble = .true.
elseif(NDensityS(iLon,iLat,iAlt,iSpecies,iBlock) <= -1.0e+300) then
write(*,*)'NDensityS(',iLon,iLat,iAlt,iSpecies,iBlock,') = ',&
NDensityS(iLon,iLat,iAlt,iSpecies,iBlock)
trouble = .true.
elseif(NDensityS(iLon,iLat,iAlt,iSpecies,iBlock) >= 1.0e+300) then
write(*,*)'NDensityS(',iLon,iLat,iAlt,iSpecies,iBlock,') = ',&
NDensityS(iLon,iLat,iAlt,iSpecies,iBlock)
trouble = .true.
elseif(NDensityS(iLon,iLat,iAlt,iSpecies,iBlock) <= 0.0) then
write(*,*)'NDensityS(',iLon,iLat,iAlt,iSpecies,iBlock,') = ',&
NDensityS(iLon,iLat,iAlt,iSpecies,iBlock)
trouble = .true.
endif
enddo !iSpecies = 1,nSpeciesTotal
if( .not. (NDensity(iLon,iLat,iAlt,iBlock) < 0.0) .and. &
.not. (NDensity(iLon,iLat,iAlt,iBlock) > 0.0) ) then
write(*,*)'NDensity(',iLon,iLat,iAlt,iBlock,') = ',&
NDensity(iLon,iLat,iAlt,iBlock)
trouble = .true.
elseif( NDensity(iLon,iLat,iAlt,iBlock) <= -1.0e+300) then
write(*,*)'NDensity(',iLon,iLat,iAlt,iBlock,') = ',&
NDensity(iLon,iLat,iAlt,iBlock)
trouble = .true.
elseif(NDensity(iLon,iLat,iAlt,iBlock) >= 1.0e+300) then
write(*,*)'NDensity(',iLon,iLat,iAlt,iBlock,') = ',&
NDensity(iLon,iLat,iAlt,iBlock)
trouble = .true.
elseif(NDensity(iLon,iLat,iAlt,iBlock) <= 0.0) then
write(*,*)'NDensity(',iLon,iLat,iAlt,iBlock,') = ',&
NDensity(iLon,iLat,iAlt,iBlock)
trouble = .true.
endif
enddo
enddo
enddo
! -------------------------------------------------------------------------------
if(trouble) then
write(*,*) 'trouble found!!'
write(*,*) 'Stop GITM'
stop
endif
! -------------------------------------------------------------------------------
if (UseBarriers) call MPI_BARRIER(iCommGITM,iError)
if (iDebugLevel > 4) write(*,*) "=====> mean major mass", iblock
!write(*,*) '==> calc_rates: Before NDensityS Statements.'
where(NDensityS(:,:,:,:,iBlock) < 1.0e3)
NDensityS(:,:,:,:,iBlock) = 1.0e3
end where
!write(*,*) '==> calc_rates: Before IDensityS Statements.'
! where(IDensityS < 1.0e3)
! IDensityS = 1.0e3
! end where
!write(*,*) '==> calc_rates: Before Ne Set.'
Ne = IDensityS(:,:,:,ie_,iBlock)
! We add 1 because this is in the denominator a lot, and the corners
! don't have anything. Total number density.
mnd = NDensity(:,:,:,iBlock)+1.0
!write(*,*) '==> calc_rates: Before MeanMajorMass Calculation.'
MeanIonMass = 0.0
MeanMajorMass = 0.0
do iSpecies = 1, nSpecies
MeanMajorMass = MeanMajorMass + &
Mass(iSpecies) * &
NDensityS(:,:,:,iSpecies,iBlock)/mnd
enddo
! MMM_3D(1:nLons,1:nLats,1:nAlts,iBlock) = &
! MeanMajorMass(1:nLons,1:nLats,1:nAlts)/AMU
! Once again, in the corners, the meanmajormass is 0.
where (MeanMajorMass == 0) &
MeanMajorMass = Mass(1)
!write(*,*) '==> calc_rates: Before MeanIonMass Calculation.'
do iIon = 1, nIons-1
MeanIonMass = MeanIonMass + &
MassI(iIon) * IDensityS(:,:,:,iIon,iBlock) / Ne
enddo
! -------------------------------------------------------------------------------
TempUnit = MeanMajorMass / Boltzmanns_Constant
! -------------------------------------------------------------------------------
!write(*,*) '==> calc_rates: Before Mixing Ratio Calculation.'
! Mixing Ratios needed for Kt and Km calculations below
! Temperature Arrays needed for Kt amd Km calculation below
do iAlt = 0, nAlts+1
! Mixing Ratios
po(:,:,iAlt) = NDensityS(:,:,iAlt,iO_,iBlock)/mnd(:,:,iAlt)
pco(:,:,iAlt) = NDensityS(:,:,iAlt,iCO_,iBlock)/mnd(:,:,iAlt)
pn2(:,:,iAlt) = NDensityS(:,:,iAlt,iN2_,iBlock)/mnd(:,:,iAlt)
pco2(:,:,iAlt) = NDensityS(:,:,iAlt,iCO2_,iBlock)/mnd(:,:,iAlt)
! Temperature Based Arrays
ttot(:,:,iAlt) = Temperature(:,:,iAlt,iBlock) * &
TempUnit(:,:,iAlt)
tt(:,:,iAlt) = ttot(:,:,iAlt)**0.69
enddo
! -------------------------------------------------------------------------------
!\
! These are needed for the Euv Heating and other thermodynamics:
!/
!write(*,*) '==> calc_rates: Before Kt, Km Calculation (after vtgcm2d).'
! -------------------------------------------------------------------------------
!---- KT, KM FORMULATION TAKEN FROM VTGCM INPUT DATASET CODE
! ( BANKS AND KOCKARTS )
!---- KM=((PO*3.9)+(PN2*3.42))*TT*1.E-06 +(PCO*COKM)+(PCO2*CO2KM)
!---- KT=((PO*75.9)+(PN2*56.))*TT +(PCO*COKT)+(PCO2*CO2KT)
! -------------------------------------------------------------------------------
do iLon = -1,nLons+2
do iLat = -1,nLats+2
do iAlt = 0, nAlts+1
! co2 factors:
is = int((ttot(iLon,iLat,iAlt)-173.3)/100.)
if (is <= 1) is = 1
if (is >= 7) is = 7
rrco2 = RGAS*AMU/Mass(iCO2_)
if (ttot(iLon,iLat,iAlt) < 500.) &
crn = 1.64-(ttot(iLon,iLat,iAlt)-500.)*2.5e-4
co2km(iLon,iLat,iAlt)=cmrf(is)+(cmrf(is+1)-cmrf(is))* &
(ttot(iLon,iLat,iAlt)- (is*100.+73.3))*0.01
co2km(iLon,iLat,iAlt)=co2km(iLon,iLat,iAlt)*1.e-06
cpco2=3.5*RGAS*AMU/Mass(iCO2_)
co2kt(iLon,iLat,iAlt)=(cpco2-rrco2)*co2km(iLon,iLat,iAlt)*crn
!
! co factors:
is=int(ttot(iLon,iLat,iAlt)/100.)
if (is <= 1) is=1
if (is >= 7) is=7
if (ttot(iLon,iLat,iAlt) > 400.) prco=0.72
if (ttot(iLon,iLat,iAlt) > 300. .and. ttot(iLon,iLat,iAlt) < 400.) &
prco = 0.73-(ttot(iLon,iLat,iAlt)-350.)*1.5e-04
if (ttot(iLon,iLat,iAlt) < 300.) &
prco=0.75-(ttot(iLon,iLat,iAlt)-250.)*2.6e-04
cokm(iLon,iLat,iAlt)=com(is)+(com(is+1)-com(is))* &
(ttot(iLon,iLat,iAlt)-100.*is)*0.01
cokm(iLon,iLat,iAlt)=cokm(iLon,iLat,iAlt)*1.65e-04
cokt(iLon,iLat,iAlt)=(3.5*RGAS* &
cokm(iLon,iLat,iAlt))/(28.*prco)
!
! Total mixture kt and km formulation : B&K (1973) from vtgcm2d code
! *******Modified the units (cgs to mks for Mars GIM code)*******
! *******Corrected units (2/20/07); S. W. Bougher *******
kmmix(iLon,iLat,iAlt) = ((pn2(iLon,iLat,iAlt)*3.42+ &
po(iLon,iLat,iAlt)*3.9)*1.0e-06*tt(iLon,iLat,iAlt) + &
pco2(iLon,iLat,iAlt)*co2km(iLon,iLat,iAlt)+ &
pco(iLon,iLat,iAlt)*cokm(iLon,iLat,iAlt))*0.10
ktmix(iLon,iLat,iAlt) = ((po(iLon,iLat,iAlt)*75.9+ &
pn2(iLon,iLat,iAlt)*56.)*tt(iLon,iLat,iAlt)+ &
pco(iLon,iLat,iAlt)*cokt(iLon,iLat,iAlt)+ &
pco2(iLon,iLat,iAlt)*co2kt(iLon,iLat,iAlt))* &
1.0E-05
enddo
enddo
enddo
! -------------------------------------------------------------------------------
if (iDebugLevel > 4) write(*,*) "=====> Before cp and kappatemp", iblock
do iAlt = 0, nAlts+1
! ---------------------------------------------------------------------------
!
! if (Is1D .and. UseKappa1DCorrection) then
! KappaTemp(:,:,iAlt,iBlock) = KappaTemp0 * &
! (Temperature(1:nLons,1:nLats,iAlt,iBlock) * &
! TempUnit(1:nLons,1:nLats,iAlt) / &
! Kappa1DCorrectionFactor)**Kappa1dCorrectionPower
! else
! KappaTemp(:,:,iAlt,iBlock) = KappaTemp0 * &
! (Temperature(1:nLons,1:nLats,iAlt,iBlock) * &
! TempUnit(1:nLons,1:nLats,iAlt))**0.75
! endif
KappaTemp(:,:,iAlt,iBlock) = ktmix(1:nLons,1:nLats,iAlt)
! ---------------------------------------------------------------------------
! This adds the eddy turbulent conduction Term (scaled by Prandtl number)
! Simplified from Yue's
! Prandtl = 10. Small to Start!
!
do iLat = 1, nLats
do iLon = 1, nLons
KappaTemp(iLon,iLat,iAlt,iBlock) = &
KappaTemp(iLon,iLat,iAlt,iBlock) + &
EddyDiffusionCoef * cp(iLon,iLat,iAlt,iBlock) * &
Rho(iLon,iLat,iAlt,iBlock)/10.
enddo
enddo
! ---------------------------------------------------------------------------
! Earth GITM formulation for Molecular Viscosity (mks)
! ViscCoef(:,:,iAlt) = 4.5e-5 * &
! (Temperature(1:nLons,1:nLats,iAlt,iBlock)*&
! TempUnit(1:nLons,1:nLats,iAlt)/ 1000.)**(-0.71)
! MTGCM formulation for Molecular Viscosity requires cgs to mks conversion
ViscCoef(:,:,iAlt) = kmmix(1:nLons,1:nLats,iAlt)
! Visc_3D(:,:,iAlt,iBlock) = kmmix(1:nLons,1:nLats,iAlt)
! ---------------------------------------------------------------------------
enddo
call end_timing("calc_rates")
if (UseBarriers) call MPI_BARRIER(iCommGITM,iError)
if (iDebugLevel > 4) write(*,*) "=====> Done with calc_rates"
end subroutine calc_rates
subroutine calc_collisions(iBlock)
use ModGITM
use ModRates
use ModConstants
use GITM_planet
use ModInputs
implicit none
integer, intent(in) :: iBlock
real, dimension(nLons, nLats, nAlts) :: Tn, Ti, e2
integer :: iError
real, dimension(-1:nLons+2, -1:nLats+2, -1:nAlts+2) :: &
Ne, mnd, Te
!\
! Need to get the neutral, ion, and electron temperature
!/
Tn = Temperature(1:nLons,1:nLats,1:nAlts,iBlock)*&
TempUnit(1:nLons,1:nLats,1:nAlts)
Ti = ITemperature(1:nLons,1:nLats,1:nAlts,iBlock)
mnd = NDensity(:,:,:,iBlock)+1.0
Ne = IDensityS(:,:,:,ie_,iBlock)
!\
! -----------------------------------------------------------
! Collision Frequencies
! -----------------------------------------------------------
!/
e_gyro = &
Element_Charge * B0(:,:,:,iMag_,iBlock) / Mass_Electron
e2 = Element_Charge * Element_Charge
!
! Ion Neutral Collision Frequency (From Kelley, 1989, pp 460):
!
if (UseBarriers) call MPI_BARRIER(iCommGITM,iError)
if (iDebugLevel > 4) write(*,*) "=====> vin",iblock
Collisions(:,:,:,iVIN_) = 2.6e-15 * (mnd + Ne)/sqrt(MeanMajorMass/AMU)
!
! Electron Neutral Collision Frequency
!
if (UseBarriers) call MPI_BARRIER(iCommGITM,iError)
if (iDebugLevel > 4) write(*,*) "=====> ven", iblock
Te = eTemperature(:,:,:,iBlock)
where(te == 0.0) te = 1000.0
Collisions(:,:,:,iVEN_) = 5.4e-16 * (mnd)*sqrt(Te)
!!!!
!!!! Electron Ion Collision Frequency
!!!!
!!!
!!! if (UseBarriers) call MPI_BARRIER(iCommGITM,iError)
!!! if (iDebugLevel > 4) write(*,*) "=====> vei", iblock
!!!
!!! tmp = (34.0 + 4.18*log((TE**3.0)/(Ne*1.0e-6)))
!!! Collisions(:,:,:,iVEI_) = tmp*Ne*TE**(-3.0/2.0) * 1.0e-6
!! Collisions(:,:,:,VEI) = (34.0 + 4.18*log((TE**3.0)/(Ne*1.0e-6))) &
!! * Ne * TE**(-3.0/2.0) * 1.0e-6
!
i_gyro = Element_Charge * B0(:,:,:,iMag_,iBlock) / MeanIonMass
end subroutine calc_collisions
subroutine calc_viscosity(iBlock)
use ModGITM
implicit none
integer, intent(in) :: iBlock
call calc_rate(iBlock)
end subroutine calc_viscosity