-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathratesND_mhd.f90
2936 lines (2728 loc) · 107 KB
/
ratesND_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
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
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!------------------------------------------------------------------------------!
! 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: !
!------------------------------------------------------------------------------!
!!--------------------------------------------------------------------
!! Computes the rates of change of the conserved variables
!! (forces, energy etc)
!! This is the core of the SPH algorithm
!!--------------------------------------------------------------------
subroutine get_rates
! USE dimen_mhd
use debug, only:trace
use loguns, only:iprint
use bound, only:pext
use artvi
use eos
use hterms
use kernels, only:radkern2
use linklist
use options
use part
use rates
use timestep
use xsph
use fmagarray
use derivB
use get_neighbour_lists
use grutils, only:metric_diag,dot_product_gr
use getBeulerpots, only:compute_rmatrix,exactlinear
use matrixcorr, only:dxdx,idxdx,jdxdx,ndxdx
use utils, only:cross_product3D
use resistivity, only:etafunc
use externf, only:external_forces,pequil
use dust, only:get_tstop
!
!--define local variables
!
implicit none
integer :: i,j,n,k
integer :: icell,iprev,nneigh,nlistdim
integer, allocatable, dimension(:) :: listneigh
integer :: idone,nclumped
!
! (particle properties - local copies and composites)
!
real :: rij,rij2
real :: rhoi,rho1i,rho2i,rho21i,rhoj,rho1j,rho2j,rho21j,rhoav1,rhoij
real :: pmassi,pmassj,projvi,projvj
real :: Prho2i,Prho2j,prterm,pri,prj
real :: hi,hi1,hj,hj1,hi21,hj21
real :: hfacwabi,hfacgrkerni
real, dimension(ndim) :: xi, dx
real, dimension(ndimV) :: fexternal !!,dri,drj
real, dimension(ntotal) :: h1
integer :: itypei,itypej
!
!--gr terms
!
real :: sqrtgi,sqrtgj,dens1i,v2i,v2j
real, dimension(ndimV) :: gdiagi,gdiagj
!
! (velocity)
!
real, dimension(ndimV) :: veli,velj,dvel,fmean
real, dimension(ndimV) :: dr
real :: dvdotr
!
! (mhd)
!
real, dimension(ndimB) :: Brhoi,Brhoj,Bi,Bj,dB
real, dimension(ndimB) :: faniso,fmagi,fmagj,Bevoli,dBevol
real, dimension(ndimB) :: curlBi,forcei,forcej,dBevoldti
real :: fiso,B2i,B2j
real :: valfven2i,valfven2j
real :: BidotdB,BjdotdB,Brho2i,Brho2j,BdotBexti
real :: projBrhoi,projBrhoj,projBi,projBj,projdB,projBconst
real, dimension(:,:), allocatable :: curlBsym
real, dimension(:), allocatable :: divBsym
real, dimension(:,:,:), allocatable :: dveldx
real, dimension(6) :: rmatrix
real :: denom,ddenom,etai,etaj
!
! (artificial viscosity quantities)
!
real :: vsig,vsigi,vsigj,vsigav
real :: spsoundi,spsoundj,alphai,alphaui,alphaBi
!! real :: rhoi5,rhoj5
real :: vsig2i,vsig2j,vsigproji,vsigprojj,vsigB,vsigu
!! real :: vsigii,vsigjj
real :: prneti,prnetj
!
! (av switch)
!
real :: source,tdecay1,sourcedivB,sourceJ,sourceB,sourceu
real :: graduterm, graddivvmag,curr2
real, dimension(:), allocatable :: del2u
!
! (alternative forms)
!
real, dimension(:), allocatable :: phi
real :: phii,phii1,phii_on_phij,phij_on_phii,uui,uuj
!
! (one fluid dust)
!
real, dimension(ndimV) :: deltavi,deltavj
real :: dustfraci(ndust),dustfracj(ndust)
real :: rhogrhodonrhoi(ndust), rhogrhodonrhoj(ndust)
real :: deltav2i,deltav2j
real :: dtstop,projdeltavi,projdeltavj
real :: rhogasi,rhogasj,projdvgas
real :: rhodusti(ndust),rhodustj(ndust)
real, dimension(ndimV) :: vgasi,vgasj,dvgas,fextrai,fextraj
real :: esum,tstop,ratio,dvmax
!
! (kernel related quantities)
!
real :: q2i,q2j
real :: wab,wabi,wabj
real :: grkern,grkerni,grkernj
real :: gradhi,gradhni
real :: grkernalti,grkernaltj,grgrkernalti,grgrkernaltj
!
! (time step criteria)
!
real :: vsigdtc,zero,fhmax, fonh, forcemag, ts_min, h_on_csts_max
integer :: ierr
!
!--div B correction
!
real :: gradpsiterm,vsig2,vsigmax !!,dtcourant2
real :: stressmax,stressterm
logical, parameter :: itiming = .false.
real :: t1,t2,t3,t4,t5
real, allocatable :: fprev(:,:)
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' Entering subroutine get_rates'
if (itiming) call cpu_time(t1)
!
!--allocate memory for local arrays
!
nlistdim = ntotal
allocate ( listneigh(nlistdim),STAT=ierr )
if (ierr.ne.0) write(iprint,*) ' Error allocating neighbour list, ierr = ',ierr
allocate ( phi(ntotal), del2u(ntotal), STAT=ierr )
if (ierr.ne.0) write(iprint,*) ' Error allocating phi, ierr = ',ierr
if (imhd.lt.0) then
allocate( curlBsym(ndimV,ntotal), STAT=ierr )
if (ierr.ne.0) write(iprint,*) ' Error allocating curlBsym, ierr = ',ierr
allocate( divBsym(ntotal), STAT=ierr )
if (ierr.ne.0) write(iprint,*) ' Error allocating divBsym, ierr = ',ierr
endif
if (imhd.ne.0 .and. iuse_exact_derivs.gt.0) then
allocate( dveldx(ndim,ndimV,ntotal), STAT=ierr )
if (ierr.ne.0) write(iprint,*) ' Error allocating dveldx, ierr = ',ierr
endif
listneigh = 0
!
!--initialise quantities
!
dtcourant = 1.e6
dtav = huge(dtav)
dtvisc = huge(dtvisc)
ts_min = huge(ts_min)
h_on_csts_max = 0.
zero = 1.e-10
vsigmax = 0.
dr(:) = 0.
nclumped = 0
do i=1,ntotal ! using ntotal just makes sure they are zero for ghosts
force(:,i) = 0.0
dudt(i) = 0.0
dendt(i) = 0.0
dBevoldt(:,i) = 0.0
daldt(:,i) = 0.0
dpsidt(i) = 0.0
gradpsi(:,i) = 0.0
fmag(:,i) = 0.0
divB(i) = 0.0
if (imhd.gt.0 .and. iambipolar.eq.0 .and. iresist.ne.4) curlB(:,i) = 0.0
if (allocated(curlBsym)) curlBsym(:,i) = 0.
if (allocated(divBsym)) divBsym(i) = 0.
if (allocated(dveldx)) then
dveldx(:,:,i) = 0.
dxdx(:,i) = 0.
endif
if (allocated(del2v)) del2v(:,i) = 0.
xsphterm(:,i) = 0.0
del2u(i) = 0.0
graddivv(:,i) = 0.0
h1(i) = 1./hh(i)
if (icty.ge.1) then
drhodt(i) = 0.
gradh(i) = 1.
gradhn(i) = 0.
dhdt(i) = 0.
endif
if (onef_dust) then
ddustevoldt(:,i) = 0.
ddeltavdt(:,i) = 0.
endif
enddo
!
!--calculate maximum neg stress for instability correction
!
stressmax = 0.
if (imhd.ne.0 .and. (imagforce.eq.2 .or. imagforce.eq.7)) then
do i=1,ntotal
call metric_diag(x(:,i),gdiagi(:),sqrtgi,ndim,ndimV,geom)
B2i = dot_product_gr(Bfield(:,i),Bfield(:,i),gdiagi(:))
if (imagforce.eq.7) then ! vec potential force
!stressterm = max(1.5*B2i - pr(i),0.)
!stressmax = max(stressterm,stressmax)
else
stressterm = max(0.5*B2i - pr(i),0.)
stressmax = max(stressterm,stressmax,maxval(Bconst(:)))
endif
enddo
if (stressmax.gt.tiny(stressmax)) write(iprint,*) 'stress correction = ',stressmax
endif
!
!--set MHD quantities to zero if mhd not set
!
if (imhd.eq.0) then ! these quantities are still used if mhd off
Bi(:) = 0.
Bj(:) = 0.
Brhoi(:) = 0.
Brhoj(:) = 0.
Brho2i = 0.
Brho2j = 0.
valfven2i = 0.
valfven2j = 0.
projBi = 0.
projBj = 0.
projBrhoi = 0.
projBrhoj = 0.
alphaBi = 0.
endif
!
!--skip the whole neighbour thing if it is doing nothing
!
if (iprterm.lt.0 .and. iav.eq.0 .and. imhd.eq.0 .and. iener.lt.3 .and. iprterm.ne.-2) then
write(iprint,*) 'skipping rates'
goto 666
endif
!
! set alternative forms for the SPH equations here
! phi can be any scalar variable
!
select case(iprterm)
case(1) ! this gives the (P_a + P_b) / (rho_a rho_b) form
phi(1:ntotal) = rho(1:ntotal)
case(2) ! this gives the HK89 form SQRT(Pa Pb)/rhoa rhob
phi(1:ntotal) = sqrt(pr(1:ntotal))/rho(1:ntotal)
case(3)
phi(1:ntotal) = 1./rho(1:ntotal)
case(4)
phi(1:ntotal) = 1./rho(1:ntotal)**gamma
case(5)
phi(1:ntotal) = rho(1:ntotal)**gamma
case(6)
phi(1:ntotal) = 1./en(1:ntotal)
case(7)
phi(1:ntotal) = 1./uu(1:ntotal)
case(8)
phi(1:ntotal) = 1./pr(1:ntotal)
case(9)
phi(1:ntotal) = sqrt(rho(1:ntotal))
case default ! this gives the usual continuity, momentum and induction eqns
phi(1:ntotal) = 1.0
end select
if (itiming) call cpu_time(t2)
!
!--Loop over all the link-list cells
!
loop_over_cells: do icell=1,ncellsloop ! step through all cells
!!print*,'> doing cell ',icell
!
!--get the list of neighbours for this cell
! (common to all particles in the cell)
!
call get_neighbour_list(icell,listneigh,nneigh)
i = ifirstincell(icell) ! start with first particle in cell
idone = -1
if (i.ne.-1) iprev = i
loop_over_cell_particles: do while (i.ne.-1) ! loop over home cell particles
!!print*,'Doing particle ',i,'of',npart,x(:,i),rho(i),hh(i)
idone = idone + 1
xi(:) = x(:,i)
itypei = itype(i)
rhoi = rho(i)
rho2i = rhoi*rhoi
!! rhoi5 = sqrt(rhoi)
rho1i = 1./rhoi
rho21i = rho1i*rho1i
dens1i = 1./dens(i)
pri = max(pr(i) - pext,0.)
! print *, 'PRI ', pr(i), pext, itype(i)
prneti = pri - pequil(iexternal_force,xi(:),rhoi)
pmassi = pmass(i)
Prho2i = pri*rho21i
spsoundi = spsound(i)
uui = uu(i)
veli(:) = vel(:,i)
v2i = dot_product(veli(:),veli(:))
alphai = alpha(1,i)
alphaui = alpha(2,i)
alphaBi = alpha(3,i)
phii = phi(i)
phii1 = 1./phii
sqrtgi = sqrtg(i)
! one fluid dust definitions
if (onef_dust) then
dustfraci = dustfrac(:,i)
if (use_smoothed_rhodust) then
rhodusti = rhodust(:,i)
rhogasi = rhogas(i)
else
rhodusti = dustfraci*rhoi
rhogasi = (1. - sum(dustfraci))*rhoi
endif
deltavi(:) = deltav(:,i)
deltav2i = dot_product(deltavi,deltavi)
rhogrhodonrhoi(:) = rhogasi*rhodusti(:)*rho1i
else
rhogasi = rhoi
rhodusti = 0.
deltav2i = 0.
endif
! mhd definitions
if (imhd.ne.0) then
Bi(:) = Bfield(:,i)
Brhoi(:) = Bi(:)*rho1i
BdotBexti = dot_product(Bi(:),Bconst(:))
! if (geom(1:6).ne.'cartes') then
call metric_diag(x(:,i),gdiagi(:),sqrtgi,ndim,ndimV,geom)
B2i = dot_product_gr(Bi,Bi,gdiagi)
! else
! B2i = dot_product(Bi,Bi)
! endif
Brho2i = B2i*rho21i
valfven2i = B2i*rho1i
if (imhd.lt.0) Bevoli(:) = Bevol(:,i)
if (iresist.eq.3) then
etai = etafunc(x(1,i),etamhd)
elseif (iresist > 0) then
etai = etamhd
endif
endif
gradhi = gradh(i)
gradhni = gradhn(i)
hi = hh(i)
if (hi.le.0.) then
write(iprint,*) ' rates: h <= 0 particle',i,hi
call quit
endif
hi1 = h1(i)
hi21 = hi1*hi1
hfacwabi = hi1**ndim
hfacgrkerni = hfacwabi*hi1
forcei(:) = 0.
fextrai(:) = 0.
dBevoldti(:) = 0.
!
!--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(:) = xi(:) - x(:,j)
!print*,' ... neighbour, h=',j,hh(j),rho(j),x(:,j)
hj = hh(j)
hj1 = h1(j) !!1./hj
hj21 = hj1*hj1
rij2 = dot_product(dx,dx)
q2i = rij2*hi21
q2j = rij2*hj21
!----------------------------------------------------------------------------
! do pairwise interaction if either particle is within range of the other
!----------------------------------------------------------------------------
if ((q2i.lt.radkern2).or.(q2j.lt.radkern2)) then ! if < 2h
rij = sqrt(rij2)
if (rij.le.epsilon(rij) .and. itype(j).eq.itypei) then
nclumped = nclumped + 1
if (rij.lt.tiny(rij)) then
dr = 0.
if (itypei /= 2) then
write(iprint,*) 'rates: dx = 0 i,j,dx,hi,hj=',i,j,dx,hi,hj,' type=',itypei
call quit
endif
endif
elseif (rij.le.epsilon(rij)) then
dr = 0. ! can happen with gas/dust if particles on top of each other
else
dr(1:ndim) = dx(1:ndim)/rij ! unit vector
endif
!do idim=1,ndim
! dri(idim) = dot_product(1./gradmatrix(idim,1:ndim,i),dr(1:ndim))
! drj(idim) = dot_product(1./gradmatrix(idim,1:ndim,j),dr(1:ndim))
!enddo
itypej = itype(j)
if ((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
! if (itype(j).eq.itypei &
! .or.(itype(j).eq.itypebnd .or. itype(j).eq.itypebnd2) &
! .or.(itypei .eq.itypebnd .or. itype(j).eq.itypebnd2)) then
call rates_core
elseif (idust.eq.2 .and. idrag_nature.gt.0) then !-- drag step if required
call drag_forces
endif
else ! if outside 2h
! PRINT*,'outside 2h, not calculated, r/h=',sqrt(q2i),sqrt(q2j)
endif ! r < 2h or not
endif ! j.ne.i
enddo loop_over_neighbours
!
!--add contributions to particle i from summation over j
! (forcei is forces on gas only, fextra is terms that apply to total fluid)
!
force(:,i) = force(:,i) + fextrai(:) + forcei(:)
dBevoldt(:,i) = dBevoldt(:,i) + dBevoldti(:)
if (idust.eq.1) ddeltavdt(:,i) = ddeltavdt(:,i) - rhoi/rhogasi*forcei(:)
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 (nclumped.gt.0) write(iprint,*) ' WARNING: clumping on ',nclumped,' pairs'
666 continue
if (itiming) call cpu_time(t3)
!----------------------------------------------------------------------------
! calculate gravitational force on all the particles
!----------------------------------------------------------------------------
if (igravity.ne.0) then
fprev = force
select case(igravity)
case(1) ! fixed plummer
call direct_sum_poisson(x(1:ndim,1:npart),pmass(1:npart), &
potengrav,force(1:ndim,1:npart),hsoft,npart)
case(2) ! fixed cubic spline
phi(1:npart) = hsoft
call direct_sum_poisson_soft(x(1:ndim,1:npart),pmass(1:npart), &
phi(1:npart),poten(1:npart),force(1:ndim,1:npart),potengrav,npart)
case(3) ! adaptive cubic spline (no extra term)
call direct_sum_poisson_soft(x(1:ndim,1:npart),pmass(1:npart), &
hh(1:npart),poten(1:npart),force(1:ndim,1:npart),potengrav,npart)
case(4:) ! adaptive cubic spline with energy conservation
call direct_sum_poisson_soft(x(1:ndim,1:npart),pmass(1:npart), &
hh(1:npart),poten(1:npart),force(1:ndim,1:npart),potengrav,npart)
case default
stop 'unknown igravity setting in forces'
end select
!
!--add v.fgrav term if evolving total energy
!
!if (iener==3) then
! do i=1,npart
! dendt(i) = dendt(i) + dot_product(vel(:,i),force(:,i)-fprev(:,i))
! enddo
!endif
endif
if (itiming) call cpu_time(t4)
if (trace) write(iprint,*) 'Finished main rates loop'
!----------------------------------------------------------------------------
! loop over the particles again, subtracting external forces and source terms
!----------------------------------------------------------------------------
!
!--calculate maximum vsig over all the particles for use in the hyperbolic cleaning
!
if (imhd.ne.0 .and. idivBzero.ge.2) then
vsig2max = vsigmax**2
endif
fhmax = 0.0
fmean(:) = 0.
dtdrag = huge(dtdrag)
dtforce= huge(dtforce)
esum = 0.
ratio = 0.
if (idust.eq.2 .and. h_on_csts_max > 1.) then
print*,' WARNING: violating h < cs*ts resolution criterion by factor of ',h_on_csts_max
endif
do i=1,npart
rhoi = rho(i)
rho1i = 1./rhoi
if (ivisc.gt.0 .and. idust.ne.1 .and. idust.ne.4) then
esum = esum + pmass(i)*dot_product(vel(:,i),force(:,i))
endif
fexternal(:) = 0.
!
!--Dust
!
if (idust.eq.2 .and. idrag_nature.ne.0 .and. (Kdrag > 0. .or. idrag_nature > 1)) then
!
!--two fluid dust: calculate drag timestep
!
dtdrag = min(dtdrag,ts_min)
elseif (idust.eq.1) then
!------------------
! one fluid dust
!------------------
if (use_smoothed_rhodust) then
rhodusti = rhodust(:,i)
rhogasi = rhogas(i)
else
rhodusti = rhoi*dustfrac(:,i)
rhogasi = (1. - sum(dustfrac(:,i)))*rhoi
endif
do k=1,ndust
tstop = get_tstop(idrag_nature,rhogasi,rhodusti(k),spsound(i),Kdrag)
dtdrag = min(dtdrag,tstop)
!
!--d/dt(deltav) : add terms that do not involve sums over particles
!
if (dustfraci(k).gt.0.) then
dtstop = 1./tstop
ddeltavdt(:,i) = ddeltavdt(:,i) - deltav(:,i)*dtstop
else
dtstop = 0.
ddeltavdt(:,i) = 0.
endif
enddo
!
!--du/dt: add thermal energy dissipation from drag
! (maybe should do this in the step routine along with the implicit drag in ddeltav/dt?)
!
if (iener.gt.0) then
deltav2i = dot_product(deltav(:,i),deltav(:,i))
dudt(i) = dudt(i) + rhodusti(1)*rho1i*deltav2i*dtstop
endif
elseif (idust.eq.3 .or. idust.eq.4) then
!--------------------------------------------
! one fluid dust in diffusion approximation
!--------------------------------------------
if (use_smoothed_rhodust) then
rhodusti(:) = rhodust(:,i)
rhogasi = rhogas(i)
else
rhodusti(:) = rhoi*dustfrac(:,i)
rhogasi = (1 - sum(dustfrac(:,i)))*rhoi
endif
!
!--compute stopping time for drag timestep
!
dustfraci = rhodusti*rho1i
do k=1,ndust
tstop = get_tstop(idrag_nature,rhogasi,rhodusti(k),spsound(i),Kdrag)
! CAUTION: Line below must be done BEFORE external forces have been applied
if (rhogasi > tiny(rhogasi)) then
deltav(:,i) = -rhoi/rhogasi*force(:,i)*tstop
else
!stop 'rhogas < 0'
deltav(:,i) = 0.
endif
ratio = max(dustfraci(k)*tstop/dtcourant,ratio)
dtdrag = min(dtdrag,0.5*hh(i)**2/(dustfraci(k)*tstop*spsound(i)**2))
enddo
if (dtdrag < 0.) then
!print*,'WARNING: dtdrag = ',dtdrag,dustfraci,dustfrac(:,i),rhodust(:,i),rhogas(i)
dtdrag = abs(dtdrag)
endif
!dvmax = maxval(abs(deltav(:,i)))
!if (dvmax > 0.) then
! dtdrag = min(dtdrag,0.1*hh(i)/dvmax)
!endif
if (idrag_nature==4) then
! this is constant ts but keeping the particles fixed
force(:,i) = 0.
vel(:,i) = 0.
! !print*,'dtdrag = ',dtdrag
endif
endif
!
!--compute JxB force from the Euler potentials (if using second derivs)
! also divide by rho for J and div B calculated the "normal" way
!
if (imhd.ne.0) then
if (iambipolar > 0) then
B2i = dot_product(Bfield(:,i),Bfield(:,i))
dtdrag = min(dtdrag,gamma_ambipolar*rho_ion*rhoi*hh(i)**2/B2i)
else
if (imagforce.eq.4) then
call cross_product3D(curlB(:,i),Bfield(:,i),fmagi(:)) ! J x B
force(:,i) = force(:,i) + fmagi(:)*rho1i ! (J x B)/rho
fmag(:,i) = fmagi(:)*rho1i
elseif (imhd.gt.0 .and. .not.(iresist==4)) then ! not for vector potential
curlB(:,i) = curlB(:,i)*rho1i
endif
endif
divB(i) = divB(i)*rho1i
if (allocated(divBsym)) then
divBsym(i) = divBsym(i)*rhoi
divB(i) = divBsym(i)
endif
endif
!
!--add external (body) forces
!
if (iexternal_force.ne.0) then
call external_forces(iexternal_force,x(1:ndim,i),fexternal(1:ndimV), &
ndim,ndimV,vel(1:ndimV,i),hh(i),spsound(i),itype(i))
force(1:ndimV,i) = force(1:ndimV,i) + fexternal(1:ndimV)
endif
!
!--add source terms (derivatives of metric) to momentum equation
!
if (allocated(sourceterms)) force(:,i) = force(:,i) + sourceterms(:,i)
!
!--make dhdt if density is not being done by summation
! (otherwise this is done in iterate_density)
!
if (icty.ge.1) then
dhdt(i) = -hh(i)/(ndim*rho(i))*drhodt(i)
endif
!
!--calculate maximum force/h for the timestep condition
! also check for errors in the force
!
! if ( any(force(:,i).gt.1.e8)) then
! write(iprint,*) 'rates: force ridiculous ',force(:,i),' particle ',i
! call quit
! endif
fmean(:) = fmean(:) + pmass(i)*force(:,i)
forcemag = sqrt(dot_product(force(:,i),force(:,i)))
fonh = forcemag/hh(i)
if (fonh.gt.fhmax .and. itype(i).ne.1) fhmax = fonh
!
!--calculate simpler estimate of vsig for divergence cleaning and
! in the dissipation switches
!
valfven2i = 0.
if (imhd.ne.0) then
! if (geom(1:6).ne.'cartes') then
call metric_diag(x(:,i),gdiagi,sqrtgi,ndim,ndimv,geom)
valfven2i = dot_product_gr(Bfield(:,i),Bfield(:,i),gdiagi)/dens(i)
! else
! valfven2i = dot_product(Bfield(:,i),Bfield(:,i))/dens(i)
! endif
endif
vsig2 = spsound(i)**2 + valfven2i ! approximate vsig only
vsig = SQRT(vsig2)
!!!dtcourant2 = min(dtcourant2,hh(i)/vsig)
if (iuse_exact_derivs.gt.0) then
call compute_rmatrix(dxdx(:,i),rmatrix,denom,ndim)
! print*,i,' dxdx= ',dxdx(:,i)
! print*,i,' rmatrix = ',rmatrix*rho1i*rho1i,'denom= ',denom*rho1i*rho1i*rho1i
! print*,i,' normal derivs, dvxdy = ',dveldx(2,1,i)*rho1i
if (abs(denom).gt.epsilon(denom)) then
ddenom = 1./denom
do k=1,ndimV
call exactlinear(dveldx(:,k,i),dveldx(:,k,i),rmatrix,ddenom)
enddo
! print*,i,' exact derivs, dvxdy = ',dveldx(2,1,i)
else
print*,'WARNING: denominator collapsed in exact linear deriv = ',denom
dveldx(:,:,i) = dveldx(:,:,i)*rho1i
endif
!read*
endif
!
!--if evolving B instead of B/rho, add the extra term from the continuity eqn
! (note that the dBevoldt term should be divided by rho)
!
select case(imhd)
case(11:19,21:) ! evolving B
dBevoldt(:,i) = sqrtg(i)*dBevoldt(:,i) + Bevol(:,i)*rho1i*drhodt(i)
if (idivBzero.ge.2) then
!gradpsi(:,i) = gradpsi(:,i)*rho1i
gradpsi(:,i) = gradpsi(:,i)*rhoi
if (nsubsteps_divB.le.0) then
dBevoldt(:,i) = dBevoldt(:,i) + gradpsi(:,i)
endif
endif
case(10,20) ! remapped B/rho or remapped B
dBevoldt(:,i) = 0.
case(1:9) ! evolving B/rho
if (iuse_exact_derivs.gt.0) then
dBevoldt(:,i) = sqrtg(i)*dBevoldt(:,i)*rho1i
!--add the B/rho dot grad v bit
!if (any(dBevoldt(:,i).gt.0.)) then
!print*,i,'dBevol/dt = ',dBevoldt(:,i)
!endif
do k=1,ndimV
dBevoldt(k,i) = dBevoldt(k,i) + dot_product(Bevol(1:ndim,i),dveldx(1:ndim,k,i))
enddo
!if (any(dBevoldt(:,i).gt.0.)) then
!print*,i,'dBevol/dt (exact) = ',dBevoldt(:,i)
!read*
!endif
else
dBevoldt(:,i) = sqrtg(i)*dBevoldt(:,i)*rho1i
if (idivBzero.ge.2) then
gradpsi(:,i) = gradpsi(:,i)*rho1i**2
endif
endif
case(-1) ! vector potential evolution, crap gauge
!
!--add the v x B term
!
call cross_product3D(vel(:,i),Bfield(:,i),curlBi)
dBevoldt(:,i) = dBevoldt(:,i)*rho1i + curlBi(:)
case(-2) ! vector potential evolution, Axel gauge
if (iuse_exact_derivs.gt.0) then
do k=1,ndim
dBevoldt(k,i) = -dot_product(Bevol(1:ndimV,i),dveldx(k,1:ndimV,i))
enddo
else
dBevoldt(:,i) = dBevoldt(:,i)*rho1i
endif
!
!--get v x Bext
!
call cross_product3D(vel(:,i),Bconst(:),curlBi)
!--add v x Bext plus the existing term, which includes dissipation
dBevoldt(:,i) = dBevoldt(:,i) + curlBi(:)
!
!--add dissipation for vector potential = -eta J
!
if (iav.ge.2) then
curlBi(:) = curlB(:,i)
!curlBi(:) = curlBsym(:,i)*rhoi
!curlB(:,i) = curlBi(:)
!curr2 = abs(dot_product(curlBsym(:,i)*rhoi,curlBsym(:,i)*rhoi))
!curr2 = abs(dot_product(curlBi,curlBsym(:,i)*rhoi))
curr2 = dot_product(curlBi,curlBi)
etai = alpha(3,i)*sqrt(valfven2i)*hh(i)
!etai = alpha(3,i)*vsig*hh(i)
!etai = alpha(3,i)*2.*sqrt(curr2/rho(i))*hh(i)**2
!etai = alpha(3,i)*sqrt(dot_product(graddivv(:,i),graddivv(:,i)))*hh(i)**2
!etai = alpha(3,i)*(sqrt(valfven2i)*hh(i) + 2.*sqrt(curr2/rho(i))*hh(i)**2)
!divB(i) = etai
dBevoldt(:,i) = dBevoldt(:,i) - etai*curlBi(:)
!if (curlBi(1).ne.0.) print*,i,curlB(1,i),curlBi(1)
!print*,i,'diss = ',-etai,curlB(:,i)
!
!--add term to thermal energy equation. This is used to construct
! dendt for the entropy/total energy equations below
!
dudt(i) = dudt(i) + etai*curr2/rho(i)
endif
case(-3) ! Generalised Euler Potentials evolution
dBevoldt(:,i) = 0.
case default
dBevoldt(:,i) = 0.
end select
!
!--calculate resistive timestep (bootstrap onto force timestep)
!
if (iresist.gt.0 .and. iresist.ne.2 .and. etamhd.gt.tiny(etamhd)) then
if (iresist.eq.3) then
etai = etafunc(x(1,i),etamhd)
elseif (iresist > 0) then
etai = etamhd
endif
dtforce = min(dtforce,hh(i)**2/etai)
endif
!
!--if using the thermal energy equation, set the energy derivative
! (note that dissipative terms are calculated in rates, but otherwise comes straight from cty)
!
if (iener.eq.3) then
! could do this in principle but does not work with
! faniso modified by subtraction of Bconst
if (damp.lt.tiny(0.)) dudt(i) = dudt(i) + pr(i)*rho1i**2*drhodt(i)
dendt(i) = dot_product(vel(:,i),fprev(:,i)) + dudt(i) !&
! + 0.5*(dot_product(Bfield(:,i),Bfield(:,i))*rho1i**2) &
! + dot_product(Bfield(:,i),dBevoldt(:,i))*rho1i
elseif (iener.eq.1) then ! entropy variable (just dissipative terms)
if (damp.lt.tiny(0.)) dendt(i) = dendt(i) + (gamma-1.)/dens(i)**(gamma-1.)*dudt(i)
elseif (iener.eq.4) then
dendt(i) = (pr(i) + en(i))*rho1i*drhodt(i) + rhoi*dudt(i)
if (iav.lt.0) stop 'iener=4 not compatible with Godunov SPH'
elseif (iener.gt.0 .and. iav.ge.0 .and. iener.ne.10 .and. iener.ne.11 .and. &
(idust.ne.1 .and. idust.ne.3 .and. idust.ne.4)) then
if (damp.lt.tiny(0.)) dudt(i) = dudt(i) + pr(i)*rho1i**2*drhodt(i)
dendt(i) = dudt(i)
else
dendt(i) = dudt(i)
endif
if (itype(i).eq.itypedust) dendt(i) = 0.
!
!--calculate time derivative of alpha (artificial dissipation coefficients)
! see Morris and Monaghan (1997) and Price and Monaghan (2004c)
!
daldt(:,i) = 0.
if (any(iavlim.ne.0)) then
tdecay1 = (avdecayconst*vsig)/hh(i) ! 1/decay time (use vsig)
!
!--artificial viscosity parameter
!
select case(iavlim(1))
case(1,2)
source = max(drhodt(i)*rho1i,0.0)
if (iavlim(1).eq.2) source = source*(2.0-alpha(1,i))
daldt(1,i) = (alphamin - alpha(1,i))*tdecay1 + avfact*source
case(3)
graddivvmag = sqrt(dot_product(graddivv(:,i),graddivv(:,i)))
!!print*,'graddivvmag = ',graddivvmag,max(drhodt(i)*rho1i,0.0)
source = hh(i)*graddivvmag*(2.0-alpha(1,i))
daldt(1,i) = (alphamin - alpha(1,i))*tdecay1 + avfact*source
end select
!
!--artificial thermal conductivity parameter
!
if (iener.gt.0 .and. iavlim(2).gt.0) then
!--this is using h*/sqrt(u)*(del^2 u) as the source
if (uu(i).gt.epsilon(uu(i))) then
sourceu = hh(i)*abs(del2u(i))/sqrt(uu(i))
!sourceu = sqrt(abs(del2u(i))) !!
else
sourceu = 0.
endif
daldt(2,i) = (alphaumin - alpha(2,i))*tdecay1 + sourceu
endif
!
!--artificial resistivity parameter
!
if (iavlim(3).ne.0 .and. imhd.ne.0) then
!--calculate source term for the resistivity parameter
sourceJ = SQRT(DOT_PRODUCT(curlB(:,i),curlB(:,i))*rho1i)
! B2i = DOT_PRODUCT(Bfield(:,i),Bfield(:,i))
! sourceJ = SQRT(vsig2*DOT_PRODUCT(curlB(:,i),curlB(:,i))/B2i)
sourcedivB = 10.*abs(divB(i))*SQRT(rho1i)
sourceB = MAX(sourceJ,sourcedivB)
!sourceB = max(sqrt(drhodt(i)*rho1i*sourceJ),0.0)
if (iavlim(3).eq.2) then
sourceB = sourceB*(2.0-alpha(3,i))
elseif (iavlim(3).eq.3) then
source = max(drhodt(i)*rho1i,0.0)*(2.-alpha(3,i))
!graddivv(:,i) = graddivv(:,i)/rho(i)
!sourceJ = sqrt(dot_product(graddivv(:,i),graddivv(:,i)))
!sourceB = sourceJ/(abs(drhodt(i)*rho1i) + sourceJ + tiny(0.))*sourceB*(2.-alpha(3,i))
sourceB = sqrt(source*sourceB)
endif
daldt(3,i) = (alphaBmin - alpha(3,i))*tdecay1 + sourceB
endif
endif
!
!--calculate time derivative of divergence correction parameter psi
!
select case(idivBzero)
case(2:7)
dpsidt(i) = -vsig2max*divB(i) - psidecayfact*psi(i)*vsigmax/hh(i)
!if (abs(psi(i)).gt.0.) print*,' idivBzero',i,vsig2max,divB(i),psi(i)
case DEFAULT
dpsidt(i) = 0.
end select
if (idust.eq.1) then
!
!--DEBUGGING: CHECK ENERGY CONSERVATION
! BY ADDING TERMS (SHOULD GIVE ZERO)
!
esum = esum + pmass(i)*(dot_product(vel(:,i),force(:,i)) &
- dot_product(vel(:,i),fexternal(:)) &
+ rhogasi*rhodusti(1)*rho1i**2*dot_product(deltav(:,i),ddeltavdt(:,i)) &
+ ((1. - 2.*dustfraci(1))*0.5*dot_product(deltav(:,i),deltav(:,i)) - uu(i))*ddustevoldt(1,i) &
+ rhogasi*rho1i*dudt(i))
elseif (idust.eq.4) then
esum = esum + pmass(i)*(dot_product(vel(:,i),force(:,i)) &
- dot_product(vel(:,i),fexternal(:)) &
- uu(i)*ddustevoldt(1,i) &
+ (1. - dustfrac(1,i))*dudt(i))
endif
enddo
if ((idust.eq.1 .or. idust.eq.4) .and. abs(esum).gt.epsilon(esum) .and. (iener.ge.2)) then
print*,' SUM (should be zero if conserving energy) = ',esum
!read*
endif
if (idust.ne.0 .and. ratio > 1.) print "(a,g8.3,a)",' WARNING: max ts/dt = ',ratio,' approximation not valid'
if (ivisc.gt.0) print*,' dEk/dt = ',esum
if (sqrt(dot_product(fmean,fmean)).gt.1.e-8 .and. mod(nsteps,100).eq.0) print*,'WARNING: fmean = ',fmean(:)
!
!--calculate timestep constraint from the forces
! dtforce is returned together with dtcourant to the main timestepping loop
!
if (fhmax.lt.0.) then
write(iprint,*) 'rates: fhmax <=0 :',fhmax
call quit
elseif (fhmax.gt.0.) then
dtforce = min(dtforce,sqrt(1./fhmax))
endif
!!print*,'dtcourant = ',dtcourant,dtcourant2,0.2*dtcourant2
!!dtcourant = 0.2*dtcourant2
!
!--set rates to zero on ghosts/fixed particles
!
do i=1,ntotal ! using ntotal just makes sure they are zero for ghosts
if (itype(i).eq.itypebnd .or. itype(i).eq.itypebnddust .or. i.gt.npart) then
force(:,i) = 0.0
drhodt(i) = 0.0
dhdt(i) = 0.0
dudt(i) = 0.0
dendt(i) = 0.0
dBevoldt(:,i) = 0.0
daldt(:,i) = 0.
dpsidt(i) = 0.0
fmag(:,i) = 0.0
divB(i) = 0.0
if (imhd.ge.0) curlB(:,i) = 0.0
xsphterm(:,i) = 0.0
gradpsi(:,i) = 0.0
endif
enddo
if (allocated(listneigh)) deallocate(listneigh)
if (allocated(phi)) deallocate(phi,del2u)
if (trace) write(iprint,*) ' Exiting subroutine get_rates'
if (itiming) then
call cpu_time(t5)
write(iprint,"(50('-'))")
write(iprint,*) 'time for intro = ',t2-t1,'s'
write(iprint,*) 'time for main = ',t3-t2,'s'
write(iprint,*) 'time for gravity = ',t4-t3,'s'
write(iprint,*) 'time for final = ',t5-t4,'s'
write(iprint,*) 'total rates time = ',t5-t1,'s'
write(iprint,"(50('-'))")
endif
!--------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------
contains
real function slope_limiter(sl,sr,ilimiter) result(s)
real, intent(in) :: sl,sr
integer, intent(in) :: ilimiter
s = 0.
select case(ilimiter)
case(6)
if (sl*sr > 0.) s = 0.5*(sl + sr)
! van Albada
!s = (sl + sl*sr)/(1. + sl*sr)
case(5)
! UMIST