-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathsetup_choi.f90
111 lines (104 loc) · 3.38 KB
/
setup_choi.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
!------------------------------------------------------------------------------!
! 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 Choi damping test
! Choi, Kim & Witta (2009), ApJS 181, 413
!------------------------------------------
subroutine setup
!
!--include relevant global variables
!
use dimen_mhd
use debug
use loguns
use bound
use options
use part
use setup_params
use eos
use uniform_distributions
!
!--define local variables
!
implicit none
integer :: i,ny
real :: massp,volume,totmass
real :: denszero,ampl,vA
real, dimension(3) :: Bzero
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' entering subroutine setup(choi)'
!
!--set boundaries
!
ibound = 3 ! boundaries
nbpts = 0 ! use ghosts not fixed
xmin(:) = 0. ! set position of boundaries
xmax(:) = 1.
ny = 6
if (ndim >= 2) xmax(2) = xmin(2) + ny*psep
if (ndim >= 3) xmax(3) = xmin(3) + ny*psep
gamma = 1.
ampl = 0.1
gamma_ambipolar = 1000.
rho_ion = 0.1
!
!--set up the uniform density grid
!
Bzero(1) = 1.
Bzero(2) = 0.
Bzero(3) = 0.
denszero = 1.0
vA = sqrt(dot_product(Bzero,Bzero)/denszero)
call set_uniform_cartesian(2,psep,xmin,xmax,adjustbound=.true.)
npart = ntotal
print*,'npart =',npart
!
!--determine particle mass
!
volume = product(xmax(:)-xmin(:))
totmass = denszero*volume
massp = totmass/float(ntotal) ! average particle mass
!
!--now assign particle properties
!
do i=1,npart
vel(:,i) = 0.
vel(3,i) = ampl*vA*sin(2.*pi*(x(1,i)-xmin(1)))
dens(i) = denszero
pmass(i) = massp
uu(i) = 1.0
Bfield(:,i) = Bzero
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