-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathsetup_hopkins_disc2D.f90
112 lines (103 loc) · 3.53 KB
/
setup_hopkins_disc2D.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
!------------------------------------------------------------------------------!
! 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 a uniform density cartesian grid of particles in ND
!----------------------------------------------------------------
subroutine setup
!
!--include relevant global variables
!
use dimen_mhd
use debug
use loguns
use bound
use options
use part
use setup_params
use eos, only:gamma,polyk
use uniform_distributions
use externf, only:eps2_soft
!
!--define local variables
!
implicit none
integer :: i
real :: massp,volume,totmass
real :: denszero,r,phi,gamm1,omega
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' entering subroutine setup(unifdis)'
!
!--set boundaries
!
ibound = 0 ! boundaries
nbpts = 0 ! use ghosts not fixed
xmin(:) = -2. ! set position of boundaries
xmax(:) = 2.
iexternal_force = 2
if (ndim /= 2) stop 'need ndim=2'
call set_uniform_cartesian(1,psep,xmin,xmax,rmin=0.5,rmax=2.)
npart = ntotal
print*,'npart =',npart
!
!--determine particle mass
!
denszero = 1.0
volume = product(xmax(:)-xmin(:))
totmass = denszero*volume
massp = totmass/float(ntotal) ! average particle mass
gamm1 = gamma - 1.
!
!--now assign particle properties
!
do i=1,ntotal
vel(:,i) = 0.
dens(i) = denszero
Bfield(:,i) = 0.
r = sqrt(dot_product(x(:,i),x(:,i)))
phi = atan2(x(2,i),x(1,i))
!
!--keplerian velocity profile balancing pressure gradient
!
omega = sqrt(1./(r**2 + eps2_soft)**1.5)
vel(1,i) = -r*sin(phi)*omega
vel(2,i) = r*cos(phi)*omega
if (gamm1 > 0.) then
uu(i) = polyk/gamm1*dens(i)**gamm1
else
uu(i) = 1.e-5
endif
pmass(i) = massp
enddo
!
!--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