{
 "cells": [
  {
   "cell_type": "code",
   "execution_count": 2,
   "metadata": {},
   "outputs": [],
   "source": [
    "using Plots, ComplexPhasePortrait, ApproxFun, SingularIntegralEquations\n",
    "gr();"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "# M3M6: Methods of Mathematical Physics\n",
    "\n",
    "$$\n",
    "\\def\\dashint{{\\int\\!\\!\\!\\!\\!\\!-\\,}}\n",
    "\\def\\infdashint{\\dashint_{\\!\\!\\!-\\infty}^{\\,\\infty}}\n",
    "\\def\\D{\\,{\\rm d}}\n",
    "\\def\\dx{\\D x}\n",
    "\\def\\dt{\\D t}\n",
    "\\def\\C{{\\mathbb C}}\n",
    "\\def\\CC{{\\cal C}}\n",
    "\\def\\HH{{\\cal H}}\n",
    "\\def\\I{{\\rm i}}\n",
    "\\def\\qqfor{\\qquad\\hbox{for}\\qquad}\n",
    "$$\n",
    "\n",
    "Dr. Sheehan Olver\n",
    "<br>\n",
    "s.olver@imperial.ac.uk\n",
    "\n",
    "<br>\n",
    "Website: https://github.com/dlfivefifty/M3M6LectureNotes\n",
    "\n",
    "# Lecture 10: Hilbert transforms\n",
    "\n",
    "1. Hilbert transform on the interval\n",
    "    - Plemelj theorem for additive jump\n",
    "2. Real / imaginary parts of the Cauchy transform\n",
    "\n",
    "The Hilbert transform is fundamental to the study of singular integrals, and is of particular use in Signal processing. We introduce it via its links to the Cauchy transform.\n",
    "\n",
    "## Hilbert transform on the interval\n",
    "\n",
    "We now investigate the additive jump that the Cauchy transform satisfies. That is, we know\n",
    "$$\n",
    "    \\CC^+ f(x) - \\CC^- f(x) = f(x)\n",
    "$$\n",
    "but what about \n",
    "$$\n",
    "    \\CC^+ f(x) + \\CC^- f(x)?\n",
    "$$\n",
    "The answer is given in terms of a principal value integral, but now the integral is on a finite domain with a singularity at a prescribed point $x$ (normally a pole):\n",
    "\n",
    "**Definition (Principal value integral on an interval)** \n",
    "The (Cauchy) principal value integral at $a < x < b$ is defined as\n",
    "$$\n",
    "\\dashint_a^b f(t) \\dt := \\lim_{\\epsilon \\rightarrow 0} \\left(\\int_a^{x - \\epsilon} + \\int_{x + \\epsilon}^b \\right) f(t) \\dt\n",
    "$$\n",
    "Typically, the singular point $x$ is obvious from the integrand.\n",
    "\n",
    "_Examples_:\n",
    "\n",
    "1. Here, the singular point is implicitly 0:\n",
    "$$\\dashint_{-1}^1  {1 \\over x} \\dx = 0$$\n",
    "3. Here, the singular point is implicitly 0:\n",
    "$$\\dashint_0^1  {1 \\over x} \\dx = \\hbox{undefined}$$\n",
    "2. Here, the singular point is also implicitly 0\n",
    "$$\\dashint_{-1}^1  \\cot x \\dx = 0$$\n",
    "2. Here, the singular point is implicitly $x$:\n",
    "$$\\dashint_{-1}^1  {1 \\over t - x} \\D t = \\log(1+x) - \\log(1-x)$$\n",
    "5. Here the singular point is also implicitly $x$:\n",
    "$$\\dashint_{-1}^1  \\cot(t-x) \\dt = 0$$\n",
    "\n",
    "Here is a plot of the integration domain, showing how it avoids the singularity:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "metadata": {},
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "log(1 - x) - log(1 + x) = -0.20067069546215122\n"
     ]
    },
    {
     "data": {
      "image/svg+xml": [
       "<?xml version=\"1.0\" encoding=\"utf-8\"?>\n",
       "<svg xmlns=\"http://www.w3.org/2000/svg\" xmlns:xlink=\"http://www.w3.org/1999/xlink\" width=\"600\" height=\"400\" viewBox=\"0 0 2400 1600\">\n",
       "<defs>\n",
       "  <clipPath id=\"clip5400\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"2000\" height=\"2000\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<defs>\n",
       "  <clipPath id=\"clip5401\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"2400\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polygon clip-path=\"url(#clip5401)\" points=\"\n",
       "0,1600 2400,1600 2400,0 0,0 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip5402\">\n",
       "    <rect x=\"480\" y=\"0\" width=\"1681\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polygon clip-path=\"url(#clip5401)\" points=\"\n",
       "176.123,1503.47 2321.26,1503.47 2321.26,125.984 176.123,125.984 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip5403\">\n",
       "    <rect x=\"176\" y=\"125\" width=\"2146\" height=\"1378\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  236.827,1503.47 236.827,125.984 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  742.759,1503.47 742.759,125.984 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  1248.69,1503.47 1248.69,125.984 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  1754.62,1503.47 1754.62,125.984 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  2260.55,1503.47 2260.55,125.984 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  176.123,1503.47 2321.26,1503.47 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  176.123,1159.1 2321.26,1159.1 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  176.123,814.729 2321.26,814.729 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  176.123,470.357 2321.26,470.357 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  176.123,125.984 2321.26,125.984 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  176.123,1503.47 2321.26,1503.47 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  176.123,1503.47 176.123,125.984 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  236.827,1503.47 236.827,1482.81 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  742.759,1503.47 742.759,1482.81 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1248.69,1503.47 1248.69,1482.81 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1754.62,1503.47 1754.62,1482.81 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2260.55,1503.47 2260.55,1482.81 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  176.123,1503.47 208.3,1503.47 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  176.123,1159.1 208.3,1159.1 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  176.123,814.729 208.3,814.729 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  176.123,470.357 208.3,470.357 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  176.123,125.984 208.3,125.984 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 236.827, 1557.47)\" x=\"236.827\" y=\"1557.47\">-1.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 742.759, 1557.47)\" x=\"742.759\" y=\"1557.47\">-0.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 1248.69, 1557.47)\" x=\"1248.69\" y=\"1557.47\">0.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 1754.62, 1557.47)\" x=\"1754.62\" y=\"1557.47\">0.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 2260.55, 1557.47)\" x=\"2260.55\" y=\"1557.47\">1.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 152.123, 1520.97)\" x=\"152.123\" y=\"1520.97\">-40</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 152.123, 1176.6)\" x=\"152.123\" y=\"1176.6\">-20</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 152.123, 832.229)\" x=\"152.123\" y=\"832.229\">0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 152.123, 487.857)\" x=\"152.123\" y=\"487.857\">20</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 152.123, 143.484)\" x=\"152.123\" y=\"143.484\">40</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:84px; text-anchor:middle;\" transform=\"rotate(0, 1248.69, 73.2)\" x=\"1248.69\" y=\"73.2\">eps = 0.05, integral = -0.20067069546215155</text>\n",
       "</g>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#009af9; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1299.28,1159.05 1299.22,1158.64 1299.09,1157.81 1298.91,1156.57 1298.67,1154.94 1298.36,1152.93 1297.99,1150.54 1297.57,1147.8 1297.08,1144.72 1296.53,1141.33 \n",
       "  1295.92,1137.64 1295.25,1133.68 1294.52,1129.46 1293.73,1125.03 1292.88,1120.39 1291.96,1115.58 1290.99,1110.61 1289.96,1105.51 1288.87,1100.3 1287.71,1095.01 \n",
       "  1286.5,1089.65 1285.23,1084.24 1283.9,1078.8 1282.51,1073.35 1281.06,1067.91 1279.55,1062.47 1277.98,1057.07 1276.36,1051.71 1274.67,1046.41 1272.93,1041.16 \n",
       "  1271.13,1035.98 1269.27,1030.88 1267.36,1025.86 1265.38,1020.93 1263.35,1016.09 1261.26,1011.34 1259.12,1006.7 1256.92,1002.15 1254.66,997.712 1252.35,993.372 \n",
       "  1249.98,989.136 1247.55,985.004 1245.08,980.976 1242.54,977.051 1239.95,973.229 1237.31,969.508 1234.61,965.886 1231.86,962.364 1229.06,958.937 1226.2,955.606 \n",
       "  1223.29,952.368 1220.33,949.222 1217.32,946.164 1214.25,943.193 1211.14,940.307 1207.97,937.504 1204.75,934.782 1201.48,932.138 1198.16,929.57 1194.8,927.076 \n",
       "  1191.38,924.654 1187.91,922.302 1184.4,920.018 1180.84,917.8 1177.23,915.645 1173.57,913.552 1169.87,911.52 1166.12,909.545 1162.33,907.627 1158.49,905.763 \n",
       "  1154.6,903.952 1150.67,902.192 1146.7,900.481 1142.68,898.819 1138.62,897.203 1134.52,895.632 1130.38,894.104 1126.19,892.619 1121.96,891.174 1117.69,889.769 \n",
       "  1113.39,888.402 1109.04,887.072 1104.65,885.778 1100.22,884.518 1095.76,883.292 1091.26,882.099 1086.72,880.937 1082.14,879.805 1077.53,878.703 1072.89,877.63 \n",
       "  1068.2,876.584 1063.49,875.566 1058.74,874.573 1053.95,873.605 1049.13,872.662 1044.28,871.743 1039.4,870.847 1034.49,869.973 1029.55,869.12 1024.58,868.289 \n",
       "  1019.57,867.478 1014.54,866.686 1009.48,865.914 1004.39,865.16 999.28,864.424 994.138,863.706 988.971,863.005 983.779,862.32 978.561,861.652 973.319,860.998 \n",
       "  968.054,860.36 962.766,859.737 957.455,859.128 952.123,858.533 946.769,857.951 941.395,857.382 936.001,856.826 930.588,856.283 925.156,855.751 919.706,855.232 \n",
       "  914.238,854.723 908.754,854.226 903.254,853.74 897.738,853.264 892.207,852.798 886.661,852.342 881.103,851.896 875.531,851.46 869.947,851.032 864.351,850.614 \n",
       "  858.744,850.204 853.126,849.803 847.499,849.41 841.863,849.026 836.218,848.649 830.565,848.279 824.905,847.918 819.239,847.563 813.567,847.216 807.889,846.876 \n",
       "  802.207,846.542 796.521,846.215 790.832,845.895 785.14,845.581 779.446,845.273 773.751,844.971 768.055,844.675 762.36,844.385 756.664,844.1 750.971,843.821 \n",
       "  745.279,843.547 739.589,843.278 733.903,843.015 728.221,842.756 722.544,842.502 716.872,842.254 711.205,842.009 705.545,841.77 699.893,841.535 694.248,841.304 \n",
       "  688.612,841.077 682.984,840.855 677.367,840.637 671.76,840.423 666.164,840.212 660.58,840.006 655.008,839.803 649.449,839.604 643.904,839.409 638.373,839.217 \n",
       "  632.857,839.029 627.357,838.844 621.872,838.662 616.405,838.484 610.955,838.308 605.523,838.136 600.109,837.967 594.715,837.801 589.341,837.638 583.988,837.478 \n",
       "  578.655,837.321 573.345,837.166 568.056,837.015 562.791,836.865 557.55,836.719 552.332,836.575 547.14,836.434 541.972,836.295 536.831,836.159 531.716,836.025 \n",
       "  526.629,835.893 521.569,835.764 516.538,835.637 511.535,835.512 506.562,835.39 501.619,835.269 496.707,835.151 491.826,835.035 486.977,834.921 482.16,834.809 \n",
       "  477.375,834.698 472.625,834.59 467.908,834.484 463.226,834.38 458.578,834.277 453.967,834.177 449.391,834.078 444.852,833.981 440.35,833.886 435.886,833.792 \n",
       "  431.46,833.7 427.073,833.61 422.725,833.521 418.417,833.434 414.149,833.349 409.921,833.265 405.735,833.183 401.591,833.103 397.488,833.023 393.428,832.946 \n",
       "  389.411,832.87 385.438,832.795 381.509,832.722 377.624,832.65 373.784,832.579 369.99,832.51 366.241,832.442 362.538,832.376 358.882,832.311 355.273,832.247 \n",
       "  351.711,832.184 348.198,832.123 344.732,832.063 341.316,832.004 337.948,831.947 334.63,831.891 331.361,831.836 328.143,831.782 324.975,831.729 321.859,831.678 \n",
       "  318.793,831.627 315.779,831.578 312.818,831.53 309.908,831.483 307.051,831.437 304.248,831.392 301.497,831.348 298.801,831.306 296.158,831.264 293.569,831.224 \n",
       "  291.035,831.184 288.556,831.146 286.132,831.108 283.763,831.072 281.45,831.037 279.193,831.002 276.992,830.969 274.848,830.936 272.76,830.905 270.729,830.875 \n",
       "  268.756,830.845 266.839,830.817 264.981,830.789 263.18,830.762 261.437,830.737 259.753,830.712 258.127,830.688 256.56,830.665 255.051,830.643 253.601,830.622 \n",
       "  252.211,830.602 250.88,830.583 249.608,830.565 248.396,830.547 247.244,830.531 246.151,830.515 245.119,830.5 244.146,830.486 243.234,830.473 242.382,830.461 \n",
       "  241.591,830.45 240.86,830.44 240.19,830.43 239.58,830.422 239.031,830.414 238.543,830.407 238.116,830.401 237.75,830.396 237.445,830.391 237.201,830.388 \n",
       "  237.018,830.385 236.895,830.384 236.834,830.383 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5403)\" style=\"stroke:#009af9; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2260.55,795.598 2260.5,795.596 2260.4,795.594 2260.25,795.591 2260.04,795.587 2259.79,795.582 2259.49,795.575 2259.14,795.568 2258.73,795.559 2258.28,795.55 \n",
       "  2257.78,795.539 2257.22,795.527 2256.62,795.515 2255.96,795.501 2255.26,795.486 2254.51,795.47 2253.7,795.453 2252.85,795.434 2251.95,795.415 2251,795.395 \n",
       "  2249.99,795.373 2248.94,795.351 2247.84,795.327 2246.69,795.302 2245.5,795.276 2244.25,795.249 2242.96,795.221 2241.61,795.191 2240.22,795.161 2238.78,795.129 \n",
       "  2237.29,795.096 2235.76,795.062 2234.18,795.027 2232.55,794.991 2230.87,794.953 2229.14,794.914 2227.37,794.874 2225.55,794.833 2223.69,794.791 2221.78,794.747 \n",
       "  2219.82,794.702 2217.82,794.656 2215.78,794.608 2213.68,794.56 2211.55,794.51 2209.36,794.458 2207.14,794.406 2204.87,794.352 2202.55,794.296 2200.19,794.24 \n",
       "  2197.79,794.182 2195.35,794.122 2192.86,794.061 2190.33,793.999 2187.76,793.935 2185.14,793.87 2182.49,793.804 2179.79,793.736 2177.05,793.666 2174.27,793.595 \n",
       "  2171.45,793.523 2168.59,793.449 2165.69,793.373 2162.75,793.296 2159.78,793.217 2156.76,793.137 2153.71,793.055 2150.62,792.971 2147.49,792.886 2144.32,792.798 \n",
       "  2141.11,792.71 2137.87,792.619 2134.6,792.527 2131.29,792.433 2127.94,792.337 2124.56,792.239 2121.14,792.139 2117.69,792.038 2114.2,791.934 2110.69,791.829 \n",
       "  2107.14,791.722 2103.55,791.612 2099.94,791.501 2096.29,791.387 2092.61,791.272 2088.9,791.154 2085.17,791.034 2081.4,790.912 2077.6,790.788 2073.77,790.661 \n",
       "  2069.92,790.532 2066.03,790.401 2062.12,790.267 2058.18,790.131 2054.22,789.993 2050.22,789.852 2046.21,789.708 2042.16,789.562 2038.1,789.414 2034,789.262 \n",
       "  2029.89,789.108 2025.75,788.951 2021.59,788.791 2017.4,788.629 2013.19,788.463 2008.97,788.295 2004.72,788.123 2000.45,787.949 1996.16,787.771 1991.85,787.59 \n",
       "  1987.52,787.406 1983.18,787.218 1978.81,787.027 1974.43,786.833 1970.03,786.635 1965.62,786.434 1961.19,786.228 1956.74,786.02 1952.28,785.807 1947.8,785.591 \n",
       "  1943.31,785.37 1938.81,785.146 1934.3,784.917 1929.77,784.684 1925.23,784.447 1920.68,784.206 1916.12,783.96 1911.55,783.71 1906.97,783.455 1902.38,783.195 \n",
       "  1897.79,782.931 1893.18,782.661 1888.57,782.386 1883.95,782.107 1879.32,781.822 1874.69,781.531 1870.06,781.235 1865.41,780.934 1860.77,780.627 1856.12,780.313 \n",
       "  1851.47,779.994 1846.81,779.669 1842.16,779.337 1837.5,778.999 1832.84,778.655 1828.18,778.303 1823.52,777.945 1818.87,777.58 1814.21,777.207 1809.56,776.827 \n",
       "  1804.9,776.44 1800.26,776.044 1795.61,775.641 1790.97,775.23 1786.33,774.81 1781.7,774.382 1777.08,773.945 1772.46,773.5 1767.84,773.045 1763.24,772.58 \n",
       "  1758.64,772.106 1754.05,771.622 1749.47,771.128 1744.9,770.624 1740.34,770.109 1735.79,769.583 1731.25,769.045 1726.73,768.496 1722.21,767.936 1717.71,767.363 \n",
       "  1713.22,766.778 1708.75,766.18 1704.29,765.569 1699.84,764.944 1695.41,764.306 1690.99,763.653 1686.6,762.986 1682.21,762.304 1677.85,761.606 1673.5,760.893 \n",
       "  1669.18,760.163 1664.87,759.417 1660.58,758.653 1656.31,757.872 1652.06,757.072 1647.83,756.254 1643.62,755.417 1639.44,754.56 1635.28,753.682 1631.14,752.784 \n",
       "  1627.02,751.864 1622.93,750.922 1618.86,749.957 1614.82,748.968 1610.8,747.956 1606.81,746.918 1602.84,745.855 1598.91,744.766 1594.99,743.649 1591.11,742.505 \n",
       "  1587.25,741.332 1583.43,740.129 1579.63,738.896 1575.86,737.631 1572.12,736.334 1568.41,735.003 1564.73,733.639 1561.09,732.238 1557.47,730.802 1553.89,729.328 \n",
       "  1550.34,727.815 1546.82,726.263 1543.34,724.669 1539.88,723.034 1536.47,721.355 1533.09,719.631 1529.74,717.861 1526.43,716.044 1523.15,714.178 1519.91,712.262 \n",
       "  1516.71,710.294 1513.54,708.273 1510.41,706.198 1507.32,704.066 1504.26,701.877 1501.25,699.628 1498.27,697.319 1495.33,694.947 1492.43,692.511 1489.57,690.01 \n",
       "  1486.75,687.441 1483.98,684.803 1481.24,682.095 1478.54,679.314 1475.88,676.46 1473.27,673.53 1470.7,670.524 1468.17,667.439 1465.68,664.275 1463.23,661.03 \n",
       "  1460.83,657.703 1458.47,654.293 1456.16,650.799 1453.89,647.22 1451.66,643.555 1449.48,639.805 1447.34,635.969 1445.25,632.047 1443.2,628.039 1441.2,623.947 \n",
       "  1439.24,619.771 1437.33,615.513 1435.47,611.174 1433.65,606.758 1431.88,602.266 1430.16,597.702 1428.48,593.071 1426.85,588.377 1425.27,583.624 1423.73,578.82 \n",
       "  1422.24,573.971 1420.8,569.084 1419.41,564.169 1418.07,559.233 1416.77,554.288 1415.53,549.343 1414.33,544.412 1413.18,539.507 1412.08,534.641 1411.03,529.828 \n",
       "  1410.03,525.085 1409.08,520.427 1408.18,515.87 1407.32,511.431 1406.52,507.129 1405.76,502.981 1405.06,499.005 1404.41,495.219 1403.8,491.641 1403.25,488.288 \n",
       "  1402.75,485.178 1402.29,482.327 1401.89,479.75 1401.54,477.461 1401.23,475.474 1400.98,473.8 1400.78,472.448 1400.63,471.427 1400.53,470.743 1400.48,470.4 \n",
       "  \n",
       "  \"/>\n",
       "<circle clip-path=\"url(#clip5403)\" style=\"fill:#000000; stroke:none; fill-opacity:1\" cx=\"1349.88\" cy=\"814.729\" r=\"18\"/>\n",
       "<circle clip-path=\"url(#clip5403)\" style=\"fill:#e26f46; stroke:none; fill-opacity:1\" cx=\"1349.88\" cy=\"814.729\" r=\"14\"/>\n",
       "<polygon clip-path=\"url(#clip5401)\" points=\"\n",
       "1985.19,390.944 2249.26,390.944 2249.26,209.504 1985.19,209.504 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1985.19,390.944 2249.26,390.944 2249.26,209.504 1985.19,209.504 1985.19,390.944 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5401)\" style=\"stroke:#009af9; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2009.19,269.984 2153.19,269.984 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:start;\" transform=\"rotate(0, 2177.19, 287.484)\" x=\"2177.19\" y=\"287.484\">f</text>\n",
       "</g>\n",
       "<circle clip-path=\"url(#clip5401)\" style=\"fill:#000000; stroke:none; fill-opacity:1\" cx=\"2093.19\" cy=\"330.464\" r=\"25\"/>\n",
       "<circle clip-path=\"url(#clip5401)\" style=\"fill:#e26f46; stroke:none; fill-opacity:1\" cx=\"2093.19\" cy=\"330.464\" r=\"21\"/>\n",
       "<g clip-path=\"url(#clip5401)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:start;\" transform=\"rotate(0, 2177.19, 347.964)\" x=\"2177.19\" y=\"347.964\">x</text>\n",
       "</g>\n",
       "</svg>\n"
      ]
     },
     "execution_count": 5,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "ε = 0.05\n",
    "x = 0.1\n",
    "Γ = Segment(-1 , x-ε) ∪ Segment(x+ε , 1)\n",
    "t = Fun(Γ)\n",
    "f = Fun(1/(t - x))\n",
    "\n",
    "@show log(1 - x) - log(1+x)\n",
    "\n",
    "plot(f; label=\"f\", title=\"eps = $(ε), integral = $(sum(f))\",ylims=(-40,40))\n",
    "scatter!([x],[0.0];label=\"x\")"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Using the principal value integral, we can define a singular integral operator called the _Hilbert transform_:\n",
    "\n",
    "**Definition (Hilbert transform)** For $a < x < b$, define\n",
    "$$\n",
    "\\HH_{[a,b]} f(x) =  {1 \\over \\pi} \\dashint_a^b {f(t) \\over {t - x}} \\D t\n",
    "$$\n",
    "\n",
    "Here's a plot of some simple hilbert transforms:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "metadata": {},
   "outputs": [
    {
     "data": {
      "image/svg+xml": [
       "<?xml version=\"1.0\" encoding=\"utf-8\"?>\n",
       "<svg xmlns=\"http://www.w3.org/2000/svg\" xmlns:xlink=\"http://www.w3.org/1999/xlink\" width=\"600\" height=\"400\" viewBox=\"0 0 2400 1600\">\n",
       "<defs>\n",
       "  <clipPath id=\"clip5800\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"2000\" height=\"2000\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<defs>\n",
       "  <clipPath id=\"clip5801\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"2400\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polygon clip-path=\"url(#clip5801)\" points=\"\n",
       "0,1600 2400,1600 2400,0 0,0 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip5802\">\n",
       "    <rect x=\"480\" y=\"0\" width=\"1681\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polygon clip-path=\"url(#clip5801)\" points=\"\n",
       "189.504,1503.47 2321.26,1503.47 2321.26,47.2441 189.504,47.2441 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip5803\">\n",
       "    <rect x=\"189\" y=\"47\" width=\"2133\" height=\"1457\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  249.837,1503.47 249.837,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  752.609,1503.47 752.609,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  1255.38,1503.47 1255.38,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  1758.15,1503.47 1758.15,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  2260.93,1503.47 2260.93,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,1324.13 2321.26,1324.13 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,1115.28 2321.26,1115.28 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,906.423 2321.26,906.423 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,697.57 2321.26,697.57 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,488.717 2321.26,488.717 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,279.864 2321.26,279.864 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,71.0116 2321.26,71.0116 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1503.47 2321.26,1503.47 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1503.47 189.504,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  249.837,1503.47 249.837,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  752.609,1503.47 752.609,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1255.38,1503.47 1255.38,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1758.15,1503.47 1758.15,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2260.93,1503.47 2260.93,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1324.13 221.48,1324.13 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1115.28 221.48,1115.28 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,906.423 221.48,906.423 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,697.57 221.48,697.57 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,488.717 221.48,488.717 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,279.864 221.48,279.864 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,71.0116 221.48,71.0116 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 249.837, 1557.47)\" x=\"249.837\" y=\"1557.47\">-1.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 752.609, 1557.47)\" x=\"752.609\" y=\"1557.47\">-0.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 1255.38, 1557.47)\" x=\"1255.38\" y=\"1557.47\">0.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 1758.15, 1557.47)\" x=\"1758.15\" y=\"1557.47\">0.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 2260.93, 1557.47)\" x=\"2260.93\" y=\"1557.47\">1.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 1341.63)\" x=\"165.504\" y=\"1341.63\">-1.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 1132.78)\" x=\"165.504\" y=\"1132.78\">-1.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 923.923)\" x=\"165.504\" y=\"923.923\">-0.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 715.07)\" x=\"165.504\" y=\"715.07\">0.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 506.217)\" x=\"165.504\" y=\"506.217\">0.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 297.364)\" x=\"165.504\" y=\"297.364\">1.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 88.5116)\" x=\"165.504\" y=\"88.5116\">1.5</text>\n",
       "</g>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#009af9; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2260.93,697.57 2260.81,680.423 2259.9,646.192 2258.06,612.148 2255.31,578.414 2251.65,545.112 2247.08,512.36 2241.61,480.272 2235.23,448.957 2227.96,418.521 \n",
       "  2219.81,389.061 2210.77,360.67 2200.86,333.432 2190.09,307.424 2178.47,282.716 2166,259.368 2152.7,237.434 2138.59,216.958 2123.67,197.976 2107.95,180.513 \n",
       "  2091.46,164.589 2074.21,150.213 2056.21,137.389 2037.47,126.109 2018.03,116.361 1997.89,108.124 1977.07,101.371 1955.59,96.0693 1933.47,92.18 1910.74,89.6591 \n",
       "  1887.4,88.4582 1863.49,88.5245 1839.03,89.8022 1814.03,92.2322 1788.53,95.7532 1762.53,100.302 1736.08,105.814 1709.18,112.224 1681.87,119.465 1654.17,127.473 \n",
       "  1626.11,136.182 1597.71,145.527 1569,155.445 1540,165.874 1510.74,176.754 1481.25,188.027 1451.55,199.637 1421.68,211.53 1391.65,223.656 1361.5,235.966 \n",
       "  1331.25,248.414 1300.93,260.958 1270.57,273.556 1240.19,286.173 1209.83,298.772 1179.52,311.324 1149.27,323.797 1119.11,336.167 1089.09,348.409 1059.21,360.502 \n",
       "  1029.51,372.427 1000.02,384.166 970.763,395.707 941.764,407.036 913.052,418.143 884.651,429.019 856.589,439.658 828.891,450.053 801.583,460.202 774.688,470.102 \n",
       "  748.232,479.752 722.238,489.152 696.731,498.303 671.734,507.207 647.269,515.867 623.359,524.287 600.026,532.471 577.291,540.425 555.174,548.154 533.696,555.664 \n",
       "  512.877,562.962 492.736,570.055 473.29,576.951 454.558,583.656 436.556,590.178 419.302,596.527 402.81,602.709 387.097,608.733 372.175,614.607 358.06,620.341 \n",
       "  344.763,625.942 332.297,631.419 320.674,636.781 309.903,642.035 299.995,647.192 290.958,652.259 282.802,657.246 275.533,662.159 269.158,667.01 263.683,671.804 \n",
       "  259.113,676.553 255.452,681.263 252.703,685.943 250.869,690.602 249.952,695.249 249.837,697.57 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5803)\" style=\"stroke:#e26f46; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2260.82,1462.26 2259.95,1460.49 2258.22,1456.95 2255.62,1451.67 2252.16,1444.66 2247.84,1435.97 2242.67,1425.64 2236.65,1413.71 2229.78,1400.24 2222.06,1385.31 \n",
       "  2213.52,1368.97 2204.15,1351.3 2193.96,1332.39 2182.97,1312.32 2171.17,1291.18 2158.58,1269.06 2145.22,1246.06 2131.09,1222.27 2116.2,1197.8 2100.57,1172.73 \n",
       "  2084.22,1147.17 2067.15,1121.22 2049.38,1094.97 2030.92,1068.51 2011.8,1041.94 1992.02,1015.35 1971.61,988.811 1950.58,962.415 1928.95,936.235 1906.75,910.345 \n",
       "  1883.98,884.813 1860.67,859.7 1836.83,835.067 1812.5,810.964 1787.68,787.441 1762.41,764.54 1736.7,742.298 1710.57,720.749 1684.06,699.919 1657.17,679.833 \n",
       "  1629.94,660.507 1602.38,641.956 1574.52,624.19 1546.39,607.215 1518.01,591.031 1489.4,575.638 1460.59,561.03 1431.61,547.201 1402.47,534.138 1373.2,521.829 \n",
       "  1343.84,510.258 1314.4,499.408 1284.9,489.261 1255.38,479.794 1225.86,470.987 1196.37,462.816 1166.93,455.257 1137.56,448.287 1108.29,441.879 1079.16,436.009 \n",
       "  1050.17,430.651 1021.36,425.78 992.752,421.37 964.371,417.396 936.24,413.833 908.385,410.656 880.828,407.842 853.595,405.367 826.707,403.209 800.19,401.345 \n",
       "  774.065,399.754 748.354,398.415 723.081,397.308 698.266,396.415 673.932,395.716 650.099,395.195 626.787,394.835 604.018,394.619 581.81,394.533 560.182,394.562 \n",
       "  539.154,394.694 518.744,394.914 498.968,395.211 479.844,395.574 461.389,395.992 443.618,396.455 426.547,396.954 410.19,397.479 394.562,398.024 379.676,398.579 \n",
       "  365.544,399.139 352.18,399.697 339.594,400.247 327.798,400.783 316.801,401.301 306.613,401.795 297.244,402.262 288.7,402.698 280.989,403.099 274.118,403.464 \n",
       "  268.093,403.788 262.919,404.07 258.601,404.308 255.141,404.5 252.544,404.645 250.812,404.742 249.945,404.791 \n",
       "  \"/>\n",
       "<polygon clip-path=\"url(#clip5801)\" points=\"\n",
       "1661.54,312.204 2249.26,312.204 2249.26,130.764 1661.54,130.764 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1661.54,312.204 2249.26,312.204 2249.26,130.764 1661.54,130.764 1661.54,312.204 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#009af9; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1685.54,191.244 1829.54,191.244 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:start;\" transform=\"rotate(0, 1853.54, 208.744)\" x=\"1853.54\" y=\"208.744\">function</text>\n",
       "</g>\n",
       "<polyline clip-path=\"url(#clip5801)\" style=\"stroke:#e26f46; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1685.54,251.724 1829.54,251.724 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip5801)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:start;\" transform=\"rotate(0, 1853.54, 269.224)\" x=\"1853.54\" y=\"269.224\">hilbert transform</text>\n",
       "</g>\n",
       "</svg>\n"
      ]
     },
     "execution_count": 7,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "x = Fun(-1 .. 1)\n",
    "f = exp(x)*sqrt(1-x^2)\n",
    "plot(f; label=\"function\")\n",
    "plot!(hilbert(f); label=\"hilbert transform\")"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "metadata": {},
   "outputs": [
    {
     "data": {
      "image/svg+xml": [
       "<?xml version=\"1.0\" encoding=\"utf-8\"?>\n",
       "<svg xmlns=\"http://www.w3.org/2000/svg\" xmlns:xlink=\"http://www.w3.org/1999/xlink\" width=\"600\" height=\"400\" viewBox=\"0 0 2400 1600\">\n",
       "<defs>\n",
       "  <clipPath id=\"clip6000\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"2000\" height=\"2000\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<defs>\n",
       "  <clipPath id=\"clip6001\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"2400\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polygon clip-path=\"url(#clip6001)\" points=\"\n",
       "0,1600 2400,1600 2400,0 0,0 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip6002\">\n",
       "    <rect x=\"480\" y=\"0\" width=\"1681\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polygon clip-path=\"url(#clip6001)\" points=\"\n",
       "189.504,1503.47 2321.26,1503.47 2321.26,47.2441 189.504,47.2441 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip6003\">\n",
       "    <rect x=\"189\" y=\"47\" width=\"2133\" height=\"1457\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  249.728,1503.47 249.728,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  752.555,1503.47 752.555,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  1255.38,1503.47 1255.38,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  1758.21,1503.47 1758.21,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  2261.04,1503.47 2261.04,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,1503.47 2321.26,1503.47 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,1139.42 2321.26,1139.42 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,775.359 2321.26,775.359 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,411.302 2321.26,411.302 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,47.2441 2321.26,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1503.47 2321.26,1503.47 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1503.47 189.504,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  249.728,1503.47 249.728,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  752.555,1503.47 752.555,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1255.38,1503.47 1255.38,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1758.21,1503.47 1758.21,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2261.04,1503.47 2261.04,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1503.47 221.48,1503.47 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1139.42 221.48,1139.42 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,775.359 221.48,775.359 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,411.302 221.48,411.302 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,47.2441 221.48,47.2441 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 249.728, 1557.47)\" x=\"249.728\" y=\"1557.47\">-1.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 752.555, 1557.47)\" x=\"752.555\" y=\"1557.47\">-0.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 1255.38, 1557.47)\" x=\"1255.38\" y=\"1557.47\">0.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 1758.21, 1557.47)\" x=\"1758.21\" y=\"1557.47\">0.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 2261.04, 1557.47)\" x=\"2261.04\" y=\"1557.47\">1.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 1520.97)\" x=\"165.504\" y=\"1520.97\">-5.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 1156.92)\" x=\"165.504\" y=\"1156.92\">-2.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 792.859)\" x=\"165.504\" y=\"792.859\">0.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 428.802)\" x=\"165.504\" y=\"428.802\">2.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 64.7441)\" x=\"165.504\" y=\"64.7441\">5.0</text>\n",
       "</g>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#009af9; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2260.92,-25430.9 2260,-7954.75 2258.17,-4456.34 2255.42,-2954.77 2251.76,-2118.83 2247.19,-1585.47 2241.71,-1215.06 2235.34,-942.432 2228.07,-733.102 2219.91,-567.101 \n",
       "  2210.87,-432.071 2200.96,-319.953 2190.19,-225.27 2178.57,-144.168 2166.1,-73.8567 2152.8,-12.268 2138.68,42.1646 2123.76,90.649 2108.05,134.131 2091.55,173.363 \n",
       "  2074.3,208.949 2056.29,241.382 2037.56,271.065 2018.11,298.335 1997.97,323.471 1977.15,346.712 1955.67,368.257 1933.55,388.28 1910.81,406.928 1887.47,424.329 \n",
       "  1863.56,440.596 1839.09,455.825 1814.09,470.103 1788.58,483.506 1762.59,496.1 1736.13,507.947 1709.23,519.1 1681.92,529.608 1654.22,539.513 1626.15,548.856 \n",
       "  1597.75,557.673 1569.03,565.995 1540.03,573.852 1510.77,581.272 1481.28,588.28 1451.58,594.899 1421.7,601.149 1391.66,607.05 1361.51,612.62 1331.26,617.876 \n",
       "  1300.93,622.832 1270.57,627.503 1240.19,631.903 1209.83,636.042 1179.51,639.933 1149.26,643.585 1119.1,647.008 1089.07,650.21 1059.19,653.199 1029.49,655.983 \n",
       "  999.994,658.567 970.732,660.957 941.73,663.157 913.015,665.173 884.611,667.007 856.546,668.663 828.845,670.14 801.534,671.442 774.636,672.568 748.177,673.516 \n",
       "  722.181,674.286 696.671,674.874 671.671,675.276 647.203,675.486 623.291,675.497 599.955,675.299 577.218,674.883 555.099,674.233 533.619,673.334 512.797,672.165 \n",
       "  492.653,670.703 473.206,668.917 454.471,666.774 436.468,664.23 419.212,661.232 402.718,657.718 387.003,653.606 372.08,648.797 357.963,643.164 344.665,636.548 \n",
       "  332.198,628.735 320.573,619.45 309.801,608.314 299.892,594.809 290.854,578.193 282.697,557.378 275.427,530.691 269.052,495.432 263.576,446.934 259.005,376.379 \n",
       "  255.344,264.871 252.595,63.2777 250.761,-408.56 249.843,-2772.08 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6003)\" style=\"stroke:#e26f46; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2260.93,508.707 2260.06,508.837 2258.33,509.095 2255.73,509.482 2252.27,509.997 2247.95,510.638 2242.78,511.403 2236.75,512.29 2229.88,513.298 2222.17,514.423 \n",
       "  2213.62,515.662 2204.25,517.012 2194.06,518.471 2183.07,520.033 2171.27,521.696 2158.68,523.455 2145.32,525.305 2131.18,527.244 2116.3,529.265 2100.67,531.364 \n",
       "  2084.31,533.537 2067.23,535.779 2049.46,538.084 2031,540.449 2011.88,542.867 1992.1,545.334 1971.69,547.845 1950.66,550.395 1929.03,552.98 1906.82,555.594 \n",
       "  1884.04,558.234 1860.73,560.893 1836.89,563.569 1812.56,566.256 1787.74,568.951 1762.46,571.649 1736.75,574.347 1710.62,577.04 1684.1,579.726 1657.21,582.4 \n",
       "  1629.98,585.059 1602.42,587.702 1574.56,590.323 1546.42,592.921 1518.04,595.494 1489.43,598.038 1460.62,600.551 1431.63,603.032 1402.49,605.479 1373.22,607.889 \n",
       "  1343.85,610.261 1314.4,612.594 1284.9,614.886 1255.38,617.137 1225.86,619.344 1196.36,621.508 1166.92,623.627 1137.55,625.701 1108.28,627.729 1079.14,629.71 \n",
       "  1050.15,631.645 1021.33,633.533 992.724,635.373 964.339,637.166 936.206,638.912 908.347,640.61 880.788,642.261 853.551,643.864 826.661,645.42 800.141,646.93 \n",
       "  774.013,648.393 748.299,649.809 723.023,651.18 698.206,652.505 673.869,653.785 650.034,655.019 626.72,656.21 603.948,657.356 581.737,658.459 560.108,659.518 \n",
       "  539.077,660.535 518.664,661.51 498.886,662.443 479.76,663.335 461.303,664.185 443.53,664.995 426.457,665.765 410.099,666.496 394.469,667.187 379.581,667.84 \n",
       "  365.448,668.454 352.083,669.03 339.496,669.568 327.698,670.069 316.7,670.532 306.511,670.959 297.14,671.349 288.595,671.703 280.884,672.021 274.012,672.303 \n",
       "  267.987,672.55 262.812,672.761 258.493,672.936 255.034,673.077 252.436,673.182 250.704,673.252 249.837,673.287 \n",
       "  \"/>\n",
       "<polygon clip-path=\"url(#clip6001)\" points=\"\n",
       "1661.54,312.204 2249.26,312.204 2249.26,130.764 1661.54,130.764 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1661.54,312.204 2249.26,312.204 2249.26,130.764 1661.54,130.764 1661.54,312.204 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#009af9; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1685.54,191.244 1829.54,191.244 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:start;\" transform=\"rotate(0, 1853.54, 208.744)\" x=\"1853.54\" y=\"208.744\">function</text>\n",
       "</g>\n",
       "<polyline clip-path=\"url(#clip6001)\" style=\"stroke:#e26f46; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1685.54,251.724 1829.54,251.724 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip6001)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:start;\" transform=\"rotate(0, 1853.54, 269.224)\" x=\"1853.54\" y=\"269.224\">hilbert transform</text>\n",
       "</g>\n",
       "</svg>\n"
      ]
     },
     "execution_count": 8,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "x = Fun(-1 .. 1)\n",
    "f = exp(x)/sqrt(1-x^2)\n",
    "plot(f; label=\"function\", ylims=(-5,5))\n",
    "plot!(hilbert(f); label=\"hilbert transform\")"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "metadata": {},
   "outputs": [
    {
     "data": {
      "image/svg+xml": [
       "<?xml version=\"1.0\" encoding=\"utf-8\"?>\n",
       "<svg xmlns=\"http://www.w3.org/2000/svg\" xmlns:xlink=\"http://www.w3.org/1999/xlink\" width=\"600\" height=\"400\" viewBox=\"0 0 2400 1600\">\n",
       "<defs>\n",
       "  <clipPath id=\"clip6200\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"2000\" height=\"2000\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<defs>\n",
       "  <clipPath id=\"clip6201\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"2400\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polygon clip-path=\"url(#clip6201)\" points=\"\n",
       "0,1600 2400,1600 2400,0 0,0 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip6202\">\n",
       "    <rect x=\"480\" y=\"0\" width=\"1681\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polygon clip-path=\"url(#clip6201)\" points=\"\n",
       "189.504,1503.47 2321.26,1503.47 2321.26,47.2441 189.504,47.2441 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip6203\">\n",
       "    <rect x=\"189\" y=\"47\" width=\"2133\" height=\"1457\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  249.736,1503.47 249.736,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  752.559,1503.47 752.559,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  1255.38,1503.47 1255.38,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  1758.2,1503.47 1758.2,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  2261.03,1503.47 2261.03,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,1503.47 2321.26,1503.47 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,1139.42 2321.26,1139.42 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,775.359 2321.26,775.359 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,411.302 2321.26,411.302 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#000000; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  189.504,47.2441 2321.26,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1503.47 2321.26,1503.47 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1503.47 189.504,47.2441 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  249.736,1503.47 249.736,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  752.559,1503.47 752.559,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1255.38,1503.47 1255.38,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1758.2,1503.47 1758.2,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2261.03,1503.47 2261.03,1481.63 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1503.47 221.48,1503.47 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,1139.42 221.48,1139.42 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,775.359 221.48,775.359 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,411.302 221.48,411.302 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  189.504,47.2441 221.48,47.2441 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 249.736, 1557.47)\" x=\"249.736\" y=\"1557.47\">-1.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 752.559, 1557.47)\" x=\"752.559\" y=\"1557.47\">-0.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 1255.38, 1557.47)\" x=\"1255.38\" y=\"1557.47\">0.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 1758.2, 1557.47)\" x=\"1758.2\" y=\"1557.47\">0.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:middle;\" transform=\"rotate(0, 2261.03, 1557.47)\" x=\"2261.03\" y=\"1557.47\">1.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 1520.97)\" x=\"165.504\" y=\"1520.97\">-5.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 1156.92)\" x=\"165.504\" y=\"1156.92\">-2.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 792.859)\" x=\"165.504\" y=\"792.859\">0.0</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 428.802)\" x=\"165.504\" y=\"428.802\">2.5</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:end;\" transform=\"rotate(0, 165.504, 64.7441)\" x=\"165.504\" y=\"64.7441\">5.0</text>\n",
       "</g>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#009af9; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  2260.91,379.56 2260,379.921 2258.16,380.642 2255.41,381.719 2251.75,383.15 2247.18,384.928 2241.7,387.048 2235.33,389.502 2228.06,392.282 2219.9,395.376 \n",
       "  2210.86,398.776 2200.96,402.468 2190.18,406.441 2178.56,410.681 2166.09,415.174 2152.79,419.906 2138.68,424.86 2123.75,430.023 2108.04,435.378 2091.55,440.908 \n",
       "  2074.29,446.598 2056.29,452.431 2037.55,458.391 2018.1,464.462 1997.96,470.628 1977.14,476.872 1955.66,483.18 1933.54,489.536 1910.8,495.926 1887.47,502.336 \n",
       "  1863.56,508.751 1839.09,515.159 1814.09,521.548 1788.58,527.905 1762.58,534.22 1736.12,540.482 1709.23,546.681 1681.92,552.808 1654.21,558.854 1626.15,564.813 \n",
       "  1597.75,570.676 1569.03,576.438 1540.03,582.093 1510.77,587.635 1481.27,593.061 1451.57,598.366 1421.69,603.548 1391.66,608.602 1361.51,613.529 1331.26,618.324 \n",
       "  1300.93,622.989 1270.57,627.52 1240.19,631.919 1209.83,636.185 1179.51,640.319 1149.26,644.321 1119.1,648.192 1089.07,651.933 1059.19,655.547 1029.49,659.033 \n",
       "  999.996,662.395 970.734,665.635 941.733,668.754 913.017,671.755 884.614,674.64 856.55,677.412 828.849,680.073 801.537,682.626 774.64,685.074 748.181,687.418 \n",
       "  722.185,689.662 696.675,691.809 671.675,693.86 647.208,695.819 623.296,697.688 599.96,699.47 577.223,701.166 555.104,702.78 533.624,704.314 512.803,705.77 \n",
       "  492.659,707.15 473.212,708.457 454.478,709.691 436.474,710.856 419.218,711.954 402.725,712.985 387.01,713.952 372.087,714.857 357.97,715.7 344.672,716.484 \n",
       "  332.205,717.209 320.58,717.878 309.808,718.49 299.899,719.048 290.862,719.551 282.705,720.002 275.435,720.401 269.059,720.748 263.584,721.045 259.013,721.291 \n",
       "  255.352,721.488 252.603,721.635 250.768,721.733 249.851,721.782 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6203)\" style=\"stroke:#e26f46; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  249.837,543.644 253.867,606.516 257.897,617.728 261.928,624.208 265.958,628.737 269.988,632.195 274.018,634.977 278.049,637.294 282.079,639.27 286.109,640.988 \n",
       "  290.139,642.502 294.169,643.852 298.2,645.066 302.23,646.167 306.26,647.172 310.29,648.095 314.321,648.946 318.351,649.734 322.381,650.467 326.411,651.151 \n",
       "  330.442,651.791 334.472,652.391 338.502,652.956 342.532,653.488 346.563,653.991 350.593,654.466 354.623,654.917 358.653,655.345 362.684,655.751 366.714,656.138 \n",
       "  370.744,656.507 374.774,656.859 378.805,657.196 382.835,657.517 386.865,657.825 390.895,658.12 394.926,658.402 398.956,658.673 402.986,658.934 407.016,659.184 \n",
       "  411.046,659.424 415.077,659.656 419.107,659.879 423.137,660.093 427.167,660.3 431.198,660.5 435.228,660.692 439.258,660.878 443.288,661.057 447.319,661.231 \n",
       "  451.349,661.398 455.379,661.56 459.409,661.717 463.44,661.869 467.47,662.016 471.5,662.158 475.53,662.296 479.561,662.429 483.591,662.559 487.621,662.685 \n",
       "  491.651,662.806 495.682,662.925 499.712,663.039 503.742,663.151 507.772,663.259 511.803,663.364 515.833,663.467 519.863,663.566 523.893,663.663 527.923,663.757 \n",
       "  531.954,663.848 535.984,663.938 540.014,664.024 544.044,664.109 548.075,664.191 552.105,664.271 556.135,664.35 560.165,664.426 564.196,664.5 568.226,664.573 \n",
       "  572.256,664.644 576.286,664.713 580.317,664.781 584.347,664.847 588.377,664.912 592.407,664.975 596.438,665.037 600.468,665.097 604.498,665.157 608.528,665.215 \n",
       "  612.559,665.271 616.589,665.327 620.619,665.382 624.649,665.435 628.679,665.488 632.71,665.54 636.74,665.59 640.77,665.64 644.8,665.689 648.831,665.738 \n",
       "  652.861,665.785 656.891,665.832 660.921,665.878 664.952,665.923 668.982,665.968 673.012,666.012 677.042,666.056 681.073,666.099 685.103,666.142 689.133,666.184 \n",
       "  693.163,666.226 697.194,666.267 701.224,666.308 705.254,666.348 709.284,666.389 713.315,666.429 717.345,666.468 721.375,666.508 725.405,666.547 729.436,666.586 \n",
       "  733.466,666.625 737.496,666.663 741.526,666.702 745.556,666.74 749.587,666.779 753.617,666.817 757.647,666.855 761.677,666.894 765.708,666.932 769.738,666.97 \n",
       "  773.768,667.008 777.798,667.047 781.829,667.085 785.859,667.124 789.889,667.163 793.919,667.202 797.95,667.241 801.98,667.28 806.01,667.319 810.04,667.359 \n",
       "  814.071,667.399 818.101,667.439 822.131,667.48 826.161,667.52 830.192,667.562 834.222,667.603 838.252,667.645 842.282,667.687 846.313,667.73 850.343,667.773 \n",
       "  854.373,667.816 858.403,667.86 862.433,667.904 866.464,667.949 870.494,667.995 874.524,668.041 878.554,668.087 882.585,668.134 886.615,668.182 890.645,668.23 \n",
       "  894.675,668.279 898.706,668.328 902.736,668.378 906.766,668.429 910.796,668.48 914.827,668.532 918.857,668.585 922.887,668.638 926.917,668.693 930.948,668.748 \n",
       "  934.978,668.804 939.008,668.86 943.038,668.918 947.069,668.976 951.099,669.035 955.129,669.095 959.159,669.156 963.19,669.218 967.22,669.28 971.25,669.344 \n",
       "  975.28,669.408 979.31,669.474 983.341,669.541 987.371,669.608 991.401,669.677 995.431,669.746 999.462,669.817 1003.49,669.889 1007.52,669.962 1011.55,670.036 \n",
       "  1015.58,670.111 1019.61,670.187 1023.64,670.264 1027.67,670.343 1031.7,670.423 1035.73,670.504 1039.76,670.586 1043.79,670.67 1047.82,670.755 1051.85,670.841 \n",
       "  1055.89,670.928 1059.92,671.017 1063.95,671.107 1067.98,671.199 1072.01,671.291 1076.04,671.386 1080.07,671.481 1084.1,671.579 1088.13,671.677 1092.16,671.778 \n",
       "  1096.19,671.879 1100.22,671.982 1104.25,672.087 1108.28,672.193 1112.31,672.301 1116.34,672.411 1120.37,672.522 1124.4,672.635 1128.43,672.749 1132.46,672.866 \n",
       "  1136.49,672.984 1140.52,673.103 1144.55,673.225 1148.58,673.348 1152.61,673.473 1156.64,673.6 1160.67,673.729 1164.7,673.859 1168.73,673.992 1172.76,674.126 \n",
       "  1176.79,674.262 1180.82,674.401 1184.85,674.541 1188.88,674.683 1192.91,674.828 1196.94,674.974 1200.97,675.123 1205,675.273 1209.03,675.426 1213.06,675.581 \n",
       "  1217.09,675.738 1221.12,675.898 1225.16,676.059 1229.19,676.223 1233.22,676.389 1237.25,676.558 1241.28,676.729 1245.31,676.902 1249.34,677.077 1253.37,677.255 \n",
       "  1257.4,677.436 1261.43,677.619 1265.46,677.804 1269.49,677.992 1273.52,678.183 1277.55,678.376 1281.58,678.572 1285.61,678.771 1289.64,678.972 1293.67,679.176 \n",
       "  1297.7,679.383 1301.73,679.592 1305.76,679.804 1309.79,680.02 1313.82,680.238 1317.85,680.459 1321.88,680.682 1325.91,680.909 1329.94,681.139 1333.97,681.372 \n",
       "  1338,681.608 1342.03,681.847 1346.06,682.09 1350.09,682.335 1354.12,682.584 1358.15,682.836 1362.18,683.091 1366.21,683.35 1370.24,683.612 1374.27,683.877 \n",
       "  1378.3,684.146 1382.33,684.419 1386.36,684.695 1390.4,684.974 1394.43,685.257 1398.46,685.544 1402.49,685.834 1406.52,686.129 1410.55,686.427 1414.58,686.728 \n",
       "  1418.61,687.034 1422.64,687.344 1426.67,687.657 1430.7,687.975 1434.73,688.296 1438.76,688.622 1442.79,688.952 1446.82,689.286 1450.85,689.624 1454.88,689.967 \n",
       "  1458.91,690.314 1462.94,690.665 1466.97,691.021 1471,691.381 1475.03,691.746 1479.06,692.115 1483.09,692.489 1487.12,692.868 1491.15,693.251 1495.18,693.64 \n",
       "  1499.21,694.033 1503.24,694.431 1507.27,694.834 1511.3,695.242 1515.33,695.656 1519.36,696.074 1523.39,696.498 1527.42,696.927 1531.45,697.361 1535.48,697.801 \n",
       "  1539.51,698.246 1543.54,698.697 1547.57,699.154 1551.6,699.616 1555.63,700.084 1559.67,700.558 1563.7,701.037 1567.73,701.523 1571.76,702.015 1575.79,702.513 \n",
       "  1579.82,703.017 1583.85,703.527 1587.88,704.044 1591.91,704.567 1595.94,705.097 1599.97,705.633 1604,706.176 1608.03,706.726 1612.06,707.282 1616.09,707.846 \n",
       "  1620.12,708.417 1624.15,708.995 1628.18,709.58 1632.21,710.172 1636.24,710.772 1640.27,711.379 1644.3,711.994 1648.33,712.616 1652.36,713.247 1656.39,713.885 \n",
       "  1660.42,714.532 1664.45,715.186 1668.48,715.849 1672.51,716.52 1676.54,717.2 1680.57,717.888 1684.6,718.585 1688.63,719.291 1692.66,720.005 1696.69,720.729 \n",
       "  1700.72,721.462 1704.75,722.205 1708.78,722.957 1712.81,723.718 1716.84,724.489 1720.87,725.27 1724.91,726.061 1728.94,726.863 1732.97,727.674 1737,728.497 \n",
       "  1741.03,729.329 1745.06,730.173 1749.09,731.028 1753.12,731.893 1757.15,732.77 1761.18,733.659 1765.21,734.559 1769.24,735.471 1773.27,736.395 1777.3,737.331 \n",
       "  1781.33,738.28 1785.36,739.241 1789.39,740.215 1793.42,741.202 1797.45,742.202 1801.48,743.216 1805.51,744.243 1809.54,745.285 1813.57,746.34 1817.6,747.41 \n",
       "  1821.63,748.494 1825.66,749.593 1829.69,750.708 1833.72,751.838 1837.75,752.983 1841.78,754.144 1845.81,755.322 1849.84,756.516 1853.87,757.727 1857.9,758.955 \n",
       "  1861.93,760.2 1865.96,761.464 1869.99,762.745 1874.02,764.045 1878.05,765.363 1882.08,766.701 1886.11,768.058 1890.14,769.436 1894.18,770.833 1898.21,772.252 \n",
       "  1902.24,773.691 1906.27,775.152 1910.3,776.635 1914.33,778.141 1918.36,779.67 1922.39,781.222 1926.42,782.798 1930.45,784.399 1934.48,786.025 1938.51,787.677 \n",
       "  1942.54,789.354 1946.57,791.059 1950.6,792.791 1954.63,794.552 1958.66,796.341 1962.69,798.159 1966.72,800.008 1970.75,801.888 1974.78,803.8 1978.81,805.744 \n",
       "  1982.84,807.722 1986.87,809.734 1990.9,811.781 1994.93,813.865 1998.96,815.986 2002.99,818.145 2007.02,820.344 2011.05,822.583 2015.08,824.863 2019.11,827.187 \n",
       "  2023.14,829.555 2027.17,831.969 2031.2,834.43 2035.23,836.939 2039.26,839.499 2043.29,842.111 2047.32,844.776 2051.35,847.497 2055.38,850.275 2059.42,853.113 \n",
       "  2063.45,856.012 2067.48,858.976 2071.51,862.006 2075.54,865.106 2079.57,868.277 2083.6,871.524 2087.63,874.848 2091.66,878.254 2095.69,881.745 2099.72,885.325 \n",
       "  2103.75,888.998 2107.78,892.769 2111.81,896.642 2115.84,900.623 2119.87,904.717 2123.9,908.931 2127.93,913.27 2131.96,917.742 2135.99,922.356 2140.02,927.119 \n",
       "  2144.05,932.041 2148.08,937.133 2152.11,942.407 2156.14,947.874 2160.17,953.55 2164.2,959.45 2168.23,965.593 2172.26,971.998 2176.29,978.69 2180.32,985.695 \n",
       "  2184.35,993.043 2188.38,1000.77 2192.41,1008.92 2196.44,1017.53 2200.47,1026.67 2204.5,1036.41 2208.53,1046.84 2212.56,1058.05 2216.59,1070.18 2220.62,1083.4 \n",
       "  2224.65,1097.94 2228.69,1114.1 2232.72,1132.29 2236.75,1153.15 2240.78,1177.61 2244.81,1207.28 2248.84,1245.1 2252.87,1297.69 2256.9,1385.77 2260.93,1856.86 \n",
       "  \n",
       "  \"/>\n",
       "<polygon clip-path=\"url(#clip6201)\" points=\"\n",
       "1661.54,312.204 2249.26,312.204 2249.26,130.764 1661.54,130.764 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#000000; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1661.54,312.204 2249.26,312.204 2249.26,130.764 1661.54,130.764 1661.54,312.204 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#009af9; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1685.54,191.244 1829.54,191.244 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:start;\" transform=\"rotate(0, 1853.54, 208.744)\" x=\"1853.54\" y=\"208.744\">function</text>\n",
       "</g>\n",
       "<polyline clip-path=\"url(#clip6201)\" style=\"stroke:#e26f46; stroke-width:4; stroke-opacity:1; fill:none\" points=\"\n",
       "  1685.54,251.724 1829.54,251.724 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip6201)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:48px; text-anchor:start;\" transform=\"rotate(0, 1853.54, 269.224)\" x=\"1853.54\" y=\"269.224\">hilbert transform</text>\n",
       "</g>\n",
       "</svg>\n"
      ]
     },
     "execution_count": 10,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "f = exp(x)\n",
    "xx = range(-0.9999; stop=0.9999, length=500)\n",
    "plot(f; ylims=(-5,5), label=\"function\")\n",
    "plot!(xx, hilbert.(f, xx); label=\"hilbert transform\")"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "## Plemelj theorem for additive jump\n",
    "\n",
    "A key result is that the additive jump of the Cauchy transform is given by the Hilbert transform:\n",
    "\n",
    "**Theorem (Plemelj on the interval III)** \n",
    "The Cauchy transform has the additive jump\n",
    "$$\n",
    "\\CC_{[a,b]}^+ f(x) + \\CC_{[a,b]}^- f(x) = - \\I \\HH f(x)\n",
    "$$\n",
    "\n",
    "\n",
    "**Sketch of Proof** We first show the result for $f(x) = 1$. Note that by deformation through half circles of radius $r$ near $x$ we have (using the notation from last lecture)\n",
    "\\begin{align*}\n",
    "\\CC_{[a,b]}^+ 1(x) + \\CC_{[a,b]}^- 1(x) &= {1 \\over 2 \\pi\\I} \\left[\\int_{\\gamma_x^+} + \\int_{\\gamma_x^-} \\right] {\\D \\zeta \\over \\zeta -x}  \\\\\n",
    "        & = {2 \\over 2 \\pi \\I} \\left[ \\int_a^{x-\\epsilon} + \\int_{x+\\epsilon}^b \\right] {\\dt \\over t - x} \\cr\n",
    "        & \\qquad + {1 \\over 2 \\pi \\I} \\left[ \\int_{\\{x + \\epsilon e^{\\I \\theta} : \\pi \\geq \\theta \\geq 0 \\}} + \\int_{\\{x + \\epsilon e^{\\I \\theta} : -\\pi \\leq \\theta \\leq 0 \\}} \\right] {\\D\\zeta \\over \\zeta - x} \\\\\n",
    "&= {1 \\over 2 \\pi \\I} \\left[ \\int_a^{x-\\epsilon} + \\int_{x+\\epsilon}^b \\right] {\\dt \\over t - x} \n",
    "\\end{align*}\n",
    "This holds true for all $\\epsilon$, including $\\epsilon \\rightarrow 0$, hence\n",
    "$$\n",
    "\\CC_{[a,b]}^+ 1(x) + \\CC_{[a,b]}^- 1(x) = {1 \\over \\pi \\I} \\dashint_a^b {\\dt \\over t -x} = -\\I \\HH 1(x)\n",
    "$$\n",
    "For general $f$, we subtract and add back in like last lecture:\n",
    "$$\n",
    "\\CC_{[a,b]}^\\pm f(x) = {1 \\over 2 \\pi \\I} \\int_a^b {f(t) - f(x) \\over t-x} \\dt + f(x) \\CC^\\pm 1(x)\n",
    "$$\n",
    "therefore\n",
    "\\begin{align*}\n",
    "\\CC_{[a,b]}^+ f(x) + \\CC_{[a,b]}^- f(x) &= {1 \\over \\pi \\I} \\int_a^b {f(t) - f(x) \\over t-x} \\dt + f(x) (\\CC^+ 1(x) + \\CC^- 1(x) ) \\\\\n",
    "    & = {1 \\over \\pi \\I} \\int_a^b {f(t) - f(x) \\over t-x} \\dt - \\I  f(x) \\HH 1(x)\\\\\n",
    "    & = {1 \\over \\pi \\I} \\int_a^b {f(t) - f(x) \\over t-x} \\dt + {f(x)  \\over \\pi \\I }\\dashint_a^b  { 1\\over t - x} \\dt \\\\\n",
    "    & = {1 \\over \\pi \\I} \\int_a^b {f(t) \\over t-x} \\dt = - \\I \\HH f(x)\n",
    "\\end{align*}\n",
    "⬛️\n",
    "\n",
    "Here we test it numerically:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 7,
   "metadata": {},
   "outputs": [
    {
     "data": {
      "text/plain": [
       "(0.0 - 0.4372398225886695im, 0.0 - 0.43723982258866934im)"
      ]
     },
     "execution_count": 7,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "f = exp(x)sqrt(1-x^2)\n",
    "cauchy(f, 0.1+0.0im)+cauchy(f, 0.1-0.0im) , -im*hilbert(f)(0.1)"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "*Examples*\n",
    "Recall  (using $\\diamond$ for the dummy variable)\n",
    "$$\n",
    "{1 \\over \\sqrt{z-1} \\sqrt{z+1}} = -2 \\I {\\cal C}\\left[{1 \\over \\sqrt{1-\\diamond^2}}\\right](z)\n",
    "$$\n",
    "Therefore:\n",
    "$$\n",
    "{\\cal H}\\left[{1 \\over \\sqrt{1-\\diamond^2}}\\right](x) = \\I (\\CC^+  + \\CC^-) \\left[{1 \\over \\sqrt{1-\\diamond^2}}\\right](x) = {1 \\over 2  \\sqrt{x-1}_+ \\sqrt{x+1}} + {1 \\over 2  \\sqrt{x-1}_- \\sqrt{x+1}} = 0\n",
    "$$"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 8,
   "metadata": {},
   "outputs": [
    {
     "data": {
      "text/plain": [
       "0.0"
      ]
     },
     "execution_count": 8,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "x = Fun()\n",
    "norm(hilbert(1/sqrt(1-x^2)))"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Similarly, recall\n",
    "$$\n",
    "\\sqrt{z-1}\\sqrt{z+1} = z  + 2 \\I {\\cal C}[\\sqrt{1-\\diamond^2}](z)\n",
    "$$\n",
    "Therefore,\n",
    "$$\n",
    "{\\cal H}\\left[{\\sqrt{1-\\diamond^2}}\\right](x) = \\I (\\CC^+  + \\CC^-) \\left[{\\sqrt{1-\\diamond^2}}\\right](x) = -x \n",
    "$$"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 9,
   "metadata": {},
   "outputs": [
    {
     "data": {
      "text/html": [
       "<?xml version=\"1.0\" encoding=\"utf-8\"?>\n",
       "<svg xmlns=\"http://www.w3.org/2000/svg\" xmlns:xlink=\"http://www.w3.org/1999/xlink\" width=\"600\" height=\"400\" viewBox=\"0 0 600 400\">\n",
       "<defs>\n",
       "  <clipPath id=\"clip00\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"600\" height=\"400\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polygon clip-path=\"url(#clip7300)\" points=\"\n",
       "0,400 600,400 600,0 0,0 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip01\">\n",
       "    <rect x=\"120\" y=\"0\" width=\"421\" height=\"400\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polygon clip-path=\"url(#clip7300)\" points=\"\n",
       "45.7332,375.813 580.315,375.813 580.315,11.811 45.7332,11.811 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip02\">\n",
       "    <rect x=\"45\" y=\"11\" width=\"536\" height=\"365\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  112.477,370.353 112.477,17.2711 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  179.326,370.353 179.326,17.2711 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  246.175,370.353 246.175,17.2711 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  313.024,370.353 313.024,17.2711 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  379.873,370.353 379.873,17.2711 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  446.722,370.353 446.722,17.2711 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  513.571,370.353 513.571,17.2711 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  53.7519,330.367 572.296,330.367 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  53.7519,284.849 572.296,284.849 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  53.7519,239.33 572.296,239.33 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  53.7519,193.812 572.296,193.812 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  53.7519,148.294 572.296,148.294 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  53.7519,102.776 572.296,102.776 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#000000; stroke-width:0.5; stroke-opacity:0.1; fill:none\" points=\"\n",
       "  53.7519,57.2576 572.296,57.2576 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  45.7332,375.813 580.315,375.813 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  45.7332,375.813 45.7332,11.811 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  112.477,375.813 112.477,370.353 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  179.326,375.813 179.326,370.353 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  246.175,375.813 246.175,370.353 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  313.024,375.813 313.024,370.353 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  379.873,375.813 379.873,370.353 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  446.722,375.813 446.722,370.353 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  513.571,375.813 513.571,370.353 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  45.7332,330.367 53.7519,330.367 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  45.7332,284.849 53.7519,284.849 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  45.7332,239.33 53.7519,239.33 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  45.7332,193.812 53.7519,193.812 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  45.7332,148.294 53.7519,148.294 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  45.7332,102.776 53.7519,102.776 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  45.7332,57.2576 53.7519,57.2576 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:middle;\" transform=\"rotate(0, 112.477, 389.613)\" x=\"112.477\" y=\"389.613\">-0.75</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:middle;\" transform=\"rotate(0, 179.326, 389.613)\" x=\"179.326\" y=\"389.613\">-0.50</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:middle;\" transform=\"rotate(0, 246.175, 389.613)\" x=\"246.175\" y=\"389.613\">-0.25</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:middle;\" transform=\"rotate(0, 313.024, 389.613)\" x=\"313.024\" y=\"389.613\">0.00</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:middle;\" transform=\"rotate(0, 379.873, 389.613)\" x=\"379.873\" y=\"389.613\">0.25</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:middle;\" transform=\"rotate(0, 446.722, 389.613)\" x=\"446.722\" y=\"389.613\">0.50</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:middle;\" transform=\"rotate(0, 513.571, 389.613)\" x=\"513.571\" y=\"389.613\">0.75</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:end;\" transform=\"rotate(0, 39.7332, 334.867)\" x=\"39.7332\" y=\"334.867\">-0.75</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:end;\" transform=\"rotate(0, 39.7332, 289.349)\" x=\"39.7332\" y=\"289.349\">-0.50</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:end;\" transform=\"rotate(0, 39.7332, 243.83)\" x=\"39.7332\" y=\"243.83\">-0.25</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:end;\" transform=\"rotate(0, 39.7332, 198.312)\" x=\"39.7332\" y=\"198.312\">0.00</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:end;\" transform=\"rotate(0, 39.7332, 152.794)\" x=\"39.7332\" y=\"152.794\">0.25</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:end;\" transform=\"rotate(0, 39.7332, 107.276)\" x=\"39.7332\" y=\"107.276\">0.50</text>\n",
       "</g>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:end;\" transform=\"rotate(0, 39.7332, 61.7576)\" x=\"39.7332\" y=\"61.7576\">0.75</text>\n",
       "</g>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  580.315,375.813 579.474,375.241 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  579.474,375.241 577.795,374.097 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  577.795,374.097 575.282,372.387 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  575.282,372.387 571.945,370.114 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  571.945,370.114 567.792,367.287 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  567.792,367.287 562.839,363.914 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  562.839,363.914 557.099,360.005 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  557.099,360.005 550.591,355.574 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  550.591,355.574 543.336,350.634 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  543.336,350.634 535.356,345.2 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  535.356,345.2 526.676,339.29 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  526.676,339.29 517.325,332.923 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  517.325,332.923 507.33,326.117 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  507.33,326.117 496.724,318.896 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  496.724,318.896 485.541,311.281 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  485.541,311.281 473.814,303.296 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  473.814,303.296 461.581,294.966 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  461.581,294.966 448.881,286.319 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  448.881,286.319 435.754,277.38 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  435.754,277.38 422.24,268.179 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  422.24,268.179 408.383,258.743 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  408.383,258.743 394.226,249.104 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  394.226,249.104 379.813,239.29 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  379.813,239.29 365.19,229.333 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  365.19,229.333 350.403,219.264 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  350.403,219.264 335.499,209.116 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  335.499,209.116 320.524,198.919 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  320.524,198.919 305.525,188.706 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  305.525,188.706 290.549,178.509 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  290.549,178.509 275.645,168.36 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  275.645,168.36 260.858,158.292 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  260.858,158.292 246.235,148.335 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  246.235,148.335 231.822,138.521 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  231.822,138.521 217.665,128.881 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  217.665,128.881 203.808,119.446 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  203.808,119.446 190.294,110.244 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  190.294,110.244 177.167,101.306 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  177.167,101.306 164.467,92.658 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  164.467,92.658 152.234,84.3287 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  152.234,84.3287 140.508,76.3439 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  140.508,76.3439 129.324,68.7288 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  129.324,68.7288 118.718,61.5071 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  118.718,61.5071 108.724,54.7018 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  108.724,54.7018 99.3718,48.3341 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  99.3718,48.3341 90.6924,42.4242 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  90.6924,42.4242 82.7125,36.9906 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  82.7125,36.9906 75.4572,32.0504 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  75.4572,32.0504 68.9494,27.6192 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  68.9494,27.6192 63.2096,23.7109 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  63.2096,23.7109 58.2558,20.3378 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  58.2558,20.3378 54.1035,17.5105 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  54.1035,17.5105 50.766,15.2379 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  50.766,15.2379 48.2536,13.5272 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  48.2536,13.5272 46.5742,12.3837 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7302)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  46.5742,12.3837 45.7332,11.811 \n",
       "  \"/>\n",
       "<polygon clip-path=\"url(#clip7300)\" points=\"\n",
       "489.799,62.931 562.315,62.931 562.315,32.691 489.799,32.691 \n",
       "  \" fill=\"#ffffff\" fill-opacity=\"1\"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#000000; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  489.799,62.931 562.315,62.931 562.315,32.691 489.799,32.691 489.799,62.931 \n",
       "  \"/>\n",
       "<polyline clip-path=\"url(#clip7300)\" style=\"stroke:#009af9; stroke-width:1; stroke-opacity:1; fill:none\" points=\"\n",
       "  495.799,47.811 531.799,47.811 \n",
       "  \"/>\n",
       "<g clip-path=\"url(#clip7300)\">\n",
       "<text style=\"fill:#000000; fill-opacity:1; font-family:Arial,Helvetica Neue,Helvetica,sans-serif; font-size:12; text-anchor:start;\" transform=\"rotate(0, 537.799, 52.311)\" x=\"537.799\" y=\"52.311\">y1</text>\n",
       "</g>\n",
       "</svg>\n"
      ]
     },
     "execution_count": 9,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "plot(hilbert(sqrt(1-x^2)))"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "One more example: since\n",
    "$$\n",
    "\\CC 1(z) = {\\log(z-1) - \\log(z+1) \\over 2 \\pi \\I},\n",
    "$$\n",
    "we have \n",
    "$$\n",
    "\\HH 1(x) = {\\log_+(x-1) + \\log_-(x-1) - 2\\log(x+1) \\over 2\\pi} = {\\log(1-x) - \\log(x+1) \\over \\pi}\n",
    "$$"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 10,
   "metadata": {},
   "outputs": [
    {
     "data": {
      "text/plain": [
       "(-0.06387546623297949, -0.06387546623297947)"
      ]
     },
     "execution_count": 10,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "hilbert(Fun(1,Legendre()), 0.1), (log(1-0.1)-log(0.1+1))/π"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "## Real / imaginary parts of the Cauchy transform\n",
    "\n",
    "Suppose $f(x)$ is real. Then we have, for example,\n",
    "$$\n",
    "\\CC f(x) = {1 \\over 2 \\pi \\I} \\int_{-1}^1 {f(t) \\over t- x} \\dt\n",
    "$$\n",
    "is strictly imaginary for $x > 1$ and $ x < -1$.  Furthermore, for all $z \\in \\C \\backslash [-1,1]$ we have the property that\n",
    "$$\n",
    "\\overline{\\CC f( z)} = \\overline{{1 \\over 2 \\pi \\I} \\int_{-1}^1 {f(x) \\over x - z} \\dx} = -{1 \\over 2 \\pi \\I} \\int_{-1}^1 {f(x) \\over x - \\bar z} \\dx = - \\CC f(\\bar z)\n",
    "$$"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 11,
   "metadata": {},
   "outputs": [
    {
     "data": {
      "text/plain": [
       "(0.0029553295890040973 + 0.02871472499617892im, -0.0029553295890040973 + 0.02871472499617892im)"
      ]
     },
     "execution_count": 11,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "f = exp(x)sqrt(1-x^2)\n",
    "cauchy(f, 10.0+im), cauchy(f, 10.0-im)"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "It follows that\n",
    "$$\n",
    "\\overline{\\CC^+} f(x) = - \\CC^- f(x)\n",
    "$$\n",
    "for the whole real line, and in particular we can find the real part of the Cauchy transform via: \n",
    "$$\n",
    "f(x) = \\CC^+ f(x) - \\CC^- f(x) = \\CC^+ f(x) + \\overline{\\CC^+ f(x)} = 2 \\Re \\CC^+ f(x)\n",
    "$$"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 12,
   "metadata": {},
   "outputs": [
    {
     "data": {
      "text/plain": [
       "(1.0996311793408589, 1.0996311793408586, 1.0996311793408586)"
      ]
     },
     "execution_count": 12,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "f(0.1), 2real(cauchy(f, 0.1+0.0im)), -2real(cauchy(f, 0.1-0.0im))"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Similarly, we can find the imaginary part via:\n",
    "\n",
    "$$\n",
    " {}- {\\cal H} f(x) = -\\I(\\CC^+ f(x) + \\CC^- f(x)) = -\\I(\\CC^+ f(x) - \\overline{\\CC^+ f(x)}) = 2 \\Im \\CC^+ f(x) = 2 \\Im \\CC^- f(x)\n",
    "$$"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 13,
   "metadata": {},
   "outputs": [
    {
     "data": {
      "text/plain": [
       "(-0.43723982258866934, -0.4372398225886695, -0.4372398225886695)"
      ]
     },
     "execution_count": 13,
     "metadata": {},
     "output_type": "execute_result"
    }
   ],
   "source": [
    "-hilbert(f, 0.1), 2imag(cauchy(f, 0.1+0.0im)), 2imag(cauchy(f, 0.1-0.0im))"
   ]
  }
 ],
 "metadata": {
  "kernelspec": {
   "display_name": "Julia 1.0.1",
   "language": "julia",
   "name": "julia-1.0"
  },
  "language_info": {
   "file_extension": ".jl",
   "mimetype": "application/julia",
   "name": "julia",
   "version": "1.0.1"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 2
}