{ "cells": [ { "cell_type": "markdown", "id": "7f6234eb", "metadata": { "editable": true }, "source": [ "\n", "\n", "# Demo - Rayleigh Benard\n", "**Mikael Mortensen** (email: `mikaem@math.uio.no`), Department of Mathematics, University of Oslo.\n", "\n", "Date: **November 21, 2019**\n", "\n", "**Summary.** Rayleigh-Benard convection arise\n", "due to temperature gradients in a fluid. The governing equations are\n", "Navier-Stokes coupled (through buoyancy) with an additional temperature\n", "equation derived from the first law of thermodynamics, using a linear\n", "correlation between density and temperature.\n", "\n", "This is a demonstration of how the Python module [shenfun](https://github.com/spectralDNS/shenfun) can be used to solve for\n", "these Rayleigh-Benard cells in a 2D channel with two walls of\n", "different temperature in one direction, and periodicity in the other direction.\n", "The solver described runs with MPI\n", "without any further considerations required from the user.\n", "Note that there is also a more physically realistic [3D solver](https://github.com/spectralDNS/shenfun/blob/master/demo/RayleighBenard.py).\n", "The solver described in this demo has been implemented in a class in the\n", "[RayleighBenard2D.py](https://github.com/spectralDNS/shenfun/blob/master/demo/RayleighBenard2D.py)\n", "module in the demo folder of shenfun. Below is an example solution, which has been run at a very high\n", "Rayleigh number (*Ra*).\n", "\n", "\n", "\n", "\n", "\n", "

Figure 1: Temperature fluctuations in the Rayleigh Benard flow. The top and bottom walls are kept at different temperatures and this sets up the Rayleigh-Benard convection. The simulation is run at Ra =1,000,000, Pr =0.7 with 256 and 512 quadrature points in x and y-directions, respectively.

\n", "" ] }, { "cell_type": "markdown", "id": "177a6e73", "metadata": { "editable": true }, "source": [ "## The Rayleigh Bénard equations\n", "\n", "\n", "The governing equations solved in domain $\\Omega=(-1, 1)\\times [0, 2\\pi)$ are" ] }, { "cell_type": "markdown", "id": "ad1449b2", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\frac{\\partial \\boldsymbol{u}}{\\partial t} + (\\boldsymbol{u} \\cdot \\nabla) \\boldsymbol{u} = - \\nabla p + \\sqrt{\\frac{Pr}{Ra}} \\nabla^2 \\boldsymbol{u} + T \\boldsymbol{i}, \\label{eq:momentum} \\tag{1}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "59b6c4c7", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " \\frac{\\partial T}{\\partial t} +\\boldsymbol{u} \\cdot \\nabla T = \\frac{1}{\\sqrt{RaPr}} \\nabla^2 T, \\label{eq:T} \\tag{2}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "5d7e6701", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " \\nabla \\cdot \\boldsymbol{u} = 0, \\label{eq:div} \\tag{3}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "2a91fb26", "metadata": { "editable": true }, "source": [ "where $\\boldsymbol{u}(x, y, t) (= u\\boldsymbol{i} + v\\boldsymbol{j})$ is the velocity vector, $p(x, y, t)$ is pressure, $T(x, y, t)$ is the temperature, and $\\boldsymbol{i}$ and\n", "$\\boldsymbol{j}$ are the unity vectors for the $x$ and $y$-directions, respectively.\n", "\n", "The equations are complemented with boundary conditions $\\boldsymbol{u}(\\pm 1, y, t) = (0, 0), \\boldsymbol{u}(x, 2 \\pi, t) = \\boldsymbol{u}(x, 0, t), T(-1, y, t) = 1, T(1, y, t) = 0, T(x, 2 \\pi, t) = T(x, 0, t)$.\n", "Note that these equations have been non-dimensionalized according to [[pandey18]](#pandey18), using dimensionless\n", "Rayleigh number $Ra=g \\alpha \\Delta T h^3/(\\nu \\kappa)$ and Prandtl number $Pr=\\nu/\\kappa$. Here\n", "$g \\boldsymbol{i}$ is the vector accelleration of gravity, $\\Delta T$ is the temperature difference between\n", "the top and bottom walls, $h$ is the hight of the channel in $x$-direction, $\\nu$ is the\n", "dynamic viscosity coefficient, $\\kappa$ is the heat transfer coefficient and $\\alpha$ is the\n", "thermal expansion coefficient. Note that the\n", "governing equations have been non-dimensionalized using the free-fall velocityscale\n", "$U=\\sqrt{g \\alpha \\Delta T h}$. See [[pandey18]](#pandey18) for more details.\n", "\n", "The governing equations contain a non-trivial coupling between velocity, pressure and temperature.\n", "This coupling can be simplified by eliminating the pressure from the equation for the wall-normal velocity\n", "component $u$. We accomplish this by taking the Laplace of the momentum equation in wall normal\n", "direction, using the pressure from the divergence of the momentum equation\n", "$\\nabla^2 p = -\\nabla \\cdot \\boldsymbol{H}+\\partial T/\\partial x$, where\n", "$\\boldsymbol{H} = (H_x, H_y) = (\\boldsymbol{u} \\cdot \\nabla) \\boldsymbol{u}$" ] }, { "cell_type": "markdown", "id": "f838dca0", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\frac{\\partial \\nabla^2 {u}}{\\partial t} = \\frac{\\partial^2 H_y}{\\partial x \\partial y} - \\frac{\\partial^2 H_x}{\\partial y\\partial y} + \\sqrt{\\frac{Pr}{Ra}} \\nabla^4 {u} + \\frac{\\partial^2 T}{\\partial y^2} . \\label{eq:rb:u} \\tag{4}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "ca6ea1f7", "metadata": { "editable": true }, "source": [ "This equation is solved with $u(\\pm 1,y,t) = \\partial u/\\partial x(\\pm 1,y,t) = 0$, where the latter follows from the\n", "divergence constraint. In summary, we now seem to have the following equations to solve:" ] }, { "cell_type": "markdown", "id": "0afef1df", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\frac{\\partial \\nabla^2 {u}}{\\partial t} = \\frac{\\partial^2 H_y}{\\partial x \\partial y} - \\frac{\\partial^2 H_x}{\\partial y\\partial y} + \\sqrt{\\frac{Pr}{Ra}} \\nabla^4 {u} + \\frac{\\partial^2 T}{\\partial y^2}, \\label{eq:rb:u2} \\tag{5}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "3d483309", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " \\frac{\\partial v}{\\partial t} + H_y = - \\frac{\\partial p}{\\partial y} + \\sqrt{\\frac{Pr}{Ra}} \\nabla^2 v, \\label{eq:v} \\tag{6}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "9f1801e4", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " \\frac{\\partial T}{\\partial t} +\\boldsymbol{u} \\cdot \\nabla T = \\frac{1}{\\sqrt{RaPr}} \\nabla^2 T, \\label{eq:T2} \\tag{7}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "2d815654", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " \\nabla \\cdot \\boldsymbol{u} = 0 \\label{eq:div2} \\tag{8}.\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "87c1dd49", "metadata": { "editable": true }, "source": [ "However, we note that Eqs. ([5](#eq:rb:u2)) and ([7](#eq:T2)) and ([8](#eq:div2)) do not depend on pressure, and,\n", "apparently, on each time step we can solve ([5](#eq:rb:u2)) for $u$, then ([8](#eq:div2)) for $v$ and finally ([7](#eq:T2)) for $T$.\n", "So what do we need ([6](#eq:v)) for? It appears to have become redundant from the elimination of the\n", "pressure from Eq. ([5](#eq:rb:u2)). It turns out that this is actually almost completely true, but\n", "([5](#eq:rb:u2)), ([7](#eq:T2)) and ([8](#eq:div2)) can only provide closure for all but one of the\n", "Fourier coefficients. To see this it is necessary to introduce some discretization and basis functions\n", "that will be used to solve the problem. To this end we use $P_N$, which is the set of all real polynomials\n", "of degree less than or equal to N and introduce the following finite-dimensional approximation spaces" ] }, { "cell_type": "markdown", "id": "469f8896", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " V_N^B(x) = \\{v \\in P_N | v(\\pm 1) = v´(\\pm 1) = 0\\}, \\label{eq:VB} \\tag{9} \n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "c8f804cb", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " V_N^D(x) = \\{v \\in P_N | v(\\pm 1) = 0\\}, \\label{eq:VD} \\tag{10} \n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "f3d0a11f", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " V_N^T(x) = \\{v \\in P_N | v(-1) = 0, v(1) = 1\\}, \\label{eq:VT} \\tag{11} \n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "eae8e9d1", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " V_N^W(x) = \\{v \\in P_N\\}, \\label{eq:VW} \\tag{12} \n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "7be2b88e", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " V_M^F(y) = \\{\\exp(\\imath l y) | l \\in [-M/2, -M/2+1, \\ldots M/2-1]\\}. \\label{eq:VF} \\tag{13}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "a5e83417", "metadata": { "editable": true }, "source": [ "Here $\\text{dim}(V_N^B) = N-4, \\text{dim}(V_N^D) = \\text{dim}(V_N^W) = \\text{dim}(V_N^T) = N-2$\n", "and $\\text{dim}(V_M^F)=M$. We note that\n", "$V_N^B, V_N^D, V_N^W$ and $V_N^T$ can be used to approximate $u, v, T$ and $p$, respectively, in the $x$-direction.\n", "Also note that for $V_M^F$ it is assumed that $M$ is an even number.\n", "\n", "We can now choose basis functions for the spaces, using Shen's composite bases for either Legendre or\n", "Chebyshev polynomials. For the Fourier space the basis functions are already given. We leave the actual choice\n", "of basis as an implementation option for later. For now we use $\\phi^B(x), \\phi^D(x), \\phi^W$ and $\\phi^T(x)$\n", "as common notation for basis functions in spaces $V_N^B, V_N^D, V_N^W$ and $V_N^T$, respectively.\n", "\n", "To get the required approximation spaces for the entire domain we use tensor products of the\n", "one-dimensional spaces in ([9](#eq:VB))-([13](#eq:VF))" ] }, { "cell_type": "markdown", "id": "38df0fa1", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " W_{BF} = V_N^B \\otimes V_M^F, \\label{eq:WBF} \\tag{14} \n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "ba2319bf", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " W_{DF} = V_N^D \\otimes V_M^F, \\label{eq:WDF} \\tag{15} \n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "97c4f4bc", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " W_{TF} = V_N^T \\otimes V_M^F, \\label{eq:WTF} \\tag{16} \n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "25acfc6e", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " W_{WF} = V_N^W \\otimes V_M^F. \\label{eq:WWF} \\tag{17}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "827f02da", "metadata": { "editable": true }, "source": [ "Space $W_{BF}$ has 2D tensor product basis functions $\\phi_k^B(x) \\exp (\\imath l y)$ and\n", "similar for the others. All in all\n", "we get the following approximations for the unknowns" ] }, { "cell_type": "markdown", "id": "fc1fc414", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " u_N(x, y, t) = \\sum_{k \\in \\boldsymbol{k}_B} \\sum_{l \\in \\boldsymbol{l}} \\hat{u}_{kl}(t) \\phi_k^B(x) \\exp(\\imath l y), \n", "\\label{_auto1} \\tag{18}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "b3e1db16", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " v_N(x, y, t) = \\sum_{k \\in \\boldsymbol{k}_D} \\sum_{l \\in \\boldsymbol{l}} \\hat{v}_{kl}(t) \\phi_k^D(x) \\exp(\\imath l y), \n", "\\label{_auto2} \\tag{19}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "bfc89552", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " p_N(x, y, t) = \\sum_{k \\in \\boldsymbol{k}_W} \\sum_{l \\in \\boldsymbol{l}} \\hat{p}_{kl}(t) \\phi_k^W(x) \\exp(\\imath l y), \n", "\\label{_auto3} \\tag{20}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "d37af4a7", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " T_N(x, y, t) = \\sum_{k \\in \\boldsymbol{k}_T} \\sum_{l \\in \\boldsymbol{l}} \\hat{T}_{kl}(t) \\phi_k^T(x) \\exp(\\imath l y),\n", "\\label{_auto4} \\tag{21}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "3d0a241c", "metadata": { "editable": true }, "source": [ "where $\\boldsymbol{k}_{x} = \\{0, 1, \\ldots \\text{dim}(V_N^x)-1\\}, \\, \\text{for} \\, x\\in(B, D, W, T)$\n", "and $\\boldsymbol{l} = \\{-M/2, -M/2+1, \\ldots, M/2-1\\}$.\n", "Note that since the problem is defined in real space we will have Hermitian symmetry. This means\n", "that $\\hat{u}_{k, l} = \\overline{\\hat{u}}_{k, -l}$, with an overbar being a complex conjugate,\n", "and similar for $\\hat{v}_{kl}, \\hat{p}_{kl}$ and\n", "$\\hat{T}_{kl}$. For this reason we can get away with\n", "solving for only the positive $l$'s, as long as we remember that the sum in the end goes over both positive\n", "and negative $l's$. This is actually automatically taken care of by the FFT provider and is\n", "not much of an additional complexity in the implementation. So from now on $\\boldsymbol{l} = \\{0, 1, \\ldots, M/2\\}$.\n", "\n", "We can now take a look at why Eq. ([6](#eq:v)) is needed. If we first solve ([5](#eq:rb:u2)) for\n", "$\\hat{u}_{kl}(t), (k, l) \\in \\boldsymbol{k}_B \\times \\boldsymbol{l}$, then we can use ([8](#eq:div2)) to\n", "solve for $\\hat{v}_{kl}(t)$. But here there is a problem. We can see this by creating the variational\n", "form required to solve ([8](#eq:div2)) by the spectral Galerkin method. To this end make $v=v_N$ in ([8](#eq:div2))\n", "a trial function, use $u=u_N$ a known function and take the weighted inner product over the\n", "domain using test function $q \\in W_{DF}$" ] }, { "cell_type": "markdown", "id": "4b48d850", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\left < \\frac{\\partial u_N}{\\partial x} + \\frac{\\partial v_N}{\\partial y}, q \\right > _w = 0.\n", "\\label{_auto5} \\tag{22}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "6505d4f0", "metadata": { "editable": true }, "source": [ "Here we are using the inner product notation" ] }, { "cell_type": "markdown", "id": "9150f8b5", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\left < a, b \\right > _w = \\int_{-1}^1 \\int_0^{2\\pi} a \\overline{b} dx_wdy_w \\left(\\approx \\sum_{i}\\sum_{j} a(x_i, y_j) \\overline{b}(x_i, y_j) w_i w_j\\right),\n", "\\label{_auto6} \\tag{23}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "d98142bf", "metadata": { "editable": true }, "source": [ "where the exact form of the\n", "weighted scalar product depends on the chosen basis; Legendre has $dx_w=dx$, Chebyshev\n", "$dx_w = dx/\\sqrt{1-x^2}$ and Fourier $dy_w=dy/2/\\pi$. The bases also have associated quadrature weights\n", "$\\{w_i \\}_{i=0}^{N-1}$ and $\\{w_j\\}_{j=0}^{M-1}$ that are used to approximate the integrals.\n", "\n", "Inserting now for the known $u_N$, the unknown $v_N$, and $q=\\phi_m^D(x) \\exp(\\imath n y)$ the\n", "continuity equation becomes" ] }, { "cell_type": "markdown", "id": "2ac1d624", "metadata": { "editable": true }, "source": [ "$$\n", "\\int_{-1}^1 \\int_{0}^{2\\pi} \\frac{\\partial}{\\partial x} \\left(\\sum_{k \\in \\boldsymbol{k}_B} \\sum_{l \\in \\boldsymbol{l}} \\hat{u}_{kl}(t) \\phi_k^B(x) \\exp(\\imath l y) \\right) \\phi_m^D(x) \\exp(-\\imath n y) dx_w dy_w + \\\\ \n", " \\int_{-1}^1 \\int_{0}^{2\\pi} \\frac{\\partial}{\\partial y} \\left(\\sum_{k \\in \\boldsymbol{k}_D} \\sum_{l \\in \\boldsymbol{l}} \\hat{v}_{kl}(t) \\phi_k^D(x) \\exp(\\imath l y) \\right) \\phi_m^D(x) \\exp(-\\imath n y) dx_w dy_w = 0.\n", "$$" ] }, { "cell_type": "markdown", "id": "1a8e691c", "metadata": { "editable": true }, "source": [ "The $x$ and $y$ domains are separable, so we can rewrite as" ] }, { "cell_type": "markdown", "id": "1cf1221d", "metadata": { "editable": true }, "source": [ "$$\n", "\\sum_{k \\in \\boldsymbol{k}_B} \\sum_{l \\in \\boldsymbol{l}} \\int_{-1}^1 \\frac{\\partial \\phi_k^B(x)}{\\partial x} \\phi_m^D(x) dx_w \\int_{0}^{2\\pi} \\exp(\\imath l y) \\exp(-\\imath n y) dy_w \\hat{u}_{kl} + \\\\ \n", " \\sum_{k \\in \\boldsymbol{k}_D} \\sum_{l \\in \\boldsymbol{l}} \\int_{-1}^1 \\phi_k^D(x) \\phi_m^D(x) dx_w \\int_{0}^{2\\pi} \\frac{\\partial \\exp(\\imath l y)}{\\partial y} \\exp(-\\imath n y) dy_w \\hat{v}_{kl} = 0.\n", "$$" ] }, { "cell_type": "markdown", "id": "4df1a6eb", "metadata": { "editable": true }, "source": [ "Now perform some exact manipulations in the Fourier direction and introduce the\n", "1D inner product notation for the $x$-direction" ] }, { "cell_type": "markdown", "id": "a3eb4492", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\left(a, b\\right)_w = \\int_{-1}^1 a(x) b(x) dx_w \\left(\\approx \\sum_{j = 0}^{N-1} a(x_j)b(x_j) w_j\\right).\n", "\\label{_auto7} \\tag{24}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "1b5d5113", "metadata": { "editable": true }, "source": [ "By also simplifying the notation using summation of repeated indices,\n", "we get the following equation" ] }, { "cell_type": "markdown", "id": "26b6a18f", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\delta_{ln} \\left(\\frac{\\partial \\phi_k^B}{\\partial x}, \\phi_m^D \\right)_w \\hat{u}_{kl}\n", " + \\imath l \\delta_{ln} \\left(\\phi_k^D, \\phi_m^D \\right)_w \\hat{v}_{kl} = 0.\n", "\\label{_auto8} \\tag{25}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "e486aa6d", "metadata": { "editable": true }, "source": [ "Now $l$ must equal $n$ and we can simplify some more" ] }, { "cell_type": "markdown", "id": "f4c027ba", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\left(\\frac{\\partial \\phi_k^B}{\\partial x}, \\phi_m^D \\right)_w \\hat{u}_{kl}\n", " + \\imath l \\left(\\phi_k^D, \\phi_m^D \\right)_w \\hat{v}_{kl} = 0. \\label{eq:div3} \\tag{26}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "fb8eb364", "metadata": { "editable": true }, "source": [ "We see that this equation can be solved for\n", "$\\hat{v}_{kl} \\text{ for } (k, l) \\in \\boldsymbol{k}_D \\times [1, 2, \\ldots, M/2]$, but try with\n", "$l=0$ and you hit division by zero, which obviously is not allowed. And this is the reason\n", "why Eq. ([6](#eq:v)) is still needed, to solve for $\\hat{v}_{k,0}$! Fortunately,\n", "since $\\exp(\\imath 0 y) = 1$, the pressure derivative $\\frac{\\partial p}{\\partial y} = 0$,\n", "and as such the pressure is still not required. When used only for\n", "Fourier coefficient 0, Eq. ([6](#eq:v)) becomes" ] }, { "cell_type": "markdown", "id": "e1966ab1", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", "\\frac{\\partial v}{\\partial t} + N_y = \\sqrt{\\frac{Pr}{Ra}} \\frac{\\partial^2 v}{\\partial x^2}. \\label{eq:vx} \\tag{27}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "69a9f741", "metadata": { "editable": true }, "source": [ "There is still one more revelation to be made from Eq. ([26](#eq:div3)). When $l=0$ we get" ] }, { "cell_type": "markdown", "id": "440cd498", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\left(\\frac{\\partial \\phi_k^B}{\\partial x}, \\phi_m^D \\right)_w \\hat{u}_{k,0} = 0,\n", "\\label{_auto9} \\tag{28}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "200607d9", "metadata": { "editable": true }, "source": [ "which is trivially satisfied if $\\hat{u}_{k,0}=0$ for $k\\in\\boldsymbol{k}_B$. Bottom line is\n", "that we only need to solve Eq. ([5](#eq:rb:u2)) for $l \\in \\boldsymbol{l}/\\{0\\}$, whereas we can use\n", "directly $\\hat{u}_{k,0}=0 \\text{ for } k \\in \\boldsymbol{k}_B$.\n", "\n", "To sum up, with the solution known at $t = t - \\Delta t$, we solve\n", "\n", ":::{table}\n", ":widths: auto\n", "\n", "| Equation | For unknown | With indices |\n", "| :------------------------------------------------------ | :-------------------------- | :------------------------------------------------------------------ |\n", "| (5) | {math}`\\hat{u}_{kl}(t)` | {math}`(k, l) \\in \\boldsymbol{k}_B \\times \\boldsymbol{l}` |\n", "| (8) | {math}`\\hat{v}_{kl}(t)` | {math}`(k, l) \\in \\boldsymbol{k}_D \\times \\boldsymbol{l}/\\{0\\}` |\n", "| (27) | {math}`\\hat{v}_{kl}(t)` | {math}`(k, l) \\in \\boldsymbol{k}_D \\times \\{0\\}` |\n", "| (7) | {math}`\\hat{T}_{kl}(t)` | {math}`(k, l) \\in \\boldsymbol{k}_T \\times \\boldsymbol{l}` |\n", "\n", ":::\n" ] }, { "cell_type": "markdown", "id": "b943e520", "metadata": { "editable": true }, "source": [ "## Temporal discretization\n", "\n", "The governing equations are integrated in time using any one of the time steppers available in\n", "shenfun. There are several possible IMEX Runge Kutta methods, see [integrators.py](https://github.com/spectralDNS/shenfun/blob/master/shenfun/utilities/integrators.py).\n", "The time steppers are used for any generic equation" ] }, { "cell_type": "markdown", "id": "b6237baa", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\frac{\\partial \\psi}{\\partial t} = \\mathcal{N} + \\mathcal{L}\\psi \\label{eq:genericpsi} \\tag{29},\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "f42997a8", "metadata": { "editable": true }, "source": [ "where $\\mathcal{N}$ and $\\mathcal{L}$ represents the nonlinear and linear contributions, respectively.\n", "The timesteppers are provided with $\\psi, \\mathcal{L}$ and $\\mathcal{N}$, and possibly some tailored\n", "linear algebra solvers, and solvers are then further assembled under the hood.\n", "\n", "All the timesteppers split one time step into one or several stages.\n", "The classes [IMEXRK222](https://shenfun.readthedocs.io/en/latest/shenfun.utilities.html#shenfun.utilities.integrators.IMEXRK222), [IMEXRK3](https://shenfun.readthedocs.io/en/latest/shenfun.utilities.html#shenfun.utilities.integrators.IMEXRK3) and [IMEXRK443](https://shenfun.readthedocs.io/en/latest/shenfun.utilities.html#shenfun.utilities.integrators.IMEXRK443)\n", "have 2, 3 and 4 steps, respectively, and the cost is proportional." ] }, { "cell_type": "markdown", "id": "1b1d02c6", "metadata": { "editable": true }, "source": [ "## Implementation\n", "\n", "To get started we need instances of the approximation spaces discussed in\n", "Eqs. ([9](#eq:VB)) - ([17](#eq:WWF)). When the spaces are created we also need\n", "to specify the family and the dimension of each space. Here we simply\n", "choose Chebyshev and Fourier with 100 and 256 quadrature points in $x$ and\n", "$y$-directions, respectively. We could replace 'Chebyshev' by 'Legendre',\n", "but the former is known to be faster due to the existence of fast transforms." ] }, { "cell_type": "code", "execution_count": 1, "id": "4917abcd", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:15.145923Z", "iopub.status.busy": "2026-06-03T08:50:15.145558Z", "iopub.status.idle": "2026-06-03T08:50:15.621801Z", "shell.execute_reply": "2026-06-03T08:50:15.621477Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "from shenfun import *\n", "\n", "N, M = 64, 128\n", "family = 'Chebyshev'\n", "VB = FunctionSpace(N, family, bc=(0, 0, 0, 0))\n", "VD = FunctionSpace(N, family, bc=(0, 0))\n", "VW = FunctionSpace(N, family)\n", "VT = FunctionSpace(N, family, bc=(1, 0))\n", "VF = FunctionSpace(M, 'F', dtype='d')" ] }, { "cell_type": "markdown", "id": "b7122e43", "metadata": { "editable": true }, "source": [ "And then we create tensor product spaces by combining these bases (like in Eqs. ([14](#eq:WBF))-([17](#eq:WWF)))." ] }, { "cell_type": "code", "execution_count": 2, "id": "32a0a207", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:15.623566Z", "iopub.status.busy": "2026-06-03T08:50:15.623319Z", "iopub.status.idle": "2026-06-03T08:50:15.980695Z", "shell.execute_reply": "2026-06-03T08:50:15.980413Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "W_BF = TensorProductSpace(comm, (VB, VF)) # Wall-normal velocity\n", "W_DF = TensorProductSpace(comm, (VD, VF)) # Streamwise velocity\n", "W_WF = TensorProductSpace(comm, (VW, VF)) # No bc\n", "W_TF = TensorProductSpace(comm, (VT, VF)) # Temperature\n", "BD = VectorSpace([W_BF, W_DF]) # Velocity vector\n", "DD = VectorSpace([W_DF, W_DF]) # Convection vector\n", "W_DFp = W_DF.get_dealiased(padding_factor=1.5)\n", "BDp = BD.get_dealiased(padding_factor=1.5)" ] }, { "cell_type": "markdown", "id": "3e9e1f8d", "metadata": { "editable": true }, "source": [ "Here the `VectorSpae` create mixed tensor product spaces by the\n", "Cartesian products `BD = W_BF` $\\times$ `W_DF` and `DD = W_DF` $\\times$ `W_DF`.\n", "These mixed space will be used to hold the velocity and convection vectors,\n", "but we will not solve the equations in a coupled manner and continue in the\n", "segregated approach outlined above.\n", "\n", "We also need containers for the computed solutions. These are created as" ] }, { "cell_type": "code", "execution_count": 3, "id": "be8b51ce", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:15.982222Z", "iopub.status.busy": "2026-06-03T08:50:15.982144Z", "iopub.status.idle": "2026-06-03T08:50:15.984084Z", "shell.execute_reply": "2026-06-03T08:50:15.983874Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "u_ = Function(BD) # Velocity vector, two components\n", "T_ = Function(W_TF) # Temperature\n", "H_ = Function(DD) # Convection vector\n", "uT_ = Function(BD) # u times T" ] }, { "cell_type": "markdown", "id": "4c90842d", "metadata": { "editable": true }, "source": [ "### Wall-normal velocity equation\n", "\n", "We implement Eq. ([5](#eq:rb:u2)) using a generic time stepper.\n", "To this end we first need to declare some test- and trial functions, as well as\n", "some model constants and the length of the time step. We store the\n", "PDE time stepper in a dictionary called `pdes` just for convenience:" ] }, { "cell_type": "code", "execution_count": 4, "id": "b3ee1794", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:15.985263Z", "iopub.status.busy": "2026-06-03T08:50:15.985195Z", "iopub.status.idle": "2026-06-03T08:50:16.011173Z", "shell.execute_reply": "2026-06-03T08:50:16.010927Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "\n", "# Specify viscosity and time step size using dimensionless Ra and Pr\n", "Ra = 1000000\n", "Pr = 0.7\n", "nu = np.sqrt(Pr/Ra)\n", "kappa = 1./np.sqrt(Pr*Ra)\n", "dt = 0.025\n", "\n", "# Choose timestepper and create instance of class\n", "\n", "PDE = IMEXRK3 # IMEX222, IMEXRK443\n", "\n", "v = TestFunction(W_BF) # The space we're solving for u in\n", "\n", "pdes = {\n", " 'u': PDE(v, # test function\n", " div(grad(u_[0])), # u\n", " lambda f: nu*div(grad(f)), # linear operator on u\n", " [Dx(Dx(H_[1], 0, 1), 1, 1)-Dx(H_[0], 1, 2), Dx(T_, 1, 2)],\n", " dt=dt,\n", " solver=chebyshev.la.Biharmonic if family == 'Chebyshev' else la.SolverGeneric1ND,\n", " latex=r\"\\frac{\\partial \\nabla^2 u}{\\partial t} = \\nu \\nabla^4 u + \\frac{\\partial^2 N_y}{\\partial x \\partial y} - \\frac{\\partial^2 N_x}{\\partial y^2}\")\n", "}\n", "pdes['u'].assemble()" ] }, { "cell_type": "markdown", "id": "7b8ac028", "metadata": { "editable": true }, "source": [ "Notice the one-to-one resemblance with ([5](#eq:rb:u2)).\n", "\n", "The right hand side depends on the convection vector $\\boldsymbol{H}$, which can be computed in many different ways.\n", "Here we will make use of\n", "the identity $(\\boldsymbol{u} \\cdot \\nabla) \\boldsymbol{u} = -\\boldsymbol{u} \\times (\\nabla \\times \\boldsymbol{u}) + 0.5 \\nabla\\boldsymbol{u} \\cdot \\boldsymbol{u}$,\n", "where $0.5 \\nabla \\boldsymbol{u} \\cdot \\boldsymbol{u}$ can be added to the eliminated pressure and as such\n", "be neglected. Compute $\\boldsymbol{H} = -\\boldsymbol{u} \\times (\\nabla \\times \\boldsymbol{u})$ by first evaluating\n", "the velocity and the curl in real space. The curl is obtained by projection of $\\nabla \\times \\boldsymbol{u}$\n", "to the no-boundary-condition space `W_TF`, followed by a backward transform to real space.\n", "The velocity is simply transformed backwards with padding." ] }, { "cell_type": "code", "execution_count": 5, "id": "1edae12b", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:16.012529Z", "iopub.status.busy": "2026-06-03T08:50:16.012446Z", "iopub.status.idle": "2026-06-03T08:50:16.020374Z", "shell.execute_reply": "2026-06-03T08:50:16.020162Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "# Get a mask for setting Nyquist frequency to zero\n", "mask = W_DF.get_mask_nyquist()\n", "Curl = Project(curl(u_), W_WF) # Instance used to compute curl\n", "\n", "def compute_convection(u, H):\n", " up = u.backward(padding_factor=1.5).v\n", " curl = Curl().backward(padding_factor=1.5)\n", " H[0] = W_DFp.forward(-curl*up[1])\n", " H[1] = W_DFp.forward(curl*up[0])\n", " H.mask_nyquist(mask)\n", " return H" ] }, { "cell_type": "markdown", "id": "4833c72d", "metadata": { "editable": true }, "source": [ "Note that the convection has a homogeneous Dirichlet boundary condition in the\n", "non-periodic direction." ] }, { "cell_type": "markdown", "id": "4bf47890", "metadata": { "editable": true }, "source": [ "### Streamwise velocity\n", "\n", "The streamwise velocity is computed using Eq. ([26](#eq:div3)) and ([27](#eq:vx)).\n", "The first part is done fastest by projecting $f=\\frac{\\partial u}{\\partial x}$\n", "to the same Dirichlet space `W_DF` used by $v$. This is most efficiently\n", "done by creating a class to do it" ] }, { "cell_type": "code", "execution_count": 6, "id": "ef83143c", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:16.021588Z", "iopub.status.busy": "2026-06-03T08:50:16.021517Z", "iopub.status.idle": "2026-06-03T08:50:16.030209Z", "shell.execute_reply": "2026-06-03T08:50:16.030005Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "f = dudx = Project(Dx(u_[0], 0, 1), W_DF)" ] }, { "cell_type": "markdown", "id": "15d1dd66", "metadata": { "editable": true }, "source": [ "Since $f$ now is in the same space as the streamwise velocity, Eq. ([26](#eq:div3))\n", "simplifies to" ] }, { "cell_type": "markdown", "id": "0cc5e0c5", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\imath l \\left(\\phi_k^D, \\phi_m^D \\right)_w \\hat{v}_{kl}\n", " = -\\left( \\phi_k^D, \\phi_m^D \\right)_w \\hat{f}_{kl}, \\label{eq:fdiv} \\tag{30} \n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "6f940ea9", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " \\hat{v}_{kl} = \\frac{\\imath \\hat{f}_{kl}}{ l}.\n", "\\label{_auto10} \\tag{31}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "93c20b46", "metadata": { "editable": true }, "source": [ "We thus compute $\\hat{v}_{kl}$ for all $k$ and $l>0$ as" ] }, { "cell_type": "code", "execution_count": 7, "id": "0a919817", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:16.031541Z", "iopub.status.busy": "2026-06-03T08:50:16.031466Z", "iopub.status.idle": "2026-06-03T08:50:16.033858Z", "shell.execute_reply": "2026-06-03T08:50:16.033642Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "K = W_BF.local_wavenumbers(scaled=True)\n", "K[1][0, 0] = 1 # to avoid division by zero. This component is computed later anyway.\n", "u_[1] = 1j*dudx()/K[1]" ] }, { "cell_type": "markdown", "id": "df413965", "metadata": { "editable": true }, "source": [ "which leaves only $\\hat{v}_{k0}$. For this we use ([27](#eq:vx)) and get" ] }, { "cell_type": "code", "execution_count": 8, "id": "115a2254", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:16.035010Z", "iopub.status.busy": "2026-06-03T08:50:16.034943Z", "iopub.status.idle": "2026-06-03T08:50:16.039814Z", "shell.execute_reply": "2026-06-03T08:50:16.039604Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "v00 = Function(VD)\n", "v0 = TestFunction(VD)\n", "h1 = Function(VD) # convection equal to H_[1, :, 0]\n", "pdes1d = {\n", " 'v0': PDE(v0,\n", " v00,\n", " lambda f: nu*div(grad(f)),\n", " -Expr(h1),\n", " dt=dt,\n", " solver=chebyshev.la.Helmholtz if family == 'Chebyshev' else la.Solver,\n", " latex=r\"\\frac{\\partial v}{\\partial t} = \\nu \\frac{\\partial^2 v}{\\partial x^2} - N_y \"),\n", "}\n", "pdes1d['v0'].assemble()" ] }, { "cell_type": "markdown", "id": "6dafd2a3", "metadata": { "editable": true }, "source": [ "A function that computes `v` for an integration stage `rk` is" ] }, { "cell_type": "code", "execution_count": 9, "id": "d52d2548", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:16.040968Z", "iopub.status.busy": "2026-06-03T08:50:16.040895Z", "iopub.status.idle": "2026-06-03T08:50:16.042681Z", "shell.execute_reply": "2026-06-03T08:50:16.042472Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "def compute_v(rk):\n", " v00[:] = u_[1, :, 0].real\n", " h1[:] = H_[1, :, 0].real\n", " u_[1] = 1j*dudx()/K[1]\n", " pdes1d['v0'].compute_rhs(rk)\n", " u_[1, :, 0] = pdes1d['v0'].solve_step(rk)" ] }, { "cell_type": "markdown", "id": "c2709468", "metadata": { "editable": true }, "source": [ "### Temperature\n", "\n", "The temperature equation ([2](#eq:T)) is implemented using a Helmholtz solver.\n", "The main difficulty with the temperature is the non-homogeneous boundary\n", "condition that requires special attention. A non-zero Dirichlet boundary\n", "condition is implemented by adding two basis functions to the\n", "basis of the function space" ] }, { "cell_type": "markdown", "id": "54e16ff6", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\phi^D_{N-2} = 0.5(1+x), \n", "\\label{_auto11} \\tag{32}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "cf5b8b0b", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " \\phi^D_{N-1} = 0.5(1-x),\n", "\\label{_auto12} \\tag{33}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "f036c207", "metadata": { "editable": true }, "source": [ "with the approximation now becoming" ] }, { "cell_type": "markdown", "id": "e4e54760", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " T_N(x, y, t) = \\sum_{k=0}^{N-1} \\sum_{l \\in \\boldsymbol{l}} \\hat{T}_{kl} \\phi^D_k(x)\\exp(\\imath l y), \n", "\\label{_auto13} \\tag{34}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "7711f3a7", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " = \\sum_{k=0}^{N-3} \\sum_{l \\in \\boldsymbol{l}} \\hat{T}_{kl} \\phi^D_k(x)\\exp(\\imath l y) + \\sum_{k=N-2}^{N-1} \\sum_{l \\in \\boldsymbol{l}} \\hat{T}_{kl} \\phi^D_k(x)\\exp(\\imath l y).\n", "\\label{_auto14} \\tag{35}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "9c0495e0", "metadata": { "editable": true }, "source": [ "The boundary condition requires" ] }, { "cell_type": "markdown", "id": "de20a9ca", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", "T_N(1, y, t) = \\sum_{k=N-2}^{N-1} \\sum_{l \\in \\boldsymbol{l}} \\hat{T}_{kl} \\phi^D_k(1)\\exp(\\imath l y), \n", "\\label{_auto15} \\tag{36}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "7a45e061", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " = \\sum_{l \\in \\boldsymbol{l}} \\hat{T}_{N-2, l} \\exp(\\imath l y), \\label{eq:TN0} \\tag{37}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "d478e8e2", "metadata": { "editable": true }, "source": [ "and" ] }, { "cell_type": "markdown", "id": "46da00c4", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", "T_N(-1, y, t) = \\sum_{k=N-2}^{N-1} \\sum_{l \\in \\boldsymbol{l}} \\hat{T}_{kl} \\phi^D_k(-1)\\exp(\\imath l y), \n", "\\label{_auto16} \\tag{38}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "08f254e3", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " = \\sum_{l \\in \\boldsymbol{l}} \\hat{T}_{N-1, l} \\exp(\\imath l y). \\label{eq:TN1} \\tag{39}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "6ba8d915", "metadata": { "editable": true }, "source": [ "We find $\\hat{T}_{N-2, l}$ and $\\hat{T}_{N-1, l}$ using orthogonality. Multiply ([37](#eq:TN0)) and\n", "([39](#eq:TN1)) by $\\exp(-\\imath m y)$ and integrate over the domain $[0, 2\\pi]$. We get" ] }, { "cell_type": "markdown", "id": "37475efa", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation}\n", " \\hat{T}_{N-2, l} = \\int_{0}^{2\\pi} T_N(1, y, t) \\exp(-\\imath l y) dy, \n", "\\label{_auto17} \\tag{40}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "11a34479", "metadata": { "editable": true }, "source": [ "\n", "\n", "\n", "$$\n", "\\begin{equation} \n", " \\hat{T}_{N-1, l} = \\int_{0}^{2\\pi} T_N(-1, y, t) \\exp(-\\imath l y) dy.\n", "\\label{_auto18} \\tag{41}\n", "\\end{equation}\n", "$$" ] }, { "cell_type": "markdown", "id": "802c82db", "metadata": { "editable": true }, "source": [ "Using this approach it is easy to see that any inhomogeneous function $T_N(\\pm 1, y, t)$\n", "of $y$ and $t$ can be used for the boundary condition, and not just a constant.\n", "However, we will not get into this here.\n", "And luckily for us all this complexity with boundary conditions will be\n", "taken care of under the hood by shenfun.\n", "\n", "A time stepper for the temperature equation is implemented as" ] }, { "cell_type": "code", "execution_count": 10, "id": "174a1ac5", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:16.044343Z", "iopub.status.busy": "2026-06-03T08:50:16.044279Z", "iopub.status.idle": "2026-06-03T08:50:16.074400Z", "shell.execute_reply": "2026-06-03T08:50:16.074152Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "uT_ = Function(BD)\n", "q = TestFunction(W_TF)\n", "pdes['T'] = PDE(q,\n", " T_,\n", " lambda f: kappa*div(grad(f)),\n", " -div(uT_),\n", " dt=dt,\n", " solver=chebyshev.la.Helmholtz if family == 'Chebyshev' else la.SolverGeneric1ND,\n", " latex=r\"\\frac{\\partial T}{\\partial t} = \\kappa \\nabla^2 T - \\nabla \\cdot \\vec{u}T\")\n", "pdes['T'].assemble()" ] }, { "cell_type": "markdown", "id": "d1745f18", "metadata": { "editable": true }, "source": [ "The `uT_` term is computed with dealiasing as" ] }, { "cell_type": "code", "execution_count": 11, "id": "62dc1037", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:16.075695Z", "iopub.status.busy": "2026-06-03T08:50:16.075621Z", "iopub.status.idle": "2026-06-03T08:50:16.077301Z", "shell.execute_reply": "2026-06-03T08:50:16.077109Z" }, "tags": [ "thebe-init" ] }, "outputs": [], "source": [ "def compute_uT(u_, T_, uT_):\n", " up = u_.backward(padding_factor=1.5)\n", " Tp = T_.backward(padding_factor=1.5)\n", " uT_ = BDp.forward(up*Tp, uT_)\n", " return uT_" ] }, { "cell_type": "markdown", "id": "7db1c217", "metadata": { "editable": true }, "source": [ "Finally all that is left is to initialize the solution and\n", "integrate it forward in time." ] }, { "cell_type": "code", "execution_count": 12, "id": "1a2cfb7b", "metadata": { "editable": true, "execution": { "iopub.execute_input": "2026-06-03T08:50:16.078416Z", "iopub.status.busy": "2026-06-03T08:50:16.078348Z", "iopub.status.idle": "2026-06-03T08:50:24.208041Z", "shell.execute_reply": "2026-06-03T08:50:24.207743Z" }, "tags": [ "thebe-init" ] }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAjgAAAGdCAYAAAAfTAk2AAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjguNCwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8fJSN1AAAACXBIWXMAAA9hAAAPYQGoP6dpAAEAAElEQVR4nOz9e7B1WVUeDj9rrX3OeWmU10uH7qbsALHkoqDBRujGwsovYiuKl1QiWMZWU4BFqVFCrDIdtD6l6gtlKlGUCBHLFKUQbBOk1BLRxit+NBiwwZjPIElpdYvdInzwtojve87ea31/rDXmGnPMMeZlXfY5zTlP1a5z9rrOvda8PPMZl1l1XdfhAhe4wAUucIELXOBTCPVpF+ACF7jABS5wgQtcYGlcEJwLXOACF7jABS7wKYcLgnOBC1zgAhe4wAU+5XBBcC5wgQtc4AIXuMCnHC4IzgUucIELXOACF/iUwwXBucAFLnCBC1zgAp9yuCA4F7jABS5wgQtc4FMOFwTnAhe4wAUucIELfMphc9oFOA20bYu//Mu/xKd/+qejqqrTLs4FLnCBC1zgAhfIQNd1+Ju/+Rs85jGPQV3HNZpzSXD+8i//EjfffPNpF+MCF7jABS5wgQtMwP3334/P+ZzPiR5zLgnOp3/6pwMA/tGN/wKbg0v9xg17FBVjhZsm76Lbnf+9a4ft23HbrvWPbfvvHR0z/O12/f7uRFyToToYy1U1ooybhV4rL7sAlTHYHilz6txcBL93YZSUTysLfzfWMdY7qiLvrrPeh9jOyy/fh/XbeBmp/F65WbnMMsrZ1FC/JYLfwb7Hyq4dk8LadWUKtPKnyinr1FwUPZcl+pNIXwKU9Sep9y9/m/bslq4XuXUyp38svaaFnL5JPW7O+470RUC6P8rph3Z1i9/5/73BjeMxnEuCQ2apzaM+E5vmCN1m6JiHDpq+d6LDdsfRdbZ+B16xDt3ta9vxf05stltGeLZue3eyBaqB6FREdPrv1QEbZOil84GH7Q8GmxLQ7zgcN1EZHA6gdlpUZm+bPHdzoN6WV/bSDsj77YWIlS+HEGhlMMlB7B3ldCzymbete0/ud1SMKMt9w2/rdjtlINj4ZR/K48pM5eXlbDLr2U6QHVbnccjKx+oVr0veOzLqz1qYU7eA/Pqfut8ig3KkjmX9ztJ+RelLCN5zUfoTXn+941kdJkxtj9bxGoL3yEGvNDUpVJ6Dee2MflIi9Ry8Y8RzWGT8YGNgd7IN3qsc06zf6PVFm43XB1V1f40c95JzSXAI7XWH2B1cYgSnQrup3P/uuE38Zdec6LTdsK1z/1fbFhUnOm2LanvQE5627QeJXU96qqND4GCD7mSLCsNApFXQJQcfDjkQAWO5hrI7DOX0jmMNjBoiL3+skyjpwGOdUe51tPIBfhlzrhUjnsH7YfuC95NDcI4ODVVweEdtO9afzQbYbr065MiyQm542b1y8/JulN8D6EqnVDV53aHfPpS/qutR0aQyDX+18q+F1CBXSjJk+7XKbw5CGibMsM3fZQ1kqXuk+hatHyEM9TfoU0R/woevnDrszpNEHQjbI6FgIHflJWgKpewTCaJvJOT0QUE5MutgiuBF+yagfPzg71zpj+gZxNqzJDZuG/VDTQ3U+crW+SY4Rxvg0mYgNoNq01TohvfaDWSnq+NMsWr7E6rtQGhaYAeg2nU9+WlrR3iI7KBugU3dkx7q3DnROdigAxzJAYzB0xqAgGmzLV4jaCBthoYtKrArwwBvWfoE2ZkyUM2a2dJz2erkxZsNWp1UpExZqkejkFC3v2Dg5O+VOllWd7zfsNn4dUj5bSq54R0KJzVDOaXqqeKQ/c6hnI7kU53fYCjz+H5cHaHyDX9zic5cxcW7VoF5MQBrA7H6X6o0yHOS0N6RvLY2mKV+p/XuZT/Csd3G+xMaFMH6EzYwAvH2md0e2T6H2ICeQdhcGdhvIHh9Izsn1ge57TP7ymxFVvZPU9U6eseyP6LyQO+PTHJzeDj2Q3XtFOocnGuCs71ug+Zog67uyUxXV2ibkdB07OlYJKdqiYb0z70iBWfX/99uK1Qt0Bpkpz7e+kTn+BjAJiA5AGxGS/uUwVKa1Sz0Aw/rXLe7vmIBY8V1gxF9oWMjZEc05GJVp4TU5Aw82jGZA5FVNlP5oPtpykcJUQh+w3hstW3Hd6WQZABeHbI6Ua/8ktwMxKbb1EBdByZcq45xE64z39ZUl+qw7ArRIZKmDXJrIaduZZkzGKnmZdfOzzJpSuTUGa2+N5FBPkeZG5Dz3j3IPoVuzQdE+g4kiW509g/Y7RHQJxsxyMM4qeFl1/ZD6e/4fRcmOxa5MVXZYf9i4wefvET6I5XkcHJzdOj3Q4cH6DY1uoAt2jjXBGd3qQaOakdqus1IZNrhXXX0V3nPVQtwIbXaAfVuOL/tUG0rVAddQHZwbQe0Heq6Qgvq/HuTFNoNgKGjPz4eb2aRG1ExaRAi5FZQN/D0hXeDqKu4ZG6gwUgjO8yvQnZOvBEHjmWFg1XKKa7ketpABKQVj6AsWgciOw9OahL+XjFUTA3pDmtGFoZ7YRgwrh0zcrpRiaj7DRq5OTwMiE3HCI5n0tVwOJB5AB2antQP9SxQMansXNFJkJzkc1rYZyXbvMHrvyg7EHesVE0pS5mRYmQmNqgp94/2K7VBcGSfAvgEV1HzAAR1gGDO/gG/Pg/fAZSbWQnS3HooVBqvP6QNdO7Wu1cJ2bH6yhxTZ1HfBIxtnY6jMuaMIfTO23YYEYcxg/qjoT3z/kiSHN/BeKP2Q93hBm1n1C8F55rgbI9qdI9gyg0jM23DOu7U+3Wqe4f2YCQ+RHiI7DR1h6rt0NUN6pMOqFvUALq2Hm6xRTV6q/UOVVza453ezNk1gWZb3eHQ+bSt+8FuAG3ZQDQUDQBrmHGikyIRuYgRG5XUxEgDkQROYDgJMCTxog7EIDbmO0qYQgkdmUTrdnxfmsmTmztPdIJJv8EkN4cbn9gQqaGJQKJ+tYfMR20g9fQbXPk9ohN2iiUkZw1H3CSxkcRBDMLAQPS1umWRm5gphSNmVrHM1RaZMfoMjXzH+pVKITi8X+GD4divsAFRkFzXDxomV3cPrT0O37XJYOp3ePDIGaHxiY/3njLIjiTC/BqZRMdCdt8EBBNjau+EnDHEGz/qdny3V0/s9gxlwsuJKZWT9UPtYR9FlYtzTXDawwo4qNAOnKJtKqBm6o14r53oO2XAUNtWo/9T2xOergGqXYW67k1e9TFQ1xW6TYfueBjQ2g7AZiQ52+HFChumq5Ak3dU1uksH5uzaK5sy0663HbpD9qPaDp0gN9Qx6WSHs/S+7D1YIzVms4DdeCWKFJvUbJcGH2VGpZUvaRKRg0imaYfIAgDVsT0KMoPWVU8UhndTH2/RbQZNcYu+ngy/2RFlzQeJfoNBbtrDjUds2k2NriH/tLFYHatj5I8GAB2RqF3nBrf6uO2fwdDvVVv4Zc8gOQDUaDAVuaYIBrNeWdciOZ7XtYFsyuef5fTJ72MRGassKUKjDGDmoKbUy5hyVyvHO1Ie6VfUmT/z5YgOjO63pttj6jlokBPBUYlgEz9N6QYQJTvDtaYSHbWsBRMvjdhoY0ju++YTl2rbort0EHW/8CYsUkUms9ThBq3ri2q0df4k5nwTnAMABzqx4WRGM0/R9kqQyR1ta9BXEEF0AKDbAs1JP7sFGjTXdkBdoavrPppk0wDHgz/FFuFsBBgbqyQ3RiSY+vsP2X7mS+SIz0B41E5pKIYbkKhs263f0TNZ0pSaI403W7XJnfFqs20+CCnl0xDIqRnkRiU2Qgnpmvg7q3YUpdf2749UkbZDi41PcqgOcXPnMAMOngcfVEkO5uTmsHbEhnzWANs3rRvqVtV2vc1929evXdO430Bl7xXMrU9yjhUVQDMXWuRmAqHh9/FQQjh4dBgfpBRToQrNKV3eWzueI1OhUQmNJN1sG0dUuZMzbN63cAVP9Cs08/dILvUnsYER5e3Reg4WxnK6H+mTHurj+ElS6Xa/BX7/mFJ0RF+ZFSZeOPEyTc+JaGLXDwF9e97CTVyCSZfmfiEmLA4NK5vri2q0hw1axPtmjnNNcHYHg4LTjCRG87mxCA7tI5JD//dkZrBcMaLToMMOlavPnOQAGMxVvSNytWUd3fFxWEml+eCwmTVgAqMpgQZOoPLIDhdoyJHdVV6AzbyHL7wRA0kSkZyJp8hNzgDEITsZpXyqz5AlhVMZEuSGKyEAAjWkM2ZL1bZDV/d+XEQU6m2LdgPXsXgkxymBx77ML6/L7d1U3sNN3/ELctMeVp6fWizCsGo7oBkHtqrt3G8Aaq/s7eHGJzlUdsN+b5KEicQmaeLU6lb0XkP5yPlUmFlc+TWzlPSRkGXhUPxGckiNSWiUAU3rQ6J9okxBoJByoHLEnCZRppKnvWutffL3UaKkZqJD7Zuo2O14f+i20T8W0THUPsA2nQNpkuOOk6ZOIK9vAtT+ScOOq7PAQFLHSVdPZBHvj/jv4+rNMIFvDzdoD/vydU2FrmB5pXNNcEjBkaTGIzcTTfpEdojo1OgJFak5VP2J5FS7um/4hxs027Z/ue1YGWSDDWbYCfOBXU6aacObXTuyQ/3rFj7RIQVnOzZq1ym11HC3WSQiC7nkJie3jAt/p2P9mbYkAqbPBC+DNygpHQgjoFNTEnSHlVNEgOHZ0zurW2f2IZKDTdO/D6niyGcj1BvUrMyHNXZHjVNtdgdjOem9a2XuyzlGF7boyQ4RnRZ92Ztru5Gg1fVInqns1CHSgKb4dHm/hZdhaqSV5WMTC/PnaNnAxQkaEJYfBrnRiE2KzIjjS0mNRbj7c/13HIsqlUEubvC3BsNt17cTKH5ZCR8OwFBvAJXctEPqghIzTM3NrdQHQhAe1h+6Y3m/CNi+izlqTobJykrgF0yMMyZeuZOu/jePE682MnHR+qPAv1FMtGh8Q11hd9SgPaywwwXByQI3RUlyU0Js6NhqN55ftX4HQSSnRYUaHdoGqDbADlVPcg4q1FsKT69HFWcwU3lggxCvmDQQ9deID5gUzu7gGvHQ2Q2+OJzo1NuuZ9EbViYakNp+cE2SHLqLNVBpmEJuZCNnHYe3L2JSk/44VpncvTescxUdCO88eedB6QkIbareNczPa4PhvVXj+xpITtdSHdqFsyb63WDPU5OEN0Nnx8iNJDZmeb1ydsCW6lWHGpUre7thHWLb16sKrC6Rmc3KQ7IksQHyVZsY8aCqzNuAVKE0SNVBXF9VGwTBMs1PEVITG8xkVGkSYrZf73xirg6GZGodiI7zyzpufaKLrXBnCZUw97812dDMMOx5aBhVp+FQj8QM12P9IQDP2b8jP53N4KNDSo1m0geK1RwNZkAKK3eM3BTlgqP+AGLShZ03cdH6I6+vdddjE62hH3Lk5qBCe6Hg5IGTGku1Sakg3AdHEp0Yyana3lzVoEPbAtUwmNR12/viyBsp6o0bOBm5yR4w+QAEiEo62sw50fHUnKHRerPuFMkB/ErNEZOiB2SRm5Skz7drnQy7V0pl8mz9rryNN8jwwcV1IKzz4CRB60Q8yZssfcMcmRMFR3LIDm74c43lVGZ5gjg7SXgoJ0+lINMoaKCn2ZdzGCC2lSM5Fbre7DUMdESeSdYeHNiC60ZNbRIZ/hUqcsiNIDbmQKaU05u1yjrkyqAQG4PQABk+NYpSYw1kVp3sCkaMatv7JAJjP0kKnhwMpakVWwyJIreqn5+ZvkGoN1KtiD2PGHi2+rEf7Nw1bKKDYf/wnS6Sq+aU9JXe9WD3TVzRM8iNVR9UWJOuwwbAbpy4GP2Re5d17dqaI6N1f51+jB76oAuCkw9ObmRnnWfi6f9KopMiObuDntx0u8qVoWv6yla3tZNWezMVIwfKICTJTWzA5OgdooeBR1RSjej4ao5o1Izk9NdoRpIDeI6CgamKEPFpmExurNwWNADRecKZMTnT1rZJ9UbpQNqDyus8ihJLbkgar5zZR5Kc/j69GjKqOOx3MjOVHFwlcaZZU1cP6g0jN7H8UBKDa/r4GwaSU22qsdzbftCpj7uB4NPMV8+lEZgpcnPUpKApc7RdEBtNUaFtoV9avPzuXtKsIn5LLIQ7xwQl1cP+f92nyjJBptQc6ke6IYiPzJVEeKp2OMb1M8qsfzMOioGTqqXiEJh6Y5mJ5fNIgSZ6AALzfb9t+M0xolOq5qT6ytjEi8NQlmPkJjaOlEy63ERpY/RHWl9b126MazfjJKtt+nu3BU37XBOcueRGHl9EcoZtZKpqabLKZwRtjYrGo8ggxCsl+Uj0x2YUfEs+OJ2r0PWgMHVuAAJIIaCZNpEcZ49mSk7vUMZJmZCSYw03hVxyo8yuAVYu2r9lDoq8w7Cc4GRZFDncdSBA0IG48jByYxEbOZB4uh4jCprJJ4jK4/5cEg2b1XH1pmblFOSG6nNrDBD1rhvrfjPUeRr4NkB1PPxWpdz9MfVoolAiqoBx5hdNEzDR6Tgwd8bIjbynHMhk+Vnd92bZEXKTVGqAScSGD2IaoclR6iR2TMkG4CZP3eEwodpWAdEJZv0DyVHVPMUfK3iWpN6AKanDM9EUrDTYMxHme7RdT85J0Wk7T+GuXH3oD3f941ySEwM9qyasK9rziJEbOY54hIcps/2zCfsjVVWm/oieLvnfkJJM5R36od1B5UxUF5mMZ6KU3ExCDWA3DBa7sXKRmWo8rg5nJHIQ4stMpFQAhp7EDF8E0Rln3GQ3BzySg2HGzcxVvMxjJASPcNhOarhJnwqD3MRm2ADvZKD7S0TuZSoGQv4FYHYg3mw5Qwp2BHRH5CAkOV0Lp4Y0wqc4SDsgB2Zm86ZBwJEwQW7MRJjU5zeVIzl8lyPPdV/uhiIQm8rNeh25qWtYZirC4uQm5qSemUMlGMho4N3JFxK573B906fGHaM4ySoRMBqxkcphbhZ3i9QCPbF15WSTvt0w4atReURHUyE9kmP5b8TA3plUUmPmOQtj/0eE3fZT7L+P9w5JTjv2j3NJjpx4aWZn97+su74PUmq5IjMdBG1XTNCuP6orAN1IWOn1RRTz8Z3Vvom8QeiTGsG5Jjh7ITLynkOn0QJAOzoc125/X/EqYrsE3tlxB1BRKXMdVt0giZHoVNt+myM5biAlcxV8nwk3cxklZXD1BtBNVYrDaJZfhaXeePsyZtiAIhcjaGReuXLt3vAHIz4AWR2INbBYCHxbtiNZ6OqeUNXHO10WtsruqQbC92YDndxol/QITeUNeO451L3/Df3fO5T236WZyi9nxFnXM0saZqYSBL4MOrmRyoo+kO1YOSL1iKs3lgIBmKSm/64TG80MlbM0TUBmEo+ypbLR62l80sOJDtVdUnMaQDW1ZvlvSPDJID2noV7z59L/bxMcp2RDJzrVTqiPQs3hJqtikmOVSZKcFLTJF6AqyxwW8VVvoZiga1ToaqbiDG28fwZKReKBDtw8xfqgvhwXPjhFSA0o2n6ZxTh6fj2aqbzEgFzF2QDdSc90PZCcJ2Q7NwgJUwevhDkOoH0xQhbeH9N5JAfbzjVsTnL831rrpiogaqoCEkpNzDSlkRvDGdMznUGU0UncxkCqdSZCAtZmSJYUniI31vpnnCTzgaIbxk83ayIzFb+A8MOhbUGnV7OBsB4HBUlurE6v3vXHUH7aeqjXZKrqNujr2Y6IzkjMfFJYjwTBu8FYb5KRT6WQ5IaVhbaZJiNJcoCRQPPr83bAw8LZ9ZPh3d62MsXGIzwWqeFktcBEBWAgNsO5ddW/eEZ2eqLT5wQjNYfqb2BqzfXfECHGXL2xHPxjv49PAvvrs50FUadksiomOYbiDRQ42jPXBmCsqyn1RuubAHss5O4N0s8u6I94X05qnDaZHywUXt2tEQyRMVwQnIkgH5up5wIAdn1n0pBk33TOD6fb1IBUYr0ZSe29fF4htUFSZlzmEV9creFqDjDuA+B8cnySNvrjBKYqQkrF4YOVBstEEZNg6Xca33k5x8gvY30jy6FVub93T24qUDoQQs7Cru7akiSDCGgfjcdVHLRdb6Yy7N6e/xDgZk1UXufU14iOLUFuStA2cGaq/tojyU/64cRMUjlZgGPwSJPh8GuZjZiD6Uj0d/6gpUFE4LnrziQ29H+M2KRIzRS1W87l+m2VI7t0DA2KlDKDSLqLCiRFz6rHQNSMyScaQKhomeexfaHiXTlFx1JzPJPVFJIDmKaqKCLPot8f+mZZz8JSlX0la/hf6T57qwJTcdphXONqHFc1A9/JyhvX2qaSc+ooLgjORKTITanCo27fjPZK6VwsTQip6Bb53WVfbmySg23l+eTw3DmeihP8nkIVB8hruDHTFBAMQtYsmzoYQIn8MhxaoxASsIS91AebRWX6O6h+LURAhxkPhV8DMIly8J0/R+EMTeWiJU36cunlpXrliHENf6TLhCM3VpklYmHdpRDkRm7n6gp95+C+F14+F6tuKTlK3N+I4zDgE5v+nNAcpeUtChzFFdK6RMJToK8T4wzfUPWGtAHU13S1H4UTQJipHIx0BzyfU24m7r7snTsWkH1kmuRwM34RyaHfoqSwcGWL+S6K4AcgrKccXK3lRDg8Tv/uBRMwFccptM34Hvm4FpQZY73n4ft8tYESXBAcDNFOrAHLqKesa0Q6cWtf22BcWLcB1EXqpL2SNXheKWV0S39NvZLyQdJbWkKQHLq2R2zWUnEsFKg3/HjNOdP9hmEbj24IHEINM5XpoCuuLWeNFqwZpPXuyHnX3a/xBwigf0fVcRd2KOSwbvSJXBbm/lxycLPIDd9m1XkqL1dtyCeHOxqrIEk75lyZm204B8rAYP3vDcIsima8v3ggXmoCdj1LuckgNv33tGoTIzZT84FJ8Jk938b9Dbma06KPJuWmjerYn/knzVSKSTJQbwxyIyOFnKM4U23oPElyAHjme43kAIibq4Cwzkh/nBIVJwXhG6ghpixT/1QzE7Pbh5HkUI438rHz+iS+dIMriz/RcnWY6lINoEA8mNjyzzdK1JsY8QEQvAGKpApmLZZ5imZnzTjLpo8Fvl+GyVvhgWHUz6DixGYFwlnSYbPR/9eUBXlcTL2JkBvPJu9+Q371T60mbg52Gchy7BSdS6v49XAfLLePdSix8tOxrTeAj9fTFqSV5Zfg5E0rrwXfFCNO4gQmlrOGzp1CbjaN6lRshdlaKotUesK0ArqJVUtI127YYqcUkUeKBEsRQe20PfSV3d1h2EfQ+6TZcd+m2aceP/3vzvsEx4rrtWywonrRCVM7qZFu5j/4b9jvzK8LvK1zMzEgfDo2Ibnh27n6FfjL8UnAxs9nxu/rRSwpdSPou2jCxn+vkQ6jOlDWw+OwHIzpdyaaR4zc0P+xsYa/y1bplwEEPmjcVE7X4PctIdsXCs4AqeIUnctIzBS/nLapUJ9o0SaDbwjr/ALzFJuxeYhVAqcaVUG+Eq7k9PeEM1X1ZarAo6rC+1beDCXcz2Yf0r5M+yUywh8Jqg+DAPfDCcxpQFEYYgweWTAc+Mb9wzlW+DX/rvhTkSrCzVS9066RGRsI09nT9TZ+nYqZO10HKHzSPF+hWilzPSpOszAn2WMCKb8bLZpJ5kNx16EwZw3GIJSzKGtOdFRMtZETHO/6SlOM9ZGVM5WybUINp9xgzrRaj0qOS08hZv5Uj90zsvxwhv2SfI6E0Ha2Vn+Pe3/DdxFlCtgmK0D3yYk7HsM3l2eaqoBEgIaApi7Lvklr/2rf1Pr7NLNj/wx7Z2OXvoJHxblrhuQrphbnYOJp5xclBEaqNzE1JxnJJSIptA7MdV6pt6qwcb2j8xu/RqScisNnKVb5+UATIy0WMtSbvtwKuaGZtifNZswImX9EqhOxrpfKs6Ei9g5Zh+ORDUV2HmeOiedvIFYvA8fDGX4asxBL9mgQ4Bj8aLiwLknzUegfM3737qcRLeM9cLMUAFO1AXRy45SZhGrTX3v4qygxtJ8+0eemHOfXT/9esj/pmlAdcYPrpgrfXSRnEVe7++cizPoZ1V+m3UipOf594NS2mJLDy9pvZ3Vv8CWSv00l8ez8vcIYb+hdysjeoE1Yqj2fPHhKcHnxzi0CAmKYlmj7oqapXIgZiWqeqgW7xtiJaZ/+WvzYkOR4DF5THDZ5JgctZFuVXgG9s48NxJHBSzMhjPevgg5GLacVZkzHCSe+4B4JmI1V3LZtlGOlqpLbCWimHTqP550QnTtdVyPCVjliCBZ7nYtIskdpnpTbtGP4tVTTFKQpzSA5w3meSiOi1+Q9fGXIXvxwJCrVaFLh/njKxIfPiKXpiLa7MimkJmWa0s4NTFeS5NQ+WfdyLw110jS1Cgdt94zZOxgVCHqG432S/SSVRyFGOSSnfxY2ybHqlkkCpqY+WBJaEYxJF4GTVSA+0ciZhOTiDDyt00UuEZmr3MwlPDKEt99WBeG76mAokENyxvtqg1zEv4cSakn/A4KmIliNNiMyJrQpC9WGlcubbQ/71TDyiSaNudBs2XL2It9toOLwDrgxBgWOurbfFS9X5BLa4LYqmtr/EBRy44GIRsw3R+xXndU1wgzf56L/nlAMY1EuNJkR5IabpR0RZX4rrSAtmklKkgz6y8lrStXRYKk/fL87ThJyY3LlgZMcM/jAf5bSuTiVmZdDIzp0DTP7s0FyANsnJ1QH+aTPUF1jKs7U6MEIspLraSRXTJS4f+moVNnl9d7TBLZy7gkOMM+Hhs7VzislNqZsKkJ4vRlFrQ9AqVmWRnK88zNUHF4m1TFaHpfDzN16JAq5ISgzNsuEACAwnVmmtKh0monS2Uc0GaOxz3p3WmccdTJWOpaYv43bZgxa0bKuidRyChahkYRH+gBE/LmssG1Jcug8z+wQAfnZhdun+9toJinVMdhQXVwZmvjHOzZyXbe/geiThIqTclBNJXUMyEZo3rd+F4farlBGcgLHYwiVT3FIz1ZxMvuqkqAKwOibGFm2JtOaikNqXFAWag/K++PHl2QvVop7gaVMSo7wxMLGS3LkMFndy+wobZvwZ2jm9dh+OVBqKo6m2Iyy+EJmKiBPeo04FwN6A1bNCGC+E1Jtsjz8rf2nhMB/gc+a6/A9Ahnkiw+USv1aC2SuqoaQ09p0YK9t3y3mY6Q6ByfMUpr5ClDIjeJHIZ2A1Zm69TuMeqWtXwbk+dv0x43kpv8dbAAXpDRGbHL9b+R13Lac/qEWJIKZkvjA6PmTadD8pURfp5GbnN/CB3PNZFVMciyVT5uspVQc9ttPAx7ZUYnN+J2WIgIQlpn9Ti9AY0Y/dDZ661OCRUKcz00h6UmRm2oXJzdVLHqHJGthnpIdWS5SJAdQOqtE1AFhkpkK0ElOZLaiESfpADqexwaljN8QlG1fkO8FadLqvTPtt3ky73AhwwekBDlqj0S1g4vQK0LMQVfWCbZNkhX6aJApBKylEvj/QYSTMlN3Ic6JQSg1w7bITX/fNLmh/72/ihrnvk8cWFL+WZYJi3xxAIXAS3OrMfPX3pUVLDGFtFkmK/7d840aJoL9/jjJiZqqgNm+OGpU60rQAiB89UtXp5b0vwHOOcEBfDPSVBVHNU/twk8pPPWGtinOn+M+f8aVkpJjJCel4oT3zTBTWbKrJDnUeC0HZMW5WA4OlnIjt2lkLNuUthByiAEdFxxrnBslvDMcFZf2r5lGdhRHaWU5BUlsODjZ0T7uuIhfjSQ37hxlEAuuZ20bCHhs9Xk50Via3FiKTY5zMb+Gdm15vdj2wNnd6ltMp3lJfvR75TpNu+smSI6aLydBcuj3maaqUhWnsI1r7XDqUkSQdcnwWQuc+lfAuSc4SyPbaTm29JKSCZgiXACM0m0N15nlDpJAHsmhfZaKw81UyftJYgOEM/LYCtC5Ic1aY2GDEX3kPvUaamTU/puLHpHQ/9XCbAG/IzeTaxF4KKoFpY7lDFrymND5vjP3ZcOY5dJAwUlrVuSUOLa/ZqgIqmtAsYExID2aH84Ekpzrc+PKNIHcePdLEBl+TLBduYel4jgzlVK+4HjFBKmFiPfH2n43KVjmOrpOaI7yiY53vwTJ0fxx9P5sfh9EJuBq55uCq7ZDtU1HOMYmT35QhObakFPAaf426qUWu9LDHFLFyWGvnnPyBEfl3Nmr66TreOWSHUrujEQjOXrlzK94muzq7jfVmTdhnuqPsX0kxrJl3jOnbBnHmEkRDeSYG0MnPt0x1UNGWafk7EkNbkUgfxxtuQ8J5qjI64NU43JmiRbZsbIUAz658a7FSE4rnYYNnxu6b7sJ103yQsFXIjdT/Wfk8VPOkfBMrkoelX67ciIj6x4ZlUREIVnavXkZLRJY6peTJDma4mepOHPXXIsgf2wyxhj2v/TDyQlIWQrnmuBU7QwZzried+3IxyFBqMYcLcz3hvnfaDPr2Ozan2X6+wDdedU/zs92yZ2Nc81UDjIleU4enAzzVH9NndwQAkdQTsbkzFCkEt8XSgcLwIhqsEjLArPBmMLk3UrOFrf8/wk5cWQUBt1bMzWS39WUD0Jyw1dh9pUFfwDzy2ssNaCtQq8gWDZFzpRFHyD/T5Eb714TiIoF1dSVSeD9KM74sjD6vX2ipIEvbTP6/+jL3ZSqOfTd85UySI4rT8p/EdDbLU9MWoA5wTVWHZH+g945Wn6oTbNa/3quCQ5BOhWXvHRJTKYQJt75W509dzDm5ilC6Yw5SXKEihPNT5F1v0giKw4eJm7Yl01VKOYnRGv48PTk0QFlWtNwykMbyr/1blqCu2QUi0oyLH8F+0K5A0guqSmFZpoFFMVPdSw3UgUwkpL6ENTtCcLsddxyAMvww7Hg+d6wgZKu768LFG/XXhkikyALsVDqKdfWzFRLQPosevsi6o0GSXRK1BxpsoqRHMtfiychDVQcIN/XZnAw1pTRKROMrHfsEW2ZD2eagtP3oQXHT7rLpyAsYpLzMKUqM9lZmQ9+fB0brt5s7M5L60CiZqqMztByOOZyuaXimGaqmIqjQZutG8fH/CTcb4wuDlfYJCZEJuTYuTVog0qKdGYl+5PnMJVwDfi+N8pzsJ5N4LfVBL9LyzbskRShzni3lYkg+TnQ/W4ATmhCkuOVTfPDAcy6bOedCgcSNcxatnfD/2W8n3o72+k4RbwnQIYaawOjmQuHg78L4zmW5FaJEZ0SNSeH5GimqvD3MVV5iqlKmYAtCa0+AsPvTvVJbadPdiaMq+ea4GiJ+OZIdpwk1Tv/IxG1cSqVTWbkJPNUv48dZ5AZuZ9/l9for4+kw7E7LhNJFcdKxsavkaMCZdh3oyQn5mgMZHckliLh9kfqQGywsY4ZzRGJ3x9ZrkEivlJw1iWyIB0ekxAEIXDMlL5Y2jPJNlHV2eTG98UYk9SZ+XAU0GK67ruyACIwzTRVQm7Wzn8TO7ZtKnNgBDCaNoB+gK9DH6pwshYpW218lHLlEh1X7mFCSN81kuOur0wOkyrORMj+KaYwSxNzDP6E2DATxnypBPrUEtMI2Ckk+jh7oBVv56wozq8VIzReJ9X2L46foxGsljrKxvZEn+LcR/ei//lv5ytBt6hcBesaBCvFumtuKmDbAaj8FbmHhurJo3UNdeXcWLSP+jt8HyWCNSDx39Y1FVr0q/32K/0qN5jpp1LtuijhcKsrG0TCmjnnmEL7+47vqNvUqJQFmAFkkcL0/QrlY62d5MjlPKEfIzfjRcJZcKn/hkQnBrX+/qFyMx5ThYPF4IdT1TXM6SgvcxNee4w4yjdNTfGDWVKZya2va0CqKeN2Ir+JC/D97JW51bOpX+T96fB76d604jg2Q33YVmiHdtmvng7XbxYRfVJxEiuNE6pti+6wHibQffm1/qnaAth0qFG5n5waX6i+8PfcNkDdjv83bLzpyz+SHK1fqrYdcBS/bwrnWsHhUNePmtMorf5rSqKzWulI+axBzNYIqVw43uxDUQNkyF/oMDbOSIIORAl/HM9TVJjYrETuyzBPAbZKE/UzSISgZikfSy8kqSAYgGIteU8RCzlQ822UqKbWOlO0rZb+NL76Yn0syP0WudFS1zvTreKHUwrpWAwg2zQly57avnSCv7nIWfnbQ+DLwv6fM+Ipqg5XdHLUHCqPzEydMvEXqTiFOXBmr5OoEWhDWaRnEaSvWCkFx7kmOJYJaUoG42qnqDctiuyGdK43i2VyNwAXPQXEZeWSGVpMzuamKiI55HAci6gKf5yexKr/bRGSY0SZxMxTQYKvJWp5RiSVFdq82MryGVhrkMm+/0o9ivrupWlKmKUksYleP0F8vMHLyK1CsMy27aYKwtGDJH/uftb/odoXU3RT6s1S5GYtTA1syFLsVBOUTlaD8wqJDs9uDPgkh2CSnBRiixYLH0HKZiyV0kmRjAqsiVfn/fZ5ZL9EIDjXBEfD1MU2PchtGtEReXPkejzeoZvas+MCYeObMrDkkBx3r1ITWKSRakmsPJLDP3I/j6YyzFN0f/57+msoM/EmJERLYFIHYgw2qbDdIj+oFX7rlEGxyLmY7pObP0mQG2B8hlM//f2JwIfZavt7+GajMV3/dBOZvQCvOG6ievNwwRg5qufWImhJ/tx3p6goDuYGoYkSnhoB2YkRnfF3+CQnmWMspuLIdqCRnK0xqBmRnnMh66K1lIy2JpWMQg36icLJ4sO82i+H0vCz9PW0Djw8xo++Cm32XgbjobPUzFP0PQbNATlFcsJojWkqTs56K1bHZZIh5frmb9+Eswf/t41ErGjJhgLkRCrEsn0WOW46p8ZhgxWRM3XgjSgEJSidOcr6E6g3UrnRCEri49/P3y5JjDtOcTIOILcXOIuWOPgvrt5UxmdlZNV3RRH2r6GrbISoUhM53iQ77jif6AB5JGe2imOAq8v8/3rbmhmNl4bnXB2JpJJKeO2sI+Uq014Izmte8xo8/vGPx6VLl3DLLbfgHe94h3nst3/7t6OqquDzBV/wBe6Y17/+9eoxV69enVbABdag0pyL610XvhTDbGVWKJoBitmHZePk0Gah2jnSL0e7vsxynNP5pPI7jL9RDFqWKUqoNxq8VWjZwObfv6B3zhmE2AzJW9BOITRzZ0ix9zj5mhGSk9upaJF8q8FMY8/IK6ujkqBoHw6N7EjVRstg612Dl4VNUqyFNz3lwbXB0BQVTT63hnoTayorkZxcwuHVAW8l6oyCLVA/VbIjVB0iOq4fZiYrL19OjOQAZSqOWWC7b9Im9zmKjkWKQ/8j5SBev5V+dgnBYXWCc9ddd+GlL30pXv7yl+Pee+/Fs5/9bDz3uc/Ffffdpx7/4z/+43jggQfc5/7778dnfdZn4Ru/8Ru94x71qEd5xz3wwAO4dOlSWeEsslFqphJ+N3JQiA0ScnXletv2la+mRGPDbJSFF2odgDRnZBEQQ83hJCfwxxH3zFZxlMgWgHVSSji4Rm7UZG7ydxXktwighS/OcIBLhYpPgTVDL8nrcdaQ/ZwMEyUNapRlWPOXiZkBLOIT5GExTFKU1M13aPVDu6eaqTS1wMozkoss9eYMVyef3ITqRWeZRRJpD3LNlRJRssPuy9tu6GYQpgLor5OZG6cwG7BmRtdU5ryIzcwxRyrL0CdYVDaVYJ2lRH8/+qM/ihe+8IV40YtehCc/+cl41atehZtvvhmvfe1r1eMvX76MG2+80X3e85734GMf+xj+xb/4F95xVVV5x914443zCtr6ZiotR04J5LlSzdHMU4FcrzgYB7lpCt9g0uM9QnKA0R8nZaoar6f4uAiHY4vkjMfbs11AN09p6k0qEVsUBZ2HlFg1CTiF0oyrHLbpIvGcBeaG9vr1Zg7prMcPgJhpKlBuFOVlzgcIyZL0k9GIjvNH2MR9SEqgqTe55qkkcl9XxnGT6tESymSCEMis0PnXjZMelejAV3OAsO/k5VZVHMi6n1BxlD4rCITIVJmnmIhcfRN10ouk8tpR/CVMKcOqBOf4+Bjvfe97cfvtt3vbb7/9drzzne/MusbP/MzP4DnPeQ4e+9jHets/8YlP4LGPfSw+53M+B8973vNw7733mte4du0aHnroIe/DMTWJkH8NBOqNRpJCdSc+4PUdtmwAZZ1VEB7ehOfHSI47p6gTSNuTTZIjP+x47X9CdOkF5iMxXiMsn3XtFLIWh1wQOe9ChmQGWChhGGGNCJwc/yv9vJDYpFSamMqj7deJjNgvzGb8d2hLTCyNPFOy2HAWlBvWnJKJKxniDsjDPys8bo3seERHqDkWyTFV8BIVJ5HIsxKRVZYfjopEN2eZ0EcFS4+k4u+Nykdlqtpu8mRrVYLzkY98BLvdDjfccIO3/YYbbsCDDz6YPP+BBx7Ar/3ar+FFL3qRt/1JT3oSXv/61+OXf/mX8aY3vQmXLl3Cl37pl+KDH/ygep1XvvKVuHz5svvcfPPNAAwfmULotsv4d7qnzHzMjxsdXkU0RsI8JbfHBh25X/P6l9tcGTJUnGhuHN7xJ5x6vbBaLaGb95tClaakg/SQE6FgIcPGfVZRtbpUPQVz/D+kUsfrQEy9keSmP0cxZSpmJTrWIj2WKcoK66YyuA691gcpOejIa2n+N8A6pHItOHV8xtp/Mei5kTLOM+poLI9YtF8VZCcgOsgjOWM54ioOAFvFMdRadU0qOW4ZfqFJvxzreTqVkSwSVbj0BsGKoirE6iYqoDcncXRdF2zT8PrXvx6f8RmfgW/4hm/wtt966634lm/5FnzRF30Rnv3sZ+MXfuEX8IQnPAGvfvWr1evceeeduHLlivvcf//9wTH1rgvMVMB0mT5YhDNCeryKxF7onPDlUoXH+x4hOfQ9RXLcsUwlCVBXAXGJ5QkxiZBhnjJ/75wZYYzkTFiXKqsMigw+Gysl1kpF+0zJa2ISXO++I7nh9+Lkg5MTy6yU+rjjMxQfuu64jf8/8fkbk5xPJcwOVVb98vRDrUlhTv+ZQ3xUoqOQHNoviXKOijNZAWy7ZDqLqQRD/i75f04fTH55kmiViBKrLtVw/fXXo2maQK358Ic/HKg6El3X4b/8l/+CO+64A4eHh9Fj67rGl3zJl5gKztHREY6OZuZ8jsCagcilH2hJCP/czpwx05o0xOidzNewWYJhTrLQ1Xp5ZaptOo7SjvN70Pe26dtqiwqNKH9X9+nIZSpyf2mEcUHRatuiq2tvJhslN4o6lFJv2mbsPLtaSaWPXjmrj7u+LFPD6wbU2w7t4XS9X5N7rbrWNnBp0eOFCjvDfZvXpsKrA0K98Y/zFUS+DlAJeP2InSsJh1bHurpf2qTd1KiP5z/vecpYZOeU6loBWEbwm4zSQT6mhM8uC+/zZd/Zjv0m0PedNXxFv9vALedAS+K4ZWWOd8Mxdb8czoYIQDsseQNzouWWapDQFt5kfWm9A3Zzn00NgP++k+GZ0CSXLd1TbVtADPn1rsi32Lvtajg8PMQtt9yCu+++29t+991341nPelb03N/93d/F//k//wcvfOELk/fpug7ve9/7cNNNNxWVj3xkVDNT6SyCrSvlhY/vwnByul+QQZkzU+YzUpJqPDZT8RyHFZuxPF+apty5ij9OjsMxQTNVeeaHmvnm1CKk3LuOFp7L/9ejX4ohfVUiKk7MiS8350tpR6s5m0pkhc4uhFJ/LfUaAWlRyA2gmqZS5CbHkZjOsVQaMx8K9Ov0ZaEyk7lKPKSZs+UpDsanZeJa22RrKQhzrjclykpzA/DKFFHBXV0uVXFi6qBBfqQfTomKljJJuVs3Gf0xJzkoz5Mlsfpimy972ctwxx134OlPfzpuu+02vO51r8N9992Hl7zkJQB689GHPvQh/OzP/qx33s/8zM/gmc98Jp7ylKcE1/zhH/5h3Hrrrfi8z/s8PPTQQ/iJn/gJvO9978NP/uRPTi5n1QI1OrRT/TUE6p3/EgM1h+zQ28Espi3PwORoGUGViywnQza7oHtJJSf4fyhXi1HFoZmIXIyzq30VxzWkuvIWuaRGSpXb8sHwGrc2g5cr9AqJlM+wqVxJQjvHH0dBtQU6RZhMm9fKBoduUwEnM6fXLYAFZrxdzetVBex8cyzNTv1zIg8kME3KTn8ayZ1kRlPadj9TH5Sbpl9s0EtwZi18Grn2aYH3CSUI/G0M0/1ii3FO7L9T0aVF12Lnyf6SFtd1So6rK+n+s2thqzikCtLCxXwxY4Fq2w79HtAd+j+8akcxrmo7dG5RzvB3FTmw79j3k+H/IQ1KDePlMz+cfnyuipWk1QnOC17wAnz0ox/FK17xCjzwwAN4ylOegre+9a0uKuqBBx4IcuJcuXIFb37zm/HjP/7j6jU//vGP4zu+4zvw4IMP4vLly3ja056G3/u938MznvGM+QVuBzWtzn+JFpIkJ+LM2W5qkcOj326pJ0tAEhjNNBWQoVonOdghkFo1U1X/W6veXFVXQNt5RMfyt0ktrFlqjvAkYKYqVRPMCVXb/2JX5G2LXeO3TKcsGP4pclZul3uhwaHt+nIuKOp6dYVJ1CY44R3qgSO7jNxajsWA7gfmTRIKYD3XZMZw8U66DVAdU53sZkXrrQqjyXTi+Z3WyuAqyCyjQM8pk3fZJftUIOxXPZIDOPuLNwlg/WfTsn3UTo75pLju6xlXkHP8vNjK4oRq26EbTOvWRKwU3IQ+jmfjRNfqawPf1cIVB1YnOADwnd/5nfjO7/xOdd/rX//6YNvly5fxyU9+0rzej/3Yj+HHfuzHFilbxZnlQtfzEvcpJCc4hxiq4gthJfgrUWZKYJEc/t1j8hb5Hmav2HRpktN2PskBPKLjIMmNUG+032uZECY5Mm42wHaFHOZAkTpXquKUoNp2wEE40tW7nvRPKU8nyA35owA9MamObYUpOyxcidzjbadkhi5n3bnQ2o5TbgbFipsaciHbXNbxK/iT8G0pkpNTP1PRpqUoddzWTHolzzl1LP89nOhIv0zuy1jvut5UhbH/dORYqjjkv7KF7y9IKg7/HkG97dBmMIHKUHL575O/27+RPMdoA23r+WLWO6CbSLTO2BRiv4g10pLG1q9Ibq8rpa1Y3pulxu2VJ9dXXv6bjjkYcwQe+0t1anXY+K3OoCR03J2n5Z+RId8UTss/8MnN+LtD52IgJDeTwsVjnWYs18QeHHdjHaweFj2cUJD/Jrcd5EZ6WYpVAO1dJdQbdw/h/6Vl+Y4/u2kfOlf+VhmyvsbiriksrbpM6WukuUpuD453YeRlJta1nm9u/dGOl9uAse8E7P6zP4fVZUaQZU4xL6cStfGCts4n2HPCs2PPh7tdADQG6MuX1NvOGxenlOlcExwOkr4oXLzkPHOfeCEa0QHmO1KVoCSvQ4zkaE7HsUYKhGGQ/bX89Vb65Sl0IiLJjaXeqIO7kptkEgrToc+BNljuC5x8lyDW8RdnMw58reyQbJNMy4FhKEfOZyq0Z+BMkszR2FIczkrOpMVU7YxcYYQlVrPmmJwDy52fT2hyrsO/A35fCsD1n+Nxo7nVhY0ba1S5cyZMZGKYEyrOxwdvX0E3yl056l0ZWd9fb32GIc1Ic8ErRNV2QSOj7MXjR5in6jGpn7NXFjYytZPN8B0AhFmqHmVVbjt2x2Z0xr28qpuqqsHx0pmrBj8cM+pHIUjuf0O9kWaLaIMd/CS4Tbj3BVnP8SA1mFrkM1tdKRysLXOUlNfN+ylytQxjJ7+U8Xg9hUCg5CTUG0LLSPs4Sy57DlNJTs1NUcDoYEl162G8ZthcWOrN0qROqxf0Pn0Fgf0fmFAy7pPRHmR/Coz+nbxNteABEGwbhvEg4YvT91MA2tbzXytBtRvHq5hZegraBkBboR6CHtpmEG3I1HYM5xxdDa4L3C9wij/QuVZwtBXA50Jfx2NaZlgtzNqaCSSvVUDgchL/yfJMyXLMwx/lrMQyURBy1ZtU/pJJM7yFVByLVJvy7p5VnCWQ6xcVHCMj4Yw09Zp6o/mqObJST/gUQhIjk2xv7BQIUzCLJIhHu0Rdk+kysjEp4Un/HJdKhxA1YUZU79zjVZOm0X8CKFNxcteas6KsIs9/iWWNOqPftfztqu24qnnpOHquCQ5HqfQFxG3JWpprT7Vh4eEy/w2PoAKmzyTnwGyQ3J9B+jYsRHICoiMIj1thPaHepLLMaueamCL3FjZG+UzPCuZ2atYsOZYJGBgHKz5o5ag3dE+V3Mh7NOEnwAQyJNusb2qA7WS80BIZZxlWaLgfnHE2n8Nc4meRHN7uZR3kGbmzfHGU/7XvAJhSomc05qahpUH9v1x00wKVpcg/dmrhPpXgm5QwrFGVDklLeY5ba3nkghyMgXkmtDmNMpU0zDKfWDMRDovkAPCJDh3Ptrl8ImygS6kxJXbf/nqiefDvpOIU5GkZr8tNZv1fU81ZWLWRndxesxgbSpu16ClgR8y5/Rpx9vy+fHKTJDPKMUVtT5AcWfe9axV06vFjpu3bB/zAi/QxEmsMrClMVVBTDuj8OoFfo+wjDRXHJZ+MqDh8jSq59M1pwf22WidxAOtvhZmqT447kJvCMfXcE5ylFhX0rrn1/1c/lv8N4CKoJFYxUVTsY8BSFaSjnEVkCKYjqEJyJNHhxIaTG+0edE2vrAly4zmh7jHr71lDKsJlqtnBcpg2yQNX7UTEHJHaGKHtap3cxI6POZOmsh77ZdfPd/faVMLEqsysZy54Oidqam4/kxUibvnenBIhSymmsWeS44RsER15b00FB3xzekrFcShYfLO/8PJjoRc1Jn8/EbZmLL+qMu2mq0jnmuDkdiBTZkFziZN0MF4Fsk1EyI4VLr6EP05/7khyJNHhH34vIK3eaAsrliDMxTP9ZXCTilfegkuu6ocjnOO9Oixm4qk2weuBjGByxxjvxVNqWMSc5pyrqTcyKqW/ryhbbIad2O+V1SI87r6CyAt1itevJdS0fSo2uXVRLlWjQfVdzB3QBvNxp5iul4D1O3OITeqcnBw8li9Of35ExTlFTOmntDKT6YxHGpeMreea4HBUbTddDjXy31gvgu6lhYe3m1qt6Iv7ZKT6AIXoxEiOlpNnKslx1zSeg0VupHqjDaCxhRfPSnSLfL6LoMCZelYOjIjzpUwZ4M7h5JYpd+5Y4VDO15yyoK17ExuQUo6juaRHIzmyfdCs1UJ1Rv1P5sAzVQn1xvW77vvp/X5LZUkdN/s+oq7EVJzoGlUSSn6ZteHGBbMfqJyjcUsTEqV/Iidot6r4BFXzIkwcQ/jZ8CSWTojFw8S1l1Ptun7mHPgkjCuIq/bKNWfyXgEBbaXgrvbDx2WZgvBHQF2vygsf3w3Eg97FtrMHpIQZaY5qE6CusZR2boWo0vfTci7uV/Bd7+ZdMw5oXc3qTF15a1JxeMSmqdRnk1JvAjkf6baj7df6BVkeLyReOb5tgGpYSbkvXzWGyE7016OsxbFlZVxIcsaxS0HP2K7/r57Pj52Yk2kuYrnBSs6RkOv8yRQc6nXpvsMSOJTdOFhpnNcnvnTDwiSHEzG3LeO3y3QRHuo+N5SXlZk5Fk+hvedawSF/GBUTxjMeylZ7sqzPPCkPTnBLQ4ZfHKW3YMdbkVWaqSqFmJIDDERPWzFcMfWkIqeKEfjwcLNJWbOJZVZtDX8i+ZyXHJBK5etgpm1dN2e269UT3Q8rMEcycqOtGC7RNpXas6VCe3MVnhyThSOyyvXaDCKbSv5Z4udyWsgpY4l6U2p+W6wvyLlXom7JY93/mho+V8VJmKpi7V9bLij5exQFv+QY6ZMWoO0jjSlcvMTR+FwTHI65EU+5qIU92lp/alEFYgkkSA79HwzMEVMVByc5vBH3+yrvM16/wDTlXa/spzvIUHEiOoUh5N3GNtdknb+ieifNI3Nmzlp9AEaJGgifhXe+JDqZYeEETb0JyqeYnHJ8b3JyobilGVid931w6sCswMN1rc489U7W8MOZGsUZc1LPqVvBRDBloliR1MTIb/G1IiSHH2P1r1FfHHkvFkXVCed2SXZSy1xMSVliuRrISQr5pC3pP3RBcMDj60UuAEsujFRo3iFp609p4J1aIPlZsrMShbD6jM0gOfK7Jl2m/HGAkeQAPmnxyIn4HlNu5LHafonsdWwSKg5f/bq/l3b/9G1WMyVwUqYk/JKDa6lfRFA/pLKBcNDUHM3pXO5vFUvsJ9UbbfAoUcRy/G6snFHR68p11yKoM5eRWcvBeG6Wd+mUbkVOUR3LVQw51OibPeQUWioFh1TDeR3SVJz+HDYBjKg47h5W35ZBCpfM9K8tuqlZLqptHypeb/uMxhQufuFkXABtvagc5FZsWn8qMFnJBH8MS1Sm1aIpFJKjmaq0hsqRS3KsLLC0X26PORHr5cg6bERMrRlmRFZ2WpVsKc9rb/5VAlrEgsSceiVzYQDjM9HzYkiikz971HydoiQlQ7VJHZdLcvisW0OumcBKmKftW6I/SLWV5Arjmnpj5caJrfE3DHylqGeEG6fet3VOyjk9aTI1VBwi+F7faag4Wn8kt1mZwufCmui471xNX+H+55rgmFFOM5UQy8eG31PrcLqmghb2nLTXJ1Sc3N8j/Qxy7MpWuDj/XzNVefsNkiOJjhaSq5GbaHkLGlGWVCpCVNXraFltM0iMFWnxcIJlpgLCuuUlfvSys1bKMUzNcXVLV2/MclmDToZqw4/ztim/iauWQb11dd6/0FRCMpfIrE2u7TBxu1/chxKzJErrTWyiSJij4vTHjX9LlgdZy4+Jm6n7+wz/iPtVbZ/sT/rh5OJh2GWug1DJmd+oQpOXMkNWFpYczTFsm9LxLGmSis4k5L5InbdCx2OOc3SeJDn9eXkmJk5uShtlPGGc0kQyfW5yEwbujbhM7KwqMdueUu/CDps5TW70d8CJTios3L82u0aCgOcgRXgskmOqRWzJBs0Ph1CSxTgn14w8LgeLmiYGFDkWD2bSfZOclBq3RPh4jOTQ/6q/Ilc/LRWHEWhJctx1Jvq6TK0Tmpmaf2+Ncs7JD3VBcMCin1ZwNJb+PW6Nj13nSdGU/0b6luRWJn3mk3dulpoQITnaDITv02bIcmZL4A7GPEkf96eRpGcKuSlxNDY7ggLn4m4TDtAxh73TUmt4nZySd4KTgJQfjvYbpYN5f5yMtkq/Z43QWIpLqXIZU31S6Irb9YxsxjMIqYWkmcroh4IyKOQmNKuJCeIS+XEWMNUtQW60c2S/KY/j6jdNCPt9iooDf4JVSmbU6NVC8pYzqeBlB1iZBcmRfji5OPcER8pdKUfjKRJuQG5SWWBnRlAVpdAv+D0lJEf7X/UvYSRHdv6eycp4Jha5yU6rr0EbPGXEASM3audBRDVY8E53tF2a0JQOjHxlYa0TCZZwKBgo1E7ck+N10qI5k3No5qlYniHv2pnm19ywcL6N31tTLL0BySC5lTfwy0iidObfkjwzS4PW8JNr+ZWsmTXXwTgXpeRvDrmx6lLMFK31m3xSyNuHpuJw3xq+LpUMgoiXp8Csn0n0+X1SPmkA0tFzBs41wak1pzeUNS4L1bZTnTVTHUvOQGz53CwpTVvINVdZpira19VyoI+THGAkOvzj9glyoyH2bLNWFEe6M9VCL02zhjIYaorXWUK2E37EN2VUckRuD5EqgOBnup6u3mhlyUWM7EybuYsoKv5dMWUTWdDMOTEfvJLkevtAvYMZNQWk+ym56vVZQNTfRqk3sTpUquIA8NqEVHGAtAOxp/QoGcA9CD8h3meXKJhyCRW/nxcmKpp8MT+cXJzhrvThCSnrV9paGsI8RZCDLDdhxDoubXs0emJmH2GRnBxTlRx8LKdjIAwVN8uzoM9NKcaZkJK/QZpW6LOCX8OasEhNzoBphWfnKFeWMznts5yL5X1ztsvyxZAiOZaJLjwuXQ9jhLJkIjZlgrOG/w1BJWpqNJi9rE14USXU+BQdlOcoiKl6aC3CafriACORZv45o0konmwvJ2Js6rkehNM9+d9M9cM59wSHE5KpIeOpp2iaptj9/IUkywbB08hYWlLhVX8FOUP3JPyQ5PTbQ98MK3R8LmIOwjmSuFsdl0XGefsjK6KfGpT6r+XCyfbtsuzuiplKqjiaD07KdCszBEv1RpsRa87DKadiulYqBDZaTuZo7M436lXMN7BWJjJLhoZrmNLGcidkQNonslLyNk1BjrqVeqfJOp6ApmxaE0OpfHvqiaLiuDGlVtapksQn+F12m1u6zwragkH8nR9OAWc9S93r3qGuDeW89u3z5sjdKUe51KzOapSaqSqWRXQJeL+9Urax75r5JTDHREhOjOhox0xBkROdZ8fWlRtvwchhULP8b1T/JLpXhr8IYfagpkX6YV7dcZ10xEylnsfeM3c2bxPPIxUBE9tuHVs6kKWUSiCMHuEgh0oTiaCCmOl63xMiuh83TwEJcqP0zVOWEeDXTG1bCqVjglVfo3WU9ZWWigOwujWQHE50tElcV0OxJNDfZZRvre3H2sJFFNUpIDYgAeSoyRoxe0f1tvVnyzNMJjn266m+OTkOujGSk5o58+9BNlp2fYvoWNv4efLepc7G3poukfwR1tovuf438vicaJ61ECXhc/OsGGaqQMXh/lVG8kbLPEX7vO8ZpKek7LnXTiqdMnrEQL2z/XDk8i8cc5zDl4B2v3rXGdutayyQsoNfe+IzyFHslmqvKRWHl4GrOG47U3Fa4VDsER02GdPyda0N73d6y/AofjgsH04uLggOBodg7jezsN12nKnMVG8M6TkmT0+dreUSiRSsxqh957Nc53SsEIEY8VJNB0uShECtMZpQwpR2WihJ8kVQlc7d+DdX7fS2CzOVeoxIEwDAVG+keYrfO0vRaZBNKudI9JTgzFOnxMDCHWmjfRGfNE1UcNdUc0oWBJWZ3gmaCW5OeZbIbzavEPrmUhVHKoLcV02qOBrJ8a4V870ROdmk6pyC9Vtiinxu/rAcnGuCIyOd1pQt3T2kejPAdcaRxQc5SmdmS83cLFIx3kjZJr7zWYiclUiSQ8e483LVpIx9qn05R4YlkkDZQb2VxivT/4Z8q4LfJYlf7HdYfGpG3bUk4MBhfollGvg2SWa5imPkwZHn5tzHOicaGZUgO4EKqflSZA0A/kzagc1U55Ac63sKqT6odLJD5ikv702E3FRb43fTtoX8cHg55iJKjCuM5Ib/Hznf9BkTavdY3yrVF8cdy0mOjNzz7pOfVHPJyaNckyqIpmLrUuXiXBOcfWKRJFXymkYnFnMydI051mdO6NysQQKwWbw14wZgkhy6rrt2gSlHH2QnzhZkVlBhntL8bzi0KKLTdDiWJCeWl6SUTGlkVs4IObToOd08Gea+4fcDIgSloGNOEZ3ca0r/K83R2EKuqqzmxzmF8PCoD5CW+A+nG/FUiqK2anUxCtGx+kyN7FiJ/wBu7vUdjvvM2WE/1dEkw0sO6JcpZVZfEhq5mYILgiMQa2RT2Gq1G2diU53kUnkuNMdjr4NZQZmKkqAcGTZlqgLC0N9C/4YYYsTGWrBOi3KJJfnrB3LfkS7mf5ODVRO0GTMka1XxxdY9qn0VB9Cdsekc068p8/vUejOvvun/e+1IC3MWCrPqhxMho9FFK89KXhz4/S6vb7QwMQBV+Z6L1X2UcuZPkWNy6lxMxXHXYSQH8ImOVK2l0zLdg64zp6we6pGoZS+6OfjhFNzifMM5A89IIKXZ/zmx8Q8Ot8VsoFYYaIzomJFWCzfeYOYcMVXlqDiBqcrdp/KOT80YnJSbWbv1hFaj6YCTGB49xdd46XjmUPidhpb/xjRvLEB6kgpLgR/OkmbbWDSVRXKszMQ5uW9iquJU5J5vEirLpyKjDwi3x0nOPhJ/qtfOTNwXIzeU/2a8Dj8mvyMzQ+yNbNBTMdtUw8b0mIpjORtrZeAqzniMMAPxqNWNkZ5B6ZOteyfdCCJdjzaOckxJD3CuCU4snj6nAygZRPms2MrGWbrOi0V0+P61UWqLj6k42v/SVMWP8UwfBqnR3o9qYlOd79K5cDSnXTczGtaf4rOTtqnUjuI0zVMWYqHipWRZey8xM5VmhtSuY+W+0e5t+QFZH/O3FJjEotdhvhISUkmzBuqA5BhmaTeRSxCf1dCO0VPVzjClLbwWICdQNesrLUfjxfvOhXxlY4RB5hCT+aSAkeT49bTyyE1wT2me4vsyFOjZgR515avmbXr8VC9TeNtPaXAHNwBqZ6HZQ7OxRLij1YEZRIevC8OPXwtLqDgEaarS8jCUevLPgRWB5Kk3KfMUnZMxAyvFkhEiUn30TLcLDpD0HHhHnXLeTak3uX5NuUqgVe45GNtEPHrEV2A7lRgE731F1XZpBH0uRF2DeAaZ6o2lyC/xPLL6kynkRlFxgnsLFUeWS1Nx3H7yx1GIPKk3MgLLwSA1U0xWSV+1iJpZ6otz7gmOaUpy+/XtsZckG1dUvYnM4twhhjkq8L9hNno6LihbzPF4IqKmKgWc8Hj/s0YUzO7Z4DYn4VR01pGp2BDR0fxyZD6JIEQ8s8XlkjIz99HU2bBi445FUiUXjhWSNTdTpchHiqiUOBfPcYzMITkl11OJeqROh6tqh21fIzn7iAqNgcqXym+jkZvexDb639SDf9hqZZ2ZyHIxSKfjDHVSCxmXa7tZamG30TOt03X6v2n/m7mIqZkER24ufHDyUCkOk1M7BWl6INQTXoo7l5UvVi7Z2dG24Dh5Dc1FaGJjNc+LqDiuGInZt+Z0XEp0chunepwaplyH+70F44QpqolEhBlmtbkdSvYaPhFoOUncvkRbySEVnOxIud3tF+9drl2Wuod2z1KUKDmxus5hrbM1HpDOhyNJjrUY576yF6fuQ+YpNbdSYR85J8MtIeWrNAkJp+EphKjE2Vg7VpIcSWw4EfLNW+P5Vnkm/R5lMiLVzDmrxBPONcFZCv6gnDfoTpmNpMiXGkFl2N3X7vBKVBy3LTIDtyKriOhYn2ll180GQZZitpidl/uG/46Nb56yEmUtOgucCOnEF+v4cyN1pkK+v/D7eGwqhHuq0mJdb4n9o2mAHaCsGeQ711rEwP9uTYy0iVA2ar2/yq23VvZi7lQ8HltaOHHNXZoY1ruuyNS6RPsMlEyL7CRM+ylnY0D3xenPT6/nF5R3pQlYX864mqlmNC7ABcFhmDLb1R20+hck1ZtJEqsSHaF9CNL3Zp9YU8XxMLHWFpkRtJV4ARdF5ee8qcBX4yV7dhA9leOcN7NF6tlqJ5DpzFBxt1/UtZiJSJqpNBXHNElahDiDLKb8BXJm13Pfj9VGYtcNTN6J4IL+GHsA3+daVNFIrgi5IfNUP0HLrL+WyrX1HY2XTIq6NPGZepz0YfPaGzNVuW1Kjilp1ur38YirZUhNDDkCQbW7SPSXhd6u207OT3OazsYSnOjEwsRnE5+6M2dzQN7MOrZfly7Da+fW3ORijkrOiOAall1YzLi5eqNFT/X3hPurhU2HZbT3ee8x9U5TdS4SpbCkWlPUQcpjlXqQIsapSA+N1KTIztKdfGyhQYIadZTIicXPnd3mF0RPXpTbxOoWrUOUmsGvtOROdr3PyAGWtb8ythcgWLJhIC48/YJmJpVJSGPh25oP5SKIZFouuswCRUniNa95DR7/+Mfj0qVLuOWWW/COd7zDPPZ3fud3UFVV8Pnf//t/e8e9+c1vxud//ufj6OgIn//5n4+3vOUts8ooZc0SmbQTTl4a6oJkVcG924yPOJfb3RezwfOOrrTTy4yoim0DFJIzsQZra22pkCqOlubcUG/oupp5SsPas6McBBmNc9dFykA0vFqoOJ5JUnnPc9WbJWbNpdfSzivJpK0RAz0jcP83quJMITwFbd4i32HARBdsT9aztl3E/0ZCNeXvUQXPneik1E916RNlwVpJbjzfm4h5am0lJ2irLFS82rZF/dDqXepdd92Fl770pXj5y1+Oe++9F89+9rPx3Oc+F/fdd1/0vA984AN44IEH3OfzPu/z3L577rkHL3jBC3DHHXfg/e9/P+644w48//nPx7vf/e61f04US/pRWCsHm2Bkp2jGnVtXtM7N6PCmqjgq+VG2BWTEGATH4/MGEZnanK/doqk43DQl1Rt3PdZRpNSaXCQdOa0OoNSJM+HLoEXx5cCLnot0ltaCr7m+NVPDVHOOL+3gk+2Are7cOhN3nwunf872BCyVD4sfF/tOmBpokEJYRpvccPPUFIU9xwxFBLBE4ZrkoF7wPGNO66XX8q4hSA7/AL4ZayRMcfPUkuOdtpzNElid4Pzoj/4oXvjCF+JFL3oRnvzkJ+NVr3oVbr75Zrz2ta+NnvfoRz8aN954o/s0zfg0X/WqV+ErvuIrcOedd+JJT3oS7rzzTnz5l385XvWqV5UVbqonfyRfQfQlLWyaoqgJjQTxRrvqLKREyYnIrjEVJ4vkuPIgWaunNsx2U3kfKaHyTKB80VQtuZ/sMEo7zspQ7RZDROqvjIRyJjJFCjkzBRBI6to5qW1LOHMv6gQelC/Mfq2etot/14gOAG/Sc1p5cUJlpKzvnRoiPoab+344q2ChxH7W9cwoxAwVxx2/8cco+V0ilV1YlmspZacV2eKnYlWCc3x8jPe+9724/fbbve2333473vnOd0bPfdrTnoabbroJX/7lX47f/u3f9vbdc889wTW/8iu/0rzmtWvX8NBDD3kfh8K1LRxE5Zsz60l3bJ36sY7REn1x/5vi5H8pEqPsnxzlkjDjWPcIy4QgrDz3ut49+Aq8xv7A94bPjEQ51pwFcUwxJ2myf/EyEMEF/K8lOXFySK917X3gNEyK8r3GfHNyFOBTWYjTc1incvBoMctZ+GxlLVzLfJl13UI1OPDHYcSG+qtAvTHMU2vCW6LH8H0sqQerNtGPfOQj2O12uOGGG7ztN9xwAx588EH1nJtuugmve93r8OY3vxm/+Iu/iCc+8Yn48i//cvze7/2eO+bBBx8suuYrX/lKXL582X1uvvlm9TirsfPBP6Y+qB1e262boMoK/eSztoJBKegwF3YwjOaIKFRxgPXkdBXCB4eTG6ne8DID8c7hLPjeeIhEo6jbJ4w7ydmocF7kqo6sCzHz1BLqTcn5ue8yXHqE0g2EF6h2ne8LJd6DjKTszwnvubSim1xaRi2DHjHl9g+/0zNPFZJ1MmnNzgF1yokSAaiqd6oeWhFVfPLFl3Fw5xjQ1OY1iU4s2V8pVrB6hagq0Zi7LthGeOITn4gnPvGJ7vttt92G+++/H//hP/wHfNmXfdmka95555142cte5r4/9NBDAcmpty12jf/W6l2HXWRENh0dtdl+opFW2w51XaHbAt3hsG0XXtsdH3GGq9G5DrTeAbt63L9axaw7oPV/d9vYM/6u6X9fV4e/hfZZkOdQ41TvVSNqTunqCth0wLZCiw41qr5VbDsAFapdh3ZTB34AMmMxlcObDdV55qmgTMPxUzrYImfN1LVYm6B6NBkVQP5e8v3Ss5C/N1UP5LH8WtllMi9o30crU0lZOdpmnGV29XCNmj0shqrtgKHO0bvlfQ1/91TnqnboD1ABU9q+nOAo7Tx6+k7/X0MOIZmj4lRbAJu+jdNV6NnU8Ot31c6feMzpa816Jvo+7R5Bf8raF++TJbmR6k2pecoq46I4S3lwrr/+ejRNEygrH/7whwMFJoZbb70VH/zgB933G2+8seiaR0dHeNSjHuV9gOl2XY7Q3p9nT+fQ1v3REnRFl2cI7POjqWpvC3BG1J5cZ2M11DfTqTQ1o5Sz5hz1hxalaze19+HliKo3rOyaeWquyWrpjoTCcOutbTLos9Eu48+gRVZZqotUdvhf7ZpR9SbVPAv8hqYi69yEScpOZGcoKBNU3bmwfAQB36mYoCvS0/pq53eTcf4shatQdOC+Zcl+yAo9j9QfNSs4U3LkvWMO+dLPZx9w4+ew6CaMdQBTWLW4h4eHuOWWW3D33Xd72++++24861nPyr7Ovffei5tuusl9v+2224Jr/sZv/EbRNQlLkpxgAF0g1bQFNXJCS2TFEv+542QDXtIKlWPSysjxkBMBM4XkpK7fX1dffZf+ug8jN1K96ctSBbOgWCcxZ8aXs0Jy/13pvbf2aBdddNM7Ll6+/iD/a0z2Dp6Zsi9lvswtR/S4hdTyqCKZ0NEpkspdS1krbMmcL2vA739E+RVyw81TtTNbkenJZiDu2IkLz+6D+Gl9VN5kS/8fiJMQvoBt0JbEpCFHvZmklgpkkUhzCZN8Brq6ieplL3sZ7rjjDjz96U/Hbbfdhte97nW477778JKXvARAbz760Ic+hJ/92Z8F0EdIPe5xj8MXfMEX4Pj4GG94wxvw5je/GW9+85vdNb/3e78XX/ZlX4Yf+ZEfwdd//dfjl37pl/D2t78dv//7vz+pjPW2Q3voS72l9ZxXmraZpgYH5drlBavQbIg3EidNM9LFZVcyD60NTxJNyK7afrevQPYMTGPCTBW7ljzXPa/WV+UC5YaHYVLnIAfwiKKwaIqBqQttuvNbdIdhgaq2Q5cY8XOlff6uNZOV9c5z6+ySvjfBtSeao2JoG6DaVOhaAJav07ZDN/RT9I45OeLEgeqlK2ezgJmRUGimsmCtGj6VnITXT8/dqhbLdNYCVp2LEZmg7ykwU2nHOfNnG56nTRhMx+KIb2QusgNaFsbqBOcFL3gBPvrRj+IVr3gFHnjgATzlKU/BW9/6Vjz2sY8FADzwwANeTpzj42N83/d9Hz70oQ/hEY94BL7gC74Av/qrv4qv/uqvdsc861nPws///M/jB37gB/CDP/iD+NzP/VzcddddeOYznzm7vLwTIWSZV4YOh8iNs6cH12/90DcRjlu1ALa9X0gqCZgVNtqyMte7Dqgr17ERkZD/L4qZHWCKzHiDI/tdHDH/n/B+ox+O59uwGWeXSTl4OEdbgbckiqwEuR3F1EzdQaj4jFXcxwvBG3W8TleQ7xhZMh2OU8Rxyk/QXWJWRbupfDPhru8P5HuotroCVLWDrwnrB1oiPeK5Vru0irQUgkSFzKm4L8v0CCpHzKXT8gaD303n+lbyw9ln5F2OSmP2W0MdlBMDQJkEyO01ej+jXRf0l5r1gdQbK+JxDvYdtVd1XXe2tc0V8NBDD+Hy5cv4R1/ycmwOHwEAaA8H/4qDCu1hjZNH1theAk6uq7C9bpDvDoYXTHWhA+qT8XPwSWDzyQ4Hf9ti83ct6pMO9fGu76icrZkRnCGpV3vY9H4ehxV2B5VTAzQPdzUsVJm5eSF/B0OK7gP2O2pRaaux8vrZgidWD0Fw1OSDnf9dZl+W+/hf7zoI92n3pZkMLfxH58tcGXybuzYbbMjLX1uBd3c4Sr38ObebsJNIERztWfAw//pk/FvvuuEvUB93aE56s0V93KG5tuul/uMW9fG2zwZ6zEaa7Q7dpQOgHtbZqmt0h41aN93fg75e7Q5EvZIdoUUouvjvc89Aeacxk1/02c7hZ6IZaGW06rFcIqXeYYhy7N9Zc9K5d+a9r20XJJNsD/v+YXeg/xgvx4kM/R3eF6+Pro4y9TGr/bdjEIP7rVtRJ0+A5gR9P8jqZtX2v7U/pxNtemh/3DzVdkMG277eSrLTbWp0hxu13qpmZCOpXdAvCpOP/H8stPIeJpqh6JlyaP2ltuagaz9afQNcnQvLRf4uvmlKew655imrz/L6r6EdNCcdqp3eDurjFtXxDvXx1pmmtrtr+M3/73/AlStXnD+t+SzjxTwHoEFsQRu2p7zUfLaVtiFz8MFWDQVVbO+xpGx7dTCcSIxmRR4U1mY95N8aOCr3oeOiC9RFnIunlHVtBHWyZQONdxyRP6hO7CWwUi7E/KxyzHxzlTEVC/ni5KiK3HePT44ArnjofjfVdlRJgsR27fR3tQ945GYmlrjGVMytd1n+OLX+Py/DqM7wa1fBB0BomjJInnY/8uPZpxqWizPWzT48YL3Mlm3vNlV2JBVvjDWbDcXgqwvio2TslJlM85y8liN9pflqTBNEwrcidLCLX9/b5i2OOZAXma249meD3vGiU8nxA1m6U7BUpyxn+sQaP7p6WFQ8FWqUXOO/e8uHaXHH4kzkvDfz2Rjb+7ql7+N9hFy2wSI6WWWJIdb+J/QNcmmGSpjgAPY71d9U9iO09QVjWY2Lk1gWYF4i2OGfjJw4QTvh9Un4BcptmmlKLUdi2xLoagwRq/Ma7rkmOG7Z9alLNgywHLG86xU0Tq0Dk5/+mhlEaDfKknufvRmdoNVgl3IgNa9jbNdW3fWvFy6g6fbxBerq0FGP141SG3bs2LnvsjQFelaOElmmWN+UuWxHbqcarUNWmG2jf1JlXgM5i27yAdtSc087mqp0wKOEfkCoukjzVPQ6BjGq2nTd9SPUMgu+MkpVHCCuaKqLFHNXUIMoxcjOeuRGLoFTT16y4VwTnCWgSehTHfbkys0x8sL3BQSIzVAczmAjtmCZGlJhixbCRRvD3BDuO0tfLnNFeHkjFFs+v9+azsUc2hpDcwc5mQtnqaiWGGLv2iIgseeYFRgwYd9SyF5EV0Bbgyw8JtIPKNeZi5LnpeW9CTA7QeXM+j9DyZnje5OFjJw41oRby7nj92m+Cu1dr1QxPSN4mBX3bMOvZCELjaFyKotI8rcd//IPHVuSwCq8Z/LUKEpZ/JINPTl7yVBxwobsJ+vj8IhNxFFRU2+WikAglAxOS5OTRVWChGgR60xjfk0xB9DcOruPd1VEMizTo+Kbp95rweSMa0BVb2KI5G4CYKryVtK/qYRzTVUvtp6fZaaSEwCNmARJ/hS1ptRUNRVFz31Csr9zT3BmLeDGZPYxAslfaDGa7I9FCQBhh8dNUcGt27CT49lBrdnbUqQmZgf2MNFMlVI/otdUjk2ttcIJDABHYOSH9qmoJ6o3CyWUW8WHYE5iv0zE3n3KQXvNmaX6Die8p7nPigZ7yw8HsAMRTgs57yOZsXgwT+VCmrHqbduv46WYqWLqllm2MwjL38ZqIzGnfe2Yh7N6A1wQnEUhzRQxyIbLF4jjDc/yx3HnGblwNKzZ+eWy+rnOxtq+KZBZp/3r2iG4nPgE6o1SPiuKSkXG4Jk7WM42U6kri9vXzBkQSqItgtlkxNHYOi+2knkuljJlFQ+YRh2UjrPqqYLonLZfTgyewugi96QfzeD8Tn9T6s1UrEBq5qrWqoqTCEtPmXFlW+qasJ9aQqVJ9lV7IJHnmuB4RCHXbmt0mppaQCHFuZ7g0pzgKTLGjM3aD4yzFO5oPN4rq0geJvs1ZEZc5DaqVMiidazrLJTzuYrTb69GQrOB55tD15KmKa7eqKGUGSaUEox5LhLPl+Vgmnsvjlx5OUpqIo7mueQih0Cu6VeTS2Ap70fyeonoS+lvZS6fkfDjWwIxM4qFoC+bsGK4ee1EHVf7yYX4UnEdq7vwM+F+KRWnxFctpd5McSLXtu1LaTzXBIdAjWJq7oQgFA+2ChDA8PqXIZX0f2x2Rp1danG5NeXXSVFCxhpFU81W+WpS+I7ssPJKJTdWebJ8bgrJTc57853Py65vIdYuUj4dk1STCVFUwfYVfSNiyHnmMWLYGZMhKx/OeN8wytK/Z7pcZwFuollonvIgTf/CTJW7+KaFRciyRWaU7VPW1wv6TqnWNP4xUtFZEjLJn4a16ucFwVkQ4cCclwtHdl4EmejP2p+NPSX6Ko1qmnSPRM0tmeHL9xY4FDNiE0ZMheqNvHa2vxJwagOzihVNGxoZsd6ZFUUlj8u61wQE55/iO/L8cLadGY2k+t9ps2lrYFkgB1bsuauKICM343GRDmswVcUX37T3aerWkgPt5H4u8uwt30XvGMPPJuXPxq+nKkErqqBr4dwTnFR+BUKueWY0TynHDA7H1j1l59Ufa8/U6LjgI5Jaucyzp4lUh5lYYXyO8zGHNFNJFYe/R/nh+2OmqZyOpFi9SXS8+0pYVu/0gWEpTCGm8jz+bGPtVgubXQI5Jo8gjX7R9aXDbEytFWaqfUxwMtsiXzE8etyKWYldP2qYUghLOZwDWDSBKkH1tdF8cRR/Ns2xWJ4fYKHACGDZgAWJc09wLGSbmLxzxr9jtlvYDoOaI+duXJtFpmUPyM7MfA8aTBNNqXPwBFnV257RyJIRAUsMlsr1LKdieeySoeEmlu77RTZjb8HHhevbHJNSDrmJwU9fnyY7U3wPYpI8MK9jr4QaKyc54fHpd3eaM3Sp3oRLhyQqesZEVZr/l0T2s8slN+I4tW5mqDhTIhL5/1P8BudOppd8PxcEZwa6BuCLVAKRwS9GmMhmzCBJjtuudGKuMz3ja83wRrtI6nKUkZzg3kLFiZlFuDqjLU7H91tmlFgUxBJQfS8iqe/PKrSID+t77nUIKSKziqLDFhacA82UPaW9L9ZHZAzWc56nVLqTDsR8rT/mh0Ph4hxTlLOSqM016pErR+GksCQiUSVBK/6Wtcerc01wlpQ+5cA2LsRYqblwqPF6ZRCNkZMcbWZmERppppqElW3wKhI5cebmYVDPV0hObKZjrbyrHb9qx5AyWa2g7nEsZQbLNS2Z76XANJU76Kw5ONVMoeUoaadywI71ASnMVVrl/iVypVSeqT7RR4uQ8VyXA469OWCX9qk5Kk4Euap0UX8rJ2dTJ2utEtm7gun7XBOcFFJmCA1dA/Wp8lDxLiMjoyQ5bvsEpWau5Le2dD0p0maiijMSlcT1lZlO7sq7WlnWVm/Wxj6Wa5DIHVQBJJ/rIo7tC7YDTlCtjj2VXkKSHK1f4JOcfWYwDhzxc2D1U5y0LJgD5yznB3KYmSjVUqW9ayjnBfdh95qDfYaIAxcEZxJynBO7Bt5CjACAuvKWbPBmG0JSBXySo3Ve1a5TP2cJwaAw00wVU0xyjlPBCAx9OLxtjCBZKs3sWeyM8PGlOu0pM+GpiHWiVph41DdrASIyJTS3FHPDlQGdeJbK/tnqzB59dNaof5WinnFyOXW5hrMeXZRSPrOCIgqhkemoP1qiHeSIAhouCM4SIBbtfDQkc07MxBQZVpIc+j+HyBSFg541FCzdkOyYLQc6TcURLUElO4Lc8Otp4eCfKurNkv47UyKXYmHiAFZTb6LnTXiXlnmq6BpiEgQsQ3IIa64XlwPLPDXJnUCZNPbXSkRtRW6VQ2bUZ7hQ5FSOs7HlYJzrZDzFsdgCz4FjQeZ6A5bzzbkgOMifLZTOGimSCoAdSWWoOIBPcixCU29b9/Guu5APxhKzk5xrlNwnRVxKEJAceY16/Lg8OKwDmKUcTURu46c6k8zSvTn9KejkzjSh+CztTzOlPUgH40mTjTZ0MJYkJxY+vpR/Q9S3KRbVE3luQRLDmepNJSIBpyLalvc5WUmZqbR9lslpgpPxVGQlJl154n1BcBIo8gMAV3FGcuMUg4hNPQyLHElO7MMht51ZG3Ohmaok0VSMcFhkKCgDIzXyGHN2VKjeaFFapdBmPsCCZqqF84/I0Oz4zTMuuAdyM/caVidvLvJYODEJJjaC6JylPkBNvZEqXyHZya2zQcqNh4vCjQkh45H+xdq3tO+NhqkqTUm/dO4JjntYbdmgq8GrFDTAyUiqjPw6WkbPEpw1P5xiJJL+ATZxkfuC83JJDttukRsrQVYKZ91mvy/IZx48l1hTWWAhzbWRjHSbQz68dfT0fFpu/ykO3iWKwHyTVPwcGX22xHPJItVLmKcKVBxLjYk6ECvHA1iH3CgRVGvh3BMchwmJ/TywQZk75WkZjYNThZkKyCA5JFkr0rUGrXKbg8IamTYjzsbmMTnXyTg3NwGglb1YOzeWICum3pzFgfhMQ2ZMzcygumao91Roif9oQdySwVbrG6auozcFOaZ6r815k4DyftaRnSkRVMLsn4N9mJmDexYquTEVx7uuFvGp3CcnpHwuUhFUa6iNFwQnA1aFl+n7+fHBytRN5UhUp+TFAUI1yevICgmNBdX2OiQrTM6oV0KJmSq1zdpXcqy2n5uhHIFRtp0mSv0srHq4BmIKGYdZ5yLEpsg0pa3gnCD0S5Ml0zyV067lMYLk5BCdZF3NmOCUmos7i+TMnVgOULPC75H0TYVKNgonnikTU07/pl5rQfUmnbeL/+/7oM7BGeiWzyaWCGHuB0a2ZANGPxwZ9hY41kmSMwOuQ1nhbZcM7EWEaYLDYkpyLTVrxbLmap33aag3ltS7RBTCPgmQd9+C51T0TGMDdwbRCVAwAARJzeQaUXxfqXxvOB8vUgdSz3ei2mu1K1/Nnv4DSh2VV0stkLNw5pLlsSYBmek19olql58qYarj+QXBKUROaDJ3NObQMhpLyMzGQE9yYh9+rF+W9Vz9+eCveeVnodBMlaPizCU5VoRByjS1BB5O5qslVY2YKTCGaESPxMQ1gHLvK+GRi0jfnGWeUtZnWmLisyR8xZr9X+JYvifMJn6nme7BCNDImQjK/i3ad0aCIkqRWsjU375svT7fBEewwoCATHAgDW7RwDkay3w4loqjkZypIJITC232kNPBzzD3WChd6HCqBCu/l0QYcEIXbNuTenOm1xqbgRKSY72z2eQmgqUG55R5Kuv9aupuJIQ8B2uSa28gzvBJTF8wb02qOdDM0MExM57ZHH/DJAr7naA9TXTgV5P7Ka8ix8E4O8VFAueb4CjoEunRJbgfThD2OzQOL6NxPfriaDBJTsT3ZkolmKS4IO+c3Cim5DETZ0ql/gF0jhW2LdUqeY8zrbxkEmRnjirMGNrS84ooW3MXVrXeS9Y9ppCbFZzsgXEAsMxTU3JXxUjOvhFrEyXL3qgEJeZgzPfNMGvFypjq99aKnlLre46KI0hONknLzCRugRzpNYf6vsD2uWtE+10QHAVLmHYscwqpRDn+DWpDjzkcGx0cNYKSgUat2CX+NrnHZnYCyeilQlOVtS12jBYmbu07j5FTc+35i2ccXpioTClfKmpkcqfuORfbfYBFmqZOcEqQIqOpDO/AfDWm9PxSJXkKTqs/SE4WCoa9Kb9Bq+tE9tfK13RBcBKY0wm49YrI0XhQh1ILb3KHqmrbuk8MvJPrmipQovhq2f7N9t/gimbhGaGP1vY5jsRWjpsgkuqMkZfSAdMi2lPXfnHPbOIcYTE/jbnkZsL5qZmtvS//Xpb5eqqS09VY1qdEKgC13laSSvmSa1FNCBWfjQnZhy2kzslRcaYg6W9a8Ft0U1X+sVNxQXAMlLw86WAnHY3JD6e/LiM3ERVH8xrPITqyXN1GH5Qlg28bxD3+VzJnmedmPv/SRhgjR6bzXYLMTFFvzorDpYagXg5pDmhQIj+KlNlhbfJnrmu1lHKTGQGT9GtjWNI8JTHXX2FymgjxnHSzyrBPUcdTq6aH18p72EEfWqjMq79DS8a8h7YclCU3Ois3I/iKiTOD3E+W+QrLtAOOC4KD6bPVAInKlOuH4y5XMItpnTpUO9WIPv09qQzpzjiXFOQg69xCM5W8bsrh2IomSN4vou5o183Fmh3ipPQGBeHgpdcPl8HIyzuTus+ZIYiKYmFBLqdhzWDrbTuqDtotJwQhuIlOxGm2FKmIR3OiwPo+M7J0qT45AqsOreUAvC9H7gCx5JiaQr5iOYsyGC+guJ17glOU6yPCSaRMqDkak8mq3dRJM9UkaLMjuqfiOEtlPQ3MNVNFr51BcqzjrH2WaapUvSlajymjXElYRFqrcwX1sCSKTb9/uvPKTQ5Yet0iiOstTay4eaqf2WaQGonIQNBtKrXM+2r3MgppXPak8uvNHL/HjMViLZWI+uVFsZKT+mQVhyAfwb7D3SOh4bF8OHP8sM49wZHoQ7r7/61MxdHzRcfP/XAAYXuuqyTBKk9YVblOTZoSqLORg3Esn8ISjohTr5Hy4E/ltMm6Xq1/tPNOy4n2NKENDHwGPosoZ5Kc1PIZJdfbF0rqyhRZ3lqctwT7nNxkZSvf1Ojq2usT3f8LrnhvqrJGtGTsXZ6lNp1lYqwQVXSWrhMp3zMr2eVSS4+cb4IjZqw5Xv05kIOk6ocjZiyTVZyI3NsvETGWydtXMBvYC1IrjGe8mqkkJ3a9GMkqUW/OUkeYHRKemlWn6lQO9uAvA/iqak7Y+b5hDgQR4qLNbHP8cHLyYE15Lmmzou/D1d8nz1w/BYsp4yh4Hgs6F2ddZ8F+e13TVPw7x5KLRZ9vgsMgZ6pTkuGlfD0oH45npprgiyM7Nn4d7n/jfkeB/40s/1KYSjZKVJyS+051Tl664z8zxGcYDDoljYH065p+D6PjWjniaXY0yAqDFknyXL2pdl3xzDVXvnf9zlRH4hQ4MRILD2uRVDyq9DQQq8dSxQn2n5U2q0B7v2v5FKVgZiuOpk4QG2b64VwQnATMAdTqVMUaSp4fTi3s4WLxzdiMI8dU5Udo+flvLFJwZgZYID3QTVBxgHhHxTsAbVZvhrlOiJwqQmbfX5JADVDqmFLncma+i3f+U0nOTHLDj8s9dk6b0XJ+LN2p50C+v6X7gdJIqtNAcTlWTDw663oZfmJmBGhGvc8yD0fA63dsiYalI6iAPRGc17zmNXj84x+PS5cu4ZZbbsE73vEO89hf/MVfxFd8xVfg7/29v4dHPepRuO222/Drv/7r3jGvf/3rUVVV8Ll69eoi5dUqvt5g4zM7b7bCGjdftqEkRDJGcsYEgiOJ6pqI/w3C3zTbcTSBEhUn1ZhiZS11KC49p+TZ5HQKp04yJ2QvJqg+ClPHr1KSk3CyPK1IF3mNosiRCXAqjkGMUgP5aSgSY5Rn/704VLwEUbWmcv2knAimoEYIzoDlA5iFTJJTYp6dRWoUlSZYcHabiiZ8GERR3XXXXXjpS1+Kl7/85bj33nvx7Gc/G8997nNx3333qcf/3u/9Hr7iK74Cb33rW/He974X/8//8//ga7/2a3Hvvfd6xz3qUY/CAw884H0uXbo0u7zWonEmjLBXKXO6BiTMVBQynlJxurp2+0YfisqZp7qm8swIrvMwHOe8Qei0/W8S8Bpi5mq5U0hOzvVUx+YF1Zs553u+DVyl2/hKoX2+v9+rW8JxHYi3lcm+AjnHJULNVxmwc3x8lMfbNrbKJs1T3u2sEPG2NXNkqeUaBncqR6f0B3OQNP/RoM3ajxw0W9afeXVQ9oUTHI1lnc71s4wpHkvDCqmPRnmuTErXnnQlk1supGIusfRZFD/6oz+KF77whXjRi14EAHjVq16FX//1X8drX/tavPKVrwyOf9WrXuV9/3f/7t/hl37pl/Arv/IreNrTnua2V1WFG2+8cZEyjuRi3FYq/aPu0KJCw9io15GQgrKr0DYdqk2FrgWw7QeSessao+FZzsvKyY3zkWCD0O7Qn5U4BYeRHCt6ag31hl87mqmy7oB27JRjzmhdE7fnavvpt1llyI6+OkVlpqv739U2QJ3rslFXXqfRbWpvUAwIdmLWX9w+SkBkolXusZBJKnWNWL0KUAHo2Llt+L95Kttv+d/Qe/IynLftok60a4DqKaFtALQV6pOuN6G3QLWpgJNuqG+DH09dA5vxd3t1ddOMa08xwpPjOC9VbmvhT28SU8+rUzkEMLU/O7Mv6zuBdP9pIeUzmLqmZ5JSlZz4+Vo7mBoqvmoLOT4+xnvf+17cfvvt3vbbb78d73znO7Ou0bYt/uZv/gaf9Vmf5W3/xCc+gcc+9rH4nM/5HDzvec8LFB6Oa9eu4aGHHvI+ARayC0t1RJu9cIWFVBzPUViZaXeb2iY3wwzbUm+8MiTUmzXJjYViJzhDMSkxgeWEh+dc3yrLvs1OvpoSZs3m4EqgVdfcdXnILtUt6aOkqFuyTJNASg3/RLC6qUWJ9IvVr0nIWFsuF24yoyT509RHd94SKqRoF9GUDjFzvUZYNk1UzVGVSqN/d+4DSh8Qew5LmKdy68gcf7fZiTmXQgY/ieWB6q9RTnJWHc4+8pGPYLfb4YYbbvC233DDDXjwwQezrvEf/+N/xN/+7d/i+c9/vtv2pCc9Ca9//evxy7/8y3jTm96ES5cu4Uu/9EvxwQ9+UL3GK1/5Sly+fNl9br75ZrdPmwUFnui1vc+CZ+usfanaJf0bfHFoEJEkh38AOHOWdxwzTUXVm9r+Hfu2wa+tEAXbJv6+fcrUKUy9Z5A6IEJoQgdkPzLP31fwHvdgAl3dkfNhBtenRNTHyWbE2H0D8sL+H4hE/78/GeNmKne8lg8nE7wu+0lVx3vHzxcbJsx/T6UOKe8wd8zax6SsUrJ5pxyMS3PBcexlvl5VIgS764JtGt70pjfhh37oh3DXXXfh0Y9+tNt+66234lu+5VvwRV/0RXj2s5+NX/iFX8ATnvAEvPrVr1avc+edd+LKlSvuc//995v3lBV/qUrqmamk0gK4hu3IC1/Ogf3PyY00TbWHFXYHunrDZ2zu/wL1pmu64LMkvOds5MTJUXGAZUiOdc2cMkzpKErPsZQltePm9Sh1Xa4oKuAO82p5TiFAZq2BZDIxzuhVvQzGmf434TXmhYiXIMsfJLImlTT76GuIscmdFt1nLQxrmKeoLo/lkalAfAdjrdxZWFG9SR2fG/ySMjudeqCDgbmrya/qg3P99dejaZpArfnwhz8cqDoSd911F174whfiv/23/4bnPOc50WPrusaXfMmXmArO0dERjo6OkuX1swzPm/lXW1856RoAJ2PjrtH/JV+cFrVbg6b3yRkqKSc23nfDNNXY6g2P6uJljf+WmDNnv6/anX7Yp/Sb0GzX9FtTfjvBtlN2daDfIn0aoudsKlTHjCjyOoVxhkuDpKpksjqm+SxYjrUPd/UjCuHnwCHrl1NNM+pPyv+mrIxj1vRgwhYxx6aQ8l+LnVft+ufTtkB9MmxnfjieP+Kx6HPq2pknLN+xFLj/jVlOpT5rUadZWZlXaAOWP06uz9hpkRjyufF8cyJ+pgAWW/191a778PAQt9xyC+6++25v+913341nPetZ5nlvetOb8O3f/u34r//1v+JrvuZrkvfpug7ve9/7cNNNNxWVz/dpYXLoDNOEq0TK7N5FMgyDBVdxnKmKERr+oXL2+2rVNKXZk+eoN7kqzRRFJ9nBFqo4JffQwiWtsEnTJl+q3qxlpkmpNwo86V+LzgOiak8w49XULXed9cxTOaGu3vETBvdsPwzxuEwFMDaljPjf5Ko1mvLGI6jUc1KTHDnoFzxDPbKsUp+rNFPJZRvcNaXpHoBMVOmBJz0F1P4yt+zhtddXbyade4YjY801pyjR5YJ5oFaPonrZy16GO+64A09/+tNx22234XWvex3uu+8+vOQlLwHQm48+9KEP4Wd/9mcB9OTmW7/1W/HjP/7juPXWW53684hHPAKXL18GAPzwD/8wbr31Vnze530eHnroIfzET/wE3ve+9+Enf/InZ5c32gFNuV4zzr598lEBm26MJNh2ACpPyfGgkRuaWReoN/nlnrK2TTdLzSmOXsk4vygKQV4vozPJeq7U2cSigxZE2wCNpfiIaKoAwjzFfRakg/FpIrccJsll26fWD16W3Bl0k7jXErk/+jLpfiYxhWJJkIrNv3etmDDQ5K4Gurrrn+EQTSWVGq7iJO9NxChiag3PYSRU9peFzTX2TNdSg9U6GFEb5yAWSeX72ORdr9p2BceWNdbVCc4LXvACfPSjH8UrXvEKPPDAA3jKU56Ct771rXjsYx8LAHjggQe8nDg/9VM/he12i+/6ru/Cd33Xd7nt3/Zt34bXv/71AICPf/zj+I7v+A48+OCDuHz5Mp72tKfh937v9/CMZzxjvR+SU09EhZLkhku0vZlqbOAVuqHyM5LD4CRWRm5o0GkPK+wOE+rNUB5NDQk7vemdbAnJKSEfvFF5jZmF6Ab7JtyHjg+2ZXR42RLwDKJDdQoYysmeCZROp2uqaNoBd5w2U46EgkvivHgEVQRzyY12XBHJSQwc3jsaTNNrwQsX1wgNi6Di5etPZgfOTVLH2p0cAPnz6O8Fz1Q/HlepZipHdiySE1Nv+LURknWZ7iBVt86Cr0pRfV2J5JSCkvzxehFbQXwpVF3XnV0tayU89NBDuHz5Mv7xP/w3qB/xSKCusL1ug+0jG2wfUWF7qcLJdRXaA2B31H+oo2gbxDuDthptjtu+0dcnfYWsT/pPc61/4fUJ0Jx0qI99b3I/j0DHGme/rduQM3G/8vnuoC9re1BhdwC0B0MHQv83cOYpKr82MC3pOJxLcrSG6hEU1jh54whmK6Lo1ow61TGUkBvZGeqOk/G64h1K9YbK3o3f+xDK8W891K3mBEDb16PmeIhMGOpTfdy60EtSBZ1CoHUsSvoBclx3f0U9a1l9C+pZ6vdPwNLkhiO68rGsT6JeVjsA3fiO6pPxPVG731ztvPaeej98tmqZqDzz4vD+2sMG7UGF7SNq7C7V2F7q39n2Ut+X0TuTfQKgP1/rWcrnVe3G50LPpM9WO/aD1vOgv9WuQ328cz4Y1XZMbBibvfNIQJlOoz1svLrMJ4VjnR77TDc5FHUaEG1cS/C6J/WmeJ2nhUmO11exforGOtlH1Sf9uCb7KD7mBe2Avf9q2wLD3+3uGt7+pz+KK1eu4FGPelS8nIv+6ocpWrZ0AiCkygnmHemHw8lE1wwdDPniDI2tdf/XaA9HZacdoqLaw8HZ87D2GykbdLhpyhGboRwaOVvTeTbbf6egDGXO0fb9VBJjbc9898XkJoKce+q+DbRvdDINjrEke2maYj5edE3LBFpS9hL/J+28FKY60NK5s8AGQqsuzV6wNAXhYOytR8f6tVxEo6e0fZFIquB8bl5nvohehvaN4SPGr6OlOTCyvPfH6/433JeMk5uziJh/oYoz7JOzJlY3UT3cMKXz8cwiETMVgNEOPXTEu0OgOa7QHgIdrc2BCt2h7ixIqg11XkSW6H9umiqJnFo67Jtfc1aUVa7EKkxV/f3tGc0k/5qIehMgd8mBmTMrymbsIlXqChUkiVV8u0Q9547sMuIkp00sMRjM9QdZgrBnR6po7ZxFT1HUEG/va8MtySHMLtzB2HSaXxHW8yAzVcWiqarjDu2mRn08PuxuUwPbwRxnkRxFybLqsjRPSf8bC2dFveHXXNK06h3HsZD6k/JTi6o3E3GuFRzZcc/KE6GSBvSzOi55korTMHVGEJfdoMbwD5kJ3DGHjNxwM4HwwZGRU2ubpqagRAWIEjWlHU7txHOurZVpDoLrpKLE3HvVr2XlxImrOKMTu7tOJEIPCN+fZZ5ac0BdU42cAsvU4z3PmQtMykEdwJgbi/qTwf+G+h+H4dZr+ZTIhSulit0fI7Ni+3nBeM6vfr+RHyel8Ch1WSunqfruiQguhWh5rYzgmdnCS6G6IGzDY6IZjIFsB3OJCwVnSKLHG38spFLCqkyBUyxzOqRZTAsM7Hh4uU2Fuu5VnE6MbnzRvNEMVbmBZzRpjQpOjmnKJDey751R7+dGV6Vyj3gzg0Ilx7rmnP1ryMH0zmIzNnI07jZAdTxsGyL0ql0/+JHzekByZL4UI7+Ses9TxtLkJlvFGRCLKpHBBeq9ZkQOemDmqTG6kpkuLXOgUV8Xea7V2P/x38rXUXPEZtf19Y7lBWs38Nbpq7b6Gly5y9h4v00xT9H/pfV63+oNv/ak3DglfVRK+VH63BykMhgvgTM27zl9aPk9imGsLi5VHOcIPKgwREacA/HwIUWHVJve72YkN56DHM3UlJw3shEG5KZiH4nZ5ofCPDmZKo4KQ8mZ7OdRot6sbOuWne84Ox0JcL9dRoeMM/wg2kZJQcCvxcEJdlG5VyJCaw0eWdfVFCpuxgyUreXD7GU4NI8WKu3LZpsIDeWO39/z3WJmqlElZL44HCzLcXBtSW5cecK6nDJPWQEFflnWaeNatvicfnMv6uXD1IfnQsFhUPNGxJh5RoegMWlScRxqoEWFeggVr1pgJ8lIzWbQNbA78DsLTm6cNA1/9kbX6b8r5CYFOmZiXY8pOUl7comKA5izCv7OyC8gCnHLs6BYABnPC3CpCMgnp2uqUQoOfHDqcXDk6o1hnqIy9NeFT6j3hNMwS+WogfQcuN+Jt5/lf6kwRElmhPL71wh/vDRPAaMaLZ2vl3pXJX4gXMnuywagrXrFeuerOM4Xx6k4HSiHk5XMT/Mj68uoOBejzDwVm1Qtod6kMsZPUcBLleskWB9sqZbyHZ82zrWC0zX1kCm4DpZpSJ+b3ib9RQIVZ/g+htxWQYfU1aNig5oUHKEAHYjr80ap+N14jclSbGKYoeaUKDklKs6UDnsuuVlEvck4J1XOVr5f1dxE1/I7fwBiMJDX1FWHEv+bNYjPPsjN3HsESlttX1O+k7Astfp/cNym8ic1huIHJAbtJZ+vbEeGCYj7x5CK0wqFSl2jTyE3VjJUB4Wox8o2nrdeMMacY6z3dVYmZATKe6OadFfIiXOh4AhYjN75suRCKg7MFq1dHxh8chqEJgSMMw/NzMX/es6DOeRmKibaXWOYo+KomKM4JW5zZhJ+ueg99P43DZwfDrYVwJz3+PMNVhkXCSQBeFm9eXj4eM7iPykLJYNvbGDImRXn1smYH46EFjnk9ok1w8LyMKIjIobCZHahggtgb6HPlMiPlCzuh6OZqXjiU2zgsrtji+F39td1ag7dRyxlA8Csy9KXzFS1Mp/RXPWmbMJ3RpScFUC+gTyCaglcEBwGbTFBAH5lZ6QhB7yRA/FKTxVRvlo+q0iRG+5UvBq5mYmibMeRxqlmTNWOLSVjU6KxVkhol9MpcdLcNhWatvOWBKCBAxug2o7ZsoPrsAGBz3jV6KmIEqGVb0ksIfnLY0oHjdi7kft4e21beCYZoBue8zCQ82U0hhwwAFxnHyM3ZJ7aHYxmRa4GU9ksB+M570klgcpkxJEbIjvimbQYiR/VWWw7R8ZHB3l5L+ZfBkD6IWl1WZqntHqVY56KPZPkMQsvhxMj42uRnLnXjUVQ5a6/FsP5Jjj1OPMZE5rNm517L9xQcXKuEWyTHdVUcqMO3hNnuWdAxVmU5Bg/NWmaWhOi3J4PkbZkg5shV8CmA7YVWnTjgAEMM2MfAbnZjIOAVG88ss3qnsNaTpgLkht5fKyel/iZjNdk/7cI8r/0161m+ZwQtGghbp6SpDTpbD/DPJUa8LiS4+7Fyzb44vA6Sz5MLXyiQwiJDXzTlKjL7hijTzWxYNqDtZbDmVJXH07oNrW6JI2F801wBIIwQsxr7N61qeFX/ra8co3HB8QGmERuSlYKB8pnuTnXnazi5KwFZJEcd5CxXblWEpl+NGbHn+G8B4SdF73fajcqhWSmakGz4G4MGa8H9WBT9WqOt6gm892gle7ZjFfeUzWxrkj81iI38rxJdTLDTCVVnLbpPLWCVBw/PJqZYpQEjVqmXu5I65m1hRq9L5IeM1PxY9CO6mM1mFeJ5HSufQ1tRHGS5BFT0kneq8vw1RtAqVuMsK+l3szFFHPVaZiq+OTrNHBBcDBftQmuJzq/Fv36VLQ9ardl52qNUPO3sciNZZJaNLHfDBXHaqSlsxBtYEk25ol8rbSeSJ+VqR1M0izCzFQ1OmEm6QcLYIiqaiPkhuVbouv125cn/blYm9zIa8ytk7ydawM6HaNFDnF4/jiGz4mLnOLr0ynmKbqnhwXMU+7amc/H1SHmi1j1bNxTcbj6WNHv3fpER16XZyvm9XkkfFVoahWTQRNnRL2R1ymtq1TuyUSnwA9Sq/cca+fAAc45wWk3FSAjSoy05rPvNQzCXue384mL5acTKjIIZxkWuUkRG36vaBK5mcn6ZiKl4kwiORn3TMJKkrZPMxaGDhy2ioPt4OfQDNsYeNZjTb2Zap6a+wyWciZeC1FzNB3D/E5iyxRoKg7aLi/r9GZUKzyfKf6u9khK3XOJpXaox/6PZvlEzke/sb7eurJnmFiBsT5z0xSAoC5TOby6PbOLSz3nfdXTFOFcQs0xQ8UTxKYEmrNxbNV4iXNNcAhyYTpgfoeQlLBz/XGkSUCqNsP1uQQryU3QqKzfliA7Jsk5AyoOsCzJyTLBTPA1WaJj0crmS/693N//PzocuwFAvCxJbrwZL/wZr3X/RRXQwra39KAxp05adZDO15yNuYpDkUOO5AB++KySlNHznWqYz5RinvL6De23r0SEpHrN/6dnUgNoMTocUw4nKq8j5sqoJVUbqsv0m7xnItQbrayA379q+8MyRB/BKuRmzsRztpqjXO+s+f+ce4ITZMxcGV5UFSMr/G9gxhLHyYaXTW5KOi86NrfCrkBywuPKfHH4eUB+VJKGkgF8KfXG89tKPFvZuXD/HOdwDKDajh2/dz4nN80o5+dEmyyBqddda0Y8V7GUZipvX0TFqbQXbSRl5M607aGIEmJr01F5/GvuT/EyJx7MSd4pkMwBmaKqsPWJuVQg6R4qUWd1me7ryqD5Mz4MMXeSKH/7WQ8pL8H5JjhBOnt40q63fQIsCds5hA4I7PPsrfB7S8UGmEds+DFmZ95Pq7xz9mWqymqgGaYqd70E0SkiNzNMU6qKk0vWYvZ1ZqaizNjAaKrqfW/6Y6utqGdsQPCcMYV6Q/dRzVPiN+b8ljlYW+7X6nrg6D3BTOX5RzWIqjiEMUTayDitqDf8/vZvzHkSeeDPJkepdCS8ZeWo/brr/HGAgOj41xqzNwd1uYb3TExCk+FcfJbUmxxMUcL5b9wn2Ynlf5qK801wQJ13BZkDx+1PVPgiCJIDjH45Grx7CofAgNgAIbkJfHkiCcRiZCeX5JyGiqMglXTtLMzWYr8jmgZ9YofDfW+I6HAlhw8IqRmvdm0A2arAEipQ9oCRqlIrjDslZioeOSTzvwDj0ho8FHocqIeBXFlOg+4pzVMxrKLOsT5PmqkAeCHjPBqQSI7nMyaIjld2Vp+DugxfrZHfo8/lYbIG02n7SGrpLIqIFZlhF85mfK4JTrupA3FDtV0nwGciwb6EWSVJnJRIB4/caE7EtbJtQKUMDJ2cpZ52Y2HQZiA5pqqSzLIxLK3eLIXYvTQVhwYPwJf4eaZbTm4sf4XZqstZIjf8GOOSc1QczUwlnY13BxUadGr+F2w7N2h7uWOYz2B7OERPMfWGm6fkAE4BCftGzCEV8FUcMlX15GQsa0B0vOtU8brc6M/Dgb3iU1dvtHqbcYkl/RnHa5apzTEHY1IrVTPsSjjXBAeAWKskcXCiY7AqkjogA3GThBKNkjRHKaqORmgk+DFEdoLGIlScNTCrga5Acs4CudE6GLOeSaIt5H4CP5XXrYDcGOrNVPPU3sjNFG4eITpTCH+QgFIx8ZFjLfni4LDqB/pt/w7rIW8RLxZXKdrD0TQV+N4I9UZD9gTOeOaTJ0Gi7nDy5uD8cHqSU7VjvdWiAPtyDn8TRN1Sb6b4KS1KbmKPM0HEU5hDcpLXNpRlTa309vPUCPVAZnk2bwDdpp61bMO5JjhdU+1txQKTCUeOd/8XEhuP1OTO1gZyUDWdTXK88i1vpkrd0z8uz1wzleTsw6m4yGl6eK78nMA5kDltthh+tyA5nr8DxrrFzVJ07Rz1Jtc8tZZzcoC5DTqz/kZVHITbLRWnBYZ3zhxoWQ4YNJWrv1qU0O5wuM5B1S/Cy0wwdC/6u4/oKdUPh0eRQifsHM5UpdRfIJxnyfoZI+oqZvjexLAYuZHHRS67hvpeahovNk8NiK7DVtcDUxr+z8S5JjgAKTeDtKus5J0+P/w+J8nSIsRGDjbW4MMHUza4miQnV8WZSXI0TFVxAN/fKQfmQHDKkScpyLooSQ6AUM1pxnfbsoFCm/EC09SbxQbP1ICxVL+u1N+igSPibMz/l461PDyaSA4A7Ph5zJHWkZvGd6r1BnxZjMI6nPLbmxtpRuDEm5uqgNDcGlyHE5vhBM/VIEO9CbCAepON0kc4keSsqeKUgqevyAYRm4ulGvLQVUNdiT2FCe03VZFKnP1UcpNDbHI6Mn5MW3mms1wlx8RpORxHVJAUaZiUln2PfjfefflMmUxTysyYkxy0bBDgYAMCYJMUub3UuXgO9kZu+PUSJCfVzmXeK/9aGNdYYqYqYBjYmwp1DWeiGu9ZeUS0bSrf54YP4oIwyPvz37FPuHbLnVKV+gz027jjMQCgAepdF9ZlR3RsciOvP1W9Wcw0NbXeTuxfp5IcS20uSfZHpthqA7Z8zLqRWuea4ADCfpshZy7heJw6p79PWrUJiA0baKrciBYiA5wYxEwlTMVJEpGFSU6WwzGQJDmTcIrRFO43GmYquXo1gTLFeiRHgSQ3lnrTX3T/6s1S5Ea7zpL1N8fZ2DueBuQDdsMWqGvmc6KZqBqF3GiDuXScXVC94cfE2qlqphJh8jnon8V4bitzAwnCZpIboUSqmJmNe3Vyw88vcI53+9ZScpTbOUtIjclrUs3xwznXBKfb9D44MvGZjKAymX2BCSvmFOrfm9XYlGojiA0nNXVGI2t3lTunIwWHEYPZKg77DR7m+OjkNs6CNVOyrqWV5ZTUG3f/iB+IZq6KXYebNXLk/H2pN0uQm5z0CPmLbC6j4nAFh66jRb8BwE4M3qRcBOSGv8NM9Sb2O9eEquJklMt61AGBUVwMYo7FRT53SyheK5pTczCF5Kzli9M1FbBVttc1qja8wMVSDQXw03mveZ/E/ox8NpZqQySFk5qqjtesrq3d8UR0PJIz/OUkZzFkRgTMdjheguSsTG5ijsa5fjjOt6sV/9ehU7F1Ppfv+XcPhY9y7kAwd5AtOd8kOnNUHGW7fD98P+CTHDLFEGI+JpzsuHK78yLl3ad5KrFkjWaiyoV0NNYInzRNmWUU1/Tus4Rpag/mVCpLrP9cUsnRnOmT59TlIePdpr7wwclF/4BDedP9v/IMPWgMOU7EEWIjSY2l4vSEpj+WiE6M5AQoMVPFkEF0ZpmqgHkk5ywm+VLMVATpx2CRHGku0QYFID7jlVi6rcwdKKaSI7U+iwEkW8UpWIOOqzm74Xr1zjbFWOTGVG8K6nLps8s2UwXn+aZWAl+EOLsMwr9GmqU8hadi/2M5NfJUyE0C+yQ53jWV99024fprRdctUG4I55rgAHAOxslVxOWCa0vOTitle51HbIioEJlpEsqNuzhGoqOSHAZScVZLADhBZl2d5GSG8MfKByw/Q9LuQwSGDwz0v5cplg04vIz0f2rGWzIgrK4MrEBu+PlLrEOVux3Qn5dWdVRzoiSoEfXmNJ2LAZikj9q/rL/8fzVfDvzzAIXQwP5/Kcfi01qKwWGOv+NUkhNL+NeM76trBkd6zU1jyIXTbmrUx7shVBzO74bMVJzcaGvpWTjXBIe8863ZLKBX9DkoUW00YsPVGklqamXA2TR9L7LlhnwvABNo4ZMcul9UxVkaEx3mslHyOxYiN/L/WPqAomUbhIojSQ6/Jic3qs+X5rsgB4LEY7OcZ6dizmARPTdWroAsizq3BxUnB5o5kbbL47wyZF9/OeVriopj1d+8Moz3ld81052HmaapLKzdjU40VQH5JEd7j9mRVAMZ50uTNC07TiT5G68zkJyLMPEyUA4cAFkRVLPvl0FuAtXGIDac1BCJAYDG6MyaeosdDfANEZ6R6HgkBwhUnL2gkOQUqThAGBYf2z8R0Vle7kzJyqOS6eQnzVKpgUJ1Tk05Y65kvpsj85vn5gxOdIzntC38chIkxy9LnoozxRQjTS+eaUoqbcbAfSrqjYCm4kh/MkngAcU0q01QI+RmimNxDGuYpmLXXC1ydQVzFV23N0lp+3o/HM3RmEdPdSy5X9fkV95zTXCA0VmvVTrzpRi3WlkNcuM5ETNyw4kNKTWbZucIzabANNVjhy3oh/aujURyxsPHgTYwUy3lhyOxNskhTBigk866pzBoWCqOp9xo/jqK6pKa8WoDwpLqzamRG3l8rCnF6memiuM5ZBaoFNbg7c6f2QTXNu25eigIvNueIDl8W3Btw0TlECM3M9SbJclN7vOf0t9mB2zMIDmBI73hdNw1UP1w2k2Netui3QD1FkDbBSHi3aZWI64snGuC020wOBn7jY3/9Y5XfBeS9yggNw4JckPEhkgN97s5aPSR/WTXoKlb7ByB6Sn1dtegrjvshktUdYsata/i7MtMRZgxCyHkqh0l14vuLxikSjsQPtOVgwHfppqq6vxBIfgdYlBwUBaADc4twGr+C5kDU9DxC5IzackSiDqotaGqjORoJpio0nYG1JtUfZekT/PHIcSeU0zFMesx8LAjN/x4k7DM7D9T7yxbSWbkhvxwkmtScQwmK+ePQ344Fz44+ZAOe8teu4zcaGap5mCnEhsiNURomtSoqQ7QI8npr1ePPjjukD2TG0KBPdlqkEuTHAtL1JugrDnPXXlGpeYpujed676vFGkS3rvgeiXqjeZvZNxLJTuKycorR6mKQ5eTAzrC/61rpcwvlukly29sT06ymoqjkRxZf7PSHSiEXSM3S6yfdprkhp9XouSUHF9qTo+ltNCuZfrhBNcfSc74PesnADjnBIfCxJMy9lIZQGXd4g6omeTmaNPrcwfNDk3V4pApNodKDTt2rXmH412jEJ2B5IB8cnpTVaDiANGcOKtEWJ0RkrO4o/lMW7f3m5SwcencCZT5Lchj1jBNFXfqK5Eb7Tivbg1qzpIqjkVyCNq6TN41DYUCSA/ea6s3JW1TI/GWuSr//uxLAbkpuc+SZHA1s+AMh2N37AI+OURcnIBQh2SG++G04GYqEgEG5+NBuWk3FwpOFuiBAxjlXqHklDqi5ZAbfowXLWWQm6ODrVNtJLE5rHdR/xvat3WmKY3ojCRnZ13KSz0f+uGshgK5dQ2Sc1b8bnKS/kmSA4RExztem/EC3gzaw0LqzapqwQxyw49XFyvk22eoOEBctUiXD2MZ2PccIqpfbz/qzXg/m/R5++kVCKKT8ieTGaMJJeRmdrvOGIOXeu6rpe9YGX2QD4Bdp5qpvNXFObm58MEpQ0kuhBhy124B4DsVZ5KbS5sTHDY7p9Rs6tZTbTZiBNxG2NnxDgHJ2e4wmqq4igOcjpkqgnJ5tpzkrEluJjtHa8cKfxwgJDqx6zgoj1MbFKaqN5M79FLH4sgxVeQcrk569SuDyE/xxVH9T5DR7ygDuOY0m/OeTj1/C4NKcgCP6NC+AOLRq8QGmE1uljJN7eW5r6jiWH1VzNGY++FUbeeiqmhl8f65jyoO2s4jOU65KfDB2cv88zWveQ0e//jH49KlS7jlllvwjne8I3r87/7u7+KWW27BpUuX8A/+wT/Af/7P/zk45s1vfjM+//M/H0dHR/j8z/98vOUtbykuV9foDsYAsu2nWchwKuYOxU3d4vBg68jNpYMTj9xs6hbXbU5w3eYYm3qHTb3DdZtrOKy33of2berxPCJFh4MSxJ2SN81OzaVz6igY4FLh0Hkz2jOi3Ih3EVUTIx28Bau+L22a6ppuPx26KIvnW8MjFA3QMXScV26n9LJrRPqI6HNR2lhW31Mhm9wUlWdh5LTLziizjGYNnov1Ma4dXabilMjNGlizfSXrzvD+Yv2Fd43BStI2QLchJUd5cEy1aTdV/72u0BZkNF692t9111146Utfipe//OW499578exnPxvPfe5zcd9996nH/9mf/Rm++qu/Gs9+9rNx77334t/+23+L7/me78Gb3/xmd8w999yDF7zgBbjjjjvw/ve/H3fccQee//zn493vfndR2bhN0ItMmOBzMaWCkXrDyc3BwQ6HB1vnb3Pp4ARN1eLTDo4HUjMSGyIyj2yOcVC1wUeSnUNGdLjy0wzEh3+PlnmhZ1A08C1EcvpzjM5t4rs/LcSkedrPfxP/XkRuZpDe2R3vRPVGkhuHukt/5Dnu3Pi9Y2Wy6pWangKIDuL8/ZWEO5eUdw1MITnBcbHri3odDLgPt4ipmn2mXneG+rkU3Ngq3D96t5BRzew21eA4z0hMPRIbAEXkBgCqrutW/aXPfOYz8cVf/MV47Wtf67Y9+clPxjd8wzfgla98ZXD893//9+OXf/mX8Sd/8idu20te8hK8//3vxz333AMAeMELXoCHHnoIv/Zrv+aO+aqv+ip85md+Jt70pjcly/TQQw/h8uXL+Id3/L/RPeoRaA+AdjOs6dL0f8kmTuGWmlOm+x6rLIZ6w01TzUHryA2ZpY42Wxw0Ozxic+KREiI2B4NueFTbBslr7WiBPOlqHA/ft22DT24PcNw2ON412HU1rp4cYNvW2LUVjk822LV9RFXX1tid1KOJqq3Q7apR4hRcKG9xTPt5Jc/X2rBxzhpJq4BlZ8NaGQPpV5gHpS9OcHxui86Z+c4wTS3SieZ20LW+LyA32v8E/px5nhalvnt1jos6oi7K9+u9K3aP5DulW2WYXhY3Tcl3kDsnUdpl9HkAybrOz4kRnxi50c49U+Qmp38x+ja1L7TcQgv9dqJ91fDe6h3b1vX/1yf9ufxvfQI0Jx2qXf+3Pu76c9sO1bZD1QKVslZV11TYbq/i//P2/xeuXLmCRz3qUdEyr6rgHB8f473vfS9uv/12b/vtt9+Od77zneo599xzT3D8V37lV+I973kPTk5OosdY17x27Roeeugh7wP45ilCstNe+IlV9ZjnxiI3pNpct7nm1JqjeosjYYaSn0durrnjSNEBMOz3TVVctZFmqorNbN1zmDB45Sg2yWO09mvNqFeo3WchA6yEquTE+q5cWf9Thdzw+qvU5fE64jjtnjNUHPV+AzTVQvuoxy/gH6XCqkep+hVBVMUBssyyMaW1bcrJTQxnktxEjjs1FSdWB4WZiltO+u+9iaodFB2u5HSeC8kQ5FJQ91btrj/ykY9gt9vhhhtu8LbfcMMNePDBB9VzHnzwQfX47XaLj3zkI9FjrGu+8pWvxOXLl93n5ptv7ncw85T2Uman8s5Qb+qmc0n8Nk1PPCS54arNUb3FIzfXAgJjfYjscJLTf3xT1cFwb8qMnLdoZz6mmK9MnBLJOTVyk9npB8j0V+DmjiXS169JbrJPbxSiQuanwSxsfaTJivvklJa1ZEDPef7BMZMyci80aKfE1uzQfLFBqe+xZ2PuF2TWNE2XkAWJ0yA3U4+fiZz+zwzYqf2//bGVG2c7JdxpdBvpiY4bqwvCxPfyiKrKL1DXdcG21PFye8k177zzTly5csV97r//fvU4La9EMhnUjM6c1BvnVLzZ4dLBiVBuetXm05pjPHJzzZGay5u/w6X6JPnhROeImbbIJweACzknUkNrW9ULNco5iaxKzt2nkrMUdPNB+rzSmW3q2BJ/hdly/gzkqjf9Pp/cOAJDm5rO+xDcMexY39Rl3C/7N4gNxoAuvy9hdkmWeYY6k4us9mgobfxZmKRHVZz1cjxsyU3kvH2qONH+RqQ+4L6A9N5iKg7t93xkN4ZDsoFVw8Svv/56NE0TKCsf/vCHAwWGcOONN6rHbzYbfPZnf3b0GOuaR0dHODo6Cra3DVA14gEu5WRaoN7UdYem7lWcpmo9cnNQtXjk5lr/O+otLtUnOGCG66P6JLj1tfZg/KI1HLfYB7Cp+x98XDU4QYNN3WLXjss30NINu8wwcRmCuMSAp4Y1GiGQ5jUycpOkzj9rMFfxHepwMl8Iu46H0yY3S6k3BDHYWcSdtre7ClXd9Uku6VyZA2qA+87qY3BMTkoAJemdqeakorAwg9yUItEOc0OSrTorn7+JyET0TDsUL9GvKGkMikLBFwgbD6/JfKVYTi6XlXowU7UYQ8a7DYBthRYdalQqO+nqqmg18VW77cPDQ9xyyy24++67ve133303nvWsZ6nn3HbbbcHxv/Ebv4GnP/3pODg4iB5jXdNCfv6OosvGwXLekHpDpqlLmxMWLTWSG1JsPr25ik9rruKoPsFRfYJPa67ioNoFH34MbSM1Ryo5zlTV7MZlHxQzFffDSYXbkvKyuqd+gakKmEZSSpfw4L990WdgdOCpMFg1aoqdG5g7TpvcRDBJvWFmKa7UVHVrfsZjhALETVWFdSlLqRORXAEMZWIRv5s5pHLCubZyYpyQEfmmXetMk5t9Y+0Q9pS1Q5ipuIozKjfwlZyaTxSqYjP66on+Xvayl+GOO+7A05/+dNx222143eteh/vuuw8veclLAPTmow996EP42Z/9WQB9xNR/+k//CS972cvw4he/GPfccw9+5md+xouO+t7v/V582Zd9GX7kR34EX//1X49f+qVfwtvf/nb8/u//fnH5uP/NYrP0mHozgKs3ZJo6bELl5vLm73BQ7ZxSc1DtcKkKVRsT9JtaBJ3yyfCDj9sGh/UOx1W/JtW2rbFpdjhuN6ibzs9urMw09wFzlqHMIKNJ1wqUnFJiE9sXXWU5Z4YPmM8+J8uxPD7A2o6YuVhSvalDolJxZ3p2LF+DzS04S0oO4J67t1yJnDlHVBwLVsK0HN+axfxJ1jZJGc8ilkQOsCPJ4vdK7D9L5GZJeWGPKk4K0XWpGvjritXo1ZuGhqh+2QZ+HR5KDgBnxkQF9CHdH/3oR/GKV7wCDzzwAJ7ylKfgrW99Kx772McCAB544AEvJ87jH/94vPWtb8W/+lf/Cj/5kz+JxzzmMfiJn/gJ/NN/+k/dMc961rPw8z//8/iBH/gB/OAP/iA+93M/F3fddRee+cxnFpVNLtXgtmf7MUzr2KV6Q6YpcijmUVLcHFVMbgwc1Vts2wYHVYtj9CrOuJQDnJnqrIGed24oeozkAHGisxS5ySnPEuCvS+tcoq8z019Bw77IzVTfGyI3lGeK/gd8hbKpgR1rA5TJu/+fcWjFVJVNZiIDeslgbvVPp0puCk3GhNiEo4TozE3MedbITeoaan1TSE4JptZjq/7SdkthrHZDv9RWqInUAKgGU1X/nSnKdXkU1ep5cM4iKA/Ok7/r36F+xCWX/6b30IbnBMVz4GTlv6mU/TWAw9bJqc1Bi+Zgh4PNDo88OnaOxZ9xeBXXba45h2LytyFTkyQ3nPycKC38anfg9pFfzknX4Gp7gGvtBn+7PcJJV+OT2yOXF+fvtge4enKAa9sG212Dk22D3Ukz5MSp+g5+yIUDDA1taqNKkYiZuR5yCYW2YnMuSgf4WJmycuIAyypoDwNyo95LITgyLFySG0lstKzd7fBsKQ8UAOxO+oci6z/PjZObFwewB/Q5A/mZUW5ityus+0vh4URuStqU+jxz8pKlXLEy+s1UfieaZFW7/n7VDkNuG3g5cfo8OH25m5OhfR6j98sBULE0bxRpdbK7hvf+wsuz8uCc+7WoosrNhBDMAGSeGr7mqDebehclNwdKT3hQ7VSSQziqT3znY/RKzsnucAglb9jK470fzpYqad0Cu4UVnaxIiuFvjvRaaKryjpsoFZ+ajX0JM2GGoypwNshNDNHIKfjkRhKbTaMwigbY7hpwm27XdqMJiz37KSpO/LeM/3urLU9UJ07TLKVhrqpafL+lIjzPILnhx3vPVKg4aynHKRVHrrPm1e3Wfzdt2xebrtk2vakKALBRFJyCx3SuCU5XD+NiM32QM6+tVVYlcuqg2eHTDo6HZRS2nmlKkhuN2HBIknOpOsHV7iBKfshMRaCkf9xklUSpNDol38PKJGdfKC1P1EdjKsl5mJGbmHrjLqP43nCzlNtVd47YNOYEZoctiOQAbV2NC88OkVUV4C3OGSu7fN85fmDZZvKzRm5mRlXNjXakaySPya3DZ5TczMJSkW/Wu0r0TZLMdvXQ0oSpCmBEB8wHJ1myEeeb4CjEZpb/TaxOOIdHP3KqqdphEUxfvSklN4SUknNQ7VxDk6uNH9Y7HAulpqlbtLtmcDZeoGecm+8hNTvZE8mZ0yGVOl2aKCU5kWgTrSz6sQt1xGsMsorvDQCn3hC5IWKzqa2HPXg+DmpO3XTOJ8f543AVB72ZVgsZt7DmQH7WlJtSTFFzlvaZ6wuw4PVWJDdBf7KAijPtHF3F4fsJUsWpBsG0xUhyAH9e61xFLhSccnQNXIU2F8ArhYiecj4Bg3pz6eDEqTePHPxuKGpKOhcTLrG8N1eFySkFaaba1Dtcazc4rLfY1k2g2lAk1WJYIedDfrTKciRn36apqBNqLE/IPiNxSlDwGlLRaQEUvxtJbojYxLN100PYDX45tSM5uwVUHGA6yZkcCbQvcrOgQrAkiurvGSA3VjoOrc6tNYlb4pqu/2I/JxAW6tEoLEmOPEdbYsnCuSY4FCKeJWnOnSGIPBxcvbnOWzNqDAm/VI1J/S4pCf0sssNVnJSZivxwCIfNDruuxq6tg0gql/xsCiY4+k0OcTQ62GX8JJYZ5Jce9Er9xRZ1VrWwikqj3KbpTPWGwMkNEZsDzQdHuWHvs9MrmXR9UnGqpv9/iooDlL/vM09uMrFP03Fx3T1lcpPKM+byMcWeX8ptILd+TjUpDupySsWR2znJAQA0QD0suumIzYWCswBKHYxlHeDqjZLYrxmS6x1GTFOUoC8Hl+qTbEXnoNrhKjghalVHY/dTGpELB/DzgUxEqoMww8JzfH5WIDmnmbirNJQ4dS11+1xys/B4Vfy8I+oNJzcuoaXFLJqe/JwMhGbX9D453FRFKk6uiXCOg+1sn5LC98KvNZmE5JjprPa9ANYm5KdFbuSxvA+O1rGZ/V7+ZDMeMh5TcQgyZZsjO07BySoKgAuCA8A3T+UdP22g8xL7DetNUWZhUm80cnOpOlavd7U7VLcDti8Omano2tIPZxXMlJrVBjbDzjzVJr00FnM4Ljhf3b6n5GdLIggNp6KwiQk3TXFy01StW3+N1mPzif3gi9bQdXqfHDJV9UuYdKYvzlgY+B16xmx4CpYgN9Y15PaidrOQSpCLyW10BgE0sTK54eeYE82FVBwgTdBLnI2z1trjl6FtdF5BH3iuCY40Ty2ygriGmpmmmCzer+a99dQbbpYCbHJTAjJTlaKpO7cmVRKFkVR7CYmMNOBSQrRPpJKflZKcWIeyj8UGp0A38xoH15LktF57I6dijdyQPw53OCZftOMdgIYSADJT1RBVZfninBnzS6IIU+r1amHHE1Sj2e1yws84S+TGuv5adS/f3zEeMm72yZr5uRU+sXXZROBcE5w1EGsADXNu5OYprt4Qck1TFrh6U0puoqHixZE74SbV3ybmTBobPGbMUpIz6lM0ScVQQnJWITdnzJ8DgDMD8+UXXL6bQb0hckOqTT/BiD/I4x1w6eAEV08OPFOVpuLkYC/O7isQG+38rN9RoBLI66+GiY9/DXJj9nspFwlW53LdBXKjTouvgWVUnGrn90nyO513QXAy4dahCvxnBieu3DWqjLol/W8AuNw33DxlIabeSPMU97/RyA1t41FU15QIKS1UfEkEkndGp5G0NTOSU272CTvrfRKbqc7GKZKzSnI44GySGwUNa29Ar95o5GaTsYjX8a6/3tHQXNq20lUcIGmmWgL7Vm1i1zpLOaaimFHM7Ge2BLnJ9f0UxIH3kd57Kc1RthByVRzZD2mkJrj2BcEpx7IrhpNNnm0a/G949JRlnipVb1Lk5tQQsO9yciOPVRuxQImpyirbwwFT6+xZz5+SMk9Z/jdAaJ7i6g2Bk5vDyASDYzc8tKbuTVWeiqPOUNdz+DQRueSp1u8VSF72fWdiCXKTjdLAlgjJKcIpqjhexm5ytldIjsybk4tzTXC6Gl4j2If/TY55Cpjve6ORG67eaBFXctFNjqWXayiavSwxU1m7k7X6lRlOlkskg9OuaZXBxFmepIuEfu7/IXIKGMPByeeGk5sD7QEPpGdT73on/C1wXDXOH2fX7ISKIzpwWReVure4o/spkZvs37EvkrNgXV2K3GSpN1rfl9Ef2uTB6FMWdOhOXSdHxZHLk2gkh8gNfb9QcFZErNIH+0REh2WeOipVbJh5ylJvOOT2WKbjxZBQb/xjEx3JBDl2CXtzNmLtnPadAYHo4UBuJg3GbGFNAC4VA/1PvjcauVFNxO0GB80xToYHtqkbHDY7Z6pSVRxupspcnyrfaXO62XAfys2pkpwV6mfRM9s3uaF9BslRVZwZaTVi0N57UV4cdk9LvZHfu6ZMuT7fBGd4N94DE/43s2/B/G+keQoACwWfZp4iWKapVciM5WhMDSnR6L1Vn0vuCagNee3ZShS5l090IGuqOA+X5HBmWazkftqhbPuG5bzp13oLyU3UB6fdAEOWb+opXRJMpuJg16Sd7433P9vR3Th175F/+yI5+47ci2ENcjPRB+esIpkXpxr/z1VmLhScCZhsngoclIfNoiJTRIdM7re0ecoiN3Il8X3CCwGdouQsOVtZahY5pW+ZeO85JGcSuTmFfjN3cLGO4wtrknmKHI1JvdHITczJHwDQbnDd5hqwPcK2rnFcNYGK07Vd2tk44zdVu6pskJ1KblIDxNpOqSWq5sp1cTIRPE1yw4+nvlEhPMUTvz2rOBbJqdpy0mPh3BOcNQQO52DMlmdw/jcscqNP8HfiFtTMgWae0lSavZihNEQqpLbic/C/ek1fvUk25rVNVSt1uouGEE/ogAGcbX8bCbeAbWieIvWmqVqn3mjkxloCxfNRazcu07e3lAlTcaq6KnY21o7NxhRykztY+Glks1Fcf8+KQliKqW0rek3NYV6/TpfwwZlsppqItAoZqjhyCQd+XLrvyi/buSc4KWTLZolKTf43BNm5zjFPEWRIuIZ9Ex89IkYnN8kGrcxYJpuqppKcuZ3yHlScNclNaee9iI9JDAlyzNUbACa54SrqSdcEbfGkq/sFaYWK44rRdKGzscTK6mGJiS8bExJ4nnbo+KrmuYltK6reBIkq4+UPIveoX1zCbLVQHc2NqOIqDjkdU1O0rnFhosqEW+CLxvwl/W+Y4yNftZiHh8+BVG+skPBrGU7I+4Dqd6PMwNVz6b1QIwa8Br03x7ql+u3IfaPhxRkkZy1yM3XQCFIDMDKahUL/mzpSl7SIKWkepm2S5By0raricGfjvhAdMJiaSsxU2SghN0uELvPrnEI+lRzszd9oZXKT6ge9a/I+UV6b9YmT/BNXMFV5Kk7EVKWRnDk41wQHKPO9CSpwjk+dcDCW0DpYC7G1pwiW343nhGz442gLbaooWGDQhEJu6sjx7a4KiY5CclazO58R0421MOOsBRn3GIUz53qxcys2iSD/G2meAuCpN7zt8UjGa+2BR3KutRsc1dtAxQH61A/STBVNTLlCJNGq5EZeM4Pk7EPF2asT9Yy2NZXcpPpCfl6Xq96UqHEToj/nmKo0kgPoDsolOPcEZx/QHIwp/80SiGUr5ts5rrWb/Sy0CUW9EeRGa8w0YHVDXh46hohOVqPObdCxQWetfnrmQFc6szkL5Oa0oTkUayka3La2bzt0nlRxrgozFYB1o1uUy5aQm9T7XJKUrEVy9l4n90xuYsSGH9PuywxYSHTkey9J/idJDm0DBt+c3JUFGM41weka1mco5qlieUz2dc2Y5I+bqWRHW0p0Ys7FEvIYqd5cazcu3wdgqzhENGarNwa54TNwDo3oBI07puLkLuOgEY61+5CJocOlKCU3Z4bYxNqfMihIB2Puf0PmKa7eHBk+OAARmxMv7cG2bjwV56DZeWYqHk2FWALKKeR2BrkpiVBL1rtTS/1/9ogNMJHc8OMi5MbqE7u2jk/4MlXtIvDTE68iRnJipipvvyBWU70r1hAyzzW6pnNrUBE2XiZVXxYviaCKQao3J10TmKUs09ReoDRwSW6IEMoPHcOPq4Qa5IjTnI6wEp9PATxsyc0KsPzeLD8coG+j9P+GpXcgkzORqmaonxWrh/xZznquE8lNJ8qQg6xzFgq8WKw8S6FmnwxMDzFPm+ctciP35Sg+FiaXf2bfqJEVb8Vwvp/3xRP65nOt4DiU5h8oACk3Y8Kxcf0pjqkRVFe7A5XceMcoqo00T/Hv5EC5BLTQcL7qMyctFvzZSuvUnJipapKKcxpYUcX5lCc3zJFf1h/pf3PEEmuSesOJjZxoXO0O3P5rOHDtk1ScT+LQORs3tehGc8xUS/t+lWQOz0Cy/q3sj7N00r2lEfUHK4iYAvy6GyM23j1YP9h/X9hsnyzA8Nd4DFNMVaqSI+9XiAuCk4kpHQZfQRzoHY2nIsfBODhHUWz4CuLaauIcu7ZCO8GfIHhWhr1Zkpsm8nx2gyxLJKdG7cuzdB8rqkqU70yRHANzyrkauZEDyVrmigkDtkzwJ2El9dNUVNp2tTvwTFXX6g0O2p48kbPxpm7LzVQ5yH1XGc8qpmxabeU0SE5RP3tK9oclyI1GzHOVGTLRE8nJMduvhoKI0JSpChiVnJr53AA+2bnIg1MCVumK/G9YnVE7FKrAdedFUM11MOb+Nzw0XKo3mmqTgrXQ5lRYif28rLOC3Ghhvm1boalbk+S462eqOGvAtMWfwho9i5KbVJWYmBRuKqwVxCU0/xsNRGS09njSNSP5qYGTXe9wvK2bwNmYoql2aBzprhpGcqAM9NZ7L3H8TpCbHJMtHaMNhEtNBtIRNmef2ACF5MY6Tk3ql2924oRmkqKdsb0IBf1XDskBwkSAHqkpMHacb4KzlmnKkM0PZ+a+0cDNUbnkhpujTroaxxnkJxcp9cZb8bnpTGKzYUkRt6DyumQjroHXTYcWOFUVJ9bpzb1P6flTfC5MzHGyn0N2JgxgMfUPGP1vNPOUBQoV599JxTkciM4ncdBnTB7MVFXd9mtTDVgyu28JuZnii2a1mehvKJg4UPkqQfqK8HAiN4mcXykzvVWnd4pqs9fIKg2GyaoodFwhOQCCjMcluHAyHjA7B17iSR42tiPjXFjrTF1rNya5yVF0CFMajrksA/oGLcnNptl5H4CZ9obvRIKa4Xw+66nqTiWsa/qV5DpAZh23QLh29LhcR1WgyNHSxCn2LAdKW5t8rWrnFsL1o65obStyQG5d/awbvS6qmOrYXkJuqG1oH1mcKe1lQsqCYufhJerkDCxBbrzdEXLTZ8m2WWMjzgn6wch9++3+10X7yIx+RlpHtMWuOeZkMznfCs7CcBFUDBRhQbD8AK62oyMj97fJXXhTqjcx1WY8JwwPP97ZtUlbaycKRb3h5KZmBKYxGuO4vR9ItmicX1Cv3ghTFd1XU3HEbHMVH5fEOWv5/ixKbpbClAy4OfePkIeYn1vg2D8hgrEnOr1687e7Q8/8vGl2OGmaMauxV2bMUrVSKiEQn1SYoONY29aUnKX8cSZhwbWflnZ4LiE3VmoMzQcxlpFbM9nTdU5VxSEoJqtsp2NANbGZUVYJXBAcBSHDLB/MNDsqyeQlEVNTnItT5IbvP2YRVURyTiIkh2BJ2QEUp2JCoNREzQz0UnaDyaoNIr2quuvbVab9ObXPOn4OinPxiPtadnQTM/KnLIZcojMhQR1vZ5sFlRsJIkPXMKqlZKYiPxygX6PK+eFwR+MBi5Bc9pyyyI01WEYWbZxEctZApF5OnWQQUr8ldf0kueHHZpKbGLFxt6g7j+TQ9XhUFX+fOY7Gi7/bQpITRE0tlDDzguBgAfOUAc7Il/S/sRL8Sd8bTbUhcnOS8KLetjW2GURHg2z43PdGqjfN4ITd/x8fAbdtPZChkeQAdVTFCaDMNHMb91kIo+ZEZ24H77Av6d/y0Vnp/nyJhjngvjhH9QlOusZFU7ljDD+cqulJ91w/rOztOSYKvr+Q5MSvh2VVnIkZmXMROGSXrpOWQuL5x8hNjKxvd41HcqI+iQx7J6gZJMffp5AcYBbROfcEZ+31J+Uq4nMW2bQS9fnrT/WvNKbaELhzMUVQ7QTxyc2HE3MudpsYuSGzFCc3lu/Eya7xyM9u+G19h9A38F3LGjcQzF5yojgWi/JIYKqKU1yenH5h5gx5coc528dnufdBvnCWskrt7lLVkxtupiI/nL/DgRcuvlsmf6fpWKz63MjklwzcNyNYiRoYB5HEzHkvg+TCqk0uikxdC5imgJDccGJDyvZOvI9Ns/NIzngfoeJQORZQQpJq1sT+TJqqtPWq5hCdc+1krJGbRdb4MTpfPpvU1sApyTRMCf40B2NObjRHY67e5K5HNcm2OyRikw7BwKjckEPdQbPrM8Mqn4NmN8yQx7xCcoajJsjKdLIjaO+y2BkyE6urQTmmqchzKHFsXusZpRBby2zpiEUiQTIwoE/c6fvhuDIZ9W/ys4qtrG6QGvokt1uDc2lZV3ROzyLcih/k0ljSNAXEyQ39zz/82D73U1ue0ThjSY+Sth09rgqP9b7HnI45Is7xFs69gpOC+tJyUptIf5MZSf5ykErqB0inYu6H4zsYl2YxzglPleoN4Gd5bqpWjTTryzX841X8GtJUlVJxSn/LqWHhvDjAQpFYGeeehQSKjfBetBz7c3GpPvHalzRTHdY71Q8HQE/y0b/OClWRCSfbNJUIRwbGAVVOVLREmRLR1dE1rOBwbOaZKti+RMK7LHJTaJoCQn9EwB8zeI6ypu6TSZKSo11/ah+4mo+h6NNSTscXq4k/DCCT/BGWCg+X4EswcHJj+dvE1BsuiwayZwqMZUv1hpumyCRF5CY26z7eNf1zVEgOkTVVomXwGtUpLRrIsRohmGiaWoLcLf2bcnO75DhnzgW1W80H7rDeBX44VV2VRx5aiKk37hif3Gizer5tTBbXheU8DVNVoc9NqVITS2hYcn4UE0xTOcEWkuwQyfGv1/sjFvfXC8MMiCgkOf0xw6kTh8wLgsNQvHq4hj10tCmkyI2q3hREUHFEQyg9mbV1qfSdaWqozURuclSu4x2AxleZpC+Oh4Vs0HvHVBVngmlqadUqFvH1qQLphwMWXZVDtqb6fAVRU8zvJpX+311jaDt+Rtz4WkarOxwXkJu5JqgpyxeY90yYpsbvumnKIjexJH+bug0CLqJL6rB3uk+VVb3XBJLTHzeN5JxrH5yziCVX/Jbk5rjdOHLD1RtriQa10UQakpVGnxo1N025xUcZuRkXRxw/ANw+Op6u0ROlsCNxq40rvgReh7lW7d/HisQaFiY35M+gfXJwJsx9sJNa8qVOUrjkJfkLe1pSajfKQDYbVl2KkJtqWN3cHdp06j6fDIWkKYbk+81tX0uSm0wfjZJ6nE1ujMAK+T+PIu2/6+SGfA85aN/aLg9J1MZHQPXNKfTJGY8bPyXFXA0f+9jHcMcdd+Dy5cu4fPky7rjjDnz84x83jz85OcH3f//346lPfSoe+chH4jGPeQy+9Vu/FX/5l3/pHfeP/tE/QlVV3uebvumbZpVVe6ilHfTUmUVq9e8UtOOlv422HMOxMFFRBBUnPMXOxUMHEyzLwNQbwDdLcQIjQUSHH0OOx66xN7u9mChMWI08k+jkOuctDWsAyen8cweIOQ7IU8+TkYAS19oDdYmTUlyqT4YoKt/HJ8fROIYs9UaQGyIwRF7oOyc2VCbN4XVRQsZRspYZQza54Y6n0hcpwyk1Vo+jdTxCbnJMU4Cf4FQGXHBiQ99lsAVN8Kj/kxnerezuY5ntXdFzYudFiI6HRN+2iDVFL8py+OZv/ma8733vw9ve9ja87W1vw/ve9z7ccccd5vGf/OQn8Yd/+If4wR/8QfzhH/4hfvEXfxF/+qd/iq/7uq8Ljn3xi1+MBx54wH1+6qd+as2fMmKlQWdqJ0vg/jeWMzFBU29kBmM1B06hegOMzsWbYfFDasCS3GwiHwC+yjM4JQP+TEaNJlijkfNzszLvzrhPbn2bETXlLjFhkCshOmsiZ9X7q4LYaPs5tCzipN7IKEhXLzUTgzfwJovpQwkLV68Ln7A0bMCUHzpGi26MqTiq426WX4qxreBZlEQu6WXIJ+xzIrFipinpdyPJDQDPL1H70DGpnGGLE9bSPkw5PkZy1AjWBdjJaj44f/Inf4K3ve1teNe73oVnPvOZAICf/umfxm233YYPfOADeOITnxicc/nyZdx9993etle/+tV4xjOegfvuuw9//+//fbf9uuuuw4033rhIWZd4kF3TufelrRbLcdI15DICYFrKeAuxNaY4sTk2/rciqEqdJb0Vw1l4o4uYEuSmNCHb8a4ZO4M6dDauUQcJr3Jz4uT9wInnGP3SrDJNcCoOkpztwa+hxDdnLiE6bhtchxMctxscNPZyJ1e7g2j7u1Qdu2ziFE0l8+HI1cq1ZUfkyuIWsrLnMtIu0zDwgRQI/YHG5HBs0cbBMXXVVP8F7SVZN6eqtRn+eKWOxFZ5LIIRLCpskBtgNMfzyWdTtU6hHN9zjW2Gf8oU36Ph8tMh+rzSrPLRJR0yb78K7rnnHly+fNmRGwC49dZbcfnyZbzzne/Mvs6VK1dQVRU+4zM+w9v+xje+Eddffz2+4Au+AN/3fd+Hv/mbvzGvce3aNTz00EPeh2Db++bNSuumc+tQWeHPEnx2mWumis1Iub+NFjElfW+kgzGRnexOT4niaITCQqYpSW4Oqtb8EPnhixtqKg7JtMWYOKuehCnnxh6/sS9WN4vITUHeiRLfnJLyLYHcfE9TQXWTmxZ6c1H+78x2LHbX99d346kY6APAW8SWb+dtRV2wcUkVJxOTyY00SaWOm4oM09T4XSedMskpAC8P2OGgcBPoO33IRA8wnxxmpkqh6H0lcmVZn9g1vP0JfxxgngCxmoLz4IMP4tGPfnSw/dGPfjQefPDBrGtcvXoV/+bf/Bt88zd/Mx71qEe57f/8n/9zPP7xj8eNN96IP/7jP8add96J97///YH6Q3jlK1+JH/7hHw62d/Xqbg4OuYnHaGYoUeqXo5mm+u2N+p37LcjMmRJyFqCtgyNDw3lY+OhXM5IbIJKvpN0A9db9JuoYpIpDs5i66fxoKup4eOi7nC3kRn2saNQtznCcW3lZmbPIjTYAyG1GHSkJxS1RdKYoTMdtAytz+LX2AKj90O+S9A0H1c6tS3VUb/G3u17lIUdjzSdMnaFGyaW2rXN+NxJyAVsNPAtujzHVv1M+Jy5IuX6G40IzldLuvX0zFg+OHiZ8ntz/dZiglCc5BUbVRhsveN9NKTNOdg02deuyu9M9W4zLhbiyF/9efXPJ5D8rLUcisgoYSU6pmlPcXf/QD/1Q4OArP+95z3v6wlTKehhdp26XODk5wTd90zehbVu85jWv8fa9+MUvxnOe8xw85SlPwTd90zfhv//3/463v/3t+MM//EP1WnfeeSeuXLniPvfff3/pzz4VTI2o8pdg0E1RXL3hEihtJ5+G7JwKddioeWi4VG84uTmqtziqt6YPDpkBDodjnHlLqDjc2diLpmKYbIqZ40ejXWsKKvExkDtDU2fHubPbDJ+GXERnfwuAfNKutZtiP5wYjuoTF11FGY1VBLP+7Fs4aM+TqzeS3Mjst/xDag5B+q557SZTxVkCaymL5rFzfXgKHItl1BT3RyRIcrMZIkWpXrnoUuGA7K4pyJMWTRr+hrLtpe3Tilotjaxy2+syRadYwfnu7/7uZMTS4x73OPzRH/0R/uqv/irY99d//de44YYbouefnJzg+c9/Pv7sz/4Mv/Vbv+WpNxq++Iu/GAcHB/jgBz+IL/7iLw72Hx0d4ejoKHqNLFiDSuYDl7Z6gvQDKJ1NcmjLMEjVxjkWyyiqBRJEyYbO1RtyLKYGTOSGZtmagnOt3Yyz8EHJobIf1jun4uyGvBDbXd+htFMWCrVmGWuoNsa95s6EY/JwNGneFNk+NkPGDJv/npHyw+G4VJ8EptxNvQOGKkoZjTfNDidNg64dErLVXZYPjoPmXMzUG43c8JwqFtniSeI2TZhDJZUsM4Y5dbfINLWkD06OspFJboJDxG/gpinPdC8SnaYS/RGOd73y3gx5cQC4nGBF/Z/sizLJTYyAmpmv2b1yckHNVQWLCc7111+P66+/PnncbbfdhitXruAP/uAP8IxnPAMA8O53vxtXrlzBs571LPM8Ijcf/OAH8du//dv47M/+7OS9/tf/+l84OTnBTTfdlP9DIjiN/B1WJ3u1PTAXAZS4Zig3HNLvhtQb3mmnVhHnA1fMPMXVG2CMhCL15shbm0s3Tx3VW9Vx+rDeYdvWOGx2+Lut38gBxUxF5YuZqYBVTVABSpOhrYU5Pgl0/sOU5MyZTAB+TpzDeodPYrk8VrmIReVEznKmaJ4ok5un1AzHtG+f7zVCblIBHVkZmqeYbxTkqDccfIkaSW5SJiqgJ9PHO+AEzXBenrPxFCQziovJDj+GAjuAwWRlkRzFDD+X5KzWnT/5yU/GV33VV+HFL34x3vWud+Fd73oXXvziF+N5z3ueF0H1pCc9CW95y1sAANvtFv/sn/0zvOc978Eb3/hG7HY7PPjgg3jwwQdxfNxHQvzf//t/8YpXvALvec978Od//ud461vfim/8xm/E0572NHzpl37p7HKvSW6sJRMIPPFYzOG4pEM2zVI8csooV8wWHwsLB8JMnM7GLPxuyAR1aZD75QeAU3m4qao/t/WuDfi+B1Y+iLUX45uLxfLG5Kg31uzU+liI7Ft7AcRYPdWWLgGgLlSbm6pBa3+WOjsJQr3hkVOaesOVAW760D5aDhW1CFodySDCU+puVBmI1auM8miLjWabqzJMW1PVG239PZ7stP8epsrQ8oZxM5eXg6kEtfjLwOuhl0Fb6xeU7fw8964tc1VOuosCrDpffeMb34inPvWpuP3223H77bfjC7/wC/FzP/dz3jEf+MAHcOXKFQDAX/zFX+CXf/mX8Rd/8Rf4h//wH+Kmm25yH4q8Ojw8xG/+5m/iK7/yK/HEJz4R3/M934Pbb78db3/729E060ZKzIFk30CYcAwoy64K5JMdjdwEuW/acVaXY66S6eJlY5fmKVJvAJ+0EJE5qHbB77HUK97ox22jckRYLYGZQHZEgcQ+FSMkBo8cv4bYMYnzpxKdknO2bY2TXYPjXcNMsWM0IfnhSPB2R/uvdocuRFzCMmmNa1KNymXOQFyKwO9GJIpL5VDhmcCtRHFTUTIglUX0pImFhazjc8h8gtBo6o1FIgmStHBSQ5DbeE4wulcJ9LxK9nFqvxGbCBlEx7v3HkjOqmtRfdZnfRbe8IY3RI/purHgj3vc47zvGm6++Wb87u/+7iLlC8qyx5n9XGmccnHkwFqKARjVG4vQtLthwcACO3XddKFzsVt+oVdhNHJD4P+fdE1IctoNjuGbqXZdPckPZwk7b6relEQLTS1XrnpjonQQjvnfJCT/tcwb212DpjCfkoymKgFFUh3VW/wtjoI1qYBhwJviD0bnG+oNHzypncnFayWOaVOTv+bcvs1UJgE3yE1sAhONBlvALJWz3hTBUm80cgNAzQt2zPwRZUQVmal2bR8ht2v1ulfyzqLJJRNkM1ihfvgb5CJb2Vx1sdjmgOQgsKK5mYgKdbLkjzOVBFnOxV5yP1bxZWeXkxXWgjX7k+oNgIDc8OywZD44qHbG6s1bbNu+UdPvmuKHMwdTIgrUBronX5ycwaMY1vPMIDlAXjj5HFCo+ElXg8IMnLlKEBtSceZMPGKomg4dQj8ED7HMxXTIMIHgZg8ZamylpTjeNS7E2L9p73BMDqpByLgYqFLIGYyy288EcsP3m6uml/YFmeYoQFdvrDX4CDLhKTd7nnS1mQxVC7TgZZMJT3OgmrQZuZGkRnvWoPsWkBwPBskpwZ4F8rOJfTsVW6Goc5dr0GBFTAG67w13MC6NpuD+N9I8BSBQbwA79f1RfeK2HVQ7XKpPPLMWAE+yleHiBL7YJ4fnBLeA38AayL1HsXqTcNq0PiYK5Xzvfis9x11XuzrP/XAs5/u5be+I+YWpuXBmmqno+ZN6A4yqwKWDEy9B3OjMr4QZM5XHSvlfYtaNRtJETLTa9pTvTS65kZMsdUHRUijtRLtnTL1x34W/lszmTglOObRt7nxFrQvea475OQbDD8ta62xcCoSFqivmKomsHDsFv+NcE5yuXifvhoSfZ8bvSC1HR5pNLkF6AiVHUW+4/40Ja6bOBsCYeQrw1Zs+1X1PZMj/RvrhSOJDOPCcloewc9bQZeRCMh9EAebUGfPciXkninPeJMhN9BoxojPRL4fKNoXo8FkpN7HmpjugnDiyDZaA109rhj3VD0wuy+AGjpqHhLeeucNbs632o3MCvw02KeC+OLLsWU66EeRmbFah3MsnLG3woe32OfOIZ4poxdQbzTQlyQ0AlxuMf2i/zAcGDGtXkQN5ZuZ8C4F6I8gNX9y1/x6+A0n4rGc+xR+nBOea4Jw2LPJSSmo0R1xNsQFGcuNnLg6rwRJr0kjnYqneHAkT1aXqxDlvcqJDBIhUHAIfUA69ziTRaWk246lqyQSsSqpzfW+E7F8yszWPn6Ly8OsqRCdFfGQ9lSRdOhr325romm1AWcK/WBqHInKTME9p6o1cuPa6zUmweC0nO5LkaCqO6WxsvMOSVeXN35ap3khyEwMfaOW5Hib0F/x6JeoNgZMbAk962h/jR1HRPk5yAD3QwipPCfQUICGxkaqNtp/O9cxcEadjvyDTf8MFwcnByukerrUH6iwypeKkfAWsnDcEy9FwSsI/7vxI5ikCORdL9Yb/Bh6Vwv/XfqMfTeDfS659FUOparDmmjsAilWcWeoN7Zuhaq1BcoCR6JSuIG35jhHZP+nqKLE56Rpc7eKZjlOQmWZno+6CGTO1L2eaEmoAKQKHykApSY5Ucc4U2KAqoSk0coC1jlWvmalMptpLrnrDwTO6A3lRVOO5o09PTMHORWju9idCnLzw9c/k5+Bg561Yz01W3jUNkpPKdJyLC4IzEbmDy84Lzx47V22GaBGdOdBCwi31ZrtrJjkYW+viUOfLwdUbjdwQUplltcSAWj6ckllMjEisobrMCadNEaRc9cbbrMzEYgMHECE5M0xWS4BCxQma7w1fumFuu5ODDxF9DTnv3TMNQE+eaS1cy/NM0XfLfAboeVRKMceXKvdcLRTbqpfmat4xf5wCU9xS6g1/V/0+X63hig7tt/KB8fuOv3Fa9IJUbzRyQ2SGL+7KPwebnVMDm4Odv3zOFH+cCSTnIooqhRXUm2vtxlVaK1IqN4LqiKWOP6q3ONmFeTs0cuNlLo7430RDxMVgxf1vCNI8BdjKE1/4EOhJDs9PclDtgDpM2gbALdugwQuXrLto5IQM6Z68rhO/ZmmKfqNPmkKyctWbXCJYN11gFjJDiWPPmu69QERb19YAG5j75QiGBVnbxg0Ax7TUh1u0ddloqVnJ/iLmqZjzvly4lpcjJwqH2oxM98+jD9VInEQ7mgzF34PK0H/PN8HwyB6+BIVWh4P7K8hVbyRi6g0htWQNbZN936bejVGkgx9OLC0IYIeKx9QbbrKsm568APCWByHwsYQi84aroMUQnQcWIEWRVUhEGALF4/GFgrMitoJEAPrgLMFnk7nOxuZK3BGkTFElfjhypsCjp4AwNJzAfW5i0JyNrQFFmzWpzpIRlKg2qRlozIl278uCKLJ/qa1eNwNE/HJSZqsFFZ2tQuZTHT5hafU0iUj2WELKeV8uXMv9OKRzKm8v3LRBkKbeEj8cYJqKM1X54aTP+shjs0xVBfe2vltrTgG6esOjQ3k/KTO6A2FWd8APtOBoaltdAiK+Xop6w3+rI9rD7zza+FF7R5sdjjbDvoMtDg+2QSJJTcnxyhTzRytf0u0CKlb2vSEz1UnXuA/HnNTxGjT1xsKUBfek/w1P7geE5imL3PDvl6qTYL8kc1ba8qgdOiGP5qA0+if72CmtMmaeivgRaP4KVkREcMtckxWVYQH/HAuhs7HwP2MmKvLF4ZmNqa0tnaqhcc/Q+P3Ku5bRU9y3DfCd9w/FQEkDoOacCtCq5/66RxQy3pd3z4Sbo9CZN+VjF9tfSupTfkCxe1nqDfe7IXBCwwMtJNEhSFXugII6Ck2Nap/ByDU3Sx0ebD1i0wz+YJcOTjyCebTZRrNly+VA1pgEXhCclcA73HEmOXaelpJDHaym4qwBmt2S/82urfOUm8ggSTNMAJ55ijsXE1LmqhRivgVrYiox0s5bS8VJrRiuhdua1zL2F5EcpQzZ+zJApIZL5Me7xsuHc5yhoJKjcUkkFWFpsxcPDe+v76deoPbFiQ2AgOzQtlgETnGYsTboL0H4I87FGrnR/D8IUsVYSsWR5aFyaGtOEbSw8H67n9UdCDO7c6Ij84G5a0feW+w3xvqfXpn1V63nv+/SwYm3NMilzQkubU68tc9iJEdzOk6GjhfgguBYWFG94URH60T5LHKpGaVUb/y8IVNd1PX8NwDcDFMDqTfaulO5K6erEQVKqHjOGi2lZGVucrqs8ye2zFRHkBOREj8/n+RMUnMWNFlJPzMJmfhvqbbG16MqQbDopDi/qTtcOjjxUi8cMiLDzVJ8oVru3yHLCeiJ4jT105Un8bsmtY9M9QYIyY16XgbJyVVxNMUzV70B9MzS0jQ1bvfTYxwZpIdfZw6pVhc6FeoN4EeHbQSxecTmZEwyOXwesTnBIwayc+ngBEebba/4MJJD188xVU3FhZPxyuCZhK/DSZAyPuY7c609iCS6Gyv1NUxTeHJ9E1LQ/G8APzy8nz3u1N8T+OXUJ7jaHpjLNGgoWrJhhoPkWpl3Z62JZZk5gMBp051izGi1Tj90Ku7P6cQz1pS/2FpGSSfkCe+obSugGR2Ne2Ifrt8D5sNAaqpGrq+2BybpLlkPjpBaCyhmnpLqjZYZnJeV/r/aHvj9zOBszZc6kesZATRw946h5npaE9/THPXGfVeOIfVpu+vXZKKo0KZuoz6HVj2NKR+aegPAU28a4fck1RsCV7e1rO5APx5QoAWB6u5hvcO2rnFcNcN7a1wZg+VqMsHVG/p9R0NkFF/vTC45EaYnccVBU/ft8bjdDM7rvdOxm2TXXb9Cg+EEXdpHnm8FR3tWlbF9YUgVh38AfTZpdaZHCeVjcfVGafSa/00OchUbQGv4+eapuSskA4XkZolZbk7rLGjBsYgUvj24heJUXHL+YmrORGXHiq7TsHdH4wi4eYqci2XiTK7Y0CApP1zNkTmkgLGtxvxwin1WCv3TNMTUEmeyG5xd6UNwaSISxD4WNm6ttxTrS+Szk2tOaY7FFrmRmd2lotNv01MANMJUl4OYeiMXdaXfdd3mJMicfd3mxN8+EL2DZoejTe+/0/vxGEkAByxhqjrfBAcYCc0eiM0xWxMnBm62ujY4QOb44VCn593T6Ny3Sv6bbP8bBjX/g+F/Q9DMU7lInaflwuGwBtqcSKgouMzKVZPIAC6vOdcXx14GIn5di8Box0lYCdc0TFrqgfYVdNZbhci7fFBt43xxNGdjYJ014ebChYYz9YYGSU5sCHIglJE4PAIHGM1UOSbdEmiZqVO+YRqkWgKkfYUkyckyV0cUJA2p0HCCzFoMhBM0ntk95quoZXXPQkE7kuoNz5z9CJctux0IzbHn2B5LKklO0ETAyAlfmqpSa5xl/+TsIy+QDZI6SR7lvgDSyTGVWXXKvrmOjt1AdKIr0Aq7tOxENP8bntwP0ELGj83byRw5iyC3seeQm4XuFZ5Xvi+W+yblh6CF2crzczLFTiI6E/1v+rrqk/Mda3ul4eLANMf+kgEnmrix7twAM86edcd9fk+a4R85xSYvBYMGizwUZ6/GSGpy2xG/h+XrwvOvWB9+XC36q2h5c8xj4nsqNJzA1RvAX5Ov3+9HkOZkdae6YDka56jX8t1IZ26ZOZsvCXJYb/HI5hiPbI69xJLXba7hus2xU3SI5JDjMU/GapVxri/OBcFZGdTRakoKX+WYf4D8tXBIDckOFS81T8nkXgyaEySgOy3mgkhOiekqF6UyO5AgN6UDsua4maPiLNBKrUGDQ1uTyCI7OSSHjptEdDIhSbjlZyHbnxVRpRFoqy3mkKASM4HMXkznN3UnQsPDZU8sXz0O2aasTLhZmErYJ1xHU29yIElOv238raXOxvwceS0N3LlY872R70Nbk4++82NS7/ogo/+NLYMiA0c4weYm0us214IUBZ/WHLucS3JhUCI5lzYUTt55UVVJFWdCP3jhZLwSegdMv6LRDPKkq4GIgzE5H19tD0BZHXMq7VG9xd86F2YfqYiSpcDLyR2MU7DUG2r02mCSkzSxBEv62ETPWyIDrHiFOeYpaxaaI93TMbuIY7HmfKzdOysT8oznxB2NgbFOcidjGnioLQLTkmWuBU4MpXMxV2+s2f94odE59aRrcDQ4F2uQjsY9mWjRFvgwLQVNWSTI7LkWQesdzbvpUaJGWbx9wSRPdy4myMgpen+xZWtoe25Wd6B/JtcSv0v/QaPCS8/8aLNzBJuUG239rPHeu7GOUfbwAdt6SOXAlgu6drLBroWfcZ5hTgDGhYKzB3jrPbXhTFKDrLi5CwGOCb3sEFDv/qWdV6Di2AtscsQUpphpimNRB1DuJzP13Ln3HrBk2LhmntLUG4vcaLlEOILjDb+cmCyenQl5xnPm9TonsSUwtrlrC/vipBSCWHg4zZ6BcdkTANG8KSlojsYxUPmjv2NuO0qYpwB4kUopcsPBc7C4ayVMThosxYebpzg052KOIE1GpTsYx8gPTSB7s5Afjp6bz0gzTzlyzdUbp8Rc80xtj9xcC7Jnkx+OzL1EPjncH4fusYaKc6HgTMCUhTY5YyXIMHEiO7zTIRXnCGPllh3vpeoEn8Al9/0wMksbyzYyYnIwzkakM+MLbMocDVb+mzXAZ6FTwyRVpDpyuX9CGLQ5Y1lxOqIRGr6NL8IqQ27tEPHxwcdUHV8F6hZb82jXVtgMt911dR+uquCgySPYMcSUkSKIGTQ3+3LnYm6a4lnBOa52B2ym708O5Lp1tJbbQbMr6wtk2Rden0oNYhDkRlMgd22NTa2vy8TrL1+jKvf+1j3dvqoVS2D4//t5ivz+MJb49KRrxvdc5xH3ps5U4BjBlAu7HjQ7fNrBsRe998jNtaH8W9WdgFJ9OJCSswE+uT3EYbPD8Y7M4PWoFtbVoirOhYKzMraegzFFcvQRHO6YIbKDfwfKTDCpUPEclK09Fc7CPRmZNWqrAQC6epOr6OwdUyJ9Mv10VsmvY8xYpXqT49egSfHe/kQkVsmSD6qSo8z2JaSjMZBnjqVoKiDf9y0GMinxGX3WekDc94H8H1j0FHcoTZEbIJzxh5E5+eHO/N3lJvyLInKu1q/I+hcjN3z7eNw8FUd1pHeBFfHMxYC+LAP3pclJfJqbBb5R3qvpB6fUv/43tUPem62ry9dtjj23g8ubv3NE+9Oaq/i05qqnOn16cxWXN3+HR26u4ZGba/g054Q8+uPIpSViKs4Uh+MLgrMHbJVw1fQ5PsnJNc/EJGc5M9PMU9ogEYPMYNyXIU8y2ReRkZ3gkqaQbEdj7bjIeSWhkG4F9EzzFGCTGy2niCtuguQA04lOkuQA2e9J1vNdWwez3W1bB21xiazGU5cOoXcnzVMA3AAT82nTTBp8Hw2m+npGeebsLMwwVcm6apmCCFT/aFkA/uH7tf4ox+9Mq8damSRk5mJeJ2TeG42gSmITIzmX6hOP+BKmhvvL3Ev+qvW9EzFPTfBpzVVHtPl3+vCQdl5OMndRVJUMG18K55rgTMk7kjzHkGhlB8vzcKRQmp+Dhw1KrOlgDPidpAxhfVhjZrRPDtbKkqypNxIaqdG25ZAcuk88gipcyFMeX0JyJClv2yog8L2ZqvGITWp9qjVy4iRzLgkiKqNxLvHBwzD5pvw2clC6YKPDgu1Cq188FNsKvpAkx50rVBwAQR3kf3PL5PYZuW+05TI005SZNTszj9iciDiZmoA7t3Nycnnzd56CSPWQL57MkxNykkP+QjysfZKKk/vTyp/GpxaWXOjQpZduK1UFOd412LZ1YKuXYeJEerQZpRZNZM3cLEj/m6UgI6iA5SJTTjX52pLkZmEVJze5XyobbAyS6GgLGqZWbpYfWTZtkOn3zWufu7bKNlON/+eHihOyzMkF/ls0q5XmKUBPBKet4zbV161oFW7rN6XMubnXB1iosr/vgClc8uPtr8cFH7PuZ6RAsMLVLfOUKycjPTKs3wsH90jMsfto+4EZk0YiNIZ5isAzZ2tL7mgkOjdvDwBTxQGWyTwPXBAcAPkkp5QMtbtqcOAtJxGc5MTMVLIS5YZlpxBN8jeAKqGMoJqKkg55EYdOjpKO2jiWIoHkJ/teS8AwTwF6By3JjUyUxpGj5qTIDi9LLJfOkiQHGM1UXMWRxIdPLvaBwLSI0P8mZp7ah8N+abp//+SELxU7JidPE+CrN4Dub8K3x9JraDlxSqCRLgCqeUomZeTgy2kAoek+lQBVXnfMSp1WozhU85RITUDljCmEchtXcfrvcRVHLo7af5mm4lwQnAFzlZyUh/eurV0klUwXb4F3tqVmqhJfgNKoCe5grA1muetQaYhlOt47MslJdAkCa7/hNFuyfAPfZ5k9cvKJ9NtDUqNtS5Ecfh/rI8vnk5lWnTEvQXI05PrE5V9vIkHK+H3cPGXO/E+jzaTUmlw/NX5aoj/21GLyE2EfwCc5VO80Z+MS5EZPSfMUh2VelIqNBuv9WmtSTQUn13ylek19skLaNQLOfXEs1Gx80daoKsEFwWHoms4cUJYyZWkd6klXBx8CJzkxM9Wl6sRFUvFF9CyyIf0Tsh2LE5VtDrnJBX8mpQOKnBX4O7t4h2yQm1ysNUgHiJinAE15iZdLEh2N5JRlmA3JjuWvU0pyyEl+19bODyeloC6uBs4EX56B/G8oCkdbk0hiHyRnkbocURiBMXrKWgZBhmNzSJIDTPNPSZmnqJwl5imJnPfFiU8scq64/yUfl3qcsMq6R2qL5vsV8//i32U4vMzdw81UVA7N0bzUT/GC4CggosM/k67T1ti1NbY7Xw73o6r0AVqGkc+BJeMujTn3iUntOSnx6Zlq+YaAfJnWxExyk7zmRBUn+1ZNF80GC1AiufEjESM5wEh0tI8FjegAy5urtkNblGYqwmTlZUHwHCTS/4arAKoPTsbsfyro/SwZ3TIHMlMwME7m+KROIzkatDqXQiy53/h/3Dwl1Ru+PcevMrpPWYMrBa3u9ffxyXWgPA0Tay3aS4sQG7M4t6aZCkCo4kxo+xcEZ4842TVuTZxjkftGPV5RcnLCxflqwSWIJbzSMJqpOlW6nbvo51lGbKCNRRAtpeLoa1aliYDvHDz44RiERtsuSU5upE2K+CxBckrSGxDmRhXy9jj6zeVfU5uR1mzgJP8bPoMmRHNJ7dtUVVqvjbpaQjJGhUEJ1RYkh8xU5GxMZqqYyUmWJWftKYJmntJyF0lSEI+IC9/zUtGplu8XT0yoLZYcLJhshLTLNbRGK4Ou4sScjUtUnAuCsxJynHQJx+3G+3BoScho2QaCVam8LMLKQLREBJU1Q+AEa0lHSB7NEvNfWjQcPjD56KYUjdRYyxKo155IfmL+N1oHsTHkdAuS6KT8ckrAiQ5Xc1Kz6jXMfSddH+EYi5haIhFggGF2SnVloygPcgZtkZupKk6qHsxWQAtAaqNmntLUG0BXH0vMNbK+WeRag7X2lDRPWZirvHG3hGII/0DLPMXDwTWnaK4gBqRHRFbxkHEA+c7GU3y5io7+VMMefj2PpCJHY4rgsFQcjejQccEyDQa5GStPmZKTOwuORzrsV7nJdRSdlI01k9ykkJXIju9XzFRZJlOjfDIbbOhQnHYI9o8PSc6UjyujoubIQSdr/SoG6YejmakAvf5caw9M02juunBTMDrD+v43hKWjpqRjao6Z2TMdeDvKSRBfuboE9Gxi6iMwmK6ace2jvDLFHWAtnyC59hRBqjgaSUi915iKMx5TQOwUsmClJujJiK8camHsvHwy2kpDytkYiNS1TJxvggOs8wQWyi1DRIeUihwzFV+ygRrZGDYYVqYp687EOiOts1ljIFg0nNdqPAlyk8raG1wuM2fNmuChral093JfzFw1pzz0kWoOkCY5FnKJOpGbs+Ro7NafYv433DylDXaaP0eOmcryPZqc9iGnTqxQ3zUfHCuDsZUTJ4Yp5in+vzTREEKCcBJ8gnNWMD9apI6vNbVWWoJiM1XBu7sgOCvCWhdH5uGgtahinaxmpuLfAT9sbxUVZY8D8RQcl66MHkMGuZHg2XlzliRI5cjJsTVbyzPkEgHZcct09/IYjeRYgwXPqZP6AKOpaymSMxfUnrS2xjHFbKW9W15fLLXBip4B9Nly7oCUctIvQqyfiOyLKSfcPBUuhWD74NB6UFq0lfTDmVqvNPMUMCo3KUISyymTOkci1u+bfnmUvbkOcy9x89RpgL+bKSHjFwQHWPQpdMrMkSKpAG3JhnA9HIvsWASId7C+Y1dbZIeei0XWsolAU4IWm3kX+MKkMvHm7jPvvxAozNJKd0/QSA1tJ0iSYxGdaKLASJQWJzk8T8kUklPqLL8GTKIgCezQaTtTTeYgY6k3+4ZpOpCDkTI4pcy9VqSS21/7/dtmULx4mgyOEjOVVh7N3zAVHh4zT832vUkQohwVjpsIrdQEXvkLy51S71N+OFS2qTj9nuCsYOknMZiptM5Wdn6WD8lWSQh4rd1kVJoTr2KS9LcErA5NdhzksEnIXSw0B7rf0rhtitkNgOnExn9zLPtuDGtm6R0vbF9L+g4AYf4ODdrihYRUB5oKO5f7pZojSU7wmzI6vpQfDrD++mxLo9REIQfCOebd4lDxjFl3zurhBKneAPBIDd82/u8ToakLUQKI+t+M/4+DNsGKdpqi3qwNmZpAJpbcF5ZYuuHh1bLXxpSnkXkOdzQGQjPVcdsknWVlx2Ql/CNoC7wtaVuf01FcbQ9wtTs09+fMEvaRu2QKubFDxDMjMzLNVHMW55TvzlrPh1BCcuw8OvHsxiPx0klOKnx8LubUpyWUxNgsWvrfEE5TvSHsK4GlTNp3KIjNoUseNxIMj3gItUVbmZwjlq2dIyc8fA6464FE8bUVpbqqR78Xauc8NYErB1NvStc/jOFAmPjID4cwx0y1KsH52Mc+hjvuuAOXL1/G5cuXcccdd+DjH/949Jxv//ZvR1VV3ufWW2/1jrl27Rr+5b/8l7j++uvxyEc+El/3dV+Hv/iLvygun55LpPgy9vUHFUf64kgzFYdGdHjnSWGsMUVEeujTwGHN1qfkD+HIJU1rRp+sgfjyC/5v1sLEU+stBVhgoCgZbPhihRrM7QrJsU1PeetTcZLDBx+ZVr+E5OTU67l+W0vVZfkb5OAN5M3yc6JtzgJkPc3xv7HASQ3fNv7P1BtRp4Os3OSPMoE4p8LDpXlKqjeH1Q6HEwiDRjLk74xNxgDm9zaYRuXCrlMTDlq4FFgZtooy1wbvpxSrEpxv/uZvxvve9z687W1vw9ve9ja8733vwx133JE876u+6qvwwAMPuM9b3/pWb/9LX/pSvOUtb8HP//zP4/d///fxiU98As973vOw2y00m1nwqfBOlmc05mGqmkRORIfIDpmpSuTlzZDHAFjPP2ZuBSyF9vv584uRxzmQ6yXx7anO0NqfGwKZVGoUB2juf5NKd8/B1/Kh/drChekwcn8/9/OxPlLNkSQn5gxaGtG2BJYyu/I0+XwWDYxRLIvM3gVieaRSWOpZywkBh5UpmKs3h55S0rpPf+yoYHMzlSRL0hxmrWguwfPfaJO8HBVHI61TiM6cBZbp98r8QSXXXEpJ5H44rnxD0MSkBVEXKZWCP/mTP8Hb3vY2vOtd78Izn/lMAMBP//RP47bbbsMHPvABPPGJTzTPPTo6wo033qjuu3LlCn7mZ34GP/dzP4fnPOc5AIA3vOENuPnmm/H2t78dX/mVX1lUzq7pkgtlLoFdW6MeKs6urd0LPG4bHNY7N0jnqCHX2g0u1Y05e9zLKsOZHVxf1uXtt7mds5bMsG46tLsKVd2ZCRkt09RU8wjdk66RdIStu2i6gTnmKSBcrFCCthEJb6oWu66vt5xENnXrzK4a4ZGDiaUK8WvTNfu2UHtrSfH7ZT1Hht7Rfwc+g2k2fXmO26Yo8lBre5MIg2b6ZQ7GDxfE2pJ2rAVOyCWsSZqWS+Wgat372NQ7N1E8bHb4u22Npm6DiSWvWxqoTCX+N1Z4ePDblP76sNrhuEAhPKq3+NtdaPZv6hYtn/iJvkWuPyVTE1BiSW6eKkHOUjuA/56oXE1dY8tuVzcd2uoMmKjuueceXL582ZEbALj11ltx+fJlvPOd74ye+zu/8zt49KMfjSc84Ql48YtfjA9/+MNu33vf+16cnJzg9ttvd9se85jH4ClPeYp53WvXruGhhx7yPknUWOzp8E6YOmvu5HjsmaDqoOHJ1ce5yeqks4kOLdnAcwwAE/KXFB4/1x+hpAF5zsVLhLiinNyUrJwtrxWoOFPNVJnncbk/tlghQdsfkBbltwah5opfj3eNYT9Xc4BRyckxVcVQmrU7NijJDpurilvRlnNR0cAemP/WD9HNTZSZIvQ55tGYT1sMuUESPEMuwHKBKeR1ij9irO+0Vg8naOapNZD1u4bIPdl+xvD60/ft4pjq77kawXnwwQfx6Ec/Otj+6Ec/Gg8++KB53nOf+1y88Y1vxG/91m/hP/7H/4j/8T/+B/7xP/7HuHbtmrvu4eEhPvMzP9M774YbbjCv+8pXvtL5AV2+fBk333yztz+eGdbeZe3TQsWBsZP1Fts0sqkS0bE6yZQfDnXQ0tHYDRxiwFgTV9vT8buRK6Zr4NlwY5lxY+QmBmsxyVKYio3If6NBC7eVuTu0j9uv+OpoYeV8H78PnUfmL+vDj5ckh35HTt4S2WlPjqpjWKX+shBdCZnBGJhnBlhjQLXWWrPbUN77kpCh2Nw8xcOxObHQ1uNzTuxsOxEWqy/MdTDWCJT0vwn2M9+b04Z0MAbm1RmZG0q2n378sg1InkP5jHDx4pb/Qz/0Q4ETsPy85z3vAQBUVTjQd12nbie84AUvwNd8zdfgKU95Cr72a78Wv/Zrv4Y//dM/xa/+6q9GyxW77p133okrV664z/3331/wi5Gt5nimrrbyZEAyT9Cgq3W6sZkUV3E45GySO4RpjsYSPDplCcRmrRohi0VSlSDmLJozwOUuwcAjKwKFprZX0NY6yqT6IB0xm859zFOYk2Qs3BZI+2WlSA4wDj5qgkBBbOia5kcQHU5y+MxZLpK4z1WutY4aKAy9jtQ1MhMQXLI4ZXFGC2srPjnPWxKdlLqjtbWYKYiDExsZMt4r2OMkz/LDmQJNVbLy35xl8D6CxokjRnKWKn+uz5qf2yh0li4Zr4p9cL77u78b3/RN3xQ95nGPexz+6I/+CH/1V38V7Pvrv/5r3HDDDdn3u+mmm/DYxz4WH/zgBwEAN954I46Pj/Gxj33MU3E+/OEP41nPepZ6jaOjIxwdHQXb+UCR5YtTA2jZ/xno2gotgKb2/XCAwSG2AY53zOdh8Mnh2Dp/hB7X2g027Q6fHplUUqWkmc4ncYjDZoer20LHyBkqz7V2E0id3gDRzks7HjOF7WYulxEzJwXHGs+ItpNyRzZ+7o9D1/V8ghL+NynEBiBrYUCNAG/bGofCf+F41zh/HH5NCU5u5L20+9G9CHQfeK+4BrDzlLnx2ZX540xBroqTa/IB/IHfXwS1bFC5VJ0Ek51L1TGuIr+9T4kok3VZIpU9l8MPDRZRZQYhkeRG/p9DPDfNWKc0n7LcBH+p/DeaeWpt9WbT7KLpD7iDsQsmcCRtnu+k7qeWX8cOmp17F3XdYTdhHlNMcK6//npcf/31yeNuu+02XLlyBX/wB3+AZzzjGQCAd7/73bhy5YpJRDR89KMfxf3334+bbroJAHDLLbfg4OAAd999N57//OcDAB544AH88R//Mf79v//3pT/HQzbJmYl+8K29Bn28azySw0EdP20/qftCkP8Nlz+JMPxNe8k1JPLD+SQOXKXZrbD2zvGuYWXd4KA5xrZtZnn4a5CdVsmAQoh1zClyE+v41OvVXeADkhoYloSamIwlTItlipUkx11TkBy+nd8DiBMbujZt50RHkpxtWw8D3w5bNAB8p9C5Sg4NULmkey2zq4xikQONFQKukZylkHLAnYog3YLhYEzQzFMpMnhQtTgezt22ddTRmBCkQjCyescyxftZ5dMh/7RtrXeYA/o9VhtYQs0pXdqkqTtsd4qzdAZWm/I8+clPxld91VfhxS9+Md71rnfhXe96F1784hfjec97nhdB9aQnPQlvectbAACf+MQn8H3f932455578Od//uf4nd/5HXzt134trr/+evyTf/JPAACXL1/GC1/4Qvzrf/2v8Zu/+Zu499578S3f8i146lOf6qKqSqCt2rwGaHbJBzqZUTU2g5INkasX1CD8VO6jDw5hqWUbtIEkFp59rd1MWq9nKuZkprWcigkWucldLVt2nKVmqux9rFwSMuGZt0/438jQUSuEXDoRS3JD+W00Px/+3TtOlN0iYkuYWEvbBpfb5zrVW4nWzgpSRH7N8PxUKDYH+eHw9P9HzsmYtvl+OPyaqd8ZXzbCJ6Up/xvAV29OI1OwU29YBBXgjxnONDoxr5KV9yylrFmKXWldW1XTfeMb34inPvWpuP3223H77bfjC7/wC/FzP/dz3jEf+MAHcOXKFQBA0zT4n//zf+Lrv/7r8YQnPAHf9m3fhic84Qm455578Omf/ununB/7sR/DN3zDN+D5z38+vvRLvxTXXXcdfuVXfgVNc7YTyHE/HGlCiZEcPns+Zn44Ux2NgWkrQc/NWHrSNYsu2aBhziwzla9GIzeSyEjESI6WA2RSNFXCREbgg4VcrDBwLM4kOYCeS4cfL8/vt+2CDy8LJzncH0dGVaUcjrlKJpdsAIZ8VAZBsfwPLOWmKBNyNBpnfJ5TTbiaH94UrJHnKnddsZiDL2V4tiKXpMMxAEg/HA5LoYlB1nM9XH2557ekYlgZEy4e3r5k2a8pgSbbtkkSHR4UQyhRalfNdf9Zn/VZeMMb3hA9puvGyv2IRzwCv/7rv5687qVLl/DqV78ar371q2eXEehVHB75tGhunLZCBwBDI+Z+OL3k3uJk1/Smo65GU7UeyeGDCffFudZu8Ej0kWXSTMVxqT7B3+JocLSLh+ouhVROkVh5JS7VJ9FQ+LWS/MnO1yI3BE4YJXnlNv5JyPDJkSuIB0sxKJ2CXKNn3L7DdvAFo8FfmqvIH4dD87eR5Eam1iccM3+treKHdrwD0CwTEUWgfFRzlE3LXJpVF+swOaG/HMGy5OKoPpndRizTaq7JNZfMaMpiv30XqFue4jAQOlKMqS6TmYqDzPVN3RX77MX6Uc20yP1vZN+nfZ9qpiJ/yxKQcriPEPF9KvnAygrOwxlLm6rkTJKQ6rC1MHKSxGMLb/IU26F6E7JiQm52XQ7LJBTLvJyzZENMFp26bpCWhVhmwU2RG67aaCtnaytqu1BHZeY0yXckEnIrBxE+WJDvDd+mKSmx2a53LZrFFpAbK7U+fVILJZaoOKWSNp/BauADD++s52YFljlwNDVgLEPZADZ1Jh5TebU6m8omncz4zTJuc+RmYZfrJmm+f5YKGSuTRI6Dccz/BhjNU0ubpnL6xTEtRvgO55hGc+uZHBNOuhrH7SZp6p2S1uR8ExwexriGHVnMCniUhzRT0cxKc9zkM2UiEzLhH+94JTHwEv41/uCzNIiIyfWzLD8cWnTzrKxPZa7ZYpAYgr3ApK74uOtGUtVnIcNvwP0vTUl1SH45JMnRTE1AfLDQyA3gp9W3UuzLMHUZOr4EpiqavL7GZPbSBINAvj9QELJ+is6pHHLiEFuTLVXntYg/YDS7az4vErLueddXnJyD8xXSlYsYgbH2yffokeru0OUVo0ninL6TE+u1xoS58PJ3FZbxfBOcBNZwOJYy7ratnYrDSY5KdNpQzYn54VyqhlTbovHLdT4IU5SEUtPLWSAy1u/Ukvlp0RQxYqMtPCk7x9RMZHZWY1ZWgjbrJHD1RFuRuT8mTnJomxaZxQcYTmS4U+iRID6c5GhOx6W+ODngCwzOwdwFPGOwJH6N3GhBB8C8dYs41s4/ZNVZqW7x/k3mARvNRTyc3HAUFkQ9CCQw8vJYDsbayu8lik3JMg3X2k1ZLiaBFNEn1XBO/31VOOZPLW+JGHFBcBIqzhokZ9fWzuGRbwuOY0RHmqqOMys0NbIjZVbsJU9b4HdSWUlloqUlJDSHs6XuvQRipimCJDbW4pOS5FjOjLOiqbTDuXKjdMop9aY/Ju6cKSOhCJwESXIDIIh04QMEfTSSA4wJAKO/PbMuNwoJy4EVQSUldlJoc3xTZCTLwwVTSI48p2RWLlUYSW60/zm8rO5KRmNCahJS+p64qj4l982+JobkvD03B46EVf45pt0cXBAcgaVITieyGndtpXZ0PKID6FUc6QgYi7CSJEfzWyE/HO5o3IjBeU1cG+yrsqz7ktRTWYSrus3yu+m3d1Fio61FpSk52rpK8v8AdTd+MiH9byS4qmKZi/rrxBQgO+KKXwOASmz491QWWjJVxXxxgOmEPaZuSL8xXp/J96E0RUFV68Q4XZa8GbXsD+aoU0uoZBxevRf5b2JKCf2vZds9GhaIJHAV58Awd80xzcgIKsvBeCqW7iNJHSazISfWS6URkeATAm1SnuM3NNVEeEFwgKwBY7KSo/jhaESHm6qAeBRGv0ZV4/xwyA7LOzveqJZm4xqkAmWF3qbWKNk3cmagnNwQJLEpSf1eFJKakQeHIqjIwVhb3oDn/eDJ0mKQJKf/fzRBhQtD6iGndC0ZritNVFr+En7fTeHzLh2IS0w3S0SDxBz5tbLMbStrLRmQq+LwiYQFOYHQSLn1nmQ/x3+vH72nm24JyXw4BX5b1jOf4lyc61w+Jw/YPqBNdvtliBps2z5tw1Jm3rP9JNaGFYlidIxLmqu4mWpnRFVxkhNTcazOltt8paPxQTMkVYsMFCVryRB8h2hLlkxXXmdHj0ZSzU+wxqGpN5LcSDs8X3+Jr8fkVuzOVHFkeebmHOIYk8eF79ryh9FmvZqpSnNQ1ExTGrkBhvDZ2vcP4Q6k0lTFfXEaRrI0FScG6Tclcak6cWWTg5Ssv0vK7Dmmw1JyNTcXTg5SxCXV1jRoDsbWAqSS3FgKikfWZ5oDUxFUqXDwXPD6tgSxthZ3lfVuiQivWF9f2m6mLBB9vgkOR+aDi5GcuQTIJR9LkZyB7XI/nFgCvYNqZ4ZNTkFqzR8/MeHGhYtzwsP9cFKNNjXzXILtx/LdaORGW1jSWohSkhxgnQRqHNHkg0K9sWfEXLkJTVXj9XTpXypERFrUbK/1iUd2JMmR94r5T6SgHc+XaShRPC0/uNJ8PTyLcS5KwsW1NhQLRZ8CUmi0TyliodhAXz/GuqI7GdN3fo68poRsNxRUsGHtfgqoHNL/xurbYg7GS6vepb9p6v2nmqckSlwqLghOoYqD/3973x4jyXWV/1VV98x47fFGjvFLcYJ5yDHxI14bzNq8Ey1YEGEFBQzGWCI4smSTBCt/8BDCQlGMFCEBgljYIiiBgCMIToIgThxIlkAwwQ6rWCYYkwC2+MVxiJLd9eKe7umu3x/Vp+rcc8+591Z1z87Odn3SaGa6q6tuVd/Hd8/5zjlYjMjMppnqpvIyG7dIYMc7zqhcq0MJAX8A0U6Y63AkK04x7fP2U9snSjh7qK0xyIlK3ktbn3EsU3H9P9PdWOQGaEiNs9iK1yTJSWmb184FrDl0Xa0uFLesAHo+HDqu+qzrqtJEx3Qe6xoEHmkiv1vN4mBVhNaea0qSw/q8CQUueUhuqEQDEfskwXvAkqeBtBipBQxju++uydyWocPRxlsbwm+1Xc4P67keOWaR9BhRDtWgCkVQSSs0fTexjZsMD5fvLSMzvLxn/pz8cPV2+ZeANLfooqVOQugJDod0GbQgOVHSMxcac2jRVJwYaCSHRyqRDodgCdI2soknNCZQBy/m+o02u62uWWVH3Hozb3Pq4CELQNeEVN7kamT/lRFTGrkhUMI7nviOHyN3SH5iQCtsPU1/470lzx+p5SMJDUfIVQW4UVkclvUGCIfR0g9ZHKUVR0vUJsXGFqyFtC1RrlIzuALjbSWFw06CNjBtFh0rbURK8jvtuS4qNvaJqFIU1gjFJqwH+lLXdhBSSFeMKLZ1D8o5PMU91cVtFRK3p6CrFcciMzs1bnqCAwR3yDGSs1xdzrxWVYTkcJcMTa6hDjfMpvPogm7+5xQ9CK/t07StmIvHGjcVx7LqUrX1pWvJxyzdjSZsLbKZSmh4e9YUIsQ1I4DuU140t4gb+qyfi+ti+AStiX3pdf5Z/l71mm4F0aw3ktxQ1IuMfuFtkonaNDeVRMrkzYkntTW2MC5bf5NCErZqMuNacdxMyroVp9LhjRfS4YQW+jYkJ5bpWCI0rrk2Sub7kdYr1xIsSJXRh7q4kK3q74uCk1iyJkp4KQpoI7wE9/2kLDoRagtb8/VAuqd2Qhy90gQny5TBpizmscRCKdYb7zNzN5VmxQFskkMgHY7EaDYMMnoSGgNoLTTWoGVqDVdF9xX0JyOSKihmDOhu6H1JbgiyOrbMvkuIWXGWAW3xkAJjbcfJSQ1/jR8vtTDauayMsWSF0ciN35aK6JAVR2sPf86hzMZWbqIYuqacD8FxSSd+99r5NVfVIouO/K5SIoS87L8JpMWruWW4x7V2WKHYEppbjvcx2V+7hkWnFNmsXhfXS0z0N1ZIzE5Hn54qWYyXbclZaYJjYkmLT6mEg1v5cAC/fINGcibTwrGS8DII6i6ORYJws78c3F10OBKciMU6qkbCtPbLha4toiGfMrRZcU1p5EaSGSm05Wb/kBXHaWssB04AfvI0PVIC8Csxc0ihL+BbZfiiKGtX8eO0hUiSG+6W0gSiMo+JFQXDNRQpmpZamzR3z3ECZoEsKVIwT2i7A+Xfd2yBsRa0mBWHgz9jjdRqSCXiVkmGFKtNWBBv52MCwhFLlq6rytQduOYOBAB0iUiqrDWu9YZDWsBl+R4NWp26QtHYaZBWnFB/06IPNWjtXWYm8JUnOKlp8XekVhWDW4BTJzkSsQk1XJ/Fj/pZRsI/KTQOuana1FFZdlE6wL1fab2RWgBJbgC3LAP/n+BkjQ5k4NXCm9WipyLRH0/aFbo/1+LkW2pkBBMAj+TIvzWSo51fIw+WboJALlXtPPS35aaSRU3bIuZe0Ij5tmJN3Z7luntgCZunZVtxJEL9NAUW2SFYFlO5oQi2kfWRkEtKYtFIUtOlFQgR1yPY/NcsDSXvcyGNTiqyiPUspU0hhFJ7aNhJ3drKExwHcjERWJjkKO4cclMRaFKUUVXVa8xCwko3EHHQSiBsZOP5jys0BubWhYTogRjkRB6LItFYOw3i0IBdpm9bIzdaUb3h3ALDyY1GZGQdpliysqoNuuh1UR2Ods5UE7Ss5QPAs/ZoUVX8/5DJXpIb0kzwH4LmqtKy0SYl/ct9awlfrAa5H6XHoUcv+TqCsZKJPIZUi4FLaPSq5hZii47WP0J9putmqF1ZBl1gHNq4ad9hbN4gN6cWSWURLy10XUPrxd4saVCofwPprtJFCvtK3Zdlxeny/AF/s76s0js9wUFAZ7ITJGcOChd3XvOirHTR8UTUpUoFCY25DodQ73xbRlLJthKkmZGsOAQ3tN2O/vJeyyetQ1xpcZMlFPh7lmtKagF4KYJQqDSdR77PJ9OTAbmIA26OGRmuXbdTkBz6vJ9heNv5TbBcYM65DaucdC/Q+dzrSk1Tt+dZfUfNuWVBRM0F1KZQYNdIQ3790OIW+x9wv0sZSTXMZur3ZD1Pq8ZZCmLjjsPaFHii9YClRPYjzXqzzNpfUmNmEbGU8HByT1nWm2UFaISwDJ3PThR5bYOVJjgqi5WkZgdIjgwXB6BacTgshlu5gZpIKossyAq7pMMhobGGtixftpGn3iY4Ye2BSujLQGxnbNWZkrs1qtMii0hqP/S+R3wEWTLbpDzzNlmN63sK5OyQsPKHhBYH15Kyrb5OnyN//Drzy/NcIKmiXk3EbIWLA82z0MTj0uqjFRjULCMWsZH9fNnYUoiWtOJouVJCaCO6tcZSKskJi/zL4KYC0PrVRCXk1v8S3M0ZapcGKwdODMvIZGz97aYtICLefXnfUohVamJWAF7kHn0foc2pHD+WJTRvcVsrTXBMJCworUmOQmrqt5jomKw4IVcVIUWMxReQ2KDvqlvg4BYmC22iTyq3xViZwNLCPQlxsbH/fh0OrpAbC35kkW/FsduwWHRKUmh04DnxUG1Ociwrjp82Pxzl4rfFXaAc91jmkyJ57pRwcQ65SNN3y+9DW4RGpU/EyT2l6W9SoJFWfh/aGLFIjrbTluRmmXWovDxLzDJq/Wifj21AYqUQCFJMrFlOQqHiXdHGkhyrIG6RUY3MOuSD/d0lK3AM8vqc5CxT97VTWHmCo4o5vYPaZzs2wQiLdFNZ5mzNVUUgES+PpLKiqQB4Cf8oo7Fc4FPcVFI/5LTZifRqxMbNa9qkHB4wcrcf2oWFdEV5XgazFXPdDZEbsspwK4L2Q+2qjolP4CkJ6jTI/hoKEZdt4e6pJkdSOMRWlvqQJEcLIW8Tzuu0W1kM5KLlRdO0cDVoRDAl6mM0G0bdU221AyFSysc1YYtZPUOuKw20UZCRVDyqSHPZLDOlgU+O0qwkhJD+Zpk1t7res2b9jAmOCRQeLt1Tmqic/72TmYD5taxs2inRe7uFlSc4JhJcVcByNTn138KKw0FVx1OEjFaKcFl4k6OLlkHLhQNUFia5u5U6HK0S+k6BExsgrLvRJnpZQFL70cJZa4I0P3cod8uyEYqg0mCFa3NI0ScnN1piP26JkWnqeQRXbIHSrudEqokkikA4zwoPEV8mFtXdaLBChLlbmnbVi44lV7elZ/fuYu21xp6EVj5DLffBNDapwuMUq0uI3KRqg7RrLwL6Ti3rDeFkZdJOgbS8a5GZEtIrIcdS2xD+nuAweFacRJLTBVyHE7LiBKOqyrwmETySipsz5W5Z26nLhH9d1PYyh08IWlZjidAuP7YgpRS1lFXCNXDrDVm9iNxYkCSnazKxOpQzwcIYy4HDn0co1X0oAkLmpNGg5c3REFqUuKWRu6mae/GtUj5J958Vfd+WxU+Kiwm0sFR1fwbzv93oKbJQ8sk5dSxosHbkciOg6SSWiWj9tA7aMPf8uuWUYCX4k25Lqx9Z36nXDrbpCEaOJWiDnDa0jKAikPXG0lRq2huZUTtlAzyY3zNZquU1Qrovym6cXmJn5yvaa1hpgiMXkeq1yKBtITpOse44lhvFisPRddKUOpZQwr9FIS1MRMA0EWbMtCr91pYiP5bVVqsQDNgTrHRN8dB6PtlqP02bulsFLFKZ4pZq7i28OFmp7gG/GKA2OcmdtZYJOcWdQGkMtPc4vKywxvNNiU4bzktpUARVqMgm/w242phF9Tdd0wFIoqO1UwORVB5JNWAaJK0siYQcS9IqKl/X3k+r8ZSWHViDRSz0Iq4pbuT03DzV9ReLnpJItd6kIi/8vFscW46VPSxuT8FuRlKtNMGx4CVXkzBIDic0QXITISpkxQm5qjh4JFUMtFPnizbX4fDJaNF8LJoeQYaLUyTVqLRLTHCha9edQL2Dl2THuEdObpp2zJIywPI8MUQiuZsqhBgxyfKyVVQVQdacskCTkVbxWFpxQu4D/hnunrIm+a67XYImNLa+cwmKoGoT0ZUCy3VrgY+XSR0l2ehw9KiuNPFnanVxa8GPaVJChKa+RjFVNxcSTrmTgFVQCuA1shzCsjd3IcSev5b/hhMeTXdFRIR+1wVfA+uFNqdzKUFbcbvsb/S+ZqVNEXqHJBNddFErT3A0K46HRJID+EQnBMtNFcJ0ljlWEkkiQnoW2r1xaBMMVRYnLFo1WINVpDDF5EmDpq04ki90Wq0pbr3h4MUpARI36z9WhJHTLsMkvmgUm1VFnBCLHtHS3Vef8zMLa+6oUFVy57zMeuNeL21x0iKpUkDfOelvZASVha05udieu4JDm4m2Sf7C11VqUTGiIxc/vrM25wBFaCwhx5J8xl1KGYQ+Y2nfrAR/0vXUltgAVxwlYgAARg5JREFUi1lYveeh9HmnNESLCCruniJoJJawTIGxjMKy3KGWlSmU9E/TQUnr6U4QzpUnOBxBV9Uy9DfGjk5zU7W14gBwIi5iGYGl0DiUD4eDnlEqIQNcN5Vsrwbe9jUxSNJzpoTzbtAxKa4pSW4Av/L2IohlNA4RzJiVLTUHjsxLY5nb5a6ZIC1Esagk040QWKjaRMlYuz232vo02taQOZ67XsdssyH7eauxYixYcqxYJKeLwFhGUnFYO+hUksOtNvI8oTEa0t90QZfPxeZDjSRZ7l4Loe+L62AkuZH9YVxbcZazpGvriBXBFwN3dVtuQu1Zp7pJQ1hpgpPVAk62i25DcpYoOgbaTYQcMgRbA9dbWIsUYHeemBWHKqLHLEyyrTySKgS58IUsBXJSdq01aZOrJDfcQqPpTehHWnEW2Sm2QWpUlktE9Fo+EprYt/p8XFQs3VOy38VIa6htKaHilphc+/55YVoNW7PBjuQakRiL8bHN3Af0A/iLUCwRW+1+ZDocCZ5bKKY70cgLvaa9B/jjT/Zba7FT7yefdHJtpmxK2mYaX0aIOrl8pLjYIjfULyxrOKFtVJ9cSyxxu6UZckPkx/Pf+tyRiq6Z31ea4KTgZJMcIM2KY4WLUyQVEHf3cKEx6XCApvijLNtgkZwuYbHaQpG6K/DDluMJ3zjR4ZOrZr2xwE3lsfBm143iWoaWgRDhlPegXVNbxFOzwbqaKJ/sxUCTnmXCpiyoqQkqU5/pYP59N3XFtpFSg4pHUAFp4uIUi2sIVJzWgkVy6pBxFi4eIoltS59YEWoWoYl9FvDrvRGsBH8hLRdBZtJtgxTLUgrUfE4tLYUx0qplMF4W5DUlyQnpvhwXotDhpLjxl4WVJzh5xIqTdpLuCf9Ih8PdVDFo0VQUKh4CZ9My4R+hrZArtd1aTpxJmbf2IbtJs/w6OtFdp7FzpPcs6w2vul1d210Y6b1QDSZCSv2bNn0wluQPcM3poSglSdq0zMLaebhlKyVpnoaUz3TZfcs+relvKNKQWwRi0SKae0qS/TbjOnwt30orSU7IfbAmyKKf66gZS1qdr0XmCPmZ0BhsjrH1NxLcrRrL4aShS36qkMA+te/L74lIAllG6o2qk+jRtd4QtHm0S9QtrSG8b8l6bFbRV4lQlXdLhwP4ZW0WERqvPMHh6Oyqsl7rAHJTtdHijB3LTrgmlbWrIR3OMqG5qbRwcaegnLJT4aZ1Qqz4IhA3a1rCYovcANrC0O2Zha1NzXshHY6mv1lm1tnYokKwSIfmnpLWG1lBfFGkToaulip8fS4wDhHyZQqM63OKnFFEdLjbitpYHR9Ppc9JnBQaS9G+VkEbiPczsphKF7EkN3R+qX+r2ubXn+qCRftXKPknR9f28Qgqj/QEyI2mYRzPqvl22ckmNaLj6YOU/iatsDEdjlacmNDFTbXSBKdQFo8dJTk74M7iHTmmD+CdzaosTgn/uoaLcx1O8DhvNzoM7gbcBdV1j4TEkYP57nPAdqFaXgsZFm6RG96GUJI8+mwwKaCyYwxFUuVFWf+kIMVSFMsG64bhhq041jkkYvlCQmhTx4yDBOVF1vQD67t1r6dniw3179DOORbiP55W55bpFLT7lgudlqckNaKShMYcsRpfksRohIZDLlA0BrXr8Lb4hV4bstx2gyHPFS+42c2NIrVBFdEPu0L539wFBNiWG5l0si2cIs+iP3G9F0EjOvS6BU2HQxbfWLRp/XdHK85KExxgF0iOAq26uAaZV2N7lqt/h2CJ8kiHwyHDxdu2T0LmZ2g7MGWyQkCvLB2yRDkDRVhvuGtKgu9A1h1Lkr/A+9alsNuKJtku4bfOeRLyvaRAumosl1PIrbVbmUsl7FpHjf7GWiSbCKWBKjCW2Ysl2uyitwPlV7jbQAqNyZrTLErtoly0HbWVOsJ6Lwa+ueDnAoz8RSL1ALAYKd5phMbVWiSbsl8R3q5gL4MzJLQ1IDUXk/ysvJZGdNok/qM5RBtvsXxhi1hxVp7gcOwmyQEat5TlpiLIsgjcBSSL82mo05gzzQi3bmiLwiJJ/7ibSrqntMETg7YgaUJjufviZnE6TjONSutNY/FyxXIy0+8iYtiurqUYAQ1NCLoIcqz+LT+Xcs+ae+pkLVROuQ72N68gnkL4JFHQrDfTsiEnXcJ0Q4Lk0CZAuq5SUkQA9L241g/KaCx1OFo0VeoiI4kN4LumAN09BaTpb4Cd60uW9SY0jpfRFrLeVElQbXLDk0EuS2Ac629yvk4l1Fq28pQ6cKH5PAUrTXDy+d0XJplZkOTEMiIvGXwQWB1O6kd4RmMnyVZLa4K2W6VIL+11idQdJ99tDpjlRaI2f5NrwphY+bm49UbLjxISMaZaLNrsQDQBvHnsgv3Lci3V5mWlto+850UyAWuZkwnLKMRqpdtfV5JfSpD+hkNbUNpYbELEdDwrnHQKFFElfwB3g0ApFwA4mcFDten0HbWrwyEUgjRqfVm6g/ln5RgMucBk/hsOTpZ55u2d0nLF2mghtT2jsimSqvX1ELnxrIpL0IFx4qSJ26kdnFDHMtEDdnoCSqeRbMVpEX210gQHaBaGHSE5odcTYYmNJWjC5ROxNsEB5IJwOxmv50STk1w021hxYgp+LXcPDfCUbMaNnojtzvPwxMknVkvUqEVBecX95pNprJxBm/Dbtgh9F6HJOVymwQ/f5q9zaMSvuUYsyZ+b5TWW6XWZ4BXEHX2RmHRH5bDOYEzYFtab0GLStkwDQeaNCoqaheUWaEqfpECLgovlbZJ9S9O4hT4jyQ0tYlrdNy1nUzw8PJ3kWOOzi6WgDbEPaXEIZL2ROZAAO+9NvQa0JDmadV1a1TnZCW2kU4tvkrg9JeJ0ESvOShOcunr2bpEcZRKUbqoYNGsI75wh/ygXGgOu6bVgzyQ1gqdukzLAtKR/QKhkgyRlTAyrTGKhkg0cls9fihq5a4p2irFdombFUVPhL5n4xAZ9F6Kl7exTrDjO8Yp7SpKbnYZWqTpFfxMCXwi4e6oLQlYffp3xXANEPwRuzZELD98syFpHfCyFEv4BvpsqdZHhVhtASV1Q5yPSLGvNuGnr1kzJk9MVKRbY1Pw7o3JYfy+kv5HiYsCdz/l8Sf1Ac5t2zWgc0/nw47juKzVZKx9zPJeWZcUh8D6Yt7GCJx/ZAV/72tdw6623Yv/+/di/fz9uvfVWfP3rXw9+Jssy9eed73xnfcz3fd/3ee/ffPPNndp40kiOBnZsqtCYg6wkdeXuRD+sFfUDxCMnsnwWJDd6FXR316tVFQ/phtYc60llTSByZu06aUDwn/p8CdabVKFsyLWyCOQioj1ztZp44uC37s1K9hez4tA5U111qeRGijDbwHKhSIRExlxg7BAOQWqcaJSOhGdqBA1oVhyN6ACum8pC49ZxheSAnw/HchNYJEe6hNXPKlbX5n9X5K/pb6R+aJlEORSgsMyNiVYNXoKTiOo4n9wsqxaVJNMEKWwPiY2bdsayGzfzhGXF0VxVKdm1JXaU4PzUT/0Ujhw5gocffhgPP/wwjhw5gltvvTX4mS996UvOz7vf/W5kWYYf+7Efc467/fbbneN+//d/v3X76sy2CslZKLpqydDExtYkOp650RQSbqI8WR3azVHBdThdC25aOwnd8hQRR0fIBLXfsuJouSws6w2BC2VDCcWWCe4alP0vRjAJMaLaBhrpkVYcKcAO7bZTFyRr4rdCZmOQAuMU/Y2EtN5wVJqzjq4p9jlOnuh61r1yoiPdVClJ2Hg/TrH0SZIjf7TjNXLj6t/cFA0Erf6UTFtAfWlRkiNdIbz9KVmWAT0zd9v+JcXFzetp5MbJObaAZTGWRZvcZ9yKE9LhAG4GdPm/tOIQFiU5O1ZY5fOf/zwefvhhPProo7juuusAAA888AAOHjyIp556Cpdeeqn6uQsuuMD5/0Mf+hC+//u/H9/0Td/kvL5v3z7v2C4o8hLTWYZBMcX2tECel7UVoshnmM5y5EVZk4wsn9XuI/f1srHC5KVZWLMNylmerHsZTwvsG6QPpo18guOzDQDNIv9/WMNaMcVoe4hBPsN0VswX21mnxFHTWRaMDhrPCuxreU5nJ5BPcALr9a5zPCuwNv8NVINhWubexBWKnKJrrCvuGAsb2SSSw2cGaf9YK6Z4cdt/ptQPu6BrJFZ67poxkDeWFHo2IUG75Z5qi0WExkNl8eJhq9r90/W0BH+cgLRxT6VshCbTAsVgVvfhQT6r+7W2oFGfHc8GGBZjbM8Kx70zmg1VzR3gL8jr+Ta28/n95iJ6bNosnsNiGr1vuQg5VhvmmtIsanIsau3dS3DvocSk9NcGy8WjERvAJTc8ozbQrXQOnQdo+pRGcjRh9Wg2bEU6aNxtoZpHHDf+bACw/8ezwukj42nRavO2Yxacf/zHf8T+/ftrcgMA3/md34n9+/fj05/+dNI5vvzlL+Ov/uqv8MY3vtF7733vex/OPfdcvOpVr8Lb3vY2HD9+3DzP1tYWjh075vwAjRrbsuRU74UtOSZkFFULxNK7WwJG6qDWbo9HHBAcUy8TGgP+gtnVikPQJkUZBZCSUwFAdNfJGT9/LWQStxK/1Qu1UntqkUl3GVaWronIQkjRG1nRZdFQ+VTXlLETTLHauEnBXHcID02NES7ZBm515DvlZWWNpfPLzN/hz3BXh+umSiGFMuFfCLy/agua5g6Wn9N0N7I8Cm9b/bcxZ3F0teJYrqdFM7untIciqDikdQSwyc0i4eHb86SSk2nhZcO3wsV5FB8XQVN/c6L4lOg9nvSv3kSK3GFSjwP4lpxU7BjBee6553Deeed5r5933nl47rnnks7xnve8B5ubm3j961/vvH7LLbfgT//0T/HJT34Sv/qrv4oPfOAD3jEc9957b60D2r9/Py6++GIAlVipccu0IzmEqKvqJISHc1gKd+pgNOj4rk4TGls5ZEIIudCmyuRNSMmFsyYsN4AkZ9uN0Fjz6QZM4pb2JiRW1IpshhZ3KzqFZzO2LDCxfgcsniTQynkjiY4sYqgtjlbl8K6IuS5jkzy5GOqEinm8wGZzbVfrYllvuoo6gXhkZL1DFz/NtX1haEpOHOlmpNwkWo03giQ5FqmhY0PkJqWwLeDrb2TbdwJtyE3KhiclckpG7RHakJtllgvRRO0aZO41q/CmLLkjdY4yL84ySE7rUXnPPfeYQmD6eeyxxwBUgmGJsizV1zW8+93vxi233IKNjQ3n9dtvvx2vfe1rcfnll+Pmm2/Gn//5n+PjH/84PvvZz6rn+aVf+iUcPXq0/nn22Wfr94qOJKeVHmeJJGfZNUY0sReflGTZhi6Qk78VXhtL+82zgnJy5lbuFhYbp3q0axK3rDf8/ECTd8MiAW57uk28oeyxMSz6/RBiod2WADmEtu6p1DQBsUlXakIa8tvoOkIkbEvoICztDR+PTvLNWd6p0Cadb5wQOMCJDt91N6HFfNFxxxTf6FiwcpJESzgIYgPY5MYqbksIlQVp7mG5ZCdesHf5rjL+/ZD1RpZh0MiNVvAV6FZok6KxtD4no/goN5OW9I9DbqyH4jf9zVNruBtun+QMW1itW2tw7rrrrmjE0jd+4zfic5/7HL785S97733lK1/B+eefH73Opz71KTz11FN4//vfHz32wIEDGA6HePrpp3HgwAHv/fX1dayvr3uvD4spSsiFhDpJ9UC3UX3Zs1nmaXKS9Dj1aeeLENfptMBsmqk7+eksU0lPaOcmOxtpWYBm4hkWU0xnOaais/P7XATSt9oGG9kYx7GB9XyC49OK/HIdDp2XD1TNJE6QEyofhF6tpfmCwBfhYTZ1nvcwm2KEdhE/O4EuUR/Jobj5xNPiRNsjiE6KO1Krx9MGTvRciwzGzfW5BiK+O96eFo5FJjndw7RAMW9XrcOZFjUx0CxEvJ/Tdz0pc6yjelab8+aOyjUMsxfrz0nXCXcvOjoccb9rxbQmXqRvI1ikZ008f97uYEZgJcmmehzb8IzKIdayqRcSv1Ow+lGINHJo7QzVHKs/x8hNG7SqKVhrwOJzCB27WVRjejibYjN/0TtuLZtiVFaasBGGzVyQA5jBN7fQsxDPebKTBOfcc8/FueeeGz3u4MGDOHr0KD7zmc/gO77jOwAA//RP/4SjR4/i+uuvj37+D/7gD3DNNdfgqquuih775JNPYjKZ4MILL4zfAEORzZCxATid5RjkM2zPclN8vBDJAVoTmzZC4xjWlAmDa1m40JgEsEVeYnv+MU10HUNIaLyIWZ+Ts628cXFxkbGcQKWwmC92fEKVLhZtwtrIxslJrRYBf+bae22xNRtUzyxRGEh9hiZjEhvTvUtyR6/xyuEWaIJLITpugb92NdikSXsjYr2R0HQ/y3BPlbMccOafDIP5qbg4nhMYt125Q3JoMdqaDeq/R+UQmxip11+LaKbk4jaeFR7JseBHHenkRrPeSPcU12/Q/131NjGSrBG3FJdIiuV2mOnzYMidaLmkQn2ubX+cznJMyxzjqT9nRjVv+TbIdCBD37W+tZFNMCmLauzN/M/JMbk1G1TWHKXifQp2TINz2WWX4Yd+6Idw++2349FHH8Wjjz6K22+/HT/yIz/iRFC98pWvxEMPPeR89tixY/izP/sz/NzP/Zx33i984Qv49V//dTz22GP4r//6L/z1X/813vCGN+Dqq6/GDTfc0KqNPDSQayFiLquF3FU7BFmPyoJ0vQC60LjIGhO/5qai+9QW3jZZXLVsxiFoJK1+z3BVyf9lttR14bbg4Nl9tfeWpTORWMTttKg4kiAXEWmpSY2+6upCCBUebIsq065f+4aTsFG5Nv+pcnu4bh9fQ2a5pxaB5c7VtDd0PA8ucHPjpOlveAp9rsOxEFrwuUsYmOfSMTQ3MetAyD3VBqmBCyFokZchpObH0r4jzT1lkRvunppMi4UlDCE3lQUpNo6hmhPGzv987pXuKh5CHitcLLGjeXDe97734YorrsChQ4dw6NAhXHnllfijP/oj55innnoKR48edV578MEHUZYlfvInf9I759raGv7mb/4GP/iDP4hLL70Ub37zm3Ho0CF8/OMfR1G0M00O85lKckK6HIlYvaBUkrNTZEjrcK5o1BYax7BoVBUhlHOBwyqNwHU43F9L0Pz9gC1klBaIkx2q2kYwLInQMqKqtPIJkljKCYr/aItRrKqyREoitNgkXOcwmfcJnsG4Ek+TYDXRpTAnG8sSF1uYziNb6hT6QrNmER2CdC9RMra2C71cTBw3L9O3Sa2bdjwfg9xFrFlv2rqnlgWv9lbLKEdvA9Chv+tJHW1ysyhkslh5XUvY3hwzqOdvyxK1Vo83f4PE84xZJEfqcoZZuot5x/LgAMA555yDP/7jPw4eU5b+IvmmN70Jb3rTm9TjL774Yhw+fHgp7QOaQTimvl3MfeDzAem5q1A4rioAnrtKunBMd9WC2J4WzmJY5ZVJT2d+HGcAgKNlqSaeIdaKKaZlXutw6N7b5MSR7ZvOctO6ELPkaINjUhZ4Aa4AnRAjaVZouDUpcXeK5pbZSfC+xl+T6KppagOe82cjc/PiuMeN1YWJXuvy/Nok9ku1YtF3GXI5phbW7JK/aDbNUASGE7mquB5Htquew5ibCtCzzAL+BmeSFRhlQ1+HM11zzr89yx0XsIXQBsM9LrxQyfFouaea7zCcj8pqq0ZSY+QmVmizC/y0BDJazic3sXpobSw65KYiaPOJ/O65W/RMbGFrNsRmPvJ0X4S1yPhfn7uuajcWXJcrZgOUp0IenL0AKhBnZUq0LDlaWQdCJtxXzevlSXFZhSYff4fh7lxlRmNCagK5FOtBKIwxJRW/vAdpVtdMmHzXqEVryHBFy3rDsVOWHP6sU91Ug2LaOclfDHznFXIP1sd7ouzuegnADhG3yI5WloPcUzyDcSxnjyx2GCpkyN1TXYpsTmc5ZrMM29PCORddi1ty+A9BCxlPceXFLAxyEQ+JgxsrmWu10TQ33IpqWW+kpkXOVcsukSLrHqXCc3kq1stYiLhf6sB1TwFhy41WD22RLMb8mtwFGkOdoqCMpyggK7C04nBLDo+ucl1Wp4iL6lSHFV8fIzmArcfh0ItUpr22bMjIKQI3CQLNrmotnzo6HI7Ya8sKWybIyYzcC/Je3HBxf8JVzy1Cw2Pp1ttEGS0TTa4cVkbkJGq8rGdAkxTPC6SLstOfh0Z0raRnIciFSqYCoHbx79yyemiwJv5FdRBVyYfqHJLkcFgkZyxCd2O1gRqXrKvDAfQw3er1qUpq6Bg3HUNDbuic1XXjRCJVDH4ysxxLy1OX9AkW9Cgqndxo9dC6lAvRdF+W1ov/VHqdwnFTpYwdWVtQFvHlBJcnWSWSU/QEJw18EWxDcrqKjglkzVnUqqNpNVL0M9aA5EX2nOvks/q+NevVTmTTTQU3q1f/z7wJiO8atWgNTVwsyY1mhZBalEWhmYS1XEwnC2rV9kQ9kubua6tl0vIihfRa2s6b62+ctrXI4CsFxpq4eBk7ZgmN5NAPYRzQA6VYcULWNWmdiCXo04hNjNxY1hut/5A1cKcK3GrQsp8vE+QWtS2VultqHCC+qZiJvE1c90XXa6PzkbovoLHmyPnUSqIqSY5lzUnFShOcYbbtLOqpJAdoH1m1LEFuCNZAtBYqvjtqwsabSLI1I0MpIbbgLppdl4OHtcu2A3B2nUCz84wp7kPiYvncrMRiUhS5DKQ+O+u4tinNUxCy4tBvSfjauqeWFXovyQ7lSTpZO/1FckVNZ1nULWCRnNqVJha+2HMli6i1YZCbHyIz/Ifgu6rC5CYGS8+V8rkUyOSg2t96ZXqbsKeWabAg9Tea5ob6QJfoqTbHh0TGTVvj0XuWJo9bcZoAANdlVR+rZJAPYcUJjj94Y9VLyZoB2CSHICOrUklOGzK0DO0Fj6QC/PBNSe5206ogQWZ1AlfbEyi8UPP3A/FQZj5Z7VT21NBz1FxRO+We0u4vNHk3O7OxebxLhqa1JiGogVH8+Kkh4rxiNelvgKZ/8wiq0EJUaSF0MScnIGS9aaO/0QgQ6XA4uBVH6ixsS07arltzwXLIecCy8BI0YhMiNwQtcqpNuQ+pcVm0sjigbxC0iuchDJ2+Hh+v1L+lC1YVuSuWmy4RfZruqxIbF6qAmZMdLZycR1KFSLUWVSk325Y1pw1WmuCsK8LUEMnheXJCZR1ChTnbWnJSkvwVuZ8avS14yQagmaz4eTX3XFu0ydEiM31aYj2565R5E0LQxIyc8LSZLGOkx5ocveKELZ9tkZdLjaDSyE3InGydI9WdtROQCxR3R9bHBBZ4i0ztRO4beV5C7RrjIeBzoqNpc9TaRInh4ZoOB3DHEuCSnLAOx7WmVtdwyY2Wlyi2iPl6vG79StP/1OU8khL7+fPKTvRxTXPDv3OP9Hbsk1z35VxfEbVbSN2ASEswd/VzkmNZcwYtnvNKE5yCxdfL6ButBgsnOdXn/Rw5KdXHLZITcmW1IUbB/C5il8TNv1JoXH9ODPgQyQkRn2VZe5pq177p2jUbz5wfK9eGtVtc855Ve9IjIReCEDHl/SrPS+dHvg9Uz5dyv+wUrIgq/qMdvyhiVgm1erzQ32iJLTXEinsuA3w8a9acWtujkBzvWE2AzCLANGwEyB1BupAkyZEgiw232qSQG6m9qY5v3MUyPDyFAC0LbYIUqusvpsmzsnSHyI3sG4vqwaQWh7fBsuq4BV/t8cO/Oz+RaENywtacXoPTCpLkePkaCmnJaB9ZpZEc+SOxSImG8CRgC0B5NWHS4QDNPUvrgrbwAs3iu4gLTSMScscmd51aunfArhtjWSYkuVkmQlVxvWyvxm4y9fm2EUiGRMHa/9pnQudpGy7blWjoz9VdHLWJ04IMEQ+5pxaNoOJuKo3kSGuORBc3BS+8qRFVWfwQcMmM1LpJYqORG4IeVp2ayyvePyVktWpCyAIq35PWH98KvDPZzTm8sHDxvXdJVyCtOJqoHWBWnbmbiqOtFYdcVdzVLYkOt+YA7SxlK01wNJEbJzkyNbemx6n/NqwaIZLTBVFhbwuhsYTVcTSLQEgE21VcbE12WrtidXSAptK43DVquTb4blFe3/pfe6ZyYeiKQiGLO4GUycIjlQmWrFjm4hDZ8YlO+3ykNH552QFpqTwZCxFhmUEGi5Aoa9yoRHW+YSBIkuOdW7iEORGQ5EZq4EKuqZRcSl1dVVrOLPqtFejl98efTagvWfqb1P4nrTec3CzDZSp1OPIa/NqpCLlFLe0XJ9aS6ND7G9kE6y3G7Y5mMj7VsZFPUObsi2Dx/IC/+6WORiSHOpcszFm95mc6BuAU5wwhRIbyuRaGL4IkqORIWbzIFHx8VmUEdqoJ5wW280pw5rqpcu9+QwiHlTYWoxCsCWwjGztZWDll52bTEPELTawWqGrxMJtitKS1S8tYzOG7pOauqsDz7QrL6jLM/CzO4cl9Ma2EBsv1wgtBcnCSy8NU689l06TkZG2RWkmcMJ3lTQb1eRZwKlZLO2X6TunYyXxsVguQ/YxHs6HnPtkQ3yXNBZOywCRvdD48W+0wm2FS5kFtW4jY0GuS3DRtatzFcsPh7vSX159kduauYeEyp5Z6TFZiUsaJSEp4Ns0VyywXsl2vaU3fqq83XxstF3hqyR3AnUeIvMqQcsowDqCZ12eAnx/ZxkpbcAjcrRHMhtvRVQUsZsmhz1rp+T2NjGy7MRlIN4xVEVcz+Vs1ulw9kh5t1TW3BO2ErB0AtV/WMZG7Rm3HKH39od2FRMjvnhp1oVkG21hxCtEP5PONmdVjZLiNFSe02x5mpbejbWNJiSX5o75VuwGZzqPuH0vIXbQsgbFGgmYitw6/1razY9enb2txtMN0dQsJzz0CNGTFStAn3VGa1UaSG9mOFKSQmy7JOfVQcGvzoOuB3DlFjpn2OyHLemN99231N5YrKzV6T4OWv8qC/J5kgIOM0mybkmOlCc5AEZY67zNXFUeKq6o+NkBy6Cf0WlvTdlNI0t+pbmRTjyRUrzcTPi97MGD3ziPJtPIV9NPct29dSBHAhqJ2OLg/lidt01xEcmKV1wv1gVThLD+HJA8xV5Ukzc57hsbGst4U2SwYBbKI20yGd4cSAe7cbttvf4ws88UopivqEgkjF4lYDhxtTMc+Y5EcwHcnaM8jLdy6sZpIy4oMtZcCfovYeCH6gtzIjYYlLrZIikaYLYS+Wy2JoSzSS5BjKHTeWNsk2W5b46rRZ2lC9bYWRFvYrpEcrsPpUviTzyNWQId7XPWs1lsU21xpgkOQ4lQtGy7QWHE4tKgqwM4VY4mJNasOP1ZLzz+Yp+3nFZMBfzKR4IOOEwNroGr3rVkZmramVbimHXabRZe3XYvkqbOgirB3fpyVKdW/lr6r1a5vHd+8ZzwDS0Qs+pX8kceEEAuVdwhuB2JC1cdjOgltsu9CKDSrGNdOkPVG09/Q35oWpZ5AZTRdohUutqCkJPykc0grTlvEvvMY5AZCVna2fgjSagOEyU1K3hsZBZpCbLTzaUEI7vszdWxRLi3t3KF5IRWxeTBmvSGkCoz5eax+FhO1E7q69eRzs4iOJDmpWGmCU5GZEKvfXsiK04bkSISOCRVXtM3Hyq7OEMnSpFSJBplVRskHBPgLcN1OYb0haM+Tfxex3aYW6i7vp94Jsh+5K6VrabvF1IXX3llqi7AtYAfgZM0G0gmMdE9xuKkP4hYYDS4hjp9DWg1T0YXscFiTrJPgz+hbvrszThBi5GMZQQUcIStOW/DvVEaoWOJfKxEoQXMJA3rYr0ZuCDIXVaz9aferbVb0LM0csT7l6odCyS7583bbrs1LKYRBWm+WUS5EWnGq19JITgyU5JMn+2zeczVxWhRrF5Kz0gSHIHcpWjZcQig3TiwBoCQ5qTlvQsUVh8XUy9hK90RIsVRotXmccHFx36FFmO+AHAtXxH0SbqM/KTgLqSFapNc0f386ifEHpZXpV2tTCizXnbWbjD1fjUDF4E4qRGRK53cI2o6LoH1+Ga4rL6yeWW+4ddDSl1XtsEhq+2e4TGhWnDban7QyCLrlwVughUU0FCZuWW1obEpyw6+hlfqQY9WOSppGiVGbMcnraK0p1qnmuu31PsMsbvHUsGgaAg3S4hMjOV2gjXVts2QRnS4h4sCKExxtMdQQsuLUfwcSAIZKG8Ty4UhyQ9Ybck9x0ITOqzpbg4g6l3xfVhOuzutrcSSx02rS0PtSeyNdCIN8ytrsk5MmKZk/sXFXAycuvEgbvzdtAgwJLbWBqVvD9N1brN6ORpgBt2o4oD9j/nzlOSUsN0GXkOnYxOQfp7mmbMIkn69lTdHqIknwvlV9z2Onva0sTAFB/TKhZjVO3DlrkZR+qgP7e5EWUen2tQoeSjeVZbXh/2tWI95mPxFcyAVsuz9TFkUtQ7N/vpmZRDUUQRUSF8tzkOWcYG0I21rwUrwGWh/TSI604liW5pS5RSc+PtHh5+lFxguirRbHclUBaSRHA1/grMy1RBxk1V5rUhtmpTCV8gXKL7ZHg80Z/IysEMmRizH/n5MbGT0lny2PcrEWTQ3Ufj6g+G5Rq2PCxYzyXCENimXFAeJCY5lMUhb548+VIJ8tf50fS9YbWQuI6wYscaR8LhZCrirtvPIz1edKdcLXyKUmcrX+r8/DrDdcf6OTUmllmtauTFd34lseQwkWu+a8sepTxeBvIBrtESE0pqTlw3kObDxIUmyRGstq45xDITcpeahSXVOh/mTpLbU+JvsSfd6KBOvSVisyKMVF3dY9ldo3Uy2FfL7R9G4WNIu8BO+vXWsArnQenPVsgiwbVDRvhpruWYnF6MvkrLbIZk1+AJbHhHLjAKjzxeR5WU9YqWULJLmR1hvunmoWMjf7o+YPplwMG9nEq/kEVINuKx9gOJthDN8qMJ5WSvphMa0ZPb8nz6ogwnelCJGTqxCaxbE5juejGWZNfgVr0KQMQH/xc/NX0HVo4PFrOvkc8niiOp6Hg/oTPT/t2dZtUMoTALp7ynUh6LoBadWwLC9WHg/XVac9+xQ31xgjDLGRTbycNOv5NjAbAPl2HS5Oz04Li5f5b6yFhNwFKTWb6jw1O+AqWBTyGVT3325BoOdA30Nzrgm2ZsPm+Rm3Ly2m8nVPT6GQG2m94VZcCepTI6NPyv40zKbOfdX3Pe9TIWurtMY6Va5F8VbL6ryR5RiVZX2f46zAcZxRn2eQTzGcVfN5l8ikLhmMOSj/kvPaLDf1lAStoK0F7vaelBk2silG8++Gg8b/RjbBqBzWucfaJPo79UbpKQYrLw6AqBVHy2MiSxpY0GoOcXIjd21ruR/JkApy8/CyBwRuxXGsL8LiIH+AauG1dDe0K2pTnTfUfgC1FcfywzvWmZahqPx3zG0lSRSfGGNWHMAlLPx5aq9Z1rGmXbpZXd5rSA9gERNuxdHIjWvx6VicNdAuTRwq9RKAEhUVINJk6aDFhoT2svAsh5b3ihATGi9biKxFJqZY5TbYM6mJrrDiSEuM9kOQx/HzauSGkKpnsSyB7jG+K5quq+ktZb+hzMxyDMk6dintrd7P2d9hS1XThsXnR0LIetOFGNF8FWpjzALXZk4F2keorbQFh2M9nzhWHILcffMv09p1T+eMd3uWq5YcwHVZzWaZSXokuSFwl8Sa0LJwsmKh6lh+VtrqvcZNtT0rqglg/rc0m44j/Y0vvJr1BqAJZ+LssvlOSLrWCDQoqh1ns+PfoO9SuS9631Liy91iTCvCd/0b2RjIG6sL9SnaMa7n29H+RNl4aTEl66BXcTxziZH2jGVa+Tb6m2biaX7TLnnoWACr3VcIoYVI+w438kmVRTWbYgtDx/VH2XRh7Lj5gmTpb+R9UiZqzWrkPZdi2lpwmZq9nINnNdYgNzpaJnO+EEtRK/8+nc/MLbrV58aNdZujtnS7Fh4ObXHSNCs0DvmmhG82+HiUfVID70/8WnxMapFAVoZm7pri4eWOUDpw7027OLlxraAb2RjHsTEnkCKT/CyfzwkA5s3eCaExh8yirUEj+1LQ7+f4iRFS+u7c7NoAZTVuLDmpWGmCM8ymGMwHUmxio46/raT0pgVJuqo4yQFQlzcAXJ+pJDdWwjyuu+AuH56bgiYKd+KwSQIt0hvZBJOswCgbYlQOq4mgLLCVz7tIxHxrP7dm4tX0Edw9pbVNA90Hd7M1x48xKtdMksPJTch6Yw1GaVbln5OuMYdwcczdLEDVn6wdECc6GuQOSpIbvtBTu6QWwtMlZVPz3kMkRyJlMXLPMfFcjXU/ZGOFSCInOdUzcMlcszufevfa9Bc3HxRdv7FWNFbGQV6Rzxe3OfGoyiNw9zONfV6eBXAtNUR2tAK8BIvcWAuOq8Ejl7W/AA/Z90vfp3Q7Wu46WmDIPRvLhq1FGmrEhh9rRXLF+pM2J8j+VL0/J28aRzBcybw/SeuNJJCuRSe8oNNnuAWMSwOo3wHwSpDQ2iIhZRBdyJAlnveDWpr1Z03ZUPD71LR4ANT+x+eTEZtTieScaOF4WmmCo0Ht/EaNqu1Z7pEcFNUOnpMcOlYjOhZkMjcp2uV6ljOLcd25+M5CjwrKMSqlALGZDKgWDT0LbpptU2uEwPUg2sLLU+jzXbbuYvLbLu+jgqhf4hzj7xgJi/j66fO06yULBNd3bXGSmEhyQnBdM651DPAn5uYeJw7Rs9xzutVMJznuc0jR2ugWhOo9XzehWVaJ5LjXboiztAzSvVbXcF1qpKfi+ilyZTh12bKi3knX9aGMmmyS5BDauKW0AAPNegO47qlQ7S0NZNEFGmtaBcOKY8BMn6DovixyExP6p4LrifiY3JpbVHmbtmdF/cy2ZgNn3pNjSFpvLMtgc7/+RESWw41sgmPYqOsB1pbzRCuOVgeRI0RyUnWgWjZ6ubnSBP3Satg8j3RtH7cQE8lpg5UmOOvZNvKMPYLcNbty8MmUrDicRTuuBUZyAJhEJwQZbi2jkWhB2zfY8oS6VXtdFw9HNeDo/O7kTsUruRWnfgYzf0GJIURunOR+bBEKgZMcPilLUKG2UbnGPutOrIC9uLvX9P3GZMXh4Ltefo0tw0Ul+xQJZqlfee4rhQRp6eRl1Ae1R1pvqnvT7z3sBoiTHO0cfKLXiConGNxNxd1+BMeSI16XCQ2le8qyFtLCzhcbLvokK854CkxA7tq8dkFLKw5gkxyJrtYbacnlruomLUIT9q3VR5JE0+8TFUEYcZGxAc+KEyA2/FpalEwb15Tb/mYh5G5HbhXEzLXaD+pNa0N0pCudt12LLqO/XfIc3xRuZBMcxxlVW/MCx6cbzvs0xrtYcVKhSSRi61Qz9/ik2k6o6Vty+FzCIS3EJETeaFHOYqUJjoTccXNsKabwMdOkePqJeV+0iE4MMiKJa26sPDJ8tyqtILFcDJ4Ggrtb6FkMgBPb661IDteCSHLDJ2GN9Yf1NyELQON29P3AU29SbV6PT6iWWZ9fl1txgLC+S7NCWNDTx7sEsmq//ZybHf3Ym6xtMlyBk5LQxETva+eQsCwHyW4FAe2e3Xb5Y0K6NYgMTfLC303nVQ0eVxNlu6oAl7xIsmPlvQL83FfV/QnrDZsP+EKTspPmz4H69EiZJhw9DuJERy5uKcSGv6a5ptzz+R2Bb3omZVbPabz9ZBUEYM7zHFruKM16Q/eREsI8RA5kwKgs6/lmbU5CR1n1XDeLUTU3zNwgjDZWnDYRu075H7Mkz8xbh6w1CIBjiQ/lMaPXtbl86I3PuN5PYqUJzkY2Rp4XtTFD+p017YQmRNNQi287ZraW+U0At87OvsFWvZDJMGtK8Cd9n6GFxtOxZNWzmNDzYSRne1ZUhG8eZWBVeJYLLuAuQE17J56p11+YcufvUTljA8Y9diy0Mdp9aubw6rrtd4sa+AINQJ1QB/m0sgbOBhgW47pvkXhW03txuMkYt5OeM7leuP7IvydrIspNkmN9Ptzn3M835MZ3K6QsSEDjmjqrGDl9Sy7ybqmCRnCfook4azjGC5M1MbYrkgMA2/UbzbPSyA7BS/5pRE/SsVJYHNvsNPfZLDT0vdD3Sd+F1p9JdxNy/UpYiS/5uAwlcLMqcXfpTwDqjZvjaov0Kau8hLTeSEE0bcpkW4fgcxi5vNF8dr6xnOQF9g9erI89gcoKvZ1Xn+ckh9+M1HhSPwpZczRyo8kjtE32wCDVGvm1rLmy/2mQm8lBC7flShOc6stQfMzz/8mK4YReC5IT0084UUaJZIe7ogCX2ABwyM2Zgy1sFiMvCokjJsyTO2g+GYQWl3XMrRBFNZmRRkcmxQL8BVdOFmRRcCMo0kOL2/jqu/j6+YAMWXGa4+1+hbyaPEfSHUpEcf7MUsgj4OpOuFmdP2tObuSuU07M8p7lcwiRnHRy47oaueVAs+JUzwXJC5LsW1pUjgRdizQRZDmSVhwAAZIDcGsOgNqiE4OmtwmRm32DSU1uzizGOHPurpabHW7N1Sx0mrtQziENyQFqohNAyHKmuaO0cehH8qX3Jz6nNecPjEn+v2iTqxG0x1EIQ6PTbmRTTLIqHw4FedB3x/sdgFpwDMCz5JC7tPrft+aEYAW1VP/r5IYshmcJ/afcUOhkx9+sVn/77ZSu8I1sihfQolxJ8pGnMWKL0Wg29EhOdUzbnDNTz4/KwfPF8MgY+l/qWM4cbNUdazMfeRNaVT/JnRyGyDHBTJ3YpBUHOfNTi8WFNBDc0rCu3BMPrTTJjRKurJk1aZLg7ZduBjpPKJxQTkYh15Q1qfKF3SZHkQlVPMsaXKczt+xISBJpEZvq2KkzKQ8TJuTUxaT6P82NoJ9LtxxwK07s+RH4giStVdW1dOE6MPP6Dy02aqK87aqnS5eB20h/J03gO2orgtKKntRcBGcW43oXHVpo4iS+2fCMmIalKzRSI9+zrDZtLanNOWwrToXKVeURZ8K8b1lJCkPjyIpAtMgNPz/PRUZWHIpiPQtjvAAAA+D/thtNodb3eGoSLWJXoiu5oU22jFaUG4rq/tptmtz3fL1fm83sShOc9WyCIisarYYxmQ4LV72d4lrgglGgibiIRctwQkPXqv73BbpkuaGO9ZLihDNYhsbOh0hO9Z5uxakwBgpgOJvWu1n+m6wQnPzJCISYmZcmi818ZFbytklGXHDs5qmxd4ypu0VpyQEQGXRuv6L7m5SF17c4aQQQJI6AzE7sRkpZz5m7pizrTai/8OdgT0oKIROTvEWym3PYUTz0zNZRZdeVi5G8Z7rfkAuB9x++2GyWo6q9bLEBgIHIpExpIngEZbXYVGPBSX2vrDdOUlCxwAAILjJ8s3NWMVI3O3SPmutEugtC4v2Q+1fCIjXVNcLEpnpdXxg1wiD7U+geXlKcqIIP2JisLT1KdKtj6Q2MI95uvqGUIB0ONYCiqbgVZzMf1ccf3T7DITnbswKD2czrewBPUdKsOVoflAJiM1qXeRH2DSbOJpsshpoHgZ5LquA6RnKq55S1suoDK05wqomMD1Z7MtVcC0RmuKtGkh2gfRgwJzWALhrli9hmPqrJzdn5yOlYoYkhZMUJPRPkqBeXuqSAyPnCB7xm4qVjNXKjLUK87S5BqxrXdH497Xf9WWVi7bJbtCxI+nXH2CiqiC7uluILNSeNgG81lMSxbodBIKs2ppMbfl91++D/zYmO70ePExsNlpaqQdUHNzCpyTYtSnJB0u65XaSc7zLQFputjDYu+mLDS5jQgkPgkS9acVoAXtSkDC6Qi4wkN3Kz45UeUb5Peg4jIeq0xlIMfui0/K50YkPtaP4Okxv/urqewyVcTZ8i0bSV9BSA45JK2SS0aeskc4n1S4oT+Pr0zHrjun/wYj1vbGUzTPK8tiI2mHp9T4aID4zmaHXtALvfkVuKk5vNfITNYuSQak3zJZ+LNp9waO7wNiRnpQnORraNIqs6tSOkI8sN7/RsUedER7oXJNkBZOSV/8it7JmafqV6f4qzipEzmdOgkx1rvlXwBlvIiuPDfya0uNTWB+gLiHQPSP9118lTWhXczk/WKH1yjpnCQ+JAOSDljjG8IPjPsd495s0zDBHHpt06gazaJV1/Nrlp7kWEwhuTs2XNsY61YLlKg5YwwCHbfHw6xFVxIViLELWft4ETA77YnFVURGc0GwLzITvJcwxmA4/ojKcFikH3umKa/s5KtSAtuXI+qO7JSPQpNguaZXIRkhMiNdQu9/20cehdR/QnbU6rUxBw6yBQu6iI7ADueNPy9ljjSGs3YX2ekmSr3K6PG5Uz8Jw448wf7/UaNEDV72YD7BtsOf2uIc1Tr+9pmZvr58ZdVEZAS6zfcYshjZ3q/nyXXWgd0iCjODeysrJkJWKlCc4gm9aFu6iQV4Vmx20t6nwR2sgLZ2c+mDWdhkfCAACKsEBPK2MANG4HAE2HEhP52XVHcztWaGKIRSRVx0wwLOehjOVavZMG4hmg+Q6oOpdfVE+KfUPWG+fcyqTmMnx7UtZM4Sk7Rbnz1SxIHBvZBGtl07eG2YvOc9R2j3xRsJ6va4kSBMcgNnSs7CP8/lN2ndbuX76fch76/oi1WO6+xo3AxiaqRck9ziU2a3MLp2apkpouJzdUNgFy4NhsAy8pTszb1pCpjXyCo9tnVBGF8101LTj7MMF4VtSLDpV0OWMQtiJZurvqb3+BkdFi3HIj5wPtvvn3wJ+DOx9UzwOwNwwhhLJcu8fZ6Qm69CfNVcXn8eoalJn5xWYjoczPsj/RZ/mzDc1Z65m7zK5nA6DcdvodWXHOzkeNtyCj4A23341mQ2zNBjixvY5JnmMfqo3z9qzAeFZg32Di9D3qd6T/lLUBtYSh1d86sSECTZtsstxwD8Kw5TxOiJGdUTnDem/BScNGNsVAFT35k6k0j9NA5y6G+vO5uyh50TJaW4wduqXgr9rpdypJEGSn4ruI2ITgPouG6EzKop4IeCI9896MHVB9zpYTBe2A6H1tYk4xZS5iBtfcZARrkSYSHXqOMsdICsEBXFJD7/NJWSOR1d++JkNiPRs4z5w/g1TISZ7Ol0JyZB+kRak61l90Zf9y39PDd2Ub6nMykkPfES2Uw2yK49MNnIktbM0XmK15lCVfdPbNSZhVHVoWC61es93T1AZpxZWENqQ50rDIhiF43mgUJ/2f5uKUfQmw+1PMdc3HJKBZeZrj6HVrHGlt1trKj5vAn3drkgOY/W5jnqOp7ndzkk39DoBDtgFgn0KwrUK1ssyLJY3g/U6uQ/JeU2Adx+faSYschitNcIZZiXVlcuSl2TnD5ztHacHYzEeOoJUvTGTaBsIF6qo2hXfmqX5gSW6sgRaaEHhiJbnIAHAWGn7P2oJD0HJfWAuu1X5JdNx7gLIDxTxSJx7tYw0wi1wB8Z0v3acsGEdWHf4cZWJCi0BqCQyr1/XMsJJEVu/ZvnF11ynuP4bQ5M5Jk+4uBTRrGI+O0wSvslBj9bn4blKOgwnvw3l1Lczgie7X88ptTVZcTnZQjB33tVXqxEqrQP8vQzwu75t/N8vaMIQQyq2kIZXY8PdSSE59/kzX3Oh5ofz+VP0f3iSktFcKjjnJIcvvsJxfTwn2oM00WXWcflegJtshWJIIAMHI19R+Z22yCSlzSsxqbGGlCc56VmKYlRhm256/mRajtbKZTAHUnY0vMLQIbWKkWmv4QOKCRYKqE1F25NX//q7csn5U92gUkDMsOX7Rman3XKp79sOw6RmFJgn3Xtx2V6/FTbzWPQBwiE7VTn/hTNkpplyTf9a9rv8MnWspZAfQF+xNuMS5uQd9h8nfC03I1e/0Z229r01OsXPw4/ii1KS34HoYfWxW7bd1X6FFSC7yGtHazGcYliWOMxf1WjbFsdkUm/mLc1dHJRqf5I3bg9zX9DeQZsEFXCsut+BKbRUdK91xZ8/nljbkJvQc3HHSbmEx77GDXqtrfwrNaVVbfKuOvUmzdXtdxhFvr0ZyNrIphuUUx2cb9Vw/yoed+92WkU9LRrxanoP6b4PYNMfqFsMY2ZOwSM8QeSvSkpVl2Z2W71EcO3YM+/fvxxf+7QIMz6peI9U9z4dB1gtOUPhiMzZeJ7QtDNbF8lH9rROE1ImBOhNnx0QMQs+FEIo+cO8nTdxrTcQp9+C0KzIpt90lplxXe4bV392eY3OsbjonWBNyynMmtL33ZSKlDwLu86veI7Kji1etRShE1CaYOdeelBlGZeGQe7LuVm0qHBcjvUbHNm0PiD0jltvqddv92LwetoamfMehPr0sxFwWi/RF+V0StPEI+H1KwhJCp1o/q3M03/2kdPsq7/vL7ncpa1BbrwHQrt8tc17ZKrdx7PgMr3jl/8PRo0dx9tlnB49faQsO0HTSjcx1bVRmzXnHYx1yKDondaCz0VhmNLLTrk26HqfNLkJ2Kj7AqnY356qP5RNbRm2xn0t1nnDipVCUBG93qO1a+817gLsTTcEiA1Be17mmM2+6lqS2z3EzKXusnkOk+j99Mm7OZ/eZVGjfm3aulD5Y/e2GAPNnSO+7/6ftrp3Xym0MM7LmVH1/WJZ1m2kO4ILV8TxfzijXCU0brVp1X/6iY5U2sDc6WXQ+0L5TzaIjoZGeNpos7XoxWH0JsOcCrT9VCCeqtEDPFYhvEKz2ymfvW3Kq9oX6HUCEZ/F+V92X3/dS+p2cc2L9znoGqVjPBlhvUWxzpS04//1vF+HsTbeTyoGrJR/ScizEdgFt0aaAoTV5hyYEgtbBUq0hVmImiVSzdCqx0dBl8bWQes3QddtYlFKfYwoWcb+1ue9lYpE+GEKI1MTIm2YFCFmUuDVO7ppTNzyadc4lOmkWUCB9oyORMo5SNBPL2LV37Y/LGJNeWyLaPOfYlu3m7d2Nfkdos7GuXrctWW3XoRTQczp2fIYLLn2mt+DEULFBMfnPf9e7F2XB2Mx2xmzrwiVMOzXAtIle3eHOfzvak8R0/M55IhNf18GQurtb9iJuWZZCzxBY/DnGsFPPeSfQtg8C9iLbxiKV+t76/Hr0PWnjf1TyxWHCXs8AbJnn1pAmyCUrgtvuZc8HBNNC0gInq89Z5DXWp5Z1va6fl/3+ZPc7QO97er9b3uY0FXS+YdbXogqCjFYvnsgCDyvMfk/28mDRqRfV1xYhX7HOs5j7TWuv+/5OEMfmnnbm/Pb1bCz2HGPYnee8LHR/fqH7bnfPsg3u9eT4XzOsBmkOgjQUyoIh72j58wFhcQv17vW55VrXJZZ3X1o7V73f+Tj2QnXeFOfTShKcr371qwCASw789y63pEePHj169OjRFsePH8f+/fuDx6wkwTnnnHMAAM8880z0AfXwcezYMVx88cV49tlnoz7QHi76Z7cY+ue3GPrn1x39s1sMy3p+ZVni+PHjuOiii6LHriTByfPKf7h///6+oy6As88+u39+HdE/u8XQP7/F0D+/7uif3WJYxvNLNUwsX93Yo0ePHj169Oixy+gJTo8ePXr06NHjtMNKEpz19XX82q/9GtbXlxEkuHron1939M9uMfTPbzH0z687+me3GHbj+a1kor8ePXr06NGjx+mNlbTg9OjRo0ePHj1Ob/QEp0ePHj169Ohx2qEnOD169OjRo0eP0w49wenRo0ePHj16nHZYSYLzrne9C5dccgk2NjZwzTXX4FOf+tRuN2lP4O/+7u/wute9DhdddBGyLMMHP/jB3W7SnsG9996Lb//2b8fm5ibOO+883HTTTXjqqad2u1l7Bvfddx+uvPLKOknYwYMH8ZGPfGS3m7Unce+99yLLMrz1rW/d7absCdxzzz3Issz5ueCCC3a7WXsG//M//4Of/umfxktf+lLs27cPr371q/H444+flGuvHMF5//vfj7e+9a34lV/5FfzLv/wLvvu7vxs33ngjnnnmmd1u2imPEydO4KqrrsLv/u7v7nZT9hwOHz6MO++8E48++igeeeQRbG9v49ChQzhx4sRuN21P4GUvexl+4zd+A4899hgee+wx/MAP/AB+9Ed/FE8++eRuN21P4Z//+Z9x//3348orr9ztpuwpvOpVr8KXvvSl+ueJJ57Y7SbtCXzta1/DDTfcgOFwiI985CP413/9V/zmb/4mXvKSl5yU669cmPh1112HAwcO4L777qtfu+yyy3DTTTfh3nvv3cWW7S1kWYaHHnoIN9100243ZU/iK1/5Cs477zwcPnwY3/M937PbzdmTOOecc/DOd74Tb3zjG3e7KXsCL7zwAg4cOIB3vetdePvb345Xv/rV+K3f+q3dbtYpj3vuuQcf/OAHceTIkd1uyp7DL/7iL+If/uEfds1LslIWnPF4jMcffxyHDh1yXj906BA+/elP71Kreqwijh49CqAp/NojHdPpFA8++CBOnDiBgwcP7nZz9gzuvPNO/PAP/zBe+9rX7nZT9hyefvppXHTRRbjkkktw880344tf/OJuN2lP4MMf/jCuvfZavOENb8B5552Hq6++Gg888MBJu/5KEZz//d//xXQ6xfnnn++8fv755+O5557bpVb1WDWUZYm7774b3/Vd34XLL798t5uzZ/DEE0/grLPOwvr6Ou644w489NBD+LZv+7bdbtaewIMPPojPfvazvZW6A6677jq8973vxUc/+lE88MADeO6553D99dfjq1/96m437ZTHF7/4Rdx333341m/9Vnz0ox/FHXfcgTe/+c1473vfe1Kuv5LVxLMsc/4vy9J7rUePncJdd92Fz33uc/j7v//73W7KnsKll16KI0eO4Otf/zo+8IEP4LbbbsPhw4d7khPBs88+i7e85S342Mc+ho2Njd1uzp7DjTfeWP99xRVX4ODBg/jmb/5mvOc978Hdd9+9iy079TGbzXDttdfiHe94BwDg6quvxpNPPon77rsPP/MzP7Pj118pC865556Loig8a83zzz/vWXV69NgJ/PzP/zw+/OEP4xOf+ARe9rKX7XZz9hTW1tbwLd/yLbj22mtx77334qqrrsJv//Zv73azTnk8/vjjeP7553HNNddgMBhgMBjg8OHD+J3f+R0MBgNMp9PdbuKewplnnokrrrgCTz/99G435ZTHhRde6G1ALrvsspMW1LNSBGdtbQ3XXHMNHnnkEef1Rx55BNdff/0utarHKqAsS9x11134i7/4C/zt3/4tLrnkkt1u0p5HWZbY2tra7Wac8njNa16DJ554AkeOHKl/rr32Wtxyyy04cuQIiqLY7SbuKWxtbeHzn/88Lrzwwt1uyimPG264wUuH8e///u94xStecVKuv3Iuqrvvvhu33norrr32Whw8eBD3338/nnnmGdxxxx273bRTHi+88AL+4z/+o/7/P//zP3HkyBGcc845ePnLX76LLTv1ceedd+JP/uRP8KEPfQibm5u1FXH//v0444wzdrl1pz5++Zd/GTfeeCMuvvhiHD9+HA8++CA++clP4uGHH97tpp3y2Nzc9LReZ555Jl760pf2GrAEvO1tb8PrXvc6vPzlL8fzzz+Pt7/97Th27Bhuu+223W7aKY9f+IVfwPXXX493vOMd+PEf/3F85jOfwf3334/777//5DSgXEH83u/9XvmKV7yiXFtbKw8cOFAePnx4t5u0J/CJT3yiBOD93HbbbbvdtFMe2nMDUP7hH/7hbjdtT+Bnf/Zn6zH7Dd/wDeVrXvOa8mMf+9huN2vP4nu/93vLt7zlLbvdjD2Bn/iJnygvvPDCcjgclhdddFH5+te/vnzyySd3u1l7Bn/5l39ZXn755eX6+nr5yle+srz//vtP2rVXLg9Ojx49evTo0eP0x0ppcHr06NGjR48eq4Ge4PTo0aNHjx49Tjv0BKdHjx49evTocdqhJzg9evTo0aNHj9MOPcHp0aNHjx49epx26AlOjx49evTo0eO0Q09wevTo0aNHjx6nHXqC06NHjx49evQ47dATnB49evTo0aPHaYee4PTo0aNHjx49Tjv0BKdHjx49evTocdqhJzg9evTo0aNHj9MO/x9uho7rKnB1AQAAAABJRU5ErkJggg==", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "%matplotlib inline\n", "\n", "import matplotlib.pyplot as plt\n", "# initialization\n", "T_b = Array(W_TF)\n", "X = W_TF.local_mesh(True)\n", "#T_b[:] = 0.5*(1-X[0]) + 0.001*np.random.randn(*T_b.shape)*(1-X[0])*(1+X[0])\n", "T_b[:] = 0.5*(1-X[0]+0.25*np.sin(np.pi*X[0]))+0.001*np.random.randn(*T_b.shape)*(1-X[0])*(1+X[0])\n", "T_ = T_b.forward(T_)\n", "T_.mask_nyquist(mask)\n", "\n", "t = 0\n", "tstep = 0\n", "end_time = 15\n", "while t < end_time-1e-8:\n", " for rk in range(PDE.steps()):\n", " compute_convection(u_, H_)\n", " compute_uT(u_, T_, uT_)\n", " pdes['u'].compute_rhs(rk)\n", " pdes['T'].compute_rhs(rk)\n", " pdes['u'].solve_step(rk)\n", " compute_v(rk)\n", " pdes['T'].solve_step(rk)\n", " t += dt\n", " tstep += 1\n", "plt.contourf(X[1], X[0], T_.backward(), 100)\n", "plt.show()" ] }, { "cell_type": "markdown", "id": "33635218", "metadata": { "editable": true }, "source": [ "A complete solver implemented in a solver class can be found in\n", "[RayleighBenard2D.py](https://github.com/spectralDNS/shenfun/blob/master/demo/RayleighBenard2D.py).\n", "Note that in the final solver it is also possible to use a $(y, t)$-dependent boundary condition\n", "for the hot wall. And the solver can also be configured to store intermediate results to\n", "an `HDF5` format that later can be visualized in, e.g., Paraview. The movie in the\n", "beginning of this demo has been created in Paraview.\n", "\n", "\n", "\n", "1. **A. Pandey, J. D. Scheel and J. Schumacher**. Turbulent Superstructures in Rayleigh-B\\'enard Convection, *Nature Communications*, 9(1), pp. 2118, [doi: 10.1038/s41467-018-04478-0](https://dx.doi.org/10.1038/s41467-018-04478-0), 2018." ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.12.3" } }, "nbformat": 4, "nbformat_minor": 5 }