diff --git a/examples/Noh_2D/getGlass.sh b/examples/Noh_2D/getGlass.sh
new file mode 100755
index 0000000000000000000000000000000000000000..ae3c977064f5e7a408aa249c5fd9089b3c52ecb1
--- /dev/null
+++ b/examples/Noh_2D/getGlass.sh
@@ -0,0 +1,2 @@
+#!/bin/bash
+wget http://virgodb.cosma.dur.ac.uk/swift-webstorage/ICs/glassPlane_128.hdf5
diff --git a/examples/Noh_2D/makeIC.py b/examples/Noh_2D/makeIC.py
new file mode 100644
index 0000000000000000000000000000000000000000..9b8604b87768b69fd7c21793e9e5580625a90a82
--- /dev/null
+++ b/examples/Noh_2D/makeIC.py
@@ -0,0 +1,98 @@
+
+###############################################################################
+ # 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 2D Sod Shock in a periodic box
+
+# Parameters
+gamma = 5./3.          # Gas adiabatic index
+gamma = 5./3.      # Gas adiabatic index
+rho0 = 1.          # Background density
+P0 = 1.e-6         # Background pressure
+fileName = "noh.hdf5" 
+
+#---------------------------------------------------
+glass = h5py.File("glassPlane_128.hdf5", "r")
+
+vol = 4.
+
+pos = glass["/PartType0/Coordinates"][:,:] * sqrt(vol)
+h = glass["/PartType0/SmoothingLength"][:] * sqrt(vol)
+numPart = size(h)
+
+# Generate extra arrays
+v = zeros((numPart, 3))
+ids = linspace(1, numPart, numPart)
+
+m = zeros(numPart)
+u = zeros(numPart)
+m[:] = rho0 * vol / numPart    
+u[:] = P0 / (rho0 * (gamma - 1))
+
+# Make radial velocities
+#r = sqrt((pos[:,0]-1.)**2 + (pos[:,1]-1.)**2)
+#theta = arctan2((pos[:,1]-1.), (pos[:,0]-1.))
+v[:,0] = -(pos[:,0] - 1.)
+v[:,1] = -(pos[:,1] - 1.)
+
+norm_v = sqrt(v[:,0]**2 + v[:,1]**2)
+v[:,0] /= norm_v
+v[:,1] /= norm_v
+
+#File
+file = h5py.File(fileName, 'w')
+
+# Header
+grp = file.create_group("/Header")
+grp.attrs["BoxSize"] = [sqrt(vol), sqrt(vol), 1.]
+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"] = 2
+
+#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=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')
+
+
+file.close()
diff --git a/examples/Noh_2D/noh.yml b/examples/Noh_2D/noh.yml
new file mode 100644
index 0000000000000000000000000000000000000000..911ebbdc1ca82f9b5cf2f70841d848008cbc5f6c
--- /dev/null
+++ b/examples/Noh_2D/noh.yml
@@ -0,0 +1,35 @@
+# 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.6   # The end time of the simulation (in internal units).
+  dt_min:     1e-7  # The minimal time-step size of the simulation (in internal units).
+  dt_max:     1e-3  # The maximal time-step size of the simulation (in internal units).
+
+# Parameters governing the snapshots
+Snapshots:
+  basename:            noh # Common part of the name of output files
+  time_first:          0.    # Time of the first output (in internal units)
+  delta_time:          5e-2  # Time difference between consecutive outputs (in internal units)
+
+# Parameters governing the conserved quantities statistics
+Statistics:
+  delta_time:          1e-5 # 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.
+  CFL_condition:         0.1      # Courant-Friedrich-Levy condition for time integration.
+
+# Parameters related to the initial conditions
+InitialConditions:
+  file_name:  ./noh.hdf5          # The file to read
+
diff --git a/examples/Noh_2D/plotSolution.py b/examples/Noh_2D/plotSolution.py
new file mode 100644
index 0000000000000000000000000000000000000000..75b53fc76fb17a934946549d25f5895d6d79714c
--- /dev/null
+++ b/examples/Noh_2D/plotSolution.py
@@ -0,0 +1,179 @@
+###############################################################################
+ # 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 Noh problem and plots the SPH answer
+ 
+
+# Parameters
+gas_gamma = 5./3.      # Polytropic index
+rho0 = 1.              # Background density
+P0 = 1.e-6             # Background pressure
+v0 = 1
+
+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("noh_%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"]
+
+x = sim["/PartType0/Coordinates"][:,0]
+y = sim["/PartType0/Coordinates"][:,1]
+vx = sim["/PartType0/Velocities"][:,0]
+vy = sim["/PartType0/Velocities"][:,1]
+u = sim["/PartType0/InternalEnergy"][:]
+S = sim["/PartType0/Entropy"][:]
+P = sim["/PartType0/Pressure"][:]
+rho = sim["/PartType0/Density"][:]
+
+r = np.sqrt((x-1)**2 + (y-1)**2)
+v = -np.sqrt(vx**2 + vy**2)
+
+N = 1000  # Number of points
+x_min = -1
+x_max = 1
+
+#x += x_min
+
+# ---------------------------------------------------------------
+# Don't touch anything after this.
+# ---------------------------------------------------------------
+x_s = np.arange(0, 2., 2./N) - 1.
+rho_s = np.ones(N) * rho0
+P_s = np.ones(N) * rho0
+v_s = np.ones(N) * v0
+
+# Shock position
+u0 = rho0 * P0 * (gas_gamma-1)
+us = 0.5 * (gas_gamma - 1) * v0
+rs = us * time
+
+# Post-shock values
+rho_s[np.abs(x_s) < rs] = rho0 * ((gas_gamma + 1) / (gas_gamma - 1))**2
+P_s[np.abs(x_s) < rs] = 0.5 * rho0 * v0**2 * (gas_gamma + 1)**2 / (gas_gamma-1)
+v_s[np.abs(x_s) < rs] = 0.
+
+# Pre-shock values
+rho_s[np.abs(x_s) >= rs] = rho0 * (1 + v0 * time/np.abs(x_s[np.abs(x_s) >=rs]))
+P_s[np.abs(x_s) >= rs] = 0
+v_s[x_s >= rs] = -v0
+v_s[x_s <= -rs] = v0
+
+# Additional arrays
+u_s = P_s / (rho_s * (gas_gamma - 1.))  #internal energy
+s_s = P_s / rho_s**gas_gamma # entropic function
+
+# Plot the interesting quantities
+figure()
+
+# Velocity profile --------------------------------
+subplot(231)
+plot(r, v, '.', color='r', ms=4.0)
+plot(x_s, v_s, '--', color='k', alpha=0.8, lw=1.2)
+xlabel("${\\rm{Radius}}~r$", labelpad=0)
+ylabel("${\\rm{Velocity}}~v_r$", labelpad=-4)
+xlim(0, 0.5)
+ylim(-1.2, 0.4)
+
+# Density profile --------------------------------
+subplot(232)
+plot(r, rho, '.', color='r', ms=4.0)
+plot(x_s, rho_s, '--', color='k', alpha=0.8, lw=1.2)
+xlabel("${\\rm{Radius}}~r$", labelpad=0)
+ylabel("${\\rm{Density}}~\\rho$", labelpad=0)
+xlim(0, 0.5)
+ylim(0.95, 19)
+
+# Pressure profile --------------------------------
+subplot(233)
+plot(r, P, '.', color='r', ms=4.0)
+plot(x_s, P_s, '--', color='k', alpha=0.8, lw=1.2)
+xlabel("${\\rm{Radius}}~r$", labelpad=0)
+ylabel("${\\rm{Pressure}}~P$", labelpad=0)
+xlim(0, 0.5)
+ylim(-0.5, 11)
+
+# Internal energy profile -------------------------
+subplot(234)
+plot(r, u, '.', color='r', ms=4.0)
+plot(x_s, u_s, '--', color='k', alpha=0.8, lw=1.2)
+xlabel("${\\rm{Radius}}~r$", labelpad=0)
+ylabel("${\\rm{Internal~Energy}}~u$", labelpad=0)
+xlim(0, 0.5)
+ylim(-0.05, 0.8)
+
+# Entropy profile ---------------------------------
+subplot(235)
+plot(r, S, '.', color='r', ms=4.0)
+plot(x_s, s_s, '--', color='k', alpha=0.8, lw=1.2)
+xlabel("${\\rm{Radius}}~r$", labelpad=0)
+ylabel("${\\rm{Entropy}}~S$", labelpad=-9)
+xlim(0, 0.5)
+ylim(-0.05, 0.2)
+
+# Information -------------------------------------
+subplot(236, frameon=False)
+
+text(-0.49, 0.9, "Noh problem with  $\\gamma=%.3f$ in 2D at $t=%.2f$"%(gas_gamma,time), 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("Noh.png", dpi=200)
diff --git a/examples/Noh_2D/run.sh b/examples/Noh_2D/run.sh
new file mode 100755
index 0000000000000000000000000000000000000000..cff200801018e04ea560bd2c3fbd84057aec2d7c
--- /dev/null
+++ b/examples/Noh_2D/run.sh
@@ -0,0 +1,19 @@
+#!/bin/bash
+
+ # Generate the initial conditions if they are not present.
+if [ ! -e glassPlane_128.hdf5 ]
+then
+    echo "Fetching initial glass file for the Noh problem..."
+    ./getGlass.sh
+fi
+if [ ! -e noh.hdf5 ]
+then
+    echo "Generating initial conditions for the Noh problem..."
+    python makeIC.py
+fi
+
+# Run SWIFT
+../swift -s -t 2 noh.yml 2>&1 | tee output.log
+
+# Plot the solution
+python plotSolution.py 12