-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathsetup_sedimentND.f90
215 lines (201 loc) · 6.4 KB
/
setup_sedimentND.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
!------------------------------------------------------------------------------!
! 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: !
!------------------------------------------------------------------------------!
!----------------------------------------------------------------
! Set up the sedimenting dust problem from Monaghan (1997b)
!----------------------------------------------------------------
subroutine setup
!
!--include relevant global variables
!
use dimen_mhd
use debug
use loguns
use bound
use options
use part
use setup_params
use externf, only:external_forces
use cons2prim, only:primitive2conservative
use uniform_distributions
use eos, only:polyk,gamma
!
!--define local variables
!
implicit none
integer :: i,ngas,npdust
real :: massp,volume,totmass,oldtolh,rhodust,masspdust,voldust
real :: denszero,fext(3),cs,gx,xmindust(ndim),xmaxdust(ndim)
real :: vdust,vgas,epsi,ts
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' entering subroutine setup(sediment)'
!
!--set boundaries
!
print*,'Setup for sedimenting dust layer '
if (iexternal_force.ne.9) stop 'need iexternal_force = 9 for this setup'
ibound(1) = 3
! boundaries
ibound(2) = 0 ! boundaries
nbpts = 0 ! use ghosts not fixed
xmin(1) = 0.
xmax(1) = 0.25
xmin(2) = 0. ! set position of boundaries
xmax(2) = 1.
xmin(2) = xmin(2) - 8.*psep
xmax(2) = xmax(2) + 8.*psep
call set_uniform_cartesian(1,psep,xmin,xmax,adjustbound=.true.)
npart = ntotal
!xmin(2) = 0.
!xmax(2) = 1.
print*,'npart =',npart
!
!--determine particle mass
!
denszero = 1.0
volume = product(xmax(:)-xmin(:))
totmass = denszero*volume
massp = totmass/float(ntotal) ! average particle mass
!
!--now assign particle properties
!
cs = 10.
polyk = cs*cs
!
!--get the actual value of the gravity force from the external force routine
!
call external_forces(9,x(:,1),fext,ndim,ndimV,vel(:,1),hh(1),cs,itypegas)
gx = -fext(2)
itype = itypegas
if (iener >= 1) then
if (gamma < 1.001) stop 'need gamma > 1 if non-isothermal'
endif
do i=1,ntotal
vel(:,i) = 0.
if (iener <= 0) then
dens(i) = denszero*exp(-gx*x(2,i)/cs**2)
pmass(i) = massp*exp(-gx*x(2,i)/cs**2)
uu(i) = 1.5*polyk ! non-isothermal
else
dens(i) = denszero
pmass(i) = massp
uu(i) = 1.5*polyk - gx*x(2,i)/(gamma - 1.) ! non-isothermal
endif
Bfield(:,i) = 0.
enddo
!
!--calculate density assuming free boundaries
! use very small htol to get it right
!
oldtolh = tolh
tolh = 1.e-12
call primitive2conservative
do i=1,ntotal
dens(i) = pmass(i)/(hh(i)/hfact)**ndim
enddo
tolh = oldtolh
!
!--now fix boundary particles
!
ibound(2) = 1
do i=1,ntotal
if (x(2,i) < 0. .or. x(2,i) > 1.) then
itype(i) = itypebnd
nbpts = nbpts + 1
endif
enddo
!
!--add dust particles
!
if (idust > 0) then
xmindust = xmin
xmaxdust = xmax
xmindust(2) = 0.6
xmaxdust(2) = 0.8
rhodust = 0.1 ! in Monaghan 97 this is rhod*thetad = 1000*0.001 = 1.0
print*,' rhodust = ',rhodust,' rhogas = ',denszero
voldust = product(xmaxdust-xmindust)
print*,' assuming Kdrag = ',Kdrag, ' gives vdust = ',-rhodust*gx/Kdrag
select case(idust)
case(2)
ngas = npart
call set_uniform_cartesian(1,psep,xmindust,xmaxdust,fill=.true.)
npart = ntotal
npdust = npart - ngas
masspdust = rhodust*voldust/real(npdust)
do i=ngas+1,npart
itype(i) = itypedust
pmass(i) = masspdust
dens(i) = rhodust
uu(i) = 0.
Bfield(:,i) = 0.
vel(:,i) = 0.
vel(2,i) = -rhodust*gx/Kdrag
enddo
case default
npdust = 0
do i=1,npart
if (x(2,i) >= xmindust(2) .and. x(2,i) <= xmaxdust(2)) then
npdust = npdust + 1
endif
enddo
masspdust = rhodust*voldust/real(npdust)
do i=1,npart
if (x(2,i) >= xmindust(2) .and. x(2,i) <= xmaxdust(2)) then
epsi = rhodust/(rhodust + dens(i))
vdust = -rhodust*gx/Kdrag
vgas = 0.
dustfrac(:,i) = 0.
dustfrac(1,i) = epsi
deltav(:,i) = 0.
deltav(2,i) = vdust - vgas
vel(2,i) = epsi*vdust + (1. - epsi)*vgas
!print*,' vy = ',vel(2,i)
!dens(i) = dens(i) !+ rhodust
pmass(i) = pmass(i) + masspdust
else
dustfrac(:,i) = 0.
deltav(:,i) = 0.
endif
enddo
end select
ts = (rhodust*denszero)/(Kdrag*(rhodust + denszero))
print*,' stopping time ts = ',ts
print*,' resolution criterion is h < ',cs*ts
if (idust.eq.2 .and. psep > cs*ts) then
print*,' warning: this is not being met '
read*
endif
endif
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' exiting subroutine setup'
return
end
!
! use this routine to modify the dump upon code restart
!
subroutine modify_dump()
implicit none
end subroutine modify_dump