{{indexmenu_n>4}}
====== Temperature ======
###
The effect of temperature can be added by calculating excited states and expectation values of excited states. The temperature dependent expectation value is then created using Boltzmann statistics.
###
###
A small example:
-- Sofar we calculated eigenstates and expectation values (or spectra) of these
-- eigenstates. At 0 K one would measure the expectation value of the lowest eigenstate
-- at finite temperature one would measure an average over several states weighted by
-- Boltzmann statistics. In this example we calculate the temperature dependent
-- x-ray absorption spectra of NiO. (Ni L23 edge 2p to 3d) within the ligand-field
-- theory approximation
-- The first part is an exact copy of example 41
Verbosity(0)
-- here we calculate the 2p to 3d x-ray absorption of NiO within the Ligand-field theory
-- approximation. The first part of the script is very much the same as calculating
-- the ground-state with the addition that we now also need a 2p core shell in the basis
-- from the previous example we know that within NiO there are 3 states close to each other
-- and then there is an energy gap of about 1 eV. We thus only need to consider the 3
-- lowest states (Npsi=3 later on)
NF=26
NB=0
IndexDn_2p={ 0, 2, 4}
IndexUp_2p={ 1, 3, 5}
IndexDn_3d={ 6, 8,10,12,14}
IndexUp_3d={ 7, 9,11,13,15}
IndexDn_Ld={16,18,20,22,24}
IndexUp_Ld={17,19,21,23,25}
-- angular momentum operators on the d-shell
OppSx_3d =NewOperator("Sx" ,NF, IndexUp_3d, IndexDn_3d)
OppSy_3d =NewOperator("Sy" ,NF, IndexUp_3d, IndexDn_3d)
OppSz_3d =NewOperator("Sz" ,NF, IndexUp_3d, IndexDn_3d)
OppSsqr_3d =NewOperator("Ssqr" ,NF, IndexUp_3d, IndexDn_3d)
OppSplus_3d=NewOperator("Splus",NF, IndexUp_3d, IndexDn_3d)
OppSmin_3d =NewOperator("Smin" ,NF, IndexUp_3d, IndexDn_3d)
OppLx_3d =NewOperator("Lx" ,NF, IndexUp_3d, IndexDn_3d)
OppLy_3d =NewOperator("Ly" ,NF, IndexUp_3d, IndexDn_3d)
OppLz_3d =NewOperator("Lz" ,NF, IndexUp_3d, IndexDn_3d)
OppLsqr_3d =NewOperator("Lsqr" ,NF, IndexUp_3d, IndexDn_3d)
OppLplus_3d=NewOperator("Lplus",NF, IndexUp_3d, IndexDn_3d)
OppLmin_3d =NewOperator("Lmin" ,NF, IndexUp_3d, IndexDn_3d)
OppJx_3d =NewOperator("Jx" ,NF, IndexUp_3d, IndexDn_3d)
OppJy_3d =NewOperator("Jy" ,NF, IndexUp_3d, IndexDn_3d)
OppJz_3d =NewOperator("Jz" ,NF, IndexUp_3d, IndexDn_3d)
OppJsqr_3d =NewOperator("Jsqr" ,NF, IndexUp_3d, IndexDn_3d)
OppJplus_3d=NewOperator("Jplus",NF, IndexUp_3d, IndexDn_3d)
OppJmin_3d =NewOperator("Jmin" ,NF, IndexUp_3d, IndexDn_3d)
Oppldots_3d=NewOperator("ldots",NF, IndexUp_3d, IndexDn_3d)
-- Angular momentum operators on the Ligand shell
OppSx_Ld =NewOperator("Sx" ,NF, IndexUp_Ld, IndexDn_Ld)
OppSy_Ld =NewOperator("Sy" ,NF, IndexUp_Ld, IndexDn_Ld)
OppSz_Ld =NewOperator("Sz" ,NF, IndexUp_Ld, IndexDn_Ld)
OppSsqr_Ld =NewOperator("Ssqr" ,NF, IndexUp_Ld, IndexDn_Ld)
OppSplus_Ld=NewOperator("Splus",NF, IndexUp_Ld, IndexDn_Ld)
OppSmin_Ld =NewOperator("Smin" ,NF, IndexUp_Ld, IndexDn_Ld)
OppLx_Ld =NewOperator("Lx" ,NF, IndexUp_Ld, IndexDn_Ld)
OppLy_Ld =NewOperator("Ly" ,NF, IndexUp_Ld, IndexDn_Ld)
OppLz_Ld =NewOperator("Lz" ,NF, IndexUp_Ld, IndexDn_Ld)
OppLsqr_Ld =NewOperator("Lsqr" ,NF, IndexUp_Ld, IndexDn_Ld)
OppLplus_Ld=NewOperator("Lplus",NF, IndexUp_Ld, IndexDn_Ld)
OppLmin_Ld =NewOperator("Lmin" ,NF, IndexUp_Ld, IndexDn_Ld)
OppJx_Ld =NewOperator("Jx" ,NF, IndexUp_Ld, IndexDn_Ld)
OppJy_Ld =NewOperator("Jy" ,NF, IndexUp_Ld, IndexDn_Ld)
OppJz_Ld =NewOperator("Jz" ,NF, IndexUp_Ld, IndexDn_Ld)
OppJsqr_Ld =NewOperator("Jsqr" ,NF, IndexUp_Ld, IndexDn_Ld)
OppJplus_Ld=NewOperator("Jplus",NF, IndexUp_Ld, IndexDn_Ld)
OppJmin_Ld =NewOperator("Jmin" ,NF, IndexUp_Ld, IndexDn_Ld)
-- total angular momentum
OppSx = OppSx_3d + OppSx_Ld
OppSy = OppSy_3d + OppSy_Ld
OppSz = OppSz_3d + OppSz_Ld
OppSsqr = OppSx * OppSx + OppSy * OppSy + OppSz * OppSz
OppLx = OppLx_3d + OppLx_Ld
OppLy = OppLy_3d + OppLy_Ld
OppLz = OppLz_3d + OppLz_Ld
OppLsqr = OppLx * OppLx + OppLy * OppLy + OppLz * OppLz
OppJx = OppJx_3d + OppJx_Ld
OppJy = OppJy_3d + OppJy_Ld
OppJz = OppJz_3d + OppJz_Ld
OppJsqr = OppJx * OppJx + OppJy * OppJy + OppJz * OppJz
-- define the coulomb operator
-- we here define the part depending on F0 seperately from the part depending on F2
-- when summing we can put in the numerical values of the slater integrals
OppF0_3d =NewOperator("U", NF, IndexUp_3d, IndexDn_3d, {1,0,0})
OppF2_3d =NewOperator("U", NF, IndexUp_3d, IndexDn_3d, {0,1,0})
OppF4_3d =NewOperator("U", NF, IndexUp_3d, IndexDn_3d, {0,0,1})
-- define onsite energies - crystal field
-- Akm = {{k1,m1,Akm1},{k2,m2,Akm2}, ... }
Akm = PotentialExpandedOnClm("Oh", 2, {0.6,-0.4})
OpptenDq_3d = NewOperator("CF", NF, IndexUp_3d, IndexDn_3d, Akm)
OpptenDq_Ld = NewOperator("CF", NF, IndexUp_Ld, IndexDn_Ld, Akm)
Akm = PotentialExpandedOnClm("Oh", 2, {1,0})
OppNeg_3d = NewOperator("CF", NF, IndexUp_3d, IndexDn_3d, Akm)
OppNeg_Ld = NewOperator("CF", NF, IndexUp_Ld, IndexDn_Ld, Akm)
Akm = PotentialExpandedOnClm("Oh", 2, {0,1})
OppNt2g_3d = NewOperator("CF", NF, IndexUp_3d, IndexDn_3d, Akm)
OppNt2g_Ld = NewOperator("CF", NF, IndexUp_Ld, IndexDn_Ld, Akm)
OppNUp_2p = NewOperator("Number", NF, IndexUp_2p, IndexUp_2p, {1,1,1})
OppNDn_2p = NewOperator("Number", NF, IndexDn_2p, IndexDn_2p, {1,1,1})
OppN_2p = OppNUp_2p + OppNDn_2p
OppNUp_3d = NewOperator("Number", NF, IndexUp_3d, IndexUp_3d, {1,1,1,1,1})
OppNDn_3d = NewOperator("Number", NF, IndexDn_3d, IndexDn_3d, {1,1,1,1,1})
OppN_3d = OppNUp_3d + OppNDn_3d
OppNUp_Ld = NewOperator("Number", NF, IndexUp_Ld, IndexUp_Ld, {1,1,1,1,1})
OppNDn_Ld = NewOperator("Number", NF, IndexDn_Ld, IndexDn_Ld, {1,1,1,1,1})
OppN_Ld = OppNUp_Ld + OppNDn_Ld
-- define L-d interaction
Akm = PotentialExpandedOnClm("Oh", 2, {1,0})
OppVeg = NewOperator("CF", NF, IndexUp_3d, IndexDn_3d, IndexUp_Ld, IndexDn_Ld,Akm) + NewOperator("CF", NF, IndexUp_Ld, IndexDn_Ld, IndexUp_3d, IndexDn_3d, Akm)
Akm = PotentialExpandedOnClm("Oh", 2, {0,1})
OppVt2g = NewOperator("CF", NF, IndexUp_3d, IndexDn_3d, IndexUp_Ld, IndexDn_Ld,Akm) + NewOperator("CF", NF, IndexUp_Ld, IndexDn_Ld, IndexUp_3d, IndexDn_3d, Akm)
-- core valence interaction
Oppcldots= NewOperator("ldots", NF, IndexUp_2p, IndexDn_2p)
OppUpdF0 = NewOperator("U", NF, IndexUp_2p, IndexDn_2p, IndexUp_3d, IndexDn_3d, {1,0}, {0,0})
OppUpdF2 = NewOperator("U", NF, IndexUp_2p, IndexDn_2p, IndexUp_3d, IndexDn_3d, {0,1}, {0,0})
OppUpdG1 = NewOperator("U", NF, IndexUp_2p, IndexDn_2p, IndexUp_3d, IndexDn_3d, {0,0}, {1,0})
OppUpdG3 = NewOperator("U", NF, IndexUp_2p, IndexDn_2p, IndexUp_3d, IndexDn_3d, {0,0}, {0,1})
-- dipole transition
t=math.sqrt(1/2)
Akm = {{1,-1,t},{1, 1,-t}}
TXASx = NewOperator("CF", NF, IndexUp_3d, IndexDn_3d, IndexUp_2p, IndexDn_2p, Akm)
Akm = {{1,-1,t*I},{1, 1,t*I}}
TXASy = NewOperator("CF", NF, IndexUp_3d, IndexDn_3d, IndexUp_2p, IndexDn_2p, Akm)
Akm = {{1,0,1}}
TXASz = NewOperator("CF", NF, IndexUp_3d, IndexDn_3d, IndexUp_2p, IndexDn_2p, Akm)
TXASr = t*(TXASx - I * TXASy)
TXASl =-t*(TXASx + I * TXASy)
-- We follow the energy definitions as introduced in the group of G.A. Sawatzky (Groningen)
-- J. Zaanen, G.A. Sawatzky, and J.W. Allen PRL 55, 418 (1985)
-- for parameters of specific materials see
-- A.E. Bockquet et al. PRB 55, 1161 (1996)
-- After some initial discussion the energies U and Delta refer to the center of a configuration
-- The L^10 d^n configuration has an energy 0
-- The L^9 d^n+1 configuration has an energy Delta
-- The L^8 d^n+2 configuration has an energy 2*Delta+Udd
--
-- If we relate this to the onsite energy of the L and d orbitals we find
-- 10 eL + n ed + n(n-1) U/2 == 0
-- 9 eL + (n+1) ed + (n+1)n U/2 == Delta
-- 8 eL + (n+2) ed + (n+1)(n+2) U/2 == 2*Delta+U
-- 3 equations with 2 unknowns, but with interdependence yield:
-- ed = (10*Delta-nd*(19+nd)*U/2)/(10+nd)
-- eL = nd*((1+nd)*Udd/2-Delta)/(10+nd)
--
-- For the final state we/they defined
-- The 2p^5 L^10 d^n+1 configuration has an energy 0
-- The 2p^5 L^9 d^n+2 configuration has an energy Delta + Udd - Upd
-- The 2p^5 L^8 d^n+3 configuration has an energy 2*Delta + 3*Udd - 2*Upd
--
-- If we relate this to the onsite energy of the p and d orbitals we find
-- 6 ep + 10 eL + n ed + n(n-1) Udd/2 + 6 n Upd == 0
-- 6 ep + 9 eL + (n+1) ed + (n+1)n Udd/2 + 6 (n+1) Upd == Delta
-- 6 ep + 8 eL + (n+2) ed + (n+1)(n+2) Udd/2 + 6 (n+2) Upd == 2*Delta+Udd
-- 5 ep + 10 eL + (n+1) ed + (n+1)(n) Udd/2 + 5 (n+1) Upd == 0
-- 5 ep + 9 eL + (n+2) ed + (n+2)(n+1) Udd/2 + 5 (n+2) Upd == Delta+Udd-Upd
-- 5 ep + 8 eL + (n+3) ed + (n+3)(n+2) Udd/2 + 5 (n+3) Upd == 2*Delta+3*Udd-2*Upd
-- 6 equations with 3 unknowns, but with interdependence yield:
-- epfinal = (10*Delta + (1+nd)*(nd*Udd/2-(10+nd)*Upd) / (16+nd)
-- edfinal = (10*Delta - nd*(31+nd)*Udd/2-90*Upd) / (16+nd)
-- eLfinal = ((1+nd)*(nd*Udd/2+6*Upd)-(6+nd)*Delta) / (16+nd)
--
--
--
-- note that ed-ep = Delta - nd * U and not Delta
-- note furthermore that ep and ed here are defined for the onsite energy if the system had
-- locally nd electrons in the d-shell. In DFT or Hartree Fock the d occupation is in the end not
-- nd and thus the onsite energy of the Kohn-Sham orbitals is not equal to ep and ed in model
-- calculations.
--
-- note furthermore that ep and eL actually should be different for most systems. We happily ignore this fact
--
-- We normally take U and Delta as experimentally determined parameters
-- number of electrons (formal valence)
nd = 8
-- parameters from experiment (core level PES)
Udd = 7.3
Upd = 8.5
Delta = 4.7
-- parameters obtained from DFT (PRB 85, 165113 (2012))
F2dd = 11.14
F4dd = 6.87
F2pd = 6.67
G1pd = 4.92
G3pd = 2.80
tenDq = 0.56
tenDqL = 1.44
Veg = 2.06
Vt2g = 1.21
zeta_3d = 0.081
zeta_2p = 11.51
Bz = 0.000001
Hz = 0.120
ed = (10*Delta-nd*(19+nd)*Udd/2)/(10+nd)
eL = nd*((1+nd)*Udd/2-Delta)/(10+nd)
epfinal = (10*Delta + (1+nd)*(nd*Udd/2-(10+nd)*Upd)) / (16+nd)
edfinal = (10*Delta - nd*(31+nd)*Udd/2-90*Upd) / (16+nd)
eLfinal = ((1+nd)*(nd*Udd/2+6*Upd) - (6+nd)*Delta) / (16+nd)
F0dd = Udd + (F2dd+F4dd) * 2/63
F0pd = Upd + (1/15)*G1pd + (3/70)*G3pd
Hamiltonian = F0dd*OppF0_3d + F2dd*OppF2_3d + F4dd*OppF4_3d + zeta_3d*Oppldots_3d + Bz*(2*OppSz_3d + OppLz_3d) + Hz * OppSz_3d + tenDq*OpptenDq_3d + tenDqL*OpptenDq_Ld + Veg * OppVeg + Vt2g * OppVt2g + ed * OppN_3d + eL * OppN_Ld
XASHamiltonian = F0dd*OppF0_3d + F2dd*OppF2_3d + F4dd*OppF4_3d + zeta_3d*Oppldots_3d + Bz*(2*OppSz_3d + OppLz_3d)+ Hz * OppSz_3d + tenDq*OpptenDq_3d + tenDqL*OpptenDq_Ld + Veg * OppVeg + Vt2g * OppVt2g + edfinal * OppN_3d + eLfinal * OppN_Ld + epfinal * OppN_2p + zeta_2p * Oppcldots + F0pd * OppUpdF0 + F2pd * OppUpdF2 + G1pd * OppUpdG1 + G3pd * OppUpdG3
-- we now can create the lowest Npsi eigenstates:
Npsi=3
-- in order to make sure we have a filling of 8 electrons we need to define some restrictions
StartRestrictions = {NF, NB, {"000000 1111111111 0000000000",8,8}, {"111111 0000000000 1111111111",16,16}}
psiList = Eigensystem(Hamiltonian, StartRestrictions, Npsi)
oppList={Hamiltonian, OppSsqr, OppLsqr, OppJsqr, OppSz_3d, OppLz_3d, Oppldots_3d, OppF2_3d, OppF4_3d, OppNeg_3d, OppNt2g_3d, OppNeg_Ld, OppNt2g_Ld, OppN_3d}
-- print of some expectation values
print(" # ");
for i = 1,#psiList do
io.write(string.format("%3i ",i))
for j = 1,#oppList do
expectationvalue = Chop(psiList[i]*oppList[j]*psiList[i])
io.write(string.format("%8.3f ",expectationvalue))
end
io.write("\n")
end
-- We calculate the x-ray absorption spectra for z, right circular and left circular polarized light for the 3 lowest eigen-states. (9 spectra in total)
XASSpectra = CreateSpectra(XASHamiltonian, {TXASz, TXASr, TXASl}, psiList, {{"Emin",-15}, {"Emax",25}, {"NE",2000}, {"Gamma",0.1}})
-- and put some additional energy broadening on it (0.4 Gaussian and energy dependent lorenzian)
XASSpectra.Broaden(0.4, {{-3.7, 0.45}, {-2.2, 0.65}, { 0.0, 0.65}, { 1.0, 2.00}, { 6 , 2.00}, { 8 , 0.80}, {13.2, 0.80}, {14.0, 0.90}, {16.0, 0.90}, {17.0, 2.00}})
-- and now we start to do things different in order to include temperature averaging
-- We have calculated three states. The energy of these states is:
Enlist = {}
for i=1,#psiList do
Enlist[i] = psiList[i] * Hamiltonian * psiList[i]
end
-- In order to calculate occupation numbers we need to calculate E^(-Energy/(kb T))
-- If the ground-state has an energy of -3.5 eV and we calculate this pre-factor at a temperature
-- of 10 Kelvin we get E^(3.5/(10*8.6*10^(-5))) = 8.3 * 10^(1763) a number so large it does
-- not fit in the computer memory. The solution is simple, we need to make sure the lowest
-- energy is zero.
for i=2,#psiList do
Enlist[i] = Enlist[i] - Enlist[1]
end
Enlist[1]=0
-- Besides the energy I would like to look at the magnetic moment, both spin and angular part
-- so also here we calculate the expectation values in a list
Szlist = {}
Lzlist = {}
for i=1,#psiList do
Szlist[i] = psiList[i] * OppSz_3d * psiList[i]
Lzlist[i] = psiList[i] * OppLz_3d * psiList[i]
end
-- We can now calcualte the temperature dependent expectation values:
p={}
print("Temperature, Total Energy, Magnetic Moment, Sz, Lz")
for T=1,1000,10 do
Z=0
SzT=0
LzT=0
ET=0
for i=1,#psiList do
p[i] = exp(-Enlist[i]/(EnergyUnits.Kelvin.value * T))
Z = Z + p[i]
SzT = SzT + p[i] * Szlist[i]
LzT = LzT + p[i] * Lzlist[i]
ET = ET + p[i] * Enlist[i]
end
SzT = SzT / Z
LzT = LzT / Z
MzT = -LzT - 2*SzT
io.write(string.format("%8i ",T))
io.write(string.format("%14.8f ",ET))
io.write(string.format("%14.8f ",MzT))
io.write(string.format("%14.8f ",SzT))
io.write(string.format("%14.8f\n",LzT))
end
-- Note that in order to see the magnetic phase transition you need to make Hz in the Hamiltonian
-- temperature dependent. This can be done using self consistent loops. (The mathematica version
-- has an example of this, will make it here at some point as well)
-- The first column shows you how much energy you need to add to the system in order to heat
-- it. (specific heat) Be aware though that most of the specific heat is due to phonons, not
-- included in this calculation. It does capture nicely the electronic contribution to the
-- specific heat. (Including the change at cross overs for excited states (important in rare-
-- earth systems) and Lambda peaks for phase transitions
-- Now the temperature dependent XAS:
-- The object XASSpectra contains 9 spectra. For 3 different polarizations and 3 different states
-- What we need to do is to sum them according to Boltzmann statistics.
T = 100
Z=0
for i=1,#psiList do
p[i] = exp(-Enlist[i]/(EnergyUnits.Kelvin.value * T))
Z = Z + p[i]
end
-- we now create the Bolzmann summ for z, right and left polarized absorption (at T=100 Kelvin)
XASSpectraT100 = Spectra.Sum(XASSpectra,{p[1]/Z,p[2]/Z,p[3]/Z, 0,0,0, 0,0,0},{0,0,0, p[1]/Z,p[2]/Z,p[3]/Z, 0,0,0},{0,0,0, 0,0,0, p[1]/Z,p[2]/Z,p[3]/Z})
-- and the same at 1000 Kelvin
T = 1000
Z=0
for i=1,#psiList do
p[i] = exp(-Enlist[i]/(EnergyUnits.Kelvin.value * T))
Z = Z + p[i]
end
-- we now create the Bolzmann summ for z, right and left polarized absorption (at T=1000 Kelvin) (note again that we still have a large magnetization here as the exchange-field is not temperature dependent)
XASSpectraT1000 = Spectra.Sum(XASSpectra,{p[1]/Z,p[2]/Z,p[3]/Z, 0,0,0, 0,0,0},{0,0,0, p[1]/Z,p[2]/Z,p[3]/Z, 0,0,0},{0,0,0, 0,0,0, p[1]/Z,p[2]/Z,p[3]/Z})
-- We can print the spectra to file
XASSpectraT100.Print({{"file","XASSpectraT100.dat"}})
XASSpectraT1000.Print({{"file","XASSpectraT1000.dat"}})
gnuplotInput = [[
set autoscale
set xtic auto
set ytic auto
set style line 1 lt 1 lw 1 lc rgb "#0000FF"
set style line 2 lt 1 lw 1 lc rgb "#FF0000"
set style line 3 lt 1 lw 1 lc rgb "#00FF00"
set xlabel "E (eV)" font "Times,12"
set ylabel "Intensity (arb. units)" font "Times,12"
set out 'XASSpecT.ps'
set size 1.0, 0.6
set terminal postscript portrait enhanced color "Times" 12
energyshift=857.6
intensityscale=48
set xrange [847:877]
plot "XASSpectraT100.dat" using ($1+energyshift):((-$3-$5-$7) * intensityscale) title 'isotropic theory T=100K' with lines ls 1,\
"XASSpectraT1000.dat" using ($1+energyshift):((-$3-$5-$7) * intensityscale) title 'isotropic theory T=1000K' with lines ls 2,\
"NiO_Experiment/XAS_L23_PRB_57_11623_1998" using 1:2 title 'isotropic experiment' with lines ls 3,\
"XASSpectraT100.dat" using ($1+energyshift):(($5-$7) * intensityscale) title 'XMCD theory T=100K' with lines ls 1,\
"XASSpectraT1000.dat" using ($1+energyshift):(($5-$7) * intensityscale) title 'XMCD theory T=1000K' with lines ls 2
]]
-- write the gnuplot script to a file
file = io.open("XASSpecT.gnuplot", "w")
file:write(gnuplotInput)
file:close()
-- and finally call gnuplot to execute the script
os.execute("gnuplot XASSpecT.gnuplot")
-- as I like pdf to view and eps to include in the manuel I transform the format
os.execute(" ps2pdf XASSpecT.ps ; ps2eps XASSpecT.ps ; mv XASSpecT.eps temp.eps ; eps2eps temp.eps XASSpecT.eps ; rm temp.eps")
###
###
The output is:
#
1 -3.503 1.999 12.000 15.095 -0.908 -0.281 -0.305 -1.042 -0.924 2.186 5.990 3.825 6.000 8.175
2 -3.395 1.999 12.000 15.160 -0.004 -0.002 -0.322 -1.043 -0.925 2.189 5.988 3.823 6.000 8.178
3 -3.286 1.999 12.000 15.211 0.903 0.278 -0.336 -1.043 -0.925 2.193 5.987 3.820 6.000 8.180
Temperature, Total Energy, Magnetic Moment, Sz, Lz
1 0.00000000 2.09639244 -0.90751752 -0.28135740
11 0.00000000 2.09639244 -0.90751752 -0.28135740
21 0.00000000 2.09639244 -0.90751752 -0.28135740
31 0.00000000 2.09639244 -0.90751752 -0.28135740
41 0.00000000 2.09639244 -0.90751752 -0.28135740
51 0.00000000 2.09639244 -0.90751752 -0.28135740
61 0.00000000 2.09639244 -0.90751752 -0.28135740
71 0.00000000 2.09639240 -0.90751750 -0.28135739
81 0.00000002 2.09639208 -0.90751736 -0.28135735
91 0.00000011 2.09639041 -0.90751664 -0.28135713
101 0.00000042 2.09638445 -0.90751406 -0.28135633
111 0.00000128 2.09636784 -0.90750687 -0.28135410
121 0.00000327 2.09632960 -0.90749031 -0.28134898
131 0.00000724 2.09625331 -0.90745727 -0.28133877
141 0.00001431 2.09611725 -0.90739835 -0.28132054
151 0.00002587 2.09589513 -0.90730217 -0.28129080
161 0.00004345 2.09555736 -0.90715590 -0.28124556
171 0.00006869 2.09507255 -0.90694596 -0.28118063
181 0.00010326 2.09440895 -0.90665860 -0.28109175
191 0.00014878 2.09353581 -0.90628050 -0.28097481
201 0.00020678 2.09242436 -0.90579920 -0.28082596
211 0.00027865 2.09104862 -0.90520345 -0.28064171
221 0.00036565 2.08938584 -0.90448341 -0.28041902
231 0.00046884 2.08741683 -0.90363076 -0.28015531
241 0.00058915 2.08512592 -0.90263871 -0.27984850
251 0.00072732 2.08250098 -0.90150202 -0.27949695
261 0.00088395 2.07953318 -0.90021685 -0.27909949
271 0.00105950 2.07621676 -0.89878071 -0.27865534
281 0.00125430 2.07254875 -0.89719232 -0.27816410
291 0.00146857 2.06852868 -0.89545148 -0.27762573
301 0.00170243 2.06415828 -0.89355892 -0.27704043
311 0.00195590 2.05944119 -0.89151624 -0.27640872
321 0.00222896 2.05438273 -0.88932572 -0.27573129
331 0.00252149 2.04898959 -0.88699027 -0.27500905
341 0.00283333 2.04326969 -0.88451332 -0.27424305
351 0.00316429 2.03723192 -0.88189871 -0.27343450
361 0.00351411 2.03088598 -0.87915065 -0.27258468
371 0.00388252 2.02424225 -0.87627363 -0.27169500
381 0.00426920 2.01731161 -0.87327236 -0.27076690
391 0.00467384 2.01010537 -0.87015173 -0.26980191
401 0.00509608 2.00263511 -0.86691677 -0.26880157
411 0.00553555 1.99491262 -0.86357258 -0.26776747
421 0.00599187 1.98694984 -0.86012432 -0.26670120
431 0.00646466 1.97875875 -0.85657719 -0.26560437
441 0.00695351 1.97035131 -0.85293636 -0.26447859
451 0.00745802 1.96173945 -0.84920700 -0.26332544
461 0.00797777 1.95293497 -0.84539423 -0.26214651
471 0.00851235 1.94394955 -0.84150309 -0.26094337
481 0.00906134 1.93479468 -0.83753856 -0.25971755
491 0.00962431 1.92548164 -0.83350554 -0.25847056
501 0.01020084 1.91602149 -0.82940880 -0.25720389
511 0.01079052 1.90642502 -0.82525301 -0.25591899
521 0.01139292 1.89670276 -0.82104276 -0.25461725
531 0.01200762 1.88686496 -0.81678245 -0.25330005
541 0.01263421 1.87692156 -0.81247641 -0.25196873
551 0.01327227 1.86688221 -0.80812881 -0.25062458
561 0.01392140 1.85675621 -0.80374369 -0.24926884
571 0.01458119 1.84655259 -0.79932494 -0.24790272
581 0.01525123 1.83628003 -0.79487633 -0.24652738
591 0.01593114 1.82594688 -0.79040147 -0.24514394
601 0.01662053 1.81556118 -0.78590385 -0.24375348
611 0.01731900 1.80513065 -0.78138681 -0.24235703
621 0.01802619 1.79466269 -0.77685356 -0.24095557
631 0.01874172 1.78416436 -0.77230714 -0.23955007
641 0.01946523 1.77364244 -0.76775051 -0.23814142
651 0.02019636 1.76310338 -0.76318645 -0.23673048
661 0.02093476 1.75255334 -0.75861763 -0.23531809
671 0.02168008 1.74199819 -0.75404658 -0.23390503
681 0.02243199 1.73144349 -0.74947573 -0.23249203
691 0.02319017 1.72089454 -0.74490736 -0.23107982
701 0.02395428 1.71035636 -0.74034365 -0.22966906
711 0.02472402 1.69983370 -0.73578666 -0.22826038
721 0.02549907 1.68933105 -0.73123833 -0.22685440
731 0.02627915 1.67885267 -0.72670050 -0.22545167
741 0.02706395 1.66840255 -0.72217491 -0.22405274
751 0.02785320 1.65798447 -0.71766318 -0.22265810
761 0.02864662 1.64760198 -0.71316687 -0.22126824
771 0.02944393 1.63725840 -0.70868740 -0.21988360
781 0.03024489 1.62695686 -0.70422614 -0.21850459
791 0.03104922 1.61670029 -0.69978434 -0.21713161
801 0.03185669 1.60649141 -0.69536319 -0.21576503
811 0.03266706 1.59633277 -0.69096380 -0.21440517
821 0.03348008 1.58622675 -0.68658719 -0.21305237
831 0.03429555 1.57617555 -0.68223432 -0.21170691
841 0.03511323 1.56618121 -0.67790606 -0.21036908
851 0.03593291 1.55624561 -0.67360325 -0.20903911
861 0.03675440 1.54637050 -0.66932663 -0.20771725
871 0.03757748 1.53655749 -0.66507689 -0.20640371
881 0.03840196 1.52680804 -0.66085468 -0.20509868
891 0.03922767 1.51712348 -0.65666057 -0.20380234
901 0.04005441 1.50750505 -0.65249509 -0.20251486
911 0.04088201 1.49795385 -0.64835873 -0.20123639
921 0.04171030 1.48847088 -0.64425191 -0.19996706
931 0.04253913 1.47905702 -0.64017502 -0.19870698
941 0.04336832 1.46971308 -0.63612841 -0.19745627
951 0.04419774 1.46043977 -0.63211237 -0.19621502
961 0.04502722 1.45123768 -0.62812719 -0.19498331
971 0.04585664 1.44210737 -0.62417308 -0.19376121
981 0.04668585 1.43304928 -0.62025025 -0.19254878
991 0.04751472 1.42406379 -0.61635886 -0.19134607
###
===== Table of contents =====
{{indexmenu>.#1|msort}}