Skip to content
Snippets Groups Projects
Commit 0f3987e4 authored by Bert Vandenbroucke's avatar Bert Vandenbroucke
Browse files

Added 2D Kelvin Helmholtz growth rate test.

parent 718e8946
No related branches found
No related tags found
1 merge request!508Evrard and other tests
#!/bin/bash
wget http://virgodb.cosma.dur.ac.uk/swift-webstorage/ICs/glassPlane_128.hdf5
# Define the system of units to use internally.
InternalUnitSystem:
UnitMass_in_cgs: 1 # Grams
UnitLength_in_cgs: 1 # Centimeters
UnitVelocity_in_cgs: 1 # Centimeters per second
UnitCurrent_in_cgs: 1 # Amperes
UnitTemp_in_cgs: 1 # Kelvin
# Parameters governing the time integration
TimeIntegration:
time_begin: 0. # The starting time of the simulation (in internal units).
time_end: 4. # The end time of the simulation (in internal units).
dt_min: 1e-6 # The minimal time-step size of the simulation (in internal units).
dt_max: 1e-2 # The maximal time-step size of the simulation (in internal units).
# Parameters governing the snapshots
Snapshots:
basename: kelvinHelmholtzGrowthRate # Common part of the name of output files
time_first: 0. # Time of the first output (in internal units)
delta_time: 0.04 # Time difference between consecutive outputs (in internal units)
# Parameters governing the conserved quantities statistics
Statistics:
delta_time: 1e-2 # Time between statistics output
# Parameters for the hydrodynamics scheme
SPH:
resolution_eta: 1.2348 # Target smoothing length in units of the mean inter-particle separation (1.2348 == 48Ngbs with the cubic spline kernel).
CFL_condition: 0.1 # Courant-Friedrich-Levy condition for time integration.
# Parameters related to the initial conditions
InitialConditions:
file_name: ./kelvinHelmholtzGrowthRate.hdf5 # The file to read
################################################################################
# This file is part of SWIFT.
# Copyright (c) 2016 Matthieu Schaller (matthieu.schaller@durham.ac.uk)
# 2017 Bert Vandenbroucke (bert.vandenbroucke@gmail.com)
#
# 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
import numpy as np
# Generates a swift IC file for the Kelvin-Helmholtz vortex in a periodic box
# Parameters
gamma = 5./3. # Gas adiabatic index
P0 = 2.5 # Pressure
rho0 = 1. # Density
d = 0.0317 # Thickness of the transition layer
B = 0.0005 # Amplitude of the seed velocity
fileOutputName = "kelvinHelmholtzGrowthRate.hdf5"
#---------------------------------------------------
glass = h5py.File("glassPlane_128.hdf5", 'r')
pos = glass["/PartType0/Coordinates"][:, :]
h = glass["/PartType0/SmoothingLength"][:]
N = len(h)
vol = 1.
# Generate extra arrays
v = np.zeros((N, 3))
ids = np.linspace(1, N, N)
m = np.ones(N) * rho0 * vol / N
u = np.ones(N) * P0 / (rho0 * (gamma - 1.))
v[pos[:, 1] <= 0.25 - d, 0] = -0.5
v[(pos[:, 1] < 0.25 + d) & (pos[:, 1] > 0.25 - d), 0] = \
-0.5 + \
0.5 * (pos[(pos[:, 1] < 0.25 + d) & (pos[:, 1] > 0.25 - d), 1] + d - 0.25) / d
v[(pos[:, 1] <= 0.75 - d) & (pos[:, 1] >= 0.25 + d), 0] = 0.5
v[(pos[:, 1] < 0.75 + d) & (pos[:, 1] > 0.75 - d), 0] = \
0.5 - \
0.5 * (pos[(pos[:, 1] < 0.75 + d) & (pos[:, 1] > 0.75 - d), 1] + d - 0.75) / d
v[pos[:, 1] >= 0.75 + d, 0] = -0.5
v[:, 1] = B * np.sin(4. * np.pi * pos[:, 0]) * \
(np.exp(-(pos[:, 1] - 0.25)**2 / 32. / d**2) + \
np.exp(-(pos[:, 1] - 0.75)**2 / 32. / d**2))
#File
fileOutput = h5py.File(fileOutputName, 'w')
# Header
grp = fileOutput.create_group("/Header")
grp.attrs["BoxSize"] = [1., 1., 1.]
grp.attrs["NumPart_Total"] = [N, 0, 0, 0, 0, 0]
grp.attrs["NumPart_Total_HighWord"] = [0, 0, 0, 0, 0, 0]
grp.attrs["NumPart_ThisFile"] = [N, 0, 0, 0, 0, 0]
grp.attrs["Time"] = 0.0
grp.attrs["NumFileOutputsPerSnapshot"] = 1
grp.attrs["MassTable"] = [0.0, 0.0, 0.0, 0.0, 0.0, 0.0]
grp.attrs["Flag_Entropy_ICs"] = [0, 0, 0, 0, 0, 0]
grp.attrs["Dimension"] = 2
#Runtime parameters
grp = fileOutput.create_group("/RuntimePars")
grp.attrs["PeriodicBoundariesOn"] = 1
#Units
grp = fileOutput.create_group("/Units")
grp.attrs["Unit length in cgs (U_L)"] = 1.
grp.attrs["Unit mass in cgs (U_M)"] = 1.
grp.attrs["Unit time in cgs (U_t)"] = 1.
grp.attrs["Unit current in cgs (U_I)"] = 1.
grp.attrs["Unit temperature in cgs (U_T)"] = 1.
#Particle group
grp = fileOutput.create_group("/PartType0")
grp.create_dataset("Coordinates", data = pos, 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')
fileOutput.close()
################################################################################
# This file is part of SWIFT.
# Copyright (c) 2016 Matthieu Schaller (matthieu.schaller@durham.ac.uk)
# 2017 Bert Vandenbroucke (bert.vandenbroucke@gmail.com)
#
# 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
import numpy as np
# Generates a swift IC file for the Kelvin-Helmholtz vortex in a periodic box
# Parameters
gamma = 5./3. # Gas adiabatic index
P0 = 2.5 # Pressure
rho0 = 1. # Density
d = 0.0317 # Thickness of the transition layer
B = 0.0005 # Amplitude of the seed velocity
N1D = 128 # Number of particles in one dimension
fileOutputName = "kelvinHelmholtzGrowthRate.hdf5"
#---------------------------------------------------
N = N1D ** 2
x = np.linspace(0., 1., N1D + 1)
x = 0.5 * (x[1:] + x[:-1])
y = x
xx, yy = np.meshgrid(x, y)
pos = np.zeros((N, 3))
pos[:, 0] = xx.reshape((N))
pos[:, 1] = yy.reshape((N))
h = np.ones(N) * 2. / N1D
vol = 1.
# Generate extra arrays
v = np.zeros((N, 3))
ids = np.linspace(1, N, N)
m = np.ones(N) * rho0 * vol / N
u = np.ones(N) * P0 / (rho0 * (gamma - 1.))
v[pos[:, 1] <= 0.25 - d, 0] = -0.5
v[(pos[:, 1] < 0.25 + d) & (pos[:, 1] > 0.25 - d), 0] = \
-0.5 + \
0.5 * (pos[(pos[:, 1] < 0.25 + d) & (pos[:, 1] > 0.25 - d), 1] + d - 0.25) / d
v[(pos[:, 1] <= 0.75 - d) & (pos[:, 1] >= 0.25 + d), 0] = 0.5
v[(pos[:, 1] < 0.75 + d) & (pos[:, 1] > 0.75 - d), 0] = \
0.5 - \
0.5 * (pos[(pos[:, 1] < 0.75 + d) & (pos[:, 1] > 0.75 - d), 1] + d - 0.75) / d
v[pos[:, 1] >= 0.75 + d, 0] = -0.5
v[:, 1] = B * np.sin(4. * np.pi * pos[:, 0]) * \
(np.exp(-(pos[:, 1] - 0.25)**2 / 32. / d**2) + \
np.exp(-(pos[:, 1] - 0.75)**2 / 32. / d**2))
#File
fileOutput = h5py.File(fileOutputName, 'w')
# Header
grp = fileOutput.create_group("/Header")
grp.attrs["BoxSize"] = [1., 1., 1.]
grp.attrs["NumPart_Total"] = [N, 0, 0, 0, 0, 0]
grp.attrs["NumPart_Total_HighWord"] = [0, 0, 0, 0, 0, 0]
grp.attrs["NumPart_ThisFile"] = [N, 0, 0, 0, 0, 0]
grp.attrs["Time"] = 0.0
grp.attrs["NumFileOutputsPerSnapshot"] = 1
grp.attrs["MassTable"] = [0.0, 0.0, 0.0, 0.0, 0.0, 0.0]
grp.attrs["Flag_Entropy_ICs"] = [0, 0, 0, 0, 0, 0]
grp.attrs["Dimension"] = 2
#Runtime parameters
grp = fileOutput.create_group("/RuntimePars")
grp.attrs["PeriodicBoundariesOn"] = 1
#Units
grp = fileOutput.create_group("/Units")
grp.attrs["Unit length in cgs (U_L)"] = 1.
grp.attrs["Unit mass in cgs (U_M)"] = 1.
grp.attrs["Unit time in cgs (U_t)"] = 1.
grp.attrs["Unit current in cgs (U_I)"] = 1.
grp.attrs["Unit temperature in cgs (U_T)"] = 1.
#Particle group
grp = fileOutput.create_group("/PartType0")
grp.create_dataset("Coordinates", data = pos, 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')
fileOutput.close()
###############################################################################
# This file is part of SWIFT.
# Copyright (c) 2017 Bert Vandenbroucke (bert.vandenbroucke@gmail.com)
#
# 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 numpy as np
import h5py
import matplotlib
matplotlib.use("Agg")
import pylab as pl
import sys
if len(sys.argv) < 2:
print "No final snapshot number provided!"
exit()
lastsnap = int(sys.argv[1])
# Read the simulation data
t = np.zeros(lastsnap + 1)
ey = np.zeros(lastsnap + 1)
for i in range(lastsnap + 1):
file = h5py.File("kelvinHelmholtzGrowthRate_{0:04d}.hdf5".format(i), 'r')
t_snap = float(file["/Header"].attrs["Time"])
vy = file["/PartType0/Velocities"][:,1]
m = file["/PartType0/Masses"][:]
ey_snap = 0.5 * m * vy**2
t[i] = t_snap
ey[i] = ey_snap.sum()
pl.semilogy(t, ey, "k.")
pl.savefig("kelvinHelmholtzGrowthRate.png")
#!/bin/bash
# Generate the initial conditions if they are not present.
if [ ! -e kelvinHelmholtzGrowthRate.hdf5 ]
then
echo "Generating initial conditions for the Kelvin-Helmholtz growth rate " \
"example..."
python makeIC.py
fi
# Run SWIFT
../swift -s -t 1 kelvinHelmholtzGrowthRate.yml 2>&1 | tee output.log
# Plot the solution
python plotSolution.py 100
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment