{
"cells": [
{
"cell_type": "markdown",
"metadata": {},
"source": [
"# Implementation of a Devito self adjoint variable density visco- acoustic isotropic modeling operator -- Linearized Ops --\n",
"\n",
"## This operator is contributed by Chevron Energy Technology Company (2020)\n",
"\n",
"This operator is based on simplfications of the systems presented in:\n",
" **Self-adjoint, energy-conserving second-order pseudoacoustic systems for VTI and TTI media for reverse time migration and full-waveform inversion** (2016)\n",
" Kenneth Bube, John Washbourne, Raymond Ergas, and Tamas Nemeth\n",
" SEG Technical Program Expanded Abstracts\n",
" https://library.seg.org/doi/10.1190/segam2016-13878451.1"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Introduction \n",
"\n",
"The goal of this tutorial set is to generate and prove correctness of modeling and inversion capability in Devito for variable density visco- acoustics using an energy conserving form of the wave equation. We describe how the linearization of the energy conserving *self adjoint* system with respect to modeling parameters allows using the same modeling system for all nonlinear and linearized forward and adjoint finite difference evolutions. There are three notebooks in this series:\n",
"\n",
"##### 1. Implementation of a Devito self adjoint variable density visco- acoustic isotropic modeling operator -- Nonlinear Ops\n",
"- Implement the nonlinear modeling operations. \n",
"- [sa_01_iso_implementation1.ipynb](sa_01_iso_implementation1.ipynb)\n",
"\n",
"##### 2. Implementation of a Devito self adjoint variable density visco- acoustic isotropic modeling operator -- Linearized Ops\n",
"- Implement the linearized (Jacobian) ```forward``` and ```adjoint``` modeling operations.\n",
"- [sa_02_iso_implementation2.ipynb](sa_02_iso_implementation2.ipynb)\n",
"\n",
"##### 3. Implementation of a Devito self adjoint variable density visco- acoustic isotropic modeling operator -- Correctness Testing\n",
"- Tests the correctness of the implemented operators.\n",
"- [sa_03_iso_correctness.ipynb](sa_03_iso_correctness.ipynb)\n",
"\n",
"There are similar series of notebooks implementing and testing operators for VTI and TTI anisotropy ([README.md](README.md)).\n",
"\n",
"Below we continue the implementation of our *self adjoint* wave equation with dissipation only Q attenuation, and linearize the modeling operator with respect to the model parameter velocity. We show how to implement finite difference evolutions to compute the action of the ```forward``` and ```adjoint``` Jacobian. \n",
"\n",
"## Outline \n",
"1. Define symbols \n",
"2. The nonlinear operator \n",
"3. The Jacobian opeator \n",
"4. Create the Devito grid and model fields \n",
"5. The simulation time range and acquistion geometry \n",
"6. Implement and run the nonlinear forward operator \n",
"7. Implement and run the Jacobian forward operator \n",
"8. Implement and run the Jacobian adjoint operator \n",
"9. References \n"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Table of symbols\n",
"\n",
"There are many more symbols for the Jacobian than we had in the previous notebook. We need to introduce terminology for the nonlinear operator, including total, background and perturbation fields for several variables. \n",
"\n",
"| Symbol | Description | Dimensionality | \n",
"|:---|:---|:---|\n",
"| $\\overleftarrow{\\partial_t}$ | shifted first derivative wrt $t$ | shifted 1/2 sample backward in time |\n",
"| $\\partial_{tt}$ | centered second derivative wrt $t$ | centered in time |\n",
"| $\\overrightarrow{\\partial_x},\\ \\overrightarrow{\\partial_y},\\ \\overrightarrow{\\partial_z}$ | + shifted first derivative wrt $x,y,z$ | shifted 1/2 sample forward in space |\n",
"| $\\overleftarrow{\\partial_x},\\ \\overleftarrow{\\partial_y},\\ \\overleftarrow{\\partial_z}$ | - shifted first derivative wrt $x,y,z$ | shifted 1/2 sample backward in space |\n",
"| $\\omega_c = 2 \\pi f$ | center angular frequency | constant |\n",
"| $b(x,y,z)$ | buoyancy $(1 / \\rho)$ | function of space |\n",
"| $Q(x,y,z)$ | Attenuation at frequency $\\omega_c$ | function of space |\n",
"| $m(x,y,z)$ | Total P wave velocity ($m_0+\\delta m$) | function of space |\n",
"| $m_0(x,y,z)$ | Background P wave velocity | function of space |\n",
"| $\\delta m(x,y,z)$ | Perturbation to P wave velocity | function of space |\n",
"| $u(t,x,y,z)$ | Total pressure wavefield ($u_0+\\delta u$)| function of time and space |\n",
"| $u_0(t,x,y,z)$ | Background pressure wavefield | function of time and space |\n",
"| $\\delta u(t,x,y,z)$ | Perturbation to pressure wavefield | function of time and space |\n",
"| $q(t,x,y,z)$ | Source wavefield | function of time, localized in space to source location |\n",
"| $r(t,x,y,z)$ | Receiver wavefield | function of time, localized in space to receiver locations |\n",
"| $F[m; q]$ | Forward nonlinear modeling operator | Nonlinear in $m$, linear in $q$: $\\quad$ maps $m \\rightarrow r$ |\n",
"| $\\nabla F[m; q]\\ \\delta m$ | Forward Jacobian modeling operator | Linearized at $[m; q]$: $\\quad$ maps $\\delta m \\rightarrow \\delta r$ |\n",
"| $\\bigl( \\nabla F[m; q] \\bigr)^\\top\\ \\delta r$ | Adjoint Jacobian modeling operator | Linearized at $[m; q]$: $\\quad$ maps $\\delta r \\rightarrow \\delta m$ |\n",
"| $\\Delta_t, \\Delta_x, \\Delta_y, \\Delta_z$ | sampling rates for $t, x, y , z$ | $t, x, y , z$ | "
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## A word about notation \n",
"\n",
"We use the arrow symbols over derivatives $\\overrightarrow{\\partial_x}$ as a shorthand notation to indicate that the derivative is taken at a shifted location. For example:\n",
"\n",
"- $\\overrightarrow{\\partial_x}\\ u(t,x,y,z)$ indicates that the $x$ derivative of $u(t,x,y,z)$ is taken at $u(t,x+\\frac{\\Delta x}{2},y,z)$.\n",
"\n",
"- $\\overleftarrow{\\partial_z}\\ u(t,x,y,z)$ indicates that the $z$ derivative of $u(t,x,y,z)$ is taken at $u(t,x,y,z-\\frac{\\Delta z}{2})$.\n",
"\n",
"- $\\overleftarrow{\\partial_t}\\ u(t,x,y,z)$ indicates that the $t$ derivative of $u(t,x,y,z)$ is taken at $u(t-\\frac{\\Delta_t}{2},x,y,z)$.\n",
"\n",
"We usually drop the $(t,x,y,z)$ notation from wavefield variables unless required for clarity of exposition, so that $u(t,x,y,z)$ becomes $u$."
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## The Nonlinear operator\n",
"\n",
"The nonlinear operator is the solution to the self adjoint scalar isotropic variable density visco- acoustic wave equation shown immediately below, and maps the velocity model vector $m$ into the receiver wavefield vector $r$.\n",
"\n",
"$$\n",
"\\frac{b}{m^2} \\left( \\frac{\\omega_c}{Q}\\overleftarrow{\\partial_t}\\ u + \\partial_{tt}\\ u \\right) =\n",
" \\overleftarrow{\\partial_x}\\left(b\\ \\overrightarrow{\\partial_x}\\ u \\right) +\n",
" \\overleftarrow{\\partial_y}\\left(b\\ \\overrightarrow{\\partial_y}\\ u \\right) +\n",
" \\overleftarrow{\\partial_z}\\left(b\\ \\overrightarrow{\\partial_z}\\ u \\right) + q\n",
"$$\n",
"\n",
"In operator notation, where the operator is nonlinear with respect to model $m$ to the left of semicolon inside the square brackets, and linear with respect to the source $q$ to the right of semicolon inside the square brackets.\n",
"\n",
"$$\n",
"F[m; q] = r\n",
"$$"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## The Jacobian operator\n",
"\n",
"In this section we linearize about a background model and take the derivative of the nonlinear operator to obtain the Jacobian forward operator. \n",
"\n",
"In operator notation, where the derivative of the modeling operator is now linear in the model perturbation vector $\\delta m$, the Jacobian operator maps a perturbation in the velocity model $\\delta m$ into a perturbation in the receiver wavefield $\\delta r$.\n",
"\n",
"$$\n",
"\\nabla F[m; q]\\ \\delta m = \\delta r\n",
"$$\n",
"\n",
"#### 1. We begin by simplifying notation \n",
"To simplify the treatment below we introduce the operator $L_t[\\cdot]$, accounting for the time derivatives inside the parentheses on the left hand side of the wave equation. \n",
"\n",
"$$\n",
"L_t[\\cdot] \\equiv \\frac{\\omega_c}{Q} \\overleftarrow{\\partial_t}[\\cdot] + \\partial_{tt}[\\cdot]\n",
"$$\n",
"\n",
"Next we re-write the wave equation using this notation.\n",
"\n",
"$$\n",
"\\frac{b}{m^2} L_t[u] =\n",
" \\overleftarrow{\\partial_x}\\left(b\\ \\overrightarrow{\\partial_x}\\ u \\right) +\n",
" \\overleftarrow{\\partial_y}\\left(b\\ \\overrightarrow{\\partial_y}\\ u \\right) +\n",
" \\overleftarrow{\\partial_z}\\left(b\\ \\overrightarrow{\\partial_z}\\ u \\right) + q\n",
"$$\n",
"\n",
"#### 2. Linearize\n",
"To linearize we treat the total model as the sum of background and perturbation models $\\left(m = m_0 + \\delta m\\right)$, and the total pressure as the sum of background and perturbation pressures $\\left(u = u_0 + \\delta u\\right)$.\n",
"\n",
"$$\n",
"\\frac{b}{(m_0+\\delta m)^2} L_t[u_0+\\delta u] =\n",
" \\overleftarrow{\\partial_x}\\left(b\\ \\overrightarrow{\\partial_x} (u_0+\\delta u) \\right) +\n",
" \\overleftarrow{\\partial_y}\\left(b\\ \\overrightarrow{\\partial_y} (u_0+\\delta u) \\right) +\n",
" \\overleftarrow{\\partial_z}\\left(b\\ \\overrightarrow{\\partial_z} (u_0+\\delta u) \\right) + q\n",
"$$\n",
"\n",
"Note that *model parameters* for this variable density isotropic visco-acoustic physics is only velocity, we do not treat perturbations to density.\n",
"\n",
"We also write the PDE for the background model, which we subtract after linearization to simplify the final expression. \n",
"\n",
"$$\n",
"\\frac{b}{m_0^2} L_t[u_0] =\n",
" \\overleftarrow{\\partial_x}\\left(b\\ \\overrightarrow{\\partial_x} u_0 \\right) +\n",
" \\overleftarrow{\\partial_y}\\left(b\\ \\overrightarrow{\\partial_y} u_0 \\right) +\n",
" \\overleftarrow{\\partial_z}\\left(b\\ \\overrightarrow{\\partial_z} u_0 \\right) + q\n",
"$$\n",
"\n",
"#### 3. Take derivative w.r.t. model parameters\n",
"Next we take the derivative with respect to velocity, keep only terms up to first order in the perturbations, subtract the background model PDE equation, and finally arrive at the following linearized equation:\n",
"\n",
"$$\n",
"\\frac{b}{m_0^2} L_t\\left[\\delta u\\right] =\n",
" \\overleftarrow{\\partial_x}\\left(b\\ \\overrightarrow{\\partial_x} \\delta u \\right) +\n",
" \\overleftarrow{\\partial_y}\\left(b\\ \\overrightarrow{\\partial_y} \\delta u \\right) +\n",
" \\overleftarrow{\\partial_z}\\left(b\\ \\overrightarrow{\\partial_z} \\delta u \\right) + \n",
" \\frac{2\\ b\\ \\delta m}{m_0^3} L_t\\left[u_0\\right]\n",
"$$\n",
"\n",
"Note that the source $q$ in the original equation above has disappeared due to subtraction of the background PDE, and has been replaced by the *Born source*: \n",
"\n",
"$$\n",
"q = \\frac{\\displaystyle 2\\ b\\ \\delta m}{\\displaystyle m_0^3} L_t\\left[u_0\\right]\n",
"$$\n",
"\n",
"This is the same equation as used for the nonlinear forward, only now in the perturbed wavefield $\\delta u$ and with the Born source."
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## The adjoint of the Jacobian operator\n",
"\n",
"In this section we introduce the adjoint of the Jacobian operator we derived above. The Jacobian adjoint operator maps a perturbation in receiver wavefield $\\delta r$ into a perturbation in velocity model $\\delta m$. In operator notation:\n",
"\n",
"$$\n",
"\\bigl( \\nabla F[m; q] \\bigr)^\\top\\ \\delta r = \\delta m\n",
"$$\n",
"\n",
"#### 1. Solve the time reversed wave equation with the receiver perturbation as source\n",
"The PDE for the adjoint of the Jacobian is solved for the perturbation to the pressure wavefield $\\delta u$ by using the same wave equation as the nonlinear forward and the Jacobian forward, with the time reversed receiver wavefield perturbation $\\widetilde{\\delta r}$ injected as source. \n",
"\n",
"Note that we use $\\widetilde{\\delta u}$ and $\\widetilde{\\delta r}$ to indicate that we solve this finite difference evolution time reversed.\n",
"\n",
"$$\n",
"\\frac{b}{m_0^2} L_t\\left[\\widetilde{\\delta u}\\right] =\n",
" \\overleftarrow{\\partial_x}\\left(b\\ \\overrightarrow{\\partial_x}\\ \\widetilde{\\delta u} \\right) +\n",
" \\overleftarrow{\\partial_y}\\left(b\\ \\overrightarrow{\\partial_y}\\ \\widetilde{\\delta u} \\right) +\n",
" \\overleftarrow{\\partial_z}\\left(b\\ \\overrightarrow{\\partial_z}\\ \\widetilde{\\delta u} \\right) + \n",
" \\widetilde{\\delta r}\n",
"$$\n",
"\n",
"#### 2. Compute zero lag correlation \n",
"\n",
"We compute the perturbation to the velocity model by zero lag correlation of the wavefield perturbation $\\widetilde{\\delta u}$ solved in step 1 as shown in the following expression: \n",
"\n",
"$$\n",
"\\delta m(x,y,z) = \\sum_t \n",
" \\left\\{ \n",
" \\widetilde{\\delta u}(t,x,y,z)\\ \n",
" \\frac{\\displaystyle 2\\ b}{\\displaystyle m_0^3} L_t\\bigl[u_0(t,x,y,z)\\bigr]\n",
" \\right\\}\n",
"$$\n",
"\n",
"Note that this correlation can be more formally derived by examining the equations for two Green's functions, one for the background model ($m_0$) and wavefield ($u_0$), and one for for the total model $(m_0 + \\delta m)$ and wavefield $(u_0 + \\delta u)$, and subtracting to derive the equation for Born scattering.\n"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Implementation\n",
"\n",
"Next we assemble the Devito objects needed to implement these linearized operators."
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Imports \n",
"\n",
"We have grouped all imports used in this notebook here for consistency."
]
},
{
"cell_type": "code",
"execution_count": 1,
"metadata": {},
"outputs": [],
"source": [
"import numpy as np\n",
"from examples.seismic import RickerSource, Receiver, TimeAxis\n",
"from devito import (Grid, Function, TimeFunction, SpaceDimension, Constant, \n",
" Eq, Operator, solve, configuration, norm)\n",
"from devito.finite_differences import Derivative\n",
"from devito.builtins import gaussian_smooth\n",
"from examples.seismic.self_adjoint import setup_w_over_q\n",
"import matplotlib as mpl\n",
"import matplotlib.pyplot as plt\n",
"from matplotlib import cm\n",
"from timeit import default_timer as timer\n",
"\n",
"# These lines force images to be displayed in the notebook, and scale up fonts \n",
"%matplotlib inline\n",
"mpl.rc('font', size=14)\n",
"\n",
"# Make white background for plots, not transparent\n",
"plt.rcParams['figure.facecolor'] = 'white'\n",
"\n",
"# We define 32 bit floating point as the precision type \n",
"dtype = np.float32\n",
"\n",
"# Set logging to debug, captures statistics on the performance of operators\n",
"# configuration['log-level'] = 'DEBUG'\n",
"configuration['log-level'] = 'INFO'"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Instantiate the Devito grid for a two dimensional problem\n",
"\n",
"We define the grid the same as in the previous notebook outlining implementation for the nonlinear forward."
]
},
{
"cell_type": "code",
"execution_count": 2,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"shape; (301, 301)\n",
"origin; (0.0, 0.0)\n",
"spacing; (10.0, 10.0)\n",
"extent; (3000.0, 3000.0)\n",
"\n",
"shape_pad; [401 401]\n",
"origin_pad; (-500.0, -500.0)\n",
"extent_pad; (4000.0, 4000.0)\n",
"\n",
"grid.shape; (401, 401)\n",
"grid.extent; (4000.0, 4000.0)\n",
"grid.spacing_map; {h_x: 10.0, h_z: 10.0}\n"
]
}
],
"source": [
"# Define dimensions for the interior of the model\n",
"nx,nz = 301,301\n",
"dx,dz = 10.0,10.0 # Grid spacing in m\n",
"shape = (nx, nz) # Number of grid points\n",
"spacing = (dx, dz) # Domain size is now 5 km by 5 km\n",
"origin = (0., 0.) # Origin of coordinate system, specified in m.\n",
"extent = tuple([s*(n-1) for s, n in zip(spacing, shape)])\n",
"\n",
"# Define dimensions for the model padded with absorbing boundaries\n",
"npad = 50 # number of points in absorbing boundary region (all sides)\n",
"nxpad,nzpad = nx + 2 * npad, nz + 2 * npad\n",
"shape_pad = np.array(shape) + 2 * npad\n",
"origin_pad = tuple([o - s*npad for o, s in zip(origin, spacing)])\n",
"extent_pad = tuple([s*(n-1) for s, n in zip(spacing, shape_pad)])\n",
"\n",
"# Define the dimensions \n",
"# Note if you do not specify dimensions, you get in order x,y,z\n",
"x = SpaceDimension(name='x', spacing=Constant(name='h_x', \n",
" value=extent_pad[0]/(shape_pad[0]-1)))\n",
"z = SpaceDimension(name='z', spacing=Constant(name='h_z', \n",
" value=extent_pad[1]/(shape_pad[1]-1)))\n",
"\n",
"# Initialize the Devito grid \n",
"grid = Grid(extent=extent_pad, shape=shape_pad, origin=origin_pad, \n",
" dimensions=(x, z), dtype=dtype)\n",
"\n",
"print(\"shape; \", shape)\n",
"print(\"origin; \", origin)\n",
"print(\"spacing; \", spacing)\n",
"print(\"extent; \", extent)\n",
"print(\"\")\n",
"print(\"shape_pad; \", shape_pad)\n",
"print(\"origin_pad; \", origin_pad)\n",
"print(\"extent_pad; \", extent_pad)\n",
"\n",
"print(\"\")\n",
"print(\"grid.shape; \", grid.shape)\n",
"print(\"grid.extent; \", grid.extent)\n",
"print(\"grid.spacing_map;\", grid.spacing_map)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Define velocity, buoyancy and $\\frac{\\omega_c}{Q}$ model parameters\n",
"\n",
"We have the following constants and fields to define:\n",
"\n",
"| Symbol | Description |\n",
"|:---:|:---|\n",
"| $$m0(x,z)$$ | Background velocity model |\n",
"| $$\\delta m(x,z)$$ | Perturbation to velocity model |\n",
"| $$b(x,z)=\\frac{1}{\\rho(x,z)}$$ | Buoyancy (reciprocal density) |\n",
"| $$\\omega_c = 2 \\pi f_c$$ | Center angular frequency |\n",
"| $$\\frac{1}{Q(x,z)}$$ | Inverse Q model used in the modeling system |\n"
]
},
{
"cell_type": "code",
"execution_count": 3,
"metadata": {},
"outputs": [
{
"name": "stderr",
"output_type": "stream",
"text": [
"Operator `WOverQ_Operator` ran in 0.01 s\n"
]
}
],
"source": [
"# NBVAL_IGNORE_OUTPUT\n",
"\n",
"# Create the velocity and buoyancy fields as in the nonlinear notebook \n",
"space_order = 8\n",
"\n",
"# Wholespace velocity\n",
"m0 = Function(name='m0', grid=grid, space_order=space_order)\n",
"m0.data[:] = 1.5\n",
"\n",
"# Perturbation to velocity: a square offset from the center of the model\n",
"dm = Function(name='dm', grid=grid, space_order=space_order)\n",
"size = 10\n",
"x0 = shape_pad[0]//2\n",
"z0 = shape_pad[1]//2\n",
"dm.data[:] = 0.0\n",
"dm.data[x0-size:x0+size, z0-size:z0+size] = 1.0\n",
"\n",
"# Constant density\n",
"b = Function(name='b', grid=grid, space_order=space_order)\n",
"b.data[:,:] = 1.0 / 1.0\n",
"\n",
"# Initialize the attenuation profile for Q=100 model\n",
"fpeak = 0.010\n",
"w = 2.0 * np.pi * fpeak\n",
"qmin = 0.1\n",
"qmax = 1000.0\n",
"wOverQ = Function(name='wOverQ', grid=grid, space_order=space_order)\n",
"setup_w_over_q(wOverQ, w, qmin, 100.0, npad)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Define the simulation time range and the acquisition geometry \n",
"\n",
"#### Simulation time range: \n",
"\n",
"In this notebook we run 3 seconds of simulation using the sample rate related to the CFL condition as implemented in ```examples/seismic/self_adjoint/utils.py```. \n",
"\n",
"We also use the convenience ```TimeRange``` as defined in ```examples/seismic/source.py```.\n",
"\n",
"#### Acquisition geometry: \n",
"\n",
"**source**:\n",
"- X coordinate: left sode of model\n",
"- Z coordinate: middle of model\n",
"- We use a 10 Hz center frequency [RickerSource](https://github.com/devitocodes/devito/blob/master/examples/seismic/source.py#L280) wavelet source as defined in ```examples/seismic/source.py```\n",
"\n",
"**receivers**:\n",
"- X coordinate: right side of model\n",
"- Z coordinate: vertical line in model\n",
"- We use a vertical line of [Receivers](https://github.com/devitocodes/devito/blob/master/examples/seismic/source.py#L80) as defined with a ```PointSource``` in ```examples/seismic/source.py```"
]
},
{
"cell_type": "code",
"execution_count": 4,
"metadata": {},
"outputs": [],
"source": [
"def compute_critical_dt(v):\n",
" \"\"\"\n",
" Determine the temporal sampling to satisfy CFL stability.\n",
" This method replicates the functionality in the Model class.\n",
" Note we add a safety factor, reducing dt by a factor 0.75 due to the\n",
" w/Q attentuation term.\n",
" Parameters\n",
" ----------\n",
" v : Function\n",
" velocity\n",
" \"\"\"\n",
" coeff = 0.38 if len(v.grid.shape) == 3 else 0.42\n",
" dt = 0.75 * v.dtype(coeff * np.min(v.grid.spacing) / (np.max(v.data)))\n",
" return v.dtype(\"%.5e\" % dt)"
]
},
{
"cell_type": "code",
"execution_count": 5,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Time min, max, dt, num; 0.000000 1200.000000 2.100000 572\n",
"time_range; TimeAxis: start=0, stop=1201.2, step=2.1, num=573\n",
"src_coordinate X; +752.0000\n",
"src_coordinate Z; +1505.0000\n",
"rec_coordinates X min/max; +2257.0000 +2257.0000\n",
"rec_coordinates Z min/max; +0.0000 +3000.0000\n"
]
}
],
"source": [
"t0 = 0.0 # Simulation time start\n",
"tn = 1200.0 # Simulation time end (1 second = 1000 msec)\n",
"dt = compute_critical_dt(m0)\n",
"time_range = TimeAxis(start=t0, stop=tn, step=dt)\n",
"print(\"Time min, max, dt, num; %10.6f %10.6f %10.6f %d\" % (t0, tn, dt, int(tn//dt) + 1))\n",
"print(\"time_range; \", time_range)\n",
"\n",
"# Source at 1/4 X, 1/2 Z, Ricker with 10 Hz center frequency\n",
"src_nl = RickerSource(name='src_nl', grid=grid, f0=fpeak, npoint=1, time_range=time_range)\n",
"src_nl.coordinates.data[0,0] = dx * 1 * nx//4\n",
"src_nl.coordinates.data[0,1] = dz * shape[1]//2\n",
"\n",
"# Receivers at 3/4 X, line in Z\n",
"rec_nl = Receiver(name='rec_nl', grid=grid, npoint=nz, time_range=time_range)\n",
"rec_nl.coordinates.data[:,0] = dx * 3 * nx//4\n",
"rec_nl.coordinates.data[:,1] = np.linspace(0.0, dz*(nz-1), nz)\n",
"\n",
"print(\"src_coordinate X; %+12.4f\" % (src_nl.coordinates.data[0,0]))\n",
"print(\"src_coordinate Z; %+12.4f\" % (src_nl.coordinates.data[0,1]))\n",
"print(\"rec_coordinates X min/max; %+12.4f %+12.4f\" % \\\n",
" (np.min(rec_nl.coordinates.data[:,0]), np.max(rec_nl.coordinates.data[:,0])))\n",
"print(\"rec_coordinates Z min/max; %+12.4f %+12.4f\" % \\\n",
" (np.min(rec_nl.coordinates.data[:,1]), np.max(rec_nl.coordinates.data[:,1])))"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Plot the model \n",
"\n",
"We plot the following ```Functions```:\n",
"- Background Velocity\n",
"- Background Density\n",
"- Velocity perturbation\n",
"- Q Model\n",
"\n",
"Each subplot also shows:\n",
"- The location of the absorbing boundary as a dotted line\n",
"- The source location as a red star\n",
"- The line of receivers as a black vertical line"
]
},
{
"cell_type": "code",
"execution_count": 6,
"metadata": {},
"outputs": [
{
"data": {
"image/png": "iVBORw0KGgoAAAANSUhEUgAAA1gAAANICAYAAADaSU4zAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjQuMSwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/Z1A+gAAAACXBIWXMAAAsTAAALEwEAmpwYAAEAAElEQVR4nOzdeVwU5R8H8M+yXMvpwe2Bmor3Bal45wFeaGalqYimQmWpaV5lgWVapnaYWvjzyrO8M8vEvCCwBDVR8BbNuMQTLwT2+f1hTA6Hciy7O8vn/Xrta9mZZ2e+Mwz75dnnGJUQQoCIiIiIiIjKzMzQARAREREREZkKVrCIiIiIiIh0hBUsIiIiIiIiHWEFi4iIiIiISEdYwSIiIiIiItIRVrCIiIiIiIh0hBUsMjpdunRBgwYNDB2GwXTp0gVdunTR+373798PlUqF/fv3K3L7RETlibnJMLmpvKxcuRIqlQpJSUmGDoVMECtYVKi8D57HH87OzujUqRO2bdtm6PAqtPT0dFhYWGDQoEFFlsnOzoaTkxPatWunx8hK7ueff0ZYWJihwyAihWBuMn5dunSRfjdmZmZwcHCAl5cXAgMDERERYejwnmjdunX44osvDB0GmQBzQwdAxi0sLAzPPPMMhBBIT0/HmjVrMGDAAGzYsOGJ/+BT+XFxcUGPHj2wY8cO3LlzB3Z2dgXK7Nq1C9euXcOwYcMMEGHhOnXqhPv378PS0lJa9vPPP2PRokWsZBFRiTA3GTd3d3fMnTsXAHD37l2cO3cOW7ZswZo1a/Dyyy9jzZo1sLCwMGiMgYGBGDx4MKysrKRl69atw4kTJzBhwgTDBUYmgRUseiJ/f3+0bdtWeh0SEgIPDw+sW7dO0Uns7t27sLW1NXQYpTZs2DD88ssv2Lp1KwIDAwusX7t27VNbufTNzMwM1tbWhg6DiEwAc5Nxc3BwKPAF3yeffIJx48Zh8eLFqFWrFj799FMDRfeIWq2GWq02aAxkuthFkErEzs4OdnZ2MDeX183nz5+PDh06wMnJCdbW1mjatCn+97//FbqNiIgIdO3aFQ4ODrC3t4e3t3eRZfMcPHgQDg4OCAgIQFZWFgDg/v37GDduHJycnGBvb49+/frhypUrUKlUshaRsLAwqFQqnDhxAoGBgahSpQqaNGkirf/mm2/QpEkTWFtbw83NDSEhIbh+/bps/7Vq1cKIESMKxJW/T3reOKP169dj9uzZqF69OqytrdGtWzecO3euwPvDw8PxzDPPQKPRoHXr1oiMjHziecjz/PPPw87ODuvWrSuwLjMzEz/++CN69uyJqlWrAgDOnDmDl19+GVWrVoW1tTVatmyJTZs2FWtfBw4cQOfOnWFrawtHR0f07dsXJ06cKFAuJSUFISEhqF69OqysrFCrVi2MGTMGmZmZsnOTNwZrxIgRWLRoEQDIuvskJSWhffv2aNasWaHxtGrVCm3atClW7ERUMTA3yRkqNz2JWq3GV199hUaNGuHrr7/GrVu3ZOvXrVuHZ599FhqNBpUrV8ZLL72EixcvFjiuBg0aICEhAV27doWNjQ2qVasmtZY9bvHixWjatKmUu1q0aIFvv/1WWp9/DFaXLl2wc+dOXLp0SZaTtFotatasiX79+hXYR05ODlxdXRVdqafywRYseqJbt24hIyMDAHD16lV8++23SE1NxfDhw2XlPv/8c/Tt2xcvv/wyVCoVtm/fjjFjxiAnJwevvfaaVG716tUICgpCw4YNMWXKFFStWhXHjx/Hzp07MXr06EJjiIiIwPPPP4/evXtj3bp1UreCESNG4IcffsCwYcPg6+uLAwcOoE+fPkUey6BBg1C7dm3MmjULDx8+BADMmjUL77//Prp27YqQkBCcP38eixYtwh9//IE//vhD1nWgJObOnQu1Wo133nkHt27dwty5czF06FD88ccfUplly5YhJCQE7dq1w/jx43Hp0iX0798flStXRo0aNZ64fRsbGwwYMADr169Heno6XFxcpHVbt27F/fv3MXToUABAYmIi2rVrBzc3N0yZMgV2dnbYunUrXnrpJaxevfqJ3Qj37dsHPz8/1K5dG2FhYXjw4AEWLVqE9u3b4/Dhw6hfvz4AIDU1Fa1bt0ZGRgaCg4PRuHFjJCcnY+vWrbh27Rrs7e0LbDskJATJycmIiIjA6tWrpeXOzs4ICgpCSEgIjh8/LqtoJSYm4ujRo/j666+f8hsgIlPG3GScuelp1Go1XnnlFbz//vuIioqSzssnn3yCd999Fy+++CJGjhyJGzdu4Ouvv0b79u3x119/wdnZWdrGrVu30KtXLwwYMAAvvfQSNm3ahKlTp6Jp06bo1auXdAxjx47Fiy++iDfffBPZ2dk4efIkoqOjERISUmhs7733Hm7duoUrV67g888/l5abmZlh2LBhmDdvHq5duyZ9cQkAu3fvRnp6eoHrjgiCqBArVqwQAAo8LC0tRXh4eIHyd+/eLbCsR48e4plnnpFe37p1Szg4OAhvb29x7949WVmtViv93LlzZ+Hl5SWEEOLHH38UVlZWIjAwUOTk5Ehl4uLiBADx5ptvyrYzYsQIAUCEhoZKy0JDQwUA8cILL8jKpqenC0tLS9GtWzfZtvOOfeHChdIyT09PERQUVOAYO3fuLDp37iy93rdvnwAgGjRoILKysqTlX375pQAg4uPjhRBCPHz4ULi4uIgWLVrIyi1btkwAkG2zKL/++muBOIUQws/PTzg4OEjnuEePHqJRo0YFznmPHj1EtWrVpHOfF/u+ffukMi1bthRVq1YVGRkZ0rIzZ84ICwsLMXDgQGlZUFCQMDMzE4cOHSoQ55O2P3bsWFHYx9CNGzeEtbW1mDx5smz59OnThYWFhbh69eqTTg0RmSjmJuPPTY+fp8Js3bpVABBffvmlEEKIS5cuCXNzczFz5kxZuXPnzgkrKysxffp02bYBiFWrVknLsrKyhJubmywnPf/886Jx48ZPjDPvfF68eFFa1qdPH+Hp6VmgbGJiogAgFi1aJFv+yiuvCGdnZ5Gdnf3EfVHFwy6C9ERfffUVIiIiEBERgTVr1qB79+54/fXX8cMPP8jK2djYAHg0e93169eRkZGB5557DufPn5e6AezevRu3b9/GtGnToNFoZO9XqVQF9r1p0yYMHDgQQUFBWLlypayv9K5duwAAb7zxhuw9b731VpHH8vrrr8te79mzBw8fPsT48eNl2w4MDISrqyt27txZ5LaeZvjw4bLJHDp27AgAuHDhAgAgNjYW6enpGDNmjKzc8OHDUalSpWLto1u3bnB3d5d1E0xLS8Nvv/2GF154ARqNBtevX8eePXvw8ssv4+7du8jIyJAePXv2xD///IMzZ84Uuv2UlBQcPXoUQUFBsm/s6tWrh379+mHXrl3Izc2FVqvF1q1b0atXr0K77hX2u32aSpUqoV+/fli3bh20Wi0AQAiBdevWoVevXnBycirxNonIdDA3lY4+ctPT5E3MlNd9fMuWLcjJycGgQYNkOcrR0RFNmzbFvn37ZO/XaDSynheWlpZo3bq1dAwA4OjoiCtXruDw4cM6iblBgwZo3bq1rLfFnTt3sH37drzyyisFuqYSsYJFT/Tss8+ie/fu6N69O4YOHYodO3agadOmGDdunNSVAQC2b98OHx8faDQaVK1aFc7Oznj33XcBQEpi58+fBwBZH/OiXL58GYMHD0a/fv3w7bffwsxMfqnm9ZF+5plnZMvr1q1b5Dbzl7106RIAwMvLS7ZcrVajXr16Zbo3Rs2aNWWvK1euDAC4ceOGbN/16tWTlTM3N0ft2rWLtQ+1Wo3BgwcjJiZG6qe+YcMG5ObmSt0Dz507ByEEwsLC4OzsLHtMmjQJwKNp3wtT1PkBgIYNG0oVtqtXr+L27dvF+r2WRFBQEP755x8puUZGRuLSpUuFTupBRBULc1Pp6CM3Pc2dO3cAQOo6nvclX4MGDQrkqbwK3+OqVatW4LxXrlxZOgYAmDp1Kuzt7dG6dWs888wzeO211wpU1EoqKCgIhw4dkq6XLVu24N69e8xJVChWsKhEzMzM0KVLF6SlpeHs2bMAgKioKAwYMAA2Njb45ptvsHPnTkRERODtt98GAKkFoiRcXV3RoUMH/Prrr4iJidFJ7Pm/mSyJolphcnNzC11e1MxEQohSx1CYvG/x8lqx1q5dCw8PD3Tt2hXAf+f+7bfflr7tzf/QdcVIV/z9/eHq6oo1a9YAANasWYNKlSohICDAwJERkbFhbpIzdG56krxJkvIqnXm/h19++aXQHPXdd9/J3l+cY2jYsCFOnz6NjRs3omvXrvjpp5/QtWtXjB07ttRxDx48GJaWlrKc1KBBA/j4+JR6m2S6WMGiEsvOzgbw37dQmzZtgrW1NXbv3o3Ro0ejd+/e6N69e4GkkfctXWEz0OVnZWWFHTt2oFGjRujduzeOHTsmW+/p6QkhhPRNUp7CZkMqiqenJwDg9OnTsuVarRZnz55FrVq1pGWVK1fGzZs3C2wj79u+ksrbd94/AnlycnIKzJr0JK1atULDhg2xbt06nD17FocPH8Yrr7wifbtXp04dAI++fcz7tjf/I+8bzKJizH9+AODUqVOwtbWFk5MTnJ2d4eDgUKzfa35P6j6oVqsxdOhQbN68GTdv3sSmTZvw0ksvlXpwNxGZNuam/xg6NxUlNzcX69atg42NDTp06ADgv/Nfs2bNQnNU+/btS7UvGxsbvPjii1i6dCmSkpIwdOhQLF68GP/880+R73lSTqpSpQr69u2LNWvWICUlBXv37mXrFRWJFSwqkezsbERERMDS0hINGzYE8Ogf4bypTPPcuHEDy5cvl73Xz88PDg4O+OSTT3D//n3ZusK+PbO3t8euXbtQs2ZN+Pn54dSpU9I6f39/AI+mYX3cwoULi30sPXr0gKWlJb766itZ7GvXrkVaWhr69u0rLXvmmWdw6NAhWdeTn376CX///Xex9/c4Hx8fODs7Y+nSpbJtfvfdd4UmyycZNmwYEhISMHnyZACQugcCj25K/Nxzz2Hp0qWFJpWrV68WuV13d3e0atUK3333nWxq4PPnz+PHH39Er169oFarYWZmhgEDBuCXX36RzUSV50nfjObd7+Xxrh2PCwoKQmZmJkJCQnDjxg3O1EREhWJuMr7clF9ubi7GjRuHxMREjBs3Dg4ODgCAgQMHQq1W48MPPyz0fOfNFlkS165dk702NzdH06ZNAeCJx2Fra4ubN28WmbeCgoJw7tw5TJgwAVqt9omz8FLFxlF59ES//vqr9M1beno6NmzYgDNnzmDatGnSh2NAQAAWLFiAHj16IDAwENevX8fSpUvh5uaG1NRUaVsODg748ssv8eqrr8LHxwdDhgxB1apVcfLkSfzzzz/YsmVLgf1XrlwZu3fvRqdOndC9e3dERkaidu3a8Pb2xsCBA6V7abRt2xYHDhyQ+nIXZ2IFJycnvP/++3j//ffh5+eH559/HhcuXMDXX3+N5s2by6bmHT16NDZt2oSePXvi5Zdfxvnz57FmzZoCfeeLy8LCArNmzUJISAiee+45DB48GElJSVixYoXU6lRcQ4cOxYwZM7B9+3Y0atQILVu2lK1fsmSJdF+pMWPG4JlnnkF6ejr++OMPJCQkPPGb1Xnz5sHPzw++vr4YM2aMNE27tbU1Pv74Y6ncnDlzEBERgS5duiAkJASNGjVCWloatmzZgq1bt8q+cX1cXteKN998E7169YK5uTkCAgKkilezZs3QvHlz/PDDD6hdu3apv8kkItPC3PSIseam27dvS13p7t27h3PnzmHLli04f/48Bg8ejI8++kgqW6dOHXzyySeYPHkyLl26hOeffx6VKlXCxYsXsX37dgwaNEh2/7Di8PPzg4uLCzp06AA3NzecO3cOCxcuRLNmzaQKeGF8fHzw/fffY8KECWjTpg3MzMwwePBgaX2vXr3g7OyMH374AV26dCkwpo1IYpjJC8nYFTYVrrW1tWjRooVYsmSJbOpaIYRYtWqVaNCggbCyshLPPPOMmD9/vli+fHmBKVCFEGLnzp2iQ4cOwsbGRtjb2wtvb2+xfPlyaX1hU7xevnxZeHp6itq1a4srV64IIR5Nvzt27FhRpUoVYWdnJ55//nlx+vRpAUB88skn0nvzpsJNSUkp9FiXLFkiGjVqJCwtLYWLi4sYM2aMbFryPPPnzxfVqlUTVlZWon379iI2NrbIqXDXr18ve+/FixcFALFixQrZ8sWLF4vatWsLKysr4ePjIw4ePFhgm8XRsWNHAUDMnj270PUXL14UI0aMEO7u7sLCwkJ4eHiIXr16iXXr1hWI/fFp1POWd+zYUWg0GmFvby969+4tjh8/XmAff//9txgxYoRwcXERlpaWolatWiI4OFhkZmYWuf3c3FwxYcIE4erqKlQqVaHXy4IFCwQA8f7775fonBCR6WFuMv7clDeVet7Dzs5O1KtXTwwbNkzs3r27yPdt27ZNdOrUSdjZ2QkbGxtRv3598dprr4kTJ07Itl3YFPBBQUGy6dW//fZb0blzZ+Hk5CTlo7feekukpaVJZQqbpv3u3bti+PDhokqVKlJOym/cuHECgFi2bNlTzwVVXCoh9DiykaicHTt2DC1btsSaNWtkXeVIuRYtWoQ333wTp0+flm5sTESkJMxNpmPy5Mn4+uuvkZaWJrWWEuXHMVikWPn7ygPAF198ATMzM3Tq1MkAEVF5+N///gdfX19WrohIEZibTFdWVhZWr16NAQMGsHJFT2TSFawuXbpApVLJHo/3pQUeDXgNDAyEo6MjHB0dERgYWGAAZHx8PDp37gyNRoNq1aoVORCT9Gvu3Lno168fPv/8cyxcuBC9e/fGqlWrMHr0aNSoUcPQ4VEZ3L17F+vXr8drr72GY8eOSffsIlIy5qSKgbnJ9KSnp2PdunUYNmwY0tPTpan+iYpi8pNcjBw5ErNnz5Ze55+edciQIbh8+bJ09/XRo0cjMDAQO3bsAPBooGaPHj3QqVMnHD58GKdOncLIkSNha2vLf/oMrF27doiIiMBHH32EO3fuoGbNmggLC8N7771n6NCojK5evYohQ4agUqVKmDJlCgYOHGjokIh0gjnJ9DE3mZ6EhAQMHToUzs7O+Pzzz/Hss88aOiQyciY9BqtLly5o0qQJvv7660LXJyYmolGjRoiKipJmJ4uKikLHjh1x6tQpeHl5YcmSJZg6dSrS0tKkRDhr1iwsWbIEV65cKdaMQERERMxJREQVg0l3EQSADRs2wMnJCY0bN8Y777yDzMxMaV1MTAzs7OzQrl07aVn79u1ha2uL6OhoqUzHjh1l3zL6+/sjOTkZSUlJejsOIiJSPuYkIiLTZ9JdBIcMGQJPT094eHjg5MmTmD59Oo4fP47du3cDAFJTU+Hs7Cz7xk+lUsHFxUW6R0ZqaiqqV68u266rq6u0rnbt2rJ14eHhCA8PBwDExR0H4FReh0dEClW16sNS3TyTlM0QOQmQ56Xjp47DqUHVcjk+IlKmm0m3cC/jnqHDMCmKq2DNmDFDdoPTwuzbtw9dunRBcHCwtKxp06aoU6cO2rRpgyNHjqBVq1blEl9wcLC0X5XKA0Dwk99ARBVOrVo/GToE0hFjz0mAPC95+LhjVOzIctsXESnPMp8Vhg7B5CiugjVhwgQMGzbsiWWKurO2j48P1Go1zp49i1atWsHNzQ1Xr16FEEL6xlAIgfT0dLi5uQEA3NzckJaWJttO3uu8MkREVDExJxERUX6Kq2A5OTnByal03e7i4+ORm5sLd3d3AICvry/u3LmDmJgYqc97TEwM7t69K7329fXF1KlT8eDBA1hbWwMAIiIi4OHhgVq1apX9gIiISLGYk4iIKD+TneTi/Pnz+PDDDxEbG4ukpCT8/PPPGDx4MFq2bCnNztSwYUP07NkTISEhiImJQUxMDEJCQtC3b194eXkBeNRn3sbGBiNGjMCJEyewZcsWfPLJJ5g4cSJnayIiomJhTiIiqjhMtoJlaWmJ3377Df7+/vDy8sK4cePg5+eHPXv2QK1WS+XWrVuH5s2bw9/fH/7+/mjevDlWr14trXd0dERERASSk5Ph4+ODsWPHYtKkSZg4caIhDouIiBSIOYmIqOIw6ftgGRonuSCiwnh7/4TY2FhDh0EVECe5IKL8lvmsQHJsiqHDMCkm24JFRERERESkb6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgFdPixYtRu3ZtWFtbw9vbG5GRkYYOiYiIKjDmJSIi48QKVjF8//33GD9+PN59910cPXoU7dq1Q69evXD58mVDh0ZERBUQ8xIRkfFiBasYFixYgBEjRmDMmDFo2LAhFi5cCHd3dyxZssTQoRERUQXEvEREZLzMDR2AsXv48CHi4uLwzjvvyJb7+fkhOjraQFEVX1iYSvo5Lk5gx47CywUHAx4e/5X99luBlJSC5dzdgZCQ/8olJwuEhxe+zYAAwNv7v7I7dgjExT09zkevRaHlvL2BgAAeUx4eU+HbLMkxESmN0vOS+Uwr6Wdtq1xoA3IKLacOt4Aq5b/vgXPGPAQ8CvnMSVbBfKml9FK4a5EbnF3oNs12mMPsiFp6nds3G8Jb+9Q4ASAnNKvQcqo4M6h/spBe85h4TIUpyTGR8rEF6ykyMjKQm5sLV1dX2XJXV1ekpqYaKCoiQ0jBnDlzeN0TGRjzEtEjqcfSmJfIKLGCpWPh4eHw8fGBj48PgHsGiyMsTFWgtYGobDYhKysLGzduNHQgRFQCj+elu1cNl5fMZ1oVaG0gKostg7YyL5FRYgXrKZycnKBWq5GWliZbnpaWBjc3twLlg4ODERsbi9jYWAA2eoqSSB+uAwCuXbtm4DiIKray5CVbZ+YlMh3Xz94AwLxExkclhCh8EAVJ2rRpg+bNmyP8sQEf9evXx8CBAzFnzpwi36dSeQAI1kOEBRV3rAxR8X0EQItH38u8b+BYlM3b+6d/v4QhKp3S5iUPH3eMih2pjxALKO5YGaLimm31KbQPtTCzNMO7WVMNHY5iLfNZgeTYQgY/U6lxkotimDhxIgIDA9G6dWu0b98e33zzDZKTk/Haa68ZOjQiPdLmeyYiQ2FeIgK02VrZM5GxYAWrGAYNGoRr165h1qxZSElJQZMmTfDzzz/D09PT0KEViS1WpHtm+K8FS3nYqkumRIl5iS1WpGtmFmaPWrAslJmX2KpruljBKqY33ngDb7zxhqHDIDIgtmARGRPmJaro2IJFxkqZVX4iMgCzfM9ERESGk9dypdQWLDJdvCKJqJjYgkVERMaDLVhkrNhFkIiKSeljsDjmiojIlCh9DBbHXJkuZV6RRGQAbMEiIiLjwRYsMlZswTJR3t7y13FxhomDTImyW7CIyLBUcfLPDuHNf4qpbJTegkWmixUsExUQIJ+SOi6O3aOorNiCRUSlp/7JQvY6x5vdo6hs2IJFxopVfiIqJs4iSERExoOzCJKx4hVJRMXEFiwiIjIebMEiY8UugiaKXQJJ95Q9BovjEokMS9sq19AhkIlR+hgsjks0XaxgmagdOwwdAZkeZbdgcVwikWFpA3IMHQKZGKW3YHFcoulSZpWfiAyAY7CIiMh4cAwWGStekURUTMpuwSIiItOi9BYsMl3sIkhExaQCIP59Vh52CSQiMi0qcxVEtoDKXJl5ieMSTRcrWERUTCLfs7JwXCIRkWkR2UL2rDQcl2i62EWQiIiIiJRHle+ZyEiwBctEBQfLX4eHGyYOMiXKnqadiAxLHS6fMS03ONtAkZCpUPo07WS6WMEyUR4e+b/OUWbzORkTTnJBRKWnSuE/waRbnOSCjBU/7YiomDhNOxERGQ9O007GilckERUTW7CIiMh4sAWLjBW7CJqob79ll0DSNWWPweK4RCLDyhnz0NAhkIlR+hgsjks0XaxgmaiUFENHQKZH2S1YHJdIZGAe/Jsj3VJ6CxbHJZou/maJqJg4BouIiIwHx2CRsTLpKzIsLAwqlUr2cHNzk9YLIRAWFgYPDw9oNBp06dIFJ0+elG3jxo0bCAwMhKOjIxwdHREYGIibN2/q+UiIjIGyW7CIDI05iUi3lN6CRabL5LsIenl5Yf/+/dJrtVot/Tx37lzMnz8fK1euhJeXFz788EP06NEDp0+fhr29PQBgyJAhuHz5Mnbt2gUAGD16NAIDA7Fjxw69HgeR4Sl7DBbHJZIxYE4i0h2lj8HiuETTZfIVLHNzc9k3hHmEEPjiiy8wbdo0DBw4EACwatUquLi4YN26dQgJCUFiYiJ27dqFqKgo+Pr6AgC+/fZbdOzYEadPn4aXl5dej4XIsJTdgsVxiWQMmJOIdEfxLVgcl2iylFnlL4ELFy7Aw8MDtWvXxuDBg3HhwgUAwMWLF5Gamgo/Pz+prEajQadOnRAdHQ0AiImJgZ2dHdq1ayeVad++PWxtbaUyxsrdXf4gKjuOwSIqq4qakwAAySr5g6iMOAaLjJVJt2C1adMGK1euRIMGDZCeno5Zs2ahXbt2OHnyJFJTUwEArq6usve4urrin3/+AQCkpqbC2dkZKtV/iUClUsHFxUV6v7EKCZEnr7AwfktCZaXsFiwiQ6vIOQkAzJdayl7nhGYZKBIyFYpvwSKTZdIVrF69eslet23bFnXq1MGqVavQtm3bctlneHg4wqUb7Nwrl30QGYayx2ARGZohchIgz0t3rzIvkelQ+hgsMl0V6oq0s7ND48aNcfbsWakPfFpamqxMWlqatM7NzQ1Xr16FEP+1/gghkJ6eXmgfegAIDg5GbGwsYmNjAdiUz4EQGQRbsIh0SR85CZDnJVtn5iUyHWzBImNVoSpYDx48wKlTp+Du7o7atWvDzc0NERERsvWRkZFS/3ZfX1/cuXMHMTExUpmYmBjcvXtX1gfeGCUnC9mDqOyUPQaL4xLJ2FSknAQAwl0rexCVleLHYHFcosky6S6C77zzDgICAlCzZk2kp6fjo48+wt27dxEUFASVSoUJEyZg9uzZaNCgAerXr49Zs2bBzs4OQ4YMAQA0bNgQPXv2REhIiNS9IiQkBH379jX62ZqkXopEOqPsFiyOSyRDq8g5CQByg7MNHQKZGKW3YHFcouky6QrWlStX8MorryAjIwPOzs5o27YtDh06BE9PTwDAlClTcP/+fYwdOxY3btxAmzZtsHv3bul+IwCwbt06vPXWW/D39wcA9OvXD19//bVBjofIsDgGi6gsmJOIdItjsMhYqcTjnblJp1QqDwDBhg6DSEdmPvZzqMGiKK2wMONpwfL2/unfcZpE+uXh445RsSMNHQaRTswymwMIACpghna6ocMpMfOZVrLXhmrBWuazAsmxvFmkLpl0CxYR6ZIKUiZTII5FJCIyLSpzFUS2gMpcmXmJYxFNl94qWFqtFvv378eBAweQlJSE+/fvw9nZGa1atYKfnx9q1Kihr1CIqFREvmdl4bhEehxzEpHyiWwhe1Yajks0XeXeafX+/fv4+OOPUaNGDfTp0we7d+/GnTt3YGlpiYsXL2LmzJmoXbs2evfujUOHDpV3OEREVIExJxGZEFW+ZyIjUe4tWPXq1YOvry/Cw8Ph5+cHCwuLAmUuXbqEdevWYdCgQZgxYwbGjBlT3mGZvIAA+esdOwwTB5kSTnJBysecZDhmO+T/cmgDcgwUCZkKTnJBxqrcK1i7du1CkyZNnljG09MT06dPx6RJk3Dp0qXyDqlC8PaWf52zY4cym8/JmCh7mnYigDnJkMyOqGWvWcGislL6NO1kusq9yv+0RPY4S0tL1KtXrxyjIaLSU/aNhokA5iQiU6L4Gw2TydL7LIIPHz7EiRMnkJ6eDq1W/o1D79699R0OERUbW7DI9DAnESkXW7DIWOm1ghUREYHAwECkp6cXWKdSqZCbm6vPcEwauwSS7il7DBbHJVJ+zEn6lduXM6aRbil9DBbHJZouvVawxo4di759++L999+Hq6srVCpO+1Je4uIMHQGZHmW3YHFcIuXHnKRfwluZnx1kvJTegsVxiaZLrxWslJQUvPvuu/D09NTnbolIJ5TdgkWUH3MSkbIpvQWLTJder8i+ffsiOjpan7skIp1RdgsWUX7MSUTKpvQWLDJdem3B+uabbzB06FDExcWhSZMmBe4/Mnz4cH2GQ0QlouwWLHYJpPyYk4iUTektWByXaLr0WsH69ddf8dtvv+Hnn3+GjY2NrL+7SqViMiMyaspuweK4RMqPOYlI2ZTegsVxiaZLr1X+d955B2+++SYyMzNx584dZGZmSo/bt2/rMxQiKjHeB4tMC3MSkbLxPlhkrPTagnXz5k289tprsLW11eduK6SwMFW+1+weRWWl7BYsovyYk/TLfKaV7HVOaJaBIiFTofQWLDJdeq3yDxw4EHv27NHnLolIZ9iCRaaFOYlI2diCRcZKry1YderUwXvvvYeDBw+iWbNmBQYUT5w4UZ/hEFGJsAWLTAtzEpGysQWLjJVeK1jLly+Hvb09oqOjC0yNq1KpmMyIjJqyZxEkyo85iUjZlD6LIJkuvVawLl68qM/dVWgcc0W6p+wWLI5LpPyYk/SLY65I15TegsVxiaaLVX4iKiaOwSIiIuPBMVhkrMr9ipw1axbu3r1brLK///47duzYUc4REVHpKLsFiwhgTiIyJUpvwSLTVe4VrPPnz6NmzZoIDg7Gjh07kJKSIq178OABjhw5gq+++gqtW7dGYGAgKleuXN4hEVGpqPI9EykPcxKR6VCZq2TPRMai3CtYK1aswP79+6FSqTB8+HBUr14d5ubm0Gg0sLW1hY+PD7777juMHj0aiYmJ6NChQ7G3ffDgQfTr1w/VqlWDSqXCypUrZeuFEAgLC4OHhwc0Gg26dOmCkydPysrcuHEDgYGBcHR0hKOjIwIDA3Hz5k1Zmfj4eHTu3BkajQbVqlXDhx9+CCE4foMqGpHvWVnCwoTsQRUTcxKR6RDZQvasNDmhWbIHmQ69THLRtGlTfPvtt1iyZAmOHz+OS5cu4f79+3ByckKLFi3g5ORUqu3euXMHTZo0wfDhwzF8+PAC6+fOnYv58+dj5cqV8PLywocffogePXrg9OnTsLe3BwAMGTIEly9fxq5duwAAo0ePRmBgoNQt5Pbt2+jRowc6deqEw4cP49SpUxg5ciRsbW0xadKkUp4RIiIyFOYkIhOhwqPv/NiARUZGJUzkay87Ozt8/fXXGDFiBIBH3xR6eHjgzTffxHvvvQcAuH//PlxcXDBv3jyEhIQgMTERjRo1QlRUFNq3bw8AiIqKQseOHXHq1Cl4eXlhyZIlmDp1KtLS0qDRaAA86sO/ZMkSXLlyBSpV0X/VKpUHgOByPe6ieHvLX8fFGSQMMikf4b9p2t83cCzK5u39E2JjYw0dBpUjY8xJAODh445RsSPL78CfQBUn7zQjvDluhspmttWnj6ZptzTDu1lTDR2OYi3zWYHk2JSnF6RiM9lpVy5evIjU1FT4+flJyzQaDTp16iTd7yQmJgZ2dnZo166dVKZ9+/awtbWVlenYsaOUyADA398fycnJSEpK0s/BlEJAgEr2ICo7TnJBVFoVPScBgPonC9mDqKw4yQUZK5OtYKWmpgIAXF1dZctdXV2ldampqXB2dpZ946dSqeDi4iIrU9g2Ht/H48LDw+Hj4wMfHx8A93R2PESGx2naiUrLUDkJkOelu1eZl8h0cJp2Mla8InUsODgYsbGx/3b/sTF0OEQ6xBYsIiV6PC/ZOjMvkelgCxYZK71McmEIbm5uAIC0tDTUrFlTWp6Wliatc3Nzw9WrVyGEkL4xFEIgPT1dViYtLU227bzXeWWMUVycSQytI6Nihv/GYCkPxyWSIVX0nAQA2la5hg6BTIyZhdmjMVgKbcHiuETTZbArMi0tDVpt+V1ItWvXhpubGyIiIqRlDx48QGRkpNS/3dfXF3fu3EFMTIxUJiYmBnfv3pWViYyMxIMHD6QyERER8PDwQK1atcot/rLasUP+ICo7ZbdgcVwiPQlzUvnTBuTIHkRlpfQWLI5LNF16rWBlZ2djypQpsLe3R7Vq1aQBuVOnTsXixYtLvL07d+7g2LFjOHbsGLRaLS5fvoxjx47h8uXLUKlUmDBhAj799FNs2bIFJ06cwIgRI2BnZ4chQ4YAABo2bIiePXsiJCQEMTExiImJQUhICPr27QsvLy8Aj6bMtbGxwYgRI3DixAls2bIFn3zyCSZOnPjU2ZqITAvHYJFpYU4iUjaOwSJjpdcrcubMmdixYwfWrFkDKysraXnr1q0L3JCxOGJjY9GyZUu0bNkS9+/fR2hoKFq2bIkPPvgAADBlyhS8/fbbGDt2LHx8fJCSkoLdu3dL9xsBgHXr1qF58+bw9/eHv78/mjdvjtWrV0vrHR0dERERgeTkZPj4+GDs2LGYNGkSJk6cWPoTQaRIym7BIsqPOYlI2ZTegkWmS6/3wXrmmWewfPlydO7cGfb29vjrr79Qp04dnD59Gm3atClwt3qlM+R9sIh0T9n3wQoIkL82ZNdZ3gfLOFS0nAQY9j5YRLqm9Ptgme2QT4VgqK6zvA+W7ul1kovk5GR4enoWWJ6Tk4OcHPbHJjJuym7B4lhEyo85iUjZlN6CxbGIpkuvXQQbN26MgwcPFlj+ww8/wDv/FF9EZGQ4BotMC3MSkbJxDBYZK722YIWGhmLYsGH4+++/kZubi40bN+LUqVNYt24ddu7cqc9QTF5wvp6J4eGGiYNMibJbsIjyY07SL3W4fJa03OBsA0VCpkLpLVhkuvRa5Q8ICMAPP/yA3bt3w8zMDDNnzsTZs2exY8cOdO/eXZ+hmDwPD5XsQVR2bMEi08KcpF+qFDPZg6is2IJFxkrvNxrOmxmJiJSGLVhkepiTiJSLLVhkrPRa5a9Tpw6uXbtWYPnNmzdRp04dfYZCRCXGFiwyLcxJRMrGFiwyVnptwUpKSkJubm6B5VlZWfjnn3/0GYrJ+/Zbvc2+TxWGsluwOC6R8mNO0q+cMQ8NHQKZGKW3YHFcounSSwVry5Yt0s87d+6Eo6Oj9Do3Nxe//fYbatWqpY9QKowU3s6AdE4FQPz7rDwFxyLyS4iKijnJQDz4N0e6pTJXQWQLqMyVmZc4FtF06aWC9eKLLwIAVCoVRo0aJVtnYWGBWrVqYf78+foIhYhKTeR7JlIm5iQi0yCyheyZyFjopYKl1T5quq1duzYOHz4MJycnfeyWiIioAOYkIhOh7I4VZML0Ogbr4sWL+twdEemUGR6Nv1JmlwaOS6T8mJOIlM3Mwgzah1rFTnLBcYmmS+/TtN+4cQO//PILLl++jIcP5RfWBx98oO9wiKjYlD3JBcclUmGYk4iUS+mTXHBcounSawXr0KFD6NOnD6ysrHD16lVUq1YNKSkpsLKyQq1atZjMdMjdXf6a/1xS2Sm7BYsoP+YkPUvO14+L/1xSGSm9BYtMl14rWJMnT8bQoUPx5ZdfwsHBAXv37oWtrS1eeeWVAgONqWxCQuSJLCyMiYzKStktWET5MSfpl/lSS9nrnNAsA0VCpkLxLVhksvRa5T9+/DjefPNNqFQqqNVqZGVlwdXVFZ9++inCwsL0GQoRlRhvNEymhTmJSNl4o2EyVnq9Ii0t//v2ytXVFZcuXQIA2NnZITk5WZ+hEFGJsQWLTAtzEpGysQWLjJVeuwi2atUKhw8fRv369dGlSxfMmDEDaWlpWLNmDZo1a6bPUExecjK7BJKuKXsMFsclUn7MSfol3PlPMOmW4sdgcVyiydJrBevjjz9GZmYmAGDWrFkYPnw43nrrLdSvXx8rVqzQZygmLzzc0BGQ6VF2CxbHJVJ+zEn6lRucbegQyMQovQWL4xJNl14rWD4+PtLPzs7O+OWXX/S5eyIqE2W3YBHlx5xEpGyKb8Eik8UrkoiKSdktWEREZFqU3oJFpkuvFazr16/j9ddfR/369VGpUiU4ODjIHkT5VcE9/II1qIJ7hg6FFD6LYHKykD2ImJOoNDTX7iGw5wZorjEvGZrSZxEU7lrZg0yHXrsIjho1CkePHkVwcDA8PDygUqme/iaq0IJwDP44j+H4C1/A19DhVHDKbsHiuETKjzmJSqPlynjU+/UiWq6KR/TENoYOp0JTegsWxyWaLr1W+X/77Td8//33ePfddzFixAgEBQXJHiV18OBB9OvXD9WqVYNKpcLKlStl60eMGAGVSiV7tG3bVlYmKysLb731FpycnGBra4t+/frhypUrsjKXL19GQEAAbG1t4eTkhHHjxuHhw4cljpdKSmAiDkEFYCIOAWCrg2EpuwWLKD/mJCoxIdDu88NQAWj3+WFAMC8ZktJbsMh06fWKdHFxgZ2dnc62d+fOHTRp0gRffvklNBpNoWW6d++OlJQU6fHzzz/L1k+YMAGbN2/G+vXrERkZidu3b6Nv377Izc0FAOTm5qJPnz7IzMxEZGQk1q9fj02bNmHSpEk6O47yEBAgfyhRR1yGIx7NqFMJD9ABlw0cUUWn7BYsovyYk/TLbIe57KFEnpF/w/rWAwCA9c0H8Iy68pR3UHlSegsWmS69T9P+wQcfYNWqVTpJar1790bv3r0BPPpmsDBWVlZwc3MrdN2tW7ewbNkyrFixAj169AAArF69Gp6entizZw/8/f2xe/dunDx5EpcuXUKNGjUAAHPnzsXo0aPx8ccfG20/fW9veVeXHTuU9y3bBByCDR59K2uDh3gbhxAFTwNHVZFxFkEyLcxJ+mV2RC17rQ3IMVAkpdfui8OwuPuoW5fF3Wy0+/xPXOpYw8BRVVycRZCMlV4rWLNmzUJSUhJcXFzg6ekJCwsL2frjx4/rfJ9RUVFwcXFBpUqV0LlzZ3z88cdwcXEBAMTFxSE7Oxt+fn5S+Ro1aqBhw4aIjo6Gv78/YmJi0LBhQymRAYC/vz+ysrIQFxeH5557TucxV0TbsB79cUa2LAtq5KVjNYA+OAuBmbIy21Efz+MV/QRZ4bEFi0wLcxI9yZD+G9Hwx3OyZTmWZjD79/tKMwHU33kOH6nmyMok9quHddtf1FeYFRpbsMhY6bWC9eKL+v3A6dmzJ1544QXUrl0bSUlJmDFjBrp27Yq4uDhYWVkhNTUVarUaTk5Osve5uroiNTUVAJCamgpXV1fZeicnJ6jVaqnM48LDwxEujabnDEPF9S66oTnS4Iq70ODRt5pWyJWVefz1PZgjHbZ4F930GmfFxhYsMi0VIScB8rx09yrzUnHtmd0FbsfSYZd+DxYPHuUl84fyf+Qff51tbY47rjbYM7uzXuOsyNiCRcZKrxWs0NBQfe4OgwcPln5u2rQpvL294enpiZ07d+KFF14ol30GBwcjODgYAKBSeZTLPopDaV0CE+CCxngDy/Ej+uIMbFH0zDp3YIEdqI/R6Id7sCyyHOmasluw8o9F3LHDMHGQ8agIOQmQ5yUPH/dy28/T5PZV1oxp6Y2dsTBhDAa8uhNeP52H5b2i439oY4FTAXWxbVlvZNsyL+mL0luw8o9FVGK3WSpcharye3h4oHr16jh79iwAwM3NDbm5ucjIyJCVS0tLk/rIu7m5IS0tTbY+IyMDubm5RfajNwZxcfKHEtyDJQbjRUyEHx5AXWiZB1BjEvwwBC+ycqV3qnzPyuLtrZI9iAytIuUkABDeWtlDCbJtLfHD9wPwy/xuyLYqPC/lWKnxy/xu2LjheVau9ExlrpI9K43ZEbXsQaaj3CtYDg4OUrKwt7cvcCNHfd7UMSMjA//88w/c3R99g+ft7Q0LCwtERERIZa5cuYLExES0a9cOAODr64vExETZNLkRERGwsrKCt7d3ucdcER2FO7KKaFzNgjmOwHDfwFZsIt8zkfIwJ1FppLRyRe4TKlgp3q6FrqPyJbKF7JnIWJR7F8GFCxfC3t4eAPD111/rdNt37tzBuXOPBqBqtVpcvnwZx44dQ5UqVVClShWEhYVh4MCBcHd3R1JSEqZPnw4XFxcMGDAAAODo6IhRo0ZhypQpcHFxQdWqVTFx4kQ0a9YM3bt3BwD4+fmhcePGGD58OObPn49r165h8uTJGDNmjFHP1qRkPkiG+b/jrbQA7sMCGmTDDIA5cuGDZMSimkFjJCJlYk6i0qgWmwKzvO5oKiBHYw7z+zkwE4BZthYesan451nDDQuosFR49J2fMhuwyISVewXr8Zs1lubGjU8SGxsrmzEpNDQUoaGhCAoKwpIlSxAfH4/vvvsON2/ehLu7O5577jn88MMPUnIFgC+++ALm5uYYNGgQ7t+/j27duuG7776DWv3omyq1Wo2dO3fijTfeQPv27aHRaDB06FB89tlnOj0W+k9HXIItcqSJLCagJ77ALrjiLmyRgw64jG/wrKHDrICUPcmF0sYlUvlgTqLS8Iz8G5b3c6SJLH75ogd6jY+AXfo9WN7PgWfk3zj8eitDh1nhKH2SC6WNS6TiUwnB25CXl0eTXAQbOgzFOY8vURO3sBGNpIksbPAQy/EjBiIBl+CIuhhv6DAroMenyNfv5ACmxtv7J8TGxho6DKqAPHzcMSp2pKHDUJy36yxGpcu3ceLFBtJEFhZ3H2LAqzvRePNp3PR0xOfnXzd0mBXOLLM5UgvWDO10Q4ejWMt8ViA5NsXQYZiUcm/BMjMzg0pVvLbbvDvVU8WWCCd8hE5YiZbSsrwJMEbgKF5EggGjq8iU3YJFBDAnUelcbeiE/TPa4+irzaVleRNgtFz+FxpvPm3A6Coupbdgkekq9wrWDz/8ICWztLQ0fPDBBxgwYAB8fX0BADExMdi2bRtmzpz5pM1QCYWFqfK9Vk5DZV8MLXLdSrSUVbxIn5Q9TTsRwJxkSOYzrWSvc0KzDBRJya3Z+XKR646+2lxW8SL9Ufo07WS6yr2C9fiNHPv164c5c+ZgzJgx0rJXX30VrVu3xrZt2/DGG2+UdzhEVGpswSLlY04iMh1swSJjpdcrcu/evbIBwHmee+457N+/X5+hEFGJsQWLTAtzEpGysQWLjJVeK1hOTk7YtGlTgeWbNm2Cs7OzPkMhohIzy/dMpGzMSUTKltdyxRYsMjbl3kXwcR9++CFGjhyJffv2Sf3dDx06hD179mDZsmX6DMXkKWnMFSmFsluwlDwukcoHc5J+KWnMFSmD0luwlDwukZ5MrxWs4cOHo0GDBvjyyy/x448/AgAaNmyI33//HW3atNFnKERUYhyDRaaFOYlI2TgGi4yV3ipY2dnZGDZsGGbPno21a9fqa7dEpDPKbsEiehxzEpHyKb0Fi0yX3qr8FhYW2L17d7HvP0JExoZjsMh0MCcRKR/HYJGx0msXwRdeeAFbtmzBO++8o8/dEpFOKLsFi2OuKD/mJCJlU3oLFsdcmS69VrBq1qyJWbNmITIyEj4+PrC1tZWtnzhxoj7DIaIS4RgsMi3MSUTKxjFYZKz0WsFauXIlKleujOPHj+P48eOydSqVislMh7y95a/j4gwTB5kSZbdgEeXHnKRfqjj5P8HCm58lVDZKb8Ei06XXCtbFixf1ubsKLSBAPq4gLo7do6is2IJFpoU5Sb/UP1nIXud4s3sUlQ1bsMhYGeyKvHPnDu7evWuo3RNRibEFi0wXcxKR8rAFi4yV3itYixYtQs2aNeHo6AgHBwd4enpi8eLF+g6DiEqMswiS6WFOIlIuziJIxkqvXQRnz56NOXPm4J133kGHDh0AAJGRkZg2bRpu376NadOm6TMck8YugaR7ym7B4rhEyo85Sb+0rXINHQKZGKW3YHFcounSawXrm2++QXh4OF555RVpWbdu3VCvXj28++67TGY6tGOHoSMg06MCIP59Vh6OS6T8mJP0SxuQY+gQyMSozFUQ2QIqc2XmJY5LNF16bVNNT0/Hs88+W2B569atkZaWps9QiKjERL5nImVjTiJSNpEtZM9ExkKvFaz69etj3bp1BZavW7cOXl5e+gyFiIgqOOYkIoVT5XsmMhJ67SIYFhaGl19+GQcPHkT79u0BAL///jsOHDiAjRs36jMUIioxZU/Tzi6BlB9zEpGyKX2ado5LNF16rWC98MIL+OOPP/D555/jp59+AgA0bNgQf/75J1q2bKnPUIioxJQ9yQXHJVJ+zElEyqb0SS44LtF06bWCBQDe3t5Ys2aNvndLRGWm7BYsosIwJxEpl9JbsMh06f2KzMrKwvLly/HOO+9g8uTJWLlyJbKySj5rypw5c/Dss8/CwcEBzs7OCAgIwIkTJ2RlhBAICwuDh4cHNBoNunTpgpMnT8rK3LhxA4GBgXB0dISjoyMCAwNx8+ZNWZn4+Hh07twZGo0G1apVw4cffgghjLu7UXCw/EFUdspuwSIqjK5yEsC89DTqcAvZg6islN6CRaZLrxWshIQE1KtXDxMnTsQff/yBQ4cOYcKECahfvz4SExNLtK39+/fjjTfeQHR0NPbu3Qtzc3N0794d169fl8rMnTsX8+fPx8KFC3H48GG4uLigR48eyMzMlMoMGTIER44cwa5du7Br1y4cOXIEgYGB0vrbt2+jR48ecHV1xeHDh/Hll1/is88+w4IFC8p+QsqRh4dK9iAqO95omEyLLnMSwLz0NKoUM9mDqKx4o2EyViqhx6+8evToARsbG6xevRoODg4AHiWKYcOGISsrC7/++mupt33nzh04Ojpi27ZtCAgIgBACHh4eePPNN/Hee+8BAO7fvw8XFxfMmzcPISEhSExMRKNGjRAVFSUNcI6KikLHjh1x6tQpeHl5YcmSJZg6dSrS0tKg0WgAALNmzcKSJUtw5coVqFRFV15UKg8Ahmk+CgtT5Xtt3N9skhLMfOznUINFYQq8vX9CbGysocOo8MozJwHGmZc8fNwxKnZkmY6rtMxnWsle54Tynj9UNrPM5ki3Z5yhnW7ocBRrmc8KJMemGDoMk6LXKv/vv/+O2bNnS4kMABwcHPDxxx8jKiqqTNvOzMyEVqtF5cqVAQAXL15Eamoq/Pz8pDIajQadOnVCdHQ0ACAmJgZ2dnZo166dVKZ9+/awtbWVlenYsaOUxADA398fycnJSEpKKlPMRMrCFiwyLeWZkwDmJaLyxhYsMlZ6neTC2tq6QD9yALh16xasra3LtO3x48ejRYsW8PX1BQCkpqYCAFxdXWXlXF1d8c8//0hlnJ2dZd/2qVQquLi4SO9PTU1F9erVC2wjb13t2rVl68LDwxEeHv7vq3tlOqay+PZbtliRril7DFb+sYjSnylVWOWZkwDjzEt3rxouL+WMeWiwfZNpUvoYrPxjEXODsw0UCemaXqv8AQEBGDNmDH7//Xfk5uYiNzcXUVFRCAkJQb9+/Uq93YkTJyIqKgqbN2+GWq3WYcQlFxwcjNjY2H+7/9gYLI6UFPmDqOyU3YLFcYmUX3nlJMB485Kts+HyEjyE/EFURkpvweK4RNOl19/ml19+iXr16qFjx46wtraGtbU1OnfujPr16+OLL74o1TbffvttrF+/Hnv37kWdOnWk5W5ubgCAtLQ0Wfm0tDRpnZubG65evSqbeUkIgfT0dFmZwrbx+D6IKgZlt2AR5VceOQlgXiLSF6W3YJHp0msFq1KlSti+fTvOnDmDLVu2YMuWLTh9+jS2bt0KR0fHEm9v/PjxUhJr0KCBbF3t2rXh5uaGiIgIadmDBw8QGRkp9W339fXFnTt3EBMTI5WJiYnB3bt3ZWUiIyPx4MEDqUxERAQ8PDxQq1atEsdMpFzKbsEiyk/XOQlgXiLSJ6W3YJHp0tsYrMuXL8PNzQ2WlpaoW7cu6tatCwB4+PAhLl++jJo1a5Zoe2PHjsXq1auxbds2VK5cWeqbbmdnBzs7O6hUKkyYMAGzZ89GgwYNUL9+fcyaNQt2dnYYMmQIAKBhw4bo2bMnQkJCpP7pISEh6Nu3L7y8vAA8mi535syZGDFiBGbMmIEzZ87gk08+QWho6BNnaiIyPcpuweK4RHqcrnMSwLxEpG9Kb8HiuETTpZcq/4YNG9C7d+9C1wkh0Lt3b2zcuLFE21y8eDEyMzPRrVs3uLu7S4958+ZJZaZMmYK3334bY8eOhY+PD1JSUrB7927Y29tLZdatW4fmzZvD398f/v7+aN68OVavXi2td3R0REREBJKTk+Hj44OxY8di0qRJmDhxYgnPApHSKbsFi+MSKU955CSAeYlI3xTfgsVxiSZLL/fB6tq1K4KCghAUFFTo+tWrV2PFihXYu3dveYeiV4a8D5a7u/w1/6GksuN9sHSF98EyrIqakwDD3gcLyfla1/gPJZUR74OlG7wPlu7ppYtgYmIiOnToUOT6du3aYfLkyfoIpcIICeGNhknXzPCoe6BCvykk+hdzkmGYL7WUveaNhqmszCzMoH2oVW4LFpksvVyRt27dwsOHRfczzcrKwu3bt/URChGVmrLHYBHlYU4iMg1KH4NFpksvFazatWvj8OHDRa7/888/OfMRkdFT5XsmUibmJCLToDJXyZ6JjIVeKlgvvPAC3nvvPaQUMhAoOTkZ77//PgYOHKiPUCqM5GQhexCVncj3rCzu7vIHVVzMSYYh3LWyB1FZiWwhe1acZJX8QSZDL2Owpk6diq1bt6J+/foYNmyYdG+QxMRErF27FjVr1sSUKVP0EUqF8e/svkT0L45LpDzMSYaRG5xt6BDI1KggTXKhRByXaLr0UsGys7PD77//junTp+P777/HjRs3ADy6yWNgYCA+/vhj2RS1RGSMOMkFmQbmJCLTwEkuyFjp7UbDjo6OWLx4MRYtWoSMjAwIIeDs7MybIhIpBie5INPBnESkfJzkgoyV3ipYeVQqFZydnfW9WyIqM2W3YHEsIhWGOYlIuZTegsWxiKZL7xUsIlIqZbdgcVwiEZFpUXoLFsclmi5lVvmJyADM8j0TEREZTl7LlVJbsMh0sQXLRAUEyF/v2GGYOMiUKLsFi4gMy2yH/F8ObUCOgSIhU6H0FiwyXaxgmShvb/lA7R07OP6EykrZY7CIyLDMjqhlr1nBorJS+hgsMl16uSK7deuGLVu2FLk+LS0NarW6yPVEZAzYgkWmgTmJyDSwBYuMlV4qWPv27cPLL7+M0NDQIssIwRYWIuPGMVhkGpiTiEwDx2CRsdJbF8ElS5bgnXfewfHjx7FmzRrY2trK1vPeI7rFLoGke8puweK4RHocc5L+5fbljGmkW0pvweK4RNOltyp///79cejQIZw8eRJt27bFhQsX9LXrCikuTv4gKjtlt2B5e6tkD6rYmJP0T3hrZQ+islJ6C5bZEbXsQaZDr1dkw4YNcfjwYdSoUQPPPvss9uzZo8/dE1GZKLsFiyg/5iQiZVN6CxaZLr1X+R0dHbFz506MGTMGvXv3xueff67vEIioVJTdgkVUGOYkIuVSegsWmS69jMHK35ddpVLhk08+QYsWLTB69Gjs3btXH2EQUZkouwWL4xIpD3MSkWlQegsWxyWaLr1UsIqajWnw4MFo0KABnn/+eX2EQURlouz7YHEsIuVhTiIyDUq/DxbHIpouvVSw9u3bhypVqhS6rkWLFoiLi8POnTv1EQoRlZqyW7CI8jAnEZkGpbdgkenSS5W/c+fOMDcvui5XtWpVDB8+vETbnDNnDp599lk4ODjA2dkZAQEBOHHihKzMiBEjoFKpZI+2bdvKymRlZeGtt96Ck5MTbG1t0a9fP1y5ckVW5vLlywgICICtrS2cnJwwbtw4PHz4sETx6ltYmEr2ICo7jsEi08CcZBjmM61kD6Ky4hgsMlaKvSL379+PN954A9HR0di7dy/Mzc3RvXt3XL9+XVaue/fuSElJkR4///yzbP2ECROwefNmrF+/HpGRkbh9+zb69u2L3NxcAEBubi769OmDzMxMREZGYv369di0aRMmTZqkt2MlMg5swSIqCnMSkf6xBYuMld5uNKxrv/76q+z16tWr4ejoiN9//x0Bj91R1MrKCm5uboVu49atW1i2bBlWrFiBHj16SNvx9PTEnj174O/vj927d+PkyZO4dOkSatSoAQCYO3cuRo8ejY8//hgODg7ldIRExkYFQPz7TESPY04i0j+VuQoiW0BlzrxExkWxLVj5ZWZmQqvVonLlyrLlUVFRcHFxQf369TFmzBikp6dL6+Li4pCdnQ0/Pz9pWY0aNdCwYUNER0cDAGJiYtCwYUMpkQGAv78/srKyEMdR81ShiHzPRFQU5iSi8ieyheyZyFgotgUrv/Hjx6NFixbw9fWVlvXs2RMvvPACateujaSkJMyYMQNdu3ZFXFwcrKyskJqaCrVaDScnJ9m2XF1dkZqaCgBITU2Fq6urbL2TkxPUarVUxhiFhfHDhuhx+cci8m+EyhNzUkE5oVmGDoFMjcI7VuQfi8i/EdNhEhWsiRMnIioqClFRUVCr1dLywYMHSz83bdoU3t7e8PT0xM6dO/HCCy+USyzh4eEIDw//99W9ctkHkWEoe5p2In0xppwEyPPS3avMS2Q6lD5NO5kuxV+Rb7/9NtavX4+9e/eiTp06Tyzr4eGB6tWr4+zZswAANzc35ObmIiMjQ1YuLS1N6iPv5uaGtLQ02fqMjAzk5uYW2o8+ODgYsbGxiI2NBWBThiMjMjac5ILoaYwtJwHyvGTrzLxEpoOTXJCxUnQFa/z48VIia9CgwVPLZ2Rk4J9//oG7uzsAwNvbGxYWFoiIiJDKXLlyBYmJiWjXrh0AwNfXF4mJibJpciMiImBlZQVvb28dHxGRMeM07URPwpxEpF+cpp2MlWK7CI4dOxarV6/Gtm3bULlyZanvuZ2dHezs7HDnzh2EhYVh4MCBcHd3R1JSEqZPnw4XFxcMGDAAAODo6IhRo0ZhypQpcHFxQdWqVTFx4kQ0a9YM3bt3BwD4+fmhcePGGD58OObPn49r165h8uTJGDNmDGdrogpG2S1YHHNF5Yk5iUj/lN6CxTFXpkuxVf7FixcjMzMT3bp1g7u7u/SYN28eAECtViM+Ph79+/dH/fr1ERQUBC8vL8TExMDe3l7azhdffIEBAwZg0KBBaN++Pezs7LBjxw6p37xarcbOnTthY2OD9u3bY9CgQRg4cKC0H6KKgy1YREVhTiLSP7ZgkbFSCSH4tW45Uak8AAQbZN/5e4pw9l4qu5mP/RxqsChMgbf3T/+O0yTSLw8fd4yKHWmQfavi5P8EC29ltjqQ8ZhlNkeaRXCGdrqhw1GsZT4rkBybYugwTIpiuwjSkwUEyOcsjYtjPZrKirMIElHpqX+ykL3O8Wb3KCobziJIxopXJBEVk7LHYBERkWlR+hgsMl2sYBFRMXEMFhERGQ+OwSJjxS6CJopdAkn3lN2CxXGJRIalbZVr6BDIxCi9BYvjEk0XK1gmascOQ0dApkfZY7A4LpHIsLQBOYYOgUyM0sdgcVyi6VLmFUlEBqDsFiwiIjItSm/BItPFChYRFRPHYBERkfHgGCwyVuwiSETFpOwWLHYJJCIyLUpvweK4RNPFChYRFZOyx2BxXCIRkWlR+hgsjks0Xcq8IonIAJTdgkVERKZF6S1YZLrYgmWigoPlr8PDDRMHmRJlt2ARkWGpw+UzpuUGZxsoEjIVSm/BItPFCpaJ8vBQ5VvC8SdUVmzBIqLSU6Xwn2DSLbZgkbHipx0RFZMq3zMREZHhqMxVsmciY8EKFhEVk8j3TEREZDgiW8ieiYwFuwiaqG+/5YcN0eM4LpHIsHLGPDR0CGRqVHj0nZ9CG7A4LtF0sYJlolJSDB0BmR5lT3LBcYlEBubBvznSLaVPcsFxiaaLv1kiKiZOckFERMaDk1yQsWIFi4iKySzfMxERkeHktVwptQWLTBe7CBJRMSm7BYvjEomITIvSW7A4LtF0sYJFRMWk7DFYHJdIRGRalD4Gi+MSTRcrWCbK3V3+mv9cUtkpuwWLiAwsOd9EM/znkspI6S1YZLpYwTJRISHyRBYWxkRGZaXsFiwiMizzpZay1zmhWQaKhEyF4luwyGQp9opctGgRmjVrBgcHBzg4OMDX1xc7d+6U1gshEBYWBg8PD2g0GnTp0gUnT56UbePGjRsIDAyEo6MjHB0dERgYiJs3b8rKxMfHo3PnztBoNKhWrRo+/PBDCMHKClVEbMEiehLmJSL9YgsWGSvFVrCqV6+OTz/9FEeOHEFsbCy6du2K559/HsePHwcAzJ07F/Pnz8fChQtx+PBhuLi4oEePHsjMzJS2MWTIEBw5cgS7du3Crl27cOTIEQQGBkrrb9++jR49esDV1RWHDx/Gl19+ic8++wwLFizQ+/ESGR5nESR6EuYlIv3iLIJkrBTbRbB///6y1x9//DGWLFmCmJgYNG3aFF988QWmTZuGgQMHAgBWrVoFFxcXrFu3DiEhIUhMTMSuXbsQFRUFX19fAMC3336Ljh074vTp0/Dy8sLatWtx7949rFq1ChqNBk2aNMGpU6ewYMECTJw4ESqV8d46PDmZ32aSrim7BYvjEqm8MS89mXBX5mcHGS/Ft2BxXKLJUmwF63G5ubnYuHEj7ty5g3bt2uHixYtITU2Fn5+fVEaj0aBTp06Ijo5GSEgIYmJiYGdnh3bt2kll2rdvD1tbW0RHR8PLywsxMTHo2LEjNBqNVMbf3x/vv/8+kpKSULt2bb0eZ0mEhxe+PCAA8Pb+7w96xw6BuLjCy4aFFW8cl7c3EBDwX9m4OIEdOwrfZnAw4OHxX9lvvxWF/qPr7i4fR5acLHhMBj+mgmOwlHRMHJdI+sS8VFBucHahy812mMPsiPq/cn2zIbwL/4fZfKaV7HVR47hUcWZQ/2Qhvda2yoU2IKfQsupwC6hS/vtcyxnzsPB/dJNVsnFkwl3LYzLwMRU6BktBx8RxiaZL0W2q8fHxsLOzg5WVFV577TVs3boVTZs2RWpqKgDA1dVVVt7V1VVal5qaCmdnZ9m3fSqVCi4uLrIyhW0jbx1RxaLsFiwifWBeItIfbY5W9kxkLBTdguXl5YVjx47h1q1b2LRpE4KCgrB//36DxhQeHo5w6Sv8ewaNhUi3NADu//usPMnJQmqVu3GDrVdUPow9L929yrxEpkNT2Rr3rz2AprK1oUMpFeGulVrlRCXmJVOi6BYsS0tL1K1bF97e3pgzZw5atGiBzz//HG5ubgCAtLQ0Wfm0tDRpnZubG65evSqbeUkIgfT0dFmZwraRt64wwcHBiI2NRWxsLAAbnRwnkXFwzPesTDduFN2NkKisjD0v2TozL5HpcKjhIHtWKlFJQNu38G6MpEwqYUJzu3bt2hUeHh5YvXo1PDw88NZbb+Hdd98FADx48AAuLi747LPPpMHEjRo1wu+//y71d4+Ojkb79u1x6tQpeHl5YcmSJZg6dSrS09Nhbf3o25HZs2dj0aJFuHLlylMHE6tUHgCCy/WYiUh5vL1/+vdLGDJ1xpaXPHzcMSp2ZPkeNBEpyjKfFUiO5cxPuqTYFqxp06YhMjISSUlJiI+Px/Tp07F//34MHToUKpUKEyZMwKeffootW7bgxIkTGDFiBOzs7DBkyBAAQMOGDdGzZ09pYHFMTAxCQkLQt29feHl5AXg0Xa6NjQ1GjBiBEydOYMuWLfjkk0+MfqYmIiLSP+YlIiICFDwGKzU1FcOGDUNqaiocHR3RrFkz/PLLL/D39wcATJkyBffv38fYsWNx48YNtGnTBrt374a9vb20jXXr1uGtt96S3tOvXz98/fXX0npHR0dERERg7Nix8PHxQeXKlTFp0iRMnDhRvwdLRERGj3mJiIgAE+siaGzYRZCICsMugmQo7CJIRPmxi6DuKbaLIBERERERkbFhBYuIiIiIiEhHWMEiIiIiIiLSEVawiIiIiIiIdIQVLCIiIiIiIh1hBYuIiIiIiEhHWMEiIiIiIiLSEVawiIiIiIiIdIQVLCIiIiIiIh1hBYuIiIiIiEhHWMEiIiIiIiLSEVawiIiIiIiIdIQVLCIiIiIiIh1hBYuIiIiIiEhHWMEiIiIiIiLSEVawiIiIiIiIdIQVLCIiIiIiIh1hBYuIiIiIiEhHWMEiIiIiIiLSEVawiIiIiIiIdIQVLCIiIiIiIh1hBYuIiIiIiEhHFFvBWrRoEZo1awYHBwc4ODjA19cXO3fulNaPGDECKpVK9mjbtq1sG1lZWXjrrbfg5OQEW1tb9OvXD1euXJGVuXz5MgICAmBrawsnJyeMGzcODx8+1MsxEhGRcjAvERERoOAKVvXq1fHpp5/iyJEjiI2NRdeuXfH888/j+PHjUpnu3bsjJSVFevz888+ybUyYMAGbN2/G+vXrERkZidu3b6Nv377Izc0FAOTm5qJPnz7IzMxEZGQk1q9fj02bNmHSpEl6PVYiIjJ+zEtERAQAKiGEMHQQulKlShXMmTMHISEhGDFiBDIyMvDTTz8VWvbWrVtwdnbGihUrMHToUADA33//DU9PT/zyyy/w9/fHL7/8gj59+uDSpUuoUaMGAGDNmjUYPXo00tPT4eDg8MR4VCoPAME6PUYiUj5v758QGxtr6DBID4wtL3n4uGNU7EjdHiQRKdoynxVIjk0xdBgmRbEtWI/Lzc3Fhg0bcOfOHbRr105aHhUVBRcXF9SvXx9jxoxBenq6tC4uLg7Z2dnw8/OTltWoUQMNGzZEdHQ0ACAmJgYNGzaUkhgA+Pv7IysrC3FxcXo4MiIiUiLmJSKiisvc0AGURXx8PHx9ffHgwQPY2dlh69ataNq0KQCgZ8+eeOGFF1C7dm0kJSVhxowZ6Nq1K+Li4mBlZYXU1FSo1Wo4OTnJtunq6orU1FQAQGpqKlxdXWXrnZycoFarpTL5hYeHIzw8HABga3sbDRoU/k2lvly9ehXOzs4GjcEY8DzwHOQxhvOQlJRk0P1T+TH2vHT7VCZ+8dmt68MuEWP4GzQ0noNHeB6M4xw8TMo26P5NkaIrWF5eXjh27Bhu3bqFTZs2ISgoCPv370eTJk0wePBgqVzTpk3h7e0NT09P7Ny5Ey+88EK5xRQcHIzgYOPpFujj48OuSOB5AHgO8vA8UHliXno6/g3yHOTheeA5MFWK7iJoaWmJunXrwtvbG3PmzEGLFi3w+eefF1rWw8MD1atXx9mzZwEAbm5uyM3NRUZGhqxcWloa3NzcpDJpaWmy9RkZGcjNzZXKEBER5WFeIiIiRVew8tNqtcjKyip0XUZGBv755x+4u7sDALy9vWFhYYGIiAipzJUrV5CYmCj1l/f19UViYqJsityIiAhYWVnB29u7HI+EiIhMAfMSEVHFo9gugtOmTUOfPn1Qo0YNZGZmYt26ddi/fz927tyJO3fuICwsDAMHDoS7uzuSkpIwffp0uLi4YMCAAQAAR0dHjBo1ClOmTIGLiwuqVq2KiRMnolmzZujevTsAwM/PD40bN8bw4cMxf/58XLt2DZMnT8aYMWOeOlOTsTCmbiGGxPPAc5CH54HKC/NS8fBvkOcgD88Dz4HJEgoVFBQkatasKSwtLYWzs7Po1q2b2LVrlxBCiHv37gk/Pz/h7OwsLCwsRM2aNUVQUJC4fPmybBsPHjwQb775pqhSpYrQaDSib9++BcpcunRJ9OnTR2g0GlGlShXx1ltviQcPHujtOImISBmYl4iISAghTOo+WERERERERIZkUmOwiIiIiIiIDIkVLCIiIiIiIh1hBcuELV68GLVr14a1tTW8vb0RGRlp6JBK7eDBg+jXrx+qVasGlUqFlStXytYLIRAWFgYPDw9oNBp06dIFJ0+elJW5ceMGAgMD4ejoCEdHRwQGBuLmzZuyMvHx8ejcuTM0Gg2qVauGDz/8EMbSi3bOnDl49tln4eDgAGdnZwQEBODEiROyMqZ+HhYtWoRmzZrBwcEBDg4O8PX1xc6dO6X1pn78RErGnGRan0XMSY8wL1GhDDHwi8rfhg0bhLm5uQgPDxcJCQnizTffFLa2tuLSpUuGDq1Udu7cKaZPny42btwoNBqNWLFihWz9J598Iuzs7MSmTZtEfHy8eOmll4S7u7u4ffu2VKZnz56iUaNGIjo6WkRHR4tGjRqJvn37Sutv3bolXF1dxUsvvSTi4+PFxo0bhZ2dnZg3b56+DvOJ/Pz8xPLly0V8fLw4fvy4eP7554Wrq6u4du2aVMbUz8O2bdvEzz//LM6ePStOnz4t3n33XWFubi7++usvIYTpHz+RUjEnmd5nEXPSI8xLVBhWsExU69atxejRo2XL6tatK6ZNm2agiHTH1tZWlsy0Wq1wc3MTs2bNkpbdu3dP2NnZiW+++UYIIURCQoIAIKKioqQykZGRAoA4deqUEEKIxYsXC3t7e3Hv3j2pzEcffSQ8PDyEVqst56MquczMTGFmZiZ+/PFHIUTFPQ+VK1cW33zzTYU9fiIlYE4y/c8i5qT/MC8RuwiaoIcPHyIuLg5+fn6y5X5+foiOjjZQVOXn4sWLSE1NlR2vRqNBp06dpOONiYmBnZ2ddLNOAGjfvj1sbW1lZTp27AiNRiOV8ff3R3JyMpKSkvRzMCWQmZkJrVaLypUrA6h45yE3NxcbNmzAnTt30K5duwp3/ERKwZxUMT6LKnpOApiX6D+sYJmgjIwM5ObmwtXVVbbc1dUVqampBoqq/OQd05OONzU1Fc7OzlCpVNJ6lUoFFxcXWZnCtvH4PozJ+PHj0aJFC/j6+gKoOOchPj4ednZ2sLKywmuvvYatW7eiadOmFeb4iZSGOQnSa1P+LKqoOQlgXqKCzA0dABGV3MSJExEVFYWoqCio1WpDh6NXXl5eOHbsGG7duoVNmzYhKCgI+/fvN3RYREQVVkXOSQDzEhXEFiwT5OTkBLVajbS0NNnytLQ0uLm5GSiq8pN3TE86Xjc3N1y9elU2444QAunp6bIyhW3j8X0Yg7fffhvr16/H3r17UadOHWl5RTkPlpaWqFu3Lry9vTFnzhy0aNECn3/+eYU5fiKlYU6C9NoUP4sqek4CmJeoIFawTJClpSW8vb0REREhWx4RESHr42sqateuDTc3N9nxPnjwAJGRkdLx+vr64s6dO4iJiZHKxMTE4O7du7IykZGRePDggVQmIiICHh4eqFWrln4O5inGjx8vJbIGDRrI1lWk8/A4rVaLrKysCnv8RMaOOcl0P4uYkwrHvEScRdBEbdiwQVhYWIilS5eKhIQEMW7cOGFrayuSkpIMHVqpZGZmiqNHj4qjR48KjUYjZs6cKY4ePSpN8fvJJ58IBwcHsXnzZhEfHy8GDRpU6DSoTZo0kaZBbdKkiWwa1Js3bwpXV1cxaNAgER8fLzZv3izs7e2NZhrUN954Q9jb24vffvtNpKSkSI/MzEypjKmfh6lTp4qDBw+KixcviuPHj4tp06YJlUolfv75ZyGE6R8/kVIxJ5neZxFz0iPMS1QYVrBM2KJFi4Snp6ewtLQUrVq1EgcOHDB0SKW2b98+AaDAIygoSAjxaDrY0NBQ4ebmJqysrESnTp1EfHy8bBvXr18XQ4cOFfb29sLe3l4MHTpU3LhxQ1bm+PHjomPHjsLKykq4ubmJsLAwo5kCtbDjByBCQ0OlMqZ+HoKCgkTNmjWFpaWlcHZ2Ft26dRO7du2S1pv68RMpGXOSaX0WMSc9wrxEhVEJwdtAExERERER6QLHYBEREREREekIK1hEREREREQ6wgoWERERERGRjrCCRUREREREpCOsYBEREREREekIK1hEREREREQ6wgoWkY5s2rQJKpVKer1y5UrY2dkZMKLiGTlyJD788MMybSM+Ph7VqlXD3bt3dRQVERGVBXMScxIZDitYZDBarRadOnVCQECAbPm9e/fg5eWF11577Ynvz8zMxPvvv49GjRpBo9HA1dUVXbp0wfr166HVassz9GIZNGgQLly4oPPtqlQqbNq0SSfbio+Px7Zt2zBhwoQybadp06Zo27YtFixYoJO4iIj0jTmpdJiTiApiBYsMxszMDCtXrsS+ffuwfPlyafnUqVORm5uL+fPnF/nemzdvwtfXF8uXL8fkyZMRGxuLqKgoBAUF4aOPPsLly5fLLe6HDx8Wq5xGo4GLi0u5xaELCxcuxMCBA+Hg4FDmbY0cORJLlixBTk6ODiIjItIv5iTDY04ikyGIDGzJkiXCwcFBJCUliT179gi1Wi0iIyOf+J7XX39d2NjYiL///rvAuvv374v79+8LIYS4fv26GD58uKhUqZKwtrYW3bp1EydOnJCV37x5s2jSpImwtLQU1atXF7NmzRJarVZa7+npKUJDQ8XIkSOFo6OjePHFF4UQQqxatUrUrFlTaDQa0adPH/H111+Lx/+kVqxYIWxtbaXXoaGhonHjxmL9+vWiTp06ws7OTvTv319cvXpVKvPnn3+KHj16iKpVqwp7e3vRvn17ER0dLYsFgPTw9PSU1v3444+iVatWwsrKStSqVUu8++67Iisrq8hzmJOTIxwdHcW2bdtkyz09PcXMmTNFUFCQsLOzE9WrVxcbNmwQN27cEIMGDRK2traibt264tdff5W9LysrS1hZWYmIiIgi90lEZOyYk5iTiMqKFSwyCn5+fqJjx46iRo0aYvLkyU8sm5ubKypXrizGjBnz1O3269dPeHl5iQMHDojjx4+LgIAAUb16dXHv3j0hhBCxsbHCzMxMfPDBB+L06dNizZo1wtbWVnz11VfSNjw9PYW9vb349NNPxdmzZ8WZM2fEoUOHhEqlErNmzRKnT58W33zzjahSpcpTk5mtra14/vnnxV9//SWio6NFzZo1RXBwsFTmt99+E999951ISEgQiYmJYuzYsaJSpUoiIyNDCCFEenq6ACCWLl0qUlJSRHp6uhBCiF27dgl7e3uxfPlyce7cObF3715Rv359MWnSpCLPzZEjRwQAceXKFdlyT09PUblyZbFo0SJx5swZMXHiRGFlZSV69eolVq1aJc6ePSteffVV4ezsLP3TkKdNmzZixowZT/29EBEZM+akR5iTiEqHFSwyChcuXBAqlUrUrVtXPHjw4Ill09LSBACxYMGCJ5Y7c+aMACAOHDggLbt586ZwcHAQS5cuFUIIMWTIEPHcc8/J3hcaGiqqVasmvfb09BR9+/aVlXnllVdE9+7dZctGjRr11GRmZWUlbt68KS2bNWuWeOaZZ4o8Bq1WK9zc3MTq1aulZQDExo0bZeU6duwoPvzwQ9myrVu3CltbW9k3n/nXq1QqkZubK1vu6ekpBg8eLL3OzMwUAMRbb70lLbt48aIAIA4fPix774ABA8SwYcOKPB4iIiVgTioccxJR8XAMFhmF5cuXQ6PR4MqVK7h48eITywohirXNxMREmJmZwdfXV1rm6OiIpk2bIiEhQSrTvn172fs6dOiAf/75B7dv35aW+fj4FNj249sFUOB1YTw9PeHo6Ci99vDwQHp6uvQ6PT0dISEhqF+/PhwdHWFvb4/09PSn9t+Pi4vDxx9/DDs7O+kxZMgQ3L17F6mpqYW+5/79+7CwsICZWcGPgWbNmkk/29nZwcbGBk2bNpWWubq6SvE+TqPR4P79+0+MlYjI2DEnPcKcRFQ65oYOgOjw4cP45JNP8OOPP2LJkiUICgpCdHQ01Gp1oeWdnZ1RqVIlJCYmlnqfj09dW5wytra2pd7X4ywsLArs4/HZpYKCgpCWlobPP/8ctWrVgpWVFbp16/bUQcxarRahoaF46aWXCqxzdnYu9D1OTk54+PAh7t27Bxsbm6fG+fiyvHOTf2as69evo1atWk+MlYjImDEnMScRlRVbsMigHjx4gOHDh2PEiBHo1asXwsPDce7cOcydO7fI95iZmWHw4MFYu3Ytrly5Uug2Hzx4gIYNG0Kr1SImJkZad/v2bcTHx6NRo0YAgIYNG+L333+XvT8qKgrVq1eHvb19kTE0bNgQhw4dki3L/7o0oqKi8NZbb6FPnz5o3Lgx7O3tkZKSIitjYWGB3Nxc2bJWrVrh1KlTqFu3boGHuXnh36O0aNECAKRvTnXhxIkTaNWqlc62R0SkT8xJcsxJRKXDChYZ1PTp0/HgwQPpXhVubm5YtGgRwsLCcPLkySLf9/HHH6NmzZpo06YNVqxYgZMnT+LcuXNYvXo1vL29kZqainr16qF///4ICQlBZGQk4uPjMWzYMDg4OGDIkCEAgEmTJuHAgQMICwvDmTNnsHbtWsyfPx9Tpkx5Ytzjxo3Dnj17MGfOHJw9exZLly7F1q1by3w+6tevjzVr1iAhIQGHDx/G4MGDYWlpKStTq1Yt/Pbbb0hNTcWNGzcAAB988AHWrVuHDz74ACdOnMCpU6ewadOmJx6Hs7MzWrVqhaioqDLHDQBJSUn4559/4Ofnp5PtERHpG3OSHHMSUSkZehAYVVwHDhwQarVa7Nu3r8C6F198UXh7e4vs7Owi33/z5k3x7rvvCi8vL2FlZSWcnZ1F586dxfr166VBsiWZEtfCwqLIKXE/++yzAvtfvny5qFGjhrC2thY9e/YUCxcuLNaUuI/LX+bYsWOidevWwtraWtSpU0d89913onHjxiI0NFQq8+OPP4q6desKc3Nz2ZS4v/76q+jQoYPQaDTC3t5eeHt7i4ULFxZ5/oQQ4ptvvhE+Pj6yZYUdr62trVixYoX0+v79+wKA2LFjh7Rs9uzZwt/f/4n7IyIyVsxJzElEuqISopijM4nI5Dx48AANGjTA6tWr0bFjx1JvJysrC/Xq1cP69esLDNAmIiIqDuYkMhXsIkhUgVlbW+O7777D9evXy7SdS5cu4b333mMiIyKiUmNOIlPBFiwiIiIiIiIdYQsWERERERGRjrCCRUREREREpCOsYBEREREREekIK1hEREREREQ6wgoWERERERGRjrCCRUREREREpCOsYBEREREREekIK1hEREREREQ6wgoWERERERGRjrCCRUREREREpCOsYBEREREREekIK1hEREREREQ6wgoWERERERGRjrCCRUREREREpCOsYBEREREREekIK1hEREREREQ6wgoWERERERGRjpgbYqcqVV0A9wGo8paU4meUsHxp3/ukAynDLlDE8uKGV5rtG3Lfhb5XPLb80c8qlZDKqVRC/hoCqn/L5b1+/Pm/n/OWF7au4LLyKlPYewwRR+H7L/oclXgf4t8y/21W+tUC/y4Tj/1c2HLkK1Pc5U/b7pO2VdLt6Gq/OjoHAoDI+1k89nMxQirOz4XttqTbaubvj127dsHY/ZeTgNLllOJ8+OnqvU86kDLsAkUsL0l4T9u+Ifdd7PeKx5aLwnPSv+VL8hkMiMd2bbgc9KS8ZOg8lXeOdLIPIf7LQ6KQnJT3rMv8UJLPf13nwrLsV0fnQPrx8XwkihdScX4uKuyS/qyvvGSQChZwD8BYABaPhVGSny3wX+hl/flp+8pHBUBdyKbMFbC8OO8p6medlxeAee6/P+fCTJ3z6EeLXKj/Xa42z4XaPAfmea/NcqHGvz8jF2r8+x4Ud3nezzkw18nyJ8dQmvcoLu7cXKhz/l2eq8W/v0aocwDVv79e5Pz7yPu5sOW5+coUtTz/e8uyD33EV9R7dRBfdg7w76lHdu5jP+f8Vzz730feW0v7c3Yp378jIwPKkJeTgPLJI7rMcfnk/deuzrcpcyNc/rS8UNS2ivpZp+X//ffLPFfKTWbqnEJzEgCYm+dCbVb+n8clWy7fV/nmFyONO/ffn3NyYJ6r/fdnSLlJVZbP8qd9dpdXvtBVfEWV0VF82f/+nJPzKCdJP+f8Vzx/Tintz2V5v77yErsIEhERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6Yi5IXbaqFFlaDRbDbFr4yL+fX7476MQV69ehbOzs74iUryizlfuvw8CAPW/D15fJVUu50v92LOlbjdtaPfv3zd0CMWi5Jyk82tS4L98VERe0iUlfwaVJfbcfM+GYDznXp3v+emMJ/bSMfr4//s3oUBeMvrYn0JfeckgFSyNRoPY2FhD7FpxfHx8eK5KgOerZHi+Sobnq2R8fHwMHUKxKDknKf2aVHL8So4dUHb8So4dUHb8So4d0F9eYhdBIiIiIiIiHWEFi4iIiIiISEcMUsEKDg42xG4VieeqZHi+Sobnq2R4vkpGKedLKXEWRsmxA8qOX8mxA8qOX8mxA8qOX8mxA/qLXyWEEE8vRkRERERERE/DLoJEREREREQ6wgoWERERERGRjui0gnXw4EH069cP1apVg0qlwsqVK59YPikpCSqVqsBj165dugzLaJX0fAGAEAJffPEFGjRoACsrK7i7u2PatGnlH6wRKOn5CgsLK/T6UqlUSE9P10/QBlSa6+vXX3+Fr68v7O3t4eTkhP79++PMmTPlH6wRKM35+uGHH9CiRQvY2NjA09MTn332WfkHagTmzJmDZ599Fg4ODnB2dkZAQABOnDjx1PfFx8ejc+fO0Gg0qFatGj788EOUtZf64sWLUbt2bVhbW8Pb2xuRkZFPLL9o0SI0bNgQGo0GXl5e+O677wqU2bx5Mxo1agQrKys0atQIW7fK75ElhEBYWBg8PDyg0WjQpUsXnDx50uCxL126FB07dkTlypVRqVIlPPfcc4iKipKVKexz0c3NrcSxl0f8K1euLPTz+sGDB2Xarz5i79KlS6GxN27cuMTH9zSl+awqzt+ePq778ohdn9d9ecSvr+u+PGI35uv+wYMHGDFiBJo1awYLCwt06dKl0HIHDhyAt7c3rK2tUadOHXzzzTcFypTq3Asd2rlzp5g+fbrYuHGj0Gg0YsWKFU8sf/HiRQFA7Nq1S6SkpEiPrKwsXYZltEp6voQQ4u233xb16tUT27ZtE+fPnxdHjhwRO3fuLP9gjUBJz1dmZqbsukpJSRGdO3cWXbp00U/ABlbS83XhwgVhZWUlJk+eLM6ePSuOHj0qevToIZ555hn9BGxgJT1fP//8s1Cr1WLRokXi/Pnz4qeffhLu7u5i4cKF+gnYgPz8/MTy5ctFfHy8OH78uHj++eeFq6uruHbtWpHvuXXrlnB1dRUvvfSSiI+PFxs3bhR2dnZi3rx5pY5jw4YNwtzcXISHh4uEhATx5ptvCltbW3Hp0qVCyy9evFjY2tqKdevWifPnz4v169cLOzs78eOPP0ploqOjhVqtFrNmzRIJCQli1qxZQq1Wi0OHDkllPvnkE2FnZyc2bdok4uPjxUsvvSTc3d3F7du3DRr7kCFDxMKFC8WRI0fEqVOnREhIiLCxsRFnzpyRyoSGhgovLy/Z52J6enqx4y7P+FesWCFsbGwKfG6XZb/6iv3atWuymJOSkoS9vb0ICwsr0fEVR0k/q4rzt6ev6748YtfndV8e8evrui+P2I35ur9z544ICQkR3377rejfv7/o3LlzgTIXLlwQNjY24s033xQJCQkiPDxcmJubi02bNkllSnvudVrBepytrW2xK1iHDx8urzAUozjn69SpU8Lc3FwkJCToJygjVpzzld/ly5eFmZmZWLt2bfkEZcSKc742btwozMzMRE5OjrRs7969AoC4evVqOUdoXIpzvl555RXx/PPPy5Z99dVXonr16kKr1ZZjdMYnMzNTmJmZyf7hzG/x4sXC3t5e3Lt3T1r20UcfCQ8Pj1Kfr9atW4vRo0fLltWtW1dMmzat0PK+vr5iwoQJsmUTJ04U7du3l16//PLLonv37rIy3bp1E4MHDxZCCKHVaoWbm5uYNWuWtP7evXvCzs5OfPPNNwaNPT+tVitcXV3FV199JS0LDQ0VjRs3LnacRSmP+FesWCFsbW11ul99xZ7fmjVrhFqtFpcvX5aWFef4Sqo4n1XF+dvT13VfHrHnV57X/eN0Fb++rvvyiD0/Y7ruHzd27NhCK1hTpkwRdevWlS0bNWqUaNu2rfS6tOfeKMZgvfDCC3BxcUH79u2xadMmQ4djtLZv3446depg165dqFOnDmrVqoWgoKAK0d1NF5YtW4bKlStj4MCBhg7FKD377LOwsLDA//73P+Tm5iIzMxOrVq3Cs88+CycnJ0OHZ3SysrJgbW0tW6bRaHDlyhVcunTJQFEZRmZmJrRaLSpXrlxkmZiYGHTs2BEajUZa5u/vj+TkZCQlJZV4nw8fPkRcXBz8/Pxky/38/BAdHV3oe4r6nf3555/Izs6W4sy/TX9/f2mbFy9eRGpqqqyMRqNBp06dityvvmIvbD8PHjwo8Hu5cOECPDw8ULt2bQwePBgXLlwoVtz6iP/+/fvw9PRE9erV0bdvXxw9erRM+9Vn7I9bunQpevbsiRo1asiWP+n4yktx/vb0cd2XV+z5ldd1XxrFjb+8r/vyjP1xxnTdF0dR131sbCyys7PLdO4NWsGys7PDvHnz8MMPP+Dnn39Gt27dMGjQIKxZs8aQYRmtCxcu4NKlS9iwYQNWrlyJ1atX49SpUwgICIBWqzV0eEYtNzcXy5cvR2BgIKysrAwdjlHy9PREREQEQkNDYWVlBUdHR8THx+Onn34ydGhGyd/fH9u3b8fu3buh1Wpx5swZzJ8/HwCQkpJi4Oj0a/z48WjRogV8fX2LLJOamgpXV1fZsrzXqampJd5nRkYGcnNzC91mUdvz9/fH8uXLcfjwYQghEBsbi//973/Izs5GRkbGE+PM22bec0n2q6/Y85sxYwbs7OzQr18/aVmbNm2wcuVK7Nq1C0uXLkVqairatWuHa9euFSv28ozfy8sLy5cvx/bt27F+/XpYW1ujffv2OHv2bKn3q6/YH3fmzBkcOHAAY8aMkS1/2vGVl+L87enjui+N0nxulNd1XxrFiV8f1315xf44Y7vui6OoY8zJyUFGRkaZzr25zqMtAScnJ0yaNEl67ePjg4yMDMydOxfDhg0zYGTGSavVIisrC6tXr0b9+vUBAKtXr4aXlxcOHz6MNm3aGDhC47Vr1y78/fffBf7w6T+pqakYNWoUhg8fjldeeQWZmZn44IMP8PLLL2Pv3r0wMzOKBm+jMWbMGJw/fx79+/dHdnY2HBwcMH78eISFhVWoczVx4kRERUUhKioKarXa0OE80fvvvy/9YyWEgKurK4KCgjB37lyj/52VNPYvv/wS3377Lfbs2QMHBwdpea9evWTl2rZtizp16mDVqlWYOHGiQeP39fWVVdLbtWuHFi1aYOHChfjqq6/KLbanKem5X7p0Kdzd3dGnTx/ZcmM9PlNibNd9cZjKdcHrXs7oMkqbNm2MolZrjNzd3WFubi5VrgCgXr16UKvVuHz5sgEjM37h4eFo164dGjVqZOhQjNaiRYtga2uLuXPnomXLlujUqRPWrFmDAwcOlGs3BKVSqVT49NNPcefOHVy6dAmpqalo3bo1AKBOnToGjk4/3n77baxfvx579+596jG7ubkhLS1NtizvdWlm83JycoJarS50m0VtT6PRYPny5bh37x6SkpJw+fJl1KpVC/b29nB2dn5inHnbzHsuyX71FXueL774AjNmzMDPP/8sXZNFsbOzQ+PGjUuUd8s7/jxqtRo+Pj5SbKXZr75jf/jwIVatWoWRI0fC3PzJ32HnP77yUpy/PX1c96VRks+N8r7uS6M0n3vlcd2XRkliN8brvjiKOkZzc3M4OTmV6dwbXQXr2LFjcHd3N3QYRql9+/bIycnB+fPnpWUXLlxAbm4uPD09DRiZcUtOTsbOnTvZevUU9+7dK9ACkfeaXVCLplarUa1aNVhaWmL9+vXw9fUt8h9GUzJ+/HipctWgQYOnlvf19UVkZKRsat6IiAh4eHigVq1aJd6/paUlvL29ERERIVseERGBdu3aPfG9FhYWqF69OtRqNTZs2IC+ffvKWlGetM3atWvDzc1NVubBgweIjIx86n7LO3YAWLBgAd5//33s3LkTHTp0eGosDx48wKlTp0qUd8sz/scJIXD8+HEptrLsV1+xb9u2DRkZGRg1atRTY8l/fOWlOH97+rjuyyt2QD/XfWmU5nOvPK778o7dGK/74ijquvfx8YGFhUXZzn2xp+AohszMTHH06FFx9OhRodFoxMyZM8XRo0elqQynTZsmunbtKpVfuXKlWLt2rUhISBCnTp0Sn332mbCwsBALFizQZVhGq6TnKzc3V7Rq1Up06tRJHDlyRBw5ckR06tRJtGnTRuTm5hrqMPSmpOcrz0cffSQcHBzE3bt39R2yQZX0fP32229CpVKJmTNnijNnzoi4uDjh7+8vatSoIe7cuWOow9Cbkp6vq1evisWLF4uEhARx9OhRMW7cOGFtbS3++OMPQx2C3rzxxhvC3t5e/Pbbb7JpdzMzM6Uy+c/XzZs3haurqxg0aJCIj48XmzdvFvb29mWept3CwkIsXbpUJCQkiHHjxglbW1uRlJQkhBAiMDBQBAYGSuVPnz4tvvvuO3HmzBnxxx9/iEGDBokqVaqIixcvSmV+//13oVarxZw5c0RiYqKYPXu2MDc3LzBdtYODg9i8ebOIj48XgwYNKtU07bqOfe7cucLCwkJ8//33st/LzZs3pTKTJk0S+/fvFxcuXBCHDh0Sffr0Efb29tJ+DRl/WFiY2LVrlzh//rw4evSoGDlypDA3N5f9TT1tv4aKPU+3bt0KzMZXkuMrjpJ+VhXnb09f1315xK7P67484tfXdV8esecxxuteCCFOnjwpjh49KgYNGiS8vb2l9+fJm6Z9/PjxIiEhQSxdulRYWFgUmKa9NOdepxWsffv2CQAFHkFBQUIIIYKCgoSnp6dUfuXKlaJhw4bCxsZG2NvbC29vb7F69WpdhmTUSnq+hBAiOTlZvPjii8LOzk44OzuLIUOGiNTUVP0HbwClOV9arVbUqlVLvP766/oP2MBKc77Wr18vWrVqJWxtbYWTk5Po27evOHnypP6DN4CSnq+rV6+Ktm3bCltbW2FjYyO6desm+2fElBV2ngCI0NBQqUxh19fx48dFx44dhZWVlXBzcxNhYWFlntJ+0aJFwtPTU1haWopWrVqJAwcOSOs6d+4sm5o3ISFBtGjRQmg0GuHg4CD69+8vTp06VWCbGzduFF5eXsLCwkI0aNBAbN68WbZeq9WK0NBQ4ebmJqysrESnTp1EfHy8wWP39PR84jUshJD+KbawsBAeHh7ihRdeKPXfuK7jnzBhgqhZs6awtLQUzs7Ows/PT0RHR5dov4aKXQghzp8/L1Qqlfj+++8L3Wdxj+9pSvPZXpy/PX1c9+URuz6v+/KIX1/XfXldN8Z83Rd1bTxu//79omXLlsLS0lLUqlVLLFmypMC+S3PuVULku5U3ERERERERlYrRjcEiIiIiIiJSKlawiIiIiIiIdIQVLCIiIiIiIh1hBYuIiIiIiEhHWMEiIiIiIiLSEVawiIiIiIiIdIQVLFKcsLAwNGnSRGfbW7lyJezs7HS2vfxmzpyJV199tdy2b0zS09Ph7OyMK1euGDoUIiJFqlWrFubNm1cu275x4wZcXV1x/vx5AMD+/fuhUqmQkZFRLvsri2XLlkGlUqFly5a4detWqbYRHx+PatWq4e7duzqOjujJWMEivenXrx+6detW6LrExESoVCrs3r1bz1EBgwYNwoULF6TXuqzApaenY/78+ZgxY4ZOtvfw4UM4Ojri2LFjOtmerrm4uGD48OEIDQ01dChERMUyYsQIqFQqqFQqWFhYwMXFBc899xwWLVqE7Oxsvcdz+PBhvPHGG9JrlUqFTZs26WTbs2fPRu/evfHMM8/oZHuPS0lJgaWlJa5du1bmbW3evBmvv/465s+fDyEE+vbti/v378vKaLVa9OvXDzVr1oS1tTXc3d0xbNgw/PPPP1KZpk2bom3btliwYEGZYyIqCVawSG9GjRqFffv2ISkpqcC6ZcuWwdPTE927d9d7XBqNBi4uLuWy7f/9739o3bo16tSpo5Pt7du3D5UrV0aLFi10sr3yMHLkSKxduxbXr183dChERMXSvXt3pKSkICkpCbt370ZAQABCQ0PRsWNHvbd+ODs7w8bGRufbvXfvHv73v/9h1KhROt82APz4449o164dqlatWqbtREREYMSIEVizZg0mTpyI/fv3Izc3Fy+++GKBCm/Xrl3xww8/4PTp09i8eTMuXLiAAQMGyMqMHDkSS5YsQU5OTpniIioJVrBIb/r06QNXV1esWLFCtjw7OxurV6/Gq6++CjMzMyQkJKBPnz6wt7eHi4sLXnnlFaSmpha5Xa1Wi48++gg1atSAlZUVmjZtiu3bt8vKJCcnY+jQoahatSpsbGzQokUL7Nu3D4C8i+DKlSsxc+ZMnDx5UvpGc+XKlXj11VfRt2/fAvutWbPmE78ZW7duHQICAmTLunTpgtdffx2TJk1ClSpV4OzsjC+//BJZWVkYO3YsKlWqhJo1a2L16tUFtrd9+3b0799fFvcvv/yCBg0awMbGBv369cOtW7ewadMm1KtXD46OjggMDJR983fw4EG0bdsWdnZ2cHR0ROvWrXHixAlpfXR0NDp37gwbGxtUq1YNr7/+Om7fvi2tF0Jg/vz5qFevHqysrFC9enVMnz5dWt+kSRN4eHhgy5YtRZ4XIiJjYmVlBTc3N1SrVg0tWrSQ/rE/cuQI5s6dK5V7+PAhpk6diurVq8PGxgbPPvssfv31V2l9Xpe73377DW3atIGNjQ18fHxw5MgRqcytW7cQGBgIFxcXWFtbo06dOvjiiy+k9Y93EaxVqxYA4KWXXoJKpUKtWrWQlJQEMzMzxMbGyo5h6dKlcHJywsOHDws9xp9//hkqlQrt27cv8jxkZWVhwIABaNWqFdLT0wEAf/zxB1q1agVra2u0bNlS2s7+/ftl7308P+X1BFm1ahVq1aoFW1tbjBw5Eg8fPsTixYtRo0YNVK1aFRMnToRWq5W2cejQIQwZMgQbNmzAyy+/DACoVKkSIiIikJWVhaCgIKm8mZkZJkyYgLZt28LT0xPt2rXDtGnTcPjwYTx48EDapp+fH65fv14gXqJyJYj0aOrUqaJmzZoiNzdXWrZ582ZhZmYmLl++LJKTk0XVqlXFlClTREJCgvjrr79E3759RevWraX3hIaGisaNG0vvX7BggbC3txdr164Vp0+fFu+//74wMzMTR48eFUIIcefOHVG3bl3Rrl07cfDgQXHu3DmxefNmsXfvXiGEECtWrBC2trZCCCHu3bsnJk2aJLy8vERKSopISUkR9+7dE9HR0UKtVovk5GRpv7t27RIWFhYiPT290GO9du2aUKlUIioqSra8c+fOwt7eXoSGhoozZ86IefPmCQCiZ8+e4osvvhBnz54VM2bMEJaWlrL9abVaUa1aNfHbb79JcZubm4tu3bqJ2NhYER0dLdzd3UW3bt1E3759xV9//SX27t0rKlWqJObNmyeEECI7O1tUqlRJTJo0SZw7d04kJiaKtWvXioSEBCGEEMePHxe2trZi3rx54syZM+LQoUOibdu2YuDAgVIc06ZNE46OjmLZsmXi7NmzIjo6WixatEh2jIMGDRLDhg172uVARGRwQUFBok+fPoWuCwgIkOWbIUOGiDZt2ogDBw6I8+fPi4ULFwoLCwtx7NgxIYQQ+/btEwDEs88+K/bu3SsSExOFn5+faNCggdBqtUIIId58803RvHlz8ccff4ikpCSxb98+8cMPP0j78PT0FJ999pkQQoj09HQBQCxdulSkpKRI+cbPz0+8/vrrsljbtm0rJkyYUORxjhs3TnTv3l22LC/eq1evilu3bokuXbqITp06iVu3bgkhhMjMzBROTk7ilVdeESdOnBC7d+8WjRo1EgDEvn37pO3cvn1bWFlZifPnzwshHuVpW1tbMWDAABEfHy927dolbG1thb+/vxgxYoRISEgQW7ZsEebm5mLTpk1F/3JK4Nq1a+Lll18Wbdq0KbCuTZs2YsaMGTrZD1FxsIJFenXmzBkBQPz666/Sst69e4uePXsKIYR4//33RdeuXWXvuX79ugAg/vjjDyFEwQqWh4eHmDlzpuw9nTt3FkOHDhVCCBEeHi7s7OzE1atXC43p8QpWYdvP07hxYzFnzhzp9csvvyyreOR39OhRAUBcuHChQGxt27aVXmu1WuHk5CQCAgKkZQ8fPhQWFhZi48aN0rI///xTVK5cWWRnZ0txAxCnTp2SykyaNEmYmZnJjvXxfx6uXbsmAIj9+/cXGnNgYKB49dVXCz2OtLQ0kZmZKaysrMSSJUuKPG4hhHj77bdFhw4dnliGiMgYPKmCNXXqVKHRaIQQQpw7d06oVCpx6dIlWZn+/ftLlZ28CsuuXbuk9VFRUQKA+Pvvv4UQjyptI0eOLDKexytYQggBQJYLhBBi48aNolKlSuL+/ftCCCESEhIEABEfH1/kdvv37y+GDx8uW5YX78mTJ0WrVq1EQECAtE0hhPjmm29E5cqVxb1796Rla9euLVDB+uGHH0TTpk2l16GhocLa2lrcvHlTWjZw4EDh5OQksrKypGWdO3cWY8eOLTLm4pgyZYqwsbERAETbtm0LzfUDBgzgl36kV+wiSHpVr149dO7cGcuXLwfwqOver7/+KvUJj4uLw8GDB2FnZyc9atSoAQDSrEePu337NpKTkwt0eejQoQMSEhIAAEePHkWzZs3g5ORUptjHjBkjdW+8fv06tm/f/sS+7Hnd8qytrQusa9asmfSzSqWCi4sLmjZtKi2zsLBA5cqVpS4awKPuF3369IG5ubm0zMrKCl5eXtJrV1dXuLm5yY7V1dVV2k6VKlUwYsQI+Pv7o0+fPliwYAEuX74slY2Li8OaNWtk5z/v3J4/fx4JCQnIysoqcrKSPBqNpsCAZCIipRFCQKVSAQCOHDkCIQQaNWok+4zcuXNngfz0+Ge8h4cHAEifw6+//jq+//57NG/eHO+88w4OHDhQ4rj69+8PS0tLqSv28uXL0bp16ydO0HT//v1C8xEA+Pv7o3r16tiyZYuszKlTp9CkSRNoNBppWZs2bQq8//HugXlq1qwJR0dH6bWrqyvq168PS0tL2bLH81xpTJ48GUePHsXu3buhVqsxbNgwCCFkZZiTSN9YwSK9GzVqFLZt24br169j5cqVqFKlivTBrNVq0adPHxw7dkz2OHv2bIExUE+TlxR1JTAwEJcuXUJUVBTWrl0LZ2dn+Pv7F1k+r5Jz48aNAussLCwKxFrYssf7pm/btq1AAnu8slXc7axYsQJ//PEHOnXqhB9//BFeXl7SGAKtVovRo0fLzv1ff/2Fs2fPlmhijevXr8PZ2bnY5YmIjFFCQoI0SZFWq4VKpcLhw4dln5GJiYnSl4Z5Hv8czstFeZ/DvXr1wqVLl/DOO+8gIyMDffr0wciRI0sUl4WFBYYPH47ly5cjJycHq1evfurkFU5OToXmIwDo27cvoqKiZONxiysnJwc7d+4skJ9Kk+dKw8nJCfXr10ePHj2wYcMG/Prrr4iKipKVYU4ifWMFi/TuxRdfhLW1NdasWYPly5dj+PDh0oduq1atcPLkSXh6eqJu3bqyh729fYFtOTg4wMPDA7///rtseVRUFBo1agQAaNmyJY4fP17s+3xYWloiNze3wPIqVarghRdewPLly7F8+XIEBQXBzKzoP6FnnnkGDg4OUktaWZw/fx7nzp1Dz549y7wtAGjevDmmTp2K/fv3o0uXLli1ahWA/85//nNft25daDQaNGzYEFZWVvjtt9+euP0TJ06gVatWOomViMgQTpw4gV27duHFF18E8CiXCCGQmppa4POxWrVqJdq2k5MTAgMDsXLlSixbtgyrVq1CVlZWoWUtLCwKzUmjR4/Gvn37sHjxYmRmZmLw4MFP3GfLli2LzEcfffQRXnvtNXTr1k12G5AGDRrgxIkTstafP//8U/beAwcOwM7ODj4+Pk/cvz7kVdbyn0vmJNI3VrBI7zQaDYYMGYKwsDCcP39e9q3b2LFjcevWLQwaNAh//PEHLly4gD179iA4OBiZmZmFbm/y5MmYN28e1q9fjzNnzuCDDz5AZGQk3nnnHQDAkCFD4OLigv79+yMyMhIXLlzAjz/+KM0imF+tWrVw6dIlHDlyBBkZGbIP6jFjxmDt2rX466+/nnrzYDMzM3Tv3r3AN2mlsX37dnTr1q3MN0S+ePEipk2bhujoaFy6dAn79u3D8ePHpcro1KlT8eeff+K1117D0aNHce7cOfz0008ICQkBANjb22P8+PGYPn06VqxYgfPnz+PPP//EkiVLpH3cu3cPcXFxOqsMEhGVt6ysLKSmpiI5ORl//fUXFixYgC5dusDb21vKJfXr18fQoUMxYsQIbNq0CRcuXEBsbCzmzZtXollTP/jgA2zbtg1nz55FYmIitmzZgjp16sDKyqrQ8rVq1cJvv/2G1NRUWQuUl5cXOnTogMmTJ+PFF1+Eg4PDE/fr7++PxMTEIu9T9fHHHyMkJATdu3fHX3/9BeBR/lSr1RgzZgwSEhKwZ88ezJ49G8B/LXPbt29Hv379in38uhITE4NFixbhr7/+wqVLl7B371688sorqFWrFjp06CCVS0pKwj///AM/Pz+9x0gVFytYZBCjR4/GjRs30K5dOzRs2FBantcaZWZmhp49e6Jx48YYO3YsrKysikw+48aNw+TJkzFlyhQ0adIEW7duxebNm9G8eXMAgK2tLQ4cOIDq1asjICAATZo0QWhoaJFdCAcOHIjevXujW7ducHZ2xvr166V1Xbp0QfXq1dGlS5di3dsqODgY33//faHfPpZEYd0DS8PGxgZnzpzBSy+9hPr16yMoKAhDhw7F1KlTATwaN3Dw4EEkJSWhc+fOaN68OaZPnw5XV1dpG3PmzMHUqVPx0UcfoWHDhhg4cCCuXLkird++fTtq1qyJjh07ljleIiJ92LNnD9zd3VGzZk1069YNP/74I8LCwnDw4EHY2tpK5VasWIGRI0diypQpaNCgAfr27YuDBw/C09Oz2PuysrLCe++9h+bNm6N9+/bIzMzEjh07iiw/f/587Nu3DzVq1EDLli1l60aNGoWHDx8W695WTZs2RevWrbFhw4Yiy8yePRtjxoxBt27d8Ndff8He3h47duzAyZMn0bJlS0yePBlhYWEA/htfXNj4K33QaDTYtGkTunbtCi8vL4waNQrNmjVDZGSkbBzZ+vXr4efnV6LfEVFZqUT+kYBEVKT79++jWrVqWLhwIYYOHVqs9/j6+uKNN95AYGBgqfaZkZEBd3d3/P3333BzcyvVNvSpdevWmDBhAoYMGWLoUIiITNqnn36KZcuW4cyZM8Uqv2vXLowfPx4JCQlQq9Wl2uf27dsxYMAApKen4++//8Zzzz2Hq1evFhhfZQyysrJQr149rF+//on3/yLSNbZgERWDVqtFeno6Zs2aBY1GI90AsTi+/fbbMg3ivX79OhYsWKCIylV6ejpefPFFvPLKK4YOhYjIZN25cwcnT57El19+ifHjxxf7fT179sTYsWNlvQ6eZtWqVYiMjERSUhJ++uknTJgwAQEBAXByckJ2djYWLlxolJUrALh06RLee+89Vq5I79iCRVQMSUlJqF27NqpXr45ly5axLzcRERnMiBEjsH79evTr1w/r168vMKOsLs2dOxeLFy9GSkoK3Nzc0KdPH3z66aeFTjxFRI+wgkVERERERKQj7CJIRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERERERDrCChYREREREZGOsIJFRERERESkI6xgERERERER6QgrWERERET0f/buPC6K+v8D+GtZrpVLlEtQUVPxToVUNNNSwQSttNJMBFOgNNM008oCzbRLq1+phd888iy1LLNM8iawFDXxtrzyAMQDQQURPr8/gNGBBTn2muH1fDz2sTszn515zyj72fd+jiEiA2GCRURERKqwePFiaDQanD592tyhmNy+ffvQvXt3ODo6QqPRYP/+/eYOyeT0XYOS/ydq8v8RMh0mWEREREQVlJ2djZiYGPTr1w/u7u7QaDR4//33yyyfm5uLKVOmwMfHBzqdDp06dcJvv/1m0JgKCgowePBgXLx4EbNnz8bSpUvh6+t73/f9888/GDt2LJo2bQp7e3vUrVsX/fv3x65duwwSV3EyU/ywtrZG/fr1MWLECJw/f94gxyhW1WtAZAzW5g6AiIiISCkyMjIwffp01K9fHx06dEB8fHy55SMiIrBmzRqMGzcOzZs3x5IlSxASEoLNmzejR48eBonpwoULOHHiBD799FNER0dX6D2LFy/G6NGjUbt2bQwfPhxNmjTB2bNnsWDBAnTv3h3fffcdnnrqKYPEFxsbiwceeAA5OTn4448/8M0332D79u04ePAgatWqZZBjlHUNwsLCMGTIENjZ2RnkOEQVwRYsMrpGjRohIiJCsftXo9jYWGg0GqSmpprkePw3IiK1qFevHs6fP4///vsPcXFx5Zb966+/sGrVKsyYMQMff/wxoqKisHnzZjRq1AiTJk0yWEzp6ekAABcXlwqV//bbb/HCCy8gODgY//77L95//31ERUVhxowZ+Pvvv1G7dm1ERkbixo0bBokvODgYw4YNw6hRo7Bo0SKMHz8ep06dwo8//ljmeyp77LKugVarhb29PTQaTeUDJ6oiJlgkGTBgAOzt7XHt2rUyy7zyyivQaDQ4fvy46QKrpAsXLiA2Ntbg/c+3bdsm6+qg1Wrh6emJZ555BkeOHDHosX755RfExsYadJ/GpsSYiUj9/v77b/Tr1w/Ozs5wcHBAz549sXPnzlLltm3bhoCAANjb2+OBBx7AV199Jf0YdS87Ozt4e3tX6Nhr1qyBlZUVoqKipHX29vYYOXIkdu/eXaFxQPeLPyIiAv7+/gCAESNGQKPRoGfPnmXu7+LFi4iOjkbr1q2xatUq6HQ62XYvLy+MHTsWly9fNnhXxmKPPfYYAODUqVMA7v7od/DgQYSFhaFOnTpo06aNVL4616CiY64uXryIUaNGwcvLC3Z2dmjZsiXmz59vwLOmmoRdBEkybNgwrF+/HmvXrsXIkSNLbc/Pz8e3336Lhx56CM2bNzdDhPodO3YMVlZ3fyu4cOECpk2bhkaNGqF9+/YGP96YMWPQpUsX5OXlYd++fYiLi8OWLVuQkpJS4Ur3fn755RfMnTtXUQlLeTGX/DciIjKFI0eOoHv37nBwcMCkSZNgb2+PBQsWoHfv3oiPj8cjjzwCoHByhL59+8LLywuxsbEoKCjA9OnT4ebmVq3j79u3Dw888ABcXV1l6zt16iRtb9SoUbXij46OxgMPPIB33nkHUVFR6N69Ozw9Pcvc58cff4zMzEysWrWqzG5zbdu2BQCj/Zj677//AgDq1q0rWz948GA0btwYM2bMwO3btwEY5xqUlJ6eji5duiA/Px+jR4+Gh4cHNm/ejNGjR+Py5cuYOnWq4U6eagQmWCQZMGAAnJ2dsWLFCr0JVnx8PNLT0/HWW2+ZIbqymbpf9cMPP4whQ4YAKPylrFmzZnjllVewZMkSvPHGG9Xa940bN+Dg4GCIMM2y/7Kw7zsRmcNbb72FnJwcJCcno1mzZgAKP7dbtGiBCRMmYM+ePQCAmJgYaDQaJCQkoH79+gCAZ599Fi1btqzW8S9evIh69eqVWl+87sKFC9WOPzAwEDY2NnjnnXcQGBiIYcOGlbk/IQS++eYbNGrUCMHBwWWWK27Vys/Pv+85VkRmZiYyMjKkMVjTp0+HTqdDaGiorFyLFi2wdu1a2TpDXwN9pk6ditzcXKSkpMDd3R0A8OKLLyIyMhIzZ87Eyy+/jNq1a1f9AlCNw5+USWJvb49BgwZh27Ztej/0ly9fDq1WKyUXubm5mDZtGpo1awY7Ozv4+Pjg1Vdfxc2bN+97rIyMDERFRcHLywv29vZo06YNFixYUKqcEAJz585F+/btodPp4Obmhj59+si6Btw7vmfbtm146KGHANztJqDRaBAbG4sFCxZAo9Fg7969pY7zf//3f9BoNFXq6leyqwMAbNq0CT169ICjoyMcHR3Rt2/fUl0WIyIiYG9vj9OnT0vJbUhICCIiIjB37lwAkHVJPH36NE6fPg2NRoPFixeXiqP4PIvdr8sFAFy5cgVDhw6Fi4sLXF1dER0djezsbFmZn376Cf3790f9+vVhZ2cHX19fTJo0CTk5ObJzKStmQP8YrIr8Hyg+3/fffx8LFizAAw88ADs7Ozz00EPYvXu3/n8QIiIUJge//fYb+vfvL30xBwA3NzdEREQgOTkZaWlpyM/Px++//44BAwZIyRUANG3aFI8//ni1Yrh165beH5js7e2l7dWNvzKOHDmCjIwMPP744+WOSTp58iQAGKxXRt++feHu7o4GDRpgyJAh8PT0xPr16+Hj4yMr99JLL8mWjXENShJCYM2aNQgJCYFGo0FGRob0CAoKwq1bt/Dnn39W6xhU87AFi2SGDRuGRYsWYdWqVZgwYYK0/ubNm1i3bh369OkDDw8PCCHw1FNPYfv27YiMjESrVq1w5MgRzJs3D4cOHcJvv/1W5od3Tk4OHn30URw9ehRjxozBAw88gHXr1iEqKgqXL1/GlClTpLJRUVH43//+h+DgYEREREAIgcTEROzYsQPdu3cvte+WLVti+vTpsm4CANCuXTv4+vrilVdewbJly9CxY0fZ+5YtW4aAgIAq/VpZsqvDihUrMGzYMPTp0wezZs1Cbm4u4uLi0L17d+zevRstWrSQ3ltQUICgoCB06tQJH330EaytrdGqVStcuHAB8fHxWLp0qVTW3d0dly5dqnR8+rpcFBsyZAh8fHwwc+ZM7N+/H3Fxcfjvv//wyy+/SGUWLVoEOzs7vPLKK3BxccGuXbvwySef4L///sOqVasAANHR0WXGrE9l/g8AhQOys7OzER0dDY1Ggw8//BADBw7EyZMnYWNjU+lrQkTqd+nSJdy8eRN+fn6lthV/1p8+fRoFBQW4desWmjZtWqqcvnWVodPpkJubW2p98Q9UJcc/3aui8VemK9y5c+cA4L7Tl2/ZsgUA0K1btwrvuzz/93//h5YtW8Le3h4NGzZEgwYN9H5HeOCBB2TLxrgGJV26dAlXr17FwoULsXDhQr1liifQIKooJlgk07NnT9SvXx8rVqyQJVg//vgjsrOzpWb3lStXYuPGjdi6datsmtmAgAAMGzYM8fHxCAoK0nuMuLg4HDx4EIsXL0Z4eDgAYPTo0QgODkZsbCwiIyNRt25dbNu2Df/73/8wevRoqXUEAF599VUIIfTu29PTE48//niZ3QSeeOIJrFy5Eh999BG0Wi2Awj7mu3fvxmeffVaha5SVlYWMjAxpDNb48eOh0WgwaNAg3LhxAy+//DIiIiJkH9QjR46En58fpk+fjhUrVkjr8/LyEBoaijlz5siO0bx5c8THx5eKvyoJlr4uF8V8fHzwyy+/SBVdvXr18O677+L3339H7969ARS2XN47jW50dDSaNWuGqVOn4qOPPkKDBg0QGBhYZsz6VPT/QLH//vsPJ06ckMYx+Pn54YknnsBvv/1WqosJEZGlqFevHs6cOVNq/cWLFwEYroWossrraXLhwgX8/PPP6Ny5s8HGWz/00EPo0qXLfcuVl3AaS0FBAQDgueeewwsvvKC3TOvWrU0ZEqkAuwiSjJWVFZ577jkkJyfLBrcuX74cDg4OePLJJwEA3333HZo3b47WrVvLmtN79OgBjUaDrVu3lnmMDRs2wN3dXfZFXKvVYvz48cjNzcXvv/8OoHD2JQCYNm1aqX1UdbrV8PBwpKamSscACluvrK2tpa6P9xMVFQV3d3d4e3sjJCQE+fn5WLlyJQICAhAfH4+rV69i6NChsuuSn5+P7t27670uo0ePrtK5VFTJLhf3evnll2XX8pVXXgEA/Pzzz9K64uSqoKBA6kf/8MMPQwiht7tlRVT0/0CxQYMGyQaJF7dMFndjISIqyd3dHbVq1cKxY8dKbTt69CiAwu7LHh4esLe3xz///FOqnL51ldG+fXv8+++/uHr1qmx9cZez8iZiqmj8lVGcMKWkpJRZZuzYscjNzcV7770nrYuLi8Njjz2GJ598Eh4eHvDz85PGrxmTMa6BvmM4OTnhzp076N27t96HvnF0ROVhgkWlFH/pXb58OYDCsTK//fYbnnzySWmChOPHj+PYsWNwd3eXPRo0aAAhRLnN6WfOnEHTpk2lFqRi9zb3A4Vd7zw9Pas9i9O9goKC4OXlhWXLlknrli9fjuDgYHh4eFRoH2+99Rbi4+Oxfft2/Pvvvzh58iQGDx4M4O6MS3369Cl1bb7//vtS18XKyqralcP9lOxyca97+7QDhf3aXV1dZdPZHjx4EP369YOjoyNq164Nd3d3qdUyMzOzSjFV9P9AsYYNG8qWi5Otkl9aiIiKabVa9O3bF+vXr5e6cgOFY0+XLFmCgIAAeHp6QqvVonfv3vjpp5+kLnRAYXL166+/ViuGp59+GgUFBbL7ZeXm5mLRokXw9/dH48aNqx1/ZTRq1AidOnXCTz/9VCpBEkJg0qRJ+P777/Haa6+hV69e0raDBw9iz549GDNmDC5evIjg4GC8/vrrlTp2VRjjGug7xtNPP41169bh77//LrW9Kj1HiNhFkEpp164d2rZti5UrV2LatGn47rvvcOfOHVlrQ0FBAVq1alVmtzpzdXu4H61Wi+effx5ffvklbty4gb///hsnT57EzJkzK7yPNm3aSN3nSiruarB48eJSg3f1sbGxgbV1xf8My2q5K2+mp+p0ucjMzMSjjz4KBwcHvPfee2jatCl0Oh3Onz+PiIgI6XyNrWQiVqysrqJERAAwY8YMbNq0CQ8//DDGjBkjTfF97do1qZcEUDgpUHG5l156CQUFBfjiiy/QunVrvV+6v/jiC1y7dk26b+TWrVtx584dAIUtQMU3u+3cuTOeeeYZTJ06FRkZGWjWrBm++eYbnDp1CvHx8QaLvzLi4uLQo0cPaXrzli1bIjU1FatWrcLRo0cxadIkfPDBB7L3HDx4EG+++Sb69OkDoPCH2IEDB1bp+JVljGtQ0vvvv49t27YhMDAQkZGRaN26Na5evYr9+/fjhx9+kE3qRFQRTLBIr2HDhmHy5MnYvXs3li9fDk9PT+mDFShsFUlOTkavXr0q3V3P19cX+/btQ35+vuyLc8nm/gceeAAbN27EpUuXypwsQZ/7xRMeHo7Zs2fjhx9+QGJiIpydnfHEE09U6hzKUtxa5O7uXmYSVhFlnUNxy03Jm0Hr6+NfESdOnJD1sc/IyMDVq1elf4OtW7ciIyMDa9askY210/fFoDL/Dyr6f4CIqDpatmyJhIQEvPHGG/jggw9QUFCAgIAALFiwQLoHFgD4+/vj119/xWuvvYZ33nkHDRo0wPTp03HkyBG93dM+/vhj2efupk2bsGnTJgCF9WdxggUA33zzDd555x0sW7YMV65cQZs2bbB+/Xo8+uijBou/Mh588EHs2bMHM2bMwLfffou0tDQUFBSgfv36SEhIQNeuXUu95+DBg/jyyy+l5fT0dIP2LimPMa5BSR4eHvjzzz/x7rvvYt26dZg/fz7q1KmDli1bYvbs2QY5BtUwgkiP//77T1hZWYkBAwYIAOKVV16RbV+yZIkAIObOnVvqvTk5OeL69evSsq+vrwgPD5eW/+///k8AEN988420Lj8/X/Tp00fY2dmJjIwMIYQQW7duFQDE6NGjSx2joKCgzP0fOXJEABBz5swp8/w6dOggHnvsMVG3bl0xcuTIsi/EPYrjWblyZZllMjMzRe3atUW3bt1Ebm5uqe3p6enS6/DwcGFnZ6d3P5MnTxYAxJUrV0ptc3NzE0899ZRs3cSJEwUAERMTI62LiYkRAMTFixdL7aN4W9++fWXX8u233xYAxKZNm4QQQvz0008CgNi6datUJj8/X/Tt21cAEIsWLapQzFX9P3Dq1CkBQMyaNavUPkueLxGRoT3xxBOiadOm5g7D6CIjI4WVlZXYvn17qW1paWmlPttHjBghpkyZYsoQiRSFLVikV/369dGjRw/89NNPAFBqZrhhw4ZhzZo1GDNmDLZv3y5NenDs2DF89913WL16NXr27Kl335GRkYiLi8PIkSOxb98+NGnSBOvWrcPmzZsxa9Ysafa4nj17IiIiAvPmzcO///4r3Y8kKSkJ7dq1w5tvvql3/w888ABcXV0xf/58ODo6wsnJCW3atJHdA2r48OF49dVXAQBhYWHVulb3cnZ2xpdffonnn38eHTp0wHPPPQdPT0+cPXsWGzduROvWrfXew6qkgIAAAIWTUDz++OOwtrZG//794eDggFGjRuH999/HqFGjEBAQgB07dsgmJKmM8+fPo1+/fggNDcXff/+NBQsWICgoSGqt7NatG+rWrYvw8HCMHTsWNjY2WLNmTal7Zd0v5pIq+n+AiMhUbt26JetSfeLECfzyyy/STKdq9sknn2Dbtm0ICwvDgQMHZC1wBw8ehFarxbfffouRI0di1apV+O2336o8yRFRjWDuDI8s19dffy0AiObNm+vdnpeXJz766CPRpk0bYWdnJ2rXri06duwo3nnnHXH58mWpXMnWCyGEuHTpkhg1apTw8PAQtra2olWrViIuLq7UMfLz88WcOXNEq1athK2trahbt67o06ePSEhIKHf/P//8s2jbtq2wsbHR29KRlpYmrK2tha+vr6wFpzwVacEqtmPHDtG3b19Ru3ZtYW9vL5o0aSLCwsJEYmKiVKa8Fqz8/Hwxfvx44enpKTQajQAgTp06JYQQ4ubNm2LkyJHCxcVFODk5iWeffVakp6dXqQXr4MGD4rnnnhPOzs7CxcVFjBo1SmRmZsrK7tq1S3Tr1k3UqlVLeHh4iJdeekkcOHCgVAtWeTFX9f8AW7CIyFS8vLzElClTRFxcnHjrrbdEnTp1hIODgzh+/Li5QzOrzz77TIwcOVIEBwcLR0dH0b17d3Ho0CFzh0Vk0TRCcJQ41TxXr16Fl5cXXnvtNdlUtEREVDONGDECW7duRWpqKuzs7BAYGIiZM2eWujF9TRMVFYUOHTqUe8sPIpJjF0GqkZYsWYLbt2/XiK4fRER0f4sWLTJ3CBbp4MGDGDp0qLnDIFIUVd8Hq2fPntBoNLJHyZvJXr16FWFhYXBxcYGLiwvCwsJKzdCWkpKCHj16QKfTwcfHB9OnT+f00Aq1ZcsWzJ07F9OnT0doaKjB7lJPRHQ/rJNIiQ4dOoQWLVqYOwwiRVF9C9aIESNk9zgqeU+goUOHShMQAMCoUaMQFhaG9evXAwCuX7+OPn364JFHHsHu3btx9OhRjBgxAg4ODpg4caLpToQMYvr06UhMTERgYCDmzZtn7nCIqIZhnURKU9UbyhPVZKpPsGrVqgUvLy+9244cOYKNGzciISEBgYGBAICvvvoK3bt3x7Fjx+Dn54fly5fj5s2bWLJkCXQ6Hdq0aYOjR49izpw5mDBhQqXvAUXmtW3bNnOHQEQ1GOskIiL1U3UXQQBYtWoV3Nzc0Lp1a7z22mvIysqStiUlJcHR0VF2U71u3brBwcEBiYmJUpnu3bvLfmUMDg7GhQsXcPr0aZOdBxERKR/rJCIi9VN1C9bQoUPh6+sLb29vHDp0CG+88QYOHDgg3W09NTUV7u7usl/8NBoNPDw8kJqaKpWpX7++bL+enp7StsaNG8u2xcXFIS4uDgBw9OhR9lsmolJOnz6NjIwMc4dBJmaOOgmQ10vJySkAvACwpYuIAEAAyIAQN8wdiKooLsGaOnXqfafV3rp1K3r27ImoqChpXdu2bdGkSRN07twZe/fuNdq0q1FRUdJxAwICsOevv4xyHCJSroBOncwdAhmIpddJgLxe0mh8AEwF4AnAxmjHJCIlyAOQBiDWzHGoj+ISrPHjx2PYsGHllmnYsKHe9QEBAdBqtThx4gQ6duwILy8vXLp0CUII6RdDIQTS09OlPvJeXl5IS0uT7ad4uax+9EREVDMot06yAeABJllENVUegHRzB6Faikuw3Nzc4ObmVqX3pqSkID8/H/Xq1QMABAYGIjs7G0lJSVKf96SkJNy4cUNaDgwMxOTJk5GTkwN7e3sAQHx8PLy9vdGoUaPqnxARESmWsuskGxS2ZCnuqwARVcsdFLZckbFohEpvnvHvv/9i+fLl6NevH9zc3HD48GFMnDgROp0Ou3fvhlarBQA8/vjjOHfunNQ/PSoqCo0aNZKmxM3MzISfnx969uyJqVOn4vjx44iIiEBMTMx9p8RlF0Ei0iegUyfs2bPH3GGQCVlCnQQUdxGMBVAfgE/Ro64xTpmILNZlAOeLHucAxEKI8+YNSWVU+7OVra0tNm/ejM8++wzZ2dlo0KABQkJCEBMTI1VkALBixQqMHTsWwcHBAIABAwbgiy++kLa7uLggPj4eY8aMQUBAAFxdXTFx4kRMmDDB5OdERETKZPF1kiNU/I2AiAAUNlxlmzuImkG1LViWgC1YRKQPW7DIXMpswaqNwgeTLCJ1ugPgWtGDLVhGx49SIiIiKvxG4AZ+MyBSmzsAeGcQk+LHKBERERVikkWkLkyuzIIfoURERHQXkywidWByZTb8+CQiIiI5JllEysbkyqz40UlERESlMckiUiYmV2bHj00iIiLSj0kWkbIwubII/MgkIiKislmDU7gTKUHxVOxkdvy4JCIiovJZA3ATgPUdc0dCRPrcsQYyNOaOgoowwSIiIqL7s74DR7dr0FrnmzsSIrpH/h0tsjNqA7AxdyhUhAkWERERVYjWOh+17a5BCyZZRJYgH1pcQ21zh0ElMMEiIiKiCtMiH7VxFdZMsojM6g60uAZXc4dBejDBIiIiokqxRj5qgy1ZRObClivLxgSLiIiIKk2LfDghC7a4be5QiGqU27BFFpzMHQaVgwkWERERVQm7CxKZ1h1ocRlu5g6D7oMJFhEREVWZNfLhiCwmWURGdgdaZLPlShGYYBEREVG1cEwWkXFxzJWyMMEiIiKiatNKSRZvRkxkSPmwZnKlMEywiIiIyCC0uANXtmQRGUw+tLjK5EpxmGARERGRwWjZXZDIINgtULmYYBEREZFBFc8uyCSLqGryeRNhRWOCRURERAZXfJ8s63wmWUSVcUer5X2uFI4JVgXNmzcPH330ES5evIjWrVvj008/Rffu3c0dFhER1VBKqJes8/PhfDkPnPeCqIKsget1AWjNHQhVBxOsCvj2228xbtw4zJs3Dw8//DDmzZuHxx9/HIcPH0bDhg3NHR4REdUwiqqX7gC4AiZZRPdjDaCOuYMgQ2CCVQFz5sxBREQEIiMjAQCff/45Nm7ciPnz52PWrFlmjo6IiGoaxdVLTLKIysfkSlWYYN3H7du3kZycjNdee022PigoCImJiWaKquKmvfuu9Lpjhw7oHxqqt1zcggW4mJoqLUeOGgXvevVKlbtw8SIW/O9/0nI9Ly9EFVXwJa3/+Wfs3bdPWg4NCYF/x473jRMAYt5+W2+55L178fOGDdIyz4nnpE9lzolIaRRbL90BcBnAc+LuuhNxwJ/R+ss/vgeo6393+Rd/4Mre0uXqdAT6Jd9dvpwM/Bqgf5+dvwKaRd1d3hUF/LNAf9lhQr68TKO/XNNIoEvc3WWeE89Jn/LOSQugrv63kTIxwbqPjIwM5Ofnw9PTU7be09MTv//+u5miIjK9/fv3Y9asWRgxYgS8vLzMHQ5RjaXoeqnkfBe3AWRWsGx2GWVt9byvrH3eLrF8q5yyJZVV7paeY/CcSjPCOe3/D+hZuzZ27NiBdu3aFa5U4jm5lPEeUiwmWAYWFxeHuLjCX0guXbpktjhKtjQQVdfgIUOQm5uL1atXY+zYseYOh4gq6N56CbhhvkCCxP3LEFXC4IVAZmYmnnnmGRw7dszc4RBJrMwdgKVzc3ODVqtFWlqabH1aWpreX/GjoqKwZ88e7NmzB+7u7qYKk8joTvzzDwDg8uXLZo6EqGarTr0EOJgoSiLjO1H0O/aJEyfMGwhRCRohBH9Suo/OnTvjwQcfvOcXQKB58+YYNGhQuYOJAwICsOevv0wRYikVHStDVFF2Oh1u374NW1tb5N4q2deCKiOgU6eiL7tEVVPVekmj8QEQC6A+AJ+iR12gNgA3AF5lPNwAeOXBxesy6tpdhhsyUBeXURvX4Ik01MY11MVl1EUG3HBZeu2afw3OaXlAOoD2Jb5uzC1jrAxRBdmNA27nA7ZaIPczc0dTDS4APO4+rnva4Kq2Ni4X/TUV/rUVvr6G2kiDJ66h9t2/uNy6yEytC6TaABkAUst4ZAC4BhQOhjxf9DgHIBZCnDf1WasauwhWwIQJExAWFoZOnTqhW7du+PLLL3HhwgW8+OKL5g6NyGTy8vJkz0RkPqyXiIC8fPkzkaVgglUBgwcPxuXLlzFjxgxcvHgRbdq0wS+//AJfX19zh1YmtliRodnY2OD27duwsbExdyhVwlZdUhMl1kvYpClsyaroBAxE92GjLWzBslHqTXnHlGjV3c9WXbVgglVBo0ePxujRo80dBpHZsAWLyLKwXqKaji1YZKk4yQURVUhxy5VSW7CIiEhdiluuFNuCRarFBIuIKoQtWEREZEnYgkWWil0EiahClD4Gi2OuiIjURfFjsIpn0rx3FkFSBbZgEVGFsAWLiIgsCVuwyFKxBUulkvfulS37d+xopkhILZTegkVEZuYTCTgDKL6N3uEF5oyGVEDxLVikWkywVOrnDRtky0ywqLrYgkVE1dI6Tr7MBIuqiS1YZKnYRZCIKsTa2lr2TEREZE7WVvJnIkvB/5JEVCFswSIiIkuSVyB/JrIU/ClapTp26GDuEEhlNBoNhBDQaJR5p3mOSyQys3NxheOvbps7EFILDQBR9KxIrSILn3UAnADksdusWjDBUqn+oaHmDoFURumTXHBcIpGZHY4G0gFkmjsQUgvFT3LxaIlxifuZYKkFuwgSUYWwiyAREZXJDD/Zc5ILslRMsIioQopbrpTagkVEREZ0x/SHLG65UmwLFqkWuwgSUYUovQWL4xKJiNRF8S1Yh4q6CNqicBwWqQYTLCKqEKWPweK4RCIidVH8GKxt0YXPLgA8ih6kCuwiSEQVovQWLCIiUhfFt2CRarEFS6XiFshnoomKjDRTJKQWSm/BIiIz67IHyANQ/GV4dYA5oyEVUHwLFqkWEyyVupiaau4QSGXYgkVE1eLsb+4ISGXYgkWWil0EiahCOIsgERFZEs4iSJaKCRYRVQhbsIiIyJKwBYssFbsIqlTkqFHmDoFURuljsDgukcjMkvyBqwCyzR0IqYXix2A9s6fwWQvABsBJjktUCyZYKuVdr565QyCVUXoLFsclEplZ1l7gCoBMcwdCaqH4FiwPjktUK3YRJKIK4RgsIiKyJByDRZZK1QlWbGwsNBqN7OHl5SVtF0IgNjYW3t7e0Ol06NmzJw4dOiTbx9WrVxEWFgYXFxe4uLggLCwM165dM/GZEJmf0luwiMyNdRJVC/sclVKpFizrEs9ERqT6/2Z+fn7Ytm2btKzV3v2Z48MPP8Ts2bOxePFi+Pn5Yfr06ejTpw+OHTsGJycnAMDQoUNx9uxZbNy4EQAwatQohIWFYf369SY9DyJzU/oYLI5LJEvAOomq7I65A7A8lRqDdafEsyX4rqiLoCMAVwA6cwZDhqT6BMva2lr2C2ExIQQ+/fRTTJkyBYMGDQIALFmyBB4eHlixYgWio6Nx5MgRbNy4EQkJCQgMDAQAfPXVV+jevTuOHTsGPz8/k54LkTkpvQWL4xLJErBOIjIcxY/BurS38Pk2Cr+RM8FSDVV3EQSAkydPwtvbG40bN8aQIUNw8uRJAMCpU6eQmpqKoKAgqaxOp8MjjzyCxMREAEBSUhIcHR3RtWtXqUy3bt3g4OAglbFUFy5elD2Iqsva2lr2TESVV1PrJACAU0egTkfAvehBVE3WVvJnIkuh6m9KnTt3xuLFi9GiRQukp6djxowZ6Nq1Kw4dOoTUohnFPD09Ze/x9PTE+fPnAQCpqalwd3eHRqORtms0Gnh4eEjvLykuLg5xcXEAgEuXLhnjtCpkwf/+J1uOefttM0VCaqH0FiwiczNHnQTI6yXghmFPqjICk+XLczX6yxFVUF6B/JnIUqg6wXr88cdly126dEGTJk2wZMkSdOnSxSjHjIqKQlRUFAAgIID3MyD10Gg0EELIvtwRUcWZo04C5PWSRuNjtOMQmZoGgCh6JrIkNapR1dHREa1bt8aJEyekPvBpaWmyMmlpadI2Ly8vXLp0CUIIabsQAunp6Xr70BOpGadpJzIs1klE1cNp2slS1agEKycnB0ePHkW9evXQuHFjeHl5IT4+XrZ9586dUv/2wMBAZGdnIykpSSqTlJSEGzduyPrAW6J6Xl6yB1F1Kb2LIMclkqWpSXUSAOB6MnA5GUgvehBVk+InuSgej1inY+EYRVINVXcRfO2119C/f380bNgQ6enpePfdd3Hjxg2Eh4dDo9Fg/PjxmDlzJlq0aIHmzZtjxowZcHR0xNChQwEALVu2RN++fREdHS31X4+OjkZoaKjFz9YUFRlp7hBIZZQ+TTvHJZK51eQ6CQCwKwBIB5Bp7kBILSo1TbslerbEDw372dlRLVSdYJ07dw7PPfccMjIy4O7uji5dumDXrl3w9fUFALz++uu4desWxowZg6tXr6Jz587YtGmTdL8RAFixYgXGjh2L4OBgAMCAAQPwxRdfmOV8iMxJ6S1YRObGOonIsBTfgkWqpeoEa9WqVeVu12g0iI2NRWxsbJllXF1dsWzZMgNHRqQ8Sm/BIjI31klEhqX4FixSLVUnWERkOEpvweJYRCIidVF8C1bxWEQtAP52qSomS7AKCgqwbds2bN++HadPn8atW7fg7u6Ojh07IigoCA0aNDBVKERUBUpvweK4RLoX6yQiA7MGcMe0h1R8C9bqotv5uADwKHqQKhh9FsFbt27hvffeQ4MGDRASEoJNmzYhOzsbtra2OHXqFKZNm4bGjRujX79+2LVrl7HDIaIqUnoLFhHAOonIaEycXAEqaMEi1TJ6C1azZs0QGBiIuLg4BAUF6f31+8yZM1ixYgUGDx6MqVOnIpK/NFfb+p9/li33Dw01UySkFkpvwSICWCeZVauvgMYAbhctb4s2ZzSkAopvwSLVMnqCtXHjRrRp06bcMr6+vnjjjTcwceJEnDlzxtgh1Qh79+2TLTPBoupiCxapAeskM6ofJV9mgkXVxBYsslRG7yJ4v4rsXra2tmjWrJkRoyGiqir+pZ8tWKRkrJOI1KO45YotWGRpTD6L4O3bt3Hw4EGkp6ejoKBAtq1fv36mDoeIKogtWKRGrJOIlIstWGSpTJpgxcfHIywsDOnp6aW2aTQa5OfzL8RQQkNCzB0CqYzSx2BxXCKVxDrJxA5FAVkAbpk7EFILxY/B6vlV4bMtAB2ADHabVQuTJlhjxoxBaGgo3n77bXh6ekKj0Zjy8DWKf8eO5g6BVEbpLVgcl0glsU4ysfMLgHQAmeYOhNRC8S1YrUuMS2SCpRomTbAuXryIN998E76+vqY8LBEZgNJbsIhKYp1EpGyKb8Ei1TL6JBf3Cg0NRWJioikPSUQGovQWLKKSWCcRKZviW7BItUzagvXll1/i+eefR3JyMtq0aVPql/Dhw4ebMhwiqgSlt2BxXCKVxDqJSNkU34K1taiLoA6Ak1kjIQMzaYL122+/YfPmzfjll19Qq1YtWX93jUbDyozIgim9BYvjEqkk1klEyqb4FqzDCwqfXQB4FD1IFUzaRfC1117Dyy+/jKysLGRnZyMrK0t6XL9+3ZShEFElWVtby56JlI51EpGyWVvJn4kshUm/KV27dg0vvvgiHBwcTHnYGmnau+/KlmPefttMkZBaKL0Fi6gk1kkmFiTky3M5ayNVT16B/JnIUpg05x80aBB+//13Ux6SiAykuPsUp7ImtWCdRKRsmhLPRJbCpC1YTZo0wVtvvYUdO3agXbt2pQYUT5gwwZThEFElKH2SC6KSWCcRKZviJ7kg1TJpgrVw4UI4OTkhMTGx1NS4Go2GlRmRBWMXQVIb1klEBmQN4I5pD6n4SS5ItUyaYJ06dcqUh6vROOaKDE3pLVgcl0glsU4ysU0aIB1AprkDIaMwcXIFqKAFa0yJcYn72dlRLTjvChFVCFuwiIjIkrAFiyyV0ROsGTNm4MaNGxUq+8cff2D9+vVGjoiIqqK45UqpLVhEAOskIjUpbrlSbAsWqZbRE6x///0XDRs2RFRUFNavX4+LFy9K23JycrB371783//9Hzp16oSwsDC4uroaOyQiqgK2YJEasE4iUg+2YJGlMnqCtWjRImzbtg0ajQbDhw9H/fr1YW1tDZ1OBwcHBwQEBOCbb77BqFGjcOTIETz88MMV3veOHTswYMAA+Pj4QKPRYPHixbLtQgjExsbC29sbOp0OPXv2xKFDh2Rlrl69irCwMLi4uMDFxQVhYWG4du2arExKSgp69OgBnU4HHx8fTJ8+HUKU6DdLpHJKb8GKeftt2YNqJtZJROqh+BasuZrCxzJN4RhFUg2TTHLRtm1bfPXVV5g/fz4OHDiAM2fO4NatW3Bzc0P79u3h5uZWpf1mZ2ejTZs2GD58OIYPH15q+4cffojZs2dj8eLF8PPzw/Tp09GnTx8cO3YMTk5OAIChQ4fi7Nmz2LhxIwBg1KhRCAsLk7qFXL9+HX369MEjjzyC3bt34+jRoxgxYgQcHBwwceLEKl4RIuVhCxapBeskInVgCxZZKpPOImhlZYX27dujffv2Btlfv3790K9fPwBARESEbJsQAp9++immTJmCQYMGAQCWLFkCDw8PrFixAtHR0Thy5Ag2btyIhIQEBAYGAgC++uordO/eHceOHYOfnx+WL1+OmzdvYsmSJdDpdGjTpg2OHj2KOXPmYMKECRZ709XkvXtly/4dO5opElILpc8iSFQS6yQT84kEnAHcKlo+vMCc0ZAKKH4WQVIt1c4ieOrUKaSmpiIoKEhap9Pp8Mgjj0j3O0lKSoKjoyO6du0qlenWrRscHBxkZbp37w6dTieVCQ4OxoULF3D69GnTnEwV/Lxhg+xBVF1swSKquppeJwEAWscBXeKAR4seRNXEFiyyVKpNsFJTUwEAnp6esvWenp7SttTUVLi7u8t+8dNoNPDw8JCV0bePe49xr7i4OAQEBCAgIACXLl0y3AkRmZnSx2ARmZO56iRAXi8BFZtBkUgJFD8Gi1RLtQmWuURFRWHPnj3Ys2cP3N3dzR0OkcGwBYtIme6tlwAHc4dDZDBswSJLZdIxWKbk5eUFAEhLS0PDhg2l9WlpadI2Ly8vXLp0CUII6RdDIQTS09NlZdLS0mT7Ll4uLmOJOnboYO4QSGWUPgaL4xLJnGp6nQQAOBdXOP7qtrkDIbVQ/BisVpGFzzoATgDyOC5RLcyWYKWlpcHd3R1WVsZpRGvcuDG8vLwQHx+Phx56CEDhPU527tyJjz76CAAQGBiI7OxsJCUlSX3ek5KScOPGDWk5MDAQkydPRk5ODuzt7QEA8fHx8Pb2RqNGjYwSuyH0Dw01dwikMkpvwSo5FpEJFt2LdZIJHI4G0gFkmjsQUgvFt2CVHIu4nwmWWpi0i2BeXh5ef/11ODk5wcfHRxqQO3nyZMybN6/S+8vOzsb+/fuxf/9+FBQU4OzZs9i/fz/Onj0LjUaD8ePH44MPPsD333+PgwcPIiIiAo6Ojhg6dCgAoGXLlujbty+io6ORlJSEpKQkREdHIzQ0FH5+fgAKp8ytVasWIiIicPDgQXz//fd4//33LX+2JiID4xgsUhvWSUTKxjFYZKlMmmBNmzYN69evx7Jly2BnZyet79SpU6kbMlbEnj170KFDB3To0AG3bt1CTEwMOnTogHfeeQcA8Prrr+PVV1/FmDFjEBAQgIsXL2LTpk3S/UYAYMWKFXjwwQcRHByM4OBgPPjgg1i6dKm03cXFBfHx8bhw4QICAgIwZswYTJw4ERMmTKj6hSBSIKW3YBGVxDqJSNkU34JFqmXSLoIrV67EwoUL0aNHD1k3jDZt2uD48eOV3l/Pnj3LvXu9RqNBbGwsYmNjyyzj6uqKZcuWlXuctm3bYseOHZWOj0hNlD4Gi+MSqSTWSUTKpvgxWIeKugjaonAcFqmGSROsCxcuwNfXt9T6O3fu4M6dO6YMhYgqSektWByXSCWxTiJSNsW3YG2LLnx2AeBR9CBVMGkXwdatW+v91e27776Dv7+/KUMhokqytraWPRMpHeskImWztpI/E1kKk35TiomJwbBhw/Dff/8hPz8fq1evxtGjR7FixQpsKDHDF1VP3AL5TDRRkZFmioTUQuktWEQlsU4ysS57gDwAxa0NqwPMGQ2pQF6B/JnIUpg0werfvz++++47zJw5E1ZWVpg2bRo6duyI9evXo3fv3qYMRfUupqaaOwRSGY1GI7s/D5HSsU4yMWcztgpaA2Cvz8ox9jUzwP41AETRM5ElMXlfn+KZkYhIWZQ+yQWRPqyTaggmV5Vn7GtmgP0rfpILUi2T9lpt0qQJLl++XGr9tWvX0KRJE1OGQkSVxC6CpDask2oQDh2tPGNfMwPsX/GTXJBqmfQj5/Tp08jPL/1XkJubi/Pnz5syFNWLHDXK3CGQyii9BYvjEqkk1kkmluQPXAWQbYZjswWr8tiCZXzP7Cl81gKwAXCS4xLVwiQJ1vfffy+93rBhA1xcXKTl/Px8bN68GY0aNTJFKDWGd7165g6BVEbpLVgcl0jFWCeZSdZe4AqATHMHQmqh+BYsD85WqlYmSbCefvppAIWD5EeOHCnbZmNjg0aNGmH27NmmCIWIqkjpLVhExVgnEamD4luwSLVMkmAVFBTOn9m4cWPs3r0bbm5upjgsERmQ0luwiIqxTqrh1DijoBrPqQJKtWDV0OtAlsekY7BOnTplysMRkQEpvQWL4xKpJNZJNZQav4Cr8ZwqoFQLltKuw3dFXQQdAbgC0JkzGDIkk8+rc/XqVfz66684e/Ysbt++Ldv2zjvvmDocIqogpbdgcVwi6cM6iUi5FD8G69LewufbKPxGzgRLNUyaYO3atQshISGws7PDpUuX4OPjg4sXL8LOzg6NGjViZWZAFy5elC3zyyVVl9JbsIhKYp1kYk4dC1sYbIuWi79cElURx2CRpTJpgjVp0iQ8//zz+Oyzz+Ds7IwtW7bAwcEBzz33XKmBxlQ9C/73P9lyzNtvmykSUgult2ARlcQ6ycQCk+XLczXmiYNUQ/EtWKRaJr3R8IEDB/Dyyy9Do9FAq9UiNzcXnp6e+OCDDxAbG2vKUIiokopbrtiCRWrBOolI2YpbrtiCRZbGpAmWra2t9NrT0xNnzpwBADg6OuLChQumDIWIKoktWKQ2rJOIlI0tWGSpTNpFsGPHjti9ezeaN2+Onj17YurUqUhLS8OyZcvQrl07U4aievW8vMwdAqmM0sdgcVwilcQ6ycSuJwN5APhlmAxE8WOw3DsWPjsCcAIAjktUC5MmWO+99x6ysrIAADNmzMDw4cMxduxYNG/eHIsWLTJlKKoXFRlp7hBIZZTegsVxiVQS6yQT2xUApAPINHcgpBaKb8F6tsS4xP0cl6gWJk2wAgICpNfu7u749ddfTXl4IuP6+WdgxozKvadLF+DTT40SjqEpvQWLqCTWSUTKpvgWLFItk98Hi6hSLl8GwsKApUuBunXNHU35ZsxA7O7dlXpL79278bBCEiylt2ARERlENoAlAMJR2LWLzEbxLVikWiZNsK5cuYK33noLmzdvRnp6OgoKCmTbr1+/bspwSAmWLAF++w345hvg1VfNHU2NpvQWLI5LpJJYJ1GV/AngSNFzLzPHUsMpvgUrvaiLoBaAMqtWKoNJE6yRI0di3759iIqKgre3NzQa9jWlcghxt/vcp58C48cD/D9jNkpvweK4RCqJdRJVmgCwtej1VgCPAeB/G7NRfAvW6qJuyi4APIoepAomnaZ98+bN+Pbbb/Hmm28iIiIC4eHhskdl7dixAwMGDICPjw80Gg0WL14s2x4REQGNRiN7dOnSRVYmNzcXY8eOhZubGxwcHDBgwACcO3dOVubs2bPo378/HBwc4ObmhldeeQW3b9+udLxUSTt3AplFo6GvXgUSEswbTw1nbW0teyZSOtZJVGn/ArhV9PpW0TKZjbWV/JnIUpj0m5KHhwccHQ3XYTk7Oxtt2rTB8OHDMXz4cL1levfujaVLl0rL9973BADGjx+PH3/8EStXrkTdunUxYcIEhIaGIjk5GVqtFvn5+QgJCUHdunWxc+dOXL58GeHh4RBC4PPPPzfYuRja+p9/li33Dw01UyTV8NlnwI0bha9v3ixsxere3awh1WRKb8EiKol1kom1+gpoDKA4F9wWbc5oqmYL7sZ/u2i5qfnCqenyCuTPRJbC5NO0v/POO1iyZIlBKrV+/fqhX79+AAp/GdTHzs4OXmWMvcjMzMTXX3+NRYsWoU+fPgCApUuXwtfXF7///juCg4OxadMmHDp0CGfOnEGDBg0AAB9++CFGjRqF9957D87OztU+D2PYu2+fbNniE6wnnwTWr5evs7Ut7CYIFD7/8gugLdHRun9/YN06U0RY42k0Gggh2I2KVIN1konVj5IvW3qC9RWAlBLrtCjsJoii50MAXi5Rpi0ACz81tdCg8J+BtRJZGpM2qs6YMQObNm2Ch4cHWrZsiXbt2skexpCQkAAPDw80b94ckZGRSE9Pl7YlJycjLy8PQUFB0roGDRqgZcuWSExMBAAkJSWhZcuWUkUGAMHBwcjNzUVycon7F1DVvfce0LAhYG9/d13JLi/3LtvbA76+he8jkyie3EKpk1wQlcQ6ico1AIAr5D9Flxzrc++yDYA6Re8jkyie3EKxk1yQapm0Bevpp5825eHQt29fDBw4EI0bN8bp06cxdepUPPbYY0hOToadnR1SU1Oh1Wrh5uYme5+npydSU1MBAKmpqfD09JRtd3Nzg1arlcrcKy4uDnFxcQCAS5cuGenMVKh1a+DgQWDkSGDDhsIugWWpVQsIDQX+9z/AwcF0MdZw7CJIalMT6iRAXi8BNwx+XqpVD8BUAMsBHMTdroH62AJoA+B5AHbGD40KKX6SC1ItkyZYMTExpjwchgwZIr1u27Yt/P394evriw0bNmDgwIFGOWZUVBSiogq7Qdx7E0tTCw0JMduxq8zBAVi1Cvjqq8Ip2XNzS5exswM+/hiIZv8LU1P6NO2qGJdIBlUT6iRAXi9pND5GO859HYoCsnB3kgglsAPwAoCdANYCuKOnjDWApwBwiLDJKX6a9p5fFT7bAtAByOB3G7WoUdOBeXt7o379+jhx4gQAwMvLC/n5+cjIyIC7u7tULi0tDd2LJlPw8vLCH3/8IdtPRkYG8vPzy+xHbwn8O3Y0dwhV16FDYSJVVoKl5HNTMKW3YCluXCKpXk2qkwAA5xcA6QAyzR1IFTRA4TemshKshqYNhwopvgWrdYlxiUywVMPoY7CcnZ2RkZEBAHBycoKzs3OZD2PLyMjA+fPnUa9ePQCAv78/bGxsEB8fL5U5d+4cjhw5gq5duwIAAgMDceTIEdk0ufHx8bCzs4O/v7/RY66RkpOB4i/xGk1hl8DiiRXy8oA9e8wXWw3GMVikBqyTqErOovR4q2L5RdvJ5DgGiyyV0VuwPv/8czg5OQEAvvjiC4PuOzs7G//88w8AoKCgAGfPnsX+/ftRp04d1KlTB7GxsRg0aBDq1auH06dP44033oCHhweeeuopAICLiwtGjhyJ119/HR4eHtKUuO3atUPv3r0BAEFBQWjdujWGDx+O2bNn4/Lly5g0aRIiIyMte7YmJdu5E7h1q3AiC09PYM6cwi6D6emF6xMSgJdeMneUNY7SW7CIANZJVEX/AMhDYWLlBGAQgDUo7PKYV7SdXQRNTvEtWKRaRk+w7r1ZY1Vu3FiePXv24NFHH5WWY2JiEBMTg/DwcMyfPx8pKSn45ptvcO3aNdSrVw+PPvoovvvuO6lyBYBPP/0U1tbWGDx4MG7duoVevXrhm2++gbZoOnCtVosNGzZg9OjR6NatG3Q6HZ5//nl89NFHBj0XuseffxZOxz5gwN2JLPr0KZwA4/vvC7dboi5d0Hv37kq9xXyj9CpP6WOwFDkukQyOdRJVyRkU9vlpi7sTWbRA4QQY+wGcNldgNZvix2BtLeoiqENh4k6qoRGi+EZDZGgBAQHY89df5g5DeUJDgUGDgBEjSm9btAhYuxYoMWEBGZ+VtbV0H6yCO/oGIlBFBXTqhD3s6kpmUDjJRSyA+gB8ih51gdoA3AB4lfFwA+CVBxevy6hrdxluyEBdXEZtXIMn0lAb11AXl1EXGXDDZem1a/41OKflFY69Kn4ocQzWPAAdAATq2ZYEYB+A0SaNiABYvXz3PlgFhm2QNi0XAB53H9c9bXBVWxuXi/6aCv/aCl9fQ22kwRPXUPvuX1xuXWSm1gVSbYAMAKllPDIAXAOAywDOFz3OAYiFEOdNfdaqZvQWLCsrqwrfmDQ/n228hPKTpxEj9CdeZHRKb8EiAlgnURWVlzwFQn/iRUan+BYsUi2jJ1jfffedVJmlpaXhnXfewVNPPYXAwMJPo6SkJKxbtw7Tpk0zdig1yrR335Utx7z9tpkiIbXgGCxSA9ZJZhRUosPM3IolukRl4RgsslRGT7DuvZHjgAEDMGvWLERGRkrrXnjhBXTq1Anr1q3D6NFsXyeyVGzBIjVgnUSkHmzBIktl9Gna77VlyxbZAOBijz76KLZt22bKUIioktiCRWrDOolI2diCRZbKpAmWm5sb1qxZU2r9mjVrZDdVJCLLw/tgkdqwTiJSNt4HiyyV0bsI3mv69OkYMWIEtm7dKvV337VrF37//Xd8/fXXpgxF9TjmigxN6S1YHJdIJbFOMrFNGuXOIkgWSfEtWGNKjEvcz3GJamHSBGv48OFo0aIFPvvsM/z0008AgJYtW+KPP/5A586dTRkKEVUSx2CR2rBOIlI2jsEiS2WyBCsvLw/Dhg3DzJkzsXz5clMdlogMROktWET3Yp1EpHyKb8Ei1TLZGCwbGxts2rSpwvcfISLLYm1tLXsmUjLWSUTKZ20lfyayFCb9pjRw4EB8//33eO2110x5WCIyAKW3YHHMFZXEOolI2fIK5M+KU3wvOBcAHkUPUgWTJlgNGzbEjBkzsHPnTgQEBMDBwUG2fcKECaYMh4gqQaPRQAjBX/xJNVgnESmbBoAoeiayJCZNsBYvXgxXV1ccOHAABw4ckG3TaDSszAwoee9e2bJ/x45mioTUgpNckNqwTjIxn0jAGcCtouXDC8wZDakAJ7kgS2XSBOvUqVOmPFyN9vOGDbJlJlhUXUrvIkhUEuskE2sdJ19mgqUu1gDumPaQnOSCLJXZhgVmZ2fjxo0b5jo8EVUSbzRMasY6iaiaTJxcAbzRMFkukydYc+fORcOGDeHi4gJnZ2f4+vpi3rx5pg6DiCqJLVikRqyTiJSLLVhkqUzaRXDmzJmYNWsWXnvtNTz88MMAgJ07d2LKlCm4fv06pkyZYspwVK1jhw7mDoFURuljsDgukUpinWRi5+IKx1/dNncgpBaKH4PVKrLwWQfACUAeu82qhUkTrC+//BJxcXF47rnnpHW9evVCs2bN8Oabb7IyM6D+oaHmDoFURuktWByXSCWxTjKxw9FAOoBMcwdCaqH4FqxHS4xL3M8ESy1M2kUwPT0dDz30UKn1nTp1QlpamilDIaJK4hgsUhvWSUTKxjFYZKlMmmA1b94cK1asKLV+xYoV8PPzM2UoRFRJSm/BIiqJdRKRsim+BYtUy6RdBGNjY/Hss89ix44d6NatGwDgjz/+wPbt27F69WpThkJElaT0MVgcl0glsU4iUjbFj8E6VNRF0BaF47BINUyaYA0cOBB//vknPvnkE/z8888AgJYtW+Kvv/5CB375IbJoSm/B4rhEKol1EpGyKb4Fa1t04bMLAI+iB6mCSRMsAPD398eyZctMfVgiqialt2AR6cM6iUi5FN+CRapl8vtg5ebmYuHChXjttdcwadIkLF68GLm5uZXez6xZs/DQQw/B2dkZ7u7u6N+/Pw4ePCgrI4RAbGwsvL29odPp0LNnTxw6dEhW5urVqwgLC4OLiwtcXFwQFhaGa9euycqkpKSgR48e0Ol08PHxwfTp0yGEqHTMphS3YIHsQVRdSm/BItLHUHUSwHrpvrrsAR7fAzxT9CCqJsW3YJFqmTTBOnz4MJo1a4YJEybgzz//xK5duzB+/Hg0b94cR44cqdS+tm3bhtGjRyMxMRFbtmyBtbU1evfujStXrkhlPvzwQ8yePRuff/45du/eDQ8PD/Tp0wdZWVlSmaFDh2Lv3r3YuHEjNm7ciL179yIsLEzafv36dfTp0weenp7YvXs3PvvsM3z00UeYM2dO9S+IEV1MTZU9iKqLswiS2hiyTgJYL92Xsz9Q1x/wKHoQVRNnESRLZdIuguPGjUOHDh2wdOlSODs7AyisKIYNG4bx48fjt99+q/C+SpZdunQpXFxc8Mcff6B///4QQuDTTz/FlClTMGjQIADAkiVL4OHhgRUrViA6OhpHjhzBxo0bkZCQgMDAQADAV199he7du+PYsWPw8/PD8uXLcfPmTSxZsgQ6nQ5t2rTB0aNHMWfOHEyYMAEajcZAV4fIsrEFi9TGkHUSwHqJjMwawB1zB2FZKtWCVXz9eB3JBEzagvXHH39g5syZUkUGAM7OznjvvfeQkJBQrX1nZWWhoKAArq6uAIBTp04hNTUVQUFBUhmdTodHHnkEiYmJAICkpCQ4Ojqia9euUplu3brBwcFBVqZ79+7Q6e5O7xIcHIwLFy7g9OnT1YqZSEnYgkVqY8w6CWC9RAbGpKCUSrVg3SnxTGREJm3Bsre3L9WPHAAyMzNhb29frX2PGzcO7du3l37xSy3qFufp6Skr5+npifPnz0tl3N3dZb/2aTQaeHh4SO9PTU1F/fr1S+2jeFvjxo1l2+Li4hAXVzjt5qVLl6p1TtUROWqU2Y5N6qT0FqySYxGjIiPNFAlZCmPWSYBl1kvAjWqfV5Ul+QNXAWSbLwRSF8WPwSoei6gFYAPgZIA5oyEDMmmC1b9/f0RGRmLBggXo0qULgMJf4qKjozFgwIAq73fChAlISEhAQkICtFrzdsSNiopCVFQUACAgwHx/KN716pnt2KROSp9FkGMRqSRj1UmA5dZLGo2P+QLJ2gtcAZBpvhBIXRQ/iyDHIqqWSbsIfvbZZ2jWrBm6d+8Oe3t72Nvbo0ePHmjevDk+/fTTKu3z1VdfxcqVK7FlyxY0adJEWu/l5QUASEtLk5VPS0uTtnl5eeHSpUuymZeEEEhPT5eV0bePe49BVBMovQWLqCRj1EkA6yUiU1F8CxaplkkTrNq1a+PHH3/E8ePH8f333+P777/HsWPH8MMPP8DFxaXS+xs3bpxUibVo0UK2rXHjxvDy8kJ8fLy0LicnBzt37pT6tgcGBiI7OxtJSUlSmaSkJNy4cUNWZufOncjJyZHKxMfHw9vbG40aNap0zERKZW1tLXsmUjpD10kA6yUiU7K2kj8TWQqTfVM6e/YsvLy8YGtri6ZNm6Jp06YAgNu3b+Ps2bNo2LBhpfY3ZswYLF26FOvWrYOrq6vUN93R0RGOjo7QaDQYP348Zs6ciRYtWqB58+aYMWMGHB0dMXToUABAy5Yt0bdvX0RHR0v906OjoxEaGgo/Pz8AhdPlTps2DREREZg6dSqOHz+O999/HzExMZypiWoUpbdgcVwi3cvQdRLAeonI1PIK5M+K811RF0FHAK4AdOUVJiUxSc6/atUq9OvXT+82IQT69euH1atXV2qf8+bNQ1ZWFnr16oV69epJj48//lgq8/rrr+PVV1/FmDFjEBAQgIsXL2LTpk1wcnKSyqxYsQIPPvgggoODERwcjAcffBBLly6Vtru4uCA+Ph4XLlxAQEAAxowZg4kTJ2LChAmVvApEylb8xU2pX+C869WTPajmMkadBLBeIjI1TYlnxbm0t/BxZW/hGEVSDZO0YMXFxWHSpEmwtbUttc3Ozg6TJ0/G/Pnz8cwzz1R4nxW5Y71Go0FsbCxiY2PLLOPq6oply5aVu5+2bdtix44dFY7NEly4eFG2zC+UVF1Kn+SCqJgx6iSA9dJ9OXUsnCK7+LJf4hdKqh7FT3JBqmWSBOvIkSN4+OGHy9zetWtXTJo0yRSh1BgL/vc/2XLM22+bKRJSC6V3ESQqxjrJTAKT5ctzFdvuQBaCk1yQpTJJF8HMzEzcvn27zO25ubm4fv26KUIhoirijYZJLVgnEalDpW40TGRCJkmwGjdujN27d5e5/a+//uLMR0QWji1YpBask4jUgS1YZKlM0kVw4MCBeOutt9CnTx/UKzEW6MKFC3j77bcRERFhilBqjHq8FwoZmNLHYHFcIhVjnWQm15OBPAD8MqxO1igcY2dCih+D5d6x8NkRgBMAcFyiWpgkwZo8eTJ++OEHNG/eHMOGDZPuDXLkyBEsX74cDRs2xOuvv26KUGqMqMhIc4dAKqP0FiyOS6RirJPMZFcAkA4g09yBkFGYOLkCVNCC9WyJcYn7OS5RLUySYDk6OuKPP/7AG2+8gW+//RZXr14FUHiTx7CwMLz33nuyKWqJyPIovQWLqBjrJCJ1UHwLFqmWyW407OLignnz5mHu3LnIyMiAEALu7u6KvacOUU2j9BYsonuxTiJSPsW3YJFqmSzBKqbRaODu7m7qwxJRNSm9BYvjEkkf1klEyqX4Fqz0oi6CWgDKrFqpDCZPsIhImZTegsVxiURE6qL4FqzVAYXPLgA8ih6kCiaZpp2IlI/3wSIiIkvC+2CRpWILlkqt//ln2XL/0FAzRUJqofQWLCIys1ZfAY0BFN/jeVu0OaMhFVB8CxapFhMsldq7b59smQkWVZfSx2ARkZnVj5IvM8GialL8GCxSLZN0EezVqxe+//77MrenpaVBq+VfB5ElYwsWqQXrJCJ1YAsWWSqTJFhbt27Fs88+i5iYmDLLCCFMEQoRVRHHYJFasE4iUgeOwSJLZbIugvPnz8drr72GAwcOYNmyZXBwcJBt571HDCs0JMTcIZDKKL0Fi+MS6V6sk8zgUBSQBeCWuQMhtVB8C1bPrwqfbQHoAGSw26xamCzBeuKJJ/Dwww/jiSeeQJcuXfDjjz+iSZMmpjp8jePfsaO5QyCVUfoYLI5LpHuxTjKD8wuAdACZ5g6E1ELxY7BalxiXyARLNUw6TXvLli2xe/duNGjQAA899BB+//13Ux6eiKpB6S1YRCWxTiJSNsW3YJFqmfw+WC4uLtiwYQMiIyPRr18/fPLJJ6YOgYiqwNraWvZMpAask4iUy9pK/kxkKUzyTalkX3aNRoP3338f7du3x6hRo7BlyxZThEFE1aD0FiyOS6RirJOI1CGvQP6sOFuLugjqADiZNRIyMJMkWGXNxjRkyBC0aNECTz75pCnCIKJq0Gg0EEIodvA/xyVSMdZJREZgDeCOaQ+pASCKnhXp8ILCZxcAHkUPUgWTJFhbt25FnTp19G5r3749kpOTsWHDBlOEQkRVpPRJLoiKsU4iMgITJ1eACia5INUySa/VHj16lDtuo27duhg+fHil9jlr1iw89NBDcHZ2hru7O/r374+DBw/KykRERECj0cgeXbp0kZXJzc3F2LFj4ebmBgcHBwwYMADnzp2TlTl79iz69+8PBwcHuLm54ZVXXsHt27crFa+pTXv3XdmDqLqU3kWQqBjrJDMJEsAwAYwpepgSh45WnrGvmQH2z0kuyFIpdljgtm3bMHr0aCQmJmLLli2wtrZG7969ceXKFVm53r174+LFi9Ljl19+kW0fP3481q5di5UrV2Lnzp24fv06QkNDkZ9f+Nean5+PkJAQZGVlYefOnVi5ciXWrFmDiRMnmuxciSwBbzRMVDbWSRbODK0rimfsa2aA/fNGw2SpFPubzm+//SZbXrp0KVxcXPDHH3+gf//+0no7Ozt4eXnp3UdmZia+/vprLFq0CH369JH24+vri99//x3BwcHYtGkTDh06hDNnzqBBgwYAgA8//BCjRo3Ce++9B2dnZyOdIZFlYQsWUdlYJxGZHluwyFIptgWrpKysLBQUFMDV1VW2PiEhAR4eHmjevDkiIyORnp4ubUtOTkZeXh6CgoKkdQ0aNEDLli2RmJgIAEhKSkLLli2ligwAgoODkZubi+TkZCOfFZHlYAsWUcWxTiIyPrZgkaVSbAtWSePGjUP79u0RGBgorevbty8GDhyIxo0b4/Tp05g6dSoee+wxJCcnw87ODqmpqdBqtXBzc5Pty9PTE6mpqQCA1NRUeHp6yra7ublBq9VKZSxRzNtvmzsEUhmlt2CVHIvIvxEyJtZJemzSAOkAMs0dCKmF4luwSo5F3K/Y+RCpBFUkWBMmTEBCQgISEhKg1d79GWPIkCHS67Zt28Lf3x++vr7YsGEDBg4caJRY4uLiEBcXBwC4dOmSUY5BZA6cRZCoYiypTgLk9RJww2jHITI1ziJIlkrxXQRfffVVrFy5Elu2bEGTJk3KLevt7Y369evjxIkTAAAvLy/k5+cjIyNDVi4tLU3qI+/l5YW0tDTZ9oyMDOTn5+vtRx8VFYU9e/Zgz549cHd3r86pEVkUpbdgEZmCpdVJgLxeAhyqeGZElkfxLVikWopOsMaNGydVZC1atLhv+YyMDJw/fx716tUDAPj7+8PGxgbx8fFSmXPnzuHIkSPo2rUrACAwMBBHjhyRTZMbHx8POzs7+Pv7G/iMiCwXx2ARlY91EpFpcQwWWSrFdhEcM2YMli5dinXr1sHV1VXqe+7o6AhHR0dkZ2cjNjYWgwYNQr169XD69Gm88cYb8PDwwFNPPQUAcHFxwciRI/H666/Dw8MDdevWxYQJE9CuXTv07t0bABAUFITWrVtj+PDhmD17Ni5fvoxJkyYhMjKSszVRjaL0FiyOuSJjYp1EZHqKb8GaWzTmygWAR9GDVEGxLVjz5s1DVlYWevXqhXr16kmPjz/+GACg1WqRkpKCJ554As2bN0d4eDj8/PyQlJQEJycnaT+ffvopnnrqKQwePBjdunWDo6Mj1q9fL/Wb12q12LBhA2rVqoVu3bph8ODBGDRokHQcopqCLVhEZWOdRGR6bMEiS6XYFiwhyr8LvE6nK3VfEn3s7Ozw+eef4/PPPy+zTMOGDfHzzz9XOkZzSt67V7bs37GjmSIhtVB6CxaRMbFOqgCfSMAZwK2i5cMLzBkNqYDiW7BItRSbYFH5ft6wQbbMBIuqi7MIElG1tI6TLzPBomriLIJkqRTbRZCITIstWEREZEnYgkWWigkWEVUIx2AREZEl4RgsslTsIqhSHTt0MHcIpDJKb8HiuEQiMzsXVzj+6ra5AyG1UHwLVqvIwmcdACcAeew2qxZMsFSqf2iouUMglbG2tkZeXh6srZX5scFxiURmdjgaSAeQae5ASC2srYC8gsJnRXq0xLjE/Uyw1EKp/yWJyMSU3oJFRETqklcgfyayFEywiKhCNBqN7JmIiMicNCWeiSyFMvv6EJHJKX2ado5LJCJSF8VP036oqIugLQrHYZFqMMEiogpRehdBjkskIjIiawB3THtIxU9ysS268NkFgEfRg1SBXQSJqEI4TTsREZXJxMkVwGnayXKxBUul4hbIZ6KJiow0UySkFkpvwSIiM+uyB8gDUNzasDrAnNGQCii+BYtUiwmWSl1MTTV3CKQySh+DRURm5uxv7ghIZRQ/BotUi10EiahC2IJFRESWhC1YZKmYYBFRhXAMFhERWRKOwSJLxS6CKhU5apS5QyCVUXoLFsclEplZkj9wFUC2uQMhtVB8C9YzewqftQBsAJzkuES1YIKlUt716pk7BFIZpY/B4rhEIjPL2gtcAZBp7kBILRQ/BsuD4xLVil0EiahClN6CRURE6qL4FixSLSZYRFQhHINFRESWhGOwyFKxiyARVYjSW7A4LpGISF0U34L1XVEXQUcArgB05gyGDIkJFhFViNLHYHFcIhGRuih+DNalvYXPt1H4jZwJlmowwVKpCxcvypb55ZKqS+ktWERkZk4dgTsAbIuWi79cElWR4luwSLWYYKnUgv/9T7Yc8/bbZoqE1ELpLVhEZGaByfLluRrzxEGqofgWLFItxU5yMXfuXLRr1w7Ozs5wdnZGYGAgNmzYIG0XQiA2Nhbe3t7Q6XTo2bMnDh06JNvH1atXERYWBhcXF7i4uCAsLAzXrl2TlUlJSUGPHj2g0+ng4+OD6dOnQwhhilMksihswSIqH+slItNiCxZZKsUmWPXr18cHH3yAvXv3Ys+ePXjsscfw5JNP4sCBAwCADz/8ELNnz8bnn3+O3bt3w8PDA3369EFWVpa0j6FDh2Lv3r3YuHEjNm7ciL179yIsLEzafv36dfTp0weenp7YvXs3PvvsM3z00UeYM2eOyc+XyNw4iyBR+VgvEZkWZxEkS6XYLoJPPPGEbPm9997D/PnzkZSUhLZt2+LTTz/FlClTMGjQIADAkiVL4OHhgRUrViA6OhpHjhzBxo0bkZCQgMDAQADAV199he7du+PYsWPw8/PD8uXLcfPmTSxZsgQ6nQ5t2rTB0aNHMWfOHEyYMAEajeV2b6jn5WXuEEhllN6CxXGJZGysl+7jejKQB4CtDepkjcIxdiak+BYs946Fz44AnACA4xLVQrEJ1r3y8/OxevVqZGdno2vXrjh16hRSU1MRFBQkldHpdHjkkUeQmJiI6OhoJCUlwdHREV27dpXKdOvWDQ4ODkhMTISfnx+SkpLQvXt36HR3p3UJDg7G22+/jdOnT6Nx48YmPc/KiIqM1Lt+/c8/Y+++fdJyaEgI/Dt21Ft22rvvypbLGseVvHcvfr6nG0zHDh3QPzRUb9m4BQtwMTVVWo4cNUrvF90LFy/KxpHV8/LiOZn5nKytrZGXlwdr67sfG0o6J45LJFNivaTHrgAgHUBmifU9vwJaR91d3hoFHF6gfx9jSnSFLGscV6tI4NG4u8uH4oBt0frLPrMH8PC/u/ydv/4JONw7As/eM44sPRlYHaB/nzXxnO5Nrkx0TtZWe5FXAFjf2x9LSf9Oz5YYl7jfgn8goUpRdIKVkpKCwMBA5OTkwNHRET/88APatm2LxMREAICnp6esvKenJ86fPw8ASE1Nhbu7u+zXPo1GAw8PD6QWfbFMTU1F/fr1S+2jeJu+iiwuLg5xcYV/WJcuXTLQmRKZn9JbsIhMwdLrJeCGQc6zyrQAXEqssy2xrNNTpixllSs53bVtOWVLdi9zROG02SU56nlfWfvkOckZ6ZzyCgpf5hXc8z4lnROplqITLD8/P+zfvx+ZmZlYs2YNwsPDsW3bNrPGFBUVhaiowl9DAgLK+MWESIHs7e2Rk5MDe3t7c4dSJfW8vKRWudq1a5s3GFItS6+XNBof8wViDaCunvUlv5A6AfCo4D7LKuek5xhllS05rNQV+r8dldynTTn75DnJGemc7G2BnNuFz9L7lHRO15MB56JWuVsnK3hgUgJFJ1i2trZo2rQpAMDf3x+7d+/GJ598grfeegsAkJaWhoYNG0rl09LS4FU0NsnLywuXLl2CEEL6tVAIgfT0dFmZtLQ02TGLl704xolqmFs3zPzLt4HUrl0boSEh5g6DVIr1UhmsAdQpY1vJ32ycUThWqyLK+uLqrOcYZZUt+U2oDvTf8LXkOuty9slzkjPSOd36U085JZ5T7kngfLTCv5XTvTRCRXO7PvbYY/D29sbSpUvh7e2NsWPH4s033wQA5OTkwMPDAx999JE0mLhVq1b4448/pP7uiYmJ6NatG44ePQo/Pz/Mnz8fkydPRnp6uvSr/cyZMzF37lycO3fuvoOJAwICsOevv4x70kSkOAGdOmHPnj3mDoNMwNLqpcIWrFgA9QH4FD3qArUBuAHwKuPhBsArDy5el1HX7jLckIG6uIzauAZPpKE2rqEuLqMuMuCGy9Jr1/xrcL6cZ/LJD4gUyxq4XtcGV7W1cbnor6nwr63w9TXURho8cQ217/7F5dZFZmpdINUGyACQWsYjA8A1ALgM4HzR4xyAWAhx3gwnq16KzZWnTJmCkJAQNGjQAFlZWVixYgW2bduGDRs2QKPRYPz48Zg5cyZatGiB5s2bY8aMGXB0dMTQoUMBAC1btkTfvn0RHR0t9U2Pjo5GaGgo/Pz8ABROlztt2jRERERg6tSpOH78ON5//33ExMRY9kxNRERkcqyXSruj1eK6vm6BRFSmO1rOO690ik2wUlNTMWzYMKSmpsLFxQXt2rXDr7/+iuDgYADA66+/jlu3bmHMmDG4evUqOnfujE2bNsHJ6W7n2BUrVmDs2LHSewYMGIAvvvhC2u7i4oL4+HiMGTMGAQEBcHV1xcSJEzFhwgTTniwREVk81kty+dAiC06lJycgovvK5x+Ooqmqi6ClYRdBItKHXQTJXEzVRbA2rkHLG14RVUs+tLiG2uwiqECKbcEiIiIiy1P8pZCIqo8tWcrEBIuIiIgMIh/WuMrkisig8vl1XXH4L0ZERETVxpYrIuNhS5ayMMEiIiKiarnD5IrI6O4wyVIMJlhERERUZXegRTac7l+QiKqNSZYyMMEiIiKiKsmHFpfhZu4wiGoUdhe0fEywiIiIqNKk+1wRkckxybJsTLCIiIioUjjmisj82F3QcjHBIiIiogornC3Q1dxhEBHYkmWpmGARERFRheTfYcsVkaXJv8Mky9IwwSIiIqL7u2ON7Iza5o6CiPS5w6/0loT/GkRERFS+OwAyNABszB0JEZXljrkDoGJMsIiIiKhsdwBcM3cQRFQhTLIsAhMsIiIi0u8OgAxzB0FElcIky+yYYBEREVFpTK6IlItJllkxwSIiIiI5JldEyscky2yYYBEREdFdTK6I1INJllkwwSIiIqJCTK6I1IdJlskxwSIiIiImV0RqxiTLpJhgERER1XScip1I/ZhkmQwTLCIiopou29wBEBGpBxMsIiKiGusOgMvmDoKITIpNWcam2ARr7ty5+Oqrr3D69GkAQOvWrTF16lSEhIQAACIiIrBkyRLZezp37oxdu3ZJy7m5uXjttdewcuVK3Lp1C7169cK8efNQv359qczZs2cxZswYbNmyBTqdDkOHDsXHH38MW1tb458kEREphvLqpTwAaVU5VSJSvDxzB6Bqik2w6tevjw8++ADNmjVDQUEBlixZgieffBLJyclo164dAKB3795YunSp9J6Slc/48ePx448/YuXKlahbty4mTJiA0NBQJCcnQ6vVIj8/HyEhIahbty527tyJy5cvIzw8HEIIfP755yY9XyIismzKqpfyAKQb4KyJSLmYZBmLRgghzB2EodSpUwezZs1CdHQ0IiIikJGRgZ9//llv2czMTLi7u2PRokV4/vnnAQD//fcffH198euvvyI4OBi//vorQkJCcObMGTRo0AAAsGzZMowaNQrp6elwdnYuN56AgADs+esvw54kESleQKdO2LNnj7nDIBOwtHpJo/EBMBUK/n2ViAzqDoAZEOK8uQNRFVV8wubn52P16tXIzs5G165dpfUJCQnw8PBA7dq10aNHD7z33nvw8PAAACQnJyMvLw9BQUFS+QYNGqBly5ZITExEcHAwkpKS0LJlS6kSA4Dg4GDk5uYiOTkZjz76qOlOkoiIFMNy6yUB4IohT5WIFE81bS0WQ9EJVkpKCgIDA5GTkwNHR0f88MMPaNu2LQCgb9++GDhwIBo3bozTp09j6tSpeOyxx5CcnAw7OzukpqZCq9XCzc1Ntk9PT0+kpqYCAFJTU+Hp6Snb7ubmBq1WK5UpKS4uDnFxcQCAo0ePIqBTJ0OfdqVcunQJ7u7uZo3BEvA68BoUs4TrUDxGh9TH0uslB4fraNHiB0OfdqVYwt+gufEaFOJ1sIxrcPr0bbMeX40UnWD5+flh//79yMzMxJo1axAeHo5t27ahTZs2GDJkiFSubdu28Pf3h6+vLzZs2ICBAwcaLaaoqChERUUZbf+VFRAQwK5I4HUAeA2K8TqQMbFeuj/+DfIaFON14DVQKytzB1Adtra2aNq0Kfz9/TFr1iy0b98en3zyid6y3t7eqF+/Pk6cOAEA8PLyQn5+PjIy5LetT0tLg5eXl1QmLU0+w1JGRgby8/OlMkRERMVYLxERkaITrJIKCgqQm5urd1tGRgbOnz+PevXqAQD8/f1hY2OD+Ph4qcy5c+dw5MgRqb98YGAgjhw5gnPnzkll4uPjYWdnB39/fyOeCRERqQHrJSKimkexXQSnTJmCkJAQNGjQAFlZWVixYgW2bduGDRs2IDs7G7GxsRg0aBDq1auH06dP44033oCHhweeeuopAICLiwtGjhyJ119/HR4eHtJ0uO3atUPv3r0BAEFBQWjdujWGDx+O2bNn4/Lly5g0aRIiIyPvO1OTpbCkbiHmxOvAa1CM14GMhfVSxfBvkNegGK8Dr4FqCYUKDw8XDRs2FLa2tsLd3V306tVLbNy4UQghxM2bN0VQUJBwd3cXNjY2omHDhiI8PFycPXtWto+cnBzx8ssvizp16gidTidCQ0NLlTlz5owICQkROp1O1KlTR4wdO1bk5OSY7DyJiEgZWC8REZEQQqjqPlhERERERETmpKoxWERERERERObEBEvF5s2bh8aNG8Pe3h7+/v7YuXOnuUOqsh07dmDAgAHw8fGBRqPB4sWLZduFEIiNjYW3tzd0Oh169uyJQ4cOycpcvXoVYWFhcHFxgYuLC8LCwnDt2jVZmZSUFPTo0QM6nQ4+Pj6YPn06LKWRd9asWXjooYfg7OwMd3d39O/fHwcPHpSVUft1mDt3Ltq1awdnZ2c4OzsjMDAQGzZskLar/fyJlIx1kro+i1gnFWK9RHqZo18iGd+qVauEtbW1iIuLE4cPHxYvv/yycHBwEGfOnDF3aFWyYcMG8cYbb4jVq1cLnU4nFi1aJNv+/vvvC0dHR7FmzRqRkpIinnnmGVGvXj1x/fp1qUzfvn1Fq1atRGJiokhMTBStWrUSoaGh0vbMzEzh6ekpnnnmGZGSkiJWr14tHB0dxccff2yq0yxXUFCQWLhwoUhJSREHDhwQTz75pPD09BSXL1+Wyqj9Oqxbt0788ssv4sSJE+LYsWPizTffFNbW1uLvv/8WQqj//ImUinWS+j6LWCcVYr1E+jDBUqlOnTqJUaNGydY1bdpUTJkyxUwRGY6Dg4OsMisoKBBeXl5ixowZ0rqbN28KR0dH8eWXXwohhDh8+LAAIBISEqQyO3fuFADE0aNHhRBCzJs3Tzg5OYmbN29KZd59913h7e0tCgoKjHxWlZeVlSWsrKzETz/9JISoudfB1dVVfPnllzX2/ImUgHWS+j+LWCfdxXqJ2EVQhW7fvo3k5GQEBQXJ1gcFBSExMdFMURnPqVOnkJqaKjtfnU6HRx55RDrfpKQkODo6SveSAYBu3brBwcFBVqZ79+7Q6XRSmeDgYFy4cAGnT582zclUQlZWFgoKCuDq6gqg5l2H/Px8rFq1CtnZ2ejatWuNO38ipWCdVDM+i2p6nQSwXqK7mGCpUEZGBvLz8+Hp6Slb7+npidTUVDNFZTzF51Te+aampsLd3R0ajUbartFo4OHhISujbx/3HsOSjBs3Du3bt0dgYCCAmnMdUlJS4OjoCDs7O7z44ov44Ycf0LZt2xpz/kRKwzoJ0rKaP4tqap0EsF6i0hR7o2GimmzChAlISEhAQkICtFqtucMxKT8/P+zfvx+ZmZlYs2YNwsPDsW3bNnOHRURUY9XkOglgvUSlsQVLhdzc3KDVapGWliZbn5aWBi8vLzNFZTzF51Te+Xp5eeHSpUuyGXeEEEhPT5eV0bePe49hCV599VWsXLkSW7ZsQZMmTaT1NeU62NraomnTpvD398esWbPQvn17fPLJJzXm/ImUhnUSpGU1fhbV9DoJYL1EpTHBUiFbW1v4+/sjPj5etj4+Pl7Wx1ctGjduDC8vL9n55uTkYOfOndL5BgYGIjs7G0lJSVKZpKQk3LhxQ1Zm586dyMnJkcrEx8fD29sbjRo1Ms3J3Me4ceOkiqxFixaybTXpOtyroKAAubm5Nfb8iSwd6yT1fhaxTtKP9RJxFkGVWrVqlbCxsRELFiwQhw8fFq+88opwcHAQp0+fNndoVZKVlSX27dsn9u3bJ3Q6nZg2bZrYt2+fNMXv+++/L5ydncXatWtFSkqKGDx4sN5pUNu0aSNNg9qmTRvZNKjXrl0Tnp6eYvDgwSIlJUWsXbtWODk5Wcw0qKNHjxZOTk5i8+bN4uLFi9IjKytLKqP26zB58mSxY8cOcerUKXHgwAExZcoUodFoxC+//CKEUP/5EykV6yT1fRaxTirEeon0YYKlYnPnzhW+vr7C1tZWdOzYUWzfvt3cIVXZ1q1bBYBSj/DwcCFE4XSwMTExwsvLS9jZ2YlHHnlEpKSkyPZx5coV8fzzzwsnJyfh5OQknn/+eXH16lVZmQMHDoju3bsLOzs74eXlJWJjYy1mClR95w9AxMTESGXUfh3Cw8NFw4YNha2trXB3dxe9evUSGzdulLar/fyJlIx1kro+i1gnFWK9RPpohOBtoImIiIiIiAyBY7CIiIiIiIgMhAkWERERERGRgTDBIiIiIiIiMhAmWERERERERAbCBIuIiIiIiMhAmGAREREREREZCBMsIiIiIiIiA2GCRWQga9asgUajkZYXL14MR0dHM0ZUMSNGjMD06dOrtY+UlBT4+Pjgxo0bBoqKiIiqg3US6yQyHyZYZDYFBQV45JFH0L9/f9n6mzdvws/PDy+++GK578/KysLbb7+NVq1aQafTwdPTEz179sTKlStRUFBgzNArZPDgwTh58qTB96vRaLBmzRqD7CslJQXr1q3D+PHjq7Wftm3bokuXLpgzZ45B4iIiMjXWSVXDOomoNCZYZDZWVlZYvHgxtm7dioULF0rrJ0+ejPz8fMyePbvM9167dg2BgYFYuHAhJk2ahD179iAhIQHh4eF49913cfbsWaPFffv27QqV0+l08PDwMFochvD5559j0KBBcHZ2rva+RowYgfnz5+POnTsGiIyIyLRYJ5kf6yRSDUFkZvPnzxfOzs7i9OnT4vfffxdarVbs3Lmz3Pe89NJLolatWuK///4rte3WrVvi1q1bQgghrly5IoYPHy5q164t7O3tRa9evcTBgwdl5deuXSvatGkjbG1tRf369cWMGTNEQUGBtN3X11fExMSIESNGCBcXF/H0008LIYRYsmSJaNiwodDpdCIkJER88cUX4t4/qUWLFgkHBwdpOSYmRrRu3VqsXLlSNGnSRDg6OoonnnhCXLp0SSrz119/iT59+oi6desKJycn0a1bN5GYmCiLBYD08PX1lbb99NNPomPHjsLOzk40atRIvPnmmyI3N7fMa3jnzh3h4uIi1q1bJ1vv6+srpk2bJsLDw4Wjo6OoX7++WLVqlbh69aoYPHiwcHBwEE2bNhW//fab7H25ubnCzs5OxMfHl3lMIiJLxzqJdRJRdTHBIosQFBQkunfvLho0aCAmTZpUbtn8/Hzh6uoqIiMj77vfAQMGCD8/P7F9+3Zx4MAB0b9/f1G/fn1x8+ZNIYQQe/bsEVZWVuKdd94Rx44dE8uWLRMODg7i//7v/6R9+Pr6CicnJ/HBBx+IEydOiOPHj4tdu3YJjUYjZsyYIY4dOya+/PJLUadOnftWZg4ODuLJJ58Uf//9t0hMTBQNGzYUUVFRUpnNmzeLb775Rhw+fFgcOXJEjBkzRtSuXVtkZGQIIYRIT08XAMSCBQvExYsXRXp6uhBCiI0bNwonJyexcOFC8c8//4gtW7aI5s2bi4kTJ5Z5bfbu3SsAiHPnzsnW+/r6CldXVzF37lxx/PhxMWHCBGFnZycef/xxsWTJEnHixAnxwgsvCHd3d+lLQ7HOnTuLqVOn3vffhYjIkrFOKsQ6iahqmGCRRTh58qTQaDSiadOmIicnp9yyaWlpAoCYM2dOueWOHz8uAIjt27dL665duyacnZ3FggULhBBCDB06VDz66KOy98XExAgfHx9p2dfXV4SGhsrKPPfcc6J3796ydSNHjrxvZWZnZyeuXbsmrZsxY4Z44IEHyjyHgoIC4eXlJZYuXSqtAyBWr14tK9e9e3cxffp02boffvhBODg4yH75LLldo9GI/Px82XpfX18xZMgQaTkrK0sAEGPHjpXWnTp1SgAQu3fvlr33qaeeEsOGDSvzfIiIlIB1kn6sk4gqhmOwyCIsXLgQOp0O586dw6lTp8otK4So0D6PHDkCKysrBAYGSutcXFzQtm1bHD58WCrTrVs32fsefvhhnD9/HtevX5fWBQQElNr3vfsFUGpZH19fX7i4uEjL3t7eSE9Pl5bT09MRHR2N5s2bw8XFBU5OTkhPT79v//3k5GS89957cHR0lB5Dhw7FjRs3kJqaqvc9t27dgo2NDaysSn8MtGvXTnrt6OiIWrVqoW3bttI6T09PKd576XQ63Lp1q9xYiYgsHeukQqyTiKrG2twBEO3evRvvv/8+fvrpJ8yfPx/h4eFITEyEVqvVW97d3R21a9fGkSNHqnzMe6eurUgZBweHKh/rXjY2NqWOce/sUuHh4UhLS8Mnn3yCRo0awc7ODr169brvIOaCggLExMTgmWeeKbXN3d1d73vc3Nxw+/Zt3Lx5E7Vq1bpvnPeuK742JWfGunLlCho1alRurERElox1EuskoupiCxaZVU5ODoYPH46IiAg8/vjjiIuLwz///IMPP/ywzPdYWVlhyJAhWL58Oc6dO6d3nzk5OWjZsiUKCgqQlJQkbbt+/TpSUlLQqlUrAEDLli3xxx9/yN6fkJCA+vXrw8nJqcwYWrZsiV27dsnWlVyuioSEBIwdOxYhISFo3bo1nJyccPHiRVkZGxsb5Ofny9Z17NgRR48eRdOmTUs9rK31/47Svn17AJB+OTWEgwcPomPHjgbbHxGRKbFOkmOdRFQ1TLDIrN544w3k5ORI96rw8vLC3LlzERsbi0OHDpX5vvfeew8NGzZE586dsWjRIhw6dAj//PMPli5dCn9/f6SmpqJZs2Z44oknEB0djZ07dyIlJQXDhg2Ds7Mzhg4dCgCYOHEitm/fjtjYWBw/fhzLly/H7Nmz8frrr5cb9yuvvILff/8ds2bNwokTJ7BgwQL88MMP1b4ezZs3x7Jly3D48GHs3r0bQ4YMga2traxMo0aNsHnzZqSmpuLq1asAgHfeeQcrVqzAO++8g4MHD+Lo0aNYs2ZNuefh7u6Ojh07IiEhodpxA8Dp06dx/vx5BAUFGWR/RESmxjpJjnUSURWZexAY1Vzbt28XWq1WbN26tdS2p59+Wvj7+4u8vLwy33/t2jXx5ptvCj8/P2FnZyfc3d1Fjx49xMqVK6VBspWZEtfGxqbMKXE/+uijUsdfuHChaNCggbC3txd9+/YVn3/+eYWmxL1XyTL79+8XnTp1Evb29qJJkybim2++Ea1btxYxMTFSmZ9++kk0bdpUWFtby6bE/e2338TDDz8sdDqdcHJyEv7+/uLzzz8v8/oJIcSXX34pAgICZOv0na+Dg4NYtGiRtHzr1i0BQKxfv15aN3PmTBEcHFzu8YiILBXrJNZJRIaiEaKCozOJSHVycnLQokULLF26FN27d6/yfnJzc9GsWTOsXLmy1ABtIiKiimCdRGrBLoJENZi9vT2++eYbXLlypVr7OXPmDN566y1WZEREVGWsk0gt2IJFRERERERkIGzBIiIiIiIiMhAmWERERERERAbCBIuIiIiIiMhAmGAREREREREZCBMsIiIiIiIiA2GCRUREREREZCBMsIiIiIiIiAyECRYREREREZGBMMEiIiIiIiIyECZYREREREREBsIEi4iIiIiIyECYYBERERERERkIEywiIiIiIiIDYYJFRERERERkIEywiIiIiIiIDIQJFhERERERkYEwwSIiIiIiIjIQa3McVKPxAHAbgKZ4zb1bS7wub1vJ15UpW539W0ocVduP5p7VGs3d5Xtf31uuZJnqvLcq+zfnsQ0dN4QofJR8Xbysr4wx36uvjDmPXdn9V2S/9xBlvL7f9vLeV5n9VGebsd5rzPjaBQdj48aNekpaFo2mKYBbxUuQf9ZW5jUqWb4q7y3vRKpxCJSxvjLh3W//5jx2hd8r7lkvoLlnWaMpWi4qryn6H6+553++BqLUsrSvUuvKLm+IMhV5jyniqGiZ4mtkkGMIAc3d3d59XbQsPZe1HnrK3G99RfZ7v31Vdj+GOK6BroH0smR1XoGQKvK6rLAr+9pU9ZJZEqzC5KonAG3Rsvae11ZlrK9oGatKrq9MGasKbDNErMZ5r0ZT+B5r68KHUl5bShzVidsKBYULd+4UPiz5taXEYaBYRdEnfUHRo7qv83H3g7qq5ct6T37Ra1HBmExdvrLXbH1GBpThJoAxRa9tcLdqrOhrm6LX1gZ6Xd6xSij+1q4tsStrC1xfXCXde0oV2VdZrw1avuiv1Dq/8AHASnsH1jb50BYta63zobUu/Gyxts6H1qpoPfKhRfHrO7BGZdbnQ4uifVZ7vfxYVXmP4uPOL3p95w6s8wuKXhc+AECTDxTtsvC5+HXJ9fl6ypS1Pr+M5aocw5jxlVXGQPHl3VMF5+Xf8/rO3eJ5RcXz7nlrVV5X5/2mqpfYRZCIiIiIiMhAmGAREREREREZCBMsIiIiIiIiA2GCRUREREREZCBMsIiIiIiIiAyECRYREREREZGBMMEiIiIiIiIyECZYREREREREBsIEi4iIiIiIyECYYBERERERERkIEywiIiIiIiIDYYJFRERERERkIEywiIiIiIiIDIQJFhERERERkYEwwSIiIiIiIjIQJlhEREREREQGwgSLiIiIiIjIQKzNcdBWrXyg050zx6Gr7dKlS3B3dzd3GFVi6NgLCoDbtwtfFz8bE6+9+Zg8fqui335sbQsf1cBrbz63bt0ydwgV0qqVK3S6H8wdRpUY/P+HAFD8ec7P9XJVJ/b8Es/mYDnXXlvi+f4sJ/aqsfj4tbj7z1GiCrb42O/DVPWSWRIsnU6HPXv2mOPQ1RYQEMDYzUTJ8Ss5dkDZ8Ss5dkDZ8QcEBJg7hAphnWQ+So5fybEDyo5fybEDyo5fybEDpquX2EWQiIiIiIjIQJhgERERERER0lCu5QAAHf5JREFUGYhZEqyoqChzHNYgGLv5KDl+JccOKDt+JccOKDt+pcSulDj1UXLsgLLjV3LsgLLjV3LsgLLjV3LsgOni1wghhEmOREREREREpHLsIkhERERERGQgTLCIiIiIiIgMxOAJVlxcHB599FHUrl0bGo0Gp0+frtD71q5di1atWsHOzg6tWrXCDz/I70kihEBsbCy8vb2h0+nQs2dPHDp0yKCx5+bmYuzYsXBzc4ODgwMGDBiAc+fKv19Xo0aNoNFoSj1CQkKkMrGxsaW2e3l5GTT2qsZfkdgs9drPmjULDz30EJydneHu7o7+/fvj4MGDsjIRERGlzq9Lly7VjnfevHlo3Lgx7O3t4e/vj507d5Zbfvv27fD394e9vT2aNGmCL7/8str7NEXs33//PYKCguDu7g4nJyd07twZP/30k6zM4sWL9f4N5OTkmD3+bdu26Y3t6NGjsnL3+/wxR+z6/u9qNBo4ODhU+vyqa8eOHRgwYAB8fHyg0WiwePHi+74nJSUFPXr0gE6ng4+PD6ZPn46SPdJNdd1LMlY9ZQrG+pw3FmN8VpqKMT5rTMFYf6+mUtn4T58+rffab9y40TQB36Mi30v0sYTrX5XYLenaz507F+3atYOzszOcnZ0RGBiIDRs2lPseo153YWCffPKJmDlzpvjkk08EAHHq1Kn7vicxMVFotVoxY8YMcfjwYTFjxgyh1WrFrl27pDLvv/++cHR0FGvWrBEpKSnimWeeEfXq1RPXr183WOwvvviiqFevnti0aZNITk4WPXr0EA8++KC4c+dOme9JT08XFy9elB579+4VGo1GLF68WCoTExMj/Pz8ZOXS09MNFnd14q9IbJZ67YOCgsTChQtFSkqKOHDggHjyySeFp6enuHz5slQmPDxc9O7dW3Z+926vilWrVglra2sRFxcnDh8+LF5++WXh4OAgzpw5o7f8yZMnRa1atcTLL78sDh8+LOLi4oS1tbVYs2ZNlfdpqthfeeUVMWvWLPHnn3+KEydOiNjYWGFlZSV27NghlVm0aJGoVauW7BpfvHjRoHFXNf6tW7cKAOLQoUOy2O79f1WRzx9zxH7t2rVS17RJkyYiIiKiUudnCBs2bBBvvPGGWL16tdDpdGLRokXlls/MzBSenp7imWeeESkpKWL16tXC0dFRfPzxx1IZU113fYxVT5mCsT7njcEYn5WmYozPGlMxxt+rKVU2/lOnTgkAYuPGjbJrn5uba5qA71GR7yUlWcr1r0rslnTt161bJ3755Rdx4sQJcezYMfHmm28Ka2tr8ffff+stb+zrbvAEq9ju3bsrXHE9++yzonfv3rJ1vXr1EkOGDBFCCFFQUCC8vLzEjBkzpO03b94Ujo6O4ssvvzRIvNeuXRM2NjZi2bJl0rqzZ88KjUYjNm7cWOH9zJgxQ7i4uIibN29K62JiYkTr1q0NEmdZqhr//WJT0rXPysoSVlZW4qeffpLWhYeHi5CQEIPEWaxTp05i1KhRsnVNmzYVU6ZM0Vv+9ddfF02bNpWtGzlypOjSpUuV91lVhjjOQw89JCZMmCAtL1q0SDg4OBgsxvJUNv7iLz2XLl0qc5/3+/wxlOpe+4SEBAFA/PHHH9K6ipyfoTk4ONz3C8+8efOEk5OT7HPw3XffFd7e3qKgoEAIYbrrXh5D1lOmYKzPeWMxxmelqRjjs8YcDPX3ai4Vib/4S/7u3btNE1Ql6PteUpKlXv+KxG7J114IIVxdXcv8rmrs624RY7CSkpIQFBQkWxccHIzExEQAwKlTp5Camioro9Pp8Mgjj0hlqis5ORl5eXmyYzRo0AAtW7as8DGEEPj6668xbNgw6HQ62baTJ0/C29sbjRs3xpAhQ3Dy5EmDxG2I+MuLTSnXHgCysrJQUFAAV1dX2fqEhAR4eHigefPmiIyMRHp6epVjvX37NpKTk0v9fw0KCioz1rL+f+/Zswd5eXlV2qepYtcnKyur1DW+desWfH19Ub9+fYSGhmLfvn0Gifle1Yk/ICAA9erVQ69evbB161bZtvt9/hiCIa79ggUL0Lp1a3Tt2rXUtvLOzxySkpLQvXt32edgcHAwLly4IHXHM8V1NyRLiNdYn/PGYIzPSlMx1meNparI36sSDBw4EB4eHujWrRvWrFlj7nAAlP295F6Wev0rEnsxS7v2+fn5WLVqFbKzs/XWmYDxr7tFJFipqanw9PSUrfP09ERqaqq0vXhdWWUMEYNWq4Wbm1uVjxEfH49Tp04hMjJStr5z585YvHgxNm7ciAULFiA1NRVdu3bF5cuXDRJ7deK/X2xKufYAMG7cOLRv3x6BgYHSur59++Kbb77B5s2bMXv2bPz111947LHHkJubW6VYMzIykJ+fX6nrUdb/7zt37iAjI6NK+zRV7CXNnTsX586dQ1hYmLTOz88PCxcuxI8//oiVK1fC3t4e3bp1w4kTJwwWe1Xjr1evHubPn4+1a9fi+++/h5+fH3r16iUbS3G/zx9zxX6vzMxMfPfdd6U+WypyfuZQ1jUt3lZeGUNed0OyhHiN9TlvDMb4rDQVY33WWKqK/L1aMkdHR3z88cf47rvv8Msvv6BXr14YPHgwli1bZu7Q9H4vKclSr39FYre0a5+SkgJHR0fY2dnhxRdfxA8//IC2bdvqLWvs625dkUJTp07Fe++9V26ZrVu3omfPntUOyNAqGrshLFiwAA899BAefPBB2frHH39cttylSxc0adIES5YswYQJE8rdp7Hjr05s92PKaz9hwgQkJCQgISEBWq1WWj9kyBDpddu2beHv7w9fX19s2LABAwcONMixa4q1a9di0qRJ+Pbbb+Hr6yutDwwMlH0Ad+3aFe3bt8fnn3+O//u//zNHqBI/Pz/4+flJy4GBgTh9+jQ++ugjdO/e3YyRVc6yZctQUFAgS2wB9ZyfIbCeKpsxP+epEP8WzcfNzQ0TJ06UlgMCApCRkYEPP/wQw4YNM1tcZX0vUYKKxm5p197Pzw/79+9HZmYm1qxZg/DwcGzbtg1t2rQxeSwVSrDGjx9/3wvVsGHDKgfh5eWFtLQ02bq0tDRplqPi57S0NNlx7i1TlorGvmvXLuTn5yMjIwPu7u6yY1TkwzE9PR0//vgj5s6de9+yjo6OaN26dYV+4TdV/GXFpoRr/+qrr2LVqlXYunUrmjRpUm5Zb29v1K9fv8qtK25ubtBqteX+fy2prP/f1tbWcHNzgxCi0vs0VezF1qxZg+HDh+Obb75B//79yy2r1WoREBBg8Bas6sR/r86dO2PVqlXS8v0+fwyhurEvWLAAgwYNQp06de5btuT5mUNZ17R4W3llqnrdzV1PVYe5P+eNwRiflaZirM8aS1WRv1el6dy5MxYtWmS241fme4mlXf/KxK6POa+9ra0tmjZtCgDw9/fH7t278cknn+Drr78uVdbY171CXQTd3NzQokWLch+1atWqchCBgYGIj4+XrYuPj5f6TTZu3BheXl6yMjk5Odi5c2eZfSsrG7u/vz9sbGxkxzh37hyOHDly32MAhVNV29nZ4bnnnrtv2ZycHBw9ehT16tW7b1lTxV9WbJZ+7ceNG4eVK1diy5YtaNGixX3PLyMjA+fPn6/QtdfH1tYW/v7+5f5/Lams/98BAQGwsbGp0j5NFTsAfPfddwgLC8PixYvx9NNP3/c4QggcOHCgyte4LIa6Tvv375fFdr/PH0OoTux//fUX/v7771LdA8tS8vzMITAwEDt37pRN1R8fHw9vb280atRIKmPI627ueqo6zP05bwzG+Kw0FWN91liqivy9Ko05r31lv5dY0vWvbOz6WNL/+4KCgjKHhBj9uld7mowSLl68KPbt2yeWL18uAIgNGzaIffv2yaZ5fOyxx2Qz8fzxxx9Cq9WKWbNmiSNHjoiZM2cKa2vrUtO0Ozs7i7Vr14qUlBQxePBgo0wV7uPjI+Lj48XevXtFz549S01/6+fnJz7//HPZ+woKCkSzZs1KzThUbOLEiWLbtm3i5MmTYteuXSIkJEQ4OTmJ06dPGyz2qsZfkdgs9dqPHj1aODk5ic2bN8umB83KyhJCFM6AM3HiRJGYmChOnToltm7dKrp06SJ8fHyqFfuqVauEjY2NWLBggTh8+LB45ZVXhIODg3TNwsLCRFhYmFS+eOrhcePGicOHD4sFCxYIGxubUtO0l7dPQ6ls7CtXrhTW1tbi008/LXOq+9jYWLFx40bx77//in379okRI0YIa2tr8eeffxo09qrE/8knn4gffvhBHD9+XBw8eFBMmTJFABBr166VylTk88ccsRcbOXKkaNasmd59VuT8DCErK0vs27dP7Nu3T+h0OjFt2jSxb98+acrqKVOmiMcee0wqf+3aNeHp6SkGDx4sUlJSxNq1a4WTk5Ns+ltTXXd9jFVPmYKxPueNwRiflaZijM8aUzHG36slx7948WKxfPlycfjwYXH06FHx0UcfCRsbGzFnzhyTx36/7yX64reU61+V2C3p2k+ePFns2LFDnDp1Shw4cEBMmTJFaDQa8csvv+iN3djX3eAJVkxMjABQ6nHvNJu+vr4iPDxc9r7Vq1cLPz8/YWNjI1q0aFHqQ6mgoEDExMQILy8vYWdnJx555BGRkpJi0NhzcnLEyy+/LOrUqSN0Op0IDQ0VZ8+elZUBIGJiYmTrtmzZIgCU+YWyOCGxsbER3t7eYuDAgeLQoUMGjb2q8VckNku99vr+n91b5ubNmyIoKEi4u7sLGxsb0bBhQxEeHl5qv1Uxd+5c4evrK2xtbUXHjh3F9u3bpW09evQQPXr0kJXftm2b6NChg7C1tRWNGjUS8+fPr9Q+Dakysffo0UPvNb63zPjx40XDhg2Fra2tcHd3F0FBQSIxMdEosVc2/g8++EA0bdpU2NvbC1dXV/Hwww+LDRs2lNrn/T5/zBG7EEJcv35dODg4iA8++EDv/ip6ftVVPAV1yUfx53h4eLjw9fWVvefAgQOie/fuws7OTnh5eYnY2NhSU9+a6rqXZKx6yhSM9TlvLMb4rDQVY3zWmIKx/l5NpbLxL168WLRs2VLUqlVLODk5CX9/f7F06VKzxH6/7yVCWO71r0rslnTtw8PDZd9FevXqJbt9hamvu0YIM92qm4iIiIiISGUsYpp2IiIiIiIiNWCCRUREREREZCBMsIiIiIiIiAyECRYREREREZGBMMEiIiIiIiIyECZYREREREREBsIEixQnNjYWbdq0Mdj+Fi9eDEdHR4Ptr6Rp06bhhRdeMNr+LUl6ejrc3d1x7tw5c4dCREREZBZMsMhkBgwYgF69eundduTIEWg0GmzatMnEUQGDBw/GyZMnpWVDJnDp6emYPXs2pk6dapD93b59Gy4uLti/f79B9mdoHh4eGD58OGJiYswdChGRYkRERCA0NNTcYRCRgTDBIpMZOXIktm7ditOnT5fa9vXXX8PX1xe9e/c2eVw6nQ4eHh5G2ff//vc/dOrUCU2aNDHI/rZu3QpXV1e0b9/eIPszhhEjRmD58uW4cuWKuUMhIqqxduzYgQEDBsDHxwcajQaLFy/WW27evHlo3Lgx7O3t4e/vj507d1Zo/ykpKRg2bBi8vb1ha2uLRo0aYfLkybh165YBz4JImZhgkcmEhITA09MTixYtkq3Py8vD0qVL8cILL8DKygqHDx9GSEgInJyc4OHhgeeeew6pqall7regoADvvvsuGjRoADs7O7Rt2xY//vijrMyFCxfw/PPPo27duqhVqxbat2+PrVu3ApB3EVy8eDGmTZuGQ4cOQaPRSJXSCy+8UOrXxYKCAjRs2BBz5swpM7YVK1agf//+snU9e/bESy+9hIkTJ6JOnTpwd3fHZ599htzcXIwZMwa1a9dGw4YNsXTp0lL7+/HHH/HEE0/I4v7111/RokUL1KpVCwMGDEBmZibWrFmDZs2awcXFBWFhYbIKb8eOHejSpQscHR3h4uKCTp064eDBg9L2xMRE9OjRA7Vq1YKPjw9eeuklXL9+XdouhMDs2bPRrFkz2NnZoX79+njjjTek7W3atIG3tze+//77Mq8LEREZV3Z2Ntq0aYPPPvsMOp1Ob5lvv/0W48aNw5tvvol9+/aha9euePzxx3H27Nly97106VL4+/vD2dkZP/zwA44ePYpZs2Zh8eLFePLJJ41wNkQKI4hMaPLkyaJhw4YiPz9fWrd27VphZWUlzp49Ky5cuCDq1q0rXn/9dXH48GHx999/i9DQUNGpUyfpPTExMaJ169bS++fMmSOcnJzE8uXLxbFjx8Tbb78trKysxL59+4QQQmRnZ4umTZuKrl27ih07doh//vlHrF27VmzZskUIIcSiRYuEg4ODEEKImzdviokTJwo/Pz9x8eJFcfHiRXHz5k2RmJgotFqtuHDhgnTcjRs3ChsbG5Genq73XC9fviw0Go1ISEiQre/Ro4dwcnISMTEx4vjx4+Ljjz8WAETfvn3Fp59+Kk6cOCGmTp0qbG1tZccrKCgQPj4+YvPmzVLc1tbWolevXmLPnj0iMTFR1KtXT/Tq1UuEhoaKv//+W2zZskXUrl1bfPzxx0IIIfLy8kTt2rXFxIkTxT///COOHDkili9fLg4fPiyEEOLAgQPCwcFBfPzxx+L48eNi165dokuXLmLQoEFSHFOmTBEuLi7i66+/FidOnBCJiYli7ty5snMcPHiwGDZs2P3+OxARkRAiPDxchISECCGEyMnJEePGjRMeHh7Czs5OdO7cWezcuVNWPjs7W4SFhQkHBwfh4eEhZs6cKUJCQkR4eLje/Ts4OIhFixaVWt+pUycxatQo2bqmTZuKKVOmlBnrzp07hVarFV999VWpbWvWrBEASsVLVNMwwSKTOn78uAAgfvvtN2ldv379RN++fYUQQrz99tvisccek73nypUrAoD4888/hRClEyxvb28xbdo02Xt69Oghnn/+eSGEEHFxccLR0VFcunRJb0z3Jlj69l+sdevWYtasWdLys88+K0s8Stq3b58AIE6ePFkqti5dukjLBQUFws3NTfTv319ad/v2bWFjYyNWr14trfvrr7+Eq6uryMvLk+IGII4ePSqVmThxorCyspKd670V9+XLlwUAsW3bNr0xh4WFiRdeeEHveaSlpYmsrCxhZ2cn5s+fX+Z5CyHEq6++Kh5++OFyyxARUaF7P6dfeeUV4eXlJX7++Wdx+PBhMWrUKOHg4CD7wS06Olo0bNhQbNq0SRw8eFAMHjxYODs7VyrBys3NFVqtVnz33Xey9aNHjxaPPPJImbF27NhR9OrVS++24jrmiy++qMBZE6kXuwiSSTVr1gw9evTAwoULARR23fvtt98wcuRIAEBycjJ27NgBR0dH6dGgQQMAwL///ltqf9evX8eFCxfQrVs32fqHH34Yhw8fBgDs27cP7dq1g5ubW7Vij4yMlLo3XrlyBT/++KMUtz7F3fLs7e1LbWvXrp30WqPRwMPDA23btpXW2djYwNXVFenp6dK6H3/8ESEhIbC2tpbW2dnZwc/PT1r29PSEl5eX7Fw9PT2l/dSpUwcREREIDg5GSEgI5syZI+sKkpycjGXLlsmuf/G1/ffff3H48GHk5uaWOVlJMZ1Ox374RESVdOPGDcyfPx8ffPABQkJC0LJlS3z55Zfw9PTE3LlzARR2/Vu4cCE++OAD9OnTB61bt8bXX38NK6vKfaXLyMhAfn4+PD09Zes9PT3L7JafkpKCvXv3YsyYMXq3F3/u31tPEdVETLDI5EaOHIl169bhypUrWLx4MerUqSONKyooKEBISAj2798ve5w4caLSMyxpNBqDxh0WFoYzZ84gISEBy5cvh7u7O4KDg8ssX5zkXL16tdQ2GxubUrHqW1dQUCAtr1u3TrpOxUpWYhXZz6JFi/Dnn3/ikUcewU8//QQ/Pz/89ttvAAqv/6hRo2TX/u+//8aJEycqNbHGlStX4O7uXuHyRERU+ENWXl6e7EdDrVaLwMBA6UfD4jKdOnWSyjg4OBj09iVlKZ7B1t/fX+/2vXv3AgAefPBBo8dCZMmYYJHJPf3007C3t8eyZcuwcOFCDB8+XEoKOnbsiEOHDsHX1xdNmzaVPZycnErty9nZGd7e3vjjjz9k6xMSEtCqVSsAQIcOHXDgwAFkZGRUKD5bW1vk5+eXWl+nTh0MHDgQCxcuxMKFCxEeHl7uL4YPPPAAnJ2dpUqxOv7991/8888/6Nu3b7X3BRRWfpMnT8a2bdvQs2dPLFmyBMDd61/y2jdt2hQ6nQ4tW7aEnZ0dNm/eXO7+Dx48iI4dOxokViIiMvyPhm5ubtBqtUhLS5OtT0tLg5eXl9733L59GwDKnDRj7ty58PPzQ+fOnQEATz31FIYMGYKHHnoIDzzwAPbs2WPAMyCyXEywyOR0Oh2GDh2K2NhY/Pvvv7JudmPGjEFmZiYGDx6MP//8EydPnsTvv/+OqKgoZGVl6d3fpEmT8PHHH2PlypU4fvw43nnnHezcuROvvfYaAGDo0KHw8PDAE088gZ07d+LkyZP46aefpFkES2rUqBHOnDmDvXv3IiMjA7m5udK2yMhILF++HH///fd9bx5sZWWF3r17IyEhobKXqJQff/wRvXr1qvYNkU+dOoUpU6YgMTERZ86cwdatW3HgwAEpGZ08eTL++usvvPjii9i3bx/++ecf/Pzzz4iOjgYAODk5Ydy4cXjjjTewaNEi/Pvvv/jrr78wf/586Rg3b95EcnKywZJBIqKa4oEHHoCtra3sR8P8/HwkJSVJn9MPPPAAbGxssHv3bqnMzZs3ZbPBVoStrS38/f0RHx8vWx8fH4+uXbvqfU9xy9T27dtLbfv666+xZcsWzJ8/X0oGDxw4gI4dO2L37t2YPn06Zs+eXakYiZSKnWTJLEaNGoX58+eja9euaNmypbS+uDXqjTfeQN++fZGTk4OGDRsiKCgIdnZ2evf1yiuvICsrC6+//jrS0tLg5+eHtWvXShWBg4MDtm/fjokTJ6J///64ffs2/Pz88Mknn+jd36BBg/D999+jV69euHbtGhYtWoSIiAgAhVOs169fH76+vhW6t1VUVBQiIiIwe/ZsaLXaSl6lu9atW4dhw4ZV+f3FatWqhePHj+OZZ55BRkYGPD098fzzz2Py5MkACseG7dixA1OnTkWPHj2Qn5+PJk2a4KmnnpL2MWvWLLi6uuLdd9/FuXPn4OnpieHDh0vbf/zxRzRs2BDdu3evdrxERDWJg4MDXnrp/9u7n1fo2jiO458n0yxkQ1JTkzGslIXESYkyReFfIKXUPDQLPyb+ACXyIwtFmhoasrKxspBsLB1JaSRNFoRZGU3UNNe9uGvqvvF0N51mmud+v+oszulc5/p2dp+uc77Xv5qZmVF1dbX8fr9WV1f19PSksbExSVJFRYVGRkZy93g8Hs3NzSmbzf6yyvX29qbb21tJPz//vr+/18XFhaqqqlRbWytJmpyc1NDQkCzLUkdHhzY2NvTw8KBgMPhlfa2trerv71coFFImk5FlWUomk9re3lYkElEsFlN3d3du/vf3d01NTUmSGhsbv9x+BPhfKnaXDaCUpNNpU1lZaWKx2B+PaW9vNzs7O3nP+fLyYlwul3l8fMz7GYXU1tZmdnd3i10GAJSM79q0u93uL9u0p1IpMzg4aMrLy01NTY2Zn583gUDABIPB3D0nJydG0qfj906D6+vrxufzGbfbbVpaWszp6el/1ppOp83s7Kypq6szZWVlRpLp6ur61DH37OzM9PT05M63trbM9PR0Pq8HKDn/GGNMMQMeUAqy2aySyaTW1tYUjUaVSCQ+NZP4zuXlpWzb1vDwcF5z39zc6OjoSKFQKK/xhfT8/KxoNKpwOOz4/wIAgK99fHzI5/MpHA7nVowKZWJiQvv7+zo/P5fH48ld39zc1MrKiq6urvT6+qpAIKCDgwM1NDQUtD6gGAhYwB9IJBLy+/3yer2KRCLq7e0tdkkAgL+Ubdu6vr6WZVlKpVJaWFjQ4eGh4vG4vF5vQWvJZDJaXl5WU1OTBgYGctfHx8fldrt1fHwsY4wWFxfV19dX0NqAYiFgAQAAlBDbtjU6Oqp4PC6Xy6Xm5mYtLS192z69GDo7O7W3t5fbyxL4mxCwAAAA4Kj6+nrd3d0VuwygKAhYAAAAAOAQ9sECAAAAAIcQsAAAAADAIQQsAAAAAHAIAQsAAAAAHELAAgAAAACHELAAAAAAwCEELAAAAABwCAELAAAAABxCwAIAAAAAh/wA96IsTwDiBnEAAAAASUVORK5CYII=\n",
"text/plain": [
"