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

Added 1D vacuum expansion test.

parent a0f62f65
Branches
Tags
1 merge request!508Evrard and other tests
###############################################################################
# This file is part of SWIFT.
# Copyright (c) 2018 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/>.
#
##############################################################################
# Generates an overdensity within a vacuum to test the vacuum handling
# capabilities of the code
import numpy as np
import h5py
fileName = "vacuum.hdf5"
numPart = 100
boxSize = 1.
gamma = 5. / 3.
coords = np.zeros((numPart, 3))
v = np.zeros((numPart, 3))
m = np.zeros(numPart)
h = np.zeros(numPart)
u = np.zeros(numPart)
ids = np.arange(numPart, dtype = 'L')
rho = np.zeros(numPart)
# first set the positions, as we try to do a reasonable volume estimate to
# set the masses
for i in range(numPart):
# we only generate particles in the range [0.25, 0.75]
coords[i,0] = 0.25 + 0.5 * (i + 0.5) / numPart
rho[i] = 1.
P = 1.
u[i] = P / (gamma - 1.) / rho[i]
m[i] = rho[i] * 0.5 / numPart
# reasonable smoothing length estimate
h[i] = 1. / numPart
#File
file = h5py.File(fileName, 'w')
# Header
grp = file.create_group("/Header")
grp.attrs["BoxSize"] = boxSize
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]
grp.attrs["Flag_Entropy_ICs"] = 0
grp.attrs["Dimension"] = 1
#Runtime parameters
grp = file.create_group("/RuntimePars")
grp.attrs["PeriodicBoundariesOn"] = 1
#Units
grp = file.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 = file.create_group("/PartType0")
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')
grp.create_dataset('Density', data=rho, dtype='f')
file.close()
###############################################################################
# This file is part of SWIFT.
# Copyright (c) 2018 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 matplotlib
matplotlib.use("Agg")
import pylab as pl
import h5py
import sys
# Parameters
gamma = 5. / 3. # Polytropic index
rhoL = 1. # Initial density in the non vacuum state
vL = 0. # Initial velocity in the non vacuum state
PL = 1. # Initial pressure in the non vacuum state
rhoR = 0. # Initial vacuum density
vR = 0. # Initial vacuum velocity
PR = 0. # Initial vacuum pressure
# Plot parameters
params = {'axes.labelsize': 10,
'axes.titlesize': 10,
'font.size': 12,
'legend.fontsize': 12,
'xtick.labelsize': 10,
'ytick.labelsize': 10,
'text.usetex': True,
'figure.figsize' : (9.90,6.45),
'figure.subplot.left' : 0.045,
'figure.subplot.right' : 0.99,
'figure.subplot.bottom' : 0.05,
'figure.subplot.top' : 0.99,
'figure.subplot.wspace' : 0.15,
'figure.subplot.hspace' : 0.12,
'lines.markersize' : 6,
'lines.linewidth' : 3.,
'text.latex.unicode': True
}
pl.rcParams.update(params)
pl.rc('font',**{'family':'sans-serif','sans-serif':['Times']})
# Read the snapshot index from the command line argument
snap = int(sys.argv[1])
# Open the file and read the relevant data
file = h5py.File("vacuum_{0:04d}.hdf5".format(snap), "r")
x = file["/PartType0/Coordinates"][:,0]
rho = file["/PartType0/Density"]
v = file["/PartType0/Velocities"][:,0]
u = file["/PartType0/InternalEnergy"]
S = file["/PartType0/Entropy"]
P = file["/PartType0/Pressure"]
time = file["/Header"].attrs["Time"][0]
scheme = file["/HydroScheme"].attrs["Scheme"]
kernel = file["/HydroScheme"].attrs["Kernel function"]
neighbours = file["/HydroScheme"].attrs["Kernel target N_ngb"][0]
eta = file["/HydroScheme"].attrs["Kernel eta"][0]
git = file["Code"].attrs["Git Revision"]
# Get the analytic solution, which is just the solution of the corresponding
# vacuum Riemann problem evaluated at the correct time
# left state sound speed (and rarefaction wave speed)
aL = np.sqrt(gamma * PL / rhoL)
# vacuum front speed
SL = vL + 2. / (gamma - 1.) * aL
# we evaluate the solution centred on 0., and shift to the correct position
# afterwards
xa = np.arange(-0.25, 0.25, 0.001)
rhoa = np.zeros(len(xa))
va = np.zeros(len(xa))
Pa = np.zeros(len(xa))
for i in range(len(xa)):
dxdt = xa[i] / time
if dxdt > vL - aL:
if dxdt < SL:
# rarefaction regime
# factor that appears in both the density and pressure expression
fac = 2. / (gamma + 1.) + \
(gamma - 1.) / (gamma + 1.) * (vL - dxdt) / aL
rhoa[i] = rhoL * fac**(2. / (gamma - 1.))
va[i] = 2. / (gamma + 1.) * (aL + 0.5 * (gamma - 1.) * vL + dxdt)
Pa[i] = PL * fac**(2. * gamma / (gamma - 1.))
else:
# vacuum regime
rhoa[i] = 0.
va[i] = 0.
Pa[i] = 0.
else:
# left state regime
rhoa[i] = rhoL
va[i] = vL
Pa[i] = PL
ua = Pa / (gamma - 1.) / rhoa
Sa = Pa / rhoa**gamma
# Plot the interesting quantities
fig, ax = pl.subplots(2, 3)
# Velocity profile
ax[0][0].plot(x, v, "r.", markersize = 4.)
ax[0][0].plot(xa + 0.75, va, "k--", alpha = 0.8, linewidth = 1.2)
ax[0][0].plot(xa + 0.25, -va[::-1], "k--", alpha = 0.8, linewidth = 1.2)
ax[0][0].set_xlabel("${\\rm{Position}}~x$", labelpad = 0)
ax[0][0].set_ylabel("${\\rm{Velocity}}~v_x$", labelpad = 0)
# Density profile
ax[0][1].plot(x, rho, "r.", markersize = 4.)
ax[0][1].plot(xa + 0.75, rhoa, "k--", alpha = 0.8, linewidth = 1.2)
ax[0][1].plot(xa + 0.25, rhoa[::-1], "k--", alpha = 0.8, linewidth = 1.2)
ax[0][1].set_xlabel("${\\rm{Position}}~x$", labelpad = 0)
ax[0][1].set_ylabel("${\\rm{Density}}~\\rho$", labelpad = 0)
# Pressure profile
ax[0][2].plot(x, P, "r.", markersize = 4.)
ax[0][2].plot(xa + 0.75, Pa, "k--", alpha = 0.8, linewidth = 1.2)
ax[0][2].plot(xa + 0.25, Pa[::-1], "k--", alpha = 0.8, linewidth = 1.2)
ax[0][2].set_xlabel("${\\rm{Position}}~x$", labelpad = 0)
ax[0][2].set_ylabel("${\\rm{Pressure}}~P$", labelpad = 0)
# Internal energy profile
ax[1][0].plot(x, u, "r.", markersize = 4.)
ax[1][0].plot(xa + 0.75, ua, "k--", alpha = 0.8, linewidth = 1.2)
ax[1][0].plot(xa + 0.25, ua[::-1], "k--", alpha = 0.8, linewidth = 1.2)
ax[1][0].set_xlabel("${\\rm{Position}}~x$", labelpad = 0)
ax[1][0].set_ylabel("${\\rm{Internal~Energy}}~u$", labelpad = 0)
# Entropy profile
ax[1][1].plot(x, S, "r.", markersize = 4.)
ax[1][1].plot(xa + 0.75, Sa, "k--", alpha = 0.8, linewidth = 1.2)
ax[1][1].plot(xa + 0.25, Sa[::-1], "k--", alpha = 0.8, linewidth = 1.2)
ax[1][1].set_xlabel("${\\rm{Position}}~x$", labelpad = 0)
ax[1][1].set_ylabel("${\\rm{Entropy}}~S$", labelpad = 0)
# Run information
ax[1][2].set_frame_on(False)
ax[1][2].text(-0.49, 0.9,
"Vacuum test with $\\gamma={0:.3f}$ in 1D at $t = {1:.2f}$".format(
gamma, time), fontsize = 10)
ax[1][2].text(-0.49, 0.8,
"Left:~~ $(P_L, \\rho_L, v_L) = ({0:.3f}, {1:.3f}, {2:.3f})$".format(
PL, rhoL, vL), fontsize = 10)
ax[1][2].text(-0.49, 0.7,
"Right: $(P_R, \\rho_R, v_R) = ({0:.3f}, {1:.3f}, {2:.3f})$".format(
PR, rhoR, vR), fontsize = 10)
ax[1][2].plot([-0.49, 0.1], [0.62, 0.62], "k-", lw = 1)
ax[1][2].text(-0.49, 0.5, "$\\textsc{{Swift}}$ {0}".format(git), fontsize = 10)
ax[1][2].text(-0.49, 0.4, scheme, fontsize = 10)
ax[1][2].text(-0.49, 0.3, kernel, fontsize = 10)
ax[1][2].text(-0.49, 0.2,
"${0:.2f}$ neighbours ($\\eta={1:.3f}$)".format(neighbours, eta),
fontsize = 10)
ax[1][2].set_xlim(-0.5, 0.5)
ax[1][2].set_ylim(0., 1.)
ax[1][2].set_xticks([])
ax[1][2].set_yticks([])
pl.tight_layout()
pl.savefig("Vacuum.png".format(snap))
#!/bin/bash
# Generate the initial conditions if they are not present.
if [ ! -e vacuum.hdf5 ]
then
echo "Generating initial conditions for the 1D vacuum expansion example..."
python makeIC.py
fi
# Run SWIFT
../swift -s -t 1 vacuum.yml 2>&1 | tee output.log
# Plot the result
python plotSolution.py 1
# 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: 0.1 # The end time of the simulation (in internal units).
dt_min: 1e-10 # 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: vacuum # Common part of the name of output files
time_first: 0. # Time of the first output (in internal units)
delta_time: 0.1 # 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: ./vacuum.hdf5 # The file to read
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment