Comparing discretization techniques

In Periodic problems and plane-wave discretizations we saw how simple 1D problems can be modelled using plane-wave basis sets. This example invites you to work out some details on these aspects yourself using a number of exercises. The solutions are given at the bottom of the page.

For this example we consider the discretization of

\[ H = - \frac12 Δ + V(x) \quad \text{with $V(x) = \cos(x)$}\]

on $[0, 2π]$ with periodic boundary conditions. The $\cos(x)$ takes the role of a lattice-periodic potential. We will be interested in the smallest eigenvalues of this discretized Hamiltonian. Of note the boundary condition matters: The spectrum we will get is different from e.g. the spectrum of $H$ on $\mathbb{R}$.

Finite differences

We approximate functions $ψ$ on $[0, 2\pi]$ by their values at grid points $x_k = 2\pi \frac{k}{N}$, $k=1, \dots, N$. The boundary conditions are imposed by $ψ(x_0) = ψ(x_N), ψ(x_{N+1}) = ψ(x_1)$. We then have

\[ \big(Hψ\big)(x_k) \approx \frac 1 2 \frac{-ψ_{k-1} + 2 ψ_k - ψ_{k+1}}{δx^2} + V(x_k) ψ(x_k)\]

with $δx = \frac{2π}{N}$.

This can be put in matrix form in the following way:

# Finite differences Hamiltonian -1/2 Delta + V on [0, 2pi] with periodic bc.
# Pass it a function V.
using LinearAlgebra

function build_finite_differences_matrix(Vfunction, N::Integer)
    δx = 2π/N

    # Finite-difference approximation to -½Δ
    T = 1/(2δx^2) * Tridiagonal(-ones(N-1), 2ones(N), -ones(N-1))
    # The type Tridiagonal is efficient, but to establish the periodic boundary conditions
    # we need to add elements not on the three diagonals, so convert to dense matrix
    T = Matrix(T)
    T[1, N] = T[N, 1] = -1 / (2δx^2)

    # Finite-difference approximation to potential: We collect all coordinates ...
    x_coords = [k * δx for k=1:N]
    V = Diagonal(Vfunction.(x_coords))  # ... and evaluate V on each of the x_coords

    T + V
end;
Exercise 1

Show that the finite-difference approximation of -½Δ is indeed an approximation of the second derivative. Obtain an estimate of the first eigenvalue of $H$. Hint: Take a look at the eigen function from LinearAlgebra.

Plane waves method

In this method, we expand states on the basis

\[ e_G(x) = \frac{1}{\sqrt{2\pi}} e^{iGx} \qquad \text{for $G=-N,\dots,N$}.\]

Exercise 2

Show that

\[ \langle e_G, e_{G'}\rangle = ∫_0^{2π} e_G^\ast(x) e_{G'}(x) d x = δ_{G, G'}\]

and (assuming $V(x) = \cos(x)$)

\[ \langle e_G, H e_{G'}\rangle = \frac 1 2 \left(|G|^2 \delta_{G,G'} + \delta_{G, G'+1} + \delta_{G, G'-1}\right).\]

What happens for a more general $V(x)$?

Exercise 3

Code this and check the first eigenvalue agrees with the finite-difference case. Compare accuracies at various basis set sizes $N$.

Using DFTK

We now use DFTK to do the same plane-wave discretization in this 1D system. To deal with a 1D case we use a 3D lattice with two lattice vectors set to zero.

using DFTK
a = 2π
lattice = a .* [[1 0 0.]; [0 0 0]; [0 0 0]];

Define Hamiltonian: Kinetic + Potential

terms = [Kinetic(),
         ExternalFromReal(r -> cos(r[1]))]  # r is a vector of size 3
model = Model(lattice; n_electrons=1, terms, spin_polarization=:spinless);  # One spinless electron

Ecut defines the number of plane waves by selecting all those $G$, which satisfy the relationship $½ |G|^2 ≤ \text{Ecut}$.

Ecut = 500
basis = PlaneWaveBasis(model; Ecut, kgrid=(1, 1, 1))
PlaneWaveBasis discretization:
    architecture         : DFTK.CPU()
    num. mpi processes   : 1
    num. julia threads   : 1
    num. DFTK  threads   : 1
    num. blas  threads   : 2
    num. fft   threads   : 1

    Ecut                 : 500.0 Ha
    fft_size             : (135, 1, 1), 135 total points
    kgrid                : MonkhorstPack([1, 1, 1])
    num.   red. kpoints  : 1
    num. irred. kpoints  : 1

    Estimated memory usage (per MPI process):
        single ψ             :   3.94KiB
        single ρ             :   1.05KiB
        peak memory (SCF)    :   44.2KiB
    
    Discretized Model(custom, 1D):
        lattice (in Bohr)    : [6.28319   , 0         , 0         ]
                               [0         , 0         , 0         ]
                               [0         , 0         , 0         ]
        unit cell volume     : 6.2832 Bohr
    
        num. electrons       : 1
        spin polarization    : spinless
        temperature          : 0 Ha
    
        terms                : Kinetic()
                               ExternalFromReal(Main.var"##354".var"#3#4"())

We now seek the ground state using the self-consistent field algorithm.

scfres = self_consistent_field(basis; tol=1e-4)
scfres.energies
Energy breakdown (in Ha):
    Kinetic             0.2090845 
    ExternalFromReal    -0.7441493

    total               -0.535064852288

On this simple linear (non-interacting) model, the SCF converges in one step. The ground state energy of is simply the lowest eigenvalue; it should match the smallest eigenvalue of $H$ computed above.

Plotting

We can also get the first eigenvector (in the plane wave basis) and plot it

using Plots

ψ_fourier = scfres.ψ[1][:, 1];    # first k-point, all G components, first eigenvector
# Transform the wave function to real space
ψ = ifft(basis, basis.kpoints[1], ψ_fourier)[:, 1, 1]
# Eigenvectors are only defined up to a phase. We fix it by imposing that psi(0) is real
ψ /= (ψ[1] / abs(ψ[1]))
plot(real(ψ); label="")

Again this should match with the result above.

Exercise 4

Look at the Fourier coefficients of ψ_fourier and compare with the result above.

The DFTK Hamiltonian

We can ask DFTK for the Hamiltonian

E, ham = energy_hamiltonian(basis, scfres.ψ, scfres.occupation; ρ=scfres.ρ)
H = ham.blocks[1]
typeof(H)
DFTK.DftHamiltonianBlock

This is an opaque data structure, which encodes the Hamiltonian. What can we do with it?

using InteractiveUtils
methodswith(typeof(H), supertypes=true)
9-element Vector{Method}:

This defines a number of methods. For instance, it can be used as a linear operator:

H * DFTK.random_orbitals(basis, basis.kpoints[1], 1)
63×1 Matrix{ComplexF64}:
  -0.01740482502094278 + 0.03861800224251265im
  -0.00946501282491465 + 0.08599349098266994im
  0.060971394012929295 + 0.06333042203064607im
   0.11934623600772912 - 0.5842911405992695im
  -0.31728373711753743 + 0.17215501068702063im
    1.0301688595488998 - 0.013499088932366149im
    -0.520487349532157 + 2.1295342413569527im
   -2.1803402643295686 + 5.009550848621473im
     5.294493447290119 + 3.856232451529862im
    3.1215352745217633 - 4.3713420326353605im
    -4.320893102210905 + 12.18928632206811im
     6.268406432138969 + 5.546716668148607im
    11.755618367668012 - 1.6108768811703689im
    1.1443533447085579 + 13.514376380283604im
    -0.383877970048948 + 6.001686819379167im
    0.6153984430806799 - 1.305092054069409im
     8.280369877595776 + 3.6342685305185998im
     3.361256270198038 + 16.655977879933562im
   -15.719683308690843 - 0.5116271064391241im
     18.32641612632942 - 6.417832068293668im
   -15.640521566599995 - 16.17647783933925im
    12.133476695963402 - 7.446354583096801im
   -11.815117130227668 - 26.68314471488191im
   -22.015127443398935 + 43.63265558540008im
    -6.141391772484941 + 15.485019942760553im
    -28.49803803568792 + 21.590340503658982im
    -20.72733327290591 + 31.535188756173557im
      -6.4720338968107 + 26.773510109074138im
    -59.21535890853072 + 14.992263956728776im
   -2.7014675695490187 + 43.49922120468371im
   -7.7631767094822495 + 41.578564138701665im
    -45.16029105815769 + 36.28492539074999im
     23.85211367283726 + 28.57714930716191im
    11.330769455606692 - 9.260428218498616im
    -18.80614081665663 + 43.41001552741439im
   -30.269594843386948 + 13.075233776151821im
     55.14394161603507 - 5.218440441054651im
     -9.31566843938927 - 37.47771934136608im
    15.701146975561294 + 16.425832231014404im
     21.55074538393716 - 33.65358031203644im
    -8.766359682432505 + 20.215928607555515im
    3.7047883484496293 - 17.544413225253983im
     25.12556774547627 - 7.587146762393372im
   -11.403753700875757 + 4.347952310271904im
   -24.340282266294786 + 5.82026770860512im
    -2.177306878231722 + 1.2627963399645887im
     1.314330100644258 + 18.051403688797595im
    -3.204930706127625 + 17.781158441651492im
    1.9092814009581354 + 21.277092943729574im
    -21.68045986771249 - 2.314719238070994im
   -1.1785683782134184 + 2.6055488333192507im
    -6.067616391877231 + 18.227646612553716im
     5.556787645116421 + 3.625866971035207im
      3.12708024638485 + 2.785105520163669im
     -5.04971433603434 - 5.460726806727035im
    -2.690873038298538 + 1.2199129464889866im
   -3.1470089273117883 - 4.191929251676254im
    1.0703703318921287 + 0.21292880628880415im
    0.4113490797538778 + 0.6033560390558571im
    0.4163695025878387 - 1.4894032162573712im
  -0.09114985627242775 - 0.008505084536585245im
   0.08939981969700417 + 0.12001167270202742im
 -0.009068704124066965 + 0.05273966639073368im

We can also get its full matrix representation:

Array(H)
63×63 Matrix{ComplexF64}:
  9.38036e-19+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im
          0.5+2.37729e-18im           0.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im
  6.84621e-18+4.62374e-18im           0.5+2.37729e-18im           2.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im   -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im
 -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im           4.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im
  9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im           8.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im
  5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          12.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im
  1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          18.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im
 -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          24.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im
  6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          32.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im
  5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          40.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im  -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im
 -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          50.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im
  1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          60.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im
 -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          72.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im
 -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          84.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im
   6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          98.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im
  1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         112.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im
  3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         128.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im
  4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         144.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im
  6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         162.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im
 -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         180.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im
 -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         200.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im
 -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         220.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im
 -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         242.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im     2.3641e-17+8.69902e-18im   7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im
  1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         264.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im   -1.24221e-17+9.18236e-18im    2.3641e-17+8.69902e-18im   7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im
  -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         288.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im  -5.05513e-18+9.40828e-18im  -1.24221e-17+9.18236e-18im    2.3641e-17+8.69902e-18im   7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im
 -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         312.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.80457e-18-1.65903e-18im  -5.05513e-18+9.40828e-18im  -1.24221e-17+9.18236e-18im    2.3641e-17+8.69902e-18im   7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im
  1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         338.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   2.94354e-19-1.42818e-17im  -6.80457e-18-1.65903e-18im  -5.05513e-18+9.40828e-18im  -1.24221e-17+9.18236e-18im    2.3641e-17+8.69902e-18im   7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im
  6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         364.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   1.25324e-17+7.3538e-18im    2.94354e-19-1.42818e-17im  -6.80457e-18-1.65903e-18im  -5.05513e-18+9.40828e-18im  -1.24221e-17+9.18236e-18im    2.3641e-17+8.69902e-18im   7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im
 -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         392.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   2.40812e-18-1.91415e-17im   1.25324e-17+7.3538e-18im    2.94354e-19-1.42818e-17im  -6.80457e-18-1.65903e-18im  -5.05513e-18+9.40828e-18im  -1.24221e-17+9.18236e-18im    2.3641e-17+8.69902e-18im   7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im
          0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         420.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im   9.85983e-18+4.45332e-18im   2.40812e-18-1.91415e-17im   1.25324e-17+7.3538e-18im    2.94354e-19-1.42818e-17im  -6.80457e-18-1.65903e-18im  -5.05513e-18+9.40828e-18im  -1.24221e-17+9.18236e-18im    2.3641e-17+8.69902e-18im   7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im
  1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         450.0+0.0im                   0.5-1.91936e-18im    4.3417e-18+1.90195e-17im   9.85983e-18+4.45332e-18im   2.40812e-18-1.91415e-17im   1.25324e-17+7.3538e-18im    2.94354e-19-1.42818e-17im  -6.80457e-18-1.65903e-18im  -5.05513e-18+9.40828e-18im  -1.24221e-17+9.18236e-18im    2.3641e-17+8.69902e-18im   7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im  -2.32476e-18-1.03486e-17im
 -2.32476e-18-1.03486e-17im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         480.5+0.0im          -6.21078e-18-8.87719e-18im    4.3417e-18+1.90195e-17im   9.85983e-18+4.45332e-18im   2.40812e-18-1.91415e-17im   1.25324e-17+7.3538e-18im    2.94354e-19-1.42818e-17im  -6.80457e-18-1.65903e-18im  -5.05513e-18+9.40828e-18im  -1.24221e-17+9.18236e-18im    2.3641e-17+8.69902e-18im   7.46908e-18-4.89191e-18im  -3.96529e-18-1.63017e-17im   -6.5294e-18-5.68244e-18im  -2.09005e-17+1.30029e-17im   1.30228e-18-1.05794e-17im   7.50769e-18+1.8538e-17im    7.79391e-18-1.28996e-18im  -4.55526e-18-1.45669e-17im  -2.07003e-17+3.24918e-17im   2.76834e-17+2.04324e-17im   1.26093e-17-3.72139e-17im   -1.9652e-17-2.7462e-17im   -5.50336e-18+2.58911e-17im   1.50263e-17+4.32004e-18im   2.08933e-18-9.22629e-18im   -2.1533e-17+1.12123e-17im  -1.02584e-17-3.88425e-20im   5.20786e-18+2.01367e-18im   1.19815e-17-3.24554e-18im  -4.04378e-18+5.31047e-18im  -5.18221e-18+5.03797e-20im
 -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im    2.3641e-17-8.69902e-18im  -1.24221e-17-9.18236e-18im  -5.05513e-18-9.40828e-18im  -6.80457e-18+1.65903e-18im   2.94354e-19+1.42818e-17im   1.25324e-17-7.3538e-18im   -4.44267e-19-4.10411e-18im   9.85983e-18-4.45332e-18im   8.49259e-18-1.66418e-17im  -6.21078e-18+8.87719e-18im         480.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im
  1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im    2.3641e-17-8.69902e-18im  -1.24221e-17-9.18236e-18im  -5.05513e-18-9.40828e-18im  -6.80457e-18+1.65903e-18im   2.94354e-19+1.42818e-17im   1.25324e-17-7.3538e-18im   -4.44267e-19-4.10411e-18im   9.85983e-18-4.45332e-18im   8.49259e-18-1.66418e-17im           0.5+2.37729e-18im         450.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im
 -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im    2.3641e-17-8.69902e-18im  -1.24221e-17-9.18236e-18im  -5.05513e-18-9.40828e-18im  -6.80457e-18+1.65903e-18im   2.94354e-19+1.42818e-17im   1.25324e-17-7.3538e-18im   -4.44267e-19-4.10411e-18im   9.85983e-18-4.45332e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         420.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im
 -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im    2.3641e-17-8.69902e-18im  -1.24221e-17-9.18236e-18im  -5.05513e-18-9.40828e-18im  -6.80457e-18+1.65903e-18im   2.94354e-19+1.42818e-17im   1.25324e-17-7.3538e-18im   -4.44267e-19-4.10411e-18im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         392.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im
  6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im    2.3641e-17-8.69902e-18im  -1.24221e-17-9.18236e-18im  -5.05513e-18-9.40828e-18im  -6.80457e-18+1.65903e-18im   2.94354e-19+1.42818e-17im   1.25324e-17-7.3538e-18im    9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         364.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im
  1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im    2.3641e-17-8.69902e-18im  -1.24221e-17-9.18236e-18im  -5.05513e-18-9.40828e-18im  -6.80457e-18+1.65903e-18im   2.94354e-19+1.42818e-17im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         338.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im
 -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im    2.3641e-17-8.69902e-18im  -1.24221e-17-9.18236e-18im  -5.05513e-18-9.40828e-18im  -6.80457e-18+1.65903e-18im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         312.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im
  -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im    2.3641e-17-8.69902e-18im  -1.24221e-17-9.18236e-18im  -5.05513e-18-9.40828e-18im  -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         288.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im
  1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im    2.3641e-17-8.69902e-18im  -1.24221e-17-9.18236e-18im   6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         264.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im
 -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im    2.3641e-17-8.69902e-18im   5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         242.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im
 -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   7.46908e-18+4.89191e-18im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         220.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im
 -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.96529e-18+1.63017e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         200.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im
 -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im   -6.5294e-18+5.68244e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         180.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im
  6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im  -2.09005e-17-1.30029e-17im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         162.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im
  4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.30228e-18+1.05794e-17im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         144.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im
  1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    7.50769e-18-1.8538e-17im    1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         128.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im
  1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   3.42423e-18+1.3351e-17im    3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im         112.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im
  6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im  -4.55526e-18+1.45669e-17im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          98.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im
 -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -7.79391e-18-1.28996e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          84.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im
 -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im   2.76834e-17-2.04324e-17im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          72.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im
  1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im    1.26093e-17+3.72139e-17im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          60.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im
 -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im   -1.9652e-17+2.7462e-17im   -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          50.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im
  5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im  -5.50336e-18-2.58911e-17im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          40.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im
  6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   1.50263e-17-4.32004e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          32.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im
 -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im   2.08933e-18+9.22629e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          24.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im
  1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   -2.1533e-17-1.12123e-17im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          18.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im
  5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im  -1.02584e-17+3.88425e-20im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im          12.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im
  9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im   5.20786e-18-2.01367e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im           8.0+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im
 -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im   1.19815e-17+3.24554e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im           4.5+0.0im                   0.5-1.91936e-18im   6.84621e-18-4.62374e-18im
  6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im  -4.04378e-18-5.31047e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im           2.0+0.0im                   0.5-1.91936e-18im
          0.5-1.91936e-18im   6.84621e-18-4.62374e-18im  -2.10986e-18-7.33721e-18im   9.53855e-19+2.05757e-17im   5.30139e-18+9.96275e-18im   1.15239e-17-1.06252e-18im  -6.09606e-18-2.01355e-17im   6.73916e-18-9.6349e-18im    5.08456e-18-5.8805e-19im   -5.68698e-18+2.3003e-19im    1.48741e-17+1.10137e-17im  -3.44396e-18-3.19427e-17im  -4.12338e-17-4.82849e-18im   6.96319e-18+3.06423e-17im   1.41552e-18-4.45966e-18im   1.06844e-17+1.07046e-17im   4.92777e-18+5.57465e-18im   6.87119e-18-1.51362e-19im  -6.53743e-18-2.77989e-18im  -2.28337e-18-2.65376e-18im  -4.23963e-18+3.77063e-18im  -1.07755e-17-3.20506e-18im   1.64216e-17-3.33489e-18im   -4.8805e-18-2.16383e-18im  -2.56782e-17-4.87146e-18im   1.35838e-17+5.69649e-18im   6.55612e-18+3.83412e-18im  -7.57355e-18-5.37897e-18im  -8.89745e-18+9.6317e-18im    1.13146e-17+6.44402e-18im  -8.98722e-18-7.74616e-18im  -5.18221e-18-5.03797e-20im   1.13146e-17-6.44402e-18im           0.0+6.78716e-18im  -7.57355e-18+5.37897e-18im   6.55612e-18-3.83412e-18im   1.35838e-17-5.69649e-18im  -2.56782e-17+4.87146e-18im   -4.8805e-18+2.16383e-18im   1.64216e-17+3.33489e-18im  -1.07755e-17+3.20506e-18im  -4.23963e-18-3.77063e-18im  -2.28337e-18+2.65376e-18im  -6.53743e-18+2.77989e-18im   6.87119e-18+1.51362e-19im   4.92777e-18-5.57465e-18im   3.61585e-19-1.24616e-17im   1.41552e-18+4.45966e-18im    6.5791e-18-3.57518e-18im  -4.12338e-17+4.82849e-18im  -3.44396e-18+3.19427e-17im   1.48741e-17-1.10137e-17im  -5.68698e-18-2.3003e-19im    5.08456e-18+5.8805e-19im    6.73916e-18+9.6349e-18im   -6.09606e-18+2.01355e-17im   1.15239e-17+1.06252e-18im   5.30139e-18-9.96275e-18im   9.53855e-19-2.05757e-17im  -2.10986e-18+7.33721e-18im   6.84621e-18+4.62374e-18im           0.5+2.37729e-18im           0.5+0.0im
Exercise 5

Compare this matrix Array(H) with the one you obtained in Exercise 3, get its eigenvectors and eigenvalues. Try to guess the ordering of $G$-vectors in DFTK.

Exercise 6

Increase the size of the problem, and compare the time spent by DFTK's internal diagonalization algorithms to a full diagonalization of Array(H). Hint: The @belapsed and @benchmark macros (from the BenchmarkTools package) are handy for this task. Note that there are some subtleties with global variables (see the BenchmarkTools docs for details). E.g. to use it to benchmark a function like eigen(H) run it as (note the $):

using BenchmarkTools
@benchmark eigen($H)

Solutions

Exercise 1

If we consider a function $f : [0, 2π] → \mathbb{R}$, to first order

\[f(x + δx) = f(x) + δx f'(x) + O(δx^2)\]

therefore after rearrangement

\[f'(x) = \frac{f(x + δx) - f(x)}{δx} + O(δx).\]

Similarly

\[f''(x) = \frac{f'(x + δx) - f'(x)}{δx} + O(δx),\]

such that overall

\[f''(x) \simeq \frac{f(x + 2δx) - f(x + δx) - f(x + δx) + f(x)}{δx^2} = \frac{f(x + 2δx) - 2f(x + δx) + f(x)}{δx^2}\]

In finite differences we consider a stick basis of vectors

\[\left\{ e_i = (0, …, 0, \underbrace{δx}_\text{$i$-th position}, 0, …, 0) \middle| i = 1, … N \right\}.\]

Keeping in mind the periodic boundary conditions (i.e. $e_0 = e_N$) projecting the Hamiltonian $H$ onto this basis thus yields the proposed structure.

We start off with $N = 100$ to obtain

Hfd = build_finite_differences_matrix(cos, 100)
L, V = eigen(Hfd)
L[1:5]
5-element Vector{Float64}:
 -0.5351640695081414
  0.3429576842587494
  0.8532206682511969
  2.053679933442362
  2.078512110605267

This is already pretty accurate (to about 4 digits) as can be estimated looking at the following convergence plot:

function fconv(N)
    L, V = eigen(build_finite_differences_matrix(cos, N))
    first(L)
end
Nrange = 10:10:100
plot(Nrange, abs.(fconv.(Nrange) .- fconv(200));
     yaxis=:log, legend=false, mark=:x, xlabel="N", ylabel="Absolute error")

Exercise 2

  • We note that

    \[\langle e_G, e_{G'}\rangle = ∫_0^{2π} e_G^\ast(x) e_{G'}(x) d x = 1/2π ∫_0^{2π} e^{i(G'-G)x} d x\]

    Since $e^{iy}$ is a periodic function with period $2\pi$, $\int_0^{2\pi} e^{i m y} = \delta_{0,m}$. Therefore if $G≠G'$ we have that $\langle e_G, e_{G'}\rangle = 0$, while $G=G'$ implies $\langle e_G, e_{G'}\rangle = 1$. In summary:

    \[\langle e_G, e_{G'}\rangle = δ_{G, G'}\]

  • Next fo $V(x) = \cos(x)$ we obtain

    \[\langle e_G, H e_{G'}\rangle = \frac 1 2 ∫_0^{2π} e_G^\ast(x) H e_{G'}(x) d x\]

    We start by applying the Hamiltonian to a plane-wave:

    \[H e_{G'}(x) = - \frac 1 2 (-|G|^2) \frac 1 {\sqrt{2π}} e^{iG'x) + cos(x) \frac 1 {\sqrt{2π}} e{iG'x}\]

    Then, using the result of the first part of the exercise and the fact that $cos(x) = \frac 1 2 \left(e{ix} + e{-ix}\right)$, we get:

    \[\begin{align*} ⟨ e_G, H e_{G'}⟩ &= \frac 1 2 G^2 δ_{G, G'} + \frac 1 {4π} \left(∫_0^{2π} e^{ix ⋅ (G'-G+1)} d x + ∫_0^{2π} e^{ix ⋅ (G'-G-1)} d x \right) \\ &= \frac 1 2 \left(|G|^2 \delta_{G,G'} + \delta_{G, G'+1} + \delta_{G, G'-1}\right) \end{align*}\]

  • In case a more general $V(x)$ was employed, this potential still has to be periodic over $[0, 2\pi]$ to fit our setting. Assuming sufficient regularity in $V$ we can employ a Fourier series:

    \[V(x) = \sum_{G=- \infty}^{\infty} \hat{V}_G e_G(x)\]

    where

    \[\hat{V}_G = \frac{1}{\sqrt{2π}} ∫_0^{2π} V(x) e^{-iGx} dx = ∫_0^{2π} V(x) e_G^\ast dx .\]

    Note that one can change of this as a change of basis from the position basis to the plane-wave basis.

    Based on this expansion

    \[\begin{align*} ⟨ e_G, V e_{G'} ⟩ &= \left\langle e_G, ∑_{G''} \hat{V}_{G''} \, e_{G'+G''} \right\rangle \\ &= \sum_{G''=-\infty}^\infty \hat{V}_{G''} ⟨ e_G, e_{G'+G''} ⟩ \\ &= \sum_{G''=-\infty}^\infty \hat{V}_{G''} \, δ_{G-G', G''} ⟩ \\ &= \hat{V}_{G-G'} \end{align*}\]

    and therefore

    \[⟨ e_G, H e_{G'} ⟩ = \frac 1 2 |G|^2 \delta_{G,G'} + \hat{V}_{G-G'},\]

    i.e. essentially the Fourier transform of $V$ determines the contribution to the matrix elements of the Hamiltonian.

Exercise 3

The Hamiltonian matrix for the plane waves method can be found this way:

# Plane waves Hamiltonian -½Δ + cos on [0, 2pi].
function build_plane_waves_matrix_cos(N::Integer)
    # Plane wave approximation to -½Δ
    Gsq = [float(i)^2 for i in -N:N]
    # Hamiltonian as derived in Exercise 2:
    1/2 * Tridiagonal(ones(2N), Gsq, ones(2N))
end;

Then we check that the first eigenvalue agrees with the finite-difference case, using $N = 10$:

Hpw_cos = build_plane_waves_matrix_cos(10)
L, V = eigen(Hpw_cos)
L[1:5]
5-element Vector{Float64}:
 -0.5350648522878154
  0.34336012839908336
  0.8536343543207997
  2.0565044112660984
  2.081227363352144

We look at the convergence plot to compare the accuracy for various numbers of plane-waves $N$:

function fconv(N)
    L = eigvals(build_plane_waves_matrix_cos(N))
    first(L)
end

Nrange = 2:10
plot(Nrange, abs.(fconv.(Nrange) .- fconv(200)); yaxis=:log, legend=false,
     ylims=(1e-16,Inf), ylabel="Absolute error", xlabel="N", mark=:x)

Notice how compared to exercise 1 the considered basis size $N$ is much smaller, indicating that plane-wave methods more quickly lead to accurate solutions than finite-difference methods.

Exercise 4

For efficiency reasons the data in Fourier space is not ordered increasingly with $G$. Therefore to plot the Fourier space representation sensibly, we need to sort by ascending values of the $G$ vectors first. For this we extract the Fourier vector of each plane-wave basis function in the index order:

coords_G_vectors = G_vectors_cart(basis, basis.kpoints[1])  # Get coordinates of first and only k-point

# Only keep first component of each vector (because the others are zero for 1D problems):
coords_Gx = [G[1] for G in coords_G_vectors]

p = plot(coords_Gx, real(ψ_fourier); label="real part", xlims=(-10, 10))
plot!(p, coords_Gx, imag(ψ_fourier); label="imaginary part")

The plot is symmetric about the zero (confirming that the orbitals are real) and only takes peaked values, which corresponds to the expected result for a cosine potential.

Exercise 5

To figure out the ordering we consider a small basis and build the Hamiltonian:

basis_small  = PlaneWaveBasis(model; Ecut=5, kgrid=(1, 1, 1))
ham_small = Hamiltonian(basis_small)
H_small = Array(ham_small.blocks[1])
H_small[abs.(H_small) .< 1e-12] .= 0  # Drop numerically zero entries
H_small
7×7 Matrix{ComplexF64}:
 0.0+0.0im          0.5+1.00424e-16im  0.0+0.0im          0.0+0.0im          0.0+0.0im          0.0+0.0im          0.5-1.00424e-16im
 0.5-1.00424e-16im  0.5+0.0im          0.5+1.00424e-16im  0.0+0.0im          0.0+0.0im          0.0+0.0im          0.0+0.0im
 0.0+0.0im          0.5-1.00424e-16im  2.0+0.0im          0.5+1.00424e-16im  0.0+0.0im          0.0+0.0im          0.0+0.0im
 0.0+0.0im          0.0+0.0im          0.5-1.00424e-16im  4.5+0.0im          0.0+0.0im          0.0+0.0im          0.0+0.0im
 0.0+0.0im          0.0+0.0im          0.0+0.0im          0.0+0.0im          4.5+0.0im          0.5+1.00424e-16im  0.0+0.0im
 0.0+0.0im          0.0+0.0im          0.0+0.0im          0.0+0.0im          0.5-1.00424e-16im  2.0+0.0im          0.5+1.00424e-16im
 0.5+1.00424e-16im  0.0+0.0im          0.0+0.0im          0.0+0.0im          0.0+0.0im          0.5-1.00424e-16im  0.5+0.0im

The equivalent version using the build_plane_waves_matrix_cos function is N=3 (both give rice to a 7×7 matrix).

Hother = build_plane_waves_matrix_cos(3)
7×7 LinearAlgebra.Tridiagonal{Float64, Vector{Float64}}:
 4.5  0.5   ⋅    ⋅    ⋅    ⋅    ⋅ 
 0.5  2.0  0.5   ⋅    ⋅    ⋅    ⋅ 
  ⋅   0.5  0.5  0.5   ⋅    ⋅    ⋅ 
  ⋅    ⋅   0.5  0.0  0.5   ⋅    ⋅ 
  ⋅    ⋅    ⋅   0.5  0.5  0.5   ⋅ 
  ⋅    ⋅    ⋅    ⋅   0.5  2.0  0.5
  ⋅    ⋅    ⋅    ⋅    ⋅   0.5  4.5

By comparing the entries we find the ordering is 0,1,2,...,-2,-1, which can also be found by inspecting

first.(G_vectors(basis_small, basis_small.kpoints[1]))
7-element Vector{Int64}:
  0
  1
  2
  3
 -3
 -2
 -1

Both matrices have the same eigenvalues:

maximum(abs, eigvals(H_small) - eigvals(Hother))
1.7763568394002505e-15

and in the eigenvectors we find the same rearrangements in the entries of the eigenvectors of both matrices, matching the DFTK ordering of is 0,1,2,...,-2,-1.

eigvecs(Hother)[:, 1]
7-element Vector{Float64}:
 -0.008461095314567344
  0.08520425570524764
 -0.42353481192309905
  0.7915641575577494
 -0.42353481192310566
  0.08520425570524905
 -0.008461095314567471
eigvecs(H_small)[:, 1]
7-element Vector{ComplexF64}:
  -0.7915641575577493 + 1.5898485476322146e-16im
  0.42353481192310877 - 1.7013307112974292e-16im
 -0.08520425570525011 + 5.1339563909648e-17im
 0.008461095314567624 - 6.797609502050578e-18im
 0.008461095314567069 + 3.3988047510250303e-18im
 -0.08520425570524359 - 1.711318796988132e-17im
   0.4235348119230965 + 0.0im

Notice, that eigenvectors are only defined up to a phase, so the sign may globally be inverted between the two eigenvectors.

Exercise 6

We benchmark the time needed for a full diagonalization (instantiation of the Array plus call of eigen) versus the time needed for running the SCF (i.e. iterative diagonalization using plane waves).

using Printf

for Ecut in 200:200:1600
   basis_time = PlaneWaveBasis(model; Ecut, kgrid=(1, 1, 1))
   t_eigen = @elapsed eigen(Array(Hamiltonian(basis_time).blocks[1]))
   t_scf   = @elapsed self_consistent_field(basis_time; tol=1e-6, callback=identity);
   @printf "%4i  eigen=%8.4fms  scf=%8.4fms\n" Ecut 1000t_eigen 1000t_scf
end
 200  eigen=  2.4049ms  scf=  4.9015ms
 400  eigen=  4.0740ms  scf=  2.7551ms
 600  eigen=  6.1433ms  scf=  2.9348ms
 800  eigen= 10.1058ms  scf=  3.1102ms
1000  eigen= 12.6920ms  scf=  3.0556ms
1200  eigen= 15.3636ms  scf=  3.2065ms
1400  eigen= 17.3888ms  scf=  3.0588ms
1600  eigen= 20.6633ms  scf=  3.3610ms