makeIC.py 3.77 KB
Newer Older
1
2
###############################################################################
 # This file is part of SWIFT.
3
 # Copyright (c) 2012 Pedro Gonnet (pedro.gonnet@durham.ac.uk),
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
 #                    Matthieu Schaller (matthieu.schaller@durham.ac.uk)
 # 
 # This program is free software: you can redistribute it and/or modify
 # it under the terms of the GNU Lesser General Public License as published
 # by the Free Software Foundation, either version 3 of the License, or
 # (at your option) any later version.
 # 
 # This program is distributed in the hope that it will be useful,
 # but WITHOUT ANY WARRANTY; without even the implied warranty of
 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 # GNU General Public License for more details.
 # 
 # You should have received a copy of the GNU Lesser General Public License
 # along with this program.  If not, see <http://www.gnu.org/licenses/>.
 # 
 ##############################################################################

import h5py
22
import random
23
24
25
26
27
28
from numpy import *

# Generates a swift IC file for the Sedov blast test in a periodic cubic box

# Parameters
periodic= 1      # 1 For periodic box
29
boxSize = 10.
30
31
32
L = 101           # Number of particles along one axis
rho = 1.          # Density
P = 1.e-5         # Pressure
33
E0= 1.e2          # Energy of the explosion
34
pert = 0.1
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
gamma = 5./3.     # Gas adiabatic index
fileName = "sedov.hdf5" 


#---------------------------------------------------
numPart = L**3
mass = boxSize**3 * rho / numPart
internalEnergy = P / ((gamma - 1.)*rho)

if L%2 == 0:
    print "Number of particles along each dimension must be odd."
    exit()

#Generate particles
coords = zeros((numPart, 3))
v      = zeros((numPart, 3))
51
52
53
54
m      = zeros(numPart)
h      = zeros(numPart)
u      = zeros(numPart)
ids    = zeros(numPart, dtype='L')
55
56
57
58
59
60
61
62
63
64
65

for i in range(L):
    for j in range(L):
        for k in range(L):
            index = i*L*L + j*L + k
            x = i * boxSize / L + boxSize / (2*L)
            y = j * boxSize / L + boxSize / (2*L)
            z = k * boxSize / L + boxSize / (2*L)
            coords[index,0] = x
            coords[index,1] = y
            coords[index,2] = z
66
            v[index,0] = 0.
67
68
69
            v[index,1] = 0.
            v[index,2] = 0.
            m[index] = mass
70
            h[index] = 1.1255 * boxSize / L
71
            u[index] = internalEnergy
72
            ids[index] = index + 1
73
            if sqrt((x - boxSize/2.)**2 + (y - boxSize/2.)**2 + (z - boxSize/2.)**2) < 2.01 * boxSize/L:
74
                u[index] = u[index] + E0 / (33. * mass)
75
                print "Particle " , index , " set to detonate."
76
77
78
79
            coords[index,0] = x + random.random() * pert * boxSize/(2.*L)
            coords[index,1] = y + random.random() * pert * boxSize/(2.*L)
            coords[index,2] = z + random.random() * pert * boxSize/(2.*L)

80
81
82
83
84
85
86
87
88

#--------------------------------------------------

#File
file = h5py.File(fileName, 'w')

# Header
grp = file.create_group("/Header")
grp.attrs["BoxSize"] = boxSize
89
90
91
92
93
94
grp.attrs["NumPart_Total"] =  [numPart, 0, 0, 0, 0, 0]
grp.attrs["NumPart_Total_HighWord"] = [0, 0, 0, 0, 0, 0]
grp.attrs["NumPart_ThisFile"] = [numPart, 0, 0, 0, 0, 0]
grp.attrs["Time"] = 0.0
grp.attrs["NumFilesPerSnapshot"] = 1
grp.attrs["MassTable"] = [0.0, 0.0, 0.0, 0.0, 0.0, 0.0]
95
grp.attrs["Flag_Entropy_ICs"] = 0
96
97
98
99
100
101
102

#Runtime parameters
grp = file.create_group("/RuntimePars")
grp.attrs["PeriodicBoundariesOn"] = periodic

#Particle group
grp = file.create_group("/PartType0")
103
104
105
106
107
108
grp.create_dataset('Coordinates', data=coords, dtype='d')
grp.create_dataset('Velocities', data=v, dtype='f')
grp.create_dataset('Masses', data=m, dtype='f')
grp.create_dataset('SmoothingLength', data=h, dtype='f')
grp.create_dataset('InternalEnergy', data=u, dtype='f')
grp.create_dataset('ParticleIDs', data=ids, dtype='L')
109
110

file.close()