Skip to content
Snippets Groups Projects
Commit 954b572d authored by Matthieu Schaller's avatar Matthieu Schaller
Browse files

Added the 'Square test' in 2D.

parent 352f3fa0
Branches
Tags
No related merge requests found
......@@ -27,7 +27,7 @@ Statistics:
SPH:
resolution_eta: 1.2348 # Target smoothing length in units of the mean inter-particle separation (1.2348 == 48Ngbs with the cubic spline kernel).
delta_neighbours: 0.1 # The tolerance for the targetted number of neighbours.
max_smoothing_length: 0.1 # Maximal smoothing length allowed (in internal units).
max_smoothing_length: 0.02 # Maximal smoothing length allowed (in internal units).
CFL_condition: 0.1 # Courant-Friedrich-Levy condition for time integration.
# Parameters related to the initial conditions
......
......@@ -19,7 +19,6 @@
import h5py
from numpy import *
import sys
# Generates a swift IC file for the Gresho-Chan vortex in a periodic box
......@@ -80,7 +79,7 @@ fileOutput = h5py.File(fileOutputName, 'w')
# Header
grp = fileOutput.create_group("/Header")
grp.attrs["BoxSize"] = boxSize
grp.attrs["BoxSize"] = [boxSize, boxSize, 0.2]
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]
......
###############################################################################
# This file is part of SWIFT.
# Copyright (c) 2016 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
from numpy import *
# Generates a swift IC file for the Square test in a periodic box
# Parameters
L = 64 # Number of particles on the side
gamma = 5./3. # Gas adiabatic index
rho0 = 4 # Gas central density
rho1 = 1 # Gas outskirt density
P0 = 2.5 # Gas central pressure
P1 = 2.5 # Gas central pressure
vx = 142.3 # Random velocity for all particles
vy = -31.4
fileOutputName = "square.hdf5"
#---------------------------------------------------
vol = 1.
numPart_out = L * L
numPart_in = L * L * rho0 / rho1 / 4
L_out = int(sqrt(numPart_out))
L_in = int(sqrt(numPart_in))
pos_out = zeros((numPart_out, 3))
for i in range(L_out):
for j in range(L_out):
index = i * L_out + j
pos_out[index, 0] = i / (float(L_out)) + 1./(2. * L_out)
pos_out[index, 1] = j / (float(L_out)) + 1./(2. * L_out)
h_out = ones(numPart_out) * (1. / L_out) * 1.2348
m_out = ones(numPart_out) * vol * rho1 / numPart_out
u_out = ones(numPart_out) * P1 / (rho1 * (gamma - 1.))
pos_in = zeros((numPart_in, 3))
for i in range(L_in):
for j in range(L_in):
index = i * L_in + j
pos_in[index, 0] = 0.25 + i / float(2. * L_in) + 1./(2. * 2. * L_in)
pos_in[index, 1] = 0.25 + j / float(2. * L_in) + 1./(2. * 2. * L_in)
h_in = ones(numPart_in) * (1. / L_in) * 1.2348
m_in = ones(numPart_in) * 0.25 * vol * rho0 / numPart_in
u_in = ones(numPart_in) * P0 / (rho0 * (gamma - 1.))
# Remove the central particles
select_out = logical_or(logical_or(pos_out[:,0] < 0.25 , pos_out[:,0] > 0.75), logical_or(pos_out[:,1] < 0.25, pos_out[:,1] > 0.75))
pos_out = pos_out[select_out, :]
h_out = h_out[select_out]
u_out = u_out[select_out]
m_out = m_out[select_out]
# Add the central region
pos = append(pos_out, pos_in, axis=0)
h = append(h_out, h_in, axis=0)
u = append(u_out, u_in)
m = append(m_out, m_in)
numPart = size(h)
ids = linspace(1, numPart, numPart)
vel = zeros((numPart, 3))
vel[:,0] = vx
vel[:,1] = vy
#File
fileOutput = h5py.File(fileOutputName, 'w')
# Header
grp = fileOutput.create_group("/Header")
grp.attrs["BoxSize"] = [vol, vol, 0.2]
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, 0, 0, 0, 0, 0]
#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")
ds = grp.create_dataset('Coordinates', (numPart, 3), 'd')
ds[()] = pos
ds = grp.create_dataset('Velocities', (numPart, 3), 'f')
ds[()] = vel
ds = grp.create_dataset('Masses', (numPart, 1), 'f')
ds[()] = m.reshape((numPart,1))
ds = grp.create_dataset('SmoothingLength', (numPart,1), 'f')
ds[()] = h.reshape((numPart,1))
ds = grp.create_dataset('InternalEnergy', (numPart,1), 'f')
ds[()] = u.reshape((numPart,1))
ds = grp.create_dataset('ParticleIDs', (numPart,1), 'L')
ds[()] = ids.reshape((numPart,1))
fileOutput.close()
###############################################################################
# This file is part of SWIFT.
# Copyright (c) 2016 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/>.
#
##############################################################################
# Computes the analytical solution of the square test
# Parameters
gas_gamma = 5./3. # Gas adiabatic index
gamma = 5./3. # Gas adiabatic index
rho0 = 4 # Gas central density
rho1 = 1 # Gas outskirt density
P0 = 2.5 # Gas central pressure
P1 = 2.5 # Gas central pressure
vx = 142.3 # Random velocity for all particles
vy = -31.4
# ---------------------------------------------------------------
# Don't touch anything after this.
# ---------------------------------------------------------------
import matplotlib
matplotlib.use("Agg")
from pylab import *
import h5py
# 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
}
rcParams.update(params)
rc('font',**{'family':'sans-serif','sans-serif':['Times']})
snap = int(sys.argv[1])
# Read the simulation data
sim = h5py.File("square_%03d.hdf5"%snap, "r")
boxSize = sim["/Header"].attrs["BoxSize"][0]
time = sim["/Header"].attrs["Time"][0]
scheme = sim["/HydroScheme"].attrs["Scheme"]
kernel = sim["/HydroScheme"].attrs["Kernel function"]
neighbours = sim["/HydroScheme"].attrs["Kernel target N_ngb"]
eta = sim["/HydroScheme"].attrs["Kernel eta"]
git = sim["Code"].attrs["Git Revision"]
# Analytical soltion
centre_x = 0.5 + time * vx
centre_y = 0.5 + time * vy
while centre_x > 1.:
centre_x -= 1.
while centre_x < 0.:
centre_x += 1.
while centre_y > 1.:
centre_y -= 1.
while centre_y < 0.:
centre_y += 1.
pos = sim["/PartType0/Coordinates"][:,:]
vel = sim["/PartType0/Velocities"][:,:]
v_norm = sqrt(vel[:,0]**2 + vel[:,1]**2)
rho = sim["/PartType0/Density"][:]
u = sim["/PartType0/InternalEnergy"][:]
S = sim["/PartType0/Entropy"][:]
P = sim["/PartType0/Pressure"][:]
x = pos[:,0] - centre_x
y = pos[:,1] - centre_y
# Box wrapping
x[x>0.5] -= 1.
x[x<-0.5] += 1.
y[y>0.5] -= 1.
y[y<-0.5] += 1.
# Azimuthal velocity profile -----------------------------
subplot(231)
scatter(x, y, c=v_norm, cmap="PuBu", edgecolors='face', s=4, vmin=-1., vmax=1.)
text(0.47, 0.47, "${\\rm{Velocity~norm}}$", ha="right", va="top", backgroundcolor="w")
plot([-0.25, 0.25], [0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
plot([-0.25, 0.25], [-0.25, -0.25], '--', color='k', alpha=0.8, lw=2)
plot([0.25, 0.25], [-0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
plot([-0.25, -0.25], [-0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
xlabel("${\\rm{Position}}~x$", labelpad=0)
ylabel("${\\rm{Position}}~y$", labelpad=-7)
xlim(-0.5, 0.5)
ylim(-0.5, 0.5)
# Radial density profile --------------------------------
subplot(232)
scatter(x, y, c=rho, cmap="PuBu", edgecolors='face', s=4, vmin=0., vmax=4.)
text(0.47, 0.47, "${\\rm{Density}}$", ha="right", va="top", backgroundcolor="w")
plot([-0.25, 0.25], [0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
plot([-0.25, 0.25], [-0.25, -0.25], '--', color='k', alpha=0.8, lw=2)
plot([0.25, 0.25], [-0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
plot([-0.25, -0.25], [-0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
xlabel("${\\rm{Position}}~x$", labelpad=0)
ylabel("${\\rm{Position}}~y$", labelpad=-7)
xlim(-0.5, 0.5)
ylim(-0.5, 0.5)
# Radial pressure profile --------------------------------
subplot(233)
scatter(x, y, c=P, cmap="PuBu", edgecolors='face', s=4, vmin=2, vmax=4)
text(0.47, 0.47, "${\\rm{Pressure}}$", ha="right", va="top", backgroundcolor="w")
plot([-0.25, 0.25], [0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
plot([-0.25, 0.25], [-0.25, -0.25], '--', color='k', alpha=0.8, lw=2)
plot([0.25, 0.25], [-0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
plot([-0.25, -0.25], [-0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
xlabel("${\\rm{Position}}~x$", labelpad=0)
ylabel("${\\rm{Position}}~y$", labelpad=-7)
xlim(-0.5, 0.5)
ylim(-0.5, 0.5)
# Internal energy profile --------------------------------
subplot(234)
scatter(x, y, c=u, cmap="PuBu", edgecolors='face', s=4, vmin=0.5, vmax=4.)
text(0.47, 0.47, "${\\rm{Internal~energy}}$", ha="right", va="top", backgroundcolor="w")
plot([-0.25, 0.25], [0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
plot([-0.25, 0.25], [-0.25, -0.25], '--', color='k', alpha=0.8, lw=2)
plot([0.25, 0.25], [-0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
plot([-0.25, -0.25], [-0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
xlabel("${\\rm{Position}}~x$", labelpad=0)
ylabel("${\\rm{Position}}~y$", labelpad=-7)
xlim(-0.5, 0.5)
ylim(-0.5, 0.5)
# Radial entropy profile --------------------------------
subplot(235)
scatter(x, y, c=S, cmap="PuBu", edgecolors='face', s=4, vmin=0., vmax=3.)
text(0.47, 0.47, "${\\rm{Entropy}}$", ha="right", va="top", backgroundcolor="w")
plot([-0.25, 0.25], [0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
plot([-0.25, 0.25], [-0.25, -0.25], '--', color='k', alpha=0.8, lw=2)
plot([0.25, 0.25], [-0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
plot([-0.25, -0.25], [-0.25, 0.25], '--', color='k', alpha=0.8, lw=2)
xlabel("${\\rm{Position}}~x$", labelpad=0)
ylabel("${\\rm{Position}}~y$", labelpad=-7)
xlim(-0.5, 0.5)
ylim(-0.5, 0.5)
# Information -------------------------------------
subplot(236, frameon=False)
text(-0.49, 0.9, "Square test with $\\gamma=%.3f$ at $t=%.2f$"%(gas_gamma,time), fontsize=10)
text(-0.49, 0.8, "Centre:~~~ $(P, \\rho) = (%.3f, %.3f)$"%(P0, rho0), fontsize=10)
text(-0.49, 0.7, "Outskirts: $(P, \\rho) = (%.3f, %.3f)$"%(P1, rho1), fontsize=10)
plot([-0.49, 0.1], [0.62, 0.62], 'k-', lw=1)
text(-0.49, 0.5, "$\\textsc{Swift}$ %s"%git, fontsize=10)
text(-0.49, 0.4, scheme, fontsize=10)
text(-0.49, 0.3, kernel, fontsize=10)
text(-0.49, 0.2, "$%.2f$ neighbours ($\\eta=%.3f$)"%(neighbours, eta), fontsize=10)
xlim(-0.5, 0.5)
ylim(0, 1)
xticks([])
yticks([])
savefig("SquareTest.png", dpi=200)
#!/bin/bash
# Generate the initial conditions if they are not present.
if [ ! -e square.hdf5 ]
then
echo "Generating initial conditions for the square test ..."
python makeIC.py
fi
# Run SWIFT
../swift -s -t 1 square.yml
# Plot the solution
python plotSolution.py 40
# 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: square # Common part of the name of output files
time_first: 0. # Time of the first output (in internal units)
delta_time: 1e-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).
delta_neighbours: 0.1 # The tolerance for the targetted number of neighbours.
max_smoothing_length: 0.02 # Maximal smoothing length allowed (in internal units).
CFL_condition: 0.1 # Courant-Friedrich-Levy condition for time integration.
# Parameters related to the initial conditions
InitialConditions:
file_name: ./square.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