{ "cells": [ { "cell_type": "markdown", "source": [ "# Fan-beam tomography: short scan\n", "\n", "This page describes fan-beam tomographic image reconstruction\n", "using the Julia package\n", "[`Sinograms.jl`](https://github.com/JuliaImageRecon/Sinograms.jl)\n", "for a \"short\" scan\n", "where the rotation (`orbit`) is 180° plus the fan angle.\n", "This case requires\n", "[Parker weighting](http://doi.org/10.1118/1.595078).\n", "\n", "This page focuses on fan-beam with an \"arc\" detector,\n", "i.e., 3rd-generation CT systems\n", "where the detector is an arc\n", "having a focal point\n", "at the X-ray source.\n", "This geometry is popular\n", "in part because it facilitates\n", "anti-scatter grids.\n", "\n", "This page was generated from a single Julia file:\n", "[06-fan-short.jl](https://github.com/JuliaImageRecon/Sinograms.jl/blob/main/docs/lit/examples/06-fan-short.jl)." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "### Setup" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "Packages needed here." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using Plots: plot, gui # these 2 must precede Sinograms for Requires to work!\n", "using Unitful: mm\n", "using Sinograms: SinoFanArc, rays, plan_fbp, fbp, sino_geom_plot!\n", "using ImageGeoms: ImageGeom, fovs, MaskCircle\n", "using ImagePhantoms: SheppLogan, shepp_logan, radon, phantom\n", "using MIRTjim: jim, prompt" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "The following line is helpful when running this file as a script;\n", "this way it will prompt user to hit a key after each figure is displayed." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "isinteractive() ? jim(:prompt, true) : prompt(:draw);" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "## Fan-beam sinogram of Shepp-Logan phantom\n", "\n", "For illustration,\n", "we start by synthesizing\n", "a fan-beam sinogram\n", "of the Shepp-Logan phantom.\n", "\n", "For completeness,\n", "we use units (from Unitful),\n", "but units are optional." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "Use `ImageGeom` to define the image geometry." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "ig = ImageGeom(MaskCircle(); dims=(128,126), deltas = (2mm,2mm) )" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Use `SinoFanArc` to define the sinogram geometry,\n", "with the `:short` option for `orbit` to make a short scan.\n", "Note that even though we specify `na = 100`\n", "we end up with `na = 67` views\n", "because of the `:short` option." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "rg = SinoFanArc( :short, ;\n", " nb = 130, d = 3.2mm, na = 100, dsd = 400mm, dod = 140mm,\n", ")" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Examine the geometry to verify the FOV.\n", "The `na=67` blue dots show the `:short` scan." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "jim(axes(ig), ig.mask; prompt=false)\n", "sino_geom_plot!(rg, ig)" ], "metadata": {}, "execution_count": null }, { "outputs": [], "cell_type": "code", "source": [ "prompt()" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Ellipse parameters for Shepp-Logan phantom:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "μ = 0.01 / mm # typical linear attenuation coefficient\n", "ob = shepp_logan(SheppLogan(); fovs = fovs(ig), u = (1, 1, μ))" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Short arc fan-beam sinogram for Shepp-Logan phantom:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "sino = radon(rays(rg), ob)\n", "jim(axes(rg), sino; title=\"Shepp-Logan 'short' sinogram\",\n", " xlabel=\"r\", ylabel=\"ϕ\")" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "## Image reconstruction via FBP\n", "Here we start with a \"plan\",\n", "which would save work if we were reconstructing many images." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "plan = plan_fbp(rg, ig)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Examine Parker weights:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "jim(axes(rg), plan.view_weight; title = \"Parker weights\",\n", " xlabel=\"r\", ylabel=\"ϕ\")" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Finally perform FBP:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "fbp_image = fbp(plan, sino);" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "A narrow color window is needed to see the soft tissue structures:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "clim = (0.9, 1.1) .* μ\n", "jim(axes(ig), fbp_image, \"FBP image for 'short' arc case\"; clim)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "For comparison, here is the ideal phantom image:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "true_image = phantom(axes(ig)..., ob, 2)\n", "jim(axes(ig)..., true_image, \"True phantom image\"; clim)" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "Here is the difference image.\n", "Better sampling would reduce the errors." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "jim(axes(ig)..., fbp_image - true_image, \"Error image\")" ], "metadata": {}, "execution_count": null }, { "cell_type": "markdown", "source": [ "---\n", "\n", "*This notebook was generated using [Literate.jl](https://github.com/fredrikekre/Literate.jl).*" ], "metadata": {} } ], "nbformat_minor": 3, "metadata": { "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.8.4" }, "kernelspec": { "name": "julia-1.8", "display_name": "Julia 1.8.4", "language": "julia" } }, "nbformat": 4 }