-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathevwrite_mhd.f90
320 lines (298 loc) · 10.4 KB
/
evwrite_mhd.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
!------------------------------------------------------------------------------!
! NDSPMHD: A Smoothed Particle (Magneto)Hydrodynamics code for (astrophysical) !
! fluid dynamics simulations in 1, 2 and 3 spatial dimensions. !
! !
! (c) 2002-2015 Daniel Price !
! !
! http://users.monash.edu.au/~dprice/ndspmhd !
! daniel.price@monash.edu -or- dprice@cantab.net (forwards to current address) !
! !
! NDSPMHD comes with ABSOLUTELY NO WARRANTY. !
! This is free software; and you are welcome to redistribute !
! it under the terms of the GNU General Public License !
! (see LICENSE file for details) and the provision that !
! this notice remains intact. If you modify this file, please !
! note section 2a) of the GPLv2 states that: !
! !
! a) You must cause the modified files to carry prominent notices !
! stating that you changed the files and the date of any change. !
! !
! ChangeLog: !
!------------------------------------------------------------------------------!
!!--------------------------------------------------------------------
!! Calculate conserved quantities etc and write to .ev file
!!--------------------------------------------------------------------
subroutine evwrite(t,etot,momtot)
use dimen_mhd, only:ndim,ndimV
use debug, only:trace
use loguns, only:iprint,ievfile
use derivB
use options
use part
use rates, only:force,potengrav,poten,dendt
use fmagarray
use timestep, only:dt
use utils, only:cross_product3D,minmaxave
use externf, only:external_potentials
!
!--define local variables
!
implicit none
integer :: i
real, intent(in) :: t
real, intent(out) :: etot,momtot
real :: ekin,etherm,emag,epot,dmomtot
real :: pmassi,rhoi,epoti!,rr
! real, dimension(ndim) :: rhat
real, dimension(ndimV) :: veli,mom,ang,angi,dmom
! real :: betai,alphai,alphatstarav,betatstarav
!
!--mhd
!
real, dimension(ndimB) :: Bi,Brhoi,fluxtot
real :: B2i,Bmagi
real :: fluxtotmag,crosshel
real :: betamhdi,betamhdmin,betamhdav
real :: fdotBi,fdotBmax,fdotBav
real :: forcemagi,force_erri,force_err_max,force_err_av
real :: divBi,divBav,divBmax,divBtot
real :: omegamhdi,omegamhdav,omegamhdmax
real :: fracdivBok
real, parameter :: omegtol = 1.E-2
real :: fmagabs,rhomax,rhomean,rhomin
real :: angtot,ekiny,emagp
!
!--one fluid dust
!
real :: totmassgas,totmassdust,dustfraci,dterm,ekindeltav
logical, save :: init = .false.
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' Entering subroutine evwrite'
ekin = 0.0
etherm = 0.0
emag = 0.0
emagp = 0.
etot = 0.0
if (igravity.ne.0) then
epot = potengrav
else
epot = 0.
endif
mom(:) = 0.0
dmom(:) = 0.
momtot = 0.0
ang(:) = 0.
ekiny = 0.
totmassdust = 0.
totmassgas = 0.
! alphatstarav = 0.
! betatstarav = 0.
!
!--mhd parameters
!
if (imhd.ne.0) then
betamhdav = 0.
betamhdmin = huge(betamhdmin)
divBmax = 0.
divBav = 0.
divBtot = 0.
FdotBmax = 0.
FdotBav = 0.
force_err_max = 0.
force_err_av = 0.
omegamhdav = 0.
omegamhdmax = 0.
fracdivBok = 0.
fluxtot(:) = 0.
fluxtotmag = 0.
crosshel = 0.
endif
!
!--should really recalculate the thermal energy from the total energy here
! (otherwise uu is from the half time step and same with Bfield)
!
! CALL conservative2primitive
do i=1,npart
pmassi = pmass(i)
rhoi = rho(i)
veli(:) = vel(:,i)
mom(:) = mom(:) + pmassi*veli(:)
dmom(:) = dmom(:) + pmassi*force(:,i)
if (ndim.eq.3) then
call cross_product3D(x(:,i),veli(:),angi(:))
ang(:) = ang(:) + pmassi*angi(:)
elseif (ndim.eq.2) then
ang(3) = ang(3) + pmassi*(x(1,i)*veli(2) - x(2,i)*veli(1))
endif
ekin = ekin + 0.5*pmassi*DOT_PRODUCT(veli,veli)
if (onef_dust) then
dustfraci = sum(dustfrac(:,i))
dterm = (1. - dustfraci)
ekindeltav = 0.5*pmassi*dustfraci*dterm*dot_product(deltav(:,i),deltav(:,i))
ekin = ekin + ekindeltav
ekiny = ekiny + ekindeltav
etherm = etherm + pmassi*uu(i)*dterm
totmassgas = totmassgas + pmassi*dterm
totmassdust = totmassdust + pmassi*dustfraci
else
if (ndim.ge.2) ekiny = ekiny + 0.5*pmassi*vel(1,i)*vel(1,i)
etherm = etherm + pmassi*uu(i)
endif
!
!--potential energy from external forces
!
call external_potentials(iexternal_force,x(:,i),epoti,ndim)
epot = epot + pmassi*epoti
! rr = dot_product(x(1:ndim,i),x(1:ndim,i))
! if (rr.gt.tiny(rr)) then
! rhat(1:ndim) = x(1:ndim,i)/rr
! else
! rhat = 0.
! endif
! alphai = dot_product(vel(1:ndim,i),rhat(1:ndim))
! betai = vel(2,i)*rhat(1) - vel(1,i)*rhat(2)
! alphatstarav = alphatstarav + alphai
! betatstarav = betatstarav + betai
!
!--mhd parameters
!
if (imhd.ne.0) then
Bi(:) = Bfield(:,i)
Brhoi(:) = Bi(:)/rhoi
B2i = DOT_PRODUCT(Bi,Bi)
Bmagi = SQRT(B2i)
forcemagi = SQRT(DOT_PRODUCT(force(:,i),force(:,i)))
divBi = abs(divB(i))
emag = emag + 0.5*pmassi*B2i/rhoi
emagp = emagp + 0.5*pmassi*dot_product(Bi(1:2),Bi(1:2))/rhoi
!
!--Plasma beta minimum/maximum/average
!
if (B2i.LT.tiny(B2i)) then
betamhdi = 0.
else
betamhdi = pr(i)/(0.5*B2i)
endif
betamhdav = betamhdav + betamhdi
if (betamhdi.LT.betamhdmin) betamhdmin = betamhdi
!
!--Maximum divergence of B
!
if (divBi.GT.divBmax) divBmax = divBi
divBav = divBav + divBi
!
!--volume integral of div B (int B.dS)
!
divBtot = divBtot + pmassi*divBi/rhoi
!
!--Max component of magnetic force in the direction of B (should be zero)
!
fmagabs = SQRT(DOT_PRODUCT(fmag(:,i),fmag(:,i)))
if (fmagabs.GT.1.e-8 .and. Bmagi.gt.1.e-8) then
fdotBi = ABS(DOT_PRODUCT(fmag(:,i),Bi(:)))/(fmagabs*Bmagi)
else
FdotBi = 0.
endif
fdotBav = fdotBav + fdotBi
if (fdotBi.GT.fdotBmax) fdotBmax = fdotBi
!
!--Compute total error in the force due to the B(div B) term
! only slight worry with this is that fmag is calculated in rates, whilst
! B has been evolved a bit further since then. A possible solution is to
! evaluate these quantities just after the call to rates.
!
if (forcemagi.GT.1.e-8 .AND. Bmagi.GT.1e-8) then
force_erri = ABS(DOT_PRODUCT(fmag(:,i),Bi(:)))/(forcemagi*Bmagi)
else
force_erri = 0.
endif
force_err_av = force_err_av + force_erri
if (force_erri.GT.force_err_max) force_err_max = force_erri
!
!--|div B| x smoothing length / |B| (see e.g. Cerqueira and Gouveia del Pino 1999)
! this quantity should be less than ~0.01.
!
if (Bmagi.lt.1e-8) then
omegamhdi = 0.
else
omegamhdi = divBi*hh(i)/Bmagi
endif
if (omegamhdi.LT.omegtol) fracdivBok = fracdivBok + 1.
if (omegamhdi.GT.omegamhdmax) omegamhdmax = omegamhdi
omegamhdav = omegamhdav + omegamhdi
!
!--Conserved magnetic flux (int B dV)
!
pmassi = pmass(i)
fluxtot(:) = fluxtot(:) + pmassi*Brhoi(:)
!
!--Conserved Cross Helicity (int v.B dV)
!
crosshel = crosshel + pmassi*DOT_PRODUCT(veli,Brhoi)
endif
enddo
etot = ekin + emag + epot
if (iprterm.ge.0 .or. iprterm.lt.-1) etot = etot + etherm
momtot = sqrt(dot_product(mom,mom))
dmomtot = sqrt(dot_product(dmom,dmom))
!print*,'dpdt = ',dmomtot
call minmaxave(rho(1:npart),rhomin,rhomax,rhomean,npart)
angtot = sqrt(dot_product(ang,ang))
!print*,' etot = ',etot,' sum of en(i) = ',sum(pmass(1:npart)*en(1:npart)),' potengrav=',potengrav,&
! sum(pmass(1:npart)*0.5*poten(1:npart))
!print*,'energy conservation=',sum(pmass(1:npart)*dendt(1:npart))
!
!--write line to .ev file
!
if (imhd.ne.0) then
fluxtotmag = SQRT(DOT_PRODUCT(fluxtot,fluxtot))
betamhdav = betamhdav/FLOAT(npart)
fracdivBok = 100.*fracdivBok/FLOAT(npart)
omegamhdav = omegamhdav/FLOAT(npart)
divBav = divBav/FLOAT(npart)
fdotBav = fdotBav/FLOAT(npart)
force_err_av = force_err_av/FLOAT(npart)
!! print*,'t=',t,' emag =',emag,' etot = ',etot, 'ekin = ',ekin,' etherm = ',etherm
if (.not.init) then
write(ievfile,"('#',24(a18,1x))") 't','ekin','etherm','emag','epot','etot','momtot',&
'angtot','rhomax','rhomean','dt','emagp','crosshel','betamhdmin','betamhdav',&
'divBav','divBmax','divBtot','fdotBav','fdotBmax','force_err_av','force_err_max',&
'omegamhdav','omegamhdmax','fracdivBok'
init = .true.
endif
write(ievfile,30) t,ekin,etherm,emag,epot,etot,momtot,angtot,rhomax,rhomean,dt, &
emagp,crosshel,betamhdmin,betamhdav, &
divBav,divBmax,divBtot, &
fdotBav,FdotBmax,force_err_av,force_err_max, &
omegamhdav,omegamhdmax,fracdivBok
30 format(24(1pe18.10,1x),1pe8.2)
else
!alphatstarav = alphatstarav/float(npart)
!betatstarav = betatstarav/float(npart)
!! print*,'t=',t,' emag =',emag,' etot = ',etot, 'ekin = ',ekin,' etherm = ',etherm
if (idust.eq.1 .or. idust.eq.3 .or. idust.eq.4) then
if (.not.init) then
write(ievfile,"('#',24(a18,1x))") 't','ekin','etherm','emag','epot','etot','momtot',&
'angtot','rhomax','rhomean','dt','ekiny','totmassgas','totmassdust'
init = .true.
endif
write(ievfile,40) t,ekin,etherm,emag,epot,etot,momtot,angtot,rhomax,rhomean,dt,ekiny,&
totmassgas,totmassdust
else
if (.not.init) then
write(ievfile,"('#',24(a18,1x))") 't','ekin','etherm','emag','epot','etot','momtot',&
'angtot','rhomax','rhomean','dt','ekiny','dmomtot'
init = .true.
endif
write(ievfile,40) t,ekin,etherm,emag,epot,etot,momtot,angtot,rhomax,rhomean,dt,ekiny,dmomtot
endif
40 format(24(1pe18.10,1x))
endif
!
!--flush the buffer so that the line is written to the file immediately
!
call flush(ievfile)
return
end subroutine evwrite