{
 "cells": [
  {
   "cell_type": "markdown",
   "source": [
    "# Gross-Pitaevskii equation in one dimension\n",
    "In this example we will use DFTK to solve\n",
    "the Gross-Pitaevskii equation, and use this opportunity to explore a few internals."
   ],
   "metadata": {}
  },
  {
   "cell_type": "markdown",
   "source": [
    "## The model\n",
    "The [Gross-Pitaevskii equation](https://en.wikipedia.org/wiki/Gross%E2%80%93Pitaevskii_equation) (GPE)\n",
    "is a simple non-linear equation used to model bosonic systems\n",
    "in a mean-field approach. Denoting by $ψ$ the effective one-particle bosonic\n",
    "wave function, the time-independent GPE reads in atomic units:\n",
    "$$\n",
    "    H ψ = \\left(-\\frac12 Δ + V + 2 C |ψ|^2\\right) ψ = μ ψ \\qquad \\|ψ\\|_{L^2} = 1\n",
    "$$\n",
    "where $C$ provides the strength of the boson-boson coupling.\n",
    "It's in particular a favorite model of applied mathematicians because it\n",
    "has a structure simpler than but similar to that of DFT, and displays\n",
    "interesting behavior (especially in higher dimensions with magnetic fields, see\n",
    "Gross-Pitaevskii equation with external magnetic field)."
   ],
   "metadata": {}
  },
  {
   "cell_type": "markdown",
   "source": [
    "We wish to model this equation in 1D using DFTK.\n",
    "First we set up the lattice. For a 1D case we supply two zero lattice vectors,"
   ],
   "metadata": {}
  },
  {
   "outputs": [],
   "cell_type": "code",
   "source": [
    "a = 10\n",
    "lattice = a .* [[1 0 0.]; [0 0 0]; [0 0 0]];"
   ],
   "metadata": {},
   "execution_count": 1
  },
  {
   "cell_type": "markdown",
   "source": [
    "which is special cased in DFTK to support 1D models.\n",
    "\n",
    "For the potential term `V` we pick a harmonic\n",
    "potential. We use the function `ExternalFromReal` which uses\n",
    "Cartesian coordinates ( see Lattices and lattice vectors)."
   ],
   "metadata": {}
  },
  {
   "outputs": [],
   "cell_type": "code",
   "source": [
    "pot(x) = (x - a/2)^2;"
   ],
   "metadata": {},
   "execution_count": 2
  },
  {
   "cell_type": "markdown",
   "source": [
    "We setup each energy term in sequence: kinetic, potential and nonlinear term.\n",
    "For the non-linearity we use the `LocalNonlinearity(f)` term of DFTK, with f(ρ) = C ρ^α.\n",
    "This object introduces an energy term $C ∫ ρ(r)^α dr$\n",
    "to the total energy functional, thus a potential term $α C ρ^{α-1}$.\n",
    "In our case we thus need the parameters"
   ],
   "metadata": {}
  },
  {
   "outputs": [],
   "cell_type": "code",
   "source": [
    "C = 1.0\n",
    "α = 2;"
   ],
   "metadata": {},
   "execution_count": 3
  },
  {
   "cell_type": "markdown",
   "source": [
    "… and with this build the model"
   ],
   "metadata": {}
  },
  {
   "outputs": [],
   "cell_type": "code",
   "source": [
    "using DFTK\n",
    "using LinearAlgebra\n",
    "\n",
    "n_electrons = 1  # Increase this for fun\n",
    "terms = [Kinetic(),\n",
    "         ExternalFromReal(r -> pot(r[1])),\n",
    "         LocalNonlinearity(ρ -> C * ρ^α),\n",
    "]\n",
    "model = Model(lattice; n_electrons, terms, spin_polarization=:spinless);  # spinless electrons"
   ],
   "metadata": {},
   "execution_count": 4
  },
  {
   "cell_type": "markdown",
   "source": [
    "We discretize using a moderate Ecut (For 1D values up to `5000` are completely fine)\n",
    "and run a direct minimization algorithm:"
   ],
   "metadata": {}
  },
  {
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "n     Energy            log10(ΔE)   log10(Δρ)   Δtime\n",
      "---   ---------------   ---------   ---------   ------\n",
      "  1   +169.2606726480                   -1.52    361ms\n",
      "  2   +126.0086591069        1.64       -0.82   1.50ms\n",
      "  3   +72.03783082690        1.73       -0.54   1.53ms\n",
      "  4   +13.64379652598        1.77       -0.36   1.50ms\n",
      "  5   +11.28133065189        0.37       -0.31    996μs\n",
      "  6   +11.05645401970       -0.65       -0.64    829μs\n",
      "  7   +10.55002381403       -0.30       -1.06    500μs\n",
      "  8   +9.450164406687        0.04       -0.54    818μs\n",
      "  9   +7.354132105128        0.32       -0.37   1.00ms\n",
      " 10   +4.731470225247        0.42       -0.21    996μs\n",
      " 11   +1.741687359323        0.48       -0.04    983μs\n",
      " 12   +1.492955275769       -0.60       -0.22    817μs\n",
      " 13   +1.323463775139       -0.77       -0.44    673μs\n",
      " 14   +1.269640343253       -1.27       -0.54    664μs\n",
      " 15   +1.227237297763       -1.37       -0.54    637μs\n",
      " 16   +1.176962610400       -1.30       -0.96    654μs\n",
      " 17   +1.150536737251       -1.58       -1.16    666μs\n",
      " 18   +1.146711574614       -2.42       -1.51    633μs\n",
      " 19   +1.145380388805       -2.88       -1.53    654μs\n",
      " 20   +1.144422364102       -3.02       -1.86    658μs\n",
      " 21   +1.144251253364       -3.77       -2.11    636μs\n",
      " 22   +1.144099904777       -3.82       -2.20    660μs\n",
      " 23   +1.144080967071       -4.72       -2.59    676μs\n",
      " 24   +1.144053337500       -4.56       -2.74    677μs\n",
      " 25   +1.144047865766       -5.26       -2.91    699μs\n",
      " 26   +1.144042505909       -5.27       -2.43    679μs\n",
      " 27   +1.144038676988       -5.42       -2.65    639μs\n",
      " 28   +1.144037765972       -6.04       -3.63    664μs\n",
      " 29   +1.144037387767       -6.42       -4.14    483μs\n",
      " 30   +1.144037052483       -6.47       -4.24    499μs\n",
      " 31   +1.144036955532       -7.01       -4.37    658μs\n",
      " 32   +1.144036882112       -7.13       -4.22    640μs\n",
      " 33   +1.144036868946       -7.88       -4.33    661μs\n",
      " 34   +1.144036862487       -8.19       -4.49    639μs\n",
      " 35   +1.144036857093       -8.27       -4.54    679μs\n",
      " 36   +1.144036854705       -8.62       -4.82    664μs\n",
      " 37   +1.144036854215       -9.31       -4.56    635μs\n",
      " 38   +1.144036853822       -9.41       -5.13    662μs\n",
      " 39   +1.144036853032       -9.10       -4.92    660μs\n",
      " 40   +1.144036852876       -9.81       -5.49    641μs\n",
      " 41   +1.144036852797      -10.10       -5.41    653μs\n",
      " 42   +1.144036852774      -10.63       -6.41    487μs\n",
      " 43   +1.144036852766      -11.10       -6.21    636μs\n",
      " 44   +1.144036852759      -11.14       -6.10    660μs\n",
      " 45   +1.144036852757      -11.68       -6.01    657μs\n",
      " 46   +1.144036852756      -12.39       -6.35    638μs\n",
      " 47   +1.144036852756      -12.31       -6.54    651μs\n",
      " 48   +1.144036852755      -12.76       -7.05    638μs\n",
      " 49   +1.144036852755      -13.34       -6.78    656μs\n",
      " 50   +1.144036852755      -13.42       -7.27    664μs\n",
      " 51   +1.144036852755      -13.45       -7.27    637μs\n",
      " 52   +1.144036852755      -14.18       -7.48    660μs\n",
      " 53   +1.144036852755      -14.29       -7.67    667μs\n",
      " 54   +1.144036852755      -15.35       -8.02    627μs\n"
     ]
    },
    {
     "output_type": "execute_result",
     "data": {
      "text/plain": "Energy breakdown (in Ha):\n    Kinetic             0.2682057 \n    ExternalFromReal    0.4707475 \n    LocalNonlinearity   0.4050836 \n\n    total               1.144036852755 "
     },
     "metadata": {},
     "execution_count": 5
    }
   ],
   "cell_type": "code",
   "source": [
    "basis = PlaneWaveBasis(model, Ecut=500, kgrid=(1, 1, 1))\n",
    "scfres = direct_minimization(basis; tol=1e-8) # This is a constrained preconditioned LBFGS\n",
    "scfres.energies"
   ],
   "metadata": {},
   "execution_count": 5
  },
  {
   "cell_type": "markdown",
   "source": [
    "## Internals\n",
    "We use the opportunity to explore some of DFTK internals.\n",
    "\n",
    "Extract the converged density and the obtained wave function:"
   ],
   "metadata": {}
  },
  {
   "outputs": [],
   "cell_type": "code",
   "source": [
    "ρ = real(scfres.ρ)[:, 1, 1, 1]  # converged density, first spin component\n",
    "ψ_fourier = scfres.ψ[1][:, 1];    # first k-point, all G components, first eigenvector"
   ],
   "metadata": {},
   "execution_count": 6
  },
  {
   "cell_type": "markdown",
   "source": [
    "Transform the wave function to real space and fix the phase:"
   ],
   "metadata": {}
  },
  {
   "outputs": [],
   "cell_type": "code",
   "source": [
    "ψ = ifft(basis, basis.kpoints[1], ψ_fourier)[:, 1, 1]\n",
    "ψ /= (ψ[div(end, 2)] / abs(ψ[div(end, 2)]));"
   ],
   "metadata": {},
   "execution_count": 7
  },
  {
   "cell_type": "markdown",
   "source": [
    "Check whether $ψ$ is normalised:"
   ],
   "metadata": {}
  },
  {
   "outputs": [],
   "cell_type": "code",
   "source": [
    "x = a * vec(first.(DFTK.r_vectors(basis)))\n",
    "N = length(x)\n",
    "dx = a / N  # real-space grid spacing\n",
    "@assert sum(abs2.(ψ)) * dx ≈ 1.0"
   ],
   "metadata": {},
   "execution_count": 8
  },
  {
   "cell_type": "markdown",
   "source": [
    "The density is simply built from ψ:"
   ],
   "metadata": {}
  },
  {
   "outputs": [
    {
     "output_type": "execute_result",
     "data": {
      "text/plain": "8.31875022839291e-16"
     },
     "metadata": {},
     "execution_count": 9
    }
   ],
   "cell_type": "code",
   "source": [
    "norm(scfres.ρ - abs2.(ψ))"
   ],
   "metadata": {},
   "execution_count": 9
  },
  {
   "cell_type": "markdown",
   "source": [
    "We summarize the ground state in a nice plot:"
   ],
   "metadata": {}
  },
  {
   "outputs": [
    {
     "output_type": "execute_result",
     "data": {
      "text/plain": "Plot{Plots.GRBackend() n=3}",
      "image/png": "iVBORw0KGgoAAAANSUhEUgAAAlgAAAGQCAIAAAD9V4nPAAAABmJLR0QA/wD/AP+gvaeTAAAgAElEQVR4nOzdd3xV9f0/8Pc5997kZpN9sxMIEEI2hJEBhK0oQh1VcVW01lFbrdZBXbU/tVptraNabEtFxT1AlmySsAlkkJCELLJubvYe997z+f1x/abITOBz7nw9/0rC4X1OHjknr3w+n/P5fATGGAEAADgq0dIXAAAAYEkIQgAAcGgIQgAAcGgIQgAAcGgIQgAAcGgIQgAAcGgIQgAAcGgIQgAAcGgIQgAAcGgIQgAAcGgWCMKCgoI1a9aM/Hi9Xi/btcBIGQwGS18CkMFgwJqIFmc0GiVJsvRVODpJkoxGI69qFgjCoqKibdu2jfz4gYEB+S4GRmhgYAC/gi1Or9fjV7DFGQwGjr+C4fIYjUaObSR0jQIAgENDEAIAgENDEAIAgENDEAIAgENDEAIAgENDEAIAgENDEAIAgENDEAIAgENTWvoCAOASGFFBGyvrZGVtopsTi/Zi8T4U4S5Y+roA7ASCEMB6tQ3SB6XSv0slRhTnLYS6kGGAtjUYj7awmDHCPRPFm8eKSnTrAFwZBCGAlfq8UvrtAePiUHHNbMWMAIGI+vuHnJyUCoVCL9H3p6W3TkhvFEr/mqVI9kXr0Cbp9fq77757cHDQ0hdiGx5//PHU1FQ5KiMIAaxOr4Hu3GMs7WDfzFdODzhPyKlEWh4pLo8UPyyXrtpieHiy4ukkNAxtT3d397fffvvvf//b0hdiA957772CggIEIYBDaB2kJVsNk8cIR5crnS6VbneMFxeGitdsNZzuYe+kKxRoGdoaJyenG2+80dJXYQO2bt0qX3H8FQlgRRr6WMYGw9wg4YNZikumoInGhXYtUVZ0s5/vNBqwNQXA6CEIAaxFj56u2Wq8PVp8KXV0TTsPFW1cpBwwsAf3YXsggFFDEAJYBYnRbbuNKX7C5Y32OYn0+Tzl8Vb2agFahQCjgyAEsAqPHTT2G9h76YrLruCqpK/nK94plr6rQRYCjAKCEMDyNteyr6vZZ/OUVzgpMMRN+Hyu4r4cY30v43RpAPYPQQhgYS0DdG+O8T+zFWOcOFSbHiA8EKv4xV4jkhBghBCEABa2Mtt4W7SQFcRt6sPTiWKvgf5WhA5SMKtnnnnm1VdfNX2s1+vT0tKampoucvySJUsKCwvNcmmXgCAEsKTPKqWabvbHKZc/NHgupUhr5yheOm6s7kazEMynp6ent7fX9PEHH3yQkJAQGBh4kePvu+++VatWmeXSLgFBCGAxPXp67KD097SRThkcubEewm/jFI8cQKMQRm39+vWVlZUffvjhk08+aTQajUbj559//txzz61Zs2ZoaMh0TElJyZtvvrlq1aozv3imd99994477iCinp4e09I5JSUlO3bsIKKtW7eWlZUR0dVXX33o0KHq6mqzfWsXgiAEsJg/HjPODxFmaWRZD+bxBLGkg22sRaMQRue111675pprjh8/7u/vbzQalyxZ8tlnn4WHh+/YsWPx4sWMMSL64IMPBgcHo6Oj169ff911151Voba2tqqqatq0aUTU3t7+0EMPEVF2dvbq1auJ6J133tm3bx8RKZXKtLS0LVu2mPs7PAeWWAOwjJIOtqZMKrxeJVN9J5HenKl4aJ9xXrBSzbPnFeR1+26jtt98f778fKx4z8SzW0TXXHONabTv66+/bmtr27x5syAIK1euTElJ2bNnz5w5c15//XXTkXfeeWdYWFh1dXVkZOTwf8/Pzx87dqxSeel8mThxYkFBAbdv5nIhCAEs44lD0qokRaCLjKdYFCrEeQvvFkuPxqPvx2b8Jk7sMON2FBPHnOeLw2tb5+XlNTQ0LFy40PRpQ0NDaWnpnDlz1q5d+/LLL0uS5Obm1t7eXldXd2YQ9vT0uLm5jeTsbm5udXV1V/gtXDkEIYAFHGpmhe3si3myt9RemSbO+d5wT4zoKVfLEzib6mf5pdOdnZ1NH7i4uMyaNeudd94Z/idXV9empqaHHnqooKAgIiKCiCIjI43Gn6zt5+fn197ebvpYrVbr9fozD+jv73dx+fEPwLa2Nn9/f1m/l5HA34kAFvDUYeNzyaKz/D2WE72ERaEiplLA5Vm4cOG2bdt6enq8vb29vb2dnJyMRmNLS4uzs3NQUBAR7dq1q6am5qz/lZqaWlVV1dPTQ0T+/v6BgYHDe0dotdqjR48mJCSYPs3PzzcNJVoWWoQA5vZDPavrpduizfRn6AtTxCnfGH4VIwbI2Q0Ldik1NfX555+fNm1aSkrK0NBQaWnpzp07J02alJKSMnXq1LCwsIGBgfHjx5/1v7y8vObMmbN9+/Zly5YR0Zo1a1auXCmKYk9PT0pKytNPPz1p0iQi6urqys/PX7RokQW+sZ9CEAKY2zNHjP9vqniFq6mNXIS7cPM48S+Fxlen4Z0ZuLRNmzap1erhTx988ME777zz5MmTarV6/Pjxpl7TTZs2mebCx8fH9/T0uLq6EtGLL74oCD/26z766KNvvvmmKQjnz59fWVn58ssv5+bmfv3118P9omvXrl2xYoW3t7eZv8FzIQgBzGp7Pes10PVRZh2VeDJRTPra8FSiwtvZnKcFm+Th4XHWV9zd3adOnXrmV0RRTExMNH3s6ek5fNjwAYsWLdqxY4dOpwsICCAihUIREBDg5eU1nIJEVFVV9cwzz8jxLYwWghDArF7JNz6RKJr5dYhQN2FphPhuibTqsvZ4ArgMw8utmSQlJZ2ZlET0l7/8xbxXdEF4KgDM53grK+ukn4+1wHP3RKL49gljv8H8ZwYgIpo2bdqtt95q6as4PwQhgPn8v+PS7+JF7guqjcREL2F6gLimHK+PwiXMnz9/8+bNclRmjC1atKikpOQixzz00EPr16+X4+wXgSAEMJNTXWxPo7TynFU8zOaJRPH1Qgn7M8HFPfLII3FxcXJU3rBhg5OTk+mV0Qt54IEHnnrqKdNCbmaDMUIAM3m7WLo3RnS33MT2mQFCgJq+Py1dF4G/gOGC+vr6DAYDEWVnZ+v1+t7e3h07diQmJt51111Hjx5dt25dWFjY/fffb3p9NDc3d8eOHe3t7YmJibfddtvwsmpff/11bm5ubGxscnJybW2taT3S995776677jId8Oqrrz744IM6nS4nJ+f222/ftm2bi4tLRkZGbGysq6vr7t27s7KyzPYt43kAMIcePa0tl34ZY+En7sFY8Z1i9I7Cxbz55psnTpwgoh9++GHlypWbN29OSEj405/+dO+997788ssJCQnffvvtY489Zjp43bp1Go0mNTX1iy++GA6511577ZlnnklMTKypqbn++us/+eQTIhoYGNi1a9fs2bNNxzz99NM9PT3l5eWmZWvWr1+/fft20z/Nnj1bpr7ZC0GLEMAcPjolzQkSI9wtvHrWjWPFxw4aT3awmDGWX8cLzqv1X380djSb7XRuMxa7pS+50L+OGzfu3XffJaKhoaEXXnihpqbG1Ld58803v/XWW0T09ttvE1F/f39WVlZUVNQHH3zg7Oz80ksv7dixIyUlhYiqq6sHBweJqLy83NnZ+eI7FJqMHz/+22+/5fUNjgSCEMAc/lEi/XWG5eezO4l0T4z4jxLpzZmWvxg4L69r75YG+812OqX3xZb6HB4sDAgIGD9+vJOTk+nj1tZW09dfeuml1atXe3l5KZVKo9HY0NDg5ubW1dU1PMswJSVl//79RNTb23vmJMKLcHNzMy3PZjYIQgDZ7Wlkeomygq2iEXb/JDH+K8OLUxVYhts6KQNCLX0J/6NQKM77scnRo0fff//9wsJCT09Pg8Hg6uoqSZKXlxcRdXZ2+vj4ENHw6tsBAQHt7e2SJImiSEQuLi79/f/L+76+vrCwMNPHLS0tI2k4coQxQgDZvVsiPRhr7kn0FxLsKmQFi5+cwkghXKnW1lY3NzfTNPk1a9bo9XoiUqvV8+bNM02Wb2pqWrt2rengqKgoLy+v0tJS06cpKSnffPON6eOOjo7hrlQiys/Pnz59ujm/EbQIAeTVOkg/1EnvZ1hR++uXMeIfjhh/NQl/B8MVmT17to+PT2Jioo+PT0hIyPCqoatXr16xYkVQUFBgYOCiRYtMLT9BEG644YaNGzeapk+8++67K1asaGlp6ejoiI2NXbly5fz584lIkqQdO3Y8++yz5vxGEIQA8lpbLl0TLo5xsvR1nGF+sPDLfspvY4k+VtJMBSuyfft2lUpFRE8//fTwF6+99toFCxaYPg4PDzdtveTs7Lx3797y8nK1Wh0REdHR0WFadzQiIiInJ8fUC3rXXXdNmDDB9B9/85vfXH/99Y888ohCoZg8efLx48c/++yzl1566eDBg8PLfK9fv3769OlRUVHm/JYRhADyWlMm/dXK3kwRBbpzvLCmzCre3wFrM5xJZ77bolKpTOlIRKIojhkzZvjjiRMnmj4e/uLWrVu3bdsWHR197Nix7du3Dy86OmHChAceeKC8vDwmJsb0FW9vbxcXlzM3u6itrX3llVfk+t4uAH0jADI60sK69DRbY3UNr19MED8+JQ0aL30kwGhNnTp18uTJbW1tM2bMKCwsNG1AYXL//fcPpyARRUVF3X777Wf+31//+tfjxo0z37USEVqEALL6T5n0iwmitbwnc4ZIDyHeR9hwWrrBvBtCgSPw9fX9xS9+MZIjx48ff+6+vuaHZwBALgNG+qxCunO89cUgERGtnCj+uwzvjgIgCAFks+G0lOwnhFt6NZkLWR4hHtAxrfmmbgNYKQQhgFw+OcVWjLPeR8xFSdeGi59VoFEIZ6uoqHj00Ud//etfb9q0ydLXYg7W+5QC2LT2QdrVKC2PtOpHbEW0+AmCEH6qoqJi1qxZkZGRaWlpDzzwwEcffWTpK5IdXpYBkMWXVdLCENHLmqYPnmtesHBnDyvrZBO8rLT/1gE9tfvF1v52s51uYVTWDTHXnvmVt95667bbbnv44YeJyNXVddWqVbfddpvZrsciEIQAsvi4QvptnFU3B4lIIdBNY8V1Fey5FAShtfj11Hu7h8y35HSoR/BZXykuLl65cqXp4ylTppSVlQ0vEGqvEIQA/DX0saI2dlWoDfzuWBEtrthlfC7FBi7VQQS7ayx7Af39/d3d3aaPu7q63N3d7TsFCWOEAHL4pIItjxSdbWHZlmn+giDQ4WZm6QsBK/Lhhx8ODQ0R0erVqxcuXGjpy5EdWoQA/H1WIf15mi3EIBER/Xys8EWVlOpvMxcMchs3blxycrJKpTIajRs3brT05cgOLUIAzqq7WU0Pm2V9y6pdyA1R4hdVDE1CGLZ8+fIDBw6sX7++sLAwPDzc0pcjOwQhAGdfVLHlkaLSdp6tRB/BWaS8FkQh/I+Hh4cjRKCJ7TysADbiyyrbW8DzZ5HCl1WYUAhERL/61a9iY2MtfRVmZWOPK4CVq+1lld0sK8hm+kVNbogSP69EixCIiFasWBEdHW3pqzArBCEAT19UsmURttQvapLiJwgCHW9FFoIjsrXnFcC6fVUtXW9r/aIm10cKX6B3FBySTT6xANapsY9OdrB5wTbWL2pyfZT4bTVahOCIEIQA3Kw/LV0dJqps86lK9Re69FTWiSwEh4MJ9QDcfFcj3T3BNmOQSCC6Nlz4roY9nmCTLVpbpFQqu7q6xo0bZ+kLsQHNzc2zZ8+WqTiCEICPHj3ta2KfzbXVICSi6yLEF48ZH0+w4W/Btnh6elZXVw8ODlr6QmyDfPMaEYQAfGyqldIDBQ+Vpa/jCswNFm7dxbT9pHGx9KU4jJCQEEtfAmCMEICT72rYdRG2/UCpRFoYKm48jXdHwbHY9nMLYCX0Em2tk64Nt/kH6rpw4bsavC8DjsXmn1sAa7BXy8Z7CUGulr6OK3ZVmLhXK/UaLH0dAGaEIATgYMNpe2gOEpGXE031E3bUo3cUHIg9PLoAFrepli0Jt5NZB0vCxY216B0FB4IgBLhSpZ2sV08JPnYShEvDhQ2nJSQhOI5RTJ+orKz84YcfvL29ly5d6uJy/ter8/LyDhw44OHhkZWVFRoayukiAazaxtPsmnDBTmKQaJyn4K4S8ltZkq/dfE8AFzPSFmFOTs6UKVNOnDjxr3/9KzMz87wzQB977LGlS5cWFBRs27btn//8J9frBLBeG2ulJWF2lRlLwgT0joLjGGmL8IUXXli1atVjjz1mNBqnTp365Zdfrlix4swDtm7d+tFHHxUVFfn5+clwnQBWqktPh5vZ3GC7GmVYEiY+e9S4KsmuvimACxnRjT44OLhjx45ly5YRkUKhuPbaazdv3nzWMV988cUdd9zR2tq6fv366upq7hcKYJ221UnpgYK7LS8oc65ZQUJxB9P1W/o6AMxiRC1CrVbLGAsODjZ9GhwcnJOTc9YxlZWVRUVFe/fujYmJufvuu99444077rjjvNXa29tPnDjx0ksvDX/llltuuciAol6v1+v1I7lOkI/ppyDYz0AYNxtP08JgMs8tavoRSJLscxsEoiwNbTqtXzFW7lPZHr1eL4poK1uYXq83Go0KheKSRyqVykv+4hpREDLGiGi4liiK5z6KQ0NDQ0NDR44cEUVx06ZNt95664oVK857lQaDYWhoqL29ffgrAwMDF3m2JUkyw5MPF2f6KSAIz8KIttaLj0820x1qzp/C4hBha51wSyQevbOZftj4pWRZ0v/hUm1EQajRaARBaGpqioyMJKLGxsbh1uGw4OBgX19f0x9K6enpnZ2dWq32vOvJ+vv7Jycnv/baayO8xKGhIWdn5xEeDDLR6/XOzs4IwrMUtDE3lXGSn5nuT0mSnJycRvJX8JVbEsmeOWZQOTmL+JmfQxRFlcq+esNtjSiKRqORVzSMqIGvVqvT09M3bdpERIyxLVu2zJs3j4gMBsPp06dNmbxw4cKSkhLT8cXFxS4uLgEBAVwuEcBqbalji0PtMyhC3YQAFyGvFe+Ogv0b6Vujf/jDH2655RadTldSUtLR0XHzzTcTUVVV1YQJE5qbm/38/G699dY33njj9ttvj4+Pf++9955//nn8xQR2b2ud9Gi8OdpnFrE4VNhSy6b62WfSAwwb6ZDvokWLdu3apVKpsrKy9u/f7+bmRkQajebjjz/28PAgIldX14MHD2ZmZkqStHbt2t///vcyXjWAFeg10JFmNltjtzmxKFTcikVHwQGMYmWZxMTExMTEM7/i4eFx6623nvnpL3/5S26XBmDddjWwVH97mzhxplkaoaiNdQzRGCdLXwqAnPASMMBl2lonLQq15yfIWUEzA4WdDWgUgp2z58cYQFZ2/KbMsEUh4tY6vC8Ddg5BCHA5KrpYn4Hi7GXHiQtZHCZsQRCCvUMQAlyOH+rZghD7n1Y50UsQBSrtRBaCPUMQAlyO7fVsQYjd5yAR0dwgYVs9ghDsGYIQYNSMjHY3Sna248SFzA8RtiMIwa45xJMMwNeRZhbiJgS5Wvo6zGJBiLi7UdLj1VGwXwhCgFHb5jD9okTkp6YoD+FwMxqFYLcQhACjtr1BWhDiQM/OghAME4I9c6CHGYCLXgPltbBM+11Z7VzzQ8TtmFYP9gtBCDA6exrZFD/BbRSrE9q8WRohv5V1YXtssFMIQoDR2V4vzXekflEiUitoeoCwpxGNQrBPjvRnLQAPOxvY+xlmDcIh49CO6r3FrWWlLadcnVwm+kQnBcbPDJlqzmuYGyzuamDXhpvznABmgiAEGIXWQaruYVPMuEVfdu3+d/L+HeUVMTUoaVbQDKMoneqo/iB/7afFX/966r3R3lHmuYy5QcIvc9AiBPuEIAQYhR310iyNqDRLg1Bi0p8P/L20reL303+dokkgov7+ficnpxkhU2+N/dn3p354bOdzKxNXXBu9yAwXM9VfqO1lTf0U6GKGswGYFYIQYBR2NbK5weZoDkpMenn/3zoGu95f/Lqz4uz9AEVBXDp+8bTglEe2/2HAMHhjzFK5r0chUEaguKdRummsY42PgiPAPQ0wCjsazBGEEmPPZf+5V9/30qxV56bgMI1bwN/m/79vyzZ9eXKD3JdERHODhZ0NmE0IdghBCDBS9b2sfZDFecsehJ+WfN0+0PnHzKdUCtXFjwx08//r/D99XPxlYXOx3Fc1N1jY2YggBDuEIAQYqe0NbG6wKMqcg8UtZZ+XfPdM+qNKUTGS4wNc/Z6c8fAfc/7SOdgl64XF+whdQ6ymB1kI9gZBCDBSO+XvF+0a6n4+59UnZjwc6BYw8v81PXjK3IjMl/f/Tb4LIyKBaE6QuBuNQrA7CEKAkdrVwOYGyRuEHxz/KD009TLmCN6bdEdLX9vOmmw5rmrY3GBhBxYdBbuDIAQYkVNdjBGN95IxCE+1V+6t3feL+Fsv4/8qRcVvU3/1j7z/DBgGuF/YsLnBAlqEYH8QhAAjsquBZcncHHzryAf3JN7u6exxef89zj8mMSDuk+Kv+F7VmaI9BYmoogtZCHYFQQgwInu0bLacQbijem+/YeDqcQuupMivku/8tmxzY08Tr6s615wgNArB3iAIAUZkdyPLku1NGSMzfpD/0UNTVorCFZ3Cz9X3+onXrClcx+vCzoUgBPuDIAS4tPJOJhKN9ZArCLdX7dG4ByQETL7yUjfELN1ff6S+u/HKS53XnCBMqwd7gyAEuDRZm4MSY+uKv74j7udcqrmpXK8bf9WnJd9wqXauaE9BKdIpDBOCHUEQAlza7kYZBwh3n85xUbkkB8bzKnjjpKW7a3KbenW8Cp5ltga9o2BXEIQAl7ZHK9cro4zYR0Vf3BV/M8eank4eS6IXfFryLceaZ5qNYUKwLwhCgEso62QiUZQ8A4RHGo8Lgjg9eArfsjfFXLe9ak/PUC/fsiZZwcLOBuxNCPYDQQhwCXvkHCD8pmzjzyYu4V7Wx8V7evCUzZU7uFcmorEegkoUyjrRKAQ7gSAEuIQ9WjZHnn7Rpt7mouaT8yIy5Si+fOLV35RtlJgscTVbI+zVIgjBTiAIAS4hW8tmaWQJwu/KNy8aO1etVMtRfLJfjLvK7aj2uBzFZwUJezFMCPYCQQhwMZXdzMhonCf/INQb9ZsrdyyNXsS98rDrJlz1TdlGOSrjxVGwJwhCgIvZ08hmy9Mc3HU6N3pMVJhniBzFTeZFzCpqPinHPIrxXoKRUXU3shDsAYIQ4GL2yrbE6ObK7deOl7E5SERqpfO8yFlbKnfJUXxWEIYJwU4gCAEuZo88U+mbepsr2qtnBo9638HRWhw1d2vVTkb8E2uWRtiD3lGwCwhCgAuq62W9BjZBhj0If6jaNTciU6VQca98lom+0U4KpxPNpdwrz0aLEOwFghDggnY3stkaUY6O0R+qdi8amyVD4fNYGDVniwwTCieNEbr0rL4XWQg2D0EIcEEy9YueaDnJiE3yncC98nktjMrac3rfoHGIb1mBKCNQRKMQ7ACCEOCCZJpBuLVy1+Kx87iXvRA/F5+JvtH76g5xrzw7CMOEYA8QhADn19RPzQNssjfnIDRIxt2ncxdGzeFb9uIWjc36oWo397KZGiGnCUEINg9BCHB+e7VShkbkPkJ4VHs83DM0wNWPc92Lygidka8r6h7q4Vs20Udo6GO6fr5VAcwNQQhwftlalilDv+iu07lZEency16ci1KdoknYX3+Yb1lRoJkBwj4ddqIA24YgBDg/OQYIDZJxX92h2WFpfMuOxJzw9F01udzLZmrEbLwvAzYOQQhwHp1DVNnFkn05B6GpX9TP1Zdv2ZFID51e0HyiV9/Ht2ymBqtvg81DEAKcR04Tmx4gqHg/HxbpFzVxUaqTA+Nz6w7yLZvqL5R1sm4936oAZoUgBDiPbK2UqeH8dFiwX9REjt5RJ5FS/IR9eHcUbBmCEOA89jbyf1Mmryk/3DPEIv2iJmkh0/J1Rdx7R2dphGwt3pcBG4YgBDhbv4EK29k0f85BmFN7MDNsJt+ao+KqckkIiD3UkMe3LN6XAVuHIAQ424FmluAjuCp51mTE9tUfTgtJ5Vl09NJDp++r57zEzMxAIa+VDRj5VgUwHwQhwNnk6Bcta61wU7nIug3vSKSFTNtff8Qg8UwtNyVNGiMcbkajEGwVghDgbDlNUkYg50cjt/5geuh0vjUvg6+Ld6hHcGFzMd+yGYECekfBdiEIAX7CINEhHUsL5D9AmB46jW/Ny5MeOo37JIpMjZDThPdlwFYhCAF+4lgri/AQfJx51mzq1bUNtE/ynciz6OVKD52eXXuAb81MjbiviRnRJgTbhCAE+ImcJpbBvTlYdzAtZJooyLHF76iNHRMhCmJVRw3Hmn5qCnIVCtuQhGCTEIQAPyHHWtv76g+nWUe/qElaaGou7wW4MUwItgtBCPA/jCi3ScrgGoQDhoGSlrKUwASONa/Q9OApBxuO8q2JvQnBdiEIAf6nrJO5KIQwN55BmNdUEOM73lXlwrHmFUoKiDvVXtk11M2xZqZG2NuI92XAJiEIAf5Hjn7RA/VHpwdP4VvzCjkpnBL8J+dpCzjWjPIQlKJQ0YVGIdgeBCHA/+Q2Mb79okR0uPGYtQUhydM7mqHBMCHYJAQhwP9kazm/Mnq6q84gGSK9wjjW5GJmyNSDDUcZ8cyt9EAhF8OEYIMQhAA/0vZT+yCbNIZnEB5oODojZCrHgrwEuQe6KF0q2qs51sT7MmCjEIQAP8rWShkaUeTaM3qoIW9aUArPivxw7x2N9xa0fUzXz7EkgDkgCAF+lKNl6Vz7RQcMg8UtpVM0iRxrcjQjZMoBrkEoCjQzUMjFWmtgaxCEAD/ivqZMvq5ovPdYq5o4cabEgLhT7ZV89+lNDxQxTAg2B0EIQETUo6eyTjbFj2cQHmk8PjUoiWNBvpwVTpN8J+TrijjWxDAh2CIEIQAR0T4dm+InOCt41jysteogJKKpQUmHG49zLDjdXzjRzvoMHEsCyA5BCEBElKOV+PaLtvW3t/S1TvSJ5liTu6mapCNcg9BZQQk+wkFs0gs2BUEIQPTjm1iz8i8AACAASURBVDI8H4cj2uMpmgRRsOpHLNp7bNdgt66vhWNNrL4NNseqn1IA89BLdKSF82a8R7T5UzVW3S9KRKIgpGgSjmrzOdbM0Ag5Wrw4CrYEQQhAR1vYOE/By4lnzTxtgZUPEJpM0STy7R3NCBQP6pgBUQi2A0EIQLm8J05UdZ5WCGKwu4ZjTZmkBiUd1R7nuNaatzOFuQsF2KQXbAeCEIBymzhPpT/amJ8alMyxoHwC3QLcndwr2qs41swIxCQKsCUIQnB0ps1407luOnFUm2+1C8qca4omke+WTOkarL4NtgRBCI6urJO5KnluxisxqbC5OCkwjldBuSUHxh9rKuRYMD0Qm/SCLUEQgqPL4b31UmnbqQBXP2/1GI41ZZUSmJCvO2FkRl4Fx3oISlGo7EajEGwDghAcHfcBwjxtQYomgWNBuXk6e2jcAkpbKzjWTAsUcjCbEGwEghAcXQ7vXemPNRUmB9pSEBJRiibhWBPXYUJs0gu2A0EIDq2pn1oGWCy/zXj1kqG4pTQhIJZXQfNIDkzgO0yYgRYh2I5RBOGnn356zz33rFq1SqvVXuSwzz///K233rriCwMwhxytlB4ocNyMt6SlNMwzxMPJnVtFs0gKjDvRclJv1PMqmOgrNPSx1kFe9QBkNNIg/Pvf/75q1arZs2e3trZmZGQMDp7/Bj948OBDDz304osv8rtCABnlNrE0rkuM5jUVJAfGcyxoHm4q13DP0JLWMl4FFQJN8xf2YZNesAUj+hVgNBpff/31d9555/bbb3/vvfc8PDy++uqrcw8bHBx84IEHXnjhBd4XCSAX7pvxHtMWptjaAKFJcmB8HtdhwjRs0gs2YkRBWFdXV1tbm5WVZfp0zpw5ubm55x724osvLlu2LDbWxkZHwGH1Gqikg6X6cwvCIeNQadupeFsbIDRJ0SQca+K9SS+GCcEWKEdyUGNjo7u7u7Ozs+lTf3//I0eOnHVMfn7++vXrDx8+fODAgYtXq6ury87Ovv7664e/8vvf/z4+/oK9Sf39/QoF1/1SYfT6+voEQRAEno0ni9vTJCZ4i9JgXx+nggUtxZGe4WxI6hviVfIn+vv7DQaDTI/DOLfI0tbyju4OJwWf1ccT3Cm/TdXW3ae2r8d3cHBQFEWVSmXpC3Foer3eaDRK0qX73tVqtSheosk3oiBUq9VDQ0PDnw4NDbm4uJx5gMFguOeee957773hsLwIHx+fiIiIn//858NfiY6OVqvVFzper9df5F/BPAwGg1qttrMgPNTGMjWM491V3F6WoomX73ZljDk5OckUhGpSR40Jr+qtTQyYzKkgxY6RCrudM7nOTrE4QRAQhBanUCiMRuNInrVLpiCNMAhDQkIGBwebm5v9/f2JqLa2NiQk5MwDKisr8/Pzb7/9diIaGBhoa2sbN27c9u3bo6Kizq3m6uoaHh5+0003jeTURCSK4ki+E5CV6adgZ0GYqzP8ZrJC5PfOaH7ziRWTr5fvdhX/j0z1EwPiCppPJGu4veyToWH7moXZwXb1/Mr9U4CREEWRMcbrpzCiKv7+/pmZmR9//DERdXR0bNq0afny5UTU1tb2zTffEFFUVNTJkye3bdu2bdu2119/3cvLa9u2baGhoVwuEUAOBokO6Xhuxqs36ktby+P8JvEqaH5JgXHHdTyHCdMDsUkv2IARtQiJ6JVXXlm2bNmuXbuKi4sXLlw4c+ZMIjp58uTPfvYzxphKpRo7dqzpyNraWoVCMfwpgHU63sbC3QWfS/flj1Rxa1nUmHBXlculD7VWCf6TX8h5TW/UqxR8+v0yNeLde41GRgq76koAezPSIExLSzt58uSBAweCgoKSk3/caC0lJaW0tPSsI6dPn3748GGe1wgggxwt55XVjjcVJgbYzI4T5+WqcgnzDDnZVh7vz+fFV381BboIJ9pZgg+SEKzXKDpYfXx8rr766uEUJCK1Wj1hwoSzDlOr1REREXyuDkA23NfaPq4rsqGtly4kOSD+ONdJFFhrDawfxnvBQeU2SRyn0uslw8lWbg0pC0rkPUyYocFu9WDtEITgiE51MVEQIj24BeHJ1rIwzxA3lSuvgpaSGDC5pKXMIHHbmzAjUMhGixCsG4IQHFG2ls3iOkCYrytOsvEBQhM3lWuIR1BpWzmvguO9BCNjNT3IQrBeCEJwRLm8lxjN1xXZ+psywxIDJhfoijkWTA8UMUwI1gxBCI4om+sroxKTiltK4/1teAbhmRICJhfoTnAsmBGIYUKwaghCcDjNA9TUzyZ7cwvC8vbKAFc/T2cPXgUtKyEgtqC5WGLcJsJnaDBMCFYNQQgOJ1srpQcKHKd45+tOJHBan9MajHH28nXxqeyo4VUwyVeo7cEmvWC9EITgcHK0LEPD884v0BXbUxASUWLA5Hx+vaMKgaYHYJNesF4IQnA4OU0sk9+bMoxYUXMxrx0brERCQCznYUIN3pcB64UgBMdi2ox3Kr/NeGs661yULn4uPrwKWoOkgLh83QlG3KIrE8OEYMUQhOBYDuhYkq/AcavYAt2JRNtfWe0s/q5+zgqnuq4GXgWn+wuF7azfwKseAE8IQnAs2VqJ71T6Al1xou2vrHauxMA4jsOErkqK9xYONaNRCNYIQQiOJVvLMrm+KZOvK4oPsMMgTPCfxDEICb2jYMUQhOBA9BIdbmYzA7i1CLW9Or1kCPUI5lXQesT7xxY281xfJlMjZmOTXrBKCEJwIHktLNpT8HLiVrBAZ2/viw4L9wodMAzo+lp4FUwPFA7omAFRCNYHQQgOZK+WZXIeIDyRYI/9okQkkBDnP6mwuYRXQW9nivAQjrehdxSsDoIQHEg27yAsbC62gz0ILyTeP7aQ6+rbGCYE64QgBEfBiPY1SemB3O75zsEuXV/LOO9IXgWtTUIA72FC7E0IVglBCI7iRDvzVQtB/LbOLWwuifObpBD4zUm0MhN8xjX2NPUM9fIqODtIzNFKSEKwNghCcBT8+0V1xfEBdrL10nkpBMVE3+iiFm7DhEGu5KESSjsQhWBdEITgKLK1nDfjLWi2q00nzivBfzLfYcJZQcJe9I6ClUEQgqPI0bLZQdyCcMAwWNlxOsZnPK+C1ikhIDafaxBmYJgQrA+CEBxCRRdjRFEe3IKwpLVs3JhItdKZV0HrFOs38VR7pd6o51VwlkbY3YggBOuCIASHsJdrc5DsegbhmVyU6nCv0JNt5bwKjvcSJEbV3chCsCIIQnAI3KfSFzaXOEIQElGCf2wBhgnBriEIwSHsbWQcN52QmFTSWhbnZ8+vjA6L948t4re+DGE2IVgfBCHYv/pe1q1nMWO4BeGp9ip/F19PZw9eBa1ZQkBsYXOJxLhFF1qEYG0QhGD/9mjZrCCRY8doYXOxXW69dF7e6jFezp41nad5FZzsLbQNsoY+ZCFYCwQh2L9sLc9+USIq0BXH+ztEv6hJfEBsAb+11gSijEAxB41CsBoIQrB/e7gOEBJRUXOJHa+1fa54/0mFOq7DhBr0joIVQRCCnWseIG0/S/DhFoQNPVoShCD3QF4FrV+CP88WIRHNDhL2YDYhWA0EIdi5PY1SRqDIcYSwQFec6EjNQSIK8wzRS3qOm/Qm+Qr1faxlgFc9gCuCIAQ7t6eR81T6wubiOAcLQiKK84vhuOioQqCZAcJeLbarB6uAIAQ7t5t3EBboihPsetOJ84r3n8S7d1RE7yhYCQQh2LPWQarrZcm+3IKwc7Crtb9t7JhIXgVtRXwA593q52CYEKwGghDs2Z5GKT1QUHAdIJzsHyMKDvfgTPAZp+3VdQ118yqY4itU97DWQV71AC6fwz3P4FD2NLLZQTxv8qLmkgTHGyCk/9ukt7illFdBpUgzA4RsDBOCFUAQgj3jP0DoSGvKnCXBf3Ih10VHMUwIVgJBCHarbZCqu1kKvwHCQeNQZUeN3W/GeyHx/pP4bkOBYUKwEghCsFt7tVK6RlDyu8eLW0odYTPeC5nsH1PWVjFkHOJVcKqfUNnN2jBMCJaGIAS7tbuRzdbwvMMLdMUOsgfhebko1RGeoaVtFbwKKkWagWFCsAIIQrBbuxpYVjDnJUbjHGmt7XPFB0wq5DqbcE6QuBu9o2BpCEKwT62DVNPDc4BQYtKJlpMOtenEueL9Oc8mzAoSdjYgCMHCEIRgn3Y1SBmBPAcIKzqq/V19vZw9uVW0QQn+nDfpneInnO7BoqNgYQhCsE+7GllWMN8BwhMOtfXSefm4eHs6e3DcpFcpUnqgsKcRw4RgSQhCsE+7GlgW/yVGHT0IiSghYDLfRUezgsVdGCYEi0IQgh3S9ZO2nyXyGyAk05oyAZM5FrRRCf6x+boTHAtmBQm7MEwIFoUgBDu0s0GarRE5LjFa391IgqBxC+BW0WYlBHAOwiRfoamfafs5lgQYHQQh2KFdjZwnThQ0FyeiOUhERKEewYxJTb06XgVFgTI04u4GDBOCxSAIwQ7tbmRzOA8QnnDMtbbPK94/Np/rJIq5wQKGCcGCEIRgb+p7Wdsgi/Pmuys9Bgj/J94/lu+0+rnBwg4ME4LlIAjB3myrZ/OCRZFfDnYMdnYMdEZ6hXOraOMSAmILuA4TTvYWevWsuhtZCJaBIAR7s6OBzeM7QKgrjvOPEQWeNW1atHdUc19r52AXr4ICUVawiEYhWAqCEOzNrkY2L4T3ACH6Rc8gCuJk/xi+exPOQ+8oWA6CEOxKcQdTCjTWg2cQ5utO4JXRsyT6T+bbOzo/RNjRICEJwSIQhGBXdtSzBVybg736vrruhgk+0Rxr2oHEwMl8N+mNcBc8VMKJdkQhWACCEOwK9wHCouaSGN/xKlHJsaYdiPEZX915uk/Pcxr8vGBhRz2CECwAQQj2w8goWyvN4b3WdoI/+kXPplKoJvpEF7eUcqyJYUKwFAQh2I8jzSzETdC48KyZr8OaMueXEDC5oJnnMOHcYHGvVtJjhRkwOwQh2I/tDWw+137RIeNQRUdVrN9EjjXtRkJAbH4TzyD0U9NYD+GgDo1CMDcEIdiPH+qkhaE8b+niltIorwi10pljTbsR5z+ptO3UkHGIY80FIcJ2LDoKZocgBDvRrae8Vpap4TtxAv2iF+SiVEd4hZW2VXCsuSBE/KEOLUIwNwQh2IndjdJ0f8GN69udBc0n4rEZ74UlBkzmuyVThkYoamftgxxLAlwaghDsxLZ6tiCE5/1skIwlLWXYdOIiEgIm5zcVcSyoVlBaoLC7Eb2jYFYIQrATP9Rxnkpf2lYe4hHk7uTGsaadSQyYfKLlpEEycqy5IETchtmEYF4IQrAHdb2sdZAl+fIMwuNNRUkBcRwL2h8PJ/cg98Dydr7DhMIPCEIwLwQh2IMtdWxBCM+tl4jouK4oMRBBeAmJAXHHufaOxvsIfQZWiS2ZwIwQhGAPtvNeYtTIjCeaT2KA8JKSAjkHoUA0PxjvjoJZIQjB5hkZba+XFoXyDMKytoog90BPZw+ONe1SUmBcUUuJxHi+3rIoFL2jYFYIQrB5h5pZqJsQ7Mp5gDARA4Qj4Onk4e/qV95eybHmolBxV4M0hFdHwVwQhGDzttRKV4Vx3j4+X1eUhAHCkUniPUzop6YJXsK+JjQKwUwQhGDzNtexxVxXVpOYVNR8Mh4DhCPDfZiQiBaHClvq0CQEM0EQgm1rGaDyTpYWyLNFWN5e6efi46324ljTjiUGxBU2F/MdJlwcJm6pRYsQzARBCLZta52UFSyquN7IedqCZE08z4p2zVvt5evizXeYcJq/UN/H6nuRhWAOCEKwbVvq2FVc3xclomNNhcmBCXxr2rcUTcKxpkKOBRUCzQ8Rt+LdUTALBCHYMInRNt4TJ4zMWNRcgk0nRiU5MOGYlmcQkmmYEL2jYBYIQrBhR1qYn1oId+e6xGjrqSB3jZezJ8eadi8pMK6wuZjvoqOLQ8XtDdiwHswBQQg27PvT0jXhnPtF87QFKYEYIBwdTyePIHdNWdspjjUDXSjaU8jBJAqQH4IQbNj3p9mSMM73cF4T3pS5HCmB8XnaAr41l4SJG0+jSQiyQxCCrWroY6d72MwAni1CvWQ42Vqe4I8BwlFL1sTnNXEOwmvChe9Po0UIshtFEObm5i5cuDApKemJJ54YHDx7D2mdTvfUU0/NmjUrNTX1oYceampq4nqdAGf7/jRbFCoquf4tV9xSGuYZgj0IL0NiQFxJa5neqOdYM8VP6DFQeSeyEOQ10t8iLS0tS5Ysufnmmz/55JP9+/c///zzZx1QVlY2MDDwwgsvrF69ur6+fvny5ZyvFOCnNtayJbwHCI83FSZjgPCyuKlcIzzDilvLONYUiK4OEzbi3VGQ2UiDcO3atampqXfffXdsbOyf//zn1atX6/U/+dMvIyPjr3/9a1ZWVlJS0uuvv75///7u7m4ZLhiAiGjQSHsapUVcV1YjoqPa/BQNZhBephRNwlFtPt+aS8KEjbUYJgR5jfT3SGFhYWpqqunjqVOntrW11dfXX+jgvLy84OBgDw9sYQNy2dnAEn0EX2eeNQcMA+XtldiD8LJN0STm8Q7C+SHiQR3rHOJbFeAnlCM8TqfTxcTEmD5WqVTu7u5NTU2RkZHnHtnQ0PDwww+/+eabFyp16tSp7777Lioqavgr77//flpa2oWO7+3tFQTOPWAwWr29vYwx6/lBfFupnB9IPT0DHGsebjoW7RVlGDD0UA/Hshz19/c7OTkpFApLX8j5RbmEnWqv0nU0uypdOJad4afaUNG3LMxa2oWDg4OiKKpUKktfiEPT6/VGo9FgMFzySFdXV1G8RJNvpEHo5eXV29tr+liSpL6+vjFjxpx7mE6nmz9//sMPP3zjjTdeqNS4cePmzZv317/+dfgrERERF3m2GWPu7u4jvE6Qj5ubm5UEISPa3GDYdrXC3V3NsWxxadm0kBRrvtkUCoU1ByERxfpNPNVblRYyjWPNn42VtjYpb5tkLd+1SqVCEFqcKQjVaj6/AUbaNRoVFVVW9uMweEVFhUKhCA0NPeuY1tbWBQsW3HTTTU8++eRFSgmC4O7uPvYM1vxggxU62sLcVTTRi3MqH9Een6JJ4lvT0UzRJB5p5Nw7ujRC2FiLJWZARiMNwhUrVmzatKmiooKI3nzzzWXLlrm5uRHR2rVrt2/fTkRtbW0LFixIT0//zW9+097e3t7ebjTyXG8JYNh3NdJ1EZxTsH2go7mvNcY3mm9ZRzNVk3RUe5xvzWBXYYKXsFeLd0dBLiMNwkmTJj377LNTpkwJDQ09cODAX/7yF9PXN2zYsH//fiLau3dvdXX1p59+Ou7/VFdXy3TR4OC+rWbXRXB+X/SI9nhSYLwoYImJKzLeZ1z7QGdzXwvfstdFiN/VoEkIchnpGCER/e53v3vwwQe7u7v9/f2Hv/j555+bPli2bNmyZcs4Xx3AOSq6WMsAm+bPuUV4VFswBRMnrpgoCMmB8XnagkVj53IsuyxCWLxFenMmWcUYNdid0f39q1arz0xBAPP7toZdFyGKvH8j5mkLMEDIRYom4SjvtdYmjRFcFHS8Fb2jIAt0BIGN+a5G4t4verqrnoiFe4bwLeuYUoOSjzQeY8Q5tJZGCOgdBZkgCMGWNPVTUTubG8y5PXioIW9aUArfmg4r2F3jonSp7KjhW3ZZhPhVFVqEIAsEIdiSr6ulq8NEZ97TbQ415k0LRhByMy04+VBDHt+aMwOFziEq6UAWAn8IQrAlX1ZJ10dybg4OGYeKmkumaBL5lnVkqUEphxuP8a0pEC2LFNAoBDkgCMFmtAxQXgtbzHuh7XzdiXHeUW4qV75lHVlKYHxJa1m/gecCeER0Q5T4VTWGCYE/BCHYjK+rpUWhossopvyMyKHGY9OCkjkXdWxqpXqS74RjTYV8y2YECrp+OtWFRiFwhiAEm/FVlXRDFP+JZIfxpowMUoOSDzdyHiYUBVoWKXyJ3lHgDUEItqF9kA42s6vCON+xzX0t7QOd433G8S0L04JTDjVwHiYkousjxa+q0DsKnCEIwTZ8UyPNDxHdePeLHmzImxqUJFrHrhr2ZOyYiAHDQEOPlm/Z2UHC6V5W1Y1GIfCEIATbsK5CumUs/7g60HBkRsgU7mVBIGF68JT99Yf5llUIdEOU+GklghB4QhCCDdD109EWdjXvflG9UX+sqXB6EIJQFjNDpu6vP8K97C1jxY/K0TsKPCEIwQZ8WildG87/fdFjTYVRXhGezh6c6wIREaUGJRe3lPbp+/mWTdcIfUYqakejELhBEIINWFch3TKO/726v+HwzJCp3MuCiVqpjvWbyH17QoHopihhXQUahcANghCsXU0Pq+xm83mvL0pEB+qPpoWkci8Lw+TqHR0nflrB0CQEXhCEYO0+qWA3RIlK3rdqVedpIzNGjYngXBfOkBYybV/9IYl3ZiX5CmoFHdIhCoEPBCFYu3UV0i1j+d+oB+qPpIVM414WzhTkHujp7FneVsG98s3jxI/ROwqcIAjBquW1sG49pWvk6BfFxAlzSAtJ3cd7EgUR3R4trKuQhhCFwAOCEKzaf8ulu8Zz346eOge7TnVUpQQm8C4MZ0sLSc2pO8i9bKSHEDtG2FyLJAQOEIRgvQwSfV4p3RbNvzmYW3coNSjZSeHEvTKcJc4/tm2gnfsSM0R05wTxv+UYJgQOEIRgvTbWShO8hHGe/IMwu/ZAZugM7mXhXKIgyNQovClK3N0oNXPe6wkcEYIQrNeH5ezOCfxv0QHDYL6uaHowBgjNJDNsRk7tAe5l3VV0dZj4KV6ZgSuGIAQr1TpIOxukG6L436IHG45O9otxd3LjXhnOa0pgYkVHdVt/O/fKd40X/4vl1uCKIQjBSq0tl64NFz1V/Ctn1x7IDEO/qPmoFKrUoOT9Dfxn1s8NFloG6HgrRgrhiiAIwUp9UCrdG8P//jRIxoONR9NDp3OvDBeRGTYzu3Y/97KiQHdPFFeXolEIVwRBCNYoR8uMjDJkmD54vKkwzCPE18Wbe2W4iBnBUwp0xb36Pu6VV04QPq2Q+gzcC4MDQRCCNVpdKv0yhv/0QSLaWZOdFZEhQ2G4GDeVa1JgXG7dIe6VQ9yE9EDx80o0CuHyIQjB6nQO0foa6bZoWfpFc+oOzgqbyb0yXNLciMxdNTlyVL43RkDvKFwJBCFYnbWnpKvCRH81/8pHtMcivMIC3fz5l4ZLSQ+dnq8r6hrq5l756jDxdA92KITLhyAE68KI3i2WfjVJljtzZ01OVjj6RS3DRameokmUo3dUIdA9E8V3i9EohMuEIATrsr2eKQSaJcNrMnqjfl/doVnh6Be1mLkRmTtrsuWofN8kcV2F1D4oR22wfwhCsC5vnZB+GyfLbXmwMS/aO8rPxUeO4jASM0NSi1tKOwY7uVfWuNDVYeIaTK6Hy4IgBCtS3c3266RbxsnUL4r3RS1MrXSeFpSSLcNya0T068niO8WShIFCGD0EIViRt4qluyeIrkr+lXv1fQcbjs4JT+dfGkZjQdTsH6p2y1F5RoDg40yb65CEMGoIQrAWvQb6sFy6X57XZPac3pcSmODl7ClHcRi56cFT6rob5NiViYgeihX/XmSUozLYNwQhWIvVJ6V5wWKkhxzT6Glr5c6FUVlyVIZRUQiKuREZP1TtkqP4LePEk510DEuPwighCMEqGCT6W5H0iDyvyTT16qo6T8/AvkvWYWFU1pbKnYz4x5VKpIdixdcL8coMjA6CEKzCZ5VStCdND5ClObilcte8yEyVQoadLGD0JvpEuyjVRc0n5Sj+q0ni1jqppgeNQhgFBCFYhdcLpccTFDIV3169e1HUXJmKw2VYEDVHpt5RDxWtnCj+rQiNQhgFBCFY3g/1zMhoYagszcEC3QlREGN8x8tRHC7Pwsg5u2tyBwwDchR/eLL4YbnUisn1MGIIQrC8P+YZn0iUZa8JIlpfvvXa6MXy1IbL5OfqGx8waac8a3AHuwo3RIl/LcTrozBSCEKwsG31TDdAPx8ry63YNdh9oOHIwrFz5CgOV2Lp+MXry7fIVHxVkviPEqlFlgYn2CEEIVjYi8eMz6eICnnag1sqd6SFTvN08pClOlyBaUFT2gc6ytoq5Cge7i5cHyW+eQKNQhgRBCFY0rZ6pu2XqzlIRBsrti1Fv6hVEgVhSfSC70/9IFP9PySJ75VIbRgphBFAEIIlvZAnY3PweFMhI4rzj5GlOlyxa8Yt3FmT3afvl6N4uLuwPFJ8HSOFMAIIQrCY72qkbj3dLFtz8JuyTcsnXC1TcbhyPi7eKZqErVU7Zar/XLL4folU14s5hXAJCEKwDCOjVUekV6cpZHpbVNurO9ZUuHjsPFmqAyc3xlz35ckNEpMlq0LchLsnii8ew5xCuAQEIVjGv0slfzUtkmfuIBF9eXL9knELXJRqmeoDF/H+kzyd3ffXH5ap/tNJim9rpOIONArhYhCEYAF9BvrjMem16XItJdOn799atWv5xCUy1QeOboi57vOT38lUfIwTPZ6gWHUYjUK4GAQhWMAr+cZZGmGqn1zNwe9PbZ0WlBLg6idTfeBoTnhaY4/2ZGu5TPUfihXz29jOBjQK4YIQhGBuld3sHyXSq9PkuvcMkvGr0u9vjFkqU33gSyEolk9Y8uXJDTLVVyvojRnir/cZ9WgWwgUgCMHcfrtfejxBEeImV3Pwh6pdIR5BWFzUhiwdv/hw47HarnqZ6i+LECM96O1iJCGcH4IQzGpTLSvrZL+VZ99BIpKY9EnxV3fG3yxTfZCDm8p1+cSrPy7+Sr5T/HWG4uXjRq0sUxbB5iEIwXx6DfTQPuNbaQon2e67bVW7fV18EgMmy3UCkMeNMdftqztU190gU/0JXsI9E8Xf7Mf8ejgPBCGYz9OHjXOChAUhcnWKSkz66MQXd8b9XKb6IB83leuyCVetK/5avlM8l6IoamNfV6ODFM6GIAQzOaBjX1Wzv8g2ZYKIQShQuAAAFolJREFUtlfv8XL2StEkyHcKkM8NMUv31u5v6NHKVN9ZQR/MUjy8X2rHAqTwUwhCMIcBI9291/jmDNHHWa5TDBmHPsj/+L7kO+Q6AcjM08njppjr/nn8Q/lOMTNAWB4hPHoQHaTwEwhCMIfHDxoTfYXro2S8374s3RDjGx3vHyvfKUBuN01aVtxSWtR8Ur5TvJyqyG1iX1ahgxT+B0EIsttcy9afZu+mydgp2j3U83nJd79MQnPQtjkrnH4Rf8u7ef9mJNf8d3cVrctSPLjPWNODKfbwIwQhyEvbTyuzDR/PUXjL1ilKRP8pWJcVkR7qESzjOcAsFo2dO2AYyKk9IN8ppvgJj8Qp7txjNCIKgYgQhCArg0Q37zT8apIiQyPXm6JEVN5eubMm+674W+Q7BZiNKIgPTbnnraP/GjAMyHeW3yeISoGeO4rBQiBCEIKsfnfQ6CzSqiQZbzOJsTcO/eO+pDu8nD3lOwuYU4omITEg9j8F6+Q7hSjQ5/OU6yowWAhECEKQzycV0uY69tk8pUwb0Jt8W7ZJKSoXj8O+g3bloSn3/FC9u7y9Ur5T+DjTV/MVD+4znmhHD6mjQxCCLHKb2CMHjF/PV4xxkvEszX0t/y389HfTHhBIzrAFs/Ny9lyZsOKNQ+9KTMYWW5Kv8OdUxbJtxmYZe2HBBiAIgb9TXeymHcYPZyvjvGXMJ4mxl/b/7YaYayO9wuQ7C1jKkugFbiq3tUVfyHqWuyaIt44Trt5i6DXIeh6waghC4KyxjxZuNr6UKsq3+7zJ5ye/1Rv1KybfIOtZwFIEElalPbK+fHNxS6msJ3p+iiLOR7hlp9GA4UJHhSAEnnT9NG+T4b4Y8c7x8t5ap9or1534+g/pvxMF3MN2y1s95rep9/1p3xv9cr5BKhD9M0PBiN22GxMqHBR+iQA3LQM0f7Ph5nHiE4ny3lddQ93PZv/54an3atwCZD0RWFxm2MwpmsSX9/9Nvin2RKQS6av5yl4DQ7vQMSEIgY/6Xpa10XBdhPBssrw3lcSkP+W+MSts5rzIWbKeCKzEb6be1zHQubboc1nP4iTSF/OUHUPstt3GIWShg0EQAgcVXWzORuP1UcKLU2RcR83kH8fWSEzCamqOQykqns/4/Ybyrbl1h2Q9kVpBGxYqGdHizYYuvaynAuuCIIQrtaeRZWwwPJMsPp8iewp+XvLtgfojz2U8jqFBh+Lj4v3irKdeO/hWYXOJrCdyVtAnWYqJY4SsjYbaXgwYOgr8NoEr8k6xdPNOw0dZyjtkfjuGiLZU7vyydMNf5j7v4eQu97nA2sT4jn82/fFn9r5U2nZK1hMpBPpHuuKWceKM74w5WmShQ0AQwmXq1tMdu43/PCnlLlXOC5Z9PvsPVbtX5699Y96LgXhBxlGlaBJ+N+2Bp3a/WNFRLfe5HosX/zNbccMOw+uFEsLQ7iEI4XIc1LHkbwwuStq/VDnWQ/YU/PLkhg/yP3pj7h+xv4SDywyb+Zup9/1ux7MFuhNyn2thiHBgqfKbamnxZkNjn9xnA0tCEMLo9BrosYPGZdsMr00T389QuCrlPZ3E2PvH/rv+1Ja3F74SgRVkgGh2eNqz6Y89m/3K3tr9cp8r0kPYvUSZFigmfaP/VymahnYLQQij8E21FP+VQddPhderlkfKfvN0DXY/sfuFktaytxe8EuDqJ/fpwFakaBJem/vC20f/9f6x/8q6GCkRKUV6LkXcdpXyg1Ipa6Mhvw1paIcQhDAiR9vEORuNz+dJH2QqPpyj8FPLfsbC5uJfbnl0rFfE6/P+6OnsIfv5wKaM9x67+qo3ytorfrfz2abeZrlPl+Aj5F6rvHmsuHiz4b59Qj16Su0LghAuYV8Tu3qr4fZc1e3RQt5y5Vz534sZMAz8/cjq53Nee3jqvfen/EIhyD4rA2yRl7Pna1nPT9Ek/XLzo+vLt8i69AwRiQL9apJYepMq0IWmrqf7c43V3Wgd2gkEIZzfkEQfn5JmrDfcsce4PEI8vmRw5URR1p0FiUhibEvlzts3PNCr712z5K20kGnyng9snCiIt02+4c0FL22u3PHA1scLm4vlPqOniv6YzAqXka8zTf3WcMMO4+5GxKHNk/lVB7BBh5vZ2lPSpxVSoq/wdKJ4TbgoCtTTI+9JJSbtOb3voxNfuCjVL2Q+Ees3Ud7zgR2J9Ap7d9Gr26v3vJj7RrR35IrJN0z2i5H1jH5q+tNUxZOJig/LpQdzjUZGt48Xb4sWItyxL6ZNGl0Q6vV6lUp15ceAtRmSKLeJfX9a+rqaqUS6LVo8eJ0ySv55EUTUPtC5rWrXN2WbfF18VibelhaSaoaTgp0RSFgQOWd2WNqmiu1/yn3Dz9X3uvGLM8NmOitk3BjaXUUPxIoPxIoHdOzDcmnqt8ZId+H6KPGqMCHBR0Ak2hCBsRG163t6eu68885t27YpFIonnnjiySefPPeYV1999eWXXzYajfPmzfvwww89PM7/gsMnn3yycePGjz/+eISX2N3dfaFScCUGjZTXyrK1LFsrZWtZzBjh6jBxeYQQ73OeR7inp8fNzU3g93R3Dnbtqz+cXbs/X3ciI3T60vGL5f4r3g709/c7OTkpFBg0vRiJSXtr9288te1ka3lm2IzMsJlTNAlO/BJxcHBQFMVz/9w3MtrTyL6tkbbWsW49mxcsZmqEDI0Q4yWISEXe9Hq90WhUq/m8tjfSIHz66aePHj26YcOGhoaG6dOnf/PNN2lpaWcecPDgwSVLlhw8eDAsLGz58uVxcXF//vOfz1sKQWgptb2srJOK2lhRO8trZSc72KQxQoZGmKUR5gSJPs4X+79cgrB9oPNka3m+ruh4U9HprrrUoOT00OkZodNdVS5XUtZxIAhHpbmvZffpfTl1B8vbKmL9JiYHxsf7TxrvM85FeUW/PS8UhGeq6ma7G9leLcttYk19LNlPSPIV4ryFOG9hgpdw8WcNRsIyQRgSErJmzZoFCxYQ0W9/+9u+vr5//vOfZx5w//33KxSKt99+m4h27dp18803NzU1nbcUglA+BolaB6llgDX1U0Mfa+ij+l5W1U3VPayyi3moKGaMMNlbiPcRknyERF9BPeLfqKMNwiHjUHNfq7ZX19ijreturOo4XdVZ06fvj/EdH+cfkxKYMMlvokrEEPXoIAgvT89Qb77uxLGmguKWsoqOao2bf6RXeIRXWKhHUJC7JsgtwNtlzMhfTh5JEJ6pY4iONLP8NlbUzk60s7JOphRorKcQ5SGEuVGYmxDsRkEugr8L+auRkSPFNwhH9Juor6+voaEhNjbW9GlsbOynn3561jGnTp1avnz58AE6na6rq8vT0/O8BYeGhtrb24c/HTNmzEV+w3b39bb32960HQOjXsOl/8joNdBZG4F268n0x4nEqNvAiMgoUY+BiKhLT8SoU88MEvXqqddAfUbqHWJdeurSU8cQ6zWQjxN5qwV/Z/J3EQJdaKyLkBlKoW5CmJvg/tMnd7CHBomIqE/fb7zArOQ+Q7+RGYlooH9AUAkGZiSinqE+RqzfMGBgxm5974BhoN8w2K3v7dH3dQ11dw52tw92DhiH/NRjAlz8gt0CQj2ClkbMifAIDT5zmdCBAWz6NlpsYEAyDAkIwlFyJZrpEzvTJ5aIjMxY091Q01VX3V1/4PQhbZ+uqa+lc6jb08ndy8nD08nDy9nDXeXqoXJzUarVSmdXpYuzQuUkOqlElVrpRESiJKgUKqVSSUROCidn8fydrs4KJyeFiogEolR3SnUnCv/xn1oHqLaH1fWyhn6qbWGHT1PzAGsbpOYB1megMSrydBY8lOThRB4qwVVJziJ5OZOCyF0lKAVyVREReTkREYlEHqoff3mqFOTy01tDrSDnEdwsHipb2sxFpVIGefvyrTmiIDSFlrv7j0v+e3p6trW1nXvM8AGmBlx7e/t5g7C8vHz9+vXbt28f/sq6desyMjIudPYX1j952rn9Qv/qCAQiuugcKff/397dxzR1tQEAf25bR4sgiLajkFVq/AC/EGFuGh2KiCGTaLTOd6QuApItm9kScB/JNmPCEv7QbTHO+LHAolm2VHFbVrdsGTq1zE1hQ0iaUYsl8llGoa1Cgfbee94/7nbTlyLW17aX9T6/v46nD/c+Sq9Pe+655wDEAaRIAbg3/Tgw42C/B3aAliCOr2ApyQMOr2BBSv4+fQyhZCzhOiUE5ISSsjCTgRgCcgZms9RMhopnII6GBIaKpwHACeAEsPJH6wvmb4tQmCkA0gH870izFOWWjdyXDt+T9o7IqBEpGZGAW0r6JdQYRbwS8FGElsA4BQAwJgH2n8/tXgq8D6ghXgnxBX0nQQqgBFDyl7APwAf06N+X0KQIiPTGo5SFo/nHkmbGMwxD0/RD42NjYyWShxT6oArh3LlzKYpyu90JCQkA4HQ6n3zyyQkxSqXS7XZzbZfLxfVMerSFCxfqdLrgh0aP/OcYDo0KLuSTZdD/AYdGwyf4dWwfdWgUhUNoh0aD+kIcExMzf/785uZm7o+3bt1KT584wS8jI8M/IC0tLTY2NiQpIoQQQuET7MhweXn5Bx980NXVde3atXPnzu3btw8ABgcHt2zZMjAwAABlZWUXLlz4+eefu7u7q6qqysvLw5g1QgghFCLBTturrKx0OBzr169PTEw8ceLE8uXLAYAQMj4+zs07Xbp06enTpysqKpxO586dO998880wZo0QQgiFSLCPT4TQoz4+8eGHH5aWls6ePTusWaGpnTp1qrCwUKPRPDwUhY3BYFiyZAn3MRQJ5Ycffpg5c+b69euFTkTUfv31V4fDUVRUFJKj/QsmzZ49e/bu3btCZyF2Fy5cMJvDvic4mtp3333X2NgodBZid+nSJZPJJHQWYvfLL7/89NNPoTrav6AQIoQQQuGDhRAhhJCoYSFECCEkagJMlqmurq6urn7Q4/aBenp6VCoVPr4qLLvdnpCQoFDg6thCGhgYUCgU/BJOSBBDQ0NSqZRbXQQJxe12MwyTlJT00Mji4uKqqqqpYwQohCzLWq3W4Avb+Ph4TAyuRCsw/C1MBz6fTyqVPnS9KBRWNE1TFIXr+wiLYRhCCLfi69TUavVDP8ELUAgRQgih6QM/WiKEEBI1LIQIIYREDQshQgghUcNCiBBCSNSCXXQ7MiwWy+eff07TdHFx8aQLKg4PD3/66addXV3r1q3bsWNH5DOMekNDQ0aj0Ww2x8fH79ixY+nSpYExNTU1DMNw7cWLF+fm5kY2x+hns9n8d67eunVrSkrKhBifz1dTU3P79u3MzMw9e/bgVNKQa21t/e233/x79Hr9hN3lzp07x22/CgBqtTpUS1+ijo6OpqYmp9O5e/du/ydVfv/9d4PBoFAoSkpK0tLSAn9wYGCgpqZmYGDg+eefz8vLC/J00+jiaW9vf+aZZwghcXFx69ata2mZuLk6IaSgoODKlSsLFix45513Dh8+LEie0e3AgQPffvutSqVyu92rV6+edDW/1157rbW11Waz2Ww2bhMuFFpNTU3V1dW2f4yOjgbG6PX6L7/8cuHChceOHXvjjTcin2TUc7lc/K/gm2++ef/99wMf+jp06JDJZOJi+vr6BMkz+vT392dnZ586derll1/+66+/+P7r16/n5eXNnTvX4/Hk5OT09PRM+EGPx7NmzRqLxTJv3rzi4mKDwRDsKcm08frrr+/bt49rv/XWWy+99NKEgPr6+tTUVK/XSwhpaGhQqVTcJlAohEZHR/n2gQMHdDpdYExMTExvb28EkxIdg8GQn58/RYDFYlEoFE6nkxBy9+5duVze398fqezESKfTVVZWBvZnZGQ0NDREPp/oxj0jSNM0ANy+fZvv37ZtW1VVFdfevXv3e++9N+EHa2trn376aZZlCSFffPHFihUrgjzjNPpGePXq1YKCAq69efPmq1evBgZs2LCB+1C2Zs0aj8fz559/RjrLaCeXy/n22NjYgxYx+eyzz44ePXrjxo1I5SU6fX19R44cqampsdvtga+aTKacnJzExEQA0Gg0Wq12wiAeCqHBwUGj0VhaWjrpq19//fVHH33kP5SNHtODxvmvXbu2efNmrj1pjeACKIriAlpbW4eGhoI642NkG2J9fX38umsqlcput5P/fdjfbrfzARKJRKlU9vb2RjpL0WhpaTlz5kxFRUXgS7m5uffv37darYWFhQcPHox8blEvPj4+MzPT5XJdvHgxIyOjqalpQoD/tQAAKpUKr4XwOXv2bFZW1pIlSwJfysrKAoDe3t69e/fq9fqIpyYiY2NjTqfTv0YEjkX7F5E5c+bIZLIgx6un0WSZGTNmcN+FAYCmaZlMxhV2nkwm4+doAIDP53viiScimqJodHZ2bt++/eOPP550ytKPP/7INcrKynJycvbv369SqSKbYJQrLCwsLCzk2pWVlQcPHvz+++/9A/BaiKQzZ87s379/0pf4DcYrKioWLVp08+bN1atXRzA1EeEWF/SvEYHveZlMxgcwDMMwTJDXxTT6Rpiamsp/qu3p6UlNTQ0M4O+Oer1eh8MROJUOPb7Ozs4NGza8/fbbZWVlU0dmZWXJ5XLcNjms1q5da7PZJnT6XwsA0NPTg9dCmNy8ebO9vf2FF16YOiwlJUWr1XZ0dEQmKxGaMWOGUqnk3/aTvuf9iwjXUKvVwRx8GhXCoqKi8+fPc+3z58/zE5FNJtPg4CAXcPnyZW7M12g0PvXUU+np6UJlG626u7vz8vJeffXVV155xb+/ubm5s7MTAMbGxvjOS5cuMQyzYMGCSGcZ7fh/ZELIxYsXly1bxv2xsbGR+4+goKDAbDbfuXOH63S73c8995xQ2Ua32traXbt2zZo1i+9pa2uzWCwA4PV6WZblO61W66SPG6FQKSoqqqurAwBCSF1dHVcjWJa9fPnyyMgIF2A0GrnLp66uLi8vL9itWh53fk/oOByORYsWbdmyZdu2bRqNpru7m+tPSkoyGo1cu6SkJD09fe/evUql8quvvhIu2ai1a9cuuVye/Q+9Xs/15+bmVldXE0IMBsPixYtffPHFrVu3xsXFnTx5UtB8o9P27ds3btyo1+tXrlyp1WqtVivXn5WVdfz4ca596NAhjUZTWlqanJz8ySefCJdsNPN4PImJiSaTyb+zpKSkvLycENLY2KjRaHQ63c6dO2fNmvXuu+8KlGYU2rRp06pVqwBg2bJl2dnZw8PDhJD29vbk5GSdTrdx48bMzMx79+4RQoaHhwGgtbWVEELTdEFBQXZ29p49e+bMmXP9+vUgTze9dp/weDz19fUMw+Tn58fHx3Odzc3NWq2WmyAHAA0NDV1dXc8++6xWqxUu06h1584d/gFhAIiNjc3IyACAtra2hIQEtVpN0/StW7esVmtcXFxOTk6QIw/okbhcrhs3bgwNDanV6rVr1/L3Ocxms1Kp5O/I/vHHHxaLZcWKFfhFJExGRkba2tpWrVrlP1+ho6ODoqi0tDSWZc1mc1tbm0wmy8zMnD9/voCpRpmWlhb+bh8ArFy5ktv3yuVy1dfXx8bGbtq0idsYjmXZpqam5cuXc3stMQxz5coVh8ORm5ubnJwc5OmmVyFECCGEImwa3SNECCGEIg8LIUIIIVHDQogQQkjUsBAihBASNSyECCGERA0LIUIIIVHDQogQQkjUsBAihBASNSyECCGERA0LIUIIIVHDQogQQkjU/gsCEhsiZg9QDgAAAABJRU5ErkJggg==",
      "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 2400 1600\">\n",
       "<defs>\n",
       "  <clipPath id=\"clip630\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"2400\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<path clip-path=\"url(#clip630)\" d=\"M0 1600 L2400 1600 L2400 0 L0 0  Z\" fill=\"#ffffff\" fill-rule=\"evenodd\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip631\">\n",
       "    <rect x=\"480\" y=\"0\" width=\"1681\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<path clip-path=\"url(#clip630)\" d=\"M156.598 1486.45 L2352.76 1486.45 L2352.76 47.2441 L156.598 47.2441  Z\" fill=\"#ffffff\" fill-rule=\"evenodd\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip632\">\n",
       "    <rect x=\"156\" y=\"47\" width=\"2197\" height=\"1440\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"218.754,1486.45 218.754,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"739.124,1486.45 739.124,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"1259.5,1486.45 1259.5,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"1779.87,1486.45 1779.87,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"2300.24,1486.45 2300.24,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"156.598,1445.72 2352.76,1445.72 \"/>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"156.598,1082.52 2352.76,1082.52 \"/>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"156.598,719.322 2352.76,719.322 \"/>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"156.598,356.126 2352.76,356.126 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,1486.45 2352.76,1486.45 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"218.754,1486.45 218.754,1467.55 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"739.124,1486.45 739.124,1467.55 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1259.5,1486.45 1259.5,1467.55 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1779.87,1486.45 1779.87,1467.55 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"2300.24,1486.45 2300.24,1467.55 \"/>\n",
       "<path clip-path=\"url(#clip630)\" d=\"M196.138 1517.37 Q192.527 1517.37 190.698 1520.93 Q188.893 1524.47 188.893 1531.6 Q188.893 1538.71 190.698 1542.27 Q192.527 1545.82 196.138 1545.82 Q199.772 1545.82 201.578 1542.27 Q203.407 1538.71 203.407 1531.6 Q203.407 1524.47 201.578 1520.93 Q199.772 1517.37 196.138 1517.37 M196.138 1513.66 Q201.948 1513.66 205.004 1518.27 Q208.082 1522.85 208.082 1531.6 Q208.082 1540.33 205.004 1544.94 Q201.948 1549.52 196.138 1549.52 Q190.328 1549.52 187.249 1544.94 Q184.194 1540.33 184.194 1531.6 Q184.194 1522.85 187.249 1518.27 Q190.328 1513.66 196.138 1513.66 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M216.3 1542.97 L221.184 1542.97 L221.184 1548.85 L216.3 1548.85 L216.3 1542.97 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M241.369 1517.37 Q237.758 1517.37 235.929 1520.93 Q234.124 1524.47 234.124 1531.6 Q234.124 1538.71 235.929 1542.27 Q237.758 1545.82 241.369 1545.82 Q245.003 1545.82 246.809 1542.27 Q248.638 1538.71 248.638 1531.6 Q248.638 1524.47 246.809 1520.93 Q245.003 1517.37 241.369 1517.37 M241.369 1513.66 Q247.179 1513.66 250.235 1518.27 Q253.314 1522.85 253.314 1531.6 Q253.314 1540.33 250.235 1544.94 Q247.179 1549.52 241.369 1549.52 Q235.559 1549.52 232.48 1544.94 Q229.425 1540.33 229.425 1531.6 Q229.425 1522.85 232.48 1518.27 Q235.559 1513.66 241.369 1513.66 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M710.861 1544.91 L727.18 1544.91 L727.18 1548.85 L705.236 1548.85 L705.236 1544.91 Q707.898 1542.16 712.481 1537.53 Q717.088 1532.88 718.268 1531.53 Q720.513 1529.01 721.393 1527.27 Q722.296 1525.51 722.296 1523.82 Q722.296 1521.07 720.351 1519.33 Q718.43 1517.6 715.328 1517.6 Q713.129 1517.6 710.676 1518.36 Q708.245 1519.13 705.467 1520.68 L705.467 1515.95 Q708.291 1514.82 710.745 1514.24 Q713.199 1513.66 715.236 1513.66 Q720.606 1513.66 723.8 1516.35 Q726.995 1519.03 726.995 1523.52 Q726.995 1525.65 726.185 1527.57 Q725.398 1529.47 723.291 1532.07 Q722.713 1532.74 719.611 1535.95 Q716.509 1539.15 710.861 1544.91 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M736.995 1542.97 L741.879 1542.97 L741.879 1548.85 L736.995 1548.85 L736.995 1542.97 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M752.111 1514.29 L770.467 1514.29 L770.467 1518.22 L756.393 1518.22 L756.393 1526.7 Q757.411 1526.35 758.43 1526.19 Q759.448 1526 760.467 1526 Q766.254 1526 769.634 1529.17 Q773.013 1532.34 773.013 1537.76 Q773.013 1543.34 769.541 1546.44 Q766.069 1549.52 759.749 1549.52 Q757.573 1549.52 755.305 1549.15 Q753.06 1548.78 750.652 1548.04 L750.652 1543.34 Q752.736 1544.47 754.958 1545.03 Q757.18 1545.58 759.657 1545.58 Q763.661 1545.58 765.999 1543.48 Q768.337 1541.37 768.337 1537.76 Q768.337 1534.15 765.999 1532.04 Q763.661 1529.94 759.657 1529.94 Q757.782 1529.94 755.907 1530.35 Q754.055 1530.77 752.111 1531.65 L752.111 1514.29 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M1226.66 1514.29 L1245.02 1514.29 L1245.02 1518.22 L1230.94 1518.22 L1230.94 1526.7 Q1231.96 1526.35 1232.98 1526.19 Q1234 1526 1235.02 1526 Q1240.8 1526 1244.18 1529.17 Q1247.56 1532.34 1247.56 1537.76 Q1247.56 1543.34 1244.09 1546.44 Q1240.62 1549.52 1234.3 1549.52 Q1232.12 1549.52 1229.85 1549.15 Q1227.61 1548.78 1225.2 1548.04 L1225.2 1543.34 Q1227.28 1544.47 1229.51 1545.03 Q1231.73 1545.58 1234.21 1545.58 Q1238.21 1545.58 1240.55 1543.48 Q1242.89 1541.37 1242.89 1537.76 Q1242.89 1534.15 1240.55 1532.04 Q1238.21 1529.94 1234.21 1529.94 Q1232.33 1529.94 1230.46 1530.35 Q1228.6 1530.77 1226.66 1531.65 L1226.66 1514.29 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M1256.78 1542.97 L1261.66 1542.97 L1261.66 1548.85 L1256.78 1548.85 L1256.78 1542.97 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M1281.84 1517.37 Q1278.23 1517.37 1276.4 1520.93 Q1274.6 1524.47 1274.6 1531.6 Q1274.6 1538.71 1276.4 1542.27 Q1278.23 1545.82 1281.84 1545.82 Q1285.48 1545.82 1287.28 1542.27 Q1289.11 1538.71 1289.11 1531.6 Q1289.11 1524.47 1287.28 1520.93 Q1285.48 1517.37 1281.84 1517.37 M1281.84 1513.66 Q1287.65 1513.66 1290.71 1518.27 Q1293.79 1522.85 1293.79 1531.6 Q1293.79 1540.33 1290.71 1544.94 Q1287.65 1549.52 1281.84 1549.52 Q1276.03 1549.52 1272.96 1544.94 Q1269.9 1540.33 1269.9 1531.6 Q1269.9 1522.85 1272.96 1518.27 Q1276.03 1513.66 1281.84 1513.66 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M1746.19 1514.29 L1768.41 1514.29 L1768.41 1516.28 L1755.86 1548.85 L1750.98 1548.85 L1762.78 1518.22 L1746.19 1518.22 L1746.19 1514.29 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M1777.53 1542.97 L1782.41 1542.97 L1782.41 1548.85 L1777.53 1548.85 L1777.53 1542.97 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M1792.64 1514.29 L1811 1514.29 L1811 1518.22 L1796.93 1518.22 L1796.93 1526.7 Q1797.94 1526.35 1798.96 1526.19 Q1799.98 1526 1801 1526 Q1806.79 1526 1810.17 1529.17 Q1813.55 1532.34 1813.55 1537.76 Q1813.55 1543.34 1810.07 1546.44 Q1806.6 1549.52 1800.28 1549.52 Q1798.11 1549.52 1795.84 1549.15 Q1793.59 1548.78 1791.19 1548.04 L1791.19 1543.34 Q1793.27 1544.47 1795.49 1545.03 Q1797.71 1545.58 1800.19 1545.58 Q1804.19 1545.58 1806.53 1543.48 Q1808.87 1541.37 1808.87 1537.76 Q1808.87 1534.15 1806.53 1532.04 Q1804.19 1529.94 1800.19 1529.94 Q1798.32 1529.94 1796.44 1530.35 Q1794.59 1530.77 1792.64 1531.65 L1792.64 1514.29 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2252.31 1544.91 L2259.95 1544.91 L2259.95 1518.55 L2251.64 1520.21 L2251.64 1515.95 L2259.9 1514.29 L2264.58 1514.29 L2264.58 1544.91 L2272.22 1544.91 L2272.22 1548.85 L2252.31 1548.85 L2252.31 1544.91 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2291.66 1517.37 Q2288.05 1517.37 2286.22 1520.93 Q2284.42 1524.47 2284.42 1531.6 Q2284.42 1538.71 2286.22 1542.27 Q2288.05 1545.82 2291.66 1545.82 Q2295.29 1545.82 2297.1 1542.27 Q2298.93 1538.71 2298.93 1531.6 Q2298.93 1524.47 2297.1 1520.93 Q2295.29 1517.37 2291.66 1517.37 M2291.66 1513.66 Q2297.47 1513.66 2300.53 1518.27 Q2303.61 1522.85 2303.61 1531.6 Q2303.61 1540.33 2300.53 1544.94 Q2297.47 1549.52 2291.66 1549.52 Q2285.85 1549.52 2282.77 1544.94 Q2279.72 1540.33 2279.72 1531.6 Q2279.72 1522.85 2282.77 1518.27 Q2285.85 1513.66 2291.66 1513.66 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2311.82 1542.97 L2316.71 1542.97 L2316.71 1548.85 L2311.82 1548.85 L2311.82 1542.97 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2336.89 1517.37 Q2333.28 1517.37 2331.45 1520.93 Q2329.65 1524.47 2329.65 1531.6 Q2329.65 1538.71 2331.45 1542.27 Q2333.28 1545.82 2336.89 1545.82 Q2340.53 1545.82 2342.33 1542.27 Q2344.16 1538.71 2344.16 1531.6 Q2344.16 1524.47 2342.33 1520.93 Q2340.53 1517.37 2336.89 1517.37 M2336.89 1513.66 Q2342.7 1513.66 2345.76 1518.27 Q2348.84 1522.85 2348.84 1531.6 Q2348.84 1540.33 2345.76 1544.94 Q2342.7 1549.52 2336.89 1549.52 Q2331.08 1549.52 2328 1544.94 Q2324.95 1540.33 2324.95 1531.6 Q2324.95 1522.85 2328 1518.27 Q2331.08 1513.66 2336.89 1513.66 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,1486.45 156.598,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,1445.72 175.496,1445.72 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,1082.52 175.496,1082.52 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,719.322 175.496,719.322 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,356.126 175.496,356.126 \"/>\n",
       "<path clip-path=\"url(#clip630)\" d=\"M63.4226 1431.51 Q59.8115 1431.51 57.9828 1435.08 Q56.1773 1438.62 56.1773 1445.75 Q56.1773 1452.86 57.9828 1456.42 Q59.8115 1459.96 63.4226 1459.96 Q67.0569 1459.96 68.8624 1456.42 Q70.6911 1452.86 70.6911 1445.75 Q70.6911 1438.62 68.8624 1435.08 Q67.0569 1431.51 63.4226 1431.51 M63.4226 1427.81 Q69.2328 1427.81 72.2883 1432.42 Q75.367 1437 75.367 1445.75 Q75.367 1454.48 72.2883 1459.08 Q69.2328 1463.67 63.4226 1463.67 Q57.6125 1463.67 54.5338 1459.08 Q51.4782 1454.48 51.4782 1445.75 Q51.4782 1437 54.5338 1432.42 Q57.6125 1427.81 63.4226 1427.81 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M83.5845 1457.12 L88.4688 1457.12 L88.4688 1463 L83.5845 1463 L83.5845 1457.12 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M108.654 1431.51 Q105.043 1431.51 103.214 1435.08 Q101.409 1438.62 101.409 1445.75 Q101.409 1452.86 103.214 1456.42 Q105.043 1459.96 108.654 1459.96 Q112.288 1459.96 114.094 1456.42 Q115.922 1452.86 115.922 1445.75 Q115.922 1438.62 114.094 1435.08 Q112.288 1431.51 108.654 1431.51 M108.654 1427.81 Q114.464 1427.81 117.52 1432.42 Q120.598 1437 120.598 1445.75 Q120.598 1454.48 117.52 1459.08 Q114.464 1463.67 108.654 1463.67 Q102.844 1463.67 99.765 1459.08 Q96.7095 1454.48 96.7095 1445.75 Q96.7095 1437 99.765 1432.42 Q102.844 1427.81 108.654 1427.81 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M65.0198 1068.32 Q61.4087 1068.32 59.58 1071.88 Q57.7745 1075.42 57.7745 1082.55 Q57.7745 1089.66 59.58 1093.22 Q61.4087 1096.77 65.0198 1096.77 Q68.6541 1096.77 70.4596 1093.22 Q72.2883 1089.66 72.2883 1082.55 Q72.2883 1075.42 70.4596 1071.88 Q68.6541 1068.32 65.0198 1068.32 M65.0198 1064.61 Q70.83 1064.61 73.8855 1069.22 Q76.9642 1073.8 76.9642 1082.55 Q76.9642 1091.28 73.8855 1095.89 Q70.83 1100.47 65.0198 1100.47 Q59.2097 1100.47 56.131 1095.89 Q53.0754 1091.28 53.0754 1082.55 Q53.0754 1073.8 56.131 1069.22 Q59.2097 1064.61 65.0198 1064.61 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M85.1818 1093.92 L90.066 1093.92 L90.066 1099.8 L85.1818 1099.8 L85.1818 1093.92 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M104.279 1095.86 L120.598 1095.86 L120.598 1099.8 L98.6539 1099.8 L98.6539 1095.86 Q101.316 1093.11 105.899 1088.48 Q110.506 1083.83 111.686 1082.48 Q113.932 1079.96 114.811 1078.23 Q115.714 1076.47 115.714 1074.78 Q115.714 1072.02 113.77 1070.29 Q111.848 1068.55 108.746 1068.55 Q106.547 1068.55 104.094 1069.31 Q101.663 1070.08 98.8854 1071.63 L98.8854 1066.91 Q101.709 1065.77 104.163 1065.19 Q106.617 1064.61 108.654 1064.61 Q114.024 1064.61 117.219 1067.3 Q120.413 1069.98 120.413 1074.48 Q120.413 1076.6 119.603 1078.53 Q118.816 1080.42 116.709 1083.02 Q116.131 1083.69 113.029 1086.91 Q109.927 1090.1 104.279 1095.86 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M62.9365 705.121 Q59.3254 705.121 57.4967 708.686 Q55.6912 712.227 55.6912 719.357 Q55.6912 726.463 57.4967 730.028 Q59.3254 733.57 62.9365 733.57 Q66.5707 733.57 68.3763 730.028 Q70.205 726.463 70.205 719.357 Q70.205 712.227 68.3763 708.686 Q66.5707 705.121 62.9365 705.121 M62.9365 701.417 Q68.7467 701.417 71.8022 706.024 Q74.8809 710.607 74.8809 719.357 Q74.8809 728.084 71.8022 732.69 Q68.7467 737.274 62.9365 737.274 Q57.1264 737.274 54.0477 732.69 Q50.9921 728.084 50.9921 719.357 Q50.9921 710.607 54.0477 706.024 Q57.1264 701.417 62.9365 701.417 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M83.0984 730.723 L87.9827 730.723 L87.9827 736.602 L83.0984 736.602 L83.0984 730.723 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M111.015 706.116 L99.2095 724.565 L111.015 724.565 L111.015 706.116 M109.788 702.042 L115.668 702.042 L115.668 724.565 L120.598 724.565 L120.598 728.454 L115.668 728.454 L115.668 736.602 L111.015 736.602 L111.015 728.454 L95.4132 728.454 L95.4132 723.94 L109.788 702.042 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M63.2606 341.924 Q59.6495 341.924 57.8208 345.489 Q56.0152 349.031 56.0152 356.16 Q56.0152 363.267 57.8208 366.831 Q59.6495 370.373 63.2606 370.373 Q66.8948 370.373 68.7004 366.831 Q70.5291 363.267 70.5291 356.16 Q70.5291 349.031 68.7004 345.489 Q66.8948 341.924 63.2606 341.924 M63.2606 338.221 Q69.0707 338.221 72.1263 342.827 Q75.205 347.41 75.205 356.16 Q75.205 364.887 72.1263 369.494 Q69.0707 374.077 63.2606 374.077 Q57.4504 374.077 54.3717 369.494 Q51.3162 364.887 51.3162 356.16 Q51.3162 347.41 54.3717 342.827 Q57.4504 338.221 63.2606 338.221 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M83.4225 367.526 L88.3067 367.526 L88.3067 373.406 L83.4225 373.406 L83.4225 367.526 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M109.071 354.262 Q105.922 354.262 104.071 356.415 Q102.242 358.568 102.242 362.318 Q102.242 366.044 104.071 368.22 Q105.922 370.373 109.071 370.373 Q112.219 370.373 114.047 368.22 Q115.899 366.044 115.899 362.318 Q115.899 358.568 114.047 356.415 Q112.219 354.262 109.071 354.262 M118.353 339.609 L118.353 343.869 Q116.594 343.035 114.788 342.596 Q113.006 342.156 111.246 342.156 Q106.617 342.156 104.163 345.281 Q101.733 348.406 101.385 354.725 Q102.751 352.711 104.811 351.646 Q106.871 350.558 109.348 350.558 Q114.557 350.558 117.566 353.73 Q120.598 356.878 120.598 362.318 Q120.598 367.642 117.45 370.859 Q114.302 374.077 109.071 374.077 Q103.075 374.077 99.9039 369.494 Q96.7326 364.887 96.7326 356.16 Q96.7326 347.966 100.621 343.105 Q104.51 338.221 111.061 338.221 Q112.82 338.221 114.603 338.568 Q116.408 338.915 118.353 339.609 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><polyline clip-path=\"url(#clip632)\" style=\"stroke:#009af9; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"218.754,1445.72 228.39,1445.72 238.027,1445.72 247.663,1445.72 257.3,1445.72 266.936,1445.72 276.573,1445.72 286.209,1445.71 295.846,1445.71 305.482,1445.71 315.119,1445.71 324.755,1445.71 334.392,1445.71 344.028,1445.71 353.665,1445.71 363.301,1445.71 372.938,1445.71 382.574,1445.7 392.211,1445.7 401.847,1445.69 411.484,1445.69 421.12,1445.68 430.757,1445.67 440.393,1445.65 450.03,1445.64 459.666,1445.61 469.303,1445.58 478.939,1445.55 488.576,1445.51 498.212,1445.45 507.849,1445.38 517.485,1445.3 527.122,1445.19 536.758,1445.06 546.395,1444.9 556.031,1444.7 565.668,1444.47 575.304,1444.18 584.941,1443.82 594.577,1443.4 604.214,1442.89 613.85,1442.28 623.487,1441.55 633.123,1440.68 642.76,1439.65 652.396,1438.43 662.033,1436.99 671.669,1435.3 681.305,1433.32 690.942,1431.02 700.578,1428.34 710.215,1425.24 719.851,1421.67 729.488,1417.56 739.124,1412.86 748.761,1407.49 758.397,1401.39 768.034,1394.48 777.67,1386.69 787.307,1377.93 796.943,1368.13 806.58,1357.2 816.216,1345.06 825.853,1331.63 835.489,1316.83 845.126,1300.58 854.762,1282.83 864.399,1263.51 874.035,1242.57 883.672,1219.98 893.308,1195.69 902.945,1169.72 912.581,1142.05 922.218,1112.71 931.854,1081.73 941.491,1049.19 951.127,1015.15 960.764,979.714 970.4,942.996 980.037,905.13 989.673,866.267 999.31,826.574 1008.95,786.229 1018.58,745.422 1028.22,704.353 1037.86,663.223 1047.49,622.241 1057.13,581.614 1066.77,541.546 1076.4,502.237 1086.04,463.879 1095.67,426.658 1105.31,390.745 1114.95,356.303 1124.58,323.481 1134.22,292.415 1143.86,263.228 1153.49,236.029 1163.13,210.915 1172.77,187.971 1182.4,167.269 1192.04,148.873 1201.68,132.833 1211.31,119.194 1220.95,107.99 1230.59,99.2491 1240.22,92.9911 1249.86,89.2307 1259.5,87.9763 1269.13,89.2307 1278.77,92.9911 1288.4,99.2491 1298.04,107.99 1307.68,119.194 1317.31,132.833 1326.95,148.873 1336.59,167.269 1346.22,187.971 1355.86,210.915 1365.5,236.029 1375.13,263.228 1384.77,292.415 1394.41,323.481 1404.04,356.303 1413.68,390.745 1423.32,426.658 1432.95,463.879 1442.59,502.237 1452.23,541.546 1461.86,581.614 1471.5,622.241 1481.13,663.223 1490.77,704.353 1500.41,745.422 1510.04,786.229 1519.68,826.574 1529.32,866.267 1538.95,905.13 1548.59,942.996 1558.23,979.714 1567.86,1015.15 1577.5,1049.19 1587.14,1081.73 1596.77,1112.71 1606.41,1142.05 1616.05,1169.72 1625.68,1195.69 1635.32,1219.98 1644.96,1242.57 1654.59,1263.51 1664.23,1282.83 1673.86,1300.58 1683.5,1316.83 1693.14,1331.63 1702.77,1345.06 1712.41,1357.2 1722.05,1368.13 1731.68,1377.93 1741.32,1386.69 1750.96,1394.48 1760.59,1401.39 1770.23,1407.49 1779.87,1412.86 1789.5,1417.56 1799.14,1421.67 1808.78,1425.24 1818.41,1428.34 1828.05,1431.02 1837.69,1433.32 1847.32,1435.3 1856.96,1436.99 1866.59,1438.43 1876.23,1439.65 1885.87,1440.68 1895.5,1441.55 1905.14,1442.28 1914.78,1442.89 1924.41,1443.4 1934.05,1443.82 1943.69,1444.18 1953.32,1444.47 1962.96,1444.7 1972.6,1444.9 1982.23,1445.06 1991.87,1445.19 2001.51,1445.3 2011.14,1445.38 2020.78,1445.45 2030.42,1445.51 2040.05,1445.55 2049.69,1445.58 2059.32,1445.61 2068.96,1445.64 2078.6,1445.65 2088.23,1445.67 2097.87,1445.68 2107.51,1445.69 2117.14,1445.69 2126.78,1445.7 2136.42,1445.7 2146.05,1445.71 2155.69,1445.71 2165.33,1445.71 2174.96,1445.71 2184.6,1445.71 2194.24,1445.71 2203.87,1445.71 2213.51,1445.71 2223.15,1445.71 2232.78,1445.71 2242.42,1445.72 2252.05,1445.72 2261.69,1445.72 2271.33,1445.72 2280.96,1445.72 2290.6,1445.72 \"/>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#e26f46; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"218.754,1445.72 228.39,1445.72 238.027,1445.72 247.663,1445.72 257.3,1445.72 266.936,1445.72 276.573,1445.72 286.209,1445.72 295.846,1445.72 305.482,1445.72 315.119,1445.72 324.755,1445.72 334.392,1445.72 344.028,1445.72 353.665,1445.72 363.301,1445.72 372.938,1445.72 382.574,1445.72 392.211,1445.72 401.847,1445.72 411.484,1445.72 421.12,1445.72 430.757,1445.72 440.393,1445.72 450.03,1445.72 459.666,1445.72 469.303,1445.72 478.939,1445.72 488.576,1445.72 498.212,1445.72 507.849,1445.72 517.485,1445.72 527.122,1445.72 536.758,1445.72 546.395,1445.72 556.031,1445.72 565.668,1445.72 575.304,1445.72 584.941,1445.72 594.577,1445.72 604.214,1445.72 613.85,1445.72 623.487,1445.72 633.123,1445.72 642.76,1445.72 652.396,1445.72 662.033,1445.72 671.669,1445.72 681.305,1445.72 690.942,1445.72 700.578,1445.72 710.215,1445.72 719.851,1445.72 729.488,1445.72 739.124,1445.72 748.761,1445.72 758.397,1445.72 768.034,1445.72 777.67,1445.72 787.307,1445.72 796.943,1445.72 806.58,1445.72 816.216,1445.72 825.853,1445.72 835.489,1445.72 845.126,1445.72 854.762,1445.72 864.399,1445.72 874.035,1445.72 883.672,1445.72 893.308,1445.72 902.945,1445.72 912.581,1445.72 922.218,1445.72 931.854,1445.72 941.491,1445.72 951.127,1445.72 960.764,1445.72 970.4,1445.72 980.037,1445.72 989.673,1445.72 999.31,1445.72 1008.95,1445.72 1018.58,1445.72 1028.22,1445.72 1037.86,1445.72 1047.49,1445.72 1057.13,1445.72 1066.77,1445.72 1076.4,1445.72 1086.04,1445.72 1095.67,1445.72 1105.31,1445.72 1114.95,1445.72 1124.58,1445.72 1134.22,1445.72 1143.86,1445.72 1153.49,1445.72 1163.13,1445.72 1172.77,1445.72 1182.4,1445.72 1192.04,1445.72 1201.68,1445.72 1211.31,1445.72 1220.95,1445.72 1230.59,1445.72 1240.22,1445.72 1249.86,1445.72 1259.5,1445.72 1269.13,1445.72 1278.77,1445.72 1288.4,1445.72 1298.04,1445.72 1307.68,1445.72 1317.31,1445.72 1326.95,1445.72 1336.59,1445.72 1346.22,1445.72 1355.86,1445.72 1365.5,1445.72 1375.13,1445.72 1384.77,1445.72 1394.41,1445.72 1404.04,1445.72 1413.68,1445.72 1423.32,1445.72 1432.95,1445.72 1442.59,1445.72 1452.23,1445.72 1461.86,1445.72 1471.5,1445.72 1481.13,1445.72 1490.77,1445.72 1500.41,1445.72 1510.04,1445.72 1519.68,1445.72 1529.32,1445.72 1538.95,1445.72 1548.59,1445.72 1558.23,1445.72 1567.86,1445.72 1577.5,1445.72 1587.14,1445.72 1596.77,1445.72 1606.41,1445.72 1616.05,1445.72 1625.68,1445.72 1635.32,1445.72 1644.96,1445.72 1654.59,1445.72 1664.23,1445.72 1673.86,1445.72 1683.5,1445.72 1693.14,1445.72 1702.77,1445.72 1712.41,1445.72 1722.05,1445.72 1731.68,1445.72 1741.32,1445.72 1750.96,1445.72 1760.59,1445.72 1770.23,1445.72 1779.87,1445.72 1789.5,1445.72 1799.14,1445.72 1808.78,1445.72 1818.41,1445.72 1828.05,1445.72 1837.69,1445.72 1847.32,1445.72 1856.96,1445.72 1866.59,1445.72 1876.23,1445.72 1885.87,1445.72 1895.5,1445.72 1905.14,1445.72 1914.78,1445.72 1924.41,1445.72 1934.05,1445.72 1943.69,1445.72 1953.32,1445.72 1962.96,1445.72 1972.6,1445.72 1982.23,1445.72 1991.87,1445.72 2001.51,1445.72 2011.14,1445.72 2020.78,1445.72 2030.42,1445.72 2040.05,1445.72 2049.69,1445.72 2059.32,1445.72 2068.96,1445.72 2078.6,1445.72 2088.23,1445.72 2097.87,1445.72 2107.51,1445.72 2117.14,1445.72 2126.78,1445.72 2136.42,1445.72 2146.05,1445.72 2155.69,1445.72 2165.33,1445.72 2174.96,1445.72 2184.6,1445.72 2194.24,1445.72 2203.87,1445.72 2213.51,1445.72 2223.15,1445.72 2232.78,1445.72 2242.42,1445.72 2252.05,1445.72 2261.69,1445.72 2271.33,1445.72 2280.96,1445.72 2290.6,1445.72 \"/>\n",
       "<polyline clip-path=\"url(#clip632)\" style=\"stroke:#3da44d; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"218.754,1445.72 228.39,1445.72 238.027,1445.72 247.663,1445.72 257.3,1445.72 266.936,1445.72 276.573,1445.72 286.209,1445.72 295.846,1445.72 305.482,1445.72 315.119,1445.72 324.755,1445.72 334.392,1445.72 344.028,1445.72 353.665,1445.72 363.301,1445.72 372.938,1445.72 382.574,1445.72 392.211,1445.72 401.847,1445.72 411.484,1445.72 421.12,1445.72 430.757,1445.72 440.393,1445.72 450.03,1445.72 459.666,1445.72 469.303,1445.72 478.939,1445.72 488.576,1445.72 498.212,1445.72 507.849,1445.72 517.485,1445.72 527.122,1445.72 536.758,1445.72 546.395,1445.72 556.031,1445.72 565.668,1445.71 575.304,1445.71 584.941,1445.71 594.577,1445.71 604.214,1445.71 613.85,1445.71 623.487,1445.71 633.123,1445.7 642.76,1445.7 652.396,1445.69 662.033,1445.67 671.669,1445.66 681.305,1445.63 690.942,1445.6 700.578,1445.55 710.215,1445.48 719.851,1445.4 729.488,1445.28 739.124,1445.12 748.761,1444.91 758.397,1444.63 768.034,1444.27 777.67,1443.8 787.307,1443.19 796.943,1442.4 806.58,1441.4 816.216,1440.14 825.853,1438.55 835.489,1436.57 845.126,1434.12 854.762,1431.11 864.399,1427.43 874.035,1422.99 883.672,1417.65 893.308,1411.29 902.945,1403.77 912.581,1394.94 922.218,1384.65 931.854,1372.76 941.491,1359.13 951.127,1343.63 960.764,1326.13 970.4,1306.55 980.037,1284.79 989.673,1260.82 999.31,1234.63 1008.95,1206.22 1018.58,1175.66 1028.22,1143.06 1037.86,1108.55 1047.49,1072.3 1057.13,1034.55 1066.77,995.534 1076.4,955.539 1086.04,914.872 1095.67,873.861 1105.31,832.845 1114.95,792.175 1124.58,752.202 1134.22,713.274 1143.86,675.732 1153.49,639.903 1163.13,606.097 1172.77,574.605 1182.4,545.694 1192.04,519.605 1201.68,496.555 1211.31,476.732 1220.95,460.294 1230.59,447.374 1240.22,438.072 1249.86,432.462 1259.5,430.587 1269.13,432.462 1278.77,438.072 1288.4,447.374 1298.04,460.294 1307.68,476.732 1317.31,496.555 1326.95,519.605 1336.59,545.694 1346.22,574.605 1355.86,606.097 1365.5,639.903 1375.13,675.732 1384.77,713.274 1394.41,752.202 1404.04,792.175 1413.68,832.845 1423.32,873.861 1432.95,914.872 1442.59,955.539 1452.23,995.534 1461.86,1034.55 1471.5,1072.3 1481.13,1108.55 1490.77,1143.06 1500.41,1175.66 1510.04,1206.22 1519.68,1234.63 1529.32,1260.82 1538.95,1284.79 1548.59,1306.55 1558.23,1326.13 1567.86,1343.63 1577.5,1359.13 1587.14,1372.76 1596.77,1384.65 1606.41,1394.94 1616.05,1403.77 1625.68,1411.29 1635.32,1417.65 1644.96,1422.99 1654.59,1427.43 1664.23,1431.11 1673.86,1434.12 1683.5,1436.57 1693.14,1438.55 1702.77,1440.14 1712.41,1441.4 1722.05,1442.4 1731.68,1443.19 1741.32,1443.8 1750.96,1444.27 1760.59,1444.63 1770.23,1444.91 1779.87,1445.12 1789.5,1445.28 1799.14,1445.4 1808.78,1445.48 1818.41,1445.55 1828.05,1445.6 1837.69,1445.63 1847.32,1445.66 1856.96,1445.67 1866.59,1445.69 1876.23,1445.7 1885.87,1445.7 1895.5,1445.71 1905.14,1445.71 1914.78,1445.71 1924.41,1445.71 1934.05,1445.71 1943.69,1445.71 1953.32,1445.71 1962.96,1445.72 1972.6,1445.72 1982.23,1445.72 1991.87,1445.72 2001.51,1445.72 2011.14,1445.72 2020.78,1445.72 2030.42,1445.72 2040.05,1445.72 2049.69,1445.72 2059.32,1445.72 2068.96,1445.72 2078.6,1445.72 2088.23,1445.72 2097.87,1445.72 2107.51,1445.72 2117.14,1445.72 2126.78,1445.72 2136.42,1445.72 2146.05,1445.72 2155.69,1445.72 2165.33,1445.72 2174.96,1445.72 2184.6,1445.72 2194.24,1445.72 2203.87,1445.72 2213.51,1445.72 2223.15,1445.72 2232.78,1445.72 2242.42,1445.72 2252.05,1445.72 2261.69,1445.72 2271.33,1445.72 2280.96,1445.72 2290.6,1445.72 \"/>\n",
       "<path clip-path=\"url(#clip630)\" d=\"M1881.72 302.578 L2279.55 302.578 L2279.55 95.2176 L1881.72 95.2176  Z\" fill=\"#ffffff\" fill-rule=\"evenodd\" fill-opacity=\"1\"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1881.72,302.578 2279.55,302.578 2279.55,95.2176 1881.72,95.2176 1881.72,302.578 \"/>\n",
       "<polyline clip-path=\"url(#clip630)\" style=\"stroke:#009af9; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1906.12,147.058 2052.53,147.058 \"/>\n",
       "<path clip-path=\"url(#clip630)\" d=\"M2092.12 142.393 Q2091.4 141.977 2090.54 141.791 Q2089.71 141.583 2088.69 141.583 Q2085.08 141.583 2083.14 143.944 Q2081.21 146.282 2081.21 150.68 L2081.21 164.338 L2076.93 164.338 L2076.93 138.412 L2081.21 138.412 L2081.21 142.44 Q2082.56 140.078 2084.71 138.944 Q2086.86 137.787 2089.94 137.787 Q2090.38 137.787 2090.91 137.856 Q2091.45 137.903 2092.09 138.018 L2092.12 142.393 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2117.72 150.31 L2117.72 152.393 L2098.14 152.393 Q2098.41 156.791 2100.77 159.106 Q2103.16 161.398 2107.4 161.398 Q2109.85 161.398 2112.14 160.796 Q2114.46 160.194 2116.72 158.99 L2116.72 163.018 Q2114.43 163.99 2112.02 164.5 Q2109.62 165.009 2107.14 165.009 Q2100.94 165.009 2097.3 161.398 Q2093.69 157.787 2093.69 151.629 Q2093.69 145.264 2097.12 141.537 Q2100.57 137.787 2106.4 137.787 Q2111.63 137.787 2114.66 141.166 Q2117.72 144.523 2117.72 150.31 M2113.46 149.06 Q2113.41 145.565 2111.49 143.481 Q2109.59 141.398 2106.45 141.398 Q2102.88 141.398 2100.73 143.412 Q2098.6 145.426 2098.27 149.083 L2113.46 149.06 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2136.49 151.305 Q2131.33 151.305 2129.34 152.486 Q2127.35 153.666 2127.35 156.514 Q2127.35 158.782 2128.83 160.125 Q2130.33 161.444 2132.9 161.444 Q2136.45 161.444 2138.58 158.944 Q2140.73 156.421 2140.73 152.254 L2140.73 151.305 L2136.49 151.305 M2144.99 149.546 L2144.99 164.338 L2140.73 164.338 L2140.73 160.402 Q2139.27 162.763 2137.09 163.898 Q2134.92 165.009 2131.77 165.009 Q2127.79 165.009 2125.43 162.787 Q2123.09 160.541 2123.09 156.791 Q2123.09 152.416 2126.01 150.194 Q2128.95 147.972 2134.76 147.972 L2140.73 147.972 L2140.73 147.555 Q2140.73 144.615 2138.78 143.018 Q2136.86 141.398 2133.37 141.398 Q2131.14 141.398 2129.04 141.93 Q2126.93 142.463 2124.99 143.527 L2124.99 139.592 Q2127.33 138.69 2129.52 138.25 Q2131.72 137.787 2133.81 137.787 Q2139.43 137.787 2142.21 140.703 Q2144.99 143.62 2144.99 149.546 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2153.76 128.319 L2158.02 128.319 L2158.02 164.338 L2153.76 164.338 L2153.76 128.319 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2177.16 128.366 Q2174.06 133.69 2172.56 138.898 Q2171.05 144.106 2171.05 149.453 Q2171.05 154.801 2172.56 160.055 Q2174.08 165.287 2177.16 170.588 L2173.46 170.588 Q2169.99 165.148 2168.25 159.893 Q2166.54 154.639 2166.54 149.453 Q2166.54 144.291 2168.25 139.06 Q2169.96 133.828 2173.46 128.366 L2177.16 128.366 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2194.41 164.916 Q2189.06 164.06 2186.58 161.814 Q2183.55 159.06 2183.55 153.435 L2183.55 138.412 L2187.86 138.412 L2187.86 153.273 Q2187.86 157.509 2189.83 159.268 Q2191.54 160.796 2194.41 161.12 L2194.41 138.412 L2198.64 138.412 L2198.64 161.097 Q2201.68 160.773 2203.23 159.245 Q2205.2 157.301 2205.2 153.25 L2205.2 138.412 L2209.5 138.412 L2209.5 153.412 Q2209.5 159.245 2206.47 161.791 Q2203.74 164.083 2198.64 164.893 L2198.64 174.199 L2194.41 174.199 L2194.41 164.916 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2216.03 128.366 L2219.73 128.366 Q2223.2 133.828 2224.92 139.06 Q2226.65 144.291 2226.65 149.453 Q2226.65 154.639 2224.92 159.893 Q2223.2 165.148 2219.73 170.588 L2216.03 170.588 Q2219.11 165.287 2220.61 160.055 Q2222.14 154.801 2222.14 149.453 Q2222.14 144.106 2220.61 138.898 Q2219.11 133.69 2216.03 128.366 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><polyline clip-path=\"url(#clip630)\" style=\"stroke:#e26f46; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1906.12,198.898 2052.53,198.898 \"/>\n",
       "<path clip-path=\"url(#clip630)\" d=\"M2076.93 190.252 L2081.19 190.252 L2081.19 216.178 L2076.93 216.178 L2076.93 190.252 M2076.93 180.159 L2081.19 180.159 L2081.19 185.553 L2076.93 185.553 L2076.93 180.159 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2110.29 195.229 Q2111.89 192.358 2114.11 190.993 Q2116.33 189.627 2119.34 189.627 Q2123.39 189.627 2125.59 192.474 Q2127.79 195.298 2127.79 200.529 L2127.79 216.178 L2123.51 216.178 L2123.51 200.668 Q2123.51 196.942 2122.19 195.136 Q2120.87 193.33 2118.16 193.33 Q2114.85 193.33 2112.93 195.53 Q2111.01 197.729 2111.01 201.525 L2111.01 216.178 L2106.72 216.178 L2106.72 200.668 Q2106.72 196.918 2105.4 195.136 Q2104.08 193.33 2101.33 193.33 Q2098.07 193.33 2096.15 195.553 Q2094.22 197.752 2094.22 201.525 L2094.22 216.178 L2089.94 216.178 L2089.94 190.252 L2094.22 190.252 L2094.22 194.28 Q2095.68 191.895 2097.72 190.761 Q2099.76 189.627 2102.56 189.627 Q2105.38 189.627 2107.35 191.062 Q2109.34 192.497 2110.29 195.229 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2148.07 203.145 Q2142.9 203.145 2140.91 204.326 Q2138.92 205.506 2138.92 208.354 Q2138.92 210.622 2140.4 211.965 Q2141.91 213.284 2144.48 213.284 Q2148.02 213.284 2150.15 210.784 Q2152.3 208.261 2152.3 204.094 L2152.3 203.145 L2148.07 203.145 M2156.56 201.386 L2156.56 216.178 L2152.3 216.178 L2152.3 212.242 Q2150.84 214.603 2148.67 215.738 Q2146.49 216.849 2143.34 216.849 Q2139.36 216.849 2137 214.627 Q2134.66 212.381 2134.66 208.631 Q2134.66 204.256 2137.58 202.034 Q2140.52 199.812 2146.33 199.812 L2152.3 199.812 L2152.3 199.395 Q2152.3 196.455 2150.36 194.858 Q2148.44 193.238 2144.94 193.238 Q2142.72 193.238 2140.61 193.77 Q2138.51 194.303 2136.56 195.367 L2136.56 191.432 Q2138.9 190.53 2141.1 190.09 Q2143.3 189.627 2145.38 189.627 Q2151.01 189.627 2153.78 192.543 Q2156.56 195.46 2156.56 201.386 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2182.39 202.914 Q2182.39 198.284 2180.47 195.738 Q2178.58 193.192 2175.13 193.192 Q2171.7 193.192 2169.78 195.738 Q2167.88 198.284 2167.88 202.914 Q2167.88 207.52 2169.78 210.066 Q2171.7 212.613 2175.13 212.613 Q2178.58 212.613 2180.47 210.066 Q2182.39 207.52 2182.39 202.914 M2186.65 212.96 Q2186.65 219.58 2183.71 222.798 Q2180.77 226.039 2174.71 226.039 Q2172.46 226.039 2170.47 225.691 Q2168.48 225.367 2166.61 224.673 L2166.61 220.529 Q2168.48 221.548 2170.31 222.034 Q2172.14 222.52 2174.04 222.52 Q2178.23 222.52 2180.31 220.321 Q2182.39 218.145 2182.39 213.724 L2182.39 211.617 Q2181.08 213.909 2179.02 215.043 Q2176.95 216.178 2174.08 216.178 Q2169.32 216.178 2166.4 212.543 Q2163.48 208.909 2163.48 202.914 Q2163.48 196.895 2166.4 193.261 Q2169.32 189.627 2174.08 189.627 Q2176.95 189.627 2179.02 190.761 Q2181.08 191.895 2182.39 194.187 L2182.39 190.252 L2186.65 190.252 L2186.65 212.96 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2205.66 180.206 Q2202.56 185.53 2201.05 190.738 Q2199.55 195.946 2199.55 201.293 Q2199.55 206.641 2201.05 211.895 Q2202.58 217.127 2205.66 222.428 L2201.95 222.428 Q2198.48 216.988 2196.75 211.733 Q2195.03 206.479 2195.03 201.293 Q2195.03 196.131 2196.75 190.9 Q2198.46 185.668 2201.95 180.206 L2205.66 180.206 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2222.9 216.756 Q2217.56 215.9 2215.08 213.654 Q2212.05 210.9 2212.05 205.275 L2212.05 190.252 L2216.35 190.252 L2216.35 205.113 Q2216.35 209.349 2218.32 211.108 Q2220.03 212.636 2222.9 212.96 L2222.9 190.252 L2227.14 190.252 L2227.14 212.937 Q2230.17 212.613 2231.72 211.085 Q2233.69 209.141 2233.69 205.09 L2233.69 190.252 L2238 190.252 L2238 205.252 Q2238 211.085 2234.96 213.631 Q2232.23 215.923 2227.14 216.733 L2227.14 226.039 L2222.9 226.039 L2222.9 216.756 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip630)\" d=\"M2244.52 180.206 L2248.23 180.206 Q2251.7 185.668 2253.41 190.9 Q2255.15 196.131 2255.15 201.293 Q2255.15 206.479 2253.41 211.733 Q2251.7 216.988 2248.23 222.428 L2244.52 222.428 Q2247.6 217.127 2249.11 211.895 Q2250.64 206.641 2250.64 201.293 Q2250.64 195.946 2249.11 190.738 Q2247.6 185.53 2244.52 180.206 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><polyline clip-path=\"url(#clip630)\" style=\"stroke:#3da44d; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1906.12,250.738 2052.53,250.738 \"/>\n",
       "<path clip-path=\"url(#clip630)\" d=\"M2079.02 246.721 Q2080.43 244.36 2083.92 242.277 Q2085.29 241.467 2089.5 241.467 Q2094.22 241.467 2097.16 245.217 Q2100.13 248.967 2100.13 255.078 Q2100.13 261.189 2097.16 264.939 Q2094.22 268.689 2089.5 268.689 Q2086.65 268.689 2084.59 267.578 Q2082.56 266.443 2081.21 264.129 L2081.21 277.879 L2076.93 277.879 L2076.93 255.309 Q2076.93 249.962 2079.02 246.721 M2095.71 255.078 Q2095.71 250.379 2093.76 247.717 Q2091.84 245.032 2088.46 245.032 Q2085.08 245.032 2083.14 247.717 Q2081.21 250.379 2081.21 255.078 Q2081.21 259.777 2083.14 262.462 Q2085.08 265.124 2088.46 265.124 Q2091.84 265.124 2093.76 262.462 Q2095.71 259.777 2095.71 255.078 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /></svg>\n"
      ],
      "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=\"clip600\">\n",
       "    <rect x=\"0\" y=\"0\" width=\"2400\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<path clip-path=\"url(#clip600)\" d=\"M0 1600 L2400 1600 L2400 0 L0 0  Z\" fill=\"#ffffff\" fill-rule=\"evenodd\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip601\">\n",
       "    <rect x=\"480\" y=\"0\" width=\"1681\" height=\"1600\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<path clip-path=\"url(#clip600)\" d=\"M156.598 1486.45 L2352.76 1486.45 L2352.76 47.2441 L156.598 47.2441  Z\" fill=\"#ffffff\" fill-rule=\"evenodd\" fill-opacity=\"1\"/>\n",
       "<defs>\n",
       "  <clipPath id=\"clip602\">\n",
       "    <rect x=\"156\" y=\"47\" width=\"2197\" height=\"1440\"/>\n",
       "  </clipPath>\n",
       "</defs>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"218.754,1486.45 218.754,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"739.124,1486.45 739.124,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"1259.5,1486.45 1259.5,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"1779.87,1486.45 1779.87,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"2300.24,1486.45 2300.24,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"156.598,1445.72 2352.76,1445.72 \"/>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"156.598,1082.52 2352.76,1082.52 \"/>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"156.598,719.322 2352.76,719.322 \"/>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:2; stroke-opacity:0.1; fill:none\" points=\"156.598,356.126 2352.76,356.126 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,1486.45 2352.76,1486.45 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"218.754,1486.45 218.754,1467.55 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"739.124,1486.45 739.124,1467.55 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1259.5,1486.45 1259.5,1467.55 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1779.87,1486.45 1779.87,1467.55 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"2300.24,1486.45 2300.24,1467.55 \"/>\n",
       "<path clip-path=\"url(#clip600)\" d=\"M196.138 1517.37 Q192.527 1517.37 190.698 1520.93 Q188.893 1524.47 188.893 1531.6 Q188.893 1538.71 190.698 1542.27 Q192.527 1545.82 196.138 1545.82 Q199.772 1545.82 201.578 1542.27 Q203.407 1538.71 203.407 1531.6 Q203.407 1524.47 201.578 1520.93 Q199.772 1517.37 196.138 1517.37 M196.138 1513.66 Q201.948 1513.66 205.004 1518.27 Q208.082 1522.85 208.082 1531.6 Q208.082 1540.33 205.004 1544.94 Q201.948 1549.52 196.138 1549.52 Q190.328 1549.52 187.249 1544.94 Q184.194 1540.33 184.194 1531.6 Q184.194 1522.85 187.249 1518.27 Q190.328 1513.66 196.138 1513.66 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M216.3 1542.97 L221.184 1542.97 L221.184 1548.85 L216.3 1548.85 L216.3 1542.97 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M241.369 1517.37 Q237.758 1517.37 235.929 1520.93 Q234.124 1524.47 234.124 1531.6 Q234.124 1538.71 235.929 1542.27 Q237.758 1545.82 241.369 1545.82 Q245.003 1545.82 246.809 1542.27 Q248.638 1538.71 248.638 1531.6 Q248.638 1524.47 246.809 1520.93 Q245.003 1517.37 241.369 1517.37 M241.369 1513.66 Q247.179 1513.66 250.235 1518.27 Q253.314 1522.85 253.314 1531.6 Q253.314 1540.33 250.235 1544.94 Q247.179 1549.52 241.369 1549.52 Q235.559 1549.52 232.48 1544.94 Q229.425 1540.33 229.425 1531.6 Q229.425 1522.85 232.48 1518.27 Q235.559 1513.66 241.369 1513.66 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M710.861 1544.91 L727.18 1544.91 L727.18 1548.85 L705.236 1548.85 L705.236 1544.91 Q707.898 1542.16 712.481 1537.53 Q717.088 1532.88 718.268 1531.53 Q720.513 1529.01 721.393 1527.27 Q722.296 1525.51 722.296 1523.82 Q722.296 1521.07 720.351 1519.33 Q718.43 1517.6 715.328 1517.6 Q713.129 1517.6 710.676 1518.36 Q708.245 1519.13 705.467 1520.68 L705.467 1515.95 Q708.291 1514.82 710.745 1514.24 Q713.199 1513.66 715.236 1513.66 Q720.606 1513.66 723.8 1516.35 Q726.995 1519.03 726.995 1523.52 Q726.995 1525.65 726.185 1527.57 Q725.398 1529.47 723.291 1532.07 Q722.713 1532.74 719.611 1535.95 Q716.509 1539.15 710.861 1544.91 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M736.995 1542.97 L741.879 1542.97 L741.879 1548.85 L736.995 1548.85 L736.995 1542.97 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M752.111 1514.29 L770.467 1514.29 L770.467 1518.22 L756.393 1518.22 L756.393 1526.7 Q757.411 1526.35 758.43 1526.19 Q759.448 1526 760.467 1526 Q766.254 1526 769.634 1529.17 Q773.013 1532.34 773.013 1537.76 Q773.013 1543.34 769.541 1546.44 Q766.069 1549.52 759.749 1549.52 Q757.573 1549.52 755.305 1549.15 Q753.06 1548.78 750.652 1548.04 L750.652 1543.34 Q752.736 1544.47 754.958 1545.03 Q757.18 1545.58 759.657 1545.58 Q763.661 1545.58 765.999 1543.48 Q768.337 1541.37 768.337 1537.76 Q768.337 1534.15 765.999 1532.04 Q763.661 1529.94 759.657 1529.94 Q757.782 1529.94 755.907 1530.35 Q754.055 1530.77 752.111 1531.65 L752.111 1514.29 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M1226.66 1514.29 L1245.02 1514.29 L1245.02 1518.22 L1230.94 1518.22 L1230.94 1526.7 Q1231.96 1526.35 1232.98 1526.19 Q1234 1526 1235.02 1526 Q1240.8 1526 1244.18 1529.17 Q1247.56 1532.34 1247.56 1537.76 Q1247.56 1543.34 1244.09 1546.44 Q1240.62 1549.52 1234.3 1549.52 Q1232.12 1549.52 1229.85 1549.15 Q1227.61 1548.78 1225.2 1548.04 L1225.2 1543.34 Q1227.28 1544.47 1229.51 1545.03 Q1231.73 1545.58 1234.21 1545.58 Q1238.21 1545.58 1240.55 1543.48 Q1242.89 1541.37 1242.89 1537.76 Q1242.89 1534.15 1240.55 1532.04 Q1238.21 1529.94 1234.21 1529.94 Q1232.33 1529.94 1230.46 1530.35 Q1228.6 1530.77 1226.66 1531.65 L1226.66 1514.29 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M1256.78 1542.97 L1261.66 1542.97 L1261.66 1548.85 L1256.78 1548.85 L1256.78 1542.97 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M1281.84 1517.37 Q1278.23 1517.37 1276.4 1520.93 Q1274.6 1524.47 1274.6 1531.6 Q1274.6 1538.71 1276.4 1542.27 Q1278.23 1545.82 1281.84 1545.82 Q1285.48 1545.82 1287.28 1542.27 Q1289.11 1538.71 1289.11 1531.6 Q1289.11 1524.47 1287.28 1520.93 Q1285.48 1517.37 1281.84 1517.37 M1281.84 1513.66 Q1287.65 1513.66 1290.71 1518.27 Q1293.79 1522.85 1293.79 1531.6 Q1293.79 1540.33 1290.71 1544.94 Q1287.65 1549.52 1281.84 1549.52 Q1276.03 1549.52 1272.96 1544.94 Q1269.9 1540.33 1269.9 1531.6 Q1269.9 1522.85 1272.96 1518.27 Q1276.03 1513.66 1281.84 1513.66 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M1746.19 1514.29 L1768.41 1514.29 L1768.41 1516.28 L1755.86 1548.85 L1750.98 1548.85 L1762.78 1518.22 L1746.19 1518.22 L1746.19 1514.29 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M1777.53 1542.97 L1782.41 1542.97 L1782.41 1548.85 L1777.53 1548.85 L1777.53 1542.97 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M1792.64 1514.29 L1811 1514.29 L1811 1518.22 L1796.93 1518.22 L1796.93 1526.7 Q1797.94 1526.35 1798.96 1526.19 Q1799.98 1526 1801 1526 Q1806.79 1526 1810.17 1529.17 Q1813.55 1532.34 1813.55 1537.76 Q1813.55 1543.34 1810.07 1546.44 Q1806.6 1549.52 1800.28 1549.52 Q1798.11 1549.52 1795.84 1549.15 Q1793.59 1548.78 1791.19 1548.04 L1791.19 1543.34 Q1793.27 1544.47 1795.49 1545.03 Q1797.71 1545.58 1800.19 1545.58 Q1804.19 1545.58 1806.53 1543.48 Q1808.87 1541.37 1808.87 1537.76 Q1808.87 1534.15 1806.53 1532.04 Q1804.19 1529.94 1800.19 1529.94 Q1798.32 1529.94 1796.44 1530.35 Q1794.59 1530.77 1792.64 1531.65 L1792.64 1514.29 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2252.31 1544.91 L2259.95 1544.91 L2259.95 1518.55 L2251.64 1520.21 L2251.64 1515.95 L2259.9 1514.29 L2264.58 1514.29 L2264.58 1544.91 L2272.22 1544.91 L2272.22 1548.85 L2252.31 1548.85 L2252.31 1544.91 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2291.66 1517.37 Q2288.05 1517.37 2286.22 1520.93 Q2284.42 1524.47 2284.42 1531.6 Q2284.42 1538.71 2286.22 1542.27 Q2288.05 1545.82 2291.66 1545.82 Q2295.29 1545.82 2297.1 1542.27 Q2298.93 1538.71 2298.93 1531.6 Q2298.93 1524.47 2297.1 1520.93 Q2295.29 1517.37 2291.66 1517.37 M2291.66 1513.66 Q2297.47 1513.66 2300.53 1518.27 Q2303.61 1522.85 2303.61 1531.6 Q2303.61 1540.33 2300.53 1544.94 Q2297.47 1549.52 2291.66 1549.52 Q2285.85 1549.52 2282.77 1544.94 Q2279.72 1540.33 2279.72 1531.6 Q2279.72 1522.85 2282.77 1518.27 Q2285.85 1513.66 2291.66 1513.66 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2311.82 1542.97 L2316.71 1542.97 L2316.71 1548.85 L2311.82 1548.85 L2311.82 1542.97 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2336.89 1517.37 Q2333.28 1517.37 2331.45 1520.93 Q2329.65 1524.47 2329.65 1531.6 Q2329.65 1538.71 2331.45 1542.27 Q2333.28 1545.82 2336.89 1545.82 Q2340.53 1545.82 2342.33 1542.27 Q2344.16 1538.71 2344.16 1531.6 Q2344.16 1524.47 2342.33 1520.93 Q2340.53 1517.37 2336.89 1517.37 M2336.89 1513.66 Q2342.7 1513.66 2345.76 1518.27 Q2348.84 1522.85 2348.84 1531.6 Q2348.84 1540.33 2345.76 1544.94 Q2342.7 1549.52 2336.89 1549.52 Q2331.08 1549.52 2328 1544.94 Q2324.95 1540.33 2324.95 1531.6 Q2324.95 1522.85 2328 1518.27 Q2331.08 1513.66 2336.89 1513.66 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,1486.45 156.598,47.2441 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,1445.72 175.496,1445.72 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,1082.52 175.496,1082.52 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,719.322 175.496,719.322 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"156.598,356.126 175.496,356.126 \"/>\n",
       "<path clip-path=\"url(#clip600)\" d=\"M63.4226 1431.51 Q59.8115 1431.51 57.9828 1435.08 Q56.1773 1438.62 56.1773 1445.75 Q56.1773 1452.86 57.9828 1456.42 Q59.8115 1459.96 63.4226 1459.96 Q67.0569 1459.96 68.8624 1456.42 Q70.6911 1452.86 70.6911 1445.75 Q70.6911 1438.62 68.8624 1435.08 Q67.0569 1431.51 63.4226 1431.51 M63.4226 1427.81 Q69.2328 1427.81 72.2883 1432.42 Q75.367 1437 75.367 1445.75 Q75.367 1454.48 72.2883 1459.08 Q69.2328 1463.67 63.4226 1463.67 Q57.6125 1463.67 54.5338 1459.08 Q51.4782 1454.48 51.4782 1445.75 Q51.4782 1437 54.5338 1432.42 Q57.6125 1427.81 63.4226 1427.81 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M83.5845 1457.12 L88.4688 1457.12 L88.4688 1463 L83.5845 1463 L83.5845 1457.12 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M108.654 1431.51 Q105.043 1431.51 103.214 1435.08 Q101.409 1438.62 101.409 1445.75 Q101.409 1452.86 103.214 1456.42 Q105.043 1459.96 108.654 1459.96 Q112.288 1459.96 114.094 1456.42 Q115.922 1452.86 115.922 1445.75 Q115.922 1438.62 114.094 1435.08 Q112.288 1431.51 108.654 1431.51 M108.654 1427.81 Q114.464 1427.81 117.52 1432.42 Q120.598 1437 120.598 1445.75 Q120.598 1454.48 117.52 1459.08 Q114.464 1463.67 108.654 1463.67 Q102.844 1463.67 99.765 1459.08 Q96.7095 1454.48 96.7095 1445.75 Q96.7095 1437 99.765 1432.42 Q102.844 1427.81 108.654 1427.81 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M65.0198 1068.32 Q61.4087 1068.32 59.58 1071.88 Q57.7745 1075.42 57.7745 1082.55 Q57.7745 1089.66 59.58 1093.22 Q61.4087 1096.77 65.0198 1096.77 Q68.6541 1096.77 70.4596 1093.22 Q72.2883 1089.66 72.2883 1082.55 Q72.2883 1075.42 70.4596 1071.88 Q68.6541 1068.32 65.0198 1068.32 M65.0198 1064.61 Q70.83 1064.61 73.8855 1069.22 Q76.9642 1073.8 76.9642 1082.55 Q76.9642 1091.28 73.8855 1095.89 Q70.83 1100.47 65.0198 1100.47 Q59.2097 1100.47 56.131 1095.89 Q53.0754 1091.28 53.0754 1082.55 Q53.0754 1073.8 56.131 1069.22 Q59.2097 1064.61 65.0198 1064.61 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M85.1818 1093.92 L90.066 1093.92 L90.066 1099.8 L85.1818 1099.8 L85.1818 1093.92 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M104.279 1095.86 L120.598 1095.86 L120.598 1099.8 L98.6539 1099.8 L98.6539 1095.86 Q101.316 1093.11 105.899 1088.48 Q110.506 1083.83 111.686 1082.48 Q113.932 1079.96 114.811 1078.23 Q115.714 1076.47 115.714 1074.78 Q115.714 1072.02 113.77 1070.29 Q111.848 1068.55 108.746 1068.55 Q106.547 1068.55 104.094 1069.31 Q101.663 1070.08 98.8854 1071.63 L98.8854 1066.91 Q101.709 1065.77 104.163 1065.19 Q106.617 1064.61 108.654 1064.61 Q114.024 1064.61 117.219 1067.3 Q120.413 1069.98 120.413 1074.48 Q120.413 1076.6 119.603 1078.53 Q118.816 1080.42 116.709 1083.02 Q116.131 1083.69 113.029 1086.91 Q109.927 1090.1 104.279 1095.86 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M62.9365 705.121 Q59.3254 705.121 57.4967 708.686 Q55.6912 712.227 55.6912 719.357 Q55.6912 726.463 57.4967 730.028 Q59.3254 733.57 62.9365 733.57 Q66.5707 733.57 68.3763 730.028 Q70.205 726.463 70.205 719.357 Q70.205 712.227 68.3763 708.686 Q66.5707 705.121 62.9365 705.121 M62.9365 701.417 Q68.7467 701.417 71.8022 706.024 Q74.8809 710.607 74.8809 719.357 Q74.8809 728.084 71.8022 732.69 Q68.7467 737.274 62.9365 737.274 Q57.1264 737.274 54.0477 732.69 Q50.9921 728.084 50.9921 719.357 Q50.9921 710.607 54.0477 706.024 Q57.1264 701.417 62.9365 701.417 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M83.0984 730.723 L87.9827 730.723 L87.9827 736.602 L83.0984 736.602 L83.0984 730.723 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M111.015 706.116 L99.2095 724.565 L111.015 724.565 L111.015 706.116 M109.788 702.042 L115.668 702.042 L115.668 724.565 L120.598 724.565 L120.598 728.454 L115.668 728.454 L115.668 736.602 L111.015 736.602 L111.015 728.454 L95.4132 728.454 L95.4132 723.94 L109.788 702.042 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M63.2606 341.924 Q59.6495 341.924 57.8208 345.489 Q56.0152 349.031 56.0152 356.16 Q56.0152 363.267 57.8208 366.831 Q59.6495 370.373 63.2606 370.373 Q66.8948 370.373 68.7004 366.831 Q70.5291 363.267 70.5291 356.16 Q70.5291 349.031 68.7004 345.489 Q66.8948 341.924 63.2606 341.924 M63.2606 338.221 Q69.0707 338.221 72.1263 342.827 Q75.205 347.41 75.205 356.16 Q75.205 364.887 72.1263 369.494 Q69.0707 374.077 63.2606 374.077 Q57.4504 374.077 54.3717 369.494 Q51.3162 364.887 51.3162 356.16 Q51.3162 347.41 54.3717 342.827 Q57.4504 338.221 63.2606 338.221 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M83.4225 367.526 L88.3067 367.526 L88.3067 373.406 L83.4225 373.406 L83.4225 367.526 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M109.071 354.262 Q105.922 354.262 104.071 356.415 Q102.242 358.568 102.242 362.318 Q102.242 366.044 104.071 368.22 Q105.922 370.373 109.071 370.373 Q112.219 370.373 114.047 368.22 Q115.899 366.044 115.899 362.318 Q115.899 358.568 114.047 356.415 Q112.219 354.262 109.071 354.262 M118.353 339.609 L118.353 343.869 Q116.594 343.035 114.788 342.596 Q113.006 342.156 111.246 342.156 Q106.617 342.156 104.163 345.281 Q101.733 348.406 101.385 354.725 Q102.751 352.711 104.811 351.646 Q106.871 350.558 109.348 350.558 Q114.557 350.558 117.566 353.73 Q120.598 356.878 120.598 362.318 Q120.598 367.642 117.45 370.859 Q114.302 374.077 109.071 374.077 Q103.075 374.077 99.9039 369.494 Q96.7326 364.887 96.7326 356.16 Q96.7326 347.966 100.621 343.105 Q104.51 338.221 111.061 338.221 Q112.82 338.221 114.603 338.568 Q116.408 338.915 118.353 339.609 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><polyline clip-path=\"url(#clip602)\" style=\"stroke:#009af9; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"218.754,1445.72 228.39,1445.72 238.027,1445.72 247.663,1445.72 257.3,1445.72 266.936,1445.72 276.573,1445.72 286.209,1445.71 295.846,1445.71 305.482,1445.71 315.119,1445.71 324.755,1445.71 334.392,1445.71 344.028,1445.71 353.665,1445.71 363.301,1445.71 372.938,1445.71 382.574,1445.7 392.211,1445.7 401.847,1445.69 411.484,1445.69 421.12,1445.68 430.757,1445.67 440.393,1445.65 450.03,1445.64 459.666,1445.61 469.303,1445.58 478.939,1445.55 488.576,1445.51 498.212,1445.45 507.849,1445.38 517.485,1445.3 527.122,1445.19 536.758,1445.06 546.395,1444.9 556.031,1444.7 565.668,1444.47 575.304,1444.18 584.941,1443.82 594.577,1443.4 604.214,1442.89 613.85,1442.28 623.487,1441.55 633.123,1440.68 642.76,1439.65 652.396,1438.43 662.033,1436.99 671.669,1435.3 681.305,1433.32 690.942,1431.02 700.578,1428.34 710.215,1425.24 719.851,1421.67 729.488,1417.56 739.124,1412.86 748.761,1407.49 758.397,1401.39 768.034,1394.48 777.67,1386.69 787.307,1377.93 796.943,1368.13 806.58,1357.2 816.216,1345.06 825.853,1331.63 835.489,1316.83 845.126,1300.58 854.762,1282.83 864.399,1263.51 874.035,1242.57 883.672,1219.98 893.308,1195.69 902.945,1169.72 912.581,1142.05 922.218,1112.71 931.854,1081.73 941.491,1049.19 951.127,1015.15 960.764,979.714 970.4,942.996 980.037,905.13 989.673,866.267 999.31,826.574 1008.95,786.229 1018.58,745.422 1028.22,704.353 1037.86,663.223 1047.49,622.241 1057.13,581.614 1066.77,541.546 1076.4,502.237 1086.04,463.879 1095.67,426.658 1105.31,390.745 1114.95,356.303 1124.58,323.481 1134.22,292.415 1143.86,263.228 1153.49,236.029 1163.13,210.915 1172.77,187.971 1182.4,167.269 1192.04,148.873 1201.68,132.833 1211.31,119.194 1220.95,107.99 1230.59,99.2491 1240.22,92.9911 1249.86,89.2307 1259.5,87.9763 1269.13,89.2307 1278.77,92.9911 1288.4,99.2491 1298.04,107.99 1307.68,119.194 1317.31,132.833 1326.95,148.873 1336.59,167.269 1346.22,187.971 1355.86,210.915 1365.5,236.029 1375.13,263.228 1384.77,292.415 1394.41,323.481 1404.04,356.303 1413.68,390.745 1423.32,426.658 1432.95,463.879 1442.59,502.237 1452.23,541.546 1461.86,581.614 1471.5,622.241 1481.13,663.223 1490.77,704.353 1500.41,745.422 1510.04,786.229 1519.68,826.574 1529.32,866.267 1538.95,905.13 1548.59,942.996 1558.23,979.714 1567.86,1015.15 1577.5,1049.19 1587.14,1081.73 1596.77,1112.71 1606.41,1142.05 1616.05,1169.72 1625.68,1195.69 1635.32,1219.98 1644.96,1242.57 1654.59,1263.51 1664.23,1282.83 1673.86,1300.58 1683.5,1316.83 1693.14,1331.63 1702.77,1345.06 1712.41,1357.2 1722.05,1368.13 1731.68,1377.93 1741.32,1386.69 1750.96,1394.48 1760.59,1401.39 1770.23,1407.49 1779.87,1412.86 1789.5,1417.56 1799.14,1421.67 1808.78,1425.24 1818.41,1428.34 1828.05,1431.02 1837.69,1433.32 1847.32,1435.3 1856.96,1436.99 1866.59,1438.43 1876.23,1439.65 1885.87,1440.68 1895.5,1441.55 1905.14,1442.28 1914.78,1442.89 1924.41,1443.4 1934.05,1443.82 1943.69,1444.18 1953.32,1444.47 1962.96,1444.7 1972.6,1444.9 1982.23,1445.06 1991.87,1445.19 2001.51,1445.3 2011.14,1445.38 2020.78,1445.45 2030.42,1445.51 2040.05,1445.55 2049.69,1445.58 2059.32,1445.61 2068.96,1445.64 2078.6,1445.65 2088.23,1445.67 2097.87,1445.68 2107.51,1445.69 2117.14,1445.69 2126.78,1445.7 2136.42,1445.7 2146.05,1445.71 2155.69,1445.71 2165.33,1445.71 2174.96,1445.71 2184.6,1445.71 2194.24,1445.71 2203.87,1445.71 2213.51,1445.71 2223.15,1445.71 2232.78,1445.71 2242.42,1445.72 2252.05,1445.72 2261.69,1445.72 2271.33,1445.72 2280.96,1445.72 2290.6,1445.72 \"/>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#e26f46; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"218.754,1445.72 228.39,1445.72 238.027,1445.72 247.663,1445.72 257.3,1445.72 266.936,1445.72 276.573,1445.72 286.209,1445.72 295.846,1445.72 305.482,1445.72 315.119,1445.72 324.755,1445.72 334.392,1445.72 344.028,1445.72 353.665,1445.72 363.301,1445.72 372.938,1445.72 382.574,1445.72 392.211,1445.72 401.847,1445.72 411.484,1445.72 421.12,1445.72 430.757,1445.72 440.393,1445.72 450.03,1445.72 459.666,1445.72 469.303,1445.72 478.939,1445.72 488.576,1445.72 498.212,1445.72 507.849,1445.72 517.485,1445.72 527.122,1445.72 536.758,1445.72 546.395,1445.72 556.031,1445.72 565.668,1445.72 575.304,1445.72 584.941,1445.72 594.577,1445.72 604.214,1445.72 613.85,1445.72 623.487,1445.72 633.123,1445.72 642.76,1445.72 652.396,1445.72 662.033,1445.72 671.669,1445.72 681.305,1445.72 690.942,1445.72 700.578,1445.72 710.215,1445.72 719.851,1445.72 729.488,1445.72 739.124,1445.72 748.761,1445.72 758.397,1445.72 768.034,1445.72 777.67,1445.72 787.307,1445.72 796.943,1445.72 806.58,1445.72 816.216,1445.72 825.853,1445.72 835.489,1445.72 845.126,1445.72 854.762,1445.72 864.399,1445.72 874.035,1445.72 883.672,1445.72 893.308,1445.72 902.945,1445.72 912.581,1445.72 922.218,1445.72 931.854,1445.72 941.491,1445.72 951.127,1445.72 960.764,1445.72 970.4,1445.72 980.037,1445.72 989.673,1445.72 999.31,1445.72 1008.95,1445.72 1018.58,1445.72 1028.22,1445.72 1037.86,1445.72 1047.49,1445.72 1057.13,1445.72 1066.77,1445.72 1076.4,1445.72 1086.04,1445.72 1095.67,1445.72 1105.31,1445.72 1114.95,1445.72 1124.58,1445.72 1134.22,1445.72 1143.86,1445.72 1153.49,1445.72 1163.13,1445.72 1172.77,1445.72 1182.4,1445.72 1192.04,1445.72 1201.68,1445.72 1211.31,1445.72 1220.95,1445.72 1230.59,1445.72 1240.22,1445.72 1249.86,1445.72 1259.5,1445.72 1269.13,1445.72 1278.77,1445.72 1288.4,1445.72 1298.04,1445.72 1307.68,1445.72 1317.31,1445.72 1326.95,1445.72 1336.59,1445.72 1346.22,1445.72 1355.86,1445.72 1365.5,1445.72 1375.13,1445.72 1384.77,1445.72 1394.41,1445.72 1404.04,1445.72 1413.68,1445.72 1423.32,1445.72 1432.95,1445.72 1442.59,1445.72 1452.23,1445.72 1461.86,1445.72 1471.5,1445.72 1481.13,1445.72 1490.77,1445.72 1500.41,1445.72 1510.04,1445.72 1519.68,1445.72 1529.32,1445.72 1538.95,1445.72 1548.59,1445.72 1558.23,1445.72 1567.86,1445.72 1577.5,1445.72 1587.14,1445.72 1596.77,1445.72 1606.41,1445.72 1616.05,1445.72 1625.68,1445.72 1635.32,1445.72 1644.96,1445.72 1654.59,1445.72 1664.23,1445.72 1673.86,1445.72 1683.5,1445.72 1693.14,1445.72 1702.77,1445.72 1712.41,1445.72 1722.05,1445.72 1731.68,1445.72 1741.32,1445.72 1750.96,1445.72 1760.59,1445.72 1770.23,1445.72 1779.87,1445.72 1789.5,1445.72 1799.14,1445.72 1808.78,1445.72 1818.41,1445.72 1828.05,1445.72 1837.69,1445.72 1847.32,1445.72 1856.96,1445.72 1866.59,1445.72 1876.23,1445.72 1885.87,1445.72 1895.5,1445.72 1905.14,1445.72 1914.78,1445.72 1924.41,1445.72 1934.05,1445.72 1943.69,1445.72 1953.32,1445.72 1962.96,1445.72 1972.6,1445.72 1982.23,1445.72 1991.87,1445.72 2001.51,1445.72 2011.14,1445.72 2020.78,1445.72 2030.42,1445.72 2040.05,1445.72 2049.69,1445.72 2059.32,1445.72 2068.96,1445.72 2078.6,1445.72 2088.23,1445.72 2097.87,1445.72 2107.51,1445.72 2117.14,1445.72 2126.78,1445.72 2136.42,1445.72 2146.05,1445.72 2155.69,1445.72 2165.33,1445.72 2174.96,1445.72 2184.6,1445.72 2194.24,1445.72 2203.87,1445.72 2213.51,1445.72 2223.15,1445.72 2232.78,1445.72 2242.42,1445.72 2252.05,1445.72 2261.69,1445.72 2271.33,1445.72 2280.96,1445.72 2290.6,1445.72 \"/>\n",
       "<polyline clip-path=\"url(#clip602)\" style=\"stroke:#3da44d; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"218.754,1445.72 228.39,1445.72 238.027,1445.72 247.663,1445.72 257.3,1445.72 266.936,1445.72 276.573,1445.72 286.209,1445.72 295.846,1445.72 305.482,1445.72 315.119,1445.72 324.755,1445.72 334.392,1445.72 344.028,1445.72 353.665,1445.72 363.301,1445.72 372.938,1445.72 382.574,1445.72 392.211,1445.72 401.847,1445.72 411.484,1445.72 421.12,1445.72 430.757,1445.72 440.393,1445.72 450.03,1445.72 459.666,1445.72 469.303,1445.72 478.939,1445.72 488.576,1445.72 498.212,1445.72 507.849,1445.72 517.485,1445.72 527.122,1445.72 536.758,1445.72 546.395,1445.72 556.031,1445.72 565.668,1445.71 575.304,1445.71 584.941,1445.71 594.577,1445.71 604.214,1445.71 613.85,1445.71 623.487,1445.71 633.123,1445.7 642.76,1445.7 652.396,1445.69 662.033,1445.67 671.669,1445.66 681.305,1445.63 690.942,1445.6 700.578,1445.55 710.215,1445.48 719.851,1445.4 729.488,1445.28 739.124,1445.12 748.761,1444.91 758.397,1444.63 768.034,1444.27 777.67,1443.8 787.307,1443.19 796.943,1442.4 806.58,1441.4 816.216,1440.14 825.853,1438.55 835.489,1436.57 845.126,1434.12 854.762,1431.11 864.399,1427.43 874.035,1422.99 883.672,1417.65 893.308,1411.29 902.945,1403.77 912.581,1394.94 922.218,1384.65 931.854,1372.76 941.491,1359.13 951.127,1343.63 960.764,1326.13 970.4,1306.55 980.037,1284.79 989.673,1260.82 999.31,1234.63 1008.95,1206.22 1018.58,1175.66 1028.22,1143.06 1037.86,1108.55 1047.49,1072.3 1057.13,1034.55 1066.77,995.534 1076.4,955.539 1086.04,914.872 1095.67,873.861 1105.31,832.845 1114.95,792.175 1124.58,752.202 1134.22,713.274 1143.86,675.732 1153.49,639.903 1163.13,606.097 1172.77,574.605 1182.4,545.694 1192.04,519.605 1201.68,496.555 1211.31,476.732 1220.95,460.294 1230.59,447.374 1240.22,438.072 1249.86,432.462 1259.5,430.587 1269.13,432.462 1278.77,438.072 1288.4,447.374 1298.04,460.294 1307.68,476.732 1317.31,496.555 1326.95,519.605 1336.59,545.694 1346.22,574.605 1355.86,606.097 1365.5,639.903 1375.13,675.732 1384.77,713.274 1394.41,752.202 1404.04,792.175 1413.68,832.845 1423.32,873.861 1432.95,914.872 1442.59,955.539 1452.23,995.534 1461.86,1034.55 1471.5,1072.3 1481.13,1108.55 1490.77,1143.06 1500.41,1175.66 1510.04,1206.22 1519.68,1234.63 1529.32,1260.82 1538.95,1284.79 1548.59,1306.55 1558.23,1326.13 1567.86,1343.63 1577.5,1359.13 1587.14,1372.76 1596.77,1384.65 1606.41,1394.94 1616.05,1403.77 1625.68,1411.29 1635.32,1417.65 1644.96,1422.99 1654.59,1427.43 1664.23,1431.11 1673.86,1434.12 1683.5,1436.57 1693.14,1438.55 1702.77,1440.14 1712.41,1441.4 1722.05,1442.4 1731.68,1443.19 1741.32,1443.8 1750.96,1444.27 1760.59,1444.63 1770.23,1444.91 1779.87,1445.12 1789.5,1445.28 1799.14,1445.4 1808.78,1445.48 1818.41,1445.55 1828.05,1445.6 1837.69,1445.63 1847.32,1445.66 1856.96,1445.67 1866.59,1445.69 1876.23,1445.7 1885.87,1445.7 1895.5,1445.71 1905.14,1445.71 1914.78,1445.71 1924.41,1445.71 1934.05,1445.71 1943.69,1445.71 1953.32,1445.71 1962.96,1445.72 1972.6,1445.72 1982.23,1445.72 1991.87,1445.72 2001.51,1445.72 2011.14,1445.72 2020.78,1445.72 2030.42,1445.72 2040.05,1445.72 2049.69,1445.72 2059.32,1445.72 2068.96,1445.72 2078.6,1445.72 2088.23,1445.72 2097.87,1445.72 2107.51,1445.72 2117.14,1445.72 2126.78,1445.72 2136.42,1445.72 2146.05,1445.72 2155.69,1445.72 2165.33,1445.72 2174.96,1445.72 2184.6,1445.72 2194.24,1445.72 2203.87,1445.72 2213.51,1445.72 2223.15,1445.72 2232.78,1445.72 2242.42,1445.72 2252.05,1445.72 2261.69,1445.72 2271.33,1445.72 2280.96,1445.72 2290.6,1445.72 \"/>\n",
       "<path clip-path=\"url(#clip600)\" d=\"M1881.72 302.578 L2279.55 302.578 L2279.55 95.2176 L1881.72 95.2176  Z\" fill=\"#ffffff\" fill-rule=\"evenodd\" fill-opacity=\"1\"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#000000; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1881.72,302.578 2279.55,302.578 2279.55,95.2176 1881.72,95.2176 1881.72,302.578 \"/>\n",
       "<polyline clip-path=\"url(#clip600)\" style=\"stroke:#009af9; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1906.12,147.058 2052.53,147.058 \"/>\n",
       "<path clip-path=\"url(#clip600)\" d=\"M2092.12 142.393 Q2091.4 141.977 2090.54 141.791 Q2089.71 141.583 2088.69 141.583 Q2085.08 141.583 2083.14 143.944 Q2081.21 146.282 2081.21 150.68 L2081.21 164.338 L2076.93 164.338 L2076.93 138.412 L2081.21 138.412 L2081.21 142.44 Q2082.56 140.078 2084.71 138.944 Q2086.86 137.787 2089.94 137.787 Q2090.38 137.787 2090.91 137.856 Q2091.45 137.903 2092.09 138.018 L2092.12 142.393 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2117.72 150.31 L2117.72 152.393 L2098.14 152.393 Q2098.41 156.791 2100.77 159.106 Q2103.16 161.398 2107.4 161.398 Q2109.85 161.398 2112.14 160.796 Q2114.46 160.194 2116.72 158.99 L2116.72 163.018 Q2114.43 163.99 2112.02 164.5 Q2109.62 165.009 2107.14 165.009 Q2100.94 165.009 2097.3 161.398 Q2093.69 157.787 2093.69 151.629 Q2093.69 145.264 2097.12 141.537 Q2100.57 137.787 2106.4 137.787 Q2111.63 137.787 2114.66 141.166 Q2117.72 144.523 2117.72 150.31 M2113.46 149.06 Q2113.41 145.565 2111.49 143.481 Q2109.59 141.398 2106.45 141.398 Q2102.88 141.398 2100.73 143.412 Q2098.6 145.426 2098.27 149.083 L2113.46 149.06 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2136.49 151.305 Q2131.33 151.305 2129.34 152.486 Q2127.35 153.666 2127.35 156.514 Q2127.35 158.782 2128.83 160.125 Q2130.33 161.444 2132.9 161.444 Q2136.45 161.444 2138.58 158.944 Q2140.73 156.421 2140.73 152.254 L2140.73 151.305 L2136.49 151.305 M2144.99 149.546 L2144.99 164.338 L2140.73 164.338 L2140.73 160.402 Q2139.27 162.763 2137.09 163.898 Q2134.92 165.009 2131.77 165.009 Q2127.79 165.009 2125.43 162.787 Q2123.09 160.541 2123.09 156.791 Q2123.09 152.416 2126.01 150.194 Q2128.95 147.972 2134.76 147.972 L2140.73 147.972 L2140.73 147.555 Q2140.73 144.615 2138.78 143.018 Q2136.86 141.398 2133.37 141.398 Q2131.14 141.398 2129.04 141.93 Q2126.93 142.463 2124.99 143.527 L2124.99 139.592 Q2127.33 138.69 2129.52 138.25 Q2131.72 137.787 2133.81 137.787 Q2139.43 137.787 2142.21 140.703 Q2144.99 143.62 2144.99 149.546 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2153.76 128.319 L2158.02 128.319 L2158.02 164.338 L2153.76 164.338 L2153.76 128.319 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2177.16 128.366 Q2174.06 133.69 2172.56 138.898 Q2171.05 144.106 2171.05 149.453 Q2171.05 154.801 2172.56 160.055 Q2174.08 165.287 2177.16 170.588 L2173.46 170.588 Q2169.99 165.148 2168.25 159.893 Q2166.54 154.639 2166.54 149.453 Q2166.54 144.291 2168.25 139.06 Q2169.96 133.828 2173.46 128.366 L2177.16 128.366 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2194.41 164.916 Q2189.06 164.06 2186.58 161.814 Q2183.55 159.06 2183.55 153.435 L2183.55 138.412 L2187.86 138.412 L2187.86 153.273 Q2187.86 157.509 2189.83 159.268 Q2191.54 160.796 2194.41 161.12 L2194.41 138.412 L2198.64 138.412 L2198.64 161.097 Q2201.68 160.773 2203.23 159.245 Q2205.2 157.301 2205.2 153.25 L2205.2 138.412 L2209.5 138.412 L2209.5 153.412 Q2209.5 159.245 2206.47 161.791 Q2203.74 164.083 2198.64 164.893 L2198.64 174.199 L2194.41 174.199 L2194.41 164.916 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2216.03 128.366 L2219.73 128.366 Q2223.2 133.828 2224.92 139.06 Q2226.65 144.291 2226.65 149.453 Q2226.65 154.639 2224.92 159.893 Q2223.2 165.148 2219.73 170.588 L2216.03 170.588 Q2219.11 165.287 2220.61 160.055 Q2222.14 154.801 2222.14 149.453 Q2222.14 144.106 2220.61 138.898 Q2219.11 133.69 2216.03 128.366 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><polyline clip-path=\"url(#clip600)\" style=\"stroke:#e26f46; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1906.12,198.898 2052.53,198.898 \"/>\n",
       "<path clip-path=\"url(#clip600)\" d=\"M2076.93 190.252 L2081.19 190.252 L2081.19 216.178 L2076.93 216.178 L2076.93 190.252 M2076.93 180.159 L2081.19 180.159 L2081.19 185.553 L2076.93 185.553 L2076.93 180.159 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2110.29 195.229 Q2111.89 192.358 2114.11 190.993 Q2116.33 189.627 2119.34 189.627 Q2123.39 189.627 2125.59 192.474 Q2127.79 195.298 2127.79 200.529 L2127.79 216.178 L2123.51 216.178 L2123.51 200.668 Q2123.51 196.942 2122.19 195.136 Q2120.87 193.33 2118.16 193.33 Q2114.85 193.33 2112.93 195.53 Q2111.01 197.729 2111.01 201.525 L2111.01 216.178 L2106.72 216.178 L2106.72 200.668 Q2106.72 196.918 2105.4 195.136 Q2104.08 193.33 2101.33 193.33 Q2098.07 193.33 2096.15 195.553 Q2094.22 197.752 2094.22 201.525 L2094.22 216.178 L2089.94 216.178 L2089.94 190.252 L2094.22 190.252 L2094.22 194.28 Q2095.68 191.895 2097.72 190.761 Q2099.76 189.627 2102.56 189.627 Q2105.38 189.627 2107.35 191.062 Q2109.34 192.497 2110.29 195.229 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2148.07 203.145 Q2142.9 203.145 2140.91 204.326 Q2138.92 205.506 2138.92 208.354 Q2138.92 210.622 2140.4 211.965 Q2141.91 213.284 2144.48 213.284 Q2148.02 213.284 2150.15 210.784 Q2152.3 208.261 2152.3 204.094 L2152.3 203.145 L2148.07 203.145 M2156.56 201.386 L2156.56 216.178 L2152.3 216.178 L2152.3 212.242 Q2150.84 214.603 2148.67 215.738 Q2146.49 216.849 2143.34 216.849 Q2139.36 216.849 2137 214.627 Q2134.66 212.381 2134.66 208.631 Q2134.66 204.256 2137.58 202.034 Q2140.52 199.812 2146.33 199.812 L2152.3 199.812 L2152.3 199.395 Q2152.3 196.455 2150.36 194.858 Q2148.44 193.238 2144.94 193.238 Q2142.72 193.238 2140.61 193.77 Q2138.51 194.303 2136.56 195.367 L2136.56 191.432 Q2138.9 190.53 2141.1 190.09 Q2143.3 189.627 2145.38 189.627 Q2151.01 189.627 2153.78 192.543 Q2156.56 195.46 2156.56 201.386 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2182.39 202.914 Q2182.39 198.284 2180.47 195.738 Q2178.58 193.192 2175.13 193.192 Q2171.7 193.192 2169.78 195.738 Q2167.88 198.284 2167.88 202.914 Q2167.88 207.52 2169.78 210.066 Q2171.7 212.613 2175.13 212.613 Q2178.58 212.613 2180.47 210.066 Q2182.39 207.52 2182.39 202.914 M2186.65 212.96 Q2186.65 219.58 2183.71 222.798 Q2180.77 226.039 2174.71 226.039 Q2172.46 226.039 2170.47 225.691 Q2168.48 225.367 2166.61 224.673 L2166.61 220.529 Q2168.48 221.548 2170.31 222.034 Q2172.14 222.52 2174.04 222.52 Q2178.23 222.52 2180.31 220.321 Q2182.39 218.145 2182.39 213.724 L2182.39 211.617 Q2181.08 213.909 2179.02 215.043 Q2176.95 216.178 2174.08 216.178 Q2169.32 216.178 2166.4 212.543 Q2163.48 208.909 2163.48 202.914 Q2163.48 196.895 2166.4 193.261 Q2169.32 189.627 2174.08 189.627 Q2176.95 189.627 2179.02 190.761 Q2181.08 191.895 2182.39 194.187 L2182.39 190.252 L2186.65 190.252 L2186.65 212.96 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2205.66 180.206 Q2202.56 185.53 2201.05 190.738 Q2199.55 195.946 2199.55 201.293 Q2199.55 206.641 2201.05 211.895 Q2202.58 217.127 2205.66 222.428 L2201.95 222.428 Q2198.48 216.988 2196.75 211.733 Q2195.03 206.479 2195.03 201.293 Q2195.03 196.131 2196.75 190.9 Q2198.46 185.668 2201.95 180.206 L2205.66 180.206 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2222.9 216.756 Q2217.56 215.9 2215.08 213.654 Q2212.05 210.9 2212.05 205.275 L2212.05 190.252 L2216.35 190.252 L2216.35 205.113 Q2216.35 209.349 2218.32 211.108 Q2220.03 212.636 2222.9 212.96 L2222.9 190.252 L2227.14 190.252 L2227.14 212.937 Q2230.17 212.613 2231.72 211.085 Q2233.69 209.141 2233.69 205.09 L2233.69 190.252 L2238 190.252 L2238 205.252 Q2238 211.085 2234.96 213.631 Q2232.23 215.923 2227.14 216.733 L2227.14 226.039 L2222.9 226.039 L2222.9 216.756 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><path clip-path=\"url(#clip600)\" d=\"M2244.52 180.206 L2248.23 180.206 Q2251.7 185.668 2253.41 190.9 Q2255.15 196.131 2255.15 201.293 Q2255.15 206.479 2253.41 211.733 Q2251.7 216.988 2248.23 222.428 L2244.52 222.428 Q2247.6 217.127 2249.11 211.895 Q2250.64 206.641 2250.64 201.293 Q2250.64 195.946 2249.11 190.738 Q2247.6 185.53 2244.52 180.206 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /><polyline clip-path=\"url(#clip600)\" style=\"stroke:#3da44d; stroke-linecap:round; stroke-linejoin:round; stroke-width:4; stroke-opacity:1; fill:none\" points=\"1906.12,250.738 2052.53,250.738 \"/>\n",
       "<path clip-path=\"url(#clip600)\" d=\"M2079.02 246.721 Q2080.43 244.36 2083.92 242.277 Q2085.29 241.467 2089.5 241.467 Q2094.22 241.467 2097.16 245.217 Q2100.13 248.967 2100.13 255.078 Q2100.13 261.189 2097.16 264.939 Q2094.22 268.689 2089.5 268.689 Q2086.65 268.689 2084.59 267.578 Q2082.56 266.443 2081.21 264.129 L2081.21 277.879 L2076.93 277.879 L2076.93 255.309 Q2076.93 249.962 2079.02 246.721 M2095.71 255.078 Q2095.71 250.379 2093.76 247.717 Q2091.84 245.032 2088.46 245.032 Q2085.08 245.032 2083.14 247.717 Q2081.21 250.379 2081.21 255.078 Q2081.21 259.777 2083.14 262.462 Q2085.08 265.124 2088.46 265.124 Q2091.84 265.124 2093.76 262.462 Q2095.71 259.777 2095.71 255.078 Z\" fill=\"#000000\" fill-rule=\"nonzero\" fill-opacity=\"1\" /></svg>\n"
      ]
     },
     "metadata": {},
     "execution_count": 10
    }
   ],
   "cell_type": "code",
   "source": [
    "using Plots\n",
    "\n",
    "p = plot(x, real.(ψ), label=\"real(ψ)\")\n",
    "plot!(p, x, imag.(ψ), label=\"imag(ψ)\")\n",
    "plot!(p, x, ρ, label=\"ρ\")"
   ],
   "metadata": {},
   "execution_count": 10
  },
  {
   "cell_type": "markdown",
   "source": [
    "The `energy_hamiltonian` function can be used to get the energy and\n",
    "effective Hamiltonian (derivative of the energy with respect to the density matrix)\n",
    "of a particular state (ψ, occupation).\n",
    "The density ρ associated to this state is precomputed\n",
    "and passed to the routine as an optimization."
   ],
   "metadata": {}
  },
  {
   "outputs": [],
   "cell_type": "code",
   "source": [
    "E, ham = energy_hamiltonian(basis, scfres.ψ, scfres.occupation; ρ=scfres.ρ)\n",
    "@assert E.total == scfres.energies.total"
   ],
   "metadata": {},
   "execution_count": 11
  },
  {
   "cell_type": "markdown",
   "source": [
    "Now the Hamiltonian contains all the blocks corresponding to $k$-points. Here, we just\n",
    "have one $k$-point:"
   ],
   "metadata": {}
  },
  {
   "outputs": [],
   "cell_type": "code",
   "source": [
    "H = ham.blocks[1];"
   ],
   "metadata": {},
   "execution_count": 12
  },
  {
   "cell_type": "markdown",
   "source": [
    "`H` can be used as a linear operator (efficiently using FFTs), or converted to a dense matrix:"
   ],
   "metadata": {}
  },
  {
   "outputs": [],
   "cell_type": "code",
   "source": [
    "ψ11 = scfres.ψ[1][:, 1] # first k-point, first eigenvector\n",
    "Hmat = Array(H)  # This is now just a plain Julia matrix,\n",
    "#                  which we can compute and store in this simple 1D example\n",
    "@assert norm(Hmat * ψ11 - H * ψ11) < 1e-10"
   ],
   "metadata": {},
   "execution_count": 13
  },
  {
   "cell_type": "markdown",
   "source": [
    "Let's check that ψ11 is indeed an eigenstate:"
   ],
   "metadata": {}
  },
  {
   "outputs": [
    {
     "output_type": "execute_result",
     "data": {
      "text/plain": "1.0563325319244079e-7"
     },
     "metadata": {},
     "execution_count": 14
    }
   ],
   "cell_type": "code",
   "source": [
    "norm(H * ψ11 - dot(ψ11, H * ψ11) * ψ11)"
   ],
   "metadata": {},
   "execution_count": 14
  },
  {
   "cell_type": "markdown",
   "source": [
    "Build a finite-differences version of the GPE operator $H$, as a sanity check:"
   ],
   "metadata": {}
  },
  {
   "outputs": [
    {
     "output_type": "execute_result",
     "data": {
      "text/plain": "0.00022340550154961012"
     },
     "metadata": {},
     "execution_count": 15
    }
   ],
   "cell_type": "code",
   "source": [
    "A = Array(Tridiagonal(-ones(N - 1), 2ones(N), -ones(N - 1)))\n",
    "A[1, end] = A[end, 1] = -1\n",
    "K = A / dx^2 / 2\n",
    "V = Diagonal(pot.(x) + C .* α .* (ρ.^(α-1)))\n",
    "H_findiff = K + V;\n",
    "maximum(abs.(H_findiff*ψ - (dot(ψ, H_findiff*ψ) / dot(ψ, ψ)) * ψ))"
   ],
   "metadata": {},
   "execution_count": 15
  }
 ],
 "nbformat_minor": 3,
 "metadata": {
  "language_info": {
   "file_extension": ".jl",
   "mimetype": "application/julia",
   "name": "julia",
   "version": "1.11.4"
  },
  "kernelspec": {
   "name": "julia-1.11",
   "display_name": "Julia 1.11.4",
   "language": "julia"
  }
 },
 "nbformat": 4
}