-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathdust_diffusion.f90
250 lines (233 loc) · 9.34 KB
/
dust_diffusion.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
!------------------------------------------------------------------------------!
! 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: !
!------------------------------------------------------------------------------!
module dustdiffusion
implicit none
contains
!!------------------------------------------------------------------------
!! Computes dust diffusion term for one fluid dust in the terminal
!! velocity approximation (assumes deltav is known)
!!
!! particle quantities are explicitly passed rather than through the
!! global modules
!!
!! assumes there has been a previous call to density
!!------------------------------------------------------------------------
subroutine dust_diffusion(npart,ntot,x,pmass,rho,hh,gradh,dustfrac,ddustevoldt,deltav,vel,pr,uu,spsound,dudt)
use dimen_mhd, only:ndim, ndimV
use debug, only:trace
use loguns, only:iprint
use bound, only:ireal
use kernels, only:interpolate_kernel,radkern2
use linklist
use options
use dust, only:get_tstop
use get_neighbour_lists
use part, only:ndust,rhodust,rhogas
integer, intent(in) :: npart,ntot
real, dimension(ndim,ntot), intent(in) :: x
real, dimension(ntot), intent(in) :: pmass,rho,hh,gradh,uu,spsound,pr
real, dimension(ndust,ntot), intent(in) :: dustfrac
real, dimension(ndimV,ntot), intent(in) :: vel
real, dimension(ndimV,ntot), intent(inout) :: deltav
real, dimension(ndust,ntot), intent(out) :: ddustevoldt
real, dimension(ntot), intent(inout) :: dudt
!
!--define local variables
!
integer :: i,j,n
integer :: icell,iprev,nneigh
integer, dimension(ntot) :: listneigh ! neighbour list
integer :: idone
!
! (particle properties - local copies)
!
real :: rij,rij2
real :: hi,hi1,hj,hj1
real :: hfacwabi,hfacwabj
real :: pmassi,pmassj,rhoi,rhoj,gradhi,rho1i,rho1j,rhogasi,rhogasj
real :: dustfraci,dustfracj,projdeltavi,projdeltavj,dvgasdotr
real :: termi,termj,term,du,disstermi,disstermj,deltav2i,deltav2j
real :: tstopi,tstopj,pdvtermi,pdvtermj,si,sj,rhodusti,rhodustj
real, dimension(ndim) :: dx
real, dimension(ndimv) :: dr,deltavi,deltavj,dvgas,vgasi,vgasj
!
! (kernel quantities)
!
real :: q2i,q2j
real :: wabi,wabj
real :: grkerni,grkernj
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' entering subroutine dust_diffusion'
!
!--initialise quantities
! (do NOT set du/dt to zero)
!
!ddustevoldt = 0.
listneigh = 0
!
! make sure deltav has been copied to ghosts
!
do i=npart+1,ntot
j = ireal(i)
deltav(:,i) = deltav(:,j)
enddo
!
!--loop over all the link-list cells
!
loop_over_cells: do icell=1,ncellsloop ! step through all cells
!
!--get the list of neighbours for this cell
! (common to all particles in the cell)
!
call get_neighbour_list(icell,listneigh,nneigh)
!
!--now loop over all particles in the current cell
!
i = ifirstincell(icell)
idone = -1 ! note density summation includes current particle
if (i.ne.-1) iprev = i
loop_over_cell_particles: do while (i.ne.-1) ! loop over home cell particles
idone = idone + 1
hi = hh(i)
hi1 = 1./hi
hfacwabi = hi1**ndim
pmassi = pmass(i)
rhoi = rho(i)
rho1i = 1./rhoi
gradhi = gradh(i)
dustfraci = dustfrac(1,i)
if (use_smoothed_rhodust) then
rhogasi = rhogas(i)
rhodusti = rhodust(1,i)
else
rhogasi = (1. - dustfraci)*rhoi
rhodusti = rhoi - rhogasi
endif
deltavi = deltav(:,i)
deltav2i = dot_product(deltavi,deltavi)
vgasi = vel(:,i) !- dustfraci*deltavi
tstopi = get_tstop(idrag_nature,rhogasi,rhodusti,spsound(i),Kdrag)
!tstopi = dustfraci*rhogasi/Kdrag
!
!--for each particle in the current cell, loop over its neighbours
!
loop_over_neighbours: do n = idone+1,nneigh
j = listneigh(n)
if ((j.ne.i) .and..not.(j.gt.npart .and. i.gt.npart)) then
! don't count particle with itself
dx(:) = x(:,i) - x(:,j)
hj = hh(j)
hj1 = 1./hj
hfacwabj = hj1**ndim
rij2 = dot_product(dx,dx)
q2i = rij2*hi1*hi1
q2j = rij2*hj1*hj1
!
!--do interaction if r/h < compact support size
!
if ((q2i.lt.radkern2).or.(q2j.lt.radkern2)) then
rij = sqrt(rij2)
pmassj = pmass(j)
rhoj = rho(j)
rho1j = 1./rhoj
dustfracj = dustfrac(1,j)
if (use_smoothed_rhodust) then
rhogasj = rhogas(j)
rhodustj = rhodust(1,j)
else
rhogasj = (1. - dustfracj)*rhoj
rhodustj = rhoj - rhogasj
endif
deltavj = deltav(:,j)
deltav2j = dot_product(deltavj,deltavj)
tstopj = get_tstop(idrag_nature,rhogasj,rhodustj,spsound(j),Kdrag)
!tstopj = dustfracj*rhogasj/Kdrag
dr = 0.
dr(1:ndim) = dx(1:ndim)/(rij + epsilon(rij)) ! unit vector
!
!--interpolate from kernel table
!
! (using hi)
call interpolate_kernel(q2i,wabi,grkerni)
wabi = wabi*hfacwabi
grkerni = grkerni*hfacwabi*hi1*gradhi
! (using hj)
call interpolate_kernel(q2j,wabj,grkernj)
wabj = wabj*hfacwabj
grkernj = grkernj*hfacwabj*hj1*gradh(j)
!
!--calculate depsilon/dt
!
projdeltavi = dot_product(deltav(:,i),dr)
projdeltavj = dot_product(deltav(:,j),dr)
select case(idustevol)
case(2)
si = sqrt(rhodusti*rhogasi)
sj = sqrt(rhodustj*rhogasj)
termi = rho1i**3*projdeltavi*grkerni
termj = rho1j**3*projdeltavj*grkernj
case(1)
si = sqrt(rhoi*dustfraci)
sj = sqrt(rhoj*dustfracj)
termi = (1. - dustfraci)*rho1i*rho1i*projdeltavi*grkerni
termj = (1. - dustfracj)*rho1j*rho1j*projdeltavj*grkernj
case default
termi = dustfraci*(1. - dustfraci)*rho1i*projdeltavi*grkerni
termj = dustfracj*(1. - dustfracj)*rho1j*projdeltavj*grkernj
si = 1.
sj = 1.
end select
term = termi + termj
ddustevoldt(1,i) = ddustevoldt(1,i) - pmassj*sj*term
ddustevoldt(1,j) = ddustevoldt(1,j) + pmassi*si*term
if (iener.gt.0) then
disstermi = 0.!dustfraci*deltav2i/tstopi
disstermj = 0.!dustfracj*deltav2j/tstopj
vgasj = vel(:,j) !- dustfracj*deltavj
dvgas = vgasi - vgasj
dvgasdotr = dot_product(dvgas,dr)
pdvtermi = pr(i)*rho1i/rhogasi*pmassj*dvgasdotr*grkerni
pdvtermj = pr(j)*rho1j/rhogasj*pmassi*dvgasdotr*grkernj
if (iener.ne.2) stop 'only thermal energy equation implemented for idust=3'
du = uu(i) - uu(j)
! dudt(i) = dudt(i) - pmassj*termi/(1. - dustfraci)*du
! dudt(j) = dudt(j) - pmassi*termj/(1. - dustfracj)*du
dudt(i) = dudt(i) - dustfraci*rho1i*pmassj*du*projdeltavi*grkerni + pdvtermi + disstermi
dudt(j) = dudt(j) - dustfracj*rho1j*pmassi*du*projdeltavj*grkernj + pdvtermj + disstermj
endif
endif
endif! .not. j>npart .and. i>npart
enddo loop_over_neighbours
iprev = i
if (iprev.ne.-1) i = ll(i)
enddo loop_over_cell_particles
enddo loop_over_cells
do i=npart+1,ntot
j = ireal(i)
ddustevoldt(:,i) = ddustevoldt(:,j)
dudt(i) = dudt(j)
enddo
return
end subroutine dust_diffusion
end module dustdiffusion