-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathdensity_sums.f90
661 lines (614 loc) · 23.5 KB
/
density_sums.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
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
!------------------------------------------------------------------------------!
! 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 density_summations
use get_neighbour_lists
implicit none
contains
!!------------------------------------------------------------------------
!! Computes the density by direct summation over the particles neighbours
!! ie. rho_a = sum_b m_b W_ab (h_a)
!!
!! Also computes the variable smoothing length terms sum_b m_b dW_ab/dh_a
!!
!! This version computes the density on all particles
!! and therefore only does each pairwise interaction once
!!------------------------------------------------------------------------
subroutine density(x,pmass,hh,vel,rho,drhodt,densn,dndt,delsqn, &
gradh,gradhn,gradsoft,gradgradh,npart,ntotal)
use dimen_mhd, only:ndim,ndimV
use debug, only:trace
use loguns, only:iprint
use kernels, only:radkern2,interpolate_kernel,interpolate_kernels_dens,interpolate_kernel_soft
use linklist, only:ll,ifirstincell,numneigh,ncellsloop
use options, only:ikernav,igravity,imhd,ikernel,ikernelalt,iprterm,onef_dust
use part, only:Bfield,uu,psi,itype,itypebnd,itypebnddust,itypegas,itypedust,dustfrac,rhogas,rhodust,ndust
use setup_params, only:hfact
use rates, only:dBevoldt
use matrixcorr, only:dxdx,ndxdx,idxdx,jdxdx
!
!--define local variables
!
implicit none
integer, intent(in) :: npart,ntotal
real, dimension(:,:), intent(in) :: x, vel
real, dimension(:), intent(in) :: pmass, hh
real, dimension(:), intent(out) :: rho,drhodt,densn,dndt,delsqn,gradh,gradhn,gradsoft,gradgradh
integer :: i,j,k,n
integer :: icell,iprev,nneigh
integer, dimension(ntotal) :: listneigh
integer :: idone
integer, parameter :: itemp = 121
!
! (particle properties - local copies)
!
real :: rij,rij2
real :: hi,hi1,hav,hav1,hj,hj1,hi21
real :: hfacwab,hfacwabi,hfacwabj
real, dimension(ndim) :: dx,xi
real, dimension(ndimV) :: veli,dvel
real, dimension(ntotal) :: rhoin
real, dimension(ntotal) :: h1,unity
real :: dvdotr,pmassi,pmassj,projBi,projBj,dustfraci,dustfracj
real, dimension(ndxdx) :: dxdxi
real, dimension(ndimV) :: dr
integer :: itypei,itypej
!
! (kernel quantities)
!
real :: q2,q2i,q2j
real :: wab,wabi,wabj,weight,wabalti,wabaltj
real :: grkern,grkerni,grkernj,grkernalti,grkernaltj,grgrkerni,grgrkernj
real :: dwdhi,dwdhj,dwaltdhi,dwaltdhj,dphidhi,dphidhj,dwdhdhi,dwdhdhj ! grad h terms
real :: wconst
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' Entering subroutine density'
!
!--initialise quantities
!
listneigh = 0
numneigh = 0
dwdhi = 0.
dwdhj = 0.
dwaltdhi = 0.
dwaltdhj = 0.
dwdhdhi = 0.
dwdhdhj = 0.
wconst = 1./hfact**ndim
dr(:) = 0.
dustfraci = 0.
do i=1,npart
if (itype(i) /= itypebnd .and. itype(i) /= itypebnddust) then
rhoin(i) = rho(i)
rho(i) = 0.
drhodt(i) = 0.
densn(i) = 0.
dndt(i) = 0.
delsqn(i) = 0.
gradh(i) = 0.
gradhn(i) = 0.
gradsoft(i) = 0.
gradgradh(i) = 0.
if (imhd.eq.5) dBevoldt(:,i) = 0.
if (imhd.eq.0) then
psi(i) = 0.
unity(i) = 0.
endif
dxdx(:,i) = 0.
if (onef_dust) then
rhogas(i) = 0.
rhodust(:,i) = 0.
endif
endif
enddo
do i=1,ntotal
h1(i) = 1./hh(i)
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
! PRINT*,'Doing particle ',i,nneigh,' neighbours',hh(i)
idone = idone + 1
pmassi = pmass(i)
if (onef_dust) dustfraci = sum(dustfrac(:,i))
xi(:) = x(:,i)
veli(:) = vel(:,i)
hi = hh(i)
hi1 = h1(i)
hfacwabi = hi1**ndim
hi21 = hi1*hi1
dxdxi(:) = 0.
itypei = itype(i)
!
!--for each particle in the current cell, loop over its neighbours
!
loop_over_neighbours: do n = idone+1,nneigh
j = listneigh(n)
!--skip particles of different type
itypej = itype(j)
if (.not. ((itypej.eq.itypei) .or. (itypei.eq.itypegas .and. itypej.eq.itypebnd) &
.or. (itypej.eq.itypegas .and. itypei.eq.itypebnd) &
.or. (itypei.eq.itypedust .and. itypej.eq.itypebnddust) &
.or. (itypej.eq.itypedust .and. itypei.eq.itypebnddust))) then
cycle loop_over_neighbours
endif
dx(:) = xi(:) - x(:,j)
hj = hh(j)
hj1 = h1(j)
rij2 = dot_product(dx,dx)
q2i = rij2*hi21
q2j = rij2*hj1*hj1
! PRINT*,' neighbour,r/h,dx,hi,hj ',j,SQRT(q2i),dx,hi,hj
!
!--do interaction if r/h < compact support size
! don't calculate interactions between ghost particles
!
if ((q2i.LT.radkern2).OR.(q2j.LT.radkern2) &
.AND. (i.LE.npart .OR. j.LE.npart)) then
!if (i.eq.itemp .or. j.eq.itemp) then
! print*,' neighbour,r/hi,r/hj,hi,hj:',i,j,sqrt(q2i),sqrt(q2j),hi,hj,rho(itemp)
!endif
if (i.LE.npart) numneigh(i) = numneigh(i) + 1
if (j.LE.npart .and. j.ne.i) numneigh(j) = numneigh(j) + 1
rij = sqrt(rij2)
dr(1:ndim) = dx(1:ndim)/(rij + epsilon(rij))
hfacwabj = hj1**ndim
!
!--weight self contribution by 1/2
!
if (j.EQ.i) then
weight = 0.5
else
weight = 1.0
endif
pmassj = pmass(j)
!
!--interpolate from kernel table
! (use either average h or average kernel gradient)
!
if (ikernav.EQ.1) then
! (using average h)
hav = 0.5*(hi + hj)
hav1 = 1./hav
hfacwab = hav1**ndim
q2 = rij2*hav1*hav1
call interpolate_kernel(q2,wab,grkern)
wab = wab*hfacwab
grkern = grkern*hfacwab*hav1
wabi = wab
wabj = wab
grkerni = grkern
grkernj = grkern
else
if (igravity.ne.0 .and. ikernav.eq.3 .and. ndim.eq.3) then
call interpolate_kernel_soft(q2i,wabi,grkerni,dphidhi)
gradsoft(i) = gradsoft(i) + weight*pmassj*dphidhi*hi21
call interpolate_kernel_soft(q2j,wabj,grkernj,dphidhj)
gradsoft(j) = gradsoft(j) + weight*pmassi*dphidhj*hj1*hj1
endif
call interpolate_kernels_dens(q2i,wabi,grkerni,grgrkerni,wabalti,grkernalti)
call interpolate_kernels_dens(q2j,wabj,grkernj,grgrkernj,wabaltj,grkernaltj)
! (using hi)
wabi = wabi*hfacwabi
wabalti = wabalti*hfacwabi
grkerni = grkerni*hfacwabi*hi1
grgrkerni = grgrkerni*hfacwabi*hi1*hi1
grkernalti = grkernalti*hfacwabi*hi1
! (using hj)
wabj = wabj*hfacwabj
wabaltj = wabaltj*hfacwabj
grkernj = grkernj*hfacwabj*hj1
grgrkernj = grgrkernj*hfacwabj*hj1*hj1
grkernaltj = grkernaltj*hfacwabj*hj1
! (calculate average)
if (ikernav.eq.2) then
wab = 0.5*(wabi + wabj)
wabi = wab
wabj = wab
grkern = 0.5*(grkerni + grkernj)
grkerni = grkern
grkernj = grkern
endif
!
!--derivative w.r.t. h for grad h correction terms (and dhdrho)
!
dwdhi = -rij*grkerni*hi1 - ndim*wabi*hi1
dwdhj = -rij*grkernj*hj1 - ndim*wabj*hj1
dwaltdhi = -rij*grkernalti*hi1 - ndim*wabalti*hi1
dwaltdhj = -rij*grkernaltj*hj1 - ndim*wabaltj*hj1
dwdhdhi = ndim*(ndim+1)*wabi*hi1**2 + 2.*(ndim+1)*rij*hi1**2*grkerni &
+ rij**2*hi1**2*grgrkerni
dwdhdhj = ndim*(ndim+1)*wabj*hj1**2 + 2.*(ndim+1)*rij*hj1**2*grkernj &
+ rij**2*hj1**2*grgrkernj
endif
!
!--calculate density and number density
!
if (itypei /= itypebnd) then
rho(i) = rho(i) + pmassj*wabi*weight
densn(i) = densn(i) + wabalti*weight
!print*,delsqn(i),grgrkerni,weight
!delsqn(i) = delsqn(i) + grgrkerni*weight
if (onef_dust) then
dustfracj = sum(dustfrac(:,j))
rhodust(:,i) = rhodust(:,i) + pmassj*dustfrac(:,j)*wabi*weight
rhogas(i) = rhogas(i) + pmassj*(1. - dustfracj)*wabi*weight
endif
endif
if (itypej /= itypebnd) then
rho(j) = rho(j) + pmassi*wabj*weight
densn(j) = densn(j) + wabaltj*weight
!delsqn(j) = delsqn(j) + grgrkernj*weight
if (onef_dust) then
rhodust(:,j) = rhodust(:,j) + pmassi*dustfrac(:,i)*wabj*weight
rhogas(j) = rhogas(j) + pmassi*(1. - dustfraci)*wabj*weight
endif
endif
!
!--drhodt, dndt
!
if (i.ne.j) then
dvel(1:ndimV) = veli(1:ndimV) - vel(1:ndimV,j)
dvdotr = dot_product(dvel,dr)
drhodt(i) = drhodt(i) + pmassj*dvdotr*grkerni !+ pmassj*dvdotr*wabi*hi1
drhodt(j) = drhodt(j) + pmassi*dvdotr*grkernj !+ pmassi*dvdotr*wabj*hj1
dndt(i) = dndt(i) + dvdotr*grkernalti
dndt(j) = dndt(j) + dvdotr*grkernaltj
if (imhd.eq.5) then
projBi = dot_product(Bfield(:,i),dr)
projBj = dot_product(Bfield(:,j),dr)
dBevoldt(:,i) = dBevoldt(:,i) - pmassj*projBi*dvel(:)*grkerni
dBevoldt(:,j) = dBevoldt(:,j) - pmassi*projBj*dvel(:)*grkernj
endif
else
!drhodt(i) = drhodt(i) + pmassj*dvdotr*wabi*hi1
endif
if (ikernav.EQ.3) then
!
!--correction term for variable smoothing lengths
! this is the small bit that should be 1-gradh
! need to divide by rho once rho is known
! also do the number density version
if (itypei /= itypebnd) then
gradh(i) = gradh(i) + weight*pmassj*dwdhi
gradhn(i) = gradhn(i) + weight*dwaltdhi
gradgradh(i) = gradgradh(i) + weight*pmassj*dwdhdhi
endif
if (itypej /= itypebnd) then
gradh(j) = gradh(j) + weight*pmassi*dwdhj
gradhn(j) = gradhn(j) + weight*dwaltdhj
gradgradh(j) = gradgradh(j) + weight*pmassi*dwdhdhj
endif
endif
if (imhd.eq.0) then
if (iprterm.eq.10) then
psi(i) = psi(i) + pmassj*wabi*uu(j)
unity(i) = unity(i) + wconst*wabi/hfacwabi
if (i.ne.j) then
psi(j) = psi(j) + pmassi*wabj*uu(i)
unity(j) = unity(j) + wconst*wabj/hfacwabj
endif
elseif (iprterm.eq.12) then
psi(i) = psi(i) + wconst*wabi/hfacwabi
if (i.ne.j) then
psi(j) = psi(j) + wconst*wabj/hfacwabj
endif
endif
endif
if (i.ne.j) then
do k=1,ndxdx
dxdxi(k) = dxdxi(k) - pmass(j)*(dx(idxdx(k)))*dr(jdxdx(k))*grkerni
dxdx(k,j) = dxdx(k,j) - pmass(i)*(dx(idxdx(k)))*dr(jdxdx(k))*grkernj
enddo
endif
! ELSE
! PRINT*,' r/h > 2 '
endif
enddo loop_over_neighbours
dxdx(:,i) = dxdx(:,i) + dxdxi(:)
iprev = i
if (iprev.NE.-1) i = ll(i) ! possibly should be only IF (iprev.NE.-1)
enddo loop_over_cell_particles
enddo loop_over_cells
if (imhd.eq.5) then
do i=1,npart
dBevoldt(:,i) = dBevoldt(:,i)/rhoin(i)**2
enddo
endif
if (onef_dust) then
do i=1,npart
!where (rhodust(:,i) < 0.) rhodust(:,i) = 0.
enddo
endif
return
end subroutine density
!!------------------------------------------------------------------------
!! Computes the density by direct summation over the particles neighbours
!! ie. rho_a = sum_b m_b W_ab (h_a)
!!
!! This version computes the density only on a selected list of particles
!! given by the contents of the array ipartlist (enables iteration on
!! unconverged particles only). It is therefore slightly slower
!! since some particle pairs may be done twice, once for each particle.
!!
!! Assumes rho_a is only a function of h_a
!!
!! This version must be used for individual particle timesteps
!!------------------------------------------------------------------------
subroutine density_partial(x,pmass,hh,vel,rho,drhodt,densn,dndt,delsqn, &
gradh,gradhn,gradsoft,gradgradh,ntotal,nlist,ipartlist)
use dimen_mhd, only:ndim,ndimV
use debug, only:trace
use loguns, only:iprint
use kernels, only:radkern2,interpolate_kernels_dens,interpolate_kernel_soft
use linklist, only:iamincell,numneigh
use options, only:igravity,imhd,ikernel,ikernelalt,iprterm,onef_dust
use part, only:Bfield,uu,psi,itype,itypebnd,dustfrac,rhodust,rhogas
use rates, only:dBevoldt
use setup_params, only:hfact
use matrixcorr, only:dxdx,ndxdx,idxdx,jdxdx
!
!--define local variables
!
implicit none
integer, intent(in) :: ntotal
real, dimension(:,:), intent(in) :: x, vel
real, dimension(:), intent(in) :: pmass, hh
real, dimension(:), intent(out) :: rho,drhodt,densn,dndt,delsqn,gradh,gradhn,gradsoft,gradgradh
integer, intent(in) :: nlist
integer, intent(in), dimension(:) :: ipartlist
integer :: i,j,k,n
integer :: icell,ipart,nneigh !!,minneigh,minpart
integer, dimension(ntotal) :: listneigh
integer :: icellprev
!
! (particle properties - local copies)
!
real :: rij,rij2
real :: hi,hi1,hi2,hi21
real :: hfacwabi,hfacgrkerni,pmassj,dustfracj
real, dimension(ndim) :: dx,xi
real, dimension(ndimV) :: veli,dvel
real, dimension(nlist) :: rhoin
real :: dvdotr,projBi
real, dimension(ndxdx) :: dxdxi
real, dimension(ndimV) :: dr
integer :: itypei
!
! (kernel quantities)
!
real :: q2i
real :: wabi,wabalti,grkerni,grgrkerni,grkernalti
real :: dwdhi,dwaltdhi,dphidhi,dwdhdhi ! grad h terms
real :: wconst,unityi
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' Entering subroutine density_partial'
!
!--initialise quantities
!
listneigh = 0
wconst = 1./hfact**ndim
dr(:) = 0.
do ipart=1,nlist
i = ipartlist(ipart)
rhoin(ipart) = rho(i)
rho(i) = 0.
drhodt(i) = 0.
densn(i) = 0.
dndt(i) = 0.
delsqn(i) = 0.
gradh(i) = 0.
gradhn(i) = 0.
gradsoft(i) = 0.
gradgradh(i) = 0.
numneigh(i) = 0
if (imhd.eq.5) dBevoldt(:,i) = 0.
if (imhd.eq.0 .and. iprterm.eq.10) psi(i) = 0.
if (onef_dust) then
rhodust(:,i) = 0.
rhogas(i) = 0.
endif
dxdx(:,i) = 0.
enddo
icellprev = 0
!
!--Loop over all the particles in the density list
!
loop_over_particles: do ipart=1,nlist ! step through all cells
i = ipartlist(ipart)
!
!--find cell of current particle
!
icell = iamincell(i)
! PRINT*,' particle ',i,' cell = ',icell
!
!--if different to previous cell used, get the list of neighbours for this cell
! (common to all particles in the cell)
!
if (icell.NE.icellprev) then
call get_neighbour_list_partial(icell,listneigh,nneigh)
endif
icellprev = icell
!!!PRINT*,'Doing particle ',i,nneigh,' neighbours',hh(i)
hi = hh(i)
hi2 = hi*hi
hi1 = 1./hi
hi21 = hi1*hi1
unityi = 0.
hfacwabi = hi1**ndim
hfacgrkerni = hfacwabi*hi1
xi = x(:,i)
veli(:) = vel(:,i)
dxdxi(:) = 0.
itypei = itype(i)
if (itypei.eq.itypebnd) cycle loop_over_particles
!
!--loop over current particle's neighbours
!
loop_over_neighbours: do n = 1,nneigh
j = listneigh(n)
!--skip particles of different type
if (itype(j).ne.itypei .and. itype(j).ne.itypebnd) cycle loop_over_neighbours
dx(:) = xi(:) - x(:,j)
!
!--calculate averages of smoothing length if using this averaging
!
rij2 = dot_product(dx,dx)
q2i = rij2*hi21
!
!--do interaction if r/h < compact support size
!
if (q2i.LT.radkern2) then
rij = sqrt(rij2)
dr(1:ndim) = dx(1:ndim)/(rij + epsilon(rij))
!!!if (i.eq.416) PRINT*,' neighbour,r/h,hi ',j,SQRT(q2i),hi
numneigh(i) = numneigh(i) + 1
pmassj = pmass(j)
!
!--interpolate from kernel table (using hi)
!
if (igravity.ne.0) then
call interpolate_kernel_soft(q2i,wabi,grkerni,dphidhi)
gradsoft(i) = gradsoft(i) + pmassj*dphidhi/hi2
if (ikernel.ne.ikernelalt) then
call interpolate_kernels_dens(q2i,wabi,grkerni,grgrkerni,wabalti,grkernalti)
else
wabalti = wabi
grkernalti = grkerni
endif
else
call interpolate_kernels_dens(q2i,wabi,grkerni,grgrkerni,wabalti,grkernalti)
endif
wabi = wabi*hfacwabi
wabalti = wabalti*hfacwabi
grkerni = grkerni*hfacgrkerni
grgrkerni = grgrkerni*hfacwabi*hi1*hi1
grkernalti = grkernalti*hfacgrkerni
!
!--derivative w.r.t. h for grad h correction terms (and dhdrho)
!
dwdhi = -rij*grkerni*hi1 - ndim*wabi*hi1
dwaltdhi = -rij*grkernalti*hi1 - ndim*wabalti*hi1
dwdhdhi = ndim*(ndim+1)*wabi*hi1**2 + 2.*(ndim+1)*rij*hi1**2*grkerni &
+ rij**2*hi1**2*grgrkerni
!
!--calculate density and number density
!
rho(i) = rho(i) + pmassj*wabi
densn(i) = densn(i) + wabalti
delsqn(i) = delsqn(i) + grgrkerni
if (onef_dust) then
dustfracj = sum(dustfrac(:,j))
rhodust(:,i) = rhodust(:,i) + pmassj*dustfrac(:,j)*wabi
rhogas(i) = rhogas(i) + pmassj*(1. - dustfracj)*wabi
endif
!
!--drhodt, dndt
!
if (i.ne.j) then
dvel(:) = veli(:) - vel(:,j)
dvdotr = dot_product(dvel,dr)
drhodt(i) = drhodt(i) + pmassj*dvdotr*grkerni !+ pmassj*dvdotr*wabi*hi1
dndt(i) = dndt(i) + dvdotr*grkernalti
if (imhd.eq.5) then
projBi = dot_product(Bfield(:,i),dr)
dBevoldt(:,i) = dBevoldt(:,i) - pmassj*projBi*dvel(:)*grkerni
endif
else
!drhodt(i) = drhodt(i) + pmassj*dvdotr*grkerni + pmassj*dvdotr*wabi*hi1
endif
!
!--correction term for variable smoothing lengths
! this is the small bit that should be 1-gradh
! need to divide by rho once rho is known
gradh(i) = gradh(i) + pmassj*dwdhi
gradhn(i) = gradhn(i) + dwaltdhi
gradgradh(i) = gradgradh(i) + pmassj*dwdhdhi
if (imhd.eq.0 .and. iprterm.eq.10) then
psi(i) = psi(i) + pmassj*wabi*uu(j)
unityi = unityi + wconst*wabi/hfacwabi
endif
if (i.ne.j) then
do k=1,ndxdx
dxdxi(k) = dxdxi(k) - pmass(j)*(dx(idxdx(k)))*dr(jdxdx(k))*grkerni
enddo
endif
endif
enddo loop_over_neighbours
dxdx(:,i) = dxdx(:,i) + dxdxi(:)
if (imhd.eq.0) then
!psi(i) = abs(uu(i) - psi(i)/unityi)
!if (psi(i)/uu(i).lt.0.01) psi(i) = 0.
else
psi(i) = 0.
endif
!psi(i) = abs(psi(i) - uu(i)*unityi)
enddo loop_over_particles
!print*,'finished density_partial, rho, gradh, h =',rho(1),gradh(1),hh(1),numneigh(1)
!print*,'maximum number of neighbours = ',MAXVAL(numneigh),MAXLOC(numneigh),rho(MAXLOC(numneigh))
!print*,'minimum number of neighbours = ',MINVAL(numneigh(1:npart)), &
! MINLOC(numneigh(1:npart)),rho(MINLOC(numneigh(1:npart)))
!minneigh = 100000
!minpart = 1
!do i=1,nlist
! j = ipartlist(i)
! if (numneigh(j).lt.minneigh) then
! minneigh = numneigh(j)
! minpart = j
! endif
!enddo
!print*,'minimum number of neighbours = ',minneigh,minpart,rho(minpart)
if (imhd.eq.5) then
do i=1,nlist
j = ipartlist(i)
dBevoldt(:,j) = dBevoldt(:,j)/rhoin(i)**2
enddo
endif
if (onef_dust) then
do i=1,nlist
j = ipartlist(i)
!where (rhodust(:,j) < 0.) rhodust(:,j) = 0.
enddo
endif
!do i=1,nlist
! j = ipartlist(i)
! psi(j) = abs(uu(j) - psi(j)) !/psi(j)
!enddo
return
end subroutine density_partial
end module density_summations