{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Fe-C diffusion and elastodiffusivity\n", "Taking data from R.G.A. Veiga, M. Perez, C. Becquart, E. Clouet and C. Domain, Acta mater. **59** (2011) p. 6963 [doi:10.1016/j.actamat.2011.07.048](http://dx.doi.org/10.1016/j.actamat.2011.07.048)\n", "\n", "Fe in the body-centered cubic phase, $a_0 = 0.28553\\text{ nm}$; C sit at octahedral sites, where the transition states between octahedral sites are represented by tetrahedral sites. The data is obtained from an EAM potential, where $C_{11} = 243\\text{ GPa}$, $C_{12}=145\\text{ GPa}$, and $C_{44} = 116\\text{ GPa}$. The tetrahedral transition state is 0.816 eV above the octahedral site, and the attempt frequency is taken as 10 THz ($10^{13}\\text{ Hz}$).\n", "\n", "The dipole tensors can be separated into *parallel* and *perpendicular* components; the parallel direction points towards the closest Fe atoms for the C, while the perpendicular components lie in the interstitial plane. For the octahedral, the parallel component is 8.03 eV, and the perpendicular is 3.40 eV; for the tetrahedral transition state, the parallel component is 4.87 eV, and the perpendicular is 6.66 eV." ] }, { "cell_type": "code", "execution_count": 1, "metadata": { "collapsed": false }, "outputs": [], "source": [ "import sys\n", "sys.path.extend(['../'])\n", "import numpy as np\n", "import matplotlib.pyplot as plt\n", "plt.style.use('seaborn-whitegrid')\n", "%matplotlib inline\n", "import onsager.crystal as crystal\n", "import onsager.OnsagerCalc as onsager\n", "from scipy.constants import physical_constants\n", "kB = physical_constants['Boltzmann constant in eV/K'][0]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Create BCC lattice (lattice constant in nm)." ] }, { "cell_type": "code", "execution_count": 2, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "#Lattice:\n", " a1 = [-0.142765 0.142765 0.142765]\n", " a2 = [ 0.142765 -0.142765 0.142765]\n", " a3 = [ 0.142765 0.142765 -0.142765]\n", "#Basis:\n", " (Fe) 0.0 = [ 0. 0. 0.]\n" ] } ], "source": [ "a0 = 0.28553\n", "Fe = crystal.Crystal.BCC(a0, \"Fe\")\n", "print(Fe)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Elastic constants converted from GPa ($10^9$ J/m$^3$) to eV/(atomic volume)." ] }, { "cell_type": "code", "execution_count": 3, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "S11 = 0.1023 S12 = -0.0382 S44 = 0.1187\n" ] } ], "source": [ "stressconv = 1e9*1e-27*Fe.volume/physical_constants['electron volt'][0]\n", "c11, c12, c44 = 243*stressconv, 145*stressconv, 116*stressconv\n", "s11, s12, s44 = (c11+c12)/((c11-c12)*(c11+2*c12)), -c12/((c11-c12)*(c11+2*c12)), 1/c44\n", "print('S11 = {:.4f} S12 = {:.4f} S44 = {:.4f}'.format(s11, s12, s44))\n", "stensor = np.zeros((3,3,3,3))\n", "for a in range(3):\n", " for b in range(3):\n", " for c in range(3):\n", " for d in range(3):\n", " if a==b and b==c and c==d: stensor[a,b,c,d] = s11\n", " elif a==b and c==d: stensor[a,b,c,d] = s12\n", " elif (a==d and b==c) or (a==c and b==d): stensor[a,b,c,d] = s44/4" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Add carbon interstitial sites at octahedral sites in the lattice. This code (1) gets the set of symmetric Wyckoff positions corresponding to the single site $[00\\frac12]$ (first translated into unit cell coordinates), and then adds that new basis to our Fe crystal to generate a *new* crystal structure, that we name \"FeC\"." ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "#Lattice:\n", " a1 = [-0.142765 0.142765 0.142765]\n", " a2 = [ 0.142765 -0.142765 0.142765]\n", " a3 = [ 0.142765 0.142765 -0.142765]\n", "#Basis:\n", " (Fe) 0.0 = [ 0. 0. 0.]\n", " (C) 1.0 = [ 0.5 0. 0.5]\n", " (C) 1.1 = [ 0.5 0.5 0. ]\n", " (C) 1.2 = [ 0. 0.5 0.5]\n" ] } ], "source": [ "uoct = np.dot(Fe.invlatt, np.array([0, 0, 0.5*a0]))\n", "FeC = Fe.addbasis(Fe.Wyckoffpos(uoct), [\"C\"])\n", "print(FeC)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next, we construct a *diffuser* based on our interstitial. We need to create a `sitelist` (which will be the Wyckoff positions) and a `jumpnetwork` for the transitions between the sites. There are tags that correspond to the unique states and transitions in the diffuser. " ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Diffuser for atom 1 (C)\n", "#Lattice:\n", " a1 = [-0.142765 0.142765 0.142765]\n", " a2 = [ 0.142765 -0.142765 0.142765]\n", " a3 = [ 0.142765 0.142765 -0.142765]\n", "#Basis:\n", " (Fe) 0.0 = [ 0. 0. 0.]\n", " (C) 1.0 = [ 0.5 0. 0.5]\n", " (C) 1.1 = [ 0.5 0.5 0. ]\n", " (C) 1.2 = [ 0. 0.5 0.5]\n", "states:\n", "i:+0.500,+0.000,+0.500\n", "transitions:\n", "i:+0.500,+0.000,+0.500^i:+0.000,-0.500,+0.500\n", "\n" ] } ], "source": [ "chem = 1 # 1 is the index corresponding to our C atom in the crystal\n", "sitelist = FeC.sitelist(chem)\n", "jumpnetwork = FeC.jumpnetwork(chem, 0.6*a0) # 0.6*a0 is the cutoff distance for finding jumps\n", "FeCdiffuser = onsager.Interstitial(FeC, chem, sitelist, jumpnetwork)\n", "print(FeCdiffuser)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next, we assemble our data: the energies, prefactors, and dipoles for the C atom in Fe, matched to the *representative* states: these are the first states in the lists, which are also identified by the tags above.\n", "\n", "*A note about units:* If $\\nu_0$ is in THz, and $a_0$ is in nm, then $a_0^2\\nu_0 = 10^{-2}\\text{ cm}^2/\\text{s}$. Thus, we multiply by `Dconv` = $10^{-2}$ so that our diffusivity is output in cm2/s." ] }, { "cell_type": "code", "execution_count": 6, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "dipole: [[ 3.4 0. 0. ]\n", " [ 0. 8.03 0. ]\n", " [ 0. 0. 3.4 ]]\n", "pre: [ 1.]\n", "dipoleT: [[ 6.66 0. 0. ]\n", " [ 0. 6.66 0. ]\n", " [ 0. 0. 4.87]]\n", "ene: [ 0.]\n", "preT: [ 0.1]\n", "eneT: [ 0.816]\n" ] } ], "source": [ "Dconv = 1e-2\n", "vu0 = 10*Dconv\n", "Etrans = 0.816\n", "dipoledict = {'Poctpara': 8.03, 'Poctperp': 3.40, \n", " 'Ptetpara': 4.87, 'Ptetperp': 6.66}\n", "FeCthermodict = {'pre': np.ones(len(sitelist)), 'ene': np.zeros(len(sitelist)),\n", " 'preT': vu0*np.ones(len(jumpnetwork)), \n", " 'eneT': Etrans*np.ones(len(jumpnetwork))}\n", "# now to construct the site and transition dipole tensors; we use a \"direction\"--either\n", "# the site position or the jump direction--to determine the parallel and perpendicular\n", "# directions.\n", "for dipname, Pname, direction in zip(('dipole', 'dipoleT'), ('Poct', 'Ptet'), \n", " (np.dot(FeC.lattice, FeC.basis[chem][sitelist[0][0]]),\n", " jumpnetwork[0][0][1])):\n", " # identify the non-zero index in our direction:\n", " paraindex = [n for n in range(3) if not np.isclose(direction[n], 0)][0]\n", " Ppara, Pperp = dipoledict[Pname + 'para'], dipoledict[Pname + 'perp']\n", " FeCthermodict[dipname] = np.diag([Ppara if i==paraindex else Pperp \n", " for i in range(3)])\n", "for k,v in FeCthermodict.items():\n", " print('{}: {}'.format(k, v))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We look at the diffusivity $D$, the elastodiffusivity $d$, and the activation volume tensor $V$ over a range of temperatures from 300K to 1200K.\n", "\n", "First, we calculate all of the pieces, including the diffusivity prefactor and activation barrier. As we only *have* one barrier, we compute the barrier at $k_\\text{B}T = 1$." ] }, { "cell_type": "code", "execution_count": 7, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "D0: 1.3588e-03 cm^2/s\n", "Eact: 0.816 eV\n" ] } ], "source": [ "Trange = np.linspace(300, 1200, 91)\n", "Tlabels = Trange[0::30]\n", "Dlist, dDlist, Vlist = [], [], []\n", "for T in Trange:\n", " beta = 1./(kB*T)\n", " D, dD = FeCdiffuser.elastodiffusion(FeCthermodict['pre'], \n", " beta*FeCthermodict['ene'], \n", " [beta*FeCthermodict['dipole']],\n", " FeCthermodict['preT'], \n", " beta*FeCthermodict['eneT'], \n", " [beta*FeCthermodict['dipoleT']])\n", " Dlist.append(D[0,0])\n", " dDlist.append([dD[0,0,0,0], dD[0,0,1,1], dD[0,1,0,1]])\n", " Vtensor = (kB*T/(D[0,0]))*np.tensordot(dD, stensor, axes=((2,3),(0,1)))\n", " Vlist.append([np.trace(np.trace(Vtensor))/3, \n", " Vtensor[0,0,0,0], Vtensor[0,0,1,1], Vtensor[0,1,0,1]])\n", "D0 = FeCdiffuser.diffusivity(FeCthermodict['pre'], \n", " np.zeros_like(FeCthermodict['ene']), \n", " FeCthermodict['preT'], \n", " np.zeros_like(FeCthermodict['eneT']))\n", "D, dbeta = FeCdiffuser.diffusivity(FeCthermodict['pre'], \n", " FeCthermodict['ene'], \n", " FeCthermodict['preT'], \n", " FeCthermodict['eneT'],\n", " CalcDeriv=True)\n", "print('D0: {:.4e} cm^2/s\\nEact: {:.3f} eV'.format(D0[0,0], dbeta[0,0]/D[0,0]))" ] }, { "cell_type": "code", "execution_count": 8, "metadata": { "collapsed": false }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAaEAAAEzCAYAAACYBryKAAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAIABJREFUeJzs3Xd4FFUXwOHfBhJ6lxZFBZELUkIMzdBDD70TOtJRVEBA\niqJIUTGgIF2pAgpIFQVBihB6pCsXLCBNPytVCJB8f9wNBgiQsrszm5z3efYhmZ2dPZOEnNyZe89x\nxMTEIIQQQljBx+oAhBBCpF6ShIQQQlhGkpAQQgjLSBISQghhGUlCQgghLCNJSAghhGUkCQkhhLCM\nJCEhhBCWSWt1AEIIUEpFOz88rLUu5fy8gdb6izj7FAG+AQ4COYAgIAbIorW+4umYhXAFGQkJYR/t\ngKrxPaGUKgB8BXwHNAJqAc09F5oQ7iEjISHcQCn1JXAROIL5Y+9V4BPgMPAQ0A0oprU+E+dl57XW\nf8dzrDzABuA0ZnR0FbiqlPrLvWchhPtJEhLCxZRSjwJfaK0nOT8vDwwFemutzzu3cUcCutexsgLr\ngH+AenLZTaQ0koSEcL06wPQ4n9cG9sYmIKcfE3CcTMAXQEmgiNb6outCFMIeJAkJ4WJa65l3bKqN\nGc3E3WdyAg71PvAX8BswGghzSYBC2IhMTBDCjZRSWYDy3JGEEugCUBPoA7RWSrV2ZWxC2IEkISHc\nqwZwGdiVhNcO1Fr/qrVeiZnUMFkpld+l0QlhMUlCQrhXLeBrrXX0A/e82804H/d1fj7LJVEJYROS\nhIRwr7vuByWF1vpP4DmgjlKqd7KjEsImZGKCEC7mXFjaGXgYKASEONf6jNda/5vAw8TcuUFrvVQp\n9RnwjlJqvaviFcJKjpiYu37WhRAeFl+ZngS8piqwESnbI7yYXI4Twj6yK6VyJmRHpVQOIJeb4xHC\n7eRynBD2EAPMx5T5KZWA/VcBwcRz2U4IbyKX44QQQlhGLscJIYSwjFyOcwNnwcq3tNbVlVKlgYnA\nDeAa0FFr/btSqjvQA7gOjNZar1FK5QIWAumBs0AXrfVVpdQ5rXV+57GLAiuBXlrrTfeJwQ+YjZmd\ndR54HnMP4X3ne67XWo9USjmAKUAAcBXoprX+SSm1CeiptT6mlMoMfI4pyvmOS79YLqaUegXT6sAX\nc17fAHOAaEyvnuec+70G1Md8LV7SWu9VSs0GFmmtv1JKpcF8L37XWj/v+TMR4t6UUj7ATEBhfrd0\nwQwq5uBlP+syEnIxpdRAzA9HOuem94DntNYhwHJgsFIqL2bx4TNAXWCsUsoXeA1YoLWuCuwHejqP\nEeM8dnHnMTrcLwE5dQcuaq2fcb7XB8BUoI3WujJQ3pkgmwDptNbBwBBg/B3nkwX4EvMDa/cEVBV4\nxnku1YBHMecz1Pk19VFKNVZKBQJVtNblMfXYptxxnLTAYuAHSUDCphoCMVrrSsAIYAJe+rMuScj1\nfgCaxvm8tdb6kPPjtJjRRjlgm9b6htb6AnAcMxKpBKx17vslpm4YgEMpVQpYBrTUWu9OQBxPOY+B\n1vo4UBaTbE44n1/nPP6t99Ra78J064yVA1gPzNBax60KbVd1gMNKqRWYG/efA09rrbc6n/8SU8Gg\nEqZBHFrrU0Aa5ygUzCh0GbBPaz3Mk8ELkVDOUk49nJ8+BvyKl/6sSxJyMa31cszwOPbz3wCUUsGY\nFe8TgKyYS2SxLgLZgCxxtsduw7l9NmY4nSOBoewHGjjfu4LzWJce8J4AN51DfYCPgSjMoktv8BAm\nibYAegMLuP1n/F7nHPdrPRHICDzi7mCFSA6tdbRSag7mZ/YzwBHnaa/5WZck5AHO6sdTgFBn+ZUL\nmEQUKyvwt3N7Fue2LJhGZmAuxzUGOgBzlVIPJeBtZwEXlVIbMdeDD2D608TKEs97AvjEqXM2yPna\nLkqpygl4T6v9CaxzjjCPYUad2eI8H/ec7/z6x36t38eU2imllGrn/pCFSDqtdWegCPAhkCHOU17z\nsy5JyM2UUu0xI6BqWuuTzs27gUpKKT+lVDagKKbtcwTmlz5APSB2aH1Ja31aa70Pc29nYQLeuizm\nkl8IsAI4BkQppQo6JyPUcR5/OxDqjLUCcCjOMY44G6l1BD5OYPKz0jbMPTaUUv6YpPu1814R/Pc1\n3Q7UVko5nF1QHVrr2FbZR5xJuB2mPE5Rj56BEAmglGrvnIQD5o+tm8Beb/xZlyTkRs7LWu8DmYHl\nSqmNSqkRzkt0EzG/NDdgbiZGYRqXtVFKbQUqYBIOxFmQqLUeD0QrpV59wNsfB3orpbYDI4F+QC9M\nAtsJfKu13oOZ6HBNKRUBhAMvxfOeuzCdQhOS/CyjtV4D7FNK7cbMIOwNDADecJ6fL7BUa/0t5j/o\nDmAJpl8P3H7OP2NGgouVUuk9dxZCJMgyIFAptQVz/+cFzB+7XvezLotVhRBCWCZFrRNSSgVgRhg/\nAXO01lssDkkIIcR9pLTLceWAc5jZaUcsjkUIIcQD2H4kdEf1gXhX98fZfSumDXJeYCAw2NPxCiGE\nSDhbj4TiqT4Q7+p+pdRIpdRCoDSQBjMFMY3nIxZCCJEYtk5C3F194M7V/WWcH7+mtW4LnAQmAW9j\n7g0JIYSwMVtfjtNaL1dKPRZn052VBm4opW4trtRa78BMRbyvyMhImRIohBCJFBQU5HjwXolj6yQU\nj/ut7k+UoKCgB+/kKW++yaRZs+g7ahS0bQuOpH+fz549i7+/f6JeM2nSJPr27Zvk9/SkpJzfvdjx\nvF15fnYk5+d5rvo5j4yMdEE0d7P75bg7RXDv1f3e69VXoVUrCA+HatXg8GGrIxJCCI/wtiR05+r+\nfhbH4zqPPAJ79kDr1lC9OvTvDxcuWB2VEEK4le0vxznrrQU7P47BlGJJmdKkgT59oEULGDwYihWD\nd9+FNm2SdYlOCCHsyttGQqlDnjwwezYsWQJvvw01a8L331sdlRBCuJwkIZsoV67c3RuDg2HvXmjc\nGKpUgVdegcuXPfPeqUBqPW+Rutj951ySkE2UL18+/ifSpoUXXoCDB+HUKXjqKVixAlxYePae753C\npdbzFqmL3X/OJQl5i/z5YcECmDMHhgyBBg3gp58e+DIhhLAzSULepnp1OHAAKleGcuXgzTfh2jWr\noxJCiCSRJOSN/PzM/aHISPMoWRI2bLA6KiGESDRJQt7sscfM/aHwcOjeHcLC8Pn1V6ujEkKIBJMk\nlBI0bAhHjkChQuSuVQsmToSbN62OSgghHkiSUEqRMSOMHs2fy5bB8uXmftGePVZHJYQQ92X7igmJ\noZQqBrwI3AQma62/szgkj7vx5JOwcSN8/DE0agTNmsHo0ZA9u9WhCSHEXVLaSKg3cAZzXiesDcVC\nDgd06ADffQc3bpi1RQsXunRtkRBCuILtk5BSqrxSapPzY4dSaqpSartSaqNSqtAduz+GaWq3FOjk\n6VhtJ0cOmD4dli2Dd96BWrXg2DGroxJCiFtsnYQS2d57EfA7cAX4C5CKn7EqVDDlf0JDTSmg11+H\nq1etjkoIIWx/Tyi2vfd85+e3tfdWSt1q7w2glArCJK0oYND9Dnz27Fk3hWytixcv3vvc2rTBp3Jl\nso0Yge9TT/HP2LFEVa7s2QCT6b7nlwLI+Xm3lH5+7mDrJJSE9t6RJPAynN26H7rKAzs7+vvDF1/A\n6tU89PzzpvJCeDjkzeu5IJPBjp0rXUnOz7ul5PM7d+6cW45r68tx8XBZe+9Ur2FDM3HB399UXJgx\nA6LlSymE8CxvS0Ipsr33xIkWVd3JlMlMWNiwwRRGrVQJDqWIL6kQwkt4WxJKke29AwKgUycYOhSu\nX7cggFKlYNs26NgRQkJMV1c39C0SQog72fqeEKSO9t67j46j3/icfD2nK1WrwqJFpiycR/n4QK9e\n0KQJ9O8PJUrA5MlmRp0QQriJt42EUiSfnD4MjuzBlaer0bDhDcqWNUt7LJEvn1nYOn26aabXqhXI\nbB8hhJtIErKBAS0HcLDPQTgZw6JFTzN+/DEGDoQ+feDffy0KqnZtc3+oSBFzvXDKFCmKKoRwOUlC\nNlH88eJsnr+Zl19+mf79K9GuXTh//BFN+fLw/fcWBZUhA4waBZs3m2uEFSuahnpCCOEikoRsxOFw\n0LFjR/bu3cv27V/y46/lCGm5nSpVYPZsC0u/FS8OW7ZAt26m9M+gQTJxQQjhEpKEbOjRRx/lq6++\nolzjsszYXpOePRcRHh5D+/Zw4YJFQfn4mCR06JC5R1SihFn0KoQQySBJyKZ8fHyYOmAq+9/fz4YN\n75M3byNiYi4RFGQ6elsmb17TJmL6dOjbF9q0AenmKoRIIklCNlekSBG2bdtGSEgFNmwoRK3am6lb\nL5r33rO4M0PsxIWCBU3FhenTpeKCECLRJAl5gbRp0zJs2DDWrVvHylPt8GlTmI8+OUvjxvDnnxYG\nljEjjB1rmujNnWvq0B05YmFAQghvk6KSkFLqRaXUbKXUNqVUL6vjcbXAwEC+X/Q9j2V5iO+qFOBy\nhrUEBsLWrRYHVrKkqbjQvj1UqwbDhlk4t1wI4U1SVBLSWr8P9AAOa62nWR2PO2TNlJXdo3czrdo0\nfo7sg1ITaNEimlGjLF7G4+MDvXubKdzHjplSQBs3WhiQEMIb2D4JJbKzKkAYYFW9AY/pHtqdg/sP\nUqjQUdKmrcCSJX9Tuzb89pvF31J/f1iyBMaPh86dzeOPP6yNSQhhW7ZOQonsrLpQKZUDqKy1/sqS\ngD0sc+bMTJ8+nQ8/fIM//gjgdIaR1K6fkXXrrI4M0yriyBHInt1M554/3+KZFEIIO7J1EuK/zqqx\nbuusCtzqrKq1bqu1/htI4/EoLVavXj0OHtxHTOGPSZOxMx06RDFkiEUVuePKkgXeew8+/9yMjGrX\nhh9/tDgoIYSd2LqKdmI7qzpf82xCjp0SW/BuHrSZxY8v5s03S7Fo0XK++qogU6ee55FHLK755u8P\nK1eS6cMPyVy2LJd79+ZSjx7g65voQ6X09slyft4tpZ+fO9g6CcXDZZ1VU2oL3latWtG6dWu6du3O\ngQO1qNugG7NnZqJxY6sjA0aOhC5dyNq7N1kbNYKZM6FcuUQdIiW3TwY5P2+Xks9P2nsbKbKzqqvl\nz5+fNWtWM+L1DJxvlY+wSaN4/oXrXLtmdWSYxa1ffgmvvAKNG5t2ERcvWh2VEMIi3paEUmRnVXdw\nOBz06N6DTT3WkLb4eGb+VZugoGv88IPVkQEOB4SFweHDcOmSKZC6erXVUQkhLGD7JKS1PumcDYfW\nOkZr3VtrXdH5OGZ1fHZXpVQV/nj3V/oHPMOJE68TGPgvixbZZJZarlwwa5apttC/v2mg56YhvxDC\nnmyfhETy+fn6MXbgGHbtas8jj3SjW7eztG9/2T5FDapXh4MH4cknzSLXGTOkDp0QqYQkoVSkePHi\nHDw4h75957Lky08oUGeudQ3z7pQhA4webaoszJplyv8cPWp1VEIIN5MklMr4+vry1ltDmbs4E9f9\nPiEw8AKTJ1+yOqz/lCwJERHQsiVUqgRvvglRUVZHJYRwE0lCqVSbGm34dfUyWracxksvnaNWrdP2\naZaaJo3pVbRvH+zeDU8/DTt2WB2VEMINJAmlYhkyZGD+/EGsXHmWHTt28nCBU2ze/j+rw/pPgQKw\nahW89ho0bw7PPYdDpnMLkaJIEhKEhlbl7Nna+Fd9m+qLy9F2yEL7lHlzOMysuSNHICqKPNWrw8qV\nVkclhHARSUICgKxZs/Ld8g/o8Wg3Fl1/icerT+P3369aHdZ/cuSAmTP5e+JEGDgQWrSQ6dxCpACS\nhMRtpvcfzv4e23D8UpBHHjnHp5/aa4ZaVHCwmc6tlEznFiIFkCQk7hJQpAg//1ibZ589Q1jYQzRq\ntIbr129YHdZ/0qc307m//ho++ghCQkwjPSGE10lRSUgpVUYpNVkpNV8pFWB1PN7M4XAwdWolNm68\nwaZNxche42VWbbVZp9RSpWD7dmjaFIKDYcwYG/SvEEIkRopKQkAQUAx4GDhlcSwpQrVq+fjtt4Lk\nz5+DJvUzM2jQQqLtdPkrTRp48UWIjIRt2yAoyEzrFkJ4BdsnoUS29/4WqA28BTTwdKwpVcaMDn74\ndAThbxRgwoR6PPXU+5w8+YvVYd3uscdgzRpTnbtRI+jXzxRHFULYmq2TUCLbey8C3gRigD+AnJ6P\nOGXr1y8/Bw9m5cKFthQpEsmYyZPsNSpyOKBtW1Od+6+/TPUFW/Q6F0Lci62TEIlr7x0GTAHmA32B\njz0baupQrFgafvwxL40aVWb4nsXkefFJjpw4YnVYt3voIVOZe/p06NULOnSAP/6wOiohRDwcMbZZ\nlRg/Z3vvRVrrYKXUTGCp1nqd87kTQKHEdleNjIyMyZ8/v8tjtYOLFy+SJUuWB+/oAosW32Tw2veI\nKbiAaU+Pon79ULe/Z2LPz3HlClnGjSPD8uVceO01/m3a1IyYbMqT3z8ryPl5r3PnzhEUFOTy/zzS\n3juF8WR74QEvQcPQcOo1Gkz/j3fz9dcjmDZtHNmzZ3fbeybp/KZPh27dyNGtGznWrIFp08w9JBtK\nye2hQc7Pm0l7b0Pae9tMkSJweF8eWrSow9q1oyhatDXr16+3Oqy7lS0Le/dC5cpmBt3EiXDzptVR\nCZHqeVsSkvbeNpQhA8ye7cvkyXm5enU1Lbos4+lBQVy2TVluJ19fGDrUtIpYuhQqVjSTGIQQlrH9\n5Tit9UngVntvoLe1EYl7adcOgoL8aNx6FGd2ziAgoCIffzyNChUqWB3a7ZSCzZth5kzT1bV3bxg2\nDNKle+BLhRCu5W0jIWFzRYvCvh25qF94CFevfkODBq8wbNgwouzWmM7HB3r2hP37TS26wEBTfUEI\n4VGShITLZcxoOnSPHp2VmJiNfPFFbsqWL8uKbSusDu1uDz8My5fDyJGmMnffviA9i4TwGElCwm06\ndYJvvvHh6tUX8cn/Os0+b0690fWIum6zUZHDYRLQ4cNw+TIUL26qLwgh3E6SkHCr4sVhzx4HxXM2\npeD674n4bR8lwkrw448/Wh3a3XLmNEO42bPhhRfMTa7ff7c6KiFSNElCwu0yZ4b582FwzyL4LTpD\n+YzjqVChAtOmTcOWi6Vr1DD3ifLlgxIl4OOPsU+rWSFSFklCwiMcDujRA9Z/lYYd2xtQs+YxZs6c\nR2hoKGfOnLE6vLtlygTh4fD55zBuHNSvDydPWh2VECmOJCHhUYGBputCVFQOYBtK1aVI9yK8+tGr\n9hwVxS5yrVjRLHKdNEk6uQrhQpKEhMdly2bWinbs6MOiRS/So/oCPp30Ka1ateIPOxYa9fU164i2\nbYPFi6FSJfjuO6ujEiJFSFFJSCnV2tlVdYJSKpPV8Yh7czhML7qVK2HZB01oEHKEhx9+nICAAD7/\n/HOrw4tf0aKwZQu0bw9Vq8Kbb4Ld1j8J4WVSVBICGgKdMO0cOlkci0iAChXM5bnvvvMlMnIcEycu\npW/fvjTu05jTv5+2Ory7+fhAnz7w7bewcyeUKSOdXIVIBtsnoUR2Vv0A+BBoBOTydKwiaR56CL74\nAurUgb59n2HixCP8nPFnCr5dkPdXvm91ePErUMBMWojt5DpggFljJIRIFFsnoUR2Vl0I5AO6AVsB\nm/WfFvfj4wPDh8O8edCjR0ZaZT/IkNLDGLh2IC/2f5F///3X6hDvFtvJ9dAh+O03KFUKvv7a6qiE\n8Cq2TkIkrrNqW+A8MAdoByz0aKTCJWrWNJPR1q2DvQtf53C/c5w7fY6goCD27t1rdXjxy53brCWa\nNAm6dIGuXeHvv62OSgivYOsq2lrr5c7OqrGyYhJNrBtKqVuN7bTWm4BNCTn22bNnXReojVy8eNHr\nz83hML/T33orKzVqZGPq1A/45Zel1K1bl7CwMF5++WV8fX2tDvNupUvjWL+erGPHkv6ppzg/ahRX\n69VL1CFSwvfvfuT8xJ1snYTiIZ1VHyAldXacOhVq14Znn83D66/3Yf/+xrTq2pJS75Qiok8EJYqX\nsDrE+M2ZA1u3krNbN1i71oyQ8uVL0EtT0vcvPnJ+3ks6qxrSWTWVadrUdFiYMQMGD36YOZOW8HzR\n56lWtRoTJkwg2q4LRytXhgMH4Mknzb2iOXOk9I8Q8fC2JCSdVVOhwoVhxw7Tc65hwzx0qDWanTt3\nsnTpUkJCQjhx4oTVIcYvfXoYM8bc4Jo40Uz/+/lnq6MSwlZsn4S01ieds+HQWsdorXtrrSs6H8es\njk94RoYMpsB1jx6XqVwZ9u0rzDfffENoaChlypeh16Re9h0VBQbCrl0QEmLKAL3/Pty8aXVUQtiC\n7ZOQEHGFhV1h3ToYPBgGDkxDv36DWLRqEXN/mkv+Afk5+NNBq0OMn6+vWVMUEQGffSalf4RwkiQk\nvM7TT5tp3MeOQfXqULxALX4f8ztFshahUrdKLF261OoQ700p2LwZOnaU0j9CIElIeKmcOWHVKqhb\n11TO2bsrM1vf2Mr6sesZNmwY7dq142+7rtXx8YHevU3pnx07/qvULUQqJElIeK3YKgtz5kCbNqbt\nT7ly5dm3bx+5cuWiVKlSrF271uow761AAdNGfNAg069o0CAcdqwMIYQbJXidkFIqY1LeQGt9JSmv\nEyKhatc2NURbtDA1RWfPzsjEiRNp3LgxzSY1I+irIFa/uZpMmWxYWN3hMG3Ea9WCF14gd82aJqtW\nrWp1ZEJ4RGJGQheT8LjgymCFuJdHH4WtW00FnbJl4cgRqFGjBgenHyTv33kJCAhg+/btVod5b3ny\nwCefcH7ECJOUeveGC/LfR6R8iUlCDqAFEJLAR0vna4TwiHTpYNo0GDIEqlWDTz+Fx/I+xqLZixg3\nbhzNmzfnlVde4dq1a1aHek/XateGw4fNFO4SJczlOiFSsMSU7YkANmmt/0nIzkqpHICN//QUKVXn\nzlC6NDRvbpbnvP02NG3alIoVK9KjRw9UQ8U7w9+hVZVWVocav+zZTYmIjRuhe3dYtAjee8/0vBAi\nhUnwSEhrXTluAlJKPaGUyur8uIZSarJSqmOc/f/WWld2bbhCJEzp0mbCmdamMvevv0KePHlYvnw5\n9RrXo80Xbagzqg5Xo65aHeq9hYSYNhF585pR0SefSOkfkeIkaXacUqoDoIGySqniwGogAHhPKTXU\nhfElJJbqSqmZzo+fUUrNUUrNjk2QIvXKkQNWrzZricqUMTXoHA4HU5+bSkTHCA6dPUSVkCocP37c\n6lDvLWNGCA83fdBHjYImTUCqNIsUJKlTtIcAz2mtvwY6A8e01pWANkBPF8X2QEqpJ4Cn+a/pXQ/n\n4yNnLCKV8/GB11+H6dPN7+/Jk81g4pmnnuH0B6fp0LoDzzzzDJMnT7Zv2R+A8uVNH/TSpc3jo49k\nVCRShKQmoULAF86P62NGQmBGR7mTE1Bi2nlrrX/UWofH2ZRGax0F/IrpsioEYJbhbN9uklGXLvDv\nv+Dj40Pfvn2JiIhg3rx51KlThx9P/Gh1qPeWLh288Ybp3jptmrnO+NNPVkclRLIkNQmdAsoopcoA\nRYHYKTyNgCT/L05sO2+lVPY7DnFZKeUH5MckIiFuia3Gfe0aVKwIJ0+a7UopIiIieKbqM6gJig/n\nfkiMnUcZJUuaE6lXD8qVM5MWpCCq8FJJTULjgCWY2W9faa13KqVex7RXeCMZ8SSqnXc8M/VmAtMx\nl+Q+TkYcIoXKlAkWLoQOHcwVrg0bzPa0adMycvhINrTdwPvvvk+LFi34/fffrQ32ftKmhZdfNslo\n2TIpiCq8liOpf/EppQKBx4C1WuurziZzl7XWyWo052znvUhrHeyccLBUa73O+dwJoFBSu6nGioyM\njMmfP39yDmFbFy9eJEuWLA/e0Uu58vwiIvx4/vkc9OhxiV69LuNwrmq7evUq7777LsuWLeOtt96i\ndu3aLnm/hEjS+UVHk3H+fLKMG8flbt249Nxzpmq3DcnPp/c6d+4cQUFBLl/7meAkpJQ6DnwFrAc2\naq3dspz7jiQUDuzQWi91PveL1vrR5L5HZGRkTFBQUHIPY0spub0wuP78Tp2CZs2gUCHTryhuZZ+t\nW7fS8dmOZG2Ylc+HfE6B3AVc9r73kqzz++UX6NnTzEefNcv0MbIZ+fn0XpGRkW5JQom5HFcf+A7o\nApxQSkUopUYopYKVUu4qhCrtvIVbFShgyv1kzAjPPAM/xrmjWblyZXbs2sH1tNcp+HZBxi8fb12g\nCfHoo/DFF9CvnykvPmQIXLXxOighSNxi1WNa68la68aYGXCDna8PB/6nlFqulOqtlCrswviknbdw\nu/TpzcChZ08IDjbduGPly5mP7975jteCXmPU+6N48cUXuXLFxjV5HQ7Tq+jgQTh+3EznjoiwOioh\n7i0mJibZjyJFimQrUqRIsyJFikwtUqTIcVcc052PvXv3xqRUZ86csToEt3L3+X3zTUxM/vwxMW+/\nHRMTHX37c3/++WdMWFhYjFIqZteuXW55f5ef39Kl5oReeCEm5uJF1x47CeTn03s5f2+6/PexSy6j\naa3Pa62Xaa17a62fdMUxhbBC5cqm3tySJRAWBpcv//dczpw5WbhwISNHjqRhw4YMeW0Il/+9fO+D\n2UHz5qYg6j//mKnd69dbHZEQt3lgElJKpVFK9VdKva+UahJne1OlVFf3hieE5xUoAN98Yy7TBQfD\nzz/f/nyrVq3Yv38/q06tovCzhTly5Ig1gSZUzpwwdy5MnQrdukHXriYpCWEDCRkJTQP8gaNAQ6XU\nx0qpdFrr5cBYt0YnhEUyZIDZs83v62eegU2bbn8+f/78HProEK9Wf5Vq1aoRHh7OTbsvGK1b14yK\n0qc3BVFXrbI6IiESlIT2a61f1lpP1Vp3BV4HRiqlUuZkeCGcHA544QWzuDUsDCZNur1cm4+PD316\n9GH37t2sWrWK6tWr85Pdy+hkyWIK6C1YAAMGmBOz86JckeIlJAldU0pVUEqNUUpl1lr/AAwHWgHp\n3RueENYLCTGFCT780FzNurMnXsGCBdm0aRONGzcmoEMAHd7rYO9iqGDahx84AA8/bO4VSZsIYZEH\nJiGt9YcMOf44AAAgAElEQVSAH6Y46WXntuta64+A5u4NTwh7KFjQzHQ+f94kpV/vqEzo4+PDgAED\n+Hjcx3x28jPy9c/H4Z8OWxNsQmXMCO++ay7LSZsIYZEEzY7TWn+jtZ6rtb7zT6WDSqlH73y4IU4h\nLJc5MyxeDLVrm7qhkZF379M4uDF/jP2D0tlKE1I5hE8//dTzgSZW7MnEtomYNUtGRcJjktrUrpFS\n6n/AWeDnOI8Tzn+FSJF8fGDECFO4um5dcxXrThnTZ+SrN77i8+WfM2LECMLCwvjrr788H2xixLaJ\nWL/e3DOqUwdOnLA6KpEKJHWd0HuY6tblgVJxHiWd/wqRojVrZtr6DBkCw4dDfLeAypUrx759+8ib\nNy+lSpVi/sr5ng80sQICzEKpGjWgbFn44IP4T04IF0lSFW2l1EUgSGt9zPUhJZ5SqjrQVmvdPb7P\n75SQAqaHDh0iKirK5bGKpPPz8yNXrly2KhD5v/+Z9aAPPQTz55tLdvHZsHED9ZfXp/WN1kwZN4XM\n99jRVgUwjx41c9R9fEwn1yJFkn1IW52fG6Tk83NXAdO0SXzdEkwDu3ddGEuS3NniO56W30kSFRVF\nSq207a0i47sJY7E8eUxPoj59TKO8Vavgscfu3q9mSE1+C/qN/v36ExAQwJw5c6hcubLnA06MokXN\nqt3Jk82q3cGDTXHUtEn9tSHE3ZJ6Oe41YIhS6lul1FKl1OK4j+QGlZwW3/G0/BYpyJw5c1iwYIHV\nYdwmXTozfbtzZ7OwdceO+PfLni07s2bNYsKECbRu3ZqBAwdy1e5VrtOkMYuldu+GtWvNCR6SYvbC\ndZKahGYBN4FjwEXM1O24jyRLRovvO4eJLh82CuvlypXLlh1PHQ4zSJg5Exo1Mpfm7qVRo0YcOHCA\nH078QJ6+eVi4aaHnAk2qQoXMkK9HDzNHfeRIkMvVwgWSOq6uBFTSWn/rymCcYlt8x/43vq3Ft1Lq\nVovvO153582t+97sOivrIbxSVFQUV65cse33LzAQFi9OS+fOOdm9+18GD76Izz3+1Js0YRJDPxlK\n+7Xtmbl+JnN7zeXff/+17bkBUL8+PoGBZB88mDSffMI/48dzvVTC5yJdvHjR3ueXTCn9/NwhqUno\nGJDBlYHE0lovd3ZXjZUVOB/n8xtKKZ87W3xrrTve7/M7Pejm4blz5xIWsPAoPz8/MmbMaOubv/7+\nsHcvNG2ahX79sjB3rlkXGp95L8+jr+5Lz1d70rJlS8aNG4dSyrMBJ5a/vxkVLVhA7k6d4Nlnzbz1\n9A8uoJKSb9xDyj4/d/1OTOrluFHAXKXUi841Q6FxH64MELgAxK1Td1cCEsJucuc2U7jTp4dq1eB+\n/3/LqrJEfhpJp06daNKkCRMnTrR/2R+HA9q3N83zfvjBLHLdvt3qqIQXSmoSWgwUAiYAK4DP4zxW\nuya0W6TFt/BK6dLBvHnQsCFUqGB+X9+Lw+GgT58+rFq1ikWLFlG9bnV2fHePGQ52kjevab40ejS0\naAEvvnh7EyYhHiBJSUhr7XOfRxoXxygtvoXXcjjg1Vfh7bfN+s8vv7z//oUKFWLr1q0UrFqQkKEh\nzJkzh6Ss5fO45s3NrLm//oJSpWDjRqsjEl7ClhP+tdYngWDnxzFAb2sj8l5RUVHMnDmThQsXUrBg\nQSpVqkR0dDR//fUXv/76K/369eOJJ56wOswUr00bePRR87t6xAjo1eve+6ZNm5Y5w+bQ70A/OnTo\nwIoVK5gxYwZ58uTxXMBJkSuXmRa4Zg106gT168M770DWrFZHJmwsqbXj0iilXldK9YizbYdS6lWl\nlEtahgvX8PPzo2fPnly9epXnn3+eXr160adPH4YPH85zzz1H165duXTpktVhpgrBwbBtG0yYAAMH\nPrgaTkBAAHv27KFYsWKUKlWKZcuWeSbQ5Kpf3zTPi442zfMeNPwTqVpSE8Y4oAtwMs62ac5tI5Mb\nlHCtffv2ER0dfVcFiGLFipE9e3Y2bNhgUWSpzxNPmMWsu3ZBq1bw77/33z9dunSMHTuWZcuW0Xtq\nb554+QlO/nby/i+yg2zZYMYMU5G7Tx8zMrJ7EVdhiaQmoTAgTGu9LnaD1nou0BmTiFKsEiVK4HA4\nPPooUaJEsmKOiIigTJky+Pr63vXc1atX+eOPP5J1fJE4OXOaYtXp0pn7RAn58gcHB3Nw8UEypMnA\nE+8+wUerPnJ/oK5Qs6a5V5QtG5QsSXoZFYk7JDUJZcRUSrjTn0C2pIdjf4cPHyYmJsajj8OHk9cc\nLSIigooVK961/fLly5w+fZqCBQsm6/gi8dKlM7dPqlY1NecS0hU8b468HH77MGPKjWHECyPo27cv\nV65ccX+wyZU5M0ycCJ98QtbRo6F1a2kpLm5JahLaBLyjlMoZu0EplQMYDWx2QVzCRc6fP8+RI0eo\nVKnSXc+tWbOG7Nmz27+QZgrl4wNjx5pZzZUqwZ49CXvdoJaDOLTvEH///TeBgYHs3LnTvYG6SuXK\n/G/9elPhVVqKC6ekJqG+mHVCZ5RS3yuljmIa3BVyPidsYvv27eTJk4fChQvftj0qKoq5c+cyYsQI\n/Pz8LIpOgLllMnUqhIbChg0JK/6eI0cOPv74Y0aPHk2TJk3oNLwTl/71ggkmGTKYGXOrVsGbb0LT\npvdfyStSvKSuEzqFaWDXAvgImIKp91Zaay2dVW0kIiKC4ODg27b9+eef9O/fn44dO1KrVi2LIhNx\nNW4Mq1fDyy9nZ+7chL+uRYsW7Nu3j43nN1KuZjkOeUuF63Ll4NtvzYgoIADmzJFRUSqV4HVCSqmD\nQFWt9d8AWusoYI3zEd/+OYHNWmvptGqBo0ePsnbtWtatW0dgYCDTpk0jOjqaq1evcuHCBfr27Wv/\nGmWpTIUKsHTpn3TsmIfffjPTuB0JqAWfP39+fpn4C7NnzyYkJISBAwcyYMAA0qRx9bpxF0uXzoyG\nmjUz9ec+/dTMqCtQwOrIhAclZrFqCaCuUur8A/c0sgPFEx9S4sTtoqqUCgE6YYqrvqm19pI/C12v\naNGiFC1alJdeesnqUEQiFC58g4gIqFvXXKUKD+eeVbjjcjgcPPvss4SEhNC5c2dWrVrFrNmzKPJk\n8ruhul1goOlX9M478PTTpgRQ9+4Jy8DC6yX2ctwCbq8Td7/Hx64LM37xdFHNoLXuBIwBarv7/YVw\nh4cfNg1NIyOhQwe4fj3hr3388cfZuHEjzZo3o/j44rQd39b+xVABfH1h2DDYvNl0CKxZM2FTBoXX\nS/BISGvtkUoISqnywFta6+pKKQfmflMAcBXoprW+9ZOptf4RCFdKzXN+vkYplREzOWKwJ+IVwh1y\n5IB168xs5qZNTY3QDAlsnuLj40P/fv0psqsIrT9pzcMvPszeV/by8MMPuzdoVyhe3FTjnjDB3Dca\nMQKeey5hw0HhlWz1nU1uV1WlVC7gfeA1rbWswBReLUMG+OwzyJ7dXJ67cCFxr29QvgG/j/2dxg81\nJjAwkEWLFnlHMdS0ac0NsYgIc5+oalU4dszqqISb2CoJ8V9X1Vi3dVUFbnVV1Vq31Vr/49wv9n/W\neMAfGKuUauaZkIVwH19f0w6iRAnTVTuxxS0yps/ItBHT+OKLL3jzzTdp3bo1f/75p3uCdTWlYMsW\naNnSFN4bNw5u3rQ6KuFitqqindyuqs77QQkiLXi9k93beyfXvdpDDx0Kb7+dheDg9Cxa9Cf58yfu\nPo+/vz+rV6/mnXfe4amnn6L9kPYMaDTAVWEnWJLaX7doQZoyZcj+8ss4Fi3in/Bwbth0Zqe09068\nJCUhpVQWzIglL6CBdVrra64MzMltXVWlvbd38ob23slxv/bQEyeaCWQtW+Zj40Z4/PHEH3/69OlU\nWFOB/jP6c37vecaPH09WD7ZaSHL7a39/U4J8xgzytGwJ/frBoEFmqGgj0t478RJ9OU4pVQw4gmku\n1xVYAvyslGp63xcmjXRVFSKOQYPgpZdMy/CkTh7rUr8Lv3z8Cw6Hg4CAALZs2eLSGN3Gx8c0Yvr2\nW9i6FcqXh/37rY5KJFNS7gmNBhpprQO11kWBXMAIYFTc/kIuIl1VhbjDCy/A4MEmEf3wQ9KOkSVL\nFmbOnMmkSZNo27Yt/Qf0559L/zz4hXbw6KOmR9ELL0Dt2vDaaxAVZXVUIomSkoROa61v/fmhtb6k\ntZ6JWa9TVylVMjkBaa1POmfDobWO0Vr31lpXdD5kiowQQO/epm149eqgddKP06BBAw4cOMCev/eQ\n97W8fPy125f3uYbDAZ07m5HQ/v0QFAR791odlUiCpCSh/8W30XlPqBfSilsIj+jeHUaONLPmvvsu\n6cd56KGH2PLhFrqr7nT8qiMd3+jI9cSskLWSvz+sXAlDhpiOrq+8AlevWh2VSISkJKF7TgzQWv8P\n8JIxfeqzYsUKqlWrJp1UU5AuXeCtt6BWLTh+POnH8fHx4YOeH7Cn6x7O7TxHcHAw33//vesCdSeH\nA9q2hYMHzfXJwEDTvlZ4haQkoTpKqepKqXtNS/GCLlupU2hoKJcuXaJcuXJWhyJcqEMHeOMNk4h+\n+SV5xwoqEsRXX3xF165dqVKlChMmTPCOsj8AefPC0qVmeNisGQwYAN7Q9C+VS0oSUsAK4G+l1BdK\nqZeUUk/Fed5LfmJTnx07dlC4cGGPTskVntGtm5k1V6NG8tvzOBwOevXqxY4dO/jss88o2aYkEYcj\nXBOoJ7RsaVqKnztn2kR8843VEYn7SEoSmoyZERcKRAJhwEGl1C9KqY+Au/tIC0vExMQwZcoUFi1a\ndOsRt4vq6NGj73pNfNuEd3jpJXOvvlYtcEVRhMKFC7NlyxaeKP0EDRs1ZNasWd5R9gfgoYdg4UJ4\n910IC4O+feGSFzT9S4WSkoTe1lrf0Fp/o7V+VWtdHsgNDMCMgoLv/3LvVqKEuQTtyUeJEkmLdcSI\nEcTExBAWFkZYWBh79+6lSpUqXL9+nXnz5t22PiS+bcL7DB0KDRpAnTpwPqFNV+4jTZo0rBq6ii0r\ntzBx4kQaNWrEr7/+mvwDe0rjxmZUdOGCaaD39ddWRyTukOgk5Gxmd+e2v7XWS7TW3YEPXBKZTR0+\nbBpAevJx+HDi4zx69CgrV66kUydTyejUqVP4+vpSsmRJfH196dixI/ny5bu1f3zbhPdxOGDsWNMg\nr2FD100UK1myJLt37yYgIIDSpUsz+5PZrjmwJ+TMCXPnwuTJZqjYs2fiq8EKt3FHAdOlbjimSKSd\nO3cSEBBA5syZAdi+fTsVK1bk2rVr3jP9ViSJw2FK/OTPbxqWuuoKmp+fH6NGjWLFihU8v+l5Cg4o\nyM/nfnbNwT0hNNT8RRcdbUZF69ZZHZHADUlIa33A1ccUiZctWzZy584NmEttq1atIjg4mNWrV3NT\nKhGneD4+MGcO/PyzacnjShUqVODkOyfJ6psVNVaxdt1a176BO2XLBjNnmsZ5PXuaLP2PrCqxkt1a\nOSSac7r4TOfHTyulFiilZiulclsdm5Xq169PhgwZWLNmDStWrCA0NJTjx4/j5+dH+vTprQ5PeECG\nDGYd54IF5mqUKz2U7SEOvHWAT+p9Qo/uPejTpw+XL1927Zu4U61a5l5Rhgzmpuvnn1sdUaplq1YO\niRVPe+90mIoNdYBngFUWhWa52Esn9xPfTCevmf0kEiRPHvP7tVo1U3KtenXXHr9ZvWaEHAzhxRdf\nJCAggHnz5hEc7CVzk7JkMfeJWraErl1NA7333zf3kITH2G4kpJQqr5Ta5PzYoZSaqpTarpTaqJQq\nFHdfrfWPWuvwOJ/vAIpjZurt82jgXmbBggX88ssvTJ8+nd9///2e24T3K1YMFi2CNm3g6FHXHz97\n9uzMnTuXcePG0bRtU6oMr8JVbyqdU62aqbaQK5cZFS1fbnVEqYqtkpAL2nuXAfZi1jC94MnYvU27\ndu3YsmULPXv2vHXvKL5tImUICTHlferXT3x31oRq2rQpm7ds5uLZi5QrV44DB7zo9nCmTPDee7Bk\niSlR3qYNyB9iHmGrJETy23tnBWYB72P6HAkhnLp0gSZNTHUFd111LfZYMb796FsGDBhAzZo1GTt2\nLDdu3HDPm7lDxYpw4AAUKAClSsHixe77YgkAHHa7B+Bs771Iax3snHCwVGu9zvncCaBQcrurRkZG\nxuTPn/+++5w7d46goKDkvI1wscjISLZv387Fixfp3Lmz1eG4xcWLF8mSJcuDd0yia9egYcPcdOp0\nmXbt3FtX7cyZM/Tr148r/15h6FtDCS4e7PbzcyXfb78le//+3HjySc6PGUN0Aq4OeNP5JZbzd6LD\n1ce1+8QEae8tbpOa23u7ypIlUKVKdho1yo5S7nsff39/vvnmG4a+P5RWX7eixYEWjG8z3nu+d/7+\nULMmvm+8QYbatWH8eFOt23Hv38PS3jvx7HY57k7S3lsIFytWzBSabtfO/Q1JfXx8eKvfW6xpvoYv\nznxByMshnD592r1v6krp05sSFGvWmJtqjRvD2bNWR5Wi2D0JSXtvIdygVy/zh/5rr3nm/eqVrccf\nb/9Bw0cb8vTTTzN//nzvWg5Qpozp3Fq6tHnMnSv3ilzEdpfjtNYncRZB1VrHIJ1ahXA5hwM++sj8\nPq1Tx/Xrh+KT3i89A14aQPuw9nTo0IEVK1Ywbdo075mJmS7df72KOnc264pmzIBHHrE6Mq9m95GQ\nEMJNcuc2iahTJ/jrL8+9b2BgIHv37uWJJ56gcNvCjF041nNv7gqlS8OePaZKbGCgKQEko6IkkyQk\nRCpWt675w/655zz7vunTp+edd95hzAtjmDp2Kl26dOG8K3pPeIqvr7mWuXEjTJtmhpMnT1odlVeS\nJCREKjd2LGzaBN9/7/n3fq7hcxzZfgRfX18CAgLYtGmT54NIjpIlYedOcz2zTBkyzptnqnSLBJMk\nJEQqlyED9O5tCgZYIUuWLMyYMYMpU6bQvn17OvTvwF8XPHh9MLnSpoUhQ2DLFjIuXgw1a5ry5SJB\nJAkJIejd2xQHcFdJn4QIDQ3l4MGD7L6xm/yv52fueheX/na3p57ijxUroF49KFsWPvhARkUJIEko\nFVmxYgXVqlVjw4YNVocibCZPHmje3NzesFKuXLnQEzU9i/ak6+quDH1tqHc1YUybFgYOhIgIUzW2\nenX44Qero7I1SUKpSGhoKJcuXaJcuXJWhyJsqF8/09ng2jWrI4GJPSby86Cf2b93PxUqVODIkSNW\nh5Q4SsE335hifRUqmBYRMiqKlyShVGTHjh0ULlyYrFmzWh2KsKHixU3Nzk8+sToSo8AjBVizZg29\nevWiWrVqhIeHe1dX4DRpTGbfsQM++wyqVIFjx6yOyna8PgnF7azq/DyPUmqPlTHZRUxMDFOmTGHR\nokW3HpUrV771/OjRo2/bf8mSJSxYsIDhw4d713924TL9+pkSaXZZ9uJwOOjevTu7du1i6Zql5Oub\nj+M/Hrc6rMR58knYvBlat4bgYAgPB/n/dYtXJ6F4OqsCDAROuO1NS5Qwy809+ShRIkmhjhgxgpiY\nGMLCwggLC2Pv3r1UqVKF69evM2/ePLZs2XJr3z179vDkk0/Srl07smTJwrx581z1FRNepE4duHHD\nTNm2k0KFCrFl3RZaPNaC4ArBzJw507vK/vj4QN++sHu3aXVbqZJ7Ogx6IdsloeR0VlVK9QIWAO5r\n63j4sPkz0ZOPw4cTHebRo0dZuXIlnTp1AuDUqVP4+vpSsmRJfH196dixI/ny5bu1/+nTp9m4cSMA\nBQoU4MyZM675egmv4nDASy+Z0ZDd+Pn6MXXwVDZt2sTUqVNp0KCB91W7L1QIvv4aOnSAypXh7bdN\n1k/FbJWEktFZNVYtoCdQTinV3ENh29LOnTsJCAggc+bMAGzfvp2KFSty7dq1eGcbNWnShJ49ewJw\n+PBhKlSo4NF4hX20b2+q0mhtdSTxK1GiBDt37iQoKIiAwACGzRpmdUiJ4+MDffqYL/L69eYSnbdN\nvHAhWyUhkt5ZFef25lrr3sAurfVnHorZlrJly3arMOT169dZtWoVwcHBrF69Ot77PQ6Hg0yZMnHy\n5EmioqKoWbOmp0MWNpEhA/Tsad3i1YTw8/Nj5MiRzFo8i3ePvMtjAx7jx7M/Wh1W4jz+uElC3bpB\ntWowZkyqHBXZqoq21nq5s7NqrKxA3IJSN5RSdzW201p3vN/n8TmbwnuC1K9fn8jISNasWcOVK1cI\nDQ3l+PHjFCtWjPTp08f7muvXr7N48WLGjBnj4WgTLioqiitXrqTY79/FixdtcW7NmvlQrVoennvu\nN3LmTNi9l+++S8vQodk4ftyXbNmiyZYtmvz5bzJr1t+39nH1+T1d+Gkie0TSZmYbAp8NZEr7KYSE\nhLjs+ImVpPNr0IA0gYFkGzgQn08/5Z/x47lRrJh7ArQhWyWheEhn1STy8/Nj1KhR993nzhu7K1as\noHfv3vj5+bF+/Xpq1arlzhCTRDqreoa/PzRtCqtX52fIkAfvHxMDzz4LPXqYpTH//OPD33/DlSvg\n75/h1n6x53fuHMyZY4oLBATct1npg2PFn8PvHubrr7/m2WefpW7duoSHh9+6FO1JSf7++fub2SCz\nZ5OnTRt44QV45RVTKNUmpLOqdFZ1qQULFvDLL78wffp0fv/9dyIiIhg7diy1atXimWee4Z9//nnw\nQUSKNmSIuTSXEA4HfPmlubL00ENQuLCpXFO1avz737gB585By5bm92+XLub1yZm5XKNGDQ4ePEhU\nVBQBAQFs27Yt6QezgsNhMvm+fWZtUblysH+/1VG5nd1HQsuBWs7OqgBdrAwmJWnXrh3t2rW79Xnu\n3Ln59ttvLYxI2E2RIuaRUIkZzRQoABMnmo9/+MEkoBEjTJGBscloL5QtWzZmz57NypUrCR0VSqWn\nK7HstWX3vARtS488YtqJz5sHtWubSQxDh4Kfn9WRuYXtRkJa65PO2XBorWO01r211hWdD1luLEQK\nU7jwf0to3njDNcds3Lgxu2bsgh+gbNmy7Pe2EYXDYboN7t8PkZFmWJlC/0i0XRISQqRe9/pj/4MP\n4JdfEnesYo8WY82naxg0aBC1a9dm9OjR3PC22Wf+/rBqlSmKWq8eDB9uj+J+LiRJSAhhazdvwokT\nppN2u3bmlklCORwOOnToQGRkJJs3b6ZUw1J8uedLt8XqFg6HWby1f79ZuF6mDOzda3VULiNJSAhh\na2nSwLvvwk8/QenS0LCh6Ru3eXPCj1GgQAHWrVtHmZplqP9ZfVqMa8GNm142KsqfH5YvN/eH6tc3\n/6aAUZEkISGEV8iWzVyV+uknU/Xmx0SuTfXx8WHegHmsbbmWjT9vpEb9Gpw6dco9wbqLwwFhYXDw\noClp8fTT5maaF5MkJITwKn5+5p59165Je33toNr8b+L/qFu1LkFBQcydO9e7iqEC5M0LS5eaKYWN\nGsHgwXDVfSUz3UmSkBAixYiOhgULHnyVKm3atAwZMoT169cTHh5Os2bN+PW3Xz0TpKs4HNCqlRkV\n/fyzuWm2c6fVUSWaJCEhRIpx/jwsXGha+Eyd+uBkFBAQwJ49e3ii6BMUGFWAhUsXeiZQV8qTBxYv\nhpEjTZmLgQPh33+tjirBJAkJIVKMHDnMOs+lS03bnoQko3Tp0vHu2HdZ0nQJr73yGp06deL8+fP3\nfoFdtWxpRkWnTpkZHNu3Wx1RgkgSEkKkOOXK3Z6MEtKyvElIE/bv30+mTJkoVaoUX3/9tfsDdbXc\nuc3JjhkDzZtD//627+Lq9UkobntvpVSAUmqLUmq2UuoeVauEEKlFbDLq+MC6+kbmzJmZMmUKM2bM\noFOXTlQcVJE/zv/h3iDdoXlzOHTI3DeyeR1Ir05C8bT3LgecA24AqbdL1D2sWLGCatWqsWHDBo+/\nd1RUFJMnT6ZixYq0b9+eadOmMWXKFEaNGsXzzz/Pj4mdbytEIsRX1y462lT/jk+dOnXYuXcnZ2+e\nxX+kPx+t+8i9AbrDQw9BeDjkymV1JPdluySUnPbewDagO/A2MNBzUXuH0NBQLl26RLly5Tz+3n5+\nfvTs2ZOrV6/y/PPP06tXL/r06cPw4cN57rnn6Nq1K5cuXfJ4XCL1mj8fQkLuPaHskTyP8HP4z/R9\nqi/93u3H8OHDiYqK8myQqYCtklAy2nvH/p1TGkgD/OP8V8SxY8cOChcuTNasWS15/3379hEdHU1Q\nUNBt24sVK0b27NktGaGJ1KtdO1MNp1Ur0wPpXh22w7uGc2z+MQ4cOED58uU5dEg6yriS3Vo5xLb3\nnu/8/Lb23kqpW+2973hd7KD6BDAJiAJcVI/Xe8XExDB16lRy5MgBwJYtW6hcufJt+wwbNozRo0c/\ncBvAhQsXGJuAOvudO3dGKXXX9oiICMqUKYNvPI26rl69yh9/eOG1d+G10qY1C17btYMpU8yoqEED\n02IiU6bb982XLx+rVq1i9uzZhISE0P/l/gzoPwA/35TZXsGTbJWEktveW2u9A9iRkPdKaovh1ze/\nzhtbTH4bUXUEr1d7/a7ngXi3J/V18W1PiBEjRpA3b17CwsIACA8Pp0+fPree37p1K1rr214T37ZY\nWbNmTVASupeIiAjq169/1/bLly9z+vRpChYs+MBjSHtv72bX82vTBurVc7BkSUb++ecy95qhXbdu\nXYoXL07b8LZMenYSS3sv5fHHH7/1vF3Pz85slYTiYbv23q9Xe/2+CeFezyXndUlx9OhRVq5cSUSE\n6Qd46tQpfH19KVWqFGD+s8TExJAxY8Zbr4lvm6ucP3+eI0eOxJvE1qxZQ/bs2e8apcVH2nt7Nzuf\nn78/vPYaQLYH7OePXqR5d+K7NG7cmFGjRtGjRw8cDoetzy+53NXe2+5JKAJoACyV9t6Js3PnTgIC\nAsicOTMA27dvp2LFily7dg2Hw8HWrVsJDQ1l5syZt14Tu+3DDz+M95jnz5/nrbfeeuB7x3c5bvv2\n7SHRj08AAA5KSURBVOTJk4fChQvftj0qKoq5c+cyYsQI/FJo50jh/U6cgMce+2+WXdo0aXml3ys0\nrtuYjh07smLFCj788EMciWkvKwD7JyFp751E2bJlI3fu3ABcv36dVatW0bx5c1avXk2JEiUoWrTo\nbfsfPXr01rZ7FXPMli1bki/HRUREEBwcfNu2P//8kxEjRtCxY0dq1aqVpOMK4W4xMeZyXebMMG6c\nKdEWq1ixYmzfvp0xY8ZQvElxmrVtxkcvfSTJKBFsl4S01ieBW+29gd7WRuSd6tevT2RkJGvWrOHK\nlSuEhoZy/PhxihUrxrFjx4iOjubAgQOcOXOGDRs2cOXKlVvbTp8+zYYNG6hZs2ay4zh69Chr165l\n3bp1BAYGMm3aNKKjo7l69SoXLlygb9++8U5iEMIuHA7Ytg0+/BBCQ6FOHVOQIPaqm6+vLyNGjKBQ\npUJ0W9ON9QPWs2nAJgo/XPj+BxZGTExMqnvs3bs35kESsk9KUL169QRts4O9e/fGTJs2LWbChAlW\nh+I2Z86csToEt/L28zt/PibmlVdiYnLmjIlZsODu5w8fPRxT/pXyMfkK5ov5/PPPPR+gGzl/J7r8\n97Gt1gkJz5oxYwZ///03y5Ytu+82IYSRNSuMHQt79pgu23fKkSUHO8fu5JPZn/D888/To0cPLl68\n6PlAvYgkoVSsR48e7Nu3j2bNmt13mxDidoUKQZEi936+atWqHDhwgOjoaAICAljz9RrPBedlJAkJ\nIYSL/P67D7HLhLJmzcqHH37Ie++/R5PlTej5silbJW4nSUgIIVxk1y4/SpaEUaP+6yvXqGEjTr56\nkr9O/kVQUBCRkZHWBmkzkoSEEMJFGjS4yp49sH8/PPUULFtmpnj75/Vn8eLFDBs2jHr16jFy5Eiu\nX79udbi2IElICCFcqFAh00zvo49MBYbQUNM2wuFw0LZtW/bt28e2Hdvw7+XP57s+tzpcy0kSEkII\nNwgJMSOioUPBJ85v2ocffpi1a9ZSq3gtGi1vRIu3WxAd7ZJqZF5JkpAQQrhJ2rQQX0lEHx8fFvZf\nyPrW6zm0+xA1atTgxIkTHo/PDrw+Cd3R3ruYUmqaUmqyUuopq2MTQoh70RpqBNbgu8XfUa9ePcqW\nLcvs2bPvWTYrpbJd2Z7EiKe9d2/gDOCP6S2UZH5+fjKLxWZu3rxpdQhCuMRff0GNGlC9Orz9dhoG\nDRpEvXr16NChA0tXLOWtCW9RslBJq8P0CNslIaVUeeAtrXV1pZQDmAIEAFeBblrrn2L31Vr/CIQr\npeY5Nz0GvAYEAZ2AqUmNo2RJ8wPw3nvvkSdPnqQexuMuXrxIlixZHryjEMIyOXPC99+bGnSlSsHg\nwfDiiyXZvXs3HUd2pPyQ8sxvNZ/mzZtbHarb2SoJOdt7dwAuOTfdau/tTE7jgSZKqZFAYaCP1vqf\nOIf4DbgC/MV/Lb+TJXfu3Pzvf/9zxaE84sqVK/wbu0AhBbpy5QqPPfbYg3cUwuayZDElgJ59Fl54\nAWbNguXL/fhk1Cds376dzp07s3z5ciZNmnSrO3JKZKskRNLbe8eaDszEtPcedL83Smj3w+rVqydo\nP7tI6SOh2PNLqd0rU3pnTjm/u2XKZCp0r1+fDh+fKM6ejeHxxx/niy++MC0iihcnPDycqlWruilq\na9kqCbmgvXck5jLcA6XU7ocpubMjyPl5Ozm/e+vc+e5ts2bNYsOGDYS9Hkbe/XnZ+MpG8uSw5vaA\nuzqr2n12nNvaewshhDcoW7Ymez/by783/+WRUY+w4KsFVofkUnZPQhFAKIC09xZCpEbt2kH/Po+x\n6YUfGVByAP2792fIkCFcu3bN6tBcwu5JaDlwzdneOxzo9//27j1GyuqM4/h3sQpUbpZy64WAtjyY\ncLHQsGIvWLlsLdh/mrYRo4Ct1ZVAsalNbINtpJSERlJDwBIagQoUI5ZLLAQp1ri22uo2doOhT5c0\ntShoUxGwxdUK9I9zgNllZ3cH3svM+PskhJ13zrznPJzZeTjnPfOenNsjIpKpRx+FsWNh3Djo9eoS\nnm9oYt++fUyYMIGmpqa8m3fByuqaEGh7bxGRQj17hnvQ3XwzLFgA69YNYtWqLbz88jomT57MzAUz\nWXr3Urpf0r3zk5Whch8JiYgIMHw4bNsGy5bB8eM1zJ49m8bGRjYd3ETt9Fqam5vzbuJ5URISEakg\nM2bA9Onh56FDh3Jo+SFuveFWJk6cyMqVKyvutj9KQiIiFaxbt27Mnz+fhoZnWLt2LXV1dRw4cCDv\nZnWZkpCISBXYvXskI0c+x1WfmsYVP7qC+gfrK2KLCCUhEZEqMGcODBjQjbVrvstNg7aypnkNw74z\njJaWlryb1iElIRGRKtC7N9x/P+zZA/sbvsSVDYep/cCssv8+kZKQiEgVGT0ann4a7pr3Qbq/toi+\nffvm3aQOKQmJiFSZmhq45RZYvz7vlnSu7L6sWioz+wIw091vM7NvA1cBnwTWu/vP822diIh0pKJH\nQm13VnX3B4BvAXuVgEREyl/ZjYQucGdVgBuBX2fZZhEROT9lNRKKO6uuJo5sKNhZFbiHsLMqZnaf\nmW00s36xXOEuqp9z9yeyarOIiJy/chsJne/OqoX3qbgo7UaKiEgyasrtPkNxZ9Vfufs1ZrYa2Ozu\nu+Jz/wAuv9CN7RobG8sraBGRCjB+/PiazkuVptxGQm2lsrNqGv+QIiJSurK6JtQO7awqIlLFyn0k\ntAWYGndWBZiTZ2NERCRZZXdNqKs6W75tZrcRvjP0P2Cxu//GzPoDG4EewEFgjru3FCl7GfA3zo6+\ntrj78ozCSzS+WH4AYWQ5yt3fNbMewHpgIGHac5a7v1ENscVjrxD6D+BZd/9BFrHFupN8b94FfJ2w\n+GaHuy/Ks+9Sjm+nu98Xz1Et/TcXmAWcBBbFstXUf+fEF8/R5f4r9+m4jrS7fBvAzAYB84CJwBeB\nJWZ2MXAvsMHdJwEvArd3UHYcsNHdr4t/MktAURLx3RHLTwN2Ed70p9UDTe7+ecJqxIWpR3RWqrHF\nLzE3FvRdZh9gUVLvzeHAje5+NWHL+zozG0W+fQfpxTfNzEZVUf/1J7xPrwamAA/G01RL/7UbX6n9\nV8lJqNXybeDTBc9NAJ5x9/fc/RjQTMj6Z14D7ASmFik7BhgPjDezp8zsETMbnEVQBZKIb3L8+UT8\n+XB7549lp6QQQzFpxzYe+JiZPWlmj5vZiNQiaV8S8U0B/kn4IMDdTxGmz1uKlM1SWvFdTIivKvov\njm7GxsVUQ4A3256fCu6/DuIrqf8qOQn1AY4WPH7PzLoVee4toC9hpd3RDo4B/Cce3wfc6+7XAtuA\nrEdCScWHu+9x9zdp/aXePm3K9km09R1LO7aDwE/c/TpgCWHqI0uJxOfuJ9z9MICZ/RT4s7vvJ9++\ng/TjO0QV9B+Au5+MU1bPApvbOUfF9h+0iu8PnI2vpP6r5CTU0fLtY7Tu2D6ELF34mt4FxwrL9gaO\nAL8DnorHthBujJqlJOI70uachRcAOyubprRjawS2A7j774GPJNPsLkssPjPrbmYbgEuBue2cP+u+\na1s/JBffnfH5F6iS/gNw9xXAYGCSmV1L+DCviv6DM/ENIcQ3iRL7r5KTUEfLt/8EfNbMLjGzvsBI\nYG98zfRY5nqgAXi+SNlfAF+JZacQPtiylFR8hQpHC2fOH/9uWzZNacf2Q2BBPP9YwrRPlpKMbzvw\norvfGaesWp2f7PuuVf0pxVcV/WdmI8zssXjsBGGq8USbshXbf0XiO0mJ/VcNq+PGxENzCP9Ize7+\nuJl9A7id8OG02N23mtlAYB3QC/g3YQuIt4uUHQY8FM/9X8IKktczCi/R+ArO+XdgpIfVcT1j2SHA\nO7Hsv6oktn6EKYBehBU+c9399Eqd1CUVH1BHWJH0XCx7inAhuYmc+g4yie+vwAYqvP/iZ8tCwgf+\nScLqvx/n+buXUXwl/f5VbBISEZHKV8nTcSIiUuGUhEREJDdKQiIikhslIRERyY2SkIiI5EZJSERE\ncqMkJCIiuVESEhGR3CgJiWTAzH72fq5fpJhy31lVpKLFvVjqCbdFWfB+q1+kM7ptj0gGzOzJeGv7\njsqcvpPxXncfEx/PcPcdBWVGEG54+Zf43LvnU7+ZNQCfIdyvrbe7Hy8tIpFkaDpOpLzcBExq7wkz\n+zjwBPAS8OWuJqAibuDsXeJFcqPpOJEOmNlOwiZeLxH+07YQ2ES4vf2HgW8CV7r7qwlVeTRu0te2\nHQOB3wKvEEZALRdSibsfMbPDnZcUSZeSkEgRZjYU2OHuy+PjWuD7QL27H43HegLXxZ0pT89t1wDv\nuPumhNrRB9hF2Ejs+sKpMzPrAXytzUsSrV8kTUpCIsXVAasKHk8DXjidgKK97v5wF85V03mRdl0K\n7ABGAyPc/a3CJ+OI6Jcp1i+SKl0TEinC3Ve3ue4yjTAiKSyzorPzmNlc4BNmdo+ZDS6xGQ8A/YDX\ngcUlvjaJ+kVSpZGQSBeYWW+gFvheqa+NiarTZFXEMcL28rXAFjPb6u6PZFi/SKo0EhLpmsmEbd7/\nmHG9d7v7a+6+jbAgYoWZDcm4DSKpURIS6ZqpwB53P9lpyWSdKPh5Xnz8UMZtEEmNkpBI15xzPShr\n7v4GMBeoM7P6PNsikhRdExIpIn45dDbwUeBywlLsgcAyd387gyacczsTd99sZo8BS81st7vvz6Ad\nIqlREhIpwt0PAIviwztyqP+iIse/mnVbRNKi6TiR8tLPzD6UdiVmdhnQP+16RDqjkZBI+TgFPEy4\nRdCYlOvaDlxDO1N+IlnSXbRFRCQ3mo4TEZHcKAmJiEhulIRERCQ3SkIiIpIbJSEREcmNkpCIiORG\nSUhERHKjJCQiIrlREhIRkdz8HxTBeNyArqcBAAAAAElFTkSuQmCC\n", "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "D, dD = np.array(Dlist), np.array(dDlist)\n", "d11_T = np.vstack((Trange, dD[:,0])).T\n", "d11pos = np.array([(T,d) for T,d in d11_T if d>=0])\n", "d11neg = np.array([(T,d) for T,d in d11_T if d<0])\n", "fig, ax1 = plt.subplots()\n", "ax1.plot(1./Trange, D, 'k', label='$D$')\n", "# ax1.plot(1./Trange, dD[:,0], 'b', label='$d_{11}$')\n", "ax1.plot(1./d11pos[:,0], d11pos[:,1], 'b', label='$d_{11}$')\n", "ax1.plot(1./d11neg[:,0], -d11neg[:,1], 'b--')\n", "ax1.plot(1./Trange, dD[:,1], 'r', label='$d_{12}$')\n", "ax1.plot(1./Trange, dD[:,2], 'g-.', label='$d_{44} = D$')\n", "ax1.set_yscale('log')\n", "ax1.set_ylabel('$D$ [cm$^2$/s]', fontsize='x-large')\n", "ax1.set_xlabel('$T^{-1}$ [K$^{-1}$]', fontsize='x-large')\n", "ax1.legend(bbox_to_anchor=(0.15,0.15,0.2,0.4), ncol=1, \n", " shadow=True, frameon=True, fontsize='x-large')\n", "ax2 = ax1.twiny()\n", "ax2.set_xlim(ax1.get_xlim())\n", "ax2.set_xticks([1./t for t in Tlabels])\n", "ax2.set_xticklabels([\"{:.0f}K\".format(t) for t in Tlabels])\n", "ax2.set_xlabel('$T$ [K]', fontsize='x-large')\n", "ax2.grid(False)\n", "ax2.tick_params(axis='x', top='on', direction='in', length=6)\n", "plt.show()\n", "# plt.savefig('Fe-C-diffusivity.pdf', transparent=True, format='pdf')" ] }, { "cell_type": "code", "execution_count": 9, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "(430.0, 420.0)" ] }, "execution_count": 9, "metadata": {}, "output_type": "execute_result" } ], "source": [ "d11pos[0,0], d11neg[-1,0]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Activation volume. We plot the isotropic value (change in diffusivity with respect to pressure), but also the $V_{xxxx}$, $V_{xxyy}$, and $V_{xyxy}$ terms. Interestingly, the $V_{xxxx}$ term is negative---which indicates that diffusivity along the $[100]$ direction *increases* with compressive stress in the $[100]$ direction." ] }, { "cell_type": "code", "execution_count": 10, "metadata": { "collapsed": false }, "outputs": [ { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAaAAAAEzCAYAAAB3xNe0AAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAIABJREFUeJzt3XmcXFWd//9X9ZruTncSQgJESAZUPgRDcIzKOuCwBEQH\n+X6dUQf8qoyIoqKAgyiO+BtH3EYQwR2GATTqyDozCAgB1BAETdBAED+ALAE6C2TrTu9L/f44t7pv\nV1cv1anqW139fj4e53HvPefcW+f07a5P3+3cVDqdRkREZLJVJN0AERGZnhSAREQkEQpAIiKSCAUg\nERFJhAKQiIgkQgFIREQSoQAkIiKJUAASEZFEVCXdAJHpzsz6o9n17r40Wn67u98Rq3Mg8BvgUWAO\nsAxIA43u3j7ZbRYpBB0BiZSGM4BjcxWY2X7A3cCfgFOBE4F3Tl7TRIpDR0AiBWZmdwKtwOOEf/I+\nD/wMWA/sCZwFLHb3l2Kr7XT37Tm2NR9YCbxIOCrqBDrNbFtxeyFSfApAIgVkZguBO9z9qmj5MOBi\n4Bx33xnlkRV8RtpWE/BLYAfwVp1qk3KjACRSWCcBP4gtLwfWZIJP5C/j2E4DcAdwCHCgu7cWroki\npUEBSKSA3P3qrKzlhKOYeJ3vjGNT3wK2AZuBS4F/LEgDRUqIbkIQKRIzawQOIysAjVMLcALwUeDd\nZvbuQrZNpBQoAIkUz/FAG/DwBNa90N03uft/E25g+I6Z7VPQ1okkTAFIpHhOBO519/4xaw7XF5s/\nN1q+tiCtEikRCkAixTPs+s9EuPtW4GPASWZ2zm63SqRE6CYEkQKKHhr9APAq4ADguOhZnsvdvWOc\nm0lnZ7j7TWZ2M/B1M7unUO0VSVIqnR72uy4ikyjX0DvjWOdY4D40FI9MYToFJ1IaZpvZHuOpaGZz\ngLlFbo9I0ekUnEjy0sCPCEP3LB1H/f8BjiTHqTqRqUSn4EREJBE6BSciIonQKbgCiwaf/Kq7/62Z\nvR64EugFuoD3ufvLZvYh4GygB7jU3X9hZnOBnwAzgGbgTHfvNLON7r5PtO2DgP8GPuLu94/w+TXA\nfxLuwNoJfJxwveBb0efd4+5fNLMU8F3gUKATOMvdnzGz+4EPu/uTZjYTuJ0wuObXC/7DKjAz+wzh\ndQXVhL79BrgO6Ce8a+djUb1LgLcRfh7nufsaM/tP4KfufreZVRL2xcvu/vHJ74lIbmZWAVwNGOF7\n5UzCgcR1TMHfcx0BFZCZXUj45aiNsq4APubuxwG3AheZ2V6EBwuPAE4GvmJm1cAlwAp3Pxb4I/Dh\naBvpaNuvi7bx/0YKPpEPAa3ufkT0Od8Gvge8x93/BjgsCoynAbXufiTwWeDyrL40AncSflmnQvA5\nFjgi6s9bgIWEPl0c/UwrzOwdZvbXwDHufhhhfLXvZm2nCvg58LSCj5SgvwPS7n408AXgm0zh33MF\noMJ6Gvg/seV3u/tj0XwV4UjjzcAD7t7r7i3AU4SjkKOBu6K6dxLGAQNImdlS4BbgH9z9d2O04eBo\nfdz9KeBNhEDzXFT+y2jbA5/n7g8T3rCZMQe4B/ihu8dHdi5lJwHrzew2wkX624E3uPuqqPxOwsgE\nRxNe7oa7vwBURkefEI4+bwH+4O6fm8zGi4xHNDTT2dHiImATU/j3XAGogNz9VsJhcWZ5M4CZHUl4\nkv2bQBPh1FhGKzALaIzlZ/KI8v+TcBg9ZxzN+CPw9uhzD4+2s2uMzwPoiw7vAX4MdBMeppwq9iQE\n0b8HzgFWMPT3e6R+x3/WVwL1wL7FbqzIRLl7v5ldR/h9vRlIxYqn1O+5AlCRRaMYfxc4JRpSpYUQ\nhDKagO1RfmOU10h4CRmEU3DvAP4fcL2Z7TnGR14LtJrZfYTzv+sI75bJaMzxeQAVsTHLPh2te6aZ\n/c04u5q0rcAvoyPLJwlHm7Ni5fF+Z//8Mz/rbxGGz1lqZmcUv8kiE+PuHwAOBK4B6mJFU+r3XAGo\niMzsvYQjn7e4+/NR9u+Ao82sxsxmAQcRXtW8mvClD/BWIHNIvcvdX3T3PxCu5/xkjI99E+EU33HA\nbcCTQLeZ7R/deHBStO0HgVOidh4OPBbbxuPRC9DeB/x4HEGvFDxAuKaGmS0gBN17o2tDMPgzfRBY\nbmap6O2lKXfPvN768SgIn0EY8uagSe2ByBjM7L3RzTYQ/snqA9ZM1d9zBaAiiU5nfQuYCdxqZveZ\n2Rei03JXEr4wVxIuHnYTXjr2HjNbBRxOCDYQe9jQ3S8H+s3s86N89FPAOWb2IPBF4HzgI4TA9RDw\niLv/nnBDQ5eZrQYuA87L8XkPE97uOVbQS5y7/wL4g5n9jnCn4DnAp4B/jfpYDdzk7o8Q/kB/C9xI\neN8ODO33s4SjwJ+b2YzJ64XImG4B/trMfk243vMJwj+5U/L3XA+iiohIInQEJCIiiSi5B1FHekAy\nq848wjWTJdHpK8zsRcL1DoDfJn17oYiIjK7kAhCxBySjUQUuj/IAMLPlwFeB+bG8VwNr3f0dk91Y\nERGZmFI8BZf9gOQbs8r7gOOBbbG8ZcC+0YX+283swElpqYiITFgpBqDsBzV7Yw9I4u73uvt2hj58\n1Qx8Obr1+CuEBylFRKSEleIpuNEekIyL3763lmgEAndfHT0HMszatWt1y5+IyAQsW7YsNXat/JRi\nAFpNGErmphwPSMbFfxhfIDwJ/+9mdiiwYaSNL1u2bKSixFx11VWce+65u72d5uZmFizIGXuL/tmT\nYSL9G0kp9ruQ/StF5dy/Uu1boX7P165dW4DWDFeKAehW4MTooSoIw8GcDzzl7rfH6sWPZr5KeGI/\nM/T4ByalpSIiMmElF4DcPU14ij3uyRz1DojN7yAagFNERKaGUrwJQUREpgEFIBERSYQCUAl485vf\nPC0/O0nTtd8yvZT677kCUAk47LDDpuVnJ2m69luml1L/PVcAEhGRRCgAiYhIIhSAREQkEQpAIiKS\nCAUgERFJhAKQiIgkQgFIREQSoQAkIiKJUAASEZFEKACJiEgiFIBERCQRCkAiIpIIBSAREUmEApCI\niCRCAUhERBKhACQiIolQABIRkUQoAImISCIUgEREJBEKQCIikoiqpBsgImUqnYb+/pDSaejshPb2\noXmZ+ezlzHyuvNHKs/Pj09HKctXJnh+lrG7bNpg1a3zr5VrOlcZTZ7R1ZsyACy+E2tqkfxNGpAAk\n00Pmj7O3d/TU1zd6Xl/f4HJ8mj2fnXp7w+ePVN7XR9POnVBfH5az68aXM1+48fzxTLPXy5VGK8sn\npdPh515RAakUpFLsU1kZljN5uebjy9nTscpHm0607mjzseXazs6w7/JYZ9Q0njq51qmsHFyuqwvT\nEqYAJPlJp8OXaWcndHUNT93dQ+czy5n5eOrpGT4fn2bNz921K/yRZcp6e4dOs/MyKbNcUQFVVeGP\ntLo693z2clXV0Pn4dLT5eIrnZ74k4p81YwZUVtJXUwNz5gzWya4bX858Gcfz49Psupkvplzrjjc/\nHgBGy49/gcZsbG5mwYIFCf3iFteO5mbqy7RvxaQAVA76+qCtjYpNm6CtLaT29sFpPHV0DE7j852d\nw6fZqasrTCsqwpdmbe3wVFMzdLm6ejC/pib3ckPD4Hz2NJZ2tbRQu/feg3lVVcPn43nxYFJVFdpd\nwtqam5mlLzGZRhSAktbSAi+9BDt2wM6dYdrSEuZbWqC1NUwz85m0a9dg6u6G+nrm1dVBY2P4Qm9o\nCKcEGhrCoXhmOZP22CMEkfr6MK2rG0y1tWGaCTLx+dra8GWegK7mZtAXtEjZKLkAZGYp4LvAoUAn\ncJa7P5NVZx6wGlji7t1mNgP4MTAfaAHe7+5bJ7flE/SJT8BDD8Hs2eEiZiY1NYW0aFGYNjaG1NQE\nM2eG+ZkzBwNLKsXmMj7FISLlp+QCEHAaUOvuR5rZYcDlUR4AZrYc+Coh2GScAzzq7l80s3cDnwfO\nm8Q2T9x11yXdAhGRRJTiSfGjgbsA3P1h4I1Z5X3A8cC2XOsAdwInFLmNIiKym0oxADUBO2PLvWY2\n0E53v9fdtwOpEdZpjZZFRKSEleIpuBagMbZc4e79OeqlR1inEdgx0sabm5t3u4GlqrW1Vf2bwtS/\nqauc+1ZMpRiAVgNvB24ys8OBx0aoFz8CWg2cAqyJpqtG2ng5X6RvLvObENS/qa2c+1fOfQPYuHFj\nUbZbigHoVuBEM1sdLZ9pZucDT7n77bF68SOg7wHXm9kqoAs4fXKaKiIiE1VyAcjd04S72uKezFHv\ngNh8B/CuIjdNREQKqBRvQhARkWmg5I6ApptvfOMbPPHEEzQ0NFBfX099fT0NDQ3U1dVRX19PXV3d\nkPkZM2YMTLOTiMhUogCUsCOOOILZs2fT1tZGW1sb7e3t7Nixg/b2djo6OgammdTZ2Tkwjaeuri6q\nqqqora0dlmbMmEFNTc3Ack1NzcByZj6Tqqurh8xnljPzI6WqqqqBaWY+vpydqqurqaysHFiurKwk\nVeIj94pIYSkAJeyoo47iqKOO2u3tpNNpnnvuOebOnUtXV9dAUOru7qazs5Pu7m66uroG8jLTnp6e\ngfnMcma+s7OTnp6egbzMfDz19vYOm8+e9vX1DSxnUia/t7eXvr4++vr6qKysHBaUMtNMgKqpqRlS\nLzOfK1VVVVFRUZGzLJOfXR7PH2l+pPLR8saTduzYwZ577jmwnEqlhpSPtpxrfrS8fKaZFF/OVSe7\nfnaeSDYFoDKRSqWora2lqWlqPoObTqcHgtFI0+bmZubOnTsQsOLBK1fq7e2lv78/Z1l/f/+wssxy\nvCx7Pj7t6ekZUp5Op4eUp9PpYeWjLXd0dFBTUzOsTrxef3//sOVMXq56o+XFy8bKzy7Prptdlp0y\nxgpso6Xx1itEyrR1vGXd3d3U1tbu1jbiednlE1mur6/n+9//PvX19UX9290dCkBSElKp1MBpu5FU\nVlaW9bMW5fwsSTqd5qWXXmKfffYZMZCNJ+VTd6Ip097xlgFs2bKFPffcc7e2MVL5RJdnzJhBbQm/\nDRUUgERkEmSOYCorK5NuSlGU8z8PxaTbsEVEJBEKQCIikggFIBERSYQCkIiIJEIBSEREEqEAJCIi\niVAAEhGRRCgAiYhIIsZ8ENXMJjSOg7u3T2Q9ERGZHsYzEkLrBLabHue2RURkmhpPkEgB7wS2jXOb\nc4EbJ9wiERGZFsYTgFYD97v7jvFs0MzmAA/uVqtERKTsjRmA3P1vsvPMbDFwEHA3MB94zt3TUf3t\nwLB1RERE4vK6C87MmszsDuBxwmm2vYBvAX8wMw0FKyIi45bvbdiXAbXAvkBHlHcusAu4ooDtEhGR\nMpdvAHobcKG7N2cy3P154OPACYVsmIiIlLd8A9BMBo98srejh1pFRGTc8g0adwGXmFnm5oW0mc0D\nvgHcU9CWiYhIWcs3AJ0LLCQ8E1QPrAQ2ALOATxa2aSIiUs7yGq3A3TcDR5nZ3wIHR+s/AdyTuQ1b\nRERkPCY6XM6fgL/ElvczM9x9QwHaJCIi00BeAcjMTgWuIQy3E5cijP9WWaB2iYhImcv3COgKwo0I\nV5L7briCMLMU8F3gUKATOMvdn4mVfwg4G+gBLnX3X0RDAD0JPBZVu9XdrypWG0VEZPfkG4DmAV9y\n9yeL0ZiY04Badz/SzA4DLo/yMLO9CDdDvIFwI8QDZnZ3tPwTd9fNECIiU0C+d8HdCJxajIZkOZpw\npIW7Pwy8MVb2ZuABd+919xbgKWApsAxYZma/MrP/MrO9J6GdIiIyQfkeAV0CrDOz04FngP54obu/\nq0DtagJ2xpZ7zazC3ftzlO0i3Ab+BLDG3e+L2ncV8A8Fao+IiBRYvgHoWqCPcK2laNeAgBagMbac\nCT6ZsqZYWSOwA/gdkHkL663Av+bacHNzc67sstDa2qr+TWHq39RVzn0rpnwD0NHA0e7+SDEaE7Ma\neDtwk5kdzuCNBRACzZfMrAaoI7wWYj1wA3Az4TThCcDaXBtesKB8B+1ubm5W/6Yw9W/qKue+AWzc\nuLEo2803AD1J+NIvtluBE81sdbR8ppmdDzzl7reb2ZXAA4Tbvy92924z+wxwrZmdA7QBZ01CO0VE\nZILyDUBfAq43s6uAZ4HeeKG731GIRkWjKpyTlf1krPw/gP/IWuc54LhCfL6IiBRfvgHo59H0mznK\n9CCqiIiMW75jwemVCyIiUhAKKCIikoh8x4J7mXCqLSd3n7/bLRIRkWkh32tA/5xj/VcDHwA+W4gG\niYjI9JDvNaDrc+Wb2e8JL6TLWS4iIpKtUNeA1gGHFWhbIiIyDeR7DejgHNmNwKeApwvSIhERmRby\nvQa0nnATQior/wXCdSAREZFxyTcA7Z+1nAa6gc3R6AUiIiLjku9NCM8XqyEiIjK9jBmARrjuk5O7\n/2n3miMiItPFeI6ARrruk01jwYmIyLiNJwBlX/cRERHZbWMGoFzXfczseOB1hOeIngDudffe7Hoi\nIiIjyfc5oL0JL4tbBjxHCEALgT+b2QnuvqXgLRQRkbKU70gI3wL6gP3d/UB3fw1wALADuLzQjRMR\nkfKVbwA6GfiEu7+UyXD3F4ELgFMK2TARESlv+QagDnK/jqEf3QEnIiJ5yDcA/RK43Mz2ymRE14Uu\ni8pERETGJd+heC4E7gOeN7PM3XGLgEeB0wvZMBERKW/5DsWzxcwOBd4KLAY6gSfcfWUxGiciIuUr\n39uw1wMrgJ+5++3FaZKIiEwH+V4Dug44DXjKzB40s4+Z2bzCN2v6eGb7M/zp5T+xvWM76bQGFBeR\n6SPfU3DfAL5hZvsD7wb+Cfimmd0HrHD3HxWhjWXtO7/7Drc/dTubdm2iq7eLvWbuxfyG+ezVEKbz\n6ucxr2HewHTP+j2ZWzeXPev3pKm2iVRqrCH6RERKU743IQDg7s8CXzWzbwMfBr4AnAgoAOXpspMu\n47KTLgOgo6eDTbs2saVty5DU3NrMus3reKX9FV5pf4Wt7Vt5pf0VOno7mDNjDnvU7cEedXtQl6pj\nwZwFzK6dzZy6OcyZMYfZM2Yza8asMK2dxawZswamNZU1CfdeRKazvAOQmTUAfwe8i/BgajNhhIQV\nhW3a9FNXXcf+c/Zn/znjG/+1u6+b7R3b2dqxlW0d23j6paepqK9ge8d2tndu5/mdz7Nu8zp2dO5g\nZ9dOdnbuZEfnDlq6WtjZtZPKVCVNtU001TbRWNsYpjWNNNY2MrN6Jo21jTTWNDKzZuZAaqhpGJyv\nbqC+up6GmgYaqhtoqGmgqmJC/9OIyDSU700INxOCTitwI3Ccuz9UjIbJ2Goqa9hr5l7sNTM8lnVA\n1QEsWLBgXOum02k6eztp6WoZSK3drbR2tQ5Md3XvorW7lc1tm/nL9r/Q2t1KW3cbbT1ttHW3sat7\n18B8W08b7T3tVKYqaagJgSme6qrqqKuuG5yPlrOnM6pmjJhatrXQVts2sFxbVUttZS01lTU6FSky\nBeX772ob8E7gHnfvK0J7ZJKkUqnwxV9dNxDAdlc6naarr4v2nnbautvo6O2gvaedjp5oGlvu6O0Y\nMm3taqWjt4PO3s6BaTx19HTQ1tlGb6qXrt4uuvq66OztpKu3i57+Hmoqa0JQqqyltioEpcx8Jkhl\n8jNlmfmRUnVFdZhWVo+4XF1ZPa5pVUXVkPnKCg0cIpLvTQjvK1ZD4swsBXwXOJTwrNFZ7v5MrPxD\nwNlAD3Cpu//CzOYCPwFmEE4LnununZPRXglSqdTA0ckedXsUfPvNzc05j/D60/0DQSkz7e7rHpbX\n09czsNzT3zNQp7uvm+6+bnr6ewbK2nvaB+Z7+nro7u+mp69nYL3uvsHlsaa9/b3D8lKkqK4Mwaiq\noorqimoqqKCmqmYgWGXK4nUy85UVlcPKB8pSlcPmM/VHW65MVQ7Jz+Rll090+vKOl+mY0UFlRSUV\nqYqB/LHmK1IVOsItU6V6wv40oNbdjzSzwwgjbZ8GEA0DdC7wBqAeeMDM7gYuIdyJd4OZXQR8BLgi\nkdbLpKpIVQwczU0Vff19ITDFAtRLG19i7ry59Pb3DpTlqteXDnnxssx8PGXXy+Rlljt7Owfz+/sG\nygbWy17OqpPvtLunm1RFiv50/0D+eObTpEmRCgEpKzDlkxfPz5Wy6w4rj62baU8mdXV20VDfMGyd\n7Hq5UiqVyll/d/Mbahp4z5L3lPR12VJt2dHAXQDu/rCZvTFW9mbggegFeC1m9hThSOlo4NKozp3R\nvAKQlKTKivBffi21A3k99T0smDW+a3hT0UhHsGNJp9P0p/sHUl+6b3A+Cm5DymMBLLOcJj2Qnytl\nb2NYefQ58bakGZzfum0rs2bPGvJ58fWz+xD/3HQ6PaR+vG5vujfndkaqH8+vrazlnYvfqQA0AU3A\nzthyr5lVuHt/jrJWYBbQGMvP5InIFJdKpcJpuRIecH+iwXW6y/cuuErg80Czu/8wyvstcAfhWkx/\ngdrVQggoGRWxbbcQglBGE7A9tk5XNN2Ra8PNzc0FamLpaW1tVf+mMPVv6irnvhVTvkdA/064C+7s\nWN73CQ+i1gL/UqB2rQbeDtxkZocDj8XKfgd8ycxqgDrgIGB9tM7bgOsJg6WuyrXhcv4vpdz/C1P/\nprZy7l859w1g48aNRdluvmPB/SPwj+4+8O4fd78e+ABwZgHbdSvQZWarCe8aOt/Mzjezt7v7ZuBK\n4AFgJXCxu3cTrvm8x8xWAYcD3y5ge0REpMDyPQKqJ1xfybaVAl5zcfc0cE5W9pOx8v8A/iNrnS2E\nIx8REZkC8j0Cuh/4upkNPORhZnMIRx+/KmC7RESkzOV7BHQu4bTXS2b2HJAivBH1aeDUgrZMRETK\nWr4jIbxgZocAJwAHA92EU2N3F/AOOBERmQbGPAVnZqeYWXVmnhB8AP5EOPKpAE6OykRERMZlPEdA\ntwN7A1ui+ZGkoYSfFJNJ8dhjj9Hd3V207RfrdtBSof5NrpqaGg455JCkmzFtjRmA3L0i17xILt3d\n3SxbtizpZoiMy9q1a5NuwrQ2oaF4ogFBa7Pz3X3DbrdIRGQSXXHFFcybN48zzjgj6aZMO/kOxXMq\ncA0wN6sohU7BicgUNH/+fLZs2ZJ0M6alfI+AriCMUn0l0FH45oiIyHSRbwCaB3zJ3Z8cs6aIiMgo\n8r2p4Eb0wKmIiBRAvgHoEuCzZvaImd1kZj+Pp2I0UGSytLS08OUvf5nFixfzwQ9+kHvvvXegrK+v\njwsuuIAjjzySr371qwm2cvrQ/ih/+Z6CuxboI4x+oGtAUlaampq4+OKLWblyJcuXL+f4448fKKus\nrGTp0qVccMEF7Lvvvgm2cvrQ/ih/+Qago4Gj3f2RYjRGpBQsXLiQ5557bkjeI488wsKFC/VllwDt\nj/KV7ym4JwkvgRMpW4sWLWLDhsFH2nbt2sXDDz/McccdV5DtP/jgg6xYsaIg25oOxtofl1566ZD6\nN954IytWrOBf/uVf6Ovrm9S2Sn7yDUBfAq43s0+a2anROHEDqRgNlPKwZMkSUqnUpKUlS5ZMuK2L\nFi3i+eefH1i+7rrr+MAHPgDAXXfdNeb6Y9V505vexG233Tbh9hXEkiWQSk1umuA+GWl/9PT0cMMN\nN/DrX/96oOz3v/89r33taznjjDNobGzkhhtu2O0flRRPvqfgMjcafDNHmR5ElRGtX78+6SaM26JF\ni3jhhRcAuPvuuznmmGOoq6ujubmZe+65h5NPPnnEdcdTp7q6mrq6hE8klMH+AHjf+97HypUrB+q+\n+OKLPPvss7z+9a9nv/3245lnnkmkzTI++b6OQWPBSdlbtGgR3d3drFmzhs2bN7N8+XIgDLT6+OOP\nc9NNN/G2t70NgJ/97GcDX5BnnHHGsDqbNm1ixYoVLF26lPb2dt7znvck2bUpaaT9kctpp51Ge3s7\nEP7pKdRpUymOCQUUMzvezD5hZueZ2UlmNqEx5URK0cKFC4Fwqic+PthJJ53E/Pnz+fu//3vq6uq4\n5ppreMMb3sBxxx3HIYccwve///1hdfbff38+/elPc+qpp3Lrrbcm1aUpbaT9kUsqlaKhoYHnn3+e\n7u5uTjjhhFHrS7LyHQtub+A24A3Ac4QAthD4s5md4O4aUEmmvJqaGhYvXsxFF11ERcXQ/9HS6TQA\nHR0dPPbYY5x22mkAzJ07l3Xr1g2rU1dXx5o1a2htbaWzs3MSe1E+RtsfufT09PDzn/+cL3/5y5PQ\nOtkd+R4BfQvoBfZ39wPd/TXAAcAO4PJCN04kKbfccgv77bffsPxUKkVvby/r16/nwAMPpLm5GQjv\nuVm8ePFAvUydq6++mqeffpqTTjqJ+vr6gfqZICXjM9L+gOE/y9tuu41zzjmHmpoa7rnnnslonkxQ\nvgHoZOAT7v5SJsPdXwQuAHQXnJS9gw46iF/84hccfPDBfOITn2DdunXccccdrF+/nnPPPReAxYsX\nD9R59atfzQsvvMBvfvMbXvva13L//ffzm9/8hmeffZZVq1Yl3Jupb8WKFWzYsIEf/OAHvPzyy6xe\nvZqvfOUrnHjiiRxxxBHs2LEj6SbKKPK9dtNBuNstWz+6A06mgYsvvnjI8tlnnz1qneOOO27gQvgx\nxxwzkP/AAw8UqYXTyxlnnDHkutC8efN45BE9Jz9V5HsE9Evg8uiFdMDAdaHLojIREZFxyfcI6ELg\nPuB5M8s8GbYIeBTQ6wRFRGTc8n0OaIuZHQq8FVgMdAJPuPvK0dcUEREZKt/bsO8D/q+73w7cHsuf\nB9zl7ssK3D4RESlTYwYgM3sLcHC0eCxwtpntyqq2GHh1YZsmIiLlbDxHQFuBfwZSUfo44Z1AGWlg\nF/CpgrdORETK1pgByN0fIzxsipndTzgFt71YDTKzGcCPgflAC/B+d9+aVecS4G1AD3C+u//ezP4a\n+F/CKyMAvufuNxarnSIisnvyvQnhb3Plm1kNsMzdf1uANp0DPOruXzSzdwOfB86LfdZfA8e4+2Fm\nth9wM/BmwvBAl7l7rpG6RUSkxOR7E8KbgR8ASxj+DFE63+2N4Gjga9H8nYQAlF1+N4C7v2BmlWY2\nF1gGHGhQWN6aAAASsUlEQVRmpwFPAZ9097YCtEdERIog34BxJeGa0HuA64EPAfsBnwPOyvfDzeyf\ngPMZHF0hBWwCdkbLrUBT1mpNwCux5VZgFvAwcLW7/8HMLgb+P8JzSyIiUoLyDUCHAoe5+6Nm9kng\nZXf/qZltIZwmy+uai7tfC1wbzzOzm4HGaLGRMNBpXEusHEJA2gHc5u6ZwHUrIViKiEiJyjcA9RIC\nAISL/UuBlcD9wBUFatNqwsCma6Jp9oiNq4GvmdllhKOvlLtvM7OHzOzj7r4GOB5Ym2vjmdGIy1Fr\na2tZ96/YWlpa+Pa3v82PfvQjjjzySE4//XSOP/54APr6+rjwwgt56KGHOPXUU/nMZz6TcGvL32Tt\nj9bWVtrb23frb0d/exOTbwB6GPiomX0GWAecSngNwxKgu0Bt+h5wvZmtArqA0wHM7GvAje6+Jir7\nLeGU3Uej9T4CfMfMugin8YaPEgksWLCgQM0sPc3NzYn3b+PGjYl+/u5oamri4osvZuXKlSxfvnzg\nyw6gsrKSpUuXcsEFF7Dvvvsm2MrpYzz741Of+hSvetWrdutzGhsb6ejo2K2/nVL42yumYv1d5xuA\nPgvcAWwBrgE+bWbPAHsBVxeiQe7eAbwrR/5FsfkvAl/MKv8jcFQh2iDT28KFC3nuueeG5D3yyCMs\nXLhQwScBo+2P3Q0+kqy8RsN2998D+wM3uPsO4E3AN4DMzQQiU96iRYvYsGHDwPKuXbt4+OGHB16r\ncOmllw5bJ1eeFMZY+wPgc5/73LD1cuVJaRnPUDyPAsdmHj51912EkQ9w903Ad7Pq7wH8yt2XFr65\nMlUtWQKPPz55n/e618H69RNbd9GiRaxdO3gJ8brrruODH/wgPT09/PSnP+XXv/71wJdbrrypYMl3\nl/D4y5O4Q4DXzXsd6z+a/04ZaX9krFq1Cncfsk6uPCk94zkFtwQ42cx2jlkzmA28buJNknI00WCQ\nhEWLFvHCCy8AcPfdd3PMMcdQV1cHwPve9z5Wrhwc/L26unpY3lQwkUCQlNH2R2trK+l0mvr6+oH6\nufKkNI33GtCKPLerF97LlLVo0SK6u7tZs2YNmzdvZvny5Uk3aVobbX+sWrWKU045hauvvnpY3jXX\nXJNEcyUP4xkLLt+3popMaQsXLgTCqZ4rr9TjZEkbaX/8+c9/5qCDDhpSN56XTuv/4FJXiKFzRMpK\nTU0Nixcv5qKLLqKiQv9/JW2k/eHupNNp1q1bx0svvcTKlStpb2+nv7+fdevW8eKLL7Jy5UpOOOGE\nBFsvoxl3ADKzy9xdr1yQaeGWW24ZsSzXf9b6b7u4cu2Pd7zjHQPzV1111bBAkytPSks+/96918wW\nFa0lIlPAihUr2LBhAz/4wQ94+eWXR8yTyfPDH/6Q7du3DwlSufKk9KTG+5+bmf2UMBLCX9z9f4va\nqiJZu3Ztetmy8n1reCk8jb127VrK+Wcs5WXt2rW4O1u2bOG8884be4URlMLfXjFFf9epQm83ryMg\nd78C2GpmnzMz3eMoIiITNu4A5O590fRBwkjT/xy9H0hERCRvE7rFx91bo/HY9jWzc81MtwqJiEhe\nxh04zOzU7Dx3vwX4A3CTme1fyIaJiEh5y+c5oI+a2d7Aa4HXRNMDgNqofBHhtdgiIiJjyicAnQgc\nBDwNPAU8EE2fJtwZ11X45omISLnKJwBd6e565YKIiBREPgHoX4vWCikbNTU1Q4bOFyllfX19STdh\nWht3AIpeQCcyqkMOOQSAK664gvnz5xd0262trTQ2NhZ0m6VE/ZPpRoORSlHMmzePLVu2FHSb7e3t\ndHR0FHSbpUT9S868efOSbsK0pAAkRXHGGWcUfJvlPtyJ+ifTjR4gFRGRRCgAiYhIIhSAREQkEQpA\nIiKSCAUgERFJhAKQiIgkQgFIREQSoQAkIiKJUAASEZFEKACJiEgiSm4oHjObAfwYmA+0AO939605\n6r0GuNXdD4mW5wI/AWYAzcCZ7t45aQ0XEZG8lOIR0DnAo+5+DPAj4PPZFczsvcBPgbmx7EuAFe5+\nLPBH4COT0FYREZmgUgxARwN3RfN3AifkqLMNOGaM9Y4vSutERKQgEj0FZ2b/BJwPpKOsFLAJ2Bkt\ntwJN2eu5+x3R+vHsxqz1ZhW+xSIiUiiJBiB3vxa4Np5nZjcTggnRdLwvwmuJ6neNtl5zc/OE2joV\ntLa2qn9TmPo3dZVz34qp5G5CAFYDpwBroumqUeqmcqx3A/DWkdYr5/eRlPv7VtS/qa2c+1fOfQPY\nuHFjUbZbigHoe8D1ZraKcDRzOoCZfQ240d3XxOqmY/OXRut9CHgls56IiJSmkgtA7t4BvCtH/kU5\n8hbE5rcQjnxERGQKKMW74EREZBpQABIRkUQoAImISCIUgEREJBEKQCIikggFIBERSYQCkIiIJEIB\nSEREElFyD6JON9/4BmzbBgcfDIsXgxnMnJl0q0REik8BKGFHHAErV8L//i98/evw9NMwd24IRAce\nGNJrXhPS/vtDbW3SLRYRKQwFoIQddVRIGX19sGEDuMOTT8JTT8EvfxmmL74I8+fDAQeEYPRXfxXS\nokUhpVIjfYqISOlRACoxlZUhuOy/P5x88tCy3l544QV45hl49ll47rlw9PTcc/D887Bp0z7Mmwf7\n7RfSvvuG9KpXwYIFg9O6uiR6JiIylALQFFJVNRicctmwYSOwgBdeCEdLmfTQQ9DcDC+9BJs2wYwZ\nsM8+Ie2992Daa69whJWZzpunU34iUjwKQGWkqioc4SxcOHKddBq2bw8BadOmoemxx2DLFti8OUxf\neSUcLc2bB3vuOTjNpLlzh6Y99ghJQUtExkMBaJpJpQYDxZIlo9dNp2HHDnj55ZBeeWUwbd0arlNt\n3Rru4tu6NQS2rVtDINxjD5gzZzDNnj04nT0bZs0aOp9JTU1QUzM5PwsRSZYCkIwolRoMIAceOL51\n0mloawtBaceOEJQyKbP87LNhfufOwWkmtbSEANbUNJgaG8O0qmo2e+0VblNvbBxMmeWZM4emhoYw\nrdJvuUhJ0p+mFFQqNRgARjsVOJJ0Gjo6oLU1BKTW1hCUWlpgw4YuqqrqaW0N+c8/z8D8rl0h8LW2\nhumuXYPzVVWDAam+Pkzjqb5+MD8zn0l1dYPT7Pl4mjFDdyGK5EsBSEpKKjX45b/XXkPLmps7WLBg\nTl7bS6ehqysEokzatQva20Nqaxs+bWkJ18Ta2kIw7OgIZdnTzs7B8q6ucO0rHpCyp/FUWzt8ubt7\n5sCNH5nyzHyuVFMzfL6mRoFQpg4FIClrqdTgl/zcucX7nP7+EIQyASkTlDo7BwNVZj6z3NU1WKej\nA7ZuraClZWh+Zj6e190d5jPT+HxPD1RXDwaj7DRSfnV17uXq6vHPj5SqqsJ0+/YqOjqG5mXXqdDg\nYNOKApBIAVRUDB79TFRzcwsLFuzeOEzp9NBg1N09NFDF8+JlPT25y+L57e2D85n8XMvZqbc3TDs6\n9iCdHlzOVSeVGgxG8Wn2fK7lkfIyqbIy9/xIZSPl5ZrfsaOW+fOH5mfXy7WcTyrHI1sFIJEykkoN\nnpYrNc3NW1iwYMGodfr6QjDKBKT4NDt/pOXs1NMzuN349uPzmeXOzqFluerEyzJp164Gqqtzl421\nPJ7U3x/2bUVF7uCUK7+xER58sLTHllQAEpGSkfnyLMUAOprm5m1jBtfdkU6HINTfP3KAys6rqSnt\n4AMKQCIiJS+VGgzO1dVJt6ZwdMlPREQSoQAkIiKJUAASEZFEKACJiEgiFIBERCQRCkAiIpKIkrsN\n28xmAD8G5gMtwPvdfWuOeq8BbnX3Q6LlOcCTwGNRlVvd/arJabWIiOSr5AIQcA7wqLt/0czeDXwe\nOC9ewczeC3wSiI/u9QbgJ+7+yUlrqYiITFgpnoI7Grgrmr8TOCFHnW3AMVl5y4BlZvYrM/svM9u7\niG0UEZHdlOgRkJn9E3A+kI6yUsAmYGe03Ao0Za/n7ndE68eznwDWuPt9ZnY6cBXwD8VpuYiI7K5E\nA5C7XwtcG88zs5uBxmixEdgxzs3dD7RH87cC/5qr0tq1a/Nv6BSycePGpJtQVOrf1FbO/SvnvhVL\nKV4DWg2cAqyJpqtGqRsfoPwa4GbgRsJpu2GRZtmyZWU4oLmIyNRUigHoe8D1ZrYK6AJOBzCzrwE3\nuvuaWN10bP4i4D/N7BygDThrktorIiITkEqn02PXEhERKbBSPAIak5mlgO8ChwKdwFnu/kys/EPA\n2UAPcKm7/8LM5gI/AWYAzcCZ7t45Qt3EnikqZN+i+vMIpzWXuHv3eJ+zmqr9i/JeJOw/gN+6++cm\npXMU/HfzfODdhCP9O9z938pp/2X17053/2K0jXLZfx8D3g/0A/8W1U1s/xW7b9E28tp3pXgb9nic\nBtS6+5HAZ4HLMwVmthdwLnAEcDLwFTOrBi4BVrj7scAfgQ+PUjfzTNFxUZrMB1oL0bePRPWXA78k\n/LJnZJ6zOgb4EeE5q8lU1P6Z2auBtbF9N2lfXpFC/W7uD/yjux8OHAmcZGZLKI/9l6t/y81sSRnt\nv7mE39PDCdekvxdtJsn9V9S+TWTfTdUANPCskLs/DLwxVvZm4AF373X3FuApQsTPfr7oxBHqLiXZ\nZ4oK0bfjo/m+aH5bru0z8nNWxVTs/i0D9jWz+8zsdjM7sGg9ya0Q/TsB2ED4IsDd04SzFZ0j1J1M\nxepfNaF/ZbH/oqOaQ929H9gH2J69fSZ//xW7b3nvu6kagJoYfFYIoNfMKkYoawVmEW7p3jlKHsCu\nKP8J4BJ3fwvw34RniiZLofqGu9/r7tsZerdgE2M8Z1Vkxe5fM/Bldz8O+ArhdMdkKkj/3L3P3bcB\nmNm/A4+4+9OUyf4bpX8bKYP9B+Du/dGpqt8CN+XYxmTvv2L07UEG+5b3vpuqAaiFwWeFACqiaJwp\ni+/UJkKEjq/TGMuL1808d3Q/8Kso71bg9QVs+1gK0bfsZ6fid5qMVbfYit2/tcD/ALj7amBBYZo9\nbgXrn5nVmtkKoAH4WI7tT+n9l9W/j0blayiT/Qfg7t8B9gaONbO3EL7Mk9p/xejbPoS+HcsE9t1U\nDUCZZ4Uws8MZvFkA4HfA0WZWY2azgIOA9dE6b4vqvJXwfNHvR6h7DfDOqG7OZ4qKqFB9i4sfIQxs\nn7GfsyqGYvfvC0RjB5rZoYRTPZOpkP37H+CP7v7R6DTVkO0z9fdfrv6Vxf4zswMtPFQP4VRxZzSN\n153s/VfMvvUzgX03JW/Djt3NsTTKOpPwQ3rK3W83sw8CHyZ8MV3q7reZ2XzgemAm8Apwurt3jFD3\nrxgcoaGNcLfI5qnWt9g2nwEO8nAXXF1Udx+i56zcfctk9C1qS7H7N5tw6D+TcDfPx9w9c1dO0RWq\nf8BJhLuPHorqpgkXjh+lDPYfI/fvz8AKpvj+i75bPk/4wu8n3OX3pST//iahb3n/7U3JACQiIlPf\nVD0FJyIiU5wCkIiIJEIBSEREEqEAJCIiiVAAEhGRRCgAiYhIIhSAREQkEQpAIiKSCAUgkSIzsyum\n8+eLjGRKvpBOZCqI3qdyDmG4k/Om2+eLjEVD8YgUmZndFw1RP1qdzKjE6919abT8dne/I1bnQMLg\nleuisu6JfL6ZrQKOIoy/1uju7fn1SKQwdApOpHScARybq8DM9gPuBh4HTh1v8BnB3zE42rtIYnQK\nTmQEZnYn4SVcjxP+Wfs88DPCMPV7AmcBi939pQJ95M7oBXvZ7ZgPrAReJBz5dO7Oh7j7DjPbNnZN\nkeJSABLJwcwWAne4+1XR8mHAxcA57r4zyqsDjoveKpk5l50Cutz9ZwVqRxPwS8KLwN4aP11mZjOA\nd2WtUtDPFykmBSCR3E4CfhBbXg6syQSfyHp3/9E4tpUau0pODcAdwCHAge7eGi+MjoRuKOLnixSV\nrgGJ5ODuV2ddZ1lOOBKJ1/nOWNsxs48BrzGzz5rZ3nk241vAbGAzcGme6xbi80WKSkdAImMws0bg\nMODT+a4bBakxA9UIWgivhD8MuNXMbnP3/5rEzxcpKh0BiYzteMKr2R+e5M+90N03uft/E25++I6Z\n7TPJbRApGgUgkbGdCNzr7v1j1iysvtj8udHytZPcBpGiUQASGduw6z+Tzd23Ah8DTjKzc5Jsi0ih\n6BqQSA7Rg58fAF4FHEC43Xo+cLm7d0xCE4YNUeLuN5nZzcDXzewed396EtohUjQKQCI5uPsLwL9F\nix9J4PMrR8j/h8lui0ix6BScSOmYbWZ7FPtDzGwOMLfYnyMyFh0BiZSGNPAjwrA/S4v8Wf8DHEmO\n03wik0mjYYuISCJ0Ck5ERBKhACQiIolQABIRkUQoAImISCIUgEREJBEKQCIikggFIBERSYQCkIiI\nJEIBSEREEvH/AwS9lS+dZiSEAAAAAElFTkSuQmCC\n", "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "V = np.array(Vlist)\n", "fig, ax1 = plt.subplots()\n", "ax1.plot(1./Trange, V[:,0], 'k', label='$V_{\\\\rm{total}}$')\n", "ax1.plot(1./Trange, V[:,1], 'b', label='$V_{11}$')\n", "ax1.plot(1./Trange, V[:,2], 'r', label='$V_{12}$')\n", "ax1.plot(1./Trange, 2*V[:,3], 'g', label='$V_{44}$')\n", "ax1.set_yscale('linear')\n", "ax1.set_ylabel('$V$ [atomic volume]', fontsize='x-large')\n", "ax1.set_xlabel('$T^{-1}$ [K$^{-1}$]', fontsize='x-large')\n", "ax1.legend(bbox_to_anchor=(0.3,0.3,0.5,0.2), ncol=2, \n", " shadow=True, frameon=True, fontsize='x-large')\n", "ax2 = ax1.twiny()\n", "ax2.set_xlim(ax1.get_xlim())\n", "ax2.set_xticks([1./t for t in Tlabels])\n", "ax2.set_xticklabels([\"{:.0f}K\".format(t) for t in Tlabels])\n", "ax2.set_xlabel('$T$ [K]', fontsize='x-large')\n", "ax2.grid(False)\n", "ax2.tick_params(axis='x', top='on', direction='in', length=6)\n", "plt.show()\n", "# plt.savefig('Fe-C-activation-volume.pdf', transparent=True, format='pdf')" ] }, { "cell_type": "code", "execution_count": 11, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Total volume: 0.0921, 1.0722A^3\n", "V_xxxx: -0.1175, -1.3681A^3\n", "V_xxyy: 0.1048, 1.2202A^3\n", "V_xyxy: 0.0061, 0.0714A^3\n" ] } ], "source": [ "print('Total volume: {v[0]:.4f}, {V[0]:.4f}A^3\\nV_xxxx: {v[1]:.4f}, {V[1]:.4f}A^3\\nV_xxyy: {v[2]:.4f}, {V[2]:.4f}A^3\\nV_xyxy: {v[3]:.4f}, {V[3]:.4f}A^3'.format(v=V[-1,:], V=V[-1,:]*1e3*Fe.volume))" ] }, { "cell_type": "code", "execution_count": 12, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Spherical average uniaxial activation volume: -0.0237 -0.2757A^3\n" ] } ], "source": [ "Vsph = 0.2*(3*V[-1,1] + 2*V[-1,2] + 4*V[-1,3]) # (3V11 + 2V12 + 2V44)/5\n", "print('Spherical average uniaxial activation volume: {:.4f} {:.4f}A^3'.format(Vsph, Vsph*1e3*Fe.volume))" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.4.5" } }, "nbformat": 4, "nbformat_minor": 0 }