-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathsetup_wave_x_ND_dust.f90
278 lines (257 loc) · 7.59 KB
/
setup_wave_x_ND_dust.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
!----------------------------------------------------------------
! Set up a travelling soundwave in the x direction
! perturbs particles from a uniform density grid
! should work in 1, 2 and 3 dimensions
!----------------------------------------------------------------
subroutine setup
!
!--include relevant global variables
!
use dimen_mhd
use debug
use loguns
use bound
use eos
use options
use part
use setup_params
use uniform_distributions, only:set_uniform_cartesian
!
!--define local variables
!
implicit none
integer :: i
integer, parameter :: itsmax = 100
integer :: ntypes
real, parameter :: tol = 1.e-8
integer :: its,iwave,ngas,npdust,jtype
real, dimension(ndimV) :: Bzero
real :: massp,masspdust
real :: ampl,wk,xlambda,dxmax,denom
real :: dxi,dxprev,xmassfrac,func,fderiv
real :: spsoundi,valfven2i,vamplx,vamply,vamplz
real :: vfast,vslow,vcrap,vwave,term,dens1
real :: denszero,uuzero,przero,Rzero,denszerodust
real :: dust_to_gas_ratio,mgastot
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' Entering subroutine setup'
10 format(/,'-------------- ',a,' ----------------')
iwave = 0 ! preset wave parameters
!
!--setup parameters (could read in from a file)
!
Bzero = 0.
if (iwave.EQ.1) then ! MHD slow wave
ampl = 0.006
denszero = 1.0
uuzero = 4.5
if (imhd.ne.0) Bzero(1:3) = sqrt(2.)
write(iprint,10) 'MHD slow wave'
elseif (iwave.EQ.2) then ! MHD fast wave
ampl = 0.0055
denszero = 1.0
uuzero = 0.3
if (imhd.ne.0) Bzero(1:3) = 0.5
write(iprint,10) 'MHD fast wave'
else ! Sound wave
!ampl = 0.1
!ampl = 1.e-4
ampl = 1.e-6
denszero = 1.0
if (abs(gamma-1.).gt.1e-3) then
uuzero = 1.0/((gamma-1.)*gamma)
else
uuzero = 1.0
endif
Rzero = -1.0 ! negative stress parameter
open(unit=20,err=30,file=trim(rootname)//'.rstress',status='old')
read(20,*) Rzero
close(unit=20)
print*,'Read stress file: R parameter = ',Rzero
30 continue
Bzero(1) = sqrt(2.*(1.-Rzero))
write(iprint,10) 'sound wave'
endif
xlambda = 1.0
wk = 2.0*pi/xlambda ! wave number
!
!--set boundaries
!
ibound = 3 ! periodic boundaries
nbpts = 0 ! use ghosts not fixed
xmin(:) = 0. ! set position of boundaries
xmax(1) = 1.0
if (ndim.GE.2) then
xmax(2:ndim) = 11.*psep ! would need to adjust this depending on grid setup
endif
!
!--initially set up a uniform density grid (also determines npart)
! (the call to set_uniform_cartesian means this works in 1,2 and 3D)
!
if (onef_dust) then
ntypes = 1
else
ntypes = 2
endif
ngas = 0
npdust = 0
do jtype=1,ntypes
call set_uniform_cartesian(2,psep,xmin,xmax,adjustbound=.true.)
if (jtype.eq.1) then
ngas = npart
itype(1:ngas) = itypegas
elseif (jtype.eq.2) then
npdust = npart - ngas
itype(ngas+1:ngas+npdust) = itypedust
endif
enddo
! npart = INT((xmax(1)-xmin(1))/psep)
mgastot = denszero*product(xmax-xmin)
massp = mgastot/real(ngas) ! average particle mass
masspdust = 0.
dust_to_gas_ratio = 1.
if (npdust.gt.0) masspdust = dust_to_gas_ratio*mgastot/real(npdust) ! average particle mass
denszerodust = dust_to_gas_ratio*denszero
if (ntypes.gt.1) print*,' ngas = ',ngas,' npdust = ',npdust
!
!--allocate memory here
!
! call alloc(npart)
!
!--setup uniform density grid of particles
!
do i=1,npart
! x(1,i) = xmin(1) + (i-1)*psep + 0.5*psep
vel(:,i) = 0.
if (itype(i).eq.itypedust) then
dens(i) = denszerodust
pmass(i) = masspdust
uu(i) = 0.
else
dens(i) = denszero
pmass(i) = massp
uu(i) = uuzero
endif
if (imhd.GT.0) then
Bfield(:,i) = Bzero
else
Bfield(:,i) = 0.
endif
if (onef_dust) then
dustfrac(1,i) = dust_to_gas_ratio/(1. + dust_to_gas_ratio)
deltav(:,i) = 0.
pmass(i) = pmass(i)/(1. - dustfrac(1,i))
endif
enddo
ntotal = npart
!
!--get sound speed from equation of state (want average sound speed, so
! before the density is perturbed)
!
call equation_of_state1(przero,spsoundi,uuzero,denszero)
print*,' gamma = ',gamma
print*,' pr = ',przero,' cs = ',spsoundi,' u = ',uu(1),' dens = ',dens(1)
!
!--work out MHD wave speeds
!
dens1 = 1./denszero
valfven2i = dot_product(Bzero,Bzero)*dens1
vfast = sqrt(0.5*(spsoundi**2 + valfven2i &
+ sqrt((spsoundi**2 + valfven2i)**2 &
- 4.*(spsoundi*Bzero(1))**2*dens1)))
vslow = sqrt(0.5*(spsoundi**2 + valfven2i &
- sqrt((spsoundi**2 + valfven2i)**2 &
- 4.*(spsoundi*Bzero(1))**2*dens1)))
vcrap = sqrt(spsoundi**2 + valfven2i)
!----------------------------
! set the wave speed to use
!----------------------------
if (iwave.EQ.1) then
vwave = vslow
if (imhd.le.0) write(iprint,*) 'Error: can''t have slow wave if no mhd'
elseif (iwave.EQ.2) then
vwave = vfast
else
vwave = spsoundi
endif
if (vwave.le.0.) then
write(iprint,*) 'Error in setup: vwave = ',vwave
stop
endif
dxmax = xmax(1) - xmin(1)
denom = dxmax - ampl/wk*(COS(wk*dxmax)-1.0)
write (iprint,*) 'Wave number = ',wk
write (iprint,*) 'Amplitude = ',ampl
!
!--now perturb the particles appropriately
!
do i=1,npart
dxi = x(1,i)-xmin(1)
dxprev = dxmax*2.
xmassfrac = dxi/dxmax ! current mass fraction(for uniform density)
! determines where particle should be
!
!--Use rootfinder on the integrated density perturbation
! to find the new position of the particle
!
its = 0
do while ((abs(dxi-dxprev).gt.tol).and.(its.lt.itsmax))
dxprev = dxi
func = xmassfrac*denom - (dxi - ampl/wk*(COS(wk*dxi)-1.0))
fderiv = -1.0 - ampl*SIN(wk*dxi)
dxi = dxi - func/fderiv ! Newton-Raphson iteration
its = its + 1
! PRINT*,'iteration',its,'dxi =',dxi
enddo
if (its.GE.itsmax) then
write(iprint,*) 'Error: soundwave - too many iterations'
call quit
endif
! if (idebug(1:5).EQ.'sound') then
! write(*,99002) i,its,dxi-dxprev,x(1,i),xmin(1)+dxi
!99002 FORMAT('Particle',i5,' converged in ',i4, &
! ' iterations, error in x =',1(1pe10.2,1x),/, &
! 'previous x = ',1(0pf8.5,1x), &
! 'moved to x = ',1(0pf8.5,1x))
! endif
x(1,i) = xmin(1)+dxi
!
!--multiply by appropriate wave speed
!
term = vwave**2 - Bfield(1,i)**2/dens(i)
vamplx = vwave*ampl
vamply = -vamplx*Bfield(1,i)*Bfield(2,i)/(dens(i)*term)
vamplz = -vamplx*Bfield(1,i)*Bfield(3,i)/(dens(i)*term)
vel(1,i) = vamplx*SIN(wk*dxi)
! vel(2,i) = vamply*SIN(wk*dxi)
! vel(3,i) = vamplz*SIN(wk*dxi)
!
!--perturb internal energy if not using a polytropic equation of state
! (do this before density is perturbed)
!
if (itype(i).eq.itypegas) then
uu(i) = uu(i) + przero/dens(i)*ampl*SIN(wk*dxi) ! if not polytropic
endif
!
!--perturb density if not using summation
!
! Bfield(2,i) = Bfield(2,i) + vwave*Bfield(2,i)*vamplx/term*SIN(wk*dxi)
! Bfield(3,i) = Bfield(3,i) + vwave*Bfield(3,i)*vamplx/term*SIN(wk*dxi)
dens(i) = dens(i)*(1.+ampl*SIN(wk*dxi))
enddo
write(iprint,*) ' Wave set: Amplitude = ',ampl,' Wavelength = ',xlambda,' k = ',wk
write(iprint,*) ' sound speed = ',spsoundi
write(iprint,*) ' alfven speed = ',sqrt(valfven2i)
write(iprint,*) ' fast speed = ',vfast
write(iprint,*) ' slow speed = ',vslow
write(iprint,*) ' wave speed = ',vwave,iwave
return
end subroutine setup
!
! use this routine to modify the dump upon code restart
!
subroutine modify_dump()
implicit none
end subroutine modify_dump