Skip to content

Instantly share code, notes, and snippets.

@ketch
Last active February 21, 2017 21:31
Show Gist options
  • Save ketch/ad1b6459fb3debad86a3760edcbd81de to your computer and use it in GitHub Desktop.
Save ketch/ad1b6459fb3debad86a3760edcbd81de to your computer and use it in GitHub Desktop.
Comparison of some Riemann solvers for the variable speed-limit LWR traffic model
all: traffic_fwave_old.so traffic_fwave.so traffic_qwave.so
traffic_fwave_old.so: rp1_traffic_vc_tracer_fwave_old.f90
f2py -c rp1_traffic_vc_tracer_fwave_old.f90 -m traffic_fwave_old
traffic_fwave.so: rp1_traffic_vc_tracer_fwave.f90
f2py -c rp1_traffic_vc_tracer_fwave.f90 -m traffic_fwave
traffic_qwave.so: rp1_traffic_vc_tracer_qwave.f90
f2py -c rp1_traffic_vc_tracer_qwave.f90 -m traffic_qwave
clean:
rm *.so
! Riemann solver for the variable speed limit LWR traffic model
! with a tracer:
! q_t + (u_max q(1-q))_x = 0
! p_t + u_max(1-q) p_x = 0
! Here the speed limit u_max (stored in aux(1,i) may vary with x and t.
! waves: 2
! equations: 2
! aux fields: 1
! Conserved quantities:
! 1 q
! 2 p
! Auxiliary variables:
! 1 u_max
! See http://www.clawpack.org/riemann.html for a detailed explanation
! of the Riemann solver API.
subroutine rp1(maxm,num_eqn,num_waves,num_aux,num_ghost,num_cells, &
ql,qr,auxl,auxr,fwave,s,amdq,apdq)
implicit none
! Inputs
integer, intent(in) :: maxm, num_eqn, num_waves, num_aux, num_ghost, num_cells
double precision, intent(in), dimension(num_eqn, 1-num_ghost:maxm+num_ghost) :: ql,qr
double precision, intent(in), dimension(num_aux, 1-num_ghost:maxm+num_ghost) :: auxl,auxr
! Outputs
double precision, intent(out) :: s(num_waves,1-num_ghost:num_cells+num_ghost)
double precision, intent(out) :: fwave(num_eqn,num_waves,1-num_ghost:num_cells+num_ghost)
double precision, intent(out), dimension(num_eqn, 1-num_ghost:maxm+num_ghost) :: amdq,apdq
! Locals
double precision :: fim1, fi, sim1, si, f0
integer :: i
do i = 2-num_ghost, num_cells+num_ghost
! Compute the wave, speeds, and flux difference
! compute characteristic speed in each cell:
sim1 = auxr(1,i-1)*(1.d0 - 2.d0*qr(1,i-1))
si = auxl(1,i )*(1.d0 - 2.d0*ql(1,i ))
s(2,i) = auxl(1,i) * (1.d0 - ql(1,i))
! compute flux in each cell and flux difference:
fim1 = auxr(1,i-1)*qr(1,i-1)*(1.d0 - qr(1,i-1))
fi = auxl(1,i )*ql(1,i )*(1.d0 - ql(1,i ))
fwave(1,1,i) = fi - fim1
fwave(2,1,i) = 0.d0
fwave(1,2,i) = 0.d0
fwave(2,2,i) = s(2,i)*(ql(2,i) - qr(2,i-1))
apdq(2,i) = fwave(2,2,i)
amdq(2,i) = 0.d0 ! Traffic always moves right
if (sim1 .lt. 0.d0 .and. si .le. 0.d0) then
! left-going
s(1,i) = sim1
amdq(1,i) = fwave(1,1,i)
apdq(1,i) = 0.d0
else if (sim1 .ge. 0.d0 .and. si .gt. 0.d0) then
! right-going
s(1,i) = si
amdq(1,i) = 0.d0
apdq(1,i) = fwave(1,1,i)
else if (sim1 .lt. 0.d0 .and. si .gt. 0.d0) then
! transonic rarefaction
s(1,i) = 0.5d0*(sim1 + si)
! entropy fix: (perhaps doesn't work for all cases!!!)
! This assumes the flux in the transonic case should
! correspond to q=0.5 on the side with the smaller umax value.
f0 = dmin1(auxr(1,i-1),auxl(1,i))*0.25d0
! split fwave between amdq and apdq:
amdq(1,i) = f0 - fim1
apdq(1,i) = fi - f0
else
! transonic shock
s(1,i) = 0.5d0*(sim1 + si)
if (fi-fim1 .lt. 0.d0) then
amdq(1,i) = fwave(1,1,i)
apdq(1,i) = 0.d0
else if (fi-fim1 .gt. 0.d0) then
amdq(1,i) = 0.d0
apdq(1,i) = fwave(1,1,i)
else
amdq(1,i) = 0.5d0 * fwave(1,1,i)
apdq(1,i) = 0.5d0 * fwave(1,1,i)
endif
endif
enddo
return
end
! Riemann solver for the variable speed limit LWR traffic model
! with a tracer:
! q_t + (u_max q(1-q))_x = 0
! p_t + u_max(1-q) p_x = 0
! Here the speed limit u_max (stored in aux(1,i) may vary with x and t.
! waves: 2
! equations: 2
! aux fields: 1
! Conserved quantities:
! 1 q
! 2 p
! Auxiliary variables:
! 1 u_max
! See http://www.clawpack.org/riemann.html for a detailed explanation
! of the Riemann solver API.
subroutine rp1(maxm,num_eqn,num_waves,num_aux,num_ghost,num_cells, &
ql,qr,auxl,auxr,fwave,s,amdq,apdq)
implicit none
! Inputs
integer, intent(in) :: maxm, num_eqn, num_waves, num_aux, num_ghost, num_cells
double precision, intent(in), dimension(num_eqn, 1-num_ghost:maxm+num_ghost) :: ql,qr
double precision, intent(in), dimension(num_aux, 1-num_ghost:maxm+num_ghost) :: auxl,auxr
! Outputs
double precision, intent(out) :: s(num_waves,1-num_ghost:num_cells+num_ghost)
double precision, intent(out) :: fwave(num_eqn,num_waves,1-num_ghost:num_cells+num_ghost)
double precision, intent(out), dimension(num_eqn, 1-num_ghost:maxm+num_ghost) :: amdq,apdq
! Locals
double precision :: fim1, fi, sim1, si, f0
integer :: i
do i = 2-num_ghost, num_cells+num_ghost
! Compute the wave, speeds, and flux difference
! compute characteristic speed in each cell:
sim1 = auxr(1,i-1)*(1.d0 - 2.d0*qr(1,i-1))
si = auxl(1,i )*(1.d0 - 2.d0*ql(1,i ))
s(2,i) = auxl(1,i) * (1.d0 - ql(1,i))
! compute flux in each cell and flux difference:
fim1 = auxr(1,i-1)*qr(1,i-1)*(1.d0 - qr(1,i-1))
fi = auxl(1,i )*ql(1,i )*(1.d0 - ql(1,i ))
fwave(1,1,i) = fi - fim1
fwave(2,1,i) = 0.d0
fwave(1,2,i) = 0.d0
fwave(2,2,i) = s(2,i)*(ql(2,i) - qr(2,i-1))
apdq(2,i) = fwave(2,2,i)
amdq(2,i) = 0.d0 ! Traffic always moves right
if (sim1 .lt. 0.d0 .and. si .le. 0.d0) then
! left-going
s(1,i) = sim1
amdq(1,i) = fwave(1,1,i)
apdq(1,i) = 0.d0
else if (sim1 .ge. 0.d0 .and. si .gt. 0.d0) then
! right-going
s(1,i) = si
amdq(1,i) = 0.d0
apdq(1,i) = fwave(1,1,i)
else if (sim1 .lt. 0.d0 .and. si .gt. 0.d0) then
! transonic rarefaction
s(1,i) = 0.5d0*(sim1 + si)
! entropy fix: (perhaps doesn't work for all cases!!!)
! This assumes the flux in the transonic case should
! correspond to q=0.5 on the side with the smaller umax value.
f0 = dmin1(auxr(1,i-1),auxl(1,i))*0.25d0
! split fwave between amdq and apdq:
amdq(1,i) = f0 - fim1
apdq(1,i) = fi - f0
else
! transonic shock
s(1,i) = 0.5d0*(sim1 + si)
if (s(1,i) .lt. 0.d0) then
amdq(1,i) = fwave(1,1,i)
apdq(1,i) = 0.d0
else if (s(1,i) .gt. 0.d0) then
amdq(1,i) = 0.d0
apdq(1,i) = fwave(1,1,i)
else
amdq(1,i) = 0.5d0 * fwave(1,1,i)
apdq(1,i) = 0.5d0 * fwave(1,1,i)
endif
endif
enddo
return
end
! Riemann solver for the variable speed limit LWR traffic model
! with a tracer:
! q_t + (u_max q(1-q))_x = 0
! p_t + u_max(1-q) p_x = 0
! Here the speed limit u_max (stored in aux(1,i) may vary with x and t.
! waves: 2
! equations: 2
! aux fields: 1
! Conserved quantities:
! 1 q
! 2 p
! Auxiliary variables:
! 1 u_max
! See http://www.clawpack.org/riemann.html for a detailed explanation
! of the Riemann solver API.
subroutine rp1(maxm,num_eqn,num_waves,num_aux,num_ghost,num_cells, &
ql,qr,auxl,auxr,wave,s,amdq,apdq)
implicit none
! Inputs
integer, intent(in) :: maxm, num_eqn, num_waves, num_aux, num_ghost, num_cells
double precision, intent(in), dimension(num_eqn, 1-num_ghost:maxm+num_ghost) :: ql,qr
double precision, intent(in), dimension(num_aux, 1-num_ghost:maxm+num_ghost) :: auxl,auxr
! Outputs
double precision, intent(out) :: s(num_waves,1-num_ghost:num_cells+num_ghost)
double precision, intent(out) :: wave(num_eqn,num_waves,1-num_ghost:num_cells+num_ghost)
double precision, intent(out), dimension(num_eqn, 1-num_ghost:maxm+num_ghost) :: amdq,apdq
! Locals
double precision :: f_l, f_r, s_l, s_r, f0, q_l, q_r, v_l, v_r
integer :: i
do i = 2-num_ghost, num_cells+num_ghost
q_r = ql(1,i)
q_l = qr(1,i-1)
v_r = auxl(1,i)
v_l = auxr(1,i-1)
! compute characteristic speed in each cell:
s_l = v_l*(1.d0 - 2.d0*q_l)
s_r = v_r*(1.d0 - 2.d0*q_r)
s(2,i) = v_r * (1.d0 - q_r)
! compute flux in each cell and flux difference:
f_l = v_l*q_l*(1.d0 - q_l)
f_r = v_r*q_r*(1.d0 - q_r)
wave(1,1,i) = q_r - q_l ! Trying this even though it's not right
wave(2,1,i) = 0.d0
wave(1,2,i) = 0.d0
wave(2,2,i) = s(2,i)*(ql(2,i) - qr(2,i-1))
if (q_r .ne. q_l) then
s(1,i) = (f_r-f_l)/(q_r - q_l)
else
s(1,i) = 0.d0
endif
apdq(2,i) = wave(2,2,i)
amdq(2,i) = 0.d0 ! Traffic always moves right
if ((f_l .ge. 0.25d0*v_r) .and. (s_r .gt. 0.d0)) then
! left-going shock, right-going rarefaction
f0 = 0.25d0*v_r
elseif ((f_r .ge. 0.25d0*v_l) .and. (s_l .lt. 0.d0)) then
! right-going shock, left-going rarefaction
f0 = 0.25d0*v_l
elseif ((s_r .le. 0.d0) .and. (f_l .gt. f_r)) then
! left-going shock
f0 = f_r
elseif ((s_l .ge. 0.d0) .and. (f_r .gt. f_l)) then
! right-going shock
f0 = f_l
elseif ((s_l .le. 0.d0) .and. (s_r .ge. 0.d0)) then
! Transonic rarefaction
if (v_r .le. v_l) then
f0 = 0.25d0*v_r
else
f0 = 0.25d0*v_l
endif
elseif (f_l .le. f_r) then
! left-going rarefaction
f0 = f_r
else
! right-going rarefaction
f0 = f_l
endif
amdq(1,i) = f0 - f_l
apdq(1,i) = f_r - f0
enddo
return
end
Display the source blob
Display the rendered blob
Raw
{
"cells": [
{
"cell_type": "code",
"execution_count": 1,
"metadata": {
"collapsed": false,
"deletable": true,
"editable": true
},
"outputs": [],
"source": [
"%matplotlib inline\n",
"from clawpack import pyclaw\n",
"from clawpack import riemann\n",
"import matplotlib.pyplot as plt\n",
"import numpy as np\n",
"from clawpack.visclaw import ianimate\n",
"from utils import riemann_tools\n",
"\n",
"\n",
"import traffic_fwave_old\n",
"import traffic_fwave\n",
"import traffic_qwave"
]
},
{
"cell_type": "code",
"execution_count": 2,
"metadata": {
"collapsed": false,
"deletable": true,
"editable": true
},
"outputs": [],
"source": [
"def test_solver(q_l, q_r, v_l, v_r, riemann_solver):\n",
" solver = pyclaw.ClawSolver1D(riemann_solver)\n",
" solver.num_waves = 2\n",
" solver.num_eqn = 2\n",
"\n",
" solver.bc_lower[0] = pyclaw.BC.extrap\n",
" solver.bc_upper[0] = pyclaw.BC.extrap\n",
" solver.aux_bc_lower[0] = pyclaw.BC.extrap\n",
" solver.aux_bc_upper[0] = pyclaw.BC.extrap\n",
"\n",
" x = pyclaw.Dimension(-1.0,1.0,50,name='x')\n",
" domain = pyclaw.Domain(x)\n",
" num_eqn = 2\n",
" num_aux = 1\n",
" state = pyclaw.State(domain,num_eqn,num_aux)\n",
"\n",
" grid = state.grid\n",
" xc=grid.p_centers[0]\n",
"\n",
" state.q[0,:] = q_l*(xc<0) + q_r*(xc>=0.)\n",
" state.q[1,:] = 0.0\n",
" state.aux[0,:] = v_l*(xc<0) + v_r*(xc>=0.) # Speed limit\n",
"\n",
" claw = pyclaw.Controller()\n",
" claw.tfinal = 1.0\n",
" claw.solution = pyclaw.Solution(state,domain)\n",
" claw.solver = solver\n",
" claw.keep_copy = True\n",
" claw.num_output_times = 10\n",
" claw.verbosity = 0\n",
"\n",
" claw.run()\n",
" return xc, claw.frames[5].q[0,:]"
]
},
{
"cell_type": "code",
"execution_count": 3,
"metadata": {
"collapsed": true,
"deletable": true,
"editable": true
},
"outputs": [],
"source": [
"def riemann_traffic_vc_exact(q_l,q_r,v_l,v_r):\n",
" r\"\"\"Exact solution to the Riemann problem for the LWR traffic model with variable speed limit.\n",
" Inputs:\n",
" - q_l, q_r : traffic density for left and right states\n",
" - v_l, v_r : speed limits for left and right states\n",
" \"\"\"\n",
" f = lambda q, v: v*q*(1-q)\n",
" c = lambda q, v: v*(1.-2.*q)\n",
" \n",
" f_l = f(q_l,v_l)\n",
" f_r = f(q_r,v_r)\n",
" c_l = c(q_l,v_l)\n",
" c_r = c(q_r,v_r)\n",
" \n",
" if (f_r <= f_l) and (q_r > 0.5):\n",
" # Left-going shock (1)\n",
" q_star = (1. + np.sqrt(1.-4*f_r/v_l))/2.\n",
" states = [q_l, q_star, q_r]\n",
" shock_speed = (f_r-f_l)/(q_star-q_l)\n",
" speeds = [shock_speed,0]\n",
" def reval(xi):\n",
" q = np.zeros((1,len(xi)))\n",
" q[0,:] = (xi<=shock_speed)*q_l \\\n",
" + (shock_speed<xi)*(xi<=0)*q_star \\\n",
" + (xi>=0)*q_r\n",
" return q\n",
" \n",
" elif (f_r >= f_l) and (q_l < 0.5):\n",
" # Right-going shock (2)\n",
" q_star = (1. - np.sqrt(1.-4*f_l/v_r))/2.\n",
" states = [q_l, q_star, q_r]\n",
" shock_speed = (f_r-f_l)/(q_r-q_star)\n",
" speeds = [0,shock_speed]\n",
" def reval(xi):\n",
" q = np.zeros((1,len(xi)))\n",
" q[0,:] = (xi<=0)*q_l \\\n",
" + (0<xi)*(xi<=shock_speed)*q_star \\\n",
" + (xi>=shock_speed)*q_r\n",
" return q\n",
" \n",
" elif (f_l >= v_r/4.) and (q_r <= 0.5):\n",
" # Left-going shock and right-going rarefaction (3)\n",
" q_star = (1. + np.sqrt(1.-v_r/v_l))/2.\n",
" shock_speed = -(f_l-v_r/4.)/(q_star-q_l)\n",
" states = [q_l, q_star, q_r]\n",
" speeds = [shock_speed,(0,c_r)]\n",
" def reval(xi):\n",
" q = np.zeros((1,len(xi)))\n",
" q[0,:] = (xi<=shock_speed)*q_l \\\n",
" + (shock_speed<xi)*(xi<=0)*q_star \\\n",
" + (0<xi)*(xi<c_r)*(1.-xi/v_r)/2. \\\n",
" + (xi>=c_r)*q_r\n",
" return q\n",
" \n",
" elif (f_r >= v_l/4.) and (q_l >= 0.5):\n",
" # Left-going rarefaction and right-going shock (4)\n",
" q_star = (1. - np.sqrt(1.-v_l/v_r))/2.\n",
" shock_speed = (f_r-v_l/4.)/(q_r-q_star)\n",
" states = [q_l,q_star,q_r]\n",
" speeds = [(c_l, 0), shock_speed]\n",
" def reval(xi):\n",
" q = np.zeros((1,len(xi)))\n",
" q[0,:] = (xi<=c_l)*q_l \\\n",
" + (c_l<xi)*(xi<=0)*(1.-xi/v_l)/2. \\\n",
" + (0<xi)*(xi<=shock_speed)*q_star \\\n",
" + (shock_speed<xi)*q_r\n",
" return q\n",
" \n",
" elif (f_l<f_r<=v_l/4.) and (q_l>0.5) and (q_r>=0.5):\n",
" # Left-going rarefaction (6)\n",
" q_star = (1. + np.sqrt(1.-4*f_r/v_l))/2.\n",
" states = [q_l, q_star, q_r]\n",
" c_star = c(q_star,v_l)\n",
" speeds = [(c_l,c_star),0]\n",
" def reval(xi):\n",
" q = np.zeros((1,len(xi)))\n",
" q[0,:] = (xi<=c_l)*q_l \\\n",
" + (c_l<xi)*(xi<=c_star)*(1.-xi/v_l)/2. \\\n",
" + (c_star<xi)*(xi<=0)*q_star \\\n",
" + (0<=xi)*q_r\n",
" return q\n",
" \n",
" elif (f_r<f_l<=v_r/4.) and (q_l<=0.5) and (q_r<0.5):\n",
" # Right-going rarefaction (7)\n",
" q_star = (1. - np.sqrt(1.-4*f_l/v_r))/2.\n",
" states = [q_l, q_star, q_r]\n",
" c_star = c(q_star,v_r)\n",
" speeds = [0,(c_star,c_r)] \n",
" def reval(xi):\n",
" q = np.zeros((1,len(xi)))\n",
" q[0,:] = (xi<=0)*q_l \\\n",
" + (0<xi)*(xi<=c_star)*q_star \\\n",
" + (c_star<xi)*(xi<=c_r)*(1.-xi/v_r)/2. \\\n",
" + (c_r<=xi)*q_r\n",
" return q\n",
" \n",
" elif (q_l>=0.5) and (q_r<=0.5) and (f_l<=v_r/4.) and (f_r<=v_l/4.):\n",
" # Transonic rarefaction (5)\n",
" if v_r < v_l:\n",
" # q* on left side (5a)\n",
" q_star = (1. + np.sqrt(1.-v_r/v_l))/2.\n",
" states = [q_l, q_star, q_r]\n",
" c_star = c(q_star, v_l)\n",
" speeds = [(c_l,c_star),(0,c_r)]\n",
" def reval(xi):\n",
" q = np.zeros((1,len(xi)))\n",
" q[0,:] = (xi<=c_l)*q_l \\\n",
" + (c_l<xi)*(xi<=c_star)*(1.-xi/v_l)/2. \\\n",
" + (c_star<xi)*(xi<=0)*q_star \\\n",
" + (0<xi)*(xi<=c_r)*(1.-xi/v_r)/2. \\\n",
" + (c_r<xi)*q_r\n",
" return q\n",
" \n",
" elif v_r >= v_l:\n",
" # q* on right side (5b)\n",
" q_star = (1. - np.sqrt(1.-v_l/v_r))/2.\n",
" states = [q_l, q_star, q_r]\n",
" c_star = c(q_star, v_r)\n",
" speeds = [(c_l,0),(c_star,c_r)]\n",
" def reval(xi):\n",
" q = np.zeros((1,len(xi)))\n",
" q[0,:] = (xi<=c_l)*q_l \\\n",
" + (c_l<xi)*(xi<=0)*(1.-xi/v_l)/2. \\\n",
" + (0<xi)*(xi<=c_star)*q_star \\\n",
" + (c_star<xi)*(xi<=c_r)*(1.-xi/v_r)/2. \\\n",
" + (c_r<xi)*q_r\n",
" return q\n",
" \n",
" else: # v_r == v_l\n",
" states = [q_l, q_r]\n",
" speeds = [(c_l,c_r),]\n",
" def reval(xi):\n",
" q = np.zeros((1,len(xi)))\n",
" q[0,:] = (xi<=c_l)*q_l \\\n",
" + (c_l<xi)*(xi<=c_r)*(1.-xi/v_l)/2. \\\n",
" + (c_r<xi)*q_r\n",
"\n",
" else:\n",
" print f_l, f_r\n",
" raise Exception('Unhandled state!')\n",
" \n",
" return np.array([states]), speeds, reval"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"# Old method generates unphysical values"
]
},
{
"cell_type": "code",
"execution_count": 4,
"metadata": {
"collapsed": false,
"deletable": true,
"editable": true
},
"outputs": [
{
"data": {
"text/plain": [
"<matplotlib.legend.Legend at 0x1141f0b10>"
]
},
"execution_count": 4,
"metadata": {},
"output_type": "execute_result"
},
{
"data": {
"image/png": "iVBORw0KGgoAAAANSUhEUgAAAiEAAAEICAYAAACJXCTkAAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAIABJREFUeJzsnXmcU+X1/99nssxkVgZQqiwuCCgqAoILClJRq9attSrY\nfq22rl+t1qXWVouI4tJSqRZ/tfi1xS4qaLXVinWnUpcWVEoLlEWqMKDADLNPZpJJnt8f92YmM5Nk\nEmYmN5Oc9+uV1yTPfXLvSeDefO455zlHjDEoiqIoiqKkmzynDVAURVEUJTdREaIoiqIoiiOoCFEU\nRVEUxRFUhCiKoiiK4ggqQhRFURRFcQQVIYqiKIqiOIKKECVrEJGXReSbTtvRHSJyoIgYEXE7bYui\nKIqTqAhRYiIiPxCRlzuNbYozNjO91sXGGHOGMeaJvj6OiFwqIn/r6+MoSjYgIp+IyCkJtk8XkYo0\n2TJQRJ4XkUYR+VRELk4wd46IBEWkIepxcNT28SLygYg02X/HR20TEXlARKrsxwMiIn39+fojKkKU\neLwNTBERF4CI7Ad4gAmdxg6x5yqKomQ6jwABYAjwdeAXInJ4gvlLjDHFUY8tACLiBf4E/A4oB54A\n/mSPA1wJnAccBYwDzgau6osP1N9REaLEYyWW6Iio+6nAW8CGTmMfG2N2AIjIQyKyTUTq7DuDqfb4\n/iLiF5GBkZ2LyAQRqRQRj/36WyKyXkSqReQVETnAHhcRWSAiu+z9/ktEjohlsIgsF5HL7eeXisjf\nRGS+vc//isgZcd43UkT2iMjEKHt3i8j0GHMPAx4FjrfvjGoS2HKfiPzDtvtP0Z+/09zL7M9eLyJb\nROSqqG3TRaRCRG62v4PPROSyqO359mfcKiI7ReRREfHFOo6iJIt9fn5o/59cIiJPi8g9ceaOFJE3\n7Tv+ShH5vYgMsLf9FhgBvGifL7d2em8R8DKwf5S3Yf8++kxFwPnAj4wxDcaYvwEvAP+zF7ubDriB\nnxljWowxDwMCnGxv/ybwU2NMhTFmO/BT4NIefoSsREWIEhNjTAD4OzDNHpoGrAD+1mks2guyEkug\nDASeBJ4RkQJbpLyHdQGIcDHwrDEmKCLnAj8EvgrsYx/nKXveafZxRgNlwIVAVZIf41gs0TQY+DHw\neCyXqDHmY+D7wO9EpBD4NfCEMWZ5jLnrgauB9+w7owEJjn8J8C1gP6AVeDjOvF3AWUApcBmwICKI\nbL6A9dmHAt8GHhGRcnvb/VjfzXgsr9RQYHYCmxQlIfbd/B+B32Kdy8/Q8dzt8hbgPmB/4DBgODAH\nwBjzP8BW4Gz7fPlx9BuNMY3AGcCOKG/Djhg23SYiNfEeSX600UCrMWZj1Ng/gUSekLPtG5S1InJN\n1PjhwBrTse/Jmqh9HW7vO9nj5CwqQpRE/JV2wTEVSxys6DT218hkY8zvjDFVxphWY8xPgXxgjL35\nSWAWWN4NYKY9BtaP+n3GmPXGmFbgXmC87Q0JAiXAoYDYcz5L0v5PjTGPGWNCWO7S/bDcsF0wxjwG\nbMYSXvsBtyd5jET81hjzb/tC+yPgwkgoq9OxXzLGfGws/gq8ivXdRggCc40xQWPMMqABGGN/j1cC\nNxpj9hhj6rG+u4zI0VH6LcdheUF/Zv+fexbrBiMmxpjNxpjXbI/AbuBB4KTeNMgYc78xZkC8R5K7\nKQbqOo3VYl1fYrEUS1TtA1wBzBaRWVH7qk2wr87ba4FizQvpiooQJRFvAyfaYYR9jDGbgHexckUG\nAkcQ5QkRkVvssEKtfXdShuWFAPgDVghjPywRE8YSNAAHAA9F3dXswbq7GmqMeRNYiBXL3SUii0Sk\nNEn7P488McY02U+LE8x/zP5MPzfGtNifaWqUm3htkseNsC3q+adYF/bBnSeJyBki8r59x1UDnNlp\nXpUtziI02Z9jH6AQ+CDqu/uLPa4oe8v+wPZOd/mfxpssIkPscM12EanDypPo8v88A2jA8jZGUwrU\nx5psjFlnjNlhjAkZY94FHgK+luS+Om8vBRo6facKKkKUxLyHJSSuAN4BMMbUATvssR3GmP+C9WMN\n3IoVLim3705qscQExphqrDv8i7BCMU9HnZDbgKs63d347BMfY8zDxpijgbFYLtXv9fYHFZFi4GfA\n48CcSP6GMWZFlJs44k5N9kIyPOr5CCyPRmWn4+ZjCbT5wBD7e1uG/b11QyXgBw6P+t7KjDGJhJai\ndMdnwNBOd+0jEsy/F+ucONIYUwp8g47/f7s7X7o9n0Tkh9JxlUqHR3fvt9kIuEVkVNTYUUCyNxeG\n9s+1FhjX6TsaF7Wvtfa+9+Y4OYWKECUuxhg/sAq4iXavBVh5ITfRMR+kBCvvYTfWiT6brncKT2Ll\nSXyN9lAMWImePxA7S11EykTkAvv5ZBE5VqwE1kagGcuL0ts8BKwyxlwOvGTbFI+dwDBpz4SPxzdE\nZKydZzIXKwcm1GmOFytstRtoFSt59rRkDDbGhLG8NwtEZF8AERkqIl9K5v2KEof3sM7l60XEIyJf\nBY5JML8E686/VkSG0vUmYSdwcJd3ddw+SETK4k0wxtzbaZVKh0cyH8oOiz4HzBWRIhE5ATgXK/el\nCyJyroiUi8UxwPVYK2IAlgMhrO8oX0Sus8fftP/+BrjJPh/3B24GFidjZ66hIkTpjr8C+2IJjwgr\n7LFoEfIKVihgI5brtpmO4QiwMtFHAZ8bY9qStowxzwMPAE/b7tx/YyWrgSVkHgOq7f1WAT/pjQ8W\nwU6MPR2IJJ7dBEwUka/HecubWHc1n4tIZZw5YF3cFmOFhQqwLmIdsPM4rseKP1djeYleSMH872Pl\nsrxvf3ev056HoygpYyelfxVrNcceLO/lcwnechcwEcvz+VKMufcBd9ghw1tiHO8/WInoW+w5fbI6\nxuZ/AR9WMvhTwDXGmLXQHnqNmjsT69yqxxIVD0TqENnf0XlYN1U1WAno59njAL8EXgT+hXU9e8ke\nUzohGqJSlN5HRJYDvzPG/J/TtihKTxGRxUCFMeYOp21Rsgv1hCiKoiiK4ggqQhRFURRFcQQNxyiK\noiiK4gjqCVEURVEUxRG0lXgf88knnxi3W79mRemOQCBQefDBB2d0oTU9n3uHxpZWPqlqImwMh+1X\nSkswxNY9TYwYWIjHnYc/EKIpEMIfCOEPhgjbHvt8dx6hsGHEwEKK8rv/d/h0TxOB1jCj9tXSOekm\n2fNZz6Y+xu12M2zYMKfNUJSMZ8uWLXGrcmYKej73Do+8tZlf/rWCYw8eyGPHWjUAq6jkrYparj7p\nwA5zW0NhNu1q4KevbuD19bu4/uRDOHXygV13GoOfvrOalZ/s4W/fP7n7yUqvkuz5rOEYRVEUJa2M\nGVJCXXOQC45uF3RTRg7m6pNGdpnrduVR3RTgw601XH/yIfzu71t59+NE5Xna8Xld+AOd6wMqmYSK\nEEVRFCWtPPPBNgYXe/nioft2O/fdjyu57smPWHjxBG46bQwLL57AdU9+lJQQKfS4aFIRktGoCFEU\nRVHSRlVDC2+s38V544ficXX/E7SmopaFF09gykirJ96UkYNZePEE1lR0bmLblUKvy8opCesq0ExF\nc0IURVGUtPHH1TtoDRsumDS8+8kQM0QzZeTgNlGSCJ/X+olrbg1R6N37n7tgMEhFRQXNzc17vY9s\npaCggGHDhuHxePbq/SpCFEVRlLRgjOGZVds4algZY75Q0ufHK/S6AGgK9EyEVFRUUFJSwoEHHkjH\nxrm5jTGGqqoqKioqOOigg/ZqHxqOiUJEfiUiu0Tk33G2i4g8LCKbRWSNiExMt42KonSPnsuZyb+3\n1/Gfz+v5WpJekJ7is0VIT5NTm5ubGTRokAqQTogIgwYN6pGHSEVIRxZjdVONxxlYXWBHAVcCv0iD\nTYqipM5i9FzOOJ75YBtedx7njOvLRrntRHtCeooKkNj09HtRERKFMeZtrNbV8TgX+I2xeB8YICL7\nJdpnRQWEw71ppdJTAgHroWQvfXEuKz2jORjiT6t38KXDv0BZ4d7lD6RKuwhpTcvxlNRREZIaQ4Ft\nUa8r7LEOiMiVIrJKRFbt3AnPPps2+5Qk+OQT6/Hqq9Zj+XJLlLz6qsOGKekkqXMZOp7PlZXJ1adQ\nuvL6+p3U+jvWBulrInkg2VArxOVyMX78+LbH/fff32v7Xr16NcuWLeu1/aWCJqb2AcaYRcAigFde\nqTC33Qbnnw8ul8OGKR1oabH+BoPgdre/VpRoos/niooKXeu5lzyzqoL9ywo44ZDuV7X0FhFPSGMW\niBCfz8fq1av7ZN+rV69m1apVnHnmmX2y/0SoJyQ1tgPRGVXD7LG4DB0K//kPPPlkn9ql9ABjIE/P\nhFwj5XNZ2Xs+q/WzYtNuzj96GK689OVWZHs4pra2ljFjxrBhwwYAZs2axWOPPQbANddcw6RJkzj8\n8MO58847296zcuVKpkyZwlFHHcUxxxxDbW0ts2fPZsmSJYwfP54lS5ak9TOoJyQ1XgCuE5GngWOB\nWmPMZ4neMGAATJgAd90FM2fCXi6lVnoZ0+l+NhwGzTvLKVI+l5W957kPtxM28LU0hmKgvU5Ib4Zj\n+uo60fma1Bm/38/48ePbXv/gBz/goosuYuHChVx66aXccMMNVFdXc8UVVwAwb948Bg4cSCgUYsaM\nGaxZs4ZDDz2Uiy66iCVLljB58mTq6uooLCxk7ty5rFq1ioULF/bNh0uAipAoROQpYDowWEQqgDsB\nD4Ax5lFgGXAmsBloAi5LZr9z58LZZ8MTT8Dll/eF5UqqRLwfkaThQAC8XmdtUnqPvjqXldSJ1AY5\n5qCBHDCoqH3DmqXwxlyorYCyYTBjNoy7MPZOUpkbRaGn91bHOE28cMypp57KM888w7XXXss///nP\ntvGlS5eyaNEiWltb+eyzz1i3bh0iwn777cfkyZMBKC0tTZv98VAREoUxZlY32w1wbar7/fKX4dhj\n4e674X/+B/Lz99pEpZcIhzuKkJYW/XfJJvrqXM41pi+ZTlVzVZfxQQWDWH7R8qTmlnrK2V71fa79\n4iHtc39/HFWtjTBQYKAdFfvobgb960GWf/39jvtNYW5n2uqEBHtPhHTnsUg34XCY9evXU1hYSHV1\nNcOGDeO///0v8+fPZ+XKlZSXl3PppZdmbLVXjYSnARFLgGzdCo8/7rQ1CrSLkAgqQhSlK7FERbzx\neHPrgtWc4v03ZxesgfUvwr//YImKWPttbYRVv4YPFtuPJxLP7YZ8dx55kr05IQALFizgsMMO48kn\nn+Syyy4jGAxSV1dHUVERZWVl7Ny5k5dffhmAMWPG8Nlnn7Fy5UoA6uvraW1tpaSkhPr6ekfsV09I\nmjjlFJg6Fe65By67DHw+py3KbYzpGNtVEaLkCs3N8OijsHt3EpNHx9905bz3MES5BRKke3yr8EE+\n/HPUwH7xu+e++/ptHQcSzL39dhg50rqmxsrVEBEKve6sCMd0zgk5/fTTueyyy/i///s//vGPf1BS\nUsK0adO45557uOuuu5gwYQKHHnoow4cP54QTTgDA6/WyZMkSvvOd7+D3+/H5fLz++ut88Ytf5P77\n72f8+PFtuSbpQkVImhCxBMhJJ1kXgBtvdNqi3KbzihgVIUqu8NJLyV9/jlgcf9t7w65M+phXJRAS\nPZl7773W3+OPh8MOiz3H53VlRZ2QUCj2Z1i/fn3b8wcffLDt+eLFi2POnzx5Mu+/3zWMFfGOpBsV\nIWlk2jTLI3LffXDFFVBc7LRFuUuscExBgXP2KEq6aGiw/o4fD1/7WuK5TyfYNmPXE1GvhDf2vSTu\n3JN3PtGhvHeiuTN2/abD60Rz998fduxo/0yxKPS6ssITkq2oCEkzd99tqfaFC+G227qfr/QNsURI\nWZlz9ihKuhk3zgpnJOLpJ+Jv+9n3onr+1e/kyOfiz33o1o79AY9MuN8JSc+NiJBEZEs4JlvRxNQ0\nc9xx1mqZH/8YamudtiZ30ZwQJVdJZXXHoIJB3Y+HWuHZbzEoFLtJVnn+wL3bbwpzE30myxOSvYmp\n/R31hDjA3Llw9NHws59BVCE7JY1E54QYoyJEUWKx/KLlVNRXcMZzZzD7+NlcMPqCrpPeugc+/RvL\nv/JLOGomAKGwYeoDb3LIkBJ+M/OYmPtNxQaAO5Z9m79/9i6vnfprGHEcAJPnd//+Qq+L+mYVIZmK\nekIcYOJE+MpX4MEHYU+iPp9Kn6FLdJVcJ9nKnxurNwIwpnxM143/WQZ/WwBHX9omQADe/biSHbXN\nXDip9yqklhQNpj4vzypYZpPMZ/B5siMxNVtREeIQd90F9fUwPwklr/Q+0SJERCumKko8NlRvQBAO\nGXBIxw17/gvPXw37HQWnP8Cjf/2Ydz+2ugw/s6qCMp+H4nw3j/71416xo7RoPxrz8gjVbO2yrdtw\nTFA9IZmKihCHOPJIuOgiePjhJNfrK71K55wQ0N4xSm6QasXPTdWbGFE6gkJPYftgsBmWXgICXPgb\n8BQwblgZ1z35Ea+t3clf1n7OsQcN5Kal/2TcsN7J+C4utHJAGmq7ipBE+LzurPCEVFRUcO655zJq\n1ChGjhzJDTfcQCAQYPny5Zx11lkx33PggQdSWVmZZktTQ0WIg8yZA34/PPCA05bkHp1zQhRFic2G\nPRsYXd6patnLt8Lna+Ari6D8QACmjBzMwosncMOSjwi0hnl/SxULL57AlJGDe8WOEm8JAPV129rG\nkrlxSPcS3WiPUIR3P67skUfIGMNXv/pVzjvvPDZt2sTGjRtpaGjg9u6WN/UDVIQ4yJgxVi+ZRx7p\nfpmZ0rt0zglRlFwjmR/wpmAT2+q3dRQhq5+ED5+AE2+CMad3mD9l5GBG7mM1qfvmlAN7TYBAuwhp\nqO96sewuHOMPhjBputuIeIQiQuTdjyu57smPeuQRevPNNykoKOCyy6w+iy6XiwULFvCrX/2Kpqam\ntnlVVVWcdtppHH744Vx++eVp+8w9QVfHOMzs2fD731sFzH7+c6etyR2iPSEi6g1RcodU/q9vqtmE\nwbSLkM//DX++EQ6cCl/sehf+7seVbPi8gQE+D7//+1aOHzmo9zwhHtsT0pRa/NrndWEMNAfDbQ3t\nesJdL65l3Y66hHP2Lcnnksf/wZDSfHbWtXDIvsU89PomHnp9U8z5Y/cv5c6zD4+7v7Vr13L00Ud3\nGCstLWXEiBFs3ry53ba77uLEE09k9uzZvPTSSzzeD5qV6b2gwxx8MHzrW7BokdXgTkkP4bDmgChK\nd7StjHn+epgzABadBC4vfO1X4Op4Dxu54x89pJiD9ili4cUTOngEekrEE1IXbICAdfefVDjGYwmP\ndNYKKfN5GFKaz/aaZoaU5lPm86TluG+//Tbf+MY3APjyl79MeXl5Wo7bE9QTkgHccQcsXmz1llm0\nyGlrcoPOvWMUJddI5gd8w6ZlFIfD7F9tL4sNt0JrC2xZDuMu7DB3TUUtCy+ewL3L1lNe6G3LEVlT\nUdsr3pC2cExeHtRth8Gj2rYlDsdYP3NNgRCxy56lRiKPRYSIILv+5EP43d+3csMpo3r0HYwdO5Zn\nn322w1hdXR1bt27lkEMO4dVXX93rfTuNXoYzgOHD4aqr4Ne/ho97ZzWb0g2aE6Io3bPps5WMDgTo\noFdCLfDG3C5zrz5pJFNGDqa6McgA+85/ysjBXH3SyF6xpS0xtVOtkO4ozLc8If5gepJTIwJk4cUT\nuOm0Mb3iEZoxYwZNTU385jdWX51QKMTNN9/MpZdeSmFh+6qladOm8eSTTwLw8ssvU11d3bMPkwb0\nMpwh/OAH4HZb1VSVvkc9IUqukmxOiDGGjS7DqECw68YEIqDWH2RAYe8X3SnyWAmv9XlieUJo9+Z0\nl5gK0NiSnnBMxCMU8XxEe4T2FhHh+eef55lnnmHUqFGMHj2agoIC7o20Eba58847efvttzn88MN5\n7rnnGDFiRI8+SzrQcEyGsN9+cO21sGCBJUgOPdRpi7IbzQlRlMTsaNxBQ14eYwKBrhvLYldCDbSG\naWhpZUBh7+dAuPPcFLoLqc+rT8kT4vNYP3PpqhUSy/MzZeTgHoekhg8fzosvvthlfPr06UyfPh2A\nQYMG9bvQjN4LZhDf/z74fFb9EKVvUU+Ikut0J8I37NkAwOhwpxPF44MZs2O+p9ZveU36QoSAFZKp\n9xa2iZBk64QA2kk3Q9HLcAaxzz5www2wZAn8619OW5PdaE6IkqskG46JlGsfdcr94LYbK5UNh7Mf\n7pKUGqHWb3lN+iIcA5YIacj3dfGEJBOOaUpTToiSGnoZzjBuvhlKS7W7bl+jnhBFSUxbufYJ34CB\nh8CYM+HGf8cVIADVTbYnpI+WpJZ4S6h3e9tyQpIhUhvEn8Ylukry6GU4wxg40BIizz8PH3zgtDXZ\nSzjc/lwLlSm5SDLhmLYiZY27oWifbvdZY4uQ8j70hNTluSxPiDFJhmPal+gqmYeKkAzku9+1xMjs\n2GFXpReICA9XzwsoKkrW0aFcezgETZVJiZDqpkg4pu88IQ1iINgE/vblp0mFY1SEZCQqQjKQ0lK4\n9VZYtgzee89pa7KTcNh6uFxatj1bEZHTRWSDiGwWkdtibD9ARN4QkTUislxEYi/5yDKS+b/eoVx7\n0x4wYSjet9v31Tb1bWJqsaeYemOHVZJcIZPvzkMkfatjlNRQEZKhXHcd7Lsv/OhHTluSnUQuxG5d\npJ6ViIgLeAQ4AxgLzBKRsZ2mzQd+Y4wZB8wF7kuvlZlLW7n2gWOsUAwk7Qlx5wnF+X1zYpV6S6kP\nNWMA6rYnFY4REQo96e2k2xeICDfffHPb6/nz5zMnC5ZSqgjJUIqK4Lbb4I03YPlyp63JPsJhCIUs\nT4h6QbKSY4DNxpgtxpgA8DRwbqc5Y4E37edvxdie1ST6Ad+wZwPFnmL2L9ofGndZg8nkhPiDDCj0\nIH1UhKfYW0zIhPGLdPCEdHcOF+a78QfTmJi6ZiksOMLqt7PgCOt1D8nPz+e5556jsrJ3evFkCipC\nMpirr4b997e8IfpD2bsYYz3UE5K1DAW2Rb2usMei+SfwVfv5V4ASEenSXkRErhSRVSKyKht+AJIK\nx1RvYnT5aEtMNNiekCTDMX3ZrK2tdLvbm1rpdq+LxpY0eULWLIUXr4fabYCx/r54fY+FiNvt5sor\nr2TBggVdtu3evZvzzz+fyZMnM3nyZN555x0AjjzySGpqajDGMGjQoLay75dccgmvvfZaj+zpLfQS\nnMH4fHD77VYl1ddeg9NOc9qi7CGSE+J2a+XUHOYWYKGIXAq8DWwHuvxSGWMWAYsAKioqsv52wBjD\nxuqNfPngL1sDKYZj+qpGCESJkJIhDEkyHAPg681wzMu3wecJCjlVrLT660QT9MOfroMPnoj9ni8c\nCWfc3+2hr732WsaNG8ett97aYfyGG27gxhtv5MQTT2Tr1q186UtfYv369Zxwwgm88847HHDAARx8\n8MGsWLGCSy65hPfee49f/OIX3R4vHagIyXC+/W348Y8tb8ipp+oPZm+hq2Oynu3A8KjXw+yxNowx\nO7A9ISJSDJxvjKlJm4UOE+9asqNxBw3BBisfBKxwTJ4bCgZ0u8+apiD7DyjoRSs7UuKxO+mWDEkt\nHON1pS8c01mAdDeeAqWlpVxyySU8/PDD+Hy+tvHXX3+ddevWtb2uq6ujoaGBqVOn8vbbb3PAAQdw\nzTXXsGjRIrZv3055eTlFRUU9tqc3UBGS4eTnWwLk8svhz3+Gs8922qLsQHNCsp6VwCgROQhLfMwE\nLo6eICKDgT3GmDDwA+BXabcyA2kr1x6pEdJg1whJorpfTVOAw/Yr7TPbIp6QuqJy2P6fpN9X6HXT\n2FvFyrrzWCw4wg7FdKJsOFz2Uo8P/93vfpeJEydy2WWXtY2Fw2Hef/99Cgo6CsBp06bxyCOPsHXr\nVubNm8fzzz/Ps88+y9SpU3tsR2+hOSH9gEsugZEjrboh0UW2lL0nkhMCmheSjRhjWoHrgFeA9cBS\nY8xaEZkrIufY06YDG0RkIzAEmOeIsWmmO9G9sXqjVa59wChrIMlCZWAlppb30fJcaBchDb4yqNuO\nS5ILsfi8rvQt0Z0x2+qvE02CfjupMnDgQC688EIef/zxtrHTTjuNn//8522vV69eDVhN7yorK9m0\naRMHH3wwJ554IvPnz2fatGm9YktvoCIkiiTqCowQkbdE5CO7tsCZ6bDL47Ga2q1eDc89l44jZj+R\nnBCwvl8l+zDGLDPGjDbGjDTGzLPHZhtjXrCfP2uMGWXPudwY03N/eRawsXojw0uGU+gptAYadyUl\nQlpaQzQFQn1WIwSickK8hWBCDPTuBJILx6Rtie64C63+OmXDAem2387ecPPNN3dYJfPwww+zatUq\nxo0bx9ixY3n00Ufbth177LGMHm15taZOncr27ds58cQTe82WnqL3gDZRdQVOxcqkXykiLxhj1kVN\nuwPrjuoXds2BZcCB6bBv1iy4916rp8xXvqK5DD0l2hPi8YDf76w9ipJu4uWEbKze2J4PAlY4Zp9D\nu91fe6GyNCSmei1Pw5CCCmD/bt+XVhECluDoRdEB0NDQ0PZ8yJAhNDU1tb0ePHgwS5Ysifm+3/72\nt23Pp0yZQjjD3OnqCWknmboCBogEPMuAHekyzuWCu+6CdeusLrtKz4jkhBijnhAlt0jkNWgKNrG1\nbiujyke1T04yHFPdx9VSAfJd+XjzvNS7rPtnS4R0j8/j1gZ2GYqKkHaSqSswB/iGiFRgeUG+E2tH\nfVVX4PzzYdw4KzTTqudTj4h4QsJhK99O80IUpb1c+5hy2xPSUmet6kiqeZ3VN6avmtdFKPYWU2+7\ncfYpsBY8JRWOCYYwmoWecagISY1ZwGJjzDDgTOC3ItLlOzTGLDLGTDLGTBo8eHCvHTwvD+bOhU2b\nIMrDpuwFxrR7QoyxViEpSi4RKxzToVw7pFSorMZveUL6slgZ2KXbwy3gLWFIfpKeEK8LY6ClNbNC\nEYqKkGi6rSsAfBtYCmCMeQ8oAHpPZSTBOefApElWaCYQSOeRs4twuN0ToiJEUSw6lGuHlAqV1fRx\nB90IxZ7Q33INAAAgAElEQVRi6oP1UDaUfZMMxxRpJ92MRUVIO211BUTEi1VX4IVOc7YCMwBE5DAs\nEbI7nUaKwN13w6efwq+0qsFeE8kJURGi5BqJIhIdyrVDan1j7JyQvg7HlHhLqA/UQ9mwtpyQ7sMx\nVry1sUXj2JmGihCbJOsK3AxcISL/BJ4CLjUOBBm/9CWYMgXuuQeam9N99OwgOidERYiitJdrb0tK\nBWiwRUgS4ZjqpiAel1Do7dule20ipHQo++R3dlbHxmfb5A+qJyTT0HS8KIwxy7ASTqPHZkc9Xwec\nkG67OiNiCZCTT4Zf/hJuuMFpi/ofkZwQsL5PFSFKrtE5J6RLuXaARjuxvrD7qHOt3+ob01cddCOU\neEtoCDTAkOEM9O4m39WM5ZSOT6GGYzIWFSH9lC9+0Xrcd59V0j1D2gD0G6JzQkBFiJI7xPPddinX\nDlY4xjcQXN3/VFQ3BhnQx0mpEB2OsRYvDivdjjEjE77H1yZC+j4cM33JdKqaq7qMDyoYxPKLlvf5\n8fsbGo7px9x9N+zcCY884rQl/Y+IJySyQkZFiJLrdCnXDlY4JolQDECNP9Dn+SBgiZDmUDPBki8A\nMLys+5BMJCckHaXbYwmQROPJMm/ePEaPHs2JJ57IrFmzmD9/foftP/nJT3j44YcBuPHGGzn55JMB\nePPNN/n6178OwDXXXMOkSZM4/PDDufPOOwH4y1/+wgUXXNC2n+XLl3PWWWcB8Oqrr3L88cczceJE\nLrjggg4F03oL9YT0Y044AU4/3eqye801UFLitEX9i0g4RkWIkot0jpp0KdcOVjgm2b4xTUGGDyzs\nfmIPKfYUA1DvG8BAYERZ9ytkejMc88A/HuA/e5JvnhfNZX+5LOb4oQMP5fvHfD/u+z744AOefvpp\nVq9eTWtrKxMnTuToo4/uMGfq1Kn89Kc/5frrr2fVqlW0tLQQDAZZsWJFW6+YefPmMXDgQEKhEDNm\nzGDNmjWccsopXHnllTQ2NlJUVMSSJUuYOXMmlZWV3HPPPbz++usUFRXxwAMP8OCDDzJ7du/0wImg\nnpB+zty5UFUFDz3ktCX9j4gnJBRSEaIoXcq1Q9J9Y8ASIX3ZvC5CW+n2fEuMDC/d3u3qGJ/HTkzt\npzkhK1as4Ctf+QqFhYWUlpZyzjnndJlz9NFH88EHH1BXV0d+fj7HH388q1atYsWKFW1dc5cuXcrE\niROZMGECa9euZd26dbjdbk4//XRefPFFWltbeemllzj33HN5//33WbduHSeccALjx4/niSee4NNP\nP+31z6aekH7O5Mlw7rkwfz5cey2UlzttUf/BGKvyrCamKrlErB/sSLn2Lx/85Y4bGnanFI7py74x\nEUq9VueMBhOkOrAPI8q2dfOOaE9Iz3NCEnksAI584si42359+q97fPwIfr+f8ePHA3D11Vdz9dVX\nc9BBB7F48WKmTJnCuHHjeOutt9i8eTOHHXYY//3vf5k/fz4rV66kvLycSy+9lGZ7eeXMmTNZuHAh\nAwcOZNKkSZSUlGCM4dRTT+Wpp57qNZtjoZ6QLGDuXKithQcfdNqS/okx2hBQyW0212zuWK4dIOiH\nQH1SnpDmYIjmYLjPq6WCVbYdoC5Qx66WoUl5QiI5IU39dInutGnT+OMf/4jf76e+vp4XX3wRn8/H\n6tWrWb16NVdffTVghWTmz5/PtGnTmDp1Ko8++igTJkxARKirq6OoqIiysjJ27tzJyy+/3Lb/k046\niQ8//JDHHnuMmTNnAnDcccfxzjvvsHnzZgAaGxvZuHFjr382FSFZwLhxcOGF8LOfQS+2qskJ+ng1\noaJkLNH/9zdUx1oZk0q11PQUKoOocEygnl0tw5LKCSnw5CGSnnDMoIJBKY0nw8SJE7nooos46qij\nOOOMM5g8eXLMeVOnTuWzzz7j+OOPZ8iQIRQUFLSFYo466igmTJjAoYceysUXX8wJJ7RXm3C5XJx1\n1lm8/PLLbUmp++yzD4sXL2bWrFmMGzeO448/nv/8Z+9yYRKh4ZgsYc4cePZZK0n1xz922pr+g/az\nUnKNWP/nN+7ZSJGniKHFUT07U+gbU52mku0AJR5LhDQEG2gIDGP8oOV8ZgwQ/45CRCj0uGhs6XsR\n0lfLcG+//XZuv/12AObMmRNzzowZMwgGg22vO3suFi9eHHf/CxcuZOHChR3GTj75ZFauXLl3BieJ\nekKyhMMOg4svhoUL4fPPnbZGUZT+xMbqjR3LtUOUJySJ5nW2JyQtIqSTJ6QkvwFXa2237/N53fiD\nWrY901ARkkXceafV1O6++5y2RFGUTCeiNyLl2juEYiCqb0z31VLbmtf5+j4cU+gpJE/yqAvUsbvF\n8twUtCRTK8SVNRVT58yZwy233OK0Gb2CipAs4pBD4NJL4dFHoSK55pI5j+aEKLlG53BMzHLtkFLf\nmBq/nRNS1PeekDzJo8hTREOggd0tVuPzdIgQB9qE9Qt6+r2oCMkyfvQj6yIzb57TlvQP9Lqi5Dox\ny7WDVajMWwIeX7f7aAvHpMETAtYy3fpAPbsClickv6X7Zbo+r2uvE1MLCgqoqqpSIdIJYwxVVVUU\nFCTu3ZMITUzNMg44AK64AhYtgltvhYMOctoiRVEykYgXMGa5drALlXUfigErHON151HgSc99bYm3\nhPpgPdWBIQRD7hQ8IXuXEzJs2DAqKirYvXv3Xr0/mykoKGDYsGF7/X4VIVnID38Ijz9u9Zb51a+c\ntiazEVFviJLbxCzXDqn1jbGrpfZ1B90IxZ5i6gP1GHGxvX5/CgLdx599HjdVDYG9Op7H4+EgvaPr\nEzQck4UMHQr/+7/wxBPQB7VlFKVfICKni8gGEdksIrfF2D5CRN4SkY9EZI2InOmEnemms+iOWa4d\nUuobU90USFsoBqI66QLbaocm7Qnx99NiZdmMipAs5bbboKAA7rrLaUsyG/WCZCci4gIeAc4AxgKz\nRGRsp2l3AEuNMROAmcD/S6+VzhMp1z6qfFTXjan0jfEH07I8N0KJt4SGgNXRdWvtMAqSyAnJptUx\n2YSKkCxl333hO9+Bp56CtWudtkZR0s4xwGZjzBZjTAB4Gji30xwDlNrPy4AdabTPcUTilGsHCLVC\n054UwjGBtIuQ+kA9IrCtbigFgR0QDid8T6HX3W8b2GUzKkKymO99D4qLrfohSmx0iW7WMhSIvj2u\nsMeimQN8Q0QqgGXAd2LtSESuFJFVIrKqMgv6IkR7/2KWawdoqgJMih100xuOaQg2YAiztXYYeSbY\nXtckDpHEVF3hklmoCMliBg2CG2+EP/wBPvrIaWsUJeOYBSw2xgwDzgR+KyJdronGmEXGmEnGmEmD\nBye3WqS/ELNcO0QVKutehBhjqPEHKUujJ6TYU4zBYDyNbKuzV2bUJs4L8XldhA20tCb2mCjpRUVI\nlnPjjVBeDrNnO21JZqI3RVnLdmB41Oth9lg03waWAhhj3gMKgOxSGXEY89B03jr6SJ7e8DSNwUbG\n/WYcRz5xJNOXTLcmpFCozB8MEWgNp9UTUuq1omjGW8+2WltA1SVeIVPotVpla15IZqEiJMsZMABu\nuQX+/Gf4+9+dtkZR0sZKYJSIHCQiXqzE0xc6zdkKzAAQkcOwREhOFILwlFXFHK9qtscb7bBTKn1j\nfOnNCQFLhGytjXhCkhUh2j8mk1ARkgNcfz0MHqzekFhoTkh2YoxpBa4DXgHWY62CWSsic0XkHHva\nzcAVIvJP4CngUpMDCQNJfcIU+saks4NuhGJvMQBhTz3VzeW05hUmEY6xymJpcmpmocXKcoDiYmvJ\n7i23wIoVMHWq0xZlDtn/k5O7GGOWYSWcRo/Njnq+Djgh3Xb1Cxp2gcsLBWXdTq1t66Cb3sRUsDwh\nIDR7h1Fcm3iZbqFHwzGZiHpCcoRrroEvfAHuuEN/eBVF6YbGSisUk4SrsLpNhKTPE1LqsXJCwl6r\nVojfOxTqEntCNCckM1ERkiMUFsLtt8Pbb8MbbzhtjaIoTpF0OCbZvjF+KxyTzsTUSDjGeOoAaM4f\n1m1OiM8WIf6g5oRkEipCcogrroDhw9s77SqaE6LkJsHaQTHHBxXY4yn2jQEoS2NialtOiNcq3e73\nDoOGndDaEvc9hXZOiHpCMgvNCckh8vOtcMxVV8HLL8OZOdEpIzEqxpRcZMMNyzn1wdsxw//Ba197\nreuExkr4wrik9lXTFMDncVFg51ykA0+eB5/bh/HY4RiPvUKmbgcMjN1oTsMxmYl6QnKMyy6Dgw9W\nb4ii5DoBd1W75yMaY6Bxd/LhmKb09o2JUOIpifKERGqFxM8LiYgQXR2TWagIyTE8Hmup7ocfwh//\n6LQ1iqKkm8jNR4u7ikG+GCLEXw3hYNLhmOqmYFpXxkQo8ZYQ9lgipNnbfdXUSDimUeuEZBQqQnKQ\nr38dxoyxxEg3PZ+yHs0JUXKVgCeOJySFQmUAtf5AWguVRSj2FmNsEdLksT0hCZbpFnjyEFFPSKah\nIiQKETldRDaIyGYRuS3OnAtFZJ2IrBWRJ9NtY2/gdsOcOfDvf8PSpU5b4yyhkAoRJQeRMEH3ntie\nkBQKlUHEE+JAOMbbHo4JuwrBNzBhOEZE8HlcmhOSYagIsRERF/AIcAYwFpglImM7zRkF/AA4wRhz\nOPDdtBvaS1x4IRxxhCVGWnPYO9naCt70e5IVxTGMAVdRLUZCDPbFEBop9I2BSE6IU+GYhvaBsu6X\n6VqddFWEZBIqQto5BthsjNlijAkATwPndppzBfCIMaYawBiTuHd0BpOXB3PnwoYN8PvfO22NcwSD\n1qohRckl3KVWj5iehmOMMdQ0BZxLTLXDMcZgi5DuO+n6NScko1AR0s5QIDqgWGGPRTMaGC0i74jI\n+yJyeqwdiciVIrJKRFZVVlb2kbk957zzYOJES4wEg05b4wzBIBQUOG2FoqQXd5l1XYobjpE8KBzY\n7X4aAyFaw4Zyx8IxdYCdaZuMJ8TjVk9IhqEiJDXcwChgOjALeExEBnSeZIxZZIyZZIyZNHhw5nYG\nF7EEyJYtsHix09Y4g3pClFzEXZbAE9KwCwoHQV73dT+qG+3mdT5nwjHktSKeFssTUjoUWmqhpT7u\ne3xeF/6gipBMQkVIO9uB4VGvh9lj0VQALxhjgsaY/wIbsURJv+XMM+G44+Duu6ElfrHBrEVFiJJr\nGBMVjonpCalMYWVM+vvGRIg0sXMV2qKjLJllupoTkmmoCGlnJTBKRA4SES8wE3ih05w/YnlBEJHB\nWOGZLek0srcRsQTItm3w2GNOW5N+WltVhCi5h7usCgm7KfWWdt2YSt8YBzroRogvQuKHZFSEZB4q\nQmyMMa3AdcArwHpgqTFmrYjMFZFz7GmvAFUisg54C/ieMabKGYt7jxkz4KSTYN48aGpy2pr0Egio\nCFFyD3dpFd7QQCTW+vQU+sZUN0Wa1zlQJ8Rj9Y/JK6xvD8cA1CUSIW5NTM0wVIREYYxZZowZbYwZ\naYyZZ4/NNsa8YD83xpibjDFjjTFHGmOedtbi3iHiDfn8c/jFL5y2Jn0Yo+EYJfeIhGPyg3G8HSmE\nY2rscExZJoRjSvazEmrVE9KvUBGiADB1Kpx6Ktx/PzQ0dD8/G4jUR1ERouQa7rIqvK0x8kECjRBs\nTD4c42BiaiSU1CZCXG5LiCTICfGpCMk4VIQobdx9N1RWwsMPO21JegiFrLL1KkKUXMNdWkV+KM7K\nGEi+UJk/SJHXhded/p+SYq8djvHVtzfjLBuWsHS75QlpxWj3zoxBRYjSxrHHwllnwU9+AjU1TlvT\n90Q8IZ70e5IVxTEMYVyle2J7QlLsG1PdFHAkKRVihGPAygtJ2EnXTdhAS2uON83KIFSEKB2YO9cS\nIAsWOG1J39PaCi6X9o5RcosW6shzt5IfU4Sk1jem1qG+MQAFrgIIu3H5Opdu397eKrgTPo9V+0Sb\n2GUOKkKUDkyYAOefb4mQqn6/7icxoZB6QZTcw59nndgxPSF7EY5xSoSICHnBkvbVMWCJkFBLu0en\nE4VeS4Q0acGyjEFFiNKFu+6yklPnz3fakr4l4glRlFyi2RYhBbFyQtrCMfsktS8nwzEArmAJrsK6\n9oFIrZA4y3R93ognRJfpZgoqQpQuHH44zJplJaju3Om0NX1HKARut9NWKEp6SegJadwF+WXgTi5b\nu7YpyACfc+7EvGAJrsKocEykVkicZbqFXuuE1xUymUNCESIiA0Tkf3vjQCKSLyJLRGSziPxdRA6M\nM+8TEfmXiKwWkVVR4wNF5DUR2WT/Le8Nu5TY3HknNDfDAw84bUnfoZ6Q7EZETheRDfY157YY2xfY\n15nVIrJRRHIgHbvdE5LfGiPvo2EXFCfnBTHGOBqOAcgLFncNx0DcZbpt4RgVIRlDd56QAUCviBDg\n20C1MeYQYAGQ6Ofti8aY8caYSVFjtwFvGGNGAW/Yr1NChUvyjB4Nl1wC/+//wfbE3bH7La2t6gnJ\nVkTEBTwCnAGMBWaJyNjoOcaYG+3rzHjg58Bz6bc0/fjzqjCtbrzhWCXbky9UVt/SSihsKHcwHJMX\nLO24OqZwELgL4i7TbQ/HqAjJFLoTIfcDI+07hZ/08FjnAk/Yz58FZkjMmsFJvf8J4Ly9sGGViPxe\nRE5O8dg5yezZVsji3nudtqRv0HBMVnMMsNkYs8UYEwCexrqGxGMW8FRaLHOY5rwqWusHIrEu/6n0\njWm0q6U6GI5xtRbj8kWJEJGEy3SLNByTcXQnQm4DPrbvFr7XeaOIrIhyZ0Y/Tomxr6HANmjr01IL\nxAhKYoBXReQDEbkyanyIMeYz+/nnwJBubI/FaKwLzXXAOhH5oYjsvxf7yQkOOgi+/W2rsd2nnzpt\nTe8TDGo4Jotpu97YVNhjXRCRA4CDgDfjbL9SRFaJyKrKytirLvoT/rwqWmtjXXpJqW9MjT/SN8ZJ\nT0in1TFgL9ONlxNinfCNmpiaMfQoMdUYMzXizuz0eL0Huz3RGDMRy416rYhMi3FcgyVWUrU3ZIz5\nszHmq8A04GBgq4gc0wN7s5o77oC8PKuaarahS3QVm5nAs8aYmLfHxphFxphJxphJgwcn5yXIZJrz\nKmmtG9h1Q2sAmmtSKFQW6aDrpCekBFeBn5CJEhWRWiEx0HBM5tEjEZKiJ2Q7MNx+nxsoA7pUojDG\nbLf/7gKex3KrAuwUkf3s9+8H7Iphz6/t4y8TkeFR9lwdNadMRK4CXgBGAd8C1uz9t5DdDBsGV10F\nixfD5s1OW9O7aGJqVtN2vbEZZo/FYiY5EooBaHZV0Vo3qGuRvqbI8twkwzF2B10nl+jmBa2qqS10\nKlhW/xmEgl3ma2Jq5tFdRLweKIm30RgzNYVjvQB8E3gP+BrwpulUwF9EioA8Y0y9/fw0YG6n999v\n//1TDHsu6zQ0vtP+fwccDzwDXGKM2ZSC/XvFBx98wPDhw7uf2A8YNcppC5Rs5uOPP+7N3a0ERonI\nQVjiYyZwcedJInIoUI51Xcp6jDE05+2xwjGdIzIpFiqr9WeAJ8QWIf5wPdY6CuxlusYSIgNGdJhf\n4NY6IZlGQk+IMaYKeEdE/t0LiamPA4NEZDNwE/bqFhHZX0SW2XOGAH8TkX8C/wBeMsb8xd52P3Cq\niGwCTrFfp8pSYIwx5rZ0CBBFUZzBzju7DngFWA8sNcasFZG5InJO1NSZwNOdb4iylbpAHWEJ0lrX\nC31j7MRUJ+uEuFptT4iJSk7ds8X6+7NxsOAIWLO0bVNenuDzaCfdTKLbtQHGmC53D3uDMaYZuCDG\n+A7gTPv5FuCoOO+vAmb00IYXevL+veHoo4/Oio6NlZVWouqZZ8KSJU5b0zM2boQ33oAdO+CYY+Ds\ns522SAHYsmVLr+7PGLMMWNZpbHan13N69aAZTlWzFQFvrR3cNRyTYt+YGn+Aknw3bpdzNS/zWq1O\nus3YImTNUvj7o/ZWYy3VffF66+W4CwG7k66Wbc8YtGKqkhSDB8N3vwtLl8I//+m0NT3HGM0HUXKP\nKn9EhPRC35imIGUOhmLAqhMCkXAM8MZcaG3uOCnot8ZtfF6XJqZmECpClKS56SYoK7OqqfZ3QiFr\n1Y+i5BJtnpCY4Zjd4PaBtzipfdU0BRxdngvgjoRjIp6QOEtzo8cLvS6aNCckY9DLsJI05eVw883w\npz/BqlXdz89kjFERouQebZ6QeCKkaB+6xmliU93kbMl2sMq2AzRHckIiZds7EzXu87o1JySD0Muw\nkhI33ACDBsGPfuS0JT0jFNJwjJJ7VPmrEOMi1DCgq9ZIoW8MWKtjnFyeC1ZOiAkL/ogImTEbPL6O\nkzw+a9ymSMMxGYWKECUlSkvh1lvhL3+Bd95x2pq9JxxWT4iSe1Q1V5EfLgcTq2R78n1jwArHOLky\nBiBP8gg3F9Fi7Doh4y6Esx+GfLuyRNkw67WdlAqRcIyKkExBL8NKylx7LQwZ0r+9ISpClFykyl9F\nQShOyfYU+saEw4Zaf5Byh8MxAKGmEpqpax8YdyGcMsd6/u3XOggQiIRjNCckU9DLsJIyRUXwgx/A\nW29Zj/5IOKzhGCX3qPJXURC2hEaHcEw4bHlCklwZU9/cSthAmcPhGLBESJsnJMKAA62/NVu7zC/U\nOiEZhYoQZa+46ioYOtTyhvTHMijqCVFykarmOJ4QfzWYUAp9YyLN65z1hIhA2F/cXickQvkB1t/q\nrp03dYluZqGXYWWvKCiwmtu98w688orT1qSOihAl1zDGxA/HpFyozPmS7RFCTaXtq2MilNmtMmq6\nipBIsbJsKCKZDehlWNlrvvUtOPDA/ukN0XCMkms0BBsIhAPkhy0R0iEck2KhsognpMyXKeGYTiLE\nUwDFX4grQkJhQyAUTpOFSiJUhCh7jddrCZBVq+CFtBfE7xnqCVFyjUiNkILQwK4bG3dbf5MMx9Q2\nWZ6QjA3HgBWSiRmOsbqVaEgmM9DLsNIjLrkEDjkEZs+2ftj7CypClFyj0m81qMuPGY6JiJDk6oRE\nPCFO1wkB2xNCQ9fwyoARcT0hgCanZgh6GVZ6hNsNc+bAmjXwhz84bU3yhMNWCMndbQtHRckOIiXb\nY+aENOwCcYGvPKl91diekDKH64SAJUIMYZpamzpuGHAA1G6HUMfluCpCMgsVIUqPmTkTxo61esqE\n+sl5HREh+flOW6Io6aE9HBMjJyRSsj1J92CtP0hpgRtXXnIl3vsKEQg3WYXJ6gMxVsiYENRt7zBc\nqOGYjEJFiNJjXC646y5Yvx6eesppa5IjHLYeKkKUXKGquYo8ycMbiuHtiIiQJKluClBe5HwoBixP\nCEBdoK7jhgEjrL+dQjLtnhAtWJYJqAiJQkROF5ENIrJZRG5LMO98ETEiMimd9mUyX/0qHHWUFZoJ\nBp22pntCIfWEKLlFlb+K8vxyhBjLwlLsG1PTFHS8ZHuEkN9qYtcQ6FywLHatEJ+GYzIKFSE2IuIC\nHgHOAMYCs0RkbIx5JcANwN/Ta2Fmk5cHd98NH38Mv/mN09Z0jzHQ2qoiRMkdqpqrGORrzwfpGI7Z\ni74xGZCUaoVjSoEY4ZiyYSB5Xaqmak5IZqEipJ1jgM3GmC3GmADwNHBujHl3Aw8Azek0rj9w1llw\nzDEwdy60tDhtTWJCIQ3HKLnFHv8eBhXESEo1JqW+MWAVK8uEQmXQHo6pD3YSIS4PlA7rGo7xWDkh\nGo7JDFSEtDMU2Bb1usIea0NEJgLDjTEvJdqRiFwpIqtEZFVlZWXvW5qhiFjekK1b4fHHnbYmMeGw\nJURUhCi5QsQT0qWwYEs9tDYnXagMoLrR+Q66EUJNVjimiycErLyQOOEYf1A9IZmAipAkEZE84EHg\n5u7mGmMWGWMmGWMmDR6c/N1FNnDqqXDiiTBvHvj9TlsTH01MVXIJYwyV/srYnpAUC5WFwoa65tbM\nCcf446yOAWuFjIZjMhoVIe1sB4ZHvR5mj0UoAY4AlovIJ8BxwAuanNqRiDdkxw549FGnrYmPekKU\nXKIx2EhLqCV2TkiKhcpqM6hvDIBp9eKmoGtiKljJqfWfQWt7fNjnURGSSagIaWclMEpEDhIRLzAT\naCtGboypNcYMNsYcaIw5EHgfOMcYs8oZczOX6dNhxgy4/35obHTamthoYqqSS0QKlQ32De4ajmnr\nG5OcCKlp66CbGZ4QgHyKuy7RBXuZroGa9kh7Xp5Q4MnDrzkhGYGKEBtjTCtwHfAKsB5YaoxZKyJz\nReQcZ63rf9x9N+zaBQsXOm1JbFSEZD/JLLkXkQtFZJ2IrBWRJ9NtY7qIFCrrjXBMpINuWYZ4QgDy\nKYkfjgGo+aTDcKHXrZ6QDEGLVkdhjFkGLOs0NjvO3OnpsKm/cvzxcOaZ8OMfwzXXQGmp0xZ1RTvp\nZi9RS+5PxUoyXykiLxhj1kXNGQX8ADjBGFMtIslnZvYzIp6QxOGY5PLXMskTEiGfEhqCccIx0CUv\nxOdxacXUDEE9IUqfMXcu7NkDP/uZ05YoOUgyS+6vAB4xxlQDGGN2pdnGtNHmCfHF6RvjK7eWtCZB\npG9MJqyOiQipgniekJIvQJ6nywqZonyXekIyBBUhSp9x9NFw3nnw059aYkRR0ki3S+6B0cBoEXlH\nRN4XkdPTZl2aqWquQhAG5A/omhPSuDulQmXVTZmVmAoJwjF5LhgwvEutEJ/XTZMu0c0IVIQofcrc\nuVBfbwmRTEKkU8VIJRdxA6OA6cAs4DERGdB5UjbU/anyV1FeUI47L0YEPsW+MbVNAUSgtCBzRIiX\n4tgiBKyQTOdluh4XTS2amJoJqAhR+pQjj4QLL4SHHoLdu522RskhultyD5Z35AVjTNAY819gI5Yo\n6UA21P2p9FcysGBgh7E2EZ5i35jqpiBlPg95DnfQhSTCMWAlp1Z3bWKn4ZjMQEWI0ufMmWMVLnvg\nAactUXKIhEvubf6I5QVBRAZjhWe2pNPIdBHdN6ZrOCbFvjH+YEYlpQLkU0ogHKAlFKNfxIAR0FQJ\nLalmTboAABv0SURBVO2Jqz6vSyumZggqQpQ+59BD4RvfgEcegc8+c9oaJRdIcsn9K0CViKwD3gK+\nZ4ypcsbivqXKX8VgXwwvTrAZWmpTCsfUNAUoy4Ck1Gi8JCrdbq+QqW1PEbI8IRqOyQRUhChpYfZs\nCAbhvvuctkTJFYwxy4wxo40xI40x8+yx2caYF+znxhhzkzFmrDHmSGPM085a3Hfsae7avE6E9uW5\nKYRjapqClGdIUmp0OAbilW4/0PobFZLROiGZg4oQJS2MHAnf+hb88pdWgztFUdJDU7AJf6s/9vLc\nFAuVAdT4AxnRNyaafFuExC7dPsL6G7VCxufVOiGZgooQJW3ccYf1d948Z+0Aq1CZro5RcoHO1VI7\n5ISk2DcGoKYxmFHLcwG8JoEnpGgf8BR29IR4XLSGDYHWcLpMVOKgIkRJGyNGwBVXwK9+BVscTv8L\nh8Gt9YKVHCBWtdQ2UuwbEwyFqW9pZYAvMzwh7b1jLBFSF4zRP0bE8oZ08oQA6g3JAFSEKGnlhz+0\nfvznznXWDhUhSq4Qr29Mh5yQftpBN0K+nZgaMxwDXURIUb518jcFNTnVaVSEKGll//3hf/8Xfvtb\n2LDBOTu0b4ySK1T6rQJrMZfoNu4GTxF4i5LaV00GVksF8BqrOVXCgmXV7clohbYnRJNTnUdFiJJ2\nvv998Pms+iFOoZ4QJVeIlGwvLyjvujHFQmWR5nWZkpgaCcd48OESV+KCZS214K8BrAZ2AE0tKkKc\nRkWIknb23Reuvx6WLIF//csZG0IhFSFKblDlr2JA/gA8eR29F23hmFRWxtiekExZotuGEYq93ZRu\nh7aQTKHXDsdorRDHURGiOMItt0BJCdx5pzPHV0+IkitEV0vtQop9Y2oiOSEZkpgaTYmnhPpgPBFi\nL9O1V8hEElO1iZ3zqAhRHGHgQLjpJnj+efjww/QfX3NClFyhyl/VISm1Q07I3oZjijLDExK9zL7E\nWxI/MbU84gmx8kIKdXVMxqAiRHGM734XysutaqrpRj0hSq5Q1VzFQN/ALuOHh562eqp8sBgWHAFr\nlna7r5qmIK48oSQ/s04eYywREjcc4yuH/LKocIwmpmYKKkIUxygrg1tvhZdegvfeS++xVYQouUJn\nTwjArCOWcpa5oX2gdhu8eH23QqS6KcAAnwfJwEp/Jd4S6gIx6oREGDCiSzjGrzkhjqMiRHGU666D\nffZJvzdEwzFKLtAUbKKptalDTogxcO+MuXho7jg56Ic3EhfwqfEHKcugpNQu4ZhgnHAMWCGZtnBM\nJDFVPSFOoyJEcZTiYrjtNnj9dfjrX9N3XF0do+QCbdVSO3lCRpRVxH5DbZxxmxrbE5JpGAPFngSr\nY8BaIVPzKRjTvkRXRYjjqAhRHOeaa2C//eBHP+qUNNeHqCdEyQUi1VIH+wZ3GN9aOyz2G8rijNtY\nHXQzZ2VMtCek1FtKY7CRUDiOsBgwAoJN0FiJK0/Id+fh19UxjqMiRHEcnw9uvx1WrLA8IulAc0KU\nXCBe35gfvjEbQ6e8Do8PZiSOi9Y0ZVY4Jppir126PV5IprxjrZCifLfWCckAVIQoGcHll8Pw4Van\n3XR4Q1SEKLlArL4xxsDLm08FDOSXAgJlw+Hsh2HchQn3V9MUyChPSITI6hjopnQ7tIkQn8el4ZgM\nQC/DSkaQn28lp15xhbVa5qyz+vZ4xnR05SpKNhLxhHReonvC8L9bfpCZT8JBU5PaV6A1TGMglFE5\nIZ0TUyGBJ6RTwbJCr0vLtmcA6glRMoZvfhMOPtgSI+Fw3x3HGBUhSm5Q5a+iLL+sS8n2E0e8RwgP\nDJuU9L5q/JFCZZnnCQGrYiok8ITkF0PhoA61QrRiqvOoCFEyBo/Hamr30UdWJdW+IiJC8vP77hiK\nkgnsad7TZWWMMTB1xHt8JhOsPJAkqY100M0gT0iE6HBM4loh7ct0fV6X1gnJAFSEKBnFxRfDoYda\nPWVCfXSTEg5bDxUhSrZT5e/aN8aNn8lDP6Qi77iU9lXd1rwuczwhMcMx8Uq3g5WcWt3exE5zQpxH\nRYiSUbhcljdk7Vqry25fEA6rJ0TJDSr9lQwu6Lg8d7h8gNcVpCJvSkr7ausbk6GrY7pNTAUrL6R2\nG4TDtidERYjTqAhRMo4LLoAjj7TESGsfeEs1HKPkCrE66B6UZ/VIqMg7NqV91diekLIMDccUeYqA\n7kTIARAKQMPnFOrqmIxARYiSceTlwdy5sGkT/O53vb9/YzQco2Q/za3NNAYbu4iQA/Pe5V87x9Is\nXZvaJSKSmFqeQYmp0eEYd56bQnch9cEEIiRSK6T6UysxVXNCHEdFSBQicrqIbBCRzSJyW4ztN4nI\nOhFZIyJviMgBTtiZC5x7Lhx9NNx1FwQCvbtvzQnJDZI4ny8Vkd0istp+XO6EnX1FzJLtoVYOyPsH\nK7Yen/L+apqCuPOEIm/mlhpO2EkXOtQK8XndWjE1A1ARYiMiLuAR4AxgLDBLRMZ2mvYRMMkYMw54\nFvhxeq3MHUTg7rvhk0/g17/u3X1rTkj2k+T5DLDEGDPefvxfWo3sY9oKlUV7Qnb+i3xpYMXWKSkv\nUa9uCjKg0JuRHXQjBQ67FSFlw62/1Z9S5HURDBmCoT6sB6B0i4qQdo4BNhtjthhjAsDTwLnRE4wx\nbxljmuyX7wOJGy0oPeL00+H44y0x0tzc/fxk0XBMTtDt+ZztxKqWyqdWPsiKT1P3hNT6AxmXlNpZ\nD5V4SxKvjvEUQMl+ULMVn1eb2GUCKkLaGQpsi3pdYY/F49vAy7E2iMiVIrJKRFZVVlb2oom5hQjc\ncw9s3w6LFvXefjUxNSdI9nw+3w6vPisiw2PtqL+ezzH7xvz/9u4+Oqr6TOD498nLEPJOEgRKEsCA\ndlHUKNJKa41Abe22VFd6ULS+1B6qXcS+2cVDfcNatT0UW7GHWu2KdvGNY7dodQXFuGtBQcVqwSIJ\nIoSXkEwgRDLJTJLf/nEnyWQyk5mEmbkzN8/nnJxk7vzuPc89mXvz5Pfy3L2baeoqZ3/LQLe20I4c\n9zEqyZKQYHmuvIHrhEDP03SzXVbBcJ0XYi9NQoZARK4CpgO/CvW+MeZhY8x0Y8z0kpKSUE1UlGbN\ngqoq+MUvoLU1YvOo6JwQ5fc8MNE/vLoBWB2qUapez909IUVZ/gmoxsDezXzcZS3NHeyoylGPj4KR\nyTMpFawChwDt7db3PFde+LLt3QrLeyamgvaE2E2TkF77gcD/hEr92/oQkTnAUmCuMaY9QbENa3ff\nDfX18NBDsTle95yQzOT+p06dmIjXszHGHXANPwKck6DYEqLR00i+Kx9Xuj9xcNfC8QY+7hxckbJu\nR1uTbzhm9Gjre0OD9T03M3fgOSFgrZA5tp/sDGsiidYKsZcmIb22AlNEZJKIuIDLgXWBDUSkEvg9\nVgJy2IYYh6UvfhG+8hW4/35oiXB/iUb3c2OScH6dip1orudxAS/nAh8mML6461cjZO8mAD7uHFyR\nsm5HW5NvOOakk6zvh/1343xXPi3eFsxAj+IunACmk8IOayftCbGXJiF+xpgOYBHwMtbN6BljzHYR\nWSYic/3NfgXkAs/6l/StC3M4FWPLloHbDb/5zYkfa6D7k3KGKK/nxSKyXUT+DiwGrrUn2vhwe9z9\nJ6VmF3PYnAIMLglv83Xi8XVSmEQl26E3Camvt77nufLoNJ14Ojzhd/I/Tbew7QCgc0LslmF3AMnE\nGPMi8GLQttsDfp6T8KAUADNmwNy5sHw5LFoEhYVDP1Y8n9CrkkcU1/OtwK2JjitRmtqaOLXo1N4N\nezdB+XnQMPguwGaP/+F1SdYTMmaM9b27JyTXlQtYVVOzM7ND7+QvWJbrOQCM1+EYm2lPiEoZy5bB\n0aPw61+f2HG6h2OUcrI+PSHHDsKRPVB+3qB6Ale9Xsum2saeku2FI11sqm1k1eu1sQ94CIKHY6J6\nfkx+KUg62a3WFCEdjrGXJiEqZZx5JsybBytWwImslOzq0iREOVt7ZzstvpbeOSH++SBMGFx9kDNK\nC1i0Zhtv1FgzPw8cbWXRmm2cUVoQy3CHrF8Skul/ku5AK2TSMyB/PFmfWiu4W7Vqqq00CVEp5a67\n4Phx+FXIxdHR0TkhyumaPE1AQKGyvW9CZg6MPbOnTTSJ+MyKElYuqGTFhl0APLixhpULKplZkRxL\nlfPzweWCTz+1lvB394RErBUyagKZLVYS4tE5IbbSJESllKlTYcECePBBOHRoaMfQnhDldN2FykpG\n+pOFTzZD2bmQnjHoJHxmRQmnfSYfgG9NL02aBASs6zhwXkh3EjJg1VSAwgmkHfP3hOhwjK00CVEp\n5447rIfa3Xff0PbXnhDldH2eG+M5CvX/gPKhLc3d+M96tnzcxMTibP687QCbapOramzgkExUc0IA\nRk1AWg6Sl9GhE1NtpkmISjlTpsA118CqVVBXN/j9dXWMcrpGj5UoFGcVw74tgOk3HySa3sBNtY38\n+39twwAr5p/FygWVLFqzLakSkZBJiC9CEuJfpluR2cRxHY6xlSYhKiXddpuVTNxzz+D31dUxyum6\nh2OKRhZZk1LTMmD89EEf5509R8hMF86fUkJl+aieOSLv1zXHOuQhC6wVMiJ9BK40V+SekEJrme6k\njCYdjrGZJiEqJU2cCN/9Ljz6KOzZM7h9tSdEOZ3b4yYvM48R6SOs+SDjzgKXVTdjMMOROSMyONbW\nwU2zpvRsm1lRwg0XVMQ65CELVSskqtLtwIT0Bh2OsZkmISplLV0KaWnWs2UGQ3tClNP1lGz3tcGB\ndwe9NBesKqm//99aPjepiBmTiuIQZWyEK90+oNyxkO6iTA5rT4jNNAlRKWv8eLjxRli9Gnbtin6/\nTr3nKIdze9zW03P3vwOd3pCTUiMl4s++U0f9sXYWz54ycEObhSrdHnFOSFoaFJQxzmhPiN00CVEp\nbckSGDHCqh8Src5O6x6klFP19ITs3WxtKO99cm40wzHeji5WVddydnkhMyuKI+9go37DMdE8SRdg\n1ATGdtXT6tOJqXbSW7FKaWPGWM+SWbMGtm+Pbp+ODk1ClLO5PW6rRsjezTD6XyB7cMMpf95Wx/6j\nHm6aPQVJ8rHLUKXbo0pCCssZ3XlIh2NsprdilfJ++lPIzYU774yufWcnpKfHNSSlbOPt9HLMe4zi\nEUXW8tww80HC5RYdnV089FotZ5QWUHXK6DhGGhuhkpCIxcoACieQ19kM7VG0VXGjSYhKecXF8IMf\nwNq18N57kdv7fDoxVTlXU5u/ZLuvDdqPDbpI2bq/H2BvUyuLLpyc9L0gAKP9eVJDg/UPRtQ9If4V\nMqO8Qyy9rGJCkxDlCD/6ERQWwu23R26rc0KUk/VUS20+aG0I6gkZaE5IZ5dh5Ws1fHZsHl+eOiZe\nIcZUZiYUFVlL75uarCSkrbMNX6dv4B0LJwJQ0qFJiJ30VqwcobAQfvITeP552LJl4LY+nyYhyrm6\nC5UVN9ZAQRkUlEa974sfHGR3w3FumpX8c0ECBQ7J5GbmAtFXTR1n6vF1avEgu+itWDnG4sXW0Mxt\ntw3crqtLkxDlXD09Ifs/gPLw9UGCc4yuLsPKjTVMPimXi08fG88QY25Iz4/JKcGXlkWpNOjkVBvp\nrVg5Rl6etWR3/Xp4443w7bxenROinKunJ+TTwyEnpYYbjlm/o56d9S0sunAyaWmpdYEE1grJd1lP\n/I04OVWE49njKROtFWInTUKUo3z/+zB2LPzsZ+Fvtj4fZGQkNi6lEsXtcZOT5iLLmKgnpRpjeHDj\nLiYWZ/P1M8bFOcLYC6wVkuuyhmOOeY9F3K8tp5QyaaBVH2JnG70VK0fJzoZbb4Wbb4aNG2H27P5t\n2ts1CVHO5fa4KSEdRhbB6FPDtgvsDXxt52G2HzjGL+edQUZ66v1v2t0T8nhOFX/4H6snaOGGhT3v\nF2cVUz2/ut9+7blllMpW9mpPiG1S79OmVAQLF0JpqTU3JFRviNdrzahXyoncbW6KvR5rPkgU447G\nGH77ag2lo0ZyaeX4BEQYe91JSHuGO+T73UNUgVa9Xkt7i5t8aeW0P0yAFafz0YZHWfV6bTxDVUE0\nCVGOk5VlDcds3gwvvdT//bY27QlRzlP1dBXTVk9jy6EtvJsB07wfMG31NKqerurTLjgxf6Omkff2\nHeXGqgoyU7AXBHqTkMGY5a1m4uENAAgGmvdR9rclzPJWxzY4NSC9FStHuu46uP9+q27IxRf3/Yew\nvd1KVJRyklD/7Q+0vduDr9YwriCLeedEv5Q32YyJoqTJpv2bOHfcuWSmWd2gC+tW4J70mX7tivet\n4BXft8nIdAFWj8ksbzWn/GMFNNdBQSkfnf5DNrqquOGCipiex3CkSYhyJJfLSkCuuw7+8he45JLe\n99raYORI+2JTiSEiXwV+A6QDjxhj7gvT7jJgLXCuMebtgY65w72Df3vtCrZ8+//6bJ/xxPl4uo72\naz8yrTDubVt9rexr2TdQ2P3kz6hln7eAN3cLW/Y0cec3pvLOJ0d4v645Jf+wRtMT8r1XvkeeK48L\nyy5kdvls3GFWALnTBc/Py6nJqaS97Hz+5H2Oh2iHIoGiMqvRgQcoTHuYG3irz752fg7i2XYox1x+\n1vJzTubkfvsE0yREOdZVV8G991pzQ+bO7a0N4vFYK2iUc4lIOvAQ8GWgDtgqIuuMMTuC2uUBN0PQ\nX5MBeLqO8tiaXf22Jbrtvz5+JUfMfloYuKcD4IlnvD0/7/wYfE05vND0Lh++MJKS3BFMLMlh0Zpt\nrFxQGfFYySiaJGTByAfZ4dvAhtpq1tWuG7DtczkzqfBu5zO7t9Bc2r+3BOBoVyuP/m4NkubCpGVB\nWiaezPC/r9V/+hBDd+Ijtnxmhto2FscMR0w0z3VWQ1ZXV2dKS1O3mzPVPfkkLFgATz0F8+db2668\nEqZOhaVL7Y1N9bV79+53Tj755OmxOJaInAfcaYz5iv/1rQDGmHuD2j0AbABuAX4SqSdk/d/Xmx+/\n9+NYhHjCKtvaKPN1MMHXQXlHB7ecVBK2bcuHITuBAPji5BJ2HDzGygWVzKwIf4xkZgzk5ED5fVVk\nFvRPynzNxey8udp6ke4j57NbmXTL9xIb5DCz/KzlXHTmRRFnRmtPiHK0+fPhnnvgjjvgssusCale\nr3XDUo42Hggco6gDPhfYQETOBsqMMX8VkVvCHUhEFgILAV5+72UAPr/37D5t3ix/N2wg8Wp79vZr\nAWgF/glw0lNh25YentrndXY2fP3rsO1gI6/tPMziWZNTNgEBa87XypXwwgvVYdtMvbT7p0xgJrvC\ntoRxbz2AkU5I6+DQ9P8I2+7UmktIEx9p0kGadPLBpFfCtv3c3nP6vH6r/J2wbZPh8xXYNlbHDEWT\nEOVoaWmwbJmVgKxZA1dfbc0Jyc21OzJlJxFJA34NXBuprTHmYeBhsHpCAP5w2+o+baatnhZ2/3i1\n/cFv+3bl/fd/PhlynkNxl6H6j5P6bd9U28hjW2tYPGsyf3prL5+vKE7pROQ737G+ojVtdfj31v9u\ndkC78EnI2rvvDjpm+N/XI7c9FnXbZPh8BbaN1TFD0SREOd6ll0JlJdx1F1xxhbU6RpMQx9sPlAW8\nLvVv65YHnA5U+x/UNhZYJyJzIw3JJKuHS39I2d+WMJLe+R8eXOz7Qv+hmE21jT1zQGZWlPD5iuI+\nr4eD4qzikCuHirOKbYhm+NIkRDmeiNUb8o1vwGOPWWXb8/LsjkrF2VZgiohMwko+LgcWdL9pjGkG\nev7aikg1UcwJAWtFQKht4VYPJKrtRlcVs75wX5+lpPv8S0lPCWr7fl1zn4RjZkUJKxdU8n5d87BJ\nQkJVUA1lMMlKMnwO4tE2FscMRyemxplOTE0OxsB558GBA1Y11XvvhQsusDsqFSiWE1MBRORrwANY\nS3T/aIy5R0SWAW8bY9YFta0miiREr2elohPt9aw9IQEi1RUQkRHA48A5gBuYb4zZk+g41eCJwN13\nw0UXWfNEdDjG+YwxLwIvBm27PUzbqkTEpJTqKzVr9MZBQF2Bi4GpwBUiMjWo2fXAEWPMZGAFcH9i\no1QnYs4cOP982L+/t2aIUkop++ituNcMoMYYs9sY4wWeAr4Z1OabQPdU4LXAbJEonhClkoII/Pzn\n0NEBzzxjdzRKKaU0CekVqq5A8CMle9oYYzqAZqDf7CQRWSgib4vI242NjXEKVw3Fl75krZY5JXim\nnlJKqYTTOSFxEFhXoK6uTmf+JpnnnrM7AqWUUqA9IYEi1RXo00ZEMoACiOLBDUoppZTqR5OQXj11\nBUTEhVVXIPgpR+uAa/w/zwM2Gl3jrJRSSg2JDsf4GWM6RGQR8DK9dQW2B9UVeBR4QkRqgCasREUp\npZRSQ6BJSIBIdQWMMW3AtxIdl1JKKeVEOhyjlFJKKVtoEqKUUkopW2gSopRSSilbaBKilFJKKVvo\nxNQ483q9jbt37/7E7jiGqqGhoWT06NGOK/uq55WUJtgdQCSpfD2n+GcjLKeeF6T8uUV1PYuWuVAD\nEZG3jTExe7x6stDzUsONUz8bTj0vcPa5ddPhGKWUUkrZQpMQpZRSStlCkxAVycN2BxAnel5quHHq\nZ8Op5wXOPjdA54QopZRSyibaE6KUUkopW2gSopRSSilbaBKi+hCRIhHZICK7/N9HhWhzlohsFpHt\nIvK+iMy3I9ZoiMhXRWSniNSIyJIQ748Qkaf9778lIhMTH+XgRXFePxKRHf7fz6sikvQ1OFRs6bUs\nExMf5eAN92tZkxAVbAnwqjFmCvCq/3WwVuBqY8xpwFeBB0SkMIExRkVE0oGHgIuBqcAVIjI1qNn1\nwBFjzGRgBXB/YqMcvCjPaxsw3RhzBrAW+GVio1RJQK/lJKfXsiYhqr9vAqv9P68GLgluYIz5yBiz\ny//zAeAwMDphEUZvBlBjjNltjPECT2GdX6DA810LzBYRSWCMQxHxvIwxrxljWv0v3wRKExyjsp9e\ny3otJz1NQlSwMcaYg/6fDwFjBmosIjMAF1Ab78CGYDywL+B1nX9byDbGmA6gGShOSHRDF815Bboe\neCmuEalkpNeyXstJT58dMwyJyCvA2BBvLQ18YYwxIhJ2DbeIjAOeAK4xxnTFNkoVCyJyFTAduMDu\nWFTs6bU8fDj1WtYkZBgyxswJ956I1IvIOGPMQf+N6XCYdvnAX4Glxpg34xTqidoPlAW8LvVvC9Wm\nTkQygALAnZjwhiya80JE5mD9MbrAGNOeoNhUAum1HLKNXsspRIdjVLB1wDX+n68B/hLcQERcwJ+B\nx40xaxMY22BtBaaIyCR/zJdjnV+gwPOdB2w0yV/BL+J5iUgl8HtgrjEm5B8f5Xh6Leu1nPyMMfql\nXz1fWGOorwK7gFeAIv/26cAj/p+vAnzAewFfZ9kde5jz+RrwEdY491L/tmVYFzRAFvAsUANsAU62\nO+YYndcrQH3A72ed3THrV8I/I3otJ0HcMTgvR1/LWrZdKaWUUrbQ4RillFJK2UKTEKWUUkrZQpMQ\npZRSStlCkxCllFJK2UKTEKWUUkrZQpMQpZRSStlCkxCllFJK2eL/AaIngtND0nj5AAAAAElFTkSu\nQmCC\n",
"text/plain": [
"<matplotlib.figure.Figure at 0x1140f9350>"
]
},
"metadata": {},
"output_type": "display_data"
}
],
"source": [
"q_l = 0.4\n",
"q_r = 0.4\n",
"v_l = 1.0\n",
"v_r = 0.1\n",
"t = 0.5\n",
"\n",
"states, speeds, reval = riemann_traffic_vc_exact(q_l,q_r,v_l,v_r)\n",
"riemann_tools.plot_riemann(states, speeds, reval, t);\n",
"x1, q1 = test_solver(q_l,q_r,v_l,v_r,traffic_fwave_old)\n",
"plt.plot(x1,q1,'-x')\n",
"x2, q2 = test_solver(q_l,q_r,v_l,v_r,traffic_fwave)\n",
"plt.plot(x2,q2,'-o')\n",
"x3, q3 = test_solver(q_l,q_r,v_l,v_r,traffic_qwave)\n",
"plt.plot(x3,q3,'-s')\n",
"plt.legend(['Exact','Old','New','q-wave'])"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"# q-wave method is more accurate"
]
},
{
"cell_type": "code",
"execution_count": 15,
"metadata": {
"collapsed": false,
"deletable": true,
"editable": true
},
"outputs": [
{
"data": {
"text/plain": [
"<matplotlib.legend.Legend at 0x1143cd050>"
]
},
"execution_count": 15,
"metadata": {},
"output_type": "execute_result"
},
{
"data": {
"image/png": "iVBORw0KGgoAAAANSUhEUgAAAtAAAAEWCAYAAABPDqCoAAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAIABJREFUeJzs3Xl81PW1+P/XyUoCCYGwBQib7JvsS1lkExFUEFAWK+KC\ntbW1+63t/V5rvbf32t5ata33Zw3WfUEWFVRcERSVVZRVdoSwE7aQfTLv3x9nIiFkQkImM5PkPB+P\neSTzmc+8P++Jkpw5c97nLc45jDHGGGOMMeUTEeoJGGOMMcYYU51YAG2MMcYYY0wFWABtjDHGGGNM\nBVgAbYwxxhhjTAVYAG2MMcYYY0wFWABtjDHGGGNMBVgAbYwxplYSkaUicluo53EpItJGRJyIRIV6\nLsYYZQG0McaYShOR34rI0hLHdvo5Nj24syudc+5a59xzVX0dEZktIiur+jrGmOCxANoYY0wgfAJ8\nT0QiAUQkBYgGepc41t53rjHGVFsWQBtjjAmEtWjA3Mt3fxjwMbC9xLHdzrlDACLyuIgcEJGzIrJe\nRIb5jjcXkRwRaVg0uIj0FpETIhLtu3+HiGwTkVMi8p6ItPYdFxF5VESO+cbdJCLdS5uwiCwXkbt8\n388WkZUi8hffmHtF5Fo/z7tCRE6KSJ9i8z0uIiNKObcL8CQwWETOicjpMubyPyKyxjfvN4u//hLn\n3u577ZkiskdEflDssREiki4iv/T9DA6LyO3FHo/1vcb9InJURJ4UkbjSrmOM8c8CaGOMMZXmnMsH\nVgPDfYeGA58CK0scK559XosG1w2Bl4H5IlLHF2B/AUwpdu5MYIFzrkBEJgK/AyYDjX3XecV33ljf\ndToC9YGbgYxyvoyBaMDfCPgz8LSISCmvdTfwG+BFEYkHngGec84tL+XcbcA9wBfOuXrOuaQyrj8L\nuANIATzA3/ycdwy4DkgEbgceLQrmfZqhr70FcCfwhIg08D32MPqz6YV+GtACeKCMORljSmEBtDHG\nmEBZwflgeRga2H5a4tiKopOdcy865zKccx7n3CNALNDJ9/DLwAzQrDIw3XcMNCD9H+fcNuecB/hv\noJcvC10AJACdAfGdc7ic8//WOZfmnCsEnkMD2aalneicSwN2oW8aUoB/L+c1yvKCc26zcy4L+A/g\n5qLylxLXfts5t9upFcD76M+2SAHwkHOuwDn3DnAO6OT7Od4N/Nw5d9I5l4n+7MKiJt2Y6sQCaGOM\nMYHyCTDUV3rQ2Dm3E/gcrY1uCHSnWAZaRH7lK0U44yttqI9mfwEWomUPKWgA7kWDcYDWwOMictr3\nvJOAAC2cc8uAfwBPAMdE5CkRSSzn/I8UfeOcy/Z9W6+M89N8r+nvzrk832sa5ivVOCciW8p53SIH\nin3/LVoS06jkSSJyrYis8pWRnAbGlzgvw/fGoki273U0BuKB9cV+du/6jhtjKsACaGOMMYHyBRoE\nzwE+A3DOnQUO+Y4dcs7tBQ00gX9DSywa+EobzqCBMM65U2hmdRpavvGqc875rnMA+IFzLqnYLc45\n97nvuX9zzvUFuqLlCr8O9AsVkXrAY8DTwINF9crOuU99pRr1nHPdfKc7f+OUkFrs+1ZoJvlEievG\nom8u/gI09f3c3sH3c7uEE0AO0K3Yz62+c66sNwnGmFJYAG2MMSYgnHM5wDrgF5zPFoPWQf+CC+uf\nE9A63+NAlIg8gNb0FvcyWhc8lfPlG6CL8n4rIt0ARKS+iNzk+76/iAz0LTbMAnLR7HWgPQ6sc87d\nBbztm5M/R4GWIhJziTG/LyJdfXXVD6E134UlzolBS12OAx7fQsex5Zmwc86LZs0fFZEmACLSQkSu\nKc/zjTHnWQBtjDEmkFYATdCgucinvmPFA+j30PKBHWi5Qi4XljAALAY6AEecc18XHXTOvQ78CXhV\nRM4Cm4GijhmJaJB4yjduBvC/gXhhRXyLGMcBP/Qd+gXQR0Ru8fOUZcAW4IiInPBzDsALwLNoKUkd\n4L6SJ/jqlu8DXkNf40z051Rev0Frt1f5fnYfcr7u3BhTTnL+EzFjjDHGhIKILAdedM7NDfVcjDGX\nZhloY4wxxhhjKsACaGOMMcYYYyrASjiMMcYYY4ypAMtAG2OMMcYYUwFRoZ6AMdXZvn37XFSU/TMy\nxhhjaoL8/PwT7dq1u+TmQvaX35hKiIqKomXLlqGehjHGGGMCYM+ePd+W5zwr4TA1ioj8S0SOichm\nP4+LiPxNRHaJyEYR6VPssdtEZKfvdlvwZm2MMcaY6sQCaFPTPItucODPtejGDB2Au4H/D8C3De/v\ngYHAAOD3ItKgSmdqjDHGmGrJAmhTozjnPgFOlnHKROB5p1YBSSKSAlwDfOCcO+mcOwV8QNmBuDGm\nlsjP15sxpuZyDtavL//5FkCb2qYFF24XnO475u/4RUTkbhFZJyLrTpwoa1feilmyBN5/HzyegA0Z\nlhYtgm++CfUsqpZzMG8e5OWFeiZVyzl4++1Qz6Lq7dunN2NMzXPmDPzf/0GfPtCvX/mfZwG0MRXk\nnHvKOdfPOdevUaNGARv3+uthwAANpBcuhADG5mFl8mTIzoZXXoGMjFDPpmqI6OucN6/m/ncEfZ0J\nCXDgwKXPNcaYcOEcrFwJs2dDSgrce68ef+KJ8o9hAbSpbQ4CqcXut/Qd83c8qJKS4MYb9bZpEyxY\nAF9/HexZVL0+fWD6dP247PXXa+bH49HRcOut8OmnsGNHqGdTdYYP19dojDHh7sQJ+OtfoVs3GDZM\nPxGdNQvWroUvv4Qf/aj8Y1kAbWqbxcAsXzeOQcAZ59xh4D1grIg08C0eHOs7FhIRETByJEydCnFx\nGki/9x4UFIRqRoEnAmPHwrhxsHgxLF+uWYGaRETfDB08CKtXh3o2VWfgwJr9+owx1ZfXCx9+CNOm\nQfPm8MtfQv368PTTcOgQPPmklm6IVGxc6wNtahQReQUYATQSkXS0s0Y0gHPuSeAdYDywC8gGbvc9\ndlJE/hNY6xvqIedcWYsRg6ZjR71lZsI772iN9LBh0KRJqGcWGHFx+kbh8GF49VXo2VOzAzXJyJH6\nScLSpXDttaGeTeBdcQWsWwd9+4LtK2SMCQeHDsEzz2igvHcvNGigGea77oLu3Ss/vrialvIxJojS\n09NdsDdScU4/Mj92DNq103KImuTrr2HrVhg9uua8SSiyfz989plmQiJq2Od/587BRx/BxImhnkng\nFZXgdOwY2nkYY8rm8cC770Jami5wLiyEESNgzhxdl1KnzqXH2LNnz/p27dpdcjmh5QqMqWZEtO4U\nYPduLe+Ij4cxYyAmJrRzC4Qrr9Qs9LJlcPo0XHcdxMaGelaB0aqVZkGee06D6Pj4UM8ocOrV0/9O\nGRmQnBzq2RhjapN9+zTT/MwzWjLXtCn86ldw553QoUPVXNMy0MZUQigy0KXJytLsX34+DBmiq4pr\ngrw8zSIkJmpGuqI1auHK49FylTFjoFmzUM8mcJzT1zVjRqhnEliWgTYm/OTn6/qZtDT44AM9Nm6c\nZpuvu04Xcl+O8magLYA2phLCJYAu4hx8/rnWE7dufXkLI8LRsWOake7SRTPUNcXixdC+PXTtGuqZ\nBM6mTfq1R4/QziOQLIA2Jnxs3w5z5+onecePQ2oq3HGH3lq1qvz45Q2ga1gVnjG1m4hmoKdO1frh\nRYs0g5ubG+qZVU6TJtr2LipK+0cfPhzqGQXGDTfAyZNaF11T9OihQbTlZowxgZKTAy+8oOWLnTvD\nY4/B0KG6sH7vXnjwwcAEzxVhGWhjKiHcMtClycnRFj65uTB4MIT5dC/JOVixQvt5TpigXTyquy1b\nYM8e3UynJjhxQnt8X3NNqGcSGJaBNiY0Nm7UEo0XX9Q1MVdcoV00Zs+uuvI3W0RojAE0wLz+eg08\nV63SW4sWMGhQ9SzvENFV1fn5mn2Ii9N+0tXxtRTp1g0aNtQ/EjNmQGRkqGdUOY0a6Ru2c+d0caEx\nxpRXZqaupZg7F9as0cXxU6ZobfNVV4VPByPLQBtTCdUhA12a9HQNpGNidCFbde4GkZGhC0g6dqz+\nLf2ysuC117QEJyEh1LOpHI9Ht6SfNi3UM6k8y0AbU7Wc090A09K0TC8rSxMLc+bA978f3M4+tojQ\nmCCorgF0kdxc7d6RkwP9++vCw+pq+3bdinXoUF1UUl15vTBvnm6WU43/1wJ0d8JGjfRj1+rMAmhj\nqsapU/rJW1qarp2Ij9c33XPmhO5TUltEaIy5pDp1tI54yhTtdLFwIaxcWT0XgHXqpOUP+/ZpFjcr\nK9QzujwREfo6Nm7UW3VmW3wbY0oqWsdy6626tfZ99+mnoU8+qQvE//UvXa8T7mV5loE2phKqewa6\nNIcPayu8qCgt76hbN9QzqjiPR+ujo6K0L2i41MxV1KpV2gv7qqtCPZPLd+CArpIv2vynOrIMtDGV\nd+yYtp6bO1f/TSUmwi23aLa5d+9Qz+48W0RojLksKSmakc7P1/KOrCz95VadPoaPitIWcadPw/z5\n0LYtDBgQ6llV3KBBWpry+uswaVL4Z2RKk5qqbwTy82vGTpnGmPLzenWNSloavPmmJjeGDIHf/Q5u\nuql6r7+xDLQxlVATM9Cl+fJLbbPWuLHW5la3jO6uXbpAZfBgaNMm1LOpuOPHYelSLe243N21Qik3\nVz8RmDw51DO5PJaBNqZi0tN1W+2nn4Zvv9VFgLfdpi3ounQJ9ezKZhloY0zA9Omjt2PHNBsaGalb\na1eXThHt2+vtiy+0Jvfaa/Xjw+qicWO4+WZ46SXNRCclhXpGFVOnjs758OGas828MeZCHo9u3JWW\npm/4vV79O/GnP+nvrdjYUM8wsCwDbUwl1JYMdEkFBbq1dmYm9OxZvTJzhYX6yx00kK5OPZed05KU\n/v21LKU6cU7bU82cGeqZVJxloI3xb88ezTQ/88z5N8m33w533gnt2oV6dhVnbeyMCYLaGkAX9/XX\nsHOnbgQyYkT1Ke/IzNRAOjVVSzuqk/ff1/Zw1a3v9fbturlK376hnknFWABtzIXy8uCNNzTb/NFH\n+nv/2mt1QeCECboOpbqyNnbGmKC48krd+OPKK/UX6qJFungv3CUkaFlESopmRnftCvWMym/sWP14\n9KOPQj2TiunUSYNRrzfUMzHGXI5t2+AXv9DdbKdP19+bDz2kdc5vvQUTJ1bv4LkiLANtahQRGQc8\nDkQCc51zD5d4/FFgpO9uPNDEOZfke6wQ2OR7bL9z7oZLXc8y0BfzeODjj+HMGejaVW/Vwdq1sHu3\ntr2rLjXGu3frAs+pU6tPh44zZ7TX+IQJoZ5J+VkG2tRm2dlaOpaWBp99pgHypEmabR4zpvp86lhe\ntojQ1DoiEgk8AVwNpANrRWSxc25r0TnOuZ8XO/8nQPHukznOuV7Bmm9NFRUFV1+t32/ZAgsWaEA6\ncmR41xv376+lBe+9py3XqsPHkFdcAQ0aaG/VGTOqxyKd+vU12D99uvq8UTGmNtqwQYPml16Cs2f1\nDeSf/6zdNJo0CfXsQi/M/zwYUyEDgF3OuT0AIvIqMBHY6uf8GcDvgzS3WqlbN72dOgWLF+tCshEj\ntF46HBXV8WVlabeRZs20bV84a9hQF+a98gpcd522iwp348bpduUzZoR6JsaY4s6e1d8laWmwfr12\n0Jk6VbPNw4ZVn0+6gqGGJd5NLdcCOFDsfrrv2EVEpDXQFlhW7HAdEVknIqtEZJK/i4jI3b7z1p04\ncSIQ867xGjSAG2/U+rgNGzQrvWnTpZ8XKnXrapP/Nm30j8k334R6RmWLiYFZs2D58upRyx0RofXQ\n27aFeibGGOe0xecdd+iakHvu0U/h/vY3OHQIXnhBdxK14PlCFkCb2mo6sMA5V1jsWGvnXD9gJvCY\niJS6955z7innXD/nXL9GjRoFY641RlH/6KlTdUOQBQu0ZMLjCfXMSpeaqlnS7GwNpDMyQj0j/0R0\nB8n9+2HNmlDP5tL69NE3U7YMx5jQyMiAxx6DHj3ge9+D117TT7NWr9buSj/5iSY/TOksgDY1yUEg\ntdj9lr5jpZkOvFL8gHPuoO/rHmA5F9ZHmwDr3FkD6cGDdfX2woW641446tNHV5yvX6+lHfn5oZ6R\nf6NGaUb63XdDPZNLGz1a+4kbY4LD69VF3jNnaieNn/9cP3FLS9MezmlpMGCAZZvLw2qgTU2yFugg\nIm3RwHk6mk2+gIh0BhoAXxQ71gDIds7liUgjYAjw56DMupZLTNQV3V4vfPIJnDihi+N6h9nbFxFt\nH5eTo/Xcyclazx2Of2h69dK2Uq++qq36wnWVfNOm8PnnmuGPjw/1bIypuY4cgWefhblztXtPUpLW\nNc+Zo5thmYoL01+rxlScc84D/Bh4D9gGvOac2yIiD4lI8ZZ004FX3YU9HLsA60Tka+Bj4OHi3TtM\n1YuI0IB06lTt0bxwoWZRCwpCPbMLxcXpHDt31gB1y5ZQz6h0rVvD+PHaoSMnJ9Sz8e+663T7X2NM\nYBUWwjvvwOTJWo72299q1vmFF7S2+e9/t+C5MqwPtDGVYH2gq9a5c7pZSEEBDB2qXTHCzcaNGkSP\nHh2erZ08Hq3fHjtWM77haOVK/QPfunWoZ1I66wNtqpP9++Ff/9LbgQPQuLG2nrvrLl28a8pmW3kb\nEwQWQAeHcxpkHT2qnTH6XfJXW3A5p7W8p09rRjXc+jE7p2UnnTpp5jwcvfyy1mWGIwugTbgrKIAl\nS7REo2j9w9VXa4nGDTfoughTPhZAGxMEFkAH3969upgvLk53wQqnYDUvT8sREhM1Ix1u9dGffqqd\nUL73vVDP5GJ792qd5uDBoZ7JxSyANuFq1y4Nmp99VhMMLVpoO7o77tBkg6k424nQGFMjtW2rt+zs\n87sGDh6sfzhCLTZW6w2PHdONQrp0gSuvDPWszhs2DDZv1iA/3LbSbttWt1Pv109bHBpjSpebC4sW\naceM5cv1TfGECZptHjcu/HdQrSksA21MJVgGOvSKNgE4fBhatgyvFkxbt2o/1auugubNQz2b8w4d\n0lZW06eH1/bq2dnw/vvalSWcWAbahIMtWzRofuEFOHlS33TedRfMnh1ev1+qO8tAG2NqBZHzJQn7\n92v3jjp1tIQiLi60c+vaVbPQK1Zou7YJE0I/J9A/tpMmwfPP646L9eqFekYqPl5vx4/rwidjarus\nLP00Ky0NVq3ST2cmT9bAedSo8G1RWRtYAG2MqTFatdJbTo5278jJgYED9VioiGh7voICLZ2Ii9OO\nGKHOktetq9t/z5unGfJwKIEBXfj0yivhu6DQmKrmnK7zSEvTfwuZmbr495FH9N+sbYAbHqyEw5hK\nsBKO8Oacbmt94ACkpGimOtSB68mTWqbQoQP07RvauRR5+219k9GjR6hnorZs0fZ74VI/biUcJhhO\nn9ZuNGlp8NVX+mb75ps12zxkSOh/d9UW5S3hsOS/MabGEtEM9NSpWi+4cKG2c8vODt2cGjbU2uN6\n9TS7dOBA6OZSZMIEzXJ98kmoZ6K6ddP6ca831DMxpmoVtei87TYtrbr3Xj3+xBO6VuHZZ7UHvgXP\n4ccy0MZUgmWgq5+8PC3vyM7WDHDbtqGdz8qV+ody/PjQ1yJ/843eJk4M/R/skydh9Wq49trQzgMs\nA20C78QJXYMwdy5s26a7r86cqZ00+vQJ/b+/2swWERpjTCliYzVYBVi3TmsNmzYNXZZn6FAtV3jn\nHW0/NW5c6BYGde4MDRroKv+ZM0PbDqthQ60bz8zU4MKY6s7r1Q2X0tLg9df1/+/Bg+Hpp7VUI9Rv\noE3FWAbamEqwDHTNcOSIZoKjonRzllD9ITt9Wntbt22r7fhCJScHXn1VV/vXrx+6eRQWwvz5WvIS\nSpaBNpVx6BA884wGynv36pvUWbO0trl791DPzpRkOxEaEwQWQNcs+fmaITp3ThewdegQmnns3q2L\nHwcNCl2JiderweugQdC6dWjmALq5SlJS6P5bgAXQpuI8Ht1SOy1NF+kWFsLIkVqiceON2mrThKeA\nLCIUkSQR+VEgJiQisSIyT0R2ichqEWnj57x9IrJJRL4SkXXFjjcUkQ9EZKfva4NAzMsYY4rExGgJ\nxdSp2n91wQLdcCTYi9muuAJmzNDM+Lx5cPZscK8PWkYybZrWRG/YEPzrF+nfX4Noy/WY6mDfPviP\n/9BttK+/Xuv4f/UrfRO2bJn+u7bguWYoMwPtC3Lfcs5V+kMGXyDe0zl3j4hMB250zk0r5bx9QD/n\n3IkSx/8MnHTOPSwi9wMNnHO/qeAcGjjnTl3+qzDmQpaBrvmOH4dPP9WActQoSEwM7vULCzWT5Zwu\nqAvFzoFr12ot8qhRwb82wMGDsHOn9tMOBctAm7Lk52t3n7Q0+OADPTZunGabr7vOtqavbgLVxu5h\n4ApfNvh/KzmnicBzvu8XAKNFKrRkp/jznwMuZ7PXdSLykoiMquC1jTG1VOPGWgt83XW6E9iCBZqV\nDZbISG0zd9VV2obviy+Cd+0i/ftrn+gFC0KTCW7RQrsW5OUF/9rG+LN9O/z619Cype7ouW0bPPCA\nZqHfeUdLNSx4rrkqlYEWkU+B0tZH/8o592GJczcD45xz6b77u4GBpWSa9wKnAAf80zn3lO/4aedc\nku97AU4V3S8vEYkErgXuALoALwDPOucOVWQcY4pYBrp22rRJ/3g2aKBZ0WBmhfft0yC6X7/g1wVn\nZMBbb+nH0DExwb12Xh4sWaLlNcFmGWhTJCdH30impeknU1FRWqoxZ47uMBqKT4hMYAVlIxXn3DDn\nXK9Sbh9e+tl+DXXO9UED3XtFZHgp13VogF3R+RY6595yzk0GhgPtgP0iEsL17ibQRGSciGz31dvf\nX8rjs0XkuO+Tla9E5K5ij93mq7PfKSK3BXfmprro0UMDud694c03YdEiOBWk4rA2bTSAPX1aO2UE\n67oAycl67Vde0T7NwRQbq1sYHzwY3OsaA7BxI/zkJ7rZyaxZcPgwPPywboS0aFHoyqtM6FSqy2dF\nMtDAQSAVSBeRKKA+kFHyic65g76vx0TkdWAA8AlwVERSnHOHRSQFOFbKfJ4BegOHgB8AS3wPPemc\ne9J3Tn1gOjAbyEez0Rsr8rpN+PJ9yvAEcDWQDqwVkcXOua0lTp3nnPtxiec2BH4P9EPfoK33Pdfq\n5k2pGjbU8o7CQl1sePo0dOmiO+lVtf79dSOY997TGswJE4LTtzkmRgOIhQv1DcQVV1T9NYuMGKFb\nHc+cGbxrmtorM1PfpM6dq11xYmJgyhTNNl91Vej6tZvwcKlft5mUHiADmoGuwLUWA7cBXwBTgWWu\nRP2IiNQFIpxzmb7vxwIPlXj+w76vb5Yyn9tLHOpVYvwXgcHAfGCWc25nBeZvgkhE/lRykWhpx0ox\nANjlnNvje86raP18yQC6NNcAHzjnTvqe+wEwDnjF3xPWr19PampqOYY2xgTKLbeEegamNsrP109f\nXvH7F8HUBLt37y7XeWW+f3LOZQCficjmACwifBpIFpFdwC+A+wFEpLmIvOM7pymwUkS+BtYAbzvn\n3vU99jBwtYjsBMb47lfUa0An59z9FjyHvatLOVaeTX1bAAeK3U/3HStpiohsFJEFIlIUAZfruSJy\nt4isK95m0RhjjDG1xyU/8HPOBeTDMudcLnBTKccPAeN93+8BrvTz/AxgdCXnsLgyzzdVT0R+CPwI\naCcixUtrEoDPAnSZJcArzrk8EfkB2tWl3A26fAtbnwJdRGibEZmyeL2wfLnWDHfooBu0VLUNG3SR\n45gxWjdc1fbt0363N98cnO3QMzO1ZOaGG6r+WmCLCGsq5+CTT3RB4IIFulC1b18t0ZgxI/gtK014\n2LNnT7nOC0LFnDEV8jKwFPgffJ9S+GQWlVZcQlGtfZGWvmPf8b0ZKzIX+HOx544o8dzl5Zm0Mf4U\n9Y8GDcQWLNCtwkePrroWV717Q69e2pM2K0vro6uya0abNtqR5LnndNvtqt4oIiFBf3YnT2odujEV\nceyY/r86d67+m0xMhDvv1K21e/cO9exMdWFbeZsaxbdAdQf6acVBYC0w0zm3pdg5Kc65w77vbwR+\n45wb5FtEuB7o4zv1S6BvWYG7tbEzlyMzEz76SLf7HT4cmjSpumvl5GhP2oYNdRFeVWaICwq0PnTc\nuKp9TaDZw1df1UxhVbMMdPXn9eobyrQ07Zzj8cCQIZptvukmiI8P9QxNuChvGzvLQJsaxTnnEZEf\nA+8BkcC/nHNbROQhYJ2vjOc+EbkB8AAn0Y4sOOdOish/okE3wEPlzHobUyEJCTBpkv5RX7lSP0Zu\n1w769Ln0cysqLk47Bxw+rAFnz55V1yUkOhpuvVUDlM6d9VZVRKBrV9i8GbpXeq9cU1Olp8Mzz8DT\nT8O332orxvvu02xzly6hnp2pziwDbUwlWAbaBMru3Vq7XLeulndUVcnFxo0adI4eDU2bVs01AFas\n0NcweHDVXQO0rd2MGVWbWbcMdPXi8cDbb2u2eelSfaM6ZowGzZMmaU9xY/wpbwbaAmhjKsECaBNo\nWVla3pGfrx8xp6QE/hrOne9bPWFC1QUUmzbpRhPjx1fN+ADHj+sbj7Fjq+4aFkBXD3v2aF3zs8/q\nJy4pKXD77Vrf3K5dqGdnqgsr4TDGmGqobl3tLuEcfPaZ3lq10o1TApVlFdGFjXl5Wh+dkKAZ6UBn\ncXv00MWFL78M06ZVzU5tjRtDdra+8ahbN/Djm/CWlwdvvKHZ5o8+0kW748drbfP48cHZXMjUTpaB\nNqYSLANtguHbb2HtWu1uMWZM4LtcHDsGy5ZpTWhVtNk7dw7mz9c2d1UR5Ho8ujPitGmBHxssAx2O\ntm3ToPn55yEjA1q31kzz7beD/Uo2lWEZaGOMqSFat9ZbTo52EsjN1driQAUKTZpo+7mtW7WLxlVX\nQfPmgRkbtG3frFm6iHHkyMCODZplbNVKP8K3j+prruxsfSOWlqafzERHw8SJmm0eM8a21jbBZRlo\nYyrBMtAmFJzTjUvS06FFCxg0KHDlF0WbSxw7pvXRgW7vtWQJtG1bNZ0zXn4ZZgZk668LWQY6tDZs\n0KD5pZdQ/cYHAAAgAElEQVTg7Fn97zBnjr4pq+p2iab2sUWExgSBBdAm1NLTNZiOidEsXFxcYMYt\nKNBOBnFxukAvkPXRn32mnRGGDQvcmKClLgcOwNChgR3XAujgO3tWPw1JS4P167VsaepUDZyHDQvO\njpemdrISDmOMqQVattRbbi58+KGWeQwYoCUflREdrS2/Tp6EefN0G/K+fQMz5yFDtIb1zTd1wWSg\ngqHWrfXNRH5+1e68aKqGc7BqlQbN8+ZpyUbPnvD3v8Mtt+iCVGPChWWgjakEy0CbcOOcLjjcvx+a\nNdNgNRAB6vbt8OWXOl6rVpUfD+DIEa3pnjEjcN0ScnK09+/kyYEZDywDXdUyMuCFF7QF3ZYtutB0\nxgzNNgey+4wx5VHeDLSV3BtjTA0iohnoqVPhiiu0O8Wbb2qbt8ro1EmDmv374bXXtLNGZTVrpvN8\n4QX9yD4Q4uIgMVGDcxO+vF7t/DJzpi4q/fnPNXBOS9Mezmlp+v+xBc8mXFkG2phKsAy0qQ7y87VH\nblYW9O6tgXVleDzaPzoyEq69tvLdD7xeDcq/973AZLed0/rZQC0otAx04Bw5ohudzJ2ru28mJen2\n73fdpeUaxoSaLSI0JggsgDbVzZdfwt69ugHJ0KGVC35Pn4b33tOuGgMGVH5u776ru8cFohf1tm1a\nztGnT+XHsgC6cgoL9f+TuXO1C4vHA8OHa4nGlCmBW/hqTCBYCYcxxpiL9OmjQUvnzvD667qLW2bm\n5Y2VlKSblyQna4/nvXsrN7dx43Qx5McfV24c0E1htm/X7LYJjf374cEH9Q3WhAmwcqWWanzzDaxY\nAd//vgXPpvqyDLQxlWAZaFPdFRRoece5c/oRemWyrF98oa3kxo/XOuTLtXMnbNyoCwErUwN7+rS2\nzJsw4fLHAMtAV0RBgWaZ587VTxQArr5as8033GDdUUz4szZ2xhhjLik6WjO/AF9/DQsWQMOGMGJE\nxcs7Bg/WUo5339XM7/jxWiddUR06aMuy55/XhYuXG3QlJWk99JkzUL/+5Y1hymfXLg2an30Wjh7V\nDX7+3/+DO+6ANm1CPTtjAs8y0MZUgmWgTU2UkaEfsQOMGqWBaEVlZmo7uZYtdXHg5cjL08WAEyde\nfg/gogWK06df3vPBMtD+5ObCokXaMWP5cn2zdN11uiBw3LjAtSY0JpisBtrUSiIyTkS2i8guEbm/\nlMd/ISJbRWSjiHwkIq2LPVYoIl/5bouDO3NjwkdyspZP3HCD9pResEAX5VVEQgLcfLO2KHvlFS3L\nqKjYWLjtNt0gZs+eij8fNIveoYPW3ZrA2LwZfvYzzTLfcouW7fzxj1rz/MYbGkRb8GxqOstAmxpD\nRCKBHcDVQDqwFpjhnNta7JyRwGrnXLaI/BAY4Zyb5nvsnHOuXkWuaRloU1ts2aJBdFISjBxZ8dKM\ntWu1bdk111xeNvnDD/Xa/S6ZFyrdSy9pW7vLqam2DLS2QJw3T7PNq1ZpWc2NN2pt88iRlW9laEy4\nsAy0qY0GALucc3ucc/nAq8DE4ic45z52zmX77q4CLPo1phy6ddNNT/r21Y1ZFi3Sbb7Lq39/zUiv\nXq1ZyoKCil1/zBgN0j74oGLPKzJqVGC6e9QmzsG6dfCDH2h7wTvv1IWZjzwCBw9q55XRoy14NrWT\n/W9vapIWwIFi99N9x/y5E1ha7H4dEVknIqtEZJK/J4nI3b7z1p04caJyMzammmnQQMs7Jk6EDRu0\nvGPTpvI9NyJCa2OvvlqD6E8+qdi1+/SB9u21prmiH56mpMCpU9ob2pTt9Gn4v//Tn3f//rpT5OTJ\n2oZu61b4xS+gUaNQz9KY0LIA2tRKIvJ9oB/wv8UOt3bO9QNmAo+JSKn7tTnnnnLO9XPO9Wtkf0VM\nLRUZqdnHqVO1k8eCBbpZhsdz6efWrQs33QTt2mkWsyL1yW3bwtix2u0hN7dic77+enjrrYo9p7Zw\nTgPk227TuvV779XjTzwBhw7pz3vIENta25giFkCbmuQgkFrsfkvfsQuIyBjg34EbnHN5Rcedcwd9\nX/cAy4HeVTlZY2qKzp01kB40SAPUBQvg+PFLP69lS+2OkZOjCw3L+4FOUpJuwjFvXvmuUyQmBpo1\ngwMHLn1ubXH8OPz1r1qiM2yYbq4za5aWbmzYAD/60eV1YTGmprMA2tQka4EOItJWRGKA6cAF3TRE\npDfwTzR4PlbseAMRifV93wgYAmzFGFNu9evDpEn6cf+WLbBwIXz11aWf17u3BtJffqkBXH7+pZ8T\nHa2B3sqV5xf5lcewYfDpp+U/vybyenVR5rRp2knjl7/UIPlf/9Js85NPaq27McY/azRjagznnEdE\nfgy8B0QC/3LObRGRh4B1zrnFaMlGPWC+6GeR+51zNwBdgH+KiBd9Y/lw8e4dxpjyi4jQjVhA29ct\nWKBlG2PGaOBbGhEtzcjJ0Z3sGjTQ7g5llQyIaCeI5ct1QeOgQeWb38CB2kmivOfXFIcOwTPPwNNP\n67brDRtqhvmuu6B791DPzpjqxdrYGVMJ1sbOmPI5d063DC8o0Cxw06Zln3/kiHbN6NGjfMHd119r\nZ4jx48s3n3nzYMqU8vUrrs5t7Dwe3RkyLQ3efhsKC/WNyZw5+uajTp1Qz9CY8FLeNnYWQBtTCRZA\nG1MxRYvVjh7VLZ4v1dd540YtBxk16tJB94EDWp4xffqlW6sVBfQTJ5Z9HlTPAHrfPs00P/OMvrFo\n2hRmz9ZWdB06hHp2xoSv8gbQVsJhjDEmaEQ0Aw1aRrBgAcTFaXlHbOzF5/fsqVnojz/WwHvCBP9Z\n09RU7bTx3HNa3xsf738e9erpOCdO1JyWbPn5sHixZpuL+mWPGwd//7vuDuivfMYYU3GWgTamEiwD\nbUzlZWfrorb8fPje97SNWmny8rQMISFBA25/9dEez/lNPlJS/F/XOe3+MXNm2fML9wz09u0wd66+\ncTh+XN9I3Hkn3H47tGoV6tkZU71YCYcxQWABtDGB4xx88YUudktNhQEDSg+Sjx/X8ovOnaFXL//j\nLVkCV1wBXbv6P2fTJr1uz57+zwnHADonR7P3aWlathIVpdn3OXN0MWZFt1o3xijbytsYY0y1IqIZ\n6KlTNXO8cKEGwSU3TGncWOucY2Lg5Zc14C7N9ddrd46VK/1fs0cP2Ly54jsbhsrGjfCTn2iWftYs\nOHwYHn5Y678XLYJrr7Xg2ZhgsAy0MZVgGWhjqlZOjmabc3K0/VzJkgTndEvw48e1A0dpdc9btsDu\n3XDDDaVfIyMD1q7VeuHShDoDnZmpJSlpaTrP2FjtIHLXXXDVVZdeMGmMKT9bRGiMMabai4vTBXDO\nwZo1ektJ0Uy1iN6uukrb4739ti4MvOaaC0s/unXTnscvvKD1ziUztMnJWl997pwuLgwHzmmwnJam\nddpZWfo6HntMd2FMTg71DI2p3SwDbUwlWAbamOA7eFBrpaOj4eqrL8w6nzwJ77+vrdpK7qaXna39\nn6dO1YWIxXk8WjIybdrF1wtmBvrUKXjxRQ2cN23S1zZ9utY2DxxY9sYyxpjKs0WExgSBBdDGhE5e\nnpZ3ZGdrsNy27fnHduyA9ethyJALyz68Xg2ihw7VhYrFrV6tmd327S88XtUBdFEZSlqaLgzMy9PX\nM2cOzJgBiYlVc11jzMWshMMYY0yNFht7fufBdes0YG7aVIPjjh31tnKlbts9fryWZ0REaFD6zjua\nrb7yyvPjDRyoixJLBtBV5dgxbT03d64G6YmJ2n7urrugd+/gzMEYc3ksA21MJVgG2pjwcuQIfPaZ\n1jmPGaNBs8ejAXNkpHapKFp0t2qVdvgYMeL88w8c0A1ehg8/fyyQGWivVzc5SUuDN9/UuQ0dqkHz\nTTeVvfmLMabqWQbaGGNMrdOsmXaoyM+HZct0YWCvXtqB48wZmD9ftxAfOBAGDdLgeNEiuPFGrS9O\nTdXAOj9f2+QFSnq6bqv99NPw7bdaKnLffRo4d+kSuOsYY4LDAmhjjDE1TkzM+bZ0GzZobXGjRprl\n3btXO1sMHKhZ5QYN4PnntUNHdLT2j16yRAPxyvB4tDNIWhosXarZ5zFj4M9/hokTS9+63BhTPVgA\nbYwxpkbr3Vtvx4/DG29oCceECbB1q7aKu/Za7b7x4ouaiU5K0rZ3hw7531a8LHv2aF3zs8/qRicp\nKXD//Vrf3K5dwF+eMSYELIA2xhhTKzRuDJMna2Z42TI4exa6d9etsL1euPVWLefo31/rol95RbPS\n5ZGXp8F5Wpp2BomI0IWLc+bo1yj7a2tMjWL/pI0xxtQqUVEwdqx+v2mT7nIYG6tlHqmpsHOndujo\n21e7e5TVRm7bNg2an39edzRs3Rr+8z9h9myw9cXG1FwWQJsaR0TGAY8DkcBc59zDJR6PBZ4H+gIZ\nwDTn3D7fY78F7gQKgfucc+8FcerGmCDr0UNvJ0/C8uUaUEdF6eLDxEQ4cUJb3RXfvTA7WxcjpqVp\nx4/oaK1pnjNHa5xta21jaj4LoE2NIiKRwBPA1UA6sFZEFjvnthY77U7glHOuvYhMB/4ETBORrsB0\noBvQHPhQRDo65wqD+yqMMcHWsKGWdxQWwscfa0/pU6d0E5YVK2DUKF2MmJYGL72k5R8dO8L//i/M\nmgVNmoT6FRhjgskCaFPTDAB2Oef2AIjIq8BEoHgAPRF40Pf9AuAfIiK+46865/KAvSKyyzfeF0Ga\nuzEmxIr6R48ZA5s3a+u5BQvgr3/V0o46dXQr8DlzYNgw21rbmNrKPmgyNU0L4ECx++m+Y6We45zz\nAGeA5HI+FxG5W0TWici6EydOBHDqxphw0r07PPKIlmfk5cHf/66dOV54QTdaseDZmNrLMtDGVJBz\n7ingKdCdCEM8HWNMFXvkEXj8cQuYjTHnWQba1DQHgdRi91v6jpV6johEAfXRxYTlea4xppaJjrbg\n2RhzIQugTU2zFuggIm1FJAZdFLi4xDmLgdt8308FljnnnO/4dBGJFZG2QAdgTZDmbYwxxphqwko4\nTI3inPOIyI+B99A2dv9yzm0RkYeAdc65xcDTwAu+RYIn0SAb33mvoQsOPcC91oHDGGOMMSWJJt6M\nMZdj3759Lsq2GDPGGGNqhPz8fNeuXbtLVmjYX35jKiEqKoqWtt2YMcYYUyPs2bPny/KcZzXQxhhj\njDHGVIAF0MYYY4wxxlSABdDGGGOMMcZUgAXQxhhjjDHGVIAtIjQmwAa8MIwc7+mLjsdFJLHm1k/D\nduwR80aQkZtx0fHkOsksn7Y8bMc2xhhjgs0CaGMCrLQAt6zj4TJ2aQFuWcfDZWwLzo0xxgSbBdDG\nVMKeozv45SuzAXD49vpt4v/8yX+agNfJd+c7JzgiwPnuf3dMv+IiLjhOa/9jD/7dj3HeCLzeKAq9\nUXhdJF5vFF5vJIXeKH1+UdWWEyiar68VfNRA/2N3nfMPnYsT8EaCi9AxXKTvawS4CMT31fm+Ft2i\nrvM/dpdRG3GeGPDEQGEM4o0BTywURutXF+n/yUDU/f6D865dy3xqmSq7dXNCAlx7LUyeDN2721bQ\nxhhTk1gAbUwleCILONHwUImj/gO+M432lHvs0rc48v9PNqb9sgueV/pXwQFepMRjQk4Zc4kc8s/v\n4u1Ai7rtljIfd4WReAticJ4YXIHevJ4YXEEszhNT5i+xw8mvU5hVX2/n6n/3vfPElHt+nR4fQXT9\ni4P0gjPJbP/p8jKfu3o1PPggtG+vgfSNN8KAARBhq0+MMaZas50IjamE3bvTnXMXbqRy42c9/J7/\nWv/VOK8X5/XgdV68hYU4byFerwecV793egyv13dc7zvn5QcH7/A79mMxv8IV5kNhPnjPf5Wir75b\nhLcA8RYQ4c0jwltAhCsg0hVwR/P9fsfesCedSPHiBbxo4F0oRcE4FIoG4V6gACGLOuRIHXIklpyI\nOH7UosDv2HfH3oKrUw9i6yKx8bhIocCbT4HLp8CbR4G3gAKX5ztWoMdcPgXefDzefL4++7nfsf2J\njYijXmQidaPqUy+qPglRSdSLSqReZH09FqnH60bW5z93+P+Zv9pvk9/Hvv0W3nhDb8ePnz/evLkG\n0pMnw/DhYBtZGmNM+NizZ8/6du3a9bvUeRZAG1MJ6enpruROhD2e8x9Ab7rNf8BVHqEae+OtX5OX\nl0POubPkZJ0lL+s0edmZFORk4sk5iyfnHN48vZF/DsnPIqLgHJGebKI82Xy/2WH/8957YeCe7WI5\nHVGfrMj6ZEclkR/bkMI6DXHxyUTWa0x0QiPikppSt2EzEhs2Y+iSUX7HfnfKu5zJO6O3/DOczTt7\nwf3vvi92v8DrP9gv6aHvPUTLhJa0rNeSJvFNiIy4+NOHwkL47DNYtAhefx32F3u5yclwww0aTI8Z\nA3XqlPvSxhhjqkB5A2jLfRgTYHERSX47ZYTz2EkR8Zz2Zpd6XCIiqBNXlzpxdWnQOKXCYyc/052M\niItrQJILHd9cO5+cM8coOHuMwqwTSHYGkbknic07RXzBKZrk7qP+6bPES17pg7dt5fe6J9evJrFZ\nO9q06Ep8vfqXnKdzjhxPDmfzzwfad75/p9/zH/j8ge++j4qIokW9FrSs1/K7oLplgt76DGrJ8OH1\nePRR+PJLDaYXLgR+NIK19TNYmwG/nXd+3IaxyayYvvyS8zXGGBMaFkAbE2CVbScXqrFvafUyo/KX\n03Hzo3AmHeq3ZEf3n7MsZkSlx36q5c9J/ex+4sj/7lgOMRwY8jAdB44t1xg5WZmcyThC5skj5Jw+\nRv7ZY3gyj5Oc+yoZURdnfpM9hfT4+Pbv7p8igYzIJmTGNiW/bnNc/ZZEN2xF3SZtSG5+BclNU4mI\njCQ+Op746Hia1W12yTktnbyU9HPppGf6br7vN2ds5kzemQvOTYpNOh9UT2nJH29ryYNflL4A8mRe\nBtdfr5np66+HRo3K9SMyxhgTJFbCYUwllFbCYS725IrdfoPze666onKDP9odzhy46LCnbjN2Dv8b\nWcf2UXByPxFn06mTfZjEvKM0KjxKgly4bDLfRXIiohGnopuQHZeCp14LIpJSuSv7Wb+XLqts5mz+\nWQ5mHiw1wD507hAe5ynzZW2Zsw5XEEtEBFx1lQbTkyaB/e9mjDFVx2qgjQkCC6BDb8cHT/vPbl/t\nv/zi7OkMMg7t4eyRveSe+Bbv6QNEnztIfM4RGhQcpbHLIEq8jEhtUWqGu2Ghl7+1+T2N2nSjWatO\nRFZgNaDH6+FY9jGuWXiN33MiXBSRp9tzfHNXsvZ0I2dvN/LSOzCgb8x3HT06dCj3JY0xxpSDBdDG\nBIEF0KFXVdltT0E+GUcP0GRuX1/zP//yXDSHI1M4GdeavKQriGrckfqtutGsXQ8Sk5L9Pq+shZtz\nesxhS8YWtpzYypl8rXt3nihy0zuQs7cbOfu60iKyG5OGdmDq5Gh69ryw17RtMGOMMRVnAbQxQWAB\ndC3gp0SksF5zdg57nLPpW/Ee30Hs2b00ytlHivcIUeL97rwTJHE0JpVz9driktsTl9Llu6x1r5d6\n+71sUXmIc45DWYfYcmILXx/dysqdW/g2dwuF0ZkAeAuiyT3QkeiT3eid0o0bv9eNiUPb0e+lPpcc\n2xhjzIUsgDYmCCyArvkqWiJSkJ/H4X3byNi3mbwj24k4uYvErH00LThAAzK/Oy/fRTG2VdPSF0B6\nHctv3+x3Ts459p5KZ/4nW/h461b2Zm8lssVWIuN1fFcQg0Tn+32+BdDGGFM6C6CNCQILoGu+QJaI\nnDp+mKN7NpF5cBue4zsYdPjFUjd4dMCa5InQrCdJbfvQqkt/4uom+B23sBA++9zLS28f4KPNW8mq\nt4VG1z7n9/wPpn5Qri4jxhhT21gAbUwQWABtKsVPeUg+0eQSQyJZABQ6IT2yBcfrdsLTuBt1W/em\neecBJDe9+P895+Crr2DWRv/11QCp9dowuPkABqYMZECzASTVqXwvcWOMqe5sIxVjjAlzO7r774/d\nYfTtHD6wkyPb15J74CvqZGyhReZGUjI/gj3Ax3CMhhyOa092w67EtOhF0479ad62C717R8JG/9c9\n/PKvOdd9Neld3uK1Ha8hCJ0bdmZgykAGpgykT5M+xEfHV/0PwBhjqinLQBtTCZaBNpVxOeUhZzKO\ncuCbNZzbt4HIY5tolLmD1ML93y1czHJ12B/TjjkpWZyJvPj3e93CeOo8v5pPPwWvFBDXdgsJ3VbR\nfPBqvE2/xisFREVE0bNRTwalDGJgykB6NO5BdER0lf4sjDEmHFgJhzFBYAG0CQe5OVmk79jAyV3r\ncIc3knjmGzrlb6GU3dPJi2nA2dnLcVFtWLxYtxX/8EMoKACJyaFuhw20H72KxO6rORWzDYcjLiqO\nvk37fhdQ3/PBPdYizxhTI1kAbUwQWABtwpV7MKnM/tWHacyhhB4UNO9HXOthbNozkDfeiGXpUsjO\n1nMi656h89VraTdiFblNV3M4b98lr2sdPowx1ZkF0MaUICLjgMeBSGCuc+7hUs65GXgQbYTwtXNu\nZlljWgBtwpafBYoFcY1Zn3obMYfW0uLcZpqimeQcF8Pe2E6catCbdO9gPvryKuYvasHp0+ef27bn\nEQZOXcPmVv/u97IbZ21EpLTeIsYYE/4sgDamGBGJBHYAVwPpwFpghnNua7FzOgCvAaOcc6dEpIlz\n7lhZ41oAbcJVeftXHzmwi4ObVlCwbzUNTn5F24JdxEghAOnSjD1RPdl4uj9vrBzOx2v64fVG0f1Z\n/x0+msY3ZWTqSEa2Gkn/pv2JjrTaaWNM9WEBtDHFiMhg4EHn3DW++78FcM79T7Fz/gzscM7NLe+4\nFkCbcHW5/atzs8+xb9PnnN6xktgj60nN2kwjNA2d5eqwxduZO9uf9Pv8ES1Gs+rIZ+QW5pIQncDQ\nlkMZ1WoUQ5sPpV5MvYC/TmOMCSQLoI0pRkSmAuOcc3f57t8KDHTO/bjYOW+gWeohaJnHg865d0sZ\n627gboANGzb07dWrVxBegTGh4bxeDn+7g0ObV1D47SqST33Nja0K/J6/50ebGDs+lytvWEVuy2V8\nfnQ5p/JOER0RzYCUAYxKHcXI1JE0jm8cxFdhjDHlYwG0McWUM4B+CygAbgZaAp8APZxzp0sZErAM\ntKmdRjzdtfQtyD2F/PirLry14xreWjORs1lNGTWmkMGTv0Y6LmPViWUcyNS67J6NejKy1UhGpY6i\nbf22VjdtjAkLtpGKMRc6CKQWu9/Sd6y4dGC1c64A2CsiO4AOaL20McZn+dkIOLP/ouN5Uoes5DVM\nHfwhnkG/YW1uV5Z+O4pFf5rE5j2/YMjQXzJl6m7ir1zGhrMf8/iXj/P4l4/TJrENI1NHMqrVKH72\n8c+sRZ4xJuxZBtrUCiIShZZnjEYD57XATOfclmLnjEMXFt4mIo2ADUAv59zFf819LANtaqOyFii2\nHzWbnRtWcPKrJTQ6tIIO3t0AfOtpynvHhrJk6wQ+WDuO7j3qMm7qURoPXs7W/GWsObwGj/OUeV1r\nkWeMqWqWgTamGOecR0R+DLyH1jf/yzm3RUQeAtY55xb7HhsrIluBQuDXZQXPxtRWy2JGMGrIwxcs\nUDzgW6DYMTKSTv1GQb9RABxN382+L94gZvcH3Bq5hLubLyRrdCwrzvbj7c/H8uyjN5KYPI3rpmTS\nesRKnkz/txC/OmOMuTTLQBtTCZaBNqb8crPPsX31UrI3vUObjE9JkeMAfJ3XnqUHR/Lmhhs496t7\n/T7/dwN/x4R2E0iMSQzWlI0xtYwtIjQmCCyANubyOK+XfdvWcnjNmyTu/5gu3m1EiqNH21ZlPi82\nMpaxrccypeMU+jTpY4sPjTEBZSUcxhhjwpZERNC220DadhsIwKnjh9n1+Rtw+v/8PqfDqnnUH7aQ\nj/a/zZI9S2hbvy1TOkzh+iuup2GdhsGaujHGWAbamMqwDLQxgVVWi7wVd+nGofH1sxl82/tE9VnI\nQb4iKiKK0a1GM7nDZAalDCJCIoI9bWNMDWElHMYEgQXQxgTYo93hzIGLDudE1ONvMQd4480IVq06\nfzy+1S56zFyIp9MS8uQMLeq1YHKHyUxqP4km8U2COHFjTE1gAbQxQWABtDGBVVqLPI+LIEq8rG40\nhX73PMWRo1G88QYsWgQrVkBhIUh0Hol9P6Lt9YvwtFhNBJEMTx3GlA5TGNpiKFERVrFojLk0C6CN\nCQILoI0JrCdX7GZU/vILWuRt7/oztm/8nBuyFvJlveF0vfdV6sTVBeDECViyRIPpDz6AvDyIabKf\nBsMX0XjkG1A3g4bRTZjaeRKTO07mlrdvsY1ajDF+WQBtTBBYAG1M8Kx66Q8M2vlXtsT0IPVHb5KY\nlHzB45mZsHSpBtNvvw3ncgpIuPITGl61kHo9VyIA4v9vnm3UYoyxANqYICgtgC4oKCA9PZ3c3NwQ\nzSp81alTh5YtWxIdHR3qqZhqat1bT9Fz7f2kR6aScNebNG7eptTzcnPhww81mH7zTcjkCEnDXqfp\njf67fFgAbYyxANqYICgtgN67dy8JCQkkJydbj9pinHNkZGSQmZlJ27ZtQz0dU41t+uRN2n10N5mS\nQN6MBbTu1KvM8z0e+PRTDaaX9+vh97x10zcRGxvo2RpjqhPrA21MiOTm5tKmTRsLnksQEZKTkzl+\n/Hiop2KquR7DJ7KrfjJJr99CnVcm8M11z9O532i/50dFwciReuvxnP9xGzeG666DyZNh3DioV68K\nJm8Mpdf67+j+c5bFjOCeq66o1NgDXhhGjvf0RcfjIpJYc+unYTt2VY8f6LGtWaYxVcCC59LZz8UE\nSvsrh5I3613OSQKtl0znq49erfSYmZnwyitw000aTE+aBM8/DydPBmDCxhTz0v6ZTDn0GD0aCj3a\nptKjoTDl0GO8tH9mpccuLUgs63i4jF3V4wd6bMtAG2OMqZZatOtCxj0fkp42ie6f/JC1Z47Sf/JP\ny+HTfs4AACAASURBVHxOcp1kv104du2C11/XUo8vvtDa6TffhMhI+OlP4ZFHquqVmHD0j3/Aiy8G\ndkyHl6S6Bzk9K7vUx097s3nkl3PwipdCcXjFi1cchRGOQvHixFEY4TtW9NgF53mhjHXts/80vnIv\noFkVjl3V45dz7Adu+ke5hrMaaGMqobQa6G3bttGlS5cQzUhFRkbSo8f5Ws/p06dz//33B2Tsr776\nikOHDjF+/OX9MguHn4+pWbIyT7P7iSn0zF3HF21+yKBZ/41EXPoD1kU7F/H7z3/P89c+T+8mvS94\n7NAhDZ7nz4ePP4akJDh1qqpegQlHqWO1zOIPvR6lVf109p9pye+/+jlv5Izg7Br/ZRZJiSfo0HI7\nVzTbRZtG+2hb/1va1EunddwhUqOOES959GjbKoiv5LwYb+VivvwI/58iVnbsqh6/vGPPv2qxLSI0\npqqFawBdr149zp07VyVjP/vss6xbt45//KN879JLCoefj6l5CvLz+OqJ79P/zPusbjSZfvekERlV\n9oes2QXZjJ4/mpGpI/nvYf9d6jmnT0ODBpCYCGfOVMXMTbj60YSn+Uu/+4mX85v65LgYvmj1EJ5G\nAyk4sws5t4eY7H0k5KXTwHOIpt4j1JesC8bJdHEciWjGqajmZMakkhffmt/Ue8XvdW9r9hBRRBEp\n0URKNFESTZSUuF/i8UiJ8n2N5ic7B/odO61T5TrNzNnufxFuZceu6vHLO3aTJraI0JiQq6qS38t5\n33vmzBkGDBjA4sWL6dSpEzNmzGDUqFHMmTOHH/7wh6xdu5acnBymTp3KH/7wBwDWrl3LT3/6U7Ky\nsoiNjeWDDz7ggQceICcnh5UrV/Lb3/6WadOmBfjVGVNx0TGx9L3vVb6Yex+DD7/Il48ep+uP5323\n4Upp4qPjmdBuAq/vfJ3fDPgN9WPrX3ROUSLb662qmZtw9dtOj14QPAPEST6jDtwPxXabz3PRHI1s\nwunY5nxTtwcuqTUxjdqSkNKeJqkdSWzQmA7FPhE5kXMCXvMfQP/qmhsrN/Gd/h8aNKhyQ7O9Cseu\n6vHLOfaePeUbzgJoY2qgnJwcevU639qrKND9xz/+wezZs/npT3/KqVOnmDNnDgB//OMfadiwIYWF\nhYwePZqNGzfSuXNnpk2bxrx58+jfvz9nz54lPj6ehx56qFIZaGOqSkRkJIN/8ASrXm7GoB1/Yeuj\n19Dih29Qv0Ejv8+5qeNNzNs+jyW7l/D9rt+/6PHISP1aWFhVszbhqkVCeqnHHbC+z5+o17Qdyakd\nSW6aSqvISMpTlHE8+zh3vn9nQOdZUlxEkt9uE+E8dlWPH+ixLYA2pgqFqkIqLi6Or7766qLjV199\nNfPnz+fee+/l66+//u74a6+9xlNPPYXH4+Hw4cNs3boVESElJYX+/fsDkJiYGLT5G1MZg2b+B+ve\nbkbPNb/h4N9HkXfnmzRpUXrv8U4NO9GzUU/m75jPLV1uuahTjGWga6+MvAY0rnNxCxapn0q/G+6p\n8HhHs45y1/t3cTT7KIkxiZzNP3vROcl1kkt5ZsUEop1cKMau6vEDPbYF0MbUIl6vl23bthEfH8+p\nU6do2bIle/fu5S9/+Qtr166lQYMGzJ4923ZRNNVevwlz2JzYhDYf/oBzaWP4dvp8WnfuU+q5UztO\n5YHPH+DLY1/St2nfCx4rCqAtA1277N26lpTYs3idEFFs+/ccYjjQ/ed0rOB4R7KOcOd7d3Ii5wT/\nvPqfFy1aNdWP9YE2phZ59NFH6dKlCy+//DK33347BQUFnD17lrp161K/fn2OHj3K0qVLAejUqROH\nDx9m7dq1AGRmZuLxeEhISCAzMzOUL8OYcuk+bCJHJi8kCg9Jr17HN2s+KPW8a9pcQ73oeizYseCi\nx4pKOCwDXXucOn6Y2Pm3cMYl8Mv1v8ObmAoI1E/lwJCHWRYzokLjHTp3iNvfvZ2TuScteK5BLANt\nqh0R+UVZjzvn/hqsuYSrkjXQ48aN4/bbb2fu3LmsWbOGhIQEhg8fzn/913/xhz/8gd69e9O5c2dS\nU1MZMmQIADExMcybN4+f/OQn5OTkEBcXx4cffsjIkSN5+OGH6dWrly0iNGGv/ZVDOJjwLjkvTKb1\n2zP56uzf6DVmxgXnxEfHc12761i0cxG/6f8bkuqcr4m0Eo7aJT8vl0NpN9Hee5Jr3lnAivUj+ONr\n/0Z8vD7e0Xcrr/TMdO58704y8zN56uqn6NHYfycIU71YGztT7YjIy0B/YLHv0PXAGnxrj51zfwjW\nXMK1jV04s5+PCYWTxw5y4qlJXFGwk42Nr6d3wYYLtlCeTxtePfQzft3v18zqNuv/b+++w6Os0saP\nf+8EUghJCISiJHQiLfQOASyUfUWwgIqrCMgKigsi7G/19XUBV9ZFWRHEnyxKc1eUIi5go5e8KJAg\nAQkIBKQkIBJKYgppc94/ZhISSDKTTJJJuT/XlWvmOc+ZM/dhJuGeM+ecJ89js6dFZ2XdTKhV5WMs\nFiIWjqb71Y1EdnmL3o9MICMDbtwAT8+it3c+8TzPbH6G5IxkFg9aTNs6bUs+aFXiTp92bBs7/VOg\nKqIgoLMxZpoxZhrQBWhkjJlVlsmzUqriqF2vIQ2nbOGse2M6xW+AhPOAgYTzBO95mZGcoX1d62LC\nWweWdBS6ati36k26X93I9w3H0vWBCTmvd3E+NJ1NPMvYTWNJyUzho0EfafJcCWkCrSqi+kDuzTnT\nbWVKKVUgH99aNK2ZeVu5N+mEHJnHyJCRnEk8w4FLB/Kc14WEld+Pu9bR7ae3OVijDz3GWa/Znp1A\nZ8+Dd9TPCT8z7ttxpGels2TQElrX0W/cKiNNoFVF9DGwX0RmishMYB+w3KURKaUqBEmMy/9EQiyD\nmwzGt7ova06syXNKFxJWbudORNF4xyTOujcm5LmVuLm7Y8zNbUiLckGs09dPM27TODJNJksGL+Gu\n2neVTtDK5TSBVhWOMWY2MBa4ZvsZa4x507VRKaUqBP+gAsu9q3kztPlQtpzdwrUb13JO6RSOyivh\n6mXk01FkUg3vp1fj42tdQJr9Wos4nkDHXIth7KaxGGNYOngpLQNallLUqjzQBFpVSMaYH4wx820/\nB10dj1KqYjjRbiqpeOQpS8WDE+2mAtYrE2ZYMthwakPOeb0aYeWUmZHOuX+OpL7lEr/+bgl3NL45\nWpz9Wjs6feP41eOM2zQOd3Fn6ZClNK/VvBQiVuWJJtBKKaWqjO0eAzjf5+/gH4zBOrR4wrdXzt6+\nLQNa0rFuR9aeWJuzmFBHoCunA4ufIzTtIFEdZ9Gqx6A854qygPCnqz8xfvN4qrtXZ9mQZTTzb1YK\n0aryRhNopVxo0a5TfHcqPk/Zd6fiWbTrlNNtx8bGMnz4cFq2bEnz5s2ZMmUK6enp7Ny5k6FDh+b7\nmCZNmhAfH5/vOaUqg4n9mxMy8BmYegSZeZ2DPn1pnrifx9rUyKkzImQEZxLPEHkpEtBFhJXRvjVz\n6XF5LXvrj6L7Q3+87byjCwijr0TzzKZn8KrmxfLBy2ns17gUolXlkSbQqsoQkSEiclxEYkTk5ULq\nPSIiRkTs7gPprPZB/ryw8mBOEv3dqXheWHmQ9kH+TrVrjOHhhx/mwQcf5OTJk5w4cYKkpCReffXV\nkghbqUqj9tDX8eYGx9fe3AFzcJPB+HrcXEyoiwgrl+g9X9H5yN845N2dbn9YmG+d7A9LhY1A/3j5\nR/6w6Q/UrF6TZYOXEewXXArRqvJKr0SoqgQRcQfeBwYCsUCEiGwwxhy9pZ4vMAXrzh5Om7UxmqMX\nEgutU8/Xk9FL9lPfz5NLiWm0qFeT+VtPMn/ryXzrt7nTjxkPFL6n6Pbt2/Hy8mLs2LEAuLu7M2/e\nPJo2bcrdd9+dU+/KlSuMGjWKuLg4evXqddv+t0pVdo1bdyEiYAidflnLL+en0yC4BV7VvBjWfBir\nj6/m2o1ruLkFAJpAVwZxp6NpuOVZ4tzvpNnEVbhXyz8NsjeF49DlQ0zcMhF/T3+WDl7KnTXvLKWI\nVXmlI9CqqugOxBhjThtj0oHPgOH51PsrMAe4UVaB+XtXp76fJ3HXb1DfzxN/7+pOtxkdHU2XLl3y\nlPn5+dGoUSNiYmJyymbNmkXfvn2Jjo7moYce4ty5c04/t1IVTdBDsxAM59bNyCkb0XIEGZYM1ses\n10WElcRvCVfJ/PdjAFR/chW+/rULrFvYIsKDvx5kwpYJBHgFsHzIck2eqygdgVZVRUPgfK7jWKBH\n7goi0hkINsZ8JSJ/KqghEXkWeBbg4MGD3Hop79zsjRTDzWkbk+9pwb/3nWPKfS3p3TzQ7uNKwu7d\nu1m3bh0A999/PwEBAWXyvEqVJ3c0vou99R+m26U1nDsRRaOQjrQIaEGnep1Ye3Itbm5PA6Ij0BVY\nVmYmpxc9RpusCxwftIJ2zQr/21zQCHTkL5E8v+156tWox5JBS6jvo9fwqqp0BFopQETcgHeAafbq\nGmMWG2O6GmO6BgY6l+hmJ88Ln+jES4PuYuETnfLMiS6uNm3acOBA3qupJSYmcu7cOVq0aOFU20pV\nRi0fmUEaHsRv+EtO2ciQkZxNPItH8whAp3BUZBEfTaZD6n5+aPfftOvzgN36+S0ijPglgue3PU8D\nnwYsG7xMk+cqThNoVVXEAblXeATZyrL5Au2AnSJyBugJbCjthYSHYxNY+ESnnBHn3s0DWfhEJw7H\nJjjV7r333ktKSgoff/wxAFlZWUybNo0xY8ZQo8bN3Qb69evHypUrAfjmm2+4du1avu0pVdnVqR/E\n4eAn6Zy0i5NR4QAMbDwQPw8/PLtbFxPqFI6Kaf8X79Hzl0/YF/gIPUZOd+gxWVlw1/wB1H87lNAV\n1p9xm8aRmpnK9RvXqVujbilHrco7TaBVVREBtBSRpiLiATwO5FwpwRiTYIwJNMY0McY0AfYCw4wx\nkaUZ1MT+zW+brtG7eSAT+zu3Cb+I8MUXX7BmzRpatmxJSEgIXl5e/O1vf8tTb8aMGezevZu2bduy\nbt06GjVq5NTzKlWRtRnxKtepSeq31rnQ2YsJq7XZirvvVR2BroB+2reZjlEz+dGzE10mLHL4cRYL\nVPe/ku+5a2k60KB0DrSqIowxmSLyArAJcAeWGmOiReR1INIYs6HwFiqe4OBgNm7ceFv5gAEDGDBg\nAAB16tRh8+bNZRyZUuWTX6067G3xB3rGzCN6z1e07XM/I0JG8O9j/yag73qyssa6OkRVBBfPHqfu\nN+O55FaPRhPWUK26h/0H2ei3DcoeHYFWVYYx5mtjTIgxprkxZrat7C/5Jc/GmAGlPfqslCp/Oj48\nnV+pTbUdr2MsFprXao7lfGcC+q8lM0uHoCuK5N+uk/LxY1QnAzPqU/xrF23KhX7boOzRBFoppZSy\n8apRkzPt/shdmT9xaNtnAJioEXg2OMeRhAgXR6ccYcnK4sQHT9Ak8wxn736fRiEdi9xGWmZ6KUSm\nKhNNoJVSSqlcOg9/gfNyJ7W+f5OszEzcTg4iM8mPLZfWuDo05YB9S6fRKWUPEa3+RGj/h4v8eGMM\nH5ycXQqRqcpEE2illFIql2rVPfi123SaWM7xw1eLcTOeXN8zjIiEbVxJzX9hmSofIjf+k15xy9hf\n+wF6PPZKsdpY+dNKtv26jqwb3vmer+NVx5kQVSWhiwiVUkqpW3QaPIaYA+/TMOpdPKqP49rOkQQO\n/jfrT61nXLtxrg5P2SzadYp70ncScmQeJiGWLsZwzq0hB9q+SveCrsNdiL0X9/J2xNt0rTWA5Q/N\np9Vdbhw7VgqBqwpPR6CVUkqpW7i5u5PS97+501ziwY7vknaxGXfV6MLaE2uxGF1hVl7ck76T4D0v\nQ8J5BIMI1DOXudfyXZHbOv/beabvmk4TvyY81+RNMG63XYlQqWz61lDK1Q6vhnntYGYt6+3h1U43\nKSJMm3bzoopz585l5syZTrerVFUS2v9hoj1CmdDwA3xqJNCv1gjO/3aefRf3uTo0ZRNyZB7e5F3w\n50U6IUfmFamd5IxkJm+fjDGG9+55Dy+pCeS9EqFSuWkCrZQrHV4NGydDwnnAWG83TnY6ifb09GTd\nunXExzt3SXClqjJxc8N94Ezqul3n/w2fQ8caA6nlWYu1J9a6OjRlYxJi8z9RUHk+LMbCK+Gv8HPC\nz8ztP5dgv+CcfaB1BFoVROdAK1WavnkZfvmx4POxEZCVlrcsIxXWvwAHVuT/mAah8Lu/F/q01apV\n49lnn2XevHnMnp13Nfnly5eZOHEi586dA+Ddd9+lT58+hIaGEh4ejr+/P4GBgcybN4/Ro0czevRo\nnnrqKQYOHGi3u0pVNq263cf2tb2Y0nIZB5KnMKz5MFYeW0l8ajyB3oH2G1Cl5sqlWPxwpzqZt5/0\nD3K4nQ8OfcCO8zv4c7c/0+vOXsDNfaA1gVYF0beGUq50a/Jsr7wIJk2axCeffEJCQkKe8ilTpjB1\n6lQiIiL4/PPPGT9+PAB9+vRhz549REdH06xZM8LDwwH4/vvv6d27t9PxKFVRfXRiFr6kknn4dR4J\neYRMk8l/Yv7j6rCqtEuxp0j652CMMaTdMhaYigcn2k11qJ0tZ7ew6NAiHmzxIL9v/fuc8uwRaJ3C\noQqiI9BKlSY7I8XMa2ebvnEL/2AY+5VTT+3n58fo0aNZsGAB3t43t2PaunUrR48ezTlOTEwkKSmJ\nsLAwdu/eTePGjXnuuedYvHgxcXFxBAQE4OPj41QsSlVkcVd7sPbqfQxjFb8l/4mu9bvy+YnPGddu\nHG6i41BlLe70MeRfw6hj+Y3lLRYwoEGadc5zQiz4B3G+3VS2ewwgxE47x68e59X/fZUOdTvwWs/X\nEJGcczoCrezRt4ZSrnTvX6D6LXuNVve2lpeAF198kSVLlpCcnJxTZrFY2Lt3L1FRUURFRREXF0fN\nmjXp168f4eHhhIeHM2DAAOrWrcvatWsJCwsrkViUqqjc3OB/1v8Vdyyc/nwGI0NGEpsUy96Le10d\nWpVz9qcfqP7xf1HDpHDxwdU8+9SThAx8BqYegZnXYeoRQgY+w8T+zQtt5+qNq0zePhlfD1/mDZiH\nh7tHnvPZCbSOQKuCaAKtlCu1fxQeWGAdcUastw8ssJaXgNq1a/Poo4+yZMmSnLJBgwbx3nvv5RxH\nRUUBEBwcTHx8PCdPnqRZs2b07duXuXPn0q9fvxKJRamKyt0dTp5vzW6v4XSO30irzDt0MaELxBza\ng99nw3HDwrVH/0PLTsX725RhyWDazmnEp8Yz/+751K1R97Y6uohQ2aNvDaVcrf2jeUZPSip5zjZt\n2rQ8u3EsWLCAyMhI2rdvT5s2bVi0aFHOuR49ehASYv3iMywsjLi4OPr27Vui8ShV0WQnUckhM8nE\nnasb/8rw5sPZcW4H8am6001Z+CliK/W+GEEanqQ++SVN23Qrdltz9s8h8lIks/rMol1gu3zr6BQO\nZY/OgVaqEkpKSsq5X79+fVJSUnKOAwMDWbVqVb6P+9e//pVzv3fv3lgsesEIpbK/xq/m3ZiohqPo\ndWEFKYxmhW0x4fjQ8a4NsJI7Er6eZlv/wBW3OlQfu4HgRi2L3dbq46tZdXwVY9uOZWizoQXW00WE\nyh79bKWUUkoVInsU0mKBNiNeIxEf/Lf/f7o16KZXJixlUVs/peXWZ7jk3gDvCZto4ETyHPlLJG/u\ne5O+DfsypfOUQuvqCLSyR98aSimlVCFyJ9D+tesS3XQsHVL30dsthLikOL6/8L1rA6ykDnz1EW3D\nJ3G2elNqT9pCYINGxW7rQtIFpu2aRpBvEHP6zcHdrfChZV1EqOzRBFoppZQqRHYSlf21fscRL3OZ\nALrvX0+AZ4AuJiwFEevm02n/dE56tuGOP27Cv079YreVkpHClB1TyMjKYME9C/Dz8LP7GF1EqOzR\nt4ZSSilViNwj0ADePr6cbvM87TOO0sujFTvO7+ByymXXBVjJ7P10Nt0O/4Uj3l1oOuUbfP1rF7st\nYwyv7XmN41ePM6ffHJr6N3XocTqFQ9mjbw2llFKqELeOQAN0Gj6ZOKnP0OP7yDJZfBHzhWuCq0SM\nxcL3y1+m5/G3OOjTl7te3Ii3j69TbX7444dsPruZqV2mEhbk+J72uohQ2aMJtFJKKVWIW0egATw8\nvbjY6SXC0s7QuloQn5/4XBcTOsFYLOz9cAq9znxAhP8gQl/8Ak+vGk61uePcDt47+B73N7ufMW3H\nFOmxOgKt7NFt7JRyoQGrBnDlxpXbyut41WHnYzvLPiCl1G3yS6ABOt//B05HfcDwC2f4e71qfHfh\nO/o21H3Ti8qSlUXEB+PpFb+OfXUepNvzS3Fzcug35loML4e/TNs6bZnZa2aey3Q7FJMuIlR26Gcr\npVwov+S5sHKlVNnLbwoHgJu7O4l9XuHR5Av44sWa42vKPrgKLjMjnQMLnqBH/Dr2Nvg93Sctczp5\nTkhLYPKOyXhX8+bdu9/Fq5pXkdvQRYTKHh2BVqoUzdk/h5+u/lSsx479dmy+5a1qt+LP3f9s9/Gz\nZ89mxYoV1KtXj+DgYLp06cL06dNzzr/99tt4enoyefJkpk6dyqFDh9i+fTvbt29nyZIlfPLJJzz3\n3HNERESQmprKiBEjmDVrFt9++y1LlixhzRprsrBz507mzp3Ll19+yebNm5kxYwZpaWk0b96cZcuW\nUbNmzWL1X6nyoqARaIAOdz/K8b3z+a/rV1kru/g15Vfq1ahXtgFWUOlpNziyYCTdknfzfeOJ9Hz6\nTcTJjDXTksn0XdP5JfkXlg5eSgOfBsVqR6dwKHv0raFUJXTgwAE+++wzoqKi+Prrr4mIiLitTlhY\nGOHh4QBERkaSlJRERkYG4eHh9OvXD7Am4ZGRkRw+fJhdu3Zx+PBh7rvvPvbt20dycjIAq1at4vHH\nHyc+Pp433niDrVu38sMPP9C1a1feeeedsuu0UqWkoBFoAHFzw3LvDEb/Fm9dTHhSFxM64kZKEsfm\nDaVz8m72tpxGr7FznE6eAf4R+Q/2XtzLaz1fo2O9jsVuRxcRKnt0BFqpUmRvpDh0RWiB55YNWVbs\n5w0PD+ehhx6iRg3rIpxhw4bdVqdLly4cOHCAxMREPD096dy5M5GRkYSHh7NgwQIAVq9ezeLFi8nM\nzOTixYscPXqU9u3bM2TIEDZu3MiIESP46quveOutt9i1axdHjx6lT58+AKSnp9OrV69i90Gp8qKw\nEWiANj2HcGhnJ7qknmfN8TWMDx1v90IdVVlS4jXOLnyA0LQj7A+dQc8RLxW7rYLWkcz/YT4PtXyo\n2O3qCLSyRxNopaqI1NRUOna0jshMnDiRiRMn0rRpU5YvX07v3r1p3749O3bsICYmhtatW/Pzzz8z\nd+5cIiIiCAgIYMyYMdy4cQOAxx9/nIULF1K7dm26du2Kr68vxhgGDhzIp59+6spuKlXi7CXQAD6/\nm8WoTY8w3duD7y58V6Qt0yqzRbtOcU/6TkKOzIOEWCy+d5DwWxZ3mcv80O0tug991qn2S2sdiS4i\nVPZoAq2UC9XxqlPgLhzO6NevH2PGjOGVV14hMzOTjRs3MmHCBKKiovLUCwsLY+7cuSxdupTQ0FBe\neuklunTpgoiQmJiIj48P/v7+XLp0iW+++YYBAwYA0L9/f8aNG8eHH37I448/DkDPnj2ZNGkSMTEx\ntGjRguTkZOLi4ggJCXGqL0q5WmFTOLK16NCH0QfqAxae3/Z8nnPebrXY/1S4UzGU5o49pdn2J+ee\n4H1LCtQWqB1sLQyshp9pyh4HkucMSwYJaQkkpCVwPe0619Ou57lfWnQRobJHE2ilXKi0tqrr3Lkz\njz32GB06dKBevXp069Yt33phYWHMnj2bXr164ePjg5eXF2Fh1pGzDh060KlTJ1q1akVwcHDO1AwA\nd3d3hg4dyvLly1mxYgUAdevWZfny5YwaNYq0tDQA3njjDU2gVYWXnUQdOwbbthVc7zf3/IeoUy3X\n2bK1kOzbAYWNtJZm25u3ZuQcm9tq5C3JsmSQkXqFzLR4TNo1LBlXuO6ekm/biZLBWxu/JsVynWRL\nAsmW66Rk5bpvK79hkgqM253qhfarsNfKnuho660m0KogYsztvxJKVUYiMgSYD7gDHxlj/n7L+ZeA\n8UAmcBkYZ4w5W1ibsbGxJigoKE/ZsWPHaN26dUmG7rSZM2dSs2bNPLtwuEp5/PdRqjAvvgjz59uv\n1255wWsaVOGyUnzJSvInM6kWWUn+ZCXVunk/uZbt2Fqe/WNJ86bd8vYFtnlkzI9OxzVhAixa5HQz\nqgI5ffr0gWbNmnW1V09HoFWVICLuwPvAQCAWiBCRDcaYo7mqHQS6GmNSROQ54C3gsbKPVilVnowd\nCzExkJpaeL1fCzk36Zpz0w3eD6jlkrZfuHYdyTXOlok7N4wHN4wnacaTdONJmsWTNONFmsWL9Cxv\n0i3e1tusGpwN+brAtmuvX49bmj+S5o8YO+mIj+2nft7iwv7N77mn8Cbt8fKC8eOda0NVXppAq6qi\nOxBjjDkNICKfAcOBnATaGLMjV/29wJNlGmEpmjlzpqtDUKrC6tABvvzSfr3QFQWf69N5IZgsjDEY\nY4HctxbbFAyThbFYMMZ2H2NdzWYMXHu/wLY7NZxirWMsgMFYrLe5y3LuZz8vBjFgTBZQcJLbv+9y\nvH3r4ONfB7+Aunh4Fu2iJKErCm5717pmRWorPwNWFbyOxJkpHErZowm0qioaAudzHccCPQqp/wzw\nTX4nRORZ4FmAgwcPcusUDgBjTJEvHVsV6JQxVVWF9hvuXAMrCk6ge4x0cmpWIUluq673Otd2KSut\ndSRK2aPT45W6hYg8CXQF3s7vvDFmsTGmqzGma2Bg4G3nvby8uHLliiaLtzDGcOXKFby8in5ZXaUq\nAm+3/KdCFFReFLXcahSpvLy0XdCOQs7uNKSUq+kItKoq4oDgXMdBtrI8ROQ+4FWgvzEmrThP9viw\npgAABztJREFUFBQURGxsLJcvXy5WoJWZl5dXviP2SlUGzm5VV5jfN1qZZz9l/IM40W4q2z0GlOu2\ndYRYVVa6C4eqEkSkGnACuBdr4hwBPGGMic5VpxOwFhhijDnpSLv57cKhlFJKqYrJ0V04dAqHqhKM\nMZnAC8Am4Biw2hgTLSKvi0j2da7fBmoCa0QkSkQ2uChcpZRSSpVjOoVDVRnGmK+5Zbm5MeYvue7f\nV+ZBKaWUUqrC0RFopZRSSimlikBHoJVyQnp6evzp06cLvVqhoy5fvhxYt27d+JJoqzzTflYu2s/K\nRftZuWg/i6WxI5V0EaFS5YSIRBpj7C5cqOi0n5WL9rNy0X5WLtrP0qNTOJRSSimllCoCTaCVUkop\npZQqAk2glSo/Frs6gDKi/axctJ+Vi/azctF+lhKdA62UUkoppVQR6Ai0UkoppZRSRaAJtFJKKaWU\nUkWgCbRSLiIitUVki4ictN0GFFLXT0RiRWRhWcZYEhzpp4g0FpEfbJdQjxaRia6I1RkO9rOjiHxv\n6+NhEXnMFbE6w9H3rYh8KyLXReTLso7RGSIyRESOi0iMiLycz3lPEVllO79PRJqUfZTOc6Cf/Wy/\nk5kiMsIVMZYEB/r5kogctf0+bhMRh/YALm8c6OdEEfnR9jf2f0WkjSvidJa9fuaq94iIGBEpta3t\nNIFWynVeBrYZY1oC22zHBfkrsLtMoip5jvTzItDLGNMR6AG8LCJ3lmGMJcGRfqYAo40xbYEhwLsi\nUqsMYywJjr5v3waeKrOoSoCIuAPvA78D2gCj8kk0ngGuGWNaAPOAOWUbpfMc7Oc5YAywsmyjKzkO\n9vMg0NUY0x5YC7xVtlE6z8F+rjTGhNr+xr4FvFPGYTrNwX4iIr7AFGBfacajCbRSrjMcWGG7vwJ4\nML9KItIFqA9sLqO4Sprdfhpj0o0xabZDTyrm3yZH+nnCGHPSdv8C8CtQt8wiLBkOvW+NMduA38oq\nqBLSHYgxxpw2xqQDn2Htb265+78WuFdEpAxjLAl2+2mMOWOMOQxYXBFgCXGknzuMMSm2w71AUBnH\nWBIc6WdirkMfoCLuIOHI7ydYB5zmADdKM5iK+J+UUpVFfWPMRdv9X7AmyXmIiBvwD2B6WQZWwuz2\nE0BEgkXkMHAemGNLMCsSh/qZTUS6Ax7AqdIOrIQVqZ8VTEOs779ssbayfOsYYzKBBKBOmURXchzp\nZ2VQ1H4+A3xTqhGVDof6KSKTROQU1hHoyWUUW0my208R6QwEG2O+Ku1gqpX2EyhVlYnIVqBBPqde\nzX1gjDEikt+IwPPA18aY2PI8yFUC/cQYcx5ob5u68R8RWWuMuVTy0RZfSfTT1s4dwL+Ap40x5W6E\nr6T6qVRFISJPAl2B/q6OpbQYY94H3heRJ4D/AZ52cUglyjbg9A7WqUelThNopUqRMea+gs6JyCUR\nucMYc9GWUP2aT7VeQJiIPA/UBDxEJMkYU9h86TJXAv3M3dYFETkChGH9irzcKIl+iogf8BXwqjFm\nbymF6pSSfD0rmDggONdxkK0svzqxIlIN8AeulE14JcaRflYGDvVTRO7D+uGwf66pZBVJUV/Pz4AP\nSjWi0mGvn75AO2CnbcCpAbBBRIYZYyJLOhidwqGU62zg5gjA08D6WysYY35vjGlkjGmCdRrHx+Ut\neXaA3X6KSJCIeNvuBwB9geNlFmHJcKSfHsAXWF/HcvXhoAjs9rMCiwBaikhT22v1ONb+5pa7/yOA\n7abiXZHMkX5WBnb7KSKdgH8Cw4wxFfXDoCP9bJnr8H7gZBnGV1IK7acxJsEYE2iMaWL7P3Mv1te1\nxJNn0ARaKVf6OzBQRE4C99mOEZGuIvKRSyMrWY70szWwT0QOAbuAucaYH10SbfE50s9HgX7AGNt2\nUlEi0tE14RabQ+9bEQkH1mBdZBcrIoNdEm0R2OY0vwBsAo4Bq40x0SLyuogMs1VbAtQRkRjgJQrf\nPadccqSfItJNRGKBkcA/RSTadREXj4Ov59tYv91bY/t9rHAfJBzs5wti3T4zCuv7tsJN33Cwn2VG\nL+WtlFJKKaVUEegItFJKKaWUUkWgCbRSSimllFJFoAm0UkoppZRSRaAJtFJKKaWUUkWgCbRSSiml\nlFJFoAm0UkoppZRSRaAJtFJKKaWUUkWgCbRSSilVANsFRQ6LiJeI+NguRtHO1XEppVxLL6SilFJK\nFUJE3gC8AG8g1hjzpotDUkq5mCbQSimlVCFExAOIAG4AvY0xWS4OSSnlYjqFQymllCpcHaAm4It1\nJFopVcXpCLRSSilVCBHZAHwGNAXuMMa84OKQlFIuVs3VASillFLllYiMBjKMMStFxB34TkTuMcZs\nd3VsSinX0RFopZRSSimlikDnQCullFJKKVUEmkArpZRSSilVBJpAK6WUUkopVQSaQCullFJKKVUE\nmkArpZRSSilVBJpAK6WUUkopVQSaQCullFJKKVUE/wcwyf9mUh8iqQAAAABJRU5ErkJggg==\n",
"text/plain": [
"<matplotlib.figure.Figure at 0x1157eef10>"
]
},
"metadata": {},
"output_type": "display_data"
}
],
"source": [
"q_l = 0.7\n",
"q_r = 0.5\n",
"v_l = 1.0\n",
"v_r = 1.5\n",
"t = 0.5\n",
"\n",
"states, speeds, reval = riemann_traffic_vc_exact(q_l,q_r,v_l,v_r)\n",
"riemann_tools.plot_riemann(states, speeds, reval, t,layout='vertical');\n",
"x1, q1 = test_solver(q_l,q_r,v_l,v_r,traffic_fwave_old)\n",
"plt.plot(x1,q1,'-x')\n",
"x2, q2 = test_solver(q_l,q_r,v_l,v_r,traffic_fwave)\n",
"plt.plot(x2,q2,'-o')\n",
"x3, q3 = test_solver(q_l,q_r,v_l,v_r,traffic_qwave)\n",
"plt.plot(x3,q3,'-s')\n",
"plt.legend(['Exact','Old','New','q-wave'])"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"# q-wave method has no entropy glitch"
]
},
{
"cell_type": "code",
"execution_count": 16,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/plain": [
"<matplotlib.legend.Legend at 0x115a651d0>"
]
},
"execution_count": 16,
"metadata": {},
"output_type": "execute_result"
},
{
"data": {
"image/png": "iVBORw0KGgoAAAANSUhEUgAAAtUAAAEWCAYAAACpJ2vsAAAABHNCSVQICAgIfAhkiAAAAAlwSFlz\nAAALEgAACxIB0t1+/AAAIABJREFUeJzs3Xd81/W1+PHXySJhQ5gSVthTtuwpeykIIipIh62t1a57\nb/vrvWptr7WtvV6t3lpa24qDKShDNjJEQEARBGRFRtgkbBKy3r8/zjcQQhIyvivfnOfjkQfJ5/sZ\n728Svjmf9/e8zxHnHMYYY4wxxpjiCwv0AIwxxhhjjCntLKg2xhhjjDGmhCyoNsYYY4wxpoQsqDbG\nGGOMMaaELKg2xhhjjDGmhCyoNsYYY4wxpoQsqDbGGFMmichSEZka6HHciYg0EhEnIhGBHosxJn8W\nVBtjjCkxEfmliCzNte1APtsm+Xd0eXPODXfOveXr64jIYyLyia+vY4wJLAuqjTHGeMN6oKeIhAOI\nSF0gEuiYa1tTz77GGBNSLKg2xhjjDVvRILqD5+s+wMfAvlzbDjnnTgCIyCsickxELonIdhHp49l+\nl4ikiEj17JOLSEcROScikZ6vvyUie0XkvIgsF5GGnu0iIi+LyBnPeXeJSNu8Biwia0XkO57PHxOR\nT0TkJc85vxGR4fkc10REkkWkU47xnhWR/nns2wp4A+ghIldE5EIBY/mdiHzmGfeHOZ9/rn2neZ77\nZRFJEJHv5Xisv4gkisjPPN+DkyIyLcfj5TzP8aiInBaRN0QkJq/rGGOKxoJqY4wxJeacSwO2AH09\nm/oCG4BPcm3LOUu9FQ24qwPvAXNFJNoTdG8CxufYdzIwzzmXLiJjgf8HjANqeq4z07PfEM91mgNV\ngIlAUiGfxj3oTUAN4A/AmyIieTzXQ8B/AO+ISHngn8Bbzrm1eey7F/g+sMk5V9E5V7WA608BvgXU\nBTKAV/PZ7wwwCqgMTANezg7wPeqgz70e8G3gdRGp5nnsRfR70wF916Ae8EwBYzLGFJIF1cYYY7xl\nHTcD6D5osLsh17Z12Ts7595xziU55zKcc38CygEtPA+/BzwEOvsMTPJsAw1Sf+ec2+ucywBeADp4\nZqvTgUpAS0A8+5ws5PiPOOf+5pzLBN5Cg9vaee3onPsbcBC9kagL/KqQ1yjI2865r5xzV4H/AiZm\np87kuvYS59whp9YBK9DvbbZ04HnnXLpz7iPgCtDC8318HPiJcy7ZOXcZ/d4FRY67MaWdBdXGGGO8\nZT3Q25O2UNM5dwD4FM21rg60JcdMtYj83JPGcNGTFlEFnSUGeB9NmaiLBuVZaIAO0BB4RUQueI5L\nBgSo55xbA7wGvA6cEZHpIlK5kOM/lf2Jc+6a59OKBez/N89z+rNz7rrnOfXxpHlcEZHdhbxutmM5\nPj+CptPUyL2TiAwXkc2eFJQLwIhc+yV5bjayXfM8j5pAeWB7ju/dMs92Y0wJWVBtjDHGWzahgfF3\ngY0AzrlLwAnPthPOuW9Ag0/g39H0jGqetIiLaHCMc+48OgP7IJr6Mcs55zzXOQZ8zzlXNcdHjHPu\nU8+xrzrnOgOt0VSHf/P2ExWRisD/Am8Cz2XnPzvnNnjSPCo659p4dnf5nSeX+jk+b4DOOJ/Ldd1y\n6A3HS0Btz/ftIzzftzs4B6QAbXJ836o45wq6cTDGFJIF1cYYY7zCOZcCbAN+ys1ZZdC86p9yaz51\nJTRv+CwQISLPoDnCOb2H5hk/wM3UD9CFf78UkTYAIlJFRCZ4Pu8qIvd4FjReBVLRWW5vewXY5pz7\nDrDEM6b8nAbiRCTqDud8RERae/K0n0dzyDNz7ROFpsmcBTI8iymHFGbAzrksdHb9ZRGpBSAi9URk\naGGON8YUzIJqY4wx3rQOqIUG0tk2eLblDKqXo6kH+9FUh1RuTX8AWAg0A045577M3uicWwD8Hpgl\nIpeAr4DsSh2V0cDxvOe8ScAfvfHEsnkWSg4DnvBs+inQSUQezueQNcBu4JSInMtnH4C3gX+haSjR\nwFO5d/DkQT8FzEGf42T0+1RY/4Hmgm/2fO9WcTOP3RhTAnLz3TRjjDHGBIKIrAXecc79PdBjMcYU\nj81UG2OMMcYYU0IWVBtjjDHGGFNClv5hjDHGGGNMCdlMtTHGGGOMMSUUEegBGFOaHT582EVE2H8j\nY4wxJhSkpaWdi4+PL1ZDJIsGjCmBiIgI4uLiAj0MY4wxxnhBQkLCkeIea+kfJqSIyD9E5IyIfJXP\n4yIir4rIQRHZKSKdcjw2VUQOeD6m+m/UxhhjjCntLKg2oeZfaFOG/AxHm0k0Ax4H/gLgaTH8LHAP\n0A14VkSq+XSkxhhjjAkZFlSbkOKcWw8kF7DLWGCGU5uBqiJSFxgKrHTOJTvnzgMrKTg4N8aUEVeu\n6IcxxhTEgmpT1tTj1lbIiZ5t+W2/jYg8LiLbRGTbuXMFdRw2RXH4MCxYAOnpgR5Jwb74ArZvD/Qo\n7mzZMrh6NdCjuLMVKwI9gjvbuFE/yqJVqyA1NdCj8L/z52HevECPwv+ysmD5cti8OdAj8b89e+D9\n90t2DguqjSki59x051wX51yXGjVqBHo4IaNRIxg2DD78MLgDmI4dITxcxxnMZf4HD4Y5c/SPZDBz\nDi5dCvQoTF6SkiAjA6KjAz0S/9q9Gz7+GMaPD/RI/Mc5WLtWbyS6dYPu3QM9Iv+4ehUWLboZTJf0\nZ25BtSlrjgP1c3wd59mW33bjRzEx8MADEBcHM2fCkWKvwfatDh2gRw946y24di3Qo8lbeDhMmACz\nZwd6JAUbOBDWrAn0KExeli+HoUMDPQr/WrFCb/LGjQORQI/GPzZvhrlzoVUrmDgRqpWB1URffKGB\n9Pr1MGSIBtOtW5f8vBZUm7JmITDFUwWkO3DROXcSWA4MEZFqngWKQzzbTAA0bAgPPQSJifpin5IS\n6BHdrlYteOQR+OAD+OabQI8mbxUrQp8+sHRpoEeSv8hIyMwM/hn1smbTJp2tLCuBZXo6vPMOtGyp\nN8xlwY4detNdt64G07VrB3pEvnXxor7D+P77UKmSBtLDh0O5ct67htWpNiFFRGYC/YEaIpKIVvSI\nBHDOvQF8BIwADgLXgGmex5JF5DfAVs+pnnfOFbTg0fhBr176x27xYqhSBQYMCK4/8hERMHkyrF4N\nJ07oeINNXBwkJ8OWLXDPPYEeTd769IFPPoG+fQM9EgN6E3vyZNkJLk+d0hnqiRPLRqrL11/Dl1/C\n3XfDgw8GejS+5Rx89plO0FSuDCNG6I28r4gL5qRAY4JcYmKis+Yv/nHqlOY5duoELVoEejS3+/pr\n2LVLZz/CgvA9wI8/1gC7WbNAjyRv8+Zp6k8wWu55z6qspELMmQP33+/b4CNYbN+ury0jRwZ6JL53\n+LDeXLdooSlsoSwpCdat06C6WzeoX//Ox2RLSEjYHh8f36U41w3Cl35jjLldnTqaEnLlir5lGWyL\n21q21Ny8t97StxmDzYABsHOn/rEJRo0aBW8aTVmyezc0bx76AbVzsHCh3gCHekB96pS+Zp46pTPT\noRpQZ2XBhg2a3rFrF4wdq5McRQmoS8rSP4wxpUrnzvpHYelSTb8YOjR4UkKqVIGpU/VFvW1bXfgT\nTMaNgxkz9OYkKirQo7lVly46W924caBHUnZlZmowMmlSoEfiWykpGmSOGKFrI0LV+fOwcqU+x4kT\ng+d10ttOnoRPP9Xn16uXppMFigXVxphSJzwcRo3SXOHZs3XVdvv2gR6VCgvTqhsbN2qu9aBBgR7R\nTSIaUM+aBY8+Gnx/ZGNitJpK+fKBHknZtHix/r8KZYcPawD2yCN6Ux6Krl7VSYdKlTSlKhjT0Uoq\nM1PTO86f13cxg6VaS4j+ShljyoLq1XVW7auvtATfkCEQGxvoUalevfQP+Lvv6luuwfIHPCpK3+5e\nsED/EAWTQYN0Zm306ECPpOw5elT/P1WsGOiR+M7GjdrIZvLkQI/EN65fv/kO3tixoZnCc+QIbNum\nEyt9++rvbDAJkpd5Y4wpvrZtoU0bXcF//boGjeHhgR6V5gnXrq2B9fDhwfNWc2wstGunzR769w/0\naG6Kjtafn3PBMetUVjin1VdCNdjMyoL58/U1IthSsrwhM1M7qKal6etMqFUwSU/XWvZXrkCDBsEz\nK50XC6qNMSFBRPOrr1zR3Nz4eOjaNdCj0pSGKVP0rfW4OO3IGAyaNdP0mZ07gyd1BrTs35YtZaej\nWzAItjQlb7p4UQPqceN0zUMocU7byF+8qK99lSoFekTedeCAvj5FRurNf+XKgR7RnVlQbYwJKRUr\narrFgQPw3nv6YnzXXYEdk4imNGzfrsH1yJHBMdNyzz3w0Uf6FmqwVIasX1+DauMf58/rDGcoNv7I\nLnM5dWpo5RVnv7Nw8qTeDAVLyps3pKbqrHRKit74l7ZW8RZUG2NCUrNm0LSpLmbZvFkXYAW64kXn\nzlrWasYMDfyD4W3aESM0H3306ODJp73rLm2mE+ibobJg6VJdvBpqVq/W/18TJgR6JN61dSskJEDv\n3oGtcuFtX30F+/bpz2zAgNK7WDmE7t2MMeZWIjpTPWKE1qTdsCHQI9KV6pMna9m9o0cDPRo1caKm\nzGRmBnokqkcPrdBgfGvLFm2MEQzvmnhLRoa+QxUfH5wdTotr1y6t2lO9ut6Q16sX6BGV3JUrsGiR\nvhZGROis9MiRpTegBpupNsaUAdHRWlrq6FGdle3ePbD1kCMj4eGHtdLFiROBzx8OD9cZvTlzgmPW\nUkTfVUhLC/y7C6EqNVVbNwdr6/riOHNGZ94nTCjdgVlOBw9q2ljbtqFRP9w5+PxzrYxUoYJWbCpX\nLtCj8h6bqTbGlBkNGmjQeOqUBpDXrgV2PIMH6+Ki99/XPzaBVKGClqj66KPAjiPboEH6Fr7xjYUL\nYcyYQI/Ce3bs0Jn3KVNCI6BOTNSZ6YsXdWa6TZtAj6hkLlyADz7QRaPVq+us9LBhoRVQg81UG2PK\noB49tIPfkiWaRzxoUODeAm/TRnOH//UvnU0P5Ar+evW0jfnmzYGfPa9QQZtYGO/bu1fXG4RKHeMl\nSzStKhTqm589qwv14uJK/8y0c/pacuIEVK2q61qCpV6/r4T40zPGmLxFRsJ99+lbxrNna+vzli0D\nM5Zq1bRCwbx5Oo7mzQMzDtDyemvXavWUZs0CNw7Q8oNffBE8ZQhDQWamzuoGQ5pPSaWm6jtOgwdD\n3bqBHk3JXLoEy5frLG5pbyl+9qyuX3FOb8579Aj0iPzHgmpjTJlWq5bOCH3+ub7dOnx4YOrZhoXp\nH9P16+H4cV0BHyj9++vbtNWqQY0agRtHkyZ6o2FBtfcsWRIarciPHtXA7aGHSveMe0qK5oFHR2st\n7WBoWlUcWVn68zh3DmrW1AmLUCpjWFgWVBtjDNCpk84SL12qs0TDhgXmj0LfvnDokC6onDgxcH9k\n778f3n5b8zkDmfcYG6t/qAMZ3IeKxES9YSztTUI2b9aZ3YcfDvRIii89XV9rIDjKfRbX8eP68wgL\n0zJ/NWsGekSBZUG1McZ4hIVpSacLF/Rt5ZYtNdD2tyZNNEf07bf1D24gAkoRncGfOVMXfwXq7eh+\n/XSB07hxgbl+qHBOa7aX5kDUOViwAFq0CHzOf3FlZcGKFbpIetiw0rmoMiNDf5cuXND1IMHcNtzf\nLKg2IUVEhgGvAOHA351zL+Z6/GUg+4318kAt51xVz2OZwC7PY0edcyG0Nt4URdWqGlDu3q1B5b33\n+n8GpkIFzbP+8EOtuRuIVuJRUbr4a8GCwAW1YWH6BzsjI/QXOfnSxx/DwIGBHkXxXb6sVXLGjtW0\npNLGOf0ZJCVpDnjVqoEeUdF9842W94uI0BSx0vgcfM1eokzIEJFw4HVgMJAIbBWRhc65Pdn7OOd+\nkmP/HwE5szVTnHMBmJc0wapNG2jdGlat0pmlkSP9G9iJaG7i1q36VvHw4f67drbq1TWg//jjwOV5\nDxigiyfvvTcw1y/tLlzQ39/Suphv/35dsDplSunM0/30Uzh2TAPR0tYOPi1Nq5FcvQqNGmkpPJuV\nzp8F1SaUdAMOOucSAERkFjAW2JPP/g8Bz/ppbKaUEtGZpatXdfFew4b+b5jRtavmLma3N/d3jnPT\nppCcDF9+CXff7d9rg86IXbjg/+uGio8+Kr3l2dau1UD6wQcDPZKi++ILvSHo0QN69gz0aIpm3z5t\nHR4VpTe1FSsGekSlQym85zMmX/WAYzm+TvRsu42INAQaA2tybI4WkW0isllE7svvIiLyuGe/befO\nnfPGuE0pUKGCLhysUUNTQhIT/Xv9evU0MJo7VwNsf+vWDU6e1Bm3QGjVCvbkd3ts8rV1q9ZkL20z\nvJmZWo2nfn1dvFua7N2rZTorVNCbgQYNAj2iwrl2DRYv1jSb9HSdlR492gLqorCZalNWTQLmOecy\nc2xr6Jw7LiLxwBoR2eWcO5T7QOfcdGA6QGJiYoD74Bl/a9JEP9at0w5uo0b5b+Y4KgoeeQSWLdPA\nuls3/1w327BhekMxapT/K0i0aaPl9Vq39u91S7Pr1+HIEW0qVJokJcGiRdpuvEKFQI+m8L75Bj77\nTBc4l6aZ9Z07tS59TIw2woqJCfSISi8Lqk0oOQ7Uz/F1nGdbXiYBP8y5wTl33PNvgoisRfOtbwuq\njQGtSnH9us7s1Kzp39m0YcNg1y6tijF2rH9zHB98UNNQHn3U/+X+KlfWts2BqCNeGi1cqL8fpcmu\nXVpScurU0pO7e/Kk3mQ3blx6gunLlzVXOiND10yMHx/oEYWGUvaGkDEF2go0E5HGIhKFBs4Lc+8k\nIi2BasCmHNuqiUg5z+c1gF7kn4ttDKAz1OPHa3WOmTM1GPCXdu2gTx946y3/tvMOC9MZxNmz/XfN\nbAMH6oJJc2dff61BXmmqf7xsmTZDue++0hFQJyfr/4P9+zWY9vdai6JyTtOB5s+HTZt04fP48YHv\nnBpKbKbahAznXIaIPAksR0vq/cM5t1tEnge2OeeyA+xJwCznXM7UjVbAX0UkC73ZfDFn1RBjChIX\np53dtmyBbds0PcIfb1vHxuqM8dy5upixSRPfXxP0ufXvr935Ro70zzVBK69kZmqt39KWI+xPWVna\nIXTy5ECPpHDS0jQ4HThQ1w4EuytX9AagcmW9wQz238XkZJ1Jz8rS/PquXQM9otAlt8YVxpiiSExM\ndHFxcYEehgkiGRkabMbEaNUQf824ZVdJ8Gcayldf6dvIPXr475pnz+qCxX79/HfN5cv136FD/XfN\nkli8WN/FKA1pMsePaxrCxImB7dxZGNeva2nLyEj9XQjmuunOaSm/U6e0rne/fqW3Bbq/JSQkbI+P\nj+9SnGOD+FfCGGNKn4gIzWM9e1arF9x9t38W1/Xvr29Dz57tv9mztm11BmzfPu1y5w81a+r31uTt\nxAl9J6E0BNRbt+qixEcfDfRICpaRoTPTGRm6niE6OtAjyt/p07Bxo37eowf06hXY8ZQ1FlQbY4wP\n1KypKSE7dmhwPWyY7zuQNW+uDT5mzPBf57l+/bTjYmys/9qpx8dr/rq/0l1Kk7Vrgz/tw7mbnUKH\nDQv0aPKXlaWNny5d0plpf1e8KazMTFi/XtM8atfWnPRgT0kJVRZUG2OMD3XooKvrly3TYGL4cN/+\nwatUSSsnLFigs8dt2vjuWtnuuw/eftt/jWk6ddJauhZU3+rjj/2bFlMcV6/qGoBRo/x3E1ZUzsGG\nDZo6ce+92lU0GB07piX8wsI03SdYv59liQXVxhjjY2FhMGKEloObO1dnlDt29N31RGDcONi8GVas\ngCFDfHet7Os99BC89562kvZHHnlMjAZopamOsS9duqQL6IJ5od+hQ5ryEYhyjIX12Wdab7pPn+Bs\nOpOeru9GXLqkC6THjSsdlVLKCnuDwBhj/KRKFZ3NjY7WEnxnzvj2et27ayfCd97RP8a+FBkJY8bo\nDLI/DBoEq1f751qlgb8rsRTV+vVw9Kh2BQ3GgHrnTl2PUKOG/h+9665Aj+hWhw7p/62lS7Xp0/jx\nWsLPAurgYjPVxhjjZ61aade11at1dnHUKN9VEqhfXzvqzZyp1Ujq1vXNdUBzuDt21GoOAwf67jqg\naSZpafpWfVkPLLZv1+97MObRZmXpuzMdO+o7NMFm/34tP9i+ffA1brl+Xf8vXbum+efWoCX4WVBt\njDEBIKL5mteuaf5zXJzvStNFR2taxpIlUKcOdO7sm+uA5jknJ+sCzQ4dfHcd0O/X5s3+LekXbNLS\nICFBK74Em/PntevnAw8E3yK/o0e1AUqzZjp7Hkz27IG9e/XGccAAS3EqTYLwvtYYY8qO8uU1IKpT\nR2eTjx3z3bVGjtS33hcu1BleX+naVUt7HT3qu2uA5g+fOOHbawS7hQs17SbY7N6tCycfeyy4Auoz\nZzTN49gxnZnu1CnQI1JXr2p98ez0qfHj/ddEyniPzVQbY0wQaNxYPzZs0M6Mo0b5ph5uhw6aLzpj\nhgbz5ct7/xqgJchmzdIygpUr++YaoM8lMVFn+sua/fuhQYPga5qyYoX+zMeNC/RIbrp4UcdVo4Y2\nmgmWlKEdOzRfunx5Tc8Ktp+lKRoLqo0xJoj06aNv6S9erKW8+vXzfgBQqxY8/DDMmaOpE40be/f8\n2SZO1OD9kUd8lzPevbvO7j3wgG/OH6yysmDbtuCqSZ2errPAfftqsB8Mrl3TxX3ly2uQHwyLJC9e\n1AoemZnaHMpypUOHBdXGGBNkoqI0ADhxQlNCunbV3E9viojQgGz1ar2OLzqvhYVpYD17tgbxviCi\ns3vXr5etWb6lS7VMY7A4dUpngidODI6Og2lp+j0SgdGj9f9UIDmn5foSE3UWf8QIrZhjQovlVBtj\nTJC66y4NfM+f18D0yhXvX2PQIO2GOHeuzn56W/nyWglk8WLvnztbWSuvd+qUBq6+7tBZWNu368eU\nKYEPqDMzNZhevFjTKcaMCWxAnZSkC5Hnz9f/z+PH67gsoA5NBc5Ui0hVYLJz7v9KeiERKQfMADoD\nScCDzrnDeex3GLgMZAIZzrkunu3VgdlAI+AwMNE5d76k4zLGmGDXrZtW7FiyRIOWwYO9mxLSsuXN\n9ub336/1tL2pbl1NMfn0U+jZ07vnBg3cU1K8f95gtXp1cKR9OKfBa/36ga+R7ZwujExO1v8f3v4d\nLoqsLP1dP3NGU7jGjAmOtBPje3eaqa4K/MBL1/o2cN451xR4Gfh9AfsOcM51yA6oPX4BrHbONQNW\ne74uEhGpVtRjjDEmGISH6x/nTp10AeBXX3n3/FWq6EzjqlVazsvb2rTRWcSvv/b+uUHrIG/f7ptz\nB5N16zRnOdAL7VJStDX9Pff4vnTinWzcqOsD2rXT3PpABdQnT2p+/wcfaLrWuHHQv78F1GXJnYLq\nF4EmIrJDRP5YwmuNBd7yfD4PGCRSpJeFnMe/BdxXjDFsE5F3RWRgEa9tjDFBoUYNbQmelaX51snJ\n3jt3WJi+PZ2c7Jt0ij59YN8+33SSjI/X9tKh7PJluHBBZ4YD6fBhTWmYPFkXvQbK9u2aFtWggZbH\nq1nT/2PIzNQGLe+/r79/48bpR+3a/h+LCbw7LVT8BdDWOZfnfaiIbADyqkD5c+fcqlzb6gHHAJxz\nGSJyEYgFzuXazwErRMQBf3XOTfdsr+2cO+n5/BRQnF/Z5sBw4EngdRF5G/iXc66MVzo1xpQ27dvr\nzNyyZZCRoQufvDUj1quXBk7vvqvBijcrd4wZo23TJ0zwfv5tzZpw9mxggit/WLxYFwIG0saNkJoa\n2PST3bth1y5NiQpUF8QjR7T6Sni4vnNQvXpgxmGCS4kWKjrn+njSNHJ/5A6oi6K3c64TGvz+UET6\n5nFdhwbfRR1vpnNusXNuHNAXiAeOiki3EozXBBkRGSYi+0TkoIjcliYkIo+JyFnPOzA7ROQ7OR6b\nKiIHPB9T/TtyY4pGBIYP17J78+Z5N/2hUSOdcXv3Xe/OLItoB7vZs73fgKZPH63zHYq++ELLrwUq\nlSArS3/HqlfXhaGBkJCgvzeZmfo75O2KOHeSng7Ll+us9Nmz+v/jvvssoDY3lWj+oYgz1ceB+kCi\niEQAVdAFi7dwzh33/HtGRBYA3YD1wGkRqeucOykidYHbXuZF5J9AR+AE8D1gkeehN5xzb3j2qQJM\nAh4D0oBvATuL8rxN8BKRcOB1YDCQCGwVkYXOuT25dp3tnHsy17HVgWeBLuhN23bPsbYg1gS1ypV1\nxm7fPnjvPa22UadOyc8bE6N51tmL0byVOxsZCWPHer++dFiYfmRk+K4udiCkp8OBA4Gbpb54UdM9\nfLGItTBOnNCbpfj4wMxMHzigM+MREdo2PJg6RJrgcqeXncvkHTQDOlNdhGstBKYCm4AHgDWeGecb\nRKQCEOacu+z5fAjwfK7jX/T8+2Ee45mWa9MtfwJE5B2gBzAXmOKcO1CE8Rs/EpHfO+f+407b8tAN\nOOicS/AcMwvNx88dVOdlKLDSOZfsOXYlMAyYmd8B27dvp36gExyNMWVCoFIdsk3L/RfWmBB06NCh\nYh9bYPqHcy4J2CgiX3lhoeKbQKyIHAR+iqd6h4jcJSIfefapDXwiIl8CnwFLnHPLPI+9CAwWkQPA\nvZ6vi2oO0MI59wsLqIPe4Dy2DS/EcTdy9z0SPdtyGy8iO0VknohkR8WFOlZEHheRbSKyrRDjMcYY\nY0wZcMc3yJxzXlmO4JxLBSbksf0EMMLzeQJwdz7HJwElyuRyzi0syfHG90TkCbSMY7yI5EzLqQRs\n9NJlFgEznXPXReR7aDWZgYU92LN4djpAYmJi7jdcjAkqhw9rzdxevaBhw5KfLz1dy5f16eO9VtQr\nVkCLFt4ZH2hJwKwsLePnDcuX679Dh3rnfIXlnFZ48feiwIwM37ewz0tqqjZuiYrS77W/UniuXNEa\n12lp+jvTsqV/rmuCU0JCQrGPDaGsMxMi3gOWAr/j1lrkl7PTMu4gO3c/W5xn2w2eG7Rsfwf+kOPY\n/rmOXVuYQRsTrBo10o+NG7VN8qhRmitdXJGR2nJ85UrNde3eveRjHDJEF6BVq6b54SXVqpXma3sr\nqA6UpUsr9JZEAAAgAElEQVRh2DD/XvPMGb3uhAnaVMcfMjJurWLjj3bzzsHnn+tNZ8WKGsQHupW5\nKf0sqDZBxTl3EbgIPFTMU2wFmolIYzRIngTcMs+TveDV8+UYILvVxXLghRxNgoYAvyzmOIwJKr16\n6UzckiXa3rp//5I1EBk8WEubzZ+vC9hKWvl/wgTt6PjII96ZoaxSRWs6B0sr76I6fVpvYPxZWWLH\nDjh2TBen+qOTQ1aWNhu6fFmD2ooVfX/NCxe0gU5GhjZS6tzZ99c0ZYcF1SakeGqgP4kGyOHAP5xz\nu0XkeWCbJwXoKREZA2QAyWglGJxzySLyGzQwB3i+kLPjxpQKUVEaAJ88qV0ZO3eG5s2Lf742beCu\nu+Ctt7RpTEmqIoSF6UK8WbM0sC6pAQNg4UJ9vqXRqlX+Tfv46CNtWDJ6tO+v5ZwGtmfO6M1ZNR/3\nOnYONm/Wd1aqVtWW6qFUHcYED7F8UGOKLzEx0cXFxQV6GMYUy7ZtcOiQvuVekoA4K0vTLe6+u2RB\nOsCpU7B1q3eCu/nztY5wWIk6Mvg/p3rDBs0v91bOekFSUzV/esgQ75RhvJMtWzTlom9fqFvXt9c6\ne/Zm3fJ77oF6eS1ZNyaXhISE7fHx8V2Kc6zdqxljTBnVpQt07KizlJGRGjQW523/sDBN31i/Ho4f\n11ni4qpTB5o00cWVPXsW/zygTXHWr9dUl9LiyhVIStKFoL529KgGnQ89pD9/X/ryS/j6a+jWTQNc\nX8nKgk8+udlZ0xs3VcYUlgXVxhhThoWH66xwcrIuFmzTRtufF0ffvjrzPXOmNiopbve/1q01MNq7\nVxcdFldsrAaopcnixXqD4mubN8OlS7ro1Jf27dNc7fbtfVtn+/hxfU5hYdC7t/4uGuNvdv9mjDGG\n6tW19TNoUFzcYLRJExgzBt5+G86dK/54evfWTnanTxf/HNnjOXiwZOfwlx07oG1b37Yid07TYipV\n0pQPXzlyRG/Srl3TYLokN0f5yciA1av1+Rw9qm3D779fZ6iNCQSbqTbGGHNDu3Ya2K1YodVCRowo\nepBXoQJMnaoLBePjiz/zPXo0vPuuLoIsbhnADh1g3jxo2rR4x/tLerrO6vpyNvfyZc19HzvWd4sD\nT5/Wms8NG/ruuXzzjZbDCw/X1J7SWuHFhB4Lqo0xxtxCRPOrL1/WgDQ+Hrp2Lfo5xo7VRYdLl8Lw\nwvRDzeMckyZpYF2SMm8VKmiusj9KthXXokU6w+8r+/frTPiUKb7JMb5wQW/EatbUYNrbJfnS0mDN\nGp35btRIZ6X9UfbPmKKwoNoYY0yeKlXSAOnAAXjvPV2AWNSKDV27ar7rjBl6rqI29oiI0MVm8+YV\nP9d40CBtLuLLoLUkEhJ0gWZJmvIUZO1aDaQnTvT+ua9d05umChXggQe8H7Dv2wdffaXlIAcMCO4b\nI2MsqDbGGFOgZs00fWLdOti0SbsyFqX7XL16OuM8Z44GRkUtbValilYqWbUK7r23aMeCjjU9XfOJ\ng2120zn9nvpiwWBmJsydqzc2TZp499xpaRpMh4XpzYo3q4ekpOis9PXrWqJx/HjvndsYX7Kg2hhj\nzB2JaP5qaqqmKtSpo10aCysqSpu6LFumM9fduhXt+o0ba4WSzz/XTnhF1bOnlukrypj9Yfly39S/\nTkrSn9OECTqL7C2ZmTrm69e1hbo3Z9d37tRFpdHR+u5CdLT3zm2MP1hQbYwxptCio3Xm8OhRrRLS\no4fmuBbWsGGwaxd88IHmXBdl5rhzZ1i5UpuHFOWaoGkrn35atGN87exZff41anj3vLt2aUrJ1Kne\nm5l3TittnD+vNwGVK3vnvJcv68LG9HQtuzdunHfOa0wgWFBtjDGmyBo00I9PP9XFiCNHQvnyhTu2\nXTttbz5jhubhFmUmdfBgTSOpVk3TQooiLg6OHYP69Yt2nK+sWOH9VuTLlml5xLFjvXM+52DjRn13\nYdAg79wAOKfdPI8d07z94cN933zGGH+woNoYY0yx9eyps4yLF+vs5cCBhZsdjY3VdJDi5PxOmKAB\n+cMP60LGwurWTWsaB0NQvXGjpqJ4ayY5LU3rQg8c6L123Nlt7Hv10rrhJZWcrB0uMzM1R76oFWWM\nCXYWVBtjjCmRyEhtunH6NMyapa3PW7a883Hh4bqAce1anQktbBc8Ea1kMWuWBtaFDUxFtPpIampg\n83WvXoUzZ7yX3338uC7smzix6NVV8rJ7t1bc6Ny55LWmndN3M06d0hn00aN929zGmECyjorGGGO8\nonZteOghLbM2a5a2wS6M/v114ePs2ZCVVbhjYmI0FWTRoqKNcdAgzQ0OpEWLNLj0hq1bNYf60UdL\nHlAfPKg/N+c0mC5Jw5zTp/VdgQ8+0Hchxo/Xyi8WUJtQZjPVxhhjvKpTJ+1kmF1ybejQO9cvbt5c\nFxPOmFH4jn+1a2u5v08+KXx6QkyMzlQHys6d0Lp10dJW8uIcfPihNuYZNqxk5zp+HDZs0OA3u1V9\ncWRm6nmSkvRnc999vmk0Y0ywsqDaGGOM14WF6eLF8+d1YWGrVnD33QUfU6mSVqz44AMNstu0ufN1\nWrXS/OQ9ezRYLYzOnTVfuEuXwu3vLRkZOs6SBK6g6SPz5mm98NjY4p/n3Dmdta9Xr2RdEBMT4bPP\n9Gfep0/JxmRMaWb3kCakiMgwEdknIgdF5Bd5PP5TEdkjIjtFZLWINMzxWKaI7PB8LPTvyI0JTdWq\naRAZEaFdGc+dK3h/Ec3PvnxZq2MURq9euqDu1KnC7d+oERw5Urh9vckbaR+HDul5Hnmk+MHr5cu6\nQHTHDs3D7t276AF1erqWN3z/fZ3pvv9+nZm2gNqUZTZTbUKGiIQDrwODgURgq4gsdM7tybHbF0AX\n59w1EXkC+AOQvRQnxTnXwa+DNqaMaNNGZ5JXrtT0ixEjCk6B6N5dS669847Oot6p5Nro0brv+PGF\na0hSq5bm/dauXbTnUVyHD+s1S9KIJbtyRnFnulNSNCUnOlrrQRcnv/nQIQ3GIyOhX7+ilzU0JpTZ\nTLUJJd2Ag865BOdcGjALuKVaq3PuY+fcNc+Xm4E4P4/RmDJLBIYM0bJv8+drykBB6tfXOtYzZ8LJ\nk3c+/6RJhV/s2Lu35mL7g3N6reJW+8jK0udVp44u9iuq9HSd3V6xQlNyRowoWkB9/boG4/Pn6yz3\n+PHamtwCamNuZTPVJpTUA47l+DoRuKeA/b8NLM3xdbSIbAMygBedcx/kdZCIPA48DvDFF18QF2dx\nuTFFUbGiph0cPKgpIf365V9bOToapkyBJUs0qOzcOf/zRkRoGsK8eXr+gojo/unpvm88smJF8VuR\nnz+vCxLHj9ec86LIytJ3Bq5c0cWMRZ0l37MH9u7VqiIDBni33bkxociCalMmicgjQBegX47NDZ1z\nx0UkHlgjIrucc4dyH+ucmw5MB0hMTHR+GbAxIahpU/1Ytw42b9aFd/mVhRs5UtMOFi7UVI/8coCr\nVIF77tFgcvDggq8/cKDWdy5uwFsYSUka3NasWfRjd++GffuK3m7cOf2enjmj34PCVFLJdvWqtg1P\nS9NFoOPHF33cxpRVFlSbUHIcyNkrLc6z7RYici/wK6Cfc+569nbn3HHPvwkishboCNwWVBtjvKtf\nP82zXrxY84779Ml7vw4dbrY3nzAh/7boDRtq977t2wue2a5USWdxfWn5cq3dXVQrVmiHynHjinbc\n5s1w9Kh+T/v3L/xxO3ZovnT58hqIe6OJjDFljeVUm1CyFWgmIo1FJAqYBNxSxUNEOgJ/BcY4587k\n2F5NRMp5Pq8B9AJyLnA0xvhQdLTOijZqpDnUCQl571erlnZR/OAD+Oab/M/XsaOmThS0D0Dbtto8\nxRc2bdIFl0WZZU5P1wWXLVvqsYW1Y4fmXdetq6kvhVmAefGippbMn683GOPHw/DhFlAbU1w2U21C\nhnMuQ0SeBJYD4cA/nHO7ReR5YJtzbiHwR6AiMFf0L91R59wYoBXwVxHJQm82X8xVNcQY4wf16+vM\n7ubNWkt61KjbZ6QjImDyZK2xfOJE/gsA771XS8dVqwZVq+a9T4sWWhauXTvvPo+UFB1bjx6FP+bU\nKZ2hnjix8G3Uv/4avvxSa4AXpqW4c7pA9PhxnQkfMcL3OeXGlBXinKWEGlNciYmJzhYqGuMbGRm6\nQLFCBW0vnteM79df60zz+PF5d+9zTtNFJk/OP3hcvVq7QOaXe7x8uf5blNzr2bM1daOwAev27Vri\nb8SIwu1/+DBs2aI3BR0KUQg0KUlL8mVlac65vWwZk7eEhITt8fHxxWoNZekfxhhjglJEhLYsb98e\nZs3SShS5tWypZfpmzNB0htxEdAZ31iwNsPPSvz+sXeu9cX/1lQa7hQmondNyd+HhhQuoT53SgP3U\nKX1eBQXUWVlaym/+fB3TmDF682EBtTG+YUG1McaYoFarlqaEpKZqcHzhwq2PV6miZfdWrco78I6O\n1lnmRYvyPn92zebMzJKPNTNTA9jCzB6npMDbb+vM8Z32z273/vXXmh5SUL71yZMaSH/4ITRrpjPm\n/foVr9mLMabwLP3DmBKw9A9j/CsrC5Yt0xne4cNvT/nYuFGD70GDbj/266/h7Nm8q4skJ+tiv4ED\nb3+sKOkfH36o165YseD9Dh+GTz/VALmgzpJXr2rjlUqVtCpHXikuoMH8unUafNetq7ncRW09bowp\nWfqHLVQ0xhhTaoSFaZrEhQs6c9uy5a2zvL16acD67ruaHpEzYG3ZUoPn3bu1bXpO1atrQFoSR49q\nXvadAurswH/y5Pz3ye5imJ0Ck18qyZEjuqAzIgL69i1aTWpjjHdZUG2MMabUqVpV25Lv2aMl+AYN\n0jQR0LJ8tWtrt8Zhw25uB+jZU+thV6+uM7o5NWsG+/dD8+ZFH49zsGGDlvvLT1aWpmW0aaONVfKS\nmakz8WlpOhOfVxWQ9HRt0HLlCjRooOkdNittTOBZUG2MMabUat1aA9TVqzXIHDVKZ21jYuDRR7V6\nSFzcrbPZo0bpTPb9999arq99ey2vV5ygetUqLeGXn4sXYcECvWaVKrc/7pye4+JFTTPJqyX5gQNa\n6SQyUhdXFrVtuTHGtyyoNsYYU6qJaEB79aoGrvXr32y6MmqUlqtbskRbnWebNEkXCU6ZcmuecsWK\ncPly0QLW8+d19ji/hitff62LF3NfCzSY/uQTXVw4aBDExt76eGqqtlJPSbm56NAYE5ys+ocxxpiQ\nUKGCti+vVUtTP44d0+2dO+vHjBkapIJWwhg3TpvD5DRwoKZWFMXSpZqqkZfVq7VG9AMP3B5Qb92q\neeHx8bpgMWdA/dVXmiqyZg0MGKCl8Nq3L9q4jDH+ZTPVxhhjQkp8vH6sX68NUkaPhjp1tCzfnDla\n/aNBA+0o2KOHVvfIruwRGalNZ5wrXJ7y5s3Qrdvt+2Zk6LV69IDGjW99bNcuXSzZtat+ZLtyRQP6\n9HTNu7ZZaWNKFwuqjTHGhKS+fbWKxpIlOgvcr58uJFy5UluId++uwXVyslbQ6OIpotWzp1bo6N27\n4POnpkJi4u01o8+c0dnrCRNuzdk+eBA+/1wD5kmTdJtzuu3IEZ1pHzoUoqK89z0wxviPBdXGGGNC\nVrlyOuN7/LhWCenaVes979mj6RX336+LGNesgYQEneGuU0eD6jtZuBDuu+/WbTt2aNrJlCk3Z68T\nEzVvulkzTfMALQm4bp1W++jYUdNTjDGlmwXVxhhjQl69epr+8dlnmpYxYoSW1HvrLc1XHjhQK39U\nq6YfDRro7HHDhnmfb88eaNr01lnljz7SxYqjR+vXZ89qsF6vntbMBti0SWfJq1bVhZMFNX4xxpQu\n9t/ZGGNMmdGtG3TqpAFwdLSW3Zs/H+6+W2e0Z8zQpixduuj2vILqzEz48ksN0kHTQObMgSFDdJb7\n0iXN065eXWemz52DDz7QVI/u3TXP2hgTeiyoNsYYU6ZERMCYMRrszpmjVTVOndIUkUmTNE3k0Uc1\n6E5Juf34xYu1VB9oF8UNGzTAzsjQQDw6WtNCNm3SEn81amhXxPxajBtjQoMF1cYYY8qkGjU0GP7y\nSw2oW7TQIHjIEM2XHjpUG7LkbBGemKjNWypV0sofly/rbPTSpfp4x47wxRcaePfuDTVrBua5GWP8\nz+6bjTHGlGl33605z6dP6yz20qWaxrF1q6Z2OKf7OQdr12oVkfnztVpHZqamdoSH60z16dO6+PH+\n+y2gNqassaDahBwRGSYi+0TkoIj8Io/Hy4nIbM/jW0SkUY7HfunZvk9Ehvpz3MaYwAkL0wYuQ4Zo\ni/ONG3WhYpUqsH+/7rNmDdxzD/zrXxpEr1+v5fOioqBXL83Jzu7kaIwpeyz9w4QUEQkHXgcGA4nA\nVhFZ6Jzbk2O3bwPnnXNNRWQS8HvgQRFpDUwC2gB3AatEpLlzLtO/z8IYEyiVK2te9ddfa7fFnTsh\nLU0reCQlaSB9/To0b64LDjt2tCDaGKMsqDahphtw0DmXACAis4CxQM6geizwnOfzecBrIiKe7bOc\nc9eBb0TkoOd8m/w0dmNMkGjZEv7zP3Uh4zPPaDDdtKnOZk+YABUrBnqExphgY+kfJtTUA47l+DrR\nsy3PfZxzGcBFILaQxyIij4vINhHZdu7cOS8O3RgTTEQ013rJEmjXTgPsadMsoDbG5M2CamOKyDk3\n3TnXxTnXpUaNGoEejjHGx5o21SYxxhhTEAuqTag5DtTP8XWcZ1ue+4hIBFAFSCrkscYYY4wxt7Gg\n2oSarUAzEWksIlHowsOFufZZCEz1fP4AsMY55zzbJ3mqgzQGmgGf+WncxhhjjCnFbKGiCSnOuQwR\neRJYDoQD/3DO7RaR54FtzrmFwJvA256FiMlo4I1nvznoosYM4IdW+cMYY4wxhSEuu6q9MabIDh8+\n7CIi7N7UGGOMCQXXr193TZo0KVYmh0UDxpRAREQEcXFxgR6GMcYYY7xg2bJlKU2aNCnWsZZTbYwx\nxhhjTAlZUG2MMcYYY0wJWVBtjDHGGGMMUK1atbPFPdZyqo0xxhhjQlB6ejqJiYmkpqYGeihBJzo6\nmri4OCIjI2/ZXrNmzWK3Srag2hgv6/Z2H1KyLty2PSasKp89uqFMnSMYxgDQf3Z/klKTbtseGx3L\n2gfX+u0cxhjjT4mJiVSqVIlGjRohIoEeTtBwzpGUlERiYiKNGzf22nktqDbGy/IKAAvaHsrnCIYx\nAHkGwwVt98U5LCg3xvhbamqqBdR5EBFiY2M5e7bYmR55sqDamBL45vQ+fjl78q0ba+S//6N/6lu4\nE5eCczzyp355bnfkePEu4PiHXhqo+zvxHAcgN453iG6snf85Rr9wP86FkeXCcJ6PLBdGlgvHOSEr\nKxxHGDTN/xwtH38ZssLBhYMLg6ww/deFebYLzoUTNTT/czTruxMyoiE9Gpeu/5IeA1kR4Hk+Mf+V\nf1DetIDxFVe5cjBuHDzxBNx1l/fPb4wpHSygzpsvvi8WVBtTAunhGSRWzR0s5b/+9/Z98xP85zhe\ntTBpZ/kff7raqdu23f4S5yjoZepK7X039gJwuU5ws7VVeL7nKN/zTbKAzBK8wEZ/++E8t7vMcLLS\nyuHSogs8/lytRWRcqEnGxZqkX6hB1rXK5PXdaPFKfyKr3P6zS78Yy76n1962fc8eePFFmDgRnn4a\nunUr1NMxxhhTDBZUG1MCDSu14W/9dt+y7f6N7fLdP/e++QmVcxR0/PR+e0o8hpfvXklmRjpZWem4\nzDSyMtPJytTPXWb29nSePP8f+Z7jT5cHE5ZxjfDMa4RlXiU86xrhmSlEkEpEVgqRLpUouc6YhpXz\nPcfrp85wRSJJDivPhbAYLoeX52pENKnh0aRGRJEeHcWnJOd7fP3v/b9bvo4KK0e1yBpUi6rp+ahB\n9ahavHMk75uhyCpJHDhw67ajR+H//g8WLID33tOP7t3hqafggQcg19ocY4zxifDwcNq1u/k6PmnS\nJH7xi1945dw7duzgxIkTjBgxwivnKykLqo0pgXLl4LaGihvz37/Qb/OHyjl8PIYOXeoV7hxv5R9U\nD/jRy4U8R/7BfVS971Pt6jlqpiQRdT2ZmGvnqZiZSNWsi5SX6wC0a9wg3+Ofv96FS5Wqc61yRdIq\nRJMaDclpFziXco5TKQfZdWkTV9KvFDi8I1GraVm9JXdVuAsRoWlTGDgQjhyB11+Hv/0NNm/Wj5//\nHH7wA3j8cahZs3BP3xhjiiMmJoYdO3b45Nw7duxg27ZtFlQbE6piwqrmW62irJ0jGMYAUDWsPBey\nruW53Rvn6P7Ir/M9LuXqZS6cOwFrH8h3n5HHFxIlGTe+znLCWalOUmRdrpavR0blXmRWr8cTl/6R\n7zl+/PGPAagUWYnm1ZvToloLWlZvSfPqzfnN75ry7LPleOcd+HNGf6RiErOAWR/dPN4WTBpj/OXi\nxYt069aNhQsX0qJFCx566CEGDhzId7/7XZ544gm2bt1KSkoKDzzwAL/+tb6+bt26laeffpqrV69S\nrlw5Vq5cyTPPPENKSgqffPIJv/zlL3nwwQcD+rzEOXfnvYwxeUpMTHRxt01Vm2DzxrpDDExbS/Ov\nXoaLiVAljv1tf8KaqP58v18Tv5yj/z/bkhR2e550bJZj9aM7OHvyMMmJB7h6+hAZSYeJuHSUCteO\nE5t+kpoumTBxBc52/0/750mOTmP/hf3sO7+P/ef3k5KRAkC4hNO4SmNaVG/BkoQl+Z7j+cq7GDMG\nwvNPQTfGlCJ79+6lVatWAPhqveKdwsjc6R/ZwW92UPz000/zr3/9i2XLlgGQnJxM9erVyczMZNCg\nQbz66qu0bNmSli1bMnv2bLp27cqlS5coX74877zzDtu2beO1114r1thzfn+yJSQkbI+Pj+9SnPPZ\nTLUxJuRp0NsEBn/7xrbmng9/nWN63E+ov/EXxJB2Y1sKURzr9SLhERHUqd+UOvXzzom5nnqNM4mH\nYOPkPB8HGPzhY1xz5UiMbMiFik3JqDGSK3FxXIqN4UTWOfad38fWU1sLHOO4cdCoETz5JHz721C1\n8G8GGGNMnvJL/xg8eDBz587lhz/8IV9++eWN7XPmzGH69OlkZGRw8uRJ9uzZg4hQt25dunbtCkDl\nyvmvcQkkC6qNMcYP1kT1Z2CvF2+Z6T7mmem+U2BeLro89Zu2I3aDy3u2O9PxWbtfk3V6NxUv7qfJ\nhY3EXvgIDurjyVTmRLl4rlRpzncrnMn3Oq2em8bZrb34zz/34JlnWzF1ShhPPQUtW5bgiRtjgkKw\nJSZkZWWxd+9eypcvz/nz54mLi+Obb77hpZdeYuvWrVSrVo3HHnusVHWDtKDaGGP8wBuz5QXNdjfP\ncV6ApNOJnNy/nSvHdhJ2Zg9Vrxyk/ekPIT7/wt8tO1wivNErMOEVMi5VZ9Hu7sz8Tk+61OjJT79X\nk6FDISzMGtkYY0ru5ZdfplWrVrzwwgtMmzaNTZs2cenSJSpUqECVKlU4ffo0S5cupX///rRo0YKT\nJ0+ydetWunbtyuXLl4mJiaFSpUpcvnw50E/lBguqjTGmlCjKbHds7Thia8cBY29sy8rMhHc65Hv+\ndz9fz47oJmyOvYsvK0ews8tGqvb4iFPAj7Y1J2phT8be3ZOkmJJ3qDTGlA0pKSl06HDzdWfYsGFM\nmzaNv//973z22WdUqlSJvn378tvf/pZf//rXdOzYkZYtW1K/fn169eoFQFRUFLNnz+ZHP/oRKSkp\nxMTEsGrVKgYMGMCLL75Ihw4dbKGiMaWdLVQ0pU1+CyarZ2bx+9ReVEv+kkbphygn6WQBm6Ni+Si6\nHpvDIjhd+TyEZxZ4/l1Td/lo5MaYosprIZ65yRYqGmOMKbbCpJCkXU9l/+4tJO/bSNTJ7Xz/8m5+\n605y7aKwpVw0T9W14tbGGJObBdXGGFOGFCaFJKpcNM079YNO/W4cd/7sSY7u2kD0N1uAFfmef9qH\nUxnfbgID6w+kfGTh64AbY0xpZ+kfxpSApX+YsqhdAd0l62ZkcDIigkgndCzXlDGtJjG83X1EhUf5\ncYTGGLD0jzux9A9jjDEBlV93yei0cnRYMZH7WnzMyVpHWVfxa/7zy9/wu8+fp1NmDXrWHMDQHtMY\nueLhfLtkfvboBn88BWOM8ToLqo0xxhTJww3ey7O75OrI/tzdpgmvvPJfLPhDBh1bfkL3vjPJbPY5\n2yufZcPFefxj8SxSIvL+05NXoG2MMaWFBdWmzBCRYcArQDjwd+fci3nsMxF4DnDAl865/FvYGVNG\n3anmds+ecOxYBP/3f/2ZPr0/yckQWf4ivSb8japdl3O24ql8z+2yspCwMJ+O3xhjfMFeuUyZICLh\nwOvAcKA18JCItM61TzPgl0Av51wb4Md+H6gxIaJ+ffjd7+DYMZg+HVrEV2HtWz9nwZMrCzzu5PPN\n2PLnqexYPYuUq8HT1MEYU3yJiYmMHTuWZs2a0aRJE55++mnS0tJYu3Yto0aNyvOYRo0ace7cOT+P\ntGQsqDZlRTfgoHMuwTmXBswiZ1cM9V3gdefceQDnXP79nI0xhVK+PHz3u7BzJ6xeDWPGFLz/o/Ur\nsyN9HQ03PoH8oTFf/n4wW2b/npNH9vlnwMaUUW+sO8Snh24NYj89dI431h0q0Xmdc4wbN4777ruP\nAwcOsH//fq5cucKvfvWrEp03GFlQbcqKesCxHF8nerbl1BxoLiIbRWSzJ13kNiLyuIhsE5Ftpe0u\n2phAEYGBA+HDDwver06t9rwWW4lBDRvxo/rtOOOO0m3vC9T9ZzcOP9+OTX/9EXs2L+Mva75m/8o3\n4eW28FxVeLkt+1e+WeIAwJiyqn1cFZ5874sbgfWnh87x5Htf0D6uSonOu2bNGqKjo5k2bRoA4eHh\nvPzyy/zjH//g2rWbC56TkpIYMmQIbdq04Tvf+Q6lsTqd5VQbc1ME0AzoD8QB60WknXPultVTzrnp\nwFajFt8AACAASURBVHTQknr+HqQxpV1MWNU8FyVmXqvKgifeYvz3DlB7yBw+Ob+IzXWjaFy+N71S\nYrn3+Dd0OfEukSdn0MhFEUUGSJYefPEY9Tf+goG9XgSa+PcJGVMK/HrRbvacuFTgPrUqlWPKm59R\nu3I5Tl+6TtNaFXll1QFeWXUgz/1b31WZZ0e3KfCcu3fvpnPnzrdsq1y5Mg0aNODgwYM3x/frX9O7\nd2+eeeYZlixZwptvvlnIZxY8LKg2ZcVxoH6Or+M823JKBLY459KBb0RkPxpkb/XPEI0pG3KXzcvK\ngmXL4JVXYG8KvPO/zeB/f8Wg4T9h8LcWszdqDu9c+4IFDSswNO6HdL5YmT+d/SfJ4eG3nTs28WXW\n8u3bthtj7qxKTCS1K5fj+IVU6lWNpkpMpN+uvX79eubPnw/AyJEjqVatmt+u7S0WVJuyYivQTEQa\no8H0JCB3ZY8PgIeAf4pIDTQdJMGvozSmDAoLgxEj9GPvXnj1VZgxA1YvLc/qpRNp2nQCE3/8JWmt\n5rD4yEfMz0qDPAJqgKQwYfNfvk+tPtOIb3uPn5+JMcHrTjPKcDPl46mBTXlny1GevrcZPZvUKNF1\nW7duzbx5827ZdunSJY4ePUrTpk1ZsSL/Dq2ljeVUmzLBOZcBPAksB/YCc5xzu0XkeRHJXjq1HEgS\nkT3Ax8C/OeeSAjNiY8qmVq3gL3+BxET44x+hYUM4eFB44ckO/PW+F7hnxyqmNvppgefodGoO8fOG\ncPA3ndg88785f/akn0ZvTOmVHVC/NrkjPx3Sgtcmd7wlx7q4Bg0axLVr15gxYwYAmZmZ/OxnP+Ox\nxx6jfPnyN/br27cv7733HgBLly7l/PnzJbpuIFhQbcoM59xHzrnmzrkmzrn/9mx7xjm30PO5c879\n1DnX2jnXzjk3K7AjNqbsqlYNfv5zOHgQ3n8f+vaFy5fhL/9TjX8fMK3AY68+uZvNLf4DgO77/kCF\n19rwxR9H8sWKd0hPu+6P4RtT6uxMvMhrkzvemJnu2aQGr03uyM7EiyU6r4iwYMEC5s6dS7NmzWje\nvDnR0dG88MILt+z37LPPsn79etq0acP8+fNp0KBBia4bCFIaV1caEywSExNdXFxcoIdhTJmwY4em\nhrz3HjT7a7t893u609NMbjmZ8pHl+Wb3Fk6v/ydNTy+lBhdIpjL7aw2nZu/HWHm+dp6dIddE9fc0\nuDGmdNu7dy+tWrUK9DCCVl7fn4SEhO3x8fFdinM+C6qNKQELqo3xvzNnYNDS/INqgOrR1flW22/x\nYIsHiY6IJj3tOrs3LCDri3dpe/lToiSDk9Qg1p0nSjJvHJdCFMd6vUjzwbbY0ZR+FlQXzNtBtaV/\nGGOMKVVq1YLY6Ng8H8u4GEvDT96mljTnpW0vMXz+cN7d+y4uXOgwaBKdfr6Iaz/azZZWv6QGF24J\nqAFiSNOZa2OMKSKbqTamBGym2pjgsHmzpobMnQsZGbqtw6ht1HngNRLZTu3ytXm8/ePc3/R+IsM9\nZcKeqwrk9TdQ4Lnb62gbU9rYTHXBbKbaGGOMyaV7d821PnwYfvUrqFEDdizuwrLH/sn56X8jPbkO\nv9n8G0Z/MJoFBxaQkZUBVfK5Ic5vuzHGFMDqVBtjjAkZ9erBb3+rgfXMmfDKK8LOT7tz/NN7qNrx\nE2Iee51nrjzD33f9nXOx5bhW/fYKA1WI5pMAjN0YU7rZTLUxxpiQExMD3/qWVgxZuxbuv1+49GUf\ntj49kyOvvELi4WiuubQ8j71IKhfOnfLvgI0xpZ4F1cYYY0KWCPTrB/Pnw6FD8LOfCWEJA9n587kF\nHndi+gNcT73mp1EaE7pEhJ/97Gc3vn7ppZd47rnnAjcgH7Kg2hhjTJnQqBG89JJ2a3z9tYL//LVK\n28XOv0zFZWX5Z3DGBIOdc+DltrqI9+W2+nUJlStXjvnz53PuXMk6M5YGFlQbY4wpUypWhB/8oOB9\nRtZrT8b1dWz613/4Z1DGBNrOObDoKbh4DHD676KnShxYR0RE8Pjjj/Pyy7eXqjx79izjx4+na9eu\ndO3alY0bNwLQrl07Lly4gHOO2NjYGy3Op0yZwsqVK0s0Hl+yhYrGGGNMLgmXy/F43dp0TvmA/bNj\neOzB5wM9JGNKZukv4NSu/B9P3AqZ12/dlp4CHz4J29/K+5g67WD4i3e89A9/+EPat2/Pv//7v9+y\n/emnn+YnP/kJvXv35ujRowwdOpS9e/fSq1cvNm7cSMOGDYmPj2fDhg1MmTKFTZs28Ze//OWO1wsU\nC6qNMcaUSbHRsSSlJt22PSYrlqt/WUJy/CwOTPgT21MX8PZrX/GTHs8yqvPdARipMX6QO6C+0/Yi\nqFy5MlOmTOHVV18lJibmxvZVq1axZ8+eG19funSJK1eu0KdPH9avX0/Dhg154oknmD59OsePH6da\ntWpUqFChxOPxFQuqjTHGlElrH1yb72OZU2Dx4qn89W/9Gdd/OAuq7+OXXz3Cr1f04fvtfsi3RrZB\nxH9jNabE7jSj/HJbT+pHLlXqw7QlJb78j3/8Yzp16sS0adNubMvKymLz5s1ER0ffsm/fvn15/fXX\nOXr0KP/93//NggULmDdvHn369CnxOHzJcqqNMcaYXMLDYexY+GhxQ7rePZdZRy/x6LlMrlXZwf8m\nTaLtfz7Nb97YR79Z/Wn3VrvbPvrP7h/op2BM0Qx6BiJjbt0WGaPbvaB69epMnDiRN99888a2IUOG\n8Oc///nG1zt27ACgfv36nDt3jgMHDhAfH0/v3r156aWX6Nu3r1fG4isWVBtjjDEF6Du4DVeGv8mP\nL53mb6eFtDXfxcV9xpyYB0i+fnv6CJBnWokxQa39RBj9qs5MI/rv6Fd1u5f87Gc/u6UKyKuvvsq2\nbdto3749rVu35o033rjx2D333EPz5s0B6NOnD8ePH6d3795eG4sviHMu0GMwptRKTEx0cXHW0tiY\nsmDb4ul02fZvbK08hK/qTOe1LW+T1ulv+e6/a2oBi8KM8YO9e/fSqlWrQA8jaOX1/UlISNgeHx/f\npTjns5lqY4wxphC6jHqcTQ2/T9dLK2iZ+ju2v/JUgfunpvppYMaYoGBBtTHGGFNI3af+jq1VhtHj\n6F/ZtuivBe7boEkK//VfcOKEnwZnjAkoC6qNMcaYQpKwMO7+wVvsjmpH+23/r8B9q//bKF5f/z4N\nG2cweTJs2eKnQRpjAsKCamOMMaYIospFE/e99zkdVovqGXm3MS8XVpHmdepS71vP0fi58Sze+zHd\nuzu6d4eZMyE93c+DNsb4nNWpNqWOiPy0oMedc//jr7EYY/5/e3ceHVWV9nv8+1TmEEggjDIjRFDG\nMImMMrzyemkcGgRUWhyuoi9K061Lba/Da8sSFUVtua3SKqiNoqggFxBklEbBICAgKOMrRAEhBCKQ\nOfv+USUdJZUUVKikyO+zVq3UOWefzXOelYQnu/bZp2pKTK7H8dGzmT/jcuLJxVNszepsotnXcxKt\nBt7Csr3LeGbt88T88R5yd6ey6Z8TuP76jtx7L9x5J9xxB9SpU3HXISLlRyPVEo66AHcCDX2vsUAq\nUN33EhE55xq2uISYuOq/KqgB4sgjZcsUzIwBTQcwb9iHPHzpwzRo8z0XPjyaSx6cwOGiPTz8MDRu\nDLfcAl9/XTHXICLlRyPVEo4aAanOuZ8BzOwxYL5z7sYKjUpEqpwoP+tRu2Pp/FJrR3miuO6i6xjS\nYghvbn2TNyLfoPWk5VTb9Xsy63xKWo1MbtwIbPz3+cmxyaU+8VFEKh8V1RKO6gF5xbbzfPtEREIr\nsVHJj3Z2js1P9iO71VBS+l1PUu36xEfFM7bDWIanDOeVTa/wvr1PpCsosduMnAwyM6FmzXMcv4hP\nv1n9Snxokf7AC5ymf0g4ehP40swe841SrwWmV2hEIlIlbW87gWyif7Uvhyg2J/SkZt4Bum35b6r9\n7WI2TRrIlx+9yLEjh0iOS+Yv3f/C3Kvnltp3o0Zw113w7bfn8gpEvPw9BbQ8ng46ceJEUlJS6NWr\nF6NGjWLy5Mm/Ov7MM8/w4osvAjBhwgT69+8PwLJly7jhhhsAuPPOO+nSpQuXXHIJjz76KACffPIJ\nw4cPP9XPihUrGDJkCACLFy+mR48epKamMnz4cI4fPx70dZRFI9USdpxzE81sIdDbt+tm59yGioxJ\nRKqmZdH96N9zEilbpsCxdEhsxN62E/g8uh939G7Ozs2fc2jNLJoe+IQLvn6YvI2PsTG+K/mtr6J1\nv5Gl9n3yJPz9797XFVfA+PHerx4Nh8lZeOrLp/j2yNn9hXbzJzeXuL91rdbc3+3+Us/96quvePfd\nd9m4cSMFBQWkpqbSuXPnX7Xp3bs3zz77LPfccw/r1q0jNzeX/Px8Vq1aRZ8+fQBvYV6rVi0KCwsZ\nMGAAmzZtYuDAgdx+++2cOHGCatWqMWvWLEaOHMnhw4d54oknWLJkCdWqVeOpp57iueee45FHHjmr\n6w+UimoJS8659cD6io5DRKq2sX0vBC6EQbee2pfiewG07NCLlh164YqK2L7xM46sfZdmBz+l/oYH\nyV3/CLRo4LfvcfMeJnPpaD58JYVFi2DRIkhJgbvvhptuguq6LVvCwKpVq7jmmmuIj48HYOjQoae1\n6dy5M1999RVZWVnExMSQmprKunXrWLVq1akR7Pfee49XX32VgoIC9u/fz9atW2nfvj2DBw9m3rx5\nDBs2jPnz5/P000+zcuVKtm7dSs+ePQHIy8ujR48e5/xaVVRLlWFmg4EXgAjgH865SX7a/R6YDXR1\nzq0LYYgicp4yj4eU1H6Q2o+iwkK+/WoZR9e9B6zxe87aY4vI7jiHEe/3oMb2m/hoymVs327cfTc8\n9BDceiuMGwctWoTsMiSMlTWi3G5GO7/H3hj8RrnGkp2dTceOHQEYO3YsY8eOpXnz5kyfPp3LLruM\n9u3bs3z5cnbu3EmbNm3Ys2cPkydPJi0tjZo1azJmzBhycnIAGDlyJC+99BK1atWiS5cuVK9eHecc\ngwYN4p133inXuMuiD5GkSjCzCGAq8J/AxcAoM7u4hHbVgfF452mLiJQ7T0QErbsN4tK7ppFcUFhi\nm+SCQibXuo272o7l+xM7WVZrLKkvXstfZn5Ez765ZGXBlCnQsiVcfTUsXw7OhfhCRALQp08f5syZ\nQ3Z2Nj///DPz5s0jLi6OjRs3snHjRsaOHQt4p4BMnjyZPn360Lt3b15++WU6deqEmZGVlUW1atVI\nTEzk4MGDLFy48FT/ffv2Zf369UybNo2RI71Tqi699FJWr17Nzp07AThx4gTbt28/59eqkWqpKroB\nO51zuwHM7F3gKmDrb9r9FXgKuC+04YlIVbQiywPH9p62vxAjYt99dHaxdErsyboLu7PMfcvHeY9Q\n687nefj+Uez84Do+eKsWc+fC3LnQrh3ccw/ccAPExVXAxUhYS45N9rv6RzBSU1MZMWIEHTp0oG7d\nunTt2rXEdr1792bixIn06NGDatWqERsbS+/e3lunOnToQKdOnWjdujWNGzc+Na0DICIigiFDhjB9\n+nRmzJgBQJ06dZg+fTqjRo0iNzcXgCeeeIKUlJTT/+FyZE5/2koVYGbDgMHOudt826OB7s65ccXa\npAIPOed+b2YrgHvLmv6Rnp7uGjVqdA4jF5Hz2fZPX6Px6geIK7ZKaDbRfN9jIvlx9cje8B4pmStI\n4jhZLo6Panfm07pRfJ2zi5iIGAY2+B2ffr+EvIijp/VdMzqZz0atCOHVSGWzbds22rRpU9Fh/Mpj\njz1GQkIC9957b0WHUmJ+du/e/VWLFi26nE1/mv4hApiZB3gO+HMAbW83s3Vmtu7w4cPnPjgROW8t\ni+7Hvp6TILExYJDYmH09J7EidgDt+lxFt/H/pNpfdrOp72t8W/NyrsnYwNvblvN2+lG6n4xncfpH\nJRbUAJl5GYwYAZ9/rqkhIqGgkWqpEsysB/CYc+4K3/aDAM65J33bicAu4JeFLOsDR4ChpY1Wa6Ra\nREIpLzeHbavnkPf1B7Q+uoq8iFz6NfX/O2jLmM0AdOniXZLvuusgOtpvcznPVMaR6spEI9UiZycN\naGVmzc0sGhgJfPzLQefcMedcbedcM+dcM7y35JdaUIuIhFp0TCwd+o+k64T3iXpgF/t6TC21/d2P\n/A/JybBuHYweDU2bwuOPw8GDIQpYKpwGT0t2LvKiolqqBOdcATAOWARsA95zzn1jZo+b2emLZoqI\nVHKxcdXoOOj6Utssb/E7Br50NfdPW0TbdgUcOACPPgpNmsCYMbBeq/2f12JjY8nIyFBh/RvOOTIy\nMoiNjS3XfjX9QyQImv4hIhWttPWF7zhynDk1YjkYGUmtoihS6cuOhfezYHb9U/Ose/XyTg25+mqI\n1Jpg55X8/HzS09NPreks/xYbG0ujRo2Iior61f5gpn/ox0dERCSMJXniOVp0ssT9t93xBZd+9j7r\nd8xkTeT3LIlfQuSVnzJmQC2SfhzO9Bdv51//imFT3i4+e2sFf02dQg2XjiU2YnvbCSyL7ud7aqSE\no6ioKJo3b17RYVQZGqkWCYJGqkWkor28chf981aQsmUKHEsHPwVxzsnjLF3+Kkv3fcjnMZmc8Hho\nlVtAh5+bMb/GQbKj80/rO8kTz6rRehaWVB3BjFSrqBYJgopqEQlHR7J+4o2lk1hyZCXpkXmltn0y\neTNXXgke3YUlVYBW/xAREZGA1apRlz9f8xwLblnHtMtfLrXtNdfkkpICL7wAWVkhClAkDKmoFhER\nqaLMjEub9Cy1zZePtOS/rxjF8tnv0bxZDuPHw86dIQpQJIyoqBYRERG//tAkifntNzDmqj/yzT3N\n6VF4Hfff9i5XXXWSJUv0tEaRX2j1DxEREfHrz53/zOzts3k49ntikz0MTt7AXy9YRqPcP/LJ3F7c\n9fK1dOh/FfmtDnAFZd8wKXK+0o2KIkHQjYoicj7oN6sfGTkZp+1Pjk1mxYgVOOdY/9N6Ptj+AYu/\nX0xuYS4Nc2MZkXWYYScyoSiWwU3qkFXCUJ1WEJFwotU/RCqIimoRqWqO5R5jwZ4FfLD9A77L/I5I\nF8klmdX5ulam33M237Q5hBGKnD2t/iEiIiIhkRiTyKjWo3j/d+/zzv96h6svupoddXNLPWftnFf4\n+diREEUoUjE0Ui0SBI1Ui4jAyfyTdJ/Z3e/xzXv2kusi+Sa2CwWth3BRn+tITK4XwghFAqPHlIuI\niEiFiY+KL/V4lxPDGVx4nNuL1tHk6/9D/sZH2RTXkZyWQ2jZZwS16jYM+MmQIpWVRqpFgqCRahER\nr3Yz2vk9Fm1x5Lls8g5fQMKO9gyNOMKYWv+iacQBCp3xbUx78hIa0ubIEmL59xMes4lmX89JpAy6\nNRSXIKKRahEREalYybHJflcQWXDtApbtW8Z7W+axIXkxs62IN3d1psb2doxMSmdYnZX8V73DZCTW\nP/389CmsQEW1VH4aqRYJgkaqRUTOzE8nf+Kjbxcw8+uPOeLZQVFBJMe/7kONzsv8nvNx19dp1roz\n5tH6CnJuaUk9kQqiolpE5OxtO/wdU5fPY1XGfIriDvttt3nPXtKtAen1+pPU+VpSUi/HExERwkil\nqlBRLVJBVFSLiASvoKiATm918nv8pcjR1Ny3iNbZG4i2Qg5Rk93JfYlvfxUXXXol0TGxIYxWzmea\nUy0iIiJhK9JTejlyZ+b/o2Hk5dx4+QRSsvYSs30h7Q4vJH75HLKWx7OpxmV4Lv4daz2pDLA0rSAi\nFUJFtYiIiFRqP2/pwQ8dFjP5+w+xvAQ6Nu7DyE53kLwnA8+2hbQ8+i9qrllCW+chAgfm+xT+2D4a\nr36A/j0nASqq5dxSUS0iIiIVrrTVQ2b+11M8/7c85mxYQ2y7paSlLmPDFwuIdLH0atuL/2g6lXr7\n87hv52MciTj9ZsbkfVNYkj+ayKjoUFyKVFGaUy0SBM2pFhEJncOHYdo0mPr3Ao4mfEWNzktI6rqU\niMRDRFoUBS7f77mr92SwM6ELhS0G0LT7UOo2bB7CyCVc6EZFkQqiolpEJPTy8+HDD+HFF+HzL4qI\na7GJpK5LSB48w+85r2W2pdnRNdTlCAB7PM04WLcnCZdcQauug4iJ9T4VUk92rNpUVItUEBXVIiIV\nKy3NW1zPmgUXTfP/VMdu9bvRrV5XGuUkkLj9W5J+WEVKzhairYCTLobt8Z3IbXY5cZ4CWn3zAnF6\nsmOVpKJapIKoqBYRqRz274f/WOy/qL6o5kV8l/kdAHGRcaTWTaVDzXbUPZRDg12baZLxBY3cAfo1\nbkhG5OlrYCcXOVbcvOWcxS+Vg5bUExERkSqtQYPSj0e/O5vnxx3FNVzH2v1rSTuQxuofVwOQUCuB\nzm0GkRLZmIz/ebvE8zM8xhf/mEBc0640uqQntS9oWt6XIGFORbWIiIicF/ytIFJwLJl33oF33kmi\nW7eBjB8/kPuGQVbhYdIOpPHlgS9JO5DGyqyVpfbfdd90ItNfh9VwkGR+rNaGnLodSGjenSbtepJY\ns7bmZFdhmv4hEgRN/xARqfz27oWpU70rh2Rmevc1aAB33QV33AF16nj3HThxgEGzB/ntp2ViS+q5\n6tQ+nk+jYxl0OLaPrvk/nhqh3GcXMKpxJMdKeIJ6kieeVaPXlu+FSbnTnGqRCqKiWkQkfJw8CW+/\nDS+8AFu3evfFxMD118P48dChA7Sb4X9edq+Gvdh9dDc/nvjx1L4oTxT1I2pSPy+SxieO82FClt/z\np2X3xpIaE1O7GdXrNad2wwupkZSMeX69tna3t3qTXXT0tPPjPEl8OXrVGV61nAkV1SIBMLPBwAtA\nBPAP59yk3xz/E3AbUAAcAm5xzn1fWp8qqkVEwo9zsHSpt7ieP9+7DdC3L2Tc7L+o3nzTZgBO5p9k\n97Hd7Dq6i13Hdnm/Ht3FD8d/KPXfvS3zZ2oX5ZNUWETNoiKSCguJLoim0JLJiapHTrWGFNVoyO3M\nKzOGspRHYR5sH5UhhjPtQzcqipTBzCKAqcAgIB1IM7OPnXNbizXbAHRxzp00szuBp4ERoY9WRETO\nJTMYOND72rkTXnoJXn8dVq6Ei65NJirx9HnZ1T3JLF36y1Y80JYE2tIB6ACQBLk1TnL/j939/ruv\n1ayBo+TBzJiigyQV/kit3DUQ4//Jj//32XY4ooBoIAYsDizWG5MnHo/FYxEJZCecXkQCZBcdZeY/\nF2MWgVkk5okEi8Rjkd59nkg8FgWeiBIL0V/6WLDAe8zMTu13vveG51Q7f+cvWpzr9xp/2zaUfbRq\nGVCXJdJItVQJZtYDeMw5d4Vv+0EA59yTftp3Al5yzvUsrV+NVIuInB+ysmD6dPjb37yF9tlqO93/\nSPeWmzcSEf8zEdUziUw4SkSNI96vv2xXzyQm8TDx7T4/+wAkKHN7z9VItUgZGgL7im2nA/6HE+BW\nYGFJB8zsduB2gA0bNqCiWkQk/NWoAffcA+PGwYIFMHOm96bGvLyyzy3up1KO9b88AkjyvYrJ9r0O\n+fpo578wT35/GXjycZ58iMjzfj21nUNEVBZRUVkc6DnRbx8X7xgIVohRhJn3BUWYFZ7aNitiQ5ON\nfvvonN4GKz7qbt73xcatSWu43e/5XX9o5fdYcWkNd1SKPgKhkWqpEsxsGDDYOXebb3s00N05N66E\ntjcC44C+zrlSP1vSSLWIiBTXb1a/Epf1S45NZsWIFQH1UdrNkoHOqa4MfVSGGM60D82pFinbD0Dj\nYtuNfPt+xcwGAg8RQEEtIiLyW4EWzqWJ8yT5vbFOKi8V1VJVpAGtzKw53mJ6JHB98Qa+edSv4B3R\nLu0TPBERkXOmPJbNK4/CPNg+KkMM5dVHIDT9Q6oMM7sSeB7vknqvO+cmmtnjwDrn3MdmtgRoB+z3\nnbLXOTe0tD41/UNEROT8oekfIgFwzi0AFvxm3yPF3g8MeVAiIiJyXvCU3UREREREREqjolpERERE\nJEgqqkVEREREgqQ51SJByMrKOr579+7vKjqO88WhQ4dq16lT53BFx3E+UC7Ll/JZvpTP8qV8lp+c\nnJyLzvZcrf4hEgQzW+ecO6u7hOV0ymf5US7Ll/JZvpTP8qV8lp9gcqnpHyIiIiIiQVJRLSIiIiIS\nJBXVIsF5taIDOM8on+VHuSxfymf5Uj7Ll/JZfs46l5pTLSIiIiISJI1Ui4iIiIgESUW1iIiIiEiQ\nVFSLnAEzq2Vmn5rZDt/Xmn7aNTGzxWa2zcy2mlmz0EYaHgLNp69tDTNLN7OXQhljuAgkl2bW0cy+\nMLNvzGyTmY2oiFgrMzMbbGbfmdlOM3ughOMxZjbLd3ytfrZLF0A+/+T7HbnJzJaaWdOKiDMclJXL\nYu1+b2bOzLTEXikCyaeZXef7/vzGzGaW1aeKapEz8wCw1DnXCljq2y7Jm8Azzrk2QDfgpxDFF24C\nzSfAX4HPQhJVeAoklyeBPzjnLgEGA8+bWVIIY6zUzCwCmAr8J3AxMMrMLv5Ns1uBTOdcS2AK8FRo\nowwfAeZzA9DFOdcemA08Hdoow0OAucTMqgPjgbWhjTC8BJJPM2sFPAj09P3O/GNZ/aqoFjkzVwEz\nfO9nAFf/toHvBzPSOfcpgHPuuHPuZOhCDCtl5hPAzDoD9YDFIYorHJWZS+fcdufcDt/7H/H+sVcn\nZBFWft2Anc653c65POBdvHktrnieZwMDzMxCGGM4KTOfzrnlxX4/rgEahTjGcBHI9yZ4Bx+eAnJC\nGVwYCiSf/xuY6pzLBHDOlTk4pqJa5MzUc87t970/gLfQ+60U4KiZfWhmG8zsGd9fxXK6MvNpZh7g\nWeDeUAYWhgL53jzFzLoB0cCucx1YGGkI7Cu2ne7bV2Ib51wBcAxIDkl04SeQfBZ3K7DwnEYUvsrM\npZmlAo2dc/NDGViYCuR7MwVIMbPVZrbGzAaX1WlkOQYocl4wsyVA/RIOPVR8wznnzKykNSkjB/gf\nBwAAAytJREFUgd5AJ2AvMAsYA7xWvpGGh3LI513AAudcelUfECyHXP7STwPgLeAm51xR+UYpcubM\n7EagC9C3omMJR77Bh+fw/l8j5SMSaAX0w/sJymdm1s45d7S0E0SkGOfcQH/HzOygmTVwzu33FSYl\nfRyUDmx0zu32nTMHuJQqWlSXQz57AL3N7C4gAYg2s+POudLmX5+XyiGXmFkNYD7wkHNuzTkKNVz9\nADQutt3It6+kNulmFgkkAhmhCS/sBJJPzGwg3j8M+zrnckMUW7gpK5fVgbbACt/gQ33gYzMb6pxb\nF7Iow0cg35vpwFrnXD6wx8y24y2y0/x1qukfImfmY+Am3/ubgLkltEkDkszsl7mq/YGtIYgtHJWZ\nT+fcDc65Js65ZningLxZFQvqAJSZSzOLBj7Cm8PZIYwtXKQBrcysuS9XI/HmtbjieR4GLHN6ipo/\nZebTzDoBrwBDA5mzWoWVmkvn3DHnXG3nXDPf78o1eHOqgrpkgfysz8E7So2Z1cY7HWR3aZ2qqBY5\nM5OAQWa2Axjo28bMupjZPwCcc4V4i7+lZrYZMGBaBcVb2ZWZTwlYILm8DugDjDGzjb5Xx4oJt/Lx\nzZEeBywCtgHvOee+MbPHzWyor9lrQLKZ7QT+ROkr1lRpAebzGbyfQL3v+378bWEjBJxLCVCA+VwE\nZJjZVmA5cJ9zrtRPpfSYchERERGRIGmkWkREREQkSCqqRURERESCpKJaRERERCRIKqpFRERERIKk\nolpEREREJEgqqkVEREREgqSiWkREREQkSCqqRURE/DCzrma2ycxizayamX1jZm0rOi4RqXz08BcR\nEZFSmNkTQCwQB6Q7556s4JBEpBJSUS0iIlIKM4sG0oAc4DLnXGEFhyQilZCmf4iIiJQuGUgAquMd\nsRYROY1GqkVEREphZh8D7wLNgQbOuXEVHJKIVEKRFR2AiIhIZWVmfwDynXMzzSwC+NzM+jvnllV0\nbCJSuWikWkREREQkSJpTLSIiIiISJBXVIiIiIiJBUlEtIiIiIhIkFdUiIiIiIkFSUS0iIiIiEiQV\n1SIiIiIiQVJRLSIiIiISpP8PJV+rOd/F24QAAAAASUVORK5CYII=\n",
"text/plain": [
"<matplotlib.figure.Figure at 0x115a59f50>"
]
},
"metadata": {},
"output_type": "display_data"
}
],
"source": [
"q_l = 0.7\n",
"q_r = 0.2\n",
"v_l = 1.0\n",
"v_r = 1.0\n",
"t = 0.5\n",
"\n",
"states, speeds, reval = riemann_traffic_vc_exact(q_l,q_r,v_l,v_r)\n",
"riemann_tools.plot_riemann(states, speeds, reval, t,layout='vertical');\n",
"x1, q1 = test_solver(q_l,q_r,v_l,v_r,traffic_fwave_old)\n",
"plt.plot(x1,q1,'-x')\n",
"x2, q2 = test_solver(q_l,q_r,v_l,v_r,traffic_fwave)\n",
"plt.plot(x2,q2,'-o')\n",
"x3, q3 = test_solver(q_l,q_r,v_l,v_r,traffic_qwave)\n",
"plt.plot(x3,q3,'-s')\n",
"plt.legend(['Exact','Old','New','q-wave'])"
]
},
{
"cell_type": "code",
"execution_count": null,
"metadata": {
"collapsed": true
},
"outputs": [],
"source": []
}
],
"metadata": {
"kernelspec": {
"display_name": "Python 2",
"language": "python",
"name": "python2"
},
"language_info": {
"codemirror_mode": {
"name": "ipython",
"version": 2
},
"file_extension": ".py",
"mimetype": "text/x-python",
"name": "python",
"nbconvert_exporter": "python",
"pygments_lexer": "ipython2",
"version": "2.7.13"
}
},
"nbformat": 4,
"nbformat_minor": 2
}
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment