{ "cells": [ { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "using LinearAlgebra" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# Classical vs. Modified Gram–Schmidt" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This is an example adapted from Trefethen and Bau, *Numerical Linear Algebra*, lecture 9. We construct a square matrix $A$ with exponentially graded singular values between $2^{-1}$ and $2^{-80}$, perform QR factorization in various ways, and check the accuracy and orthogonality of the results." ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "80×80 Array{Float64,2}:\n", " 0.00231288 -0.00188713 0.00115339 … -0.000777549 -0.00196706 \n", " 0.00447612 0.0119956 0.00404125 0.00753016 0.00083697 \n", " 0.000730346 0.00228595 0.00128298 0.00071391 0.000510262\n", " -0.00165406 0.00106365 0.000736025 4.42455e-5 0.00118117 \n", " 0.003155 -0.00739048 0.000132767 -0.0014877 -0.00334916 \n", " 0.000346554 -0.00636389 -0.00172838 … -0.00366515 -0.000830444\n", " -0.00235921 -0.00687301 -0.000762344 -0.00675667 0.000155346\n", " 0.000472659 0.00683781 0.00318707 0.00217628 0.00242981 \n", " -0.00434494 0.00826401 0.00196047 0.00511342 0.00391427 \n", " -0.00131739 0.0088279 0.0003912 0.0055408 0.00253646 \n", " -0.000332517 -0.00588532 -0.000676342 … -0.00485406 0.000115519\n", " -0.000461368 -0.000437043 -0.00198776 0.000728049 0.00106017 \n", " 0.00210663 0.00154606 9.68214e-5 0.00219168 -0.000360623\n", " ⋮ ⋱ \n", " -0.000297086 -0.000820249 -0.000581252 -0.00193116 0.000442198\n", " -0.00257615 0.0033895 -0.000220947 0.000960119 0.00199674 \n", " 0.00276878 0.00638691 0.00230822 … -0.00269343 0.00112297 \n", " -0.00217747 0.00454191 0.000260071 0.00577961 0.00127744 \n", " -0.00365117 -0.00727022 -0.00288899 -0.00559812 2.08879e-5 \n", " -0.00160264 0.00883925 0.000150251 0.00276509 0.00288517 \n", " -0.000872615 0.00200195 0.000797861 -0.00144595 0.00204208 \n", " -0.00141865 -0.00393479 -0.00054564 … -0.00428038 -0.00132189 \n", " -0.000173487 0.00315454 -0.00112636 0.00162983 0.000267716\n", " 0.00138451 0.00430689 0.000437716 0.000179362 0.000557466\n", " -0.000243029 -0.0106028 -0.00260152 -0.00219126 -0.0019942 \n", " -0.00246574 -0.00295037 -0.000836444 -0.00417088 0.000334032" ] }, "execution_count": 2, "metadata": {}, "output_type": "execute_result" } ], "source": [ "U = qr(randn(80,80)).Q # U = random 80×80 orthogonal matrix\n", "V = qr(randn(80,80)).Q # V = random 80×80 orthogonal matrix\n", "Σ = Diagonal(2.0 .^ (-1:-1:-80)) # Σ = diagonal matrix with entries 2^-1, 2^-2, ..., 2^-80\n", "A = U * Σ * V'" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now, let's implement classical and modified Gram–Schmidt algorithms to return the (reduced) QR factorizations of a matrix A:" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "mgs (generic function with 1 method)" ] }, "execution_count": 3, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# Classical Gram–Schmidt (Trefethen algorithm 7.1), implemented in the simplest way\n", "# (We could make it faster by unrolling loops to avoid temporaries arrays etc.)\n", "function clgs(A)\n", " m,n = size(A)\n", " T = float(eltype(A))\n", " Q = similar(A, T)\n", " R = zeros(T,n,n)\n", " @views for j = 1:n\n", " aⱼ = A[:,j]\n", " vⱼ = copy(aⱼ) # use copy so that modifying vⱼ doesn't change aⱼ\n", " for i = 1:j-1\n", " qᵢ = Q[:,i]\n", " R[i,j] = qᵢ'aⱼ\n", " vⱼ .-= R[i,j] .* qᵢ\n", " end\n", " R[j,j] = norm(vⱼ)\n", " Q[:,j] .= vⱼ ./ R[j,j]\n", " end\n", " return Q, R\n", "end\n", "\n", "# Modified Gram–Schmidt (Trefethen algorithm 8.1)\n", "function mgs(A)\n", " m,n = size(A)\n", " T = float(eltype(A))\n", " Q = similar(A, T)\n", " R = zeros(T,n,n)\n", " @views for j = 1:n\n", " aⱼ = A[:,j]\n", " vⱼ = copy(aⱼ)\n", " for i = 1:j-1\n", " qᵢ = Q[:,i]\n", " R[i,j] = qᵢ'vⱼ # ⟵ NOTICE: mgs has vⱼ, clgs has aⱼ\n", " vⱼ .-= R[i,j] .* qᵢ\n", " end\n", " R[j,j] = norm(vⱼ)\n", " Q[:,j] .= vⱼ ./ R[j,j]\n", " end\n", " return Q, R\n", "end" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next, use these on our matrix $A$ from above:" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [], "source": [ "QC, RC = clgs(A);\n", "QM, RM = mgs(A);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Finally, we plot the diagonal elements of $R$, because $r_{jj} = \\Vert v_j \\Vert$, the norm of the projection of $a_j$ on to the space orthogonal to the previous columns." ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAkwAAAG0CAYAAADATXgqAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAAPYQAAD2EBqD+naQAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4xLjEsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy8QZhcZAAAgAElEQVR4nOzde1iUZcLH8e8wAmp5SE3TxENtZWZpgFYmhdha2radSyvTcmvt5HawzROIoGnWZr2rtmvualaalmm75bulAmpab0ZiJdVmYWpq5NbiKUGG+/1jnGEGZoYZmAMDv891zaU8c88z98Awz4/7aDHGGERERETEq5hIV0BERESkvlNgEhEREamBApOIiIhIDRSYRERERGqgwCQiIiJSAwUmERERkRooMImIiIjUoEmkK9AQVFRUsHfvXlq0aIHFYol0dURERMQPxhgOHTpEp06diInx3YakwBQEe/fuJSEhIdLVEBERkVrYvXs3nTt39llGgSkIWrRoAdi/4S1btoxwbURERMQfBw8eJCEhwXkd90WBKQgc3XAtW7ZUYBIREYky/gyn0aBvERERkRooMImIiIjUQIFJREREpAYKTCIiIiI1UGASERERqYECk4iIiEgNFJhEREREaqDAJCIiIlIDBSYRERGRGigwiYiIiNRAgUlERESkBgpMIiIiIjVQYKrn1q1bx5EjRyJdDRERkUZNgakeKygoYOjQofTt25fPP/880tURERFptBSY6rGjR4/Stm1bvvjiC/r27cuCBQswxkS6WiIiIo2OAlM91r9/fwoKCrjyyis5duwY99xzD3fccQeHDh2KdNVEREQaFQWmeq59+/asXr2aGTNmYLVaWbJkCcnJyRQUFES6aiIiIo2GAlMUiImJYfz48eTl5dG5c2f+/e9/c/HFF/PCCy+oi05ERCQMFJiiyIABAygoKODqq6+mtLSU+++/n2HDhlFSUhLpqomIiDRoCkxRpm3btvzzn//kmWeeoUmTJixfvpzExETy8/MjXTUREZEGS4EpClksFh577DE2btxI165d+fbbb+nfvz9//vOf1UUnIiISAgpMUeziiy9m69atXHfddZSVlTF27FhuvPFGfv7550hXTUREpEFRYIpCNhvk5cHSpbBt2ym8/vqbPP/888TGxrJy5UoSExP56KOPIl1NERGRBkOBKZpkZlI4PJtu3WDgQLjtNvu/3btbuOKDEjaPGsUZZ5zBzp07ufTSS3n22WfVRSciIhIECkxRpPArKz1fy2DUnmy343ftyabnaxk0P5TAJ598wk033UR5eTmPPfYY1157LT/99FOEaiwiItIwKDBFCZsNrnw/nXSyyCaDydhD02SyySKDDLK4alM6J5/ciuXLlzNv3jzi4+P55z//SZ8+fdi8eXOEX4GIiEj0shj12dTZwYMHadWqFSUlJbRs2TIkz5GXZ+9+A3tIyiaDUuKIp4x0sphGOgC5uZCaai9XUFDALbfcwtdff43VamX69Ok8/vjjxMTEYLPBxo2wbx907AgpKWC1hqTqIiIi9VIg12+1MEWJffsq/z+NdGdYKiXOGZaqluvTpw/5+fncfvvt2Gw2xo8fz9XnnMP714+vNg6qWzcoHJ4NmZnhekkiIiJRQ4GpHnOdDffDD5XHJ5PtDEvxlDm758DeWuSqRYsWvPzyyyxYsICmTZvyrx07uHXVUwzcM9qtnGMcVOFXamYSERGpqkmkKyAeZGZS+JWVK99PZ8+eysNWK7xrG8QgcpzdcI7uOQuwKCGdlJTqp7NYLIwePZrk5ItITr6FveVf8Cp/53J2sp73mMyTznFQizalU2RT95yIiIgrBaZ6yDkbDty62xxhaR1pzuPTSMcCZJHBsEvBak33eE6An3/uRXn5FuABKniJ9eQwkFgexFSOg9ptH9vkGAclIiIiCkz1jmM23CggmwwAZ0uSIyxdaV0HtsrHLEpIZ9il0PMcm8fzOQZ3FxYCnAQsAgbSnFHkYugN/EB/52Ncx0GJiIiIAlO9s3Ej7NlT2bJkX0JgmvtsOBvMng0dOrjOcKvSsuSlWw/sY6AGkkNH4GYsbMcAvwbSgQw6dlR/nIiIiCsN+gZ2795NamoqPXv25IILLuD111+PWF38nQ3XoQMMH27vOvM03sjbIpeOMU9p5LGELLZzmD4kAQbIIj5+EGeeuTckr01ERCRaKTABTZo04bnnnqOwsJC1a9fyyCOPcOTIkYjUxXWWWyCz4Vxn1K1b532RS0c3X+XaTc0p4GOu4yZOBkpL15OY2Jt33303tC9UREQkiigwAR07dqRPnz4AtG/fnjZt2kRsO5GUFOjcGdJPhJt0smhKqTP8pJNNQgKVs+E87C93xRU4u+FySCWbDI4RTzYZ5JDqttClQ37C6ywbOpbeHTpw4MABrrrqKiZMmEB5eXlYX7+IiEh9FBWBacOGDVxzzTV06tQJi8XCqlWrqpWZN28e3bt3p2nTpiQlJbFx48ZaPdfHH39MRUUFCQkJda12rVit8O6Ayu1OXGfDZZBFFhn869JsZzdcTV1vuaS5desNItd5zsmTYckS++rgRUUw9J3n+XDnTu677z4AZs6cSWpqKrt37w7fN0BERKQeiorAdOTIEXr37s2cOXM83r9s2TIefvhhJk2axNatW0lJSWHIkCHs2rXLWSYpKYlevXpVu+3dWzle5z//+Q933nkn8+fP91mf0tJSDh486HYLpp7n2CgclsXCzu6tQIsS0ikcluWcDedrfzlH6xTgtVtv0KDq46CaNm3KvHnzWLZsGS1atGDTpk306dOHd955J6ivUUREJKqYKAOYlStXuh3r16+fGTNmjNuxHj16mPHjx/t93mPHjpmUlBSzePHiGstOmTLFYB8l7XYrKSnx+/n8UV5uTG6uMUuW2P8tL3e/PzfXGLDfJpNlDJhjxBkDZjJZzmOTyXIrk06WSUiofr6qduzYYZKSkpyv77HHHjNlZWVBfY0iIiKRUlJS4vf1OypamHwpKysjPz+fwYMHux0fPHgwmzdv9uscxhhGjRpFWloaI0aMqLH8hAkTKCkpcd5C1WVltdpbf7zNhvM1ow5wtjLV1K3nzZlnnsmmTZt46KGHAPjTn/5ESkoK3333XZBeoYiISHSI+sB04MABbDYbHTp0cDveoUMH9u/f79c5Nm3axLJly1i1ahV9+vShT58+fPbZZ17Lx8fH07JlS7dbJPiaUTfQZfsUV1W79WoSHx/P//zP//Dmm2/SunVr/u///o8+ffqwatUqt5l5eXn2LkIREZGGqMEsXGmxWNy+NsZUO+bNgAEDqKioCEW1Qsoxo+6uPfZB4lX3l8sjjc6dYdEiKC72scilH66//nouvPBCbr31Vj766COuv/56Tj75Dxw+PAtOtGh17mwfsN7zHBtkZgb1tYqIiERS1Aemdu3aYbVaq7UmFRcXV2t1amgcM+p6vlZ9Rp1zf7kB0HNQ4AHJk27durFx40Z+f8FAFn21mcOHnwc2AcuAM7hrj70uuQOz2L/UNaAF5elFREQiJuq75OLi4khKSmLNmjVux9esWUP//v29PKrh8HdGnT/86WKzWuNYe2QTt3AbbQD4GLiQDM5ztnKl5aZz2232daG6dYPC4dlqcRIRkagWFS1Mhw8fZseOHc6vi4qKKCgooE2bNnTp0oVHH32UESNGkJyczCWXXML8+fPZtWsXY8aMiWCtwyQzk57ATpdNdgPuevOy75ynLjbHXnfLeZWHOI18nmUzB8mikANAS467ndrR6lQ4LIueQXnBIiIiERDyOXtBkJub63Ea/8iRI51l5s6da7p27Wri4uJMYmKiWb9+fdjqF8i0xPpo+zD35Qcct/QTyxBsH5blLLtkiXErc4hYM97lZ9IHzP38oVbLGIiIiIRTINdvizHGRCytNRAHDx6kVatWlJSURGzGXG3ZbPZus1F7st2WIXAMHM8gi0UJ6RQV2cci5eXZu9qgcoHMUuLIpYwbac5RjnIyMBcrd2Jzm6mXm2tfHkFERKQ+COT6HfVjmKRuHF1s06hcMdyx75z963R274Y//9k+tslm87zX3Say+JqjdKEbh4GR2BhNDNMY53wu13WjREREookCUyPna/FLKzbnViqPPFK5se/9/7EvY5BLqtvMvBfI4ht2kg5YgL9TQXu6A18C7utGiYiIRBMFpkbO1+KXKWx026PO4eJfcgD4sFlatfM1AVJJxbCWkziZYn4glt60afMyKSmhfCUiIiKho8DUyDkWv6zaxZZOFoPIYR1p1Tb2HUgeGWQxr206a9fCkiWQM7ByzNMgcoFBHOFrunMGxynjp5/u5J577ubIkSORfcEiIiK1oEHfQRDNg77Bvk6SY/HLbJetVByDuteRxiBynC1PHgdye12awMZ1rQczb3suFcbQs2dPli9fznnnnRfeFykiIlJFINdvBaYgiPbA5C3sgD00WbExgRnO7rqmlDrvX7LEvjmwg83jelCQk7OeW24Zzn/+s4/4+GbMmTOH0aPv8nv7GhERkWBTYAqzqA9MJ7iGnR9+sA/0BvflA7y2MHnjFsZ+BEYA7wLwm64XsGT4UFrMmOGzLtpiRUREQkHLCkitWK328DN8ODz0kPexTfavs0lIoMaB3IVfWen5Wgaj9mQDpwKrgRlYiOHt7z7lgrkvsm3btsoHZGZSODybbt3s6z1pixUREakPFJjEI8fGvllU39g3gyyyyOBfl2b7bPWx2eDK9yvXd7IPHI9hMsfZQAUtacnOQ//hoosu4i9/+QvGmCoBq5Jzi5Wv1MwkIiLhFxV7yUlkODf2fT8dXMY2LUpIZ9il1Lixr+uimMCJ0DTN2a13kPuBkZSWvsN9991HTk4umzbN5+4TZQHnquOO4LZoUzpFNnXPiYhIeGkMUxA0lDFM3tR2PNHSpfYuNYdjxHsYOF7B7bfPZtmy8ZSXlwNnAsuYzOraj5sSERHxg8YwSVC5jm1KTfW/dcfXopiVi2HG8LvfPcbGjRtp164r8A3Qn2m05hixzsdMc1nuYMUK+552Nt8NXCIiIkGjwCQh42tRTMfA8c6d7cGnqOhiHnlkK3AdUAaM5TaOs/9EaFrLIOd558zRQHAREQkvBSYJGV8Dx3NJJYsM7v9PNldcYe+6mzTpFGJi3uQBziIWWAl05HTmkswgctxCE2gguIiIhI8Ck4SUc+B453S344596C46sS+dw3sVVzCHr3mOZKA7sJMHKeA+fkWaS2hyHQh+1aZ0dc+JiEhIadB3EDT0Qd/B4DpwvH17GDUKRu2p7KpzzIZzbMVypXUdNlsJ8DvgDQAuoR1vc4CTNBBcRESCQCt9h5kCU2Dy8uxjkMD3KuKzZ8OOHYa5c18AHgHK6AwsA5Jq2KJFRESkJpolJ/Xavn2V/59GujMsVZ0N16ED3HSTBbgf+JBTaMMe4DLgOcqYxFRnWdcZeSIiIsGmwCRh599yA5VrPtln2r3Nd/zEeZyPDRgPfEImjzLeOdNu6VItNyAiIqGhwCRh589yA4596lxn2j1NFtvZBswHmvK/wGs8xZDi0c6ZdlpuQEREQkGBScIu0H3q3GfaWYB7gI841dKOvcDfyv4OTAcqAC03ICIiwadB30GgQd+1kJlJ4VdWrnw/nT0u+9QlJMC/Ls2271NXpYWo6ky7O+88TNu9KXxGwYkSv+ZhLmQ2s+z7ziWkU1SkfedERMQzzZILMwWm2qvtPnWuM+2u4QbWsZKjwGnAIEbxKgsBLTcgIiLeBXL9bhKmOol45NinLlCuM+3+yZtsJY47OM524FVeAroC6ezbp+YlERGpO41hkqhUdaZdH46zkVh+B4ABpgJXEBe3z+PjRUREAqHAJFHJ00y7NpRxGlm8CsQRB+Rx3329ee+99yJdXRERiXIKTBKVfM20+5IsPqOMc1qfxo8//shVV13FpEmTKC8vj3CtRUQkWikwSdTytrHvooR0yodlUfDAaO677z6MMTz55JOkpaWxx3VKnoiIiJ80S+6EJk2a0KtXLwCSk5NZsGCB34/VLLnIqmmm3fLly/nd737HoUOHaNu2LYsXL2bo0KGRq7CIiNQLWlagFtq1a8eBAwdq9VgFpvrvm2++4dZbbyU/Px+Axx9/nOnTpxMbGxvhmomISKRo812RKs4880w2bdrEgw8+BMDTTz9Nnz6X8e2330W4ZiIiEg2iIjBt2LCBa665hk6dOmGxWFi1alW1MvPmzaN79+40bdqUpKQkNm7cGNBzHDx4kKSkJAYMGMD69euDVXWpLzIz+WbULFat+h/gTaAVhYUf8qtfXciclNu075yIiPgUFYHpyJEj9O7dmzlz5ni8f9myZTz88MNMmjSJrVu3kpKSwpAhQ9i1a5ezTFJSEr169ap227t3LwA7d+4kPz+fv/zlL9x5550cPHjQa31KS0s5ePCg203qt8KvrPR8LYNRe7KB64GtQD+M+ZmH3l/KiKXvUlZWFuFaiohIvWWiDGBWrlzpdqxfv35mzJgxbsd69Ohhxo8fX6vnuOqqq8yWLVu83j9lyhSDfXVEt1tJSUmtnk9Co7zcmNxcY155xZhTTzVmMlnGgJlMlgFjJpBhHnX5+SUnJ5tvvvkm0tUWEZEwKSkp8fv6HRUtTL6UlZWRn5/P4MGD3Y4PHjyYzZs3+3WOn3/+mdLSUgD27NlDYWEhZ5xxhtfyEyZMoKSkxHnbvXt37V+ABF9mJoXDs+nWzb7f3B13wI8/2tdoWkca2WRwjHieJIvmZAH/AE7h448/5sILL+SNN96I8AsQEZH6JuoD04EDB7DZbHTo0MHteIcOHdi/f79f5/jiiy9ITk6md+/e/OY3v+H555+nTZs2XsvHx8fTsmVLt5vUH+7db5Umk80gcijHSjxllBJ3YsHLa4ACzj67PwcPHuTmm2/mgQce4NixYxGpv4iI1D9RH5gcLBaL29fGmGrHvOnfvz+fffYZ27Zto6CggOuuuy4UVZQwsNngyvfTSSeLbDKYjD00TT6xhco60miCjVLiiKfMeT90Ye7cPJ544gnAPomgf//+fP311xF6JSIiUp9EfWBq164dVqu1WmtScXFxtVYnafg2boQ9e+zdb47QdIx4Z1gaRA7pZNGUUuf96WSTkAADB8Yyc+ZMVq9eTbt27di6dSuJiYm89tprkX5ZIiISYVEfmOLi4khKSmLNmjVux9esWUP//v0jVCuJlH37Kv8/jXRnS1I5VmdYct13LoMssshg9SXZbNwIS5dCs2ZDyM8vICUlhcOHDzN8+HDuvfdefvnllwi9KhERibQmka6APw4fPsyOHTucXxcVFVFQUECbNm3o0qULjz76KCNGjCA5OZlLLrmE+fPns2vXLsaMGRPBWkskdOxY+f/JZDvHKsVTxjrSnGHJ4ZSWNr5uncZ7/2vjseWVxzt3Pp13+qfxegpMf/99XnzxRT788EOWL19Ojx49wvRqRESk3gj9pL26y83N9TiNf+TIkc4yc+fONV27djVxcXEmMTHRrF+/Pmz1C2RaooRWebkxnTsbk15lCQHXJQVOPdW+1EBurjGf3+peznFzPH77sCzz3nvvmfbt2xvAnHTSSWbx4sWRfpkiIhIEgVy/oyIw1XcKTPXL9mH2sJPuIwQZUxmuqq7P5Pg6nSyTkGAvt3fvXjNw4EBnWL/rrrvM4cOHI/xKRUSkLhrVOkwiVfU8x0bhsCwWdnbvfluUkE7hsCx6nmMDfA8Qt3+dzu7d9nIdO3ZkzZo1ZGZmYrFYWLhwIf369WP79u2ReIkiIhJmFmOMiXQlol0gux1L+Nhs9rCzb599bFNKClitlfcvXQq33Vb59THinWOemlLqPL5kCQwfXlkuNzeX2267jf3799OsWTPmzJnDXXfd5fcyFiIiUj8Ecv1WC5M0WFYrpKbaw05qqntYAt8DxCvXZ3IvBzBw4EC2bdvGr3/9a3755RdGjx7NnXfeSUnJYfLy7EEsL88e2EREpGFQYJJGKyUFOneG9BOLWnpbnyklpfpj27dvz7/+9S+mT59OjMXCK6+8Qrt2yQwc+Cm33WbfkqVbNygcng2ZmeF+aSIiEmQKTNJoWa3w7oBsssggw4/1maq2GsXExDBx4kQWDbyL04Hy8q+AfsBfAcNde7Lp+VoGhV9Zqz+5iIhElahYh0kkVJwDxN9Phz2Vx72vz2QPWT3PsUFmJjYbTPz337iVU/mSp1hNKTCGnsxlHJ+Rgf3ci9ZBcbHnsVQi4q6m8Yd1LV9f1Jd615d61Hshn7PXCGhZgehXXm5fl2nJEv/XZzLGXtZx3yQyzdNgmpxYeuB0mpnRjHF7PNiXMtg+LMuYKVMi+IpF6qcVK+y/I1V/Z1asCE75QFX9bCgvD855Q13vaKuHL6H6GRijdZjCToGpYfFnfabOnY1Zu9aYBx90/6A5Rpz5AEzCidAUB+ZKrjZQUS10fXZLVsg+BETCIdgXshUrjLFYTLU/MiwW+63qRTzQ8rWpTyjCRG3qHYrQUF/qUVMdQxnoFJjCTIGpYXFtNXKEpGPEGQNmHanVWp0ct7WkOcv+B8zZ9HAudNmDngZ+dp5vLWnGaq3+IfDZLVmmaOQUhSip94J9IXP8oeLpd8txEXcsJFub8rV5faEIY7WpdyhCQ7DrUV8CXaAUmMJMgalhWbKkequRIwhVbXWqGpbWkuYMWhVgBjPExGA1gOkK5v9OlKl6jilM8Xhc3XdSH4XiQub6h4qvW25u7coHIpRhzN96T55sL/v6676/18uX1y6oBLMeYEzbtsENUqEOxA4KTGGmwNSw+GphmkxWtdBUNSxVfezd/N50w97SFAsGnjWTmOrXOaqOmWqMwt0FoHrU/PyhuJBV/UPF223JktqVD0Qow5i/9XbcqrZE13S/vy1Pwa5HoEEqkj8DVwpMYabA1LD4s4Fv1SBVNei4hqZ1pJr/grkBi3F00cFvzWOMr3YOT8/nuqddY1NfBqQ2lnr4E8ZCdSELVQvT7NmBh0t/w8SDDwYeWv2td21v/rbyhboevuoHxkyd6vvnEspA7EqBKcwUmBoebxv4uoYm1666qh+ia9d6Cl0V5ip+Y+KcoSnB5NHE7RyeWrRqewGKduEYv6B6uJ/fnzAW7AuZI6S98ooxp57q+TU6XqenMUzeykNgrS+uYXH27MBCQCCh1Z96ByOU1PRHVjjqUZfvn1qYGigFpgZoyhSzfViWx64HR6tRTcHGW+j6HfeZX2EPTVYw07Eam8s5vAWxuv4lFU2C3e0TiXEUwew686cejpmbtXm+QMJYMC9knkKat9fna5acvxd9X+epWo9AuqACDa2B1ru2N8f4I2/vhXDVozbfv5oCncYwRSkFpobL9aLnudXIR9eZl9BltRrzDy4zw3G0NGHO5CxTDG4z7Vyfw98LUEMR6otyqMdReHvOUA/Qrc1rrO3stLpeyLyFNE+3hITA1mGqKeyceqq9RcvXYObaXPQDCa3+hsVg3AJdxyoSN2+z8jwFOs2Si2IKTI2Ht1YjX4OzvS+KOdXAfANNDWDaEmc24D7TzvFcjW0MU7C6feq6zszkyYHXI5AgUNNMIsfxqut9+XsB8mcWVW1CYV0vZDWFNHAPNTW99+vSnRbooOpAbjWFVke9/X2f1fZW088lXPUI9H3meK9Vfa/4CtCBUmAKMwWmRsRLq1FCQgDT/6udY5tpSzsDmBgwAxlowGag8c6SC0YLU7DWmQmkHv4Egap1AM8ziR5/PDh/+dc0jqe24dTbhcyfVrT6NAvNn9vs2XULrf4GyEDGY9Vm1lowxjVVfV7HezdY3Xqe/gjSSt8NiAJT4xOMX2DXc3x+2wRzTdfepnIW3a8N7A8siEU5T92fvrp9aur+CPSiHEjLkLeLUKRmHtXlIl6X8FL19+D11yMzcLw2P/dAL+K1PW+gXZTeWu2qBlFHV2Kg79maQmig9Sgv9xyeaxukwj30QIEpzBSYJBgqKirMggV/N/HxzQxg2rQ5zbz33rpIVyssAvnA9dUq4/phHkh3WqAtQ56ChzGhad0Ixc31Ih7qMUmhHjheVShmfzmCQV3OW9vxdoGO3/Ln/R7seji+7/4EKV/vybpMXKgtBaYwU2CSYNq+fbvp2bOnAYzFYjFTpkwx5fV8AFMgLW6eWiUCWUG46tfBuBjWpvXA0wUkGlqYPF3EQz0mKVQDx70J1uyvqvWoy3kDXWbB39AQ6Pgjf9elClY3mOt5pk71/j7z9LsdjrXOFJjCTIFJgu3IkSPm7rvvNo4uuoEDB5q9e/dGuloeBTILLdAZTVX/6nR01QUjLLheDP1tGZo82fcFJJxr27iu91Xb56s6WL22g2sjMXC8JnWd/RXIEgSBvvZQCPa6VKHiqzXZ359BMCkwhZkCk4TKyy+/bE466SQDmFNPPdW8++67bvdHepuMQLphajNGqOrFJlgtOMEcx+PtexLq0FRT61BtXk9t30/BHjge7JXLa1oU01OY8FWPQMfbhWuWa7DWpQo1T98/X3UM5fdPgSnMFJgklL744gtzwQUXGEcX3aRJk8zx9HSPs/XCuVlvIN0wtRkj5OlCG6wxQlUvhsHuIgrl2jb+zu6rqeUumBehYA4cD1U3UW0GMwfr3OEMJJF+LwQqlOPZ/KHAFGYKTBJqR48eNb///e+No4su6dSuZjfuC1tC6JchqM16N7UdI+Tpw7Iu56mpOy3YFz1/ZpAFOpPI3/WjfM2iCsVFPByrMgdjH71QtmiFurUsELX9PY2EcO0Z540CU5gpMEm4LF261LRo0cIAphnNzTuEbrNefy74/n7Q1bZlKNDBwnW9IISri6immUQJCZ7XYQq0Lr5eT7C7c0PZyhLMffRC2Y0d6S5yTyIdSGqiFqZGRoFJwumVV742kGgcrU2PgzlErDEEbyuVYHYp1baFydeg20DHaYRj37m6qGml72Ct9+UrpAVjAHAoAmdd9vOTyAeSmoRrzzhvFJjCTIFJwsn+F+MxAw86Q9MlYP5NrNsHTU1dUN7UZXC2tw+62szgCcb6M5Ea1FqfBbO1xpNgh7/6fsGv7yIdSPwRyTFgCkxhpsAk4eR6AbmJYabVidB0CphselUb1+RoPfBnMHhdB2f7+qAL9qBbf7oMI2KXg7wAACAASURBVDWGpL6KVGtNXVq06nuXUjSoT4PSfdUxEr+/Cky18Oyzz5qePXuac8891zz00EOmoqLC78cqMEk4OS56jgHeD/CIgb7O1qaHwUwgw+2Dx9/B4MGatu/tgy4SY4SCWT7aRaK1pq4tWmphCo76NCjdm0j8PgZy/W6C8OOPPzJnzhy2b99ObGwsl112GR9++CGXXHJJpKsGwJEj3u+zWqFpU//KxsRAs2a1K3v0qP1XzBOLBZo3r13ZX36Bigrv9TjppNqVPXYMbLbglG3e3F5vgNJSKC8PTtlmzezfZ4CyMjh+3L+yNhusvHgm574xk/HMYC7jgTIgA/gzz3GMvmTxAIeYy7P8kZk8gb3sKxvHU3jQ/r4B+3vH8f/jx6GoyHsdajJzJrRvD6edBmlpEB9vP15ebv9eAFx5JWzfDps2wf799rIDB1a+h13LehIXB7Gxld+HY8eql+nb1/5vbGzla/NU9q234PHHYe/eymOdO8Ps2TBkiPc6xMba6wH29+Mvv3gv26RJ5ffBGPvvRjDKBvJ771p23z7v5Vw5yjnOa7O5/8wuvdT+fajpM8Jmg4ce8vx54Dg2dixccUXlz6rqZ0RyMnTq5P5zcmWx2H9uKSn6jPBV1vV376ef4PTT7d8zm833+6fqZ0RZmfey8fH293GgZV1/7x2/v1D5O+v6ex9Roc9v9V9xcbHp0qWL+fnnn80vv/xi+vbta3bs2OH340PdwuTrr6qhQ93LNm/uvezll7uXbdfOe9nkZPeyXbt6L9uzp3vZnj29l+3a1b1scrL3su3auZe9/HLvZZs3dy87dKjv75urm27yXfbw4cqyI0f6LltcXFn2/vt9ly0qqiw7bpzvsp9/Xll2yhTfZeMYYADTCszt/NHvv8rnzPHvL3l/bsuXV553+XLfZRcurCz79tu+y86ZU1m2ppaHWbMqy370kX/19mfs1rhxlectKvJd9v77K8sWF/suO3JkZdnDh32Xvekm9/ewr7KunxGBttb4+ow480z3Ovj6jAjkFshnhONn5mgl0WeEXU2fER99VFl21iz/3gvG1PwZ8fbblWUXLvRdtrafEcEWyPU7JsJ5zS8bNmzgmmuuoVOnTlgsFlatWlWtzLx58+jevTtNmzYlKSmJjRs3+n3+U089lXHjxtGlSxc6derEFVdcwZlnnhnMlyASNmW8DlxCCfAqPv7clmqMiXQNQislpbJ1w5uEBHs58N1C8M038OabwatbbcTEwBtvwA03RLYe0jhYjKn/HxH/+7//y6ZNm0hMTOTGG29k5cqVXHfddc77ly1bxogRI5g3bx6XXnopf/3rX1mwYAGFhYV06dIFgKSkJEo9tPW/9957NGvWjGHDhrF06VKaNWvGkCFDyMrK4rLLLvOrfgcPHqRVq1aUlJTQsmXL4LxoF+qSC7xsQ2tuj4uzN6Xv2wft2sFFF1U2kzts2ABDh9r/P4FMbExlFrFALO05jWLeA87kj3+E1FR7t8pJJ1Vvbn/rLbj9du91cXX66TBrFlx7rftxb83t3l6bo7k9GF1yDq5dZ65lXb9PvqxeDZ4+AlzP+8Yb8Ic/uHcXdeoETz9t/57Uty45gCVLfP98V6ywBxCbDbp08d4VBvZwVVRkfw5Pv/e1+V77+ozw1DXo+pHbmD8jAinr2s0WSNlwdMl5EsouuYCu36Fr6AoNwKxcudLtWL9+/cyYMWPcjvXo0cOMHz/er3MuX77c3O/Sdj5r1izz1FNPeS1/7NgxU1JS4rzt3r3b7yY9kYBMmeL3FihVB4PbZ8u9Y5rR3AAmjngDr7mdI5Bp+wkJdds+or4I1qyrUE/PDyV/BgAHY7B1NExpl8atwXXJ+VJWVkZ+fj6DBw92Oz548GA2b97s1zkSEhLYvHkzx44dw2azkZeXxznnnOO1/IwZM2jVqpXzlpCQUKfXIOJN4VdWer6Wwag92W7H79qTTc/XMij8qrKZyWqFdwdkk0UGGWQxjXRgKGP4PefTijJKgWHA74Ff+P57KLgxm+03Z5KXB0uXQl6e/S/pG26AnTshN9feIpGba29JuPlme+vU8OH2f6u2clVls1Ht3JHWsWPdy9ls9pYlTy2pjmMPP1w/Xq8n3n6+rl1bgQ4Q98Rqheeft/+/aleg4+vnnqv5fSRSL4QhwAUVVVqYvv/+ewOYTZs2uZWbPn26Ofvss/0+78SJE02PHj1Mz549a1xWQC1MEg6Ov84nu7UY1bAFiocWqclkmeNgbqOrAU7cLjD3MdYYMBkW93WbgrHiszGhW026roLR6tEYproH8zVGw5T2UGpsy1dEkwa9DpO3wLR582a3ctOmTTPnnHNOWOqkdZgkFFwvWI6QdIw4t/Dk7YLl+ICePNn98bdxp4H2BjAngbmWG0PSpVTfu6vqupBfY1hMMdjdaY01NNTXPxzErlF1ybVr1w6r1cr+/fvdjhcXF9OhQ4cI1Uqk7ly7OqaRTilxxFNGKXEnutuql3OwWu1dZj17Vj4+nSxeZTFF/MxA4AjwFiuAu098ZWcMTDbZfHdXZq26lKKhu+qGG+wDtk8/3f14587+zboKRrdefRfs7jTHe9Lf7tyG4M034aabYM8e9+Pff28/HulZhhKYqA9McXFxJCUlsWbNGrfja9asoX///hGqlUjduV5sJ5PtDEvxlDGZbI/lfJ3DEbq6cZy3ieUyBmK/7i0E+gHbnc+VRQY/HbQSwOocThs3Vr9AuDIGdu+mVucOJn/G8XiTkmIPV96m6Fss7tPzo1Vdg2VjFg1/OEhgomKl78OHD7Njxw7n10VFRRQUFNCmTRu6dOnCo48+yogRI0hOTuaSSy5h/vz57Nq1izFjxkSw1iJ147go37XHHmDSTwzknkw22WRgARYlpPu8KDvO8f33MMlUhq7mlHEZA+lCV9ayiP0UAn35DVeSxSrnc/X0c+Cvq2AMFg4XR6tHbR73/PP2VgKLxf2i2NAGM99wg32JhI0b7T+zjh3t76uG8NpCKZA/HGrzHpQICEMXYZ3l5ua6DFatvI10WRp37ty5pmvXriYuLs4kJiaa9evXh61+GsMkobJ9WOUAb9cxEP7uDWeMfaxEupeB45PJMg/zR/Nrl9+rXvQ2cKjWg5Ybw4Boh8Y+mFm8awzj3BqCQK7fUbFwZX0X6oUrpRHLzKTwKytXvp/u9tdqQgL869Jsep5jg8xM3+fIzoaMDJ5pmcXjByvHPmVYsplq7C1XE8lmNsdJByqANrQltkMO339/gXtLQna2vQ/Bx3PabNCtm71Vy9Oni2PvL8eCh9HOZlPri1SXl2ffJ7EmublqYYqkQK7fCkxBoMAkoVani3JmJlit2Camu53jwAH49OZsBpLDQPIoJY6PKOMG4jlAKXExTfjzC3O55557sFgszuBFVhakp/t8SsdgV/DcXaXxL9LQNbY/HKKVAlOYKTBJtCocbl8A03V81B/IYNjJbVl3+D8ADBs2jL+ecQYtn3zSr7Dk8Oab9kGvVVvGnntOYUkaB/3hUP8pMIWZApNEpRMtRhWZWWy4PN3Z8nTZ+mzIzOBPZ57JhG++wQb8Clh+331cOG9eQE+h7ipp7PSHQ/2mwBRmCkwSlU501XlsMToxVumD6dMZVl7OLuxLeDz77LPcf//99i46EfGL/nCovxSYwkyBSRqkEy1QP8XGctfx4/zjxOEbb7yRBQsW0Lp164hWT0SkrgK5fkf9wpUiEgIuA7zblJWxaupUngVirVZWrFhBYmIiW7ZsiXQtRUTCRoFJRNx5mA1nycjgkawsNtlsdGvdmqKiIi699FKee+451EgtIo2BuuSCQF1y0qB4GtvkOAb898gRRn/9NW+e2Ajr2h49+Ptvf0ubp54K3vM5+LHuk4hIbalLTkRqLzOzenixWu2tTkDrmTN54403mDNnDnFWK299+SUX/vWvfPjhh7V7Pse5s7PdjztaujQ6VkTqgajYS05EIswRoE6EJkt6Og/89BOX2Gzc0qYN3/z0EykpKTz55JM89thjxMQE8LdYlXOTnh7QIpkiIuGgLrkgUJecNBqOIBMXB2VlkJXFwT/8gXvvvZdly5YBcPXVV7No0SLatWtX53MrLIlIKGlZgTBTYJJGJT7eHmji4qC0FABjDC+++CJjx46ltLSU008/nddee40BAwbU+dwiIqGiMUwiEhrZ2ZWBpqzMOe7IYrFw77338tFHH3HOOefw/fffk5qayowZM6ioqKjTuUVE6gMFJhHxj+u4otJS+79VBmtfcMEFfPzxx9xxxx3YbDYmTpzIkCFDKC4utg8m9xaCBg2q8dwiIpGkQd8iUjNPg7A9DdYGTj75ZBafcQZp113HA+++y3vvvUefPn1YctVVpC5cCHl59r0hHEsFDBoEOTmQllbjuUVEIkWBSURqZrN5HoTt+NpmcztsadKEu1atot+DD3JLTg6FhYUMeuklMrp3Z3JODs6FArKzK8PSunV+nVtEJBI06DsINOhbxIMTrVJHJk/moe+/Z+HChQCknXIKr/78M6dpNpyIRJhmyYWZApOIFy5LBbxcVsZ9cXEcKSujPfAqcIVmw4lIBGmWnIjUD+npzllvI+Li+HjbNs7v0IFiYDAwuayMcm17IiJRQIFJREKnylIBPR54gP/74Qd+n5yMAaYDaVOn8v3jj0e6piIiPikwiUjtBbJUQFoa5OTQLC2Nv2zZwtKlS2nRogUbgT7PPMP/jhgRzpqLiAREgUlEas/bxrmelgpISXGGJrKzGTZsGJ988gkXXnghB4Chr7zCE088wfHjx8P+MkREaqLAJCK1l55efZFJb0sFZGbav87Kci4V8Ktf/YoPPviABx98EIBZs2Zx+eWXs2vXrjC/EBER3zRLLgg0S04avSBsnLtixQpGjx5NSUkJp5xyCosWLeK3v/1tiCosIqJlBcJOgUmEoGyc++233zJs2DC2bNkCwCOPPMLMmTOJi4sLZk1FRAAtKyAi4RakjXPPOOMM3n//fR5++GEAZs+ezYABAygqKqos5GugeXZ25ZYrIiJBpMAkInXjx6a8gYiLi2P27NmsWrWK1q1bs2XLFi688ELefPNNewFvA80d9bBaq59URKSOtJeciNReAJvyBuraa6+loKCAYcOG8eGHH3LjjTfy4IMP8swzzxBf9fye6iEiEkRqYRKR2vO1Ka/LbLhaycyk6+LFbNiwgT/+8Y8AzJkzh/79+7PjP/+B1FR7SIqPV1gSkZBrdIHp+uuv55RTTuGmm24K6D4R8SAz03tISU+v23iiE11vsTNn8tRTT/HOO+/Qtm1bPvnkExKff57l7dtXjpmKi1NYEpGQanSBaezYsSxevDjg+0QkzKqs8TR06FAKRo1iAHAIuHX5cu4rK+OX2Ng6DTQXEfFHowtMAwcOpEWLFgHfJyIR4Bqa4uPp/Kc/kTtlChMvuwwL8Bfg4nPP5auxY+s00FxEpCb1KjBt2LCBa665hk6dOmGxWFi1alW1MvPmzaN79+40bdqUpKQkNm7cGIGaikjYpKe7db01sVqZvmED/7rzTk499VQ+/fRTkl54gVd79vQ+e05LDYhIHdWrwHTkyBF69+7NnDlzPN6/bNkyHn74YSZNmsTWrVtJSUlhyJAhbtsoJCUl0atXr2q3vXv3Bq2epaWlHDx40O0mIiFSdY2nnBzIymLwSy+xbds2UlNTOXL8OHcUFvK7jh05euyY+2O11ICIBIOppwCzcuVKt2P9+vUzY8aMcTvWo0cPM378+IDOnZuba2688caA73OYMmWKAardSkpKAqqHiNQgK8sYsP/r6WtjTHl5uZkyZYqxWCwGMOe1b2+2b9/usayIiKuSkhK/r9/1qoXJl7KyMvLz8xk8eLDb8cGDB7N58+aw1mXChAmUlJQ4b7t37w7r84s0Ct7WeKqyMKbVaiUzM5O1a9dy2skns724mL7nncciLTUgIkEUNQtXHjhwAJvNRocOHdyOd+jQgf379/t9niuvvJJPPvmEI0eO0LlzZ1auXEnfvn1rvM9VfHw88fHxdXtBIuKbrzWeHPe7SEtLo2DHDu7o2JG1xnAXkLtjB3MPH+bkk092P0dmpr2bzlOYys62n1vjnkTERdQEJgeLxeL2tTGm2jFf3n333VrdJyJh5iuweGk16jB/Pu8awwyrlQybjcWLF/PR6tUsz8nh/PPPryzo2F4lJwdycyuPu7ZqiYi4iJouuXbt2mG1Wqu1JhUXF1drdRKRRuhE2InJymJSeTm5d99NJ+DLAwfol5jIiy++iDHG/TF5eZWz6rS9ioj4EDWBKS4ujqSkJNasWeN2fM2aNfTv3z9CtRKResFD2Lnsb3+j4IknGAIcKy/n3nvv5fbbb+fQ5MmVZV3WeFJYEhFf6lWX3OHDh9mxY4fz66KiIgoKCmjTpg1dunTh0UcfZcSIESQnJ3PJJZcwf/58du3axZgxYyJYaxGJOC/jnU6dOZO3mzXjmcWLmfjttyxdupQtwPL77uNCR9lp07S9iojUyGKqtVFHTl5eHgMHDqx2fOTIkSxatAiwL1w5a9Ys9u3bR69evZg9ezaXXXZZmGvq7uDBg7Rq1YqSkhJatmwZ0bqIiGebY2MZVl7Obuwt1rNnz+a+AwewTJlSucaTWphEGpVArt/1KjBFKwUmkXruRJfdT7Gx3HX8OP84cfgmYMHEibSaPl1jmEQaoUCu31EzhklEpFZcglCbsjJWTZ3KbCAWeAO4cOlSPv74Y49rPImIOCgwiUjD5aHVyJKRwcOpqWwCurVuTVFREf379+f555/HTJ5sL1tljScREXXJBYG65ETqqRoWqPzvkSOM/vpr3nzzTQCuu+46/v73v3PKKaeEt54iEhEawxRmCkwi0csYw9y5c3nssccoKyuja9euvPbaa1x88cWRrpqIhFjIxzA99dRTAHz66accP368NqcQEakXLBYLDz74IB+MGsWZbdrw3XffkZKSwjPPPENFRYW9UHa2tkoRaeRqtQ7TgAEDAMjMzOSLL74gNjaWXr16cf7553P++efTt29frb4tIlElsXNnPvnpJ+7t1Ytln3/O448/Tl5eHi+dfz5tZ87UdikijVxQuuSOHj3K559/zmeffcZnn33Gxo0bGTp0KNmNZKaJuuREGojsbExGBvOvuYY/vPcepaWldAaWjh7NgAULIl07EQmykI9heuqpp3jiiSf49NNPOffcc4mNja1WJikpifz8/EBPHZUUmEQakBMz67Y1acIt5eX8G7DGxJA9bRpPPPEEMTEx7mVtNnXXiUSpkIxhstlsrFixgkOHDrl1yV1wwQVccMEF3HbbbcyYMYO3336bH374gQ8//LBur0JEJBLS0yEujt7l5eTHxnJH797YKiqYOHEiQ4YMobi42F7OsWSB1RrZ+opIWATUwtSsWTO2b9/OGWec4XZcXXJqYRJpMBxB6MR2KWbqVBZu3cqDq1bxC9CxY0eWXHUVqQsXalVwkSgXslly/fr1o6ioqNrx5s2b069fP0aPHs1zzz1Hfn4+q1evDqzWIiKR5rrQZWkpZGVhmTKFuxMT2fLgg5wL7Nu3j0ELF5I1cCC2iRMjXWMRCZOAAtPYsWOZOHEiu3fvrrGsuuREJKp42kvOZbuU89q3Z0tsLHcBFcCU3FwGn3UW+/fv93wujWsSaVACCkw333wzW7Zs4bzzzuOOO+5gwYIF5OfnU1ZWVq2sp4HgIiL1ls3muYvNEZpycjjp+HH+HhfHYqB5TAw5RUX0/tWvWLt2bWV5jW0SaZACGsP03XffUVBQwLZt25z/7ty5E6vVSo8ePfj0009DWdd6S2OYRBq4qq1P2dl8mZHBLSedxGdHjmCxWJg4cSKZTZrQZOpUjW0SiRJh3Rrl0KFDFBQU8Omnn/LAAw/U5VRRS4FJpAHz1FV34vgvGRk83KkT8/fuBeAyYMm4cZz+9NORqauIBER7yYWZApNIA1bDBr7YbLw2bRr32GwcBtq1a8fLL7/MVVddFdA5NOZJJPxCvpeciEijkZnpvXstPR2sVobZbHwSG0sf4MCBAwwZMoQJEyZU7rVptdpbqaoutaLxTiJRQ4FJRKS2XLrrzior44P0dBwDE2bOnElqaqp9VrHLbDtnaPLW1Sci9ZK65IJAXXIijZCPsU1vZGQwOj6eg6WltGnThkWLFnHNNddUWxRTYUkkstQlJyISat6WIbDZuCktja2//z3Jycn89NNP/Pa3v+Wxxx6jrLzc3v1WVmYPTQpLIlFDgUlEpDa8jW2yWiEnhzPatWPTpk08/PDDADz77LOkZGWx02arbGFqJNtHiTQETSJdARGRBsURojIyiANmz55N6g8/MGrpUj4CLmzalL8vWcL1n39u755zfYyI1FtqYRIRCTbXQd7x8Vy7dCkFwMWdO/PfY8e44YYbGPvjj5Redpnn2XOg7VVE6hkFJhGRUEhPr+x6i4mha1YWG779lscffxyAP//5z1y6YwffAOTkuD9Wyw2I1DsKTCIioZCdXTm4u6ICsO+xOWvWLN5++23atm1L/t69JMbH83peXs3LDWRmeh/zpNYokZBTYBIRCTbX0FNaWm0NpquvvpqCggIGDBjAwdJSbgHuz8jgWFyc97WZtPilSGQZqbOSkhIDmJKSkkhXRUQiLSvLGLD/W8Px48ePmwkTJhjAAKY3mK9iY/0/t7fnEhG/BHL91iw5EZFg8rY+k+Nrm815qEmTJjz55JNcvmcPI15+mW1A4vHj/PWmm7j9jTeqn9tlBh7TpmnxS5EwanQrfV9//fXk5eUxaNAg3vDwgXT06FHOPfdcbr75Zp555hm/zqmVvkWk1k50qe0dN47bP/6YvLw8AEYnJvI/GzfSvHnz6o+Jj68cH1VaGt76ijQgWunbh7Fjx7J48WKv90+fPp2LLroojDUSkUbLZaxTp6efZu3atWRkZGCxWPjbJ59wUffufPHFF9Uf4whLWvxSJGwaXWAaOHAgLVq08Hjf119/zZdffsnQoUPDXCsRaZSqdN9ZrVamTp3K2hEj6BAXx+fFxSQnJ/PSSy/ZyzsCVmqqx8HkIhI69SowbdiwgWuuuYZOnTphsVhYtWpVtTLz5s2je/fuNG3alKSkJDZu3Bi05x83bhwzZswI2vlERHzysr1K2q9+xbayMq444wyOHj3KqFGjGHnhhRxxrAyelmb/13WBTIUmkZCqV4HpyJEj9O7dmzlz5ni8f9myZTz88MNMmjSJrVu3kpKSwpAhQ9i1a5ezTFJSEr169ap227t3r8/nfuuttzj77LM5++yza6xnaWkpBw8edLuJiARNejodsrL417ffMm3QIGJiYlhcUEAy8NkDD7iHLEdochlMLiIhEPI5e7UEmJUrV7od69evnxkzZozbsR49epjx48cHdO7c3Fxz4403uh0bP3686dy5s+natatp27atadmypZk6darHx0+ZMsU5Ddj1pmUFRCSoTiwbsL5JE9PpxOdM06ZNzYsvvmgqKioiXTuRqBfIsgL1qoXJl7KyMvLz8xk8eLDb8cGDB7N58+Y6n3/GjBns3r2bnTt38swzz3DPPfeQ4Wj+rmLChAmUlJQ4b7t3767z84uIVHNie5XLysspiI3lqquu4tixY9xzzz3cfvvtHDp0KNI1FGk0oiYwHThwAJvNRocOHdyOd+jQgf379/t9niuvvJKbb76Z1atX07lzZ7Zs2RJwXeLj42nZsqXbTUQk6FxmxJ16/DjvXHwxM2fOxGq1snTpUpKSkigoKIh0LUUahagJTA4Wi8Xta2NMtWO+vPvuu/z4448cPXqUPXv20Ldv32plRo0a5fcaTCIiIeFhe5WYzEyeKCtjw4YNJCQk8PXXX3PxxRczb948TONaUk8k7KImMLVr1w6r1VqtNam4uLhaq5OISFTztAGvy4y4/uvWsXXrVq655hpKS0t54IEHuOWWWygpKYlsvUUasKgJTHFxcSQlJbFmzRq342vWrKF///4RqpWISAj42l7lxIy4tm3b8tZbb/GnP/2JJk2a8MYbb5CYmMjHH39sL5uZ6X2pgexs+/0i4rd6tZfc4cOH2bFjh/ProqIiCgoKaNOmDV26dOHRRx9lxIgRJCcnc8kllzB//nx27drFmDFjIlhrEZEg8xVmXEKUxWLh0Ucf5dJLL+XWW2/l22+/pX///jz99NOMjYnB4pi44hq8XFuvRMR/IZ+zF4Dc3FyP0/VHjhzpLDN37lzTtWtXExcXZxITE8369esjV+ETApmWKCISVFOmGJOVZX766Sdz/fXXOz83r7vuOvNTSooxYF+ewBjnMgXOr0UauUCu341u891Q0Oa7IhIxLi1GZvJk5syZw7hx4ygrK6MrsCwpiYvy8yv3nvPU1SfSSAVy/a5XXXIiIhIgR/jJyMACPJSeTv/CQm75y1/4FhiwbRszrVYeKSsjJi5OYUmklqJm0LeIiHjhuqdcfDxJf/kLn0ycyM0330x5eTnjbDZ+a7Hwn7Iy7TknUksKTCIiDcGJVcEdC122mj6dZb168QIQ36QJ7xhDn5Yt2aSNekVqRYFJRKQhcFkVnLIyGDQIy5QpjMnK4sMtWzjrrLPYc/Agl8fEMDMjg4qpUysfqyUIRGqkwCQiEu08rApOTg6kpUF6On369CE/P5/bbrsNW0UFE4CrX3mFH3/80f54q9X++KqhyXFeqzXsL0mkvlFgEhGJZr5WBc/JcYagFi1a8Morr7BgwQKaNm3Kv3bsoE+fPqxfv959DJQjNHk6r0gjpllyIiLRzNeq4I77T7BYLIwePZqLLrqIW265hS+++IK0tDQyMzOZOHEiVrCHpGnTtASBSBVahykItA6TiESbI0eO8OCDD7Jo0SIABg0axCuvvMJpXbtWjoUqLY1sJUVCLJDrt7rkREQaoZNOOomFCxfy0ksv0bx5c9atW0efs85inevAcc2mE3FSYBIRacTup1OvUQAAIABJREFUvPNO8vPz6dW+PT8cPsyvLRYynniC8ilTPA8EF2mkFJhERBq5Hq+/zkfFxdyTlIQxhuzsbAbl5rJ33DiFJpETFJhERBo7m41mWVnM//hjlixZwsknn8yGDRvovWgR/xoxwm3guEhjpcAkItLYZWY6Z8MNHz6cTz75hD59+nDgwAGGvPwyE0pLKS8vdy8/aJDnBS210KU0UFpWIBocOeL9PqsVmjb1r2xMDDRrVruyR4+CtwmVFgs0b167sr/8AhUV3utx0km1K3vsmO+/igMp27y5vd5gnzXkeuGoS9lmzezfZ7APsD1+PDhlmzatXGgwkLLHj9vLexMfD02aBF62vNz3bKu4OIiNDbyszWb/2XkTG2svH2jZigr7ey0YZZs0sX8vwP47cfRocMoG8ntfi8+Is846iw8++IAJAwYwPz+f/5k5ky15eSxatIjOnTtDXh6sX+/+2KNHYcYM+5IEkye7P48+I2pXVp8Rdq6/95FkpM5KSkoMYEpKSkLzBPaPFs+3oUPdyzZv7r3s5Ze7l23XznvZ5GT3sl27ei/bs6d72Z49vZft2tW9bHKy97Lt2rmXvfxy72WbN3cvO3So7++bq5tu8l328OHKsiNH+i5bXFxZ9v77fZctKqosO26c77Kff15ZdsoU32U/+qiy7KxZvsvm5laWnTPHd9m3364su3Ch77LLl1eWXb7cd9mFCyvLvv2277Jz5lSWzc31XXbWrMqyH33ku+yUKZVlP//cd9lx4yrLFhX5Lnv//ZVli4t9lx05srLs4cO+y950k3Hjq2zVz4jYWO9lu3Vz/174+oxo0cL+b1aWvWzr1t7LBvIZ0bq1e1l9RtjpMyIkArl+q0tORETsdu70exuUirPPrlwdPD4e/vtf74UdW7D447//9d7d54m6ACVMtHBlEIR84Up1yQVeVs3tgZetD83t6pLzr2xdu+RmzqzsOhs/vvJrxxYpHs5b2ro1k8vLmXfi6+S+fXlp+XK6nXOO/b0QGws//1z9+R3ndl01/OhRGDrU3q13+eXwzjuV5a++2n48LQ3WrXOfpeeor7dz6zPCv7LR/BkRZAFdv0PWztWIhLxLTkQk2LKy7N0dcXHuXWs1lF0JpnXTpoYT/6705xyOxzvud3ydlub7uOt5vZ3DV71FahDI9VuBKQgUmEQkKjkCSVyc9zIegkoRmItatjSAAczYsWPNsYwM/0JT1XDl7binulUtm5rq+/lcx2OJeKDAFGYKTCISdfxpYfLWipOWZkrBPJaQ4AxNSR07mh0XX+y5vCO8eAtoVY/7qptrWW/1U+uT+EmDvkVExLvsbPvYoKws+9gRx+Dtqit622zuY48cUlKIS0vjmbvv5p///Cdt2rQhf98+Ej/8kNd79XIf7+N4ro0bKzf1dd2nLjvb/figQd7rVrUsVK+74/nS0qrX23G/BolLbYQhwDV4amESkagRglaZXbt2mUsvvdTZ2nRf377ml19+8X+skuO44+u0NM918zaGqWqLVNXzBuE1SsMUyPVbC1eGkc1m47iv2QhSb8TGxmL1c3q1SFTx1mrk+LoW26AkJCSQm5tLRkYGM2fO5IUtW/igeXOWGcPZaWmQk+P+nOnp9sUvc3LcW4JSUuz/5uTYW4Kq1jE11f0cUNka5Wh1iouzz65ztDQ5yrq2qrmeNzPTPgvMW2uUzaYWKQG00ndYGGPYv38///W1TonUO61bt+a0007D4pgCLNIQ+Lr4ewoNfoqNjWXGjBlcfvnljBgyhAJjSAL+2ro1t3np1nP717VujqDiUFPIy8mp3t3nGqimTbMf93QOq9U9WDm4BizX+nkKV5mZ9i7HlJTq31+FrgZDgSkMHGGpffv2NG/eXBfges4Yw9GjRykuLgagY8eOEa6RSPS4assWtgG3WSysN4bb33yT3DZteP7oUZq7rsEWSHCrKWzk5VWGoaotS46wFONlyK6jtcuf1ihv4WrjRntoq8pT6HK8HrVqRZ+QdxA2Ar76QMvLy01hYaE5cOBABGomdXHgwAFTWFhoysvLI10VkejgMkbo+PHjJv3yy43lxLimXr16mcLCwpA9n8fjVddzCqSstyULHOVTUz0/3p91ojzdN2WK57FXjuOelkjQ0gl1pmUFwszXN/yXX34xhYWF5ujRoxGomdTF0aNHTWFhoX3wqoj45iUgrBk50nQ4EZqaN29uFi1aFLznnDLF+wBubwHGV6gJZMmCQNaU8qTq+b0NeK9pIHzVcKV1qQKiwBRm/gQmXXSjj352IgHwcbHe9/jjZlD37s5ZdCNHjjSHXTesDbZghR1fq4v7s6ZUTQEmNdX37L5AW69CuS5VAw1jCkxhpsDUMOlnJxI85eXlJjs728TExBjAnHvuueazzz4LzZPVdHG3Wj0vlOnPkgXejrk+LpDlDfxdtNPfVc69dRl6KhtICGqgi4QqMPlw3XXXmdatW5sbb7zR7fiXX35pevfu7bw1bdrUrFy50q9zNrbAVFRUZACzdevWkD/XwoULTatWrYJ2vtzcXAOYn3/+ucayDfFnJxJpeXl5plOnTgYwTZs2NQsWLDAVFRXhq0BdQ423cFWb1iFv4cif1qtAWtH8DVeeXruv41EeloxRYPIpJyfH/OMf/6gWmFwdOnTItG3b1u8m43AFpvJyY3JzjVmyxP5vpMYihzMwHT161Pzwww9BO58Ck0jkFRcXm6uuusrZRXf77bebgwcPhv6JfYUdb+Wrdpt5Cgg1jTMKZDNhT2Vd769Nl6E/4cpxn+P1egtYVb8f3sJSIK1XEezuU2CqQW5urs/A9Oqrr5pbbrnF7/OFIzCtWGFM586V722wf71iRZ1OWyvhDEzBpsAkUj/YbDYzc+ZMY7VaDWDOPvtsU1BQELonrE2XUtX7vIUJf2ay+RNgqgav2rZSBRquqpYPtBuxrt/vCHb3RW1gWr9+vfnNb35jOnbsaACPXWJz58413bp1M/Hx8SYxMdFs2LAh4OepKTBde+21ZkUASSTUgWnFCmMsFvewBPZjFkvoQpPjA+3MM880cXFxJiEhwUybNq1aYCovLzd333236datm2natKk5++yzzXPPPed2rtzcXNO3b1/TvHlz06pVK9O/f3+zc+dOY4wxBQUFJjU11Zx88smmRYsWJjEx0WzZssUY47lL7q233jJJSUkmPj7etG3b1lx//fXO+15++WWTlJRkTj75ZNOhQwczfPhwtxYqBSaR+mXTpk0m4cQmvvHx8eaFF14ITRddoK0YwbyI+zP+yDUcudbFn1ly/nQZ+uoCNKb25/DW3edad39ar3x1aXr7ngah5SlqA9Pq1avNpEmTzIoVKzwGptdee83ExsaaF1980RQWFpo//OEP5qSTTjLfffeds0xiYqI577zzqt2+//57ZxlfgamkpMS0a9fO50Xy2LFjpqSkxHnbvXt3yAJTeXn1lqWqoSkhITTdc3/84x/NKaecYhYtWmR27NhhNm7caF588cVqgamsrMxkZGSYjz76yHz77bfmlVdeMc2bNzfLli0zxhhz/Phx06pVKzNu3DizY8cOU1hYaBYtWuT8uZ133nnmjjvuMF988YX597//bZYvX+78S7NqYHr77beN1Wo1GRkZprCw0BQUFJjp06c77//b3/5mVq9ebb755hvzwQcfmIsvvtgMGTLEef//t3fncVFX+//AX8OwL4JIIqggZqICMoioKBGioeTuLbW8Ct+0e9VAicott4smltfUbqa5ZK5X6oua3TQ1E3CNdXAjccEd5KqxKKtwfn/4c76O7OtnZnw9H4/PQ+Zzznw+78NnxnlzzvmcYcJEpHnu3bsnhg0bphqiGzNmjMjJyZE2qMYaJqrtvJ+qzldT71Vthgxr6h1qaC9VdedsyFpYzfCdgFqbMD2rsoSpV69eYsqUKWr7unTpImbPnl2nY1eXMG3dulWMHz++2ucvXLhQ9cZ+dmuKhOno0aqTpWe3o0frdfgq5eXlCSMjI7Fhw4YKZbUZkps2bZrqd3z//n0BQMTExFRa18LCosq1WZ5PmLy9vWu8Ps+Kj48XAER+fr4QggkTkaYqLy8XK1asEPr6+gKA6Nixo0hMTJQ6rIZp6qGm2g4ZPq1b2fyjppgH9Xzd2vZePV+3svJGHqbTyYSpuLhYyOVysXv3brV606dPF76+vnU6dnUJ09ChQ8W+ffuqfX5z9jDt3Fm7hGnnznodvkq///67ACCuXr1aoayyhGnt2rXC09NT2NjYCDMzM2FgYCC8vLxU5cHBwcLIyEgMHTpUrFq1Sty5c0dVtnDhQqGvry8GDBggIiMjxeXLl1VlzydMJiYm4ttvv60y7uTkZDF8+HDh4OAgzM3NhampqQAgzp8/L4RgwkSk6U6fPi0cHR0FAGFgYCBWr17dvHfRNaamnMxc12SsPr1DtR16qyoZe/ZY9VmWoab9jUAnE6bbt28LAOLEiRNq9T799FPRuXPnWh83ICBA2NjYCBMTE9G2bVsRHx+vKsvJyRGtW7cWxcXFdYq1KecwSdXDdObMmVonTFFRUcLY2FisWbNGJCcni0uXLom//e1vwt3dXe15ycnJYunSpcLb21uYm5uLU6dOqcouXrwovvjiC/H6668LQ0NDVWL8fMJkbW1dZcL08OFDYWNjI9555x0RFxcn0tLSxMGDB9ViZcJEpPkePHggRo0apeq9HzVqlHjw4IHUYWmWxlhDqaoeqaf1a7OW07Pq0jtU156k2kw0rwedTphOnjypVm/JkiXC2dm5ucNT05QJ09M5TJVN+kYTzmEqLCwUJiYmtRqSCwkJEf7PTcwbMGBAhYTpWX369BGhoaGVlo0bN04MGzZMCFExYfLz86tySC4xMVEAEDdu3FDt27ZtGxMmIi1UXl4uvvzyS2FoaCgACEdHR3H69Gmpw9JOzbFA5fO9QDXNP6pL3aqWWmgEdUmYqvj6Zs1jY2MDuVyOrKwstf3Z2dmwtbWVKKqmJ5cDq1c/+VkmUy97+njVqif1GpOxsTFmzZqFmTNnYuvWrbhy5QpOnz6NTZs2VajbqVMnJCYm4uDBg0hPT8f8+fORkJCgKs/IyMCcOXNw6tQpXL9+HYcOHUJ6ejq6du2KwsJChISEICYmBtevX8eJEyeQkJCArl27VhrXwoUL8e9//xsLFy5EWloazp49i88//xwA4ODgAENDQ/zrX//C1atXsW/fPixevLhxfzFE1CxkMhlCQ0Nx8uRJdOzYEdevX4ePjw9WrFgBIYTU4WmXRYuA+fMrL5s//0n5U2VlQERExfrz5z/ZX1ZW8RiLFwMLFjwpLy5+8u9vvwH+/lWf18+v5rrz5z/Z/9tv6sdesODJOZtbo6VpjQyofNL31KlT1fZ17dq1zpO+G5tU6zC1b9+06zCVlZWJJUuWCEdHR2FgYCAcHBzE0qVLK/QwFRUVieDgYGFpaSmsrKzE1KlTxezZs1U9TFlZWWLkyJHCzs5OGBoaCkdHR7FgwQJRVlYmiouLxbhx40T79u2FoaGhsLe3FyEhIarfV2XLCkRHRwuFQiEMDQ2FjY2NGD16tKps586dqmUnvL29xb59+9jDRKTlcnJyxFtvvaUaohs6dKi4d++e1GGREHXrkWrovKua9tdDXXqYZEJoTqr+8OFDXL58GQDg4eGBL774Av3794e1tTUcHBwQFRWFCRMmYN26dfD29sb69euxYcMGnD9/Ho6OjpLFnZeXB0tLS+Tm5qJFixZqZUVFRcjIyICTkxOMjY0bdJ6yMuDYMSAzE7CzA159tfF7luj/NOa1I6KGEULgm2++QVhYGIqLi9GuXTvs2rUL/fr1kzq0F9uiRU8+iCrrSVq8+MkH19MerLrUrU/9eqju8/t5GpUwxcTEoH///hX2BwUF4bvvvgMAfP311/j888+RmZkJV1dXrFy5Er6+vs0cqbrmSpioefHaEWkepVKJMWPG4NKlS5DL5ViyZAlmzpwJPT2tmWFCGkRrEyZtxYRJN/HaEWmm/Px8TJkyBTt37gQADB48GFu3bsVLL70kcWSkbeqSMDElJyIirWJhYYHt27djw4YNMDY2xi+//AKFQoHY2FipQyMdxoSJiIi0jkwmw+TJkxEfH48uXbrgzp078Pf3x+LFi1FW2Z1cRA3EhImIiLSWm5sbEhMTMXHiRJSXl2PBggUYNGhQhSVoiBqKCRMREWk1MzMzbNmyBd999x1MTU1x5MgRKBQKHDlyROrQSIcwYSIiIp0QFBSEhIQEuLi44O7du3j99dexcOFCDtFRo2DCREREOqNbt26Ij4/H5MmTIYRAREQEBgwYgDt37kgdGmk5JkxERKRTTE1NsWHDBmzfvh3m5uaIjY2FQqHAwYMHpQ6NtBgTJpKcn58fwsLCVI87dOiAVatWqR5nZWXh9ddfh5mZGaysrAA8uUNm7969DTpvcHAwRo4c2aBjEJHmGj9+PJKSkuDu7o7//ve/GDx4MObOnYvHjx9LHRppISZMpHESEhLwt7/9TfV45cqVyMzMhFKpRHp6OgAgMzMTgYGBUoVIRFqic+fOOH36NKZOnQoAiIyMRP/+/XHr1i2JIyNtw4RJGyxaVPU3My9e3ODv0tE0L730EkxNTVWPr1y5Ak9PT7zyyito3bo1AKBNmzYwMjKSKkQi0iLGxsb4+uuvERUVBQsLCxw/fhwKhQI///yz1KGRFmHCpA3kcmDBgopJ0+LFT/Y30Tfw+vn5ITQ0FGFhYWjZsiVsbW2xfv16PHr0CP/zP/8DCwsLvPzyyzhw4IDqObGxsejVqxeMjIxgZ2eH2bNnq3V/P3r0CBMnToS5uTns7OywYsWKCud9dkiuQ4cOiI6OxtatWyGTyRAcHAyg4pDc7du3MXbsWLRs2RKtWrXCiBEjcO3aNVV5WVkZwsPDYWVlhVatWmHmzJngtwIRvVjGjBmDlJQUeHp64v79+xg6dCg+/vhjlJaWSh0aaQEmTNpg/nwgIkI9aXqaLEVEVP5Nzo1ky5YtsLGxQXx8PEJDQzF16lS89dZb6Nu3L5KTkzFo0CBMmDABBQUFuH37Nt544w14eXkhNTUVa9euxaZNm7BkyRLV8T7++GMcPXoUe/bswaFDhxATE4OkpKQqz5+QkIDBgwdjzJgxyMzMxOrVqyvUKSgoQP/+/WFubo64uDgcP34c5ubmGDx4MEpKSgAAK1aswLfffotNmzbh+PHjePDgAfbs2dP4vzAi0mgvv/wyTpw4genTpwMA/vnPf8LX1xfXr1+XODLSeIIaLDc3VwAQubm5FcoKCwvFhQsXRGFhYcNPFBEhBCCEoeGTfyMiGn7Marz22mvCx8dH9fjx48fCzMxMTJgwQbUvMzNTABCnTp0Sc+fOFc7OzqK8vFxVvmbNGmFubi7KyspEfn6+MDQ0FLt27VKV379/X5iYmIgZM2ao9jk6OoqVK1eqHo8YMUIEBQWpxQZA7NmzRwghxKZNmyqct7i4WJiYmIiDBw8KIYSws7MTy5YtU5WXlpaKdu3aiREjRlTZ/ka9dkSkcXbv3i2srKwEAGFlZSX27t0rdUjUzKr7/H4ee5i0yfz5gKEhUFLy5N8m7Fl6qnv37qqf5XI5WrVqBTc3N9U+W1tbAEB2djbS0tLg7e0NmUymKu/Xrx8ePnyIW7du4cqVKygpKYG3t7eq3NraGs7Ozg2KMSkpCZcvX4aFhQXMzc1hbm4Oa2trFBUV4cqVK8jNzUVmZqbaefX19dGzZ88GnZeItNuoUaOQkpKCXr16IScnByNHjkRYWJiqZ5roWUyYtMnixf+XLJWUVD0RvBEZGBioPZbJZGr7niZH5eXlEEKoJUsAVPOEZDJZk80ZKi8vh6enJ5RKpdqWnp6Od955p0nOSUS6oUOHDjh27Bg+/PBDAMDq1avRr18/XL16VeLISNMwYdIWz85ZKi6uOKdJA3Tr1g0nT55US4xOnjwJCwsLtG3bFp06dYKBgQFOnz6tKv/zzz9VSwXUV48ePXDp0iW0bt0anTp1UtssLS1haWkJOzs7tfM+fvy42rlTRPTiMDQ0xD//+U/89NNPsLa2RmJiIjw8PPC///u/UodGGoQJkzaobIJ3ZRPBJTZt2jTcvHkToaGh+OOPP/Djjz9i4cKFCA8Ph56eHszNzTFp0iR8/PHHOHLkCM6dO4fg4GDo6TXsZTh+/HjY2NhgxIgROHbsGDIyMhAbG4sZM2ao1lqZMWMGli1bhj179uCPP/7AtGnTkJOT0xjNJiIdMXToUCiVSvTr1w95eXl46623MG3aNBQVFUkdGmkAJkzaoKys8rvhniZNGvLFkm3btsX+/fsRHx8Pd3d3TJkyBZMmTcK8efNUdZYvXw5fX18MHz4cAwcOhI+PDzw9PRt0XlNTU8TFxcHBwQGjR49G165d8e6776KwsBAtWrQAAHz44YeYOHEigoOD4e3tDQsLC4waNapB5yUi3dO+fXscPXoUs2fPBgCsXbsWffr0aXBPOGk/mWiqiSUvkLy8PFhaWiI3N1f1Af1UUVERMjIy4OTkBGNjY4kipPrgtSN6sf3yyy+YMGEC7t27B3Nzc3zzzTecF6ljqvv8fh57mIiIiCoxePBgpKamwtfXFw8fPsT48ePx3nvvoaCgQOrQSAJMmIiIiKpgb2+PI0eOYP78+ZDJZNi4cSN69+6NtLQ0qUOjZsaEiYiIqBr6+vqIiIjAoUOHYGtri3PnzqFnz57YunWr1KFRM2LCREREVAsDBw6EUqmEv78/CgoKEBQUhODgYDx69Ejq0KgZMGEiIiKqpTZt2uDQoUOIiIiAnp4etmzZAi8vL5w7d07q0KiJMWEiIiKqA7lcjvnz5+PIkSOws7NDWloaevXqhU2bNjXZNxqQ9JgwERER1YOfnx+USiUCAgJQWFiIyZMnY8KECcjPz5c6NGoCTJiIiIjqqXXr1jhw4ACWLl0KuVyOHTt2oGfPnkhNTZU6NGpkTJiIiIgaQE9PD3PmzEFMTAzatWuH9PR09O7dG+vWreMQnQ55oRKmmzdvws/PD926dUP37t3xww8/qJX/5z//gbOzM1555RVs3LhRoiiJiEgb+fj4ICUlBUOGDEFxcTGmTp2KcePGIS8vT+rQqBG8UAmTvr4+Vq1ahQsXLuDXX3/FBx98oLod9PHjxwgPD8dvv/2G5ORkfPbZZ3jw4IHEEeuG48ePIzIyUuowiIianI2NDfbt24fly5dDX18f33//PXr06IGkpCSpQ6MGeqESJjs7OygUCgBPxp2tra1VSVF8fDxcXFzQtm1bWFhY4I033sDBgwelDFdn+Pj4YM6cOVKHQUTULPT09PDRRx/h2LFjcHBwwJUrV9C3b1989dVXHKLTYhqVMMXFxWHYsGGwt7eHTCbD3r17K9T5+uuvVV+G6unpiWPHjtXrXImJiSgvL0f79u0BAHfu3EHbtm1V5e3atcPt27fr1xBSM2rUKCiVSqnDICJqVn369EFKSgpGjBiBkpIShIaG4s0330ROTo7UoVE9aFTC9OjRI7i7u+Orr76qtDwqKgphYWH45JNPkJKSgldffRWBgYG4ceOGqo6npydcXV0rbHfu3FHVuX//PiZOnIj169er9lWW9ctkskrjKC4uRl5entqmqyIjI+Hl5QULCwu0bt0aI0eOxMWLF+t0jLS0NHTt2rWJIiQi0lzW1tbYs2cPVq1aBQMDA+zevRseHh6Ij4+XOjSqI32pA3hWYGAgAgMDqyz/4osvMGnSJEyePBkAsGrVKhw8eBBr165VzZGpaZy4uLgYo0aNwpw5c9C3b1/V/rZt26r1KN26dQu9e/eu9BiRkZH4xz/+Uet2abPY2Fi8//778PLywuPHj/HJJ58gICAAFy5cgJmZGTw9PVFcXFzheYcOHYK9vT0ePXoEfX19GBkZSRA9EZH0ZDIZZsyYgX79+mHs2LG4evUqfHx88NlnnyEsLKzKP85JwwgNBUDs2bNH9bi4uFjI5XKxe/dutXrTp08Xvr6+tTpmeXm5GDdunFi4cGGFstLSUtGpUydx69YtkZeXJzp16iTu3btX6XGKiopEbm6uart586YAIHJzcyvULSwsFBcuXBCFhYW1ilHTZWdnCwAiNja2VvVPnz4txo4d28RRNQ1du3ZEJL2cnBzx5ptvCgACgBg2bJi4f/++1GG9sHJzc6v8/H6eRvUwVefevXsoKyuDra2t2n5bW1tkZWXV6hgnTpxAVFQUunfvrpoftW3bNri5uUFfXx8rVqxA//79UV5ejpkzZ6JVq1aVHsfIyKhBPSZCCBQUFNT7+Q1hamraoL9mcnNzATzpZq6Ns2fPws3Nrd7nIyLSJZaWlvj++++xbt06fPDBB/jpp5+gUCiwa9cutVEP0jxakzA99fyHvRCi1gmAj48PysvLqywfPnw4hg8f3qD4aqOgoADm5uZNfp7KPHz4EGZmZvV6rhAC4eHh8PHxgaura62ec/bsWQwYMKBe5yMi0kUymQxTp06Ft7c3xowZg0uXLsHX1xdLly7FRx99BD09jZpeTP+f1lwVGxsbyOXyCr1J2dnZFXqdqGmEhITgzJkz+Pe//13r5xw7dgxeXl5NGBURkXZSKBRISkrCO++8g7KyMsyaNQtDhw7Ff//7X6lDo0poTQ+ToaEhPD09cfjwYYwaNUq1//DhwxgxYoSEkdWdqakpHj58KNm56yM0NBT79u1DXFwc2rVrV2P9kpIS9OnTB4MHD4adnV29zklEpOssLCywfft29O/fH6GhoThw4IBqiO7VV1+VOjx6hkYlTA8fPsTly5dVjzMyMqBUKmFtbQ0HBweEh4djwoQJ6NmzJ7y9vbF+/XrcuHEDU6ZMkTDqupPJZPUeFmtuQgiEhoZiz549iImJgZOTU62eZ2hoiOTk5CaOjohI+8lkMkyePBm9e/fGmDFj8Mcff8DPzw8RERGYM2cOh+g0hEZdhcRDEX1dAAAdDElEQVTERHh4eMDDwwMAEB4eDg8PDyxYsAAAMHbsWKxatQoRERFQKBSIi4vD/v374ejoKGXYOu3999/H9u3bsXPnTlhYWCArKwtZWVkoLCyUOjQiIp3i5uaGxMREBAUFoby8HPPmzcPgwYNx9+5dqUMjADIhuE57Q+Xl5cHS0hK5ublo0aKFWllRUREyMjJUq5Nrm6om1G/evBnBwcHNG0wz0/ZrR0Taa8uWLZg2bRoKCgrQpk0b7NixA/7+/lKHpXOq+/x+nkb1MJHmEUJUuul6skREJKWgoCAkJCTAxcUFWVlZGDhwIBYuXIiysjKpQ3thMWEiIiLSQN26dUN8fDwmTZoEIQQiIiIwcOBAta/6oubDhImIiEhDmZqaYuPGjdi+fTvMzc0RExMDhUKBgwcPSh3aC4cJExERkYYbP348kpKS4O7ujv/+978YPHgw5s6di8ePH0sd2guDCRMREZEW6Ny5M06dOqVaSicyMhL9+/fHrVu3JI7sxcCEiYiISEuYmJhg7dq1iIqKgoWFBY4fPw6FQoGff/5Z6tB0HhMmIiIiLTNmzBgkJyejR48euH//PoYOHYqZM2eitLRU6tB0FhMmIiIiLdSpUyecPHkSoaGhAIDly5fD19cX169flzgy3cSEiYiISEsZGRnhyy+/RHR0NCwtLXH69Gl4eHjgxx9/lDo0ncOEiYiISMuNHj0aKSkp8PLywp9//omRI0ciLCwMJSUlUoemM5gwERER6QAnJyccP34cH3zwAQBg9erV6NevH65evSpxZLqBCRMREZGOMDQ0xBdffIEff/wRLVu2VH2pfXR0tNShaT0mTERERDpm+PDhUCqV8Pb2Rl5eHt58802EhISgqKhI6tC0FhMmqpcTJ05AJpNp9EZE9CJzcHBAbGwsZs2aBQBYs2YN+vbti0uXLkkcmXZiwkT1kpiYCCGERm9ERC86AwMDLFu2DPv374eNjQ1SUlLg6emJXbt2SR2a1mHCRPXCHhwiIu0RGBgIpVIJX19f5Ofn4+2338bf//53FBYWSh2a1mDCRHV29uxZuLq6Sh0GERHVQdu2bXHkyBHMmzcPMpkM69evR58+fXDx4kWpQ9MKTJiozo4dOwZfX1+pwyAiojrS19fH4sWLcejQIdja2uLMmTPw9PTE9u3bpQ5N4zFhojorLS2Fvr6+2r5169bBzc0NJiYmsLS0hL+/v0TRERFRTQYOHAilUgl/f388evQIEyZMwLvvvotHjx5JHZrG0q+5CjWV6l6XcjlgbFy7unp6gIlJzXXNzOoWX2Xu3r0LOzs7tX3R0dGYPXu2qns3Pz8f165da/jJiIioybRp0waHDh3C0qVLsWjRImzevBm///47vv/+e7i4uEgdnsZhD5OEzM2r3v7yF/W6rVtXXTcwUL1uhw6V16uPNWvWYMeOHarHv/zyCwICAtTqpKenw8HBAQEBAXBwcICLiwuGDBlSvxMSEVGzkcvlmD9/Po4cOQI7OztcuHABXl5e+Pbbb3m38XOYMFGVSkpK0KdPH0RFRan25ebmwsrKSq3ee++9B7lcDmtra5ibm+PKlStNGte1a9fQs2fPJj0HEdGLxM/PD0qlEgEBASgsLMSkSZMwYcIE5OfnSx2axuCQnIQePqy6TC5Xf5ydXXVdvefS3sYaDTM0NISnpyfy8vLw6NEj6OnpweTZsT88mc80btw4eHl5YcOGDbCyskLHjh0bJwAiImo2rVu3xoEDB/DZZ59h/vz52LFjBxISEvD999/D3d1d6vAkxx4mCZmZVb09O3+pprrP5TBV1quvgIAAHDp0CEeOHMGAAQPUyvbs2YPLly9j/fr16NmzJzp16qS2RtO1a9fg7u6O4OBgdOvWDVOnTsXevXvRu3dvuLi4qFacHTZsGDw9PeHq6ordu3ernr9hwwa4ubnB3d0ds2fPVu0vLS1FUFAQunbtirFjx7LrmIioEejp6WHOnDmIiYlB27ZtkZ6ejt69e+Obb7554f+fZcJENRo+fDj27duHq1evVug9KikpQWZmJrZt24Zr167h3Llz+Oabb1BaWqqqk5aWhjlz5uDs2bOIiYnBiRMn8PvvvyM0NBRfffUVAGDLli1ISkrCiRMnMHfuXAghcObMGXz55Zc4fvw4UlNTMXPmzArHvHDhAu7evYvjx483zy+DiOgF4OPjA6VSiSFDhqC4uBhTpkzBuHHjkJeXJ3VokmHCRDVydXVFWloa9J4f+wMwbtw4hISEYO7cuejcuTMGDhyIuLg4GBgYqOo4OzvD2dkZcrkcXbt2xcCBAwEA3bt3V91Nt3LlSri7u8PX1xc3btxAVlYWYmJiMHbsWFhaWgIArK2t1Y7ZpUsXyGQyeHh48K48IqJGZmNjg3379mH58uXQ19fH999/jx49eiA5OVnq0CTBhIlqpW/fvpVOtNbX18eKFStw8+ZNlJSUICsrS+2uOgAwMjJS/aynp6d6rKenh7KyMhw9ehQnTpzA6dOnkZqaCgcHBxQXF1cbz7PHlMvlKCsra0jziIioEnp6evjoo49w7NgxODg44MqVK/D29sa//vWvF26IjgkT1cqcOXPQq1evJjl2Xl4eWrVqBRMTE8THxyM9PR0A4O/vj6ioKOTm5gIAHjx40CTnJyKi6vXp0wcpKSkYPnw4SkpKMH36dLz55pvIycmROrRm80IlTDdv3oSfnx+6deuG7t2744cfflCV5efnw8vLCwqFAm5ubtiwYYOEkWqel156qdIhucYwaNAg5ObmQqFQYM2aNXBzcwPwZChwxowZ6NevHxQKBZYvX94k5ycioppZW1tj7969WLlyJQwMDLB79254eHggPj5e6tCahUy8QH1qmZmZuHv3LhQKBbKzs9GjRw9cvHgRZmZmKCsrQ3FxMUxNTVFQUABXV1ckJCSgVatWNR43Ly8PlpaWyM3NRYsWLdTKioqKkJGRAScnJxg/f+sbaTReOyKiyiUkJGDs2LHIyMiAgYEBPvvsM4SFhandJa0Nqvv8ft4L1cNkZ2cHhUIB4Ml6E9bW1qphHrlcDlNTUwBPPijLyspeuPFZIiKi2vDy8kJycjL+8pe/oLS0FOHh4RgxYoROT53QqIQpLi4Ow4YNg729PWQyGfbu3Vuhztdff636i9/T0xPHjh2r17kSExNRXl6O9u3bq/bl5OTA3d0d7dq1w8yZM2FjY1PvthAREekyKysr/PDDD/jqq69gaGiIn376CQqFAqdOnZI6tCahUQnTo0eP4O7urlqb53lRUVEICwvDJ598gpSUFLz66qsIDAzEjRs3VHWeLn74/Hbnzh1Vnfv372PixIlYv3692vGtrKyQmpqKjIwM7Ny5E3fv3m2ahhIREekAmUyG999/H6dPn0anTp1w8+ZNvPrqq/j8889RXl4udXiNSmPnMMlkMuzZswcjR45U7evduzd69OiBtWvXqvZ17doVI0eORGRkZK2OW1xcjNdffx3vvfceJkyYUGW9qVOnwt/fH2+99Valx3j2tve8vDy0b9+ec5h0DK8dEVHt5eXl4e9//zt27doFAAgMDMTWrVs1erRGJ+cwlZSUICkpCQEBAWr7AwICcPLkyVodQwiB4OBg+Pv7V0iW7t69q1rBNC8vD3FxcXB2dq70OJGRkbC0tFRtzw7rERERvYhatGiBnTt3Yv369TA2NsaBAwegUCjqPXVG02hNwnTv3j2UlZXB1tZWbb+trS2ysrJqdYwTJ04gKioKe/fuhUKhgEKhwNmzZwEAt27dgq+vL9zd3eHj44OQkBB079690uPMmTMHubm5qu3mzZsNaxwREZEOkMlkeO+99/D777/D2dkZt2/fhp+fHz799FOtH6LTlzqAunr+lkUhRK1vY/Tx8anygnl6ekKpVNbqOEZGRmorTRMREdH/6d69OxITEzFt2jRs27YN8+bNQ2xsLLZt21ah40NbaE0Pk42NDeRyeYXepOzsbK395RMREekqc3NzbN26FZs3b4apqSkOHz4MhUKBo0ePSh1avWhNwmRoaAhPT08cPnxYbf/hw4fRt29fiaIiIiKi6gQHByMhIQEuLi7IysrCwIED8Y9//EPrvgNUoxKmhw8fQqlUqobGMjIyoFQqVcsGhIeHY+PGjfj222+RlpaGDz74ADdu3MCUKVOkDJuIiIiq0a1bN8THx2PSpEkoLy/HokWL8PrrryMzM1Pq0GpNoxKmxMREeHh4wMPDA8CTBMnDwwMLFiwAAIwdOxarVq1CREQEFAoF4uLisH//fjg6OkoZNhEREdXA1NQUGzduxPbt22FmZoajR49CoVBUGDnSVBqVMPn5+UEIUWH77rvvVHWmTZuGa9euobi4GElJSfD19ZUu4BfEunXr4ObmBhMTE1haWsLf31/qkIiISEuNHz8eSUlJcHd3R3Z2NgYNGoR58+bh8ePHUodWLa27S06nPHpUdZlcDjy7WGJ1dfX0ABOTmuuamdUtPgDR0dGYPXs21q9fjz59+iA/Px/Xrl2r83GIiIiecnZ2xqlTpxAeHo5169bh008/RVxcHHbu3Il27dpJHV6lNKqH6YVjbl719pe/qNdt3brquoGB6nU7dKi8Xj2kp6fDwcEBAQEBcHBwgIuLC4YMGVK/9tbg2rVr6NmzZ5Mcm4iINIuJiQnWrl2LqKgoWFhY4NixY1AoFNi/f7/UoVWKCRNV67333oNcLoe1tTXMzc1x5coVqUMiIiIdMmbMGCQnJ6NHjx64f/8+hgwZgpkzZ6K0tFTq0NQwYZLSw4dVb9HR6nWzs6uue+CAet1r1yqvV0elpaUYN24cvLy8EB8fD6VSiY4dO/7/U1yDu7s7goOD0a1bN0ydOhV79+5F79694eLigkuXLgEAhg0bpvpC5N27d6uOvWHDBri5ucHd3R2zZ89WO2dQUBC6du2KsWPHQkO/6pCIiBpRp06dcPLkSYSEhAAAli9fjtdee011l7xGENRgubm5AoDIzc2tUFZYWCguXLggCgsLJYisYaKiooSjo2OlZRkZGcLAwED88ccf4vHjx6JLly7io48+EkIIsXbtWjF9+nQhhBD3798XQgiRk5MjnJ2dRXl5uUhNTRWurq4iJydHrc7TY6alpYny8nLx2muvibi4uCZuZdW0+doREWmr6OhoYWlpKQCIli1bih9//LHJzlXd5/fz2MNEVSopKUFmZia2bduGa9eu4dy5c/jmm29U3aTOzs5wdnaGXC5H165dMXDgQABPlsR/OjF85cqVcHd3h6+vL27cuIGsrCzExMRg7NixsLS0BABYW1urzuns7IwuXbpAJpPBw8ODE8yJiF4wo0ePRkpKCry8vPDnn39ixIgRCA8PR0lJiaRxMWGiKo0bNw4hISGYO3cuOnfujIEDByIuLg4GBgYAoPZ9enp6eqrHenp6KCsrw9GjR3HixAmcPn0aqampcHBwQHFxcbXnfPaYcrlc61aCJSKihnNycsLx48cRHh4O4Mkf3z4+PnhU3R3jTYwJE1VJX18fK1aswM2bN1FSUoKsrCzs2LGj1s/Py8tDq1atYGJigvj4eKSnpwMA/P39ERUVhdzcXADAgwcPmiR+IiLSXoaGhlixYgX27duHli1bwt3dHWb1WB6nsXAdJmoygwYNwpo1a6BQKODu7g43NzcAgKurK2bMmIF+/fpBX18fgYGBiIyMlDhaIiLSRMOGDUNqaipatWolaRwyIXgbUkPl5eXB0tISubm5aNGihVpZUVERMjIy4OTkBONnF6IkjcdrR0Sk26r7/H4eh+SIiIiIasCEiYiIiKgGTJiIiIiIasCEiYiIiKgGTJiaCefWax9eMyIieooJUxN7ushjQUGBxJFQXT29Zk+vIRERvbi4DlMTk8vlsLKyQnZ2NgDA1NQUMplM4qioOkIIFBQUIDs7G1ZWVpDL5VKHREREEmPC1AzatGkDAKqkibSDlZWV6toREdGLjQlTM5DJZLCzs0Pr1q1VX1xLms3AwIA9S0REpMKEqRnJ5XJ+CBMREWkhTvomIiIiqgETJiIiIqIaMGEiIiIiqgHnMDWCpwsc5uXlSRwJERER1dbTz+3aLFTMhKkR5OfnAwDat28vcSRERERUV/n5+bC0tKy2jkzw+x8arLy8HHfu3IGFhUWjL0qZl5eH9u3b4+bNm2jRokWjHltTsI26gW3UDWyjbmAba0cIgfz8fNjb20NPr/pZSuxhagR6enpo165dk56jRYsWOvuif4pt1A1so25gG3UD21izmnqWnuKkbyIiIqIaMGEiIiIiqoF80aJFi6QOgqonl8vh5+cHfX3dHUFlG3UD26gb2EbdwDY2Lk76JiIiIqoBh+SIiIiIasCEiYiIiKgGTJiIiIiIasCEiYiIiKgGTJg02Ndffw0nJycYGxvD09MTx44dkzqkBomLi8OwYcNgb28PmUyGvXv3qpULIbBo0SLY29vDxMQEfn5+OH/+vETR1l1kZCS8vLxgYWGB1q1bY+TIkbh48aJaneLiYoSGhsLGxgZmZmYYPnw4bt26JVHEdbd27Vp0795dtVCct7c3Dhw4oCrX9vZVJjIyEjKZDGFhYap92t7ORYsWQSaTqW1t2rRRlWv7e/Gp27dv469//StatWoFU1NTKBQKJCUlqcp1oZ0dOnSocC1lMhnef/99ANr/WgWAx48fY968eXBycoKJiQk6duyIiIgIlJeXq+o0y7UUpJF27dolDAwMxIYNG8SFCxfEjBkzhJmZmbh+/brUodXb/v37xSeffCKio6MFALFnzx618mXLlgkLCwsRHR0tzp49K8aOHSvs7OxEXl6eRBHXzaBBg8TmzZvFuXPnhFKpFEOGDBEODg7i4cOHqjpTpkwRbdu2FYcPHxbJycmif//+wt3dXTx+/FjCyGtv37594ueffxYXL14UFy9eFHPnzhUGBgbi3LlzQgjtb9/z4uPjRYcOHUT37t3FjBkzVPu1vZ0LFy4ULi4uIjMzU7VlZ2eryrX9vSiEEA8ePBCOjo4iODhY/P777yIjI0P8+uuv4vLly6o6utDO7Oxstet4+PBhAUAcPXpUCKH9r1UhhFiyZIlo1aqV+M9//iMyMjLEDz/8IMzNzcWqVatUdZrjWjJh0lC9evUSU6ZMUdvXpUsXMXv2bIkialzPJ0zl5eWiTZs2YtmyZap9RUVFwtLSUqxbt06KEBssOztbABCxsbFCCCFycnKEgYGB2LVrl6rO7du3hZ6envjll1+kCrPBWrZsKTZu3Khz7cvPzxevvPKKOHz4sHjttddUCZMutHPhwoXC3d290jJdeS/OmjVL+Pj4VFmuK+183owZM8TLL78sysvLdeK1KoQQQ4YMEe+++67avtGjR4u//vWvQojmu5YcktNAJSUlSEpKQkBAgNr+gIAAnDx5UqKomlZGRgaysrLU2mxkZITXXntNa9ucm5sLALC2tgYAJCUlobS0VK2N9vb2cHV11co2lpWVYdeuXXj06BG8vb11rn3vv/8+hgwZgoEDB6rt15V2Xrp0Cfb29nBycsK4ceNw9epVALrzXty3bx969uyJt956C61bt4aHhwc2bNigKteVdj6rpKQE27dvx7vvvguZTKYzr1UfHx8cOXIE6enpAIDU1FQcP34cb7zxBoDmu5a6u/ynFrt37x7Kyspga2urtt/W1hZZWVkSRdW0nrarsjZfv35dipAaRAiB8PBw+Pj4wNXVFcCTNhoaGqJly5ZqdbXtup49exbe3t4oKiqCubk59uzZg27dukGpVOpE+wBg165dSE5ORkJCQoUyXbiOvXv3xtatW9G5c2fcvXsXS5YsQd++fXH+/HmdeS9evXoVa9euRXh4OObOnYv4+HhMnz4dRkZGmDhxos6081l79+5FTk4OgoODAejGaxUAZs2ahdzcXHTp0gVyuRxlZWX49NNP8fbbbwNovs8PJkwaTCaTqT0WQlTYp2t0pc0hISE4c+YMjh8/XmNdbWujs7MzlEolcnJyEB0djaCgIMTGxlZZX9vad/PmTcyYMQOHDh2CsbFxrZ+nTe0MDAxU/ezm5gZvb2+8/PLL2LJlC/r06QNA+9+L5eXl6NmzJ5YuXQoA8PDwwPnz57F27VpMnDhRVU/b2/msTZs2ITAwEPb29tXW07Y2RkVFYfv27di5cydcXFygVCoRFhYGe3t7BAUFqeo19bXkkJwGsrGxgVwur/AXQHZ2doUMWlc8vUNHF9ocGhqKffv24ejRo2jXrp1qf5s2bVBSUoI///xTrb62tdHQ0BCdOnVCz549ERkZCXd3d6xevVpn2peUlITs7Gx4enpCX18f+vr6iI2NxZdffgl9fX3Y2trqRDufZWZmBjc3N1y6dEln3ot2dnbo1q2b2r6uXbvixo0bAHTr/xwAuH79On799VdMnjxZtU9X3pMff/wxZs+ejXHjxsHNzQ0TJkzABx98gMjISADNdy2ZMGkgQ0NDeHp64vDhw2r7Dx8+jL59+0oUVdNycnJCmzZt1NpcUlKC2NhYrWmzEAIhISHYvXs3fvvtNzg5OamVe3p6wsDAQK2NmZmZOHfunNa0sTJCCBQXF+tM+wYMGICzZ89CqVSqtp49e2L8+PGqn3Whnc8qLi5GWloa7OzsdOK9CAD9+vWrsKxHeno6HB0dAejG/znP2rx5M1q3bo0hQ4ao9unKe7KgoAB6eurpilwuVy0r0GzXstGmj1OjerqswKZNm8SFCxdEWFiYMDMzE9euXZM6tHrLz88XKSkpIiUlRQAQX3zxhUhJSVEtlbBs2TJhaWkpdu/eLc6ePSvefvttrbrFd+rUqcLS0lLExMSo3eZbUFCgqjNlyhTRrl078euvv4rk5GTh7++vVbf4zpkzR8TFxYmMjAxx5swZMXfuXKGnpycOHTokhND+9lXl2bvkhND+dn744YciJiZGXL16VZw+fVoMHTpUWFhYqP5/0fb3ohBPloTQ19cXn376qbh06ZLYsWOHMDU1Fdu3b1fV0YV2CiFEWVmZcHBwELNmzapQpu2vVSGECAoKEm3btlUtK7B7925hY2MjZs6cqarTHNeSCZMGW7NmjXB0dBSGhoaiR48eqtvTtdXRo0cFgApbUFCQEOLJraELFy4Ubdq0EUZGRsLX11ecPXtW2qDroLK2ARCbN29W1SksLBQhISHC2tpamJiYiKFDh4obN25IF3Qdvfvuu6rX5EsvvSQGDBigSpaE0P72VeX5hEnb2/l0jRoDAwNhb28vRo8eLc6fP68q1/b34lM//fSTcHV1FUZGRqJLly5i/fr1auW60s6DBw8KAOLixYsVyrT9tSqEEHl5eWLGjBnCwcFBGBsbi44dO4pPPvlEFBcXq+o0x7WUCSFE4/VXEREREekezmEiIiIiqgETJiIiIqIaMGEiIiIiqgETJiIiIqIaMGEiIiIiqgETJiIiIqIaMGEiIiIiqgETJiIiIqIaMGEiIiIiqgETJiIiIqIaMGEiIqqFjz/+GEOHDpU6DCKSCL9LjoioFnJyciCXy2FhYSF1KEQkASZMRERERDXgkBwRUQ3u3bsHmUyG8+fPSx0KEUmECRMRUQ1SU1NhZGQEZ2dnqUMhIokwYSIiqkFqaipcXFygr68vdShEJBEmTERENVAqlVAoFFKHQUQSYsJERFSD1NRUuLu7Sx0GEUmICRMRUTVKS0uRlpbGhInoBceEiYioGhcuXEBpaSkTJqIXHBMmIqJqKJVKODo6wsrKSupQiEhCTJiIiKqRkJCAXr16SR0GEUmMCRMRUSWKioqQnJyM6OhoDBo0SOpwiEhiTJiIiCqxatUqDBs2DCNHjsTEiROlDoeIJMbvkiMiIiKqAXuYiIiIiGrAhImIiIioBkyYiIiIiGrAhImIiIioBkyYiIiIiGrAhImIiIioBkyYiIiIiGrAhImIiIioBkyYiIiIiGrAhImIiIioBv8Pi7ILIb+BbuYAAAAASUVORK5CYII=", "text/plain": [ "Figure(PyObject
)" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/plain": [ "PyObject Text(0.5, 28.0, '$j$')" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "using PyPlot\n", "\n", "n = size(A,2)\n", "semilogy(diag(RC), \"bo\")\n", "semilogy(diag(RM), \"rx\")\n", "semilogy(2.0 .^ -(1:n), \"k-\")\n", "semilogy(fill(sqrt(eps()), n), \"b--\")\n", "semilogy(fill(eps(), n), \"r--\")\n", "legend([\"classical\",\"modified\",L\"2^{-j}\", L\"\\sqrt{\\epsilon_\\mathrm{mach}}\", L\"\\epsilon_\\mathrm{mach}\"], loc=\"lower left\")\n", "ylabel(L\"r_{jj}\")\n", "xlabel(L\"j\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note that the errors in classical Gram–Schmidt are at the level of $\\sqrt{\\epsilon_\\mathrm{mach}}$, whereas they are at the level of $\\epsilon_\\mathrm{mach}$ for modified Gram–Schmidt. Clearly, the errors are accumulating too fast for stability in the former algorithm." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Arbitrary-precision Gram–Schmidt" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "It turns out that the \"errors\" in the modifed Gram–Schmidt $r_{jj}$ are mainly due to errors in constructing $A$ in the first place (which make `A` less singular than we thought: its actual singular values don't go down to $2^{-80}$), not errors during the Gram–Schmidt process. We can see this by comparing to a computation in **arbitrary precision** arithmetic via `big` in Julia. (The default is about 80 decimal places, which is more than sufficient for us here.) In arbitrary precision, even classical Gram–Schmidt is fine." ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAkwAAAG0CAYAAADATXgqAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAAPYQAAD2EBqD+naQAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4xLjEsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy8QZhcZAAAgAElEQVR4nOzdeVxUVf8H8M+AgOwKCMgqJZJLLqC4YAlZipYKammZuJallEsuoYEI9qCZayg/0lLaXCqhnvTJ7REtUVGSFhcEHgpcCTAWwUFm7u+PcUYGZmCQgWHg8369eOnce+bec4mYr+d8z/eIBEEQQERERERqGei6A0REREQtHQMmIiIionowYCIiIiKqBwMmIiIionowYCIiIiKqBwMmIiIionowYCIiIiKqRztdd6A1kEqluHHjBiwtLSESiXTdHSIiItKAIAgoLS2Fk5MTDAzqHkNiwKQFN27cgKurq667QURERI8gLy8PLi4udbZhwKQFlpaWAGTfcCsrKx33hoiIiDRRUlICV1dXxed4XRgwaYF8Gs7KyooBExERkZ7RJJ2GSd9ERERE9WDARERERFQPBkxERERE9WAOExG1CBKJBPfv39d1N4ioFTEyMoKhoaFWrsWA6YHg4GAkJydj+PDh+Oabb3TdHaI2QxAE3Lp1C//884+uu0JErVCHDh3g6OjY6DqJDJgeePvttzFz5kwkJCTouitEbYo8WLK3t4eZmRmLvxKRVgiCgPLycuTn5wMAOnfu3KjrMWB6ICAgAMnJybruBlGbIpFIFMGSra2trrtDRK2MqakpACA/Px/29vaNmp7Ti6TvkydPYsyYMXBycoJIJEJSUlKtNtu2bYOHhwfat28PHx8f/PTTTzroKRE1hDxnyczMTMc9IaLWSv77pbE5knoxwnT37l306dMHM2bMwIQJE2qd37t3LxYsWIBt27bBz88P8fHxGDVqFC5dugQ3NzcAgI+PD8Rica33Hj58GE5OTg3qj1gsVrpWSUlJA5+IiKrjNBwRNRVt/X7Ri4Bp1KhRGDVqlNrzGzZswKxZszB79mwAwKZNm3Do0CHExcUhJiYGAJCWlqa1/sTExGDVqlVau15DSKQCUnOKkF96D/aW7eHrYQNDA37YEBERNSW9CJjqUllZibS0NLz77rtKx0eMGIGUlJQmuWdYWBgWLVqkeC3fi6ZJRUYis6AcM1wDce1OheKwS0dT7Mz7EZ52ZkBkZNP2gYiIqI3SixymuhQUFEAikcDBwUHpuIODA27duqXxdUaOHIkXX3wRBw8ehIuLC86dO6e2rYmJiWLfuObaPy6zoByeW9dh8dl92D93CC6uGon9c4dg8dl98Ny6DpkF5U3eB6KWTCIVcDq7EN+lX8fp7EJIpILO+vLnn39CJBIhPT29ye+1a9cudOjQQWvXS05Ohkgk0lmZh+b83qni7++PBQsW1NmmS5cu2LRpUzP1SD+pyzdubFtd0vsRJrmac5SCIDRo3vLQoUPa7pLWSKQCZrgGYnFwCYIS44F+zkB4OLwTYuGdGI+k4DlY7xaIZKnA6Tlqk3784yZWH7hca/T1vee7I7BX45YSt3STJk3C6NGjm/2+Fy5cwJo1a3Dy5EkUFRXB0dERTz75JObMmYMXXnhB8fv322+/xQcffIArV65AKpXCzc0NgYGBWL9+vcrrurq64ubNm7Czs2vOx2mQc+fOwdzcXPFaJBIhMTERQUFBOuxVy3Lz5k107NhR6211Se9HmOzs7GBoaFhrNCk/P7/WqJO+Ss0pwrU7FXDbFANERQEREYCJiezPqCi4boxBXlEFUnOKdN1Vomb34x838eaXv+AJR0ul0dcnHC3x5pe/4Mc/buq6i03K1NQU9vb2zXrP7777DoMGDUJZWRkSEhJw6dIlfP311wgKCsJ7772H4uJiAMDRo0cxefJkTJw4EampqUhLS8P777+PyspKtdc2NDSEo6Mj2rVr3n/PN2QFVadOnRq8srMxK7SaswJ+Xf9tGsLR0REmJiZab6tLeh8wGRsbw8fHB0eOHFE6fuTIEQwZMkRHvdKu/NJ7AAAvB0sgPBwwNgYqK2V/hofDy9FSqR1RWyGRClh94DKGP2GPj6f2h7dbR5ibtIO3W0d8PLU/hj9hj/cPXm6S6TmpVIq1a9eia9euMDExgZubG95//33V/ZRIMGvWLHh4eMDU1BReXl7YvHmzUpvk5GT4+vrC3NwcHTp0gJ+fH/766y8AwK+//oqAgABYWlrCysoKPj4+OH/+PADVU3Lff/89+vfvj/bt28POzg7jx49XnPviiy/Qv39/WFpawtHREa+88oqisJ8m7t69i1mzZuH555/HgQMHMGLECDz++OPw9fXF7Nmz8euvv8La2hoA8MMPP2Do0KFYsmQJvLy80K1bNwQFBeGjjz5Se/2aU3LyKcJjx46hf//+MDMzw5AhQ5CRkVFnP5ctW4Zu3brBzMwMjz32GMLDw5UCj8jISPTt2xeffvopHnvsMZiYmEAQZD8nVVVVCA0NRYcOHWBra4v33ntPcQ5QnpLr0qULANluESKRSPFa3fV//PFHDB06VHHtF154AdnZ2bWef9++ffD390f79u3x8ccfw8rKqtYuFP/+979hbm6O0tJSld8Df39/hIaG1vssq1evxvTp02FtbY3XXnsNAHD9+nVMmjQJHTt2hK2tLcaNG4c///xT6fqffvopevbsCRMTE3Tu3BmhoaGKc9Wn2SorKxEaGorOnTujffv26NKli2JBVs22APD777/jmWeegampKWxtbfH666+jrKxMcX769OkICgrChx9+iM6dO8PW1hbz5s1r8sBSLwKmsrIypKenK/4HysnJQXp6OnJzcwEAixYtwo4dO/Dpp5/i8uXLWLhwIXJzc/HGG2/osttaY2/ZHgCQcbsUiI5+GCxVVgLR0ci4VarUThMtKd+D6FHJR1/nBnSFQY3paAMDEd7079pko69hYWFYu3YtwsPDcenSJXz11VdqR7WlUilcXFywb98+XLp0CREREVi+fDn27dsHQPYBHRQUhGHDhuG3337D6dOn8frrryumtaZMmaLIrZQvcjEyMlJ5rwMHDmD8+PF4/vnnceHCBUWgIVdZWYno6Gj8+uuvSEpKQk5ODqZPn67xcx8+fBiFhYVYunSp2jbyfjs6OuLixYv4448/NL6+OitWrMD69etx/vx5tGvXDjNnzqyzvaWlJXbt2oVLly5h8+bN2L59OzZu3KjUJisrC/v27cO3336rlDOVkJCAdu3a4ezZs9iyZQs2btyIHTt2qLyPPN91586duHnzplL+q6rr3717F4sWLcK5c+dw7NgxGBgYIDg4GFKpVOm6y5Ytw9tvv43Lly8jODgYkydPxs6dO5Xa7Ny5ExMnToSlpaXa74Mmz7Ju3Tr06tULaWlpCA8PR3l5OQICAmBhYYGTJ0/i559/hoWFBQIDAxUjUHFxcZg3bx5ef/11/P777/j+++/RtWtXlX3YsmULvv/+e+zbtw8ZGRn44osvFIFlTeXl5QgMDETHjh1x7tw5fP311zh69KhSMAYAx48fR3Z2No4fP46EhATs2rULu3btUvt90ApBDxw/flwAUOtr2rRpijZbt24V3N3dBWNjY8Hb21s4ceJEs/WvuLhYACAUFxc3yfWrJFLBb80xITF4jiAAghAVJTsRFSUIgJAYPEcYuvaYUCWR1n+xlSuFq/OWCH5rjgnuy35QfPmtOSZcnbdEEFaubJJnIFKloqJCuHTpklBRUfFI70+6cE1wX/aDUHbvvsrzpffuC+7LfhCSLlxrTDdrKSkpEUxMTITt27erPJ+TkyMAEC5cuKD2GnPnzhUmTJggCIIgFBYWCgCE5ORklW0tLS2FXbt2qTy3c+dOwdraWvF68ODBwpQpUzR9FCE1NVUAIJSWlgqC8PD37Z07d1S2X7NmjQBAKCoqUrqGubm54uvf//63IAiCUFZWJowePVoAILi7uwuTJk0SPvnkE+HevXtq+1Pzeyfvz9GjRxVtDhw4IABo0M/NBx98IPj4+Cher1y5UjAyMhLy8/OV2g0bNkzo3r27IJU+/H26bNkyoXv37orX7u7uwsaNGxWvAQiJiYlK11F3/Zry8/MFAMLvv/8uCMLD59+0aZNSu7NnzwqGhobC9evXBUEQhL///lswMjJS+zPTkGcJCgpSet8nn3wieHl5Kb1PLBYLpqamwqFDhwRBEAQnJydhxYoVau9d/Xvy1ltvCc8884zS9dS1/fjjj4WOHTsKZWVlivMHDhwQDAwMhFu3bgmCIAjTpk0T3N3dhaqqKkWbF198UZg0aZLK69f1e6Yhn996McLk7+8PQRBqfVWPJufOnYs///wTYrEYaWlpePrpp3XXYS2RjwL98NsNrLv8PYIeJHinhYSiTFyFtJBQJAXPQVBiPD7N/VGjhG+utqPWRGn0VYVHGX3VxOXLlyEWizF8+HCN3/N///d/6N+/Pzp16gQLCwts375dMUpuY2OD6dOnY+TIkRgzZgw2b96Mmzcf5l4tWrQIs2fPxrPPPos1a9YoTeHUlJ6eXme/Lly4gHHjxsHd3R2Wlpbw9/cHAEVfHkXv3r0VswB3795FVVUVAMDc3BwHDhxAVlYW3nvvPVhYWOCdd96Br68vyssb9rumd+/eir/L9wSrayrxm2++wdChQ+Ho6AgLCwuEh4fXekZ3d3d06tSp1nsHDRqktGho8ODByMzMhEQiaVCfVV0/Ozsbr7zyCh577DFYWVnBw8MDQO3vf/VRQQDw9fVFz5498dlnnwEAPv/8c7i5udX7WafJs9S8V1paGrKysmBpaQkLCwtYWFjAxsYG9+7dQ3Z2NvLz83Hjxg2Nf/6nT5+O9PR0eHl54e2338bhw4fVtr18+TL69OmjlFTv5+cHqVSqNA3bs2dPpW1OOnfu3KCp5UehFwFTmxMZiczQpRi27jhe3n4G8/ek42xWPuKeCcHN4gr8NGUeeq08hAlxKVg/6CVkzlsiq8NUD/lqO3mQ5Z0QK8v3SIhVBGMz3QI5PUd6w9fDBi4dTbHteBakNX5upVIBcclZcLUxha+HjVbvK9+fSlP79u3DwoULMXPmTBw+fBjp6emYMWOGUoLtzp07cfr0aQwZMgR79+5Ft27dcObMGQCyfJiLFy/i+eefx3//+1/06NEDiYmJDe7b3bt3MWLECFhYWOCLL77AuXPnFNfRNNnX09MTAJQ+vExMTNC1a1e1UzKPP/44Zs+ejR07duCXX37BpUuXsHfvXo3uJ1d9ClIeANScxpI7c+YMJk+ejFGjRuGHH37AhQsXsGLFilrPWP1DuSmouv6YMWNQWFiI7du34+zZszh79iyA2t9/Ve+dPXu2Ylpu586dmDFjhlaqWNe8l1QqhY+PjyIIln9dvXoVr7zySoN//r29vZGTk4Po6GhUVFTgpZdewsSJE1W2FepY4V79eM0paZFIpPbnQVsYMLVAqkaBnv5qKzpbm+LN/36GgV07YfPkvtj92iAkLw6AZ+wHGhWt5Go7am0MDUR47/nuOHYlH69/fh5pf92Rjb7+dQevf34ex67kY8Xo7lovt+Hp6QlTU1McO3ZMo/Y//fQThgwZgrlz56Jfv37o2rWrylGifv36ISwsDCkpKejVqxe++uorxblu3bph4cKFOHz4MMaPH18rn0Wud+/eavt15coVFBQUYM2aNXjqqafwxBNPNPhf5SNGjICNjQ3Wrl3boPfJdenSBWZmZrh79+4jvV8Tp06dgru7O1asWIH+/fvD09NTkUCvCXmgWv21p6en2o1bjYyMNBp9KiwsxOXLl/Hee+9h+PDh6N69O+7cuaNxv1599VXk5uZiy5YtuHjxIqZNm1bvexr6LIAswMnMzIS9vb0iEJZ/WVtbw9LSEl26dNH45x8ArKysMGnSJGzfvh179+7Ft99+i6Ki2p81PXr0UIxUyp06dQoGBgbo1q2bxvdrCgyYWhhNRoGW9hiHF3o7YfDjtg36IOBqO2qNAnt1RtwUb1y5VYoJcSmK0deM26WIm+LdJHWY2rdvj2XLlmHp0qX47LPPkJ2djTNnzuCTTz5R2b5r1644f/48Dh06hKtXryI8PFwpOTgnJwdhYWE4ffo0/vrrLxw+fBhXr15F9+7dUVFRgdDQUCQnJ+Ovv/7CqVOncO7cOXTv3l3lvVauXIndu3dj5cqVuHz5Mn7//Xd88MEHAAA3NzcYGxvjo48+wv/+9z98//33iI6ObtCzW1hYYMeOHThw4ACef/55HDp0CP/73//w22+/Ke4j/zCOjIzE0qVLkZycjJycHFy4cAEzZ87E/fv38dxzzzXovg3RtWtX5ObmYs+ePcjOzsaWLVvUjsipkpeXh0WLFiEjIwO7d+/GRx99hPnz56ttLw8ebt26VWcAJF9x9vHHHyMrKwv//e9/lXaNqE/Hjh0xfvx4LFmyBCNGjICLi4vWnwWQLTKws7PDuHHj8NNPPyEnJwcnTpzA/Pnzce3aNQCy/7br16/Hli1bkJmZiV9++UXt6seNGzdiz549uHLlCq5evYqvv/4ajo6OKguuTpkyBe3bt8e0adPwxx9/4Pjx43jrrbcwdepUnZcKajWFK1sLpVGgfs6y0Z/Vq2WBTVQUXENCkReXgtScIgx+3LZB166e7+GdEFt7tV1IqFI7In0R2Ksznuvh2Kz7LIaHh6Ndu3aIiIjAjRs30LlzZ7Urc9944w2kp6dj0qRJEIlEePnllzF37lz85z//ASDbTf3KlStISEhAYWGhYon2nDlzUFVVhcLCQoSEhOD27duKMgHq9rP09/fH119/jejoaKxZswZWVlaKPJdOnTph165dWL58ObZs2QJvb298+OGHGDt2bIOePTg4GCkpKVi7di1CQkJQVFQEa2tr9O/fH3v27MELL7wAABg2bBi2bt2q6HvHjh3Rr18/HD58GF5eXg26Z0OMGzcOCxcuRGhoKMRiMZ5//nmEh4cjUsPto0JCQlBRUQFfX18YGhrirbfewuuvv662/fr167Fo0SJs374dzs7OtZbfyxkYGGDPnj14++230atXL3h5eWHLli2KPDJNzJo1C1999VW9qwQf9VkA2c/jyZMnsWzZMowfPx6lpaVwdnbG8OHDFTtbTJs2Dffu3cPGjRuxePFi2NnZqZ1ms7CwwNq1a5GZmQlDQ0MMGDAABw8ehIFB7TEbMzMzHDp0CPPnz8eAAQNgZmaGCRMmYMOGDRo9b1MSPchQp0YoKSmBtbU1iouLG71Nynfp1zF/TzourhoJc5N2sikzeWAjFqNMXIVeKw9h8+S+GNfXud7rVd+s187cBMv2/4bFZ/fJKoZHRclGmqKjgYgIWcXwQS8heXEAK4ZTs7h37x5ycnLg4eGB9u0ZqBPV58svv8T8+fNx48YNGBsb19nW398fffv2bfPbuNT1e6Yhn98cYWphtDYKpGaz3mXn9iHov5/hyhM+uBsSCi9xFTJCQpF34TqCEuPR08kKhgbPNNnzERFRw5WXlyMnJwcxMTGYM2dOvcESaR8DphZGvuond0EYvFWMAuVduA7XQS/Vu+pHkTgeXAK3TTHwcrBExu1SmI1cCQA47tQTa+NSFO1dB72Enk5WGq22IyKi5vXBBx/g/fffx9NPP42wsDBdd6dNYsDUwhgaiLAz70d4Pkjwdn2EUaC6NuvFlTQkBc/BlwNfxJcT+qCgTFwt34MjS0RELVFkZKTGOVhyycnJTdKXtooBUwvkaWeGzHlL8KFrIK5pOApUPVepoFRcb+L4tbgUGIhEGuVBERERtXUMmFqiyEh4AjhRLQhSOwqkJlcJAEzXvA/YW9QuHyCWVeFl+QAiIiLNMGBqwQwNRPWWDlCVq/TtL9fw99L30P3nL1Ho6wdblg8gIiJqFBau1GPqily+euQzvPPzl0j16Avb1FOQrloFiMWK6t55C8OaZLsIIiKi1oojTHpMXZFLg8pKFPr6wTf1FNYPnQK7Z0MwgeUDiIiIHhlHmPRYXVud2I56FpfeXIyP/F7Gyu8vPtJmvUTUMvj7+2PBggWK1126dFEqRnjr1i0899xzMDc3V2w3IRKJkJSU1Kj7Tp8+HUFBQY26BlFrwREmPVZnkUtDQ1QsWw7EpSD8+e6wszRh+QCiVuLcuXNKO8xv3LgRN2/eRHp6OqytrQEAN2/eRMeOHZu8L4IgYMeOHfj0009x8eJFSKVSuLu749lnn8Vbb72Frl27AgDu3r2LqKgofP3117hx4wYsLS3Rs2dPLF68WLGVClFLxoBJj2la5HK6nwe3OqHWKzISMDSU/ezXFB0NSCSyNq1Ip06dlF5nZ2fDx8cHnp6eimOOjo5N3g9BEPDKK68gKSkJy5cvx8aNG2Fvb4+cnBwcOXIEq1evxq5duwDI9tNLTU1FbGwsevTogcLCQqSkpKCwsLDJ+0mkFQI1WnFxsQBAKC4ubvZ7X523RBAAITF4jnD+zyKh9N594fyfRUJi8BxBAGTnm0iVRCqkZBUISReuCSlZBUKVRNpk96LWqaKiQrh06ZJQUVHx6BeJihIEQPanJse1ZNiwYUJoaKgwf/58oUOHDoK9vb0QHx8vlJWVCdOnTxcsLCyExx57TDh48KDS+5KTk4UBAwYIxsbGgqOjo7Bs2TLh/v37ivNlZWXC1KlTBXNzc8HR0VH48MMPhWHDhgnz589XtHF3dxc2btyo+DsAxde0adMEQRAEAEJiYqLiPdeuXRNeeukloUOHDoKNjY0wduxYIScnR3G+qqpKWLhwoWBtbS3Y2NgIS5YsEUJCQoRx48ap/R7s3r1bACB89913Ks9LpQ9/J1hbWwu7du2q/xtLpGV1/Z5pyOc3R5j03KMUuWw0NbWfXDqayqqU25m1un/RUwsmH1mKiHj4+sEoq2LUtYkkJCRg6dKlSE1Nxd69e/Hmm28iKSkJwcHBihGXqVOnIjc3F2ZmZrh+/TpGjx6N6dOn47PPPsOVK1fw2muvoX379ooqzkuWLMHx48eRmJgIR0dHLF++HGlpaejbt6/KPpw7dw4hISGwsrLC5s2bYWpqWqtNeXk5AgIC8NRTT+HkyZNo164dVq9ejcDAQPz2228wNjbG+vXr8emnn+KTTz5Bjx49sH79eiQmJuKZZ9RP4e/evRteXl4YO3asyvMi0cORbUdHRxw8eBDjx4+HpaVlA77LRC1EU0RzbY0uR5jkmnO0p/qoVtpfRULZvftC2l/NM6pFrYtWRpjk5CNKxsZNOrIkN2zYMGHo0KGK11VVVYK5ubkwdepUxbGbN28KAITTp08LgiAIy5cvF7y8vJRGXrZu3SpYWFgIEolEKC0tFYyNjYU9e/YozhcWFgqmpqZqR5gEQRDGjRunGFmSQ7URpk8++aTWfcVisWBqaiocOnRIEARB6Ny5s7BmzRrF+fv37wsuLi51jjA98cQTwtixY5WOzZ8/XzA3NxfMzc0FZ2dnxfETJ04ILi4ugpGRkdC/f39hwYIFws8//6z22kTaoq0RJq6SayXkRS7H9XXG4MdtmyxnSV3tJ++EWAQ92P9uplsgJFKhSe5PpFaNlaJNObIk17t3b8XfDQ0NYWtriyeffFJxzMHBAQCQn58PALh8+TIGDx6sNPLi5+eHsrIyXLt2DdnZ2aisrMTgwYMV521sbODl5dWofqalpSErKwuWlpawsLCAhYUFbGxscO/ePWRnZ6O4uBg3b95Uum+7du3Qv3//eq9d/VkAYMWKFUhPT0dERATKysoUx59++mn873//w7FjxzBhwgRcvHgRTz31FKKjoxv1bETNhVNy1CDqaj/J96lzenUe8uJPY+ORq/DravdgVR4TzqkZREcrrxSNjm7yoMnIyEjptUgkUjomDyakUikAWZJ0zQBDEARFW/nftU0qlcLHxwdffvllrXM1E8gbwtPTE1euXKl1vU6dOsHe3r5WeyMjIzz11FN46qmn8O6772L16tWIiorCsmXLYGxs/Mj9IGoOHGEiBYlUwOnsQnyXfh2nswtVjhKprf1kYIDM26VYuDcdABB7PAsvbz+DYeuOIzN0KXOaqGlVz1mqVtUeLWz0okePHkhJSVEKjFJSUmBpaQlnZ2d07doVRkZGOHPmjOL8nTt3cPXq1Ubd19vbG5mZmbC3t0fXrl2VvqytrWFtbY3OnTsr3beqqgppaWl1Xvfll19GRkYGvvvuu0fqV48ePVBVVYV797ivJbV8HGGiBiVx11X7yXPrOrz1wh282zMIO6cPgLWZEXIXhMEzMV5WLFMXz0atn6oEb1WJ4C3A3LlzsWnTJrz11lsIDQ1FRkYGVq5ciUWLFsHAwAAWFhaYNWsWlixZAltbWzg4OGDFihUwMGjcv22nTJmCdevWYdy4cYiKioKLiwtyc3Oxf/9+LFmyBC4uLpg/fz7WrFkDT09PdO/eHRs2bMA///xT53UnT56M/fv3Y/LkyQgLC8PIkSPh4OCAv/76C3v37oWhoaGirb+/P15++WX0798ftra2uHTpEpYvX46AgABYWVk16vmImgMDJlK5gW/G7VKVwY6q2k+SFe/h08CZeO3ILkz+YQfKxFV4+l+jYfj+ang/yGta7xaIZKnA6TnSPolE9Wo4+WuJpPn7pIazszMOHjyIJUuWoE+fPrCxscGsWbPw3nvvKdqsW7cOZWVlGDt2LCwtLfHOO++guLi4Ufc1MzPDyZMnsWzZMowfPx6lpaVwdnbG8OHDFcHKO++8g5s3b2L69OkwMDDAzJkzERwcXOe9RSIR9u7di+3bt2Pnzp344IMPcP/+fbi4uGD48OHYsGGDou3IkSORkJCA5cuXo7y8HE5OTnjhhRcQIQ9qiVo6bWejt0UtYZXco6qSSAW/NccUK9wUK4serDhKDJ4jDF17TGnVXc3aT8cu3xbcl/0g7B4zW3aNGiuVzv9ZJLgv+0FIySrQ0VNSS6XVVXJERCqwDhNpRX1J3K4hociLS0FqThEGP24LQH3tp4/8JsPHrSM8/2+D0kolL3EVgIf5T0RERPqGAVMbVyuJWx4s1Qh2TmX9jfzSe7L96CJWwtNAhBNSAak5RTiV9Tdij2dj86R+8Pw8RTYFUm2lUkZIKICH+U9ERPaVBbUAACAASURBVET6hgFTG1fnBr4BAbjo0Qewfw6xx7MV76meDD44MhK+HjZISr+Ba4uWo38de9r5etjo6jGJiIgahQFTG1ffBr4Dk5OxcFgh/HZtRvfOViqTwQ0NRLIA6kGCt2tIKLzEVcgICUXehesISoxHTycrGBqo32KBiIioJWPA1MapC3YuvToPp49dxfwTX2D+iS+Az7sB4eHwTohVufJNJ3vaUashNFHBRiIibf1+EQn8TdVoJSUlsLa2RnFxsX7WE1FThwkA9uYfwcCcX4Hk5IdTdVFRSAsJxYS4FOx+bZAiGRyQFb9MzSl6mO/ESt9UB4lEgqtXr8Le3h62trb1v4GIqIEKCwuRn5+Pbt26KdUGAxr2+c0RJgIiI+EJKJK480vvIfN2GWKPZ6FX3IeASTvAxESjlW/yPe2INGFoaIgOHToo9lozMzOrtXUIEdGjEAQB5eXlyM/PR4cOHWoFSw3FgIkUqgc7p7MLEXs8S3UyOFe+kRY5OjoCeLhBLRGRNnXo0EHxe6YxGDCRSvUlg3PlG2mLSCRC586dYW9vj/v37+u6O0TUihgZGTV6ZEmOAROpxJVv1NwMDQ219ouNiEjbGDCRWtpc+cZkcCIi0mdcJacFer9Krh6NCnbUrMCrXvwSkZHavScREZEGuEqOtKoxK98yC8rhuXUdFgeXwG1TDLwcLFUWv1R4xACLiIioKTFgoiYjkQqY4RqIxcElCEqMl23uW0fxS+ARAiwiIqJmYKDrDlDrlZpThGt3KuC2KUa2yi4iQlbPKSICiIqC68YY5BVVIDWnCMDDACspeA6CEuPhnRALc5N28E6IRdCDAGumWyAkUs4iExFR82LARE1GXtTSy8FSVpJAXsdJXvzS0VKpXUMDLCIioubCgImajLyoZcbtUln9pprFL2+VKrVraIBFRETUXBgwUZOpXvxSPkoEsVgxepS3MAyuNqaK4pf1BViXbpQAADJvl+F0diGn5oiIqNkw6ZuaTEOLX6qtLh4QAERE4PSxq8CgyYg9noXY41lcOUdERM2GARM1qfqKXz5ua4rT2YWKekuf5v6IbjUCrIsefTAwORnzT3yBQY/Zolfch1w5R0REzYoBEzWtyEh4AjhRoxDlwM+2IBvA02bDcG37GUXziIzbMO49EEUlFVggD7Dsn8PCYYWYf+ILDMz5FXiwck5daQIiIiJtY8BEzaJm8cvMwgo19Zb+jS7HzmLo3Kew+7VBOJX1N2KPZ8Nv12bg824PV85VVspWzoWEIi8uBak5RY9cXJOIiKg+TPqupry8HO7u7li8eLGuu9KqaVJvaZb7KPh62MDTQbYyrntnK66cIyIinWHAVM3777+PgQMH6robrV5D6i01tDQBERFRU2DA9EBmZiauXLmC0aNH67orrV5D6i01tDQBERFRU9CLgOnkyZMYM2YMnJycIBKJkJSUVKvNtm3b4OHhgfbt28PHxwc//fRTg+6xePFixMTEaKvLVIeGjBrJSxPIp+rSQkJRJq5CWkioYkrv09wfmfBNRERNSi+Svu/evYs+ffpgxowZmDBhQq3ze/fuxYIFC7Bt2zb4+fkhPj4eo0aNwqVLl+Dm5gYA8PHxgVgsrvXew4cP49y5c+jWrRu6deuGlJSUWm1qEovFStcqKSlpxNO1PWrrLUVHy0aNLlyH66CXFKNG9ZUm8LQz09WjEBFRWyHoGQBCYmKi0jFfX1/hjTfeUDr2xBNPCO+++65G13z33XcFFxcXwd3dXbC1tRWsrKyEVatWqW2/cuVKAUCtr+Li4oY/UBt1dd4SQQCExOA5wvk/i4TSe/eF838WCYnBcwQBkJ2voUoiFVKyCoSkC9eElKwCoUoi1UHPiYiotSguLtb481skCIJe7S8hEomQmJiIoKAgAEBlZSXMzMzw9ddfIzg4WNFu/vz5SE9Px4kTJxp0/V27duGPP/7Ahx9+qLaNqhEmV1dXFBcXw8rKqoFP1EZFRiKzoBwzXANx7U6F4rCrjSk+zW1Y9W5JjRpPvh42nKIjIqJ6lZSUwNraWqPPb72YkqtLQUEBJBIJHBwclI47ODjg1q1bTXJPExMTmJiYNMm12ww1BS1lwc4zGl9DVdDFLVOIiEjb9D5gkhOJlEcUBEGodUwT06dP11KPSBM1C1o2RGZBuZril9wyhYiItEsvVsnVxc7ODoaGhrVGk/Lz82uNOlHroUnxy5lugZBI9WrGmYiIWii9D5iMjY3h4+ODI0eOKB0/cuQIhgwZoqNeUVNrSPFLIiKixtKLKbmysjJkZWUpXufk5CA9PR02NjZwc3PDokWLMHXqVPTv3x+DBw/Gxx9/jNzcXLzxxhs67DU1pVrFL1evVi5+Ka5SakdERNQYehEwnT9/HgEBAYrXixYtAgBMmzYNu3btwqRJk1BYWIioqCjcvHkTvXr1wsGDB+Hu7q6rLlMTq1780jshtnbxy5BQpXZERESNoRcBk7+/P+qrfjB37lzMnTu3mXpEutbQ4pdERESNoRcBE1FN8i1TPB8keLuGhMJLXIWMkFDkXbiOoMR49HSy0rxEARERUR0YMJHe4pYpRETUXPSu0ndL1JBKoaR9rPRNRESPok1V+iZqTPFLIiIiTeh9HSYiIiKipsaAiYiIiKgeDJiIiIiI6sGAiYiIiKgeDJiIiIiI6sGAiYiIiKgeLCtAbQprNhER0aNgwERtQ2QkMgvKMcM1ENfuVCgOu3Q0lW2xYmcGREbqrn9ERNSiMWCiNiGzoByeW9dhcXAJ3DbFwMvBEhm3S5G7IAyeifHInLcEnrruJBERtVgMmKjVkk+/3SquwL8cR2B5cAmCE+OBfs5AeDi8E2Lh/WDz3vVugUiWCpyeIyIilRgwUeujZvrtXz4T8fT1P2AbEQGsXg1UVgJRUXANCUVeXApSc4q4xQoREanEVXLU6iim387uw/65Q/DBxN4AgLDzX8M29RSkBoayYMnYGAgPh5ejJQAgv/SeLrtNREQtGAMmalUkUgEzXAORFDwHQYnx8E6IhWtHM7x1ajfGJ32My937w0AqgWBsLAuaoqORcasUAGBv2V7HvSciopaKU3LUqqTmFOHanQq4bYqR5SpFRGDQ6tUYXFmJy937o/vl81g/dAqG7NqMwV9tAyIikHfhOlwHvQRfDxtdd5+IiFooBkzUqsin1bwcLIHwcGD1aogqKyE1MET3y+exP+h1fOQ1Fi5Fd2EcEoq8C9cRlBiPnk5WMDR4Rse9JyKilooBE7Uq8mm1jNul8E6IVeQqGVRWotDXD2v6vwiUirHs298BAK6DXkJPJytZHSYiIiI1GDBRq+LrYQOXjqbIXRAG78R4ICpKNtIUHQ3biAiEOX+Nf/lMxPLRPeBoJa/0zZElIiKqGwMmalUMDUSyyt0P6iu5hoTCS1yFjAfTb8GJ8ejlZA3Pfs/puqtERKRHGDBRq+NpZ4bMeUvwoWsgrsWlKI7XNf3GPeaIiKguIkEQBF13Qt+VlJTA2toaxcXFsLKy0nV36AGNgiDuMUdE1GY15PObI0zUahkaiOqt3M095oiISBMsXEltlqoil+Ym7eCdEIugBzlQM90CIZFyEJaIqK1jwERtllKRy6goICICMDGR/RkVBdeNMcgrqkBqTpGuu0pERDrGgInarFpFLuXbpXCPOSIiqoEBE7VZ1YtcIjr6YbDEPeaIiKgGBkzUZlUvcimfhoNYrJiey1sYBlcbU433mJNIBZzOLsR36ddxOruQuU9ERK0IV8lRm1VfkUuN95hjaQIiolaPARO1aY9S5LImliYgImr9GDBR2xYZCU8AJ1QWuax/jzl5aYLFwSUISowH+jkD4eHwToiF94ORq/VugUiWCqwcTkSkxxgwEUF9kcv6qoUrlSbo5yzLhVq9WpY4HhUF15BQ5MWlIDWnqN4imkRE1HIxYCJSRcO8pFqlCeTBkrw0gbgKAHAq62/uU0dEpMcYMBGpoGleUvXSBN4JscqlCQICcNGjD2D/HGKPZyuuzWRwIiL9w7ICRDU0ZMuUOksTJCdj4M7NWHhmD755YzAurhqJ/XOHYPHZffDcug6ZBeW6flQiItIQAyaiGjTdMmXXqRz88NsNrLv8vSKQSgsJRZm4CqmvzsPmYa8CAOaf+AL9P99aK+ia4ToSp7IKWLeJiEgPiARB4G/pRiopKYG1tTWKi4thZWWl6+5QI32Xfh3z96Tj4qqRMDdpJwuW5FNtYWG4lH8Xo638Fe0X/PwlTEyMAABi8X1sGjpFcW5v/hEMzPkVSE5+OFXn74+zHn0wyf45pftyqo6IqHk15PObOUxENdSVl1T4n6PokXoKbw29iU5rV2OCjwsybg9B7oIwBCXG4/S0t7F5cl9k3i5D7PEs9Ir7EKgRdGX2HICBW9fhraEFaLcyArOfekwpP+r0tLeRn36dCeJERC0IAyaiGqrnJXknxsum5cLDIY2Kgu3KlUj16It3fv4S0qPdYDAkQrnmUo9xSO7thNScIsQez1IZdCVf/RtpL8zGOz/sQO4BV5g/GwPvaePhnZyMzcNexUbHEcCedAAcdSIiaimYw0RUg3zLlJp5SZ8/G4L1Q6fANycdhb5+MFi5UmVuU2pOkdpk8NwF7+K1I7sAEbDjuelw27RGdo3kZACAw4PRrZ3TBzBBnIioBeEIE5EK6rZMgd/LGP2kI7rbWwDp51TWXMovvad2n7r9z8+A5Hwe3vn3DmTOWwKckI06VRoa4eDYmZicGI/r/1SgZHJfBHwdz2rhREQtBAMmIlVUbJlSUCpG9IHLqHh3BVCz5lJ0NDJCQgE8zIGqK+ga6mmHgRdlAZfUyBjG9yvh62GD64vC8M6GGEh9vwbus1o4EVFLwYCJqA7Vt0yRSAXsTPmzVm4ToqOBiAjkXbgO10EvwdfDRvZmFUGXnbkJlu3/DTev3pNNw0VF4d9jZiLrrWV4Z0MMkoLnoJOhEYzvqx65IiIi3WDARKQhddNsGSGhyLtwHUGJ8ejpZFVr096a+9TVvIZlxX185PcynDuaYnJivKxRHSNXRETU/BgwETWAumk210EvoaeTlWw12yNe43bJgxEkf3/g+HH1I1dERNTsGDARNYSKabaH9ZKeqfftaq+xaR0Gn/hCNuq0MUajkSsiImo+DJgAZGRkYNKkSUqvd+/ejaCgIB32ilqymtNsjb5Glw6NHrkiIqKmw61RaigrK0OXLl3w119/wdzcXKP3cGsU0haJypErlhIgImoK3BqlEb7//nsMHz5c42CJSJu0MXJFRETapxeVvk+ePIkxY8bAyckJIpEISUlJtdps27YNHh4eaN++PXx8fPDTTz890r327dunND1HREREpBcjTHfv3kWfPn0wY8YMTJgwodb5vXv3YsGCBdi2bRv8/PwQHx+PUaNG4dKlS3BzcwMA+Pj4QCwW13rv4cOH4eTkBEA2NHfq1Cns2bOnzv6IxWKla5WUlDTm8YiIiKiF07scJpFIhMTERKWE7IEDB8Lb2xtxcXGKY927d0dQUBBiYmI0vvbnn3+OQ4cO4YsvvqizXWRkJFatWlXrOHOYiIiI9EdDcpj0YkquLpWVlUhLS8OIESOUjo8YMQIpKSlq3qWaptNxYWFhKC4uVnzl5eU16D5ERESkX/RiSq4uBQUFkEgkcHBwUDru4OCAW7duaXyd4uJipKam4ttvv623rYmJCUxMTBrcVyIiItJPeh8wyYlEykuvBUGodawu1tbWuH37tra7RURERK2A3gdMdnZ2MDQ0rDWalJ+fX2vUiai1Yd0mIqLmofcBk7GxMXx8fHDkyBEEBwcrjh85cgTjxo3TYc+ImlBkJDILyjHDNRDX7lQoDrt0NJVt7mtnBkRG6q5/REStjF4ETGVlZcjKylK8zsnJQXp6OmxsbODm5oZFixZh6tSp6N+/PwYPHoyPP/4Yubm5eOONN3TYa6Kmk1lQDs+t67A4uARum2Lg5WCJjNulyF0QBs/EeGTOWwJPXXeSiKgV0YuA6fz58wgICFC8XrRoEQBg2rRp2LVrFyZNmoTCwkJERUXh5s2b6NWrFw4ePAh3d3dddZmoyUikAma4BmJxcAmCEuOBfs5AeDi8E2LhnRiPpOA5WO8WiGSpwOk5IiIt0bs6TC0R95Kj5nQ6uxAvbz+D/XOHwDshFoiIAIyNgcpKICoKaSGhmBCXgt2vDeI2K0REdWhTdZiI2pr80nsAAC8HSyA8/GGwZGwMhIfDy9FSqR0RETUeAyYiPWNv2R4AkHG7FIiOfhgsVVYC0dHIuFWq1I6IiBqPARORnvH1sIFLR1PkLgiTTcdFRQFisezPiAjkLQyDq40pfD1sdN1VIqJWQy+SvonoIUMDkax0wIMEb9eQUHiJq5AREoq8C9cRlBgPBysT/PCbF2szERFpCZO+tYBJ39TsVNRhWvDzlzAxMQIAiMX3sWnoFACszUREpE5DPr85wkSkjyIj4QngRPVK37eOYHDCZiQFz4Hbpk14jbWZiIi0hgETkR4zNBBh8OO2kEgFDOs+FouD77E2ExFRE2DSN1ErkJpThGt3KuC2KUaR/A0TE0VSuOvGGOQVVSA1p0jXXSUi0ksMmIhaAdZmIiJqWgyYiFoB1mYiImpaDJiIWgHWZiIialpM+iZqBTSpzdTTyQqGBs/ouqtERHqJARORHpBULx+gphilp50ZMuctwYeugbgWl6I47jroJfR0spLVYSIiokfCwpVawMKV1GRUFKgE6i5GqS640iToIiJqS1i4kqiVyCwoh+fWdVgcXAK3TTHw0qAYpbw2k8IjBF1ERKSMARNRCyWRCpjhGojFwSWNKkb5KEEXEREp4yo5ohZKG8Uo5UFXUvAcBCXGwzshFuYm7eCdEIugB0HXTLdASKScmSciqgsDJqIWShvFKDUNunadysF36ddxOruQwRMRkQqckiNqoaoXo/ROiK1djDIkVKmdKrWCrtWrH15HIoHp2n8BVv6IPnBZ8R7mNhER1caAiaiFql6M0jsxXjZCFB4uq+QdEYG8C9fhOuilOotR1hV0Ff7nKHqknsJbQ2+i09rVmODjwtwmIiI1OCVH1ELJi1HKc43SQkJRJq5CWkioIifp09wf60z4VlcBXLpqFWxTTyHVoy/e+flLvHr0M+Y2ERHVgSNMRC1YY4tRqqsA/u2zISg4chXv/PwlCn39YLtyJfD++7IRqKgouIaEIi8uBak5RcolCoiI2igGTEQtWWQkPAGcUFl0UrNtTtQFXfB7GaOfdER3ewsg/ZxyQrm4CkDdCeVERG0JAyYiPVCrGGVDqAi6CkrFiD5wGRXvrgDqSCgvKBXju/TrrAxORG0eAyaiNqJ60CWRCtiZ8qfahPJ2nyfBcGI0V88RET3AgImoDVKX25QREop2nyehT+Yv2P9NOCRHj7IyOBERGDARtVnqcpsMJ0Zj/zfh6JP5i2y67hG2YyEiam0YMBG1VXXkNkmOHpUFSxERD4tdcvUcEbVhrMNE1MbJc5vG9XWGnaUJgMZtx0JE1BoxYCIiheqVwREdXXv13K1SpXZERG0FAyYiUlBXGVy+cW/ewjC42pjWuR0LEVFrxBwmIlKoa/Vc3oXrCEqMR08nK42LZhIRtRYMmIhISWO3YyEiao1EgiBwd81GKikpgbW1NYqLi2FlZaXr7hBphUTldiwsJUBErUdDPr85wkREKqnajqUpgygGaETUkjFgIqL6RUYis6AcM1wDce1OheKwVrZLacprExFpCQMmIqpXZkE5PLeuw+LgErhtitHqdilNeW0iIm15pLICa9euBQD89ttvuH//vlY7REQti0QqYIZrIJKC5yAoMR7eCbEwN2kH74RYBD1YTTfTLRASacPTIZvy2kRE2vRII0xDhw4FAERGRuLy5cswMjJCr1698OSTT+LJJ5/EgAED4ODgoNWOEpFupOYU4dqdCrhtigH6OWt1u5SmvDYRkTY90giTn58fAGD//v24fPkyzpw5gwULFsDe3h5Hjx7F6NGjER4ertWOEpFuyLdBaYrtUpry2kRE2qSVKTkzMzP4+vpi1qxZ2LRpE9LS0nDw4EGtdpSIdKMpt0vhVixEpC80DpgkEgm+/fZblJaWKk3J9e7dG71798Yrr7yCmJgY/PDDD7h9+zbOnDnTZJ0moubTlNulcCsWItIXGucwGRoa4tVXX8XFixeVpuQAoLy8HH/88Qd+//13HD16FCtXrsTo0aMRHR3dNL0mombTlNulcCsWItIXDUr69vX1RU5ODh577DGl4/IpOV9fX8UxHx8fBkxErURTbpfCrViISB80aGuUb7/9Fh988AG++eYbuLq61tn2/v37MDIyanQH9QG3RqG2oqHVuFW1B6DyGqz0TUTNrSGf3w0KmAwMZClPFhYWGDt2LPz9/dGvXz88+eSTMDY2blyv9RgDJqIaVFTvXvDzlzAxkf0jSiy+j01DpwBgRW8i0p0m20suJycH6enp+PXXX5Geno6YmBj8+eefMDQ0xBNPPIHffvutUR0notZBVfXuP948goE7NwMAzs6Yj9dWjWRFbyLSGw0qK+Du7o5x48YhIiIC+/fvR3Z2Nv755x8cO3YMc+bMaao+alVwcDA6duyIiRMn1jr3ww8/wMvLC56entixY4cOekek/1RV725vZIjfrhUr2vx+rRjtjQxZ0ZuI9EaDpuRag+PHj6OsrAwJCQn45ptvFMerqqrQo0cPHD9+HFZWVvD29sbZs2dhY1P/cmZOyRE9dDq7EC9vP4P9c4fAOyEWiIiA1MgYBvcrcWNRGAQAzhtiFMcQFYW0kFBMiEvB7tcGsaI3ETWbhnx+P1LhSn0WEBAAS0vLWsdTU1PRs2dPODs7w9LSEqNHj8ahQ4d00EMi/aaqerfB/UqIDdvB+l9RsP5XFMSG7WTBEit6E5GeaFEB08mTJzFmzBg4OTlBJBIhKSmpVptt27bBw8MD7du3h4+PD3766Set3PvGjRtwdnZWvHZxccH169e1cm2itkRV9W6pkTFMJFUoXh6B4uURMJFUQWrEit5EpD9aVMB09+5d9OnTB7GxsSrP7927FwsWLMCKFStw4cIFPPXUUxg1ahRyc3MVbXx8fNCrV69aXzdu3Kjz3qpmJkUi1UuaxWIxSkpKlL6ISEZV9W7h3j1sf246nDbEwHlDDHY8Nx3CvXus6E1EeqNBq+Sa2qhRozBq1Ci15zds2IBZs2Zh9uzZAIBNmzbh0KFDiIuLQ0xMDAAgLS3tke7t7OysNKJ07do1DBw4UGXbmJgYrFq16pHuQ9TaqazefV+C3i7WijZPulij4r6EFb2JSG+0qICpLpWVlUhLS8O7776rdHzEiBFISUlR8y7N+fr64o8//sD169dhZWWFgwcPIiIiQmXbsLAwLFq0SPG6pKSk3kKeRK1ZzaKTA21Na1XvXnA1H788EwIAEF/Nx6SVshxBVvQmIn2gNwFTQUEBJBIJHBwclI47ODjg1q1bGl9n5MiR+OWXX3D37l24uLggMTERAwYMQLt27bB+/XoEBARAKpVi6dKlsLVVvVrHxMQEJiYmjXoeolZBRYFKAHDpOAw7C3/EyfITOPva27JA6rVBSpW+NytV9ObIEhG1bHoTMMnVzCsSBEFtrpEqda18Gzt2LMaOHfvIfSNqa1QVqKxZjFJVmQCWDiAifdOikr7rYmdnB0NDw1qjSfn5+bVGnYio6akqUGlu0o7FKImoVdKbgMnY2Bg+Pj44cuSI0vEjR45gyJAhOuoVUduVmlOEa3cq4LYpRrHaDSYmipVxrhtjkFdUgdScIl13lYio0VrUlFxZWRmysrIUr+V719nY2MDNzQ2LFi3C1KlT0b9/fwwePBgff/wxcnNz8cYbb+iw10RtU60ClatXy+oqyYtRiquU2hER6bMWFTCdP38eAQEBitfylWjTpk3Drl27MGnSJBQWFiIqKgo3b95Er169cPDgQbi7u+uqy0RtVvUCld4JsQ+DJXkxypBQpXZERPqsRQVM/v7+KgtIVjd37lzMnTu3mXpEROpUL1DpnRgvm5YLD5dV946IQN6F63Ad9FKTFaOsWcpAttpO8wUgREQN0aICJiLSHyoLVIqrmr4YpdpSBqay/tiZAZGR2r0nEbV5DJiI6JF52pnVKlAJNG0xSk1KGXhq/a5E1NYxYCKiRxcZCU8AJ1ROj2m/GKW8lMHi4BIEJcYD/ZyB8HB4J8TC+8FI13q3QCRLBU7PEZFWMWAiokYzNBA1SzFKpVIG/ZxlJQzkq/OiouAaEoq8uBSk5hSxOCYRaZXe1GEiorZJIhVwOrsQ36Vfx6msAgDVShnIV+XJSxk4WgJgKQMi0j6OMBFRy6QmuRsA/nhzMQbm/MpSBkTUbBgwEVGLpCq5+9LNEqRMn4+BJ76QNWrmUgZE1HYxYCKiFkddcveAz7diwINgafOwVzFk6jx0b45SBkTU5jFgIqIWp87kbn9/nPXog432z2Hj/51WvKcpSxkQETFgIqIWp8596o4fR09xFbDyEEIDHoeng2WTljIgIgIYMBFRC6TpPnV+XTsplQ/gdilE1FQYMBFRi9Pgfeq4XQoRNTEGTETU4jR0nzpul0JETY0BExG1SJruU8ftUoioOTBgIqKWScN96rhdChE1BwZMRNSi1bdPXZ0r6sLD4SWuUmpXHZPEiUhTDJiISK9puqJOabsUJokTUQMxYCIivabJijqXgS9CKgj4Lv067C3bw/bvcnTbxiRxItIcAyYi0mt1ragz/+o7BCXG42ZxBabceenhe6z9sd/zGJPEiUhjDJiISO+pW1G3zKknnriShoAbFzFw7hB4OVji21+u4e+l76FP5i8o9PWDLZPEiUgDDJiISP+pWFFnZ26CZR1N0dnaVDaSlBCr2MC3+89fYn/Q69g4eBJOpo+CSMMkcSJquwx03QEiIm2Rr6gb19cZBgaih+UGoqJk5QZMTNA97kOsOBt3rwAAIABJREFUHzoF4ndXIOjfnz4MluRJ4rdKAdRIEieiNo8BExG1SrXKDTwIigRjYyS+MBMma97HOz9/ictvLgbEYkVQlbcwDK42pg+3XSEiAqfkiKiVUlduQFRZie/2R8A29RTWD50Cu1fnwbWObVeIiAAGTETUStVVbsA2IgK/enpj21OvQPL9Raz8/iKA2tuuEBHJMWAiolZJkw18/1N8HIWL3lW77QoRkRwDJiJqterbwLebnRnA0gFEpAGRIAiCrjuh70pKSmBtbY3i4mJYWVnpujtEVAP3jCMiVRry+c0RJiJq9erbwJeIqD4sK0BERERUDwZMRERERPXglBwRUTXMdyIiVRgwEREBQGQkMgvKMcM1ENfuVCgOu3Q0lZUnsDMDIiN11z8i0ikGTEREADILyuG5dR0WB5fAbVMMvBwskXG7FLkLwuCZGI/MeUvgqetOEpHOMIeJiNo8iVTADNdAJAXPQVBiPLwTYmFu0g7eCbEIelD4cqZbICRSVmEhaqsYMBFRm5eaU4RrdyrgtilGsQkvTExkf0ZFwXVjDPKKKpCaU6TrrhKRjjBgIqI2L7/0HgDAy8FStt+csbFis16Eh8PL0RIAcCrrb3yXfh2nsws52kTUxjCHiYjaPHvL9gCAjNul8E6IfRgsVVYCAQG46NEHsH8OscezFe9hMjhR28IRJiJq83w9bODS0RS5C8IU03AQi2V/Jidj4M7NWHhmD755YzAurhqJ/XOHYPHZffDcug6ZBeW67j4RNQMGTETU5hkaiLAz70dFgndaSCjKxFVIfXUeNg97FQAw/8QX6P/5ViaDE7VR3HxXC7j5LlEroKYOEwDszT+CgTm/AsnJD6fqoqKQFhKKCXEp2P3aIO5VR6SHuPkuEVFDRUbCE8CJapW+M2+XIfZ4FnrFfQiYtJOtnKueDC6uAiBLBmdlcKLWjQETEVE1hgYixWjR6exCxB7PYjI4ETGHiYhIHSaDE5EcAyYiIjWYDE5Eckz61gImfRO1YkwGJ2q1mPRNRKQtjUgGl1cQJyL91+am5IKDg9GxY0dMnDhR6XhpaSkGDBiAvn374sknn8T27dt11EMiaonkyeDj+jrDr6sdAFllcERHKyeDR0cj41YpgIcVxIlI/7W5gOntt9/GZ599Vuu4mZkZTpw4gfT0dJw9exYxMTEoLCzUQQ+JqKWrMxk8IgJ5C8PgamMKXw8bpfdJpAJOZxdyPzoiPdTmpuQCAgKQnJxc67ihoSHMzMwAAPfu3YNEIgHTu4hIFXkyuOeDBG/XkFB4iauQERKKvAvXEZQYj55OVjA0eEb2BjV5UCxBQKQ/WtQI08mTJzFmzBg4OTlBJBIhKSmpVptt27bBw8MD7du3h4+PD3766Set3f+ff/5Bnz594OLigqVLl8LOzk5r1yai1sXTzgyZ85bgw4EvYUJcCnqtPIQJcSlYP+glZM5bIguCHsgsKIfn1nVYfHYf9s8dwhIERHqoRY0w3b17F3369MGMGTMwYcKEWuf37t2LBQsWYNu2bfDz80N8fDxGjRqFS5cuwc3NDQDg4+MDsVhc672HDx+Gk5NTnffv0KEDfv31V9y+fRvjx4/HxIkT4eDgoJ2HI6LWRUUy+MNK388omkmkAma4BmJxcAmCEuOBfs5AeDi8E2Lh/WCEar1bIJKlAiuEE7VgLSpgGjVqFEaNGqX2/IYNGzBr1izMnj0bALBp0yYcOnQIcXFxiImJAQCkpaU1uh8ODg7o3bs3Tp48iRdffLHWebFYrBSUlZSUNPqeRKSfqlcGVyU1pwjX7lTAbVOMLFiKiABWr1aUIHANCUVeXApSc4pYgoCoBWtRU3J1qaysRFpaGkaMGKF0fMSIEUhJSWn09W/fvq0IfEpKSnDy5El4eXmpbBsTEwNra2vFl6ura6PvT0StR/Xk7lNZBQAALwdLIDz84Wo6eQkCR0sALEFA1NK1qBGmuhQUFEAikdSaInNwcMCtW7c0vs7IkSPxyy+/4O7/t3fnYVFW7R/Av8OwKaKJ5IICkhruCyCKYLlkYmpJvS1elJpaueASvr6KJiiuLZqZSWqLVpb93hTfMlPJXUlNkVxDLFxTCRdA0wFm7t8fwzzMsA0oMszw/VwXl/E8Z57nHIeRu3Pu5z63b6NJkyaIj49H586dcfHiRYwYMQIiAhFBREQE2rdvX+w1oqKiEBkZqXyflZXFoImISi1yeXz0v/VFLguXIBgSAYAlCIiqOqsJmAxUKtM1fhEpcqw0W7ZsKfa4v78/kpOTy3QNJycnODk5lfmeRFQ9KMndYVnwWjwfvg1ccfJyFhKHTUCXXV/pG8XG6meaZs/WlyA4cgmeXV8oUoKAiKoWqwmY3N3doVari8wmpaenMzGbiCyupOTuzl9+hM75wdIHj7+Mbq+MRavSShAQUZVkNQGTo6Mj/P39kZCQgLCwMOV4QkICnnnmGQv2jIjITHJ3jx444NMB79fvg/c//kV5jWfXF9DGozaa1auBX/64VuhJOz4xR1SVVKmA6datWzhz5ozyfVpaGpKTk+Hm5gYvLy9ERkbilVdeQUBAAIKCgrBixQqcP38eo0aNsmCviYgKkraV5G5DsOToCOzYgTaaPCBmCyJ6NkOLBq6o7+qMLl8swR8AHqv5OC6u3K9ciwUtiaqeKhUwHTp0CD179lS+NyRWDx06FKtWrcKLL76Ia9euITY2FpcvX0bbtm2xadMmeHt7W6rLREQACpK2U65mw2/10hKTu4ObP6yUD0i9dqdIzlPK1WycnxiFFvHL9QUwLTYiIjJWpQKmHj16mN2OZMyYMRgzZkwl9YiIqGyM95fzi19uNrn7XgtaaostlMnlO6IHrUoFTERE1qq8+8uVu6Al96MjsigGTEREFUTZX84zFBfjCgrqGpK7jfeXKzXnacYM+GryTNoVV7KAy3dElYcBExFRRSnj/nJA2XOe6rs6cz86oiqAARMRUQUzt78cUL6cJ+5HR2R5DJiIiCygPDlP5V2+I6KKx4CJiMhCyprzVJ7lOyJ6MBgwERFZShlznspbsoCIKh4DJiIiCzOX81TekgVEVPEYMBERWYHylCwgoorHgImIqIoyqer9yjgE+rhhF2C2ZAERVTwGTEREVY2Zqt5BrOpNVOkYMBERVTGs6k1U9dhZugNERFTAUNV7Q9gbGBS/HH6rl8LFyR5+q5diUH7S93CvUGh1pW9UTkQViwETEVEVYlLVOzZWX9XbyUn/Z2wsPN+fjwvX7+Bg2nVLd5WoWmHARERUhRSp6m0oUGmo6t3Q1aQdEVUOBkxERFWIcVVvzJ5dtKr3lWyTdkRUORgwERFVIcZVvQ3LcNBolOW5C29GwdOtBqt6E1UyPiVHRFSFsKo3UdXEgImIqIphVW+iqkclInw29T5lZWWhTp06yMzMRO3atS3dHSKyEdpiN+VVWbpbRDajPL+/OcNERFRFmduUtywYdBFVDAZMRES2yMz2Ki24vQpRuTBgIiKyQdxehahiMWAiIrIBxktv7i5OmOIZin+HZWFQ/HKgU2Ngxgz4rV4Kv/yn7xZ6hWKnTrg8R1RGDJiIiKxZCUtvANDoIWegRw99Pac5c/TFL2Nj4TkkAhfiEnEw7fp950gRVRcMmIiIrFhxS2+f7P0TeTNj0WXvGv3SW2Ki6fYqmjwA3F6FqDxY6ZuIyEppdYJXPUOxIewNDIpfDr/VS+HiZI+wjZ9j0t41WDtwJHad/pvbqxBVAAZMRERW6mDadVy8cQdei+crW6fAyQleixdgZZ9hgAAjE1bh/MSp3F6F6D4xYCIislKGJTXfBq7AjBkFs0iOjujx6MN4aeMnWBgSjvX9X8UtTR4OD4lQZqM+O7+ZCd9E5cAcJiIiK2VYUku5mg2/1UtNlt5anPgVB16dgA/r9wF+TsXin1MBcHsVonvFgImIyEoF+rihSd0aOD8xCn7xy/VLbjNmALNnA9HRuFzXF018nfH2cx2QcUtjVOm7+I17WRWcqGQMmIiIrJTaTqWv2p1fW8lzSAR8NXlIGRKBC0cuYVD8cv1sUvPepV+IVcGJzGLARERkxVq410Tq2Ml4zzMUF+MSlePlWXpjVXAi8xgwERFZs5kz0QLArmKX04pfejNmKE3AquBEpWPARERkA9R2qnuq2m1SmqBTY1YFJyoBywoQEVVjpZUmwIwZ8G3oatKOqLpiwEREVI0ZlybA7NmsCk5UAgZMRETVmHFpAkRH60sTsCo4URHMYSIiqsbKWpqgLAnkRLaMARMRUTVXEaUJiGydSkTE0p2wdllZWahTpw4yMzNRu3ZtS3eHiOiesNI3VTfl+f3NGSYiIgJw76UJiKoDJn0TERERmcEZJiIiKhWX6ogYMBERUUlK2JQ3+tD/4ZnsP1A3tDcODBlvGkjNnQNotdysl2wOAyYiIipWSZvy1gidiXqnDuHg37cxWBOotJ+e9B1eS1ilr+FEZGMYMBERVTNlWWIrbVNenDqEfV7tEZyWjH05+/DQ/FhkTouGR8IqLOoejtZhIxF6D/e8l7ZElYUBExFRdVHCEluTujX0xSvdaypLaaVtyruyzzBsCXsN+79ehknvzwc+WgiXnBzoZs3CSe9+iN90Cn1aN9QHOWbu2axeDWVZz++zD3BHCwz3Mt8/ospW7QKmsLAw7Ny5E71798Z3331nci4tLQ3Dhw/H1atXoVarsX//fri4uFiop0REFaukJbbzE6PQIn45UsdORov8tkU25c0PlnQOjpjr9y+s6fMows8NxsQD/4U6f/85u+hojD53A8/FJeJg2nUENatX4j1r9u2DFr8fRlyvIXhbsx8AMO7IX5i0dw0+arEN9fs/gdoL5pq0PT3m37j2xzXmTJFFVLuyAuPHj8cXX3xR7Llhw4YhNjYWJ0+exK5du+Dk5FTJvSMiejAMS2wbwt7AoPjl8Fu9FC5O9vBbvRSD8rdFGe4VCq1OX8u4pE157XJzMG7fNxAA4/Z9A3Wu6Wa9vg1dAegDrtLu2fL3wwCAHK0O340KwtGYJ/Ft6DAc9OmIDqlJcEzcV6Ttj8euYPDK/ZiwNhmDV+7HZ6HD9TNfanXl/4VS9SPV0I4dO+S5554zOXb8+HHp3bv3PV0vMzNTAEhmZmZFdI+IqMIlnskQ7ykb5fC56yKxsSKAiKOj/s/YWDl09rp4T9koiWcyREQkT6uT4AXbJD7sDaWNiMi5iVNFADnq6y8CiHbWLP0N8q95MTJKuU5p91zRZ5h8M2CkCCDnJk6VxDMZ8l5IuAggJ1sFiACiy297euxkWZh/Lmnkm3Lrbq5ciowSAWRh93D56dhflvprJStXnt/fVWqGaffu3Rg4cCA8PDygUqmwYcOGIm2WLVsGHx8fODs7w9/fH3v27KmQe6empqJWrVp4+umn4efnh3nz5lXIdYmIqoIiS2yGWSFHR2DGDJOZIaBgU17D7NPhIRG4pcnDlfGT8UvT9miXchgZgcGwi47W32DGDOhmzULjRfPxVtJ3CPRxK/GehmU9z/fnY2FIOLwWL0BgKw9M2rsGOdEz8c9PW6FR20OVkwNxdMSrnqE4Onw8FoaEo9Mn78Oltgs8Fs3X50yNmIC5m04pM2NED0qVymG6ffs2OnTogFdffRXPPfdckfPffvstJk6ciGXLliE4OBjLly9Hv379cPLkSXh5eQEA/P39odFoirx269at8PDwKPHeubm52LNnD5KTk1G/fn2Ehoaic+fO6NOnT5G2Go3G5B5ZWVn3MlwiokpjvMTmt3ppQbCUv5SWMiTCpB1Q8qa8UT7todMB1xq1QuNzN+Db0BUpV7IR590PrbufxtPN9U+1lXRPu5z8Zb2RXfBh8GC8eeC/UOfmQKO2x4kR49H20yVw0uZB6+AIdU4OwjZ+Bo/35iEqv62d4To6HWYdjUewQzBW7UuDu6sTc5vowamEGa97AkDi4+NNjgUGBsqoUaNMjrVs2VKmTp1armsXtySXmJgoffv2Vb5/55135J133in29TExMQKgyBeX5Iioqippic2wVBYf9oaEvL1N8rS6Yl+beCZDNhy5KIlnMiRPq5Ofjv0lwQu2ifeUjcpXyNvbTJbHzC3rfTNwpKzsM8xk+c2wHPdeSLgknsmQk6P/LQLIukGvK22VpcRevZS2xv1YYWhnGCNRCcqzJGc1AZNGoxG1Wi3r1683aTd+/Hh57LHHynXt4gKm3Nxc6dixo1y/fl20Wq0MGDBAfvjhh2Jff/fuXcnMzFS+Lly4wICJiKq802MnK8HRobPXJfturhw6e10JaE6PnVyu6xUXSBU+ljKm6D0P/HlNFj/+sj6oMQpsMgKDRQA54NNROs/ZKjf/yZFV+9KU3KbiAr29Xu2Z20T3rDwBU5VakitNRkYGtFotGjRoYHK8QYMGuHLlSpmv07dvXyQlJeH27dto0qQJ4uPj0blzZ9jb22PevHl47LHHICJ48sknMWDAgGKv4eTkxCfoiMjqlLTE5tn1BbTxqK2vc1QOajsVgprV039T0jYqKVfh2L4LrmfdwUSje05R61Nof2/pj9tDIuCryYMmKBi/3biDwNQkvLh5NTpk5wAAJqj0RSulRw+oZswAAGinv4XP9vyJ1xJW4VCzTgj45H3gi4/gkpMD6dEDj9SthTk/noSrswMybmm4VEf3zWoCJgOVyrTaq4gUOVaaLVu2lHiuX79+6Nev3z33jYioSps5Ey0A7Cq2knav+7p0yTWefkDTbQcQMqY7vnmtq3LPLk4HkdpqMkYaB2/Oj8FzZF98eu4nvKBWoflLHVHf1RnumduxCOE4OWwCRufnTK07fBFz/f6FrDu5eKZdA+DTE0qOlKpXLwyKjsYff99C+I3BSh+5dQvdD6sJmNzd3aFWq4vMJqWnpxeZdSIiopKZzAxVgNK2UfHLf8puoXc/7DTe4mTWLLPBm6fhBnEL0fr4Zaz/8RSeM5qlAoBmSxag+cbPTZLYU69m44fu4Zi0Zw26NauH9ssXmt26hcgcqwmYHB0d4e/vj4SEBISFhSnHExIS8Mwzz1iwZ0RE1Vtp26ggNhaeQyJwwaj6t7GyBm+hbRuhT+uGSnCVka3B7B9PIfCrZcCi+fpZoxkzoIuNRYuYGHgMGImFIeGYtHoJ8M3HJW/dko/715E5VSpgunXrFs6cOaN8n5aWhuTkZLi5ucHLywuRkZF45ZVXEBAQgKCgIKxYsQLnz5/HqFGjLNhrIqLqraRtVJQaT5o8k3b3yji40uoEutmz4ZGwCrpZs2CXn9t0IHwsEhNOY9LGT/BJn2EQR0eoStm6pTz76xkwuKqeqlTAdOjQIfTs2VP5PjIyEgAwdOhQrFq1Ci+++CKuXbuG2NhYXL58GW3btsWmTZvg7e1tqS4TEVV791Lj6X6p7VTo2dwNi+6G46R3PyW36UDaNXwYrM9bGnbztBIsKVu3/CcKQEHwVp799UpMbD/0f3gm+w/U6/dE0WTy2bOZZG4jqlTA1KNHD4iUXq11zJgxGDNmTCX1iIiIzAn0cUOTujVwfmIU/OKXK8tjmD0biI7GhSOX4Nn1BQT6uFXofZsvKzm3KaSFO+p9vqagLz17AtHRuHnjDuAYjNSrt7AvNQNTPEOxstWOknOvvEKxUydQ26lKDK5qhM5EvVOHcA1APePAKH/8TDK3EQ+6xkF1wL3kiKi6q+gaT+VhXPtp7+m/lcKVyj53IgV72RUqdKnUeMovglnS/nrmCn8aCm6enTBVNhy5qBTnLFI8Myam5IKasbH681RpbLIOExERVV0VXeOpPAonjjcoZqludchg3A1JwaS9a/DirTMYNasvjo6ahKC9a7AwJBxtPngXof5NC5YTtVq0/WwJgADsO/M39p35Gxdv3EGTRfOAzNQiie3bQwZjU/RMTPpgARoufQ9O2jys7DMMnoWfyFOrgeho6ERwIHyskgfVdeS/oNq5s+hs1MyZwJ490IWE4MCQ8aZ5U6wpVakYMBER0f17gDWeyqukpTrnHuEY2KERHv3oPaC2C4JycrAwJBz7w8egTuQ009yrPXvgtH07xoWE40MU1HJKem0SArZv1wc++e03h43Ee2uSIMGDMeHA/8EpNxc6R0ccCB+DeWuSEBfuh9C2jfQXmDEDqVez0SImBokJp/Fh8GCM2/cNgvbuLH4we/YA27cj6dwNDNYEKoeLrSk1c6a+X/kJ8CaYS3X/KmHGy+ZxSY6IqOoxLNW9u/mUeE/ZKL+mXdOfyF920zk6SvCCbfLNgJEigH4ZTaTItisX34yS7aeuKst3e73bmyzfregzTJ5btq9geS//uHbWLBmx6qDJHn0/HftLmk7dqCztGfbQiw97QxZ2Dy9++5f8+12KjCp9+xfDsmPhJb+SjpdH/lJicdvhWPNSok3sJWdNGDAREVVdG45cFO8pG+XW3dyC4CE/UDHsX/deSLi8n5BSZK+7v/PPGwKbE/m5Siv7DJM8rU7JVTrq62+aN5V/n8xu3YvkQY1YdVDfTq1W+qLV6mTEqoNK/lWeQ0EwprQ3E4wZ3/fcRDO5VOWVf90VfYbZ1EbH5fn9rRIx81gamZWVlYU6deogMzMTtWvXtnR3iIjIyC9/XMPglfuRmLsPHkZFLpWn+Dp1Rfcn3yryum/TE9DF6yFg/nwgJwc6OzXsdFqsHTASU9sMwmdDO+PopZvoPOJ5BJ8/imuBwah3YG/BBXr3BrZvxz6v9sj43ybUd3U27QdQsAQYG4uPQgbj3S0pSHlvEJy0edCo7eH77w34T19fjOnZHHByKlg2jIrCpewcBDsEY0b/VnB3dUJ9V2fcuJ2D+s/2R8AfR6BR2xfkUi2aV7AsWEhZ6kptPn4ZJ8dMRuSeNfgrMgp15sXqq6cvmq+vnr7s3RKvX5WV5/c3c5iIiMimBfq4YXrSd0WKXGLGDOhE4BkTg+nu36H1sneRcUuD1Ku3sHTHGbSNew94Z74SpNjl5OBaYLC+ztPNOxi++lcAwETPNvCqVxOe/Z4wvfG2bcgKfgzBiXtweOl7+N/AERi37xt47F2jP18ocNN0T8E4AZy0eRBHRzjl5GDcvm/wDgajz7oVaFEox6pxfo7VbKMcq3H7vsFTfxyBqNVw0uaVnEsFlLlop1YnmPPjKbQcMQG6Jx6FR0wMsHSh2erpNueBz3dVA1ySIyKq2lJHR8rC7uEyYtVBk7IHI1YdlIXdwyV1dKTSNvFMhnhP2ajkChXOKbrwZpR4T9ko727+Xfae/rtgmc14aUxEtFqdDFyyWxZ11y/v3VXbK6UNjMss5Gl1BUtbRrlU5pb7DDlWSSPflMx/ciSu91CT40quVo8eEh/2RpHlO0MpiFMt/eXwuety626uHD5XtBSE4e/j8Lnr+hcaSi84OoqImJRfuCcWLLXAHKZKxoCJiKjq++nYXxK8YJtJDk7I29tME6fFNIAxqeUkouQSGXKYDNdtOnVjkWBs4JLd4j1lowxcske0+UFGnspOlvUaIt5TNir3TTyToSSMH3qko3LdPK1Ofm3WSQmCDAGJcf8M57X5OU+ZQd1N+1eo/tQnu/8oUq+quKBwfdjrEjBnq6xPuijvbv69xBwwiY2V7Lu54j1lo2w4crHgL6o8QdCDTFY3g3WYiIiICim8gW9J+TolbbuSciUbcd790Lr7aTzdvOB1oW0bIS7cD3MKlTFQ26nQvnEdbLixHXb5y2nqnBzcydXCvZYjov93HEHN3HEg7RrUotPXg4p7V7mu2k6Fh57siX2btPi1SRsg7RraNamDdYcvYq7fv5B1JxfPtGsAfHoCdjk5yFPZofYve3ApMgpzHILRJu06ggzLj9HRmLR3DRZOA2bnlzKYtHcNDrw6AV186pnUlboWGIz5/s/j72wN3vw2WRnP9aDucDmyv8hSYubNO4BDsOnWN/n1pgCYljkorvq5UR+V743bFVcmwQKY9F0BmPRNRGR7Nh+/jDk/njLJ7/F0q4HpT7UqNsHZOHk6I1uD2T+eKjHR/OPeQ7Eg4HmT1yvJ3YV8tOMM3t2SUuT44hc7YNDGz/WFMB0cYZebg6xu3WG3fTvaxmzBBy91xDMdG+sTtkdPRtfzR9Ht/DFlQ+KFIeH4MHgwPn7ZTynaaUhsXz/odUT6Po13nmsHn4drIa9nTwSdPYqMwGC4GyW2S8+eUO3ciU/6DMOrmz8zDT7zk94Lj73EIMhw3igR/kEHS0z6JiIiuk9lnZEyMK44/r/kSwUJ3sa/+PP/HBUdjduaPGT+eyr6tm6IKeuPIun8Deh0Ajuj6+t0giPnb6BJXWe8/VwHZNzSKMFY4FfLgPxgTKa/hZWhw/FawipcjIoGHPUzPoaE7bt9h2K9gx32RD8FVU4OdA6O+DB4MDo3rYuL+UU7JT+x/VSrADy7YQXSQm7D87Wu6PzlR8DZo9jn1R4jnpyOL89eR6tGtZFyJRsX6vpiEHYi7OZpk7+nLms+gt327dD17Am76GhoZ8+BOjdH/31phTUNwVJ+tXXMnl11CnE+sIXBaoQ5TEREZCzxTIa8HzxYnzhejIuRUfJ+8GAlL6mkPKgRqw5K06kbTfKsSsqx+unYX0rxy7jeQ+TmPzmyal+akq9lSPI2LrhpXLTTOJfKsDeezihXaen2VJP8L0MO2K6XI0z26DNcY9fLERK8YJuS7J6rsit+X0FDrlLh/fwM3z/A3CbmMBEREVlQoI8bJg8cieMNXbGimFmjmPZhSGmQjXE+bgBKzoPydKtRpBxASTlWD7s6Y+fzb0AEUN/NRYdZW5XX/JS1Ay3iFposj72WnzO0MCQc6v6vwuPGP/gweDB83F3w7IYV+uU5w2zPjBkYqsnDu1tSENGzGVo0cFXqPg29cQfvh+ViUvxyRB78L1Q5OVjz1AhMbxyKubu/Vsok2OfPXrX66F38ckuD9In/QafVS+G1eAHQq1fxy3e9elWZ3CbmMFUA5jAREVFhm49fxujLN+L9AAAQsUlEQVQ1Sejdsj5G92hekDi+8wy2/Z5etC4SylZE0vj6xeVYRYW2RF0XJyWX6uY0fcJ3kSCjZ09g504ll8mgfm0n/PS/mah3cJ9JPtHhIRF4Li4R37zWFUHN6kGrEzz+7g60bOiKFa8EwK6Gs7K01z32J4RvXY3R27/Q176KjlaCnV+atkfQ2aNKYc1DzToh4I8jRftnHDRt3/5AcpuYw0RERGRh5Zk1MjDOgyrL9c3lWGl1gtVOamwIewNPT38Ldkav123bju//NRp1s/7BmpFdkJ51F/M2/Y6oQ9/pg6XCFdGPXIJn1xcQmD8rdjDtOi7euIMlgzvBbu4cJfdIlZODt5e9ieDzR7EwJBzdwsciCMDmsJE4mZCCyD1rkKeyUwpr/tXWH4s8WqN12EiEGg8wPyjS5eUBe/bCLj/3Sqa/BXVZ34QKxICJiIjoASlv4nh5mQuw1HYqeCxegNFrkrDxy0NFZ7p8ByIu3A/Bzd0BAG0/W4IW8cuxIewNeA6JgK8mDylDInDhyCUMil+ONh61obbrBQBIz76rf82nS4DYmUqAdWrMZATHvQfNYz3wYdBgNM++a1ItPNfuEhx27YTWQV9mYYCfJ14voVr45rCRuBA5Da/l5uhnpHJzsDJ0eKlbvTwoDJiIiIgeoPLMGj0I5ZnpauFeE6ljJ+M9z1BcNG7b9QW08ait3y4lX31XZ4zb9w0cCy333Zw0FQuPXcak3WswTtcY9V/rqsxG/V/6z3DYpV8G7LbqAwR9vQx20dGYFXkXwQ7BOJh2Xfm7Km7/umvTovHaovlYNCYXqOT96xgwERER2bgyz3TNnIkWAHYVm0vVy6RpoI8bThWz3Bfo44bJA4ajpqMadVWCQB83bDz6l1JmYUPYG9jQ9QVM9HFTgqzG0dEYFxKO9Jc6AtAvJV6InIbIPWugmzULHvmJ3y4L50FXxxmRMTH4JNIB2sK1nx4gBkxERETVQHlmusrStrTlPreaDni78wto37gOOly4CVdnB6hFh7UDRyLKdyDinmpVEOjMmIFLN+9AfeCsUi38YNp13P5Hg78io5RgycAuOhqXMu/i1oGzJjNSDxqfkqsAfEqOiIiqq5Ke1nuqbUP8eOyKyXFnBzsseqEDnmrnoRzT6QSvf3kIv1/JUopzpl69haU7zuDErL5wcSo6t3NLk2dSzfxe8Sk5IiIiqhSlLff9J7SVcvxsxj9YvO001iddQoPaNUySz38+lY56Lo4I/+SAybVXJ54tdruYlCvZAGC6f90DxoCJiIiI7ktJS3iFj/s2rFUk+bxeLUeoAHTyeghjejaHbwNXnLychVc+PYB3tqSgqXvNIjNScTvPwNOthlLioDJwSa4CcEmOiIiobIyLc7q7OGHK+qMFxS+NErg3Hb2MMV8nwdnBDl+O6KLsX1da4c/y4pIcERERVUnGs06//HGtoPhloafdnmrfCJOv+eLdLSl4/uNflOOlFf58kBgwERERkUUYil/6NnAt9vzQbk2L7F9XkYU/y4MBExEREVmEIWk75Wo2/LzqFjlvSO4Obv6wRYt/AjDZVoaIiIio0gT6uKFJ3RpYtuMMdDrTlGpLJXeXhAETERERWYTaToW3+rfCtt/T8fqXh3D43A3c0uTh8LkbeP3LQ9j2ezqmGxe5tCA+JVcB+JQcERHRvSup+OX0p1o90ORuPiVHREREVqPMe91ZEAMmIiIisrjy7HVnCcxhIiIiIjKDARMRERGRGQyYiIiIiMxgwERERERkBgMmIiIiIjMYMBERERGZwYCJiIiIyAwGTERERERmMGAiIiIiMoOVviuAYTu+rKwsC/eEiIiIysrwe7ss2+oyYKoA2dnZAABPT08L94SIiIjKKzs7G3Xq1Cm1jUrKElZRqXQ6Hf766y+4urpCparYjQKzsrLg6emJCxcumN1J2VpxjLaBY7QNHKNt4BjLRkSQnZ0NDw8P2NmVnqXEGaYKYGdnhyZNmjzQe9SuXdtmf+gNOEbbwDHaBo7RNnCM5pmbWTJg0jcRERGRGQyYiIiIiMxQz5w5c6alO0GlU6vV6NGjB+ztbXcFlWO0DRyjbeAYbQPHWLGY9E1ERERkBpfkiIiIiMxgwERERERkBgMmIiIiIjMYMBERERGZwYCpClu2bBl8fHzg7OwMf39/7Nmzx9Jdui+7d+/GwIED4eHhAZVKhQ0bNpicFxHMnDkTHh4eqFGjBnr06IETJ05YqLflN3/+fHTu3Bmurq6oX78+Bg0ahJSUFJM2Go0G48aNg7u7O1xcXPD000/j4sWLFupx+cXFxaF9+/ZKobigoCD89NNPynlrH19x5s+fD5VKhYkTJyrHrH2cM2fOhEqlMvlq2LChct7aP4sGly5dwssvv4x69eqhZs2a6NixIw4fPqyct4VxNm3atMh7qVKpMHbsWADW/7MKAHl5eXjrrbfg4+ODGjVq4JFHHkFsbCx0Op3SplLeS6Eqae3ateLg4CArV66UkydPyoQJE8TFxUXOnTtn6a7ds02bNsn06dNl3bp1AkDi4+NNzi9YsEBcXV1l3bp1cuzYMXnxxRelUaNGkpWVZaEel0/fvn3l888/l+PHj0tycrL0799fvLy85NatW0qbUaNGSePGjSUhIUGSkpKkZ8+e0qFDB8nLy7Ngz8vu+++/lx9//FFSUlIkJSVFpk2bJg4ODnL8+HERsf7xFXbw4EFp2rSptG/fXiZMmKAct/ZxxsTESJs2beTy5cvKV3p6unLe2j+LIiLXr18Xb29vGTZsmBw4cEDS0tLk559/ljNnzihtbGGc6enpJu9jQkKCAJAdO3aIiPX/rIqIzJkzR+rVqycbN26UtLQ0+e9//yu1atWSxYsXK20q471kwFRFBQYGyqhRo0yOtWzZUqZOnWqhHlWswgGTTqeThg0byoIFC5Rjd+/elTp16sjHH39siS7et/T0dAEgu3btEhGRmzdvioODg6xdu1Zpc+nSJbGzs5PNmzdbqpv3rW7duvLJJ5/Y3Piys7OlRYsWkpCQII8//rgSMNnCOGNiYqRDhw7FnrOVz+KUKVMkJCSkxPO2Ms7CJkyYIM2aNROdTmcTP6siIv3795fhw4ebHHv22Wfl5ZdfFpHKey+5JFcF5eTk4PDhw3jyySdNjj/55JNITEy0UK8erLS0NFy5csVkzE5OTnj88cetdsyZmZkAADc3NwDA4cOHkZubazJGDw8PtG3b1irHqNVqsXbtWty+fRtBQUE2N76xY8eif//+eOKJJ0yO28o4U1NT4eHhAR8fH7z00kv4888/AdjOZ/H7779HQEAAnn/+edSvXx+dOnXCypUrlfO2Mk5jOTk5+OqrrzB8+HCoVCqb+VkNCQnBtm3bcPr0aQDAb7/9hr179+Kpp54CUHnvpe2W/7RiGRkZ0Gq1aNCggcnxBg0a4MqVKxbq1YNlGFdxYz537pwlunRfRASRkZEICQlB27ZtAejH6OjoiLp165q0tbb39dixYwgKCsLdu3dRq1YtxMfHo3Xr1khOTraJ8QHA2rVrkZSUhF9//bXIOVt4H7t06YIvvvgCjz76KK5evYo5c+agW7duOHHihM18Fv/880/ExcUhMjIS06ZNw8GDBzF+/Hg4OTlhyJAhNjNOYxs2bMDNmzcxbNgwALbxswoAU6ZMQWZmJlq2bAm1Wg2tVou5c+di8ODBACrv9wcDpipMpVKZfC8iRY7ZGlsZc0REBI4ePYq9e/eabWttY/T19UVycjJu3ryJdevWYejQodi1a1eJ7a1tfBcuXMCECROwdetWODs7l/l11jTOfv36Kf/drl07BAUFoVmzZli9ejW6du0KwPo/izqdDgEBAZg3bx4AoFOnTjhx4gTi4uIwZMgQpZ21j9PYp59+in79+sHDw6PUdtY2xm+//RZfffUVvv76a7Rp0wbJycmYOHEiPDw8MHToUKXdg34vuSRXBbm7u0OtVhf5P4D09PQiEbStMDyhYwtjHjduHL7//nvs2LEDTZo0UY43bNgQOTk5uHHjhkl7axujo6MjmjdvjoCAAMyfPx8dOnTABx98YDPjO3z4MNLT0+Hv7w97e3vY29tj165dWLJkCezt7dGgQQObGKcxFxcXtGvXDqmpqTbzWWzUqBFat25tcqxVq1Y4f/48ANv6NwcAzp07h59//hkjR45UjtnKZ3Ly5MmYOnUqXnrpJbRr1w6vvPIK3nzzTcyfPx9A5b2XDJiqIEdHR/j7+yMhIcHkeEJCArp162ahXj1YPj4+aNiwocmYc3JysGvXLqsZs4ggIiIC69evx/bt2+Hj42Ny3t/fHw4ODiZjvHz5Mo4fP241YyyOiECj0djM+Hr37o1jx44hOTlZ+QoICEB4eLjy37YwTmMajQanTp1Co0aNbOKzCADBwcFFynqcPn0a3t7eAGzj3xxjn3/+OerXr4/+/fsrx2zlM/nPP//Azs40XFGr1UpZgUp7LyssfZwqlKGswKeffionT56UiRMniouLi5w9e9bSXbtn2dnZcuTIETly5IgAkEWLFsmRI0eUUgkLFiyQOnXqyPr16+XYsWMyePBgq3rEd/To0VKnTh3ZuXOnyWO+//zzj9Jm1KhR0qRJE/n5558lKSlJevXqZVWP+EZFRcnu3bslLS1Njh49KtOmTRM7OzvZunWriFj/+Epi/JSciPWPc9KkSbJz5075888/Zf/+/TJgwABxdXVV/n2x9s+iiL4khL29vcydO1dSU1NlzZo1UrNmTfnqq6+UNrYwThERrVYrXl5eMmXKlCLnrP1nVURk6NCh0rhxY6WswPr168Xd3V3+85//KG0q471kwFSFffTRR+Lt7S2Ojo7i5+enPJ5urXbs2CEAinwNHTpURPSPhsbExEjDhg3FyclJHnvsMTl27JhlO10OxY0NgHz++edKmzt37khERIS4ublJjRo1ZMCAAXL+/HnLdbqchg8frvxMPvzww9K7d28lWBKx/vGVpHDAZO3jNNSocXBwEA8PD3n22WflxIkTynlr/ywa/PDDD9K2bVtxcnKSli1byooVK0zO28o4t2zZIgAkJSWlyDlr/1kVEcnKypIJEyaIl5eXODs7yyOPPCLTp08XjUajtKmM91IlIlJx81VEREREtoc5TERERERmMGAiIiIiMoMBExEREZEZDJiIiIiIzGDARERERGQGAyYiIiIiMxgwEREREZnBgImIiIjIDAZMRERERGYwYCIiIiIygwETEVEZTJ48GQMGDLB0N4jIQriXHBFRGdy8eRNqtRqurq6W7goRWQADJiIiIiIzuCRHRGRGRkYGVCoVTpw4YemuEJGFMGAiIjLjt99+g5OTE3x9fS3dFSKyEAZMRERm/Pbbb2jTpg3s7e0t3RUishAGTEREZiQnJ6Njx46W7gYRWRADJiIiM3777Td06NDB0t0gIgtiwEREVIrc3FycOnWKARNRNceAiYioFCdPnkRubi4DJqJqjgETEVEpkpOT4e3tjYceesjSXSEiC2LARERUil9//RWBgYGW7gYRWRgDJiKiYty9exdJSUlYt24d+vbta+nuEJGFMWAiIirG4sWLMXDgQAwaNAhDhgyxdHeIyMK4lxwRERGRGZxhIiIiIjKDARMRERGRGQyYiIiIiMxgwERERERkBgMmIiIiIjMYMBERERGZwYCJiIiIyAwGTERERERmMGAiIiIiMoMBExEREZEZ/w/mz7liqty40AAAAABJRU5ErkJggg==", "text/plain": [ "Figure(PyObject
)" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/plain": [ "PyObject Text(0.5, 28.0, '$j$')" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "QCbig, RCbig = clgs(big.(A))\n", "\n", "semilogy(Float64.(diag(RCbig)), \"o\", mfc=\"none\") # convert back to Float64 for plotting, since PyPlot doesn't understand BigFloat\n", "semilogy(diag(RM), \"rx\")\n", "legend([\"classical GS in arbitrary precision\",\"modified GS\"], loc=\"upper right\")\n", "ylabel(L\"r_{jj}\")\n", "xlabel(L\"j\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Orthogonality of $Q$" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can also look at how close $Q$ is to being orthogonal by computing $\\Vert Q^*Q - I\\Vert_1$:" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "(52.22409280066077, 4.203714197941792, 3.363088553594440047513746051409831884824413129653818425769182742811618036868859e-37)" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" } ], "source": [ "opnorm(QC'*QC - I, 1), opnorm(QM'*QM - I, 1), opnorm(QCbig'*QCbig - I, 1)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Unfortunately, both classical and modified Gram–Schmidt can yield very non-orthogonal $Q$ matrices! In general, they are *both* unstable as a way to compute $Q$, but the modified Gram–Schmidt algorithm is still useful for least-squares algorithms because $R$ is computed in a backwards-stable way, and the algorithm can be tweaked to solve the least-squares problem $Ax \\approx b$ as explained in Trefethen chapter 19. (The key is to not compute $Q^* b$ via $Q$, but rather by adding $b$ as an additional column of $A$...)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# Householder QR" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The `qr` function in Julia (and Matlab and SciPy), in contrast, defaults to using a Householder QR algorithm. This gives an orthogonal $Q$ matrix and an accurate $R$:" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAkwAAAG0CAYAAADATXgqAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAAPYQAAD2EBqD+naQAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4xLjEsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy8QZhcZAAAgAElEQVR4nOzdfVyUZb4/8M8wMoDAgCaMioIm1WporvjIjCVmmp5SsVpFdDOsLdmGJkBNXCYFxcQhWWn1uLmunWKJ7UG39bBHbTczMMssVk3rhIefKKGQASOPo3D//hgYGWZ4GBiYGebzfr142dz3Pfdc0/qKz173dX2/IkEQBBARERFRu1xsPQAiIiIie8fARERERNQJBiYiIiKiTjAwEREREXWCgYmIiIioEwxMRERERJ1gYCIiIiLqxABbD6A/aGpqwo8//ghvb2+IRCJbD4eIiIi6QBAE3Lx5E8OHD4eLS8dzSAxMVvDjjz9i5MiRth4GERERdcOVK1cwYsSIDq9hYLICb29vAPp/4VKp1MajISIioq7QarUYOXKk4fd4RxiYrKDlMZxUKmVgIiIicjBdWU7DRd9EREREnWBgIiIiIuoEAxMRERFRJ7iGiYiIel1jYyNu3bpl62GQk3F1dYVYLLbKvRiYmkVEROD48eN4+OGH8f7779t6OERE/YIgCLh27RoqKyttPRRyUr6+vhg6dGiP6yQyMDWLjY1FdHQ03nrrLVsPhYio32gJS/7+/hg4cCCL+1KfEQQBtbW1KCsrAwAMGzasR/djYGoWHh6O48eP23oYRET9RmNjoyEs3XXXXbYeDjkhDw8PAEBZWRn8/f179HjOIRZ9nzhxAo8//jiGDx8OkUiEQ4cOmVyze/dujB49Gu7u7ggNDcVnn31mg5ESEVGLljVLAwcOtPFIyJm1/P3r6Ro6h5hhqqmpwQMPPIBnnnkGTzzxhMn5nJwcqFQq7N69G3K5HHv37sX8+fNx4cIFBAYGAgBCQ0PR0NBg8t6jR49i+PDhFo2noaHB6F5ardbCb0RE5Dz4GI5syVp//xwiMM2fPx/z589v9/zrr7+O1atX49lnnwUAZGRk4MiRI9izZw+2bdsGADhz5ozVxrNt2zZs3rzZavfripLKOlTU6No9P8hTggBfjz4cERERkfNwiMDUEZ1OhzNnzuCVV14xOj537lycPHmyVz5zw4YNiIuLM7xu6UXTW0oq63Bw4XNoEIBMeaTJeWV+NtxEQMRHbzI0ERER9QKHWMPUkZ9++gmNjY2QyWRGx2UyGa5du9bl+8ybNw9PPfUUcnNzMWLECJw+fbrda93c3Ax94/qif1xFjQ4NAhCfl4UvGk/isFJh+Pmi8STi87LQIKDDGSgiIkdTUlmH8yVV7f6UVNbZeog2t2nTJkycOLFH9zh+/DhEIlGHpR8OHDgAX1/fHn2Oo3P4GaYWbZ9RCoJg0XPLI0eOWHtIVpUpj8SK6UGQaVIhk7oDSUlASgqgScX1hERkisMwz9aDJCKykpLKOsxJ/xR1txrbvcbDVYyP4x+y+sz6qlWrUFlZabLB6Pjx4wgPD0dFRYXThwdr+/bbb7F582Z88skn0Gq1CAwMxLJly7BhwwajTQOjRo3C5cuXAQDu7u4ICgrC6tWrkZCQ0Otr5Rw+MA0ZMgRisdhkNqmsrMxk1snRlavW6cOSWg1s2QLodEByMsqjY4HMPFsPj4jIaipqdKi71YiMpRMR7O9lcr6wrBqqnAJU1Oi4FMFB3Lp1C66uribHT506hTlz5mDOnDn47//+b8hkMnz55ZeIj4/Hv/71L3zyySeQSCSG65OTk/Hcc8+hvr4eH3/8MdasWQOpVIrnn3++V8fv8I/kJBIJQkNDcezYMaPjx44dQ1hYmI1G1YuSkgCJRB+WJBL9ayKifirY3wshAT4mP+ZClC188MEHuP/+++Hm5oZRo0YhPT3d6Ly5Uji+vr44cOAAAP063BdffBHDhg2Du7s7Ro0aZdisBABVVVX4zW9+A39/f0ilUsyePRv//ve/Tcbx9ttvY9SoUfDx8cGyZctw8+ZNw7mGhgbExsbC398f7u7uUCgUHS47AfSP4AIDAzFw4EBERETgxo0bJtf8/e9/R2hoKNzd3XH33Xdj8+bNuH37ttF3/8///E8sWrQInp6e2LJli8k9BEHA6tWrMXbsWHz44YeYOnUqgoKC8NRTT+Hvf/87Pv/8c+zcudPoPd7e3hg6dChGjRqFZ599FhMmTMDRo0c7/D7W4BCBqbq6GgUFBSgoKAAAFBUVoaCgAMXFxQCAuLg47Nu3D/v378fFixfx8ssvo7i4GC+88IIth907UlLuhCWdTv/aAlwTQERkHWfOnMGvfvUrLFu2DOfOncOmTZuQlJRkCENdsWvXLnz00Uf461//iu+//x7vvPMORo0aBUAfJv7jP/4D165dQ25uLs6cOYNJkybh4Ycfxs8//2y4x6VLl3Do0CEcPnwYhw8fxqefforXXnvNcH7dunX44IMP8NZbb+Hrr79GcHAw5s2bZ3SP1r744gtER0cjJiYGBQUFCA8PNwk7R44cwYoVKxAbG4sLFy5g7969OHDgALZu3Wp03auvvopFixbh3LlziI6ONvmsgoICXLhwAXFxcXBxMY4kDzzwAObMmYPs7Gyz4xQEAcePH8fFixfNzlxZneAAPvnkEwGAyc/TTz9tuOYPf/iDEBQUJEgkEmHSpEnCp59+2mfjq6qqEgAIVVVVvXL/c1crhaD1h4VrCYmCAAhCcrL+RHKyIADCtYREIWj9YeHc1coO73O1olbInBklaBRRQtD6wyY/GkWUkDkzSrhaUdsr34OInEtdXZ1w4cIFoa6uzuL3tvx3r73/rnV2vieefvppQSwWC56enkY/7u7uAgChoqJCEARBWL58ufDII48YvXft2rXCuHHjDK8BCAcPHjS6xsfHR/jzn/8sCIIgKJVKYfbs2UJTU5PJOP75z38KUqlUqK+vNzo+ZswYYe/evYIgCMKrr74qDBw4UNBqtUZjmDZtmiAIglBdXS24uroKWVlZhvM6nU4YPny4kJaWJgjCnd+xLd8rMjJSePTRR40+c+nSpYKPj4/h9cyZM4XU1FSja95++21h2LBhRt9dpVKZfK/W3n33XQGA8M0335g9HxsbK3h4eBhet/ye9/T0FFxdXQUAgru7u5Cfn9/uZ3T099CS398OsYZp1qxZ0P+7b19MTAxiYmL6aES9q23NpcKyaijzsyHLy8LFNQkofep5+JdUAdGx8NPWQ6ZJhVIRBSgVHd639W67FdODUK5aZzjnl5EGWV4W0hVRXBNARE4vPDwce/bsMTr2xRdfYMWKFYbXFy9exKJFi4yukcvlyMjIQGNjY5facKxatQqPPPII7rvvPjz66KN47LHHMHfuXAD6Gazq6mqTtjJ1dXW4dOmS4fWoUaPg7e1teD1s2DBD/7RLly7h1q1bkMvlhvOurq6YOnUqLl68aHZMFy9eREREhNGxGTNm4H/+538Mr8+cOYPTp08bzSg1Njaivr4etbW1hoXakydP7vTfQUcEQTBavwQAa9euxapVq1BeXo6NGzdi9uzZfbIExyECkzNpr+aSSmhCuiIKOFcK8bMqZCii9CfEYVAqouAm0hev7Ax32xERdc7T0xPBwcFGx65evWr0WjCzG7vt/7kXiUQmx1q36Jg0aRKKiorwj3/8Ax9//DF+9atfYc6cOXj//ffR1NSEYcOGme1z2nqXXtvHUSKRCE1NTUbjsWQneWcTFADQ1NSEzZs3Y8mSJSbn3N3dDf/s6enZ4X3uueceAMCFCxfMlkf47rvvcO+99xodGzJkCIKDgxEcHIwPPvgAwcHBmD59OubMmdPpuHuCgcnOtDcLVLZqMobt0mDsHg3SFVHGO0eUCosqfXO3HRE5isKyaouO96Vx48YhL8/4v5knT57Evffea5hd8vPzQ2lpqeH8Dz/8gNraWqP3SKVSLF26FEuXLsWTTz6JRx99FD///DMmTZqEa9euYcCAAYZ1TZYKDg6GRCJBXl4eli9fDkAf2L766iuoVKp2v9epU6eMjrV9PWnSJHz//fcmodJSv/zlL/GLX/wCO3fuxLJly4zWMf373//Gxx9/jDfeeKPd9w8aNAhKpRIJCQn45ptverW0AAOTHTI/C7QL2KMxzAIdbt450m1JSXfCUstuu5Iq630JIqIeGOQpgYerGKqcgnav8XAVd2lmvbfEx8djypQpSElJwdKlS/H555/jjTfewO7duw3XzJ49G2+88QamT5+OpqYmrF+/3mhGaOfOnRg2bBgmTpwIFxcXvPfeexg6dCh8fX0xZ84czJgxA4sXL8b27dtx33334ccff0Rubi4WL17cpcddnp6eWLNmDdauXYvBgwcjMDAQaWlpqK2txerVq82+JzY2FmFhYUhLS8PixYtx9OhRo8dxAKBWq/HYY49h5MiReOqpp+Di4oKzZ8/i3LlzZnfDtUckEmHfvn2YO3cunnjiCWzYsAFDhw7FF198gfj4eMybN6/TcgG//e1vsX37dnzwwQd48sknu/zZlmJgslO9PgtkbrdddGzP70tEZAUBvh74OP4hu+6hOWnSJPz1r3+FWq1GSkoKhg0bhuTkZKxatcpwTXp6Op555hk8+OCDGD58OH7/+98b9Tb18vLC9u3b8cMPP0AsFmPKlCnIzc01zLTk5uZi48aNiI6ORnl5OYYOHYoHH3zQojqDr732GpqamrBy5UrcvHkTkydPxpEjRzBo0CCz10+fPh379u3Dq6++ik2bNmHOnDn43e9+h5RWu7LnzZuHw4cPIzk5GWlpaXB1dcUvfvELQ09XS8jlcpw6dQqbN2/G/PnzDbv3XnzxRezcubPTtWB+fn5YuXIlNm3ahCVLlpjstrMWkdCVh5XUIa1WCx8fH1RVVfW4Tcr5kio8lpmHw0qFfgbJze1OsGloMD3fDnPNelsKvf1Dexxj92iA5OQ7a5jUalxPSMQ0cVin9yYi6or6+noUFRVh9OjRRutaiDrS1NSE1atX48iRI/j0008N65y6q6O/h5b8/uYMkz3r5ixQR816lfnZGJuXhVNB41Hbzd12REREvcXFxQV/+tOfkJmZic8++6zHgclaGJjslF9GGqBJNZkF8tPWA+KOt0+2u3D8Zj0GZm8AAOSPnIDMA60qvVq4246IiKi3uLi44KWXXrL1MIwwMNmhlppL1xMS9WuWujEL1O7C8eJzuLgmAZnSWaY9mizcbUdEROQsGJjszCBPCdxEQLoiCpniMOMF3mZmgdpbqwQAny+PwWyJGNI2C8cbmxeOB/d0px0REZGTYGCyMwG+Hoj46E1U1OjMF5BsNQvU0VolALgUux5XRcBvJRKIWD6AiIio2xiY7FCAr0eXHou1t1apsKwal2LXIz4vC/mBE+6EJZYPICIi6hYGJgdnbq2SX0YaFudloXyKHPLT+biekAjZjq0WLRwnIiKiOxiY+oG2RS5lOh3yAydAfjof6YoojFkeg2CWDyAiIuq23imHSX0vKcnw2E2QSPBN0Hj9wnF5JFQ5BXgsMw+PZeZhmjgM6SwfQERkF2bNmmXU023UqFHIyMgwvL527RoeeeQReHp6GhruikQiHDp0qEefu2rVKixevLhH93A2nGHqL1oVuRTpdPj1zLvxryefB3IKWD6AiMhBnD59Gp6enobXO3fuRGlpKQoKCuDjo9/VXFpa2m5bE2sSBAH79u3D/v378e2336KpqQlBQUGYM2cOlEqlofFuTU0NkpOT8d577+HHH3+Et7c37r//fiQkJOCxxx7r9XH2FQamfsBckUupWo0ZukZAHMbyAUTkmDZtAsRi/X/X2kpJARob9df0I35+fkavL126hNDQUKNq10OHDu31cQiCgOXLl+PQoUNITEzEzp074e/vj6KiIhw7dgxbtmzBgQMHAAAvvPACvvzyS7zxxhsYN24cbty4gZMnT+LGjRu9Ps4+JVCPVVVVCQCEqqqqPv3cc1crBY0iShAA4VpConDuaqXh51pCoiAAgkYRJZy7WmnVz71aUWv0WW1/rlbUWvXziMgx1dXVCRcuXBDq6uq6d4PkZEEA9H925biVPPTQQ8KLL74ovPTSS4Kvr6/g7+8v7N27V6iurhZWrVoleHl5CXfffbeQm5tr9L7jx48LU6ZMESQSiTB06FBh/fr1wq1btwznq6urhZUrVwqenp7C0KFDBY1GIzz00EPCSy+9ZLgmKChI2Llzp+GfARh+nn76aUEQBAGAcPDgQcN7rl69KvzqV78SfH19hcGDBwsLFy4UioqKDOdv374tvPzyy4KPj48wePBgYe3atcKvf/1rYdGiRe3+O8jOzhYACH/729/Mnm9qajL8s4+Pj3DgwIHO/8XaSEd/Dy35/c0ZJgdmaZFLa+is9pMyPxtuIiDiozf5yI+IeqZlZkmtvvO6ebevYUa9l7z11ltYt24dvvzyS+Tk5GDNmjU4dOgQIiIiDDMuK1euRHFxMQYOHIiSkhIsWLAAq1atwn/913/hu+++w3PPPQd3d3dsap4FW7t2LT755BMcPHgQQ4cORWJiIs6cOYOJEyeaHcPp06fx61//GlKpFL///e/h4WH639Ta2lqEh4dj5syZOHHiBAYMGIAtW7bg0UcfxdmzZyGRSJCeno79+/fjT3/6E8aNG4f09HQcPHgQs2fPbvf7Z2dn47777sPChQvNnheJRIZ/Hjp0KHJzc7FkyRJ4e3tb8G/ZwfRGmnM2tpphEoS+n+2x1awWETmeHs8wtWiZUZJIenVmqcVDDz0kKBQKw+vbt28Lnp6ewsqVKw3HSktLBQDC559/LgiCICQmJgr33Xef0czLH/7wB8HLy0tobGwUbt68KUgkEuHdd981nL9x44bg4eHR7gyTIAjCokWLDDNLLdBqhulPf/qTyec2NDQIHh4ewpEjRwRBEIRhw4YJr732muH8rVu3hBEjRnQ4w/SLX/xCWLhwodGxl156SfD09BQ8PT2FgIAAw/FPP/1UGDFihODq6ipMnjxZUKlUQl5eXrv37mucYSIAXS9yaU3maj9pN6oh1aQa+tSNaW7P0oKLzImo25KS7rR3aulY0MsmTJhg+GexWIy77roL48ePNxyTyWQAgLKyMgDAxYsXMWPGDKOZF7lcjurqaly9ehUVFRXQ6XSYMWOG4fzgwYNx33339WicZ86cQWFhocnMTn19PS5duoSqqiqUlpYafe6AAQMwefJkCILQ4b1bfxcA2LhxI1588UV8+OGHSE1NNRx/8MEH8X//9384deoU8vPz8a9//Qu///3vsXnzZiT1wf9WfYWBibqlde0nYcsWSJtrP0VJZwEAVDkFhmv5mI6IeqTVLmBDx4Je/kXs6upq9FokEhkdawkTTU1NAPSLpNsGjJZAIhKJOg0n3dXU1ITQ0FBkZWWZnGu7gNwS99xzD7777juT+/n5+cHf39/keldXV8ycORMzZ87EK6+8gi1btiA5ORnr16+HRNI/StiwDhOhpLIO50uq2v0pqawz/8bm2k8inQ63RS6QF5/FP7THAQAZSyfisFKBLxpPIj4vCw0CTJoEExF1qvWapYYG/Z9qtf64HRk3bhxOnjxpFIxOnjwJb29vBAQEIDg4GK6urjh16pThfEVFBf73f/+3R587adIk/PDDD/D390dwcLDRj4+PD3x8fDBs2DCjz719+zbOnDnT4X0jIyPx/fff429/+1u3xjVu3Djcvn0b9fX13Xq/PeIMk5Pr0SLu5v/X1ySRYIBOh2r5gxi7RwOlohTBSgVC9u8CNKm4npCITHGY+WbCRETtMbfA29xCcDsQExODjIwMKJVKvPjii/j+++/x6quvIi4uDi4uLvDy8sLq1auxdu1a3HXXXZDJZNi4cSNcXHo2bxEVFYUdO3Zg0aJFSE5OxogRI1BcXIwPP/wQa9euxYgRI/DSSy/htddewz333IOxY8fi9ddfR2VlZYf3XbZsGT788EMsW7YMGzZswLx58yCTyXD58mXk5ORALBYbrp01axYiIyMxefJk3HXXXbhw4QISExMRHh4OqVTao+9nTxiYnFx7DXwBfX0nWV4W0hVRqKjRGQWm1rWfLkTH4siy3yI+LwvV8gcRn5eFprvf00+dJyejPDrWeAcfEVFXNDaa3w3X8rqxse/H1I6AgADk5uZi7dq1eOCBBzB48GCsXr0av/vd7wzX7NixA9XV1Vi4cCG8vb0RHx+PqqqqHn3uwIEDceLECaxfvx5LlizBzZs3ERAQgIcfftgQVuLj41FaWopVq1bBxcUF0dHRiIiI6PCzRSIRcnJy8Oabb+LPf/4z0tLScOvWLYwYMQIPP/wwXn/9dcO18+bNw1tvvYXExETU1tZi+PDheOyxx6BuCbX9hEjorQerTkSr1cLHxwdVVVUOl6bPl1Thscw8fNF4ErI2xS+hVuN6QiKmicNwWKlASIAPzpdUGcLR9YRElKvWobCsGqqcAvxDexxj92hwW+SCAUKTfr1BQ4PhM1ruQUTOob6+HkVFRRg9ejTc3d1tPRxyUh39PbTk9zfXMBEA/SJuw9oANzfDNHjrGSfAuPbTNHEYHsvMMyzwni+dhfzACRggNEFovTiTiIjIwfGRHN3RZutuiTIBhUU/AwAKW5UJGPfHdGjrbmO/hyv8vd2MZ5iKz0KbmATp1mTDLJWfth4Qh9nqWxEREfUYAxPd0Wbr7tWJ03Fp5HhAHmlUJgC4sxj8vo/eRLC/F5T52Rjb8pguJh4oqQKiY+GnrYdMkwqlIgpQKmz0xYiIiHqGgYkAmDbwvb52I6ZpUjHt8lkAwJhd2xHs72W4tvVicFu0aCEiIupLDEwEZX42ZC2zQ9GxQEkVCpfH4NKpy4jPy9Iv8P5LEGQ7tupnodqUCgjw9UDER2+iokZnvnSAUsFK30ROjHuLyJas9fePgcnJdTg71FyXSX7lHKZrUoFdmnZLBdiiRQsR2beWyti1tbVmG8cS9YXa2loAptXbLcXA5OTamx1qWcg9Ztd2jBw9GJD5GvdxKulZ7RAi6v/EYjF8fX0N/dYGDhxo0j6EqLcIgoDa2lqUlZXB19fXqNhmdzAwUYezQ8H+XgjI1Jj2cYqO7eNREpEjGjp0KIA7TWqJ+pqvr6/h72FPMDBRh9ouBmepACKyhEgkwrBhw+Dv749bt27ZejjkZFxdXXs8s9SCgYnaZW4xuKWlAkoq6zpsusvF4ETOQSwWW+0XF5EtMDCRWdYoFdCdxr4MWEREZI8YmMgsa5QKsLSxb3cCFhERUV9gYKJ2WaNUQKY8EiumB0GmSYVM6n5nHVSbWk6A5QGLiIiorzAwUa8rV63ThyW1+k6vOjO1nADLAhYREVFfcbH1AMhJJCXdKUvQUsupHeWqdfpdeWo14Oam/zM52WjGiYiIqC8xMFHfaNPYFykpHV9vQcAiIiLqbQxM1Ov8MtIMs0RoaDDMHvllpLX/pjYB6/rajSgsqwagr0J+vqTK8FNSWddH34SIiJwV1zBRr+pOLafWxTJLlAm4OnE6pmlScenUZUAeCVVOgdH9uXOOiIh6GwMT9Zqu1HKSiASU3WzA+ZIqFJZVGwLW/ynXozY6FoVFP+PSyPGYdvks4vOyAABjdm1HsL8Xd84REVGfYWCiXtNeLaeym/X49jdx0EGE9LDlwIHThnMqoQn5gRNw+utiZLQErOaaTPF5WZBfOYeRowfr+9tx5xwREfURBibqVeZqOZ0vAb4RRGbrLQ2sPYG787NxKnA8MpZOBACocgowZtd2XP9LEKZrUgGZb4elCYiIiKyNi75bqa2tRVBQEBISEmw9lH4vUx6J6wmJkGlSEbJ/F0ICfBCyfxfuzkzTzxrJIxHs74Vgfy8AQLC/F2Q7tnLnHBER2QQDUytbt27FtGnTbD0Mp2FxvSVLSxMQERFZCQNTsx9++AHfffcdFixYYOuhOJcu1lvqVmkCIiIiK3GIwHTixAk8/vjjGD58OEQiEQ4dOmRyze7duzF69Gi4u7sjNDQUn332mUWfkZCQgG3btllryNRVXZg1UuZnQ9a8wPt8dKy+/lJ0rOGRnjI/2wYDJyIiZ+IQi75ramrwwAMP4JlnnsETTzxhcj4nJwcqlQq7d++GXC7H3r17MX/+fFy4cAGBgYEAgNDQUDQ0NJi89+jRozh9+jTuvfde3HvvvTh58mSn42loaDC6l1ar7cG3c16t6y0Zesap1fDT1gPiMABdK03gJtJfR0RE1FtEgiAIth6EJUQiEQ4ePIjFixcbjk2bNg2TJk3Cnj17DMfGjh2LxYsXd2nWaMOGDXjnnXcgFotRXV2NW7duIT4+Hmq12uz1mzZtwubNm02OV1VVQSqVduNbOZfzJVU4suy3iG8paNlqzZJfRhpkmlSkK6Iw790/ICTAByWVdaio0bV7v0GeEtZgIiIii2m1Wvj4+HTp97fDByadToeBAwfivffeQ0REhOG6l156CQUFBfj0008tuv+BAwdw/vx5aDSadq8xN8M0cuRIBqYuKqmsw8GFz6FB0O+Wa6ur1bsZpIiIqCcsCUwO8UiuIz/99BMaGxshk8mMjstkMly7dq1XPtPNzQ1ubm69cm9n0F5BSwOlotOwY63QRURE1BUOH5haiEQio9eCIJgc64pVq1ZZaUTUEXMFLS1RUaNDgwCzxS/ZMoWIiKzN4QPTkCFDIBaLTWaTysrKTGadqH/JlEdixfQgyDSpkEnd7ywcZ8sUIiKyMocoK9ARiUSC0NBQHDt2zOj4sWPHEBYWZqNRUV+xuPglERFRNzhEYKqurkZBQQEKCgoAAEVFRSgoKEBxcTEAIC4uDvv27cP+/ftx8eJFvPzyyyguLsYLL7xgy2FTX+li8UsiIqLucohHcl999RXCw8MNr+Pi4gAATz/9NA4cOIClS5fixo0bSE5ORmlpKUJCQpCbm4ugoCBbDZn6krnil9Gxth4VERH1Iw4RmGbNmoXOqh/ExMQgJiamj0ZE9qIrxS+JiIh6yiECE5E5yvxsyFqKX0bHAiVVQHQs/LT1+pYpiihAqbD1MImIqB9gYCKHxJYpRETUlxyu0rc9sqRSKFkPK30TEVFPOFWlb+fQvWoAACAASURBVHJePS1+SURE1FUOUVaAiIiIyJYYmIiIiIg6wcBERERE1AkGJiIiIqJOMDARERERdYKBiYiIiKgTLCtAToE1m4iIqCcYmKjfK6msw8GFz6FBADLlkSbnlfnZcBMBER+9ydBERERmMTBRv1dRo0ODAMTnZWHF9CCUq9YZzvllpEGWl4V0RRQqanQMTEREZBYDE/U7bR+/FZZVI1MeiQXjh2GsJhUeEjGkW5OBlBRAk4rrCYnIFIdhng3HTERE9o2BifqVjh6/zZfOQlbgUchTUyBotkOk0wHJySiPjjVu3ktERNQGd8lRv9L68dsXjSdxWKlAxtKJAIB/aI9DXnwWt0Uu+rAkkQBJSTYeMREROQIGJup3MuWRuJ6QCJkmFSH7dyHY3wvK/GyM3aNBtfxBDBCa0CSRADqd/rEcERFRJ/hIjvqlctU6yKTugFqNcVu2IESnQ7X8QXjln0C6Igrz3v0DQvbvAtRq+GnrAXGYrYdMRER2jDNM1H8lJQESCVx0OtwWucAr/wQurklApjwShWXVOB8da5iJUuZn23q0RERkxzjDRP1XSgqg00GQSDBAp0N+4ARESWcBAFQ5BfprxGFQKqLgJtIXryQiIjKHgYn6Jb+MNECTCiQnQ5SUBO1GNeSpKfiH9jjmS2chY+lEBPt76S9WKljpm4iIOsTARP2OMj8bsrwsXE9I1JcMKKkCYuLhp2vEWE0qlIpSBCsVCAnwMdRsavlpi0GKiIgABibqZwZ5SuAmAtIVUcgUhxnXV2rz+I0tU4iIqKsYmKhfCfD1QMRHb6KiRme+cnerx2/nS6rYMoWIiLqEgYn6nQBfjy4HnEx5JFZMD4JMk6ovQ5CUxJYpRERkgoGJnF7rmk3YskVf0JItU4iIqBXWYSICDDWbwJYpRERkBgMTEWCo2QS2TCEiIjP4SI6cXuuaTYY1TBa0TGkpTdAeliYgInJ8DEzk1MzWbIqOhZ+2Xt8yRREFKBXtvp+lCYiInAMDEzktS2o2taeiRsfSBEREToCBiZyWJTWbOsLSBERE/R8DEzm19mo2WdoyhaUJiIj6NwYmojYsWZdkJCnpTlhqKU1QUtVHoyYiot7EwETUhiXrkoy0KU2g3ahG4ZPPAwAKy6pNPoe754iIHAcDE5EZlq5LaluaQLtRDWlqCvze+RsQmQpVToHR/bl7jojIsTAwEbWjq+uSzJUmKHzyefi98zfIi88iKzsR5YdyEezvBYC754iIHBEDE1FH2qxLKlEmoLDoZwD6x2xSjwGQiATzpQkiU5GVnQh58Vlo398L6dZkwyzVxTUJyJTOwhg+qiMicggMTEQdabMu6erE6bg0cjwgj7zzmC1sOQD9TJNEJOD+P+6Etu4WVDkFKD+Uqw9LqSmAZjug0+GLoAnIO1cKyGH2UV3LPfy93QzHGaKIiGyLgYmoHW3XJV1fuxHTNKmYdvksAGDMru1mH7P5e7sZwk6wv5d+Zqk5LDVJJMgbOR7xeVmGe0g9BuDb38Rh0uXzkBefRboiCtEHThvGwfVORES2x+a7RGYo87Mha17gfT46FudLqvD58hikK6IA6HfQzfjLboQE+CBk/y7DtebKELSepXLR6XfWXVyTYLiHv7e7ISyVT5EjUx6JjKUTcVipwBeNJxGfl4UGAR32qyMiot7FGSaiNjpsmdIciORXzmG6JhXYpemwSKW5Wap4TSoujk9AuiIK8ZpU+O3SwEWnQ7X8Qfjln4DSLRvBSgVC9u9itXAiIjvBwETURnstUwrLqqHKKcCYXdsxcvRgQObbYZFKs7vnlsfg0qnLiN+jQa4iCk3Ns04N4gH4fzkfwS8jDfGaVDTd/R6rhRMR2REGJiIz2muZAujXJQVkaowWgyMlBYiONVzT1Vkql+Z1TW46Hfwy0lCuWgffnWlwY7VwIiK7wsBEZKG2j9mQkgKo1fDT1gPiMACdz1ItGD8MY/OygORkXIiOxZFlv0W8JhWen+fBrfG2YeapbRAjIiLbYGAisoC5x2yIjoWfth4yTSqUiihAqQDQ/iyVMj8bY1vdo7CsGpnySPy64f/BL/8E8gMnwOfkCf0apjZBjIiIbIOBiaiLOnzMJg6DUhEFN5H+OkvvoczPht/pfOQHTsDXQSG4/2YDzrcTxIiIqO+JBEEQbD0IW/v++++xdOlSo9fZ2dlYvHhxl96v1Wrh4+ODqqoqSKXS3hom2YGSyroOt/d3pcBk23uU3azHt7+Jg04QmS1LwDpMRES9w5Lf3wxMbVRXV2PUqFG4fPkyPD09u/QeBibqKWsEMSIisowlv7/5SK6Njz76CA8//HCXwxKRNXS0K4+IiGzPISp9nzhxAo8//jiGDx8OkUiEQ4cOmVyze/dujB49Gu7u7ggNDcVnn33Wrc/661//avR4joiIiMghZphqamrwwAMP4JlnnsETTzxhcj4nJwcqlQq7d++GXC7H3r17MX/+fFy4cAGBgYEAgNDQUDQ0NJi89+jRoxg+fDgA/dRcfn4+3n333Q7H09DQYHQvrVbbk69HREREds7h1jCJRCIcPHjQaEH2tGnTMGnSJOzZs8dwbOzYsVi8eDG2bdvW5Xu//fbbOHLkCN55550Or9u0aRM2b95scpxrmIiIiByHJWuYHOKRXEd0Oh3OnDmDuXPnGh2fO3cuTp48adG9uvo4bsOGDaiqqjL8XLlyxaLPISIiIsfiEI/kOvLTTz+hsbERMpnM6LhMJsO1a9e6fJ+qqip8+eWX+OCDDzq91s3NDW5ubhaPlYiIiByTwwemFiKRyOi1IAgmxzri4+OD69evW3tYRERE1A84fGAaMmQIxGKxyWxSWVmZyawTUX/Buk1ERH3L4QOTRCJBaGgojh07hoiICMPxY8eOYdGiRTYcGVHvKKmsw8GFz6FBACuDExH1EYcITNXV1SgsLDS8LioqQkFBAQYPHozAwEDExcVh5cqVmDx5MmbMmIE//vGPKC4uxgsvvGDDURP1jooaHRoEID4vCyumB6Fctc5wzi8jDbK8LKQrolBRo2NgIiKyEocITF999RXCw8MNr+Pi4gAATz/9NA4cOIClS5fixo0bSE5ORmlpKUJCQpCbm4ugoCBbDZmoV2XKI7FiehBkmlTIpO5AUhKQkgJoUnE9IRGZ4jDMs/UgiYj6EYcITLNmzUJn5aJiYmIQExPTRyMisr1y1Tp9WFKrgS1bAJ0OSE5GeXQskJln6+EREfUrDl+HicipJSUBEok+LEkk+tdERGR1DExEjiwl5U5Y0un0r4mIyOoYmIgclF9Gmv5xXHIy0NCg/1Ot1h8nIiKrcog1TERkTJmfDVleFq4nJOrXLJVUoeyp32BYqRZjNalQKqJQuHSi0XtYm4mIqPsYmIgczCBPCdxEQLoiCpniMMMCb1VeFr4RuSBXEQWx0ARVToHhPazNRETUMwxMRA4mwNcDER+9iYoanaF0QGFZNS7lZyM+LwsX1ySgNHYtDnvr+x2yNhMRUc8xMBE5oABfD5Pgo2quzTRWk4qxw6SszUREZEUMTET9CGszERH1Du6SI+pvWJuJiMjqGJiI+hvWZiIisjoGJqJ+hLWZiIh6B9cwEdmxkso6VNTo2j3furaSudpMiI6Fn7YesubaTFAq+mroRET9CgMTkZ0qqazDwYXPoUEAMuWRJudb11ZqrzYTAEAcBqUiChKRgLKbDThfUmV0Hxa0JCLqHAMTkZ2qqNGhQQDi87KwYnoQylXrDOfa1lYKCfAxqc0EAGU36/Htb+KggwjpYcuBA6eNPoMFLYmIuoaBiciOZTbXVpJpUvXlAjqorWSuNtP5EuAbQdSl0MXARETUPgYmIjvX09pKloQuIiIyj7vkiBxBD2srlavWGXbMwc3NsJOu9YwTERG1j4GJyBFYo7ZSm9BVokxAYVk1AH0vuvMlVUY/JZV1Vv4SRESOi4/kiOycX0YaoEnVzxC1PE5Tq+GnrQfEYV2/UZvQdXXidFwaOR6QR0KVU2B0KReDExEZY2AismPWqq3UNnRdX7sR0zSpmHb5LABgzK7tCPb3MlzLxeBERMYYmIjsVFdqK7mJ9Nd1xFzoKlweg0unLiM+LwvxeVm4/pcgyHZs5WJwIqJ2iARBEGw9CEen1Wrh4+ODqqoqSKVSWw+H+hFLKn239/7Oil/Kr5zD9Mtn76xvSk7G+ehYPJaZh8NKBUICfKzyXYiI7I0lv785w0Rkx8zVVrL0/eYKWhaWVUOVU4Axu7Zj5OjBgMzXeAdem2rgRETOjoGJqJ/rKHQF+3shIFNjtBhcu1GNwiefBwDDLrrW2EqFiJwRAxORE2u7GFy7UQ1pagr83vkbEJnK3XNERM0YmIiclNnF4E8+D793/gZ58VlkZSei/FAud88REYGBicgpdbgDLzIVWdmJkBefhfb9vZBuTebuOSJyegxMRE6os8Xg5Ydy9WEpNQXQbLe4fx0RUX/DwETkpDpbDC7dmnwnLHH3HBE5OfaSIyLzrNG/joion2BgIiITfhlpgFqt3z3X0KD/U63WHycickJ8JEdERqzVv46IqD9hYCIiA2v1ryMi6m/YS84K2EuO+pOe9q8jInIU7CVHRN1mbvdc6xBVUaMzCVQ9CVEMaETkCBiYiKhDJZV1OLjwOTQIQKY80uR8T9ql9Oa9iYisiYGJiDpUUaNDgwDE52VhxfQglKvWGc71tF1Kb96biMiaGJiIqFOZ8kismB4EmSYVMqm7voilldql9Oa9iYishYGJiLqkXLVOH2jUamDLFqu2S+nNexMRWUO3Cldu374dAHD27FncunXLqgMiIjuWlHSn8ndLuxRHuDcRUQ91a4ZJodAXrdu0aRMuXrwIV1dXhISEYPz48Rg/fjymTJkCmUxm1YESkR0w1y4lOtb+701E1EPdmmGSy+UAgA8//BAXL17EqVOnoFKp4O/vj48//hgLFixAEv/fIVG/0pvtUtiKhYjsXbdmmLZv347169fj7NmzGDt2LAYOHIipU6di6tSphmtCQ0ORwmadRP1Cb7ZLYSsWInIEXQ5MjY2NOHToEObOndulR3KnTp3qtUETUd/pTrsUc8Uoy27WQ1t3GwAg9XCFv7cbym7WQyIS2IqFiOyeRa1RPDw88O233+Luu+82Ol5bW4vz58/j3LlzOHfuHD777DMsWLDAaWaY2BqF+jtLqnG3V4xSlZeFRpF+FYBYaEKGIspwTpmfDYlIwP1/3Al/b7d2701EZE291hpl6tSpKCoqMglMfCRH1L+Za5fSHnPFKAvLqnEpPxvxeVkA9LNVGUsnItjfy6hApb+3G0ICfHrzqxARdYtFi75jY2ORmJiIK1eudHqtvT6Si4iIwKBBg/Dkk0+anDt8+DDuu+8+3HPPPdi3b58NRkfUP2TKI3E9IREyTSpC9u9CsL+XyTXB/l4I2b8LspYClWZaoxAR2QuLZpieeuopAMD999+PhQsXYtasWfjlL3+J8ePHQyIxXmPg6upqvVFaUWxsLKKjo/HWW28ZHb99+zbi4uLwySefQCqVYtKkSViyZAkGDx5so5ESObbWxSjHbdmCEJ0O1xMSAQDxmlQ03f0eC1QSkcOwaIapqKgIBw8eREJCAmpra7Ft2zZMnToVXl5emDBhQm+N0arCw8Ph7e1tcvzLL7/E/fffj4CAAHh7e2PBggU4cuSIDUZI1I80F6N00enQIB6ActU6lKvWoUE8AC4sUElEDsSiwBQUFIRFixZBrVbjww8/xKVLl1BZWYl//vOfeP7553s8mBMnTuDxxx/H8OHDIRKJcOjQIZNrdu/ejdGjR8Pd3R2hoaH47LPPevy5APDjjz8iICDA8HrEiBEoKSmxyr2JnFZzMcomiQRujbfhl5EGv4w0uDXeRlPrApVERHaux73kvL29MXPmTMycObPHg6mpqcEDDzyAZ555Bk888YTJ+ZycHKhUKuzevRtyuRx79+7F/PnzceHCBQQGBgLQLzZvaGgwee/Ro0cxfPjwdj/b3GZBkUhk9tqGhgajz9BqtZ1+NyJn45eRBmhSgeRkXIiOxZFlv0W8JhWAftH3vHf/gJD9u/QFKrX1gDjMxiMmImqfXTXfnT9/PubPn9/u+ddffx2rV6/Gs88+CwDIyMjAkSNHsGfPHmzbtg0AcObMmW59dkBAgNGM0tWrVzFt2jSz127btg2bN2/u1ucQOYO2xSgLy6pNriksq2aBSiJyGHYVmDqi0+lw5swZvPLKK0bH586di5MnT/b4/lOnTsX58+dRUlICqVSK3NxcqNVqs9du2LABcXFxhtdarRYjR47s8RiIHE17BSrNFaNUCU1Ib669JBaaoMop0L+BBSqJyAE4TGD66aef0NjYaNLUVyaT4dq1a12+z7x58/D111+jpqYGI0aMwMGDBzFlyhQMGDAA6enpCA8PR1NTE9atW4e77rrL7D3c3Nzg5uZm9hyRs2ivQCUAIGw5lPnZiD/5F0MxyrJVk40qfR9uXaBSqWCBSiKyaw4TmFq0XVckCEK7a43M6Wjn28KFC7Fw4cJuj43ImZgrUNnCfDFKFqQkIsflMIFpyJAhEIvFJrNJZWVlJrNORNQ3MuWRWDE9CDJNqr7mUlKSftdbSzFKcRjm2XqQRERWYFFZAVuSSCQIDQ3FsWPHjI4fO3YMYWHcXUNkK+WqdUByMqBWA25u+j+Tk41mnIiIHJ1dBabq6moUFBSgoEC/GLSoqAgFBQUoLi4GAMTFxWHfvn3Yv38/Ll68iJdffhnFxcV44YUXbDlsImouUAkWoySifsquHsl99dVXCA8PN7xu2Yn29NNP48CBA1i6dClu3LiB5ORklJaWIiQkBLm5uQgKCrLVkIkIMBSoROtilNGxth4VEZHV2FVgmjVrltkCkq3FxMQgJiamj0ZERJ1pXaDSsIapF4tRmitl0Bp32xFRb7CrwEREjqVtgUqUVPVqMcoOSxk0j8dNBER89CZDExFZFQMTEXXLIE8J3EQwKVAJoNeKUXa1lEFFjY6BiYisioGJiLolwNcDER+9iYoanfnSAb1UjJKlDIjIFhiYiKjbAnw9bDKTU65apw9LajWwZYt+oXlysv6xYOuZLiIiK2FgIiK7ZG5xd0sT38KyagxSJiCgJSy1lDIoqbLFUInICTAwEZHd6Wxx96XY9Rh65RwCWMqAiPoIAxMR2Z32FncXllXjUux6xOdlAQCuJyRCtmNrr5cyICJiYCIiu2RucbdfRhoWN4eldEUUxiyPQXAvlzIgIgIYmIjIjrVd3C3T6fBF0ATkjRyvf1SXU3Dn4l4qZUBEBDAwEZG9S0q6sxNOIsGIglMYU/QzkFOAjKUTEezvdefaXiplQETEwERE9q1Nn7qATA0qmhd3B/t7ISTAx2hHXUWNzmR3HUMUEfUUAxMR2a2u9KljuxQi6gsMTERkl7rap47tUoioLzAwEZHdsaRPXUWNju1SiKjXMTARkd2xpE9dy3oltkshot7EwEREdqlbfera7Khrr12KubYrrXGROBG1xcBERP1Hmx115tqlcJE4EXUHAxMR9Qvt7agTl2oB6Syjxr0ti8QjfjkCtes3GN2Di8SJyBwGJiJyeOZ21JU99RsM3Pcepu/RQKkohar1G+SRmF58DvLM7dB6SyDdmsxF4kTUIQYmInJoHe2oUwZO0AejK2fxy1Wvw9/bzdDAV158FvmBEyBPTQE027lInIg6xMBERA6tvR11hWXVUAFYMH4Ypu/RAO/tNWrge3FNAqKks/B/O5fApZNF4kREDExE5PA62lHXuPF3wDCpUQPfdEUUxjyngjJ2/Z2w1M4icSIiAHCx9QCIiHpdUpIhFDVJJMiUR+K+NzMQ37zuCQ0N+sXiarV+8TgRURucYSKi/q9VuQEXnQ5Z2YkYW3xWP9O0PAbB7bRdISJqwcBERP1a23ID2o1qyFNTkB84QV+HKafgzsVt2q4QEbVgYCKifstsA9+YeAy82QB5ZhqU+dkYs2s7gv29Wr1JwUrfRGSCgYmI+qUOG/gOfBBKRQncRMCU0YMZjoioUyJBEARbD8LRabVa+Pj4oKqqClKp1NbDIaJm7BlHRB2x5Pc3Z5iIqN/qVgNfIiIzWFaAiIiIqBOcYSIiAh/fEVHHGJiIyOmVVNbh4MLn0CBAX2qgDWV+NtxEQMRHbzI0ETkpBiYicnoVNTo0CEB8XhZWTA9CuWqd4ZxfRhpkeVlIV0ShokbHwETkpBiYiIign1laMT0IMk0qZFJ3fTuVlBRAk4rrCYnIFIcZNfclIufCwERE1KxctU4flpob9UKnA5KT9UUvW9dxIiKnw8BERNRaUtKdsCSRoESZgMKinwEAhWXVJpdzMTiRc2BgIiJqrVWjXuh0uDpxOi6NHA/II6Fq3XcOXAxO5EwYmIiImrVt1Ht97UZM06Ri2uWzAGDUd46LwYmcCwMTERHMN+otXB6DS6cuIz4vC/F5Wbj+lyDIdmzlYnAiJ8TAREROr8NGvc11meRXzmG6JhXYpeFicCInxMBERE4vwNcDER+9iYoandFsUWFZNVQ5BRizaztGjh4MyHy5GJzISTEwERGh40a9wf5eCMjUcDE4kRNjYCIi6gQXgxMRAxMRUQe4GJyIAAYmIqJ2cTE4EbVgYCIiakd3FoMjKQkoqbLZmImodzAwERF1wNLF4EhJAaJj+3iURNTbnC4wRURE4Pjx43j44Yfx/vvvG47fvHkTs2fPxq1bt9DY2IjY2Fg899xzNhwpEdm7tovBkZICqNXw09YD4jDDdSWVdaio0bV7H5YgILJ/TheYYmNjER0djbfeesvo+MCBA/Hpp59i4MCBqK2tRUhICJYsWYK77rrLRiMlIntmbjE4omPhp62HTJMKpSIKUCpQUlmHgwufQ4MAZDave2p7H5YgILJ/TheYwsPDcfz4cZPjYrEYAwcOBADU19ejsbERgiD08eiIyBF0uBhcHAalIgpuIv11FTU6NAhAfF4WVkwPQrlqneFSliAgchx2FZhOnDiBHTt24MyZMygtLcXBgwexePFio2t2796NHTt2oLS0FPfffz8yMjIwc+ZMq3x+ZWUlHnroIfzwww/YsWMHhgwZYpX7ElH/0t5icAOlwvCYraJGh0x5JFZMD4JMkwqZ1P3O4zuWICByGHYVmGpqavDAAw/gmWeewRNPPGFyPicnByqVCrt374ZcLsfevXsxf/58XLhwAYGBgQCA0NBQNDQ0mLz36NGjGD58eIef7+vri3//+9+4fv06lixZgieffBIymczkuoaGBqPP0Gq1ln5VInJwHS0GN6dctU4fltRqYMsWliAgcjB2FZjmz5+P+fPnt3v+9ddfx+rVq/Hss88CADIyMnDkyBHs2bMH27ZtAwCcOXOmx+OQyWSYMGECTpw4gaeeesrk/LZt27B58+Yefw4R9S/mFne39JkrLKvGIGUCAlrCEksQEDkUuwpMHdHpdDhz5gxeeeUVo+Nz587FyZMne3z/69evw8PDA1KpFFqtFidOnMCaNWvMXrthwwbExcUZXmu1WowcObLHYyAix9XZ4u5Lsesx9Mo5BLAEAZFDcpjA9NNPP6GxsdHkEZlMJsO1a9e6fJ958+bh66+/Rk1NDUaMGIGDBw9iypQpuHr1KlavXg1BECAIAl588UVMmDDB7D3c3Nzg5ubWo+9DRP1Le4u7C8uqcSl2PeLzsgAA1xMS77RRMVOCgIjsk8MEphYikcjotSAIJsc6cuTIEbPHQ0NDUVBQYPYcEVFXmFvc7ZeRhsXNYSldEYUxy2MQ3E4JAiKyXw4TmIYMGQKxWGwym1RWVmZ2YTYRkS20Xdwt0+nwRdAE5I0cr39Ul9Pq/5g1lyCQiASU3WzA+TbrmVjQksh+OExgkkgkCA0NxbFjxxAREWE4fuzYMSxatMiGIyMiaiMp6c5OOIkEIwpOYUzRz0BOATKWTkSwvxfKbtbj29/EQQcR0sOWAwdOG92CBS2J7ItdBabq6moUFhYaXhcVFaGgoACDBw9GYGAg4uLisHLlSkyePBkzZszAH//4RxQXF+OFF16w4aiJiNpISTHqLxeQqUFF8+LuYH8vhAT44HwJ8I0gYkFLIgdhV4Hpq6++Qnh4uOF1y060p59+GgcOHMDSpUtx48YNJCcno7S0FCEhIcjNzUVQUJCthkxEZKSr/eUA82ue2itoyX50RLZlV4Fp1qxZnbYjiYmJQUxMTB+NiIio67raX661rhS0ZD86Ituzq8BEROSoLOkvZ6LNmqe2BS3Zj47I9hiYiIiswJL+cibarHkyV9CS/eiIbIuBiYjISiztLwdYtuaJ/eiIbMfF1gMgInJWyvxsyJpniM5Hx+J8SRXOR8fqq4FrUqHMzzZ9U1LSnZmolsd3RNTrOMNERGQDlqx5Mtod14XHd0RkfQxMREQ2YMmap5bAZMnjOyKyLgYmIiIbsWTNU3dKFhCR9TAwERHZmbZFKstu1kMiErpXsoCIrIKBiYjIjrRbpDJsOQD9TJNEJOD+P+6Ev7eb/lxHJQuIyCoYmIiI7EhXi1T6e7shJMDHhiMlci4MTEREdoZFKonsDwMTEZEdYpFKIvvCwpVERPaKRSqJ7AYDExGRvTJXpJKIbIKBiYjIDvllpOkfxyUnAw0N+j/Vav1xIupzXMNERGRnWKSSyP4wMBER2RFLeswRUd8RCYIg2HoQjk6r1cLHxwdVVVWQSqW2Hg4RObi2lb7bYpFKIuuw5Pc3Z5iIiOyMJT3miKhvMDAREfVDnKUisi4GJiKifqbdfnTNlPnZcBMBER+9ydBE1EUMTERE/UxX+9FV1OgYmIi6iIGJiMiBmXv0VlhWjUx5JBaMH4ax7EdHZBUMTEREDqqzR2+550pRFTQB09mPjqjHGJiIiBxUe4/eCsuqcSl2PeKbH71NLf0OLq370ZVU2XjkRI6HrVGIiBxYpjwS1xMSIdOkImT/LoQE+GDGX3Yjd41hRgAAIABJREFUPi8LF9ckAMCdsMR+dETdxhkmIiIHV65ap1+n1PzoTabTIV0RhQXQzz5dT0iEbMdWfVhSq+GnrQfEYbYeNpFD4QwTEVF/kJRkmEVqkujbpozdo0G6IgqfL4/B+ZIqnI+ONcxGKfOzbTxgIsfCGSYiov4gJUX/yE0igYtOB8WVc/p+dPJIIKfgznXsR0fULQxMREQOzi8jDdCkAsnJhvIB09RqSBfMRSaAjKUTEezvdecNSoVRpW/tKxtR1wijek2t7+0hBqSvbe2jb0NknxiYiIgcmDI/G7LmdUrl0bH6HXDRsfDT1mOsJhVKRSmClQqEBPiYfX9JZR0+OHUFsZ++jXdOXTYqT6DMz0Z8XhZ2PbQST1TWscglOTUGJiIiBzXIUwI3EfSP3sRhxrWVuvjoraJGh9enL8Uj42SI36MxlCdoqQh+cU0CXpfOwmxWBScnx8BEROSgAnw9EPHRm6io0Zmv2t3m0VtHGjf+DhgmhUythmyXxlDkspFFLokAMDARETm0AF8P6838JCXdqQjOIpdERlhWgIiI9FrttGORSyJjDExERKTfaadW63faNTTo/1Sr9ceJiI/kiIicXctOO0NZAsDwp0ythlIRBSgVNhwhke0xMBEROTmx0ISLaxL0C7xbr1mKjoW4VAvx2RLbDY7ITjAwERE5sUGeEuwN/zUybjWa3w0nnQWPcDGeYlVwcnIMTERETizA1wMfxz+Eihqdybmym/XQ1t2G1MMVFTU6k2u6WrKAqD9gYCIicnLmShOUVNbhqf/8HHW3Gtt9n9sAF+xZEQp/bzej4wxS1B8xMBERkYmKGh3qbjWa9KEru1mPNe98jYbbTWi43YToA6dN3uvhKsbH8Q8xNFG/wsBERERmqfKyMKPxJGQ77jTePV8CrDn+Np6oLcIHA0djVOYOo0Al3roFR86WoOKFGQjw9UBJZZ3Zx30tWs9GWXItUV9jYCIichKWBpJGkQtkmlRA6n6n3ACAKVe+xcjis5gSWAsff687jX1TUoA9GuQqogyfNyf90y491gMEw8xVezhzRbbEwERE5AS6El7aBpJMeSRWTA+CTK3WX5CUBL+MNIQUn0X5FDnkp/NxPSMN2LFVH5bUalxPSESmOAzz0LPHem3XRxWWVUOVU4AKNgEmG2FgIiJyAhU1Ojz/yX9h3oQAfaPdNto+SmtRrloHmdRdXwV8yxbIdDqkK6IwZtd2XIpdj3hNKtCqWW+5mWa9wa1noaB/rNdwuwnxc+9F+tH/RcbSiQCA/6dci3kTAvD9cyqocgrg7+1meJ9fRhpUnxexgCbZjNO1RomIiMCgQYPw5JNPmpwrKipCeHg4xo0bh/Hjx6OmpsYGIyQi6h2NIheM3aNByP5dCAnwufOzfxfG7tGgUdTOr4SkJEN/uSaJBJnySAD6Gaimlr5zLc1621DlZZltr6LKy8KyV56BKi8Lwf5eCPb3MoxvjnI5VHlZhmsbZj4EmSYVjSIXFJZV43xJFc6XVOFf313HxZi1KHpxreFY65+Syjrr/IsjghPOMMXGxiI6OhpvvfWWyblVq1Zhy5YtmDlzJn7++We4ubmZuQMRkWMy94jN3KO01grLquGXkQZZc1hy0emgzM/GlUfuhTI/Gy5tm/VGxxq9v6N1UH7FZzEl8KbR+J4XrsAr/wSmBE4AAGg3qiHNO2G4RpVTYPhnZX42ZudlIV0RhUwzRTe55omsyekCU3h4OI4fP25y/Ntvv4WrqytmzpwJABg8eHAfj4yIqPe1fcTW3qO0QZ4SeLiKcSl2PRa3hBJ5JJT52YjPy0J+9DnIi8/i/5TrUbt+gz5UqdUQl2oB6SzDfbq6DqpctQ7K/Gx45Z9AtfxByPNPoOluf7g0PwKUBw9B/IHfY8H4YWjc+DuIt27B2LwsfLzst8gMmm+yToprnsja7OqR3IkTJ/D44/+/vXsPi7LO+wf+hoFBBPGQgoJCaqWWooWUOBPSSRe11s3fXsriad2eJ61QEtLCE3lAKzCUX7LldtjSxKvLw2OurWJ5gixdEcXDY2om4irkKicVBofP88fM3M7NaTwAw0zv13XNZXPfN/d8v86MvPsen4e/vz9cXFywadOmWtesXLkS3bt3R6tWrRASEoK9e/c2ymufOnUK3t7eeOGFF/DYY48hKSmpUe5LRNTiWHWx1deVFtDOEz9Wf4+4rDUojE/AsIwPsCVGj0dXpWBfUDB0+UeQHRiMp1s/iZFpWXhCMxgp+mj0SU/G6/syUFRWidNF5QCAfX96BaUJc00hzcMDfslJSNFHI/vDdUjRR8MvOQkP9/BVXuuXdZtRqXGDq1UX4KVpbyj379vDF33Sk5Gij0b5zLcA3BonZXlYhyeixtCiWpiuXbuG/v37489//jNGjx5d6/y6desQGxuLlStXQqfT4cMPP0RkZCSOHz+OwMBAAEBISAgqKytr/ez27dvh7+9f72tXVVVh7969yM3Nha+vL373u98hNDQUzz33XK1rKysrVa9RWlp6N9UlIrKPhQtvhaV6utIAwMfdFViwAH5z58JPOdoWpWNHojyrHR58PAxbYvTKrDfL2CaN0aia+Ra7LhcebmE44W7q0jO6m0JQ3JXrSNNF4fX9X8HVYEClxg37/vQKei1eBA/jTRjdtdBYdQFaX1ut1UIj1ei1KlXVogWYZgRqFi9C7JELOG0eUG6N6znR3WhRgSkyMhKRkZH1nl+2bBn+8pe/4KWXXgIApKamYtu2bUhPT8eSJUsAAAcPHryr1+7atStCQ0PRrVs3AMDw4cORm5tbZ2BasmQJ3n777bt6HSIie9IsXgSkJ6MwPgG/xs6stysNAJCYWOc9fJaaFrL0BuCLW7PeUscMwAPmWWzPQr18wH/vXgPXKlMo8qgyhaAUROH1fRmmAOSuhUeVAZ1GDUef/CO1ugBTEoEZGldlzJSrwYDQ88fQJ3stYvQXldlzF4pvYP2olzFt9xfYqo9WjXmy4NgmuivSQgGQjRs3Ks8rKytFo9HIhg0bVNdNmzZNwsPD7+jeO3fulNGjR6uOVVVVyYABA+TKlStiNBpl5MiR8vXXX9f58xUVFVJSUqI8zp8/LwCkpKTkjspBRNRcCq5el+VDxosAkqyPlqBZW5RHsj5aBJDlQ8ZLwdXrd3zvvIJieV8XJZfiE+p83aJQnQggx6fGS15BsZyJmSkCSFZgsKo8ludZgcHy0Oyt8u2JQtmYU6CUTwCRBQtEREyvBSj3try25fi+SdMlaNYW2ZhTIHkFxcpjY06BBM3aInkFxUr5rM/XfNzN3wc5jpKSktv+/d2iWpgacvnyZRiNRvj5+amO+/n54dKlS7d9n2HDhiEnJwfXrl1D165dsXHjRoSGhsLNzQ1JSUkIDw+HiGDo0KEYOXJknffw8PDgDDoicigB7TwxaVA3FIYmYFjsTPVsuBg9ClODMEkD+Nxlq0t9s+ECRo8EDmQjOzAYbWfPMa2rtOIdVBz+Ebo9u1GuC8ewjA8wDICvIRLl+9tBl70HB2QffHpH4uiFEpyx3CwiQrn3r7EzsfqHc6YB6IHB0JnXg/IzGPB9t35o7+kOoPYaUNbrOd3NYp702+UwgcnCxcVF9VxEah1ryLZt2+o9Z6tLkIjIkfksXQwfAH51nbTaL+5u1LtkwXffoVwXjmj9TGyxur7VUxGAmwbeTz55K9AsW2r6c+FC+BhvhRiNVOPE1HjTgpsXSgCYZsGl6aLQv2s7HD1/BYMu/i805vFR+4KCEZeerOqqs9zXLzkJRvPWLfWtRG7BmXZkzWECU8eOHaHRaGq1JhUVFdVqdSIiouZX35IFv9Sx+nd946MAqFqo2ntp8eFTE5BaZax9DwAvBUUipmAtNFbjowAgRR+NuKw1KF3ZA1i8wLSeU9JCnJgajzSfCPQ0z+Czfh3rUGS96OXpGtfWdX1LwM2Lm5bDBCatVouQkBBkZmbiD3/4g3I8MzMTv//97+1YMiIiUsydeyssWZYsMLcK3Y2Adp7YETekziBQVFaBLiuS0ce8HMGvsTOhWbwIcenJODE1HqXhc+GTtBCS/A58zOs5pZkHtlsPBo9dl6vqeqvZVVfXwHHrve6KyipQeuMmfDzdlb3vrFkHlaYKNexebHotKjCVl5fj9OnTyvOzZ88iNzcXHTp0QGBgIGbMmIHx48dj4MCBCAsLw0cffYT8/HxMmTLFjqUmIiKg7lXBC9+YjdN/euWe7hvQzrPuX/ILVwDpyaqlD47OnoOUvIuIMx+HVgsX85IFPVe8g1SYApD1/nWDHvTFWL/nlK43S1ddRmEmfjhVhPvT3oOPp5sy4w9AvZsG18USVADYDDU1Nx22sBWk2L3Y9FpUYPrXv/6Fp556Snk+Y8YMAMDEiRPx2WefYcyYMfjPf/6DBQsW4OLFi+jbty+2bt2KoKAgexWZiOg3r8FVwZOTcOaHc/CMGIf2XtrGfWGj0RSKaiy8maaLwvB+XdDtm+3wNo9r8qgyoNeqVJz8r1gAUELFGRdXDPpsOWL0l1XjnWKy12JQ1hpk66OVa6fu+kK1ObB16Hq+5Ay+btsT96e9pwos1psaA6gz1Fgvv1BfELvdIFVzkDs1nhYVmCIiIiAiDV7zyiuv4JVX7u3/VoiIqPFYVgX3sawKbpmFZ559F5echP8K7wGfdr9r3BeuYxyUJbxtzbuIuH1ZSnj78su3MDg9GVvzLirh7eo1A9J0UZhQ+YtplfHUIOC9xeiU+i7istYo450sMwotmwN38NICmsFK6Ol0/hh65h9BaGAF2loHloULgfRkbDUPMgdMmw6HGb+Hn9Uge8s6VtYtWg/4et9WkLJuvaKm1aICExEROaa6VwWHafadTyvVrLemVFd461lUjn3ZwRh8Pg9xWWuU8Hb1mmkBzU7mZQ90yUmoNi9NkKKPhnb0fwOZPyn3VmYCJichRh8NxOjr3BcP7y1WNjWuOcjcsvxCIYCbs+corUM1W7T6BrRVglTc0IeQsv2nBvfLa8iF4hvKwHVHGcDeEjEwERHRvbvNWW/Noa7wFmsVdny+N+1BGrgyxbSKuLkV6mTyKHiYxztppBqGxLdVrVGAaSag174sxGWtQXWPr5TNgXuueAdnps1CnHk9KBgMWDFkPJbVHGRu3j4mLjkJKw6cx+hNHyotWr+G6pCmi1KtkRWbtQZjv8uHURuIB2L0SuvVheIb6LUqtdb2L5ZAZAlBtzOA/XYGg3MGHgMTERE5m3rCm2rZAw8P+BgMKE2Yi57/72XETJsFD+NNVGu18DBvu6LLP6JqjQJMC196Z+/BTRdXuFltDpwK1NrrbtmgMarWIcAUWHqueAcnVnXBtPRkVPutg6vBYGrhOpCNGI+1qrFUoeePoVP+EYQGlinH6tv+JTZrDc5kr0WsLkoJQXUNYK+vlaq+wMMZeCau9i4AERFRczhdVI6jk6eh2rzpcLVWi/xX4tBrVappDFN8AlwrK1EYnwCduZvNJ2mhaaZfUTlistfCLzkJv4bq4CbVykzAmOy1AEzdatZ73cVkr0W39q0BmAZjW4LKA77eMM6eg0qNmxKuoqOScGJqPOKy1qBT6rsATOHMUg5d/hHluNviRZi2+wv8MGm6KayNGYDUMQNgdHFFXNYaZBRm4kaVEQfOXlHKPeiz5TC6uCrdfX0D2qoG4Z8uKsfRCyWqh2UtqqvXDHh55+f4pnQXtsToaz2+Kd2Fl3d+brNr0NGxhYmIiJyaZSB47LpcxGSvRV9zt5uHwYCSweHQ5R/BiiHjMXr2HACmsUUrDpzHtN1fKGObIt9/Fx7Gm0pL0PLwcej3YQq6rHgPcenJODPhFEYd3a+sSG5ZD+r8X84AQ+eoynO6qBy6l8fAw3gTRnctNOZwtSNxvmmwutU2L3V19ynHp70BrMvFA77eaO+lxVsR4wAAceZZf7HArc2L9dH4W8Q4jDWHpDtZa6r0RpUy4L20vSd8Fi+4dVEdA9udFVuYiIjIqVkWv/zR+L3SknTq3GWlJalcF47Rmz5UupMC2nli9KYPURifgAdHR5q66Yw3cdPFFbr8I0jRR+P9sLGY/NkBRPpEIDswGD2P7kd2YDAifSIwMi1LOd7t0A9Yu2422ntpleDWadRwZaB5zxkblFXJDYlv46Mh0ah2N7WAGd1N3X3nr1xHmi5K1TKmkWr0WpVaq47DMj7Ar6E6xGWtwZllLyphqeeKd1RdZpauurihDwEAUscMwCeTBsLD7VYssMzMi12XizRdFFL00fBJWojS2ebtb8wD2wvjE5BmHpvlzNjCRERETi8gLRlITlIPBjfP4POeNw/eacnqTYPbed6a7WbuZnMzGFSbBVt4lT+HX/e3gf/jYdhiHn90uqgc0UjCge+WIOxANmC+/6E9S9HK3M0W/fRbSB0zAD6TBmL5yy6I27Mag/Lz4Gq1zUtM9lqkIAqv78tQuu9czWOs+mSvVe2XF9DO01TPA9mmMVZVt8ZYbfH1rjW+KDZrDZ4tDUCKT4TSXVh5sxqpYwag16pUbDtyAfenvWe6dl0utInzkZIIxCUtBJLfUba++bWOrW+ccZA4AxMRETm/eha5VJ7XteyBuQVF+bmFC+E9bx76frJCfZ//nwwA6FTHyxZu/Ac6fbJC2V+vlcEAPP00Cj/fAKRlKWOKev3PRygfmQ9d9h5lm5fW7yQhLu1dDMrPU1q2lAVBs9YgOzBYtV+epbzlunB4Z+9RwtWXX76FTsbna22wbOlm+zJwGxCzVzke9uVK+Jm72awHiHfr0BqxVgPbodXiQkw8Tp+9AuDWDD3r9aPqoxoknpgIaDR1z6ZcuND03jQ0C7O5CN2zkpISASAlJSX2LgoRETWGBQtEANOft3O8hryCYgmatUU25hRIXkGxGLVaEUCMWq3kFRTLxpwCCZq1RfIKihu8743wISKAlOnCJa+gWPnZZH20CKCcF/P95emnRQBJ1kdLXkGxXIpPMB0H5FJ8Qp33sD6X3a2f8txSPktdjk+Nl+zAfqrXWz5kvATN2lLr8b4uSpL10fLQ7K3y7YlC5T55BcVyfGq8vK+Lsln32/27vhd38vubLUxEREQ13U2LlJWGBppvG/sq0sxT/5WZavW8XqunIgA3DbyffFK15UmsefuX+1ppoP1hn6m7TqOB63ffqVYotx7A7pechNU/nFO1UqXoo+GuccW05CR0WpFsajmqQ0z2WvTJWgMAqoHt09KT8XAXH7wUFKksrHm6qBxnzPcHAN+pgxtc/Vyp87x5t57XbN1rCZostv2GsIWJiIhqKrh6XWnhsbTgWD8vuHr9ru/be843EjRri9JKVKFxEwEkKzBYgmZtkd5zvlHubylHWZheaeWytEJtzCkwXWc+VqFxU8poaaUqSZirtEQl66NVLUmWFqnlQ8Yrr2dpkSoK1Sl1FRGlxci69UrF0qJkaS1rwpYlizv5/c3A1AgYmIiIqJYm7Gq66zBm1TVo6TK0/NxNd9O541Pj5fjUeFW4qtCHy/Gp8apuRlXXoD5ceYm8gmLleFZgcK0uQ0ugUgWmBQtE5s+/dZ1Wa3pe39+R5fp7xMDUzBiYiIiolqb8hX83YaxGC87yIeOVYGNpObJ+Xmlutapv3JUl+CjjpMyvaXluCViW0CUaTd0tTJZymcdf1RyP1ZRjmxiYmhkDExERNas7DWM1Q4bluVUL1bcnLslDs7eqBoNbuvosgcq6q0/VUlQjjFm6+5R71TMovVZYqlG+Ml14o3dpWmNgamYMTERE1GLV1yITEVHruHVXX1mYvsGQUqtrrUZ3375J05WWJuv7WMY2GRtoSSq4el2WDxmv6tZrKLzdLc6SIyIiIpP6Zvzt3HlrnSMz6wU+vefORV9AWeDTb948wKdVrfucLipHp9R34We1sObajAQMOmdeO8onwrSwpWYwYsyrmltvXny5fyjksUG4GROPAPM9r14zYNmgMXjuYT882EqD6g/+Fx7m64dlfICet7FpcGPj1ihERETOLDGx/qn5c+eqF4VsaDmFBQtU4cqydMKZabPgl5yEFH00erxu2uol7NwRZAcG46Mh0fhkUqiy7UqaLgrZgcFwk2pl8+EvDhTgCc1gPJuyW9nwVynO7Dnwbe+l2tS47ycrVAtqNhe2MBEREZFJQytq1whRAe088WP19/Ax7883LHamacuYGD0KU4OgS07CAdkHn96ROHqhBJU3q/FN6S70yT+irGbeKfVdxCUnYXi/Loj0iajVYtQp9V2lxct6faZOpRWAZnCT/BXUh4GJiIiI7oqPu6t6fz4Lczeej1WLlLL4ZR37+fWZN0+1L57lej/z9UpYM//pN28eYvTRquubGgMTERER3Z07aJHSSDUK4xPgV0d3X2FpBTT7zt7T9U2NgYmIiIiaXKo+Gs/G6NUtUWa/xs5EqiYLz97D9U2NgYmIiIhapNNF5Xd0vCkxMBEREVGLYr15cX1Umxc3AwYmIiIiaha322IU0M4TO+KG4Oo1Q733au+lbbY1mAAGJiIiImpid9NiFNDOs1kDkS0MTERERNSkGmoxKiqrQOmNm/DxdMfVa4Za1zR3S1J9GJiIiIioydXVYnSh+Ab++Nd9uFFlrOenTC1PO+KG2D00MTARERGRXVy9ZsCNKiNSxwyoc7uT03bYM64+DExERERkVw/4eqNvQFt7F6NB3HyXiIiIyAYGJiIiIiIbGJiIiIiIbGBgIiIiIrKBg76JiIjIrlrSnnH1YWAiIiIiu2iJe8bVh4GJiIiI7KIl7hlXHwYmIiIispuWtmdcfTjom4iIiMgGBiYiIiIiGxiYiIiIiGxgYCIiIiKygYGJiIiIyAYGJiIiIiIbGJiIiIiIbGBgIiIiIrKBgYmIiIjIBq703QhEBABQWlpq55IQERHR7bL83rb8Hm8IA1MjKCsrAwB069bNziUhIiKiO1VWVoa2bds2eI2L3E6sogZVV1fj3//+N9q0aQMXF5dGvXdpaSm6deuG8+fPw8fHp1Hv3VKwjs6BdXQOrKNzYB1vj4igrKwM/v7+cHVteJQSW5gagaurK7p27dqkr+Hj4+O0H3oL1tE5sI7OgXV0DqyjbbZaliw46JuIiIjIBgYmIiIiIhs0iYmJifYuBDVMo9EgIiICbm7O24PKOjoH1tE5sI7OgXVsXBz0TURERGQDu+SIiIiIbGBgIiIiIrKBgYmIiIjIBgYmIiIiIhsYmFqwlStXonv37mjVqhVCQkKwd+9eexfpnuzZswfPP/88/P394eLigk2bNqnOiwgSExPh7+8PT09PRERE4NixY3Yq7Z1bsmQJQkND0aZNG/j6+mLUqFE4efKk6prKykrExMSgY8eO8PLywgsvvICCggI7lfjOpaenIzg4WFkoLiwsDN98841y3tHrV5clS5bAxcUFsbGxyjFHr2diYiJcXFxUj86dOyvnHf27aHHhwgWMGzcO9913H1q3bo0BAwbg4MGDynlnqOf9999f6710cXHBq6++CsDxP6sAcPPmTcyZMwfdu3eHp6cnevTogQULFqC6ulq5plneS6EWKSMjQ9zd3WXVqlVy/PhxmT59unh5ecm5c+fsXbS7tnXrVpk9e7asX79eAMjGjRtV55cuXSpt2rSR9evXS15enowZM0a6dOkipaWldirxnRk2bJh8+umncvToUcnNzZURI0ZIYGCglJeXK9dMmTJFAgICJDMzU3JycuSpp56S/v37y82bN+1Y8tu3efNm+cc//iEnT56UkydPSkJCgri7u8vRo0dFxPHrV9P+/fvl/vvvl+DgYJk+fbpy3NHrOX/+fHnkkUfk4sWLyqOoqEg57+jfRRGRK1euSFBQkEyaNEl+/PFHOXv2rOzYsUNOnz6tXOMM9SwqKlK9j5mZmQJAdu7cKSKO/1kVEVm0aJHcd999smXLFjl79qx89dVX4u3tLampqco1zfFeMjC1UI8//rhMmTJFdax3797y5ptv2qlEjatmYKqurpbOnTvL0qVLlWMVFRXStm1b+etf/2qPIt6zoqIiASC7d+8WEZHi4mJxd3eXjIwM5ZoLFy6Iq6ur/POf/7RXMe9Z+/bt5W9/+5vT1a+srEwefPBByczMlCFDhiiByRnqOX/+fOnfv3+d55zluzhr1izR6/X1nneWetY0ffp06dmzp1RXVzvFZ1VEZMSIETJ58mTVsRdffFHGjRsnIs33XrJLrgUyGAw4ePAghg4dqjo+dOhQfP/993YqVdM6e/YsLl26pKqzh4cHhgwZ4rB1LikpAQB06NABAHDw4EFUVVWp6ujv74++ffs6ZB2NRiMyMjJw7do1hIWFOV39Xn31VYwYMQLPPvus6riz1PPUqVPw9/dH9+7dMXbsWPz8888AnOe7uHnzZgwcOBB//OMf4evri0cffRSrVq1SzjtLPa0ZDAasXr0akydPhouLi9N8VvV6Pb799lv89NNPAIDDhw8jKysLw4cPB9B876XzLv/pwC5fvgyj0Qg/Pz/VcT8/P1y6dMlOpWpalnrVVedz587Zo0j3REQwY8YM6PV69O3bF4CpjlqtFu3bt1dd62jva15eHsLCwlBRUQFvb29s3LgRDz/8MHJzc52ifgCQkZGBnJwcHDhwoNY5Z3gfn3jiCXz++ed46KGHUFhYiEWLFmHw4ME4duyY03wXf/75Z6Snp2PGjBlISEjA/v37MW3aNHh4eGDChAlOU09rmzZtQnFxMSZNmgTAOT6rADBr1iyUlJSgd+/e0Gg0MBqNWLx4MaKiogA03+8PBqYWzMXFRfVcRGodczbOUufXXnsNR44cQVZWls1rHa2OvXr1Qm5uLoqLi7F+/XpMnDgRu3fvrvd6R6vf+fPnMX36dGzfvh2tWrW67Z9zpHpGRkYq/92vXz+EhYWhZ8+e+Pvf/45BgwYBcPzvYnV1NQYOHIikpCQAwKOPPopjx44hPT0dEyZMUK5z9Hpa+/i/GxFTAAAFVklEQVTjjxEZGQl/f/8Gr3O0Oq5btw6rV6/Gl19+iUceeQS5ubmIjY2Fv78/Jk6cqFzX1O8lu+RaoI4dO0Kj0dT6P4CioqJaCdpZWGboOEOdY2JisHnzZuzcuRNdu3ZVjnfu3BkGgwFXr15VXe9oddRqtXjggQcwcOBALFmyBP3798fy5cudpn4HDx5EUVERQkJC4ObmBjc3N+zevRsrVqyAm5sb/Pz8nKKe1ry8vNCvXz+cOnXKab6LXbp0wcMPP6w61qdPH+Tn5wNwrn9zAODcuXPYsWMHXnrpJeWYs3wn33jjDbz55psYO3Ys+vXrh/Hjx+P111/HkiVLADTfe8nA1AJptVqEhIQgMzNTdTwzMxODBw+2U6maVvfu3dG5c2dVnQ0GA3bv3u0wdRYRvPbaa9iwYQO+++47dO/eXXU+JCQE7u7uqjpevHgRR48edZg61kVEUFlZ6TT1e+aZZ5CXl4fc3FzlMXDgQERHRyv/7Qz1tFZZWYkTJ06gS5cuTvFdBACdTldrWY+ffvoJQUFBAJzj3xxrn376KXx9fTFixAjlmLN8J69fvw5XV3Vc0Wg0yrICzfZeNtrwcWpUlmUFPv74Yzl+/LjExsaKl5eX/PLLL/Yu2l0rKyuTQ4cOyaFDhwSALFu2TA4dOqQslbB06VJp27atbNiwQfLy8iQqKsqhpvhOnTpV2rZtK7t27VJN871+/bpyzZQpU6Rr166yY8cOycnJkaefftqhpvi+9dZbsmfPHjl79qwcOXJEEhISxNXVVbZv3y4ijl+/+ljPkhNx/HrGxcXJrl275Oeff5YffvhBRo4cKW3atFH+fXH076KIaUkINzc3Wbx4sZw6dUrWrFkjrVu3ltWrVyvXOEM9RUSMRqMEBgbKrFmzap1z9M+qiMjEiRMlICBAWVZgw4YN0rFjR5k5c6ZyTXO8lwxMLdgHH3wgQUFBotVq5bHHHlOmpzuqnTt3CoBaj4kTJ4qIaWro/PnzpXPnzuLh4SHh4eGSl5dn30LfgbrqBkA+/fRT5ZobN27Ia6+9Jh06dBBPT08ZOXKk5Ofn26/Qd2jy5MnKZ7JTp07yzDPPKGFJxPHrV5+agcnR62lZo8bd3V38/f3lxRdflGPHjinnHf27aPH1119L3759xcPDQ3r37i0fffSR6ryz1HPbtm0CQE6ePFnrnKN/VkVESktLZfr06RIYGCitWrWSHj16yOzZs6WyslK5pjneSxcRkcZrryIiIiJyPhzDRERERGQDAxMRERGRDQxMRERERDYwMBERERHZwMBEREREZAMDExEREZENDExERERENjAwEREREdnAwERERERkAwMTERERkQ0MTEREt+GNN97AyJEj7V0MIrIT7iVHRHQbiouLodFo0KZNG3sXhYjsgIGJiIiIyAZ2yRER2XD58mW4uLjg2LFj9i4KEdkJAxMRkQ2HDx+Gh4cHevXqZe+iEJGdMDAREdlw+PBhPPLII3Bzc7N3UYjIThiYiIhsyM3NxYABA+xdDCKyIwYmIiIbDh8+jP79+9u7GERkRwxMREQNqKqqwokTJxiYiH7jGJiIiBpw/PhxVFVVMTAR/cYxMBERNSA3NxdBQUFo166dvYtCRHbEwERE1IADBw7g8ccft3cxiMjOGJiIiOpQUVGBnJwcrF+/HsOGDbN3cYjIzhiYiIjqkJqaiueffx6jRo3ChAkT7F0cIrIz7iVHREREZANbmIiIiIhsYGAiIiIisoGBiYiIiMgGBiYiIiIiGxiYiIiIiGxgYCIiIiKygYGJiIiIyAYGJiIiIiIbGJiIiIiIbGBgIiIiIrLh/wBZf9Om+0iqpwAAAABJRU5ErkJggg==", "text/plain": [ "Figure(PyObject
)" ] }, "metadata": {}, "output_type": "display_data" }, { "name": "stdout", "output_type": "stream", "text": [ "L1 NORM of Householder Q'Q - I: 7.997590220285011e-15\n" ] } ], "source": [ "QRH = qr(A)\n", "QH, RH = QRH.Q, QRH.R\n", "semilogy(abs.(diag(RH)), \"s\", mfc=\"none\") # need abs because householder has arbitrary phase factor in R\n", "semilogy(diag(RM), \"rx\")\n", "legend([\"Householder QR\",\"modified GS\"], loc=\"upper right\")\n", "ylabel(L\"r_{jj}\")\n", "xlabel(L\"j\")\n", "\n", "println(\"L1 NORM of Householder Q'Q - I: \", opnorm(QH'*QH - I, 1))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "On top of that, Householder is faster than Gram–Schmidt, so what's not to love? But modified Gram–Schmidt is still useful to know since it is the basis of the Arnoldi and GMRES iterative algorithms for large (usually sparse) systems." ] } ], "metadata": { "@webio": { "lastCommId": null, "lastKernelId": null }, "kernelspec": { "display_name": "Julia 1.3.0", "language": "julia", "name": "julia-1.3" }, "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.3.1" } }, "nbformat": 4, "nbformat_minor": 1 }