{ "cells": [ { "cell_type": "markdown", "id": "bf560305-a936-46c2-82eb-ec0e7cb6e218", "metadata": {}, "source": [ "# wsssss introduction\n", "\n", "This notebook will show some of the possibilities of `wsssss`. The documentation is available [here.](https://wsssss.readthedocs.io/latest/wsssss.html)\n", "\n", "First we create a grid of MESA inlists, then run that grid with MESA and gyre. We first ensure that either the `MESA_DIR` environment variable is set or that we pass `mesa_dir` to `MesaGrid`, as it needs to know which version of MESA to create a grid for. We also tell `MesaGrid` to copy the basic files and folders needed for MESA (rn, re, mk, etc.). Next we set some inputs for the grid. We create a grid with 2 initial masses and 2 different overshooting, resulting in 4 sets of inlists." ] }, { "cell_type": "code", "execution_count": 1, "id": "623384ca-7847-4009-a179-3e39ac9d0f77", "metadata": {}, "outputs": [], "source": [ "import os\n", "\n", "wsssss_dir = \"/home/walter/Github/wsssss\"\n", "if 'MESA_DIR' not in os.environ.keys():\n", " os.environ['MESA_DIR'] = '/home/walter/Software/MESA/mesa-24.08.1'\n", "mesa_dir = os.environ['MESA_DIR']\n" ] }, { "cell_type": "code", "execution_count": 2, "id": "initial_id", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:45:18.689441467Z", "start_time": "2026-03-19T15:45:17.621793228Z" } }, "outputs": [], "source": [ "from wsssss.inlists import create_grid as cg\n", "\n", "# If mesa_dir is not set it will use the MESA_DIR environment variable\n", "grid = cg.MesaGrid(mesa_dir=mesa_dir, add_base_workdir=True)\n", "\n", "grid.star_job['history_columns_file'] = 'history_test.list'\n", "# Tell MesaGrid where the file is located:\n", "grid.add_file(os.path.join(f'{wsssss_dir}/tests/make_data', grid.star_job['history_columns_file']))\n", "\n", "grid.star_job['pgstar_flag'] = False\n", "grid.star_job['create_pre_main_sequence_model'] = False\n", "\n", "grid.kap['use_Type2_opacities'] = True\n", "grid.kap['Zbase'] = 0.02\n", "\n", "grid.controls['initial_mass'] = [1, 2]\n", "grid.controls['initial_z'] = 0.02\n", "grid.controls['max_model_number'] = 10\n", "grid.controls['mesh_delta_coeff'] = 5\n", "\n", "grid.controls['history_interval'] = 1\n", "grid.controls['photo_interval'] = 8\n", "\n", "grid.controls['write_profiles_flag'] = True\n", "grid.controls['profile_interval'] = 9\n", "grid.controls['write_pulse_data_with_profile'] = True\n", "grid.controls['pulse_data_format'] = 'GYRE'\n", "\n", "grid.controls[f'{cg.non_mesa_key_start}group_unpack'].append([{'overshoot_scheme(1)': 'exponential',\n", " 'overshoot_zone_type(1)': 'nonburn',\n", " 'overshoot_zone_loc(1)': 'shell',\n", " 'overshoot_bdy_loc(1)': 'bottom',\n", " 'overshoot_f(1)': 0.02,\n", " 'overshoot_f0(1)': 0.001},\n", " {'overshoot_scheme(1)': 'exponential',\n", " 'overshoot_zone_type(1)': 'nonburn',\n", " 'overshoot_zone_loc(1)': 'shell',\n", " 'overshoot_bdy_loc(1)': 'bottom',\n", " 'overshoot_f(1)': 0.03,\n", " 'overshoot_f0(1)': 0.01}])\n", "\n", "\n", "def name_function(unpacked):\n", " return f'm{unpacked[\"controls\"][\"initial_mass\"]:.3f}_fOV{unpacked[\"controls\"][\"overshoot_f(1)\"]:.3f}'\n", "\n", "\n", "grid.set_name_function(name_function)\n", "\n", "# If we want, we can apply an arbritraty function on each unpacked inlist with grid.set_inlist_finalize_function\n", "\n", "def finalize_function(unpacked):\n", " if unpacked[\"controls\"]['overshoot_f(1)'] == 0.03:\n", " unpacked[\"controls\"][\"overshoot_D_min\"] = 1e-2 \n", " # We could also set the history file output name to something other than the default:\n", " # unpacked[\"controls\"][\"star_history_name\"] = f'm{unpacked[\"controls\"][\"initial_mass\"]:.3f}_fOV{unpacked[\"controls\"][\"overshoot_f(1)\"]:.3f}.data'\n", " return unpacked\n", "\n", "grid.set_inlist_finalize_function(finalize_function)" ] }, { "cell_type": "code", "execution_count": 3, "id": "22723f2b074a2bd9", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:45:18.691069506Z", "start_time": "2026-03-19T15:45:18.685900082Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "0000\n", " inlist\n", " star_job:\n", " read_extra_star_job_inlist(5): True\n", " extra_star_job_inlist_name(5): inlist_project\n", " eos:\n", " read_extra_eos_inlist(5): True\n", " extra_eos_inlist_name(5): inlist_project\n", " kap:\n", " read_extra_kap_inlist(5): True\n", " extra_kap_inlist_name(5): inlist_project\n", " controls:\n", " read_extra_controls_inlist(5): True\n", " extra_controls_inlist_name(5): inlist_project\n", " pgstar:\n", " read_extra_pgstar_inlist(5): True\n", " extra_pgstar_inlist_name(5): inlist_project\n", " star_job\n", " history_columns_file: history_test.list\n", " pgstar_flag: False\n", " create_pre_main_sequence_model: False\n", " eos\n", " kap\n", " use_Type2_opacities: True\n", " Zbase: 0.02\n", " controls\n", " initial_mass: 1\n", " initial_z: 0.02\n", " max_model_number: 10\n", " mesh_delta_coeff: 5\n", " history_interval: 1\n", " photo_interval: 8\n", " write_profiles_flag: True\n", " profile_interval: 9\n", " write_pulse_data_with_profile: True\n", " pulse_data_format: GYRE\n", " overshoot_scheme(1): exponential\n", " overshoot_zone_type(1): nonburn\n", " overshoot_zone_loc(1): shell\n", " overshoot_bdy_loc(1): bottom\n", " overshoot_f(1): 0.02\n", " overshoot_f0(1): 0.001\n", " pgstar\n", "\n", "0001\n", " inlist\n", " star_job:\n", " read_extra_star_job_inlist(5): True\n", " extra_star_job_inlist_name(5): inlist_project\n", " eos:\n", " read_extra_eos_inlist(5): True\n", " extra_eos_inlist_name(5): inlist_project\n", " kap:\n", " read_extra_kap_inlist(5): True\n", " extra_kap_inlist_name(5): inlist_project\n", " controls:\n", " read_extra_controls_inlist(5): True\n", " extra_controls_inlist_name(5): inlist_project\n", " pgstar:\n", " read_extra_pgstar_inlist(5): True\n", " extra_pgstar_inlist_name(5): inlist_project\n", " star_job\n", " history_columns_file: history_test.list\n", " pgstar_flag: False\n", " create_pre_main_sequence_model: False\n", " eos\n", " kap\n", " use_Type2_opacities: True\n", " Zbase: 0.02\n", " controls\n", " initial_mass: 2\n", " initial_z: 0.02\n", " max_model_number: 10\n", " mesh_delta_coeff: 5\n", " history_interval: 1\n", " photo_interval: 8\n", " write_profiles_flag: True\n", " profile_interval: 9\n", " write_pulse_data_with_profile: True\n", " pulse_data_format: GYRE\n", " overshoot_scheme(1): exponential\n", " overshoot_zone_type(1): nonburn\n", " overshoot_zone_loc(1): shell\n", " overshoot_bdy_loc(1): bottom\n", " overshoot_f(1): 0.02\n", " overshoot_f0(1): 0.001\n", " pgstar\n", "\n", "0002\n", " inlist\n", " star_job:\n", " read_extra_star_job_inlist(5): True\n", " extra_star_job_inlist_name(5): inlist_project\n", " eos:\n", " read_extra_eos_inlist(5): True\n", " extra_eos_inlist_name(5): inlist_project\n", " kap:\n", " read_extra_kap_inlist(5): True\n", " extra_kap_inlist_name(5): inlist_project\n", " controls:\n", " read_extra_controls_inlist(5): True\n", " extra_controls_inlist_name(5): inlist_project\n", " pgstar:\n", " read_extra_pgstar_inlist(5): True\n", " extra_pgstar_inlist_name(5): inlist_project\n", " star_job\n", " history_columns_file: history_test.list\n", " pgstar_flag: False\n", " create_pre_main_sequence_model: False\n", " eos\n", " kap\n", " use_Type2_opacities: True\n", " Zbase: 0.02\n", " controls\n", " initial_mass: 1\n", " initial_z: 0.02\n", " max_model_number: 10\n", " mesh_delta_coeff: 5\n", " history_interval: 1\n", " photo_interval: 8\n", " write_profiles_flag: True\n", " profile_interval: 9\n", " write_pulse_data_with_profile: True\n", " pulse_data_format: GYRE\n", " overshoot_scheme(1): exponential\n", " overshoot_zone_type(1): nonburn\n", " overshoot_zone_loc(1): shell\n", " overshoot_bdy_loc(1): bottom\n", " overshoot_f(1): 0.03\n", " overshoot_f0(1): 0.01\n", " overshoot_D_min: 0.01\n", " pgstar\n", "\n", "0003\n", " inlist\n", " star_job:\n", " read_extra_star_job_inlist(5): True\n", " extra_star_job_inlist_name(5): inlist_project\n", " eos:\n", " read_extra_eos_inlist(5): True\n", " extra_eos_inlist_name(5): inlist_project\n", " kap:\n", " read_extra_kap_inlist(5): True\n", " extra_kap_inlist_name(5): inlist_project\n", " controls:\n", " read_extra_controls_inlist(5): True\n", " extra_controls_inlist_name(5): inlist_project\n", " pgstar:\n", " read_extra_pgstar_inlist(5): True\n", " extra_pgstar_inlist_name(5): inlist_project\n", " star_job\n", " history_columns_file: history_test.list\n", " pgstar_flag: False\n", " create_pre_main_sequence_model: False\n", " eos\n", " kap\n", " use_Type2_opacities: True\n", " Zbase: 0.02\n", " controls\n", " initial_mass: 2\n", " initial_z: 0.02\n", " max_model_number: 10\n", " mesh_delta_coeff: 5\n", " history_interval: 1\n", " photo_interval: 8\n", " write_profiles_flag: True\n", " profile_interval: 9\n", " write_pulse_data_with_profile: True\n", " pulse_data_format: GYRE\n", " overshoot_scheme(1): exponential\n", " overshoot_zone_type(1): nonburn\n", " overshoot_zone_loc(1): shell\n", " overshoot_bdy_loc(1): bottom\n", " overshoot_f(1): 0.03\n", " overshoot_f0(1): 0.01\n", " overshoot_D_min: 0.01\n", " pgstar\n" ] } ], "source": [ "grid.unpack_inlists()\n", "\n", "for i, inlist in enumerate(grid.unpacked):\n", " print(f'{i:04}')\n", " for namelist in inlist.keys():\n", " if isinstance(namelist, str):\n", " print(f' {namelist}')\n", " for key, value in inlist[namelist].items():\n", " if not isinstance(value, dict):\n", " if not key.startswith(cg.non_mesa_key_start):\n", " print(f' {key}: {value}')\n", " else:\n", " print(f' {key}:')\n", " for kkey, vvalue in value.items():\n", " if not key.startswith(cg.non_mesa_key_start):\n", " print(f' {kkey}: {vvalue}')\n", " print()\n" ] }, { "cell_type": "code", "execution_count": 4, "id": "48450caff6f15486", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:45:18.832519842Z", "start_time": "2026-03-19T15:45:18.686405636Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Cannot copy calling file /tmp/ipykernel_409779/2927543967.py to grid directory /home/walter/Github/wsssss/examples/example_wsssss.\n" ] } ], "source": [ "grid_dir = os.path.abspath('example_wsssss')\n", "grid.create_grid(grid_dir, rm_dir=False)" ] }, { "cell_type": "code", "execution_count": 5, "id": "f709b068bca739b8", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:45:18.850812023Z", "start_time": "2026-03-19T15:45:18.837039244Z" }, "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "['m2.000_fOV0.020',\n", " 'star',\n", " 'm1.000_fOV0.020',\n", " 'out_m2.000_fOV0.030',\n", " 'm1.000_fOV0.030',\n", " 'out_m1.000_fOV0.030',\n", " 'out_m1.000_fOV0.020',\n", " 'out_m2.000_fOV0.020',\n", " 'm2.000_fOV0.030']" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "os.listdir(grid_dir)" ] }, { "cell_type": "markdown", "id": "2a7d5a21d48090cb", "metadata": { "collapsed": false }, "source": [ "# Run `mesa-go`\n", "\n", "First we compile `star` and copy it to the grid directory as we will copy it instead of compile it for each run of MESA. Next, run `mesa-go` ([docs](https://wsssss.readthedocs.io/latest/cmd.html#mesa-grid-overseer)), the MESA grid runner. We tell it that each run has a non-default name with `--sub-dirs m*`, and that before starting MESA to run `cp ../star ./` to copy the `star` executable we created earlier, and after running MESA to run the `gyre-driver` ([docs](https://wsssss.readthedocs.io/latest/cmd.html#gyre-driver)) on all the profiles in LOGS and to be lenient with which version of gyre is used.\n", "```\n", "cd /path/to/example_wsssss\n", "cd m1.000_fOV0.020; ./mk; cp star ../; cd ..\n", "mesa-go --sub-dirs m* --cmd-pre-each 'cp ../star ./' --cmd-post-each 'gyre-driver 01 MESA LOGS/*.GYRE --lenient'\n", "check-grid --no-slurm --list-all\n", "```\n", "Once `mesa-go` is finished running we check if the grid ran successfully with `check-grid` ([docs](https://wsssss.readthedocs.io/latest/cmd.html#mesa-grid-check)). This should show the following:\n", "```\n", "\n", " termination_code count\n", "--------------------------------------------\n", "max_model_number 4\n", "--------------------------------------------\n", "m1.000_fOV0.021 pre-CHeX 9 max_model_number\n", "m1.000_fOV0.030 pre-CHeX 9 max_model_number\n", "m2.000_fOV0.021 pre-CHeX 9 max_model_number\n", "m2.000_fOV0.030 pre-CHeX 9 max_model_number\n", "--------------------------------------------\n", "```\n" ] }, { "cell_type": "markdown", "id": "3a00b989-05d4-4c08-89d6-b213742217e4", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:48:36.993403096Z", "start_time": "2026-03-19T15:48:36.950260680Z" }, "collapsed": false }, "source": [ "# Inlists\n", "\n", "We can also evaluate an inlist and see which options it sets as well as compare different inlists. We also see that `overshoot_d_min` is not set in `m1.000_fOV0.020/inlist` but is set to 0.01 in `m1.000_fOV0.030/inlist` as expected from `grid.set_inlist_finalize_function`." ] }, { "cell_type": "code", "execution_count": 6, "id": "9823610e6e6c6ef9", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:48:47.762896358Z", "start_time": "2026-03-19T15:48:47.686897137Z" }, "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "{'star_job': {'history_columns_file': 'history_test.list',\n", " 'pgstar_flag': False,\n", " 'create_pre_main_sequence_model': False},\n", " 'eos': {},\n", " 'kap': {'use_type2_opacities': True, 'zbase': 0.02},\n", " 'controls': {'initial_mass': 1,\n", " 'initial_z': 0.02,\n", " 'max_model_number': 10,\n", " 'mesh_delta_coeff': 5,\n", " 'history_interval': 1,\n", " 'photo_interval': 8,\n", " 'write_profiles_flag': True,\n", " 'profile_interval': 9,\n", " 'write_pulse_data_with_profile': True,\n", " 'pulse_data_format': 'GYRE',\n", " 'overshoot_scheme(1)': 'exponential',\n", " 'overshoot_zone_type(1)': 'nonburn',\n", " 'overshoot_zone_loc(1)': 'shell',\n", " 'overshoot_bdy_loc(1)': 'bottom',\n", " 'overshoot_f(1)': 0.02,\n", " 'overshoot_f0(1)': 0.001},\n", " 'pgstar': {}}" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from wsssss.inlists import inlists as inl\n", "inl.evaluate_inlist(f'{grid_dir}/m1.000_fOV0.020/inlist')" ] }, { "cell_type": "code", "execution_count": 7, "id": "e3314a144c43f7bf", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:06.028742297Z", "start_time": "2026-03-19T15:49:06.003724252Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "left : /home/walter/Github/wsssss/examples/example_wsssss/m1.000_fOV0.020/inlist\n", "right: /home/walter/Github/wsssss/examples/example_wsssss/m1.000_fOV0.030/inlist\n", "\n", "####### star_job differences #######\n", "\n", "####### controls differences #######\n", "overshoot_f(1): (0.02, 0.03)\n", "overshoot_f0(1): (0.001, 0.01)\n", "overshoot_d_min: (None, 0.01)\n" ] } ], "source": [ "inl.compare_inlist(f'{grid_dir}/m1.000_fOV0.020/inlist', f'{grid_dir}/m1.000_fOV0.030/inlist')" ] }, { "cell_type": "markdown", "id": "4fd76915-4c46-460a-9473-02bb69a92f01", "metadata": {}, "source": [ "However, comparing MESA version strings should be done carefully:" ] }, { "cell_type": "code", "execution_count": 8, "id": "54957a4d760ab2f8", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:12.110476371Z", "start_time": "2026-03-19T15:49:12.098518157Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "True True\n", "False True\n" ] } ], "source": [ "print(inl.compare_version('11701', 'r24.03.1', '<'), '11701' < 'r24.03.1')\n", "print(inl.compare_version('11701', '8118', '<'), '11701' < '8118')" ] }, { "cell_type": "code", "execution_count": 9, "id": "40e9f3c9b098f771", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:13.118102228Z", "start_time": "2026-03-19T15:49:13.106149810Z" }, "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "'r24.08.1'" ] }, "execution_count": 9, "metadata": {}, "output_type": "execute_result" } ], "source": [ "inl.get_mesa_version(mesa_dir)" ] }, { "cell_type": "markdown", "id": "4d4f1861-6b75-46ea-86b2-aeec7c38730c", "metadata": { "ExecuteTime": { "start_time": "2026-03-19T15:45:19.016508514Z" }, "collapsed": false }, "source": [ "# MESA histories and profiles\n", "\n", "`wsssss` can also read MESA histories and profiles, as well as gyre summaries and mode files." ] }, { "cell_type": "code", "execution_count": 10, "id": "d104c7d9e26d02ac", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:22.585335761Z", "start_time": "2026-03-19T15:49:22.537389212Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "['model_number', 'elapsed_time', 'star_age', 'star_mass', 'effective_T', 'photosphere_L', 'photosphere_r', 'luminosity', 'radius', 'gravity', 'center_T', 'center_Rho', 'center_P', 'center_degeneracy', 'center_mu', 'center_ye', 'center_h1', 'center_he3', 'center_he4', 'center_c12', 'center_n14', 'center_o16', 'center_ne20', 'center_mg24', 'surface_h1', 'surface_he3', 'surface_he4', 'surface_c12', 'surface_n14', 'surface_o16', 'surface_ne20', 'surface_mg24', 'mass_conv_core', 'he_core_mass', 'he_core_radius', 'co_core_mass', 'co_core_radius', 'log_LH', 'log_LHe', 'log_LZ', 'log_extra_L', 'log_Lneu', 'log_abs_Lgrav', 'pp', 'cno', 'tri_alpha', 'c_alpha', 'n_alpha', 'o_alpha', 'delta_nu', 'delta_Pg', 'nu_max', 'log_L', 'mix_type_1', 'mix_qtop_1', 'mix_type_2', 'mix_qtop_2', 'mix_type_3', 'mix_qtop_3', 'mix_type_4', 'mix_qtop_4', 'mix_type_5', 'mix_qtop_5', 'mix_type_6', 'mix_qtop_6', 'mix_type_7', 'mix_qtop_7', 'mix_type_8', 'mix_qtop_8', 'mix_type_9', 'mix_qtop_9', 'mix_type_10', 'mix_qtop_10', 'mix_type_11', 'mix_qtop_11', 'mix_type_12', 'mix_qtop_12', 'mix_type_13', 'mix_qtop_13', 'mix_type_14', 'mix_qtop_14', 'mix_type_15', 'mix_qtop_15', 'mix_type_16', 'mix_qtop_16', 'mix_type_17', 'mix_qtop_17', 'mix_type_18', 'mix_qtop_18', 'mix_type_19', 'mix_qtop_19', 'mix_type_20', 'mix_qtop_20', 'burn_type_1', 'burn_qtop_1', 'burn_type_2', 'burn_qtop_2', 'burn_type_3', 'burn_qtop_3', 'burn_type_4', 'burn_qtop_4', 'burn_type_5', 'burn_qtop_5', 'burn_type_6', 'burn_qtop_6', 'burn_type_7', 'burn_qtop_7', 'burn_type_8', 'burn_qtop_8', 'burn_type_9', 'burn_qtop_9', 'burn_type_10', 'burn_qtop_10', 'burn_type_11', 'burn_qtop_11', 'burn_type_12', 'burn_qtop_12', 'burn_type_13', 'burn_qtop_13', 'burn_type_14', 'burn_qtop_14', 'burn_type_15', 'burn_qtop_15', 'burn_type_16', 'burn_qtop_16', 'burn_type_17', 'burn_qtop_17', 'burn_type_18', 'burn_qtop_18', 'burn_type_19', 'burn_qtop_19', 'burn_type_20', 'burn_qtop_20', 'burn_type_21', 'burn_qtop_21', 'burn_type_22', 'burn_qtop_22', 'burn_type_23', 'burn_qtop_23', 'burn_type_24', 'burn_qtop_24', 'burn_type_25', 'burn_qtop_25', 'burn_type_26', 'burn_qtop_26', 'burn_type_27', 'burn_qtop_27', 'burn_type_28', 'burn_qtop_28', 'burn_type_29', 'burn_qtop_29', 'burn_type_30', 'burn_qtop_30']\n" ] } ], "source": [ "from wsssss import load_data as ld\n", "import numpy as np\n", "\n", "hist = ld.History(f'{grid_dir}/m1.000_fOV0.020/LOGS/history.data')\n", "print(hist.columns)" ] }, { "cell_type": "markdown", "id": "645048b2-21a3-4a85-9d66-bdfc21bde49a", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:23.205345940Z", "start_time": "2026-03-19T15:49:23.185721425Z" }, "collapsed": false }, "source": [ "If we don't need all these columns we can choose which columns to keep in `History` (or and data class in `load_data`) :" ] }, { "cell_type": "code", "execution_count": 11, "id": "ffafa34943dd7548", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:23.795247653Z", "start_time": "2026-03-19T15:49:23.791528949Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "['model_number', 'star_age', 'effective_T', 'luminosity', 'nu_max', 'delta_nu', 'delta_Pg']\n" ] } ], "source": [ "hist = ld.History(f'{grid_dir}/m1.000_fOV0.020/LOGS/history.data', keep_columns=['model_number', 'star_age', 'effective_T', 'luminosity', 'nu_max', 'delta_nu', 'delta_Pg'])\n", "print(hist.columns)" ] }, { "cell_type": "markdown", "id": "2e9607c0-6525-42f0-8133-a09f36204af7", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:24.855423495Z", "start_time": "2026-03-19T15:49:24.844769330Z" }, "collapsed": false }, "source": [ "We can access the data in two ways, through attributes of `hist.data` or through the `hist.get` method:" ] }, { "cell_type": "code", "execution_count": 12, "id": "d73c9617ebb10b4a", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:25.287931422Z", "start_time": "2026-03-19T15:49:25.245656617Z" }, "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "True" ] }, "execution_count": 12, "metadata": {}, "output_type": "execute_result" } ], "source": [ "hist.data.star_age\n", "hist.get('star_age')\n", "np.all(hist.get('star_age') == hist.data.star_age)" ] }, { "cell_type": "markdown", "id": "d08015ea-1422-4736-ae6d-6e4dc257297c", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:26.531892130Z", "start_time": "2026-03-19T15:49:26.504754035Z" }, "collapsed": false }, "source": [ "We can also index directly on `hist`, this will return a new `History` with a correctly set profile index, containing only profiles in the new `History` file." ] }, { "cell_type": "code", "execution_count": 13, "id": "d584b3c09ac48ef9", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:26.532289725Z", "start_time": "2026-03-19T15:49:26.510682776Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "10\n", "[( 1, 100000. , 5616.90892797, 0.70925563, 168.08043363, 0., 3945.14980227)\n", " ( 2, 220000. , 5614.14501653, 0.70804627, 168.05210824, 0., 3945.08631765)\n", " ( 3, 364000. , 5604.89019466, 0.70420514, 167.92708726, 0., 3943.76604529)\n", " ( 4, 536800. , 5580.62037896, 0.69460378, 167.59231705, 0., 3938.01161344)\n", " ( 5, 744160. , 5575.69560644, 0.6946222 , 167.152432 , 0., 3925.75774334)\n", " ( 6, 992992. , 5577.13854459, 0.69672805, 166.90750047, 0., 3917.43838629)\n", " ( 7, 1291590.4 , 5578.83351588, 0.69861458, 166.72173235, 0., 3911.01712141)\n", " ( 8, 1649908.48 , 5579.89506039, 0.70000223, 166.57075811, 0., 3905.86420042)\n", " ( 9, 2079890.176 , 5580.69086712, 0.7011625 , 166.43786954, 0., 3901.34770013)\n", " (10, 2595868.2112, 5581.41662209, 0.7022522 , 166.31199997, 0., 3897.06714391)]\n", "[[ 1 2 1]\n", " [ 9 1 2]\n", " [10 3 3]]\n" ] } ], "source": [ "print(len(hist))\n", "print(hist.data)\n", "print(hist.index)" ] }, { "cell_type": "code", "execution_count": 14, "id": "b9cfb21cbfcf658", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:26.532536573Z", "start_time": "2026-03-19T15:49:26.517056585Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "5\n", "[(1, 100000., 5616.90892797, 0.70925563, 168.08043363, 0., 3945.14980227)\n", " (2, 220000., 5614.14501653, 0.70804627, 168.05210824, 0., 3945.08631765)\n", " (3, 364000., 5604.89019466, 0.70420514, 167.92708726, 0., 3943.76604529)\n", " (4, 536800., 5580.62037896, 0.69460378, 167.59231705, 0., 3938.01161344)\n", " (5, 744160., 5575.69560644, 0.6946222 , 167.152432 , 0., 3925.75774334)]\n", "[[1 2 1]]\n" ] } ], "source": [ "new_hist = hist[:5]\n", "print(len(new_hist))\n", "print(new_hist.data)\n", "print(new_hist.index)" ] }, { "cell_type": "markdown", "id": "a6bb3a66-9363-4c1d-a0b8-7684ba2235da", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:26.596025349Z", "start_time": "2026-03-19T15:49:26.520719361Z" }, "collapsed": false }, "source": [ "The same is true for `gyre` summaries:" ] }, { "cell_type": "code", "execution_count": 15, "id": "a5e13055e00e18bf", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:26.811554667Z", "start_time": "2026-03-19T15:49:26.669945070Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "[(1., 1.82764478e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3309.71129536, 0., 0, 0, 19, 19)\n", " (1., 1.57116548e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3471.5380448 , 0., 0, 0, 20, 20)\n", " (1., 1.38862698e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3634.47656824, 0., 0, 0, 21, 21)\n", " (1., 1.24263570e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3797.7041817 , 0., 0, 0, 22, 22)\n", " (1., 1.12116935e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3961.31386135, 0., 0, 0, 23, 23)\n", " (1., 1.03091459e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4124.99185735, 0., 0, 0, 24, 24)\n", " (1., 9.63169050e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4288.82138142, 0., 0, 0, 25, 25)\n", " (1., 9.15110781e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4453.22548755, 0., 0, 0, 26, 26)\n", " (1., 8.77471803e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4617.86112542, 0., 0, 0, 27, 27)\n", " (1., 8.44784618e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4782.93365092, 0., 0, 0, 28, 28)\n", " (1., 2.01195016e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3223.42320367, 0., 1, 0, 17, 18)\n", " (1., 1.69127179e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3385.16288116, 0., 1, 0, 18, 19)\n", " (1., 1.48022252e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3547.55896439, 0., 1, 0, 19, 20)\n", " (1., 1.31330986e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3710.78091331, 0., 1, 0, 20, 21)\n", " (1., 1.18119987e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3874.44966317, 0., 1, 0, 21, 22)\n", " (1., 1.07390212e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4038.06751587, 0., 1, 0, 22, 23)\n", " (1., 9.93904450e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4202.03953297, 0., 1, 0, 23, 24)\n", " (1., 9.38204390e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4366.1774459 , 0., 1, 0, 24, 25)\n", " (1., 8.94148813e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4530.80566806, 0., 1, 0, 25, 26)\n", " (1., 8.60259786e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4695.84074713, 0., 1, 0, 26, 27)]\n", "\n", "[(1., 1.82764478e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3309.71129536, 0., 0, 0, 19, 19)\n", " (1., 1.57116548e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3471.5380448 , 0., 0, 0, 20, 20)\n", " (1., 1.38862698e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3634.47656824, 0., 0, 0, 21, 21)\n", " (1., 1.24263570e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3797.7041817 , 0., 0, 0, 22, 22)\n", " (1., 1.12116935e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 3961.31386135, 0., 0, 0, 23, 23)\n", " (1., 1.03091459e-08, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4124.99185735, 0., 0, 0, 24, 24)\n", " (1., 9.63169050e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4288.82138142, 0., 0, 0, 25, 25)\n", " (1., 9.15110781e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4453.22548755, 0., 0, 0, 26, 26)\n", " (1., 8.77471803e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4617.86112542, 0., 0, 0, 27, 27)\n", " (1., 8.44784618e-09, 2.71503057e+33, 1.98840987e+33, 6.18702395e+10, 4782.93365092, 0., 0, 0, 28, 28)]\n" ] } ], "source": [ "gs = ld.GyreSummary(f'{hist.LOGS}/../gyre_out/profile1.data.GYRE.sgyre_l')\n", "print(gs.data)\n", "print()\n", "print(gs[gs.data.l == 0].data)" ] }, { "cell_type": "markdown", "id": "d9d7773d-2b3d-4709-84e4-3a38c364d39a", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:26.957571817Z", "start_time": "2026-03-19T15:49:26.825009723Z" }, "collapsed": false }, "source": [ "Convenience function to load all MESA profiles or `gyre` summaries are also available:" ] }, { "cell_type": "code", "execution_count": 16, "id": "b89b158090faa72b", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:27.045068889Z", "start_time": "2026-03-19T15:49:26.957052265Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "[MESA profile data file at /home/walter/Github/wsssss/examples/example_wsssss/m1.000_fOV0.020/LOGS/profile1.data\n", "{'model_number': 1, 'num_zones': 3828, 'star_mass': 1.0, 'star_age': 100000.0, 'Teff': 5616.908927965037, 'photosphere_L': 0.7092556343734885, 'center_h1': 0.6969903923887202, 'center_he4': 0.2825119329953253, 'date': '20260320'}, MESA profile data file at /home/walter/Github/wsssss/examples/example_wsssss/m1.000_fOV0.020/LOGS/profile2.data\n", "{'model_number': 9, 'num_zones': 217, 'star_mass': 1.0, 'star_age': 2079890.176, 'Teff': 5580.690867117586, 'photosphere_L': 0.7011624968709207, 'center_h1': 0.6968676556622809, 'center_he4': 0.2826257295568369, 'date': '20260320'}, MESA profile data file at /home/walter/Github/wsssss/examples/example_wsssss/m1.000_fOV0.020/LOGS/profile3.data\n", "{'model_number': 10, 'num_zones': 217, 'star_mass': 1.0, 'star_age': 2595868.2112000003, 'Teff': 5581.416622090657, 'photosphere_L': 0.70225220496579, 'center_h1': 0.696835018765734, 'center_he4': 0.28265520404943156, 'date': '20260320'}]\n" ] }, { "data": { "text/plain": [ "[GyreSummary at /home/walter/Github/wsssss/examples/example_wsssss/m1.000_fOV0.020/gyre_out/profile1.data.GYRE.sgyre_l,\n", " GyreSummary at /home/walter/Github/wsssss/examples/example_wsssss/m1.000_fOV0.020/gyre_out/profile2.data.GYRE.sgyre_l,\n", " GyreSummary at /home/walter/Github/wsssss/examples/example_wsssss/m1.000_fOV0.020/gyre_out/profile3.data.GYRE.sgyre_l]" ] }, "execution_count": 16, "metadata": {}, "output_type": "execute_result" } ], "source": [ "profs = ld.load_profs(hist)\n", "print(profs)\n", "ld.load_gss(hist)" ] }, { "cell_type": "markdown", "id": "b0599a0b-a107-4000-93f6-79803729dadb", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:27.254336Z", "start_time": "2026-03-19T15:49:27.131643892Z" }, "collapsed": false }, "source": [ "We can also load a `gyre` summary from a MESA profile:" ] }, { "cell_type": "code", "execution_count": 17, "id": "43ae2bdccf47faba", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:27.305939203Z", "start_time": "2026-03-19T15:49:27.270015467Z" }, "collapsed": false }, "outputs": [], "source": [ "gs = ld.load_gs_from_profile(profs[0])" ] }, { "cell_type": "markdown", "id": "b69c5313-b38b-4c50-88eb-679b02e10ad2", "metadata": {}, "source": [ "Various convenience functions are available in `functions`, for example calculating $\\Delta\\nu$ from a gyre summary. The full list is available in the [documentation.](https://wsssss.readthedocs.io/latest/wsssss.html#module-wsssss.functions)" ] }, { "cell_type": "code", "execution_count": 18, "id": "32402f79-d83b-46b5-ba6a-85181c828e84", "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "168.08043363435422 163.70286603873035\n" ] } ], "source": [ "from wsssss import functions as uf\n", "print(hist.data.delta_nu[hist.get_profile_index(1)[0]], uf.calc_deltanu(gs, hist))" ] }, { "cell_type": "markdown", "id": "05211e85-0a0a-49d1-8020-1593c24aff07", "metadata": {}, "source": [ "# Constants\n", "\n", "Over the course of MESA's development the physical constant definitions have changed. `wsssss` knows about this and can get the constants used by the version of MESA used by a `History`. The constants in `wsssss` are generated from `$MESA_DIR/const/public/const_def.f90` for pre- and post-15140 as the definitions changed in 15140." ] }, { "cell_type": "code", "execution_count": 19, "id": "c62b7b07cd0ab75e", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:27.635333437Z", "start_time": "2026-03-19T15:49:27.621309923Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "6.6743e-08\n" ] } ], "source": [ "c = uf.get_constants(hist)\n", "print(c.standard_cgrav)" ] }, { "cell_type": "code", "execution_count": 20, "id": "8df850e999e09e58", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:27.803128819Z", "start_time": "2026-03-19T15:49:27.794100918Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "6.6743e-08 6.67428e-08\n", "False\n" ] } ], "source": [ "from wsssss.constants import pre15140\n", "print(c.standard_cgrav, pre15140.standard_cgrav)\n", "print(c.standard_cgrav == pre15140.standard_cgrav)" ] }, { "cell_type": "code", "execution_count": 21, "id": "6a0c1ad3dcace050", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:28.141694826Z", "start_time": "2026-03-19T15:49:28.137360877Z" }, "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "['version', 'max_extra_inlists', 'pi', 'pi2', 'pi4', 'eulercon', 'eulernum', 'ln2', 'ln3', 'lnpi', 'ln10', 'iln10', 'a2rad', 'rad2a', 'one_third', 'two_thirds', 'four_thirds', 'five_thirds', 'one_sixth', 'four_thirds_pi', 'ln4pi3', 'two_13', 'four_13', 'sqrt2', 'sqrt_2_div_3', 'avo', 'amu', 'clight', 'qe', 'kerg', 'boltzm', 'planck_h', 'hbar', 'cgas', 'ev2erg', 'mev_to_ergs', 'mev_amu', 'mev2gr', 'qconv', 'kev', 'boltz_sigma', 'crad', 'au', 'pc', 'dayyer', 'secday', 'secyer', 'ly', 'mn', 'mp', 'me', 'rbohr', 'fine', 'hion', 'sige', 'weinberg_theta', 'num_neu_fam', 'standard_cgrav', 'mu_sun', 'mu_earth', 'mu_jupiter', 'agesun', 'msun', 'rsun', 'lsun', 'teffsun', 'loggsun', 'mbolsun', 'm_earth', 'r_earth', 'r_earth_polar', 'm_jupiter', 'r_jupiter', 'r_jupiter_polar', 'semimajor_axis_jupiter', 'arg_not_provided', 'missing_value', 'crystallized', 'no_mixing', 'convective_mixing', 'overshoot_mixing', 'semiconvective_mixing', 'thermohaline_mixing', 'rotation_mixing', 'rayleigh_taylor_mixing', 'minimum_mixing', 'anonymous_mixing', 'leftover_convective_mixing', 'phase_separation_mixing', 'number_of_mixing_types']\n" ] } ], "source": [ "print(list(c.__dict__.keys())[8:])" ] }, { "cell_type": "markdown", "id": "a6d79951-e147-4456-ad7c-74432ed5e100", "metadata": { "ExecuteTime": { "end_time": "2026-03-19T15:49:28.366895737Z", "start_time": "2026-03-19T15:49:28.298920047Z" }, "collapsed": false }, "source": [ "# Plotting\n", "There are many plotting functions avialable in `wsssss.plotting.plotting` ([docs](https://wsssss.readthedocs.io/latest/wsssss.plotting.html)), some of them are exhibited here:" ] }, { "cell_type": "code", "execution_count": 22, "id": "9d445329772dd1e3", "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "(
,\n", " )" ] }, "execution_count": 22, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAjcAAAG2CAYAAACDLKdOAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjguMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/H5lhTAAAACXBIWXMAAA9hAAAPYQGoP6dpAAAtG0lEQVR4nO3dfVxUdd7/8fcwwICieJeDN5iWea9o3hBamxXpdkM/97q2zLrU2NpuVk0kLcyU7Ea0XV0frZjpttleV6bWZrnVWkZqaRQmUlreVFqiCWoqo5AgM+f3BzUsKyqD6IEvr+fjcR6P5jvfc87nw5mad+ecmXFYlmUJAADAEEF2FwAAAFCTCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCi2hpsPP/xQCQkJat26tRwOh954442zrrN27Vpdfvnlcrlc6tixoxYvXnze6wQAAHWHreGmsLBQMTExSk9Pr9L83bt366abbtI111yjnJwcJSUl6Z577tG77757nisFAAB1haO2/HCmw+HQihUrNGzYsNPOeeSRR/T2229r69at/rHbb79dR48e1apVqy5AlQAAoLYLtruAQGRmZio+Pr7C2NChQ5WUlHTadYqLi1VcXOx/7PP5dPjwYTVv3lwOh+N8lQoAAGqQZVk6duyYWrduraCgM194qlPhJi8vT263u8KY2+2Wx+PRTz/9pPDw8FPWSUtL0/Tp0y9UiQAA4DzKzc1V27ZtzzinToWb6pg8ebKSk5P9jwsKCtSuXTu9tetPio6IkSQdcmzTByEpamK117Un0xSsBgHv58ugl/VFyGL1OnmXuvvuDHj9UhXpg5DJOur4TteenKkWVteAt0Ef5eijDH2Uo49y9FGGPsrVhT5yj3+umy+ZqEaNGp11W3Uq3ERFRSk/P7/CWH5+vho3blzpWRtJcrlccrlcp4xHR8ToskZX6gfHRr0V+pgutnrp1pKVcoWd/Y/2nz52ztSOkMWKP5mqgaEpAa9frGN6NfQWlTj2KrFklVq7+ge8DfooRx9l6KMcfZSjjzL0Ua6u9VGVW0rq1PfcxMXFKSMjo8LY6tWrFRcXV63t/eDYqGWhN+siq1vZAVX1DuhHIdN11clUDfRW/4AedHyl4SVvqbVVvRcmfZShjzL0UY4+ytFHGfooZ0of/8nWcHP8+HHl5OQoJydHUtlHvXNycrRnzx5JZZeURo0a5Z9///33a9euXXr44Ye1fft2zZ8/X8uXL9eECRMC3vchxzYjDqgpL0z6KEcfZeijHH2Uo48y9HFmtoabzz77TH369FGfPn0kScnJyerTp4+mTZsmSdq/f78/6EhShw4d9Pbbb2v16tWKiYnR7Nmz9de//lVDhw4NeN8fhKTU+QNqyguTPsrRRxn6KEcf5eijDH2cXa35npsLxePxKDIyUn880k33hK2tswfUlBcmfZSjjzL0UY4+ytFHmfrcx9fH1ium5fUqKChQ48aNzzi3Tt1zU5OuPZlWZw7of6qrL8z/RB/l6KMMfZSjj3L0UYY+qq7ehpvqfMytNhxQU16Y9FGOPsrQRzn6KEcfZegjMPX2stTnB1brskZXVnm92nBATXlh0ke52tbHgF3TdbKoJOBt7LxoqXa2elmd9t+pTgdvD3j9UkeRPunwuI6Hfa/Y3U+o6U+dA97GkfAd+rTDNEWcuFhX7H5cwVbg/wNDH2Xooxx9lLO7jx+8W/X/Lp9cpctS9Tbc/GNzqtoF9avSOnYfUMmMF6ZEH/+utvVx2dZEbbzi4YC3AQAXQrFHmhkpwk1lfgk3KQWS68x/G6Be6vRiqppHdKjS3L0b9ir3wz2K/lU7tR105q9DPxtviVfbln2looNF6nZ7N0W0rtqZsOM/HNNXS79Sg4saqOvwbgprEqbG0dX/l/vrlV9rx+vb1fm/uuiyWy6r0jqlJ0r16exPdGzvMV0x8Qo1ubRptfcvSUe/PaJP/vSJGrVtpNiHrlBwWNW+b7U6tZ9OdXuqbu2nw/Eow/GQfvhqv/7fiHsJN5Xxh5ted8nlDLW7HEhq2LKhuo/sIWeo84zzavKN9BeBvqH+5xvp2WquqkB6q24IOJtDXx7U1y/vkeOb5jWyPTv95tX/UqM21f+7fLH4C+Us3Kze9/ZRr7t6nXHuyaKTen/Cah3ZdURD5g5Ri+4XVXu//+7Qlwf1XtJ7anpJU8X/+XqFNAipsZqrKtDeAq25qjgeZer78fh+8/e6Zsxgwk1l/PfcrP5CjRrWzJsCqufYvmNacevrdpeB/3CuwcBOvKYAc53QCc3UTMJNZX4JN3k//yYV7HX4m8MqORb4Daw4P0IbhapZx2Z2l3FOeE0BZjpWeEwx1/eqUripUz+cCfPU9TdS1D68pgAzNfBU/YMb9fZ7bgAAgJkINwAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjEK4AQAARiHcAAAAoxBuAACAUQg3AADAKIQbAABgFMINAAAwCuEGAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjEK4AQAARiHcAAAAoxBuAACAUQg3AADAKIQbAABgFMINAAAwCuEGAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjGJ7uElPT1f79u0VFham2NhYZWVlnXH+3Llz1blzZ4WHhys6OloTJkzQiRMnLlC1AACgtrM13CxbtkzJyclKTU1Vdna2YmJiNHToUB04cKDS+UuWLFFKSopSU1O1bds2vfDCC1q2bJkeffTRC1w5AACorWwNN3PmzNHvf/97JSYmqlu3blqwYIEaNGigv/3tb5XO//jjjzVo0CDdcccdat++vYYMGaIRI0ac9WwPAACoP2wLNyUlJdq0aZPi4+PLiwkKUnx8vDIzMytdZ+DAgdq0aZM/zOzatUvvvPOObrzxxtPup7i4WB6Pp8ICAADMFWzXjg8dOiSv1yu3211h3O12a/v27ZWuc8cdd+jQoUO68sorZVmWSktLdf/995/xslRaWpqmT59eo7UDAIDay/YbigOxdu1azZgxQ/Pnz1d2drZef/11vf3223ryySdPu87kyZNVUFDgX3Jzcy9gxQAA4EKz7cxNixYt5HQ6lZ+fX2E8Pz9fUVFRla4zdepUjRw5Uvfcc48kqWfPniosLNS9996rKVOmKCjo1KzmcrnkcrlqvgEAAFAr2XbmJjQ0VH379lVGRoZ/zOfzKSMjQ3FxcZWuU1RUdEqAcTqdkiTLss5fsQAAoM6w7cyNJCUnJ2v06NHq16+fBgwYoLlz56qwsFCJiYmSpFGjRqlNmzZKS0uTJCUkJGjOnDnq06ePYmNj9c0332jq1KlKSEjwhxwAAFC/2Rpuhg8froMHD2ratGnKy8tT7969tWrVKv9Nxnv27Klwpuaxxx6Tw+HQY489pn379umiiy5SQkKCnn76abtaAAAAtYzDqmfXczwejyIjI5WXn6/GjRvbXQ4AAKgCj8ejKLdbBQUFZ33/rlOflgIAADgbwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjEK4AQAARiHcAAAAoxBuAACAUQg3AADAKIQbAABgFMINAAAwCuEGAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjEK4AQAARiHcAAAAoxBuAACAUQg3AADAKIQbAABgFMINAAAwCuEGAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjEK4AQAARiHcAAAAoxBuAACAUWwPN+np6Wrfvr3CwsIUGxurrKysM84/evSoxowZo1atWsnlcqlTp0565513LlC1AACgtgu2c+fLli1TcnKyFixYoNjYWM2dO1dDhw7Vjh071LJly1Pml5SU6Prrr1fLli312muvqU2bNvr+++/VpEmTC188AAColRyWZVl27Tw2Nlb9+/fXvHnzJEk+n0/R0dEaN26cUlJSTpm/YMEC/fGPf9T27dsVEhJSrX16PB5FRkYqLz9fjRs3Pqf6AQDAheHxeBTldqugoOCs79+2XZYqKSnRpk2bFB8fX15MUJDi4+OVmZlZ6TorV65UXFycxowZI7fbrR49emjGjBnyer2n3U9xcbE8Hk+FBQAAmMu2cHPo0CF5vV653e4K4263W3l5eZWus2vXLr322mvyer165513NHXqVM2ePVtPPfXUafeTlpamyMhI/xIdHV2jfQAAgNrF9huKA+Hz+dSyZUstXLhQffv21fDhwzVlyhQtWLDgtOtMnjxZBQUF/iU3N/cCVgwAAC40224obtGihZxOp/Lz8yuM5+fnKyoqqtJ1WrVqpZCQEDmdTv9Y165dlZeXp5KSEoWGhp6yjsvlksvlqtniAQBArWXbmZvQ0FD17dtXGRkZ/jGfz6eMjAzFxcVVus6gQYP0zTffyOfz+cd27typVq1aVRpsAABA/WPrZank5GQtWrRIL730krZt26YHHnhAhYWFSkxMlCSNGjVKkydP9s9/4IEHdPjwYY0fP147d+7U22+/rRkzZmjMmDF2tQAAAGoZW7/nZvjw4Tp48KCmTZumvLw89e7dW6tWrfLfZLxnzx4FBZXnr+joaL377ruaMGGCevXqpTZt2mj8+PF65JFH7GoBAADUMrZ+z40d+J4bAADqnjrxPTcAAADnA+EGAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCjVCjfZ2dnasmWL//Gbb76pYcOG6dFHH1VJSUmNFQcAABCoaoWb++67Tzt37pQk7dq1S7fffrsaNGigV199VQ8//HCNFggAABCIaoWbnTt3qnfv3pKkV199Vb/61a+0ZMkSLV68WP/4xz9qsj4AAICAVCvcWJYln88nSXr//fd14403SpKio6N16NChmqsOAAAgQNUKN/369dNTTz2l//3f/9W6det00003SZJ2794tt9tdowUCAAAEolrhZu7cucrOztbYsWM1ZcoUdezYUZL02muvaeDAgTVaIAAAQCAclmVZNbWxEydOyOl0KiQkpKY2WeM8Ho8iIyOVl5+vxo0b210OAACoAo/Hoyi3WwUFBWd9/67WmZvc3Fzt3bvX/zgrK0tJSUn6+9//XquDDQAAMF+1ws0dd9yhNWvWSJLy8vJ0/fXXKysrS1OmTNETTzxRowUCAAAEolrhZuvWrRowYIAkafny5erRo4c+/vhjvfzyy1q8eHFN1gcAABCQaoWbkydPyuVySSr7KPgtt9wiSerSpYv2799fc9UBAAAEqFrhpnv37lqwYIE++ugjrV69Wr/+9a8lST/88IOaN29eowUCAAAEolrhZtasWXr++ec1ePBgjRgxQjExMZKklStX+i9XAQAA2CG4OisNHjxYhw4dksfjUdOmTf3j9957rxo0aFBjxQEAAASqWuFGkpxOZ4VgI0nt27c/13oAAADOSbXDzWuvvably5drz549KikpqfBcdnb2ORcGAABQHdW65+bZZ59VYmKi3G63Nm/erAEDBqh58+batWuXbrjhhpquEQAAoMqqFW7mz5+vhQsX6i9/+YtCQ0P18MMPa/Xq1XrwwQdVUFBQ0zUCAABUWbXCzZ49e/w/kBkeHq5jx45JkkaOHKlXXnml5qoDAAAIULXCTVRUlA4fPixJateunT755BNJ0u7du1WDv8MJAAAQsGqFm2uvvVYrV66UJCUmJmrChAm6/vrrNXz4cP3mN7+p0QIBAAAC4bCqcarF5/PJ5/MpOLjsw1bLli3Thg0bdNlll+n++++v1b8M7vF4FBkZqbz8/LP+ZDoAAKgdPB6PotxuFRQUnPX9u1rhRpJOnDihL774QgcOHJDP5yvfoMOhhISE6mzygiDcAABQ9wQSbqr1PTerVq3SyJEj9eOPP57ynMPhkNfrrc5mAQAAzlm17rkZN26cbrvtNu3fv99/ieqXhWADAADsVK1wk5+fr+TkZLnd7pquBwAA4JxUK9z89re/1dq1a2u4FAAAgHNXrRuKi4qKdOutt+qiiy5Sz549T/l01IMPPlhjBdY0bigGAKDuOe83FL/yyit67733FBYWprVr18rhcPifczgctTrcAAAAs1Ur3EyZMkXTp09XSkqKgoKqdWULAADgvKhWMikpKdHw4cMJNgAAoNapVjoZPXq0li1bVtO1AAAAnLNqXZbyer165pln9O6776pXr16n3FA8Z86cGikOAAAgUNUKN1u2bFGfPn0kSVu3bq3w3L/fXAwAAHChVSvcrFmzpqbrAAAAqBHcEQwAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjEK4AQAARiHcAAAAoxBuAACAUQg3AADAKLUi3KSnp6t9+/YKCwtTbGyssrKyqrTe0qVL5XA4NGzYsPNbIAAAqDNsDzfLli1TcnKyUlNTlZ2drZiYGA0dOlQHDhw443rfffedJk6cqKuuuuoCVQoAAOoC28PNnDlz9Pvf/16JiYnq1q2bFixYoAYNGuhvf/vbadfxer268847NX36dF1yySUXsFoAAFDb2RpuSkpKtGnTJsXHx/vHgoKCFB8fr8zMzNOu98QTT6hly5a6++67z7qP4uJieTyeCgsAADCXreHm0KFD8nq9crvdFcbdbrfy8vIqXWf9+vV64YUXtGjRoirtIy0tTZGRkf4lOjr6nOsGAAC1l+2XpQJx7NgxjRw5UosWLVKLFi2qtM7kyZNVUFDgX3Jzc89zlQAAwE7Bdu68RYsWcjqdys/PrzCen5+vqKioU+Z/++23+u6775SQkOAf8/l8kqTg4GDt2LFDl156aYV1XC6XXC7XeageAADURraeuQkNDVXfvn2VkZHhH/P5fMrIyFBcXNwp87t06aItW7YoJyfHv9xyyy265pprlJOTwyUnAABg75kbSUpOTtbo0aPVr18/DRgwQHPnzlVhYaESExMlSaNGjVKbNm2UlpamsLAw9ejRo8L6TZo0kaRTxgEAQP1ke7gZPny4Dh48qGnTpikvL0+9e/fWqlWr/DcZ79mzR0FBderWIAAAYCOHZVmW3UVcSB6PR5GRkcrLz1fjxo3tLgcAAFSBx+NRlNutgoKCs75/c0oEAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjEK4AQAARiHcAAAAoxBuAACAUQg3AADAKIQbAABgFMINAAAwCuEGAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYJRguwuwy49FJSoJLrG7DAAAUAXHiqr+nl1vw80/d/6osIaEGwAA6oIThceqPJfLUgAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjEK4AQAARiHcAAAAoxBuAACAUQg3AADAKIQbAABgFMINAAAwCuEGAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEapFeEmPT1d7du3V1hYmGJjY5WVlXXauYsWLdJVV12lpk2bqmnTpoqPjz/jfAAAUL/YHm6WLVum5ORkpaamKjs7WzExMRo6dKgOHDhQ6fy1a9dqxIgRWrNmjTIzMxUdHa0hQ4Zo3759F7hyAABQGzksy7LsLCA2Nlb9+/fXvHnzJEk+n0/R0dEaN26cUlJSzrq+1+tV06ZNNW/ePI0aNeqs8z0ejyIjIzXz3RyFNWx0zvUDAIDz70ThMaUM7a2CggI1btz4jHNtPXNTUlKiTZs2KT4+3j8WFBSk+Ph4ZWZmVmkbRUVFOnnypJo1a1bp88XFxfJ4PBUWAABgLlvDzaFDh+T1euV2uyuMu91u5eXlVWkbjzzyiFq3bl0hIP27tLQ0RUZG+pfo6OhzrhsAANRett9zcy5mzpyppUuXasWKFQoLC6t0zuTJk1VQUOBfcnNzL3CVAADgQgq2c+ctWrSQ0+lUfn5+hfH8/HxFRUWdcd0//elPmjlzpt5//3316tXrtPNcLpdcLleN1AsAAGo/W8/chIaGqm/fvsrIyPCP+Xw+ZWRkKC4u7rTrPfPMM3ryySe1atUq9evX70KUCgAA6ghbz9xIUnJyskaPHq1+/fppwIABmjt3rgoLC5WYmChJGjVqlNq0aaO0tDRJ0qxZszRt2jQtWbJE7du399+bExERoYiICNv6AAAAtYPt4Wb48OE6ePCgpk2bpry8PPXu3VurVq3y32S8Z88eBQWVn2B67rnnVFJSot/+9rcVtpOamqrHH3/8QpYOAABqIdu/5+ZC43tuAACoe+rM99wAAADUNMINAAAwCuEGAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjEK4AQAARiHcAAAAoxBuAACAUQg3AADAKIQbAABgFMINAAAwCuEGAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYBTCDQAAMArhBgAAGIVwAwAAjEK4AQAARiHcAAAAoxBuAACAUQg3AADAKIQbAABgFMINAAAwCuEGAAAYhXADAACMQrgBAABGIdwAAACjEG4AAIBRgu0uwC6lDbaptGEDu8uo9xzehnKe6GB3GQAAg9TbcHO81widbGx3FZCkhtueV1BxlN1lQIRNAGaot+Em/JsZahjUvUpzi93LVNJqiUL33yFX/vCA92UF/aSiSx6XL+x7Ndj1hJxFnQLehrfBThVdMk1BJy5Wg12Py+ELD3gbta2P8H0PqLDzeBV2vS/g7eD8IWzWDgRNoPrqbbhxFnVQsKPHWef91HaeSlotUdj3ExS+d2zA+7Gcx3WsW6J8YXvV6MuXFXw8JuBtlEZ8rqJLpstZ1FWNvnpRDm9EwNuorX04j3eX5Sys8jZqW0AzJWg22PW4rJACFXa9j7BZixA0UZPqU2Cut+GmKn5qO08nLv7zOQcCb4OdavTl36sdCI51HyVnUadzCja1tY9A/kWrrQEtULW5D+em96scNmtjQCNoAmdWlwNzqVVU5bmEm9OozYEgEPRRhj7KnamPqobN2hzQAlGb+yBo0geBuaJiT9Xn1ttws//oTwourfw/HK5OCxV+8Tz99NVYHd2ZKKnql00kScGFioh7QM7wb3R8/fM6erRjwNtwNtmiiG73yVvQUUcz0/VjqSPgbdDHz+jDjz7K1f4+WlZpG65OCxXeasnPfdwb0P4llffh2qvj6/+qo0d7BrwJZ5Mtiug+Xd6Czjqe+ZwOlzYMeBv08bPz2EfQ3n/KEVy112jIxa8r7JLlOrHrNp38/r8CrkHOEwrvNUPOhrkq+nyKfMc6BryJoEbfqEHM0/IWRuunLx7VSd9uSdOrtG69DTcnvZa8pb5Txht2WaTwbuk6vnWMCrffI+nUOWfiCC5Ukyv+oKDG3+rwh8+p9Ej3gLcR3HSLGg96QCcLLtXR9emySsMD3gZ9lKGPcvRRjj7K0Ee5etFHQXSV+wi7ZPnPffw+oP1LP/dx5R8U1PAHHf5woUqPBB7QgptuUdOYmTpZ0ElH18+XVdpQvhBvldevFV/il56ervbt2yssLEyxsbHKyso64/xXX31VXbp0UVhYmHr27Kl33nmnRupo2GWRInqkn/MBDY78Vkc+fK76B/RXD6i04FL/AQ0UfZShj3L0UY4+ytBHOfooZ0oftoebZcuWKTk5WampqcrOzlZMTIyGDh2qAwcOVDr/448/1ogRI3T33Xdr8+bNGjZsmIYNG6atW7eeUx2mHFD6KEMf5eijHH2UoY9y9FHOlD4kyWFZllWtNWtIbGys+vfvr3nz5kmSfD6foqOjNW7cOKWkpJwyf/jw4SosLNRbb73lH7viiivUu3dvLViw4Kz783g8ioyM1L1vvKCgk30lmXNA6aMMfZSjj3L0UYY+ytFHubrQhy9kkxYOu1sFBQVq3PjM38Jr6z03JSUl2rRpkyZPnuwfCwoKUnx8vDIzMytdJzMzU8nJyRXGhg4dqjfeeKPS+cXFxSouLvY/LigokCT5grdIksLavaWQdit0ZPNvdGLP5VLIpoB6cDhPqGHP2fI6ftCRjx6S73hJwNsIarRLDWJmq+iH1jq+5V5Zju1SSECboI+f0Uc5+ihHH2Xooxx9lKsrfXiDdkiSqnROxrLRvn37LEnWxx9/XGF80qRJ1oABAypdJyQkxFqyZEmFsfT0dKtly5aVzk9NTbUksbCwsLCwsBiw5ObmnjVfGP9pqcmTJ1c40+Pz+XT48GE1b95cDofDxsrOjcfjUXR0tHJzc896eg7nF8ei9uBY1B4ci9rFhONhWZaOHTum1q1bn3WureGmRYsWcjqdys/PrzCen5+vqKjKv0ExKioqoPkul0sul6vCWJMmTapfdC3TuHHjOvtCNQ3HovbgWNQeHIvapa4fj8jIyCrNs/XTUqGhoerbt68yMjL8Yz6fTxkZGYqLi6t0nbi4uArzJWn16tWnnQ8AAOoX2y9LJScna/To0erXr58GDBiguXPnqrCwUImJiZKkUaNGqU2bNkpLS5MkjR8/XldffbVmz56tm266SUuXLtVnn32mhQsX2tkGAACoJWwPN8OHD9fBgwc1bdo05eXlqXfv3lq1apXcbrckac+ePQoKKj/BNHDgQC1ZskSPPfaYHn30UV122WV644031KPH2X/h2yQul0upqamnXHLDhcexqD04FrUHx6J2qW/Hw/bvuQEAAKhJtn9DMQAAQE0i3AAAAKMQbgAAgFEINwAAwCiEmzomLS1N/fv3V6NGjdSyZUsNGzZMO3bssLusem/mzJlyOBxKSkqyu5R6a9++ffqf//kfNW/eXOHh4erZs6c+++wzu8uqd7xer6ZOnaoOHTooPDxcl156qZ588smq/R4QzsmHH36ohIQEtW7dWg6H45TfXLQsS9OmTVOrVq0UHh6u+Ph4ff311/YUe54RbuqYdevWacyYMfrkk0+0evVqnTx5UkOGDFFhYaHdpdVbGzdu1PPPP69evXrZXUq9deTIEQ0aNEghISH617/+pa+++kqzZ89W06ZN7S6t3pk1a5aee+45zZs3T9u2bdOsWbP0zDPP6C9/+YvdpRmvsLBQMTExSk9Pr/T5Z555Rs8++6wWLFigTz/9VA0bNtTQoUN14sSJC1zp+cdHweu4gwcPqmXLllq3bp1+9atf2V1OvXP8+HFdfvnlmj9/vp566in17t1bc+fOtbuseiclJUUbNmzQRx99ZHcp9d7NN98st9utF154wT/23//93woPD9f//d//2VhZ/eJwOLRixQoNGzZMUtlZm9atW+uhhx7SxIkTJUkFBQVyu91avHixbr/9dhurrXmcuanjCgoKJEnNmjWzuZL6acyYMbrpppsUHx9vdyn12sqVK9WvXz/deuutatmypfr06aNFixbZXVa9NHDgQGVkZGjnzp2SpM8//1zr16/XDTfcYHNl9dvu3buVl5dX4b9VkZGRio2NVWZmpo2VnR+2f0Mxqs/n8ykpKUmDBg2qd9/QXBssXbpU2dnZ2rhxo92l1Hu7du3Sc889p+TkZD366KPauHGjHnzwQYWGhmr06NF2l1evpKSkyOPxqEuXLnI6nfJ6vXr66ad155132l1avZaXlydJ/m///4Xb7fY/ZxLCTR02ZswYbd26VevXr7e7lHonNzdX48eP1+rVqxUWFmZ3OfWez+dTv379NGPGDElSnz59tHXrVi1YsIBwc4EtX75cL7/8spYsWaLu3bsrJydHSUlJat26NccCFwyXpeqosWPH6q233tKaNWvUtm1bu8updzZt2qQDBw7o8ssvV3BwsIKDg7Vu3To9++yzCg4OltfrtbvEeqVVq1bq1q1bhbGuXbtqz549NlVUf02aNEkpKSm6/fbb1bNnT40cOVITJkzw//gx7BEVFSVJys/PrzCen5/vf84khJs6xrIsjR07VitWrNAHH3ygDh062F1SvXTddddpy5YtysnJ8S/9+vXTnXfeqZycHDmdTrtLrFcGDRp0ylci7Ny5UxdffLFNFdVfRUVFFX7sWJKcTqd8Pp9NFUGSOnTooKioKGVkZPjHPB6PPv30U8XFxdlY2fnBZak6ZsyYMVqyZInefPNNNWrUyH+tNDIyUuHh4TZXV380atTolPucGjZsqObNm3P/kw0mTJiggQMHasaMGbrtttuUlZWlhQsXauHChXaXVu8kJCTo6aefVrt27dS9e3dt3rxZc+bM0e9+9zu7SzPe8ePH9c033/gf7969Wzk5OWrWrJnatWunpKQkPfXUU7rsssvUoUMHTZ06Va1bt/Z/osooFuoUSZUuL774ot2l1XtXX321NX78eLvLqLf++c9/Wj169LBcLpfVpUsXa+HChXaXVC95PB5r/PjxVrt27aywsDDrkksusaZMmWIVFxfbXZrx1qxZU+n7w+jRoy3Lsiyfz2dNnTrVcrvdlsvlsq677jprx44d9hZ9nvA9NwAAwCjccwMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCiEGwAAYBTCDQAAMArhBsAFN3jwYCUlJVV5/uLFi9WkSZPzVs+5+O677+RwOJSTk2N3KQB+RrgBAABGIdwAQC1UUlJidwlAnUW4ASCp7FLRuHHjlJSUpKZNm8rtdmvRokUqLCxUYmKiGjVqpI4dO+pf//pXhfXWrVunAQMGyOVyqVWrVkpJSVFpaan/+cLCQo0aNUoRERFq1aqVZs+efcq+i4uLNXHiRLVp00YNGzZUbGys1q5dW+Xaf7k09Prrr+uaa65RgwYNFBMTo8zMTP+cxx9/XL17966w3ty5c9W+fXv/47vuukvDhg3TjBkz5Ha71aRJEz3xxBMqLS3VpEmT1KxZM7Vt21YvvvjiKTVs375dAwcOVFhYmHr06KF169ZVeH7r1q264YYbFBERIbfbrZEjR+rQoUP+5wcPHqyxY8cqKSlJLVq00NChQ6vcP4CKCDcA/F566SW1aNFCWVlZGjdunB544AHdeuutGjhwoLKzszVkyBCNHDlSRUVFkqR9+/bpxhtvVP/+/fX555/rueee0wsvvKCnnnrKv81JkyZp3bp1evPNN/Xee+9p7dq1ys7OrrDfsWPHKjMzU0uXLtUXX3yhW2+9Vb/+9a/19ddfB1T/lClTNHHiROXk5KhTp04aMWJEhaBVFR988IF++OEHffjhh5ozZ45SU1N18803q2nTpvr00091//3367777tPevXsrrDdp0iQ99NBD2rx5s+Li4pSQkKAff/xRknT06FFde+216tOnjz777DOtWrVK+fn5uu222yps46WXXlJoaKg2bNigBQsWBFQ3gH9j9y93Aqgdrr76auvKK6/0Py4tLbUaNmxojRw50j+2f/9+S5KVmZlpWZZlPfroo1bnzp0tn8/nn5Oenm5FRERYXq/XOnbsmBUaGmotX77c//yPP/5ohYeH+39B/fvvv7ecTqe1b9++CvVcd9111uTJky3LsqwXX3zRioyMPG3tu3fvtiRZf/3rX/1jX375pSXJ2rZtm2VZlpWammrFxMRUWO/Pf/6zdfHFF/sfjx492rr44ostr9frH+vcubN11VVXnfJ3eeWVVyrse+bMmf45J0+etNq2bWvNmjXLsizLevLJJ60hQ4ZU2Hdubq4lyf+rzFdffbXVp0+f0/YIoOqCbU1WAGqVXr16+f/Z6XSqefPm6tmzp3/M7XZLkg4cOCBJ2rZtm+Li4uRwOPxzBg0apOPHj2vv3r06cuSISkpKFBsb63++WbNm6ty5s//xli1b5PV61alTpwq1FBcXq3nz5tWuv1WrVv5au3TpUuVtdO/eXUFB5Se13W63evTo4X/8y9/ll7/BL+Li4vz/HBwcrH79+mnbtm2SpM8//1xr1qxRRETEKfv79ttv/b337du3ynUCOD3CDQC/kJCQCo8dDkeFsV9CjM/nq7F9Hj9+XE6nU5s2bZLT6azwXGVh4EzOVGtQUJAsy6ow/+TJk2fcxi/bqWwskL/B8ePHlZCQoFmzZp3y3C8hTJIaNmxY5W0COD3uuQFQbV27dlVmZmaF0LBhwwY1atRIbdu21aWXXqqQkBB9+umn/uePHDminTt3+h/36dNHXq9XBw4cUMeOHSssUVFRNVbrRRddpLy8vAq11uR303zyySf+fy4tLdWmTZvUtWtXSdLll1+uL7/8Uu3btz+lRwINUPMINwCq7Q9/+INyc3M1btw4bd++XW+++aZSU1OVnJysoKAgRURE6O6779akSZP0wQcfaOvWrbrrrrsqXPbp1KmT7rzzTo0aNUqvv/66du/eraysLKWlpentt9+usVoHDx6sgwcP6plnntG3336r9PT0Uz75dS7S09O1YsUKbd++XWPGjNGRI0f0u9/9TpI0ZswYHT58WCNGjNDGjRv17bff6t1331ViYqK8Xm+N1QCgDOEGQLW1adNG77zzjrKyshQTE6P7779fd999tx577DH/nD/+8Y+66qqrlJCQoPj4eF155ZWn3Fvy4osvatSoUXrooYfUuXNnDRs2TBs3blS7du1qrNauXbtq/vz5Sk9PV0xMjLKysjRx4sQa2/7MmTM1c+ZMxcTEaP369Vq5cqVatGghSWrdurU2bNggr9erIUOGqGfPnkpKSlKTJk0qBD0ANcNh/edFaAAAgDqM/2UAAABGIdwAAACjEG4AAIBRCDcAAMAohBsAAGAUwg0AADAK4QYAABiFcAMAAIxCuAEAAEYh3AAAAKMQbgAAgFEINwAAwCj/HxIMp4EdeUz0AAAAAElFTkSuQmCC", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAiIAAAGwCAYAAACU8g7/AAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjguMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/H5lhTAAAACXBIWXMAAA9hAAAPYQGoP6dpAABZ2klEQVR4nO3deXwTZeIG8GdyJ22S3iflKOUGOQUFuRTFi1Vx1fVAVFR2BVdlvdBVvFZZ9eexLup6666Ksp67KiqsHAoe3CiH3OVoC72StrmT+f3xJmlSWmhL00mb5/v55DOTyUzyZiidp+81kizLMoiIiIgUoFK6AERERJS4GESIiIhIMQwiREREpBgGESIiIlIMgwgREREphkGEiIiIFMMgQkRERIrRKF2AYwkEAjh06BDMZjMkSVK6OERERNQMsiyjpqYGeXl5UKmOXecR10Hk0KFDKCgoULoYRERE1Ar79+9Hly5djrlPXAcRs9kMQHwRi8WicGmIiIioOex2OwoKCsLX8WOJ6yASao6xWCwMIkRERB1Mc7pVsLMqERERKYZBhIiIiBTDIEJERESKYRAhIiIixTCIEBERkWIYRIiIiEgxDCJERESkGAYRIiIiUgyDCBERESmGQYSIiIgUwyBCREREimEQiSDLMp5bukPpYhARESWMuL7pXXvz+mV8uP4gbj6jV4uPLbE58ZfPtsKoVSPHasD5J+WhT87x7zpIRESUyDpEEBn31/8hKyMVqSadeCTpkGrSIi1Jh7QkHTLNemSZDcg065GkV0OnVjXrjn8NOT1+GLTq8PPtpTVYV1yFnpnJGNkjLbzd7vKipNqFoqxkqCRxd8GKWg/2VtThid8Oxu4jdbjyle/xyezTkJ9ibJNzQERE1Bl1iCBS6fCiuqy22ftrVBJMOjWS9JropU4Dk16DJJ0aJp0GSXqxLEg14vzBeXB6/TDpRBBZ/usRPPHlNsw4rQdeXrkba/ZV4qYJRVhXXIW7P9iEUwrTUVHrwfriKqyaewacXj9STTr0y7WgX64F/167H1sO2RlEiIiIjqFDBJFFvz8FXpURlQ4Pqh0eVNZ5UO3worLOg4o6Nw7b3ThS60a1wwsA8AVk2F0+2F2+Zr3/8G6p4SBiDNaIvPbtHsybMgAnd0/Dmf1zcPqTy3DThCI8/fWvmH/xSRjWNRUb91fj9/9aC0DUpnh8Aewoq0GJzYXNB23oG2ya2bi/Ck5vIBiE1OGlSauGRs1uOkRElLg6RBDpl2uFxWI57n4+fwAOrx8Otx91Hl/90uNDndsfvfT44XCLZdc0UWsR2TRTZnehW5oJAJCs18AXkBEIyNh9pA4D8kRZCjOTwsHF6fWjzO7CC8t3waRTY8EVw1CQZkJFrRuX/uN7uH2BRsus16ga1NgEl03W6Bz/9dY2TREREbW3DhFEmkujVsGiVsFi0LbqeKfXD2OwaaZnZjI2HrDhzP4GbC2xI9dqgEolIddqQHGFA72yzdhaUhMOLi6vH6cUpmP+xSdFvecLy3ehT44ZdW4fHB4/6oLhxx+QAQBuXwBunweVdSfwxRtoadNUeKlTN/m6QaOGSsVw01AgIEMK9hMiIqKW61RB5ESpVRK6p4takDln9cYtC9fj4w0HsbOsFo9cNBAAcMukXrj53fUY2ysDOw7XIkkfrBFp0NE15M/n9T9qmyzL8PgDETU2/qig4vA0rNFp3uuhWpeWNk01hyQBJm3TQeV4QUa83vmapt5avRdev4wbxhUqXRQiog6JQSTCkIIUDClIASBqRP4z+zQcqXEjLUkXvmCO7ZWJhTdacaDKicM1Lny+uVQc2zUFPbOSm/U5kiRBr1FDr1EjNUnXZuVvbdPU8V4HAFkG6jx+1Hn8ONJmJW5+01TS8YJPRADSadov3NR5/Egxta4Gbt4nP2PlznIk6zXIsxpx1zl90SMjqY1LSEQU3xhEjkGSJGRZDFHbPt14CGv2ViJZr8GSrWV4+rIhAIC+OcfvwxJrJ9o01ZhAQIbT64ejYWBpdpA5ukYn1k1TWrXUzBqaY7+erNcgSS+WTYWbWrcPXVJFH6M1eyvx7NId8AdknDsoF1eO6gpJklBqc+G9n/bDbNBgXO8MABKKspJR5fDirrP74qz+2Xjvp/2Y++EmLLzx1LY7EUREHQCDSAudPygXfXPMqHX7MHN8T1iNbXfRj0cqlSRqJPQaAPo2ec9YN015/TJsTi9sTm+blBcAdGoVkvTqcDC5//z+GF2UgTq3D0k6DWrdPtz63ga8ff0oZFsMuPrVH9EjIwnDu6Vi2qs/YPbpRTBq1bjhrbW4dEQBirKSUef2IVmvgSRJGNI1BS+t3N1m5SUi6igYRFpIpZLQO5szpp6IeGiaqj1uU5V4zekVTVMefwAeRwBVwSHiHr8IPbVuH5L0Gmw6UI1B+VZ0SxdNK1OG5OGH3RUIyDL651lwwZB8AMC3O8uRHOxXVOv24YufS7C1xI7/birB+SflAQCKK+rw9o/FSNZpkGyor5URy/owFNqm7eD9bIgosTGIUKcRi6Ypnz8g+sa4RTipDQaU/sEh3KFajRqXN2pUkcPtg16rhs3pRVpE2BI1K+K/XZ3Hh56ZyUhL0uH+Kf0xNNg/6Ykvt+M/m0qaXUadRhUMJWok67UwGzQw6zUwG0SQMRu0SNZrYAk912uD28V6aD8GGiJSAoMI0TFo1CpYjaomm+Dq3H4k6dUozEzC3I82Y+P+aqSadFi09gAWXDEMBq0Kzy3dCYfHB39AxtJth3Fy8HYBdW4/zh6Yg1xr/ey7P+2pRGFmMmac1iMi+IjwUxvxvDaiKcrjC6Ay3M/G2ervatCqkKzX1geWo0JLMNQEnycHn5sjnoeamoiImotBhOgEXD+2B7ItBiTpNXj56hF4ecVuOLx+3H9+//BND2eM7YFrXvsJaUk6ZFv04RqSULNOpJN7pIWDyvF4g/1saoNNSTWu+mWt24saly/i4UVt8LUad/B58LVQ85PLG4DL60Z5rbvV50MlAWaDFhajBhaDVjxC60bx3GrUhNctxujXk3RqBhmiBCPJsiwrXYim2O12WK1W2Gy2Zs2sShSPDttdyLIYUO3wYOrzq/Dx7DGwGLQ4WO1ErsWg+ERxPn+gPqS4fMF1EVzsLl8wsHij9ol+LkKPL3Div0pUEiJCSlNhpj7ImIN9aMwRfWn0mvieWXjxzyU4rVcmkvX8O5A6r5Zcv/k/gSjGHv5sK6rqPPD4Arjz7L7hPizxckNEjVqFFJMOKabWdxyWZRluXwB2lxd2py+49IqJ9ZzeRrfbnF7URLzm8QcQkIFqhzd836jWUKskJOnUMBu0USOdxFw0odASvd1i1CLFpEWqSSvOhVEbs8n2nvvfTpxSmN6qY2VZjuuQRdQaDCJEMfbc5UOVLkLMSZIEg1YNg1aNrFYMKgsHmWAwsTUzzNQ4vfX9aIKT7/nbaGZhs16DlCQtUoMhLdUk1q1GEVhSk6K3Z5r1jc6u3JDN6Q2HUX9AxuEaF9KSdNBrjn2sLMsYOO9LFGUlQ6tWYVL/bNw4tlDxGjWiE8UgQkSKiwoyDSYRbK5AQIbD60etK7KTr69BJ19/I9tEaKl2iLt6h+afqXGL/jT7K5vfATgjWYe8FCPyrEaxTDEgL8WI/BQDBhekAhCzFKtUEo7UuHHdGz+hR0YSdh2pxU0TinDeSbnw+QN4eeUebNxfjTP7Z0OnUeHM/tnw+AMw6jT4ZPZpqHF5cVVwrprJA3Jadb6I4gWDCBF1CiqVFB65cyL8ATEhXpXDg2qHB1V1XlQ7vWLd4UGVwxsOLaH1yjoP3L4Ayms9KK/1YNMBW9R7dkk14tu7Tg/ek0o0+by4fBcuHpaPa8b0QGWdB+f/bSXOHpiD17/bi91HavGXiwbi/TUH8OzSX7Fx3lmorvMiLUnUpJgNWvTPNaPM7jqh70oUDxhEiIgiqFUS0pJ0UfO/HI8sy6h2eHGw2okSmwuHqp3iEVxPCQ7/rnJ4wu/780EbLju5AACQlqSDTqMKT3L31KVDkJ6sx4zTemDBNzuh16hR6fCgxuXDKyt3w+nx4+sth/H78T0BAEu3lqGizoOMZB2yzAZkWfRIT9JDzWYb6gAYRIiITpAkSUhN0iE1SYeB+dYm96us84Q7BVuMWlTWeQCIodx1Hj/Meg1c3kB4WLfL5w8Hl6o6D/JTjMg06yHLwH9uHoNcqxGlNhdmvbMOLm8g6rPUKqk+mJj1yLKElnpkB8NKtsWA9IibehIpgUGEiKid1Lp9yLaIezZdM7o7Hvt8K648pRu+2FyCa8d0h0ol4bReGfjn6r2YOb4n/rF8V/juzpV1HvTKNodvFwAANS4vPt5wEOcMzEWVw4PyWjfK7G5U1LrhD8gos4vnxyJJQHqSHtkWvQgqZgOyLXpkBoNLdnCZadZz9l2KCc4jQkSkkD3ldViztxLdM5JwcncxkZ3L68ff/7cT28tqMCjfil8O2fCPaSPwysrdqHJ4cMfkvsd9X58/gMo6D8rsbhyucYWXh2vcOGwPLd04EgwszZWepENuiiHcGTc/pb5Tbn6KERnJeo7iIQAtu34ziBARxZGfD9qQZdHDrNdi3qc/46QuKbjqlG5wef3w+gMwt+G9lPwBGZV1HhFSQmHF7kZZeOnGEbsLR2rd8PqPf6nQqiXkWiNHC4WCihg5lGs1HjWbMHVOnNCMiKiDcnj8eOg/W1Dr9mFEt1RcPrIrAISHN7cltUpCZrDZZUBe0/sFAjKqHKKGpdTuxMHqiA651U4cqnah1O6C1y+juNKB4kpHk++VYtJG1KiIwNI1zYRu6Unolm5iUElArBEhIqIT5vMHUFbjDoeTgxEhJfS8phmTzGWa9eieLoJJ/TIJ3TJMbXpnbYot1ogQEVG70qhVyA82xzTF7vKiJCKYhJb7KhzYV1GHKocXR2rcOFLjxk97q446Pi1Jh27pJhFMGixTTFpOf99BMYgQEVG7sBi0sORow3embsjm8GJfZR32VjhQXCGW+4LLIzVuVNaJyePWF1c38t4adM9IalCTIpYZyTqGlDjGphkiIop7tW4fiiOCiVjWYV+FAyW2Y88wm6RTi2CSYUKPjCQMyrdiUJcU5FkNDCgx0i6jZlasWIEnnngCa9euRUlJCT766CNceOGF4devueYavPnmm1HHTJ48GYsXL272ZzCIEBHR8bi8fhRXOrC3XASTUEDZW1GHg9VONHWVy0jWhUPJSflWnNTF2up7HVG0dukjUldXh8GDB+O6667D1KlTG93n7LPPxuuvvx5+rtfrW/txREREjTJo1eidbUbv7KObfNw+Pw5UOUUNSrkDv5bVYPNBG7aX1qC81oNvth/BN9uPhPfPtugxKD8FJ3WxYlAXK07KtyI9mdeuWGp1EDnnnHNwzjnnHHMfvV6PnJzm3xnS7XbD7a6fBdBut7e2eERERNBr1OiZmYyemclR211eP7aV1mDzgWpsOmDDpgM27DhcE5yNtgxLtpaF981PMWJIQQqGdk3BkIIUDMy3tvlQ6kQW086qy5YtQ1ZWFlJTU3H66afjkUceQXp6epP7P/bYY3jwwQdjWSQiIiIYtGoMKRDBIsTh8WHLITs2HbBh80EbNh2oxu5y0bxzsNqJzzaXAAA0Kgl9c804uXsaxvTMwMjCNA4tPgFt0llVkqSj+ogsXLgQJpMJPXr0wK5du3DPPfcgOTkZq1evhlrdeJJsrEakoKCAfUSIiEgRNS4vNh+wYf3+amzYX431xdUor42+f49aJWFQvhVjitIxpmcGhnVLTfgak3af4r2xINLQ7t270bNnTyxZsgRnnHFGs96XnVWJiCieyLKMg9VOrC+uxurdFVi1sxx7K6JnktVpVBjRLRVjijIwumc6BuVbE+4Ox3E5oVlhYSEyMjKwc+fOZgcRIiKieCJJErqkmtAl1YQpg8W8+AernVi1sxyrdlXgu53lOFzjxqpdFVi1qwIAYNZrMKowHaN7pmNMUQZ6Zydz2HCEdgsiBw4cQEVFBXJzc9vrI4mIiGIuP8WIS0YU4JIRBZBlGbuO1IZDyepdFbC7fFiytb4DbEayPhhK0jG6ZwYK0kwt+rziCgdyUwzQdpJallYHkdraWuzcuTP8fM+ePdiwYQPS0tKQlpaGBx98EBdffDFycnKwa9cu3HnnnSgqKsLkyZPbpOBERETxRpIkFGWZUZRlxtWndoc/IGPLITu+21WO73aW46e9lSivdePTjYfw6cZDAICCNCNOK8rAlaO6YWC+9Zjvv7/Sgemv/4iXpg1Hr0aGK3dEre4jsmzZMkycOPGo7dOnT8cLL7yACy+8EOvXr0d1dTXy8vJw1lln4eGHH0Z2dnazP4N9RIiIqDNx+/zYUFyN73aJ/iUb9lfDF6i/DI/vnYlZE4swskfaUcdW1Xlw+cvf4+5z+mJCn6z2LHaLtXtn1VhhECEios6s1u3DT3sr8cn6g/jPphL4g6FkYL4FV4zshktGdIFWrYLT48dVr/6Ay0d2xW+Hd1G41MfXkut352hgIiIi6oCS9RpM7JOFZ343FN/8aQKuGNUVOrUKPx+045DNCa1aBZ8/gJvfXY/T+2bht8O7YOfhWny47gAOVDmO/wEdAGtEiIiI4khlnQefbyrB70YWQK2ScM9HP0OrlvDgbwbgn9/vwxur9uKCwfn4blc55pzZG6cUNj1RqFLicvguERERHV9akg5XndoNAODxBdA93YTrxxZiS4kdb67ai4/+MAZWkxa/G1mAm99Zj/d/f6rCJT4xDCJERERxSqdRYeb4ngCAZduP4JrR3WE1ienkS2wudIbpSNhHhIiIqAPonW3GttIayLKMUpsLd/57I34/vieqHR6880MxPtlwED5/QOlithiDCBERUQcwqV8WsswGnPX0Clzyj1X4w4SeKMxMwtTnV6HU5sSuw7W496OflS5mi7GzKhERUQciyzJkGVCpJFzx8ve4ZnR3nDUgBwBw1783IS/FiNmnF0GtkmBzeMNNOe2Jw3eJiIg6KUmSoFJJ8AfEDfhCIcTrD2BtcRV8gQDUKgn7Kx2Y9PRyPLtkB+K4zoFBhIiIqKOpdfugksR9btbuq0SNy4u7P9iMYV1TMOfM3qis82D66z/izsl9UFHnxhNfble6yE1iECEiIupgXl25B19vKcPfLh+Kr7aU4exnVqLK4cGjFw2CyxvAdW/8hIuHdcElIwrw4G8GYOnWw9hfGZ8ToHH4LhERUQdzw7geuPmd9fhhTyWqHB4UZibh71cMBQDMfmcd/AEZ20tr4PL6YXd6UVHnjtuhvgwiREREHYxJp8HLV4/A8l+PwB+QMaFPJtQqCXM/3AydRoWPZ43BB+sO4KynV6DK4cGsiUXokmpSutiNYhAhIiLqgFQqCRP71t+FNxCQMSDPgktGiKnhLx1RgN8MzoPbG1Bk5ExzsY8IERFRJ6BSSZh2ancYtOrwNoNWDatJi/JaN3aU1ShYuqYxiBAREXVi//fVdox6dCluensdvHE48yqDCBERUSd2/WmFsBq12HG4Fv/6fp/SxTkKgwgREVEnZjVpMefM3gCAfyzfDY8vvmpFGESIiIg6uUtGdEG2RY9SuwufbjykdHGiMIgQERF1cnqNGlef2h0AsGjNfmUL0wCDCBERUQK4aGg+JAn4YU8lDlU7lS5OGIMIERFRAshLMWJoQQoAYMWvR5QtTAQGESIiogQxrncmAGDFDgYRIiIiamdje4kgsnpXBWRZVrg0AoMIERFRghiYb4FGJaHK4cUhm0vp4gBgECEiIkoYeo0aRVnJAIBfDtoULo3AIEJERJRA+udaAAA7DtcqXBKBQYSIiCiBdEkzAQAOVMXHEF4GESIiogTSJdUIADgYJ3OJMIgQERElkC4pIogcqHIoXBKBQYSIiCiBZFn0AICKWo/CJREYRIiIiBKI1agDANhdXgQCys8lwiBCRESUQKxGLQBAloEal0/h0jCIEBERJRSdRgWTTg0AqHYq3zzDIEJERJRgLAZRK2J3skaEiIiI2plBKy7/bp9f4ZIwiBARESUcvUY0zXh8AYVLwiBCRESUcHSaUI0IgwgRERG1M72GTTNERESkEL2WNSJERESkkFAfEQYRIiIianc6NWtEiIiISCGSpHQJ6jGIEBERkWIYRIiIiBKVzJveERERUTtj0wwRERERGESIiIhIQQwiRERECUr5HiIMIkRERAlHQvx0EmEQISIiIsUwiBARESWoOBi9yyBCRESUcOKnZYZBhIiIiJTDIEJERESKYRAhIiJKUHIcdBJhECEiIkowcdRFhEGEiIiIlMMgQkRElKCUb5hhECEiIko4UhzdfrfVQWTFihWYMmUK8vLyIEkSPv7446jXZVnG/fffj9zcXBiNRkyaNAk7duw40fISERFRJ9LqIFJXV4fBgwdjwYIFjb7++OOP429/+xtefPFF/PDDD0hKSsLkyZPhcrlaXVgiIiJqO3EwaAaa1h54zjnn4Jxzzmn0NVmW8cwzz+DPf/4zLrjgAgDAW2+9hezsbHz88cf43e9+19qPJSIiok4kJn1E9uzZg9LSUkyaNCm8zWq1YtSoUVi9enWTx7ndbtjt9qgHERERta346SESoyBSWloKAMjOzo7anp2dHX6tMY899hisVmv4UVBQEIviERERUZyIq1Ezc+fOhc1mCz/279+vdJGIiIg6rTjoIhKbIJKTkwMAKCsri9peVlYWfq0xer0eFosl6kFERERtK45G78YmiPTo0QM5OTlYunRpeJvdbscPP/yAU089NRYfSURERB1Qq0fN1NbWYufOneHne/bswYYNG5CWloauXbvi1ltvxSOPPIJevXqhR48euO+++5CXl4cLL7ywLcpNREREJygebnrX6iCyZs0aTJw4Mfx8zpw5AIDp06fjjTfewJ133om6ujrceOONqK6uxmmnnYbFixfDYDCceKmJiIio1eKoZab1QWTChAnHTFKSJOGhhx7CQw891NqPICIiok4urkbNEBERUWJhECEiIiLFMIgQERElmE5x910iIiKiE8UgQkRElKDiYPQugwgREVGiiZ+GGQYRIiIiUhCDCBERESmGQYSIiChByXFw/10GESIiokQTR51EGESIiIhIMQwiRERECYrDd4mIiKjdSXHUNsMgQkRERIphECEiIkpQcdAywyBCRESUaOLonncMIkRERKQcBhEiIiJSDIMIERFRguLwXSIiImp3cdRFhEGEiIiIlMMgQkRElKB40zsiIiJqdxy+S0RERAQGESIiIlIQgwgREVGC4vBdIiIiane8+y4RERERGESIiIhIQQwiEbwBLz7a8ZHSxSAiIoqpeBq+q1G6APHE4/fgn1v/iYt6XdSi44rtxfjT8j8BACpdlTBpTNCr9ehq6YqnJjwVi6ISERF1CgwiEXwBH9SSOvx8t2031pSuQRdzF5yaeyqkYITcY9uDX6t+xYjsESh3lqNPWh8smrIIADBn2RxcWHQhxnUZp8h3oPa14fAGWPQWFFoLlS4KEVGH1CGCyPVfXg+LxQKj1gijxgidSgeNSlP/kDTRz1UaaFXao15Tq9TiNUkbtW+yNhkDMgYgIAfCQWRt2Vo89sNjuHbgtfhk5yf4vuR7zBk+B6sPrcbTa5/GRb0uwvwf52NN2Rp8c+k34bI2DDPUuS0/sBw9U3q2OIg8u+5ZfHvwW/gCPlS7q5FhzAAA3Dz0ZoZYImo3chyM3+0QQeTnip+hdsTu4t47tTc++M0H8Mv+cIh4d9u7+OOwP2Jcl3E4s9uZOPPfZ+KWobfgX1v/hXtG3YMhWUMwqeskXPRpdDNOQA5AJYmuN8+sfQarS1ZDI4kQpJaCD5UaKkkFjaSBSlKJgCRpoFKpovYJrasklQhSEeuRS7Wkjlo/6niV+KzI7c0+vkG5G5YpdIwUTw2O7cgv+6GRxH8jf8CPFQdWwOaxYXTeaGSZsgCIvkc/lvwInVqHnik94Q/4ccuwW3DLsFtQbC/GLd/cEq5Ro87vg18/wMW9L1a6GJTg4ulXdocIIo+PfRwqowoOnwNOnxMevwe+gA8+2SeWkY/GtgV88MreqOf+gD+8b5fkLgCCtRkqEUQqnBXIS8oDAOjUOpg0Jjh8DhxxHEFestiebkyHTqWLKqtP9kGjEqf1QM0BbKnY0l6nSVHHCjSh0NVYOAoFscbCTmPHNwxAxwpLDUNXYwGqyeOPEwSzTFkw68zwB/zh4HnnijuRYcxAr9RemPHlDLx45ovIS8rDbd/chvzkfKQaUvHYj4/hnO7n4IaTbgAQDDKqDvHfkNrIy5tfbnEQqXRVYubXMwEAtZ5aaNVa6NV6pOpT8dJZL8WimETtpkP8BhxXMA4WiyXmnxNZm9EzpSc2HtmIotQilNSWAADMOjOKUoqw6cgmTOo2CZuObDrqPSIvTNefdD1+U/Qb+AN++OXgI2I9IAdEKIpYD8iBqP2itjU4Puq9AgH45OC+wZAVWg/vHzwmIIt9w+vBzzje8QE5cMxzF5AD8Aa8MfiXiT/zx87HeYXnhZvzDtUewm7bbjw5/klIkgSH14EPd3yIs7ufjRpPDeaOmgsAqHJVhcMuEP3zcsRxBDO+nHFUOGsY2hoNcI2EqKjasyaOb7h+VBCL2CcUxhoGt6iQGPEeDWvcGgt1iVqTFnKw9iB+LPkRmaZMjMkbEz4fJbUl2FK5BYMzB8PutqOHtUe41uyBVQ9gePZwTOk5RcmiUycRBy0zHSOItJeAHIBZawYAXD/oevxp2Z/w3aHvsNe+F/eMugcA8IfBf8DtK27Hhzs+RJohDQaN4aj3CDXv9E3ri75pfdv3S8SQLMtHhaCGASoq1ASaEaAaBKGG4eh4x4fXm3jPo8Jc4DihLBjAGj0+4li9Wg+gvhat3FmOvOS88IUkJykHe+17ccR5BPnJ+eFzmJOUEw4eAKKaA91+N/bY97Tjv6jyGgtCTQaa0D6RzY2N1GA1FnraKjw1Fd6Ot49WrUVuUi78AX94Rstfq37FHcvvwLUDr8VXe7/CN8Xf4L5T78MvFb/gz9/+GRcVXYTFexZjdclqLL90efjnJPJnhhLHodpD4dr4thE/fwQwiEToYu6C5854DoC4YPzr3H+h0lUJq94arj7vYu6ChecthNPnxG7bbpS7yqPe46ExDyHNkNbuZW8PkiSJzr/QAPw9CKC+Fq3AXIBd1bvg8XugU+uw/vB69ErthR6WHthWtS3ciXnD4Q0YkTMifHzkRSXTmInXJ7/eZDhqbgBrLGiFQmFjAa/hZzX1OZEhrTnhMLK2ryk+2QfIAJrepVPINmVjySVLopri3t/+Pq4fdD2m9JyCKYVTMPmDybjNcxve3vI2bh12K8YXjEeNpwbj3hsXHV4D/nCt2j9/+Se+2PtFODyFO+k37LTf2PNgQNKqtOFO/c05PnJbY+8XtS5pEr7Wqy1UOCswa+ksfHRB55znikHkGCRJQroxPWrbLxW/4Jl1z8Cis2CffR8eGv1Q1OuRf/1S56dWiV/+qYZUTOs/DdO+mIYcUw7qfHW4ZdgtMGgMOLv72bjq86uQbkyHy+cK16YA0c2Beo0+KqR0Fo3VpB0rADWsjToqZDUIPU0GqYj9wrVuxwhRDWvmjhX6Ij/3eGEsIAeQpE0CED2qrtxZHv59oVapkW5Mh81jw2Hn4fB2s84Mq84adTH3yb5wB+ndtt3YXL65Pf85W0yn0kGr1oqlSivW1WI98jWNWgOdSlf/WnAZuX681/RqPQwaA3RqHQzq+qVeo4derYdOpeuQwWjlwZUYmz82Ju8dBy0zDCItNTBjIJ4/43k4fU5YdJYO+UNNbSfUZAcAV/a7ElN7TYXD60CaIS38s3HDoBswfcB0qKDCXSvvQndL9/AxPVN64i+n/aW9i92uWJNWzy/X12YUWguxuXwzhmUPQ7WrGuXOcmSbstHD0gM/V/yMotQi7LXtRa23Nvo9ImpEruh3BSYUTAiHooYd9r0Bb6PLhq8ftX+wc7/X7230/Zp67pN9R31nT8ADT8CDOtTF/gQfhwRJBJIGASXqodEfFWIMagNMWhNMGhOStEkwaU1I0ohlw+2xCDvL9y/HVf2vAiDCbKWrEhnGjKiaspaKp0sXg0gr6NQitRM1ZNSIuW4izf9xPipdlSh3liM3KRcn55wctX8Xc5f2LiYpJCAHkKxNBgBc1f8q3PK/W7DpyCbste/F7SNuh0alwYxBM3DrN7fif8X/g0FtgFVvjXqPyOa8Xqm90Cu1V7t/j6aEmhC9AW/9wy+WHr9HLAMeeP31y4avhUZFhrc1WDbcFn4/vwfugBtunxtuf/3D5XNBDv7dL0OGy++Cy++CHfaYnAONpIFRa6wPJ8GlUWuMem7SmMJBJkmThCRtEtKN6cgwZiDdkA6tWgtAzPi9pUJ0XN5WuQ13LL8DvVJ74YjjCOadOg9FqUWtKmcc5RBIcjzMZtIEu90Oq9UKm83WLqNmiGJBlmVUu6th1BiP6txMiU2WZVS4KmDRWcJ/3IR+Jds9dtjddty98m68fd7b4WP22/cjxZACs86sSJk7GlmW4Qv44PK76gNKRFhx+V3w+D1w+VxHBRiP3yOCi88Fh8+BOm8dHD4HHF7xiHzu8rvatNwzT5qJ2UNnY9XBVfh096eYP3Y+rvzsStw6/FacnHMy1h9ej39s+gdenPQi3H43vtjzBQJyAGd3PxsmrQkAMHflXPgDfowrGIex+WOjQu0vh2w472/fYs6ZvfHHM9o+zLbk+s0aEaIYkyQJqYZUpYtBcUiSpPCsuiF7bHvw6I+PIsOYge2V2zF35Nyo1wssBe1ZxA5PkiTRKVethRmxC2/+gL8+pEQs67x1IrT46sIBJjLU1Hnr4PQ6UeOtQbmzHJXOSvhkH0bmjAQgZm+e0GUCbG4bqtxV4RrVvKQ8HHYcBgDcvux2dDF3QVFKER5c/SCMGiPmnToPj572KLZVbsPTa5/G/635P3S3dMe4LuNwVf+rMCDPihyLgcN3iYgoWmFKIRacsQB2tx1phrSoeWcofqlVaph15hOuqQrIAdjcNqToUwAAG49sxKyhs+AP+KFVaSHLMiRJwqpDq3BSxkn4peIX2D123DXyLgBAtbsaP5b+GO6n0i+9H9x+N54Y9wS6Wbrh+5LvoZE08Adk1LjiY94nBhEiojijV+uRacpUuhikAJWkiqpBffvct8NhdGDGQDz242PoZumGt7e+jefPeB5rytZE9Tur8dRE3a+q2lWNPbY9GJI1BBqVJjwR3obiStR5/O30rY6t9V1uiYiIKKYia8QeGfMIJnWdBKPGiNcmv4bu1u7oau6K7VXbwzUp/9n9H0womBA+ZuXBlTgl75SjbiXx6YZDABDuyKsk1ogQERF1AJIkYWTuSIzEyPC2k3NOxtrDa3H9V9cjSZsEi84SNZ/VigMrMLFgIgDgk52fYETOCOQn58MfD51DglgjQkRE1EEt3L4Q0/pNw2uTX0NXc1dcVFR/R3hvwIvvS77HmPwxAIBFvy7C5f+9HNsrtytV3EYxiBAREXVQWaYs3LHiDlz+38tR7a7G7/r+LvxapbMS5xeeD6veKm5K6veiyl2FG766Ac5AhYKljsZ5RIiIiBKA3WPHjC9nYFvlNmRpTsKuzVfgljN64bYze7f9Z7Xg+s0aESIiogRg0Vnw+LjHoVFpcNi3CWrjXqWLBIBBhIiIKGH0sPbA+YXnAwA01vUKl0ZgECEiIkogk7tPBgBokrfHweBdBhEiIqKEMjhzMABApa2GO1B7nL1jj0GEiIgogZh1ZhikNACA3VeicGkYRIiIiBKOQSXuxOsK2BUuCYMIERFRwtFJ4uZ8bj+DCBEREbUzjWQAAHhll8IlYRAhIiJKOBKk4Jry42YYRIiIiBKU8jEkxkHkgQcegCRJUY++ffvG8iOJiIjoOOKpRkQT6w8YMGAAlixZUv+Bmph/JBERER2TCCJyIgQRjUaDnJycZu3rdrvhdrvDz+125XvzEhERdT7BIKJ8Dol9H5EdO3YgLy8PhYWFuPLKK1FcXNzkvo899hisVmv4UVBQEOviERERJRxJip+mmZgGkVGjRuGNN97A4sWL8cILL2DPnj0YO3YsampqGt1/7ty5sNls4cf+/ftjWTwiIqIElSBNM+ecc054/aSTTsKoUaPQrVs3vP/++5gxY8ZR++v1euj1+lgWiYiIKOFJ4TXlg0i7Dt9NSUlB7969sXPnzvb8WCIiIorQLS0JANA316xwSdo5iNTW1mLXrl3Izc1tz48lIiKiCKlJovUhy6x8K0RMg8jtt9+O5cuXY+/evVi1ahUuuugiqNVqXH755bH8WCIiImqGTt9H5MCBA7j88stRUVGBzMxMnHbaafj++++RmZkZy48lIiKiY1BJoh5CjoPxuzENIgsXLozl2xMREVErSHE0aob3miEiIkowoXlE4qFGhEGEiIgowbBGhIiIiBTHGhEiIiJqd+GmGdaIEBERUXtTBS//DCJERETU7thZlYiIiBTHGhEiIiJqd+FRM6wRISIiovbGzqpERESkGNaIEBERkWLMOjMyjZlI0iYpXRRIcjzEoSbY7XZYrVbYbDZYLBali0NERETN0JLrN2tEiIiISDEMIkRERKQYBhEiIiJSDIMIERERKYZBhIiIiBTDIEJERESKYRAhIiIixTCIEBERkWIYRIiIiEgxDCJERESkGAYRIiIiUgyDCBERESmGQYSIiIgUwyBCREREimEQISIiIsUwiBAREZFiGESIiIhIMQwiREREpBgGESIiIlIMgwgREREpJjGDyOZ/A4fWt+17yjLgsrftexIREXVyGqULoIj9PwC6pNYd+/7VwJFfAY0eMOcAZ/0FyCgC/B7gpfHAH9s44BAREXViHSOIHN4K+LMBnRnQJ4sQcCJ8bkCtFeuOSmDdW4DPBQy8GMjoJba77MCvXwLGVCCrL6DSiODhrAamPAt0HQX8+DLw1b3AFe+JIKLWiZoRnwvQGk+sjERERAmgYwSRV88E9FL9c5VWBJJQMNFbAIMF0CWLAKDWARqDCCzhhwHoNwVI6RoMDXrA5wHenAKM+r0IGe9cCkz/D5CcLbb3vwCwHwS+uAMYfTMw4jpxrEYnypFWCDirxLrPIwLIoulAXTmgUgNXfSgCz46vgW+fFiGl6Axg7J8ASTr6exIRESWYjhFETBmA5AB8TvE84BUBIBQCmit7YEQQ0QEH1wDWAmDYNPH6wcuAXz4WNSDpRcDYOWJ72c8iuACiNmXFk4AxBTi4DjjzIbHd7xblmfKsqEVZdA2wbxWQMwj4bA4wYwlgSgf+NRXIHQL0mgQc3iaCjsYQEZwaWaq1DC5KO7xN/FwQEVGb6hhB5JYNgMUC+H2Apxbw1ImluxZw28XDZRfbfK762gmfWwSE0Lo5R7yfzy1qNewlQFJG/efoLeI9ao8A1vz67UmZIrgAgN8LDJ0G6ExA6c/1TTB+D5DVX4QQALDkA65qoHQTkD8CMGeL7f0vAPZ/L4LI/h+A//yxGSdAOnZQ0TYVZI5xTJNLw9HHqjvGj0lMvXMJcOvmlh/nrgV2LhE/J8lZQNaA+hq1si0iYGb1a9uyEhF1IB3rCqPWiJoIY8qJvY/fK2o4sgcASx4AvE4RNHZ8CYy5RdTArH0DCPhFc8reb4GupwSPdQPZ/UXNyrlPAp/fDty4TIQfjSHiM4K1LiqtWA9x2wGtSawbrKLGxBcRliKXYbKoDQrVCLU3Sd3KUKMXF+BmByNjMFQZ649V65SvDfJ7RR+hSB6HKN/xylZTAnxxFzDiWqByjxitNSPY92jHV6LTNIMIESWwjhVE2kp6EaA3AykFwCl/AF4/F5BUQPfTgMKJ4uLS/TTRN8WQIi4WoVE2Pnd9M03ByeJ9dv1P9CuJ7EQb6hCbPxw4sl30E0nOAtb9U3RuBYABF4pHY2S56YBy1PJ4+zhb8F7BZWR4kv2At0482p0UDDOGiFATEVjCtUGGBtsiwkz42GO9xzHCj9dRHx5rD4uRU7oksT75UaDHWNGJ+Ys7AdtBoOcEoKYUOO//RMhN6QpMuFsc/85lwK5vgIFTxWumNPHzI8viZ08VHFEfCIgQk5RZX4NCRNQJJWYQOfvR+vUR1wLDrxHrkRefCXOBifcCkIGXJgCZwf4BV38i+nqEXPiCuKAkZ4lOqCGDfguk9RQXuWkfAj/8QzQpXfRi/cicY5EkcazWcPx9YyEQiG7WCi29xwo1LQg6kQHJ6xTBx+sUr3mdAORgQWQRBLwOoN0qhILhJ3cIcN0XgNdVX9u14gkxumrkDUDFLhEsbl4DLP+r6IM09SVg40Jg1XP1QUSjF0HDWQVU7AQseeK9vA5gzWtA/wuBkg0ikEz+C2A/BLz7OxGYj2wHJt4D9D2vvb48EVG7Sswg0lBj1ev/vlZcFGvLxIXH2kVsT+sRvV9KQf16wcj69R7jIvbpKi4wHYlKBaiMygxDluWIYBIMLF5X9DIcilwNlo3tH3w0fI+o1xoJP4FgrZDXUX8eSjYCI2eK9fSeIlx6XcCelSJwAsCAi4Cv54l1n1P0E3ppAlD+K3Dy9fXNfD6X6G908gzRx+nVM8XPycqnRDgecZ2oWXnlTKDPueLndNl8MSJLaxJlOmppBLRJR29TqWP/70ZE1AoMIk257J+i3wdHrLQ/Saofdt1eosJPMJyEJh72OuuDiN4MOCvFuscBBHyinGqNCEGAaKYJ1WR5nUDP04FL3gCqi4E3zhOjsYypIuCkdhf7afTivQCg7BcRWADRwVoC4K4RQ9Q3vQdU7m7591Prjw4uOlNwmSSGvusi1hvbrk2qb6YMPbQm/v8gohPCIHIsbJtPHMcKPz5nfdPMiOuAr+8HRv8R+PkDYNjV4tiBFwNLHwJO+T2w9k3R5wQINs0E11O6AgOmAj+9Aoy7Q9SkhPoeRda6mNJETVxWXzHqxusSAQgA+l0AOCvE+3qdwWaryGWDbSH+4AgyV3Vbn7iI0NIgoBwz4ETu2wkDzu7lollt1I1Kl4Qo7jGIEB1PZt/6+WL6ngck5wDFq0Q/oN5ni+2nzgY2vQ9s+UTUgFQXi+0GqxhlFTJqJvD2pcDoW0Stiyai5iQUWEbNBL66Dxh5I7DtM7EMXZTPfKD55ZblYwSWOrH01EU/vKF1R/1Qea+jfsi8x1G/n/iQ+o7MbdqXWTq61ib0PBxoTMEQE9onuZH9GwSjUI1QrEOObX99zVlL/fqV+PkypIgpAYrOqG9a2/CuaPpTqu8YUQxIsizLx99NGXa7HVarFTabDRaLReniEDWtZBPgsonajOWPA93GHP+v4dB/PUkSc+TUHQEsuWJb+U6geLWYvbf7mNiWvTUCAVFTFBViIsJLKMg0FWI8xwg8iPWvpIiQE1kjY7CIuYT05uC6Nbg018/eHF63imVTtaarF4iReKf8QTSrrX1T/Hz0mwLkniT2kWXRzKbR1/dBA4DF94hmuq6nAFs+Fs1qF78sXvvbUOCm79u32ZKoFVpy/WaNCFFbMFiAbf8VnU77ni+aao4n8q9ytaY+hADiRooZRW1fzraiUkUPa28roVqcqHDiCIaXyCDT1DZHg2AUsV94Hp7IWpwjJ1Zetb4+wCRlAOPuFJMVuuyiKS4QAP71WzGRYY+xwEczxciqnEHA53cAVXvEEO26chFQzrgfcNvEUO6BU8U0AgsiOsH7vaLM+1aJqQHSe0aXJxCoHwJO1EEwiBC1hdTuYpgtnRhJCtZUmNr+vQP+poONu1bUXIRmaY6csfmo7TUiKAGi303dEfGo3AW4gredcNlErUnFTjEPz6k3ie0nXy/6FmlNwIEfgRuXi+/88U31TXMum+gjVPYLsPU/QO5gsV2WRUfo5Y+LoeJv/xa46gNRa/bzh8C3T4lAlNYDOO9p9nGjDoNBhIgSg0odbFoxn/h7Bfz1ocRlF+HBUQ7kDRWvh4KIo7z+1hKA6J/ic4uOrLlD6mvFUnuI/UPH7v9BNM3tWy1mbgaCo7cMwHlPiRq02jKg+AfAmCY6UM9cIZoGP54FbF4EDL0SqNwLOCrEdlO6+O4duRMwdUoMIkRELaVSiyHYoXtLNRQKIpY8cU8qd43oi7L1v8Dgy0QzXKhDMyAmtOs3JXisHTjnCTFqatXfge+eBc57Urxnavf6ez9JkuiHUvaLqDUxpYntPcaJWwkAotblwxvqP0etE4HElB4MJxkRz9OBpPTo56Z09kehmGMQISJqa93HiJqQpAzRZPfWhWJ7lxFAv9+Idb0ZWHStqCU5vFVMbgfUhxhAzPz83Ajg9Hujt4f3s4haEk/EkCX7ASA5U6yrdeIO43Xloo+M3yNuHVBT0vzvojPX16gkZQBJWWIm6eTs+qU5R6zrklnjQi3GIEJE1NZOnVW/ftKlwKBLxHrkRfrSt0QfEl0SsOgaILO32D7s6vqaFl2SGH216xvR9yOymScUTHKHivWfXgHMueJ+Vld/IvaJvJ+VxyGaaRwVosnIUSnW68ojtlcGXwuuy37AUyMe1fuO/721pqNDSnLO0dt4DyWKwOG7RERKeOcycQfw8h3iAj71Hy07vqYsWCNiFJ1YN70vRtwMmHr0aJrWCATE+9VVRAeY2sOif0ptWcT64foOvM1lTAvWpmQD5jzAmg9Y8kUNTmjdwN/7HVVLrt8MIkRESnDZxf2H9BZxI8yO3qThrgXqDotQUlPadGCpO1x/O4Pj0VuC4SQUUrqIR2jdkqfM/bDouDiPCBFRvDNYRJ+RzkKfLB5phcfeLxAQs86GQkpNGWA/KB62g4DtgOjn4rKJkUlH7MCRrU2/nyk9Ipg0DCv5orlKrW3b70ptikGEiIjaj0oV7PSaIZqmmuKuDYaTA/UhxX4guAw+99bVNxuVbmr8fSSVaAIKBZOUruK2DZn9gMw+IjyRohhEiIgo/uiTRVDI7NP467IMOKuaDin2A4D9UPRIoYNrjn6flK4ilGT1rV9m9G77WYOpSQwiRETU8UhScFhxmpgyvzGBgOhga9tfH1Kq9orh0oe3iv4q1cXisePLyDcXASWrX7D2pG8woPSJzay/jSnbIj6/o/cdagZ2ViUiosTkqBSB5EgwmBzZLpaO8iYOaBBQcgaJe/6kdm/bwFCxS4yq+v23HfZOy3E1ambBggV44oknUFpaisGDB+O5557DyJEjj38gGESIiEgBdeXAkW3BcLLt+AHFmAbkDxOhJH84kDesflK5lpJl4O1LgJE3Ar3Pav13UFjcjJp57733MGfOHLz44osYNWoUnnnmGUyePBnbt29HVlZWLD+aiIiodZIygKTTxN2PI0UGlMNbxdT8pZvFKKCdS8QjxNo1OpzkDm5ex9hfF4t5ZXqdeUJfwVdVBW9xMfx2OyStFsahQ6HSx+d0/TGtERk1ahROPvlk/P3vfwcABAIBFBQU4Oabb8bdd9993ONZI0JERHHN5wbKfgYOrgs+1or5YdDg0iqpRHNOKJzkDROjhiKHFntdwD/GiplxLXnwlpbCs3cfkk4Z1eJiVS18D1Xvvovk0yfCV3YYzg0b0OOjD8NhRA4EAEmCFKM+KHFRI+LxeLB27VrMnTs3vE2lUmHSpElYvXp1o8e43W643e7wc7vdHqviERERnTiNvr7WI8RlF7UlB9cGH+vFKJ7DW8Rj/b+CxxqA/BHAFQvFvYdW/Q3oc66YFh+A99AhVC9adNwgIns88BQXQ9u1Kzy7dsHQrx/UKSkw9O2DrFtuAQDsnjoVri1bYBo6FHIggKp33kX1e+8h/foZsF5wQSzOTLPFLIiUl5fD7/cjOzs7ant2dja2bdvW6DGPPfYYHnzwwVgViYiIKPYMFnEX5B7j6rfVlIoak0Pr6gOKywaYUkUIqS4GNi4EZq4QISZ/GNQpKfDbbAAAWZZR/f4i1Hz9NdRWKzJvng1d9+7wFBfjwKzZMJ08Ar7KKjjXrUOvFcuhTrHCX22DLMvw7NoFf3kFtDki4By46SbULlsOAHDv3tPeZ+cocTV8d+7cuZgzZ074ud1uR0FBgYIlIiIiagPmHKDvueIBiE6plbsBe/BOyElZwO/eFv1IDm0UQcRqhb+6GgBQ89XXqFm6BF2e+xucGzdi/+zZKPz0U1S++RbSb7wR1innw1dejp1nngVZlqG2WlH344/Ye8ml8JWXI2ncWGhzcwEAfpsNktEI8+mnI/2GGxQ4GdFiFkQyMjKgVqtRVlYWtb2srAw5OTmNHqPX66GP0840REQn6vBTTyN53FiYRrR8avdD99wL17atkNQaqNNSkX3HHdAXFcWglNQUWZYhu92Q3W4EXG7IHjdklwsBtwey24WAywU5tO52Q3a5xbrTiYDDGVzWQQ49dzgQcLuR8XsvzBMnimHBgJizBBBBJNhFoW71KqRc/FuojEYknXIK1GYLPHv2wFtWCmtBFwCAJiMDmowMyA4H1Ckp0GZno8e/F0H2elF8w42w/fczWM8/D11ffx2SWg1JGx9T38csiOh0OgwfPhxLly7FhRdeCEB0Vl26dClmz54dq48lIopbfrsNgbq64+4XcDpFR0KdDpBlSGo1fKWlyLn3XpiGD0fV+++j7NFH0fW119qh1MqQZRnw+SB7PJC9XgQ8HiC4lL1eyB4vZK8nuPSG9zvustFt4n0CHrcIEq5gkHC7EXC7wttkjycm3/XAH26CcdgwmM84HekzZgDdTgUASBoN4BM3CJR0OsguZ/gYyaCH7HZD36MQzl9+gXHIEAScTvgrKuC32USzTrA2RdJqkfvAPBRfNwPJE8ZDnRxf09rHtGlmzpw5mD59OkaMGIGRI0fimWeeQV1dHa699tpYfixRu/BXV6Nm2TKkBIN2S3j274f9s88h+3xQm5ORPGECdN26tX0hKa5IGi1kvx8AULtyJY78/e+QNFrounZF7oMPQNLpUPHKK7Av/hKarCwEHA4kjx2L9BnXAVoNZJ84Vl/UC1WV7x7zs2S/H7LPB9nrA3xesR56eL3iIu+N2O71QfYFt4efR2zzeiO2RW5vYlvke/u8QNR+DbY1ERLimkoFyWCASq+HpNeHl5LBAJVOB8lgqN9uMkJlMkFlNEFlDK6bjJCMEduD+/gqK6FJS6v/HLUass8H8+mn48gzzyL5jDPg3b8f3oOHoCsqQuq0q7B/xgx49u2Dt3g/tAUFCNTVQZObC5XVAlmWIUkSdN27I2PWLLh//RWmYcOUO2+NiGkQueyyy3DkyBHcf//9KC0txZAhQ7B48eKjOrASdUT+2lrYPvyoVUHEu38/apYsQcasm+ArK8PeK69Cjw8/gJbz68RU+K/s4EW6qXXZ5wMaXfcDfl+j67I/uJ/Xd9S6+YwzYOjbF5JGA9nrQ8DlQum8B9D9vYVQZ2Sg9P55qP7gA5jPPBPVH36Ewk8+hqTV4tBdd4nPByCpNah65x3UfLkYjp/WIOOPN4vv5PNhx7jxouxeLxAMAIjfSbNbR6WCpNNB0mrFI7R+1DL6dZVOBwSXje6v1UUfo9dDZTBA0umhMgTDhd4Qva7XtVuzhiYjA/6aGiSdcgq8lxzCwVtuhcpiRsGCv0Ol00GVlYXu770H967d0HXvhuLp10CdkgJJklD05ZdR75Vy8dR2KXNLxbyz6uzZs9kUQ52SFPxLJaTmm29Qt3IlNDm5SLvqSqhMJsg+HypefhnODRthnjwZrq1bkXPvPUCwnd88caI49uslcP3yi2JBRPb7xV/qXm/0xTh0cfb6Ii66R78m+4PPQxfkiHWxrz/4vt7gRb3BepPHB/drdP0YQcIf/MzgttBrCNZGtDdtTk4wiKgBvw/eAwegzcuDJlPMvpl8+kTULP4S+j59oS8qCl/kDIMHw19VBUBU0xsGDoBx0CB4ivfDX14e3u632Zr33TQaSBqNuOgG16HVQNLUP4/aFrGfpNUEj4/eJmm10du1wfeI2qY9eru2/vOaDhURS7U6Nv84ca77u++E11MuntpomKh8658AgPIFC2AcMiT8cwVA3F/nwDqg/5SYl7W14mrUDFE8kmU5+i/igB9qi0X8dRsMIjVLlqDyzbeQfd+fUbdiBQ7+6XYUvPA8Kt98E96yMuQ/9X+o/uBD2D/9FDn33gNJq0HAZodryxZ49h+Ae/duGPqJjmq2T/8D957d9RdjXzAEhKqxIy/UoefeBhdkn6/R15rat9P99dxSGo3ovKdWA1ptxLoGkjr4mlYDBNehUYuLqVotwoUmuF/kesTxuu7dw58j+3xQmS3w26rDHx+orYVkMkKTlQXP/v3h6nRfSUk4lEgaNfQ9eyLplFOg79MHey6aCvPkydCkpaFo6RIE6urCF300CBvhAJAAN1DrlPxeMWtrZl8xdXztYTFzqyQBKjUs550L78GDsE45H5q8vPrjDm0APpkN3LhMqZI3C4NIG5D9fuz//R/gt9sAfwC6wh7IufdeqK1WpYsWF+RAINgO7BVt1aG24Mg26lCHs4jX69ulvdFtzt7ofRH5esP3DXZEC7dbN+NiHX4t4nkkXWEhen7+GRDRkcz++edInzkTht69oe/VC5VvvwN/bR3qvvsO2XPnQpWUhNQrr0D5888DEH/Feg4eRNW778K9azeSx40Lj/Evf+kleHbubN9/pMYEL8b1f8VGX1zFBU4NhP7qjdy3BRfrFl/gG6xLmuDxmmB40Ggj1iPKHbEeVU61ut0u0JJa9PPQZGVCk5mFildfhXHwYFS8+hpy7r8fui750BcVoeTeP0Obn4e671YhaexYcXBE8NWkpiLl0ktQ8fIryL7rzvDPDnVSai3QbTSweZGYyXXoNECtAXYvBwrHQ1dQAF3DqS4CAeDzO8Q8Jf+9BTj15vBonHjDIBJh928uQI9/LxI91VtAUqvhWLsWfb5fDajVKH3oYVS88gqy/vSnGJX02ORAQHT4Cj3cbtHT3BPsCOZxh18LuN0R2+tfC3g8oqd4Y9s9XjGELfQe4feL2B4RBJSqDo+VcJt9xIUh4PZAZTKK7ZIEdXISZKcDkFRiKuXQsaFzoVaLDooPPww5EEDx1dNRs2wZzBMmIO+xR2H75NPoi7Sm4QVYE7xoh6q/Iy7I2ogLb9RzbXj9qHAR8Zdz+LlK1b4nNgGYTj4ZqqQkSJKELgv+jupFi1C7YiVyH5gH45AhAIC8v86H46c1kFSS6OioE1MaZMycCXVqavi90q+7Dt6SEiW+BinBYAVOvj56W+F4UZsZCtIrnwKyB4qb5W1aKOYkufSfonnm3d8BN/0AaFp2fWsPDCIRJIMBfrtdjMOWZbh+2QJ/dTVMw4ZCZTIBABzr16PyzbegNpuhTk+DeeJEGAcPhtpigb+uDprUVBj69YVzw8ao95ZlWQw9c7kQcLogOx3BdacYKuZ0IeB0hNdllzO4zVm/7nJCDm4Lr7tc9Rd/txuBYA1BXJOk+g5nGg2g09a3R2sj25Sj16ENVTc3sa82Yl9Nw30jqqcbXqxDVddqdX3bdaMX64gLPoJBJBgsTMOGwf7FYhiHDoVn507IXh/U6elInjgRFa++iuy774btk08hB29hIGm09YFGpULWHbej5P55SB43DsZBg2AcNEiZfxuKqcipulVGI9KuvvqofWoWL4bfZoO/pha2Tz9Bt9dfBwDoe/aM2k9lMEDfo0dsC0zxL7I279A6oNsYMWPr8seBaz4TYSStp2jekQNNv4+CGEQiqK1W+G02aDIyUHLffZBdbmgLuuDwX+ej2z//CajVOHTX3Sh48UWoDHrsu3o6DL17iyBitaLipZche72oW7UK+f/3JADAtW0b9l15FQIul2I1A6Ge4JJOB0mnE73HIx96veg1rtNBpdM3sV0HKeq1iPfR64M9z0OPYE/1UAjQNQgMnaTTmaTRQJ2SAgBIvepKlD3yCPZNmwZJpUbeE09AUqmQevnvIPu8KH3gQZjPOB2aYBW6OiUFyePqp382Dh6MlKkXwXuoBLou+Up8HYoTxuEj4NywAZqcHHR/+2028VLz/eY5ABKwbD4wfDpgDf4uWf9PoOdEQGsA1r4JrHkVMOcBk+bVT6KmoIQMIn67Hb7ycvirquCvroZh4EBos7PD0+l6DhyEa9Nm9PjkY0jBuxNWLVoE44ABMA0dCn2h+CvEMuX88H0A1BYL1Gmp0PfqBfeOHXBu2gxDv37QFRaGh/+FqdVQGY2QjAaoDEYxVMxoFOPLQ+sGg3jdaKpfNxihMjbYNzRWvUHYCIUHaLXsoBYjklaL7m+Lm1ep9HrkPvzwUfvIPh+Sx41H6hVXwP7Z59AXFgIAtNlZyJw9K2rftOnTY19oinva7CxoJ5+ldDGoIzIGm+6GTgPSgzVojkrgu2eBGV8BP38IrHsLuPIDwFsHvHsFcP0SQGdSrszopEHEb7fDU7wf3kMH4T14CN6DB+E9FFwePIhAbW3U/vlPPQXtueeEa0QkrRaanOzwBVzfrx9qly+HafhwBByO8HEqna4+iKRYYejTR1StDx6MPVMvRnJwbv+iZd9A9nihMhqgMhgYDhJJIICKl1+Gr6Ic+p5FyH3sUaVLRESdXXb/+vXvngFGzQSSs4CV/wf89nUx8gaZosakYgeQO1ipkgLoBEHEV1EBx7p1cG3aDNev2+H+dQd8zejApTKboU5NhTo1BVKSSIOiRsQG04gRcO/YiYDHA5VOB+++fdDm5MIwYADcO3fC8dNP0GRloebrJTCNHCneLxhiANGjPePGG1D2+OPo8vTT0ER0MKPEojIakcfwQURKGXcHoDGIdWc1kNFLrNdVAKU/i/4jCuuQQcRTXAz7F4thX7wY7q1bG91HnZEBbX4edPn50ObnQ5uXJ5b5+dDm5oY7n0YdkyLChNpiQerll2P/dTOg79sXznXrUPDSP6DS69Hl78+h4uVXIGk1sE6dCu+hQwCA1EsvjWrLTbnsMlgvuCA2J4CIiKg59Ob69V5nAt8+DfQYD3x1LzDmj6Izq8IkWY7fmYzsdjusVitsNhvMZjNqly9H+fMvwLVpU/1OkgR9URGMQ4bA0L8f9MF5HNQWS4s/z71zJwIuN4wDBwAAvIcOwW+zQd+rV3ikhK+8HCqTCQGnE4funovU310G8xlntMn3JSIiipmAH1j7BnBgDVB0BjDw4uhRN20o8vptOc71uEMEkapDh1Dzl7+gdslS8YJajaRRo2A+52yYJ01q16YPx/r1qPjHS5AMBpjPOAOW889jfw8iIqIInS6IbDx/CrQ7dkDSapF69TSkX3cdNOnpShePiIio4zm8Fdi7EvC6RPNMDLQkiHSIPiLubdtgyMxEwfMLwrMPEhERUSuUbBTTvxecErMg0hIdZg7n/CefYAghIiI6UTnBmZuPbFO2HEEdIojo+vRB0ujRSheDiIio47MGb5DnqgbctcfctT10iCCSFJyrg4iIiE6QwQLog9NN2A8qWxZ0kCCizc9TughERESdhylNLJ1VypYDHSSIaLKzlS4CERFR5xGayIxNM82jTlZ+5jciIqJOQxeccdVTo2w50EGCiKTVKl0EIiKiziNUI+KpU7Yc6CBBBJoOMd0JERFRx6DRi6XPpWw50EGCCGtEiIiI2lL83JqkYwQRtVrpIhAREXU+cXCXlw4RRGJ1d0AiIqKEFEfX1Y4RRIiIiKhTYhAhIiJKOKwRaREpjqqQiIiIqO10iCBCREREMcDOqs3EGhEiIqK2E0fX1Y4RRIiIiCgGWCPSPHGU3IiIiDq++LmudowgQkRERG2PfUSaR221Kl0EIiIiioGOEUQsFqWLQERE1HnEUZeHDhFEiIiIKBbYNENERETtjjUiREREpDR2ViUiIqJExiBCRESUaMKdVVkjQkRERAmMQYSIiCjhsLMqERERKY2dVYmIiKjdcUIzIiIiUh5rRIiIiCiBMYgQERElnGDTDPuIEBERUSJjECEiIko07KxKREREymPTDBERESUwBhEiIqKEw86qRERERAwiRERECYd33yUiIiJiECEiIiIFMYgQERElHHZWJSIiImIQISIiSjjhiVVZI0JEREQJLGZBpHv37pAkKeoxf/78WH0cERERNVv89BHRxPLNH3roIdxwww3h52azOZYfR0RERB1MTIOI2WxGTk5OLD+CiIiIOrCY9hGZP38+0tPTMXToUDzxxBPw+XzH3N/tdsNut0c9iIiIqI3F0cyqMasR+eMf/4hhw4YhLS0Nq1atwty5c1FSUoKnnnqqyWMee+wxPPjgg7EqEhEREcUZSZab31Pl7rvvxl//+tdj7rN161b07dv3qO2vvfYaZs6cidraWuj1+kaPdbvdcLvd4ed2ux0FBQWw2WywWCzNLSYREREdy6d/BNa9CUz8MzD+jjZ/e7vdDqvV2qzrd4tqRP70pz/hmmuuOeY+hYWFjW4fNWoUfD4f9u7diz59+jS6j16vbzKkEBERUefToiCSmZmJzMzMVn3Qhg0boFKpkJWV1arjiYiIqK110j4iq1evxg8//ICJEyfCbDZj9erVuO2223DVVVchNTU1Fh9JREREzRXurKq8mAQRvV6PhQsX4oEHHoDb7UaPHj1w2223Yc6cObH4OCIiIuqgYhJEhg0bhu+//z4Wb01EREQnLH5mVuW9ZoiIiEgxDCJEREQJizUiRERE1N7iqLMqgwgREVGiYh8RIiIian+sESEiIiJiECEiIko4cXT3XQYRIiIiUgyDCBERUaJiZ1UiIiJqf+ysSkRERIpjjQgRERG1N05oRkRERMQgQkRElLjYWZWIiIjaH5tmiIiISHGsESEiIqL2xs6qREREpBiVBtAYxVJhkizHQU+VJtjtdlitVthsNlgsFqWLQ0RERM3Qkus3a0SIiIhIMQwiREREpBgGESIiIlIMgwgREREphkGEiIiIFMMgQkRERIphECEiIiLFMIgQERGRYhhEiIiISDEMIkRERKQYBhEiIiJSDIMIERERKYZBhIiIiBTDIEJERESK0ShdgGORZRmAuJ0wERERdQyh63boOn4scR1EKioqAAAFBQUKl4SIiIhaqqamBlar9Zj7xHUQSUtLAwAUFxcf94tQ27Lb7SgoKMD+/fthsViULk7C4HlXDs+9MnjelRHr8y7LMmpqapCXl3fcfeM6iKhUoguL1WrlD6hCLBYLz70CeN6Vw3OvDJ53ZcTyvDe3AoGdVYmIiEgxDCJERESkmLgOInq9HvPmzYNer1e6KAmH514ZPO/K4blXBs+7MuLpvEtyc8bWEBEREcVAXNeIEBERUefGIEJERESKYRAhIiIixTCIEBERkWIUDyILFixA9+7dYTAYMGrUKPz444/H3H/RokXo27cvDAYDBg0ahM8//7ydStr5tOTcv/zyyxg7dixSU1ORmpqKSZMmHfffihrX0p/5kIULF0KSJFx44YWxLWAn1dLzXl1djVmzZiE3Nxd6vR69e/fm75tWaum5f+aZZ9CnTx8YjUYUFBTgtttug8vlaqfSdg4rVqzAlClTkJeXB0mS8PHHHx/3mGXLlmHYsGHQ6/UoKirCG2+8EfNyAgBkBS1cuFDW6XTya6+9Jv/yyy/yDTfcIKekpMhlZWWN7v/dd9/JarVafvzxx+UtW7bIf/7zn2WtVitv3ry5nUve8bX03F9xxRXyggUL5PXr18tbt26Vr7nmGtlqtcoHDhxo55J3bC097yF79uyR8/Pz5bFjx8oXXHBB+xS2E2npeXe73fKIESPkc889V/7222/lPXv2yMuWLZM3bNjQziXv+Fp67t9++21Zr9fLb7/9trxnzx75yy+/lHNzc+XbbrutnUvesX3++efyvffeK3/44YcyAPmjjz465v67d++WTSaTPGfOHHnLli3yc889J6vVannx4sUxL6uiQWTkyJHyrFmzws/9fr+cl5cnP/bYY43uf+mll8rnnXde1LZRo0bJM2fOjGk5O6OWnvuGfD6fbDab5TfffDNWReyUWnPefT6fPHr0aPmVV16Rp0+fziDSCi097y+88IJcWFgoezye9ipip9XScz9r1iz59NNPj9o2Z84cecyYMTEtZ2fWnCBy5513ygMGDIjadtlll8mTJ0+OYckExZpmPB4P1q5di0mTJoW3qVQqTJo0CatXr270mNWrV0ftDwCTJ09ucn9qXGvOfUMOhwNerzd8Y0I6vtae94ceeghZWVmYMWNGexSz02nNef/0009x6qmnYtasWcjOzsbAgQPx6KOPwu/3t1exO4XWnPvRo0dj7dq14eab3bt34/PPP8e5557bLmVOVEpeXxW76V15eTn8fj+ys7OjtmdnZ2Pbtm2NHlNaWtro/qWlpTErZ2fUmnPf0F133YW8vLyjfnCpaa05799++y1effVVbNiwoR1K2Dm15rzv3r0b//vf/3DllVfi888/x86dO3HTTTfB6/Vi3rx57VHsTqE15/6KK65AeXk5TjvtNMiyDJ/Ph9///ve455572qPICaup66vdbofT6YTRaIzZZyveWZU6nvnz52PhwoX46KOPYDAYlC5Op1VTU4Np06bh5ZdfRkZGhtLFSSiBQABZWVl46aWXMHz4cFx22WW499578eKLLypdtE5v2bJlePTRR/H8889j3bp1+PDDD/HZZ5/h4YcfVrpoFCOK1YhkZGRArVajrKwsantZWRlycnIaPSYnJ6dF+1PjWnPuQ5588knMnz8fS5YswUknnRTLYnY6LT3vu3btwt69ezFlypTwtkAgAADQaDTYvn07evbsGdtCdwKt+XnPzc2FVquFWq0Ob+vXrx9KS0vh8Xig0+liWubOojXn/r777sO0adNw/fXXAwAGDRqEuro63Hjjjbj33nuhUvHv51ho6vpqsVhiWhsCKFgjotPpMHz4cCxdujS8LRAIYOnSpTj11FMbPebUU0+N2h8Avv766yb3p8a15twDwOOPP46HH34YixcvxogRI9qjqJ1KS8973759sXnzZmzYsCH8+M1vfoOJEydiw4YNKCgoaM/id1it+XkfM2YMdu7cGQ5+APDrr78iNzeXIaQFWnPuHQ7HUWEjFAhl3hotZhS9vsa8O+wxLFy4UNbr9fIbb7whb9myRb7xxhvllJQUubS0VJZlWZ42bZp89913h/f/7rvvZI1GIz/55JPy1q1b5Xnz5nH4biu19NzPnz9f1ul08r///W+5pKQk/KipqVHqK3RILT3vDXHUTOu09LwXFxfLZrNZnj17trx9+3b5v//9r5yVlSU/8sgjSn2FDqul537evHmy2WyW3333XXn37t3yV199Jffs2VO+9NJLlfoKHVJNTY28fv16ef369TIA+amnnpLXr18v79u3T5ZlWb777rvladOmhfcPDd+944475K1bt8oLFixIjOG7sizLzz33nNy1a1dZp9PJI0eOlL///vvwa+PHj5enT58etf/7778v9+7dW9bpdPKAAQPkzz77rJ1L3Hm05Nx369ZNBnDUY968ee1f8A6upT/zkRhEWq+l533VqlXyqFGjZL1eLxcWFsp/+ctfZJ/P186l7hxacu69Xq/8wAMPyD179pQNBoNcUFAg33TTTXJVVVX7F7wD++abbxr9nR0619OnT5fHjx9/1DFDhgyRdTqdXFhYKL/++uvtUlZJllnXRURERMpgrx8iIiJSDIMIERERKYZBhIiIiBTDIEJERESKYRAhIiIixTCIEBERkWIYRIiIiEgxDCJERESkGAYRIiIiUgyDCBERESmGQYSIiIgUwyBCRK0yYcIE3Hzzzbj11luRmpqK7OxsvPzyy6irq8O1114Ls9mMoqIifPHFFwAAv9+PGTNmoEePHjAajejTpw+effbZqPdctmwZRo4ciaSkJKSkpGDMmDHYt28fAGDjxo2YOHEizGYzLBYLhg8fjjVr1rT79yaitsUgQkSt9uabbyIjIwM//vgjbr75ZvzhD3/AJZdcgtGjR2PdunU466yzMG3aNDgcDgQCAXTp0gWLFi3Cli1bcP/99+Oee+7B+++/DwDw+Xy48MILMX78eGzatAmrV6/GjTfeCEmSAABXXnklunTpgp9++glr167F3XffDa1Wq+TXJ6I2wLvvElGrTJgwAX6/HytXrgQgajysViumTp2Kt956CwBQWlqK3NxcrF69GqeccspR7zF79myUlpbi3//+NyorK5Geno5ly5Zh/PjxR+1rsVjw3HPPYfr06bH9YkTUrlgjQkStdtJJJ4XX1Wo10tPTMWjQoPC27OxsAMDhw4cBAAsWLMDw4cORmZmJ5ORkvPTSSyguLgYApKWl4ZprrsHkyZMxZcoUPPvssygpKQm/15w5c3D99ddj0qRJmD9/Pnbt2tUeX5GIYoxBhIharWHTiCRJUdtCzSqBQAALFy7E7bffjhkzZuCrr77Chg0bcO2118Lj8YT3f/3117F69WqMHj0a7733Hnr37o3vv/8eAPDAAw/gl19+wXnnnYf//e9/6N+/Pz766KN2+JZEFEsMIkTULr777juMHj0aN910E4YOHYqioqJGazWGDh2KuXPnYtWqVRg4cCDeeeed8Gu9e/fGbbfdhq+++gpTp07F66+/3p5fgYhigEGEiNpFr169sGbNGnz55Zf49ddfcd999+Gnn34Kv75nzx7MnTsXq1evxr59+/DVV19hx44d6NevH5xOJ2bPno1ly5Zh3759+O677/DTTz+hX79+Cn4jImoLGqULQESJYebMmVi/fj0uu+wySJKEyy+/HDfddFN4eK/JZMK2bdvw5ptvoqKiArm5uZg1axZmzpwJn8+HiooKXH311SgrK0NGRgamTp2KBx98UOFvRUQniqNmiIiISDFsmiEiIiLFMIgQERGRYhhEiIiISDEMIkRERKQYBhEiIiJSDIMIERERKYZBhIiIiBTDIEJERESKYRAhIiIixTCIEBERkWIYRIiIiEgx/w/QPjtkC79SGgAAAABJRU5ErkJggg==", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAkgAAAG0CAYAAADJpthQAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjguMywgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy/H5lhTAAAACXBIWXMAAA9hAAAPYQGoP6dpAABcYElEQVR4nO3de1hU1f4/8PcwwIDIDIICkoCkpqKiImajpqUkGpmlXSxUStPwoCmelPiVlmbeupiXUtOT2jl6TC09HUkNMTATFUEU1LxFisKAN2YA5eLM+v3hl32cAZX7zOD79Tzz1Ky1Zu3PmkHmw9prry0TQggQERERkcTG3AEQERERWRomSEREREQmmCARERERmWCCRERERGSCCRIRERGRCSZIRERERCaYIBERERGZsDV3ANbIYDAgOzsbzs7OkMlk5g6HiIiIqkAIgYKCAnh5ecHG5v5zREyQaiA7Oxve3t7mDoOIiIhqICsrC61atbpvGyZINeDs7AzgzhusVCrNHA0RERFVhU6ng7e3t/Q9fj9MkGqg/LSaUqlkgkRERGRlqrI8hou0iYiIiEwwQSIiIiIywQSJiIiIyATXIBGR2ej1epSVlZk7DCJqROzt7R94CX9VMEEiogYnhIBGo0F+fr65QyGiRsbGxgZ+fn6wt7evVT9MkIiowZUnR+7u7mjSpAk3XCWiOlG+kXNOTg58fHxq9buFCRIRNSi9Xi8lR25ubuYOh4gamRYtWiA7Oxu3b9+GnZ1djfvhIm0ialDla46aNGli5kiIqDEqP7Wm1+tr1Q8TJCIyC55WI6L6UFe/W5ggEREREZlggkRERERkggkSERERkQkmSERE1XTp0iW89NJLmDlzprlDeejxs7AsjenzYIJERFRNUVFRaNeuHbZs2VKn/X711Vdo3bo1HBwc0KtXLxw+fLhO+2+M6uOz2LdvH4YOHQovLy/IZDJs3769zvpu7Orj85g/fz569uwJZ2dnuLu744UXXsDp06frrP97YYJERFYr82oRFu76A5P/fRQLd/2BzKtF9X5MrVaLhIQE9O3bF15eXnXW7/fff49p06bhww8/RGpqKrp27YqQkBDk5eXV2THqzbXzQHZaxce18/V62Pr6LIqKitC1a1d89dVXddbnw6C+Po/ExERERkbi4MGDiIuLQ1lZGQYNGoSiovr9986NIonIKm0+koX3fjgOmUwGIQRkMhlWJZ7HwhEBeDnIu96Ou3fvXjz55JNITExEnz596qzfL774AuPHj8ebb74JAFi5ciViY2Px7bff4r333quz49S5a+eBZYH3rp+cCri1qZdD19dnMWTIEAwZMqTO+ntY1NfnsWvXLqPn69atg7u7O1JSUtCvX786O44pziARkdXJvFqE9344DoMA9AZh9N/oH47jr3qcSfrtt9/wxBNPYMeOHXjxxRcr1M+bNw9Nmza97+PixYtGryktLUVKSgqCg4OlMhsbGwQHByMpKanexlInSgpqV18L9fFZUM011Oeh1WoBAK6urnU+hrtxBomIrM7mI1l3NoMTokKdTCbD90eyED24Q70c+8iRI+jTpw8cHBwQGFhx5iQiIgKvvPLKffswPf1w9epV6PV6eHh4GJV7eHjgjz/+qH3QjVR9fBaNSVZBFj468BFm956NVs6t6v14DfF5GAwGTJ06FX369EHnzp1rFe+DMEEiIqtz6cYtiEqSIwAQQuDSjVv1duy//voLV65cwezZsyutd3V1rfe/bOkOfhb3t/LYSvyp/RMrjq3AJ30/qffjNcTnERkZiYyMDOzfv79W/VQFT7ERkdVp1czxnrcTkMlkaNXMsd6OnZubC3t7e7z00kuV1tfkNELz5s0hl8uRm5tb4Vienp71NhZrVx+fRWMigwwOtg6QoWFu61Pfn8ekSZOwY8cO/Prrr2jVqv5nxDiDRERW55Ugb6xKrPwKKSEEXq3HRdp2dnb4/PPPYWNT+d+XNTmNYG9vjx49eiA+Ph4vvPACgDunEuLj4zFp0qQ6ibsxqo/PojGZ23dugx6vvj4PIQQmT56Mbdu2ISEhAX5+fnUS74MwQSIiq+PX3AkLRwQg2uQqNiEEFo4IQOvmTvVy3PXr16OoqAgKhQIHDx7E1atX8dxzzxm1qelphGnTpiE8PBxBQUF4/PHH8eWXX6KoqEi6qs1iKZxrV19D9flZFBYW4ty5c9LzzMxMpKWlwdXVFT4+PrWOvTGqz88jMjISGzduxH/+8x84OztDo9EAAFQqFRwd62+2mAkSEVmll4O80bO1K74/koVLN26hVTNHvBrkXW/JUXFxMX788UesX78eEyZMQOvWrfHdd9/VWf+vvvoqrly5glmzZkGj0aBbt27YtWtXhYXbFsetzZ1L+Su7Wk3hXC+X+Ffrs7hdDBgMldfZ2AC2DhWKjxw5gqefflp6Pm3aNABAeHg41q1bV9vwG536/rexYsUKAMBTTz1lVL527Vq88cYbdXYcUzJxr5WOdE86nQ4qlQparRZKpdLc4RBZleLiYmRmZsLPzw8ODhW/nIjqzO1iIO/U/du4d6w0SSLrdb/fMdX5/uYibSIiapzuNXNU3Tb0UGKCRERERGSCCRIRERGRCSZIRGQWXP5IRPWhrn63MEEiogZlZ2cHALh586aZIyGixqi0tBQAIJfLa9UPL/MnogYll8vh4uKCvLw8AECTJk3uuSs2Ua2UlgC3HzCbUFwCGDhX0FgYDAZcuXIFTZo0ga1t7VIcJkhE1ODKb59RniQR1Qt9GVBw5f5tdLaA3K5h4qEGYWNjAx8fn1r/4cUEiYganEwmQ8uWLeHu7o6ysjJzh0ON2Q0HoKyo8jo7J6AZd8ZubOzt7e95u5PqYIJERGYjl8trvU6A6L5aPmbuCMhK8cQrERERkQmLTZAWLFgAmUyGqVOnGpUnJSVhwIABcHJyglKpRL9+/XDr1i2p/vr16wgLC4NSqYSLiwvGjRuHwsJCoz6OHz+OJ598Eg4ODvD29saiRYsaYkhERERkJSwyQUpOTsaqVasQEBBgVJ6UlITBgwdj0KBBOHz4MJKTkzFp0iSjc41hYWE4ceIE4uLisGPHDuzbtw8TJkyQ6nU6HQYNGgRfX1+kpKTg008/xUcffYRvvvmmwcZHREREls3iblZbWFiIwMBAfP3115g7dy66deuGL7/8EgDwxBNP4JlnnsHHH39c6WtPnToFf39/JCcnIygoCACwa9cuPPvss7h06RK8vLywYsUKvP/++9BoNLC3twcAvPfee9i+fTv++OOPSvstKSlBSUmJ9Fyn08Hb25s3qyUiIrIiVn2z2sjISISGhiI4ONioPC8vD4cOHYK7uzt69+4NDw8P9O/fH/v375faJCUlwcXFRUqOACA4OBg2NjY4dOiQ1KZfv35ScgQAISEhOH36NG7cuFFpTPPnz4dKpZIe3t7edTlkIiIisjAWlSBt2rQJqampmD9/foW6P//8EwDw0UcfYfz48di1axcCAwMxcOBAnD17FgCg0Wjg7u5u9DpbW1u4urpCo9FIbTw8PIzalD8vb2MqJiYGWq1WemRlZdVuoERERGTRLOYy/6ysLEyZMgVxcXFwcHCoUG8wGAAAb7/9Nt58800AQPfu3REfH49vv/220qSqrigUCigUinrrnyzXrVI9fjt7BUWlt9G1lQsebdHU3CEREVEDsJgEKSUlBXl5eQgMDJTK9Ho99u3bh+XLl+P06dMAAH9/f6PXdezYERcvXgRwZ3de0515b9++jevXr0s793p6eiI3N9eoTfnz8jZEALDt6CXM3H4ChSW3pbJnOnpg8chuaKqwmH86RERUDyzmFNvAgQORnp6OtLQ06REUFISwsDCkpaXh0UcfhZeXl5QolTtz5gx8fX0BAGq1Gvn5+UhJSZHq9+7dC4PBgF69eklt9u3bZ7R7b1xcHNq3b49mzZo1wEjJGhw4dxXTvj9mlBwBwN4/chG1Kc08QRERUYOxmD+DnZ2d0blzZ6MyJycnuLm5SeXTp0/Hhx9+iK5du6Jbt25Yv349/vjjD2zduhXAndmkwYMHY/z48Vi5ciXKysowadIkjBw5El5eXgCA119/HbNnz8a4ceMQHR2NjIwMLFmyBIsXL27YAZNFW5F4HjYyQG9yjadeAHGncvHnlUKebiMiasQsJkGqiqlTp6K4uBhRUVG4fv06unbtiri4OLRp00Zqs2HDBkyaNAkDBw6EjY0NRowYgaVLl0r1KpUKv/zyCyIjI9GjRw80b94cs2bNMtoriejoxRsVkqO7pWXlM0EiImrELG4fJGtQnX0UyDr1/GQPrhSU3LN+1egeCOnENWtERNbEqvdBIrIEw7s/Arms8jonhRz92rVo2ICIiKhBMUEiqsTb/dvgkWZNIJf9L0uS28ggAzDn+c5wtOcd6ImIGjMmSESVcHWyx/bIPnjrST+0cFbAyV6O3m3c8K+3emFEj1bmDo+IiOoZ1yDVANcgERERWR+uQSIiIiKqBSZIRERERCaYIBERERGZYIJEREREZIIJEhEREZEJJkhEREREJpggEREREZlggkRERERkggkSERERkQkmSEREREQmbM0dABER0T1dOw+UFFRep3AG3No0bDz00GCCRERElunaeWBZ4P3bTE5lkkT1gqfYiIjIMt1r5qi6bYhqgAkSERERkQkmSEREREQmmCARERERmWCCRERERGSCCRIRERGRCSZIRERkmRTOddOGqAa4DxIREVkmtzZ39jniRpFkBkyQiIjIcjEBIjPhKTYiIiIiE0yQiIiIiEwwQSIiIiIywQSJiIiIyAQTJCIiIiITFpsgLViwADKZDFOnTq1QJ4TAkCFDIJPJsH37dqO6ixcvIjQ0FE2aNIG7uzumT5+O27dvG7VJSEhAYGAgFAoF2rZti3Xr1tXfQIiIiMjqWGSClJycjFWrViEgIKDS+i+//BIymaxCuV6vR2hoKEpLS3HgwAGsX78e69atw6xZs6Q2mZmZCA0NxdNPP420tDRMnToVb731Fnbv3l1v4yHLUFBchqzrN1Fcpjd3KEREZOEsbh+kwsJChIWFYfXq1Zg7d26F+rS0NHz++ec4cuQIWrZsaVT3yy+/4OTJk9izZw88PDzQrVs3fPzxx4iOjsZHH30Ee3t7rFy5En5+fvj8888BAB07dsT+/fuxePFihISEVBpTSUkJSkpKpOc6na4OR0z1LVdXjDn/PYldGRrohUATezlGP+GLaYMeg8JWbu7wiIjIAlncDFJkZCRCQ0MRHBxcoe7mzZt4/fXX8dVXX8HT07NCfVJSErp06QIPDw+pLCQkBDqdDidOnJDamPYdEhKCpKSke8Y0f/58qFQq6eHt7V3T4VED0xWXYcSKA9h14k5yBAA3S/VY/dufiNxwFOL/yoiIiO5mUQnSpk2bkJqaivnz51daHxUVhd69e2PYsGGV1ms0GqPkCID0XKPR3LeNTqfDrVu3Ku03JiYGWq1WemRlZVVrXGQ+3x/OwuX8W9AbjBMhgwD2nMpFWla+eQIjIiKLZjGn2LKysjBlyhTExcXBwcGhQv1PP/2EvXv34ujRow0em0KhgEKhaPDjUu3tOZWLe00SyW1k2PtHHrr7NGvYoIiIyOJZzAxSSkoK8vLyEBgYCFtbW9ja2iIxMRFLly6Fra0t4uLicP78ebi4uEj1ADBixAg89dRTAABPT0/k5uYa9Vv+vPyU3L3aKJVKODo61vMoqaEZHnAK7UH1RET0cLKYGaSBAwciPT3dqOzNN99Ehw4dEB0djebNm+Ptt982qu/SpQsWL16MoUOHAgDUajU++eQT5OXlwd3dHQAQFxcHpVIJf39/qc3PP/9s1E9cXBzUanV9DY3M6OkO7ki5cAOGSvIgvUHgqfbuDR8UERFZPItJkJydndG5c2ejMicnJ7i5uUnllS3M9vHxgZ+fHwBg0KBB8Pf3x+jRo7Fo0SJoNBp88MEHiIyMlE6RRUREYPny5ZgxYwbGjh2LvXv3YvPmzYiNja3nEZI5vP64D9Yf+AtXC0uN1iHZyIDebZojyJen14iIqCKLOcVWF+RyOXbs2AG5XA61Wo1Ro0ZhzJgxmDNnjtTGz88PsbGxiIuLQ9euXfH5559jzZo197zEn6ybSxN7bI3ojb5tm6N85yx7uQ1ee9wH34zpUel+WkRERDLB65yrTafTQaVSQavVQqlUmjscqqK8gmJcLyrFIy6OcHawM3c4RETUwKrz/W0xp9iI6pu7swPcnSteIUlERGSqUZ1iIyIiIqoLTJCIiIiITDBBIiIiIjLBBImIiIjIBBMkIiIiIhNMkIiIiIhMMEEiIiIiMsEEiYiIiMgEEyQiIiIiE9xJm4iIGta180BJQeV1CmfArU3DxkNUCSZIRETUcK6dB5YF3r/N5FQmSWR2PMVGREQN514zR9VtQ1TPmCARERERmWCCRERERGSCCRIRERGRCSZIRERERCaYIBERERGZYIJEREQNR+FcN22I6hn3QSIioobj1ubOPkfcKJIsHBMkIiJqWEyAyArwFBsRERGRCSZIRERERCaYIBERERGZYIJEREREZIIJEhEREZEJJkhEREREJniZPzUaF6/dxLFL+XBSyNG7TXM42MnNHRIREVkpJkhk9W6W3sa7W45jZ3oOxP+VKR1sMffFLni+q5dZYyMiIuvEU2xk9aZtPoZdGf9LjgBAV3wbUzYdxcE/r5ktLiIisl4WmyAtWLAAMpkMU6dOBQBcv34dkydPRvv27eHo6AgfHx+888470Gq1Rq+7ePEiQkND0aRJE7i7u2P69Om4ffu2UZuEhAQEBgZCoVCgbdu2WLduXQONiuraX1eLsCtDA4OoWGcDYEXC+QaPiYiIrJ9FnmJLTk7GqlWrEBAQIJVlZ2cjOzsbn332Gfz9/XHhwgVEREQgOzsbW7duBQDo9XqEhobC09MTBw4cQE5ODsaMGQM7OzvMmzcPAJCZmYnQ0FBERERgw4YNiI+Px1tvvYWWLVsiJCTELOOlmjuadeOedXoBHPnregNGQ0REjYVMCFHJ397mU1hYiMDAQHz99deYO3cuunXrhi+//LLStlu2bMGoUaNQVFQEW1tb7Ny5E8899xyys7Ph4eEBAFi5ciWio6Nx5coV2NvbIzo6GrGxscjIyJD6GTlyJPLz87Fr164qxajT6aBSqaDVaqFUKms9Zqq5XRkaRPwr5Z71LZraI/mDZxowIiIislTV+f62uFNskZGRCA0NRXBw8APblg/Q1vbORFhSUhK6dOkiJUcAEBISAp1OhxMnTkhtTPsOCQlBUlLSPY9TUlICnU5n9CDL0O+x5mhiX/nVanIZ8EL3Rxo4IiIiagwsKkHatGkTUlNTMX/+/Ae2vXr1Kj7++GNMmDBBKtNoNEbJEQDpuUajuW8bnU6HW7duVXqs+fPnQ6VSSQ9vb+9qjYvqTxN7W8x+vhMAwEb2v3K5jQyeKkdE9Oddw4mIqPosJkHKysrClClTsGHDBjg4ONy3rU6nQ2hoKPz9/fHRRx/Ve2wxMTHQarXSIysrq96PSVX3cpA3vhv7OB73c4XC1gbNmtghXN0a/5nUB25NFeYOj4iIrJDFLNJOSUlBXl4eAgMDpTK9Xo99+/Zh+fLlKCkpgVwuR0FBAQYPHgxnZ2ds27YNdnZ2UntPT08cPnzYqN/c3Fyprvy/5WV3t1EqlXB0dKw0NoVCAYWCX7SWrN9jLdDvsRbmDoOIiBoJi5lBGjhwINLT05GWliY9goKCEBYWhrS0NMjlcuh0OgwaNAj29vb46aefKsw0qdVqpKenIy8vTyqLi4uDUqmEv7+/1CY+Pt7odXFxcVCr1fU/SCIiIrIKFjOD5OzsjM6dOxuVOTk5wc3NDZ07d5aSo5s3b+Jf//qX0WLpFi1aQC6XY9CgQfD398fo0aOxaNEiaDQafPDBB4iMjJRmgCIiIrB8+XLMmDEDY8eOxd69e7F582bExsY2+JiJiIjIMllMgvQgqampOHToEACgbdu2RnWZmZlo3bo15HI5duzYgYkTJ0KtVsPJyQnh4eGYM2eO1NbPzw+xsbGIiorCkiVL0KpVK6xZs4Z7IBEREZHE4vZBsgbcB4mIiMj6VOf7u1YzSGVlZdBoNLh58yZatGgBV1fX2nRHREREZBGqvUi7oKAAK1asQP/+/aFUKtG6dWt07NgRLVq0gK+vL8aPH4/k5OT6iJWIiIioQVQrQfriiy/QunVrrF27FsHBwdi+fTvS0tJw5swZJCUl4cMPP8Tt27cxaNAgDB48GGfPnq2vuImIiIjqTbXWIL322mv44IMP0KlTp/u2Ky4uxrp162Bvb4+xY8fWOkhLwzVIRERE1qc6399cpF0DTJCIiIisT4PcrHbAgAGYPXt2hfIbN25gwIABNe2WiIiIyOxqfBVbQkIC0tPTcfToUWzYsAFOTk4AgNLSUiQmJtZZgEREZKWunQdKCiqvUzgDbryZNFmuWl3mv2fPHrz99tt44okn8N///hetW7euo7CIiMiqXTsPLAu8f5vJqUySyGLV6l5sLVu2RGJiIrp06YKePXsiISGhjsIiIiKrdq+Zo+q2ITKTGidIMpkMwJ073W/cuBFTpkzB4MGD8fXXX9dZcERERETmUONTbKYXv33wwQfo2LEjwsPDax0UERERkTnVOEHKzMxEixYtjMpGjBiBDh064MiRI7UOjIiIiMhcapwg+fr6VlreqVOnB24kSURERGTJqp0g6XS6KrXjBopERERkraqdILm4uEgLtCsjhIBMJoNer69VYEREZMUUznXThshMqp0g/frrr9L/CyHw7LPPYs2aNXjkkUfqNDAiIrJibm3u7HPEjSLJSlU7Qerfv7/Rc7lcjieeeAKPPvponQVFRESNABMgsmK12iiSiIiIqDFigkRERERkok4SpPst2iYiIiKyNtVegzR8+HCj58XFxYiIiICTk5NR+Y8//li7yIiIiIjMpNoJkkqlMno+atSoOguGiIiIyBJUO0Fau3ZtfcRBREREZDFqfKuRrKwseHt712UsRNXyywkN1v7+F87mFcDd2QGv9fLBaz29YSvntQdERFQ7MiGEqMkLbWxs4Orqiq5du6Jbt27So7S0FEuXLsX69evrOlaLodPpoFKpoNVqeUsVM1my5ywW7zkDuQzQC6D8MoEBHd3xzeggyG144QARERmrzvd3jWeQMjMzcfToUaSlpeHo0aPYvHkzsrOzAfA+bFS/Llwrwpd7zgC4kxwBQHmWH38qDzszcvBcgJd5giMiokahxgmSr68vfH198cILL0hlSUlJCA8Px5w5c+oiNqJK/ZSWDRvZ/5Kju9nIgG2pl5kgERFRrdTpYg21Wo0lS5bgs88+q8tuiYzoisvuufeWQQDaW2UNHBERETU2NU6QSktLKy1v164dTpw4UeOAiB6kSysX3DZUvnROLpOhm7dLwwZERESNTo1PsTVt2hT+/v7o3r07unXrhu7du8PLywvLli1DcHBwXcZIZCSkkwdaqhyQV1AC/V2JkkwGyG1kGK32NWN0RETUGNR4Bmnv3r0YP3487OzssGHDBgwePBiPPfYYli1bBr1ej1mzZmHLli34448/atT/ggULIJPJMHXqVKmsuLgYkZGRcHNzQ9OmTTFixAjk5uYave7ixYsIDQ1FkyZN4O7ujunTp+P27dtGbRISEhAYGAiFQoG2bdti3bp1NYqRzENhK8fG8U/Ax7UJgP9dwaZytMM/3giCr5vTvV9MRERUBTWeQerbty/69u0rPTcYDDh9+jTS0tKQlpaGw4cPY/Xq1cjLy4Ner69W38nJyVi1ahUCAgKMyqOiohAbG4stW7ZApVJh0qRJGD58OH7//XcAgF6vR2hoKDw9PXHgwAHk5ORgzJgxsLOzw7x58wDcufouNDQUERER2LBhA+Lj4/HWW2+hZcuWCAkJqenbQQ3Mr7kT4qf1x8E/r+H8lUK0cHbA0x1aQGErN3doRETUGIhquHDhQnWai6ysLKHRaKr1moKCAtGuXTsRFxcn+vfvL6ZMmSKEECI/P1/Y2dmJLVu2SG1PnTolAIikpCQhhBA///yzsLGxMTrmihUrhFKpFCUlJUIIIWbMmCE6depkdMxXX31VhISEVDlGrVYrAAitVlutsREREZH5VOf7u1qn2Hr27Im3334bycnJ92yj1WqxevVqdO7cGT/++CM8PDyqlbBFRkYiNDS0wjqmlJQUlJWVGZV36NABPj4+SEpKAnBnm4EuXboYHTMkJAQ6nU5aOJ6UlFSh75CQEKmPypSUlECn0xk9iIiIqPGq1im2kydP4pNPPsEzzzwDBwcH9OjRA15eXnBwcMCNGzdw8uRJnDhxAoGBgVi0aBGeffbZagWzadMmpKamVpqAaTQa2Nvbw8XFxajcw8MDGo1GamOakJU/f1AbnU6HW7duwdHRscKx58+fj9mzZ1drLERERGS9qjWD5Obmhi+++AI5OTlYvnw52rVrh6tXr+Ls2bMAgLCwMKSkpCApKanayVFWVhamTJmCDRs2wMHBoVqvrW8xMTHQarXSIysry9whERERUT2q0SJtR0dHvPTSS3jppZfqLJCUlBTk5eUhMDBQKtPr9di3bx+WL1+O3bt3o7S0FPn5+UazSLm5ufD09AQAeHp64vDhw0b9ll/ldncb0yvfcnNzoVQqK509AgCFQgGFQlHrMRIREZF1sJjbng8cOBDp6enSVXBpaWkICgpCWFiY9P92dnaIj4+XXnP69GlcvHgRarUawJ2dvNPT05GXlye1iYuLg1KphL+/v9Tm7j7K25T3QURERFTjy/zrmrOzMzp37mxU5uTkBDc3N6l83LhxmDZtGlxdXaFUKjF58mSo1Wo88cQTAIBBgwbB398fo0ePxqJFi6DRaPDBBx8gMjJSmgGKiIjA8uXLMWPGDIwdOxZ79+7F5s2bERsb27ADJiIiIotlMQlSVSxevBg2NjYYMWIESkpKEBISgq+//lqql8vl2LFjByZOnAi1Wg0nJ6cKN8/18/NDbGwsoqKisGTJErRq1Qpr1qzhHkhEREQkkQkhKr+pFd2TTqeDSqWCVquFUqk0dzhERERUBdX5/raYNUhERERElqLGp9jCw8Mxbtw49OvXry7jISKixuDaeaCkoPI6hTPg1qZh4yGqphonSFqtFsHBwfD19cWbb76J8PBwPPLII3UZGxERWaNr54FlgfdvMzmVSRJZtBqfYtu+fTsuX76MiRMn4vvvv0fr1q0xZMgQbN26FWVlZXUZIxERWZN7zRxVtw2RGdVqDVKLFi0wbdo0HDt2DIcOHULbtm0xevRoeHl5ISoqStphm4iIiMia1Mki7ZycHMTFxSEuLg5yuRzPPvss0tPT4e/vj8WLF9fFIYiIiIgaTI0TpLKyMvzwww947rnn4Ovriy1btmDq1KnIzs7G+vXrsWfPHmzevNloDyIiIiIia1DjRdotW7aEwWDAa6+9hsOHD6Nbt24V2jz99NNG900jIiIisgY1TpAWL16Ml19+GQ4ODvds4+LigszMzJoegoiIiMgsanyK7dKlS9i4cWOF8m+//RYLFy6sVVBERGTFFM5104bIjGp8q5HWrVtj48aN6N27t1H5oUOHMHLkyEY9c8RbjRARPQA3iiQLVJ3v7xqfYtNoNGjZsmWF8hYtWiAnJ6em3RIRUWPABIisXI1PsXl7e+P333+vUP7777/Dy8urVkERERERmVONZ5DGjx+PqVOnoqysDAMGDAAAxMfHY8aMGfj73/9eZwESERERNbQaJ0jTp0/HtWvX8Le//Q2lpaUAAAcHB0RHRyMmJqbOAiQiIiJqaDVepF2usLAQp06dgqOjI9q1aweFQlFXsVksLtImIiKyPg2ySLtc06ZN0bNnz9p2Q0RERGQxapUgxcfHIz4+Hnl5eTAYDEZ13377ba0CIyIiIjKXGidIs2fPxpw5cxAUFISWLVtCJpPVZVxED6Q3CMSd1CA2XYObpbcR5OuKV3t6w9XJ3tyhERGRlavxGqSWLVti0aJFGD16dF3HZPG4Bsn8Sm7r8db6I/jt7FXYyACDAGQyQOlgh+/ffgIdPPm5EBGRsep8f9d4H6TS0tIKu2gTNZQ1v2Vi/7mrAO4kRwAgBFBYXIbIDamo5bUHRET0kKtxgvTWW29Vei82oobwr4MXUFkOpBfA+StFSMvKb/CYiIio8ajxGqTi4mJ888032LNnDwICAmBnZ2dU/8UXX9Q6OKJ7ySsouW99jrYY3RsoFiIianxqnCAdP34c3bp1AwBkZGQY1XHBNtU372aOuHDtJu51Is3XrUmDxkNERI1LjROkX3/9tS7jIKqWN3q3xkf/PVmhXG4jQ2cvJTp5qcwQFRERNRY1XoNEZE6j1a0xIvARAICtjQy2NndmLT2VDlj+eqA5QyMiokagVhtF/vbbb1i1ahXOnz+PrVu34pFHHsE///lP+Pn5oW/fvnUVI1EFchsZPnu5K0Y94YvY4zm4WaZHkG8zPNulJRzs5OYOj4iIrFyNE6QffvgBo0ePRlhYGI4ePYqSkjuLZrVaLebNm4eff/65zoIkqoxMJkN3n2bo7tPM3KEQEVEjU+NTbHPnzsXKlSuxevVqoyvY+vTpg9TU1DoJjoiIiMgcapwgnT59Gv369atQrlKpkJ+fX5uYiIiIiMyqxgmSp6cnzp07V6F8//79ePTRR2vU54oVKxAQEAClUgmlUgm1Wo2dO3dK9RqNBqNHj4anpyecnJwQGBiIH374waiP69evIywsDEqlEi4uLhg3bhwKCwuN2hw/fhxPPvkkHBwc4O3tjUWLFtUoXiIiImqcapwgjR8/HlOmTMGhQ4cgk8mQnZ2NDRs24N1338XEiRNr1GerVq2wYMECpKSk4MiRIxgwYACGDRuGEydOAADGjBmD06dP46effkJ6ejqGDx+OV155BUePHpX6CAsLw4kTJxAXF4cdO3Zg3759mDBhglSv0+kwaNAg+Pr6IiUlBZ9++ik++ugjfPPNNzV9K4iIiKixETVkMBjE3LlzhZOTk5DJZEImkwkHBwfxwQcf1LTLSjVr1kysWbNGCCGEk5OT+O6774zqXV1dxerVq4UQQpw8eVIAEMnJyVL9zp07hUwmE5cvXxZCCPH111+LZs2aiZKSEqlNdHS0aN++fZVj0mq1AoDQarU1HhcRERE1rOp8f9d4Bkkmk+H999/H9evXkZGRgYMHD+LKlSv4+OOP6yRx0+v12LRpE4qKiqBWqwEAvXv3xvfff4/r16/DYDBg06ZNKC4uxlNPPQUASEpKgouLC4KCgqR+goODYWNjg0OHDklt+vXrB3t7e6lNSEgITp8+jRs3blQaS0lJCXQ6ndGDiIiIGq8aX+Y/Z86cCmW7du2S/n/WrFk16jc9PR1qtRrFxcVo2rQptm3bBn9/fwDA5s2b8eqrr8LNzQ22trZo0qQJtm3bhrZt2wK4s0bJ3d3dqD9bW1u4urpCo9FIbfz8/IzaeHh4SHXNmlW8ZHz+/PmYPXt2jcZDRERE1qfGCdK2bduMnpeVlSEzMxO2trZo06ZNjROk9u3bIy0tDVqtFlu3bkV4eDgSExPh7++PmTNnIj8/H3v27EHz5s2xfft2vPLKK/jtt9/QpUuXmg7lgWJiYjBt2jTpuU6ng7e3d70dj4iIiMyrxgnS3Qujy+l0Orzxxht48cUXaxyQvb29NCPUo0cPJCcnY8mSJZgxYwaWL1+OjIwMdOrUCQDQtWtX/Pbbb/jqq6+wcuVKeHp6Ii8vz6i/27dv4/r16/D09ARw5+q73Nxcozblz8vbmFIoFFAoFDUeExEREVmXOr0Xm1KpxOzZszFz5sw669NgMKCkpAQ3b94EANjYGIcsl8thMBgAAGq1Gvn5+UhJSZHq9+7dC4PBgF69eklt9u3bh7KyMqlNXFwc2rdvX+npNSIiqoVr54HstMof186bMzKi+6rVvdgqo9VqodVqa/TamJgYDBkyBD4+PigoKMDGjRuRkJCA3bt3o0OHDmjbti3efvttfPbZZ3Bzc8P27duly/kBoGPHjhg8eDDGjx+PlStXoqysDJMmTcLIkSPh5eUFAHj99dcxe/ZsjBs3DtHR0cjIyMCSJUuwePHiOnsPiIgIdxKgZQ+4efTkVMCtTcPEQ1QNNU6Qli5davRcCIGcnBz885//xJAhQ2rUZ15eHsaMGYOcnByoVCoEBARg9+7deOaZZwAAP//8M9577z0MHToUhYWFaNu2LdavX49nn31W6mPDhg2YNGkSBg4cCBsbG4wYMcIoVpVKhV9++QWRkZHo0aMHmjdvjlmzZhntlURERHWgpKBu2hCZgUwIIWryQtMrwWxsbNCiRQsMGDAAMTExcHZ2rpMALZFOp4NKpYJWq4VSqTR3OERElik7Dfim//3bTEgEvLo1RDRE1fr+rvEMUmZmZk1fSkRERGTR6nSRNhEREVFjUOMZpLv3BXqQL774oqaHISIiImpwtdoH6ejRoygrK0P79u0BAGfOnIFcLkdg4P+uWpDJZLWPkoiIiKgB1ThBGjp0KJydnbF+/Xpp/6AbN27gzTffxJNPPom///3vdRYkERFZIUUVLtapShsiM6jxVWyPPPIIfvnlF2lX63IZGRkYNGgQsrOz6yRAS8Sr2IiIquja+Xtfyq9w5h5I1KAa5Co2nU6HK1euVCi/cuUKCgq4rwUREYEJEFmtGl/F9uKLL+LNN9/Ejz/+iEuXLuHSpUv44YcfMG7cOAwfPrwuYyQiIiJqUDWeQVq5ciXeffddvP7669J9zWxtbTFu3Dh8+umndRYgERERUUOr8RqkckVFRTh//s4NB9u0aQMnJ6c6CcyScQ0SERGR9WmQNUjlnJycEBAQUNtuiIiIiCxGrXbS/u233zBq1Cio1WpcvnwZAPDPf/4T+/fvr5PgiIiIiMyhxgnSDz/8gJCQEDg6OuLo0aMoKSkBAGi1WsybN6/OAiQiIiJqaDVOkObOnYuVK1di9erVsLOzk8r79OmD1NTUOgmOqLr+ulqEj346gUGLE/H8sv1Yve9PFJbcNndYRERkZWq8Bun06dPo169fhXKVSoX8/PzaxERUIykXriNszSGU6QX0hjvXHqRf1mLzkSxsjegNVRO7B/RARER0R41nkDw9PXHu3LkK5fv378ejjz5aq6CIqstgEIj6/hhKbxuk5AgABIA/rxRi6d6z5guOiIisTo0TpPHjx2PKlCk4dOgQZDIZsrOzsWHDBrz77ruYOHFiXcZI9EBpl/Jx8fpNGCrZtEIvgC1Hsho+KCIislo1PsX23nvvwWAwYODAgbh58yb69esHhUKBd999F5MnT67LGIke6Hph6X3rdcW3YTAI2NjIGigiIiKyZjVOkGQyGd5//31Mnz4d586dQ2FhIfz9/dG0adO6jI+oStp73vuO4DIAfs2dmBwREVGV1egUW1lZGQYOHIizZ8/C3t4e/v7+ePzxx5kckdl4uzZBcEd3yGUVkyABYEI/rosjIqKqq1GCZGdnh+PHj9d1LES18vnL3dDNxwUAIJfJYCO7M3v0dv9H8WpPb7PGRkRE1qXG92KLioqCQqHAggUL6jomi8d7sVkuIQSOXLiBw5nXobC1QUgnT3i7NjF3WEREZAEa5F5st2/fxrfffos9e/agR48eFW5S+8UXX9S0a6Iak8lk6NnaFT1bu5o7FCIismI1TpAyMjIQGBgIADhz5oxRnaySdSBERERE1qLaCdKff/4JPz8//Prrr/URDxEREZHZVXuRdrt27XDlyhXp+auvvorc3Nw6DYqIiIjInKqdIJmu6f75559RVFRUZwERERERmVuNbzVCRERE1FhVO0GSyWQVFmFzUTYRERE1JtVepC2EwBtvvAGFQgEAKC4uRkRERIXL/H/88ce6iZCIiIiogVV7Bik8PBzu7u5QqVRQqVQYNWoUvLy8pOflj5pYsWIFAgICoFQqoVQqoVarsXPnTqM2SUlJGDBgAJycnKBUKtGvXz/cunVLqr9+/TrCwsKgVCrh4uKCcePGobCw0KiP48eP48knn4SDgwO8vb2xaNGiGsVLREREjVO1Z5DWrl1bH3EAAFq1aoUFCxagXbt2EEJg/fr1GDZsGI4ePYpOnTohKSkJgwcPRkxMDJYtWwZbW1scO3YMNjb/y/PCwsKQk5ODuLg4lJWV4c0338SECROwceNGAHd20Rw0aBCCg4OxcuVKpKenY+zYsXBxccGECRPqbWxERERkPWp8q5GG4urqik8//RTjxo3DE088gWeeeQYff/xxpW1PnToFf39/JCcnIygoCACwa9cuPPvss7h06RK8vLywYsUKvP/++9BoNLC3twcAvPfee9i+fTv++OOPKsXEW40QEdXStfNASUHldQpnwK1Nw8ZDD4UGudVIfdPr9diyZQuKioqgVquRl5eHQ4cOISwsDL1798b58+fRoUMHfPLJJ+jbty+AO6ffXFxcpOQIAIKDg2FjY4NDhw7hxRdfRFJSEvr16yclRwAQEhKChQsX4saNG2jWrFmFWEpKSlBSUiI91+l09ThyIqJG7tp5YFng/dtMTmWSRGZlcZf5p6eno2nTplAoFIiIiMC2bdvg7++PP//8EwDw0UcfYfz48di1axcCAwMxcOBAnD17FgCg0Wjg7u5u1J+trS1cXV2h0WikNh4eHkZtyp+XtzE1f/58o/VV3t68MzwRUY3da+aoum2I6pHFJUjt27dHWloaDh06hIkTJyI8PBwnT56EwWAAALz99tt488030b17dyxevBjt27fHt99+W68xxcTEQKvVSo+srKx6PR4RERGZl8WdYrO3t0fbtm0BAD169EBycjKWLFmC9957DwDg7+9v1L5jx464ePEiAMDT0xN5eXlG9bdv38b169fh6ekptTG9NUr58/I2phQKhbStARERETV+FjeDZMpgMKCkpAStW7eGl5cXTp8+bVR/5swZ+Pr6AgDUajXy8/ORkpIi1e/duxcGgwG9evWS2uzbtw9lZWVSm7i4OLRv377S9UdERET08LGoBCkmJgb79u3DX3/9hfT0dMTExCAhIQFhYWGQyWSYPn06li5diq1bt+LcuXOYOXMm/vjjD4wbNw7AndmkwYMHY/z48Th8+DB+//13TJo0CSNHjoSXlxcA4PXXX4e9vT3GjRuHEydO4Pvvv8eSJUswbdo0cw6diIiILIhFnWLLy8vDmDFjkJOTA5VKhYCAAOzevRvPPPMMAGDq1KkoLi5GVFQUrl+/jq5duyIuLg5t2vzvSocNGzZg0qRJGDhwIGxsbDBixAgsXbpUqlepVPjll18QGRmJHj16oHnz5pg1axb3QCIiIiKJxe+DZIm4DxIRUS3wMn8yk0axDxIRETVSbm3uJEDcKJIsGBMkIiJqeEyAyMJZ1CJtIiIiIkvABImIiIjIBBMkIiIiIhNMkIiIiIhMMEEiIiIiMsEEiR46pbcNOJNbgKzrN8FtwIiIqDK8zJ8eGkIIrPktE18lnEP+zTv34uvY0hkfDe2EXo+6mTk6IiKyJJxBoofGsr3n8MnPp6TkCABOawow6h+HcCwr33yBERGRxWGCRA8FXXEZvvr1XIVyg7jzWBp/1gxRERGRpWKCRA+FI39dR8ltQ6V1eoNA4pkrXI9EREQSJkj0UJDJZA+ob6BAiIjIKjBBoodCz9aucLSTV1onl8kwoIP7A5MoIiJ6eDBBoodCU4Utpga3q1BuIwNs5TJMGfiYGaIiIiJLxcv86aExod+jaOpgi2XxZ6HRlQAAAn2b4YNQf/h7Kc0cHRERWRImSPTQkMlkCOvli5E9fZCjvQUHOzmaN1WYOywiIrJATJDooSO3kaFVsybmDoOIiCwY1yARERERmWCCRERERGSCCRIRERGRCSZIRERERCaYIBERERGZYIJEREREZIIJEhEREZEJJkhEREREJrhRJBERWb5r54GSgsrrFM6AW5uGjYcaPSZIRERk2a6dB5YF3r/N5FQmSVSneIqNiIgs271mjqrbhqgamCARERERmWCCRERERGTCohKkFStWICAgAEqlEkqlEmq1Gjt37qzQTgiBIUOGQCaTYfv27UZ1Fy9eRGhoKJo0aQJ3d3dMnz4dt2/fNmqTkJCAwMBAKBQKtG3bFuvWravHUREREZG1sagEqVWrVliwYAFSUlJw5MgRDBgwAMOGDcOJEyeM2n355ZeQyWQVXq/X6xEaGorS0lIcOHAA69evx7p16zBr1iypTWZmJkJDQ/H0008jLS0NU6dOxVtvvYXdu3fX+/iIiIjIOsiEEMLcQdyPq6srPv30U4wbNw4AkJaWhueeew5HjhxBy5YtsW3bNrzwwgsAgJ07d+K5555DdnY2PDw8AAArV65EdHQ0rly5Ant7e0RHRyM2NhYZGRnSMUaOHIn8/Hzs2rWr0hhKSkpQUlIiPdfpdPD29oZWq4VSqaynkRMREQAgOw34pv/920xIBLy6NUQ0ZMV0Oh1UKlWVvr8tagbpbnq9Hps2bUJRURHUajUA4ObNm3j99dfx1VdfwdPTs8JrkpKS0KVLFyk5AoCQkBDodDppFiopKQnBwcFGrwsJCUFSUtI9Y5k/fz5UKpX08Pb2roshEhFRVSic66YNUTVY3D5I6enpUKvVKC4uRtOmTbFt2zb4+/sDAKKiotC7d28MGzas0tdqNBqj5AiA9Fyj0dy3jU6nw61bt+Do6Fih35iYGEybNk16Xj6DREREDcCtzZ19jrhRJDUgi0uQ2rdvj7S0NGi1WmzduhXh4eFITEzEuXPnsHfvXhw9erTBY1IoFFAoFA1+XCIi+j9MgKiBWVyCZG9vj7Zt2wIAevTogeTkZCxZsgSOjo44f/48XFxcjNqPGDECTz75JBISEuDp6YnDhw8b1efm5gKAdErO09NTKru7jVKprHT2iIiIiB4+FrsGqZzBYEBJSQnee+89HD9+HGlpadIDABYvXoy1a9cCANRqNdLT05GXlye9Pi4uDkqlUjpNp1arER8fb3SMuLg4aZ0TERERkUXNIMXExGDIkCHw8fFBQUEBNm7ciISEBOzevRuenp6VLsz28fGBn58fAGDQoEHw9/fH6NGjsWjRImg0GnzwwQeIjIyUTpFFRERg+fLlmDFjBsaOHYu9e/di8+bNiI2NbdCxEhERkeWyqAQpLy8PY8aMQU5ODlQqFQICArB7924888wzVXq9XC7Hjh07MHHiRKjVajg5OSE8PBxz5syR2vj5+SE2NhZRUVFYsmQJWrVqhTVr1iAkJKS+hkVERERWxuL3QbJE1dlHgYiIiCxDo9gHicgSCSFQXKYH/64gImrcLOoUG5GlEkLgnwcvYPW+P5F14xaa2MvxUo9WiAp+DM2c7M0dHhER1THOIBFVwZz/nsSs/5xA1o1bAICbpXpsOHQRI1YeQEFxmZmjIyKiusYEiegB/rpahLUH/qpQrjcIZF4twvfJWQ0fFBER1SsmSEQPsPuEBjayyuuEAHYcy2nYgIiIqN4xQSJ6gNLbBshk98iQABTf1jdgNERE1BCYIBE9wBNt3KA3VH7VmtxGhr5tmzdwREREVN+YIBE9QJBvMzzu5wq5ySSSjQxQ2NogvHdrs8RFRET1hwkS0QPIZDKsCQ/CoE6euPtM26MtmmLThCfg7drEfMEREVG94D5IRFWgdLDDilE9kJ1/C+fyCuHqZI9OXsr7rk0iIiLrxQSJqBq8XBzh5eJo7jCIiKie8RQbERERkQkmSEREREQmmCARERERmWCCRERERGSCCRIRERGRCSZIRERERCaYIBERERGZYIJEREREZIIbRRIRUeN37TxQUlB5ncIZcGvTsPGQxWOCREREjdu188CywPu3mZzKJImM8BQbERE1bveaOapuG3qoMEEiIiIiMsEEiYiIiMgEEyQiIiIiE0yQiIiIiEwwQSIiIiIywQSJiIgaN4Vz3bShhwr3QSIiosbNrc2dfY64USRVAxMkIiJq/JgAUTVZ1Cm2FStWICAgAEqlEkqlEmq1Gjt37gQAXL9+HZMnT0b79u3h6OgIHx8fvPPOO9BqtUZ9XLx4EaGhoWjSpAnc3d0xffp03L5926hNQkICAgMDoVAo0LZtW6xbt66hhkhERERWwKJmkFq1aoUFCxagXbt2EEJg/fr1GDZsGI4ePQohBLKzs/HZZ5/B398fFy5cQEREBLKzs7F161YAgF6vR2hoKDw9PXHgwAHk5ORgzJgxsLOzw7x58wAAmZmZCA0NRUREBDZs2ID4+Hi89dZbaNmyJUJCQsw5fCIiIrIQMiGEMHcQ9+Pq6opPP/0U48aNq1C3ZcsWjBo1CkVFRbC1tcXOnTvx3HPPITs7Gx4eHgCAlStXIjo6GleuXIG9vT2io6MRGxuLjIwMqZ+RI0ciPz8fu3btqjSGkpISlJSUSM91Oh28vb2h1WqhVCrreMRERERUH3Q6HVQqVZW+vy3qFNvd9Ho9Nm3ahKKiIqjV6krblA/Q1vbORFhSUhK6dOkiJUcAEBISAp1OhxMnTkhtgoODjfoJCQlBUlLSPWOZP38+VCqV9PD29q7t8IiIiMiCWVyClJ6ejqZNm0KhUCAiIgLbtm2Dv79/hXZXr17Fxx9/jAkTJkhlGo3GKDkCID3XaDT3baPT6XDr1q1KY4qJiYFWq5UeWVlZtRojERERWTaLWoMEAO3bt0daWhq0Wi22bt2K8PBwJCYmGiVJOp0OoaGh8Pf3x0cffVTvMSkUCigUino/Dj18tDfLcPxyPhS2cnT3cYGd3OL+ZiEieihZXIJkb2+Ptm3bAgB69OiB5ORkLFmyBKtWrQIAFBQUYPDgwXB2dsa2bdtgZ2cnvdbT0xOHDx826i83N1eqK/9vedndbZRKJRwdHettXER3u603YNHu01j3+18o1RsAAG5O9vjw+U54vquXmaMjIiKL/3PVYDBIC6R1Oh0GDRoEe3t7/PTTT3BwcDBqq1arkZ6ejry8PKksLi4OSqVSmoFSq9WIj483el1cXNw91zkR1YeFu/7A6n1/SskRAFwrKsWUfx9F4pkrZoyMiIgAC0uQYmJisG/fPvz1119IT09HTEwMEhISEBYWJiVHRUVF+Mc//gGdTgeNRgONRgO9Xg8AGDRoEPz9/TF69GgcO3YMu3fvxgcffIDIyEjpFFlERAT+/PNPzJgxA3/88Qe+/vprbN68GVFRUeYcOj1E8m+WYt2Bv1DZ5aMyGbAs/myDx0RERMYs6hRbXl4exowZg5ycHKhUKgQEBGD37t145plnkJCQgEOHDgGAdAquXGZmJlq3bg25XI4dO3Zg4sSJUKvVcHJyQnh4OObMmSO19fPzQ2xsLKKiorBkyRK0atUKa9as4R5I1GCOXdKiTF/57hoGAaRcuAG9QUBuI2vgyIiIqJzF74NkiaqzjwKRqaTz1/Da6oP3rLeTy3Bm7hDIZEyQiIjqUqPYB4moserh2wzNmthVWie3kWFwJ08mR0REZsYEiaiB2dva4MOhnQAAd59Fk9vI0MRejqhnHjNTZEREVM6i1iARPSxe6P4IlI62WLLnLI5d0kozR9MGPYZHWzQ1d3hERA89JkhEZjKggwcGdPBAmd4AuUwGGy7KJiKyGEyQiMyMu2cTEVke/mYmIiIiMsEEiYiIiMgEEyQiIiIiE0yQiIiIiEwwQSIiIiIywQSJiIiIyAQTJCIiIiITTJCIiIiITDBBIiIiIjLBBImIiIjIBBMkIiIiIhNMkIiIiIhMMEEiIiIiMsEEiYiIiMgEEyQiIiIiE0yQiIiIiEwwQSIiIiIywQSJiIiIyAQTJCIiIiITTJCIiIiITDBBIiIiIjLBBImIiIjIBBMkIiIiIhNMkIiIiIhMMEEiIiIiMmFRCdKKFSsQEBAApVIJpVIJtVqNnTt3SvXFxcWIjIyEm5sbmjZtihEjRiA3N9eoj4sXLyI0NBRNmjSBu7s7pk+fjtu3bxu1SUhIQGBgIBQKBdq2bYt169Y1xPCIiIjISlhUgtSqVSssWLAAKSkpOHLkCAYMGIBhw4bhxIkTAICoqCj897//xZYtW5CYmIjs7GwMHz5cer1er0doaChKS0tx4MABrF+/HuvWrcOsWbOkNpmZmQgNDcXTTz+NtLQ0TJ06FW+99RZ2797d4OMlIiIiyyQTQghzB3E/rq6u+PTTT/HSSy+hRYsW2LhxI1566SUAwB9//IGOHTsiKSkJTzzxBHbu3InnnnsO2dnZ8PDwAACsXLkS0dHRuHLlCuzt7REdHY3Y2FhkZGRIxxg5ciTy8/Oxa9euKsWk0+mgUqmg1WqhVCrrftBERERU56rz/W3bQDFVm16vx5YtW1BUVAS1Wo2UlBSUlZUhODhYatOhQwf4+PhICVJSUhK6dOkiJUcAEBISgokTJ+LEiRPo3r07kpKSjPoobzN16tR7xlJSUoKSkhLpuVarBXDnjSYiIiLrUP69XZW5IYtLkNLT06FWq1FcXIymTZti27Zt8Pf3R1paGuzt7eHi4mLU3sPDAxqNBgCg0WiMkqPy+vK6+7XR6XS4desWHB0dK8Q0f/58zJ49u0K5t7d3jcdJRERE5lFQUACVSnXfNhaXILVv3x5paWnQarXYunUrwsPDkZiYaNaYYmJiMG3aNOl5fn4+fH19cfHixQe+wdZKp9PB29sbWVlZjfI0Isdn/Rr7GBv7+IDGP0aOz/IIIVBQUAAvL68HtrW4BMne3h5t27YFAPTo0QPJyclYsmQJXn31VZSWliI/P99oFik3Nxeenp4AAE9PTxw+fNiov/Kr3O5uY3rlW25uLpRKZaWzRwCgUCigUCgqlKtUKqv5oaip8isKGyuOz/o19jE29vEBjX+MHJ9lqerEhkVdxVYZg8GAkpIS9OjRA3Z2doiPj5fqTp8+jYsXL0KtVgMA1Go10tPTkZeXJ7WJi4uDUqmEv7+/1ObuPsrblPdBREREZFEzSDExMRgyZAh8fHxQUFCAjRs3IiEhAbt374ZKpcK4ceMwbdo0uLq6QqlUYvLkyVCr1XjiiScAAIMGDYK/vz9Gjx6NRYsWQaPR4IMPPkBkZKQ0AxQREYHly5djxowZGDt2LPbu3YvNmzcjNjbWnEMnIiIiC2JRCVJeXh7GjBmDnJwcqFQqBAQEYPfu3XjmmWcAAIsXL4aNjQ1GjBiBkpIShISE4Ouvv5ZeL5fLsWPHDkycOBFqtRpOTk4IDw/HnDlzpDZ+fn6IjY1FVFQUlixZglatWmHNmjUICQmpcpwKhQIffvhhpafdGovGPkaOz/o19jE29vEBjX+MHJ91s/h9kIiIiIgamsWvQSIiIiJqaEyQiIiIiEwwQSIiIiIywQSJiIiIyAQTpBr46quv0Lp1azg4OKBXr14VNqe0FvPnz0fPnj3h7OwMd3d3vPDCCzh9+rRRm+LiYkRGRsLNzQ1NmzbFiBEjKmy0aS0WLFgAmUxmdN+9xjC+y5cvY9SoUXBzc4OjoyO6dOmCI0eOSPVCCMyaNQstW7aEo6MjgoODcfbsWTNGXHV6vR4zZ86En58fHB0d0aZNG3z88cdG91GytvHt27cPQ4cOhZeXF2QyGbZv325UX5XxXL9+HWFhYVAqlXBxccG4ceNQWFjYgKO4t/uNr6ysDNHR0ejSpQucnJzg5eWFMWPGIDs726gPax2fqYiICMhkMnz55ZdG5ZY8PqBqYzx16hSef/55qFQqODk5oWfPnrh48aJU3xh+tzJBqqbvv/8e06ZNw4cffojU1FR07doVISEhRptTWovExERERkbi4MGDiIuLQ1lZGQYNGoSioiKpTVRUFP773/9iy5YtSExMRHZ2NoYPH27GqGsmOTkZq1atQkBAgFG5tY/vxo0b6NOnD+zs7LBz506cPHkSn3/+OZo1aya1WbRoEZYuXYqVK1fi0KFDcHJyQkhICIqLi80YedUsXLgQK1aswPLly3Hq1CksXLgQixYtwrJly6Q21ja+oqIidO3aFV999VWl9VUZT1hYGE6cOIG4uDjs2LED+/btw4QJExpqCPd1v/HdvHkTqampmDlzJlJTU/Hjjz/i9OnTeP75543aWev47rZt2zYcPHiw0ltaWPL4gAeP8fz58+jbty86dOiAhIQEHD9+HDNnzoSDg4PUxtp/twIABFXL448/LiIjI6Xner1eeHl5ifnz55sxqrqRl5cnAIjExEQhhBD5+fnCzs5ObNmyRWpz6tQpAUAkJSWZK8xqKygoEO3atRNxcXGif//+YsqUKUKIxjG+6Oho0bdv33vWGwwG4enpKT799FOpLD8/XygUCvHvf/+7IUKsldDQUDF27FijsuHDh4uwsDAhhPWPD4DYtm2b9Lwq4zl58qQAIJKTk6U2O3fuFDKZTFy+fLnBYq8K0/FV5vDhwwKAuHDhghCicYzv0qVL4pFHHhEZGRnC19dXLF68WKqzpvEJUfkYX331VTFq1Kh7vqYx/G4VQgjOIFVDaWkpUlJSEBwcLJXZ2NggODgYSUlJZoysbmi1WgCAq6srACAlJQVlZWVG4+3QoQN8fHysaryRkZEIDQ01GgfQOMb3008/ISgoCC+//DLc3d3RvXt3rF69WqrPzMyERqMxGqNKpUKvXr2sYoy9e/dGfHw8zpw5AwA4duwY9u/fjyFDhgCw/vGZqsp4kpKS4OLigqCgIKlNcHAwbGxscOjQoQaPuba0Wi1kMpl0j01rH5/BYMDo0aMxffp0dOrUqUJ9YxhfbGwsHnvsMYSEhMDd3R29evUyOg3XGH63AjzFVi1Xr16FXq+Hh4eHUbmHhwc0Go2ZoqobBoMBU6dORZ8+fdC5c2cAgEajgb29vdHNgQHrGu+mTZuQmpqK+fPnV6hrDOP7888/sWLFCrRr1w67d+/GxIkT8c4772D9+vUAII3DWn9m33vvPYwcORIdOnSAnZ0dunfvjqlTpyIsLAyA9Y/PVFXGo9Fo4O7ublRva2sLV1dXqxtzcXExoqOj8dprr0k3O7X28S1cuBC2trZ45513Kq239vHl5eWhsLAQCxYswODBg/HLL7/gxRdfxPDhw5GYmAigcfxuBSzsViNkPpGRkcjIyMD+/fvNHUqdycrKwpQpUxAXF2d0brwxMRgMCAoKwrx58wAA3bt3R0ZGBlauXInw8HAzR1d7mzdvxoYNG7Bx40Z06tQJaWlpmDp1Kry8vBrF+B5mZWVleOWVVyCEwIoVK8wdTp1ISUnBkiVLkJqaCplMZu5w6oXBYAAADBs2DFFRUQCAbt264cCBA1i5ciX69+9vzvDqFGeQqqF58+aQy+UVVuLn5ubC09PTTFHV3qRJk7Bjxw78+uuvaNWqlVTu6emJ0tJS5OfnG7W3lvGmpKQgLy8PgYGBsLW1ha2tLRITE7F06VLY2trCw8PDqscHAC1btoS/v79RWceOHaWrScrHYa0/s9OnT5dmkbp06YLRo0cjKipKmhG09vGZqsp4PD09K1wUcvv2bVy/ft1qxlyeHF24cAFxcXHS7BFg3eP77bffkJeXBx8fH+l3zoULF/D3v/8drVu3BmDd4wPufA/a2to+8PeOtf9uBZggVYu9vT169OiB+Ph4qcxgMCA+Ph5qtdqMkdWMEAKTJk3Ctm3bsHfvXvj5+RnV9+jRA3Z2dkbjPX36NC5evGgV4x04cCDS09ORlpYmPYKCghAWFib9vzWPDwD69OlTYWuGM2fOwNfXF8CdmzN7enoajVGn0+HQoUNWMcabN2/Cxsb415RcLpf+irX28ZmqynjUajXy8/ORkpIitdm7dy8MBgN69erV4DFXV3lydPbsWezZswdubm5G9dY8vtGjR+P48eNGv3O8vLwwffp07N69G4B1jw+48z3Ys2fP+/7esfbvDom5V4lbm02bNgmFQiHWrVsnTp48KSZMmCBcXFyERqMxd2jVNnHiRKFSqURCQoLIycmRHjdv3pTaRERECB8fH7F3715x5MgRoVarhVqtNmPUtXP3VWxCWP/4Dh8+LGxtbcUnn3wizp49KzZs2CCaNGki/vWvf0ltFixYIFxcXMR//vMfcfz4cTFs2DDh5+cnbt26ZcbIqyY8PFw88sgjYseOHSIzM1P8+OOPonnz5mLGjBlSG2sbX0FBgTh69Kg4evSoACC++OILcfToUekqrqqMZ/DgwaJ79+7i0KFDYv/+/aJdu3bitddeM9eQjNxvfKWlpeL5558XrVq1EmlpaUa/d0pKSqQ+rHV8lTG9ik0Iyx6fEA8e448//ijs7OzEN998I86ePSuWLVsm5HK5+O2336Q+rP13qxBCMEGqgWXLlgkfHx9hb28vHn/8cXHw4EFzh1QjACp9rF27Vmpz69Yt8be//U00a9ZMNGnSRLz44osiJyfHfEHXkmmC1BjG99///ld07txZKBQK0aFDB/HNN98Y1RsMBjFz5kzh4eEhFAqFGDhwoDh9+rSZoq0enU4npkyZInx8fISDg4N49NFHxfvvv2/0ZWpt4/v1118r/XcXHh4uhKjaeK5duyZee+010bRpU6FUKsWbb74pCgoKzDCaiu43vszMzHv+3vn111+lPqx1fJWpLEGy5PEJUbUx/uMf/xBt27YVDg4OomvXrmL79u1GfTSG360yIe7akpaIiIiIuAaJiIiIyBQTJCIiIiITTJCIiIiITDBBIiIiIjLBBImIiIjIBBMkIiIiIhNMkIiIiIhMMEEiIiIiMsEEiYiIiMgEEyQiIiIiE0yQiIiszFNPPYWpU6eaOwyiRo0JElEVvPHGG5DJZBUe586dM3dojdZHH31U4f3u0KGDUZv58+ejZ8+ecHZ2hru7O1544QWcPn26ysdYsGABZDJZhWRDr9dj5syZ8PPzg6OjI9q0aYOPP/4Y1nLryn379mHo0KHw8vKCTCbD9u3bK7Rp3bp1pT/TkZGRRu0uX76MUaNGwc3NDY6OjujSpQuOHDly3+N/9dVXaN26NRwcHNCrVy8cPnzYqH7FihUICAiAUqmEUqmEWq3Gzp07az3uN954Ay+88EKF8oSEBMhkMuTn59f6GPTwYIJEVEWDBw9GTk6O0cPPz69Cu9LSUjNEZ/meeuoprFu3rlqv6dSpk9H7vX//fqP6xMREREZG4uDBg4iLi0NZWRkGDRqEoqKiB/adnJyMVatWISAgoELdwoULsWLFCixfvhynTp3CwoULsWjRIixbtqxa8ZtLUVERunbtiq+++uqebZKTk43e27i4OADAyy+/LLW5ceMG+vTpAzs7O+zcuRMnT57E559/jmbNmt2z3++//x7Tpk3Dhx9+iNTUVHTt2hUhISHIy8uT2rRq1QoLFixASkoKjhw5ggEDBmDYsGE4ceJEHYyeqI4IInqg8PBwMWzYsErr+vfvLyIjI8WUKVOEm5ubeOqpp4QQQuj1ejFv3jzRunVr4eDgIAICAsSWLVuMXltYWChGjx4tnJychKenp/jss89E//79xZQpU6Q2vr6+YvHixUav69q1q/jwww+rfJz+/fuLyZMni+nTp4tmzZoJDw8P6fXl9Hq9WLhwoWjTpo2wt7cX3t7eYu7cuWL9+vXC1dVVFBcXG7UfNmyYGDVqVNXewP+LYe3atVVu/+GHH4quXbtWub0QQuTl5QkAIjEx8b7tCgoKRLt27URcXFyF91sIIUJDQ8XYsWONyoYPHy7CwsLu2Wf//v3FpEmTxJQpU4SLi4twd3cX33zzjSgsLBRvvPGGaNq0qWjTpo34+eefjV5XXFwsJk+eLFq0aCEUCoXo06ePOHz4sFRflZ+R+wEgtm3b9sB2U6ZMEW3atBEGg0Eqi46OFn379q3Scco9/vjjIjIyUnqu1+uFl5eXmD9//n1f16xZM7FmzZp71v/222/C1tZW3Lp1SyrLzMwUAMRff/0lhLj3v9Nff/1VABA3btwwep3po3///lUfKDV6nEEiqgPr16+Hvb09fv/9d6xcuRLAndM/3333HVauXIkTJ04gKioKo0aNQmJiovS66dOnIzExEf/5z3/wyy+/ICEhAampqdU6dlWOUx6jk5MTDh06hEWLFmHOnDnSrAEAxMTEYMGCBZg5cyZOnjyJjRs3wsPDAy+//DL0ej1++uknqW1eXh5iY2MxduzYmrxdVXb27Fl4eXnh0UcfRVhYGC5evHjf9lqtFgDg6up633aRkZEIDQ1FcHBwpfW9e/dGfHw8zpw5AwA4duwY9u/fjyFDhty33/Xr16N58+Y4fPgwJk+ejIkTJ+Lll19G7969kZqaikGDBmH06NG4efOm9JoZM2bghx9+wPr165Gamoq2bdsiJCQE169fB1A3PyMPUlpain/9618YO3YsZDKZVP7TTz8hKCgIL7/8Mtzd3dG9e3esXr36vv2kpKQYva82NjYIDg5GUlJSpa/R6/XYtGkTioqKoFar79l3WloaOnbsCAcHB6ns6NGjaNasGXx9faszXHh7exvNnh09ehRubm7o169ftfqhRs7cGRqRNQgPDxdyuVw4OTlJj5deekkIcWfmoHv37kbti4uLRZMmTcSBAweMyseNGydee+01IcSdWQx7e3uxefNmqf7atWvC0dGxyjNIVTlOeYymMwE9e/YU0dHRQgghdDqdUCgUYvXq1ZWOf+LEiWLIkCHS888//1w8+uijRrMND1LdGaSff/5ZbN68WRw7dkzs2rVLqNVq4ePjI3Q6XaXt9Xq9CA0NFX369Llvv//+979F586dpZmIymZj9Hq9iI6OFjKZTNja2gqZTCbmzZv3wPHd/R7fvn1bODk5idGjR0tlOTk5AoBISkoSQtyZHbKzsxMbNmyQ2pSWlgovLy+xaNGiKv+M3A+qMIP0/fffC7lcLi5fvmxUrlAohEKhEDExMSI1NVWsWrVKODg4iHXr1lXaz+XLlwWACj+P06dPF48//rhR2fHjx4WTk5OQy+VCpVKJ2NjY+8b41ltviTFjxhiVzZo1S5qxFaLyf6dOTk7CwcHBaAbpbrdu3RK9evUSzz33nNDr9feNgR4utuZNz4isx9NPP40VK1ZIz52cnKT/79Gjh1Hbc+fO4ebNm3jmmWeMyktLS9G9e3cAwPnz51FaWopevXpJ9a6urmjfvn2VY6rKccqZrrVp2bKltC7k1KlTKCkpwcCBAys9zvjx49GzZ09cvnwZjzzyCNatWyctXL+XefPmYd68edLzW7du4eDBg5g0aZJUdvLkSfj4+FT6+rtnawICAtCrVy/4+vpi8+bNGDduXIX2kZGRyMjIqLBO6W5ZWVmYMmUK4uLijGYiTG3evBkbNmzAxo0b0alTJ6SlpWHq1Knw8vJCeHj4PV9393ssl8vh5uaGLl26SGUeHh4AIL3v58+fR1lZGfr06SO1sbOzw+OPP45Tp07Vyc9IVfzjH//AkCFD4OXlZVRuMBgQFBQkfY7du3dHRkYGVq5ced/3oSrat2+PtLQ0aLVabN26FeHh4UhMTIS/v3+l7dPS0vD6668blR09ehTdunUzKjP9dwoAhw4dwqhRoyrtd+zYsSgoKEBcXBxsbHhShf6HCRJRFTk5OaFt27b3rLtbYWEhACA2NhaPPPKIUZ1CoajWcW1sbCpcPVVWVlbt49jZ2Rk9l8lkMBgMAABHR8f7xtC9e3d07doV3333HQYNGoQTJ04gNjb2vq+JiIjAK6+8Ij0PCwvDiBEjMHz4cKnM9Av5flxcXPDYY49VeuXgpEmTsGPHDuzbtw+tWrW6Zx8pKSnIy8tDYGCgVKbX67Fv3z4sX74cJSUlkMvlmD59Ot577z2MHDkSANClSxdcuHAB8+fPv29iUNl7fHdZeUJZ/r5bggsXLmDPnj348ccfK9S1bNmyQsLSsWNH/PDDD5X21bx5c8jlcuTm5hqV5+bmwtPT06jM3t5e+vfUo0cPJCcnY8mSJVi1alWFfvV6PTIyMiok/ampqRgxYoRRWWX/Ti9dulRpvHPnzsXu3btx+PBhODs7V9qGHl5Ml4nqgb+/PxQKBS5evIi2bdsaPby9vQEAbdq0gZ2dHQ4dOiS97saNG9K6l3ItWrRATk6O9Fyn0yEzM7PKx6mKdu3awdHREfHx8fds89Zbb2HdunVYu3YtgoODH9i/q6urUTyOjo5wd3c3KrO1rfrfaIWFhTh//jxatmwplQkhMGnSJGzbtg179+6t9KrCuw0cOBDp6elIS0uTHkFBQQgLC0NaWhrkcjkA4ObNmxVmE+RyeZ0nNm3atJHWrpUrKytDcnIy/P39q/wzUhtr166Fu7s7QkNDK9T16dOnwrYJZ86cueeaH3t7e/To0cPo58hgMCA+Pv6+64vK25WUlFRad/r0aRQXFxsl1ElJSbh8+XKFGaSq+uGHHzBnzhxs3rwZbdq0qVEf1LhxBomoHjg7O+Pdd99FVFQUDAYD+vbtC61Wi99//x1KpRLh4eFo2rQpxo0bh+nTp8PNzQ3u7u54//33K3wxDxgwAOvWrcPQoUPh4uKCWbNmSV/kVTlOVTg4OCA6OhozZsyAvb09+vTpgytXruDEiRPS6azXX38d7777LlavXo3vvvuubt+wSrz77rsYOnQofH19kZ2djQ8//BByuRyvvfaa1CYyMhIbN27Ef/7zHzg7O0Oj0QAAVCqVNCu2fPlybNu2DfHx8XB2dkbnzp2NjuPk5AQ3Nzej8qFDh+KTTz6Bj48POnXqhKNHj+KLL76o80XpTk5OmDhxIqZPnw5XV1f4+Phg0aJFuHnzJsaNG1flnxFThYWFRjNtmZmZSEtLk45RzmAwYO3atQgPD680WY2KikLv3r0xb948vPLKKzh8+DC++eYbfPPNN1Kbu99fAJg2bRrCw8MRFBSExx9/HF9++SWKiorw5ptvSq+JiYnBkCFD4OPjg4KCAmzcuBEJCQnYvXt3peNJS0sDACxbtgzvvPMOzp07h3feeQdAzbbVyMjIwJgxYxAdHY1OnTpJPzf29vYPXOBPDw8mSET15OOPP0aLFi0wf/58/Pnnn3BxcUFgYCD+3//7f1KbTz/9FIWFhRg6dCicnZ3x97//XboSq1xMTAwyMzPx3HPPQaVS4eOPP5ZmkKp6nKqYOXMmbG1tMWvWLGRnZ6Nly5aIiIiQ6lUqFUaMGIHY2NhKN+Ora5cuXcJrr72Ga9euoUWLFujbty8OHjyIFi1aSG3K15o89dRTRq9du3Yt3njjDQDA1atXcf78+Wode9myZZg5cyb+9re/IS8vD15eXnj77bcxa9asWo2pMgsWLIDBYMDo0aNRUFCAoKAg7N69W9prqCo/I6aOHDmCp59+Wno+bdo0AEB4eLjRXlR79uzBxYsX75n49ezZE9u2bUNMTAzmzJkDPz8/fPnllwgLC5PamL6/r776Kq5cuYJZs2ZBo9GgW7du2LVrl7T+CrizBmvMmDHIycmBSqVCQEAAdu/eXWEtXbm0tDSEhITgzz//RJcuXeDv74/Zs2dj4sSJWLp0Kf75z3/e9/2o7P25efMm5s6di7lz50rl/fv3R0JCQrX6osZLJkwXNxCRWT311FPo1q0bvvzyS3OHUsHAgQPRqVMnLF261Nyh0EMkJCQEPXv2NEpmiOob1yAR0QPduHED27ZtQ0JCQoVbURDVt2PHjhldDUjUEHiKjYgeqHv37rhx4wYWLlxY55eYE92PRqNBbm4uEyRqcDzFRkRERGSCp9iIiIiITDBBIiIiIjLBBImIiIjIBBMkIiIiIhNMkIiIiIhMMEEiIiIiMsEEiYiIiMgEEyQiIiIiE0yQiIiIiEwwQSIiIiIy8f8B2Zq6999ZBWsAAAAASUVORK5CYII=", "text/plain": [ "
" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "from wsssss.plotting import plotting as pl\n", "\n", "hist = ld.History(f'{grid_dir}/m1.000_fOV0.020/LOGS/history.data')\n", "profs = ld.load_profs(hist)\n", "gss = ld.load_gss(hist)\n", "\n", "pl.make_kipp(hist)\n", "pl.make_structural(profs[-1])\n", "pl.make_echelle(gss[-1], hist)" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "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.11.2" } }, "nbformat": 4, "nbformat_minor": 5 }