{ "metadata": { "name": "" }, "nbformat": 3, "nbformat_minor": 0, "worksheets": [ { "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "## Managing Sparse Data\n", "\n", "### Introduction\n", "\n", "The aim of this example is demonstrate how to take data sampled at arbitrary random locations and put it onto a defined grid. \n", "It also contains some specific suggestions for processing EN3 ocean profiles data. \n", "\n", "We take our example data from the EN3 datasheet of vertical ocean profiles\n", "observations. These data are produced by mobile submerging ocean buoys.\n", "[General information](http://www.metoffice.gov.uk/hadobs/en3/) on this data\n", "and an account of the [file format](http://www.metoffice.gov.uk/hadobs/en3/en3_file_formats.html)\n", "is provided by the Met Office.\n", "\n", "From our point of view, there are two key aspects of this data that we need to deal with:\n", "\n", " 1. each measurement value can also have a related 'quality control' (\"QC\") indicator that needs to be taken into account.\n", " 2. measurements occur at arbitrary times and locations, so the time, depth, latitude and longitude are provided by additional _measured_ values for each data point, instead of following a planned sampling pattern.\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Task summary \n", "\n", " 1. Read in potential temperature data from an EN3 observations file, taking the values from the POTM_CORRECTED (\"corrected pot. temp\") netCDF variable.\n", " 1. Also read the other information relating to these datapoints, including locations, depths and times.\n", " 1. Apply a quality threshold to the related quality control flag variable ('POTM_CORRECTED_QC'), \n", " such that only data values with a quality figure of less than 3 are included in\n", " the calculations.\n", " 1. Create a global 3d spatial grid with regular longitude, latitude and depth (down to 500m).\n", " 1. Collect the observations that fall within each grid cell.\n", " 1. Calculate mean, standard deviation and count (number) of datapoints within each gridcell: each observation point contained in a particular cell contributes to that cell's statistical result.\n", " 1. Implement a 'count' threshold on the data, so that each gridcell must contain >3 observation points to produce a result (otherwise it is set to 'missing data').\n", " 1. Preserve an indication of the original observation times in the gridded results.\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Import required modules, including Iris" ] }, { "cell_type": "code", "collapsed": false, "input": [ "import matplotlib.pyplot as plt\n", "import numpy as np\n", "\n", "import iris\n", "import iris.analysis\n", "import iris.quickplot as qplt" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 1 }, { "cell_type": "markdown", "metadata": {}, "source": [ "Check the Iris version." ] }, { "cell_type": "code", "collapsed": false, "input": [ "print('Iris version: {}'.format(iris.__version__))" ], "language": "python", "metadata": {}, "outputs": [ { "output_type": "stream", "stream": "stdout", "text": [ "Iris version: 1.6.0-dev\n" ] } ], "prompt_number": 2 }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Load data\n", "\n", "Define the data file." ] }, { "cell_type": "code", "collapsed": false, "input": [ "obs_filepath = './example_data/EN3_v2a_Profiles_195001.nc'" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 3 }, { "cell_type": "markdown", "metadata": {}, "source": [ "We first define some helper functions to load data from netCDf file variables and mask any points of poor\n", "quality.\n", "\n", "Each data variable has an associated variable containing its 'quality control' flags (\"QC\").\n", "On loading, we will check the related QC variable, where available, and mask individual datapoints according to a QC threshold value.\n", "\n", "In practical terms, the QC values are numeric but coded as string. Typical values are 1='good', 4='poor' and 0='missing QC data'. We also have to deal with possible missing datapoints in both the data\n", "variable and the QC variable: we mask any missing datapoints, but ignore\n", "missing QC data (as some QC variables in the file have _all_ values = 'missing').\n" ] }, { "cell_type": "code", "collapsed": false, "input": [ "def ensure_full_mask(array):\n", " \"\"\"Return MaskedArray version of data with a full mask array.\"\"\"\n", " full_mask_array = np.ma.getmaskarray(array)\n", " return np.ma.MaskedArray(array, mask=full_mask_array)\n", "\n", "def load_en3_variable(var_name):\n", " \"\"\"\n", " Load a named data variable from the EN3 file as an Iris cube with a masked\n", " data array.\n", "\n", " \"\"\"\n", " # Load cube from nc variable\n", " cube = iris.load_cube(obs_filepath, var_name)\n", " # Force data to be a masked array, with a fully expanded mask array\n", " cube.data = ensure_full_mask(cube.data)\n", " # Also infer a data mask if there is a '_fillvalue' attribute. This\n", " # _ought_ to be automatic with netCDF4-python, but some EN3 files seem to\n", " # have a mis-spelling here (should be '_FillValue', with capitals).\n", " if '_fillvalue' in cube.attributes:\n", " mdi = cube.attributes['_fillvalue']\n", " cube.data.mask |= cube.data == mdi\n", " cube.data.set_fill_value(mdi)\n", " return cube\n", "\n", "def load_en3_with_quality_mask(main_var_name, qc_var_name, qc_max_valid=9):\n", " \"\"\"\n", " Load a named data variable from the EN3 file as an Iris cube.\n", "\n", " Generates a masked data cube, in which data is also masked where the\n", " related QC value exceeds a given threshold.\n", "\n", " Args:\n", "\n", " * main_var_name (string):\n", " the data variable long_name within the netCDF file.\n", " * qc_var_name (string):\n", " name of the related QC variable\n", "\n", " Kwargs:\n", " * qc_max_valid (int):\n", " Threshold value for QC data. QC values larger than this result in the\n", " datapoint being masked.\n", "\n", " \"\"\"\n", " # Load main and QC data values\n", " cube = load_en3_variable(main_var_name)\n", " qc_data = load_en3_variable(qc_var_name).data\n", " # Missing QC data equates to okay\n", " qc_data.set_fill_value(qc_max_valid)\n", " qc_data = qc_data.filled()\n", " # Convert to numbers\n", " qc_data = np.array(qc_data, dtype=int)\n", " # data QC figure means invalid if *too large*\n", " cube.data.mask |= qc_data > qc_max_valid\n", " return cube" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 4 }, { "cell_type": "markdown", "metadata": {}, "source": [ "**NOTE:** This approach is still somewhat simplified, as the original files\n", "contain extra QC data that may need to be considered (see the file format specification, linked above). We are ignoring that here, for simplicity. " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Load the required data into cube variables." ] }, { "cell_type": "code", "collapsed": false, "input": [ "# Get the main data (potential temperatures), applying quality levels\n", "potm = load_en3_with_quality_mask('corrected pot. temp',\n", " 'quality on pot. temperature',\n", " qc_max_valid=2)\n", "\n", "# Get depth and location info\n", "depth = load_en3_with_quality_mask('corrected depth', 'quality on depth')\n", "longitude = load_en3_with_quality_mask(\n", " 'longitude of the station, best estimated value',\n", " 'quality on position (latitude and longitude)',\n", " qc_max_valid=1)\n", "latitude = load_en3_with_quality_mask(\n", " 'latitude of the station, best estimated value',\n", " 'quality on position (latitude and longitude)',\n", " qc_max_valid=1)\n", "\n", "# Get time reference + convert to units string\n", "reftime = iris.load_cube(obs_filepath, 'date of reference for julian days')\n", "reftime_str = ''.join(reftime.data)\n", "assert len(reftime_str) == 14\n", "ref_unit_str = 'days since {:4s}-{:2s}-{:2s} {:2s}:{:2s}:{:2s}'.format(\n", " reftime_str[:4], reftime_str[4:6], reftime_str[6:8],\n", " reftime_str[8:10], reftime_str[10:12], reftime_str[12:14])\n", "\n", "# Get time data\n", "time = load_en3_with_quality_mask(\n", " 'julian day (utc) of the location relative to reference_date_time',\n", " 'quality on date and time')\n", "time.units = ref_unit_str\n", "\n", "# Show some results\n", "print '\\npotm cube:\\n', potm\n", "print '\\ndepth cube:\\n', depth\n", "print '\\nlongitude cube:\\n', longitude\n", "print '\\ntime cube:\\n', time" ], "language": "python", "metadata": {}, "outputs": [ { "output_type": "stream", "stream": "stdout", "text": [ "\n", "potm cube:\n", "corrected pot. temp / (degree_celsius) (-- : 2581; -- : 55)\n", " Attributes:\n", " _fillvalue: 99999.0\n", " c_format: %9.3f\n", " comment: corrected value\n", " fortran_format: f9.3\n", " resolution: 0.001\n", "\n", "depth cube:\n", "corrected depth / (metre) (-- : 2581; -- : 55)\n", " Attributes:\n", " _fillvalue: 99999.0\n", " c_format: %7.1f\n", " fortran_format: f7.1\n", " resolution: 0.1\n", "\n", "longitude cube:\n", "longitude of the station, best estimated value / (degrees) (-- : 2581)\n", " Attributes:\n", " _fillvalue: 99999.0\n", "\n", "time cube:\n", "julian day (utc) of the location relative to reference_date_time / (days since 1950-01-01 00:00:00) (-- : 2581)\n", " Attributes:\n", " _fillvalue: 99999.0\n", " conventions: relative julian days with decimal part (as parts of day)\n" ] } ], "prompt_number": 5 }, { "cell_type": "markdown", "metadata": {}, "source": [ "**NOTE:**\n", "The data here is actually dimensioned by `[<\"profile id\">, <\"depth sample number\">]`. It turns out that these dimensions actually have no practical meaning for our purposes, except that position and time information only depend on the first of these. These dimensions will be lost when data is put onto the grid.\n", "\n", "**NOTE:**\n", "You would not normally see a '\\_fillvalue' attribute. This is due to a mis-spelling in the netCDF source datafiles -- see above." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Perform the regridding operation\n", "We want to place our data onto a regular grid.\n", "\n", "To do this, for each datapoint we work out which gridcell it belongs in, by comparing its measured latitude, longitude and depth values to the gridcell boundaries. \n", "We then calculate a result for each gridcell, which is a statistical combination of all the datapoints that fall within that cell. \n", "\n", "Generally there is less data than fills the whole grid (in this case, for instance, all gridcells over land will be empty). Thus, many cells will have no result, which is why we call this a 'sparse' arrangement.\n", "\n", "#### Define the parameters of the desired global grid" ] }, { "cell_type": "code", "collapsed": false, "input": [ "z_min, z_max, dz = 0.0, 500.0, 50.0\n", "y_min, y_max, dy = -90.0, 90.0, 2.0\n", "x_min, x_max, dx = -180.0, 180.0, 3.0\n", "nz = int((z_max - z_min) / dz)\n", "ny = int((y_max - y_min) / dy)\n", "nx = int((x_max - x_min) / dx)" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 6 }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Define a function for the main regridding operation\n", "\n", "Because of the sparse nature of the source data, we can use the \"cube.aggregated_by\" method. \n", "This collects data into categories based on attached categorical coordinate values (see Iris documentation for this method, under [Iris.cube.Cube](http://scitools.org.uk/iris/docs/latest/iris/iris/cube.html)).\n", "\n", "The basis of this operation is explained more fully in the \"Coordinate Categorisation\" example (see: \n", "[Coordinate Categorisation](http://nbviewer.ipython.org/urls/raw.github.com/SciTools/iris_example_code/master/coord_categorisation.ipynb)).\n", "\n", "In outline, we will do the following:\n", "\n", " * calculate gridcell indices for each data point from the associated values of depth, longitude and latitude\n", " * attach these index-value arrays to the data as auxiliary coordinates\n", " * categorise the data over the index coordinates, performing an aggregation (mean, std-dev or count) for each result \n", " * expand the results, which are still sparse, into a full cube on the required grid" ] }, { "cell_type": "code", "collapsed": false, "input": [ "def aggregate_to_grid(phenom, aggregator=iris.analysis.MEAN, **agg_kwargs):\n", " # Regrid a data cube onto the grid previously specified (by 'nx', 'dx',\n", " # 'x_min', 'x_max' and the equivalents in y and z).\n", " #\n", " # For this, each datapoint is assigned to a gridcell 'box' according to its\n", " # corresponding coordinate values (depth, longitude, latitude).\n", " #\n", " # N.B. agg_kwargs is passed to the aggregator operation.\n", "\n", " # Make a 1-D 'flattened' version of the source data cube\n", " # (as 1D coordinates are required by the 'aggregate_by' operation)\n", " cube_flat = iris.cube.Cube(phenom.data.flat[...])\n", " cube_flat.metadata = phenom.metadata\n", "\n", " # Make coordinates containing gridcell indices for each grid coord value\n", " def add_index_coord(coord_points, coord_name, start, step):\n", " # Calculate gridcell indices of all the points.\n", " cell_indices = np.floor((coord_points - start) / step)\n", " # Make these all integers, for eventual use as array indices.\n", " cell_indices = np.array(cell_indices, dtype=int)\n", " # Add these as a \"categorical coord\" to aggregate by.\n", " cube_flat.add_aux_coord(\n", " iris.coords.AuxCoord(cell_indices, long_name=coord_name),\n", " 0)\n", "\n", " # Add a coordinate containing gridcell indices in the z-dimension (depth).\n", " add_index_coord(depth.data.flat[:], 'i_depth', z_min, dz)\n", " # Repeat for lats and lons -- except these need broadcasting to 2d first\n", " lons_2d, _ = np.broadcast_arrays(longitude.data[:, None], phenom.data)\n", " lats_2d, _ = np.broadcast_arrays(latitude.data[:, None], phenom.data)\n", " add_index_coord(lons_2d.flat[:], 'i_lon', x_min, dx)\n", " add_index_coord(lats_2d.flat[:], 'i_lat', y_min, dy)\n", "\n", " # Aggregate the data to get a statistical result for each 'inhabited' cell.\n", " result_cells = cube_flat.aggregated_by(('i_depth', 'i_lat', 'i_lon'),\n", " aggregator=aggregator,\n", " **agg_kwargs)\n", "\n", " # Make a full-grid result cube.\n", " # N.B. metadata comes from aggregation (includes units + cell_method)\n", " full_data_empty = np.ma.MaskedArray(np.zeros((nz, ny, nx),\n", " dtype=result_cells.data.dtype),\n", " mask=True)\n", " result_cube = iris.cube.Cube(full_data_empty)\n", " result_cube.metadata = result_cells.metadata\n", "\n", " # Assign aggregation results to the appropriate gridcells ...\n", " i_z, i_y, i_x = [result_cells.coord(coord_name).points\n", " for coord_name in ('i_depth', 'i_lat', 'i_lon')]\n", " # ... but discarding results that lie outside the target grid.\n", " i_ok = np.where((i_z >= 0) & (i_z < nz) &\n", " (i_y >= 0) & (i_y < ny) &\n", " (i_x >= 0) & (i_x < nx))\n", " i_z, i_y, i_x = i_z[i_ok], i_y[i_ok], i_x[i_ok]\n", " result_cube.data[[i_z, i_y, i_x]] = result_cells.data[i_ok]\n", "\n", " # Add DimCoords defining the grid.\n", " result_cube.add_dim_coord(iris.coords.DimCoord.from_regular(\n", " z_min - 0.5 * dz, dz, nz, with_bounds=True,\n", " standard_name='depth', units='metres'), 0)\n", " result_cube.add_dim_coord(iris.coords.DimCoord.from_regular(\n", " y_min - 0.5 * dy, dy, ny, with_bounds=True,\n", " standard_name='latitude', units='degrees'), 1)\n", " result_cube.add_dim_coord(iris.coords.DimCoord.from_regular(\n", " x_min - 0.5 * dx, dx, nx, with_bounds=True,\n", " standard_name='longitude', units='degrees'), 2)\n", "\n", " return result_cube" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 7 }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Make and output the results\n", "\n", "First regrid the data to form the main required result cubes." ] }, { "cell_type": "code", "collapsed": false, "input": [ "# Regrid the main data to get a gridded mean, std-dev and count.\n", "potm_mean = aggregate_to_grid(potm)\n", "potm_std_dev = aggregate_to_grid(potm, iris.analysis.STD_DEV)\n", "potm_counts = aggregate_to_grid(potm,\n", " iris.analysis.COUNT,\n", " function=lambda points: np.ones(points.shape))\n", " # NOTE: 'function' arg looks a bit odd, as must return an *array* of bool" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 8 }, { "cell_type": "markdown", "metadata": {}, "source": [ "Additionally, mask data according to the \"`count > 3`\" requirement (task point \\#7)." ] }, { "cell_type": "code", "collapsed": false, "input": [ "# Add occupancy 'count > 3' as an extra condition for valid data on the others\n", "invalid_counts = potm_counts.data <= 3\n", "potm_mean.data.mask |= invalid_counts\n", "potm_std_dev.data.mask |= invalid_counts" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 9 }, { "cell_type": "markdown", "metadata": {}, "source": [ "Also add a bounded 'time' coord from the minimum and maximum times of the values contributing to each cell (task point \\#8).\n", "\n", "**NOTE:** In principle, the 'aggregate\\_by' method should be able to _automatically_ calculate these time bounds, if the time values were attached to the data as an auxiliary coordinate. However, at present it does not calculate these correctly." ] }, { "cell_type": "code", "collapsed": false, "input": [ "# Expand the times to 2d so we can use the same code to get min+max times.\n", "time_data2d, _ = np.broadcast_arrays(time.data[:, None], potm.data)\n", "time2d = iris.cube.Cube(time_data2d)\n", "time2d.metadata = time.metadata\n", "\n", "# Add a bounded time coordinate to all the output cubes.\n", "time_min = aggregate_to_grid(time2d, iris.analysis.MIN)\n", "time_max = aggregate_to_grid(time2d, iris.analysis.MAX)\n", "time_centres = 0.5 * (time_min.data + time_max.data)\n", "time_bounds = np.concatenate((time_min.data[..., None],\n", " time_max.data[..., None]),\n", " axis=-1)\n", "time_coord = iris.coords.AuxCoord(time_centres,\n", " bounds=time_bounds,\n", " units=time.units,\n", " standard_name='time')\n", "potm_mean.add_aux_coord(time_coord, (0, 1, 2))\n", "potm_std_dev.add_aux_coord(time_coord, (0, 1, 2))\n", "potm_counts.add_aux_coord(time_coord, (0, 1, 2))" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 10 }, { "cell_type": "markdown", "metadata": {}, "source": [ "Save results to a file." ] }, { "cell_type": "code", "collapsed": false, "input": [ "# Save result cubes to netCDF\n", "iris.save((potm_mean, potm_std_dev, potm_counts), 'temp_sparse_regrid.nc')" ], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 11 }, { "cell_type": "markdown", "metadata": {}, "source": [ "Plot, and print, a 2d slice of the cube to give an idea of the results." ] }, { "cell_type": "code", "collapsed": false, "input": [ "# Plot the mean result at the top level\n", "plt.figure(figsize=(12, 8))\n", "qplt.pcolormesh(potm_mean[0])\n", "plt.gca().coastlines()\n", "plt.show()\n", "\n", "print potm_mean[0]" ], "language": "python", "metadata": {}, "outputs": [ { "metadata": {}, "output_type": "display_data", "png": "iVBORw0KGgoAAAANSUhEUgAAAqwAAAHVCAYAAAA0B6mkAAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAIABJREFUeJzs3Xd8VFX6+PFPei+T3iuQEAIh1BAIRaRKIm0VxIAgiqss\nKuAqrBRRXAUBFVcXWQGB3aUE6dVIDUgngAkhhVTSM+l9Zs7vD365X7MUQcFEPe/X676STLn3OZNk\n5rnnnnMePSGEQJIkSZIkSZJaKf2WDkCSJEmSJEmS7kUmrJIkSZIkSVKrJhNWSZIkSZIkqVWTCask\nSZIkSZLUqsmEVZIkSZIkSWrVZMIqSZIkSZIktWoyYZUkSfoFFi5cSHR0dEuHIUmS9LsmE1ZJkn4V\n//nPf+jWrRvW1taoVCoef/xxjhw50tJh4ePjw+HDh3/28/X09B5iNP9n3bp1REREPNBz9PX1uXHj\nxiOJR5IkqSXJhFWSpEdu+fLlzJ49m/fff5/y8nKKi4t5/fXXOXDgwAPvS6PR3HabVqv92bHp6enx\ne6qf8ntqiyRJUhOZsEqS9EiVl5ezYMEC1qxZw+DBg9HT08PAwIAnnniCDz/8EIC6ujpeeOEF7Ozs\nsLe3Z9q0adTX1wNw9OhRPDw8WLJkCe7u7kyZMoV33nmHsWPHEh0djUql4uuvv6akpITx48djZ2eH\ng4MDs2bNQqfTKXGsWLECX19frKysCAwM5OLFi0RHR5OVlUVkZCRWVlZ89NFHABw+fJjQ0FCsra0J\nDAxsllhfv36d7t27Y21tzeDBgykuLr5r25ti//vf/46zszMuLi589dVXyv1qtZqxY8diY2ODs7Mz\nf/vb3xBCcO3aNf785z/z/fffY2VlhZ2d3U++zn379gUgJCQEKysrtm7dCsDmzZsJDAzE2tqaLl26\ncO7cOeU5Pj4+fPTRR3Tu3BkrKyuef/55CgoKGDZsGFZWVvTp0we1Wg1ARkYG+vr6rF69Gk9PT+zs\n7Hjvvfd+Mi5JkqSHQkiSJD1C+/fvFxYWFvd8zKxZs0Tfvn1FWVmZKCsrE/379xezZs0SQghx5MgR\nYWhoKBYsWCC0Wq2oq6sTCxYsEKampmL//v1CCCHq6urEoEGDxCuvvCLq6+uFWq0WPXv2FCtWrBBC\nCLF27Vrh7e0tfvjhByGEEBkZGSIrK0sIIYSPj4/47rvvlFhSU1OFra2tiI2NFUIIcfToUWFjYyNy\nc3OFEEKEhISIOXPmCJ1OJ86ePStsbGxEdHT0HdvVFPvcuXOFTqcTZ86cEZaWliI+Pl4IIcSYMWPE\n008/LWpra0Vubq4ICgoSK1euFEIIsW7dOtGnT58Heq319PREWlqa8vOJEyeEo6OjuHz5shBCiH//\n+9/C1dVV1NXVKW3v3bu3UKvV4ubNm8LFxUWEhoaKxMREUV9fLwYNGiTmzp0rhBAiPT1d6OnpicmT\nJ4uGhgaRkpIiXFxcxK5dux4oRkmSpJ9D9rBKkvRIlZSU/GQP4aZNm5g/fz42NjbY2Ngwf/58/v3v\nfyv3GxkZ8fbbb6Ovr4+JiQkAffr0YejQoQDk5+dz/PhxPvroI4yNjVGpVLz66qts2bIFgDVr1jBn\nzhw6dOgAgLe3N56enneMZePGjURGRjJw4EAA+vXrR1hYGLt37yY5OZmkpCQWLFiAnp4e3bt3Z9So\nUfe8DG9gYMD8+fPR09OjR48ejBw5kq1bt1JbW8uuXbtYvHgxpqamuLq6Mnv2bKXd99rn/frqq694\n6aWX6NSpEwDPPPMM1tbWHD9+XHnMK6+8gkqlws3NjYiICHr16kX79u0xNjZm5MiRXL58udk+582b\nh5GREW3atGHq1Kls3rz5F8cpSZL0UwxbOgBJkn7f7O3tlcvKd1NQUICXl5fys6enJ4WFhc32YWjY\n/O3KxcVF+T4nJ4fGxkZcXV2V23Q6HR4eHsCthNbPz+++4s3JyWHr1q3s3r1buU2j0dC/f3+Kioqw\ns7NTkmYADw8PMjIy7rq/Oz2+sLAQtVqNRqO5rd0FBQX3Fef9tmXLli2sXLlSua2xsbHZMAZnZ2fl\nexMTk2Y/GxsbK0Mzfhx/E3d3d06fPv3Q4pUkSbob2cMqSdIj1atXL4B7TrBydnYmMzNT+Tk7Oxsn\nJ6e7Pv5/Z+a7uLhgaWmJWq2mtLSU0tJSysvLSUhIAMDNze2us+f/d1+urq5MmTJF2U9paSmVlZW8\n9dZbODo6olarqauraxbrvdzp8c7Oztjb22NgYHBbu5sS8Yex+oCrqysLFy5s1paqqirGjx9/1+f8\nVM9uTk5Os+9/fOIgSZL0qMiEVZKkR8rGxoZFixYxdepUvv32W3Q6HY2Njezfv58333wTgKeffpr3\n3nuPsrIyysvLeffdd3nmmWfuus//Tar8/f3p3r07c+fOpbq6GoDMzExOnjwJwOTJk/nwww9JTEwE\nbk0gako07ezsSE9PV/YVHR3N9u3bOXLkCEIIGhsbOXnyJLm5ubRr146AgADee+89dDod58+fZ+fO\nnfdMLrVarfL4M2fOsGvXLsaOHYupqSlRUVHMmzePuro68vLyWL58uZJM2tnZkZeXR2Nj432/1v/b\nlqlTp/LFF19w6dIl4NbktkOHDlFVVXXf+/xfixcvpqGhgdTUVNasWcNTTz31s/clSZJ0v2TCKknS\nIzdz5kyWLFnC3LlzsbW1xcnJiRUrVvDEE08At5KgNm3a4Ofnh6+vL/7+/rz//vvK8/83IdTT07vt\ntq1bt5Kbm4u3tzfW1tZERkaSlZUFwKRJk5g+fboy+3348OGUlJQA8MYbbzBv3jxsbW1Zvnw5bdu2\n5b///S9z587FxsYGFxcXJeGEW7PuDx48iK2tLXPnzv3JogEuLi6Ym5vj5uZGVFQUy5cvJyQkBIAv\nv/yShoYGnJ2dCQkJYcSIEUyfPh2Axx9/HD8/P+zt7ZXe5vfff5/hw4ff9Vhvv/02Tz/9NCqVipiY\nGPr27cvSpUuZNGkSVlZWeHt7s2rVqnsm2D++706vc1hYmHKC8NJLLxEZGXnP9kuSJD0MeuJhjOyX\nJEmSbnP06FGio6N/ctjAb0FGRgZ+fn5oNBr09WVfhyRJvy75riNJkiRJkiS1ajJhlSRJeoQeVenW\nlvB7aoskSb8tckiAJEmSJEmS1Kq1yDqsdnZ2lJaWtsShJUmSJEmSpJ9BpVL95Lraj0qL9LDq6ek9\nlCoukiRJkiRJ0q+jJfM3OYZVkiRJkiRJatVkwipJkiRJkiS1ajJhlSRJkiRJklo1mbBKkiRJkiRJ\nrZpMWCVJkiRJkqRWTSaskiRJkiRJUqsmE1ZJkiRJkiSpVZMJqyRJkiRJktSqyYRVkiRJkiRJatVa\npDSrJEm/H1qtlnPnzpGenk59ff1tm06nw9DQ8Cc3IyOjZj/b2dnh6+uLo6Mjenp6Ld1MSZIkqQXJ\n0qySdAfV1dUkJCSQlZVFbm4uZWVllJWVUVlZiUajQavVKl8BfH190Wg0lJeXK1tlZSVGRkaYmJhg\namra7GvT9z++3dTUFCsrK2xsbLCxscHOzg5/f3+sra1/9fZXVVWRkJDAlStXuHr1KomJiTQ0NGBg\nYIChoSEGBgYYGBgghODcuXO4uLgQFBSktO3Hm76+vvJ63e/W2NhIcXExGRkZ1NbW4uPjg4+PD76+\nvrd9ValUMqGVJEn6FbRk/iYTVukPR6vVUlJSwtWrVzly5AinT58GwNLSks6dO9OmTRvmzp2Lvb09\nPj4+uLu7o1KpsLW1xcrKSukB/HHSlpaWhrGxsZJs2tjYYGVlhUajob6+nrq6Ourq6pRex7v9XFlZ\nqSS8xcXFpKWlUVtbS9u2bblw4QKWlpa/qO0NDQ1cuXKFa9eu4eHhQWBgIC4uLkrCt2HDBt58803K\nyspo3749HTt2pFOnTnTo0AFTU1O0Wq2yaTQadDodoaGheHp6PnAsQggaGxuVBLVp+/HPTScBN27c\nIDU1lbS0NFJTU0lOTkatVgNgbW3NqVOn6NChwy96bSRJkqR7kwmrJD1kGo2G1157jZKSErKysqiu\nrqauro6SkhLKysqwsbEhICCAAQMG0Lt3b4yMjCgvL+fs2bNcuHCBWbNmMWzYsF8l1qqqKq5fv05R\nURGFhYUUFRWRk5PD1atXuXz5MsXFxXh5eTFp0iRyc3PJy8ujoqKCyspKZauqqkKr1SKEUDYALy8v\nhg4dir6+PufOnePq1av4+/vToUMHcnNzSUpKor6+noCAABwdHdm3bx8ABgYGTJs2jalTpxIaGvrQ\n26yvr/+L3wPc3d1ZtGgR3bt3Jzg4GD09PbRaLWq1muLiYoqKiiguLqa4uBi1Wo2ZmZly4qFSqVCp\nVLRv3x4DA4OH1CpJkqTfN5mwStIvcPLkScrKypSeOT09PZycnOjbty9mZmZMnTqV6OhorKyssLOz\nw97evtUkKXdK3F5//XXMzc0xNTWlurqaDz74gO7duzNkyBA8PT1xdXXF2toaKysrZbO0tMTIyAi4\n9f/V1GN6/fp1Dh48CECPHj0IDQ29rZe2pKSE5ORkSkpKKCgoYNeuXezatUu5v6CgACcnp4fa7srK\nSo4fPw6gDIdo2oyMjCgsLCQ7O1vZsrKysLa2JjAwEIC3334bgEmTJpGfn09mZiZFRUWUlZVha2uL\ng4MDDg4OODo64uDggJ2dHbW1tZSWllJWVkZeXh4XLlzA0NCQkpKSFhl2IUm/hvr6enJzc7l58yZq\ntZouXbrg4eHR0mFJv1EyYZWkB6TVavnmm29Ys2YNBw4cACAyMhJjY2Nqamo4ffo0gwYNwszMjCtX\nrlBSUsJTTz1Fly5dCA0NxdHR8bZL0BqNBgMDAyUJtLCwQF//4S6kodVqOX78OJcvX+bq1ausWbMG\nuDUcwdvbGzs7O1JTU6mvryckJAQXFxcaGhp49913ad++/UON5bdq+vTp/OMf/+CZZ56he/futG3b\nFh8fH5ycnFCpVBgaNp9LmpOTo/QuN22ZmZn4+vqir69PbGwsLi4uLdSaR6euro6CggLKy8vx8fHB\n2tqaxsZGzp49i1arxdXVFTc3NywsLFo6VABKS0v5+uuvlZ7x6upqvLy8AGjfvj0JCQn88MMPGBoa\nYmlpSb9+/XjhhRdaOOrWJy8vj++++47Y2FiOHj1KXl4eLi4uuLu7Y2Njw/nz51GpVAwcOJDRo0fz\n+OOPyzHg0n2TCaskPYDc3FzGjBmDRqPhjTfeoG3btkqy9+PHrF27lvj4eDIzM0lNTaWiogJTU1N0\nOh21tbX3dazDhw8zYMCAXxxzQ0MDGzZs4IMPPsDa2ppevXphb2/Pe++9h06nIzw8nIiICLp160a3\nbt3w9vaWHyIPQK1WExYWRkpKinKbi4sLnp6e3Lhxg169etGxY0c6duxIcHAwAQEBGBsbt2DED66h\noYHs7Gw8PDwwMTEBoLy8nOvXr9PQ0EBGRgbR0dEA2NjYUFNTg5OTE9bW1mRmZmJnZ0d1dTU+Pj6Y\nm5uTl5dHXl4ehoaGeHt707FjRzp06ECXLl1+teEwP3by5En69OnT7Lbg4GAqKyvJzMy87fF9+vTh\nxIkTv1Z4rV5hYSHOzs7Kz0uWLCEyMpJ27do1O/HW6XRcvXqV2NhYvv76axoaGnj55ZeZNGkSNjY2\nLRG69BsiE1ZJegDJyckMHDiQ/v37M27cOLp37w7cmsRjaWl51x6j2tpasrKyOHToEAsXLkStVmNi\nYsKgQYMICwtDX1+fvLw86uvrCQ8PJzw8nDZt2twzcRRCUF9fT21trbJVVFRw9OhRkpOTlTGUSUlJ\ndOrUiblz59K3b1/09PTYu3cvI0aMwNPTk8LCQurr63FycuLmzZu39RJK91ZUVERQUBDFxcXKbSNH\njmTmzJl069YNMzOzFozu58nMzOSdd97hxo0bpKenk5+fj4ODA1ZWVjz22GNoNBo2btxI+/btMTMz\no7a2litXrqDRaACwtbXF2dmZWbNmYWtrS2lpKc7OzgQEBKCvr09NTQ3V1dXk5uaSkJDAF198QWFh\nIXZ2dpSUlLRo2w8fPszKlSspKCjAwcGBwMBAxowZg7u7O87OzsrwF+n/6HQ6Nm/eTHFxMd999x2F\nhYV8880397x6IIQgLi6Ozz77jEOHDhETE8PAgQMRQqDT6VrN0Cmp9ZAJq/S7oNPpSE9PR6VSNevt\nvJO1a9cyZcoUfHx8cHV1xdHRET8/P/z8/PD398ff3582bdoob5j5+fnk5+djamqKSqXCyMiIr776\nit27d5OYmKh8wLZv357ExMR7HjsvL4/p06fzzTffKLeNGTOGmJiYB27z/yazvr6+mJub06tXL0JD\nQ5WxlF5eXnh4eJCZmUl6ejrp6elKInLjxg3S0tKoqKhg8uTJyjAB6cGUlZUxYcIE1Go1q1evJjg4\nuKVD+kVyc3MJCwujd+/eLF68GC8vLwwNDTl48CBpaWno6+vTs2fPZpPiysrK+O677/juu++UYRCW\nlpbY2NhgbW2NoaEh5eXlAJibm2Nubo6FhQXm5uZ06NCB6OhoOfTkd0Cn0zF37lw+//xz3N3dCQ8P\nJygoCDc3N9zd3XF3d8fV1bXZUnzjxo1j8eLFfPXVVxw7doyqqiocHBzo0qULvXr1olevXvTs2VOO\n9/6Dkwmr9Js3Y8YMVq5cedvtNTU1d+zd2rNnD5GRkQCYmZnx1FNPUVJSgqWlJWq1muTkZCwtLVm9\nejVhYWFERkayZ88e2rVrR0FBAVZWVkRERPDGG28QGhrKggUL+Mc//sG8efMICQnB1tYWW1tbzM3N\nqauro6amhoyMDFavXs3u3btpbGzEyMiILl260KFDB+bMmUObNm0euN2NjY2MHz+ebdu2sWnTJp5+\n+mngVi/w+vXrSU9PJyMjg/T0dEpKSvD09MTX1xdfX1/8/PyafbW3t5fDAH6BZ599FmNjY1atWvW7\n6YFLS0sjKioKY2NjZsyYwfjx4zE1Nf3J5+3evZuPPvqIH374gYaGBry9vQkJCWHw4MGUlJQwYcKE\nZpePpd8njUbD1atXOXXqFKmpqdy8eVOZgJWXl4eRkZFyMmNpaUlmZiavvPIK06dPx8bGhry8POLi\n4ti3bx979+5FrVazf/9+hg4d2tJNk1qITFil34SmdTd/PPavsLCQt956i0OHDlFSUqJcRqqtrWX0\n6NHExMTcMQnTaDT4+/uTlZXF4MGDMTExoaSkhCtXruDk5ERAQACxsbE0NjZy/vx5ioqKiI6OZsSI\nEXh5ebFp0yaSk5N57rnn6N27N66urgghiImJIT09nbKyMsrLy6mursbMzAxzc3McHBzw9/dn/fr1\nzWJ5//33ee211+77srEQgvXr17N3714uX75Mbm4uUVFRzJkzR+nVGzt2LNu2bSMoKIgnn3ySqKgo\nQkNDlbGH0sM3atQooqOjGT16dEuH8lDpdDoOHTrExx9/zKVLlxg6dCghISF0796diIiIZo9Vq9Uc\nP36cF154gVWrVtGzZ0+ysrIIDw8HbiX1Wq2WAwcOoFKpsLKyIj4+viWaJbUib731Fh9++CFeXl54\ne3ujVqspLS1FrVaj0+mwtbXFxsaGoKAg/vOf/2Bubt7SIUstRCasUqtRVVXFgQMH2L59OxcuXEBP\nTw99fX3KysrIz89n2LBh7NmzR3n85s2bGTduHAD+/v5UVVVRUFAAwL///W+eeeaZBzq+VqslJSWF\nlJQUSktLuXHjBs899xw+Pj5kZmaye/dusrOzWbduHYWFhc2eO3LkSLZv337XfTc0NDBnzhw2bNjA\na6+9RlJSEhs2bFDu/6m/yYMHDxIbG0tcXBw1NTW88cYbhISEEBgYeFuPXmZmJkePHuXKlStcvnyZ\nc+fOUVFRcdd9x8TEMGbMmHseX7q3mTNncvz4cZYvX07fvn1bOpxHIjk5mb/+9a/s3LkTuJXMNp0Q\nZmdn07lzZzp37sz06dMZNWoUAJs2bWL8+PEAhISEoNFoSEhIUPb5431If0w7d+6kpKREWZ84NTWV\nc+fOkZaWRmJiorLWc5PNmzfz1FNPtWDEUkuRCavUooQQHDlyhM8++4zY2Fh69erFqFGj6NOnD3p6\nesqEoFWrVuHi4sInn3yCkZGR8iEnhCAjI4PExES8vLxo06YNOTk5bNiwAY1GQ9u2bZk0adIDLREl\nhCAlJYWrV69iZGREbm4uQUFB9OnThx49enDhwgXatWuHTqejvLwcPT09Bg8ejK2trbKPDh06MG7c\nOCwsLAgKCiItLQ0vLy8CAgLYsWOH0qOq0+morq7GyspKOfadPsDt7OwoLS1Vfj5//jxdu3b9ybbs\n2rWLJ598Erg1c71jx45YWVmh1Wqpra2lrq6OTz/9lJCQkPt+faRbv7dt27Zx9uxZLl26RGpqqjKb\nfO/evQwfPryFI3w0rl+/zmeffcb+/fupqanh8ccfx8LCgnPnzjF8+HAWLVp023OEEMoJoL6+PkFB\nQfc1tOCPSAhBTk4OZ86c4fTp02RmZlJaWkppaSnm5uYsWrTooawc0hoVFBQwbNgw6urqePbZZ+nU\nqRPBwcF4eXnx3HPPNTvBd3R05LPPPqOoqEh5D27q4NDT08PPz4+ePXvKNV9/Z2TCKrWYprGkBgYG\nzJ8/n759+1JTU8OlS5c4efIk33//Pfr6+tja2tKpUydKSko4e/YsBgYGyuUjLy8v/P39GT58OMHB\nwSQkJNCxY0cGDx5MRESE0nv44osvkpaWhomJCW3atLltCZX169ezcuVKjI2NuX79ulIqVavVUlZW\nRkpKCuPHj8fW1pZ9+/YRHx+PhYUFNjY2FBYWUlNTw4svvkhwcDBCCNatW8elS5coLCxk2LBhXLx4\nEQATExOsra1xdHREpVJRW1tLZWWlUj2qpqYGOzs73N3dlUkKjY2NHD9+nMzMTPr06UO/fv2U2dc/\nRQjB6dOn2bhxI1u2bMHPz48hQ4bQvn37ZmVZGxoa0Gg0zTYDAwPMzMwwMzPDwsKCjh07EhoaKpMN\nbo2PjoqK4tixY5iZmTFt2jReffVVTE1N/zDjgVNTUzly5AgajQYrKyvGjh0r/zbuQafTkZaWxsWL\nF7l48SIpKSmUl5dTUVFBRUWF8r2VlRU9e/bEx8cHuFVco6SkhHPnzuHh4UF8fHyr//tKTk7m8uXL\nzVYw0dfXp0+fPlRUVBAXF8eJEyf4/vvvqaqqUjoUZs+ezd///vd7tk+r1fLvf/+bRYsW0a9fP5yc\nnJpV2dNqtVy/fp0zZ85gYmJCWFgY0dHRREZGPvS1rR8mtVrNlStXSE9PR61WY2BggKGhIUZGRhgb\nG+Pt7U3btm1xd3dv1e14lGTCKrWYffv28eSTT+Lh4UF1dTWurq64uLgQHBxM7969CQ8Px8rKipkz\nZxIbG8uSJUsYMmQIOp2OrKwsMjMzycrKIikpSem1HD16NKdOnSIuLo633nqLadOmMWTIELKysggO\nDqaqqgobGxtOnjzZLJZFixbx9ddf4+zsjLOzM++9916z+vAXLlzg6aefJjU1lfLycuLi4jhz5gz5\n+fkMGDCAtLQ0Pv74Y3bu3Env3r2Jj49nwoQJ3Lx5Ew8PD4YMGYK5uTnGxsbK6gOlpaUUFhZSWFhI\ndXV1syWqDA0N0dfXx9ramm7dujFlyhQ6der0i5Z6aWho4OTJkxw8eJD09PRmFZ6MjY0xMjLC0NAQ\nQ0NDDAwMmvXCVlRUcPnyZa5du0aHDh0ICwtj5MiRPPbYYz87nt+DiooKYmNjefXVV9m5cyddunRp\n6ZBalJ5ewR1vF+KPN8mqadLQqVOnuHDhAvHx8djb29OlSxe6dOlCQEAAhoaG5Ofnk5eXR25uLrm5\nueTk5JCWloa5uTlt2rRRVi3p1KkTI0aMaDWT+kpLSzl27BjHjx+noKBAuXRfWVlJeno6/fr1w9zc\nXDnpraur4/jx49ja2tKnTx8iIiIIDw/Hzs4OnU6HVqt9qOPshRCkp6dz4sQJ3n33XZ5++mkWL178\n0Pb/MOXk5NCjRw/atGlDmzZtsLKyQk9PTykqU1dXR0ZGBsnJyVRUVODv70+7du1o164dr776Kvr6\n+sTFxSkrwzg7O6NSqVr9ic2Dkgmr1Kr96U9/IiYmhhEjRuDt7U1DQwP19fXU19dTWFjIyZMnaWho\nAGD58uWUlJSQlpbGhQsXKC4uZsWKFejp6VFVVUV1dTXHjx+npqaG7777rtlxEhIS+O9//4uVlRXb\ntm3Dw8ODMWPGIISgrKyMmJgYBgwYwKhRoxgyZAgeHh5ERERgZmbGqlWrGDlyJL6+vsTHx/Ppp5/i\n7OyslD5t+tBqWjGgqUfFxMQEf39//Pz8cHBwUJb4aaqYVV1dzZUrV/jXv/5Ft27d+Oqrr1riV9BM\nTU0NFy9eZPv27WzatInPP/9cuSz8R5SQkMBzzz3HtWvX+Pjjj5k6dWpLh9SiZMIKc+fOZcuWLZSW\nltK7d2969+5Nt27dCA0NxcbGhjNnzrBz5052795NRkYGbdu2JTAwUNnatWt3x6tALUkIwZ49e1iz\nZg21tbUUFBSQmppKeHg4/fv3x8PDo1m5Zm9v7xZdCaKpR/rMmTMcPnyYS5cucejQIcLCwlospnt5\n+eWX2bFjB/PmzWP06NE4OzvT2NhIRkYGxcXF6OvrK1t1dTVpaWls3ryZCxcucPbsWSZOnEhcXFyz\nfRoZGeHm5kZmZiYODg7cuHFDGXr2WyUTVulX09jYyJEjR5Qz8erqaqqqqqiqqsLW1hYfHx/8/Pww\nNzfHyMgIIyMjqqqqyM/PJysri/LyckxMTDA2NqasrIzXX3+92f779evHjRs3UKvVtGvXDnt7e2Ux\nf0tLSywtLXF3d+fll1++55n8hQsX2LFjBykpKRgYGGBra4urqyuzZ8/m6NGjvPPOO1y5cgU7Ozus\nrKyUmfpffPEFo0eP5tKlS1RUVODp6Ym3tzceHh44ODigUqkwNzdXFsbWaDRUVVUpCWxFRQW1tbXY\n29vj5OREUVERKSkpJCUlMWnSJD755JNH/Su6b/X19fzzn/9k165dxMfH85///IchQ4a0dFi/qpyc\nHHJzc/kbO0RMAAAgAElEQVTLX/7C9evX2bp1K4MGDWrpsFrUHzlhbWhoYP/+/UybNo0pU6bw3nvv\nKSet+/fvZ9u2bezZswdnZ2eioqJ48skn6dKlS6tfIL+6upqRI0eSk5PDnDlzcHJywtbWli5durSa\nim0NDQ1s27aN3bt3c/bsWYqKiujWrRvdu3enX79+9O3bt1WfVNfW1nLgwAG2bt3Kvn37sLW1JS8v\nD3d3d2XIg06na7ZZWFjw1Vdf0b59e9LS0oiNjeXs2bOcPXuWnJwcysrKlP1PmDCBdevW/eaLwsiE\nVfrVJCYm0rNnT2praxk7dix2dnZYWlpiampKZWWlsph905jKxsZG6uvr0dfXZ/jw4URHRzNgwADl\nQ2Dp0qUUFBRgaGiIjY0NwcHBBAcH4+Pjc8cxPjqdjtLSUurq6jAzM7trgYHKykpiYmIwNDSkY8eO\nGBsbo1arUavV9OzZE2dnZxITE1m0aBFxcXHk5eVhYmKijDXy8PDA09MTOzs7hgwZQm1tLSUlJcrk\nidraWmpqaqipqSEwMBBPT0+sra2Veus//PAD9fX1dOjQAVNTUxobG5Xk29LSEnt7e7y8vFrN5Z6T\nJ08yZswY3n//faZMmdLS4fxqfvz69+7dm379+rXaS46/lj9iwpqbm8uwYcO4cuUKERERTJw4kXHj\nxilri06bNo3c3Fyef/55IiMj8fPza+mQFWVlZSQkJJCVlUV2drayRnTT8KCGhgYKCwtZunQpbdu2\n5ezZs62q57fJZ599xl/+8hcAHBwcGDNmDAEBAdy4cQMXFxfefPNNcnNzycvLo6CgAGNj42Y9wtbW\n1tja2t715KGgoIDExEQcHR1xdHTE3t7+kSV/dXV1ZGdn4+XlRUNDAwcOHKCgoICKigpCQ0PvO/lu\nqoT4expXLhNW6Vf12WefMXPmTM6dO3dfM9N1Oh07duxg1qxZZGRk4O/vz+TJkwkNDWXIkCE/2TtR\nWlrKqVOnyMrK4uWXXwZuzTCtq6ujc+fOHDlypNk+Pv74YxYuXMiAAQMwMTHhhx9+QKvVYm5uzsWL\nFwkODubbb78lJCSEv/zlL4wdOxZfX1+0Wq0y1mjPnj3MmDGDwMBAdu7cSUVFBZ999hmbNm2irq4O\nd3d3SktLqays5ODBgxw7doyjR4+SmJhITU0N7du3Jz09nQkTJrB9+3YqKioICwtTeqbz8vKwsbFh\n9OjRjBs3rlXM8L9+/TrDhg0jOjqahQsXtppk+mFKTk5m1apVwK21fG/evKkMH8nOzm7h6FoHPa87\n3y6yft04fi06nY6VK1fy2muvAf+3PJ1Wq+WLL75g4cKFzJo1i9mzZz/UsacVFRVYWlr+4sk3Tz75\nJMnJyXTq1AkPDw9MTEzQaDQ0NjbS2NiIsbGxsg6qo6MjTz31VKvtpWtsbCQ3N5fU1FQSEhJITk7G\n19eXtWvXkp+fj4mJCW5ubsrl9h9Pdi0vLyckJIQVK1YoY25/PP52zpw5xMTEYGZmRlFREWq1Wpnk\nZWZmxtChQxk7diyRkZEP7bJ7U/nsJvb29kpVxTfeeIMlS5Y8lOP8lrRo/iZaQAsdVhJCfP7558LT\n01OcPn36no+rrKwUmzdvFmPGjBHm5ubC1NRUGBoaCh8fHzFx4kQxdOhQAYiePXuK8PBw4ebmJoYM\nGSI2btwoDh8+LKKiokRgYKCYP3++ePHFF4VKpRITJ04Ubdu2FQYGBsLMzEzY2dmJFStWKMc8deqU\nGDhwoHBzcxOAcHJyahbTzJkzBSBGjBgh/vOf/whPT09RV1en3H/s2DExevRoERgYKBwdHcWqVauE\nRqMR8fHxwtHRUSxcuFBkZGQInU4nhBDi9ddfF0FBQWLo0KFi8ODB4uDBgyI7O1u5v6ioSDz//PPC\n1tZWuLi4iF69eono6GjRq1cv4eHhIQwMDAQgAFFZWfmwfkW/SH5+vujRo4d44oknRG5ubkuH81Bd\nvXpVfPPNN6J3794CEFFRUWLTpk0iLi5O1NTUtHR4rQaed95+j86fPy/8/f1Ft27dxIYNG5T3g3Pn\nzomuXbuKvn37isTExId6zH379gk/Pz9hZmYmjI2Nhaenp+jRo4cYOXKkOHLkiPK4V199VURFRYn8\n/Px77i86OlqsWbPmocbY2mRlZYlTp07d8zEajUaMHz9edO7cWQQEBAgvLy/h6OgoLCwshL6+vgDE\n5s2blcdrtVqh1WqFTqcTxcXFYs2aNWL48OGiTZs2oqKi4hfHXFlZKT755BPRrl070aZNGxEYGCiC\ngoJE586dRYcOHcSUKVPEgQMHxNmzZ0VaWpooLS0VWq32Fx+3tWvJ/E32sP6BrF27lvnz53P8+HF8\nfX3v+BitVsuaNWuYP3++Mkh8woQJjB8/nu7du2NpaQnA6dOnWbp0KSNGjMDf3x8vLy++//57NmzY\nQHZ2Nq+88gqdO3dm+fLlpKenk52dzebNm+nXrx9qtZpFixaxbds2rl+/rlRNiY2NVcYfent788QT\nT/CPf/xDiU2j0XD27FnWr1/P0aNHad++PcbGxgwaNIhr166xadMm3n33Xbp160ZAQAAmJiY0Njbi\n7+/P0qVLlbKpTYQQLFu2jJs3b7JkyZK79r68+OKLXL58GXNzczIyMigoKKB///6Ym5tz5coVsrKy\nCA0NJSoqinHjxt31tf21HDt2jMGDBxMTE6OUv/2t+3GBiiZt27YlOTm5hSJqvf4IPaxarZaYmBj+\n9a9/4eHhwZo1a5QrCgkJCfTu3ZtPPvmEiRMnPtQrDTk5OXTt2pX169czePBg6uvrlRUGvvnmG86f\nP8+RI0eoq6ujf//+FBcX4+bmpgxdOn36NGVlZcrEVUtLSwoKCkhLS+PGjRst/t7RWgkh0Gg0P9lD\nXlxcTGhoKG+++SbTp0//RccMDQ1FCMH06dNxdXWlpKQEtVqtfP3x901fq6urUalU2NnZYWdnh729\nPW5ubvztb3/D29v7F8XTWvwhhwS8/fbbWFtbExoaSpcuXe46lrElNL0JeXh4tPrB+PeroaEBFxcX\nbGxsGDZsmDLe78eLOmu1WmbOnMnp06eJiIhgxYoVfP7554wYMQILCwu0Wi1GRkZYW1vf93E1Gg3D\nhw8nLi6OzZs3N0ugnnnmGXx8fHj//feV2yoqKjh48CC7d+9m3759+Pv7M3/+fHr06IFGo+HkyZNs\n3rwZnU5HQEAAy5YtQ6vVYmhoyIABA5g9ezaPPfaY8iGl0+lwd3fn+PHjtG3b9me9djExMVy6dIl2\n7drRtm1bZeH/JlVVVZw+fZpt27YRExPDuHHjGDlyJJmZmRQXF/Pss8/i5ub2s479oE6fPk2/fv2Y\nOHEiPXr0QKvVMnbsWBwcHH6V4z8Kqamp2Nvbc+HCBa5evcqWLVs4ffo0f/rTn9i8efPvcujDL/FH\nSFgLCgpwcXFhxIgRzJkzh169eil/B4WFhbRr144LFy7g7+//UI5XWlrKp59+ymeffcacOXOYOXNm\ns/tPnDjB2LFjeeWVVygoKGDz5s107twZa2trtm/fjpubG+PHjyc8PBxHR0dl4mpxcTG7d+/m008/\nBX662p50b3v37mXUqFGsWLGCP//5z8pwjcTERIYPH46npydDhw6lW7duZGZmkpSURGFhIf7+/rRv\n3x5vb29iY2PJzc0lISGB9PR0Ll++fN/5SWNjo1LStimJPX/+PF988QUzZswgMjKSTp06/abfs/6Q\nCWv37t3JyMigqKgIuFVWsCUqYlRUVHDmzBmqq6vR09Nj586d7NixA1NTU9RqNd7e3kydOpWePXtS\nUFCAmZkZVlZW9O7du9WOI7qbhoYGLl68yKlTpzh58iTHjh1j2rRpFBUVcfnyZRISEnB0dGT9+vWU\nlpayY8cO1q5de9t+HvRPprKykqysrGZrqgJkZGTQtWtX1q5dy8CBA28bxK7RaNi5cycffPABGRkZ\nCCHo1asXffr0Yf/+/VhZWbFy5Uq8vb0pKipix44dLFu2jJ49e/Lll18qA91nzpzJihUrADAwMECj\n0TxQ/A8iOzubIUOG4OzsjLe3N0ZGRmzfvp1Zs2YRERFBUFDQIz05Ky0tZdWqVZw4cYLKykrS0tLI\nzc0lLCyMN998k+HDh7eaWcX3q6nC2KJFi5g+fToqlYqLFy8yYcIEpZdfenC5ubmsXLmSqqoqZQ3g\nkpIScnJyyMnJwcnJiT59+tCnTx/CwsJa3XI8c+fO5ejRoyQnJ6PT6di4cSM9evTA3NycZcuW8ckn\nnxAeHs4rr7zC4MGDf3aSsH79el577TWefPJJ5s6de8cT3zlz5rB69WqmTZvGkiVLmDVrFu+//z4v\nvvgia9euRaVSUVxcfNdjCCFQq9XY29v/rBil/3P16lVefPFFysrKGDp0KGFhYbzxxhssWLAADw8P\nDhw4QHx8PH5+fgQGBuLk5ERqaipJSUnKMmFBQUHKONm33nrrF3c4NC2NuHfvXurr6xk+fDiPP/64\nsgrP/261tbXY2toqE8zCwsJo06bNQ3qFfp78/HyCgoIoLS394yWsTzzxBJaWlnTo0IEBAwbQp0+f\nR35cIQQnT57k1KlTpKenc/r0aVJSUujatatS8WjIkCE89dRTeHh4UFtbS2JiIm+//TY5OTkEBARQ\nW1vLvn37cHFxYcCAAdja2lJQUICenp6yiHC7du1wdHTE0tISKysrLCwslDfLmpoazMzMHskZlk6n\no7CwEGdnZ4QQzJs3DwMDA5599lnatWvX7HEASUlJrFq1irZt2xISEkKnTp2azT5t6rmEW7+zzp07\nY2RkREpKilKitKnih4uLC3v37sXR0fGBYt6wYQNr167l/Pnz9O7dm0GDBjFw4EA6dux422QGrVar\n9Hi7urpy4cKF295IampqeO6558jOzmb79u24uLig0+no27cvJ0+eZOLEiURERODi4oK3tzd+fn4P\ntNSKRqNBT0/vgXrer127xrJly/jhhx9ITEzE3Nycjh07MmXKFP70pz89khOf+Ph4Bg4cSHBwMAUF\nBVy/fl25LzU19aH1PD1KOp2O+fPns3jxYuzs7AgMDOTUqVP06NGDiooKbt68yZgxY+54UiX9HyEE\naWlp5OXlKQUqtmzZwj//+U8mTZqEn58f9fX1NDQ0oFKp8PT0xMPDg5s3bxIXF8exY8c4deoUcOvk\nb9myZS3cotvt3r2b6dOnU11dTU1NDXV1dc0+VE+cOPGzP2Oee+45zM3N+fzzz+/5uOXLl3Pq1Cne\neecdZXjQ4sWLqauro7i4WJYo/RXpdDrOnj3L4cOHOXr0KKNHj+all15q6bAQQpCcnMzevXs5duwY\n+vr6ysozP95MTU0pKyujqKiIgoICYmNjeeaZZ1p0WcW6ujp69epFfHz8Hy9h/bUPm5yczKRJkygu\nLiYyMhIfHx969Ojxs9axKyws5Nq1a2RnZ1NWVoazszNarZaUlBSuX79OcnIyarVaqTiir6+Pl5cX\n+fn5VFVVYWFhQZcuXejcuTP+/v74+vri6+uLt7f3PZe/EP+/xrWRkREODg7KDPf9+/crtdQNDQ3p\n2bMnVlZWpKSkMGzYMNauXcuoUaOYO3cuhw8f5t1332Xy5MnMmzePQ4cOcfjwYSwtLfH29ubxxx9v\nlgRqNBr2799PVFRUs1gMDQ1v66nMysrC09PzgV7LJuXl5Xz77bfExsaye/duGhsbiY+Px83NjaKi\nIpYtW8aKFSvo0KEDgwYNYvXq1QwcOJDo6GhGjBjRLLndv38/Y8aMYebMmbz77rvo6emh1Wo5cOAA\ncXFxFBQUkJ+fT2ZmJjdu3MDQ0BBra2tlaRUHBwdcXV0JDAxEpVLR2NiIWq2muLiYrVu3oq+vz4QJ\nE3B1dSUqKkop33g/hBDcvHmTM2fO8PHHH6PVatm/f/9DX6Zm3rx5NDQ08OGHHyq35eXl8fHHH6NW\nq1m6dOl9lZVtSenp6cryQ4GBgWg0GlJTU5X7e/bsyYkTJ1pN1aHW5MqVK2zcuJHz589z8eJFrK2t\n8fT0VCq59ezZk3feeee+/l/j4uKIiIgAbo3ntrCwoE2bNjz99NOttkdQp9NRV1dHdXU1NjY2t73H\nHz16lMOHD+Pq6qqUX26avf6/J6OZmZl06tSJxMRE3N3d73rMU6dO8eqrr/LJJ58wcuRITp061eK9\nYtLvw+XLlxk8eDDJycnY2NjQ0NBAeno6VVVVaDQatFotKpUKNzc3rK2tlQ6xpvua1jyvqalhzJgx\nnD17FicnJ5ydnZWvdnZ2yjKWdXV1WFtbExQURFBQEO3bt1fmmvwhhwT8WodtaGjg4MGDvPPOOwwb\nNoyFCxf+6uNSS0tLyc7OxtXVFQcHBwoLC7lw4QKXL18mPT1d2bKzs3FxcWHw4MFERUXx2GOPYWFh\nwZUrV1i5ciWHDh2itrYWuFXz2NTUlP79+zN8+HDCwsJo164dpqamfPXVVwghlLGLarWaJUuW8Omn\nnxIWFkZ4eDiffvopenp6qFQqpk6dSl1dHdevX+e7777D3d2dbt26KVVfAgIC+OCDD/j6669ZsGAB\nQgjMzc3p2rUr3bt3v69kq66uTknyKyoqMDc3JzAwkDZt2igfJl27duXixYsYGBgQEBDA+fPnWb58\nOQsWLCAkJISNGzdy9uxZkpKSWLNmDb6+vtTV1WFsbKy0DWDWrFnExsaSk5OjVNRqKiX7v0vCaDQa\nTExMmDZtGiUlJRQXFytntdXV1fj5+dG1a1fc3Nyws7MjIiICQ0NDdu3aRW5uLt988w2vvPLKz1pG\nSgjBs88+q9SDb3pDeBj27NmjLLfl7+/P9evXSUxMpGvXrpw6dQq1Wq1MoGvNNBoN33//PXv27OGb\nb75Bo9EwZswYxowZQ8+ePf+w9bzvRAjB4cOHWbp0KVeuXGHq1KmEh4fTpUsXnJycftG+L168yOzZ\ns1Gr1YwaNYqkpCT27dvHgAEDCAoKavbh5+npSdu2bVvtOL2UlBTatWtHhw4dCA8Pb1aStbKykm7d\nuhEeHs6wYcMIDw8nISGBkJAQunbtyuDBg1m8ePEd2/bVV19x4sQJAGxtbfnrX/96xwRYkh7UzZs3\niY6O5ty5c9jY2FBUVKSsHd5Uxru0tJSbN28qn8+VlZU0NDRgZmbG888/z5w5c8jKyiIiIoI1a9YQ\nGhpKQUGBUpq8pKQEIyMjTE1NMTExobS0lMTERBITE0lOTsbV1ZWwsDD++9///vESVjs7Ox5//HFm\nz55N9+7dH9q+z58/z+TJkyktLVWqGIWGhvLss88yZcqUVt0bo9VqSU1NZd++fezcuZPvv/8eExMT\nLCwsmDFjBqNHj6Zdu3bo6ekpdZ8fpD0NDQ0YGxuj0Wj44Ycf8PLyum08pVar5cKFC1y5coWkpCSu\nXbtGUlIS2dnZqFQqNBoNGo0GW1tbzMzMKCkpYcmSJVhaWipr6pWVlVFeXk55eTmlpaVcvXqVzMxM\nVCoVFhYWNDY2Ym5uTmNjI3l5ecoHXUNDA9euXSM5OVnptbxx4wbr1q1j69atpKenY21tTU1NDW3b\ntiUpKYk9e/Zw9OhRdu/eTXx8/G1tjoyMJDMzk7KyMry8vFCr1Xz++ef0798fuDXrt3379lRWVt72\n3ISEBKZOncpTTz11W0WvJk0fXFOmTKFHjx6EhYXd16D6efPmsXXrVnJzc5Ue44c9trWmpoaTJ0+S\nlZWFjY0NVVVV/Pe//8XR0ZGNGzc+1GM9LBqNhpycHI4dO8a+ffv49ttv8fX1Zfjw4YwaNYrQ0NBW\nmwi1lMbGRmJiYli6dCl1dXXMnj2bCRMmPNSa8HdSXl6ulDZt+tArLCwkNTUVAwMDRowYQWRkJP36\n9btrLEIIvv32W2UogpOT032fhOh0Oq5evcqJEyc4ceIEMTExuLi4kJ2dTUFBgVIDvry8nCNHjnDw\n4EHMzMxYuHAhx48f59VXX2XNmjVMnjy5WZvOnDnDiRMn2LVrF/n5+URFRREQEICXlxfz5s3Dzs6O\nf/7zn7etvbxmzRreeecd/vSnP3HgwAGSkpL44IMPmD179s9/kSXpRxobG5X5Pne7MtxUrdHKygoz\nMzMKCgr44IMPWL9+PY8//jiHDx/G3t6epKSk+34v1Wg0pKWlKfNe/nAJ6xtvvMG2bdtIT0/HxsaG\nefPm8eyzz/7inoC2bdtSVlbGE088wbhx4wgLC2v1lz7vRghBeXm5Utu+JWk0GnJzcyktLaWxsZHX\nX3+duLg4IiMj0Wq1mJmZKZfVraysqK2tRaVSAbeSM39/fyZOnIhKpcLV1ZXz58+zbt06HnvsMdq3\nb09ubi6WlpYMGzaMxx577Lbj/+1vf+Nf//oXISEhxMfHo9FomDhxIs7OzkpvYlOVlR979913aWho\n4K233sLc3JzJkyejUqlYsWIFtbW1PPbYYwwePJh33nnntud+9NFHLF++nL179xIaGnrb/UII5s6d\ny8GDB6moqMDZ2Zm8vDyEEIwePZo333zzrn/PGzZsYMuWLRw7dozAwEAGDx5Mx44d8ff3p127dg+0\nEsO9CCGYMWMGGzduJCIigscee4yXXnqpRSuv5OXlsWnTJnJycpRx15mZmWRlZZGfn4+TkxNhYWEM\nHz6cYcOG4erq2mKxtmbl5eX861//4tNPP8XX15fZs2czfPjwFu91FkKQkJDA7t272bNnDwkJCQwc\nOJDIyEhGjRrV7IqMEEKJ197enqqqKtzc3AgPD+fLL79sdtWhoaGB8+fPKwnqyZMncXJyIiIigrS0\nNI4ePYqRkRFdunQhJSUFCwsLDA0NMTMzIyIigiFDhpCTk8P777+PVqtlwYIFTJs27Z5jyG/cuMGO\nHTvYvn07V69eZcCAARw6dIiNGzcyatSo2x4fGxvLiBEjWLlyJdOmTePatWsEBAQ8xFdXkn6emzdv\nsm/fPrp3706nTp1+9vvEH7pwQHV1tVi3bp2YOHGisLGxESNHjhQ7d+4UDQ0NP2vfpaWlYt++fWLO\nnDnC3d1dDBw4UCQlJT2s0CUhxIcffigcHBwEIHx9fUXHjh2Fn5+fcHZ2FhYWFsLAwEAEBgYKJycn\n8dxzz4ng4GBhZGQkLly40Gw/FRUV4r333hMvvfSSeOutt8STTz4pVCqVaN++vejevbtwc3MTBgYG\nwsDAQOjr64u5c+cKIYTQ6XTi008/FYAYN26cmD9/vqiurr6v2KOiokRUVJQ4duyY6NSpk5g0aZJS\nKOB/bdiwQTg4ONwWd5P169eLzp07i3379omNGzcKQGRnZ4vLly+L6dOni8DAwJ+Mq76+Xhw5ckT8\n7W9/E2PGjBGdO3cWtra24sMPPxSNjY331aZ7WbFihQgODhbl5eW/eF+/VG1trfjzn/8sbG1txeTJ\nk8XSpUvFRx99JNavXy+OHTsm0tPTf/b/fWvDqjtvD0tsbKxwdnYW48ePF+fOnXt4O34ECgsLxbp1\n68TIkSOFq6ur2L17d7P7jxw5IiwsLERRUZE4e/as6N27t9DX1xfPP/+8EEKIjz/+WPTv319YWlqK\n0NBQMWPGDLF169ZmC/JnZGSItWvXig8++EBs3rxZaDSau8ZTXV39swp9FBQUiDVr1oizZ8/e8f76\n+nrx2muvia5du4qCggIREBAgnn766Qc+jiS1Zi2UNt46dosc9C4NLi8vF6tXrxbdu3cXgBg1apSY\nNGmS+Mtf/iLefvttsWzZMpGRkXHfx6mvrxcrV64UDg4OYurUqWL37t2ivr7+YTXjD6+6ulpcu3ZN\nxMfHi5SUFJGXlycqKiqUD4vS0lIxZ84cMXv2bJGdnX1f+9RoNOLChQvi+++/F1lZWaK+vl40NjaK\n2traZlVEjh49KqZNm3bXZPNuampqxIwZM4SdnZ1Yu3btXZ9fWloqVq9eLby8vMSBAwfu+JiYmBjh\n4OAggoODha2trZgxY4bS9jNnzghLS0tx+PDhB4pPCCHS0tLEoEGDROfOncXSpUvF559/LrZs2fLA\nVVSKioqEvr6+WLp06QPH8LClpKSIzp07i6eeekqUlZW1dDiP3KNKWHU6nVi+fLlwcXH5WX9bD5NW\nqxUZGRmisLDwvp+zY8cO0bFjx2a3LVmyRPj4+IjLly8rlfi+/PJLcfDgQbF27Vrh4uIi9u3b16r/\nbpKTk0W3bt1EZGSkSElJEZ06dRIvvPDCT1a5kqTfmpZMWFt00tWPlyn6sUOHDjFkyBDat2/PG2+8\nQUVFhbKMzdatW4mKiuKtt96646WWuro6Ll68yOnTpzl9+rQysWnQoEHk5eWhr6/P1q1bZUURCYAt\nW7awevVqhg4dikqlIicnh0uXLnH48GEGDRrEs88+y5NPPnnHsT5CCK5evQrcWmaraUmvL774ggUL\nFvDPf/6T0aNH3/XYeXl5HDhwQBnH+9xzzyn/D0IINm/ezLlz50hPT+fQoUMUFRVhZmb2QO37+OOP\nWbVqFdeuXXug5z0MarWanTt3sm3bNuLi4li8eDEvv/zyLxqDmp6eztGjR8nLy6OsrIzRo0crk+2E\nEJw5cwatVkvXrl3RaDTs3buXb775hry8PCwsLNDX12fXrl2PfCKM3pd3vl28+PP3WVtbywsvvEBi\nYiLbt29v8co5Tk5Oyjrajo6OBAcH06FDB+Vrhw4dUKlU1NfXk5ycTFxcHMuWLcPT05MBAwbwww8/\noFKpGDJkCGvWrMHMzIzPP/+cuXPnsu3/sXfmcTWlbwD/3tteKJTKmomyjGzJEkplbcTYMskYk6xj\nrA3ZzWDsyxB+jLFFSFH2ncqefamoRLRIaN87vz8ad6QblbKe7+fTp3vf8y7POZ3ufc7zPouXF3p6\netSpU4effvrprf9HHxsfHx8cHR0xNzcnNjaW0NBQevXqxaZNm0R/a5Evjq8yS4CJiQk3b97kxo0b\nmJiY5FNec3NzOXbsGPHx8Tg4OOQb+/TpU0aMGMGePXtwd3fHwcGBjIwMNmzYwP79+/H398fY2JhW\nrVd/9UUAACAASURBVFrRqlUrmjZtyrBhwzh//jxKSkooKirSr18//vnnnw992iKfGBcvXqRr1670\n6NGDnJwcpFIp1atXx8jICDs7uxL5Pp85cwYHBwf8/f1lKZkKo3///jx79gxTU1MuXLhAUlISa9eu\nzReEePnyZQYPHkzPnj2ZM2dOseWZMGECZ8+exdvbu8yqbaWmphISEsLt27e5ceMGt27dIjY2lgcP\nHmBjY0OfPn2wtbV9b7/ca9eu0axVZ9CxAnUDUFCBSHdQ1maiU3u8vb1RVlZGQ0ODoKAgFBQUMDc3\np3fv3tSpU4eEhAR69uzJrVu3+Pbbb0vn5AuhtBXWBw8e0KdPH4yNjfn7779LlFEiOzubHTt2cPXq\nVe7evcvjx4+xtrbG3t4+X6WoohIXF8fatWtxc3Ojfv36uLq6cufOHW7fvs2dO3e4c+cOampqJCUl\nYWBgQHJyMoIgkJqaSr9+/bC0tOTp06d4eXkRFxfHqVOn0NXVLfZ5fWy2b9+Op6cnT548YezYsfTt\n2/eTDu4VEXkfvkof1kqVKgnt27cX9PT0hAoVKggKCgpC3bp1hVWrVhXon5OTIzRo0EDmH1m7dm1h\n/vz5si2ie/fuCVWrVhW0tLQEX1/fAuOLu20s8uVz8OBBQVtbWzAyMhLMzc2F6tWrC3Z2dsKOHTuE\n1NTUt47Nzc0V5syZI/Tq1Uvo2bOnYGdnJ3Ts2FHmygIUyf909uzZgp6enuDm5iYkJCQIW7ZsEfT0\n9AQ7OzvBw8NDmDRpkqCrqyt4eHiU+B7OyMgQZsyYIejo6Aje3t4lmqMwHj16JFhbWwtqampCo0aN\nBHt7e2HevHnC/v37hcDAQCE5OblU15s7d66Aqr7sGlO5nUDPbIFWe4WpU6cKV65ckV2n5OTkAlvI\n2dnZgrOzs6Cvr1/q1+JNmC3/pyT4+voKVapUEZYuXVri++Dhw4dC27ZthXbt2gkLFy4U9u/fL1y6\ndEmYNWuWoK2tLaxdu7Zkwgl5vskqKiqy/5vMzExh3759QnBwsPDkyRMhIyND2LlzpyCRSIQlS5YI\nL168KPFaIiIiH5ePpDYKgiAIH6226MaNG+nevTv37t2TVYW6fv06Xbt2ZcuWLbRq1YoffviBVq1a\nkZGRQVhYGM+ePctXNeoVdevW5cGDB2zfvp1hw4Zx4sQJTE1NMTIyIiMjg7i4ODQ1NdHR0WHfvn3s\n3r2b69evc+nSpVJNqSXyeXDw4EEGDx7M6tWrGT16NLdv32bDhg2oqKjw999/M3ToUDp06ECnTp3o\n2bNnAcvkvHnz8PT0ZNq0aUilUhQUFFBVVSUiIgITExNGjBhRpOpVM2bM4LvvvmPWrFm4urpibW3N\npk2biImJYfPmzVSqVImbN2++V+YMZWVlZs+eTevWrfntt9/kRjaXhFf5KseOHcvhw4c/SJni3377\njal/vFZt6OUVyMkA/R7MmdMjX1951cukUimdOnVi7969nD59utSuRVmRmZnJ9OnT8fDwYO/evbRu\n3brEc82fP5+6devy999/I5VKycnJwcvLCx8fH1RVVWUZPUqCqqoqlSpVIj4+Hn19fczMzEhLSyM5\nORmJRIKVlRUtW7ZEQ0ODly9ffrZZW0REvkZCQkJYs2YNzZs3l7kAfSw+ucIBqampXLlyBX9/f9av\nX092djYKCgpUr16dgICAd84dHh6Ou7u7LJ+nmpoaOjo6xMTEEBISwuDBg1m6dCmQl+ahrLZJRT4N\nBEHgyJEjXL16leTkZK5cuUJERAR//fUXjx8/Zt++fezdu1fWPyEhgSVLliCVSrl79y7Xr1/H2dmZ\nLl260KhRIwBmzpzJX3/9RY0aNUhOTiYhIYHExERZ5S8dHR1q1aqFvb09/fr1o2bNmu+U88WLF3h6\nejJ9+nR27dqFhYVFqV6H5ORk9PT0CA4OLpUSkcHBwXTq1IlHjx6VgnRFR9ItBp6fh8v9ITcTDIZC\nkzUIewqmaMnNzeX27duyNEj+/v5oa2uzYsUKWR7eMpPzd/ntwox3j83JycHd3Z1Zs2ZhYmLChg0b\n0NbWLpEc2dnZeHh4MH36dIYPH8748ePZunUrCxYsQFtbG1dXV2xtbd87FVbv3r05c+YM5ubmhISE\nkJWVRWhoKKGhoZw8eRJ/f38MDQ0ZOHCgWP1JROQz4ty5c5ibmwN5ZYo3bdr09fmwFmXZ7OxsIiMj\nycjIoG7dusyZM4dZs2YBYGZmxoABA2jQoAGmpqZvfWp/5R8bGBjI9u3bUVZWZuXKlezcuZPvvvuu\ntE5L5BOka9euREZG0r17d8qXL0/VqlVxcHCgXbt2XLp0CQBfX18aNGjA1q1bWbNmDS1btuTChQu0\natWKffv2oaCgQE5ODseOHcPGxgbIU2xfFTJ4lX/2VXL07OxsTp8+zc6dO/n7778ZN26c7CHpFa92\nCdatW0efPn1kFi5fX1/Gjx/PrVu3ih1g9S7+/PNPdu/ezd69e0tcQhfygsUcHBxo27Ytf/zxR6H9\n5CpthRhihSlFW1vyKvbmeANICgJlHchOwKBmVapXry5LPh8cHMzFixfR0dGhXbt2sp9vvvmmxIEw\nkkLc3oWf5fQtRGHlTV0t4rXXaS8hyAf8FtDWuDLz5s2TlUQtLunp6bi7u7No0SJ0dXWZPXs2mpqa\nODo6Uq1aNaZNm0b79u1LNSjoyZMn+Pv7s2/fPjQ0NFi3rhBHXhERkc+K5ORk7O3tqV+/PkuWLBEV\n1qKwatUqxo4di0QiyVfHfurUqXIDUgICAujduzfPnj2jWrVqZGVlMXz4cLKzs/nuu+8wMzMTozi/\nQNLT09m1axfh4eHMnz+flJSUAlHh7u7upKWloaSkxIgRIyhfvjz29vb079+fgIAAWanGmjVrkpWV\nRffu3bGysip2dLmjoyPbtm3jp59+4vvvv8fS0pIKFSpw/fr1fMUIbGxs6NSpE02bNqV3794EBQWV\nuvU/NzeXBQsWsHz5ctzc3OjTp0+xxj958oQFCxbg7u7Ozz//zIIFC956PcpCYZVHeno6jx8/JjIy\nksePHxMbG4uRkRFmZmbo6emVfOI3KI7CKnf8djmNt6Lhrg/c3QOR56G2JbQcSe7GziX6bHrx4gUr\nVqxgzZo1mJqaMnHiRCwtLTlx4gQ//PADK1as4IcffhA/90RERIpFfHw8ZmZmhIeHiwprcXn27Bl+\nfn7UqVOHRo0aFfgA9vX1ZejQoaxZswY7OztCQ0PR09MrUt17kc+X8+fPM3jwYAwMDDA1NaVJkybv\nVMyCgoKoU6cOUqkURUVFjI2NCQkJAfIyCZiZmXHjxg02btxIz549adeuHTNmzGD//v1UrlyZevXq\nUa9ePdlW+5EjR6hVqxaOjo7UrFmTuLg4tm7dysGDB7l48SJNmjTBysqKBg0aEBwczKJFi0hJScHe\n3p4HDx7QokULVq1aVWbX6PLlyzg4ONCpUydWrVpVZOWlY8eOlCtXjjVr1hRJEfxQCuuHolQV1sw0\n2DEZzmwB427Q4Huo2wVUyuXNWcLrYWdnJys/Wr9+fVl7o0aNmDp1Kv379y/ZxCIiIp8c0dHRODs7\nExUVhZqaGmfPni2VeXNzc7l27RrHjh3jxIkTpKenU7NmTSpWrIibm5uosJYmHh4euLq64u7uTtu2\nbctsHZFPiw0bNjB16lRWrlxJ3759iz0+JSWFHj16kJiYiLW1NWPHjkVXVxdBENDQ0CAtLQ0zMzP8\n/PwoX748hw8fJjs7m+DgYEJCQoiKiiIjIwNra2vu37+Pp6cnnTp1YvLkybK646mpqfj7++Pn50dI\nSAj379/n9u3b5ObmynUdKCuSkpJo164dY8aMyVdL/W1s3bqVTZs2ceLEiSL1FxXWN8a/rrAusQMF\nJei0HtQrFZyzBNfj4sWL9OnTh/v37xcovfvHH38QFBTE9u3yzLwiIiKfEllZWURGRhIeHk54eDhJ\nSUmMHj06X4n21NRUevTowbfffkt4eDi6uroldsPJysriyJEjBAcHc/nyZU6cOIGOjg4dO3bExsYG\nTU1NHj16xKNHj5g2bdpHU1g/WpaAsmTNmjUyq1dZIQgCR48eZfPmzSxdurRUtx5Fis+2bduYMWMG\n/v7+1K1bt0RzaGhocPz48QLtEokET09Ptm/fztSpU1FRUaFx48YoKytjZWVFp06d5M63cOFC/ve/\n/9GxY0fGjBlDzZo1UVdXp3HjxnTu3FnWLzc3lydPnrx3ntLiUL58ebZs2YKVlRWpqanvrKkOcPr0\naczMzIq+iLzp0osn56dEURXTQnn22uvUDLAYDlUqIfz6nvP+y19//cXPP/9cQFmFPOv47t27S2ch\nERGRMiEuLo5ff/0Vb29v9PT0+Oabb6hduzZ3794lKSmJWbNm8ejRI9zc3Pjnn3/o1KkT33//Pb16\n9eLOnTvFXi8pKYn169ezbNky2a5kly5dWLx4caFxDtOmTXvf0ywx7xca+gmSlpaGh4cHUqmU1atX\nv3tACbh+/TpSqZQuXbrg4eHBjRs3ZMeWLVtGq1atqFixInZ2duzcuZP4+PgykUMkj9u3bzNu3DiO\nHj1aYmX1Xdja2rJt2zYaNGgAQMuWLfH393/rmAoVKuDi4sKhQ4d49OgRR44cYeHChUyePDlfP6lU\nSo0aNWTuKh/q6dXExISTJ0+yZ88eNDU10dLSomLFilSpUoUJEybw/PlzWd+goCD27dvHpEmTPohs\nXzzahvAsrFSnHDJkCKtXr8bT07PAPXT16lUkEgnnz5+Xe39lZ2cTFBTEuXPnOHr0KL6+vvn+/iIi\nImWLl5cXJiYmVK9enRcvXvDw4UNOnTrFP//8g5eXF6tXr6ZHjx40bdqU9PR0XF1defz4MY6OjqxZ\ns6bQohtZWVmcOnUKb29vWexPSEgIv/76K7Vq1eLixYvs3bsXf39/li1bxuDBg98rKLcs+WJcAhIS\nEvJlCtDR0UFFRYXq1aszZ84cVFRUCA0NxczMTKZ0FBdBEPDx8WHYsGE8ffqUpUuXMmbMGKRSKWlp\nady5c4fhw4djZWXF2LFjOX78OB4eHpw9e5YaNWrQtm1batWqRUJCAunp6TRs2BBTU1O+/fbbfKZ+\nkeLh6OiIiYkJv/322wdb8/Lly/Tq1YuwsLBi/e127NiBl5cXnp6eco9fvXoVR0dHOnXqxLJlyz5I\ncIwgCCQmJpKbmwvk/S8tWLAALy8vHBwcEAQBf39/fvjhB1xcXIo8r2SenMZCLKxCYVH1XxCSv157\nc3IJvIyEXstLzcIKeT7czs7OVK1alf3798vuzdjYWNauXcumTZuws7OjU6dOPHjwgNDQUK5cucK1\na9fQ19dHW1sbDQ0NFBQUuHz5Ms7OzsyfP18M0hIRKSPi4+P55ZdfuHr1Khs3bqRNmzZy+x06dIjw\n8HAMDQ2ZPHkySkpKTJgwgU6dOhEbG0vVqlWJiopizZo1NG7cGHV1dXx9fTly5AiGhoaoqKgQHR1N\n9+7d2bZtG0OHDmX48OHFVk6/ytKsZbHs33//zc6dO7l69SoNGjTg1KlT7Ny5k/nz51OuXDlq167N\nkSNHWL9+fbFrUx89elRmXZo2bRoDBgzgwoULXLlyhS1btnD58mWMjIwwNjZm1KhRtG/fXjY2Ozub\nGzduEBAQQHR0NFpaWigpKXH79m0CAwMJDw+nQYMGWFtbM27cuM+yPOHHpGnTpvzzzz/5ou7LmiNH\njtCvXz+uX79O7dq1izwuKCiIrl278uDBgwJKwM2bN7GxsWHevHksX76cwYMHM2HChNIWvcgEBQXh\n4+ODmpoaFStWxN7eXpa6600kI+U0ykv3KiqsedzyhQA3GHGkoOuEIiAIcHM/lNeBb1oBRS/rOn36\ndE6fPs2ZM2cK5Fd98eIFTk5OREdHExcXx7Nnz9DU1ERXV5eqVatiZGSEqqoq0dHReHl5MXjwYBYv\nXiwqrCIib5CRkcHTp0+JjY0lPj4eVVVVNDU1qVChApqammhoaHDq1Cn+/vtvwsPD0dfXp3r16vz5\n559UrlxZNo+joyN+fn7cvHnzrek5nz59iouLCydOnGDx4sXY29sDYGFhwcOHD3nx4gWqqqoMGTKE\nEydOoKurS/fu3bG1tZVlnPHz88PNzY0pU6bI4iqKi6iwfkBatmzJyJEjSUxMxMvLi8aNGzNkyBCZ\nXCdPnkQqlfLLL7/kGzdkyBB8fHwIDAykZs2aGBkZIQgCzZs3p3///nTr1q3QL/N3kZKSwo0bN9ix\nYwfbtm3j559/xsXF5b0qHH1NGBkZ0aFDBxo3boyVlVWZ+i5DXu3w8ePH4+3tXejTcGEIgkD16tVx\ndHSUObQD3Lt3DxsbGxYtWoS9vT2RkZFYWFgwevRoxo0bVxanUarIVVgLulICIHyYuLIyIz09Xfa/\nnpqaKreq1rtISEigZs2aJM54BGpvZC6JvQU7x0DS07yfgX9DY7u3Kqy5ubkcPHiQEydO4O3tzaVL\nl9764Kurq4uhoSHm5uYYGRnJimDcu3ePrKwsypcvT//+/cXCKiIi/xIVFcWUKVO4cOECsbGxpKSk\noKOjg66uLpUrVyYjI4PExERZIZnExESaNGnC0KFDadq0KTExMTg5OTF06FDKly9PXFwckZGR3Llz\nh9u3b/P7778zffp0uWvv27ePoUOHMmDAAGbOnEn58uXJyclh8uTJHD58mOvXr6OgoIAgCGX+cPkx\n9bcvMujqbVhZWfHHH39gamqKubk58+bNw93dHUVFRRQVFYmKiqJbt275FNbU1FS0tLRQV1cnKysL\niUTCvXv3Su3G0NDQoE2bNrRp04ZJkyYxf/586tevT/v27bGwsKBatWooKCigpKSEpqYm9erVE5XZ\n11i1ahWXL18mMDCQRYsWce/ePZSUlMpsvfv37yORSHj06BGtW7cu1n0gkUjw9vbm4MGDshKwDx48\n4NatW8yePVv21FyjRg3OnDmDlZUVWVlZH9TdQaRwFi1ahKurK4qKimRlZSEIAvfv38fQ0LBY82hq\natK1a1d2rrCC1j/DN20gNhjunYKbe+G7mdBuGERchrW9wLgDUL7Q+Xbt2sXUqVMZOHCgzLryNurU\nqUPlypXR0tLC2dm5WLKLiHxNZGVlsXLlSubNm8fw4cPx9vZGV1eXihUrvrVCnDzlcciQIURHR5Oe\nno62tjZNmzbll19+oWbNmlSrVq3AHElJSbi4uHDkyBE8PT1lWY/S0tLo2bMnOTk5nD59WpYP+0vf\nCfnqLKwSScQbLQIQy4sX9YiNjcXGxoaHDx+Sk5ODsnYAZHpCxg6QVgeUITcWFFqD6gaQ/PcFIiSW\nrpyxsbGcOHECf39/4uPjyc7OJjMzk5cvX/Lo0SOCg4NRV1cv3UW/AKytrXFwcMDJyalM1zl37hxD\nhgyhf//+zJhRhHqbcoiPj2fevHm0atUKOzs7uRb6J0+eYGlpyZgxYwpY/YuCRJImb2W5fQUh/x6+\npDB3VXnD5VlTvzALq7u7OzNnzsTPz4+kpCQ6duzIY0V9MBsAKhrQ5mf49wusKNv32dnZKI09Dhc2\nw+ProN8QajQFy2F5qa7Cz0FCNFz3Ad26CL4zC53rl19+oXbt2m91IUlOTmbv3r1UrFgRCwsLfvzx\nRw4cOEBGRkaxr4WIyNeAn58fo0aNomrVqqxcuRIjI6MyXzMsLIwRI0Zw7949oqOjZQU/Xs8hf/Xq\nVczNzWWllZs3b17mcr1CdAn4oGtHyGnNYc2aY8ycOZPhw4czceJETE1NuRdeAZR7QuoSUOkKSh0g\nzRNyjoL6FVBoJpuhtBXWt2Fvb4+enh5Lly4tduWlL53AwEC+++47hg8fTmhoKBkZGbRu3Zq2bdvS\ntGnTUrW8RkdH06JFCzZu3EjHjh1Lbd43iYiIoG3btixcuBAHB4dijf1gCisgFJKn9EuhSZMmLF26\nFCsrK7y8vPjf//7HsTQDyEqD63tg6XNQzAtyKqq/6Zv4+flh0d8fEg9Ayvm8xgbhENKCe0HnC82C\n0axZM1avXk2rVq3kHr958ybW1ta0bNmSmzdvMmLECCIjI1FSUmLFihUlE1ZE5AslJiYGFxcXzpw5\nw7Jly+jVq9cHsV5mZGQwatQoYmNjWbZsGbVq1ZL7nZWTk8Pu3bvx9fXl4cOHBAQElLlsr/iY+tsX\nl9aqeKQCGwArNm/ezJgxYwgKCqJGjRp5CojWZVCfCqoDIfMgpEwFhSagEZdPWf3QLF++nFu3btGy\nZUuuXr360eT4FDE1NeXPP/8kJyeHjh070qNHD8LCwnB2dkZHRwcnJyeuXbtWKmvp6+szZswYDh06\nVCrzFYaBgQHTpk1jyZIlZbqOSOHExcVx584dvv32WwB69+7N0aNHwfF/edbVxj1kympJCQwMzCt4\nkZsMlYdC5SGg+i0oaIH+bPr37y/XGnr27FliY2PfGnQ4ePBglixZwv79+/Hz8+P8+fOEhobyww8/\nvJfMIiJfEjk5OaxevRoTExOqVatGUFAQvXv3/mBb7ba2tjx79oyNGzdSp06dQg0sCgoK2Nvb06tX\nL3R0dD6IbJ8CX7SFVSLPGBF6ATgPxAF7gJrAj1DjGcQuBf3poNUTlLQh4rVxOQ/y3AIy5d9AH9LC\nCnn+MZs3b2b8+PE8efIENTW1DyvAZ0h0dDRbt25l6dKl9O3bl+7du2NmZvbWyMy3kZycTJs2bRgz\nZkyZuyD8+OOPXLhwgcmTJ/PDDz8U+e8tWlhLB0EQGD9+PLGxsaxfv14WaCXpNBFOrwbnHdDY7u2T\nxMhpy37t9YWl8DICsv/KyxLw0BFebAf9uZCbhZHmdpycnPjtt9949OgR/fr14+JtKWQ+gmqLoeJ/\nZVeFN55jK1WqxP379/NFJ4uIiPzHlStXGDFiBKqqqqxZs4aGDRt+0PWDg4OxsrIiMjLynTunycnJ\nODs7c+3aNaytrXFzc/tAUoouAWW3zusKa2YQJPwPkjeDtAJkPcprV2sCmQ+gnDlUmwvqTf4bU9gX\nc3iZiVxsrKysMDU1Zdy4cejr639scT4Lnj59yvLlyzl37hxXr17FwsKCAQMGUKVKFdTU1ChXrhxV\nqlRBR0enUKf6Q4cOMWrUKLp168bKlSvL/Ak8NzeXo0ePsmLFCkJCQnBzc6NTp07v/GCTfCO/vSj3\nsOQt3gfCV1LhUxAEMjIykEgkXLx4kQ4dOqCpqcnly5cxNDRkzZo1XLt2jX379rFp0yZZBTOJvAqJ\n71JYA1fDk4ugvBmi5sALT1BtCAmHQMcZ1YSVaGlpMWnSJDw8POjUqRNzdnYGIRvKWcBr9+CbCuu3\n337Lxo0badGixXtfExGRL4mEhASmT5/Ozp07WbBgAYMGDfoowUuzZs0iISGBZcuWvbPvn3/+yblz\n55gwYQKNGjX6oA+iYpaAsiI3ERI3QcJKyArNa5Mog5I+VJoGmrYgUQTVGqD0eeY+dXNzY/HixXz7\n7bf88ccfjBwpL7+QyOtUqVKFefPystqnpKSwa9cuPDw8SExMJC0tjaSkJCIjI3F0dMxXLe3ixYvs\n2bOHAwcOkJaWxurVq+nSpcsHkflVZbUuXbpw6NAhJk6cyIMHD6hTpw6mpqaMHz9etl0tUjqcOnWK\nUaNGERoaioaGBr/++ivt2rXjzJkzBAYGoqCgwIgRIwA4c+YM9vb2XLp0iZo1a5ZsQSM78JsNBreB\nXNAwBb3JkHYdkv0pX748ffr0ISIigt69ezNhwgTmHPj3gUXIguQLeX6vSvpERLSjVq1a/P7777KE\n40+fPi2dCyMi8oXg4+PDiBEjMDIyIjAw8KNWeIqPjy/S+pmZmZw5c4YqVapgaWlZ9oJ9QnzZFla1\n9pD+b/lMpbqg3hn0nUHNJH/HwoxUn7CFNTg4mLS0NBQUFNDU1OTly5dYWlry4MGDEm9xi/yHq6sr\nampqsgwAOTk5GBgYIJVK2bFjB2ZmZh894C0lJYWQkBCOHTvGggULuHXrVoHUKKKFtfhERETg4uLC\npUuXWLFiBT169CA0NBRbW1tsbW1RU1Pj8uXLXL9+nZ9++omFCxcikUhYsmQJM2bMoFevXri321pw\n4ndZWAH2OEJSe9BoBcFtQKUOaP8MabdpWiOQZs2a4eDgQIcOHYiLi0O3xQFI2A/JJ0DZMG+nKPMx\nJOyhSZMmZGdns2rVKlq0aCFmFREReQMXFxcCAgJ49uwZqampXLp0iWrVqhEQEMClS5fQ0tKiW7du\n6OnpAXluA6tWrWLJkiVUqlSpVGXZvXs3Gzdu5MCBA3KPp6amMnjwYA4fPkyDBg0YOHDgRzFQiRbW\nskJvN+REg5IxSP/NsfOFuHoOGDCAq1evUqVKFczMzPDx8eGHH36gfv36LFq0CEdHx48t4mdNhQoV\niIqKkr1/VarSwsJC9qDwsdHQ0KBZs2Y0a9aM6OhoBg4cSP/+/XFyciq2fJIf308WuWVYCxFBmPR+\na5UlN27coF27dkycOJHNmzfLlLxXPqDjx4+nSpUqREVFkZOTQ+fOnbG0tMTW1pYJEyYwcuRImjRp\nApWPQoNOb19MECA+FDIT4d5+SH0KsTcgOwt0hkKTZ/B4MkSOAa2euLi48OTJE0aPHo1UKuXx48dA\nR9DqATXWgtJ/wRf3dt4nPj6eunXrin6rIiKFsGjRItnrBQsW0KhRI5o1a8adO3fo168fly9fZtq0\nabi7u2NlZcXhw4fZtGkTDRs2ZOLEiaUqi4WFBU5OTmRnZ6OoWFA1O3LkCJGRkYSGhn5VgVav80Vb\nWL9kZs2aha+vL3p6erRv357JkycDeTXuO3bsSFBQkOjT+h7cvXuXjh07Ehoami/Aaf78+Tx8+JA1\na9Z8ROkKkpiYiLu7Ox4eHrRo0YKlS/OSnUokm+T2F4Sf8r0vTGEVthRt/S9FYU1OTqZOnTocO3aM\nRo0aydofP35Mo0aNeP78OTNnzmTlypUIgkB6ejp3797lm2/+M2V7eXkxa9YsTp06hba2doE1l9S+\ngAAAIABJREFUoqOjmTNnDr6+vkgkElRVVTE3N6dSpUqkp6ezYMECypUrV6iMWVlZHD16lLZt2+bL\nzSgiIvJ+REVFcenSJVq0aCHbrTpx4gQDBw6kffv2nDx5khYtWhAZGcnNmzdLff1GjRqxYcMGzMzM\ngLzYhePHj3Pnzh22bdvGjz/+yK+//lrq6xYHMehKpNgIgsDy5csZP348BgYG3Lx5k/Ll8woZ9OzZ\nk7S0NAYPHoydnZ24FVhC+vXrh6GhIX/++aesLTQ0FHNzc2JiYj7JqiJnzpxh+vTp+Pn5AaLCWhJ6\n9OhB69atZQ+BgiDg7OxMXFwcPj4+9OrVi3r16jFv3jy51WwEQWD48OHs2LEDFRUV3N3d6dTpP2tr\nz549qVy5MhMnTqRevXqf5H0kIiLyH3fu3MHMzIzFixczbNgwWSXCOnXqlOo648ePR0tLS+aKtnnz\nZn7//XdsbW1p0KABgwYN+ugZgUSXgC8Qyc+FHyuN9D8SiYRx48bJfvfv31/m+7J582a8vb35/fff\nCQ0NZdq0aXljehUij/f7y/MlsnLlSkxMTKhZsybDhg1DKpUSFBT0yZbFTUtL49ixY6VauWgc8wu0\nXaSlnJ4dSm3ND00bTsleP918hEeXAnjsNoC+YWF4enqya9cuAFLO/EXNw27EXLuIz/GjeE3uSpsK\ndwpOKIHdy+aSO2YoSd370iMrGcXUZwB4+gXKaofLq2wmIiLy6fH8+XMAtLS0kEqlVKpUibQ0eSkD\n3w9bW1umTJnCjBkzyMnJYceOHcycOZMff3xPn60vBFFh/cwZO3Ysjo6OqKr+VwdTU1OTwYMHExAQ\nIFqy3wNdXV1OnjyJs7Mz27Ztw8nJiTlz5rBs2bJPxir24sULLl26REBAAKtWrUJBQYH4+HiCgoKo\nX78+UET/xYI714xdWlBZLQypU4rc9hxdjSLP8SmQfu8xCpoaXK//E23Uy9O7d2+WLFlC+/btaaRw\nk5z4l2RFRCGtUI57mu355sx49NrLrz6VtWELEjVVpO3NyX0QQc6J0zgtWMb69etFZVVE5DMhLi6O\nAQMGYGdnJ4sV0dDQICYmJp/bUGnQrl07wsLC2Lt3L/Pnz0dFRYWePXuW6hqfM6LC+gUgz08uLi6O\n+vXrM3fuXKpWrUq7du1AqJsvV6PIu2nYsCEBAQGsWbOG3bt3s2DBArp37/6xxSIiIoIRI0Zw9uxZ\nmjdvTkJCAtra2mhoaKClpYWhoeHHFvGzpOZcJ2rOzSsCcVawLPBgUr6XFRqdWpNy9DwSFWUyX6SS\nGp0AwMu70STdf0pOZg6ZGc/JOXEKITOTNOMmkJuLQrs2bN++HQsLiw9+XiIiIiVj+/btWFpasnbt\nWurXr8/ly5dxdHRk/fr1pV6SW1lZGRcXFwYNGsTChQtxdnYuNBf414iosH6h+Pr64uLigq6uLiNH\njsz74s2SQM1B0GgxKBYe1CGSH6lUyqhRoxg1alSZzC8IAh4eHly7do3ExESZVVxDQwN9fX3ZT0ZG\nBjExMYSFhbF+/XomTJiAr68vmzdvZvXq1UyfPp1BgwZx69YtlJXfr0yoCAWU1czQSMIa9IasbMp1\nb49idV3OOm0lJz0bBTUltOrrUcGoCgpqygiZYSgOGUzmzD9Q27cbaeNGSCQSLNTlmLJFREQ+WZ49\ne0bdunVRV1enUaNGxMTE8OOPP/LHH39w8+ZNTExM3j1JMXBxcWH48OFiQKUcRIX1IyBZ+kaDvL/C\nm/kZAR7Ln094cz7AwcGBdevWERwcTL169XB1deWHqQHwcBNEeYPxFKg9DBTUCsoDCOPffg4ipUNK\nSgoTJ04kICAAR0dHDA0NSUtLY/z4vD+AsbExlSpVQklJCXV1dfT19alatSoBAQEYGxsTFRWFq6sr\nhw8fxtTUFMiLNI2IiCAiIgJIQ14uN4kkK3/DOPklh0X+IysyBrKy0ejcmvI9LMl++hzrS5PQqFWZ\n2wuPolGjIlVt6hG5/zb3PcLJ2eGJ0uCBKDQp3S80ERGRD0d4eDjW1tYApKeno6amRoUKFZgxYwbj\nx4/n2LFjpeoiJpVKRWW1EMQsAWWExLmQA/XltJWBwvqK+Ph45s+fz/r160mo5QDa34L/ZFCpALnZ\nYDMJzEeCYn6fOlFhLVtu3brFvHnzOHjwIOXKlaNz586kpqby4sULDA0NKVeuXL4cgQAxMTHo6v5X\nkS01NZXevXtjZmbG7NmzCQ8PJzU1lfbt25OYmEj9+vW5ffs+YAm0fm2mKQUF6lBQYZV6yPdLhc/P\nN7U0aMexfO+zXyaT8CwH5Vp6RPSbSnLATYTUdCrYtmF9vxF07doVDY2v7zqJiHwp5Obmoqenx+XL\nl6lRowbVqlXj+PHjDBkyhMqVK3Po0CH27t37SbiJfSjELAEiZUblypVZtGgRkydPRtvUEcJ9oao5\nvAwDRQlcWA8SKbQf87FF/Wo4c+YM1tbW6Ovro6SkRIUKFahTpw4GBgZUrFiRJUuWcP36ddLS0ggJ\nCcHf35+kpCRZ2jLIyxdYo0YN7O3tmTIlTwF9lQv0+fPnpKSkoKGhgURiDJwDqgElLBkqko943/OE\njVxJTkIKCpU1yYp5jkRZES17G2r8bxISqZQ+tPrYYoqIiLwngYGBaGtrU6tWrX93rPLa1NTUGDhw\nIIaGhmRny7MuiZQFosL6BSDZKacx582GytDzILy8D9tb5FXXUa4AlZvBpX9EhfUDsXfvXr7//nsA\nLC0tGTZsGObm5rItJS8vL+7cucOWLVtQVVWlcePGNG7cuMA8AQEB5ObmMmjQIFnEuUQS8UavMKAR\neS4BNxEV1qJjh6ec1rySx9kvklDSrUi9nVORmJki5OaSEfoYhQoaSMQACRGRz5r/clc/AxYBXZBI\nNpGU1IeXL1+ycuVKZs+eja2tLfb29h9P0K8QUWEtK1oUfkgY+u7hknVyGnXltBkUUR6AFhIEByPm\nGLpw//59XF1dsbS0JDY2luRhKeL2ZRkiCAILFy6UJaM/cOAA3bp1K9AvMzOTFi1ayD32Oo8ePQJg\n0KBBKCoq/htJmgBIgOZAA2AVYAwMAmYCAwElYFvBCU8VLCYA4v3wOhV5iS99SehpxuwbOZwevYXQ\n0BmMHDmS2bNni6mq3kAQBGJjY7l9+zZNmzYVS8SKfEbE8kpZfZVjWk1NjaysLBQUFAgICMDW1vZj\nCvhVIpoDvkKcnZ3x8fGhZs2aLFu2jPr16+fL4ypS+nzzzTdMnjwZdXV1goKCClVIO3fuzOnTpzE3\nN2ft2rX5fIUEQWD06NHs2bOHCRMmEBcXx/Xr17lw4QL+/v7AGfI+YE+Sp6zaArnAZaAG4vNp6aCp\nqcnSpUu5evUq9+7dIyQkhA4dOpCYmPixRStTQkJC8PHxITIyUu7xly9fcvjwYbZu3YqVlRVSqRR9\nfX06duzInj173nt9QRAIDQ1l7ty51KtXj2bNmvHLL78wbNgwfH1933t+EZE80oGFQHfARtYaGxuL\ntrY2U6ZMYe3ataxatepjCfjVIgZdlRFyLaT/UmIL6ws5bQaFTFDAJeDftR3yfjdq1IitW7fSpEmT\ndwsj8t5s3rwZHR0dLC0t31kq19bWloMHD+Zra9u2Lebm5ixdupQePXqwdevWAg8ZeS4B0cB1oA3w\nBJhN3o3TA3hV715+ab+CFtavF/kuAeBL3wJtgiAwcuRIbt26xeHDhylX7stLGbd7925GjBhBVlYW\ntWrVomfPnjRo0ABbW1vKlSvH0qVL+eOPP2jWrBmVK1emS5cuTJ06FU1NTTw8PGjcuPF75ZM8duwY\nAwcOREFBge7duzN48GAyMzMJDAwkNjaWM2fOcP78+VI8Y5GvEUEQkEptgXhgRL5jCQm90NXVRU1N\njcTERPr168f27ds/ipwfEzHoSuSDY2hoSHh4uKiwfiAGDRr0zj6Sqv++yF4FymMh8wggAAIBZ+8R\ncP4ZqP3D7t0/snt3NWCGnFn0//0B0AR+KrKMkqv53+9vZiO3X6xc3xT4WZ6rwVeARCLBzc2Nzp07\ns2vXLn7++S11md+Du3fvsnbtWkJCQpg5cyZt2rQpk3Xe5PTp04wYMYL9+/fTuHFjvL29Wbt2Lb//\n/jvdunXjwIED+Pj4MG3aNDp37oy6ujrlypVjxYoV9O7dm6ZNm5Z4bXd3dw4cOMCxY8fw8PDAxsZG\n5u8tCAL169cnJCSE48ePl9bpinxFSCR/v/YuG/AAHgDjCvQtX748rVu3ZtKkSSgqKpZqCWyRoiEq\nrGVEUayoxR2v8rzglmNmTIXC52hQ+PwaGhokJyeXRDSRUqZ3797ExMTA06cgZABZoDEaKm6GtHOQ\nEwE5NyD7DqS4APMBebWlD8mdXxBGyG1/nTeV1a8dNVKL1V8qlaKurk7FihXLRJ6UlBRMTU1xcXHh\n6NGjHD16FAcHBxYsWED16tXLZM3ExEQuX76MjY0NJ0+epGXLlgB0796dAQMGANCvXz8AnJycGDRo\nEBMnTgRAQUGBKVOmMH369BKvv3z5clatWsXMmTOZOXMm9erVkx07dOgQAwcOJCsri5SUFJycnEq8\njogIpAL/I89LMhAoX6DHxYsXOXXqFIaGhqxZswZFRVF9+tCIPqxfKRERERgYGHxsMb5qwsLCGDdu\nHN7e3v/m8ZOSp6xOAOX2kDQZkgZCuhugBuquUCkMGAKIlaw+NWJiYtDT03uvOZ4/f46/vz9ZWfkL\nO7x8+ZKsrCxmzZrFunXraNmyJVFRUYwZM4bc3Nz3WrMwFi1ahI2NDX/++ScdOnSQtZcvX56YmBjs\n7OxYu3YtAwcOZMyYMfTp04dZs2bh4eFBw4YN2bBhA46Ojty/f7/Ya0dERDBu3DhOnjzJwIEDZcrq\n8+fPsba2Ztq0aTg6OpKQkEB8fDzLly8vtfMW+dqIBxYAVYBfkKesAuzYsYNJkybx+PFjvvvuO+7c\nuQNAQkICt2/f/lDCftWIjwifKKacldPaqFhzaGXGFGh7qZz3hRoWFibWm/9AyKyXggDBV+HsIfDZ\nCwmR0GQQjI/E9eo5UNADrTOQNA1SloH6CNC6Dgo13pgx/UOfQpGYw8QCbdNY/BEkKR2yktNJCY9D\n4xsdlMq9PSgxLS2NBw8eULVq1bf2K4ykpCR8fHwYN24c+vr6xMfH06NHD+bPn0+FChXQ0dHB2NiY\nX3/9lVmzZuHs7ExUVBTVqlXDy8uLvn0L+ta+LwEBARw6dIguXboUOKarq8ucOXNwc3PD0NCQFStW\nUKlSJdnxXr16MXXqVHbu3ImxsTHDhg3Dzc2tyH6s6urqGBgYsGfPHsaM+S/lnlQq5eTJkwAy/0Gx\nKpBIyXlIXoBqZ8CavCwr8lmx4hx5iu0MYDFHjthgbKzJo0ePqF69Ovfu3fsQAn/ViEFXnyjyFFaF\nQiKpLtK+QJs8ZRXyFNbw8HDMzMx4+vTpewVCiBQNmcL6v1nguxGse4OiLdSyAIV/nxkfBsBGK5Aq\nQ60fodF8UKrwKqNKfnYUvpbgV8rCy13++wJtoch/+PkcFdaUlBRWr17N4sWLqVixIo8ePUJFRQWJ\nRIJUKsXAwIDatWuTlpZGTk4OSkpKPH78mPr16+Pu7l5omcaoqCjOnj3LhQsXyM3NpU6dOhgaGhIT\nE4Orqyt169ZlyZIltGjRgrt37zJ69Gjs7e0ZOjTPPyg2NpbZs2fj7e3N1atXqVq1Kn5+fvTt2xdb\nW1scHBywtrYulTKRbm5uzJkzh5CQECpUKNztqChs27aNIUOGYGxszNGjR6lSpUqRxj18+BAjIyOe\nPn2aTylNT0+nadOmLFu2TK4yLSLyJpmZmRw8eBBra2vKly9PWloaW7ZsYfjwScBo8oJU/6VR1/9e\nCwI8Xw1RwcAGIIr/0v1lcfz4eRISErC3t+fZs2dfxcPTx9TfRIX1E6UsFVZXV1eysrJYvPjzUyY+\nRyRXgXOH4c8RsPkiVKoC++V0vJ4B2cmg8lq+SlFh/eD07t2bkJAQdu7cScOGDcnNzeXFi7wUHTk5\nOYSHhxMREYGGhgaKioqy7fvOnTsXmov14cOHmJiY0L59e1q3bo2SkhJhYWGEhYWRk5PD/PnzMTMz\nIzc3lwkTJhAeHo6joyPz5s3j5MmT+Xxjp0yZQlBQEDt37kRZWZmoqCh27drF//73PwwMDBg7diw2\nNjYoKCiU6PzDwsJo2bIlV65coVatWiWa402Sk5NxdXVl8+bNtG/fnvXr16Ovr//WMTk5OaiqqnL3\n7l3q1q2b79i+ffsYOHAgxsbGuLi40KdPn1KRU+TLQhAEPD09cXV1leVRbdOmDQcOHKBly5bs328D\nGOUf9LrC+nIHPJ0NGU5ACtAXeHUvzkNVdT7ffPMNzZo1Y/Xq1fmqEX6piAqrSAHKUmGdMWMGL1++\n5K+//no/IUWKhOQqsH05XPUD+9HQooN8hfWWnDZRYf3gnD9/HgcHBx48eFBqcx48eJAVK1Zw5MiR\nt/abNWsWx48f5+rVq8THxzN16lT279/Pzp07ZdH2qamp9OvXj9jYWNauXYuuri7Vq1cnMzOTjRs3\nsm7dOhISEti+fTtmZmbFlnXnzp0sX76cc+fOvdVaGxsby/3798nOziY7OxtDQ0Nq16791rkzMzMZ\nM2YMmZmZbNiw4Z2yrF+/nnnz5hEQEEC1atUKzLV582bmzp2Ll5cXKSkptG9f8LNQ5Ovk5s2bTJs2\njTt37rBu3Tqsra05ePAgjx8/xtLSEiMjIyQSOYGqryus4TZQeTg86kJeMZZUoCVgAngSHR3w3n7r\nnxtiWqvPnAH8I7d9GwXT28xgaoG2kDef8ACoU+T163BHTqv8qjIGBJPZ34y470bg+9dIACKoJ7ev\nSPGQ3C3kwLMYUFGDU3vyfo4/BXQ+pGgfDU/sCrT1pWhJ3uX9rwD8ztz3kulNFvFrvvexms9JUUko\ntflTU1Px9PRETU1+/tvXUVZWJkElBOXyOexV7EvrpVJyWpTHurMZfX7RZoBLFVTUpUzdJ7B71TN6\n9GvDk/BM5s6dy5QpUxg2bBhDhw5l9+7dtG7dmrNnz9KqVatiydurVy8mT57Mxo0bGTBgQD6rcVpa\nGqtWreKff/4hOjqaBg0aoKSkhIKCArdu3aJp06YMHz6c7t27o6SkJPf8bG1tcXV1LZIszs7OxMfH\n07FjR06dOoWu7n8p1ZSVlVFQUODhw4f07dtXFvzyLsutyJdPcHAwNjY2TJgwgfnz59OgQV7KnHdV\nECyAUnVI2A34kGeJ9SYv9dV+wPerU1Y/NqLC+oliSCg7eXfuTvnKKqQlF0xOr18pCgCpVgVyU9Le\nT0CRwnkWA9/Vg6T/lB4FK0to1wZp2zYotwJJ62ckqWu/c6pueBdsHAkH6VV68haT/hSsWnSMdgXa\nXlI2KZ7Kmkq1K5D2MoOQkBCMjY1LNEd2djaLFi0iMDCQixcvYmFhwaZNm946JiEhAS0tLaLvJVOn\npRaKSnn+5W1/qEYH83RWTnhC92p3qG6ozJhl1eg7Woe+o3X4Z2gLpk6dyrVr1xg9ejTt27enb9++\nTJs2Ta7S+Ir09HSuX7/OkydPZH65lSpV4vz589SsWRMnJyecnJxo3rw5FStWlB0zNTVly5YtNGvW\nLJ/bQXp6Ol5eXixfvhwnJycsLS1p164dCQkJPH78GC0tLe7cucP169eLFdU/efJk0tLS6Nq1K+fP\nn8+nQPfs2RMnJyfmzp3Lo0ePaN26NUeOHCnx303k0+L+/ftkZmbSsGHDYo07d+4cHTt2ZNKkSW/t\nJwhdC7S9nj4y5/T35IwZDcnRoKwKWySgPgwYVix5REoHMeLmKyTrZghK9d6+dSdSQsKDwUI/n7Kq\n6umO2j5P1A77oDJtUqkExYiUHUpqimhW0yAuLq7Ec7i5ubFlyxbs7e05evQo27ZtQ0tL661jQkND\nGTVqFNbOtXh0M4nAff+59ejVVGauZ2223a7HgN90mdTzAc6t75EQn826det4+PAhd+/excLCgh49\nevD06VPi4uKoWbNmvjVyc3NZtWoVfn5+aGpqMnDgQDZv3oylpSVWVlYMHjwYCwsL/Pz+8y1ZvXo1\nLi4ufP/993h5eeHt7U2LFi0K+MiqqqoyYMAA/Pz8CA4Opk+fPoSFhSEIAi1btkRPTw8HBwciIiLo\n379/sa7nrFmzMDAwoFGjRmzevFnWXqlSJYyNjXFycmLhwoWoqKjk5TQW+awJCgqiS5cumJubY2lp\nyalTp4o81t3dnQkTJtC6dev3lkM4dxapmRnseA5NbMDvLf5YImWOaGEtJvK3KeX7721DXiBA0Z/8\nh+BWhF6WRZ5PhlRK9oMnZD98gmKtau/uL1J0sjLzfntegW/qg6oaigbP5HbdjH2BtkHsLEvpRIpA\nUmwqz8MTS+T/CeDt7c2KFSs4duzYO1PHHX7NSVloJqChqUjnUQakJ2fz+E4Spt3zbznqVFXCuq8W\nFj01WeMahUPDYA7svcCaNWvQ1tamYcOG+Pr6YmpqioWFBb179+bAgQOyYJDIyEhGjx6NoaEhmZmZ\nKCsr4+vrS0ZGBjVq1ODGjRtYWVnRunVrVFRUsLOzo3HjxsW+Brq6ugwYMEBWYOB9kUgktGnThj17\n9jBz5kxZ5ThBEHjy5AkWFhYsX778X79E8YHwc0UQBObMmcNff/3FjBkzGDlyJEOGDCE9PZ179+5x\n6dIlbt68SdeuXWnfvr3soUkQBJYtW8bKlStRVlbGz8+PRo2KlwZSLsnJSJo0BfXy0GkIeP4JXZxl\nh7VznxQYUkX6tECbAQX94a9gKnfJWGrKbRcRg66KTVH96uQrqxAiR2ENk6PwqlG0Lfv4QnxVoyiY\nD/Il/1l4no5fjERJEZ0FYwnBpEhricgnIiKCbdu2sfr4RpQ0lFHXq0CHv/vJjlehaB9gICqsH4Oz\nb3xxuPYMpWIVJXzWlczC2rx5c2bPns133333zr6H34iqc6x2jokb67J19iO+G66PzcC8FFCKhQRc\nHt8Wh8+fKjRv3pyjR49y/Phxjh07Ru/evalatSodO3akefPmLFq0CMj7Ym/dujXPnz9n4cKFWFpa\noqWlRU5ODgYGBvTu3Ztly5Z9kkpft27dGDx4cIGcs7t37+ann34iKirqvVNwiXwc7t27x7lz5zh7\n9ixXrlxh//79VK1aFVtbW5SV84qkXLx4kbZt21KvXj08PT0JDw8nOjqaSpUqERwcTP369bl+/Tom\nJialdv8uWbKEsLAw1gxfDSE34Dd78AkGoPK3BZVV+PIV1o+pv4kuAV8p6lYtSDkUQG6q6Mv6PgQG\nBsp8/6JOhxFzNoJ2qwpG0Yt8HuTmCpw/mMivy98s1lB0KlasSFBQUInGmlhqMev7IHSqq1C/1btT\n5Fg7aPP48WOCgoJITU1FSUmJsWPHUqNGDRQUFNiwYQPu7u7s27cPyPuyadu2Lffv38fW1lbmpqCg\noMCYMWPYtm1bieQua9zc3Lh9+7bcvKvm5uZYW1tTp06dMqv6JVJ2BAYG0rZtW06ePEn58uU5evQo\nVatW5cqVKxw8eJC9e/fStGlTHjx4wK5du/j999/ZuHEjEomEcuXKAXkFJVRUVJg5cyZJSUmlIldu\nbi5bt26lc+fOkJsLa2eDmVWpzC1SMkSXgK8UjW7tSNp1lKg+E0ndfQR19YJBWiJvJyUlhRYtWgB5\nioCimhK9zo1GUbXwQBeRTxuJBIybqzOx631O+bx8p9/pmwQGBnLjxg1GjRpVovUnbq1P6vN0NLUL\n3kPXTiUQcTuVCpUVSU/N5dmTTB7cTiUhIYGQkBBOnz6NkVH+jCO1a9fGx8cHW1tbrl27RvXq1Zk+\nfTo3btygWbNmnDlzRlah6tatW/Tv3/+TtK6eP3+e6dOny81zqa+vz9q1a2ncuLFYCOUzRCqVIggC\n9vb2tG3bVpZ838TEhCdPnsitHjdkyBDs7e1JSUlBWVkZIyMjEhMTGTJkCE5OTuzateu97mNBEPjr\nr79QUlLCzs4OhkyBl8/gz0/zge5rQVRYv1IkUil6G2YR3e83+vfvj69v0VINieTRhlMIarlU6NCE\nxFPXGTFiBN+ObEPF+rrvHvw29r3xISuvcEpKIWPl5awvauGVSgWbsuV7m6AoZ/0MOes8UC+YdD5J\nTp1u4xz5JQ2jFOSnJ6pHRL73TyiYbaEwVxn5/CeTRCJh+fG6jLG6R2BgIDY2NrJjb27fv6ILeQEh\nR44cYd68eUycOBGt7xdxikX5+tUgUs7o/OnrpFIJmtpKJMZnUa6iIlLpf/eD25gI0pJzqNtMA/UK\nClSsooS/93PKlSvHunXrZHla38TMzIyOHTty4MABhg0bxoULFzh+/DgA3vda8W2rvMo9dbu9YMkv\nj7G1tf3kKkh1796d+fPno6uri6WlZYGt//j4+K+iytCXyKuk+3PnziU2NpagoCB8fX0ZN24choaG\nDBo0CBMTE5o3by4bY2JiwpYtW9iyZQuzZs1i5syZKCsrs27dOlq1asX69etlFeKKiiAISCQSnj59\nyrBhwwgNDWXnzp15iu+lk/DbclB9d2o6kbJD9GEtJvLS9wB0xD/f+zf94l7HnMB3riMvf2Uk8rcp\nxxcpOEt+MYKsuBc8qf8TZ8+eFVPBFIM2/Be1el6St00Um6yEqmrBp/oKezILTiBHQZSriIoKawFK\nX2EFE0LyvzcxwcLCAgMDA27fvo2JiQmGI/airFrQgteFUyxduhQ3NzecnZ355ZdfuFyuU74+aSk5\nqEQ/Ijkxl5QkAakCqKlLWLVenacPUpFIoElXXbIzc7m2K5wHdzP+z955h0V1fH/4XQQEO1JEQMWG\ngopIERELYheMQePXKNZojBrFYI3G3nsvP6NYsGKvqNgL2HvBFlABO00BpezO7w8iCWFBkMUFve/z\n7JMwO3fmc929d889c+Yc3scraNdHD+v6RWjbS48DvtFsXxrFi6fJNHAvjr1rUab8FMEZqBi0AAAg\nAElEQVS7t4lZpq8COHjwIL///jtXrlyhevXUvMv95itwbFkcLe1/zunKiXfM7JFEcHAwRYsWzWy4\ndCgUCiIjIzE0zLvcwh8r8504cYLg4GDmzp2LtbU1ZmZmvHjxgps3b9KhQwfevXuXtkwsUbBISkrC\nzc2NCxcuYGZmxvLlywkLC2PmzJk8ePAgQ4neoKAgjI2NqVevHkuWLOF//0vdN3Dnzh1cXFwIDw/P\ntPLcfzl+/Dj/+9//GDBgACtXrqRnz55MmDAh7XhbW1uWLl2qkswDH8lJ/vb8hFQ4QEJtaBnqMXr0\naLy8vDh06FC+XA7M7xh7tSf+2iN0dO6pW4qEinB1deXDhw9cvHgRV1dX1q5di1XyK9oNNmPTpMcc\nWP4MhVygoSFDJi+Bvr4+p06dSksj9SQ4gYB1r7l4KIYXoYkkJwnKlJVRrKQGxYrL+PBeEP1GgdNP\npbH7riyJ8Slc3vOCIiU1+WWSEQ7NitHbKYQ9q6LRLiyjbS893Lqnvp7cT+T84Th8Jr2molVhNDU/\nfRtv1aoVy5Yto1atWmhpaeHh4UGDthmr/Ng1KY6ra+20BP3NmjVj4sSJBAYGMmrUKLS0tGjcuDGL\nFy+mWLFivHnzhhEjRrB//35KlCjB69evefz4scpKun5ES0uLUaNGMWrUKA4fPsycOXN4/PgxYWFh\nlClThpo1azJz5kzJWC3AaGtrc+TIESIiIjAwMGDnzp30798fExOTDBvtAOrXr8/s2bORy+VERUWl\ntdeoUYMKFSpw+fJlnJ2dPzlvcnIygwYNwtvbm5s3b7Jjxw7q16+f9v6NGzd4+fIltra2qjlRic9G\nMlhVRMrb/xh6JeyUd8wDlHljlXmyYIDS4zcPsuHOumVUmD2AsiM6A8rLvUoox2R4J2459Fe3jG8G\nEfWfa610zrypyriHebq/hyz492bEk1i4pOBaPYrVI0Mob1mExVfsKKanhVAIWmr6U6RIER4/fszK\nlSu5cuUKfjvv0KZPGX5bXoly1XQorqdJeVl4hnkv/6vKnOMPqSnmKvMIgLl7y9O+ygN+GmuU7pgK\n1QpToVphfvi1NInvFQTJHDKM+99VHJlMxt69ezl8+DDnz5+nbt26PH+yh7IVtDMcu3r1ahYsWIAQ\nAh8fH1q3bk2ZMmUICAjA3NycpUuXYm5ujpaWFvr6+ri4uBAWFsbOnTvp169fni/Nt2zZMnUjjMRX\niaamJiNHjmTTpk0EBQVlmZ6qTZs2zJs3j59+Su+VLFasGMnJyek7B6beN5JTIPQZhDyDQS03E3f4\nPHGmxdk8ujkyWQvMcOfpvw679DiRChYJ2fbWSuQdUkhADrmGVYa2Wm8z7gi+XMIm0zHqce2T85yg\nvtL2JgSl+1uZsQrKDdbrZNR0AUcAkiJeE9zIC2PvjpQZ2F4yWD/Fin+MphvhUH8WvNsESvd8JGZz\nTGWPj8rqO2RW4l7Z/VTZmMoiS6KUtFXNZB4lvDXKaPgkFMoY75WkRGQRErI9j35UXIa2mNIZ51GW\n1q1wJh9ECoUytBX5T1q51y8V2Bu/wbBsIT68Fzg2LUJCnIJylbVobjuX69evs2XLFtzc3LCysqJr\n165KN4vkBCEE7u7utGzZkg5ek5T2efwfQ/sjnwo7GjNmDFOnpqbis7W1RS6X065dOyZOnPhJXe/e\nvcPf358dO3awbds2ANauXUubNm2wsLDgwYMHeRoeIPH1EhkZSZ06dfDw8GDUqFHZKn1atmxZzpw5\nQ5UqqfHgz549w97enoCAAGrWrJnWb+c0GaNXwOMXYGoA5mXh7NtypLyKptI5H3QsU2+2B/gnHd2H\nD4LRv7xDSxs2r5Qy6oAUEiCRD9A2NaTakbncrtULfc/mFNCqml+Uey9g93VYchJ8umVirEp8FUwb\nHke3IXoMm2vIk4dJ3L38gaIlNHh8P5nTp09jZGTE3bt3VWqoxcbGEhERQcWKqq9K17x5c27fvs2e\nPXu4evUqkLr0OXbs2E+GGOzcuZOePXuio6MDQL9+/ejRowevXr1CJpNJxqrEZ/OxWtngwYOzZaxG\nRUXx9u1bSpYsiY+PD7Nnz+b58+cMHjw4zVh9//49CxcuZNE82DwB6lpB9Dt4/Bw8q/5OUVd7ZDIZ\nKVGxvNt3lvUJCXx4DyH3U/DfnohNXS2WbZPy++YHJINVIg2dSiboezbjmkl79IsUY9CgQbRr1y7T\n3cffKu/evWNrIIzaDZ3sYMNP4GLx6eMkCh5CCHZvTOTSmSS23jIDoEJVbSpUTfUoN3IDa9ZlNcRn\n07FjRxwcHGjVqpWS0hO5o3HjxjRu3JiUlBTi4+MpXrx4tlNCde7cmcePHxMREUG3bt1o0KABANu2\nbaNq1Ry45SUk/sOoUaMoXLgw9vb2HD58OC1tYGbo6elRunRpqlevjp2dHatWraJ+/fpp3+XAwECa\nNGmClZUV0/vBrI1w7CqkpEDJYhD1biA1EgN5fzeEpx4j0LGuSrBxCjq6MipW1cT/WlFMy2dcgZFQ\nD1JIQFbczrgB6VpNywxtuQ0JeEHGXI/BSkIPAGpxK93fJzJJtZPTkIB/I4TAIngf1+ae5t6Gq9QZ\n2ogakzqgoan8wl1DP6XtXwv//nz2bk3i9/4fqFpCsKYHWP17I7tRxmMBKSTgXxSUkABdkcDJQ0ks\nn5FA5GsFc9eWwLyucs/hfzMMqIrBgwezbt06fHx8qNfhF6V9PjckQNXExMRQqVIlzp07J2Ubkcg1\nc+bM4datW6xb9+mHwY/pqD4Sn1iIN28E+/fBhLGCN39Xxq5oAm9iQLcwNLGFoFvwTEOfog1siD95\nhbILh1CqS6t0IQH/Rtl9ygDVFCkoSKjTfpMM1qxQYrDGV83ohSgalUl1lbLZO0dlBut7MibyL67k\n4nj0nzyOH1Fe7jXjBRejZO3/34Zt1L1XnBm8l0J6xWm8pY/Sub4Fg1UIwd6tyfze/wN+R4vSYlNG\nwwmAOdn8Xq9Vko1BmXGpo6Qts9h/ZVmIYpW0KbPblM2dWVYjZc8tlZS0Pcvk+EZ5f+0nJmT8900o\nojyHYsJ/rrXkZEG/ju94Gqqg3zBd2nXWRlMzdTxT3qhebCYoFAqOHz9Ox44dcXV1xc3NDT09PYyN\njbGxsUFXN3s5IYUQhISEEBISgrm5ORUrVsxWZoGcMGnSJO7fv59vK2VJFCzevHmDpaUlx48fz3LT\nlTLiEwvh3lpBXByMHS/DvCJ07SxITAQzM7hxHYzKQPnysGwFbNkwlTZt2mBjk7mT6Y3STcySwfql\nkUICJLKkdHUj3Pf1xMd0Gu9CXlO80rcXn/bmlYIBXRKIeiNYu7cI1raFYJO6VUnkFTNGJaBQgP/l\nkmhpqS/Nm4aGBs2aNePRo0fs37+fAwcOEB8fz4sXL/jrr79o3bo1PXr0yDLJf0REBH379uXq1atU\nq1aNp0+f8uLFCxYtWkTnzp2znWs1KxQKBUuXLuXUqVO5HktCAsDAwIDp06fTo0cPAgMDs/1w9pHB\n3jL69RWEhUHTZjKMjQXW1nDyJEyYLKNHz9R+GhoyRo8erXL9EnmD5GHNCsnDmsaJEUdQJMupOz81\nH17YgVsUKVuSkpbGrNf1UqrhayAyMhKHekZ810mLYRMKp3najIcpc10ieVj/TWYe1v+u1qdk0k+Z\nppfZm1tZIYPseFhTUgT1zGPYcrQEVapnPNEv6WHNitevX7Nz507GjRvHwYMHleaIVCgU1K1bl2bN\nmjFp0iS0tVPDNq5du0avXr0IDw/HysqKpk2b0qdPH+RyOYaGhpkaBwkJCejq6mbI1Xzq1Cn69u3L\n/ft5Exoh8W0ihMDT05MPHz6wffv2bMdYxyemXrfb/ASTJwr2H5TRwElw/y8ZuroZ77tFC8s/Oabk\nYf0HKSQgv3I8m0ZFZmR36TM04zzKqgwlFM0YJ1g4UUkVJZRXGXpAxtiyRCXWzwXqZmiLCIddtafz\n/fVRFDUrxfpSw9Esok1KQhJjh41mzJgxX2XRge7du2NoaMjcuXNVO/AgJf9WGcOJlceWZsYHJW0Z\nQ66VG4iZpc50VfIdVnZdKEuFqsxYVpbpKQ8M1pRMjG3NEllfkzt27GDevHkEBmasCpcfWbp0KYsX\nL6ZTp064urri7OycttyflJRE6dKlefLkCfr6GT+g8PBwHj16xLhx43j06BEaGhrExsbi5OREmzZt\naNu2LSVLlmTu3LksWbKEpKQkChcuTOPGjWnbti0uLi5oa2vTuHFjevbsma2UWBISOSExMREzMzOu\nXLmSVpQju6xZs4YVK1agr69P5cqVWbRoUR6p/LZQp/0mJeKRyBZFzfSwHNCISyP3IJPJKFHZANet\nP9H+1mimTZvGmzf5w/P0KVJS/rGOkpOTkcuzfro+duwYAwcOzGtZEvmEv/7665M7k/MT/fv3Z86c\nOSQnJzN06FCMjIzo0aMHN2/e5OnTp2kFDZRhZmaGi4sLp0+f5tmzZ4SHhxMeHk7//v25c+cOzs7O\nVKpUiRcvXhAcHMyHDx+4dOkS+/fv5/jx47Rs2RJ7e3s6dOjA2LFjv+yJS3wTjBkzhuTkZEqUyHla\nqR9++AErKyusra2ZOXNmHqiT+NJIMawS2cb69+ZsrzaZl+dCKW1tSuS1cIwbVsHd3Z0NGzbg7e2t\nbolKEULQvHlzzp49S1JSEvr6+mhpafHmzRsMDQ0ZP348vXr1YteuXSQnJ/Pjjz9SqFDqslLp0qV5\n/fp1nuTClMh/1KlThwMHDqhbRrbR0NDA3d0dd3d3pk2bxrNnz9iwYQOtWrXi/fv3zJgxAzu77Ffd\nK1myJB4eHnh4eCCXy3n//n26cqfVqlVDoVCkeVvzmuTkZI4ePcrWrVuJjIzE1tYWZ2dnmjdvnudz\nS6ifU6dOMWDAgM8quVu8eHFWr16dB6ok1IUUEpAVUkhAGh9jXR+tv8jN2UepOcSVS8N347q9N9P1\nfqJp06ZcunQJc3NzpXryCoVCwYoVK9i5cydXr16lcOHCNGzYkJcvX9KwYcM0z4+hoSFLlizB09OT\nV69ekZycTNmyZbl79y4dOnQgOjqamjVrkpCQgI2NDStWrEAmkzFq1ChSUlKYPXt2+omHZRL+kN0Y\n1oISEqAsAaiy9F1fQUiAEILffvuNe/fucfjw4Uz7FQQSExP58OFDnpRJ9ff3x83NjcTExLS42O3b\nt/Pdd9+l/f25dGdlur/fv3zLduNhlP+fA4YNqhJzM4znh+/gfn8qfrpfd3YSCbh16xYdOnRg/Pjx\neHp6qluOBFKWgPzL7YxNR7waKO3anDOfP0/FjB++sg9G6aJIiHLDqXjNjGmX7LiSrU0jJ5mXoa0w\niSxlCKLrAGY9m8XcEXPp3LYj+37YyKvN7gwYMIBx48bh6+v7yfFVRUxMDN27dycqKophw4ZRr149\n3r59y7lz5zA2NmbRokWYmZkhhODt27c8f/4cDQ2NdBVUrK2tuXz5MhEREVhZWfHu3TuaNWvG7Nmz\nGTFiBD/99BPOzs6MHj0aPT0Vlv9yVtIWqaTtlpI2gJVKbhjzlHwXlOVsVeYYi1fSltnd4RXwo5oe\nOJXvMcxATm9sz549w9vbm5CQEI4cOZJjWfmNwoUL55kH9GPIxMKFCxk+fDj+/v507NiRsLAwzMzM\nVDqXbpkSaJXUxXry95SwSL1uz3RYyv0FR2CUZLB+7dSqVQtDQ0PKli376c4SXz1SDKtEjpDJZIwc\nOZKTJ09y8OBBvL29ad68OQ4ODpw+fZrFixfnuYbY2FiGDBlC5cqVqVChAsePH+f777/H2NgYCwsL\nevToQcuWLdm/fz+XLl3i8uXLJCUlMXz4cKXjlSxZEiur1EINxYsXZ926dWkxT1WrVuWHH35g0iTl\ntdwlCjZJSUkEBQXh6OhI5cqVOXHiBKVKZczaIfEPhoaGmJubM2LECKKjo5k1axaAyo3Vj1Tq1YC/\nfP5xCNjM/IF7cwOkrARfOUIIrl69yu3bt6lXr5665UjkAyQPa35AmWdM2dJpzTxXki3i4+PZsGED\nQgjc3Nx4+fIl69ev58SJE7Ro0YJdu3Yxffp0HB2VrXHnjpiYGBo0aEC9evW4du1aljtHZTIZFSpk\nDI34FC9fvkwXs9quXTsmT56cvYPrZBIqcK0AhMB8QyQkJNCvXz/8/PwwNTVl/vz5/PDDD+qWVWA4\nefIk5ubmlC5dOq2tf//+LF++XOVzVfNqyiH7yVT/rTm6ZUtRvEoZak1qR/v27RkxYgTfffedalc/\nJPIFy5YtY8qUKfz+++8UKZIxzaPEt4cUw5oVyuIUM1tlm5qL88mtwVrzy/5bvn79GiMjIxo0aEBA\nQAAJCQk0btyYX375he7duzN//nwuXLiAv7+/ylNdzZ49m+vXr7Nhw4Y8S6P19OlTbG1tuXv3LkZG\nRnz48IFy5cpx7tw5qlT5e036dxkfUkDn3498WYU95neDdVY2U1V9pHc+P59P4OXlxdOnT9m0aZP0\nY/iZTJw4kQkTJrBp0ybc3d2xt7dn7NixdO3aVeVzTZ48mQULFnD27FksLS0RQrBx40Z27drF+fPn\n2bhxIy4uLmn9Y2NjCQwMRFNTk6ZNm6ZtopQoGDx69Ih69eoRGBgolfrNZ0hprSQKFIaGhjx69Ihi\nxYpRpEgRZsyYga+vL7Nnz2b16tUMGzaMyMhIOnbsSEhIiErnfvHiBTY2Nnma87V8+fL06tUrLYRA\nR0eHH3/8kW3btnH8+HFkMhmymaA7F/68nmcyJPKIgwcPsnXrVnx8fCRjNRcMHjyY0qVL4+joSPHi\nxdm2bRve3t7cu3dP5XONHTuWmTNn0rBhQ/bt24dMJqNr167s2LGDtWvX0rlzZzw8PFi+fDm//fYb\n5ubmTJ06FS8vLxYsWKByPRJ5hxCCX375hdGjR0vGqkQ6JINV4rOoXLkyBw8eZO/evQQHB/O///2P\n+fPnM3nyZDQ1NTl16hRFihTBy0u1VbBKlSpFRESESsdUxvjx4zlx4gQnT54EoGjRogghCA8Pp2TJ\nkjQ3T+03sODvz/mmGD16NL/++isbN25UmkxfIvuUKlWK/v37p4XLWFtb88cffzB06NA8ma9Pnz5M\nnDiR7du3s2LFCvbs2QNA8+bNuXPnDj/88APnz59HS0uLU6dOERQUhLm5OevXr88TPRJ5w6ZNm4iK\nilL5b4dEwUcKCQDolYm3TtnvWW5DApQtvSqLJM7HIQHKWLBgASdPniQ+Pp4ffviBX375hVWrVrFg\nwQJ8fX2Vlo78HHbs2MGCBQs4cyYXWRmyosI/n8+B99AjEvYbwmS7NnTt2pXOnTsDIB8hY85FMCsO\nnjX+PkAKCcjXfAz1uHfvHgYGBuqW81UQExODtbU169ato0mTJiQmJlKnTh169OjByJEjVT7fypUr\n8fb2RiaTkZSURGRkZIYcnc+fP6dFixbcvn2bTZs2MX78eBITE6lQoQI2NjY4ODhQo0YNqlevnuZh\nT0lJ4e3bt6SkpCCXy3n79i3Pnj3j2bNnmJiY0KRJE5Wfi0RGIiIisLOzY8+ePTg6OhIdHU3hwoWl\nlZB8hJTWKr+iLAclgK+KPyxltd1DgVEFxyioV68emzZtwsfHB1dXVzw8POjZsyfv3r3D3d2d33//\nXSVPzGvWrKFXr14qUPxp3HRhZilo8QpK3LjB1q1bATh79iwNZ6fmo2zduvU/B5wtwKVplX3Xc5ID\ntgBw9+5drKysJGNVhZQqVYoFCxYwcOBALl26RJEiRThy5AiVKlXCy8sLXV1dlc7XrVs3+vbtC6Su\neigUinTv37lzBzc3N/r27culS5cAuH37NuHh4Tx+/JirV69y4MABZs+ezcOHDzEyMiI5OZnXr19T\nrFgxNDU1KVSoEMWLF8fU1JSyZcty+vRppk6dSs+ePVV6LhLpSUpKomPHjnh5eeHo6MiDBw9o3rw5\nGhoaXL9+PU9yCksULCQPK2TuYc2samduDFZlnqzM5ilABmtYWBgODg40a9aMsLAwqlevzqxZsyhZ\nsiSPHz/G0dGRadOm0bt378+eIzExEX19fZ49e/ZZpfqyRYWMn89rORS5H0fRoqnZ6MeNG8fkyZNx\ndXXl2LFj/3RskIXBejaff5bKNhhm9Tg7I5+fz394+PAh7du3p3///gwYMEDdcr4qhBB0796dlJQU\nNm3ahEwmS/OyDhgwINfFBP7L4cOHWbZsGf7+/gQGBvL+/Xu0tbV5+/YtPXv2ZM6cOXh6erJjxw7G\njBnDhQsXlN4vUlJSePLkCTo6OpQpUwZNTeVf+Pv379OgQQPOnj0rxVTmEa9evaJPnz4A7N69O81z\nP3z4cMaOHUt4eHje3fMlcoQ67bevx2Btp+QHV1kydGUPaZldB5LBmiPi4+MZN24c8+bNw8TEhGrV\nqnH48GG0tLQIDg6mffv2NGrUiOXLl6OhkfPw6UePHtGiRQuVb+RKhxKDFYAn/3wW7969o3z58sTE\nxKT/HksGa74kMjKSunXrMmjQILy8vD7ruyeRNe/fv6dIkSJs3bqVjh07cu3aNQYNGsSHDx/w9/fH\nyEhZebTPR6FQUKVKFUJDQ2nQoAHJyckkJiYyZswYOnToAKQa0t7e3uzYsYM2bdrQoEED2rdvn/bg\nmRPGjBmDv78/3t7edO3aNU83fX5LCCHYtm0bXl5e9OrVi/Hjx6Ojo0P//v0pVKgQdevWZdu2bezb\nt0/dUiX+RjJYVUFuDNbMSkMC7FGxzsmZ3OjGFhwD4FPs2rWL3r17U6JECQYOHMiwYcMA0sIDbG1t\nmT9/fo7HvXr1Kr179+batWuqlpyGEIKbN29y9+5dIiIiKFSoEDo6OlhYWODi4pKWHufs2bM8evQo\n/TKhQRY/Ym++ns+3oDFy5Eiio6P5888/1S3lq6Zy5cps3boVOzs7IPVa+vXXX5HL5axYsSJP5hRC\nfNJ4vHr1KkFBQezevRuAQ4cOZepNzYzk5GT279/PtGnTKFOmDGvWrMHQ0PCzdUukrnoMGjSIsLAw\nVq9enS5vd/Hixbl37x5du3Zl4MCBaQ8hEupHSmsl8VXh4eFB586defLkCY0bN05rL168OLt37yYg\nIOCzfsBevXqVLlG5qtm9ezfVqlXDxsaGLl26MHz4cHbt2sX169cZMWIEFSpUYMuWLezevRsjIyMp\npq2AcPjwYbp166ZuGV89HTp0oGvXrty9exdI/WGbPHkyO3bsIDRUWY3g3JMdT6etrS0DBw7k0KFD\nyGQyunXrRlxcxtLVWaGlpYWHhweBgYHUrFkTBwcHXr9+/bmyVc7x48fR1NTE3t6eq1evqltOBhQK\nBbt27eLYsWN07NiR+vXr4+TkRLNmzbh+/XqGIjNVqlRh9+7d3L59G3d3dzWplshvfPMGq1sQdLyk\nbhVfH4sWLSIuLi6t7vhH9PT02LFjB6NGjSIhISFHY7579y7P6qPPnj0bDw8PQkJCKFOmTFq5yUeP\nHpGQkECfPn34888/WbhwIatWrcLZ2ZnFixfnr1hsiQz4+voSGRlJ7dq11S3lq2fWrFmMGDGCxo0b\nExQUBIC+vj52dnZ5kps1p2hqauLn58eWLVvS9OUUbW1tOnToQFRUFA8fPlSxwpwhhGDmzJmMHz+e\ndevWYWtry88//4ynpydyeWZxZl8WIQQXL16kQ4cO9O7dm86dO9OiRQtmzZrF/fv3GTZsGFpaWhmO\nGz16NAMHDsTb2zvP7vkSBY+CGRLgouSpWtlSfzZCAv4MhbB46FkOqhyDP2uDhzEYfLxG/hsS4KRk\n7upK5lG28x/AMpP2rygkIDu4ubnh7OzM6NGjs31MbGwsFhYWHDp0iDp16mT7OIVCwaVLlzh37hwP\nHjygcOHCREZG8uOPP9KmTRuEEFhbW3P79m2mT59O7969MTQ0JCYmhmvXrvH06VO2bNlCYGAgNjY2\nWFpaUuHPP1lJ6tepI2AI9M0iDdRbJd/FEu/z+Wc+Pgvv1cR8rh2Qy+UUKVKEy5cvU6tWLXXL+WZY\nuXIl/v7+7Nq1C4Dvv/8eT09POnbsqGZlqeU+fX19CQoK+uxY5vDwcFxcXOjRowdjxoxJ8/LevXuX\nwoULU7lyZVVKTseZM2eYNWsWd+7c4eXLl5iYmFC7dm1evnzJsmXLqFWrFrVr12bZsmU4OzvnmY6P\nPHz4kOHDh3P69Gnq1KmDs7MzQ4YMYc2aNdy5c4fTp08Dqdkdhg4dihAiWzHEQghWrFhBt27dPivm\nWCLvkGJYc4oyg1Unk76H/jNPJhtjEuTgeBNKFYKz78C+GKyvCtWvqthgBdiU/3/s85rVq1ezefNm\njhzJWeb9qVOn8urVKxYuXJit/sHBwbRt2xZtbW1cXFywsrIiMTGR0NBQbty4kaN8rlFRUVy7do2f\nf/4Zp9BQqgFBwDngLbBcC9pqQPFs7seQDNa85c2bN5QvX57Q0FDKlCmjbjnfDA8fPsTCwgK5XI6G\nhgYLFizA39+fQ4cOqXXDm5+fHyNGjMDf358aNWp8+oAsePHiBS1btqR27drUq1cPPz8/Ll26RJs2\nbdi+fbuKFGdkxowZXLp0iWnTpmFqakrhwoUzeCj79u1L+fLlGTNmTJ7pSElJYc6cOcyZM4dRo0bx\n448/cuTIEXr16sWiRYvw8fGhf//+2NnZYWdnJ21S+4qQDNackgcG60eSFbA/Ctrfh5/LwKTrzzE2\nNv6ng2Sw5hohBBoaGowZMyatSk52uXXrFs2bNycwMDBbnozhw4ejqanJtGnT0m6acXFxODs706NH\nD4YMGZJj/deuXaO7rS23SfWutgB8gTcyiBCwSAs8slG6XDJY855x48axc+dOrl+/nuONNhKfx7Fj\nx/D29ubmzZtAqnHj7OyMo6MjEydORE9P74vqSU5OxtPTk4CAAPbu3UujRo1UMm5MTAwLFy4kJCQE\nd3d3rKyscHd3z7N43cjISJo3b87PP/9M//79M+338OFD6tevj7+/f4aQLFXwsTknwyoAACAASURB\nVHTqvXv38PX1xdzcHEjNo2piYoKWlhZ//PEHAwcOVPncEupH2nSVj9DSAA8DOGgFK19CjRo1OH78\nuLplfVXIZDIcHBwyVKjJDrVq1WLChAm4ubmRkpJVeodUQkNDqVmzZron/PXr12NiYoK3t3eG/oEy\nWYbX0f+8Im1tcf27/zZgKBACJJHqaX0swCsZDD7A/5KgdiI4JUKvJNglB0XBsPW+CiZOnEhCQgJ3\n7txRt5RvhipVqvDmzZu0HzVNTU22b99OVFQUTZo04c2bN19Uz8mTJwkODubs2bMqM1YhtWjCx/jR\njh07YmRklOPNXNkhIiKCYcOGUbVqVVq0aEG/fv2y7F+1alVWrVpF27Zt+e2337J1n8wJc+fO5cqV\nK/j7+6cZq48fP6Zdu3bEx8ejUCiYPHkysbGxKp1XQqJguhwyq0B1NxuWQDbzYSr8/cHNDV1dXTZs\n2ICr698myjnJ2lAFfn5+ODk5kZKSwuDBg3NkvPbr149169axceNGevTokWm/gIAAjh07xvLly9O1\nh4eHY2BgkKtlqi6AHaAH6JJaFEoIMAZ0UmAmkAxUV0BzoNn165w+fZqhEycyolAhJk2aRN9spORR\nKwXEi5oVMpkMmUym8uT1EplTvnx5ChUqxL1797C0TA3aL1euHOvXr8fb25shQ4bg6+v7RbQIIbh2\n7RomJibUrJlZbWvV8d/KW7nl2rVrtGrVii5dunDjxg3KlSuXrePatWuHvr4+DRs2ZNy4cSrLrqJQ\nKJg/fz5Hjx5Nd892c3OjW7dubNmyBT8/P+bNmydtSJVQOZKHNROcnZ2pX78+o0ePlvI35gEVK1bk\n7Nmz3L59GwsLC65fv56j45cvX86QIUNwdXVl4cKFJCUlcfr0af7880/69u1L1apV6devH35+funy\nJR48eJBVq1bRrFmzz9aeCBwBngCxpG64qg6YAyeBqYAJqYbsEaABpHlHtm/fjpGREb/++isDBgxA\noVCwdu1ali5dip+fX44zJ0h8GkdHR8nD+gWRyWR89913LFu2LJ0BJ5PJGDVqFHv37iUmJibPdcTE\nxNCkSRMWLlyIra1tns9XunRp5HK5ytJdBQcH06ZNG5YvX878+fOzbax+pEGDBvTp04eRI0eqRM9H\ndHR0CA4OTtcWHR1N8+bNefDgAZMmTWLDhg2UKlUqV/Ps37+fVatWcf/+fW7dupWrsSS+EoQayPW0\nlih/SRRItm7dKkqUKCHc3NzEzJkzxa1bt7J1XGxsrDhw4IBwcHAQRYsWFfb29qJXr15i4cKF4vr1\n60Iul2c4pnz58gIQUVFRSsc8CxleR/7z8gYBiKp//1cG4gKI8yBa/d3271c7EFZWVgIQpUqVEo6O\njkJXV1eMGTNG3LhxI13fP/74I1f/lhIZGThwoJg/f766ZXxTPHnyRDg5OYnGjRuLmJiYdO/169dP\n/Pzzz3k2t0KhEMePHxfDhw8XzZs3FykpKXk2139p1aqV2LlzZ67HCQoKEmXLlhW+vr65Guft27ei\nfPny4tKlS7nW9JHz588LIyMj8fTp07S26dOni2LFigkjIyOxZs2aXM8RGhoq9PT0ROfOnQUgSpQo\nkesxJVSDmsxGIYQQBXPTVWYVhaRqQgWWmJgYDh48yLlz59i8eTN//PFHtuuQf/jwgejoaMqWLfvJ\nvm3btuX777+nd+/eSt8PzGSJ3vlf39ekpCTmzJnDqlWr0NfXp0iRIuzatYvH+vrIgUBSPawKoBDQ\nHbCpV4+DBw9StGhRtLS0SElJQVNTk/fv3zNlyhQ0NTXR1NSkX79+UgUdFdO9e3eaNGlCr1691C3l\nm0Iul9O1a1csLCyYOHFiWntsbCxOTk40atSI//u//1P5vAsXLmTJkiXY2toyevToL5qDd/LkycTF\nxTFz5szPOj4lJYW5c+cyZ84c1q5di5ubW641DR48GFNTU0aMGJHrsT4yffp0Dh8+zLFjx9Kq/6Wk\npKChoaGSTBBz587lwYMHrFixgt27dzNjxgzOnz+f63Elco+UJSCnSAbrV83du3cZOnQoN27coE+f\nPnh4eFCjRg2VxCE2bNiQYcOG0a5dO6XvZ8dgzYyrmRybCNSOj6dIkSLZ1imhGpKSkqhevTp+fn55\nsmNaImtCQkKoW7cujx49SrdEHBcXh6mpKceOHcPe3j5Xc3z48IFLly6xbt06Ll68yPPnz1mzZo1a\nKiQdPnyY6dOnc/Lkyc86vlOnTrx584ZVq1ZRsWJFlWiaMGECCQkJacVQYmNjmTJlCt7e3piYmHzW\nmHK5nKZNm9K8eXP++OMPlej8N7/99hsVKlTA29ub5ORkDAwMuH37do7DIiRUj5QlQELiX1hZWXHw\n4EECAgKIjo6mW7dumJmZMW3aNJ49e5bWTwhBYGAga9eu5f79+58c98qVK9y9excbG5u0tuTkZIQQ\neVoZpjBIxqqauHfvXlpWCokvT6VKlWjYsCF79+5N116sWDFmzpzJtGnTcjX+2rVrMTExYejQoZQt\nW5ZVq1Zx/vx5tZXzrFu3LleuXMnx/eTx48eMHj2aQ4cOsWbNGpUZqwBHjhxJi9k/ceIEdnZ27N+/\nn0mTJn32mIUKFWL9+vX4+Pjw008/8fbtW1XJ5caNGyxcuJCSJVOr/Jw/f56iRYtiamqqsjkkCib5\n2sOaXEq5x0ors9wGkof1qyU4OJiZM2eyf/9+TExMaNGiBSdOnCAhIYE6depw6tQpqlatyurVq6lU\nqZLSMapVq8aDBw+oWbMmurq6lC9fnsOHD6OhoUFcXBx169aly/nzKNuekRsPK6RuwPovltIu2jzn\n8OHDjBo1Kl/WV/9WWL9+PfPnz+fIkSPo6/9TDu748eMMGzbssz+b9+/fU716dbZs2YKTk5Oq5OYa\na2trpk6dStu2bT/Z982bN9SvX5/IyEi6devGoEGDVF4pq0aNGhQuXBgTExNu3rzJ4sWLqV+/Pubm\n5rx79y5XS/jv3r3j119/5fnz5xw6dCgtPCA3+Pj40KdPH34GUoCtQBugBjBeumeqHXV6WAtkWqvk\nTNLKZaxILPG1YGlpydq1a1EoFJw5c4bTp08zdepUWrRogYaGBnK5nAULFtCwYUMOHz6sNIXNnj17\n2LlzJ9bW1mhraxMdHc306dMpVaoU2traHD16lEGDBtGkSRN8fX1zfPNNzqS9xGecr0TuiYuLY/Dg\nwbn24knkDk9PT27fvo29vT1r166lcePGAAQGBvL8+fMcjfUxG8j+/fvZuHEjDRs2zHfe8wULFtCz\nZ08aNWqU5iXMjA0bNuDg4PBZ95vscvPmTc6cOcPjx4/ZvHkzxYsXB6BEiRI8e/YMMzOzzx67ePHi\nrF69mqZNm7J8+XKVFAsoW7Ys5oApcAqoQ6qxKiFRID2smaEVIz19fets2rSJIUOGsGvXrs/yukRF\nRaGvr09sbCwlSuTM1LyQiYc1s1EkD2vekZycjLa2No6OjtJmjXzC/v376du3L127dmX69OncunWL\nOnXqkJiY+Mn4dLlczsSJE1m0aBGWlpa0adOGzp07U6VKlS+kPmf0798fDQ0Nli5dmmW/iRMnkpSU\nxNSpU7+Qsn9wcXHh999/p1WrVrkea8+ePSxdupSAgIBcj5WQkIBR0aLEk2qoagEtSY1fnCbdM9WO\nFMMqIaEiunTpgo+PD+3atWPgwIEsXLgwRzleY2NjKVKkyBdJNB8hk2V4SagGLS0tNm3axOPHj3n3\n7p265UgA7u7u3Lhxg4sXL9K3b1+sra2xtbXlzJkzmR5z7tw5hg4dStOmTTlz5gx37tzh3LlzjB07\nNt8aqwBjxozBz8/vk4UELC0tM+Q0/VI0adKEY8eOqWSshg0bcu7cOZKTM1tnyj5FihTB8+//rwxc\nB2YDC0jdjBUWFpbrOSQKJvnaw/pWN/Mf8Hxfh11CrURERLBu3ToiIiLYs2cPFhYWuLq60r59e44c\nOYKrqyu1atXKcNzcuXMJCQn5pGckJ7xUYohmVizRVPIgqJTvv/+eNm3a0LdvX3VLkfibuLg4WrZs\niYODAwYGBkRGRjJ//nwg9YHR19eXly9fcuvWLS5cuEDv3r2pWrUqXbt2RVOz4ESxVa5cmQMHDlC9\nevVM+1y5coWff/5ZLTHWly5dokOHDty8eZNSpUrh5+dHUFAQderUoWfPnjkez8bGhv/7v/+jXr16\nuda2c+dOBg4cyMuXLzE2NqZcuXKkpKRgaGhImTJlWLt2ba7nkPg8pBhWCQkVY2pqyujRo4HUnIFn\nzpxh9+7dtGvXjkePHgEoveiio6OlPKhfGZ+KI5T4shQrVowDBw7g4uLCrVu3qFevHteuXeP8+fPM\nmTMHS0tLHB0d6dSpE+vWrct1xSR1YWdnx+XLl7M0WAsVKkR8fPwXVPUPDg4ONGzYkBUrVjBy5EgW\nLlxIoUKFOH/+PDVr1sTOzi5HpaNLly6tsnMpWrQor169oly5cowfP57g4GBmz56d9v6sWbMoWbIk\nhQsXVsl8EgUDKSRA4qunRIkSuLm5sXLlSh4+fMiJEyeA1PKu/yUqKoq4uDiVzX3hwgWakbqcJUgt\n3araauMSWVGrVi0phjUfUqpUKU6dOkXr1q0JCgrC09OT48eP8/vvv7Nv3z7Gjh1Lly5dlBqriYmJ\nPHz4MN/Xqndzc2PFihVZ6ly3bh3ff//9F1SVHn19fXR0dAB4/fo1P/zwAzKZDE9PT+zs7Lhw4UK2\nx3rz5g0GBgYq0dWyZUsGDBjAkydPqF27NrNmzSIxMZGJEydSpEgRjI2NKV++vMrK4EoUEL5USa1/\nk91pY3XI9CUhkRsePXokzMzMxMaNG9O1165dWxw9ejTX4ycmJoqpU6cKPT090QBEBRDz/i6/Oh3E\nFRDhSl7iu/SvKyh/SWSP+/fvi6JFiwo/Pz+RnJysbjkS/+Lhw4cCEL/++muW/RQKhYiKihK3bt0S\nCxcuFDVq1BCACAoK+kJKP4+UlBRRs2ZNsWXLlkz7VK1aVdy8efMLqkrP8OHDxeDBg8X9+/dFmTJl\nRFJSkhAi9d/cx8dHWFhYZKu0bUJCgjAwMBBhYWG51qRQKIQQQiQlJYnVq1eLKlWqiK5duwofHx8R\nHx8vAgMDBSDq1q0rpkyZkuv5JHKGmsxGIYQQ+drD+uhD5i8JidzwMb6se/fuaZ7WpKQkbty48dmb\nOQ4fPkznzp0xMDCgUqVKBAYGMnPmTM4CWlWqsKRSJWrUqMFma2uaacLhGmDa8l+v77Ie/yUQ+1nK\nvl0sLCzYtGkTkydPxtXVlcePH6tbksTfGBkZ0ahRIw4dOkRCQkK69xITEzl06BD9+vWjevXqVKxY\nEQ8PD27fvs3kyZOpVq1ajrN4fGkKFSqEj48PgwYNynTVJiUlRa1FRdq3b8++ffvw9fWlS5cuaGml\nJoeUyWT06tULQ0ND1q9fn+UYDx48oEGDBrRo0eKzK2d9JCQkBA0NDZKSktDS0qJXr14kJydz7do1\nfvnlF+bNm8f//vc/AC5evMjx48dzNZ9EAUMdVnJ2p83MuyR5mCRUhZubmyhWrFja3+7u7mL27Nk5\nHmf9+vXCzMxMTJ8+XRw9elTcvHlTKBQKoVAoBCD09PREUFCQkMvlQggh7jVI9bZ2KYvoboJoro/4\noyoi2R3xtDnCUQ9hXQJR7G+v7MeXtfT9/yzkcrmYMmWKMDMzkzyt+QxPT0/Rtm1bMX36dNGyZUth\nZGQktLS0hKOjo5g7d64ICgpK87oJIURycrLQ19cX4eHhalSdfb777juxevXqDO0vXrwQZcuWFX/9\n9ZcaVKUybtw44e3tLcqXLy+uXbuW4f0rV64IAwMD0bt3b7F3795078XExIjffvtN6Ovri0WLFqX7\njD6XLl26CEBcuXIlrc3S0jLdPfDfr/+ukEnkPWoyG1PnVsukksEqkU9YsmSJcHZ2FhYWFuLOnTvi\nzp07wtDQUBw6dChH47Rp00Zs3749Q/v+/fsFIKpUqZK23CaEEKIl4qAdYoUVYk1NxAFbREtDhJ4W\nQl8b0ckEcb0xYup/btANpO9/rnBwcBC//vqr+PDhg7qlSPzNhw8fxIQJE4S3t7fYsWOHCA8Pz3IZ\net++fcLR0VElBtKXYOfOnaJJkybp2mJiYoSJiYlwcXERly5dUpMyIaZMmSKcnJxErVq1Mu1z+/Zt\nsWzZMlGpUiUxcuRIIUSqfgcHB9GtWzfx8uVLlWhZs2aNAISNjY149+5dWvuOHTtEx44dhZaWlgBE\n165dRa1atUTjxo3T31MlvgiSwZoJ/pDpS0JCFcTExKQ91evr6wtfX18RGBgoGjduLHbs2CGESI2p\nCgwMFA0bNhSAsLe3F+vXr0/7wQwICBCAePDgQYbxZ8+eLTw8PISurq6IjIzMUotCoRAvXrwQjx49\nSte+Z88e0a1bN2FoaCj69eunojP/Nnn16pVo27at6NChQ4ExeCTS4+XlJaZPn65uGdkmISFBlCxZ\nMp1hN2PGDFGnTh0xZcoUoa2tLdasWaMWbX379hXOzs7CyclJ3L59O8u+vr6+wtXVVTx+/FjY29uL\nAQMGqOwaWrlypQCEoaGhSExMTPeeXC4XGhoaon///sLc3FzEx8eL+Ph44eLiIiwtLcWFCxdUokEi\ne0gGayZIBqvEl0ChUIiEhARx6tQpYWxsLA4cOCAOHTqUugRvbS0cHR1FtWrVhJubmzA1NRUeHh7C\nxsZGdOnSRcycOVMYGxuLdevWpS33/5vo6GhRrVo10aNHj1zrjI6OFrq6uuLSpUuSsZULEhMThYOD\ng9qMBInPJzk5WZQrV07cuHFD3VJyROfOncXEiROFEEJMnz5dAKJJkyaiYcOGYsGCBcLIyEgtnlZT\nU1Nx9OhR4eHhIUqVKiVevXqVad9FixYJW1tbYWRkJGbNmqX0fvc5zJ07VwDCwMAg0/nPnj0r9PT0\nxLRp04RcLhetW7cWOjo6ol27dsLd3V0lOiSyh2SwZoJksEp8aXr27CkAYWpqKgYOHCgmTZok/Pz8\nROvWrUXFihXTdqgOGTJELFiwQHTo0EEEBwd/UX3m5uaiUaNG6eK8JHKGj4+P6Nq1q7plSOSQbdu2\niQYNGqhbRo4JDQ0V5ubm4vfffxd6enri6dOnIjY2VhQrVkzEx8eLtWvXCmdn5y+qKSoqSgBi7Nix\nQqFQiD59+mS56/79+/fixx9/VGl2huTkZAGINm3aZPCsfkShUIgNGzYICwsLUbRoUWFlZSUaNGgg\nHj9+LO7fvy8qVKigMj0Sn0YyWDNBMlglvjRyuVxs3rxZeHl5CS8vL6GtrS2aNGkiALF//34hhBA/\n//yz6Nu3r9o0pqSkiBUrVghDQ0Nx8OBBtekoyISEhAh9fX1x//59dUuRyCYKhULUrVtXaax4QSA0\nNFQ4Ozunpbl68eKFMDAwEEKkGm6VKlUSgYGBX0zPxw2hJUqUEMOGDRMHDhwQrq6uX2z+j7x+/TrL\n97dt2ybMzc3FuHHjhJmZmRg4cGCadzclJUXo6OiI+Pj4LyFVQkgGq4REvsXPzy9tw9Pu3buFEEIE\nBwcLGxsbcfjwYbVqCwwMFBoaGtLN+jPx8fERxsbGYu3atdImrAJAQECAsLS0VNlStLpRKBRCV1dX\nxMXFCSGE6Nixo1i3bt0X1fDu3TsRGRkpqlWrJtavXy+MjIyEn5/fF9vMdPnyZXHx4sUs+9jZ2QlA\nWFhYiG7duglAvHjxQpw9e1bY29uLMmXKZIj7l8g7JINVQiIfEx4eLmJiYtK1nT17VhgYGKh9l2rN\nmjVFuXLlxN27d9Wqo6DycYOdpaWliIiIyJM5kpOTRUhIiLh37564e/euCA0NzZN5vmbkcrmoW7fu\nV5fGyMbGRhw8eFAkJSUJJycnsW/fPrXo2Lhxo6hatarYv3+/sLKyEmXKlBFXr17N0zmvX78u9PX1\nhbGxsdiwYUOm/VatWpWWNuvUqVOievXqQgghGjRoIDQ0NISxsbGUqu4Lok77TVPFaV0lJL46TE1N\nM7Q5Ozujr69PcHAw1tbWalCVyq1bt1i3bh0uLi5069aN4sWLc+fOHeLi4nBwcMDV1RVnZ2c0NaVL\nXRn169fn5MmTzJgxAycnJ1asWEGrVq1yPe6TJ09Yv349/v7+XL16FSMjI3R1dZHJZLx69Yrx48cz\nePBgFZzBt4Gvry8ymYwff/xR3VJUyvjx4+nevTtyuZx69erh4uKiFh1dunTh7NmznD17ljt37rBj\nxw5cXV1ZunQpXbp0Ufl84eHh/Pjjj8yaNYvy5cvz888/U7FiRZycnJDJZOn69u7dO+3/nZ2dSUpK\nomfPnlSqVInw8HCGDRsm3d++EWR/W8xfdlKZLN/XgZaQyIoPHz5QsmRJQkNDc13dRRVcvnyZo0eP\ncv/+fRwcHDA1NeX8+fMcPnyY8PBwfvrpJyZPnpxWyUYiIwEBAfzyyy94eHgwe/ZsChUqlKFPbGws\nDx8+JDo6mrJly1KjRo10P7ABAQHMmDGDmzdv0qlTJ9q3b4+Tk1O6akZPnz7FxcWFxYsX4+bm9kXO\nrSATHR1NjRo12LNnDw4ODuqWo1KEEBw9ehQrKyulD8aqIjk5mZYtW1KoUCHKlSuHubk5Li4uNGrU\nKK3PpEmTeP/+PdOnTwfgzp07NGvWjE2bNtGkSROVaTlx4gQ9e/Zk0KBBDBs2DIVCwaJFi1i2bBkA\nmzZtwt7ePtPjY2JiWLRoEYsXL8bX15fWrVurTJvEp1Gr/aYOt66appWQUCkVKlQoEEvxISEhwt3d\nXTRs2FC8ePFC3XLyNVFRUaJp06aiWbNmaWmGXr16JRYsWCAsLS1F0aJFRe3atUWzZs2Eubm50NXV\nFZaWlmLgwIGibt26wtzcXGzevPmTMbGbNm0S9vb2aTuj5XK5uH37tpgyZUqWqYW+NRQKhfjuu+/E\n4MGD1S2lQKNQKIS+vr4YM2aMWLlypRg9erQoU6aM8Pf3T+uze/duYWFhkW55fc+ePWlL8LlFLpeL\nYcOGiXLlyikNfWjVqpUAxIkTJ7I1npTaTz2o036TDFYJic9AoVCISpUqFZjUUnK5XHh5eQlnZ2cp\nhvITJCUlifnz54ty5cqJ8uXLi5IlSwpPT09x+vTpdD+SCoVCxMbGivPnz4tp06aJI0eOZDumWS6X\ni++//17UrVtXNGvWTGhoaAhTU9PUcr1duuTVqRU4ZsyYIRwdHTNNeSSRfQICAoShoaG4fPmyEEKI\noKAgYWBgkM5ArFq1aroCAgqFQmhra4uEhIRczz916lRhZ2cnoqKiMrz3f//3f6JUqVJi586duZ5H\nIm+RDFYJiQKGQqEQVatWFZMmTVK3lGwTFxcnfvrpJ9GmTRt1SykQJCUliXv37om3b9/myfgJCQli\n9+7dYvv27SI2NlbI5XIRGxsrzMzMRKdOnb55D9L27duFqampePr0qbqlfDXs2rVLGBoaij179ggh\nhDh27JgwNDQUvr6+QojUPM/9+/dP++7FxMQIXV3dLEvlZoft27cLExMTpRsbT5w4IczMzMS1a9dy\nNYfEl0EyWCUkChjnz58XZcqUybOd5XnFxzKR0rJz/uX9+/fC3t5ebN26Vd1S1Mb58+eFoaFhgVnB\nKEgEBQUJMzOztGX527dvCwsLCzFixAgRExMjqlWrlpbfWS6XixIlSog3b97keJ6wsDDh4+MjvLy8\nhKGhYaZZB1q2bCnmzZv3+Sck8UVRp/2moZ7IWQmJgsugQYNo1qwZS5YsyRcbrnKCrq4uLVq0YO/e\nveqWIpEJOjo61KxZk+joaHVLUQu7d+/G3d2dNWvWYGtrq245Xx1OTk6sXbuWgQMHEhERQY0aNTh/\n/jz79+9nzZo1DBw4EG9vby5fvoyGhgZOTk5MmDAhR3M8fPiQWrVqcfz4ccqVK8eJEyeoU6eO0r5e\nXl7MnTuXFy9eqODsJL5q1GElq2laCYlcExYWJoBPVmfJz2zevFnY2tqKO3fuqFuKhBI+hpvkdR7M\n/EZCQoIYOnSoKFeu3CeTyUvknhkzZggjI6O0Cn5PnjwR5cqVS8uNqqenJyZOnChiYmKEgYGBePjw\nYZbj3b17V3Tu3Fls3bpV9O/fX/Ts2TPbWho2bCj8/Pyy3T85OTlDwZT69euL33777ZsPpclr1Gm/\nSR5WCYlskpKSQseOHZHJZOjr66tbzmfj4eGBk5MTAwYMULcUCSXcu3ePd+/eYWNjo24pX4ygoCBs\nbGwIDw/nypUrX136qvzIyJEj2bp1K3369OH9+/eUL1+ehw8fEhwczMqVK9HQ0GDLli38+OOPaGtr\nZ7kq8/btW1q3bo2ZmRnLly8nLi6OYcOGZVtL48aNmTBhAo8ePUpr8/PzY+PGjSQkJKTrGxISgp6e\nXrrrIywsjKCgIJYvX/7V5eqV+AfJYJWQyAYBAQFYWFigo6NDcHBwhuTWBYnChQujr69P8eLF1S1F\nQgkxMTHExcUV6O9YdhFCMH/+fNq3b8+0adPYsmULhoaG6pb1zdC4cWPq1KnD1q1bgdR7g6GhIe7u\n7hw6dIjXr19ja2tLr169+O677zId5/r165QtW5ZZs2Zx/PhxfH19qVGjRrY0pKSkEBUVRVhYGA0b\nNkQulzNixAjGjRvH+vXr0dPTo3Tp0lhaWlKtWjWsrKyIi4ujb9++aWPcunWLJk2aEBwczLlz53L3\njyKRb5HKQ0hIZIOjR4/SunVrlixZUuANiV27drFkyRLOnDmDEKLAn8/XRkREBA0aNFC3jDxHCMGv\nv/7KuXPnuHDhAhUqVFC3pG+Sbt26sXz5cjw9PdNVjLK3t2fs2LFMmTKFhQsXUqVKlUzHePnyJcbG\nxp81/9atW7l8+TJ79+6lYsWK7Nu3j9mzZ6Ovr09YWBhJSUkkJSURHR1N6dKlSUlJoWfPngwdOjRt\nDGtra06cOEHTpk0pVarUZ+mQyP9IBquExCeIjY1l9erVBAYGfhXGXa1aemgazwAAIABJREFUtTA3\nN6d+/fokJCRQrlw5AgICkMvlVK1a9as4x4LMtm3baNeunbpl5Dnjxo3j4sWLnD59WvL2qxEPDw/W\nrFlDly5d2LBhA9ra2mnveXl58eDBAwIDA+ncuXOmY5ibmxMaGvpZ81+/fp3vv/8+rZrWL7/8AsDx\n48epVKkSOjo6FCpUiBs3bhASEoKdnV2GhxszMzMAQkNDefDgQYY54uLimDdvHjKZjCZNmlC3bt10\n5ylRMJBCAiQkPsGkSZOoW7cu1apVU7cUlVClShWuXLlCTEwMb9++JSkpCTs7O1xcXKhcuTLe3t5c\nvXpVKp+sBkJDQzlx4gSdOnVSt5Q8Zf369WzZsgV/f3/JWFUzOjo67N27l7dv3zJt2rQM7585c4bI\nyEhevXqV6RjW1tZERUVx9uzZHM0dGxvL1atX0x6SY2NjOXPmDKampuzZs4d169bRrFkzNDQ0qFOn\nDseOHcvUE7948WIAnj9/ntY2cuRIrKysMDQ05N69e8TGxuLl5YWJiQlPnjzJkVaJfIA6dnqpaVoJ\nic+idOnSaWU6v0Z69+4tdu/eLRQKhbhx44YYO3asMDc3FzVq1BC7d+9Wt7xvioiICKGvry/CwsLU\nLSXPiIiIEBUqVBBnz55VtxSJfxEaGipKly4tnj9/nq79zZs3on379qJy5cpZHr9582ZhYmKS7dKq\nMTExomLFiqJfv35pO/63bNki3NzcxPPnz4WlpaWwsbERQ4cOFcWLFxcDBgzIMnuEQqEQ169fT5cl\nYOzYsaJUqVLil19+SVfhr0+fPmLhwoXpjpfL5SIgIEBs3bpVXLhwQco2kAnqtN9kfwv4oshkMsl7\nI1EgEEJQrlw55s6dS5UqVbCzs1O3pC+CQqHg2LFjeHp6cvz4cQwMDDA0NKRQoULqlvbVM3bsWO7d\nu8fWrVu/uvCMrVu30r9/f4YOHcro0aPVLUfiPwwZMoS4uDj+/PPPdO2JiYmUL1+eNWvW0LhxY4oW\nLar0+IULF3Lx4kU2btz4yblu375Np06duHPnTlrbsGHDKFWqFGPGjCE+Pp5ly5Yxffp0Fi9ejKen\nZ47PRwjBmTNn2LZtG0FBQXTp0gUhBH/99RdyuTztPAMCAhgxYgQymYwqVapw48YNdHV1GTZsGB07\ndkRH5//Zu++oKo73f+DvSweRjhTpiIgoYBBFRRDBgr0gtgTsDXtNYhJ7b4kxKoldsYIVGwJ2jAJR\nQJCudOlN+uXO74/8wvn4tVEu7AWf1zk5J7K7M+9FhIfZ2RmZevfdWnFav3FRJXPULSH1xufzmYWF\nBQPwVX7dysjIMABMUVGRqaqqMjc3N/bXX3+xe/fusdu3bzN/f39WVVXFdcxWpbS0lHXt2pX9/vvv\nXEcRqr179zJFRcXaveyJ6CkoKGA6OjosICDgg2MXL15knTt3ZhoaGuzNmzcfvf7u3busZ8+e7N27\nd6yoqOizo5QPHz5kPXr0qP1zRkYGU1NTY69evWKZmZls5cqVzNvbm+nq6rIbN258MbtAIGB+fn7s\n0KFD7Pjx4yw0NLR2S9m8vLza7+EAmKqqKrOwsGDl5eXM39+faWpqsosXL9bmrampYbdu3WKOjo5M\nRkaGWVhYsO+//54lJiZ+McenjB07likrK7MnT540uA1RwOXPQSpYCfmC6upq1rFjRzZhwgSuozS7\ne/fusZcvXzKBQMBSU1PZ0aNH2cSJE1nfvn2Zs7Mz69GjB+vYsWPtlAIiHImJiUxTU5OdPXuW6yhC\nIRAImJWVFbtz5w7XUcgXBAQEsHbt2rGIiIiPHv/pp5+Yh4fHR49lZWWxrl27MmlpaSYlJcXc3NzY\nokWLmLW1Nfv1119ZcHAwq6mpYYwxtm3bNjZ//vzaa6Ojo5mmpiaLjo5mM2fOZOLi4mzHjh0sKCiI\nqamp1V73v3744QempqbGevXqxXx8fBgAxuPx3itM8/PzWVxcHGvfvj0bMmQIq6ioYPn5+WzMmDHM\n0NCQaWhofHYaQ0VFBXvy5AlbunQpU1NTY87OzuzIkSP13t56w4YNtbmmTJlSW0y3NFSwEiKiysrK\n2MSJE1m3bt1YQUEB13FE0s2bN5m5uTkbNmwYKy4u5jpOqxEREcE0NTXZ6dOnuY7SaCdOnGDW1tas\nsrKS6yikDs6ePct0dXVZXl7eB8d++eUXtnr16s9eLxAIWHFxMXN3d2cbNmxgN27cYOPHj2cGBgZs\n/vz5bPLkyQwAO3/+/HvXLV26lMnKyrK5c+cySUlJNmTIEBYTE8Pk5OQ+Oq973bp1tUXgkSNHmLu7\nO/vuu+/YnDlz2IABAxgAtmzZMsbYv3NxJSUlmaurK9PS0mLnz59nAQEBzMfHp86fl/LycnbmzBnm\n6urKFBUVmZ2dXb2uLy4uZkOHDmUAmJOTE6uurq7ztaKCy/qN5rAS8hn/LYR96dIlyMrKch1HZPH5\nfHh6eiIkJATXr1+HlpYW15FahaioKDg6OsLPzw89evTgOk6DvHv3Ds7Ozli4cCEmTZrEdRxSR8uW\nLcODBw9w5MgRdO3atfbjI0aMwJQpUzBmzJh6t5mUlIQpU6Zg5MiRSE5OxsaNG6GgoPDeOcbGxujR\nowcMDAwQHh6OTp06ISAgAImJiSgoKPhgOarU1FTIycl9dPfBiIgISElJoVOnTgD+XTUgNja2dh61\ntrY2BgwYABsbGzg4ONRrznhFRQXu3LmDxYsXY+jQodizZ88X5/ifPn0akydPhpKSEgoLC3H//n3Y\n29vXuU9RQHNYCRFRjo6OrWKEqzkIBAK2efNmpqGhwVasWMEOHz7MwsLCuI7V4l26dInp6em12BH+\nFStWsLFjx9Jc5xampqaGeXl5MTU1NbZu3braKT/29vbs3r17TdJneXk5U1RUZIsXL2arVq1iR48e\nfW/u6dWrV4XWV1VVFTt06BBbunQpMzMzY127dmWHDx9m5eXl9WqnsLCQOTk5sSVLlnz0eHl5OVu9\nejXLzMxko0aNYra2tkxJSYnxeDyWkpLCysvL2cGDB1vMlCou6zcqWAn5jFOnTrE+ffpwHaNFiY6O\nZj/99BObOnUqa9euHS1fJASurq7swIEDXMeot+rqaqapqcliYmK4jkIaKCMjg9nY2LCFCxeytLQ0\n1q1bN3bv3r0mKbAiIyOZgoICc3R0ZEuXLmWMMfbkyRMmKyvL1NXV2cOHD4XST05ODnNycqpdwksg\nELDbt28zFxcXpq+vzwoLC+vVXnx8PFNWVmZ79uxhWVlZjLF/Bzv69+/PXFxcGACmqanJpKSkmJGR\nEVu5cmXttdHR0QwAW79+vVDuralRwUqIiDp16hTT1dVtsRPkubZp0ya2aNEirmO0eBcuXGC2trZc\nx6i3Z8+eMXNzc65jkEYqKChgTk5OTFtbmykpKTExMTEmLi7OOnXqxP76668P5ro2dG6mQCBgFy5c\nYDt27HhvJYKgoCAWHBzcqHv4XytXrmQA2OHDhz84Nn36dLZ8+fJ6t3n37l327bffMkVFRTZnzhxm\nZmbGALAOHTqwOXPmME1NTaasrMx8fHzee9ogEAiYtLQ0A9AinuZxWb/RHFZCPiMuLg4zZ85EVVUV\nLl68SHMz6+nZs2eYNm0aXr58yXWUFi09PR0dO3bE48ePYWVlxXWcOsvIyIC5uTnWrl2LRYsWcR2H\nNFJAQACGDBmCa9euoX///nj8+DG2bduG0NBQLFmyBG3btkVOTg62bNmCnj17wsnJCUOGDEHPnj25\njv6evn374tGjR3j16lXt/Nb/JCYmws7ODhkZGfVeBzk+Ph4LFy7ErVu3kJCQgLy8PJw+fRopKSno\n378/pk6d+tE1bG/evIkhQ4YAAMaNGwczMzOsXbtWJNdh5rJ+o4KVkC9gjGHjxo04fPgwIiMjaSvJ\neqipqUG7du0QGRkJbW1truO0aAcPHsT169dx7do1rqPUS0hICHr06AGBQCCSP4BJ3QUEBODQoUO4\nffs2Nm/ejP79+6Njx444fvw4YmJiUFJSgnfv3mH69OmorKxEUFAQDh06hCtXrqBXr14i8/efkJCA\n2NhYDB069INj6enpMDU1RXh4OIyNjb/YlkAggJjYv7vch4eHw87ODlevXkW/fv3qdb/nz5/H3bt3\nUVVVhSNHjsDGxgYrV65EaWkpvvvuu9o+uEYFKyEtgKurK+zs7LB48WKuo7QYAoEAampqiIqKotHp\nRiooKIC+vj4yMzM/udOQKDp06BBu3boFHx8frqMQIfH398exY8cQHByMd+/eoby8HHJycrC1tYWr\nqysGDhxY++/97NmzmDVrFng8HgwMDDBhwgR4enq+tzpAamoqtLW1P/mWfVVVFbKysqCrq9uk91VU\nVARTU1O4ublhxowZ0NfXh6Ki4gfnVVRUwMvLC8eOHUNcXByGDRuGb7/9FsOGDcPUqVMRGhoKHx+f\nD0Zv68rX1xexsbFISEjA0aNHISsri7CwMJiZmTX2FhuNVgkgpAV4+fIlU1dXpzff6yE8PJyZmJhw\nHaPVcHZ2Zr6+vlzHqJd9+/axmTNnch2DNJH09HSWm5vLMjIymLe3Nxs1ahRTVVVl6urqzNLSktnZ\n2bGePXuy4cOHs/Xr1zNdXV2mpKTEIiMjWUpKCps9ezaTlZVl8vLyzNzcnEVFRTFnZ2emqqrKRo4c\nyWbNmsUAMGlp6do+y8vL2bFjx9jcuXPZoEGDhLaCwI4dO1jbtm2Zk5MTU1RUZG3atGF9+vSpfTnr\nP5MnT2bOzs7s7t277O3bt8zLy4vp6emx69evM4FAwPbt28fMzMwave6wQCBg+/fvr10lYcmSJczS\n0pLTFQW4rN9ohJWQeti3bx8CAgJw+fJlrqO0CNu3b0dsbCwOHz7MdZQWr7y8HCYmJjh16hT69evH\ndZw6e/ToESZPnozY2Fjak/0rwRhDRkYGsrOzUVxcDElJSYSGhsLb2xsKCgpwc3ODu7s7Bg0ahPv3\n77937dGjRzF16tQP2jx69Cjc3d1x7do1bN26FfHx8QCAvLw8qKmpISQkBAYGBo3KnZKSgvDwcEhJ\nScHCwgLq6urYtGkTTp48iQULFqBt27bo3LkzFi9eDA8PD8ydO7f22ps3b2L+/PkIDg5Gu3btMGzY\nMMjJycHDwwPDhg1rVK7s7OzadV43bdqEa9euNbrNhqIpAYS0EOXl5TA0NERAQAC6dOnCdRyR9s8/\n/2DQoEG4e/cufa6EoLi4GDo6OkhKSoKamhrXcepl8ODBmDhxIjw8PLiOQkTI3bt3ERsbi86dO8PC\nwgKrVq3C5cuXIS0tjXXr1sHKygodOnSofW9g3bp1uHDhAiorK5GdnY3Jkydj1qxZTf4i4pkzZ/Dk\nyRPk5+fXfl/7/vvvoaGh8d5569evh7e3N/z8/NCmTRssXrwY/v7+WLx4MdauXdvg/m/cuPHBfNua\nmhpO5rXSlABCWpDNmzezb7/9lusYIi0sLIxpaGiwS5cucR2lVfHw8GBbtmzhOka9+fn5MRMTE5aY\nmMh1FCLCfHx8mKGhYe0j8Ozs7Npj+fn5DAAbPnw4MzQ0ZG/fvuUw6acdPHiQtWvXjs2cOZNt2LCh\n9l48PT0bvPlHZmYmGz58eO1Ws927d2c5OTlCTl43XNZvNMJKSD0VFRWhQ4cOePDggUhMghc1ISEh\nGDZsGLy8vDBq1Ciu47Qq/23V+uzZs0Y//mxOjDF8//33SElJwZkzZ7iOQ0RYVVUVrl27hqSkJCxd\nurT2RazKykosWLAAxsbGmDRpUpO/gPUx0dHR2LlzJ7y8vCApKfnJ8yIjI/Ho0SPExcUhJiYGt27d\nAgDExMTA1NS0URlCQ0NhY2MD4N/vB507d25Ue/VFUwIIaWG2b9+Op0+fwtfXl+soIqWkpASmpqY4\ncOAARo4cyXWcVmnUqFFwd3dv0F7uXAoKCsKGDRtw9+5drqOQr8yJEyfg4eGB7t27Y8+ePbCzs2tQ\nO05OTggLC8OBAwcwceLEOl8XHx8PJSUlqKurN6jf/8vHxwcPHjzAkiVLYGhoKJQ264oKVkJamPLy\ncnTs2BE+Pj4ityg2l3788Uekp6fj+PHjXEdptRwcHLBmzRr079+f6yj1cuPGDezdu7d2tKk++Hw+\ngoODYW9v3wTJSGtWUlKCtLQ09OvXD2VlZZCTk0P//v1hZmaGnj17IisrC//88w/69OmDcePGfbKd\ngoICGBgYYObMmYiIiIC/v38z3oXo4LJ+E42VaAlpYWRlZbF27VosXrwYNTU1XMcRCQ8ePMChQ4ew\nadMmrqO0WkFBQYiLi0O3bt24jlJv2traSE9Pb9C1O3bsgIODA63lSurFy8sL+vr6sLKygpiYGLZs\n2YKXL1+if//+yM7OxooVK3Dy5Enk5uZi4sSJSElJ+WRb3t7eGDRoENasWYOgoCAkJyc3450QgEZY\nCWkwgUAABwcHuLm5YcGCBVzH4ZSPjw/mzJmDM2fOYMCAAVzHaZViY2Nhb2+PM2fOtLjRVQDIycmB\nmZkZcnNz633t4MGDYW1tjUOHDiE5OZmWxyJfdP/+fYwdOxaPHz+Gvr4+pKWlP7rzVFVVFU6cOIHN\nmzcjNjb2k3NT58+fDz6fj4MHD2LDhg1ISEj4Kp8k0QgrIS2QmJgYvLy8sGjRIlRWVnIdhzN8Ph8L\nFy6El5cXFatNaNGiRfjxxx9bZLEKANLS0qiqqqrXNYmJifD09ERUVBQWL14MCwsLXL9+vYkSkpao\ntLQUwcHB4PP5tR+rqKjApEmT4OTkBFNTU8jIyHxym9QTJ05g5syZuHz58mdfpIqOjq6dN963b18k\nJSUJ90bIF1HBSkgjaGlpQU5ODsXFxVxH4czbt2/BGMPYsWO5jtJqpaSkIDQ0FHPmzOE6SoNJSUmh\noqICAoGgztdERkZi//79uHLlCtTV1dGrVy9ERkY2YUrSUiQnJ6NXr15QVlbGlClTYGNjg1evXgEA\ndu7ciW+++Qbe3t5fbOfBgwf48ccfYWFh8dnzxMXFUV1dDQDQ0dFBbGwsNm7ciCtXrjT+ZkidUMFK\nSCMoKyvX7nhSXl7OdRxOyMrK1nvkjNSPv78/XFxcIC0tzXWUBpORkYGenh6eP39ep/Pj4uLw008/\n4ffff8c333wDADAyMqotSsjXKzw8HBYWFkhKSkJgYCBiY2MxZ84cDBw4EFu3boWXlxfc3d0hISHx\n2Xbu3buHgIAAfP/991/sU1ZWFu/evQMAdOjQAYMHD0ZISAhmzJghckVrYmIitLS0EBERwXUUoaKC\nlZBG2rFjB6SkpNCrVy+EhoZyHafZycnJoaysjOsYrdqdO3fQp08frmM0SlpaGvLz82FiYvLFc+/d\nu4c+ffpg/vz58PT0BPDvWq7Z2dm4cuUKioqKmjouEVGJiYkYM2YM9uzZg6ysLPTt2xc8Hg+zZ8/G\n3r17UVBQgP79+0NZWfmLbWloaKC6uhp79+794oCDgoLCe/OvT5w4gStXruCPP/4QuRdN8/Ly8Pbt\nW6xcuZLrKMLVzBsVsP//khdLT0/nomtCmoRAIGDHjh1jmpqabN68eSw/P5/rSM3m3bt3TEpKiusY\nrda7d+9Y27ZtW/zX1LRp09iKFSs+eozP57Po6Gi2f/9+1q1bN6asrMyCgoLeOycsLKx21yA+n98c\nkYmISUlJYUOHDmW//PKLUNscOXIkU1NTY9u2bWPV1dUfnFNRUcGUlJQ+WreEhoayzp07Cy2PMJSU\nlDAVFRUmLS3NKioqhNo2R2UjY4wxzkZYo6KiuOqaEKHj8Xjw8PBAVFQUBAIBzMzMsGfPnq9imkB0\ndHSz77byNYmNjYWhoWGdRoxEUU1NDTZt2lQ7V/BjvL290bdvX9y/fx/bt29HUlISHB0d3zunY8eO\nAICAgIDa3Y/I1yE3Nxeurq7o2bMngoKChLopia6uLi5fvozg4GDcunULLi4ueP369XvnBAcHw9jY\nGNra2h9cb2VlhczMTGRkZAgtU2PJy8vj8ePHqKyshJOTk1DaTE5OhoeHh1DaaijOCtb/+82IkNZA\nRUUFBw4cgL+/Px4+fAhjY2Ps3bu3Va8iEBER8cUXFkjDlZSUQF5enusY9Zabm4stW7bAyMgIt27d\nwv3796GkpPTRc8+fP4+NGzfi7NmzcHZ2/uh5gYGB0NHRabGrJJCG4fP5GDduHB4/fozMzEwcP368\ndk6zMJmYmMDf3x/dunVDt27dMHz48NqXrLy8vD5ZrImLi2PixImYMmWK0DM1RqdOnXDhwgVs3Lix\n0W2Vl5cjKCgIJ06cEEKyhuOsYP3SZGhCWjILCwtcvHgR169fx9mzZ7FhwwauIzWZly9fokuXLlzH\naLUMDQ1b1CLljDGsW7cOJiYmiI+Px6VLl/Dw4cOPjk4B/85JfPbsGaZOnfrZdq9cuYKffvrpk8sT\nkdYpKioKKSkpePz4MebOndukOwtKSEhg+/btyM7OhpSUFAwNDeHi4oKHDx9i8uTJn7xu9+7dCAkJ\nwdu3b5ssW0O4urqiX79+jWrj0KFDkJeXx7Rp0zhfpYSqRkKaULdu3fDnn39i9OjRePHiBY4ePSq0\n/aRFQU5ODs6dO4erV69yHaXVateuHQoKCpCYmAhjY2Ou43yWQCDAokWLEBwcjOjoaGhpaX3xmjZt\n2gDAF1dASEhIwPjx44WSk7QcERER0NPTg5GREfbv398sfUpJScHX1xdRUVEIDQ3FhQsXPvuUQ1pa\nGqqqqsjPz4empmazZGwujx49Qp8+fVBRUcH5hh20SgAhTaxLly54+fIleDweLl68yHUcofr7779h\nYWGB7t27cx2l1ZKRkcGOHTswfPhwkX47Pi0tDQMHDsTz588RFBRUp2IV+LcgB4A3b9588pyMjAy8\nfPmy0aNFpGW5f/8+li9fztkTKnNzc3h4eNRpSo6xsTH8/f1b3S6eu3btwpQpU7Bs2TKhTC9oDCpY\nCWkG0tLSsLOzQ3x8PNdRhCopKQkdOnTgOkarN2/ePHTv3h2//vor11E+6s2bN3B0dISDgwOCgoKg\nqKhY52vFxMQwbtw4HDp06JPnXLhwAcOHD2/R69CS+klJSYGnpyd2794NOzs7ruN80datW7F//37Y\n2dnhwYMHXMcRGlVVVUybNg3jx4+vfRrCFSpYCWkmxsbGSExM5DqG0JSWluLQoUPo27cv11G+CvPn\nz8fx48dRU1PDdZT35Obm4ptvvsGCBQvw888/Q0pKqt5tWFpafvYta29v78/OISQtH5/Px7Zt2+Du\n7o7evXvD3Nwc7u7umDRpEtfR6qRbt2549eoVPD09MXHiRCxcuLB2o4GmFhISAl9f39qXxForKlgJ\naSZKSkqtZgtXxhhmz56Nb775Bm5ublzH+SrY2NhAXV0dfn5+XEd5T2JiInR0dLBw4cIGXV9dXY19\n+/bBxcXlo8fj4uKQmppKqwO0Mnw+H1evXsWePXswdepUKCgowM/PD46Ojli/fj1ycnKwcuXKFvWS\nnbi4OCZNmoTIyEhER0fD1ta2WfrdtWsXXF1d0b1791a9ExwVrIQ0E3l5+Wb7jbupBQUFISwsDAcO\nHGhRP1BaMh6Ph2XLlmHDhg0iNcratWtXlJSUwN/fv0HXb9q0CTIyMnB1df3ocR8fH4wdO5ZWlmlF\ngoOD0alTJ2zbtg0REREwMjJCamoqHj58iKlTp8LZ2ZnzF3waQ0VFBcePH0dGRgZu3LjR5P399ttv\n6NChAyorK+Hg4CByv9QKCxWshDQjgUDAdQShCAsLg4uLC+Tk5LiO8lUZN24c2rRpgz179nAdpZac\nnBw2bdqEnTt31vvaly9f4sCBA9ixY8cnf/Hx8fH5ZDFLWo60tDScPHkSU6ZMwYgRI7B79248fvwY\nR48exc8//wxVVVWuIwpV+/btcfbsWSxfvrzJX8TS0NBAYGAgJCUloaKigkmTJmH16tWt5ufNf6hg\nJaSZvH79GoaGhlzHaBTGGDZs2IBt27Zh0KBBXMf56vB4PBw7dgy7d+9GaGgo13FqDR48GHfu3EFM\nTEy9rgsLC4OTkxMcHBw+ejwpKQlpaWk0T7oV0NfXh7u7O9q1a4eYmBiMGDGC60hNbsCAASgvL8c3\n33yDPXv2gM/nN1lfenp6CA4OhqOjI8rLy7F582bs2rWryfrjAj1jIaSZJCYmtsg36svKyvD06VO8\nevUK+fn5OHToEP755x/o6+tzHe2rZGhoCFdXV9y/f19klhNTUVGBsbFxvedot2vXDpmZmZ887ufn\nh2HDhtFWrC0cn8+HQCBA3759sXXrVoiJfR1jZTweD0eOHEFISAiWLl0KGxubJl3xoG3btjhw4ABm\nz56NS5cu4bvvvmuyvrhABSshzSQhIQH29vZcx6izJ0+e4IcffkBoaCgsLCzQtWtXiIuLw8vLi4pV\njikqKqK0tJTrGLVqamqQlpZWu6ZqXdnb2+O777775KYI/v7+cHd3F1ZMwhGBQAB5eXlcvnz5qylW\n/+Po6AhHR0fExMTg2bNnzbJEl5WVFaysrJq8n+ZGBSshzSQ+Ph7Tpk3jOsZnvX37Fr/++isCAwOR\nnp6O7du34/r165yvv0feJyMjg5KSEq5j1OLxeKisrERubi4MDAzqfF2bNm0wdOhQ3Lt374OCtaSk\nBA8fPuR8/3LSeF5eXrC3t4eKigrXUTgzd+5cDBs2DAoKCpgxYwbXcVqkr+tXHUI4FBMTA1NTU65j\nfNK5c+dgaWmJ8vJy7Nq1C2/evMG3335LxaoIMjc3R3h4ONcxIBAIcOPGDQwcOBB9+vSBpaVlvdvQ\n0tL66C5X58+fh6Oj41dd5LQWMTExDfraaE1sbGzw8OFD7NixA6tXr251O2I1Bx7j4LPG4/HoL4t8\nVXJzc2FiYoL8/HyRXAZq165d8PLygre3N2xsbLiOQ74gMzMT5ubmyM3N5eQRK5/Px+HDh7Fnzx60\nadMGnp6emDx5coN2onr58iUcHBywePFizJ07F3Fxcdi/fz+uX7/+y1kGAAAgAElEQVSOc+fOYeDA\ngU1wB6Q5vXr1Cg4ODkhISICCggLXcTiVl5eHfv36YcKECVi9ejXXceqNy/qNRlgJaQb/ja6KYrF6\n5coV7NmzB4GBgVSsthBaWlpQVVXFixcvOOl//PjxOHPmDLy8vBAaGopp06Y1eNvULl264MGDB3jz\n5g0MDQ0xefJk9OzZE5GRkVSsthIvXrxATk4OVqxYwXUUzqmqqsLf3x/bt29HYWEh13FaFBphJaQZ\nnD17FpcuXcK5c+e4jvKBPn364Mcff8TQoUO5jkLqYdWqVZCUlMTGjRubvW8jIyNMmzYNP/30k1Db\nLSoqgpycHCQlJYXaLuHW//6iTj/7/zV06FBMmzYNY8eO5TpKvdAIKyGt3Nu3b6Gpqcl1jI9KSUlB\nx44duY5B6snR0RGPHj3ipO9hw4ahvLxc6O0qKipSsdoKLVq0CABEYt61qFBRUWk1Ox82FypYCWkG\nolywmpmZ1XvBd8K93r17IzIyEomJic3ed1BQEIYNG9bs/ZKWadeuXejXrx+Cg4O5jiIySktLm3Qj\ngdaIClZCmhhjDIGBgejatSvXUT4qISGBRlhbIAUFBSxduhTTpk1DWVlZs/VbWlqKxMRE9OjRo9n6\nJC2buLg4vv32W9y/f5/rKCJj3rx52LJlCxWt9UAFKyFN7M6dOygrK8OQIUO4jvKBBw8eoLq6ukXu\nwEWA77//Hrq6uhg3blyzzSvz8/ODhYUF7T5F6sXFxQXBwcG4cOEC11FEgrOzM3JyckRqPWVRRwUr\nIU3s2rVrcHd3F7kdXi5duoSxY8fiwIEDVHy0UOLi4jh27BjS0tJw6dKlJu2LMYYdO3ZgwYIF2LNn\nT5P2RVofbW1tbNy4Eb/++isEAgHXcThXVFSEysrKr36Zr/oQrZ+ghLRCampqIrWNJvDvnNrZs2fj\n5s2bNBexhZOQkMCOHTuwatUqVFVVNVk/W7Zsgbe3N0JCQtC7d+8m64e0XpMnT0ZaWhqioqK4jsK5\ny5cvY/DgwTRYUA9UsBLSxMTFxUWuYF28eDGmT5+O7t27cx2FCMHAgQNhYmKCvXv3Nkn7sbGx2L17\nN65fvw59ff0m6YO0fmJiYlBSUmrWOdeiys9nCkb1uwJk8t7/j3wSFayENLFLly6J1Bqnly9fRkhI\nCH755ReuoxAh+u2337B161a8fv1aqO2WlpbCzc0N69atQ/v27YXaNvn6uLi4wNvbm+sYnMvKAQx0\nuU7RskhwHYCQ1q6iogKqqqpcxwAAnDt3DgsWLMCVK1cgKyvLdRwiRCYmJli1ahXc3d1x7949oTxq\nLCwsxIgRI9C9e3fMmzdPCCnJ127JkiXo3Lkzli1b1iSj9VlZWQgICMDr169RUlKC4uJiWFlZYcqU\nKQ3ejU3YBAIBXsUDRvSwol5ohJWQJmZqaorY2FhOM1RVVWHRokVYsWIFAgIC0KtXL07zkKaxdOlS\nvH79WiijrBkZGbC3t4eVlRX++usvkdxWmLQ8GhoaGD16NK5du9bott68eQNfX1/cuHEDq1evRrdu\n3dCpUydcvHgRZWVlUFZWhrm5OS5fvgxjY2MMGDAAW7Zs4XxL1KtXr8JAF9CjBxb1QiOshDQxDQ0N\npKenc5rhzz//xPPnzxEeHg5lZWVOs5CmExkZifT09EaPJD179gyurq6YM2cOfvjhBypWiVApKSk1\neJenoqIiBAUF4ejRo3j27Bm6dOkCHo+H7t274/fff4etrS0kJN4vbebPn4+oqCikpKTgzJkzUFZW\nxunTpzFx4sT3zhMIBCgrK4O8vHyD7+1LYmNj8f3332PLlouA1ugm66c1ooKVkCb26NEjeHh4cNb/\n27dvsW3bNvj6+lKx2sopKChASkoKioqK9b6WMYbQ0FB4eXnBz88Pf/zxR4vb55yIPsYYYmJiPigW\n64LP58PS0hIdO3bE6NGjce7cuTpPbTI3N4e5uTlcXFwgJyeHqVOnws7ODrq6/04kFQgEmDFjBlJS\nUhAQEFDvbHXx6tUr9O3bF2vXrsWoUaOapI/WjApWQppQamoqMjMzYWNj0yTtBwUFwc7ODlJSUh89\nnp6eDjs7O8ycOZN2JvoKGBkZQU1NDS9evIC9vf17x2pqahAYGIgXL14gNzcXYmJi6Ny5M6KiohAb\nG4vg4GC0bdsW7u7uiI6OhoqKCkd3QVozd3d3vHz58pPTkgoLC3Hnzh1kZWUhJSUFJiYmUFVVRUZG\nBp4/fw5DQ0P4+/s3KsOBAwdgaGiIoUOHIjg4GKWlpVi2bBl8fHyadBOVzZs3Y8mSJZg/f36T9dGa\n0RxWQprQiRMnMGrUqCZZa6+mpgZOTk7o2bMn1q5di/v377+321FycjImTZoEDw8PWhHgK6KmpgYH\nBwecPn0aAoEA3t7eGD9+PAwMDPDjjz8iMzMTCgoKkJeXx40bNyAhIQFdXV34+fkhISEBa9asoWKV\nNJn09HR06dIF5eXl7328vLwcf/31F3r27IlDhw4hKioKSkpKePToEY4fP47o6Gioq6vj3Llzjc7A\n4/GwcuVKSElJITAwEA4ODlBRUYGOjg6sra0RERHR6D7e48/Do9083L99CvM7/wT48/79TwQ8efIE\nZ86cabad8hqDxzhIyePxWsQnh5DGqK6uhoGBAW7evAkLCwuht88Yw4wZM3DkyBH07t0bOTk5qKqq\nQmVlJczNzfH8+XPMnj0b69atg6SkpND7J6Lr1KlT2LlzJzIyMqCvr48FCxagR48e6NSpE9fRyFeu\noKAAx44dw+bNmzF27Fjo6upCUlISBw4cgKmpKdzc3DBt2rRmyXLw4EGsXLkSXbt2xaxZszBlyhRI\nSUnB1dVVuEtv+fMwZSfQzRhY9L/TVgdyXweNGTMGly5dwrx587BkyZIvjjBzWb9RwUpIEzl79iy8\nvLxw9+7dJu3n2rVrcHNzQ0pKCt6+fYs2bdrg+fPn6Nevn8gsp0WaX2VlJdLS0mBoaChy2wITkpWV\nhUOHDqG0tBSlpaVwcXHB4MGDmz1HdnY2ampqYGJiAnFxcZibmyMvLw+nTp0S2lSuf/7gYdBqIPIg\noPm/Dy9EoGB9/fo1jIyMYGVlhTdv3uDhw4fo0qXLJ8+ngpWQVsjZ2RkzZ87E+PHjm7Sf6upqyMrK\n4tGjR7C1tW3SvgghdZeAj68M3wGpzZyE1NWVK1dw8eJF3L17FzNmzBDKdKrRfXgYaA3M/b+7YItA\nwQoAt27dgoeHB7p06YLevXtjw4YNnzyXClZCWpnk5GRYW1sjLS0NMjIyTd7fmTNnMG/ePBQUFDR5\nX4SQuqGCteXKyclBp06d8OLFi9qVBBqiuroaykpSSP0bUFb6Pwd1RacOunDhAtzc3Gr/HBgYiP79\n+39wHpf1Gz0nIkTIkpKSMGDAAEyePLlZilXg3xcWjIyMmqUvQghp7dTV1TFo0CDcvn27Ue3ExMRA\nVfkjxaqIGTdu3HvL2K1atar2/xljEAgEXMR6DxWshAhZYWEheDweduzY0Sz9vXz5EqtWrcLJkyeb\npT9CCPkaeHh44KeffkJQUNBHC7bs7GwcPHjwvZ2zGGPw9/eHh4cHRo8ejWHDhsGpT3Ombrjly5fD\nysoKERERyMjIwJAhQ7Bjxw7Y2tpCSUkJgYGBnOajdVgJETI1NTUUFxc36Zv5fD4foaGhOHnyJM6d\nO4fdu3ejc+fOTdYfIYR8bQYNGoS//voL8+bNQ1FREaytrdGjRw8MHz4cRUVFmDBhArKyslBZWQlP\nT0/k5+fD1dUVhYWFmD17NrS1tbFy5UrYtu/N9a3UCWMM4uLi6Nq1KxISEuDr64t//vkHs2bNwrZt\n2zhf7o7msBIiZIwx6Onpwd/fH2ZmZo1u78cff4Sfnx/U1dVRXFyM9PR05ObmokOHDpg0aRImT54M\nQ0NDISQnhBDyMYmJiYiMjMTDhw9x5coV5Ofn4/z580hJScHRo0eRk5ODkpISzJgxA2vWrGlxK3OE\nh4dj/vz56N69O/bs2VP78SdPniAsLAy7d+9GQkICxMXF6aUrQloLPp8PKysrbN++HUOGDGl0e1u2\nbMHTp0/h6ekJRUVFtG/fHhoaGh/sl00IIaTpMcbA5/Nrn6IxxhAUFARZWVn07t0yRlP/r6FDhyIz\nMxMXL16EgYEBACAsLAwuLi6wtbXF/PnzMXDgQE7rN/qJR4iQeXt7Q01NDS4uLkJp782bN+jXrx8G\nDBgglPYIIYQ0HI/He2/KF4/Hg5OTE4eJGm/Tpk3Yt28funXrhps3b8LW1hbPnj2Dubk5rl69ynU8\nAFSwEiJ0f//9N8aMGQMer/5b7/H5fPj4+MDS0hJmZmbg8/l4+/YtNDQ0miApIYSQhsrOzsb69eth\nbGyMMWPGQFNTE0+fPsXr16+hp6cHR0dHriN+VlxcHGbOnAkHBwcsXLgQO3fuhJycHE6ePAktLS2s\nWbMGN2/e5DpmrZY1yYKQFiA8PByWlpYNuvbXX3/F+vXr4eTkBF1dXcjLy6OgoAATJkwQckpCCCGN\nERMTg9OnTyM6OhrW1tZQVVXF8uXLcePGDYwcOZLreF80ceJE9OnTB/v27YOxsTEGDx6M33//HT17\n9oSbmxtWrFgBa2trrmPWohFWQoRIIBAgMjIS2traYIzVe5S1TZs26NWrF3777TfcvXsXAwYMaLa1\nXAkhhNRdr169oK+vj969e+OPP/5AZWUl2rZtC4FAAElJSdTU1EBcXJzrmB+1fPlyVFRUYN26dejX\nrx9kZWXh6OiI6dOnw8PDA5s2bcKyZcu4jvkeeumKECGKiorCsGHDUFBQAA0NDfj7+0NfX7/O18fE\nxMDZ2Rl5eXlgjCExMRHt27dvwsSEEEIa6t69e1iyZAmeP39e+7GzZ89i9erViI2NFcmXYx8+fIjx\n48cjOjoaSkpKqKiogLOzM/r164cNGzaguroaUlJSH72WdroipJX4+++/wePxoKOjg6ysrE/+o/8U\nU1NTLFu2DBUVFejVqxeqqqqaKCkhhJDGsrW1RUZGBoKDgwEAfn5+mDhxIlavXi1yxSpjDN7e3hgz\nZgwOHz4MJSUlvHv3Dm5ubtDV1cX69evB4/Hq/XOruYjWZ5OQFk5BQQHV1dXQ19cHn8+HlpZWva5f\nt24dLl++jLi4OJiYmDRRSkIIIcIgIyMDLy8vjBs3Dj169MDjx4/h7u4OBwcHrqO9JywsDMuXL0dB\nQQFu3boFa2trVFdXY8KECYiPj0dERITIrx1LBSshQsTn8zFmzBhUVFRAV1e33tffvn0b+/bto2KV\nEEJaiFGjRkFZWRmJiYnYuHEjzM3NuY5U69q1a/jll1+QlZWFn3/+GTNnzoSEhAS8vLywYcMGZGRk\noGPHjpCWluY66hfRHFZChMjZ2RmamppISEjAunXrMGjQoDpdV1hYiIKCAnTr1g0JCQlQU1Nr4qSE\nEEJaq5qaGqxduxbHjh3DX3/9BWdn59opCqmpqbCyskJAQAC6detWr3Zp4wBCWgGBQIDQ0FDExcWh\nXbt2dbqmqqoK69evx969eyEvL4/JkydTsUoIIaRR5s2bh9jYWISGhn6wjvfmzZsxderUeherXKOC\nlRAh8fPzg7Kycp2LVQDYv38/7t+/j7i4OGhqajZhOkIIIV+DwMBA+Pv7IyIiAm3btn3v2I0bN+Dr\n64vY2FiO0jUcFayECIGvry+mTZuGa9eu1eu6a9euwdbWlnayIoQQIhTt2rVDUVERVFVVYWpqinnz\n5kFFRQUnT55EbGwszpw5A2VlZa5j1ptovxJGSAvx6tUrqKurY8iQIVi+fHmdr/Py8kJQUBBcXFyw\nf/9+ZGZmNmFKQgghrV3Xrl2Rn5+PkpIS7Nq1C0+ePMH58+fRv39/REVFwcnJieuIDUIvXREiBIwx\neHp64vbt2zh79ixsbGzqfG1lZSVOnjyJzZs3Y8GCBViyZEkTJiWEEEIahl66IqSFCwgIwL179/D8\n+XMoKCjU61oJCQkoKiqioKAAbm5uTZSQEEIIabmoYCVECLy8vLB48eJ6F6tZWVlwdHSEgoICzp49\nS9uwEiLipkyZAmVlZYwZMwZmZma0qgchzYTmsBIiBJWVlQ16y3/t2rVwdnbGkydP6rxmKyGEG0lJ\nSTh+/Dhqamowa9Ys6OnpcR2JkK8GjbASIgSqqqrIy8ur8/nFxcU4evQofHx8EBsbCx6P14TpCCHC\nIBAIIC4ujlOnTuHdu3eYPXs215EI+WpQwUqIEGhoaCAqKqpO56ampsLOzg49e/bErVu3oKKi0sTp\nCCHC0KFDB6Snp0NCQgJt2rRpEdtZEtJa0CoBhAhBamoqrK2tce/ePXTu3Pmz5/722294/Pgxzp8/\n36C+IiIiICMjAwC4fPkypk+fDlVV1Qa1RQghhNQVl/UbzWElRAiUlJRQWVkJJSWlL57r5OSEp0+f\norq6ul59xMXFYcCAAejZsyc6deoEa2trBAUFoUOHDkhNTW1odEIIIUTkUcFKiBDIyMhAXV0dFy9e\nRE1NzSfPe/36Nc6ePYuUlBRISUmBz+fXuY8rV64gICAAEhISUFFRgYaGBh49eoRhw4ZBS0tLGLdB\nCCGEiCSaEkCIkMTHx8PZ2RmVlZXYsGEDpk+fDjGxf38njI2NhaenJx48eADg340GOnfujODgYLRp\n06ZO7QsEApw6dQp2dnaorq5GcXExLC0tISUl1WT3RAghhPyHy/qNClZChEhaWhpVVVUAAAMDA+zb\ntw8DBw7E0qVLcfnyZbi6umL9+vXIzs5G79694efnV69dsQghhBCuUMFKSCshEAiwceNGrFmzpvaR\nv5iYGHg8HrS1tREXF1c7Inr69GmsXr0aW7ZsgaurKyQkaNEOQgghoosKVkJamdu3b2P58uXg8/mw\ntLSEkpISNmzYAHV19ffO8/Pzw7Zt25CWlgZHR0d07doVrq6u0NXV5Sg5IYQQ8nFUsBLSCtXU1MDX\n1xc7d+5Efn4+7O3t0alTJ4wfPx76+voAgLdv32L48OEQFxdHWVkZIiMjceLECXz33XccpyeEEELe\nRwUrIa0YYwyhoaF48eIFXrx4gbNnz0JFRQVKSkqIj4+v3S3n3bt3sLS0xKxZszhOTAghhHyIClZC\nviKVlZVISkpCfn4+rKys6rxKACGEEMIlKlgJIYQQQohIo52uCCGEEEII+QQqWAkhhBBCiEijgpUQ\nQgghhIg0KlgJIYQQQohIo4KVEEIIIYSINCpYCSGEEEKISKOClRBCCCGEiDQqWAkhhBBCiEijgpUQ\nQgghhIg0KlgJIYQQQohIk+A6ACGkYSorK3H37l28efMGAwcOhJGREdeRCCGEkCbBYxxsCsvlXrSE\ntGSlpaWIi4tDYGAgDhw4gLZt2yIxMRHv3r1DcXEx2rZty3VEQgghrRSX9RuNsBIiwkpLS5GTk4O/\n//4bEydOfO+YhoYGsrOz0aNHDyxZsoSKVUIIIa0WFayEiKCXL19i8+bN8PPzQ0lJSe3Hf/rpJ9jb\n24PH4+HNmzcYPXo0VFVVOUxKCCGEND2aEkCICOHz+fD19cX8+fOxYsUKODo6YujQobh27Rp69uzJ\ndTxCRFpNTQ3ExMTA4/G4jtLilZSU4MyZMzA0NIS9vT2kpaW5jiRU1dXV0NPTw5w5c/Dzzz9DTIze\nQa8LmhJAyFeupKQEv//+O/bt2wc9PT1cvHgRffv2RXx8PHJyciAQCLiOSIjIYIwhOjoaAQEB+Pvv\nvxETE4Pk5GQUFRVBQkICqqqqaN++Pfr27YvZs2fD1NSU68hC9eLFCxgZGUFBQUHobQsEApw5cwYr\nV66EjY0NcnJyEBsbi+nTp2PlypV1eqLD5/Ph5+eH58+fo7y8HKGhoUhNTYW9vT1++OEHdOjQQei5\nP4UxhsrKSggEAjDGkJubi4yMDDx79gz5+fk4d+4cLC0tMWrUqGbLRBqGClZCOFRTU4OjR49izZo1\ncHR0xO3bt2FmZoY///wTs2bNQmZmJgAgLS2N46SEiIbMzEwMHDgQxcXFcHFxweDBg7Fs2TIYGBhA\nRUUFVVVVyMnJQWpqKi5cuIARI0bg6tWrraZoDQsLg6OjI7S0tBAYGAgdHR2htMsYw4kTJ7BmzRpo\namrCx8cHvXr1AgAkJydj8+bNMDExwfTp0zF37lwoKSlBXFwcAoEAWVlZSEtLw5s3bxAeHo6LFy/C\nwMAAzs7OUFVVxeLFi2FoaIgrV66gR48eGDJkCHR1dVFTU4NBgwbB3t4ekpKSH2SqqKjAqVOncPr0\nabx9+xby8vIwMTFB+/btISUlhZiYGBQWFkJcXBxaWlrQ0dGBmpoaEhISEBkZibi4OGRlZUFSUhLi\n4uIAAFVVVWhra0NDQwN37txBTU0NJkyYgLKyMtjb20NMTAwlJSWorq5GTU0NSkpKUFxcDMYYGGPg\n8/kQCASQlpaGjIwMlJSU0K5dO+jq6grl74F8Gk0JIIQjjDGMHz8eqamp+O2339CjRw88fPgQs2bN\ngpaWFjZv3oyOHTtCWVmZHnGSrxJjDCEhIdi/fz/Cw8MBAAkJCVi6dCnWrl37xX8XfD4fP//8M44e\nPQozMzOMGjUK2trasLa2/ugycEVFRUhISEDXrl0hJSXVJPf0KdnZ2bh9+zZCQ0MhLS0NZWVlZGVl\nISUlBWlpaZCQkACPx8OrV6/w559/IiUlBb/88gssLCxgYGAAWVlZpKamorq6GqqqqpCTkwOfz0da\nWhratGkDc3NzdO3aFa9fv0ZpaSlqamqQmpqKnJwc6OvrIyoqCgCwb98+2NrafjRjSkoKdu7cifPn\nz6OyshI1NTXg8Xi1BZuBgQFMTU0xcuRIdOrU6aNt5Ofnw9vbGyUlJeDz+bh69SpiYmJgaWkJY2Nj\nqKmpgcfjITs7G7du3UL37t0xb948GBoaoqSkBPHx8cjMzERFRQVMTU2hpqaG6upqZGZmIi0tDdnZ\n2ejQoQMsLCxgYmICLS0tSEh8fmzu7t272Lp1K168eAExMTHIy8tDSkoK4uLiaNu2LRQUFGqnDEhI\nSEBMTAxVVVUoLy9HYWEhUlNToaenh6tXr6J9+/aN+CoQfVzWb1SwEsKRqqoqKCkpIT4+Hu3bt8ej\nR48wduxYeHl5YeTIkVSkkq9KYWEh7ty5g1u3buHu3buQlpauHdGaNWsWnJycIBAIYG5uDjk5uXq1\nXV1djQsXLiA4OBgZGRl4+PAhrKys4OrqCiMjI1RVVSEoKAjHjx9Hu3btkJKSAkNDQxgYGMDAwAA6\nOjpQV1eHjIwM5OTkMHz48NoRu8aoqKhASEgILl68iGPHjqF///6wtbVFTU0NCgoKoKGhAV1dXejq\n6kIgEKCiogK2tra19//u3TuEhoYiPT0dZWVl0NHRgZSUFPLy8lBeXg4xMTHo6OigrKwMYWFhePXq\nFYyNjaGoqAgej1d7XykpKVBVVcWwYcM4mctZXFyMsLAwpKSkICcnB4wxqKqqwsHBAcbGxs2ep76q\nq6shJSWFefPmYeLEibC2toasrCzXsZoEFayEfKWGDBmCCRMmwN3dHX5+fli4cCEuX74MCwsLrqMR\n8kWPHz9GWVkZOnfuXK+RpbS0NLx48QLx8fF49eoVwsPD8erVK9jZ2WHw4MFwdnYGYwxlZWWwtrYW\nehFVXl4OPz8/XL16FZmZmeDxeOjduzc8PDxgZGSEkpISJCYmIjk5GcnJyUhLS0NOTg4qKysRERGB\n0aNHY8OGDQ3qmzGGmzdvYs+ePfj777/RqVMnODo6YtmyZdDQ0BDqfZLmExcXh4MHD+LRo0eIioqC\nubk5bG1tYWJiAiMjI+jr60NPT++DecdJSUmQk5ODpqYmR8nrhwpWQlqxiooKHDx4EP/88w86deqE\nFStWQFJSEqdPn8ayZcvwzz//QEtLCwBw5MgRfP/992jfvj3c3Nwwffp0tGvXjuM7IOR91dXVmD59\nOh49egR9fX28ePECs2fPxuTJk9GmTRuUlZXB3Ny89ikBYwzXr19HcHAw7ty5g9evX8PGxgYdOnRA\n586dYW5ujh49ekBGRobjO/uyrKwsWFtb4/Dhwxg0aFC9rmWM4eeff8b58+exfv16uLi4QFFRsYmS\nEq6Ul5cjJCQEISEhSExMRFJSElJSUpCcnAxZWVmYmppCUVERVVVVCA8PB5/PR5cuXTBz5kyMGDEC\nSkpKXN/CJ1HBSkgrJRAI0LdvX6iqqmLEiBE4c+YMOnTogNGjR2PWrFnw8/P7YDSVz+fjyZMnOH78\nOPz8/HDz5k1069aNozsg5N+XA5OSkpCcnIyEhAQcP34cKioq8PX1hYyMDNLS0rBx40YEBgaiuroa\nAKCiooJ58+ZBX18fV65cwe3bt/Htt9+ib9++sLe3/+K8QlFVWloKd3d3ODg4YOHChZ88JyUlBWVl\nZSgqKkJqaioSExPx7NkzpKWlITAwkEZTv0KMMWRlZSEuLg4lJSUQCATo168fpKSkcPPmTRw6dAj3\n79+HgYEBNDQ0oK6uDm1tbVhaWuK7774TiWliVLAS0kpVVFRAQUEBlZWV4PF4KCwshK2tLUxNTdG+\nfXvs37//s9efOHEC27dvR2hoaIsYfSL1wxjD69ev4e3tjSNHjsDe3h7i4uLIy8vDqlWr0Lt379rz\nLl26hGvXriEyMhKysrJQVVWFh4cHfHx8EB0djczMTJSWlqKiogJiYmKQlZWFpqYmNDU1oaSkhC5d\numD16tVfnFtXWlqK169fIyUlBYmJibh58yYePnwINTU1GBoaQkdHByNHjsSIESM++mb3f3n9/Pxw\n4sQJZGVloVevXliwYIHQ3mjnQlRUFNatW4cbN27A2toap06d+uDNcMYYNm7ciC1btkBXVxdt2rRB\n27ZtoaenB0NDQ3Tq1AkjRoyAvLw8R3dBRF1FRQViYmKQk5OD7OxsZGRk4Pz585CUlES/fv3g7OwM\nR0dHzopXKlgJaaXy8vKgr6+Pd+/e1X4sNTUV33zzDfh8PvJlsDgAACAASURBVAICAmBtbf3J6xlj\n+Pbbb/HmzRtcvnwZ6urqzRGb/I/q6mpkZWVBXV29wYunM8aQk5OD5ORkREVFITg4GJGRkYiPj4eU\nlBQGDRqEOXPmICwsrLYIXLNmDUaOHAljY2Pcu3cPKSkpmD17Nrp3746qqio8efIEZ86cwdSpU9Gn\nTx9oa2tDXl4e0tLSEAgEKCsrQ2ZmJt6+fYvi4mJcuHABDx48gKOjIxISEpCRkYHy8nLw+XwA/46i\nVldXgzEGAwOD2jl3zs7OGDBgAJSVlYX2OW2JRo4ciTt37mDIkCHYtGnTB8tkVVRUYNu2bbhw4QKN\noBKhqqysxJ07dxAaGooLFy6Ax+Ohe/fu0NPTA/Dvkzx5eXnIycmhpqYGAKCoqAgVFRXo6enByMhI\naNMMqGAlpJWIjIzEpUuX8PTpU5SWlqKoqAgWFhY4fvz4e+cdOXIEf/zxB7KyshAcHFz7jedjBAIB\nli5ditzcXJw6daqpb+GrVFpaCl9fX4SHh0NeXh6WlpYYPXo0cnJyMHjwYKSlpaGwsBAzZ87E8OHD\nISUlheTkZCQlJSEjIwNv374FYwwKCgqYM2cOpKWlcfv2bURERODVq1dISUmBjIwMDAwM0LFjR/Tu\n3RtWVlbo0KEDNDU1PzpaUlBQgM2bN6O6uhpdunTBhAkTGj0y9+LFC4SGhsLU1BS6urqQlZWtXS5J\nTEwMkpKSkJWVpV1/PiI2NhYbN27EqVOn8PTpU/To0aP22NWrV+Hp6QlLS0vs3bv3o0tmESIMAoEA\nT548QUxMDFJTU2t3dnv37h3KysogISEBxhiKioqQl5dX+6RER0cHPXv2RLdu3TBgwAB07ty5Qf1T\nwUrI/7d06VKEhoZCXV0d7du3h7KyMhQUFKCqqgoDAwNIS0sjPT0d3bt3h4GBAddx3+Pp6YkrV65g\n3LhxcHR0RNu2bZGVlYX+/ft/8OLU6dOn8eeff8LFxQV3797FrVu3Ptt2dHQ0Ro8ejdjY2Ka8ha/S\nmzdvMHjwYBgZGcHR0RFlZWXw9fVFfn4+0tPTMWLECFy+fBm5ublYsmRJ7dvi/41c6OjoQENDA2Ji\nYkhLS8OKFSugr6+PoUOHwtraGmZmZtDX10fbtm25vlXSCPfu3cOQIUPg5eWF7777DsC/L9ccP34c\na9euhY+PD+zs7DhOSciHampq8PLlSzx79gzPnz/HpUuX4OHhga1bt9a7LSpYyVettLQUGRkZMDY2\nhrq6OvLz82FgYICamhqkpaV98mslJiamzrvX5OXlIS8vDx07dhRm9Pfs2LEDW7duhY2NDTIyMpCU\nlITS0lIAwK5du5CQkICoqCgkJycjLy8Pt2/fRkZGBlavXv3FQrSqqgomJib4888/6/1mMvlQVVUV\nHjx4gPPnz8PHxwdr1qzBokWLao//t+B6eXk5tLS06vU4raamRihrdBLREhkZiVWrViEwMBDKysqQ\nkpJCdnY2bGxs8Ndff31yoXxCRElYWBjc3NywdOlSeHp61vt6KlhJkxIIBCgoKKjTHtBcOH78OKZM\nmQJJSUm0adMGYmJiKC8vR1VVFWRkZKCiogItLS0YGBjAxMQEUlJSePToEbZt21b79jxjDHl5ecjM\nzERhYSG6dOlSO+cuNjYWLi4uKCwsxJAhQ2Bvb1+7a0pRUREGDRoEbW3tz2Z8+vRp7a5Tn5OSkoLw\n8HDo6upCQ0MDS5cuhaysLCoqKqCsrAxFRUUIBALk5ubi2bNnKCkpweHDh9G/f/8vfp4OHjyIixcv\n4vbt2yLxtmhLkpubCz8/P7x48QIhISEIDw+Hubk5Ro4ciRkzZtDSYaTOqqurkZubi4qKCujq6rbY\n1Q7I1yckJARDh/6/9u48LMpyfwP4zSL7sMwAI8IgNKhsMxMKKridUvPkUmqUop5MLT1pZeupU9cv\nS9vTckntmGWlXrabC6ZmGmomaiqbAoHsKDPAADMgODPM7w/zPXBEk0TnFe7Pdc3F8K7fecXh5pnn\nfZ4xWLlyJe6///6/dAwGVrphLly40OpGEZVKhcTERDzwwAMIDw/H2bNn0aNHDygUCkRHRyMkJATB\nwcEoLi7G6tWrcezYsSveFFRYWAi9Xn/dQy4ZjUbcdttt2LFjB5RKJerr6+Hq6go3Nzfo9Xrk5eXh\nwIED2LRpE3r06AG1Wi2sKysrQ2FhIfLy8uDo6IgePXpAIpEgKysLEokE9fX1aGpqwqpVqzBx4kR8\n8sknyMjIQFFRETw9PeHs7Izdu3cjNjYWbm5ucHd3h1wuR2BgIKKiotC9e3esXLkS27ZtQ2NjI+Lj\n49HY2IiamhrcdtttSEhIwIwZMyCVSmEymZCbm4ujR4/iwIEDOHjwIIqKiuDu7o7u3bvDz88Pcrkc\nQUFBCAkJEfoTXelO67auk1QqhVarFfU4fWKzZs0avPDCCxgxYgTi4uLQr18/xMbGXjaANxFRZ1NW\nVobs7GycPn0aCxcuxMcff4xx48b95eMxsNJVmUwm5OTkoKysDFKpFD169MCFCxdQUVEBg8EgBCBn\nZ2dYLBbU1dWhrKwMer0eEokEW7ZsgV6vx3fffYeysjLhuNXV1fD29oZSqURBQcFl53Vzc8PUqVOF\nfaKiouDn54f8/HwcO3ZMGARZo9Fg06ZN8PDwQGNjI7p169buj0Q3btyIJ598EuPGjYNUKoVGo8G0\nadNatSQ2NjZi48aNqK6uRn19Pby9vREUFISePXsiLCysVeun2WxGeXm5MA/01erR6XQ4duwYmpqa\nYDQace7cOZSUlCArKwtnzpzBAw88gOeeew5WqxUHDx6Ep6cnPD098fvvv2PXrl3YvHkzgItdG0JC\nQhAXF4fBgwdjyJAhCAsLa/c0kv/LarXi5MmTePTRRxEaGopNmzZd1/G6EovFIrSAffjhh5gzZ46N\nKyIiunGsVivS09ORnp6O/fv3Y/PmzVCpVAgODsacOXOEofL+KgZWAgA0NDSgqqoKZ86cQUZGBtLS\n0oT5n4ODg6FQKFBdXY3y8nI4OztDLpdDIpGgpKQEhYWFMJlMcHBwgIeHBwIDAyGVSoU7BfV6Pcxm\nM9zd3dHU1CSMC+ro6AgHBwch0JnNZpjNZkgkEowaNQoDBgxAaGgompubkZWVBZ1Oh7CwMKhUKgwa\nNAhWqxXz58/Hli1bEB4ejsOHD8PDwwP33HMPYmNj4e/vj549e0KtVrf66OxS/9T6+no4OzvDxcUF\nVVVVSE1NRU1NDb788kt4enpi6dKlop+mVK/XA7g4jEhH3V1tsViwePFirFmzBuXl5ZBKpXj55Zcx\ne/Zsdge4Rs3NzUhJScH69euxbt06rFy5EnPnzrV1WUREN0RtbS2effZZ7N69GwkJCVCr1Zg9e3aH\ndgdkYO2ktFotPv74Y+Tl5cFgMGDWrFkwmUwoKChAZmYmDh8+jIKCAmH8Q3t7e+FueJVKBZVKhX79\n+gkfgV/Npet5tTDT1NQkBERnZ2dYrVZYLBaYzWbhq6OjI1xcXODk5NSuYHT8+HGUlpbijjvugFar\nxdatW3Hq1ClotVrk5eWhpKQEffr0EVoxy8vLIZPJhEH1z58/j6qqKvj4+EChUCAoKAh79uzByy+/\njH/961/XXMetzGw24/jx40hJScGGDRvg7e2NDz74AEql8rpbabsCk8mEjIwMHD58GIcOHcL+/fvh\n7e2NadOmYcKECejVq5etSyQi6nB6vR7vv/8+/vOf/2Ds2LFYunTpDRuVhIG1k5owYQJqa2sxefJk\nmEwmrFu3TpgtJjw8HPHx8ejduzccHR2FoNhZVVdXIzc3F66urnB3d0ePHj0uC2EWiwVarRbFxcUo\nLCyEm5vbdfW1ESudTofTp0/j9OnTyMzMFMbpLC0thVKpxLBhwzB69GjcfffdbE1tw/nz51FRUYGC\nggLhGh4/fhyZmZlC3+CEhAQkJCQgPDyc15CIOi2DwYAhQ4ZArVbjxRdfvOGjVTCwdiL19fVISUlB\ncXExPv30U+h0OuTn59u6LLoBLvUjBi62cDc2NsJgMKC+vh4ODg5wcnJCt27dYDab0dTUhLy8PLz1\n1ls4e/YsIiIiEBERAZVKhYiICISEhAhTOXZVWq0WaWlpKCgoQFFREcrLy1FeXo7S0lIYDAahJb6p\nqQndu3dHcHAwIiIiEBUVhb59+0Kj0fBGKiLqUv7v//4PRUVF+Oyzz27KH+e2zG8cj+MPKSkpmDNn\nDqZPn47x48fDzs4OtbW10Ov10Gq10Gq1UCgUSExMFPp7Njc3C1MtHjlyBIcOHUJqaipiY2PRq1cv\nDBs27LrvoCfxWrZsWavuCv7+/ggMDIS7uzuam5tx4cIFXLhwAd26dYOTkxNkMhmWL1+OESNGsNXv\nD1qtFgcOHMDmzZuRnJyMvn37CtOCDh48GAEBAVAoFPDy8oKTkxNcXFzg5eXF60dEBCAtLQ2JiYld\n4j2RgfUPTk5OyMnJQXFxMUaPHg0nJyd4eXnB29sb/v7+8Pf3x3fffYcXX3wREokENTU10Ol0kEql\niIyMRFxcHObPn4+hQ4fCy8vL1i+HboJ58+ZBp9MhNzcXxcXFyM3NBQCMGzcOTzzxxC3XxcNqtaK8\nvBwA0KNHjw59A7RYLMIUgXl5eThx4gT279+Ps2fPYvDgwRg5ciRWrFjR5eerJyJqj5kzZ+KRRx5B\nYGAghg8fbutybqguEVitViu++eYblJWVwc/PD/7+/oiIiEBQUJCwTVxcHHx9fREQEID9+/dDoVC0\neZwTJ07A3t4e3t7e8PX1ve65venWVFxcjH379sFiscDNzQ39+/fHsGHDsHPnTjz//PMAcEvcLLZ/\n/368/fbbKCoqQkFBgdBR32AwoFevXhgzZgymTJmCqKiodh+7oKAAb775Jo4ePYqcnBzIZDL06tUL\nSqUS0dHRePTRR6FSqTgrFBHRXzR+/Hi4uroiKSkJDz30EPz8/CCTyRAbGyv6EXbaq0v0Yf3oo48w\ne/ZsTJw4Ec7OztBqtTh58iRCQkIwa9YsPPTQQ3B1dUV2djYWL16MLVu2wMPDA+PGjcNLL70EuVx+\n02qlW0O/fv2gUCgwcOBABAUFwWAwoK6uDq6urhg5cqTob/aprq7G/Pnz8fPPP2PBggWIi4tDaGio\n0Ae0rq4Op06dwrfffov169dj0aJFeOSRR9p1jsTERMjlcsyYMQPh4eH8446I6AYwm83YvHkzVq1a\nhYyMDFRVVQG4OJ1wdHR0h56LN13dYDqdDu+88w7WrVuHhx56CHfddRfUajW2bduGl156CS4uLigq\nKhIChtVqRXZ2NtasWYMPP/wQFosFEokEEokEUqkUvr6+kMlk8PDwEIaI8vT0hFqtFoKMmMMKXT+5\nXI6xY8figQcewF133XXT/r2vZfiyP5OamoqpU6di9OjReOONN/40SGZnZyM+Ph4lJSXtCp3x8fEY\nMWIEIiIiAFzsFqDX61FVVYXKykoYDAaYTCY0NjbC19cXH3300V9+TUREXUFhYSEee+wxZGRkwGQy\n4cKFCzAajZDL5VAqlQgMDIRcLoebmxvmzZvX4Q1uDKw3SW5uLjZt2oQff/wRWVlZ8Pb2hre3N6xW\nK1JTU1tNYdpSU1MTDAYDDAYDqqurUVlZicrKSmHaz6amJuj1epw8eRLHjh0DAPTv3x9JSUmYMmXK\nzXyJdJOcPn0a27Ztw8qVKzF58mRMmzYN0dHRNyy4Wq1WbN68GY899hi6d++OuXPnYtasWX96vrq6\nOhQUFAiP5ORkZGZm4oMPPkBiYuI1nbe0tBTDhw/HZ599hvj4+GuuOTU1FcuXLxf+r9vb28PHxwcy\nmUwYg/fSeL9JSUn46KOPMGPGDHYRICJqQ1FREWJjY/Hss89i0qRJwg29Eonkpt0zwcDaiVitVmRl\nZeH111/HF198gbS0tE7Xj4T+Kzs7GytXrkRycjKqqqqE1ndfX18EBwejT58+6NOnDzQaDQIDA694\nnPT0dGg0GgCAp6cn3Nzc4OrqCg8PD7i7u6O+vh5WqxUrVqzAsWPH8O6776KsrKzV7GH/a+/evRg+\nfDgiIyOhVCoREhKCQYMGYfz48Zf9cdbQ0ICioiIUFhYiPz8f2dnZyM7ORkZGBuzs7BAbG4vVq1e3\n2be7PaxWK4xGo9DKWllZCa1Wi++++w5btmyBRCJBXV3ddZ2DiKiz0ev1mDRpEjQaDd59912b1cHA\nKjJr167F0qVLMWjQIERGRsLFxQWNjY3Q6/XCw8/PT5j6rLi4GFlZWcLsVWfOnMGAAQMwa9YsTJo0\niS1GXYDVahU+7r70KCgoQE5ODnJycnDixAm4uLhgwIABUKvViI6ORmRkJORyOTw9PWGxWPD6669j\n/fr1KCgogEwmQ21tLcxms3AOlUoFX19f5OfnY+HChZg+ffpVa9qxYweefvpp6HQ6xMTEwNXVFU5O\nTkJY1Wq1OHv2LM6dO4f6+nr07NkTISEhwsQW4eHhiI6ORmBg4GUtuZeCp8ViQXNzM5qbm2EymVBR\nUYHS0lKUlpaipKQEJSUlKC4uxtmzZ1FdXY2amho4OTnB19e31UOpVCIyMhL9+vXjjFRERC2YTCY4\nOTkBAH788UfceeedHTYNeHsxsIrM4sWL8dxzz7VaFhsbi0GDBiEkJAQ+Pj4oLy/HL7/8gszMTPTs\n2RNRUVGIiorCgAEDoNFo0K1bNxtVT2JktVpx5swZpKamIjMzE5mZmTh16hSqqqpQV1cHDw8P+Pj4\nwN/fXwixgYGBCA0NRVRUFLp37w6j0Qi9Xg+r1Yrhw4dfc/eD8vJyZGRkCN1XKisrUV5ejoaGBqEV\n19XVFVarFc7OzvDy8oKnpyc8PT0hk8kQHBwMiUQCi8WCtWvXYvv27Thy5AiMRiO6desGe3t72Nvb\nw8HBAXK5HEFBQcIjODgYCoUCgYGBkEql8Pb2Ft54iYjo2uzZswc//fQTtm7dit69e2Px4sVQKpU3\nvQ4GVhHLycnBihUrkJqaitOnT8Pf3x8//fQTQkNDbV0adRIWiwUGgwF6vR4VFRUoKSkRWijz8/OR\nmZmJs2fPIjw8HCqVCmq1Gv3798fgwYOvGFr1ej0yMjKQnp6O33//HUVFRcJH/vb29ggPD0dYWBjc\n3d2FwGlnZ4empibU1dWhrq4OtbW1qKysRHFxMZycnODh4YHQ0FDMnz8fAwYMaLPllYiIbpympiYs\nWrQIH374ITQaDUaNGoXhw4dDo9FctYtYR2FgvUVUVVUhNjYWS5YswcSJE21dDnUhRqMRWVlZQgjd\ns2cP3Nzc8NZbbyE+Ph779u3DoUOHkJ6ejvT0dOj1eqhUKqhUKvTp00eYPSokJAQymaxd577U3eHc\nuXMIDw+32UdRRER0UUNDA/bs2YNdu3Zh3759KCkpQVxcHPr374+4uDjExMQgODi4w0MsA+stYtSo\nUYiKisKSJUvYskQ21dzcjO+//x5PPPEEamtrERcXh2HDhkGj0UClUiE0NJTBkoioi9Dr9Th8+DCO\nHj2Ko0ePIi0tTZhSXqlUYtSoUXjyySevO7t0ycDa0NCAwsJChIWFCf09dTod0tLScPvtt8PX1xfA\nxebvlJQUFBYWora2FjU1NaitrUVtbS3kcjmGDh2KoUOHwtvb+4bWXF9fD19fX7i7u8Pf3x9+fn7w\n9fXFwIEDMXv2bE7HSjZhNBrR3NwsDPhPREQEAI2NjSgsLEReXh5efvllODs7Q6FQCEMLymQyaDQa\n/O1vfxMaOMxmM3Q6HfR6PQIDAy/LNl0ysIaFhQl3Fffq1Qtubm44ffo0VCoVMjIyEBERgZ49e2L3\n7t2IiopCREQEvLy84O3tDS8vL3h5eaGkpAQ///wzDh8+jJ49e2LgwIGIjo6GXC5v9ZDJZB3SImo2\nm4XheHQ6HXQ6HbZs2YLvv/8eMpkMUqlU+CGQSqWQSCRwcHAQbki59LW5uRnTp09HSEjI9V9MIiIi\noqswGo349ddfW41kU1lZiZSUFNTW1sLT0xPnzp2DXq+HTCaDt7c3SktL4eLiAqVSCaVSibCwMCxa\ntKjrBdbvv/8e9957L86fP4/Tp0+jtrYWCQkJcHZ2RlNTEw4ePIjCwkKMHj0aAQEBVz2eyWRCRkYG\nfv31V+Tk5KCioqLVo7GxEbfddhvCwsKgVCoREBAgBN9Lkwd4eXnBw8MDer1eGOrn3Llz0Gq1wg0n\n//twd3eHXC6Hr68v6urqhB+C6upqVFVV4dy5cyguLhZuoikrK0NDQwOAi3f8DR8+/GZcbiIiIqLL\nWK1WZGZmorm5Wcgzl/q9Wq1WaLVa5OfnIz8/H3l5eVi4cGHXC6w387QGg6HVBddqtULXgpqaGuG5\nwWCAj48PAgIC0L17dwQEBMDPzw8mkwlGoxFGoxH19fXCc4PBINzV7efnh9DQULi7u6OsrAwlJSW4\ncOECgoKCEBgYiICAAOHRq1cv3HPPPewHS0RERLeMLtkl4Fa86epKzGYzysrKUFBQgPr6egQFBUGh\nUMDHx4ehlIiIiDoFBlYiIiIiEjVb5jeOe0NEREREosbASkRERESixsBKRERERKLGwEpEREREotax\nk8y2A++eJyIiIrp1+Pj42OzcNgmsHCGAiIiIiK4VuwQQERERkagxsBIRERGRqDGwEhEREZGoMbAS\nERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIR\nERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDGwEhER\nEZGoMbASERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDGwEhEREZGoMbASERER\nkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGR\nqDGwEhEREZGoMbASERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDGwEhEREZGo\nMbASERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagx\nsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDGw\nEhEREZGoMbASERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDGwEhEREZGoMbAS\nERERkagxsBIRERGRqDGwEhEREZGoMbASERERkagxsBIRERGRqDleaYVUKoVer7+ZtRARERFRF+Tj\n44Pq6uorrrezWq3WNlfY2QF45Y/vWubabi2eO7ZzWcvlbS37K/uIqY6rsGvx3OFPDul4le1arm+5\nz4085pUux80+Zlvb/dU6rrbseuu4Icds8d/U0dL6KwB7B/N/V3e7uNyhxfrWzy9u69hymf0f+6DF\nslbP/9injfVtbddy2ysf0/IX9hF/HX9W05/v0zG1d5462rGP5Y99zC3qsDT/d7259VcAsLu0e4tl\nbT63tLGs5XZ/tr6tba/1PO2p41rPc711dNR5/uyc7bnG13OejqyjjeOYWqy/9ONpsly+rOW2LQ9p\n+p+vV1p/rctaLv+zY/6Vff5KHa8AuEIkBcAuAUREREQkcgysRERERCRqDKxEREREJGoMrEREREQk\nagysRERERCRqDKxEREREJGoMrEREREQkagysRERERCRqDKxEREREJGoMrEREREQkagysRERERCRq\nDKxEREREJGoMrEREREQkagysRERERCRqDKxEREREJGoMrEREREQkagysRERERCRqDKxEREREJGoM\nrEREREQkagysRERERCRqDKxEREREJGoMrEREREQkagysRERERCRqDKxEREREJGoMrEREREQkagys\nRERERCRqDKxEREREJGoMrEREREQkagysRERERCRqDKxEREREJGoMrEREREQkagysRERERCRqdlar\n1drmCju7m10LEREREXVBPj4+qK6uvuJ6xyutuEKOJSIiIiK6qdglgIiIiIhEjYGViIiIiESNgZWI\niIiIRK1DA+vOnTuhUqkQGRmJt99+uyMPTX/CYrEgJiYG48aNs3UpXcaCBQvQu3dvhIeHIzExEQ0N\nDbYuqdOaOXMm5HI5VCqVsOzpp59GZGQkIiMjMXbsWFRVVdmwws6presOACtWrIBGo4FKpcJzzz1n\no+o6r5KSEgwdOhQqlQp9+vTBO++8AwCorq7GyJEjoVarMWrUKNTU1Ni40s7lStf9kiVLlsDe3v6q\nNwbRjdNhgbWpqQmPPvoodu7cifT0dHzzzTc4ceJERx2e/sSyZcsQGRnJ0R1ukry8PKxfvx6ZmZnI\nzs6Gg4MDNm3aZOuyOq0ZM2Zg586drZaNGzcOmZmZOHXqFKKjo/Haa6/ZqLrOq63rnpycjF27duG3\n335DRkYGXnjhBRtV13k5OTlh1apVyMjIwG+//Ya1a9ciLS0NCxYswJgxY5Ceno67774bCxYssHWp\nncqVrjtwMcz++OOP6Nmzp42r7Lo6LLCmpqYiKioKgYGBcHR0xKRJk5CcnNxRh6erKC0txY4dO/Dw\nww9zdIebRCqVolu3bqivr4fZbEZDQwPfyG6gIUOGwMfHp9WyO+64A/b2F9/CBg0ahLKyMluU1qm1\ndd3Xrl2L559/Ho6OFweZkclktiitU5PL5YiOjgYAeHh4QK1Wo6ysDDt27MA//vEPAMC0adP4O7aD\ntXXdy8vLAVz8ROd/W1zp5uqwwFpaWgqFQiF8HxQUhNLS0o46PF3FU089hXfffVf45U03nlQqxTPP\nPIPg4GD06NED3t7eGDFihK3L6rLWrFmDe++919ZldAnZ2dnYtWsXbr/9dsTHx+PQoUO2LqlTKyws\nxNGjRzF48GDodDrhDwRfX19otVobV9d5tbzuW7ZskYjUfQAAB0hJREFUQVBQENRqta3L6tI6LOHw\no2jb2L59O/z9/RETE8PW1ZsoPz8fS5cuRWFhIcrLy2E0GrFx40Zbl9Ulvf7663BycsLUqVNtXUqX\n0NzcDIPBgJMnT2L58uWYPHky33tuEKPRiMTERCxbtgyenp62LqfLMBqNuP/++7Fs2TI4ODjgjTfe\nwKuvviqs58+7bXRYYA0KCkJJSYnwfUlJSasWV7oxDh06hK1btyI0NBRJSUnYu3cvHnzwQVuX1ekd\nOXIECQkJkMlkcHR0xMSJE3Hw4EFbl9XlfPbZZ0hOTuYfCzeRQqHAxIkTAQBxcXFwcnJCRUWFjavq\nfEwmE+677z5MnToV48ePBwD4+fmhsrISAKDT6eDv72/LEjulS9d9ypQpGD9+PPLz81FYWAiNRoPQ\n0FCUlpaiX79+bN22gQ4LrHFxccjMzERZWRlMJhO++uor3H333R11eLqCN954AyUlJSgoKMAXX3yB\nO++8E59//rmty+r0wsLCcPjwYZw/fx5WqxV79uxBWFiYrcvqUnbu3Il33nkHW7duhYuLi63L6TLG\njBmDvXv3AgByc3PR0NDA4NTBrFYrZs2ahcjISDz11FPC8tGjR2PDhg0AgA0bNmD06NG2KrFTauu6\nq1QqVFRUoKCgAAUFBQgKCsLx48f5M28DHRZYXVxcsHr1aowaNQoajQYTJ05E3759O+rwdI3YNePm\niIuLQ2JiItRqNcLDw9HU1IR58+bZuqxOKykpCQkJCcjJyYFCocAnn3yCxx9/HEajESNHjkRMTAzm\nzp1r6zI7nUvXPTc3FwqFAuvWrcNjjz2GM2fOIDo6GhMnTsSnn37K/vMd7JdffsGGDRuwb98+xMTE\nICYmBjt37sSrr76K5ORkqNVq/PDDD1i4cKGtS+1U2rruP/zwQ6tt+DvWduys7IxBRERERCLGP4uJ\niIiISNQYWImIiIhI1BhYiYiIiEjUGFiJiIiISNQYWImIiIhI1BhYiajT8/DwEJ7v2LEDffr0aTXR\nSUt79+7FPffcg759+yIuLg5LlixBc3Nzm9vu3r0bffv2hVqthkqlwq5duwAABoNBGBYnJiYGfn5+\nwriOq1atgkajgVqtRmxsLH777bc2jz1z5kzI5XKoVKpWy1955RUEBQW1Gu6IiKiz47BWRNTpSSQS\nGAwG/PTTT/jnP/+J3bt3IzQ09LLtVq9eje3bt2P58uVQKpVoaGjAsmXLcOLECXz11VeXbZ+eno6A\ngAD4+fkhKysLw4cPx9mzZy8bqzE2NhZLly7F4MGDYTQahQC9bds2LF68GCkpKZcd+8CBA/Dw8MCD\nDz6IjIwMYfmrr74KiUSCp59++novCxHRLYMtrETUJezfvx+zZ89GcnJym2H1999/x9dff41t27ZB\nqVQCANzc3PDvf/8b0dHRwgxDLanVavj5+QEAoqKi0NzcjMbGxlbb5ObmQqvVYvDgwQBat/YajUYE\nBAS0We+QIUPg4+PT5jq2MxBRV+No6wKIiG60xsZGTJgwASkpKejdu3eb26xbtw4vvfQS7Ozs8PDD\nDyM1NRVjx47F+fPnsWjRIiQlJWHatGlXPMc333wDjUYDV1fXVsu/+OILTJ48udWyVatW4b333kN9\nfT0OHTrU7tezcuVKrF27Fv369cPy5cshlUrbfQwiolsJW1iJqNNzcnLCoEGDsHbt2ituk56ejgED\nBuDbb7+Fp6cnMjIy4O/vj7q6OkgkEtTX119x31OnTuGFF17AmjVrLlv35ZdfIikpqdWyuXPnIi8v\nD++99x5mzpzZrtcyb9485Ofn49SpU1AqlXjiiSfatT8R0a2IgZWIOj17e3t89dVXOHLkCN58880r\nbufg4ICcnByMGjUKAPD3v/8dVqsVFy5cuOIc4qWlpZgwYQLWr19/WVeDtLQ0mM1mxMTEtLnvpEmT\ncPTo0Xa9Fl9fX9jZ2cHOzg5z5sxp9/5ERLciBlYi6hJcXFyQnJyMjRs34pNPPrlsfXR0NFJTU9G7\nd2/s3r0bAIS7/t977z3cd999l+1TU1ODMWPG4K233kJ8fPxl6zdt2oQpU6a0WlZYWCg8T05ORkRE\nRLteh1arFZ5/++23iIqKatf+RES3IgZWIur0LrWO+vj4YOfOnXjttdewffv2VttMnz4dr732GsaP\nH4+amhqoVCpUVFTg2LFjsFqtmDdv3mXH/eCDD5Cfn4+FCxcKw0zpdDph/ddff31Zd4AlS5ZArVYj\nOjoaixcvxueffw4AKC8vx5gxY4TtkpKSkJCQgNzcXCgUCqxbtw4A8Mwzz0Cj0SAiIgLJyclYsWJF\nx1wkIiIR47BWRER/WLJkCX799Ve8//77UCgUaGpqwvbt2zFw4EAEBgbaujwioi6LgZWIqIUffvgB\ny5Ytg1arhaurK+6//348/vjjcHBwsHVpRERdFgMrEREREYka+7ASERERkagxsBIRERGRqDGwEhER\nEZGoMbASERERkagxsBIRERGRqP0/OZ8e7fY5cOUAAAAASUVORK5CYII=\n", "text": [ "" ] }, { "output_type": "stream", "stream": "stdout", "text": [ "corrected pot. temp / (degree_celsius) (latitude: 90; longitude: 120)\n", " Dimension coordinates:\n", " latitude x -\n", " longitude - x\n", " Auxiliary coordinates:\n", " time x x\n", " Scalar coordinates:\n", " depth: 25.0 metres, bound=(0.0, 50.0) metres\n", " Attributes:\n", " _fillvalue: 99999.0\n", " c_format: %9.3f\n", " comment: corrected value\n", " fortran_format: f9.3\n", " resolution: 0.001\n", " Cell methods:\n", " mean: i_depth, i_lat, i_lon\n" ] } ], "prompt_number": 12 } ], "metadata": {} } ] }