{ "metadata": { "name": "" }, "nbformat": 3, "nbformat_minor": 0, "worksheets": [ { "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Rozwi\u0105zywanie stochastycznych r\u00f3wna\u0144 r\u00f3\u017cniczkowych z CUDA\n", "\n", "R\u00f3wnania stochastyczne s\u0105 niezwykle po\u017cytecznym narz\u0119dziem w modelowaniu zar\u00f3wno proces\u00f3w fizycznych, biolgicznych czy chemicznych a nawet ekonomicznych (wycena instrument\u00f3w pochodnych).\n", "\n", "Klasycznym przyk\u0142adem problemu z jakim si\u0119 spotykami przy numerycznym rozwi\u0105zywaniu r\u00f3wna\u0144 stochastycznych jest konieczno\u015b\u0107 u\u015brednienia po wielu niezale\u017cnych od siebie realizacjach procesu losowego. M\u00f3wi\u0105c wprost musimy rozwi\u0105za\u0107 numerycznie wiele razy to samo r\u00f3wnanie r\u00f3\u017cniczkowe, za ka\u017cdym razem zmieniaj\u0105c \"seed\" generatora liczb losowych. Jest to idealny problem dla urz\u0105dzenia GPU, gdzie generacja niezale\u017cnych trajektorii wielu kopii tego samego uk\u0142adu jest w stanie wykorzysta\u0107 maksymalnie jego mo\u017cliwo\u015bci obliczeniowe.\n", "\n", "Poni\u017cej przedstawiamy implementacj\u0119 algorytmu, wg. pierwszego przyk\u0142adu z pracy: http://arxiv.org/abs/0903.3852 \n", "\n", "R\u00f3\u017anic\u0105 b\u0119dzie skorzystanie z pycuda, zamiast C. Co ciekawe, taka modyfikacja jest w stanie przy\u015bpieszy\u0107 kernel obliczniowe o ok 25%. Spowodowane jest to zastosowaniem metoprogramowania. Pewne parametry, kt\u00f3re nie zmieniaj\u0105 si\u0119 podczas wykonywania kodu s\u0105 \"klejane\" do \u017ar\u00f3d\u0142a jako konkrente warto\u015bci liczbowe, co u\u0142atwia kompilatorowi `nvcc` optymalizacje.\n", "\n", "W tym przyk\u0142adzie wykorzystamy w\u0142asny generator liczb losowych i transformacj\u0119 Boxa-Mullera (zamiast np. curand). \n", "\n", "Przyk\u0142ad ten mo\u017ce by\u0107 z \u0142atwo\u015bci\u0105 zmodyfikowany na dowolny uk\u0142ad SDE, dlatego mo\u017cna do traktowa\u0107 jako szablon dla w\u0142asnych zagadnie\u0144.\n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Struktura programu\n", "\n", "### Szablony \n", "\n", "Niezwykle pomocne w programowaniu w pycuda jest zastosowanie metaprogramowania - to jest - piszemy program pisz\u0105cy nasz kernel. Tutaj mamy najprostszy wariant, po prostu pewne parametry r\u00f3wna\u0144, wpisujemy automatycznie do tekstu j\u0105dra. W pythonie jest przydatne formatowanie \"string\u00f3w\" np.:\n", "\n", "\n" ] }, { "cell_type": "code", "collapsed": false, "input": [ "print '%(language)04d a nawiasy {} ' % {\"language\": 1234, \"number\": 2}" ], "language": "python", "metadata": {}, "outputs": [ { "output_type": "stream", "stream": "stdout", "text": [ "1234 a nawiasy {} \n" ] } ], "prompt_number": 6 }, { "cell_type": "markdown", "metadata": {}, "source": [ "lub:" ] }, { "cell_type": "code", "collapsed": false, "input": [ "print '{zmienna} a nawiasy: {{}}'.format( **{\"zmienna\": 123} )" ], "language": "python", "metadata": {}, "outputs": [ { "output_type": "stream", "stream": "stdout", "text": [ "123 a nawiasy: {}\n" ] } ], "prompt_number": 4 }, { "cell_type": "markdown", "metadata": {}, "source": [ "W pewnych bardziej zaawansowanych przypadkach, mo\u017cna zastosowa\u0107 system szablon\u00f3w np. mako templates (w projekcie http://sailfish.us.edu.pl). \n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Struktura kernela\n", "\n", "J\u0105dro:\n", "\n", " __global__ void SDE(float *cx,float *cv,unsigned int *rng_state, float ct)\n", "\n", "jest funkcj\u0105 CUDA typu `__global__`, jako parametry przyjmuje tablice `cx` i `cv`, b\u0119d\u0105ce zmiennymi stanu uk\u0142adu dw\u00f3ch r\u00f3wna\u0144 r\u00f3\u017cniczkowch:\n", "\n", "$$ \\dot x = v$$\n", "$$ \\dot v = -2\\pi \\cos(2.0\\pi x) + A \\cos(\\omega t) + F - \\gamma v$$\n", "\n", "Ponadto w wywo\u0142aniu przekazujemy czas (przez warto\u015b\u0107) oraz wska\u017anik do stanu generatora liczb losowych na GPU.\n", "\n", "Funkje dost\u0119pne dla j\u0105dra z GPU to:\n", "\n", "generator liczb losowych o rozk\u0142adzie jednostajnym:\n", "\n", " __device__ float rng_uni(unsigned int *state)\n", "\n", "transformacja Boxa-Mullera:\n", "\n", " __device__ void bm_trans(float& u1, float& u2)\n", "\n", "i wreszczcie funkcja obliczaj\u0105ca prawe strony uk\u0142adu r\u00f3wna\u0144:\n", "\n", " __device__ inline void diffEq(float &nx, float &nv, float x, float v, float t)\n", "\n", "\n", "Zauwa\u017cmy, \u017ce dla poprawienia wydajno\u015bci, ka\u017cde wywo\u0142anie kernela, powoduje wielokrotne (okre\u015blone przez parametr `spp`) wykonanie p\u0119tli iteracyjnej. \n", "\n", "\n" ] }, { "cell_type": "code", "collapsed": false, "input": [ "import pycuda.gpuarray as gpuarray\n", "\n", "import numpy\n", "from pycuda.curandom import rand as curand\n", "from pycuda.compiler import SourceModule\n", "import pycuda.driver as cuda\n", "\n", "\n", "cuda.init()\n", "device = cuda.Device(0)\n", "ctx = device.make_context()\n", "print device.name(), device.compute_capability(),device.total_memory()/1024.**3,\"GB\"\n", "\n", "\n", "blocks = 2**11\n", "block_size = 2**8\n", "N = blocks*block_size\n", "\n", "omega = 4.9\n", "spp = 100\n", "dt = 2.0*np.pi/omega/spp\n", "pars = {'samples':spp,'N':N,'dt':dt,'gam':0.9,'d0':0.001,'omega':omega,'force':0.1,'amp':4.2}\n", "\n", "rng_src = \"\"\"\n", "#define PI 3.14159265358979f\n", "/*\n", " * Return a uniformly distributed random number from the\n", " * [0;1] range.\n", " */\n", "__device__ float rng_uni(unsigned int *state)\n", "{\n", "\tunsigned int x = *state;\n", "\n", "\tx = x ^ (x >> 13);\n", "\tx = x ^ (x << 17);\n", "\tx = x ^ (x >> 5);\n", "\n", "\t*state = x;\n", "\n", "\treturn x / 4294967296.0f;\n", "}\n", "/*\n", " * Generate two normal variates given two uniform variates.\n", " */\n", "__device__ void bm_trans(float& u1, float& u2)\n", "{\n", "\tfloat r = sqrtf(-2.0f * logf(u1));\n", "\tfloat phi = 2.0f * PI * u2;\n", "\tu1 = r * cosf(phi);\n", "\tu2 = r * sinf(phi);\n", "}\n", "\n", "\"\"\"\n", "\n", "src = \"\"\"\n", " __device__ inline void diffEq(float &nx, float &nv, float x, float v, float t)\n", "{{\n", "\tnx = v;\n", "\tnv = -2.0f * PI * cosf(2.0f * PI * x) + {amp}f * cosf({omega}f * t) + {force}f - {gam}f * v;\n", "}}\n", "\n", "__global__ void SDE(float *cx,float *cv,unsigned int *rng_state, float ct)\n", " {{\n", " int idx = blockDim.x*blockIdx.x + threadIdx.x;\n", " float n1, n2; \t\n", " unsigned int lrng_state;\n", " float xim, vim, xt1, vt1, xt2, vt2,t,x,v;\n", " lrng_state = rng_state[idx]; \n", " t = ct;\n", " x = cx[idx];\n", "\t v = cv[idx]; \n", " \n", " for (int i=1;i<={samples};i++) {{\n", " \tn1 = rng_uni(&lrng_state);\n", "\t\tn2 = rng_uni(&lrng_state);\n", "\t\tbm_trans(n1, n2);\n", "\tdiffEq(xt1, vt1, x, v, t);\n", "\t\txim = x + xt1 * {dt}f;\n", "\t\tvim = v + vt1 * {dt}f + sqrtf({dt}f * {gam}f * {d0}f * 2.0f) * n1;\n", "\t\tt = ct + i * {dt}f;\n", "\tdiffEq(xt2, vt2, xim, vim, t);\n", "\t\tx += 0.5f * {dt}f * (xt1 + xt2);\n", "\t\tv += 0.5f * {dt}f * (vt1 + vt2) + sqrtf(2.0f * {dt}f * {gam}f * {d0}f) * n2;\n", " }}\n", " cx[idx] = x;\n", "\t cv[idx] = v;\n", "\n", "\t rng_state[idx] = lrng_state;;\n", " \n", " }}\n", " \"\"\".format(**pars)\n", "\n", "mod = SourceModule(rng_src + src,options=[\"--use_fast_math\"])\n", "SDE = mod.get_function(\"SDE\")\n", "\n", "print \"kernel ready for \",block_size,\"N =\",N,N/1e6" ], "language": "python", "metadata": {}, "outputs": [ { "output_type": "stream", "stream": "stdout", "text": [ "GeForce GTX 680 (3, 0) 3.99932861328 GB\n", "kernel ready for " ] }, { "output_type": "stream", "stream": "stdout", "text": [ " 256 N = 524288 0.524288\n" ] } ], "prompt_number": 7 }, { "cell_type": "code", "collapsed": false, "input": [ "print spp,N" ], "language": "python", "metadata": {}, "outputs": [ { "output_type": "stream", "stream": "stdout", "text": [ "100 524288\n" ] } ], "prompt_number": 8 }, { "cell_type": "markdown", "metadata": {}, "source": [ "Maj\u0105c gotowe j\u0105dro, mo\u017cna wykonac testowe uruchomienie:" ] }, { "cell_type": "code", "collapsed": false, "input": [ "import time\n", "\n", "x = np.zeros(N,dtype=np.float32)\n", "v = np.ones(N,dtype=np.float32)\n", "rng_state = np.array(np.random.randint(1,2147483648,size=N),dtype=np.uint32)\n", "\n", "x_g = gpuarray.to_gpu(x)\n", "v_g = gpuarray.to_gpu(v)\n", "rng_state_g = gpuarray.to_gpu(rng_state)\n", "\n", "start = time.time()\n", "for i in range(0,200000,spp):\n", " t = i * 2.0 * np.pi /omega /spp;\n", " SDE(x_g, v_g, rng_state_g, np.float32(t), block=(block_size,1,1), grid=(blocks,1))\n", "\n", "ctx.synchronize()\n", "elapsed = (time.time() - start)\n", "x=x_g.get()\n", "print elapsed,N/1e6, 200000*N/elapsed/1024.**3,\"Giter/sek\"" ], "language": "python", "metadata": {}, "outputs": [ { "output_type": "stream", "stream": "stdout", "text": [ "4.28112101555 0.524288 22.8109062195 Giter/sek\n" ] } ], "prompt_number": 9 }, { "cell_type": "markdown", "metadata": {}, "source": [ "Wynikiem dzia\u0142ania programu jest $N$ liczb okre\u015blaj\u0105cych ko\u0144cowe po\u0142o\u017cenie cz\u0105stki. Mo\u017cemy zwizualizowa\u0107 je wykorzystuj\u0105c np. hostogram:" ] }, { "cell_type": "code", "collapsed": false, "input": [ "h=histogram(x,bins=50,range=(-150, 100) ) \n", "figsize(8,4)\n", "plot(h[1][1:],h[0])" ], "language": "python", "metadata": {}, "outputs": [ { "metadata": {}, "output_type": "pyout", "prompt_number": 17, "text": [ "[]" ] }, { "metadata": {}, "output_type": "display_data", "png": "iVBORw0KGgoAAAANSUhEUgAAAfgAAAD9CAYAAABQkmV4AAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAIABJREFUeJzt3X9clfX9//HHMSgrg5QSjGPD5BihqJQRW6tohL/6hBZN\n47MQS2vDamZus9Yy3D6Kbf2YmVZruIh+4I9SWUu+lpM1a2ERZknLs0SFI1KGOioUwev7x1tP/hbh\nwHU4PO+327lBF+c6PK8rPK/zfl/v6/12WJZlISIiIgGli90BRERExPdU4EVERAKQCryIiEgAUoEX\nEREJQCrwIiIiAUgFXkREJAA1q8A3NTURHx/PDTfcAEBtbS0pKSn069ePoUOHsmvXLu9zc3JycLlc\nxMTEsHLlSu/20tJS4uLicLlcTJ482bt97969jB07FpfLRWJiIlu2bPHVsYmIiHRazSrwc+bMITY2\nFofDAcDs2bNJSUlh48aNJCcnM3v2bADKy8tZuHAh5eXlFBUVMWnSJA7eZp+VlUVubi5utxu3201R\nUREAubm5hIWF4Xa7mTJlCtOmTWuL4xQREelUTlrgq6qqeOONN5g4caK3WBcWFpKZmQlAZmYmy5Yt\nA2D58uWkp6cTHBxMVFQU0dHRlJSUUF1dTV1dHQkJCQCMGzfOu8+hr5WWlsaqVat8f5QiIiKdzEkL\n/JQpU/jDH/5Aly7fPbWmpobw8HAAwsPDqampAWDbtm04nU7v85xOJx6P56jtkZGReDweADweD717\n9wYgKCiI0NBQamtrfXBoIiIinVfQiX74+uuv07NnT+Lj4ykuLj7mcxwOh7frvi21x+8QERHxJ62Z\nTf6ELfh3332XwsJC+vTpQ3p6On//+9/JyMggPDyc7du3A1BdXU3Pnj0B0zKvrKz07l9VVYXT6SQy\nMpKqqqqjth/cZ+vWrQA0Njaye/duevTocdwD1aPtHg8//LDtGQL9oXOs8xwoD53jtn+01gkL/KxZ\ns6isrKSiooKCggJ+9KMfkZ+fT2pqKnl5eQDk5eUxevRoAFJTUykoKKChoYGKigrcbjcJCQlEREQQ\nEhJCSUkJlmWRn5/PqFGjvPscfK0lS5aQnJzc6oMSERHp7E7YRX+kg93k999/P2PGjCE3N5eoqCgW\nLVoEQGxsLGPGjCE2NpagoCDmz5/v3Wf+/PmMHz+e+vp6Ro4cyfDhwwGYMGECGRkZuFwuwsLCKCgo\n8OXxiYiIdEoOyxf9AO3A4XD4pMtCjq+4uJikpCS7YwQ0neP2ofPc9nSO215r654KvIiIiB9qbd3T\nVLUiIiIBSAVeREQkAKnAi4iIBCAVeBERkQCkAi8iIhKAVOBFREQCkAq8iIhIAFKBFxERCUAq8CIi\nIgFIBV5ERCQAqcCLiIgEIBV4kQC0cydUVNidQkTspAIvEiAaGmD5crj5ZoiKgiFDICcHmprsTiYi\ndlCBF+nALAveew/uugsiI+HRR2HoUNi8Gdatg5Ur4Uc/gq1b7U4qIu0tyO4AInLqNm2CF180D4CM\nDFi7Fvr0+e453bvDW2+Zoj9kCMydC2PH2pNXRNqf1oMX6UC++QZ+9StYtAhuucUU9ssvB4fjxPuV\nlsL//i8kJppCHxLSPnlFpOW0HrxIJ/Gvf8HgwVBXB263KdQJCScv7gCXXQYffghnnAHx8ea1RCSw\nqQUv4ucaGmDGDMjNhXnzIC2tda+3dCn87Gfmuv2vfw1BulAn4pfatAW/Z88errjiCgYPHkxsbCwP\nPPAAANnZ2TidTuLj44mPj2fFihXefXJycnC5XMTExLBy5Urv9tLSUuLi4nC5XEyePNm7fe/evYwd\nOxaXy0ViYiJbtmxp8cGIBJqPPzat9I8/ho8+an1xB7jxRigrg3/+03yvUfYigemEBb5r166sXr2a\ndevWsX79elavXs2aNWtwOBzcd999lJWVUVZWxogRIwAoLy9n4cKFlJeXU1RUxKRJk7yfPrKyssjN\nzcXtduN2uykqKgIgNzeXsLAw3G43U6ZMYdq0aW18yCL+r6kJfv97MwJ+8mRz+1t4uO9e/4ILYMUK\n090/Y4bvXldE/MdJr8GfddZZADQ0NNDU1ET37t0BjtltsHz5ctLT0wkODiYqKoro6GhKSkqorq6m\nrq6OhIQEAMaNG8eyZcsAKCwsJDMzE4C0tDRWrVrlmyMT6aA+/xyuucYU4Pffh9tua9519lMVFAQL\nF8Lzz8Nf/+r71xcRe520wO/fv5/BgwcTHh7OtddeS//+/QGYO3cugwYNYsKECezatQuAbdu24XQ6\nvfs6nU48Hs9R2yMjI/F4PAB4PB569+4NQFBQEKGhodTW1vruCEU6kI8+MiPdb74ZVq0yE9a0pfBw\nMyJ/wgQzcE9EAsdJh9d06dKFdevWsXv3boYNG0ZxcTFZWVlMnz4dgIceeoipU6eSm5vb5mGzs7O9\n3yclJZGUlNTmv1OkvezdC7feCo8/bm5/ay+JifC738FNN5lJc84+u/1+t4h8p7i4mOLiYp+9XrPH\nz4aGhnL99dfzwQcfHFZYJ06cyA033ACYlnllZaX3Z1VVVTidTiIjI6mqqjpq+8F9tm7dygUXXEBj\nYyO7d++mR48ex8xwaIEXCTTTp0O/fqbIt7c774SSEpg4EV5+uW0uCYjIiR3ZcJ3RygEyJ+yi37Fj\nh7f7vb6+njfffJP4+Hi2b9/ufc7SpUuJi4sDIDU1lYKCAhoaGqioqMDtdpOQkEBERAQhISGUlJRg\nWRb5+fmMGjXKu09eXh4AS5YsITk5uVUHJNIRrVkD+fnwzDP2FFeHw9yCt3Ej/PGP7f/7RcT3TtiC\nr66uJjMzk/3797N//34yMjJITk5m3LhxrFu3DofDQZ8+fXj22WcBiI2NZcyYMcTGxhIUFMT8+fNx\nHHi3mj9/PuPHj6e+vp6RI0cyfPhwACZMmEBGRgYul4uwsDAKCgra+JBF/EtdHWRmmuJ+/vn25Tjz\nTHjtNbjiCrj0UjPQT0Q6Lk10I2KzO+80t8W1wzCWZnnzTfOBY+1aOGRsrIi0s9bWPRV4ERv97W9w\n991m9Lw/zQ+fkwOFhVBcbKa3FZH2pwIv0kHt2AGDBplBbf7WHW5ZZta88HB4+mm704h0TlpsRqQD\nsizIyoL0dP8r7mAG3T3/PKxeDS+8YHcaEWkJteBFbPDSSzBrllnGtWtXu9Mc30cfQUoKlJfDeefZ\nnUakc1EXvUgHU1VlRqkXFZmv/m7KFPj6a3juObuTiHQuKvAiHcj+/TBsmOmW/81v7E7TPLt3Q2ws\nvPqqmfVORNqHrsGLdCBPPw3//S/cf7/dSZovNNSsbHfXXVpaVqQjUQtepJ1s3gxDhsA778DFF9ud\n5tRYFiQlwS23mMGBItL21EUv0kGkpcHgwfDQQ3YnaZlPPjHr03/yCfTsaXcakcCnAi/SAbz1Ftxx\nhxmNfuaZdqdpualTYedOWLDA7iQigU8FXsTP7dtnWu7/939w4412p2md//4XLrkEFi+GH/zA7jQi\ngU2D7ET83Pz5cMEFMHq03UlaLyQEHn3UDLhrbLQ7jYiciFrwIm3oyy/NLWb/+If5Gggsy1yLT0sz\n8+iLSNtQF72IH7vzTjj7bHjiCbuT+NaGDWZU/SefmPnqRcT3VOBF/FRpKVx/Pfz733DuuXan8b1f\n/tL0UDz/vN1JRAKTCryIH7IsuPJKmDDBPAJRXZ0ZcFdQAD/8od1pRAKPBtmJ+KGXXjKj52+7ze4k\nbeecc+CxxzTgTsRfqcCL+FhdHUybBk8+CV0C/F/YmDFmlblnn7U7iYgcSV30Ij52//1QXQ15eXYn\naR+lpeYWwM8/h9NPtzuNSOBo0y76PXv2cMUVVzB48GBiY2N54IEHAKitrSUlJYV+/foxdOhQdu3a\n5d0nJycHl8tFTEwMK1eu9G4vLS0lLi4Ol8vF5MmTvdv37t3L2LFjcblcJCYmsmXLlhYfjIjd3G74\n859h9my7k7Sfyy4ztwDm59udREQOdcIC37VrV1avXs26detYv349q1evZs2aNcyePZuUlBQ2btxI\ncnIysw+8m5WXl7Nw4ULKy8spKipi0qRJ3k8fWVlZ5Obm4na7cbvdFBUVAZCbm0tYWBhut5spU6Yw\nbdq0Nj5kkbYzZYrpnu/Vy+4k7euBB+CRR7TanIg/OekVwrPOOguAhoYGmpqa6N69O4WFhWRmZgKQ\nmZnJsmXLAFi+fDnp6ekEBwcTFRVFdHQ0JSUlVFdXU1dXR0JCAgDjxo3z7nPoa6WlpbFq1SrfH6VI\nO/jb30wL/pAOqk7jmmvMtfhXX7U7iYgcFHSyJ+zfv59LL72Uzz//nKysLPr3709NTQ3hB2a3CA8P\np6amBoBt27aRmJjo3dfpdOLxeAgODsbpdHq3R0ZG4vF4APB4PPTu3duECQoiNDSU2tpaevTocVSW\n7Oxs7/dJSUkkJSWd+hGLtIF9+0zrfc6cznkd2uEwrfjp0+HHPzb/LSKnpri4mOLiYp+93kkLfJcu\nXVi3bh27d+9m2LBhrF69+rCfOxwOHO30r/nQAi/iT155BSIjYcQIu5PY5/rr4cEHoaioc58HkZY6\nsuE6Y8aMVr1es2/iCQ0N5frrr6e0tJTw8HC2b98OQHV1NT0PLA4dGRlJZWWld5+qqiqcTieRkZFU\nVVUdtf3gPlu3bgWgsbGR3bt3H7P1LuKvmppg1iz4zW/sTmKvLl3MHQQ5OXYnERE4SYHfsWOHd4R8\nfX09b775JvHx8aSmppJ34B6gvLw8Rh9YJis1NZWCggIaGhqoqKjA7XaTkJBAREQEISEhlJSUYFkW\n+fn5jBo1yrvPwddasmQJycnJbXawIm3h1Vehe3ezAEtnN2YMeDywZo3dSUTkhF301dXVZGZmsn//\nfvbv309GRgbJycnEx8czZswYcnNziYqKYtGiRQDExsYyZswYYmNjCQoKYv78+d7u+/nz5zN+/Hjq\n6+sZOXIkw4cPB2DChAlkZGTgcrkICwujoKCgjQ9ZxHf27zfrvOfk6LozQFAQ/OpX5nz87W92pxHp\n3DTRjUgrFBbCww/Dhx+qwB+0Zw/07WsK/ODBdqcR6bg0F72ITSwLZs40195V3L/Ttau5o6AzTfYj\n4o9U4EVa6K23zLzzN95odxL/89OfwqpVZl4AEbGHCrxIC/3f/8Gvfx34C8q0xDnnwKRJ8Pvf251E\npPPSNXiRFnj7bbMU7GefmYFlcrSvvgKXC9avh0PmuRKRZtI1eBEbzJxpZm5TcT++sDAYPx4ef9zu\nJCKdk1rwIqfo/fchLQ3+85/OOS3tqfB4IC4ONm40c9WLSPOpBS/SzmbONPd6q7ifXGSk+TA0d67d\nSUQ6H7XgRU7B+vUwbBhs2gRnnml3mo7B7YYf/MCcs3POsTuNSMehFrxIO5o1C6ZOVXE/FS4XJCfD\nn/9sdxKRzkUteJFm+uwzuOoq0xLt1s3uNB3LP/8Jd94J5eWaFEikudSCF2knOTlwzz0q7i3xwx+a\nr1qERqT9qAUv0gwVFTBkCHz+OZx7rt1pOqbHH4d16+CFF+xOItIxtLbuqcCLNENWFvToYUbQS8vs\n2AHR0ebDUvfudqcR8X8q8CJtbNcuiIoyo8HPP9/uNB1beroZUX/PPXYnEfF/ugYv0sYWL4brrlNx\n94U774Q//cmsxCcibUsFXuQkXngBxo2zO0VgSEoy68WXlNidRCTwqcCLnMCmTeb2uOHD7U4SGBwO\nuOMO04oXkbala/AiJ/Db35rBYU8+aXeSwPHFF9CvH2zZAqGhdqcR8V+6Bi/SRixL3fNtoWdPSEmB\nl1+2O4lIYDtpga+srOTaa6+lf//+DBgwgCcPNGWys7NxOp3Ex8cTHx/PihUrvPvk5OTgcrmIiYlh\n5cqV3u2lpaXExcXhcrmYPHmyd/vevXsZO3YsLpeLxMREtmzZ4stjFGmRf/3LLChz2WV2Jwk8Gmwn\n0vZOWuCDg4N54okn2LBhA++99x7z5s3j008/xeFwcN9991FWVkZZWRkjRowAoLy8nIULF1JeXk5R\nURGTJk3ydjFkZWWRm5uL2+3G7XZTVFQEQG5uLmFhYbjdbqZMmcK0adPa8JBFmudg611Tq/pecjLs\n3g2lpXYnEQlcJy3wERERDB48GIBu3bpxySWX4PF4AI55bWD58uWkp6cTHBxMVFQU0dHRlJSUUF1d\nTV1dHQkJCQCMGzeOZcuWAVBYWEhmZiYAaWlprFq1yjdHJ9JCe/aY2+N+8hO7kwSmLl1g4kQNthNp\nS0Gn8uTNmzdTVlZGYmIi77zzDnPnzuWFF15gyJAhPPbYY5x77rls27aNxMRE7z5OpxOPx0NwcDBO\np9O7PTIy0vtBwePx0Lt3bxMoKIjQ0FBqa2vp0aPHYb8/Ozvb+31SUhJJSUmnerwizfL66xAfDwf+\nLKUN3HYbxMbCY49pGVkRgOLiYoqLi332es0u8F9//TU333wzc+bMoVu3bmRlZTF9+nQAHnroIaZO\nnUpubq7Pgh3LoQVepC1pcF3b69XL3BdfUGBunRPp7I5suM6YMaNVr9esUfT79u0jLS2NW2+9ldGj\nRwPQs2dPHA4HDoeDiRMnsnbtWsC0zCsrK737VlVV4XQ6iYyMpKqq6qjtB/fZunUrAI2Njezevfuo\n1rtIe/niC3j7bbjpJruTBL477oDnnrM7hUhgOmmBtyyLCRMmEBsby7333uvdXl1d7f1+6dKlxMXF\nAZCamkpBQQENDQ1UVFTgdrtJSEggIiKCkJAQSkpKsCyL/Px8Ro0a5d0nLy8PgCVLlpCcnOzTgxQ5\nFQUFkJqqZWHbw7BhsH27WWVORHzrpBPdrFmzhquvvpqBAwfiODCceNasWbzyyiusW7cOh8NBnz59\nePbZZwkPD/f+fMGCBQQFBTFnzhyGDRsGmNvkxo8fT319PSNHjvTecrd3714yMjIoKysjLCyMgoIC\noqKiDg+qiW6knQwZArNnm/nnpe3NmGF6TebNszuJiH/RanIiPrRhg2lVbtkCp51md5rOobISBg0y\nX88+2+40Iv5DM9mJ+FB+Ptx6q4p7e+rd2ywhu3ix3UlEAosKvMgBTU3w4ouQkWF3ks7n4Mx2IuI7\nKvAiB6xeDRER0L+/3Uk6n5EjzWWRTz6xO4lI4FCBFzlA977bJygIbr9drXgRX9IgOxHg66/B6YSN\nG81qZ9L+tmyBSy81g+3OOsvuNCL20yA7ER947TW4+moVdzt973uQmAiLFtmdRCQwqMCLoO55f/Gz\nn8Ezz9idQiQwqIteOr3KShg8GDwe6NrV7jSdW2Mj9OljFvsZNMjuNCL2Uhe9SCu99BLcfLOKuz8I\nCjLLyD77rN1JRDo+teClU7MsuOQSyM2FK6+0O40AVFXBwIGwdavWA5DOTS14kVb4+98hONjMpCb+\nwemEq64yi/6ISMupwEun9tRTcPfdcGAdJfETGmwn0nrqopdOa+tWiI8391+rK9i/NDVB377w6qtw\n2WV2pxGxh7roRVromWfMvPMq7v7ntNPgjjs02E6kNdSCl05pzx648EJYswb69bM7jRxLdTXExpoe\nlpAQu9OItD+14EVaYNEiMy2qirv/6tULkpPh5ZftTiLSManAS6c0b54ZXCf+7ac/NZdS1HkncupU\n4KXTWbsWvvgCRoywO4mcTHIy1NWZ/2cicmpU4KXTmTcPJk0yA7nEv3XpAnfeqcF2Ii1x0gJfWVnJ\ntddeS//+/RkwYABPPvkkALW1taSkpNCvXz+GDh3Krl27vPvk5OTgcrmIiYlh5cqV3u2lpaXExcXh\ncrmYPHmyd/vevXsZO3YsLpeLxMREtmzZ4stjFPH68ksoLDRrj0vHcNttZrW/Q95iRKQZTlrgg4OD\neeKJJ9iwYQPvvfce8+bN49NPP2X27NmkpKSwceNGkpOTmT17NgDl5eUsXLiQ8vJyioqKmDRpkncU\nYFZWFrm5ubjdbtxuN0VFRQDk5uYSFhaG2+1mypQpTJs2rQ0PWTqz3Fy46SYIC7M7iTRXz54wfDi8\n+KLdSUQ6lpMW+IiICAYPHgxAt27duOSSS/B4PBQWFpKZmQlAZmYmy5YtA2D58uWkp6cTHBxMVFQU\n0dHRlJSUUF1dTV1dHQkJCQCMGzfOu8+hr5WWlsaqVat8f6TS6TU2wtNPw1132Z1ETtVPf2q66TXY\nTqT5gk7lyZs3b6asrIwrrriCmpoawsPDAQgPD6empgaAbdu2kZiY6N3H6XTi8XgIDg7G6XR6t0dG\nRuLxeADweDz07t3bBAoKIjQ0lNraWnr06HHY78/OzvZ+n5SURFJS0qnEl07u9dchMtLcHicdS1IS\nNDTAu+9qUSAJXMXFxRQXF/vs9Zpd4L/++mvS0tKYM2cO55xzzmE/czgcONphMu9DC7zIqdKtcR2X\nw/HdYDsVeAlURzZcZ8yY0arXa9Yo+n379pGWlkZGRgajR48GTKt9+/btAFRXV9OzZ0/AtMwrKyu9\n+1ZVVeF0OomMjKSqquqo7Qf32bp1KwCNjY3s3r37qNa7SGt8+il8/DGkpdmdRFoqM9MMkKyttTuJ\nSMdw0gJvWRYTJkwgNjaWe++917s9NTWVvLw8APLy8ryFPzU1lYKCAhoaGqioqMDtdpOQkEBERAQh\nISGUlJRgWRb5+fmMGjXqqNdasmQJycnJPj9Q6dzmzzdzm59xht1JpKXOOw9SU+G55+xOItIxnHQu\n+jVr1nD11VczcOBAbzd8Tk4OCQkJjBkzhq1btxIVFcWiRYs499xzAZg1axYLFiwgKCiIOXPmMGzY\nMMDcJjd+/Hjq6+sZOXKk95a7vXv3kpGRQVlZGWFhYRQUFBAVFXV4UM1FLy1UVwff+x6sX2/WGpeO\n66OPYORI2LRJH9Yk8LW27mmxGQl48+fD6tWweLHdScQXhg2DW24x98eLBDIVeJETsCzo398Ued10\nERjeegsmTzZjKrpoLk4JYFpNTuQEVq82ReCaa+xOIr6SnAynnw4rVtidRMS/qcBLQHvqKXNrXDvc\nxSntxOGAX/4S/vAHu5OI+Dd10UvA+vRT0y2/aROcfbbdacSX9u2D6GgzruLA5JgiAUdd9CLHkZNj\nrtWquAee4GCYMkWteJETUQteAtKmTaZl9/nnEBpqdxppC19/DVFRUFICffvanUbE99SCFzmGRx6B\nn/1MxT2QdetmFqF5/HG7k4j4J7XgJeB4PBAXBxs3mtnPJHBt3w6xsfp/LYFJLXiRIzz6qJkERW/4\ngS8iwqwvMG+e3UlE/I9a8BJQvvwSLr4YPvkELrjA7jTSHv79bzPPQUUFnHWW3WlEfEcteJFDPPEE\njB2r4t6ZxMRAYiI8/7zdSUT8i1rwEjB27TKjqT/4APr0sTuNtKc1a2D8ePjsMzjtNLvTiPiGWvAi\nBzz1FNxwg4p7Z3TllXD++bB0qd1JRPyHWvASEL7+Gi66CP75T3MNXjqf114zt0e+956mJpbAoBa8\nCPDss2ZaWhX3zmvUKKitNR/yREQteAkAe/aY1vsbb8DgwXanETs9+yy8/jr89a92JxFpPbXgpdP7\ny1/g0ktV3AXGjYP334cNG+xOImI/teClQ9u3D1wueOUV+P737U4j/mDWLLOSYH6+3UlEWqe1dS/I\nh1lE2t1LL5lb41Tc5aC77jJ/E59/rkVopHM7aRf97bffTnh4OHFxcd5t2dnZOJ1O4uPjiY+PZ8WK\nFd6f5eTk4HK5iImJYeXKld7tpaWlxMXF4XK5mDx5snf73r17GTt2LC6Xi8TERLZs2eKrY5MA19Rk\nloR98EG7k4g/CQ2FrCwzol6kMztpgb/tttsoKio6bJvD4eC+++6jrKyMsrIyRowYAUB5eTkLFy6k\nvLycoqIiJk2a5O1eyMrKIjc3F7fbjdvt9r5mbm4uYWFhuN1upkyZwrRp03x9jBKgXn0VevSAa6+1\nO4n4m8mTYckSqKqyO4mIfU5a4K+66iq6d+9+1PZjXRdYvnw56enpBAcHExUVRXR0NCUlJVRXV1NX\nV0dCQgIA48aNY9myZQAUFhaSmZkJQFpaGqtWrWrVAUnnYFnmWuuDD+qeZznaeeeZBYcefdTuJCL2\nafE1+Llz5/LCCy8wZMgQHnvsMc4991y2bdtGYmKi9zlOpxOPx0NwcDBOp9O7PTIyEo/HA4DH46F3\n794mTFAQoaGh1NbW0qNHj6N+Z3Z2tvf7pKQkkpKSWhpfOrgVK0yRv/56u5OIv5o6FQYMgF//Gnr2\ntDuNyMkVFxdTXFzss9drUYHPyspi+vTpADz00ENMnTqV3Nxcn4U6nkMLvHRelgUzZ5o3brXe5Xgu\nuABuucUsQJSTY3cakZM7suE6Y8aMVr1ei+6D79mzJw6HA4fDwcSJE1m7di1gWuaVlZXe51VVVeF0\nOomMjKTqkIthB7cf3Gfr1q0ANDY2snv37mO23kUO+uc/4Ysv4Oab7U4i/u5Xv4I//Ql27rQ7iUj7\na1GBr66u9n6/dOlS7wj71NRUCgoKaGhooKKiArfbTUJCAhEREYSEhFBSUoJlWeTn5zNq1CjvPnl5\neQAsWbKE5OTk1h6TBLhZs+D++7VqmJxcVBSkpsLcuXYnEWl/J53oJj09nX/84x/s2LGD8PBwZsyY\nQXFxMevWrcPhcNCnTx+effZZwsPDAZg1axYLFiwgKCiIOXPmMGzYMMDcJjd+/Hjq6+sZOXIkTz75\nJGBuk8vIyKCsrIywsDAKCgqIioo6OqgmuhGgtBRGjzb3OJ9+ut1ppCP47DO46irYtAm6dbM7jUjz\ntbbuaSY76VBuvtm8WR8ylYLISY0dC5dfDr/4hd1JRJpPBV46jU8/NSvGbdoEZ59tdxrpSNavh+HD\nTc/PmWfanUakebTYjHQajzwCP/+5irucuoEDYcgQWLDA7iQi7UcteOkQNm+Gyy4zLbBzz7U7jXRE\nJSUwZgz85z8QHGx3GpGTUwteOoVHH4U771Rxl5a74gro1w9efNHuJCLtQy148Xvbt0NsrLkGf+Bm\nDZEWKS42HxQ//VS3WYr/UwteAt4f/wg/+YmKu7TeNdfA+efD4sV2JxFpe2rBi1/buROio+HDD+F7\n37M7jQTDsgUQAAAVFklEQVSCN94wEyWtWwdd1MQRP6YWvAS0efPMTGQq7uIrI0aYQXYHFrQUCVhq\nwYvf+uYb6NMH3n4bYmLsTiOB5G9/M634jz5SK178l1rwErCee85cM1VxF18bOdLMp7Bokd1JRNqO\nWvDil/buhb59obAQLr3U7jQSiN58E+6+GzZsgKAWLZwt0rbUgpeA9MILMGCAiru0neuug4gIeOkl\nu5OItA214MXvfPON6ZZfuBB+8AO700gge/ttGD/erDin2e3E36gFLwHn0UfhyitV3KXtXX21uRT0\nl7/YnUTE99SCF79SVQWDBum+d2k/771n5qjfuBG6drU7jch31IKXgPLAA5CVpeIu7Scx0aw299xz\ndicR8S214MVvrF0LN95orod262Z3GulMPvwQ/ud/zEpzZ51ldxoRQy14CQiWBffeCzNnqrhL+7v0\nUvj+9+Hpp+1OIuI7asGLXygogD/8Ad5/XzOLiT0++QSSk00r/pxz7E4j0g4t+Ntvv53w8HDi4uK8\n22pra0lJSaFfv34MHTqUXbt2eX+Wk5ODy+UiJiaGlStXereXlpYSFxeHy+Vi8uTJ3u179+5l7Nix\nuFwuEhMT2bJlS4sPRjqm+nqYNg2eeELFXewzYAD86Ecwd67dSUR846Rvp7fddhtFRUWHbZs9ezYp\nKSls3LiR5ORkZs+eDUB5eTkLFy6kvLycoqIiJk2a5P30kZWVRW5uLm63G7fb7X3N3NxcwsLCcLvd\nTJkyhWnTpvn6GMXPPfYYJCSYW5ZE7JSdbT5oHtJmEemwTlrgr7rqKrp3737YtsLCQjIzMwHIzMxk\n2YFlmZYvX056ejrBwcFERUURHR1NSUkJ1dXV1NXVkZCQAMC4ceO8+xz6Wmlpaaxatcp3Ryd+b9s2\ns977739vdxIRuPhiuP56U+RFOroWzcBcU1NDeHg4AOHh4dTU1ACwbds2EhMTvc9zOp14PB6Cg4Nx\nOp3e7ZGRkXg8HgA8Hg+9e/c2YYKCCA0Npba2lh49ehz1e7Ozs73fJyUlkZSU1JL44kcefBDuuMOs\nGifiD6ZPh8svh5//HMLC7E4jnUlxcTHFxcU+e71WL7HgcDhwOBy+yHJShxZ46fhKS+H//T/497/t\nTiLynYsugrQ0M+jzwNVHkXZxZMN1xowZrXq9Fg1pCg8PZ/v27QBUV1fTs2dPwLTMKysrvc+rqqrC\n6XQSGRlJVVXVUdsP7rN161YAGhsb2b179zFb7xJYDt4W99vfQkiI3WlEDveb38Cf/gQHOhpFOqQW\nFfjU1FTy8vIAyMvLY/To0d7tBQUFNDQ0UFFRgdvtJiEhgYiICEJCQigpKcGyLPLz8xk1atRRr7Vk\nyRKSk5N9cVzi5159Ferq4Lbb7E4icrQLL4R77jEPkQ7LOolbbrnF6tWrlxUcHGw5nU5rwYIF1ldf\nfWUlJydbLpfLSklJsXbu3Ol9/syZM62+fftaF198sVVUVOTd/sEHH1gDBgyw+vbta91zzz3e7Xv2\n7LF+/OMfW9HR0dYVV1xhVVRUHDNHM6JKB1Ffb1lRUZb197/bnUTk+OrrLeviiy3rtdfsTiKdVWvr\nnia6kXY3c6a5/v7aa3YnETmxf/wDfvIT2LABQkPtTiOdTWvrngq8tKtPPoFrr4UPPtCCMtIx3HEH\nnH46zJtndxLpbFTgpcPYt8+s3JWVBRMn2p1GpHl27oT+/WHJEvjBD+xOI52JFpuRDiMnB3r2hAkT\n7E4i0nzdu5vJmO68Exoa7E4j0nxqwUu7WLcOhg41y3IeMueRSIdgWXDDDWbFuQcftDuNdBbqohe/\n19BgZgabOhXGjbM7jUjLbN1qlpV9913o18/uNNIZqIte/N7vfgdRUZCRYXcSkZa78EJ46CHTVa+2\nhnQEasFLm3r/ffif/4GPPoKICLvTiLROU9N3A0Vvv93uNBLo1EUvfmvPHtOlOX063HKL3WlEfOPg\neJKPP4YDa26JtAkVePFb06bB55/D4sXQTusRibSLadOgshJeftnuJBLIVODFL/3rX3DTTbB+PZx/\nvt1pRHzr228hLg6eegpGjLA7jQQqDbITv/PttzB+vHnzU3GXQHTWWfDMM+Za/M6ddqcROTa14MXn\npkyBmhp1X0rgu/de2LgR/vpXOO00u9NIoFEXvfiV4mKzOMf69RAWZncakba1bx8MG2ZG1s+aZXca\nCTTqohe/sXkzpKfDggUq7tI5BAfDwoWmt2rxYrvTiBxOLXjxif/+F6680qy89fOf251GpH2VlZlb\n51atgoED7U4jgUJd9GK7piYYNcrMMf/007olTjqnl182M929/z706GF3GgkEKvBiu6lTzeQfRUWm\ny1Kks/rFL8z4kzfegKAgu9NIR6dr8GKrP//ZjCBevFjFXWT2bDNP/a9/bXcSkVYW+KioKAYOHEh8\nfDwJCQkA1NbWkpKSQr9+/Rg6dCi7du3yPj8nJweXy0VMTAwrV670bi8tLSUuLg6Xy8XkyZNbE0na\nUXGxWTrz9dfVJSkCptVeUABLlpivInZqVYF3OBwUFxdTVlbG2rVrAZg9ezYpKSls3LiR5ORkZs+e\nDUB5eTkLFy6kvLycoqIiJk2a5O16yMrKIjc3F7fbjdvtpqioqJWHJW3tP/8x88u//LKWzhQ5VFgY\nLF0K99xjLl2J2KXVXfRHXh8oLCwkMzMTgMzMTJYtWwbA8uXLSU9PJzg4mKioKKKjoykpKaG6upq6\nujpvD8C4ceO8+4h/2rnTrBA3YwYkJ9udRsT/DBoE8+bBjTfCjh12p5HOqtUt+Ouuu44hQ4bw3HPP\nAVBTU0P4gSWWwsPDqampAWDbtm04nU7vvk6nE4/Hc9T2yMhIPB5Pa2JJG9q3D8aMMZN7/PSndqcR\n8V9jxsDYsfDjH5vpm0XaW6vGeb7zzjv06tWLL7/8kpSUFGJiYg77ucPhwOHDe6ays7O93yclJZGU\nlOSz15bmufdec53xscfsTiLi/2bONOsyJCebwajnnWd3IvFnxcXFFBcX++z1WlXge/XqBcD555/P\njTfeyNq1awkPD2f79u1ERERQXV1Nz549AdMyr6ys9O5bVVWF0+kkMjKSqqqqw7ZHRkYe8/cdWuCl\nfVmWebMqLoZ339UtQCLNcdppkJdnRtX/8IfmVtKoKLtTib86suE6Y8aMVr1ei7vov/32W+rq6gD4\n5ptvWLlyJXFxcaSmppKXlwdAXl4eo0ePBiA1NZWCggIaGhqoqKjA7XaTkJBAREQEISEhlJSUYFkW\n+fn53n3EPzQ1wd13m1vh3nwTQkPtTiTScXTpYm6fu+suU+Q18E7aS4vbYTU1Ndx4440ANDY28pOf\n/IShQ4cyZMgQxowZQ25uLlFRUSxatAiA2NhYxowZQ2xsLEFBQcyfP9/bfT9//nzGjx9PfX09I0eO\nZPjw4T44NPGF+nq49VYzsO7tt1XcRVrqnnugVy8zpe0rr2iAqrQ9zWQnx7VzJ6Smmilon38ezjjD\n7kQiHd8//mEG4P3xj2ZxJpHj0Ux20ia2bjXdiQkJ8NJLKu4ivnLNNWZRmmnTNFhV2pYKvBzl44/N\nynATJpg3oC76KxHxqQED4J134C9/gfvug/377U4kgUhd9HKY4mJz7+6cOWamOhFpOzt3mpUYIyJg\nwQLo1s3uROJP1EUvPrN4sbk2+MorKu4i7aF7d1i5EkJC4LLLNMJefEsFXti503QTTpliboP70Y/s\nTiTSeXTtalZlfPhhSEkxU9yqs1J8QQW+E9u3D+bOhYsvhm++gdJSM4e2iLS///1fM4nUggWQlmY+\neIu0hgp8J2RZZonXuDgzfeaqVfDss3BgCQERsYnLZYp8794QHw//+pfdiaQj0yC7Tuajj2DqVPB4\nzAj5ESPAh8sFiIiPLF8Od95p/r3+4he6m6Uz0iA7aZbt2+GOO8wsWjfdBOvXw8iRKu4i/mrUKHj/\nfSgsNP9Wv/jC7kTS0ajAB7i6Ovjtb819t+eeC599BpMmQXCw3clE5GQuvNDcunrZZeaS2uOPw549\ndqeSjkIFPkDt3QtPPmmu6X32GaxdC3/4gynyItJxBAWZlRxXrTLrQbhcZtR9Y6PdycTfqcAHmKYm\nePFFiIkxS1MWFZmpZi+6yO5kItIaAwbAsmVmvoqXX4bYWFi4ULPgyfFpkF2AsCx44w144AE4+2yz\nPOU119idSkTagmWZFv0DD5gP9TNnwvDhGlMTaFpb91TgA8C778L998NXX8GsWWYFOP1DFwl8lgVL\nl8JvfgPnnQfZ2ZCUpBH3gUIFvhOyLNiwAV591Tx27zb/sMeNg9NOszudiLS3pibIz4dHHzWTVmVk\nmIfLZXcyaQ0V+E7CsqCszBT0JUugvt7MdpWWBt//vgq7iJj3iXXr4IUXzHX6vn3NB/8xY6BHD7vT\nyalSgQ9gTU1m9Purr8Jrr5lut4NF/fLL1Q0vIse3b59ZyOaFF8xg25QUU+yHD4fTT7c7nTSHCnyA\n2bTJLPjy1lvw97+bZSRvugluvhkGDlRRF5FTt2uXGX2fl2dms7zsMrjiCkhMNF8vuMDuhHIsKvAd\n3FdfmUL+1lumsNfXw3XXmU/byckQGWl3QhEJJDt3mhny3nsPSkrM17PP/q7YJyaaefDPOsvupBIw\nBb6oqIh7772XpqYmJk6cyLRp0w77eUcu8I2NUFkJFRWHPz79FP7zH7jqKlPQr7sO+ve3r5VeXFxM\nUlKSPb+8k9A5bh86z81nWfD554cX/E8+gV69zL32/fubr7GxZn6Nc84x++kct73W1r0gH2Zpsaam\nJu6++27eeustIiMjufzyy0lNTeWSSy6xO9pxNTWZbq8dO0wr/ODXg99/+SVs3mwKucdjVmrr08dM\nONOnj7kOdvfdpqvMX66H6R9s29M5bh86z83ncEB0tHnceqvZ1tho3rvKy80dO2++CXPmmFkxzzvP\nFPsdO4oZMSKJ884z284///CvXbvae1ziJwV+7dq1REdHExUVBcAtt9zC8uXLW13g9++HhobDH/v2\nmT/eI78e/H7fPlO4jyzahxbvr74yt6aFhJg/5LAw8zj0+4sugrFjTTG/8EI44wwfnCgRkXYQFGRu\nsXO5zKI3BzU1wZYtpug/84y5e8ftNsvafvnld42bHTtMw+Xg++GJHt26mffHrl2P/fWMMzT2qKX8\nosB7PB569+7t/W+n00lJSckpvcbEiWZt80OLeWOj+eM4/fTDH0FB3z2Cg4/+73PP/e6PLyLCTBF5\n5B9ljx66NU1EOpfTTjONl4sugtJSePjhYz/PssxCV0c2kA4+3G5zKeCrr8x9+3v3mkV0jvy6Z49p\ndJ12mnlvPtGjSxfzcDha9vXg9wfzH3osR25rjiM/lBzrQ8rBbcf72lp+cQ3+1VdfpaioiOeeew6A\nF198kZKSEubOnet9jkMf4UREpJPp8NfgIyMjqays9P53ZWUlTqfzsOf4wecQERGRDsMvZiweMmQI\nbrebzZs309DQwMKFC0lNTbU7loiISIflFy34oKAgnnrqKYYNG0ZTUxMTJkzw6xH0IiIi/s4vWvAA\nI0aM4LPPPiMnJ4cXX3yR0047jQ8//ND7882bN3PmmWcSHx9PfHw8kyZN8v6stLSUuLg4XC4XkydP\ntiN+h7J48WL69+9/1DkGyMnJweVyERMTw8qVK73bdY5bJzs7G6fT6f37XbFihfdnxzvncuqKioqI\niYnB5XLxyCOP2B0nYERFRTFw4EDi4+NJSEgAoLa2lpSUFPr168fQoUPZtWuXzSk7nttvv53w8HDi\n4uK82050Xk/5vcLyM59++qn12WefWUlJSVZpaal3e0VFhTVgwIBj7nP55ZdbJSUllmVZ1ogRI6wV\nK1a0S9aO6njneMOGDdagQYOshoYGq6Kiwurbt6+1f/9+y7J0jlsrOzvbeuyxx47afqxz3tTUZEPC\njq+xsdHq27evVVFRYTU0NFiDBg2yysvL7Y4VEKKioqyvvvrqsG2//OUvrUceecSyLMuaPXu2NW3a\nNDuidWhvv/229eGHHx5W2453XlvyXuE3LfiDYmJi6NevX7OfX11dTV1dnfdT5bhx41i2bFlbxQsI\nxzvHy5cvJz09neDgYKKiooiOjqakpETn2EesYwwUPdY5X7t2rQ3pOr5D59MIDg72zqchvnHk329h\nYSGZmZkAZGZm6j2hBa666iq6d+9+2LbjndeWvFf4XYE/kYqKCuLj40lKSmLNmjWAuYf+0BH3kZGR\neDweuyJ2aNu2bTvsXDqdTjwez1HbdY5bZu7cuQwaNIgJEyZ4u92Od87l1B1rPg2dS99wOBxcd911\nDBkyxHs7c01NDeHh4QCEh4dTU1NjZ8SAcbzz2pL3ClsG2aWkpLB9+/ajts+aNYsbbrjhmPtccMEF\nVFZW0r17dz788ENGjx7Nhg0b2jpqh9WScyytc7xzPnPmTLKyspg+fToADz30EFOnTiU3N/eYr6M5\nH1pG563tvPPOO/Tq1Ysvv/ySlJQUYmJiDvu5w+HQ+W8DJzuvJzvnthT4N99885T3Of300zn9wKTt\nl156KX379sXtdhMZGUlVVZX3eVVVVURqCbYWneMj5yOoqqrC6XTqHDdTc8/5xIkTvR+yjnXOdW5b\npjnzaUjL9OrVC4Dzzz+fG2+8kbVr1xIeHs727duJiIigurqanj172pwyMBzvvLbkvcKvu+gPveaz\nY8cOmpqaANi0aRNut5uLLrqIXr16ERISQklJCZZlkZ+fz+jRo+2K3OEceo5TU1MpKCigoaGBiooK\n3G43CQkJRERE6By3UnV1tff7pUuXekfNHu+cy6nTfBpt49tvv6Wurg6Ab775hpUrVxIXF0dqaip5\neXkA5OXl6T3BR453Xlv0XtEWIwNb47XXXrOcTqfVtWtXKzw83Bo+fLhlWZa1ZMkSq3///tbgwYOt\nSy+91Hr99de9+3zwwQfWgAEDrL59+1r33HOPXdE7jOOdY8uyrJkzZ1p9+/a1Lr74YquoqMi7Xee4\ndTIyMqy4uDhr4MCB1qhRo6zt27d7f3a8cy6n7o033rD69etn9e3b15o1a5bdcQLCpk2brEGDBlmD\nBg2y+vfv7z2vX331lZWcnGy5XC4rJSXF2rlzp81JO55bbrnF6tWrlxUcHGw5nU5rwYIFJzyvp/pe\n4Rdz0YuIiIhv+XUXvYiIiLSMCryIiEgAUoEXEREJQCrwIiIiAUgFXkREJACpwIuIiASg/w9aeEcW\n2fB7dgAAAABJRU5ErkJggg==\n", "text": [ "" ] } ], "prompt_number": 17 }, { "cell_type": "code", "collapsed": false, "input": [], "language": "python", "metadata": {}, "outputs": [] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Dane referencyjne dla walidacji\n", "\n", "W tablicy `hist_ref` znajduj\u0105 si\u0119 dane referencyjne dla cel\u00f3w walidacji. Mo\u017cemy sprawdzi\u0107 czy program dzia\u0142a tak jak ten w pracy referencyjnej:" ] }, { "cell_type": "code", "collapsed": false, "input": [ "hist_ref = (array([ 46, 72, 134, 224, 341, 588, 917, 1504, 2235,\\\n", " 3319, 4692, 6620, 8788, 11700, 15139, 18702, 22881, 26195,\\\n", " 29852, 32700, 35289, 36232, 36541, 35561, 33386, 30638, 27267,\\\n", " 23533, 19229, 16002, 12646, 9501, 7111, 5079, 3405, 2313,\\\n", " 1515, 958, 573, 370, 213, 103, 81, 28, 15,\\\n", " 7, 3, 2, 0, 0]),\\\n", " array([-150., -145., -140., -135., -130., -125., -120., -115., -110.,\\\n", " -105., -100., -95., -90., -85., -80., -75., -70., -65.,\\\n", " -60., -55., -50., -45., -40., -35., -30., -25., -20.,\\\n", " -15., -10., -5., 0., 5., 10., 15., 20., 25.,\\\n", " 30., 35., 40., 45., 50., 55., 60., 65., 70.,\\\n", " 75., 80., 85., 90., 95., 100.]) )\n", "\n", "\n", "\n", "hist(x,bins=50,range=(-150, 100) ) \n", "plot((hist_ref[1][1:]+hist_ref[1][:-1])/2.0,hist_ref[0],'r')" ], "language": "python", "metadata": {}, "outputs": [ { "metadata": {}, "output_type": "pyout", "prompt_number": 18, "text": [ "[]" ] }, { "metadata": {}, "output_type": "display_data", "png": "iVBORw0KGgoAAAANSUhEUgAAAfgAAAD9CAYAAABQkmV4AAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAIABJREFUeJzt3X1cVGX+//HXGNidYkAKxliYjhGKSqvIbuVShLff0LI0\ntgRTd39f3G+rrrVWmy1um9q2bpnFbttiIbuFrpW6lWRZdLcbJunWiuVUqDAipXgzmYLA+f0xOIGI\nIAycmeH9fDzOw/Hic858zlHOZ64z51yXxTAMAxEREfErXcxOQERERDxPBV5ERMQPqcCLiIj4IRV4\nERERP6QCLyIi4odU4EVERPxQiwp8TU0NsbGx3HjjjQBUVFSQlJTEgAEDGDVqFIcOHXLHLl68GJvN\nRlRUFBs3bnS3FxYWEhMTg81mY/bs2e72yspKpkyZgs1mIz4+nt27d3tq30RERDqtFhX4ZcuWER0d\njcViAWDJkiUkJSWxc+dOEhMTWbJkCQBFRUWsWrWKoqIi8vLymDVrFicfs09PTycrKwu73Y7dbicv\nLw+ArKwsQkNDsdvtzJ07l/nz57fHfoqIiHQqzRb40tJSXnvtNWbOnOku1uvXryctLQ2AtLQ01q5d\nC8C6detISUkhMDCQyMhI+vfvT0FBAWVlZTidTuLi4gBITU11r1N/W5MmTWLTpk2e30sREZFOptkC\nP3fuXB599FG6dPk+tLy8nLCwMADCwsIoLy8HYO/evVitVnec1WrF4XA0ao+IiMDhcADgcDjo06cP\nAAEBAfTo0YOKigoP7JqIiEjnFXCmH77yyiv06tWL2NhY8vPzTxtjsVjcl+7bU0e8h4iIiDdpy2jy\nZ+zB/+tf/2L9+vX07duXlJQU3nrrLaZOnUpYWBj79u0DoKysjF69egGunnlJSYl7/dLSUqxWKxER\nEZSWljZqP7nOnj17AKiurubw4cOEhIQ0uaNa2m/5zW9+Y3oO/r7oGOs4+8uiY9z+S1udscAvWrSI\nkpISiouLyc3N5frrrycnJ4fk5GSys7MByM7OZuLEiQAkJyeTm5tLVVUVxcXF2O124uLiCA8PJygo\niIKCAgzDICcnhwkTJrjXObmtNWvWkJiY2OadEhER6ezOeIn+VCcvk997771MnjyZrKwsIiMjWb16\nNQDR0dFMnjyZ6OhoAgICyMzMdK+TmZnJtGnTOHbsGOPGjWPMmDEAzJgxg6lTp2Kz2QgNDSU3N9eT\n+yciItIpWQxPXAfoABaLxSOXLKRp+fn5JCQkmJ2GX9Mx7hg6zu1Px7j9tbXuqcCLiIh4obbWPQ1V\nKyIi4odU4EVERPyQCryIiIgfUoEXERHxQyrwIiIifkgFXkRExA+pwIuIiPghFXgRERE/pAIvIiLi\nh1TgRURE/JAKvIiIiB9SgRcREfFDKvAiPigoKIQAi4UYi4VUi4XfWyxca7FgqbcEBYWYnaaImOis\n5oMXEZMcPw6ffgoffwxbt/Km8yADuYBSrGwllp0M4HmeZROJ3MOj7KcnTqfF7KxFxETqwYt4s//8\nB264AYKD4ac/hQ8/hOhofgmEs48oPieFXH7Db4mmiAOEsp2B/IynUXkX6dw0H7yIlwkKCiHAeZCH\ngFuA3wDPAZWNIk//+xDDJ/yZ/6UL/yb+448hNrYdsxWR9qL54EX8SU0NKc6D7KAXFtKJZj9PY1CJ\nAQ2Wpn3KYK7hfZ4BGDMG5syBI0faP3cR8Soq8CLe4l//grg47gBG8zo/J5MKQlu1KYMurADYvh2c\nTrjySli1CnQVTKTTOGOBP378OCNGjGDo0KFER0dz3333AZCRkYHVaiU2NpbY2Fg2bNjgXmfx4sXY\nbDaioqLYuHGju72wsJCYmBhsNhuzZ892t1dWVjJlyhRsNhvx8fHs3r3b0/so4t3KyiA1FSZPhnnz\nGAn8h6Ge2fbFF0NWFqxeDQsWwJ//7Jntioj3M5px9OhRwzAM48SJE8aIESOM9957z8jIyDCWLl3a\nKHb79u3GkCFDjKqqKqO4uNjo16+fUVtbaxiGYQwfPtwoKCgwDMMwxo4da2zYsMEwDMN46qmnjPT0\ndMMwDCM3N9eYMmXKafNoQaoivuf55w0jNNQw5s83DKfTMAyj7hq80czSkpiABtf1bWB8DUZMvbbu\n3YNNPgAi0pS21r1mL9FfcMEFAFRVVVFTU0NwcPDJDwaNYtetW0dKSgqBgYFERkbSv39/CgoKKCsr\nw+l0EhcXB0Bqaipr164FYP369aSlpQEwadIkNm3a1KYPLCI+o6QE7roL3nwTliyBbt08/AbV1P/e\n3o7BXHJYRRQX8C1g4HQe9PB7ioi3aLbA19bWMnToUMLCwrjuuusYOHAgAMuXL2fIkCHMmDGDQ4cO\nAbB3716sVqt7XavVisPhaNQeERGBw+EAwOFw0KdPHwACAgLo0aMHFRUVnttDEW81Zw783//BUA9d\njm+Bv3MHBYxgOXd12HuKiDmaHeimS5cubNu2jcOHDzN69Gjy8/NJT0/nwQcfBGDBggXMmzePrKys\ndk82IyPD/TohIYGEhIR2f08RTwoKCsHpPMh44DEg5qWXqFy4sENz+D+epJAf8BP+zvMd+s4icib5\n+fnk5+d7bHstHsmuR48ejB8/ni1btjQorDNnzuTGG28EXD3zkpIS989KS0uxWq1ERERQWlraqP3k\nOnv27OGSSy6hurqaw4cPExJy+iE26xd4EV/kdB7kAr7lSQbyU56hkqTTRLXvEDVH6cYUVvEGSWxu\n13cSkbNxasd1YRs//J/xEv3+/fvdl9+PHTvGG2+8QWxsLPv27XPHvPzyy8TExACQnJxMbm4uVVVV\nFBcXY7fbiYuLIzw8nKCgIAoKCjAMg5ycHCZMmOBeJzs7G4A1a9aQmJjYph0S8XYLeIh/8SPePG1x\n7xj/YSgZZJALUNl4CB0R8X1n7MGXlZWRlpZGbW0ttbW1TJ06lcTERFJTU9m2bRsWi4W+ffvy9NNP\nAxAdHc3kyZOJjo4mICCAzMxMLBZXbyQzM5Np06Zx7Ngxxo0bx5gxYwCYMWMGU6dOxWazERoaSm5u\nbjvvsoh5BgIzyCKGT81OhUxmcQP/xw/uvRcee8zsdETEwzRUrUhHqa3lvXPO4Xky+TPpZwi00Nxo\ndZ6KCcZCxaWXwpNPQt1XbSLiHTRUrYivePZZugJ/4WdmZ+J2EOD552HmTKh3n4yI+D714EU6wjff\nwMCBxH7zDds6qHfe0hjDMGDRIsjLg7feggDNIi3iDdSDF/EFv/oV3H4728zOoynz50PXrvDQQ2Zn\nIiIeoh68SHt75x244w4oKsISFERH9s6bjwnENeIdhANbgRuBLfUiuncP5sgRDT4l0tHUgxfxZlVV\nkJ4Oy5ZB9+5mZ3Ma3w9nuw+D+8niCeKxUONu13C2Ir5JBV6kPf3hD3D55XDTTWZn0iLPMY0Aqrmd\nv5udioi0kS7Ri7SXr76C4cNhyxbo2xegblwIb7pE3zgmnn+zhluI4jO+pTvuG/FEpEPpEr2It/rl\nL2HePHdx9xUf8kPe4nruZ5HZqYhIG6gHL+JBJyeTSQCeBaKAxgPBencPHqA3e/mEwcTzIV9i0++e\niAnUgxfxIk7nQbpQzR8Zynxyqaw3H3vzxdZ7lHEJf+BuljLP7FREpJVU4EU8bCo5HOc8VjPZ7FTa\n5DHmMpDtJk6JIyJtoUv0Ih50ocXC50RwC2soIP40EeZffj+bmBtZzyNM4MqqKggMbGZbIuJJukQv\n4kXuBt7j2iaKu+/5JzeyB1yT0YiIT1EPXsRTHA4OWK38gGJ2E9lEkPf0zlsaE4WFHRdfDNu3Q69e\nzWxPRDxFPXgRb/HAAzwDZyjuvukzcA21++tfm52KiJwF9eBFPGHrVhg7lh7l5Rw5Y4/Yu3rnLZ5x\n7uBBiIqCV1+FH/ygmXgR8QT14EXMZhiuAW0yMjhidi7t5aKLXDPNzZ7t2l8R8Xoq8CJt9c9/Qnk5\nzJxpdibta/p0+O47eOEFszMRkRbQJXqRtjhxAgYNcs0WN2ZMC8aa977L72czpezVwAvAAOB4vQhN\nKSviee16if748eOMGDGCoUOHEh0dzX333QdARUUFSUlJDBgwgFGjRnHo0CH3OosXL8ZmsxEVFcXG\njRvd7YWFhcTExGCz2Zg9e7a7vbKykilTpmCz2YiPj2f37t2t3hmRDvfnP0NkJIwZY3Ym7ej7KWU/\nwOBjkpnOk9QfoU9Tyop4nzMW+PPOO4+3336bbdu28cknn/D222/z/vvvs2TJEpKSkti5cyeJiYks\nWbIEgKKiIlatWkVRURF5eXnMmjXL/ekjPT2drKws7HY7drudvLw8ALKysggNDcVutzN37lzmz5/f\nzrss4iEHD7q+l/7DH8zOpEMt5j7u4VECOGF2KiJyBs1+B3/BBRcAUFVVRU1NDcHBwaxfv560tDQA\n0tLSWLt2LQDr1q0jJSWFwMBAIiMj6d+/PwUFBZSVleF0OomLiwMgNTXVvU79bU2aNIlNmzZ5fi9F\n2sPvfuea5z0mxuxMOlQB8XzF5aSg7+JFvFlAcwG1tbVcddVVfPnll6SnpzNw4EDKy8sJCwsDICws\njPLycgD27t1LfPz3I3hZrVYcDgeBgYFYrVZ3e0REBA6HAwCHw0GfPn1cyQQE0KNHDyoqKggJCWmU\nS0ZGhvt1QkICCQkJZ7/HIp7w5ZeQne0a/KUTWsx9LGM2f+MODN2rK+IR+fn55Ofne2x7zRb4Ll26\nsG3bNg4fPszo0aN5++23G/zcYrHU3VjU/uoXeJGOdnIqWIBngBLgt+HhpuZklje5gaNcyATWsZab\nzE5HxC+c2nFduHBhm7bX4o/ePXr0YPz48RQWFhIWFsa+ffsAKCsro1fd8JURERGUlJS41yktLcVq\ntRIREUFpaWmj9pPr7NmzB4Dq6moOHz582t67iNlcxd2gJ+VM4iIy+ZqGU8F2pqc8LCzifu5nEZ1r\nv0V8xxkL/P79+913yB87dow33niD2NhYkpOTyc7OBiA7O5uJEycCkJycTG5uLlVVVRQXF2O324mL\niyM8PJygoCAKCgowDIOcnBwmTJjgXufkttasWUNiYmK77ayIJ/w/nmYNt7CfnmanYqp1TOBCjpKI\n7psR8UZnvERfVlZGWloatbW11NbWMnXqVBITE4mNjWXy5MlkZWURGRnJ6tWrAYiOjmby5MlER0cT\nEBBAZmam+/J9ZmYm06ZN49ixY4wbN44xdY8VzZgxg6lTp2Kz2QgNDSU3N7edd1mk9bpSySwySeIN\ns1MxnUEXlnAv97NIJV7EC2mgG5EWslgs3MFKUlnJqCYLvD8OdNN0TAAnsGNjCrsp0O+niEdpLHqR\nDjSXx3icOWan4TWqCeRR7uE+sxMRkUZU4EVa6BqgG9+ygbFmp+JVVjCdeID//tfsVESkHhV4kRaa\nAyxjtp77PsVxzudxgEceMTsVEalH38GLtERxMQcuv5zLcHKUbmcI7FzfwZ8UhIXDoaHw0UfQt28z\n2xORltB38CId4cknWQHNFPfO6wjA//t/8OijZqciInXUgxdpjtMJkZFcVlHBnjb3dr2r5+3JGKO8\nHKKioKgIOukIfyKepB68SHt77jm4/nr2mJ2Ht+vVC+64Ax57zOxMRAT14EXOrLYWBgyA7Gws11xD\n23u73tfz9lgP3jBgzx6IjYUvvoDg4GbWEZEzUQ9epD29+qqrUP3oR2Zn4hsuvRSSk+Gpp8zORKTT\nUw9e5EwSE2H6dLj99rphl9WDP71AoBqAK4B3gEjgeL2I7t2DOXKkopntiMhJ6sGLtJdPPoHPPoNb\nbzU7Ex9QzckZ9T7H4CPGk0KWuw0M91S7ItIxVOBFmrJsGcyaBV27mp2Jz3mCX/ALnkBTyYqY54yz\nyYl0Wl9/DS+9BDt3mp2JT3qTGziP41zLe7zHSLPTEemU1IMXAYKCQrBYLO5lQVgYzxw6hKVXL3eb\ntJxBF5ZzF3ex3OxURDot3WQnAg1uoOtKJbuIJIk32M6g+lHoJruWx3TDyS4iiWUrJVyK+1E6EWkR\n3WQn4mG3kct/GXRKcZez9S3dyWEq6fzJ7FREOiX14EX4vgffhRq2M5Cf8xRvkXhqFOrBn11MP77g\n3/yQS9nDcS7Q77DIWVAPXsSDJvEih7iIt7je7FT8wpf0p4ARpPCC2amIdDoq8CJuBr/mYX7HA7h6\npOIJy7mr7pE5EelIzRb4kpISrrvuOgYOHMigQYN44gnXL2pGRgZWq5XY2FhiY2PZsGGDe53Fixdj\ns9mIiopi48aN7vbCwkJiYmKw2WzMnj3b3V5ZWcmUKVOw2WzEx8eze/duT+6jSIv8D69gYOFVxpud\nil95g6S6R+ZEpEMZzSgrKzO2bt1qGIZhOJ1OY8CAAUZRUZGRkZFhLF26tFH89u3bjSFDhhhVVVVG\ncXGx0a9fP6O2ttYwDMMYPny4UVBQYBiGYYwdO9bYsGGDYRiG8dRTTxnp6emGYRhGbm6uMWXKlEbb\nbUGqIq0GGP9mhHELqw0wmlg4w89aGuOJbfhezCyeNFbrd1jkrLS17jXbgw8PD2fo0KEAdOvWjSuv\nvBKHw3Hyw0Gj+HXr1pGSkkJgYCCRkZH079+fgoICysrKcDqdxMXFAZCamsratWsBWL9+PWlpaQBM\nmjSJTZs2tf4Ti0grJAI9OMxL3Gx2Kn5pJamuWxb3aNJdkY5yViPZ7dq1i61btxIfH88HH3zA8uXL\nWblyJcOGDWPp0qVcdNFF7N27l/j4ePc6VqsVh8NBYGAgVqvV3R4REeH+oOBwOOjTp48roYAAevTo\nQUVFBSEhIQ3ePyMjw/06ISGBhISEs91fkdN6AFjMfdRyjtmp+KVv6c5KYM6f/gSLF5udjohXys/P\nJz8/32Pba3GB//bbb7nllltYtmwZ3bp1Iz09nQcffBCABQsWMG/ePLKysjyW2OnUL/AiHvP++1wK\nvECK2Zn4tSeBOX/9Kzz4IJx/vtnpiHidUzuuCxcubNP2WnQX/YkTJ5g0aRJ33HEHEydOBKBXvSE8\nZ86cyebNmwFXz7ykpMS9bmlpKVarlYiICEpLSxu1n1xnT92lu+rqag4fPtyo9y7Sbh5+mEeAagLN\nzsSvfQkQFwcv6JE5kY7QbIE3DIMZM2YQHR3NnDlz3O1lZWXu1y+//DIxMTEAJCcnk5ubS1VVFcXF\nxdjtduLi4ggPDycoKIiCggIMwyAnJ4cJEya418nOzgZgzZo1JCaeOsCISDspLIRPP+U5s/PoLH7x\nC3jiibr78kSkPTU7kt3777/PyJEjGTx4sHvCjUWLFvHCCy+wbds2LBYLffv25emnnyYsLMz98xUr\nVhAQEMCyZcsYPXo04HpMbtq0aRw7doxx48a5H7mrrKxk6tSpbN26ldDQUHJzc4mMjGyYqEayk/Zw\n880wciSWuXPpmFHfzB9dzswYo6YGoqPhL3+BkZplTuRM2lr3NFStdF7//S8kJsJXX2Hp1g0V+PaP\nMQwDnnwS3nkH/vGPZuJFOjcVeJHWuv12iImBe+9tMJtc01Tg2xpjGAY4nXDZZbBtG1x6aTPriHRe\nKvAirfHFFxAfD199BUFBKvAdFOP+HZ4zx3UnvR6ZE2mSCrxIa8ycCRERUPcYigp8x8S4f4ftdrj6\naigtha5dm1lPpHNSgRc5W3v2wNChriITGgqowHdMTCBQ7f5bPvA4sLZeRPfuwRw5UtHMdkQ6B00X\nK3K2Hn3U1YOvK+7SUapxfQhwLdlkkcaEBm1O50EzExTxK+rBS+eybx9ceSXs2AHh4e5m9eA7PqY7\nR9jDpdiws5+e7hj9nou4qAcvcjb++leYPLlBcRdzOAnin9zIT3je7FRE/JIKvHQetbXw7LOuy/Pi\nFbJJI41ss9MQ8Usq8NJ5vPsuXHABDBtmdiZS522uoyffMIhPzU5FxO+owEvnsWIFTJ8OdUMui/lq\nOYccpqoXL9IOdJOd+L2goBAszoPsBmzA/iYjdZOdGTED+Jx8EuhDCTUE6vdcpI5ushNphtN5kCk8\nzSZuZn+9R7IaLmKWnVzBbi5jNK+bnYqIX1GBl05hOitYwXSz05AmPMc0XaYX8TBdohe/N9Bi4Q16\ncyl7qCGgiSg9B29mzEUcpJi+9OUwB/V7LgLoEr1Is+7E9ThW08VdzHaIYDYyiilmJyLiR1Tgxb+d\nOMFU4FnuNDsTaYbrmXgR8RQVePFvr77KTsDOALMzkWa8zmj6Anz+udmpiPgFFXjxbytWsMLsHKRF\nagjgbwDZutlOxBN0k534r7IyiI6m26FDHPWCG8k69n18M2YQFj61WmHXLjjnnGa2J+Lf2v0mu5KS\nEq677joGDhzIoEGDeOKJJwCoqKggKSmJAQMGMGrUKA4dOuReZ/HixdhsNqKioti4caO7vbCwkJiY\nGGw2G7Nnz3a3V1ZWMmXKFGw2G/Hx8ezevbvVOyTilpMDN9/MUbPzkBb7L0DPnvD222anIuLzmi3w\ngYGBPPbYY2zfvp0PP/yQp556ih07drBkyRKSkpLYuXMniYmJLFmyBICioiJWrVpFUVEReXl5zJo1\ny/0JJD09naysLOx2O3a7nby8PACysrIIDQ3Fbrczd+5c5s+f3467LJ2CYXw/NK34lrQ0XaYX8YBm\nC3x4eDhDhw4FoFu3blx55ZU4HA7Wr19PWprrnte0tDTWrl0LwLp160hJSSEwMJDIyEj69+9PQUEB\nZWVlOJ1O4uLiAEhNTXWvU39bkyZNYtOmTZ7fU+lc/v1vV5H/0Y/MzkTO1k9+Av/8Jxw5YnYmIj7t\nrB4M3rVrF1u3bmXEiBGUl5cTFhYGQFhYGOXl5QDs3buX+Ph49zpWqxWHw0FgYCBWq9XdHhERgcPh\nAMDhcNCnTx9XQgEB9OjRg4qKCkJCQhq8f0ZGhvt1QkICCQkJZ5O+dCaaWMZ39ewJCQmwZo2uwEin\nkp+fT35+vse21+IC/+233zJp0iSWLVtG9+7dG/zMYrFg6YATaf0CL9Kkb7+FF1+EoiKzM5HWSkuD\nxx9XgZdO5dSO68KFC9u0vRY9JnfixAkmTZrE1KlTmThxIuDqte/btw+AsrIyevXqBbh65iUlJe51\nS0tLsVqtREREUFpa2qj95Dp79uwBoLq6msOHDzfqvYu02Jo1cM010Lu32ZlIa40f7/qAVlxsdiYi\nPqvZAm8YBjNmzCA6Opo5c+a425OTk8muuxEmOzvbXfiTk5PJzc2lqqqK4uJi7HY7cXFxhIeHExQU\nREFBAYZhkJOTw4QJExpta82aNSQmJnp8R6UT0c11PizAdUXw3HNZvn8/v7n8cvcVwpNLUJA+/Iu0\nRLPPwb///vuMHDmSwYMHuy/DL168mLi4OCZPnsyePXuIjIxk9erVXHTRRQAsWrSIFStWEBAQwLJl\nyxg9ejTgekxu2rRpHDt2jHHjxrkfuausrGTq1Kls3bqV0NBQcnNziYyMbJionoOXlti5E669FkpK\noGtXgLr/t97yrLc35eLdMbF8zMvchA07J+jaIEbnAukM2lr3NNCN+Jf774fKSli61N2kAu+7MRtJ\nYhVTyGJmgxidC6QzUIEXOam6Gi67DF5/HQYNcjerwPtuzNW8z0pSuYLPqSbQHaNzgXQGmi5W5KRX\nXwWrtUFxF9/2AddQTF+mkmN2KiI+RwVe/MfSpVDvRlDxD7/lQX7NwwRwwuxURHyKCrz4tKCgECwW\nC8MtFva89x6BP/lJo7uuxbe9y48poQ8/4XmzUxHxKfoOXnzaye/XnyeFLQzjj8w7XRTe8z2zN+Xi\nOzEJvM1f+BlXsoMaAnUukE5B38FLp9eHPYxiI39tcKe1+JN8EiijNym8YHYqIj5DPXjxaRaLhT/w\nSwDuZmlTUZjdA/XOXHwr5no2kcksotlJjc4F0gm0te6d1WQzIt4mCJjGc8Sy1exUpJ29xfV8Q09u\nY6fZqYj4BF2iF582E9jIKEq41OxUpN1ZWMhvWABQU2N2MiJeTwVefNeJE8wGlp72xjrxR29yAxUA\nq1ebnYqI11OBF9+1Zg3FQCHDzM5EOoyFhQAPPaRevEgzVODFNxkGLF3a5G114r82AgQFuaYFFpEm\nqcCLb3r3XXA6ecXsPMQcv/mNqxdfW2t2JiJeSwVefNPSpTB3brMPXYmfGjMGLrgAXnzR7ExEvJae\ngxff8/nnrjnfd+3CcuGFeMMz2i2P8aZcfDfGMAx45RW47z74z3+gi/oq4n80kp10Po8/Dv/7v64e\nnHRe48fDuefCyy+bnYmIV1IPXnzL/v0wYADs2AFhYV4213tLYrwpF9+NcZ8L1q+HBx+ErVtBEwuJ\nn9FIdtK5/OlPcPPNEBZmdiZimoAGswRuA+7r0oUN9SK6dw/myJGKDs9MxJuoBy++4/hxiIyEt96C\n6GgA9eAVwxRyuYvlXMP7dT9zxeh8Ib6u3b+Dnz59OmFhYcTExLjbMjIysFqtxMbGEhsby4YN3392\nXrx4MTabjaioKDZu3OhuLywsJCYmBpvNxuzZs93tlZWVTJkyBZvNRnx8PLt37271zoif+/vf4aqr\n3MVdBOAf3EovvuZa3jM7FRGv0myBv/POO8nLy2vQZrFY+OUvf8nWrVvZunUrY8eOBaCoqIhVq1ZR\nVFREXl4es2bNcn/6SE9PJysrC7vdjt1ud28zKyuL0NBQ7HY7c+fOZf78+Z7eR/EHhgF//CPM07C0\n0lAt5/AI87mfRWanIuJVmi3w1157LcHBwY3aT3fZYN26daSkpBAYGEhkZCT9+/enoKCAsrIynE4n\ncXFxAKSmprJ27VoA1q9fT1paGgCTJk1i06ZNbdoh8R9BQSFYLBYsFgvju3Rha1ERlhtucLdZdFOV\n1MlhKoP4L1dRaHYqIl6j1TfZLV++nJUrVzJs2DCWLl3KRRddxN69e4mPj3fHWK1WHA4HgYGBWK1W\nd3tERAQOhwMAh8NBnz59XMkEBNCjRw8qKioICQlp9J4ZGRnu1wkJCSQkJLQ2ffEBTudBTn7XOp+R\n/J50IOWUKBV5gSrO5Q/czX0s5lY0hK34pvz8fPLz8z22vVYV+PT0dB588EEAFixYwLx588jKyvJY\nUk2pX+C1F3U7AAAZoklEQVSl84jn31gp5R/canYq4sWe4afczyKi2MFnZicj0gqndlwXLlzYpu21\naqCbXr16uS+Rzpw5k82bNwOunnlJSYk7rrS0FKvVSkREBKWlpY3aT66zZ88eAKqrqzl8+PBpe+/S\nec3nEf7A3dToqU45g++4kCf4BfN5xOxURLxCqwp8WVmZ+/XLL7/svsM+OTmZ3NxcqqqqKC4uxm63\nExcXR3h4OEFBQRQUFGAYBjk5OUyYMMG9TnZ2NgBr1qwhMTGxrfskfiSKHcTzIc9yp9mpiA94ip9z\nI//kMrMTEfECzXaJUlJSeOedd9i/fz99+vRh4cKF5Ofns23bNiwWC3379uXpp58GIDo6msmTJxMd\nHU1AQACZmZnuG6EyMzOZNm0ax44dY9y4cYwZMwaAGTNmMHXqVGw2G6GhoeTm5rbj7oqvuYdHeZL/\n4zjnm52K+IBDBPMMP+Vu9eJFNNCNeC+rxcJ/CMGGnYM09bWNdw/C4t25+GdML8rZQTghZWUQHt7M\ntkS8lyabEb81B8gm7QzFXaSxrwnjeXBNSiTSiakHL97p4EEqQkIYwh5K6XOGQO/pObYsxpty8d+Y\nS7GwOyQEvvgCTjOOh4gvUA9e/FNmJuuhmeIucnp7AG68EZ580uxUREyjHrx4n2PHoG9fosvL2eFF\nvUL14H0rxtixA0aOhK++gm7dmokX8T7qwYv/ee45iItjh9l5iG+LioIf/xieecbsTERMoR68eJfq\narjiCli5Ess11+BtvUL14H0nxjAM+PhjSE6GL7+Ec89tZh0R76IevPiXF1+E3r3h6qvNzkT8wVVX\nQUwMrFxpdiYiHU49ePEehgE/+AEsXAg33lg3SJJ39QrVg/eVmECgGoCrgZXAFe4Wl+7dgzlypKKZ\n7YiYRz148R9vvAGVlTB+vNmZiM+rxvUhwOADDL7kBqbxF3cbGHWzFYr4L/XgxXskJkJqKqSlAagH\nrxiPxcTzb3K5jQHspIpz3TE6p4g3Uw9e/MOWLWC3Q8qp872LtN2H/JDtDGQ6K8xORaTDqAcv3uHW\nW1031s2Z425SD14xnowZxke8zE305wsqOQ/14MXbtbXuqcCLKYKCQtzfgV4J5AOXA0cbRXp30Tj7\nGG/KpfPFrCOZN7mB5fwCFXjxdirw4pPq987XcyNvcx2P8ctTozC7IHg+xpty6XwxQ9nKq4ynP19w\njAt1ThGvpu/gxaf9mHwGsp2n+LnZqUgnsI1Y/s0P+V/+bHYqIu1OPXgxhcViwUINBYxgKfNYxW2n\ni8LsHp/nY7wpl84ZM4hPeYMk+lPOtzqniBdTD1581hRWAbCaySZnIp3Jf4nhHX6sa0bi99SDF1Oc\na7HwGZFM4zne5cdNRHlHj8+zMd6US+eNiWIH7xBNr8OHISiomW2JmEM9ePFJPwc+JeYMxV2k/XzG\nlWwEWL7c7FRE2k2zBX769OmEhYURExPjbquoqCApKYkBAwYwatQoDh065P7Z4sWLsdlsREVFsXHj\nRnd7YWEhMTEx2Gw2Zs+e7W6vrKxkypQp2Gw24uPj2b17t6f2TbzVwYPcC9zLErMzkU7stwCPPw6H\nD5udiki7aLbA33nnneTl5TVoW7JkCUlJSezcuZPExESWLHGdqIuKili1ahVFRUXk5eUxa9Ys9+WF\n9PR0srKysNvt2O129zazsrIIDQ3Fbrczd+5c5s+f7+l9FG/z8MO8DOwg2uxMpBOzA/zP/8Bjj5md\niki7aLbAX3vttQQHBzdoW79+PWl144WnpaWxdu1aANatW0dKSgqBgYFERkbSv39/CgoKKCsrw+l0\nEhcXB0Bqaqp7nfrbmjRpEps2bfLc3on32bULnn2W35idhwjAggXw5JNQoVnlxP8EtGal8vJywsLC\nAAgLC6O8vByAvXv3Eh8f746zWq04HA4CAwOxWq3u9oiICBwOBwAOh4M+ffq4kgkIoEePHlRUVBAS\nEtLofTMyMtyvExISSEhIaE36YqZf/xruuovyhQvNzkQELr8cbroJ/vhH+N3vzM5GOrn8/Hzy8/M9\ntr1WFfj6LBZL3ahk7a9+gRcftGULvP02PP20a853EW/wwANw1VWueRAuvtjsbKQTO7XjurCN58lW\n3UUfFhbGvn37ACgrK6NXr16Aq2deUlLijistLcVqtRIREUFpaWmj9pPr7NmzB4Dq6moOHz582t67\n+DjDgHvugYwM6NbN7GxEvnfZZTBlCixaZHYmIh7VqgKfnJxMdnY2ANnZ2UycONHdnpubS1VVFcXF\nxdjtduLi4ggPDycoKIiCggIMwyAnJ4cJEyY02taaNWtITEz0xH6Jt3ntNSgvh+nTzc5EpE6A+wpk\nrz/9iX2PPcaIur+fXIKC1NkQH2Y047bbbjN69+5tBAYGGlar1VixYoVx4MABIzEx0bDZbEZSUpJx\n8OBBd/zDDz9s9OvXz7jiiiuMvLw8d/uWLVuMQYMGGf369TPuuusud/vx48eNW2+91ejfv78xYsQI\no7i4+LR5tCBV8VYnThhGdLRhrF/vbgIMV7f+TIs/xnhTLoqp//dbWWUUEWWcy7EGMSJmaev/P41k\nJ+3vr3+Fv/3N9f173f0a/jnXe0tivCkXxZwas5pbKaYv8/m9O0bnHTGLposVr1N/rveLgP8CE4Et\njSK9+2TfPjHelItiTo25mG/4hMHcxMsUEI8KvJhJQ9WK13EVdwOoJYubWMMv2ELdFVH3IuJ99tOT\nu1jOc0zjPI6ZnY5Im6jAS7uZRSaXsZtfuS93ini/F7mF/zCEh1hgdioibaJL9OJxFouFIWzlDZL4\nEf/iC2yni8LbL9e2T4w35aKYpmJC2c+nxDCJffxL5x0xiS7Ri9e5EFjFFGazrIniLuLdDnAxP+cp\nngP47juTsxFpHRV48bhM4H2u4QV+YnYqIq32Mje7bgx94AGzUxFpFRV48azsbIYBv+AJszMRabO7\nAHJz4f33zU5F5KypwIvnfPYZ3H03k4HvuNDsbETarALgT3+CO+/UpXrxOSrw4hnHjrnG8374Ybab\nnYuIJ02YACNGwP33m52JyFnRXfTiGbNmwYEDkJuLpUsXvOFOaO+M8aZcFNOSGMMwXPPFx8TA88/D\nj3/czDoinqG76MV8L74Ir78Of/mLeyhaEb8SEuIacvmOO1yTJon4APXgpW2Ki12XL199FYYPBzrz\nOPMtifGmXBTTfEwgUO3+WwaQANxQr7V792COHKloZjsiZ089eDFPZSXcdhvce6+7uIv4l2rqD7G8\nkBq+ZRyPMNfddnLeBRFvox68nJX6E8k8h2tQm1tPG+lNvTBvivGmXBTTmpiLOMgWhvEAvyOXFDQh\njbQX9eClQ52cSOZufk8MsaTxLQ0nkdGJTvzbIYK5mZd4gl8wiE/NTkekSQFmJyC+ZzyvMIfHiedD\nPe8undInDGEOj/MSN6Mvp8RbqQcvZ2UgsILpTOJFSuljdjoipnme29nAWHIAamvNTkekERV4ablv\nvmE98Ev+SAHxZmcjYrq7+QMXAfzud2anItJImwp8ZGQkgwcPJjY2lri4OAAqKipISkpiwIABjBo1\nikOHDrnjFy9ejM1mIyoqio0bN7rbCwsLiYmJwWazMXv27LakJO2lqgpuuYVVwN+5w+xsRLzCCboy\nGVxjQLz2mtnpiDTQpgJvsVjIz89n69atbN68GYAlS5aQlJTEzp07SUxMZMmSJQAUFRWxatUqioqK\nyMvLY9asWe67A9PT08nKysJut2O328nLy2vjbolHGQb8/OcQHMyvzc5FxMvsA1i1yjVe/Zdfmp2O\niFubL9Gfegv/+vXrSUtLAyAtLY21a9cCsG7dOlJSUggMDCQyMpL+/ftTUFBAWVkZTqfTfQUgNTXV\nvY54iSeegIICyMnRPfIip3P11bBgAdx8syalEa/R5h78DTfcwLBhw3jmmWcAKC8vJywsDICwsDDK\n64Z13Lt3L1ar1b2u1WrF4XA0ao+IiMDhcLQlLfGkvDxYsgTWr4fu3c3ORsR7/fznMHQopKTAiRNm\nZyPStsfkPvjgA3r37s0333xDUlISUVFRDX5usVjqhi31jIyMDPfrhIQEEhISPLZtOY3PPoPUVHjp\nJYiMNDsbES8V4D7PBQIvAc6uXbkDqEVD2UrL5efnk5+f77HttanA9+7dG4CePXty0003sXnzZsLC\nwti3bx/h4eGUlZXRq1cvwNUzLykpca9bWlqK1WolIiKC0tLSBu0RERGnfb/6BV7a2c6dMG6cq/d+\nzTVmZyPixU4OZwsngFs4zquM5y9E8lOewek8x9TsxHec2nFduHBhm7bX6kv03333HU6nE4CjR4+y\nceNGYmJiSE5OJjs7G4Ds7GwmTpwIQHJyMrm5uVRVVVFcXIzdbicuLo7w8HCCgoIoKCjAMAxycnLc\n60jHCgoKwWKxMNJiofyKK/hpcTGWGTPcV2I8eTVGxF9Vch4TWEcUn/E4c8xORzqxVvfgy8vLuemm\nmwCorq7m9ttvZ9SoUQwbNozJkyeTlZVFZGQkq1evBiA6OprJkycTHR1NQEAAmZmZ7oKRmZnJtGnT\nOHbsGOPGjWPMmDEe2DU5W07nQe5gJUuZx+38nTdJOk2UirxIc47SjXG8xltczyJwPYmiD8jSwTTZ\njLgYBhldupBGJP/DKxQxsIlA758MxLtjvCkXxbR3TCj7yacngx56CB54oJntiDSkyWak7SorITWV\nsUA8H56huIvI2TjAxdwAsHIlPPaY2elIJ6PJZjq7AwfgppugVy+uA44RZnZGIn6lHGDTJhg5Ei68\nEH72M7NTkk5CPfjOzG6HH/7QtaxezTGz8xHxV336wJtvwm9/C3/7m9nZSCeh7+A7q/feg1tvdZ1w\n6noUrpsezf/e0r9jvCkXxXRMTCCuR+ngSmAT8Fvgz/Ui9Ky8nI6+g5ezU1Xlutnnlltc3wvqcqFI\nOzv5nLzBDgyuYwczuYrXGMMllAIGTudBk3MUf6QC30kEBYUwxGJh27nn8s+HH6b3119jGT1az7iL\ndLDPiSKeD/mQeLYSSwrPm52S+Cldou8Mqqv5dWAgs+nJr/g92aRx+ufZvemypr/GeFMuijH73/wq\nCllJKjso4pavv4aePZt5L+lMdIlezmzHDvjRj/gx8AMKyWYaGqxGxDt8zA/4AYUUAwwZ4prUScRD\nVOD9VU0NLF3qejRn+nRGA6X0MTsrETlFJefxK3DNKT93rmte+cOHzU5L/IAu0fujL75wnSS6dIFn\nn4XLL9cd8l4T4025KMZ7/s1dd9pfCDwKjAfuBtbUW1N32nc+ba17KvB+ICgoBKfzIEOAe4CxwEPA\nMk49rXjLCa8zx3hTLorx1n/z69nEIu7nQo7yWx5kDbdgcI7OgZ2MvoPv7AyDeOdBNnIDr3IJ/+ER\nLucgj2Ng1D2a0/zJRUS8yVskEs+H3M0f+CV/5FNimAyur95EWkgF3ledOOEaESs2lj8Cf+d2+lLM\no/yKw1xkdnYi0mYWXmcMP+TfzGMpcwEGD3Z9V69CLy2gS/S+xumEZ56Bxx+Hfv3gnnuwjB+P91yO\nVIwu0Sumvd7HyMuDjAzXTXgLFrhGowzQlCL+St/B+7mgoBCOOg9yHZAKJAOv47oRp7BBpLeczBSj\nAq+Y9nmf74e8TQIWAFcA/wBeAP4FdNONeH5F38H7s+3b+bXzILuJ4BGuopDHGUA5t2FQqO/XRTqZ\n74e8fQODkRj8CDtlPMTTRFPMpfzaeRC2bYNO2BmSxtSD9zZffw0vvOAaJ768nCUOBzn8t5k52r2p\nt6IY9eAV0/HvYxDDp6QwhPsiI+G88yAlxbXYbM1sV7yVLtH7sKDuwYR+e4jhwDAgDhgKrAdWAm8D\ntYD3nKgU4xsne8V4V0xH5uK6jB8PpACTge+A94B365bybhdxRJPb+AQVeF+ydy9s2QIffQQffcT+\n11+nkkv4iOFsYRgfMZwPuJqjdKu3kjedqBTjWyd7xXhHjHm5WKglis8Yybtcy3uM5F0CKSX81ltd\no1xeey3ExLgGxRKv4zcFPi8vjzlz5lBTU8PMmTOZP39+g5/7TIGvrobiYvj8c/jsM9efJ1/X1MDw\n4e7lkgkTKPOak5Bi/P1krxizYrwpF4PLCGQkNYwErgXCADtQDOyq+/Pr87vx4scfwWWXwfnnn3ZL\n+fn5JCQkNJOPtIVfFPiamhquuOIK3nzzTSIiIhg+fDgvvPACV155pTvG1AJvGK7H0w4ccC0VFd+/\nPvn3khJXEf/qK+jdG664gsy38vn0RCWfAZ8DZaffeDNv7k0nB8X418leMfo3h1D2058viGQXkeyi\nL8VE8hdG22ywZw8EB0NkJFxyievcVrdkvPUWGXff7Wq/+GJdBWgHba17XvEA5ebNm+nfvz+RkZEA\n3Hbbbaxbt65BgW+SYUBtravnfOpSVQXffedajh37/nX95ehR1zOlR46c+c+uXSE01LWEhPDiO+9S\nXn2CA0AF4AA+w/VJ+PiuXbBr18kEz5C8pdXHTETEEw5wMQe4mALi67WuALsdC9B73z4i9+2jN3AJ\n0Ltu+crSBT75BMrKXOfInj1d58fgYAgJcS0nX5/8s0cPCAyEc85xPb9/uj8DA103CZ5//vd/6ln/\nVvGKo+ZwOOjT5/uZzqxWKwUFBY3iXBOmQH/gP7iS74rrRrQTuB4iqb+cwHWDScPFwncY7r8fBQ4D\nR+othzmHI9TUvXYtVSdOuD4M7NlTL6OWfFIWEfE1rkfyDGBv3dKIcQ45n34KuG7tC9u7l5C9ewkG\nQuqWYLoQQm3da7gI13n7nHp/du1yDlcNjnF9hVld7Rql8/hx13LsmGuxWBoW/K5dXR8GunQ582Kx\nuJamXtdf4Myvz6T+z1vyuin33gtXX918XAt5RYG3tGTH6/kCuLDV79aSyx0tHQayJXk3F+OJbSjG\nt2K8KRfFdEyMN+XiyRiXE0Bp3dJQbfMr19a4nt1vzrffuhZ/9sorHt2cVxT4iIgISkpK3H8vKSnB\narU2iPGCWwVERER8hlfcFTFs2DDsdju7du2iqqqKVatWkZycbHZaIiIiPssrevABAQE8+eSTjB49\nmpqaGmbMmNGyG+xERETktLyiBw8wduxYPv/8cxYvXszf/vY3zjnnHD7++GP3z3ft2sX5559PbGws\nsbGxzJo1y/2zwsJCYmJisNlszJ4924z0fco//vEPBg4c2OgYAyxevBibzUZUVBQbN250t+sYt01G\nRgZWq9X9/3fDhg3unzV1zOXs5eXlERUVhc1m45FHHjE7Hb8RGRnJ4MGDiY2NJS4uDoCKigqSkpIY\nMGAAo0aN4tChQyZn6XumT59OWFgYMTEx7rYzHdezPlcYXmbHjh3G559/biQkJBiFhYXu9uLiYmPQ\noEGnXWf48OFGQUGBYRiGMXbsWGPDhg0dkquvauoYb9++3RgyZIhRVVVlFBcXG/369TNqa2sNw9Ax\nbquMjAxj6dKljdpPd8xrampMyND3VVdXG/369TOKi4uNqqoqY8iQIUZRUZHZafmFyMhI48CBAw3a\n7rnnHuORRx4xDMMwlixZYsyfP9+M1Hzau+++a3z88ccNaltTx7U15wqv6cGfFBUVxYABA1ocX1ZW\nhtPpdH+qTE1NZe3ate2Vnl9o6hivW7eOlJQUAgMDiYyMpH///hQUFOgYe4hxmhtFT3fMN2/ebEJ2\nvq/+eBqBgYHu8TTEM079/7t+/XrS0tIASEtL0zmhFa699lqCg4MbtDV1XFtzrvC6An8mxcXFxMbG\nkpCQwPvvvw+4nqGvf8d9REQEDofDrBR92t69exscS6vVisPhaNSuY9w6y5cvZ8iQIcyYMcN92a2p\nYy5n73TjaehYeobFYuGGG25g2LBhPPPMMwCUl5cTFhYGQFhYGOXl5Wam6DeaOq6tOVeYcpNdUlIS\n+/bta9S+aNEibrzxxtOuc8kll1BSUkJwcDAff/wxEydOZPv27e2dqs9qzTGWtmnqmD/88MOkp6fz\n4IMPArBgwQLmzZtHVlbWabdztuNCiIuOW/v54IMP6N27N9988w1JSUlERUU1+LnFYtHxbwfNHdfm\njrkpBf6NN94463W6du1K165dAbjqqqvo168fdrudiIgISku/H16htLSUiIgIj+Xqq1pzjE8dj6C0\ntBSr1apj3EItPeYzZ850f8g63THXsW2dloynIa3Tu3dvAHr27MlNN93E5s2bCQsLY9++fYSHh1NW\nVkavXr1MztI/NHVcW3Ou8OpL9PW/89m/fz81Na4R5r766ivsdjuXX345vXv3JigoiIKCAgzDICcn\nh4kTJ5qVss+pf4yTk5PJzc2lqqqK4uJi7HY7cXFxhIeH6xi3UVnZ91MNvfzyy+67Zps65nL2NJ5G\n+/juu+9wOp0AHD16lI0bNxITE0NycjLZ2dkAZGdn65zgIU0d11adK9rjzsC2eOmllwyr1Wqcd955\nRlhYmDFmzBjDMAxjzZo1xsCBA42hQ4caV111lfHKK6+419myZYsxaNAgo1+/fsZdd91lVuo+o6lj\nbBiG8fDDDxv9+vUzrrjiCiMvL8/drmPcNlOnTjViYmKMwYMHGxMmTDD27dvn/llTx1zO3muvvWYM\nGDDA6Nevn7Fo0SKz0/ELX331lTFkyBBjyJAhxsCBA93H9cCBA0ZiYqJhs9mMpKQk4+DBgyZn6ntu\nu+02o3fv3kZgYKBhtVqNFStWnPG4nu25wiumixURERHP8upL9CIiItI6KvAiIiJ+SAVeRETED6nA\ni4iI+CEVeBERET+kAi8iIuKH/j+DTO63HeU90wAAAABJRU5ErkJggg==\n", "text": [ "" ] } ], "prompt_number": 18 }, { "cell_type": "code", "collapsed": false, "input": [], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 18 }, { "cell_type": "code", "collapsed": false, "input": [], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 18 }, { "cell_type": "code", "collapsed": false, "input": [], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 18 }, { "cell_type": "code", "collapsed": false, "input": [], "language": "python", "metadata": {}, "outputs": [], "prompt_number": 18 }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Zadania\n", "\n", "### ElementwiseKernel\n", "\n", "Napisa\u0107 program ca\u0142kuj\u0105cy r\u00f3wnanie stochastyczne, wykorzystuj\u0105ce `ElementwiseKernel`" ] }, { "cell_type": "code", "collapsed": false, "input": [ "from pycuda.elementwise import ElementwiseKernel\n", "SDE2 = ElementwiseKernel(\n", " \"float *x,unsigned int *rng_state\",\n", " \"x[i] = x[i]+0.01*f(x[i]) + sqrtf(0.01*0.1)*rng_uni( &(rng_state[i]) ) \",\n", " \"SDE\",\n", " preamble=\"\"\"\n", " __device__ float f(float x)\n", " { \n", " return sin(x)+0.1;\n", " }\n", " \"\"\"+rng_src)\n" ], "language": "python", "metadata": {}, "outputs": [] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Wydajno\u015b\u0107\n", "\n", "Zbada\u0107 jak zmienia si\u0119 wydajno\u015b\u0107 programu w zale\u017cno\u015bci od wielko\u015bci bloku" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### R\u00f3wnanie dyfuzji 2d\n", "\n", "Napisa\u0107 program rozwi\u0105zuj\u0105cy r\u00f3wnanie dyfuzji 2d metod\u0105 symulacji odpowiedniego SDE." ] }, { "cell_type": "code", "collapsed": false, "input": [], "language": "python", "metadata": {}, "outputs": [] } ], "metadata": {} } ] }