diff --git a/week10/nmap10.ipynb b/week10/nmap10.ipynb
new file mode 100644
index 0000000000000000000000000000000000000000..2fb021b2be9b8ce28418602575216e3325c2b9cb
--- /dev/null
+++ b/week10/nmap10.ipynb
@@ -0,0 +1,791 @@
+{
+ "cells": [
+  {
+   "cell_type": "markdown",
+   "id": "incomplete-medline",
+   "metadata": {},
+   "source": [
+    "<h1 style=\"text-align: center; vertical-align: middle;\">Numerical Methods in Accelerator Physics</h1>\n",
+    "<h2 style=\"text-align: center; vertical-align: middle;\">Python examples -- Week 10</h2>"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "executive-television",
+   "metadata": {},
+   "source": [
+    "<h2>Run this first!</h2>\n",
+    "\n",
+    "Imports and modules:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "biological-product",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "from config10 import (np, plt, tqdm, trange, c, epsilon_0,\n",
+    "                    beta, gamma, Machine, track_one_turn,\n",
+    "                    charge, mass, emittance, hamiltonian, U,\n",
+    "                    plot_hamiltonian, plot_rf_overview,\n",
+    "                    plot_dist)\n",
+    "%matplotlib inline"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "a11d88ee-9907-42d9-bd4d-9115f8e7873e",
+   "metadata": {},
+   "source": [
+    "If the progress bar by `tqdm` (`trange`) in this document does not work, run this:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "9bb4244f-11d1-4a18-b695-cf5e4de54b88",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "!jupyter nbextension enable --py widgetsnbextension"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "21fcb2eb-7252-46c5-b101-9ca9e539c248",
+   "metadata": {},
+   "source": [
+    "<h2>Example: Longitudinal Space Charge in CERN PS</h2>\n",
+    "\n",
+    "Consider CERN PS for illustration once again:\n",
+    "- non-accelerating rf bucket at $V_{rf}=25$ kV\n",
+    "- injection energy $\\gamma=3.13$\n",
+    "- bunch of $N_p=1\\times 10^{12}$ particles\n",
+    "- rms bunch length of $\\sigma_z=2$ m\n",
+    "- vacuum pipe of $\\approx 10$ cm\n",
+    "- beam radius of $\\approx 1$ cm"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "e9db48a4-c79a-4b0b-9b63-33f502efb440",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "m = Machine(gamma_ref=3.13, phi_s=0, voltage=25000)\n",
+    "#m = Machine(gamma_ref=7, phi_s=np.pi, voltage=25000)\n",
+    "# adjust here later!"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "dff4b590-3874-4164-a7ad-5afa105d6cc5",
+   "metadata": {},
+   "source": [
+    "The rf harmonic determines the width of an rf bucket:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "97dc0aed-769d-4cd8-b20a-37a45af03fd0",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "m.harmonic"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "676a2bf0-03d2-4779-9bc1-92b9bda738e6",
+   "metadata": {},
+   "source": [
+    "Can you determine the location of the left and right unstable fix points, which limit the rf bucket?"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "e3f7175c-addc-4ca6-92e2-0750d02b71b3",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "rfbucket_z_right = m.circumference/m.harmonic/2\n",
+    "rfbucket_z_left = -m.circumference/m.harmonic/2"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "a40193ab-713c-44d9-9ecf-7a8ba1b8b5ff",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plot_rf_overview(m);"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "2d93e85a-4b73-4ce1-8b6c-c384524ee1e2",
+   "metadata": {},
+   "source": [
+    "We initialize a Gaussian bunch distribution with small rms bunch length $\\sigma_z=2$ m (so that the small-amplitude approximation applies):"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "02622e47-ed37-4362-be90-4c61c51a4840",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "sigma_z = 2"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "e8af2e54-5b29-473e-b43f-65ff7d60db5d",
+   "metadata": {},
+   "source": [
+    "The matched rms momentum spread $\\sigma_{\\Delta p}$ (remember, $\\sigma_{\\Delta p}$ and $\\sigma_z$ are linked via equal Hamiltonian values, the equilibrium condition):"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "1d2c4ca9-73bb-491f-b19d-9a0590a70b90",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "sigma_deltap = np.sqrt(\n",
+    "    2 * m.p0() / np.abs(m.eta(0)) * \n",
+    "    charge * m.voltage * np.pi * m.harmonic / (beta(gamma(m.p0())) * c * m.circumference**2)\n",
+    ") * sigma_z"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "a99772e3-3b32-4607-9b23-309cf03fcf8f",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "sigma_deltap / m.p0()"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "04751471-2160-4308-ab06-a724d9270fe0",
+   "metadata": {},
+   "source": [
+    "Number of simulated macro-particles (not too low such as to well resolve the line charge density of the bunch):"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "ebe73b63-47e2-4b6f-8f79-f5b1c22fe6e6",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "N = 100000"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "5383dd25-88e6-4c2b-931a-10c66c82dd82",
+   "metadata": {},
+   "source": [
+    "Generation of macro-particles by initialising the phase-space distributio"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "51c17b2c-63f0-456c-bddb-dc306d96f213",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "np.random.seed(12345)\n",
+    "\n",
+    "z = np.random.normal(loc=0, scale=sigma_z, size=N)\n",
+    "deltap = np.random.normal(loc=0, scale=sigma_deltap, size=N)"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "bb77203a-dd91-477a-8da4-234880e2c9a2",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plot_hamiltonian(m, dpmax=0.005)\n",
+    "plt.scatter(z[::N // 1000], deltap[::N // 1000] / m.p0(), marker='.', s=1);"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "31f3c420-9ccd-4be9-95e5-34f8ed67f539",
+   "metadata": {},
+   "source": [
+    "Next, we define a regular 1D grid for discretisation -- we will slice up the beam into bins (the slices):"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "341788c9-04ae-4824-8814-4ed0e620c2e7",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "N_slices = 100\n",
+    "\n",
+    "slice_boundaries = np.linspace(rfbucket_z_left, rfbucket_z_right, N_slices + 1, endpoint=True)"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "4e40b5a4-2ebe-488c-8984-58cf1b16a55a",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "slice_centres = (slice_boundaries[1:] + slice_boundaries[:-1]) / 2"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "1e3dc474-8d39-426d-81e5-cbf572169365",
+   "metadata": {},
+   "source": [
+    "And we are now in the position to compute a discretised line charge density $\\lambda(z)$ from the $z$ distribution of the macro-particles:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "4d690835-26a7-48c0-bebe-29b2dd4b716b",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "N_per_slice, slice_boundaries = np.histogram(z, bins=slice_boundaries)"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "e16706b9-9c91-43be-8cdb-a71329c3848c",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plt.step(slice_centres, N_per_slice, where='mid')\n",
+    "plt.xlabel('$z$ [m]')\n",
+    "plt.ylabel('Number of particles / slice');"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "49262b41-dcdd-45f4-a865-807dbf6ce9fc",
+   "metadata": {},
+   "source": [
+    "The bunch has a total intensity of $N_p=1\\times 10^{12}$ real particles. What total charge does a macro-particle carry?"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "8383e875-c944-4ac4-969b-2052fb1e93d3",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "charge_per_mp = 1e12 / N * charge"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "a4d4e0af-0e25-482d-91fe-a21eec19e92b",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "lmbda = charge_per_mp * N_per_slice / np.diff(slice_boundaries)"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "944e2c0d-c3c7-4047-9068-ecd3718fcd1f",
+   "metadata": {},
+   "source": [
+    "We evaluate the geometry factor $g$:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "0f7d1e55-e107-4419-b07d-e90bba346393",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "r_pipe = 0.1\n",
+    "r_beam = 0.01\n",
+    "\n",
+    "g = 1 + 2 * np.log(r_pipe / r_beam)"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "73961e95-9322-4b3f-adde-964af5516893",
+   "metadata": {},
+   "source": [
+    "If the space-charge field is given by $E_z\\propto \\frac{d\\lambda}{dz}$, the space-charge potential is given by the same expression but with $V_\\mathrm{sc}\\propto -\\lambda$!"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "fd354881-f501-407e-b20e-5dd18dba7266",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "V_sc = g / (4 * np.pi * epsilon_0 * m.gamma_ref**2) * lmbda"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "091df43e-12da-447b-bc70-538acc2cb5ba",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plt.plot(slice_centres, V_sc)\n",
+    "plt.xlabel('$z$ [m]')\n",
+    "plt.ylabel('$|V_{sc}|$ [V]');"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "788e9ebf-3e67-4e9d-8204-6366251602cd",
+   "metadata": {},
+   "source": [
+    "Compute the derivative of the line charge density via second-order (symmetric) finite difference,\n",
+    "\n",
+    "$\\lambda'(z)\\approx\\frac{\\lambda(z+\\Delta z)-\\lambda(z-\\Delta z)}{2\\Delta z}$"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "4cd27fda-e134-47e1-b0b7-a169d9ed04e0",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "lmbda_prime = np.gradient(lmbda, slice_centres)"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "dd01468b-f396-4433-933f-cfc45bf2af7b",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plt.step(slice_centres, 1e9 * lmbda_prime)\n",
+    "plt.xlabel('$z$ [m]')\n",
+    "plt.ylabel(\"$\\lambda'(z)$ [nC/m${}^2$]\");"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "d62d4c74-7a04-4ff2-bdec-53adf47e0e3f",
+   "metadata": {},
+   "source": [
+    "The longitudinal space-charge field is thus given by:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "2e3fbc84-b42d-48a5-8102-37b3f31df138",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "E_z = -g / (4 * np.pi * epsilon_0 * m.gamma_ref**2) * lmbda_prime"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "d83dc894-6819-495e-81df-74183fef8c07",
+   "metadata": {},
+   "source": [
+    "<h3>Comparison with external phase focusing from rf</h3>\n",
+    "Compute the effective external RF field as a mean via the rf wave kick per circumference:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "f5058fdb-f929-4b8c-8a79-9250edd11795",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "E_rf = m.voltage / m.circumference * np.sin(\n",
+    "    -m.harmonic * slice_centres * 2 * np.pi / m.circumference + m.phi_s)\n",
+    "# E_rf = -np.gradient(U(slice_centres, m) / charge * c, slice_centres)"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "a08b4c37-b9cc-4b73-a62d-65862fa70443",
+   "metadata": {},
+   "source": [
+    "The external phase focusing by the rf wave and the space-charge field sum up:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "97ab1cc6-ea52-41b6-add1-555f5f396810",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plt.plot(slice_centres, E_rf, c='k', label='RF')\n",
+    "plt.plot(slice_centres, E_z, c='r', label='SC')\n",
+    "plt.plot(slice_centres, E_rf + E_z, c='lightblue', label='RF+SC')\n",
+    "\n",
+    "plt.legend()\n",
+    "plt.xlabel('$z$ [m]')\n",
+    "plt.ylabel('mean $E$ [V/m]');"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "df150e5e-e918-4118-84b7-c4e57d635e90",
+   "metadata": {},
+   "source": [
+    "<h2>Microwave instability</h2>\n",
+    "We implement the space-charge kick for the tracking:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "4b6d8e70-8369-4212-ab81-6263997935c2",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "def space_charge_kick(z_n, deltap_n, machine, charge_per_mp=charge_per_mp, N_slices=N_slices, g=g):\n",
+    "    m = machine\n",
+    "    rfbucket_z_right = m.circumference / m.harmonic / 2\n",
+    "    rfbucket_z_left = -rfbucket_z_right\n",
+    "    \n",
+    "    # slicing\n",
+    "    slice_boundaries = np.linspace(rfbucket_z_left, rfbucket_z_right, N_slices + 1, endpoint=True)\n",
+    "    N_per_slice, _ = np.histogram(z_n, bins=slice_boundaries)\n",
+    "    # find slice index for each particle:\n",
+    "    slice_idx_per_mp = np.floor((z_n - slice_boundaries[0]) / \n",
+    "                                (slice_boundaries[1] - slice_boundaries[0])).astype(int)\n",
+    "    \n",
+    "    # get line charge density derivative\n",
+    "    lmbda = charge_per_mp * N_per_slice / np.diff(slice_boundaries)\n",
+    "    lmbda_prime = np.gradient(lmbda, slice_boundaries[:-1])\n",
+    "    \n",
+    "    # space-charge field per slice\n",
+    "    E_z = -g / (4 * np.pi * epsilon_0 * m.gamma_ref**2) * lmbda_prime\n",
+    "    \n",
+    "    # momentum update\n",
+    "    deltap_n1 = deltap_n - charge * E_z[slice_idx_per_mp] * m.circumference / (beta(m.gamma_ref) * c)\n",
+    "    return deltap_n1"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "d540691a-ae7f-44cb-bb22-56713e1c80e5",
+   "metadata": {},
+   "source": [
+    "The kicks which the particles receive look as follows:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "a3ed02d7-f0a9-43bb-9b2d-f973ded1d6f9",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plt.scatter(z, (space_charge_kick(z, deltap, m) - deltap) / m.p0())\n",
+    "plt.xlabel('$z$ [m]')\n",
+    "plt.ylabel('$(\\Delta p_{n+1} - \\Delta p_{n}) / p_0$');"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "417f9c34-7ffc-4621-80db-dac75feffbdc",
+   "metadata": {},
+   "source": [
+    "We implement the one-turn tracking:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "42b21b20-196d-4e42-8413-eb85e88ff46b",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "def track_one_turn(z_n, deltap_n, machine, **kwargs):\n",
+    "    m = machine\n",
+    "    # half drift\n",
+    "    z_nhalf = z_n - m.eta(deltap_n) * deltap_n / m.p0() * m.circumference / 2\n",
+    "    # rf kick\n",
+    "    amplitude = charge * m.voltage / (beta(gamma(m.p0())) * c)\n",
+    "    phi = m.phi_s - m.harmonic * 2 * np.pi * z_nhalf / m.circumference\n",
+    "    \n",
+    "    m.update_gamma_ref()\n",
+    "    deltap_n1 = deltap_n + amplitude * (np.sin(phi) - np.sin(m.phi_s))\n",
+    "    # space-charge kick\n",
+    "    deltap_n1 = space_charge_kick(z, deltap_n1, m, **kwargs)\n",
+    "    # half drift\n",
+    "    z_n1 = z_nhalf - m.eta(deltap_n1) * deltap_n1 / m.p0() * m.circumference / 2\n",
+    "    return z_n1, deltap_n1"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "5bb0ac97-a935-45eb-a229-8c617e6ea8d4",
+   "metadata": {},
+   "source": [
+    "We simulate just above transition, continuing from the end of the previous part."
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "b257377d-e733-44e5-9d2f-612d0dd1b8b6",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "assert m.gamma_ref > 1 / np.sqrt(m.alpha_c), 'Initialise the machine above transition energy!'"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "b722f9ae-96fb-4f17-8f29-2a0b2717d070",
+   "metadata": {},
+   "source": [
+    "We simulate for a bit more than 1 synchrotron period above transition at 10x higher intensity than the initial nominal parameters:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "145673db-5edd-437d-98f8-6ec40f6978bc",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "n_turns = 3000\n",
+    "\n",
+    "intensity_factor = 10"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "02b88c5d-1717-4cde-920c-7d3977897bd8",
+   "metadata": {},
+   "source": [
+    "... and record the longitudinal emittance:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "54c56f03-cfc7-4154-8118-591ab6689d1b",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "epsn_z = np.zeros(n_turns, dtype=np.float64)\n",
+    "\n",
+    "epsn_z[0] = emittance(z, deltap)"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "87b981c1-1556-4358-ac3e-86620acd02e9",
+   "metadata": {},
+   "source": [
+    "... as well as the bunch profiles (just as we did in the tomography lecture 09):"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "ef5a773d-ec22-4018-b291-d61b2a4dc0e6",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "profiles = [np.histogram(z, bins=slice_boundaries)[0]]"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "ad4e1857-1473-40f2-8cf3-12d11f093865",
+   "metadata": {},
+   "source": [
+    "Let's track! (May take a couple of minutes!)"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "0f881e19-a9a8-4df5-9ffe-1dfad5027c53",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "for i_turn in trange(1, n_turns):\n",
+    "    z, deltap = track_one_turn(z, deltap, m, \n",
+    "                               charge_per_mp=charge_per_mp * intensity_factor, \n",
+    "                               N_slices=400)\n",
+    "    \n",
+    "    # record:\n",
+    "    epsn_z[i_turn] = emittance(z, deltap)\n",
+    "    profiles += [np.histogram(z, bins=slice_boundaries)[0]]"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "0b8a936b-d4b4-4cbe-b37b-836259f7336c",
+   "metadata": {},
+   "source": [
+    "Let us have a look at the final bunch profile:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "24759fb3-3e84-416b-af13-78e73d0ab289",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "N_per_slice, slice_boundaries = np.histogram(z, bins=slice_boundaries)\n",
+    "\n",
+    "plt.step(slice_centres, N_per_slice, where='mid')\n",
+    "plt.xlabel('$z$ [m]')\n",
+    "plt.ylabel('Number of particles / slice');"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "3892214c-3772-4602-b984-25bca4cd47ce",
+   "metadata": {},
+   "source": [
+    "Let's plot the emittance evolution:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "5bbfe4d1-b709-44a7-8306-aae2b50cd08f",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plt.plot(np.arange(n_turns), 100 * (epsn_z - epsn_z[0]) / epsn_z[0])\n",
+    "\n",
+    "plt.xlabel('Turns')\n",
+    "plt.ylabel('$\\Delta \\epsilon_z/\\epsilon_{z0}$ [%]');"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "5a90979d-e8eb-4e1f-86be-fbb245f4bb31",
+   "metadata": {},
+   "source": [
+    "Look into phase space and a couple of particles:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "2f965ba3-1ab8-41e9-aacc-55da67153600",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plot_hamiltonian(m, dpmax=0.005)\n",
+    "plt.scatter(z[::N // 1000], deltap[::N // 1000] / m.p0(), marker='.', s=1);"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "16f81006-fec6-4b88-9ca6-8c593749da23",
+   "metadata": {},
+   "source": [
+    "A density plot of the longitudinal phase space shows the damage:\n"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "a39f0636-d15d-4039-9009-735da0b6cda7",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plt.hist2d(z, deltap/m.p0(), bins=100)\n",
+    "\n",
+    "plt.xlabel('$z$ [m]')\n",
+    "plt.ylabel('$\\delta$');"
+   ]
+  },
+  {
+   "cell_type": "markdown",
+   "id": "d372e6d2-f86d-4ec8-885a-cbb66cd8a687",
+   "metadata": {},
+   "source": [
+    "And the bunch profiles over time:"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "577fd2c8-7de9-4f3e-a1c1-6a4040fe2587",
+   "metadata": {},
+   "outputs": [],
+   "source": [
+    "plt.imshow(profiles, origin='lower', cmap='jet', extent=(rfbucket_z_left, rfbucket_z_right, 0, n_turns))\n",
+    "plt.gca().set_aspect(np.diff(plt.xlim()) / np.diff(plt.ylim()))\n",
+    "plt.xlabel('$z$ [m]')\n",
+    "plt.ylabel('Turn')\n",
+    "plt.colorbar(label='Density');"
+   ]
+  },
+  {
+   "cell_type": "code",
+   "execution_count": null,
+   "id": "04e1687f-8412-475d-91e4-f1f2102f3360",
+   "metadata": {},
+   "outputs": [],
+   "source": []
+  }
+ ],
+ "metadata": {
+  "kernelspec": {
+   "display_name": "Python 3 (ipykernel)",
+   "language": "python",
+   "name": "python3"
+  },
+  "language_info": {
+   "codemirror_mode": {
+    "name": "ipython",
+    "version": 3
+   },
+   "file_extension": ".py",
+   "mimetype": "text/x-python",
+   "name": "python",
+   "nbconvert_exporter": "python",
+   "pygments_lexer": "ipython3",
+   "version": "3.9.20"
+  }
+ },
+ "nbformat": 4,
+ "nbformat_minor": 5
+}