diff --git a/examples/Supplementary code for Proc R Soc paper 'Integral equation methods for acoustic scattering by fractals'.ipynb b/examples/Supplementary code for Proc R Soc paper 'Integral equation methods for acoustic scattering by fractals'.ipynb index b1f1c0c..4e433f2 100644 --- a/examples/Supplementary code for Proc R Soc paper 'Integral equation methods for acoustic scattering by fractals'.ipynb +++ b/examples/Supplementary code for Proc R Soc paper 'Integral equation methods for acoustic scattering by fractals'.ipynb @@ -6,9 +6,9 @@ "source": [ "# Supplementary code for *Integral equation methods for acoustic scattering by fractals*\n", "\n", - "This notebook contains sample code which can be used to reproduce the field plots in Section 5(a) of the Proc. R. Soc. paper **Integral equation methods for acoustic scattering by fractals**, by *A. M. Caetano, S. N. Chandler-Wilde, X. Claeys, A. Gibbs, D. P. Hewett and A. Moiola*.\n", + "This notebook contains sample code which can be used to reproduce the 2D field plots in Section 5(a) of the Proc. R. Soc. paper **Integral equation methods for acoustic scattering by fractals**, by *A. M. Caetano, S. N. Chandler-Wilde, X. Claeys, A. Gibbs, D. P. Hewett and A. Moiola*.\n", "\n", - "All results in the above paper were produced using the open-source Julia package [IFSintegrals](https://github.com/AndrewGibbs/IFSintegrals). A technical understanding of this package is **not** required to understand the key ideas in this notebook. Throughout this notebook, we will make frequent references to equation numbers of the associated paper." + "All results in the above paper were produced using the open-source Julia package [IFSintegrals](https://github.com/AndrewGibbs/IFSintegrals). A technical understanding of this package is **not** required to understand the key ideas in this notebook. Throughout this notebook, we will make frequent references to equation and other numbers of the associated paper." ] }, { @@ -49,12 +49,13 @@ "cell_type": "markdown", "metadata": {}, "source": [ - "The following values of ```eg_index``` correspond to different scattering configurations, each in $\\mathbb{R}^2$:\n", - "1. The Koch curve (see Example 2.3)\n", - "2. Cantor set with contraction factor $\\rho=1/3$\n", - "3. Cantor dust with contraction factor $\\rho=1/3$\n", - "4. A Koch snowflake, constructed as the union of three Koch curves (see Fig. 8(b))\n", - "5. The union of a Cantor set with contraction factor $\\rho=1/3$ and a Koch curve.\n", + "The following values of ```eg_index``` correspond to different scattering configurations, each in $\\mathbb{R}^2$. Plots for 1,3,4 and 5 are included in the paper:\n", + "1. The Koch curve (see Example 2.3 and Figure 5(a))\n", + "2. Cantor set $C$ with contraction factor $\\rho=1/3$ (see Example 2.2)\n", + "3. Cantor dust (i.e. $C\\times C$) with contraction factor $\\rho=1/3$\n", + "4. A Koch snowflake boundary, constructed as the union of three Koch curves (see Fig. 8(b))\n", + "5. A solid Koch snowflake (see Example 2.4 and Fig 8(a))\n", + "6. The union of a Cantor set with contraction factor $\\rho=1/3$ and a Koch curve.\n", "\n", "Now we define our scatterer $\\Gamma$ and plot it." ] @@ -67,18 +68,20 @@ "source": [ "if eg_index == 1 # Koch curve\n", " Γ = KochCurve()\n", - "elseif eg_index == 2 # Cantor dust\n", + "elseif eg_index == 2 # Cantor set\n", " Γ = CantorSet() + [0,0]\n", - "elseif eg_index == 3 # Cantor set\n", + "elseif eg_index == 3 # Cantor dust\n", " Γ = CantorDust()\n", - "elseif eg_index == 4 # Koch snowflake, as a union of Koch curves\n", + "elseif eg_index == 4 # Koch snowflake boundary, as a union of Koch curves\n", " Koch_height = 0.288675134594813\n", " γ = KochCurve()+[-1/2,Koch_height]\n", " γ₁ = γ\n", " γ₂ = rotation2(2π/3)*γ\n", " γ₃ = rotation2(4π/3)*γ\n", " Γ = UnionInvariantMeasure([γ₁,γ₂,γ₃])\n", - "elseif eg_index == 5 # Cantor set union Koch curve\n", + "elseif eg_index == 5\n", + " Γ = KochFlake()\n", + "elseif eg_index == 6 # Cantor set union Koch curve\n", " γ₁ = CantorSet() + [0,-0.2] # embed Cantor set in ℝ² and translate\n", " γ₂ = KochCurve()\n", " Γ = UnionInvariantMeasure([γ₁, γ₂])\n", @@ -92,7 +95,7 @@ "cell_type": "markdown", "metadata": {}, "source": [ - "We note that for cases ```eg_index```$=4,5$, we are translating the fractals by adding a vector, and ```eg_index```$=4$ also has a rotation. By default, the Cantor set is defined as a subset of $\\mathbb{R}$. By adding $[0,0]$ in line 4 (similarly in line 15), the Cantor set is embedded in $\\mathbb{R}^2$." + "We note that for cases ```eg_index```$=4,6$, we are translating the fractals by adding a vector, and ```eg_index```$=4$ also has a rotation. By default, the Cantor set is defined as a subset of $\\mathbb{R}$. By adding $[0,0]$ in line 4 (similarly in line 17), the Cantor set is embedded in $\\mathbb{R}^2$." ] }, { @@ -143,15 +146,15 @@ "metadata": {}, "outputs": [], "source": [ - "h = 0.05 # max mesh width of BEM elements\n", - "hQ = 0.01; # max mesh width for quadrature in BEM" + "h = 0.05 # max element diameter of the IEM elements\n", + "hQ = 0.01; # max element diameter for IEM quadrature in " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ - "Define the boundary integral operator $A$ of (3.15), and its discretised stucture $A_h$, which contains essential information about the mesh, Galerkin matrix (4.16), etc. Then solve the system to obtain the discrete solution $\\phi_h$." + "Define the boundary integral operator $A$ of (3.15), and its discretised stucture $A_h$, which contains essential information about the mesh, Galerkin matrix (4.16), etc. Then solve the system to obtain the discrete solution $\\phi_N$." ] }, { @@ -169,7 +172,7 @@ "cell_type": "markdown", "metadata": {}, "source": [ - "The final line above represents the Galerkin system (4.15). Now we can construct an approximation to the total field using (3.14) and (3.10):" + "The final line above solves the Galerkin system (4.15). Now we can construct an approximation to the total field using (3.14) and (3.10):" ] }, { @@ -178,7 +181,7 @@ "metadata": {}, "outputs": [], "source": [ - "𝘈ϕ_N = SingleLayerPotentialHelmholtz(ϕ_N,k; h_quad = hQ)# returns function\n", + "𝘈ϕ_N = SingleLayerPotentialHelmholtz(ϕ_N,k; h_quad = hQ)# returns potential as function\n", "u_N(x) = uⁱ(x) + 𝘈ϕ_N(x) # total field" ] }, @@ -227,7 +230,7 @@ "outputs": [], "source": [ "heatmap(x₁,x₂,transpose(real(u_N.(X))), aspect_ratio = 1, \n", - " title=\"Total field\", label=false, c = :jet)\n", + " title=\"Rela part of total field\", label=false, c = :jet)\n", "plot!(Γ,color = \"black\", markersize=0.5, label=false)" ] } diff --git a/src/fractal_presets.jl b/src/fractal_presets.jl index c496456..7c9b26e 100644 --- a/src/fractal_presets.jl +++ b/src/fractal_presets.jl @@ -309,7 +309,10 @@ function KochCurve() false true true true; false false true true] - return InvariantMeasure(IFS, connectedness = touching_singularities); + G = get_group_operations2D(GroupGenerator2D(1,[π/2],[0.5,0])) + return InvariantMeasure(IFS, + connectedness = touching_singularities, + symmetry_group = G); end