{"metadata":{"kernelspec":{"display_name":"Python 3 (ipykernel)","language":"python","name":"python3"},"language_info":{"name":"python","version":"3.10.16","mimetype":"text/x-python","codemirror_mode":{"name":"ipython","version":3},"pygments_lexer":"ipython3","nbconvert_exporter":"python","file_extension":".py"}},"nbformat_minor":5,"nbformat":4,"cells":[{"id":"ce2efe0c-f86c-4cbd-ab3f-2aae9c5a574d","cell_type":"markdown","source":"## Thermal Expansion \nCalculate the thermal expansion for a Morse Pair potential using the [LAMMPS](https://www.lammps.org/) molecular dynamics\nsimulation code. In the following three methods to calculate the thermal expansion are introduced and compared for a \nMorse Pair Potential for Aluminium. \n\nAs a first step the potential is defined for the [LAMMPS](https://www.lammps.org/) molecular dynamics simulation code \nby specifying the `pair_style` and `pair_coeff` commands for the [Morse Pair Potential](https://docs.lammps.org/pair_morse.html)\nas well as the Aluminium bulk structure: ","metadata":{}},{"id":"7a96bac5-5afe-4924-90e5-04f6e6b2bedb","cell_type":"code","source":"from ase.build import bulk\nimport pandas\n\npotential_dataframe = pandas.DataFrame(\n {\n \"Config\": [\n [\"pair_style morse/smooth/linear 9.0\", \"pair_coeff * * 0.5 1.8 2.95\"]\n ],\n \"Filename\": [[]],\n \"Model\": [\"Morse\"],\n \"Name\": [\"Morse\"],\n \"Species\": [[\"Al\"]],\n }\n)\n\nstructure = bulk(\"Al\", cubic=True)","metadata":{"trusted":true},"outputs":[],"execution_count":1},{"id":"fe980834-e174-464a-8c07-04db9889a8c6","cell_type":"markdown","source":"The `pandas.DataFrame` based format to specify interatomic potentials is the same `pylammpsmpi` uses to interface with \nthe [NIST database for interatomic potentials](https://www.ctcms.nist.gov/potentials). In comparison to just providing\nthe `pair_style` and `pair_coeff` commands, this extended format enables referencing specific files for the interatomic\npotentials `\"Filename\": [[]],` as well as the atomic species `\"Species\": [[\"Al\"]],` to enable consistency checks if the \ninteratomic potential implements all the interactions to simulate a given atomic structure. ","metadata":{}},{"id":"b825e555-5e3d-43c2-8c51-50207ead60b5","cell_type":"markdown","source":"Finally, the last step of the preparation before starting the actual calculation is optimizing the interatomic structure. \nWhile for the Morse potential used in this example this is not necessary, it is essential for extending this example to\nother interactomic potentials. For the structure optimization the `optimize_positions_and_volume()` function is imported\nand applied on the `ase.atoms.Atoms` bulk structure for Aluminium:","metadata":{}},{"id":"7bd4ed1a-bffe-40f9-99f8-706797418877","cell_type":"code","source":"from atomistics.workflows import optimize_positions_and_volume\n\ntask_dict = optimize_positions_and_volume(structure=structure)\ntask_dict","metadata":{"trusted":true},"outputs":[{"execution_count":2,"output_type":"execute_result","data":{"text/plain":"{'optimize_positions_and_volume': Atoms(symbols='Al4', pbc=True, cell=[4.05, 4.05, 4.05])}"},"metadata":{}}],"execution_count":2},{"id":"8c4ff286010f6006","cell_type":"markdown","source":"It returns a `task_dict` with a single task, the optimization of the positions and the volume of the Aluminium structure.\nThis task is executed with the [LAMMPS](https://www.lammps.org/) molecular dynamics simulation code using the \n`evaluate_with_lammpslib()` function:","metadata":{}},{"id":"ecf124099126f2d9","cell_type":"code","source":"from atomistics.calculators import evaluate_with_lammpslib\n\nresult_dict = evaluate_with_lammpslib(\n task_dict=task_dict,\n potential_dataframe=potential_dataframe,\n)\nstructure_opt = result_dict[\"structure_with_optimized_positions_and_volume\"]\nstructure_opt","metadata":{"trusted":true},"outputs":[{"name":"stderr","output_type":"stream","text":"/srv/conda/envs/notebook/lib/python3.10/site-packages/atomistics/calculators/__init__.py:63: UserWarning: calc_static_with_qe(), evaluate_with_qe() and optimize_positions_and_volume_with_qe() are not available as the import of the module named 'pwtools' failed.\n raise_warning(module_list=quantum_espresso_function, import_error=e)\n"},{"execution_count":3,"output_type":"execute_result","data":{"text/plain":"Atoms(symbols='Al4', pbc=True, cell=[4.047310585424964, 4.047310585424964, 4.047310585424964])"},"metadata":{}}],"execution_count":3},{"id":"6750cc4c-106f-4066-80ad-860bf6980732","cell_type":"markdown","source":"The `result_dict` just contains a single element, the `ase.atoms.Atoms` structure object with optimized positions and \nvolume. After this step the preparation is completed and the three different approximations can be compared in the following.","metadata":{}},{"id":"6c120581-efd6-4204-8413-75ee81065db1","cell_type":"markdown","source":"### Equation of State \nThe first approximation to calculate the thermal expansion is based on the Equation of State derived by [Moruzzi, V. L. et al.](https://link.aps.org/doi/10.1103/PhysRevB.37.790).\nSo in analogy to the previous example of calculating the elastic properties from the Equation of State, the `EnergyVolumeCurveWorkflow`\nis initialized with the default parameters: ","metadata":{}},{"id":"b69b6ee2-b526-4913-a6c7-36018e8960af","cell_type":"code","source":"from atomistics.workflows import (\n analyse_results_for_energy_volume_curve, \n get_tasks_for_energy_volume_curve,\n get_thermal_properties_for_energy_volume_curve,\n)\n\ntask_dict = get_tasks_for_energy_volume_curve(\n structure=structure_opt.copy(),\n num_points=11,\n vol_range=0.05,\n axes=[\"x\", \"y\", \"z\"],\n)\ntask_dict","metadata":{"trusted":true},"outputs":[{"execution_count":4,"output_type":"execute_result","data":{"text/plain":"{'calc_energy': {0.95: Atoms(symbols='Al4', pbc=True, cell=[3.9786988461213992, 3.9786988461213992, 3.9786988461213992]),\n 0.96: Atoms(symbols='Al4', pbc=True, cell=[3.992610493736228, 3.992610493736228, 3.992610493736228]),\n 0.97: Atoms(symbols='Al4', pbc=True, cell=[4.00642586504517, 4.00642586504517, 4.00642586504517]),\n 0.98: Atoms(symbols='Al4', pbc=True, cell=[4.020146608667117, 4.020146608667117, 4.020146608667117]),\n 0.99: Atoms(symbols='Al4', pbc=True, cell=[4.033774328510742, 4.033774328510742, 4.033774328510742]),\n 1.0: Atoms(symbols='Al4', pbc=True, cell=[4.047310585424964, 4.047310585424964, 4.047310585424964]),\n 1.01: Atoms(symbols='Al4', pbc=True, cell=[4.060756898772644, 4.060756898772644, 4.060756898772644]),\n 1.02: Atoms(symbols='Al4', pbc=True, cell=[4.074114747931804, 4.074114747931804, 4.074114747931804]),\n 1.03: Atoms(symbols='Al4', pbc=True, cell=[4.087385573728375, 4.087385573728375, 4.087385573728375]),\n 1.04: Atoms(symbols='Al4', pbc=True, cell=[4.100570779804249, 4.100570779804249, 4.100570779804249]),\n 1.05: Atoms(symbols='Al4', pbc=True, cell=[4.113671733924125, 4.113671733924125, 4.113671733924125])}}"},"metadata":{}}],"execution_count":4},{"id":"8f5d1e8d-0204-4dca-9298-878b9b2f6406","cell_type":"markdown","source":"After the initialization the `generate_structures()` function is called to generate the atomistic structures which are\nthen in the second step evaluated with the [LAMMPS](https://www.lammps.org/) molecular dynamics simulation code to derive\nthe equilibrium properties:","metadata":{}},{"id":"7d1f126e-4fd0-41c5-986b-91d3b5910e3e","cell_type":"code","source":"result_dict = evaluate_with_lammpslib(\n task_dict=task_dict, potential_dataframe=potential_dataframe\n)\nresult_dict","metadata":{"trusted":true},"outputs":[{"execution_count":5,"output_type":"execute_result","data":{"text/plain":"{'energy': {0.95: -14.609207927145922,\n 0.96: -14.656740101454448,\n 0.97: -14.692359030099395,\n 0.98: -14.71688372487553,\n 0.99: -14.731079276327012,\n 1.0: -14.735659820057943,\n 1.01: -14.731295089579728,\n 1.02: -14.718611862249293,\n 1.03: -14.69819671584233,\n 1.04: -14.670598736769113,\n 1.05: -14.636332030744796}}"},"metadata":{}}],"execution_count":5},{"id":"fb679eb1-338f-4485-a953-791e147fe632","cell_type":"markdown","source":"While in the previous example the fit of the energy volume curve was used directly, here the output of the fit, in\nparticular the derived equilibrium properties are the input for the Debye model as defined by [Moruzzi, V. L. et al.](https://link.aps.org/doi/10.1103/PhysRevB.37.790):","metadata":{}},{"id":"2e0f7aab-6744-4b6f-a454-38c28833a3ac","cell_type":"code","source":"structure_opt.get_volume()","metadata":{"trusted":true},"outputs":[{"execution_count":6,"output_type":"execute_result","data":{"text/plain":"66.29787349319821"},"metadata":{}}],"execution_count":6},{"id":"11c8b18d-64ff-4c93-b646-668b00eb1cf8","cell_type":"code","source":"fit_dict = analyse_results_for_energy_volume_curve(\n output_dict=result_dict,\n task_dict=task_dict,\n)\nfit_dict","metadata":{"trusted":true},"outputs":[{"execution_count":7,"output_type":"execute_result","data":{"text/plain":"{'b_prime_eq': 6.218602473906521,\n 'bulkmodul_eq': 218.25686441091622,\n 'volume_eq': 66.29753112311258,\n 'energy_eq': -14.735788882589588,\n 'fit_dict': {'fit_type': 'polynomial',\n 'least_square_error': 1.2103551168057661e-08,\n 'poly_fit': array([-3.72876215e-04, 8.44360952e-02, -6.27903075e+00, 1.39077950e+02]),\n 'fit_order': 3},\n 'energy': [-14.609207927145922,\n -14.656740101454448,\n -14.692359030099395,\n -14.71688372487553,\n -14.731079276327012,\n -14.735659820057943,\n -14.731295089579728,\n -14.718611862249293,\n -14.69819671584233,\n -14.670598736769113,\n -14.636332030744796],\n 'volume': [62.98297981853827,\n 63.645958553470244,\n 64.30893728840229,\n 64.97191602333424,\n 65.63489475826624,\n 66.29787349319821,\n 66.96085222813018,\n 67.62383096306218,\n 68.28680969799419,\n 68.94978843292616,\n 69.61276716785807]}"},"metadata":{}}],"execution_count":7},{"id":"9c7cd51c-6058-4d1d-8948-56d29c3b13e7","cell_type":"code","source":"import numpy as np\n\nthermal_properties_dict = get_thermal_properties_for_energy_volume_curve(\n fit_dict=fit_dict,\n masses=structure.get_masses(),\n temperatures=np.arange(1, 1500, 50),\n output_keys=[\"temperatures\", \"volumes\"],\n)\ntemperatures_ev, volume_ev = (\n thermal_properties_dict[\"temperatures\"],\n thermal_properties_dict[\"volumes\"],\n)","metadata":{"trusted":true},"outputs":[],"execution_count":8},{"id":"35ab7b86-0688-4520-ad47-ea54b4bfde86","cell_type":"markdown","source":"The output of the Debye model provides the change of the temperature specific optimal volume `volume_ev`\nwhich can be plotted over the temperature `temperatures_ev` to determine the thermal expansion. ","metadata":{}},{"id":"88ccd1f0-98c5-4e13-ab2c-febe5d3f235b","cell_type":"markdown","source":"### Quasi-Harmonic Approximation \nWhile the [Moruzzi, V. L. et al.](https://link.aps.org/doi/10.1103/PhysRevB.37.790) approach based on the Einstein crystal\nis limited to a single frequency, the quasi-harmonic model includes the volume dependent free energy. Inside the \n`atomistics` package the harmonic and quasi-harmonic model are implemented based on an interface to the [Phonopy](https://phonopy.github.io/phonopy/)\nframework. Still the user interface is still structured in the same three steps of (1) generating structures, (2) evaluating \nthese structures and (3) fitting the corresponding model. Starting with the initialization of the `QuasiHarmonicWorkflow`\nwhich combines the `PhonopyWorkflow` with the `EnergyVolumeCurveWorkflow`:","metadata":{}},{"id":"493663b9-ea0c-4234-87ef-8f70774794f4","cell_type":"code","source":"from atomistics.workflows import (\n get_tasks_for_quasi_harmonic_approximation,\n analyse_results_for_quasi_harmonic_approximation,\n get_thermal_properties_for_quasi_harmonic_approximation,\n)\nfrom phonopy.units import VaspToTHz\n\ntask_dict, qh_dict = get_tasks_for_quasi_harmonic_approximation(\n structure=structure_opt.copy(),\n num_points=11,\n vol_range=0.10,\n # fit_type='birchmurnaghan',\n interaction_range=10,\n factor=VaspToTHz,\n displacement=0.01,\n)\ntask_dict","metadata":{"trusted":true},"outputs":[{"execution_count":9,"output_type":"execute_result","data":{"text/plain":"{'calc_energy': {0.9: Atoms(symbols='Al108', pbc=True, cell=[11.7229062192894, 11.7229062192894, 11.7229062192894]),\n 0.92: Atoms(symbols='Al108', pbc=True, cell=[11.80910715486485, 11.80910715486485, 11.80910715486485]),\n 0.94: Atoms(symbols='Al108', pbc=True, cell=[11.894067681419225, 11.894067681419225, 11.894067681419225]),\n 0.96: Atoms(symbols='Al108', pbc=True, cell=[11.977831481208684, 11.977831481208684, 11.977831481208684]),\n 0.98: Atoms(symbols='Al108', pbc=True, cell=[12.060439826001351, 12.060439826001351, 12.060439826001351]),\n 1.0: Atoms(symbols='Al108', pbc=True, cell=[12.141931756274893, 12.141931756274893, 12.141931756274893]),\n 1.02: Atoms(symbols='Al108', pbc=True, cell=[12.222344243795412, 12.222344243795412, 12.222344243795412]),\n 1.04: Atoms(symbols='Al108', pbc=True, cell=[12.301712339412747, 12.301712339412747, 12.301712339412747]),\n 1.06: Atoms(symbols='Al108', pbc=True, cell=[12.38006930767338, 12.38006930767338, 12.38006930767338]),\n 1.08: Atoms(symbols='Al108', pbc=True, cell=[12.457446749652004, 12.457446749652004, 12.457446749652004]),\n 1.1: Atoms(symbols='Al108', pbc=True, cell=[12.533874715230777, 12.533874715230777, 12.533874715230777])},\n 'calc_forces': {(0.9,\n 0): Atoms(symbols='Al108', pbc=True, cell=[11.7229062192894, 11.7229062192894, 11.7229062192894]),\n (0.92,\n 0): Atoms(symbols='Al108', pbc=True, cell=[11.80910715486485, 11.80910715486485, 11.80910715486485]),\n (0.94,\n 0): Atoms(symbols='Al108', pbc=True, cell=[11.894067681419225, 11.894067681419225, 11.894067681419225]),\n (0.96,\n 0): Atoms(symbols='Al108', pbc=True, cell=[11.977831481208684, 11.977831481208684, 11.977831481208684]),\n (0.98,\n 0): Atoms(symbols='Al108', pbc=True, cell=[12.060439826001351, 12.060439826001351, 12.060439826001351]),\n (1.0,\n 0): Atoms(symbols='Al108', pbc=True, cell=[12.141931756274893, 12.141931756274893, 12.141931756274893]),\n (1.02,\n 0): Atoms(symbols='Al108', pbc=True, cell=[12.222344243795412, 12.222344243795412, 12.222344243795412]),\n (1.04,\n 0): Atoms(symbols='Al108', pbc=True, cell=[12.301712339412747, 12.301712339412747, 12.301712339412747]),\n (1.06,\n 0): Atoms(symbols='Al108', pbc=True, cell=[12.38006930767338, 12.38006930767338, 12.38006930767338]),\n (1.08,\n 0): Atoms(symbols='Al108', pbc=True, cell=[12.457446749652004, 12.457446749652004, 12.457446749652004]),\n (1.1,\n 0): Atoms(symbols='Al108', pbc=True, cell=[12.533874715230777, 12.533874715230777, 12.533874715230777])}}"},"metadata":{}}],"execution_count":9},{"id":"9dcd4a1e-7122-4f57-93c1-bd9267084f70","cell_type":"markdown","source":"In contrast to the previous workflows which only used the `calc_energy` function of the simulation codes the `PhonopyWorkflow`\nand correspondingly also the `QuasiHarmonicWorkflow` require the calculation of the forces `calc_forces` in addition to\nthe calculation of the energy. Still the general steps of the workflow remain the same: ","metadata":{}},{"id":"2e96e588-e279-4d6c-8f40-2eafa982933b","cell_type":"code","source":"result_dict = evaluate_with_lammpslib(\n task_dict=task_dict,\n potential_dataframe=potential_dataframe,\n)","metadata":{"trusted":true},"outputs":[],"execution_count":10},{"id":"8fa40f79-f919-47df-ab07-c9a9dcd04b3d","cell_type":"markdown","source":"The `structure_dict` is evaluated with the [LAMMPS](https://www.lammps.org/) molecular dynamics simulation code to \ncalculate the corresponding energies and forces. The output is not plotted here as the forces for the 108 atom cells \nresult in 3x108 outputs per cell. Still the structure of the `result_dict` again follows the labels of the `structure_dict`\nas explained before. Finally, in the third step the individual free energy curves at the different temperatures are \nfitted to determine the equilibrium volume at the given temperature using the `analyse_structures()` \nand `get_thermal_properties()` functions:","metadata":{}},{"id":"371977fd-cd4e-469f-8955-17a2946c8629","cell_type":"code","source":"eng_internal_dict, phonopy_collect_dict = analyse_results_for_quasi_harmonic_approximation(\n qh_dict=qh_dict,\n output_dict=result_dict,\n dos_mesh=20,\n number_of_snapshots=None,\n)\nthermal_properties_dict_qm = get_thermal_properties_for_quasi_harmonic_approximation(\n eng_internal_dict=eng_internal_dict,\n task_dict=task_dict,\n qh_dict=qh_dict,\n fit_type=\"polynomial\",\n fit_order=3,\n temperatures=np.arange(1, 1500, 50),\n output_keys=[\"free_energy\", \"temperatures\", \"volumes\"],\n quantum_mechanical=True,\n)\ntemperatures_qh_qm, volume_qh_qm = (\n thermal_properties_dict_qm[\"temperatures\"],\n thermal_properties_dict_qm[\"volumes\"],\n)","metadata":{"trusted":true},"outputs":[],"execution_count":11},{"id":"6c7145b4-9a55-4212-a34f-a1300f7b440f","cell_type":"markdown","source":"The optimal volume at the different `temperatures` is stored in the `volume_qh_qm` in analogy to the previous section. Here the extension `_qm` indicates that the quantum-mechanical harmonic oszillator is used. ","metadata":{}},{"id":"70002fc3-2436-43ed-8a4c-0b3c1f9a3812","cell_type":"code","source":"thermal_properties_dict_cl = get_thermal_properties_for_quasi_harmonic_approximation(\n eng_internal_dict=eng_internal_dict,\n task_dict=task_dict,\n qh_dict=qh_dict,\n fit_type=\"polynomial\",\n fit_order=3,\n temperatures=np.arange(1, 1500, 50),\n output_keys=[\"free_energy\", \"temperatures\", \"volumes\"],\n quantum_mechanical=False,\n)\ntemperatures_qh_cl, volume_qh_cl = (\n thermal_properties_dict_cl[\"temperatures\"],\n thermal_properties_dict_cl[\"volumes\"],\n)","metadata":{"trusted":true},"outputs":[],"execution_count":12},{"id":"bb0db978-365f-43af-9a20-6ebb58fb8da9","cell_type":"markdown","source":"For the classical harmonic oszillator the resulting volumes are stored as `volume_qh_cl`. ","metadata":{}},{"id":"eb795fbd-0477-492a-b883-9cb31b58d3e2","cell_type":"markdown","source":"### Molecular Dynamics\nFinally, the third and most commonly used method to determine the volume expansion is using a molecular dynamics \ncalculation. While the `atomistics` package already includes a `LangevinWorkflow` at this point we use the [Nose-Hoover\nthermostat implemented in LAMMPS](https://docs.lammps.org/fix_nh.html) directly via the LAMMPS calculator interface. ","metadata":{}},{"id":"746294bd5691f89e","cell_type":"code","source":"from atomistics.calculators import calc_molecular_dynamics_thermal_expansion_with_lammpslib\n\nstructure_md = structure_opt.copy().repeat(11)\nresult_dict = calc_molecular_dynamics_thermal_expansion_with_lammpslib(\n structure=structure_md, # atomistic structure\n potential_dataframe=potential_dataframe, # interatomic potential defined as pandas.DataFrame\n Tstart=15, # temperature to for initial velocity distribution\n Tstop=1500, # final temperature\n Tstep=5, # temperature step\n Tdamp=0.1, # temperature damping of the thermostat\n run=100, # number of MD steps for each temperature\n thermo=100, # print out from the thermostat\n timestep=0.001, # time step for molecular dynamics\n Pstart=0.0, # initial pressure\n Pstop=0.0, # final pressure\n Pdamp=1.0, # barostat damping\n seed=4928459, # random seed\n dist=\"gaussian\", # Gaussian velocity distribution\n)\ntemperature_md_lst, volume_md_lst = result_dict[\"temperatures\"], result_dict[\"volumes\"]","metadata":{"trusted":true},"outputs":[{"name":"stderr","output_type":"stream","text":"100%|██████████| 298/298 [06:26<00:00, 1.30s/it]\n"}],"execution_count":13},{"id":"d2efeb52-ee54-4eb0-878a-184f353941bf","cell_type":"markdown","source":"The `calc_molecular_dynamics_thermal_expansion_with_lammpslib()` function defines a loop over a vector of temperatures in\n5K steps. For each step 100 molecular dynamics steps are executed before the temperature is again increased by 5K. For \n~280 steps with the Morse Pair Potential this takes approximately 5 minutes on a single core. These simulations can be \nfurther accelerated by adding the `cores` parameter. The increase in computational cost is on the one hand related to \nthe large number of force and energy calls and on the other hand to the size of the atomistic structure, as these \nsimulations are typically executed with >5000 atoms rather than the 4 or 108 atoms in the other approximations. The \nvolume for the individual temperatures is stored in the `volume_md_lst` list. ","metadata":{}},{"id":"eff137a4-61fc-4cbe-8c60-b0dc534a5f3f","cell_type":"markdown","source":"### Summary\nTo visually compare the thermal expansion predicted by the three different approximations, the [matplotlib](https://matplotlib.org)\nis used to plot the volume over the temperature:","metadata":{}},{"id":"da8f641d-c5e6-4c10-8aeb-c891109e2e6d","cell_type":"code","source":"import matplotlib.pyplot as plt\n\nplt.plot(\n temperature_md_lst,\n np.array(volume_md_lst) / len(structure_md) * len(structure_opt),\n label=\"Molecular Dynamics\",\n color=\"C2\",\n)\nplt.plot(temperatures_qh_qm, volume_qh_qm, label=\"Quasi-Harmonic (qm)\", color=\"C3\")\nplt.plot(temperatures_qh_cl, volume_qh_cl, label=\"Quasi-Harmonic (classic)\", color=\"C0\")\nplt.plot(temperatures_ev, volume_ev, label=\"Moruzzi Model\", color=\"C1\")\nplt.axhline(structure_opt.get_volume(), linestyle=\"--\", color=\"red\")\nplt.legend()\nplt.xlabel(\"Temperature (K)\")\nplt.ylabel(\"Volume ($\\AA^3$)\")","metadata":{"trusted":true},"outputs":[{"execution_count":14,"output_type":"execute_result","data":{"text/plain":"Text(0, 0.5, 'Volume ($\\\\AA^3$)')"},"metadata":{}},{"output_type":"display_data","data":{"text/plain":"
","image/png":"iVBORw0KGgoAAAANSUhEUgAAAkUAAAGwCAYAAACnyRH2AAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjMsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvZiW1igAAAAlwSFlzAAAPYQAAD2EBqD+naQAAs6pJREFUeJzs3Xd4jecbwPHvOSc7klhJJGSRBCEIsXeNWLW3UruqSlG7ahWl1OjQUqNVtWrWHrU3FXslRohEgkhkj/P+/njr/JoGNcLJuD/XlevyvO9znnO/aZNz55kaRVEUhBBCCCFyOa2xAxBCCCGEyAokKRJCCCGEQJIiIYQQQghAkiIhhBBCCECSIiGEEEIIQJIiIYQQQghAkiIhhBBCCABMjB1AdqHX67l79y42NjZoNBpjhyOEEEKIF6AoCo8fP8bZ2Rmt9vl9QZIUvaC7d+/i4uJi7DCEEEII8Qpu375NkSJFnltHkqIXZGNjA6jfVFtbWyNHI4QQQogXERMTg4uLi+Fz/HkkKXpBT4bMbG1tJSkSQgghspkXmfoiE62FEEIIIZCkSAghhBACkKRICCGEEAKQOUWZLi0tjZSUFGOHIUSOZ2Zm9p/La4UQ4mVIUpRJFEUhPDycR48eGTsUIXIFrVaLh4cHZmZmxg5FCJFDSFKUSZ4kRA4ODlhZWckGj0K8QU82Uw0LC8PV1VV+3oQQmUKSokyQlpZmSIgKFChg7HCEyBXs7e25e/cuqampmJqaGjscIUQOIAPymeDJHCIrKysjRyJE7vFk2CwtLc3IkQghcgpJijKRdOEL8fbIz5sQIrNJUiSEEEIIgSRFQgghhBCAJEUik+3duxeNRvNWtyYYP3485cqVe2vvlx1pNBrWr19v7DCEECJLk6Qol+vevTsajYZ+/fpluNe/f380Gg3du3d/+4FlEXXq1EGj0aDRaDA3N6dw4cK8++67rF271tihvZSwsDAaN25s7DCEEDlJWipc3W7sKDKVJEUCFxcXVqxYQUJCguFaYmIiy5cvx9XV1YiRvT3P24W8T58+hIWFERQUxJo1a/Dx8aFjx4707dv3LUb4egoVKoS5ubmxwxBC5BQpCbCqK/zWHk4tMXY0mUaSojdEURTiU+KN8qUoykvFWr58eVxdXdP1fqxduxYXFxf8/PzS1U1KSmLgwIE4ODhgYWFBjRo1OHHixHPbP3z4MLVq1cLS0hIXFxcGDhxIXFxcujaHDx+Oi4sL5ubmeHl5sXDhQgCWLFlC3rx507W3fv365648OnHiBA0aNKBgwYLY2dlRu3Zt/vrrr3R1NBoNP/zwAy1atMDa2povvvjime1ZWVlRqFAhXFxcqFKlCtOmTePHH39kwYIF7Nq1C4B33nmHAQMGpHvdgwcPMDc3588//wTA3d2dKVOm0LNnT2xsbHB1dWX+/PnpXjNixAi8vb2xsrKiaNGijB07Nl3C9mSocNGiRbi6upInTx4+/PBD0tLSmD59OoUKFcLBwYHJkydneN5/Dp/duXOHjh07kj9/fqytrfH39+fYsWMAnDlzhrp162JjY4OtrS0VKlTg5MmTz/z+CCFymYRHsLQ1XNkCOnOwKmjsiDKNbN74hiSkJlD5t8pGee9jnY9hZfpyeyb16NGDxYsX06VLFwAWLVpEz5492bt3b7p6w4cPZ82aNfz888+4ubkxffp0AgICCAoKIn/+/BnaPXfuHAEBAUyaNImFCxcSGRnJgAEDGDBgAIsXLwagW7duHDlyhLlz51K2bFlu3LjB/fv3X+3hgcePH/P+++8zd+5cAGbOnEmTJk24du0aNjY2hnrjxo1j6tSpzJo1C51O91Lv8f777zN06FDWrl1L/fr16d27NwMGDGDmzJmGHplly5bh7OxM3bp1Da+bOXMmkyZNYvTo0fz+++98+OGH1KpVixIlSgBgY2PDkiVLcHZ25ty5c/Tp0wcbGxuGDx9uaCM4OJitW7eybds2goODadu2LTdu3MDb25t9+/Zx+PBhevbsSb169ahSpUqG2GNjY6lduzaFCxdm48aNFCpUiL/++gu9Xg9Aly5d8PPzY968eeh0OgIDA2VzRCGEKiYMfm0DERfA3BY6rQD36saOKtNIUiQA6Nq1K6NGjeLmzZtoNBoOHTrEihUr0iVFcXFxzJs3jyVLlhjmpyxYsICdO3eycOFChg0blqHdr776is6dO/PJJ58A4OXlxdy5c6lduzbz5s0jJCSEVatWsXPnTurXrw9A0aJFX+tZ3nnnnXTlH3/8kXz58rFv3z6aNWtmuN65c2d69uz5Su+h1Wrx9vbm5s2bALRp04aPP/6YDRs20L59ewAWL15smLP1RJMmTejfvz+g9grNmjWLvXv3GpKizz77zFDX3d2doUOHsnLlynRJkV6vZ9GiRdjY2ODj40PdunW5cuUKW7ZsQavVUrx4caZNm8bevXufmhT99ttvREZGcuLECUMi6+npabgfEhLCsGHDDDF5eXm90vdICJHDPAiGpS3hUQjkcYT31kAhX2NHlakkKXpDLE0sOdb5mNHe+2UVLFiQpk2b8vPPP6MoCk2bNqVgwfRdosHBwaSkpFC9+v//KjA1NaVSpUpcunTpqe2eOnWKoKAgli1bZrimKAp6vZ4bN25w7tw5dDodtWvXfumYnyUiIoLPP/+cP//8k3v37pGWlkZ8fDwhISHp6vn7+7/W+yiKYkh4zM3Nee+991i0aBHt27cnMDCQM2fOZFjxVaZMGcO/NRoNhQoVIiIiwnDt999/Z/bs2QQFBREbG0tqaiq2trbp2nB3d0/X4+Xo6IhOp0t3Yryjo2O6dv8pMDAQPz+/p/bsAQwZMoTevXuzdOlS6tevT7t27ShWrNiLfVOEEDnT3dPwa1uIvw/5i8J7ayG/h7GjynSSFL0hGo3mpYewjK1nz56GeTHfffddhvtP5ir9ez7PP5ODf9Pr9XzwwQcMHDgwwz1XV1eCgoKeG5NWq80wR+p5k6JBXVEXGRnJ7NmzcXNzw9zcnKpVq5KcnJyunrW19XPbeZ60tDSuXbtGxYoVDdd69+5NuXLluHPnDosWLaJevXq4ubmle92/h6E0Go1h2Oro0aN07NiRCRMmEBAQgJ2dHStWrGDmzJn/2cbz2v03S8vnJ83jx4+nc+fObN68ma1btzJu3DhWrFhBq1atnvs6IUQOdX0vrOgCybFQqIzaQ5THwdhRvREy0VoYNGrUiOTkZJKTkwkICMhw39PTEzMzMw4ePGi4lpKSwsmTJylZsuRT2yxfvjwXLlzA09Mzw5eZmRm+vr7o9Xr27dv31Nfb29vz+PHjdBOzAwMDn/scBw4cYODAgTRp0oRSpUphbm7+WnOUnubnn38mKiqKNm3aGK75+vri7+/PggUL+O233156aO7QoUO4ubkxZswY/P398fLy4tatW5kaN6i9VYGBgTx8+PCZdby9vRk8eDA7duygdevWhvlfQohc5sI6WNZOTYg8akH3zTk2IQJJisQ/6HQ6Ll26xKVLl5468dja2poPP/yQYcOGsW3bNi5evEifPn2Ij4+nV69eT21zxIgRHDlyhI8++ojAwECuXbvGxo0b+fjjjwF1KOj999+nZ8+erF+/nhs3brB3715WrVoFQOXKlbGysmL06NEEBQXx22+/sWTJkuc+h6enJ0uXLuXSpUscO3aMLl26/GfvyPPEx8cTHh7OnTt3OHbsGCNGjKBfv358+OGH6SZRg9pb9OWXX5KWlvbSPSuenp6EhISwYsUKgoODmTt3LuvWrXvluJ+lU6dOFCpUiJYtW3Lo0CGuX7/OmjVrOHLkCAkJCQwYMIC9e/dy69YtDh06xIkTJ56Z9AohcrDjC2B1D0hLBp8W0OV3sLD979dlY5IUiXRsbW0zzGH5py+//JI2bdrQtWtXypcvT1BQENu3bydfvnxPrV+mTBn27dvHtWvXqFmzJn5+fowdOxYnJydDnXnz5tG2bVv69+9PiRIl6NOnj6FnKH/+/Pz6669s2bIFX19fli9fzvjx45/7DIsWLSIqKgo/Pz+6du1q2ELgVS1YsAAnJyeKFStGq1atuHjxIitXruT777/PULdTp06YmJjQuXNnLCwsXup9WrRoweDBgxkwYADlypXj8OHDjB079pXjfhYzMzN27NiBg4MDTZo0wdfXly+//BKdTodOp+PBgwd069YNb29v2rdvT+PGjZkwYUKmxyGEyKIUBfZMhS2fAgr494S2i8Ek5+91plFedlObXComJgY7Ozuio6MzJA2JiYncuHEDDw+Pl/4gFDnL7du3cXd358SJE5QvX97Y4eRo8nMnxBugT4Mtw+CkulcctUdCnZHwnL3hsrrnfX7/m0y0FiITpKSkEBYWxsiRI6lSpYokREKI7Cc1Cdb2gYsbAA00+Qoq9TF2VG+VJEVCZIJDhw5Rt25dvL29+f33340djhBCvJzEGFjZBW7sB60ptJ4PpVsbO6q3TpIiITJBnTp1Xvp4FSGEyBJiI2BZWwg7A2Z5oOMyKFrH2FEZhSRFQgghRG4VdROWtoKH19UzzN77HZz9/vNlOZUkRUIIIURudDdQ3YMoLgLsXKHrOijo+Z8vy8kkKRJCCCFym+A/YWVXdVNGx9LqHkS2Tv/9uhwuS+5TFBoaynvvvUeBAgWwsrKiXLlynDp1ynD/3r17dO/eHWdnZ6ysrGjUqBHXrl17bptLlixBo9Fk+EpMTHzTjyOEEEJkHWdW/n+Xavea0GOLJER/y3I9RVFRUVSvXp26deuydetWHBwcCA4OJm/evIB6zlbLli0xNTVlw4YN2Nra8vXXX1O/fn0uXrz43POsbG1tuXLlSrprsr+JEEKIXEFR4NAc2DVOLZduAy3n5YpNGV9UluspmjZtGi4uLixevJhKlSrh7u5OvXr1DKd0X7t2jaNHjzJv3jwqVqxI8eLF+f7774mNjWX58uXPbfvJqeT//BJZx5IlSwzJb26UWc9fq1Ytfvvtt9cP6D98+umnTz3oVwiRBen1sG3k/xOiqgOg9U+SEP1LlkuKNm7ciL+/P+3atcPBwQE/Pz8WLFhguJ+UlASk7+HR6XQZDip9mtjYWNzc3ChSpAjNmjXj9OnTz6yblJRETExMuq+c6vbt2/Tq1QtnZ2fMzMxwc3Nj0KBBPHjw4K3G0aFDB65evfrcOuPHj6dcuXIZrt+8eRONRvOfh8VmZS/y/P9l06ZNhIeH07Fjx0yK6tmGDx/O4sWLuXHjxht/LyHEa0hJhN97wLEf1HLDLyBgMmizXApgdFnuO3L9+nXmzZuHl5cX27dvp1+/fgwcOJBffvkFgBIlSuDm5saoUaOIiooiOTmZL7/8kvDwcMLCwp7ZbokSJViyZAkbN25k+fLlWFhYUL169WfORZo6dSp2dnaGLxcXlzfyvMZ2/fp1/P39uXr1KsuXLycoKIgffviB3bt3U7Vq1eeepJ7ZLC0tX+uMsteVnJxstPeGzHn+uXPn0qNHD7Rv4Zedg4MDDRs25Icffnjj7yWEeEUJj+DXNnBx/d+bMv4E1T42dlRZl5LFmJqaKlWrVk137eOPP1aqVKliKJ88eVIpW7asAig6nU4JCAhQGjdurDRu3PiF3yctLU0pW7as8vHHHz/1fmJiohIdHW34un37tgIo0dHRGeomJCQoFy9eVBISEl74/bOKRo0aKUWKFFHi4+PTXQ8LC1OsrKyUfv36Ga4Byrp169LVs7OzUxYvXmwoDx8+XPHy8lIsLS0VDw8P5bPPPlOSk5MN9wMDA5U6deooefLkUWxsbJTy5csrJ06cUBRFURYvXqzY2dk9N95x48YpZcuWzXD9xo0bCqCcPn1aURRFSU1NVXr27Km4u7srFhYWire3tzJ79ux0r3n//feVFi1aKFOmTFGcnJwUNzc3QzsrV65UatSooVhYWCj+/v7KlStXlOPHjysVKlRQrK2tlYCAACUiIsLQVlpamjJhwgSlcOHCipmZmVK2bFll69atGeJbs2aNUqdOHcXS0lIpU6aMcvjwYUOdpz3/hg0blAoVKijm5uZKgQIFlFatWj3zexMZGaloNBrl/Pnz6a5fvXpVqVmzpmJubq6ULFlS2bFjR7r/lq/6zIqiKEuWLFFcXFyeGdOblJ1/7oR4K6JDFeW7KooyzlZRJhdWlOA9xo7IKKKjo5/5+f1vWW6itZOTEz4+PumulSxZkjVr1hjKFSpUIDAwkOjoaJKTk7G3t6dy5cr4+/u/8PtotVoqVqz4zJ4ic3NzzM1ffaxVURSUhIRXfv3r0FhaonmBw/sePnzI9u3bmTx5MpaWlunuFSpUiC5duhhOg3+R9gBsbGxYsmQJzs7OnDt3jj59+mBjY8Pw4cMB6NKlC35+fsybNw+dTkdgYCCmpqYv/5D/Qa/XU6RIEVatWkXBggU5fPgwffv2xcnJifbt2xvq7d69G1tbW3bu3JluR+px48Yxe/ZsXF1d6dmzJ506dcLW1pY5c+ZgZWVF+/bt+fzzz5k3bx4Ac+bMYebMmfz444/4+fmxaNEimjdvzoULF/Dy8jK0O2bMGGbMmIGXlxdjxoyhU6dOBAUFYWKS8Udx8+bNtG7dmjFjxrB06VKSk5PZvHnzM5/54MGDWFlZUbJkyXTfh9atW1OwYEGOHj1KTEwMn3zyyVNf/7LPDFCpUiVu377NrVu3cHNz++//MEKItyPistpDFHMH8jiqS+6dyhg7qiwvyyVF1atXz7BC7OrVq0/9hWtnZweok69PnjzJpEmTXvh9FEUhMDAQX1/f1wv4We0nJHClfIU30vZ/Kf7XKTRWVv9Z79q1ayiKku5D9J9KlixJVFQUkZGRLzys89lnnxn+7e7uztChQ1m5cqUhKQoJCWHYsGGUKFECIF3C8KLOnTtHnjx50l1T/nXEhqmpKRMmTDCUPTw8OHz4MKtWrUqXFFlbW/PTTz9hZmYGqHOTQJ1EHBAQAMCgQYPo1KkTu3fvpnr16gD06tWLJUuWGNqZMWMGI0aMMMzlmTZtGnv27GH27Nl89913hnqffvopTZs2BWDChAmUKlWKoKAgw/fjnyZPnkzHjh3TPUfZsmWf+X25efMmjo6O6YbOdu3axaVLl7h58yZFihQBYMqUKTRu3DjD61/2mQEKFy5seG9JioTIIm4dgeUdIDEaCnjBe2sgn/x8vogsN6do8ODBHD16lClTphAUFMRvv/3G/Pnz+eijjwx1Vq9ezd69e7l+/TobNmygQYMGtGzZkoYNGxrqdOvWjVGjRhnKEyZMYPv27Vy/fp3AwEB69epFYGAg/fr1e6vPl508STSeJAwv4vfff6dGjRoUKlSIPHnyMHbsWEJCQgz3hwwZQu/evalfvz5ffvklwcHBT20nJCSEPHnyGL6mTJliuFe8eHECAwPTfW3ZsiVDGz/88AP+/v7Y29uTJ08eFixYkC4WAF9f36c+X5ky//+LytHR0VD3n9ciIiIAiImJ4e7du4bk4Ynq1atz6dKlZ7br5KTuC/KknX8LDAykXr16T733NAkJCRm2mLh06RKurq6GhAigatWqT339yzzzE096GOPj4184TiHEG3RxI/zSQk2IilSEntslIXoJWa6nqGLFiqxbt45Ro0YxceJEPDw8mD17Nl26dDHUCQsLY8iQIdy7dw8nJye6devG2LFj07UTEhKS7i/mR48e0bdvX8LDw7Gzs8PPz4/9+/dTqVKlN/IcGktLiv916r8rvqH3fhGenp5oNBouXrxIy5YtM9y/fPky9vb2hmXiGo0mQ49MSkqK4d9Hjx419GwEBARgZ2fHihUrmDlzpqHO+PHj6dy5M5s3b2br1q2MGzeOFStW0KpVq3TtOjs7p1tJlj9/fsO/zczM8PRMvxX9v4efVq1axeDBg5k5cyZVq1bFxsaGr776imPHjqWr96x9rf45pPdk6PDf1/R6fbrX/HuIUVGUDNee1u6/23ni30Oa/6VgwYJERUVliOHfnjUU+irP/GQivr29/UvFKoR4A44vgC3DAAW8G0PbRWD236MG4v+yXFIE0KxZM5o1a/bM+wMHDvzP/VH27t2brjxr1ixmzZqVGeG9EI1G80JDWMZUoEABGjRowPfff8/gwYPTfQiHh4ezbNmydD109vb26Vb4Xbt2LV0PwaFDh3Bzc2PMmDGGa7du3crwvt7e3nh7ezN48GA6derE4sWLMyRFJiYmGRKfl3HgwAGqVatG//79Ddee1Sv1umxtbXF2dubgwYPUqlXLcP3w4cOvlXSXKVOG3bt306NHjxeq7+fnR3h4OFFRUeTLlw8AHx8fQkJCuHv3Ls7OzgAcOXLklWP6t/Pnz2NqakqpUqUyrU0hxEtSFPhzEhz4+w/Q8u9D069BlyU/4rO0LDd8Jt6ub7/9lqSkJAICAti/fz+3b99m27ZtNGjQAG9vbz7//HND3XfeeYdvv/2Wv/76i5MnT9KvX790PQmenp6EhISwYsUKgoODmTt3LuvWrTPcT0hIYMCAAezdu5dbt25x6NAhTpw48cw5Ta/D09OTkydPsn37dq5evcrYsWM5ceJEpr/PE8OGDWPatGmsXLmSK1euMHLkSAIDAxk0aNArtzlu3DiWL1/OuHHjuHTpEufOnWP69OnPrO/n54e9vT2HDh0yXKtfvz7FixenW7dunDlzhgMHDqRLWl/XgQMHqFmz5kv3agkhMklqMqz/8P8JUZ1R8O4cSYhekSRFuZyXlxcnTpygaNGitG/fHjc3Nxo3boy3tzeHDh1KN6F55syZuLi4UKtWLTp37synn36K1T96w1q0aMHgwYMZMGAA5cqV4/Dhw+mGNXU6HQ8ePKBbt254e3vTvn17GjdunG4icWbp168frVu3pkOHDlSuXJkHDx6k6zXKbAMHDmTo0KEMHToUX19ftm3bxsaNG19pIvkTderUYfXq1WzcuJFy5crxzjvvZBj++yedTkfPnj1ZtmyZ4ZpWq2XdunUkJSVRqVIlevfuzeTJk185pn9bvnw5ffr0ybT2hBAvITEafmsHZ5aDRqcmQ3VGwguuFhYZaZSnTToQGcTExGBnZ0d0dDS2trbp7iUmJnLjxg08PDxyxFlq48aN4+uvv2bHjh3PnJQrsqZ79+5RqlQpTp069dzVYBqNhnXr1j11LtmL2rx5M8OGDePs2bNP3VLgTctpP3dCvJToUPVQ14gLYGoN7ZaAd8P/fFlu9LzP73+T/jWRwYQJE3B3d+fYsWNUrlz5reyOLDKHo6MjCxcuJCQk5I0vkY+Li2Px4sVGSYiEyNXCz6sJ0eO76h5EnVeBczljR5UjyG8z8VQvOrlXZD0tWrR4K+/zz/2ehBBvSfAeWNUNkmKgYHHoslqW3GciSYqEyIVk1FyIbChwOWwcAPpUcKsOHZeBZT5jR5WjSFIkhBBCZGWKAvtnwJ4v1HLpNtByHpi8+lFU4ukkKRJCCCGyqrQU2DwE/vpFLVf/BOqNA5nr+UZIUiSEEEJkRUmPYXV3CNoFGi00ng6VZAuMN0mSIiGEECKreRyurjALPwsmluqRHSWaGDuqHE+SIiGEECIribgMy9pC9G2wKqguuS9SwdhR5QqSFAkhhBBZxc2DsKKzult1AU/o8jvk9zB2VLmGJEVCCCFEVnDud/Ucs7RkcKkMnVaAVX5jR5WryPR1kWUsWbKEvHnzGjsMo8ms569Vqxa//fbbW33PF9G9e/cXPlYkIiICe3t7QkND32xQQmQFiqIe6Lqml5oQlWwO3TZIQmQEkhQJbt++Ta9evXB2dsbMzAw3NzcGDRrEgwcP3mocHTp04OrVq8+tM378eMqVK5fh+s2bN9FoNAQGBr6Z4N6CF3n+/7Jp0ybCw8Pp2LFjJkWVeebMmcOSJUteqK6DgwNdu3Zl3LhxbzYoIYwtLQU2fgy7J6rlKh9Bu5/B1NK4ceVSkhTlctevX8ff35+rV6+yfPlygoKC+OGHH9i9ezdVq1bl4cOHby0WS0tLHBwc3tr7/VtycrLR3hsy5/nnzp1Ljx49suR5dXZ2di/VK9WjRw+WLVtGVFTUmwtKCGNKjFYnVJ9eqi65bzIDGk2RPYiMSL7zb4iiKMQnpxrl62WOcPjoo48wMzNjx44d1K5dG1dXVxo3bsyuXbsIDQ1lzJgxhroajYb169ene33evHnT/fU/YsQIvL29sbKyomjRoowdO5aUlBTD/TNnzlC3bl1sbGywtbWlQoUKnDx5EsjcoZy0tDR69eqFh4cHlpaWFC9enDlz5qSr82Q4Z+rUqTg7O+Pt7W3ocVq1ahU1a9bE0tKSihUrcvXqVU6cOIG/vz958uShUaNGREZGGtrS6/VMnDiRIkWKYG5uTrly5di2bZvh/pN2165dS926dbGysqJs2bIcOXLEUOdpz79x40b8/f2xsLCgYMGCtG7d+pnPfP/+fXbt2kXz5s3TXX/06BF9+/bF0dERCwsLSpcuzaZNm57aRnBwMC1atMDR0ZE8efJQsWJFdu3ala7O999/j5eXFxYWFjg6OtK2bVvDvd9//x1fX18sLS0pUKAA9evXJy4uLt33+5/fs2nTpuHp6Ym5uTmurq5MnjzZcN/X15dChQqxbt26Zz6zENnWoxBYGADX96qn3HdaIXsQZQEy0foNSUhJw+fz7UZ574sTA7Ay++//tA8fPmT79u1MnjwZS8v0XbWFChWiS5curFy5ku+//x6NRvNC721jY8OSJUtwdnbm3Llz9OnTBxsbG4YPHw5Aly5d8PPzY968eeh0OgIDAzE1NX35h/wPer2eIkWKsGrVKgoWLMjhw4fp27cvTk5O6Q4y3b17N7a2tuzcuTNdMjlu3Dhmz56Nq6srPXv2pFOnTtja2jJnzhysrKxo3749n3/+OfPmzQPUoaGZM2fy448/4ufnx6JFi2jevDkXLlzAy8vL0O6YMWOYMWMGXl5ejBkzhk6dOhEUFPTUk+Y3b95M69atGTNmDEuXLiU5OZnNmzc/85kPHjyIlZUVJUuWTPd9aNy4MY8fP+bXX3+lWLFiXLx4EZ1O99Q2YmNjadKkCV988QUWFhb8/PPPvPvuu1y5cgVXV1dOnjzJwIEDWbp0KdWqVePhw4ccOHAAgLCwMDp16sT06dNp1aoVjx8/5sCBA89M0keNGsWCBQuYNWsWNWrUICwsjMuXL6erU6lSJQ4cOEDPnj2f+dxCZDuhf8HyjhB7D2ycoPNKcCpr7KgEkhTlateuXUNRlHQfov9UsmRJoqKiiIyMfOFhnc8++8zwb3d3d4YOHcrKlSsNSVFISAjDhg2jRIkSAOkShhd17tw58uTJk+7avz94TU1NmTBhgqHs4eHB4cOHWbVqVbqkyNramp9++gkzMzNA7dEB+PTTTwkICABg0KBBdOrUid27d1O9enUAevXqla6HbMaMGYwYMcIwl2fatGns2bOH2bNn89133xnqffrppzRt2hSACRMmUKpUKYKCggzfj3+aPHkyHTt2TPccZcs++xfnzZs3cXR0TDd0tmvXLo4fP86lS5fw9vYGoGjRos9so2zZsune44svvmDdunVs3LiRAQMGEBISgrW1Nc2aNcPGxgY3Nzf8/PwANSlKTU2ldevWuLmpp3b7+vo+9X0eP37MnDlz+Pbbb3n//fcBKFasGDVq1EhXr3Dhwpw+ffqZ8QqR7VzeDL/3gtQEcCyt7kFkV9jYUYm/SVL0hlia6rg4McBo750ZniQaTxKGF/H7778ze/ZsgoKCiI2NJTU1FVtbW8P9IUOG0Lt3b5YuXUr9+vVp164dxYoVy9BOSEgIPj4+hvLo0aMZPXo0AMWLF2fjxo3p6oeGhlKnTp1013744Qd++uknbt26RUJCAsnJyRkmafv6+j71+cqUKWP4t6Ojo6HuP69FREQAEBMTw927dw0J0xPVq1fnzJkzz2zXyckJUFdaPS0pCgwMpE+fF+9OT0hIwMLCIkMbRYoUMSRE/yUuLo4JEyawadMm7t69S2pqKgkJCYSEhADQoEED3NzcKFq0KI0aNaJRo0a0atXKMBxYr149fH19CQgIoGHDhrRt25Z8+TKe4n3p0iWSkpKoV6/ec+OxtLQkPj7+Bb8DQmRhigJH58H20YACnvWh3RIwtzF2ZOIfZE7RG6LRaLAyMzHK14sOdXl6eqLRaLh48eJT71++fBl7e3vDPBeNRpOhR+af84WOHj1Kx44dady4MZs2beL06dOMGTMm3QTm8ePHc+HCBZo2bcqff/6Jj4/PU+eMODs7ExgYaPjq16+f4Z6ZmRmenp7pvp70TDyxatUqBg8eTM+ePdmxYweBgYH06NEjw2Rqa2vrpz77P4f0nnw//31Nr9ene82/v++KomS49rR2/93OE/8e0vwvBQsWzDAp+WXbGDZsGGvWrGHy5MkcOHCAwMBAfH19Dd83Gxsb/vrrL5YvX46TkxOff/45ZcuW5dGjR+h0Onbu3MnWrVvx8fHhm2++oXjx4ty4ceOVn+3hw4fY29u/1DMIkeWkpcLW4bB9FKCAf0/otFISoixIkqJcrECBAjRo0IDvv/+ehISEdPfCw8NZtmwZ3bt3N1yzt7cnLCzMUL527Vq6v+IPHTqEm5sbY8aMwd/fHy8vL27dupXhfb29vRk8eDA7duygdevWLF68OEMdExOTdElP/vwvt1/HgQMHqFatGv3798fPzw9PT0+Cg4Nfqo0XZWtri7OzMwcPHkx3/fDhw88cmnwRZcqUYffu3S9c38/Pj/Dw8HSJUZkyZbhz584LL/U/cOAA3bt3p1WrVoaJzk+GFJ8wMTGhfv36TJ8+nbNnz3Lz5k3+/PNPQE30qlevzoQJEzh9+jRmZmZPTXq9vLywtLT8z+c7f/68YXhOiGwpKRZWdoHj89Vyg0nQ9GvQyUBNViT/VXK5b7/9lmrVqhEQEMAXX3yBh4cHFy5cYNiwYXh7e/P5558b6r7zzjt8++23VKlSBb1ez4gRI9L1fHh6ehISEsKKFSuoWLEimzdvTveBmJCQwLBhw2jbti0eHh7cuXOHEydO0KZNm0x/Lk9PT3755Re2b9+Oh4cHS5cu5cSJE3h4vJnt8ocNG8a4ceMoVqwY5cqVY/HixQQGBrJs2bJXbnPcuHHUq1ePYsWK0bFjR1JTU9m6dathfta/+fn5YW9vz6FDh2jWrBkAtWvXplatWrRp04avv/4aT09PLl++jEajoVGjRhna8PT0ZO3atbz77rtoNBrGjh2bridr06ZNXL9+nVq1apEvXz62bNmCXq+nePHiHDt2jN27d9OwYUMcHBw4duwYkZGRT00MLSwsGDFiBMOHD8fMzIzq1asTGRnJhQsX6NWrFwDx8fGcOnWKKVOmvPL3UAijigmD39r/fairBbSeDz4tjB2VeA7pKcrlvLy8OHHiBEWLFqV9+/a4ubnRuHFjvL29OXToULoJzTNnzsTFxYVatWrRuXNnPv30U6ysrAz3W7RoweDBgxkwYADlypXj8OHDjB071nBfp9Px4MEDunXrhre3N+3bt6dx48bpJhJnln79+tG6dWs6dOhA5cqVefDgAf3798/093li4MCBDB06lKFDh+Lr68u2bdvYuHHjK00kf6JOnTqsXr2ajRs3Uq5cOd555x2OHTv2zPo6nY6ePXtmSMTWrFlDxYoV6dSpEz4+PgwfPpy0tLSntjFr1izy5ctHtWrVePfddwkICKB8+fKG+3nz5mXt2rW88847lCxZkh9++IHly5dTqlQpbG1t2b9/P02aNMHb25vPPvuMmTNn0rhx46e+19ixYxk6dCiff/45JUuWpEOHDoZ5WgAbNmzA1dWVmjVrvsy3TYisIfw8/FRPTYisCsL7myQhygY0ystsapOLxcTEYGdnR3R0dLqJwwCJiYncuHEDDw+PDBNds6Nx48bx9ddfs2PHDqpWrWrscMRLuHfvHqVKleLUqVMZ5lllN5UqVeKTTz6hc+fOT72f037uRA4StAtWdYfkx1DQG7qshnzuxo4q13re5/e/SU+RyGDChAnMnTuXY8eOPXMSsMiaHB0dWbhwoWG1WHYVERFB27Zt6dSpk7FDEeLlnFwMy9qrCZF7Tei1QxKibER6il5QbuopEiI7kJ87kaXo02Dn53DkW7VctjO8OwdMXnxLE/FmvExPkUy0FkIIIV5HUiys7QNXtqjlumOg1jB4we1RRNYhSZEQQgjxqqJDYXkHCD8HOnNoNQ9KZ/6KWvF2SFIkhBBCvIq7p+G3jhAbDtb20HE5uFQ0dlTiNUhSJIQQQrysixthbV/1DDMHH/WU+3zZe8WnkKRICCGEeHGKAgdnwe6/91fzbABtF4HF8yfwiuxBkiIhhBDiRaQmw6bBEPirWq70AQRMkSM7chDZp0jkCBqNhvXr1xs1hr1796LRaHj06NELv8bd3Z3Zs2e/sZiEEJkk/iEsbaUmRBotNJkBTaZLQpTDSFKUy3Xv3h2NRpPuFPon+vfvj0ajSXcobFYVFhb2zOMkIOc8pxDCCO4HwU/14dZBMLOBzquhUh9jRyXeAEmKBC4uLqxYsYKEhATDtcTERJYvX46rq+trta0oCqmpqa8b4n8qVKgQ5ubmz63zJp9TCJFD3divnmH2MBjsXNUdqr3qGzsq8YZIUiQoX748rq6urF271nBt7dq1uLi44Ofnl65uUlISAwcOxMHBAQsLC2rUqMGJEycM958MIW3fvh1/f3/Mzc05cOAA3bt3p2XLluna+uSTT6hTpw4AN2/eRKPRZPh6cr9OnTpPvX/z5k3gxYbPMvM5AbZs2YK3tzeWlpbUrVvXEMs/HT58mFq1amFpaYmLiwsDBw4kLi7uuXEKIbKIv5aqQ2aJj6BIReizGxx9jB1VjqZXjHu0lCRFb4qiQHKccb5e4eSWHj16sHjxYkN50aJF9OzZM0O94cOHs2bNGn7++Wf++usvPD09CQgI4OHDhxnqTZ06lUuXLlGmTJn/fH8XFxfCwsIMX6dPn6ZAgQLUqlULUJOXf95v3bo1xYsXx9HR0SjPefv2bVq3bk2TJk0IDAykd+/ejBw5Ml0b586dIyAggNatW3P27FlWrlzJwYMHGTBgwEvFLIR4y/R69ciOjQNAn6puxvj+H5DHwdiR5Rh6RU9obCh7b+9lwdkFDN83nFYbWtF8fXOjxiUzxN6UlHiY4myc9x59F8ysX+olXbt2ZdSoUYYem0OHDrFixQr27t1rqBMXF8e8efNYsmSJYf7OggUL2LlzJwsXLmTYsGGGuhMnTqRBgwYv/P46nY5ChQoB6pBWy5YtqVq1KuPHjwcgf/78hrqzZs3izz//5NixY1haWhrlOefNm0fRokWZNWsWGo2G4sWLc+7cOaZNm2Zo56uvvqJz58588sknAHh5eTF37lxq167NvHnz5LwuIbKipFhY9wFc3qSWa4+EOiPlyI7X8CjxEdceXeNq1FWuRV3j2qNrBEUFEZ8a/9T6cSlxWJu+3GdYZsmSSVFoaCgjRoxg69atJCQk4O3tzcKFC6lQoQIA9+7dY8SIEezYsYNHjx5Rq1YtvvnmG7y8vJ7b7po1axg7dizBwcEUK1aMyZMn06pVq7fxSFlewYIFadq0KT///DOKotC0aVMKFiyYrk5wcDApKSlUr17dcM3U1JRKlSpx6dKldHX9/f1fOZZevXrx+PFjdu7ciVabvjNz69atjBw5kj/++ANvb++XbjuznvPSpUtUqVIFzT9+UVatWjVdO6dOnSIoKIhly5YZrimKgl6v58aNG5QsWfKl4xdCvEGPbsPyTnDvHOjMoMX3UKadsaPKlqISo2i/qT3hceHPreds7Yyfox9eeb3wyueFV14vrEys3lKUGWW5pCgqKorq1atTt25dtm7dioODA8HBweTNmxdQP1RatmyJqakpGzZswNbWlq+//pr69etz8eJFrK2fnl0eOXKEDh06MGnSJFq1asW6deto3749Bw8epHLlypn/IKZWao+NMZi+2v9QPXv2NAztfPfddxnuK38Py2n+9ReToigZrv37v4NWqzW8/omUlJQM7/HFF1+wbds2jh8/jo2NTbp7Fy9epGPHjnz55Zc0bNjwBZ8qo8x4zn8/y9Po9Xo++OADBg4cmOGeTOwWIou5fRxWdIG4CPXIjg7LwPUNfDbkEpcfXv7PhAigl28v2hdv/xYiejFZLimaNm0aLi4u6eZ9uLu7G/597do1jh49yvnz5ylVqhQA33//PQ4ODixfvpzevXs/td3Zs2fToEEDRo0aBcCoUaPYt28fs2fPZvny5RnqJyUlkZSUZCjHxMS83INoNC89hGVsjRo1Ijk5GYCAgIAM9z09PTEzM+PgwYN07twZUBObkydPGoaInsXe3p7z58+nuxYYGIipqamhvGbNGiZOnMjWrVspVqxYuroPHjzg3XffpXXr1gwePPhVHs8gM57Tx8cnw8Tuo0ePpiuXL1+eCxcu4Onp+VrxCiHesDMrYePHkJYEjr7Q6TfIK3+4vI6KhSriZO1EWFzYM+uUsy9HQ7dX/wP3TchySdHGjRsJCAigXbt27Nu3j8KFC9O/f3/69FH3hHiSqPxzPoZOpzN8iD0rKTpy5EiGD9OAgIBnbpw3depUJkyYkAlPlH3odDrD8JBOp8tw39ramg8//JBhw4aRP39+XF1dmT59OvHx8fTq1eu5bb/zzjt89dVX/PLLL1StWpVff/2V8+fPG1Z9nT9/nm7dujFixAhKlSpFeLj6F4aZmRn58+endevWWFpaMn78eMM9UJOtp8X6pp+zX79+zJw5kyFDhvDBBx9w6tQplixZkq6dESNGUKVKFT766CP69OmDtbU1ly5dYufOnXzzzTcvFbMQ4g3Q6+HPSXDwa7Vcohm0+hHM8xg3rmwqPiWeU/dOcTTsKIdCDxkSop6le1K5UGVcbF1wsnbCRJvlUg+DLBfZ9evXmTdvHkOGDGH06NEcP36cgQMHYm5uTrdu3ShRogRubm6MGjWKH3/8EWtra77++mvCw8MJC3t2RhoeHp5hpZKjo2O6D9h/GjVqFEOGDDGUY2JicHFxyZyHzMJsbZ9/fs+XX36JXq+na9euPH78GH9/f7Zv306+fPme+7qAgADGjh3L8OHDSUxMpGfPnnTr1o1z584BcPLkSeLj4/niiy/44osvDK+rXbs2e/fuZf/+/UD6XkOAGzduZLj2Np7T1dWVNWvWMHjwYL7//nsqVarElClT0q1kK1OmDPv27WPMmDHUrFkTRVEoVqwYHTp0eOl4hRCZ7N8TqmsMgXfGglYWZb+Mfbf3serqKvbf2Z/hnpnWjInVJ9K0aFMjRPZqNMqLTI54i8zMzPD39+fw4cOGawMHDuTEiRMcOXIEUCew9urVizNnzqDT6ahfv75hQu6WLVue2e7PP/9Mp06dDNeWLVtGr169SExM/M+4YmJisLOzIzo6OsMHamJiIjdu3MDDw0NWFAnxlsjPnXhl6SZUm0Pzb6Cs/LHysrYHb2X1T59yrLjmqavz5tWfR43CNYwQWXrP+/z+tyyXEjs5OeHjk35zrJIlSxISEmIoV6hQgcDAQB49ekRYWBjbtm3jwYMHeHh4PLPdQoUKZegVioiIeOl9boQQQmRjt4/DgrpqQmTtAN03S0L0CpLvhOI84nuGrtNT83zGvpXCeQpT1r6sESJ7PVlu+Kx69epcuXIl3bWrV6/i5uaWoa6dnR2gTr4+efIkkyZNema7VatWZefOnenmFe3YsYNq1aplUuRCCCGytDMr/p5Qnfz3hOrlkDfnT4vITBcfXOTS8h8p+tMuLBL1xJtB2t9TM7v5dKO6c3XKO5bHwiR79t5muaRo8ODBVKtWjSlTptC+fXuOHz/O/PnzmT9/vqHO6tWrsbe3x9XVlXPnzjFo0CBatmyZbpl2t27dKFy4MFOnTgVg0KBB1KpVi2nTptGiRQs2bNjArl27OHjw4Ft/RiGEEG+RXg9/ToSDs9SyTKh+JYev7OT8qEHUvKj2DF0uDH+8VxSf0nUZU7w9rrbZf8VelkuKKlasyLp16xg1ahQTJ07Ew8OD2bNn06VLF0OdsLAwhgwZwr1793BycqJbt26MHTs2XTshISHpNv6rVq0aK1as4LPPPmPs2LEUK1aMlStXvpk9ioQQQmQNSbGwti9c2ayWa34KdcfIhOqXFH/yJGaDRlLzgUKaBh52rkflj0fi+MCCasUK/ncD2USWm2idVb3IRGt3d/eXPnZCCPFqEhISuHnzpky0Fs/2KOTvCdXn1QnVLb6FMllno8Ds4MaDIE5OGYbPlstoFQjPC98019Gm+RT+POnE7ssR/PBeBRqVLmTsUJ/pZSZaZ7meouzoyQaE8fHxkhQJ8ZY82YDzZfepErlEyDFY2QXiItUJ1Z2WQ5FXP34oN0q+eZO7H3an9I0HAOwpo2FxfS2xKd5MX2tCdHwEZjotD+OSjRxp5pGkKBPodDry5s1LREQEAFZWVhmOiBBCZB69Xk9kZCRWVlaYmMivMfEvp3+FTYPVCdWFfKHTCrArYuyosqzAiECOhx8nNjmW2JRYYpMf47rvGjV/v0b+ZD2xFrCzgycN3p/IneNp/Hb0HgkoeDnkYU5HP3ycn9/7kp3Ib5NM8uSE9yeJkRDizdJqtbi6usofIOL/0lJh51g4+r1aLvmuOqE6mx259DYpisLAPwcSlRQFQJ54hQ+26ql8VZ1Zc95Nw7fNtNTwacC41bFcDFOPvHqviitjmvhgaZazemolKcokGo0GJycnHBwcnnrQqRAic5mZmaVbTCFyufiH8HtPuL5HLdcZBbWGy4Tq/6DRaBhcYTDrg9bDibP03ZhI/lhI1cLy2lr+qKQhOboyy3a5gBKDiUki33aqSqNSzsYO/Y2QpCiT6XQ6meMghBBvU+QVWN4RHl4HUyu1d8inubGjyjZauDWl2tprPPztOACPC9nyReM4rhe0JvFua1IflwZAZ30Nc6dVKFZjAUmKhBBCiKzl6nb4vRckPwY7V/WE+0K+xo4q20i6do3QYcNJunwZgAOVrPmxVhzxyZ4k3miPkmoHpGLusB3T/AcpZO2Ak7WTcYN+gyQpEkIIkf0oChyaDbsmAAq4VYf2v4B1ztkzJzMoikJ0UjSRCZFEJkRyP+E+kfGR3I+LxGHzCfzXXcIkVSHGEuY11XLSM5XU+01JvF8T0FAknynzulTHt0gLYz/KWyFJkRBCiOwlJUE9ruPcarXs3xMaTQMTM+PGlQUkpiYy7cQ0rjy8YkiCUvWp6erkj1Hov1lPmZvqZOq/imqY11RLlFlBUkO6khSvLhzqVMmFsc18sDLLPalC7nlSIYQQ2V90qLr/0N3ToDWBxtOgYm9jR5Vl3Iq5xe9Xf3/m/WoX9fTeridPIiSZwC/1tDh06krn5AC+2RVOUrIeO0tTprXxpVHpnDtM9iySFAkhhMgebp9QE6LYe2CZXx0u86hp7KiyFO983vxY/0cuR10mMl4dMtt+czvWCQo9d+gN55YFOcE37+q4m9cK/V5T4h/dBaBK0fzM6lAOJ7vcuRGxJEVCCCGyvsDf4I9B6oaMDqXUCdX53I0dVZaj0WioVrga1QpXA+DAnQOUuJFC8YU7KPgY0jSwtrqGK++WIeGBNQk36pOWYouJVsOQht58UKsYOm3u3ftLkiIhhBBZV1oq7BoHR75Vy3LC/Qt7/PgBR0d9QLMTau9QWD61dyjUxZZmutEcvnYbRQH3AlbM7uhHOZe8xg04C5CkSAghRNaUEKVuyBj8p1quPRJqj5ANGZ8jTZ/GuqB1XDy2hXqLztIsTE2IdpbT8Es9LQHFPiD+bGkWXbgNqJOpP2vqg7W5pAMgSZEQQoisKOIyrOgMD4PVDRlbzoNSLY0dVZZ24f4FJh2eQLGtF+i4T4+JHh5ZwQ9NtJzy1JISVZXlu9zQ6+PJZ2XKl23KEFAq655ubwySFAkhhMhaLm2CdR9AcizYuUDH38CpjLGjytKS0pIYtvJ9eq2Pw0ftBCKmUgnCBjSnllkewg6YcfWe+pFf29uer9qWwcHWwogRZ02SFAkhhMga9HrY9yXsm6aW3WtCuyWyIeN/UBSFyDWrmfRjHFbJkGKmQzO4F5W6f8LOi/f4Yu05HsYlY26iZXSTknSr6iYHKT+DJEVCCCGMLzFG7R26skUtV/4QGk4Cnalx48qiHiU+Ijg6mJu3zmI3ZzmFT93GCrhcGL57F8JZSv3V1Vn/VyQAJZ1smdOxHN6ONsYNPIuTpEgIIYRx3Q+CFZ3g/lXQmUOzWeDXxdhRZVm7bu1i8N7BVLim54MtevLGq6far66hZUNVDSlJLiTc6MD65Eg0QJ9aRRna0BtzEzms/L9IUiSEEMJ4rm6HNb0hKQZsnKHDr1CkgrGjyrISUxPZeXEDH25Oo+5ZdWVZSEH49l0dNxx1JD+oQ3JkPUCHxuQRdi4b+KTBEkmIXpAkRUIIId4+RYEDM+HPLwAFXKqoO1TbOBo7sixr7+29/L5iPO1X38M+BvTAH5U1rKqlJUmfn8RbHUhLcAfAxOYMFk7rqVSkHOY6c2OGna1IUiSEEOLtSoqFDf3h4ga1LAe6PlOKPoXgR8F0XtuWzvv0fPT3RowP8plwqX8D8voWp/4dezYeMyctBazMtIxo4k4n/waYmYw2cvTZjyRFQggh3p6H12FFF4i4CFpTaPIV+PcwdlRZ0p6QPYw+OBrHW4+ZtimNIg/U6zvLaSg65nNaezRn9LpzbL9wDwB/t3zM6lAOl/xWRow6e5OkSAghxNsR/Ces7gGJjyCPI7RfCq6VjR1VlnTl4RVG7x1O4/1xtD6koFMgKo+Gk90rUrZFT9LiShAw+wD3Y5Mw1WkY3EDOLcsMkhQJIYR4sxQFDn+jnmGm6KFwBXVCta2zsSPLclL0KVy4f4GRS99jzB9pFAtXr9s2aYzX2LGUtbLhi82XWH78JABeDnmY1aEcpQvbGTHqnEOSIiGEEG9Ocjz8MRDOrVbL5d6DpjPBVHZT/rc9IXsYfWAktQ/HMm2vHrM00OexosiEidg1bcqpW1EM+ekAtx7EA9CrhgfDAopjYSoryzKLJEVCCCHejKibsPI9CD8HWhNo9CVU7A2ym/JTXb10kKE/x1AqRC3HVShOma9/RFPQnpk7rvDdniD0CjjbWTCjXVmqecpO35lNkiIhhBCZL2g3rOmlnnRvVVA9rsOjprGjypIeJDzg8A8TqLh4J5bJkGgKBYd/Son3ehIcGcsn3x/ifGgMAK38CjO+eSnsLGWn7zdBkiIhhBCZR1Hg4Cz4c5I6f8i5PHRYCnZFjB1ZlqIoCufvn2fjkUUU+2EHZa/rAbhUBDZ0dGN22zYsOXyTL7deJilVj52lKZNblaZZGZmH9SZJUiSEECJzJD2G9f3h0ka17NcVmsyQ+UP/cD/hPpuvb2b9tXU4H7xGj516rJMg1UTDw26Nqf7hp1ROs+XjZec5GHQfgFp/n2rvKKfav3GSFAkhhHh994NgZReIvPz3/kPToUIPmT/0t203tjFs/zAA7GIVPtiqxz9I3YjxnrstVz9qROM6vTkarGfs+oPEJKZiYaplTJOSvFdFTrV/WyQpEkII8Xoub1FPuE+KARsn9bgOl0rGjsroopOiiUmKIT41Xk2IFIVqlxR6bddjk6ge4rqqppaNVeJIi9zKkp+KcPeeOjxWtogdX3coRzH7PEZ+itxFkiIhhBCvRq+HfV/Cvmlq2bUqtPtZzi9DHSZrtKYRSWlJANjEK/TerqfqZbV36LojfNdMx20HDamx3iSGtSU21RadVsPH73jyUV1PTHVaYz5CriRJkRBCiJeX8AjW9oFrO9Rypb7QcLKcX/a3G9E3DAlRpSt6+mzTYxev9g6tq6ZhbTUtqRpzUu+1IPFhBQCK2lvzdftylHPJa8TIczdJioQQQrycexdhRWeIugEmFtBsNpTrZOyospSE1ASsExR67tBT86LaO6QUdeHxsG74FrEmLeQxG47YERurfgz3qO7O8IASWJrJRozGJEmREEKIF3d+LWz4CFLiwc5VXW7vXM7YUWU5YTs3MfOnNPLHgqLVULB3HwoO+IhkjY6vd17llwPXURQonNeSr9qVoVox2YgxK5CkSAghxH9LS4Xd49UzzACK1oE2i8C6gDGjynLSHj/m3tSplFn7BwCh+aHyd0ux8avA+dBohqwK5Oq9WADa+xdhbDMfbCxkI8asQpIiIYQQzxf3AH7vATf2qeXqg+Cdz0EnHyH/FHvgIKGfjUZ/LxI9sLmShhW1tKwvWohFu67xzZ/XSNUrFMxjzpetfanvIxPSsxr5P1oIIcSz3TkFq7pBzB0wtYaW30GpVsaOKktJi4nh9pRJJKzfBEB4XnVl2RUXDWlJ9gz45QZnQ6MBaOJbiC9a+pLfWiakZ0VZcr1faGgo7733HgUKFMDKyopy5cpx6tQpw/3Y2FgGDBhAkSJFsLS0pGTJksybN++5bS5ZsgSNRpPhKzEx8U0/jhBCZD+KAicXweJGakJUwBP67JaECHiY+JCDoQeZf3Y+M7/pwvH61UhYv0ntHfLXMKyXjstFtCQ/rE7KraGcDY3G1sKEOR3L8V3n8pIQZWFZrqcoKiqK6tWrU7duXbZu3YqDgwPBwcHkzZvXUGfw4MHs2bOHX3/9FXd3d3bs2EH//v1xdnamRYsWz2zb1taWK1eupLtmYSHbpgshRDopCbBpCJz5TS2XaAYt54GFrXHjMqLY5FimHp/KifAThMWFYZWo8P4uPU3OqSvLwvLB2vbOJPp40MTcl2NnvblyTz3PrJa3PdPblKGQnXzeZHVZLimaNm0aLi4uLF682HDN3d09XZ0jR47w/vvvU6dOHQD69u3Ljz/+yMmTJ5+bFGk0GgoVKvQmwhZCiJzh4Q1Y1RXCz4FGC/XGqXOIcvkxE9ceXWNjsHqmm1+Qng+26skfC3pgS0UN+5sUIaBEC+xT32XSpkvEJqViZaZjdJOSdKnsKsd0ZBNZbvhs48aN+Pv7065dOxwcHPDz82PBggXp6tSoUYONGzcSGhqKoijs2bOHq1evEhAQ8Ny2Y2NjcXNzo0iRIjRr1ozTp08/s25SUhIxMTHpvoQQIke7uh3m11YTIquC0HU91Pgk1ydEAGXtyzKyxACGb7Ng1Go1IbqbD8Z11fFLfR3X42KZszmZEWvOEZuUir9bPrYOqinnlmUzGkVRFGMH8U9PhrOGDBlCu3btOH78OJ988gk//vgj3bp1AyA5OZk+ffrwyy+/YGJiglar5aeffqJr167PbPfo0aMEBQXh6+tLTEwMc+bMYcuWLZw5cwYvL68M9cePH8+ECRMyXI+OjsbWNvd2IQshciB9Guz9EvZPV8tFKqrHddgVNm5cWUBsciydt3Qm/8lg+mz7f+/Qk5VlFlZ2FExtxJVrZUlK0WFmomVoA2961yyKTivJUFYQExODnZ3dC31+Z7mkyMzMDH9/fw4fPmy4NnDgQE6cOMGRI0cAmDFjBgsWLGDGjBm4ubmxf/9+Ro0axbp166hfv/4LvY9er6d8+fLUqlWLuXPnZriflJREUlKSoRwTE4OLi4skRUKInCX+IazpDcG71XLFPhAwJdce13Eu8hx/XP8Dd1t3ClkXIioihPtTp1PrgvpRGZofFr5ryXsdJlHYogTf7Ypi58V7gHqI64x2ZfFytDHmI4h/eZmkKMvNKXJycsLHxyfdtZIlS7JmzRoAEhISGD16NOvWraNp06YAlClThsDAQGbMmPHCSZFWq6VixYpcu3btqffNzc0xNzd/jScRQogsLvQvWPU+RIeAiSW8OxvKdjR2VEb18Z8f8yDxAQD+V/X03abHJw70Gvijkoak7q1YWnsCW89F0OPX80TFp2Cq0zConhf9ahfDRA5xzdayXFJUvXr1DCvErl69ipubGwApKSmkpKSg1ab/H0+n06HX61/4fRRFITAwEF9f39cPWgghsptTP8OWTyEtGfJ5QIdfoVBpY0dldKY6U/LEK/TY+f8zy+4W1HK2Ty16tpuAiZKXQcvPsvlcGAA+TrbMbF+Wkk4ygpATZLmkaPDgwVSrVo0pU6bQvn17jh8/zvz585k/fz6gLquvXbs2w4YNw9LSEjc3N/bt28cvv/zC119/bWinW7duFC5cmKlTpwIwYcIEqlSpgpeXFzExMcydO5fAwEC+++47ozynEEIYRUqCmgyd/lUtezeGVj+AZV6jhpUVKIpC5zAPSiy+Q954tXdoY2UNLacup55TGbadD+ez9fu4H5uMiVbDR3U9+aiuJ2Ym0juUU2S5pKhixYqsW7eOUaNGMXHiRDw8PJg9ezZdunQx1FmxYgWjRo2iS5cuPHz4EDc3NyZPnky/fv0MdUJCQtL1Jj169Ii+ffsSHh6OnZ0dfn5+7N+/n0qVKr3V5xNCCKOJuqUutw87oy63rzsGagwBbe7+UFcUhYPnNhE5aQpVzj0C4E4BdVfqcrXbYmvlyicrTrM+8C4AxR1tmNm+LKUL2xkxavEmZLmJ1lnVy0zUEkKILOfaTljbBxKiwDI/tF0Ixd4xdlRGpSgK+2/v49CiKTRYf5s8iZCmgZ21bPilUjypJhoGeP/Awj8TiHichFYD/WoXY1B9L8xNdMYOX7ygbD3RWgghRCbSp8HeqbD/K7Xs7Aftl0JeF+PGZWTHwo7x0+4vqfPbFVpdV/sGrjvC6d7V6df2Sw5vH8L5K8WZeukRAMXsrZnZvhzlXPIaL2jxxklSJIQQOVVsJKzp9f/T7Sv2/nu5fe5eWZucmsTarz6g364krJIhRQera2j5o7KGNP0xYvdv5s6l90mITkSjgd41PBjasDgWptI7lNO9dlKUkpJCeHg48fHx2Nvbkz9//syISwghxOsIOQqru8PjMDC1gnfnQpl2xo7K6JJv3SLss7F0P6HuQ3elMPzQREdoQQ1KmjlEtGbFpYJAIu4FrJjRriz+7vK5llu8UlIUGxvLsmXLWL58OcePH0+3yWGRIkVo2LAhffv2pWLFipkWqBBCiBegKHDkO9j5OShpUNBbHS5zKGHsyIwiKS2Jk+EnOXh7P5pVm2m44wHmqZBspuX3umZsKJdCSftS1M3ThXWHrQiPUT/PuldzZ3ij4liZyYBKbvLS/7VnzZrF5MmTcXd3p3nz5owcOZLChQtjaWnJw4cPOX/+PAcOHKBBgwZUqVKFb7755qnHaAghhMhkidGwvj9c3qSWS7dRe4jM8xg3rrcsPC6cuX/N5fDdwzxIfEDh+wofbk7DW108xjk3DT821hCRL5WyBargmNCfeQfDgCTcClgxvU0ZKhctYNRnEMbx0qvP2rVrx+eff/6fmx4mJSWxcOFCzMzM6N2792sFmRXI6jMhRJYWfg5WdYOH10FrCo2mqnOIctlhpClpKdRZVYeY5Bh0aQotjiq0OaTHNA3izWFvS3fe6T+ZR0nRnL2VxvKDEB6dhEYDPap5MCygOJZmMncoJ3mjq89Wr179QvXMzc3p37//yzYvhBDiZZ1eBpuHQGoi2Lmoh7kWqWDsqIxCp9VRzbkal49u5cMtaXiox5JxylPDggAt79VuTVFbX77YdJHVp+4A4F7Aiq/alaWizB3K9WSwVAghsquUBNgyDE4vVcueDaD1fLDKfR/uqfpU+u7sS+Dt4/Q4ZkmPQ3o0eoixhCX1tdyp6sGXVcaQGFOMhrP2cS9G7R3qWd2DTxtK75BQvVRSFBUVhaIo5M+fn8jISPbv30/x4sUpXVrOyxFCiLfqQbB6mOu9c3/vTj0aagzNtbtT6xU9CcdPMmNLGk5RsQAcLqnht8ZWdK76IZM8OjF1yzXW/HUCAI+C1nzVtoysLBPpvHBS9NNPPzF16lT0ej3Dhw9n2bJllClThnHjxjFw4ED69u37JuMUQgjxxKU/1AnVSTFgVVDdnbpoHWNH9dZdjbrKvMB5RD+4S60NN/j8RAoAD2zgpwAtp7y0rGz2C+GRBWg697Chd6hXdXXfIekdEv/2whOty5Yty7Fjx4iPj8fV1ZUbN25gb29PTEwMtWrVIjAw8A2Halwy0VoIYXRpKbB7Ahz+Ri27VIF2i8HW2bhxGUF0UjRt/2iLy1936bVDT361c4jtfhp+q6MlwUJD2fxVKRD3ARsCwwEoWtCar9qVoYKb9A7lJm9korVOp8PCwgILCws8PT2xt7cH1FPrNblsdYMQQrx10Xfg955w+5harjoA6o8HnalRw3rbElMTufDgAp/8/j69duipfEX9u17vUgizMZ/wXuXq9DezZf+VKEavO0fE43A0GuhTsyhDGnjLrtTiuV44KTIxMSExMRELCwv27dtnuP748eM3EpgQQoi/XdsJa/tCwkMwt4MW34JPc2NH9dY8THzIt6e/5fDdw9yLC6dmYApf/6lXD3DVgkOfDyjY/0O05uY8jEtm2OoLbPj7RPui9tZ81bYsFdzyGfkpRHbwwknRn3/+ibm5el6OnZ2d4XpCQgILFy7M/MiEECK3S0uFPV/AwVlq2akctFsC+T2MGdVb99O5n1h9dTWODxVGb9NT+pbaO/TIvQAlZnxDvtJ+KIrC5rNhfL7hPA/iktH+3Ts0WHqHxEt44aQoT56n74jq4OCAg4NDpgUkhBACiLkLv/eCkMNquWIfCJicKw9z7VisPfFLfqPVviTMUiHZBIoMGU6Jbl3RmJgQ8TiRsevPs/2CuimRt2MevmpblrJyor14SbJPkRBCZDVBu9Xhsvj7YGYDzedC6dbGjuqNik6KZn3QelL0KZjrzDHXmWNhYkGem5EUnL2KDtfUM8nOumtY0EjLrp49UBSFNafuMHHTRaITUjDRauhf15OP6hbD3ER6h8TLy/Sk6MGDB5w5c4bAwECGDBmS2c0LIUTOpU+DvV/C/q8ABQr5qrtTFyhm7MjeuI3BG5lxcoahbJqi0PaQnuZHFXQKxFrAL/W07PXV4F+oIncfJTB63Tn2XokEoHRhW6a3KYuPs6wOFq/uhZfkBwUFMXbsWPLmzcuUKVPIly8f165dIzAw0JAEnTlzhrt376IoCtbW1jlqErYsyRdCvFGPw2FNb7h5QC1X6AGNvgRTC+PG9ZaExobS7o92PE5+TPlbJny4U4NdZAIAF8vYsb65PfctU0nT6/Gz/IjNJ8yJTUrFzETLJ/W96FuzKCa63LlxpXi+N7Ikv0uXLrz33nt4eHhQqlQpHj9+TFxcHHZ2dvj4+FC6dGm2bt3KwoULqVevHi4uLq/9IEIIkStc36cmRHERYGqtDpf5tjV2VG9VPvN8lNAWptIfF6h9PhUAE0dHCn0+lpL16tEGCHkQz8i1Z1l56gGQSnnXvExvWwZPBxujxi5yjhdOiu7fv0/p0qUpWrQoERERjBgxgv79+1O4cGFDnUWLFlGpUiVJiIQQ4kXo09Shsr1fAgo4lIL2P0NBL2NH9tYkpyXTfmM7nA8F8cFuPbYJoAeC63nSdNpydHnyoNcr/HzkJtO3XSEhJQ0LUy3DAkrQvZo7Oq3skycyzwsnRXPmzKFfv37Y29vzww8/MGfOHC5cuMD06dPx9vZ+kzEKIUTOExsBa/vA9b1q2a8rNJ4OZlZGDettSU5LZn3QeuZvnUifbXrK3FRnciS7FSJ5eG9q12iGzjwPwZGxjPj9LCdvRQFQpWh+prUpg1sBa2OGL3KoF06KmjVrRrNmzQzlHj16MG/ePGrVqkWbNm0YN27cGwlQCCFynJsH1eX2seFgagVNv4ZynYwd1VvzZ8ifTD30BZX33WPmQb26zF4HhT4eiEOv3mhMTUlN0zNvbzCzdl0lOVVPHnMTRjUpQaeKrmild0i8Ia88K02n0zFgwAAuXbqETqejRIkS6PV60tLSMjM+IYTIOfRpsG86/PyumhDZl4A+e3JVQgSwYeMMhswLo8teNSE656bh09467rethcbUlPOh0bT47hDTtl0mOVVPbW97tg+uRZfKbpIQiTfqhVef/ZeLFy8yePBgTp8+zfDhw/noo4+wtLTMjKazBFl9JoR4LY/D1eGyG/vVctlO0HQmmOWOYaA0fRpX75whYu4cCm4+jlaBx38vs9/nqwGNhu/fWcDhC3mZv/86aXoFO0tTxjbzoU35wnLGpnhlL/P5nWlJ0RObNm3i008/JTo6mrCwsMxs2qgkKRJCvLKg3bDuA4iL/Hu4bCaU62zsqN64/Xf288OZHzh3/xzlr+npvUNPwZi/75XS8HN9LZVLNKCSUyWsUn2Zs+0+1+/HAdC0jBPj3y2FvU3u28FbZK43siT/iZCQEFxdXZ95v1mzZgQEBPDtt98CEBoamm6FmhBC5BppqbBnMhz8Wi07loa2i8E+5y9OORd5jo92f0TeWIXBO/VUvaz+/X0vLywI0HK2qDp7Y2ylyXy9M5hfj94CwNHWnEktStOwVCFjhS5ysZdOiipWrEjz5s3p06cPlSpVemqd+Ph4rK2tKV26NB988AEff/zxawcqhBDZSvQddTL17aNq2b8nBEwB05wzreB58pnnpcFferrs1WOVBGka2FRJw+81tCSZqUNh77tPp+ncI4RFJwLQqZILIxuXxM7S1Jihi1zspZOiS5cuMWXKFBo1aoSpqSn+/v44OztjYWFBVFQUFy9e5MKFC/j7+/PVV1/RuHHjNxG3EEJkXZe3wIb+kBAF5rbqZoylWhk7qjfuUeIj9tzeQ+Dh9fj/fII+d9XeoaBCML+xjiIVajLGrSHlCtRk9o7bfLv1LpCIWwErprb2pVqxgsZ9AJHrvfKcosTERLZs2cKBAwe4efMmCQkJFCxYED8/PwICAihdunRmx2pUMqdICPGfUpNh1zg4+r1advZTh8vyexg3rjcoPiWesYfGsuPWDsyTFdoe1NPsuHpeWbwZrKitZXt5DQc6H8LWzJaNZ+4y4Y+LPIxLRquB3jWLMri+N5ZmcoCreDPe6JyiJywsLGjdujWtW+fsk5uFEOKFPLwOv/eEu6fVcpWPoP54MDEzalhv2uRjk9lxawd+QXp67dDjEK1eP1pcw+IGWqJsNPgW9CUuwYzBy0/y5+UIAEoUsmF62zKUKZLXeMEL8S+vnBQJIYT42/m18McgSIoBy3zQch4Uzx1TB+pZlafYunWGidSRtrCwoZa/vLToNDo6eXfAMa0NDWftVw9w1WkZWM+TD2oXw1QOcBVZjCRFQgjxqlISYNsoOLVYLbtUgbYLwa6IceN6g1L1qUQnRaNPSyVh1Tqcv1uAU5xCmgY2V9Kw+u+J1OOrjqds3gBGrz3H8ZvXAKjglo9pbXzlAFeRZUlSJIQQryLyKqzuDhEXAA3UHAJ1RoMu5/5a1St62mxsg/5KMH23puEZrl6/6gwLGum45aghj2kePvUbSujtMoz6+QDJaXqszHSMaFSCrlVkR2qRteXcn14hhHgTFAVO/wpbh0NKPFjbQ6sfwbOesSPLVLtu7eJo2FGqOVfDwcoBe0t7LFKgxpprNDmpoFUg3hyW1zHhr6oFeNerBe2823HvoTWj1p7l6r2rANQtbs8XrXwpnDd3bEUgsrdM39E6p5LVZ0IIEqNh02A4v0Yte9SC1j+BjaNx48pktx/fpsnaJumuVbimTqR+siO18k41vMZPxdTBAYDYpFS+2naZX47eQlGggLUZ45qX4t0yTnJEhzCql/n8fu1ZbgcOHOC9996jatWqhIaGArB06VIOHjz4uk0LIUTWceck/FBTTYg0Oqg3Drquz3EJEYCztTMB7gEAFIhWGPZ7GiN+VxOie3lhcnsteb4cZ0iIdl+6R4Ov9/HzETUhaluhCLuG1KZ5WWdJiES28lrDZ2vWrKFr16506dKF06dPk5SUBMDjx4+ZMmUKW7ZsyZQghRDCaPR6ODRbPa5Dnwp5XaHNInCpaOzI3hidVsdX1abyySV3Hv70A6bJCqna/+9IPeGdL3G1dSXycRLj/7jA5rPqOZeu+a2Y0sqXGl6yCaPInl4rKfriiy/44Ycf6NatGytWrDBcr1atGhMnTnzt4IQQwqgeh6sHuV7fq5ZLtYZ3Z4OFnTGjeuNuH9xB2IQJ2Nx+iClw0QV+CtBxx17t9THTmrHqxG0mb7lEdEIKOq2G3jU9+KSebMIosrfXSoquXLlCrVq1Mly3tbXl0aNHr9O0EEIY17WdsK4fxN9XT7ZvPB383oMcPByUGhVFxFcziF27FhsgxhKWvqMlrVFN1tT5GksTS86H3WPq5pscDj4LQOnCtnzZugylC+fsRFHkDq81p8jJyYmgoKAM1w8ePEjRokVfud3Q0FDee+89ChQogJWVFeXKlePUqVOG+7GxsQwYMIAiRYpgaWlJyZIlmTdv3n+2u2bNGnx8fDA3N8fHx4d169a9coxCiBwqNQm2jYZlbdWEyNEX+u6F8l1zZEIUnxLPmsur+fPb0Vxp2IDotWsB2FVOwyd9dTi07cA39b7FVGvBvH3BtP3+NIeDH2BhqmVMk5Ks719dEiKRY7xWT9EHH3zAoEGDWLRoERqNhrt373LkyBE+/fRTPv/881dqMyoqiurVq1O3bl22bt2Kg4MDwcHB5M2b11Bn8ODB7Nmzh19//RV3d3d27NhB//79cXZ2pkWLFk9t98iRI3To0IFJkybRqlUr1q1bR/v27Tl48CCVK1d+pViFEDnM/SBY0xPCzqjlSh9Ag4lgamHcuN6gnzdMoNC8jTip62S46QALAnSklirK/uZrMNWacub2I0auPcelMHXpWU2vgkxu6YtrASsjRi5E5nvtJfljxoxh1qxZJCYmAmBubs6nn37KpEmTXqm9kSNHcujQIQ4cOPDMOqVLl6ZDhw6MHTvWcK1ChQo0adLkme/boUMHYmJi2Lp1q+Fao0aNyJcvH8uXL89QPykpyTBxHNQlfS4uLrIkX4icKnA5bB4KKXFgmR9afp+jj+rQx8Vxd+5sYpYuQ6NXSDCDVTW1bPXXoNdqaF6sOaMqTmDG9iv8cuQmegXyWpkytqkPrcsXllVlItt4q0vyJ0+ezP379zl+/DhHjx4lMjLylRMigI0bN+Lv70+7du1wcHDAz8+PBQsWpKtTo0YNNm7cSGhoKIqisGfPHq5evUpAQMAz2z1y5AgNGzZMdy0gIIDDhw8/tf7UqVOxs7MzfLm4uLzyMwkhsrDEGFjbF9b3UxMi95rw4aEcmxAdu3uUXqNLcaiOP49//hWNXuFocQ2D++jYXEmL/u8dpz1NW9Pg630sOawmRC3LObNrSG3aVCgiCZHIsbLc5o0WFmo39ZAhQ2jXrh3Hjx/nk08+4ccff6Rbt24AJCcn06dPH3755RdMTEzQarX89NNPdO3a9ZntmpmZsWTJEjp37my49ttvv9GjR490PUJPSE+RELlA6Cn4vRdE3VD3Hqo7CmoMAW3OXEGVGHKL3z9sTIVg9df+vbzq4a356tTDp4APPgV8sDfzZO6OMLaeV8/wcMlvyeSWvtTytjdi5EK8upfpKXrtYz4SExM5e/YsERER6PX6dPeaN2/+0u3p9Xr8/f2ZMmUKAH5+fly4cIF58+YZkqK5c+dy9OhRNm7ciJubG/v376d///44OTlRv379Z7b9779uFEV55l885ubmmJubv3T8QohsQJ8Gh+b8f+8hOxdosxBcc+b8wgcx4ayd8D4Vt4dQIRVStbChioa11bSsbLMOr3xe6PUKy46H8NHWczxOSkWn1dCnZlEG1fOSZfYi13itpGjbtm1069aN+/fvZ7in0WhIS0t76TadnJzw8fFJd61kyZKsWaNuq5+QkMDo0aNZt24dTZs2BaBMmTIEBgYyY8aMZyZFhQoVIjw8PN21iIgIHB1z3m60QojniA5V9x66+fe8RZ+W6t5DlvmMGVWm2x2yGzOtGU4X7nFr/FhqRKnXz7tp+ClAy0ctpvBZMfUP1yvhjxm19ix/hTwCoKxLXqa28sXHWXrFRe7yWknRgAEDaNeuHZ9//nmmJRfVq1fnypUr6a5dvXoVNzc3AFJSUkhJSUGrTT8dSqfTZeip+qeqVauyc+dOBg8ebLi2Y8cOqlWrlilxCyGygYsbYONASHwEptbQZDqU65LjltrHJscyacMg3t+lp+AVBWfgYR51z6FrFRxZ/u4KHKwcSExJY+7ua8zff51UvUIecxOGBRTnvSpu6OQ0e5ELvVZSFBERwZAhQzK1t2Xw4MFUq1aNKVOm0L59e44fP878+fOZP38+oG4MWbt2bYYNG4alpSVubm7s27ePX375ha+//trQTrdu3ShcuDBTp04FYNCgQdSqVYtp06bRokULNmzYwK5du+SMNiFyg6RY2DYSTi9Vy87loc1PUKCYceN6A5TkZO4vWMCs+WlYpECaBrb6azjXohSf1h5LGfsyABwKus/odee49SAegIBSjoxvXgonOznNXuRer5UUtW3blr1791KsWOb9YqlYsSLr1q1j1KhRTJw4EQ8PD2bPnk2XLl0MdVasWMGoUaPo0qULDx8+xM3NjcmTJ9OvXz9DnZCQkHS9SdWqVWPFihV89tlnjB07lmLFirFy5UrZo0iInC70L1jTGx4GAxqoMRjqjgadqbEjy3RxR48RPmkSycHBWACXisCad/PTvulwhhZ7F61Gy4PYJCZvvsTa0+rGRIVsLZjQohQBpQoZN3ghsoDXWn0WHx9Pu3btsLe3x9fXF1PT9L9kBg4c+NoBZhUvM3tdCJEF/HsytW1haPUjeNQ0dmSZLuVeBBHTpxOzeTMAj6zg13e07C+t4Ujno+Qxy4OiKKz5K5TJmy8SFZ+CRgPdqrjxaUBxbCxyXoIoxBNvbfXZb7/9xvbt27G0tGTv3r3pVnJpNJoclRQJIbKRDJOpW0Cz2WCV36hhZTYlJYWHy5YROnsmpomp6DWwvbyGlbW0WOQtwKwqY8ljlofgyFg+W3eeI9cfAFCikA1TW/vi55qzJpcL8bpeq6eoUKFCDBw4kJEjR2aY+JzTSE+RENlELplMHX/yJOETJ5F09SoAV51hYYCOuKKODCo/iObFmpOYksa8vcHM2xtMcpoeC1Mtg+p507umB6a6nP07W4gn3lpPUXJyMh06dMjxCZEQIhtIjlMnU//1i1p29lP3Hsphk6lT798n4qsZRG/YAECarRXzayRy1M+S/Z0OYmGiboB7OPg+n607z/X7cQDU9rbni5alcckv55UJ8SyvlRS9//77rFy5ktGjR2dWPEII8fJC/4K1feBBEDl1MnVCYix3li4k+Yef0cYloGjgr8oF+K7yI2KttJQvWAoLE4sME6ntbcwZ964PTX2d5HgOIf7DayVFaWlpTJ8+ne3bt1OmTJkME63/uUReCCEynT4NDs+FP79QJ1PbOEPr+TluMvX+TT+QMH0urhEKWiC4EPwUoCPYORpLEyv6lHyPbj7vs/JECFO2XCY6QZ1I/V5lN4Y1Ko6tTKQW4oW8VlJ07tw5/Pz8ADh//ny6e/IXiRDijYq6Bev6QcjfhzqXbA7vzslRk6lT7kUQ8dVX2G/aBMBjC9jR2J67dX2oYFOYljYuNCvWjKgYc/osOceJm+q21SWdbJnSqrRMpBbiJb1WUrRnz57MikMIIV6MosDZlbBlGCTFgFkeaDwtR02mVpKTebh0Kfe/+x59fDxoNOwoBytqaVnSYQHF8xcHIDEljW/+VHekTklTsDTVMaSBNz2qu2MiE6mFeGmvfSCsEEK8NfEPYdNguLheLbtUVvceyu9h1LAyU+zBQ9yeOA5C1DlBoW7WzG8Al+yTAIhKUnuD9l+N5LP15wl5qO5IXb+kA+Obl6JIPplILcSreq2kaOLEic+9//nnn79O80II8X/Be2D9h/A4DLQmUGckVB8Mupzxt11KaCjXJn2GZu9RQN2AcVldLft9E1E0Gsx15lRwrEAhM28+Xn6aP87cBdQdqcc3L0VAKUeZtiDEa3qtfYqezCd6IiUlhRs3bmBiYkKxYsX466+/XjvArEL2KRLCSFISYNcEODZPLRfwUidTFy5v3LgyiT4pibNzJqFZuhazFIU0DeyqaMrV1uUpWrg0xfMXp2T+krjYuLHq5F2mb7vM48RUtBp4v5o7QxsWJ495zkgMhXgT3to+RadPn37qm3fv3p1WrVq9TtNCCAHh52BNH4i8pJYr9oYGk8Asew4RpehTCH4UTH6L/OQzz0fi/oPcmzIV89u3AbjgCqe7VGBg25k4Wv//oO3zodG0X3aMM3eiAShd2JaprcrgW8TOKM8hRE6V6X9e2NraMnHiRJo1a0bXrl0zu3khRG6gT4Mj38LuSaBPAWsHaPEdeDc0dmSvZdyhcfxx/Q8coxS679JTIUjtqI+xNWFRXT2HS2qoWdgWBysHAB4npjBzx1V+OXITvQJ5zE34tKE3Xau6o9PKUJkQme2N9Lk+evSI6OjoN9G0ECKne3RbXWp/66BaLt4Ums8F64LGjSsT+Fp7Yb0vjXePKZimQaoWNlfSsKaaQqK5ulrsQOgB7ifc51hQChP/uEjEY3WCdbMyToxt5oOjrYUxH0GIHO21kqK5c+emKyuKQlhYGEuXLqVRo0avFZgQIpdRFDi3GjZ/CknR6rlljb8Ev67Zfqm9oijEbNqM/1dLKBeh9g6dcdewuKGWuwX+/2xtvNrgalmBoStusP9qJADuBayY2KI0tbztjRK7ELnJayVFs2bNSlfWarXY29vz/vvvM2rUqNcKTAiRiyREweahcH6NWi5SUZ1Mnb+ocePKBIkXLxL+xWQS/l54Em9vw/H2pfje5kS6ZE/R69BENebLLQ9JTo3ETKflwzrF+LBOMSxMdcYKX4hc5bWSohs3bmRWHEKI3CpoF2z4GB7fBY1OXWpfY0i2X2qf+vAhkbPn8Gj1alAUFAtz/qhuxsry8aSYnAT+nxDZUx3lQQcWX7kPQA3PgkxsUYqi9nmMFL0QuVP2/q0jhMi+kuNgx1g4uVAtF/CEVvOhSAXjxvWalJQUwn/9mXvfzMU8PgWAe9W8mVT+FhE2CRS1K0atIrW4G3uXkKhHXLpagutRPkA89jbmjG3mw7tl5PBWIYzhpZOiIUOGvHBdORBWCPFUIcdg3QcQ9Xdvc6UPoP74bLvU/om4I0cI/WISacE3MAduOMLiBjouu1wHoLFHY8ZXHY+5zpLfjoew48xlYhNT0WigWxU3hgbI4a1CGNNLJ0VP25voaeSvHCFEBqlJsGeKerK9ogfbItDyOyhax9iRvZbkO6FETJvG4507AYixhBW1tZTpMQTvx9exjg2liUcT2nm343xoDJ9tOM2Z248A8C1sx+RWpSlTJK/xHkAIAbxCUiSHwAohXkn4OVj7AURcUMtlO6uryyyy7waE+oQEIufP58HChWiSU0jTwI7yGlbV1BJnqeHrMj0NfyBGJ6QwbuMFlh69haKAjbkJwxoVp0tlN9lzSIgs4rXnFD169IiFCxdy6dIlNBoNPj4+9OzZEzu77PuLTgiRidJS4fAc2DNV3YjRqiC8OwdKNjN2ZK9MURT2LPkC83kryB+jRwOcd9OwuL6W2w5qguNh54FGo0FRFNadDmXKlkvcj00GoEU5Z8Y0KYmD7DkkRJbyWknRyZMnCQgIwNLSkkqVKqEoCl9//TWTJ09mx44dlC+fM84mEkK8ogfB6tyhOyfUcolm0Gw25Mm+e+4kXLjAlc+H43RBnScUaQu/1NNyrLgm3RL7Co4VuHrvMZ+tP8/xGw8BKGZvzaQWpanmmf03ohQiJ3qtA2Fr1qyJp6cnCxYswMREza9SU1Pp3bs3169fZ//+/ZkWqLHJgbBCvARFgRM/wc7PISUezG2h8XQo2zHbbsSYGhlJxOzZRK9dB4pCkglsqKJlYxUNyab/fyZ3W3c+8B3I2atFWHTwBql6BQtTLQPredG7RlHMTLRGfAohcp+X+fx+raTI0tKS06dPU6JEiXTXL168iL+/P/Hx8a/adJYjSZEQLyg6FDZ8BNf/nn/oUVs9tyyvi3HjekX65GQe/vwzD374EX1cHADRtcsS8l4tJl2fl65uQ7cA6hUYwqRNlwiLTlSv+Tjy+bs+FMmXvVfWCZFdvczn92sNn9na2hISEpIhKbp9+zY2Njav07QQIrtRFDi7CrYMU4/pMLGEBhOgYh/QZr/eEUVReLxrFxHTvyLl71Psw1yt+bZ2IteKXIDrF9LV1ycXYP0BZ9bEqSt0i+SzZELzUtQr6ZihbSFE1vRaSVGHDh3o1asXM2bMoFq1amg0Gg4ePMiwYcPo1KlTZsUohMjqYiNg02C4vEktF/aHVj9CQU/jxvWKEq9c4d6UqcQfOwZAvJ0F1ztV43ePCK49uoyp1pRqztUA0KfpuHrDnaCbHqCYYqrT8GHtYvSv6ynHcwiRzbxWUjRjxgw0Gg3dunUjNTUVAFNTUz788EO+/PLLTAlQCJGFKYp6XtmWYZDwELSmUHsE1BicLY/pSH34kMg5c9WjOfR6FDNTNvjrWVM1hSSz/fBIrZeiT2F4xeEEh1kwbsMFQh6qUwVqehVkQnM5nkOI7OqV5hQFBgZSrlw5Qzk+Pp7g4GAURcHT0xMrq5w3di5zioT4l9hI2DwYLv2hlgv5QssfoFBp48b1CpTkZB7+9huR336HEhsLgPad6kzwu8E5swgqFaqEh50HYXFh3I29iza1AHlje7P70gMACtlaMLaZD018C8nGtUJkMW98orVWq8XPz4/evXvTuXPnXLEnkSRFQvzD+bXqqfYJD0FrArWGQc2hoMteR1QoikLsvn1EfDmN5Js3AfVojiX1dVxyVZObwnkKs+rdVdia2ZKUmsZPB27wzZ/XSEzRo9Nq6FndnUH1vcljnv16xoTIDd54UnTkyBEWLVrEqlWrSElJoXXr1vTq1Yu6deu+ctBZnSRFQgBx92HzELi4QS07loaW88CpjHHjegVJ165xb9p04g4eBCDZzoq1dS1YVzwaRxsnwuPCMdGY8EvjX/C192X/1UjGbbzAjfvqCrRK7vmZ2LIUJQrJ7wMhsrK3tiQ/ISGBVatWsXjxYg4cOIC7uzs9e/bk/fffp0iRIq/abJYkSZHI9S6sU3uH4h+ovUM1h0LNT8HEzNiRvZTUBw+I/OYbHq1S5w1hasr+ajYsLB9NgoUGM60Zf7b/k5sxN7HQWZBH68oXmy6y9Xw4AAXzmDOmaQlalissQ2VCZANvLSn6p+DgYBYvXswvv/xCWFgYDRo0YMuWLZnRdJYgSZHIteLuw5ZP1aQIwKEUtPwenMsZNayXpU9K4uEvv6Tbb+hqmQJsbpSPEyZ3SFVSaeDWgLZebalWuBrJqXp+Onidb3YHkZCShk6r4f2q7nzSwEtOshciGzFKUgQQGxvLsmXLGD16NI8ePSItLS2zmjY6SYpErnRxA2waAvH3QaODmkOg1vBs1TukKAqPt20jYsZMUkJDATAtWZzR5YMM84YAzLRmrGm+Bnc7dw5cU4fKrkeqyVNF93xMbFGakk7ysy9EdvPWNm98Yt++fSxatIg1a9ag0+lo3749vXr1yoymhRDGEPcAtg5Tl9sD2JeEVvPA2c+4cb2khLNnuTf1SxJOqxsqmjg4YD94MJ1S53EnXoONqQ1f1f6KvOZ5cc7jTGKSBR8t+4vN58IAdahsdJMStPKToTIhcoNXTopu377NkiVLWLJkCTdu3KBatWp88803tG/fHmtr68yMUQjxNl36Q92IMS5S7R2q8Ym695CJubEje2Epd+8SMWs2MX+o2wUoFmbEtW9IVOuanNYkcufoXQCG+g+leuHqJKfqWXToBnN3XyM+OQ2tBt6v5s7gBt4yVCZELvJKSVGDBg3Ys2cP9vb2dOvWjZ49e1K8ePHMjk0I8TbFRsLW4XBhrVq2L6HOHSpcwbhxvQR9XBz3f/qJh4sWoyQlAbDXV8Py2mlE2WyDY9sMdd1s3Wjj3YZDQff5fMN5gv8eKvN3U4fKfJxlqEyI3OaVkiJLS0vWrFlDs2bN0OlkG3shsjVFgXOrYesIdd8hjQ6qD4TaI8HUwtjRvRAlLY3odeuImDOHtMj7AKSWKc6Y8kGEOJvgYedBYZ0F5ibmWJhYYKmzpI5TKz767S82n30yVGbGqMYlaV1ehsqEyK1eKSnauHFjZschhDCG6DvqROpr29WyY2lo8W22mjsUd/Qo976cRtLlywAohR3Z1tSRxfkvgEZDPZc6zK4721D/yQaMI5cFkZDyGK0GulVVh8rsLGWoTIjcLEseXR0aGsp7771HgQIFsLKyoly5cpw6dcpwX6PRPPXrq6++emabS5YseeprEhMT38YjCZG16PVwYiF8V0VNiHRmUPcz6Ls32yRESdeuEfLBB4R076EmRHms+bOFG10632dxgYtoNFoauTfisyqfGV6z90oEjWYf4KvtV0hISaOiez42fVyT8c1LSUIkhMic1WeZKSoqiurVq1O3bl22bt2Kg4MDwcHB5M2b11AnLCws3Wu2bt1Kr169aNOmzXPbtrW15cqVK+muWVhkj+EBITLNg2DYOBBuqTs5U6QiNP8WHEoYN64XlBoZSeTcb3i0Zg3o9eh1Ws5Ud+RbvwgeW4UCGpoXa05v39542HkAcPthPBM3XWTnxXsA2Nuoq8pkA0YhxD9luaRo2rRpuLi4sHjxYsM1d3f3dHUKFSqUrrxhwwbq1q1L0aJFn9u2RqPJ8Fohco20VDj6PeyZDKmJYGoF9T6HSn1Bm/XnBurj43mweDEPFi5CiVdPpQ8uZ8/cSg8JKxAJaPCw82Bw+cHUdVWPHEpMSWPe3mB+2BdMUqoeE62GHtXdGVjPCxtZVSaE+JcslxRt3LiRgIAA2rVrx759+yhcuDD9+/enT58+T61/7949Nm/ezM8///yfbcfGxuLm5kZaWhrlypVj0qRJ+Pk9faggKSmJpL9Xr4C6+ZMQ2da9C7DhI7ir7teDR21oPhfyuRs1rBfxZBJ15Jy5pEZGAmBRtgymH/dmVPAQQO3pyWOah40t1fmOiqKw4+I9Jm26yJ2oBACqFSvAhOal8HK0McpzCCGyviyXFF2/fp158+YxZMgQRo8ezfHjxxk4cCDm5uZ069YtQ/2ff/4ZGxsbWrdu/dx2S5QowZIlS/D19SUmJoY5c+ZQvXp1zpw5g5eXV4b6U6dOZcKECZn2XEIYRWoSHJipfulTwdwOAiaD33uQDYaNYg8cJGL6dJKuXQPAtEgR7IcMJqq6DxOPTUlXt4pTFQCuR8Yy4Y+L7LuqJlBOdhZ81tSHJr6FZKhMCPFcmXrMR2YwMzPD39+fw4cPG64NHDiQEydOcOTIkQz1S5QoQYMGDfjmm29e6n30ej3ly5enVq1azJ07N8P9p/UUubi4yDEfIvu4cxI2DIDIS2q5eFNoOhNsnYwb1wtIvHyZiOlfEff37wGtnR02fXvws/c9tobuIiopCgBTrSmLAhZRPH9x9GmmfLsniJ8OXCclTcFMp6VPLQ8+quuJlVmW+/tPCPGWvPVjPjKTk5MTPj4+6a6VLFmSNWvWZKh74MABrly5wsqVK1/6fbRaLRUrVuTa33+B/pu5uTnm5tlnB18hDJLjYM8Udf6QogergtDkKyjVKsv3DqWEhxM5Zy7R69eDoqAxNSXfe+9xv31tupwYyf3r6h5EplpTfAr40KN0D8ral2XT2TAmb75EeIy6mrROcXvGvVsKj4Kyu74Q4sVluaSoevXqGVaIXb16FTc3twx1Fy5cSIUKFShbtuxLv4+iKAQGBuLr6/vKsQqR5VzdAZuHQnSIWi7TARp9CVb5jRvXf0iLjePBwp94uHgJyt/bZNg2aYL9kMEk2NvQZ01jHqc8BmBe/XlUKlQJM50Zl8Nj6LzgGEeuPwDAJb8lnzcrRf2SDjJUJoR4aVkuKRo8eDDVqlVjypQptG/fnuPHjzN//nzmz5+frl5MTAyrV69m5syZT22nW7duFC5cmKlTpwIwYcIEqlSpgpeXFzExMcydO5fAwEC+++67N/5MQrxxj+/BthFwYZ1atnOBpl+Dd0PjxvUflORkolat5v7335P28CEAlv4VcBw+nHivwsw69xOr9qwiWZ8MgI2pDVWdqhKbqGfKrgssPXqLNL2CuYmW/nU8+aB2USxMs/5KOiFE1pTlkqKKFSuybt06Ro0axcSJE/Hw8GD27Nl06dIlXb0VK1agKAqdOnV6ajshISFotf/fm/LRo0f07duX8PBw7Ozs8PPzY//+/VSqVOmNPo8Qb5ReD38tgZ3jISkaNFqo0h/qjgazrDt0pOj1PN62jYjZc0gJUXu1zNzdcfh0KHnq1eOP638weW0f4lPVpfde+bxoVrQZTdybsepkKF9tv8LDODVRaly6EKOblMQlv5XRnkcIkTNkuYnWWdXLTNQS4q2IuAR/fAK3j6plZz94dw44vfxw8tsUd/QoEV/NIPHCBQB09gWx/2gAedu0RmNqyubrmxl5YCQAPgV8GFR+EFWdqvJXyCPGb7zAudBoADwd8jD+3VLU8CpotGcRQmR92XqitRDiP6QkwP4ZcGgO6FPALA+881mW34Qx8fJlImZ+TdyBAwBora2xfL8Tu6taEZR0lps7N3Ir5hbRSWrSU6lQJRY0XMD9x8kMXXWGtadDAbAxN+GTBt50q+qGqS5LnlQkhMimJCkSIju5vhc2DYaH19Vy8SbqyjK7IkYN63lSQkOJnDuX6I1/gKKAqSn5OnQgsWszPjw1mpCLIRle42rjyodlP2bB/hvM3X2NuOQ0ANr7F2F4oxIUzCMrQ4UQmU+SIiGyg7gHsGMMnFmulm2coPF0KPlull1mnxoVxYMf5xO1bBlKSgqgriijb2cm3lnIgb3vAeBs7Uwb7za42brhbuuOq60rx4IfM/zXi1y/r/YOlXXJy4TmpSjnktdYjyOEyAUkKRIiK1MUNRHaPgYSHgIaqNRHHS6zsDN2dE+lT0zk4S9LebBgAfrH6jJ6qypVcBg6FEvf0ow8MJIDoeoQWo3CNRhTeQxFbNSerlsP4vh42QV2XVIPbi2Yx4wRjUrQpnwRtNqsmfwJIXIOSYqEyKoeBMOmT+DGfrXsWFqdSF3E36hhPYuSmkr0+vVEfvMtqffUpMa8RAnshw4hwteZ9eHHObbnJ3aH7Da85rt636HVaIlPTuX7PcHMP3Cd5L8Pbu1ezZ2B9b2wlYNbhRBviSRFQmQ1KQlw4Gs4NBvSksHEEuqMhKofgS7rJQiKovB4x04i58wh+bo618nU2Rn7TwZxqUJB+hweS8T1iHSvqexUmaEVhqJBw4bAUKZuuWzYjbqGZ0HGN/fB00EObhVCvF2SFAmRlVzZCluHw6O/Jx8Xq6eeV5bfw7hxPUPc4cNEzJpN4rlzAOjs7CjwYT/yde6M1syMk6e/JSJeTYgcLB3oWKIj1QpXo1SBUpwPjabdD0c4eUs9x6xIPks+a+pDQClH2Y1aCGEUkhQJkRU8vAHbRsLVbWrZtjA0mgolm2fJidQJZ88S8fUs4o+qeyRprKwo0L07+Xv2QJcnD5ceXGLvnb0cDD1oeM239b6lZIGS3I9NYuSas6w8eRtFAUtTHR/VLUbvmrIbtRDCuCQpEsKYUhLV/YYOfg2piaA1hWoDoNawLLkjdVJQEJFz5vB45y4ANKam5O3UkYIffIBJgQL8de8v5h+dz6HQQ4bXmGhNmFZzGp55i/PTgevM2X2Nx4mpALQo58zIxiVwsrM0yvMIIcQ/SVIkhLFc3QFbh0HUTbXsURuazAB7b6OG9TQpoaFEfvsd0Rs2qEeLaLXYtWiB/YCPMC1cGIANQRv47NBnhtfUc61HeYfy1Chcg5B7eWg0ez/BkXEAlC5sy/h3S+HvnrUPqhVC5C6SFAnxtkXdgu2j4fImtWzjBAFToFSrLDdUlvrgAfd/+JFHK1YY9hqyaVAf+0GDMPf0TFdXIf2JQSn6FGo6tuGL9RfZffkSAAWszRgWUJx2/i7oZIm9ECKLkaRIiLclNQkOz4X9MyE1AbQmUOVDqD0CzLPWSqu0x495uHgxD5b8jBKvHspqVbUKDoMHY1mmTLq6h0MPs/TSUo6FHTNcU9LMuHS1OA3/3EdKmmJYYv9xPS/sLLPeCjohhABJioR4O4J2wZbh8DBYLbvXVI/ncChp3Lj+RZ+YSNRvy3nw44+kRatnkFmULo3DkMFYV6sGQGhsKCfCT3Au8hynI09zLeqa4fWuNm44Ky05ebEw1+P0gEJtb3vGNvPB0yGPMR5JCCFemCRFQrxJ0Xdg2yi4tFEt53GEhpPBt22WGirTJyfzaPVqHvzwI6mRkQCYFSuG/aCB2DRoYFgify3qGu3/aE+qkmp4rYnWhI7FO1I6TzPm/xnFzjvRgB73AlaMbebDOyUcZIm9ECJbkKRIiDchOQ4OzVVXlqUmgEYHlfupmzBa2Bo7OgMlJYVH69dzf948Uu+GAerGiwU/+gi7Fs3RmJhw4f4FVl9dTZqSRlhcGKlKKrZmtrT2ak0Z+zIUtijFgr0RzDt9EwBrMx0f1/OiR3V3zE1kib0QIvuQpEiIzKTXw/nfYec4eHxXveZaTR0qK1TauLH9g5KWRsymTUR+9z0pIepGkSYODhT8sB9527RBY2ZmqPvVya84de9UuteXcyjHR2U/YcH+6wzce5aEFPUU+7YVijA8oDgOthZv72GEECKTSFIkRGa5fULdgDH0pFrO6woNJoFPiywzVKbo9Tzevp3Ib741HMmhK1CAgn37kLdDB7QWFmy7sY3PDn1GUloShfMU5l6ceo5Z91LdsTO3Q4sOk4SK1Ju5j9BHCQBUcMvHuHd9KFMkr7EeTQghXpskRUK8rug7sGs8nFutls3yQM2hUKU/mGaNHhNFUYj9808i535D0pUrgHokR/7evcjfpQtaKytD3YsPL5KUlgSok6oBClgUYKDfQK6ExzPxj4scv3kLACc7C0Y2LkHzss4yb0gIke1JUiTEq0qOU+cMHZqrzhtCA35d4J2xYFPI2NEBajIUd/AgkXPmknj+PADaPHnI36M7+d9/H12e/68Ii0+JZ8T+Eey9szddG1WdqjKs/EQ+W3eJVafUozksTLV8UKsY/WoXw9JM5g0JIXIGSYqEeFl6PZxbpfYOPVYnJ+NWXd2A0bmcMSNLJ+7oMSLnziXhr78A9Xyy/F27UqBHd3R58xrqxSbHcjriNIvOL+LkvZPp2lD0OhLu16DVN2eITVJXnDUvqx7N4ZxXjuYQQuQskhQJ8TJCjqnzhu6qiQZ53aDhpCx1cGv8yZNEfvvd/w9rNTcnX6dOFOjTG5MCBf5fLyWe3y7/xuLzi4lJjgHA2tSab9/5Fg87D/Zdecis7SEciEoEUvEtbMe4d33kaA4hRI4lSZEQL+JRiNozdH6NWjazgVpDofKHWWbeUNzx49z/7nvij/29s7SpKfnataPABx9g6uhgqJeYmkjjtY25n3DfcM3Z2plSBUvRy7cXJimuDFp2kYNB6n17G3OGBxSnTfkiaOVoDiFEDiZJkRDPk/BIPZrjyHfqKfZooHxXqPsZ2DgaOzoA4o4d5/533xF//Lh6wdSUvK1bU7BvH8Nhrf807cQ0Q0JUOE9h+pfrT1OPpkQnpDFr51V+O36ANL2CmU5L75oe9K/rSR5z+VUhhMj55DedEE+THA/Hf4SDsyHxkXrNvaY6b8ipzPNe+VYoikL8sePc//Zb4k+q84A0pqbYtW1DwT59MHV2zvCanbd2MvXYVENCVMi6EH+0/ANF0bH40E3m7L7G40R13lBAKUfGNPHBtYBVhnaEECKnkqRIiH9KS4G/foF90yE2XL1mX0JdUVaiqdHnDSmKQvzRo0R+9x0JJ9UNFTWmpuRt15YCffpg6uSU4TVp+jTmnZnHj2d/THf9+3e+Z9+Vh0zecokb9+MAKOlky+fNfKha7H/t3XlclNX+wPHPsA2LMCiI7CjuC25grrmkiWk/U8sttzbL2y2XyquWXbVSq9vtamXd9JrLLZe6qJmpuSsmqaG4i7iCIinKKtvAnN8fo4Mji5DoDPh9v17zgjnPd5453wFmvjzPec7xKLIfIYSo6qQoEgKMV5QdWwXbPoCUc8Y2XSB0exuaDwIby152rpQiKyqKq1/MK7yazMEB94ED8Rj9EvbexU8BsD1+O3MPzOVMmnEh2uGNh9PUsylp6S5MX3WNX0+fAsCzmgMTwxvyTGgAtjJuSAjxkJKiSDzclIK4TbD1ffjjiLHNpSZ0ngihz4Gd1sLdU9z4dQ/J8+aRffAgcLMYGjTIWAzVKjquSSnF4eTDrDy5kp/O/gSAq4MrUx6ZQvtaPfl08ylW7IvHoLJwsLXhxUfr8GrXurg62j/Q3IQQwtpIUSQeXheiYOsMiI8y3te6QYex0O4voK1W+mPvs1uTLiZ/MY/sQ4cA46X17oMH4fHiS2ZXk91yOuU0a8+uZVv8Ni6kXzC1P9/0eUY1eZGI35PptnQHGTfnG+od4s2UJxoTUEPGDQkhBEhRJB5GSUdg63vGI0QAdo7wyMvQaQI4W3YOHmUwkLF5C9e+/pqc48cB0Dg6Un3wYGq8+AL2XkWLIYD9SfsZs3kMeYY8AJzsnHgs8DEGNRjE1WRfBsw7wIVrWQA083Pj3T5NaBss44aEEOJ2UhSJh8e1M7B9lnEVewCNrfHy+i6TwK3o1VoPktLrSVv3M9cWLDAt1KpxcqL64MF4vPgCdjVrFvu4Xy/9yrKTy/g96XfyDHmE1QpjYIOBdA3oyvmr+bz/43GizhoHZMt8Q0IIUTopikTVl3wafp0Dh5aDwXjqiGZPQ7d3wKOuRbtmyMkhNSKCawsXkp9oXDLExs2NGsOHUX3ECOyqVy/xsUopPvn9E06nngagrXdb5vWYR3oWvP9TLCv2G9cpc7Cz4eVHgxnTta7MNySEEKWQd0hRdV0+BJGfwvEfAWVsq/c4dH8XfFpYtGsFmZmkLF/O9cVLKLh2DQBbT088nhuF+5AhZgu13ikxM5GFRxay8+JO/sj6A4AA1wA+6TyX/+xK4Mvtp7mRVwDAk819mPxEI/yry7ghIYS4GymKRNWiFFz41VgMndla2N6gF3R6AwLbWq5vQH5KCteXLiXlu2UY0o3rjdn7+lLjpRdxHzAAG8eSlwxJz0vn8wOf87+4/5F/84iXncaOxh5NeFQ3nj5z93IpNRuAFv463n1S1ikTQojykKJIVA0GA8T9YiyGLt5c7kJjYzxN1mkC1Gpq0e7pk5K4vmgRKd//gMo2Fi4OwcF4vDwaXZ8+aOxLvhxeX6Bnw/kNzD8833RVWVuftjzf9Hk0OcH845ezfBR/BQAfnSOTejWibwtfGTckhBDlJEWRqNwK8o2TLu7+F1wxXq2FrRZaDTNeXl+jjkW7l3fhAtf+8x9S1/wIej0Ajk2b4vHKy7j26IHGxuau+/g0+lO+PfEtAD4uPrzf8X18tSF8tPEk6w4bJ3J0drDlL13q8tKjwTg5WHaiSSGEqKykKBKVkz4HYr6FXz+D1Jtz8ji4QpsXoN2r4Fr8DM8PStaBg1xf9A0ZW7YaT+kBzmFheIwZg0vHDmjKsFzI6ZTT7Li4g+0J2wF4os4TjG3+N777LZmFu3eSl29Ao4FBoQG82bMBXm4ln3oTQghxd1IUicolJw32L4TfvoIbxlNGOHsYJ1xsMxqc3C3WNVVQQMaWrVxftIjsmBhTu0uXzni+/DLOoaFl2s+JayeYf3g+W+K3FO5b2eB4ozv9vjhIcqZxLqL2wR5MfbIxTX11FZqHEEI8rKQoEpXDlZNwYAkc/A5y04xtugDo8Dq0GgEOlru6ypCVReqq1VxfsgR9QgJgXKTV7am+eDz3HNp69cq0n10Xd/HJ759wLs249poGDZ39O1NNH8buw7VYetI4uDrY04W3ezeme2OvMh1xEkIIUTZWWRRdunSJSZMmsWHDBrKzs2nQoAELFy4k9OZ/2iV9EHz88cdMnDixxP1GRETw7rvvcubMGerWrcvMmTPp37//fclBVIC8LDi+BqKXQMJvhe2eDYyDp0MGgq3l1uvSX7lCynfLSFmxAkOasVCz1elwHzqEGsOGlTjhYnF+PvszkyMnA8YrynrW7km47yj+G5nJutirgAGdkz3je9RnWNsgHOzuPhZJCCFE+VhdUZSSkkLHjh3p1q0bGzZswMvLizNnzuDu7m6KuXz5stljNmzYwIsvvsjTTz9d4n6joqIYPHgw77//Pv3792f16tUMGjSI3bt307atZS/TFndIOmIshA5/X3hUSGNrvKw+dJRxrqEyDFC+X3Lj4ri2eDHpa39C3Rw8bR8YSI1RI3Hv3x8b57IdtdIX6NmasJXFRxdz7NoxwDjf0LyuS1m46zKjN56nwKCws9Ewsn1txnavh7uzw33LSwghHnYapW6OArUSkydP5tdffyUyMrLMj+nXrx8ZGRls3bq1xJjBgweTnp7Ohg0bTG29evWievXqLF++/K7PkZ6ejk6nIy0tDTc3tzL3TZRRbiYcjTCeIrsUXdjuHgitR0LL4eDmY7HuKaXI2ruXa998w41dhb+bTq1aUeP553Dt3h2Nbdmu+tpyYQvLTy4nLiWOlNyUwucw2PG8/1IW7U4k8+airY83qcWUJxoRXNOyC9QKIURlVZ7Pb6s7UrR27VrCw8MZOHAgO3fuxM/Pj1dffZXRo0cXG//HH3/w888/s2TJklL3GxUVxYQJE8zawsPDmTNnTrHxubm55Obmmu6n35xoT1SwxIMQvRiO/A/yMo1tNnbQqA+0HgXB3Sx6VMiQm0vGxo1cW7yE3BMnjI0aDa49elDj+edxbt2qXPsbs3kMvyb+arpf06kmT9XthyGjFT/8lsfnsfGAcdHWd3o3oX1dWbRVCCEeFKsris6ePctXX33FG2+8wdtvv82+ffsYO3YsWq2WkSNHFolfsmQJrq6uDBgwoNT9JiUlUatWLbO2WrVqkZSUVGz87NmzmTFjxp9PRJQsJw2O/GA8RZZ0uLC9RrCxEGr5LFQrfjX4B0WfmEjKipWk/vADBSnGozkaJyfc+/enxqiROAQFlWk/BmVgW/w24lLi0NppzQqiOV3n4FLQnNnrYzl00Xia0EfnyMTwhvRr6SeTLwohxANmdUWRwWAgLCyMWbNmAdCqVSuOHTvGV199VWxR9M033zBs2DAcS1ke4ZY7B2grpUoctD1lyhTeeOMN0/309HQCAgLKk4q4XdoliNsEcZvhzDbIN87qjK0DNO5rHCtU+1Gw4NVUSimyoqK4vmwZmdu2G2fJBuy8vak+ZAjugweVukBrcY4lH2PCjglF2qeHzeOHXS5sPGacfdvFwZa/dK3Li51k8kUhhLAUqyuKfHx8aNKkiVlb48aNiYiIKBIbGRlJbGwsK1euvOt+vb29ixwVunLlSpGjR7dotVq0Wm05ei7MGArg4n449YuxEPrjiPl2z4bGQqj5EHCx7CmigsxM0lavIWX5cvLOnjW1O7drR/Vhz+LarRsau7L/qRiUgYsZF4lNieW3xN/MtjVwa4lTxtNMXpaJviADGw0MbhPIhMfr4+Uqky8KIYQlWV1R1LFjR2JjY83aTp06RVAxpytuXabfosXdVzxv3749mzdvNhtXtGnTJjp06HDvnRZGN64ZF2E99Yvxa3bKbRs14B8G9cOh/uPGVeotPMdOblwc15ctI+3HtaisLABsnJ3R9etH9WeHlnl+oVtir8cybc8005Vkt1PKlkA1nLjDIaRn5wOKLg1q8nbvxjT0dq2IdIQQQtwjqyuKJkyYQIcOHZg1axaDBg1i3759zJ8/n/nz55vFpaen88MPP/DPf/6z2P2MHDkSPz8/Zs+eDcC4cePo3LkzH330EU899RQ//vgjW7ZsYffu3fc9pypLKeOYoLhNcGoTXPodlKFwu6MO6vUwFkL1uoOLp+X6epPS68nYuo2U774ja/9+U7tD3bpUH/Ysur5PYVvN5U/tO/JSpFlB1LhGY5p6NMMptz0/7bfh+PUcIJ9G3q683bsxnRuUfR4jIYQQ95/VFUVt2rRh9erVTJkyhffee486deowZ84chg0bZha3YsUKlFIMHTq02P3Ex8djc9tVSx06dGDFihVMnTqVd999l7p167Jy5UqZo6islIL0REiOhauxkHTUeDQow3zOKGo1Mx4Jqh8O/m3A1jp+xfR/XCE14n+krlhJ/pWby4PY2uL62GNUHzYM57aP/OnZoa9mXeXHMz+y/ETh1A5aWy2Tmn/Fhxti2X/eeMSspquWt3o24JnQAGxlELUQQlgdq5unyFo9NPMUGQogNd5Y+Fw9CcmnjF+vnoK8jKLx9s4Q3PVmIdQTdP4PvMslMeTmkrl9O6mrVnFj96+mgdO2Hh64D3yG6oMHY+/z5+c+up5znXkH5xERF0GBKgDAxd6Fdp5PkprYmR0nja+Xo70NLz8azCtd6uKitY4iUQghHhaVep4icR8pBfm5xvmA8jIhNwOunzUWPFdPGo8CJcdBfk7xj9fYGi+br9nQeAvqAEGdwN56Bggrpcg5eoy01atI+3m9afkNAKfWrak+dCiu4T2xcbi3maF/PvszM3+bSYbeWPi0rNmSXoFPc+psMMt2XkRfkIFGA8+09ufNng3x1lnPaySEEKJ4UhRZmqEAfv/G+L1SgCr6tbRtZoXOjcKvuRl33L9ZCN08olEqWy141jcWPp4NC4ugGnXBzjqXmci/epW0tT+RtmY1uXGnTe123t7onnoKXb+n0NapUyHPte7sOqZETgGM44YmtJrIsXMefBxxmrRs44Kwj9b3ZMoTjWniW4WPKgohRBUjRZGlGQpg/VsP/nntncHBxbjS/K2i51YBVL022Fj/XDkqL4+MHTtIW7WazMhIKDAWfBqtFtcePdAN6I9Lu3ZlXn7jbhYcXsCKkyu4km0ck6QUDAv4lL8tO0XCdWNbI29XpvRuTBcZRC2EEJWOFEWWprGBJk/dunPzMvVivkLJ2+wcwKGa8aatZix2bt13cLnZdtt9B5dKUfQURylF7okTpK5aTfq6dRSkppq2ObVogW7AANye6IVtBY/7up5znS9jviRfGdckczeEUpD8FONWHALAy1XLWz0b8nSovwyiFkKISkqKIkuztYNBSy3dC6uXe/YsGZu3kL5+Pbm3zWNl5+VlPD3Wvx/a4OAKf97Y67H8cOoH1p1dR77Kx8+hNf55r7IlNhkAZwdbxnSpy0uP1sHZQf6chBCiMpN3cWGVlFLkHDtOxubNZGzZQt6ZM6ZtGgcHXHt0R9e/Py4dOlTY6bFbsvOz2XhuI/879T8OJxvXZjPku+CQ9ixxyc05aUiWmaiFEKIKkqJIWA1VUEBWdDQZm7eQsXUL+Ym3zYFkb49Lu3a49uiBW69wbHW6Cn/+5Oxk5h+ez7oz60xXldkqR/wNIzhzoR4pecZB74818mLyE41oUEtmohZCiKpEiiJhUYa8PG7s2UPGli1kbt1mWpEejKvSV+vcGdcePajWtQu2rvenCMnSZxF5KZK3dhYOePdzCaCB/TD2HvPiaHoeoGjq68Y7vRvToZ7lZ+YWQghR8aQoEg9cQeYNbuzaaSyEdu7CcOOGaZutTke1xx7D9fEeuHTogI3j/T01FXs9llEbR3FDX9iHnl7jOH6qAWsvZwB5+OoceSu8If1a+mEjg6iFEKLKkqJI3HcFGRlkHzxI1u/RZEVHk3P4MEqvN223q1UL1+7dce35OM5hYeVakb68LmVeYlv8Ns6knuFi5kX2Xt5r2vak/xhOnmpExM4cIANXRzte61aPUR1q42hfOa/WE0IIUXZSFIkKl3/1KlnR0aYiKDc21rTExi0OQUG49nwc18cfx7FZMzS3rVNX0bL0Wfxw6gc2ntvI0WtHi2w36HVUy3iWFSeDUCoHe1sNI9vX5rVu9ajuYp2TVQohhKh4UhSJe6KUQh8fbyqAsqJ/R38hvkicfWAgzqGhOIeF4hwain1Q0J9egLW8FhxZwH+O/AcADRpCa4US5h1GDXt/fjvuxoaYbP4oMA6i7tvCl4nhDQmo4fxA+iaEEMJ6SFEkykzp9egvXSIvIYG8c+fIOnCQrOjfKbiabB6o0aBt2NBUBDm1DsW+lpdlOg1EJUaZvl/bby0+LgF8+1s8H22LIzUrC4C2dWrwdu/GtAhwt1AvhRBCWJoURcKMISvLWPTEx6OPv/k1IZ68+AT0iYlFToMBaOztcQwJKSyCWrWq8Bmly8ugDCRkJHAs+RjHrh0DQCkNH27bzLG4uiRczwagvlc1pvRuRLeGXg/syJUQQgjrJEVRFaWUAr0eQ3a28ZaVjSE7C3X7/aws9JcTTcVPXkJ80aM+d9A4OuIQGIh9YABOTZviHBaGY0jIfb9KrDxWxa3ik98/ISMvw9SWf6MOhuS+bMzyAbLxctXyZs8GPN3aHzvb+zeeSQghROUhRZGFqfx8Lr4+1ri6qFIo1M0NCtTNr6b7Cm5uV7dtV3l5N4ufLFRWtqkQurVAannZ6nTYBwaaih+HgEAcggKxDwjArmZNqz6isvLkSj7Y+wEAWlstPvahZFx+jPNJ1QBwubksx4uyLIcQQog7yKeCpSlF5vbt9/c57OywcXIy3TQuztg4OWPj5ISdlxcOgQHGAiggEIfAgPsyW/T9dCH9Alvjt7I9fjuHrhoXaB1SdwzXEzuw6mAiSoGdjYahjwQytnt9arpqLdxjIYQQ1kiKIkuztcXng/eN399a9f7W97eOyGgoPDpze7txAxp7e2ycbxY8Tk7Ggsf5tiLIoepdVq6UYk/iHpYeX8qexD2F7QVO1Da8zNKNfuTmJwLQO8SbieGNqOPpYqnuCiGEqASkKLIwjY0N7s88Y+luVCqZeZm8tOkl0wBqgEe8OqK98Ti7j7pwNKcAMPBInRpMeaIRrQKrW66zQgghKg0pikSlkW/IJy4ljg3nNphdUdajxhT2H6rFpdRsoIAGtaox+Qm5okwIIUT5SFEkrFZKTgqHrx7m0NVDxFyN4WjyUbLzjZfSKwUFN+rjkDqINSddgWy83Rx54/EGPB3qj62sUSaEEKKcpCgSVukf+//B0uNLi7RXs69GkEMnLie048IfjmQDro52/KVrXV7oWEfWKBNCCPGnSVEkrNKVrCtF2t5r8yVbY1xY9/tlABxsbRjZPoi/yhplQgghKoAURcIqKKVIvJHIiWsnOH7tOGm5aaZthnwXXG8MZOJ3GeQXpAPQr6Uvb/aUNcqEEEJUHCmKhEVFnIpg1elVnE87T3peutk2VeCAfUYvspLbk6TXAIpH63syqVcjmvlVrrmUhBBCWD8pioTFGJSB6VHTTfftbOyo716fhtWbkn41hB1HnEjJMq61FuKnY/ITjehYz9NCvRVCCFHVSVEkLGZnwk7T9w42DuwZEsWm48l8simWC9eyAANBHs5MDG9I72Y+2MgVZUIIIe4jKYrEA5WWm8aWC1tYd3Ydv//xOwCBroEMDZrB0//ey9FLxlNontW0jOtejyGPBGIvC7YKIYR4AKQoEg/EzoSdRMRFEHkpknxDPgA2Ght6+77ChXOtmLrvD8C4YOsrXeryYqc6uGjl11MIIcSDI5864r5SSvH5wc9ZcGSBqa1B9Qa08/w/Yk/XZ/mWVOA69rYahrcL4rVu9fCoJgu2CiGEePCkKBL3RXx6PBvObWDDuQ2cSTsDwNBGQ+nh15+fowv4el08+YZUNBp4qoVcXi+EEMLypCgSFcqgDEzcOZFNFzaZ2rS2Wsa3nMyVy8157uezZOUVANClQU3+1qshTX3l8nohhBCWJ0WRqDBn086y6OgiNl3YhI3GhnY+7egR8ARXkxry6eqLXL8RB0ALfx2TejWig1xeL4QQwopIUSTuSb4hn50JO1keu5y9l/ea2t8K/RuO2Z2ZsyaOS6lnAQiu6cLEng3p1cxbVq8XQghhdaQoEn+KUop1Z9fx2cHPSLqRBBivJuvs14UGDk+z9BeIu3IYAG83R8b3qM8zof7YyeX1QgghrJQURaLc0nLTmLl3JhvObQCgurY6A+oPoL5jbxbuvMpP8akA6JzsebVrXUZ1qC2r1wshhLB6UhSJMtMb9Hwf+z1fHfqKtNw0bDW2jGkxhvaezzBn01nmxBrHDDna2/BCxzq80qUuOid7C/daCCGEKBspikSZzfxtJhFxEQDUc6/HmCZT2RBtw0eH9qIU2NpoGNImgHHd6+Pl5mjh3gohhBDlI0WRuKuc/BwWHFnAqrhVAIxtPpX4C014bXEC+gIFwJPNfXizZ0PqeLpYsqtCCCHEn2aVo14vXbrE8OHD8fDwwNnZmZYtWxIdHW0Wc+LECfr27YtOp8PV1ZV27doRHx9f4j4XL16MRqMpcsvJybnf6VRqB/44QL8f+zH/8HwMBQ7ULRjPv1br+O9v8egLFI/W92Td65344tnWUhAJIYSo1KzuSFFKSgodO3akW7dubNiwAS8vL86cOYO7u7sp5syZM3Tq1IkXX3yRGTNmoNPpOHHiBI6OpZ+ycXNzIzY21qztbo95mF3OvMzr214nLScbbWYfspMfJSYXoIAWAe5MCm8ocw0JIYSoMqyuKProo48ICAhg0aJFprbatWubxbzzzjv07t2bjz/+2NQWHBx8131rNBq8vb0rrK9VVUJ6Av898V9Wx60lPbkphuu9yMwzLsFR36sab/ZsSHjTWjLXkBBCiCrF6k6frV27lrCwMAYOHIiXlxetWrViwYLCxUQNBgM///wzDRo0IDw8HC8vL9q2bcuaNWvuuu/MzEyCgoLw9/fnySef5ODBgyXG5ubmkp6ebnaryjLzMvnx9I+M2TyGPqv/j6V7j5Ec+xdykwagz3PGz92JTwa2YOP4zjL5ohBCiCpJo5RSlu7E7W6dznrjjTcYOHAg+/btY/z48Xz99deMHDmSpKQkfHx8cHZ25oMPPqBbt25s3LiRt99+m+3bt9OlS5di9/vbb79x+vRpQkJCSE9PZ+7cuaxfv55Dhw5Rv379IvHTp09nxowZRdrT0tJwc3Or2KQtbHXcambvm02WPpuCGw3JvRKOIdcXAA8XB157rB7Ptg1EaydzDQkhhKhc0tPT0el0Zfr8trqiyMHBgbCwMPbs2WNqGzt2LPv37ycqKorExET8/PwYOnQoy5YtM8X07dsXFxcXli9fXqbnMRgMtG7dms6dO/PZZ58V2Z6bm0tubq7pfnp6OgEBAVWqKNIX6Plo/0esjF1JflYQmutPkZVhLIZctXaM7hzMC53qUE1rdWdZhRBCiDIpT1FkdZ92Pj4+NGnSxKytcePGREQY58fx9PTEzs6u2Jjdu3eX+XlsbGxo06YNcXFxxW7XarVotdpy9t76GZSBbfHb+D72e45eO0pqejXyrj5HfmYjALR2NozqUJu/dKlLdRcHC/dWCCGEeHCsrijq2LFjkSvETp06RVBQEGA8ktSmTZtSY8pCKUVMTAwhISH33ulKYn/Sfj747QPOpp3FkOdB7tXe5Ke3BDTY2mgYFObP2O718dE5WbqrQgghxANndUXRhAkT6NChA7NmzWLQoEHs27eP+fPnM3/+fFPMxIkTGTx4MJ07dzaNKfrpp5/YsWOHKWbkyJH4+fkxe/ZsAGbMmEG7du2oX78+6enpfPbZZ8TExDBv3rwHneIDpzfomXdwHt8c/YYCvSvq+iByUlphUMbB0n2a+/Dm4w0IrlnNwj0VQgghLMfqiqI2bdqwevVqpkyZwnvvvUedOnWYM2cOw4YNM8X079+ff//738yePZuxY8fSsGFDIiIi6NSpkykmPj4eG5vCi+tSU1N5+eWXSUpKQqfT0apVK3bt2sUjjzzyQPN70PIK8piwYwI7zv9O3rUnMKR2osBgfF26NqzJWz0b0sxPZ+FeCiGEEJZndQOtrVV5BmpZg2vZ19iRsIPvT6zlYFwN9Nc7oQzGK/va1K7OxPBGPFKnhmU7KYQQQtxnlXqgtbg3Ofk5LDm2hP8cXkL61ZbkXeuDKjAuv9HU142J4Q3p0qCmzDMkhBBC3EGKoiribNpZfjrzE+tObyA+0Y+85NdQ+cbTYgE1HJjcqxlPNPPGxkaKISGEEKI4UhRVcnqDnoVHFvLvQ/PJSW1C7tWhKL1xPTJfd0fG92jAgFZ+2Nla3eTlQgghhFWRoqgSS89LZ9jPwzl9yZHcq69iyPUBwMPFntcfq89QmYVaCCGEKDMpiiohfYGebQnbeGfTEq5f7oUh2zg/k5ujHa90qctzHWrjIrNQCyGEEOUin5yVSHx6PEeSjzB9+1KuX+pAQZZxmgIHO3ipU11e6VwXnbO9hXsphBBCVE5SFFUCf9z4gw/2fsDWuFhyr/akIHMkALY2iiGP+DOueyO8XB0t3EshhBCicpOiyEolZyfjaOuIva09r6yfxtFTQeRn9AJAo1E8E+rP+B4N8XOXJTmEEEKIiiBFkZVRSrH70m7e3PkmN7IcyU3uQX5aH8AGDdC9qTtv92ohS3IIIYQQFUyKIiuilGLugbksOPg9edd6ok9pw60fUZtgB97v25ZG3tY/m7YQQghRGUlRZEWu38hjZ0wNbpz5GyjjgOl6Pnqm9mlJ13q1Lds5IYQQooqToqi8btwA22Lm/rG1BUdH87iS2NiAU+FYoPRrqSz+9RxL91zgRp4LThRQzTWJt59oSf8WLcxiycqCkpar02jA2fnPxWZng8FQcp9dXP5cbE4OFBRUTKyzs7HfALm5kJ9fMbFOTsafCUBeHuj1FRPr6Fj4u1KeWL3eGF8SrRbs7Mofm59vfC1K4uAA9vbljy0oMP7sSmJvb4wvb6zBYPxdq4hYOzvjawHGv4msrIqJLc/f/T28R5QrVt4jjN/Le0T5Y6vye0RZKVEmaWlpClBpxreQorfevc0f4OxcfBwo1aWLKezk5XR1zVlXcmxYmPl+g4JKjm3SxDy2SZOSY4OCzGPDwkqO9fQ0j+3SpeRYZ2fz2N69S46989fvmWdKj83MLIwdNar02CtXCmNffbX02HPnCmPfeqv02KNHC2OnTSs9dt++wtiPPy49dvv2wtgvvig9dt26wthFi0qP/f77wtjvvy89dtGiwth160qP/eKLwtjt20uP/fjjwth9+0qPnTatMPbo0dJj33qrMPbcudJjX321MPbKldJjR40qjM3MLD32mWeUmdJi/+R7hFLK+DdYUqy8RxTe5D3CeJP3COPt5nuE6fM7LU3djaz9YGF1a7ogy5EJIYQQlqdRSilLd6IySE9PR6fTkZaYiJtbMYOd7+HQ+KWLV/HWOWFbXHUkh8YLyaFxIzk0Xv5YOX1mJO8Rfy5W3iOMKul7hOnzOy2t+M/v20hRVEbleVGFEEIIYR3K8/ktp8+EEEIIIZCiSAghhBACkKJICCGEEAKQokgIIYQQApCiSAghhBACkKJICCGEEAKQokgIIYQQApCiSAghhBACkKJICCGEEAKQokgIIYQQApCiSAghhBACkKJICCGEEAKQokgIIYQQAgA7S3egslBKAcbVdoUQQghROdz63L71OV4aKYrKKCMjA4CAgAAL90QIIYQQ5ZWRkYFOpys1RqPKUjoJDAYDiYmJuLq6otFoKnTf6enpBAQEkJCQgJubW4Xu29o8TLnCw5Xvw5QrPFz5Pky5wsOV78OQq1KKjIwMfH19sbEpfdSQHCkqIxsbG/z9/e/rc7i5uVXZX8o7PUy5wsOV78OUKzxc+T5MucLDlW9Vz/VuR4hukYHWQgghhBBIUSSEEEIIAUhRZBW0Wi3Tpk1Dq9Vauiv33cOUKzxc+T5MucLDle/DlCs8XPk+TLmWhQy0FkIIIYRAjhQJIYQQQgBSFAkhhBBCAFIUCSGEEEIAUhQJIYQQQgBSFFncl19+SZ06dXB0dCQ0NJTIyEhLd6ncZs+eTZs2bXB1dcXLy4t+/foRGxtrFqOUYvr06fj6+uLk5ETXrl05duyYWUxubi6vv/46np6euLi40LdvXy5evPggUym32bNno9FoGD9+vKmtquV66dIlhg8fjoeHB87OzrRs2ZLo6GjT9qqSb35+PlOnTqVOnTo4OTkRHBzMe++9h8FgMMVU5lx37drF//3f/+Hr64tGo2HNmjVm2ysqt5SUFEaMGIFOp0On0zFixAhSU1Pvc3ZFlZavXq9n0qRJhISE4OLigq+vLyNHjiQxMdFsH5Ul37v9bG/3yiuvoNFomDNnjll7Zcn1vlPCYlasWKHs7e3VggUL1PHjx9W4ceOUi4uLunDhgqW7Vi7h4eFq0aJF6ujRoyomJkb16dNHBQYGqszMTFPMhx9+qFxdXVVERIQ6cuSIGjx4sPLx8VHp6emmmDFjxig/Pz+1efNmdeDAAdWtWzfVokULlZ+fb4m07mrfvn2qdu3aqnnz5mrcuHGm9qqU6/Xr11VQUJB67rnn1N69e9W5c+fUli1b1OnTp00xVSXfDz74QHl4eKh169apc+fOqR9++EFVq1ZNzZkzxxRTmXNdv369euedd1RERIQC1OrVq822V1RuvXr1Us2aNVN79uxRe/bsUc2aNVNPPvnkg0rTpLR8U1NTVY8ePdTKlSvVyZMnVVRUlGrbtq0KDQ0120dlyfduP9tbVq9erVq0aKF8fX3Vv/71L7NtlSXX+02KIgt65JFH1JgxY8zaGjVqpCZPnmyhHlWMK1euKEDt3LlTKaWUwWBQ3t7e6sMPPzTF5OTkKJ1Op/79738rpYxvUvb29mrFihWmmEuXLikbGxu1cePGB5tAGWRkZKj69eurzZs3qy5dupiKoqqW66RJk1SnTp1K3F6V8u3Tp4964YUXzNoGDBighg8frpSqWrne+cFZUbkdP35cAeq3334zxURFRSlAnTx58j5nVbLSCoVb9u3bpwDTP6WVNd+Scr148aLy8/NTR48eVUFBQWZFUWXN9X6Q02cWkpeXR3R0ND179jRr79mzJ3v27LFQrypGWloaADVq1ADg3LlzJCUlmeWq1Wrp0qWLKdfo6Gj0er1ZjK+vL82aNbPK1+Ovf/0rffr0oUePHmbtVS3XtWvXEhYWxsCBA/Hy8qJVq1YsWLDAtL0q5dupUye2bt3KqVOnADh06BC7d++md+/eQNXK9U4VlVtUVBQ6nY62bduaYtq1a4dOp7Pq/MH4vqXRaHB3dweqVr4Gg4ERI0YwceJEmjZtWmR7Vcr1XsmCsBaSnJxMQUEBtWrVMmuvVasWSUlJFurVvVNK8cYbb9CpUyeaNWsGYMqnuFwvXLhginFwcKB69epFYqzt9VixYgUHDhxg//79RbZVtVzPnj3LV199xRtvvMHbb7/Nvn37GDt2LFqtlpEjR1apfCdNmkRaWhqNGjXC1taWgoICZs6cydChQ4Gq97O9XUXllpSUhJeXV5H9e3l5WXX+OTk5TJ48mWeffda0KGpVyvejjz7Czs6OsWPHFru9KuV6r6QosjCNRmN2XylVpK0yee211zh8+DC7d+8usu3P5Gptr0dCQgLjxo1j06ZNODo6lhhXFXIF43+YYWFhzJo1C4BWrVpx7NgxvvrqK0aOHGmKqwr5rly5km+//ZZly5bRtGlTYmJiGD9+PL6+vowaNcoUVxVyLUlF5FZcvDXnr9frGTJkCAaDgS+//PKu8ZUt3+joaObOncuBAwfK3afKlmtFkNNnFuLp6YmtrW2RCvvKlStF/lurLF5//XXWrl3L9u3b8ff3N7V7e3sDlJqrt7c3eXl5pKSklBhjDaKjo7ly5QqhoaHY2dlhZ2fHzp07+eyzz7CzszP1tSrkCuDj40OTJk3M2ho3bkx8fDxQtX62EydOZPLkyQwZMoSQkBBGjBjBhAkTmD17NlC1cr1TReXm7e3NH3/8UWT/V69etcr89Xo9gwYN4ty5c2zevNl0lAiqTr6RkZFcuXKFwMBA03vWhQsXePPNN6lduzZQdXKtCFIUWYiDgwOhoaFs3rzZrH3z5s106NDBQr36c5RSvPbaa6xatYpt27ZRp04ds+116tTB29vbLNe8vDx27txpyjU0NBR7e3uzmMuXL3P06FGrej26d+/OkSNHiImJMd3CwsIYNmwYMTExBAcHV5lcATp27FhkeoVTp04RFBQEVK2fbVZWFjY25m+Jtra2pkvyq1Kud6qo3Nq3b09aWhr79u0zxezdu5e0tDSry/9WQRQXF8eWLVvw8PAw215V8h0xYgSHDx82e8/y9fVl4sSJ/PLLL0DVybVCPOiR3aLQrUvyFy5cqI4fP67Gjx+vXFxc1Pnz5y3dtXL5y1/+onQ6ndqxY4e6fPmy6ZaVlWWK+fDDD5VOp1OrVq1SR44cUUOHDi32cl9/f3+1ZcsWdeDAAfXYY49ZxaXMd3P71WdKVa1c9+3bp+zs7NTMmTNVXFyc+u6775Szs7P69ttvTTFVJd9Ro0YpPz8/0yX5q1atUp6enupvf/ubKaYy55qRkaEOHjyoDh48qAD16aefqoMHD5qutqqo3Hr16qWaN2+uoqKiVFRUlAoJCbHIZdul5avX61Xfvn2Vv7+/iomJMXvfys3NrXT53u1ne6c7rz5TqvLker9JUWRh8+bNU0FBQcrBwUG1bt3adBl7ZQIUe1u0aJEpxmAwqGnTpilvb2+l1WpV586d1ZEjR8z2k52drV577TVVo0YN5eTkpJ588kkVHx//gLMpvzuLoqqW608//aSaNWumtFqtatSokZo/f77Z9qqSb3p6uho3bpwKDAxUjo6OKjg4WL3zzjtmH5KVOdft27cX+3c6atQopVTF5Xbt2jU1bNgw5erqqlxdXdWwYcNUSkrKA8qyUGn5njt3rsT3re3bt5v2UVnyvdvP9k7FFUWVJdf7TaOUUg/iiJQQQgghhDWTMUVCCCGEEEhRJIQQQggBSFEkhBBCCAFIUSSEEEIIAUhRJIQQQggBSFEkhBBCCAFIUSSEEEIIAUhRJIQQQggBSFEkhBBWb+HChfTs2bNcj/niiy/o27fvfeqREFWTFEVCiDLRaDSl3p577jlLd7HCde3alfHjx1u0D7m5ufz973/n3XffNbVNnz6dli1bmsVFRkbi7u7O66+/jlKK0aNHs3//fnbv3v2AeyxE5SVFkRCiTC5fvmy6zZkzBzc3N7O2uXPnWrqLZabX6yvN80VERFCtWjUeffTREmN+/vlnwsPDGTduHJ9//jkajQatVsuzzz7L559//qefW4iHjRRFQogy8fb2Nt10Oh0ajcasbdeuXYSGhuLo6EhwcDAzZswgPz/f9HiNRsPXX3/Nk08+ibOzM40bNyYqKorTp0/TtWtXXFxcaN++PWfOnDE95tYRka+//pqAgACcnZ0ZOHAgqampZn1btGgRjRs3xtHRkUaNGvHll1+atp0/fx6NRsP3339P165dcXR05Ntvv+XatWsMHToUf39/nJ2dCQkJYfny5abHPffcc+zcuZO5c+eajoadP3+exYsX4+7ubvb8a9asQaPRFOn3N998Q3BwMFqtFqUUaWlpvPzyy3h5eeHm5sZjjz3GoUOHSn3dV6xYUeppsGXLljFgwAA+/PBDZsyYYbatb9++rFmzhuzs7FKfQwhhJEWREOKe/fLLLwwfPpyxY8dy/Phxvv76axYvXszMmTPN4t5//31GjhxJTEwMjRo14tlnn+WVV15hypQp/P777wC89tprZo85ffo033//PT/99BMbN24kJiaGv/71r6btCxYs4J133mHmzJmcOHGCWbNm8e6777JkyRKz/UyaNImxY8dy4sQJwsPDycnJITQ0lHXr1nH06FFefvllRowYwd69ewGYO3cu7du3Z/To0aajYQEBAWV+TW71OyIigpiYGAD69OlDUlIS69evJzo6mtatW9O9e3euX79e4n4iIyMJCwsrdtu8efN4/vnnWbhwIWPHji2yPSwsDL1ez759+8rcbyEeakoIIcpp0aJFSqfTme4/+uijatasWWYx//3vf5WPj4/pPqCmTp1quh8VFaUAtXDhQlPb8uXLlaOjo+n+tGnTlK2trUpISDC1bdiwQdnY2KjLly8rpZQKCAhQy5YtM3vu999/X7Vv314ppdS5c+cUoObMmXPXvHr37q3efPNN0/0uXbqocePGlZq7UkqtXr1a3f52Om3aNGVvb6+uXLliatu6datyc3NTOTk5Zo+tW7eu+vrrr4vtT0pKigLUrl27zNqnTZumHBwcirx+xalevbpavHhxqTFCCCM7SxZkQoiqITo6mv3795sdGSooKCAnJ4esrCycnZ0BaN68uWl7rVq1AAgJCTFry8nJIT09HTc3NwACAwPx9/c3xbRv3x6DwUBsbCy2trYkJCTw4osvMnr0aFNMfn4+Op3OrI93Hm0pKCjgww8/ZOXKlVy6dInc3Fxyc3NxcXG515cDgKCgIGrWrGm6Hx0dTWZmJh4eHmZx2dnZZqcM79wG4OjoWGSbv78/7u7ufPzxxzzxxBP4+PgUuw8nJyeysrL+bBpCPFSkKBJC3DODwcCMGTMYMGBAkW23f6Db29ubvr81Bqe4NoPBUOJz3YrRaDSmuAULFtC2bVuzOFtbW7P7dxY7//znP/nXv/7FnDlzCAkJwcXFhfHjx5OXl1dyooCNjQ1KKbO24gZS3/l8BoMBHx8fduzYUST2zjFKt3h4eKDRaEhJSSmyzdXVlS1bttCzZ0+6du3K9u3b8fX1LRJ3/fp1s+JMCFEyKYqEEPesdevWxMbGUq9evQrfd3x8PImJiaYP/KioKGxsbGjQoAG1atXCz8+Ps2fPMmzYsHLtNzIykqeeeorhw4cDxqIlLi6Oxo0bm2IcHBwoKCgwe1zNmjXJyMjgxo0bpsLn1pih0rRu3ZqkpCTs7OyoXbt2mfro4OBAkyZNOH78eLHzFFWvXp0tW7YQHh5uKoz8/PxM28+cOUNOTg6tWrUq0/MJ8bCTgdZCiHv297//naVLlzJ9+nSOHTvGiRMnWLlyJVOnTr3nfTs6OjJq1CgOHTpEZGQkY8eOZdCgQXh7ewPGK71mz57N3LlzOXXqFEeOHGHRokV8+umnpe63Xr16bN68mT179nDixAleeeUVkpKSzGJq167N3r17OX/+PMnJyRgMBtq2bYuzszNvv/02p0+fZtmyZSxevPiuefTo0YP27dvTr18/fvnlF86fP8+ePXuYOnWqaZB5ccLDw0uda0in07Fp0yY8PT3p2rUrFy9eNG2LjIwkODiYunXr3rV/QggpioQQFSA8PJx169axefNm2rRpQ7t27fj0008JCgq6533Xq1ePAQMG0Lt3b3r27EmzZs3MLrl/6aWX+M9//sPixYsJCQmhS5cuLF68mDp16pS633fffZfWrVubjrJ4e3vTr18/s5i33noLW1tbmjRpQs2aNYmPj6dGjRp8++23rF+/3nQZ//Tp0++ah0ajYf369XTu3JkXXniBBg0aMGTIEM6fP28aX1Wc0aNHs379etLS0kqMcXNz45dffqFWrVp07dqVhIQEAJYvX2421koIUTqNuvPkuBBCWInp06ezZs2aMp2eqsoGDRpEq1atmDJlSpkfc/ToUbp3786pU6eKDDoXQhRPjhQJIYSV+8c//kG1atXK9ZjExESWLl0qBZEQ5SADrYUQwsoFBQXx+uuvl+sx5V1AVgghp8+EEEIIIQA5fSaEEEIIAUhRJIQQQggBSFEkhBBCCAFIUSSEEEIIAUhRJIQQQggBSFEkhBBCCAFIUSSEEEIIAUhRJIQQQggBwP8DXJmgZj5YqLQAAAAASUVORK5CYII="},"metadata":{}}],"execution_count":14},{"id":"03887d0e-da24-49ca-9bd8-9452fd666b3c","cell_type":"markdown","source":"Both the [Moruzzi, V. L. et al.](https://link.aps.org/doi/10.1103/PhysRevB.37.790) and the quantum mechanical version of the quasi-harmonic approach start at a larger equilibrium volume as they include the zero point vibrations, resulting in an over-prediction of the volume expansion with increasing temperature. The equilibrium volume is indicated by the dashed red line. Finally, the quasi-harmonic approach with the classical harmonic oszillator agrees very well with the thermal expansion calculated from molecular dynamics for this example of using the Morse Pair Potential. ","metadata":{}},{"id":"e7aeb2f2-12c5-492a-a821-384b030b4a68","cell_type":"code","source":"","metadata":{"trusted":true},"outputs":[],"execution_count":null}]}