{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "daa26d0c",
   "metadata": {},
   "source": [
    "# Porting walkthrough: Brunel, side by side\n",
    "\n",
    "**What you'll learn / who it's for (porting / simulation).** A NEST network ported to `brainpy.state`, step for step. We take the Brunel (2000) random balanced network with **delta synapses** and show the NEST (PyNEST) construction next to the `brainpy.state` `Simulator` port — with a parameter-mapping table and the semantic divergences a port must account for.\n",
    "\n",
    "The NEST code blocks below are shown for **reference only** (they need a live NEST install); the `brainpy.state` cells are runnable. Source: `examples/nest_like/brunel_delta.py` and NEST's `brunel_delta_nest.py`."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "id": "2ba1260e",
   "metadata": {
    "execution": {
     "iopub.execute_input": "2026-06-17T09:12:03.670630Z",
     "iopub.status.busy": "2026-06-17T09:12:03.670484Z",
     "iopub.status.idle": "2026-06-17T09:12:07.527351Z",
     "shell.execute_reply": "2026-06-17T09:12:07.526751Z"
    }
   },
   "outputs": [
    {
     "name": "stderr",
     "output_type": "stream",
     "text": [
      "An NVIDIA GPU may be present on this machine, but a CUDA-enabled jaxlib is not installed. Falling back to cpu.\n"
     ]
    }
   ],
   "source": [
    "import jax\n",
    "jax.config.update('jax_enable_x64', True)\n",
    "import brainstate\n",
    "brainstate.environ.set(precision=64)\n",
    "\n",
    "import numpy as np\n",
    "import brainunit as u\n",
    "import braintools\n",
    "import matplotlib.pyplot as plt\n",
    "from brainpy import state as bp"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "f40fff73",
   "metadata": {},
   "source": [
    "## Shared parameters\n",
    "\n",
    "Both simulators use the identical Brunel knobs. With delta synapses the postsynaptic amplitude is a membrane-voltage jump, so `J_ex = J` directly — there is no current-to-voltage (`ComputePSPnorm`) rescaling that the alpha variant needs."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "id": "8d5a52f6",
   "metadata": {
    "execution": {
     "iopub.execute_input": "2026-06-17T09:12:07.529762Z",
     "iopub.status.busy": "2026-06-17T09:12:07.529319Z",
     "iopub.status.idle": "2026-06-17T09:12:07.533525Z",
     "shell.execute_reply": "2026-06-17T09:12:07.532908Z"
    }
   },
   "outputs": [],
   "source": [
    "order = 100                          # small for a fast notebook run\n",
    "g, eta, epsilon = 5.0, 2.0, 0.1\n",
    "delay = 1.5                          # ms\n",
    "NE, NI = 4 * order, 1 * order\n",
    "CE, CI = int(epsilon * NE), int(epsilon * NI)\n",
    "N_rec = 50\n",
    "tauMem, CMem, theta, tref = 20.0, 1.0, 20.0, 2.0\n",
    "J = 0.1                              # mV (delta jump)\n",
    "J_ex = J\n",
    "J_in = -g * J_ex\n",
    "nu_th = theta / (J * CE * tauMem)\n",
    "p_rate = 1000.0 * (eta * nu_th) * CE"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "fd20ba61",
   "metadata": {},
   "source": [
    "## Step 1 — create populations and devices\n",
    "\n",
    "**NEST (reference):**\n",
    "\n",
    "```python\n",
    "import nest\n",
    "nest.SetKernelStatus({'resolution': 0.1})\n",
    "nest.SetDefaults('iaf_psc_delta', {\n",
    "    'C_m': CMem, 'tau_m': tauMem, 't_ref': tref,\n",
    "    'E_L': 0.0, 'V_reset': 0.0, 'V_th': theta})\n",
    "nodes_ex = nest.Create('iaf_psc_delta', NE)\n",
    "nodes_in = nest.Create('iaf_psc_delta', NI)\n",
    "noise = nest.Create('poisson_generator', params={'rate': p_rate})\n",
    "espikes = nest.Create('spike_recorder')\n",
    "ispikes = nest.Create('spike_recorder')\n",
    "```\n",
    "\n",
    "**brainpy.state (runnable):** there is no global kernel — a `Simulator` owns the populations, devices, and connections, and parameters carry explicit units."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "id": "951b6204",
   "metadata": {
    "execution": {
     "iopub.execute_input": "2026-06-17T09:12:07.535207Z",
     "iopub.status.busy": "2026-06-17T09:12:07.535020Z",
     "iopub.status.idle": "2026-06-17T09:12:10.298196Z",
     "shell.execute_reply": "2026-06-17T09:12:10.297210Z"
    }
   },
   "outputs": [],
   "source": [
    "npar = dict(C_m=CMem * u.pF, tau_m=tauMem * u.ms, t_ref=tref * u.ms,\n",
    "            E_L=0. * u.mV, V_reset=0. * u.mV, V_th=theta * u.mV,\n",
    "            V_initializer=braintools.init.Constant(0. * u.mV))\n",
    "\n",
    "sim = bp.Simulator(dt=0.1 * u.ms)\n",
    "ne = sim.create(bp.iaf_psc_delta, NE, params=npar)\n",
    "ni = sim.create(bp.iaf_psc_delta, NI, params=npar)\n",
    "noise = sim.create(bp.poisson_generator, rate=p_rate * u.Hz)\n",
    "esr = sim.create(bp.spike_recorder)\n",
    "isr = sim.create(bp.spike_recorder)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "f90559b6",
   "metadata": {},
   "source": [
    "## Step 2 — connect\n",
    "\n",
    "**NEST (reference):**\n",
    "\n",
    "```python\n",
    "nest.Connect(noise, nodes_ex + nodes_in,\n",
    "             syn_spec={'weight': J_ex, 'delay': delay})\n",
    "conn_ex = {'rule': 'fixed_indegree', 'indegree': CE}\n",
    "conn_in = {'rule': 'fixed_indegree', 'indegree': CI}\n",
    "nest.Connect(nodes_ex, nodes_ex + nodes_in, conn_ex,\n",
    "             {'weight': J_ex, 'delay': delay})\n",
    "nest.Connect(nodes_in, nodes_ex + nodes_in, conn_in,\n",
    "             {'weight': J_in, 'delay': delay})\n",
    "nest.Connect(nodes_ex[:N_rec], espikes)\n",
    "nest.Connect(nodes_in[:N_rec], ispikes)\n",
    "```\n",
    "\n",
    "**brainpy.state (runnable):** `connect(src, dst, weight=, delay=, rule=)`, with `ne + ni` population algebra and `fixed_indegree(K)` — the same shape as `nest.Connect`. The delta weights are voltages in `mV`."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "id": "320913b5",
   "metadata": {
    "execution": {
     "iopub.execute_input": "2026-06-17T09:12:10.300535Z",
     "iopub.status.busy": "2026-06-17T09:12:10.300106Z",
     "iopub.status.idle": "2026-06-17T09:12:14.124711Z",
     "shell.execute_reply": "2026-06-17T09:12:14.123742Z"
    }
   },
   "outputs": [],
   "source": [
    "sim.connect(noise, ne, weight=J_ex * u.mV, delay=delay * u.ms, rule=bp.all_to_all)\n",
    "sim.connect(noise, ni, weight=J_ex * u.mV, delay=delay * u.ms, rule=bp.all_to_all)\n",
    "sim.connect(ne, ne + ni, weight=J_ex * u.mV, delay=delay * u.ms,\n",
    "            rule=bp.fixed_indegree(CE), comm='sparse', allow_multapses=True, seed=1)\n",
    "sim.connect(ni, ne + ni, weight=J_in * u.mV, delay=delay * u.ms,\n",
    "            rule=bp.fixed_indegree(CI), comm='sparse', allow_multapses=True, seed=2)\n",
    "sim.connect(ne[:N_rec], esr)\n",
    "sim.connect(ni[:N_rec], isr)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "5d9d0ef9",
   "metadata": {},
   "source": [
    "## Step 3 — simulate and read back\n",
    "\n",
    "**NEST (reference):** `nest.Simulate(1000.0)`, then `nest.GetStatus(espikes, 'n_events')` → rate.\n",
    "\n",
    "**brainpy.state (runnable):** `sim.simulate(T)` returns a result; read the rate of a combined-population recorder via `esr.segments[0].population`."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 5,
   "id": "42be3aa7",
   "metadata": {
    "execution": {
     "iopub.execute_input": "2026-06-17T09:12:14.126980Z",
     "iopub.status.busy": "2026-06-17T09:12:14.126589Z",
     "iopub.status.idle": "2026-06-17T09:12:21.046792Z",
     "shell.execute_reply": "2026-06-17T09:12:21.045897Z"
    }
   },
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "excitatory rate: 60.56 spks/s\n",
      "inhibitory rate: 60.18 spks/s\n"
     ]
    },
    {
     "data": {
      "image/png": "iVBORw0KGgoAAAANSUhEUgAAAxYAAAGGCAYAAADmRxfNAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjgsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvwVt1zgAAAAlwSFlzAAAPYQAAD2EBqD+naQAAap9JREFUeJzt3XlcVdX6P/APgwwyigqICSKZs0E4ZmomXiXTyqnMQs3b5FBmpY0OWVrZoHkzb5NYOZSlddPMzNTUzJxzSHMqTcUhFTQVEdbvj76cH0dAz96cddZah8/79eJ1b4dnP/vZa6+9YXvOevARQggQERERERGVga/qAoiIiIiIyHx8sCAiIiIiojLjgwUREREREZUZHyyIiIiIiKjM+GBBRERERERlxgcLIiIiIiIqMz5YEBERERFRmfHBgoiIiIiIyowPFkREREREVGZ8sCAqJ2rWrIl+/fq5FHvgwAEEBQVh1apVV4y98cYbceONN0qvyVV33nknevXq5dacJli2bBl8fHywbNky1aWQBKNHj4aPjw+OHz/ukf3JuDaJyPvxwYLoCjIzM+Hj4+P0FR0djXbt2mHhwoWqy5Pi+eefR/PmzdGqVSuP7nf79u0YPXo0fv/9d9s5RowYgc8//xybN292X2GGmjlzJiZOnGh7+7Nnz2L06NHl+mGlrGNInjdlyhRkZmba3v7QoUMYPXo0Nm3a5LaaVODcJRX4YEHkoueffx4fffQRPvzwQwwfPhzHjh3DzTffjPnz56suza2OHTuG6dOn48EHH/T4vrdv344xY8aU6cEiJSUFTZo0wWuvvea+wgzQpk0bnDt3Dm3atHG85o4HizFjxvDBohz+crZz5068++67qsuwxR0PFmPGjOGDBZENfLAgclF6ejruvvtu3HPPPXj88cexYsUKVKhQAbNmzbrsdhcvXsSFCxc8VGXZffzxx/D390eXLl1Ul2Jbr169MHfuXJw5c0Z1KR7j6+uLoKAg+Prqf1v/+++/jcrrbgUFBTh//rzW+wsMDESFChUkVUQymXIdkHfS/ycQkaYiIyMRHBwMf39/x2u///47fHx88Oqrr2LixIlISkpCYGAgtm/f7vhI1aX/Gl/SZ+NvvPFGNGzYENu3b0e7du1QsWJFVK9eHa+88kqxOnJzczFq1ChcffXVCAwMRI0aNTB8+HDk5ubaOq4vvvgCzZs3R2hoaLHvvfPOO0hKSkJwcDCaNWuGFStWlJjDTk2ZmZno2bMnAKBdu3aOj50VjsuXX36Jzp07Iy4uDoGBgUhKSsLYsWORn59fLFeHDh3w999/Y/HixTZGwD0KCgowceJENGjQAEFBQYiJicEDDzyAkydPOmJGjRoFX19fLFmyxGnb+++/HwEBAU4f5zp48CAGDBjgOP7ExEQ89NBDjofWS+fRjTfeiAULFuCPP/5wjGXNmjUBABcuXMDIkSORmpqKiIgIhISEoHXr1li6dKljf7///juqVq0KABgzZowjx+jRox0x33//PVq3bo2QkBBERkbi1ltvxa+//up0LIVrA7Zv34677roLlSpVwg033IBp06bBx8cHGzduLDZ248aNg5+fHw4ePFjq+JaWFwB++eUX9OvXD7Vq1UJQUBBiY2Nx77334q+//nLKcfr0aQwdOhQ1a9ZEYGAgoqOj0aFDB2zYsOGKYwi4Ps99fHwwePBgzJgxAw0aNEBgYCC++eabUo/tco4fP45evXohPDwclStXxiOPPFLsoeFy+3v11Vdx/fXXo3LlyggODkZqaio+++yzYvu5dI1F4f1r1apVGDZsGKpWrYqQkBDcfvvtOHbsWLFtb7nlFnz77bdITk5GUFAQ6tevj7lz5zpi9u7dCx8fH7zxxhvF9v3jjz/Cx8enxH+0ycrKQv/+/XHVVVchMDAQ1apVw6233uq4r9asWRPbtm3D8uXLHeescA3YiRMn8Pjjj6NRo0YIDQ1FeHg40tPTna6zZcuWoWnTpgCA/v37O3IUfQdkzZo16NSpEyIiIlCxYkW0bdvWpfVohdfoJ598gqeffhqxsbEICQlB165dceDAgWLxc+bMQWpqKoKDg1GlShXcfffdxa6Jfv36ITQ0FHv27MHNN9+MsLAw9OnT54pzl0gW/yuHEBEAZGdn4/jx4xBC4OjRo5g8eTLOnDmDu+++u1jstGnTcP78edx///0IDAxEVFSU5f2dPHkSnTp1Qrdu3dCrVy989tlnGDFiBBo1aoT09HQA//zy2rVrV6xcuRL3338/6tWrhy1btuCNN97Ab7/9hi+++MLSPvPy8rB27Vo89NBDxb73/vvv44EHHsD111+PoUOHYu/evejatSuioqJQo0YNR5zdmtq0aYOHH34Yb775Jp5++mnUq1cPABz/m5mZidDQUAwbNgyhoaH4/vvvMXLkSOTk5GDChAlOuerXr4/g4GCsWrUKt99+u6UxcJcHHngAmZmZ6N+/Px5++GHs27cP//nPf7Bx40asWrUKFSpUwLPPPouvvvoKAwYMwJYtWxAWFoZFixbh3XffxdixY3HttdcC+OejGc2aNcOpU6dw//33o27dujh48CA+++wznD17FgEBAcX2/8wzzyA7Oxt//vmn45e3wofFnJwcvPfee+jduzfuu+8+nD59Gu+//z46duyIn3/+GcnJyahatSrefvttPPTQQ7j99tvRrVs3AEDjxo0BAN999x3S09NRq1YtjB49GufOncPkyZPRqlUrbNiwodgvMT179kTt2rUxbtw4CCHQo0cPDBo0CDNmzEBKSopT7IwZM3DjjTeievXqVxznS/MCwOLFi7F37170798fsbGx2LZtG9555x1s27YNP/30E3x8fAAADz74ID777DMMHjwY9evXx19//YWVK1fi119/xXXXXXfZMbQ6z7///nt8+umnGDx4MKpUqWL7l7xevXqhZs2aGD9+PH766Se8+eabOHnyJD788EOX9jdp0iR07doVffr0wYULFzB79mz07NkT8+fPR+fOna+4/yFDhqBSpUoYNWoUfv/9d0ycOBGDBw/GJ5984hS3a9cu3HHHHXjwwQfRt29fTJs2DT179sQ333yDDh06oFatWmjVqhVmzJiBRx991GnbGTNmICwsDLfeemux/Xfv3h3btm3DkCFDULNmTRw9ehSLFy/G/v37UbNmTUycOBFDhgxBaGgonnnmGQBATEwMgH8eZr744gv07NkTiYmJOHLkCP773/+ibdu22L59O+Li4lCvXj08//zzGDlyJO6//360bt0aAHD99dc7xjU9PR2pqamOfxiYNm0abrrpJqxYsQLNmjW74hi++OKL8PHxwYgRI3D06FFMnDgRaWlp2LRpE4KDgwHAce9o2rQpxo8fjyNHjmDSpElYtWoVNm7ciMjISEe+ixcvomPHjrjhhhvw6quvomLFioiNjS117hJJJYjosqZNmyYAFPsKDAwUmZmZTrH79u0TAER4eLg4evRoiXn27dvn9PrSpUsFALF06VLHa23bthUAxIcffuh4LTc3V8TGxoru3bs7Xvvoo4+Er6+vWLFihVPOqVOnCgBi1apVjtcSEhJE3759L3usu3fvFgDE5MmTnV6/cOGCiI6OFsnJySI3N9fx+jvvvCMAiLZt27qlpjlz5hQbi0Jnz54t9toDDzwgKlasKM6fP1/se9dcc41IT0+/7PHKsmLFCgFAzJgxw+n1b775ptjrW7ZsEQEBAeLf//63OHnypKhevbpo0qSJyMvLc8RkZGQIX19fsXbt2mL7KigoEEKUPI86d+4sEhISim1z8eJFp/MohBAnT54UMTEx4t5773W8duzYMQFAjBo1qliO5ORkER0dLf766y/Ha5s3bxa+vr4iIyPD8dqoUaMEANG7d+9iOXr37i3i4uJEfn6+47UNGzYIAGLatGnF4ou6XN6S5sqsWbMEAPHDDz84XouIiBCDBg267H5KG0Mr8xyA8PX1Fdu2bbvsvi6n8Hi7du3q9PrAgQMFALF582aX9nfp2Fy4cEE0bNhQ3HTTTU6vX3ptFt6/0tLSHHNOCCEeffRR4efnJ06dOuW0LQDx+eefO17Lzs4W1apVEykpKY7X/vvf/woA4tdff3Wqp0qVKiXeq06ePCkAiAkTJhT7XlENGjRwuicVOn/+vNNcE+Kfe3ZgYKB4/vnnHa+tXbu2xDlYUFAgateuLTp27Og0BmfPnhWJiYmiQ4cOl62r8BqtXr26yMnJcbz+6aefCgBi0qRJQoj/f79t2LChOHfunCNu/vz5AoAYOXKk47W+ffsKAOLJJ58str/S5i6RTPwoFJGL3nrrLSxevBiLFy/Gxx9/jHbt2uHf//6309v7hbp37+74GIldoaGhTu+GBAQEoFmzZti7d6/jtTlz5qBevXqoW7cujh8/7vi66aabAMDpoy2uKPyoSKVKlZxeX7duHY4ePYoHH3zQ6V/H+/Xrh4iICKdYd9dUqPBf8oB/PsJy/PhxtG7dGmfPnsWOHTuKxVeqVMljrTkvNWfOHERERKBDhw5OY5CamorQ0FCnMWjYsCHGjBmD9957Dx07dsTx48cxffp0x0fsCgoK8MUXX6BLly5o0qRJsX0V/uu7FX5+fo7zWFBQgBMnTuDixYto0qSJ42NAl3P48GFs2rQJ/fr1c3o3rnHjxujQoQO+/vrrYtuU1AwgIyMDhw4dchqPGTNmIDg4GN27d3fpWErKW3SunD9/HsePH0eLFi0AwOn4IiMjsWbNGhw6dMilfRVldZ63bdsW9evXt7yfSw0aNMjpv4cMGQIAxca8tP0VHZuTJ08iOzsbrVu3dum8A/98TK/onGvdujXy8/Pxxx9/OMXFxcU5vVsYHh6OjIwMbNy4EVlZWQD+efclKCgIM2bMcMQtWrQIx48fL/Gd4ODgYAQEBGDZsmVOHyl0VWBgoGMNUn5+Pv766y+EhoaiTp06Lh3/pk2bsGvXLtx1113466+/HOf877//Rvv27fHDDz+goKDginkyMjIQFhbm+O8ePXqgWrVqjnNYeL8dOHAggoKCHHGdO3dG3bp1sWDBgmI5S3qXmUgFfhSKyEXNmjVz+sWud+/eSElJweDBg3HLLbc4/cKdmJhY5v1dddVVxX5prFSpEn755RfHf+/atQu//vprqQ8xR48etbVv8X8fKSlU+EtD7dq1nV6vUKECatWq5fSarJq2bduGZ599Ft9//z1ycnKcvpednV0sXghxxV+6T5w4YXthfVRUVIkfQQL+GYPs7GxER0eX+P1Lx+CJJ57A7Nmz8fPPP2PcuHFOvxAeO3YMOTk5aNiwoa06SzN9+nS89tpr2LFjB/Ly8hyvuzJ3C+dDnTp1in2vXr16WLRoEf7++2+EhIRcNm+HDh1QrVo1zJgxA+3bt0dBQQFmzZqFW2+91ekXr8spKe+JEycwZswYzJ49u9hYF50rr7zyCvr27YsaNWogNTUVN998MzIyMorN6ZJYnefuuCcAxa/BpKQk+Pr6Flu7Vdr+5s+fjxdeeAGbNm1yWgvi6gNqfHy8038X/iPEpb/oX3311cVyXnPNNQD+Wb8TGxuLyMhIdOnSBTNnzsTYsWMB/PNgWb16dccDWlGBgYF4+eWX8dhjjyEmJgYtWrTALbfcgoyMDMTGxl6x9oKCAkyaNAlTpkzBvn37nNZnVa5c+Yrb79q1CwDQt2/fUmOys7OL/cPMpS49hz4+Prj66qsd5/By11fdunWxcuVKp9f8/f1x1VVXXbF+Ik/ggwWRTb6+vmjXrh0mTZqEXbt2oUGDBo7vFf1XwUKl/eAuafEx8M+/Kpek6C/9BQUFaNSoEV5//fUSY4uufXBF4Q9XO/8aKKsmADh16hTatm2L8PBwPP/880hKSkJQUBA2bNiAESNGlPivhCdPniz2A/xS3bp1w/Llyy3XA/zzL9Kl/WHAgoICREdHO/1LbFGX/jK6d+9exy8tW7ZssVWPFR9//DH69euH2267DU888QSio6Ph5+eH8ePHY8+ePVL2WdI14efnh7vuugvvvvsupkyZglWrVuHQoUMl/mu1lby9evXCjz/+iCeeeALJyckIDQ1FQUEBOnXq5DRXevXqhdatW2PevHn49ttvMWHCBLz88suYO3euYx1TaazO85LqdIfS7isl7W/FihXo2rUr2rRpgylTpqBatWqoUKECpk2bhpkzZ7q0P1fuS1ZkZGRgzpw5+PHHH9GoUSP873//w8CBA0vtbjZ06FB06dIFX3zxBRYtWoTnnnsO48ePx/fff19src6lxo0bh+eeew733nsvxo4di6ioKPj6+mLo0KEuvdNQGDNhwgQkJyeXGKNiHUPRd2KIVOODBVEZXLx4EQBcamta+K9Yp06dcnr90o8QWJGUlITNmzejffv2tj4Sc6n4+HgEBwdj3759Tq8nJCQA+Odf7Ir+S2JeXh727dvnWGRc1ppKi1+2bBn++usvzJ071+nvNFxaZ6GLFy/iwIED6Nq162X399prr9l+iCp6zJdKSkrCd999h1atWl3xF8qCggL069cP4eHhGDp0KMaNG4cePXo4FktXrVoV4eHh2Lp1q+UaSxvPzz77DLVq1cLcuXOdYkaNGuXS9oXzYefOncW+t2PHDlSpUsXp3YrLycjIwGuvvYavvvoKCxcuRNWqVdGxY0eXti3JyZMnsWTJEowZMwYjR450vF744HapatWqYeDAgRg4cCCOHj2K6667Di+++KLjwaK0MXD3teeqXbt2Ob0bsXv3bhQUFLi0GPzzzz9HUFAQFi1ahMDAQMfr06ZNc3udu3fvLvau4W+//QYATrV26tQJVatWxYwZM9C8eXOcPXsW99xzz2VzJyUl4bHHHsNjjz2GXbt2ITk5Ga+99ho+/vhjAJef9+3atcP777/v9PqpU6dQpUoVx39f7pwD/3ysKy0t7bI1Xs6lc1EIgd27dzsaIxS9vi5952bnzp2O71+JJ+clUSE+4hLZlJeXh2+//RYBAQGOzkWXU/hD6YcffnC8lp+fj3feecd2Db169cLBgwdL/ENW586ds9zPvEKFCmjSpAnWrVvn9HqTJk1QtWpVTJ061emjQ5mZmcUelMpSU+Evo5fmLPxX0qL/KnrhwgVMmTKlxDzbt2/H+fPnHZ1cSpOamoq0tDRbX5f7uEOvXr2Qn5/v+HhHURcvXnQ6vtdffx0//vgj3nnnHYwdOxbXX389HnroIcf6EF9fX9x222346quvip2XS8fkUiEhISV+TKyk8VyzZg1Wr17tFFexYkUAxc9HtWrVkJycjOnTpzt9b+vWrfj2229x8803l1rTpRo3bozGjRvjvffew+eff44777zTqYWzVSUdG4BifygsPz+/2NhER0cjLi7O6SNCpY2hu689V7311ltO/z158mQAuOI7LMA/Y+Pj4+P0Lunvv/9uuXucKw4dOoR58+Y5/jsnJwcffvghkpOTnT625O/vj969e+PTTz9FZmYmGjVq5PgF+/Dhw04f1Tt79myx1rpJSUkICwsrds4unbPAP8d/6byYM2dOsRaupd2HUlNTkZSUhFdffbXEf0y6tO1uaT788EOcPn3a8d+fffYZDh8+7DiHTZo0QXR0NKZOnep0XAsXLsSvv/7qUveuwuMoae4SycR3LIhctHDhQsci4aNHj2LmzJnYtWsXnnzySYSHh19x+wYNGqBFixZ46qmncOLECURFRWH27NmOdz3suOeee/Dpp5/iwQcfxNKlS9GqVSvk5+djx44d+PTTT7Fo0aISF/xezq233opnnnkGOTk5juOqUKECXnjhBTzwwAO46aabcMcdd2Dfvn2YNm1asc+jl6Wm5ORk+Pn54eWXX0Z2djYCAwNx00034frrr0elSpXQt29fPPzww/Dx8cFHH31U6i/VixcvRsWKFdGhQwdLx+4ubdu2xQMPPIDx48dj06ZN+Ne//oUKFSpg165dmDNnDiZNmoQePXrg119/xXPPPYd+/fo5/iBhZmYmkpOTMXDgQHz66acA/vkIx7fffou2bds6WpsePnwYc+bMwcqVK51aTxaVmpqKTz75BMOGDUPTpk0RGhqKLl264JZbbsHcuXNx++23o3Pnzti3bx+mTp2K+vXrO/3CFBwcjPr16+OTTz7BNddcg6ioKDRs2BANGzbEhAkTkJ6ejpYtW2LAgAGOdrMRERFOf+vCFRkZGXj88ccBwNLHoEoSHh6ONm3a4JVXXkFeXh6qV6+Ob7/9tti7W6dPn8ZVV12FHj164Nprr0VoaCi+++47rF271umvtpc2hu669vr164fp06dj3759Lr3rsG/fPnTt2hWdOnXC6tWr8fHHH+Ouu+667DtohTp37ozXX38dnTp1wl133YWjR4/irbfewtVXX+20dssdrrnmGgwYMABr165FTEwMPvjgAxw5cqTEd0cyMjLw5ptvYunSpXj55Zcdrz/11FNOY/Pbb7+hffv26NWrF+rXrw9/f3/MmzcPR44cwZ133unYLjU1FW+//TZeeOEFXH311YiOjsZNN92EW265Bc8//zz69++P66+/Hlu2bMGMGTOK3cOSkpIQGRmJqVOnIiwsDCEhIWjevDkSExPx3nvvIT09HQ0aNED//v1RvXp1HDx4EEuXLkV4eDi++uqrK45NVFQUbrjhBvTv3x9HjhzBxIkTcfXVV+O+++4D8M/99uWXX0b//v3Rtm1b9O7d29FutmbNmsXa85amtLlLJJWaZlRE5iip3WxQUJBITk4Wb7/9tlPbwcJ2s6W1Q9yzZ49IS0sTgYGBIiYmRjz99NNi8eLFJbabbdCgQbHt+/btW6x94IULF8TLL78sGjRoIAIDA0WlSpVEamqqGDNmjMjOznbEudJuVgghjhw5Ivz9/cVHH31U7HtTpkwRiYmJIjAwUDRp0kT88MMPom3btsVaO5alpnfffVfUqlVL+Pn5OY3LqlWrRIsWLURwcLCIi4sTw4cPF4sWLSqxPW3z5s3F3XfffcVjle2dd94RqampIjg4WISFhYlGjRqJ4cOHi0OHDomLFy+Kpk2biquuusqpVacQQkyaNEkAEJ988onjtT/++ENkZGSIqlWrisDAQFGrVi0xaNAgR9vYktrNnjlzRtx1110iMjJSAHDMnYKCAjFu3DiRkJAgAgMDRUpKipg/f36J8+vHH38UqampIiAgoFjr2e+++060atVKBAcHi/DwcNGlSxexfft2p+0L26QeO3as1HE6fPiw8PPzE9dcc43LY3u5vH/++ae4/fbbRWRkpIiIiBA9e/YUhw4dcqo/NzdXPPHEE+Laa68VYWFhIiQkRFx77bViypQpTrlKG0MhXJ/nAEpta9u9e3cRHBwsTp486dLxbt++XfTo0UOEhYWJSpUqicGDBzu1JL3S/t5//31Ru3ZtERgYKOrWrSumTZvmyF1Uae1mL215XNK8S0hIEJ07dxaLFi0SjRs3duxrzpw5pR5fgwYNhK+vr/jzzz8drxW2Ui1s0X38+HExaNAgUbduXRESEiIiIiJE8+bNxaeffuqUKysrS3Tu3FmEhYU5tcM+f/68eOyxx0S1atVEcHCwaNWqlVi9enWJ97Avv/xS1K9fX/j7+xdrPbtx40bRrVs3UblyZREYGCgSEhJEr169xJIlS0o9vqJjNWvWLPHUU0+J6OhoERwcLDp37iz++OOPYvGffPKJSElJEYGBgSIqKkr06dPHaXwKxygkJKTE/V1u7hLJ4iOEzRVXROS1BgwYgN9++63Uv6yts02bNuG6667Dhg0bSl1gSXo5fvw4qlWrhpEjR+K5555TXY5HxcTEICMjo9gfeTRZzZo10bBhQ8yfP9/lbVJSUhAVFVXsr9B7k2XLlqFdu3aYM2cOevToobocIim4xoKIihk1ahTWrl2LVatWqS7Fspdeegk9evTgQ4VBMjMzkZ+ff8VFu95m27ZtOHfuHEaMGKG6FKXWrVuHTZs2ISMjQ3UpRFRGXGNBRMXEx8cXWyRpitmzZ6sugVz0/fffY/v27XjxxRdx2223ubTGwJs0aNCg2N9kKU+2bt2K9evX47XXXkO1atVwxx13qC6JiMqI71gQEZESzz//PIYNG4bk5GRHdyMqPz777DP0798feXl5mDVrltNfmSYiM3GNBRERERERlRnfsSAiIiIiojLjgwUREREREZWZ1y/eLigowKFDhxAWFsY/b09EREREZIEQAqdPn0ZcXBx8fS//noTXP1gcOnQINWrUUF0GEREREZGxDhw4gKuuuuqyMV7/YBEWFgbgn8EIDw9XXA0RERERkTlycnJQo0YNx+/Ul+P1DxaFH38KDw/ngwURERERkQ2uLCng4m0iIiIiIiozpQ8Wo0ePho+Pj9NX3bp1Hd8/f/48Bg0ahMqVKyM0NBTdu3fHkSNHFFZMREREREQlUf6ORYMGDXD48GHH18qVKx3fe/TRR/HVV19hzpw5WL58OQ4dOoRu3boprJaIiIiIiEqifI2Fv78/YmNji72enZ2N999/HzNnzsRNN90EAJg2bRrq1auHn376CS1atPB0qUREREREVArl71js2rULcXFxqFWrFvr06YP9+/cDANavX4+8vDykpaU5YuvWrYv4+HisXr1aVblERERERFQCpe9YNG/eHJmZmahTpw4OHz6MMWPGoHXr1ti6dSuysrIQEBCAyMhIp21iYmKQlZVVas7c3Fzk5uY6/jsnJ0dW+URERERE9H+UPlikp6c7/n/jxo3RvHlzJCQk4NNPP0VwcLCtnOPHj8eYMWPcVSIREREREblA+UehioqMjMQ111yD3bt3IzY2FhcuXMCpU6ecYo4cOVLimoxCTz31FLKzsx1fBw4ckFw1ERERERFp9WBx5swZ7NmzB9WqVUNqaioqVKiAJUuWOL6/c+dO7N+/Hy1btiw1R2BgoOOP4fGP4hEREREReYbSj0I9/vjj6NKlCxISEnDo0CGMGjUKfn5+6N27NyIiIjBgwAAMGzYMUVFRCA8Px5AhQ9CyZUt2hCIiIiIi0ozSdyz+/PNP9O7dG3Xq1EGvXr1QuXJl/PTTT6hatSoA4I033sAtt9yC7t27o02bNoiNjcXcuXNVlmzL1KlTUbNmTUydOtUtcTLz6hBrJV5WXhk16HBMVuNVHpsO80tWrEnjJWveyqzFSqzqe5MO81yX60dWzapjTckpY/9XitNp7snML7turQgvl52dLQCI7OxsZTUkJCQIACIhIcEtcTLz6hBrJV5WXhk16HBMVuNVHpsO80tWrEnjJWveyqzFSqzqe5MO81yX68dKvOp5o7pWU47pSnE6zT2Z+WXXLZuV36X5YOEBb7/9tkhISBBvv/22W+Jk5tUh1kq8rLwyatDhmKzGqzw2HeaXrFiTxkvWvJVZi5VY1fcmHea5LtePrJpVx5qSU8b+rxSn09yTmV923bJZ+V3aRwghyvy2h8ZycnIQERGB7OxsLuQmIiIiIrLAyu/SWnWFIiIiIiIiM/HBgoiIiIiIyowPFh7kykp/HboMlOdOF94eq0sd3lyzyq5BMnPrcn5l12PKfU11TtW16pRb9nznuTYvVnZurUlf8aGYDou3C7my0t+VGF1iTcnJWP3q8OaaVcXJzq3L+ZVdjyn3NdU5VdeqU27Z853n2rxY2bk9jV2hitDpwcKVlf46dBkoz50uvD1Wlzq8uWaVXYNk5tbl/Mqux5T7muqcqmvVKbfs+c5zbV6s7Nyexq5QRbArFBERERGRPewKRUREREREHsUHCw8xaUGVjLwmLPzSYREXFwFaj3X3uZVVp465Zcernm+mHZuVWFPulVbidThfKnPqcPwy6tXh+GXHq56LWpH+wSzFdFlj4erCHKsLeEzJK6NOd8fK2reVeFl57cTrULO7zpvVGmXUqWNu2fGq55tpx2Yl1pR7pZV4Hc6Xypw6HL+MenU4ftnxqueibFy8XYQuDxYmLaiSkdeEhV86LOLiIkDrse4+t7Lq1DG37HjV8820Y7MSa8q90kq8DudLZU4djl9GvTocv+x41XNRNi7eLoKLt4mIiIiI7OHibSIiIiIi8ig+WBARERERUZnxwcID2MXEeqwJ3aZMPFfeXocpOXWJ1SW3aXlNjLUSz+tI/djqcFymdV1Smc/qNjLPg3LSV3wopsPibSsr+63EysytOlZVnKycsmpgHWbm1CVWl9ym5TUx1ko8ryP1Y6vDccm8fmTUojKf1W1kngcZ2BWqCB0eLNjFxHqsCd2mTDxX3l6HKTl1idUlt2l5TYy1Es/rSP3Y6nBcpnVdUpnP6jYyz4MM7ApVBLtCERERERHZw65QRERERETkUXywICIiIiKiMuODhQeo6rChulOE6v27u1Zv7qZhNVaXsTAp1pSOPDrklVWD7OtSVj3eNh9lxZqS07RaVderelx1iTWG9BUfiumweNvKin5XY12Jk7Ffd9eoS6zK8ZQdb9K5lZlbdazKsZV1DnQ4t7rMRZn1eNt8lBVrSk7TalVdr+px1SVWJXaFKkKHBwtVHTZUd4pQvX931+rN3TSsxuoyFibFmtKRR4e8smqQfV3Kqsfb5qOsWFNymlar6npVj6susSqxK1QR7ApFRERERGQPu0IREREREZFH8cHCQ1QvUpJRh+yFjTLjdanDarxp88hOfp3iTa1FdazucbrE6lKHDrEm1uHun5Pedj2Y+HNI9ZhZjdWS9A9mKabDGgsh1C9SklGHzP3LjtelDqvxps0jO/l1ije1FtWxusfpEqtLHTrEmliHK3EmXF+q67QaKzO/6jGzGuspXLxdhC4PFqoXKcmoQ/bCRpnxutRhNd60eWQnv07xptaiOlb3OF1idalDh1gT63D3z0lvux5M/DmkesysxnoKF28XwcXbRERERET2cPE2ERERERF5FB8siIiIiIiozPhgoYAu3Q9Ux+rQwcmkGmTFmnhsnoi3uo3M/DqMuUmxusxpV+NVj5fVeNPy2t3GhJ95qq8fHeajLrFWtjG+89PlSF/xoZgui7eLsrPi38o2psS6ms/qeMk4Jh1qkBVr4rF5It7qNjLz6zDmJsXqMqddjVc9XlbjTctrdxt3x5pyvavev051yJpjdvKqxK5QRej4YKFL9wPVsTp0cDKpBlmxJh6bJ+KtbiMzvw5jblKsLnPa1XjV42U13rS8drcx4Wee6utHh/moS6yVbXTs/HQ57ApVBLtCERERERHZw65QRERERETkUXywICIiIiKiMuODhYeY1r2ivOf05D5Un2crsTp005BZs6y8OswJ2XNfVZwOeVWfX13mosx4HY5R5XmTsX/ZsXbiZdXjzjgd5qzWpK/4UEyXxdt2OgBY2UZGbHnO6cl9qD7PVmJlj7XqmmXl1WFOyJ77quJ0yKv6/OoyF2XG63CMKs+bjP3LjrUTL6sed8bpMGc9jV2hitDlwcK07hXlPacn96H6PFuJ1aGbhsyaZeXVYU7Invuq4nTIq/r86jIXZcbrcIwqz5uM/cuOtRMvqx53xukwZz2NXaGKYFcoIiIiIiJ72BWKiIiIiIg8ig8WHqRyMZcpCyZN2r+3xuqyWNNOvGmxOuWWWYvq61qXWDvbyVr8KntRrcx4b1tErcPc1OF6tpNb5j5M+f1CO9I/mKWYLmsshFC7mMvdcdy/98bKOgeeiDctVqfcMmtRfV3rEmtnOzv5XdlGVl5PxLv7/Ki+TnSYmzpcz3Zyy9yHqjirsZ5g5OLt8ePHCwDikUcecbx27tw5MXDgQBEVFSVCQkJEt27dRFZWlqW8Oj1YqFzMZcqCSZP2762xuizWtBNvWqxOuWXWovq61iXWznayFr/KXlQrM97bFlHrMDd1uJ7t5Ja5D1N+v/AE4xZvr127Fr169UJ4eDjatWuHiRMnAgAeeughLFiwAJmZmYiIiMDgwYPh6+uLVatWuZybi7eJiIiIiOwxavH2mTNn0KdPH7z77ruoVKmS4/Xs7Gy8//77eP3113HTTTchNTUV06ZNw48//oiffvpJYcVERERERHQp5Q8WgwYNQufOnZGWlub0+vr165GXl+f0et26dREfH4/Vq1d7ukwiIiIiIroMpQ8Ws2fPxoYNGzB+/Phi38vKykJAQAAiIyOdXo+JiUFWVlapOXNzc5GTk+P0pZounRmsbKOyy4Hq8dK9Tl06cugQq0tu3eOsxnpbrabNVW+8r+pQq6waZNasQ6zseBNr0WXeaUn6io9S7N+/X0RHR4vNmzc7Xmvbtq1j8faMGTNEQEBAse2aNm0qhg8fXmreUaNGCQDFvlQu3pbVCcBKrNVt3B2nOqeM/auq0+p5Vz2eMmN1ya17nNVYb6vVtLnqjfdVHWqVVYPVeNNiZcebWIsu885TjOgKNW/ePAFA+Pn5Ob4ACB8fH+Hn5ye+++47AUCcPHnSabv4+Hjx+uuvl5r3/PnzIjs72/F14MAB5Q8WunRmsLKNyi4HqsdL9zp16cihQ6wuuXWPsxrrbbWaNle98b6qQ62yapBZsw6xsuNNrEWXeecpRnSFOn36NP744w+n1/r374+6detixIgRqFGjBqpWrYpZs2ahe/fuAICdO3eibt26WL16NVq0aOHSftgVioiIiIjIHiu/S/t7qKZiwsLC0LBhQ6fXQkJCULlyZcfrAwYMwLBhwxAVFYXw8HAMGTIELVu2dPmhgoiIiIiIPEPZg4Ur3njjDfj6+qJ79+7Izc1Fx44dMWXKFNVlERERERHRJZS3my1q2bJljj+OBwBBQUF46623cOLECfz999+YO3cuYmNj1RVYBt7YNUZmftUdF0zpNiWrVjt1WNnGhPOgQ60mxcqeWzL3ocO9SVasafVaiTXxvFndRpdjlFmLiXXrcP60JX3Fh2JWFpzI5OpKf1VxVmNl57daiwnjJmPfsmq1U4eVbUw4DzrUalKs7Lklcx863JtkxZpWr5VYE8+b1W10OUaZtZhYtw7nz5OM6ArlKbo8WHhj1xiZ+VV3XDCl25SsWu3UYWUbE86DDrWaFCt7bsnchw73JlmxptVrJdbE82Z1G12OUWYtJtatw/nzJCO6QnkKu0IREREREdlj5XdprdZYEBERERGRmfhg4UEqFkvqtIBKVj261W013sR6TFrgZkqtpuxfZrwu9aqOVb1/HWqVmVt1LMfBeqyssbWzjS5zXlvSP5ilmC5rLIRQs1jS3fssS7ysenSr22q8ifW4+9zoME9V12rK/mXG61Kv6ljV+9ehVpm5VcdyHKzHyhpbO9voMuc9iYu3i9DpwULFYkmdFlDJqke3uq3Gm1iPSQvcTKnVlP3LjNelXtWxqvevQ60yc6uO5ThYj5U1tna20WXOexIXbxfBxdtERERERPZw8TYREREREXkUHyyIiIiIiKjM+GDhIaq7RajOq7qzhurOE6o7angqPzvJWItV1SlJhw5QOtZiNd7E68Pk/DqMt93tVF+bquaeKT+nTclpBOkrPhTTZfG2qyv9rXYEMCWvlf2bkFPWeMqM9UR+ledOdm5V8031tSM7XqdarMabeH2YnF+H8ba7neprU9XcUzlGqs+P7PnnSewKVYQuDxaqu0Wozqu6s4bqzhOqO2p4Kj87yViLVdUpSYcOUDrWYjXexOvD5Pw6jLfd7VRfm6rmnik/p03JqQq7QhXBrlBERERERPawKxQREREREXkUHyyIiIiIiKjM+GDhYbp0KDGh+4IusbJyq+4SIivWlJw6xJo2Zz2R35vuTTqcMxPGydU4meNpNd7U3J7Yh0n3StW12onXnvQVH4rpsni7kNXV/1biZcS6O87EWFm5ZYytrFpl1Vve55dpc9YT+b3p3qTDOTNhnFyNkzmeVuNNze2JfZh0r1Rdq514FdgVqgjdHix06VBiQvcFXWJl5VbdJURWrCk5dYg1bc56Ir833Zt0OGcmjJOrcTLH02q8qbk9sQ+T7pWqa7UTrwK7QhXBrlBERERERPawKxQREREREXkUHyw8yJRFSqpzujuv6nGXHW/CIk7Vc9oT8bKuBZn5GavPfDPhXmp1/yr2rXo+WYnl/cs9+1A1d3U4Ji1J/2CWYjqtsXB1gY6VhTwyYlXndHde1eMuO97dx6f6/PM8eC4/Y/WZbybcS63uX8W+Vc8nK7G8f7lnH6rmrg7H5ClcvF2ETg8WpixSUp3T3XlVj7vseBMWcaqe056Il3UtyMzPWH3mmwn3Uqv7V7Fv1fPJSizvX+7Zh6q5q8MxeQoXbxfBxdtERERERPZw8TYREREREXkUHyyIiIiIiKjM+GDhAd7cIcJKvOpODyq7o8jKq0MNsnPrUouJ50WH8y27Fm89hzoclw6xutRhWqzseNXXh4ycJvweYwTpKz4U02Hxtqur+12N80S8jNyy6tU9TmZeHWqQnVuXWkw8Lzqcb9m1eOs51OG4dIjVpQ7TYmXHq74+ZOR0JU6H41aBXaGK0OHBwps7RFiJV93pQWV3FFl5dahBdm5dajHxvOhwvmXX4q3nUIfj0iFWlzpMi5Udr/r6kJHThN9jVGFXqCLYFYqIiIiIyB52hSIiIiIiIo/igwUREREREZUZHyw8RGWng7Js401dF9xdpyfiTe2gJDM3c5bvzkq6d42RVaMusa7G636edLlGrG6jU92m/dzT5bx4faco6Ss+FNNh8bYQajsdlGUbd9fiSpysGt1dpyfiZdchsxZZuZmzfHdWUjleKs+pLrGuxut+nnS5Rqxuo1PdsuJ1ySsrt8x6ZWFXqCJ0ebBQ2emgLNt4U9cFd9fpiXhTOyjJzM2c5buzku5dY2TVqEusq/G6nyddrhGr2+hUt2k/93Q5LyZ2imJXqCLYFYqIiIiIyB52hSIiIiIiIo/ig4UCOizA02UBmCkLo3QZL11qUT0fdIm1Eq+6Xo6BOTll55a5iFVWrOqcsmqQmVtVnJ141efXlVjVx2MM6R/MUkyXNRZFWVmQo0OsDnWYklO33LJqUT0fdIm1Eq+6Xo6BOTll57Zahw61qM4pqwaZuVXF2YlXfX5diVV9PCpx8XYROj5Y6LAAT5cFYKYsjNJlvHSpRfV80CXWSrzqejkG5uSUnVvmIlZZsapzyqpBZm5VcXbiVZ9fV2JVH49KXLxdBBdvExERERHZw8XbRERERETkUXywICIiIiKiMuODhQep6NBgQlcEHbqZmNbtojzUofLY3DlnVI2NCV1mVOf0VH7V4yEjr+5xMmO9vQ5T5pbqnLJzG9styuoCjjNnzohnn31WtGzZUiQlJYnExESnL93otHjb1ZX/7oxzNZeqfLJyWt3Gal5ZNbMO+fW6az6oPA7dr1WT5oMn8qseDxl5dY+TGevtdZgyt1TnlJ3bzn1MFqldoe68805RrVo1MXz4cPHGG2+IiRMnOn1ZMWXKFNGoUSMRFhYmwsLCRIsWLcTXX3/t+P65c+fEwIEDRVRUlAgJCRHdunUTWVlZlvah04OFig4NJnRF0KGbiWndLspDHSqPzZ1zRtXYmNBlRnVOT+VXPR4y8uoeJzPW2+swZW6pzik7t07doqR2hYqMjMSCBQvQqlUrK5uV6KuvvoKfnx9q164NIQSmT5+OCRMmYOPGjWjQoAEeeughLFiwAJmZmYiIiMDgwYPh6+uLVatWubwPdoUiIiIiIrLHyu/Slh8sEhMT8fXXX6NevXplKrI0UVFRmDBhAnr06IGqVati5syZ6NGjBwBgx44dqFevHlavXo0WLVq4lI8PFkRERERE9khtNzt27FiMHDkSZ8+etV1gSfLz8zF79mz8/fffaNmyJdavX4+8vDykpaU5YurWrYv4+HisXr261Dy5ubnIyclx+tKF6oVTduJVL+KyEq96EZcpC85MOiZZed29+Fl2vC65TatD9Zy0E6/6OpYdb9K5U3lPMSlW9zhZsXYWV3PxdgmSk5NFWFiYCA0NFQ0bNhQpKSlOX1b98ssvIiQkRPj5+YmIiAixYMECIYQQM2bMEAEBAcXimzZtKoYPH15qvlGjRgkAxb50WGPh6kIcV+M8ES8rt4y8MvatOtaUnLL2LyuvO+eTJ+J1yW1aHarnpJ141dex7HiTzp3Ke4pJsbrHyYq1Ovdk1eEJUhdvjx49+rJfVuXm5opdu3aJdevWiSeffFJUqVJFbNu2zfaDxfnz50V2drbj68CBA9o8WKheOGUnXvUiLivxqhdxmbLgzKRjkpXX3YufZcfrktu0OlTPSTvxqq9j2fEmnTuV9xSTYnWPkxVrZ3E1F28rkJaWhqSkJNxxxx1o3749Tp48icjISMf3ExISMHToUDz66KMu5eMaCyIiIiIie6z8Lu1vdyfr16/Hr7/+CgBo0KABUlJS7KZyUlBQgNzcXKSmpqJChQpYsmQJunfvDgDYuXMn9u/fj5YtW7plX0RERERE5B6WHyyOHj2KO++8E8uWLXO8k3Dq1Cm0a9cOs2fPRtWqVV3O9dRTTyE9PR3x8fE4ffo0Zs6ciWXLlmHRokWIiIjAgAEDMGzYMERFRSE8PBxDhgxBy5YtXe4IRUREREREnmG5K9SQIUNw+vRpbNu2DSdOnMCJEyewdetW5OTk4OGHH7aU6+jRo8jIyECdOnXQvn17rF27FosWLUKHDh0AAG+88QZuueUWdO/eHW3atEFsbCzmzp1rtWQtmNAFyKTuGjJrVr1/O/HujNPhXMkeV1e3Ud1JxKT5pfr+YcK50uHc6/DzQGa8acfnrnmjev7LqkGXWNm57W6jJasLOMLDw8XPP/9c7PU1a9aIiIgIq+mk0+Uvb7u6ut/VOBk5rexbVq0yY12NV71/O/HujNPhXMkeV1e3UXk9Wo2VGe/usbIa7y3nSodzL+PYdYo37fjcNW9Uz39ZNegSKzu33W08RWpXqNDQULFx48Zir2/YsEGEhYVZTSedLg8WJnQBMqm7hsyaVe/fTrw743Q4V7LH1dVtVHcSMWl+qb5/mHCudDj3Ovw8kBlv2vG5a96onv+yatAlVnZuu9t4itSuULfeeitOnTqFWbNmIS4uDgBw8OBB9OnTB5UqVcK8efOspJOOXaGIiIiIiOyR+pe3//Of/yAnJwc1a9ZEUlISkpKSkJiYiJycHEyePNl20UREREREZC7LXaFq1KiBDRs24LvvvsOOHTsAAPXq1UNaWprbiyMiIiIiIjNYesciLy8P/v7+2LZtGzp06IAhQ4ZgyJAhfKhwkUldKGTG6tSJwaSaVZ83XWrWYRx0mJNW43Wq2co2OtTtzq48Kq8JHea4DvNbh7wyatDlXmMnXtVYmHrv0prVBRyJiYli06ZNNpZ+qKHL4m0hzOpCITNW5jjoMGZW4005b3a28dZx0GFOWo3XqWYr2+hQtytxJlwTOsxxHea3Dnll1KDLvcZOvKqxMPXe5WlSu0K999574uabbxZ//fWXreI8TacHC5O6UMiM1akTg0k1qz5vutSswzjoMCetxutUs5VtdKjbnV15VF4TOsxxHea3Dnll1KDLvcZOvKqxMPXe5WlSu0KlpKRg9+7dyMvLQ0JCAkJCQpy+v2HDBivppGNXKCIiIiIie6z8Lm158fZtt91mty4iIiIiIvJW0t8/UczEj0KZ8rayDm+Vqx5TWW+j8mMV3vV2u6firWzjzR8dsRtvZRsdxsNqvA73Nk/W5ImPoqi6X+pw/epyPlRcszqMpydJXWNhGp0eLBI0WPDjzhpk1SkjVvX+rW6j8rhMG1NTxsAT8Va2UXnvsBtv4vjpNEd0uLd5sibZx2tlO3fXosP1q8v5UHHN6jCeniT1wcLHx0f4+vqW+qUbnR4sdHjC5TsWZvyrHt+x4DsWduKtbKPDv3hajTdx/HSaIzrc2zxZE9+xUF+nJ3LzHQv5pC7e/vLLL53+Oy8vDxs3bsT06dMxZswYDBgwwEo66bh4m4iIiIjIHiu/S1t+sCjNzJkz8cknnxR78FCNDxZERERERPZY+V3a0l/evpwWLVpgyZIl7kpHREREREQGccuDxblz5/Dmm2+ievXq7kjnlWT9eXmVf95e5THJqNPONqrPl4mxsnK7M86Uue3unKrnpyfyq5p7suaU1W14P9Qr1k68lW1U36N0Gjd316LLOGjH6gKOyMhIUalSJcdXZGSk8PPzE2FhYeLLL7+0vCBENl0Wb7u62t9qVwB357Wyf5XHJKNOO9uoPl8mxsrK7c44U+a2u3Oqnp+eyK9q7smaU1a34f1Qr1g78Va2UX2P0mnc3F2LLuPgCVK7QmVmZjp9ffjhh2LhwoXixIkTtoqVTZcHC9UdREzoMmQlXpfuJqrPl4mxsnLr3vHD3fuWkVP1/PREflVzT9acsroN74d6xdqJt7KN6nuUTuPm7lp0GQdPkNoVyjRcvE1EREREZI/0xdsrVqzA3Xffjeuvvx4HDx4EAHz00UdYuXKlnXRERERERGQ4yw8Wn3/+OTp27Ijg4GBs2LABubm5AIDs7GyMGzfO7QUSEREREZH+LD9YvPDCC5g6dSreffddVKhQwfF6q1atsGHDBrcW50106IygexcJVcehujOJp/LL3Ifq7hy61CArVvU9QWYNdrbR6fpSPdYqO9144/51yCmrBpmxqutw53lQef6N7wgFWO8KFRwcLPbt2yeEECI0NFTs2bNHCCHEnj17RGBgoNV00umyeNvKKn8rsVbi3R1nJd5dMSbEWY31VH6Z+5A1Z63E6lCDrFjV9wSZNdjZRqfrS/VYuztW9XWkev865JRVg8xY1XW48zyoPP927m2eILUrVGJioli8eLEQwvnBYvr06aJevXpW00mny4OFDp0RdO8ioeo4VHcm8VR+mftQ3Z1Dlxpkxaq+J8iswc42Ol1fqsdaZacbb9y/Djll1SAzVnUd7jwPKs+/jh2hhJDcFWr8+PH4+OOP8cEHH6BDhw74+uuv8ccff+DRRx/Fc889hyFDhlhJJx27QhERERER2WPld2l/q8mffPJJFBQUoH379jh79izatGmDwMBAPP7449o9VBARERERkWdYXrzt4+ODZ555BidOnMDWrVvx008/4dixYxg7dqyM+ryGzguV3L1fUxY+qV6gZ9qCXG8dL1mxqvcvO16H49NpEbfM/KrHz5RF13bqsLKNLtept80xlfmsxOuSV3vSP5ilmC5rLKwsyHEl1tV8Kvbr7n3KiFNdp9VYHWrw1vGSFat6/7LjdTg+q8dman7V4+dKnGnXp51tdLlOvW2OqcxnJV6XvCpIXbx95swZ8eyzz4qWLVuKpKQkkZiY6PSlG10eLHReqOTu/Zqy8En1Aj3TFuR663jJilW9f9nxOhyfTou4ZeZXPX6mLLq2U4eVbXS5Tr1tjqnMZyVel7wqSF283bt3byxfvhz33HMPqlWrBh8fH6fvP/LII1bSScfF20RERERE9khdvL1w4UIsWLAArVq1sl0gERERERF5F8uLtytVqoSoqCgZtRARERERkaEsP1iMHTsWI0eOxNmzZ2XU47VM6Tqguk6ZsTrUoXr/OsSq3r/M3Docm+rjciVeh3NlN15GParvu6rnjDd2nZI9D02rw9VY1efChOPRntUFHMnJySIsLEyEhoaKhg0bipSUFKcv3eiyeNvVVf+uxsnKq7pOmbE61KF6/zrEqt6/zNw6HJvq43IlXodzZTdeRj2q77uq54w7x0jW/q3Gyp6HptXhaqzqc2HC8aggtSvU6NGjL/ulG10eLEzpOqC6TpmxOtShev86xKrev8zcOhyb6uNyJV6Hc2U3XkY9qu+7queMN3adkj0PTavD1VjV58KE41FBalco07ArFBERERGRPVZ+l7a8xoKIiIiIiOhSfLAgIiIiIqIy44OFh6nuVCErVnbnA3fXobrzhY61qI5VPddU79/uNjrMTR3mvC7jpkMdss6bnW10GWc72+pWj9V40/LaiTe1Fq/pAFUS6Ss+FNNl8XYhK6v+TYqVtX9ZdVitQebx6VKL6ljVc031/u1uo8Pc1GHO6zJuOtQh67zZ2UaXcbazrW71WI03La+deFNrKct8UUFqVyjT6PZgobpThaxY2Z0P3F2H6s4XOtaiOlb1XFO9f7vb6DA3dZjzuoybDnXIOm92ttFlnO1sq1s9VuNNy2sn3tRadO4AVRKpXaG6d++OZs2aYcSIEU6vv/LKK1i7di3mzJljJZ107ApFRERERGSP1K5QP/zwA26++eZir6enp+OHH36wmo6IiIiIiLyA5QeLM2fOICAgoNjrFSpUQE5OjluK8laq//y76gWTOi161WGcdajX5MWbVuJV16Dy/OkwprLyqjxXKuJUHa8u9x9T5rsu90gZsariVO/fhLHUhtXPWTVt2lSMGTOm2OujRo0S1113ndV00um0xsKVxToyF//IiJdVr9VadarF1Xgd6pU9zrLrUT3WrsarPH86jKmsvCrPlYo4Vcery/3HlPmuyz1SRqyqONX7N2EsZZK6ePt///uf8Pf3FxkZGSIzM1NkZmaKe+65R/j7+4t58+ZZyjVu3DjRpEkTERoaKqpWrSpuvfVWsWPHDqeYc+fOiYEDB4qoqCgREhIiunXrJrKyslzeh04PFqr//LvqBZM6LXrVYZx1qNfkxZtW4lXXoPL86TCmsvKqPFcq4lQdry73H1Pmuy73SBmxquJU79+EsZRJ6uJtAFiwYAHGjRuHTZs2ITg4GI0bN8aoUaPQtm1bS3k6deqEO++8E02bNsXFixfx9NNPY+vWrdi+fTtCQkIAAA899BAWLFiAzMxMREREYPDgwfD19cWqVatc2gcXbxMRERER2WPld2lbDxayHDt2DNHR0Vi+fDnatGmD7OxsVK1aFTNnzkSPHj0AADt27EC9evWwevVqtGjR4oo5+WBBRERERGSP1K5Qo0ePRkFBQbHXs7Oz0bt3b6vpiuUAgKioKADA+vXrkZeXh7S0NEdM3bp1ER8fj9WrV5dpX0RERERE5D6WHyzef/993HDDDdi7d6/jtWXLlqFRo0bYs2eP7UIKCgowdOhQtGrVCg0bNgQAZGVlISAgAJGRkU6xMTExyMrKKjFPbm4ucnJynL50YEJXGB26TMmuw4TODjp0H1HdrcVOFwwd5o03zgXZdVjdzsTjLO+xKrs56VCDLtexzFq8/Rh1ONfGsLqA48SJE6Jnz54iLCxMvPPOO+Lxxx8XFSpUEE8//bTIy8uzvCCk0IMPPigSEhLEgQMHHK/NmDFDBAQEFItt2rSpGD58eIl5Ro0aJQAU+1K9eDtBYScHV/NZ2a/MeJl1yIhVFWc1VlYNrsbLPC4723AuqO9EY2c7E4+zvMe6+/4gK16H+6kO46BTrOx4d59DT9wrVZDaFarQU089JXx8fESFChXEd999ZzeNEEKIQYMGiauuukrs3bvX6fUlS5YIAOLkyZNOr8fHx4vXX3+9xFznz58X2dnZjq8DBw5o8WBhQlcYHbpMya7DhM4OOnQfUd2txU4XDB3mjTfOBdl1WN3OxOMs77EquznpUIMu17HMWrz9GHU41ypJ7wo1efJkPPnkk7jtttuwfv16+Pn5YebMmbj22mst5RFCYMiQIZg3bx6WLVuG2rVrO32/cPH2rFmz0L17dwDAzp07UbduXS7eJiIiIiKSzMrv0v5Wk3fq1Anr1q3D9OnT0aNHD5w7dw7Dhg1DixYtMGbMGAwfPtzlXIMGDcLMmTPx5ZdfIiwszLFuIiIiAsHBwYiIiMCAAQMwbNgwREVFITw8HEOGDEHLli1deqggIiIiIiLPsPyORYcOHTB9+nTExcU5vb5gwQL8+9//xuHDh13fuY9Pia9PmzYN/fr1AwCcP38ejz32GGbNmoXc3Fx07NgRU6ZMQWxsrEv74DsWRERERET2SG03u3jx4mIPFQDQuXNnbNmyxVIu8c8aj2JfhQ8VABAUFIS33noLJ06cwN9//425c+e6/FChCxVdB0zoimBixwnVsSbt37S8JnYqkVmHrHMnM783x+pSh4nHJzO36ljeY82bC7J/1ihnZxHHDz/8IPr06SNatGgh/vzzTyGEEB9++KFYsWKFnXRSWVlwIourK//dGedqLlX1yarRaqwudZgyZrLGVoe8MueNzFpUj53V/cvM782xutRh4vHJzK06lvdY8+aCzHplkdoV6rPPPhPBwcHi3//+twgMDBR79uwRQggxefJkkZ6ebr1ayXR4sFDRdcCErggmdpxQHWvS/k3La2KnEpl1yDp3MvN7c6wudZh4fDJzq47lPda8uSD7Z40MUrtCpaSk4NFHH0VGRgbCwsKwefNm1KpVCxs3bkR6enqpf7hOFa6xICIiIiKyR+oai507d6JNmzbFXo+IiMCpU6espiMiIiIiIi9g+cEiNjYWu3fvLvb6ypUrUatWLbcU5W1ULxSSGSs7t5VtTFpkpdNCYVfjdZkXqs+JDuMg8zrSLb8Oc8MT9Zh4flRfizJj7W7jzWNi0s9iu9uprt24hdolsfo5q3Hjxon69euLn376SYSFhYkVK1aIjz/+WFStWlW8+eabNj65JZcOayysLLwxLVZ2bivbyKhD1ljIHGNZ8brMC9XnRIdxkHkd6ZZfh7nhiXpMPD+qr0WZsXa38eYxUVmrnXNnZzvVtds9TtmkLt4uKCgQL7zwgggJCRE+Pj7Cx8dHBAUFiWeffdZWsbLp8GCheqGQzFjZua1sY9IiK50WCrsar8u8UH1OdBgHmdeRbvl1mBueqMfE86P6WpQZa3cbbx4Tk34W291Ode06LNQuidTF24UuXLiA3bt348yZM6hfvz5CQ0PtpJGOi7eJiIiIiOyx8ru0v92dBAQEoH79+nY3JyIiIiIiL2J58TYREREREdGl+GDhIe7uCuDObgS6dL2wEuvtXVE8ES9jbsioQUacyjpNqNFqvLtz6nD96dZBSHWszLEz6dh0qNdKrCn3BR1ircTrMme0JH3Fh2I6LN4Wwv1dAVzN585cnoh353HpWIdO8TLmhowaZMSprNOEGq3GuzunDtef7Lkvsx7V88a0cTCtXiuxptwXdIi1Eq/LnPEUqV2hTKPLg4W7uwK4sxuBLl0vrMR6e1cUT8TLmBsyapARp7JOE2q0Gu/unDpcf7p1EFIdK3PsTDo2Heq1EmvKfUGHWCvxuswZT/FIVyhTsCsUEREREZE9Vn6X5hoLIiIiIiIqMz5YEBERERFRmfHBwkNUdXBwV/cVO/GmdI3QoQ5v3L/McyCrDm88XzrUKzNWdR06jIGdLjI61KI6Vod7jtV4Xe6rVmJ1qEG3Wtx9b9GK9BUfiumyeNvVVf6ejnM1j514V2Ldnc9OrA51eOP+ZZ4DWXV44/nSoV6Zsarr0GEMrI6XLrWojtXhnmM1Xpf7qpVYHWrQrRZ331tkY1eoInR5sFDVwcFd3VfsxJvSNUKHOrxx/zLPgaw6vPF86VCvzFjVdegwBna6yOhQi+pYHe45VuN1ua9aidWhBt1qcfe9RTZ2hSqCXaGIiIiIiOxhVygiIiIiIvIoPlh4iKrF26r3r3ohmAk5VZ9LO/GqF6mpngMy8qqsU/UiWnfXqsOCTnfWqfI+ZjVWt1pUj51pNXvzeFmJVV2DsQu3AS7e9hRXF+G4O071/q3UWV5zqj6XduLdHatDrarPg8o6VZ57GbXKnE8q6lR5H7Maq1stqsdOZrzK+6CJ42UlVnUNdq45mbh4uwhdHixULd5WvX/VC8FMyKn6XNqJV71ITfUckJFXZZ2qF9G6u1YdFnS6s06V9zGrsbrVonrsTKvZm8fLSqzqGnRauC0EF2874eJtIiIiIiJ7uHibiIiIiIg8ig8WRERERERUZnyw8BCVHQO8sWOCzFhTcuoSq0sdppwLE7o/mXAspsbqUocOsVbiTfi5aCVW9zirsa7Gy9y/rJrtxss+Vm1JX/GhmC6Lt62s8Hcl1t35rMZayalDvVZiTcmpS6wudZhyLlReY6pq9Ob5ZDVWlzp0iLUS7875rUOs7nFWY12Nl7l/O9uYeO2qwK5QRejyYKGyY4A3dkyQGWtKTl1idanDlHNhQvcnE47F1Fhd6tAh1kq8CT8XrcTqHmc11tV4mfuXVbPdeNnH6knsClUEu0IREREREdnDrlBERERERORRfLDwAF0WB1mJ12WBleqFUt60aM9qrMzcOsTKzi07v7fOD9Vzw6T961CvLjXrUIOsWN4b9Lr/6jJu2pL+wSzFdFhjYXVRjpV4Wbll1qBDHSr3ryrOaqzM3DrEys4tO7+3zg/Vc8Ok/etQr9VYWbl1qEFWLO8Net1/dRk3T+Li7SJ0eLDQZXGQlXhdFlipXijlTYv2rMbKzK1DrOzcsvN76/xQPTdM2r8O9epSsw41yIrlvUGv+68u4+ZJXLxdBBdvExERERHZw8XbRERERETkUXywICIiIiKiMuODhYeUZaW/SV0mVNcruwODN3bbsBJryjwwKadu8arnuOr9y6jBpLHSLV6n3LLym5JT1v49EW91G1U/D9kVygA6LN4Womwr/a1s62qs1Xpk1CCjXpnHJasOk86FKfPApJy6xaue46r3L6MGk8ZKt3idcsvKb0pOWfv3RLzVbVyJ1WHcPIVdoYrQ5cGiLCv9Teoyobpe2R0YvLHbhpVYU+aBSTl1i1c9x1XvX0YNJo2VbvE65ZaV35ScsvbviXir26j6eciuUAZgVygiIiIiInvYFYqIiIiIiDyKDxZERERERFRmfLDwIJVdd2TmltntQHUHCpO6VrlzfnnjeKqe06pjVeeUVYPMmk25/lXPV5l16JBXZqydeJ3q8eb7iupYYztESV/xoZgui7eFUNt1R2ZuV2Ot1iszt7vzqq7T1VjVdaocT9VzWnWs6pyyarAar3IcVB+TzPGXVYcOeWXG2onXqR5vvq+ojrUzL2QxpivU8uXLxS233CKqVasmAIh58+Y5fb+goEA899xzIjY2VgQFBYn27duL3377zdI+dHqwUNl1R2Zumd0OVHegMKlrlTvnlzeOp+o5rTpWdU5ZNcis2ZTrX/V8lVmHDnllxtqJ16keb76vqI7VqUOUMV2hFi5ciFWrViE1NRXdunXDvHnzcNtttzm+//LLL2P8+PGYPn06EhMT8dxzz2HLli3Yvn07goKCXNoHu0IREREREdlj5Xdpfw/VVKL09HSkp6eX+D0hBCZOnIhnn30Wt956KwDgww8/RExMDL744gvceeedniyViIiIiIguQ9vF2/v27UNWVhbS0tIcr0VERKB58+ZYvXp1qdvl5uYiJyfH6Us1dy7UUbWY10qsLguvdKnZanx5yG1aHbrUa3UbjpteC2NVx3pjTtW16pZbdrxptZg0P+3k1pLsz2W5CpessVi1apUAIA4dOuQU17NnT9GrV69S84waNUoAKPalco2FqwtwXIlzNZeqfFZzyqrBarzMmq3Gl4fcptWhS71Wt+G46bUwVnWsN+ZUXatuuWXHm1aLSfPTTm5PMWbxdlHuerA4f/68yM7OdnwdOHBA+YOFOxfqqFrMayVWl4VXutRsNb485DatDl3qtboNx02vhbGqY70xp+padcstO960Wkyan3Zye4oxi7eL8vHxcVq8vXfvXiQlJWHjxo1ITk52xLVt2xbJycmYNGmSS3m5eJuIiIiIyB4rv0tru8YiMTERsbGxWLJkieO1nJwcrFmzBi1btlRYGRERERERXUppV6gzZ85g9+7djv/et28fNm3ahKioKMTHx2Po0KF44YUXULt2bUe72bi4OKeWtEREREREpJ7SdyzWrVuHlJQUpKSkAACGDRuGlJQUjBw5EgAwfPhwDBkyBPfffz+aNm2KM2fO4JtvvnH5b1jowqSuBKq7a6iqU/XxyIr1xv2b1pXEaqwO9epQA2PtxcrOLXsfqudzeY81dV6amNsrOkCVRPqKD8V0+MvbMjoCWMkpI6+MY1JZp+rjkRXrjfuXNfftxJt2fCbVwFh7sbJzy96H6vlc3mNNnZcm5rZ77algZFcoWXR4sDCpK4Hq7hqq6lR9PLJivXH/pnUlsRqrQ7061MBYe7Gyc8veh+r5XN5jTZ2XJubWtQNUSYzsCiULu0IREREREdnjFV2hiIiIiIjIHHywICIiIiKiMuODhQeZ0MnIlC4bpuS0Eq9Dpwp2DdKrDh3Oh4w6TNq/6lpl1WDCvVb349ClBh1iXY1XXavq4y9LvDGkr/hQTIfF24Vc7QDg7jhT9q16/7KOydV4GTmtxsoYK1m1yozVpQ4dzoeMOkzav+paZdVgwr1W9+PQpQYdYl2NV12r6uMvS7xK7ApVhE4PFiZ0MjKly4YpOa3E69Cpgl2D9KpDh/Mhow6T9q+6Vlk1mHCv1f04dKlBh1hX41XXqvr4yxKvErtCFcGuUERERERE9rArFBEREREReRQfLDxE9cIi0xZByapDZrzKek07vzqcW1POg+zcpsXKzK/D+fOWOJNylqd41WPt7ljZ9w4ZNcvKqQ3pH8xSTJc1Fq4u0nE1TodYWfuXWYfMeJX1mnZ+dTi3ppwH2blNi5WZX4fz5y1xJuUsT/Gqx9rdsbLvHVbjVY+ZLFy8XYQuDxaqFxaZtghKVh0y41XWa9r51eHcmnIeZOc2LVZmfh3On7fEmZSzPMWrHmt3x8q+d8ioWVZOmbh4uwgu3iYiIiIisoeLt4mIiIiIyKP4YEFERERERGXGBwsP0aFTgLs6LJjWUUGHGkyq16RuJDJjrcSbML9M6A6kWx06zENPdI+R3VnHG+89ptyjVY+9LrEy41Xdq7UlfcWHYros3nZ1Vb+rcVZjXY13V4zdeBnjpEMNJtUrq1ZZ46XDeTBhfqmcJ1bHXZc6dJiHssfZ7naq57ysvCbkVH08smqQGSsz3t01272OZWJXqCJ0ebDQoVOAuzosmNZRQYcaTKrXpG4kMmOtxJswv0zoDqRbHTrMQ090j5HdWccb7z2m3KNVj70usTLjVd2rPYldoYpgVygiIiIiInvYFYqIiIiIiDyKDxZERERERFRmfLDwAB06IsjqDKFLblmx3pjTU9vo0t1Dh1i726jsJmNnG9UdaHSZR56qSbd4HXKr7NCj8meADvXK6nik8h6u+t5jZIco6Ss+FNNh8baVFf6qY63k1Cm3rFhvzOmpbWTF6lKHJ8bc3edfl3OuehzsxNsZO9k16RavQ2535zShRl3qdfd42tlG1fnX5Z4mC7tCFaHDg4UOHRFkdYbQJbesWG/M6altdOnuoUOs3W1UdpOxs43qDjS6zCNP1aRbvA65VXboUfkzQId6ZXU8UnkPV33v0aVDFLtCFcGuUERERERE9rArFBEREREReRQfLDxEhwWkuiyaMrEO1bGq9281trzkthLvredah3Nh0jVndxud6tHhnNuN16kmGYuddZpXOuRWXYPs61ZL0j+YpZgOayyE0GMBqZVYmflNrEN1rOr9W40tL7mtxHvrudbhXJh0zdndRqd6dDjnduN1qsmVOF1qNTW36hpkX7eewsXbRejyYKHDAlJdFk2ZWIfqWNX7txpbXnJbiffWc63DuTDpmrO7jU716HDO7cbrVJOMxc46zSsdcquuQfZ16ylcvF0EF28TEREREdnDxdtERERERORRfLAgIiIiIqIy44OFh3hT5xjV3UNMGCMrsTrUaSfeyjY6HKMOHTxM6l5k0th621jJyOktNXrjufZEvKm1yDjfquuVdUzakL7iQzFdFm+7utrf3XEyYlXuW0ac6lgd6rQTb2UbHY5R5bmQmZ9j631jJSOnt9TojefaE/Gm1iLjfKuuV9YxycSuUEXo8mDhTZ1jVHcPMWGMrMTqUKedeCvb6HCMOnTwMKl7kUlj621jJSOnt9TojefaE/Gm1iLjfKuuV9YxycSuUEWwKxQRERERkT3sCkVERERERB7FBwsiIiIiIiozPlh4iDd2PJJVr51tVI+bDvXKystOLeqOzaQ5qMP5UtUdyUqsSfvXoQbdc3rTsXgqXofjk5Fbl+NSTvqKD8V0Wbzt6up+d8fJjJVVr51tVI+bDvXKyutKrOrx1+G8mnT+ZdSgw/ly51yVFWvS/nWoQfec3nQsnorX4fhk5NbluGRgV6gidHmw8MaOR7LqtbON6nHToV5ZedmpRd2xmTQHdThfqrojWYk1af861KB7Tm86Fk/F63B8MnLrclwysCtUEewKRURERERkD7tCERERERGRR/HBwkNUL+qRvXja6jY6LAzzxJjI3ofqWG9a3Gg1Vod6dahDh2tZt1rKe6zq68hbrwdZdcioQfX+7ea2so3Ka0Fr0j+YpZguayysLMBxNVZGTk9tYzW36vGzu43sfaiOdSVO9bmTeZ2orleHOnS4lnWrpbzHqr6OvPV6kFWHjBpU799ubivbqLwWPM3rFm//5z//EQkJCSIwMFA0a9ZMrFmzxuVtdXmwUL2oR/biaavb6LAwzBNjInsfqmO9aXGj1Vgd6tWhDh2uZd1qKe+xqq8jb70eZNUhowbV+7eb28o2Kq8FT/OqxduffPIJMjIyMHXqVDRv3hwTJ07EnDlzsHPnTkRHR19xey7eJiIiIiKyx6sWb7/++uu477770L9/f9SvXx9Tp05FxYoV8cEHH6gujYiIiIiI/o/WDxYXLlzA+vXrkZaW5njN19cXaWlpWL16tcLKiIiIiIioKK0fLI4fP478/HzExMQ4vR4TE4OsrKwSt8nNzUVOTo7Tlw5M6eZh0v5N6qahS9cL3buveKKDhuqxtRpvUg2mja3VeNXHp3unp/J6L3JnrOo5pjrWhPOow1hqTfqKjzI4ePCgACB+/PFHp9efeOIJ0axZsxK3GTVqlABQ7Ev14m0rK/1Vx5qyfyv7lpVXda2q6zXlnNrdRlaslXiTajBtbK3Gqz4+lft3V4zVWN3j3B2reo6pjjXhPOowlp7mNV2hcnNzhZ+fn5g3b57T6xkZGaJr164lbnP+/HmRnZ3t+Dpw4IAWDxamdPMwaf8mddPQpeuF7t1XPNFBQ/XYWo03qQbTxtZqvOrj073TU3m9F7kzVvUcUx1rwnnUYSw9zau6QjVv3hzNmjXD5MmTAQAFBQWIj4/H4MGD8eSTT15xe3aFIiIiIiKyx8rv0v4eqsm2YcOGoW/fvmjSpAmaNWuGiRMn4u+//0b//v1Vl0ZERERERP9H+weLO+64A8eOHcPIkSORlZWF5ORkfPPNN8UWdBMRERERkTrafxSqrPhRKCIiIiIie7zqD+QREREREZH++GBBRERERERlxgcLIiIiIiIqMz5YEBERERFRmfHBgoiIiIiIyowPFkREREREVGZ8sCAiIiIiojLT/g/klVXhn+nIyclRXAkRERERkVkKf4d25U/fef2DxenTpwEANWrUUFwJEREREZGZTp8+jYiIiMvGeP1f3i4oKMChQ4cQFhYGHx8ft+fPyclBjRo1cODAAf5lbwLAOUHFcU5QSTgv6FKcE3QpHeaEEAKnT59GXFwcfH0vv4rC69+x8PX1xVVXXSV9P+Hh4bwJkBPOCboU5wSVhPOCLsU5QZdSPSeu9E5FIS7eJiIiIiKiMuODBRERERERlRkfLMooMDAQo0aNQmBgoOpSSBOcE3QpzgkqCecFXYpzgi5l2pzw+sXbREREREQkH9+xICIiIiKiMuODBRERERERlRkfLIiIiIiIqMz4YFEGb731FmrWrImgoCA0b94cP//8s+qSSJLx48ejadOmCAsLQ3R0NG677Tbs3LnTKeb8+fMYNGgQKleujNDQUHTv3h1Hjhxxitm/fz86d+6MihUrIjo6Gk888QQuXrzoyUMhSV566SX4+Phg6NChjtc4J8qngwcP4u6770blypURHByMRo0aYd26dY7vCyEwcuRIVKtWDcHBwUhLS8OuXbuccpw4cQJ9+vRBeHg4IiMjMWDAAJw5c8bTh0JukJ+fj+eeew6JiYkIDg5GUlISxo4di6JLXDknvNsPP/yALl26IC4uDj4+Pvjiiy+cvu+u8//LL7+gdevWCAoKQo0aNfDKK6/IPrTiBNkye/ZsERAQID744AOxbds2cd9994nIyEhx5MgR1aWRBB07dhTTpk0TW7duFZs2bRI333yziI+PF2fOnHHEPPjgg6JGjRpiyZIlYt26daJFixbi+uuvd3z/4sWLomHDhiItLU1s3LhRfP3116JKlSriqaeeUnFI5EY///yzqFmzpmjcuLF45JFHHK9zTpQ/J06cEAkJCaJfv35izZo1Yu/evWLRokVi9+7djpiXXnpJREREiC+++EJs3rxZdO3aVSQmJopz5845Yjp16iSuvfZa8dNPP4kVK1aIq6++WvTu3VvFIVEZvfjii6Jy5cpi/vz5Yt++fWLOnDkiNDRUTJo0yRHDOeHdvv76a/HMM8+IuXPnCgBi3rx5Tt93x/nPzs4WMTExok+fPmLr1q1i1qxZIjg4WPz3v//11GEKIYTgg4VNzZo1E4MGDXL8d35+voiLixPjx49XWBV5ytGjRwUAsXz5ciGEEKdOnRIVKlQQc+bMccT8+uuvAoBYvXq1EOKfG4uvr6/IyspyxLz99tsiPDxc5ObmevYAyG1Onz4tateuLRYvXizatm3reLDgnCifRowYIW644YZSv19QUCBiY2PFhAkTHK+dOnVKBAYGilmzZgkhhNi+fbsAINauXeuIWbhwofDx8REHDx6UVzxJ0blzZ3Hvvfc6vdatWzfRp08fIQTnRHlz6YOFu87/lClTRKVKlZx+dowYMULUqVNH8hE540ehbLhw4QLWr1+PtLQ0x2u+vr5IS0vD6tWrFVZGnpKdnQ0AiIqKAgCsX78eeXl5TnOibt26iI+Pd8yJ1atXo1GjRoiJiXHEdOzYETk5Odi2bZsHqyd3GjRoEDp37ux07gHOifLqf//7H5o0aYKePXsiOjoaKSkpePfddx3f37dvH7KyspzmRUREBJo3b+40LyIjI9GkSRNHTFpaGnx9fbFmzRrPHQy5xfXXX48lS5bgt99+AwBs3rwZK1euRHp6OgDOifLOXed/9erVaNOmDQICAhwxHTt2xM6dO3Hy5EkPHQ3g77E9eZHjx48jPz/f6ZcBAIiJicGOHTsUVUWeUlBQgKFDh6JVq1Zo2LAhACArKwsBAQGIjIx0io2JiUFWVpYjpqQ5U/g9Ms/s2bOxYcMGrF27ttj3OCfKp7179+Ltt9/GsGHD8PTTT2Pt2rV4+OGHERAQgL59+zrOa0nnvei8iI6Odvq+v78/oqKiOC8M9OSTTyInJwd169aFn58f8vPz8eKLL6JPnz4AwDlRzrnr/GdlZSExMbFYjsLvVapUSUr9l+KDBZFFgwYNwtatW7Fy5UrVpZBCBw4cwCOPPILFixcjKChIdTmkiYKCAjRp0gTjxo0DAKSkpGDr1q2YOnUq+vbtq7g6UuHTTz/FjBkzMHPmTDRo0ACbNm3C0KFDERcXxzlBXocfhbKhSpUq8PPzK9bd5ciRI4iNjVVUFXnC4MGDMX/+fCxduhRXXXWV4/XY2FhcuHABp06dcoovOidiY2NLnDOF3yOzrF+/HkePHsV1110Hf39/+Pv7Y/ny5XjzzTfh7++PmJgYzolyqFq1aqhfv77Ta/Xq1cP+/fsB/P/zermfH7GxsTh69KjT9y9evIgTJ05wXhjoiSeewJNPPok777wTjRo1wj333INHH30U48ePB8A5Ud656/zr8vOEDxY2BAQEIDU1FUuWLHG8VlBQgCVLlqBly5YKKyNZhBAYPHgw5s2bh++//77Y242pqamoUKGC05zYuXMn9u/f75gTLVu2xJYtW5xuDosXL0Z4eHixX0RIf+3bt8eWLVuwadMmx1eTJk3Qp08fx//nnCh/WrVqVawV9W+//YaEhAQAQGJiImJjY53mRU5ODtasWeM0L06dOoX169c7Yr7//nsUFBSgefPmHjgKcqezZ8/C19f51y0/Pz8UFBQA4Jwo79x1/lu2bIkffvgBeXl5jpjFixejTp06HvsYFAC2m7Vr9uzZIjAwUGRmZort27eL+++/X0RGRjp1dyHv8dBDD4mIiAixbNkycfjwYcfX2bNnHTEPPvigiI+PF99//71Yt26daNmypWjZsqXj+4WtRf/1r3+JTZs2iW+++UZUrVqVrUW9SNGuUEJwTpRHP//8s/D39xcvvvii2LVrl5gxY4aoWLGi+Pjjjx0xL730koiMjBRffvml+OWXX8Stt95aYmvJlJQUsWbNGrFy5UpRu3ZtthY1VN++fUX16tUd7Wbnzp0rqlSpIoYPH+6I4ZzwbqdPnxYbN24UGzduFADE66+/LjZu3Cj++OMPIYR7zv+pU6dETEyMuOeee8TWrVvF7NmzRcWKFdlu1iSTJ08W8fHxIiAgQDRr1kz89NNPqksiSQCU+DVt2jRHzLlz58TAgQNFpUqVRMWKFcXtt98uDh8+7JTn999/F+np6SI4OFhUqVJFPPbYYyIvL8/DR0OyXPpgwTlRPn311VeiYcOGIjAwUNStW1e88847Tt8vKCgQzz33nIiJiRGBgYGiffv2YufOnU4xf/31l+jdu7cIDQ0V4eHhon///uL06dOePAxyk5ycHPHII4+I+Ph4ERQUJGrVqiWeeeYZp7agnBPebenSpSX+DtG3b18hhPvO/+bNm8UNN9wgAgMDRfXq1cVLL73kqUN08BGiyJ9+JCIiIiIisoFrLIiIiIiIqMz4YEFERERERGXGBwsiIiIiIiozPlgQEREREVGZ8cGCiIiIiIjKjA8WRERERERUZnywICIiIiKiMuODBRERERERlRkfLIiIqJhly5bBx8cHp06dUrL/JUuWoF69esjPz5e2jxYtWuDzzz+Xlp+IqLzhX94mIirnbrzxRiQnJ2PixImO1y5cuIATJ04gJiYGPj4+Hq8pNTUVw4YNQ58+faTtY/78+Xj00Uexc+dO+Pry39mIiMqKd1IiIiomICAAsbGxSh4qVq5ciT179qB79+5S95Oeno7Tp09j4cKFUvdDRFRe8MGCiKgc69evH5YvX45JkybBx8cHPj4++P3334t9FCozMxORkZGYP38+6tSpg4oVK6JHjx44e/Yspk+fjpo1a6JSpUp4+OGHnT6+lJubi8cffxzVq1dHSEgImjdvjmXLll22ptmzZ6NDhw4ICgpyvDZ69GgkJyfjgw8+QHx8PEJDQzFw4EDk5+fjlVdeQWxsLKKjo/Hiiy86thFCYPTo0YiPj0dgYCDi4uLw8MMPO77v5+eHm2++GbNnz3bPYBIRlXP+qgsgIiJ1Jk2ahN9++w0NGzbE888/DwCoWrUqfv/992KxZ8+exZtvvonZs2fj9OnT6NatG26//XZERkbi66+/xt69e9G9e3e0atUKd9xxBwBg8ODB2L59O2bPno24uDjMmzcPnTp1wpYtW1C7du0Sa1qxYgXuuuuuYq/v2bMHCxcuxDfffIM9e/agR48e2Lt3L6655hosX74cP/74I+69916kpaWhefPm+Pzzz/HGG29g9uzZaNCgAbKysrB582annM2aNcNLL71UxlEkIiKADxZEROVaREQEAgICULFiRcTGxl42Ni8vD2+//TaSkpIAAD169MBHH32EI0eOIDQ0FPXr10e7du2wdOlS3HHHHdi/fz+mTZuG/fv3Iy4uDgDw+OOP45tvvsG0adMwbty4Evfzxx9/OOKLKigowAcffICwsDDHvnbu3Imvv/4avr6+qFOnDl5++WUsXboUzZs3x/79+xEbG4u0tDRUqFAB8fHxaNasmVPOuLg4HDhwAAUFBVxnQURURryLEhGRSypWrOh4qACAmJgY1KxZE6GhoU6vHT16FACwZcsW5Ofn45prrkFoaKjja/ny5dizZ0+p+zl37pzTx6AK1axZE2FhYU77ql+/vtMDQdH99+zZE+fOnUOtWrVw3333Yd68ebh48aJTzuDgYBQUFCA3N9fiaBAR0aX4jgUREbmkQoUKTv/t4+NT4msFBQUAgDNnzsDPzw/r16+Hn5+fU1zRh5FLValSBSdPnizz/mvUqIGdO3fiu+++w+LFizFw4EBMmDABy5cvd2x34sQJhISEIDg4+HKHTkRELuCDBRFRORcQECDl70WkpKQgPz8fR48eRevWrS1tt337drfUEBwcjC5duqBLly4YNGgQ6tatiy1btuC6664DAGzduhUpKSlu2RcRUXnHBwsionKuZs2aWLNmDX7//XeEhoYiKirKLXmvueYa9OnTBxkZGXjttdeQkpKCY8eOYcmSJWjcuDE6d+5c4nYdO3bE9OnTy7z/zMxM5Ofno3nz5qhYsSI+/vhjBAcHIyEhwRGzYsUK/Otf/yrzvoiIiGssiIjKvccffxx+fn6oX78+qlativ3797st97Rp05CRkYHHHnsMderUwW233Ya1a9ciPj6+1G369OmDbdu2YefOnWXad2RkJN599120atUKjRs3xnfffYevvvoKlStXBgAcPHgQP/74I/r371+m/RAR0T/4l7eJiEg7TzzxBHJycvDf//5X2j5GjBiBkydP4p133pG2DyKi8oTvWBARkXaeeeYZJCQkOBZiyxAdHY2xY8dKy09EVN7wHQsiIiIiIiozvmNBRERERERlxgcLIiIiIiIqMz5YEBERERFRmfHBgoiIiIiIyowPFkREREREVGZ8sCAiIiIiojLjgwUREREREZUZHyyIiIiIiKjM+GBBRERERERlxgcLIiIiIiIqs/8HP6Yre9N8loQAAAAASUVORK5CYII=",
      "text/plain": [
       "<Figure size 800x400 with 1 Axes>"
      ]
     },
     "metadata": {},
     "output_type": "display_data"
    }
   ],
   "source": [
    "res = sim.simulate(1000. * u.ms)\n",
    "erate = res.rate(esr.segments[0].population)\n",
    "irate = res.rate(isr.segments[0].population)\n",
    "print(f'excitatory rate: {erate:.2f} spks/s')\n",
    "print(f'inhibitory rate: {irate:.2f} spks/s')\n",
    "\n",
    "spk = np.asarray(res.spikes(esr.segments[0].population))\n",
    "ts, ids = np.nonzero(spk > 0)\n",
    "plt.figure(figsize=(8, 4))\n",
    "plt.scatter(ts * 0.1, ids, s=1.0, color='k')\n",
    "plt.xlabel('time (ms)'); plt.ylabel('exc neuron')\n",
    "plt.title('Brunel (delta) — excitatory raster, brainpy.state port')\n",
    "plt.tight_layout()"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "8c642248",
   "metadata": {},
   "source": [
    "## Parameter mapping\n",
    "\n",
    "| NEST | brainpy.state | Note |\n",
    "|------|---------------|------|\n",
    "| `nest.SetKernelStatus({'resolution': 0.1})` | `Simulator(dt=0.1 * u.ms)` | per-`Simulator`, no global kernel |\n",
    "| `nest.Create('iaf_psc_delta', NE)` | `sim.create(bp.iaf_psc_delta, NE, params=...)` | returns a `NodeView` |\n",
    "| `'C_m': 1.0` (pF implied) | `C_m=1.0 * u.pF` | units are explicit |\n",
    "| `nodes_ex + nodes_in` | `ne + ni` | `NodeView` concatenation |\n",
    "| `nodes_ex[:N_rec]` | `ne[:N_rec]` | `NodeView` slice |\n",
    "| `{'rule': 'fixed_indegree', 'indegree': CE}` | `rule=bp.fixed_indegree(CE)` | same draw; `seed=` to reproduce |\n",
    "| `syn_spec={'weight': J, 'delay': d}` | `weight=J * u.mV, delay=d * u.ms` | delta weight is a voltage (`mV`) |\n",
    "| `nest.Simulate(T)` | `sim.simulate(T * u.ms)` | returns a `SimulationResult` |\n",
    "| `nest.GetStatus(sr, 'n_events')` | `res.rate(...)` / `res.n_events(...)` | read from the result |"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "0b589181",
   "metadata": {},
   "source": [
    "## Divergences this port runs into\n",
    "\n",
    "- **Delta weights are voltages.** With `iaf_psc_delta` the postsynaptic amplitude is an instantaneous membrane-voltage jump, so weights are `mV` (`J_ex = J`), not the `pA` currents the alpha variant uses.\n",
    "- **Independent PRNG streams.** The Poisson drive and the random `fixed_indegree` wiring draw from JAX's PRNG, which diverges sample-by-sample from NEST's. Parity for Brunel is therefore **distributional** — the seed-mean firing rate matches within a 5% band (category D), not spike-by-spike. See [validation status](validation-status.rst).\n",
    "- **Sparse communication.** `comm='sparse'` selects the event-driven backend so large fan-in stays memory-light; `allow_multapses=True` matches NEST's `fixed_indegree` default (repeated edges sum).\n",
    "\n",
    "For the full catalog — where learning state lives, parameter-location maps, the STDP pairing conventions, and the documented numerical bands — see the [semantic divergences](divergences/index.rst)."
   ]
  },
  {
   "cell_type": "markdown",
   "id": "76689773",
   "metadata": {},
   "source": [
    "## See also\n",
    "\n",
    "- [Connect a network](tutorials/03-connect-network.ipynb) — the same Brunel built from scratch.\n",
    "- [Semantic divergences](divergences/index.rst) — the full porting catalog.\n",
    "- [Validation status](validation-status.rst) — how Brunel parity is asserted (distributional, category D).\n",
    "- [Example gallery](../examples/nest-gallery.rst) — full-scale Brunel and other ported networks."
   ]
  }
 ],
 "metadata": {
  "kernelspec": {
   "display_name": "Python 3",
   "language": "python",
   "name": "python3"
  },
  "language_info": {
   "codemirror_mode": {
    "name": "ipython",
    "version": 3
   },
   "file_extension": ".py",
   "mimetype": "text/x-python",
   "name": "python",
   "nbconvert_exporter": "python",
   "pygments_lexer": "ipython3",
   "version": "3.13.11"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 5
}
