-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathsetup_sphericalblastND_mhd.f90
253 lines (234 loc) · 8.46 KB
/
setup_sphericalblastND_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
!------------------------------------------------------------------------------!
! 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: !
!------------------------------------------------------------------------------!
!!-------------------------------------------------------------------------!!
!! !!
!! setup for spherical adiabatic (mhd) blast waves !!
!! in 1,2 and 3 dimensions !!
!! !!
!! density is set to unity all over, whilst the pressure (or equivalently !!
!! the thermal energy) is set to some large quantity in a small circle !!
!! around the origin !!
!! !!
!! magnetic field of strength 10G in the x-direction !!
!! !! !!
!!-------------------------------------------------------------------------!!
subroutine setup
!
!--include relevant global variables
!
use dimen_mhd, only:ndim,ndimV
use debug, only:trace
use loguns, only:iprint
use bound
use eos
use kernels, only:interpolate_kernel
use options, only:damp,ibound,geom
use part
use setup_params, only:psep
use uniform_distributions
!
!--define local variables
!
implicit none
integer :: ipart
real :: denszero,przero
real :: totmass,gam1,massp
!
!--set boundaries
!
!geom = 'sphrpt'
nbpts = 0
select case(geom)
case('sphrpt')
ibound = 1 ! fixed ghosts
xmin(1) = 0.1
xmax(1) = 0.5
case default
ibound = 3 ! fixed ghosts
xmin(:) = -0.5 ! same xmin in all dimensions
xmax(:) = 0.5
xmin(2) = -0.75
xmax(2) = 0.75
end select
denszero = 1.0
przero = 0.1
gam1 = gamma - 1.
if (abs(gam1).lt.1.e-3) stop 'eos cannot be isothermal for this setup'
!
!--setup uniform density grid of particles
! (determines particle number and allocates memory)
!
call set_uniform_cartesian(2,psep,xmin,xmax,fill=.true.) ! 2 = close packed arrangement
ntotal = npart
!
!--determine particle mass
!
totmass = denszero*product(xmax(:)-xmin(:)) ! assumes cartesian boundaries
massp = totmass/float(ntotal) ! average particle mass
!
!--now assign particle properties
!
do ipart=1,ntotal
vel(:,ipart) = 0.
!--uniform density and smoothing length
dens(ipart) = denszero
pmass(ipart) = massp
uu(ipart) = przero/(gam1*denszero)
enddo
!
!--add the blast wave if damping is off
!
if (abs(damp).lt.tiny(damp)) call modify_dump
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' exiting subroutine setup'
return
end
subroutine modify_dump
use loguns, only:iprint
use part
use options, only:imhd
use timestep, only:time
use setup_params, only:pi,psep,hfact
use eos, only:gamma
use cons2prim, only:specialrelativity
implicit none
integer :: ipart
real :: prblast,pri,uui,enblast,enzero,const
real :: rblast,radius,przero,gam1,denszero
real, dimension(ndim) :: xblast, dblast
real, dimension(ndimV) :: Bzero
real :: rbuffer, exx, hsmooth
real :: q2, wab, grkern
logical, parameter :: dosedov = .false.
logical, parameter :: terryparams= .false.
logical, parameter :: delzanna_bucciantini_2002 = .true.
write(iprint,*) 'modifying dump by adding blast'
gam1 = gamma - 1.0
if (abs(gam1).lt.1.e-3) stop 'eos cannot be isothermal for this setup'
!
!--setup parameters for the problem
!
xblast(:) = 0.0 ! co-ordinates of the centre of the initial blast
bzero(:) = 0.0
const = 1./sqrt(4.*pi)
if (imhd.ne.0) then
bzero(1) = 3.0 !!sqrt(2.*pi) !10.0*const ! uniform field in bx direction
! bzero(2) = sqrt(2.*pi)
endif
if (specialrelativity) then
if (delzanna_bucciantini_2002) then
write(iprint,*) 'using setup for Del Zanna & Bucciantini special relativistic problem'
rblast = 0.4
przero = 1.0
denszero = 1.
prblast = 1000.
else
rblast = 0.04
przero = 1.0
denszero = 1.
! prblast = 1000.
enblast = 1.0
prblast = gam1*enblast/(4./3.*pi*rblast**3)
write(iprint,*) 'using setup for special relativistic problem'
endif
elseif (dosedov) then
rblast = 2.*hfact*psep ! radius of the initial blast
przero = 0.0 ! initial pressure
enblast = 1.0
enzero = 0.
prblast = gam1*enblast/(4./3.*pi*rblast**3)
!enblast = enblast/massp ! enblast is now the energy to put in a single particle
elseif (terryparams) then
rblast = 0.125 ! radius of the initial blast
przero = 1.0 ! initial pressure
prblast = 100.0 ! initial pressure within rblast
if (imhd.ne.0) then
Bzero(1) = 10. !*const
endif
else
rblast = 0.1 ! radius of the initial blast
przero = 0.1 ! initial pressure
prblast = 10.0 ! initial pressure within rblast
endif
rbuffer = rblast !+10.*psep ! radius of the smoothed front
denszero = 1.0
!
!--smoothing length for kernel smoothing
!
hsmooth = hfact*(pmass(1)/denszero)**dndim
! write(iprint,*) 'hsmooth = ',hsmooth
write(iprint,10) ndim
write(iprint,20) prblast,rblast,denszero,przero
10 format(/,1x,i1,'-dimensional adiabatic mhd blast wave problem')
20 format(/,' central pressure = ',f10.3,', blast radius = ',f6.3,/, &
' density = ',f6.3,', pressure = ',f6.3,/)
do ipart=1,ntotal
vel(:,ipart) = 0.
!--uniform density and smoothing length
dens(ipart) = denszero
dblast(:) = x(:,ipart)-xblast(:)
radius = sqrt(dot_product(dblast,dblast))
!
!--smooth energy injection using the sph kernel
!
! q2 = radius**2/hsmooth**2
! call interpolate_kernel(q2,wab,grkern)
! uui = enblast*wab/hsmooth**ndim
if (radius.le.rblast) then
pri = prblast
!uui = enblast
elseif (radius.lt.rbuffer) then ! smooth out front
exx = exp((radius-rblast)/(psep))
pri = (prblast + przero*exx)/(1.0+exx)
!uui = (enblast + enzero*exx)/(1.0+exx)
else
pri = przero
!uui = enzero
endif
uu(ipart) = pri/(gam1*denszero)
!
!--euler potentials setup (for use in other codes)
!
if (imhd.eq.-2 .and. ndim.eq.3) then
if (abs(Bzero(1)).gt.tiny(Bzero)) then ! Bx
Bevol(1,ipart) = -Bzero(1)*x(3,ipart)
Bevol(2,ipart) = x(2,ipart)
elseif (abs(Bzero(2)).gt.tiny(Bzero)) then ! By
Bevol(1,ipart) = -Bzero(2)*x(1,ipart)
Bevol(2,ipart) = x(3,ipart)
elseif (abs(Bzero(3)).gt.tiny(Bzero)) then ! Bz
Bevol(1,ipart) = -Bzero(3)*x(2,ipart)
Bevol(2,ipart) = x(1,ipart)
else
Bevol(:,ipart) = 0.
endif
endif
Bfield(:,ipart) = bzero(:)
enddo
!
!--set the constant components of the mag field which can be subtracted
!
bconst(:) = bzero(:)
time = 0.
end subroutine modify_dump