{ "cells": [ { "cell_type": "markdown", "metadata": { "collapsed": true }, "source": [ "# 04 - Full waveform inversion with Devito and Dask" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Introduction\n", "\n", "In this tutorial we show how [Devito](http://www.opesci.org/devito-public) and [scipy.optimize.minimize](https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.minimize.html) are used with [Dask](https://dask.pydata.org/en/latest/#dask) to perform [full waveform inversion](https://www.slim.eos.ubc.ca/research/inversion) (FWI) on distributed memory parallel computers." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## scipy.optimize.minimize \n", "\n", "In this tutorial we use [scipy.optimize.minimize](https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.minimize.html) to solve the FWI gradient based minimization problem rather than the simple grdient decent algorithm in the previous tutorial.\n", "\n", "```python\n", "scipy.optimize.minimize(fun, x0, args=(), method=None, jac=None, hess=None, hessp=None, bounds=None, constraints=(), tol=None, callback=None, options=None)\n", "```\n", "\n", "> Minimization of scalar function of one or more variables.\n", ">\n", "> In general, the optimization problems are of the form:\n", ">\n", "> minimize f(x) subject to\n", ">\n", "> g_i(x) >= 0, i = 1,...,m\n", "> h_j(x) = 0, j = 1,...,p\n", "> where x is a vector of one or more variables. g_i(x) are the inequality constraints. h_j(x) are the equality constrains.\n", "\n", "[scipy.optimize.minimize](https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.minimize.html) provides a wide variety of methods for solving minimization problems depending on the context. Here we are going to focus on using L-BFGS via [scipy.optimize.minimize(method=’L-BFGS-B’)](https://docs.scipy.org/doc/scipy/reference/optimize.minimize-lbfgsb.html#optimize-minimize-lbfgsb)\n", "\n", "```python\n", "scipy.optimize.minimize(fun, x0, args=(), method='L-BFGS-B', jac=None, bounds=None, tol=None, callback=None, options={'disp': None, 'maxls': 20, 'iprint': -1, 'gtol': 1e-05, 'eps': 1e-08, 'maxiter': 15000, 'ftol': 2.220446049250313e-09, 'maxcor': 10, 'maxfun': 15000})```\n", "\n", "The argument `fun` is a callable function that returns the misfit between the simulated and the observed data. If `jac` is a Boolean and is `True`, `fun` is assumed to return the gradient along with the objective function - as is our case when applying the adjoint-state method." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## What is Dask?\n", "\n", "> [Dask](https://dask.pydata.org/en/latest/#dask) is a flexible parallel computing library for analytic computing.\n", ">\n", "> Dask is composed of two components:\n", ">\n", "> * Dynamic task scheduling optimized for computation...\n", "> * “Big Data” collections like parallel arrays, dataframes, and lists that extend common interfaces like NumPy, Pandas, or Python iterators to larger-than-memory or distributed environments. These parallel collections run on top of the dynamic task schedulers.\n", ">\n", "> Dask emphasizes the following virtues:\n", "> \n", "> * Familiar: Provides parallelized NumPy array and Pandas DataFrame objects\n", "> * Flexible: Provides a task scheduling interface for more custom workloads and integration with other projects.\n", "> * Native: Enables distributed computing in Pure Python with access to the PyData stack.\n", "> * Fast: Operates with low overhead, low latency, and minimal serialization necessary for fast numerical algorithms\n", "> * Scales up: Runs resiliently on clusters with 1000s of cores\n", "> * Scales down: Trivial to set up and run on a laptop in a single process\n", "> * Responsive: Designed with interactive computing in mind it provides rapid feedback and diagnostics to aid humans\n", "\n", "**We are going to use it here to parallelise the computation of the functional and gradient as this is the vast bulk of the computational expense of FWI and it is trivially parallel over data shots.**" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Setting up (synthetic) data\n", "In a real world scenario we work with collected seismic data; for the tutorial we know what the actual solution is and we are using the workers to also generate the synthetic data." ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "#NBVAL_IGNORE_OUTPUT\n", "\n", "# Set up inversion parameters.\n", "param = {'t0': 0.,\n", " 'tn': 1000., # Simulation last 1 second (1000 ms)\n", " 'f0': 0.010, # Source peak frequency is 10Hz (0.010 kHz)\n", " 'nshots': 5, # Number of shots to create gradient from\n", " 'm_bounds': (0.08, 0.25), # Set the min and max slowness\n", " 'shape': (101, 101), # Number of grid points (nx, nz).\n", " 'spacing': (10., 10.), # Grid spacing in m. The domain size is now 1km by 1km.\n", " 'origin': (0, 0), # Need origin to define relative source and receiver locations.\n", " 'nbl': 40} # nbl thickness.\n", "\n", "import numpy as np\n", "\n", "import scipy\n", "from scipy import signal, optimize\n", "\n", "from devito import Grid\n", "\n", "from distributed import Client, LocalCluster, wait\n", "\n", "import cloudpickle as pickle\n", "\n", "# Import acoustic solver, source and receiver modules.\n", "from examples.seismic import Model, demo_model, AcquisitionGeometry, Receiver\n", "from examples.seismic.acoustic import AcousticWaveSolver\n", "from examples.seismic import AcquisitionGeometry\n", "\n", "# Import convenience function for plotting results\n", "from examples.seismic import plot_image\n", "\n", "def get_true_model():\n", " ''' Define the test phantom; in this case we are using\n", " a simple circle so we can easily see what is going on.\n", " '''\n", " return demo_model('circle-isotropic', vp=3.0, vp_background=2.5, \n", " origin=param['origin'], shape=param['shape'],\n", " spacing=param['spacing'], nbl=param['nbl'])\n", "\n", "def get_initial_model():\n", " '''The initial guess for the subsurface model.\n", " '''\n", " # Make sure both model are on the same grid\n", " grid = get_true_model().grid\n", " return demo_model('circle-isotropic', vp=2.5, vp_background=2.5, \n", " origin=param['origin'], shape=param['shape'],\n", " spacing=param['spacing'], nbl=param['nbl'],\n", " grid=grid)\n", "\n", "def wrap_model(x, astype=None):\n", " '''Wrap a flat array as a subsurface model.\n", " '''\n", " model = get_initial_model()\n", " if astype:\n", " model.vp = x.astype(astype).reshape(model.vp.data.shape)\n", " else:\n", " model.vp = x.reshape(model.vp.data.shape)\n", " return model\n", "\n", "def load_model(filename):\n", " \"\"\" Returns the current model. This is used by the\n", " worker to get the current model.\n", " \"\"\"\n", " pkl = pickle.load(open(filename, \"rb\"))\n", " \n", " return pkl['model']\n", "\n", "def dump_model(filename, model):\n", " ''' Dump model to disk.\n", " '''\n", " pickle.dump({'model':model}, open(filename, \"wb\"))\n", " \n", "def load_shot_data(shot_id, dt):\n", " ''' Load shot data from disk, resampling to the model time step.\n", " '''\n", " pkl = pickle.load(open(\"shot_%d.p\"%shot_id, \"rb\"))\n", " \n", " return pkl['geometry'].resample(dt), pkl['rec'].resample(dt)\n", "\n", "def dump_shot_data(shot_id, rec, geometry):\n", " ''' Dump shot data to disk.\n", " '''\n", " pickle.dump({'rec':rec, 'geometry': geometry}, open('shot_%d.p'%shot_id, \"wb\"))\n", " \n", "def generate_shotdata_i(param):\n", " \"\"\" Inversion crime alert! Here the worker is creating the\n", " 'observed' data using the real model. For a real case\n", " the worker would be reading seismic data from disk.\n", " \"\"\"\n", " true_model = get_true_model()\n", " shot_id = param['shot_id']\n", " \n", " src_coordinates = np.empty((1, len(param['shape'])))\n", " src_coordinates[0, :] = [30, param['shot_id']*1000./(param['nshots']-1)]\n", " \n", " # Number of receiver locations per shot.\n", " nreceivers = 101\n", "\n", " # Set up receiver data and geometry.\n", " rec_coordinates = np.empty((nreceivers, len(param['shape'])))\n", " rec_coordinates[:, 1] = np.linspace(0, true_model.domain_size[0], num=nreceivers)\n", " rec_coordinates[:, 0] = 980. # 20m from the right end\n", "\n", " # Geometry \n", " geometry = AcquisitionGeometry(true_model, rec_coordinates, src_coordinates,\n", " param['t0'], param['tn'], src_type='Ricker',\n", " f0=param['f0'])\n", " # Set up solver.\n", " solver = AcousticWaveSolver(true_model, geometry, space_order=4)\n", "\n", " # Generate synthetic receiver data from true model.\n", " true_d, _, _ = solver.forward(vp=true_model.vp)\n", "\n", " dump_shot_data(shot_id, true_d, geometry)\n", "\n", "def generate_shotdata(param):\n", " # Define work list\n", " work = [dict(param) for i in range(param['nshots'])]\n", " for i in range(param['nshots']):\n", " work[i]['shot_id'] = i\n", " generate_shotdata_i(work[i])\n", " \n", " # Map worklist to cluster\n", " futures = client.map(generate_shotdata_i, work)\n", "\n", " # Wait for all futures\n", " wait(futures)" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "Operator `Forward` run in 0.02 s\n", "Operator `Forward` run in 0.06 s\n", "Operator `Forward` run in 0.02 s\n", "Operator `Forward` run in 0.03 s\n", "Operator `Forward` run in 0.04 s\n" ] } ], "source": [ "#NBVAL_IGNORE_OUTPUT\n", "\n", "# Start Dask cluster\n", "cluster = LocalCluster(n_workers=2, death_timeout=600)\n", "client = Client(cluster)\n", "\n", "# Generate shot data.\n", "generate_shotdata(param)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Dask specifics\n", "\n", "Previously we defined a function to calculate the individual contribution to the functional and gradient for each shot, which was then used in a loop over all shots. However, when using distributed frameworks such as Dask we instead think in terms of creating a worklist which gets *mapped* onto the worker pool. The sum reduction is also performed in parallel. For now however we assume that the scipy.optimize.minimize itself is running on the *master* process; this is a reasonable simplification because the computational cost of calculating (f, g) far exceeds the other compute costs." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Because we want to be able to use standard reduction operators such as sum on (f, g) we first define it as a type so that we can define the `__add__` (and `__rand__` method)." ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "# Define a type to store the functional and gradient.\n", "class fg_pair:\n", " def __init__(self, f, g):\n", " self.f = f\n", " self.g = g\n", " \n", " def __add__(self, other):\n", " f = self.f + other.f\n", " g = self.g + other.g\n", " \n", " return fg_pair(f, g)\n", " \n", " def __radd__(self, other):\n", " if other == 0:\n", " return self\n", " else:\n", " return self.__add__(other)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Create operators for gradient based inversion\n", "To perform the inversion we are going to use [scipy.optimize.minimize(method=’L-BFGS-B’)](https://docs.scipy.org/doc/scipy/reference/optimize.minimize-lbfgsb.html#optimize-minimize-lbfgsb).\n", "\n", "First we define the functional, ```f```, and gradient, ```g```, operator (i.e. the function ```fun```) for a single shot of data. This is the work that is going to be performed by the worker on a unit of data." ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [], "source": [ "from devito import Function\n", "\n", "# Create FWI gradient kernel for a single shot\n", "def fwi_gradient_i(param):\n", " # Load the current model and the shot data for this worker.\n", " # Note, unlike the serial example the model is not passed in\n", " # as an argument. Broadcasting large datasets is considered\n", " # a programming anti-pattern and at the time of writing it\n", " # it only worked relaiably with Dask master. Therefore, the\n", " # the model is communicated via a file.\n", " model0 = load_model(param['model'])\n", " \n", " dt = model0.critical_dt\n", "\n", " geometry, rec = load_shot_data(param['shot_id'], dt)\n", " geometry.model = model0\n", " # Set up solver.\n", " solver = AcousticWaveSolver(model0, geometry, space_order=4)\n", "\n", " # Compute simulated data and full forward wavefield u0\n", " d, u0, _ = solver.forward(save=True)\n", " \n", " # Compute the data misfit (residual) and objective function\n", " residual = Receiver(name='rec', grid=model0.grid,\n", " time_range=geometry.time_axis,\n", " coordinates=geometry.rec_positions)\n", "\n", " residual.data[:] = d.data[:residual.shape[0], :] - rec.data[:residual.shape[0], :]\n", " f = .5*np.linalg.norm(residual.data.flatten())**2\n", " \n", " # Compute gradient using the adjoint-state method. Note, this\n", " # backpropagates the data misfit through the model.\n", " grad = Function(name=\"grad\", grid=model0.grid)\n", " solver.gradient(rec=residual, u=u0, grad=grad)\n", " \n", " # Copying here to avoid a (probably overzealous) destructor deleting\n", " # the gradient before Dask has had a chance to communicate it.\n", " g = np.array(grad.data[:])\n", " \n", " # return the objective functional and gradient.\n", " return fg_pair(f, g)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Define the global functional-gradient operator. This does the following:\n", "* Maps the worklist (shots) to the workers so that the invidual contributions to (f, g) are computed.\n", "* Sum individual contributions to (f, g) and returns the result." ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "def fwi_gradient(model, param):\n", " # Dump a copy of the current model for the workers\n", " # to pick up when they are ready.\n", " param['model'] = \"model_0.p\"\n", " dump_model(param['model'], wrap_model(model))\n", "\n", " # Define work list\n", " work = [dict(param) for i in range(param['nshots'])]\n", " for i in range(param['nshots']):\n", " work[i]['shot_id'] = i\n", " \n", " # Distribute worklist to workers.\n", " fgi = client.map(fwi_gradient_i, work, retries=1)\n", " \n", " # Perform reduction.\n", " fg = client.submit(sum, fgi).result()\n", " \n", " # L-BFGS in scipy expects a flat array in 64-bit floats.\n", " return fg.f, -fg.g.flatten().astype(np.float64)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## FWI with L-BFGS-B\n", "Equipped with a function to calculate the functional and gradient, we are finally ready to define the optimization function." ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [], "source": [ "from scipy import optimize\n", "\n", "# Define bounding box constraints on the solution.\n", "def apply_box_constraint(vp):\n", " # Maximum possible 'realistic' velocity is 3.5 km/sec\n", " # Minimum possible 'realistic' velocity is 2 km/sec\n", " return np.clip(vp, 2.0, 3.5)\n", "\n", "# Many optimization methods in scipy.optimize.minimize accept a callback\n", "# function that can operate on the solution after every iteration. Here\n", "# we use this to apply box constraints and to monitor the true relative\n", "# solution error.\n", "relative_error = []\n", "def fwi_callbacks(x):\n", " # Apply boundary constraint\n", " x.data[:] = apply_box_constraint(x)\n", " \n", " # Calculate true relative error\n", " true_x = get_true_model().vp.data.flatten()\n", " relative_error.append(np.linalg.norm((x-true_x)/true_x))\n", "\n", "def fwi(model, param, ftol=0.1, maxiter=5):\n", " result = optimize.minimize(fwi_gradient,\n", " model.vp.data.flatten().astype(np.float64),\n", " args=(param, ), method='L-BFGS-B', jac=True,\n", " callback=fwi_callbacks,\n", " options={'ftol':ftol,\n", " 'maxiter':maxiter,\n", " 'disp':True})\n", "\n", " return result" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We now apply our FWI function and have a look at the result." ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ " fun: 211.63605994835484\n", " hess_inv: <32761x32761 LbfgsInvHessProduct with dtype=float64>\n", " jac: array([-6.69854363e-12, -3.05454377e-11, -7.75552331e-11, ...,\n", " -1.31146288e-10, -5.26459223e-11, -1.16629224e-11])\n", " message: b'STOP: TOTAL NO. of ITERATIONS REACHED LIMIT'\n", " nfev: 6\n", " nit: 5\n", " status: 1\n", " success: False\n", " x: array([2.5, 2.5, 2.5, ..., 2.5, 2.5, 2.5])\n" ] } ], "source": [ "#NBVAL_IGNORE_OUTPUT\n", "\n", "model0 = get_initial_model()\n", "\n", "# Baby steps\n", "result = fwi(model0, param)\n", "\n", "# Print out results of optimizer.\n", "print(result)" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAUYAAAEKCAYAAABuTfznAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDMuMC4yLCBodHRwOi8vbWF0cGxvdGxpYi5vcmcvOIA7rQAAIABJREFUeJzsvX2cJFlZ5/t9OqIys6qypovudrodZrCVF/Fl9bKM7PUO9yOIsKN3kevL6oCyiq6oK64oei+gF/3gru+is4KLI2/iIq4rILO7IzC+Ir4gDKIDjODsiDBvzHTPVE9VV2VmZdS5f5xzIp44eSIysjrrbSZ+/cnOyBMnIk5ERfzieTvPI8YYWrRo0aJFgWMHPYAWLVq0OGxoibFFixYtArTE2KJFixYBWmJs0aJFiwAtMbZo0aJFgJYYW7Ro0SLAnhKjiFwrIh8XkdtF5KWR9V0R+a9u/ftF5OxejqdFixYtmmDPiFFEEuA1wFcDXwg8V0S+MOj2ncCDxpjHAb8E/OxejadFixYtmmIvJcanALcbY+4wxoyA3waeE/R5DvAbbvl3gWeIiOzhmFq0aHGAEJGrROSPReRjIvJREfmBSJ/jIvLfReRvXZ8XqHXfJiL/4D7ftlfjTPdqx8CjgU+r33cC/6KqjzFmLCIXgJPAOd1JRF4IvBCA7vKTOfPEPRpyixYtAPinW84ZYz4L4HEiZrPhZvfAu40x19Z0GQMvMcZ8SERWgFtE5GZjzMdUn+8DPmaMebaIfBbwcRF5C9AHfhy4GjBu2xuNMQ/OenrTsJfEODcYY24AbgCQs1cb/r8PHuyAWrR4uOPfyj/5xU3guxtu9hNwqm69MeYe4B63vC4it2EFJE2MBlhx2mMfeABLqP8SuNkY8wCAiNwMXAu8teHwGmMvifEu4Cr1+0rXFutzp4ikwHHg/B6OqUWLFjNCmIkoTomIllxucILN5H6ts/VJwPuDVa8GbgTuBlaAbzbG7IhITAt9dPOhNcdeEuMHgMeLyOdiCfA64HlBnxuBbwP+EvhG4I9Mm9WiRYtDhWPAYvPu54wxV0/rJCJ94G3Ai40xDwWr/yXwYeArgccCN4vInzUfwqVjz5wvxpgx8CLg3cBtwO8YYz4qIq8Uka913V4PnBSR24EfAiZCelq0aHGwEGCh4afR/kQWsKT4FmPM2yNdXgC83VjcDvwj8ESaaaFzwZ7aGI0xNwE3BW2vUMsD4F/v5RhatGhxaZhRla7fl7Ubvh64zRjzqopunwKeAfyZiJwGPh+4A7gd+CkReZTr9yzgZXMaWglHwvnSokWLg4OXGOeEa4DnA7eKyIdd28uBxwAYY14L/CTwJhG51R3+/zXGnAMQkZ/EmukAXukdMfNGS4wtWrSoxTwlRmPM+9wu6/rcjZUGY+veALxhTsOpREuMLVq0qMWcJcYjgZYYW7RoUYsZvdIPC7TE2KJFi1q0EmOLFi1aRPBII4pH2vm2aNFiRrQSY4sWLVoEmKdX+qjgkXa+LVq0mBGt86VFixYtArSqdIsWLVoEaFXpFi1atAjQSowtWrRoEaCVGFu0aNEiQCsxtmjRokUAofVKt2jRokUJAiw0ZYrxXo5k/9ASY4sWLWohAmlLjC1atGhRQAQWkoMexf6iJcYWLeaFqqfpiEtRM0mMU/clVwFvBk5jy6TeYIy5PujzI8C3uJ8p8AXAZxljHhCRTwLrQAaMmxTe2g1aYjwI+Kt+xB+YFgH83zN8qsLfR+zvLgIL3bntbgy8xBjzIRFZAW4RkZuNMXldaWPMzwM/b48tzwZ+MChh8HRf6mCv0BLjfuCIPxgtZsSY+ifrqN0PcwxkNMbcA9zjltdF5DZsbeiPVWzyXOCt8zl6c+xZ+dQWLVo8TOCJsclnlt2KnAWeBLy/Yv0ScC221KqHAd4jIreIyAtnO2JztBLjXmFe9qbdSBPz+Kse9J1x2KWoadDjn3Ytj4Jppfn9cEpEPqh+32CMuSHsJCJ9LOG92BjzUMW+ng38eaBGP9UYc5eIXA7cLCJ/b4x5b+PRNcRB3/4PP+yWEPf7odjLv3zdvpue527HdynXca9eQk1JMt3lGPYaAjT3Sp+b5hARkQUsKb7FGPP2mq7XEajRxpi73Pd9IvIO4ClAS4yHGk1tR7O21x1jWvu0dbP0mcd28yDNuu0vZf9626ZjqesXG0uVg+ZSxrDXmKONUUQEeD1wmzHmVTX9jgNfAXyralsGjjnb5DK2xOor5zOyMlpinAfqHoDY71lv+Glezqq2uvZL7bub/rvd3yxS5m7IJLbNPBwk4TZVpFf3tzsM5CjA/LzS1wDPB24VkQ+7tpcDjwEwxrzWtX0d8B5jzEW17WngHZZbSYHfMsa8a24jU2iJsUWLFvWYr1f6fW6P0/q9CXhT0HYH8KXzGUk9WmK8VDSRLMYV33X7SSPtsbbY70vFbva3XzbLS1GHZ3FyNJX2mowrPG7snjnMUuMjMO/YI+x054w6UgyXq25uTXYh8YUhELHfTXApNs06HMTdMy+Pe0yVDduqzCExkqwjsNh+q+6dKjPJQZNjOyWwxcyIPUAhGcZufB3/1WOSEKuIcDeEWDWWOlvoLPuP4aDvriYOqTqCjP2uI7+w/2488NMcSAeBVmJsMROqVKoYMXposvNk2AvapxFj+HuaZOMftli/aerlXoW/7Hes5UHd6VXXPLz2VVIowfJBSI4tMba4JISkGN7gPfXRv3crGe72IdEPZhNCvBSPet0YLrX/NM9uE0lslnHECEz/raoky2lSeZ16HXsZ7jc5ztcrfSTQEmOLFi3q0UqMLXYFbVOqUqN7QJ9JiTHcB8Ag2LbJX2lW6TJmBwvPo0nfqt9Nx9W0X53XVh9/Fu/ubqTWJjbJKoRaRFOnzEE/pS0xtpgZ+mavU6H7FMSYRvoPaEZAer/h75gKPquqqMdF5Pc8VOwmTpEm7dPIJSShpqRTZTuuWhcboz+mvnb+hefbYmp5EyLcb3V6timBDwvsKTGKyLXA9djL+jpjzM8E638I+LfYP/P9wHcYY/5pL8c0V1R5o6FwqsCkbdE/IBvuM1Afvd+qhzFGjDGHTex3uH3M6F+13ESSjP0OMSsxVkl7MVIM91HVPm0sdeQ4Dfq6+5ce2L9/+Df2yzG7b5VzZr/RSozzg4gkwGuAZwJ3Ah8QkRt1Qkrgb4CrjTGbIvK9wM8B37xXY9oTxNRm/1B4IuxTJsU1t7zmPp4UqyTPJgSpj6kJMPR8676h00fvcxrZ1Y1tGoHUnUeMxJpIezFJMdzPpd7tTaVF/zu8xoPgO9z3YSWf1vkyVzwFuN1N40FEfht4DiohpTHmj1X/v0JNGD/0qLOpeSLqqTYoSNHnHvakGCPCQfANcRIOj6sfRE2KvYr2kCz1eOsQk3BiaCI91tntYn2q2qoIdF6ouy762Po6aztyj2oCj+2nyizgv/dLnW4lxrni0cCn1e87gX9R0/87gd+PrXAJKW1SyhOPmdPwWrRo0QgtMR4MRORbgauxaYYm4BJd3gAgZ682+zeyCkyT4GK2vTGFTdFvoz3V2hY1UH1h0i4Vs1NVHbdKYtRSI2o5VK+n2QPrVMCYutu0z7Q7s6pfnVRdtY8YYqqxH+s0u2SVrRAm7cy6/2HGYR/fnLGXp3sXcJX6faVrK0FEvgr4UeArjDHDPRzP3qPKdqdV5RQ45drPYK9KD0uCXsU+5z7eBukJ0jtqUsoEMA4+4XgGTCfGXtAvZo+sOucmaPrwV9k86/Y7rW03jpTYyyH2rfv5v5N/yfm/U2hSmXYtNKHW2VL3C61Xeq74APB4EflcLCFeBzxPdxCRJwG/BlxrjLlvD8eyd6gzuntoUlwFztrmM9fcwTX8Bec4xSc5y6c/Y98jO59chnspyFETY+jJ9vvXtkqNkJR1f+0l7VEmR/0wVxF+HVnGpKZZbGK7uTOnecmbOIhix9VEGMahasfave7b2471i8dvr9tCqfGwolWl5wdjzFhEXgS8G/u+eYMx5qMi8krgg8aYG7ElEvvAf3PJJz9ljPnavRrTJaPJw6TJI+zfA85A/3+7H4Cv4g/5GV7Ka/g+ALLT9rV8N1ew01u2JOqlRigToydHKFRv7bDRxw8Jc4NJaUSTpe6vVezwIY85HKrWN3HU1KHOYaHbQ/Lbjdc8dk6a1PqurY/9G61C78wDDD58wrZ7J5t/eTV5yTQdW9M+80TrlZ4vjDE3ATcFba9Qy1+1l8ffc1Q9QBCfveIepP7yOgBrrPKf+PfczRVkJCyyCcDqqTXWwJKjf/hgkhhDSVITpB9DlZrtt4sRl36YPVHqb38+VaQYXgu/fCle7Gnq8CzEWIUqNb4XfDwxKlI8dfw8d66eKG8XC8OqI8ZLfXnsFR6BEmNbPrVFixb1mGP5VBG5SkT+WEQ+JiIfFZEfiPT5ERH5sPt8REQyETnh1l0rIh8XkdtF5KXzOsUQj7D3wJxRJTFpddYjnCcN3MfljOgwpMsmi3l7JxmysrrOZpqx3VuE3kKxjzqJcYPCkQPlGTWhPSuUpLy0otVpKNscdbt20sSkxyppyK+rsgOGttkY6tTjuraYFz88plaZw5lLffXBfh9bvcjK8Q0r7WtpWh9zoNrD48TO7bA9lfOVGMfAS4wxHxKRFeAWEblZT/wwxvw81syGiDwb+EFjzAMNJ43MBYftT3B0UOWZjJFRYJ8bZdZgs56sMHTGmxEdMrejlIxOdwTAMB0zSJfcvruQyqQRv+rN7VVr3xaq2DBJJFqdhrLNUROqJszQ/hheo3C5ytPaRI2ettykrx5PSOp1xBg6W/rbLPU36TKkywh6224/C8Ux9X0QXqPYy+OwYk5eaWPMPcA9bnldRG7DxjxXkdtzKUqoTp00Mi+0xDgLYjao8AEfMxmv6O2L7nc2tnfZMOmSkE0cJiEjYUySJqRpxkLPkuQ2wFg9XTEJRI/Jk+IG8YewilhCB5K2O4a2R21/9O2xu2oaSYbkWOU4qpMCq6TQqrHEpN1wthAUZKi/gWO9EZ3eiA6jyb+jv2bhCykmKcYcVYcJs0mMp0Tkg+r3DS4OeXK3ImeBJwHvr1i/BFwLvMg1zTppZNdoiXG3qFKFvISgnSB9SlJakk6SIUDiOmR1r+fUWKkxlrIsfNC9Wh2So5ZoQwKKEU6T2En9UIdqdj726tOaOGbsGGF7bJs6THuRhBJwX7VrNbpv/7Cd3pBOMrQvMcaQ/10rJEZv6qiLazxspAizeqXPGWOunrpLkT7wNuDFxpiHKro9G/hzY8wDjY8+J7TEOE/4ByEMp/HrAtgHKlNdyoSYjRPG4ySXMBmr9fqB7hN/0HUsnY5ZjGX0CeMgq4gnTJ2lpUaC31Uqdh2qiDcmQTaxU1b91ipzSIxhNiRtH+6ZXILv9kakZKShtBhKxBtMIgx9ip3PYcGcvdIisoAlxbcYY95e0/U6CjUaGk4amQdar3SLFi3qMV+vtACvB24zxryqpt9x7BThd6rmfNKIiHSwxHnj7Cc0Ha3EWIU6w3jVOh1gHRrdlaSTJJnbTRZVnzMSK4eME7Jx5E+kJRx/fC/96BjEAYUn23+jlkOJMZbRp+oa+G18u7YxhuP0y3UPT5VkGI6nyvZYhZgqH0qLMSdLlfOlN6TbszNXO92yffGYU6V3/DHCMcYcPbPiIKTL+U4JvAZ4PnCriHzYtb0ceAyAMea1ru3rgPcYYy76DasmjcxtZAotMTZF7EpVPaxVNjH3yTJ7l42TJPdEA+XlLHInppn95Kq1GtRYykSCWg5JMBYMrvuE+6giIe2UiYUD1YX3xPZVZ8NsokqHqLJvNvE+R9u3WeiNchtxl2FuDimp06EqHQv2P6xqcwxzVKWNMe9ze5zW703AmyLtE5NG9gItMc4DTaQYRVqbGzZmceV410mN2YTDJUkykjSj2xsyHsdf16E0mY0TdsZJQZyDhUJqDKVDTY5hW5UkWXduHiFRxggpbNf7CK9lzPZZN5a6OzrmfY5JcVFiNBzrjej2hnkolSdEL/V7wtzRUl0si85uHEZVnuz9gDC7rfiIoyXGGGb1nsa2jz34G8AaDNZWAFjvb0ECHYalB6yDfei63SHjrlaxqweWZYlTvZOcSLNxyvagY4nSkyQUqrUmRq1ux6TLKpU2hO8XPtCeFKo81dOk72kS47S/WWMCDNb5/r0hnd6QJFVE6F5qHklq27f1eEJpscnLs8qJdFBos+u0mBkxO1LoBYYiscA5oG/DOdZ6qySnc4uiDRIGErYmbI7F4RzplVTwhCxJrGreTXM1PBsnZL0hw0GX7TSD1MVcpDKZ4EDbJzVh6t/6nKZJbjGCbGJfqyPI2Hq9XWyfMedAE1W6ZHu0KnSa2sD74u8UUaP9/mP20di4dXvVOR30U/oInCv9CDvdFi1a7AqPMKZ4hJ3uJWLam99LGaHjQanR3EseOLzTX2arv8nS8iYdRiy57DpJYHfU32Nl2RrSiY/TCZg+/jFJx2Rpwk7P5wHuErV/xyS6cEqhVhHrnE3+u8rGGB4zRJ2EWHXtq9piEmNUMgw+qU0Wf6w3IknHdHqjUuxpl/LvCYRB9Pqe0OOJOeqqzukg0KrSLXYFfRU1MYbODt/XpxE7BcNBF5atWubTjnlVDXAEWMynHiuS9Oq0b89IGQ07OSGOXbhP7pQZOFXae7Bj5+HHH4YExRDzuE5TeZsSY2yfMdSFBIVqtD6nkBwhIEYDPiynN6TbG9FJ7JzoLrZdq9El55lWo2NB/lU26DrTxEGiVaVbNELVw66lqljIxphyPOGgkOpCO5WXGMdOfgQY0mGLpTwbj09GMRp0GA66k17pcUopjKfKsRGOMTy30CbnEQvV0deniW3NH6NJW9heF4IDZfKbcKYE62GCFHt996LqjXLbYodhLiGG0mIeJaBJMeaVDom6SvJOg230fvYTbaLaFpWoCw+pUhf9Nv4BSSk/LAPYGXQY0mGTpfwhG7m7MCNhk6U8JdkWS2xmS2xuLDIadNnZ8Fl3ZFK1pea7at00L3AYIqIf3Jg62MRpovelv8O2mHc71ieUAGOEWRuvWJBi1039s6Q4dAniLDnaQ5ancxZTNykTY90YNfz1Cs/1oNFKjC0aoc5uFkoBUE0GA2BtgfP9U4xWuywlhY0RrIQ4yrp53ONgYwk2evHSBjFirEOsTxWBVkk0sX1VeZNjs2Q0qoiuTvWOqcd1EmNUZUYliyiToo9XtKQ4cjbFQn0uSYukVlr356r/JnXqvr4m4fU/DB5paImxRYsWLSbQEmOLmRCzFfnvmL1Ie3QhL5e6vXYZD/Qv44FYOirv0fbbxaoENpUWZ/lr10mBejmmKscCwpsEh1c5TmLzr6ukwyY2xglV2pUqD6RFHyXgVeiEzE0DLCTGLI8hUM4tPYtIn6uWXOvMD4fxqWy90i0mUKX6NfUuhg+xr/o3wIbvxB7cmLFd53rUHuHdeDRn/ctX2SSr5laHoSrTyLvKoxx6yYn0gWpnim7T+RWDxBAAC71RiRS9LdHbFVNFg+VL40KnBi4EKiTFaeaBaeaFg0YrMbYAJg3gs24bPhCa6Hw8o+6bUq4G6BOihqEyIRnVSahVY6tC3Z1QZ3sMibEqw1CdrVITRyjtVZ1bjEinSoax9m2VX3E44WgBakmxmLeUxiV4fX6hrTP20tD9w+tzUGi90i1yzEqIfpvYDe4lFS/tORU6J8gecCo4po6H9GjiBQ7HrttiKn+sX+y89HHCY4aSYeh8gOYSo0evpq/uH5JdP2irkxj723liCIBOb5Q7WlJFgnWkOCbBJh/rxJ1hIdnrv0ETh1TVE7qb+3O3aCXGFrtCjFS8FOjb/cNyDlupAorEDV5aPOXaT7nf/iHXUlGd7S9GkCGBxdrD/cWWq/YdkxinEWMMnkAGQds0abHK0xxrC8JybM2WIR0flpMMcztix81qgclkER5eWhy6clgTEqMeqz9uVdRCrH8M+0mIHkeYGEVkGRgYYyqmJ8VxRE+3RYsW+4YjRIwicgyb2ftbgC8DhkBXRM4B/xP4NWPM7dP2c0RO9xCjSgUK50qHjhOwUuFZ4InY6hVnVLu2M8acDfrYTdXamKrXRIqskjpj+6k65jT4ffTU72nSUWw2SygZ6notueRm4xX9VD+dUd3PbEkYM5E1x8HPRNJq9GjYmbyWMc+5Xq/7o/pU/b0PEOboeKX/GPgD4GXAR4wxOwAicgJ4OvCzIvIOY8x/qdvJIbnsRxQxex8UD6EmC61OnXXLTwSeCjyOoqBVvr19eBd6I5bc1LQkzegkxawLnYIsy5J8auD2wCWX2OhNhvgA0QzeodNEn6P/1mSliS8WbD5BiMZtvB20L7gTlzJZaDW6SpWG8jULibEftttjexW6mOrnbIzOtlgHnbxj5FToMQnDQWc60fl1oellEPQNzy3EPqvT5hiMQkfgLiEiVwFvBk5jb4objDHXR/o9Dfhl7A1yzhjzFa79k8A6kAHjSEXCrzLGhDcZrtLg24C3uWJctWiJcRZoO1AdKaLWo7bxD+rjXNt18BVPfhdbLPIZTrM2tG7pbJxYEnRhI0U83UjZv8pZv0dJh+Fyl83lRbawUwX99MF8xox39nii1Da9sKqhflh1WxW5bgTbYYAtdYG2IxfLH2ABWLTfVURcZ2cMf1d4n495W6KXFNNiNjpQsimG0PPVffKOoZcWsy6jQXdyrFXhOYPg259nLA6zys66jzAC46Rp3bydaR3GwEuMMR8SkRXgFhG52RjzMd9BRFaBXwWuNcZ8SkQuD/bxdGPMOSLwpCgijwXuNMYMHcl+CfBmY8xajDhDtMQYoipMp8oxEYuvCz2qYwqJsA88zTZ/5ZP/B2/hW3klr+Av+D8Ydp0kknZJkowlNllhnRXWAVhhnUU282BjD+8z9cklPDGuJytsHl9i/fgK6xf6DPpubvVarwgb0ll/dJnP2PmHKvOG+uT9DfaFvoUlwy3XHpMU/YGXJtu9lK1NEiGahuv0B65Wi92RlxQTxqVrmdY4WewlcBmMAolxNOiwE5MY60hRX0vdd5rEeAAEaUTI0qZUMapda4y5B7jHLa+LyG3Ao4GPqW7PA95ujPmU63ffzIO20uHVIvI44AZstcHfAr6mycYtMVahLo4xJr34b7+spStPiH2sHfFK27zGKj/Gf+B2Hsd5TrJxcSXfZac3gsRKLx13sy2yyQobLLGZ28GAUn7GJZbYdGTkczx2GdI9PmTdSUzr6Qo7vSVIZfIhDAkwbPPk6Ul1AwoyBHgIS4b+s612oi+YJ8ZFyiq2G0iYFahKctS7jEqM2y5wu/A+J0kxg8V7oaHwPsfsi5oUfRq4obMtDgddGCxUvzCLnZTtwPqc9NirohAOEFnS2Mh4SkQ+qH7fYIy5IdZRRM4CTwLeH6x6ArAgIn8CrADXG2Pe7NYZ4D0iYrCOlOi+gR1XVfDrgF8xxvyKiPxN05NoibFFixa1MMhEsbYanIvY/SYgIn2sVPdiY8xDweoUeDLwDOyb8y9F5K+MMZ8AnmqMucup1zeLyN8bY94bOcS2iDwX+Dbg2a5tqm1RD6CFx7SrUSet+G9tfyRo72NjFN26u7nCKcornD9/im2XRYc0o9ffpHN8OFk90Mkues4ujMhIVG5vnS+wGHTSde2nMtbTjG1WKGXy9l5h7SzS63IJUS9vAw9gJUXct1alQ4nRS4taYlzECgYal5W93tNMGVGJS1X2c4lmbdcsj1MsX6+sdL2Kw+l6gEmhSmddhoOO/bvF5oLHHC/62upDxZxHeruDtDEiuaNvHnDOj7cBbzHGvD3S5U7gvKspfVFE3gt8KfAJY8xdYNVrEXkH8BQgRowvAL4H+I/GmH8Ukc8FfrPpGFti3A2qjOtVnrtY6AawdmGVrd4S62srNrein2ubLjBKM0b9LpvJUu58WWKTEV26jBg7IvTQKmE5y/cWIbIkIesnrI8TdsbLttGTYmyeduiBxn9vYklRE+M6hTqtHS7blMnQH+QyyvZH1We8MBmCVKdST6jSwzyA22ffhskpfp4M6zzS3lTh7YpgEwSPvBodJoyIOV7GFE4vjV7QPzRrxF4K+/jkGiTPEXqpEBEBXg/cZox5VUW3dwKvFpEUWzTzXwC/5IK1jznb5DLwLOCVwf5vAH4f+ANjzL/Pz8GYfwR+tuk4ay+viFyJDZb8P4ErsHf7R7CBkr/vY4Qe8ZhmpwslMHdzZ+OU4QCXdVtK63cGHTY3Flk6vsmmc050lKNgSZU28JKPJ0T/gHecJOkJ0q8fkzDuWnLc8KE9vYX4PN7Qvpg7aLyk6D/nXbuXFkNi1BdrG1zy3QI+bGfT/XZe6qp54SFKkqPNlpMnhUiGLLFVCssJpUV7HctlCorvNP/2dkWwZSnsC41qYkS1x5wtfrlKWjwEmFGVnoZrgOcDt4rIh13by4HHABhjXmuMuU1E3gX8HdbN/TpjzEdE5POAd1huJQV+yxjzrmD/rwe+GvghERkB7wHeZYz521kGWUmMIvJGrLfof2CZ9j7sn+0JwLXAj4rISyv0e7+Pa4HrsUmLXmeM+ZmKft8A/C7wZcaYD8b6HFrEPI8hMfqCUl5NzYkxoetSXY3SzHo2AcYppBnZOGUzW8oDkPWDO6LDopMGy2p1GYkLWs60Cuj/9UZ5CMtOulDtMIhN8Ss5WPwnbI9FRYTxJ97xEqrc28XPGDFWxTcCOO9zko7zkgRdhqWaOqG06JGRltKJFe1+i9TGLIJToWXy2sSkPz92f09UTV2MvEAnlg8A8yJGY8z7iFZim+j388DPB213YFXquu3ej3Xm/ISInMRKlS8RkS8BPoQlyd+Zdvw6ifEXjTEfibR/BHi7iHRwLB+DiCTAa4BnYm0GHxCRG3W8kuu3AvwAk56pFi1aHALM28a4XzDGnAfe6j6IyJOxQt1UVBJjBSnq9SPg9pouTwFudyyPiPw28BzK8UoAP4mVSH+kyYAPDWKSooe3JemAar+NkiySNKPTG7GYZGT9LUaDQkWz68dk44TNxKrSWpIZqrARH77j8wfGoMN+bJWZDkmyRJJaSXPnkm1W2n6IWt4O1qcUDhffZ4lcdc7bF4rudYgGfrvypi5bjq2as8lSHsZYbzI0AAAgAElEQVQUTyPmi9aOA2nRHibJ1ehRnpB2odpBpG2dWpWOxSlOc7ocsLRoVemj5Y5wQeL/BjvPLB+8tjvWYerZisi/wpLX57j+YvdvLpuy6aOBT6vfd2KNqHrf/xy4yhjzP0WkkhhF5IXACwE4USmk7h+qVCUoq56eEM+p7c4V7dk4ye1fJJAtu9kVy3Y2hS+u5G1a692VEjF6g/iIjvvVKT3sM6k/VedTea4rlGMPvVe5To323mhPjlCQpCdIv5+lsu0t9D7HxuR+H3PEmKSZn5tSKmK1xFbJxqivk33ZdMgmCDPN1ejc5FE1JzymQhP5HarSfpvYlEyN8Jyr+s0J1vlSUcP88OIm4K+AW2kwHSdEk9fALwNfD9xqjDGzHqAKLgvGq4Bvn9bXBXHeACBnr57bGCYQhtmEiIXgVHlwtcS4RnHj9smJcmdjCU6TZ4vWdaI3kyWGSYQgu0WtaU2Mi2wxCogxYVwK/m5ElJqMoEjEsBpcnw2BjZPASSZnuNTFNS2ob98mxWIsAURVEt9wnO7v0fH5FZOhm90yLEmMelqlRuh00e0jLS0O1OwcLTGGzpTQIVPlnNHEWGVPPUDYWe5HTpXuGWN+aLcbNyHGT2OzVMxKSHcBV6nfV7o2jxXgi4E/cV6mM8CNIvK1h8IBE1PRqsJ09DaaGPUsEf+Q9CnyMV4pbF25BN3zziFgH9yMhA4j1lmBRJVTHSeMhh2yNGGYdEvzd4d0Jx74cGaMnso2oRqlBnoySUT+ge+5Nj0nOicFFZdYFU4SSnxV7eHvWGIIPz6fhQjdb5s0Lab4dZ3EmKpwpo4qVaAxpDuhXttTsvGLw3DqX5UaracxxkKN/PmFpO77TJvlo/exLzh6qjTwmyLyXVjncW5jcskkpqLJ2f4/wE0i8qfBAapikDw+ADzeBVbehQ37eZ7a/gJFalbc9J8fPnBSrLsZq9S5cPtQlfbT53ysoD/re+HCuVVOPvq8kxrdlDWlCmckZIkLFxnb2sXDQYc0zch6lhg7iVek6+see4L10mOWJfn84e10bEN2PCFqIvLnor2qTYPdY+1VmBaG48eVz4GmIMZVYNWw0N9icbmYEunTiOkXhrYxQiEh2hdJoTLq7EV5oghvV4RJUtRjnEZuMWm3ygMfux77iDmH6+wXRliv9o9ihV7c9+c12bjJZf6PFGGpjQ0Nbp7ii4B3Y8N13mCM+aiIvBL4oDHmxqb7atGixcHiCBLjS4DHVWXhmYYmxHiFMeaLd7NzY8xNWCOobntFRd+n7eYY+4oqiVG/4cOsM1DYG73E6FXpTwJnepw/dZKl7mZpapr/7jDM8/9Bh/E4IRunZOOUsbc9ph1GvRGdpMOW8j6HyRBK+RtJbDIFF8eY9bfYTjPod2FVJqWiqmuhZ2349t42pBnH0iz3egO5dAo2uN1+Jy7APbHxmz7QvS5eMT8OhKVPl/pbpZlCxUS+6Ya6YkJl4kwPhQ03H2eoGofqsb8m4T3hl/V9U5VerIkqvY+S4xGVGG+nmC0wM5pc3ptE5FnGmPfs9iAPC0w8/ExevapZIroOtHK+cC9wJ1w4dZrOY0f5zednZnhvaAxerQarMwwHXbq9Tp7HEQpijIWl5POtu/bp63RHZFn85k+SjDBFl7dnehOAtt/pzDWhvdMevzinkVPrCyNAQUbFFLzCA19kUCxfFz8+Xd1vySVhC22J4bblY/rjdfIXUpGItlsmOo0wPCcM0YH4i1W3V6nSB2zeM0h+jY4QLgIfFpE/pmwCnE+4DvC9wA+LyBDrcmwarvPwRZV9EaqzWuv1OoznTmBVuD/9bEZn7IO40l3Ppb6RSozqUZKyANKMHUeUnd4wJ8xOb0SWFE4Hj9j84ISsRID++JoAtWdX/w6T6eq4ytj8Y01sBZVZIsqT7Lq9brmj+GmRWyyWSDT0JOsXgD++P5f8+gWSs3ZijdwYNlkqJMY8Ea3EvcwxOyhquUmoUUxCjkmj4fI+4IhKjL/nPrvC1EtsjAnTnviJ4C38ze2Jzt/UWkL07bq//+1LqX4SYIELG2cAWD+1wlJ/M88f6Ikuc2r0BMYJPqhZr0/GST6dEJgoIK8J0IcMLbHFCus50fXz/D/rrLBRSprr+/VZL0JhsiGLG9ssDN15xmYqJuR33nYXRr1jbHa9fFcQoz1CPz+ybtekCeRSpY63C0NMYuE4oVOqoPmCjDc3Fq03uolqr6FVZ90P4k/ebsJz9iGU54gS40eMMbfoBheT3QhTiVFEXqntgi7+8DexVbgeuYipQlBWp0OJwntR/e8BVp32ZOkIdmdjmY3V5YnM0xOkmMZYhzxcxafu9wrpYm5728olvhV30BVFgKusseqCMFd5kFXWeJRrW8ksMV52fhsuYD8PFWPnIlZx8cQYSj+BtLTQhYXuDsvLG9DfsMl2/LW6DC6ePMZ6d4U15y5fxy7rb99eRZbFrOhJSVET44iOIsVFti7aIHQrLabxRBF6OSYZhhJjuJ1HU4I7ILX6CMYx/rqI/Bs/g8/lZnwxNnxnKpoUcrhKRF7mdt4F3gH8wy4H26JFiyMGPyWwyecQ4RuBN4vIE10847/DJpRohCZn8h3AWxw5Ph24yRjzy7sa6sMNdeqQXr/qlvtqGawEsqaWQ0/mao/tfpdtN5PDT3XbGZff3t77m6RjUjf/GmAp2cxtfl7lBUpqsZYMT3GeVdY4yXlOOQ/RSc6zevECvfuwmcV89Y0L2GxjXmq86Nq9xOjPp0qV9mpn1y0vqw/AcftZPrHD8vELnDl5AQBzOTx4oscaq5znVC5Jrjk518u9QL5s7YVlx46XFHX7JouM6FrJ8+IimxtW8rRB3VKtQvvvJqq0b4/ZKg/JTJcQR1GVNsbcISLXYe2MnwKeZYyZTE5agUpidPOYPa4Hfg34c+C9IvLPjTEf2uWYH36IXcWUIlgaJmtEw6QtMuq1FBjbJ27Hlf/0KrSeE9ztDUuVBYGSzVAX1VrNaWStRICeEC/nPk49YHVjuRtLhp/BEqEnxvOUVWlPjBsUarTy+Rj3sIs/L/+c9SiTo1alT2AJ8gR25iEgl8OJywecOH0vpy+/j3Pdk244pzjPSc5xkkV34avmjY9JlIOlU1K9t4ZLDF0C2h1HjHlqsapga39OsTKoVU46iJPjIYT1Sh+NudIicitFQDfYuycB3i8iGGO+pMl+6iTGXwx+Pwh8oWs3wFc2H+7DEDqUIryp/fS5cLqad9RoKdGjqjZz6S+0gE/CiooR9KS42LUVBD0xehuith0COSGedFQCcJr7uJzPcPri/fTuBu52h7zPLYcSYyAtbjpiXL9og8fGFLOnYTLHjk8hkeJSRyzDYhfkuFuxTJkUHTFyufvcB8tX7LB8+f0ArFyxrkKHyl5o7/FO3FHzhLN0rESZWQIcDTpsbiw5CXGh7DzTs1zCp2aaDTEkwGnximnQfgjCdQ6ZmlyHxg6WOlSerTHm6fM4wMMGVapPuM6rzGcoO11SylMDoZhPlAcrq334/YYxcT2Z+Kt5SbGrwmqACVL0JPgo1qxU6KRDgNN8hssvPMDCPVgi9MSoJcb7sIQIcB7MBXjgAjyUTdYIHFNd8QWKgqkpLrfORbjsIlzm9r+yDEvH7XFKCcK16u4/wIlsQHrF3WiNzztftIfZIyO17a72NmDrbw+6NkFGLJ6wSpUOiUyjTmVugkPCR0dIlT5vjNmo6yAi/Wl96lTpb8UWq4kmj3AFrT/bZeR9eCH0MofwEoSereDhpUUtMfr5xtpj7Y+zSlFS9ZTqnwZ9/b7H4DN865kk1vvsYw6LWENfm3qVtVyVPsl5VnmwJDGevPgAC146jEmMXpV2BLX5ADxw0RLhAxTSYUiM+tLoqi/+EvsMjItYcvW8e8IR5eUDEO/lBmuzjBFVCpcl2wyvOJ+rfVsu/MeH/WhJMp/dMuhYQgTY6BUvrjryqgq/qZIUY98xKbFO7T5AgpynjVFErgLeDJzGap43GGOuj/R7Gjaz1wK28uBXuPZpVQHe6UomvBO4xRXUwpVFeDrwTcCvYysGVKLucp/ERo7fAtwC3I99NB8HfAU2Au+ldTtv0aLF0cecnS9j4CXGmA+57P23iMjNOrO/SzL7q8C1xphPuVKpjaoCGGOeISJfA3w3cI2IPMod8+PYWlXfZoy5d9og61Tp60Xk1Vhb4jXAl2CFgduA5xtjPjXL1Tiy0DbEUC0aRNb77Dk6dVcP+xq5V/3G9fti7KvmFJOqtN/nhEpXjq9P0iJgouuSkIF1viyqyDwvMS4FtkeAnq9pdUF9oKy6XgCjbIm+UOo6xaRUXfHFS40aXrLUkqPP5rhNoWLnavgDcCJTBYETigBxbX5wHu3V5Q02j1tb6jor+fnr6Yx+1s+YxM43zzNyq09EIi19x9piqnddhcNDYkOchnlOCTTG3APc45bXReQ2bFJrndn/ecDbPccYY7xlu1FVgFiOhllR+ycxxmTAze7zyECTm1Tf9KHzxT+sp7Aqst7GP8xnXdv/DjwNjp21bLNzzsWqeFvkNCO9g53KV3hgU/Xtp8TpOcSLeQjzJksXna4+oGy386rrRbfOeZq3XLtXk8NagCEp+vawQGrICTHNct31TTfgpH8uL1BcYx3ec5ldt3Aclo57U8IWS2zlWbwnE9PaZBz5iyZa9CsYVJO20Fkzzcnil6ep0weEGSXGUyKiUwfe4BJNT0BEzgJPYrLe0xOABZeKcAW43hjzZhpUBZgXDtmf4Aihyl6k7YZnlIHwXK8gTJerqP899/M1yzeRkPG/eCwfX/18AC78/ZmCGGuOlUZmvoTlQIvadmXiTFwS164nQD1TRds2fSyi+x6742vHCsQLGVSZ6MK+miQ1wS6673EG226cCz5GUhM5lEh9aehm+HQ3XRXoYSnDjp6/nY2TsnMrLFPgB6aX62Y86eVY6dew76VgH0N8ZiDGc8aYq6d1EpE+8DbgxcaYh4LVKfBk4BnY2+AvReSvZhjuJaMlxlnQ5Eb00syq4fgpq9INBx0GPUeMZ4Cvsl1/ePkXeMUDP8cfnHgq/5Vv5jPd0wBc6J0pHy9y3GMV0wGrUM50U+SokZj6t48PHFiiTNUy2KdBS5vb7sfCmPiUQxVUnoxtiY+kW9R88dIzMCE5TvUme0x7Wqqu5T5fz3lj3lUCRWQBS4pvMca8PdLlTqx3+SJwUUTeiy2beif1VQHmhiZTAls0eSBCqTEdu9komZXsvKTRh+NPvJfjT7yXb+eNvPYkPPO+93GS886rPJzc7yU8XDqTjpeaEkUSE8fKgt9T1Pg6zOOtG0qmpWHHrs0YEv9x5wpl6dn/nsC8yKxq2zov9yFUoT3mOSXQJaB5PXBbTRWAdwJPFZFURJaw6vJtqKoArnzzdUA04bWI/KKIfNEuThdo8Kdw86O/gckyhK/c7UFbtGhxtDBHr/Q1wPOBW11YDcDLcTXqjTGvNcbcJiLvAv4OW+HvdSoZxERVgIrj3AbcICIp8Ebgra6cSiM0eUe9E2vyvoXSJK9HEKZJDzEPwjjN04WNAzvWhdutqvxbX/QtvOyuX+aPLv9yznNy8o1bM4NiZ5zkGbwvCbNIRodBJdzLOOOYxznEuGbdtG3rcBiubQXmWT7VxT1PTVtojPl5bM2WsL2Rx9kY8zrgdSLy+cALgL8TkT8Hft0Y88fTtm/yJ7zSGHNtg36PHISqTxjGswFsCOtrNpnBjidGn2bsT2zX//RF38/Hr3gCAJ/g8zk/dPPe6oz/ETUvy5L5EkaoZe7yoZ3Hs67qDxao8+CmkLk273KqQp5VPOZcqRt8aDapskeG46sLAJ92zAPEvG2M+wUX9/hE9zkH/C3wQyLy3caY6+q2bUKMfyEi/8wYc+ulD/WIoyq8wt/QftbEBnAv7OSxJGrdOWwZcODe9PP4jad9L/0r72c8Thice5Rd4YkxjGssEWRSqpuSJYUVbVfn5M/FxwnG+swROo5xQbXp5VStX9Dj9GPsug8UiShSyFJrOs9KlsX4tUnSjJ3wGsxyviFJhvsJvdiHlPzqcMTmSgMgIr+EnTf9R8BPGWP+2q36WRH5+LTtK89WZalIgReIyB1YVdqXNmiUpeJhCX+zh2mmUsoZc/TUv3OqzQd6u+WNx32WDeOpipOLxsoVqjqQu1X8Mmhi8HWly7VNRnQwbt8SPsialBNmJwyqeWAhWE6xXuiQQ/xUwcUuLPh3zDJF6dRIgDfLsNktsuUUn27puthj2GmV2z4xR09KJo9GJ1MV8F0l1VbFPR5yHKG50h5/B/yYnxIY4CnTNq671eeSpeLIIbQhTVOTUsoP0RgboH2OouhVj3JMm5++fjtFTsYrKc+V1tX3UrWfyF8sGydk3dRljemy6OaXjElcun4b0u1nuajwbjaXbUjR8vJOOSeil8SWsfGBLj3YoiOipWERzO3nO+vLFSa+08PWxOgJ0ZOgzzp2mf8kNqEEmhj9SpezEYrl7ctsXkWwM1+2XO2WqsJiaZqBy3dJrzdpFqk6Cd0nos7nLxf9kulV9K863jTskwR6FPMxAt9qjHmjbhCRPzTGPKOJE6aSGI0x/+R29pvGmOcHB/hNrGepRYsWD3McJRujiPSwM0tPuXnS3tFzGXbmTCM0UY5KsUDOoPnkpgc4kqiLN/PfWiLQ8BLgBpOS3gZlCXCMlSqrVGidx1EfN+g/HidkmbUzjpyECLhZHzYd2SJbeWbrIqv3FmtdO6G7c+IBFmJTAn0wtZN4F9y1WRnDtlJSQiF7gcnsOrqvtjEuYu/kRcBXXjuBlRZPnqTIywg2F+NpbIoTn5sR9/sErB3v51m9N1z2bp+QVks9vhZOpzdiwWU83+4nxchimkBs2UttoR2gysYI0b/hpcSL7jWsV/rIlE/9bmxtlysAnUz7IeDVTXdSZ2N8GTa+aFFEHqJg3hEQnfv4iEKVDckT470UN7m2haGW9cNQV1UwLMyeGkjHeYmDzNkbh0knMsNFf8aRdtd2PON0dqE+jkKd61IKixfcR+VjvIxCxZ4IzqbsXAFlR8SS4mWOuy7rw8JlWMI7TkGAJ7DEeIVru8K1XwEPnV7gPCdZwzqxitIGi6VSq/ZUXKGwZMhS37Zv+vGmC+UkwVWB3/rvVHWtmhKj3l/orZ52rD3GUVKlXQqz60Xk+40xv7Lb/dSp0j8N/LSI/LQx5mW7PcDDBnUSQQhPjmuqv3cYhPvRUmGsXVcVdIQYIhu7gvBA1Ytd39i6Sp5vH9FleKLD6e799LS319v3jlPY9bDL8hCcvGA/PoP31gC2ssnkEh7+VLxNcgFYTBwRhjVf/PGcNAgUEqMnRUeY91/e5zwn+QynOe/SffsaMFtRiTFjiS1rd3TnmqQZm+nY5mdMu5AGySWInFCdTVKfsP5dd+8cUhwVYhSRrzTG/BFwl4h8fbi+YgriBJr8aV7uDvBUrJf6z4wxuy5k/bDAtJvaS4Br6ndKkbxWlzxYDT4E/frb6BovYSEssHGSo0HBiFnXDk7XUS6TYKfkmYaiwP3G8gqXP/4+Tp10NV909u7LKTJpBynKlhwxLnk1PDbF0F8zVVd6WjGsvLSBJ0alQj90eoFzSVHzZY1VznFSFciyRV83WSqlzUpcxqHit2W0UTIiOW6ncA4HI7Y3NH0TTzOnEZP4poXzhKja5wGS6FGyMWJzxf4R8OzIOgPMjRhfg80Y+Fb3+3tE5JnGmO9rcoCHPWI3vm/XJQxwv3VYjidDn8Hbpyk7A5wyHOtv0umVJxt5tTkkyJAcgVwS8oToyWHorI++eD3Ag6y6clgnOccpVk88CMCpE+c5edUDNl9jVZVAneXGF8OK1ZXWxNBVyz4ERxfDOs4kOQLbzo543hGgJ8EHWVUZJsu1puPS4mY+l7rrDKqbLNnXxfKQzd5Sbh6w5oCFeCB2bHnaXOuYFHkIbYseRymO0Rjz4+77BZeynyZn+5XAF/gSByLyG0DV/MQWLVo8zDDPKYH7BRH5KeDnjDFr7vejsJnDf6zJ9k2I8XbsBO9/cr+vcm2PbHipIHzba/tgCO+V7qu+WmK80rWfGtBfXafTG5EkmZ3yh7cl2i4liXGcgFOzs3Ga9wFIusXgdICzrZDXz6UrW6d5nc9weqKi4OryGqvLD7JyVdH+qAcGSFiYCiYS25ZUaX/OoSqtJUZVJdActzWkfRVsKCRDLSFCYQrYco4WKALZvRnBO6W6DBm7311GeY2YDiM2WSQlgwQy7ZQZp7aMqkYoGYZVBWNJb2MOmGlOmQPGEVOlPb7aGPNy/8MY86AreTA3YlwBbhORv8bq6E8BPigiN7oDfu3sYz7iqEtNpYnRk+MG2EsnZZXSJ669EjgLvSttKajV42t5Kv4M620G6yBJ04xsnJbtjUGd6VTVm7aOlqxkY8tI83rKnli6nMzrUOsyCItBXWpfmnXlxDr9E75tS9WyHrE03CQZ79AdUuR7VDApDL2anx5j2O3kdZ59cHZRxGrFjdO+TbwzxffdymtCdydsqRqht97LPxkJHfz1HZVSso27xQtpu9exAeDhjBiYvB/C8ghNQnGaEmKTJBZ7gKOiSiskItI1xgwBRGSRStfkJJqc7St2O7KHLcLwDX3De0+ydrKsQZ4+X4ff+PIHZ6F/9n4uX7YGPGv/cqSo7GNZ4glvTDZOSslqE1c1sNsb5fWmO8kQnbm7eOtb58tktu9x3t87JxLGLLFFl6Gjoc18jP63J1Sw0linO6LTHZIux/MfavLyDiBd6hRQe7dH8nF0ftlelyJJhP8uE+A433vYnqpr4s/VE6NvzyXPXpeF3ojttFcfgxirCx7WfKl74kKnTBO74yx9d4mjFK6j8BbgD0XEz355AfAbTTeeSozGmD8Vkc8BHm+M+QPHvKkxZn3atg8rxLyFsVgUHWajy6euUajOfuqfU6EXrnyIy5fvYxXr8OgycqRYfsF5lRooOWXSNMvV7q4jAbDFsCYS0ubDT0qxffp4+jhga8rYcZX33S35tkelPv7YIfn6Y/nj+rS5oZfcE6X/DvuH0MfV5WOL4xZeaB2s1IGS/Swvq0on336UdBj1Omz3tm2Mo1Z9Yw6XqqJaVS/TYmCTKnYobeqA8n3CUSRGY8zPisjfkufL5yeNMe9uuv3Uyysi3wW8EOsbfCxW8Xstth7DIwuhB3JA2Zao4xV1+I33SvvSBt6W6OyKJ0+ey1VSD50AIrcNOtU5dZnBtWSoSUrb0srDL2IWx27/XkLburjI5sYSO4MODBaKMWvveuxFEHuYU6DnEjP4uMuqUgzeHOCLUuljVamgOvDdHedYf5Ol/iaLy1sl4vR2xNArrYPbtVvBvyB8QJNdb2fIHOuN2OktxNVpmC4x6vMKr2uVnfGQaLBH0MYI8DfYWCvjlhujSWmD78Nm3X0IwBjzDxTzEFq0aPEwxw7H3MTS6Z/DAhH5JuCvgW8Evgl4v4h8Y9Ptm7yPhsaYkS3VAC5VuGk4uGuB67F+yNcZY34m0uebgJ9w+/xbY8zzmg19j9DEZqOlxfBtH7MxrlJ4o7X3+UronXmAFdbpMiyptZvavpY5FddJV0ma0emOlMNjmNsBtaSj7WheXfXISNgiydXm4aDLztpyMWPHS4x+2adTC+1pVVNbUnHLfu7xQnl9U4QSqXZs9YFVe1/urC6zsdphPE5IjttBdRnmMYxa4tEe6hBFobDCDpkwtiVq04ydnilmxPi/fyjZ1qnSXhWO2R7rVOkQ+6xOz0uVFpGrgDdj5y8ZbHnV64M+T8NWDvhH1/R2X0pFRD6JnYGaAeOaioQ/CnyZr0ktIp8F/AHwu03G2eTS/qmI+DnTzwT+HfDfp23kkk28BngmtrrXB0TkRmPMx1SfxwMvA65x7vTDI4mGBBm7QXXIjt8mVKehINFVysR4Bpb6W6RkZWO/8xhvssRo2CnlXez0RnZ+b40jJLQrjrG5GBN1EvZ3lu97e9ApkuzqlGnnKMhSE6bOO1kX0BxTC+tCVWLTIrWH379wvKmiFAqzwKjXZdR3gexJl0U2c9XYvxgSF6zjUbZfqlyNuUnCWju7vWHhnYZiJownMH2v1NWnDlVpPS87dm1CVXyf1es52xjH2HjCD4nICnCLiNysecHhz4wxVakPn26MOVexzuOYJ0WH88xQ/K/JJX4p8J3ArdjMFTcBr2uw3VOA240xdwCIyG8DzwH0Bfgu4DXGmAcBghPZf0y7GtMcMKHEuBr0815oP8NldUCSZHkS1cLzWhDieJzk4Ted3oilxGdS3MpDanx2Re8E0QQ4zXNbnEtSEOMGxXTGc5EPbv0A7Evfp43wJ6vTR8TS0vpv36b6aAKMXUttu9UefkeeO/0OIzdvfLTccSreyHnhx/n5Z6S5VFjEd8ZLIXgvvc/Es+2JUUuM+t6pkhirZkb5F2fsZTFNc9kHkjTMz8ZojLkHuMctr4vIbdh0YCExXireJSLvppix9800qBXj0cQrvSMivwf8njHm/hkG9mjg0+r3ndgyiBpPAHBFahLgJ4wx7wp3JCIvxDqA4MRjZhjCHiHmZYRqidGvP+U+rv1YaoO3NxOrLvsHejjo5pJcpzek49JiLSWbeTxhXzlrVlifcLyAd+BY4g09xDlZlpwfxKUZnXjXE+PYYOcEPoA1P/sgBZ1XJ6ZneyL085B1fp3LYOASjw2WCvPDamQ3/hqHzppxmhcJ84l67b9u6fy1Z1qToQ8D8gHgHj4TT7fXYbvvDjroxT3NWlUOySsN+kBBorFkIvsQjjMdM00JPCUiH1S/bzDGRLNxichZ4EnA+yOrv9x5le8GflhVAzTAe0TEAL9WtW9jzI+IyDdg/SN+HO9oehKVZ+vqv/448CKcCCoiGfArcyydmgKPB56GVTDf6+rLrOlO7uRvAJCzVzeyb84doXQ4LeRC2xh9Py/tOIEjSa0quzFYKTzCYEkK62lN04ylxAVVu55n7tQAACAASURBVBkpq6yxwjqLTpX2HmltU/SIqYmeEoZ0ikqDY6l+yDVZ5ue8TkGMnhxx31sUkqOGlhg1MV6GnUewBe6c4CSML7PSq8+Arsejr2tpnJLPGR/1u4wSKy9uBS+FqNRMOZSokLCtZNllxKiUv7FbXLcwd6O+VjFzQYwYw6dR963CPhDmjKr0uRq7Xw4R6QNvA15sjHkoWP0h4HOMMRtutsrvYXkC4KnGmLuc2e1mEfl7Y8x7o+M25m3uGDOjTuf+QSzbfpkx5oQx5gRW4rtGRH6wwb7vwk4f9LjStWncCdxojNk2xvwj8AmKC9CiRYtDAIME9XOqP00gIgtYwnpLLA2YMeYhY8yGW74JWBCRU+73Xe77PuAdBPVbRGRdRB6KfNZdXtlGqJOPnw88Uxs5jTF3iMi3Au8BfmnKvj8APF5EPhdLiNcBocf594DnAm90J/4E4I6mgz8waO9iTHLUUqNGUMDJz2DZ3FhiZ2OpPBe3Z0jSjMXlwpa4yhonOT8xPU+n0ArzLo5dsLKvBwO+SFTXqe/d4pyqzqcRvHTod7DFpMSoK8Todf6CaUlym3waZR0idl8/VXI06DBc7pC4WM0ijrMzYYsNEaqOeUC4U6eBwhEzpmwj1NcxVKPDflCWGMO42Kpz3vcA7/kc0GmirwduM8a8qqLPGeAzxhgjIk/BCnDnRWQZ61RZd8vPAkraqzFmZXKPs6PubBdinh9jzP2O8WthjBmLyIuAd2Pth28wxnxURF4JfNAYc6Nb9ywR+RjW/f4jxpjz1Xs9QGiVJrxqMdVO3+QT3mtrDfDT94aDLjvh8XpDlvqFTRHIVehV1pwHupi2B8XD7N/cIzV7xPqtiyJRmyyyubFYqO91HmbtIc4faq/6blEuYxVuWBXPM/UWAqQ6gDxEidiLEKROb0SaFDNu7HfKkE7JzlgHO9s8y0OqRk6V7vU3GYwTGC+UCU2r0tqjPnEN3be2MVa9ZAm22WfM0St9DVboulVEPuzaXo5NVIMx5rXY2MPvFRH/hr3OkeRp4B0udDAFfivmk/AQkadiZ+x5wWvFaaZTUUeMo12uy+HE4JuCtleoZQP8kPscflSFT3iEiUxD+IfFzQZZ7VqiO3fyJOdBJUaFhd6IxW6RKgHIIxt9vKInRO9Z9VPZRhR5FzdZzLPT+IQRWyyxcXHFZqpecwQVxinqc9Zxg94utiZYu+C26giFc0UXONDQdkYoV3zxH4rvKmKMEWQeJmP3nfVG1qHVg3FSnoZYzmFeXMfq+EbrjukwZClxBNuzJLnTT620HxJYVdRCjOy01LjP4TjTMM9wHWPM+5iiBhhjXk2kPouLcPnSJscRkR8HrgY+H3gjdoLTf6FwxtSi7k/wpRU6uRBPqvXIwbQ4M/8d81KukTtXruLT/DNu5eM8gU+c/HzWetb9Oh4nLPW31DQ/G7Ad8yrbZT/XuMhSA4Vk6NN0eWLcuLjCxtoKrPWaxSWG3nZ/ThtLVHuZPSnqQqrblAlR9/cOmIAYY57aGCIe9Z1BJ3fl0AN/ybzzpZjP7SXJciINj2KOtw3o8ZL6YneTrJ+wMU6sSh0jRj3mUGIM+8bUab99lbayDzAI2c6RmxL4dViP94cAjDF3u7jJRqi8zMaYI3clWrRoMX+YHWEYZoc//Bg59dsn2F6etoHGIRPaDxGqQnH8upg6FzowtMSoJbIrrSax8th1fiZ7KT+a/AfOc4rxsnMOZF07Bc3JKNpp4h0J2jbmFcIwR6HPt5h/Llhxb7C2YlVoHci9EYwxPG8vLYb21I0F7NT5MPxmi3Lgt1apwwBvLzUuUdgqpSwtVsUChhK6TuwxWGAHGwA0HidkPRcb2rVmCJuLsTBJdAMLUSg1hkk5vL0xD/puqkp7W2IY4lPn/DrAeEZjpDT76ojgd0Tk14BVlwjnO4Bfb7pxS4wa2ttchyo7l344w3nFnoRUWdXfP/X1/Osv/W+ss8Iaq2wNLaFl48RWrOsuBrNYEtYrA7ntlMItNy8GCmLcurhoVecNxxh6/rNWpf2nzgHTD9pSYENszCFg1eF1iphEPSMGyuq034knxZXyujpSDBFRpS0sOWrKy8aJDZpPhpCnqbWqMkDordbXWpPjmISlJCHrJ1wYqKmCHlXZl0Ini79XQhvjYXk6DUeOGI0xv+CmMD+EtTO+whhzc9PtD8ulPzyoIsXQzlNlJPekuEGZGPUMkk+69nvhPU99ji015kN5/P5722z2l9jqL7HYtZayNVYr50LrnIZbF630Nhx0rUNnY2GSBAfqW7dVpeP3D3PM/pXq/QjlgO1pUwW99KgdL0yWla2yy4V23OjfZYGdcWo9yJAHaFu74zC3MfrZLhnpxDWGgiC1hJmR0umOikS2/prExujXhcknwtAvfy0PiSXfGGG8fTSIUUReg/VW/7kjwsZkqNES4yyIXa06iXFiuhpFkgaw0uMnscR4lmIO9Smgv8B2f4EL/cu4oHMbplmeuXui7stgofxg6XFsVLSH5F2lzsWcAr49pSBdT5IDofA2+4tQN00wpeSsnOZ99t+DSP8YeoLPbF/ykyunTOLCcmKe6bLUWI4bHdFxCSbUWGJOuHCcYR9N7r49PJ8mGs3cIexkR4YqPgH8goh8NvA7wFuNMTPlYoSWGHePaVKjJiJ/I3tVVNv1/p5Czb7StZ/BkaPbxqe5Sq0HOI95DO19VVJUKInEpJOqB9mfqx9/2B5ehzBJax607r3XsRmdOrBdfdcRXZ3EWNXPl5egWyLHZNmSnk80keT/lwkyVLG9h9onmChJ/J7AqqRvVHtIhuHflWD9fsMAR0SVdinMrndVB64D3uCqDrwVS5KfaLKfxml4WrRo8QjFjsAgbfY5JDDG/JMx5meNMU/Czq77v4Hbmm5/eM7k4QAt4cSkMV8V0EM7OrSH2EsbucSotompulVSRkyijH1XtfmxaGipSNv0tMMgJklOSG2R42hpyi9P8/5rVTo23iisWr2NnZY5Gro8jd0xvqqi/4bQO11coKKIVodOMiyPPRwr6jscaygNVvVH/d5vyfEgJNVLgEuo/dVYqfEZwJ8AP9F0+5YYd4sq+49eHyYi9cSoE9j67XR6Le2ICR+sUD2O9bkUxB7W0Dmg22Oe1FpirBhjSIxVqnSV2SAMfZl2zNSS4yjN8nyXo7RLmmQldRooJbX12b2L31adXmKr8NjraYDhfaDJTp+rXhde44OGTch4JOA80c8FvgZb2uC3gRcaYy7WbhjgMFz2w4kmNp0qcgy31eShCTDMUB0jhfDhjz1kejx151KHuj4x4s2nN1ImpXCfdRJ0rK+2Z9Y5X/T+dg1hJ+0wVHW4h0k36nzRDpci4a0lUF9CoZRkN3xJ6nFD9ctk2kvO99lPqfEIESO2IsBvYbOEP7jbnbTEGEMaLMduiqq3eij5+D5jtc4To8/mvUoQrmMo5UeMhdmEzp060omdVx2qVOyYCh/G6uljaafCtHnkfpsqYoyNvUpiDNdXmQdSIF1g2xHj0M1hT7rxQU4GgI+L7N6Myi82P/7Q819HjDEcBqnRMDnl/ZDCGPOV89jPQV/yo4G6BzpGItPsYr79FHAWFs48xOrJtVIdZJ9IdpR12XTJJWzSh94kMVbNcZ7lLV9nh4zZ9Lz65781oaHW6WtSZ2/zy1UqedV4w/HoddPiAPP9d0u7AnxkzwRi5Gjbh+Aze/d75UBufc51CWyrcNDkaCCY9POwR0uMLVq0qMfRUqXngpYYd4s6ySzmYR0H63HrVwecPvmZPAEtFF7QER2GSYet426K3/El1s+ssLmxaKVHP7F/Q8rSSMxGVaUe67aYsyKMg4yZGUKJx597lYOoCtOkx9i5eOixVd3Vsf2nUMRRdktSY9aN76ib58Esz6XOSx70emU74zRnUJ1kfBjQEmOLSsTU6aoHNFSVvL0pkkJroTdyXs3NnBh9Xehimp81Mm6xyEqyzubxJTaPF/WmR4NOXkBrZ9DBpzXL7ZR+rHUPqCbW0LETc/TUISQgv/9Ztov9Dtuqrr9W5+vGNUHo1lM9Ya5UarWfHQPlmjAJCYmzUW6HDrWYZ16PR38fRgJqibHFzIiRZdjmpSfvkVZhHf5hCuHn6oZVNHwewSU2GSa+TGiX4bItX5BlSV5tcDxOyMapJcxxQj57QU8dDG11dU6SmIQZuwZQ9rL63/o71j+GOttu1Ti13VFLujo8KjZLBvAxjgMgG6sOXb/pZJ5G74TxYT/RxL5+zLEXjx/TYZUaW2JssSuEUleYuiulIEUdx7hK/jD5tGF+udh1kpcs0MXgM1e1zv/2DpssSegsW4lzlFkpMhsnpZKsO2DJsUoKDlVp7dyZ5vDx5xsuVxFjLBRHL4djnMW5FHPMNBmTCgCHsrdap3sb0mGJrfzvk3hi9NqBDsHy7TEnnF93WIkR5kaMInIV8GbgNJZyb3DT+HSfpwHvBHwZgrf7yqQici1wPXaG++uMMT8zn5GVcVj/DIcTWiWqk6rCBA3+QfHE6JNFrA5YXN4iYZznXARrW4Ry6VP7u0yQ4foQSWLLs/oSqTthDekwLrJJ1p3YfOsqYqx78Os8zeH+68ixrq1qnI1myxRJJ0YqztEHgEMRPTAR96htzLE6LlU2Rr88y/ntB3aYQ8xojjE2xvBDLqP2LSJyszHmY0G/PzPG/CvdICIJ8BrgmdgKox8QkRsj214y2rnSLVq0qIdXpZt8pu3KmHuMMb7cwDp2/vKjG47kKcDtxpg7jDEj7KyW58xyKk3RSoyXgvDqaXVTS1decljFBXXbDDMnzlhPtFeJvdRYlGgq1GiPWFEiLTlmpPncX61C7ww6Ltu2GltMMtTfum+oSlc9BFry6QW/6yTG0F457YFrKkFVOZr0ucSgvNU7qZ0VnaYZyXKhSqfqrxSVGkOJMYyrjEnaWsU+LE/nbDbGUyLyQfX7BmPMDbGOInIWW5fl/ZHVXy4ifwvcDfywMeajWAL9tOpzJ7bW/dxxWC790UKVZzS0IekHUdkZe6fsTKWVZN2VQR26pAWu7Kcqe6qJMcvKpKizKmeqZGhuS/Qe6oFMBoSHBKhV5rBvzGFQ51QJvfB+/TQ7WsxLHusz613bhAj9GKPLC+y4Mred3ogsKb+IpqrTHrHSqWH/w+iZno0Yzxljrp7WSUT6wNuAFxtjwqJ7HwI+xxizISJfg60///jmA750tMQYwywPX5UzIXbjOwliqV8uh+qh6x5nSvrzROdthd5bqokx9zp7jzOUpb1QOqyTDKtsiVXnrslOh6j0gvZppBgidJr4tip7Y93frMp+GbM9hvtKgbRLlmaMBp08f2PiHGB5Idbchqu2ixFjnWe+yvZa56TaD8zxeK4u/duAtxhj3h6u10RpjLlJRH7V1YW+C7hKdb3Stc0dLTE2RRMHQzi/Vz8gSoJIknI6qyKgG/d7jA7U0YRYCr0J4xVDUgtJsa60Qagyh6pt1XmHklEYvxcmb6178GNSXYwcY+OpQyjdNpEew2OkRcIJ733udocM6dDBlkcYVxFj7KVxlDDHcB0REeD1wG3GmFdV9DkDfMZV+XsK1hdyHpuY7/Ei8rlYQrwOeN58RlbGUfsTHQya2NP8d1WYRoMrXbZdOftVmhQhPeN0umQYSoG6yFWVPdGfo1aXY5JYaCOsywxURYxVCMlLQ0uJMTNFU8S8wzORrk044cN3Ot0uHVf3BUZF3KMenyZCnV+z6ZMXeuUP4ondoVwe/NJwDfB84FYR+bBreznwGABjzGuBbwS+V0TG7sjXGWMMMBaRFwHvxobrvMHZHueOlhhbtGhRDwORTGy725Ux76NUyyLa59XAqyvW3QTcNJ/RVKMlxipMk0RCNVkv+1kPHhVSgvdp2t2VHSverJ8lhepWHDstSztaMgxV5vDj+1fZEkPEJESYlBan2RghfrftxpmiJaiYFBjrH+ujswM1OWYKbPTyuMatNKO7bHMxdhlOlhjVkqL/DnM1zoqDcswcNofQHqMlxnlCk2IYdlFhW9JTzGKhOBCp6TuWMqmF5LcWaY85WmIq6TTzQBNibOKBrkKsf4wIY8QYkmwVaer9hVMDG8Tl7aR2htIwzRgud1lkkzGJCqBXx6mzvx4VsmmnBLaoRZObQ5NjKEmmQGrymROLbDobVSE52sMkLh9jl9Gwkxv1twcd2OhNSoAhKYbtMXti7OGPEVpdlvEmTpZZvKkxW5pv046T2Lc+rxhiBBr2rZw/HY7TaoLb6SJbvSFLy5tkbEGV80VLjFXOl7qxHzRaYmwxM2LEokmCoI0ibdUqazyKNTZZYp0VRnmAd8omS2xdXGRzY8nGI4IN0J4mGa4RD8sJHQ0xAtLjr5IM9blpYowRYZX01kT99dtXqdpVhBmS65jqfVQdu7FnfoHN3hKLy1s27GqwMDmWmGMqNobdOJP2C/OdEngk0BLjpaBK8qhSncfAQHLv8xXcwxfyMT7NVfwvHsuayy6xySLrF/oM1lYsGcYIsEpibJrZu05lrlL/QmKMvQTC8439jj380zz/+vc0Iq3CLNJkQ9vjTrrMem9I9/hwclZQTMrWx2liHqh7OewnDiNh7yEO+nK3aNHisKNVpVvMhGn2svDqOlV3nRUAVnmQH7z4y/zS8ov5JGfztGNrF1YZ3HsCztn+UyVGbUMMYxGrvLShVKslxZgqHfZv4lypUp8v9SGbx11bJbWGqjlMlRwHnOAc2L+X7x8zS2jVPpQYp91LB4kjVAxrXjjoS360EHtoYmpYjDC8je9OuP+WxwBw05P/L7aWl/gMp7mdx3L3hc8GYPDJE3Z6fB0xhjNctA2xjgxhkgRhkhRjfetCb5qqzU1JsZETRO2vrn/MHjmr2u3JK7TTjoEN9zcLibFupkuMGA+rA2aOcYxHBXtKjNOSSorIY4DfwOadSYCXugDOo4swVENjANxLTnR33PtF3PG4L4LUwJoUD9a9lElxTW2vnSlNJcOYdNg0/nBa6E2Vd7iqbRpizpRL6a9JMXTKNN23PkZIjJ4Ew3VVxDiLxFg1vtgLei/RqtLzQ8Okkj8G/I4x5j+LyBdiI9rP7tWY9gzT1DJPKv4BuNe13wt8hDz8I3+wfMhNbCpfnXQRElgTh0pTybBOQqwjxd0Som+blSCr1sWIcLfSa5joNtY/DNnyx/MEWud8mWUs+wHDPKcEHgns5WXOk0oCiMhvY5NKamI0wGVu+Tg299rhR0wK8fC/dUycltr6qo+XBkOiC/MfxmxcIQHGJLuqWMOqMJtpdsOmBLibB72KuKoIs8k+p7XtBnVkqK9lrEqgJsRY6NQsx9tPtKr0XNEkqeRPAO8Rke8HloGviu1IRF4IvBCAE4+Z9zhbtGhRh1aV3nc8F3iTMeYXReTLgd8UkS82xuzoTi4D8A0AcvZqs//DnIJQMtASgl+vpQh91b26HEpEPqOMlvSqbs46W2LKpMRYJR3OIin65UuxJ+4GMRvitD7T2mdBzHzhv8Pr7vvraxRGDsRw0E9liJYY54omSSW/E7gWwBjzlyLSw5aLum8PxzUf1KnTUFaX6ogrnJ0CZe9yTFWPkUOM7Ora9Xe4f/1dhWkPdhWBxPpO69ekreqc6oi/6d3f1CQQI8kqE8ms2K19dB5ow3Xmig8wPankp4BnAG8SkS/AUsX9ezim+WKaUT8kOpiMEfQ2xybxh3V2thgx1BFgDHVOpBjRjam+BtOcIU3bp7U1IcK6a9MUVRJx1T5jTrLYLKS6fRwmybG1Mc4HxphoUkkReSXwQWPMjcBLgF8XkR/Evpe+3SWkfHghlpXaf3vpYj8kgJj0Oc2jHvvtUUWK07AbVXca8U6TGGcln6oXQbhNlYNolrCcurHWjXG/0M6Vni9iSSWNMa9Qyx/DZvQ9uqh7YDRiXsgxk3Y/It9U/A5R5SFvEhYSI8MmY9htyMs81l8qWTY9XlOzQhNirDrGLGPdb3vfHFVpEbkKeDNw2u35BmPM9f9/e2cea0dVx/HPl0UWQSg8NOyFsAQEKdgSDCogilANixIF2UpYUiBEFhXEqAiKIIhoAtQqpKCVNaTWshalVJACBQoPqSBbWGwERAEDFkp//nHOhXnTuTNz35t7Z+7r75PcvDNnzp353rlzv+935mxtyk4A7iHM4H19zHsXGIxFnjOzfatRNpS6/xc5jtN0qu2uswQ41cwelLQm8ICk2an+za1+0OcBt6Xe/5aZjatMTRvcGKsgWZXK+w+f1bG3NYt01miTdBrKfWPJaK+TiC4vShxOlbnseUdSbriNNUUURd95+UWRYlpPP/wKK4pSzWwRsCim35C0kNC177FU0RMJKwlOqObMndEPX0n/kGcgaTMkUbZljq335c2CXdSgktfRul1Xm+S+4TSmjJThHnukpgrFjVzp7TJdl5L78kyxyBybUI2GrnXXkTQW2BG4N5W/IXAAsAfLGuOqkuZHReea2YzqlbkxVk8nEV0rnTaj5Ows6R9Ppy3N6fNlnT9NFXdFk+6sdj/qIlPMu27J7aJosuyz0HbbWcfuJZ01vgxE42oxNfZDHoKkNQgR4UnJdaQjFwGnmdnSsNrqEDY1sxclbQ78SdKgmT1VWl1JmnT7ji7KtGIm89L57UxwuOZYNf1w5+SZSaemmJXX6T+YMtXnJl7XziLGV8xsfF4BSSsTTHG6md2QUWQ8cHU0xQFgoqQlZjbDzF4EMLOnJc0hRJxujH1FXheXJOmIo1Umy1jzflxV/dg6fY7ZCzqJmIYbIY4kvx1VfEd1RostKtKg4HaXAQvN7MKsMma2WaL8NGCWmc2QNAZ408wWSxog9Gj5STXKhuLG6DhOPtWOfNkVOAwYlLQg5p0BbAJgZlNy3rsN8EtJS4EVCM8Y0402leDG2CuKWq7zKIpe0h2Ms6rqnfTdKxMddFJ2JHR6/G5EgJ1cj073Fe1vQrRYYXcdM7sLWObBYU75SYn0X4Dtq1GSjxtjr6nKFNuVyaqq53Xb6aQjc/qcVd89IzWBTp8pdksHjLxDe1U6qsAnkXBqIa8VsopvKK9hYCQ3/Eg0duOH1k1jTTKS76TovU00oKX4RLVOQ+j2N5P3A+xGVNspdRhEnR3YoZmm2MInkXAcx0kx+qZ2ycWN0VmWJkcu/UI/PUN0lsGN0XGqouyvyQ2x8bgxOk5VuOGNGtwYHccpYPlrlnZjdByngOVv0Rc3RsdxClj+eni7MTqOU4BHjI7jOCncGB3HcVIY3vjiOI4zBH/G6DiOk8Kr0o7jOCk8YnQcx0mx/EWMK9QtwHGcptOKGMu88pG0saQ7JD0m6a+Svp5TdoKkJZIOTOQdIenv8XXEiD5WDh4xOo5TQKVDApcAp5rZg5LWBB6QNDu9doukFYHzgNsSeesA3yesImjxvTPN7N9ViWvhEaPjOAW0qtJlXgVHMltkZg/G9BvAQmDDjKInEpZYfSmR93lgtpm9Gs1wNrD38D5TPh4xOo5TgtKNLwOS5ie2p5rZ1KyCksYS1oW+N5W/IXAAsAcwIbFrQ+D5xPYLZJvqiHFjdByngI4aX14xs/FFhSStQYgITzKz11O7LwJOM7OlYRnq3uPG6DhOAdW2SktamWCK083showi44GroykOABMlLQFeBHZPlNsImFOZsARujI7jFFBdP0YFt7sMWGhmF2aezWyzRPlpwCwzmxEbX86RNCbu3gv4diXCUrgxOo5TQKWt0rsChwGDkhbEvDOATQDMbEq7N5rZq5LOBu6PWWeZ2atVCUvixug4TgHVVaXN7C6g9INDM5uU2r4cuLwSMTm4MTqOU8DyNySwa/0YJV0u6SVJj7bZL0m/kPSkpEck7dQtLY7jjITq+jH2C93s4D2N/M6X+wBbxtexwKVd1OI4zrCpbkhgv9C1qrSZzY0dONuxH3ClmRkwT9LaktY3s0Xd0uQ4znDwVQJ7Sbte7MsYo6RjCVElwGKOzq6eN5QB4JW6RZSkn7RCf+ntJ60AW7+fXHQrnDlQ8n399Bnb0heNL3FI0VQASfPL9KxvCv2kt5+0Qn/p7SetEPS20mbWlfHITabOSSReBDZObG8U8xzHcWqlTmOcCRweW6d3AV7z54uO4zSBrlWlJV1FGNc4IOkFwjxqK8N7vdtvAiYCTwJvAkeWPHTmTB0Npp/09pNW6C+9/aQV+k9vpSg0CjuO4zgtfKJax3GcFG6MjuM4KRprjJL2lvR4HDJ4esb+VSRdE/ffW9CZvKuU0HpKXPznEUl/lLRpHToTenL1Jsp9WZJJqq2bSRmtkr6SWFzpd73WmNJSdC9sEheDeijeDxPr0Bm1+LDddphZ417AisBTwObAB4CHgW1TZY4HpsT0QcA1Dda6B7B6TB9Xl9ayemO5NYG5wDxgfFO1EoaUPgSMidsfbvK1JTRqHBfT2wLP1qj308BOwKNt9k8EbibMhrMLcG9dWnv9amrEuDPwpJk9bWZvA1cThhAm2Q+4IqavB/ZUPfOgF2o1szvM7M24OY/QZ7MuylxbgLMJq7T9r5fiUpTRegxwscWV4szsJeqjjF4DPhTTawH/6KG+oULM5gJ58xm+N2zXzOYBa0tavzfq6qWpxlhm0Zv3ypjZEuA1YN2eqGujI1K0QM9RhP/CdVGoN1aZNjazG3spLIMy13YrYCtJd0uaJ6nOURpl9J4JHBq7sN1EWA2vqfRs8amm0RdDAkcLkg4lrGexW91a2iFpBeBCYFLNUsqyEqE6vTshEp8raXsz+0+tqtpzMDDNzH4q6RPAbyRtZ2ZL6xbmvE9TI8YywwXfKyNpJUK15F89UddGRyRzaKOkzwLfAfY1s8U90pZFkd41ge2AOZKeJTxbmllTA0yZa/sCMNPM3jGzZ4AnCEZZB2X0HgVcC2Bm9wCrEiaYaCLL77Dduh9ytnnouxLwNLAZ7z/E/miqzAkMbXy5tsFadyQ8lN+yH65tqvwc6mt8KXNt9wauiOkBjVkNPgAAA/1JREFUQtVv3QbrvRmYFNPbEJ4xqsb7YSztG1++wNDGl/vq0tnz61K3gJwvbCLhv/9TwHdi3lmEiAvCf9rrCEMK7wM2b7DW24F/Agvia2aTr22qbG3GWPLailD1fwwYBA5q8rUltETfHU1zAbBXjVqvIkzz9w4h8j4KmAxMTlzbi+NnGazzPuj1y4cEOo7jpGjqM0bHcZzacGN0HMdJ4cboOI6Two3RcRwnhRuj4zhOCjfGUYSkjSU9I2mduD0mbo/t0vkmSzo8pidJ2iCx79eStq3oPPtL+l5MT5N04DCPs56kW6rQ5Ixu3BhHEWb2PHApcG7MOheYambPdul8U8zsyrg5Cdggse9oM3usolN9C7hkpAcxs5eBRZJ2HbkkZzTjxjj6+Bmwi6STgE8CF6QLSBor6W+SpktaKOl6SavHfXvGuQIH43x9q8T8cxNzSl4Q886U9I0YwY0HpktaIGk1SXNawwglHRyP96ik8xI6/ivpR5IejhNAfCRD61bAYjNbZr1iSWfHCHJFSc9K+nE8/3xJO0m6VdJTkiYn3jYDOGT4l9dZHnBjHGWY2TvANwkGeVLczmJr4BIz2wZ4HThe0qrANOCrZrY9YYjbcZLWBQ4gDG/7GPDD1DmvB+YDh5jZODN7q7UvVq/PAz4DjAMmSNo/7v4gMM/MdiDM/XhMhs5dgQfTmZLOB9YDjjSzd2P2c2Y2Dvhz/BwHEoay/SDx1vnAp9pcE8cB3BhHK/sQhnptl1PmeTO7O6Z/S4gutwaeMbMnYv4VhMlMXyPMy3iZpC8RVnUsywRgjpm9bGF6uOnxmABvA7Ni+gHCuN006wMvp/K+C6xlZpNt6NCtmfHvIGFS1Tdi9XmxpLXjvpdIVPkdJws3xlGGpHHA5wiR0sk5E4umx4K2HRsaDW1nwoTAXwSqasB4J2Fs75I9Dd5bhHHxSe4HPt5qZErQmrVoaSLd2m4de9V4TMdpixvjKCLOYH4poQr9HHA+Gc8YI5vE+QABvgbcBTwOjJW0Rcw/DLhT0hqECO0m4GRgh4zjvUGYsizNfcBukgYkrUiYj/DODj7WQmCLVN4thIalGyVlnTOPrYDMNU4cp4Ub4+jiGMJzttlx+xJgG0lZE+M+DpwgaSEwBrjUzP4HHAlcJ2mQEGlNIRjeLEmPEAz0lIzjTQOmtBpfWplmtgg4HbiDMKPMA2b2+w4+01xgx/SyFWZ2HfArwlyRq2W+M5s9gLpnJncajs+usxwS+zXOMrO8Z5CNQdLPgT+Y2e0VHGsusJ/FNWIcJwuPGJ1+4Bxg9ZEeRNJ6wIVuik4RHjE6juOk8IjRcRwnhRuj4zhOCjdGx3GcFG6MjuM4KdwYHcdxUvwfHbVLWCnevIkAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "#NBVAL_SKIP\n", "\n", "# Show what the update does to the model\n", "from examples.seismic import plot_image, plot_velocity\n", "\n", "model0.vp = result.x.astype(np.float32).reshape(model0.vp.data.shape)\n", "plot_velocity(model0)" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAASQAAAD8CAYAAADe49kaAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDMuMC4yLCBodHRwOi8vbWF0cGxvdGxpYi5vcmcvOIA7rQAAIABJREFUeJztvXuwZVd93/n59dW9unTTqBGSxUuOhAMkgCeyrQLPOLZJiD0yY6I45TBoKrZ5JLKnTOJkXGWDkwyuOKkhE2PiKWY8aQ+KYIKFsTEx5WF4jB1CpcYQJMyYlwkgwEjo1Wq31HTTfVv3rvlj75/2d//OWvucc++59+x77/pWnTrnrLP32mvvs9d3/97LUkpUVFRUjAFHlj2AioqKCkclpIqKitGgElJFRcVoUAmpoqJiNKiEVFFRMRpUQqqoqBgNdo2QzOwmM/u8mX3RzF63W8epqKg4OLDdiEMysxXgPwM/ANwDfBy4JaX02YUfrKKi4sBgtySkFwJfTCndnVLaAN4J3LxLx6qoqDgguGyX+n0G8DX5fg/wotLGZpaqMauiYnexBadSSlcD3HTTTenUqVMz7XfXXXd9IKV0064OrsVuEdJUmNmtwK0ABqwvayAVFYcE5+Gr/vnUqVPceeedM+1nZlft2qACdouQ7gWule/PbNseR0rpJHASYMWsJtRVVOwpEvDYsgcxgd0ipI8Dzzaz62mI6BXAf7dLx6qoqJgbCbiw7EFMYFcIKaX0mJm9FvgAsALcllL6zG4cq6KiYjs4XBISKaX3Ae/brf4rKip2gkNGSBUVFWNGJaSKiorRoBJSRUXFqFAJqaKiYhTYAi4uexATqIRUUXEoMU6VrWZsVFQcWjw242sYZnabmT1oZp/O/PazZpZmjfauhFRRcSjhEtLOCQm4HZjIdTOza4EfBP501lFVQqqoOJRYHCGllD4CnM789Gbg59qDzYRqQ6qoOJTYYjdTR8zsZuDelNL/Z2Yz71cJqaLi0GJmo/ZVZqalAU62yfFZmNlR4Bdo1LW5UAmpouJQYi4v26mU0o1zdP5twPWAS0fPBD5hZi9MKd0/tGMlpIqKQ4ndc/unlD4FfIt/N7OvADemlKZWhKtG7YqKQ4nFGbXN7A7gD4Hnmtk9Zvaa7Y6qSkgVFYcSi5OQUkq3TPn9uln7qoRUUXEocYgKtFVUVIwd40wdqYRUUXEokYDNZQ9iApWQKioOJaqEVFFRMSpUQqqoqBgFdjd1ZLuohFRRcSgxTpVt24GRZnatmf17M/usmX3GzH6mbf9FM7vXzD7Zvl66uOFWVFQsBgstP7Iw7ERCegz42ZTSJ8zsOHCXmX2o/e3NKaVf3vnwDi5W5th2fL6QioOB8UlI2yaklNJ9wH3t57Nm9jngGYsa2EHAPKSz11jm2JZFsEPnfPhIf5wq20JsSGZ2HfAdwMeA7wFea2Y/DtxJI0X92SKOM2bkbvZZ9eGtQvt2JsmsRLPMJMZ47NL5KxZBGJuUr88s1+1gkdYBJSQzeyLwbuAfpJQeNbNfA36J5ox/CXgT8OrMfrcCtwLMXr5pfPAbed4JHifhtJt90UQzJuktjiV3LfS8cgQ2K1nodvNeA91+/5PTAfSymdkqDRm9I6X0OwAppQfk918Hfi+3b1vg6STAitnMJS4rKioWhQMkIVlTeemtwOdSSr8i7U9r7UsAPwJMrERwELBC/6mde9rGp+hWoX0l099Qv4vCTvv28c6ics2DaZLISqb9SGEcQ5LMNGlpHql3/0lMB09l+x7gx4BPmdkn27ZfAG4xsxtozvgrwE/uaIQjgt60R+R7boJs0kyQ3E2/Gvb1/uIxhpCbAIsmhlkn5DLsUX6dctdX20qkmfu/Yj9bTJ6b/9dRhczdA+PGASOklNJ/JG/+ed/2hzNeuASzEtoilIhWaMgH+ZwjoKGnc+mpv0p/ApQmxU5Iarv7LpKgpvWVsy3lJKx4bVSiKpHTimwzTRqOGD85HTBCOkxQMoqE4uRzqf3uZOEv3T6nmkVJKYfN8Nknk06ykpQ2NPFKbSXsRHqbtw/HPOSW21b/p3iNVNrJSVA5tTDnqdPvSorjJ6XxjbAS0hREMopP5A0aMnJJaL3dRskI8pKRInejD90uJZtJRI6MNsNvsxwPOtJVzEsw83oTc9dlXltbVLGG+s5t4++6fXyPfY1fjTuAXraKior9iqqy7StMU7VUOlqnk5BWZRt9OuakC+9/2hhKmObF0+/R3lTyXs0yrni8WcY2rb3UZ0lN0t9jv7NKbUMSzJHwPqT2+v2xv6SkSkj7Bko+kYz8xlQyWidvU1LbUuw/onTTDql6pYlX8vjlft+OTasE7Xd1yu+K0tj0e26baEPL9VfCZnjP7XMk8/kI3X+rUFLS8YyXlBZHSGZ2G/DDwIMppRe0bf8SeBnNVPkS8KqU0plpfVVCCsjZCPTG9BsrSkZOPBdobsL18PssNiSFE1sMH9D20sTIobRdznM4FIIwFDFdkg5LY8uRSan/aTawS8wm1U3rO6J0LXx/feCo08P7VAfEtGPtPRYmId0OvAV4u7R9CHh9SukxM/sXwOuBn5/WUSWkDEo39iU6M6B70WjbvP0YcJSOjJSQciRQmqRKRCppuaqo7z42JbBcnznj8CxkVJqUQ+Q35I7XyZrbJo631J/2FT/PilLYhfen5Lcq25fU8Xg9S1685WNxRu2U0kfafFZt+6B8/Sjwo7P0VQlJMOQV26K58bx9rX13kjrefj/e/laa7Gs0N7b36/3kSDCqf368S3Rk5JNhQ7aLT++cC7uEXDBgCbNIPKUgwrjfrLa0aKfJqYbzYsiOpLYj/0+PsL0EYRiT+ranNqRXA785y4aVkCoqDiXmIqSrzOxO+X6yzUWdCjP7R+2B3jHL9pWQAqJko0ZsD3rU33K2pPN0T0GNTzom/Xm795nz5Pkxourg0pC/fHw+ng0Zi6p+Q54ib9fo5ByGpJWcnWQoQrqEIWkpl8qxHcwai+X/kxuzNfpeEZWfUuzUuFS3mQnpVErpxnl7N7NX0hi7X5JSmimBvhJSi1mMzNHe4CR1BHikbVsHrqG5ac8DZ9v2r9MQ0nEaG9Mx2T4X2a2qgQda+jhWmfTiORFpkKZuo9vOkjahZKNq3DwqXa7/WY3JjiHD+qw2LT1eNDbn1FvvKxKrE77aBjXsoBTaETEOtW13VTYzuwn4OeD7U0rnZ92vElIGJY9XlA50cj29ff924EXAw8BXaYgI4KG27X46wzft58tpbEkxnmmFyYhv5Ptq+KxEBX2SUvtXtEvNahh2MtJJOwtB7bT201BYwiw5ZrnjKTnptfD/YL39fka2X6WRhJT0/RqsMRnftcwieNOxULf/HcCLaVS7e4A30HjVLgc+1BQG4aMppZ+a1tehJ6TSTVx6ksen4DpwAnhW+/37gddeAb/zSP8GdSnnGI3U9FDbfqZt8zCBy9t2JahIPLlxu6rppOTb+vcN+tLUCp2aVwr+i0QVg/9yBDWE7UhEs3oCtxNE6urvJh0BufR6nOZ/vbttv5/OIaEPoyjJ5q7dOIlpoV62WzLNb91OX4eekCLmSUR1AjhO52V7FPitR8Cr1PmN7r/709e/n2tfTkyXy376Uslpjbz9xieHSkI+RiWqDeAifbKCPtkOnX9Ua5SgYPsTcEjqmRYyMavXS89Bg1tVYj0OXAk8he7B8RCdl/Qc/QeNS0k6Hpe6ZpUil4PlK44RlZAqKg4laurIvsDQM0NVJ7fDuATiT+QHaKSPCzRPUscanSF7lb4kdEG2PyftHmDpKp3vuxXeoS81qIfPpSD1DvlT3duPhH0vhX6G0iv89yHvW0k6KNmKSiqaf89JUoqobucM+25Xc4k1qmx+zaPx2u8B76ekNqp3dJrncjmohDR6lGwhbtBcZ5KQHK72nKWzMWzI72pDWqELiLxA38t2XtrPtsfQAMjL6WweMVXD+1SPnE9gVdl8QuW8chvyeyQqx1Z4H4JOxkg403LohmxFQ2QUwySUkNxqEsnICcntePqf6LGP0D1EvE1JK0ahj08pclRCGiWGonShm9QuqXi7ppH4dt5fKaHWJ0ckklwk9REacjpPn9w8lkmlJJiUlPQctGKBbxMN3n5Ouu1maPfx5mJq9PNQhn7OKB1JZh4CUsT8Px+vRrerwT/a6aBvs4thHj4e97ZBd+2GchRL12L5qIQ0ekQVwaWRaFy+QGdchulqyVB6hKp9GjB5lE6V0ye+T7DLpX2NLnxAJb2SSpFLg/DtNVQAafPPQ7FFOWj/uXIuQ678OG5FlNqUjDakPScN6n96OX2VeI0u/SeOYUX2h07FLl3rUqmX5aMWaNtX0EnmT9Oj9G86VYMccbKpeuOvnEQEnRcn9uWxL9B/2q8z+cR317S35RYU0M9r4fdIGoS2IVVOkSORSEa5/n27IeQks0hG8Xqp2g1lL2bMNSzZuZTYvMJDbvvxoqpsFRUVY0Ian+y2iJVrv0Jny30spXSjmV1Jk917Hc1SSC/fT8tpq5i9Gl4uwUQVTJ/upSe8PsVjvFMulmeN5mmvkpjbri7Q2Jf0ie+xRap+qJSUk0CG1Mh4TqrS6VhLklJOEvI2VYnidro/oW+/7rno75jnB5PSUZQe/T/VaxYlRB1Lzq4EnXPCoQbu0WJ8rr+FSUh/JaV0Sr6/Dvj9lNIbzex17fepxZnGiJXwiu3RKO1kE70sJTXNbT6ar+b9bNLYN2KZEbctXaJLbVAjraoTPuFibaZS2Y44jiEoOUSSiGSjn2N5ltXMttOOC5OF6zzyXK+Xk9Ea/QeKf1cDNrJNJG8lmJjOo168eI/kHByjIKrESAbSx26pbDfT5LYAvA34MPuEkKLEUnoi+m/+NHbz4AW6yZaL7M5FMw/F90Sbk2b5b8hxPfbpLHmPUS7aWw3bEbNITjoRZ5W09Lg5Iipd61I+nL6rAVvJQ8koGq/V46hjyRncNfJdx7NBWWobLRKzZwPvIRZBSAn4oJkl4F+3dVKukeW076dJgN+3GHLNQ0MCntXv9bWjFOAxTBojFI+RO26OkFQ10bgaV+X8PjvDbJJTyZicIxklHw2gzMUp5Sa1kpEeN2d4L01oVXejUduvixrto9pN+L4S2qNjQhEdD34Mf9jEfZzEcqESS8UBlpD+ckrpXjP7FprM3j/RH1NKqSWrHszsVuBWyC9/u0zkXP8+4del3T1vZ+nKj6i4H0nA+9L3XFQx8lt0Y7tdySeGH8ttTTG5VgMsNQByXdq0amVJIiyRh5KmjhvZXvuYxe2v28c5E9W1KDl63FgknlhSOJJUSSKOKlskU+geCiUpMdrsRsMDo2HHDjvO+Usp3du+Pwi8B3gh8ICZPQ2gfX8ws9/JlNKNKaUbx0ZIFRUHHi4hzfLaQ+xIQjKzY8CRlNLZ9vMPAv8UeC/wE8Ab2/ff3elA9wI5yUBtE6oKrNBl7Z+mMy5fTiM1HaULVoS+mJ9LbdAUB5iUOFbC5xihPe28LmSO6eqMSlRaSSAnsajUsipjj2MeGstQRLYiZ7z2dv1fkM8eyFqShGK8UU6VK3lK1d6kkqn+rva1zfAbhbbYvmcYjajWYacq2zXAe9oCTJcBv5FSer+ZfRx4l5m9hqZO2ct3eJw9h0qzHgyYc6mfaLe9v20/zmSwInQ2pKiC6YTKqSBqOEa+5wjMPUw58soFcZZQcr8fCe9+jtsx5g4FEEby0XcNLlXjtZYYjv+T2sw038+/D6lsDrcPlYJNneDjPtFDOxoOSIxSZdsRIaWU7gb+Uqb9YeAlO+l72YgTIer+akvQGJRLme1yfTs5XKBvjFYjtZJSznYS3d458spJ3ZHMZkE0TM+KWchqmh0tLvfkhusoIa0wGeIQPWzRsK/GfejHXPkxfFzr0k+MwdqiL1npfTNLKMCeI9EX70aCGqndYugJpoGMGq8S3ejQTRZ3w/sNqr/7b9Dlql0I7Rrgl1OLSkbx3OR2AsoZpkvnGn/XEIEcOcXtS+OLap4ed5PJc42qmS5moAGQqqqpJBTLA+ckp/hf5lQpraXtjgX/PV7feL7RQZJDzoC/6zhoEtJBQrwZ4k0ZY3b0KYq0X6Lxaq3S5ZtpfzF+6Dxd1PUFuhv9IpMSzyzIEYmrGpGQcnFIJduJbxcjrCNJKSIpOcH45CtJCjn1VFVcjcNSe1bJVpSLxRpaXEGPrRHfumy6buclW3KxZdHuNhocYLd/RUXFfkSVkPYHcnFI0Tjp7f4E9yfv6fb9HP0ncunJG21I8cl+NOyvx94ucoGK/j3GCkEnZUUD8Gpm2zhGVcF8W1V34rvaxWIQqLdNKycS6xvlgkNj2kjOWxmdBipNuYTrKq6367kM2duGjPp7ggVKSGZ2G836aw+mlF7Qtm0rn7USUkCcpJekLYrkWtzMy9O6fcFXFtGaQn7TRzuJTxZfhhv6Ez13Uw+lm8yK6D3TKOWcjSyqN3Ey50hTz3UjfFZv4BD5bIY2JSP1mpXc+EpGUWXLBTrqQ2Mjs70iRoZrH6oST/tf9lx7WqzKdjvwFuDt0ratfNZKSAE5e4iSgz69/bOuIuI3oJe/8gJeZ+mkIJUo/EYfss/A5JhKT9hZ4pIUOSO1j6MUs6MkULJFaaiB3vcbsn3OhqQOBCWfadn7WgtqKCI7hm3kpF4Nzbgo2+v/o949tRVuyXtJ6hoFFpjLllL6iJldF5q3lc9aCUmQ8674zR8NsRqYeIwmHgm6Im4XaIIlNYcpxiB5/y7e65O25M1ihnZHlKxypKbbxMTXWQhpmsp2hMn7fiu0xTAFVeli8vIKfWnHicCDUKPKFgkpqqe5cSsR5iSkGC+m0m70+vn5KlGNBrvLjtvKZ62ENAM01kRvuBUaAjpBs44X7feLNJKRLqX9KF1OmS4WoNCJpzd5lJScMGIJj5KKl4uFKeWUqRoT1Z5YTXGNMiGpyuLYCu0lyUQlJL8mUaLMSTwa6BjJdIXJeKOSez8XVqDnHcMQYiwY5O2Qo8J8gZFXmdmd8v1km0Q/26EK+aw5VEKqqDismJ0lT6WUbpyz9wfM7GkppftK+aw5VEISxAjh+Js+DV06WqVbNgcaqeh+4Gs0a7R5jpuviBrrXpee1DHyWm0o7pVzdRLy9hOVpvy7wp/i2scR+gGGhH5jze4ho3uMMvegyJL9S1NC/FxdjZ0W0KjXVasgqMSXq3NEGKMfP0pImvqTU91dslOJSj2rOe/bPDmAC8fup45sK5+1ElIGOTJSbxn0y5Bcols6+2vAp4Cn0qhxV7ftbizOlcGIK1z4GHRiaGS3vi5Keww18L6VuPwclIwi8eSimnOEF43wOrE8/y+mspTUSiVg9cRpBDZhbEN2rkj4JRuXv0c1Ue1Ifg10v3hOMZRB748Y6b/jEhuLwAJTR8zsDhoD9lVmdg/wBhoimjuftRJSQO5ppU9DlTY8ReQh4Ott+yrw92kI4gE6G5LeqDFOxtMb4jiUjJyQPM3kfPv5YmjXbT3nSqU671vPQyd7jNuBSYkkhiZEqJEe+lLELI4dlRJLbvzoZYsu/JzhP0pFOdugktFF+ss/+ZiUYJB2va6a0xYjykeDBUlIKaVbCj+9ZN6+Dj0h5YyxEf5E1JvpCJ0nzQ3YAK8C/tmVcNvp5oZW46erHsfpAh6fREdI8SmshOThA05E59s+/LiXkycplzB0kkYyioR0lPkJKcZW6fFU4lFvV0lK0n5yqtmQhJQzyJekEiUjJx5P6/E0npy0rPeIHy/GRcUYq5XQHvvZc9TUkfFBySgXma2/xaf7Kt2Nexr48237WeBfn4a76cjK+9Q8KJ/sl9MnJYWWqY11snN2HrednJM+ztGpFD5J/fx8H+3b17Y/yuQxIzlAnpBivp7f92v0A03jJPVtdZ6o98y/5ySk6CGM7n1HJJioml0KLx+LSjyboV3HEO8jJzrdRq/VUlEJqaKiYhQ4iPWQDgJKRllHLhgO+Rzjiu6nK/p/hsnFiteZFPk1GFGf5q4K6FM+RkVrew7RUKzHdA+hp724Knm8bYt2rrhQgSJm5OdicFxicKnHVcs4drV5RU+gHz8nIWqUebQhleK+og3JVTX18Om4VuiXZ4n/m8Yqeb/xfx0NqoQ0bqj6oRO+lL4Bzc17TL6foUsTOU/fKLrCpOHZi7vlXOJOhP4b9A3VOcOqGrF98nuIQCRCV9e8zK6TkROS27lUTYzXw8ehNrY4kVWVc7UtuuC1MJpD2zTfLnoDYTJIs0TYsRRKVC/dUeFpPiWydPh13Ay/O/n670NBkiuF9l3FAV4GaV8jkkyuEJu3azpAlBJipK8Sg/+uBlM3RrvdZ0X2g3Laxxp9QlKPlE8uTXfwKOfcU90ntRPqMXkdl3aXjqLh3W1rPrF18vu1irE5Li1oXyvyysG39z40IltTR2Kdp4g44ZWUYvUFt92pXU/JTAkSOgN4/D9y9ZNyY9lzVKP2+DBUGiLnfVuhkyaOtW3uSVNpxdUH6C877cSgXjO9UY8xmaqQUyVdhVH1T5M8XRUqTfKYrxZd6rFUR0wbyfUZCVgllag+RWkoVwIF+kGEUULSvDVvj3FRDpUmc4GM6mVTiSmG6ZRUesdG+Bw9ft5HDMNYGqoNqaKiYhSoEtL4EA29Gm8UdX4Xv12NOUr3hD3NZEiAP71j2Qz/fF6213SFXOxPDlH6yNlPSgbwlbD9NPe7nlvM3leVTcvvupSh7dPu/1kkVlXVYqxUTl3TAMg4Bj/HGPHuUq/akLRmk48pIoYPaPyXBonmBJM954aDRkhm9lyainCOZwH/I03y+9+lCWAG+IWU0vu2PcI9gor1ObvSJt2aa9DEDjnciO3brsl2eoPGGxaam18n9TH53SdeKYAwYohgoD8pNTVFo8C13o+3a9BhrHnkkzp62TTCXA345+nbabQf7U/PR71smu6SC97MnbNWaNB2tR9Flc1tQhp/pERcurZ6D6naG4+7dBw0o3ZK6fPADQBmtgLcS7Ny7auAN6eUfnkhI9xl6NMqerm8zaUmnQCeWAvNRPbStZ67Fm9GzT3TCGD1gMVkTpewosSGbJubZAqXGJQM3YalIQ3e53n6Ln7IT3ToT8I48TUpWD1YSlIuJWraS7wGMRctlhiZRkjRrqOGabUhKWnqtVUJKQY5Eo6n25TqJ43GfgQH2ob0EuBLKaWvtotG7gsoAWnkbbzJXPxWle043Q3n0drQl5bipF6jM2arsTsXD6MTRaUkHW/cPudO9ie1hyPoeV1o23JqYpSGFJHEo3SwGT5Hw3spYdjrR8VUlpizpikssb1kcI/CQPRKxuvuJOv9eaVLmFSVkW28n2NM/lcwIg44aCpbwCuAO+T7a83sx4E7gZ+dpbj3sqD/SW4i+oRyqcFVMU2vuERXZuQhmiz/K+hLUeuhvzhJoZOG/LOTkS5E6eOLsUhq39CJo544J8OH6Ca+2qlKT+0oCek1y5XQmOXpr2PXutSasuIhCNGbFl9+TP3vdOIPBWlqmoi35ySamI+WW6VWk69dqo0hICUpfCkYISHtuBKCma0Bfx34rbbp14Bvo1Hn7gPeVNjvVjO708zunKmUXEVFxeLgqSOzvPYQi5CQfgj4RErpAQB/BzCzXwd+L7dTWwLzJMDKjOUtdxP+hM0ZIXNeHU9C9UhoL2H7EI36doK+yrYSXioJQT8uxts1i7yUwKlShj/w/OntaoOrml73Ww3R2ke8DqomfjxN/4t+oFXX1UCtqqR/L0lVsT6TpqzkbEV6HY4wrGJGG10piTbGlOUkJI2c9yBQ71+z/T3uLBrEtTjbUjFCCWkRhHQLoq552cr2648An17AMXYNObc59NUUDwyME8Mnx0XgKW37lTSk5DYmZHvtN0IJwr9rPpWqAiX7TvSurdKpmLltY+CepmO8Od1H38rzSeAxeSHvAJfxofSHj3+Gy/hn9h09D5YajN2YD32bjY5fKyPEIM2cJ1RJwwlAbUQlQorXXhG9ov45jjWmz6h3LnduS8dB87IBmNkx4AeAn5Tm/9nMbqA55a+E30YLjQiON80Ryk9qjZqGRgo5TZNkqzE7ak+aJRZFb+QYR+NEEmOMtE6Rt6m0U5I81GZzFHhZegvNyjWxLNxj0gaRkLrbaR24jH+cXs9H7X+aWDbce9PYH41ZGrJL+Xn4+UfjeZRIvNCa1jjydpVscnFKGj/k5+D/dzyuv9ZlX/0d2R/5bWk4iEbtlNI5OuHA235sRyPaY8TJGwPZoB/3EtMAont4nYaU7m+/+438JCZJyfvYIj8Jh8bpbaXysz4WLbqmQZ2ar+Zryl2X/hbwDeA32vdvyFn4SySk1I7ajjBJSM3ru9N/ATwRgA37fzlHY1D3YnLIdyUr6CSpaBRWxO8uealko/FPSoIqIeU8lEek3cnS2zQEoORBiw+V0WEUemMfhzpSO0Ijm/VJpqpNfPL5NPX/dp3Gw7ZJIympdOAEUErM9WPlxqQkpOSjEptLOEo8x+XlwZxXtN9PAHYFcObp7S//Fzx2vh+9CNNnL1uwstG8AFbPZw1vaxdh7bvgyV6fpY1D+OalrlyLLh3llTFzQ4mSjapeSkgx7slJ8GLY37ESPqs3zVXCnB1siChzGQFLx0GUkCoqKvYpKiHtD+Tq3+SM3efp1IFz0r5KI4n4Phqf5EJGLL1K5rOPJRqeh9SwJ9FJQ9BIQFfQLWTpXja7Un74MHBvu0TBo/QLc7tY4oaY6J5TS60a3zRyUQsuHQPeBfwTOpEIeMIj8IQz8C1n4JsXOwnpDPBIu2lckvw8fVVOU3PUeB1XadFTcntV1FyiNOqnqkGsMRg0qnj+m3rjRocFqmxm9g+Bv0NDdZ8CXpVSivUJp+LQE1LuP1ECUlyiH+2s6prakJRA/IY8TT9XTAlMvWOxaFms++NEtE6/ZtEVNIR0BX2P3wngKUfaL/qDf/+CnJTP/CE9aSiRzgcfcztirMS30yMkzrQX6Aw84TQ84eGm+cQ3ej89fr0eoSOFuOx2JB7nVTVuQ0de8VSirU5tRTFswbfP3UPR7hX/96VjgV42M3sGzWI7z0spfdPM3kUTLH37vH0dekJyTDMub9DM1ejC3qIr1E/7WfPUNMI4pkw4NMLX4WSkHjGS86Q/AAAgAElEQVTo24miNOQE5Pl0VwJPeCIN8VwdfngKTTjr15hOSGo/0ps4+rJV7CsZuo4C19JcyEfadmedh/sntfYQfMtpWL2UT+qNQ1HjtUtUWmwtuv3VDhUfQFHo22LS2eFQO19J6IhR9EsnpsWrbJcBTzCzSzT/8tenbF/s5FAjF+AG/RvLS2mcoU827ua9nE4r8VQPTdz042iiZbyxo/E6zumch8xVMtrPT6G/OOXaFXREFAnp+2n0SLckQ0dIF9t3nb05kSKKCqXBO1F5Hohb268NJxXjEFrx8ckPwdal7pBOMDF3zdNwlEtVtY7u/RJ5KBlpak5p23j63n/cRsc5CiyIkFJK95rZLwN/CnwT+GBK6YPb6evQE5JCc7XivHOnkM9dnxBPpR+ot9Vu8zDdfPd2Lz8bbQtaGTInZGgSqc9nf7lk5uR0JbD2xLbxyswLGuZylUm8XY/bj9w/HvWbUjJbruCSD17FRDegeX9/QS6CB0tFUag95ok2/t+ln1yApKptuRinnIrm0H5ysWh+2lFyiuq2EtcoXf2O+VYducrM7pTvJ9tMCwDM7MnAzcD1NHfVb5nZ304p/dt5h1UJqaLisGJ2CelUSunGgd//GvDllNJDAGb2O8B/BVRC2i6iXUIfHmqH0DSLs3R2Yt/vIRqzzNeYTPlwxFrNKh3FEhcxq11VOF0ZxLWgJ7ibj/Y9vmg3PE8jRqgkpEWbog1JDS7TxAy3Brtop5GDuq8bwFRPzrnIzsFKK8UdP99cdw8yjbY3Vem0i1yKSCkqP6d6ubCm/4tvr/v50EuS02i8bYtNHflT4LvN7CiNyvYSmkofc6MSUouhwDa/wU/Quc61BMgGXWT2Z2lI6jvb3x9t291drTYl78fJqBQCoBHZukpJLhygZ9C6PLx84xPtwHLFgHJ6aySjXPZpzONYo2PeOAv9BH63/f5DdBfZdVHoiPNs17Z+frJuE0ySiIYD5HLWdEilxQEUmhqSO25uv1zC72jsRws0aqeUPmZmvw18giaM/49oE+fnRSUkgZJEnHPH6HvON2mU5Us0zqKvtO3fSSO/rtAsp/2Ftn2LztitZpf4hFfknuDe7pMiZpf31p1ek5c+1tWlmMt3GJKC4rb6OSdaqNi5It836IjTJalSBrO4GdeOwPpWfxUU71pPTYm/ZD+axcajp6H2I5WMpvVV8s4tHQtkx5TSG4A37LSfQ09IpSdr9JL4Cq9H5fc1GpJ5mMa4DfDfAy97Ltz/+UYAaENqHi8wHr14paC8UlusVqgTYjU2HpHXvBbW7cyiaTM+ujH95GPFNWVZDcZqtz1ysSPkkpfLEQW++H+XqiZEYVDVtegNzWE7l3xPMdJI7dGotGNC7n/SWD/NtN+k0SqubV8vOwr8yXN46l9sPV5MVhf0Y+RCDGZFvOGzf2RuxsaQ82VAWV/HEYnUxxuYYCX8NM9NPIsAOPTbvBNm1KS0OeNrD3HoJaSKikOJ+dz+e4ZKSDMiivBqy12nWXIF4P3n4aZv/888+Lkm+Dh61GY5Rq699GTWJ7BFKSi3Uez4gGHUEsmYkJjv5twjVEJqoQSjRm0nCU9L0Bi/SzTq2xU0BmyAtwJ3tzUy76bLjsi5g3PHj3WSCO0lzDURIxFNqxY3K6YNIvrCY5h1NKzpUyAY92ZVq2LEdQ5DhB/7nPdSjZrzq4Q0PkTjda5W0Qpd9LUS13m6pPnnte0foPG4Xd1u64SkFQf1vTSWaQXb4j6bNPXSLFrmc2LXbosRkWzi5+iq0jgnrXmr3yU0oRSBEIegh4s5ZBpZHwkrZ5tboV/VMuf8iEHmo8ZIjdqHnpAUOk/8pnQpCBoVzFMSVmlSGK6mkZo8HGCNhpC+QBNGk9OeoqcMJlNKoO+dVwlK567mWl0C1jSwUF37OtGjq2i7ro1cVKH27x4AjXFQ95hb+zUBLVYZ8EDN9pw2t/plRqJk6dfWu74ony8x/UGgp5bzpsWwgijRllTu0aHakMYNv4l9/sQStsfoVDfady3k7yE1T2/b7qdbgQS6SpEx0jc+wXM3uAo5Ohl9Hvt4LgBrOqHP0U9885N6lElXO/QLRudUqCFEKSjnxtc4o6M0AZHQJfh6fdtYMvIcj6fvx1Vuc5HXekqXk5c8dZ8IrXGk23nFhpKEpNiJB3XPMEKmrIRUUXEYUVW2ccOllJjCQdt2nP4qsC6EuJFbY/xOMGnfOEsnUekxcsej3X+DLmHen/hrNEJDrtTQWWD9Iqydlcb48oM+g0nxa8gg4wWlS67AXMzTUPmRY3SS0CN05SH9pe1nILXrDZynq3mkkooOVQva+f9QclYg3/VU4mc/PZfKVFXWy5ILVo/jG4XUtJ+XQTKz24AfBh5MKb2gbbsS+E3gOhqzyctTSn9mZgb8KvBSmnvnlSmlTyx+6ItBvOGijcfbjtIQjbernRU6glKjqaaP6fE08DjG/znU9HMhjEfjBqPxdgW4up3QK3pSkfWeF3bWTjwE3Qfpekqu/IiemBKSlrr0fmKpS08AfITGY3C6fXlYu+u8j/SLS7raFtNB/NquM+kUcPuRJrnm+HVIaFA7UlTZSiSjhDg67GMJ6XbgLcDbpe11wO+nlN5oZq9rv/88jWXg2e3rRTRLa79oUQPeLcSJHuvarNEvGQudLSPeoJrQqfYIl44uZ5ID/LP3c4T5wkScK9SweuVpeIIn0EZj8QdpKka6aAWd9JLL9tcCbSVDjC4Mp4SkJWydkM4A97TtZ+gTkufbtOT0Z1udLc5NTG7n1qH44aLna7XdfoXuml6kf31LpFSas9G85mMZJfHksJ+N2imlj5jZdaH5ZuDF7ee30ZSL//m2/e0ppQR81MxOhNVsRwWVSoZupiN0c8m/H6df8RW6eKVH6dtnXUvxORnTSXKhQUMStXrgoF8eV23aT/kGPMUz5rV6/hng3wD/kLwRObeI2ZCvXcXKWC3SWdwJ72N06piO57S8A98839PaHt80Vk7wwzvZQ5/wS2kmR+geKI5crrG/D5FTVNlWBraPIQhLw9IHMImd2JCuEZK5H7im/fwMmnJAjnvatlES0rxQc4h/15vPHUKP0MwrT8b10iW64q0jF3NUMu3oBNmivxKrrj8GHQed2YITp+DKU027efDUCeC59Gtql1ZsVD97SWVT16GfaK7I/0n6Rf798yNw7lKfN/0cvEa2lvvWYeihQy7uRA6iQiuukDm1XCmnXFCrv1RFH7W0tJ8lpGlIKSUzS/PsY2a3ArcC2CIGUVFRMTsOYOrIA66KmdnTgAfb9nvpyrcDPJMu1etxtDV5TwKszElmuwUVuaO0ovGG0C8BvUn3BL9EF9V9nK6u/gm6KocqUbmUE8V4H4s/xPTeuSif1Wjrph+VMHR1Ei8ud/yR5nXiq7ByDXDuW5sfnnYG0qOTKluMRFSRLnrW/OKodKSuwOfwuE7rnrOz4ZWThnILoOh/pCWfNJg1Skcxziv+r9CvWaXX3/stVRhQG9JovGlDGOEAd0JI7wV+Anhj+/670v5aM3snjTH7kbHaj2AyXyynjVwKL+iL5yth2zM0c+8aOkJ6Ep19w7eDzoMWa6VFj18uslsnkn8/R6cJRY+/q49uyzoOHDsPR+1PAXhO+jGwM/DEb8ATzwHfkFH667H2FXEZ3e2kR30iruA+av/+cROVLmutbTlbeu5/iTYh1Rj9AQKTqlp0PkTvpm4D/SwWVwd1cYEj5O8Z/y2GG4wG+zkOyczuoDFgX2Vm99BUhnsj8C4zew3wVeDl7ebvo3H5f5HmHnvVgse8K9gkHzHtruKcWcWDmjfphAk3hzyLpmibL/RxjE4yit43dUVHI7s7rNSRpTYofYJvMJw25ob0KLw4UX7M/k+OAX8z/R80ZKTmcSUjJyQlpkhIlwHH+LD9457B31dYyq1DOVTZMRKPB4CrJKSkM1FJkz5RaY30tdCObJ/zkmokuEpdOUSv7WgKkO1nQkop3VL46SWZbRPw0zsZ1F5C3ewa+KZQg7E/1Y/SuZHP00klD9NIRU5Gria5O9rd1dGIGpNufQLEVAiVBCK26Mrket+l0CGVwDTochX4A/s7vCX9GX3KeIxhCQm626khpzfYX+zFaynxaHuUgFT60Mq20FXujDlmel1WM+cb1W4nI03nidvre6zjHb1pfi3jecQQMFcFdexL44XRiW01UvtxuC0nuuNX6Huv3IPl75fox/Kdpok5VK+a9+PE4GoI9J/AOqliXGHMffN91Ybkal/Mb9PcL982xizp+R4B/oI9uXfMT6bpZr7vtc49oV7AUgJq9FrF1Le4bLhDJSNHST1ycoqSjV/bi0zW5taYrrh9zBdG9lNpaYs+6Sv5jEJKGqmENIprU1FRscfw1JFZXjOgjTf8bTP7EzP7nJn9l9sZVpWQWpRUKH/K5VS2TRpJ6SE6CelqGskolh5x1c7tKG4sVSlBc7DUvqMSQUxfiDYPfXePmwdq6qq73r9KDWqrUlwCnm82ofrlyqUocnE46l1UB4GrUHHZ8A3613KVLscvxmpNS/tQ9cmN06vyGfLpKH5ctVtF56Krl5pCFNN6vH00mtJiJaRfBd6fUvpRM1uj86HMhUNNSCW7gbb7je9kpEGH7h3SJbOfS/NPuD3BJ53boXyt+RjUt0ZnH4FOVSmtA+Y2mJzBNkZve1UP3/ZqulAAXUnFj+XEFI/piMSXs1PlUlm0Dls0aqunLZe2syHbrks/6llTMs0FOioxuNqVU+X8uqobXwksR7T+QPGQC/1v1KHg1yM6H/Zce1pgYKSZXQF8H/BKgJSS+1jmxqEmJMi7fOPN7U8+lZD85j5LY8hWb5p7bVTi9UmoNdK8b5cKPJAZOvuJEoTDDaNqs4jGUz8XP+4WnR3G7VsesO3tmtqitpxoYykRz4a0aeiSGpL1peSuxzgr268zaQAvaRR6XUpkGY3gkZD8v4v3geYgKlnHwnDe7mQZwzWGwgT2HIsbyPU0z+V/Y2Z/CbgL+JmU0rnh3SZx6AmplAag2fi03z0lTPf19IZrZTtX5aKHySeMTgBNfj9OX1qJ0pGPTw2nQxUQo3FW1UGPQ+oFTNKXmrRaSUmd0wkfyVfDJfwaaB02Ha+Sp05qN4pr37pdVIXiGHOSGnRSi3vwtH//D2MgrEqOUVV1acr/1w35HB8cTnhLTS2Zz6h9lZnp0tgn28Bmx2U0a6T+vXYV21+lSbb/J/MO61ATUi6xNgYo+nZr9OOQfPuz9IOR/YaLEcZb0s+KbO+EUFKfcq5hJQH1YG3QlwY25RWDK9Wrp0SlxKjtJVLSsWiMlp+7kq+TkJIPcp45KQ85L/+sdaIcTiox0LFESJBX2dZl+6g6DgVYan/QEXLJGzcKKWl2le1USunGgd/vAe5JKX2s/f7bNIQ0N6qXraLiMGKBXraU0v3A18zsuW3TS4DPbmdYh1pCishFUrtk41LF+bDPBbp62Q6XFrQsCfQT4KM9J0olKo3og0ztGyplqCoT1Rs1sOZQelqrvUUjxGNfUQVRqVAljRhUWEJUw6LK5rWN4jlE+078PZeaEwMv3f7lUpNKPyvSFg31uZiqC/SlJj2/nEdyT6Wmxcch/T3gHa2H7W62maFxqAkpFtTKubBj/ZyYYqC2AWhubN1O9/WbXI3Xx+nbbqLtyidkVBGiSqilXXOJqCXjfc7wnpuMakvLIY4xV1Ilp2rG9jhHYrvblWIycokw4/8bvW5Rbb1EV3HSyQ/ZLqb3bLbbRBuWq/jxYZDzBC7NlrRAQkopfRIYUutmwqEmJCj/J/Fmieq2G4rdAOo32hU0xOJhAHoDuo1GbUVeZjrGBPkYon1G3eOaGX9O2jXeyEklpi/4RFKjc5xwfg3cPrNO33Cr41Qvm5PGhnyGvudNbS85ItL2KCFdCNv4OXkun9qR4meFPyT0+NGGdCZsv0bfzqiG92g01xACvWYrme33HAe5HtJBgN7Y+q6koHEvfsMdoUkX8e2fSpM68gX6N+ImXaCjutS1hEnOc6PxS9B5qCIhaekO9Wq5xBMlB59IF2QsfoxYlsOJyytd5iZ49LJtyH4XpV1fMaUkNz9UhYJOOooPCnX5r4bfc2pclFQiIam3zeHkpvvq+UQniXrUhoy1SyOmUVjW+6iEVFFxGLGfVx05yNDAtpwYrWqFGmQ12XKNbgGN48Crnw5v+3pTxlYjtWMqCfTrIcU0EE+UjRKSpoKohKTqHPJblGxUJbko2+lvKvG4SplLRI2pLMi4YzCoqmtxe4XGIUWVLaprKmWpRFrqM/ddE6o1ZivW7c6FJrg6nZPEcvFMS409ChihgFQJyaFklDN0ax4TdJ4Zrwr5tbb9PcClr3eF0jS3ySeWBgVeoltnTe0dvm1OZXMSyRU5020hX8NbSUQNxJFcot0mV3kgd63UWB49gZv0SSlC+3RC0DFGdU0JqWRw9/1KBm/dT8l4nT7ULqTbq8qm8WiqtkWVbdkR2yNN9j/chFSaENG+EKUh6HvMnkJXGfKrwJtpltTWJck8GNIn9vnQT5xMMf0iV1IkrnbiE1gjrHXc0WDuk0kJoGRIXpV3PadYWsPHrt5AlRJV6iHsp4mq3n8MW8gFik6DSk2xn1zUtUpIagP07XOeNn/l+tdIbcV2zmWRGKFN+3ATUg651AB1DSshqdtYa2dfQyMxPUDnpXG3vmbw+/FiDSaYdMvHSX0htPuYlAQ1GTQnPeRMCNGI78f0HK9ISPp52tiVjCIpKRnFWB5V07x/R3Q+aH/+HtXtuI326dfQiV2TYnWc0Wu5Sp/ENUwgRpWPIRq5Skj7APHpBv04oyjxeLyRSgFHafLa1mnqaHtZkodp7Ex+w2sahyJ6nmLcjk8uJzLNUleC1DHmAgZVKiG0bzGpDjmZRA9hKfu9REgaVOiIapRKJd5Pbuw5z1xOQvHrqJJNjA2KdqQoIUXij+PxtguZNpXMStd9GagSUkVFxSigcWNjQiWkgPiEjSL4euZ3NVIfk989ox46ycZtP7n4lqgKaNpBLsZF23XcOQO2jyGHqPbEkh1+3q6aaCqISy+lpFsNjFR1MCfZxGvgapB77eI5DNll4vXISYi5+CT/n6P662N2iSf+L34cNeyrNBttUTGNZRnqU5WQ9gHUMAmT3hC1IfmE8xIk3g7dJNLaPiqu6428Gt59m5wKOQ05Fczf1TulRtocAUS3fK7dbUqRHCPBqcoWU2EUfg1ypKOhAzqmeJ6OSCROnDkXfM6+pA4IKP9/Xr3B/y+/P87Jb3pcJ59lu//3rQ3JzG4Dfhh4MKX0grbtXwIvo7lPvgS8KqV0xsyuAz4HfL7d/aMppZ/ahXHvCvypGW8YFW81p8snsMb+qE1J6yHpxNUb1G0qngOXe2IPjXeWNkdMI4mxV273GIoTUoJxCWlaQqsSmJYSiVKiT9ZcMmqMji5JR0eYJHwynx1rhe3cZuiEdJbuv1Li0WPlvG/R6F8Kd1gGxjIOxSwS0u3AW4C3S9uHgNenlB4zs38BvB74+fa3L6WUbljoKPcQuRtdDb2O1fDu8UCn6aoxxjpJ0fsEk6uZRCN1abL7uOZBlLg01QIm02Zyx9NJpirJLISk554b+9C5KpzQo+dKVTPvR+tXRelxhcmVZlRaVKeBP2RcbdMHkxqut6Q95/UbC0aayjadkFJKH2klH237oHz9KPCjix3W8lCa5D6ZdNK6lKQSjwdKet1qJyS/IX0ixWTZi6E9Ll3k2I76VrJNxfiYnBve2zcy25U8eI5ISNFrGCfs0LnlJJ5ISsg2JZtanIQlL2G0Ger/ofs4afn/FGOoxqCe5TDSzJGFEPergf9bvl9vZn9kZv/BzL53Af1XVFTsAjZnfO0ldmTUNrN/RLOM6TvapvuAb00pPWxm3wX8OzN7fkrp0cy+twK3Alj8cckYkpLUY5Rj81UayehKmgjuo9IOZQnJs+JVQsqVtsjZdWK7bhuD9XwsOU+Rf47BobEciRrEp4n98aaOhnHFkCRR+s2vl1YByElNQ4h9+zWIQaBeUuXysK0HPkabmOe4RftUTjLZa0/bvjVql2Bmr6Qxdr+kXT6blNJF2moTKaW7zOxLwHOAO+P+bZHwkwArZtOXRd0DlAyPMdcJ+jYgtQlpVUh1+7uYv0mfeLyGdSy4pkGFcfXXnMqjbvaSzSZOhOgOL036aNTWkIdZ7Fu5G79EZCWj/pBdS8eir6FjRvJQ+DVR0vbk4hgdv0KfmGNfpQdXVOWq27/BtgjJzG4Cfg74/pTSeWm/GjidUto0s2cBz6YpZ7mvEO0tavjVCbBJX8KBxnY07UmvBKdPY139YqiY2Za0x/HEyGhFKXM9F+1dmkTQv5Evhd8UuRt+yEZQ+k3XZVMDvEt6jmmSURxj7n9Sgo5pLJ7o7DWTvA9/0EQbnB9z2cSTw76VkMzsDuDFNEuh3AO8gcardjnwIWvWc3f3/vcB/9TM3P77Uyml07s09l1FjBOKiwS61OKF0XSBx6ieweTTOU7uqP4oEV6i8wT5cd1lH0lHvYEasxNTH+LL22M8lPabOw99n/bEHTKy57ZT6PlroGE0ypcCNaeNS/tTqUbJKabCINuuMCkNzaLSLgv7lpBSSrdkmt9a2PbdwLt3OqhlIecR8QnqL73hLtDVJXp62/Ykuly1S/TLyTpK6hbkJ3dUCXIJpT4ZlQS1aiN0pKlRyB7/BP11x3Kqh79vDbRpu07meC7xWucIJNfm+3sfKrm6LWeWlIgYs3QptEdS8/9TI8/9uD7WSD7e7xgn/kH2slVUVOxDbM34mhVmttJ62H9vu2OqqSMB0bMEfbVNVQeXgHxJauiWM4pZ4iVjdJSU4ljUWI189id1bnHGaHBfpx+kqblaakOKEs20HLetQrsjF4g5hFmejjGaWwMRfaxDT/4YVe3SoEbi50q1eN8udcZr5L+rJKuxSKX+loVdUtl+hiZT40nb7aASUgG55M+osm3STHoPgoT+cka5tIlLoT3eFDlbU+xH01NK6hl0FS29uqXaipxgcwGJcUw+hrisdYmooJuM6s72mkqz2pqGPHFuR4seUH9gKCnlUjuiUd+317pPMHktosqmHtASfKz+eSxYJCGZ2TOB/wb458D/sN1+KiEVkCOMGAHsdhB38UO32GO0H6hkE9ujTQb6k1y3j2SkEpITkLqqteys5l/lIqx1LBqd7Z/j2C8xSZg6fu9TpcpZpSUdUw4x51ClJkUumVa/K3H6GFdDm44jOggu0pFSqf+x2pAWTI7/isbzfnzahkOohBSgNzj03evqRfGJqGkYvr/GseREeO0/uvB1Iudc+Rqv5GQIXZpKVM10AUXN2dIndo5MNDE4xkMpUZXc/ioV+PaayJqTlHyfWUMIYsBhDup5i/v6yrJKzH6ummSsx48PCC3LW1L1xoo5iPIqM9NYwpNtHCEAZubJ93eZ2Yt3MqZKSBlEtSzn7p2Wp5T7s9U9HAP+IvG4Ohg9Z17FUAkI+mRUStCNT/wciahqqYRUkpBK55qDJ7LGhFbto0Qy044xFECZK4Pr+zixa5zTBv0FI+PxVW0rLcbpyBHsGKSmOb1sp1JKQ6vSfg/w183spbTFUs3s36aU/va846petoqKQwg3as/ymtpXSq9PKT0zpXQd8ArgD7ZDRlAlpCJUNVN7TVTlcgbh3FNVDau5fKboTXMPWjRUR9VMJSRfPCDGFUVEO5F/j3WcchJSLlcuoqR2Ie1qTyoZfWeZDENR8Qr1SPp+WppEVW6XkrSEiKqTm7Kv/z8lI/iYMSYDu6MSUkB0+zt5XKBJEXECiCrctP6G/nw1Xqsq4MscKfEcJU9I0XPmY9dj+BjUFuRGWT0nJSYfe87Opec3KzEokUfizkWrI78NYdbjaylc/a81RURVUx2je079unl/pXSbEmGPAbsVqZ1S+jDw4e3uf+gJaRY7hZPFeSZtPpr/BX0bDvRJy43cJSO1LvLo34/RSD4eZe0evdw6cUOesyh9qV1MPXXTvGkweZ2mkUFpUkYy82sTCVWJKke2Q+PIufxh8mGi3j+tIR6jvlVidgyl2owZY5TmDj0hKdzYmPujcgZcjVfxkhQn6LL3z9KPE4rkA31JSBd+PNL2pQQEk6786N3TMSoRKSFF1SxWGIiepJIaVUr9GPJ4qdqTM+yrwTe255Dzns26Xa6a5BGa/2GdSU+jBrv6GOMDKAaTxnspl2KyDIw1daQSUkCJlHJxSS41qD3kqcCfB75Os4qtF/n35a/P0iceXRrbvWfQEY9KQ5BP9IX+jR4lm5iDpd47bVcbUUndiBO7JLEMkdIQ1J5UstGpNDONIGNtpxg97qq49rFC84DxQv9+XH/3mC/ol5Up2ddy6tyysW+TaysqKg4mxiCpRVRCCsgZJvVzfKps0JUggUbCuel6eP+Xm+W0/Qn7KM0CAGfoaihBc1McpUk/ycUVxRIhMZ7IoXaeaIx2dbEU6BhTQXzfHHKqWUlqcpT6mmdC5IznQ6pcRE7S0gL+uTilLbr/VQMnoS+x5lTaXMDpmFAlpJEid6PHkq5bme1ULXqIbt2nPwAufLlp+yqN6ka7zUM0hHSMLhnX1bJSvlnJgwPD6R0+tpKROpf2EckoHjNHQENGZu9zHoIaChnIbTsrct5T/S16387RLH3+SNvuxKOOjDienOF/jOoa7ONVRw4rchNMJ4rmO63QSEPQEM6naAjmXPsdmht7ncbGpKuRaPyQ3uhDGedqrI4u+1KbGrVLROSInrqYf1Yi74iY3xVJYZYJkSMnNX7PG7eUS5r2BR2huUa6ErEuY6XR3DGaPzoCfJ8xkpGjSkj7FFEqcjgpaXLtCo1UFLPRn07fSJ3LN4O8NBGfvKXgRVfLVBLaKPSBfC4ZpKN6VpKYclDJSI27evztJNjGGCY//1ljlrQvvy5xHLmHgV9vL9mifZSqHQwFpS4b1cs2cuiNPnRz6yEtkCwAAApSSURBVARQT41m2B9nMiZIiSfahPz4+h69XSVJKNqF5pGE4vkpEUUpqJRjpsh5vGYhJsUsE3aWbUpkNy0mqtTuqlqswxT/hxyZz3qsvUS1IVVUVIwKlZD2MWIQ3xad7WGDzlbkUokbqHU/X9c+F40cA+ZKtqKYka9P8JzqUMo102PnbESzeNOmYeiGn1U6KsWDDe0zD4YkqWgz9Kh4VX/9f5hmLxqLquaoRu19gGifcORuNi8m7/Ygv1nP00X6RtUsVwZE+42EpKRUIp5SOsc0lCbPrPaXWX+P13LaNtFeNG0bxXYIc+ic9H/3/1eJCCbJaFoq0pikkjGNxTHLMki30SwI+WBK6QVt2y8Cf5fGkw3wCyml97W/vR54Dc35/v2U0gd2Ydy7hmjkLdmW1MsWlxKKhmfFkI1haPKWXO05r1xuIpcki1IkdIyYLmEWqaa0bUn6mUY+JaN77vPQ2FaZTki5/dVOF22FOYxREtnPEtLtwFuAt4f2N6eUflkbzOx5NPVQnk/jWPp/zOw5KaUxknEROthYKjU+Mb0kiXrNjtBfTHA3MZR5X5rwJUkwti/aMDytn1kknxwZx9/iPtp3lEg1rWQIUWItOQDi8bYrwe42ErMtF7XXmGVdto+Y2XUz9ncz8M52Se0vm9kXgRcCf7jtES4ZcdJojpQWlPftNLhRJ82snqoccsQSb/RZ7DVKVLrfvLE8jmlP2Hkn4bx5ctspTxIfNkOklAsR8L5KY5tFOh0Lxigh7SRu67Vm9sdmdpuZPbltewZdjCDAPW1bRUXFiLDIipGLxHaN2r8G/BLNef0S8Cbg1fN0YGa3ArcC2DYHsZvQp+FQXpd72tym4Fn7XigtxhvNqmZMU3/mKWMxzaZSMiTvtO/cb0NPwNz20/LkZpU0Sza0aXasXLDjLHFYY5eODlQcUkrpAf9sZr8O+EqV9wLXyqbPbNtyfZwETgKsmKXtjGMvULK3KKKNY4uGmNzgjXyOnrZYmXJoHI6h/LN4k7k6EtXFnE1lCPN6sBbR1yxjm2VSDamk0wzpQwGi8/QzRoxRZdsWIZnZ01JK97VffwT4dPv5vcBvmNmv0Bi1nw38px2PcsnIPfmiZKMS0CYdMcU6zjlCKnneZoloniUuJ2cj2a49ayfYzbyuWQzt84QRzBJxrfvOaqAfC/Zt6oiZ3QG8mGZtpnuANwAvNrMbaM7rK8BPAqSUPmNm7wI+CzwG/PR+87ANQSWNaU8XJ53o1Ym1uIe8NdshikXF6IwlKXQ7hvMhkp52M5ZixIb2H7t6lsO+VdlSSrdkmt86sP0/p1lO90CidEPOGhNTWl1DY1ty++5UtRqSkHIYizg/TyjBPDFRinltXftNGiphUeM2s2tpwoKuoeG6kymlX91OXzVSu6LiEGLBgZGPAT+bUvqEmR0H7jKzD6WUPjtvR5WQdoChgL9pEk2uCFzse8j7tRO7z6zBgMvAtHHNG5RZwk4lxLFev3mwqHNo7cn3tZ/PmtnnaMJ9KiEtC/P8ubOoX0MG7XmPN3Zst1DbvNgJiR+k6w27Z0Nqg6i/A/jYdvavhLRLGPLczCMqzyvN7IeJs4wxzmpbOyyY08t2lZndKd9PtmE7PZjZE4F3A/8gpfTodsZVCWkPMZabfyzj2EscxnOehjkejKdSSjcObWBmqzRk9I6U0u9sd0yVkCoqDiEWqbKZmdF43j+XUvqVnfQ1lnCTioqKPcYCc9m+B/gx4K+a2Sfb10u3M6YqIVVUHEIs0u2fUvqPLCgltRJSRcUhxRjtapWQKioOITyNaWyohFRRcUhRJaSKiopRYD/X1K6oqDiAqBJSRUXFKLBvy49UVFQcPOzbAm0VFRUHE1VCqqioGAWqUbuiomJUqBJSRUXFKFAlpIqKilGhSkgVFRWjQPWyVVRUjAZjjUOaWg/JzG4zswfN7NPS9ptS9+QrZvbJtv06M/um/Pa/7+bgKyoqtgcnpAXVQ1oYZpGQbgfeQrPuEgAppf/WP5vZm4BHZPsvpZRuWNQAKyoqdgf70qidUvpIu5LABNrSlS8H/upih1VRUbGb2Lcq2xR8L/BASukL0na9mf2Rmf0HM/ve0o5mdquZ3Wlmd6YdDqKiomJ+bM342kvs1Kh9C3CHfL8P+NaU0sNm9l3AvzOz5+eWRGmXUTkJsGJWOamiYg+RgI1lDyKDbROSmV0G/E3gu7wtpXQRuNh+vsvMvgQ8B7gz20lFRcVSMNbAyJ2obH8N+JOU0j3eYGZXm9lK+/lZwLOBu3c2xIqKit3AIr1sZnaTmX3ezL5oZq/b7phmcfvfAfwh8Fwzu8fMXtP+9Ar66hrA9wF/3IYB/DbwUyml09sdXEVFxe5gkW7/Vgj5X4EfAp4H3GJmz9vOuGbxst1SaH9lpu3dNKtXVlRUjBwLVNleCHwxpXQ3gJm9E7gZ+Oy8HdVI7YqKQ4gFp448A/iafL8HeNF2OhoFIW3BqfNwDji17LHsIa6inu9BxhjP98/5hy34wLlmjLNg3czUMXWy9ZIvHKMgpJTS1WZ2Z0rpxmWPZa9Qz/dgY+znm1K6aYHd3QtcK9+f2bbNjZ0GRlZUVFR8HHi2mV1vZms0Dq/3bqejUUhIFRUV+xcppcfM7LXAB4AV4LaU0me209eYCGlXdNIRo57vwcahOt+U0vuA9+20H0upZm1UVFSMA9WGVFFRMRosnZAWFXI+ZrRF7D7VFq27s2270sw+ZGZfaN+fvOxx7gSFQn7Zc7QG/0v7n/+xmX3n8ka+PRTO9xfN7F4pUPhS+e317fl+3sz+6+WMevxYKiEtMuR8H+CvpJRuEFfw64DfTyk9G/j99vt+xu1AdCWXzvGHaPIcnw3cCvzaHo1xkbidyfMFeHP7P9/Q2lVo7+lXAM9v9/nfPOezoo9lS0iPh5ynlDYADzk/DLgZeFv7+W3A31jiWHaMlNJHgJi3WDrHm4G3pwYfBU6Y2dP2ZqSLQeF8S7gZeGdK6WJK6cvAF2nu/YqAZRNSLuT8GUsay24iAR80s7vM7Na27ZqU0n3t5/uBa5YztF1F6RwP8v/+2lYNvU3U8IN8vgvFsgnpsOAvp5S+k0ZV+Wkz+z79MTWuzgPt7jwM50ijen4bcANNscI3LXc4+w/LJqSFhZyPGSmle9v3B4H30IjrD7ia0r4/uLwR7hpK53gg//eU0gMppc2U0hbw63Rq2YE8393AsglpYSHnY4WZHTOz4/4Z+EHg0zTn+RPtZj8B/O5yRrirKJ3je4Efb71t3w08IqrdvkWwg/0Izf8Mzfm+wswuN7PraYz5/2mvx7cfsNRI7UWGnI8Y1wDvaRZo4TLgN1JK7zezjwPvagvefZVm9ZZ9i7aQ34uBq8zsHuANwBvJn+P7gJfSGHfPA6/a8wHvEIXzfbGZ3UCjmn4F+EmAlNJnzOxdNPWBHgN+OqU0xkU/lo4aqV1RUTEaLFtlq6ioqHgclZAqKipGg0pIFRUVo0ElpIqKitGgElJFRcVoUAmpoqJiNKiEVFFRMRpUQqqoqBgN/n95dAC1WaldRAAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "#NBVAL_SKIP\n", "\n", "# Plot percentage error\n", "plot_image(100*np.abs(model0.vp.data-get_true_model().vp.data)/get_true_model().vp.data, vmax=15, cmap=\"hot\")" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAZkAAAEaCAYAAADDgSq4AAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADl0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uIDMuMC4yLCBodHRwOi8vbWF0cGxvdGxpYi5vcmcvOIA7rQAAIABJREFUeJzt3Xl4HfV97/H3V5sleZH3VYttYbDN4gUvYjVhSTAJW8JuJw2hYbkhTW9625KUW27btLftvc1D25AAAR6SyDYQgiHQGzAJiwkgr9iAMZsXLd5kvMqSrfV7/5gRPpZl+1g6Ix0dfV7Po8fnzMyZ+Z4jWR/9fjO/35i7IyIiEoW07i5ARERSl0JGREQio5AREZHIKGRERCQyChkREYmMQkZERCKjkBERkcgoZEQAM7vFzFaa2QEz22ZmvzOz87u7LpGeTiEjvZ6ZfR+4H/gnYARQCPwUuLo762plZhndXYNIRylkpFczszzg74HvuPsz7l7r7o3u/ry7/6WZ9TGz+81sa/h1v5n1CV97kZlVmdlfmFl12AK6NVw328y2m1l6zLGuNbN3w8dpZnaPmW0ws11m9pSZDQ7XjTUzN7PbzKwCeCVc/g0zKw+3/59mttnMLj2J/f2JmVWY2Wdm9jcxdaWb2Q/D19aY2SozKwjXTTSzl81st5l9ZGY3dMG3RVKIQkZ6u3OAbGDxMdb/DVACTAWmALOAe2PWjwTygDHAbcADZjbI3ZcBtcDFMdveAiwMH38XuAaYA4wG9gAPtDn2HGAS8CUzm0zQupoHjIo5Zqt49nc+cBpwCfC3ZjYpXP594GbgCmAA8C2gzsz6Ai+HNQ8HbgJ+GtYiEh9315e+eu0XwS/t7cdZvwG4Iub5l4DN4eOLgINARsz6aqAkfPwj4LHwcX+C0CkKn68HLol53SigEcgAxgIOjI9Z/7fAopjnuUADcOlJ7C8/Zv1y4Kbw8UfA1e289xuBN9osewi4r7u/b/rqOV/q65Xebhcw1Mwy3L2pnfWjgfKY5+Xhss9f3+Z1dUC/8PFC4C0zuwv4KrDa3Vv3VQQsNrOWmNc2E5wTalXZpo7Pn7t7nZntilkfz/62H6POAoIwbasImG1me2OWZQC/amdbkXapu0x6u7eBeoKupvZsJfhl26owXHZC7v4BQSjN5ciuMggCY667D4z5ynb3LbG7iHm8DchvfWJmOcCQk9zfsVQCxcdY/nqbffZz97vi2KcIoJCRXs7d9xF0RT1gZteYWa6ZZZrZXDP7V2ARcK+ZDTOzoeG2pSdxiIXA94ALgV/HLH8Q+EczKwII93+8q9meBq40s3PNLAv4X4B1Yn+xHgH+wcwmWOAsMxsCvACcamZfDz+TTDObGXMuR+SEFDLS67n7vxGc/L4X2EnwF/zdwLME51VWAu8C7wGrw2XxWkRwMv4Vd/8sZvm/A78FlphZDVAGzD5OjesITu4/QdCqOUBw/qe+I/tr48fAU8ASYD/wKJDj7jXAFwlO+G8l6G77F6BPnPsVwdx10zKRnsbM+gF7gQnuvqm76xE5FrVkRHoIM7sy7M7rC/xfgpbV5u6tSuT4FDIiPcfVBN1WW4EJBJcgqytCkpq6y0REJDJqyYiISGQUMiIiEpleP+J/6NChPnbs2O4uQ0SkR1m1atVn7j7sRNv1+pAZO3YsK1eu7O4yRER6FDMrP/FW6i4TEZEIKWRERCQyChkREYmMQkZERCKjkBERkcgoZEREJDIKGRERiYxCRkREIqOQERGRyChkREQkMgoZERGJjEJGREQio5AREZHIKGRERCQyChkREYmMQkZERCKjkBERkcgoZEREJDIKGRERiYxCRkREIqOQERGRyChkREQkMgoZERGJjEJGREQio5AREZHIKGRERCQyChkREYmMQkZERCKjkBERkcgoZEREJDIKGRERiYxCRkREIqOQERGRyChkREQkMgoZERGJjEJGREQio5AREZHIKGRERCQyChkREYmMQkZERCKjkBERkcgoZEREJDIKGRERiYxCRkREIqOQERGRyChkREQkMgoZERGJjEJGREQio5AREZHIKGRERCQyChkREYmMQkZERCKjkBERkcgoZEREJDIKGRERiYxCRkREIpPR3QVEwcz6Aj8FGoDX3H1BN5ckItIrRd6SMbN0M3vHzF7oxD4eM7NqM3u/nXWXm9lHZvapmd0TLv4q8LS7fxu4qqPHFRGRzumK7rLvAevbW2Fmw82sf5tlp7Sz6ePA5e28Ph14AJgLTAZuNrPJQD5QGW7W3OHKRUSkUyINGTPLB74MPHKMTeYAz5pZn3D7bwP/2XYjd18K7G7n9bOAT919o7s3AE8AVwNVBEEDx3iPZnalmT28b9++k3hHIiJyMqJuydwP/BXQ0t5Kd/818BLwpJnNA74FXH8S+x/D4RYLBOEyBngG+JqZ/Qx4/hjHft7db8/LyzuJw4mIyMmI7MS/mX0FqHb3VWZ20bG2c/d/NbMngJ8Bxe5+oLPHdvda4NbO7kdERDrnuC2Z8KT91zq47/OAq8xsM0E31sVmVtrOMS4AzgAWA/ed5DG2AAUxz/PDZSIikgSOGzLu3gz8sCM7dvcfuHu+u48FbgJecff5sduY2TTgYYLzKLcCQ8zsRydxmBXABDMbZ2ZZ4XF+25F6RUQk8eI5J7PEzP7czEaZ2YDWrwQdPxe4wd03uHsL8A2gvO1GZrYIeBs4zcyqzOw2AHdvAu4mOK+zHnjK3dclqDYREekkc/fjb2BW2c5id/fCaErqWjNmzPCVK1d2dxkiIj2Kma1y9xkn2u6EJ/7dveBE24iIiLTnhCFjZhnA7cCF4aLXgEfCrioREZFjiucS5geAvsBj4fP5wHSC4BERETmmeEKmxN2nxDxfYmZroypIRERSRzxXl7WY2djWJ+Hjdkfwi4iIxIqnJfPXwBtm9hFgwCnAbZFWJSIiKeG4IWNmacA+4FRgUrh4vbsfjLowERHp+Y4bMu7eYmYPuftUYHUX1SQiIikinnMyr5rZ1ZFXIiIiKSeeczLfBL5nZvXAQYLzMu7ug6MsTEREer4TnZMxYAqa2VhERDrgRLMwO/D/3L257VcX1SciIj1YPOdk1oRT8ouIiJyUeM7JTANWmNkGoJbD52SmR1qZiIj0ePGEzFWRVyEiIinphN1l7r4BGAacFz7eCzRGXZiIiPR88Uz1fy9wHlAM/BLIBhYC50dbmoiI9HTxnPi/DriC4HwM7r4FSNTtl0VEJIXFEzL14aXMDmBmudGWJCIiqSKekHnGzB4A8szsVmAJh29gJiIickwnPCfj7v9iZnOBBoLR///o7r+LvDIREenx4rmEmTBUFCwiInJS4ukuExER6RCFjIiIRCaukDGzLDM7JepiREQktZwwZMzsy8B7wMvh86lmtjjqwkREpOeLpyXz98BsgulkcPc1gFo1IiJyQvGETKO7722zzKMoRkREUks8lzCvN7MbgDQzGwf8GVAWbVkiIpIK4mnJ3A2cDbQAzwD1wJ9HWZSIiKSGeFoyxe7+18BfR12MiIiklnhaMg+Y2ftmdp+ZTYy8IhERSRnx3LTsAuBLQA3wCzN7x8zuibwyERHp8eIajOnuW9z9x8A3CcbM/EOURYmISGqIZzDmBDO718zWAD8HVgCFkVcmIiI9Xjwn/hcCTwBXuXtFxPV0uQ07D3DtT9/s0GsH5mTy1en5fOn0kWRlaBo4EZG24rmfzMyuKKS7pJnRr09cdzw4yifVB/juoncY2i+LG2YUcPOsQgoG68ahIiKtLLizcjsrzBa5+81m9g5HjvA3wN19elcUGLUZM2b4ypUrO/TalhZn6Sc7KS2r4JUPd+DAF04bzvySQuacOpz0NEtssSIiScLMVrn7jBNud5yQyXf3KjMrbm+9u2/oZI1JoTMhE2vL3oM8ubyCRSsq2VlTz5iBOdwyu5AbZhQwrH+fBFQqIpI8Oh0yMTv6J3f/4YmW9VSJCplWjc0tvPzBDkrLynlrwy4y040vnT6S+SVFzB43GDO1bkSk50tkyKxu2zVmZmvdfUona0wKiQ6ZWBt2HmDhsgqeXlXFvoONnDK8H/NmF/LV6fnk5WRGckwRka6QiO6yO4A7gVOBj2JW9QdWu/uNiSi0u0UZMq0ONTbz/NqtlC6rYG3lXnIy07lqymjmlRRyVv7ASI8tIhKFRITMIGAI8L+B2BH+Ne5enZAqk0BXhEys97fsY8Gycp59ZysHG5s5Kz+P+bOLuHLKaHKy0rusDhGRzkhYd1nMDgcD2a3P3X1rx8tLHl0dMq32H2pk8eotlJaV80n1AQZkZ/C1s/OZN7uQU4b37/J6RERORiLPyVwB3A/kA7uA0cAn7p4Sk2V2V8i0cndWbN5DaVk5v3t/G43NTsn4wcwvKeKLkzXIU0SSU7whE88oxH8CzgOWuPs0M7sMuKGzBUrAzJg1bjCzxg3mswOTeWplJQuXVXD3wncY2q8PN87M5+ZZheQP0iBPEel54mnJrHT3GWa2Fpjq7q6ry6LV0uK8/slOFpSV88qHwemvYJBnEReeOkyDPEWk2yWyJbPPzPoBfwR+aWbVwMHOFijHlpZmfOG04XzhtOFs2XuQRcsqeGJFJX94fAX5g3K4eVYhN84sYGg/DfIUkeQWT0umP0GoGPANIA/4lbvvjL686CVjS6Y9jc0tLFkXDPJ8e2MwyPPyM0Yxf3YhszTIU0S6WMKvLktVPSVkYn1a3TrIs5L9h5qY0DrI8+x8BmRrkKeIRC8R42T20M7EmK3/uvvgRBTa3XpiyLQ62NDM8+9uZUFZOWur9pGTmc7VU0czv6SIM8bkdXd5IpLCEhEyxx0Z6O7NHawtqfTkkIn1XlUwyPO5NcEgzyn5ecwrKeLKszTIU0QSL6HdZWZWApzq7r8MB2X2S5UbmKVKyLTad7CRxaurKF1WwafhIM/rzi5gXkkhxcP6dXd5IpIiEjkY816CcTLF7n6qmY0BnnT38xNTavdKtZBp5e4s27SbBcsqeDEc5HnO+CHBIM/TR5CZrkGeItJxibyE+TpgGrAawN23mNmATtYnETMzSsYPoWT8EHbWHB7k+Z2FqxnWvw83zSzgplmFjBmY092likgKiydk6sMBmA5gZhp63sMM69+H73zhFO6cU8zSj3dSWlbOT179lAde/ZSLJw5nXkkRcyYMI02DPEUkweIJmWfM7AEgz8xuBW4DHou2LIlCeprxhYnD+cLE4VTtqWPR8gqeXFHJ79dXUzA4h1tmFXH9jHwN8hSRhIn3xP9c4IsEly+/5O6/i7qwrpKq52Ti1dDUwpIPtlNaVk7Zxt1kphtzzxjF/JIiZo4dpEGeItKuhJz4Dy9jftHdL0tkccmkt4dMrE+raygtq+A3q6uoOdTEqSP6MW92EddOH6NBniJyhEReXfYKcI27709UcclEIXO0gw2td/Is592qfeRmBYM8583WIE8RCSQyZBYDU4ElQG3rcnf/fmeLTAYKmeN7t2ovC8oqeG7tFg41tjC1YCDzZhdy5ZTRZGdqkKdIb5XIkLmtveXu/mgHa0sqCpn47DvYyDOrqygtK2fDzlrycjK5LryT53gN8hTpdTRBZpwUMifH3SnbuJvSZeW89P52mlqcc4uDQZ6XTdYgT5HeIpGDMUU+Z2acUzyEc4qHUF1ziF+vrGLhsgr+24LVDI8Z5DlagzxFBLVk1JJJgOYW57WPqiktK+e1j3diwMUTRzC/pJALNchTJCUlvCVjZn3cvb5zZUkqSk8zLpk0gksmjaBydzDI86mVlfx+/Q4KB+dyy+xCrj87nyEa5CnS68Rz4n8W8CiQ5+6FZjYF+FN3/25XFBg1tWSi0dDUwovrtrOgrJxlm3aTlZ7G3DNHMr+kiBlFGuQp0tMl8uqyMuBG4Fl3nxYue9/dz0hIpd1MIRO9T3bUsGBZBb9ZVUVNfROnjejP/JJCrpk2hv4a5CnSIyUyZJa7+ywzeycmZNa6+5QE1dqtFDJdp66hKRjkWVbBe1taB3mOYX5JIaeP1iBPkZ4kkedkKsMuMw+nmfku8HFnC5TeJzcrgxtnFnLjzELWVu6ltKycxe9UsWh5BdMKBzJ/dhFfPmuUBnmKpJB4WjLDgf8ALg0X/R64290/i7i2LqGWTPfaV9fIb1ZXUbqsnI07axmYm8l10/OZV1LEuKF9u7s8ETkGDcaMk0ImObg7b2/cxYKyCl5aFwzyPP+UocybXcilGuQpknQS1l1mZj8Hjkoid7+9g7WJHMXMOLd4KOcWD6V6/yGeWlnJouWV3NU6yHNWITfPKmBUngZ5ivQk8XSX3RjzNBu4FqjUJcwSteYW59UPqyldVs7rH+8kzYxLwjt5XnDKUA3yFOlGkXWXmVka8Ed3P7ejxSUThUzPULm7joXLK3hqRSW7ahsoGpLLLbMKuX5GAYP7ZnV3eSK9TpQhUwwscffijhaXTBQyPUt9UzMvvr+dBWUVLN8cDPK8IhzkebYGeYp0mUSek9nD4XMyacBu4J7OlSfSMX0ygrE1V08dw8c7alhQVs4zq7fw7JqtTBzZn3klRVwzdbQGeYokiRPdftmAAmBLuKjFU+xyNLVker7a+iZ+u3YrpWXlrNu6n75Z6Vw9bQzzZxcxefSA7i5PJCUlcsR/ykwh0x6FTOpwd9ZW7aO0rJzn126lvqmF6YUDmV9SxBVnapCnSCIlMmRKgX9z93cSVVwyUcikpr11DTy9KrjXzcbPgkGe15+dz7zZRYzVIE+RTut0yJhZhrs3mdk64DRgA1ALGODuPj2RBXcXhUxqc3fe3rCL0mXlLFm3g6YW54IJQ5k3u4hLJw0nQ4M8RTokESGz2t2nh1eTHcXdN3SyxqSgkOk9qvcf4okVlSxaXsG2fYcYMaAPN80s5OZZhYzMy+7u8kR6lESEzOezLqcyhUzv09Tcwqsf7aS0rJylnwSDPC+dNJx5s4s4X4M8ReKSiEuYh5nZ94+10t1/3KHKRLpZRnoal00ewWWTR1CxKxzkubKSl9btYOyQ1jt5FjBIgzxFOu14LZltwM8IzsEcxd3/LsK6uoxaMgKHB3mWlpWzYvMesjLS+MqZo5hXUsj0Qg3yFGkrYedkEl5ZklHISFsfbt/PwmUVPLN6Cwfqm5g4sj/zS4q4ZtoY+vWJ5xZMIqlP52TipJCRY6mtb+K5NcEgzw+2BYM8r5k2hvklRUwapUGe0rslImQGu/vuhFeWZBQyciLuzprKvZSWVfDCu8Egz7OLBjG/pJC5Z2iQp/ROumlZnBQycjJaB3kuWFbBps9qGZSbyfUzCrhlVqEGeUqvopCJk0JGOqKlxXlrwy4WLCtnyQc7aA4Hec4vKeKSiRrkKalPIRMnhYx01o79h3hieTDIc/v+Q4wckM1Nswq4eVYhIwZokKekJoVMnBQykihNzS288mE1pcsqWPrxTtLTjMsmjWB+SRHnFg/RIE9JKQm7n4yIxCcjPY0vnj6SL54+kvJdtSxcFgzyfHHddsYN7cstswq57ux8DfKUXkUtGbVkJEKHGg8P8lxZHg7yPGsU80uKmFYwUIM8pcdSd1mcFDLSVdZv28+CZeUsXr2F2oZmJo0awPySQq6ZOoa+GuQpPYxCJk4KGelqB+qbeG7NFkrLKli/bT/9+mRw7bQxzCspZOJIDfKUnkEhEyeFjHQXd2d1xV4WlJXzwnvbaGhqYUbRIOaXFDH3zJH0ydAgT0leCpk4KWQkGeypbR3kWc7mXXUM7pvFtdPGcProARQMzqVwcC7D+vXRFWqSNBQycVLISDJpaXHe3PAZC8oqeHl9MMizVVZGGvmDcigYlEvB4BwKB+eGj4N/83Izu7Fy6W10CbNID5SWZlwwYRgXTBjGocZmtuw9SMXuOqp211G55yCVu+uo2F3HOxV72H+o6YjX9s/OiAmenCB8wuf5g3I0x5p0C4WMSJLKzkyneFg/iof1a3f9voONVO6uo2pPEDyVuw9SuaeOT6prePWjauqbWo7Yfnj/Pp93vRUMyiE/DKDCIbmMHJBNurriJAIKGZEeKi8nk7wxeZwxJu+odS0tzs4D9VTurqNyTxBAQRDVsXzTbp5bc5CYnjgy043RA3OObAV93hWXw+C+WRrTIx2ikBFJQWlpxogB2YwYkM2MsYOPWt/Q1MK2fQePaAFVhiH00rr97K5tOGL7vlnpFAzOJb81hAaFLaLBwfPcLP0qkfbpJ0OkF8rKSKNoSF+KhrR/e4ID9U1UtWkBVYVB9NaGz6hraD5i+yF9s8iP6YprbQkVDs5l1MBsMjUrda+lkBGRo/Trk8HEkQPaHRzq7uyqbQi74g4ecV5obeVefvfeNppi+uLSDEbl5bTbAioYlMuw/n3UFZfCFDIiclLMjKH9+jC0Xx+mFQ46an1Tcwvb9x8Kr4qL6Yrbc5DXPt7Jzpr6I7bPzkwLuuEGtTkXFJ4bGpCtS7N7MoWMiCRURnoQGvmDcqH46PWHGps/74qr3FNHxa7DFyes3LyHmvojL83Oy8kMWz9Byyc/7JIrHJzLmEE5mhkhySlkRKRLZWemc8rw/pwyvP9R69w9vDQ7tgVUR8Xug3y4rYbff1BNQ/PhS7PNYET/7M8D6PDYoKAVNEKXZnc7hYyIJA0zY2BuFgNzszgzv/1Ls6tr6o9qAVXuqePtjbtYvGYLsZOYZKWnMWZQTjBTQszFCK2hNDA3U+eDIqaQEZEeIy3NGJmXzci8bGa2c2l2fVMzW/ceimkBHT4v9P5729hT13jE9gNzM7n+7HxuO388I/N0q+woKGREJGX0yUhn3NC+jBva/qXZNYeO7Ip7p3Ivj725mcff2szVU8dw55zx7XbjSccpZESk1+ifncnk0ZlMHn340uzK3XU88sZGnlxZydOrqrh00gjuumg8Zxcd3VKSk6dZmDULs4gAuw7U84u3y/nl25vZW9fIjKJB3DmnmIsnDtctFtqhqf7jpJARkVh1DU08uaKSR97YxJa9Bzl1RD9uv7CYq6aMJitDMxe0UsjESSEjIu1pbG7hhXe38tDrG/lwew2j8rK57fxx3DSrkH59dKZBIRMnhYyIHI+789rHO3nwtQ0s27SbAdkZfOOcsXzzvLEM7denu8vrNgqZOClkRCRe71Ts4aHXN/LSB9vJSk/jurPzuf3C8cecaDSVKWTipJARkZO1YecBfr50I8+s3kJTSwtzzxzFXXOK2723T6pSyMRJISMiHbVj/yEee3MTC8sqqKlv4vxThnLHnPGcf8rQlJ9JQCETJ4WMiHTW/kONLFxWwWN/3ER1TT1njBnAHRcWM/eMkWSk6L10FDJxUsiISKLUNzWzePUWHl66kY2f1VI4OJdvXzCO62cUkJ2ZWrNFK2TipJARkURraXGWfLCDB1/fwJrKvQzpm8U3zx3L188pYmBuVneXlxAKmTgpZEQkKu7Osk27efD1Dbz20U5ys9K5aWYhf3rBOEYPzOnu8jpFIRMnhYyIdIX12/bz8NKN/HbtVgy4aupo7pxTzKkjeuaEnAqZOClkRKQrVe2p45E3NvHkikoONjZzycTh3DGnmJljB/WoK9IUMnFSyIhId9hT28Av3y7nF29vZndtA9MLB3LnnGIunTSiR0zIqZCJk0JGRLrTwYZmfr2qkoeXbqRqz0GKh/XljguLuXraaPpkJO8VaQqZOClkRCQZNDW38F/vbePB1zeyftt+Rgzow23nj+PmWYX0z87s7vKOopCJk0JGRJKJu7P0k8948LUNvL1xF/2zM5hfUsSt541leP/kuUW0QiZOChkRSVZrK/fy0NIN/O797WSmp/G16cGEnMe6vXRXUsjESSEjIslu02e1/PyNjTy9qorG5hYuP30kd84pZkrBwG6rSSETJ4WMiPQU1TWHePzNzfyqrJyaQ02cM34Id15UzIUTun5CToVMnBQyItLT1BxqZNHyCh794yZ27K9n0qgB3DlnPF8+c1SXTcipkImTQkZEeqr6pmaeW7OVh17fwIadteQPyuHbF4znhhkF5GRFe/mzQiZOChkR6elaWpzfrw8m5FxdsZdBuZn8yblj+ZNzxjKobzQTcipk4qSQEZFUsmLzbh58bQN/+LCanMx0bpxZwJ9eMI78QbkJPY5CJk4KGRFJRR9tr+HhpRt5bs0WHLjyrFHcMaeYSaMGJGT/Cpk4KWREJJVt3XuQR/+4iUXLK6hraOai04Zx55xiZo8b3Kkr0hQycVLIiEhvsLeugV+9Xc7jb21mV20DUwoG8oO5EykZP6RD+4s3ZFLz5tMiInKEgblZfPeSCbx5z8X8wzVnsKe2gdr6psiPmxH5EUREJGlkZ6bz9ZIibp5ZQFoXDOBUyIiI9EJdNWhT3WUiIhIZhYyIiERGISMiIpFRyIiISGQUMiIiEhmFjIiIREYhIyIiken108qY2U6gvIMvHwp8lsByRGLlAfu6u4gU19s/4868/yJ3H3aijXp9yHSGma2MZ+4ekY4ws4fd/fburiOV9fbPuCvev7rLRJLX891dQC/Q2z/jyN+/WjKdoJaMiMjxqSXTOQ93dwEiIslMLRkREYmMWjIiIhIZhYyIiERG95MRSTFm1hf4KdAAvObuC7q5pJTTmz/jk33vaskkkJn1NbNfmNnPzWxed9cj3cfMCszsVTP7wMzWmdn3OrGvx8ys2szeb2fd5Wb2kZl9amb3hIu/Cjzt7t8GrurocZOdmWWb2XIzWxt+xn/XiX31yM/YzNLN7B0ze6ET+4j0vStkTuBY34Bk/sGTpNAE/IW7TwZKgO+Y2eTYDcxsuJn1b7PslHb29ThweduFZpYOPADMBSYDN4fHyAcqw82aO/k+klk9cLG7TwGmApebWUnsBr3gM/4esL69Fcny3hUyJ/Y4bb4BPeAHT7qZu29z99Xh4xqCXwRj2mw2B3jWzPoAmNm3gf9sZ19Lgd3tHGYW8Km7b3T3BuAJ4GqgiuBnEVL4/7gHDoRPM8OvtpfLpuxnbGb5wJeBR46xSVK895T9AUyUY3wDkvYHT5KPmY0FpgHLYpe7+6+Bl4Anw+7VbwHXn8Sux3D4jxoIfv7GAM8AXzOzn5HiI9rD7qI1QDXwsrv3ps/4fuCvgJb2VibLe9eJ/45p78OfDfwFripfAAAE3ElEQVQH8BMz+zIp/p9b4mNm/YDfAH/u7vvbrnf3fzWzJ4CfAcUxf5l3mLvXArd2dj89gbs3A1PNbCCw2MzOcPf322yTcp+xmX0FqHb3VWZ20bG2S4b3rr+2E8jda939Vne/qzddbSLtM7NMgoBZ4O7PHGObC4AzgMXAfSd5iC1AQczz/HBZr+Pue4FXaf/cQip+xucBV5nZZoKelIvNrLTtRsnw3hUyHZOsP3iSJMzMgEeB9e7+42NsM41gaqKrCf4yHGJmPzqJw6wAJpjZODPLAm4Cftu5ynsOMxsWtmAwsxzgMuDDNtuk5Gfs7j9w93x3H0tQ0yvuPj92m2R57wqZjknKHzxJKucBXyf4C3NN+HVFm21ygRvcfYO7twDfoJ17G5nZIuBt4DQzqzKz2wDcvQm4m6DffT3wlLuvi+4tJZ1RwKtm9i7B/8mX3b3tpby9+TNOiveuuctOIPwGXERwg7IdwH3u/mj4C+N+IB14zN3/sfuqFBFJTgoZERGJjLrLREQkMgoZERGJjEJGREQio5AREZHIKGRERCQyChkREYmMQkZ6JTM7EP471sxuSfC+f9jm+VuJ3H+imdk3zewn3V2HpCaFjPR2Y4GTChkzO9HEskeEjLufe5I19SjhrS9E2qWQkd7un4ELwmlf/ns4dfz/MbMVZvaumd0BYGYXmdkbZvZb4INw2bNmtsqCuzLeHi77ZyAn3N+CcFlrq8nCfb9vZu+Z2Y0x+37NzJ42sw/NbEE499kRwm3+xYK7QX4cTn54VEvEzF5onZnXzA6Ex1xnZr83s1nhfjaaWeyN9QrC5Z+Y2X0x+5ofHm+NmT3UGijhfv/NzNYC5yTqmyGpR1P9S293D/A/3P0rAGFY7HP3mRbc7OlNM1sSbjsdOMPdN4XPv+Xuu8PJGVeY2W/c/R4zu9vdp7ZzrK8S3MFxCsE0RSvMbGm4bhpwOrAVeJNg7rM/trOPDHefFU5rdB9w6QneX1+CyRP/0swWAz8imEhyMvALDs+5N4tgtt66sK7/AmqBG4Hz3L3RzH4KzAN+Ge53mbv/xQmOL72cQkbkSF8EzjKz68LnecAEoAFYHhMwAH9mZteGjwvC7XYdZ9/nA4vCe6DsMLPXgZnA/nDfVQAW3IRrLO2HTOstA1aF25xIA/Bi+Pg9oD4MjPfavP5ld98VHv+ZsNYm4GyC0AHIIbg5GAR3fv1NHMeXXk4hI3IkA77r7i8dsTDofqpt8/xS4Bx3rzOz14DsThy3PuZxM8f+v1nfzjZNHNn1HVtHox+eoLCl9fXu3tLm3FLbSQyd4LP4hbv/oJ06DoVhKXJcOicjvV0N0D/m+UvAXRbccAwzO9XM+rbzujxgTxgwE4GSmHWNra9v4w3gxvC8zzDgQmB5At7DZoK7Q6aZWQFB19fJuszMBoddf9cQdNn9AbjOzIYDhOuLElCv9CJqyUhv9y7QHJ7Afhz4d4JupNXhyfedBL9023oRuNPM1gMfAWUx6x4G3jWz1e4+L2b5YoKT5GsJWgp/5e7bw5DqjDeBTQQXJKwHVndgH8sJur/ygVJ3XwlgZvcCS8wsDWgEvkM79yQRORZN9S8iIpFRd5mIiERGISMiIpFRyIiISGQUMiIiEhmFjIiIREYhIyIikVHIiIhIZBQyIiISmf8POPStyFJz1dkAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "#NBVAL_SKIP\n", "import matplotlib.pyplot as plt\n", "\n", "# Plot objective function decrease\n", "plt.figure()\n", "plt.loglog(relative_error)\n", "plt.xlabel('Iteration number')\n", "plt.ylabel('True relative error')\n", "plt.title('Convergence')\n", "plt.show()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This notebook is part of the tutorial \"Optimised Symbolic Finite Difference Computation with Devito\" presented at the Intel® HPC Developer Conference 2017." ] } ], "metadata": { "anaconda-cloud": {}, "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.7.2" }, "widgets": { "state": {}, "version": "1.1.2" } }, "nbformat": 4, "nbformat_minor": 1 }